Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2015-2020 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 "adjmat/AdjacencyMatrixBase.h"
23 : #include "multicolvar/MultiColvarShortcuts.h"
24 : #include "core/ActionRegister.h"
25 : #include "core/ActionShortcut.h"
26 : #include "core/PlumedMain.h"
27 : #include "core/ActionSet.h"
28 : #include "tools/IFile.h"
29 :
30 : //+PLUMEDOC MATRIX HBPAMM_MATRIX
31 : /*
32 : Adjacency matrix in which two electronegative atoms are adjacent if they are hydrogen bonded
33 :
34 : This method allows you to calculate an adjacency matrix and use all the methods discussed in the documentation for
35 : the [CONTACT_MATRIX](CONTACT_MATRIX.md) upon it to define CVs. The $i,j$ element of the matrix that is calculated
36 : by this action is one if there is a hydrogen bond connecting atom $i$ to atom $j$. Furthermore, we determine whether
37 : there is a hydrogen atom between these two atoms by using the PAMM technique that is discussed in the articles from the
38 : bibliography below and in the documentation for the [PAMM](PAMM.md) action.
39 :
40 : The example shown below illustrates how the method is used in practice
41 :
42 : ```plumed
43 : #SETTINGS INPUTFILES=regtest/pamm/rt-hbpamm/b3lyp.pamm
44 : m: HBPAMM_MATRIX GROUP=1-192:3 CLUSTERS=regtest/pamm/rt-hbpamm/b3lyp.pamm GROUPC=2-192:3,3-192:3
45 : PRINT ARG=m FILE=colvar
46 : ```
47 :
48 : This example could be used to investigate the connectivity between a collection of water molecules in liquid water. Notice, however,
49 : that the output matrix here is not symmetric.
50 :
51 : If you want to investigate whether there are hydrogen bonds between two groups of molecules you can use an input like this:
52 :
53 : ```plumed
54 : #SETTINGS INPUTFILES=regtest/pamm/rt-hbpamm/b3lyp.pamm
55 : m: HBPAMM_MATRIX GROUPA=1 GROUPB=2-192:3 CLUSTERS=regtest/pamm/rt-hbpamm/b3lyp.pamm GROUPC=2-192:3,3-192:3
56 : PRINT ARG=m FILE=colvar
57 : ```
58 :
59 : This input outputs a $1\times 63$ matrix in which the $1,i$th element tells you whether or not atom 1 donates a hydrogen bond
60 : to the $i$th element in the group of 63 atoms that was specified using the ACCEPTORS keyword.
61 :
62 : In general, it is better to use this action through the [HBPAMM_SA](HBPAMM_SA.md), [HBPAMM_SD](HBPAMM_SD.md) and [HBPAMM_SH](HBPAMM_SH.md)
63 : keywords, which can be used to calculate the number of hydrogen bonds each donor, acceptor or hydrogen atom in your system participates in.
64 :
65 : */
66 : //+ENDPLUMEDOC
67 :
68 : //+PLUMEDOC COLVAR HBPAMM_SA
69 : /*
70 : Calculate the number of hydrogen bonds each acceptor participates in using the HBPamm method
71 :
72 : This shortcut action allows you to calculate the number of hydrogen bonds each of the atoms specified using the SITES
73 : keyword accepts from its neighbours. The number of hydrogen bonds that a particular site accepts is determined by using the
74 : PAMM tehcnique that is discussed in the articles from the bibliography below and in the documentation for the [PAMM](PAMM.md) action
75 :
76 : The following example shows how you can use this action to calculate how many hydrogen bonds each of the water moecules in
77 : a box of water accepts from its neighbours.
78 :
79 : ```plumed
80 : #SETTINGS INPUTFILES=regtest/pamm/rt-hbpamm/b3lyp.pamm
81 : acceptors: HBPAMM_SA SITES=1-192:3 CLUSTERS=regtest/pamm/rt-hbpamm/b3lyp.pamm HYDROGENS=2-192:3,3-192:3 MEAN
82 : DUMPATOMS ARG=acceptors ATOMS=1-192:3 FILE=acceptors.xyz
83 : ```
84 :
85 : The output here is an xyz with five columns. As explained in the documentation for [DUMPATOMS](DUMPATOMS.md), the first four
86 : columns are the usual columns that you would expect in an xyz file. The fifth column then contains the number of hydrogen bonds
87 : that have been accepted by each of the atoms.
88 :
89 : */
90 : //+ENDPLUMEDOC
91 :
92 : //+PLUMEDOC COLVAR HBPAMM_SD
93 : /*
94 : Calculate the number of hydrogen bonds each donor participates in using the HBPamm method
95 :
96 : This shortcut action allows you to calculate the number of hydrogen bonds each of the atoms specified using the SITES
97 : keyword donates to its neighbours. The number of hydrogen bonds that a particular site donates is determined by using the
98 : PAMM tehcnique that is discussed in the articles from the bibliography below and in the documentation for the [PAMM](PAMM.md) action
99 :
100 : The following example shows how you can use this action to calculate how many hydrogen bonds each of the water moecules in
101 : a box of water donates to its neighbours.
102 :
103 : ```plumed
104 : #SETTINGS INPUTFILES=regtest/pamm/rt-hbpamm/b3lyp.pamm
105 : donors: HBPAMM_SD SITES=1-192:3 CLUSTERS=regtest/pamm/rt-hbpamm/b3lyp.pamm HYDROGENS=2-192:3,3-192:3 MEAN
106 : DUMPATOMS ARG=donors ATOMS=1-192:3 FILE=donors.xyz
107 : ```
108 :
109 : The output here is an xyz with five columns. As explained in the documentation for [DUMPATOMS](DUMPATOMS.md), the first four
110 : columns are the usual columns that you would expect in an xyz file. The fifth column then contains the number of hydrogen bonds
111 : that have been donated by each of the atoms.
112 :
113 : */
114 : //+ENDPLUMEDOC
115 :
116 : //+PLUMEDOC COLVAR HBPAMM_SH
117 : /*
118 : Calculate the number of hydrogen bonds each hydrogen participates in using the HBPamm method
119 :
120 : This shortcut action allows you to calculate the number of hydrogen bonds each of the atoms specified using the SITES
121 : keyword donates to its neighbours. The number of hydrogen bonds that a particular site donates is determined by using the
122 : PAMM tehcnique that is discussed in the articles from the bibliography below and in the documentation for the [PAMM](PAMM.md) action
123 :
124 : The following example shows how you can use this action to calculate how many hydrogen bonds each of the water moecules in
125 : a box of water donates to its neighbours.
126 :
127 : ```plumed
128 : #SETTINGS INPUTFILES=regtest/pamm/rt-hbpamm/b3lyp.pamm
129 : hyd: HBPAMM_SH SITES=1-192:3 CLUSTERS=regtest/pamm/rt-hbpamm/b3lyp.pamm HYDROGENS=2-192:3,3-192:3 MEAN
130 : DUMPATOMS ARG=hyd ATOMS=2-192:3,3-192:3 FILE=hydrogens.xyz
131 : ```
132 :
133 : The output here is an xyz with five columns. As explained in the documentation for [DUMPATOMS](DUMPATOMS.md), the first four
134 : columns are the usual columns that you would expect in an xyz file. The fifth column then contains the number of hydrogen bonds
135 : that each hydrogen atom participates in.
136 :
137 : */
138 : //+ENDPLUMEDOC
139 :
140 : namespace PLMD {
141 : namespace pamm {
142 :
143 : class HBPammMatrix : public adjmat::AdjacencyMatrixBase {
144 : private:
145 : double regulariser;
146 : Tensor incoord_to_hbcoord;
147 : std::vector<double> weight;
148 : std::vector<Vector> centers;
149 : std::vector<Tensor> kmat;
150 : public:
151 : /// Create manual
152 : static void registerKeywords( Keywords& keys );
153 : /// Constructor
154 : explicit HBPammMatrix(const ActionOptions&);
155 : ///
156 : double calculateWeight( const Vector& pos1, const Vector& pos2, const unsigned& natoms, MultiValue& myvals ) const ;
157 : };
158 :
159 : PLUMED_REGISTER_ACTION(HBPammMatrix,"HBPAMM_MATRIX")
160 :
161 31 : void HBPammMatrix::registerKeywords( Keywords& keys ) {
162 31 : adjmat::AdjacencyMatrixBase::registerKeywords( keys );
163 31 : keys.use("GROUPC");
164 31 : keys.add("compulsory","ORDER","dah","the order in which the groups are specified in the input. Can be dah (donor/acceptor/hydrogens), "
165 : "adh (acceptor/donor/hydrogens) or hda (hydrogens/donor/acceptors");
166 31 : keys.add("compulsory","CLUSTERS","the name of the file that contains the definitions of all the kernels for PAMM");
167 31 : keys.add("compulsory","REGULARISE","0.001","don't allow the denominator to be smaller then this value");
168 31 : keys.add("compulsory","GAUSS_CUTOFF","6.25","the cutoff at which to stop evaluating the kernel function is set equal to sqrt(2*x)*(max(adc)+cov(adc))");
169 31 : keys.needsAction("PAMM");
170 31 : keys.needsAction("ONES");
171 31 : keys.needsAction("MATRIX_VECTOR_PRODUCT");
172 31 : keys.addDOI("10.1063/1.4900655");
173 31 : keys.addDOI("10.1021/acs.jctc.7b00993");
174 31 : }
175 :
176 6 : HBPammMatrix::HBPammMatrix(const ActionOptions& ao):
177 : Action(ao),
178 6 : AdjacencyMatrixBase(ao) {
179 : double DP2CUTOFF;
180 12 : parse("GAUSS_CUTOFF",DP2CUTOFF);
181 : std::string sorder;
182 12 : parse("ORDER",sorder);
183 6 : if( sorder=="dah" ) {
184 2 : incoord_to_hbcoord(0,0)=1;
185 2 : incoord_to_hbcoord(0,1)=-1;
186 2 : incoord_to_hbcoord(0,2)=0;
187 2 : incoord_to_hbcoord(1,0)=1;
188 2 : incoord_to_hbcoord(1,1)=1;
189 2 : incoord_to_hbcoord(1,2)=0;
190 2 : incoord_to_hbcoord(2,0)=0;
191 2 : incoord_to_hbcoord(2,1)=0;
192 2 : incoord_to_hbcoord(2,2)=1;
193 2 : log.printf(" GROUPA is list of donor atoms \n");
194 4 : } else if( sorder=="adh" ) {
195 2 : incoord_to_hbcoord(0,0)=-1;
196 2 : incoord_to_hbcoord(0,1)=1;
197 2 : incoord_to_hbcoord(0,2)=0;
198 2 : incoord_to_hbcoord(1,0)=1;
199 2 : incoord_to_hbcoord(1,1)=1;
200 2 : incoord_to_hbcoord(1,2)=0;
201 2 : incoord_to_hbcoord(2,0)=0;
202 2 : incoord_to_hbcoord(2,1)=0;
203 2 : incoord_to_hbcoord(2,2)=1;
204 2 : log.printf(" GROUPA is list of acceptor atoms \n");
205 2 : } else if( sorder=="hda" ) {
206 2 : incoord_to_hbcoord(0,0)=-1;
207 2 : incoord_to_hbcoord(0,1)=0;
208 2 : incoord_to_hbcoord(0,2)=1;
209 2 : incoord_to_hbcoord(1,0)=1;
210 2 : incoord_to_hbcoord(1,1)=0;
211 2 : incoord_to_hbcoord(1,2)=1;
212 2 : incoord_to_hbcoord(2,0)=0;
213 2 : incoord_to_hbcoord(2,1)=1;
214 2 : incoord_to_hbcoord(2,2)=0;
215 2 : log.printf(" GROUPA is list of hydrogen atoms \n");
216 : } else {
217 0 : plumed_error();
218 : }
219 : // Read in the regularisation parameter
220 12 : parse("REGULARISE",regulariser);
221 :
222 : // Read in the kernels
223 : double sqr2pi = sqrt(2*pi);
224 : double sqrt2pi3 = sqr2pi*sqr2pi*sqr2pi;
225 : std::string fname;
226 6 : parse("CLUSTERS", fname);
227 6 : double sfmax=0, ww;
228 6 : Vector cent;
229 6 : Tensor covar;
230 6 : IFile ifile;
231 6 : ifile.open(fname);
232 6 : ifile.allowIgnoredFields();
233 : for(unsigned k=0;; ++k) {
234 144 : if( !ifile.scanField("height",ww) ) {
235 : break;
236 : }
237 66 : ifile.scanField("ptc",cent[0]);
238 66 : ifile.scanField("ssc",cent[1]);
239 66 : ifile.scanField("adc",cent[2]);
240 66 : ifile.scanField("sigma_ptc_ptc",covar[0][0]);
241 66 : ifile.scanField("sigma_ptc_ssc",covar[0][1]);
242 66 : ifile.scanField("sigma_ptc_adc",covar[0][2]);
243 66 : covar[1][0] = covar[0][1];
244 66 : ifile.scanField("sigma_ssc_ssc",covar[1][1]);
245 66 : ifile.scanField("sigma_ssc_adc",covar[1][2]);
246 66 : covar[2][0] = covar[0][2];
247 66 : covar[2][1] = covar[1][2];
248 66 : ifile.scanField("sigma_adc_adc",covar[2][2]);
249 66 : weight.push_back( ww / ( sqrt2pi3 * sqrt(covar.determinant()) ) );
250 66 : centers.push_back( cent );
251 66 : kmat.push_back( covar.inverse() );
252 :
253 66 : Vector eigval;
254 66 : Tensor eigvec;
255 66 : diagMatSym( covar, eigval, eigvec );
256 : unsigned ind_maxeval=0;
257 66 : double max_eval=eigval[0];
258 198 : for(unsigned i=1; i<3; ++i) {
259 132 : if( eigval[i]>max_eval ) {
260 : max_eval=eigval[i];
261 : ind_maxeval=i;
262 : }
263 : }
264 66 : double rcut = cent[2] + sqrt(2.0*DP2CUTOFF)*fabs(sqrt(max_eval)*eigvec(2,ind_maxeval));
265 66 : if( rcut > sfmax ) {
266 24 : sfmax = rcut;
267 : }
268 66 : ifile.scanField();
269 66 : }
270 6 : ifile.close();
271 6 : setLinkCellCutoff( false, sfmax );
272 12 : }
273 :
274 16286 : double HBPammMatrix::calculateWeight( const Vector& pos1, const Vector& pos2, const unsigned& natoms, MultiValue& myvals ) const {
275 16286 : Vector ddij, ddik, ddin, in_dists, hb_pamm_dists, hb_pamm_ders, real_ders;
276 16286 : ddin = pbcDistance( pos1, pos2 );
277 16286 : in_dists[2] = ddin.modulo();
278 16286 : if( in_dists[2]<epsilon ) {
279 : return 0;
280 : }
281 :
282 : double tot=0;
283 16286 : Vector disp, der, tmp_der;
284 1572892 : for(unsigned i=0; i<natoms; ++i) {
285 1556606 : ddij = getPosition(i,myvals);
286 1556606 : in_dists[0] = ddij.modulo();
287 1556606 : ddik = pbcDistance( pos2, getPosition(i,myvals) );
288 1556606 : in_dists[1] = ddik.modulo();
289 1556606 : if( in_dists[1]<epsilon ) {
290 8210 : continue;
291 : }
292 :
293 1548396 : hb_pamm_dists = matmul( incoord_to_hbcoord, in_dists );
294 1548396 : disp = hb_pamm_dists - centers[0];
295 1548396 : der = matmul( kmat[0], disp );
296 1548396 : double vv = weight[0]*exp( -dotProduct( disp, der ) / 2. );
297 1548396 : der *= -vv;
298 :
299 1548396 : double denom = regulariser + vv;
300 6193584 : for(unsigned j=0; j<3; ++j) {
301 4645188 : hb_pamm_ders[j] = der[j];
302 : }
303 17032356 : for(unsigned k=1; k<weight.size(); ++k) {
304 15483960 : disp = hb_pamm_dists - centers[k];
305 15483960 : tmp_der = matmul( kmat[k], disp );
306 15483960 : double tval = weight[k]*exp( -dotProduct( disp, tmp_der ) / 2. );
307 15483960 : denom += tval;
308 15483960 : hb_pamm_ders += -tmp_der*tval;
309 : }
310 1548396 : double vf = vv / denom;
311 1548396 : tot += vf;
312 1548396 : if( fabs(vf)<epsilon ) {
313 1546155 : continue;
314 : }
315 : // Now get derivatives
316 2241 : real_ders = matmul( der / denom - vf*hb_pamm_ders/denom, incoord_to_hbcoord );
317 :
318 : // And add the derivatives to the underlying atoms
319 2241 : addAtomDerivatives( 0, -(real_ders[0]/in_dists[0])*ddij - (real_ders[2]/in_dists[2])*ddin, myvals );
320 2241 : addAtomDerivatives( 1, -(real_ders[1]/in_dists[1])*ddik + (real_ders[2]/in_dists[2])*ddin, myvals );
321 2241 : addThirdAtomDerivatives( i, (real_ders[0]/in_dists[0])*ddij + (real_ders[1]/in_dists[1])*ddik, myvals );
322 4482 : addBoxDerivatives( -(real_ders[0]/in_dists[0])*Tensor( ddij, ddij )
323 4482 : -(real_ders[1]/in_dists[1])*Tensor( ddik, ddik )
324 6723 : -(real_ders[2]/in_dists[2])*Tensor( ddin, ddin ), myvals );
325 : }
326 16286 : return tot;
327 : }
328 :
329 : class HBPammShortcut : public ActionShortcut {
330 : public:
331 : static void registerKeywords( Keywords& keys );
332 : HBPammShortcut(const ActionOptions&);
333 : };
334 :
335 : PLUMED_REGISTER_ACTION(HBPammShortcut,"HBPAMM_SD")
336 : PLUMED_REGISTER_ACTION(HBPammShortcut,"HBPAMM_SA")
337 : PLUMED_REGISTER_ACTION(HBPammShortcut,"HBPAMM_SH")
338 :
339 17 : void HBPammShortcut::registerKeywords( Keywords& keys ) {
340 17 : HBPammMatrix::registerKeywords( keys );
341 17 : keys.remove("GROUP");
342 17 : keys.remove("GROUPA");
343 17 : keys.remove("GROUPB");
344 34 : keys.remove("COMPONENTS");
345 17 : keys.add("optional","SITES","The list of atoms which can be part of a hydrogen bond. When this command is used the set of atoms that can donate a "
346 : "hydrogen bond is assumed to be the same as the set of atoms that can form hydrogen bonds. The atoms involved must be specified"
347 : "as a list of labels of \\ref mcolv or labels of a \\ref multicolvarfunction actions. If you would just like to use "
348 : "the atomic positions you can use a \\ref DENSITY command to specify a group of atoms. Specifying your atomic positions using labels of "
349 : "other \\ref mcolv or \\ref multicolvarfunction commands is useful, however, as you can then exploit a much wider "
350 : "variety of functions of the contact matrix as described in \\ref contactmatrix");
351 17 : keys.add("optional","DONORS","The list of atoms which can donate a hydrogen bond. The atoms involved must be specified "
352 : "as a list of labels of \\ref mcolv or labels of a \\ref multicolvarfunction actions. If you would just like to use "
353 : "the atomic positions you can use a \\ref DENSITY command to specify a group of atoms. Specifying your atomic positions using labels of "
354 : "other \\ref mcolv or \\ref multicolvarfunction commands is useful, however, as you can then exploit a much wider "
355 : "variety of functions of the contact matrix as described in \\ref contactmatrix");
356 17 : keys.add("optional","ACCEPTORS","The list of atoms which can accept a hydrogen bond. The atoms involved must be specified "
357 : "as a list of labels of \\ref mcolv or labels of a \\ref multicolvarfunction actions. If you would just like to use "
358 : "the atomic positions you can use a \\ref DENSITY command to specify a group of atoms. Specifying your atomic positions using labels of "
359 : "other \\ref mcolv or \\ref multicolvarfunction commands is useful, however, as you can then exploit a much wider "
360 : "variety of functions of the contact matrix as described in \\ref contactmatrix");
361 17 : keys.add("compulsory","HYDROGENS","The list of hydrogen atoms that can form part of a hydrogen bond. The atoms must be specified using a comma separated list, "
362 : "an index range or by using a \\ref GROUP");
363 17 : multicolvar::MultiColvarShortcuts::shortcutKeywords( keys );
364 17 : keys.needsAction("HBPAMM_MATRIX");
365 34 : keys.setValueDescription("vector","a vector specifiying the number of hydrogen bonds each of the specified atoms participates within");
366 17 : }
367 :
368 6 : HBPammShortcut::HBPammShortcut(const ActionOptions&ao):
369 : Action(ao),
370 6 : ActionShortcut(ao) {
371 6 : std::string mwords = getShortcutLabel() + "_mat: HBPAMM_MATRIX";
372 6 : if( getName()=="HBPAMM_SD" ) {
373 : std::string site_str;
374 4 : parse("SITES",site_str);
375 2 : if( site_str.length()>0 ) {
376 4 : mwords += " GROUP=" + site_str;
377 : } else {
378 : std::string d_str;
379 0 : parse("DONORS",d_str);
380 0 : mwords += " GROUPA=" + d_str;
381 : std::string a_str;
382 0 : parse("ACCEPTORS",a_str);
383 0 : mwords += " GROUPB=" + a_str;
384 : }
385 : std::string h_str;
386 2 : parse("HYDROGENS",h_str);
387 4 : mwords += " GROUPC=" + h_str + " ORDER=dah";
388 4 : } else if( getName()=="HBPAMM_SA" ) {
389 : std::string site_str;
390 4 : parse("SITES",site_str);
391 2 : if( site_str.length()>0 ) {
392 4 : mwords += " GROUP=" + site_str;
393 : } else {
394 : std::string a_str;
395 0 : parse("ACCEPTORS",a_str);
396 0 : mwords += " GROUPA=" + a_str;
397 : std::string d_str;
398 0 : parse("DONORS",d_str);
399 0 : mwords += " GROUPB=" + d_str;
400 : }
401 : std::string h_str;
402 2 : parse("HYDROGENS",h_str);
403 4 : mwords += " GROUPC=" + h_str + " ORDER=adh";
404 2 : } else if( getName()=="HBPAMM_SH" ) {
405 : std::string h_str;
406 2 : parse("HYDROGENS",h_str);
407 4 : mwords += " GROUPA=" + h_str + " ORDER=hda";
408 : std::string site_str;
409 4 : parse("SITES",site_str);
410 2 : if( site_str.length()>0 ) {
411 2 : mwords += " GROUPB=" + site_str;
412 4 : mwords += " GROUPC=" + site_str;
413 : } else {
414 : std::string d_str;
415 0 : parse("DONORS",d_str);
416 0 : mwords += " GROUPB=" + d_str;
417 : std::string a_str;
418 0 : parse("ACCEPTORS",a_str);
419 0 : mwords += " GROUPC=" + a_str;
420 : }
421 : }
422 : std::map<std::string,std::string> keymap;
423 6 : multicolvar::MultiColvarShortcuts::readShortcutKeywords( keymap, this );
424 12 : readInputLine( mwords + " " + convertInputLineToString() );
425 6 : ActionWithValue* mb=plumed.getActionSet().selectWithLabel<ActionWithValue*>( getShortcutLabel() + "_mat");
426 6 : plumed_assert( mb );
427 : std::string nsize;
428 6 : Tools::convert( (mb->copyOutput(getShortcutLabel() + "_mat"))->getShape()[1], nsize );
429 12 : readInputLine( getShortcutLabel() + "_ones: ONES SIZE=" + nsize );
430 12 : readInputLine( getShortcutLabel() + ": MATRIX_VECTOR_PRODUCT ARG=" + getShortcutLabel() + "_mat," + getShortcutLabel() + "_ones");
431 12 : multicolvar::MultiColvarShortcuts::expandFunctions( getShortcutLabel(), getShortcutLabel(), "", keymap, this );
432 6 : }
433 :
434 : }
435 : }
|