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