LCOV - code coverage report
Current view: top level - mapping - Path.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 93 102 91.2 %
Date: 2024-10-18 13:59:31 Functions: 6 7 85.7 %

          Line data    Source code
       1             : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
       2             :    Copyright (c) 2013-2023 The plumed team
       3             :    (see the PEOPLE file at the root of the distribution for a list of names)
       4             : 
       5             :    See http://www.plumed.org for more information.
       6             : 
       7             :    This file is part of plumed, version 2.
       8             : 
       9             :    plumed is free software: you can redistribute it and/or modify
      10             :    it under the terms of the GNU Lesser General Public License as published by
      11             :    the Free Software Foundation, either version 3 of the License, or
      12             :    (at your option) any later version.
      13             : 
      14             :    plumed is distributed in the hope that it will be useful,
      15             :    but WITHOUT ANY WARRANTY; without even the implied warranty of
      16             :    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
      17             :    GNU Lesser General Public License for more details.
      18             : 
      19             :    You should have received a copy of the GNU Lesser General Public License
      20             :    along with plumed.  If not, see <http://www.gnu.org/licenses/>.
      21             : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
      22             : #include "Path.h"
      23             : #include "core/ActionRegister.h"
      24             : #include "core/ActionWithValue.h"
      25             : #include "core/ActionWithArguments.h"
      26             : #include "core/PlumedMain.h"
      27             : #include "core/ActionSet.h"
      28             : #include "tools/PDB.h"
      29             : 
      30             : //+PLUMEDOC COLVAR PATH
      31             : /*
      32             : Path collective variables with a more flexible framework for the distance metric being used.
      33             : 
      34             : The Path Collective Variables developed by Branduardi and co-workers \cite brand07 allow one
      35             : to compute the progress along a high-dimensional path and the distance from the high-dimensional
      36             : path.  The progress along the path (s) is computed using:
      37             : 
      38             : \f[
      39             : s = \frac{ \sum_{i=1}^N i \exp( -\lambda R[X - X_i] ) }{ \sum_{i=1}^N \exp( -\lambda R[X - X_i] ) }
      40             : \f]
      41             : 
      42             : while the distance from the path (z) is measured using:
      43             : 
      44             : \f[
      45             : z = -\frac{1}{\lambda} \ln\left[ \sum_{i=1}^N \exp( -\lambda R[X - X_i] ) \right]
      46             : \f]
      47             : 
      48             : In these expressions \f$N\f$ high-dimensional frames (\f$X_i\f$) are used to describe the path in the high-dimensional
      49             : space. The two expressions above are then functions of the distances from each of the high-dimensional frames \f$R[X - X_i]\f$.
      50             : Within PLUMED there are multiple ways to define the distance from a high-dimensional configuration.  You could calculate
      51             : the RMSD distance or you could calculate the amount by which a set of collective variables change.  As such this implementation
      52             : of the path CV allows one to use all the difference distance metrics that are discussed in \ref dists. This is as opposed to
      53             : the alternative implementation of path (\ref PATHMSD) which is a bit faster but which only allows one to use the RMSD distance.
      54             : 
      55             : The \f$s\f$ and \f$z\f$ variables are calculated using the above formulas by default.  However, there is an alternative method
      56             : of calculating these collective variables, which is detailed in \cite bernd-path.  This alternative method uses the tools of
      57             : geometry (as opposed to algebra, which is used in the equations above).  In this alternative formula the progress along the path
      58             : \f$s\f$ is calculated using:
      59             : 
      60             : \f[
      61             : s = i_2 + \textrm{sign}(i_2-i_1) \frac{ \sqrt{( \mathbf{v}_1\cdot\mathbf{v}_2 )^2 - |\mathbf{v}_3|^2(|\mathbf{v}_1|^2 - |\mathbf{v}_2|^2) } }{2|\mathbf{v}_3|^2} - \frac{\mathbf{v}_1\cdot\mathbf{v}_3 - |\mathbf{v}_3|^2}{2|\mathbf{v}_3|^2}
      62             : \f]
      63             : 
      64             : where \f$\mathbf{v}_1\f$ and \f$\mathbf{v}_3\f$ are the vectors connecting the current position to the closest and second closest node of the path,
      65             : respectfully and \f$i_1\f$ and \f$i_2\f$ are the projections of the closest and second closest frames of the path. \f$\mathbf{v}_2\f$, meanwhile, is the
      66             : vector connecting the closest frame to the second closest frame.  The distance from the path, \f$z\f$ is calculated using:
      67             : 
      68             : \f[
      69             : z = \sqrt{ \left[ |\mathbf{v}_1|^2 - |\mathbf{v}_2| \left( \frac{ \sqrt{( \mathbf{v}_1\cdot\mathbf{v}_2 )^2 - |\mathbf{v}_3|^2(|\mathbf{v}_1|^2 - |\mathbf{v}_2|^2) } }{2|\mathbf{v}_3|^2} - \frac{\mathbf{v}_1\cdot\mathbf{v}_3 - |\mathbf{v}_3|^2}{2|\mathbf{v}_3|^2} \right) \right]^2 }
      70             : \f]
      71             : 
      72             : The symbols here are as they were for \f$s\f$.  If you would like to use these equations to calculate \f$s\f$ and \f$z\f$ then you should use the GPATH flag.
      73             : The values of \f$s\f$ and \f$z\f$ can then be referenced using the gspath and gzpath labels.
      74             : 
      75             : \par Examples
      76             : 
      77             : In the example below the path is defined using RMSD distance from frames.
      78             : 
      79             : \plumedfile
      80             : p1: PATH REFERENCE=file.pdb TYPE=OPTIMAL LAMBDA=500.0
      81             : PRINT ARG=p1.spath,p1.zpath STRIDE=1 FILE=colvar FMT=%8.4f
      82             : \endplumedfile
      83             : 
      84             : The reference frames in the path are defined in the pdb file shown below.  In this frame
      85             : each configuration in the path is separated by a line containing just the word END.
      86             : 
      87             : \auxfile{file.pdb}
      88             : ATOM      1  CL  ALA     1      -3.171   0.295   2.045  1.00  1.00
      89             : ATOM      5  CLP ALA     1      -1.819  -0.143   1.679  1.00  1.00
      90             : ATOM      6  OL  ALA     1      -1.177  -0.889   2.401  1.00  1.00
      91             : ATOM      7  NL  ALA     1      -1.313   0.341   0.529  1.00  1.00
      92             : END
      93             : ATOM      1  CL  ALA     1      -3.175   0.365   2.024  1.00  1.00
      94             : ATOM      5  CLP ALA     1      -1.814  -0.106   1.685  1.00  1.00
      95             : ATOM      6  OL  ALA     1      -1.201  -0.849   2.425  1.00  1.00
      96             : ATOM      7  NL  ALA     1      -1.296   0.337   0.534  1.00  1.00
      97             : END
      98             : ATOM      1  CL  ALA     1      -2.990   0.383   2.277  1.00  1.00
      99             : ATOM      5  CLP ALA     1      -1.664  -0.085   1.831  1.00  1.00
     100             : ATOM      6  OL  ALA     1      -0.987  -0.835   2.533  1.00  1.00
     101             : ATOM      7  NL  ALA     1      -1.227   0.364   0.646  1.00  1.00
     102             : END
     103             : \endauxfile
     104             : 
     105             : In the example below the path is defined using the values of two torsional angles (t1 and t2).
     106             : In addition, the \f$s\f$ and \f$z\f$ are calculated using the geometric expressions described
     107             : above rather than the algebraic expressions that are used by default.
     108             : 
     109             : \plumedfile
     110             : t1: TORSION ATOMS=5,7,9,15
     111             : t2: TORSION ATOMS=7,9,15,17
     112             : pp: PATH TYPE=EUCLIDEAN REFERENCE=epath.pdb GPATH NOSPATH NOZPATH
     113             : PRINT ARG=pp.* FILE=colvar
     114             : \endplumedfile
     115             : 
     116             : Notice that the LAMBDA parameter is not required here as we are not calculating \f$s\f$ and \f$s\f$
     117             : using the algebraic formulas defined earlier.  The positions of the frames in the path are defined
     118             : in the file epath.pdb.  An extract from this file looks as shown below.
     119             : 
     120             : \auxfile{epath.pdb}
     121             : REMARK ARG=t1,t2 t1=-4.25053  t2=3.88053
     122             : END
     123             : REMARK ARG=t1,t2 t1=-4.11     t2=3.75
     124             : END
     125             : REMARK ARG=t1,t2 t1=-3.96947  t2=3.61947
     126             : END
     127             : \endauxfile
     128             : 
     129             : The remarks in this pdb file tell PLUMED the labels that are being used to define the position in the
     130             : high dimensional space and the values that these arguments have at each point on the path.
     131             : 
     132             : */
     133             : //+ENDPLUMEDOC
     134             : 
     135             : //+PLUMEDOC COLVAR GPROPERTYMAP
     136             : /*
     137             : Property maps but with a more flexible framework for the distance metric being used.
     138             : 
     139             : This colvar calculates a property map using the formalism developed by Spiwok \cite Spiwok:2011ce.
     140             : In essence if you have the value of some property, \f$X_i\f$, that it takes at a set of high-dimensional
     141             : positions then you calculate the value of the property at some arbitrary point in the high-dimensional space
     142             : using:
     143             : 
     144             : \f[
     145             : X=\frac{\sum_i X_i*\exp(-\lambda D_i(x))}{\sum_i  \exp(-\lambda D_i(x))}
     146             : \f]
     147             : 
     148             : Within PLUMED there are multiple ways to define the distance from a high-dimensional configuration, \f$D_i\f$.  You could calculate
     149             : the RMSD distance or you could calculate the amount by which a set of collective variables change.  As such this implementation
     150             : of the property map allows one to use all the different distance metric that are discussed in \ref dists. This is as opposed to
     151             : the alternative implementation \ref PROPERTYMAP which is a bit faster but which only allows one to use the RMSD distance.
     152             : 
     153             : \par Examples
     154             : 
     155             : The input shown below can be used to calculate the interpolated values of two properties called X and Y based on the values
     156             : that these properties take at a set of reference configurations and using the formula above.  For this input the distances
     157             : between the reference configurations and the instantaneous configurations are calculated using the OPTIMAL metric that is
     158             : discussed at length in the manual pages on \ref RMSD.
     159             : 
     160             : \plumedfile
     161             : p2: GPROPERTYMAP REFERENCE=allv.pdb PROPERTY=X,Y LAMBDA=69087
     162             : PRINT ARG=p2.X,p2.Y,p2.zpath STRIDE=1 FILE=colvar
     163             : \endplumedfile
     164             : 
     165             : The additional input file for this calculation, which contains the reference frames and the values of X and Y at these reference
     166             : points has the following format.
     167             : 
     168             : \auxfile{allv.pdb}
     169             : REMARK X=1 Y=2
     170             : ATOM      1  CL  ALA     1      -3.171   0.295   2.045  1.00  1.00
     171             : ATOM      5  CLP ALA     1      -1.819  -0.143   1.679  1.00  1.00
     172             : ATOM      6  OL  ALA     1      -1.177  -0.889   2.401  1.00  1.00
     173             : ATOM      7  NL  ALA     1      -1.313   0.341   0.529  1.00  1.00
     174             : ATOM      8  HL  ALA     1      -1.845   0.961  -0.011  1.00  1.00
     175             : ATOM      9  CA  ALA     1      -0.003  -0.019   0.021  1.00  1.00
     176             : ATOM     10  HA  ALA     1       0.205  -1.051   0.259  1.00  1.00
     177             : ATOM     11  CB  ALA     1       0.009   0.135  -1.509  1.00  1.00
     178             : ATOM     15  CRP ALA     1       1.121   0.799   0.663  1.00  1.00
     179             : ATOM     16  OR  ALA     1       1.723   1.669   0.043  1.00  1.00
     180             : ATOM     17  NR  ALA     1       1.423   0.519   1.941  1.00  1.00
     181             : ATOM     18  HR  ALA     1       0.873  -0.161   2.413  1.00  1.00
     182             : ATOM     19  CR  ALA     1       2.477   1.187   2.675  1.00  1.00
     183             : END
     184             : FIXED
     185             : REMARK X=2 Y=3
     186             : ATOM      1  CL  ALA     1      -3.175   0.365   2.024  1.00  1.00
     187             : ATOM      5  CLP ALA     1      -1.814  -0.106   1.685  1.00  1.00
     188             : ATOM      6  OL  ALA     1      -1.201  -0.849   2.425  1.00  1.00
     189             : ATOM      7  NL  ALA     1      -1.296   0.337   0.534  1.00  1.00
     190             : ATOM      8  HL  ALA     1      -1.807   0.951  -0.044  1.00  1.00
     191             : ATOM      9  CA  ALA     1       0.009  -0.067   0.033  1.00  1.00
     192             : ATOM     10  HA  ALA     1       0.175  -1.105   0.283  1.00  1.00
     193             : ATOM     11  CB  ALA     1       0.027   0.046  -1.501  1.00  1.00
     194             : ATOM     15  CRP ALA     1       1.149   0.725   0.654  1.00  1.00
     195             : ATOM     16  OR  ALA     1       1.835   1.491  -0.011  1.00  1.00
     196             : ATOM     17  NR  ALA     1       1.380   0.537   1.968  1.00  1.00
     197             : ATOM     18  HR  ALA     1       0.764  -0.060   2.461  1.00  1.00
     198             : ATOM     19  CR  ALA     1       2.431   1.195   2.683  1.00  1.00
     199             : END
     200             : \endauxfile
     201             : 
     202             : */
     203             : //+ENDPLUMEDOC
     204             : 
     205             : namespace PLMD {
     206             : namespace mapping {
     207             : 
     208             : PLUMED_REGISTER_ACTION(Path,"PATH")
     209             : PLUMED_REGISTER_ACTION(Path,"GPROPERTYMAP")
     210             : 
     211          14 : void Path::registerKeywords( Keywords& keys ) {
     212          14 :   ActionShortcut::registerKeywords( keys ); Path::registerInputFileKeywords( keys );
     213          28 :   keys.add("optional","PROPERTY","the property to be used in the index. This should be in the REMARK of the reference");
     214          28 :   keys.add("compulsory","LAMBDA","the lambda parameter is needed for smoothing, is in the units of plumed");
     215          28 :   keys.addOutputComponent("gspath","GPATH","scalar","the position along the path calculated using the geometric formula");
     216          28 :   keys.addOutputComponent("gzpath","GPATH","scalar","the distance from the path calculated using the geometric formula");
     217          28 :   keys.addOutputComponent("spath","default","scalar","the position along the path calculated");
     218          28 :   keys.addOutputComponent("zpath","default","scalar","the distance from the path calculated");
     219          14 : }
     220             : 
     221          32 : void Path::registerInputFileKeywords( Keywords& keys ) {
     222          64 :   keys.add("compulsory","REFERENCE","a pdb file containing the set of reference configurations");
     223          64 :   keys.add("compulsory","TYPE","OPTIMAL-FAST","the manner in which distances are calculated. More information on the different "
     224             :            "metrics that are available in PLUMED can be found in the section of the manual on "
     225             :            "\\ref dists");
     226          64 :   keys.addInputKeyword("optional","ARG","scalar","the list of arguments you would like to use in your definition of the path");
     227          64 :   keys.add("optional","COEFFICIENTS","the coefficients of the displacements along each argument that should be used when calculating the euclidean distance");
     228          64 :   keys.addFlag("NOPBC",false,"ignore the periodic boundary conditions when calculating distances");
     229          64 :   keys.addFlag("NOSPATH",false,"do not calculate the spath CV");
     230          64 :   keys.addFlag("NOZPATH",false,"do not calculate the zpath CV");
     231          64 :   keys.addFlag("GPATH",false,"calculate the trigonometric path");
     232         128 :   keys.needsAction("DRMSD"); keys.needsAction("RMSD"); keys.needsAction("LOWEST"); keys.needsAction("GPATH");
     233         128 :   keys.needsAction("EUCLIDEAN_DISTANCE"); keys.needsAction("CUSTOM"); keys.needsAction("SUM"); keys.needsAction("COMBINE");
     234          96 :   keys.needsAction("NORMALIZED_EUCLIDEAN_DISTANCE"); keys.needsAction("PDB2CONSTANT"); keys.needsAction("CONSTANT");
     235          32 : }
     236             : 
     237          10 : Path::Path( const ActionOptions& ao ):
     238             :   Action(ao),
     239          10 :   ActionShortcut(ao)
     240             : {
     241          30 :   bool nospath, nozpath, gpath; parseFlag("NOSPATH",nospath); parseFlag("NOZPATH",nozpath); parseFlag("GPATH",gpath);
     242          10 :   if( gpath ) {
     243           0 :     readInputLine( getShortcutLabel() + "_gpath: GPATH " + convertInputLineToString() );
     244           0 :     readInputLine( getShortcutLabel() + "_gspath: COMBINE ARG=" + getShortcutLabel() + "_gpath.s PERIODIC=NO");
     245           0 :     readInputLine( getShortcutLabel() + "_gzpath: COMBINE ARG=" + getShortcutLabel() + "_gpath.z PERIODIC=NO");
     246             :   }
     247          10 :   if( nospath && nozpath ) return;
     248             :   // Setup the properties
     249             :   std::vector<std::string> properties, pnames;
     250          10 :   if( getName()=="PATH") { properties.resize(1); }
     251           6 :   else { parseVector("PROPERTY",pnames); properties.resize( pnames.size() ); }
     252          20 :   std::string type, reference_data, reference; parse("REFERENCE",reference);
     253          20 :   std::vector<std::string> argnames; parseVector("ARG",argnames); parse("TYPE",type);
     254          11 :   if( type.find("DRMSD")!=std::string::npos ) readInputLine( getShortcutLabel() + "_data: DRMSD SQUARED TYPE=" + type + " REFERENCE=" + reference );
     255           9 :   else readInputFrames( reference, type, argnames, false, this, reference_data );
     256             :   // Find the shortest distance to the frames
     257          20 :   readInputLine( getShortcutLabel() + "_mindist: LOWEST ARG=" + getShortcutLabel() + "_data");
     258             :   // Now create all other parts of the calculation
     259          10 :   std::string lambda; parse("LAMBDA",lambda);
     260             :   // Now create MATHEVAL object to compute exponential functions
     261          20 :   readInputLine( getShortcutLabel() + "_weights: CUSTOM ARG=" + getShortcutLabel() + "_data," + getShortcutLabel() + "_mindist FUNC=exp(-(x-y)*" + lambda + ") PERIODIC=NO" );
     262             :   // Create denominator
     263          20 :   readInputLine( getShortcutLabel() + "_denom: SUM ARG=" + getShortcutLabel() + "_weights PERIODIC=NO");
     264             :   // Now compte zpath variable
     265          10 :   if( !nozpath ) {
     266          20 :     readInputLine( getShortcutLabel() + "_z: CUSTOM ARG=" + getShortcutLabel() + "_denom," + getShortcutLabel() + "_mindist FUNC=y-log(x)/" + lambda + " PERIODIC=NO");
     267          20 :     readInputLine( getShortcutLabel() + "_zpath: COMBINE ARG=" + getShortcutLabel() + "_z PERIODIC=NO");
     268             :   }
     269             :   // Now get coefficients for properies for spath
     270          10 :   readPropertyInformation( pnames, getShortcutLabel(), reference, this );
     271             :   // Now create COMBINE objects to compute numerator of path
     272          23 :   for(unsigned i=0; i<properties.size(); ++i) {
     273          13 :     if( pnames.size()>0 ) {
     274          12 :       readInputLine( pnames[i] + "_numer_prod: CUSTOM ARG=" + getShortcutLabel() + "_weights," + pnames[i] + "_ref FUNC=x*y PERIODIC=NO");
     275          12 :       readInputLine( pnames[i] + "_numer: SUM ARG=" + pnames[i]  + "_numer_prod PERIODIC=NO");
     276          12 :       readInputLine( getShortcutLabel() + "_" + pnames[i] + ": CUSTOM ARG=" + pnames[i]  + "_numer," + getShortcutLabel() + "_denom FUNC=x/y PERIODIC=NO");
     277           7 :     } else if( !nospath ) {
     278          14 :       readInputLine( getShortcutLabel() + "_s_prod: CUSTOM ARG=" + getShortcutLabel() + "_weights," + getShortcutLabel() + "_ind FUNC=x*y PERIODIC=NO");
     279          14 :       readInputLine( getShortcutLabel()  + "_numer: SUM ARG=" + getShortcutLabel() + "_s_prod PERIODIC=NO");
     280          14 :       readInputLine( getShortcutLabel() + "_s: CUSTOM ARG=" + getShortcutLabel() + "_numer," + getShortcutLabel() + "_denom FUNC=x/y PERIODIC=NO");
     281          14 :       readInputLine( getShortcutLabel() + "_spath: COMBINE ARG=" + getShortcutLabel() + "_s PERIODIC=NO");
     282             :     }
     283             :   }
     284          20 : }
     285             : 
     286          19 : std::string Path::fixArgumentName( const std::string& argin ) {
     287          19 :   std::string argout = argin; std::size_t dot=argin.find(".");
     288          23 :   if( dot!=std::string::npos ) argout = argin.substr(0,dot) + "_" + argin.substr(dot+1);
     289          19 :   return argout;
     290             : }
     291             : 
     292          14 : void Path::readInputFrames( const std::string& reference, const std::string& type, std::vector<std::string>& argnames, const bool& displacements, ActionShortcut* action, std::string& reference_data ) {
     293          14 :   FILE* fp=std::fopen(reference.c_str(),"r"); PDB pdb; if(!fp) action->error("could not open reference file " + reference );
     294          14 :   bool do_read=pdb.readFromFilepointer(fp,false,0.1); if( !do_read ) action->error("missing file " + reference );
     295          14 :   if( pdb.getPositions().size()!=0 && argnames.size()==0 ) {
     296           9 :     reference_data = action->getShortcutLabel() + "_data_ref";
     297          11 :     if( displacements ) action->readInputLine( action->getShortcutLabel() + "_data: RMSD DISPLACEMENT SQUARED REFERENCE=" + reference + " TYPE=" + type );
     298          14 :     else action->readInputLine( action->getShortcutLabel() + "_data: RMSD SQUARED REFERENCE=" + reference + " TYPE=" + type );
     299           5 :   } else if( pdb.getPositions().size()!=0 ) {
     300           0 :     reference_data = action->getShortcutLabel() + "_atomdata_ref";
     301           0 :     if( displacements ) action->readInputLine( action->getShortcutLabel() + "_atomdata: RMSD DISPLACEMENT SQUARED REFERENCE=" + reference + " TYPE=" + type );
     302           0 :     else action->readInputLine( action->getShortcutLabel() + "_atomdata: RMSD SQUARED REFERENCE=" + reference + " TYPE=" + type );
     303           5 :   } else if( argnames.size()==0 ) {
     304           0 :     argnames.resize( pdb.getArgumentNames().size() );
     305           0 :     for(unsigned i=0; i<argnames.size(); ++i) argnames[i] = pdb.getArgumentNames()[i];
     306             :   }
     307          14 :   std::vector<Value*> theargs; if( argnames.size()>0 ) ActionWithArguments::interpretArgumentList( argnames, action->plumed.getActionSet(), action, theargs );
     308             : 
     309          14 :   if( theargs.size()>0 ) {
     310             :     std::string instargs, refargs;
     311          20 :     for(unsigned i=0; i<theargs.size(); ++i) {
     312          15 :       std::string iargn = fixArgumentName( theargs[i]->getName() );
     313          30 :       action->readInputLine( action->getShortcutLabel() + "_ref_" + iargn + ": PDB2CONSTANT REFERENCE=" + reference + " ARG=" + theargs[i]->getName() );
     314          25 :       if( i==0 ) { instargs=" ARG1=" + theargs[i]->getName(); refargs=" ARG2=" + action->getShortcutLabel() + "_ref_" + iargn; }
     315          30 :       else { instargs +="," + theargs[i]->getName(); refargs +="," + action->getShortcutLabel() + "_ref_" + iargn; }
     316          20 :       if( pdb.getPositions().size()==0 && i==0 ) reference_data = action->getShortcutLabel() + "_ref_" + iargn;
     317          20 :       else reference_data += "," + action->getShortcutLabel() + "_ref_" + iargn;
     318             :     }
     319          10 :     std::string comname="EUCLIDEAN_DISTANCE SQUARED"; std::string coeffstr; action->parse("COEFFICIENTS",coeffstr);
     320           5 :     if( coeffstr.length()>0 ) {
     321           1 :       if( displacements ) action->error("cannot use COEFFICIENTS arguments with GEOMETRIC PATH");
     322           2 :       action->readInputLine( action->getShortcutLabel() + "_coeff: CONSTANT VALUES=" + coeffstr );
     323           2 :       action->readInputLine( action->getShortcutLabel() + "_coeff2: CUSTOM ARG=" + action->getShortcutLabel() + "_coeff FUNC=x*x PERIODIC=NO");
     324           2 :       comname = "NORMALIZED_EUCLIDEAN_DISTANCE SQUARED METRIC=" + action->getShortcutLabel() + "_coeff2";
     325           4 :     } else if( displacements ) comname = "DISPLACEMENT";
     326             : 
     327          10 :     if( pdb.getPositions().size()==0 ) action->readInputLine( action->getShortcutLabel() + "_data: " + comname + instargs + refargs );
     328           0 :     else action->readInputLine( action->getShortcutLabel() + "_argdata: " + comname + instargs + refargs );
     329             :   }
     330          14 : }
     331             : 
     332          15 : void Path::readPropertyInformation( const std::vector<std::string>& pnames, const std::string& lab, const std::string& refname, ActionShortcut* action ) {
     333          15 :   if( pnames.size()>0 ) {
     334          15 :     for(unsigned i=0; i<pnames.size(); ++i) {
     335          20 :       action->readInputLine( pnames[i] + ": CONSTANT VALUE=1");
     336          20 :       action->readInputLine( pnames[i] + "_ref: PDB2CONSTANT REFERENCE=" + refname + " ARG=" + pnames[i] );
     337             :     }
     338             :   } else {
     339          10 :     ActionWithValue* av=action->plumed.getActionSet().selectWithLabel<ActionWithValue*>( lab + "_data" );
     340          10 :     unsigned nfram = av->copyOutput(0)->getShape()[0]; std::string indices = "VALUES=1";
     341         850 :     for(unsigned i=1; i<nfram; ++i) { std::string num; Tools::convert( i+1, num ); indices += "," + num; }
     342          20 :     action->readInputLine( lab + "_ind: CONSTANT " + indices );
     343             :   }
     344          15 : }
     345             : 
     346             : }
     347             : }

Generated by: LCOV version 1.16