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, 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.
43 #include "gromacs/domdec/domdec_struct.h"
44 #include "gromacs/gmxlib/nrnb.h"
45 #include "gromacs/math/functions.h"
46 #include "gromacs/math/invertmatrix.h"
47 #include "gromacs/math/units.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/math/vecdump.h"
50 #include "gromacs/mdlib/gmx_omp_nthreads.h"
51 #include "gromacs/mdlib/mdrun.h"
52 #include "gromacs/mdlib/sim_util.h"
53 #include "gromacs/mdlib/update.h"
54 #include "gromacs/mdtypes/commrec.h"
55 #include "gromacs/mdtypes/group.h"
56 #include "gromacs/mdtypes/inputrec.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/mdtypes/state.h"
59 #include "gromacs/pbcutil/boxutilities.h"
60 #include "gromacs/pbcutil/pbc.h"
61 #include "gromacs/random/gammadistribution.h"
62 #include "gromacs/random/normaldistribution.h"
63 #include "gromacs/random/tabulatednormaldistribution.h"
64 #include "gromacs/random/threefry.h"
65 #include "gromacs/random/uniformrealdistribution.h"
66 #include "gromacs/trajectory/energy.h"
67 #include "gromacs/utility/cstringutil.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/smalloc.h"
71 #define NTROTTERPARTS 3
73 /* Suzuki-Yoshida Constants, for n=3 and n=5, for symplectic integration */
75 /* for n=3, w0 = w2 = 1/(2-2^-(1/3)), w1 = 1-2*w0 */
76 /* for n=5, w0 = w1 = w3 = w4 = 1/(4-4^-(1/3)), w1 = 1-4*w0 */
78 #define MAX_SUZUKI_YOSHIDA_NUM 5
79 #define SUZUKI_YOSHIDA_NUM 5
81 static const double sy_const_1
[] = { 1. };
82 static const double sy_const_3
[] = { 0.828981543588751, -0.657963087177502, 0.828981543588751 };
83 static const double sy_const_5
[] = { 0.2967324292201065, 0.2967324292201065, -0.186929716880426, 0.2967324292201065, 0.2967324292201065 };
85 static const double* sy_const
[] = {
95 static const double sy_const[MAX_SUZUKI_YOSHIDA_NUM+1][MAX_SUZUKI_YOSHIDA_NUM+1] = {
99 {0.828981543588751,-0.657963087177502,0.828981543588751},
101 {0.2967324292201065,0.2967324292201065,-0.186929716880426,0.2967324292201065,0.2967324292201065}
104 /* these integration routines are only referenced inside this file */
105 static void NHC_trotter(t_grpopts
*opts
, int nvar
, gmx_ekindata_t
*ekind
, real dtfull
,
106 double xi
[], double vxi
[], double scalefac
[], real
*veta
, t_extmass
*MassQ
, gmx_bool bEkinAveVel
)
109 /* general routine for both barostat and thermostat nose hoover chains */
112 double Ekin
, Efac
, reft
, kT
, nd
;
114 t_grp_tcstat
*tcstat
;
119 int mstepsi
, mstepsj
;
120 int ns
= SUZUKI_YOSHIDA_NUM
; /* set the degree of integration in the types/state.h file */
121 int nh
= opts
->nhchainlength
;
124 mstepsi
= mstepsj
= ns
;
126 /* if scalefac is NULL, we are doing the NHC of the barostat */
129 if (scalefac
== nullptr)
134 for (i
= 0; i
< nvar
; i
++)
137 /* make it easier to iterate by selecting
138 out the sub-array that corresponds to this T group */
144 iQinv
= &(MassQ
->QPinv
[i
*nh
]);
145 nd
= 1.0; /* THIS WILL CHANGE IF NOT ISOTROPIC */
146 reft
= std::max
<real
>(0, opts
->ref_t
[0]);
147 Ekin
= gmx::square(*veta
)/MassQ
->Winv
;
151 iQinv
= &(MassQ
->Qinv
[i
*nh
]);
152 tcstat
= &ekind
->tcstat
[i
];
154 reft
= std::max
<real
>(0, opts
->ref_t
[i
]);
157 Ekin
= 2*trace(tcstat
->ekinf
)*tcstat
->ekinscalef_nhc
;
161 Ekin
= 2*trace(tcstat
->ekinh
)*tcstat
->ekinscaleh_nhc
;
166 for (mi
= 0; mi
< mstepsi
; mi
++)
168 for (mj
= 0; mj
< mstepsj
; mj
++)
170 /* weighting for this step using Suzuki-Yoshida integration - fixed at 5 */
171 dt
= sy_const
[ns
][mj
] * dtfull
/ mstepsi
;
173 /* compute the thermal forces */
174 GQ
[0] = iQinv
[0]*(Ekin
- nd
*kT
);
176 for (j
= 0; j
< nh
-1; j
++)
180 /* we actually don't need to update here if we save the
181 state of the GQ, but it's easier to just recompute*/
182 GQ
[j
+1] = iQinv
[j
+1]*((gmx::square(ivxi
[j
])/iQinv
[j
])-kT
);
190 ivxi
[nh
-1] += 0.25*dt
*GQ
[nh
-1];
191 for (j
= nh
-1; j
> 0; j
--)
193 Efac
= exp(-0.125*dt
*ivxi
[j
]);
194 ivxi
[j
-1] = Efac
*(ivxi
[j
-1]*Efac
+ 0.25*dt
*GQ
[j
-1]);
197 Efac
= exp(-0.5*dt
*ivxi
[0]);
208 /* Issue - if the KE is an average of the last and the current temperatures, then we might not be
209 able to scale the kinetic energy directly with this factor. Might take more bookkeeping -- have to
210 think about this a bit more . . . */
212 GQ
[0] = iQinv
[0]*(Ekin
- nd
*kT
);
214 /* update thermostat positions */
215 for (j
= 0; j
< nh
; j
++)
217 ixi
[j
] += 0.5*dt
*ivxi
[j
];
220 for (j
= 0; j
< nh
-1; j
++)
222 Efac
= exp(-0.125*dt
*ivxi
[j
+1]);
223 ivxi
[j
] = Efac
*(ivxi
[j
]*Efac
+ 0.25*dt
*GQ
[j
]);
226 GQ
[j
+1] = iQinv
[j
+1]*((gmx::square(ivxi
[j
])/iQinv
[j
])-kT
);
233 ivxi
[nh
-1] += 0.25*dt
*GQ
[nh
-1];
240 static void boxv_trotter(t_inputrec
*ir
, real
*veta
, real dt
, tensor box
,
241 gmx_ekindata_t
*ekind
, tensor vir
, real pcorr
, t_extmass
*MassQ
)
248 tensor ekinmod
, localpres
;
250 /* The heat bath is coupled to a separate barostat, the last temperature group. In the
251 2006 Tuckerman et al paper., the order is iL_{T_baro} iL {T_part}
254 if (ir
->epct
== epctSEMIISOTROPIC
)
263 /* eta is in pure units. veta is in units of ps^-1. GW is in
264 units of ps^-2. However, eta has a reference of 1 nm^3, so care must be
265 taken to use only RATIOS of eta in updating the volume. */
267 /* we take the partial pressure tensors, modify the
268 kinetic energy tensor, and recovert to pressure */
270 if (ir
->opts
.nrdf
[0] == 0)
272 gmx_fatal(FARGS
, "Barostat is coupled to a T-group with no degrees of freedom\n");
274 /* alpha factor for phase space volume, then multiply by the ekin scaling factor. */
275 alpha
= 1.0 + DIM
/((double)ir
->opts
.nrdf
[0]);
276 alpha
*= ekind
->tcstat
[0].ekinscalef_nhc
;
277 msmul(ekind
->ekin
, alpha
, ekinmod
);
278 /* for now, we use Elr = 0, because if you want to get it right, you
279 really should be using PME. Maybe print a warning? */
281 pscal
= calc_pres(ir
->ePBC
, nwall
, box
, ekinmod
, vir
, localpres
)+pcorr
;
284 GW
= (vol
*(MassQ
->Winv
/PRESFAC
))*(DIM
*pscal
- trace(ir
->ref_p
)); /* W is in ps^2 * bar * nm^3 */
290 * This file implements temperature and pressure coupling algorithms:
291 * For now only the Weak coupling and the modified weak coupling.
293 * Furthermore computation of pressure and temperature is done here
297 real
calc_pres(int ePBC
, int nwall
, matrix box
, tensor ekin
, tensor vir
,
303 if (ePBC
== epbcNONE
|| (ePBC
== epbcXY
&& nwall
!= 2))
309 /* Uitzoeken welke ekin hier van toepassing is, zie Evans & Morris - E.
310 * Wrs. moet de druktensor gecorrigeerd worden voor de netto stroom in
314 fac
= PRESFAC
*2.0/det(box
);
315 for (n
= 0; (n
< DIM
); n
++)
317 for (m
= 0; (m
< DIM
); m
++)
319 pres
[n
][m
] = (ekin
[n
][m
] - vir
[n
][m
])*fac
;
325 pr_rvecs(debug
, 0, "PC: pres", pres
, DIM
);
326 pr_rvecs(debug
, 0, "PC: ekin", ekin
, DIM
);
327 pr_rvecs(debug
, 0, "PC: vir ", vir
, DIM
);
328 pr_rvecs(debug
, 0, "PC: box ", box
, DIM
);
331 return trace(pres
)/DIM
;
334 real
calc_temp(real ekin
, real nrdf
)
338 return (2.0*ekin
)/(nrdf
*BOLTZ
);
346 void parrinellorahman_pcoupl(FILE *fplog
, gmx_int64_t step
,
347 const t_inputrec
*ir
, real dt
, const tensor pres
,
348 tensor box
, tensor box_rel
, tensor boxv
,
349 tensor M
, matrix mu
, gmx_bool bFirstStep
)
351 /* This doesn't do any coordinate updating. It just
352 * integrates the box vector equations from the calculated
353 * acceleration due to pressure difference. We also compute
354 * the tensor M which is used in update to couple the particle
355 * coordinates to the box vectors.
357 * In Nose and Klein (Mol.Phys 50 (1983) no 5., p 1055) this is
360 * M_nk = (h') * (h' * h + h' h) * h
362 * with the dots denoting time derivatives and h is the transformation from
363 * the scaled frame to the real frame, i.e. the TRANSPOSE of the box.
364 * This also goes for the pressure and M tensors - they are transposed relative
365 * to ours. Our equation thus becomes:
368 * M_gmx = M_nk' = b * (b * b' + b * b') * b'
370 * where b is the gromacs box matrix.
371 * Our box accelerations are given by
373 * b = vol/W inv(box') * (P-ref_P) (=h')
378 real vol
= box
[XX
][XX
]*box
[YY
][YY
]*box
[ZZ
][ZZ
];
379 real atot
, arel
, change
, maxchange
, xy_pressure
;
380 tensor invbox
, pdiff
, t1
, t2
;
384 gmx::invertBoxMatrix(box
, invbox
);
388 /* Note that PRESFAC does not occur here.
389 * The pressure and compressibility always occur as a product,
390 * therefore the pressure unit drops out.
392 maxl
= std::max(box
[XX
][XX
], box
[YY
][YY
]);
393 maxl
= std::max(maxl
, box
[ZZ
][ZZ
]);
394 for (d
= 0; d
< DIM
; d
++)
396 for (n
= 0; n
< DIM
; n
++)
399 (4*M_PI
*M_PI
*ir
->compress
[d
][n
])/(3*ir
->tau_p
*ir
->tau_p
*maxl
);
403 m_sub(pres
, ir
->ref_p
, pdiff
);
405 if (ir
->epct
== epctSURFACETENSION
)
407 /* Unlike Berendsen coupling it might not be trivial to include a z
408 * pressure correction here? On the other hand we don't scale the
409 * box momentarily, but change accelerations, so it might not be crucial.
411 xy_pressure
= 0.5*(pres
[XX
][XX
]+pres
[YY
][YY
]);
412 for (d
= 0; d
< ZZ
; d
++)
414 pdiff
[d
][d
] = (xy_pressure
-(pres
[ZZ
][ZZ
]-ir
->ref_p
[d
][d
]/box
[d
][d
]));
418 tmmul(invbox
, pdiff
, t1
);
419 /* Move the off-diagonal elements of the 'force' to one side to ensure
420 * that we obey the box constraints.
422 for (d
= 0; d
< DIM
; d
++)
424 for (n
= 0; n
< d
; n
++)
426 t1
[d
][n
] += t1
[n
][d
];
433 case epctANISOTROPIC
:
434 for (d
= 0; d
< DIM
; d
++)
436 for (n
= 0; n
<= d
; n
++)
438 t1
[d
][n
] *= winv
[d
][n
]*vol
;
443 /* calculate total volume acceleration */
444 atot
= box
[XX
][XX
]*box
[YY
][YY
]*t1
[ZZ
][ZZ
]+
445 box
[XX
][XX
]*t1
[YY
][YY
]*box
[ZZ
][ZZ
]+
446 t1
[XX
][XX
]*box
[YY
][YY
]*box
[ZZ
][ZZ
];
448 /* set all RELATIVE box accelerations equal, and maintain total V
450 for (d
= 0; d
< DIM
; d
++)
452 for (n
= 0; n
<= d
; n
++)
454 t1
[d
][n
] = winv
[0][0]*vol
*arel
*box
[d
][n
];
458 case epctSEMIISOTROPIC
:
459 case epctSURFACETENSION
:
460 /* Note the correction to pdiff above for surftens. coupling */
462 /* calculate total XY volume acceleration */
463 atot
= box
[XX
][XX
]*t1
[YY
][YY
]+t1
[XX
][XX
]*box
[YY
][YY
];
464 arel
= atot
/(2*box
[XX
][XX
]*box
[YY
][YY
]);
465 /* set RELATIVE XY box accelerations equal, and maintain total V
466 * change speed. Dont change the third box vector accelerations */
467 for (d
= 0; d
< ZZ
; d
++)
469 for (n
= 0; n
<= d
; n
++)
471 t1
[d
][n
] = winv
[d
][n
]*vol
*arel
*box
[d
][n
];
474 for (n
= 0; n
< DIM
; n
++)
476 t1
[ZZ
][n
] *= winv
[d
][n
]*vol
;
480 gmx_fatal(FARGS
, "Parrinello-Rahman pressure coupling type %s "
481 "not supported yet\n", EPCOUPLTYPETYPE(ir
->epct
));
486 for (d
= 0; d
< DIM
; d
++)
488 for (n
= 0; n
<= d
; n
++)
490 boxv
[d
][n
] += dt
*t1
[d
][n
];
492 /* We do NOT update the box vectors themselves here, since
493 * we need them for shifting later. It is instead done last
494 * in the update() routine.
497 /* Calculate the change relative to diagonal elements-
498 since it's perfectly ok for the off-diagonal ones to
499 be zero it doesn't make sense to check the change relative
503 change
= fabs(dt
*boxv
[d
][n
]/box
[d
][d
]);
505 if (change
> maxchange
)
512 if (maxchange
> 0.01 && fplog
)
516 "\nStep %s Warning: Pressure scaling more than 1%%. "
517 "This may mean your system\n is not yet equilibrated. "
518 "Use of Parrinello-Rahman pressure coupling during\n"
519 "equilibration can lead to simulation instability, "
520 "and is discouraged.\n",
521 gmx_step_str(step
, buf
));
525 preserve_box_shape(ir
, box_rel
, boxv
);
527 mtmul(boxv
, box
, t1
); /* t1=boxv * b' */
528 mmul(invbox
, t1
, t2
);
529 mtmul(t2
, invbox
, M
);
531 /* Determine the scaling matrix mu for the coordinates */
532 for (d
= 0; d
< DIM
; d
++)
534 for (n
= 0; n
<= d
; n
++)
536 t1
[d
][n
] = box
[d
][n
] + dt
*boxv
[d
][n
];
539 preserve_box_shape(ir
, box_rel
, t1
);
540 /* t1 is the box at t+dt, determine mu as the relative change */
541 mmul_ur0(invbox
, t1
, mu
);
544 void berendsen_pcoupl(FILE *fplog
, gmx_int64_t step
,
545 const t_inputrec
*ir
, real dt
,
546 const tensor pres
, const matrix box
,
550 real scalar_pressure
, xy_pressure
, p_corr_z
;
554 * Calculate the scaling matrix mu
558 for (d
= 0; d
< DIM
; d
++)
560 scalar_pressure
+= pres
[d
][d
]/DIM
;
563 xy_pressure
+= pres
[d
][d
]/(DIM
-1);
566 /* Pressure is now in bar, everywhere. */
567 #define factor(d, m) (ir->compress[d][m]*dt/ir->tau_p)
569 /* mu has been changed from pow(1+...,1/3) to 1+.../3, since this is
570 * necessary for triclinic scaling
576 for (d
= 0; d
< DIM
; d
++)
578 mu
[d
][d
] = 1.0 - factor(d
, d
)*(ir
->ref_p
[d
][d
] - scalar_pressure
) /DIM
;
581 case epctSEMIISOTROPIC
:
582 for (d
= 0; d
< ZZ
; d
++)
584 mu
[d
][d
] = 1.0 - factor(d
, d
)*(ir
->ref_p
[d
][d
]-xy_pressure
)/DIM
;
587 1.0 - factor(ZZ
, ZZ
)*(ir
->ref_p
[ZZ
][ZZ
] - pres
[ZZ
][ZZ
])/DIM
;
589 case epctANISOTROPIC
:
590 for (d
= 0; d
< DIM
; d
++)
592 for (n
= 0; n
< DIM
; n
++)
594 mu
[d
][n
] = (d
== n
? 1.0 : 0.0)
595 -factor(d
, n
)*(ir
->ref_p
[d
][n
] - pres
[d
][n
])/DIM
;
599 case epctSURFACETENSION
:
600 /* ir->ref_p[0/1] is the reference surface-tension times *
601 * the number of surfaces */
602 if (ir
->compress
[ZZ
][ZZ
])
604 p_corr_z
= dt
/ir
->tau_p
*(ir
->ref_p
[ZZ
][ZZ
] - pres
[ZZ
][ZZ
]);
608 /* when the compressibity is zero, set the pressure correction *
609 * in the z-direction to zero to get the correct surface tension */
612 mu
[ZZ
][ZZ
] = 1.0 - ir
->compress
[ZZ
][ZZ
]*p_corr_z
;
613 for (d
= 0; d
< DIM
-1; d
++)
615 mu
[d
][d
] = 1.0 + factor(d
, d
)*(ir
->ref_p
[d
][d
]/(mu
[ZZ
][ZZ
]*box
[ZZ
][ZZ
])
616 - (pres
[ZZ
][ZZ
]+p_corr_z
- xy_pressure
))/(DIM
-1);
620 gmx_fatal(FARGS
, "Berendsen pressure coupling type %s not supported yet\n",
621 EPCOUPLTYPETYPE(ir
->epct
));
624 /* To fullfill the orientation restrictions on triclinic boxes
625 * we will set mu_yx, mu_zx and mu_zy to 0 and correct
626 * the other elements of mu to first order.
628 mu
[YY
][XX
] += mu
[XX
][YY
];
629 mu
[ZZ
][XX
] += mu
[XX
][ZZ
];
630 mu
[ZZ
][YY
] += mu
[YY
][ZZ
];
637 pr_rvecs(debug
, 0, "PC: pres ", pres
, 3);
638 pr_rvecs(debug
, 0, "PC: mu ", mu
, 3);
641 if (mu
[XX
][XX
] < 0.99 || mu
[XX
][XX
] > 1.01 ||
642 mu
[YY
][YY
] < 0.99 || mu
[YY
][YY
] > 1.01 ||
643 mu
[ZZ
][ZZ
] < 0.99 || mu
[ZZ
][ZZ
] > 1.01)
646 sprintf(buf
, "\nStep %s Warning: pressure scaling more than 1%%, "
648 gmx_step_str(step
, buf2
), mu
[XX
][XX
], mu
[YY
][YY
], mu
[ZZ
][ZZ
]);
651 fprintf(fplog
, "%s", buf
);
653 fprintf(stderr
, "%s", buf
);
657 void berendsen_pscale(const t_inputrec
*ir
, const matrix mu
,
658 matrix box
, matrix box_rel
,
659 int start
, int nr_atoms
,
660 rvec x
[], const unsigned short cFREEZE
[],
663 ivec
*nFreeze
= ir
->opts
.nFreeze
;
665 int nthreads gmx_unused
;
667 #ifndef __clang_analyzer__
668 // cppcheck-suppress unreadVariable
669 nthreads
= gmx_omp_nthreads_get(emntUpdate
);
672 /* Scale the positions */
673 #pragma omp parallel for num_threads(nthreads) schedule(static)
674 for (n
= start
; n
< start
+nr_atoms
; n
++)
676 // Trivial OpenMP region that does not throw
679 if (cFREEZE
== nullptr)
690 x
[n
][XX
] = mu
[XX
][XX
]*x
[n
][XX
]+mu
[YY
][XX
]*x
[n
][YY
]+mu
[ZZ
][XX
]*x
[n
][ZZ
];
694 x
[n
][YY
] = mu
[YY
][YY
]*x
[n
][YY
]+mu
[ZZ
][YY
]*x
[n
][ZZ
];
698 x
[n
][ZZ
] = mu
[ZZ
][ZZ
]*x
[n
][ZZ
];
701 /* compute final boxlengths */
702 for (d
= 0; d
< DIM
; d
++)
704 box
[d
][XX
] = mu
[XX
][XX
]*box
[d
][XX
]+mu
[YY
][XX
]*box
[d
][YY
]+mu
[ZZ
][XX
]*box
[d
][ZZ
];
705 box
[d
][YY
] = mu
[YY
][YY
]*box
[d
][YY
]+mu
[ZZ
][YY
]*box
[d
][ZZ
];
706 box
[d
][ZZ
] = mu
[ZZ
][ZZ
]*box
[d
][ZZ
];
709 preserve_box_shape(ir
, box_rel
, box
);
711 /* (un)shifting should NOT be done after this,
712 * since the box vectors might have changed
714 inc_nrnb(nrnb
, eNR_PCOUPL
, nr_atoms
);
717 void berendsen_tcoupl(t_inputrec
*ir
, gmx_ekindata_t
*ekind
, real dt
)
721 real T
, reft
= 0, lll
;
725 for (i
= 0; (i
< opts
->ngtc
); i
++)
729 T
= ekind
->tcstat
[i
].T
;
733 T
= ekind
->tcstat
[i
].Th
;
736 if ((opts
->tau_t
[i
] > 0) && (T
> 0.0))
738 reft
= std::max
<real
>(0, opts
->ref_t
[i
]);
739 lll
= std::sqrt(1.0 + (dt
/opts
->tau_t
[i
])*(reft
/T
-1.0));
740 ekind
->tcstat
[i
].lambda
= std::max
<real
>(std::min
<real
>(lll
, 1.25), 0.8);
744 ekind
->tcstat
[i
].lambda
= 1.0;
749 fprintf(debug
, "TC: group %d: T: %g, Lambda: %g\n",
750 i
, T
, ekind
->tcstat
[i
].lambda
);
755 void andersen_tcoupl(t_inputrec
*ir
, gmx_int64_t step
,
756 const t_commrec
*cr
, const t_mdatoms
*md
, t_state
*state
, real rate
, const gmx_bool
*randomize
, const real
*boltzfac
)
758 const int *gatindex
= (DOMAINDECOMP(cr
) ? cr
->dd
->gatindex
: nullptr);
761 gmx::ThreeFry2x64
<0> rng(ir
->andersen_seed
, gmx::RandomDomain::Thermostat
);
762 gmx::UniformRealDistribution
<real
> uniformDist
;
763 gmx::TabulatedNormalDistribution
<real
, 14> normalDist
;
765 /* randomize the velocities of the selected particles */
767 for (i
= 0; i
< md
->homenr
; i
++) /* now loop over the list of atoms */
769 int ng
= gatindex
? gatindex
[i
] : i
;
772 rng
.restart(step
, ng
);
776 gc
= md
->cTC
[i
]; /* assign the atom to a temperature group if there are more than one */
780 if (ir
->etc
== etcANDERSENMASSIVE
)
782 /* Randomize particle always */
787 /* Randomize particle probabilistically */
789 bRandomize
= uniformDist(rng
) < rate
;
796 scal
= std::sqrt(boltzfac
[gc
]*md
->invmass
[i
]);
800 for (d
= 0; d
< DIM
; d
++)
802 state
->v
[i
][d
] = scal
*normalDist(rng
);
810 void nosehoover_tcoupl(t_grpopts
*opts
, gmx_ekindata_t
*ekind
, real dt
,
811 double xi
[], double vxi
[], t_extmass
*MassQ
)
816 /* note that this routine does not include Nose-hoover chains yet. Should be easy to add. */
818 for (i
= 0; (i
< opts
->ngtc
); i
++)
820 reft
= std::max
<real
>(0, opts
->ref_t
[i
]);
822 vxi
[i
] += dt
*MassQ
->Qinv
[i
]*(ekind
->tcstat
[i
].Th
- reft
);
823 xi
[i
] += dt
*(oldvxi
+ vxi
[i
])*0.5;
827 void trotter_update(t_inputrec
*ir
, gmx_int64_t step
, gmx_ekindata_t
*ekind
,
828 gmx_enerdata_t
*enerd
, t_state
*state
,
829 tensor vir
, t_mdatoms
*md
,
830 t_extmass
*MassQ
, int **trotter_seqlist
, int trotter_seqno
)
833 int n
, i
, d
, ngtc
, gc
= 0, t
;
834 t_grp_tcstat
*tcstat
;
836 gmx_int64_t step_eff
;
838 double *scalefac
, dtc
;
840 rvec sumv
= {0, 0, 0};
843 if (trotter_seqno
<= ettTSEQ2
)
845 step_eff
= step
-1; /* the velocity verlet calls are actually out of order -- the first half step
846 is actually the last half step from the previous step. Thus the first half step
847 actually corresponds to the n-1 step*/
855 bCouple
= (ir
->nsttcouple
== 1 ||
856 do_per_step(step_eff
+ir
->nsttcouple
, ir
->nsttcouple
));
858 trotter_seq
= trotter_seqlist
[trotter_seqno
];
860 if ((trotter_seq
[0] == etrtSKIPALL
) || (!bCouple
))
864 dtc
= ir
->nsttcouple
*ir
->delta_t
; /* This is OK for NPT, because nsttcouple == nstpcouple is enforcesd */
865 opts
= &(ir
->opts
); /* just for ease of referencing */
868 snew(scalefac
, opts
->ngtc
);
869 for (i
= 0; i
< ngtc
; i
++)
873 /* execute the series of trotter updates specified in the trotterpart array */
875 for (i
= 0; i
< NTROTTERPARTS
; i
++)
877 /* allow for doubled intgrators by doubling dt instead of making 2 calls */
878 if ((trotter_seq
[i
] == etrtBAROV2
) || (trotter_seq
[i
] == etrtBARONHC2
) || (trotter_seq
[i
] == etrtNHC2
))
887 switch (trotter_seq
[i
])
891 boxv_trotter(ir
, &(state
->veta
), dt
, state
->box
, ekind
, vir
,
892 enerd
->term
[F_PDISPCORR
], MassQ
);
896 NHC_trotter(opts
, state
->nnhpres
, ekind
, dt
, state
->nhpres_xi
.data(),
897 state
->nhpres_vxi
.data(), nullptr, &(state
->veta
), MassQ
, FALSE
);
901 NHC_trotter(opts
, opts
->ngtc
, ekind
, dt
, state
->nosehoover_xi
.data(),
902 state
->nosehoover_vxi
.data(), scalefac
, nullptr, MassQ
, (ir
->eI
== eiVV
));
903 /* need to rescale the kinetic energies and velocities here. Could
904 scale the velocities later, but we need them scaled in order to
905 produce the correct outputs, so we'll scale them here. */
907 for (t
= 0; t
< ngtc
; t
++)
909 tcstat
= &ekind
->tcstat
[t
];
910 tcstat
->vscale_nhc
= scalefac
[t
];
911 tcstat
->ekinscaleh_nhc
*= (scalefac
[t
]*scalefac
[t
]);
912 tcstat
->ekinscalef_nhc
*= (scalefac
[t
]*scalefac
[t
]);
914 /* now that we've scaled the groupwise velocities, we can add them up to get the total */
915 /* but do we actually need the total? */
917 /* modify the velocities as well */
918 for (n
= 0; n
< md
->homenr
; n
++)
920 if (md
->cTC
) /* does this conditional need to be here? is this always true?*/
924 for (d
= 0; d
< DIM
; d
++)
926 state
->v
[n
][d
] *= scalefac
[gc
];
931 for (d
= 0; d
< DIM
; d
++)
933 sumv
[d
] += (state
->v
[n
][d
])/md
->invmass
[n
];
942 /* check for conserved momentum -- worth looking at this again eventually, but not working right now.*/
948 for (d
= 0; d
< DIM
; d
++)
950 consk
[d
] = sumv
[d
]*exp((1 + 1.0/opts
->nrdf
[0])*((1.0/DIM
)*log(det(state
->box
)/state
->vol0
)) + state
->nosehoover_xi
[0]);
952 fprintf(debug
, "Conserved kappa: %15.8f %15.8f %15.8f\n", consk
[0], consk
[1], consk
[2]);
960 extern void init_npt_masses(t_inputrec
*ir
, t_state
*state
, t_extmass
*MassQ
, gmx_bool bInit
)
962 int n
, i
, j
, d
, ngtc
, nh
;
964 real reft
, kT
, ndj
, nd
;
966 opts
= &(ir
->opts
); /* just for ease of referencing */
967 ngtc
= ir
->opts
.ngtc
;
968 nh
= state
->nhchainlength
;
974 snew(MassQ
->Qinv
, ngtc
);
976 for (i
= 0; (i
< ngtc
); i
++)
978 if ((opts
->tau_t
[i
] > 0) && (opts
->ref_t
[i
] > 0))
980 MassQ
->Qinv
[i
] = 1.0/(gmx::square(opts
->tau_t
[i
]/M_2PI
)*opts
->ref_t
[i
]);
984 MassQ
->Qinv
[i
] = 0.0;
988 else if (EI_VV(ir
->eI
))
990 /* Set pressure variables */
994 if (state
->vol0
== 0)
996 state
->vol0
= det(state
->box
);
997 /* because we start by defining a fixed
998 compressibility, we need the volume at this
999 compressibility to solve the problem. */
1003 /* units are nm^3 * ns^2 / (nm^3 * bar / kJ/mol) = kJ/mol */
1004 /* Consider evaluating eventually if this the right mass to use. All are correct, some might be more stable */
1005 MassQ
->Winv
= (PRESFAC
*trace(ir
->compress
)*BOLTZ
*opts
->ref_t
[0])/(DIM
*state
->vol0
*gmx::square(ir
->tau_p
/M_2PI
));
1006 /* An alternate mass definition, from Tuckerman et al. */
1007 /* MassQ->Winv = 1.0/(gmx::square(ir->tau_p/M_2PI)*(opts->nrdf[0]+DIM)*BOLTZ*opts->ref_t[0]); */
1008 for (d
= 0; d
< DIM
; d
++)
1010 for (n
= 0; n
< DIM
; n
++)
1012 MassQ
->Winvm
[d
][n
] = PRESFAC
*ir
->compress
[d
][n
]/(state
->vol0
*gmx::square(ir
->tau_p
/M_2PI
));
1013 /* not clear this is correct yet for the anisotropic case. Will need to reevaluate
1014 before using MTTK for anisotropic states.*/
1017 /* Allocate space for thermostat variables */
1020 snew(MassQ
->Qinv
, ngtc
*nh
);
1023 /* now, set temperature variables */
1024 for (i
= 0; i
< ngtc
; i
++)
1026 if (opts
->tau_t
[i
] > 0 && opts
->ref_t
[i
] > 0 && opts
->nrdf
[i
] > 0)
1028 reft
= std::max
<real
>(0, opts
->ref_t
[i
]);
1031 for (j
= 0; j
< nh
; j
++)
1041 MassQ
->Qinv
[i
*nh
+j
] = 1.0/(gmx::square(opts
->tau_t
[i
]/M_2PI
)*ndj
*kT
);
1046 for (j
= 0; j
< nh
; j
++)
1048 MassQ
->Qinv
[i
*nh
+j
] = 0.0;
1055 int **init_npt_vars(t_inputrec
*ir
, t_state
*state
, t_extmass
*MassQ
, gmx_bool bTrotter
)
1057 int i
, j
, nnhpres
, nh
;
1059 real bmass
, qmass
, reft
, kT
;
1062 opts
= &(ir
->opts
); /* just for ease of referencing */
1063 nnhpres
= state
->nnhpres
;
1064 nh
= state
->nhchainlength
;
1066 if (EI_VV(ir
->eI
) && (ir
->epc
== epcMTTK
) && (ir
->etc
!= etcNOSEHOOVER
))
1068 gmx_fatal(FARGS
, "Cannot do MTTK pressure coupling without Nose-Hoover temperature control");
1071 init_npt_masses(ir
, state
, MassQ
, TRUE
);
1073 /* first, initialize clear all the trotter calls */
1074 snew(trotter_seq
, ettTSEQMAX
);
1075 for (i
= 0; i
< ettTSEQMAX
; i
++)
1077 snew(trotter_seq
[i
], NTROTTERPARTS
);
1078 for (j
= 0; j
< NTROTTERPARTS
; j
++)
1080 trotter_seq
[i
][j
] = etrtNONE
;
1082 trotter_seq
[i
][0] = etrtSKIPALL
;
1087 /* no trotter calls, so we never use the values in the array.
1088 * We access them (so we need to define them, but ignore
1094 /* compute the kinetic energy by using the half step velocities or
1095 * the kinetic energies, depending on the order of the trotter calls */
1099 if (inputrecNptTrotter(ir
))
1101 /* This is the complicated version - there are 4 possible calls, depending on ordering.
1102 We start with the initial one. */
1103 /* first, a round that estimates veta. */
1104 trotter_seq
[0][0] = etrtBAROV
;
1106 /* trotter_seq[1] is etrtNHC for 1/2 step velocities - leave zero */
1108 /* The first half trotter update */
1109 trotter_seq
[2][0] = etrtBAROV
;
1110 trotter_seq
[2][1] = etrtNHC
;
1111 trotter_seq
[2][2] = etrtBARONHC
;
1113 /* The second half trotter update */
1114 trotter_seq
[3][0] = etrtBARONHC
;
1115 trotter_seq
[3][1] = etrtNHC
;
1116 trotter_seq
[3][2] = etrtBAROV
;
1118 /* trotter_seq[4] is etrtNHC for second 1/2 step velocities - leave zero */
1121 else if (inputrecNvtTrotter(ir
))
1123 /* This is the easy version - there are only two calls, both the same.
1124 Otherwise, even easier -- no calls */
1125 trotter_seq
[2][0] = etrtNHC
;
1126 trotter_seq
[3][0] = etrtNHC
;
1128 else if (inputrecNphTrotter(ir
))
1130 /* This is the complicated version - there are 4 possible calls, depending on ordering.
1131 We start with the initial one. */
1132 /* first, a round that estimates veta. */
1133 trotter_seq
[0][0] = etrtBAROV
;
1135 /* trotter_seq[1] is etrtNHC for 1/2 step velocities - leave zero */
1137 /* The first half trotter update */
1138 trotter_seq
[2][0] = etrtBAROV
;
1139 trotter_seq
[2][1] = etrtBARONHC
;
1141 /* The second half trotter update */
1142 trotter_seq
[3][0] = etrtBARONHC
;
1143 trotter_seq
[3][1] = etrtBAROV
;
1145 /* trotter_seq[4] is etrtNHC for second 1/2 step velocities - leave zero */
1148 else if (ir
->eI
== eiVVAK
)
1150 if (inputrecNptTrotter(ir
))
1152 /* This is the complicated version - there are 4 possible calls, depending on ordering.
1153 We start with the initial one. */
1154 /* first, a round that estimates veta. */
1155 trotter_seq
[0][0] = etrtBAROV
;
1157 /* The first half trotter update, part 1 -- double update, because it commutes */
1158 trotter_seq
[1][0] = etrtNHC
;
1160 /* The first half trotter update, part 2 */
1161 trotter_seq
[2][0] = etrtBAROV
;
1162 trotter_seq
[2][1] = etrtBARONHC
;
1164 /* The second half trotter update, part 1 */
1165 trotter_seq
[3][0] = etrtBARONHC
;
1166 trotter_seq
[3][1] = etrtBAROV
;
1168 /* The second half trotter update */
1169 trotter_seq
[4][0] = etrtNHC
;
1171 else if (inputrecNvtTrotter(ir
))
1173 /* This is the easy version - there is only one call, both the same.
1174 Otherwise, even easier -- no calls */
1175 trotter_seq
[1][0] = etrtNHC
;
1176 trotter_seq
[4][0] = etrtNHC
;
1178 else if (inputrecNphTrotter(ir
))
1180 /* This is the complicated version - there are 4 possible calls, depending on ordering.
1181 We start with the initial one. */
1182 /* first, a round that estimates veta. */
1183 trotter_seq
[0][0] = etrtBAROV
;
1185 /* The first half trotter update, part 1 -- leave zero */
1186 trotter_seq
[1][0] = etrtNHC
;
1188 /* The first half trotter update, part 2 */
1189 trotter_seq
[2][0] = etrtBAROV
;
1190 trotter_seq
[2][1] = etrtBARONHC
;
1192 /* The second half trotter update, part 1 */
1193 trotter_seq
[3][0] = etrtBARONHC
;
1194 trotter_seq
[3][1] = etrtBAROV
;
1196 /* The second half trotter update -- blank for now */
1204 bmass
= DIM
*DIM
; /* recommended mass parameters for isotropic barostat */
1207 snew(MassQ
->QPinv
, nnhpres
*opts
->nhchainlength
);
1209 /* barostat temperature */
1210 if ((ir
->tau_p
> 0) && (opts
->ref_t
[0] > 0))
1212 reft
= std::max
<real
>(0, opts
->ref_t
[0]);
1214 for (i
= 0; i
< nnhpres
; i
++)
1216 for (j
= 0; j
< nh
; j
++)
1226 MassQ
->QPinv
[i
*opts
->nhchainlength
+j
] = 1.0/(gmx::square(opts
->tau_t
[0]/M_2PI
)*qmass
*kT
);
1232 for (i
= 0; i
< nnhpres
; i
++)
1234 for (j
= 0; j
< nh
; j
++)
1236 MassQ
->QPinv
[i
*nh
+j
] = 0.0;
1243 static real
energyNoseHoover(const t_inputrec
*ir
, const t_state
*state
, const t_extmass
*MassQ
)
1247 int nh
= state
->nhchainlength
;
1249 for (int i
= 0; i
< ir
->opts
.ngtc
; i
++)
1251 const double *ixi
= &state
->nosehoover_xi
[i
*nh
];
1252 const double *ivxi
= &state
->nosehoover_vxi
[i
*nh
];
1253 const double *iQinv
= &(MassQ
->Qinv
[i
*nh
]);
1255 int nd
= ir
->opts
.nrdf
[i
];
1256 real reft
= std::max
<real
>(ir
->opts
.ref_t
[i
], 0);
1257 real kT
= BOLTZ
* reft
;
1261 if (inputrecNvtTrotter(ir
))
1263 /* contribution from the thermal momenta of the NH chain */
1264 for (int j
= 0; j
< nh
; j
++)
1268 energy
+= 0.5*gmx::square(ivxi
[j
])/iQinv
[j
];
1269 /* contribution from the thermal variable of the NH chain */
1279 energy
+= ndj
*ixi
[j
]*kT
;
1283 else /* Other non Trotter temperature NH control -- no chains yet. */
1285 energy
+= 0.5*BOLTZ
*nd
*gmx::square(ivxi
[0])/iQinv
[0];
1286 energy
+= nd
*ixi
[0]*kT
;
1294 /* Returns the energy from the barostat thermostat chain */
1295 static real
energyPressureMTTK(const t_inputrec
*ir
, const t_state
*state
, const t_extmass
*MassQ
)
1299 int nh
= state
->nhchainlength
;
1301 for (int i
= 0; i
< state
->nnhpres
; i
++)
1303 /* note -- assumes only one degree of freedom that is thermostatted in barostat */
1304 real reft
= std::max
<real
>(ir
->opts
.ref_t
[0], 0.0); /* using 'System' temperature */
1305 real kT
= BOLTZ
* reft
;
1307 for (int j
= 0; j
< nh
; j
++)
1309 double iQinv
= MassQ
->QPinv
[i
*nh
+ j
];
1312 energy
+= 0.5*gmx::square(state
->nhpres_vxi
[i
*nh
+ j
]/iQinv
);
1313 /* contribution from the thermal variable of the NH chain */
1314 energy
+= state
->nhpres_xi
[i
*nh
+ j
]*kT
;
1318 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
]);
1326 /* Returns the energy accumulated by the V-rescale or Berendsen thermostat */
1327 static real
energyVrescale(const t_inputrec
*ir
, const t_state
*state
)
1330 for (int i
= 0; i
< ir
->opts
.ngtc
; i
++)
1332 energy
+= state
->therm_integral
[i
];
1338 real
NPT_energy(const t_inputrec
*ir
, const t_state
*state
, const t_extmass
*MassQ
)
1342 if (ir
->epc
!= epcNO
)
1344 /* Compute the contribution of the pressure to the conserved quantity*/
1346 real vol
= det(state
->box
);
1350 case epcPARRINELLORAHMAN
:
1354 /* contribution from the pressure momenta */
1355 energyNPT
+= 0.5*gmx::square(state
->veta
)/MassQ
->Winv
;
1357 /* contribution from the PV term */
1358 energyNPT
+= vol
*trace(ir
->ref_p
)/(DIM
*PRESFAC
);
1360 if (ir
->epc
== epcMTTK
)
1362 /* contribution from the MTTK chain */
1363 energyNPT
+= energyPressureMTTK(ir
, state
, MassQ
);
1367 // Not supported, excluded in integratorHasConservedEnergyQuantity()
1371 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().");
1381 energyNPT
+= energyVrescale(ir
, state
);
1384 energyNPT
+= energyNoseHoover(ir
, state
, MassQ
);
1387 case etcANDERSENMASSIVE
:
1388 // Not supported, excluded in integratorHasConservedEnergyQuantity()
1391 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().");
1398 static real
vrescale_sumnoises(real nn
,
1399 gmx::ThreeFry2x64
<> *rng
,
1400 gmx::NormalDistribution
<real
> *normalDist
)
1403 * Returns the sum of nn independent gaussian noises squared
1404 * (i.e. equivalent to summing the square of the return values
1405 * of nn calls to a normal distribution).
1407 const real ndeg_tol
= 0.0001;
1409 gmx::GammaDistribution
<real
> gammaDist(0.5*nn
, 1.0);
1411 if (nn
< 2 + ndeg_tol
)
1416 nn_int
= (int)(nn
+ 0.5);
1418 if (nn
- nn_int
< -ndeg_tol
|| nn
- nn_int
> ndeg_tol
)
1420 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);
1424 for (i
= 0; i
< nn_int
; i
++)
1426 gauss
= (*normalDist
)(*rng
);
1432 /* Use a gamma distribution for any real nn > 2 */
1433 r
= 2.0*gammaDist(*rng
);
1439 static real
vrescale_resamplekin(real kk
, real sigma
, real ndeg
, real taut
,
1440 gmx_int64_t step
, gmx_int64_t seed
)
1443 * Generates a new value for the kinetic energy,
1444 * according to Bussi et al JCP (2007), Eq. (A7)
1445 * kk: present value of the kinetic energy of the atoms to be thermalized (in arbitrary units)
1446 * sigma: target average value of the kinetic energy (ndeg k_b T/2) (in the same units as kk)
1447 * ndeg: number of degrees of freedom of the atoms to be thermalized
1448 * taut: relaxation time of the thermostat, in units of 'how often this routine is called'
1450 real factor
, rr
, ekin_new
;
1451 gmx::ThreeFry2x64
<64> rng(seed
, gmx::RandomDomain::Thermostat
);
1452 gmx::NormalDistribution
<real
> normalDist
;
1456 factor
= exp(-1.0/taut
);
1463 rng
.restart(step
, 0);
1465 rr
= normalDist(rng
);
1469 (1.0 - factor
)*(sigma
*(vrescale_sumnoises(ndeg
-1, &rng
, &normalDist
) + rr
*rr
)/ndeg
- kk
) +
1470 2.0*rr
*std::sqrt(kk
*sigma
/ndeg
*(1.0 - factor
)*factor
);
1475 void vrescale_tcoupl(t_inputrec
*ir
, gmx_int64_t step
,
1476 gmx_ekindata_t
*ekind
, real dt
,
1477 double therm_integral
[])
1481 real Ek
, Ek_ref1
, Ek_ref
, Ek_new
;
1485 for (i
= 0; (i
< opts
->ngtc
); i
++)
1489 Ek
= trace(ekind
->tcstat
[i
].ekinf
);
1493 Ek
= trace(ekind
->tcstat
[i
].ekinh
);
1496 if (opts
->tau_t
[i
] >= 0 && opts
->nrdf
[i
] > 0 && Ek
> 0)
1498 Ek_ref1
= 0.5*opts
->ref_t
[i
]*BOLTZ
;
1499 Ek_ref
= Ek_ref1
*opts
->nrdf
[i
];
1501 Ek_new
= vrescale_resamplekin(Ek
, Ek_ref
, opts
->nrdf
[i
],
1505 /* Analytically Ek_new>=0, but we check for rounding errors */
1508 ekind
->tcstat
[i
].lambda
= 0.0;
1512 ekind
->tcstat
[i
].lambda
= std::sqrt(Ek_new
/Ek
);
1515 therm_integral
[i
] -= Ek_new
- Ek
;
1519 fprintf(debug
, "TC: group %d: Ekr %g, Ek %g, Ek_new %g, Lambda: %g\n",
1520 i
, Ek_ref
, Ek
, Ek_new
, ekind
->tcstat
[i
].lambda
);
1525 ekind
->tcstat
[i
].lambda
= 1.0;
1530 void rescale_velocities(gmx_ekindata_t
*ekind
, t_mdatoms
*mdatoms
,
1531 int start
, int end
, rvec v
[])
1534 t_grp_tcstat
*tcstat
;
1535 unsigned short *cACC
, *cTC
;
1540 tcstat
= ekind
->tcstat
;
1545 gstat
= ekind
->grpstat
;
1546 cACC
= mdatoms
->cACC
;
1550 for (n
= start
; n
< end
; n
++)
1560 /* Only scale the velocity component relative to the COM velocity */
1561 rvec_sub(v
[n
], gstat
[ga
].u
, vrel
);
1562 lg
= tcstat
[gt
].lambda
;
1563 for (d
= 0; d
< DIM
; d
++)
1565 v
[n
][d
] = gstat
[ga
].u
[d
] + lg
*vrel
[d
];
1572 for (n
= start
; n
< end
; n
++)
1578 lg
= tcstat
[gt
].lambda
;
1579 for (d
= 0; d
< DIM
; d
++)
1588 /* set target temperatures if we are annealing */
1589 void update_annealing_target_temp(t_inputrec
*ir
, real t
, gmx_update_t
*upd
)
1591 int i
, j
, n
, npoints
;
1592 real pert
, thist
= 0, x
;
1594 for (i
= 0; i
< ir
->opts
.ngtc
; i
++)
1596 npoints
= ir
->opts
.anneal_npoints
[i
];
1597 switch (ir
->opts
.annealing
[i
])
1602 /* calculate time modulo the period */
1603 pert
= ir
->opts
.anneal_time
[i
][npoints
-1];
1604 n
= static_cast<int>(t
/ pert
);
1605 thist
= t
- n
*pert
; /* modulo time */
1606 /* Make sure rounding didn't get us outside the interval */
1607 if (fabs(thist
-pert
) < GMX_REAL_EPS
*100)
1616 gmx_fatal(FARGS
, "Death horror in update_annealing_target_temp (i=%d/%d npoints=%d)", i
, ir
->opts
.ngtc
, npoints
);
1618 /* We are doing annealing for this group if we got here,
1619 * and we have the (relative) time as thist.
1620 * calculate target temp */
1622 while ((j
< npoints
-1) && (thist
> (ir
->opts
.anneal_time
[i
][j
+1])))
1628 /* Found our position between points j and j+1.
1629 * Interpolate: x is the amount from j+1, (1-x) from point j
1630 * First treat possible jumps in temperature as a special case.
1632 if ((ir
->opts
.anneal_time
[i
][j
+1]-ir
->opts
.anneal_time
[i
][j
]) < GMX_REAL_EPS
*100)
1634 ir
->opts
.ref_t
[i
] = ir
->opts
.anneal_temp
[i
][j
+1];
1638 x
= ((thist
-ir
->opts
.anneal_time
[i
][j
])/
1639 (ir
->opts
.anneal_time
[i
][j
+1]-ir
->opts
.anneal_time
[i
][j
]));
1640 ir
->opts
.ref_t
[i
] = x
*ir
->opts
.anneal_temp
[i
][j
+1]+(1-x
)*ir
->opts
.anneal_temp
[i
][j
];
1645 ir
->opts
.ref_t
[i
] = ir
->opts
.anneal_temp
[i
][npoints
-1];
1649 update_temperature_constants(upd
, ir
);