Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2017 of Pipolo Silvio and Fabio Pietrucci.
3 :
4 : The piv module is free software: you can redistribute it and/or modify
5 : it under the terms of the GNU Lesser General Public License as published by
6 : the Free Software Foundation, either version 3 of the License, or
7 : (at your option) any later version.
8 :
9 : The piv module is distributed in the hope that it will be useful,
10 : but WITHOUT ANY WARRANTY; without even the implied warranty of
11 : MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 : GNU Lesser General Public License for more details.
13 :
14 : You should have received a copy of the GNU Lesser General Public License
15 : along with plumed. If not, see <http://www.gnu.org/licenses/>.
16 : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
17 : #include "colvar/Colvar.h"
18 : #include "core/ActionRegister.h"
19 : #include "core/PlumedMain.h"
20 : #include "core/ActionWithVirtualAtom.h"
21 : #include "tools/NeighborList.h"
22 : #include "tools/SwitchingFunction.h"
23 : #include "tools/PDB.h"
24 : #include "tools/Pbc.h"
25 : #include "tools/Communicator.h"
26 : #include "tools/Stopwatch.h"
27 : #include "core/ActionSet.h"
28 :
29 : #include <string>
30 : #include <cmath>
31 : #include <iostream>
32 :
33 : using namespace std;
34 :
35 : namespace PLMD
36 : {
37 : namespace piv
38 : {
39 :
40 : //+PLUMEDOC PIVMOD_COLVAR PIV
41 : /*
42 : Calculates the PIV-distance.
43 :
44 : PIV distance is the squared Cartesian distance between the PIV \cite gallet2013structural \cite pipolo2017navigating
45 : associated to the configuration of the system during the dynamics and a reference configuration provided
46 : as input (PDB file format).
47 : PIV can be used together with \ref FUNCPATHMSD to define a path in the PIV space.
48 :
49 : \par Examples
50 :
51 : The following example calculates PIV-distances from three reference configurations in Ref1.pdb, Ref2.pdb and Ref3.pdb
52 : and prints the results in a file named colvar.
53 : Three atoms (PIVATOMS=3) with names (pdb file) A B and C are used to construct the PIV and all PIV blocks (AA, BB, CC, AB, AC, BC) are considered.
54 : SFACTOR is a scaling factor that multiplies the contribution to the PIV-distance given by the single PIV block.
55 : NLIST sets the use of neighbor lists for calculating atom-atom distances.
56 : The SWITCH keyword specifies the parameters of the switching function that transforms atom-atom distances.
57 : SORT=1 means that the PIV block elements are sorted (SORT=0 no sorting.)
58 : Values for SORT, SFACTOR and the neighbor list parameters have to be specified for each block.
59 : The order is the following: AA,BB,CC,AB,AC,BC. If ONLYDIRECT (ONLYCROSS) is used the order is AA,BB,CC (AB,AC,BC).
60 : The sorting operation within each PIV block is performed using the counting sort algorithm, PRECISION specifies the size of the counting array.
61 :
62 : \plumedfile
63 : PIV ...
64 : LABEL=Pivd1
65 : PRECISION=1000
66 : NLIST
67 : REF_FILE=Ref1.pdb
68 : PIVATOMS=3
69 : ATOMTYPES=A,B,C
70 : SFACTOR=0.3,0.5,1.0,0.2,0.2,0.2
71 : SORT=1,1,1,1,1,1
72 : SWITCH1={RATIONAL R_0=0.6 MM=12 NN=4}
73 : SWITCH2={RATIONAL R_0=0.4 MM=10 NN=5}
74 : SWITCH3={RATIONAL R_0=0.4 MM=10 NN=5}
75 : SWITCH4={RATIONAL R_0=0.5 MM=12 NN=6}
76 : SWITCH5={RATIONAL R_0=0.5 MM=12 NN=6}
77 : SWITCH6={RATIONAL R_0=0.5 MM=12 NN=6}
78 : NL_CUTOFF=0.8,0.6,0.6,0.7,0.7,0.7
79 : NL_STRIDE=10,10,10,10,10,10
80 : NL_SKIN=0.1,0.1,0.1,0.1,0.1,0.1
81 : ... PIV
82 : PIV ...
83 : LABEL=Pivd2
84 : PRECISION=1000
85 : NLIST
86 : REF_FILE=Ref2.pdb
87 : PIVATOMS=3
88 : ATOMTYPES=A,B,C
89 : SFACTOR=0.3,0.5,1.0,0.2,0.2,0.2
90 : SORT=1,1,1,1,1,1
91 : SWITCH1={RATIONAL R_0=0.6 MM=12 NN=4}
92 : SWITCH2={RATIONAL R_0=0.4 MM=10 NN=5}
93 : SWITCH3={RATIONAL R_0=0.4 MM=10 NN=5}
94 : SWITCH4={RATIONAL R_0=0.5 MM=12 NN=6}
95 : SWITCH5={RATIONAL R_0=0.5 MM=12 NN=6}
96 : SWITCH6={RATIONAL R_0=0.5 MM=12 NN=6}
97 : NL_CUTOFF=0.8,0.6,0.6,0.7,0.7,0.7
98 : NL_STRIDE=10,10,10,10,10,10
99 : NL_SKIN=0.1,0.1,0.1,0.1,0.1,0.1
100 : ... PIV
101 : PIV ...
102 : LABEL=Pivd3
103 : PRECISION=1000
104 : NLIST
105 : REF_FILE=Ref3.pdb
106 : PIVATOMS=3
107 : ATOMTYPES=A,B,C
108 : SFACTOR=0.3,0.5,1.0,0.2,0.2,0.2
109 : SORT=1,1,1,1,1,1
110 : SWITCH1={RATIONAL R_0=0.6 MM=12 NN=4}
111 : SWITCH2={RATIONAL R_0=0.4 MM=10 NN=5}
112 : SWITCH3={RATIONAL R_0=0.4 MM=10 NN=5}
113 : SWITCH4={RATIONAL R_0=0.5 MM=12 NN=6}
114 : SWITCH5={RATIONAL R_0=0.5 MM=12 NN=6}
115 : SWITCH6={RATIONAL R_0=0.5 MM=12 NN=6}
116 : NL_CUTOFF=0.8,0.6,0.6,0.7,0.7,0.7
117 : NL_STRIDE=10,10,10,10,10,10
118 : NL_SKIN=0.1,0.1,0.1,0.1,0.1,0.1
119 : ... PIV
120 :
121 : PRINT ARG=Pivd1,Pivd2,Pivd3 FILE=colvar
122 : \endplumedfile
123 :
124 : WARNING:
125 : Both the "CRYST" and "ATOM" lines of the PDB files must conform precisely to the official pdb format, including the width of each alphanumerical field:
126 :
127 : \verbatim
128 : CRYST1 31.028 36.957 23.143 89.93 92.31 89.99 P 1 1
129 : ATOM 1 OW1 wate 1 15.630 19.750 1.520 1.00 0.00
130 : \endverbatim
131 :
132 : In each pdb frame, atoms must be numbered in the same order and with the same element symbol as in the input of the MD program.
133 :
134 : The following example calculates the PIV-distances from two reference configurations Ref1.pdb and Ref2.pdb
135 : and uses PIV-distances to define a Path Collective Variable (\ref FUNCPATHMSD) with only two references (Ref1.pdb and Ref2.pdb).
136 : With the VOLUME keyword one scales the atom-atom distances by the cubic root of the ratio between the specified value and the box volume of the initial step of the trajectory file.
137 :
138 : \plumedfile
139 : PIV ...
140 : LABEL=c1
141 : PRECISION=1000
142 : VOLUME=12.15
143 : NLIST
144 : REF_FILE=Ref1.pdb
145 : PIVATOMS=2
146 : ATOMTYPES=A,B
147 : ONLYDIRECT
148 : SFACTOR=1.0,0.2
149 : SORT=1,1
150 : SWITCH1={RATIONAL R_0=0.6 MM=12 NN=4}
151 : SWITCH2={RATIONAL R_0=0.5 MM=10 NN=5}
152 : NL_CUTOFF=1.2,1.2
153 : NL_STRIDE=10,10
154 : NL_SKIN=0.1,0.1
155 : ... PIV
156 : PIV ...
157 : LABEL=c2
158 : PRECISION=1000
159 : VOLUME=12.15
160 : NLIST
161 : REF_FILE=Ref2.pdb
162 : PIVATOMS=2
163 : ATOMTYPES=A,B
164 : ONLYDIRECT
165 : SFACTOR=1.0,0.2
166 : SORT=1,1
167 : SWITCH1={RATIONAL R_0=0.6 MM=12 NN=4}
168 : SWITCH2={RATIONAL R_0=0.4 MM=10 NN=5}
169 : NL_CUTOFF=1.2,1.2
170 : NL_STRIDE=10,10
171 : NL_SKIN=0.1,0.1
172 : ... PIV
173 :
174 : p1: FUNCPATHMSD ARG=c1,c2 LAMBDA=0.180338
175 : METAD ARG=p1.s,p1.z SIGMA=0.01,0.2 HEIGHT=0.8 PACE=500 LABEL=res
176 : PRINT ARG=c1,c2,p1.s,p1.z,res.bias STRIDE=500 FILE=colvar FMT=%15.6f
177 : \endplumedfile
178 :
179 : When using PIV please cite \cite pipolo2017navigating .
180 :
181 : (See also \ref PRINT)
182 :
183 : */
184 : //+ENDPLUMEDOC
185 :
186 : class PIV : public Colvar
187 : {
188 : private:
189 : bool pbc, serial, timer;
190 : ForwardDecl<Stopwatch> stopwatch_fwd;
191 : Stopwatch& stopwatch=*stopwatch_fwd;
192 : int updatePIV;
193 : size_t Nprec;
194 : unsigned Natm,Nlist,NLsize;
195 : double Fvol,Vol0,m_PIVdistance;
196 : std::string ref_file;
197 : std::unique_ptr<NeighborList> nlall;
198 : std::vector<SwitchingFunction> sfs;
199 : std::vector<std:: vector<double> > rPIV;
200 : std::vector<double> scaling,r00;
201 : std::vector<double> nl_skin;
202 : std::vector<double> fmass;
203 : std::vector<bool> dosort;
204 : std::vector<Vector> compos;
205 : std::vector<string> sw;
206 : std::vector<std::unique_ptr<NeighborList>> nl;
207 : std::vector<std::unique_ptr<NeighborList>> nlcom;
208 : std::vector<Vector> m_deriv;
209 : Tensor m_virial;
210 : bool Svol,cross,direct,doneigh,test,CompDer,com;
211 :
212 : /// Local structure, used to store data that should be
213 : /// shared across multiple PIV instances
214 6 : struct SharedData {
215 : int prev_stp=-1;
216 : int init_stp=1;
217 : std:: vector<std:: vector<Vector> > prev_pos;
218 : std:: vector<std:: vector<double> > cPIV;
219 : std:: vector<std:: vector<int> > Atom0;
220 : std:: vector<std:: vector<int> > Atom1;
221 : };
222 : /// Owning pointer. Will only be allocated by the first PIV instance
223 : std::unique_ptr<SharedData> sharedData_unique;
224 : /// Raw pointer. Will have the same value for all PIV instances
225 : SharedData* sharedData=nullptr;
226 :
227 : public:
228 : static void registerKeywords( Keywords& keys );
229 : explicit PIV(const ActionOptions&);
230 : // active methods:
231 : virtual void calculate();
232 0 : void checkFieldsAllowed() {}
233 : };
234 :
235 : PLUMED_REGISTER_ACTION(PIV,"PIV")
236 :
237 14 : void PIV::registerKeywords( Keywords& keys )
238 : {
239 14 : Colvar::registerKeywords( keys );
240 28 : keys.add("numbered","SWITCH","The switching functions parameter."
241 : "You should specify a Switching function for all PIV blocks."
242 : "Details of the various switching "
243 : "functions you can use are provided on switchingfunction.");
244 28 : keys.add("compulsory","PRECISION","the precision for approximating reals with integers in sorting.");
245 28 : keys.add("compulsory","REF_FILE","PDB file name that contains the ith reference structure.");
246 28 : keys.add("compulsory","PIVATOMS","Number of atoms to use for PIV.");
247 28 : keys.add("compulsory","SORT","Whether to sort or not the PIV block.");
248 28 : keys.add("compulsory","ATOMTYPES","The atom types to use for PIV.");
249 28 : keys.add("optional","SFACTOR","Scale the PIV-distance by such block-specific factor");
250 28 : keys.add("optional","VOLUME","Scale atom-atom distances by the cubic root of the cell volume. The input volume is used to scale the R_0 value of the switching function. ");
251 28 : keys.add("optional","UPDATEPIV","Frequency (in steps) at which the PIV is updated.");
252 28 : keys.addFlag("TEST",false,"Print the actual and reference PIV and exit");
253 28 : keys.addFlag("COM",false,"Use centers of mass of groups of atoms instead of atoms as specified in the Pdb file");
254 28 : keys.addFlag("ONLYCROSS",false,"Use only cross-terms (A-B, A-C, B-C, ...) in PIV");
255 28 : keys.addFlag("ONLYDIRECT",false,"Use only direct-terms (A-A, B-B, C-C, ...) in PIV");
256 28 : keys.addFlag("DERIVATIVES",false,"Activate the calculation of the PIV for every class (needed for numerical derivatives).");
257 28 : keys.addFlag("NLIST",false,"Use a neighbor list for distance calculations.");
258 28 : keys.addFlag("SERIAL",false,"Perform the calculation in serial - for debug purpose");
259 28 : keys.addFlag("TIMER",false,"Perform timing analysis on heavy loops.");
260 28 : keys.add("optional","NL_CUTOFF","Neighbor lists cutoff.");
261 28 : keys.add("optional","NL_STRIDE","Update neighbor lists every NL_STRIDE steps.");
262 28 : keys.add("optional","NL_SKIN","The maximum atom displacement tolerated for the neighbor lists update.");
263 28 : keys.reset_style("SWITCH","compulsory");
264 28 : keys.setValueDescription("scalar","the PIV-distance");
265 14 : }
266 :
267 12 : PIV::PIV(const ActionOptions&ao):
268 : PLUMED_COLVAR_INIT(ao),
269 12 : pbc(true),
270 12 : serial(false),
271 12 : timer(false),
272 12 : updatePIV(1),
273 12 : Nprec(1000),
274 12 : Natm(1),
275 12 : Nlist(1),
276 12 : NLsize(1),
277 12 : Fvol(1.),
278 12 : Vol0(0.),
279 12 : m_PIVdistance(0.),
280 12 : rPIV(std:: vector<std:: vector<double> >(Nlist)),
281 12 : scaling(std:: vector<double>(Nlist)),
282 12 : r00(std:: vector<double>(Nlist)),
283 12 : nl_skin(std:: vector<double>(Nlist)),
284 12 : fmass(std:: vector<double>(Nlist)),
285 12 : dosort(std:: vector<bool>(Nlist)),
286 12 : compos(std:: vector<Vector>(NLsize)),
287 12 : sw(std:: vector<string>(Nlist)),
288 12 : nl(Nlist),
289 12 : nlcom(NLsize),
290 12 : m_deriv(std:: vector<Vector>(1)),
291 12 : Svol(false),
292 12 : cross(true),
293 12 : direct(true),
294 12 : doneigh(false),
295 12 : test(false),
296 12 : CompDer(false),
297 24 : com(false)
298 : {
299 12 : log << "Starting PIV Constructor\n";
300 :
301 : {
302 : // look for another PIV instance previously allocated
303 12 : auto* previous=plumed.getActionSet().selectLatest<PIV*>(this);
304 :
305 : // Uncommenting the following line, it is possible to force
306 : // a separate object per instance of the PIB object.
307 : // Results are unaffected, but performance is worse.
308 : // I think this is the expected behavior. GB
309 : // previous=nullptr;
310 :
311 12 : if(!previous) {
312 : // if not found, allocate the shared data struct
313 6 : sharedData_unique=Tools::make_unique<SharedData>();
314 : // then set the raw pointer
315 6 : sharedData=sharedData_unique.get();
316 : } else {
317 : // if found, use the previous raw pointer
318 6 : sharedData=previous->sharedData;
319 6 : log << "(a previous PIV action was found)\n";
320 : }
321 : }
322 :
323 : // Precision on the real-to-integer transformation for the sorting
324 12 : parse("PRECISION",Nprec);
325 12 : if(Nprec<2) error("Precision must be => 2");
326 :
327 : // PBC
328 12 : bool nopbc=!pbc;
329 12 : parseFlag("NOPBC",nopbc);
330 12 : pbc=!nopbc;
331 12 : if(pbc) {
332 12 : log << "Using Periodic Boundary Conditions\n";
333 : } else {
334 0 : log << "Isolated System (NO PBC)\n";
335 : }
336 :
337 : // SERIAL/PARALLEL
338 12 : parseFlag("SERIAL",serial);
339 12 : if(serial) {
340 0 : log << "Serial PIV construction\n";
341 : } else {
342 12 : log << "Parallel PIV construction\n";
343 : }
344 :
345 : // Derivatives
346 12 : parseFlag("DERIVATIVES",CompDer);
347 12 : if(CompDer) log << "Computing Derivatives\n";
348 :
349 : // Timing
350 12 : parseFlag("TIMER",timer);
351 12 : if(timer) {
352 1 : log << "Timing analysis\n";
353 1 : stopwatch.start();
354 1 : stopwatch.pause();
355 : }
356 :
357 : // Test
358 12 : parseFlag("TEST",test);
359 :
360 : // UPDATEPIV
361 24 : if(keywords.exists("UPDATEPIV")) {
362 24 : parse("UPDATEPIV",updatePIV);
363 : }
364 :
365 : // Test
366 12 : parseFlag("COM",com);
367 12 : if(com) log << "Building PIV using COMs\n";
368 :
369 : // Volume Scaling
370 12 : parse("VOLUME",Vol0);
371 12 : if (Vol0>0) {
372 12 : Svol=true;
373 : }
374 :
375 : // PIV direct and cross blocks
376 12 : bool oc=false,od=false;
377 12 : parseFlag("ONLYCROSS",oc);
378 12 : parseFlag("ONLYDIRECT",od);
379 12 : if (oc&&od) {
380 0 : error("ONLYCROSS and ONLYDIRECT are incompatible options!");
381 : }
382 12 : if(oc) {
383 4 : direct=false;
384 4 : log << "Using only CROSS-PIV blocks\n";
385 : }
386 12 : if(od) {
387 4 : cross=false;
388 4 : log << "Using only DIRECT-PIV blocks\n";
389 : }
390 :
391 : // Atoms for PIV
392 12 : parse("PIVATOMS",Natm);
393 12 : std:: vector<string> atype(Natm);
394 12 : parseVector("ATOMTYPES",atype);
395 : //if(atype.size()!=getNumberOfArguments() && atype.size()!=0) error("not enough values for ATOMTYPES");
396 :
397 : // Reference PDB file
398 12 : parse("REF_FILE",ref_file);
399 12 : PDB mypdb;
400 12 : FILE* fp=fopen(ref_file.c_str(),"r");
401 12 : if (fp!=NULL) {
402 12 : log<<"Opening PDB file with reference frame: "<<ref_file.c_str()<<"\n";
403 12 : mypdb.readFromFilepointer(fp,usingNaturalUnits(),0.1/getUnits().getLength());
404 12 : fclose (fp);
405 : } else {
406 0 : error("Error in reference PDB file");
407 : }
408 :
409 : // Build COM/Atom lists of AtomNumbers (this might be done in PBC.cpp)
410 : // Atomlist or Plist used to build pair lists
411 12 : std:: vector<std:: vector<AtomNumber> > Plist(Natm);
412 : // Atomlist used to build list of atoms for each COM
413 12 : std:: vector<std:: vector<AtomNumber> > comatm(1);
414 : // NLsize is the number of atoms in the pdb cell
415 12 : NLsize=mypdb.getAtomNumbers().size();
416 : // In the following P stands for Point (either an Atom or a COM)
417 : unsigned resnum=0;
418 : // Presind (array size: number of residues) contains the contains the residue number
419 : // this is because the residue numbers may not always be ordered from 1 to resnum
420 : std:: vector<unsigned> Presind;
421 : // Build Presind
422 19404 : for (unsigned i=0; i<mypdb.getAtomNumbers().size(); i++) {
423 19392 : unsigned rind=mypdb.getResidueNumber(mypdb.getAtomNumbers()[i]);
424 : bool oldres=false;
425 7705008 : for (unsigned j=0; j<Presind.size(); j++) {
426 7685616 : if(rind==Presind[j]) {
427 : oldres=true;
428 : }
429 : }
430 19392 : if(!oldres) {
431 4848 : Presind.push_back(rind);
432 : }
433 : }
434 12 : resnum=Presind.size();
435 :
436 : // Pind0 is the atom/COM used in Nlists (for COM Pind0 is the first atom in the pdb belonging to that COM)
437 : unsigned Pind0size;
438 12 : if(com) {
439 : Pind0size=resnum;
440 : } else {
441 12 : Pind0size=NLsize;
442 : }
443 12 : std:: vector<unsigned> Pind0(Pind0size);
444 : // If COM resize important arrays
445 12 : comatm.resize(NLsize);
446 12 : if(com) {
447 0 : nlcom.resize(NLsize);
448 0 : compos.resize(NLsize);
449 0 : fmass.resize(NLsize,0.);
450 : }
451 12 : log << "Total COM/Atoms: " << Natm*resnum << " \n";
452 : // Build lists of Atoms/COMs for NLists
453 : // comatm filled also for non_COM calculation for analysis purposes
454 36 : for (unsigned j=0; j<Natm; j++) {
455 : unsigned oind;
456 38808 : for (unsigned i=0; i<Pind0.size(); i++) {
457 38784 : Pind0[i]=0;
458 : }
459 38808 : for (unsigned i=0; i<mypdb.getAtomNumbers().size(); i++) {
460 : // Residue/Atom AtomNumber: used to build NL for COMS/Atoms pairs.
461 38784 : AtomNumber anum=mypdb.getAtomNumbers()[i];
462 : // ResidueName/Atomname associated to atom
463 38784 : string rname=mypdb.getResidueName(anum);
464 38784 : string aname=mypdb.getAtomName(anum);
465 : // Index associated to residue/atom: used to separate COM-lists
466 38784 : unsigned rind=mypdb.getResidueNumber(anum);
467 : unsigned aind=anum.index();
468 : // This builds lists for NL
469 : string Pname;
470 : unsigned Pind;
471 38784 : if(com) {
472 : Pname=rname;
473 0 : for(unsigned l=0; l<resnum; l++) {
474 0 : if(rind==Presind[l]) {
475 : Pind=l;
476 : }
477 : }
478 : } else {
479 : Pname=aname;
480 : Pind=aind;
481 : }
482 38784 : if(Pname==atype[j]) {
483 14544 : if(Pind0[Pind]==0) {
484 : // adding the atomnumber to the atom/COM list for pairs
485 14544 : Plist[j].push_back(anum);
486 14544 : Pind0[Pind]=aind+1;
487 : oind=Pind;
488 : }
489 : // adding the atomnumber to list of atoms for every COM/Atoms
490 14544 : comatm[Pind0[Pind]-1].push_back(anum);
491 : }
492 : }
493 : // Output Lists
494 24 : log << " Groups of type " << j << ": " << Plist[j].size() << " \n";
495 : string gname;
496 : unsigned gsize;
497 24 : if(com) {
498 0 : gname=mypdb.getResidueName(comatm[Pind0[oind]-1][0]);
499 0 : gsize=comatm[Pind0[oind]-1].size();
500 : } else {
501 48 : gname=mypdb.getAtomName(comatm[Pind0[oind]-1][0]);
502 : gsize=1;
503 : }
504 24 : log.printf(" %6s %3s %13s %10i %6s\n", "type ", gname.c_str()," containing ",gsize," atoms");
505 : }
506 :
507 : // This is to build the list with all the atoms
508 : std:: vector<AtomNumber> listall;
509 19404 : for (unsigned i=0; i<mypdb.getAtomNumbers().size(); i++) {
510 19392 : listall.push_back(mypdb.getAtomNumbers()[i]);
511 : }
512 :
513 : // PIV blocks and Neighbour Lists
514 12 : Nlist=0;
515 : // Direct adds the A-A ad B-B blocks (N)
516 12 : if(direct) {
517 8 : Nlist=Nlist+unsigned(Natm);
518 : }
519 : // Cross adds the A-B blocks (N*(N-1)/2)
520 12 : if(cross) {
521 8 : Nlist=Nlist+unsigned(double(Natm*(Natm-1))/2.);
522 : }
523 : // Resize vectors according to Nlist
524 12 : rPIV.resize(Nlist);
525 :
526 : // PIV scaled option
527 12 : scaling.resize(Nlist);
528 36 : for(unsigned j=0; j<Nlist; j++) {
529 24 : scaling[j]=1.;
530 : }
531 24 : if(keywords.exists("SFACTOR")) {
532 24 : parseVector("SFACTOR",scaling);
533 : //if(scaling.size()!=getNumberOfArguments() && scaling.size()!=0) error("not enough values for SFACTOR");
534 : }
535 : // Neighbour Lists option
536 12 : parseFlag("NLIST",doneigh);
537 12 : nl.resize(Nlist);
538 12 : nl_skin.resize(Nlist);
539 12 : if(doneigh) {
540 12 : std:: vector<double> nl_cut(Nlist,0.);
541 12 : std:: vector<int> nl_st(Nlist,0);
542 12 : parseVector("NL_CUTOFF",nl_cut);
543 : //if(nl_cut.size()!=getNumberOfArguments() && nl_cut.size()!=0) error("not enough values for NL_CUTOFF");
544 12 : parseVector("NL_STRIDE",nl_st);
545 : //if(nl_st.size()!=getNumberOfArguments() && nl_st.size()!=0) error("not enough values for NL_STRIDE");
546 12 : parseVector("NL_SKIN",nl_skin);
547 : //if(nl_skin.size()!=getNumberOfArguments() && nl_skin.size()!=0) error("not enough values for NL_SKIN");
548 36 : for (unsigned j=0; j<Nlist; j++) {
549 24 : if(nl_cut[j]<=0.0) error("NL_CUTOFF should be explicitly specified and positive");
550 24 : if(nl_st[j]<=0) error("NL_STRIDE should be explicitly specified and positive");
551 24 : if(nl_skin[j]<=0.) error("NL_SKIN should be explicitly specified and positive");
552 24 : nl_cut[j]=nl_cut[j]+nl_skin[j];
553 : }
554 12 : log << "Creating Neighbor Lists \n";
555 : // WARNING: is nl_cut meaningful here?
556 24 : nlall= Tools::make_unique<NeighborList>(listall,true,pbc,getPbc(),comm,nl_cut[0],nl_st[0]);
557 12 : if(com) {
558 : //Build lists of Atoms for every COM
559 0 : for (unsigned i=0; i<compos.size(); i++) {
560 : // WARNING: is nl_cut meaningful here?
561 0 : nlcom[i] = Tools::make_unique<NeighborList>(comatm[i],true,pbc,getPbc(),comm,nl_cut[0],nl_st[0]);
562 : }
563 : }
564 : unsigned ncnt=0;
565 : // Direct blocks AA, BB, CC, ...
566 12 : if(direct) {
567 24 : for (unsigned j=0; j<Natm; j++) {
568 16 : nl[ncnt]= Tools::make_unique<NeighborList>(Plist[j],true,pbc,getPbc(),comm,nl_cut[j],nl_st[j]);
569 16 : ncnt+=1;
570 : }
571 : }
572 : // Cross blocks AB, AC, BC, ...
573 12 : if(cross) {
574 24 : for (unsigned j=0; j<Natm; j++) {
575 24 : for (unsigned i=j+1; i<Natm; i++) {
576 16 : nl[ncnt]= Tools::make_unique<NeighborList>(Plist[i],Plist[j],true,false,pbc,getPbc(),comm,nl_cut[ncnt],nl_st[ncnt]);
577 8 : ncnt+=1;
578 : }
579 : }
580 : }
581 : } else {
582 0 : log << "WARNING: Neighbor List not activated this has not been tested!! \n";
583 0 : nlall= Tools::make_unique<NeighborList>(listall,true,pbc,getPbc(),comm);
584 0 : for (unsigned j=0; j<Nlist; j++) {
585 0 : nl[j]= Tools::make_unique<NeighborList>(Plist[j],Plist[j],true,true,pbc,getPbc(),comm);
586 : }
587 : }
588 : // Output Nlist
589 12 : log << "Total Nlists: " << Nlist << " \n";
590 36 : for (unsigned j=0; j<Nlist; j++) {
591 24 : log << " list " << j+1 << " size " << nl[j]->size() << " \n";
592 : }
593 : // Calculate COM masses once and for all from lists
594 12 : if(com) {
595 0 : for(unsigned j=0; j<compos.size(); j++) {
596 : double commass=0.;
597 0 : for(unsigned i=0; i<nlcom[j]->getFullAtomList().size(); i++) {
598 0 : unsigned andx=nlcom[j]->getFullAtomList()[i].index();
599 0 : commass+=mypdb.getOccupancy()[andx];
600 : }
601 0 : for(unsigned i=0; i<nlcom[j]->getFullAtomList().size(); i++) {
602 0 : unsigned andx=nlcom[j]->getFullAtomList()[i].index();
603 0 : if(commass>0.) {
604 0 : fmass[andx]=mypdb.getOccupancy()[andx]/commass;
605 : } else {
606 0 : fmass[andx]=1.;
607 : }
608 : }
609 : }
610 : }
611 :
612 : // Sorting
613 12 : dosort.resize(Nlist);
614 12 : std:: vector<int> ynsort(Nlist);
615 12 : parseVector("SORT",ynsort);
616 36 : for (unsigned i=0; i<Nlist; i++) {
617 24 : if(ynsort[i]==0||CompDer) {
618 : dosort[i]=false;
619 : } else {
620 : dosort[i]=true;
621 : }
622 : }
623 :
624 : //build box vectors and correct for pbc
625 12 : log << "Building the box from PDB data ... \n";
626 12 : Tensor Box=mypdb.getBoxVec();
627 12 : log << " Done! A,B,C vectors in Cartesian space: \n";
628 12 : log.printf(" A: %12.6f%12.6f%12.6f\n", Box[0][0],Box[0][1],Box[0][2]);
629 12 : log.printf(" B: %12.6f%12.6f%12.6f\n", Box[1][0],Box[1][1],Box[1][2]);
630 12 : log.printf(" C: %12.6f%12.6f%12.6f\n", Box[2][0],Box[2][1],Box[2][2]);
631 12 : log << "Changing the PBC according to the new box \n";
632 12 : Pbc mypbc;
633 12 : mypbc.setBox(Box);
634 12 : log << "The box volume is " << mypbc.getBox().determinant() << " \n";
635 :
636 : //Compute scaling factor
637 12 : if(Svol) {
638 12 : Fvol=cbrt(Vol0/mypbc.getBox().determinant());
639 12 : log << "Scaling atom distances by " << Fvol << " \n";
640 : } else {
641 0 : log << "Using unscaled atom distances \n";
642 : }
643 :
644 12 : r00.resize(Nlist);
645 12 : sw.resize(Nlist);
646 36 : for (unsigned j=0; j<Nlist; j++) {
647 48 : if( !parseNumbered( "SWITCH", j+1, sw[j] ) ) break;
648 : }
649 12 : if(CompDer) {
650 : // Set switching function parameters here only if computing derivatives
651 : // now set at the beginning of the dynamics to solve the r0 issue
652 6 : log << "Switching Function Parameters \n";
653 6 : sfs.resize(Nlist);
654 : std::string errors;
655 18 : for (unsigned j=0; j<Nlist; j++) {
656 12 : if(Svol) {
657 : double r0;
658 : std::string old_r0;
659 12 : vector<string> data=Tools::getWords(sw[j]);
660 : data.erase(data.begin());
661 12 : Tools::parse(data,"R_0",old_r0);
662 12 : Tools::convert(old_r0,r0);
663 12 : r0*=Fvol;
664 12 : std::string new_r0; Tools::convert(r0,new_r0);
665 12 : std::size_t pos = sw[j].find("R_0");
666 12 : sw[j].replace(pos+4,old_r0.size(),new_r0);
667 12 : }
668 12 : sfs[j].set(sw[j],errors);
669 : std::string num;
670 12 : Tools::convert(j+1, num);
671 12 : if( errors.length()!=0 ) error("problem reading SWITCH" + num + " keyword : " + errors );
672 12 : r00[j]=sfs[j].get_r0();
673 12 : log << " Swf: " << j << " r0=" << (sfs[j].description()).c_str() << " \n";
674 : }
675 : }
676 :
677 : // build COMs from positions if requested
678 12 : if(com) {
679 0 : for(unsigned j=0; j<compos.size(); j++) {
680 0 : compos[j][0]=0.;
681 0 : compos[j][1]=0.;
682 0 : compos[j][2]=0.;
683 0 : for(unsigned i=0; i<nlcom[j]->getFullAtomList().size(); i++) {
684 0 : unsigned andx=nlcom[j]->getFullAtomList()[i].index();
685 0 : compos[j]+=fmass[andx]*mypdb.getPositions()[andx];
686 : }
687 : }
688 : }
689 : // build the rPIV distances (transformation and sorting is done afterwards)
690 12 : if(CompDer) {
691 6 : log << " PIV | block | Size | Zeros | Ones |" << " \n";
692 : }
693 36 : for(unsigned j=0; j<Nlist; j++) {
694 11516328 : for(unsigned i=0; i<nl[j]->size(); i++) {
695 11516304 : unsigned i0=(nl[j]->getClosePairAtomNumber(i).first).index();
696 11516304 : unsigned i1=(nl[j]->getClosePairAtomNumber(i).second).index();
697 : //calculate/get COM position of centers i0 and i1
698 11516304 : Vector Pos0,Pos1;
699 11516304 : if(com) {
700 : //if(pbc) makeWhole();
701 0 : Pos0=compos[i0];
702 0 : Pos1=compos[i1];
703 : } else {
704 11516304 : Pos0=mypdb.getPositions()[i0];
705 11516304 : Pos1=mypdb.getPositions()[i1];
706 : }
707 11516304 : Vector ddist;
708 11516304 : if(pbc) {
709 11516304 : ddist=mypbc.distance(Pos0,Pos1);
710 : } else {
711 0 : ddist=delta(Pos0,Pos1);
712 : }
713 11516304 : double df=0.;
714 : // Transformation and sorting done at the first timestep to solve the r0 definition issue
715 11516304 : if(CompDer) {
716 1104 : rPIV[j].push_back(sfs[j].calculate(ddist.modulo()*Fvol, df));
717 : } else {
718 11515200 : rPIV[j].push_back(ddist.modulo()*Fvol);
719 : }
720 : }
721 24 : if(CompDer) {
722 12 : if(dosort[j]) {
723 0 : std::sort(rPIV[j].begin(),rPIV[j].end());
724 : }
725 : int lmt0=0;
726 : int lmt1=0;
727 1116 : for(unsigned i=0; i<rPIV[j].size(); i++) {
728 1104 : if(int(rPIV[j][i]*double(Nprec-1))==0) {
729 0 : lmt0+=1;
730 : }
731 1104 : if(int(rPIV[j][i]*double(Nprec-1))==1) {
732 0 : lmt1+=1;
733 : }
734 : }
735 12 : log.printf(" |%10i|%15zu|%15i|%15i|\n", j, rPIV[j].size(), lmt0, lmt1);
736 : }
737 : }
738 :
739 12 : checkRead();
740 : // From the plumed manual on how to build-up a new Colvar
741 24 : addValueWithDerivatives();
742 12 : requestAtoms(nlall->getFullAtomList());
743 12 : setNotPeriodic();
744 : // getValue()->setPeridodicity(false);
745 : // set size of derivative vector
746 12 : m_deriv.resize(getNumberOfAtoms());
747 12 : }
748 :
749 327 : void PIV::calculate()
750 : {
751 :
752 : // Local variables
753 :
754 327 : if(sharedData_unique) {
755 : // This is executed by the first PIV instance.
756 : // We initialize variables with the correct Nlist.
757 321 : sharedData_unique->prev_pos.resize(Nlist);
758 321 : sharedData_unique->cPIV.resize(Nlist);
759 321 : sharedData_unique->Atom0.resize(Nlist);
760 321 : sharedData_unique->Atom1.resize(Nlist);
761 : }
762 :
763 : // create references to minimize the impact of the code below
764 327 : auto & prev_stp(sharedData->prev_stp);
765 : auto & init_stp(sharedData->init_stp);
766 : auto & prev_pos(sharedData->prev_pos);
767 : auto & cPIV(sharedData->cPIV);
768 : auto & Atom0(sharedData->Atom0);
769 : auto & Atom1(sharedData->Atom1);
770 :
771 327 : std:: vector<std:: vector<int> > A0(Nprec);
772 327 : std:: vector<std:: vector<int> > A1(Nprec);
773 : size_t stride=1;
774 : unsigned rank=0;
775 :
776 327 : if(!serial) {
777 327 : stride=comm.Get_size();
778 327 : rank=comm.Get_rank();
779 : } else {
780 : stride=1;
781 : rank=0;
782 : }
783 :
784 : // Transform (and sort) the rPIV before starting the dynamics
785 327 : if (((prev_stp==-1) || (init_stp==1)) &&!CompDer) {
786 6 : if(prev_stp!=-1) {init_stp=0;}
787 : // Calculate the volume scaling factor
788 6 : if(Svol) {
789 6 : Fvol=cbrt(Vol0/getBox().determinant());
790 : }
791 : //Set switching function parameters
792 6 : log << "\n";
793 6 : log << "REFERENCE PDB # " << prev_stp+2 << " \n";
794 : // Set switching function parameters here only if computing derivatives
795 : // now set at the beginning of the dynamics to solve the r0 issue
796 6 : log << "Switching Function Parameters \n";
797 6 : sfs.resize(Nlist);
798 : std::string errors;
799 18 : for (unsigned j=0; j<Nlist; j++) {
800 12 : if(Svol) {
801 : double r0;
802 : std::string old_r0;
803 12 : vector<string> data=Tools::getWords(sw[j]);
804 : data.erase(data.begin());
805 12 : Tools::parse(data,"R_0",old_r0);
806 12 : Tools::convert(old_r0,r0);
807 12 : r0*=Fvol;
808 12 : std::string new_r0; Tools::convert(r0,new_r0);
809 12 : std::size_t pos = sw[j].find("R_0");
810 12 : sw[j].replace(pos+4,old_r0.size(),new_r0);
811 12 : }
812 12 : sfs[j].set(sw[j],errors);
813 : std::string num;
814 12 : Tools::convert(j+1, num);
815 12 : if( errors.length()!=0 ) error("problem reading SWITCH" + num + " keyword : " + errors );
816 12 : r00[j]=sfs[j].get_r0();
817 12 : log << " Swf: " << j << " r0=" << (sfs[j].description()).c_str() << " \n";
818 : }
819 : //Transform and sort
820 6 : log << "Building Reference PIV Vector \n";
821 6 : log << " PIV | block | Size | Zeros | Ones |" << " \n";
822 6 : double df=0.;
823 18 : for (unsigned j=0; j<Nlist; j++) {
824 11515212 : for (unsigned i=0; i<rPIV[j].size(); i++) {
825 11515200 : rPIV[j][i]=sfs[j].calculate(rPIV[j][i], df);
826 : }
827 12 : if(dosort[j]) {
828 12 : std::sort(rPIV[j].begin(),rPIV[j].end());
829 : }
830 : int lmt0=0;
831 : int lmt1=0;
832 11515212 : for(unsigned i=0; i<rPIV[j].size(); i++) {
833 11515200 : if(int(rPIV[j][i]*double(Nprec-1))==0) {
834 26 : lmt0+=1;
835 : }
836 11515200 : if(int(rPIV[j][i]*double(Nprec-1))==1) {
837 63358 : lmt1+=1;
838 : }
839 : }
840 12 : log.printf(" |%10i|%15zu|%15i|%15i|\n", j, rPIV[j].size(), lmt0, lmt1);
841 : }
842 6 : log << "\n";
843 : }
844 : // Do the sorting only once per timestep to avoid building the PIV N times for N rPIV PDB structures!
845 327 : if ((getStep()>prev_stp&&getStep()%updatePIV==0)||CompDer) {
846 324 : if (CompDer) log << " Step " << getStep() << " Computing Derivatives NON-SORTED PIV \n";
847 : //
848 : // build COMs from positions if requested
849 324 : if(com) {
850 0 : if(pbc) makeWhole();
851 0 : for(unsigned j=0; j<compos.size(); j++) {
852 0 : compos[j][0]=0.;
853 0 : compos[j][1]=0.;
854 0 : compos[j][2]=0.;
855 0 : for(unsigned i=0; i<nlcom[j]->getFullAtomList().size(); i++) {
856 0 : unsigned andx=nlcom[j]->getFullAtomList()[i].index();
857 0 : compos[j]+=fmass[andx]*getPosition(andx);
858 : }
859 : }
860 : }
861 : // update neighbor lists when an atom moves out of the Neighbor list skin
862 324 : if (doneigh) {
863 : bool doupdate=false;
864 : // For the first step build previous positions = actual positions
865 324 : if (prev_stp==-1) {
866 6 : bool docom=com;
867 18 : for (unsigned j=0; j<Nlist; j++) {
868 9708 : for (unsigned i=0; i<nl[j]->getFullAtomList().size(); i++) {
869 9696 : Vector Pos;
870 9696 : if(docom) {
871 0 : Pos=compos[i];
872 : } else {
873 9696 : Pos=getPosition(nl[j]->getFullAtomList()[i].index());
874 : }
875 9696 : prev_pos[j].push_back(Pos);
876 : }
877 : }
878 : doupdate=true;
879 : }
880 : // Decide whether to update lists based on atom displacement, every stride
881 324 : std:: vector<std:: vector<Vector> > tmp_pos(Nlist);
882 324 : if (getStep() % nlall->getStride() ==0) {
883 324 : bool docom=com;
884 972 : for (unsigned j=0; j<Nlist; j++) {
885 20520 : for (unsigned i=0; i<nl[j]->getFullAtomList().size(); i++) {
886 19872 : Vector Pos;
887 19872 : if(docom) {
888 0 : Pos=compos[i];
889 : } else {
890 19872 : Pos=getPosition(nl[j]->getFullAtomList()[i].index());
891 : }
892 19872 : tmp_pos[j].push_back(Pos);
893 19872 : if (pbcDistance(tmp_pos[j][i],prev_pos[j][i]).modulo()>=nl_skin[j]) {
894 : doupdate=true;
895 : }
896 : }
897 : }
898 : }
899 : // Update Nlists if needed
900 324 : if (doupdate==true) {
901 18 : for (unsigned j=0; j<Nlist; j++) {
902 9708 : for (unsigned i=0; i<nl[j]->getFullAtomList().size(); i++) {
903 9696 : prev_pos[j][i]=tmp_pos[j][i];
904 : }
905 12 : nl[j]->update(prev_pos[j]);
906 12 : log << " Step " << getStep() << " Neighbour lists updated " << nl[j]->size() << " \n";
907 : }
908 : }
909 324 : }
910 : // Calculate the volume scaling factor
911 324 : if(Svol) {
912 324 : Fvol=cbrt(Vol0/getBox().determinant());
913 : }
914 324 : Vector ddist;
915 : // Global to local variables
916 324 : bool doserial=serial;
917 : // Build "Nlist" PIV blocks
918 972 : for(unsigned j=0; j<Nlist; j++) {
919 648 : if(dosort[j]) {
920 : // from global to local variables to speedup the for loop with if statements
921 6 : bool docom=com;
922 6 : bool dopbc=pbc;
923 : // Vectors collecting occupancies: OrdVec one rank, OrdVecAll all ranks
924 6 : std:: vector<int> OrdVec(Nprec,0);
925 6 : cPIV[j].resize(0);
926 6 : Atom0[j].resize(0);
927 6 : Atom1[j].resize(0);
928 : // Building distances for the PIV vector at time t
929 6 : if(timer) stopwatch.start("1 Build cPIV");
930 5757606 : for(unsigned i=rank; i<nl[j]->size(); i+=stride) {
931 5757600 : unsigned i0=(nl[j]->getClosePairAtomNumber(i).first).index();
932 5757600 : unsigned i1=(nl[j]->getClosePairAtomNumber(i).second).index();
933 5757600 : Vector Pos0,Pos1;
934 5757600 : if(docom) {
935 0 : Pos0=compos[i0];
936 0 : Pos1=compos[i1];
937 : } else {
938 5757600 : Pos0=getPosition(i0);
939 5757600 : Pos1=getPosition(i1);
940 : }
941 5757600 : if(dopbc) {
942 5757600 : ddist=pbcDistance(Pos0,Pos1);
943 : } else {
944 0 : ddist=delta(Pos0,Pos1);
945 : }
946 5757600 : double df=0.;
947 : //Integer sorting ... faster!
948 : //Transforming distances with the Switching function + real to integer transformation
949 5757600 : int Vint=int(sfs[j].calculate(ddist.modulo()*Fvol, df)*double(Nprec-1)+0.5);
950 : //Integer transformed distance values as index of the Ordering Vector OrdVec
951 5757600 : OrdVec[Vint]+=1;
952 : //Keeps track of atom indices for force and virial calculations
953 5757600 : A0[Vint].push_back(i0);
954 5757600 : A1[Vint].push_back(i1);
955 : }
956 6 : if(timer) stopwatch.stop("1 Build cPIV");
957 6 : if(timer) stopwatch.start("2 Sort cPIV");
958 6 : if(!doserial && comm.initialized()) {
959 : // Vectors keeping track of the dimension and the starting-position of the rank-specific pair vector in the big pair vector.
960 0 : std:: vector<int> Vdim(stride,0);
961 0 : std:: vector<int> Vpos(stride,0);
962 : // Vectors collecting occupancies: OrdVec one rank, OrdVecAll all ranks
963 0 : std:: vector<int> OrdVecAll(stride*Nprec);
964 : // Big vectors containing all Atom indexes for every occupancy (Atom0O(Nprec,n) and Atom1O(Nprec,n) matrices in one vector)
965 : std:: vector<int> Atom0F;
966 : std:: vector<int> Atom1F;
967 : // Vector used to reconstruct arrays
968 0 : std:: vector<unsigned> k(stride,0);
969 : // Zeros might be many, this slows down a lot due to MPI communication
970 : // Avoid passing the zeros (i=1) for atom indices
971 0 : for(unsigned i=1; i<Nprec; i++) {
972 : // Building long vectors with all atom indexes for occupancies ordered from i=1 to i=Nprec-1
973 : // Can this be avoided ???
974 0 : Atom0F.insert(Atom0F.end(),A0[i].begin(),A0[i].end());
975 0 : Atom1F.insert(Atom1F.end(),A1[i].begin(),A1[i].end());
976 0 : A0[i].resize(0);
977 0 : A1[i].resize(0);
978 : }
979 : // Resize partial arrays to fill up for the next PIV block
980 0 : A0[0].resize(0);
981 0 : A1[0].resize(0);
982 0 : A0[Nprec-1].resize(0);
983 0 : A1[Nprec-1].resize(0);
984 : // Avoid passing the zeros (i=1) for atom indices
985 0 : OrdVec[0]=0;
986 0 : OrdVec[Nprec-1]=0;
987 :
988 : // Wait for all ranks before communication of Vectors
989 0 : comm.Barrier();
990 :
991 : // pass the array sizes before passing the arrays
992 0 : int dim=Atom0F.size();
993 : // Vdim and Vpos keep track of the dimension and the starting-position of the rank-specific pair vector in the big pair vector.
994 0 : comm.Allgather(&dim,1,&Vdim[0],1);
995 :
996 : // TO BE IMPROVED: the following may be done by the rank 0 (now every rank does it)
997 : int Fdim=0;
998 0 : for(unsigned i=1; i<stride; i++) {
999 0 : Vpos[i]=Vpos[i-1]+Vdim[i-1];
1000 0 : Fdim+=Vdim[i];
1001 : }
1002 0 : Fdim+=Vdim[0];
1003 : // build big vectors for atom pairs on all ranks for all ranks
1004 0 : std:: vector<int> Atom0FAll(Fdim);
1005 0 : std:: vector<int> Atom1FAll(Fdim);
1006 : // TO BE IMPROVED: Allgathers may be substituted by gathers by proc 0
1007 : // Moreover vectors are gathered head-to-tail and assembled later-on in a serial step.
1008 : // Gather the full Ordering Vector (occupancies). This is what we need to build the PIV
1009 0 : comm.Allgather(&OrdVec[0],Nprec,&OrdVecAll[0],Nprec);
1010 : // Gather the vectors of atom pairs to keep track of the idexes for the forces
1011 0 : comm.Allgatherv(Atom0F.data(),Atom0F.size(),&Atom0FAll[0],&Vdim[0],&Vpos[0]);
1012 0 : comm.Allgatherv(Atom1F.data(),Atom1F.size(),&Atom1FAll[0],&Vdim[0],&Vpos[0]);
1013 :
1014 : // Reconstruct the full vectors from collections of Allgathered parts (this is a serial step)
1015 : // This is the tricky serial step, to assemble together PIV and atom-pair info from head-tail big vectors
1016 : // Loop before on l and then on i would be better but the allgather should be modified
1017 : // Loop on blocks
1018 : //for(unsigned m=0;m<Nlist;m++) {
1019 : // Loop on Ordering Vector size excluding zeros (i=1)
1020 0 : if(timer) stopwatch.stop("2 Sort cPIV");
1021 0 : if(timer) stopwatch.start("3 Reconstruct cPIV");
1022 0 : for(unsigned i=1; i<Nprec; i++) {
1023 : // Loop on the ranks
1024 0 : for(unsigned l=0; l<stride; l++) {
1025 : // Loop on the number of head-to-tail pieces
1026 0 : for(unsigned m=0; m<OrdVecAll[i+l*Nprec]; m++) {
1027 : // cPIV is the current PIV at time t
1028 0 : cPIV[j].push_back(double(i)/double(Nprec-1));
1029 0 : Atom0[j].push_back(Atom0FAll[k[l]+Vpos[l]]);
1030 0 : Atom1[j].push_back(Atom1FAll[k[l]+Vpos[l]]);
1031 0 : k[l]+=1;
1032 : }
1033 : }
1034 : }
1035 0 : if(timer) stopwatch.stop("3 Reconstruct cPIV");
1036 : } else {
1037 6000000 : for(unsigned i=1; i<Nprec; i++) {
1038 11757594 : for(unsigned m=0; m<OrdVec[i]; m++) {
1039 5757600 : cPIV[j].push_back(double(i)/double(Nprec-1));
1040 5757600 : Atom0[j].push_back(A0[i][m]);
1041 5757600 : Atom1[j].push_back(A1[i][m]);
1042 : }
1043 : }
1044 : }
1045 : }
1046 : }
1047 : }
1048 327 : Vector distance;
1049 327 : double dfunc=0.;
1050 : // Calculate volume scaling factor
1051 327 : if(Svol) {
1052 327 : Fvol=cbrt(Vol0/getBox().determinant());
1053 : }
1054 :
1055 : // This test may be run by specifying the TEST keyword as input, it pritnts rPIV and cPIV and quits
1056 327 : if(test) {
1057 : unsigned limit=0;
1058 0 : for(unsigned j=0; j<Nlist; j++) {
1059 0 : if(dosort[j]) {
1060 0 : limit = cPIV[j].size();
1061 : } else {
1062 0 : limit = rPIV[j].size();
1063 : }
1064 0 : log.printf("PIV Block: %6i %12s %6i \n", j, " Size:", limit);
1065 0 : log.printf("%6s%6s%12s%12s%36s\n"," i"," j", " c-PIV "," r-PIV "," i-j distance vector ");
1066 0 : for(unsigned i=0; i<limit; i++) {
1067 : unsigned i0=0;
1068 : unsigned i1=0;
1069 0 : if(dosort[j]) {
1070 0 : i0=Atom0[j][i];
1071 0 : i1=Atom1[j][i];
1072 : } else {
1073 0 : i0=(nl[j]->getClosePairAtomNumber(i).first).index();
1074 0 : i1=(nl[j]->getClosePairAtomNumber(i).second).index();
1075 : }
1076 0 : Vector Pos0,Pos1;
1077 0 : if(com) {
1078 0 : Pos0=compos[i0];
1079 0 : Pos1=compos[i1];
1080 : } else {
1081 0 : Pos0=getPosition(i0);
1082 0 : Pos1=getPosition(i1);
1083 : }
1084 0 : if(pbc) {
1085 0 : distance=pbcDistance(Pos0,Pos1);
1086 : } else {
1087 0 : distance=delta(Pos0,Pos1);
1088 : }
1089 0 : dfunc=0.;
1090 : double cP,rP;
1091 0 : if(dosort[j]) {
1092 0 : cP = cPIV[j][i];
1093 0 : rP = rPIV[j][rPIV[j].size()-cPIV[j].size()+i];
1094 : } else {
1095 0 : double dm=distance.modulo();
1096 0 : cP = sfs[j].calculate(dm*Fvol, dfunc);
1097 0 : rP = rPIV[j][i];
1098 : }
1099 0 : log.printf("%6i%6i%12.6f%12.6f%12.6f%12.6f%12.6f\n",i0,i1,cP,rP,distance[0],distance[1],distance[2]);
1100 : }
1101 : }
1102 0 : log.printf("This was a test, now exit \n");
1103 0 : exit();
1104 : }
1105 :
1106 327 : if(timer) stopwatch.start("4 Build For Derivatives");
1107 : // non-global variables Nder and Scalevol defined to speedup if structures in cycles
1108 327 : bool Nder=CompDer;
1109 327 : bool Scalevol=Svol;
1110 327 : if(getStep()%updatePIV==0) {
1111 : // set to zero PIVdistance, derivatives and virial when they are calculated
1112 29799 : for(unsigned j=0; j<m_deriv.size(); j++) {
1113 117888 : for(unsigned k=0; k<3; k++) {m_deriv[j][k]=0.;}
1114 : }
1115 1308 : for(unsigned j=0; j<3; j++) {
1116 3924 : for(unsigned k=0; k<3; k++) {
1117 2943 : m_virial[j][k]=0.;
1118 : }
1119 : }
1120 327 : m_PIVdistance=0.;
1121 : // Re-compute atomic distances for derivatives and compute PIV-PIV distance
1122 981 : for(unsigned j=0; j<Nlist; j++) {
1123 : unsigned limit=0;
1124 : // dosorting definition is to speedup if structure in cycles with non-global variables
1125 654 : bool dosorting=dosort[j];
1126 654 : bool docom=com;
1127 654 : bool dopbc=pbc;
1128 654 : if(dosorting) {
1129 12 : limit = cPIV[j].size();
1130 : } else {
1131 642 : limit = rPIV[j].size();
1132 : }
1133 11574918 : for(unsigned i=rank; i<limit; i+=stride) {
1134 : unsigned i0=0;
1135 : unsigned i1=0;
1136 11574264 : if(dosorting) {
1137 11515200 : i0=Atom0[j][i];
1138 11515200 : i1=Atom1[j][i];
1139 : } else {
1140 59064 : i0=(nl[j]->getClosePairAtomNumber(i).first).index();
1141 59064 : i1=(nl[j]->getClosePairAtomNumber(i).second).index();
1142 : }
1143 11574264 : Vector Pos0,Pos1;
1144 11574264 : if(docom) {
1145 0 : Pos0=compos[i0];
1146 0 : Pos1=compos[i1];
1147 : } else {
1148 11574264 : Pos0=getPosition(i0);
1149 11574264 : Pos1=getPosition(i1);
1150 : }
1151 11574264 : if(dopbc) {
1152 11574264 : distance=pbcDistance(Pos0,Pos1);
1153 : } else {
1154 0 : distance=delta(Pos0,Pos1);
1155 : }
1156 11574264 : dfunc=0.;
1157 : // this is needed for dfunc and dervatives
1158 11574264 : double dm=distance.modulo();
1159 11574264 : double tPIV = sfs[j].calculate(dm*Fvol, dfunc);
1160 : // PIV distance
1161 : double coord=0.;
1162 11574264 : if(!dosorting||Nder) {
1163 59064 : coord = tPIV - rPIV[j][i];
1164 : } else {
1165 11515200 : coord = cPIV[j][i] - rPIV[j][rPIV[j].size()-cPIV[j].size()+i];
1166 : }
1167 : // Calculate derivatives, virial, and variable=sum_j (scaling[j] *(cPIV-rPIV)_j^2)
1168 : // WARNING: dfunc=dswf/(Fvol*dm) (this may change in future Plumed versions)
1169 11574264 : double tmp = 2.*scaling[j]*coord*Fvol*Fvol*dfunc;
1170 11574264 : Vector tmpder = tmp*distance;
1171 : // 0.5*(x_i-x_k)*f_ik (force on atom k due to atom i)
1172 11574264 : if(docom) {
1173 0 : Vector dist;
1174 0 : for(unsigned k=0; k<nlcom[i0]->getFullAtomList().size(); k++) {
1175 0 : unsigned x0=nlcom[i0]->getFullAtomList()[k].index();
1176 0 : m_deriv[x0] -= tmpder*fmass[x0];
1177 0 : for(unsigned l=0; l<3; l++) {
1178 0 : dist[l]=0.;
1179 : }
1180 0 : Vector P0=getPosition(x0);
1181 0 : for(unsigned l=0; l<nlcom[i0]->getFullAtomList().size(); l++) {
1182 0 : unsigned x1=nlcom[i0]->getFullAtomList()[l].index();
1183 0 : Vector P1=getPosition(x1);
1184 0 : if(dopbc) {
1185 0 : dist+=pbcDistance(P0,P1);
1186 : } else {
1187 0 : dist+=delta(P0,P1);
1188 : }
1189 : }
1190 0 : for(unsigned l=0; l<nlcom[i1]->getFullAtomList().size(); l++) {
1191 0 : unsigned x1=nlcom[i1]->getFullAtomList()[l].index();
1192 0 : Vector P1=getPosition(x1);
1193 0 : if(dopbc) {
1194 0 : dist+=pbcDistance(P0,P1);
1195 : } else {
1196 0 : dist+=delta(P0,P1);
1197 : }
1198 : }
1199 0 : m_virial -= 0.25*fmass[x0]*Tensor(dist,tmpder);
1200 : }
1201 0 : for(unsigned k=0; k<nlcom[i1]->getFullAtomList().size(); k++) {
1202 0 : unsigned x1=nlcom[i1]->getFullAtomList()[k].index();
1203 0 : m_deriv[x1] += tmpder*fmass[x1];
1204 0 : for(unsigned l=0; l<3; l++) {
1205 0 : dist[l]=0.;
1206 : }
1207 0 : Vector P1=getPosition(x1);
1208 0 : for(unsigned l=0; l<nlcom[i1]->getFullAtomList().size(); l++) {
1209 0 : unsigned x0=nlcom[i1]->getFullAtomList()[l].index();
1210 0 : Vector P0=getPosition(x0);
1211 0 : if(dopbc) {
1212 0 : dist+=pbcDistance(P1,P0);
1213 : } else {
1214 0 : dist+=delta(P1,P0);
1215 : }
1216 : }
1217 0 : for(unsigned l=0; l<nlcom[i0]->getFullAtomList().size(); l++) {
1218 0 : unsigned x0=nlcom[i0]->getFullAtomList()[l].index();
1219 0 : Vector P0=getPosition(x0);
1220 0 : if(dopbc) {
1221 0 : dist+=pbcDistance(P1,P0);
1222 : } else {
1223 0 : dist+=delta(P1,P0);
1224 : }
1225 : }
1226 0 : m_virial += 0.25*fmass[x1]*Tensor(dist,tmpder);
1227 : }
1228 : } else {
1229 11574264 : m_deriv[i0] -= tmpder;
1230 11574264 : m_deriv[i1] += tmpder;
1231 11574264 : m_virial -= tmp*Tensor(distance,distance);
1232 : }
1233 11574264 : if(Scalevol) {
1234 11574264 : m_virial+=1./3.*tmp*dm*dm*Tensor::identity();
1235 : }
1236 11574264 : m_PIVdistance += scaling[j]*coord*coord;
1237 : }
1238 : }
1239 :
1240 327 : if (!serial && comm.initialized()) {
1241 0 : comm.Barrier();
1242 0 : comm.Sum(&m_PIVdistance,1);
1243 0 : if(!m_deriv.empty()) comm.Sum(&m_deriv[0][0],3*m_deriv.size());
1244 0 : comm.Sum(&m_virial[0][0],9);
1245 : }
1246 : }
1247 327 : prev_stp=getStep();
1248 :
1249 : //Timing
1250 327 : if(timer) stopwatch.stop("4 Build For Derivatives");
1251 327 : if(timer) {
1252 1 : log.printf("Timings for action %s with label %s \n", getName().c_str(), getLabel().c_str() );
1253 1 : log<<stopwatch;
1254 : }
1255 :
1256 : // Update derivatives, virial, and variable (PIV-distance^2)
1257 29799 : for(unsigned i=0; i<m_deriv.size(); ++i) setAtomsDerivatives(i,m_deriv[i]);
1258 327 : setValue (m_PIVdistance);
1259 327 : setBoxDerivatives (m_virial);
1260 327 : }
1261 : //Close Namespaces at the very beginning
1262 : }
1263 : }
1264 :
|