LCOV - code coverage report
Current view: top level - logmfd - LogMFD.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 298 334 89.2 %
Date: 2024-10-11 08:09:47 Functions: 14 15 93.3 %

          Line data    Source code
       1             : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
       2             : Copyright (c) 2019
       3             : National Institute of Advanced Industrial Science and Technology (AIST), Japan.
       4             : This file contains module for LogMFD method proposed by Tetsuya Morishita(AIST).
       5             : 
       6             : The LogMFD module is free software: you can redistribute it and/or modify
       7             : it under the terms of the GNU Lesser General Public License as published by
       8             : the Free Software Foundation, either version 3 of the License, or
       9             : (at your option) any later version.
      10             : 
      11             : The LogMFD module is distributed in the hope that it will be useful,
      12             : but WITHOUT ANY WARRANTY; without even the implied warranty of
      13             : MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
      14             : GNU Lesser General Public License for more details.
      15             : 
      16             : You should have received a copy of the GNU Lesser General Public License
      17             : along with plumed.  If not, see <http://www.gnu.org/licenses/>.
      18             : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
      19             : 
      20             : //+PLUMEDOC LOGMFDMOD_BIAS LOGMFD
      21             : /*
      22             : Used to perform LogMFD, LogPD, and TAMD/d-AFED.
      23             : 
      24             : \section LogMFD LogMFD
      25             : 
      26             : Consider a physical system of \f$N_q\f$ particles, for which the Hamiltonian is given as
      27             : 
      28             : \f[
      29             :   {H_{\rm MD}}\left( {\bf{\Gamma}} \right) = \sum\limits_{j = 1}^{{N_q}} {\frac{{{\bf{p}}_j^2}}{{2{m_j}}}}  + \Phi \left( {\bf{q}} \right)
      30             : \f]
      31             : 
      32             : where \f${\bf q}_j\f$, \f${\bf p}_j\f$ (\f$\bf{\Gamma}={\bf q},{\bf p}\f$), and \f$m_j\f$ are the position, momentum, and mass of particle \f$j\f$ respectively,
      33             : and \f$\Phi\f$ is the potential energy function for \f${\bf q}\f$.
      34             : The free energy \f$F({\bf X})\f$ as a function of a set of \f$N\f$ collective variables (CVs) is given as
      35             : 
      36             : \f{eqnarray*}{
      37             :   F\left( {{\bf X}} \right) &=&  - {k_B}T\log \int {\exp \left[ { - \beta {H_{\rm MD}}} \right]\prod\limits_{i = 1}^N {\delta \left( {{s_i}\left( {{\bf q}} \right) - {X_i}} \right)} d{\bf{\Gamma }}} \\
      38             :   &\simeq&  - {k_B}T\log \int {\exp \left[ { - \beta \left\{ {{H_{\rm MD}} + \sum\limits_i^N {\frac{{{K_i}}}{2}{{\left( {{s_i}\left( {{\bf q}} \right) - {X_i}} \right)}^2}} } \right\}} \right]} d{\bf{\Gamma }}
      39             : \f}
      40             : 
      41             : where \f$\bf{s}\f$ are the CVs, \f$k_B\f$ is Boltzmann constant, \f$\beta=k_BT\f$,
      42             : and \f$K_i\f$ is the spring constant which is large enough to invoke
      43             : 
      44             : \f[
      45             :  \delta \left( x \right) = \lim_{k \to \infty } \sqrt {\beta k/2\pi} \exp \left( -\beta kx^2/2 \right)
      46             : \f]
      47             : 
      48             : In mean-force dynamics, \f${\bf X}\f$ are treated as fictitious dynamical variables, which are associated with the following Hamiltonian,
      49             : 
      50             : \f[
      51             :  {H_{\rm log}} = \sum\limits_{i = 1}^N {\frac{{P_{{X_i}}^2}}{{2{M_i}}} + \Psi \left( {{\bf X}} \right)}
      52             : \f]
      53             : 
      54             : where \f${P_{{X_i}}}\f$  and \f$M_i\f$ are the momentum and mass of \f$X_i\f$, respectively, and \f$\Psi\f$ is the potential function for \f${\bf X}\f$.
      55             : We assume that \f$\Psi\f$ is a functional of \f$F({\bf X})\f$; \f$Ψ[F({\bf X})]\f$. The simplest form of \f$\Psi\f$ is \f$\Psi = F({\bf X})\f$,
      56             : which corresponds to TAMD/d-AFED \cite AbramsJ2008, \cite Maragliano2006 (or the extended Lagrangian dynamics in the limit of the adiabatic decoupling between \f$\bf{q}\f$ and \f${\bf X}\f$).
      57             :  In the case of LogMFD, the following form of \f$\Psi\f$ is introduced \cite MorishitaLogMFD;
      58             : 
      59             : 
      60             : \f[
      61             :   {\Psi _{\rm log}}\left( {{\bf X}} \right) = \gamma \log \left[ {\alpha F\left( {{\bf X}} \right) + 1} \right]
      62             : \f]
      63             : 
      64             : where \f$\alpha\f$ (ALPHA) and \f$\gamma\f$ (GAMMA) are positive parameters. The logarithmic form of \f$\Psi_{\rm log}\f$ ensures the dynamics of \f${\bf X}\f$ on a much smoother energy surface [i.e., \f$\Psi_{\rm log}({\bf X})\f$] than \f$F({\bf X})\f$, thus enhancing the sampling in the \f${\bf X}\f$-space. The parameters \f$\alpha\f$ and \f$\gamma\f$ determine the degree of flatness of \f$\Psi_{\rm log}\f$, but adjusting only \f$\alpha\f$ is normally sufficient to have a relatively flat surface (with keeping the relation \f$\gamma=1/\alpha\f$).
      65             : 
      66             : The equation of motion for \f$X_i\f$ in LogMFD (no thermostat) is
      67             : 
      68             : \f[
      69             :  {M_i}{\ddot X_i} =  - \left( {\frac{{\alpha \gamma }}{{\alpha F + 1}}} \right)\frac{{\partial F}}{{\partial {X_i}}}
      70             : \f]
      71             : 
      72             : where \f$-\partial F/\partial X_i\f$  is evaluated as a canonical average under the condition that \f${\bf X}\f$ is fixed;
      73             : 
      74             : \f{eqnarray*}{
      75             :  - \frac{{\partial F}}{{\partial {X_i}}} &\simeq& \frac{1}{Z}\int {{K_i}\left( {{s_i}\left( {{\bf q}} \right) - {X_i}} \right)\exp \left[ { - \beta \left\{ {{H_{\rm MD}} + \sum\limits_i^N {\frac{{{K_i}}}{2}{{\left( {{s_i}\left( {{\bf q}} \right) - {X_i}} \right)}^2}} } \right\}} \right]} d{\bf{\Gamma }}\\
      76             :  &\equiv& {\left\langle {{K_i}\left( {{s_i}\left( {{\bf q}} \right) - {X_i}} \right)} \right\rangle _{{\bf X}}}
      77             : \f}
      78             : 
      79             : where
      80             : 
      81             : \f[
      82             :  Z = \int {\exp \left[ { - \beta \left\{ {{H_{\rm MD}} + \sum\limits_i^N {\frac{{{K_i}}}{2}{{\left( {{s_i}\left( {{\bf q}} \right) - {X_i}} \right)}^2}} } \right\}} \right]} d{\bf{\Gamma }}
      83             : \f]
      84             : 
      85             : The mean-force (MF) is practically evaluated by performing a shot-time canonical MD run each time \f${\bf X}\f$ is updated according to the equation of motion for \f${\bf X}\f$.
      86             : 
      87             : If the canonical average for the MF is effectively converged, the dynamical variables \f${\bf q}\f$ and \f${\bf X}\f$ are decoupled and they evolve adiabatically, which can be exploited for the on-the-fly evaluation of \f$F({\bf X})\f$. I.e., \f$H_{\rm log}\f$ should be a constant of motion in this case, thus \f$F({\bf X})\f$ can be evaluated each time \f${\bf X}\f$ is updated as
      88             : 
      89             : 
      90             : \f[
      91             :  F\left( {{{\bf X}}\left( t \right)} \right) = \frac{1}{\alpha} \left[
      92             :   \exp \frac{1}{\gamma} \left\{ \left( H_{\rm log} - \sum_i \frac{P_{X_i}^2}{2M_i} \right) \right\} - 1 \right]
      93             : \f]
      94             : 
      95             : 
      96             : This means that \f$F({\bf X})\f$ can be constructed without post processing (on-the-fly free energy reconstruction). Note that the on-the-fly free energy reconstruction is also possible in TAMD/d-AFED if the Hamiltonian-like conserved quantity is available (e.g., the Nose-Hoover type dynamics).
      97             : 
      98             : 
      99             : 
     100             : \section LogPD LogPD
     101             : 
     102             : 
     103             : The accuracy in the MF is critical to the on-the-fly free energy reconstruction. To improve the evaluation of the MF, parallel-dynamics (PD) is incorporated into LogMFD, leading to logarithmic parallel-dynamics (LogPD) \cite MorishitaLogPD.
     104             : 
     105             : 
     106             : In PD, the MF is evaluated by a non-equilibrium path-ensemble based on the Crooks-Jarzynski non-equilibrium work relation. To this end, multiple replicas of the MD system which run in parallel are introduced. The CVs [\f${\bf s}({\bf q})\f$] in each replica is restrained to the same value of \f${\bf X}(t)\f$. A canonical MD run with \f$N_m\f$ steps is performed in each replica, then the MF on \f$X_i\f$ is evaluated using the MD trajectories from all replicas.
     107             : The MF is practically calculated as
     108             : 
     109             : 
     110             : \f[
     111             :  - \frac{{\partial F}}{{\partial {X_i}}} = \sum\limits_{k = 1}^{{N_r}} {{W_k}} \sum\limits_{n = 1}^{{N_m}} {\frac{1}{{{N_m}}}{K_i}\left[ {{s_i}\left( {{{\bf q}}_n^k} \right) - {X_i}} \right]}
     112             : \f]
     113             : 
     114             : 
     115             : 
     116             : where \f${\bf q}^k_n\f$  indicates the \f${\bf q}\f$-configuration at time step \f$n\f$ in the \f$N_m\f$-step shot-time MD run in the \f$k\f$th replica among \f$N_r\f$ replicas. \f$W_k\f$ is given as
     117             : 
     118             : \f[
     119             :  {W_k} = \frac{{{e^{ - \beta {w_k}\left( t \right)}}}}{{\sum\limits_{k=1}^{{N_r}} {{e^{ - \beta {w_k}\left( t \right)}}} }}
     120             : \f]
     121             : 
     122             : 
     123             : where
     124             : 
     125             : 
     126             : \f[
     127             :  {w_k}\left( t \right) = \int\limits_{{X_0}}^{X\left( t \right)} {\sum\limits_{i=1}^N {\frac{{\partial H_{\rm MD}^k}}{{\partial {X_i}}}d{X_i}} }
     128             : \f]
     129             : 
     130             : \f[
     131             :  H_{\rm MD}^k\left( {{\bf{\Gamma }},{{\bf X}}} \right) = {H_{\rm MD}}\left( {{{\bf{\Gamma }}^k}} \right) + \sum\limits_{i = 1}^N {\frac{{{K_i}}}{2}{{\left( {s_i^k - {X_i}} \right)}^2}}
     132             : \f]
     133             : 
     134             : and \f$s^k_i\f$ is the \f$i\f$th CV in the \f$k\f$th replica.
     135             : 
     136             : \f$W_k\f$ comes from the Crooks-Jarzynski non-equilibrium work relation by which we can evaluate an equilibrium ensemble average from a set of non-equilibrium trajectories. Note that, to avoid possible numerical errors in the exponential function, the following form of \f$W_k\f$ is instead used in PLUMED,
     137             : 
     138             : \f[
     139             :  {W_k}\left( t \right) = \frac{{\exp \left[ { - \beta \left\{ {{w_k}\left( t \right) - {w_{\min }}\left( t \right)} \right\}} \right]}}{{\sum\nolimits_k {\exp \left[ { - \beta \left\{ {{w_k}\left( t \right) - {w_{\min }}\left( t \right)} \right\}} \right]} }}
     140             : \f]
     141             : 
     142             : where
     143             : 
     144             : \f[
     145             :  {w_{\min }}\left( t \right) = {\rm Min}_k\left[ {{w_k}\left( t \right)} \right]
     146             : \f]
     147             : 
     148             : 
     149             : With the MF evaluated using the PD approach, reconstructing free energy profiles can be performed more efficiently (requiring less elapsed computing time) in LogPD than with a single MD system in LogMFD. In the case that there exists more than one stable state separated by high energy barriers in the hidden subspace orthogonal to the CV-subspace, LogPD is particularly of use to incorporate all the contributions from such hidden states with appropriate weights (in the limit \f$N_r\to\infty\f$ ).
     150             : 
     151             : Note that LogPD calculations should always be initiated with an equilibrium \f${\bf q}\f$-configuration in each replica, because the Crooks-Jarzynski non-equilibrium work relation is invoked. Also note that LogPD is currently available only with Gromacs, while LogMFD can be performed with LAMMPS, Gromacs, and NAMD.
     152             : 
     153             : \section Thermostat Using LogMFD/PD with a thermostat
     154             : 
     155             : Introducing a thermostat on \f${\bf X}\f$ is often recommended in LogMFD/PD to maintain the adiabatic decoupling between \f${\bf q}\f$ and \f${\bf X}\f$. In the framework of the LogMFD approach, the Nose-Hoover type thermostat and the Gaussian isokinetic (velocity scaling) thermostat can be used to control the kinetic energy of \f${\bf X}\f$.
     156             : 
     157             : \subsection Nose-Hoover Nose-Hoover LogMFD/PD
     158             : 
     159             : The equation of motion for \f$X_i\f$ coupled to a Nose-Hoover thermostat variable \f$\eta\f$ (single heat bath) is
     160             : 
     161             : \f[
     162             :  {M_i}{\ddot X_i} =  - \left( {\frac{{\alpha \gamma }}{{\alpha F + 1}}} \right)\frac{{\partial F}}{{\partial {X_i}}} - {M_i}{\dot X_i}\dot \eta
     163             : \f]
     164             : 
     165             : The equation of motion for \f$\eta\f$ is
     166             : 
     167             : \f[
     168             :  Q\ddot \eta  = \sum\limits_{i = 1}^N {\frac{{P_{{X_i}}^2}}{{{M_i}}} - N{k_B}T}
     169             : \f]
     170             : 
     171             : where \f$N\f$ is the number of the CVs. Since the following pseudo-Hamiltonian is a constant of motion in Nose-Hoover LogMFD/PD,
     172             : 
     173             : \f[
     174             :  H_{\rm log}^{\rm NH} = \sum\limits_{i = 1}^N {\frac{{P_{{X_i}}^2}}{{2{M_i}}} + \gamma \log \left[ {\alpha F\left( {{\bf X}} \right) + 1} \right] + \frac{1}{2}Q{{\dot \eta }^2} + N{k_B}T\eta}
     175             : \f]
     176             : 
     177             : \f$F({\bf X}(t))\f$ is obtained at each MFD step as
     178             : 
     179             : \f[
     180             :  F\left( {{{\bf X}}\left( t \right)} \right) = \frac{1}{\alpha }\left[ {\exp \left\{ {{{ \frac{1}{\gamma} \left( {H_{\rm log}^{\rm NH} - \sum_i {\frac{{P_{{X_i}}^2}}{{2{M_i}}}}  - \frac{1}{2}Q{{\dot \eta }^2} - N{k_B}T\eta} \right)}  }} \right\} - 1} \right]
     181             : \f]
     182             : 
     183             : 
     184             : 
     185             : \subsection VS Velocity scaling LogMFD/PD
     186             : 
     187             : The velocity scaling algorithm (which is equivalent to the Gaussian isokinetic dynamics in the limit \f$\Delta t\to 0\f$) can also be employed to control the velocity of \f${\bf X}\f$, \f$\bf{V}_x\f$.
     188             : 
     189             : The following algorithm is introduced to perform isokinetic LogMFD calculations \cite MorishitaVsLogMFD;
     190             : 
     191             : \f{eqnarray*}{
     192             : {V_{{X_i}}}\left( {{t_{n + 1}}} \right) &=& V_{X_i}^\prime \left( {{t_n}} \right) + \Delta t\left[ { - \left( {\frac{{\alpha \gamma }}{{\alpha F\left( {{t_n}} \right) + 1}}} \right)\frac{{\partial F\left( {{t_n}} \right)}}{{\partial {X_i}}}} \right]\\
     193             : S\left( {{t_{n + 1}}} \right) &=& \sqrt {\frac{{N{k_B}T}}{{\sum\limits_i {{M_i}V_{{X_i}}^2\left( {{t_{n + 1}}} \right)} }}} \\
     194             : {V_{{X_i}}}^\prime \left( {{t_{n + 1}}} \right) &=& S\left( {{t_{n + 1}}} \right){V_{{X_i}}}\left( {{t_{n + 1}}} \right)\\
     195             : {X_i}\left( {{t_{n + 1}}} \right) &=& {X_i}\left( {{t_n}} \right) + \Delta t V_{X_i}^\prime \left( {{t_{n + 1}}} \right)\\
     196             : {\Psi_{\rm log}}\left( {{t_{n + 1}}} \right) &=& N{k_B}T\log S\left( {{t_{n + 1}}} \right) + {\Psi_{\rm log}}\left( {{t_n}} \right)\\
     197             : F\left( {{t_{n + 1}}} \right) &=& \frac{1}{\alpha} \left[
     198             :     \exp \left\{ \Psi_{\rm log} \left( t_{n+1} \right) / \gamma \right\} - 1 \right]
     199             : \f}
     200             : 
     201             : Note that \f$V_{X_i}^\prime\left( {{t_0}} \right)\f$ is assumed to be initially given, which meets the following relation,
     202             : 
     203             : \f[
     204             :   \sum\limits_{i = 1}^N M_i V_{X_i}^{\prime 2} \left( t_0 \right)  = N{k_B}{T}
     205             : \f]
     206             : 
     207             : It should be stressed that, in the same way as in the NVE and Nose-Hoover LogMFD/PD, \f$F({\bf X}(t))\f$ can be evaluated at each MFD step (on-the-fly free energy reconstruction) in Velocity Scaling LogMFD/PD.
     208             : 
     209             : 
     210             : \par Examples
     211             : \section Examples Examples
     212             : 
     213             : \subsection Example-LoGMFD Example of LogMFD
     214             : 
     215             : The following input file tells plumed to restrain collective variables
     216             : to the fictitious dynamical variables in LogMFD/PD.
     217             : 
     218             : plumed.dat
     219             : \plumedfile
     220             : UNITS TIME=fs LENGTH=1.0 ENERGY=kcal/mol MASS=1.0 CHARGE=1.0
     221             : phi: TORSION ATOMS=5,7,9,15
     222             : psi: TORSION ATOMS=7,9,15,17
     223             : 
     224             : # LogMFD
     225             : LOGMFD ...
     226             : LABEL=logmfd
     227             : ARG=phi,psi
     228             : KAPPA=100.0,100.0
     229             : DELTA_T=0.5
     230             : INTERVAL=500
     231             : TEMP=300.0
     232             : FLOG=5.0
     233             : MFICT=5000000.0,5000000.0
     234             : VFICT=3.5e-4,3.5e-4
     235             : ALPHA=4.0
     236             : THERMOSTAT=NVE
     237             : VETA=0.0
     238             : META=20000.0
     239             : FICT_MAX=3.1,3.1
     240             : FICT_MIN=-3.1,-3.1
     241             : ... LOGMFD
     242             : \endplumedfile
     243             : 
     244             : To submit this simulation with Gromacs, use the following command line
     245             : to execute a LogMFD run with Gromacs-MD.
     246             : Here TOPO/topol0.tpr is an input file
     247             : which contains atomic coordinates and Gromacs parameters.
     248             : 
     249             : \verbatim
     250             : gmx_mpi mdrun -s TOPO/topol0.tpr -plumed
     251             : \endverbatim
     252             : 
     253             : This command will output files named logmfd.out and replica.out.
     254             : 
     255             : The output file logmfd.out records free energy and all fictitious dynamical variables at every MFD step.
     256             : 
     257             : logmfd.out
     258             : 
     259             : \verbatim
     260             : # LogMFD
     261             : # CVs : phi psi
     262             : # Mass for CV particles : 5000000.000000000 5000000.000000000
     263             : # Mass for thermostat   :   20000.000000000
     264             : # 1:iter_md, 2:Flog, 3:2*Ekin/gkb[K], 4:eta, 5:Veta,
     265             : # 6:phi_fict(t), 7:phi_vfict(t), 8:phi_force(t),
     266             : # 9:psi_fict(t), 10:psi_vfict(t), 11:psi_force(t),
     267             :        1       4.99918574     308.24149708       0.00000000       0.00000000      -2.85605938       0.00035002       5.19074544       2.79216364       0.00035000      -0.53762989
     268             :        2       4.99836196     308.26124159       0.00000000       0.00000000      -2.85588436       0.00035005       4.71247605       2.79233863       0.00035000      -0.00532474
     269             :        3       4.99743572     308.28344595       0.00000000       0.00000000      -2.85570932       0.00035007       5.34358230       2.79251363       0.00035000      -0.05119816
     270             : ...
     271             : \endverbatim
     272             : 
     273             : The output file replica.out records all collective variables at every MFD step.
     274             : 
     275             : replica.out
     276             : 
     277             : \verbatim
     278             : # Replica No. 0 of 1.
     279             : # 1:iter_md, 2:work, 3:weight,
     280             : # 4:phi(q)
     281             : # 5:psi(q)
     282             :        1   -8.142952e-04     1.000000e+00       -2.80432694       2.78661234
     283             :        2   -1.638105e-03     1.000000e+00       -2.80893462       2.79211039
     284             :        3   -2.564398e-03     1.000000e+00       -2.80244854       2.79182665
     285             : ...
     286             : \endverbatim
     287             : 
     288             : \subsection Example-LogPD Example of LogPD
     289             : 
     290             : Use the following command line to execute a LogPD run using two MD replicas (note that only Gromacs is currently available for LogPD).
     291             : Here TOPO/topol0.tpr and TOPO/topol1.tpr are input files
     292             : which contain atomic coordinates of each replica and Gromacs parameters.
     293             : 
     294             : \verbatim
     295             : mpirun -np 2 gmx_mpi mdrun -s TOPO/topol -plumed -multi 2
     296             : \endverbatim
     297             : 
     298             : This command will output files named logmfd.out, replica.out.0 and replica.out.1.
     299             : 
     300             : The output file logmfd.out records free energy and all fictitious dynamical variables at every MFD step.
     301             : 
     302             : logmfd.out
     303             : 
     304             : \verbatim
     305             : # LogPD, replica parallel of LogMFD
     306             : # number of replica : 2
     307             : # CVs : phi psi
     308             : # Mass for CV particles : 5000000.000000000 5000000.000000000
     309             : # Mass for thermostat   :   20000.000000000
     310             : # 1:iter_md, 2:Flog, 3:2*Ekin/gkb[K], 4:eta, 5:Veta,
     311             : # 6:phi_fict(t), 7:phi_vfict(t), 8:phi_force(t),
     312             : # 9:psi_fict(t), 10:psi_vfict(t), 11:psi_force(t),
     313             :        1       5.00224715     308.16814691       0.00000000       0.00000000      -0.95937173       0.00034994     -12.91277494       0.78923967       0.00035000       0.07353010
     314             :        2       5.00476934     308.10774854       0.00000000       0.00000000      -0.95919679       0.00034989     -11.20093553       0.78941467       0.00034999      -3.21098229
     315             :        3       5.00702463     308.05376594       0.00000000       0.00000000      -0.95902187       0.00034983     -10.81712171       0.78958965       0.00034998      -2.07196718
     316             : ...
     317             : \endverbatim
     318             : 
     319             : 
     320             : The output file replica.out.0 records all collective variables of replica No.0 at every MFD step.
     321             : 
     322             : replica.out.0
     323             : 
     324             : \verbatim
     325             : # Replica No. 0 of 2.
     326             : # 1:iter_md, 2:work, 3:weight,
     327             : # 4:phi(q)
     328             : # 5:psi(q)
     329             :        1    1.843110e-03     5.003389e-01       -1.10929125       0.83348865
     330             :        2    3.466179e-03     5.010942e-01       -1.05020764       0.78731283
     331             :        3    4.927870e-03     5.017619e-01       -1.04968867       0.79635198
     332             : ...
     333             : \endverbatim
     334             : 
     335             : The output file replica.out.1 records all collective variables of replica No.1 at every MFD step.
     336             : 
     337             : replica.out.1
     338             : 
     339             : \verbatim
     340             : # Replica No. 1 of 2.
     341             : # 1:iter_md, 2:work, 3:weight,
     342             : # 4:phi(q)
     343             : # 5:psi(q)
     344             :        1    2.651173e-03     4.996611e-01       -1.06802968       0.74605205
     345             :        2    6.075530e-03     4.989058e-01       -1.09264741       0.72681448
     346             :        3    9.129358e-03     4.982381e-01       -1.08517238       0.74084241
     347             : ...
     348             : \endverbatim
     349             : 
     350             : */
     351             : //+ENDPLUMEDOC
     352             : 
     353             : #include "bias/Bias.h"
     354             : #include "core/ActionRegister.h"
     355             : #include "core/Atoms.h"
     356             : #include "core/PlumedMain.h"
     357             : 
     358             : #include <iostream>
     359             : 
     360             : using namespace std;
     361             : using namespace PLMD;
     362             : using namespace bias;
     363             : 
     364             : namespace PLMD {
     365             : namespace logmfd {
     366             : /**
     367             :    \brief class for LogMFD parameters, variables and subroutines.
     368             :  */
     369             : class LogMFD : public Bias {
     370             :   bool firsttime;               ///< flag that indicates first MFD step or not.
     371             :   int    interval;              ///< input parameter, period of MD steps when fictitious dynamical variables are evolved.
     372             :   double delta_t;               ///< input parameter, one time step of MFD when fictitious dynamical variables are evolved.
     373             :   string thermostat;            ///< input parameter, type of thermostat for canonical dyanamics.
     374             :   double kbt;                   ///< k_B*temperature
     375             : 
     376             :   int    TAMD;                  ///< input parameter, perform TAMD instead of LogMFD.
     377             :   double alpha;                 ///< input parameter, alpha parameter for LogMFD.
     378             :   double gamma;                 ///< input parameter, gamma parameter for LogMFD.
     379             :   std::vector<double> kappa;    ///< input parameter, strength of the harmonic restraining potential.
     380             : 
     381             :   std::vector<double> fict_max; ///< input parameter, maximum of each fictitous dynamical variable.
     382             :   std::vector<double> fict_min; ///< input parameter, minimum of each fictitous dynamical variable.
     383             : 
     384             :   std::vector<double>  fict;    ///< current values of each fictitous dynamical variable.
     385             :   std::vector<double> vfict;    ///< current velocity of each fictitous dynamical variable.
     386             :   std::vector<double> mfict;    ///< mass of each fictitous dynamical variable.
     387             : 
     388             :   double xeta;                  ///< current eta variable of thermostat.
     389             :   double veta;                  ///< current velocity of eta variable of thermostat.
     390             :   double meta;                  ///< mass of eta variable of thermostat.
     391             : 
     392             :   double flog;                  ///< current free energy
     393             : 
     394             :   double hlog;                  ///< value invariant
     395             :   double phivs;                 ///< potential used in VS method
     396             :   double work;                  ///< current works done by fictitious dynamical variables in this replica.
     397             :   double weight;                ///< current weight of this replica.
     398             : 
     399             :   std::vector<double> ffict;    ///< current force of each fictitous dynamical variable.
     400             :   std::vector<double> fict_ave; ///< averaged values of each collective variable.
     401             : 
     402             :   std::vector<Value*>  fictValue; ///< pointers to fictitious dynamical variables
     403             :   std::vector<Value*> vfictValue; ///< pointers to velocity of fictitious dynamical variables
     404             : 
     405             : public:
     406             :   static void registerKeywords(Keywords& keys);
     407             : 
     408             :   explicit LogMFD(const ActionOptions&);
     409             :   void calculate();
     410             :   void update();
     411             :   void updateNVE();
     412             :   void updateNVT();
     413             :   void updateVS();
     414             :   void calcMeanForce();
     415             :   double calcEkin();
     416             :   double calcFlog();
     417             :   double calcClog();
     418             : 
     419             : private:
     420             :   double sgn( double x ) {
     421          55 :     return x>0.0 ? 1.0 : x<0.0 ? -1.0 : 0.0;
     422             :   }
     423             : };
     424             : 
     425       10425 : PLUMED_REGISTER_ACTION(LogMFD,"LOGMFD")
     426             : 
     427             : /**
     428             :    \brief instruction of parameters for Plumed manual.
     429             : */
     430           4 : void LogMFD::registerKeywords(Keywords& keys) {
     431           4 :   Bias::registerKeywords(keys);
     432           4 :   keys.use("ARG");
     433           8 :   keys.add("compulsory","INTERVAL",
     434             :            "Period of MD steps (\\f$N_m\\f$) to update fictitious dynamical variables." );
     435           8 :   keys.add("compulsory","DELTA_T",
     436             :            "Time step for the fictitious dynamical variables (MFD step)." );
     437           8 :   keys.add("compulsory","THERMOSTAT",
     438             :            "Type of thermostat for the fictitious dynamical variables. NVE, NVT, VS are available." );
     439           8 :   keys.add("optional","TEMP",
     440             :            "Temperature of the fictitious dynamical variables in LogMFD/PD thermostat. "
     441             :            "If not provided or provided as 0, it will be taken from the temperature of the MD system." );
     442             : 
     443           8 :   keys.add("optional","TAMD",
     444             :            "When TAMD=1, TAMD/d-AFED calculations can be performed instead of LogMFD. Otherwise, the LogMFD protocol is switched on (default)." );
     445             : 
     446           8 :   keys.add("optional","ALPHA",
     447             :            "Alpha parameter for LogMFD. "
     448             :            "If not provided or provided as 0, it will be taken as 1/gamma. "
     449             :            "If gamma is also not provided, Alpha is set as 4, which is a sensible value when the unit of kcal/mol is used." );
     450           8 :   keys.add("optional","GAMMA",
     451             :            "Gamma parameter for LogMFD. "
     452             :            "If not provided or provided as 0, it will be taken as 1/alpha. "
     453             :            "If alpha is also not provided, Gamma is set as 0.25, which is a sensible value when the unit of kcal/mol is used." );
     454           8 :   keys.add("compulsory","KAPPA",
     455             :            "Spring constant of the harmonic restraining potential for the fictitious dynamical variables." );
     456             : 
     457           8 :   keys.add("compulsory","FICT_MAX",
     458             :            "Maximum values reachable for the fictitious dynamical variables. The variables will elastically bounce back at the boundary (mirror boundary)." );
     459           8 :   keys.add("compulsory","FICT_MIN",
     460             :            "Minimum values reachable for the fictitious dynamical variables. The variables will elastically bounce back at the boundary (mirror boundary)." );
     461             : 
     462           8 :   keys.add("optional","FICT",
     463             :            "The initial values of the fictitious dynamical variables. "
     464             :            "If not provided, they are set equal to their corresponding CVs for the initial atomic configuration." );
     465           8 :   keys.add("optional","VFICT",
     466             :            "The initial velocities of the fictitious dynamical variables. "
     467             :            "If not provided, they will be taken as 0." );
     468           8 :   keys.add("optional","MFICT",
     469             :            "Masses of each fictitious dynamical variable. "
     470             :            "If not provided, they will be taken as 10000." );
     471             : 
     472           8 :   keys.add("optional","XETA",
     473             :            "The initial eta variable of the Nose-Hoover thermostat "
     474             :            "for the fictitious dynamical variables. "
     475             :            "If not provided, it will be taken as 0." );
     476           8 :   keys.add("optional","VETA",
     477             :            "The initial velocity of eta variable. "
     478             :            "If not provided, it will be taken as 0." );
     479           8 :   keys.add("optional","META",
     480             :            "Mass of eta variable. "
     481             :            "If not provided, it will be taken as \\f$N*kb*T*100*100\\f$." );
     482             : 
     483           8 :   keys.add("compulsory","FLOG",
     484             :            "The initial free energy value in the LogMFD/PD run."
     485             :            "The origin of the free energy profile is adjusted by FLOG to "
     486             :            "realize \\f$F({\\bf X}(t)) > 0\\f$ at any \\f${\\bf X}(t)\\f$, "
     487             :            "resulting in enhanced barrier-crossing. "
     488             :            "(The value of \\f$H_{\\rm log}\\f$ is automatically "
     489             :            "set according to FLOG).");
     490             : 
     491           8 :   keys.add("optional","WORK",
     492             :            "The initial value of the work done by fictitious dynamical "
     493             :            "variables in each replica. "
     494             :            "If not provided, it will be taken as 0.");
     495             : 
     496           4 :   componentsAreNotOptional(keys);
     497           8 :   keys.addOutputComponent("_fict","default",
     498             :                           "For example, the fictitious collective variable for LogMFD is specified as "
     499             :                           "ARG=dist12 and LABEL=logmfd in LOGMFD section in Plumed input file, "
     500             :                           "the associated fictitious dynamical variable can be specified as "
     501             :                           "PRINT ARG=dist12,logmfd.dist12_fict FILE=COLVAR");
     502           8 :   keys.addOutputComponent("_vfict","default",
     503             :                           "For example, the fictitious collective variable for LogMFD is specified as "
     504             :                           "ARG=dist12 and LABEL=logmfd in LOGMFD section in Plumed input file, the "
     505             :                           "velocity of the associated fictitious dynamical variable can be specified as "
     506             :                           "PRINT ARG=dist12,logmfd.dist12_vfict FILE=COLVAR");
     507           4 : }
     508             : 
     509             : 
     510             : /**
     511             :    \brief constructor of LogMFD class
     512             :    \details This constructor initializes all parameters and variables,
     513             :    reads input file and set value of parameters and initial value of variables,
     514             :    and writes messages as Plumed log.
     515             : */
     516           3 : LogMFD::LogMFD( const ActionOptions& ao ):
     517             :   PLUMED_BIAS_INIT(ao),
     518           3 :   firsttime(true),
     519           3 :   interval(10),
     520           3 :   delta_t(1.0),
     521           6 :   thermostat("NVE"),
     522           3 :   kbt(-1.0),
     523           3 :   TAMD(0),
     524           3 :   alpha(0.0),
     525           3 :   gamma(0.0),
     526           3 :   kappa(getNumberOfArguments(),0.0),
     527           3 :   fict_max(getNumberOfArguments(),0.0),
     528           3 :   fict_min(getNumberOfArguments(),0.0),
     529           3 :   fict (getNumberOfArguments(),0.0),
     530           3 :   vfict(getNumberOfArguments(),0.0),
     531           3 :   mfict(getNumberOfArguments(),10000.0),
     532           3 :   xeta(0.0),
     533           3 :   veta(0.0),
     534           3 :   meta(0.0),
     535           3 :   flog(0.0),
     536           3 :   hlog(0.0),
     537           3 :   phivs(0.0),
     538           3 :   work(0.0),
     539           3 :   weight(0.0),
     540           3 :   ffict(getNumberOfArguments(),0.0),
     541           3 :   fict_ave(getNumberOfArguments(),0.0),
     542           3 :   fictValue(getNumberOfArguments(),NULL),
     543           6 :   vfictValue(getNumberOfArguments(),NULL)
     544             : {
     545           3 :   parse("INTERVAL",interval);
     546           3 :   parse("DELTA_T",delta_t);
     547           3 :   parse("THERMOSTAT",thermostat);
     548           3 :   parse("TEMP",kbt); // read as temperature
     549             : 
     550           3 :   parse("TAMD",TAMD);
     551           3 :   parse("ALPHA",alpha);
     552           3 :   parse("GAMMA",gamma);
     553           3 :   parseVector("KAPPA",kappa);
     554             : 
     555           3 :   parseVector("FICT_MAX",fict_max);
     556           3 :   parseVector("FICT_MIN",fict_min);
     557             : 
     558           3 :   parseVector("FICT",fict);
     559           3 :   parseVector("VFICT",vfict);
     560           3 :   parseVector("MFICT",mfict);
     561             : 
     562           3 :   parse("XETA",xeta);
     563           3 :   parse("VETA",veta);
     564           3 :   parse("META",meta);
     565             : 
     566           3 :   parse("FLOG",flog);
     567             : 
     568             :   // read initial value of work for each replica of LogPD
     569           3 :   if( multi_sim_comm.Get_size()>1 ) {
     570           0 :     vector<double> vwork(multi_sim_comm.Get_size(),0.0);
     571           0 :     parseVector("WORK",vwork);
     572             :     // initial work of this replica
     573           0 :     work = vwork[multi_sim_comm.Get_rank()];
     574             :   }
     575             :   else {
     576           3 :     work = 0.0;
     577             :   }
     578             : 
     579           3 :   if( kbt>=0.0 ) {
     580           3 :     kbt *= plumed.getAtoms().getKBoltzmann();
     581             :   }
     582             :   else {
     583           0 :     kbt = plumed.getAtoms().getKbT();
     584             :   }
     585             : 
     586           3 :   if( meta == 0.0 ) {
     587           2 :     const double nkt = getNumberOfArguments()*kbt;
     588           2 :     meta = nkt*100.0*100.0;
     589             :   }
     590             : 
     591           3 :   if(alpha == 0.0 && gamma == 0.0) {
     592           0 :     alpha = 4.0;
     593           0 :     gamma = 1 / alpha;
     594             :   }
     595           3 :   else if(alpha != 0.0 && gamma == 0.0) {
     596           3 :     gamma = 1 / alpha;
     597             :   }
     598           0 :   else if(alpha == 0.0 && gamma != 0.0) {
     599           0 :     alpha = 1 / gamma;
     600             :   }
     601             : 
     602           3 :   checkRead();
     603             : 
     604             :   // output messaages to Plumed's log file
     605           3 :   if( multi_sim_comm.Get_size()>1 ) {
     606           0 :     if( TAMD ) {
     607           0 :       log.printf("TAMD-PD, replica parallel of TAMD, no logarithmic flattening.\n");
     608             :     }
     609             :     else {
     610           0 :       log.printf("LogPD, replica parallel of LogMFD.\n");
     611             :     }
     612           0 :     log.printf("number of replica : %d.\n", multi_sim_comm.Get_size() );
     613             :   }
     614             :   else {
     615           3 :     if( TAMD ) {
     616           0 :       log.printf("TAMD, no logarithmic flattening.\n");
     617             :     }
     618             :     else {
     619           3 :       log.printf("LogMFD, logarithmic flattening.\n");
     620             :     }
     621             :   }
     622             : 
     623           3 :   log.printf("  with harmonic force constant      ");
     624           6 :   for(unsigned i=0; i<kappa.size(); i++) log.printf(" %f",kappa[i]);
     625           3 :   log.printf("\n");
     626             : 
     627           3 :   log.printf("  with interval of cv(ideal) update ");
     628           3 :   log.printf(" %d", interval);
     629           3 :   log.printf("\n");
     630             : 
     631           3 :   log.printf("  with time step of cv(ideal) update ");
     632           3 :   log.printf(" %f", delta_t);
     633           3 :   log.printf("\n");
     634             : 
     635           3 :   if( !TAMD ) {
     636           3 :     log.printf("  with alpha, gamma                 ");
     637           3 :     log.printf(" %f %f", alpha, gamma);
     638           3 :     log.printf("\n");
     639             :   }
     640             : 
     641           3 :   log.printf("  with Thermostat for cv(ideal)     ");
     642           3 :   log.printf(" %s", thermostat.c_str());
     643           3 :   log.printf("\n");
     644             : 
     645           3 :   log.printf("  with initial free energy          ");
     646           3 :   log.printf(" %f", flog);
     647           3 :   log.printf("\n");
     648             : 
     649           3 :   log.printf("  with mass of cv(ideal)");
     650           6 :   for(unsigned i=0; i<mfict.size(); i++) log.printf(" %f", mfict[i]);
     651           3 :   log.printf("\n");
     652             : 
     653           3 :   log.printf("  with initial value of cv(ideal)");
     654           6 :   for(unsigned i=0; i<fict.size(); i++) log.printf(" %f", fict[i]);
     655           3 :   log.printf("\n");
     656             : 
     657           3 :   log.printf("  with initial velocity of cv(ideal)");
     658           6 :   for(unsigned i=0; i<vfict.size(); i++) log.printf(" %f", vfict[i]);
     659           3 :   log.printf("\n");
     660             : 
     661           3 :   log.printf("  with maximum value of cv(ideal)    ");
     662           6 :   for(unsigned i=0; i<fict_max.size(); i++) log.printf(" %f",fict_max[i]);
     663           3 :   log.printf("\n");
     664             : 
     665           3 :   log.printf("  with minimum value of cv(ideal)    ");
     666           6 :   for(unsigned i=0; i<fict_min.size(); i++) log.printf(" %f",fict_min[i]);
     667           3 :   log.printf("\n");
     668             : 
     669           3 :   log.printf("  and kbt                           ");
     670           3 :   log.printf(" %f",kbt);
     671           3 :   log.printf("\n");
     672             : 
     673             :   // setup Value* variables
     674           6 :   for(unsigned i=0; i<getNumberOfArguments(); i++) {
     675           3 :     std::string comp = getPntrToArgument(i)->getName()+"_fict";
     676           3 :     addComponentWithDerivatives(comp);
     677             : 
     678           3 :     if(getPntrToArgument(i)->isPeriodic()) {
     679             :       std::string a,b;
     680           0 :       getPntrToArgument(i)->getDomain(a,b);
     681           0 :       componentIsPeriodic(comp,a,b);
     682             :     }
     683             :     else {
     684           3 :       componentIsNotPeriodic(comp);
     685             :     }
     686           3 :     fictValue[i] = getPntrToComponent(comp);
     687             : 
     688           3 :     comp = getPntrToArgument(i)->getName()+"_vfict";
     689           3 :     addComponent(comp);
     690             : 
     691           3 :     componentIsNotPeriodic(comp);
     692           3 :     vfictValue[i] = getPntrToComponent(comp);
     693             :   }
     694           3 : }
     695             : 
     696             : /**
     697             :    \brief calculate forces for fictitious variables at every MD steps.
     698             :    \details This function calculates initial values of fictitious variables
     699             :    and write header messages to LogMFD log files at the first MFD step,
     700             :    calculates restraining fources comes from difference between the fictitious variable
     701             :    and collective variable at every MD steps.
     702             : */
     703        1500 : void LogMFD::calculate() {
     704        1500 :   if( firsttime ) {
     705           3 :     firsttime = false;
     706             : 
     707             :     // set initial values of fictitious variables if they were not specified.
     708           6 :     for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     709           3 :       if( fict[i] != 0.0 ) continue;
     710             : 
     711             :       // use the collective variables as the initial of the fictitious variable.
     712           0 :       fict[i] = getArgument(i);
     713             : 
     714             :       // average values of fictitious variables by all replica.
     715           0 :       if( multi_sim_comm.Get_size()>1 ) {
     716           0 :         multi_sim_comm.Sum(fict[i]);
     717           0 :         fict[i] /= multi_sim_comm.Get_size();
     718             :       }
     719             :     }
     720             : 
     721             :     // initialize accumulation value to zero
     722           6 :     for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     723           3 :       fict_ave[i] = 0.0;
     724             :     }
     725             : 
     726             :     // calculate invariant for NVE
     727           3 :     if(thermostat == "NVE") {
     728             :       // kinetic energy
     729           1 :       const double ekin = calcEkin();
     730             :       // potential energy
     731           1 :       const double pot = TAMD ? flog : sgn(flog)*gamma * std::log1p( alpha*fabs(flog) );
     732             :       // invariant
     733           1 :       hlog = pot + ekin;
     734             :     }
     735           2 :     else if(thermostat == "NVT") {
     736           1 :       const double nkt = getNumberOfArguments()*kbt;
     737             :       // kinetic energy
     738           1 :       const double ekin = calcEkin();
     739             :       // bath energy
     740           1 :       const double ekin_bath = 0.5*veta*veta*meta + xeta*nkt;
     741             :       // potential energy
     742           1 :       const double pot = TAMD ? flog : sgn(flog)*gamma * std::log1p( alpha*fabs(flog) );
     743             :       // invariant
     744           1 :       hlog = pot + ekin + ekin_bath;
     745             :     }
     746           1 :     else if(thermostat == "VS") {
     747             :       // initial velocities
     748           2 :       for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     749           1 :         vfict[i] = sqrt(kbt/mfict[i]);
     750             :       }
     751             :       // initial VS potential
     752           1 :       phivs = TAMD ? flog : sgn(flog)* gamma*std::log1p( alpha*fabs(flog) );
     753             : 
     754             :       // invariant
     755           1 :       hlog = 0.0;
     756             :     }
     757             : 
     758           3 :     weight = 1.0; // for replica parallel
     759             : 
     760             :     // open LogMFD's log file
     761           3 :     if( multi_sim_comm.Get_rank()==0 && comm.Get_rank()==0 ) {
     762           3 :       FILE *outlog = std::fopen("logmfd.out", "w");
     763             : 
     764             :       // output messages to LogMFD's log file
     765           3 :       if( multi_sim_comm.Get_size()>1 ) {
     766             :         fprintf(outlog, "# LogPD, replica parallel of LogMFD\n");
     767           0 :         fprintf(outlog, "# number of replica : %d\n", multi_sim_comm.Get_size() );
     768             :       }
     769             :       else {
     770             :         fprintf(outlog, "# LogMFD\n");
     771             :       }
     772             : 
     773             :       fprintf(outlog, "# CVs :");
     774           6 :       for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     775             :         fprintf(outlog, " %s",  getPntrToArgument(i)->getName().c_str() );
     776             :       }
     777             :       fprintf(outlog, "\n");
     778             : 
     779             :       fprintf(outlog, "# Mass for CV particles :");
     780           6 :       for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     781           3 :         fprintf(outlog, "%18.9f", mfict[i]);
     782             :       }
     783             :       fprintf(outlog, "\n");
     784             : 
     785             :       fprintf(outlog, "# Mass for thermostat   :");
     786           3 :       fprintf(outlog, "%18.9f", meta);
     787             :       fprintf(outlog, "\n");
     788             :       fprintf(outlog, "# 1:iter_md, 2:Flog, 3:2*Ekin/gkb[K], 4:eta, 5:Veta,\n");
     789             : 
     790           6 :       for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     791           3 :         fprintf(outlog, "# %u:%s_fict(t), %u:%s_vfict(t), %u:%s_force(t),\n",
     792             :                 6+i*3, getPntrToArgument(i)->getName().c_str(),
     793             :                 7+i*3, getPntrToArgument(i)->getName().c_str(),
     794           3 :                 8+i*3, getPntrToArgument(i)->getName().c_str() );
     795             :       }
     796             : 
     797           3 :       fclose(outlog);
     798             :     }
     799             : 
     800           3 :     if( comm.Get_rank()==0 ) {
     801             :       // the number of replica is added to file name to distingwish replica.
     802           3 :       FILE *outlog2 = fopen("replica.out", "w");
     803           6 :       fprintf(outlog2, "# Replica No. %d of %d.\n",
     804           3 :               multi_sim_comm.Get_rank(), multi_sim_comm.Get_size() );
     805             : 
     806             :       fprintf(outlog2, "# 1:iter_md, 2:work, 3:weight,\n");
     807           6 :       for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     808           3 :         fprintf(outlog2, "# %u:%s(q)\n",
     809             :                 4+i, getPntrToArgument(i)->getName().c_str() );
     810             :       }
     811           3 :       fclose(outlog2);
     812             :     }
     813             : 
     814             :     // output messages to Plumed's log file
     815             :     //    log.printf("LOGMFD thermostat parameters Xeta Veta Meta");
     816             :     //    log.printf(" %f %f %f", xeta, veta, meta);
     817             :     //    log.printf("\n");
     818             :     //    log.printf("# 1:iter_md, 2:Flog, 3:2*Ekin/gkb[K], 4:eta, 5:Veta,");
     819             :     //    log.printf("# 6:X1(t), 7:V1(t), 8:F1(t), 9:X2(t), 10:V2(t), 11:F2(t), ...");
     820             : 
     821             :   } // firsttime
     822             : 
     823             :   // calculate force for fictitious variable
     824             :   double ene=0.0;
     825        3000 :   for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     826             :     // difference between fictitious variable and collective variable.
     827        1500 :     const double diff = difference(i,fict[i],getArgument(i));
     828             :     // restraining force.
     829        1500 :     const double f = -kappa[i]*diff;
     830             :     setOutputForce(i,f);
     831             : 
     832             :     // restraining energy.
     833        1500 :     ene += 0.5*kappa[i]*diff*diff;
     834             : 
     835             :     // accumulate force, later it will be averaged.
     836        1500 :     ffict[i] += -f;
     837             : 
     838             :     // accumulate varience of collective variable, later it will be averaged.
     839        1500 :     fict_ave[i] += diff;
     840             :   }
     841             : 
     842             :   setBias(ene);
     843        3000 :   for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     844             :     // correct fict so that it is inside [min:max].
     845        1500 :     fict[i] = fictValue[i]->bringBackInPbc(fict[i]);
     846        1500 :     fictValue[i]->set(fict[i]);
     847        1500 :     vfictValue[i]->set(vfict[i]);
     848             :   }
     849        1500 : } // calculate
     850             : 
     851             : /**
     852             :    \brief update fictitious variables.
     853             :    \details This function manages evolution of fictitious variables.
     854             :    This function calculates mean force, updates fictitious variables by one MFD step,
     855             :    bounces back variables, updates free energy, and record logs.
     856             : */
     857        1500 : void LogMFD::update() {
     858        1500 :   if(getStep() == 0 || getStep()%interval != 0 ) return;
     859             : 
     860             :   // calc mean force for fictitious variables
     861          15 :   calcMeanForce();
     862             : 
     863             :   // update fictitious variables
     864          15 :   if(thermostat == "NVE") {
     865           5 :     updateNVE();
     866             :   }
     867          10 :   else if(thermostat == "NVT") {
     868           5 :     updateNVT();
     869             :   }
     870           5 :   else if(thermostat == "VS") {
     871           5 :     updateVS();
     872             :   }
     873             : 
     874             :   // bounce back variables if they are beyond their min and max.
     875          30 :   for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     876          15 :     if(fict[i] > fict_max[i]) {
     877           0 :       fict[i] = fict_max[i] - (fict[i] - fict_max[i]);
     878           0 :       vfict[i] *= -1.0;
     879             :     }
     880          15 :     if(fict[i] < fict_min[i]) {
     881           0 :       fict[i] = fict_min[i] + (fict_min[i] - fict[i]);
     882           0 :       vfict[i] *= -1.0;
     883             :     }
     884             :   }
     885             : 
     886             :   // update free energy
     887          15 :   flog = calcFlog();
     888             : 
     889             :   // record log for fictitious variables
     890          15 :   if( multi_sim_comm.Get_rank()==0 && comm.Get_rank()==0 ) {
     891          15 :     FILE *outlog = std::fopen("logmfd.out", "a");
     892             : 
     893          15 :     const double ekin = calcEkin();
     894          15 :     const double temp = 2.0*ekin/getNumberOfArguments()/plumed.getAtoms().getKBoltzmann();
     895             : 
     896          15 :     fprintf(outlog, "%*d", 8, (int)getStep()/interval);
     897          15 :     fprintf(outlog, "%17.8f", flog);
     898             :     fprintf(outlog, "%17.8f", temp);
     899          15 :     fprintf(outlog, "%17.8f", xeta);
     900          15 :     fprintf(outlog, "%17.8f", veta);
     901          30 :     for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     902          15 :       fprintf(outlog, "%17.8f", fict[i]);
     903          15 :       fprintf(outlog, "%17.8f", vfict[i]);
     904          15 :       fprintf(outlog, "%17.8f", ffict[i]);
     905             :     }
     906             :     fprintf(outlog," \n");
     907          15 :     fclose(outlog);
     908             :   }
     909             : 
     910             :   // record log for collective variables
     911          15 :   if( comm.Get_rank()==0 ) {
     912             :     // the number of replica is added to file name to distingwish replica.
     913          15 :     FILE *outlog2 = fopen("replica.out", "a");
     914          15 :     fprintf(outlog2, "%*d", 8, (int)getStep()/interval);
     915          15 :     fprintf(outlog2, "%16.6e ", work);
     916          15 :     fprintf(outlog2, "%16.6e ", weight);
     917          30 :     for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     918          15 :       fprintf(outlog2, "%17.8f", fict_ave[i]);
     919             :     }
     920             :     fprintf(outlog2," \n");
     921          15 :     fclose(outlog2);
     922             :   }
     923             : 
     924             :   // reset mean force
     925          30 :   for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     926          15 :     ffict[i] = 0.0;
     927          15 :     fict_ave[i] = 0.0;
     928             :   }
     929             : 
     930             : } // update
     931             : 
     932             : /**
     933             :    \brief update fictitious variables by NVE mechanics.
     934             :    \details This function updates ficitious variables by one NVE-MFD step using mean forces
     935             :    and flattening coefficient and free energy.
     936             :  */
     937           5 : void LogMFD::updateNVE() {
     938           5 :   const double dt = delta_t;
     939             : 
     940             :   // get latest free energy and flattening coefficient
     941           5 :   flog = calcFlog();
     942           5 :   const double clog = calcClog();
     943             : 
     944             :   // update all ficitious variables by one MFD step
     945          10 :   for( unsigned i=0; i<getNumberOfArguments(); ++i ) {
     946             :     // update velocity (full step)
     947           5 :     vfict[i]+=clog*ffict[i]*dt/mfict[i];
     948             :     // update position (full step)
     949           5 :     fict[i]+=vfict[i]*dt;
     950             :   }
     951           5 : } // updateNVE
     952             : 
     953             : /**
     954             :    \brief update fictitious variables by NVT mechanics.
     955             :    \details This function updates ficitious variables by one NVT-MFD step using mean forces
     956             :    and flattening coefficient and free energy.
     957             :  */
     958           5 : void LogMFD::updateNVT() {
     959           5 :   const double dt = delta_t;
     960           5 :   const double nkt = getNumberOfArguments()*kbt;
     961             : 
     962             :   // backup vfict
     963           5 :   std::vector<double> vfict_backup(getNumberOfArguments());
     964          10 :   for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     965           5 :     vfict_backup[i] = vfict[i];
     966             :   }
     967             : 
     968             :   const int niter=5;
     969          30 :   for(unsigned j=0; j<niter; ++j) {
     970             :     // get latest free energy and flattening coefficient
     971          25 :     flog = calcFlog();
     972          25 :     const double clog = calcClog();
     973             : 
     974             :     // restore vfict from vfict_backup
     975          50 :     for(unsigned l=0; l<getNumberOfArguments(); ++l) {
     976          25 :       vfict[l] = vfict_backup[l];
     977             :     }
     978             : 
     979             :     // evolve vfict from vfict_backup by dt/2
     980          50 :     for(unsigned m=0; m<getNumberOfArguments(); ++m) {
     981          25 :       vfict[m] *= exp(-0.25*dt*veta);
     982          25 :       vfict[m] += 0.5*dt*clog*ffict[m]/mfict[m];
     983          25 :       vfict[m] *= exp(-0.25*dt*veta);
     984             :     }
     985             :   }
     986             : 
     987             :   // get latest free energy and flattening coefficient
     988           5 :   flog = calcFlog();
     989           5 :   const double clog = calcClog();
     990             : 
     991             :   // evolve vfict by dt/2, and evolve fict by dt
     992          10 :   for(unsigned i=0; i<getNumberOfArguments(); ++i) {
     993           5 :     vfict[i] *= exp(-0.25*dt*veta);
     994           5 :     vfict[i] += 0.5*dt*clog*ffict[i]/mfict[i];
     995           5 :     vfict[i] *= exp(-0.25*dt*veta);
     996           5 :     fict[i]  += dt*vfict[i];
     997             :   }
     998             : 
     999             :   // evolve xeta and veta by dt
    1000           5 :   xeta += 0.5*dt*veta;
    1001           5 :   const double ekin = calcEkin();
    1002           5 :   veta += dt*(2.0*ekin-nkt)/meta;
    1003           5 :   xeta += 0.5*dt*veta;
    1004           5 : } // updateNVT
    1005             : 
    1006             : /**
    1007             :    \brief update fictitious variables by VS mechanics.
    1008             :    \details This function updates ficitious variables by one VS-MFD step using mean forces
    1009             :    and flattening coefficient and free energy.
    1010             :  */
    1011           5 : void LogMFD::updateVS() {
    1012           5 :   const double dt = delta_t;
    1013           5 :   const double nkt = getNumberOfArguments()*kbt;
    1014             : 
    1015             :   // get latest free energy and flattening coefficient
    1016           5 :   flog = calcFlog();
    1017           5 :   const double clog = calcClog();
    1018             : 
    1019          10 :   for(unsigned i=0; i<getNumberOfArguments(); ++i) {
    1020             :     // update velocity (full step)
    1021           5 :     vfict[i] += clog*ffict[i]*dt/mfict[i];
    1022             :   }
    1023             : 
    1024           5 :   const double ekin = calcEkin();
    1025           5 :   const double svs = sqrt(nkt/ekin/2);
    1026             : 
    1027          10 :   for(unsigned i=0; i<getNumberOfArguments(); ++i) {
    1028             :     // update position (full step)
    1029           5 :     vfict[i] *= svs;
    1030           5 :     fict[i] += vfict[i]*dt;
    1031             :   }
    1032             : 
    1033             :   // evolve VS potential
    1034           5 :   phivs += nkt*std::log(svs);
    1035           5 : } // updateVS
    1036             : 
    1037             : /**
    1038             :    \brief calculate mean force for fictitious variables.
    1039             :    \details This function calculates mean forces by averaging forces accumulated during one MFD step,
    1040             :    update work variables done by fictitious variables by one MFD step,
    1041             :    calculate weight variable of this replica for LogPD.
    1042             : */
    1043          15 : void LogMFD::calcMeanForce() {
    1044             :   // cale temporal mean force for each CV
    1045          30 :   for(unsigned i=0; i<getNumberOfArguments(); ++i) {
    1046          30 :     ffict[i] /= interval;
    1047             :     // average of diff (getArgument(i)-fict[i])
    1048          15 :     fict_ave[i] /= interval;
    1049             :     // average of getArgument(i)
    1050          15 :     fict_ave[i] += fict[i];
    1051             : 
    1052             :     // correct fict_ave so that it is inside [min:max].
    1053          15 :     fict_ave[i] = fictValue[i]->bringBackInPbc(fict_ave[i]);
    1054             :   }
    1055             : 
    1056             :   // accumulate work, it was initialized as 0.0
    1057          30 :   for(unsigned i=0; i<getNumberOfArguments(); ++i) {
    1058          15 :     work -= ffict[i] * vfict[i] * delta_t; // modified sign
    1059             :   }
    1060             : 
    1061             :   // for replica parallel
    1062          15 :   if( multi_sim_comm.Get_size()>1 ) {
    1063             :     // find the minimum work among all replicas
    1064           0 :     double work_min = work;
    1065           0 :     multi_sim_comm.Min(work_min);
    1066             : 
    1067             :     // weight of this replica.
    1068             :     // here, work is reduced by work_min to avoid all exp(-work/kbt)s disconverge
    1069           0 :     if( kbt == 0.0 ) {
    1070           0 :       weight = work==work_min ? 1.0 : 0.0;
    1071             :     }
    1072             :     else {
    1073           0 :       weight = exp(-(work-work_min)/kbt);
    1074             :     }
    1075             : 
    1076             :     // normalize the weight
    1077           0 :     double sum_weight = weight;
    1078           0 :     multi_sim_comm.Sum(sum_weight);
    1079           0 :     weight /= sum_weight;
    1080             : 
    1081             :     // weighting force of this replica
    1082           0 :     for(unsigned i=0; i<getNumberOfArguments(); ++i) {
    1083           0 :       ffict[i] *= weight;
    1084             :     }
    1085             : 
    1086             :     // mean forces of all replica.
    1087           0 :     for(unsigned i=0; i<getNumberOfArguments(); ++i) {
    1088           0 :       multi_sim_comm.Sum(ffict[i]);
    1089             :     }
    1090             :     // now, mean force is obtained.
    1091             :   }
    1092          15 : } // calcMeanForce
    1093             : 
    1094             : /**
    1095             :    \brief calculate kinetic energy of fictitious variables.
    1096             :    \retval kinetic energy.
    1097             :    \details This function calculates sum of kinetic energy of all fictitious variables.
    1098             :  */
    1099          82 : double LogMFD::calcEkin() {
    1100             :   double ekin=0.0;
    1101         164 :   for(unsigned i=0; i<getNumberOfArguments(); ++i) {
    1102          82 :     ekin += mfict[i]*vfict[i]*vfict[i]*0.5;
    1103             :   }
    1104          82 :   return ekin;
    1105             : } // calcEkin
    1106             : 
    1107             : /**
    1108             :    \brief calculate free energy of fictitious variables.
    1109             :    \retval free energy.
    1110             :    \details This function calculates free energy by using invariant of canonical mechanics.
    1111             :  */
    1112          55 : double LogMFD::calcFlog() {
    1113          55 :   const double nkt = getNumberOfArguments()*kbt;
    1114          55 :   const double ekin = calcEkin();
    1115             :   double pot;
    1116             : 
    1117          55 :   if (thermostat == "NVE") {
    1118          10 :     pot = hlog - ekin;
    1119             :   }
    1120          45 :   else if (thermostat == "NVT") {
    1121          35 :     const double ekin_bath = 0.5*veta*veta*meta+xeta*nkt;
    1122          35 :     pot = hlog - ekin - ekin_bath;
    1123             :   }
    1124          10 :   else if (thermostat == "VS") {
    1125          10 :     pot = phivs;
    1126             :   }
    1127             :   else {
    1128             :     pot = 0.0; // never occurs
    1129             :   }
    1130             : 
    1131         110 :   return TAMD ? pot : sgn(pot)*expm1(fabs(pot)/gamma)/alpha;
    1132             : } // calcFlog
    1133             : 
    1134             : /**
    1135             :    \brief calculate coefficient for flattening.
    1136             :    \retval flattering coefficient.
    1137             :    \details This function returns 1.0 for TAMD, flattening coefficient for LogMFD.
    1138             :  */
    1139          40 : double LogMFD::calcClog() {
    1140          40 :   return TAMD ? 1.0 : alpha*gamma/(alpha*fabs(flog)+1.0);
    1141             : } // calcClog
    1142             : 
    1143             : } // logmfd
    1144             : } // PLMD

Generated by: LCOV version 1.15