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