f4561d03217086bcda32924c2a31c66951284e10
[alexxy/gromacs.git] / src / gromacs / mdlib / md_support.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37
38 #include "gmxpre.h"
39
40 #include "md_support.h"
41
42 #include <climits>
43 #include <cmath>
44
45 #include <algorithm>
46
47 #include "gromacs/domdec/domdec.h"
48 #include "gromacs/gmxlib/network.h"
49 #include "gromacs/gmxlib/nrnb.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/mdlib/mdrun.h"
52 #include "gromacs/mdlib/sim_util.h"
53 #include "gromacs/mdlib/simulationsignal.h"
54 #include "gromacs/mdlib/tgroup.h"
55 #include "gromacs/mdlib/update.h"
56 #include "gromacs/mdlib/vcm.h"
57 #include "gromacs/mdtypes/commrec.h"
58 #include "gromacs/mdtypes/df_history.h"
59 #include "gromacs/mdtypes/energyhistory.h"
60 #include "gromacs/mdtypes/group.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/mdtypes/state.h"
64 #include "gromacs/pbcutil/pbc.h"
65 #include "gromacs/timing/wallcycle.h"
66 #include "gromacs/topology/mtop_util.h"
67 #include "gromacs/trajectory/trajectoryframe.h"
68 #include "gromacs/utility/arrayref.h"
69 #include "gromacs/utility/cstringutil.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/gmxassert.h"
72 #include "gromacs/utility/logger.h"
73 #include "gromacs/utility/smalloc.h"
74 #include "gromacs/utility/snprintf.h"
75
76 // TODO move this to multi-sim module
77 bool multisim_int_all_are_equal(const gmx_multisim_t *ms,
78                                 gmx_int64_t           value)
79 {
80     bool         allValuesAreEqual = true;
81     gmx_int64_t *buf;
82
83     GMX_RELEASE_ASSERT(ms, "Invalid use of multi-simulation pointer");
84
85     snew(buf, ms->nsim);
86     /* send our value to all other master ranks, receive all of theirs */
87     buf[ms->sim] = value;
88     gmx_sumli_sim(ms->nsim, buf, ms);
89
90     for (int s = 0; s < ms->nsim; s++)
91     {
92         if (buf[s] != value)
93         {
94             allValuesAreEqual = false;
95             break;
96         }
97     }
98
99     sfree(buf);
100
101     return allValuesAreEqual;
102 }
103
104 int multisim_min(const gmx_multisim_t *ms, int nmin, int n)
105 {
106     int     *buf;
107     gmx_bool bPos, bEqual;
108     int      s, d;
109
110     snew(buf, ms->nsim);
111     buf[ms->sim] = n;
112     gmx_sumi_sim(ms->nsim, buf, ms);
113     bPos   = TRUE;
114     bEqual = TRUE;
115     for (s = 0; s < ms->nsim; s++)
116     {
117         bPos   = bPos   && (buf[s] > 0);
118         bEqual = bEqual && (buf[s] == buf[0]);
119     }
120     if (bPos)
121     {
122         if (bEqual)
123         {
124             nmin = std::min(nmin, buf[0]);
125         }
126         else
127         {
128             /* Find the least common multiple */
129             for (d = 2; d < nmin; d++)
130             {
131                 s = 0;
132                 while (s < ms->nsim && d % buf[s] == 0)
133                 {
134                     s++;
135                 }
136                 if (s == ms->nsim)
137                 {
138                     /* We found the LCM and it is less than nmin */
139                     nmin = d;
140                     break;
141                 }
142             }
143         }
144     }
145     sfree(buf);
146
147     return nmin;
148 }
149
150 /* TODO Specialize this routine into init-time and loop-time versions?
151    e.g. bReadEkin is only true when restoring from checkpoint */
152 void compute_globals(FILE *fplog, gmx_global_stat *gstat, t_commrec *cr, t_inputrec *ir,
153                      t_forcerec *fr, gmx_ekindata_t *ekind,
154                      t_state *state, t_mdatoms *mdatoms,
155                      t_nrnb *nrnb, t_vcm *vcm, gmx_wallcycle_t wcycle,
156                      gmx_enerdata_t *enerd, tensor force_vir, tensor shake_vir, tensor total_vir,
157                      tensor pres, rvec mu_tot, gmx::Constraints *constr,
158                      gmx::SimulationSignaller *signalCoordinator,
159                      matrix box, int *totalNumberOfBondedInteractions,
160                      gmx_bool *bSumEkinhOld, int flags)
161 {
162     tensor   corr_vir, corr_pres;
163     gmx_bool bEner, bPres, bTemp;
164     gmx_bool bStopCM, bGStat,
165              bReadEkin, bEkinAveVel, bScaleEkin, bConstrain;
166     real     prescorr, enercorr, dvdlcorr, dvdl_ekin;
167
168     /* translate CGLO flags to gmx_booleans */
169     bStopCM       = flags & CGLO_STOPCM;
170     bGStat        = flags & CGLO_GSTAT;
171     bReadEkin     = (flags & CGLO_READEKIN);
172     bScaleEkin    = (flags & CGLO_SCALEEKIN);
173     bEner         = flags & CGLO_ENERGY;
174     bTemp         = flags & CGLO_TEMPERATURE;
175     bPres         = (flags & CGLO_PRESSURE);
176     bConstrain    = (flags & CGLO_CONSTRAINT);
177
178     /* we calculate a full state kinetic energy either with full-step velocity verlet
179        or half step where we need the pressure */
180
181     bEkinAveVel = (ir->eI == eiVV || (ir->eI == eiVVAK && bPres) || bReadEkin);
182
183     /* in initalization, it sums the shake virial in vv, and to
184        sums ekinh_old in leapfrog (or if we are calculating ekinh_old) for other reasons */
185
186     /* ########## Kinetic energy  ############## */
187
188     if (bTemp)
189     {
190         /* Non-equilibrium MD: this is parallellized, but only does communication
191          * when there really is NEMD.
192          */
193
194         if (PAR(cr) && (ekind->bNEMD))
195         {
196             accumulate_u(cr, &(ir->opts), ekind);
197         }
198         if (!bReadEkin)
199         {
200             calc_ke_part(state, &(ir->opts), mdatoms, ekind, nrnb, bEkinAveVel);
201         }
202     }
203
204     /* Calculate center of mass velocity if necessary, also parallellized */
205     if (bStopCM)
206     {
207         calc_vcm_grp(0, mdatoms->homenr, mdatoms,
208                      as_rvec_array(state->x.data()), as_rvec_array(state->v.data()), vcm);
209     }
210
211     if (bTemp || bStopCM || bPres || bEner || bConstrain)
212     {
213         if (!bGStat)
214         {
215             /* We will not sum ekinh_old,
216              * so signal that we still have to do it.
217              */
218             *bSumEkinhOld = TRUE;
219
220         }
221         else
222         {
223             gmx::ArrayRef<real> signalBuffer = signalCoordinator->getCommunicationBuffer();
224             if (PAR(cr))
225             {
226                 wallcycle_start(wcycle, ewcMoveE);
227                 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
228                             ir, ekind, constr, bStopCM ? vcm : nullptr,
229                             signalBuffer.size(), signalBuffer.data(),
230                             totalNumberOfBondedInteractions,
231                             *bSumEkinhOld, flags);
232                 wallcycle_stop(wcycle, ewcMoveE);
233             }
234             signalCoordinator->finalizeSignals();
235             *bSumEkinhOld = FALSE;
236         }
237     }
238
239     /* Do center of mass motion removal */
240     if (bStopCM)
241     {
242         check_cm_grp(fplog, vcm, ir, 1);
243         /* At initialization, do not pass x with acceleration-correction mode
244          * to avoid (incorrect) correction of the initial coordinates.
245          */
246         rvec *xPtr = nullptr;
247         if (vcm->mode == ecmANGULAR || (vcm->mode == ecmLINEAR_ACCELERATION_CORRECTION && !(flags & CGLO_INITIALIZATION)))
248         {
249             xPtr = as_rvec_array(state->x.data());
250         }
251         do_stopcm_grp(mdatoms->homenr, mdatoms->cVCM,
252                       xPtr, as_rvec_array(state->v.data()), *vcm);
253         inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
254     }
255
256     if (bEner)
257     {
258         /* Calculate the amplitude of the cosine velocity profile */
259         ekind->cosacc.vcos = ekind->cosacc.mvcos/mdatoms->tmass;
260     }
261
262     if (bTemp)
263     {
264         /* Sum the kinetic energies of the groups & calc temp */
265         /* compute full step kinetic energies if vv, or if vv-avek and we are computing the pressure with inputrecNptTrotter */
266         /* three maincase:  VV with AveVel (md-vv), vv with AveEkin (md-vv-avek), leap with AveEkin (md).
267            Leap with AveVel is not supported; it's not clear that it will actually work.
268            bEkinAveVel: If TRUE, we simply multiply ekin by ekinscale to get a full step kinetic energy.
269            If FALSE, we average ekinh_old and ekinh*ekinscale_nhc to get an averaged half step kinetic energy.
270          */
271         enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, &dvdl_ekin,
272                                        bEkinAveVel, bScaleEkin);
273         enerd->dvdl_lin[efptMASS] = (double) dvdl_ekin;
274
275         enerd->term[F_EKIN] = trace(ekind->ekin);
276     }
277
278     /* ##########  Long range energy information ###### */
279
280     if (bEner || bPres || bConstrain)
281     {
282         calc_dispcorr(ir, fr, box, state->lambda[efptVDW],
283                       corr_pres, corr_vir, &prescorr, &enercorr, &dvdlcorr);
284     }
285
286     if (bEner)
287     {
288         enerd->term[F_DISPCORR]  = enercorr;
289         enerd->term[F_EPOT]     += enercorr;
290         enerd->term[F_DVDL_VDW] += dvdlcorr;
291     }
292
293     /* ########## Now pressure ############## */
294     if (bPres || bConstrain)
295     {
296
297         m_add(force_vir, shake_vir, total_vir);
298
299         /* Calculate pressure and apply LR correction if PPPM is used.
300          * Use the box from last timestep since we already called update().
301          */
302
303         enerd->term[F_PRES] = calc_pres(fr->ePBC, ir->nwall, box, ekind->ekin, total_vir, pres);
304
305         /* Calculate long range corrections to pressure and energy */
306         /* this adds to enerd->term[F_PRES] and enerd->term[F_ETOT],
307            and computes enerd->term[F_DISPCORR].  Also modifies the
308            total_vir and pres tesors */
309
310         m_add(total_vir, corr_vir, total_vir);
311         m_add(pres, corr_pres, pres);
312         enerd->term[F_PDISPCORR] = prescorr;
313         enerd->term[F_PRES]     += prescorr;
314     }
315 }
316
317 /* check whether an 'nst'-style parameter p is a multiple of nst, and
318    set it to be one if not, with a warning. */
319 static void check_nst_param(const gmx::MDLogger &mdlog,
320                             const char *desc_nst, int nst,
321                             const char *desc_p, int *p)
322 {
323     if (*p > 0 && *p % nst != 0)
324     {
325         /* Round up to the next multiple of nst */
326         *p = ((*p)/nst + 1)*nst;
327         GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
328                 "NOTE: %s changes %s to %d", desc_nst, desc_p, *p);
329     }
330 }
331
332 void setCurrentLambdasRerun(gmx_int64_t step, const t_lambda *fepvals,
333                             const t_trxframe *rerun_fr, const double *lam0,
334                             t_state *globalState)
335 {
336     GMX_RELEASE_ASSERT(globalState != nullptr, "setCurrentLambdasGlobalRerun should be called with a valid state object");
337
338     if (rerun_fr->bLambda)
339     {
340         if (fepvals->delta_lambda == 0)
341         {
342             globalState->lambda[efptFEP] = rerun_fr->lambda;
343         }
344         else
345         {
346             /* find out between which two value of lambda we should be */
347             real frac      = step*fepvals->delta_lambda;
348             int  fep_state = static_cast<int>(std::floor(frac*fepvals->n_lambda));
349             /* interpolate between this state and the next */
350             /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
351             frac           = frac*fepvals->n_lambda - fep_state;
352             for (int i = 0; i < efptNR; i++)
353             {
354                 globalState->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
355                     frac*(fepvals->all_lambda[i][fep_state+1] - fepvals->all_lambda[i][fep_state]);
356             }
357         }
358     }
359     else if (rerun_fr->bFepState)
360     {
361         globalState->fep_state = rerun_fr->fep_state;
362         for (int i = 0; i < efptNR; i++)
363         {
364             globalState->lambda[i] = fepvals->all_lambda[i][globalState->fep_state];
365         }
366     }
367 }
368
369 void setCurrentLambdasLocal(gmx_int64_t step, const t_lambda *fepvals,
370                             const double *lam0, t_state *state)
371 /* find the current lambdas.  If rerunning, we either read in a state, or a lambda value,
372    requiring different logic. */
373 {
374     if (fepvals->delta_lambda != 0)
375     {
376         /* find out between which two value of lambda we should be */
377         real frac = step*fepvals->delta_lambda;
378         if (fepvals->n_lambda > 0)
379         {
380             int fep_state = static_cast<int>(std::floor(frac*fepvals->n_lambda));
381             /* interpolate between this state and the next */
382             /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
383             frac          = frac*fepvals->n_lambda - fep_state;
384             for (int i = 0; i < efptNR; i++)
385             {
386                 state->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
387                     frac*(fepvals->all_lambda[i][fep_state + 1] - fepvals->all_lambda[i][fep_state]);
388             }
389         }
390         else
391         {
392             for (int i = 0; i < efptNR; i++)
393             {
394                 state->lambda[i] = lam0[i] + frac;
395             }
396         }
397     }
398     else
399     {
400         /* if < 0, fep_state was never defined, and we should not set lambda from the state */
401         if (state->fep_state > -1)
402         {
403             for (int i = 0; i < efptNR; i++)
404             {
405                 state->lambda[i] = fepvals->all_lambda[i][state->fep_state];
406             }
407         }
408     }
409 }
410
411 static void min_zero(int *n, int i)
412 {
413     if (i > 0 && (*n == 0 || i < *n))
414     {
415         *n = i;
416     }
417 }
418
419 static int lcd4(int i1, int i2, int i3, int i4)
420 {
421     int nst;
422
423     nst = 0;
424     min_zero(&nst, i1);
425     min_zero(&nst, i2);
426     min_zero(&nst, i3);
427     min_zero(&nst, i4);
428     if (nst == 0)
429     {
430         gmx_incons("All 4 inputs for determining nstglobalcomm are <= 0");
431     }
432
433     while (nst > 1 && ((i1 > 0 && i1 % nst != 0)  ||
434                        (i2 > 0 && i2 % nst != 0)  ||
435                        (i3 > 0 && i3 % nst != 0)  ||
436                        (i4 > 0 && i4 % nst != 0)))
437     {
438         nst--;
439     }
440
441     return nst;
442 }
443
444 int check_nstglobalcomm(const gmx::MDLogger &mdlog, int nstglobalcomm, t_inputrec *ir)
445 {
446     if (!EI_DYNAMICS(ir->eI))
447     {
448         nstglobalcomm = 1;
449     }
450
451     if (nstglobalcomm == -1)
452     {
453         // Set up the default behaviour
454         if (!(ir->nstcalcenergy > 0 ||
455               ir->nstlist > 0 ||
456               ir->etc != etcNO ||
457               ir->epc != epcNO))
458         {
459             /* The user didn't choose the period for anything
460                important, so we just make sure we can send signals and
461                write output suitably. */
462             nstglobalcomm = 10;
463             if (ir->nstenergy > 0 && ir->nstenergy < nstglobalcomm)
464             {
465                 nstglobalcomm = ir->nstenergy;
466             }
467         }
468         else
469         {
470             /* The user has made a choice (perhaps implicitly), so we
471              * ensure that we do timely intra-simulation communication
472              * for (possibly) each of the four parts that care.
473              *
474              * TODO Does the Verlet scheme (+ DD) need any
475              * communication at nstlist steps? Is the use of nstlist
476              * here a leftover of the twin-range scheme? Can we remove
477              * nstlist when we remove the group scheme?
478              */
479             nstglobalcomm = lcd4(ir->nstcalcenergy,
480                                  ir->nstlist,
481                                  ir->etc != etcNO ? ir->nsttcouple : 0,
482                                  ir->epc != epcNO ? ir->nstpcouple : 0);
483         }
484     }
485     else
486     {
487         // Check that the user's choice of mdrun -gcom will work
488         if (ir->nstlist > 0 &&
489             nstglobalcomm > ir->nstlist && nstglobalcomm % ir->nstlist != 0)
490         {
491             nstglobalcomm = (nstglobalcomm / ir->nstlist)*ir->nstlist;
492             GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
493                     "WARNING: nstglobalcomm is larger than nstlist, but not a multiple, setting it to %d",
494                     nstglobalcomm);
495         }
496         if (ir->nstcalcenergy > 0)
497         {
498             check_nst_param(mdlog, "-gcom", nstglobalcomm,
499                             "nstcalcenergy", &ir->nstcalcenergy);
500         }
501         if (ir->etc != etcNO && ir->nsttcouple > 0)
502         {
503             check_nst_param(mdlog, "-gcom", nstglobalcomm,
504                             "nsttcouple", &ir->nsttcouple);
505         }
506         if (ir->epc != epcNO && ir->nstpcouple > 0)
507         {
508             check_nst_param(mdlog, "-gcom", nstglobalcomm,
509                             "nstpcouple", &ir->nstpcouple);
510         }
511
512         check_nst_param(mdlog, "-gcom", nstglobalcomm,
513                         "nstenergy", &ir->nstenergy);
514
515         check_nst_param(mdlog, "-gcom", nstglobalcomm,
516                         "nstlog", &ir->nstlog);
517     }
518
519     if (ir->comm_mode != ecmNO && ir->nstcomm < nstglobalcomm)
520     {
521         GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
522                 "WARNING: Changing nstcomm from %d to %d",
523                 ir->nstcomm, nstglobalcomm);
524         ir->nstcomm = nstglobalcomm;
525     }
526
527     GMX_LOG(mdlog.info).appendTextFormatted(
528             "Intra-simulation communication will occur every %d steps.\n", nstglobalcomm);
529     return nstglobalcomm;
530 }
531
532 void rerun_parallel_comm(t_commrec *cr, t_trxframe *fr,
533                          gmx_bool *bLastStep)
534 {
535     rvec    *xp, *vp;
536
537     if (MASTER(cr) && *bLastStep)
538     {
539         fr->natoms = -1;
540     }
541     xp = fr->x;
542     vp = fr->v;
543     gmx_bcast(sizeof(*fr), fr, cr);
544     fr->x = xp;
545     fr->v = vp;
546
547     *bLastStep = (fr->natoms < 0);
548
549 }
550
551 // TODO Most of this logic seems to belong in the respective modules
552 void set_state_entries(t_state *state, const t_inputrec *ir)
553 {
554     /* The entries in the state in the tpx file might not correspond
555      * with what is needed, so we correct this here.
556      */
557     state->flags = 0;
558     if (ir->efep != efepNO || ir->bExpanded)
559     {
560         state->flags |= (1<<estLAMBDA);
561         state->flags |= (1<<estFEPSTATE);
562     }
563     state->flags |= (1<<estX);
564     GMX_RELEASE_ASSERT(state->x.size() >= static_cast<unsigned int>(state->natoms), "We should start a run with an initialized state->x");
565     if (EI_DYNAMICS(ir->eI))
566     {
567         state->flags |= (1<<estV);
568     }
569
570     state->nnhpres = 0;
571     if (ir->ePBC != epbcNONE)
572     {
573         state->flags |= (1<<estBOX);
574         if (inputrecPreserveShape(ir))
575         {
576             state->flags |= (1<<estBOX_REL);
577         }
578         if ((ir->epc == epcPARRINELLORAHMAN) || (ir->epc == epcMTTK))
579         {
580             state->flags |= (1<<estBOXV);
581             state->flags |= (1<<estPRES_PREV);
582         }
583         if (inputrecNptTrotter(ir) || (inputrecNphTrotter(ir)))
584         {
585             state->nnhpres = 1;
586             state->flags  |= (1<<estNHPRES_XI);
587             state->flags  |= (1<<estNHPRES_VXI);
588             state->flags  |= (1<<estSVIR_PREV);
589             state->flags  |= (1<<estFVIR_PREV);
590             state->flags  |= (1<<estVETA);
591             state->flags  |= (1<<estVOL0);
592         }
593         if (ir->epc == epcBERENDSEN)
594         {
595             state->flags  |= (1<<estBAROS_INT);
596         }
597     }
598
599     if (ir->etc == etcNOSEHOOVER)
600     {
601         state->flags |= (1<<estNH_XI);
602         state->flags |= (1<<estNH_VXI);
603     }
604
605     if (ir->etc == etcVRESCALE || ir->etc == etcBERENDSEN)
606     {
607         state->flags |= (1<<estTHERM_INT);
608     }
609
610     init_gtc_state(state, state->ngtc, state->nnhpres, ir->opts.nhchainlength); /* allocate the space for nose-hoover chains */
611     init_ekinstate(&state->ekinstate, ir);
612
613     if (ir->bExpanded)
614     {
615         snew(state->dfhist, 1);
616         init_df_history(state->dfhist, ir->fepvals->n_lambda);
617     }
618 }