LCOV - code coverage report
Current view: top level - cltools - Driver.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 543 585 92.8 %
Date: 2018-12-19 07:49:13 Functions: 10 16 62.5 %

          Line data    Source code
       1             : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
       2             :    Copyright (c) 2012-2018 The plumed team
       3             :    (see the PEOPLE file at the root of the distribution for a list of names)
       4             : 
       5             :    See http://www.plumed.org for more information.
       6             : 
       7             :    This file is part of plumed, version 2.
       8             : 
       9             :    plumed is free software: you can redistribute it and/or modify
      10             :    it under the terms of the GNU Lesser General Public License as published by
      11             :    the Free Software Foundation, either version 3 of the License, or
      12             :    (at your option) any later version.
      13             : 
      14             :    plumed is distributed in the hope that it will be useful,
      15             :    but WITHOUT ANY WARRANTY; without even the implied warranty of
      16             :    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
      17             :    GNU Lesser General Public License for more details.
      18             : 
      19             :    You should have received a copy of the GNU Lesser General Public License
      20             :    along with plumed.  If not, see <http://www.gnu.org/licenses/>.
      21             : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
      22             : #include "CLTool.h"
      23             : #include "CLToolRegister.h"
      24             : #include "tools/Tools.h"
      25             : #include "wrapper/Plumed.h"
      26             : #include "tools/Communicator.h"
      27             : #include "tools/Random.h"
      28             : #include "tools/Pbc.h"
      29             : #include <cstdio>
      30             : #include <cstring>
      31             : #include <vector>
      32             : #include <map>
      33             : #include "tools/Units.h"
      34             : #include "tools/PDB.h"
      35             : #include "tools/FileBase.h"
      36             : #include "tools/IFile.h"
      37             : 
      38             : // when using molfile plugin
      39             : #ifdef __PLUMED_HAS_MOLFILE_PLUGINS
      40             : #ifndef __PLUMED_HAS_EXTERNAL_MOLFILE_PLUGINS
      41             : /* Use the internal ones. Alternatively:
      42             :  *    ifeq (,$(findstring __PLUMED_HAS_EXTERNAL_MOLFILE_PLUGINS,$(CPPFLAGS)))
      43             :  *    CPPFLAGS+=-I../molfile
      44             :  */
      45             : #include "molfile/libmolfile_plugin.h"
      46             : #include "molfile/molfile_plugin.h"
      47             : using namespace PLMD::molfile;
      48             : #else
      49             : #include <libmolfile_plugin.h>
      50             : #include <molfile_plugin.h>
      51             : #endif
      52             : #endif
      53             : 
      54             : #ifdef __PLUMED_HAS_XDRFILE
      55             : #include <xdrfile/xdrfile_trr.h>
      56             : #include <xdrfile/xdrfile_xtc.h>
      57             : #endif
      58             : 
      59             : using namespace std;
      60             : 
      61             : namespace PLMD {
      62             : namespace cltools {
      63             : 
      64             : //+PLUMEDOC TOOLS driver-float
      65             : /*
      66             : Equivalent to \ref driver, but using single precision reals.
      67             : 
      68             : The purpose of this tool is just to test what PLUMED does when linked from
      69             : a single precision code.
      70             : 
      71             : \par Examples
      72             : 
      73             : \verbatim
      74             : plumed driver-float --plumed plumed.dat --ixyz trajectory.xyz
      75             : \endverbatim
      76             : 
      77             : See also examples in \ref driver
      78             : 
      79             : */
      80             : //+ENDPLUMEDOC
      81             : //
      82             : 
      83             : 
      84             : //+PLUMEDOC TOOLS driver
      85             : /*
      86             : driver is a tool that allows one to to use plumed to post-process an existing trajectory.
      87             : 
      88             : The input to driver is specified using the command line arguments described below.
      89             : 
      90             : In addition, you can use the special \subpage READ command inside your plumed input
      91             : to read in colvar files that were generated during your MD simulation.  The values
      92             : read in can then be treated like calculated colvars.
      93             : 
      94             : \warning
      95             : Notice that by default the driver has no knowledge about the masses and charges
      96             : of your atoms! Thus, if you want to compute quantities depending charges (e.g. \ref DHENERGY)
      97             : or masses (e.g. \ref COM) you should pass the proper information to the driver.
      98             : You can do it either with the --pdb option or with the --mc option. The latter
      99             : will read a file produced by \ref DUMPMASSCHARGE .
     100             : 
     101             : 
     102             : \par Examples
     103             : 
     104             : The following command tells plumed to postprocess the trajectory contained in `trajectory.xyz`
     105             :  by performing the actions described in the input file `plumed.dat`.  If an action that takes the
     106             : stride keyword is given a stride equal to \f$n\f$ then it will be performed only on every \f$n\f$th
     107             : frames in the trajectory file.
     108             : \verbatim
     109             : plumed driver --plumed plumed.dat --ixyz trajectory.xyz
     110             : \endverbatim
     111             : 
     112             : Notice that `xyz` files are expected to be in internal PLUMED units, that is by default nm.
     113             : You can change this behavior by using the `--length-units` option:
     114             : \verbatim
     115             : plumed driver --plumed plumed.dat --ixyz trajectory.xyz --length-units A
     116             : \endverbatim
     117             : The strings accepted by the `--length-units` options are the same ones accepted by the \ref UNITS action.
     118             : Other file formats typically have their default coordinates (e.g., `gro` files are always in nm)
     119             : and it thus should not be necessary to use the `--length-units` option. Additionally,
     120             : consider that the units used by the `driver` might be different by the units used in the PLUMED input
     121             : file `plumed.dat`. For instance consider the command:
     122             : \verbatim
     123             : plumed driver --plumed plumed.dat --ixyz trajectory.xyz --length-units A
     124             : \endverbatim
     125             : where `plumed.dat` is
     126             : \verbatim
     127             : # no explicit UNITS action here
     128             : d: DISTANCE ATOMS=1,2
     129             : PRINT ARG=d FILE=colvar
     130             : \endverbatim
     131             : In this case, the driver reads the `xyz` file assuming it to contain coordinates in Angstrom units.
     132             : However, the resulting `colvar` file contains a distance expressed in nm.
     133             : 
     134             : The following command tells plumed to postprocess the trajectory contained in trajectory.xyz.
     135             :  by performing the actions described in the input file plumed.dat.
     136             : \verbatim
     137             : plumed driver --plumed plumed.dat --ixyz trajectory.xyz --trajectory-stride 100 --timestep 0.001
     138             : \endverbatim
     139             : Here though
     140             : `--trajectory-stride` is set equal to the frequency with which frames were output during the trajectory
     141             : and the `--timestep` is equal to the simulation timestep.  As such the `STRIDE` parameters in the `plumed.dat`
     142             : files are referred to the original timestep and any files output resemble those that would have been generated
     143             : had we run the calculation we are running with driver when the MD simulation was running.
     144             : 
     145             : PLUMED can read natively xyz files (in PLUMED units) and gro files (in nm). In addition,
     146             : PLUMED includes by default support for a
     147             : subset of the trajectory file formats supported by VMD, e.g. xtc and dcd:
     148             : 
     149             : \verbatim
     150             : plumed driver --plumed plumed.dat --pdb diala.pdb --mf_xtc traj.xtc --trajectory-stride 100 --timestep 0.001
     151             : \endverbatim
     152             : 
     153             : where `--mf_` prefixes the extension of one of the accepted molfile plugin format.
     154             : If PLUMED has been \ref Installation "installed" with full molfile support, other formats will be available.
     155             : Just type `plumed driver --help` to see which plugins are available.
     156             : 
     157             : Molfile plugin require periodic cell to be triangular (i.e. first vector oriented along x and
     158             : second vector in xy plane). This is true for many MD codes. However, it could be false
     159             : if you rotate the coordinates in your trajectory before reading them in the driver.
     160             : Also notice that some formats (e.g. amber crd) do not specify atom number. In this case you can use
     161             : the `--natoms` option:
     162             : \verbatim
     163             : plumed driver --plumed plumed.dat --imf_crd trajectory.crd --natoms 128
     164             : \endverbatim
     165             : 
     166             : Check the available molfile plugins and limitations at [this link](http://www.ks.uiuc.edu/Research/vmd/plugins/molfile/).
     167             : 
     168             : Additionally, you can use the xdrfile implementation of xtc and trr. To this aim, just
     169             : download and install properly the xdrfile library (see [this link](http://www.gromacs.org/Developer_Zone/Programming_Guide/XTC_Library)).
     170             : If the xdrfile library is installed properly the PLUMED configure script should be able to
     171             : detect it and enable it.
     172             : Notice that the xdrfile implementation of xtc and trr
     173             : is more robust than the molfile one, since it provides support for generic cell shapes.
     174             : In addition, it allows \ref DUMPATOMS to write compressed xtc files.
     175             : 
     176             : 
     177             : */
     178             : //+ENDPLUMEDOC
     179             : //
     180             : 
     181             : #ifdef __PLUMED_HAS_MOLFILE_PLUGINS
     182        2523 : static vector<molfile_plugin_t *> plugins;
     183        2523 : static map <string, unsigned> pluginmap;
     184       15138 : static int register_cb(void *v, vmdplugin_t *p) {
     185             :   //const char *key = p->name;
     186       15138 :   std::pair<std::map<string,unsigned>::iterator,bool> ret;
     187       15138 :   ret = pluginmap.insert ( std::pair<string,unsigned>(string(p->name),plugins.size()) );
     188       15138 :   if (ret.second==false) {
     189             :     //cerr<<"MOLFILE: found duplicate plugin for "<<key<<" : not inserted "<<endl;
     190             :   } else {
     191             :     //cerr<<"MOLFILE: loading plugin "<<key<<" number "<<plugins.size()-1<<endl;
     192       15138 :     plugins.push_back(reinterpret_cast<molfile_plugin_t *>(p));
     193             :   }
     194       15138 :   return VMDPLUGIN_SUCCESS;
     195             : }
     196             : #endif
     197             : 
     198             : template<typename real>
     199         576 : class Driver : public CLTool {
     200             : public:
     201             :   static void registerKeywords( Keywords& keys );
     202             :   explicit Driver(const CLToolOptions& co );
     203             :   int main(FILE* in,FILE*out,Communicator& pc);
     204             :   void evaluateNumericalDerivatives( const long int& step, Plumed& p, const std::vector<real>& coordinates,
     205             :                                      const std::vector<real>& masses, const std::vector<real>& charges,
     206             :                                      std::vector<real>& cell, const double& base, std::vector<real>& numder );
     207             :   string description()const;
     208             : };
     209             : 
     210             : template<typename real>
     211        1682 : void Driver<real>::registerKeywords( Keywords& keys ) {
     212        1682 :   CLTool::registerKeywords( keys ); keys.isDriver();
     213        1682 :   keys.addFlag("--help-debug",false,"print special options that can be used to create regtests");
     214        1682 :   keys.add("compulsory","--plumed","plumed.dat","specify the name of the plumed input file");
     215        1682 :   keys.add("compulsory","--timestep","1.0","the timestep that was used in the calculation that produced this trajectory in picoseconds");
     216        1682 :   keys.add("compulsory","--trajectory-stride","1","the frequency with which frames were output to this trajectory during the simulation"
     217             : #ifdef __PLUMED_HAS_XDRFILE
     218             :            " (0 means that the number of the step is read from the trajectory file,"
     219             :            " currently working only for xtc/trr files read with --ixtc/--trr)"
     220             : #endif
     221        3364 :           );
     222        1682 :   keys.add("compulsory","--multi","0","set number of replicas for multi environment (needs mpi)");
     223        1682 :   keys.addFlag("--noatoms",false,"don't read in a trajectory.  Just use colvar files as specified in plumed.dat");
     224        1682 :   keys.add("atoms","--ixyz","the trajectory in xyz format");
     225        1682 :   keys.add("atoms","--igro","the trajectory in gro format");
     226             : #ifdef __PLUMED_HAS_XDRFILE
     227        1682 :   keys.add("atoms","--ixtc","the trajectory in xtc format (xdrfile implementation)");
     228        1682 :   keys.add("atoms","--itrr","the trajectory in trr format (xdrfile implementation)");
     229             : #endif
     230        1682 :   keys.add("optional","--length-units","units for length, either as a string or a number");
     231        1682 :   keys.add("optional","--mass-units","units for mass in pdb and mc file, either as a string or a number");
     232        1682 :   keys.add("optional","--charge-units","units for charge in pdb and mc file, either as a string or a number");
     233        1682 :   keys.add("optional","--dump-forces","dump the forces on a file");
     234        1682 :   keys.add("optional","--dump-forces-fmt","( default=%%f ) the format to use to dump the forces");
     235        1682 :   keys.addFlag("--dump-full-virial",false,"with --dump-forces, it dumps the 9 components of the virial");
     236        1682 :   keys.add("optional","--pdb","provides a pdb with masses and charges");
     237        1682 :   keys.add("optional","--mc","provides a file with masses and charges as produced with DUMPMASSCHARGE");
     238        1682 :   keys.add("optional","--box","comma-separated box dimensions (3 for orthorombic, 9 for generic)");
     239        1682 :   keys.add("optional","--natoms","provides number of atoms - only used if file format does not contain number of atoms");
     240        1682 :   keys.add("optional","--debug-forces","output a file containing the forces due to the bias evaluated using numerical derivatives "
     241        3364 :            "and using the analytical derivatives implemented in plumed");
     242        1682 :   keys.add("hidden","--debug-float","[yes/no] turns on the single precision version (to check float interface)");
     243        1682 :   keys.add("hidden","--debug-dd","[yes/no] use a fake domain decomposition");
     244        1682 :   keys.add("hidden","--debug-pd","[yes/no] use a fake particle decomposition");
     245        1682 :   keys.add("hidden","--debug-grex","use a fake gromacs-like replica exchange, specify exchange stride");
     246        1682 :   keys.add("hidden","--debug-grex-log","log file for debug=grex");
     247             : #ifdef __PLUMED_HAS_MOLFILE_PLUGINS
     248        1682 :   MOLFILE_INIT_ALL
     249        1682 :   MOLFILE_REGISTER_ALL(NULL, register_cb)
     250       16820 :   for(int i=0; i<plugins.size(); i++) {
     251       15138 :     string kk="--mf_"+string(plugins[i]->name);
     252       30276 :     string mm=" molfile: the trajectory in "+string(plugins[i]->name)+" format " ;
     253             :     //cerr<<"REGISTERING "<<kk<<mm<<endl;
     254       15138 :     keys.add("atoms",kk,mm);
     255             :   }
     256             : #endif
     257        1682 : }
     258             : template<typename real>
     259         288 : Driver<real>::Driver(const CLToolOptions& co ):
     260         288 :   CLTool(co)
     261             : {
     262         288 :   inputdata=commandline;
     263         288 : }
     264             : template<typename real>
     265           0 : string Driver<real>::description()const { return "analyze trajectories with plumed"; }
     266             : 
     267             : template<typename real>
     268         288 : int Driver<real>::main(FILE* in,FILE*out,Communicator& pc) {
     269             : 
     270         288 :   Units units;
     271         576 :   PDB pdb;
     272             : 
     273             : // Parse everything
     274         288 :   bool printhelpdebug; parseFlag("--help-debug",printhelpdebug);
     275         288 :   if( printhelpdebug ) {
     276           0 :     fprintf(out,"%s",
     277             :             "Additional options for debug (only to be used in regtest):\n"
     278             :             "  [--debug-float yes]     : turns on the single precision version (to check float interface)\n"
     279             :             "  [--debug-dd yes]        : use a fake domain decomposition\n"
     280             :             "  [--debug-pd yes]        : use a fake particle decomposition\n"
     281             :            );
     282           0 :     return 0;
     283             :   }
     284             :   // Are we reading trajectory data
     285         288 :   bool noatoms; parseFlag("--noatoms",noatoms);
     286             : 
     287         576 :   std::string fakein;
     288             : 
     289         288 :   bool debug_float=false;
     290         288 :   fakein="";
     291         288 :   if(parse("--debug-float",fakein)) {
     292           0 :     if(fakein=="yes") debug_float=true;
     293           0 :     else if(fakein=="no") debug_float=false;
     294           0 :     else error("--debug-float should have argument yes or no");
     295             :   }
     296             : 
     297         288 :   if(debug_float && sizeof(real)!=sizeof(float)) {
     298           0 :     CLTool* cl=cltoolRegister().create(CLToolOptions("driver-float"));    //new Driver<float>(*this);
     299           0 :     cl->setInputData(this->getInputData());
     300           0 :     int ret=cl->main(in,out,pc);
     301           0 :     delete cl;
     302           0 :     return ret;
     303             :   }
     304             : 
     305         288 :   bool debug_pd=false;
     306         288 :   fakein="";
     307         288 :   if(parse("--debug-pd",fakein)) {
     308           4 :     if(fakein=="yes") debug_pd=true;
     309           0 :     else if(fakein=="no") debug_pd=false;
     310           0 :     else error("--debug-pd should have argument yes or no");
     311             :   }
     312         288 :   if(debug_pd) fprintf(out,"DEBUGGING PARTICLE DECOMPOSITION\n");
     313             : 
     314         288 :   bool debug_dd=false;
     315         288 :   fakein="";
     316         288 :   if(parse("--debug-dd",fakein)) {
     317          30 :     if(fakein=="yes") debug_dd=true;
     318           0 :     else if(fakein=="no") debug_dd=false;
     319           0 :     else error("--debug-dd should have argument yes or no");
     320             :   }
     321         288 :   if(debug_dd) fprintf(out,"DEBUGGING DOMAIN DECOMPOSITION\n");
     322             : 
     323         288 :   if( debug_pd || debug_dd ) {
     324          34 :     if(noatoms) error("cannot debug without atoms");
     325             :   }
     326             : 
     327             : // set up for multi replica driver:
     328         288 :   int multi=0;
     329         288 :   parse("--multi",multi);
     330         576 :   Communicator intracomm;
     331         576 :   Communicator intercomm;
     332         288 :   if(multi) {
     333          63 :     int ntot=pc.Get_size();
     334          63 :     int nintra=ntot/multi;
     335          63 :     if(multi*nintra!=ntot) error("invalid number of processes for multi environment");
     336          63 :     pc.Split(pc.Get_rank()/nintra,pc.Get_rank(),intracomm);
     337          63 :     pc.Split(pc.Get_rank()%nintra,pc.Get_rank(),intercomm);
     338             :   } else {
     339         225 :     intracomm.Set_comm(pc.Get_comm());
     340             :   }
     341             : 
     342             : // set up for debug replica exchange:
     343         288 :   bool debug_grex=parse("--debug-grex",fakein);
     344         288 :   int  grex_stride=0;
     345         288 :   FILE*grex_log=NULL;
     346         288 :   if(debug_grex) {
     347          18 :     if(noatoms) error("must have atoms to debug_grex");
     348          18 :     if(multi<2)  error("--debug_grex needs --multi with at least two replicas");
     349          18 :     Tools::convert(fakein,grex_stride);
     350          18 :     string n; Tools::convert(intercomm.Get_rank(),n);
     351          36 :     string file;
     352          18 :     parse("--debug-grex-log",file);
     353          18 :     if(file.length()>0) {
     354          18 :       file+="."+n;
     355          18 :       grex_log=fopen(file.c_str(),"w");
     356          18 :     }
     357             :   }
     358             : 
     359             : // Read the plumed input file name
     360         576 :   string plumedFile; parse("--plumed",plumedFile);
     361             : // the timestep
     362         288 :   double t; parse("--timestep",t);
     363         288 :   real timestep=real(t);
     364             : // the stride
     365         288 :   unsigned stride; parse("--trajectory-stride",stride);
     366             : // are we writing forces
     367         576 :   string dumpforces(""), debugforces(""), dumpforcesFmt("%f");;
     368         288 :   bool dumpfullvirial=false;
     369         288 :   if(!noatoms) {
     370         287 :     parse("--dump-forces",dumpforces);
     371         287 :     parse("--debug-forces",debugforces);
     372             :   }
     373         288 :   if(dumpforces!="" || debugforces!="" ) parse("--dump-forces-fmt",dumpforcesFmt);
     374         288 :   if(dumpforces!="") parseFlag("--dump-full-virial",dumpfullvirial);
     375         288 :   if( debugforces!="" && (debug_dd || debug_pd) ) error("cannot debug forces and domain/particle decomposition at same time");
     376         288 :   if( debugforces!="" && sizeof(real)!=sizeof(double) ) error("cannot debug forces in single precision mode");
     377             : 
     378         576 :   string trajectory_fmt;
     379             : 
     380         288 :   bool use_molfile=false;
     381             : #ifdef __PLUMED_HAS_MOLFILE_PLUGINS
     382         288 :   molfile_plugin_t *api=NULL;
     383         288 :   void *h_in=NULL;
     384             :   molfile_timestep_t ts_in; // this is the structure that has the timestep
     385         288 :   ts_in.coords=NULL;
     386         288 :   ts_in.A=-1; // we use this to check whether cell is provided or not
     387             : #endif
     388             : 
     389             : // Read in an xyz file
     390         576 :   string trajectoryFile(""), pdbfile(""), mcfile("");
     391         576 :   bool pbc_cli_given=false; vector<double> pbc_cli_box(9,0.0);
     392         288 :   int command_line_natoms=-1;
     393             : 
     394         288 :   if(!noatoms) {
     395         287 :     std::string traj_xyz; parse("--ixyz",traj_xyz);
     396         574 :     std::string traj_gro; parse("--igro",traj_gro);
     397         574 :     std::string traj_xtc;
     398         574 :     std::string traj_trr;
     399             : #ifdef __PLUMED_HAS_XDRFILE
     400         287 :     parse("--ixtc",traj_xtc);
     401         287 :     parse("--itrr",traj_trr);
     402             : #endif
     403             : #ifdef __PLUMED_HAS_MOLFILE_PLUGINS
     404        2870 :     for(int i=0; i<plugins.size(); i++) {
     405        2583 :       string molfile_key="--mf_"+string(plugins[i]->name);
     406        5166 :       string traj_molfile;
     407        2583 :       parse(molfile_key,traj_molfile);
     408        2583 :       if(traj_molfile.length()>0) {
     409          14 :         fprintf(out,"\nDRIVER: Found molfile format trajectory %s with name %s\n",plugins[i]->name,traj_molfile.c_str());
     410          14 :         trajectoryFile=traj_molfile;
     411          14 :         trajectory_fmt=string(plugins[i]->name);
     412          14 :         use_molfile=true;
     413          14 :         api = plugins[i];
     414             :       }
     415             :     }
     416             : #endif
     417             :     { // check that only one fmt is specified
     418         287 :       int nn=0;
     419         287 :       if(traj_xyz.length()>0) nn++;
     420         287 :       if(traj_gro.length()>0) nn++;
     421         287 :       if(traj_xtc.length()>0) nn++;
     422         287 :       if(traj_trr.length()>0) nn++;
     423         287 :       if(nn>1) {
     424           0 :         fprintf(stderr,"ERROR: cannot provide more than one trajectory file\n");
     425           0 :         if(grex_log)fclose(grex_log);
     426           0 :         return 1;
     427             :       }
     428             :     }
     429         287 :     if(traj_xyz.length()>0 && trajectoryFile.length()==0) {
     430         247 :       trajectoryFile=traj_xyz;
     431         247 :       trajectory_fmt="xyz";
     432             :     }
     433         287 :     if(traj_gro.length()>0 && trajectoryFile.length()==0) {
     434          22 :       trajectoryFile=traj_gro;
     435          22 :       trajectory_fmt="gro";
     436             :     }
     437         287 :     if(traj_xtc.length()>0 && trajectoryFile.length()==0) {
     438           1 :       trajectoryFile=traj_xtc;
     439           1 :       trajectory_fmt="xdr-xtc";
     440             :     }
     441         287 :     if(traj_trr.length()>0 && trajectoryFile.length()==0) {
     442           3 :       trajectoryFile=traj_trr;
     443           3 :       trajectory_fmt="xdr-trr";
     444             :     }
     445         287 :     if(trajectoryFile.length()==0) {
     446           0 :       fprintf(stderr,"ERROR: missing trajectory data\n");
     447           0 :       if(grex_log)fclose(grex_log);
     448           0 :       return 1;
     449             :     }
     450         574 :     string lengthUnits(""); parse("--length-units",lengthUnits);
     451         287 :     if(lengthUnits.length()>0) units.setLength(lengthUnits);
     452         574 :     string chargeUnits(""); parse("--charge-units",chargeUnits);
     453         287 :     if(chargeUnits.length()>0) units.setCharge(chargeUnits);
     454         574 :     string massUnits(""); parse("--mass-units",massUnits);
     455         287 :     if(massUnits.length()>0) units.setMass(massUnits);
     456             : 
     457         287 :     parse("--pdb",pdbfile);
     458         287 :     if(pdbfile.length()>0) {
     459          37 :       bool check=pdb.read(pdbfile,false,1.0);
     460          37 :       if(!check) error("error reading pdb file");
     461             :     }
     462             : 
     463         287 :     parse("--mc",mcfile);
     464             : 
     465         574 :     string pbc_cli_list; parse("--box",pbc_cli_list);
     466         287 :     if(pbc_cli_list.length()>0) {
     467          14 :       pbc_cli_given=true;
     468          14 :       vector<string> words=Tools::getWords(pbc_cli_list,",");
     469          14 :       if(words.size()==3) {
     470          14 :         for(int i=0; i<3; i++) sscanf(words[i].c_str(),"%100lf",&(pbc_cli_box[4*i]));
     471           0 :       } else if(words.size()==9) {
     472           0 :         for(int i=0; i<9; i++) sscanf(words[i].c_str(),"%100lf",&(pbc_cli_box[i]));
     473             :       } else {
     474           0 :         string msg="ERROR: cannot parse command-line box "+pbc_cli_list;
     475           0 :         fprintf(stderr,"%s\n",msg.c_str());
     476           0 :         return 1;
     477          14 :       }
     478             : 
     479             :     }
     480             : 
     481         574 :     parse("--natoms",command_line_natoms);
     482             : 
     483             :   }
     484             : 
     485         288 :   if( debug_dd && debug_pd ) error("cannot use debug-dd and debug-pd at the same time");
     486         288 :   if(debug_pd || debug_dd) {
     487          34 :     if( !Communicator::initialized() ) error("needs mpi for debug-pd");
     488             :   }
     489             : 
     490         576 :   Plumed p;
     491         288 :   int rr=sizeof(real);
     492         288 :   p.cmd("setRealPrecision",&rr);
     493         288 :   int checknatoms=-1;
     494         288 :   long int step=0;
     495         288 :   if(Communicator::initialized()) {
     496         115 :     if(multi) {
     497          63 :       if(intracomm.Get_rank()==0) p.cmd("GREX setMPIIntercomm",&intercomm.Get_comm());
     498          63 :       p.cmd("GREX setMPIIntracomm",&intracomm.Get_comm());
     499          63 :       p.cmd("GREX init");
     500             :     }
     501         115 :     p.cmd("setMPIComm",&intracomm.Get_comm());
     502             :   }
     503         288 :   p.cmd("setMDLengthUnits",&units.getLength());
     504         288 :   p.cmd("setMDChargeUnits",&units.getCharge());
     505         288 :   p.cmd("setMDMassUnits",&units.getMass());
     506         288 :   p.cmd("setMDEngine","driver");
     507         288 :   p.cmd("setTimestep",&timestep);
     508         288 :   p.cmd("setPlumedDat",plumedFile.c_str());
     509         288 :   p.cmd("setLog",out);
     510             : 
     511         288 :   if(multi) {
     512          63 :     string n;
     513          63 :     Tools::convert(intercomm.Get_rank(),n);
     514          63 :     trajectoryFile=FileBase::appendSuffix(trajectoryFile,"."+n);
     515             :   }
     516             : 
     517             : 
     518             :   int natoms;
     519             : 
     520         576 :   FILE* fp=NULL; FILE* fp_forces=NULL; OFile fp_dforces;
     521             : #ifdef __PLUMED_HAS_XDRFILE
     522         288 :   XDRFILE* xd=NULL;
     523             : #endif
     524         288 :   if(!noatoms) {
     525         287 :     if (trajectoryFile=="-")
     526           0 :       fp=in;
     527             :     else {
     528         287 :       if(use_molfile==true) {
     529             : #ifdef __PLUMED_HAS_MOLFILE_PLUGINS
     530          14 :         h_in = api->open_file_read(trajectoryFile.c_str(), trajectory_fmt.c_str(), &natoms);
     531          14 :         if(natoms==MOLFILE_NUMATOMS_UNKNOWN) {
     532           2 :           if(command_line_natoms>=0) natoms=command_line_natoms;
     533           0 :           else error("this file format does not provide number of atoms; use --natoms on the command line");
     534             :         }
     535          14 :         ts_in.coords = new float [3*natoms];
     536             : #endif
     537         273 :       } else if(trajectory_fmt=="xdr-xtc" || trajectory_fmt=="xdr-trr") {
     538             : #ifdef __PLUMED_HAS_XDRFILE
     539           4 :         xd=xdrfile_open(trajectoryFile.c_str(),"r");
     540           4 :         if(!xd) {
     541           0 :           string msg="ERROR: Error opening trajectory file "+trajectoryFile;
     542           0 :           fprintf(stderr,"%s\n",msg.c_str());
     543           0 :           return 1;
     544             :         }
     545           4 :         if(trajectory_fmt=="xdr-xtc") read_xtc_natoms(&trajectoryFile[0],&natoms);
     546           4 :         if(trajectory_fmt=="xdr-trr") read_trr_natoms(&trajectoryFile[0],&natoms);
     547             : #endif
     548             :       } else {
     549         269 :         fp=fopen(trajectoryFile.c_str(),"r");
     550         269 :         if(!fp) {
     551           0 :           string msg="ERROR: Error opening trajectory file "+trajectoryFile;
     552           0 :           fprintf(stderr,"%s\n",msg.c_str());
     553             : // cppcheck detects a false positive here. I suppress it:
     554             : // cppcheck-suppress memleak
     555           0 :           return 1;
     556             :         }
     557             :       }
     558             :     }
     559         287 :     if(dumpforces.length()>0) {
     560         183 :       if(Communicator::initialized() && pc.Get_size()>1) {
     561          87 :         string n;
     562          87 :         Tools::convert(pc.Get_rank(),n);
     563          87 :         dumpforces+="."+n;
     564             :       }
     565         183 :       fp_forces=fopen(dumpforces.c_str(),"w");
     566             :     }
     567         287 :     if(debugforces.length()>0) {
     568           3 :       if(Communicator::initialized() && pc.Get_size()>1) {
     569           0 :         string n;
     570           0 :         Tools::convert(pc.Get_rank(),n);
     571           0 :         debugforces+="."+n;
     572             :       }
     573           3 :       fp_dforces.open(debugforces);
     574             :     }
     575             :   }
     576             : 
     577         576 :   std::string line;
     578         576 :   std::vector<real> coordinates;
     579         576 :   std::vector<real> forces;
     580         576 :   std::vector<real> masses;
     581         576 :   std::vector<real> charges;
     582         576 :   std::vector<real> cell;
     583         576 :   std::vector<real> virial;
     584         576 :   std::vector<real> numder;
     585             : 
     586             : // variables to test particle decomposition
     587             :   int pd_nlocal;
     588             :   int pd_start;
     589             : // variables to test random decomposition (=domain decomposition)
     590         576 :   std::vector<int>  dd_gatindex;
     591         576 :   std::vector<int>  dd_g2l;
     592         576 :   std::vector<real> dd_masses;
     593         576 :   std::vector<real> dd_charges;
     594         576 :   std::vector<real> dd_forces;
     595         576 :   std::vector<real> dd_coordinates;
     596             :   int dd_nlocal;
     597             : // random stream to choose decompositions
     598         576 :   Random rnd;
     599             : 
     600             :   while(true) {
     601       19546 :     if(!noatoms) {
     602       19000 :       if(use_molfile==true) {
     603             : #ifdef __PLUMED_HAS_MOLFILE_PLUGINS
     604             :         int rc;
     605         186 :         rc = api->read_next_timestep(h_in, natoms, &ts_in);
     606             :         //if(rc==MOLFILE_SUCCESS){
     607             :         //       printf(" read this one :success \n");
     608             :         //}
     609         186 :         if(rc==MOLFILE_EOF) {
     610             :           //printf(" read this one :eof or error \n");
     611         302 :           break;
     612             :         }
     613             : #endif
     614       18814 :       } else if(trajectory_fmt=="xyz" || trajectory_fmt=="gro") {
     615       18762 :         if(!Tools::getline(fp,line)) break;
     616             :       }
     617             :     }
     618             : 
     619       19263 :     bool first_step=false;
     620       19263 :     if(!noatoms) {
     621       18717 :       if(use_molfile==false && (trajectory_fmt=="xyz" || trajectory_fmt=="gro")) {
     622       18493 :         if(trajectory_fmt=="gro") if(!Tools::getline(fp,line)) error("premature end of trajectory file");
     623       18493 :         sscanf(line.c_str(),"%100d",&natoms);
     624             :       }
     625             :     }
     626       19263 :     if(checknatoms<0 && !noatoms) {
     627         287 :       pd_nlocal=natoms;
     628         287 :       pd_start=0;
     629         287 :       first_step=true;
     630         287 :       masses.assign(natoms,NAN);
     631         287 :       charges.assign(natoms,NAN);
     632             : //case pdb: structure
     633         287 :       if(pdbfile.length()>0) {
     634        3271 :         for(unsigned i=0; i<pdb.size(); ++i) {
     635        3234 :           AtomNumber an=pdb.getAtomNumbers()[i];
     636        3234 :           unsigned index=an.index();
     637        3234 :           if( index>=unsigned(natoms) ) error("atom index in pdb exceeds the number of atoms in trajectory");
     638        3234 :           masses[index]=pdb.getOccupancy()[i];
     639        3234 :           charges[index]=pdb.getBeta()[i];
     640             :         }
     641             :       }
     642         287 :       if(mcfile.length()>0) {
     643           4 :         IFile ifile;
     644           4 :         ifile.open(mcfile);
     645             :         int index; double mass; double charge;
     646         440 :         while(ifile.scanField("index",index).scanField("mass",mass).scanField("charge",charge).scanField()) {
     647         432 :           masses[index]=mass;
     648         432 :           charges[index]=charge;
     649           4 :         }
     650         287 :       }
     651       18976 :     } else if( checknatoms<0 && noatoms ) {
     652           1 :       natoms=0;
     653             :     }
     654       19263 :     if( checknatoms<0 ) {
     655         288 :       checknatoms=natoms;
     656         288 :       p.cmd("setNatoms",&natoms);
     657         288 :       p.cmd("init");
     658             :     }
     659       19263 :     if(checknatoms!=natoms) {
     660           0 :       std::string stepstr; Tools::convert(step,stepstr);
     661           0 :       error("number of atoms in frame " + stepstr + " does not match number of atoms in first frame");
     662             :     }
     663             : 
     664       19263 :     coordinates.assign(3*natoms,real(0.0));
     665       19263 :     forces.assign(3*natoms,real(0.0));
     666       19263 :     cell.assign(9,real(0.0));
     667       19263 :     virial.assign(9,real(0.0));
     668             : 
     669       19263 :     if( first_step || rnd.U01()>0.5) {
     670        9819 :       if(debug_pd) {
     671          16 :         int npe=intracomm.Get_size();
     672          16 :         vector<int> loc(npe,0);
     673          32 :         vector<int> start(npe,0);
     674          64 :         for(int i=0; i<npe-1; i++) {
     675          48 :           int cc=(natoms*2*rnd.U01())/npe;
     676          48 :           if(start[i]+cc>natoms) cc=natoms-start[i];
     677          48 :           loc[i]=cc;
     678          48 :           start[i+1]=start[i]+loc[i];
     679             :         }
     680          16 :         loc[npe-1]=natoms-start[npe-1];
     681          16 :         intracomm.Bcast(loc,0);
     682          16 :         intracomm.Bcast(start,0);
     683          16 :         pd_nlocal=loc[intracomm.Get_rank()];
     684          16 :         pd_start=start[intracomm.Get_rank()];
     685          16 :         if(intracomm.Get_rank()==0) {
     686           4 :           fprintf(out,"\nDRIVER: Reassigning particle decomposition\n");
     687           4 :           fprintf(out,"DRIVER: "); for(int i=0; i<npe; i++) fprintf(out,"%d ",loc[i]); printf("\n");
     688           4 :           fprintf(out,"DRIVER: "); for(int i=0; i<npe; i++) fprintf(out,"%d ",start[i]); printf("\n");
     689             :         }
     690          16 :         p.cmd("setAtomsNlocal",&pd_nlocal);
     691          32 :         p.cmd("setAtomsContiguous",&pd_start);
     692        9803 :       } else if(debug_dd) {
     693         134 :         int npe=intracomm.Get_size();
     694         134 :         int rank=intracomm.Get_rank();
     695         134 :         dd_charges.assign(natoms,0.0);
     696         134 :         dd_masses.assign(natoms,0.0);
     697         134 :         dd_gatindex.assign(natoms,-1);
     698         134 :         dd_g2l.assign(natoms,-1);
     699         134 :         dd_coordinates.assign(3*natoms,0.0);
     700         134 :         dd_forces.assign(3*natoms,0.0);
     701         134 :         dd_nlocal=0;
     702       13196 :         for(int i=0; i<natoms; ++i) {
     703       13062 :           double r=rnd.U01()*npe;
     704       13062 :           int n; for(n=0; n<npe; n++) if(n+1>r)break;
     705       13062 :           plumed_assert(n<npe);
     706       13062 :           if(n==rank) {
     707        5883 :             dd_gatindex[dd_nlocal]=i;
     708        5883 :             dd_g2l[i]=dd_nlocal;
     709        5883 :             dd_charges[dd_nlocal]=charges[i];
     710        5883 :             dd_masses[dd_nlocal]=masses[i];
     711        5883 :             dd_nlocal++;
     712             :           }
     713             :         }
     714         134 :         if(intracomm.Get_rank()==0) {
     715          61 :           fprintf(out,"\nDRIVER: Reassigning particle decomposition\n");
     716             :         }
     717         134 :         p.cmd("setAtomsNlocal",&dd_nlocal);
     718         134 :         p.cmd("setAtomsGatindex",&dd_gatindex[0]);
     719             :       }
     720             :     }
     721             : 
     722       19263 :     int plumedStopCondition=0;
     723       19263 :     if(!noatoms) {
     724       18717 :       if(use_molfile) {
     725             : #ifdef __PLUMED_HAS_MOLFILE_PLUGINS
     726         172 :         if(pbc_cli_given==false) {
     727         172 :           if(ts_in.A>0.0) { // this is negative if molfile does not provide box
     728             :             // info on the cell: convert using pbcset.tcl from pbctools in vmd distribution
     729         170 :             real cosBC=cos(ts_in.alpha*pi/180.);
     730             :             //double sinBC=sin(ts_in.alpha*pi/180.);
     731         170 :             real cosAC=cos(ts_in.beta*pi/180.);
     732         170 :             real cosAB=cos(ts_in.gamma*pi/180.);
     733         170 :             real sinAB=sin(ts_in.gamma*pi/180.);
     734         170 :             real Ax=ts_in.A;
     735         170 :             real Bx=ts_in.B*cosAB;
     736         170 :             real By=ts_in.B*sinAB;
     737         170 :             real Cx=ts_in.C*cosAC;
     738         170 :             real Cy=(ts_in.C*ts_in.B*cosBC-Cx*Bx)/By;
     739         170 :             real Cz=sqrt(ts_in.C*ts_in.C-Cx*Cx-Cy*Cy);
     740         170 :             cell[0]=Ax/10.; cell[1]=0.; cell[2]=0.;
     741         170 :             cell[3]=Bx/10.; cell[4]=By/10.; cell[5]=0.;
     742         170 :             cell[6]=Cx/10.; cell[7]=Cy/10.; cell[8]=Cz/10.;
     743             :           } else {
     744           2 :             cell[0]=0.0; cell[1]=0.0; cell[2]=0.0;
     745           2 :             cell[3]=0.0; cell[4]=0.0; cell[5]=0.0;
     746           2 :             cell[6]=0.0; cell[7]=0.0; cell[8]=0.0;
     747             :           }
     748             :         } else {
     749           0 :           for(unsigned i=0; i<9; i++)cell[i]=pbc_cli_box[i];
     750             :         }
     751             :         // info on coords
     752             :         // the order is xyzxyz...
     753     1012447 :         for(unsigned i=0; i<3*natoms; i++) {
     754     1012275 :           coordinates[i]=real(ts_in.coords[i]/10.); //convert to nm
     755             :           //cerr<<"COOR "<<coordinates[i]<<endl;
     756             :         }
     757             : #endif
     758       18545 :       } else if(trajectory_fmt=="xdr-xtc" || trajectory_fmt=="xdr-trr") {
     759             : #ifdef __PLUMED_HAS_XDRFILE
     760             :         int localstep;
     761             :         float time;
     762             :         matrix box;
     763          52 :         rvec* pos=new rvec[natoms];
     764             :         float prec,lambda;
     765          52 :         int ret=exdrOK;
     766          52 :         if(trajectory_fmt=="xdr-xtc") ret=read_xtc(xd,natoms,&localstep,&time,box,pos,&prec);
     767          52 :         if(trajectory_fmt=="xdr-trr") ret=read_trr(xd,natoms,&localstep,&time,&lambda,box,pos,NULL,NULL);
     768          52 :         if(stride==0) step=localstep;
     769          56 :         if(ret==exdrENDOFFILE) break;
     770          51 :         if(ret!=exdrOK) break;
     771          48 :         for(unsigned i=0; i<3; i++) for(unsigned j=0; j<3; j++) cell[3*i+j]=box[i][j];
     772        7176 :         for(unsigned i=0; i<natoms; i++) for(unsigned j=0; j<3; j++)
     773        7128 :             coordinates[3*i+j]=real(pos[i][j]);
     774          48 :         delete [] pos;
     775             : #endif
     776             :       } else {
     777       18493 :         if(trajectory_fmt=="xyz") {
     778       18015 :           if(!Tools::getline(fp,line)) error("premature end of trajectory file");
     779             : 
     780       18015 :           std::vector<double> celld(9,0.0);
     781       18015 :           if(pbc_cli_given==false) {
     782       17815 :             std::vector<std::string> words;
     783       17815 :             words=Tools::getWords(line);
     784       17815 :             if(words.size()==3) {
     785       17357 :               sscanf(line.c_str(),"%100lf %100lf %100lf",&celld[0],&celld[4],&celld[8]);
     786         458 :             } else if(words.size()==9) {
     787         458 :               sscanf(line.c_str(),"%100lf %100lf %100lf %100lf %100lf %100lf %100lf %100lf %100lf",
     788             :                      &celld[0], &celld[1], &celld[2],
     789             :                      &celld[3], &celld[4], &celld[5],
     790         458 :                      &celld[6], &celld[7], &celld[8]);
     791           0 :             } else error("needed box in second line of xyz file");
     792             :           } else {                      // from command line
     793         200 :             celld=pbc_cli_box;
     794             :           }
     795       18015 :           for(unsigned i=0; i<9; i++)cell[i]=real(celld[i]);
     796             :         }
     797       18493 :         int ddist=0;
     798             :         // Read coordinates
     799     1373437 :         for(int i=0; i<natoms; i++) {
     800     1354944 :           bool ok=Tools::getline(fp,line);
     801     1354944 :           if(!ok) error("premature end of trajectory file");
     802             :           double cc[3];
     803     1354944 :           if(trajectory_fmt=="xyz") {
     804             :             char dummy[1000];
     805     1314026 :             int ret=std::sscanf(line.c_str(),"%999s %100lf %100lf %100lf",dummy,&cc[0],&cc[1],&cc[2]);
     806     1314026 :             if(ret!=4) error("cannot read line"+line);
     807       40918 :           } else if(trajectory_fmt=="gro") {
     808             :             // do the gromacs way
     809       40918 :             if(!i) {
     810             :               //
     811             :               // calculate the distance between dots (as in gromacs gmxlib/confio.c, routine get_w_conf )
     812             :               //
     813             :               const char      *p1, *p2, *p3;
     814         478 :               p1 = strchr(line.c_str(), '.');
     815         478 :               if (p1 == NULL) error("seems there are no coordinates in the gro file");
     816         478 :               p2 = strchr(&p1[1], '.');
     817         478 :               if (p2 == NULL) error("seems there is only one coordinates in the gro file");
     818         478 :               ddist = p2 - p1;
     819         478 :               p3 = strchr(&p2[1], '.');
     820         478 :               if (p3 == NULL)error("seems there are only two coordinates in the gro file");
     821         478 :               if (p3 - p2 != ddist)error("not uniform spacing in fields in the gro file");
     822             :             }
     823       40918 :             Tools::convert(line.substr(20,ddist),cc[0]);
     824       40918 :             Tools::convert(line.substr(20+ddist,ddist),cc[1]);
     825       40918 :             Tools::convert(line.substr(20+ddist+ddist,ddist),cc[2]);
     826           0 :           } else plumed_error();
     827     1354944 :           if(!debug_pd || ( i>=pd_start && i<pd_start+pd_nlocal) ) {
     828     1353324 :             coordinates[3*i]=real(cc[0]);
     829     1353324 :             coordinates[3*i+1]=real(cc[1]);
     830     1353324 :             coordinates[3*i+2]=real(cc[2]);
     831             :           }
     832             :         }
     833       18493 :         if(trajectory_fmt=="gro") {
     834         478 :           if(!Tools::getline(fp,line)) error("premature end of trajectory file");
     835         478 :           std::vector<string> words=Tools::getWords(line);
     836         478 :           if(words.size()<3) error("cannot understand box format");
     837         478 :           Tools::convert(words[0],cell[0]);
     838         478 :           Tools::convert(words[1],cell[4]);
     839         478 :           Tools::convert(words[2],cell[8]);
     840         478 :           if(words.size()>3) Tools::convert(words[3],cell[1]);
     841         478 :           if(words.size()>4) Tools::convert(words[4],cell[2]);
     842         478 :           if(words.size()>5) Tools::convert(words[5],cell[3]);
     843         478 :           if(words.size()>6) Tools::convert(words[6],cell[5]);
     844         478 :           if(words.size()>7) Tools::convert(words[7],cell[6]);
     845         478 :           if(words.size()>8) Tools::convert(words[8],cell[7]);
     846             :         }
     847             : 
     848             :       }
     849             : 
     850       18713 :       p.cmd("setStepLong",&step);
     851       18713 :       p.cmd("setStopFlag",&plumedStopCondition);
     852             : 
     853       18713 :       if(debug_dd) {
     854       10996 :         for(int i=0; i<dd_nlocal; ++i) {
     855       10764 :           int kk=dd_gatindex[i];
     856       10764 :           dd_coordinates[3*i+0]=coordinates[3*kk+0];
     857       10764 :           dd_coordinates[3*i+1]=coordinates[3*kk+1];
     858       10764 :           dd_coordinates[3*i+2]=coordinates[3*kk+2];
     859             :         }
     860         232 :         p.cmd("setForces",&dd_forces[0]);
     861         232 :         p.cmd("setPositions",&dd_coordinates[0]);
     862         232 :         p.cmd("setMasses",&dd_masses[0]);
     863         232 :         p.cmd("setCharges",&dd_charges[0]);
     864             :       } else {
     865             : // this is required to avoid troubles when the last domain
     866             : // contains zero atoms
     867             : // Basically, for empty domains we pass null pointers
     868             : #define fix_pd(xx) (pd_nlocal!=0?&xx:NULL)
     869       18481 :         p.cmd("setForces",fix_pd(forces[3*pd_start]));
     870       18481 :         p.cmd("setPositions",fix_pd(coordinates[3*pd_start]));
     871       18481 :         p.cmd("setMasses",fix_pd(masses[pd_start]));
     872       18481 :         p.cmd("setCharges",fix_pd(charges[pd_start]));
     873             :       }
     874       18713 :       p.cmd("setBox",&cell[0]);
     875       18713 :       p.cmd("setVirial",&virial[0]);
     876             :     } else {
     877         546 :       p.cmd("setStepLong",&step);
     878         546 :       p.cmd("setStopFlag",&plumedStopCondition);
     879             :     }
     880       19259 :     p.cmd("calc");
     881             : 
     882             : // this is necessary as only processor zero is adding to the virial:
     883       19259 :     intracomm.Bcast(virial,0);
     884       19259 :     if(debug_pd) intracomm.Sum(forces);
     885       19259 :     if(debug_dd) {
     886       10996 :       for(int i=0; i<dd_nlocal; i++) {
     887       10764 :         forces[3*dd_gatindex[i]+0]=dd_forces[3*i+0];
     888       10764 :         forces[3*dd_gatindex[i]+1]=dd_forces[3*i+1];
     889       10764 :         forces[3*dd_gatindex[i]+2]=dd_forces[3*i+2];
     890             :       }
     891         232 :       dd_forces.assign(3*natoms,0.0);
     892         232 :       intracomm.Sum(forces);
     893             :     }
     894       19259 :     if(debug_grex &&step%grex_stride==0) {
     895          72 :       p.cmd("GREX savePositions");
     896          72 :       if(intracomm.Get_rank()>0) {
     897          36 :         p.cmd("GREX prepare");
     898             :       } else {
     899          36 :         int r=intercomm.Get_rank();
     900          36 :         int n=intercomm.Get_size();
     901          36 :         int partner=r+(2*((r+step/grex_stride)%2))-1;
     902          36 :         if(partner<0)partner=0;
     903          36 :         if(partner>=n) partner=n-1;
     904          36 :         p.cmd("GREX setPartner",&partner);
     905          36 :         p.cmd("GREX calculate");
     906          36 :         p.cmd("GREX shareAllDeltaBias");
     907         144 :         for(int i=0; i<n; i++) {
     908         108 :           string s; Tools::convert(i,s);
     909         108 :           real a=NAN; s="GREX getDeltaBias "+s; p.cmd(s.c_str(),&a);
     910         108 :           if(grex_log) fprintf(grex_log," %f",a);
     911             :         }
     912          36 :         if(grex_log) fprintf(grex_log,"\n");
     913             :       }
     914             :     }
     915             : 
     916             : 
     917       19259 :     if(fp_forces) {
     918       10599 :       fprintf(fp_forces,"%d\n",natoms);
     919       10599 :       string fmtv=dumpforcesFmt+" "+dumpforcesFmt+" "+dumpforcesFmt+" "+dumpforcesFmt+" "+dumpforcesFmt+" "+dumpforcesFmt+" "+dumpforcesFmt+" "+dumpforcesFmt+" "+dumpforcesFmt+"\n";
     920       21198 :       string fmt=dumpforcesFmt+" "+dumpforcesFmt+" "+dumpforcesFmt+"\n";
     921       10599 :       if(dumpfullvirial) {
     922         117 :         fprintf(fp_forces,fmtv.c_str(),virial[0],virial[1],virial[2],virial[3],virial[4],virial[5],virial[6],virial[7],virial[8]);
     923             :       } else {
     924       10482 :         fprintf(fp_forces,fmt.c_str(),virial[0],virial[4],virial[8]);
     925             :       }
     926       10599 :       fmt="X "+fmt;
     927      477029 :       for(int i=0; i<natoms; i++)
     928      477029 :         fprintf(fp_forces,fmt.c_str(),forces[3*i],forces[3*i+1],forces[3*i+2]);
     929             :     }
     930       19259 :     if(debugforces.length()>0) {
     931             :       // Now call the routine to work out the derivatives numerically
     932           3 :       numder.assign(3*natoms+9,real(0.0)); real base=0;
     933           3 :       p.cmd("getBias",&base);
     934           3 :       if( fabs(base)<epsilon ) printf("WARNING: bias for configuration appears to be zero so debugging forces is trivial");
     935           3 :       evaluateNumericalDerivatives( step, p, coordinates, masses, charges, cell, base, numder );
     936             : 
     937             :       // And output everything to a file
     938           3 :       fp_dforces.fmtField(" " + dumpforcesFmt);
     939          36 :       for(int i=0; i<3*natoms; ++i) {
     940          33 :         fp_dforces.printField("parameter",(int)i);
     941          33 :         fp_dforces.printField("analytical",forces[i]);
     942          33 :         fp_dforces.printField("numerical",-numder[i]);
     943          33 :         fp_dforces.printField();
     944             :       }
     945             :       // And print the virial
     946          30 :       for(int i=0; i<9; ++i) {
     947          27 :         fp_dforces.printField("parameter",3*natoms+i );
     948          27 :         fp_dforces.printField("analytical",virial[i] );
     949          27 :         fp_dforces.printField("numerical",-numder[3*natoms+i]);
     950          27 :         fp_dforces.printField();
     951             :       }
     952             :     }
     953             : 
     954       19259 :     if(noatoms && plumedStopCondition) break;
     955             : 
     956       19258 :     step+=stride;
     957             :   }
     958         288 :   p.cmd("runFinalJobs");
     959             : 
     960         288 :   if(fp_forces) fclose(fp_forces);
     961         288 :   if(debugforces.length()>0) fp_dforces.close();
     962         288 :   if(fp && fp!=in)fclose(fp);
     963             : #ifdef __PLUMED_HAS_XDRFILE
     964         288 :   if(xd) xdrfile_close(xd);
     965             : #endif
     966             : #ifdef __PLUMED_HAS_MOLFILE_PLUGINS
     967         288 :   if(h_in) api->close_file_read(h_in);
     968         288 :   if(ts_in.coords) delete [] ts_in.coords;
     969             : #endif
     970         288 :   if(grex_log) fclose(grex_log);
     971             : 
     972         576 :   return 0;
     973             : }
     974             : 
     975             : template<typename real>
     976           3 : void Driver<real>::evaluateNumericalDerivatives( const long int& step, Plumed& p, const std::vector<real>& coordinates,
     977             :     const std::vector<real>& masses, const std::vector<real>& charges,
     978             :     std::vector<real>& cell, const double& base, std::vector<real>& numder ) {
     979             : 
     980           3 :   int natoms = coordinates.size() / 3; real delta = sqrt(epsilon);
     981           3 :   std::vector<Vector> pos(natoms); real bias=0;
     982           6 :   std::vector<real> fake_forces( 3*natoms ), fake_virial(9);
     983          14 :   for(int i=0; i<natoms; ++i) {
     984          11 :     for(unsigned j=0; j<3; ++j) pos[i][j]=coordinates[3*i+j];
     985             :   }
     986             : 
     987          14 :   for(int i=0; i<natoms; ++i) {
     988          44 :     for(unsigned j=0; j<3; ++j) {
     989          33 :       pos[i][j]=pos[i][j]+delta;
     990          33 :       p.cmd("setStepLong",&step);
     991          33 :       p.cmd("setPositions",&pos[0][0]);
     992          33 :       p.cmd("setForces",&fake_forces[0]);
     993          33 :       p.cmd("setMasses",&masses[0]);
     994          33 :       p.cmd("setCharges",&charges[0]);
     995          33 :       p.cmd("setBox",&cell[0]);
     996          33 :       p.cmd("setVirial",&fake_virial[0]);
     997          33 :       p.cmd("prepareCalc");
     998          33 :       p.cmd("performCalcNoUpdate");
     999          33 :       p.cmd("getBias",&bias);
    1000          33 :       pos[i][j]=coordinates[3*i+j];
    1001          33 :       numder[3*i+j] = (bias - base) / delta;
    1002             :     }
    1003             :   }
    1004             : 
    1005             :   // Create the cell
    1006           3 :   Tensor box( cell[0], cell[1], cell[2], cell[3], cell[4], cell[5], cell[6], cell[7], cell[8] );
    1007             :   // And the virial
    1008           6 :   Pbc pbc; pbc.setBox( box ); Tensor nvirial;
    1009          30 :   for(unsigned i=0; i<3; i++) for(unsigned k=0; k<3; k++) {
    1010          27 :       double arg0=box(i,k);
    1011          27 :       for(int j=0; j<natoms; ++j) pos[j]=pbc.realToScaled( pos[j] );
    1012          27 :       cell[3*i+k]=box(i,k)=box(i,k)+delta; pbc.setBox(box);
    1013          27 :       for(int j=0; j<natoms; j++) pos[j]=pbc.scaledToReal( pos[j] );
    1014          27 :       p.cmd("setStepLong",&step);
    1015          27 :       p.cmd("setPositions",&pos[0][0]);
    1016          27 :       p.cmd("setForces",&fake_forces[0]);
    1017          27 :       p.cmd("setMasses",&masses[0]);
    1018          27 :       p.cmd("setCharges",&charges[0]);
    1019          27 :       p.cmd("setBox",&cell[0]);
    1020          27 :       p.cmd("setVirial",&fake_virial[0]);
    1021          27 :       p.cmd("prepareCalc");
    1022          27 :       p.cmd("performCalcNoUpdate");
    1023          27 :       p.cmd("getBias",&bias);
    1024          27 :       cell[3*i+k]=box(i,k)=arg0; pbc.setBox(box);
    1025          27 :       for(int j=0; j<natoms; j++) for(unsigned n=0; n<3; ++n) pos[j][n]=coordinates[3*j+n];
    1026          27 :       nvirial(i,k) = ( bias - base ) / delta;
    1027             :     }
    1028           3 :   nvirial=-matmul(box.transpose(),nvirial);
    1029           6 :   for(unsigned i=0; i<3; i++) for(unsigned k=0; k<3; k++)  numder[3*natoms+3*i+k] = nvirial(i,k);
    1030             : 
    1031           3 : }
    1032             : 
    1033             : }
    1034        2523 : }

Generated by: LCOV version 1.13