LCOV - code coverage report
Current view: top level - piv - PIV.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 438 618 70.9 %
Date: 2025-03-25 09:33:27 Functions: 3 5 60.0 %

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

Generated by: LCOV version 1.16