2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gromacs/domdec/domdec_struct.h"
46 #include "gromacs/gmxlib/nrnb.h"
47 #include "gromacs/math/functions.h"
48 #include "gromacs/math/invertmatrix.h"
49 #include "gromacs/math/units.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/math/vecdump.h"
52 #include "gromacs/mdlib/expanded.h"
53 #include "gromacs/mdlib/gmx_omp_nthreads.h"
54 #include "gromacs/mdlib/sim_util.h"
55 #include "gromacs/mdlib/update.h"
56 #include "gromacs/mdtypes/commrec.h"
57 #include "gromacs/mdtypes/group.h"
58 #include "gromacs/mdtypes/inputrec.h"
59 #include "gromacs/mdtypes/md_enums.h"
60 #include "gromacs/mdtypes/state.h"
61 #include "gromacs/pbcutil/boxutilities.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/random/gammadistribution.h"
64 #include "gromacs/random/normaldistribution.h"
65 #include "gromacs/random/tabulatednormaldistribution.h"
66 #include "gromacs/random/threefry.h"
67 #include "gromacs/random/uniformrealdistribution.h"
68 #include "gromacs/utility/cstringutil.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/smalloc.h"
72 #define NTROTTERPARTS 3
74 /* Suzuki-Yoshida Constants, for n=3 and n=5, for symplectic integration */
76 /* for n=3, w0 = w2 = 1/(2-2^-(1/3)), w1 = 1-2*w0 */
77 /* for n=5, w0 = w1 = w3 = w4 = 1/(4-4^-(1/3)), w1 = 1-4*w0 */
79 #define MAX_SUZUKI_YOSHIDA_NUM 5
80 #define SUZUKI_YOSHIDA_NUM 5
82 static const double sy_const_1[] = { 1. };
83 static const double sy_const_3[] = { 0.828981543588751, -0.657963087177502, 0.828981543588751 };
84 static const double sy_const_5[] = { 0.2967324292201065, 0.2967324292201065, -0.186929716880426, 0.2967324292201065, 0.2967324292201065 };
86 static const double* sy_const[] = {
96 static const double sy_const[MAX_SUZUKI_YOSHIDA_NUM+1][MAX_SUZUKI_YOSHIDA_NUM+1] = {
100 {0.828981543588751,-0.657963087177502,0.828981543588751},
102 {0.2967324292201065,0.2967324292201065,-0.186929716880426,0.2967324292201065,0.2967324292201065}
105 /* these integration routines are only referenced inside this file */
106 static void NHC_trotter(t_grpopts *opts, int nvar, gmx_ekindata_t *ekind, real dtfull,
107 double xi[], double vxi[], double scalefac[], real *veta, t_extmass *MassQ, gmx_bool bEkinAveVel)
110 /* general routine for both barostat and thermostat nose hoover chains */
113 double Ekin, Efac, reft, kT, nd;
115 t_grp_tcstat *tcstat;
120 int mstepsi, mstepsj;
121 int ns = SUZUKI_YOSHIDA_NUM; /* set the degree of integration in the types/state.h file */
122 int nh = opts->nhchainlength;
125 mstepsi = mstepsj = ns;
127 /* if scalefac is NULL, we are doing the NHC of the barostat */
130 if (scalefac == nullptr)
135 for (i = 0; i < nvar; i++)
138 /* make it easier to iterate by selecting
139 out the sub-array that corresponds to this T group */
145 iQinv = &(MassQ->QPinv[i*nh]);
146 nd = 1.0; /* THIS WILL CHANGE IF NOT ISOTROPIC */
147 reft = std::max<real>(0, opts->ref_t[0]);
148 Ekin = gmx::square(*veta)/MassQ->Winv;
152 iQinv = &(MassQ->Qinv[i*nh]);
153 tcstat = &ekind->tcstat[i];
155 reft = std::max<real>(0, opts->ref_t[i]);
158 Ekin = 2*trace(tcstat->ekinf)*tcstat->ekinscalef_nhc;
162 Ekin = 2*trace(tcstat->ekinh)*tcstat->ekinscaleh_nhc;
167 for (mi = 0; mi < mstepsi; mi++)
169 for (mj = 0; mj < mstepsj; mj++)
171 /* weighting for this step using Suzuki-Yoshida integration - fixed at 5 */
172 dt = sy_const[ns][mj] * dtfull / mstepsi;
174 /* compute the thermal forces */
175 GQ[0] = iQinv[0]*(Ekin - nd*kT);
177 for (j = 0; j < nh-1; j++)
181 /* we actually don't need to update here if we save the
182 state of the GQ, but it's easier to just recompute*/
183 GQ[j+1] = iQinv[j+1]*((gmx::square(ivxi[j])/iQinv[j])-kT);
191 ivxi[nh-1] += 0.25*dt*GQ[nh-1];
192 for (j = nh-1; j > 0; j--)
194 Efac = exp(-0.125*dt*ivxi[j]);
195 ivxi[j-1] = Efac*(ivxi[j-1]*Efac + 0.25*dt*GQ[j-1]);
198 Efac = exp(-0.5*dt*ivxi[0]);
209 /* Issue - if the KE is an average of the last and the current temperatures, then we might not be
210 able to scale the kinetic energy directly with this factor. Might take more bookkeeping -- have to
211 think about this a bit more . . . */
213 GQ[0] = iQinv[0]*(Ekin - nd*kT);
215 /* update thermostat positions */
216 for (j = 0; j < nh; j++)
218 ixi[j] += 0.5*dt*ivxi[j];
221 for (j = 0; j < nh-1; j++)
223 Efac = exp(-0.125*dt*ivxi[j+1]);
224 ivxi[j] = Efac*(ivxi[j]*Efac + 0.25*dt*GQ[j]);
227 GQ[j+1] = iQinv[j+1]*((gmx::square(ivxi[j])/iQinv[j])-kT);
234 ivxi[nh-1] += 0.25*dt*GQ[nh-1];
241 static void boxv_trotter(t_inputrec *ir, real *veta, real dt, tensor box,
242 gmx_ekindata_t *ekind, tensor vir, real pcorr, t_extmass *MassQ)
249 tensor ekinmod, localpres;
251 /* The heat bath is coupled to a separate barostat, the last temperature group. In the
252 2006 Tuckerman et al paper., the order is iL_{T_baro} iL {T_part}
255 if (ir->epct == epctSEMIISOTROPIC)
264 /* eta is in pure units. veta is in units of ps^-1. GW is in
265 units of ps^-2. However, eta has a reference of 1 nm^3, so care must be
266 taken to use only RATIOS of eta in updating the volume. */
268 /* we take the partial pressure tensors, modify the
269 kinetic energy tensor, and recovert to pressure */
271 if (ir->opts.nrdf[0] == 0)
273 gmx_fatal(FARGS, "Barostat is coupled to a T-group with no degrees of freedom\n");
275 /* alpha factor for phase space volume, then multiply by the ekin scaling factor. */
276 alpha = 1.0 + DIM/((double)ir->opts.nrdf[0]);
277 alpha *= ekind->tcstat[0].ekinscalef_nhc;
278 msmul(ekind->ekin, alpha, ekinmod);
279 /* for now, we use Elr = 0, because if you want to get it right, you
280 really should be using PME. Maybe print a warning? */
282 pscal = calc_pres(ir->ePBC, nwall, box, ekinmod, vir, localpres)+pcorr;
285 GW = (vol*(MassQ->Winv/PRESFAC))*(DIM*pscal - trace(ir->ref_p)); /* W is in ps^2 * bar * nm^3 */
291 * This file implements temperature and pressure coupling algorithms:
292 * For now only the Weak coupling and the modified weak coupling.
294 * Furthermore computation of pressure and temperature is done here
298 real calc_pres(int ePBC, int nwall, matrix box, tensor ekin, tensor vir,
304 if (ePBC == epbcNONE || (ePBC == epbcXY && nwall != 2))
310 /* Uitzoeken welke ekin hier van toepassing is, zie Evans & Morris - E.
311 * Wrs. moet de druktensor gecorrigeerd worden voor de netto stroom in
315 fac = PRESFAC*2.0/det(box);
316 for (n = 0; (n < DIM); n++)
318 for (m = 0; (m < DIM); m++)
320 pres[n][m] = (ekin[n][m] - vir[n][m])*fac;
326 pr_rvecs(debug, 0, "PC: pres", pres, DIM);
327 pr_rvecs(debug, 0, "PC: ekin", ekin, DIM);
328 pr_rvecs(debug, 0, "PC: vir ", vir, DIM);
329 pr_rvecs(debug, 0, "PC: box ", box, DIM);
332 return trace(pres)/DIM;
335 real calc_temp(real ekin, real nrdf)
339 return (2.0*ekin)/(nrdf*BOLTZ);
347 /*! \brief Sets 1/mass for Parrinello-Rahman in wInv; NOTE: PRESFAC is not included, so not in GROMACS units! */
348 static void calcParrinelloRahmanInvMass(const t_inputrec *ir, const matrix box,
353 /* TODO: See if we can make the mass independent of the box size */
354 maxBoxLength = std::max(box[XX][XX], box[YY][YY]);
355 maxBoxLength = std::max(maxBoxLength, box[ZZ][ZZ]);
357 for (int d = 0; d < DIM; d++)
359 for (int n = 0; n < DIM; n++)
361 wInv[d][n] = (4*M_PI*M_PI*ir->compress[d][n])/(3*ir->tau_p*ir->tau_p*maxBoxLength);
366 void parrinellorahman_pcoupl(FILE *fplog, gmx_int64_t step,
367 const t_inputrec *ir, real dt, const tensor pres,
368 tensor box, tensor box_rel, tensor boxv,
369 tensor M, matrix mu, gmx_bool bFirstStep)
371 /* This doesn't do any coordinate updating. It just
372 * integrates the box vector equations from the calculated
373 * acceleration due to pressure difference. We also compute
374 * the tensor M which is used in update to couple the particle
375 * coordinates to the box vectors.
377 * In Nose and Klein (Mol.Phys 50 (1983) no 5., p 1055) this is
380 * M_nk = (h') * (h' * h + h' h) * h
382 * with the dots denoting time derivatives and h is the transformation from
383 * the scaled frame to the real frame, i.e. the TRANSPOSE of the box.
384 * This also goes for the pressure and M tensors - they are transposed relative
385 * to ours. Our equation thus becomes:
388 * M_gmx = M_nk' = b * (b * b' + b * b') * b'
390 * where b is the gromacs box matrix.
391 * Our box accelerations are given by
393 * b = vol/W inv(box') * (P-ref_P) (=h')
396 real vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
397 real atot, arel, change;
398 tensor invbox, pdiff, t1, t2;
400 gmx::invertBoxMatrix(box, invbox);
404 /* Note that PRESFAC does not occur here.
405 * The pressure and compressibility always occur as a product,
406 * therefore the pressure unit drops out.
409 calcParrinelloRahmanInvMass(ir, box, winv);
411 m_sub(pres, ir->ref_p, pdiff);
413 if (ir->epct == epctSURFACETENSION)
415 /* Unlike Berendsen coupling it might not be trivial to include a z
416 * pressure correction here? On the other hand we don't scale the
417 * box momentarily, but change accelerations, so it might not be crucial.
419 real xy_pressure = 0.5*(pres[XX][XX]+pres[YY][YY]);
420 for (int d = 0; d < ZZ; d++)
422 pdiff[d][d] = (xy_pressure-(pres[ZZ][ZZ]-ir->ref_p[d][d]/box[d][d]));
426 tmmul(invbox, pdiff, t1);
427 /* Move the off-diagonal elements of the 'force' to one side to ensure
428 * that we obey the box constraints.
430 for (int d = 0; d < DIM; d++)
432 for (int n = 0; n < d; n++)
434 t1[d][n] += t1[n][d];
441 case epctANISOTROPIC:
442 for (int d = 0; d < DIM; d++)
444 for (int n = 0; n <= d; n++)
446 t1[d][n] *= winv[d][n]*vol;
451 /* calculate total volume acceleration */
452 atot = box[XX][XX]*box[YY][YY]*t1[ZZ][ZZ]+
453 box[XX][XX]*t1[YY][YY]*box[ZZ][ZZ]+
454 t1[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
456 /* set all RELATIVE box accelerations equal, and maintain total V
458 for (int d = 0; d < DIM; d++)
460 for (int n = 0; n <= d; n++)
462 t1[d][n] = winv[0][0]*vol*arel*box[d][n];
466 case epctSEMIISOTROPIC:
467 case epctSURFACETENSION:
468 /* Note the correction to pdiff above for surftens. coupling */
470 /* calculate total XY volume acceleration */
471 atot = box[XX][XX]*t1[YY][YY]+t1[XX][XX]*box[YY][YY];
472 arel = atot/(2*box[XX][XX]*box[YY][YY]);
473 /* set RELATIVE XY box accelerations equal, and maintain total V
474 * change speed. Dont change the third box vector accelerations */
475 for (int d = 0; d < ZZ; d++)
477 for (int n = 0; n <= d; n++)
479 t1[d][n] = winv[d][n]*vol*arel*box[d][n];
482 for (int n = 0; n < DIM; n++)
484 t1[ZZ][n] *= winv[ZZ][n]*vol;
488 gmx_fatal(FARGS, "Parrinello-Rahman pressure coupling type %s "
489 "not supported yet\n", EPCOUPLTYPETYPE(ir->epct));
493 for (int d = 0; d < DIM; d++)
495 for (int n = 0; n <= d; n++)
497 boxv[d][n] += dt*t1[d][n];
499 /* We do NOT update the box vectors themselves here, since
500 * we need them for shifting later. It is instead done last
501 * in the update() routine.
504 /* Calculate the change relative to diagonal elements-
505 since it's perfectly ok for the off-diagonal ones to
506 be zero it doesn't make sense to check the change relative
510 change = std::fabs(dt*boxv[d][n]/box[d][d]);
512 if (change > maxchange)
519 if (maxchange > 0.01 && fplog)
523 "\nStep %s Warning: Pressure scaling more than 1%%. "
524 "This may mean your system\n is not yet equilibrated. "
525 "Use of Parrinello-Rahman pressure coupling during\n"
526 "equilibration can lead to simulation instability, "
527 "and is discouraged.\n",
528 gmx_step_str(step, buf));
532 preserve_box_shape(ir, box_rel, boxv);
534 mtmul(boxv, box, t1); /* t1=boxv * b' */
535 mmul(invbox, t1, t2);
536 mtmul(t2, invbox, M);
538 /* Determine the scaling matrix mu for the coordinates */
539 for (int d = 0; d < DIM; d++)
541 for (int n = 0; n <= d; n++)
543 t1[d][n] = box[d][n] + dt*boxv[d][n];
546 preserve_box_shape(ir, box_rel, t1);
547 /* t1 is the box at t+dt, determine mu as the relative change */
548 mmul_ur0(invbox, t1, mu);
551 void berendsen_pcoupl(FILE *fplog, gmx_int64_t step,
552 const t_inputrec *ir, real dt,
553 const tensor pres, const matrix box,
554 const matrix force_vir, const matrix constraint_vir,
555 matrix mu, double *baros_integral)
558 real scalar_pressure, xy_pressure, p_corr_z;
562 * Calculate the scaling matrix mu
566 for (d = 0; d < DIM; d++)
568 scalar_pressure += pres[d][d]/DIM;
571 xy_pressure += pres[d][d]/(DIM-1);
574 /* Pressure is now in bar, everywhere. */
575 #define factor(d, m) (ir->compress[d][m]*dt/ir->tau_p)
577 /* mu has been changed from pow(1+...,1/3) to 1+.../3, since this is
578 * necessary for triclinic scaling
584 for (d = 0; d < DIM; d++)
586 mu[d][d] = 1.0 - factor(d, d)*(ir->ref_p[d][d] - scalar_pressure) /DIM;
589 case epctSEMIISOTROPIC:
590 for (d = 0; d < ZZ; d++)
592 mu[d][d] = 1.0 - factor(d, d)*(ir->ref_p[d][d]-xy_pressure)/DIM;
595 1.0 - factor(ZZ, ZZ)*(ir->ref_p[ZZ][ZZ] - pres[ZZ][ZZ])/DIM;
597 case epctANISOTROPIC:
598 for (d = 0; d < DIM; d++)
600 for (n = 0; n < DIM; n++)
602 mu[d][n] = (d == n ? 1.0 : 0.0)
603 -factor(d, n)*(ir->ref_p[d][n] - pres[d][n])/DIM;
607 case epctSURFACETENSION:
608 /* ir->ref_p[0/1] is the reference surface-tension times *
609 * the number of surfaces */
610 if (ir->compress[ZZ][ZZ])
612 p_corr_z = dt/ir->tau_p*(ir->ref_p[ZZ][ZZ] - pres[ZZ][ZZ]);
616 /* when the compressibity is zero, set the pressure correction *
617 * in the z-direction to zero to get the correct surface tension */
620 mu[ZZ][ZZ] = 1.0 - ir->compress[ZZ][ZZ]*p_corr_z;
621 for (d = 0; d < DIM-1; d++)
623 mu[d][d] = 1.0 + factor(d, d)*(ir->ref_p[d][d]/(mu[ZZ][ZZ]*box[ZZ][ZZ])
624 - (pres[ZZ][ZZ]+p_corr_z - xy_pressure))/(DIM-1);
628 gmx_fatal(FARGS, "Berendsen pressure coupling type %s not supported yet\n",
629 EPCOUPLTYPETYPE(ir->epct));
631 /* To fullfill the orientation restrictions on triclinic boxes
632 * we will set mu_yx, mu_zx and mu_zy to 0 and correct
633 * the other elements of mu to first order.
635 mu[YY][XX] += mu[XX][YY];
636 mu[ZZ][XX] += mu[XX][ZZ];
637 mu[ZZ][YY] += mu[YY][ZZ];
642 /* Keep track of the work the barostat applies on the system.
643 * Without constraints force_vir tells us how Epot changes when scaling.
644 * With constraints constraint_vir gives us the constraint contribution
645 * to both Epot and Ekin. Although we are not scaling velocities, scaling
646 * the coordinates leads to scaling of distances involved in constraints.
647 * This in turn changes the angular momentum (even if the constrained
648 * distances are corrected at the next step). The kinetic component
649 * of the constraint virial captures the angular momentum change.
651 for (int d = 0; d < DIM; d++)
653 for (int n = 0; n <= d; n++)
655 *baros_integral -= 2*(mu[d][n] - (n == d ? 1 : 0))*(force_vir[d][n] + constraint_vir[d][n]);
661 pr_rvecs(debug, 0, "PC: pres ", pres, 3);
662 pr_rvecs(debug, 0, "PC: mu ", mu, 3);
665 if (mu[XX][XX] < 0.99 || mu[XX][XX] > 1.01 ||
666 mu[YY][YY] < 0.99 || mu[YY][YY] > 1.01 ||
667 mu[ZZ][ZZ] < 0.99 || mu[ZZ][ZZ] > 1.01)
670 sprintf(buf, "\nStep %s Warning: pressure scaling more than 1%%, "
672 gmx_step_str(step, buf2), mu[XX][XX], mu[YY][YY], mu[ZZ][ZZ]);
675 fprintf(fplog, "%s", buf);
677 fprintf(stderr, "%s", buf);
681 void berendsen_pscale(const t_inputrec *ir, const matrix mu,
682 matrix box, matrix box_rel,
683 int start, int nr_atoms,
684 rvec x[], const unsigned short cFREEZE[],
687 ivec *nFreeze = ir->opts.nFreeze;
689 int nthreads gmx_unused;
691 #ifndef __clang_analyzer__
692 // cppcheck-suppress unreadVariable
693 nthreads = gmx_omp_nthreads_get(emntUpdate);
696 /* Scale the positions */
697 #pragma omp parallel for num_threads(nthreads) schedule(static)
698 for (n = start; n < start+nr_atoms; n++)
700 // Trivial OpenMP region that does not throw
703 if (cFREEZE == nullptr)
714 x[n][XX] = mu[XX][XX]*x[n][XX]+mu[YY][XX]*x[n][YY]+mu[ZZ][XX]*x[n][ZZ];
718 x[n][YY] = mu[YY][YY]*x[n][YY]+mu[ZZ][YY]*x[n][ZZ];
722 x[n][ZZ] = mu[ZZ][ZZ]*x[n][ZZ];
725 /* compute final boxlengths */
726 for (d = 0; d < DIM; d++)
728 box[d][XX] = mu[XX][XX]*box[d][XX]+mu[YY][XX]*box[d][YY]+mu[ZZ][XX]*box[d][ZZ];
729 box[d][YY] = mu[YY][YY]*box[d][YY]+mu[ZZ][YY]*box[d][ZZ];
730 box[d][ZZ] = mu[ZZ][ZZ]*box[d][ZZ];
733 preserve_box_shape(ir, box_rel, box);
735 /* (un)shifting should NOT be done after this,
736 * since the box vectors might have changed
738 inc_nrnb(nrnb, eNR_PCOUPL, nr_atoms);
741 void berendsen_tcoupl(const t_inputrec *ir, const gmx_ekindata_t *ekind, real dt,
742 std::vector<double> &therm_integral)
744 const t_grpopts *opts = &ir->opts;
746 for (int i = 0; (i < opts->ngtc); i++)
752 Ek = trace(ekind->tcstat[i].ekinf);
753 T = ekind->tcstat[i].T;
757 Ek = trace(ekind->tcstat[i].ekinh);
758 T = ekind->tcstat[i].Th;
761 if ((opts->tau_t[i] > 0) && (T > 0.0))
763 real reft = std::max<real>(0, opts->ref_t[i]);
764 real lll = std::sqrt(1.0 + (dt/opts->tau_t[i])*(reft/T-1.0));
765 ekind->tcstat[i].lambda = std::max<real>(std::min<real>(lll, 1.25), 0.8);
769 ekind->tcstat[i].lambda = 1.0;
772 /* Keep track of the amount of energy we are adding to the system */
773 therm_integral[i] -= (gmx::square(ekind->tcstat[i].lambda) - 1)*Ek;
777 fprintf(debug, "TC: group %d: T: %g, Lambda: %g\n",
778 i, T, ekind->tcstat[i].lambda);
783 void andersen_tcoupl(t_inputrec *ir, gmx_int64_t step,
784 const t_commrec *cr, const t_mdatoms *md, t_state *state, real rate, const gmx_bool *randomize, const real *boltzfac)
786 const int *gatindex = (DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
789 gmx::ThreeFry2x64<0> rng(ir->andersen_seed, gmx::RandomDomain::Thermostat);
790 gmx::UniformRealDistribution<real> uniformDist;
791 gmx::TabulatedNormalDistribution<real, 14> normalDist;
793 /* randomize the velocities of the selected particles */
795 for (i = 0; i < md->homenr; i++) /* now loop over the list of atoms */
797 int ng = gatindex ? gatindex[i] : i;
800 rng.restart(step, ng);
804 gc = md->cTC[i]; /* assign the atom to a temperature group if there are more than one */
808 if (ir->etc == etcANDERSENMASSIVE)
810 /* Randomize particle always */
815 /* Randomize particle probabilistically */
817 bRandomize = uniformDist(rng) < rate;
824 scal = std::sqrt(boltzfac[gc]*md->invmass[i]);
828 for (d = 0; d < DIM; d++)
830 state->v[i][d] = scal*normalDist(rng);
838 void nosehoover_tcoupl(t_grpopts *opts, gmx_ekindata_t *ekind, real dt,
839 double xi[], double vxi[], t_extmass *MassQ)
844 /* note that this routine does not include Nose-hoover chains yet. Should be easy to add. */
846 for (i = 0; (i < opts->ngtc); i++)
848 reft = std::max<real>(0, opts->ref_t[i]);
850 vxi[i] += dt*MassQ->Qinv[i]*(ekind->tcstat[i].Th - reft);
851 xi[i] += dt*(oldvxi + vxi[i])*0.5;
855 void trotter_update(t_inputrec *ir, gmx_int64_t step, gmx_ekindata_t *ekind,
856 gmx_enerdata_t *enerd, t_state *state,
857 tensor vir, t_mdatoms *md,
858 t_extmass *MassQ, int **trotter_seqlist, int trotter_seqno)
861 int n, i, d, ngtc, gc = 0, t;
862 t_grp_tcstat *tcstat;
864 gmx_int64_t step_eff;
866 double *scalefac, dtc;
868 rvec sumv = {0, 0, 0};
871 if (trotter_seqno <= ettTSEQ2)
873 step_eff = step-1; /* the velocity verlet calls are actually out of order -- the first half step
874 is actually the last half step from the previous step. Thus the first half step
875 actually corresponds to the n-1 step*/
883 bCouple = (ir->nsttcouple == 1 ||
884 do_per_step(step_eff+ir->nsttcouple, ir->nsttcouple));
886 trotter_seq = trotter_seqlist[trotter_seqno];
888 if ((trotter_seq[0] == etrtSKIPALL) || (!bCouple))
892 dtc = ir->nsttcouple*ir->delta_t; /* This is OK for NPT, because nsttcouple == nstpcouple is enforcesd */
893 opts = &(ir->opts); /* just for ease of referencing */
896 snew(scalefac, opts->ngtc);
897 for (i = 0; i < ngtc; i++)
901 /* execute the series of trotter updates specified in the trotterpart array */
903 for (i = 0; i < NTROTTERPARTS; i++)
905 /* allow for doubled intgrators by doubling dt instead of making 2 calls */
906 if ((trotter_seq[i] == etrtBAROV2) || (trotter_seq[i] == etrtBARONHC2) || (trotter_seq[i] == etrtNHC2))
915 switch (trotter_seq[i])
919 boxv_trotter(ir, &(state->veta), dt, state->box, ekind, vir,
920 enerd->term[F_PDISPCORR], MassQ);
924 NHC_trotter(opts, state->nnhpres, ekind, dt, state->nhpres_xi.data(),
925 state->nhpres_vxi.data(), nullptr, &(state->veta), MassQ, FALSE);
929 NHC_trotter(opts, opts->ngtc, ekind, dt, state->nosehoover_xi.data(),
930 state->nosehoover_vxi.data(), scalefac, nullptr, MassQ, (ir->eI == eiVV));
931 /* need to rescale the kinetic energies and velocities here. Could
932 scale the velocities later, but we need them scaled in order to
933 produce the correct outputs, so we'll scale them here. */
935 for (t = 0; t < ngtc; t++)
937 tcstat = &ekind->tcstat[t];
938 tcstat->vscale_nhc = scalefac[t];
939 tcstat->ekinscaleh_nhc *= (scalefac[t]*scalefac[t]);
940 tcstat->ekinscalef_nhc *= (scalefac[t]*scalefac[t]);
942 /* now that we've scaled the groupwise velocities, we can add them up to get the total */
943 /* but do we actually need the total? */
945 /* modify the velocities as well */
946 for (n = 0; n < md->homenr; n++)
948 if (md->cTC) /* does this conditional need to be here? is this always true?*/
952 for (d = 0; d < DIM; d++)
954 state->v[n][d] *= scalefac[gc];
959 for (d = 0; d < DIM; d++)
961 sumv[d] += (state->v[n][d])/md->invmass[n];
970 /* check for conserved momentum -- worth looking at this again eventually, but not working right now.*/
975 extern void init_npt_masses(t_inputrec *ir, t_state *state, t_extmass *MassQ, gmx_bool bInit)
977 int n, i, j, d, ngtc, nh;
979 real reft, kT, ndj, nd;
981 opts = &(ir->opts); /* just for ease of referencing */
982 ngtc = ir->opts.ngtc;
983 nh = state->nhchainlength;
989 snew(MassQ->Qinv, ngtc);
991 for (i = 0; (i < ngtc); i++)
993 if ((opts->tau_t[i] > 0) && (opts->ref_t[i] > 0))
995 MassQ->Qinv[i] = 1.0/(gmx::square(opts->tau_t[i]/M_2PI)*opts->ref_t[i]);
999 MassQ->Qinv[i] = 0.0;
1003 else if (EI_VV(ir->eI))
1005 /* Set pressure variables */
1009 if (state->vol0 == 0)
1011 state->vol0 = det(state->box);
1012 /* because we start by defining a fixed
1013 compressibility, we need the volume at this
1014 compressibility to solve the problem. */
1018 /* units are nm^3 * ns^2 / (nm^3 * bar / kJ/mol) = kJ/mol */
1019 /* Consider evaluating eventually if this the right mass to use. All are correct, some might be more stable */
1020 MassQ->Winv = (PRESFAC*trace(ir->compress)*BOLTZ*opts->ref_t[0])/(DIM*state->vol0*gmx::square(ir->tau_p/M_2PI));
1021 /* An alternate mass definition, from Tuckerman et al. */
1022 /* MassQ->Winv = 1.0/(gmx::square(ir->tau_p/M_2PI)*(opts->nrdf[0]+DIM)*BOLTZ*opts->ref_t[0]); */
1023 for (d = 0; d < DIM; d++)
1025 for (n = 0; n < DIM; n++)
1027 MassQ->Winvm[d][n] = PRESFAC*ir->compress[d][n]/(state->vol0*gmx::square(ir->tau_p/M_2PI));
1028 /* not clear this is correct yet for the anisotropic case. Will need to reevaluate
1029 before using MTTK for anisotropic states.*/
1032 /* Allocate space for thermostat variables */
1035 snew(MassQ->Qinv, ngtc*nh);
1038 /* now, set temperature variables */
1039 for (i = 0; i < ngtc; i++)
1041 if (opts->tau_t[i] > 0 && opts->ref_t[i] > 0 && opts->nrdf[i] > 0)
1043 reft = std::max<real>(0, opts->ref_t[i]);
1046 for (j = 0; j < nh; j++)
1056 MassQ->Qinv[i*nh+j] = 1.0/(gmx::square(opts->tau_t[i]/M_2PI)*ndj*kT);
1061 for (j = 0; j < nh; j++)
1063 MassQ->Qinv[i*nh+j] = 0.0;
1070 int **init_npt_vars(t_inputrec *ir, t_state *state, t_extmass *MassQ, gmx_bool bTrotter)
1072 int i, j, nnhpres, nh;
1074 real bmass, qmass, reft, kT;
1077 opts = &(ir->opts); /* just for ease of referencing */
1078 nnhpres = state->nnhpres;
1079 nh = state->nhchainlength;
1081 if (EI_VV(ir->eI) && (ir->epc == epcMTTK) && (ir->etc != etcNOSEHOOVER))
1083 gmx_fatal(FARGS, "Cannot do MTTK pressure coupling without Nose-Hoover temperature control");
1086 init_npt_masses(ir, state, MassQ, TRUE);
1088 /* first, initialize clear all the trotter calls */
1089 snew(trotter_seq, ettTSEQMAX);
1090 for (i = 0; i < ettTSEQMAX; i++)
1092 snew(trotter_seq[i], NTROTTERPARTS);
1093 for (j = 0; j < NTROTTERPARTS; j++)
1095 trotter_seq[i][j] = etrtNONE;
1097 trotter_seq[i][0] = etrtSKIPALL;
1102 /* no trotter calls, so we never use the values in the array.
1103 * We access them (so we need to define them, but ignore
1109 /* compute the kinetic energy by using the half step velocities or
1110 * the kinetic energies, depending on the order of the trotter calls */
1114 if (inputrecNptTrotter(ir))
1116 /* This is the complicated version - there are 4 possible calls, depending on ordering.
1117 We start with the initial one. */
1118 /* first, a round that estimates veta. */
1119 trotter_seq[0][0] = etrtBAROV;
1121 /* trotter_seq[1] is etrtNHC for 1/2 step velocities - leave zero */
1123 /* The first half trotter update */
1124 trotter_seq[2][0] = etrtBAROV;
1125 trotter_seq[2][1] = etrtNHC;
1126 trotter_seq[2][2] = etrtBARONHC;
1128 /* The second half trotter update */
1129 trotter_seq[3][0] = etrtBARONHC;
1130 trotter_seq[3][1] = etrtNHC;
1131 trotter_seq[3][2] = etrtBAROV;
1133 /* trotter_seq[4] is etrtNHC for second 1/2 step velocities - leave zero */
1136 else if (inputrecNvtTrotter(ir))
1138 /* This is the easy version - there are only two calls, both the same.
1139 Otherwise, even easier -- no calls */
1140 trotter_seq[2][0] = etrtNHC;
1141 trotter_seq[3][0] = etrtNHC;
1143 else if (inputrecNphTrotter(ir))
1145 /* This is the complicated version - there are 4 possible calls, depending on ordering.
1146 We start with the initial one. */
1147 /* first, a round that estimates veta. */
1148 trotter_seq[0][0] = etrtBAROV;
1150 /* trotter_seq[1] is etrtNHC for 1/2 step velocities - leave zero */
1152 /* The first half trotter update */
1153 trotter_seq[2][0] = etrtBAROV;
1154 trotter_seq[2][1] = etrtBARONHC;
1156 /* The second half trotter update */
1157 trotter_seq[3][0] = etrtBARONHC;
1158 trotter_seq[3][1] = etrtBAROV;
1160 /* trotter_seq[4] is etrtNHC for second 1/2 step velocities - leave zero */
1163 else if (ir->eI == eiVVAK)
1165 if (inputrecNptTrotter(ir))
1167 /* This is the complicated version - there are 4 possible calls, depending on ordering.
1168 We start with the initial one. */
1169 /* first, a round that estimates veta. */
1170 trotter_seq[0][0] = etrtBAROV;
1172 /* The first half trotter update, part 1 -- double update, because it commutes */
1173 trotter_seq[1][0] = etrtNHC;
1175 /* The first half trotter update, part 2 */
1176 trotter_seq[2][0] = etrtBAROV;
1177 trotter_seq[2][1] = etrtBARONHC;
1179 /* The second half trotter update, part 1 */
1180 trotter_seq[3][0] = etrtBARONHC;
1181 trotter_seq[3][1] = etrtBAROV;
1183 /* The second half trotter update */
1184 trotter_seq[4][0] = etrtNHC;
1186 else if (inputrecNvtTrotter(ir))
1188 /* This is the easy version - there is only one call, both the same.
1189 Otherwise, even easier -- no calls */
1190 trotter_seq[1][0] = etrtNHC;
1191 trotter_seq[4][0] = etrtNHC;
1193 else if (inputrecNphTrotter(ir))
1195 /* This is the complicated version - there are 4 possible calls, depending on ordering.
1196 We start with the initial one. */
1197 /* first, a round that estimates veta. */
1198 trotter_seq[0][0] = etrtBAROV;
1200 /* The first half trotter update, part 1 -- leave zero */
1201 trotter_seq[1][0] = etrtNHC;
1203 /* The first half trotter update, part 2 */
1204 trotter_seq[2][0] = etrtBAROV;
1205 trotter_seq[2][1] = etrtBARONHC;
1207 /* The second half trotter update, part 1 */
1208 trotter_seq[3][0] = etrtBARONHC;
1209 trotter_seq[3][1] = etrtBAROV;
1211 /* The second half trotter update -- blank for now */
1219 bmass = DIM*DIM; /* recommended mass parameters for isotropic barostat */
1222 snew(MassQ->QPinv, nnhpres*opts->nhchainlength);
1224 /* barostat temperature */
1225 if ((ir->tau_p > 0) && (opts->ref_t[0] > 0))
1227 reft = std::max<real>(0, opts->ref_t[0]);
1229 for (i = 0; i < nnhpres; i++)
1231 for (j = 0; j < nh; j++)
1241 MassQ->QPinv[i*opts->nhchainlength+j] = 1.0/(gmx::square(opts->tau_t[0]/M_2PI)*qmass*kT);
1247 for (i = 0; i < nnhpres; i++)
1249 for (j = 0; j < nh; j++)
1251 MassQ->QPinv[i*nh+j] = 0.0;
1258 static real energyNoseHoover(const t_inputrec *ir, const t_state *state, const t_extmass *MassQ)
1262 int nh = state->nhchainlength;
1264 for (int i = 0; i < ir->opts.ngtc; i++)
1266 const double *ixi = &state->nosehoover_xi[i*nh];
1267 const double *ivxi = &state->nosehoover_vxi[i*nh];
1268 const double *iQinv = &(MassQ->Qinv[i*nh]);
1270 int nd = ir->opts.nrdf[i];
1271 real reft = std::max<real>(ir->opts.ref_t[i], 0);
1272 real kT = BOLTZ * reft;
1276 if (inputrecNvtTrotter(ir))
1278 /* contribution from the thermal momenta of the NH chain */
1279 for (int j = 0; j < nh; j++)
1283 energy += 0.5*gmx::square(ivxi[j])/iQinv[j];
1284 /* contribution from the thermal variable of the NH chain */
1294 energy += ndj*ixi[j]*kT;
1298 else /* Other non Trotter temperature NH control -- no chains yet. */
1300 energy += 0.5*BOLTZ*nd*gmx::square(ivxi[0])/iQinv[0];
1301 energy += nd*ixi[0]*kT;
1309 /* Returns the energy from the barostat thermostat chain */
1310 static real energyPressureMTTK(const t_inputrec *ir, const t_state *state, const t_extmass *MassQ)
1314 int nh = state->nhchainlength;
1316 for (int i = 0; i < state->nnhpres; i++)
1318 /* note -- assumes only one degree of freedom that is thermostatted in barostat */
1319 real reft = std::max<real>(ir->opts.ref_t[0], 0.0); /* using 'System' temperature */
1320 real kT = BOLTZ * reft;
1322 for (int j = 0; j < nh; j++)
1324 double iQinv = MassQ->QPinv[i*nh + j];
1327 energy += 0.5*gmx::square(state->nhpres_vxi[i*nh + j]/iQinv);
1328 /* contribution from the thermal variable of the NH chain */
1329 energy += state->nhpres_xi[i*nh + j]*kT;
1333 fprintf(debug, "P-T-group: %10d Chain %4d ThermV: %15.8f ThermX: %15.8f", i, j, state->nhpres_vxi[i*nh + j], state->nhpres_xi[i*nh + j]);
1341 /* Returns the energy accumulated by the V-rescale or Berendsen thermostat */
1342 static real energyVrescale(const t_inputrec *ir, const t_state *state)
1345 for (int i = 0; i < ir->opts.ngtc; i++)
1347 energy += state->therm_integral[i];
1353 real NPT_energy(const t_inputrec *ir, const t_state *state, const t_extmass *MassQ)
1357 if (ir->epc != epcNO)
1359 /* Compute the contribution of the pressure to the conserved quantity*/
1361 real vol = det(state->box);
1365 case epcPARRINELLORAHMAN:
1367 /* contribution from the pressure momenta */
1369 calcParrinelloRahmanInvMass(ir, state->box, invMass);
1370 for (int d = 0; d < DIM; d++)
1372 for (int n = 0; n <= d; n++)
1374 if (invMass[d][n] > 0)
1376 energyNPT += 0.5*gmx::square(state->boxv[d][n])/(invMass[d][n]*PRESFAC);
1381 /* Contribution from the PV term.
1382 * Not that with non-zero off-diagonal reference pressures,
1383 * i.e. applied shear stresses, there are additional terms.
1384 * We don't support this here, since that requires keeping
1385 * track of unwrapped box diagonal elements. This case is
1386 * excluded in integratorHasConservedEnergyQuantity().
1388 energyNPT += vol*trace(ir->ref_p)/(DIM*PRESFAC);
1392 /* contribution from the pressure momenta */
1393 energyNPT += 0.5*gmx::square(state->veta)/MassQ->Winv;
1395 /* contribution from the PV term */
1396 energyNPT += vol*trace(ir->ref_p)/(DIM*PRESFAC);
1398 if (ir->epc == epcMTTK)
1400 /* contribution from the MTTK chain */
1401 energyNPT += energyPressureMTTK(ir, state, MassQ);
1405 energyNPT += state->baros_integral;
1408 GMX_RELEASE_ASSERT(false, "Conserved energy quantity for pressure coupling is not handled. A case should be added with either the conserved quantity added or nothing added and an exclusion added to integratorHasConservedEnergyQuantity().");
1418 energyNPT += energyVrescale(ir, state);
1421 energyNPT += energyNoseHoover(ir, state, MassQ);
1424 case etcANDERSENMASSIVE:
1425 // Not supported, excluded in integratorHasConservedEnergyQuantity()
1428 GMX_RELEASE_ASSERT(false, "Conserved energy quantity for temperature coupling is not handled. A case should be added with either the conserved quantity added or nothing added and an exclusion added to integratorHasConservedEnergyQuantity().");
1435 static real vrescale_sumnoises(real nn,
1436 gmx::ThreeFry2x64<> *rng,
1437 gmx::NormalDistribution<real> *normalDist)
1440 * Returns the sum of nn independent gaussian noises squared
1441 * (i.e. equivalent to summing the square of the return values
1442 * of nn calls to a normal distribution).
1444 const real ndeg_tol = 0.0001;
1446 gmx::GammaDistribution<real> gammaDist(0.5*nn, 1.0);
1448 if (nn < 2 + ndeg_tol)
1453 nn_int = (int)(nn + 0.5);
1455 if (nn - nn_int < -ndeg_tol || nn - nn_int > ndeg_tol)
1457 gmx_fatal(FARGS, "The v-rescale thermostat was called with a group with #DOF=%f, but for #DOF<3 only integer #DOF are supported", nn + 1);
1461 for (i = 0; i < nn_int; i++)
1463 gauss = (*normalDist)(*rng);
1469 /* Use a gamma distribution for any real nn > 2 */
1470 r = 2.0*gammaDist(*rng);
1476 static real vrescale_resamplekin(real kk, real sigma, real ndeg, real taut,
1477 gmx_int64_t step, gmx_int64_t seed)
1480 * Generates a new value for the kinetic energy,
1481 * according to Bussi et al JCP (2007), Eq. (A7)
1482 * kk: present value of the kinetic energy of the atoms to be thermalized (in arbitrary units)
1483 * sigma: target average value of the kinetic energy (ndeg k_b T/2) (in the same units as kk)
1484 * ndeg: number of degrees of freedom of the atoms to be thermalized
1485 * taut: relaxation time of the thermostat, in units of 'how often this routine is called'
1487 real factor, rr, ekin_new;
1488 gmx::ThreeFry2x64<64> rng(seed, gmx::RandomDomain::Thermostat);
1489 gmx::NormalDistribution<real> normalDist;
1493 factor = exp(-1.0/taut);
1500 rng.restart(step, 0);
1502 rr = normalDist(rng);
1506 (1.0 - factor)*(sigma*(vrescale_sumnoises(ndeg-1, &rng, &normalDist) + rr*rr)/ndeg - kk) +
1507 2.0*rr*std::sqrt(kk*sigma/ndeg*(1.0 - factor)*factor);
1512 void vrescale_tcoupl(t_inputrec *ir, gmx_int64_t step,
1513 gmx_ekindata_t *ekind, real dt,
1514 double therm_integral[])
1518 real Ek, Ek_ref1, Ek_ref, Ek_new;
1522 for (i = 0; (i < opts->ngtc); i++)
1526 Ek = trace(ekind->tcstat[i].ekinf);
1530 Ek = trace(ekind->tcstat[i].ekinh);
1533 if (opts->tau_t[i] >= 0 && opts->nrdf[i] > 0 && Ek > 0)
1535 Ek_ref1 = 0.5*opts->ref_t[i]*BOLTZ;
1536 Ek_ref = Ek_ref1*opts->nrdf[i];
1538 Ek_new = vrescale_resamplekin(Ek, Ek_ref, opts->nrdf[i],
1542 /* Analytically Ek_new>=0, but we check for rounding errors */
1545 ekind->tcstat[i].lambda = 0.0;
1549 ekind->tcstat[i].lambda = std::sqrt(Ek_new/Ek);
1552 therm_integral[i] -= Ek_new - Ek;
1556 fprintf(debug, "TC: group %d: Ekr %g, Ek %g, Ek_new %g, Lambda: %g\n",
1557 i, Ek_ref, Ek, Ek_new, ekind->tcstat[i].lambda);
1562 ekind->tcstat[i].lambda = 1.0;
1567 void rescale_velocities(gmx_ekindata_t *ekind, t_mdatoms *mdatoms,
1568 int start, int end, rvec v[])
1571 t_grp_tcstat *tcstat;
1572 unsigned short *cACC, *cTC;
1577 tcstat = ekind->tcstat;
1582 gstat = ekind->grpstat;
1583 cACC = mdatoms->cACC;
1587 for (n = start; n < end; n++)
1597 /* Only scale the velocity component relative to the COM velocity */
1598 rvec_sub(v[n], gstat[ga].u, vrel);
1599 lg = tcstat[gt].lambda;
1600 for (d = 0; d < DIM; d++)
1602 v[n][d] = gstat[ga].u[d] + lg*vrel[d];
1609 for (n = start; n < end; n++)
1615 lg = tcstat[gt].lambda;
1616 for (d = 0; d < DIM; d++)
1625 /* set target temperatures if we are annealing */
1626 void update_annealing_target_temp(t_inputrec *ir, real t, gmx_update_t *upd)
1628 int i, j, n, npoints;
1629 real pert, thist = 0, x;
1631 for (i = 0; i < ir->opts.ngtc; i++)
1633 npoints = ir->opts.anneal_npoints[i];
1634 switch (ir->opts.annealing[i])
1639 /* calculate time modulo the period */
1640 pert = ir->opts.anneal_time[i][npoints-1];
1641 n = static_cast<int>(t / pert);
1642 thist = t - n*pert; /* modulo time */
1643 /* Make sure rounding didn't get us outside the interval */
1644 if (std::fabs(thist-pert) < GMX_REAL_EPS*100)
1653 gmx_fatal(FARGS, "Death horror in update_annealing_target_temp (i=%d/%d npoints=%d)", i, ir->opts.ngtc, npoints);
1655 /* We are doing annealing for this group if we got here,
1656 * and we have the (relative) time as thist.
1657 * calculate target temp */
1659 while ((j < npoints-1) && (thist > (ir->opts.anneal_time[i][j+1])))
1665 /* Found our position between points j and j+1.
1666 * Interpolate: x is the amount from j+1, (1-x) from point j
1667 * First treat possible jumps in temperature as a special case.
1669 if ((ir->opts.anneal_time[i][j+1]-ir->opts.anneal_time[i][j]) < GMX_REAL_EPS*100)
1671 ir->opts.ref_t[i] = ir->opts.anneal_temp[i][j+1];
1675 x = ((thist-ir->opts.anneal_time[i][j])/
1676 (ir->opts.anneal_time[i][j+1]-ir->opts.anneal_time[i][j]));
1677 ir->opts.ref_t[i] = x*ir->opts.anneal_temp[i][j+1]+(1-x)*ir->opts.anneal_temp[i][j];
1682 ir->opts.ref_t[i] = ir->opts.anneal_temp[i][npoints-1];
1686 update_temperature_constants(upd, ir);