Turn SystemAtomIterator into proper iterator
[gromacs.git] / src / gromacs / gmxpreprocess / readir.cpp
blob3d92b2fe21fec9da248cd02e7dda6570f46e52d3
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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.
37 #include "gmxpre.h"
39 #include "readir.h"
41 #include <cctype>
42 #include <climits>
43 #include <cmath>
44 #include <cstdlib>
46 #include <algorithm>
47 #include <string>
49 #include "gromacs/awh/read-params.h"
50 #include "gromacs/fileio/readinp.h"
51 #include "gromacs/fileio/warninp.h"
52 #include "gromacs/gmxlib/chargegroup.h"
53 #include "gromacs/gmxlib/network.h"
54 #include "gromacs/gmxpreprocess/keyvaluetreemdpwriter.h"
55 #include "gromacs/gmxpreprocess/toputil.h"
56 #include "gromacs/math/functions.h"
57 #include "gromacs/math/units.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/mdlib/calc_verletbuf.h"
60 #include "gromacs/mdrunutility/mdmodules.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/mdtypes/pull-params.h"
64 #include "gromacs/options/options.h"
65 #include "gromacs/options/treesupport.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/topology/block.h"
68 #include "gromacs/topology/ifunc.h"
69 #include "gromacs/topology/index.h"
70 #include "gromacs/topology/mtop_util.h"
71 #include "gromacs/topology/symtab.h"
72 #include "gromacs/topology/topology.h"
73 #include "gromacs/utility/cstringutil.h"
74 #include "gromacs/utility/exceptions.h"
75 #include "gromacs/utility/fatalerror.h"
76 #include "gromacs/utility/filestream.h"
77 #include "gromacs/utility/gmxassert.h"
78 #include "gromacs/utility/ikeyvaluetreeerror.h"
79 #include "gromacs/utility/keyvaluetree.h"
80 #include "gromacs/utility/keyvaluetreebuilder.h"
81 #include "gromacs/utility/keyvaluetreetransform.h"
82 #include "gromacs/utility/smalloc.h"
83 #include "gromacs/utility/strconvert.h"
84 #include "gromacs/utility/stringcompare.h"
85 #include "gromacs/utility/stringutil.h"
86 #include "gromacs/utility/textwriter.h"
88 #define MAXPTR 254
89 #define NOGID 255
91 /* Resource parameters
92 * Do not change any of these until you read the instruction
93 * in readinp.h. Some cpp's do not take spaces after the backslash
94 * (like the c-shell), which will give you a very weird compiler
95 * message.
98 typedef struct t_inputrec_strings
100 char tcgrps[STRLEN], tau_t[STRLEN], ref_t[STRLEN],
101 acc[STRLEN], accgrps[STRLEN], freeze[STRLEN], frdim[STRLEN],
102 energy[STRLEN], user1[STRLEN], user2[STRLEN], vcm[STRLEN], x_compressed_groups[STRLEN],
103 couple_moltype[STRLEN], orirefitgrp[STRLEN], egptable[STRLEN], egpexcl[STRLEN],
104 wall_atomtype[STRLEN], wall_density[STRLEN], deform[STRLEN], QMMM[STRLEN],
105 imd_grp[STRLEN];
106 char fep_lambda[efptNR][STRLEN];
107 char lambda_weights[STRLEN];
108 char **pull_grp;
109 char **rot_grp;
110 char anneal[STRLEN], anneal_npoints[STRLEN],
111 anneal_time[STRLEN], anneal_temp[STRLEN];
112 char QMmethod[STRLEN], QMbasis[STRLEN], QMcharge[STRLEN], QMmult[STRLEN],
113 bSH[STRLEN], CASorbitals[STRLEN], CASelectrons[STRLEN], SAon[STRLEN],
114 SAoff[STRLEN], SAsteps[STRLEN];
116 } gmx_inputrec_strings;
118 static gmx_inputrec_strings *is = nullptr;
120 void init_inputrec_strings()
122 if (is)
124 gmx_incons("Attempted to call init_inputrec_strings before calling done_inputrec_strings. Only one inputrec (i.e. .mdp file) can be parsed at a time.");
126 snew(is, 1);
129 void done_inputrec_strings()
131 sfree(is);
132 is = nullptr;
136 enum {
137 egrptpALL, /* All particles have to be a member of a group. */
138 egrptpALL_GENREST, /* A rest group with name is generated for particles *
139 * that are not part of any group. */
140 egrptpPART, /* As egrptpALL_GENREST, but no name is generated *
141 * for the rest group. */
142 egrptpONE /* Merge all selected groups into one group, *
143 * make a rest group for the remaining particles. */
146 static const char *constraints[eshNR+1] = {
147 "none", "h-bonds", "all-bonds", "h-angles", "all-angles", nullptr
150 static const char *couple_lam[ecouplamNR+1] = {
151 "vdw-q", "vdw", "q", "none", nullptr
154 static void GetSimTemps(int ntemps, t_simtemp *simtemp, double *temperature_lambdas)
157 int i;
159 for (i = 0; i < ntemps; i++)
161 /* simple linear scaling -- allows more control */
162 if (simtemp->eSimTempScale == esimtempLINEAR)
164 simtemp->temperatures[i] = simtemp->simtemp_low + (simtemp->simtemp_high-simtemp->simtemp_low)*temperature_lambdas[i];
166 else if (simtemp->eSimTempScale == esimtempGEOMETRIC) /* should give roughly equal acceptance for constant heat capacity . . . */
168 simtemp->temperatures[i] = simtemp->simtemp_low * std::pow(simtemp->simtemp_high/simtemp->simtemp_low, static_cast<real>((1.0*i)/(ntemps-1)));
170 else if (simtemp->eSimTempScale == esimtempEXPONENTIAL)
172 simtemp->temperatures[i] = simtemp->simtemp_low + (simtemp->simtemp_high-simtemp->simtemp_low)*(std::expm1(temperature_lambdas[i])/std::expm1(1.0));
174 else
176 char errorstr[128];
177 sprintf(errorstr, "eSimTempScale=%d not defined", simtemp->eSimTempScale);
178 gmx_fatal(FARGS, "%s", errorstr);
185 static void _low_check(bool b, const char *s, warninp_t wi)
187 if (b)
189 warning_error(wi, s);
193 static void check_nst(const char *desc_nst, int nst,
194 const char *desc_p, int *p,
195 warninp_t wi)
197 char buf[STRLEN];
199 if (*p > 0 && *p % nst != 0)
201 /* Round up to the next multiple of nst */
202 *p = ((*p)/nst + 1)*nst;
203 sprintf(buf, "%s should be a multiple of %s, changing %s to %d\n",
204 desc_p, desc_nst, desc_p, *p);
205 warning(wi, buf);
209 static bool ir_NVE(const t_inputrec *ir)
211 return (EI_MD(ir->eI) && ir->etc == etcNO);
214 static int lcd(int n1, int n2)
216 int d, i;
218 d = 1;
219 for (i = 2; (i <= n1 && i <= n2); i++)
221 if (n1 % i == 0 && n2 % i == 0)
223 d = i;
227 return d;
230 static void process_interaction_modifier(const t_inputrec *ir, int *eintmod)
232 if (*eintmod == eintmodPOTSHIFT_VERLET)
234 if (ir->cutoff_scheme == ecutsVERLET)
236 *eintmod = eintmodPOTSHIFT;
238 else
240 *eintmod = eintmodNONE;
245 void check_ir(const char *mdparin, t_inputrec *ir, t_gromppopts *opts,
246 warninp_t wi)
247 /* Check internal consistency.
248 * NOTE: index groups are not set here yet, don't check things
249 * like temperature coupling group options here, but in triple_check
252 /* Strange macro: first one fills the err_buf, and then one can check
253 * the condition, which will print the message and increase the error
254 * counter.
256 #define CHECK(b) _low_check(b, err_buf, wi)
257 char err_buf[256], warn_buf[STRLEN];
258 int i, j;
259 real dt_pcoupl;
260 t_lambda *fep = ir->fepvals;
261 t_expanded *expand = ir->expandedvals;
263 set_warning_line(wi, mdparin, -1);
265 if (ir->coulombtype == eelRF_NEC_UNSUPPORTED)
267 sprintf(warn_buf, "%s electrostatics is no longer supported",
268 eel_names[eelRF_NEC_UNSUPPORTED]);
269 warning_error(wi, warn_buf);
272 /* BASIC CUT-OFF STUFF */
273 if (ir->rcoulomb < 0)
275 warning_error(wi, "rcoulomb should be >= 0");
277 if (ir->rvdw < 0)
279 warning_error(wi, "rvdw should be >= 0");
281 if (ir->rlist < 0 &&
282 !(ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0))
284 warning_error(wi, "rlist should be >= 0");
286 sprintf(err_buf, "nstlist can not be smaller than 0. (If you were trying to use the heuristic neighbour-list update scheme for efficient buffering for improved energy conservation, please use the Verlet cut-off scheme instead.)");
287 CHECK(ir->nstlist < 0);
289 process_interaction_modifier(ir, &ir->coulomb_modifier);
290 process_interaction_modifier(ir, &ir->vdw_modifier);
292 if (ir->cutoff_scheme == ecutsGROUP)
294 warning_note(wi,
295 "The group cutoff scheme is deprecated since GROMACS 5.0 and will be removed in a future "
296 "release when all interaction forms are supported for the verlet scheme. The verlet "
297 "scheme already scales better, and it is compatible with GPUs and other accelerators.");
299 if (ir->rlist > 0 && ir->rlist < ir->rcoulomb)
301 gmx_fatal(FARGS, "rcoulomb must not be greater than rlist (twin-range schemes are not supported)");
303 if (ir->rlist > 0 && ir->rlist < ir->rvdw)
305 gmx_fatal(FARGS, "rvdw must not be greater than rlist (twin-range schemes are not supported)");
308 if (ir->rlist == 0 && ir->ePBC != epbcNONE)
310 warning_error(wi, "Can not have an infinite cut-off with PBC");
314 if (ir->cutoff_scheme == ecutsVERLET)
316 real rc_max;
318 /* Normal Verlet type neighbor-list, currently only limited feature support */
319 if (inputrec2nboundeddim(ir) < 3)
321 warning_error(wi, "With Verlet lists only full pbc or pbc=xy with walls is supported");
324 // We don't (yet) have general Verlet kernels for rcoulomb!=rvdw
325 if (ir->rcoulomb != ir->rvdw)
327 // Since we have PME coulomb + LJ cut-off kernels with rcoulomb>rvdw
328 // for PME load balancing, we can support this exception.
329 bool bUsesPmeTwinRangeKernel = (EEL_PME_EWALD(ir->coulombtype) &&
330 ir->vdwtype == evdwCUT &&
331 ir->rcoulomb > ir->rvdw);
332 if (!bUsesPmeTwinRangeKernel)
334 warning_error(wi, "With Verlet lists rcoulomb!=rvdw is not supported (except for rcoulomb>rvdw with PME electrostatics)");
338 if (ir->vdwtype == evdwSHIFT || ir->vdwtype == evdwSWITCH)
340 if (ir->vdw_modifier == eintmodNONE ||
341 ir->vdw_modifier == eintmodPOTSHIFT)
343 ir->vdw_modifier = (ir->vdwtype == evdwSHIFT ? eintmodFORCESWITCH : eintmodPOTSWITCH);
345 sprintf(warn_buf, "Replacing vdwtype=%s by the equivalent combination of vdwtype=%s and vdw_modifier=%s", evdw_names[ir->vdwtype], evdw_names[evdwCUT], eintmod_names[ir->vdw_modifier]);
346 warning_note(wi, warn_buf);
348 ir->vdwtype = evdwCUT;
350 else
352 sprintf(warn_buf, "Unsupported combination of vdwtype=%s and vdw_modifier=%s", evdw_names[ir->vdwtype], eintmod_names[ir->vdw_modifier]);
353 warning_error(wi, warn_buf);
357 if (!(ir->vdwtype == evdwCUT || ir->vdwtype == evdwPME))
359 warning_error(wi, "With Verlet lists only cut-off and PME LJ interactions are supported");
361 if (!(ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype) ||
362 EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD))
364 warning_error(wi, "With Verlet lists only cut-off, reaction-field, PME and Ewald electrostatics are supported");
366 if (!(ir->coulomb_modifier == eintmodNONE ||
367 ir->coulomb_modifier == eintmodPOTSHIFT))
369 sprintf(warn_buf, "coulomb_modifier=%s is not supported with the Verlet cut-off scheme", eintmod_names[ir->coulomb_modifier]);
370 warning_error(wi, warn_buf);
373 if (EEL_USER(ir->coulombtype))
375 sprintf(warn_buf, "Coulomb type %s is not supported with the verlet scheme", eel_names[ir->coulombtype]);
376 warning_error(wi, warn_buf);
379 if (ir->nstlist <= 0)
381 warning_error(wi, "With Verlet lists nstlist should be larger than 0");
384 if (ir->nstlist < 10)
386 warning_note(wi, "With Verlet lists the optimal nstlist is >= 10, with GPUs >= 20. Note that with the Verlet scheme, nstlist has no effect on the accuracy of your simulation.");
389 rc_max = std::max(ir->rvdw, ir->rcoulomb);
391 if (ir->verletbuf_tol <= 0)
393 if (ir->verletbuf_tol == 0)
395 warning_error(wi, "Can not have Verlet buffer tolerance of exactly 0");
398 if (ir->rlist < rc_max)
400 warning_error(wi, "With verlet lists rlist can not be smaller than rvdw or rcoulomb");
403 if (ir->rlist == rc_max && ir->nstlist > 1)
405 warning_note(wi, "rlist is equal to rvdw and/or rcoulomb: there is no explicit Verlet buffer. The cluster pair list does have a buffering effect, but choosing a larger rlist might be necessary for good energy conservation.");
408 else
410 if (ir->rlist > rc_max)
412 warning_note(wi, "You have set rlist larger than the interaction cut-off, but you also have verlet-buffer-tolerance > 0. Will set rlist using verlet-buffer-tolerance.");
415 if (ir->nstlist == 1)
417 /* No buffer required */
418 ir->rlist = rc_max;
420 else
422 if (EI_DYNAMICS(ir->eI))
424 if (inputrec2nboundeddim(ir) < 3)
426 warning_error(wi, "The box volume is required for calculating rlist from the energy drift with verlet-buffer-tolerance > 0. You are using at least one unbounded dimension, so no volume can be computed. Either use a finite box, or set rlist yourself together with verlet-buffer-tolerance = -1.");
428 /* Set rlist temporarily so we can continue processing */
429 ir->rlist = rc_max;
431 else
433 /* Set the buffer to 5% of the cut-off */
434 ir->rlist = (1.0 + verlet_buffer_ratio_nodynamics)*rc_max;
440 /* GENERAL INTEGRATOR STUFF */
441 if (!EI_MD(ir->eI))
443 if (ir->etc != etcNO)
445 if (EI_RANDOM(ir->eI))
447 sprintf(warn_buf, "Setting tcoupl from '%s' to 'no'. %s handles temperature coupling implicitly. See the documentation for more information on which parameters affect temperature for %s.", etcoupl_names[ir->etc], ei_names[ir->eI], ei_names[ir->eI]);
449 else
451 sprintf(warn_buf, "Setting tcoupl from '%s' to 'no'. Temperature coupling does not apply to %s.", etcoupl_names[ir->etc], ei_names[ir->eI]);
453 warning_note(wi, warn_buf);
455 ir->etc = etcNO;
457 if (ir->eI == eiVVAK)
459 sprintf(warn_buf, "Integrator method %s is implemented primarily for validation purposes; for molecular dynamics, you should probably be using %s or %s", ei_names[eiVVAK], ei_names[eiMD], ei_names[eiVV]);
460 warning_note(wi, warn_buf);
462 if (!EI_DYNAMICS(ir->eI))
464 if (ir->epc != epcNO)
466 sprintf(warn_buf, "Setting pcoupl from '%s' to 'no'. Pressure coupling does not apply to %s.", epcoupl_names[ir->epc], ei_names[ir->eI]);
467 warning_note(wi, warn_buf);
469 ir->epc = epcNO;
471 if (EI_DYNAMICS(ir->eI))
473 if (ir->nstcalcenergy < 0)
475 ir->nstcalcenergy = ir_optimal_nstcalcenergy(ir);
476 if (ir->nstenergy != 0 && ir->nstenergy < ir->nstcalcenergy)
478 /* nstcalcenergy larger than nstener does not make sense.
479 * We ideally want nstcalcenergy=nstener.
481 if (ir->nstlist > 0)
483 ir->nstcalcenergy = lcd(ir->nstenergy, ir->nstlist);
485 else
487 ir->nstcalcenergy = ir->nstenergy;
491 else if ( (ir->nstenergy > 0 && ir->nstcalcenergy > ir->nstenergy) ||
492 (ir->efep != efepNO && ir->fepvals->nstdhdl > 0 &&
493 (ir->nstcalcenergy > ir->fepvals->nstdhdl) ) )
496 const char *nsten = "nstenergy";
497 const char *nstdh = "nstdhdl";
498 const char *min_name = nsten;
499 int min_nst = ir->nstenergy;
501 /* find the smallest of ( nstenergy, nstdhdl ) */
502 if (ir->efep != efepNO && ir->fepvals->nstdhdl > 0 &&
503 (ir->nstenergy == 0 || ir->fepvals->nstdhdl < ir->nstenergy))
505 min_nst = ir->fepvals->nstdhdl;
506 min_name = nstdh;
508 /* If the user sets nstenergy small, we should respect that */
509 sprintf(warn_buf,
510 "Setting nstcalcenergy (%d) equal to %s (%d)",
511 ir->nstcalcenergy, min_name, min_nst);
512 warning_note(wi, warn_buf);
513 ir->nstcalcenergy = min_nst;
516 if (ir->epc != epcNO)
518 if (ir->nstpcouple < 0)
520 ir->nstpcouple = ir_optimal_nstpcouple(ir);
524 if (ir->nstcalcenergy > 0)
526 if (ir->efep != efepNO)
528 /* nstdhdl should be a multiple of nstcalcenergy */
529 check_nst("nstcalcenergy", ir->nstcalcenergy,
530 "nstdhdl", &ir->fepvals->nstdhdl, wi);
532 if (ir->bExpanded)
534 /* nstexpanded should be a multiple of nstcalcenergy */
535 check_nst("nstcalcenergy", ir->nstcalcenergy,
536 "nstexpanded", &ir->expandedvals->nstexpanded, wi);
538 /* for storing exact averages nstenergy should be
539 * a multiple of nstcalcenergy
541 check_nst("nstcalcenergy", ir->nstcalcenergy,
542 "nstenergy", &ir->nstenergy, wi);
546 if (ir->nsteps == 0 && !ir->bContinuation)
548 warning_note(wi, "For a correct single-point energy evaluation with nsteps = 0, use continuation = yes to avoid constraining the input coordinates.");
551 /* LD STUFF */
552 if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
553 ir->bContinuation && ir->ld_seed != -1)
555 warning_note(wi, "You are doing a continuation with SD or BD, make sure that ld_seed is different from the previous run (using ld_seed=-1 will ensure this)");
558 /* TPI STUFF */
559 if (EI_TPI(ir->eI))
561 sprintf(err_buf, "TPI only works with pbc = %s", epbc_names[epbcXYZ]);
562 CHECK(ir->ePBC != epbcXYZ);
563 sprintf(err_buf, "TPI only works with ns = %s", ens_names[ensGRID]);
564 CHECK(ir->ns_type != ensGRID);
565 sprintf(err_buf, "with TPI nstlist should be larger than zero");
566 CHECK(ir->nstlist <= 0);
567 sprintf(err_buf, "TPI does not work with full electrostatics other than PME");
568 CHECK(EEL_FULL(ir->coulombtype) && !EEL_PME(ir->coulombtype));
569 sprintf(err_buf, "TPI does not work (yet) with the Verlet cut-off scheme");
570 CHECK(ir->cutoff_scheme == ecutsVERLET);
573 /* SHAKE / LINCS */
574 if ( (opts->nshake > 0) && (opts->bMorse) )
576 sprintf(warn_buf,
577 "Using morse bond-potentials while constraining bonds is useless");
578 warning(wi, warn_buf);
581 if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
582 ir->bContinuation && ir->ld_seed != -1)
584 warning_note(wi, "You are doing a continuation with SD or BD, make sure that ld_seed is different from the previous run (using ld_seed=-1 will ensure this)");
586 /* verify simulated tempering options */
588 if (ir->bSimTemp)
590 bool bAllTempZero = TRUE;
591 for (i = 0; i < fep->n_lambda; i++)
593 sprintf(err_buf, "Entry %d for %s must be between 0 and 1, instead is %g", i, efpt_names[efptTEMPERATURE], fep->all_lambda[efptTEMPERATURE][i]);
594 CHECK((fep->all_lambda[efptTEMPERATURE][i] < 0) || (fep->all_lambda[efptTEMPERATURE][i] > 1));
595 if (fep->all_lambda[efptTEMPERATURE][i] > 0)
597 bAllTempZero = FALSE;
600 sprintf(err_buf, "if simulated tempering is on, temperature-lambdas may not be all zero");
601 CHECK(bAllTempZero == TRUE);
603 sprintf(err_buf, "Simulated tempering is currently only compatible with md-vv");
604 CHECK(ir->eI != eiVV);
606 /* check compatability of the temperature coupling with simulated tempering */
608 if (ir->etc == etcNOSEHOOVER)
610 sprintf(warn_buf, "Nose-Hoover based temperature control such as [%s] my not be entirelyconsistent with simulated tempering", etcoupl_names[ir->etc]);
611 warning_note(wi, warn_buf);
614 /* check that the temperatures make sense */
616 sprintf(err_buf, "Higher simulated tempering temperature (%g) must be >= than the simulated tempering lower temperature (%g)", ir->simtempvals->simtemp_high, ir->simtempvals->simtemp_low);
617 CHECK(ir->simtempvals->simtemp_high <= ir->simtempvals->simtemp_low);
619 sprintf(err_buf, "Higher simulated tempering temperature (%g) must be >= zero", ir->simtempvals->simtemp_high);
620 CHECK(ir->simtempvals->simtemp_high <= 0);
622 sprintf(err_buf, "Lower simulated tempering temperature (%g) must be >= zero", ir->simtempvals->simtemp_low);
623 CHECK(ir->simtempvals->simtemp_low <= 0);
626 /* verify free energy options */
628 if (ir->efep != efepNO)
630 fep = ir->fepvals;
631 sprintf(err_buf, "The soft-core power is %d and can only be 1 or 2",
632 fep->sc_power);
633 CHECK(fep->sc_alpha != 0 && fep->sc_power != 1 && fep->sc_power != 2);
635 sprintf(err_buf, "The soft-core sc-r-power is %d and can only be 6 or 48",
636 static_cast<int>(fep->sc_r_power));
637 CHECK(fep->sc_alpha != 0 && fep->sc_r_power != 6.0 && fep->sc_r_power != 48.0);
639 sprintf(err_buf, "Can't use positive delta-lambda (%g) if initial state/lambda does not start at zero", fep->delta_lambda);
640 CHECK(fep->delta_lambda > 0 && ((fep->init_fep_state > 0) || (fep->init_lambda > 0)));
642 sprintf(err_buf, "Can't use positive delta-lambda (%g) with expanded ensemble simulations", fep->delta_lambda);
643 CHECK(fep->delta_lambda > 0 && (ir->efep == efepEXPANDED));
645 sprintf(err_buf, "Can only use expanded ensemble with md-vv (for now)");
646 CHECK(!(EI_VV(ir->eI)) && (ir->efep == efepEXPANDED));
648 sprintf(err_buf, "Free-energy not implemented for Ewald");
649 CHECK(ir->coulombtype == eelEWALD);
651 /* check validty of lambda inputs */
652 if (fep->n_lambda == 0)
654 /* Clear output in case of no states:*/
655 sprintf(err_buf, "init-lambda-state set to %d: no lambda states are defined.", fep->init_fep_state);
656 CHECK((fep->init_fep_state >= 0) && (fep->n_lambda == 0));
658 else
660 sprintf(err_buf, "initial thermodynamic state %d does not exist, only goes to %d", fep->init_fep_state, fep->n_lambda-1);
661 CHECK((fep->init_fep_state >= fep->n_lambda));
664 sprintf(err_buf, "Lambda state must be set, either with init-lambda-state or with init-lambda");
665 CHECK((fep->init_fep_state < 0) && (fep->init_lambda < 0));
667 sprintf(err_buf, "init-lambda=%g while init-lambda-state=%d. Lambda state must be set either with init-lambda-state or with init-lambda, but not both",
668 fep->init_lambda, fep->init_fep_state);
669 CHECK((fep->init_fep_state >= 0) && (fep->init_lambda >= 0));
673 if ((fep->init_lambda >= 0) && (fep->delta_lambda == 0))
675 int n_lambda_terms;
676 n_lambda_terms = 0;
677 for (i = 0; i < efptNR; i++)
679 if (fep->separate_dvdl[i])
681 n_lambda_terms++;
684 if (n_lambda_terms > 1)
686 sprintf(warn_buf, "If lambda vector states (fep-lambdas, coul-lambdas etc.) are set, don't use init-lambda to set lambda state (except for slow growth). Use init-lambda-state instead.");
687 warning(wi, warn_buf);
690 if (n_lambda_terms < 2 && fep->n_lambda > 0)
692 warning_note(wi,
693 "init-lambda is deprecated for setting lambda state (except for slow growth). Use init-lambda-state instead.");
697 for (j = 0; j < efptNR; j++)
699 for (i = 0; i < fep->n_lambda; i++)
701 sprintf(err_buf, "Entry %d for %s must be between 0 and 1, instead is %g", i, efpt_names[j], fep->all_lambda[j][i]);
702 CHECK((fep->all_lambda[j][i] < 0) || (fep->all_lambda[j][i] > 1));
706 if ((fep->sc_alpha > 0) && (!fep->bScCoul))
708 for (i = 0; i < fep->n_lambda; i++)
710 sprintf(err_buf, "For state %d, vdw-lambdas (%f) is changing with vdw softcore, while coul-lambdas (%f) is nonzero without coulomb softcore: this will lead to crashes, and is not supported.", i, fep->all_lambda[efptVDW][i],
711 fep->all_lambda[efptCOUL][i]);
712 CHECK((fep->sc_alpha > 0) &&
713 (((fep->all_lambda[efptCOUL][i] > 0.0) &&
714 (fep->all_lambda[efptCOUL][i] < 1.0)) &&
715 ((fep->all_lambda[efptVDW][i] > 0.0) &&
716 (fep->all_lambda[efptVDW][i] < 1.0))));
720 if ((fep->bScCoul) && (EEL_PME(ir->coulombtype)))
722 real sigma, lambda, r_sc;
724 sigma = 0.34;
725 /* Maximum estimate for A and B charges equal with lambda power 1 */
726 lambda = 0.5;
727 r_sc = std::pow(lambda*fep->sc_alpha*std::pow(sigma/ir->rcoulomb, fep->sc_r_power) + 1.0, 1.0/fep->sc_r_power);
728 sprintf(warn_buf, "With PME there is a minor soft core effect present at the cut-off, proportional to (LJsigma/rcoulomb)^%g. This could have a minor effect on energy conservation, but usually other effects dominate. With a common sigma value of %g nm the fraction of the particle-particle potential at the cut-off at lambda=%g is around %.1e, while ewald-rtol is %.1e.",
729 fep->sc_r_power,
730 sigma, lambda, r_sc - 1.0, ir->ewald_rtol);
731 warning_note(wi, warn_buf);
734 /* Free Energy Checks -- In an ideal world, slow growth and FEP would
735 be treated differently, but that's the next step */
737 for (i = 0; i < efptNR; i++)
739 for (j = 0; j < fep->n_lambda; j++)
741 sprintf(err_buf, "%s[%d] must be between 0 and 1", efpt_names[i], j);
742 CHECK((fep->all_lambda[i][j] < 0) || (fep->all_lambda[i][j] > 1));
747 if ((ir->bSimTemp) || (ir->efep == efepEXPANDED))
749 fep = ir->fepvals;
751 /* checking equilibration of weights inputs for validity */
753 sprintf(err_buf, "weight-equil-number-all-lambda (%d) is ignored if lmc-weights-equil is not equal to %s",
754 expand->equil_n_at_lam, elmceq_names[elmceqNUMATLAM]);
755 CHECK((expand->equil_n_at_lam > 0) && (expand->elmceq != elmceqNUMATLAM));
757 sprintf(err_buf, "weight-equil-number-samples (%d) is ignored if lmc-weights-equil is not equal to %s",
758 expand->equil_samples, elmceq_names[elmceqSAMPLES]);
759 CHECK((expand->equil_samples > 0) && (expand->elmceq != elmceqSAMPLES));
761 sprintf(err_buf, "weight-equil-number-steps (%d) is ignored if lmc-weights-equil is not equal to %s",
762 expand->equil_steps, elmceq_names[elmceqSTEPS]);
763 CHECK((expand->equil_steps > 0) && (expand->elmceq != elmceqSTEPS));
765 sprintf(err_buf, "weight-equil-wl-delta (%d) is ignored if lmc-weights-equil is not equal to %s",
766 expand->equil_samples, elmceq_names[elmceqWLDELTA]);
767 CHECK((expand->equil_wl_delta > 0) && (expand->elmceq != elmceqWLDELTA));
769 sprintf(err_buf, "weight-equil-count-ratio (%f) is ignored if lmc-weights-equil is not equal to %s",
770 expand->equil_ratio, elmceq_names[elmceqRATIO]);
771 CHECK((expand->equil_ratio > 0) && (expand->elmceq != elmceqRATIO));
773 sprintf(err_buf, "weight-equil-number-all-lambda (%d) must be a positive integer if lmc-weights-equil=%s",
774 expand->equil_n_at_lam, elmceq_names[elmceqNUMATLAM]);
775 CHECK((expand->equil_n_at_lam <= 0) && (expand->elmceq == elmceqNUMATLAM));
777 sprintf(err_buf, "weight-equil-number-samples (%d) must be a positive integer if lmc-weights-equil=%s",
778 expand->equil_samples, elmceq_names[elmceqSAMPLES]);
779 CHECK((expand->equil_samples <= 0) && (expand->elmceq == elmceqSAMPLES));
781 sprintf(err_buf, "weight-equil-number-steps (%d) must be a positive integer if lmc-weights-equil=%s",
782 expand->equil_steps, elmceq_names[elmceqSTEPS]);
783 CHECK((expand->equil_steps <= 0) && (expand->elmceq == elmceqSTEPS));
785 sprintf(err_buf, "weight-equil-wl-delta (%f) must be > 0 if lmc-weights-equil=%s",
786 expand->equil_wl_delta, elmceq_names[elmceqWLDELTA]);
787 CHECK((expand->equil_wl_delta <= 0) && (expand->elmceq == elmceqWLDELTA));
789 sprintf(err_buf, "weight-equil-count-ratio (%f) must be > 0 if lmc-weights-equil=%s",
790 expand->equil_ratio, elmceq_names[elmceqRATIO]);
791 CHECK((expand->equil_ratio <= 0) && (expand->elmceq == elmceqRATIO));
793 sprintf(err_buf, "lmc-weights-equil=%s only possible when lmc-stats = %s or lmc-stats %s",
794 elmceq_names[elmceqWLDELTA], elamstats_names[elamstatsWL], elamstats_names[elamstatsWWL]);
795 CHECK((expand->elmceq == elmceqWLDELTA) && (!EWL(expand->elamstats)));
797 sprintf(err_buf, "lmc-repeats (%d) must be greater than 0", expand->lmc_repeats);
798 CHECK((expand->lmc_repeats <= 0));
799 sprintf(err_buf, "minimum-var-min (%d) must be greater than 0", expand->minvarmin);
800 CHECK((expand->minvarmin <= 0));
801 sprintf(err_buf, "weight-c-range (%d) must be greater or equal to 0", expand->c_range);
802 CHECK((expand->c_range < 0));
803 sprintf(err_buf, "init-lambda-state (%d) must be zero if lmc-forced-nstart (%d)> 0 and lmc-move != 'no'",
804 fep->init_fep_state, expand->lmc_forced_nstart);
805 CHECK((fep->init_fep_state != 0) && (expand->lmc_forced_nstart > 0) && (expand->elmcmove != elmcmoveNO));
806 sprintf(err_buf, "lmc-forced-nstart (%d) must not be negative", expand->lmc_forced_nstart);
807 CHECK((expand->lmc_forced_nstart < 0));
808 sprintf(err_buf, "init-lambda-state (%d) must be in the interval [0,number of lambdas)", fep->init_fep_state);
809 CHECK((fep->init_fep_state < 0) || (fep->init_fep_state >= fep->n_lambda));
811 sprintf(err_buf, "init-wl-delta (%f) must be greater than or equal to 0", expand->init_wl_delta);
812 CHECK((expand->init_wl_delta < 0));
813 sprintf(err_buf, "wl-ratio (%f) must be between 0 and 1", expand->wl_ratio);
814 CHECK((expand->wl_ratio <= 0) || (expand->wl_ratio >= 1));
815 sprintf(err_buf, "wl-scale (%f) must be between 0 and 1", expand->wl_scale);
816 CHECK((expand->wl_scale <= 0) || (expand->wl_scale >= 1));
818 /* if there is no temperature control, we need to specify an MC temperature */
819 if (!integratorHasReferenceTemperature(ir) && (expand->elmcmove != elmcmoveNO) && (expand->mc_temp <= 0.0))
821 sprintf(err_buf, "If there is no temperature control, and lmc-mcmove!='no', mc_temp must be set to a positive number");
822 warning_error(wi, err_buf);
824 if (expand->nstTij > 0)
826 sprintf(err_buf, "nstlog must be non-zero");
827 CHECK(ir->nstlog == 0);
828 sprintf(err_buf, "nst-transition-matrix (%d) must be an integer multiple of nstlog (%d)",
829 expand->nstTij, ir->nstlog);
830 CHECK((expand->nstTij % ir->nstlog) != 0);
834 /* PBC/WALLS */
835 sprintf(err_buf, "walls only work with pbc=%s", epbc_names[epbcXY]);
836 CHECK(ir->nwall && ir->ePBC != epbcXY);
838 /* VACUUM STUFF */
839 if (ir->ePBC != epbcXYZ && ir->nwall != 2)
841 if (ir->ePBC == epbcNONE)
843 if (ir->epc != epcNO)
845 warning(wi, "Turning off pressure coupling for vacuum system");
846 ir->epc = epcNO;
849 else
851 sprintf(err_buf, "Can not have pressure coupling with pbc=%s",
852 epbc_names[ir->ePBC]);
853 CHECK(ir->epc != epcNO);
855 sprintf(err_buf, "Can not have Ewald with pbc=%s", epbc_names[ir->ePBC]);
856 CHECK(EEL_FULL(ir->coulombtype));
858 sprintf(err_buf, "Can not have dispersion correction with pbc=%s",
859 epbc_names[ir->ePBC]);
860 CHECK(ir->eDispCorr != edispcNO);
863 if (ir->rlist == 0.0)
865 sprintf(err_buf, "can only have neighborlist cut-off zero (=infinite)\n"
866 "with coulombtype = %s or coulombtype = %s\n"
867 "without periodic boundary conditions (pbc = %s) and\n"
868 "rcoulomb and rvdw set to zero",
869 eel_names[eelCUT], eel_names[eelUSER], epbc_names[epbcNONE]);
870 CHECK(((ir->coulombtype != eelCUT) && (ir->coulombtype != eelUSER)) ||
871 (ir->ePBC != epbcNONE) ||
872 (ir->rcoulomb != 0.0) || (ir->rvdw != 0.0));
874 if (ir->nstlist > 0)
876 warning_note(wi, "Simulating without cut-offs can be (slightly) faster with nstlist=0, nstype=simple and only one MPI rank");
880 /* COMM STUFF */
881 if (ir->nstcomm == 0)
883 ir->comm_mode = ecmNO;
885 if (ir->comm_mode != ecmNO)
887 if (ir->nstcomm < 0)
889 warning(wi, "If you want to remove the rotation around the center of mass, you should set comm_mode = Angular instead of setting nstcomm < 0. nstcomm is modified to its absolute value");
890 ir->nstcomm = abs(ir->nstcomm);
893 if (ir->nstcalcenergy > 0 && ir->nstcomm < ir->nstcalcenergy)
895 warning_note(wi, "nstcomm < nstcalcenergy defeats the purpose of nstcalcenergy, setting nstcomm to nstcalcenergy");
896 ir->nstcomm = ir->nstcalcenergy;
899 if (ir->comm_mode == ecmANGULAR)
901 sprintf(err_buf, "Can not remove the rotation around the center of mass with periodic molecules");
902 CHECK(ir->bPeriodicMols);
903 if (ir->ePBC != epbcNONE)
905 warning(wi, "Removing the rotation around the center of mass in a periodic system, this can lead to artifacts. Only use this on a single (cluster of) molecules. This cluster should not cross periodic boundaries.");
910 if (EI_STATE_VELOCITY(ir->eI) && !EI_SD(ir->eI) && ir->ePBC == epbcNONE && ir->comm_mode != ecmANGULAR)
912 sprintf(warn_buf, "Tumbling and flying ice-cubes: We are not removing rotation around center of mass in a non-periodic system. You should probably set comm_mode = ANGULAR or use integrator = %s.", ei_names[eiSD1]);
913 warning_note(wi, warn_buf);
916 /* TEMPERATURE COUPLING */
917 if (ir->etc == etcYES)
919 ir->etc = etcBERENDSEN;
920 warning_note(wi, "Old option for temperature coupling given: "
921 "changing \"yes\" to \"Berendsen\"\n");
924 if ((ir->etc == etcNOSEHOOVER) || (ir->epc == epcMTTK))
926 if (ir->opts.nhchainlength < 1)
928 sprintf(warn_buf, "number of Nose-Hoover chains (currently %d) cannot be less than 1,reset to 1\n", ir->opts.nhchainlength);
929 ir->opts.nhchainlength = 1;
930 warning(wi, warn_buf);
933 if (ir->etc == etcNOSEHOOVER && !EI_VV(ir->eI) && ir->opts.nhchainlength > 1)
935 warning_note(wi, "leapfrog does not yet support Nose-Hoover chains, nhchainlength reset to 1");
936 ir->opts.nhchainlength = 1;
939 else
941 ir->opts.nhchainlength = 0;
944 if (ir->eI == eiVVAK)
946 sprintf(err_buf, "%s implemented primarily for validation, and requires nsttcouple = 1 and nstpcouple = 1.",
947 ei_names[eiVVAK]);
948 CHECK((ir->nsttcouple != 1) || (ir->nstpcouple != 1));
951 if (ETC_ANDERSEN(ir->etc))
953 sprintf(err_buf, "%s temperature control not supported for integrator %s.", etcoupl_names[ir->etc], ei_names[ir->eI]);
954 CHECK(!(EI_VV(ir->eI)));
956 if (ir->nstcomm > 0 && (ir->etc == etcANDERSEN))
958 sprintf(warn_buf, "Center of mass removal not necessary for %s. All velocities of coupled groups are rerandomized periodically, so flying ice cube errors will not occur.", etcoupl_names[ir->etc]);
959 warning_note(wi, warn_buf);
962 sprintf(err_buf, "nstcomm must be 1, not %d for %s, as velocities of atoms in coupled groups are randomized every time step", ir->nstcomm, etcoupl_names[ir->etc]);
963 CHECK(ir->nstcomm > 1 && (ir->etc == etcANDERSEN));
966 if (ir->etc == etcBERENDSEN)
968 sprintf(warn_buf, "The %s thermostat does not generate the correct kinetic energy distribution. You might want to consider using the %s thermostat.",
969 ETCOUPLTYPE(ir->etc), ETCOUPLTYPE(etcVRESCALE));
970 warning_note(wi, warn_buf);
973 if ((ir->etc == etcNOSEHOOVER || ETC_ANDERSEN(ir->etc))
974 && ir->epc == epcBERENDSEN)
976 sprintf(warn_buf, "Using Berendsen pressure coupling invalidates the "
977 "true ensemble for the thermostat");
978 warning(wi, warn_buf);
981 /* PRESSURE COUPLING */
982 if (ir->epc == epcISOTROPIC)
984 ir->epc = epcBERENDSEN;
985 warning_note(wi, "Old option for pressure coupling given: "
986 "changing \"Isotropic\" to \"Berendsen\"\n");
989 if (ir->epc != epcNO)
991 dt_pcoupl = ir->nstpcouple*ir->delta_t;
993 sprintf(err_buf, "tau-p must be > 0 instead of %g\n", ir->tau_p);
994 CHECK(ir->tau_p <= 0);
996 if (ir->tau_p/dt_pcoupl < pcouple_min_integration_steps(ir->epc) - 10*GMX_REAL_EPS)
998 sprintf(warn_buf, "For proper integration of the %s barostat, tau-p (%g) should be at least %d times larger than nstpcouple*dt (%g)",
999 EPCOUPLTYPE(ir->epc), ir->tau_p, pcouple_min_integration_steps(ir->epc), dt_pcoupl);
1000 warning(wi, warn_buf);
1003 sprintf(err_buf, "compressibility must be > 0 when using pressure"
1004 " coupling %s\n", EPCOUPLTYPE(ir->epc));
1005 CHECK(ir->compress[XX][XX] < 0 || ir->compress[YY][YY] < 0 ||
1006 ir->compress[ZZ][ZZ] < 0 ||
1007 (trace(ir->compress) == 0 && ir->compress[YY][XX] <= 0 &&
1008 ir->compress[ZZ][XX] <= 0 && ir->compress[ZZ][YY] <= 0));
1010 if (epcPARRINELLORAHMAN == ir->epc && opts->bGenVel)
1012 sprintf(warn_buf,
1013 "You are generating velocities so I am assuming you "
1014 "are equilibrating a system. You are using "
1015 "%s pressure coupling, but this can be "
1016 "unstable for equilibration. If your system crashes, try "
1017 "equilibrating first with Berendsen pressure coupling. If "
1018 "you are not equilibrating the system, you can probably "
1019 "ignore this warning.",
1020 epcoupl_names[ir->epc]);
1021 warning(wi, warn_buf);
1025 if (EI_VV(ir->eI))
1027 if (ir->epc > epcNO)
1029 if ((ir->epc != epcBERENDSEN) && (ir->epc != epcMTTK))
1031 warning_error(wi, "for md-vv and md-vv-avek, can only use Berendsen and Martyna-Tuckerman-Tobias-Klein (MTTK) equations for pressure control; MTTK is equivalent to Parrinello-Rahman.");
1035 else
1037 if (ir->epc == epcMTTK)
1039 warning_error(wi, "MTTK pressure coupling requires a Velocity-verlet integrator");
1043 /* ELECTROSTATICS */
1044 /* More checks are in triple check (grompp.c) */
1046 if (ir->coulombtype == eelSWITCH)
1048 sprintf(warn_buf, "coulombtype = %s is only for testing purposes and can lead to serious "
1049 "artifacts, advice: use coulombtype = %s",
1050 eel_names[ir->coulombtype],
1051 eel_names[eelRF_ZERO]);
1052 warning(wi, warn_buf);
1055 if (EEL_RF(ir->coulombtype) && ir->epsilon_rf == 1 && ir->epsilon_r != 1)
1057 sprintf(warn_buf, "epsilon-r = %g and epsilon-rf = 1 with reaction field, proceeding assuming old format and exchanging epsilon-r and epsilon-rf", ir->epsilon_r);
1058 warning(wi, warn_buf);
1059 ir->epsilon_rf = ir->epsilon_r;
1060 ir->epsilon_r = 1.0;
1063 if (ir->epsilon_r == 0)
1065 sprintf(err_buf,
1066 "It is pointless to use long-range electrostatics with infinite relative permittivity."
1067 "Since you are effectively turning of electrostatics, a plain cutoff will be much faster.");
1068 CHECK(EEL_FULL(ir->coulombtype));
1071 if (getenv("GMX_DO_GALACTIC_DYNAMICS") == nullptr)
1073 sprintf(err_buf, "epsilon-r must be >= 0 instead of %g\n", ir->epsilon_r);
1074 CHECK(ir->epsilon_r < 0);
1077 if (EEL_RF(ir->coulombtype))
1079 /* reaction field (at the cut-off) */
1081 if (ir->coulombtype == eelRF_ZERO && ir->epsilon_rf != 0)
1083 sprintf(warn_buf, "With coulombtype = %s, epsilon-rf must be 0, assuming you meant epsilon_rf=0",
1084 eel_names[ir->coulombtype]);
1085 warning(wi, warn_buf);
1086 ir->epsilon_rf = 0.0;
1089 sprintf(err_buf, "epsilon-rf must be >= epsilon-r");
1090 CHECK((ir->epsilon_rf < ir->epsilon_r && ir->epsilon_rf != 0) ||
1091 (ir->epsilon_r == 0));
1092 if (ir->epsilon_rf == ir->epsilon_r)
1094 sprintf(warn_buf, "Using epsilon-rf = epsilon-r with %s does not make sense",
1095 eel_names[ir->coulombtype]);
1096 warning(wi, warn_buf);
1099 /* Allow rlist>rcoulomb for tabulated long range stuff. This just
1100 * means the interaction is zero outside rcoulomb, but it helps to
1101 * provide accurate energy conservation.
1103 if (ir_coulomb_might_be_zero_at_cutoff(ir))
1105 if (ir_coulomb_switched(ir))
1107 sprintf(err_buf,
1108 "With coulombtype = %s rcoulomb_switch must be < rcoulomb. Or, better: Use the potential modifier options!",
1109 eel_names[ir->coulombtype]);
1110 CHECK(ir->rcoulomb_switch >= ir->rcoulomb);
1113 else if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
1115 if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE)
1117 sprintf(err_buf, "With coulombtype = %s, rcoulomb should be >= rlist unless you use a potential modifier",
1118 eel_names[ir->coulombtype]);
1119 CHECK(ir->rlist > ir->rcoulomb);
1123 if (ir->coulombtype == eelSWITCH || ir->coulombtype == eelSHIFT)
1125 sprintf(err_buf,
1126 "Explicit switch/shift coulomb interactions cannot be used in combination with a secondary coulomb-modifier.");
1127 CHECK( ir->coulomb_modifier != eintmodNONE);
1129 if (ir->vdwtype == evdwSWITCH || ir->vdwtype == evdwSHIFT)
1131 sprintf(err_buf,
1132 "Explicit switch/shift vdw interactions cannot be used in combination with a secondary vdw-modifier.");
1133 CHECK( ir->vdw_modifier != eintmodNONE);
1136 if (ir->coulombtype == eelSWITCH || ir->coulombtype == eelSHIFT ||
1137 ir->vdwtype == evdwSWITCH || ir->vdwtype == evdwSHIFT)
1139 sprintf(warn_buf,
1140 "The switch/shift interaction settings are just for compatibility; you will get better "
1141 "performance from applying potential modifiers to your interactions!\n");
1142 warning_note(wi, warn_buf);
1145 if (ir->coulombtype == eelPMESWITCH || ir->coulomb_modifier == eintmodPOTSWITCH)
1147 if (ir->rcoulomb_switch/ir->rcoulomb < 0.9499)
1149 real percentage = 100*(ir->rcoulomb-ir->rcoulomb_switch)/ir->rcoulomb;
1150 sprintf(warn_buf, "The switching range should be 5%% or less (currently %.2f%% using a switching range of %4f-%4f) for accurate electrostatic energies, energy conservation will be good regardless, since ewald_rtol = %g.",
1151 percentage, ir->rcoulomb_switch, ir->rcoulomb, ir->ewald_rtol);
1152 warning(wi, warn_buf);
1156 if (ir->vdwtype == evdwSWITCH || ir->vdw_modifier == eintmodPOTSWITCH)
1158 if (ir->rvdw_switch == 0)
1160 sprintf(warn_buf, "rvdw-switch is equal 0 even though you are using a switched Lennard-Jones potential. This suggests it was not set in the mdp, which can lead to large energy errors. In GROMACS, 0.05 to 0.1 nm is often a reasonable vdw switching range.");
1161 warning(wi, warn_buf);
1165 if (EEL_FULL(ir->coulombtype))
1167 if (ir->coulombtype == eelPMESWITCH || ir->coulombtype == eelPMEUSER ||
1168 ir->coulombtype == eelPMEUSERSWITCH)
1170 sprintf(err_buf, "With coulombtype = %s, rcoulomb must be <= rlist",
1171 eel_names[ir->coulombtype]);
1172 CHECK(ir->rcoulomb > ir->rlist);
1174 else if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE)
1176 if (ir->coulombtype == eelPME || ir->coulombtype == eelP3M_AD)
1178 sprintf(err_buf,
1179 "With coulombtype = %s (without modifier), rcoulomb must be equal to rlist.\n"
1180 "For optimal energy conservation,consider using\n"
1181 "a potential modifier.", eel_names[ir->coulombtype]);
1182 CHECK(ir->rcoulomb != ir->rlist);
1187 if (EEL_PME(ir->coulombtype) || EVDW_PME(ir->vdwtype))
1189 // TODO: Move these checks into the ewald module with the options class
1190 int orderMin = 3;
1191 int orderMax = (ir->coulombtype == eelP3M_AD ? 8 : 12);
1193 if (ir->pme_order < orderMin || ir->pme_order > orderMax)
1195 sprintf(warn_buf, "With coulombtype = %s, you should have %d <= pme-order <= %d", eel_names[ir->coulombtype], orderMin, orderMax);
1196 warning_error(wi, warn_buf);
1200 if (ir->nwall == 2 && EEL_FULL(ir->coulombtype))
1202 if (ir->ewald_geometry == eewg3D)
1204 sprintf(warn_buf, "With pbc=%s you should use ewald-geometry=%s",
1205 epbc_names[ir->ePBC], eewg_names[eewg3DC]);
1206 warning(wi, warn_buf);
1208 /* This check avoids extra pbc coding for exclusion corrections */
1209 sprintf(err_buf, "wall-ewald-zfac should be >= 2");
1210 CHECK(ir->wall_ewald_zfac < 2);
1212 if ((ir->ewald_geometry == eewg3DC) && (ir->ePBC != epbcXY) &&
1213 EEL_FULL(ir->coulombtype))
1215 sprintf(warn_buf, "With %s and ewald_geometry = %s you should use pbc = %s",
1216 eel_names[ir->coulombtype], eewg_names[eewg3DC], epbc_names[epbcXY]);
1217 warning(wi, warn_buf);
1219 if ((ir->epsilon_surface != 0) && EEL_FULL(ir->coulombtype))
1221 if (ir->cutoff_scheme == ecutsVERLET)
1223 sprintf(warn_buf, "Since molecules/charge groups are broken using the Verlet scheme, you can not use a dipole correction to the %s electrostatics.",
1224 eel_names[ir->coulombtype]);
1225 warning(wi, warn_buf);
1227 else
1229 sprintf(warn_buf, "Dipole corrections to %s electrostatics only work if all charge groups that can cross PBC boundaries are dipoles. If this is not the case set epsilon_surface to 0",
1230 eel_names[ir->coulombtype]);
1231 warning_note(wi, warn_buf);
1235 if (ir_vdw_switched(ir))
1237 sprintf(err_buf, "With switched vdw forces or potentials, rvdw-switch must be < rvdw");
1238 CHECK(ir->rvdw_switch >= ir->rvdw);
1240 if (ir->rvdw_switch < 0.5*ir->rvdw)
1242 sprintf(warn_buf, "You are applying a switch function to vdw forces or potentials from %g to %g nm, which is more than half the interaction range, whereas switch functions are intended to act only close to the cut-off.",
1243 ir->rvdw_switch, ir->rvdw);
1244 warning_note(wi, warn_buf);
1247 else if (ir->vdwtype == evdwCUT || ir->vdwtype == evdwPME)
1249 if (ir->cutoff_scheme == ecutsGROUP && ir->vdw_modifier == eintmodNONE)
1251 sprintf(err_buf, "With vdwtype = %s, rvdw must be >= rlist unless you use a potential modifier", evdw_names[ir->vdwtype]);
1252 CHECK(ir->rlist > ir->rvdw);
1256 if (ir->vdwtype == evdwPME)
1258 if (!(ir->vdw_modifier == eintmodNONE || ir->vdw_modifier == eintmodPOTSHIFT))
1260 sprintf(err_buf, "With vdwtype = %s, the only supported modifiers are %s and %s",
1261 evdw_names[ir->vdwtype],
1262 eintmod_names[eintmodPOTSHIFT],
1263 eintmod_names[eintmodNONE]);
1264 warning_error(wi, err_buf);
1268 if (ir->cutoff_scheme == ecutsGROUP)
1270 if (((ir->coulomb_modifier != eintmodNONE && ir->rcoulomb == ir->rlist) ||
1271 (ir->vdw_modifier != eintmodNONE && ir->rvdw == ir->rlist)))
1273 warning_note(wi, "With exact cut-offs, rlist should be "
1274 "larger than rcoulomb and rvdw, so that there "
1275 "is a buffer region for particle motion "
1276 "between neighborsearch steps");
1279 if (ir_coulomb_is_zero_at_cutoff(ir) && ir->rlist <= ir->rcoulomb)
1281 sprintf(warn_buf, "For energy conservation with switch/shift potentials, rlist should be 0.1 to 0.3 nm larger than rcoulomb.");
1282 warning_note(wi, warn_buf);
1284 if (ir_vdw_switched(ir) && (ir->rlist <= ir->rvdw))
1286 sprintf(warn_buf, "For energy conservation with switch/shift potentials, rlist should be 0.1 to 0.3 nm larger than rvdw.");
1287 warning_note(wi, warn_buf);
1291 if (ir->vdwtype == evdwUSER && ir->eDispCorr != edispcNO)
1293 warning_note(wi, "You have selected user tables with dispersion correction, the dispersion will be corrected to -C6/r^6 beyond rvdw_switch (the tabulated interaction between rvdw_switch and rvdw will not be double counted). Make sure that you really want dispersion correction to -C6/r^6.");
1296 if (ir->eI == eiLBFGS && (ir->coulombtype == eelCUT || ir->vdwtype == evdwCUT)
1297 && ir->rvdw != 0)
1299 warning(wi, "For efficient BFGS minimization, use switch/shift/pme instead of cut-off.");
1302 if (ir->eI == eiLBFGS && ir->nbfgscorr <= 0)
1304 warning(wi, "Using L-BFGS with nbfgscorr<=0 just gets you steepest descent.");
1307 /* ENERGY CONSERVATION */
1308 if (ir_NVE(ir) && ir->cutoff_scheme == ecutsGROUP)
1310 if (!ir_vdw_might_be_zero_at_cutoff(ir) && ir->rvdw > 0 && ir->vdw_modifier == eintmodNONE)
1312 sprintf(warn_buf, "You are using a cut-off for VdW interactions with NVE, for good energy conservation use vdwtype = %s (possibly with DispCorr)",
1313 evdw_names[evdwSHIFT]);
1314 warning_note(wi, warn_buf);
1316 if (!ir_coulomb_might_be_zero_at_cutoff(ir) && ir->rcoulomb > 0)
1318 sprintf(warn_buf, "You are using a cut-off for electrostatics with NVE, for good energy conservation use coulombtype = %s or %s",
1319 eel_names[eelPMESWITCH], eel_names[eelRF_ZERO]);
1320 warning_note(wi, warn_buf);
1324 /* IMPLICIT SOLVENT */
1325 if (ir->coulombtype == eelGB_NOTUSED)
1327 sprintf(warn_buf, "Invalid option %s for coulombtype",
1328 eel_names[ir->coulombtype]);
1329 warning_error(wi, warn_buf);
1332 if (ir->bQMMM)
1334 if (ir->cutoff_scheme != ecutsGROUP)
1336 warning_error(wi, "QMMM is currently only supported with cutoff-scheme=group");
1338 if (!EI_DYNAMICS(ir->eI))
1340 char buf[STRLEN];
1341 sprintf(buf, "QMMM is only supported with dynamics, not with integrator %s", ei_names[ir->eI]);
1342 warning_error(wi, buf);
1346 if (ir->bAdress)
1348 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
1352 /* interpret a number of doubles from a string and put them in an array,
1353 after allocating space for them.
1354 str = the input string
1355 n = the (pre-allocated) number of doubles read
1356 r = the output array of doubles. */
1357 static void parse_n_real(char *str, int *n, real **r, warninp_t wi)
1359 auto values = gmx::splitString(str);
1360 *n = values.size();
1362 snew(*r, *n);
1363 for (int i = 0; i < *n; i++)
1367 (*r)[i] = gmx::fromString<real>(values[i]);
1369 catch (gmx::GromacsException &)
1371 warning_error(wi, "Invalid value " + values[i] + " in string in mdp file. Expected a real number.");
1377 static void do_fep_params(t_inputrec *ir, char fep_lambda[][STRLEN], char weights[STRLEN], warninp_t wi)
1380 int i, j, max_n_lambda, nweights, nfep[efptNR];
1381 t_lambda *fep = ir->fepvals;
1382 t_expanded *expand = ir->expandedvals;
1383 real **count_fep_lambdas;
1384 bool bOneLambda = TRUE;
1386 snew(count_fep_lambdas, efptNR);
1388 /* FEP input processing */
1389 /* first, identify the number of lambda values for each type.
1390 All that are nonzero must have the same number */
1392 for (i = 0; i < efptNR; i++)
1394 parse_n_real(fep_lambda[i], &(nfep[i]), &(count_fep_lambdas[i]), wi);
1397 /* now, determine the number of components. All must be either zero, or equal. */
1399 max_n_lambda = 0;
1400 for (i = 0; i < efptNR; i++)
1402 if (nfep[i] > max_n_lambda)
1404 max_n_lambda = nfep[i]; /* here's a nonzero one. All of them
1405 must have the same number if its not zero.*/
1406 break;
1410 for (i = 0; i < efptNR; i++)
1412 if (nfep[i] == 0)
1414 ir->fepvals->separate_dvdl[i] = FALSE;
1416 else if (nfep[i] == max_n_lambda)
1418 if (i != efptTEMPERATURE) /* we treat this differently -- not really a reason to compute the derivative with
1419 respect to the temperature currently */
1421 ir->fepvals->separate_dvdl[i] = TRUE;
1424 else
1426 gmx_fatal(FARGS, "Number of lambdas (%d) for FEP type %s not equal to number of other types (%d)",
1427 nfep[i], efpt_names[i], max_n_lambda);
1430 /* we don't print out dhdl if the temperature is changing, since we can't correctly define dhdl in this case */
1431 ir->fepvals->separate_dvdl[efptTEMPERATURE] = FALSE;
1433 /* the number of lambdas is the number we've read in, which is either zero
1434 or the same for all */
1435 fep->n_lambda = max_n_lambda;
1437 /* allocate space for the array of lambda values */
1438 snew(fep->all_lambda, efptNR);
1439 /* if init_lambda is defined, we need to set lambda */
1440 if ((fep->init_lambda > 0) && (fep->n_lambda == 0))
1442 ir->fepvals->separate_dvdl[efptFEP] = TRUE;
1444 /* otherwise allocate the space for all of the lambdas, and transfer the data */
1445 for (i = 0; i < efptNR; i++)
1447 snew(fep->all_lambda[i], fep->n_lambda);
1448 if (nfep[i] > 0) /* if it's zero, then the count_fep_lambda arrays
1449 are zero */
1451 for (j = 0; j < fep->n_lambda; j++)
1453 fep->all_lambda[i][j] = static_cast<double>(count_fep_lambdas[i][j]);
1455 sfree(count_fep_lambdas[i]);
1458 sfree(count_fep_lambdas);
1460 /* "fep-vals" is either zero or the full number. If zero, we'll need to define fep-lambdas for internal
1461 bookkeeping -- for now, init_lambda */
1463 if ((nfep[efptFEP] == 0) && (fep->init_lambda >= 0))
1465 for (i = 0; i < fep->n_lambda; i++)
1467 fep->all_lambda[efptFEP][i] = fep->init_lambda;
1471 /* check to see if only a single component lambda is defined, and soft core is defined.
1472 In this case, turn on coulomb soft core */
1474 if (max_n_lambda == 0)
1476 bOneLambda = TRUE;
1478 else
1480 for (i = 0; i < efptNR; i++)
1482 if ((nfep[i] != 0) && (i != efptFEP))
1484 bOneLambda = FALSE;
1488 if ((bOneLambda) && (fep->sc_alpha > 0))
1490 fep->bScCoul = TRUE;
1493 /* Fill in the others with the efptFEP if they are not explicitly
1494 specified (i.e. nfep[i] == 0). This means if fep is not defined,
1495 they are all zero. */
1497 for (i = 0; i < efptNR; i++)
1499 if ((nfep[i] == 0) && (i != efptFEP))
1501 for (j = 0; j < fep->n_lambda; j++)
1503 fep->all_lambda[i][j] = fep->all_lambda[efptFEP][j];
1509 /* make it easier if sc_r_power = 48 by increasing it to the 4th power, to be in the right scale. */
1510 if (fep->sc_r_power == 48)
1512 if (fep->sc_alpha > 0.1)
1514 gmx_fatal(FARGS, "sc_alpha (%f) for sc_r_power = 48 should usually be between 0.001 and 0.004", fep->sc_alpha);
1518 /* now read in the weights */
1519 parse_n_real(weights, &nweights, &(expand->init_lambda_weights), wi);
1520 if (nweights == 0)
1522 snew(expand->init_lambda_weights, fep->n_lambda); /* initialize to zero */
1524 else if (nweights != fep->n_lambda)
1526 gmx_fatal(FARGS, "Number of weights (%d) is not equal to number of lambda values (%d)",
1527 nweights, fep->n_lambda);
1529 if ((expand->nstexpanded < 0) && (ir->efep != efepNO))
1531 expand->nstexpanded = fep->nstdhdl;
1532 /* if you don't specify nstexpanded when doing expanded ensemble free energy calcs, it is set to nstdhdl */
1534 if ((expand->nstexpanded < 0) && ir->bSimTemp)
1536 expand->nstexpanded = 2*static_cast<int>(ir->opts.tau_t[0]/ir->delta_t);
1537 /* if you don't specify nstexpanded when doing expanded ensemble simulated tempering, it is set to
1538 2*tau_t just to be careful so it's not to frequent */
1543 static void do_simtemp_params(t_inputrec *ir)
1546 snew(ir->simtempvals->temperatures, ir->fepvals->n_lambda);
1547 GetSimTemps(ir->fepvals->n_lambda, ir->simtempvals, ir->fepvals->all_lambda[efptTEMPERATURE]);
1550 static void
1551 convertYesNos(warninp_t /*wi*/, gmx::ArrayRef<const std::string> inputs, const char * /*name*/, gmx_bool *outputs)
1553 int i = 0;
1554 for (const auto &input : inputs)
1556 outputs[i] = (gmx_strncasecmp(input.c_str(), "Y", 1) == 0);
1557 ++i;
1561 template <typename T> void
1562 convertInts(warninp_t wi, gmx::ArrayRef<const std::string> inputs, const char *name, T *outputs)
1564 int i = 0;
1565 for (const auto &input : inputs)
1569 outputs[i] = gmx::fromStdString<T>(input);
1571 catch (gmx::GromacsException &)
1573 auto message = gmx::formatString("Invalid value for mdp option %s. %s should only consist of integers separated by spaces.",
1574 name, name);
1575 warning_error(wi, message);
1577 ++i;
1581 static void
1582 convertReals(warninp_t wi, gmx::ArrayRef<const std::string> inputs, const char *name, real *outputs)
1584 int i = 0;
1585 for (const auto &input : inputs)
1589 outputs[i] = gmx::fromString<real>(input);
1591 catch (gmx::GromacsException &)
1593 auto message = gmx::formatString("Invalid value for mdp option %s. %s should only consist of real numbers separated by spaces.",
1594 name, name);
1595 warning_error(wi, message);
1597 ++i;
1601 static void
1602 convertRvecs(warninp_t wi, gmx::ArrayRef<const std::string> inputs, const char *name, rvec *outputs)
1604 int i = 0, d = 0;
1605 for (const auto &input : inputs)
1609 outputs[i][d] = gmx::fromString<real>(input);
1611 catch (gmx::GromacsException &)
1613 auto message = gmx::formatString("Invalid value for mdp option %s. %s should only consist of real numbers separated by spaces.",
1614 name, name);
1615 warning_error(wi, message);
1617 ++d;
1618 if (d == DIM)
1620 d = 0;
1621 ++i;
1626 static void do_wall_params(t_inputrec *ir,
1627 char *wall_atomtype, char *wall_density,
1628 t_gromppopts *opts,
1629 warninp_t wi)
1631 opts->wall_atomtype[0] = nullptr;
1632 opts->wall_atomtype[1] = nullptr;
1634 ir->wall_atomtype[0] = -1;
1635 ir->wall_atomtype[1] = -1;
1636 ir->wall_density[0] = 0;
1637 ir->wall_density[1] = 0;
1639 if (ir->nwall > 0)
1641 auto wallAtomTypes = gmx::splitString(wall_atomtype);
1642 if (wallAtomTypes.size() != size_t(ir->nwall))
1644 gmx_fatal(FARGS, "Expected %d elements for wall_atomtype, found %zu",
1645 ir->nwall, wallAtomTypes.size());
1647 for (int i = 0; i < ir->nwall; i++)
1649 opts->wall_atomtype[i] = gmx_strdup(wallAtomTypes[i].c_str());
1652 if (ir->wall_type == ewt93 || ir->wall_type == ewt104)
1654 auto wallDensity = gmx::splitString(wall_density);
1655 if (wallDensity.size() != size_t(ir->nwall))
1657 gmx_fatal(FARGS, "Expected %d elements for wall-density, found %zu", ir->nwall, wallDensity.size());
1659 convertReals(wi, wallDensity, "wall-density", ir->wall_density);
1660 for (int i = 0; i < ir->nwall; i++)
1662 if (ir->wall_density[i] <= 0)
1664 gmx_fatal(FARGS, "wall-density[%d] = %f\n", i, ir->wall_density[i]);
1671 static void add_wall_energrps(gmx_groups_t *groups, int nwall, t_symtab *symtab)
1673 int i;
1674 t_grps *grps;
1675 char str[STRLEN];
1677 if (nwall > 0)
1679 srenew(groups->grpname, groups->ngrpname+nwall);
1680 grps = &(groups->grps[egcENER]);
1681 srenew(grps->nm_ind, grps->nr+nwall);
1682 for (i = 0; i < nwall; i++)
1684 sprintf(str, "wall%d", i);
1685 groups->grpname[groups->ngrpname] = put_symtab(symtab, str);
1686 grps->nm_ind[grps->nr++] = groups->ngrpname++;
1691 static void read_expandedparams(std::vector<t_inpfile> *inp,
1692 t_expanded *expand, warninp_t wi)
1694 /* read expanded ensemble parameters */
1695 printStringNewline(inp, "expanded ensemble variables");
1696 expand->nstexpanded = get_eint(inp, "nstexpanded", -1, wi);
1697 expand->elamstats = get_eeenum(inp, "lmc-stats", elamstats_names, wi);
1698 expand->elmcmove = get_eeenum(inp, "lmc-move", elmcmove_names, wi);
1699 expand->elmceq = get_eeenum(inp, "lmc-weights-equil", elmceq_names, wi);
1700 expand->equil_n_at_lam = get_eint(inp, "weight-equil-number-all-lambda", -1, wi);
1701 expand->equil_samples = get_eint(inp, "weight-equil-number-samples", -1, wi);
1702 expand->equil_steps = get_eint(inp, "weight-equil-number-steps", -1, wi);
1703 expand->equil_wl_delta = get_ereal(inp, "weight-equil-wl-delta", -1, wi);
1704 expand->equil_ratio = get_ereal(inp, "weight-equil-count-ratio", -1, wi);
1705 printStringNewline(inp, "Seed for Monte Carlo in lambda space");
1706 expand->lmc_seed = get_eint(inp, "lmc-seed", -1, wi);
1707 expand->mc_temp = get_ereal(inp, "mc-temperature", -1, wi);
1708 expand->lmc_repeats = get_eint(inp, "lmc-repeats", 1, wi);
1709 expand->gibbsdeltalam = get_eint(inp, "lmc-gibbsdelta", -1, wi);
1710 expand->lmc_forced_nstart = get_eint(inp, "lmc-forced-nstart", 0, wi);
1711 expand->bSymmetrizedTMatrix = (get_eeenum(inp, "symmetrized-transition-matrix", yesno_names, wi) != 0);
1712 expand->nstTij = get_eint(inp, "nst-transition-matrix", -1, wi);
1713 expand->minvarmin = get_eint(inp, "mininum-var-min", 100, wi); /*default is reasonable */
1714 expand->c_range = get_eint(inp, "weight-c-range", 0, wi); /* default is just C=0 */
1715 expand->wl_scale = get_ereal(inp, "wl-scale", 0.8, wi);
1716 expand->wl_ratio = get_ereal(inp, "wl-ratio", 0.8, wi);
1717 expand->init_wl_delta = get_ereal(inp, "init-wl-delta", 1.0, wi);
1718 expand->bWLoneovert = (get_eeenum(inp, "wl-oneovert", yesno_names, wi) != 0);
1721 /*! \brief Return whether an end state with the given coupling-lambda
1722 * value describes fully-interacting VDW.
1724 * \param[in] couple_lambda_value Enumeration ecouplam value describing the end state
1725 * \return Whether VDW is on (i.e. the user chose vdw or vdw-q in the .mdp file)
1727 static bool couple_lambda_has_vdw_on(int couple_lambda_value)
1729 return (couple_lambda_value == ecouplamVDW ||
1730 couple_lambda_value == ecouplamVDWQ);
1733 namespace
1736 class MdpErrorHandler : public gmx::IKeyValueTreeErrorHandler
1738 public:
1739 explicit MdpErrorHandler(warninp_t wi)
1740 : wi_(wi), mapping_(nullptr)
1744 void setBackMapping(const gmx::IKeyValueTreeBackMapping &mapping)
1746 mapping_ = &mapping;
1749 bool onError(gmx::UserInputError *ex, const gmx::KeyValueTreePath &context) override
1751 ex->prependContext(gmx::formatString("Error in mdp option \"%s\":",
1752 getOptionName(context).c_str()));
1753 std::string message = gmx::formatExceptionMessageToString(*ex);
1754 warning_error(wi_, message.c_str());
1755 return true;
1758 private:
1759 std::string getOptionName(const gmx::KeyValueTreePath &context)
1761 if (mapping_ != nullptr)
1763 gmx::KeyValueTreePath path = mapping_->originalPath(context);
1764 GMX_ASSERT(path.size() == 1, "Inconsistent mapping back to mdp options");
1765 return path[0];
1767 GMX_ASSERT(context.size() == 1, "Inconsistent context for mdp option parsing");
1768 return context[0];
1771 warninp_t wi_;
1772 const gmx::IKeyValueTreeBackMapping *mapping_;
1775 } // namespace
1777 void get_ir(const char *mdparin, const char *mdparout,
1778 gmx::MDModules *mdModules, t_inputrec *ir, t_gromppopts *opts,
1779 WriteMdpHeader writeMdpHeader, warninp_t wi)
1781 char *dumstr[2];
1782 double dumdub[2][6];
1783 int i, j, m;
1784 char warn_buf[STRLEN];
1785 t_lambda *fep = ir->fepvals;
1786 t_expanded *expand = ir->expandedvals;
1788 const char *no_names[] = { "no", nullptr };
1790 init_inputrec_strings();
1791 gmx::TextInputFile stream(mdparin);
1792 std::vector<t_inpfile> inp = read_inpfile(&stream, mdparin, wi);
1794 snew(dumstr[0], STRLEN);
1795 snew(dumstr[1], STRLEN);
1797 if (-1 == search_einp(inp, "cutoff-scheme"))
1799 sprintf(warn_buf,
1800 "%s did not specify a value for the .mdp option "
1801 "\"cutoff-scheme\". Probably it was first intended for use "
1802 "with GROMACS before 4.6. In 4.6, the Verlet scheme was "
1803 "introduced, but the group scheme was still the default. "
1804 "The default is now the Verlet scheme, so you will observe "
1805 "different behaviour.", mdparin);
1806 warning_note(wi, warn_buf);
1809 /* ignore the following deprecated commands */
1810 replace_inp_entry(inp, "title", nullptr);
1811 replace_inp_entry(inp, "cpp", nullptr);
1812 replace_inp_entry(inp, "domain-decomposition", nullptr);
1813 replace_inp_entry(inp, "andersen-seed", nullptr);
1814 replace_inp_entry(inp, "dihre", nullptr);
1815 replace_inp_entry(inp, "dihre-fc", nullptr);
1816 replace_inp_entry(inp, "dihre-tau", nullptr);
1817 replace_inp_entry(inp, "nstdihreout", nullptr);
1818 replace_inp_entry(inp, "nstcheckpoint", nullptr);
1819 replace_inp_entry(inp, "optimize-fft", nullptr);
1820 replace_inp_entry(inp, "adress_type", nullptr);
1821 replace_inp_entry(inp, "adress_const_wf", nullptr);
1822 replace_inp_entry(inp, "adress_ex_width", nullptr);
1823 replace_inp_entry(inp, "adress_hy_width", nullptr);
1824 replace_inp_entry(inp, "adress_ex_forcecap", nullptr);
1825 replace_inp_entry(inp, "adress_interface_correction", nullptr);
1826 replace_inp_entry(inp, "adress_site", nullptr);
1827 replace_inp_entry(inp, "adress_reference_coords", nullptr);
1828 replace_inp_entry(inp, "adress_tf_grp_names", nullptr);
1829 replace_inp_entry(inp, "adress_cg_grp_names", nullptr);
1830 replace_inp_entry(inp, "adress_do_hybridpairs", nullptr);
1831 replace_inp_entry(inp, "rlistlong", nullptr);
1832 replace_inp_entry(inp, "nstcalclr", nullptr);
1833 replace_inp_entry(inp, "pull-print-com2", nullptr);
1834 replace_inp_entry(inp, "gb-algorithm", nullptr);
1835 replace_inp_entry(inp, "nstgbradii", nullptr);
1836 replace_inp_entry(inp, "rgbradii", nullptr);
1837 replace_inp_entry(inp, "gb-epsilon-solvent", nullptr);
1838 replace_inp_entry(inp, "gb-saltconc", nullptr);
1839 replace_inp_entry(inp, "gb-obc-alpha", nullptr);
1840 replace_inp_entry(inp, "gb-obc-beta", nullptr);
1841 replace_inp_entry(inp, "gb-obc-gamma", nullptr);
1842 replace_inp_entry(inp, "gb-dielectric-offset", nullptr);
1843 replace_inp_entry(inp, "sa-algorithm", nullptr);
1844 replace_inp_entry(inp, "sa-surface-tension", nullptr);
1846 /* replace the following commands with the clearer new versions*/
1847 replace_inp_entry(inp, "unconstrained-start", "continuation");
1848 replace_inp_entry(inp, "foreign-lambda", "fep-lambdas");
1849 replace_inp_entry(inp, "verlet-buffer-drift", "verlet-buffer-tolerance");
1850 replace_inp_entry(inp, "nstxtcout", "nstxout-compressed");
1851 replace_inp_entry(inp, "xtc-grps", "compressed-x-grps");
1852 replace_inp_entry(inp, "xtc-precision", "compressed-x-precision");
1853 replace_inp_entry(inp, "pull-print-com1", "pull-print-com");
1855 printStringNewline(&inp, "VARIOUS PREPROCESSING OPTIONS");
1856 printStringNoNewline(&inp, "Preprocessor information: use cpp syntax.");
1857 printStringNoNewline(&inp, "e.g.: -I/home/joe/doe -I/home/mary/roe");
1858 setStringEntry(&inp, "include", opts->include, nullptr);
1859 printStringNoNewline(&inp, "e.g.: -DPOSRES -DFLEXIBLE (note these variable names are case sensitive)");
1860 setStringEntry(&inp, "define", opts->define, nullptr);
1862 printStringNewline(&inp, "RUN CONTROL PARAMETERS");
1863 ir->eI = get_eeenum(&inp, "integrator", ei_names, wi);
1864 printStringNoNewline(&inp, "Start time and timestep in ps");
1865 ir->init_t = get_ereal(&inp, "tinit", 0.0, wi);
1866 ir->delta_t = get_ereal(&inp, "dt", 0.001, wi);
1867 ir->nsteps = get_eint64(&inp, "nsteps", 0, wi);
1868 printStringNoNewline(&inp, "For exact run continuation or redoing part of a run");
1869 ir->init_step = get_eint64(&inp, "init-step", 0, wi);
1870 printStringNoNewline(&inp, "Part index is updated automatically on checkpointing (keeps files separate)");
1871 ir->simulation_part = get_eint(&inp, "simulation-part", 1, wi);
1872 printStringNoNewline(&inp, "mode for center of mass motion removal");
1873 ir->comm_mode = get_eeenum(&inp, "comm-mode", ecm_names, wi);
1874 printStringNoNewline(&inp, "number of steps for center of mass motion removal");
1875 ir->nstcomm = get_eint(&inp, "nstcomm", 100, wi);
1876 printStringNoNewline(&inp, "group(s) for center of mass motion removal");
1877 setStringEntry(&inp, "comm-grps", is->vcm, nullptr);
1879 printStringNewline(&inp, "LANGEVIN DYNAMICS OPTIONS");
1880 printStringNoNewline(&inp, "Friction coefficient (amu/ps) and random seed");
1881 ir->bd_fric = get_ereal(&inp, "bd-fric", 0.0, wi);
1882 ir->ld_seed = get_eint64(&inp, "ld-seed", -1, wi);
1884 /* Em stuff */
1885 printStringNewline(&inp, "ENERGY MINIMIZATION OPTIONS");
1886 printStringNoNewline(&inp, "Force tolerance and initial step-size");
1887 ir->em_tol = get_ereal(&inp, "emtol", 10.0, wi);
1888 ir->em_stepsize = get_ereal(&inp, "emstep", 0.01, wi);
1889 printStringNoNewline(&inp, "Max number of iterations in relax-shells");
1890 ir->niter = get_eint(&inp, "niter", 20, wi);
1891 printStringNoNewline(&inp, "Step size (ps^2) for minimization of flexible constraints");
1892 ir->fc_stepsize = get_ereal(&inp, "fcstep", 0, wi);
1893 printStringNoNewline(&inp, "Frequency of steepest descents steps when doing CG");
1894 ir->nstcgsteep = get_eint(&inp, "nstcgsteep", 1000, wi);
1895 ir->nbfgscorr = get_eint(&inp, "nbfgscorr", 10, wi);
1897 printStringNewline(&inp, "TEST PARTICLE INSERTION OPTIONS");
1898 ir->rtpi = get_ereal(&inp, "rtpi", 0.05, wi);
1900 /* Output options */
1901 printStringNewline(&inp, "OUTPUT CONTROL OPTIONS");
1902 printStringNoNewline(&inp, "Output frequency for coords (x), velocities (v) and forces (f)");
1903 ir->nstxout = get_eint(&inp, "nstxout", 0, wi);
1904 ir->nstvout = get_eint(&inp, "nstvout", 0, wi);
1905 ir->nstfout = get_eint(&inp, "nstfout", 0, wi);
1906 printStringNoNewline(&inp, "Output frequency for energies to log file and energy file");
1907 ir->nstlog = get_eint(&inp, "nstlog", 1000, wi);
1908 ir->nstcalcenergy = get_eint(&inp, "nstcalcenergy", 100, wi);
1909 ir->nstenergy = get_eint(&inp, "nstenergy", 1000, wi);
1910 printStringNoNewline(&inp, "Output frequency and precision for .xtc file");
1911 ir->nstxout_compressed = get_eint(&inp, "nstxout-compressed", 0, wi);
1912 ir->x_compression_precision = get_ereal(&inp, "compressed-x-precision", 1000.0, wi);
1913 printStringNoNewline(&inp, "This selects the subset of atoms for the compressed");
1914 printStringNoNewline(&inp, "trajectory file. You can select multiple groups. By");
1915 printStringNoNewline(&inp, "default, all atoms will be written.");
1916 setStringEntry(&inp, "compressed-x-grps", is->x_compressed_groups, nullptr);
1917 printStringNoNewline(&inp, "Selection of energy groups");
1918 setStringEntry(&inp, "energygrps", is->energy, nullptr);
1920 /* Neighbor searching */
1921 printStringNewline(&inp, "NEIGHBORSEARCHING PARAMETERS");
1922 printStringNoNewline(&inp, "cut-off scheme (Verlet: particle based cut-offs, group: using charge groups)");
1923 ir->cutoff_scheme = get_eeenum(&inp, "cutoff-scheme", ecutscheme_names, wi);
1924 printStringNoNewline(&inp, "nblist update frequency");
1925 ir->nstlist = get_eint(&inp, "nstlist", 10, wi);
1926 printStringNoNewline(&inp, "ns algorithm (simple or grid)");
1927 ir->ns_type = get_eeenum(&inp, "ns-type", ens_names, wi);
1928 printStringNoNewline(&inp, "Periodic boundary conditions: xyz, no, xy");
1929 ir->ePBC = get_eeenum(&inp, "pbc", epbc_names, wi);
1930 ir->bPeriodicMols = get_eeenum(&inp, "periodic-molecules", yesno_names, wi) != 0;
1931 printStringNoNewline(&inp, "Allowed energy error due to the Verlet buffer in kJ/mol/ps per atom,");
1932 printStringNoNewline(&inp, "a value of -1 means: use rlist");
1933 ir->verletbuf_tol = get_ereal(&inp, "verlet-buffer-tolerance", 0.005, wi);
1934 printStringNoNewline(&inp, "nblist cut-off");
1935 ir->rlist = get_ereal(&inp, "rlist", 1.0, wi);
1936 printStringNoNewline(&inp, "long-range cut-off for switched potentials");
1938 /* Electrostatics */
1939 printStringNewline(&inp, "OPTIONS FOR ELECTROSTATICS AND VDW");
1940 printStringNoNewline(&inp, "Method for doing electrostatics");
1941 ir->coulombtype = get_eeenum(&inp, "coulombtype", eel_names, wi);
1942 ir->coulomb_modifier = get_eeenum(&inp, "coulomb-modifier", eintmod_names, wi);
1943 printStringNoNewline(&inp, "cut-off lengths");
1944 ir->rcoulomb_switch = get_ereal(&inp, "rcoulomb-switch", 0.0, wi);
1945 ir->rcoulomb = get_ereal(&inp, "rcoulomb", 1.0, wi);
1946 printStringNoNewline(&inp, "Relative dielectric constant for the medium and the reaction field");
1947 ir->epsilon_r = get_ereal(&inp, "epsilon-r", 1.0, wi);
1948 ir->epsilon_rf = get_ereal(&inp, "epsilon-rf", 0.0, wi);
1949 printStringNoNewline(&inp, "Method for doing Van der Waals");
1950 ir->vdwtype = get_eeenum(&inp, "vdw-type", evdw_names, wi);
1951 ir->vdw_modifier = get_eeenum(&inp, "vdw-modifier", eintmod_names, wi);
1952 printStringNoNewline(&inp, "cut-off lengths");
1953 ir->rvdw_switch = get_ereal(&inp, "rvdw-switch", 0.0, wi);
1954 ir->rvdw = get_ereal(&inp, "rvdw", 1.0, wi);
1955 printStringNoNewline(&inp, "Apply long range dispersion corrections for Energy and Pressure");
1956 ir->eDispCorr = get_eeenum(&inp, "DispCorr", edispc_names, wi);
1957 printStringNoNewline(&inp, "Extension of the potential lookup tables beyond the cut-off");
1958 ir->tabext = get_ereal(&inp, "table-extension", 1.0, wi);
1959 printStringNoNewline(&inp, "Separate tables between energy group pairs");
1960 setStringEntry(&inp, "energygrp-table", is->egptable, nullptr);
1961 printStringNoNewline(&inp, "Spacing for the PME/PPPM FFT grid");
1962 ir->fourier_spacing = get_ereal(&inp, "fourierspacing", 0.12, wi);
1963 printStringNoNewline(&inp, "FFT grid size, when a value is 0 fourierspacing will be used");
1964 ir->nkx = get_eint(&inp, "fourier-nx", 0, wi);
1965 ir->nky = get_eint(&inp, "fourier-ny", 0, wi);
1966 ir->nkz = get_eint(&inp, "fourier-nz", 0, wi);
1967 printStringNoNewline(&inp, "EWALD/PME/PPPM parameters");
1968 ir->pme_order = get_eint(&inp, "pme-order", 4, wi);
1969 ir->ewald_rtol = get_ereal(&inp, "ewald-rtol", 0.00001, wi);
1970 ir->ewald_rtol_lj = get_ereal(&inp, "ewald-rtol-lj", 0.001, wi);
1971 ir->ljpme_combination_rule = get_eeenum(&inp, "lj-pme-comb-rule", eljpme_names, wi);
1972 ir->ewald_geometry = get_eeenum(&inp, "ewald-geometry", eewg_names, wi);
1973 ir->epsilon_surface = get_ereal(&inp, "epsilon-surface", 0.0, wi);
1975 /* Implicit solvation is no longer supported, but we need grompp
1976 to be able to refuse old .mdp files that would have built a tpr
1977 to run it. Thus, only "no" is accepted. */
1978 ir->implicit_solvent = (get_eeenum(&inp, "implicit-solvent", no_names, wi) != 0);
1980 /* Coupling stuff */
1981 printStringNewline(&inp, "OPTIONS FOR WEAK COUPLING ALGORITHMS");
1982 printStringNoNewline(&inp, "Temperature coupling");
1983 ir->etc = get_eeenum(&inp, "tcoupl", etcoupl_names, wi);
1984 ir->nsttcouple = get_eint(&inp, "nsttcouple", -1, wi);
1985 ir->opts.nhchainlength = get_eint(&inp, "nh-chain-length", 10, wi);
1986 ir->bPrintNHChains = (get_eeenum(&inp, "print-nose-hoover-chain-variables", yesno_names, wi) != 0);
1987 printStringNoNewline(&inp, "Groups to couple separately");
1988 setStringEntry(&inp, "tc-grps", is->tcgrps, nullptr);
1989 printStringNoNewline(&inp, "Time constant (ps) and reference temperature (K)");
1990 setStringEntry(&inp, "tau-t", is->tau_t, nullptr);
1991 setStringEntry(&inp, "ref-t", is->ref_t, nullptr);
1992 printStringNoNewline(&inp, "pressure coupling");
1993 ir->epc = get_eeenum(&inp, "pcoupl", epcoupl_names, wi);
1994 ir->epct = get_eeenum(&inp, "pcoupltype", epcoupltype_names, wi);
1995 ir->nstpcouple = get_eint(&inp, "nstpcouple", -1, wi);
1996 printStringNoNewline(&inp, "Time constant (ps), compressibility (1/bar) and reference P (bar)");
1997 ir->tau_p = get_ereal(&inp, "tau-p", 1.0, wi);
1998 setStringEntry(&inp, "compressibility", dumstr[0], nullptr);
1999 setStringEntry(&inp, "ref-p", dumstr[1], nullptr);
2000 printStringNoNewline(&inp, "Scaling of reference coordinates, No, All or COM");
2001 ir->refcoord_scaling = get_eeenum(&inp, "refcoord-scaling", erefscaling_names, wi);
2003 /* QMMM */
2004 printStringNewline(&inp, "OPTIONS FOR QMMM calculations");
2005 ir->bQMMM = (get_eeenum(&inp, "QMMM", yesno_names, wi) != 0);
2006 printStringNoNewline(&inp, "Groups treated Quantum Mechanically");
2007 setStringEntry(&inp, "QMMM-grps", is->QMMM, nullptr);
2008 printStringNoNewline(&inp, "QM method");
2009 setStringEntry(&inp, "QMmethod", is->QMmethod, nullptr);
2010 printStringNoNewline(&inp, "QMMM scheme");
2011 ir->QMMMscheme = get_eeenum(&inp, "QMMMscheme", eQMMMscheme_names, wi);
2012 printStringNoNewline(&inp, "QM basisset");
2013 setStringEntry(&inp, "QMbasis", is->QMbasis, nullptr);
2014 printStringNoNewline(&inp, "QM charge");
2015 setStringEntry(&inp, "QMcharge", is->QMcharge, nullptr);
2016 printStringNoNewline(&inp, "QM multiplicity");
2017 setStringEntry(&inp, "QMmult", is->QMmult, nullptr);
2018 printStringNoNewline(&inp, "Surface Hopping");
2019 setStringEntry(&inp, "SH", is->bSH, nullptr);
2020 printStringNoNewline(&inp, "CAS space options");
2021 setStringEntry(&inp, "CASorbitals", is->CASorbitals, nullptr);
2022 setStringEntry(&inp, "CASelectrons", is->CASelectrons, nullptr);
2023 setStringEntry(&inp, "SAon", is->SAon, nullptr);
2024 setStringEntry(&inp, "SAoff", is->SAoff, nullptr);
2025 setStringEntry(&inp, "SAsteps", is->SAsteps, nullptr);
2026 printStringNoNewline(&inp, "Scale factor for MM charges");
2027 ir->scalefactor = get_ereal(&inp, "MMChargeScaleFactor", 1.0, wi);
2029 /* Simulated annealing */
2030 printStringNewline(&inp, "SIMULATED ANNEALING");
2031 printStringNoNewline(&inp, "Type of annealing for each temperature group (no/single/periodic)");
2032 setStringEntry(&inp, "annealing", is->anneal, nullptr);
2033 printStringNoNewline(&inp, "Number of time points to use for specifying annealing in each group");
2034 setStringEntry(&inp, "annealing-npoints", is->anneal_npoints, nullptr);
2035 printStringNoNewline(&inp, "List of times at the annealing points for each group");
2036 setStringEntry(&inp, "annealing-time", is->anneal_time, nullptr);
2037 printStringNoNewline(&inp, "Temp. at each annealing point, for each group.");
2038 setStringEntry(&inp, "annealing-temp", is->anneal_temp, nullptr);
2040 /* Startup run */
2041 printStringNewline(&inp, "GENERATE VELOCITIES FOR STARTUP RUN");
2042 opts->bGenVel = (get_eeenum(&inp, "gen-vel", yesno_names, wi) != 0);
2043 opts->tempi = get_ereal(&inp, "gen-temp", 300.0, wi);
2044 opts->seed = get_eint(&inp, "gen-seed", -1, wi);
2046 /* Shake stuff */
2047 printStringNewline(&inp, "OPTIONS FOR BONDS");
2048 opts->nshake = get_eeenum(&inp, "constraints", constraints, wi);
2049 printStringNoNewline(&inp, "Type of constraint algorithm");
2050 ir->eConstrAlg = get_eeenum(&inp, "constraint-algorithm", econstr_names, wi);
2051 printStringNoNewline(&inp, "Do not constrain the start configuration");
2052 ir->bContinuation = (get_eeenum(&inp, "continuation", yesno_names, wi) != 0);
2053 printStringNoNewline(&inp, "Use successive overrelaxation to reduce the number of shake iterations");
2054 ir->bShakeSOR = (get_eeenum(&inp, "Shake-SOR", yesno_names, wi) != 0);
2055 printStringNoNewline(&inp, "Relative tolerance of shake");
2056 ir->shake_tol = get_ereal(&inp, "shake-tol", 0.0001, wi);
2057 printStringNoNewline(&inp, "Highest order in the expansion of the constraint coupling matrix");
2058 ir->nProjOrder = get_eint(&inp, "lincs-order", 4, wi);
2059 printStringNoNewline(&inp, "Number of iterations in the final step of LINCS. 1 is fine for");
2060 printStringNoNewline(&inp, "normal simulations, but use 2 to conserve energy in NVE runs.");
2061 printStringNoNewline(&inp, "For energy minimization with constraints it should be 4 to 8.");
2062 ir->nLincsIter = get_eint(&inp, "lincs-iter", 1, wi);
2063 printStringNoNewline(&inp, "Lincs will write a warning to the stderr if in one step a bond");
2064 printStringNoNewline(&inp, "rotates over more degrees than");
2065 ir->LincsWarnAngle = get_ereal(&inp, "lincs-warnangle", 30.0, wi);
2066 printStringNoNewline(&inp, "Convert harmonic bonds to morse potentials");
2067 opts->bMorse = (get_eeenum(&inp, "morse", yesno_names, wi) != 0);
2069 /* Energy group exclusions */
2070 printStringNewline(&inp, "ENERGY GROUP EXCLUSIONS");
2071 printStringNoNewline(&inp, "Pairs of energy groups for which all non-bonded interactions are excluded");
2072 setStringEntry(&inp, "energygrp-excl", is->egpexcl, nullptr);
2074 /* Walls */
2075 printStringNewline(&inp, "WALLS");
2076 printStringNoNewline(&inp, "Number of walls, type, atom types, densities and box-z scale factor for Ewald");
2077 ir->nwall = get_eint(&inp, "nwall", 0, wi);
2078 ir->wall_type = get_eeenum(&inp, "wall-type", ewt_names, wi);
2079 ir->wall_r_linpot = get_ereal(&inp, "wall-r-linpot", -1, wi);
2080 setStringEntry(&inp, "wall-atomtype", is->wall_atomtype, nullptr);
2081 setStringEntry(&inp, "wall-density", is->wall_density, nullptr);
2082 ir->wall_ewald_zfac = get_ereal(&inp, "wall-ewald-zfac", 3, wi);
2084 /* COM pulling */
2085 printStringNewline(&inp, "COM PULLING");
2086 ir->bPull = (get_eeenum(&inp, "pull", yesno_names, wi) != 0);
2087 if (ir->bPull)
2089 snew(ir->pull, 1);
2090 is->pull_grp = read_pullparams(&inp, ir->pull, wi);
2093 /* AWH biasing
2094 NOTE: needs COM pulling input */
2095 printStringNewline(&inp, "AWH biasing");
2096 ir->bDoAwh = (get_eeenum(&inp, "awh", yesno_names, wi) != 0);
2097 if (ir->bDoAwh)
2099 if (ir->bPull)
2101 ir->awhParams = gmx::readAndCheckAwhParams(&inp, ir, wi);
2103 else
2105 gmx_fatal(FARGS, "AWH biasing is only compatible with COM pulling turned on");
2109 /* Enforced rotation */
2110 printStringNewline(&inp, "ENFORCED ROTATION");
2111 printStringNoNewline(&inp, "Enforced rotation: No or Yes");
2112 ir->bRot = (get_eeenum(&inp, "rotation", yesno_names, wi) != 0);
2113 if (ir->bRot)
2115 snew(ir->rot, 1);
2116 is->rot_grp = read_rotparams(&inp, ir->rot, wi);
2119 /* Interactive MD */
2120 ir->bIMD = FALSE;
2121 printStringNewline(&inp, "Group to display and/or manipulate in interactive MD session");
2122 setStringEntry(&inp, "IMD-group", is->imd_grp, nullptr);
2123 if (is->imd_grp[0] != '\0')
2125 snew(ir->imd, 1);
2126 ir->bIMD = TRUE;
2129 /* Refinement */
2130 printStringNewline(&inp, "NMR refinement stuff");
2131 printStringNoNewline(&inp, "Distance restraints type: No, Simple or Ensemble");
2132 ir->eDisre = get_eeenum(&inp, "disre", edisre_names, wi);
2133 printStringNoNewline(&inp, "Force weighting of pairs in one distance restraint: Conservative or Equal");
2134 ir->eDisreWeighting = get_eeenum(&inp, "disre-weighting", edisreweighting_names, wi);
2135 printStringNoNewline(&inp, "Use sqrt of the time averaged times the instantaneous violation");
2136 ir->bDisreMixed = (get_eeenum(&inp, "disre-mixed", yesno_names, wi) != 0);
2137 ir->dr_fc = get_ereal(&inp, "disre-fc", 1000.0, wi);
2138 ir->dr_tau = get_ereal(&inp, "disre-tau", 0.0, wi);
2139 printStringNoNewline(&inp, "Output frequency for pair distances to energy file");
2140 ir->nstdisreout = get_eint(&inp, "nstdisreout", 100, wi);
2141 printStringNoNewline(&inp, "Orientation restraints: No or Yes");
2142 opts->bOrire = (get_eeenum(&inp, "orire", yesno_names, wi) != 0);
2143 printStringNoNewline(&inp, "Orientation restraints force constant and tau for time averaging");
2144 ir->orires_fc = get_ereal(&inp, "orire-fc", 0.0, wi);
2145 ir->orires_tau = get_ereal(&inp, "orire-tau", 0.0, wi);
2146 setStringEntry(&inp, "orire-fitgrp", is->orirefitgrp, nullptr);
2147 printStringNoNewline(&inp, "Output frequency for trace(SD) and S to energy file");
2148 ir->nstorireout = get_eint(&inp, "nstorireout", 100, wi);
2150 /* free energy variables */
2151 printStringNewline(&inp, "Free energy variables");
2152 ir->efep = get_eeenum(&inp, "free-energy", efep_names, wi);
2153 setStringEntry(&inp, "couple-moltype", is->couple_moltype, nullptr);
2154 opts->couple_lam0 = get_eeenum(&inp, "couple-lambda0", couple_lam, wi);
2155 opts->couple_lam1 = get_eeenum(&inp, "couple-lambda1", couple_lam, wi);
2156 opts->bCoupleIntra = (get_eeenum(&inp, "couple-intramol", yesno_names, wi) != 0);
2158 fep->init_lambda = get_ereal(&inp, "init-lambda", -1, wi); /* start with -1 so
2159 we can recognize if
2160 it was not entered */
2161 fep->init_fep_state = get_eint(&inp, "init-lambda-state", -1, wi);
2162 fep->delta_lambda = get_ereal(&inp, "delta-lambda", 0.0, wi);
2163 fep->nstdhdl = get_eint(&inp, "nstdhdl", 50, wi);
2164 setStringEntry(&inp, "fep-lambdas", is->fep_lambda[efptFEP], nullptr);
2165 setStringEntry(&inp, "mass-lambdas", is->fep_lambda[efptMASS], nullptr);
2166 setStringEntry(&inp, "coul-lambdas", is->fep_lambda[efptCOUL], nullptr);
2167 setStringEntry(&inp, "vdw-lambdas", is->fep_lambda[efptVDW], nullptr);
2168 setStringEntry(&inp, "bonded-lambdas", is->fep_lambda[efptBONDED], nullptr);
2169 setStringEntry(&inp, "restraint-lambdas", is->fep_lambda[efptRESTRAINT], nullptr);
2170 setStringEntry(&inp, "temperature-lambdas", is->fep_lambda[efptTEMPERATURE], nullptr);
2171 fep->lambda_neighbors = get_eint(&inp, "calc-lambda-neighbors", 1, wi);
2172 setStringEntry(&inp, "init-lambda-weights", is->lambda_weights, nullptr);
2173 fep->edHdLPrintEnergy = get_eeenum(&inp, "dhdl-print-energy", edHdLPrintEnergy_names, wi);
2174 fep->sc_alpha = get_ereal(&inp, "sc-alpha", 0.0, wi);
2175 fep->sc_power = get_eint(&inp, "sc-power", 1, wi);
2176 fep->sc_r_power = get_ereal(&inp, "sc-r-power", 6.0, wi);
2177 fep->sc_sigma = get_ereal(&inp, "sc-sigma", 0.3, wi);
2178 fep->bScCoul = (get_eeenum(&inp, "sc-coul", yesno_names, wi) != 0);
2179 fep->dh_hist_size = get_eint(&inp, "dh_hist_size", 0, wi);
2180 fep->dh_hist_spacing = get_ereal(&inp, "dh_hist_spacing", 0.1, wi);
2181 fep->separate_dhdl_file = get_eeenum(&inp, "separate-dhdl-file", separate_dhdl_file_names, wi);
2182 fep->dhdl_derivatives = get_eeenum(&inp, "dhdl-derivatives", dhdl_derivatives_names, wi);
2183 fep->dh_hist_size = get_eint(&inp, "dh_hist_size", 0, wi);
2184 fep->dh_hist_spacing = get_ereal(&inp, "dh_hist_spacing", 0.1, wi);
2186 /* Non-equilibrium MD stuff */
2187 printStringNewline(&inp, "Non-equilibrium MD stuff");
2188 setStringEntry(&inp, "acc-grps", is->accgrps, nullptr);
2189 setStringEntry(&inp, "accelerate", is->acc, nullptr);
2190 setStringEntry(&inp, "freezegrps", is->freeze, nullptr);
2191 setStringEntry(&inp, "freezedim", is->frdim, nullptr);
2192 ir->cos_accel = get_ereal(&inp, "cos-acceleration", 0, wi);
2193 setStringEntry(&inp, "deform", is->deform, nullptr);
2195 /* simulated tempering variables */
2196 printStringNewline(&inp, "simulated tempering variables");
2197 ir->bSimTemp = (get_eeenum(&inp, "simulated-tempering", yesno_names, wi) != 0);
2198 ir->simtempvals->eSimTempScale = get_eeenum(&inp, "simulated-tempering-scaling", esimtemp_names, wi);
2199 ir->simtempvals->simtemp_low = get_ereal(&inp, "sim-temp-low", 300.0, wi);
2200 ir->simtempvals->simtemp_high = get_ereal(&inp, "sim-temp-high", 300.0, wi);
2202 /* expanded ensemble variables */
2203 if (ir->efep == efepEXPANDED || ir->bSimTemp)
2205 read_expandedparams(&inp, expand, wi);
2208 /* Electric fields */
2210 gmx::KeyValueTreeObject convertedValues = flatKeyValueTreeFromInpFile(inp);
2211 gmx::KeyValueTreeTransformer transform;
2212 transform.rules()->addRule()
2213 .keyMatchType("/", gmx::StringCompareType::CaseAndDashInsensitive);
2214 mdModules->initMdpTransform(transform.rules());
2215 for (const auto &path : transform.mappedPaths())
2217 GMX_ASSERT(path.size() == 1, "Inconsistent mapping back to mdp options");
2218 mark_einp_set(inp, path[0].c_str());
2220 MdpErrorHandler errorHandler(wi);
2221 auto result
2222 = transform.transform(convertedValues, &errorHandler);
2223 ir->params = new gmx::KeyValueTreeObject(result.object());
2224 mdModules->adjustInputrecBasedOnModules(ir);
2225 errorHandler.setBackMapping(result.backMapping());
2226 mdModules->assignOptionsToModules(*ir->params, &errorHandler);
2229 /* Ion/water position swapping ("computational electrophysiology") */
2230 printStringNewline(&inp, "Ion/water position swapping for computational electrophysiology setups");
2231 printStringNoNewline(&inp, "Swap positions along direction: no, X, Y, Z");
2232 ir->eSwapCoords = get_eeenum(&inp, "swapcoords", eSwapTypes_names, wi);
2233 if (ir->eSwapCoords != eswapNO)
2235 char buf[STRLEN];
2236 int nIonTypes;
2239 snew(ir->swap, 1);
2240 printStringNoNewline(&inp, "Swap attempt frequency");
2241 ir->swap->nstswap = get_eint(&inp, "swap-frequency", 1, wi);
2242 printStringNoNewline(&inp, "Number of ion types to be controlled");
2243 nIonTypes = get_eint(&inp, "iontypes", 1, wi);
2244 if (nIonTypes < 1)
2246 warning_error(wi, "You need to provide at least one ion type for position exchanges.");
2248 ir->swap->ngrp = nIonTypes + eSwapFixedGrpNR;
2249 snew(ir->swap->grp, ir->swap->ngrp);
2250 for (i = 0; i < ir->swap->ngrp; i++)
2252 snew(ir->swap->grp[i].molname, STRLEN);
2254 printStringNoNewline(&inp, "Two index groups that contain the compartment-partitioning atoms");
2255 setStringEntry(&inp, "split-group0", ir->swap->grp[eGrpSplit0].molname, nullptr);
2256 setStringEntry(&inp, "split-group1", ir->swap->grp[eGrpSplit1].molname, nullptr);
2257 printStringNoNewline(&inp, "Use center of mass of split groups (yes/no), otherwise center of geometry is used");
2258 ir->swap->massw_split[0] = (get_eeenum(&inp, "massw-split0", yesno_names, wi) != 0);
2259 ir->swap->massw_split[1] = (get_eeenum(&inp, "massw-split1", yesno_names, wi) != 0);
2261 printStringNoNewline(&inp, "Name of solvent molecules");
2262 setStringEntry(&inp, "solvent-group", ir->swap->grp[eGrpSolvent].molname, nullptr);
2264 printStringNoNewline(&inp, "Split cylinder: radius, upper and lower extension (nm) (this will define the channels)");
2265 printStringNoNewline(&inp, "Note that the split cylinder settings do not have an influence on the swapping protocol,");
2266 printStringNoNewline(&inp, "however, if correctly defined, the permeation events are recorded per channel");
2267 ir->swap->cyl0r = get_ereal(&inp, "cyl0-r", 2.0, wi);
2268 ir->swap->cyl0u = get_ereal(&inp, "cyl0-up", 1.0, wi);
2269 ir->swap->cyl0l = get_ereal(&inp, "cyl0-down", 1.0, wi);
2270 ir->swap->cyl1r = get_ereal(&inp, "cyl1-r", 2.0, wi);
2271 ir->swap->cyl1u = get_ereal(&inp, "cyl1-up", 1.0, wi);
2272 ir->swap->cyl1l = get_ereal(&inp, "cyl1-down", 1.0, wi);
2274 printStringNoNewline(&inp, "Average the number of ions per compartment over these many swap attempt steps");
2275 ir->swap->nAverage = get_eint(&inp, "coupl-steps", 10, wi);
2277 printStringNoNewline(&inp, "Names of the ion types that can be exchanged with solvent molecules,");
2278 printStringNoNewline(&inp, "and the requested number of ions of this type in compartments A and B");
2279 printStringNoNewline(&inp, "-1 means fix the numbers as found in step 0");
2280 for (i = 0; i < nIonTypes; i++)
2282 int ig = eSwapFixedGrpNR + i;
2284 sprintf(buf, "iontype%d-name", i);
2285 setStringEntry(&inp, buf, ir->swap->grp[ig].molname, nullptr);
2286 sprintf(buf, "iontype%d-in-A", i);
2287 ir->swap->grp[ig].nmolReq[0] = get_eint(&inp, buf, -1, wi);
2288 sprintf(buf, "iontype%d-in-B", i);
2289 ir->swap->grp[ig].nmolReq[1] = get_eint(&inp, buf, -1, wi);
2292 printStringNoNewline(&inp, "By default (i.e. bulk offset = 0.0), ion/water exchanges happen between layers");
2293 printStringNoNewline(&inp, "at maximum distance (= bulk concentration) to the split group layers. However,");
2294 printStringNoNewline(&inp, "an offset b (-1.0 < b < +1.0) can be specified to offset the bulk layer from the middle at 0.0");
2295 printStringNoNewline(&inp, "towards one of the compartment-partitioning layers (at +/- 1.0).");
2296 ir->swap->bulkOffset[0] = get_ereal(&inp, "bulk-offsetA", 0.0, wi);
2297 ir->swap->bulkOffset[1] = get_ereal(&inp, "bulk-offsetB", 0.0, wi);
2298 if (!(ir->swap->bulkOffset[0] > -1.0 && ir->swap->bulkOffset[0] < 1.0)
2299 || !(ir->swap->bulkOffset[1] > -1.0 && ir->swap->bulkOffset[1] < 1.0) )
2301 warning_error(wi, "Bulk layer offsets must be > -1.0 and < 1.0 !");
2304 printStringNoNewline(&inp, "Start to swap ions if threshold difference to requested count is reached");
2305 ir->swap->threshold = get_ereal(&inp, "threshold", 1.0, wi);
2308 /* AdResS is no longer supported, but we need grompp to be able to
2309 refuse to process old .mdp files that used it. */
2310 ir->bAdress = (get_eeenum(&inp, "adress", no_names, wi) != 0);
2312 /* User defined thingies */
2313 printStringNewline(&inp, "User defined thingies");
2314 setStringEntry(&inp, "user1-grps", is->user1, nullptr);
2315 setStringEntry(&inp, "user2-grps", is->user2, nullptr);
2316 ir->userint1 = get_eint(&inp, "userint1", 0, wi);
2317 ir->userint2 = get_eint(&inp, "userint2", 0, wi);
2318 ir->userint3 = get_eint(&inp, "userint3", 0, wi);
2319 ir->userint4 = get_eint(&inp, "userint4", 0, wi);
2320 ir->userreal1 = get_ereal(&inp, "userreal1", 0, wi);
2321 ir->userreal2 = get_ereal(&inp, "userreal2", 0, wi);
2322 ir->userreal3 = get_ereal(&inp, "userreal3", 0, wi);
2323 ir->userreal4 = get_ereal(&inp, "userreal4", 0, wi);
2324 #undef CTYPE
2327 gmx::TextOutputFile stream(mdparout);
2328 write_inpfile(&stream, mdparout, &inp, FALSE, writeMdpHeader, wi);
2330 // Transform module data into a flat key-value tree for output.
2331 gmx::KeyValueTreeBuilder builder;
2332 gmx::KeyValueTreeObjectBuilder builderObject = builder.rootObject();
2333 mdModules->buildMdpOutput(&builderObject);
2335 gmx::TextWriter writer(&stream);
2336 writeKeyValueTreeAsMdp(&writer, builder.build());
2338 stream.close();
2341 /* Process options if necessary */
2342 for (m = 0; m < 2; m++)
2344 for (i = 0; i < 2*DIM; i++)
2346 dumdub[m][i] = 0.0;
2348 if (ir->epc)
2350 switch (ir->epct)
2352 case epctISOTROPIC:
2353 if (sscanf(dumstr[m], "%lf", &(dumdub[m][XX])) != 1)
2355 warning_error(wi, "Pressure coupling incorrect number of values (I need exactly 1)");
2357 dumdub[m][YY] = dumdub[m][ZZ] = dumdub[m][XX];
2358 break;
2359 case epctSEMIISOTROPIC:
2360 case epctSURFACETENSION:
2361 if (sscanf(dumstr[m], "%lf%lf", &(dumdub[m][XX]), &(dumdub[m][ZZ])) != 2)
2363 warning_error(wi, "Pressure coupling incorrect number of values (I need exactly 2)");
2365 dumdub[m][YY] = dumdub[m][XX];
2366 break;
2367 case epctANISOTROPIC:
2368 if (sscanf(dumstr[m], "%lf%lf%lf%lf%lf%lf",
2369 &(dumdub[m][XX]), &(dumdub[m][YY]), &(dumdub[m][ZZ]),
2370 &(dumdub[m][3]), &(dumdub[m][4]), &(dumdub[m][5])) != 6)
2372 warning_error(wi, "Pressure coupling incorrect number of values (I need exactly 6)");
2374 break;
2375 default:
2376 gmx_fatal(FARGS, "Pressure coupling type %s not implemented yet",
2377 epcoupltype_names[ir->epct]);
2381 clear_mat(ir->ref_p);
2382 clear_mat(ir->compress);
2383 for (i = 0; i < DIM; i++)
2385 ir->ref_p[i][i] = dumdub[1][i];
2386 ir->compress[i][i] = dumdub[0][i];
2388 if (ir->epct == epctANISOTROPIC)
2390 ir->ref_p[XX][YY] = dumdub[1][3];
2391 ir->ref_p[XX][ZZ] = dumdub[1][4];
2392 ir->ref_p[YY][ZZ] = dumdub[1][5];
2393 if (ir->ref_p[XX][YY] != 0 && ir->ref_p[XX][ZZ] != 0 && ir->ref_p[YY][ZZ] != 0)
2395 warning(wi, "All off-diagonal reference pressures are non-zero. Are you sure you want to apply a threefold shear stress?\n");
2397 ir->compress[XX][YY] = dumdub[0][3];
2398 ir->compress[XX][ZZ] = dumdub[0][4];
2399 ir->compress[YY][ZZ] = dumdub[0][5];
2400 for (i = 0; i < DIM; i++)
2402 for (m = 0; m < i; m++)
2404 ir->ref_p[i][m] = ir->ref_p[m][i];
2405 ir->compress[i][m] = ir->compress[m][i];
2410 if (ir->comm_mode == ecmNO)
2412 ir->nstcomm = 0;
2415 opts->couple_moltype = nullptr;
2416 if (strlen(is->couple_moltype) > 0)
2418 if (ir->efep != efepNO)
2420 opts->couple_moltype = gmx_strdup(is->couple_moltype);
2421 if (opts->couple_lam0 == opts->couple_lam1)
2423 warning(wi, "The lambda=0 and lambda=1 states for coupling are identical");
2425 if (ir->eI == eiMD && (opts->couple_lam0 == ecouplamNONE ||
2426 opts->couple_lam1 == ecouplamNONE))
2428 warning(wi, "For proper sampling of the (nearly) decoupled state, stochastic dynamics should be used");
2431 else
2433 warning_note(wi, "Free energy is turned off, so we will not decouple the molecule listed in your input.");
2436 /* FREE ENERGY AND EXPANDED ENSEMBLE OPTIONS */
2437 if (ir->efep != efepNO)
2439 if (fep->delta_lambda > 0)
2441 ir->efep = efepSLOWGROWTH;
2445 if (fep->edHdLPrintEnergy == edHdLPrintEnergyYES)
2447 fep->edHdLPrintEnergy = edHdLPrintEnergyTOTAL;
2448 warning_note(wi, "Old option for dhdl-print-energy given: "
2449 "changing \"yes\" to \"total\"\n");
2452 if (ir->bSimTemp && (fep->edHdLPrintEnergy == edHdLPrintEnergyNO))
2454 /* always print out the energy to dhdl if we are doing
2455 expanded ensemble, since we need the total energy for
2456 analysis if the temperature is changing. In some
2457 conditions one may only want the potential energy, so
2458 we will allow that if the appropriate mdp setting has
2459 been enabled. Otherwise, total it is:
2461 fep->edHdLPrintEnergy = edHdLPrintEnergyTOTAL;
2464 if ((ir->efep != efepNO) || ir->bSimTemp)
2466 ir->bExpanded = FALSE;
2467 if ((ir->efep == efepEXPANDED) || ir->bSimTemp)
2469 ir->bExpanded = TRUE;
2471 do_fep_params(ir, is->fep_lambda, is->lambda_weights, wi);
2472 if (ir->bSimTemp) /* done after fep params */
2474 do_simtemp_params(ir);
2477 /* Because sc-coul (=FALSE by default) only acts on the lambda state
2478 * setup and not on the old way of specifying the free-energy setup,
2479 * we should check for using soft-core when not needed, since that
2480 * can complicate the sampling significantly.
2481 * Note that we only check for the automated coupling setup.
2482 * If the (advanced) user does FEP through manual topology changes,
2483 * this check will not be triggered.
2485 if (ir->efep != efepNO && ir->fepvals->n_lambda == 0 &&
2486 ir->fepvals->sc_alpha != 0 &&
2487 (couple_lambda_has_vdw_on(opts->couple_lam0) &&
2488 couple_lambda_has_vdw_on(opts->couple_lam1)))
2490 warning(wi, "You are using soft-core interactions while the Van der Waals interactions are not decoupled (note that the sc-coul option is only active when using lambda states). Although this will not lead to errors, you will need much more sampling than without soft-core interactions. Consider using sc-alpha=0.");
2493 else
2495 ir->fepvals->n_lambda = 0;
2498 /* WALL PARAMETERS */
2500 do_wall_params(ir, is->wall_atomtype, is->wall_density, opts, wi);
2502 /* ORIENTATION RESTRAINT PARAMETERS */
2504 if (opts->bOrire && gmx::splitString(is->orirefitgrp).size() != 1)
2506 warning_error(wi, "ERROR: Need one orientation restraint fit group\n");
2509 /* DEFORMATION PARAMETERS */
2511 clear_mat(ir->deform);
2512 for (i = 0; i < 6; i++)
2514 dumdub[0][i] = 0;
2517 double gmx_unused canary;
2518 int ndeform = sscanf(is->deform, "%lf %lf %lf %lf %lf %lf %lf",
2519 &(dumdub[0][0]), &(dumdub[0][1]), &(dumdub[0][2]),
2520 &(dumdub[0][3]), &(dumdub[0][4]), &(dumdub[0][5]), &canary);
2522 if (strlen(is->deform) > 0 && ndeform != 6)
2524 warning_error(wi, gmx::formatString("Cannot parse exactly 6 box deformation velocities from string '%s'", is->deform).c_str());
2526 for (i = 0; i < 3; i++)
2528 ir->deform[i][i] = dumdub[0][i];
2530 ir->deform[YY][XX] = dumdub[0][3];
2531 ir->deform[ZZ][XX] = dumdub[0][4];
2532 ir->deform[ZZ][YY] = dumdub[0][5];
2533 if (ir->epc != epcNO)
2535 for (i = 0; i < 3; i++)
2537 for (j = 0; j <= i; j++)
2539 if (ir->deform[i][j] != 0 && ir->compress[i][j] != 0)
2541 warning_error(wi, "A box element has deform set and compressibility > 0");
2545 for (i = 0; i < 3; i++)
2547 for (j = 0; j < i; j++)
2549 if (ir->deform[i][j] != 0)
2551 for (m = j; m < DIM; m++)
2553 if (ir->compress[m][j] != 0)
2555 sprintf(warn_buf, "An off-diagonal box element has deform set while compressibility > 0 for the same component of another box vector, this might lead to spurious periodicity effects.");
2556 warning(wi, warn_buf);
2564 /* Ion/water position swapping checks */
2565 if (ir->eSwapCoords != eswapNO)
2567 if (ir->swap->nstswap < 1)
2569 warning_error(wi, "swap_frequency must be 1 or larger when ion swapping is requested");
2571 if (ir->swap->nAverage < 1)
2573 warning_error(wi, "coupl_steps must be 1 or larger.\n");
2575 if (ir->swap->threshold < 1.0)
2577 warning_error(wi, "Ion count threshold must be at least 1.\n");
2581 sfree(dumstr[0]);
2582 sfree(dumstr[1]);
2585 static int search_QMstring(const char *s, int ng, const char *gn[])
2587 /* same as normal search_string, but this one searches QM strings */
2588 int i;
2590 for (i = 0; (i < ng); i++)
2592 if (gmx_strcasecmp(s, gn[i]) == 0)
2594 return i;
2598 gmx_fatal(FARGS, "this QM method or basisset (%s) is not implemented\n!", s);
2599 } /* search_QMstring */
2601 /* We would like gn to be const as well, but C doesn't allow this */
2602 /* TODO this is utility functionality (search for the index of a
2603 string in a collection), so should be refactored and located more
2604 centrally. */
2605 int search_string(const char *s, int ng, char *gn[])
2607 int i;
2609 for (i = 0; (i < ng); i++)
2611 if (gmx_strcasecmp(s, gn[i]) == 0)
2613 return i;
2617 gmx_fatal(FARGS,
2618 "Group %s referenced in the .mdp file was not found in the index file.\n"
2619 "Group names must match either [moleculetype] names or custom index group\n"
2620 "names, in which case you must supply an index file to the '-n' option\n"
2621 "of grompp.",
2625 static bool do_numbering(int natoms, gmx_groups_t *groups,
2626 gmx::ArrayRef<std::string> groupsFromMdpFile,
2627 t_blocka *block, char *gnames[],
2628 int gtype, int restnm,
2629 int grptp, bool bVerbose,
2630 warninp_t wi)
2632 unsigned short *cbuf;
2633 t_grps *grps = &(groups->grps[gtype]);
2634 int j, gid, aj, ognr, ntot = 0;
2635 const char *title;
2636 bool bRest;
2637 char warn_buf[STRLEN];
2639 title = gtypes[gtype];
2641 snew(cbuf, natoms);
2642 /* Mark all id's as not set */
2643 for (int i = 0; (i < natoms); i++)
2645 cbuf[i] = NOGID;
2648 snew(grps->nm_ind, groupsFromMdpFile.size()+1); /* +1 for possible rest group */
2649 for (int i = 0; i != groupsFromMdpFile.size(); ++i)
2651 /* Lookup the group name in the block structure */
2652 gid = search_string(groupsFromMdpFile[i].c_str(), block->nr, gnames);
2653 if ((grptp != egrptpONE) || (i == 0))
2655 grps->nm_ind[grps->nr++] = gid;
2658 /* Now go over the atoms in the group */
2659 for (j = block->index[gid]; (j < block->index[gid+1]); j++)
2662 aj = block->a[j];
2664 /* Range checking */
2665 if ((aj < 0) || (aj >= natoms))
2667 gmx_fatal(FARGS, "Invalid atom number %d in indexfile", aj);
2669 /* Lookup up the old group number */
2670 ognr = cbuf[aj];
2671 if (ognr != NOGID)
2673 gmx_fatal(FARGS, "Atom %d in multiple %s groups (%d and %d)",
2674 aj+1, title, ognr+1, i+1);
2676 else
2678 /* Store the group number in buffer */
2679 if (grptp == egrptpONE)
2681 cbuf[aj] = 0;
2683 else
2685 cbuf[aj] = i;
2687 ntot++;
2692 /* Now check whether we have done all atoms */
2693 bRest = FALSE;
2694 if (ntot != natoms)
2696 if (grptp == egrptpALL)
2698 gmx_fatal(FARGS, "%d atoms are not part of any of the %s groups",
2699 natoms-ntot, title);
2701 else if (grptp == egrptpPART)
2703 sprintf(warn_buf, "%d atoms are not part of any of the %s groups",
2704 natoms-ntot, title);
2705 warning_note(wi, warn_buf);
2707 /* Assign all atoms currently unassigned to a rest group */
2708 for (j = 0; (j < natoms); j++)
2710 if (cbuf[j] == NOGID)
2712 cbuf[j] = grps->nr;
2713 bRest = TRUE;
2716 if (grptp != egrptpPART)
2718 if (bVerbose)
2720 fprintf(stderr,
2721 "Making dummy/rest group for %s containing %d elements\n",
2722 title, natoms-ntot);
2724 /* Add group name "rest" */
2725 grps->nm_ind[grps->nr] = restnm;
2727 /* Assign the rest name to all atoms not currently assigned to a group */
2728 for (j = 0; (j < natoms); j++)
2730 if (cbuf[j] == NOGID)
2732 cbuf[j] = grps->nr;
2735 grps->nr++;
2739 if (grps->nr == 1 && (ntot == 0 || ntot == natoms))
2741 /* All atoms are part of one (or no) group, no index required */
2742 groups->ngrpnr[gtype] = 0;
2743 groups->grpnr[gtype] = nullptr;
2745 else
2747 groups->ngrpnr[gtype] = natoms;
2748 snew(groups->grpnr[gtype], natoms);
2749 for (j = 0; (j < natoms); j++)
2751 groups->grpnr[gtype][j] = cbuf[j];
2755 sfree(cbuf);
2757 return (bRest && grptp == egrptpPART);
2760 static void calc_nrdf(const gmx_mtop_t *mtop, t_inputrec *ir, char **gnames)
2762 t_grpopts *opts;
2763 pull_params_t *pull;
2764 int natoms, ai, aj, i, j, d, g, imin, jmin;
2765 int *nrdf2, *na_vcm, na_tot;
2766 double *nrdf_tc, *nrdf_vcm, nrdf_uc, *nrdf_vcm_sub;
2767 ivec *dof_vcm;
2768 int mol, ftype, as;
2770 /* Calculate nrdf.
2771 * First calc 3xnr-atoms for each group
2772 * then subtract half a degree of freedom for each constraint
2774 * Only atoms and nuclei contribute to the degrees of freedom...
2777 opts = &ir->opts;
2779 const gmx_groups_t &groups = mtop->groups;
2780 natoms = mtop->natoms;
2782 /* Allocate one more for a possible rest group */
2783 /* We need to sum degrees of freedom into doubles,
2784 * since floats give too low nrdf's above 3 million atoms.
2786 snew(nrdf_tc, groups.grps[egcTC].nr+1);
2787 snew(nrdf_vcm, groups.grps[egcVCM].nr+1);
2788 snew(dof_vcm, groups.grps[egcVCM].nr+1);
2789 snew(na_vcm, groups.grps[egcVCM].nr+1);
2790 snew(nrdf_vcm_sub, groups.grps[egcVCM].nr+1);
2792 for (i = 0; i < groups.grps[egcTC].nr; i++)
2794 nrdf_tc[i] = 0;
2796 for (i = 0; i < groups.grps[egcVCM].nr+1; i++)
2798 nrdf_vcm[i] = 0;
2799 clear_ivec(dof_vcm[i]);
2800 na_vcm[i] = 0;
2801 nrdf_vcm_sub[i] = 0;
2804 snew(nrdf2, natoms);
2805 for (const AtomProxy &atomP : AtomRange(*mtop))
2807 const t_atom &local = atomP.atom();
2808 int i = atomP.globalAtomNumber();
2809 nrdf2[i] = 0;
2810 if (local.ptype == eptAtom || local.ptype == eptNucleus)
2812 g = getGroupType(groups, egcFREEZE, i);
2813 for (int d = 0; d < DIM; d++)
2815 if (opts->nFreeze[g][d] == 0)
2817 /* Add one DOF for particle i (counted as 2*1) */
2818 nrdf2[i] += 2;
2819 /* VCM group i has dim d as a DOF */
2820 dof_vcm[getGroupType(groups, egcVCM, i)][d] = 1;
2823 nrdf_tc [getGroupType(groups, egcTC, i)] += 0.5*nrdf2[i];
2824 nrdf_vcm[getGroupType(groups, egcVCM, i)] += 0.5*nrdf2[i];
2828 as = 0;
2829 for (const gmx_molblock_t &molb : mtop->molblock)
2831 const gmx_moltype_t &molt = mtop->moltype[molb.type];
2832 const t_atom *atom = molt.atoms.atom;
2833 for (mol = 0; mol < molb.nmol; mol++)
2835 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
2837 gmx::ArrayRef<const int> ia = molt.ilist[ftype].iatoms;
2838 for (i = 0; i < molt.ilist[ftype].size(); )
2840 /* Subtract degrees of freedom for the constraints,
2841 * if the particles still have degrees of freedom left.
2842 * If one of the particles is a vsite or a shell, then all
2843 * constraint motion will go there, but since they do not
2844 * contribute to the constraints the degrees of freedom do not
2845 * change.
2847 ai = as + ia[i + 1];
2848 aj = as + ia[i + 2];
2849 if (((atom[ia[i + 1]].ptype == eptNucleus) ||
2850 (atom[ia[i + 1]].ptype == eptAtom)) &&
2851 ((atom[ia[i + 2]].ptype == eptNucleus) ||
2852 (atom[ia[i + 2]].ptype == eptAtom)))
2854 if (nrdf2[ai] > 0)
2856 jmin = 1;
2858 else
2860 jmin = 2;
2862 if (nrdf2[aj] > 0)
2864 imin = 1;
2866 else
2868 imin = 2;
2870 imin = std::min(imin, nrdf2[ai]);
2871 jmin = std::min(jmin, nrdf2[aj]);
2872 nrdf2[ai] -= imin;
2873 nrdf2[aj] -= jmin;
2874 nrdf_tc [getGroupType(groups, egcTC, ai)] -= 0.5*imin;
2875 nrdf_tc [getGroupType(groups, egcTC, aj)] -= 0.5*jmin;
2876 nrdf_vcm[getGroupType(groups, egcVCM, ai)] -= 0.5*imin;
2877 nrdf_vcm[getGroupType(groups, egcVCM, aj)] -= 0.5*jmin;
2879 i += interaction_function[ftype].nratoms+1;
2882 gmx::ArrayRef<const int> ia = molt.ilist[F_SETTLE].iatoms;
2883 for (i = 0; i < molt.ilist[F_SETTLE].size(); )
2885 /* Subtract 1 dof from every atom in the SETTLE */
2886 for (j = 0; j < 3; j++)
2888 ai = as + ia[i + 1 + j];
2889 imin = std::min(2, nrdf2[ai]);
2890 nrdf2[ai] -= imin;
2891 nrdf_tc [getGroupType(groups, egcTC, ai)] -= 0.5*imin;
2892 nrdf_vcm[getGroupType(groups, egcVCM, ai)] -= 0.5*imin;
2894 i += 4;
2896 as += molt.atoms.nr;
2900 if (ir->bPull)
2902 /* Correct nrdf for the COM constraints.
2903 * We correct using the TC and VCM group of the first atom
2904 * in the reference and pull group. If atoms in one pull group
2905 * belong to different TC or VCM groups it is anyhow difficult
2906 * to determine the optimal nrdf assignment.
2908 pull = ir->pull;
2910 for (i = 0; i < pull->ncoord; i++)
2912 if (pull->coord[i].eType != epullCONSTRAINT)
2914 continue;
2917 imin = 1;
2919 for (j = 0; j < 2; j++)
2921 const t_pull_group *pgrp;
2923 pgrp = &pull->group[pull->coord[i].group[j]];
2925 if (pgrp->nat > 0)
2927 /* Subtract 1/2 dof from each group */
2928 ai = pgrp->ind[0];
2929 nrdf_tc [getGroupType(groups, egcTC, ai)] -= 0.5*imin;
2930 nrdf_vcm[getGroupType(groups, egcVCM, ai)] -= 0.5*imin;
2931 if (nrdf_tc[getGroupType(groups, egcTC, ai)] < 0)
2933 gmx_fatal(FARGS, "Center of mass pulling constraints caused the number of degrees of freedom for temperature coupling group %s to be negative", gnames[groups.grps[egcTC].nm_ind[getGroupType(groups, egcTC, ai)]]);
2936 else
2938 /* We need to subtract the whole DOF from group j=1 */
2939 imin += 1;
2945 if (ir->nstcomm != 0)
2947 int ndim_rm_vcm;
2949 /* We remove COM motion up to dim ndof_com() */
2950 ndim_rm_vcm = ndof_com(ir);
2952 /* Subtract ndim_rm_vcm (or less with frozen dimensions) from
2953 * the number of degrees of freedom in each vcm group when COM
2954 * translation is removed and 6 when rotation is removed as well.
2956 for (j = 0; j < groups.grps[egcVCM].nr+1; j++)
2958 switch (ir->comm_mode)
2960 case ecmLINEAR:
2961 case ecmLINEAR_ACCELERATION_CORRECTION:
2962 nrdf_vcm_sub[j] = 0;
2963 for (d = 0; d < ndim_rm_vcm; d++)
2965 if (dof_vcm[j][d])
2967 nrdf_vcm_sub[j]++;
2970 break;
2971 case ecmANGULAR:
2972 nrdf_vcm_sub[j] = 6;
2973 break;
2974 default:
2975 gmx_incons("Checking comm_mode");
2979 for (i = 0; i < groups.grps[egcTC].nr; i++)
2981 /* Count the number of atoms of TC group i for every VCM group */
2982 for (j = 0; j < groups.grps[egcVCM].nr+1; j++)
2984 na_vcm[j] = 0;
2986 na_tot = 0;
2987 for (ai = 0; ai < natoms; ai++)
2989 if (getGroupType(groups, egcTC, ai) == i)
2991 na_vcm[getGroupType(groups, egcVCM, ai)]++;
2992 na_tot++;
2995 /* Correct for VCM removal according to the fraction of each VCM
2996 * group present in this TC group.
2998 nrdf_uc = nrdf_tc[i];
2999 nrdf_tc[i] = 0;
3000 for (j = 0; j < groups.grps[egcVCM].nr+1; j++)
3002 if (nrdf_vcm[j] > nrdf_vcm_sub[j])
3004 nrdf_tc[i] += nrdf_uc*(static_cast<double>(na_vcm[j])/static_cast<double>(na_tot))*
3005 (nrdf_vcm[j] - nrdf_vcm_sub[j])/nrdf_vcm[j];
3010 for (i = 0; (i < groups.grps[egcTC].nr); i++)
3012 opts->nrdf[i] = nrdf_tc[i];
3013 if (opts->nrdf[i] < 0)
3015 opts->nrdf[i] = 0;
3017 fprintf(stderr,
3018 "Number of degrees of freedom in T-Coupling group %s is %.2f\n",
3019 gnames[groups.grps[egcTC].nm_ind[i]], opts->nrdf[i]);
3022 sfree(nrdf2);
3023 sfree(nrdf_tc);
3024 sfree(nrdf_vcm);
3025 sfree(dof_vcm);
3026 sfree(na_vcm);
3027 sfree(nrdf_vcm_sub);
3030 static bool do_egp_flag(t_inputrec *ir, gmx_groups_t *groups,
3031 const char *option, const char *val, int flag)
3033 /* The maximum number of energy group pairs would be MAXPTR*(MAXPTR+1)/2.
3034 * But since this is much larger than STRLEN, such a line can not be parsed.
3035 * The real maximum is the number of names that fit in a string: STRLEN/2.
3037 #define EGP_MAX (STRLEN/2)
3038 int j, k, nr;
3039 char ***gnames;
3040 bool bSet;
3042 gnames = groups->grpname;
3044 auto names = gmx::splitString(val);
3045 if (names.size() % 2 != 0)
3047 gmx_fatal(FARGS, "The number of groups for %s is odd", option);
3049 nr = groups->grps[egcENER].nr;
3050 bSet = FALSE;
3051 for (size_t i = 0; i < names.size() / 2; i++)
3053 j = 0;
3054 while ((j < nr) &&
3055 gmx_strcasecmp(names[2*i].c_str(), *(gnames[groups->grps[egcENER].nm_ind[j]])))
3057 j++;
3059 if (j == nr)
3061 gmx_fatal(FARGS, "%s in %s is not an energy group\n",
3062 names[2*i].c_str(), option);
3064 k = 0;
3065 while ((k < nr) &&
3066 gmx_strcasecmp(names[2*i+1].c_str(), *(gnames[groups->grps[egcENER].nm_ind[k]])))
3068 k++;
3070 if (k == nr)
3072 gmx_fatal(FARGS, "%s in %s is not an energy group\n",
3073 names[2*i+1].c_str(), option);
3075 if ((j < nr) && (k < nr))
3077 ir->opts.egp_flags[nr*j+k] |= flag;
3078 ir->opts.egp_flags[nr*k+j] |= flag;
3079 bSet = TRUE;
3083 return bSet;
3087 static void make_swap_groups(
3088 t_swapcoords *swap,
3089 t_blocka *grps,
3090 char **gnames)
3092 int ig = -1, i = 0, gind;
3093 t_swapGroup *swapg;
3096 /* Just a quick check here, more thorough checks are in mdrun */
3097 if (strcmp(swap->grp[eGrpSplit0].molname, swap->grp[eGrpSplit1].molname) == 0)
3099 gmx_fatal(FARGS, "The split groups can not both be '%s'.", swap->grp[eGrpSplit0].molname);
3102 /* Get the index atoms of the split0, split1, solvent, and swap groups */
3103 for (ig = 0; ig < swap->ngrp; ig++)
3105 swapg = &swap->grp[ig];
3106 gind = search_string(swap->grp[ig].molname, grps->nr, gnames);
3107 swapg->nat = grps->index[gind+1] - grps->index[gind];
3109 if (swapg->nat > 0)
3111 fprintf(stderr, "%s group '%s' contains %d atoms.\n",
3112 ig < 3 ? eSwapFixedGrp_names[ig] : "Swap",
3113 swap->grp[ig].molname, swapg->nat);
3114 snew(swapg->ind, swapg->nat);
3115 for (i = 0; i < swapg->nat; i++)
3117 swapg->ind[i] = grps->a[grps->index[gind]+i];
3120 else
3122 gmx_fatal(FARGS, "Swap group %s does not contain any atoms.", swap->grp[ig].molname);
3128 static void make_IMD_group(t_IMD *IMDgroup, char *IMDgname, t_blocka *grps, char **gnames)
3130 int ig, i;
3133 ig = search_string(IMDgname, grps->nr, gnames);
3134 IMDgroup->nat = grps->index[ig+1] - grps->index[ig];
3136 if (IMDgroup->nat > 0)
3138 fprintf(stderr, "Group '%s' with %d atoms can be activated for interactive molecular dynamics (IMD).\n",
3139 IMDgname, IMDgroup->nat);
3140 snew(IMDgroup->ind, IMDgroup->nat);
3141 for (i = 0; i < IMDgroup->nat; i++)
3143 IMDgroup->ind[i] = grps->a[grps->index[ig]+i];
3148 void do_index(const char* mdparin, const char *ndx,
3149 gmx_mtop_t *mtop,
3150 bool bVerbose,
3151 t_inputrec *ir,
3152 warninp_t wi)
3154 t_blocka *grps;
3155 gmx_groups_t *groups;
3156 int natoms;
3157 t_symtab *symtab;
3158 t_atoms atoms_all;
3159 char warnbuf[STRLEN], **gnames;
3160 int nr;
3161 real tau_min;
3162 int nstcmin;
3163 int i, j, k, restnm;
3164 bool bExcl, bTable, bAnneal, bRest;
3165 char warn_buf[STRLEN];
3167 if (bVerbose)
3169 fprintf(stderr, "processing index file...\n");
3171 if (ndx == nullptr)
3173 snew(grps, 1);
3174 snew(grps->index, 1);
3175 snew(gnames, 1);
3176 atoms_all = gmx_mtop_global_atoms(mtop);
3177 analyse(&atoms_all, grps, &gnames, FALSE, TRUE);
3178 done_atom(&atoms_all);
3180 else
3182 grps = init_index(ndx, &gnames);
3185 groups = &mtop->groups;
3186 natoms = mtop->natoms;
3187 symtab = &mtop->symtab;
3189 snew(groups->grpname, grps->nr+1);
3191 for (i = 0; (i < grps->nr); i++)
3193 groups->grpname[i] = put_symtab(symtab, gnames[i]);
3195 groups->grpname[i] = put_symtab(symtab, "rest");
3196 restnm = i;
3197 srenew(gnames, grps->nr+1);
3198 gnames[restnm] = *(groups->grpname[i]);
3199 groups->ngrpname = grps->nr+1;
3201 set_warning_line(wi, mdparin, -1);
3203 auto temperatureCouplingTauValues = gmx::splitString(is->tau_t);
3204 auto temperatureCouplingReferenceValues = gmx::splitString(is->ref_t);
3205 auto temperatureCouplingGroupNames = gmx::splitString(is->tcgrps);
3206 if (temperatureCouplingTauValues.size() != temperatureCouplingGroupNames.size() ||
3207 temperatureCouplingReferenceValues.size() != temperatureCouplingGroupNames.size())
3209 gmx_fatal(FARGS, "Invalid T coupling input: %zu groups, %zu ref-t values and "
3210 "%zu tau-t values",
3211 temperatureCouplingGroupNames.size(),
3212 temperatureCouplingReferenceValues.size(),
3213 temperatureCouplingTauValues.size());
3216 const bool useReferenceTemperature = integratorHasReferenceTemperature(ir);
3217 do_numbering(natoms, groups, temperatureCouplingGroupNames, grps, gnames, egcTC,
3218 restnm, useReferenceTemperature ? egrptpALL : egrptpALL_GENREST, bVerbose, wi);
3219 nr = groups->grps[egcTC].nr;
3220 ir->opts.ngtc = nr;
3221 snew(ir->opts.nrdf, nr);
3222 snew(ir->opts.tau_t, nr);
3223 snew(ir->opts.ref_t, nr);
3224 if (ir->eI == eiBD && ir->bd_fric == 0)
3226 fprintf(stderr, "bd-fric=0, so tau-t will be used as the inverse friction constant(s)\n");
3229 if (useReferenceTemperature)
3231 if (size_t(nr) != temperatureCouplingReferenceValues.size())
3233 gmx_fatal(FARGS, "Not enough ref-t and tau-t values!");
3236 tau_min = 1e20;
3237 convertReals(wi, temperatureCouplingTauValues, "tau-t", ir->opts.tau_t);
3238 for (i = 0; (i < nr); i++)
3240 if ((ir->eI == eiBD) && ir->opts.tau_t[i] <= 0)
3242 sprintf(warn_buf, "With integrator %s tau-t should be larger than 0", ei_names[ir->eI]);
3243 warning_error(wi, warn_buf);
3246 if (ir->etc != etcVRESCALE && ir->opts.tau_t[i] == 0)
3248 warning_note(wi, "tau-t = -1 is the value to signal that a group should not have temperature coupling. Treating your use of tau-t = 0 as if you used -1.");
3251 if (ir->opts.tau_t[i] >= 0)
3253 tau_min = std::min(tau_min, ir->opts.tau_t[i]);
3256 if (ir->etc != etcNO && ir->nsttcouple == -1)
3258 ir->nsttcouple = ir_optimal_nsttcouple(ir);
3261 if (EI_VV(ir->eI))
3263 if ((ir->etc == etcNOSEHOOVER) && (ir->epc == epcBERENDSEN))
3265 gmx_fatal(FARGS, "Cannot do Nose-Hoover temperature with Berendsen pressure control with md-vv; use either vrescale temperature with berendsen pressure or Nose-Hoover temperature with MTTK pressure");
3267 if (ir->epc == epcMTTK)
3269 if (ir->etc != etcNOSEHOOVER)
3271 gmx_fatal(FARGS, "Cannot do MTTK pressure coupling without Nose-Hoover temperature control");
3273 else
3275 if (ir->nstpcouple != ir->nsttcouple)
3277 int mincouple = std::min(ir->nstpcouple, ir->nsttcouple);
3278 ir->nstpcouple = ir->nsttcouple = mincouple;
3279 sprintf(warn_buf, "for current Trotter decomposition methods with vv, nsttcouple and nstpcouple must be equal. Both have been reset to min(nsttcouple,nstpcouple) = %d", mincouple);
3280 warning_note(wi, warn_buf);
3285 /* velocity verlet with averaged kinetic energy KE = 0.5*(v(t+1/2) - v(t-1/2)) is implemented
3286 primarily for testing purposes, and does not work with temperature coupling other than 1 */
3288 if (ETC_ANDERSEN(ir->etc))
3290 if (ir->nsttcouple != 1)
3292 ir->nsttcouple = 1;
3293 sprintf(warn_buf, "Andersen temperature control methods assume nsttcouple = 1; there is no need for larger nsttcouple > 1, since no global parameters are computed. nsttcouple has been reset to 1");
3294 warning_note(wi, warn_buf);
3297 nstcmin = tcouple_min_integration_steps(ir->etc);
3298 if (nstcmin > 1)
3300 if (tau_min/(ir->delta_t*ir->nsttcouple) < nstcmin - 10*GMX_REAL_EPS)
3302 sprintf(warn_buf, "For proper integration of the %s thermostat, tau-t (%g) should be at least %d times larger than nsttcouple*dt (%g)",
3303 ETCOUPLTYPE(ir->etc),
3304 tau_min, nstcmin,
3305 ir->nsttcouple*ir->delta_t);
3306 warning(wi, warn_buf);
3309 convertReals(wi, temperatureCouplingReferenceValues, "ref-t", ir->opts.ref_t);
3310 for (i = 0; (i < nr); i++)
3312 if (ir->opts.ref_t[i] < 0)
3314 gmx_fatal(FARGS, "ref-t for group %d negative", i);
3317 /* set the lambda mc temperature to the md integrator temperature (which should be defined
3318 if we are in this conditional) if mc_temp is negative */
3319 if (ir->expandedvals->mc_temp < 0)
3321 ir->expandedvals->mc_temp = ir->opts.ref_t[0]; /*for now, set to the first reft */
3325 /* Simulated annealing for each group. There are nr groups */
3326 auto simulatedAnnealingGroupNames = gmx::splitString(is->anneal);
3327 if (simulatedAnnealingGroupNames.size() == 1 &&
3328 gmx_strncasecmp(simulatedAnnealingGroupNames[0].c_str(), "N", 1) == 0)
3330 simulatedAnnealingGroupNames.resize(0);
3332 if (!simulatedAnnealingGroupNames.empty() &&
3333 simulatedAnnealingGroupNames.size() != size_t(nr))
3335 gmx_fatal(FARGS, "Not enough annealing values: %zu (for %d groups)\n",
3336 simulatedAnnealingGroupNames.size(), nr);
3338 else
3340 snew(ir->opts.annealing, nr);
3341 snew(ir->opts.anneal_npoints, nr);
3342 snew(ir->opts.anneal_time, nr);
3343 snew(ir->opts.anneal_temp, nr);
3344 for (i = 0; i < nr; i++)
3346 ir->opts.annealing[i] = eannNO;
3347 ir->opts.anneal_npoints[i] = 0;
3348 ir->opts.anneal_time[i] = nullptr;
3349 ir->opts.anneal_temp[i] = nullptr;
3351 if (!simulatedAnnealingGroupNames.empty())
3353 bAnneal = FALSE;
3354 for (i = 0; i < nr; i++)
3356 if (gmx_strncasecmp(simulatedAnnealingGroupNames[i].c_str(), "N", 1) == 0)
3358 ir->opts.annealing[i] = eannNO;
3360 else if (gmx_strncasecmp(simulatedAnnealingGroupNames[i].c_str(), "S", 1) == 0)
3362 ir->opts.annealing[i] = eannSINGLE;
3363 bAnneal = TRUE;
3365 else if (gmx_strncasecmp(simulatedAnnealingGroupNames[i].c_str(), "P", 1) == 0)
3367 ir->opts.annealing[i] = eannPERIODIC;
3368 bAnneal = TRUE;
3371 if (bAnneal)
3373 /* Read the other fields too */
3374 auto simulatedAnnealingPoints = gmx::splitString(is->anneal_npoints);
3375 if (simulatedAnnealingPoints.size() != simulatedAnnealingGroupNames.size())
3377 gmx_fatal(FARGS, "Found %zu annealing-npoints values for %zu groups\n",
3378 simulatedAnnealingPoints.size(), simulatedAnnealingGroupNames.size());
3380 convertInts(wi, simulatedAnnealingPoints, "annealing points", ir->opts.anneal_npoints);
3381 for (k = 0, i = 0; i < nr; i++)
3383 if (ir->opts.anneal_npoints[i] == 1)
3385 gmx_fatal(FARGS, "Please specify at least a start and an end point for annealing\n");
3387 snew(ir->opts.anneal_time[i], ir->opts.anneal_npoints[i]);
3388 snew(ir->opts.anneal_temp[i], ir->opts.anneal_npoints[i]);
3389 k += ir->opts.anneal_npoints[i];
3392 auto simulatedAnnealingTimes = gmx::splitString(is->anneal_time);
3393 if (simulatedAnnealingTimes.size() != size_t(k))
3395 gmx_fatal(FARGS, "Found %zu annealing-time values, wanted %d\n",
3396 simulatedAnnealingTimes.size(), k);
3398 auto simulatedAnnealingTemperatures = gmx::splitString(is->anneal_temp);
3399 if (simulatedAnnealingTemperatures.size() != size_t(k))
3401 gmx_fatal(FARGS, "Found %zu annealing-temp values, wanted %d\n",
3402 simulatedAnnealingTemperatures.size(), k);
3405 convertReals(wi, simulatedAnnealingTimes, "anneal-time", ir->opts.anneal_time[i]);
3406 convertReals(wi, simulatedAnnealingTemperatures, "anneal-temp", ir->opts.anneal_temp[i]);
3407 for (i = 0, k = 0; i < nr; i++)
3409 for (j = 0; j < ir->opts.anneal_npoints[i]; j++)
3411 if (j == 0)
3413 if (ir->opts.anneal_time[i][0] > (ir->init_t+GMX_REAL_EPS))
3415 gmx_fatal(FARGS, "First time point for annealing > init_t.\n");
3418 else
3420 /* j>0 */
3421 if (ir->opts.anneal_time[i][j] < ir->opts.anneal_time[i][j-1])
3423 gmx_fatal(FARGS, "Annealing timepoints out of order: t=%f comes after t=%f\n",
3424 ir->opts.anneal_time[i][j], ir->opts.anneal_time[i][j-1]);
3427 if (ir->opts.anneal_temp[i][j] < 0)
3429 gmx_fatal(FARGS, "Found negative temperature in annealing: %f\n", ir->opts.anneal_temp[i][j]);
3431 k++;
3434 /* Print out some summary information, to make sure we got it right */
3435 for (i = 0, k = 0; i < nr; i++)
3437 if (ir->opts.annealing[i] != eannNO)
3439 j = groups->grps[egcTC].nm_ind[i];
3440 fprintf(stderr, "Simulated annealing for group %s: %s, %d timepoints\n",
3441 *(groups->grpname[j]), eann_names[ir->opts.annealing[i]],
3442 ir->opts.anneal_npoints[i]);
3443 fprintf(stderr, "Time (ps) Temperature (K)\n");
3444 /* All terms except the last one */
3445 for (j = 0; j < (ir->opts.anneal_npoints[i]-1); j++)
3447 fprintf(stderr, "%9.1f %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3450 /* Finally the last one */
3451 j = ir->opts.anneal_npoints[i]-1;
3452 if (ir->opts.annealing[i] == eannSINGLE)
3454 fprintf(stderr, "%9.1f- %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3456 else
3458 fprintf(stderr, "%9.1f %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3459 if (std::fabs(ir->opts.anneal_temp[i][j]-ir->opts.anneal_temp[i][0]) > GMX_REAL_EPS)
3461 warning_note(wi, "There is a temperature jump when your annealing loops back.\n");
3470 if (ir->bPull)
3472 make_pull_groups(ir->pull, is->pull_grp, grps, gnames);
3474 make_pull_coords(ir->pull);
3477 if (ir->bRot)
3479 make_rotation_groups(ir->rot, is->rot_grp, grps, gnames);
3482 if (ir->eSwapCoords != eswapNO)
3484 make_swap_groups(ir->swap, grps, gnames);
3487 /* Make indices for IMD session */
3488 if (ir->bIMD)
3490 make_IMD_group(ir->imd, is->imd_grp, grps, gnames);
3493 auto accelerations = gmx::splitString(is->acc);
3494 auto accelerationGroupNames = gmx::splitString(is->accgrps);
3495 if (accelerationGroupNames.size() * DIM != accelerations.size())
3497 gmx_fatal(FARGS, "Invalid Acceleration input: %zu groups and %zu acc. values",
3498 accelerationGroupNames.size(), accelerations.size());
3500 do_numbering(natoms, groups, accelerationGroupNames, grps, gnames, egcACC,
3501 restnm, egrptpALL_GENREST, bVerbose, wi);
3502 nr = groups->grps[egcACC].nr;
3503 snew(ir->opts.acc, nr);
3504 ir->opts.ngacc = nr;
3506 convertRvecs(wi, accelerations, "anneal-time", ir->opts.acc);
3508 auto freezeDims = gmx::splitString(is->frdim);
3509 auto freezeGroupNames = gmx::splitString(is->freeze);
3510 if (freezeDims.size() != DIM * freezeGroupNames.size())
3512 gmx_fatal(FARGS, "Invalid Freezing input: %zu groups and %zu freeze values",
3513 freezeGroupNames.size(), freezeDims.size());
3515 do_numbering(natoms, groups, freezeGroupNames, grps, gnames, egcFREEZE,
3516 restnm, egrptpALL_GENREST, bVerbose, wi);
3517 nr = groups->grps[egcFREEZE].nr;
3518 ir->opts.ngfrz = nr;
3519 snew(ir->opts.nFreeze, nr);
3520 for (i = k = 0; (size_t(i) < freezeGroupNames.size()); i++)
3522 for (j = 0; (j < DIM); j++, k++)
3524 ir->opts.nFreeze[i][j] = static_cast<int>(gmx_strncasecmp(freezeDims[k].c_str(), "Y", 1) == 0);
3525 if (!ir->opts.nFreeze[i][j])
3527 if (gmx_strncasecmp(freezeDims[k].c_str(), "N", 1) != 0)
3529 sprintf(warnbuf, "Please use Y(ES) or N(O) for freezedim only "
3530 "(not %s)", freezeDims[k].c_str());
3531 warning(wi, warn_buf);
3536 for (; (i < nr); i++)
3538 for (j = 0; (j < DIM); j++)
3540 ir->opts.nFreeze[i][j] = 0;
3544 auto energyGroupNames = gmx::splitString(is->energy);
3545 do_numbering(natoms, groups, energyGroupNames, grps, gnames, egcENER,
3546 restnm, egrptpALL_GENREST, bVerbose, wi);
3547 add_wall_energrps(groups, ir->nwall, symtab);
3548 ir->opts.ngener = groups->grps[egcENER].nr;
3549 auto vcmGroupNames = gmx::splitString(is->vcm);
3550 bRest =
3551 do_numbering(natoms, groups, vcmGroupNames, grps, gnames, egcVCM,
3552 restnm, vcmGroupNames.empty() ? egrptpALL_GENREST : egrptpPART, bVerbose, wi);
3553 if (bRest)
3555 warning(wi, "Some atoms are not part of any center of mass motion removal group.\n"
3556 "This may lead to artifacts.\n"
3557 "In most cases one should use one group for the whole system.");
3560 /* Now we have filled the freeze struct, so we can calculate NRDF */
3561 calc_nrdf(mtop, ir, gnames);
3563 auto user1GroupNames = gmx::splitString(is->user1);
3564 do_numbering(natoms, groups, user1GroupNames, grps, gnames, egcUser1,
3565 restnm, egrptpALL_GENREST, bVerbose, wi);
3566 auto user2GroupNames = gmx::splitString(is->user2);
3567 do_numbering(natoms, groups, user2GroupNames, grps, gnames, egcUser2,
3568 restnm, egrptpALL_GENREST, bVerbose, wi);
3569 auto compressedXGroupNames = gmx::splitString(is->x_compressed_groups);
3570 do_numbering(natoms, groups, compressedXGroupNames, grps, gnames, egcCompressedX,
3571 restnm, egrptpONE, bVerbose, wi);
3572 auto orirefFitGroupNames = gmx::splitString(is->orirefitgrp);
3573 do_numbering(natoms, groups, orirefFitGroupNames, grps, gnames, egcORFIT,
3574 restnm, egrptpALL_GENREST, bVerbose, wi);
3576 /* QMMM input processing */
3577 auto qmGroupNames = gmx::splitString(is->QMMM);
3578 auto qmMethods = gmx::splitString(is->QMmethod);
3579 auto qmBasisSets = gmx::splitString(is->QMbasis);
3580 if (ir->eI != eiMimic)
3582 if (qmMethods.size() != qmGroupNames.size() ||
3583 qmBasisSets.size() != qmGroupNames.size())
3585 gmx_fatal(FARGS, "Invalid QMMM input: %zu groups %zu basissets"
3586 " and %zu methods\n", qmGroupNames.size(),
3587 qmBasisSets.size(), qmMethods.size());
3589 /* group rest, if any, is always MM! */
3590 do_numbering(natoms, groups, qmGroupNames, grps, gnames, egcQMMM,
3591 restnm, egrptpALL_GENREST, bVerbose, wi);
3592 nr = qmGroupNames.size(); /*atoms->grps[egcQMMM].nr;*/
3593 ir->opts.ngQM = qmGroupNames.size();
3594 snew(ir->opts.QMmethod, nr);
3595 snew(ir->opts.QMbasis, nr);
3596 for (i = 0; i < nr; i++)
3598 /* input consists of strings: RHF CASSCF PM3 .. These need to be
3599 * converted to the corresponding enum in names.c
3601 ir->opts.QMmethod[i] = search_QMstring(qmMethods[i].c_str(),
3602 eQMmethodNR,
3603 eQMmethod_names);
3604 ir->opts.QMbasis[i] = search_QMstring(qmBasisSets[i].c_str(),
3605 eQMbasisNR,
3606 eQMbasis_names);
3609 auto qmMultiplicities = gmx::splitString(is->QMmult);
3610 auto qmCharges = gmx::splitString(is->QMcharge);
3611 auto qmbSH = gmx::splitString(is->bSH);
3612 snew(ir->opts.QMmult, nr);
3613 snew(ir->opts.QMcharge, nr);
3614 snew(ir->opts.bSH, nr);
3615 convertInts(wi, qmMultiplicities, "QMmult", ir->opts.QMmult);
3616 convertInts(wi, qmCharges, "QMcharge", ir->opts.QMcharge);
3617 convertYesNos(wi, qmbSH, "bSH", ir->opts.bSH);
3619 auto CASelectrons = gmx::splitString(is->CASelectrons);
3620 auto CASorbitals = gmx::splitString(is->CASorbitals);
3621 snew(ir->opts.CASelectrons, nr);
3622 snew(ir->opts.CASorbitals, nr);
3623 convertInts(wi, CASelectrons, "CASelectrons", ir->opts.CASelectrons);
3624 convertInts(wi, CASorbitals, "CASOrbitals", ir->opts.CASorbitals);
3626 auto SAon = gmx::splitString(is->SAon);
3627 auto SAoff = gmx::splitString(is->SAoff);
3628 auto SAsteps = gmx::splitString(is->SAsteps);
3629 snew(ir->opts.SAon, nr);
3630 snew(ir->opts.SAoff, nr);
3631 snew(ir->opts.SAsteps, nr);
3632 convertInts(wi, SAon, "SAon", ir->opts.SAon);
3633 convertInts(wi, SAoff, "SAoff", ir->opts.SAoff);
3634 convertInts(wi, SAsteps, "SAsteps", ir->opts.SAsteps);
3636 else
3638 /* MiMiC */
3639 if (qmGroupNames.size() > 1)
3641 gmx_fatal(FARGS, "Currently, having more than one QM group in MiMiC is not supported");
3643 /* group rest, if any, is always MM! */
3644 do_numbering(natoms, groups, qmGroupNames, grps, gnames, egcQMMM,
3645 restnm, egrptpALL_GENREST, bVerbose, wi);
3647 ir->opts.ngQM = qmGroupNames.size();
3650 /* end of QMMM input */
3652 if (bVerbose)
3654 for (i = 0; (i < egcNR); i++)
3656 fprintf(stderr, "%-16s has %d element(s):", gtypes[i], groups->grps[i].nr);
3657 for (j = 0; (j < groups->grps[i].nr); j++)
3659 fprintf(stderr, " %s", *(groups->grpname[groups->grps[i].nm_ind[j]]));
3661 fprintf(stderr, "\n");
3665 nr = groups->grps[egcENER].nr;
3666 snew(ir->opts.egp_flags, nr*nr);
3668 bExcl = do_egp_flag(ir, groups, "energygrp-excl", is->egpexcl, EGP_EXCL);
3669 if (bExcl && ir->cutoff_scheme == ecutsVERLET)
3671 warning_error(wi, "Energy group exclusions are not (yet) implemented for the Verlet scheme");
3673 if (bExcl && EEL_FULL(ir->coulombtype))
3675 warning(wi, "Can not exclude the lattice Coulomb energy between energy groups");
3678 bTable = do_egp_flag(ir, groups, "energygrp-table", is->egptable, EGP_TABLE);
3679 if (bTable && !(ir->vdwtype == evdwUSER) &&
3680 !(ir->coulombtype == eelUSER) && !(ir->coulombtype == eelPMEUSER) &&
3681 !(ir->coulombtype == eelPMEUSERSWITCH))
3683 gmx_fatal(FARGS, "Can only have energy group pair tables in combination with user tables for VdW and/or Coulomb");
3686 for (i = 0; (i < grps->nr); i++)
3688 sfree(gnames[i]);
3690 sfree(gnames);
3691 done_blocka(grps);
3692 sfree(grps);
3698 static void check_disre(const gmx_mtop_t *mtop)
3700 if (gmx_mtop_ftype_count(mtop, F_DISRES) > 0)
3702 const gmx_ffparams_t &ffparams = mtop->ffparams;
3703 int ndouble = 0;
3704 int old_label = -1;
3705 for (int i = 0; i < ffparams.numTypes(); i++)
3707 int ftype = ffparams.functype[i];
3708 if (ftype == F_DISRES)
3710 int label = ffparams.iparams[i].disres.label;
3711 if (label == old_label)
3713 fprintf(stderr, "Distance restraint index %d occurs twice\n", label);
3714 ndouble++;
3716 old_label = label;
3719 if (ndouble > 0)
3721 gmx_fatal(FARGS, "Found %d double distance restraint indices,\n"
3722 "probably the parameters for multiple pairs in one restraint "
3723 "are not identical\n", ndouble);
3728 static bool absolute_reference(t_inputrec *ir, gmx_mtop_t *sys,
3729 bool posres_only,
3730 ivec AbsRef)
3732 int d, g, i;
3733 gmx_mtop_ilistloop_t iloop;
3734 int nmol;
3735 t_iparams *pr;
3737 clear_ivec(AbsRef);
3739 if (!posres_only)
3741 /* Check the COM */
3742 for (d = 0; d < DIM; d++)
3744 AbsRef[d] = (d < ndof_com(ir) ? 0 : 1);
3746 /* Check for freeze groups */
3747 for (g = 0; g < ir->opts.ngfrz; g++)
3749 for (d = 0; d < DIM; d++)
3751 if (ir->opts.nFreeze[g][d] != 0)
3753 AbsRef[d] = 1;
3759 /* Check for position restraints */
3760 iloop = gmx_mtop_ilistloop_init(sys);
3761 while (const InteractionLists *ilist = gmx_mtop_ilistloop_next(iloop, &nmol))
3763 if (nmol > 0 &&
3764 (AbsRef[XX] == 0 || AbsRef[YY] == 0 || AbsRef[ZZ] == 0))
3766 for (i = 0; i < (*ilist)[F_POSRES].size(); i += 2)
3768 pr = &sys->ffparams.iparams[(*ilist)[F_POSRES].iatoms[i]];
3769 for (d = 0; d < DIM; d++)
3771 if (pr->posres.fcA[d] != 0)
3773 AbsRef[d] = 1;
3777 for (i = 0; i < (*ilist)[F_FBPOSRES].size(); i += 2)
3779 /* Check for flat-bottom posres */
3780 pr = &sys->ffparams.iparams[(*ilist)[F_FBPOSRES].iatoms[i]];
3781 if (pr->fbposres.k != 0)
3783 switch (pr->fbposres.geom)
3785 case efbposresSPHERE:
3786 AbsRef[XX] = AbsRef[YY] = AbsRef[ZZ] = 1;
3787 break;
3788 case efbposresCYLINDERX:
3789 AbsRef[YY] = AbsRef[ZZ] = 1;
3790 break;
3791 case efbposresCYLINDERY:
3792 AbsRef[XX] = AbsRef[ZZ] = 1;
3793 break;
3794 case efbposresCYLINDER:
3795 /* efbposres is a synonym for efbposresCYLINDERZ for backwards compatibility */
3796 case efbposresCYLINDERZ:
3797 AbsRef[XX] = AbsRef[YY] = 1;
3798 break;
3799 case efbposresX: /* d=XX */
3800 case efbposresY: /* d=YY */
3801 case efbposresZ: /* d=ZZ */
3802 d = pr->fbposres.geom - efbposresX;
3803 AbsRef[d] = 1;
3804 break;
3805 default:
3806 gmx_fatal(FARGS, " Invalid geometry for flat-bottom position restraint.\n"
3807 "Expected nr between 1 and %d. Found %d\n", efbposresNR-1,
3808 pr->fbposres.geom);
3815 return (AbsRef[XX] != 0 && AbsRef[YY] != 0 && AbsRef[ZZ] != 0);
3818 static void
3819 check_combination_rule_differences(const gmx_mtop_t *mtop, int state,
3820 bool *bC6ParametersWorkWithGeometricRules,
3821 bool *bC6ParametersWorkWithLBRules,
3822 bool *bLBRulesPossible)
3824 int ntypes, tpi, tpj;
3825 int *typecount;
3826 real tol;
3827 double c6i, c6j, c12i, c12j;
3828 double c6, c6_geometric, c6_LB;
3829 double sigmai, sigmaj, epsi, epsj;
3830 bool bCanDoLBRules, bCanDoGeometricRules;
3831 const char *ptr;
3833 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
3834 * force-field floating point parameters.
3836 tol = 1e-5;
3837 ptr = getenv("GMX_LJCOMB_TOL");
3838 if (ptr != nullptr)
3840 double dbl;
3841 double gmx_unused canary;
3843 if (sscanf(ptr, "%lf%lf", &dbl, &canary) != 1)
3845 gmx_fatal(FARGS, "Could not parse a single floating-point number from GMX_LJCOMB_TOL (%s)", ptr);
3847 tol = dbl;
3850 *bC6ParametersWorkWithLBRules = TRUE;
3851 *bC6ParametersWorkWithGeometricRules = TRUE;
3852 bCanDoLBRules = TRUE;
3853 ntypes = mtop->ffparams.atnr;
3854 snew(typecount, ntypes);
3855 gmx_mtop_count_atomtypes(mtop, state, typecount);
3856 *bLBRulesPossible = TRUE;
3857 for (tpi = 0; tpi < ntypes; ++tpi)
3859 c6i = mtop->ffparams.iparams[(ntypes + 1) * tpi].lj.c6;
3860 c12i = mtop->ffparams.iparams[(ntypes + 1) * tpi].lj.c12;
3861 for (tpj = tpi; tpj < ntypes; ++tpj)
3863 c6j = mtop->ffparams.iparams[(ntypes + 1) * tpj].lj.c6;
3864 c12j = mtop->ffparams.iparams[(ntypes + 1) * tpj].lj.c12;
3865 c6 = mtop->ffparams.iparams[ntypes * tpi + tpj].lj.c6;
3866 c6_geometric = std::sqrt(c6i * c6j);
3867 if (!gmx_numzero(c6_geometric))
3869 if (!gmx_numzero(c12i) && !gmx_numzero(c12j))
3871 sigmai = gmx::sixthroot(c12i / c6i);
3872 sigmaj = gmx::sixthroot(c12j / c6j);
3873 epsi = c6i * c6i /(4.0 * c12i);
3874 epsj = c6j * c6j /(4.0 * c12j);
3875 c6_LB = 4.0 * std::sqrt(epsi * epsj) * gmx::power6(0.5 * (sigmai + sigmaj));
3877 else
3879 *bLBRulesPossible = FALSE;
3880 c6_LB = c6_geometric;
3882 bCanDoLBRules = gmx_within_tol(c6_LB, c6, tol);
3885 if (!bCanDoLBRules)
3887 *bC6ParametersWorkWithLBRules = FALSE;
3890 bCanDoGeometricRules = gmx_within_tol(c6_geometric, c6, tol);
3892 if (!bCanDoGeometricRules)
3894 *bC6ParametersWorkWithGeometricRules = FALSE;
3898 sfree(typecount);
3901 static void
3902 check_combination_rules(const t_inputrec *ir, const gmx_mtop_t *mtop,
3903 warninp_t wi)
3905 bool bLBRulesPossible, bC6ParametersWorkWithGeometricRules, bC6ParametersWorkWithLBRules;
3907 check_combination_rule_differences(mtop, 0,
3908 &bC6ParametersWorkWithGeometricRules,
3909 &bC6ParametersWorkWithLBRules,
3910 &bLBRulesPossible);
3911 if (ir->ljpme_combination_rule == eljpmeLB)
3913 if (!bC6ParametersWorkWithLBRules || !bLBRulesPossible)
3915 warning(wi, "You are using arithmetic-geometric combination rules "
3916 "in LJ-PME, but your non-bonded C6 parameters do not "
3917 "follow these rules.");
3920 else
3922 if (!bC6ParametersWorkWithGeometricRules)
3924 if (ir->eDispCorr != edispcNO)
3926 warning_note(wi, "You are using geometric combination rules in "
3927 "LJ-PME, but your non-bonded C6 parameters do "
3928 "not follow these rules. "
3929 "This will introduce very small errors in the forces and energies in "
3930 "your simulations. Dispersion correction will correct total energy "
3931 "and/or pressure for isotropic systems, but not forces or surface tensions.");
3933 else
3935 warning_note(wi, "You are using geometric combination rules in "
3936 "LJ-PME, but your non-bonded C6 parameters do "
3937 "not follow these rules. "
3938 "This will introduce very small errors in the forces and energies in "
3939 "your simulations. If your system is homogeneous, consider using dispersion correction "
3940 "for the total energy and pressure.");
3946 void triple_check(const char *mdparin, t_inputrec *ir, gmx_mtop_t *sys,
3947 warninp_t wi)
3949 char err_buf[STRLEN];
3950 int i, m, c, nmol;
3951 bool bCharge, bAcc;
3952 real *mgrp, mt;
3953 rvec acc;
3954 gmx_mtop_atomloop_block_t aloopb;
3955 ivec AbsRef;
3956 char warn_buf[STRLEN];
3958 set_warning_line(wi, mdparin, -1);
3960 if (ir->cutoff_scheme == ecutsVERLET &&
3961 ir->verletbuf_tol > 0 &&
3962 ir->nstlist > 1 &&
3963 ((EI_MD(ir->eI) || EI_SD(ir->eI)) &&
3964 (ir->etc == etcVRESCALE || ir->etc == etcBERENDSEN)))
3966 /* Check if a too small Verlet buffer might potentially
3967 * cause more drift than the thermostat can couple off.
3969 /* Temperature error fraction for warning and suggestion */
3970 const real T_error_warn = 0.002;
3971 const real T_error_suggest = 0.001;
3972 /* For safety: 2 DOF per atom (typical with constraints) */
3973 const real nrdf_at = 2;
3974 real T, tau, max_T_error;
3975 int i;
3977 T = 0;
3978 tau = 0;
3979 for (i = 0; i < ir->opts.ngtc; i++)
3981 T = std::max(T, ir->opts.ref_t[i]);
3982 tau = std::max(tau, ir->opts.tau_t[i]);
3984 if (T > 0)
3986 /* This is a worst case estimate of the temperature error,
3987 * assuming perfect buffer estimation and no cancelation
3988 * of errors. The factor 0.5 is because energy distributes
3989 * equally over Ekin and Epot.
3991 max_T_error = 0.5*tau*ir->verletbuf_tol/(nrdf_at*BOLTZ*T);
3992 if (max_T_error > T_error_warn)
3994 sprintf(warn_buf, "With a verlet-buffer-tolerance of %g kJ/mol/ps, a reference temperature of %g and a tau_t of %g, your temperature might be off by up to %.1f%%. To ensure the error is below %.1f%%, decrease verlet-buffer-tolerance to %.0e or decrease tau_t.",
3995 ir->verletbuf_tol, T, tau,
3996 100*max_T_error,
3997 100*T_error_suggest,
3998 ir->verletbuf_tol*T_error_suggest/max_T_error);
3999 warning(wi, warn_buf);
4004 if (ETC_ANDERSEN(ir->etc))
4006 int i;
4008 for (i = 0; i < ir->opts.ngtc; i++)
4010 sprintf(err_buf, "all tau_t must currently be equal using Andersen temperature control, violated for group %d", i);
4011 CHECK(ir->opts.tau_t[0] != ir->opts.tau_t[i]);
4012 sprintf(err_buf, "all tau_t must be positive using Andersen temperature control, tau_t[%d]=%10.6f",
4013 i, ir->opts.tau_t[i]);
4014 CHECK(ir->opts.tau_t[i] < 0);
4017 if (ir->etc == etcANDERSENMASSIVE && ir->comm_mode != ecmNO)
4019 for (i = 0; i < ir->opts.ngtc; i++)
4021 int nsteps = gmx::roundToInt(ir->opts.tau_t[i]/ir->delta_t);
4022 sprintf(err_buf, "tau_t/delta_t for group %d for temperature control method %s must be a multiple of nstcomm (%d), as velocities of atoms in coupled groups are randomized every time step. The input tau_t (%8.3f) leads to %d steps per randomization", i, etcoupl_names[ir->etc], ir->nstcomm, ir->opts.tau_t[i], nsteps);
4023 CHECK(nsteps % ir->nstcomm != 0);
4028 if (EI_DYNAMICS(ir->eI) && !EI_SD(ir->eI) && ir->eI != eiBD &&
4029 ir->comm_mode == ecmNO &&
4030 !(absolute_reference(ir, sys, FALSE, AbsRef) || ir->nsteps <= 10) &&
4031 !ETC_ANDERSEN(ir->etc))
4033 warning(wi, "You are not using center of mass motion removal (mdp option comm-mode), numerical rounding errors can lead to build up of kinetic energy of the center of mass");
4036 if (ir->epc == epcPARRINELLORAHMAN &&
4037 ir->etc == etcNOSEHOOVER)
4039 real tau_t_max = 0;
4040 for (int g = 0; g < ir->opts.ngtc; g++)
4042 tau_t_max = std::max(tau_t_max, ir->opts.tau_t[g]);
4044 if (ir->tau_p < 1.9*tau_t_max)
4046 std::string message =
4047 gmx::formatString("With %s T-coupling and %s p-coupling, "
4048 "%s (%g) should be at least twice as large as %s (%g) to avoid resonances",
4049 etcoupl_names[ir->etc],
4050 epcoupl_names[ir->epc],
4051 "tau-p", ir->tau_p,
4052 "tau-t", tau_t_max);
4053 warning(wi, message.c_str());
4057 /* Check for pressure coupling with absolute position restraints */
4058 if (ir->epc != epcNO && ir->refcoord_scaling == erscNO)
4060 absolute_reference(ir, sys, TRUE, AbsRef);
4062 for (m = 0; m < DIM; m++)
4064 if (AbsRef[m] && norm2(ir->compress[m]) > 0)
4066 warning(wi, "You are using pressure coupling with absolute position restraints, this will give artifacts. Use the refcoord_scaling option.");
4067 break;
4073 bCharge = FALSE;
4074 aloopb = gmx_mtop_atomloop_block_init(sys);
4075 const t_atom *atom;
4076 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
4078 if (atom->q != 0 || atom->qB != 0)
4080 bCharge = TRUE;
4084 if (!bCharge)
4086 if (EEL_FULL(ir->coulombtype))
4088 sprintf(err_buf,
4089 "You are using full electrostatics treatment %s for a system without charges.\n"
4090 "This costs a lot of performance for just processing zeros, consider using %s instead.\n",
4091 EELTYPE(ir->coulombtype), EELTYPE(eelCUT));
4092 warning(wi, err_buf);
4095 else
4097 if (ir->coulombtype == eelCUT && ir->rcoulomb > 0)
4099 sprintf(err_buf,
4100 "You are using a plain Coulomb cut-off, which might produce artifacts.\n"
4101 "You might want to consider using %s electrostatics.\n",
4102 EELTYPE(eelPME));
4103 warning_note(wi, err_buf);
4107 /* Check if combination rules used in LJ-PME are the same as in the force field */
4108 if (EVDW_PME(ir->vdwtype))
4110 check_combination_rules(ir, sys, wi);
4113 /* Generalized reaction field */
4114 if (ir->opts.ngtc == 0)
4116 sprintf(err_buf, "No temperature coupling while using coulombtype %s",
4117 eel_names[eelGRF]);
4118 CHECK(ir->coulombtype == eelGRF);
4120 else
4122 sprintf(err_buf, "When using coulombtype = %s"
4123 " ref-t for temperature coupling should be > 0",
4124 eel_names[eelGRF]);
4125 CHECK((ir->coulombtype == eelGRF) && (ir->opts.ref_t[0] <= 0));
4128 bAcc = FALSE;
4129 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4131 for (m = 0; (m < DIM); m++)
4133 if (fabs(ir->opts.acc[i][m]) > 1e-6)
4135 bAcc = TRUE;
4139 if (bAcc)
4141 clear_rvec(acc);
4142 snew(mgrp, sys->groups.grps[egcACC].nr);
4143 for (const AtomProxy &atomP : AtomRange(*sys))
4145 const t_atom &local = atomP.atom();
4146 int i = atomP.globalAtomNumber();
4147 mgrp[getGroupType(sys->groups, egcACC, i)] += local.m;
4149 mt = 0.0;
4150 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4152 for (m = 0; (m < DIM); m++)
4154 acc[m] += ir->opts.acc[i][m]*mgrp[i];
4156 mt += mgrp[i];
4158 for (m = 0; (m < DIM); m++)
4160 if (fabs(acc[m]) > 1e-6)
4162 const char *dim[DIM] = { "X", "Y", "Z" };
4163 fprintf(stderr,
4164 "Net Acceleration in %s direction, will %s be corrected\n",
4165 dim[m], ir->nstcomm != 0 ? "" : "not");
4166 if (ir->nstcomm != 0 && m < ndof_com(ir))
4168 acc[m] /= mt;
4169 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4171 ir->opts.acc[i][m] -= acc[m];
4176 sfree(mgrp);
4179 if (ir->efep != efepNO && ir->fepvals->sc_alpha != 0 &&
4180 !gmx_within_tol(sys->ffparams.reppow, 12.0, 10*GMX_DOUBLE_EPS))
4182 gmx_fatal(FARGS, "Soft-core interactions are only supported with VdW repulsion power 12");
4185 if (ir->bPull)
4187 bool bWarned;
4189 bWarned = FALSE;
4190 for (i = 0; i < ir->pull->ncoord && !bWarned; i++)
4192 if (ir->pull->coord[i].group[0] == 0 ||
4193 ir->pull->coord[i].group[1] == 0)
4195 absolute_reference(ir, sys, FALSE, AbsRef);
4196 for (m = 0; m < DIM; m++)
4198 if (ir->pull->coord[i].dim[m] && !AbsRef[m])
4200 warning(wi, "You are using an absolute reference for pulling, but the rest of the system does not have an absolute reference. This will lead to artifacts.");
4201 bWarned = TRUE;
4202 break;
4208 for (i = 0; i < 3; i++)
4210 for (m = 0; m <= i; m++)
4212 if ((ir->epc != epcNO && ir->compress[i][m] != 0) ||
4213 ir->deform[i][m] != 0)
4215 for (c = 0; c < ir->pull->ncoord; c++)
4217 if (ir->pull->coord[c].eGeom == epullgDIRPBC &&
4218 ir->pull->coord[c].vec[m] != 0)
4220 gmx_fatal(FARGS, "Can not have dynamic box while using pull geometry '%s' (dim %c)", EPULLGEOM(ir->pull->coord[c].eGeom), 'x'+m);
4228 check_disre(sys);
4231 void double_check(t_inputrec *ir, matrix box,
4232 bool bHasNormalConstraints,
4233 bool bHasAnyConstraints,
4234 warninp_t wi)
4236 real min_size;
4237 char warn_buf[STRLEN];
4238 const char *ptr;
4240 ptr = check_box(ir->ePBC, box);
4241 if (ptr)
4243 warning_error(wi, ptr);
4246 if (bHasNormalConstraints && ir->eConstrAlg == econtSHAKE)
4248 if (ir->shake_tol <= 0.0)
4250 sprintf(warn_buf, "ERROR: shake-tol must be > 0 instead of %g\n",
4251 ir->shake_tol);
4252 warning_error(wi, warn_buf);
4256 if ( (ir->eConstrAlg == econtLINCS) && bHasNormalConstraints)
4258 /* If we have Lincs constraints: */
4259 if (ir->eI == eiMD && ir->etc == etcNO &&
4260 ir->eConstrAlg == econtLINCS && ir->nLincsIter == 1)
4262 sprintf(warn_buf, "For energy conservation with LINCS, lincs_iter should be 2 or larger.\n");
4263 warning_note(wi, warn_buf);
4266 if ((ir->eI == eiCG || ir->eI == eiLBFGS) && (ir->nProjOrder < 8))
4268 sprintf(warn_buf, "For accurate %s with LINCS constraints, lincs-order should be 8 or more.", ei_names[ir->eI]);
4269 warning_note(wi, warn_buf);
4271 if (ir->epc == epcMTTK)
4273 warning_error(wi, "MTTK not compatible with lincs -- use shake instead.");
4277 if (bHasAnyConstraints && ir->epc == epcMTTK)
4279 warning_error(wi, "Constraints are not implemented with MTTK pressure control.");
4282 if (ir->LincsWarnAngle > 90.0)
4284 sprintf(warn_buf, "lincs-warnangle can not be larger than 90 degrees, setting it to 90.\n");
4285 warning(wi, warn_buf);
4286 ir->LincsWarnAngle = 90.0;
4289 if (ir->ePBC != epbcNONE)
4291 if (ir->nstlist == 0)
4293 warning(wi, "With nstlist=0 atoms are only put into the box at step 0, therefore drifting atoms might cause the simulation to crash.");
4295 if (ir->ns_type == ensGRID)
4297 if (gmx::square(ir->rlist) >= max_cutoff2(ir->ePBC, box))
4299 sprintf(warn_buf, "ERROR: The cut-off length is longer than half the shortest box vector or longer than the smallest box diagonal element. Increase the box size or decrease rlist.\n");
4300 warning_error(wi, warn_buf);
4303 else
4305 min_size = std::min(box[XX][XX], std::min(box[YY][YY], box[ZZ][ZZ]));
4306 if (2*ir->rlist >= min_size)
4308 sprintf(warn_buf, "ERROR: One of the box lengths is smaller than twice the cut-off length. Increase the box size or decrease rlist.");
4309 warning_error(wi, warn_buf);
4310 if (TRICLINIC(box))
4312 fprintf(stderr, "Grid search might allow larger cut-off's than simple search with triclinic boxes.");
4319 void check_chargegroup_radii(const gmx_mtop_t *mtop, const t_inputrec *ir,
4320 rvec *x,
4321 warninp_t wi)
4323 real rvdw1, rvdw2, rcoul1, rcoul2;
4324 char warn_buf[STRLEN];
4326 calc_chargegroup_radii(mtop, x, &rvdw1, &rvdw2, &rcoul1, &rcoul2);
4328 if (rvdw1 > 0)
4330 printf("Largest charge group radii for Van der Waals: %5.3f, %5.3f nm\n",
4331 rvdw1, rvdw2);
4333 if (rcoul1 > 0)
4335 printf("Largest charge group radii for Coulomb: %5.3f, %5.3f nm\n",
4336 rcoul1, rcoul2);
4339 if (ir->rlist > 0)
4341 if (rvdw1 + rvdw2 > ir->rlist ||
4342 rcoul1 + rcoul2 > ir->rlist)
4344 sprintf(warn_buf,
4345 "The sum of the two largest charge group radii (%f) "
4346 "is larger than rlist (%f)\n",
4347 std::max(rvdw1+rvdw2, rcoul1+rcoul2), ir->rlist);
4348 warning(wi, warn_buf);
4350 else
4352 /* Here we do not use the zero at cut-off macro,
4353 * since user defined interactions might purposely
4354 * not be zero at the cut-off.
4356 if (ir_vdw_is_zero_at_cutoff(ir) &&
4357 rvdw1 + rvdw2 > ir->rlist - ir->rvdw)
4359 sprintf(warn_buf, "The sum of the two largest charge group "
4360 "radii (%f) is larger than rlist (%f) - rvdw (%f).\n"
4361 "With exact cut-offs, better performance can be "
4362 "obtained with cutoff-scheme = %s, because it "
4363 "does not use charge groups at all.",
4364 rvdw1+rvdw2,
4365 ir->rlist, ir->rvdw,
4366 ecutscheme_names[ecutsVERLET]);
4367 if (ir_NVE(ir))
4369 warning(wi, warn_buf);
4371 else
4373 warning_note(wi, warn_buf);
4376 if (ir_coulomb_is_zero_at_cutoff(ir) &&
4377 rcoul1 + rcoul2 > ir->rlist - ir->rcoulomb)
4379 sprintf(warn_buf, "The sum of the two largest charge group radii (%f) is larger than rlist (%f) - rcoulomb (%f).\n"
4380 "With exact cut-offs, better performance can be obtained with cutoff-scheme = %s, because it does not use charge groups at all.",
4381 rcoul1+rcoul2,
4382 ir->rlist, ir->rcoulomb,
4383 ecutscheme_names[ecutsVERLET]);
4384 if (ir_NVE(ir))
4386 warning(wi, warn_buf);
4388 else
4390 warning_note(wi, warn_buf);