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 "PathBase.h" 23 : #include "core/ActionRegister.h" 24 : 25 : //+PLUMEDOC COLVAR GPROPERTYMAP 26 : /* 27 : Property maps but with a more flexible framework for the distance metric being used. 28 : 29 : This colvar calculates a property map using the formalism developed by Spiwok \cite Spiwok:2011ce. 30 : In essence if you have the value of some property, \f$X_i\f$, that it takes at a set of high-dimensional 31 : positions then you calculate the value of the property at some arbitrary point in the high-dimensional space 32 : using: 33 : 34 : \f[ 35 : X=\frac{\sum_i X_i*\exp(-\lambda D_i(x))}{\sum_i \exp(-\lambda D_i(x))} 36 : \f] 37 : 38 : Within PLUMED there are multiple ways to define the distance from a high-dimensional configuration, \f$D_i\f$. You could calculate 39 : the RMSD distance or you could calculate the amount by which a set of collective variables change. As such this implementation 40 : of the property map allows one to use all the different distance metric that are discussed in \ref dists. This is as opposed to 41 : the alternative implementation \ref PROPERTYMAP which is a bit faster but which only allows one to use the RMSD distance. 42 : 43 : \par Examples 44 : 45 : The input shown below can be used to calculate the interpolated values of two properties called X and Y based on the values 46 : that these properties take at a set of reference configurations and using the formula above. For this input the distances 47 : between the reference configurations and the instantaneous configurations are calculated using the OPTIMAL metric that is 48 : discussed at length in the manual pages on \ref RMSD. 49 : 50 : \plumedfile 51 : p2: GPROPERTYMAP REFERENCE=allv.pdb PROPERTY=X,Y LAMBDA=69087 52 : PRINT ARG=p2.X,p2.Y,p2.zpath STRIDE=1 FILE=colvar 53 : \endplumedfile 54 : 55 : The additional input file for this calculation, which contains the reference frames and the values of X and Y at these reference 56 : points has the following format. 57 : 58 : \auxfile{allv.pdb} 59 : REMARK X=1 Y=2 60 : ATOM 1 CL ALA 1 -3.171 0.295 2.045 1.00 1.00 61 : ATOM 5 CLP ALA 1 -1.819 -0.143 1.679 1.00 1.00 62 : ATOM 6 OL ALA 1 -1.177 -0.889 2.401 1.00 1.00 63 : ATOM 7 NL ALA 1 -1.313 0.341 0.529 1.00 1.00 64 : ATOM 8 HL ALA 1 -1.845 0.961 -0.011 1.00 1.00 65 : ATOM 9 CA ALA 1 -0.003 -0.019 0.021 1.00 1.00 66 : ATOM 10 HA ALA 1 0.205 -1.051 0.259 1.00 1.00 67 : ATOM 11 CB ALA 1 0.009 0.135 -1.509 1.00 1.00 68 : ATOM 15 CRP ALA 1 1.121 0.799 0.663 1.00 1.00 69 : ATOM 16 OR ALA 1 1.723 1.669 0.043 1.00 1.00 70 : ATOM 17 NR ALA 1 1.423 0.519 1.941 1.00 1.00 71 : ATOM 18 HR ALA 1 0.873 -0.161 2.413 1.00 1.00 72 : ATOM 19 CR ALA 1 2.477 1.187 2.675 1.00 1.00 73 : END 74 : FIXED 75 : REMARK X=2 Y=3 76 : ATOM 1 CL ALA 1 -3.175 0.365 2.024 1.00 1.00 77 : ATOM 5 CLP ALA 1 -1.814 -0.106 1.685 1.00 1.00 78 : ATOM 6 OL ALA 1 -1.201 -0.849 2.425 1.00 1.00 79 : ATOM 7 NL ALA 1 -1.296 0.337 0.534 1.00 1.00 80 : ATOM 8 HL ALA 1 -1.807 0.951 -0.044 1.00 1.00 81 : ATOM 9 CA ALA 1 0.009 -0.067 0.033 1.00 1.00 82 : ATOM 10 HA ALA 1 0.175 -1.105 0.283 1.00 1.00 83 : ATOM 11 CB ALA 1 0.027 0.046 -1.501 1.00 1.00 84 : ATOM 15 CRP ALA 1 1.149 0.725 0.654 1.00 1.00 85 : ATOM 16 OR ALA 1 1.835 1.491 -0.011 1.00 1.00 86 : ATOM 17 NR ALA 1 1.380 0.537 1.968 1.00 1.00 87 : ATOM 18 HR ALA 1 0.764 -0.060 2.461 1.00 1.00 88 : ATOM 19 CR ALA 1 2.431 1.195 2.683 1.00 1.00 89 : END 90 : \endauxfile 91 : 92 : */ 93 : //+ENDPLUMEDOC 94 : 95 : namespace PLMD { 96 : namespace mapping { 97 : 98 : class PropertyMap : public PathBase { 99 : public: 100 : static void registerKeywords( Keywords& keys ); 101 : explicit PropertyMap(const ActionOptions&); 102 : }; 103 : 104 10425 : PLUMED_REGISTER_ACTION(PropertyMap,"GPROPERTYMAP") 105 : 106 4 : void PropertyMap::registerKeywords( Keywords& keys ) { 107 4 : PathBase::registerKeywords( keys ); 108 4 : ActionWithValue::useCustomisableComponents( keys ); 109 8 : keys.addFlag("NOMAPPING",false,"do not calculate the position on the manifold"); 110 4 : } 111 : 112 3 : PropertyMap::PropertyMap(const ActionOptions& ao): 113 : Action(ao), 114 3 : PathBase(ao) 115 : { 116 6 : bool nos; parseFlag("NOMAPPING",nos); 117 : 118 : std::string empty; 119 3 : if(!nos) { 120 9 : for(std::map<std::string,std::vector<double> >::iterator it=property.begin(); it!=property.end(); ++it) { 121 18 : empty="LABEL="+it->first; addVessel( "SPATH", empty, 0 ); 122 : } 123 : } 124 3 : readVesselKeywords(); 125 3 : checkRead(); 126 3 : } 127 : 128 : } 129 : }