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, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
46 #include "gromacs/math/units.h"
49 #include "gromacs/topology/index.h"
50 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/pbcutil/pbc.h"
58 #include "gromacs/topology/mtop_util.h"
59 #include "chargegroup.h"
61 #include "calc_verletbuf.h"
63 #include "gromacs/topology/block.h"
64 #include "gromacs/topology/symtab.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/smalloc.h"
71 /* Resource parameters
72 * Do not change any of these until you read the instruction
73 * in readinp.h. Some cpp's do not take spaces after the backslash
74 * (like the c-shell), which will give you a very weird compiler
78 typedef struct t_inputrec_strings
80 char tcgrps[STRLEN], tau_t[STRLEN], ref_t[STRLEN],
81 acc[STRLEN], accgrps[STRLEN], freeze[STRLEN], frdim[STRLEN],
82 energy[STRLEN], user1[STRLEN], user2[STRLEN], vcm[STRLEN], x_compressed_groups[STRLEN],
83 couple_moltype[STRLEN], orirefitgrp[STRLEN], egptable[STRLEN], egpexcl[STRLEN],
84 wall_atomtype[STRLEN], wall_density[STRLEN], deform[STRLEN], QMMM[STRLEN],
86 char fep_lambda[efptNR][STRLEN];
87 char lambda_weights[STRLEN];
90 char anneal[STRLEN], anneal_npoints[STRLEN],
91 anneal_time[STRLEN], anneal_temp[STRLEN];
92 char QMmethod[STRLEN], QMbasis[STRLEN], QMcharge[STRLEN], QMmult[STRLEN],
93 bSH[STRLEN], CASorbitals[STRLEN], CASelectrons[STRLEN], SAon[STRLEN],
94 SAoff[STRLEN], SAsteps[STRLEN], bTS[STRLEN], bOPT[STRLEN];
95 char efield_x[STRLEN], efield_xt[STRLEN], efield_y[STRLEN],
96 efield_yt[STRLEN], efield_z[STRLEN], efield_zt[STRLEN];
98 } gmx_inputrec_strings;
100 static gmx_inputrec_strings *is = NULL;
102 void init_inputrec_strings()
106 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.");
111 void done_inputrec_strings()
117 static char swapgrp[STRLEN], splitgrp0[STRLEN], splitgrp1[STRLEN], solgrp[STRLEN];
120 egrptpALL, /* All particles have to be a member of a group. */
121 egrptpALL_GENREST, /* A rest group with name is generated for particles *
122 * that are not part of any group. */
123 egrptpPART, /* As egrptpALL_GENREST, but no name is generated *
124 * for the rest group. */
125 egrptpONE /* Merge all selected groups into one group, *
126 * make a rest group for the remaining particles. */
129 static const char *constraints[eshNR+1] = {
130 "none", "h-bonds", "all-bonds", "h-angles", "all-angles", NULL
133 static const char *couple_lam[ecouplamNR+1] = {
134 "vdw-q", "vdw", "q", "none", NULL
137 void init_ir(t_inputrec *ir, t_gromppopts *opts)
139 snew(opts->include, STRLEN);
140 snew(opts->define, STRLEN);
141 snew(ir->fepvals, 1);
142 snew(ir->expandedvals, 1);
143 snew(ir->simtempvals, 1);
146 static void GetSimTemps(int ntemps, t_simtemp *simtemp, double *temperature_lambdas)
151 for (i = 0; i < ntemps; i++)
153 /* simple linear scaling -- allows more control */
154 if (simtemp->eSimTempScale == esimtempLINEAR)
156 simtemp->temperatures[i] = simtemp->simtemp_low + (simtemp->simtemp_high-simtemp->simtemp_low)*temperature_lambdas[i];
158 else if (simtemp->eSimTempScale == esimtempGEOMETRIC) /* should give roughly equal acceptance for constant heat capacity . . . */
160 simtemp->temperatures[i] = simtemp->simtemp_low * pow(simtemp->simtemp_high/simtemp->simtemp_low, (1.0*i)/(ntemps-1));
162 else if (simtemp->eSimTempScale == esimtempEXPONENTIAL)
164 simtemp->temperatures[i] = simtemp->simtemp_low + (simtemp->simtemp_high-simtemp->simtemp_low)*((exp(temperature_lambdas[i])-1)/(exp(1.0)-1));
169 sprintf(errorstr, "eSimTempScale=%d not defined", simtemp->eSimTempScale);
170 gmx_fatal(FARGS, errorstr);
177 static void _low_check(gmx_bool b, char *s, warninp_t wi)
181 warning_error(wi, s);
185 static void check_nst(const char *desc_nst, int nst,
186 const char *desc_p, int *p,
191 if (*p > 0 && *p % nst != 0)
193 /* Round up to the next multiple of nst */
194 *p = ((*p)/nst + 1)*nst;
195 sprintf(buf, "%s should be a multiple of %s, changing %s to %d\n",
196 desc_p, desc_nst, desc_p, *p);
201 static gmx_bool ir_NVE(const t_inputrec *ir)
203 return ((ir->eI == eiMD || EI_VV(ir->eI)) && ir->etc == etcNO);
206 static int lcd(int n1, int n2)
211 for (i = 2; (i <= n1 && i <= n2); i++)
213 if (n1 % i == 0 && n2 % i == 0)
222 static void process_interaction_modifier(const t_inputrec *ir, int *eintmod)
224 if (*eintmod == eintmodPOTSHIFT_VERLET)
226 if (ir->cutoff_scheme == ecutsVERLET)
228 *eintmod = eintmodPOTSHIFT;
232 *eintmod = eintmodNONE;
237 void check_ir(const char *mdparin, t_inputrec *ir, t_gromppopts *opts,
239 /* Check internal consistency.
240 * NOTE: index groups are not set here yet, don't check things
241 * like temperature coupling group options here, but in triple_check
244 /* Strange macro: first one fills the err_buf, and then one can check
245 * the condition, which will print the message and increase the error
248 #define CHECK(b) _low_check(b, err_buf, wi)
249 char err_buf[256], warn_buf[STRLEN];
255 t_lambda *fep = ir->fepvals;
256 t_expanded *expand = ir->expandedvals;
258 set_warning_line(wi, mdparin, -1);
260 /* BASIC CUT-OFF STUFF */
261 if (ir->rcoulomb < 0)
263 warning_error(wi, "rcoulomb should be >= 0");
267 warning_error(wi, "rvdw should be >= 0");
270 !(ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0))
272 warning_error(wi, "rlist should be >= 0");
275 process_interaction_modifier(ir, &ir->coulomb_modifier);
276 process_interaction_modifier(ir, &ir->vdw_modifier);
278 if (ir->cutoff_scheme == ecutsGROUP)
281 "The group cutoff scheme is deprecated in Gromacs 5.0 and will be removed in a future "
282 "release when all interaction forms are supported for the verlet scheme. The verlet "
283 "scheme already scales better, and it is compatible with GPUs and other accelerators.");
285 /* BASIC CUT-OFF STUFF */
286 if (ir->rlist == 0 ||
287 !((ir_coulomb_might_be_zero_at_cutoff(ir) && ir->rcoulomb > ir->rlist) ||
288 (ir_vdw_might_be_zero_at_cutoff(ir) && ir->rvdw > ir->rlist)))
290 /* No switched potential and/or no twin-range:
291 * we can set the long-range cut-off to the maximum of the other cut-offs.
293 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
295 else if (ir->rlistlong < 0)
297 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
298 sprintf(warn_buf, "rlistlong was not set, setting it to %g (no buffer)",
300 warning(wi, warn_buf);
302 if (ir->rlistlong == 0 && ir->ePBC != epbcNONE)
304 warning_error(wi, "Can not have an infinite cut-off with PBC");
306 if (ir->rlistlong > 0 && (ir->rlist == 0 || ir->rlistlong < ir->rlist))
308 warning_error(wi, "rlistlong can not be shorter than rlist");
310 if (IR_TWINRANGE(*ir) && ir->nstlist <= 0)
312 warning_error(wi, "Can not have nstlist<=0 with twin-range interactions");
316 if (ir->rlistlong == ir->rlist)
320 else if (ir->rlistlong > ir->rlist && ir->nstcalclr == 0)
322 warning_error(wi, "With different cutoffs for electrostatics and VdW, nstcalclr must be -1 or a positive number");
325 if (ir->cutoff_scheme == ecutsVERLET)
329 /* Normal Verlet type neighbor-list, currently only limited feature support */
330 if (inputrec2nboundeddim(ir) < 3)
332 warning_error(wi, "With Verlet lists only full pbc or pbc=xy with walls is supported");
334 if (ir->rcoulomb != ir->rvdw)
336 warning_error(wi, "With Verlet lists rcoulomb!=rvdw is not supported");
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;
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 ||
362 (EEL_RF(ir->coulombtype) && ir->coulombtype != eelRF_NEC) ||
363 EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD))
365 warning_error(wi, "With Verlet lists only cut-off, reaction-field, PME and Ewald electrostatics are supported");
367 if (!(ir->coulomb_modifier == eintmodNONE ||
368 ir->coulomb_modifier == eintmodPOTSHIFT))
370 sprintf(warn_buf, "coulomb_modifier=%s is not supported with the Verlet cut-off scheme", eintmod_names[ir->coulomb_modifier]);
371 warning_error(wi, warn_buf);
374 if (ir->nstlist <= 0)
376 warning_error(wi, "With Verlet lists nstlist should be larger than 0");
379 if (ir->nstlist < 10)
381 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.");
384 rc_max = max(ir->rvdw, ir->rcoulomb);
386 if (ir->verletbuf_tol <= 0)
388 if (ir->verletbuf_tol == 0)
390 warning_error(wi, "Can not have Verlet buffer tolerance of exactly 0");
393 if (ir->rlist < rc_max)
395 warning_error(wi, "With verlet lists rlist can not be smaller than rvdw or rcoulomb");
398 if (ir->rlist == rc_max && ir->nstlist > 1)
400 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.");
405 if (ir->rlist > rc_max)
407 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.");
410 if (ir->nstlist == 1)
412 /* No buffer required */
417 if (EI_DYNAMICS(ir->eI))
419 if (inputrec2nboundeddim(ir) < 3)
421 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.");
423 /* Set rlist temporarily so we can continue processing */
428 /* Set the buffer to 5% of the cut-off */
429 ir->rlist = (1.0 + verlet_buffer_ratio_nodynamics)*rc_max;
434 /* No twin-range calculations with Verlet lists */
435 ir->rlistlong = ir->rlist;
438 if (ir->nstcalclr == -1)
440 /* if rlist=rlistlong, this will later be changed to nstcalclr=0 */
441 ir->nstcalclr = ir->nstlist;
443 else if (ir->nstcalclr > 0)
445 if (ir->nstlist > 0 && (ir->nstlist % ir->nstcalclr != 0))
447 warning_error(wi, "nstlist must be evenly divisible by nstcalclr. Use nstcalclr = -1 to automatically follow nstlist");
450 else if (ir->nstcalclr < -1)
452 warning_error(wi, "nstcalclr must be a positive number (divisor of nstcalclr), or -1 to follow nstlist.");
455 if (EEL_PME(ir->coulombtype) && ir->rcoulomb > ir->rvdw && ir->nstcalclr > 1)
457 warning_error(wi, "When used with PME, the long-range component of twin-range interactions must be updated every step (nstcalclr)");
460 /* GENERAL INTEGRATOR STUFF */
461 if (!(ir->eI == eiMD || EI_VV(ir->eI)))
465 if (ir->eI == eiVVAK)
467 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]);
468 warning_note(wi, warn_buf);
470 if (!EI_DYNAMICS(ir->eI))
474 if (EI_DYNAMICS(ir->eI))
476 if (ir->nstcalcenergy < 0)
478 ir->nstcalcenergy = ir_optimal_nstcalcenergy(ir);
479 if (ir->nstenergy != 0 && ir->nstenergy < ir->nstcalcenergy)
481 /* nstcalcenergy larger than nstener does not make sense.
482 * We ideally want nstcalcenergy=nstener.
486 ir->nstcalcenergy = lcd(ir->nstenergy, ir->nstlist);
490 ir->nstcalcenergy = ir->nstenergy;
494 else if ( (ir->nstenergy > 0 && ir->nstcalcenergy > ir->nstenergy) ||
495 (ir->efep != efepNO && ir->fepvals->nstdhdl > 0 &&
496 (ir->nstcalcenergy > ir->fepvals->nstdhdl) ) )
499 const char *nsten = "nstenergy";
500 const char *nstdh = "nstdhdl";
501 const char *min_name = nsten;
502 int min_nst = ir->nstenergy;
504 /* find the smallest of ( nstenergy, nstdhdl ) */
505 if (ir->efep != efepNO && ir->fepvals->nstdhdl > 0 &&
506 (ir->nstenergy == 0 || ir->fepvals->nstdhdl < ir->nstenergy))
508 min_nst = ir->fepvals->nstdhdl;
511 /* If the user sets nstenergy small, we should respect that */
513 "Setting nstcalcenergy (%d) equal to %s (%d)",
514 ir->nstcalcenergy, min_name, min_nst);
515 warning_note(wi, warn_buf);
516 ir->nstcalcenergy = min_nst;
519 if (ir->epc != epcNO)
521 if (ir->nstpcouple < 0)
523 ir->nstpcouple = ir_optimal_nstpcouple(ir);
526 if (IR_TWINRANGE(*ir))
528 check_nst("nstlist", ir->nstlist,
529 "nstcalcenergy", &ir->nstcalcenergy, wi);
530 if (ir->epc != epcNO)
532 check_nst("nstlist", ir->nstlist,
533 "nstpcouple", &ir->nstpcouple, wi);
537 if (ir->nstcalcenergy > 0)
539 if (ir->efep != efepNO)
541 /* nstdhdl should be a multiple of nstcalcenergy */
542 check_nst("nstcalcenergy", ir->nstcalcenergy,
543 "nstdhdl", &ir->fepvals->nstdhdl, wi);
544 /* nstexpanded should be a multiple of nstcalcenergy */
545 check_nst("nstcalcenergy", ir->nstcalcenergy,
546 "nstexpanded", &ir->expandedvals->nstexpanded, wi);
548 /* for storing exact averages nstenergy should be
549 * a multiple of nstcalcenergy
551 check_nst("nstcalcenergy", ir->nstcalcenergy,
552 "nstenergy", &ir->nstenergy, wi);
557 if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
558 ir->bContinuation && ir->ld_seed != -1)
560 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)");
566 sprintf(err_buf, "TPI only works with pbc = %s", epbc_names[epbcXYZ]);
567 CHECK(ir->ePBC != epbcXYZ);
568 sprintf(err_buf, "TPI only works with ns = %s", ens_names[ensGRID]);
569 CHECK(ir->ns_type != ensGRID);
570 sprintf(err_buf, "with TPI nstlist should be larger than zero");
571 CHECK(ir->nstlist <= 0);
572 sprintf(err_buf, "TPI does not work with full electrostatics other than PME");
573 CHECK(EEL_FULL(ir->coulombtype) && !EEL_PME(ir->coulombtype));
577 if ( (opts->nshake > 0) && (opts->bMorse) )
580 "Using morse bond-potentials while constraining bonds is useless");
581 warning(wi, warn_buf);
584 if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
585 ir->bContinuation && ir->ld_seed != -1)
587 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)");
589 /* verify simulated tempering options */
593 gmx_bool bAllTempZero = TRUE;
594 for (i = 0; i < fep->n_lambda; i++)
596 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]);
597 CHECK((fep->all_lambda[efptTEMPERATURE][i] < 0) || (fep->all_lambda[efptTEMPERATURE][i] > 1));
598 if (fep->all_lambda[efptTEMPERATURE][i] > 0)
600 bAllTempZero = FALSE;
603 sprintf(err_buf, "if simulated tempering is on, temperature-lambdas may not be all zero");
604 CHECK(bAllTempZero == TRUE);
606 sprintf(err_buf, "Simulated tempering is currently only compatible with md-vv");
607 CHECK(ir->eI != eiVV);
609 /* check compatability of the temperature coupling with simulated tempering */
611 if (ir->etc == etcNOSEHOOVER)
613 sprintf(warn_buf, "Nose-Hoover based temperature control such as [%s] my not be entirelyconsistent with simulated tempering", etcoupl_names[ir->etc]);
614 warning_note(wi, warn_buf);
617 /* check that the temperatures make sense */
619 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);
620 CHECK(ir->simtempvals->simtemp_high <= ir->simtempvals->simtemp_low);
622 sprintf(err_buf, "Higher simulated tempering temperature (%g) must be >= zero", ir->simtempvals->simtemp_high);
623 CHECK(ir->simtempvals->simtemp_high <= 0);
625 sprintf(err_buf, "Lower simulated tempering temperature (%g) must be >= zero", ir->simtempvals->simtemp_low);
626 CHECK(ir->simtempvals->simtemp_low <= 0);
629 /* verify free energy options */
631 if (ir->efep != efepNO)
634 sprintf(err_buf, "The soft-core power is %d and can only be 1 or 2",
636 CHECK(fep->sc_alpha != 0 && fep->sc_power != 1 && fep->sc_power != 2);
638 sprintf(err_buf, "The soft-core sc-r-power is %d and can only be 6 or 48",
639 (int)fep->sc_r_power);
640 CHECK(fep->sc_alpha != 0 && fep->sc_r_power != 6.0 && fep->sc_r_power != 48.0);
642 sprintf(err_buf, "Can't use postive delta-lambda (%g) if initial state/lambda does not start at zero", fep->delta_lambda);
643 CHECK(fep->delta_lambda > 0 && ((fep->init_fep_state > 0) || (fep->init_lambda > 0)));
645 sprintf(err_buf, "Can't use postive delta-lambda (%g) with expanded ensemble simulations", fep->delta_lambda);
646 CHECK(fep->delta_lambda > 0 && (ir->efep == efepEXPANDED));
648 sprintf(err_buf, "Can only use expanded ensemble with md-vv for now; should be supported for other integrators in 5.0");
649 CHECK(!(EI_VV(ir->eI)) && (ir->efep == efepEXPANDED));
651 sprintf(err_buf, "Free-energy not implemented for Ewald");
652 CHECK(ir->coulombtype == eelEWALD);
654 /* check validty of lambda inputs */
655 if (fep->n_lambda == 0)
657 /* Clear output in case of no states:*/
658 sprintf(err_buf, "init-lambda-state set to %d: no lambda states are defined.", fep->init_fep_state);
659 CHECK((fep->init_fep_state >= 0) && (fep->n_lambda == 0));
663 sprintf(err_buf, "initial thermodynamic state %d does not exist, only goes to %d", fep->init_fep_state, fep->n_lambda-1);
664 CHECK((fep->init_fep_state >= fep->n_lambda));
667 sprintf(err_buf, "Lambda state must be set, either with init-lambda-state or with init-lambda");
668 CHECK((fep->init_fep_state < 0) && (fep->init_lambda < 0));
670 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",
671 fep->init_lambda, fep->init_fep_state);
672 CHECK((fep->init_fep_state >= 0) && (fep->init_lambda >= 0));
676 if ((fep->init_lambda >= 0) && (fep->delta_lambda == 0))
680 for (i = 0; i < efptNR; i++)
682 if (fep->separate_dvdl[i])
687 if (n_lambda_terms > 1)
689 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.");
690 warning(wi, warn_buf);
693 if (n_lambda_terms < 2 && fep->n_lambda > 0)
696 "init-lambda is deprecated for setting lambda state (except for slow growth). Use init-lambda-state instead.");
700 for (j = 0; j < efptNR; j++)
702 for (i = 0; i < fep->n_lambda; i++)
704 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]);
705 CHECK((fep->all_lambda[j][i] < 0) || (fep->all_lambda[j][i] > 1));
709 if ((fep->sc_alpha > 0) && (!fep->bScCoul))
711 for (i = 0; i < fep->n_lambda; i++)
713 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],
714 fep->all_lambda[efptCOUL][i]);
715 CHECK((fep->sc_alpha > 0) &&
716 (((fep->all_lambda[efptCOUL][i] > 0.0) &&
717 (fep->all_lambda[efptCOUL][i] < 1.0)) &&
718 ((fep->all_lambda[efptVDW][i] > 0.0) &&
719 (fep->all_lambda[efptVDW][i] < 1.0))));
723 if ((fep->bScCoul) && (EEL_PME(ir->coulombtype)))
725 real sigma, lambda, r_sc;
728 /* Maximum estimate for A and B charges equal with lambda power 1 */
730 r_sc = pow(lambda*fep->sc_alpha*pow(sigma/ir->rcoulomb, fep->sc_r_power) + 1.0, 1.0/fep->sc_r_power);
731 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.",
733 sigma, lambda, r_sc - 1.0, ir->ewald_rtol);
734 warning_note(wi, warn_buf);
737 /* Free Energy Checks -- In an ideal world, slow growth and FEP would
738 be treated differently, but that's the next step */
740 for (i = 0; i < efptNR; i++)
742 for (j = 0; j < fep->n_lambda; j++)
744 sprintf(err_buf, "%s[%d] must be between 0 and 1", efpt_names[i], j);
745 CHECK((fep->all_lambda[i][j] < 0) || (fep->all_lambda[i][j] > 1));
750 if ((ir->bSimTemp) || (ir->efep == efepEXPANDED))
753 expand = ir->expandedvals;
755 /* checking equilibration of weights inputs for validity */
757 sprintf(err_buf, "weight-equil-number-all-lambda (%d) is ignored if lmc-weights-equil is not equal to %s",
758 expand->equil_n_at_lam, elmceq_names[elmceqNUMATLAM]);
759 CHECK((expand->equil_n_at_lam > 0) && (expand->elmceq != elmceqNUMATLAM));
761 sprintf(err_buf, "weight-equil-number-samples (%d) is ignored if lmc-weights-equil is not equal to %s",
762 expand->equil_samples, elmceq_names[elmceqSAMPLES]);
763 CHECK((expand->equil_samples > 0) && (expand->elmceq != elmceqSAMPLES));
765 sprintf(err_buf, "weight-equil-number-steps (%d) is ignored if lmc-weights-equil is not equal to %s",
766 expand->equil_steps, elmceq_names[elmceqSTEPS]);
767 CHECK((expand->equil_steps > 0) && (expand->elmceq != elmceqSTEPS));
769 sprintf(err_buf, "weight-equil-wl-delta (%d) is ignored if lmc-weights-equil is not equal to %s",
770 expand->equil_samples, elmceq_names[elmceqWLDELTA]);
771 CHECK((expand->equil_wl_delta > 0) && (expand->elmceq != elmceqWLDELTA));
773 sprintf(err_buf, "weight-equil-count-ratio (%f) is ignored if lmc-weights-equil is not equal to %s",
774 expand->equil_ratio, elmceq_names[elmceqRATIO]);
775 CHECK((expand->equil_ratio > 0) && (expand->elmceq != elmceqRATIO));
777 sprintf(err_buf, "weight-equil-number-all-lambda (%d) must be a positive integer if lmc-weights-equil=%s",
778 expand->equil_n_at_lam, elmceq_names[elmceqNUMATLAM]);
779 CHECK((expand->equil_n_at_lam <= 0) && (expand->elmceq == elmceqNUMATLAM));
781 sprintf(err_buf, "weight-equil-number-samples (%d) must be a positive integer if lmc-weights-equil=%s",
782 expand->equil_samples, elmceq_names[elmceqSAMPLES]);
783 CHECK((expand->equil_samples <= 0) && (expand->elmceq == elmceqSAMPLES));
785 sprintf(err_buf, "weight-equil-number-steps (%d) must be a positive integer if lmc-weights-equil=%s",
786 expand->equil_steps, elmceq_names[elmceqSTEPS]);
787 CHECK((expand->equil_steps <= 0) && (expand->elmceq == elmceqSTEPS));
789 sprintf(err_buf, "weight-equil-wl-delta (%f) must be > 0 if lmc-weights-equil=%s",
790 expand->equil_wl_delta, elmceq_names[elmceqWLDELTA]);
791 CHECK((expand->equil_wl_delta <= 0) && (expand->elmceq == elmceqWLDELTA));
793 sprintf(err_buf, "weight-equil-count-ratio (%f) must be > 0 if lmc-weights-equil=%s",
794 expand->equil_ratio, elmceq_names[elmceqRATIO]);
795 CHECK((expand->equil_ratio <= 0) && (expand->elmceq == elmceqRATIO));
797 sprintf(err_buf, "lmc-weights-equil=%s only possible when lmc-stats = %s or lmc-stats %s",
798 elmceq_names[elmceqWLDELTA], elamstats_names[elamstatsWL], elamstats_names[elamstatsWWL]);
799 CHECK((expand->elmceq == elmceqWLDELTA) && (!EWL(expand->elamstats)));
801 sprintf(err_buf, "lmc-repeats (%d) must be greater than 0", expand->lmc_repeats);
802 CHECK((expand->lmc_repeats <= 0));
803 sprintf(err_buf, "minimum-var-min (%d) must be greater than 0", expand->minvarmin);
804 CHECK((expand->minvarmin <= 0));
805 sprintf(err_buf, "weight-c-range (%d) must be greater or equal to 0", expand->c_range);
806 CHECK((expand->c_range < 0));
807 sprintf(err_buf, "init-lambda-state (%d) must be zero if lmc-forced-nstart (%d)> 0 and lmc-move != 'no'",
808 fep->init_fep_state, expand->lmc_forced_nstart);
809 CHECK((fep->init_fep_state != 0) && (expand->lmc_forced_nstart > 0) && (expand->elmcmove != elmcmoveNO));
810 sprintf(err_buf, "lmc-forced-nstart (%d) must not be negative", expand->lmc_forced_nstart);
811 CHECK((expand->lmc_forced_nstart < 0));
812 sprintf(err_buf, "init-lambda-state (%d) must be in the interval [0,number of lambdas)", fep->init_fep_state);
813 CHECK((fep->init_fep_state < 0) || (fep->init_fep_state >= fep->n_lambda));
815 sprintf(err_buf, "init-wl-delta (%f) must be greater than or equal to 0", expand->init_wl_delta);
816 CHECK((expand->init_wl_delta < 0));
817 sprintf(err_buf, "wl-ratio (%f) must be between 0 and 1", expand->wl_ratio);
818 CHECK((expand->wl_ratio <= 0) || (expand->wl_ratio >= 1));
819 sprintf(err_buf, "wl-scale (%f) must be between 0 and 1", expand->wl_scale);
820 CHECK((expand->wl_scale <= 0) || (expand->wl_scale >= 1));
822 /* if there is no temperature control, we need to specify an MC temperature */
823 sprintf(err_buf, "If there is no temperature control, and lmc-mcmove!= 'no',mc_temperature must be set to a positive number");
824 if (expand->nstTij > 0)
826 sprintf(err_buf, "nst-transition-matrix (%d) must be an integer multiple of nstlog (%d)",
827 expand->nstTij, ir->nstlog);
828 CHECK((mod(expand->nstTij, ir->nstlog) != 0));
833 sprintf(err_buf, "walls only work with pbc=%s", epbc_names[epbcXY]);
834 CHECK(ir->nwall && ir->ePBC != epbcXY);
837 if (ir->ePBC != epbcXYZ && ir->nwall != 2)
839 if (ir->ePBC == epbcNONE)
841 if (ir->epc != epcNO)
843 warning(wi, "Turning off pressure coupling for vacuum system");
849 sprintf(err_buf, "Can not have pressure coupling with pbc=%s",
850 epbc_names[ir->ePBC]);
851 CHECK(ir->epc != epcNO);
853 sprintf(err_buf, "Can not have Ewald with pbc=%s", epbc_names[ir->ePBC]);
854 CHECK(EEL_FULL(ir->coulombtype));
856 sprintf(err_buf, "Can not have dispersion correction with pbc=%s",
857 epbc_names[ir->ePBC]);
858 CHECK(ir->eDispCorr != edispcNO);
861 if (ir->rlist == 0.0)
863 sprintf(err_buf, "can only have neighborlist cut-off zero (=infinite)\n"
864 "with coulombtype = %s or coulombtype = %s\n"
865 "without periodic boundary conditions (pbc = %s) and\n"
866 "rcoulomb and rvdw set to zero",
867 eel_names[eelCUT], eel_names[eelUSER], epbc_names[epbcNONE]);
868 CHECK(((ir->coulombtype != eelCUT) && (ir->coulombtype != eelUSER)) ||
869 (ir->ePBC != epbcNONE) ||
870 (ir->rcoulomb != 0.0) || (ir->rvdw != 0.0));
874 warning_error(wi, "Can not have heuristic neighborlist updates without cut-off");
878 warning_note(wi, "Simulating without cut-offs can be (slightly) faster with nstlist=0, nstype=simple and only one MPI rank");
883 if (ir->nstcomm == 0)
885 ir->comm_mode = ecmNO;
887 if (ir->comm_mode != ecmNO)
891 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");
892 ir->nstcomm = abs(ir->nstcomm);
895 if (ir->nstcalcenergy > 0 && ir->nstcomm < ir->nstcalcenergy)
897 warning_note(wi, "nstcomm < nstcalcenergy defeats the purpose of nstcalcenergy, setting nstcomm to nstcalcenergy");
898 ir->nstcomm = ir->nstcalcenergy;
901 if (ir->comm_mode == ecmANGULAR)
903 sprintf(err_buf, "Can not remove the rotation around the center of mass with periodic molecules");
904 CHECK(ir->bPeriodicMols);
905 if (ir->ePBC != epbcNONE)
907 warning(wi, "Removing the rotation around the center of mass in a periodic system (this is not a problem when you have only one molecule).");
912 if (EI_STATE_VELOCITY(ir->eI) && ir->ePBC == epbcNONE && ir->comm_mode != ecmANGULAR)
914 warning_note(wi, "Tumbling and or flying ice-cubes: We are not removing rotation around center of mass in a non-periodic system. You should probably set comm_mode = ANGULAR.");
917 sprintf(err_buf, "Twin-range neighbour searching (NS) with simple NS"
918 " algorithm not implemented");
919 CHECK(((ir->rcoulomb > ir->rlist) || (ir->rvdw > ir->rlist))
920 && (ir->ns_type == ensSIMPLE));
922 /* TEMPERATURE COUPLING */
923 if (ir->etc == etcYES)
925 ir->etc = etcBERENDSEN;
926 warning_note(wi, "Old option for temperature coupling given: "
927 "changing \"yes\" to \"Berendsen\"\n");
930 if ((ir->etc == etcNOSEHOOVER) || (ir->epc == epcMTTK))
932 if (ir->opts.nhchainlength < 1)
934 sprintf(warn_buf, "number of Nose-Hoover chains (currently %d) cannot be less than 1,reset to 1\n", ir->opts.nhchainlength);
935 ir->opts.nhchainlength = 1;
936 warning(wi, warn_buf);
939 if (ir->etc == etcNOSEHOOVER && !EI_VV(ir->eI) && ir->opts.nhchainlength > 1)
941 warning_note(wi, "leapfrog does not yet support Nose-Hoover chains, nhchainlength reset to 1");
942 ir->opts.nhchainlength = 1;
947 ir->opts.nhchainlength = 0;
950 if (ir->eI == eiVVAK)
952 sprintf(err_buf, "%s implemented primarily for validation, and requires nsttcouple = 1 and nstpcouple = 1.",
954 CHECK((ir->nsttcouple != 1) || (ir->nstpcouple != 1));
957 if (ETC_ANDERSEN(ir->etc))
959 sprintf(err_buf, "%s temperature control not supported for integrator %s.", etcoupl_names[ir->etc], ei_names[ir->eI]);
960 CHECK(!(EI_VV(ir->eI)));
962 if (ir->nstcomm > 0 && (ir->etc == etcANDERSEN))
964 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]);
965 warning_note(wi, warn_buf);
968 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]);
969 CHECK(ir->nstcomm > 1 && (ir->etc == etcANDERSEN));
972 if (ir->etc == etcBERENDSEN)
974 sprintf(warn_buf, "The %s thermostat does not generate the correct kinetic energy distribution. You might want to consider using the %s thermostat.",
975 ETCOUPLTYPE(ir->etc), ETCOUPLTYPE(etcVRESCALE));
976 warning_note(wi, warn_buf);
979 if ((ir->etc == etcNOSEHOOVER || ETC_ANDERSEN(ir->etc))
980 && ir->epc == epcBERENDSEN)
982 sprintf(warn_buf, "Using Berendsen pressure coupling invalidates the "
983 "true ensemble for the thermostat");
984 warning(wi, warn_buf);
987 /* PRESSURE COUPLING */
988 if (ir->epc == epcISOTROPIC)
990 ir->epc = epcBERENDSEN;
991 warning_note(wi, "Old option for pressure coupling given: "
992 "changing \"Isotropic\" to \"Berendsen\"\n");
995 if (ir->epc != epcNO)
997 dt_pcoupl = ir->nstpcouple*ir->delta_t;
999 sprintf(err_buf, "tau-p must be > 0 instead of %g\n", ir->tau_p);
1000 CHECK(ir->tau_p <= 0);
1002 if (ir->tau_p/dt_pcoupl < pcouple_min_integration_steps(ir->epc))
1004 sprintf(warn_buf, "For proper integration of the %s barostat, tau-p (%g) should be at least %d times larger than nstpcouple*dt (%g)",
1005 EPCOUPLTYPE(ir->epc), ir->tau_p, pcouple_min_integration_steps(ir->epc), dt_pcoupl);
1006 warning(wi, warn_buf);
1009 sprintf(err_buf, "compressibility must be > 0 when using pressure"
1010 " coupling %s\n", EPCOUPLTYPE(ir->epc));
1011 CHECK(ir->compress[XX][XX] < 0 || ir->compress[YY][YY] < 0 ||
1012 ir->compress[ZZ][ZZ] < 0 ||
1013 (trace(ir->compress) == 0 && ir->compress[YY][XX] <= 0 &&
1014 ir->compress[ZZ][XX] <= 0 && ir->compress[ZZ][YY] <= 0));
1016 if (epcPARRINELLORAHMAN == ir->epc && opts->bGenVel)
1019 "You are generating velocities so I am assuming you "
1020 "are equilibrating a system. You are using "
1021 "%s pressure coupling, but this can be "
1022 "unstable for equilibration. If your system crashes, try "
1023 "equilibrating first with Berendsen pressure coupling. If "
1024 "you are not equilibrating the system, you can probably "
1025 "ignore this warning.",
1026 epcoupl_names[ir->epc]);
1027 warning(wi, warn_buf);
1033 if (ir->epc > epcNO)
1035 if ((ir->epc != epcBERENDSEN) && (ir->epc != epcMTTK))
1037 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.");
1043 if (ir->epc == epcMTTK)
1045 warning_error(wi, "MTTK pressure coupling requires a Velocity-verlet integrator");
1049 /* ELECTROSTATICS */
1050 /* More checks are in triple check (grompp.c) */
1052 if (ir->coulombtype == eelSWITCH)
1054 sprintf(warn_buf, "coulombtype = %s is only for testing purposes and can lead to serious "
1055 "artifacts, advice: use coulombtype = %s",
1056 eel_names[ir->coulombtype],
1057 eel_names[eelRF_ZERO]);
1058 warning(wi, warn_buf);
1061 if (ir->epsilon_r != 1 && ir->implicit_solvent == eisGBSA)
1063 sprintf(warn_buf, "epsilon-r = %g with GB implicit solvent, will use this value for inner dielectric", ir->epsilon_r);
1064 warning_note(wi, warn_buf);
1067 if (EEL_RF(ir->coulombtype) && ir->epsilon_rf == 1 && ir->epsilon_r != 1)
1069 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);
1070 warning(wi, warn_buf);
1071 ir->epsilon_rf = ir->epsilon_r;
1072 ir->epsilon_r = 1.0;
1075 if (getenv("GALACTIC_DYNAMICS") == NULL)
1077 sprintf(err_buf, "epsilon-r must be >= 0 instead of %g\n", ir->epsilon_r);
1078 CHECK(ir->epsilon_r < 0);
1081 if (EEL_RF(ir->coulombtype))
1083 /* reaction field (at the cut-off) */
1085 if (ir->coulombtype == eelRF_ZERO)
1087 sprintf(warn_buf, "With coulombtype = %s, epsilon-rf must be 0, assuming you meant epsilon_rf=0",
1088 eel_names[ir->coulombtype]);
1089 CHECK(ir->epsilon_rf != 0);
1090 ir->epsilon_rf = 0.0;
1093 sprintf(err_buf, "epsilon-rf must be >= epsilon-r");
1094 CHECK((ir->epsilon_rf < ir->epsilon_r && ir->epsilon_rf != 0) ||
1095 (ir->epsilon_r == 0));
1096 if (ir->epsilon_rf == ir->epsilon_r)
1098 sprintf(warn_buf, "Using epsilon-rf = epsilon-r with %s does not make sense",
1099 eel_names[ir->coulombtype]);
1100 warning(wi, warn_buf);
1103 /* Allow rlist>rcoulomb for tabulated long range stuff. This just
1104 * means the interaction is zero outside rcoulomb, but it helps to
1105 * provide accurate energy conservation.
1107 if (ir_coulomb_might_be_zero_at_cutoff(ir))
1109 if (ir_coulomb_switched(ir))
1112 "With coulombtype = %s rcoulomb_switch must be < rcoulomb. Or, better: Use the potential modifier options!",
1113 eel_names[ir->coulombtype]);
1114 CHECK(ir->rcoulomb_switch >= ir->rcoulomb);
1117 else if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
1119 if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE)
1121 sprintf(err_buf, "With coulombtype = %s, rcoulomb should be >= rlist unless you use a potential modifier",
1122 eel_names[ir->coulombtype]);
1123 CHECK(ir->rlist > ir->rcoulomb);
1127 if (ir->coulombtype == eelSWITCH || ir->coulombtype == eelSHIFT ||
1128 ir->vdwtype == evdwSWITCH || ir->vdwtype == evdwSHIFT)
1131 "The switch/shift interaction settings are just for compatibility; you will get better "
1132 "performance from applying potential modifiers to your interactions!\n");
1133 warning_note(wi, warn_buf);
1136 if (ir->coulombtype == eelPMESWITCH || ir->coulomb_modifier == eintmodPOTSWITCH)
1138 if (ir->rcoulomb_switch/ir->rcoulomb < 0.9499)
1140 real percentage = 100*(ir->rcoulomb-ir->rcoulomb_switch)/ir->rcoulomb;
1141 sprintf(warn_buf, "The switching range for 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.",
1142 percentage, ir->rcoulomb_switch, ir->rcoulomb, ir->ewald_rtol);
1143 warning(wi, warn_buf);
1147 if (ir->vdwtype == evdwSWITCH || ir->vdw_modifier == eintmodPOTSWITCH)
1149 if (ir->rvdw_switch == 0)
1151 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.");
1152 warning(wi, warn_buf);
1156 if (EEL_FULL(ir->coulombtype))
1158 if (ir->coulombtype == eelPMESWITCH || ir->coulombtype == eelPMEUSER ||
1159 ir->coulombtype == eelPMEUSERSWITCH)
1161 sprintf(err_buf, "With coulombtype = %s, rcoulomb must be <= rlist",
1162 eel_names[ir->coulombtype]);
1163 CHECK(ir->rcoulomb > ir->rlist);
1165 else if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE)
1167 if (ir->coulombtype == eelPME || ir->coulombtype == eelP3M_AD)
1170 "With coulombtype = %s (without modifier), rcoulomb must be equal to rlist,\n"
1171 "or rlistlong if nstcalclr=1. For optimal energy conservation,consider using\n"
1172 "a potential modifier.", eel_names[ir->coulombtype]);
1173 if (ir->nstcalclr == 1)
1175 CHECK(ir->rcoulomb != ir->rlist && ir->rcoulomb != ir->rlistlong);
1179 CHECK(ir->rcoulomb != ir->rlist);
1185 if (EEL_PME(ir->coulombtype) || EVDW_PME(ir->vdwtype))
1187 if (ir->pme_order < 3)
1189 warning_error(wi, "pme-order can not be smaller than 3");
1193 if (ir->nwall == 2 && EEL_FULL(ir->coulombtype))
1195 if (ir->ewald_geometry == eewg3D)
1197 sprintf(warn_buf, "With pbc=%s you should use ewald-geometry=%s",
1198 epbc_names[ir->ePBC], eewg_names[eewg3DC]);
1199 warning(wi, warn_buf);
1201 /* This check avoids extra pbc coding for exclusion corrections */
1202 sprintf(err_buf, "wall-ewald-zfac should be >= 2");
1203 CHECK(ir->wall_ewald_zfac < 2);
1206 if (ir_vdw_switched(ir))
1208 sprintf(err_buf, "With switched vdw forces or potentials, rvdw-switch must be < rvdw");
1209 CHECK(ir->rvdw_switch >= ir->rvdw);
1211 if (ir->rvdw_switch < 0.5*ir->rvdw)
1213 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.",
1214 ir->rvdw_switch, ir->rvdw);
1215 warning_note(wi, warn_buf);
1218 else if (ir->vdwtype == evdwCUT || ir->vdwtype == evdwPME)
1220 if (ir->cutoff_scheme == ecutsGROUP && ir->vdw_modifier == eintmodNONE)
1222 sprintf(err_buf, "With vdwtype = %s, rvdw must be >= rlist unless you use a potential modifier", evdw_names[ir->vdwtype]);
1223 CHECK(ir->rlist > ir->rvdw);
1227 if (ir->vdwtype == evdwPME)
1229 if (!(ir->vdw_modifier == eintmodNONE || ir->vdw_modifier == eintmodPOTSHIFT))
1231 sprintf(err_buf, "With vdwtype = %s, the only supported modifiers are %s a\
1233 evdw_names[ir->vdwtype],
1234 eintmod_names[eintmodPOTSHIFT],
1235 eintmod_names[eintmodNONE]);
1239 if (ir->cutoff_scheme == ecutsGROUP)
1241 if (((ir->coulomb_modifier != eintmodNONE && ir->rcoulomb == ir->rlist) ||
1242 (ir->vdw_modifier != eintmodNONE && ir->rvdw == ir->rlist)) &&
1245 warning_note(wi, "With exact cut-offs, rlist should be "
1246 "larger than rcoulomb and rvdw, so that there "
1247 "is a buffer region for particle motion "
1248 "between neighborsearch steps");
1251 if (ir_coulomb_is_zero_at_cutoff(ir) && ir->rlistlong <= ir->rcoulomb)
1253 sprintf(warn_buf, "For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rcoulomb.",
1254 IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
1255 warning_note(wi, warn_buf);
1257 if (ir_vdw_switched(ir) && (ir->rlistlong <= ir->rvdw))
1259 sprintf(warn_buf, "For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rvdw.",
1260 IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
1261 warning_note(wi, warn_buf);
1265 if (ir->vdwtype == evdwUSER && ir->eDispCorr != edispcNO)
1267 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.");
1270 if (ir->nstlist == -1)
1272 sprintf(err_buf, "With nstlist=-1 rvdw and rcoulomb should be smaller than rlist to account for diffusion and possibly charge-group radii");
1273 CHECK(ir->rvdw >= ir->rlist || ir->rcoulomb >= ir->rlist);
1275 sprintf(err_buf, "nstlist can not be smaller than -1");
1276 CHECK(ir->nstlist < -1);
1278 if (ir->eI == eiLBFGS && (ir->coulombtype == eelCUT || ir->vdwtype == evdwCUT)
1281 warning(wi, "For efficient BFGS minimization, use switch/shift/pme instead of cut-off.");
1284 if (ir->eI == eiLBFGS && ir->nbfgscorr <= 0)
1286 warning(wi, "Using L-BFGS with nbfgscorr<=0 just gets you steepest descent.");
1289 /* ENERGY CONSERVATION */
1290 if (ir_NVE(ir) && ir->cutoff_scheme == ecutsGROUP)
1292 if (!ir_vdw_might_be_zero_at_cutoff(ir) && ir->rvdw > 0 && ir->vdw_modifier == eintmodNONE)
1294 sprintf(warn_buf, "You are using a cut-off for VdW interactions with NVE, for good energy conservation use vdwtype = %s (possibly with DispCorr)",
1295 evdw_names[evdwSHIFT]);
1296 warning_note(wi, warn_buf);
1298 if (!ir_coulomb_might_be_zero_at_cutoff(ir) && ir->rcoulomb > 0)
1300 sprintf(warn_buf, "You are using a cut-off for electrostatics with NVE, for good energy conservation use coulombtype = %s or %s",
1301 eel_names[eelPMESWITCH], eel_names[eelRF_ZERO]);
1302 warning_note(wi, warn_buf);
1306 if (EI_VV(ir->eI) && IR_TWINRANGE(*ir) && ir->nstlist > 1)
1308 sprintf(warn_buf, "Twin-range multiple time stepping does not work with integrator %s.", ei_names[ir->eI]);
1309 warning_error(wi, warn_buf);
1312 /* IMPLICIT SOLVENT */
1313 if (ir->coulombtype == eelGB_NOTUSED)
1315 ir->coulombtype = eelCUT;
1316 ir->implicit_solvent = eisGBSA;
1317 fprintf(stderr, "Note: Old option for generalized born electrostatics given:\n"
1318 "Changing coulombtype from \"generalized-born\" to \"cut-off\" and instead\n"
1319 "setting implicit-solvent value to \"GBSA\" in input section.\n");
1322 if (ir->sa_algorithm == esaSTILL)
1324 sprintf(err_buf, "Still SA algorithm not available yet, use %s or %s instead\n", esa_names[esaAPPROX], esa_names[esaNO]);
1325 CHECK(ir->sa_algorithm == esaSTILL);
1328 if (ir->implicit_solvent == eisGBSA)
1330 sprintf(err_buf, "With GBSA implicit solvent, rgbradii must be equal to rlist.");
1331 CHECK(ir->rgbradii != ir->rlist);
1333 if (ir->coulombtype != eelCUT)
1335 sprintf(err_buf, "With GBSA, coulombtype must be equal to %s\n", eel_names[eelCUT]);
1336 CHECK(ir->coulombtype != eelCUT);
1338 if (ir->vdwtype != evdwCUT)
1340 sprintf(err_buf, "With GBSA, vdw-type must be equal to %s\n", evdw_names[evdwCUT]);
1341 CHECK(ir->vdwtype != evdwCUT);
1343 if (ir->nstgbradii < 1)
1345 sprintf(warn_buf, "Using GBSA with nstgbradii<1, setting nstgbradii=1");
1346 warning_note(wi, warn_buf);
1349 if (ir->sa_algorithm == esaNO)
1351 sprintf(warn_buf, "No SA (non-polar) calculation requested together with GB. Are you sure this is what you want?\n");
1352 warning_note(wi, warn_buf);
1354 if (ir->sa_surface_tension < 0 && ir->sa_algorithm != esaNO)
1356 sprintf(warn_buf, "Value of sa_surface_tension is < 0. Changing it to 2.05016 or 2.25936 kJ/nm^2/mol for Still and HCT/OBC respectively\n");
1357 warning_note(wi, warn_buf);
1359 if (ir->gb_algorithm == egbSTILL)
1361 ir->sa_surface_tension = 0.0049 * CAL2JOULE * 100;
1365 ir->sa_surface_tension = 0.0054 * CAL2JOULE * 100;
1368 if (ir->sa_surface_tension == 0 && ir->sa_algorithm != esaNO)
1370 sprintf(err_buf, "Surface tension set to 0 while SA-calculation requested\n");
1371 CHECK(ir->sa_surface_tension == 0 && ir->sa_algorithm != esaNO);
1378 if (ir->cutoff_scheme != ecutsGROUP)
1380 warning_error(wi, "AdresS simulation supports only cutoff-scheme=group");
1384 warning_error(wi, "AdresS simulation supports only stochastic dynamics");
1386 if (ir->epc != epcNO)
1388 warning_error(wi, "AdresS simulation does not support pressure coupling");
1390 if (EEL_FULL(ir->coulombtype))
1392 warning_error(wi, "AdresS simulation does not support long-range electrostatics");
1397 /* count the number of text elemets separated by whitespace in a string.
1398 str = the input string
1399 maxptr = the maximum number of allowed elements
1400 ptr = the output array of pointers to the first character of each element
1401 returns: the number of elements. */
1402 int str_nelem(const char *str, int maxptr, char *ptr[])
1407 copy0 = strdup(str);
1410 while (*copy != '\0')
1414 gmx_fatal(FARGS, "Too many groups on line: '%s' (max is %d)",
1422 while ((*copy != '\0') && !isspace(*copy))
1441 /* interpret a number of doubles from a string and put them in an array,
1442 after allocating space for them.
1443 str = the input string
1444 n = the (pre-allocated) number of doubles read
1445 r = the output array of doubles. */
1446 static void parse_n_real(char *str, int *n, real **r)
1451 *n = str_nelem(str, MAXPTR, ptr);
1454 for (i = 0; i < *n; i++)
1456 (*r)[i] = strtod(ptr[i], NULL);
1460 static void do_fep_params(t_inputrec *ir, char fep_lambda[][STRLEN], char weights[STRLEN])
1463 int i, j, max_n_lambda, nweights, nfep[efptNR];
1464 t_lambda *fep = ir->fepvals;
1465 t_expanded *expand = ir->expandedvals;
1466 real **count_fep_lambdas;
1467 gmx_bool bOneLambda = TRUE;
1469 snew(count_fep_lambdas, efptNR);
1471 /* FEP input processing */
1472 /* first, identify the number of lambda values for each type.
1473 All that are nonzero must have the same number */
1475 for (i = 0; i < efptNR; i++)
1477 parse_n_real(fep_lambda[i], &(nfep[i]), &(count_fep_lambdas[i]));
1480 /* now, determine the number of components. All must be either zero, or equal. */
1483 for (i = 0; i < efptNR; i++)
1485 if (nfep[i] > max_n_lambda)
1487 max_n_lambda = nfep[i]; /* here's a nonzero one. All of them
1488 must have the same number if its not zero.*/
1493 for (i = 0; i < efptNR; i++)
1497 ir->fepvals->separate_dvdl[i] = FALSE;
1499 else if (nfep[i] == max_n_lambda)
1501 if (i != efptTEMPERATURE) /* we treat this differently -- not really a reason to compute the derivative with
1502 respect to the temperature currently */
1504 ir->fepvals->separate_dvdl[i] = TRUE;
1509 gmx_fatal(FARGS, "Number of lambdas (%d) for FEP type %s not equal to number of other types (%d)",
1510 nfep[i], efpt_names[i], max_n_lambda);
1513 /* we don't print out dhdl if the temperature is changing, since we can't correctly define dhdl in this case */
1514 ir->fepvals->separate_dvdl[efptTEMPERATURE] = FALSE;
1516 /* the number of lambdas is the number we've read in, which is either zero
1517 or the same for all */
1518 fep->n_lambda = max_n_lambda;
1520 /* allocate space for the array of lambda values */
1521 snew(fep->all_lambda, efptNR);
1522 /* if init_lambda is defined, we need to set lambda */
1523 if ((fep->init_lambda > 0) && (fep->n_lambda == 0))
1525 ir->fepvals->separate_dvdl[efptFEP] = TRUE;
1527 /* otherwise allocate the space for all of the lambdas, and transfer the data */
1528 for (i = 0; i < efptNR; i++)
1530 snew(fep->all_lambda[i], fep->n_lambda);
1531 if (nfep[i] > 0) /* if it's zero, then the count_fep_lambda arrays
1534 for (j = 0; j < fep->n_lambda; j++)
1536 fep->all_lambda[i][j] = (double)count_fep_lambdas[i][j];
1538 sfree(count_fep_lambdas[i]);
1541 sfree(count_fep_lambdas);
1543 /* "fep-vals" is either zero or the full number. If zero, we'll need to define fep-lambdas for internal
1544 bookkeeping -- for now, init_lambda */
1546 if ((nfep[efptFEP] == 0) && (fep->init_lambda >= 0))
1548 for (i = 0; i < fep->n_lambda; i++)
1550 fep->all_lambda[efptFEP][i] = fep->init_lambda;
1554 /* check to see if only a single component lambda is defined, and soft core is defined.
1555 In this case, turn on coulomb soft core */
1557 if (max_n_lambda == 0)
1563 for (i = 0; i < efptNR; i++)
1565 if ((nfep[i] != 0) && (i != efptFEP))
1571 if ((bOneLambda) && (fep->sc_alpha > 0))
1573 fep->bScCoul = TRUE;
1576 /* Fill in the others with the efptFEP if they are not explicitly
1577 specified (i.e. nfep[i] == 0). This means if fep is not defined,
1578 they are all zero. */
1580 for (i = 0; i < efptNR; i++)
1582 if ((nfep[i] == 0) && (i != efptFEP))
1584 for (j = 0; j < fep->n_lambda; j++)
1586 fep->all_lambda[i][j] = fep->all_lambda[efptFEP][j];
1592 /* make it easier if sc_r_power = 48 by increasing it to the 4th power, to be in the right scale. */
1593 if (fep->sc_r_power == 48)
1595 if (fep->sc_alpha > 0.1)
1597 gmx_fatal(FARGS, "sc_alpha (%f) for sc_r_power = 48 should usually be between 0.001 and 0.004", fep->sc_alpha);
1601 expand = ir->expandedvals;
1602 /* now read in the weights */
1603 parse_n_real(weights, &nweights, &(expand->init_lambda_weights));
1606 snew(expand->init_lambda_weights, fep->n_lambda); /* initialize to zero */
1608 else if (nweights != fep->n_lambda)
1610 gmx_fatal(FARGS, "Number of weights (%d) is not equal to number of lambda values (%d)",
1611 nweights, fep->n_lambda);
1613 if ((expand->nstexpanded < 0) && (ir->efep != efepNO))
1615 expand->nstexpanded = fep->nstdhdl;
1616 /* if you don't specify nstexpanded when doing expanded ensemble free energy calcs, it is set to nstdhdl */
1618 if ((expand->nstexpanded < 0) && ir->bSimTemp)
1620 expand->nstexpanded = 2*(int)(ir->opts.tau_t[0]/ir->delta_t);
1621 /* if you don't specify nstexpanded when doing expanded ensemble simulated tempering, it is set to
1622 2*tau_t just to be careful so it's not to frequent */
1627 static void do_simtemp_params(t_inputrec *ir)
1630 snew(ir->simtempvals->temperatures, ir->fepvals->n_lambda);
1631 GetSimTemps(ir->fepvals->n_lambda, ir->simtempvals, ir->fepvals->all_lambda[efptTEMPERATURE]);
1636 static void do_wall_params(t_inputrec *ir,
1637 char *wall_atomtype, char *wall_density,
1641 char *names[MAXPTR];
1644 opts->wall_atomtype[0] = NULL;
1645 opts->wall_atomtype[1] = NULL;
1647 ir->wall_atomtype[0] = -1;
1648 ir->wall_atomtype[1] = -1;
1649 ir->wall_density[0] = 0;
1650 ir->wall_density[1] = 0;
1654 nstr = str_nelem(wall_atomtype, MAXPTR, names);
1655 if (nstr != ir->nwall)
1657 gmx_fatal(FARGS, "Expected %d elements for wall_atomtype, found %d",
1660 for (i = 0; i < ir->nwall; i++)
1662 opts->wall_atomtype[i] = strdup(names[i]);
1665 if (ir->wall_type == ewt93 || ir->wall_type == ewt104)
1667 nstr = str_nelem(wall_density, MAXPTR, names);
1668 if (nstr != ir->nwall)
1670 gmx_fatal(FARGS, "Expected %d elements for wall-density, found %d", ir->nwall, nstr);
1672 for (i = 0; i < ir->nwall; i++)
1674 sscanf(names[i], "%lf", &dbl);
1677 gmx_fatal(FARGS, "wall-density[%d] = %f\n", i, dbl);
1679 ir->wall_density[i] = dbl;
1685 static void add_wall_energrps(gmx_groups_t *groups, int nwall, t_symtab *symtab)
1693 srenew(groups->grpname, groups->ngrpname+nwall);
1694 grps = &(groups->grps[egcENER]);
1695 srenew(grps->nm_ind, grps->nr+nwall);
1696 for (i = 0; i < nwall; i++)
1698 sprintf(str, "wall%d", i);
1699 groups->grpname[groups->ngrpname] = put_symtab(symtab, str);
1700 grps->nm_ind[grps->nr++] = groups->ngrpname++;
1705 void read_expandedparams(int *ninp_p, t_inpfile **inp_p,
1706 t_expanded *expand, warninp_t wi)
1708 int ninp, nerror = 0;
1714 /* read expanded ensemble parameters */
1715 CCTYPE ("expanded ensemble variables");
1716 ITYPE ("nstexpanded", expand->nstexpanded, -1);
1717 EETYPE("lmc-stats", expand->elamstats, elamstats_names);
1718 EETYPE("lmc-move", expand->elmcmove, elmcmove_names);
1719 EETYPE("lmc-weights-equil", expand->elmceq, elmceq_names);
1720 ITYPE ("weight-equil-number-all-lambda", expand->equil_n_at_lam, -1);
1721 ITYPE ("weight-equil-number-samples", expand->equil_samples, -1);
1722 ITYPE ("weight-equil-number-steps", expand->equil_steps, -1);
1723 RTYPE ("weight-equil-wl-delta", expand->equil_wl_delta, -1);
1724 RTYPE ("weight-equil-count-ratio", expand->equil_ratio, -1);
1725 CCTYPE("Seed for Monte Carlo in lambda space");
1726 ITYPE ("lmc-seed", expand->lmc_seed, -1);
1727 RTYPE ("mc-temperature", expand->mc_temp, -1);
1728 ITYPE ("lmc-repeats", expand->lmc_repeats, 1);
1729 ITYPE ("lmc-gibbsdelta", expand->gibbsdeltalam, -1);
1730 ITYPE ("lmc-forced-nstart", expand->lmc_forced_nstart, 0);
1731 EETYPE("symmetrized-transition-matrix", expand->bSymmetrizedTMatrix, yesno_names);
1732 ITYPE("nst-transition-matrix", expand->nstTij, -1);
1733 ITYPE ("mininum-var-min", expand->minvarmin, 100); /*default is reasonable */
1734 ITYPE ("weight-c-range", expand->c_range, 0); /* default is just C=0 */
1735 RTYPE ("wl-scale", expand->wl_scale, 0.8);
1736 RTYPE ("wl-ratio", expand->wl_ratio, 0.8);
1737 RTYPE ("init-wl-delta", expand->init_wl_delta, 1.0);
1738 EETYPE("wl-oneovert", expand->bWLoneovert, yesno_names);
1746 void get_ir(const char *mdparin, const char *mdparout,
1747 t_inputrec *ir, t_gromppopts *opts,
1751 double dumdub[2][6];
1755 char warn_buf[STRLEN];
1756 t_lambda *fep = ir->fepvals;
1757 t_expanded *expand = ir->expandedvals;
1759 init_inputrec_strings();
1760 inp = read_inpfile(mdparin, &ninp, wi);
1762 snew(dumstr[0], STRLEN);
1763 snew(dumstr[1], STRLEN);
1765 if (-1 == search_einp(ninp, inp, "cutoff-scheme"))
1768 "%s did not specify a value for the .mdp option "
1769 "\"cutoff-scheme\". Probably it was first intended for use "
1770 "with GROMACS before 4.6. In 4.6, the Verlet scheme was "
1771 "introduced, but the group scheme was still the default. "
1772 "The default is now the Verlet scheme, so you will observe "
1773 "different behaviour.", mdparin);
1774 warning_note(wi, warn_buf);
1777 /* remove the following deprecated commands */
1780 REM_TYPE("domain-decomposition");
1781 REM_TYPE("andersen-seed");
1783 REM_TYPE("dihre-fc");
1784 REM_TYPE("dihre-tau");
1785 REM_TYPE("nstdihreout");
1786 REM_TYPE("nstcheckpoint");
1788 /* replace the following commands with the clearer new versions*/
1789 REPL_TYPE("unconstrained-start", "continuation");
1790 REPL_TYPE("foreign-lambda", "fep-lambdas");
1791 REPL_TYPE("verlet-buffer-drift", "verlet-buffer-tolerance");
1792 REPL_TYPE("nstxtcout", "nstxout-compressed");
1793 REPL_TYPE("xtc-grps", "compressed-x-grps");
1794 REPL_TYPE("xtc-precision", "compressed-x-precision");
1796 CCTYPE ("VARIOUS PREPROCESSING OPTIONS");
1797 CTYPE ("Preprocessor information: use cpp syntax.");
1798 CTYPE ("e.g.: -I/home/joe/doe -I/home/mary/roe");
1799 STYPE ("include", opts->include, NULL);
1800 CTYPE ("e.g.: -DPOSRES -DFLEXIBLE (note these variable names are case sensitive)");
1801 STYPE ("define", opts->define, NULL);
1803 CCTYPE ("RUN CONTROL PARAMETERS");
1804 EETYPE("integrator", ir->eI, ei_names);
1805 CTYPE ("Start time and timestep in ps");
1806 RTYPE ("tinit", ir->init_t, 0.0);
1807 RTYPE ("dt", ir->delta_t, 0.001);
1808 STEPTYPE ("nsteps", ir->nsteps, 0);
1809 CTYPE ("For exact run continuation or redoing part of a run");
1810 STEPTYPE ("init-step", ir->init_step, 0);
1811 CTYPE ("Part index is updated automatically on checkpointing (keeps files separate)");
1812 ITYPE ("simulation-part", ir->simulation_part, 1);
1813 CTYPE ("mode for center of mass motion removal");
1814 EETYPE("comm-mode", ir->comm_mode, ecm_names);
1815 CTYPE ("number of steps for center of mass motion removal");
1816 ITYPE ("nstcomm", ir->nstcomm, 100);
1817 CTYPE ("group(s) for center of mass motion removal");
1818 STYPE ("comm-grps", is->vcm, NULL);
1820 CCTYPE ("LANGEVIN DYNAMICS OPTIONS");
1821 CTYPE ("Friction coefficient (amu/ps) and random seed");
1822 RTYPE ("bd-fric", ir->bd_fric, 0.0);
1823 STEPTYPE ("ld-seed", ir->ld_seed, -1);
1826 CCTYPE ("ENERGY MINIMIZATION OPTIONS");
1827 CTYPE ("Force tolerance and initial step-size");
1828 RTYPE ("emtol", ir->em_tol, 10.0);
1829 RTYPE ("emstep", ir->em_stepsize, 0.01);
1830 CTYPE ("Max number of iterations in relax-shells");
1831 ITYPE ("niter", ir->niter, 20);
1832 CTYPE ("Step size (ps^2) for minimization of flexible constraints");
1833 RTYPE ("fcstep", ir->fc_stepsize, 0);
1834 CTYPE ("Frequency of steepest descents steps when doing CG");
1835 ITYPE ("nstcgsteep", ir->nstcgsteep, 1000);
1836 ITYPE ("nbfgscorr", ir->nbfgscorr, 10);
1838 CCTYPE ("TEST PARTICLE INSERTION OPTIONS");
1839 RTYPE ("rtpi", ir->rtpi, 0.05);
1841 /* Output options */
1842 CCTYPE ("OUTPUT CONTROL OPTIONS");
1843 CTYPE ("Output frequency for coords (x), velocities (v) and forces (f)");
1844 ITYPE ("nstxout", ir->nstxout, 0);
1845 ITYPE ("nstvout", ir->nstvout, 0);
1846 ITYPE ("nstfout", ir->nstfout, 0);
1847 ir->nstcheckpoint = 1000;
1848 CTYPE ("Output frequency for energies to log file and energy file");
1849 ITYPE ("nstlog", ir->nstlog, 1000);
1850 ITYPE ("nstcalcenergy", ir->nstcalcenergy, 100);
1851 ITYPE ("nstenergy", ir->nstenergy, 1000);
1852 CTYPE ("Output frequency and precision for .xtc file");
1853 ITYPE ("nstxout-compressed", ir->nstxout_compressed, 0);
1854 RTYPE ("compressed-x-precision", ir->x_compression_precision, 1000.0);
1855 CTYPE ("This selects the subset of atoms for the compressed");
1856 CTYPE ("trajectory file. You can select multiple groups. By");
1857 CTYPE ("default, all atoms will be written.");
1858 STYPE ("compressed-x-grps", is->x_compressed_groups, NULL);
1859 CTYPE ("Selection of energy groups");
1860 STYPE ("energygrps", is->energy, NULL);
1862 /* Neighbor searching */
1863 CCTYPE ("NEIGHBORSEARCHING PARAMETERS");
1864 CTYPE ("cut-off scheme (Verlet: particle based cut-offs, group: using charge groups)");
1865 EETYPE("cutoff-scheme", ir->cutoff_scheme, ecutscheme_names);
1866 CTYPE ("nblist update frequency");
1867 ITYPE ("nstlist", ir->nstlist, 10);
1868 CTYPE ("ns algorithm (simple or grid)");
1869 EETYPE("ns-type", ir->ns_type, ens_names);
1870 /* set ndelta to the optimal value of 2 */
1872 CTYPE ("Periodic boundary conditions: xyz, no, xy");
1873 EETYPE("pbc", ir->ePBC, epbc_names);
1874 EETYPE("periodic-molecules", ir->bPeriodicMols, yesno_names);
1875 CTYPE ("Allowed energy error due to the Verlet buffer in kJ/mol/ps per atom,");
1876 CTYPE ("a value of -1 means: use rlist");
1877 RTYPE("verlet-buffer-tolerance", ir->verletbuf_tol, 0.005);
1878 CTYPE ("nblist cut-off");
1879 RTYPE ("rlist", ir->rlist, 1.0);
1880 CTYPE ("long-range cut-off for switched potentials");
1881 RTYPE ("rlistlong", ir->rlistlong, -1);
1882 ITYPE ("nstcalclr", ir->nstcalclr, -1);
1884 /* Electrostatics */
1885 CCTYPE ("OPTIONS FOR ELECTROSTATICS AND VDW");
1886 CTYPE ("Method for doing electrostatics");
1887 EETYPE("coulombtype", ir->coulombtype, eel_names);
1888 EETYPE("coulomb-modifier", ir->coulomb_modifier, eintmod_names);
1889 CTYPE ("cut-off lengths");
1890 RTYPE ("rcoulomb-switch", ir->rcoulomb_switch, 0.0);
1891 RTYPE ("rcoulomb", ir->rcoulomb, 1.0);
1892 CTYPE ("Relative dielectric constant for the medium and the reaction field");
1893 RTYPE ("epsilon-r", ir->epsilon_r, 1.0);
1894 RTYPE ("epsilon-rf", ir->epsilon_rf, 0.0);
1895 CTYPE ("Method for doing Van der Waals");
1896 EETYPE("vdw-type", ir->vdwtype, evdw_names);
1897 EETYPE("vdw-modifier", ir->vdw_modifier, eintmod_names);
1898 CTYPE ("cut-off lengths");
1899 RTYPE ("rvdw-switch", ir->rvdw_switch, 0.0);
1900 RTYPE ("rvdw", ir->rvdw, 1.0);
1901 CTYPE ("Apply long range dispersion corrections for Energy and Pressure");
1902 EETYPE("DispCorr", ir->eDispCorr, edispc_names);
1903 CTYPE ("Extension of the potential lookup tables beyond the cut-off");
1904 RTYPE ("table-extension", ir->tabext, 1.0);
1905 CTYPE ("Separate tables between energy group pairs");
1906 STYPE ("energygrp-table", is->egptable, NULL);
1907 CTYPE ("Spacing for the PME/PPPM FFT grid");
1908 RTYPE ("fourierspacing", ir->fourier_spacing, 0.12);
1909 CTYPE ("FFT grid size, when a value is 0 fourierspacing will be used");
1910 ITYPE ("fourier-nx", ir->nkx, 0);
1911 ITYPE ("fourier-ny", ir->nky, 0);
1912 ITYPE ("fourier-nz", ir->nkz, 0);
1913 CTYPE ("EWALD/PME/PPPM parameters");
1914 ITYPE ("pme-order", ir->pme_order, 4);
1915 RTYPE ("ewald-rtol", ir->ewald_rtol, 0.00001);
1916 RTYPE ("ewald-rtol-lj", ir->ewald_rtol_lj, 0.001);
1917 EETYPE("lj-pme-comb-rule", ir->ljpme_combination_rule, eljpme_names);
1918 EETYPE("ewald-geometry", ir->ewald_geometry, eewg_names);
1919 RTYPE ("epsilon-surface", ir->epsilon_surface, 0.0);
1920 EETYPE("optimize-fft", ir->bOptFFT, yesno_names);
1922 CCTYPE("IMPLICIT SOLVENT ALGORITHM");
1923 EETYPE("implicit-solvent", ir->implicit_solvent, eis_names);
1925 CCTYPE ("GENERALIZED BORN ELECTROSTATICS");
1926 CTYPE ("Algorithm for calculating Born radii");
1927 EETYPE("gb-algorithm", ir->gb_algorithm, egb_names);
1928 CTYPE ("Frequency of calculating the Born radii inside rlist");
1929 ITYPE ("nstgbradii", ir->nstgbradii, 1);
1930 CTYPE ("Cutoff for Born radii calculation; the contribution from atoms");
1931 CTYPE ("between rlist and rgbradii is updated every nstlist steps");
1932 RTYPE ("rgbradii", ir->rgbradii, 1.0);
1933 CTYPE ("Dielectric coefficient of the implicit solvent");
1934 RTYPE ("gb-epsilon-solvent", ir->gb_epsilon_solvent, 80.0);
1935 CTYPE ("Salt concentration in M for Generalized Born models");
1936 RTYPE ("gb-saltconc", ir->gb_saltconc, 0.0);
1937 CTYPE ("Scaling factors used in the OBC GB model. Default values are OBC(II)");
1938 RTYPE ("gb-obc-alpha", ir->gb_obc_alpha, 1.0);
1939 RTYPE ("gb-obc-beta", ir->gb_obc_beta, 0.8);
1940 RTYPE ("gb-obc-gamma", ir->gb_obc_gamma, 4.85);
1941 RTYPE ("gb-dielectric-offset", ir->gb_dielectric_offset, 0.009);
1942 EETYPE("sa-algorithm", ir->sa_algorithm, esa_names);
1943 CTYPE ("Surface tension (kJ/mol/nm^2) for the SA (nonpolar surface) part of GBSA");
1944 CTYPE ("The value -1 will set default value for Still/HCT/OBC GB-models.");
1945 RTYPE ("sa-surface-tension", ir->sa_surface_tension, -1);
1947 /* Coupling stuff */
1948 CCTYPE ("OPTIONS FOR WEAK COUPLING ALGORITHMS");
1949 CTYPE ("Temperature coupling");
1950 EETYPE("tcoupl", ir->etc, etcoupl_names);
1951 ITYPE ("nsttcouple", ir->nsttcouple, -1);
1952 ITYPE("nh-chain-length", ir->opts.nhchainlength, 10);
1953 EETYPE("print-nose-hoover-chain-variables", ir->bPrintNHChains, yesno_names);
1954 CTYPE ("Groups to couple separately");
1955 STYPE ("tc-grps", is->tcgrps, NULL);
1956 CTYPE ("Time constant (ps) and reference temperature (K)");
1957 STYPE ("tau-t", is->tau_t, NULL);
1958 STYPE ("ref-t", is->ref_t, NULL);
1959 CTYPE ("pressure coupling");
1960 EETYPE("pcoupl", ir->epc, epcoupl_names);
1961 EETYPE("pcoupltype", ir->epct, epcoupltype_names);
1962 ITYPE ("nstpcouple", ir->nstpcouple, -1);
1963 CTYPE ("Time constant (ps), compressibility (1/bar) and reference P (bar)");
1964 RTYPE ("tau-p", ir->tau_p, 1.0);
1965 STYPE ("compressibility", dumstr[0], NULL);
1966 STYPE ("ref-p", dumstr[1], NULL);
1967 CTYPE ("Scaling of reference coordinates, No, All or COM");
1968 EETYPE ("refcoord-scaling", ir->refcoord_scaling, erefscaling_names);
1971 CCTYPE ("OPTIONS FOR QMMM calculations");
1972 EETYPE("QMMM", ir->bQMMM, yesno_names);
1973 CTYPE ("Groups treated Quantum Mechanically");
1974 STYPE ("QMMM-grps", is->QMMM, NULL);
1975 CTYPE ("QM method");
1976 STYPE("QMmethod", is->QMmethod, NULL);
1977 CTYPE ("QMMM scheme");
1978 EETYPE("QMMMscheme", ir->QMMMscheme, eQMMMscheme_names);
1979 CTYPE ("QM basisset");
1980 STYPE("QMbasis", is->QMbasis, NULL);
1981 CTYPE ("QM charge");
1982 STYPE ("QMcharge", is->QMcharge, NULL);
1983 CTYPE ("QM multiplicity");
1984 STYPE ("QMmult", is->QMmult, NULL);
1985 CTYPE ("Surface Hopping");
1986 STYPE ("SH", is->bSH, NULL);
1987 CTYPE ("CAS space options");
1988 STYPE ("CASorbitals", is->CASorbitals, NULL);
1989 STYPE ("CASelectrons", is->CASelectrons, NULL);
1990 STYPE ("SAon", is->SAon, NULL);
1991 STYPE ("SAoff", is->SAoff, NULL);
1992 STYPE ("SAsteps", is->SAsteps, NULL);
1993 CTYPE ("Scale factor for MM charges");
1994 RTYPE ("MMChargeScaleFactor", ir->scalefactor, 1.0);
1995 CTYPE ("Optimization of QM subsystem");
1996 STYPE ("bOPT", is->bOPT, NULL);
1997 STYPE ("bTS", is->bTS, NULL);
1999 /* Simulated annealing */
2000 CCTYPE("SIMULATED ANNEALING");
2001 CTYPE ("Type of annealing for each temperature group (no/single/periodic)");
2002 STYPE ("annealing", is->anneal, NULL);
2003 CTYPE ("Number of time points to use for specifying annealing in each group");
2004 STYPE ("annealing-npoints", is->anneal_npoints, NULL);
2005 CTYPE ("List of times at the annealing points for each group");
2006 STYPE ("annealing-time", is->anneal_time, NULL);
2007 CTYPE ("Temp. at each annealing point, for each group.");
2008 STYPE ("annealing-temp", is->anneal_temp, NULL);
2011 CCTYPE ("GENERATE VELOCITIES FOR STARTUP RUN");
2012 EETYPE("gen-vel", opts->bGenVel, yesno_names);
2013 RTYPE ("gen-temp", opts->tempi, 300.0);
2014 ITYPE ("gen-seed", opts->seed, -1);
2017 CCTYPE ("OPTIONS FOR BONDS");
2018 EETYPE("constraints", opts->nshake, constraints);
2019 CTYPE ("Type of constraint algorithm");
2020 EETYPE("constraint-algorithm", ir->eConstrAlg, econstr_names);
2021 CTYPE ("Do not constrain the start configuration");
2022 EETYPE("continuation", ir->bContinuation, yesno_names);
2023 CTYPE ("Use successive overrelaxation to reduce the number of shake iterations");
2024 EETYPE("Shake-SOR", ir->bShakeSOR, yesno_names);
2025 CTYPE ("Relative tolerance of shake");
2026 RTYPE ("shake-tol", ir->shake_tol, 0.0001);
2027 CTYPE ("Highest order in the expansion of the constraint coupling matrix");
2028 ITYPE ("lincs-order", ir->nProjOrder, 4);
2029 CTYPE ("Number of iterations in the final step of LINCS. 1 is fine for");
2030 CTYPE ("normal simulations, but use 2 to conserve energy in NVE runs.");
2031 CTYPE ("For energy minimization with constraints it should be 4 to 8.");
2032 ITYPE ("lincs-iter", ir->nLincsIter, 1);
2033 CTYPE ("Lincs will write a warning to the stderr if in one step a bond");
2034 CTYPE ("rotates over more degrees than");
2035 RTYPE ("lincs-warnangle", ir->LincsWarnAngle, 30.0);
2036 CTYPE ("Convert harmonic bonds to morse potentials");
2037 EETYPE("morse", opts->bMorse, yesno_names);
2039 /* Energy group exclusions */
2040 CCTYPE ("ENERGY GROUP EXCLUSIONS");
2041 CTYPE ("Pairs of energy groups for which all non-bonded interactions are excluded");
2042 STYPE ("energygrp-excl", is->egpexcl, NULL);
2046 CTYPE ("Number of walls, type, atom types, densities and box-z scale factor for Ewald");
2047 ITYPE ("nwall", ir->nwall, 0);
2048 EETYPE("wall-type", ir->wall_type, ewt_names);
2049 RTYPE ("wall-r-linpot", ir->wall_r_linpot, -1);
2050 STYPE ("wall-atomtype", is->wall_atomtype, NULL);
2051 STYPE ("wall-density", is->wall_density, NULL);
2052 RTYPE ("wall-ewald-zfac", ir->wall_ewald_zfac, 3);
2055 CCTYPE("COM PULLING");
2056 CTYPE("Pull type: no, umbrella, constraint or constant-force");
2057 EETYPE("pull", ir->ePull, epull_names);
2058 if (ir->ePull != epullNO)
2061 is->pull_grp = read_pullparams(&ninp, &inp, ir->pull, &opts->pull_start, wi);
2064 /* Enforced rotation */
2065 CCTYPE("ENFORCED ROTATION");
2066 CTYPE("Enforced rotation: No or Yes");
2067 EETYPE("rotation", ir->bRot, yesno_names);
2071 is->rot_grp = read_rotparams(&ninp, &inp, ir->rot, wi);
2074 /* Interactive MD */
2076 CCTYPE("Group to display and/or manipulate in interactive MD session");
2077 STYPE ("IMD-group", is->imd_grp, NULL);
2078 if (is->imd_grp[0] != '\0')
2085 CCTYPE("NMR refinement stuff");
2086 CTYPE ("Distance restraints type: No, Simple or Ensemble");
2087 EETYPE("disre", ir->eDisre, edisre_names);
2088 CTYPE ("Force weighting of pairs in one distance restraint: Conservative or Equal");
2089 EETYPE("disre-weighting", ir->eDisreWeighting, edisreweighting_names);
2090 CTYPE ("Use sqrt of the time averaged times the instantaneous violation");
2091 EETYPE("disre-mixed", ir->bDisreMixed, yesno_names);
2092 RTYPE ("disre-fc", ir->dr_fc, 1000.0);
2093 RTYPE ("disre-tau", ir->dr_tau, 0.0);
2094 CTYPE ("Output frequency for pair distances to energy file");
2095 ITYPE ("nstdisreout", ir->nstdisreout, 100);
2096 CTYPE ("Orientation restraints: No or Yes");
2097 EETYPE("orire", opts->bOrire, yesno_names);
2098 CTYPE ("Orientation restraints force constant and tau for time averaging");
2099 RTYPE ("orire-fc", ir->orires_fc, 0.0);
2100 RTYPE ("orire-tau", ir->orires_tau, 0.0);
2101 STYPE ("orire-fitgrp", is->orirefitgrp, NULL);
2102 CTYPE ("Output frequency for trace(SD) and S to energy file");
2103 ITYPE ("nstorireout", ir->nstorireout, 100);
2105 /* free energy variables */
2106 CCTYPE ("Free energy variables");
2107 EETYPE("free-energy", ir->efep, efep_names);
2108 STYPE ("couple-moltype", is->couple_moltype, NULL);
2109 EETYPE("couple-lambda0", opts->couple_lam0, couple_lam);
2110 EETYPE("couple-lambda1", opts->couple_lam1, couple_lam);
2111 EETYPE("couple-intramol", opts->bCoupleIntra, yesno_names);
2113 RTYPE ("init-lambda", fep->init_lambda, -1); /* start with -1 so
2115 it was not entered */
2116 ITYPE ("init-lambda-state", fep->init_fep_state, -1);
2117 RTYPE ("delta-lambda", fep->delta_lambda, 0.0);
2118 ITYPE ("nstdhdl", fep->nstdhdl, 50);
2119 STYPE ("fep-lambdas", is->fep_lambda[efptFEP], NULL);
2120 STYPE ("mass-lambdas", is->fep_lambda[efptMASS], NULL);
2121 STYPE ("coul-lambdas", is->fep_lambda[efptCOUL], NULL);
2122 STYPE ("vdw-lambdas", is->fep_lambda[efptVDW], NULL);
2123 STYPE ("bonded-lambdas", is->fep_lambda[efptBONDED], NULL);
2124 STYPE ("restraint-lambdas", is->fep_lambda[efptRESTRAINT], NULL);
2125 STYPE ("temperature-lambdas", is->fep_lambda[efptTEMPERATURE], NULL);
2126 ITYPE ("calc-lambda-neighbors", fep->lambda_neighbors, 1);
2127 STYPE ("init-lambda-weights", is->lambda_weights, NULL);
2128 EETYPE("dhdl-print-energy", fep->bPrintEnergy, yesno_names);
2129 RTYPE ("sc-alpha", fep->sc_alpha, 0.0);
2130 ITYPE ("sc-power", fep->sc_power, 1);
2131 RTYPE ("sc-r-power", fep->sc_r_power, 6.0);
2132 RTYPE ("sc-sigma", fep->sc_sigma, 0.3);
2133 EETYPE("sc-coul", fep->bScCoul, yesno_names);
2134 ITYPE ("dh_hist_size", fep->dh_hist_size, 0);
2135 RTYPE ("dh_hist_spacing", fep->dh_hist_spacing, 0.1);
2136 EETYPE("separate-dhdl-file", fep->separate_dhdl_file,
2137 separate_dhdl_file_names);
2138 EETYPE("dhdl-derivatives", fep->dhdl_derivatives, dhdl_derivatives_names);
2139 ITYPE ("dh_hist_size", fep->dh_hist_size, 0);
2140 RTYPE ("dh_hist_spacing", fep->dh_hist_spacing, 0.1);
2142 /* Non-equilibrium MD stuff */
2143 CCTYPE("Non-equilibrium MD stuff");
2144 STYPE ("acc-grps", is->accgrps, NULL);
2145 STYPE ("accelerate", is->acc, NULL);
2146 STYPE ("freezegrps", is->freeze, NULL);
2147 STYPE ("freezedim", is->frdim, NULL);
2148 RTYPE ("cos-acceleration", ir->cos_accel, 0);
2149 STYPE ("deform", is->deform, NULL);
2151 /* simulated tempering variables */
2152 CCTYPE("simulated tempering variables");
2153 EETYPE("simulated-tempering", ir->bSimTemp, yesno_names);
2154 EETYPE("simulated-tempering-scaling", ir->simtempvals->eSimTempScale, esimtemp_names);
2155 RTYPE("sim-temp-low", ir->simtempvals->simtemp_low, 300.0);
2156 RTYPE("sim-temp-high", ir->simtempvals->simtemp_high, 300.0);
2158 /* expanded ensemble variables */
2159 if (ir->efep == efepEXPANDED || ir->bSimTemp)
2161 read_expandedparams(&ninp, &inp, expand, wi);
2164 /* Electric fields */
2165 CCTYPE("Electric fields");
2166 CTYPE ("Format is number of terms (int) and for all terms an amplitude (real)");
2167 CTYPE ("and a phase angle (real)");
2168 STYPE ("E-x", is->efield_x, NULL);
2169 STYPE ("E-xt", is->efield_xt, NULL);
2170 STYPE ("E-y", is->efield_y, NULL);
2171 STYPE ("E-yt", is->efield_yt, NULL);
2172 STYPE ("E-z", is->efield_z, NULL);
2173 STYPE ("E-zt", is->efield_zt, NULL);
2175 CCTYPE("Ion/water position swapping for computational electrophysiology setups");
2176 CTYPE("Swap positions along direction: no, X, Y, Z");
2177 EETYPE("swapcoords", ir->eSwapCoords, eSwapTypes_names);
2178 if (ir->eSwapCoords != eswapNO)
2181 CTYPE("Swap attempt frequency");
2182 ITYPE("swap-frequency", ir->swap->nstswap, 1);
2183 CTYPE("Two index groups that contain the compartment-partitioning atoms");
2184 STYPE("split-group0", splitgrp0, NULL);
2185 STYPE("split-group1", splitgrp1, NULL);
2186 CTYPE("Use center of mass of split groups (yes/no), otherwise center of geometry is used");
2187 EETYPE("massw-split0", ir->swap->massw_split[0], yesno_names);
2188 EETYPE("massw-split1", ir->swap->massw_split[1], yesno_names);
2190 CTYPE("Group name of ions that can be exchanged with solvent molecules");
2191 STYPE("swap-group", swapgrp, NULL);
2192 CTYPE("Group name of solvent molecules");
2193 STYPE("solvent-group", solgrp, NULL);
2195 CTYPE("Split cylinder: radius, upper and lower extension (nm) (this will define the channels)");
2196 CTYPE("Note that the split cylinder settings do not have an influence on the swapping protocol,");
2197 CTYPE("however, if correctly defined, the ion permeation events are counted per channel");
2198 RTYPE("cyl0-r", ir->swap->cyl0r, 2.0);
2199 RTYPE("cyl0-up", ir->swap->cyl0u, 1.0);
2200 RTYPE("cyl0-down", ir->swap->cyl0l, 1.0);
2201 RTYPE("cyl1-r", ir->swap->cyl1r, 2.0);
2202 RTYPE("cyl1-up", ir->swap->cyl1u, 1.0);
2203 RTYPE("cyl1-down", ir->swap->cyl1l, 1.0);
2205 CTYPE("Average the number of ions per compartment over these many swap attempt steps");
2206 ITYPE("coupl-steps", ir->swap->nAverage, 10);
2207 CTYPE("Requested number of anions and cations for each of the two compartments");
2208 CTYPE("-1 means fix the numbers as found in time step 0");
2209 ITYPE("anionsA", ir->swap->nanions[0], -1);
2210 ITYPE("cationsA", ir->swap->ncations[0], -1);
2211 ITYPE("anionsB", ir->swap->nanions[1], -1);
2212 ITYPE("cationsB", ir->swap->ncations[1], -1);
2213 CTYPE("Start to swap ions if threshold difference to requested count is reached");
2214 RTYPE("threshold", ir->swap->threshold, 1.0);
2217 /* AdResS defined thingies */
2218 CCTYPE ("AdResS parameters");
2219 EETYPE("adress", ir->bAdress, yesno_names);
2222 snew(ir->adress, 1);
2223 read_adressparams(&ninp, &inp, ir->adress, wi);
2226 /* User defined thingies */
2227 CCTYPE ("User defined thingies");
2228 STYPE ("user1-grps", is->user1, NULL);
2229 STYPE ("user2-grps", is->user2, NULL);
2230 ITYPE ("userint1", ir->userint1, 0);
2231 ITYPE ("userint2", ir->userint2, 0);
2232 ITYPE ("userint3", ir->userint3, 0);
2233 ITYPE ("userint4", ir->userint4, 0);
2234 RTYPE ("userreal1", ir->userreal1, 0);
2235 RTYPE ("userreal2", ir->userreal2, 0);
2236 RTYPE ("userreal3", ir->userreal3, 0);
2237 RTYPE ("userreal4", ir->userreal4, 0);
2240 write_inpfile(mdparout, ninp, inp, FALSE, wi);
2241 for (i = 0; (i < ninp); i++)
2244 sfree(inp[i].value);
2248 /* Process options if necessary */
2249 for (m = 0; m < 2; m++)
2251 for (i = 0; i < 2*DIM; i++)
2260 if (sscanf(dumstr[m], "%lf", &(dumdub[m][XX])) != 1)
2262 warning_error(wi, "Pressure coupling not enough values (I need 1)");
2264 dumdub[m][YY] = dumdub[m][ZZ] = dumdub[m][XX];
2266 case epctSEMIISOTROPIC:
2267 case epctSURFACETENSION:
2268 if (sscanf(dumstr[m], "%lf%lf",
2269 &(dumdub[m][XX]), &(dumdub[m][ZZ])) != 2)
2271 warning_error(wi, "Pressure coupling not enough values (I need 2)");
2273 dumdub[m][YY] = dumdub[m][XX];
2275 case epctANISOTROPIC:
2276 if (sscanf(dumstr[m], "%lf%lf%lf%lf%lf%lf",
2277 &(dumdub[m][XX]), &(dumdub[m][YY]), &(dumdub[m][ZZ]),
2278 &(dumdub[m][3]), &(dumdub[m][4]), &(dumdub[m][5])) != 6)
2280 warning_error(wi, "Pressure coupling not enough values (I need 6)");
2284 gmx_fatal(FARGS, "Pressure coupling type %s not implemented yet",
2285 epcoupltype_names[ir->epct]);
2289 clear_mat(ir->ref_p);
2290 clear_mat(ir->compress);
2291 for (i = 0; i < DIM; i++)
2293 ir->ref_p[i][i] = dumdub[1][i];
2294 ir->compress[i][i] = dumdub[0][i];
2296 if (ir->epct == epctANISOTROPIC)
2298 ir->ref_p[XX][YY] = dumdub[1][3];
2299 ir->ref_p[XX][ZZ] = dumdub[1][4];
2300 ir->ref_p[YY][ZZ] = dumdub[1][5];
2301 if (ir->ref_p[XX][YY] != 0 && ir->ref_p[XX][ZZ] != 0 && ir->ref_p[YY][ZZ] != 0)
2303 warning(wi, "All off-diagonal reference pressures are non-zero. Are you sure you want to apply a threefold shear stress?\n");
2305 ir->compress[XX][YY] = dumdub[0][3];
2306 ir->compress[XX][ZZ] = dumdub[0][4];
2307 ir->compress[YY][ZZ] = dumdub[0][5];
2308 for (i = 0; i < DIM; i++)
2310 for (m = 0; m < i; m++)
2312 ir->ref_p[i][m] = ir->ref_p[m][i];
2313 ir->compress[i][m] = ir->compress[m][i];
2318 if (ir->comm_mode == ecmNO)
2323 opts->couple_moltype = NULL;
2324 if (strlen(is->couple_moltype) > 0)
2326 if (ir->efep != efepNO)
2328 opts->couple_moltype = strdup(is->couple_moltype);
2329 if (opts->couple_lam0 == opts->couple_lam1)
2331 warning(wi, "The lambda=0 and lambda=1 states for coupling are identical");
2333 if (ir->eI == eiMD && (opts->couple_lam0 == ecouplamNONE ||
2334 opts->couple_lam1 == ecouplamNONE))
2336 warning(wi, "For proper sampling of the (nearly) decoupled state, stochastic dynamics should be used");
2341 warning(wi, "Can not couple a molecule with free_energy = no");
2344 /* FREE ENERGY AND EXPANDED ENSEMBLE OPTIONS */
2345 if (ir->efep != efepNO)
2347 if (fep->delta_lambda > 0)
2349 ir->efep = efepSLOWGROWTH;
2355 fep->bPrintEnergy = TRUE;
2356 /* always print out the energy to dhdl if we are doing expanded ensemble, since we need the total energy
2357 if the temperature is changing. */
2360 if ((ir->efep != efepNO) || ir->bSimTemp)
2362 ir->bExpanded = FALSE;
2363 if ((ir->efep == efepEXPANDED) || ir->bSimTemp)
2365 ir->bExpanded = TRUE;
2367 do_fep_params(ir, is->fep_lambda, is->lambda_weights);
2368 if (ir->bSimTemp) /* done after fep params */
2370 do_simtemp_params(ir);
2375 ir->fepvals->n_lambda = 0;
2378 /* WALL PARAMETERS */
2380 do_wall_params(ir, is->wall_atomtype, is->wall_density, opts);
2382 /* ORIENTATION RESTRAINT PARAMETERS */
2384 if (opts->bOrire && str_nelem(is->orirefitgrp, MAXPTR, NULL) != 1)
2386 warning_error(wi, "ERROR: Need one orientation restraint fit group\n");
2389 /* DEFORMATION PARAMETERS */
2391 clear_mat(ir->deform);
2392 for (i = 0; i < 6; i++)
2396 m = sscanf(is->deform, "%lf %lf %lf %lf %lf %lf",
2397 &(dumdub[0][0]), &(dumdub[0][1]), &(dumdub[0][2]),
2398 &(dumdub[0][3]), &(dumdub[0][4]), &(dumdub[0][5]));
2399 for (i = 0; i < 3; i++)
2401 ir->deform[i][i] = dumdub[0][i];
2403 ir->deform[YY][XX] = dumdub[0][3];
2404 ir->deform[ZZ][XX] = dumdub[0][4];
2405 ir->deform[ZZ][YY] = dumdub[0][5];
2406 if (ir->epc != epcNO)
2408 for (i = 0; i < 3; i++)
2410 for (j = 0; j <= i; j++)
2412 if (ir->deform[i][j] != 0 && ir->compress[i][j] != 0)
2414 warning_error(wi, "A box element has deform set and compressibility > 0");
2418 for (i = 0; i < 3; i++)
2420 for (j = 0; j < i; j++)
2422 if (ir->deform[i][j] != 0)
2424 for (m = j; m < DIM; m++)
2426 if (ir->compress[m][j] != 0)
2428 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.");
2429 warning(wi, warn_buf);
2437 /* Ion/water position swapping checks */
2438 if (ir->eSwapCoords != eswapNO)
2440 if (ir->swap->nstswap < 1)
2442 warning_error(wi, "swap_frequency must be 1 or larger when ion swapping is requested");
2444 if (ir->swap->nAverage < 1)
2446 warning_error(wi, "coupl_steps must be 1 or larger.\n");
2448 if (ir->swap->threshold < 1.0)
2450 warning_error(wi, "Ion count threshold must be at least 1.\n");
2458 static int search_QMstring(const char *s, int ng, const char *gn[])
2460 /* same as normal search_string, but this one searches QM strings */
2463 for (i = 0; (i < ng); i++)
2465 if (gmx_strcasecmp(s, gn[i]) == 0)
2471 gmx_fatal(FARGS, "this QM method or basisset (%s) is not implemented\n!", s);
2475 } /* search_QMstring */
2477 /* We would like gn to be const as well, but C doesn't allow this */
2478 int search_string(const char *s, int ng, char *gn[])
2482 for (i = 0; (i < ng); i++)
2484 if (gmx_strcasecmp(s, gn[i]) == 0)
2491 "Group %s referenced in the .mdp file was not found in the index file.\n"
2492 "Group names must match either [moleculetype] names or custom index group\n"
2493 "names, in which case you must supply an index file to the '-n' option\n"
2500 static gmx_bool do_numbering(int natoms, gmx_groups_t *groups, int ng, char *ptrs[],
2501 t_blocka *block, char *gnames[],
2502 int gtype, int restnm,
2503 int grptp, gmx_bool bVerbose,
2506 unsigned short *cbuf;
2507 t_grps *grps = &(groups->grps[gtype]);
2508 int i, j, gid, aj, ognr, ntot = 0;
2511 char warn_buf[STRLEN];
2515 fprintf(debug, "Starting numbering %d groups of type %d\n", ng, gtype);
2518 title = gtypes[gtype];
2521 /* Mark all id's as not set */
2522 for (i = 0; (i < natoms); i++)
2527 snew(grps->nm_ind, ng+1); /* +1 for possible rest group */
2528 for (i = 0; (i < ng); i++)
2530 /* Lookup the group name in the block structure */
2531 gid = search_string(ptrs[i], block->nr, gnames);
2532 if ((grptp != egrptpONE) || (i == 0))
2534 grps->nm_ind[grps->nr++] = gid;
2538 fprintf(debug, "Found gid %d for group %s\n", gid, ptrs[i]);
2541 /* Now go over the atoms in the group */
2542 for (j = block->index[gid]; (j < block->index[gid+1]); j++)
2547 /* Range checking */
2548 if ((aj < 0) || (aj >= natoms))
2550 gmx_fatal(FARGS, "Invalid atom number %d in indexfile", aj);
2552 /* Lookup up the old group number */
2556 gmx_fatal(FARGS, "Atom %d in multiple %s groups (%d and %d)",
2557 aj+1, title, ognr+1, i+1);
2561 /* Store the group number in buffer */
2562 if (grptp == egrptpONE)
2575 /* Now check whether we have done all atoms */
2579 if (grptp == egrptpALL)
2581 gmx_fatal(FARGS, "%d atoms are not part of any of the %s groups",
2582 natoms-ntot, title);
2584 else if (grptp == egrptpPART)
2586 sprintf(warn_buf, "%d atoms are not part of any of the %s groups",
2587 natoms-ntot, title);
2588 warning_note(wi, warn_buf);
2590 /* Assign all atoms currently unassigned to a rest group */
2591 for (j = 0; (j < natoms); j++)
2593 if (cbuf[j] == NOGID)
2599 if (grptp != egrptpPART)
2604 "Making dummy/rest group for %s containing %d elements\n",
2605 title, natoms-ntot);
2607 /* Add group name "rest" */
2608 grps->nm_ind[grps->nr] = restnm;
2610 /* Assign the rest name to all atoms not currently assigned to a group */
2611 for (j = 0; (j < natoms); j++)
2613 if (cbuf[j] == NOGID)
2622 if (grps->nr == 1 && (ntot == 0 || ntot == natoms))
2624 /* All atoms are part of one (or no) group, no index required */
2625 groups->ngrpnr[gtype] = 0;
2626 groups->grpnr[gtype] = NULL;
2630 groups->ngrpnr[gtype] = natoms;
2631 snew(groups->grpnr[gtype], natoms);
2632 for (j = 0; (j < natoms); j++)
2634 groups->grpnr[gtype][j] = cbuf[j];
2640 return (bRest && grptp == egrptpPART);
2643 static void calc_nrdf(gmx_mtop_t *mtop, t_inputrec *ir, char **gnames)
2646 gmx_groups_t *groups;
2648 int natoms, ai, aj, i, j, d, g, imin, jmin;
2650 int *nrdf2, *na_vcm, na_tot;
2651 double *nrdf_tc, *nrdf_vcm, nrdf_uc, n_sub = 0;
2652 gmx_mtop_atomloop_all_t aloop;
2654 int mb, mol, ftype, as;
2655 gmx_molblock_t *molb;
2656 gmx_moltype_t *molt;
2659 * First calc 3xnr-atoms for each group
2660 * then subtract half a degree of freedom for each constraint
2662 * Only atoms and nuclei contribute to the degrees of freedom...
2667 groups = &mtop->groups;
2668 natoms = mtop->natoms;
2670 /* Allocate one more for a possible rest group */
2671 /* We need to sum degrees of freedom into doubles,
2672 * since floats give too low nrdf's above 3 million atoms.
2674 snew(nrdf_tc, groups->grps[egcTC].nr+1);
2675 snew(nrdf_vcm, groups->grps[egcVCM].nr+1);
2676 snew(na_vcm, groups->grps[egcVCM].nr+1);
2678 for (i = 0; i < groups->grps[egcTC].nr; i++)
2682 for (i = 0; i < groups->grps[egcVCM].nr+1; i++)
2687 snew(nrdf2, natoms);
2688 aloop = gmx_mtop_atomloop_all_init(mtop);
2689 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
2692 if (atom->ptype == eptAtom || atom->ptype == eptNucleus)
2694 g = ggrpnr(groups, egcFREEZE, i);
2695 /* Double count nrdf for particle i */
2696 for (d = 0; d < DIM; d++)
2698 if (opts->nFreeze[g][d] == 0)
2703 nrdf_tc [ggrpnr(groups, egcTC, i)] += 0.5*nrdf2[i];
2704 nrdf_vcm[ggrpnr(groups, egcVCM, i)] += 0.5*nrdf2[i];
2709 for (mb = 0; mb < mtop->nmolblock; mb++)
2711 molb = &mtop->molblock[mb];
2712 molt = &mtop->moltype[molb->type];
2713 atom = molt->atoms.atom;
2714 for (mol = 0; mol < molb->nmol; mol++)
2716 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
2718 ia = molt->ilist[ftype].iatoms;
2719 for (i = 0; i < molt->ilist[ftype].nr; )
2721 /* Subtract degrees of freedom for the constraints,
2722 * if the particles still have degrees of freedom left.
2723 * If one of the particles is a vsite or a shell, then all
2724 * constraint motion will go there, but since they do not
2725 * contribute to the constraints the degrees of freedom do not
2730 if (((atom[ia[1]].ptype == eptNucleus) ||
2731 (atom[ia[1]].ptype == eptAtom)) &&
2732 ((atom[ia[2]].ptype == eptNucleus) ||
2733 (atom[ia[2]].ptype == eptAtom)))
2751 imin = min(imin, nrdf2[ai]);
2752 jmin = min(jmin, nrdf2[aj]);
2755 nrdf_tc [ggrpnr(groups, egcTC, ai)] -= 0.5*imin;
2756 nrdf_tc [ggrpnr(groups, egcTC, aj)] -= 0.5*jmin;
2757 nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2758 nrdf_vcm[ggrpnr(groups, egcVCM, aj)] -= 0.5*jmin;
2760 ia += interaction_function[ftype].nratoms+1;
2761 i += interaction_function[ftype].nratoms+1;
2764 ia = molt->ilist[F_SETTLE].iatoms;
2765 for (i = 0; i < molt->ilist[F_SETTLE].nr; )
2767 /* Subtract 1 dof from every atom in the SETTLE */
2768 for (j = 0; j < 3; j++)
2771 imin = min(2, nrdf2[ai]);
2773 nrdf_tc [ggrpnr(groups, egcTC, ai)] -= 0.5*imin;
2774 nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2779 as += molt->atoms.nr;
2783 if (ir->ePull == epullCONSTRAINT)
2785 /* Correct nrdf for the COM constraints.
2786 * We correct using the TC and VCM group of the first atom
2787 * in the reference and pull group. If atoms in one pull group
2788 * belong to different TC or VCM groups it is anyhow difficult
2789 * to determine the optimal nrdf assignment.
2793 for (i = 0; i < pull->ncoord; i++)
2797 for (j = 0; j < 2; j++)
2799 const t_pull_group *pgrp;
2801 pgrp = &pull->group[pull->coord[i].group[j]];
2805 /* Subtract 1/2 dof from each group */
2807 nrdf_tc [ggrpnr(groups, egcTC, ai)] -= 0.5*imin;
2808 nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2809 if (nrdf_tc[ggrpnr(groups, egcTC, ai)] < 0)
2811 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[ggrpnr(groups, egcTC, ai)]]);
2816 /* We need to subtract the whole DOF from group j=1 */
2823 if (ir->nstcomm != 0)
2825 /* Subtract 3 from the number of degrees of freedom in each vcm group
2826 * when com translation is removed and 6 when rotation is removed
2829 switch (ir->comm_mode)
2832 n_sub = ndof_com(ir);
2839 gmx_incons("Checking comm_mode");
2842 for (i = 0; i < groups->grps[egcTC].nr; i++)
2844 /* Count the number of atoms of TC group i for every VCM group */
2845 for (j = 0; j < groups->grps[egcVCM].nr+1; j++)
2850 for (ai = 0; ai < natoms; ai++)
2852 if (ggrpnr(groups, egcTC, ai) == i)
2854 na_vcm[ggrpnr(groups, egcVCM, ai)]++;
2858 /* Correct for VCM removal according to the fraction of each VCM
2859 * group present in this TC group.
2861 nrdf_uc = nrdf_tc[i];
2864 fprintf(debug, "T-group[%d] nrdf_uc = %g, n_sub = %g\n",
2868 for (j = 0; j < groups->grps[egcVCM].nr+1; j++)
2870 if (nrdf_vcm[j] > n_sub)
2872 nrdf_tc[i] += nrdf_uc*((double)na_vcm[j]/(double)na_tot)*
2873 (nrdf_vcm[j] - n_sub)/nrdf_vcm[j];
2877 fprintf(debug, " nrdf_vcm[%d] = %g, nrdf = %g\n",
2878 j, nrdf_vcm[j], nrdf_tc[i]);
2883 for (i = 0; (i < groups->grps[egcTC].nr); i++)
2885 opts->nrdf[i] = nrdf_tc[i];
2886 if (opts->nrdf[i] < 0)
2891 "Number of degrees of freedom in T-Coupling group %s is %.2f\n",
2892 gnames[groups->grps[egcTC].nm_ind[i]], opts->nrdf[i]);
2901 static void decode_cos(char *s, t_cosines *cosine)
2904 char format[STRLEN], f1[STRLEN];
2916 sscanf(t, "%d", &(cosine->n));
2923 snew(cosine->a, cosine->n);
2924 snew(cosine->phi, cosine->n);
2926 sprintf(format, "%%*d");
2927 for (i = 0; (i < cosine->n); i++)
2930 strcat(f1, "%lf%lf");
2931 if (sscanf(t, f1, &a, &phi) < 2)
2933 gmx_fatal(FARGS, "Invalid input for electric field shift: '%s'", t);
2936 cosine->phi[i] = phi;
2937 strcat(format, "%*lf%*lf");
2944 static gmx_bool do_egp_flag(t_inputrec *ir, gmx_groups_t *groups,
2945 const char *option, const char *val, int flag)
2947 /* The maximum number of energy group pairs would be MAXPTR*(MAXPTR+1)/2.
2948 * But since this is much larger than STRLEN, such a line can not be parsed.
2949 * The real maximum is the number of names that fit in a string: STRLEN/2.
2951 #define EGP_MAX (STRLEN/2)
2952 int nelem, i, j, k, nr;
2953 char *names[EGP_MAX];
2957 gnames = groups->grpname;
2959 nelem = str_nelem(val, EGP_MAX, names);
2962 gmx_fatal(FARGS, "The number of groups for %s is odd", option);
2964 nr = groups->grps[egcENER].nr;
2966 for (i = 0; i < nelem/2; i++)
2970 gmx_strcasecmp(names[2*i], *(gnames[groups->grps[egcENER].nm_ind[j]])))
2976 gmx_fatal(FARGS, "%s in %s is not an energy group\n",
2977 names[2*i], option);
2981 gmx_strcasecmp(names[2*i+1], *(gnames[groups->grps[egcENER].nm_ind[k]])))
2987 gmx_fatal(FARGS, "%s in %s is not an energy group\n",
2988 names[2*i+1], option);
2990 if ((j < nr) && (k < nr))
2992 ir->opts.egp_flags[nr*j+k] |= flag;
2993 ir->opts.egp_flags[nr*k+j] |= flag;
3002 static void make_swap_groups(
3011 int ig = -1, i = 0, j;
3015 /* Just a quick check here, more thorough checks are in mdrun */
3016 if (strcmp(splitg0name, splitg1name) == 0)
3018 gmx_fatal(FARGS, "The split groups can not both be '%s'.", splitg0name);
3021 /* First get the swap group index atoms */
3022 ig = search_string(swapgname, grps->nr, gnames);
3023 swap->nat = grps->index[ig+1] - grps->index[ig];
3026 fprintf(stderr, "Swap group '%s' contains %d atoms.\n", swapgname, swap->nat);
3027 snew(swap->ind, swap->nat);
3028 for (i = 0; i < swap->nat; i++)
3030 swap->ind[i] = grps->a[grps->index[ig]+i];
3035 gmx_fatal(FARGS, "You defined an empty group of atoms for swapping.");
3038 /* Now do so for the split groups */
3039 for (j = 0; j < 2; j++)
3043 splitg = splitg0name;
3047 splitg = splitg1name;
3050 ig = search_string(splitg, grps->nr, gnames);
3051 swap->nat_split[j] = grps->index[ig+1] - grps->index[ig];
3052 if (swap->nat_split[j] > 0)
3054 fprintf(stderr, "Split group %d '%s' contains %d atom%s.\n",
3055 j, splitg, swap->nat_split[j], (swap->nat_split[j] > 1) ? "s" : "");
3056 snew(swap->ind_split[j], swap->nat_split[j]);
3057 for (i = 0; i < swap->nat_split[j]; i++)
3059 swap->ind_split[j][i] = grps->a[grps->index[ig]+i];
3064 gmx_fatal(FARGS, "Split group %d has to contain at least 1 atom!", j);
3068 /* Now get the solvent group index atoms */
3069 ig = search_string(solgname, grps->nr, gnames);
3070 swap->nat_sol = grps->index[ig+1] - grps->index[ig];
3071 if (swap->nat_sol > 0)
3073 fprintf(stderr, "Solvent group '%s' contains %d atoms.\n", solgname, swap->nat_sol);
3074 snew(swap->ind_sol, swap->nat_sol);
3075 for (i = 0; i < swap->nat_sol; i++)
3077 swap->ind_sol[i] = grps->a[grps->index[ig]+i];
3082 gmx_fatal(FARGS, "You defined an empty group of solvent. Cannot exchange ions.");
3087 void make_IMD_group(t_IMD *IMDgroup, char *IMDgname, t_blocka *grps, char **gnames)
3092 ig = search_string(IMDgname, grps->nr, gnames);
3093 IMDgroup->nat = grps->index[ig+1] - grps->index[ig];
3095 if (IMDgroup->nat > 0)
3097 fprintf(stderr, "Group '%s' with %d atoms can be activated for interactive molecular dynamics (IMD).\n",
3098 IMDgname, IMDgroup->nat);
3099 snew(IMDgroup->ind, IMDgroup->nat);
3100 for (i = 0; i < IMDgroup->nat; i++)
3102 IMDgroup->ind[i] = grps->a[grps->index[ig]+i];
3108 void do_index(const char* mdparin, const char *ndx,
3111 t_inputrec *ir, rvec *v,
3115 gmx_groups_t *groups;
3119 char warnbuf[STRLEN], **gnames;
3120 int nr, ntcg, ntau_t, nref_t, nacc, nofg, nSA, nSA_points, nSA_time, nSA_temp;
3123 int nacg, nfreeze, nfrdim, nenergy, nvcm, nuser;
3124 char *ptr1[MAXPTR], *ptr2[MAXPTR], *ptr3[MAXPTR];
3125 int i, j, k, restnm;
3127 gmx_bool bExcl, bTable, bSetTCpar, bAnneal, bRest;
3128 int nQMmethod, nQMbasis, nQMcharge, nQMmult, nbSH, nCASorb, nCASelec,
3129 nSAon, nSAoff, nSAsteps, nQMg, nbOPT, nbTS;
3130 char warn_buf[STRLEN];
3134 fprintf(stderr, "processing index file...\n");
3140 snew(grps->index, 1);
3142 atoms_all = gmx_mtop_global_atoms(mtop);
3143 analyse(&atoms_all, grps, &gnames, FALSE, TRUE);
3144 free_t_atoms(&atoms_all, FALSE);
3148 grps = init_index(ndx, &gnames);
3151 groups = &mtop->groups;
3152 natoms = mtop->natoms;
3153 symtab = &mtop->symtab;
3155 snew(groups->grpname, grps->nr+1);
3157 for (i = 0; (i < grps->nr); i++)
3159 groups->grpname[i] = put_symtab(symtab, gnames[i]);
3161 groups->grpname[i] = put_symtab(symtab, "rest");
3163 srenew(gnames, grps->nr+1);
3164 gnames[restnm] = *(groups->grpname[i]);
3165 groups->ngrpname = grps->nr+1;
3167 set_warning_line(wi, mdparin, -1);
3169 ntau_t = str_nelem(is->tau_t, MAXPTR, ptr1);
3170 nref_t = str_nelem(is->ref_t, MAXPTR, ptr2);
3171 ntcg = str_nelem(is->tcgrps, MAXPTR, ptr3);
3172 if ((ntau_t != ntcg) || (nref_t != ntcg))
3174 gmx_fatal(FARGS, "Invalid T coupling input: %d groups, %d ref-t values and "
3175 "%d tau-t values", ntcg, nref_t, ntau_t);
3178 bSetTCpar = (ir->etc || EI_SD(ir->eI) || ir->eI == eiBD || EI_TPI(ir->eI));
3179 do_numbering(natoms, groups, ntcg, ptr3, grps, gnames, egcTC,
3180 restnm, bSetTCpar ? egrptpALL : egrptpALL_GENREST, bVerbose, wi);
3181 nr = groups->grps[egcTC].nr;
3183 snew(ir->opts.nrdf, nr);
3184 snew(ir->opts.tau_t, nr);
3185 snew(ir->opts.ref_t, nr);
3186 if (ir->eI == eiBD && ir->bd_fric == 0)
3188 fprintf(stderr, "bd-fric=0, so tau-t will be used as the inverse friction constant(s)\n");
3195 gmx_fatal(FARGS, "Not enough ref-t and tau-t values!");
3199 for (i = 0; (i < nr); i++)
3201 ir->opts.tau_t[i] = strtod(ptr1[i], NULL);
3202 if ((ir->eI == eiBD || ir->eI == eiSD2) && ir->opts.tau_t[i] <= 0)
3204 sprintf(warn_buf, "With integrator %s tau-t should be larger than 0", ei_names[ir->eI]);
3205 warning_error(wi, warn_buf);
3208 if (ir->etc != etcVRESCALE && ir->opts.tau_t[i] == 0)
3210 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.");
3213 if (ir->opts.tau_t[i] >= 0)
3215 tau_min = min(tau_min, ir->opts.tau_t[i]);
3218 if (ir->etc != etcNO && ir->nsttcouple == -1)
3220 ir->nsttcouple = ir_optimal_nsttcouple(ir);
3225 if ((ir->etc == etcNOSEHOOVER) && (ir->epc == epcBERENDSEN))
3227 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");
3229 if ((ir->epc == epcMTTK) && (ir->etc > etcNO))
3231 if (ir->nstpcouple != ir->nsttcouple)
3233 int mincouple = min(ir->nstpcouple, ir->nsttcouple);
3234 ir->nstpcouple = ir->nsttcouple = mincouple;
3235 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);
3236 warning_note(wi, warn_buf);
3240 /* velocity verlet with averaged kinetic energy KE = 0.5*(v(t+1/2) - v(t-1/2)) is implemented
3241 primarily for testing purposes, and does not work with temperature coupling other than 1 */
3243 if (ETC_ANDERSEN(ir->etc))
3245 if (ir->nsttcouple != 1)
3248 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");
3249 warning_note(wi, warn_buf);
3252 nstcmin = tcouple_min_integration_steps(ir->etc);
3255 if (tau_min/(ir->delta_t*ir->nsttcouple) < nstcmin)
3257 sprintf(warn_buf, "For proper integration of the %s thermostat, tau-t (%g) should be at least %d times larger than nsttcouple*dt (%g)",
3258 ETCOUPLTYPE(ir->etc),
3260 ir->nsttcouple*ir->delta_t);
3261 warning(wi, warn_buf);
3264 for (i = 0; (i < nr); i++)
3266 ir->opts.ref_t[i] = strtod(ptr2[i], NULL);
3267 if (ir->opts.ref_t[i] < 0)
3269 gmx_fatal(FARGS, "ref-t for group %d negative", i);
3272 /* set the lambda mc temperature to the md integrator temperature (which should be defined
3273 if we are in this conditional) if mc_temp is negative */
3274 if (ir->expandedvals->mc_temp < 0)
3276 ir->expandedvals->mc_temp = ir->opts.ref_t[0]; /*for now, set to the first reft */
3280 /* Simulated annealing for each group. There are nr groups */
3281 nSA = str_nelem(is->anneal, MAXPTR, ptr1);
3282 if (nSA == 1 && (ptr1[0][0] == 'n' || ptr1[0][0] == 'N'))
3286 if (nSA > 0 && nSA != nr)
3288 gmx_fatal(FARGS, "Not enough annealing values: %d (for %d groups)\n", nSA, nr);
3292 snew(ir->opts.annealing, nr);
3293 snew(ir->opts.anneal_npoints, nr);
3294 snew(ir->opts.anneal_time, nr);
3295 snew(ir->opts.anneal_temp, nr);
3296 for (i = 0; i < nr; i++)
3298 ir->opts.annealing[i] = eannNO;
3299 ir->opts.anneal_npoints[i] = 0;
3300 ir->opts.anneal_time[i] = NULL;
3301 ir->opts.anneal_temp[i] = NULL;
3306 for (i = 0; i < nr; i++)
3308 if (ptr1[i][0] == 'n' || ptr1[i][0] == 'N')
3310 ir->opts.annealing[i] = eannNO;
3312 else if (ptr1[i][0] == 's' || ptr1[i][0] == 'S')
3314 ir->opts.annealing[i] = eannSINGLE;
3317 else if (ptr1[i][0] == 'p' || ptr1[i][0] == 'P')
3319 ir->opts.annealing[i] = eannPERIODIC;
3325 /* Read the other fields too */
3326 nSA_points = str_nelem(is->anneal_npoints, MAXPTR, ptr1);
3327 if (nSA_points != nSA)
3329 gmx_fatal(FARGS, "Found %d annealing-npoints values for %d groups\n", nSA_points, nSA);
3331 for (k = 0, i = 0; i < nr; i++)
3333 ir->opts.anneal_npoints[i] = strtol(ptr1[i], NULL, 10);
3334 if (ir->opts.anneal_npoints[i] == 1)
3336 gmx_fatal(FARGS, "Please specify at least a start and an end point for annealing\n");
3338 snew(ir->opts.anneal_time[i], ir->opts.anneal_npoints[i]);
3339 snew(ir->opts.anneal_temp[i], ir->opts.anneal_npoints[i]);
3340 k += ir->opts.anneal_npoints[i];
3343 nSA_time = str_nelem(is->anneal_time, MAXPTR, ptr1);
3346 gmx_fatal(FARGS, "Found %d annealing-time values, wanter %d\n", nSA_time, k);
3348 nSA_temp = str_nelem(is->anneal_temp, MAXPTR, ptr2);
3351 gmx_fatal(FARGS, "Found %d annealing-temp values, wanted %d\n", nSA_temp, k);
3354 for (i = 0, k = 0; i < nr; i++)
3357 for (j = 0; j < ir->opts.anneal_npoints[i]; j++)
3359 ir->opts.anneal_time[i][j] = strtod(ptr1[k], NULL);
3360 ir->opts.anneal_temp[i][j] = strtod(ptr2[k], NULL);
3363 if (ir->opts.anneal_time[i][0] > (ir->init_t+GMX_REAL_EPS))
3365 gmx_fatal(FARGS, "First time point for annealing > init_t.\n");
3371 if (ir->opts.anneal_time[i][j] < ir->opts.anneal_time[i][j-1])
3373 gmx_fatal(FARGS, "Annealing timepoints out of order: t=%f comes after t=%f\n",
3374 ir->opts.anneal_time[i][j], ir->opts.anneal_time[i][j-1]);
3377 if (ir->opts.anneal_temp[i][j] < 0)
3379 gmx_fatal(FARGS, "Found negative temperature in annealing: %f\n", ir->opts.anneal_temp[i][j]);
3384 /* Print out some summary information, to make sure we got it right */
3385 for (i = 0, k = 0; i < nr; i++)
3387 if (ir->opts.annealing[i] != eannNO)
3389 j = groups->grps[egcTC].nm_ind[i];
3390 fprintf(stderr, "Simulated annealing for group %s: %s, %d timepoints\n",
3391 *(groups->grpname[j]), eann_names[ir->opts.annealing[i]],
3392 ir->opts.anneal_npoints[i]);
3393 fprintf(stderr, "Time (ps) Temperature (K)\n");
3394 /* All terms except the last one */
3395 for (j = 0; j < (ir->opts.anneal_npoints[i]-1); j++)
3397 fprintf(stderr, "%9.1f %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3400 /* Finally the last one */
3401 j = ir->opts.anneal_npoints[i]-1;
3402 if (ir->opts.annealing[i] == eannSINGLE)
3404 fprintf(stderr, "%9.1f- %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3408 fprintf(stderr, "%9.1f %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3409 if (fabs(ir->opts.anneal_temp[i][j]-ir->opts.anneal_temp[i][0]) > GMX_REAL_EPS)
3411 warning_note(wi, "There is a temperature jump when your annealing loops back.\n");
3420 if (ir->ePull != epullNO)
3422 make_pull_groups(ir->pull, is->pull_grp, grps, gnames);
3424 make_pull_coords(ir->pull);
3429 make_rotation_groups(ir->rot, is->rot_grp, grps, gnames);
3432 if (ir->eSwapCoords != eswapNO)
3434 make_swap_groups(ir->swap, swapgrp, splitgrp0, splitgrp1, solgrp, grps, gnames);
3437 /* Make indices for IMD session */
3440 make_IMD_group(ir->imd, is->imd_grp, grps, gnames);
3443 nacc = str_nelem(is->acc, MAXPTR, ptr1);
3444 nacg = str_nelem(is->accgrps, MAXPTR, ptr2);
3445 if (nacg*DIM != nacc)
3447 gmx_fatal(FARGS, "Invalid Acceleration input: %d groups and %d acc. values",
3450 do_numbering(natoms, groups, nacg, ptr2, grps, gnames, egcACC,
3451 restnm, egrptpALL_GENREST, bVerbose, wi);
3452 nr = groups->grps[egcACC].nr;
3453 snew(ir->opts.acc, nr);
3454 ir->opts.ngacc = nr;
3456 for (i = k = 0; (i < nacg); i++)
3458 for (j = 0; (j < DIM); j++, k++)
3460 ir->opts.acc[i][j] = strtod(ptr1[k], NULL);
3463 for (; (i < nr); i++)
3465 for (j = 0; (j < DIM); j++)
3467 ir->opts.acc[i][j] = 0;
3471 nfrdim = str_nelem(is->frdim, MAXPTR, ptr1);
3472 nfreeze = str_nelem(is->freeze, MAXPTR, ptr2);
3473 if (nfrdim != DIM*nfreeze)
3475 gmx_fatal(FARGS, "Invalid Freezing input: %d groups and %d freeze values",
3478 do_numbering(natoms, groups, nfreeze, ptr2, grps, gnames, egcFREEZE,
3479 restnm, egrptpALL_GENREST, bVerbose, wi);
3480 nr = groups->grps[egcFREEZE].nr;
3481 ir->opts.ngfrz = nr;
3482 snew(ir->opts.nFreeze, nr);
3483 for (i = k = 0; (i < nfreeze); i++)
3485 for (j = 0; (j < DIM); j++, k++)
3487 ir->opts.nFreeze[i][j] = (gmx_strncasecmp(ptr1[k], "Y", 1) == 0);
3488 if (!ir->opts.nFreeze[i][j])
3490 if (gmx_strncasecmp(ptr1[k], "N", 1) != 0)
3492 sprintf(warnbuf, "Please use Y(ES) or N(O) for freezedim only "
3493 "(not %s)", ptr1[k]);
3494 warning(wi, warn_buf);
3499 for (; (i < nr); i++)
3501 for (j = 0; (j < DIM); j++)
3503 ir->opts.nFreeze[i][j] = 0;
3507 nenergy = str_nelem(is->energy, MAXPTR, ptr1);
3508 do_numbering(natoms, groups, nenergy, ptr1, grps, gnames, egcENER,
3509 restnm, egrptpALL_GENREST, bVerbose, wi);
3510 add_wall_energrps(groups, ir->nwall, symtab);
3511 ir->opts.ngener = groups->grps[egcENER].nr;
3512 nvcm = str_nelem(is->vcm, MAXPTR, ptr1);
3514 do_numbering(natoms, groups, nvcm, ptr1, grps, gnames, egcVCM,
3515 restnm, nvcm == 0 ? egrptpALL_GENREST : egrptpPART, bVerbose, wi);
3518 warning(wi, "Some atoms are not part of any center of mass motion removal group.\n"
3519 "This may lead to artifacts.\n"
3520 "In most cases one should use one group for the whole system.");
3523 /* Now we have filled the freeze struct, so we can calculate NRDF */
3524 calc_nrdf(mtop, ir, gnames);
3530 /* Must check per group! */
3531 for (i = 0; (i < ir->opts.ngtc); i++)
3533 ntot += ir->opts.nrdf[i];
3535 if (ntot != (DIM*natoms))
3537 fac = sqrt(ntot/(DIM*natoms));
3540 fprintf(stderr, "Scaling velocities by a factor of %.3f to account for constraints\n"
3541 "and removal of center of mass motion\n", fac);
3543 for (i = 0; (i < natoms); i++)
3545 svmul(fac, v[i], v[i]);
3550 nuser = str_nelem(is->user1, MAXPTR, ptr1);
3551 do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcUser1,
3552 restnm, egrptpALL_GENREST, bVerbose, wi);
3553 nuser = str_nelem(is->user2, MAXPTR, ptr1);
3554 do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcUser2,
3555 restnm, egrptpALL_GENREST, bVerbose, wi);
3556 nuser = str_nelem(is->x_compressed_groups, MAXPTR, ptr1);
3557 do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcCompressedX,
3558 restnm, egrptpONE, bVerbose, wi);
3559 nofg = str_nelem(is->orirefitgrp, MAXPTR, ptr1);
3560 do_numbering(natoms, groups, nofg, ptr1, grps, gnames, egcORFIT,
3561 restnm, egrptpALL_GENREST, bVerbose, wi);
3563 /* QMMM input processing */
3564 nQMg = str_nelem(is->QMMM, MAXPTR, ptr1);
3565 nQMmethod = str_nelem(is->QMmethod, MAXPTR, ptr2);
3566 nQMbasis = str_nelem(is->QMbasis, MAXPTR, ptr3);
3567 if ((nQMmethod != nQMg) || (nQMbasis != nQMg))
3569 gmx_fatal(FARGS, "Invalid QMMM input: %d groups %d basissets"
3570 " and %d methods\n", nQMg, nQMbasis, nQMmethod);
3572 /* group rest, if any, is always MM! */
3573 do_numbering(natoms, groups, nQMg, ptr1, grps, gnames, egcQMMM,
3574 restnm, egrptpALL_GENREST, bVerbose, wi);
3575 nr = nQMg; /*atoms->grps[egcQMMM].nr;*/
3576 ir->opts.ngQM = nQMg;
3577 snew(ir->opts.QMmethod, nr);
3578 snew(ir->opts.QMbasis, nr);
3579 for (i = 0; i < nr; i++)
3581 /* input consists of strings: RHF CASSCF PM3 .. These need to be
3582 * converted to the corresponding enum in names.c
3584 ir->opts.QMmethod[i] = search_QMstring(ptr2[i], eQMmethodNR,
3586 ir->opts.QMbasis[i] = search_QMstring(ptr3[i], eQMbasisNR,
3590 nQMmult = str_nelem(is->QMmult, MAXPTR, ptr1);
3591 nQMcharge = str_nelem(is->QMcharge, MAXPTR, ptr2);
3592 nbSH = str_nelem(is->bSH, MAXPTR, ptr3);
3593 snew(ir->opts.QMmult, nr);
3594 snew(ir->opts.QMcharge, nr);
3595 snew(ir->opts.bSH, nr);
3597 for (i = 0; i < nr; i++)
3599 ir->opts.QMmult[i] = strtol(ptr1[i], NULL, 10);
3600 ir->opts.QMcharge[i] = strtol(ptr2[i], NULL, 10);
3601 ir->opts.bSH[i] = (gmx_strncasecmp(ptr3[i], "Y", 1) == 0);
3604 nCASelec = str_nelem(is->CASelectrons, MAXPTR, ptr1);
3605 nCASorb = str_nelem(is->CASorbitals, MAXPTR, ptr2);
3606 snew(ir->opts.CASelectrons, nr);
3607 snew(ir->opts.CASorbitals, nr);
3608 for (i = 0; i < nr; i++)
3610 ir->opts.CASelectrons[i] = strtol(ptr1[i], NULL, 10);
3611 ir->opts.CASorbitals[i] = strtol(ptr2[i], NULL, 10);
3613 /* special optimization options */
3615 nbOPT = str_nelem(is->bOPT, MAXPTR, ptr1);
3616 nbTS = str_nelem(is->bTS, MAXPTR, ptr2);
3617 snew(ir->opts.bOPT, nr);
3618 snew(ir->opts.bTS, nr);
3619 for (i = 0; i < nr; i++)
3621 ir->opts.bOPT[i] = (gmx_strncasecmp(ptr1[i], "Y", 1) == 0);
3622 ir->opts.bTS[i] = (gmx_strncasecmp(ptr2[i], "Y", 1) == 0);
3624 nSAon = str_nelem(is->SAon, MAXPTR, ptr1);
3625 nSAoff = str_nelem(is->SAoff, MAXPTR, ptr2);
3626 nSAsteps = str_nelem(is->SAsteps, MAXPTR, ptr3);
3627 snew(ir->opts.SAon, nr);
3628 snew(ir->opts.SAoff, nr);
3629 snew(ir->opts.SAsteps, nr);
3631 for (i = 0; i < nr; i++)
3633 ir->opts.SAon[i] = strtod(ptr1[i], NULL);
3634 ir->opts.SAoff[i] = strtod(ptr2[i], NULL);
3635 ir->opts.SAsteps[i] = strtol(ptr3[i], NULL, 10);
3637 /* end of QMMM input */
3641 for (i = 0; (i < egcNR); i++)
3643 fprintf(stderr, "%-16s has %d element(s):", gtypes[i], groups->grps[i].nr);
3644 for (j = 0; (j < groups->grps[i].nr); j++)
3646 fprintf(stderr, " %s", *(groups->grpname[groups->grps[i].nm_ind[j]]));
3648 fprintf(stderr, "\n");
3652 nr = groups->grps[egcENER].nr;
3653 snew(ir->opts.egp_flags, nr*nr);
3655 bExcl = do_egp_flag(ir, groups, "energygrp-excl", is->egpexcl, EGP_EXCL);
3656 if (bExcl && ir->cutoff_scheme == ecutsVERLET)
3658 warning_error(wi, "Energy group exclusions are not (yet) implemented for the Verlet scheme");
3660 if (bExcl && EEL_FULL(ir->coulombtype))
3662 warning(wi, "Can not exclude the lattice Coulomb energy between energy groups");
3665 bTable = do_egp_flag(ir, groups, "energygrp-table", is->egptable, EGP_TABLE);
3666 if (bTable && !(ir->vdwtype == evdwUSER) &&
3667 !(ir->coulombtype == eelUSER) && !(ir->coulombtype == eelPMEUSER) &&
3668 !(ir->coulombtype == eelPMEUSERSWITCH))
3670 gmx_fatal(FARGS, "Can only have energy group pair tables in combination with user tables for VdW and/or Coulomb");
3673 decode_cos(is->efield_x, &(ir->ex[XX]));
3674 decode_cos(is->efield_xt, &(ir->et[XX]));
3675 decode_cos(is->efield_y, &(ir->ex[YY]));
3676 decode_cos(is->efield_yt, &(ir->et[YY]));
3677 decode_cos(is->efield_z, &(ir->ex[ZZ]));
3678 decode_cos(is->efield_zt, &(ir->et[ZZ]));
3682 do_adress_index(ir->adress, groups, gnames, &(ir->opts), wi);
3685 for (i = 0; (i < grps->nr); i++)
3697 static void check_disre(gmx_mtop_t *mtop)
3699 gmx_ffparams_t *ffparams;
3700 t_functype *functype;
3702 int i, ndouble, ftype;
3703 int label, old_label;
3705 if (gmx_mtop_ftype_count(mtop, F_DISRES) > 0)
3707 ffparams = &mtop->ffparams;
3708 functype = ffparams->functype;
3709 ip = ffparams->iparams;
3712 for (i = 0; i < ffparams->ntypes; i++)
3714 ftype = functype[i];
3715 if (ftype == F_DISRES)
3717 label = ip[i].disres.label;
3718 if (label == old_label)
3720 fprintf(stderr, "Distance restraint index %d occurs twice\n", label);
3728 gmx_fatal(FARGS, "Found %d double distance restraint indices,\n"
3729 "probably the parameters for multiple pairs in one restraint "
3730 "are not identical\n", ndouble);
3735 static gmx_bool absolute_reference(t_inputrec *ir, gmx_mtop_t *sys,
3736 gmx_bool posres_only,
3740 gmx_mtop_ilistloop_t iloop;
3750 for (d = 0; d < DIM; d++)
3752 AbsRef[d] = (d < ndof_com(ir) ? 0 : 1);
3754 /* Check for freeze groups */
3755 for (g = 0; g < ir->opts.ngfrz; g++)
3757 for (d = 0; d < DIM; d++)
3759 if (ir->opts.nFreeze[g][d] != 0)
3767 /* Check for position restraints */
3768 iloop = gmx_mtop_ilistloop_init(sys);
3769 while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
3772 (AbsRef[XX] == 0 || AbsRef[YY] == 0 || AbsRef[ZZ] == 0))
3774 for (i = 0; i < ilist[F_POSRES].nr; i += 2)
3776 pr = &sys->ffparams.iparams[ilist[F_POSRES].iatoms[i]];
3777 for (d = 0; d < DIM; d++)
3779 if (pr->posres.fcA[d] != 0)
3785 for (i = 0; i < ilist[F_FBPOSRES].nr; i += 2)
3787 /* Check for flat-bottom posres */
3788 pr = &sys->ffparams.iparams[ilist[F_FBPOSRES].iatoms[i]];
3789 if (pr->fbposres.k != 0)
3791 switch (pr->fbposres.geom)
3793 case efbposresSPHERE:
3794 AbsRef[XX] = AbsRef[YY] = AbsRef[ZZ] = 1;
3796 case efbposresCYLINDER:
3797 AbsRef[XX] = AbsRef[YY] = 1;
3799 case efbposresX: /* d=XX */
3800 case efbposresY: /* d=YY */
3801 case efbposresZ: /* d=ZZ */
3802 d = pr->fbposres.geom - efbposresX;
3806 gmx_fatal(FARGS, " Invalid geometry for flat-bottom position restraint.\n"
3807 "Expected nr between 1 and %d. Found %d\n", efbposresNR-1,
3815 return (AbsRef[XX] != 0 && AbsRef[YY] != 0 && AbsRef[ZZ] != 0);
3819 check_combination_rule_differences(const gmx_mtop_t *mtop, int state,
3820 gmx_bool *bC6ParametersWorkWithGeometricRules,
3821 gmx_bool *bC6ParametersWorkWithLBRules,
3822 gmx_bool *bLBRulesPossible)
3824 int ntypes, tpi, tpj, thisLBdiff, thisgeomdiff;
3827 double geometricdiff, LBdiff;
3828 double c6i, c6j, c12i, c12j;
3829 double c6, c6_geometric, c6_LB;
3830 double sigmai, sigmaj, epsi, epsj;
3831 gmx_bool bCanDoLBRules, bCanDoGeometricRules;
3834 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
3835 * force-field floating point parameters.
3838 ptr = getenv("GMX_LJCOMB_TOL");
3843 sscanf(ptr, "%lf", &dbl);
3847 *bC6ParametersWorkWithLBRules = TRUE;
3848 *bC6ParametersWorkWithGeometricRules = TRUE;
3849 bCanDoLBRules = TRUE;
3850 bCanDoGeometricRules = TRUE;
3851 ntypes = mtop->ffparams.atnr;
3852 snew(typecount, ntypes);
3853 gmx_mtop_count_atomtypes(mtop, state, typecount);
3854 geometricdiff = LBdiff = 0.0;
3855 *bLBRulesPossible = TRUE;
3856 for (tpi = 0; tpi < ntypes; ++tpi)
3858 c6i = mtop->ffparams.iparams[(ntypes + 1) * tpi].lj.c6;
3859 c12i = mtop->ffparams.iparams[(ntypes + 1) * tpi].lj.c12;
3860 for (tpj = tpi; tpj < ntypes; ++tpj)
3862 c6j = mtop->ffparams.iparams[(ntypes + 1) * tpj].lj.c6;
3863 c12j = mtop->ffparams.iparams[(ntypes + 1) * tpj].lj.c12;
3864 c6 = mtop->ffparams.iparams[ntypes * tpi + tpj].lj.c6;
3865 c6_geometric = sqrt(c6i * c6j);
3866 if (!gmx_numzero(c6_geometric))
3868 if (!gmx_numzero(c12i) && !gmx_numzero(c12j))
3870 sigmai = pow(c12i / c6i, 1.0/6.0);
3871 sigmaj = pow(c12j / c6j, 1.0/6.0);
3872 epsi = c6i * c6i /(4.0 * c12i);
3873 epsj = c6j * c6j /(4.0 * c12j);
3874 c6_LB = 4.0 * pow(epsi * epsj, 1.0/2.0) * pow(0.5 * (sigmai + sigmaj), 6);
3878 *bLBRulesPossible = FALSE;
3879 c6_LB = c6_geometric;
3881 bCanDoLBRules = gmx_within_tol(c6_LB, c6, tol);
3884 if (FALSE == bCanDoLBRules)
3886 *bC6ParametersWorkWithLBRules = FALSE;
3889 bCanDoGeometricRules = gmx_within_tol(c6_geometric, c6, tol);
3891 if (FALSE == bCanDoGeometricRules)
3893 *bC6ParametersWorkWithGeometricRules = FALSE;
3901 check_combination_rules(const t_inputrec *ir, const gmx_mtop_t *mtop,
3905 gmx_bool bLBRulesPossible, bC6ParametersWorkWithGeometricRules, bC6ParametersWorkWithLBRules;
3907 check_combination_rule_differences(mtop, 0,
3908 &bC6ParametersWorkWithGeometricRules,
3909 &bC6ParametersWorkWithLBRules,
3911 if (ir->ljpme_combination_rule == eljpmeLB)
3913 if (FALSE == bC6ParametersWorkWithLBRules || FALSE == 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.");
3922 if (FALSE == 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.");
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,
3949 char err_buf[STRLEN];
3950 int i, m, c, nmol, npct;
3951 gmx_bool bCharge, bAcc;
3952 real gdt_max, *mgrp, mt;
3954 gmx_mtop_atomloop_block_t aloopb;
3955 gmx_mtop_atomloop_all_t aloop;
3958 char warn_buf[STRLEN];
3960 set_warning_line(wi, mdparin, -1);
3962 if (ir->cutoff_scheme == ecutsVERLET &&
3963 ir->verletbuf_tol > 0 &&
3965 ((EI_MD(ir->eI) || EI_SD(ir->eI)) &&
3966 (ir->etc == etcVRESCALE || ir->etc == etcBERENDSEN)))
3968 /* Check if a too small Verlet buffer might potentially
3969 * cause more drift than the thermostat can couple off.
3971 /* Temperature error fraction for warning and suggestion */
3972 const real T_error_warn = 0.002;
3973 const real T_error_suggest = 0.001;
3974 /* For safety: 2 DOF per atom (typical with constraints) */
3975 const real nrdf_at = 2;
3976 real T, tau, max_T_error;
3981 for (i = 0; i < ir->opts.ngtc; i++)
3983 T = max(T, ir->opts.ref_t[i]);
3984 tau = max(tau, ir->opts.tau_t[i]);
3988 /* This is a worst case estimate of the temperature error,
3989 * assuming perfect buffer estimation and no cancelation
3990 * of errors. The factor 0.5 is because energy distributes
3991 * equally over Ekin and Epot.
3993 max_T_error = 0.5*tau*ir->verletbuf_tol/(nrdf_at*BOLTZ*T);
3994 if (max_T_error > T_error_warn)
3996 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.",
3997 ir->verletbuf_tol, T, tau,
3999 100*T_error_suggest,
4000 ir->verletbuf_tol*T_error_suggest/max_T_error);
4001 warning(wi, warn_buf);
4006 if (ETC_ANDERSEN(ir->etc))
4010 for (i = 0; i < ir->opts.ngtc; i++)
4012 sprintf(err_buf, "all tau_t must currently be equal using Andersen temperature control, violated for group %d", i);
4013 CHECK(ir->opts.tau_t[0] != ir->opts.tau_t[i]);
4014 sprintf(err_buf, "all tau_t must be postive using Andersen temperature control, tau_t[%d]=%10.6f",
4015 i, ir->opts.tau_t[i]);
4016 CHECK(ir->opts.tau_t[i] < 0);
4019 for (i = 0; i < ir->opts.ngtc; i++)
4021 int nsteps = (int)(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) && (ir->etc == etcANDERSENMASSIVE));
4027 if (EI_DYNAMICS(ir->eI) && !EI_SD(ir->eI) && ir->eI != eiBD &&
4028 ir->comm_mode == ecmNO &&
4029 !(absolute_reference(ir, sys, FALSE, AbsRef) || ir->nsteps <= 10) &&
4030 !ETC_ANDERSEN(ir->etc))
4032 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");
4035 /* Check for pressure coupling with absolute position restraints */
4036 if (ir->epc != epcNO && ir->refcoord_scaling == erscNO)
4038 absolute_reference(ir, sys, TRUE, AbsRef);
4040 for (m = 0; m < DIM; m++)
4042 if (AbsRef[m] && norm2(ir->compress[m]) > 0)
4044 warning(wi, "You are using pressure coupling with absolute position restraints, this will give artifacts. Use the refcoord_scaling option.");
4052 aloopb = gmx_mtop_atomloop_block_init(sys);
4053 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
4055 if (atom->q != 0 || atom->qB != 0)
4063 if (EEL_FULL(ir->coulombtype))
4066 "You are using full electrostatics treatment %s for a system without charges.\n"
4067 "This costs a lot of performance for just processing zeros, consider using %s instead.\n",
4068 EELTYPE(ir->coulombtype), EELTYPE(eelCUT));
4069 warning(wi, err_buf);
4074 if (ir->coulombtype == eelCUT && ir->rcoulomb > 0 && !ir->implicit_solvent)
4077 "You are using a plain Coulomb cut-off, which might produce artifacts.\n"
4078 "You might want to consider using %s electrostatics.\n",
4080 warning_note(wi, err_buf);
4084 /* Check if combination rules used in LJ-PME are the same as in the force field */
4085 if (EVDW_PME(ir->vdwtype))
4087 check_combination_rules(ir, sys, wi);
4090 /* Generalized reaction field */
4091 if (ir->opts.ngtc == 0)
4093 sprintf(err_buf, "No temperature coupling while using coulombtype %s",
4095 CHECK(ir->coulombtype == eelGRF);
4099 sprintf(err_buf, "When using coulombtype = %s"
4100 " ref-t for temperature coupling should be > 0",
4102 CHECK((ir->coulombtype == eelGRF) && (ir->opts.ref_t[0] <= 0));
4105 if (ir->eI == eiSD1 &&
4106 (gmx_mtop_ftype_count(sys, F_CONSTR) > 0 ||
4107 gmx_mtop_ftype_count(sys, F_SETTLE) > 0))
4109 sprintf(warn_buf, "With constraints integrator %s is less accurate, consider using %s instead", ei_names[ir->eI], ei_names[eiSD2]);
4110 warning_note(wi, warn_buf);
4114 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4116 for (m = 0; (m < DIM); m++)
4118 if (fabs(ir->opts.acc[i][m]) > 1e-6)
4127 snew(mgrp, sys->groups.grps[egcACC].nr);
4128 aloop = gmx_mtop_atomloop_all_init(sys);
4129 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
4131 mgrp[ggrpnr(&sys->groups, egcACC, i)] += atom->m;
4134 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4136 for (m = 0; (m < DIM); m++)
4138 acc[m] += ir->opts.acc[i][m]*mgrp[i];
4142 for (m = 0; (m < DIM); m++)
4144 if (fabs(acc[m]) > 1e-6)
4146 const char *dim[DIM] = { "X", "Y", "Z" };
4148 "Net Acceleration in %s direction, will %s be corrected\n",
4149 dim[m], ir->nstcomm != 0 ? "" : "not");
4150 if (ir->nstcomm != 0 && m < ndof_com(ir))
4153 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4155 ir->opts.acc[i][m] -= acc[m];
4163 if (ir->efep != efepNO && ir->fepvals->sc_alpha != 0 &&
4164 !gmx_within_tol(sys->ffparams.reppow, 12.0, 10*GMX_DOUBLE_EPS))
4166 gmx_fatal(FARGS, "Soft-core interactions are only supported with VdW repulsion power 12");
4169 if (ir->ePull != epullNO)
4171 gmx_bool bPullAbsoluteRef;
4173 bPullAbsoluteRef = FALSE;
4174 for (i = 0; i < ir->pull->ncoord; i++)
4176 bPullAbsoluteRef = bPullAbsoluteRef ||
4177 ir->pull->coord[i].group[0] == 0 ||
4178 ir->pull->coord[i].group[1] == 0;
4180 if (bPullAbsoluteRef)
4182 absolute_reference(ir, sys, FALSE, AbsRef);
4183 for (m = 0; m < DIM; m++)
4185 if (ir->pull->dim[m] && !AbsRef[m])
4187 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.");
4193 if (ir->pull->eGeom == epullgDIRPBC)
4195 for (i = 0; i < 3; i++)
4197 for (m = 0; m <= i; m++)
4199 if ((ir->epc != epcNO && ir->compress[i][m] != 0) ||
4200 ir->deform[i][m] != 0)
4202 for (c = 0; c < ir->pull->ncoord; c++)
4204 if (ir->pull->coord[c].vec[m] != 0)
4206 gmx_fatal(FARGS, "Can not have dynamic box while using pull geometry '%s' (dim %c)", EPULLGEOM(ir->pull->eGeom), 'x'+m);
4218 void double_check(t_inputrec *ir, matrix box, gmx_bool bConstr, warninp_t wi)
4222 char warn_buf[STRLEN];
4225 ptr = check_box(ir->ePBC, box);
4228 warning_error(wi, ptr);
4231 if (bConstr && ir->eConstrAlg == econtSHAKE)
4233 if (ir->shake_tol <= 0.0)
4235 sprintf(warn_buf, "ERROR: shake-tol must be > 0 instead of %g\n",
4237 warning_error(wi, warn_buf);
4240 if (IR_TWINRANGE(*ir) && ir->nstlist > 1)
4242 sprintf(warn_buf, "With twin-range cut-off's and SHAKE the virial and the pressure are incorrect.");
4243 if (ir->epc == epcNO)
4245 warning(wi, warn_buf);
4249 warning_error(wi, warn_buf);
4254 if ( (ir->eConstrAlg == econtLINCS) && bConstr)
4256 /* If we have Lincs constraints: */
4257 if (ir->eI == eiMD && ir->etc == etcNO &&
4258 ir->eConstrAlg == econtLINCS && ir->nLincsIter == 1)
4260 sprintf(warn_buf, "For energy conservation with LINCS, lincs_iter should be 2 or larger.\n");
4261 warning_note(wi, warn_buf);
4264 if ((ir->eI == eiCG || ir->eI == eiLBFGS) && (ir->nProjOrder < 8))
4266 sprintf(warn_buf, "For accurate %s with LINCS constraints, lincs-order should be 8 or more.", ei_names[ir->eI]);
4267 warning_note(wi, warn_buf);
4269 if (ir->epc == epcMTTK)
4271 warning_error(wi, "MTTK not compatible with lincs -- use shake instead.");
4275 if (bConstr && ir->epc == epcMTTK)
4277 warning_note(wi, "MTTK with constraints is deprecated, and will be removed in GROMACS 5.1");
4280 if (ir->LincsWarnAngle > 90.0)
4282 sprintf(warn_buf, "lincs-warnangle can not be larger than 90 degrees, setting it to 90.\n");
4283 warning(wi, warn_buf);
4284 ir->LincsWarnAngle = 90.0;
4287 if (ir->ePBC != epbcNONE)
4289 if (ir->nstlist == 0)
4291 warning(wi, "With nstlist=0 atoms are only put into the box at step 0, therefore drifting atoms might cause the simulation to crash.");
4293 bTWIN = (ir->rlistlong > ir->rlist);
4294 if (ir->ns_type == ensGRID)
4296 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
4298 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 %s.\n",
4299 bTWIN ? (ir->rcoulomb == ir->rlistlong ? "rcoulomb" : "rvdw") : "rlist");
4300 warning_error(wi, warn_buf);
4305 min_size = min(box[XX][XX], min(box[YY][YY], box[ZZ][ZZ]));
4306 if (2*ir->rlistlong >= 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);
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,
4323 real rvdw1, rvdw2, rcoul1, rcoul2;
4324 char warn_buf[STRLEN];
4326 calc_chargegroup_radii(mtop, x, &rvdw1, &rvdw2, &rcoul1, &rcoul2);
4330 printf("Largest charge group radii for Van der Waals: %5.3f, %5.3f nm\n",
4335 printf("Largest charge group radii for Coulomb: %5.3f, %5.3f nm\n",
4341 if (rvdw1 + rvdw2 > ir->rlist ||
4342 rcoul1 + rcoul2 > ir->rlist)
4345 "The sum of the two largest charge group radii (%f) "
4346 "is larger than rlist (%f)\n",
4347 max(rvdw1+rvdw2, rcoul1+rcoul2), ir->rlist);
4348 warning(wi, warn_buf);
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->rlistlong - ir->rvdw)
4359 sprintf(warn_buf, "The sum of the two largest charge group "
4360 "radii (%f) is larger than %s (%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.",
4365 ir->rlistlong > ir->rlist ? "rlistlong" : "rlist",
4366 ir->rlistlong, ir->rvdw,
4367 ecutscheme_names[ecutsVERLET]);
4370 warning(wi, warn_buf);
4374 warning_note(wi, warn_buf);
4377 if (ir_coulomb_is_zero_at_cutoff(ir) &&
4378 rcoul1 + rcoul2 > ir->rlistlong - ir->rcoulomb)
4380 sprintf(warn_buf, "The sum of the two largest charge group radii (%f) is larger than %s (%f) - rcoulomb (%f).\n"
4381 "With exact cut-offs, better performance can be obtained with cutoff-scheme = %s, because it does not use charge groups at all.",
4383 ir->rlistlong > ir->rlist ? "rlistlong" : "rlist",
4384 ir->rlistlong, ir->rcoulomb,
4385 ecutscheme_names[ecutsVERLET]);
4388 warning(wi, warn_buf);
4392 warning_note(wi, warn_buf);