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, 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.
46 #include "gromacs/domdec/domdec_struct.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/gmxlib/gmx_omp_nthreads.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/listed-forces/disre.h"
52 #include "gromacs/listed-forces/orires.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/math/vecdump.h"
56 #include "gromacs/mdlib/constr.h"
57 #include "gromacs/mdlib/force.h"
58 #include "gromacs/mdlib/mdrun.h"
59 #include "gromacs/mdlib/sim_util.h"
60 #include "gromacs/mdlib/tgroup.h"
61 #include "gromacs/mdtypes/commrec.h"
62 #include "gromacs/mdtypes/group.h"
63 #include "gromacs/mdtypes/inputrec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/pbcutil/boxutilities.h"
66 #include "gromacs/pbcutil/mshift.h"
67 #include "gromacs/pbcutil/pbc.h"
68 #include "gromacs/pulling/pull.h"
69 #include "gromacs/random/random.h"
70 #include "gromacs/timing/wallcycle.h"
71 #include "gromacs/utility/exceptions.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/futil.h"
74 #include "gromacs/utility/gmxassert.h"
75 #include "gromacs/utility/gmxomp.h"
76 #include "gromacs/utility/smalloc.h"
78 /*For debugging, start at v(-dt/2) for velolcity verlet -- uncomment next line */
79 /*#define STARTFROMDT2*/
103 gmx_sd_sigma_t
*sdsig
;
106 /* andersen temperature control stuff */
107 gmx_bool
*randomize_group
;
111 typedef struct gmx_update
114 /* xprime for constraint algorithms */
118 /* Variables for the deform algorithm */
119 gmx_int64_t deformref_step
;
120 matrix deformref_box
;
124 static void do_update_md(int start
, int nrend
, double dt
,
125 t_grp_tcstat
*tcstat
,
127 gmx_bool bNEMD
, t_grp_acc
*gstat
, rvec accel
[],
130 unsigned short ptype
[], unsigned short cFREEZE
[],
131 unsigned short cACC
[], unsigned short cTC
[],
132 rvec x
[], rvec xprime
[], rvec v
[],
134 gmx_bool bNH
, gmx_bool bPR
)
137 int gf
= 0, ga
= 0, gt
= 0;
139 real vn
, vv
, va
, vb
, vnrel
;
145 /* Update with coupling to extended ensembles, used for
146 * Nose-Hoover and Parrinello-Rahman coupling
147 * Nose-Hoover uses the reversible leap-frog integrator from
148 * Holian et al. Phys Rev E 52(3) : 2338, 1995
150 for (n
= start
; n
< nrend
; n
++)
165 lg
= tcstat
[gt
].lambda
;
170 rvec_sub(v
[n
], gstat
[ga
].u
, vrel
);
172 for (d
= 0; d
< DIM
; d
++)
174 if ((ptype
[n
] != eptVSite
) && (ptype
[n
] != eptShell
) && !nFreeze
[gf
][d
])
176 vnrel
= (lg
*vrel
[d
] + dt
*(imass
*f
[n
][d
] - 0.5*vxi
*vrel
[d
]
177 - iprod(M
[d
], vrel
)))/(1 + 0.5*vxi
*dt
);
178 /* do not scale the mean velocities u */
179 vn
= gstat
[ga
].u
[d
] + accel
[ga
][d
]*dt
+ vnrel
;
181 xprime
[n
][d
] = x
[n
][d
]+vn
*dt
;
186 xprime
[n
][d
] = x
[n
][d
];
191 else if (cFREEZE
!= NULL
||
192 nFreeze
[0][XX
] || nFreeze
[0][YY
] || nFreeze
[0][ZZ
] ||
195 /* Update with Berendsen/v-rescale coupling and freeze or NEMD */
196 for (n
= start
; n
< nrend
; n
++)
198 w_dt
= invmass
[n
]*dt
;
211 lg
= tcstat
[gt
].lambda
;
213 for (d
= 0; d
< DIM
; d
++)
216 if ((ptype
[n
] != eptVSite
) && (ptype
[n
] != eptShell
) && !nFreeze
[gf
][d
])
218 vv
= lg
*vn
+ f
[n
][d
]*w_dt
;
220 /* do not scale the mean velocities u */
222 va
= vv
+ accel
[ga
][d
]*dt
;
223 vb
= va
+ (1.0-lg
)*u
;
225 xprime
[n
][d
] = x
[n
][d
]+vb
*dt
;
230 xprime
[n
][d
] = x
[n
][d
];
237 /* Plain update with Berendsen/v-rescale coupling */
238 for (n
= start
; n
< nrend
; n
++)
240 if ((ptype
[n
] != eptVSite
) && (ptype
[n
] != eptShell
))
242 w_dt
= invmass
[n
]*dt
;
247 lg
= tcstat
[gt
].lambda
;
249 for (d
= 0; d
< DIM
; d
++)
251 vn
= lg
*v
[n
][d
] + f
[n
][d
]*w_dt
;
253 xprime
[n
][d
] = x
[n
][d
] + vn
*dt
;
258 for (d
= 0; d
< DIM
; d
++)
261 xprime
[n
][d
] = x
[n
][d
];
268 static void do_update_vv_vel(int start
, int nrend
, double dt
,
269 rvec accel
[], ivec nFreeze
[], real invmass
[],
270 unsigned short ptype
[], unsigned short cFREEZE
[],
271 unsigned short cACC
[], rvec v
[], rvec f
[],
272 gmx_bool bExtended
, real veta
, real alpha
)
281 g
= 0.25*dt
*veta
*alpha
;
283 mv2
= series_sinhx(g
);
290 for (n
= start
; n
< nrend
; n
++)
292 w_dt
= invmass
[n
]*dt
;
302 for (d
= 0; d
< DIM
; d
++)
304 if ((ptype
[n
] != eptVSite
) && (ptype
[n
] != eptShell
) && !nFreeze
[gf
][d
])
306 v
[n
][d
] = mv1
*(mv1
*v
[n
][d
] + 0.5*(w_dt
*mv2
*f
[n
][d
]))+0.5*accel
[ga
][d
]*dt
;
314 } /* do_update_vv_vel */
316 static void do_update_vv_pos(int start
, int nrend
, double dt
,
318 unsigned short ptype
[], unsigned short cFREEZE
[],
319 rvec x
[], rvec xprime
[], rvec v
[],
320 gmx_bool bExtended
, real veta
)
326 /* Would it make more sense if Parrinello-Rahman was put here? */
331 mr2
= series_sinhx(g
);
339 for (n
= start
; n
< nrend
; n
++)
347 for (d
= 0; d
< DIM
; d
++)
349 if ((ptype
[n
] != eptVSite
) && (ptype
[n
] != eptShell
) && !nFreeze
[gf
][d
])
351 xprime
[n
][d
] = mr1
*(mr1
*x
[n
][d
]+mr2
*dt
*v
[n
][d
]);
355 xprime
[n
][d
] = x
[n
][d
];
359 } /* do_update_vv_pos */
361 static void do_update_visc(int start
, int nrend
, double dt
,
362 t_grp_tcstat
*tcstat
,
365 unsigned short ptype
[], unsigned short cTC
[],
366 rvec x
[], rvec xprime
[], rvec v
[],
367 rvec f
[], matrix M
, matrix box
, real
368 cos_accel
, real vcos
,
369 gmx_bool bNH
, gmx_bool bPR
)
374 real lg
, vxi
= 0, vv
;
379 fac
= 2*M_PI
/(box
[ZZ
][ZZ
]);
383 /* Update with coupling to extended ensembles, used for
384 * Nose-Hoover and Parrinello-Rahman coupling
386 for (n
= start
; n
< nrend
; n
++)
393 lg
= tcstat
[gt
].lambda
;
394 cosz
= cos(fac
*x
[n
][ZZ
]);
396 copy_rvec(v
[n
], vrel
);
404 for (d
= 0; d
< DIM
; d
++)
406 if ((ptype
[n
] != eptVSite
) && (ptype
[n
] != eptShell
))
408 vn
= (lg
*vrel
[d
] + dt
*(imass
*f
[n
][d
] - 0.5*vxi
*vrel
[d
]
409 - iprod(M
[d
], vrel
)))/(1 + 0.5*vxi
*dt
);
412 vn
+= vc
+ dt
*cosz
*cos_accel
;
415 xprime
[n
][d
] = x
[n
][d
]+vn
*dt
;
419 xprime
[n
][d
] = x
[n
][d
];
426 /* Classic version of update, used with berendsen coupling */
427 for (n
= start
; n
< nrend
; n
++)
429 w_dt
= invmass
[n
]*dt
;
434 lg
= tcstat
[gt
].lambda
;
435 cosz
= cos(fac
*x
[n
][ZZ
]);
437 for (d
= 0; d
< DIM
; d
++)
441 if ((ptype
[n
] != eptVSite
) && (ptype
[n
] != eptShell
))
446 /* Do not scale the cosine velocity profile */
447 vv
= vc
+ lg
*(vn
- vc
+ f
[n
][d
]*w_dt
);
448 /* Add the cosine accelaration profile */
449 vv
+= dt
*cosz
*cos_accel
;
453 vv
= lg
*(vn
+ f
[n
][d
]*w_dt
);
456 xprime
[n
][d
] = x
[n
][d
]+vv
*dt
;
461 xprime
[n
][d
] = x
[n
][d
];
468 static gmx_stochd_t
*init_stochd(t_inputrec
*ir
)
477 ngtc
= ir
->opts
.ngtc
;
481 snew(sd
->bd_rf
, ngtc
);
483 else if (EI_SD(ir
->eI
))
486 snew(sd
->sdsig
, ngtc
);
489 for (n
= 0; n
< ngtc
; n
++)
491 if (ir
->opts
.tau_t
[n
] > 0)
493 sdc
[n
].gdt
= ir
->delta_t
/ir
->opts
.tau_t
[n
];
494 sdc
[n
].eph
= exp(sdc
[n
].gdt
/2);
495 sdc
[n
].emh
= exp(-sdc
[n
].gdt
/2);
496 sdc
[n
].em
= exp(-sdc
[n
].gdt
);
500 /* No friction and noise on this group */
506 if (sdc
[n
].gdt
>= 0.05)
508 sdc
[n
].b
= sdc
[n
].gdt
*(sdc
[n
].eph
*sdc
[n
].eph
- 1)
509 - 4*(sdc
[n
].eph
- 1)*(sdc
[n
].eph
- 1);
510 sdc
[n
].c
= sdc
[n
].gdt
- 3 + 4*sdc
[n
].emh
- sdc
[n
].em
;
511 sdc
[n
].d
= 2 - sdc
[n
].eph
- sdc
[n
].emh
;
516 /* Seventh order expansions for small y */
517 sdc
[n
].b
= y
*y
*y
*y
*(1/3.0+y
*(1/3.0+y
*(17/90.0+y
*7/9.0)));
518 sdc
[n
].c
= y
*y
*y
*(2/3.0+y
*(-1/2.0+y
*(7/30.0+y
*(-1/12.0+y
*31/1260.0))));
519 sdc
[n
].d
= y
*y
*(-1+y
*y
*(-1/12.0-y
*y
/360.0));
523 fprintf(debug
, "SD const tc-grp %d: b %g c %g d %g\n",
524 n
, sdc
[n
].b
, sdc
[n
].c
, sdc
[n
].d
);
528 else if (ETC_ANDERSEN(ir
->etc
))
537 snew(sd
->randomize_group
, ngtc
);
538 snew(sd
->boltzfac
, ngtc
);
540 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
541 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
543 for (n
= 0; n
< ngtc
; n
++)
545 reft
= std::max
<real
>(0, opts
->ref_t
[n
]);
546 if ((opts
->tau_t
[n
] > 0) && (reft
> 0)) /* tau_t or ref_t = 0 means that no randomization is done */
548 sd
->randomize_group
[n
] = TRUE
;
549 sd
->boltzfac
[n
] = BOLTZ
*opts
->ref_t
[n
];
553 sd
->randomize_group
[n
] = FALSE
;
560 gmx_update_t
init_update(t_inputrec
*ir
)
566 if (ir
->eI
== eiBD
|| EI_SD(ir
->eI
) || ir
->etc
== etcVRESCALE
|| ETC_ANDERSEN(ir
->etc
))
568 upd
->sd
= init_stochd(ir
);
577 static void do_update_sd1(gmx_stochd_t
*sd
,
578 int start
, int nrend
, double dt
,
579 rvec accel
[], ivec nFreeze
[],
580 real invmass
[], unsigned short ptype
[],
581 unsigned short cFREEZE
[], unsigned short cACC
[],
582 unsigned short cTC
[],
583 rvec x
[], rvec xprime
[], rvec v
[], rvec f
[],
584 int ngtc
, real ref_t
[],
586 gmx_bool bFirstHalfConstr
,
587 gmx_int64_t step
, int seed
, int* gatindex
)
592 int gf
= 0, ga
= 0, gt
= 0;
599 for (n
= 0; n
< ngtc
; n
++)
602 /* The mass is encounted for later, since this differs per atom */
603 sig
[n
].V
= sqrt(kT
*(1 - sdc
[n
].em
*sdc
[n
].em
));
608 for (n
= start
; n
< nrend
; n
++)
611 int ng
= gatindex
? gatindex
[n
] : n
;
613 ism
= sqrt(invmass
[n
]);
627 gmx_rng_cycle_3gaussian_table(step
, ng
, seed
, RND_SEED_UPDATE
, rnd
);
629 for (d
= 0; d
< DIM
; d
++)
631 if ((ptype
[n
] != eptVSite
) && (ptype
[n
] != eptShell
) && !nFreeze
[gf
][d
])
635 sd_V
= ism
*sig
[gt
].V
*rnd
[d
];
636 vn
= v
[n
][d
] + (invmass
[n
]*f
[n
][d
] + accel
[ga
][d
])*dt
;
637 v
[n
][d
] = vn
*sdc
[gt
].em
+ sd_V
;
638 /* Here we include half of the friction+noise
639 * update of v into the integration of x.
641 xprime
[n
][d
] = x
[n
][d
] + 0.5*(vn
+ v
[n
][d
])*dt
;
646 xprime
[n
][d
] = x
[n
][d
];
653 /* We do have constraints */
654 if (bFirstHalfConstr
)
656 /* First update without friction and noise */
659 for (n
= start
; n
< nrend
; n
++)
672 for (d
= 0; d
< DIM
; d
++)
674 if ((ptype
[n
] != eptVSite
) && (ptype
[n
] != eptShell
) && !nFreeze
[gf
][d
])
676 v
[n
][d
] = v
[n
][d
] + (im
*f
[n
][d
] + accel
[ga
][d
])*dt
;
677 xprime
[n
][d
] = x
[n
][d
] + v
[n
][d
]*dt
;
682 xprime
[n
][d
] = x
[n
][d
];
689 /* Update friction and noise only */
690 for (n
= start
; n
< nrend
; n
++)
693 int ng
= gatindex
? gatindex
[n
] : n
;
695 ism
= sqrt(invmass
[n
]);
705 gmx_rng_cycle_3gaussian_table(step
, ng
, seed
, RND_SEED_UPDATE
, rnd
);
707 for (d
= 0; d
< DIM
; d
++)
709 if ((ptype
[n
] != eptVSite
) && (ptype
[n
] != eptShell
) && !nFreeze
[gf
][d
])
713 sd_V
= ism
*sig
[gt
].V
*rnd
[d
];
715 v
[n
][d
] = vn
*sdc
[gt
].em
+ sd_V
;
716 /* Add the friction and noise contribution only */
717 xprime
[n
][d
] = xprime
[n
][d
] + 0.5*(v
[n
][d
] - vn
)*dt
;
725 static void check_sd2_work_data_allocation(gmx_stochd_t
*sd
, int nrend
)
727 if (nrend
> sd
->sd_V_nalloc
)
729 sd
->sd_V_nalloc
= over_alloc_dd(nrend
);
730 srenew(sd
->sd_V
, sd
->sd_V_nalloc
);
734 static void do_update_sd2_Tconsts(gmx_stochd_t
*sd
,
739 /* This is separated from the update below, because it is single threaded */
748 for (gt
= 0; gt
< ngtc
; gt
++)
750 kT
= BOLTZ
*ref_t
[gt
];
751 /* The mass is encounted for later, since this differs per atom */
752 sig
[gt
].V
= sqrt(kT
*(1-sdc
[gt
].em
));
753 sig
[gt
].X
= sqrt(kT
*sqr(tau_t
[gt
])*sdc
[gt
].c
);
754 sig
[gt
].Yv
= sqrt(kT
*sdc
[gt
].b
/sdc
[gt
].c
);
755 sig
[gt
].Yx
= sqrt(kT
*sqr(tau_t
[gt
])*sdc
[gt
].b
/(1-sdc
[gt
].em
));
759 static void do_update_sd2(gmx_stochd_t
*sd
,
761 int start
, int nrend
,
762 rvec accel
[], ivec nFreeze
[],
763 real invmass
[], unsigned short ptype
[],
764 unsigned short cFREEZE
[], unsigned short cACC
[],
765 unsigned short cTC
[],
766 rvec x
[], rvec xprime
[], rvec v
[], rvec f
[],
769 gmx_bool bFirstHalf
, gmx_int64_t step
, int seed
,
774 /* The random part of the velocity update, generated in the first
775 * half of the update, needs to be remembered for the second half.
778 int gf
= 0, ga
= 0, gt
= 0;
779 real vn
= 0, Vmh
, Xmh
;
787 for (n
= start
; n
< nrend
; n
++)
789 real rnd
[6], rndi
[3];
790 ng
= gatindex
? gatindex
[n
] : n
;
791 ism
= sqrt(invmass
[n
]);
805 gmx_rng_cycle_6gaussian_table(step
*2+(bFirstHalf
? 1 : 2), ng
, seed
, RND_SEED_UPDATE
, rnd
);
808 gmx_rng_cycle_3gaussian_table(step
*2, ng
, seed
, RND_SEED_UPDATE
, rndi
);
810 for (d
= 0; d
< DIM
; d
++)
816 if ((ptype
[n
] != eptVSite
) && (ptype
[n
] != eptShell
) && !nFreeze
[gf
][d
])
822 sd_X
[n
][d
] = ism
*sig
[gt
].X
*rndi
[d
];
824 Vmh
= sd_X
[n
][d
]*sdc
[gt
].d
/(tau_t
[gt
]*sdc
[gt
].c
)
825 + ism
*sig
[gt
].Yv
*rnd
[d
*2];
826 sd_V
[n
][d
] = ism
*sig
[gt
].V
*rnd
[d
*2+1];
828 v
[n
][d
] = vn
*sdc
[gt
].em
829 + (invmass
[n
]*f
[n
][d
] + accel
[ga
][d
])*tau_t
[gt
]*(1 - sdc
[gt
].em
)
830 + sd_V
[n
][d
] - sdc
[gt
].em
*Vmh
;
832 xprime
[n
][d
] = x
[n
][d
] + v
[n
][d
]*tau_t
[gt
]*(sdc
[gt
].eph
- sdc
[gt
].emh
);
836 /* Correct the velocities for the constraints.
837 * This operation introduces some inaccuracy,
838 * since the velocity is determined from differences in coordinates.
841 (xprime
[n
][d
] - x
[n
][d
])/(tau_t
[gt
]*(sdc
[gt
].eph
- sdc
[gt
].emh
));
843 Xmh
= sd_V
[n
][d
]*tau_t
[gt
]*sdc
[gt
].d
/(sdc
[gt
].em
-1)
844 + ism
*sig
[gt
].Yx
*rnd
[d
*2];
845 sd_X
[n
][d
] = ism
*sig
[gt
].X
*rnd
[d
*2+1];
847 xprime
[n
][d
] += sd_X
[n
][d
] - Xmh
;
856 xprime
[n
][d
] = x
[n
][d
];
863 static void do_update_bd_Tconsts(double dt
, real friction_coefficient
,
864 int ngtc
, const real ref_t
[],
867 /* This is separated from the update below, because it is single threaded */
870 if (friction_coefficient
!= 0)
872 for (gt
= 0; gt
< ngtc
; gt
++)
874 rf
[gt
] = sqrt(2.0*BOLTZ
*ref_t
[gt
]/(friction_coefficient
*dt
));
879 for (gt
= 0; gt
< ngtc
; gt
++)
881 rf
[gt
] = sqrt(2.0*BOLTZ
*ref_t
[gt
]);
886 static void do_update_bd(int start
, int nrend
, double dt
,
888 real invmass
[], unsigned short ptype
[],
889 unsigned short cFREEZE
[], unsigned short cTC
[],
890 rvec x
[], rvec xprime
[], rvec v
[],
891 rvec f
[], real friction_coefficient
,
892 real
*rf
, gmx_int64_t step
, int seed
,
895 /* note -- these appear to be full step velocities . . . */
901 if (friction_coefficient
!= 0)
903 invfr
= 1.0/friction_coefficient
;
906 for (n
= start
; (n
< nrend
); n
++)
909 int ng
= gatindex
? gatindex
[n
] : n
;
919 gmx_rng_cycle_3gaussian_table(step
, ng
, seed
, RND_SEED_UPDATE
, rnd
);
920 for (d
= 0; (d
< DIM
); d
++)
922 if ((ptype
[n
] != eptVSite
) && (ptype
[n
] != eptShell
) && !nFreeze
[gf
][d
])
924 if (friction_coefficient
!= 0)
926 vn
= invfr
*f
[n
][d
] + rf
[gt
]*rnd
[d
];
930 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
931 vn
= 0.5*invmass
[n
]*f
[n
][d
]*dt
932 + sqrt(0.5*invmass
[n
])*rf
[gt
]*rnd
[d
];
936 xprime
[n
][d
] = x
[n
][d
]+vn
*dt
;
941 xprime
[n
][d
] = x
[n
][d
];
947 static void dump_it_all(FILE gmx_unused
*fp
, const char gmx_unused
*title
,
948 int gmx_unused natoms
, rvec gmx_unused x
[], rvec gmx_unused xp
[],
949 rvec gmx_unused v
[], rvec gmx_unused f
[])
954 fprintf(fp
, "%s\n", title
);
955 pr_rvecs(fp
, 0, "x", x
, natoms
);
956 pr_rvecs(fp
, 0, "xp", xp
, natoms
);
957 pr_rvecs(fp
, 0, "v", v
, natoms
);
958 pr_rvecs(fp
, 0, "f", f
, natoms
);
963 static void calc_ke_part_normal(rvec v
[], t_grpopts
*opts
, t_mdatoms
*md
,
964 gmx_ekindata_t
*ekind
, t_nrnb
*nrnb
, gmx_bool bEkinAveVel
)
967 t_grp_tcstat
*tcstat
= ekind
->tcstat
;
968 t_grp_acc
*grpstat
= ekind
->grpstat
;
971 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
972 an option, but not supported now.
973 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
976 /* group velocities are calculated in update_ekindata and
977 * accumulated in acumulate_groups.
978 * Now the partial global and groups ekin.
980 for (g
= 0; (g
< opts
->ngtc
); g
++)
982 copy_mat(tcstat
[g
].ekinh
, tcstat
[g
].ekinh_old
);
985 clear_mat(tcstat
[g
].ekinf
);
986 tcstat
[g
].ekinscalef_nhc
= 1.0; /* need to clear this -- logic is complicated! */
990 clear_mat(tcstat
[g
].ekinh
);
993 ekind
->dekindl_old
= ekind
->dekindl
;
994 nthread
= gmx_omp_nthreads_get(emntUpdate
);
996 #pragma omp parallel for num_threads(nthread) schedule(static)
997 for (thread
= 0; thread
< nthread
; thread
++)
999 // This OpenMP only loops over arrays and does not call any functions
1000 // or memory allocation. It should not be able to throw, so for now
1001 // we do not need a try/catch wrapper.
1002 int start_t
, end_t
, n
;
1010 start_t
= ((thread
+0)*md
->homenr
)/nthread
;
1011 end_t
= ((thread
+1)*md
->homenr
)/nthread
;
1013 ekin_sum
= ekind
->ekin_work
[thread
];
1014 dekindl_sum
= ekind
->dekindl_work
[thread
];
1016 for (gt
= 0; gt
< opts
->ngtc
; gt
++)
1018 clear_mat(ekin_sum
[gt
]);
1024 for (n
= start_t
; n
< end_t
; n
++)
1034 hm
= 0.5*md
->massT
[n
];
1036 for (d
= 0; (d
< DIM
); d
++)
1038 v_corrt
[d
] = v
[n
][d
] - grpstat
[ga
].u
[d
];
1040 for (d
= 0; (d
< DIM
); d
++)
1042 for (m
= 0; (m
< DIM
); m
++)
1044 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1045 ekin_sum
[gt
][m
][d
] += hm
*v_corrt
[m
]*v_corrt
[d
];
1048 if (md
->nMassPerturbed
&& md
->bPerturbed
[n
])
1051 0.5*(md
->massB
[n
] - md
->massA
[n
])*iprod(v_corrt
, v_corrt
);
1057 for (thread
= 0; thread
< nthread
; thread
++)
1059 for (g
= 0; g
< opts
->ngtc
; g
++)
1063 m_add(tcstat
[g
].ekinf
, ekind
->ekin_work
[thread
][g
],
1068 m_add(tcstat
[g
].ekinh
, ekind
->ekin_work
[thread
][g
],
1073 ekind
->dekindl
+= *ekind
->dekindl_work
[thread
];
1076 inc_nrnb(nrnb
, eNR_EKIN
, md
->homenr
);
1079 static void calc_ke_part_visc(matrix box
, rvec x
[], rvec v
[],
1080 t_grpopts
*opts
, t_mdatoms
*md
,
1081 gmx_ekindata_t
*ekind
,
1082 t_nrnb
*nrnb
, gmx_bool bEkinAveVel
)
1084 int start
= 0, homenr
= md
->homenr
;
1085 int g
, d
, n
, m
, gt
= 0;
1088 t_grp_tcstat
*tcstat
= ekind
->tcstat
;
1089 t_cos_acc
*cosacc
= &(ekind
->cosacc
);
1094 for (g
= 0; g
< opts
->ngtc
; g
++)
1096 copy_mat(ekind
->tcstat
[g
].ekinh
, ekind
->tcstat
[g
].ekinh_old
);
1097 clear_mat(ekind
->tcstat
[g
].ekinh
);
1099 ekind
->dekindl_old
= ekind
->dekindl
;
1101 fac
= 2*M_PI
/box
[ZZ
][ZZ
];
1104 for (n
= start
; n
< start
+homenr
; n
++)
1110 hm
= 0.5*md
->massT
[n
];
1112 /* Note that the times of x and v differ by half a step */
1113 /* MRS -- would have to be changed for VV */
1114 cosz
= cos(fac
*x
[n
][ZZ
]);
1115 /* Calculate the amplitude of the new velocity profile */
1116 mvcos
+= 2*cosz
*md
->massT
[n
]*v
[n
][XX
];
1118 copy_rvec(v
[n
], v_corrt
);
1119 /* Subtract the profile for the kinetic energy */
1120 v_corrt
[XX
] -= cosz
*cosacc
->vcos
;
1121 for (d
= 0; (d
< DIM
); d
++)
1123 for (m
= 0; (m
< DIM
); m
++)
1125 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1128 tcstat
[gt
].ekinf
[m
][d
] += hm
*v_corrt
[m
]*v_corrt
[d
];
1132 tcstat
[gt
].ekinh
[m
][d
] += hm
*v_corrt
[m
]*v_corrt
[d
];
1136 if (md
->nPerturbed
&& md
->bPerturbed
[n
])
1138 /* The minus sign here might be confusing.
1139 * The kinetic contribution from dH/dl doesn't come from
1140 * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
1141 * where p are the momenta. The difference is only a minus sign.
1143 dekindl
-= 0.5*(md
->massB
[n
] - md
->massA
[n
])*iprod(v_corrt
, v_corrt
);
1146 ekind
->dekindl
= dekindl
;
1147 cosacc
->mvcos
= mvcos
;
1149 inc_nrnb(nrnb
, eNR_EKIN
, homenr
);
1152 void calc_ke_part(t_state
*state
, t_grpopts
*opts
, t_mdatoms
*md
,
1153 gmx_ekindata_t
*ekind
, t_nrnb
*nrnb
, gmx_bool bEkinAveVel
)
1155 if (ekind
->cosacc
.cos_accel
== 0)
1157 calc_ke_part_normal(state
->v
, opts
, md
, ekind
, nrnb
, bEkinAveVel
);
1161 calc_ke_part_visc(state
->box
, state
->x
, state
->v
, opts
, md
, ekind
, nrnb
, bEkinAveVel
);
1165 extern void init_ekinstate(ekinstate_t
*ekinstate
, const t_inputrec
*ir
)
1167 ekinstate
->ekin_n
= ir
->opts
.ngtc
;
1168 snew(ekinstate
->ekinh
, ekinstate
->ekin_n
);
1169 snew(ekinstate
->ekinf
, ekinstate
->ekin_n
);
1170 snew(ekinstate
->ekinh_old
, ekinstate
->ekin_n
);
1171 snew(ekinstate
->ekinscalef_nhc
, ekinstate
->ekin_n
);
1172 snew(ekinstate
->ekinscaleh_nhc
, ekinstate
->ekin_n
);
1173 snew(ekinstate
->vscale_nhc
, ekinstate
->ekin_n
);
1174 ekinstate
->dekindl
= 0;
1175 ekinstate
->mvcos
= 0;
1178 void update_ekinstate(ekinstate_t
*ekinstate
, gmx_ekindata_t
*ekind
)
1182 for (i
= 0; i
< ekinstate
->ekin_n
; i
++)
1184 copy_mat(ekind
->tcstat
[i
].ekinh
, ekinstate
->ekinh
[i
]);
1185 copy_mat(ekind
->tcstat
[i
].ekinf
, ekinstate
->ekinf
[i
]);
1186 copy_mat(ekind
->tcstat
[i
].ekinh_old
, ekinstate
->ekinh_old
[i
]);
1187 ekinstate
->ekinscalef_nhc
[i
] = ekind
->tcstat
[i
].ekinscalef_nhc
;
1188 ekinstate
->ekinscaleh_nhc
[i
] = ekind
->tcstat
[i
].ekinscaleh_nhc
;
1189 ekinstate
->vscale_nhc
[i
] = ekind
->tcstat
[i
].vscale_nhc
;
1192 copy_mat(ekind
->ekin
, ekinstate
->ekin_total
);
1193 ekinstate
->dekindl
= ekind
->dekindl
;
1194 ekinstate
->mvcos
= ekind
->cosacc
.mvcos
;
1198 void restore_ekinstate_from_state(t_commrec
*cr
,
1199 gmx_ekindata_t
*ekind
, const ekinstate_t
*ekinstate
)
1205 for (i
= 0; i
< ekinstate
->ekin_n
; i
++)
1207 copy_mat(ekinstate
->ekinh
[i
], ekind
->tcstat
[i
].ekinh
);
1208 copy_mat(ekinstate
->ekinf
[i
], ekind
->tcstat
[i
].ekinf
);
1209 copy_mat(ekinstate
->ekinh_old
[i
], ekind
->tcstat
[i
].ekinh_old
);
1210 ekind
->tcstat
[i
].ekinscalef_nhc
= ekinstate
->ekinscalef_nhc
[i
];
1211 ekind
->tcstat
[i
].ekinscaleh_nhc
= ekinstate
->ekinscaleh_nhc
[i
];
1212 ekind
->tcstat
[i
].vscale_nhc
= ekinstate
->vscale_nhc
[i
];
1215 copy_mat(ekinstate
->ekin_total
, ekind
->ekin
);
1217 ekind
->dekindl
= ekinstate
->dekindl
;
1218 ekind
->cosacc
.mvcos
= ekinstate
->mvcos
;
1219 n
= ekinstate
->ekin_n
;
1224 gmx_bcast(sizeof(n
), &n
, cr
);
1225 for (i
= 0; i
< n
; i
++)
1227 gmx_bcast(DIM
*DIM
*sizeof(ekind
->tcstat
[i
].ekinh
[0][0]),
1228 ekind
->tcstat
[i
].ekinh
[0], cr
);
1229 gmx_bcast(DIM
*DIM
*sizeof(ekind
->tcstat
[i
].ekinf
[0][0]),
1230 ekind
->tcstat
[i
].ekinf
[0], cr
);
1231 gmx_bcast(DIM
*DIM
*sizeof(ekind
->tcstat
[i
].ekinh_old
[0][0]),
1232 ekind
->tcstat
[i
].ekinh_old
[0], cr
);
1234 gmx_bcast(sizeof(ekind
->tcstat
[i
].ekinscalef_nhc
),
1235 &(ekind
->tcstat
[i
].ekinscalef_nhc
), cr
);
1236 gmx_bcast(sizeof(ekind
->tcstat
[i
].ekinscaleh_nhc
),
1237 &(ekind
->tcstat
[i
].ekinscaleh_nhc
), cr
);
1238 gmx_bcast(sizeof(ekind
->tcstat
[i
].vscale_nhc
),
1239 &(ekind
->tcstat
[i
].vscale_nhc
), cr
);
1241 gmx_bcast(DIM
*DIM
*sizeof(ekind
->ekin
[0][0]),
1242 ekind
->ekin
[0], cr
);
1244 gmx_bcast(sizeof(ekind
->dekindl
), &ekind
->dekindl
, cr
);
1245 gmx_bcast(sizeof(ekind
->cosacc
.mvcos
), &ekind
->cosacc
.mvcos
, cr
);
1249 void set_deform_reference_box(gmx_update_t upd
, gmx_int64_t step
, matrix box
)
1251 upd
->deformref_step
= step
;
1252 copy_mat(box
, upd
->deformref_box
);
1255 static void deform(gmx_update_t upd
,
1256 int start
, int homenr
, rvec x
[], matrix box
,
1257 const t_inputrec
*ir
, gmx_int64_t step
)
1259 matrix bnew
, invbox
, mu
;
1263 elapsed_time
= (step
+ 1 - upd
->deformref_step
)*ir
->delta_t
;
1264 copy_mat(box
, bnew
);
1265 for (i
= 0; i
< DIM
; i
++)
1267 for (j
= 0; j
< DIM
; j
++)
1269 if (ir
->deform
[i
][j
] != 0)
1272 upd
->deformref_box
[i
][j
] + elapsed_time
*ir
->deform
[i
][j
];
1276 /* We correct the off-diagonal elements,
1277 * which can grow indefinitely during shearing,
1278 * so the shifts do not get messed up.
1280 for (i
= 1; i
< DIM
; i
++)
1282 for (j
= i
-1; j
>= 0; j
--)
1284 while (bnew
[i
][j
] - box
[i
][j
] > 0.5*bnew
[j
][j
])
1286 rvec_dec(bnew
[i
], bnew
[j
]);
1288 while (bnew
[i
][j
] - box
[i
][j
] < -0.5*bnew
[j
][j
])
1290 rvec_inc(bnew
[i
], bnew
[j
]);
1294 m_inv_ur0(box
, invbox
);
1295 copy_mat(bnew
, box
);
1296 mmul_ur0(box
, invbox
, mu
);
1298 for (i
= start
; i
< start
+homenr
; i
++)
1300 x
[i
][XX
] = mu
[XX
][XX
]*x
[i
][XX
]+mu
[YY
][XX
]*x
[i
][YY
]+mu
[ZZ
][XX
]*x
[i
][ZZ
];
1301 x
[i
][YY
] = mu
[YY
][YY
]*x
[i
][YY
]+mu
[ZZ
][YY
]*x
[i
][ZZ
];
1302 x
[i
][ZZ
] = mu
[ZZ
][ZZ
]*x
[i
][ZZ
];
1306 void update_tcouple(gmx_int64_t step
,
1307 t_inputrec
*inputrec
,
1309 gmx_ekindata_t
*ekind
,
1314 gmx_bool bTCouple
= FALSE
;
1318 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1319 if (inputrec
->etc
!= etcNO
&&
1320 !(inputrecNvtTrotter(inputrec
) || inputrecNptTrotter(inputrec
) || inputrecNphTrotter(inputrec
)))
1322 /* We should only couple after a step where energies were determined (for leapfrog versions)
1323 or the step energies are determined, for velocity verlet versions */
1325 if (EI_VV(inputrec
->eI
))
1333 bTCouple
= (inputrec
->nsttcouple
== 1 ||
1334 do_per_step(step
+inputrec
->nsttcouple
-offset
,
1335 inputrec
->nsttcouple
));
1340 dttc
= inputrec
->nsttcouple
*inputrec
->delta_t
;
1342 switch (inputrec
->etc
)
1347 berendsen_tcoupl(inputrec
, ekind
, dttc
);
1350 nosehoover_tcoupl(&(inputrec
->opts
), ekind
, dttc
,
1351 state
->nosehoover_xi
, state
->nosehoover_vxi
, MassQ
);
1354 vrescale_tcoupl(inputrec
, step
, ekind
, dttc
,
1355 state
->therm_integral
);
1358 /* rescale in place here */
1359 if (EI_VV(inputrec
->eI
))
1361 rescale_velocities(ekind
, md
, 0, md
->homenr
, state
->v
);
1366 /* Set the T scaling lambda to 1 to have no scaling */
1367 for (i
= 0; (i
< inputrec
->opts
.ngtc
); i
++)
1369 ekind
->tcstat
[i
].lambda
= 1.0;
1374 void update_pcouple(FILE *fplog
,
1376 t_inputrec
*inputrec
,
1382 gmx_bool bPCouple
= FALSE
;
1386 /* if using Trotter pressure, we do this in coupling.c, so we leave it false. */
1387 if (inputrec
->epc
!= epcNO
&& (!(inputrecNptTrotter(inputrec
) || inputrecNphTrotter(inputrec
))))
1389 /* We should only couple after a step where energies were determined */
1390 bPCouple
= (inputrec
->nstpcouple
== 1 ||
1391 do_per_step(step
+inputrec
->nstpcouple
-1,
1392 inputrec
->nstpcouple
));
1395 clear_mat(pcoupl_mu
);
1396 for (i
= 0; i
< DIM
; i
++)
1398 pcoupl_mu
[i
][i
] = 1.0;
1405 dtpc
= inputrec
->nstpcouple
*inputrec
->delta_t
;
1407 switch (inputrec
->epc
)
1409 /* We can always pcoupl, even if we did not sum the energies
1410 * the previous step, since state->pres_prev is only updated
1411 * when the energies have been summed.
1415 case (epcBERENDSEN
):
1418 berendsen_pcoupl(fplog
, step
, inputrec
, dtpc
, state
->pres_prev
, state
->box
,
1422 case (epcPARRINELLORAHMAN
):
1423 parrinellorahman_pcoupl(fplog
, step
, inputrec
, dtpc
, state
->pres_prev
,
1424 state
->box
, state
->box_rel
, state
->boxv
,
1425 M
, pcoupl_mu
, bInitStep
);
1433 static rvec
*get_xprime(const t_state
*state
, gmx_update_t upd
)
1435 if (state
->nalloc
> upd
->xp_nalloc
)
1437 upd
->xp_nalloc
= state
->nalloc
;
1438 srenew(upd
->xp
, upd
->xp_nalloc
);
1444 void update_constraints(FILE *fplog
,
1446 real
*dvdlambda
, /* the contribution to be added to the bonded interactions */
1447 t_inputrec
*inputrec
, /* input record and box stuff */
1452 rvec force
[], /* forces on home particles */
1457 gmx_wallcycle_t wcycle
,
1459 gmx_constr_t constr
,
1460 gmx_bool bFirstHalf
,
1463 gmx_bool bLastStep
, bLog
= FALSE
, bEner
= FALSE
, bDoConstr
= FALSE
;
1465 int start
, homenr
, nrend
, i
, m
;
1467 rvec
*xprime
= NULL
;
1474 if (bFirstHalf
&& !EI_VV(inputrec
->eI
))
1479 /* for now, SD update is here -- though it really seems like it
1480 should be reformulated as a velocity verlet method, since it has two parts */
1483 homenr
= md
->homenr
;
1484 nrend
= start
+homenr
;
1486 dt
= inputrec
->delta_t
;
1490 * APPLY CONSTRAINTS:
1493 * When doing PR pressure coupling we have to constrain the
1494 * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1495 * it is enough to do this once though, since the relative velocities
1496 * after this will be normal to the bond vector
1501 /* clear out constraints before applying */
1502 clear_mat(vir_part
);
1504 xprime
= get_xprime(state
, upd
);
1506 bLastStep
= (step
== inputrec
->init_step
+inputrec
->nsteps
);
1507 bLog
= (do_per_step(step
, inputrec
->nstlog
) || bLastStep
|| (step
< 0));
1508 bEner
= (do_per_step(step
, inputrec
->nstenergy
) || bLastStep
);
1509 /* Constrain the coordinates xprime */
1510 wallcycle_start(wcycle
, ewcCONSTR
);
1511 if (EI_VV(inputrec
->eI
) && bFirstHalf
)
1513 constrain(NULL
, bLog
, bEner
, constr
, idef
,
1514 inputrec
, cr
, step
, 1, 1.0, md
,
1515 state
->x
, state
->v
, state
->v
,
1516 bMolPBC
, state
->box
,
1517 state
->lambda
[efptBONDED
], dvdlambda
,
1518 NULL
, bCalcVir
? &vir_con
: NULL
, nrnb
, econqVeloc
);
1522 constrain(NULL
, bLog
, bEner
, constr
, idef
,
1523 inputrec
, cr
, step
, 1, 1.0, md
,
1524 state
->x
, xprime
, NULL
,
1525 bMolPBC
, state
->box
,
1526 state
->lambda
[efptBONDED
], dvdlambda
,
1527 state
->v
, bCalcVir
? &vir_con
: NULL
, nrnb
, econqCoord
);
1529 wallcycle_stop(wcycle
, ewcCONSTR
);
1533 dump_it_all(fplog
, "After Shake",
1534 state
->natoms
, state
->x
, xprime
, state
->v
, force
);
1538 if (inputrec
->eI
== eiSD2
)
1540 /* A correction factor eph is needed for the SD constraint force */
1541 /* Here we can, unfortunately, not have proper corrections
1542 * for different friction constants, so we use the first one.
1544 for (i
= 0; i
< DIM
; i
++)
1546 for (m
= 0; m
< DIM
; m
++)
1548 vir_part
[i
][m
] += upd
->sd
->sdc
[0].eph
*vir_con
[i
][m
];
1554 m_add(vir_part
, vir_con
, vir_part
);
1558 pr_rvecs(debug
, 0, "constraint virial", vir_part
, DIM
);
1565 if (inputrec
->eI
== eiSD1
&& bDoConstr
&& !bFirstHalf
)
1567 wallcycle_start(wcycle
, ewcUPDATE
);
1568 xprime
= get_xprime(state
, upd
);
1570 nth
= gmx_omp_nthreads_get(emntUpdate
);
1572 #pragma omp parallel for num_threads(nth) schedule(static)
1573 for (th
= 0; th
< nth
; th
++)
1577 int start_th
, end_th
;
1579 start_th
= start
+ ((nrend
-start
)* th
)/nth
;
1580 end_th
= start
+ ((nrend
-start
)*(th
+1))/nth
;
1582 /* The second part of the SD integration */
1583 do_update_sd1(upd
->sd
,
1584 start_th
, end_th
, dt
,
1585 inputrec
->opts
.acc
, inputrec
->opts
.nFreeze
,
1586 md
->invmass
, md
->ptype
,
1587 md
->cFREEZE
, md
->cACC
, md
->cTC
,
1588 state
->x
, xprime
, state
->v
, force
,
1589 inputrec
->opts
.ngtc
, inputrec
->opts
.ref_t
,
1591 step
, inputrec
->ld_seed
,
1592 DOMAINDECOMP(cr
) ? cr
->dd
->gatindex
: NULL
);
1594 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
1596 inc_nrnb(nrnb
, eNR_UPDATE
, homenr
);
1597 wallcycle_stop(wcycle
, ewcUPDATE
);
1601 /* Constrain the coordinates xprime for half a time step */
1602 wallcycle_start(wcycle
, ewcCONSTR
);
1604 constrain(NULL
, bLog
, bEner
, constr
, idef
,
1605 inputrec
, cr
, step
, 1, 0.5, md
,
1606 state
->x
, xprime
, NULL
,
1607 bMolPBC
, state
->box
,
1608 state
->lambda
[efptBONDED
], dvdlambda
,
1609 state
->v
, NULL
, nrnb
, econqCoord
);
1611 wallcycle_stop(wcycle
, ewcCONSTR
);
1615 if ((inputrec
->eI
== eiSD2
) && !(bFirstHalf
))
1617 wallcycle_start(wcycle
, ewcUPDATE
);
1618 xprime
= get_xprime(state
, upd
);
1620 nth
= gmx_omp_nthreads_get(emntUpdate
);
1622 #pragma omp parallel for num_threads(nth) schedule(static)
1623 for (th
= 0; th
< nth
; th
++)
1627 int start_th
, end_th
;
1629 start_th
= start
+ ((nrend
-start
)* th
)/nth
;
1630 end_th
= start
+ ((nrend
-start
)*(th
+1))/nth
;
1632 /* The second part of the SD integration */
1633 do_update_sd2(upd
->sd
,
1634 FALSE
, start_th
, end_th
,
1635 inputrec
->opts
.acc
, inputrec
->opts
.nFreeze
,
1636 md
->invmass
, md
->ptype
,
1637 md
->cFREEZE
, md
->cACC
, md
->cTC
,
1638 state
->x
, xprime
, state
->v
, force
, state
->sd_X
,
1639 inputrec
->opts
.tau_t
,
1640 FALSE
, step
, inputrec
->ld_seed
,
1641 DOMAINDECOMP(cr
) ? cr
->dd
->gatindex
: NULL
);
1643 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
1645 inc_nrnb(nrnb
, eNR_UPDATE
, homenr
);
1646 wallcycle_stop(wcycle
, ewcUPDATE
);
1650 /* Constrain the coordinates xprime */
1651 wallcycle_start(wcycle
, ewcCONSTR
);
1652 constrain(NULL
, bLog
, bEner
, constr
, idef
,
1653 inputrec
, cr
, step
, 1, 1.0, md
,
1654 state
->x
, xprime
, NULL
,
1655 bMolPBC
, state
->box
,
1656 state
->lambda
[efptBONDED
], dvdlambda
,
1657 NULL
, NULL
, nrnb
, econqCoord
);
1658 wallcycle_stop(wcycle
, ewcCONSTR
);
1663 /* We must always unshift after updating coordinates; if we did not shake
1664 x was shifted in do_force */
1666 if (!(bFirstHalf
)) /* in the first half of vv, no shift. */
1668 /* NOTE This part of the update actually does not belong with
1669 * the constraints, since we also call it without constraints.
1670 * But currently we always integrate to a temporary buffer and
1671 * then copy the results back here.
1673 wallcycle_start_nocount(wcycle
, ewcUPDATE
);
1675 if (graph
&& (graph
->nnodes
> 0))
1677 unshift_x(graph
, state
->box
, state
->x
, upd
->xp
);
1678 if (TRICLINIC(state
->box
))
1680 inc_nrnb(nrnb
, eNR_SHIFTX
, 2*graph
->nnodes
);
1684 inc_nrnb(nrnb
, eNR_SHIFTX
, graph
->nnodes
);
1689 #ifndef __clang_analyzer__
1690 // cppcheck-suppress unreadVariable
1691 nth
= gmx_omp_nthreads_get(emntUpdate
);
1693 #pragma omp parallel for num_threads(nth) schedule(static)
1694 for (i
= start
; i
< nrend
; i
++)
1696 // Trivial statement, does not throw
1697 copy_rvec(upd
->xp
[i
], state
->x
[i
]);
1700 wallcycle_stop(wcycle
, ewcUPDATE
);
1702 dump_it_all(fplog
, "After unshift",
1703 state
->natoms
, state
->x
, upd
->xp
, state
->v
, force
);
1705 /* ############# END the update of velocities and positions ######### */
1708 void update_box(FILE *fplog
,
1710 t_inputrec
*inputrec
, /* input record and box stuff */
1713 rvec force
[], /* forces on home particles */
1719 int start
, homenr
, i
, n
, m
;
1722 homenr
= md
->homenr
;
1724 dt
= inputrec
->delta_t
;
1728 /* now update boxes */
1729 switch (inputrec
->epc
)
1733 case (epcBERENDSEN
):
1734 /* We should only scale after a step where the pressure (kinetic
1735 * energy and virial) was calculated. This happens after the
1736 * coordinate update, whereas the current routine is called before
1737 * that, so we scale when step % nstpcouple = 1 instead of 0.
1739 if (inputrec
->nstpcouple
== 1 || (step
% inputrec
->nstpcouple
== 1))
1741 berendsen_pscale(inputrec
, pcoupl_mu
, state
->box
, state
->box_rel
,
1742 start
, homenr
, state
->x
, md
->cFREEZE
, nrnb
);
1745 case (epcPARRINELLORAHMAN
):
1746 /* The box velocities were updated in do_pr_pcoupl in the update
1747 * iteration, but we dont change the box vectors until we get here
1748 * since we need to be able to shift/unshift above.
1750 for (i
= 0; i
< DIM
; i
++)
1752 for (m
= 0; m
<= i
; m
++)
1754 state
->box
[i
][m
] += dt
*state
->boxv
[i
][m
];
1757 preserve_box_shape(inputrec
, state
->box_rel
, state
->box
);
1759 /* Scale the coordinates */
1760 for (n
= start
; (n
< start
+homenr
); n
++)
1762 tmvmul_ur0(pcoupl_mu
, state
->x
[n
], state
->x
[n
]);
1766 switch (inputrec
->epct
)
1768 case (epctISOTROPIC
):
1769 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1770 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1771 Side length scales as exp(veta*dt) */
1773 msmul(state
->box
, exp(state
->veta
*dt
), state
->box
);
1775 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1776 o If we assume isotropic scaling, and box length scaling
1777 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1778 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1779 determinant of B is L^DIM det(M), and the determinant
1780 of dB/dt is (dL/dT)^DIM det (M). veta will be
1781 (det(dB/dT)/det(B))^(1/3). Then since M =
1782 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1784 msmul(state
->box
, state
->veta
, state
->boxv
);
1794 if (inputrecDeform(inputrec
))
1796 deform(upd
, start
, homenr
, state
->x
, state
->box
, inputrec
, step
);
1799 dump_it_all(fplog
, "After update",
1800 state
->natoms
, state
->x
, upd
->xp
, state
->v
, force
);
1803 void update_coords(FILE *fplog
,
1805 t_inputrec
*inputrec
, /* input record and box stuff */
1808 rvec
*f
, /* forces on home particles */
1810 gmx_ekindata_t
*ekind
,
1815 t_commrec
*cr
, /* these shouldn't be here -- need to think about it */
1816 gmx_constr_t constr
)
1818 gmx_bool bNH
, bPR
, bDoConstr
= FALSE
;
1820 int start
, homenr
, nrend
;
1824 bDoConstr
= (NULL
!= constr
);
1826 /* Running the velocity half does nothing except for velocity verlet */
1827 if ((UpdatePart
== etrtVELOCITY1
|| UpdatePart
== etrtVELOCITY2
) &&
1828 !EI_VV(inputrec
->eI
))
1830 gmx_incons("update_coords called for velocity without VV integrator");
1834 homenr
= md
->homenr
;
1835 nrend
= start
+homenr
;
1837 xprime
= get_xprime(state
, upd
);
1839 dt
= inputrec
->delta_t
;
1841 /* We need to update the NMR restraint history when time averaging is used */
1842 if (state
->flags
& (1<<estDISRE_RM3TAV
))
1844 update_disres_history(fcd
, &state
->hist
);
1846 if (state
->flags
& (1<<estORIRE_DTAV
))
1848 update_orires_history(fcd
, &state
->hist
);
1852 bNH
= inputrec
->etc
== etcNOSEHOOVER
;
1853 bPR
= ((inputrec
->epc
== epcPARRINELLORAHMAN
) || (inputrec
->epc
== epcMTTK
));
1855 /* ############# START The update of velocities and positions ######### */
1857 dump_it_all(fplog
, "Before update",
1858 state
->natoms
, state
->x
, xprime
, state
->v
, f
);
1860 if (inputrec
->eI
== eiSD2
)
1862 check_sd2_work_data_allocation(upd
->sd
, nrend
);
1864 do_update_sd2_Tconsts(upd
->sd
,
1865 inputrec
->opts
.ngtc
,
1866 inputrec
->opts
.tau_t
,
1867 inputrec
->opts
.ref_t
);
1869 if (inputrec
->eI
== eiBD
)
1871 do_update_bd_Tconsts(dt
, inputrec
->bd_fric
,
1872 inputrec
->opts
.ngtc
, inputrec
->opts
.ref_t
,
1876 nth
= gmx_omp_nthreads_get(emntUpdate
);
1878 #pragma omp parallel for num_threads(nth) schedule(static) private(alpha)
1879 for (th
= 0; th
< nth
; th
++)
1883 int start_th
, end_th
;
1885 start_th
= start
+ ((nrend
-start
)* th
)/nth
;
1886 end_th
= start
+ ((nrend
-start
)*(th
+1))/nth
;
1888 switch (inputrec
->eI
)
1891 if (ekind
->cosacc
.cos_accel
== 0)
1893 do_update_md(start_th
, end_th
, dt
,
1894 ekind
->tcstat
, state
->nosehoover_vxi
,
1895 ekind
->bNEMD
, ekind
->grpstat
, inputrec
->opts
.acc
,
1896 inputrec
->opts
.nFreeze
,
1897 md
->invmass
, md
->ptype
,
1898 md
->cFREEZE
, md
->cACC
, md
->cTC
,
1899 state
->x
, xprime
, state
->v
, f
, M
,
1904 do_update_visc(start_th
, end_th
, dt
,
1905 ekind
->tcstat
, state
->nosehoover_vxi
,
1906 md
->invmass
, md
->ptype
,
1907 md
->cTC
, state
->x
, xprime
, state
->v
, f
, M
,
1909 ekind
->cosacc
.cos_accel
,
1915 /* With constraints, the SD1 update is done in 2 parts */
1916 do_update_sd1(upd
->sd
,
1917 start_th
, end_th
, dt
,
1918 inputrec
->opts
.acc
, inputrec
->opts
.nFreeze
,
1919 md
->invmass
, md
->ptype
,
1920 md
->cFREEZE
, md
->cACC
, md
->cTC
,
1921 state
->x
, xprime
, state
->v
, f
,
1922 inputrec
->opts
.ngtc
, inputrec
->opts
.ref_t
,
1924 step
, inputrec
->ld_seed
, DOMAINDECOMP(cr
) ? cr
->dd
->gatindex
: NULL
);
1927 /* The SD2 update is always done in 2 parts,
1928 * because an extra constraint step is needed
1930 do_update_sd2(upd
->sd
,
1931 bInitStep
, start_th
, end_th
,
1932 inputrec
->opts
.acc
, inputrec
->opts
.nFreeze
,
1933 md
->invmass
, md
->ptype
,
1934 md
->cFREEZE
, md
->cACC
, md
->cTC
,
1935 state
->x
, xprime
, state
->v
, f
, state
->sd_X
,
1936 inputrec
->opts
.tau_t
,
1937 TRUE
, step
, inputrec
->ld_seed
,
1938 DOMAINDECOMP(cr
) ? cr
->dd
->gatindex
: NULL
);
1941 do_update_bd(start_th
, end_th
, dt
,
1942 inputrec
->opts
.nFreeze
, md
->invmass
, md
->ptype
,
1943 md
->cFREEZE
, md
->cTC
,
1944 state
->x
, xprime
, state
->v
, f
,
1947 step
, inputrec
->ld_seed
, DOMAINDECOMP(cr
) ? cr
->dd
->gatindex
: NULL
);
1951 alpha
= 1.0 + DIM
/((double)inputrec
->opts
.nrdf
[0]); /* assuming barostat coupled to group 0. */
1956 do_update_vv_vel(start_th
, end_th
, dt
,
1957 inputrec
->opts
.acc
, inputrec
->opts
.nFreeze
,
1958 md
->invmass
, md
->ptype
,
1959 md
->cFREEZE
, md
->cACC
,
1961 (bNH
|| bPR
), state
->veta
, alpha
);
1964 do_update_vv_pos(start_th
, end_th
, dt
,
1965 inputrec
->opts
.nFreeze
,
1966 md
->ptype
, md
->cFREEZE
,
1967 state
->x
, xprime
, state
->v
,
1968 (bNH
|| bPR
), state
->veta
);
1973 gmx_fatal(FARGS
, "Don't know how to update coordinates");
1977 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
1983 void correct_ekin(FILE *log
, int start
, int end
, rvec v
[], rvec vcm
, real mass
[],
1984 real tmass
, tensor ekin
)
1987 * This is a debugging routine. It should not be called for production code
1989 * The kinetic energy should calculated according to:
1990 * Ekin = 1/2 m (v-vcm)^2
1991 * However the correction is not always applied, since vcm may not be
1992 * known in time and we compute
1993 * Ekin' = 1/2 m v^2 instead
1994 * This can be corrected afterwards by computing
1995 * Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
1997 * Ekin = Ekin' - m v vcm + 1/2 m vcm^2
2004 /* Local particles */
2007 /* Processor dependent part. */
2009 for (i
= start
; (i
< end
); i
++)
2013 for (j
= 0; (j
< DIM
); j
++)
2019 svmul(1/tmass
, vcm
, vcm
);
2020 svmul(0.5, vcm
, hvcm
);
2022 for (j
= 0; (j
< DIM
); j
++)
2024 for (k
= 0; (k
< DIM
); k
++)
2026 dekin
[j
][k
] += vcm
[k
]*(tm
*hvcm
[j
]-mv
[j
]);
2029 pr_rvecs(log
, 0, "dekin", dekin
, DIM
);
2030 pr_rvecs(log
, 0, " ekin", ekin
, DIM
);
2031 fprintf(log
, "dekin = %g, ekin = %g vcm = (%8.4f %8.4f %8.4f)\n",
2032 trace(dekin
), trace(ekin
), vcm
[XX
], vcm
[YY
], vcm
[ZZ
]);
2033 fprintf(log
, "mv = (%8.4f %8.4f %8.4f)\n",
2034 mv
[XX
], mv
[YY
], mv
[ZZ
]);
2037 extern gmx_bool
update_randomize_velocities(t_inputrec
*ir
, gmx_int64_t step
, const t_commrec
*cr
,
2038 t_mdatoms
*md
, t_state
*state
, gmx_update_t upd
, gmx_constr_t constr
)
2041 real rate
= (ir
->delta_t
)/ir
->opts
.tau_t
[0];
2043 if (ir
->etc
== etcANDERSEN
&& constr
!= NULL
)
2045 /* Currently, Andersen thermostat does not support constrained
2046 systems. Functionality exists in the andersen_tcoupl
2047 function in GROMACS 4.5.7 to allow this combination. That
2048 code could be ported to the current random-number
2049 generation approach, but has not yet been done because of
2050 lack of time and resources. */
2051 gmx_fatal(FARGS
, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
2054 /* proceed with andersen if 1) it's fixed probability per
2055 particle andersen or 2) it's massive andersen and it's tau_t/dt */
2056 if ((ir
->etc
== etcANDERSEN
) || do_per_step(step
, (int)(1.0/rate
)))
2058 andersen_tcoupl(ir
, step
, cr
, md
, state
, rate
,
2059 upd
->sd
->randomize_group
, upd
->sd
->boltzfac
);