diff --git a/src/bias/MetaD.cpp b/src/bias/MetaD.cpp index d6a97038baae02e79f6b0e99efcff42b590cab2c..5b0995647acc4722133ce7367e56a796d6415bfc 100644 --- a/src/bias/MetaD.cpp +++ b/src/bias/MetaD.cpp @@ -285,6 +285,7 @@ void MetaD::registerKeywords(Keywords& keys){ Bias::registerKeywords(keys); componentsAreNotOptional(keys); keys.addOutputComponent("bias","default","the instantaneous value of the bias potential"); + keys.addOutputComponent("acc","ACCELERATION","the metadynamics acceleration factor"); keys.use("ARG"); keys.add("compulsory","SIGMA","the widths of the Gaussian hills"); keys.add("compulsory","PACE","the frequency for hill addition"); @@ -535,10 +536,13 @@ isFirstStep(true) if(mw_n_>1){ + if(walkers_mpi) error("MPI version of multiple walkers is not compatible with filesystem version of multiple walkers"); log.printf(" %d multiple walkers active\n",mw_n_); log.printf(" walker id %d\n",mw_id_); log.printf(" reading stride %d\n",mw_rstride_); log.printf(" directory with hills files %s\n",mw_dir_.c_str()); + } else { + if(walkers_mpi) log.printf(" Multiple walkers active using MPI communnication\n"); } addComponent("bias"); componentIsNotPeriodic("bias"); @@ -584,7 +588,6 @@ isFirstStep(true) } } - if(walkers_mpi && mw_n_>1) error("MPI version of multiple walkers is not compatible with filesystem version of multiple walkers"); // creating vector of ifile* for hills reading // open all files at the beginning and read Gaussians if restarting @@ -624,7 +627,7 @@ isFirstStep(true) log<<" Bibliography "<<plumed.cite("Laio and Parrinello, PNAS 99, 12562 (2002)"); if(welltemp_) log<<plumed.cite( "Barducci, Bussi, and Parrinello, Phys. Rev. Lett. 100, 020603 (2008)"); - if(mw_n_>1) log<<plumed.cite( + if(mw_n_>1||walkers_mpi) log<<plumed.cite( "Raiteri, Laio, Gervasio, Micheletti, and Parrinello, J. Phys. Chem. B 110, 3533 (2006)"); if(adaptive_!=FlexibleBin::none) log<<plumed.cite( "Branduardi, Bussi, and Parrinello, J. Chem. Theory Comput. 8, 2247 (2012)"); diff --git a/user-doc/tutorials/belfast-8.txt b/user-doc/tutorials/belfast-8.txt index ee9e9e8c47edff52ededb0e6a093b84c09193d10..e71146c5dca733bf01c6d6471d90e60a46938471 100644 --- a/user-doc/tutorials/belfast-8.txt +++ b/user-doc/tutorials/belfast-8.txt @@ -294,10 +294,7 @@ INCLUDE FILE=plumed-common.dat METAD ... ARG=cv2,cv3 SIGMA=0.3,0.3 HEIGHT=0.2 PACE=100 LABEL=mw GRID_MIN=-pi,-pi GRID_MAX=pi,pi GRID_BIN=200,200 -WALKERS_N=4 -WALKERS_ID=# -WALKERS_DIR=./ -WALKERS_RSTRIDE=500 +WALKERS_MPI ... METAD PRINT ARG=cv1,cv2,cv3,cv4 STRIDE=1000 FILE=COLVAR @@ -309,6 +306,10 @@ and the simulation can be run in a similar way without doing exchanges: mpirun -np 4 mdrun_mpi -s topol -plumed plumed -multi 4 >& log & \endverbatim +alternatively Multiple Walkers can be run as independent simulations sharing via the file system their own biasing potential, +this is usefull because it provides a parallelisation that do not needs a parallel code. In this case the walkers read with +a given frequency the gaussian deposited by the others and add them to their own. + \Reference This tutorial is freely inspired to the work of Biarnes et al.