Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2012-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 "Colvar.h"
23 : #include "tools/Communicator.h"
24 : #include "tools/NeighborList.h"
25 : #include "core/ActionRegister.h"
26 : #include "tools/SwitchingFunction.h"
27 :
28 : namespace PLMD {
29 : namespace colvar {
30 :
31 : //+PLUMEDOC COLVAR CONTACTMAP
32 : /*
33 : Calculate the distances between a number of pairs of atoms and transform each distance by a switching function.
34 :
35 : The transformed distance can be compared with a reference value in order to calculate the squared distance
36 : between two contact maps. Each distance can also be weighted for a given value. CONTACTMAP can be used together
37 : with \ref FUNCPATHMSD to define a path in the contactmap space.
38 :
39 : The individual contact map distances related to each contact can be accessed as components
40 : named `cm.contact-1`, `cm.contact-2`, etc, assuming that the label of the CONTACTMAP is `cm`.
41 :
42 : \par Examples
43 :
44 : The following example calculates switching functions based on the distances between atoms
45 : 1 and 2, 3 and 4 and 4 and 5. The values of these three switching functions are then output
46 : to a file named colvar.
47 :
48 : \plumedfile
49 : CONTACTMAP ATOMS1=1,2 ATOMS2=3,4 ATOMS3=4,5 ATOMS4=5,6 SWITCH={RATIONAL R_0=1.5} LABEL=f1
50 : PRINT ARG=f1.* FILE=colvar
51 : \endplumedfile
52 :
53 : The following example calculates the difference of the current contact map with respect
54 : to a reference provided. In this case REFERENCE is the fraction of contact that is formed
55 : (i.e. the distance between two atoms transformed with the SWITCH), while R_0 is the contact
56 : distance. WEIGHT gives the relative weight of each contact to the final distance measure.
57 :
58 : \plumedfile
59 : CONTACTMAP ...
60 : ATOMS1=1,2 REFERENCE1=0.1 WEIGHT1=0.5
61 : ATOMS2=3,4 REFERENCE2=0.5 WEIGHT2=1.0
62 : ATOMS3=4,5 REFERENCE3=0.25 WEIGHT3=1.0
63 : ATOMS4=5,6 REFERENCE4=0.0 WEIGHT4=0.5
64 : SWITCH={RATIONAL R_0=1.5}
65 : LABEL=cmap
66 : CMDIST
67 : ... CONTACTMAP
68 :
69 : PRINT ARG=cmap FILE=colvar
70 : \endplumedfile
71 :
72 : The next example calculates calculates fraction of native contacts (Q)
73 : for Trp-cage mini-protein. R_0 is the distance at which the switch function is guaranteed to
74 : be 1.0 – it doesn't really matter for Q and should be something very small, like 1 A.
75 : REF is the reference distance for the contact, e.g. the distance from a crystal structure.
76 : LAMBDA is the tolerance for the distance – if set to 1.0, the contact would have to have exactly
77 : the reference value to be formed; instead for lambda values of 1.5–1.8 are usually used to allow some slack.
78 : BETA is the softness of the switch function, default is 50nm.
79 : WEIGHT is the 1/(number of contacts) giving equal weight to each contact.
80 :
81 : When using native contact Q switch function, please cite \cite best2013
82 :
83 : \plumedfile
84 : # The full (much-longer) example available in regtest/basic/rt72/
85 :
86 : CONTACTMAP ...
87 : ATOMS1=1,67 SWITCH1={Q R_0=0.01 BETA=50.0 LAMBDA=1.5 REF=0.4059} WEIGHT1=0.003597
88 : ATOMS2=1,68 SWITCH2={Q R_0=0.01 BETA=50.0 LAMBDA=1.5 REF=0.4039} WEIGHT2=0.003597
89 : ATOMS3=1,69 SWITCH3={Q R_0=0.01 BETA=50.0 LAMBDA=1.5 REF=0.3215} WEIGHT3=0.003597
90 : ATOMS4=5,61 SWITCH4={Q R_0=0.01 BETA=50.0 LAMBDA=1.5 REF=0.4277} WEIGHT4=0.003597
91 : ATOMS5=5,67 SWITCH5={Q R_0=0.01 BETA=50.0 LAMBDA=1.5 REF=0.3851} WEIGHT5=0.003597
92 : ATOMS6=5,68 SWITCH6={Q R_0=0.01 BETA=50.0 LAMBDA=1.5 REF=0.3811} WEIGHT6=0.003597
93 : ATOMS7=5,69 SWITCH7={Q R_0=0.01 BETA=50.0 LAMBDA=1.5 REF=0.3133} WEIGHT7=0.003597
94 : LABEL=cmap
95 : SUM
96 : ... CONTACTMAP
97 :
98 : PRINT ARG=cmap FILE=colvar
99 : \endplumedfile
100 : (See also \ref switchingfunction)
101 :
102 : */
103 : //+ENDPLUMEDOC
104 :
105 : class ContactMap : public Colvar {
106 : private:
107 : bool pbc, serial, docomp, dosum, docmdist;
108 : std::unique_ptr<NeighborList> nl;
109 : std::vector<SwitchingFunction> sfs;
110 : std::vector<double> reference, weight;
111 : public:
112 : static void registerKeywords( Keywords& keys );
113 : explicit ContactMap(const ActionOptions&);
114 : // active methods:
115 : void calculate() override;
116 0 : void checkFieldsAllowed() override {}
117 : };
118 :
119 : PLUMED_REGISTER_ACTION(ContactMap,"CONTACTMAP")
120 :
121 10 : void ContactMap::registerKeywords( Keywords& keys ) {
122 10 : Colvar::registerKeywords( keys );
123 20 : keys.add("numbered","ATOMS","the atoms involved in each of the contacts you wish to calculate. "
124 : "Keywords like ATOMS1, ATOMS2, ATOMS3,... should be listed and one contact will be "
125 : "calculated for each ATOM keyword you specify.");
126 20 : keys.reset_style("ATOMS","atoms");
127 20 : keys.add("numbered","SWITCH","The switching functions to use for each of the contacts in your map. "
128 : "You can either specify a global switching function using SWITCH or one "
129 : "switching function for each contact. Details of the various switching "
130 : "functions you can use are provided on \\ref switchingfunction.");
131 20 : keys.add("numbered","REFERENCE","A reference value for a given contact, by default is 0.0 "
132 : "You can either specify a global reference value using REFERENCE or one "
133 : "reference value for each contact.");
134 20 : keys.add("numbered","WEIGHT","A weight value for a given contact, by default is 1.0 "
135 : "You can either specify a global weight value using WEIGHT or one "
136 : "weight value for each contact.");
137 20 : keys.reset_style("SWITCH","compulsory");
138 20 : keys.addFlag("SUM",false,"calculate the sum of all the contacts in the input");
139 20 : keys.addFlag("CMDIST",false,"calculate the distance with respect to the provided reference contact map");
140 20 : keys.addFlag("SERIAL",false,"Perform the calculation in serial - for debug purpose");
141 20 : keys.addOutputComponent("contact","default","scalar","By not using SUM or CMDIST each contact will be stored in a component");
142 20 : keys.setValueDescription("scalar","the sum of all the switching function on all the distances");
143 10 : }
144 :
145 8 : ContactMap::ContactMap(const ActionOptions&ao):
146 : PLUMED_COLVAR_INIT(ao),
147 8 : pbc(true),
148 8 : serial(false),
149 8 : docomp(true),
150 8 : dosum(false),
151 8 : docmdist(false)
152 : {
153 8 : parseFlag("SERIAL",serial);
154 8 : parseFlag("SUM",dosum);
155 8 : parseFlag("CMDIST",docmdist);
156 8 : if(docmdist==true&&dosum==true) error("You cannot use SUM and CMDIST together");
157 8 : bool nopbc=!pbc;
158 8 : parseFlag("NOPBC",nopbc);
159 8 : pbc=!nopbc;
160 :
161 : // Read in the atoms
162 : std::vector<AtomNumber> t, ga_lista, gb_lista;
163 8 : for(int i=1;; ++i ) {
164 1186 : parseAtomList("ATOMS", i, t );
165 593 : if( t.empty() ) break;
166 :
167 585 : if( t.size()!=2 ) {
168 0 : std::string ss; Tools::convert(i,ss);
169 0 : error("ATOMS" + ss + " keyword has the wrong number of atoms");
170 : }
171 585 : ga_lista.push_back(t[0]); gb_lista.push_back(t[1]);
172 585 : t.resize(0);
173 :
174 : // Add a value for this contact
175 585 : std::string num; Tools::convert(i,num);
176 595 : if(!dosum&&!docmdist) {addComponentWithDerivatives("contact-"+num); componentIsNotPeriodic("contact-"+num);}
177 585 : }
178 : // Create neighbour lists
179 16 : nl=Tools::make_unique<NeighborList>(ga_lista,gb_lista,serial,true,pbc,getPbc(),comm);
180 :
181 : // Read in switching functions
182 8 : std::string errors; sfs.resize( ga_lista.size() ); unsigned nswitch=0;
183 593 : for(unsigned i=0; i<ga_lista.size(); ++i) {
184 585 : std::string num, sw1; Tools::convert(i+1, num);
185 1170 : if( !parseNumbered( "SWITCH", i+1, sw1 ) ) break;
186 585 : nswitch++; sfs[i].set(sw1,errors);
187 585 : if( errors.length()!=0 ) error("problem reading SWITCH" + num + " keyword : " + errors );
188 : }
189 8 : if( nswitch==0 ) {
190 0 : std::string sw; parse("SWITCH",sw);
191 0 : if(sw.length()==0) error("no switching function specified use SWITCH keyword");
192 0 : for(unsigned i=0; i<ga_lista.size(); ++i) {
193 0 : sfs[i].set(sw,errors);
194 0 : if( errors.length()!=0 ) error("problem reading SWITCH keyword : " + errors );
195 : }
196 8 : } else if( nswitch!=sfs.size() ) {
197 0 : std::string num; Tools::convert(nswitch+1, num);
198 0 : error("missing SWITCH" + num + " keyword");
199 : }
200 :
201 : // Read in reference values
202 : nswitch=0;
203 8 : reference.resize(ga_lista.size(), 0.);
204 18 : for(unsigned i=0; i<ga_lista.size(); ++i) {
205 32 : if( !parseNumbered( "REFERENCE", i+1, reference[i] ) ) break;
206 10 : nswitch++;
207 : }
208 8 : if( nswitch==0 ) {
209 6 : parse("REFERENCE",reference[0]);
210 575 : for(unsigned i=1; i<ga_lista.size(); ++i) {
211 569 : reference[i]=reference[0];
212 569 : nswitch++;
213 : }
214 : }
215 8 : if(nswitch == 0 && docmdist) error("with CMDIST one must use REFERENCE to setup the reference contact map");
216 :
217 : // Read in weight values
218 : nswitch=0;
219 8 : weight.resize(ga_lista.size(), 1.0);
220 588 : for(unsigned i=0; i<ga_lista.size(); ++i) {
221 1162 : if( !parseNumbered( "WEIGHT", i+1, weight[i] ) ) break;
222 580 : nswitch++;
223 : }
224 8 : if( nswitch==0 ) {
225 1 : parse("WEIGHT",weight[0]);
226 5 : for(unsigned i=1; i<ga_lista.size(); ++i) {
227 4 : weight[i]=weight[0];
228 : }
229 : nswitch = ga_lista.size();
230 : }
231 :
232 : // Output details of all contacts
233 593 : for(unsigned i=0; i<sfs.size(); ++i) {
234 585 : log.printf(" The %uth contact is calculated from atoms : %d %d. Inflection point of switching function is at %s. Reference contact value is %f\n",
235 1170 : i+1, ga_lista[i].serial(), gb_lista[i].serial(), ( sfs[i].description() ).c_str(), reference[i] );
236 : }
237 :
238 8 : if(dosum) {
239 10 : addValueWithDerivatives(); setNotPeriodic();
240 5 : log.printf(" colvar is sum of all contacts in contact map\n");
241 : }
242 8 : if(docmdist) {
243 4 : addValueWithDerivatives(); setNotPeriodic();
244 2 : log.printf(" colvar is distance between the contact map matrix and the provided reference matrix\n");
245 : }
246 :
247 8 : if(dosum || docmdist) {
248 7 : docomp=false;
249 : } else {
250 1 : serial=true;
251 1 : docomp=true;
252 : }
253 :
254 : // Set up if it is just a list of contacts
255 8 : requestAtoms(nl->getFullAtomList());
256 8 : checkRead();
257 8 : }
258 :
259 2015 : void ContactMap::calculate() {
260 :
261 2015 : double ncoord=0.;
262 2015 : Tensor virial;
263 2015 : std::vector<Vector> deriv(getNumberOfAtoms());
264 :
265 : unsigned stride;
266 : unsigned rank;
267 2015 : if(serial) {
268 : // when using components the parallelisation do not work
269 : stride=1;
270 : rank=0;
271 : } else {
272 2010 : stride=comm.Get_size();
273 2010 : rank=comm.Get_rank();
274 : }
275 :
276 : // sum over close pairs
277 148635 : for(unsigned i=rank; i<nl->size(); i+=stride) {
278 146620 : Vector distance;
279 146620 : unsigned i0=nl->getClosePair(i).first;
280 146620 : unsigned i1=nl->getClosePair(i).second;
281 146620 : if(pbc) {
282 146620 : distance=pbcDistance(getPosition(i0),getPosition(i1));
283 : } else {
284 0 : distance=delta(getPosition(i0),getPosition(i1));
285 : }
286 :
287 146620 : double dfunc=0.;
288 146620 : double coord = weight[i]*(sfs[i].calculate(distance.modulo(), dfunc) - reference[i]);
289 146620 : Vector tmpder = weight[i]*dfunc*distance;
290 146620 : Tensor tmpvir = weight[i]*dfunc*Tensor(distance,distance);
291 146620 : if(!docmdist) {
292 146595 : deriv[i0] -= tmpder;
293 146595 : deriv[i1] += tmpder;
294 146595 : virial -= tmpvir;
295 146595 : ncoord += coord;
296 : } else {
297 25 : tmpder *= 2.*coord;
298 25 : tmpvir *= 2.*coord;
299 25 : deriv[i0] -= tmpder;
300 25 : deriv[i1] += tmpder;
301 25 : virial -= tmpvir;
302 25 : ncoord += coord*coord;
303 : }
304 :
305 146620 : if(docomp) {
306 25 : Value* val=getPntrToComponent( i );
307 25 : setAtomsDerivatives( val, i0, deriv[i0] );
308 25 : setAtomsDerivatives( val, i1, deriv[i1] );
309 25 : setBoxDerivatives( val, -tmpvir );
310 : val->set(coord);
311 : }
312 : }
313 :
314 2015 : if(!serial) {
315 2010 : comm.Sum(&ncoord,1);
316 2010 : if(!deriv.empty()) comm.Sum(&deriv[0][0],3*deriv.size());
317 2010 : comm.Sum(&virial[0][0],9);
318 : }
319 :
320 2015 : if( !docomp ) {
321 295200 : for(unsigned i=0; i<deriv.size(); ++i) setAtomsDerivatives(i,deriv[i]);
322 2010 : setValue (ncoord);
323 2010 : setBoxDerivatives (virial);
324 : }
325 2015 : }
326 :
327 : }
328 : }
|