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.
45 #include "gromacs/legacyheaders/typedefs.h"
46 #include "gromacs/math/units.h"
47 #include "gromacs/legacyheaders/names.h"
48 #include "gromacs/legacyheaders/macros.h"
49 #include "gromacs/topology/index.h"
50 #include "gromacs/utility/cstringutil.h"
51 #include "gromacs/legacyheaders/readinp.h"
52 #include "gromacs/legacyheaders/warninp.h"
55 #include "gromacs/legacyheaders/network.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/pbcutil/pbc.h"
58 #include "gromacs/topology/mtop_util.h"
59 #include "gromacs/legacyheaders/chargegroup.h"
60 #include "gromacs/legacyheaders/inputrec.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);
556 if (ir->nsteps == 0 && !ir->bContinuation)
558 warning_note(wi, "For a correct single-point energy evaluation with nsteps = 0, use continuation = yes to avoid constraining the input coordinates.");
562 if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
563 ir->bContinuation && ir->ld_seed != -1)
565 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)");
571 sprintf(err_buf, "TPI only works with pbc = %s", epbc_names[epbcXYZ]);
572 CHECK(ir->ePBC != epbcXYZ);
573 sprintf(err_buf, "TPI only works with ns = %s", ens_names[ensGRID]);
574 CHECK(ir->ns_type != ensGRID);
575 sprintf(err_buf, "with TPI nstlist should be larger than zero");
576 CHECK(ir->nstlist <= 0);
577 sprintf(err_buf, "TPI does not work with full electrostatics other than PME");
578 CHECK(EEL_FULL(ir->coulombtype) && !EEL_PME(ir->coulombtype));
579 sprintf(err_buf, "TPI does not work (yet) with the Verlet cut-off scheme");
580 CHECK(ir->cutoff_scheme == ecutsVERLET);
584 if ( (opts->nshake > 0) && (opts->bMorse) )
587 "Using morse bond-potentials while constraining bonds is useless");
588 warning(wi, warn_buf);
591 if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
592 ir->bContinuation && ir->ld_seed != -1)
594 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)");
596 /* verify simulated tempering options */
600 gmx_bool bAllTempZero = TRUE;
601 for (i = 0; i < fep->n_lambda; i++)
603 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]);
604 CHECK((fep->all_lambda[efptTEMPERATURE][i] < 0) || (fep->all_lambda[efptTEMPERATURE][i] > 1));
605 if (fep->all_lambda[efptTEMPERATURE][i] > 0)
607 bAllTempZero = FALSE;
610 sprintf(err_buf, "if simulated tempering is on, temperature-lambdas may not be all zero");
611 CHECK(bAllTempZero == TRUE);
613 sprintf(err_buf, "Simulated tempering is currently only compatible with md-vv");
614 CHECK(ir->eI != eiVV);
616 /* check compatability of the temperature coupling with simulated tempering */
618 if (ir->etc == etcNOSEHOOVER)
620 sprintf(warn_buf, "Nose-Hoover based temperature control such as [%s] my not be entirelyconsistent with simulated tempering", etcoupl_names[ir->etc]);
621 warning_note(wi, warn_buf);
624 /* check that the temperatures make sense */
626 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);
627 CHECK(ir->simtempvals->simtemp_high <= ir->simtempvals->simtemp_low);
629 sprintf(err_buf, "Higher simulated tempering temperature (%g) must be >= zero", ir->simtempvals->simtemp_high);
630 CHECK(ir->simtempvals->simtemp_high <= 0);
632 sprintf(err_buf, "Lower simulated tempering temperature (%g) must be >= zero", ir->simtempvals->simtemp_low);
633 CHECK(ir->simtempvals->simtemp_low <= 0);
636 /* verify free energy options */
638 if (ir->efep != efepNO)
641 sprintf(err_buf, "The soft-core power is %d and can only be 1 or 2",
643 CHECK(fep->sc_alpha != 0 && fep->sc_power != 1 && fep->sc_power != 2);
645 sprintf(err_buf, "The soft-core sc-r-power is %d and can only be 6 or 48",
646 (int)fep->sc_r_power);
647 CHECK(fep->sc_alpha != 0 && fep->sc_r_power != 6.0 && fep->sc_r_power != 48.0);
649 sprintf(err_buf, "Can't use postive delta-lambda (%g) if initial state/lambda does not start at zero", fep->delta_lambda);
650 CHECK(fep->delta_lambda > 0 && ((fep->init_fep_state > 0) || (fep->init_lambda > 0)));
652 sprintf(err_buf, "Can't use postive delta-lambda (%g) with expanded ensemble simulations", fep->delta_lambda);
653 CHECK(fep->delta_lambda > 0 && (ir->efep == efepEXPANDED));
655 sprintf(err_buf, "Can only use expanded ensemble with md-vv for now; should be supported for other integrators in 5.0");
656 CHECK(!(EI_VV(ir->eI)) && (ir->efep == efepEXPANDED));
658 sprintf(err_buf, "Free-energy not implemented for Ewald");
659 CHECK(ir->coulombtype == eelEWALD);
661 /* check validty of lambda inputs */
662 if (fep->n_lambda == 0)
664 /* Clear output in case of no states:*/
665 sprintf(err_buf, "init-lambda-state set to %d: no lambda states are defined.", fep->init_fep_state);
666 CHECK((fep->init_fep_state >= 0) && (fep->n_lambda == 0));
670 sprintf(err_buf, "initial thermodynamic state %d does not exist, only goes to %d", fep->init_fep_state, fep->n_lambda-1);
671 CHECK((fep->init_fep_state >= fep->n_lambda));
674 sprintf(err_buf, "Lambda state must be set, either with init-lambda-state or with init-lambda");
675 CHECK((fep->init_fep_state < 0) && (fep->init_lambda < 0));
677 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",
678 fep->init_lambda, fep->init_fep_state);
679 CHECK((fep->init_fep_state >= 0) && (fep->init_lambda >= 0));
683 if ((fep->init_lambda >= 0) && (fep->delta_lambda == 0))
687 for (i = 0; i < efptNR; i++)
689 if (fep->separate_dvdl[i])
694 if (n_lambda_terms > 1)
696 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.");
697 warning(wi, warn_buf);
700 if (n_lambda_terms < 2 && fep->n_lambda > 0)
703 "init-lambda is deprecated for setting lambda state (except for slow growth). Use init-lambda-state instead.");
707 for (j = 0; j < efptNR; j++)
709 for (i = 0; i < fep->n_lambda; i++)
711 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]);
712 CHECK((fep->all_lambda[j][i] < 0) || (fep->all_lambda[j][i] > 1));
716 if ((fep->sc_alpha > 0) && (!fep->bScCoul))
718 for (i = 0; i < fep->n_lambda; i++)
720 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],
721 fep->all_lambda[efptCOUL][i]);
722 CHECK((fep->sc_alpha > 0) &&
723 (((fep->all_lambda[efptCOUL][i] > 0.0) &&
724 (fep->all_lambda[efptCOUL][i] < 1.0)) &&
725 ((fep->all_lambda[efptVDW][i] > 0.0) &&
726 (fep->all_lambda[efptVDW][i] < 1.0))));
730 if ((fep->bScCoul) && (EEL_PME(ir->coulombtype)))
732 real sigma, lambda, r_sc;
735 /* Maximum estimate for A and B charges equal with lambda power 1 */
737 r_sc = pow(lambda*fep->sc_alpha*pow(sigma/ir->rcoulomb, fep->sc_r_power) + 1.0, 1.0/fep->sc_r_power);
738 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.",
740 sigma, lambda, r_sc - 1.0, ir->ewald_rtol);
741 warning_note(wi, warn_buf);
744 /* Free Energy Checks -- In an ideal world, slow growth and FEP would
745 be treated differently, but that's the next step */
747 for (i = 0; i < efptNR; i++)
749 for (j = 0; j < fep->n_lambda; j++)
751 sprintf(err_buf, "%s[%d] must be between 0 and 1", efpt_names[i], j);
752 CHECK((fep->all_lambda[i][j] < 0) || (fep->all_lambda[i][j] > 1));
757 if ((ir->bSimTemp) || (ir->efep == efepEXPANDED))
760 expand = ir->expandedvals;
762 /* checking equilibration of weights inputs for validity */
764 sprintf(err_buf, "weight-equil-number-all-lambda (%d) is ignored if lmc-weights-equil is not equal to %s",
765 expand->equil_n_at_lam, elmceq_names[elmceqNUMATLAM]);
766 CHECK((expand->equil_n_at_lam > 0) && (expand->elmceq != elmceqNUMATLAM));
768 sprintf(err_buf, "weight-equil-number-samples (%d) is ignored if lmc-weights-equil is not equal to %s",
769 expand->equil_samples, elmceq_names[elmceqSAMPLES]);
770 CHECK((expand->equil_samples > 0) && (expand->elmceq != elmceqSAMPLES));
772 sprintf(err_buf, "weight-equil-number-steps (%d) is ignored if lmc-weights-equil is not equal to %s",
773 expand->equil_steps, elmceq_names[elmceqSTEPS]);
774 CHECK((expand->equil_steps > 0) && (expand->elmceq != elmceqSTEPS));
776 sprintf(err_buf, "weight-equil-wl-delta (%d) is ignored if lmc-weights-equil is not equal to %s",
777 expand->equil_samples, elmceq_names[elmceqWLDELTA]);
778 CHECK((expand->equil_wl_delta > 0) && (expand->elmceq != elmceqWLDELTA));
780 sprintf(err_buf, "weight-equil-count-ratio (%f) is ignored if lmc-weights-equil is not equal to %s",
781 expand->equil_ratio, elmceq_names[elmceqRATIO]);
782 CHECK((expand->equil_ratio > 0) && (expand->elmceq != elmceqRATIO));
784 sprintf(err_buf, "weight-equil-number-all-lambda (%d) must be a positive integer if lmc-weights-equil=%s",
785 expand->equil_n_at_lam, elmceq_names[elmceqNUMATLAM]);
786 CHECK((expand->equil_n_at_lam <= 0) && (expand->elmceq == elmceqNUMATLAM));
788 sprintf(err_buf, "weight-equil-number-samples (%d) must be a positive integer if lmc-weights-equil=%s",
789 expand->equil_samples, elmceq_names[elmceqSAMPLES]);
790 CHECK((expand->equil_samples <= 0) && (expand->elmceq == elmceqSAMPLES));
792 sprintf(err_buf, "weight-equil-number-steps (%d) must be a positive integer if lmc-weights-equil=%s",
793 expand->equil_steps, elmceq_names[elmceqSTEPS]);
794 CHECK((expand->equil_steps <= 0) && (expand->elmceq == elmceqSTEPS));
796 sprintf(err_buf, "weight-equil-wl-delta (%f) must be > 0 if lmc-weights-equil=%s",
797 expand->equil_wl_delta, elmceq_names[elmceqWLDELTA]);
798 CHECK((expand->equil_wl_delta <= 0) && (expand->elmceq == elmceqWLDELTA));
800 sprintf(err_buf, "weight-equil-count-ratio (%f) must be > 0 if lmc-weights-equil=%s",
801 expand->equil_ratio, elmceq_names[elmceqRATIO]);
802 CHECK((expand->equil_ratio <= 0) && (expand->elmceq == elmceqRATIO));
804 sprintf(err_buf, "lmc-weights-equil=%s only possible when lmc-stats = %s or lmc-stats %s",
805 elmceq_names[elmceqWLDELTA], elamstats_names[elamstatsWL], elamstats_names[elamstatsWWL]);
806 CHECK((expand->elmceq == elmceqWLDELTA) && (!EWL(expand->elamstats)));
808 sprintf(err_buf, "lmc-repeats (%d) must be greater than 0", expand->lmc_repeats);
809 CHECK((expand->lmc_repeats <= 0));
810 sprintf(err_buf, "minimum-var-min (%d) must be greater than 0", expand->minvarmin);
811 CHECK((expand->minvarmin <= 0));
812 sprintf(err_buf, "weight-c-range (%d) must be greater or equal to 0", expand->c_range);
813 CHECK((expand->c_range < 0));
814 sprintf(err_buf, "init-lambda-state (%d) must be zero if lmc-forced-nstart (%d)> 0 and lmc-move != 'no'",
815 fep->init_fep_state, expand->lmc_forced_nstart);
816 CHECK((fep->init_fep_state != 0) && (expand->lmc_forced_nstart > 0) && (expand->elmcmove != elmcmoveNO));
817 sprintf(err_buf, "lmc-forced-nstart (%d) must not be negative", expand->lmc_forced_nstart);
818 CHECK((expand->lmc_forced_nstart < 0));
819 sprintf(err_buf, "init-lambda-state (%d) must be in the interval [0,number of lambdas)", fep->init_fep_state);
820 CHECK((fep->init_fep_state < 0) || (fep->init_fep_state >= fep->n_lambda));
822 sprintf(err_buf, "init-wl-delta (%f) must be greater than or equal to 0", expand->init_wl_delta);
823 CHECK((expand->init_wl_delta < 0));
824 sprintf(err_buf, "wl-ratio (%f) must be between 0 and 1", expand->wl_ratio);
825 CHECK((expand->wl_ratio <= 0) || (expand->wl_ratio >= 1));
826 sprintf(err_buf, "wl-scale (%f) must be between 0 and 1", expand->wl_scale);
827 CHECK((expand->wl_scale <= 0) || (expand->wl_scale >= 1));
829 /* if there is no temperature control, we need to specify an MC temperature */
830 sprintf(err_buf, "If there is no temperature control, and lmc-mcmove!= 'no',mc_temperature must be set to a positive number");
831 if (expand->nstTij > 0)
833 sprintf(err_buf, "nst-transition-matrix (%d) must be an integer multiple of nstlog (%d)",
834 expand->nstTij, ir->nstlog);
835 CHECK((mod(expand->nstTij, ir->nstlog) != 0));
840 sprintf(err_buf, "walls only work with pbc=%s", epbc_names[epbcXY]);
841 CHECK(ir->nwall && ir->ePBC != epbcXY);
844 if (ir->ePBC != epbcXYZ && ir->nwall != 2)
846 if (ir->ePBC == epbcNONE)
848 if (ir->epc != epcNO)
850 warning(wi, "Turning off pressure coupling for vacuum system");
856 sprintf(err_buf, "Can not have pressure coupling with pbc=%s",
857 epbc_names[ir->ePBC]);
858 CHECK(ir->epc != epcNO);
860 sprintf(err_buf, "Can not have Ewald with pbc=%s", epbc_names[ir->ePBC]);
861 CHECK(EEL_FULL(ir->coulombtype));
863 sprintf(err_buf, "Can not have dispersion correction with pbc=%s",
864 epbc_names[ir->ePBC]);
865 CHECK(ir->eDispCorr != edispcNO);
868 if (ir->rlist == 0.0)
870 sprintf(err_buf, "can only have neighborlist cut-off zero (=infinite)\n"
871 "with coulombtype = %s or coulombtype = %s\n"
872 "without periodic boundary conditions (pbc = %s) and\n"
873 "rcoulomb and rvdw set to zero",
874 eel_names[eelCUT], eel_names[eelUSER], epbc_names[epbcNONE]);
875 CHECK(((ir->coulombtype != eelCUT) && (ir->coulombtype != eelUSER)) ||
876 (ir->ePBC != epbcNONE) ||
877 (ir->rcoulomb != 0.0) || (ir->rvdw != 0.0));
881 warning_error(wi, "Can not have heuristic neighborlist updates without cut-off");
885 warning_note(wi, "Simulating without cut-offs can be (slightly) faster with nstlist=0, nstype=simple and only one MPI rank");
890 if (ir->nstcomm == 0)
892 ir->comm_mode = ecmNO;
894 if (ir->comm_mode != ecmNO)
898 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");
899 ir->nstcomm = abs(ir->nstcomm);
902 if (ir->nstcalcenergy > 0 && ir->nstcomm < ir->nstcalcenergy)
904 warning_note(wi, "nstcomm < nstcalcenergy defeats the purpose of nstcalcenergy, setting nstcomm to nstcalcenergy");
905 ir->nstcomm = ir->nstcalcenergy;
908 if (ir->comm_mode == ecmANGULAR)
910 sprintf(err_buf, "Can not remove the rotation around the center of mass with periodic molecules");
911 CHECK(ir->bPeriodicMols);
912 if (ir->ePBC != epbcNONE)
914 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).");
919 if (EI_STATE_VELOCITY(ir->eI) && ir->ePBC == epbcNONE && ir->comm_mode != ecmANGULAR)
921 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.");
924 sprintf(err_buf, "Twin-range neighbour searching (NS) with simple NS"
925 " algorithm not implemented");
926 CHECK(((ir->rcoulomb > ir->rlist) || (ir->rvdw > ir->rlist))
927 && (ir->ns_type == ensSIMPLE));
929 /* TEMPERATURE COUPLING */
930 if (ir->etc == etcYES)
932 ir->etc = etcBERENDSEN;
933 warning_note(wi, "Old option for temperature coupling given: "
934 "changing \"yes\" to \"Berendsen\"\n");
937 if ((ir->etc == etcNOSEHOOVER) || (ir->epc == epcMTTK))
939 if (ir->opts.nhchainlength < 1)
941 sprintf(warn_buf, "number of Nose-Hoover chains (currently %d) cannot be less than 1,reset to 1\n", ir->opts.nhchainlength);
942 ir->opts.nhchainlength = 1;
943 warning(wi, warn_buf);
946 if (ir->etc == etcNOSEHOOVER && !EI_VV(ir->eI) && ir->opts.nhchainlength > 1)
948 warning_note(wi, "leapfrog does not yet support Nose-Hoover chains, nhchainlength reset to 1");
949 ir->opts.nhchainlength = 1;
954 ir->opts.nhchainlength = 0;
957 if (ir->eI == eiVVAK)
959 sprintf(err_buf, "%s implemented primarily for validation, and requires nsttcouple = 1 and nstpcouple = 1.",
961 CHECK((ir->nsttcouple != 1) || (ir->nstpcouple != 1));
964 if (ETC_ANDERSEN(ir->etc))
966 sprintf(err_buf, "%s temperature control not supported for integrator %s.", etcoupl_names[ir->etc], ei_names[ir->eI]);
967 CHECK(!(EI_VV(ir->eI)));
969 if (ir->nstcomm > 0 && (ir->etc == etcANDERSEN))
971 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]);
972 warning_note(wi, warn_buf);
975 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]);
976 CHECK(ir->nstcomm > 1 && (ir->etc == etcANDERSEN));
979 if (ir->etc == etcBERENDSEN)
981 sprintf(warn_buf, "The %s thermostat does not generate the correct kinetic energy distribution. You might want to consider using the %s thermostat.",
982 ETCOUPLTYPE(ir->etc), ETCOUPLTYPE(etcVRESCALE));
983 warning_note(wi, warn_buf);
986 if ((ir->etc == etcNOSEHOOVER || ETC_ANDERSEN(ir->etc))
987 && ir->epc == epcBERENDSEN)
989 sprintf(warn_buf, "Using Berendsen pressure coupling invalidates the "
990 "true ensemble for the thermostat");
991 warning(wi, warn_buf);
994 /* PRESSURE COUPLING */
995 if (ir->epc == epcISOTROPIC)
997 ir->epc = epcBERENDSEN;
998 warning_note(wi, "Old option for pressure coupling given: "
999 "changing \"Isotropic\" to \"Berendsen\"\n");
1002 if (ir->epc != epcNO)
1004 dt_pcoupl = ir->nstpcouple*ir->delta_t;
1006 sprintf(err_buf, "tau-p must be > 0 instead of %g\n", ir->tau_p);
1007 CHECK(ir->tau_p <= 0);
1009 if (ir->tau_p/dt_pcoupl < pcouple_min_integration_steps(ir->epc))
1011 sprintf(warn_buf, "For proper integration of the %s barostat, tau-p (%g) should be at least %d times larger than nstpcouple*dt (%g)",
1012 EPCOUPLTYPE(ir->epc), ir->tau_p, pcouple_min_integration_steps(ir->epc), dt_pcoupl);
1013 warning(wi, warn_buf);
1016 sprintf(err_buf, "compressibility must be > 0 when using pressure"
1017 " coupling %s\n", EPCOUPLTYPE(ir->epc));
1018 CHECK(ir->compress[XX][XX] < 0 || ir->compress[YY][YY] < 0 ||
1019 ir->compress[ZZ][ZZ] < 0 ||
1020 (trace(ir->compress) == 0 && ir->compress[YY][XX] <= 0 &&
1021 ir->compress[ZZ][XX] <= 0 && ir->compress[ZZ][YY] <= 0));
1023 if (epcPARRINELLORAHMAN == ir->epc && opts->bGenVel)
1026 "You are generating velocities so I am assuming you "
1027 "are equilibrating a system. You are using "
1028 "%s pressure coupling, but this can be "
1029 "unstable for equilibration. If your system crashes, try "
1030 "equilibrating first with Berendsen pressure coupling. If "
1031 "you are not equilibrating the system, you can probably "
1032 "ignore this warning.",
1033 epcoupl_names[ir->epc]);
1034 warning(wi, warn_buf);
1040 if (ir->epc > epcNO)
1042 if ((ir->epc != epcBERENDSEN) && (ir->epc != epcMTTK))
1044 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.");
1050 if (ir->epc == epcMTTK)
1052 warning_error(wi, "MTTK pressure coupling requires a Velocity-verlet integrator");
1056 /* ELECTROSTATICS */
1057 /* More checks are in triple check (grompp.c) */
1059 if (ir->coulombtype == eelSWITCH)
1061 sprintf(warn_buf, "coulombtype = %s is only for testing purposes and can lead to serious "
1062 "artifacts, advice: use coulombtype = %s",
1063 eel_names[ir->coulombtype],
1064 eel_names[eelRF_ZERO]);
1065 warning(wi, warn_buf);
1068 if (ir->epsilon_r != 1 && ir->implicit_solvent == eisGBSA)
1070 sprintf(warn_buf, "epsilon-r = %g with GB implicit solvent, will use this value for inner dielectric", ir->epsilon_r);
1071 warning_note(wi, warn_buf);
1074 if (EEL_RF(ir->coulombtype) && ir->epsilon_rf == 1 && ir->epsilon_r != 1)
1076 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);
1077 warning(wi, warn_buf);
1078 ir->epsilon_rf = ir->epsilon_r;
1079 ir->epsilon_r = 1.0;
1082 if (getenv("GMX_DO_GALACTIC_DYNAMICS") == NULL)
1084 sprintf(err_buf, "epsilon-r must be >= 0 instead of %g\n", ir->epsilon_r);
1085 CHECK(ir->epsilon_r < 0);
1088 if (EEL_RF(ir->coulombtype))
1090 /* reaction field (at the cut-off) */
1092 if (ir->coulombtype == eelRF_ZERO)
1094 sprintf(warn_buf, "With coulombtype = %s, epsilon-rf must be 0, assuming you meant epsilon_rf=0",
1095 eel_names[ir->coulombtype]);
1096 CHECK(ir->epsilon_rf != 0);
1097 ir->epsilon_rf = 0.0;
1100 sprintf(err_buf, "epsilon-rf must be >= epsilon-r");
1101 CHECK((ir->epsilon_rf < ir->epsilon_r && ir->epsilon_rf != 0) ||
1102 (ir->epsilon_r == 0));
1103 if (ir->epsilon_rf == ir->epsilon_r)
1105 sprintf(warn_buf, "Using epsilon-rf = epsilon-r with %s does not make sense",
1106 eel_names[ir->coulombtype]);
1107 warning(wi, warn_buf);
1110 /* Allow rlist>rcoulomb for tabulated long range stuff. This just
1111 * means the interaction is zero outside rcoulomb, but it helps to
1112 * provide accurate energy conservation.
1114 if (ir_coulomb_might_be_zero_at_cutoff(ir))
1116 if (ir_coulomb_switched(ir))
1119 "With coulombtype = %s rcoulomb_switch must be < rcoulomb. Or, better: Use the potential modifier options!",
1120 eel_names[ir->coulombtype]);
1121 CHECK(ir->rcoulomb_switch >= ir->rcoulomb);
1124 else if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
1126 if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE)
1128 sprintf(err_buf, "With coulombtype = %s, rcoulomb should be >= rlist unless you use a potential modifier",
1129 eel_names[ir->coulombtype]);
1130 CHECK(ir->rlist > ir->rcoulomb);
1134 if (ir->coulombtype == eelSWITCH || ir->coulombtype == eelSHIFT)
1137 "Explicit switch/shift coulomb interactions cannot be used in combination with a secondary coulomb-modifier.");
1138 CHECK( ir->coulomb_modifier != eintmodNONE);
1140 if (ir->vdwtype == evdwSWITCH || ir->vdwtype == evdwSHIFT)
1143 "Explicit switch/shift vdw interactions cannot be used in combination with a secondary vdw-modifier.");
1144 CHECK( ir->vdw_modifier != eintmodNONE);
1147 if (ir->coulombtype == eelSWITCH || ir->coulombtype == eelSHIFT ||
1148 ir->vdwtype == evdwSWITCH || ir->vdwtype == evdwSHIFT)
1151 "The switch/shift interaction settings are just for compatibility; you will get better "
1152 "performance from applying potential modifiers to your interactions!\n");
1153 warning_note(wi, warn_buf);
1156 if (ir->coulombtype == eelPMESWITCH || ir->coulomb_modifier == eintmodPOTSWITCH)
1158 if (ir->rcoulomb_switch/ir->rcoulomb < 0.9499)
1160 real percentage = 100*(ir->rcoulomb-ir->rcoulomb_switch)/ir->rcoulomb;
1161 sprintf(warn_buf, "The switching range should be 5%% or less (currently %.2f%% using a switching range of %4f-%4f) for accurate electrostatic energies, energy conservation will be good regardless, since ewald_rtol = %g.",
1162 percentage, ir->rcoulomb_switch, ir->rcoulomb, ir->ewald_rtol);
1163 warning(wi, warn_buf);
1167 if (ir->vdwtype == evdwSWITCH || ir->vdw_modifier == eintmodPOTSWITCH)
1169 if (ir->rvdw_switch == 0)
1171 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.");
1172 warning(wi, warn_buf);
1176 if (EEL_FULL(ir->coulombtype))
1178 if (ir->coulombtype == eelPMESWITCH || ir->coulombtype == eelPMEUSER ||
1179 ir->coulombtype == eelPMEUSERSWITCH)
1181 sprintf(err_buf, "With coulombtype = %s, rcoulomb must be <= rlist",
1182 eel_names[ir->coulombtype]);
1183 CHECK(ir->rcoulomb > ir->rlist);
1185 else if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE)
1187 if (ir->coulombtype == eelPME || ir->coulombtype == eelP3M_AD)
1190 "With coulombtype = %s (without modifier), rcoulomb must be equal to rlist,\n"
1191 "or rlistlong if nstcalclr=1. For optimal energy conservation,consider using\n"
1192 "a potential modifier.", eel_names[ir->coulombtype]);
1193 if (ir->nstcalclr == 1)
1195 CHECK(ir->rcoulomb != ir->rlist && ir->rcoulomb != ir->rlistlong);
1199 CHECK(ir->rcoulomb != ir->rlist);
1205 if (EEL_PME(ir->coulombtype) || EVDW_PME(ir->vdwtype))
1207 if (ir->pme_order < 3)
1209 warning_error(wi, "pme-order can not be smaller than 3");
1213 if (ir->nwall == 2 && EEL_FULL(ir->coulombtype))
1215 if (ir->ewald_geometry == eewg3D)
1217 sprintf(warn_buf, "With pbc=%s you should use ewald-geometry=%s",
1218 epbc_names[ir->ePBC], eewg_names[eewg3DC]);
1219 warning(wi, warn_buf);
1221 /* This check avoids extra pbc coding for exclusion corrections */
1222 sprintf(err_buf, "wall-ewald-zfac should be >= 2");
1223 CHECK(ir->wall_ewald_zfac < 2);
1226 if (ir_vdw_switched(ir))
1228 sprintf(err_buf, "With switched vdw forces or potentials, rvdw-switch must be < rvdw");
1229 CHECK(ir->rvdw_switch >= ir->rvdw);
1231 if (ir->rvdw_switch < 0.5*ir->rvdw)
1233 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.",
1234 ir->rvdw_switch, ir->rvdw);
1235 warning_note(wi, warn_buf);
1238 else if (ir->vdwtype == evdwCUT || ir->vdwtype == evdwPME)
1240 if (ir->cutoff_scheme == ecutsGROUP && ir->vdw_modifier == eintmodNONE)
1242 sprintf(err_buf, "With vdwtype = %s, rvdw must be >= rlist unless you use a potential modifier", evdw_names[ir->vdwtype]);
1243 CHECK(ir->rlist > ir->rvdw);
1247 if (ir->vdwtype == evdwPME)
1249 if (!(ir->vdw_modifier == eintmodNONE || ir->vdw_modifier == eintmodPOTSHIFT))
1251 sprintf(err_buf, "With vdwtype = %s, the only supported modifiers are %s a\
1253 evdw_names[ir->vdwtype],
1254 eintmod_names[eintmodPOTSHIFT],
1255 eintmod_names[eintmodNONE]);
1259 if (ir->cutoff_scheme == ecutsGROUP)
1261 if (((ir->coulomb_modifier != eintmodNONE && ir->rcoulomb == ir->rlist) ||
1262 (ir->vdw_modifier != eintmodNONE && ir->rvdw == ir->rlist)) &&
1265 warning_note(wi, "With exact cut-offs, rlist should be "
1266 "larger than rcoulomb and rvdw, so that there "
1267 "is a buffer region for particle motion "
1268 "between neighborsearch steps");
1271 if (ir_coulomb_is_zero_at_cutoff(ir) && ir->rlistlong <= ir->rcoulomb)
1273 sprintf(warn_buf, "For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rcoulomb.",
1274 IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
1275 warning_note(wi, warn_buf);
1277 if (ir_vdw_switched(ir) && (ir->rlistlong <= ir->rvdw))
1279 sprintf(warn_buf, "For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rvdw.",
1280 IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
1281 warning_note(wi, warn_buf);
1285 if (ir->vdwtype == evdwUSER && ir->eDispCorr != edispcNO)
1287 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.");
1290 if (ir->nstlist == -1)
1292 sprintf(err_buf, "With nstlist=-1 rvdw and rcoulomb should be smaller than rlist to account for diffusion and possibly charge-group radii");
1293 CHECK(ir->rvdw >= ir->rlist || ir->rcoulomb >= ir->rlist);
1295 sprintf(err_buf, "nstlist can not be smaller than -1");
1296 CHECK(ir->nstlist < -1);
1298 if (ir->eI == eiLBFGS && (ir->coulombtype == eelCUT || ir->vdwtype == evdwCUT)
1301 warning(wi, "For efficient BFGS minimization, use switch/shift/pme instead of cut-off.");
1304 if (ir->eI == eiLBFGS && ir->nbfgscorr <= 0)
1306 warning(wi, "Using L-BFGS with nbfgscorr<=0 just gets you steepest descent.");
1309 /* ENERGY CONSERVATION */
1310 if (ir_NVE(ir) && ir->cutoff_scheme == ecutsGROUP)
1312 if (!ir_vdw_might_be_zero_at_cutoff(ir) && ir->rvdw > 0 && ir->vdw_modifier == eintmodNONE)
1314 sprintf(warn_buf, "You are using a cut-off for VdW interactions with NVE, for good energy conservation use vdwtype = %s (possibly with DispCorr)",
1315 evdw_names[evdwSHIFT]);
1316 warning_note(wi, warn_buf);
1318 if (!ir_coulomb_might_be_zero_at_cutoff(ir) && ir->rcoulomb > 0)
1320 sprintf(warn_buf, "You are using a cut-off for electrostatics with NVE, for good energy conservation use coulombtype = %s or %s",
1321 eel_names[eelPMESWITCH], eel_names[eelRF_ZERO]);
1322 warning_note(wi, warn_buf);
1326 if (EI_VV(ir->eI) && IR_TWINRANGE(*ir) && ir->nstlist > 1)
1328 sprintf(warn_buf, "Twin-range multiple time stepping does not work with integrator %s.", ei_names[ir->eI]);
1329 warning_error(wi, warn_buf);
1332 /* IMPLICIT SOLVENT */
1333 if (ir->coulombtype == eelGB_NOTUSED)
1335 ir->coulombtype = eelCUT;
1336 ir->implicit_solvent = eisGBSA;
1337 fprintf(stderr, "Note: Old option for generalized born electrostatics given:\n"
1338 "Changing coulombtype from \"generalized-born\" to \"cut-off\" and instead\n"
1339 "setting implicit-solvent value to \"GBSA\" in input section.\n");
1342 if (ir->sa_algorithm == esaSTILL)
1344 sprintf(err_buf, "Still SA algorithm not available yet, use %s or %s instead\n", esa_names[esaAPPROX], esa_names[esaNO]);
1345 CHECK(ir->sa_algorithm == esaSTILL);
1348 if (ir->implicit_solvent == eisGBSA)
1350 sprintf(err_buf, "With GBSA implicit solvent, rgbradii must be equal to rlist.");
1351 CHECK(ir->rgbradii != ir->rlist);
1353 if (ir->coulombtype != eelCUT)
1355 sprintf(err_buf, "With GBSA, coulombtype must be equal to %s\n", eel_names[eelCUT]);
1356 CHECK(ir->coulombtype != eelCUT);
1358 if (ir->vdwtype != evdwCUT)
1360 sprintf(err_buf, "With GBSA, vdw-type must be equal to %s\n", evdw_names[evdwCUT]);
1361 CHECK(ir->vdwtype != evdwCUT);
1363 if (ir->nstgbradii < 1)
1365 sprintf(warn_buf, "Using GBSA with nstgbradii<1, setting nstgbradii=1");
1366 warning_note(wi, warn_buf);
1369 if (ir->sa_algorithm == esaNO)
1371 sprintf(warn_buf, "No SA (non-polar) calculation requested together with GB. Are you sure this is what you want?\n");
1372 warning_note(wi, warn_buf);
1374 if (ir->sa_surface_tension < 0 && ir->sa_algorithm != esaNO)
1376 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");
1377 warning_note(wi, warn_buf);
1379 if (ir->gb_algorithm == egbSTILL)
1381 ir->sa_surface_tension = 0.0049 * CAL2JOULE * 100;
1385 ir->sa_surface_tension = 0.0054 * CAL2JOULE * 100;
1388 if (ir->sa_surface_tension == 0 && ir->sa_algorithm != esaNO)
1390 sprintf(err_buf, "Surface tension set to 0 while SA-calculation requested\n");
1391 CHECK(ir->sa_surface_tension == 0 && ir->sa_algorithm != esaNO);
1398 if (ir->cutoff_scheme != ecutsGROUP)
1400 warning_error(wi, "AdresS simulation supports only cutoff-scheme=group");
1404 warning_error(wi, "AdresS simulation supports only stochastic dynamics");
1406 if (ir->epc != epcNO)
1408 warning_error(wi, "AdresS simulation does not support pressure coupling");
1410 if (EEL_FULL(ir->coulombtype))
1412 warning_error(wi, "AdresS simulation does not support long-range electrostatics");
1417 /* count the number of text elemets separated by whitespace in a string.
1418 str = the input string
1419 maxptr = the maximum number of allowed elements
1420 ptr = the output array of pointers to the first character of each element
1421 returns: the number of elements. */
1422 int str_nelem(const char *str, int maxptr, char *ptr[])
1427 copy0 = gmx_strdup(str);
1430 while (*copy != '\0')
1434 gmx_fatal(FARGS, "Too many groups on line: '%s' (max is %d)",
1442 while ((*copy != '\0') && !isspace(*copy))
1461 /* interpret a number of doubles from a string and put them in an array,
1462 after allocating space for them.
1463 str = the input string
1464 n = the (pre-allocated) number of doubles read
1465 r = the output array of doubles. */
1466 static void parse_n_real(char *str, int *n, real **r)
1471 *n = str_nelem(str, MAXPTR, ptr);
1474 for (i = 0; i < *n; i++)
1476 (*r)[i] = strtod(ptr[i], NULL);
1480 static void do_fep_params(t_inputrec *ir, char fep_lambda[][STRLEN], char weights[STRLEN])
1483 int i, j, max_n_lambda, nweights, nfep[efptNR];
1484 t_lambda *fep = ir->fepvals;
1485 t_expanded *expand = ir->expandedvals;
1486 real **count_fep_lambdas;
1487 gmx_bool bOneLambda = TRUE;
1489 snew(count_fep_lambdas, efptNR);
1491 /* FEP input processing */
1492 /* first, identify the number of lambda values for each type.
1493 All that are nonzero must have the same number */
1495 for (i = 0; i < efptNR; i++)
1497 parse_n_real(fep_lambda[i], &(nfep[i]), &(count_fep_lambdas[i]));
1500 /* now, determine the number of components. All must be either zero, or equal. */
1503 for (i = 0; i < efptNR; i++)
1505 if (nfep[i] > max_n_lambda)
1507 max_n_lambda = nfep[i]; /* here's a nonzero one. All of them
1508 must have the same number if its not zero.*/
1513 for (i = 0; i < efptNR; i++)
1517 ir->fepvals->separate_dvdl[i] = FALSE;
1519 else if (nfep[i] == max_n_lambda)
1521 if (i != efptTEMPERATURE) /* we treat this differently -- not really a reason to compute the derivative with
1522 respect to the temperature currently */
1524 ir->fepvals->separate_dvdl[i] = TRUE;
1529 gmx_fatal(FARGS, "Number of lambdas (%d) for FEP type %s not equal to number of other types (%d)",
1530 nfep[i], efpt_names[i], max_n_lambda);
1533 /* we don't print out dhdl if the temperature is changing, since we can't correctly define dhdl in this case */
1534 ir->fepvals->separate_dvdl[efptTEMPERATURE] = FALSE;
1536 /* the number of lambdas is the number we've read in, which is either zero
1537 or the same for all */
1538 fep->n_lambda = max_n_lambda;
1540 /* allocate space for the array of lambda values */
1541 snew(fep->all_lambda, efptNR);
1542 /* if init_lambda is defined, we need to set lambda */
1543 if ((fep->init_lambda > 0) && (fep->n_lambda == 0))
1545 ir->fepvals->separate_dvdl[efptFEP] = TRUE;
1547 /* otherwise allocate the space for all of the lambdas, and transfer the data */
1548 for (i = 0; i < efptNR; i++)
1550 snew(fep->all_lambda[i], fep->n_lambda);
1551 if (nfep[i] > 0) /* if it's zero, then the count_fep_lambda arrays
1554 for (j = 0; j < fep->n_lambda; j++)
1556 fep->all_lambda[i][j] = (double)count_fep_lambdas[i][j];
1558 sfree(count_fep_lambdas[i]);
1561 sfree(count_fep_lambdas);
1563 /* "fep-vals" is either zero or the full number. If zero, we'll need to define fep-lambdas for internal
1564 bookkeeping -- for now, init_lambda */
1566 if ((nfep[efptFEP] == 0) && (fep->init_lambda >= 0))
1568 for (i = 0; i < fep->n_lambda; i++)
1570 fep->all_lambda[efptFEP][i] = fep->init_lambda;
1574 /* check to see if only a single component lambda is defined, and soft core is defined.
1575 In this case, turn on coulomb soft core */
1577 if (max_n_lambda == 0)
1583 for (i = 0; i < efptNR; i++)
1585 if ((nfep[i] != 0) && (i != efptFEP))
1591 if ((bOneLambda) && (fep->sc_alpha > 0))
1593 fep->bScCoul = TRUE;
1596 /* Fill in the others with the efptFEP if they are not explicitly
1597 specified (i.e. nfep[i] == 0). This means if fep is not defined,
1598 they are all zero. */
1600 for (i = 0; i < efptNR; i++)
1602 if ((nfep[i] == 0) && (i != efptFEP))
1604 for (j = 0; j < fep->n_lambda; j++)
1606 fep->all_lambda[i][j] = fep->all_lambda[efptFEP][j];
1612 /* make it easier if sc_r_power = 48 by increasing it to the 4th power, to be in the right scale. */
1613 if (fep->sc_r_power == 48)
1615 if (fep->sc_alpha > 0.1)
1617 gmx_fatal(FARGS, "sc_alpha (%f) for sc_r_power = 48 should usually be between 0.001 and 0.004", fep->sc_alpha);
1621 expand = ir->expandedvals;
1622 /* now read in the weights */
1623 parse_n_real(weights, &nweights, &(expand->init_lambda_weights));
1626 snew(expand->init_lambda_weights, fep->n_lambda); /* initialize to zero */
1628 else if (nweights != fep->n_lambda)
1630 gmx_fatal(FARGS, "Number of weights (%d) is not equal to number of lambda values (%d)",
1631 nweights, fep->n_lambda);
1633 if ((expand->nstexpanded < 0) && (ir->efep != efepNO))
1635 expand->nstexpanded = fep->nstdhdl;
1636 /* if you don't specify nstexpanded when doing expanded ensemble free energy calcs, it is set to nstdhdl */
1638 if ((expand->nstexpanded < 0) && ir->bSimTemp)
1640 expand->nstexpanded = 2*(int)(ir->opts.tau_t[0]/ir->delta_t);
1641 /* if you don't specify nstexpanded when doing expanded ensemble simulated tempering, it is set to
1642 2*tau_t just to be careful so it's not to frequent */
1647 static void do_simtemp_params(t_inputrec *ir)
1650 snew(ir->simtempvals->temperatures, ir->fepvals->n_lambda);
1651 GetSimTemps(ir->fepvals->n_lambda, ir->simtempvals, ir->fepvals->all_lambda[efptTEMPERATURE]);
1656 static void do_wall_params(t_inputrec *ir,
1657 char *wall_atomtype, char *wall_density,
1661 char *names[MAXPTR];
1664 opts->wall_atomtype[0] = NULL;
1665 opts->wall_atomtype[1] = NULL;
1667 ir->wall_atomtype[0] = -1;
1668 ir->wall_atomtype[1] = -1;
1669 ir->wall_density[0] = 0;
1670 ir->wall_density[1] = 0;
1674 nstr = str_nelem(wall_atomtype, MAXPTR, names);
1675 if (nstr != ir->nwall)
1677 gmx_fatal(FARGS, "Expected %d elements for wall_atomtype, found %d",
1680 for (i = 0; i < ir->nwall; i++)
1682 opts->wall_atomtype[i] = gmx_strdup(names[i]);
1685 if (ir->wall_type == ewt93 || ir->wall_type == ewt104)
1687 nstr = str_nelem(wall_density, MAXPTR, names);
1688 if (nstr != ir->nwall)
1690 gmx_fatal(FARGS, "Expected %d elements for wall-density, found %d", ir->nwall, nstr);
1692 for (i = 0; i < ir->nwall; i++)
1694 sscanf(names[i], "%lf", &dbl);
1697 gmx_fatal(FARGS, "wall-density[%d] = %f\n", i, dbl);
1699 ir->wall_density[i] = dbl;
1705 static void add_wall_energrps(gmx_groups_t *groups, int nwall, t_symtab *symtab)
1713 srenew(groups->grpname, groups->ngrpname+nwall);
1714 grps = &(groups->grps[egcENER]);
1715 srenew(grps->nm_ind, grps->nr+nwall);
1716 for (i = 0; i < nwall; i++)
1718 sprintf(str, "wall%d", i);
1719 groups->grpname[groups->ngrpname] = put_symtab(symtab, str);
1720 grps->nm_ind[grps->nr++] = groups->ngrpname++;
1725 void read_expandedparams(int *ninp_p, t_inpfile **inp_p,
1726 t_expanded *expand, warninp_t wi)
1728 int ninp, nerror = 0;
1734 /* read expanded ensemble parameters */
1735 CCTYPE ("expanded ensemble variables");
1736 ITYPE ("nstexpanded", expand->nstexpanded, -1);
1737 EETYPE("lmc-stats", expand->elamstats, elamstats_names);
1738 EETYPE("lmc-move", expand->elmcmove, elmcmove_names);
1739 EETYPE("lmc-weights-equil", expand->elmceq, elmceq_names);
1740 ITYPE ("weight-equil-number-all-lambda", expand->equil_n_at_lam, -1);
1741 ITYPE ("weight-equil-number-samples", expand->equil_samples, -1);
1742 ITYPE ("weight-equil-number-steps", expand->equil_steps, -1);
1743 RTYPE ("weight-equil-wl-delta", expand->equil_wl_delta, -1);
1744 RTYPE ("weight-equil-count-ratio", expand->equil_ratio, -1);
1745 CCTYPE("Seed for Monte Carlo in lambda space");
1746 ITYPE ("lmc-seed", expand->lmc_seed, -1);
1747 RTYPE ("mc-temperature", expand->mc_temp, -1);
1748 ITYPE ("lmc-repeats", expand->lmc_repeats, 1);
1749 ITYPE ("lmc-gibbsdelta", expand->gibbsdeltalam, -1);
1750 ITYPE ("lmc-forced-nstart", expand->lmc_forced_nstart, 0);
1751 EETYPE("symmetrized-transition-matrix", expand->bSymmetrizedTMatrix, yesno_names);
1752 ITYPE("nst-transition-matrix", expand->nstTij, -1);
1753 ITYPE ("mininum-var-min", expand->minvarmin, 100); /*default is reasonable */
1754 ITYPE ("weight-c-range", expand->c_range, 0); /* default is just C=0 */
1755 RTYPE ("wl-scale", expand->wl_scale, 0.8);
1756 RTYPE ("wl-ratio", expand->wl_ratio, 0.8);
1757 RTYPE ("init-wl-delta", expand->init_wl_delta, 1.0);
1758 EETYPE("wl-oneovert", expand->bWLoneovert, yesno_names);
1766 void get_ir(const char *mdparin, const char *mdparout,
1767 t_inputrec *ir, t_gromppopts *opts,
1771 double dumdub[2][6];
1775 char warn_buf[STRLEN];
1776 t_lambda *fep = ir->fepvals;
1777 t_expanded *expand = ir->expandedvals;
1779 init_inputrec_strings();
1780 inp = read_inpfile(mdparin, &ninp, wi);
1782 snew(dumstr[0], STRLEN);
1783 snew(dumstr[1], STRLEN);
1785 if (-1 == search_einp(ninp, inp, "cutoff-scheme"))
1788 "%s did not specify a value for the .mdp option "
1789 "\"cutoff-scheme\". Probably it was first intended for use "
1790 "with GROMACS before 4.6. In 4.6, the Verlet scheme was "
1791 "introduced, but the group scheme was still the default. "
1792 "The default is now the Verlet scheme, so you will observe "
1793 "different behaviour.", mdparin);
1794 warning_note(wi, warn_buf);
1797 /* ignore the following deprecated commands */
1800 REM_TYPE("domain-decomposition");
1801 REM_TYPE("andersen-seed");
1803 REM_TYPE("dihre-fc");
1804 REM_TYPE("dihre-tau");
1805 REM_TYPE("nstdihreout");
1806 REM_TYPE("nstcheckpoint");
1807 REM_TYPE("optimize-fft");
1809 /* replace the following commands with the clearer new versions*/
1810 REPL_TYPE("unconstrained-start", "continuation");
1811 REPL_TYPE("foreign-lambda", "fep-lambdas");
1812 REPL_TYPE("verlet-buffer-drift", "verlet-buffer-tolerance");
1813 REPL_TYPE("nstxtcout", "nstxout-compressed");
1814 REPL_TYPE("xtc-grps", "compressed-x-grps");
1815 REPL_TYPE("xtc-precision", "compressed-x-precision");
1817 CCTYPE ("VARIOUS PREPROCESSING OPTIONS");
1818 CTYPE ("Preprocessor information: use cpp syntax.");
1819 CTYPE ("e.g.: -I/home/joe/doe -I/home/mary/roe");
1820 STYPE ("include", opts->include, NULL);
1821 CTYPE ("e.g.: -DPOSRES -DFLEXIBLE (note these variable names are case sensitive)");
1822 STYPE ("define", opts->define, NULL);
1824 CCTYPE ("RUN CONTROL PARAMETERS");
1825 EETYPE("integrator", ir->eI, ei_names);
1826 CTYPE ("Start time and timestep in ps");
1827 RTYPE ("tinit", ir->init_t, 0.0);
1828 RTYPE ("dt", ir->delta_t, 0.001);
1829 STEPTYPE ("nsteps", ir->nsteps, 0);
1830 CTYPE ("For exact run continuation or redoing part of a run");
1831 STEPTYPE ("init-step", ir->init_step, 0);
1832 CTYPE ("Part index is updated automatically on checkpointing (keeps files separate)");
1833 ITYPE ("simulation-part", ir->simulation_part, 1);
1834 CTYPE ("mode for center of mass motion removal");
1835 EETYPE("comm-mode", ir->comm_mode, ecm_names);
1836 CTYPE ("number of steps for center of mass motion removal");
1837 ITYPE ("nstcomm", ir->nstcomm, 100);
1838 CTYPE ("group(s) for center of mass motion removal");
1839 STYPE ("comm-grps", is->vcm, NULL);
1841 CCTYPE ("LANGEVIN DYNAMICS OPTIONS");
1842 CTYPE ("Friction coefficient (amu/ps) and random seed");
1843 RTYPE ("bd-fric", ir->bd_fric, 0.0);
1844 STEPTYPE ("ld-seed", ir->ld_seed, -1);
1847 CCTYPE ("ENERGY MINIMIZATION OPTIONS");
1848 CTYPE ("Force tolerance and initial step-size");
1849 RTYPE ("emtol", ir->em_tol, 10.0);
1850 RTYPE ("emstep", ir->em_stepsize, 0.01);
1851 CTYPE ("Max number of iterations in relax-shells");
1852 ITYPE ("niter", ir->niter, 20);
1853 CTYPE ("Step size (ps^2) for minimization of flexible constraints");
1854 RTYPE ("fcstep", ir->fc_stepsize, 0);
1855 CTYPE ("Frequency of steepest descents steps when doing CG");
1856 ITYPE ("nstcgsteep", ir->nstcgsteep, 1000);
1857 ITYPE ("nbfgscorr", ir->nbfgscorr, 10);
1859 CCTYPE ("TEST PARTICLE INSERTION OPTIONS");
1860 RTYPE ("rtpi", ir->rtpi, 0.05);
1862 /* Output options */
1863 CCTYPE ("OUTPUT CONTROL OPTIONS");
1864 CTYPE ("Output frequency for coords (x), velocities (v) and forces (f)");
1865 ITYPE ("nstxout", ir->nstxout, 0);
1866 ITYPE ("nstvout", ir->nstvout, 0);
1867 ITYPE ("nstfout", ir->nstfout, 0);
1868 CTYPE ("Output frequency for energies to log file and energy file");
1869 ITYPE ("nstlog", ir->nstlog, 1000);
1870 ITYPE ("nstcalcenergy", ir->nstcalcenergy, 100);
1871 ITYPE ("nstenergy", ir->nstenergy, 1000);
1872 CTYPE ("Output frequency and precision for .xtc file");
1873 ITYPE ("nstxout-compressed", ir->nstxout_compressed, 0);
1874 RTYPE ("compressed-x-precision", ir->x_compression_precision, 1000.0);
1875 CTYPE ("This selects the subset of atoms for the compressed");
1876 CTYPE ("trajectory file. You can select multiple groups. By");
1877 CTYPE ("default, all atoms will be written.");
1878 STYPE ("compressed-x-grps", is->x_compressed_groups, NULL);
1879 CTYPE ("Selection of energy groups");
1880 STYPE ("energygrps", is->energy, NULL);
1882 /* Neighbor searching */
1883 CCTYPE ("NEIGHBORSEARCHING PARAMETERS");
1884 CTYPE ("cut-off scheme (Verlet: particle based cut-offs, group: using charge groups)");
1885 EETYPE("cutoff-scheme", ir->cutoff_scheme, ecutscheme_names);
1886 CTYPE ("nblist update frequency");
1887 ITYPE ("nstlist", ir->nstlist, 10);
1888 CTYPE ("ns algorithm (simple or grid)");
1889 EETYPE("ns-type", ir->ns_type, ens_names);
1890 CTYPE ("Periodic boundary conditions: xyz, no, xy");
1891 EETYPE("pbc", ir->ePBC, epbc_names);
1892 EETYPE("periodic-molecules", ir->bPeriodicMols, yesno_names);
1893 CTYPE ("Allowed energy error due to the Verlet buffer in kJ/mol/ps per atom,");
1894 CTYPE ("a value of -1 means: use rlist");
1895 RTYPE("verlet-buffer-tolerance", ir->verletbuf_tol, 0.005);
1896 CTYPE ("nblist cut-off");
1897 RTYPE ("rlist", ir->rlist, 1.0);
1898 CTYPE ("long-range cut-off for switched potentials");
1899 RTYPE ("rlistlong", ir->rlistlong, -1);
1900 ITYPE ("nstcalclr", ir->nstcalclr, -1);
1902 /* Electrostatics */
1903 CCTYPE ("OPTIONS FOR ELECTROSTATICS AND VDW");
1904 CTYPE ("Method for doing electrostatics");
1905 EETYPE("coulombtype", ir->coulombtype, eel_names);
1906 EETYPE("coulomb-modifier", ir->coulomb_modifier, eintmod_names);
1907 CTYPE ("cut-off lengths");
1908 RTYPE ("rcoulomb-switch", ir->rcoulomb_switch, 0.0);
1909 RTYPE ("rcoulomb", ir->rcoulomb, 1.0);
1910 CTYPE ("Relative dielectric constant for the medium and the reaction field");
1911 RTYPE ("epsilon-r", ir->epsilon_r, 1.0);
1912 RTYPE ("epsilon-rf", ir->epsilon_rf, 0.0);
1913 CTYPE ("Method for doing Van der Waals");
1914 EETYPE("vdw-type", ir->vdwtype, evdw_names);
1915 EETYPE("vdw-modifier", ir->vdw_modifier, eintmod_names);
1916 CTYPE ("cut-off lengths");
1917 RTYPE ("rvdw-switch", ir->rvdw_switch, 0.0);
1918 RTYPE ("rvdw", ir->rvdw, 1.0);
1919 CTYPE ("Apply long range dispersion corrections for Energy and Pressure");
1920 EETYPE("DispCorr", ir->eDispCorr, edispc_names);
1921 CTYPE ("Extension of the potential lookup tables beyond the cut-off");
1922 RTYPE ("table-extension", ir->tabext, 1.0);
1923 CTYPE ("Separate tables between energy group pairs");
1924 STYPE ("energygrp-table", is->egptable, NULL);
1925 CTYPE ("Spacing for the PME/PPPM FFT grid");
1926 RTYPE ("fourierspacing", ir->fourier_spacing, 0.12);
1927 CTYPE ("FFT grid size, when a value is 0 fourierspacing will be used");
1928 ITYPE ("fourier-nx", ir->nkx, 0);
1929 ITYPE ("fourier-ny", ir->nky, 0);
1930 ITYPE ("fourier-nz", ir->nkz, 0);
1931 CTYPE ("EWALD/PME/PPPM parameters");
1932 ITYPE ("pme-order", ir->pme_order, 4);
1933 RTYPE ("ewald-rtol", ir->ewald_rtol, 0.00001);
1934 RTYPE ("ewald-rtol-lj", ir->ewald_rtol_lj, 0.001);
1935 EETYPE("lj-pme-comb-rule", ir->ljpme_combination_rule, eljpme_names);
1936 EETYPE("ewald-geometry", ir->ewald_geometry, eewg_names);
1937 RTYPE ("epsilon-surface", ir->epsilon_surface, 0.0);
1939 CCTYPE("IMPLICIT SOLVENT ALGORITHM");
1940 EETYPE("implicit-solvent", ir->implicit_solvent, eis_names);
1942 CCTYPE ("GENERALIZED BORN ELECTROSTATICS");
1943 CTYPE ("Algorithm for calculating Born radii");
1944 EETYPE("gb-algorithm", ir->gb_algorithm, egb_names);
1945 CTYPE ("Frequency of calculating the Born radii inside rlist");
1946 ITYPE ("nstgbradii", ir->nstgbradii, 1);
1947 CTYPE ("Cutoff for Born radii calculation; the contribution from atoms");
1948 CTYPE ("between rlist and rgbradii is updated every nstlist steps");
1949 RTYPE ("rgbradii", ir->rgbradii, 1.0);
1950 CTYPE ("Dielectric coefficient of the implicit solvent");
1951 RTYPE ("gb-epsilon-solvent", ir->gb_epsilon_solvent, 80.0);
1952 CTYPE ("Salt concentration in M for Generalized Born models");
1953 RTYPE ("gb-saltconc", ir->gb_saltconc, 0.0);
1954 CTYPE ("Scaling factors used in the OBC GB model. Default values are OBC(II)");
1955 RTYPE ("gb-obc-alpha", ir->gb_obc_alpha, 1.0);
1956 RTYPE ("gb-obc-beta", ir->gb_obc_beta, 0.8);
1957 RTYPE ("gb-obc-gamma", ir->gb_obc_gamma, 4.85);
1958 RTYPE ("gb-dielectric-offset", ir->gb_dielectric_offset, 0.009);
1959 EETYPE("sa-algorithm", ir->sa_algorithm, esa_names);
1960 CTYPE ("Surface tension (kJ/mol/nm^2) for the SA (nonpolar surface) part of GBSA");
1961 CTYPE ("The value -1 will set default value for Still/HCT/OBC GB-models.");
1962 RTYPE ("sa-surface-tension", ir->sa_surface_tension, -1);
1964 /* Coupling stuff */
1965 CCTYPE ("OPTIONS FOR WEAK COUPLING ALGORITHMS");
1966 CTYPE ("Temperature coupling");
1967 EETYPE("tcoupl", ir->etc, etcoupl_names);
1968 ITYPE ("nsttcouple", ir->nsttcouple, -1);
1969 ITYPE("nh-chain-length", ir->opts.nhchainlength, 10);
1970 EETYPE("print-nose-hoover-chain-variables", ir->bPrintNHChains, yesno_names);
1971 CTYPE ("Groups to couple separately");
1972 STYPE ("tc-grps", is->tcgrps, NULL);
1973 CTYPE ("Time constant (ps) and reference temperature (K)");
1974 STYPE ("tau-t", is->tau_t, NULL);
1975 STYPE ("ref-t", is->ref_t, NULL);
1976 CTYPE ("pressure coupling");
1977 EETYPE("pcoupl", ir->epc, epcoupl_names);
1978 EETYPE("pcoupltype", ir->epct, epcoupltype_names);
1979 ITYPE ("nstpcouple", ir->nstpcouple, -1);
1980 CTYPE ("Time constant (ps), compressibility (1/bar) and reference P (bar)");
1981 RTYPE ("tau-p", ir->tau_p, 1.0);
1982 STYPE ("compressibility", dumstr[0], NULL);
1983 STYPE ("ref-p", dumstr[1], NULL);
1984 CTYPE ("Scaling of reference coordinates, No, All or COM");
1985 EETYPE ("refcoord-scaling", ir->refcoord_scaling, erefscaling_names);
1988 CCTYPE ("OPTIONS FOR QMMM calculations");
1989 EETYPE("QMMM", ir->bQMMM, yesno_names);
1990 CTYPE ("Groups treated Quantum Mechanically");
1991 STYPE ("QMMM-grps", is->QMMM, NULL);
1992 CTYPE ("QM method");
1993 STYPE("QMmethod", is->QMmethod, NULL);
1994 CTYPE ("QMMM scheme");
1995 EETYPE("QMMMscheme", ir->QMMMscheme, eQMMMscheme_names);
1996 CTYPE ("QM basisset");
1997 STYPE("QMbasis", is->QMbasis, NULL);
1998 CTYPE ("QM charge");
1999 STYPE ("QMcharge", is->QMcharge, NULL);
2000 CTYPE ("QM multiplicity");
2001 STYPE ("QMmult", is->QMmult, NULL);
2002 CTYPE ("Surface Hopping");
2003 STYPE ("SH", is->bSH, NULL);
2004 CTYPE ("CAS space options");
2005 STYPE ("CASorbitals", is->CASorbitals, NULL);
2006 STYPE ("CASelectrons", is->CASelectrons, NULL);
2007 STYPE ("SAon", is->SAon, NULL);
2008 STYPE ("SAoff", is->SAoff, NULL);
2009 STYPE ("SAsteps", is->SAsteps, NULL);
2010 CTYPE ("Scale factor for MM charges");
2011 RTYPE ("MMChargeScaleFactor", ir->scalefactor, 1.0);
2012 CTYPE ("Optimization of QM subsystem");
2013 STYPE ("bOPT", is->bOPT, NULL);
2014 STYPE ("bTS", is->bTS, NULL);
2016 /* Simulated annealing */
2017 CCTYPE("SIMULATED ANNEALING");
2018 CTYPE ("Type of annealing for each temperature group (no/single/periodic)");
2019 STYPE ("annealing", is->anneal, NULL);
2020 CTYPE ("Number of time points to use for specifying annealing in each group");
2021 STYPE ("annealing-npoints", is->anneal_npoints, NULL);
2022 CTYPE ("List of times at the annealing points for each group");
2023 STYPE ("annealing-time", is->anneal_time, NULL);
2024 CTYPE ("Temp. at each annealing point, for each group.");
2025 STYPE ("annealing-temp", is->anneal_temp, NULL);
2028 CCTYPE ("GENERATE VELOCITIES FOR STARTUP RUN");
2029 EETYPE("gen-vel", opts->bGenVel, yesno_names);
2030 RTYPE ("gen-temp", opts->tempi, 300.0);
2031 ITYPE ("gen-seed", opts->seed, -1);
2034 CCTYPE ("OPTIONS FOR BONDS");
2035 EETYPE("constraints", opts->nshake, constraints);
2036 CTYPE ("Type of constraint algorithm");
2037 EETYPE("constraint-algorithm", ir->eConstrAlg, econstr_names);
2038 CTYPE ("Do not constrain the start configuration");
2039 EETYPE("continuation", ir->bContinuation, yesno_names);
2040 CTYPE ("Use successive overrelaxation to reduce the number of shake iterations");
2041 EETYPE("Shake-SOR", ir->bShakeSOR, yesno_names);
2042 CTYPE ("Relative tolerance of shake");
2043 RTYPE ("shake-tol", ir->shake_tol, 0.0001);
2044 CTYPE ("Highest order in the expansion of the constraint coupling matrix");
2045 ITYPE ("lincs-order", ir->nProjOrder, 4);
2046 CTYPE ("Number of iterations in the final step of LINCS. 1 is fine for");
2047 CTYPE ("normal simulations, but use 2 to conserve energy in NVE runs.");
2048 CTYPE ("For energy minimization with constraints it should be 4 to 8.");
2049 ITYPE ("lincs-iter", ir->nLincsIter, 1);
2050 CTYPE ("Lincs will write a warning to the stderr if in one step a bond");
2051 CTYPE ("rotates over more degrees than");
2052 RTYPE ("lincs-warnangle", ir->LincsWarnAngle, 30.0);
2053 CTYPE ("Convert harmonic bonds to morse potentials");
2054 EETYPE("morse", opts->bMorse, yesno_names);
2056 /* Energy group exclusions */
2057 CCTYPE ("ENERGY GROUP EXCLUSIONS");
2058 CTYPE ("Pairs of energy groups for which all non-bonded interactions are excluded");
2059 STYPE ("energygrp-excl", is->egpexcl, NULL);
2063 CTYPE ("Number of walls, type, atom types, densities and box-z scale factor for Ewald");
2064 ITYPE ("nwall", ir->nwall, 0);
2065 EETYPE("wall-type", ir->wall_type, ewt_names);
2066 RTYPE ("wall-r-linpot", ir->wall_r_linpot, -1);
2067 STYPE ("wall-atomtype", is->wall_atomtype, NULL);
2068 STYPE ("wall-density", is->wall_density, NULL);
2069 RTYPE ("wall-ewald-zfac", ir->wall_ewald_zfac, 3);
2072 CCTYPE("COM PULLING");
2073 CTYPE("Pull type: no, umbrella, constraint or constant-force");
2074 EETYPE("pull", ir->ePull, epull_names);
2075 if (ir->ePull != epullNO)
2078 is->pull_grp = read_pullparams(&ninp, &inp, ir->pull, &opts->pull_start, wi);
2081 /* Enforced rotation */
2082 CCTYPE("ENFORCED ROTATION");
2083 CTYPE("Enforced rotation: No or Yes");
2084 EETYPE("rotation", ir->bRot, yesno_names);
2088 is->rot_grp = read_rotparams(&ninp, &inp, ir->rot, wi);
2091 /* Interactive MD */
2093 CCTYPE("Group to display and/or manipulate in interactive MD session");
2094 STYPE ("IMD-group", is->imd_grp, NULL);
2095 if (is->imd_grp[0] != '\0')
2102 CCTYPE("NMR refinement stuff");
2103 CTYPE ("Distance restraints type: No, Simple or Ensemble");
2104 EETYPE("disre", ir->eDisre, edisre_names);
2105 CTYPE ("Force weighting of pairs in one distance restraint: Conservative or Equal");
2106 EETYPE("disre-weighting", ir->eDisreWeighting, edisreweighting_names);
2107 CTYPE ("Use sqrt of the time averaged times the instantaneous violation");
2108 EETYPE("disre-mixed", ir->bDisreMixed, yesno_names);
2109 RTYPE ("disre-fc", ir->dr_fc, 1000.0);
2110 RTYPE ("disre-tau", ir->dr_tau, 0.0);
2111 CTYPE ("Output frequency for pair distances to energy file");
2112 ITYPE ("nstdisreout", ir->nstdisreout, 100);
2113 CTYPE ("Orientation restraints: No or Yes");
2114 EETYPE("orire", opts->bOrire, yesno_names);
2115 CTYPE ("Orientation restraints force constant and tau for time averaging");
2116 RTYPE ("orire-fc", ir->orires_fc, 0.0);
2117 RTYPE ("orire-tau", ir->orires_tau, 0.0);
2118 STYPE ("orire-fitgrp", is->orirefitgrp, NULL);
2119 CTYPE ("Output frequency for trace(SD) and S to energy file");
2120 ITYPE ("nstorireout", ir->nstorireout, 100);
2122 /* free energy variables */
2123 CCTYPE ("Free energy variables");
2124 EETYPE("free-energy", ir->efep, efep_names);
2125 STYPE ("couple-moltype", is->couple_moltype, NULL);
2126 EETYPE("couple-lambda0", opts->couple_lam0, couple_lam);
2127 EETYPE("couple-lambda1", opts->couple_lam1, couple_lam);
2128 EETYPE("couple-intramol", opts->bCoupleIntra, yesno_names);
2130 RTYPE ("init-lambda", fep->init_lambda, -1); /* start with -1 so
2132 it was not entered */
2133 ITYPE ("init-lambda-state", fep->init_fep_state, -1);
2134 RTYPE ("delta-lambda", fep->delta_lambda, 0.0);
2135 ITYPE ("nstdhdl", fep->nstdhdl, 50);
2136 STYPE ("fep-lambdas", is->fep_lambda[efptFEP], NULL);
2137 STYPE ("mass-lambdas", is->fep_lambda[efptMASS], NULL);
2138 STYPE ("coul-lambdas", is->fep_lambda[efptCOUL], NULL);
2139 STYPE ("vdw-lambdas", is->fep_lambda[efptVDW], NULL);
2140 STYPE ("bonded-lambdas", is->fep_lambda[efptBONDED], NULL);
2141 STYPE ("restraint-lambdas", is->fep_lambda[efptRESTRAINT], NULL);
2142 STYPE ("temperature-lambdas", is->fep_lambda[efptTEMPERATURE], NULL);
2143 ITYPE ("calc-lambda-neighbors", fep->lambda_neighbors, 1);
2144 STYPE ("init-lambda-weights", is->lambda_weights, NULL);
2145 EETYPE("dhdl-print-energy", fep->edHdLPrintEnergy, edHdLPrintEnergy_names);
2146 RTYPE ("sc-alpha", fep->sc_alpha, 0.0);
2147 ITYPE ("sc-power", fep->sc_power, 1);
2148 RTYPE ("sc-r-power", fep->sc_r_power, 6.0);
2149 RTYPE ("sc-sigma", fep->sc_sigma, 0.3);
2150 EETYPE("sc-coul", fep->bScCoul, yesno_names);
2151 ITYPE ("dh_hist_size", fep->dh_hist_size, 0);
2152 RTYPE ("dh_hist_spacing", fep->dh_hist_spacing, 0.1);
2153 EETYPE("separate-dhdl-file", fep->separate_dhdl_file,
2154 separate_dhdl_file_names);
2155 EETYPE("dhdl-derivatives", fep->dhdl_derivatives, dhdl_derivatives_names);
2156 ITYPE ("dh_hist_size", fep->dh_hist_size, 0);
2157 RTYPE ("dh_hist_spacing", fep->dh_hist_spacing, 0.1);
2159 /* Non-equilibrium MD stuff */
2160 CCTYPE("Non-equilibrium MD stuff");
2161 STYPE ("acc-grps", is->accgrps, NULL);
2162 STYPE ("accelerate", is->acc, NULL);
2163 STYPE ("freezegrps", is->freeze, NULL);
2164 STYPE ("freezedim", is->frdim, NULL);
2165 RTYPE ("cos-acceleration", ir->cos_accel, 0);
2166 STYPE ("deform", is->deform, NULL);
2168 /* simulated tempering variables */
2169 CCTYPE("simulated tempering variables");
2170 EETYPE("simulated-tempering", ir->bSimTemp, yesno_names);
2171 EETYPE("simulated-tempering-scaling", ir->simtempvals->eSimTempScale, esimtemp_names);
2172 RTYPE("sim-temp-low", ir->simtempvals->simtemp_low, 300.0);
2173 RTYPE("sim-temp-high", ir->simtempvals->simtemp_high, 300.0);
2175 /* expanded ensemble variables */
2176 if (ir->efep == efepEXPANDED || ir->bSimTemp)
2178 read_expandedparams(&ninp, &inp, expand, wi);
2181 /* Electric fields */
2182 CCTYPE("Electric fields");
2183 CTYPE ("Format is number of terms (int) and for all terms an amplitude (real)");
2184 CTYPE ("and a phase angle (real)");
2185 STYPE ("E-x", is->efield_x, NULL);
2186 STYPE ("E-xt", is->efield_xt, NULL);
2187 STYPE ("E-y", is->efield_y, NULL);
2188 STYPE ("E-yt", is->efield_yt, NULL);
2189 STYPE ("E-z", is->efield_z, NULL);
2190 STYPE ("E-zt", is->efield_zt, NULL);
2192 CCTYPE("Ion/water position swapping for computational electrophysiology setups");
2193 CTYPE("Swap positions along direction: no, X, Y, Z");
2194 EETYPE("swapcoords", ir->eSwapCoords, eSwapTypes_names);
2195 if (ir->eSwapCoords != eswapNO)
2198 CTYPE("Swap attempt frequency");
2199 ITYPE("swap-frequency", ir->swap->nstswap, 1);
2200 CTYPE("Two index groups that contain the compartment-partitioning atoms");
2201 STYPE("split-group0", splitgrp0, NULL);
2202 STYPE("split-group1", splitgrp1, NULL);
2203 CTYPE("Use center of mass of split groups (yes/no), otherwise center of geometry is used");
2204 EETYPE("massw-split0", ir->swap->massw_split[0], yesno_names);
2205 EETYPE("massw-split1", ir->swap->massw_split[1], yesno_names);
2207 CTYPE("Group name of ions that can be exchanged with solvent molecules");
2208 STYPE("swap-group", swapgrp, NULL);
2209 CTYPE("Group name of solvent molecules");
2210 STYPE("solvent-group", solgrp, NULL);
2212 CTYPE("Split cylinder: radius, upper and lower extension (nm) (this will define the channels)");
2213 CTYPE("Note that the split cylinder settings do not have an influence on the swapping protocol,");
2214 CTYPE("however, if correctly defined, the ion permeation events are counted per channel");
2215 RTYPE("cyl0-r", ir->swap->cyl0r, 2.0);
2216 RTYPE("cyl0-up", ir->swap->cyl0u, 1.0);
2217 RTYPE("cyl0-down", ir->swap->cyl0l, 1.0);
2218 RTYPE("cyl1-r", ir->swap->cyl1r, 2.0);
2219 RTYPE("cyl1-up", ir->swap->cyl1u, 1.0);
2220 RTYPE("cyl1-down", ir->swap->cyl1l, 1.0);
2222 CTYPE("Average the number of ions per compartment over these many swap attempt steps");
2223 ITYPE("coupl-steps", ir->swap->nAverage, 10);
2224 CTYPE("Requested number of anions and cations for each of the two compartments");
2225 CTYPE("-1 means fix the numbers as found in time step 0");
2226 ITYPE("anionsA", ir->swap->nanions[0], -1);
2227 ITYPE("cationsA", ir->swap->ncations[0], -1);
2228 ITYPE("anionsB", ir->swap->nanions[1], -1);
2229 ITYPE("cationsB", ir->swap->ncations[1], -1);
2230 CTYPE("Start to swap ions if threshold difference to requested count is reached");
2231 RTYPE("threshold", ir->swap->threshold, 1.0);
2234 /* AdResS defined thingies */
2235 CCTYPE ("AdResS parameters");
2236 EETYPE("adress", ir->bAdress, yesno_names);
2239 snew(ir->adress, 1);
2240 read_adressparams(&ninp, &inp, ir->adress, wi);
2243 /* User defined thingies */
2244 CCTYPE ("User defined thingies");
2245 STYPE ("user1-grps", is->user1, NULL);
2246 STYPE ("user2-grps", is->user2, NULL);
2247 ITYPE ("userint1", ir->userint1, 0);
2248 ITYPE ("userint2", ir->userint2, 0);
2249 ITYPE ("userint3", ir->userint3, 0);
2250 ITYPE ("userint4", ir->userint4, 0);
2251 RTYPE ("userreal1", ir->userreal1, 0);
2252 RTYPE ("userreal2", ir->userreal2, 0);
2253 RTYPE ("userreal3", ir->userreal3, 0);
2254 RTYPE ("userreal4", ir->userreal4, 0);
2257 write_inpfile(mdparout, ninp, inp, FALSE, wi);
2258 for (i = 0; (i < ninp); i++)
2261 sfree(inp[i].value);
2265 /* Process options if necessary */
2266 for (m = 0; m < 2; m++)
2268 for (i = 0; i < 2*DIM; i++)
2277 if (sscanf(dumstr[m], "%lf", &(dumdub[m][XX])) != 1)
2279 warning_error(wi, "Pressure coupling not enough values (I need 1)");
2281 dumdub[m][YY] = dumdub[m][ZZ] = dumdub[m][XX];
2283 case epctSEMIISOTROPIC:
2284 case epctSURFACETENSION:
2285 if (sscanf(dumstr[m], "%lf%lf",
2286 &(dumdub[m][XX]), &(dumdub[m][ZZ])) != 2)
2288 warning_error(wi, "Pressure coupling not enough values (I need 2)");
2290 dumdub[m][YY] = dumdub[m][XX];
2292 case epctANISOTROPIC:
2293 if (sscanf(dumstr[m], "%lf%lf%lf%lf%lf%lf",
2294 &(dumdub[m][XX]), &(dumdub[m][YY]), &(dumdub[m][ZZ]),
2295 &(dumdub[m][3]), &(dumdub[m][4]), &(dumdub[m][5])) != 6)
2297 warning_error(wi, "Pressure coupling not enough values (I need 6)");
2301 gmx_fatal(FARGS, "Pressure coupling type %s not implemented yet",
2302 epcoupltype_names[ir->epct]);
2306 clear_mat(ir->ref_p);
2307 clear_mat(ir->compress);
2308 for (i = 0; i < DIM; i++)
2310 ir->ref_p[i][i] = dumdub[1][i];
2311 ir->compress[i][i] = dumdub[0][i];
2313 if (ir->epct == epctANISOTROPIC)
2315 ir->ref_p[XX][YY] = dumdub[1][3];
2316 ir->ref_p[XX][ZZ] = dumdub[1][4];
2317 ir->ref_p[YY][ZZ] = dumdub[1][5];
2318 if (ir->ref_p[XX][YY] != 0 && ir->ref_p[XX][ZZ] != 0 && ir->ref_p[YY][ZZ] != 0)
2320 warning(wi, "All off-diagonal reference pressures are non-zero. Are you sure you want to apply a threefold shear stress?\n");
2322 ir->compress[XX][YY] = dumdub[0][3];
2323 ir->compress[XX][ZZ] = dumdub[0][4];
2324 ir->compress[YY][ZZ] = dumdub[0][5];
2325 for (i = 0; i < DIM; i++)
2327 for (m = 0; m < i; m++)
2329 ir->ref_p[i][m] = ir->ref_p[m][i];
2330 ir->compress[i][m] = ir->compress[m][i];
2335 if (ir->comm_mode == ecmNO)
2340 opts->couple_moltype = NULL;
2341 if (strlen(is->couple_moltype) > 0)
2343 if (ir->efep != efepNO)
2345 opts->couple_moltype = gmx_strdup(is->couple_moltype);
2346 if (opts->couple_lam0 == opts->couple_lam1)
2348 warning(wi, "The lambda=0 and lambda=1 states for coupling are identical");
2350 if (ir->eI == eiMD && (opts->couple_lam0 == ecouplamNONE ||
2351 opts->couple_lam1 == ecouplamNONE))
2353 warning(wi, "For proper sampling of the (nearly) decoupled state, stochastic dynamics should be used");
2358 warning_note(wi, "Free energy is turned off, so we will not decouple the molecule listed in your input.");
2361 /* FREE ENERGY AND EXPANDED ENSEMBLE OPTIONS */
2362 if (ir->efep != efepNO)
2364 if (fep->delta_lambda > 0)
2366 ir->efep = efepSLOWGROWTH;
2370 if (fep->edHdLPrintEnergy == edHdLPrintEnergyYES)
2372 fep->edHdLPrintEnergy = edHdLPrintEnergyTOTAL;
2373 warning_note(wi, "Old option for dhdl-print-energy given: "
2374 "changing \"yes\" to \"total\"\n");
2377 if (ir->bSimTemp && (fep->edHdLPrintEnergy == edHdLPrintEnergyNO))
2379 /* always print out the energy to dhdl if we are doing
2380 expanded ensemble, since we need the total energy for
2381 analysis if the temperature is changing. In some
2382 conditions one may only want the potential energy, so
2383 we will allow that if the appropriate mdp setting has
2384 been enabled. Otherwise, total it is:
2386 fep->edHdLPrintEnergy = edHdLPrintEnergyTOTAL;
2389 if ((ir->efep != efepNO) || ir->bSimTemp)
2391 ir->bExpanded = FALSE;
2392 if ((ir->efep == efepEXPANDED) || ir->bSimTemp)
2394 ir->bExpanded = TRUE;
2396 do_fep_params(ir, is->fep_lambda, is->lambda_weights);
2397 if (ir->bSimTemp) /* done after fep params */
2399 do_simtemp_params(ir);
2404 ir->fepvals->n_lambda = 0;
2407 /* WALL PARAMETERS */
2409 do_wall_params(ir, is->wall_atomtype, is->wall_density, opts);
2411 /* ORIENTATION RESTRAINT PARAMETERS */
2413 if (opts->bOrire && str_nelem(is->orirefitgrp, MAXPTR, NULL) != 1)
2415 warning_error(wi, "ERROR: Need one orientation restraint fit group\n");
2418 /* DEFORMATION PARAMETERS */
2420 clear_mat(ir->deform);
2421 for (i = 0; i < 6; i++)
2425 m = sscanf(is->deform, "%lf %lf %lf %lf %lf %lf",
2426 &(dumdub[0][0]), &(dumdub[0][1]), &(dumdub[0][2]),
2427 &(dumdub[0][3]), &(dumdub[0][4]), &(dumdub[0][5]));
2428 for (i = 0; i < 3; i++)
2430 ir->deform[i][i] = dumdub[0][i];
2432 ir->deform[YY][XX] = dumdub[0][3];
2433 ir->deform[ZZ][XX] = dumdub[0][4];
2434 ir->deform[ZZ][YY] = dumdub[0][5];
2435 if (ir->epc != epcNO)
2437 for (i = 0; i < 3; i++)
2439 for (j = 0; j <= i; j++)
2441 if (ir->deform[i][j] != 0 && ir->compress[i][j] != 0)
2443 warning_error(wi, "A box element has deform set and compressibility > 0");
2447 for (i = 0; i < 3; i++)
2449 for (j = 0; j < i; j++)
2451 if (ir->deform[i][j] != 0)
2453 for (m = j; m < DIM; m++)
2455 if (ir->compress[m][j] != 0)
2457 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.");
2458 warning(wi, warn_buf);
2466 /* Ion/water position swapping checks */
2467 if (ir->eSwapCoords != eswapNO)
2469 if (ir->swap->nstswap < 1)
2471 warning_error(wi, "swap_frequency must be 1 or larger when ion swapping is requested");
2473 if (ir->swap->nAverage < 1)
2475 warning_error(wi, "coupl_steps must be 1 or larger.\n");
2477 if (ir->swap->threshold < 1.0)
2479 warning_error(wi, "Ion count threshold must be at least 1.\n");
2487 static int search_QMstring(const char *s, int ng, const char *gn[])
2489 /* same as normal search_string, but this one searches QM strings */
2492 for (i = 0; (i < ng); i++)
2494 if (gmx_strcasecmp(s, gn[i]) == 0)
2500 gmx_fatal(FARGS, "this QM method or basisset (%s) is not implemented\n!", s);
2504 } /* search_QMstring */
2506 /* We would like gn to be const as well, but C doesn't allow this */
2507 int search_string(const char *s, int ng, char *gn[])
2511 for (i = 0; (i < ng); i++)
2513 if (gmx_strcasecmp(s, gn[i]) == 0)
2520 "Group %s referenced in the .mdp file was not found in the index file.\n"
2521 "Group names must match either [moleculetype] names or custom index group\n"
2522 "names, in which case you must supply an index file to the '-n' option\n"
2529 static gmx_bool do_numbering(int natoms, gmx_groups_t *groups, int ng, char *ptrs[],
2530 t_blocka *block, char *gnames[],
2531 int gtype, int restnm,
2532 int grptp, gmx_bool bVerbose,
2535 unsigned short *cbuf;
2536 t_grps *grps = &(groups->grps[gtype]);
2537 int i, j, gid, aj, ognr, ntot = 0;
2540 char warn_buf[STRLEN];
2544 fprintf(debug, "Starting numbering %d groups of type %d\n", ng, gtype);
2547 title = gtypes[gtype];
2550 /* Mark all id's as not set */
2551 for (i = 0; (i < natoms); i++)
2556 snew(grps->nm_ind, ng+1); /* +1 for possible rest group */
2557 for (i = 0; (i < ng); i++)
2559 /* Lookup the group name in the block structure */
2560 gid = search_string(ptrs[i], block->nr, gnames);
2561 if ((grptp != egrptpONE) || (i == 0))
2563 grps->nm_ind[grps->nr++] = gid;
2567 fprintf(debug, "Found gid %d for group %s\n", gid, ptrs[i]);
2570 /* Now go over the atoms in the group */
2571 for (j = block->index[gid]; (j < block->index[gid+1]); j++)
2576 /* Range checking */
2577 if ((aj < 0) || (aj >= natoms))
2579 gmx_fatal(FARGS, "Invalid atom number %d in indexfile", aj);
2581 /* Lookup up the old group number */
2585 gmx_fatal(FARGS, "Atom %d in multiple %s groups (%d and %d)",
2586 aj+1, title, ognr+1, i+1);
2590 /* Store the group number in buffer */
2591 if (grptp == egrptpONE)
2604 /* Now check whether we have done all atoms */
2608 if (grptp == egrptpALL)
2610 gmx_fatal(FARGS, "%d atoms are not part of any of the %s groups",
2611 natoms-ntot, title);
2613 else if (grptp == egrptpPART)
2615 sprintf(warn_buf, "%d atoms are not part of any of the %s groups",
2616 natoms-ntot, title);
2617 warning_note(wi, warn_buf);
2619 /* Assign all atoms currently unassigned to a rest group */
2620 for (j = 0; (j < natoms); j++)
2622 if (cbuf[j] == NOGID)
2628 if (grptp != egrptpPART)
2633 "Making dummy/rest group for %s containing %d elements\n",
2634 title, natoms-ntot);
2636 /* Add group name "rest" */
2637 grps->nm_ind[grps->nr] = restnm;
2639 /* Assign the rest name to all atoms not currently assigned to a group */
2640 for (j = 0; (j < natoms); j++)
2642 if (cbuf[j] == NOGID)
2651 if (grps->nr == 1 && (ntot == 0 || ntot == natoms))
2653 /* All atoms are part of one (or no) group, no index required */
2654 groups->ngrpnr[gtype] = 0;
2655 groups->grpnr[gtype] = NULL;
2659 groups->ngrpnr[gtype] = natoms;
2660 snew(groups->grpnr[gtype], natoms);
2661 for (j = 0; (j < natoms); j++)
2663 groups->grpnr[gtype][j] = cbuf[j];
2669 return (bRest && grptp == egrptpPART);
2672 static void calc_nrdf(gmx_mtop_t *mtop, t_inputrec *ir, char **gnames)
2675 gmx_groups_t *groups;
2677 int natoms, ai, aj, i, j, d, g, imin, jmin;
2679 int *nrdf2, *na_vcm, na_tot;
2680 double *nrdf_tc, *nrdf_vcm, nrdf_uc, n_sub = 0;
2681 gmx_mtop_atomloop_all_t aloop;
2683 int mb, mol, ftype, as;
2684 gmx_molblock_t *molb;
2685 gmx_moltype_t *molt;
2688 * First calc 3xnr-atoms for each group
2689 * then subtract half a degree of freedom for each constraint
2691 * Only atoms and nuclei contribute to the degrees of freedom...
2696 groups = &mtop->groups;
2697 natoms = mtop->natoms;
2699 /* Allocate one more for a possible rest group */
2700 /* We need to sum degrees of freedom into doubles,
2701 * since floats give too low nrdf's above 3 million atoms.
2703 snew(nrdf_tc, groups->grps[egcTC].nr+1);
2704 snew(nrdf_vcm, groups->grps[egcVCM].nr+1);
2705 snew(na_vcm, groups->grps[egcVCM].nr+1);
2707 for (i = 0; i < groups->grps[egcTC].nr; i++)
2711 for (i = 0; i < groups->grps[egcVCM].nr+1; i++)
2716 snew(nrdf2, natoms);
2717 aloop = gmx_mtop_atomloop_all_init(mtop);
2718 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
2721 if (atom->ptype == eptAtom || atom->ptype == eptNucleus)
2723 g = ggrpnr(groups, egcFREEZE, i);
2724 /* Double count nrdf for particle i */
2725 for (d = 0; d < DIM; d++)
2727 if (opts->nFreeze[g][d] == 0)
2732 nrdf_tc [ggrpnr(groups, egcTC, i)] += 0.5*nrdf2[i];
2733 nrdf_vcm[ggrpnr(groups, egcVCM, i)] += 0.5*nrdf2[i];
2738 for (mb = 0; mb < mtop->nmolblock; mb++)
2740 molb = &mtop->molblock[mb];
2741 molt = &mtop->moltype[molb->type];
2742 atom = molt->atoms.atom;
2743 for (mol = 0; mol < molb->nmol; mol++)
2745 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
2747 ia = molt->ilist[ftype].iatoms;
2748 for (i = 0; i < molt->ilist[ftype].nr; )
2750 /* Subtract degrees of freedom for the constraints,
2751 * if the particles still have degrees of freedom left.
2752 * If one of the particles is a vsite or a shell, then all
2753 * constraint motion will go there, but since they do not
2754 * contribute to the constraints the degrees of freedom do not
2759 if (((atom[ia[1]].ptype == eptNucleus) ||
2760 (atom[ia[1]].ptype == eptAtom)) &&
2761 ((atom[ia[2]].ptype == eptNucleus) ||
2762 (atom[ia[2]].ptype == eptAtom)))
2780 imin = min(imin, nrdf2[ai]);
2781 jmin = min(jmin, nrdf2[aj]);
2784 nrdf_tc [ggrpnr(groups, egcTC, ai)] -= 0.5*imin;
2785 nrdf_tc [ggrpnr(groups, egcTC, aj)] -= 0.5*jmin;
2786 nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2787 nrdf_vcm[ggrpnr(groups, egcVCM, aj)] -= 0.5*jmin;
2789 ia += interaction_function[ftype].nratoms+1;
2790 i += interaction_function[ftype].nratoms+1;
2793 ia = molt->ilist[F_SETTLE].iatoms;
2794 for (i = 0; i < molt->ilist[F_SETTLE].nr; )
2796 /* Subtract 1 dof from every atom in the SETTLE */
2797 for (j = 0; j < 3; j++)
2800 imin = min(2, nrdf2[ai]);
2802 nrdf_tc [ggrpnr(groups, egcTC, ai)] -= 0.5*imin;
2803 nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2808 as += molt->atoms.nr;
2812 if (ir->ePull == epullCONSTRAINT)
2814 /* Correct nrdf for the COM constraints.
2815 * We correct using the TC and VCM group of the first atom
2816 * in the reference and pull group. If atoms in one pull group
2817 * belong to different TC or VCM groups it is anyhow difficult
2818 * to determine the optimal nrdf assignment.
2822 for (i = 0; i < pull->ncoord; i++)
2826 for (j = 0; j < 2; j++)
2828 const t_pull_group *pgrp;
2830 pgrp = &pull->group[pull->coord[i].group[j]];
2834 /* Subtract 1/2 dof from each group */
2836 nrdf_tc [ggrpnr(groups, egcTC, ai)] -= 0.5*imin;
2837 nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2838 if (nrdf_tc[ggrpnr(groups, egcTC, ai)] < 0)
2840 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)]]);
2845 /* We need to subtract the whole DOF from group j=1 */
2852 if (ir->nstcomm != 0)
2854 /* Subtract 3 from the number of degrees of freedom in each vcm group
2855 * when com translation is removed and 6 when rotation is removed
2858 switch (ir->comm_mode)
2861 n_sub = ndof_com(ir);
2868 gmx_incons("Checking comm_mode");
2871 for (i = 0; i < groups->grps[egcTC].nr; i++)
2873 /* Count the number of atoms of TC group i for every VCM group */
2874 for (j = 0; j < groups->grps[egcVCM].nr+1; j++)
2879 for (ai = 0; ai < natoms; ai++)
2881 if (ggrpnr(groups, egcTC, ai) == i)
2883 na_vcm[ggrpnr(groups, egcVCM, ai)]++;
2887 /* Correct for VCM removal according to the fraction of each VCM
2888 * group present in this TC group.
2890 nrdf_uc = nrdf_tc[i];
2893 fprintf(debug, "T-group[%d] nrdf_uc = %g, n_sub = %g\n",
2897 for (j = 0; j < groups->grps[egcVCM].nr+1; j++)
2899 if (nrdf_vcm[j] > n_sub)
2901 nrdf_tc[i] += nrdf_uc*((double)na_vcm[j]/(double)na_tot)*
2902 (nrdf_vcm[j] - n_sub)/nrdf_vcm[j];
2906 fprintf(debug, " nrdf_vcm[%d] = %g, nrdf = %g\n",
2907 j, nrdf_vcm[j], nrdf_tc[i]);
2912 for (i = 0; (i < groups->grps[egcTC].nr); i++)
2914 opts->nrdf[i] = nrdf_tc[i];
2915 if (opts->nrdf[i] < 0)
2920 "Number of degrees of freedom in T-Coupling group %s is %.2f\n",
2921 gnames[groups->grps[egcTC].nm_ind[i]], opts->nrdf[i]);
2930 static void decode_cos(char *s, t_cosines *cosine)
2933 char format[STRLEN], f1[STRLEN];
2945 sscanf(t, "%d", &(cosine->n));
2952 snew(cosine->a, cosine->n);
2953 snew(cosine->phi, cosine->n);
2955 sprintf(format, "%%*d");
2956 for (i = 0; (i < cosine->n); i++)
2959 strcat(f1, "%lf%lf");
2960 if (sscanf(t, f1, &a, &phi) < 2)
2962 gmx_fatal(FARGS, "Invalid input for electric field shift: '%s'", t);
2965 cosine->phi[i] = phi;
2966 strcat(format, "%*lf%*lf");
2973 static gmx_bool do_egp_flag(t_inputrec *ir, gmx_groups_t *groups,
2974 const char *option, const char *val, int flag)
2976 /* The maximum number of energy group pairs would be MAXPTR*(MAXPTR+1)/2.
2977 * But since this is much larger than STRLEN, such a line can not be parsed.
2978 * The real maximum is the number of names that fit in a string: STRLEN/2.
2980 #define EGP_MAX (STRLEN/2)
2981 int nelem, i, j, k, nr;
2982 char *names[EGP_MAX];
2986 gnames = groups->grpname;
2988 nelem = str_nelem(val, EGP_MAX, names);
2991 gmx_fatal(FARGS, "The number of groups for %s is odd", option);
2993 nr = groups->grps[egcENER].nr;
2995 for (i = 0; i < nelem/2; i++)
2999 gmx_strcasecmp(names[2*i], *(gnames[groups->grps[egcENER].nm_ind[j]])))
3005 gmx_fatal(FARGS, "%s in %s is not an energy group\n",
3006 names[2*i], option);
3010 gmx_strcasecmp(names[2*i+1], *(gnames[groups->grps[egcENER].nm_ind[k]])))
3016 gmx_fatal(FARGS, "%s in %s is not an energy group\n",
3017 names[2*i+1], option);
3019 if ((j < nr) && (k < nr))
3021 ir->opts.egp_flags[nr*j+k] |= flag;
3022 ir->opts.egp_flags[nr*k+j] |= flag;
3031 static void make_swap_groups(
3040 int ig = -1, i = 0, j;
3044 /* Just a quick check here, more thorough checks are in mdrun */
3045 if (strcmp(splitg0name, splitg1name) == 0)
3047 gmx_fatal(FARGS, "The split groups can not both be '%s'.", splitg0name);
3050 /* First get the swap group index atoms */
3051 ig = search_string(swapgname, grps->nr, gnames);
3052 swap->nat = grps->index[ig+1] - grps->index[ig];
3055 fprintf(stderr, "Swap group '%s' contains %d atoms.\n", swapgname, swap->nat);
3056 snew(swap->ind, swap->nat);
3057 for (i = 0; i < swap->nat; i++)
3059 swap->ind[i] = grps->a[grps->index[ig]+i];
3064 gmx_fatal(FARGS, "You defined an empty group of atoms for swapping.");
3067 /* Now do so for the split groups */
3068 for (j = 0; j < 2; j++)
3072 splitg = splitg0name;
3076 splitg = splitg1name;
3079 ig = search_string(splitg, grps->nr, gnames);
3080 swap->nat_split[j] = grps->index[ig+1] - grps->index[ig];
3081 if (swap->nat_split[j] > 0)
3083 fprintf(stderr, "Split group %d '%s' contains %d atom%s.\n",
3084 j, splitg, swap->nat_split[j], (swap->nat_split[j] > 1) ? "s" : "");
3085 snew(swap->ind_split[j], swap->nat_split[j]);
3086 for (i = 0; i < swap->nat_split[j]; i++)
3088 swap->ind_split[j][i] = grps->a[grps->index[ig]+i];
3093 gmx_fatal(FARGS, "Split group %d has to contain at least 1 atom!", j);
3097 /* Now get the solvent group index atoms */
3098 ig = search_string(solgname, grps->nr, gnames);
3099 swap->nat_sol = grps->index[ig+1] - grps->index[ig];
3100 if (swap->nat_sol > 0)
3102 fprintf(stderr, "Solvent group '%s' contains %d atoms.\n", solgname, swap->nat_sol);
3103 snew(swap->ind_sol, swap->nat_sol);
3104 for (i = 0; i < swap->nat_sol; i++)
3106 swap->ind_sol[i] = grps->a[grps->index[ig]+i];
3111 gmx_fatal(FARGS, "You defined an empty group of solvent. Cannot exchange ions.");
3116 void make_IMD_group(t_IMD *IMDgroup, char *IMDgname, t_blocka *grps, char **gnames)
3121 ig = search_string(IMDgname, grps->nr, gnames);
3122 IMDgroup->nat = grps->index[ig+1] - grps->index[ig];
3124 if (IMDgroup->nat > 0)
3126 fprintf(stderr, "Group '%s' with %d atoms can be activated for interactive molecular dynamics (IMD).\n",
3127 IMDgname, IMDgroup->nat);
3128 snew(IMDgroup->ind, IMDgroup->nat);
3129 for (i = 0; i < IMDgroup->nat; i++)
3131 IMDgroup->ind[i] = grps->a[grps->index[ig]+i];
3137 void do_index(const char* mdparin, const char *ndx,
3140 t_inputrec *ir, rvec *v,
3144 gmx_groups_t *groups;
3148 char warnbuf[STRLEN], **gnames;
3149 int nr, ntcg, ntau_t, nref_t, nacc, nofg, nSA, nSA_points, nSA_time, nSA_temp;
3152 int nacg, nfreeze, nfrdim, nenergy, nvcm, nuser;
3153 char *ptr1[MAXPTR], *ptr2[MAXPTR], *ptr3[MAXPTR];
3154 int i, j, k, restnm;
3156 gmx_bool bExcl, bTable, bSetTCpar, bAnneal, bRest;
3157 int nQMmethod, nQMbasis, nQMcharge, nQMmult, nbSH, nCASorb, nCASelec,
3158 nSAon, nSAoff, nSAsteps, nQMg, nbOPT, nbTS;
3159 char warn_buf[STRLEN];
3163 fprintf(stderr, "processing index file...\n");
3169 snew(grps->index, 1);
3171 atoms_all = gmx_mtop_global_atoms(mtop);
3172 analyse(&atoms_all, grps, &gnames, FALSE, TRUE);
3173 free_t_atoms(&atoms_all, FALSE);
3177 grps = init_index(ndx, &gnames);
3180 groups = &mtop->groups;
3181 natoms = mtop->natoms;
3182 symtab = &mtop->symtab;
3184 snew(groups->grpname, grps->nr+1);
3186 for (i = 0; (i < grps->nr); i++)
3188 groups->grpname[i] = put_symtab(symtab, gnames[i]);
3190 groups->grpname[i] = put_symtab(symtab, "rest");
3192 srenew(gnames, grps->nr+1);
3193 gnames[restnm] = *(groups->grpname[i]);
3194 groups->ngrpname = grps->nr+1;
3196 set_warning_line(wi, mdparin, -1);
3198 ntau_t = str_nelem(is->tau_t, MAXPTR, ptr1);
3199 nref_t = str_nelem(is->ref_t, MAXPTR, ptr2);
3200 ntcg = str_nelem(is->tcgrps, MAXPTR, ptr3);
3201 if ((ntau_t != ntcg) || (nref_t != ntcg))
3203 gmx_fatal(FARGS, "Invalid T coupling input: %d groups, %d ref-t values and "
3204 "%d tau-t values", ntcg, nref_t, ntau_t);
3207 bSetTCpar = (ir->etc || EI_SD(ir->eI) || ir->eI == eiBD || EI_TPI(ir->eI));
3208 do_numbering(natoms, groups, ntcg, ptr3, grps, gnames, egcTC,
3209 restnm, bSetTCpar ? egrptpALL : egrptpALL_GENREST, bVerbose, wi);
3210 nr = groups->grps[egcTC].nr;
3212 snew(ir->opts.nrdf, nr);
3213 snew(ir->opts.tau_t, nr);
3214 snew(ir->opts.ref_t, nr);
3215 if (ir->eI == eiBD && ir->bd_fric == 0)
3217 fprintf(stderr, "bd-fric=0, so tau-t will be used as the inverse friction constant(s)\n");
3224 gmx_fatal(FARGS, "Not enough ref-t and tau-t values!");
3228 for (i = 0; (i < nr); i++)
3230 ir->opts.tau_t[i] = strtod(ptr1[i], NULL);
3231 if ((ir->eI == eiBD || ir->eI == eiSD2) && ir->opts.tau_t[i] <= 0)
3233 sprintf(warn_buf, "With integrator %s tau-t should be larger than 0", ei_names[ir->eI]);
3234 warning_error(wi, warn_buf);
3237 if (ir->etc != etcVRESCALE && ir->opts.tau_t[i] == 0)
3239 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.");
3242 if (ir->opts.tau_t[i] >= 0)
3244 tau_min = min(tau_min, ir->opts.tau_t[i]);
3247 if (ir->etc != etcNO && ir->nsttcouple == -1)
3249 ir->nsttcouple = ir_optimal_nsttcouple(ir);
3254 if ((ir->etc == etcNOSEHOOVER) && (ir->epc == epcBERENDSEN))
3256 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");
3258 if ((ir->epc == epcMTTK) && (ir->etc > etcNO))
3260 if (ir->nstpcouple != ir->nsttcouple)
3262 int mincouple = min(ir->nstpcouple, ir->nsttcouple);
3263 ir->nstpcouple = ir->nsttcouple = mincouple;
3264 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);
3265 warning_note(wi, warn_buf);
3269 /* velocity verlet with averaged kinetic energy KE = 0.5*(v(t+1/2) - v(t-1/2)) is implemented
3270 primarily for testing purposes, and does not work with temperature coupling other than 1 */
3272 if (ETC_ANDERSEN(ir->etc))
3274 if (ir->nsttcouple != 1)
3277 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");
3278 warning_note(wi, warn_buf);
3281 nstcmin = tcouple_min_integration_steps(ir->etc);
3284 if (tau_min/(ir->delta_t*ir->nsttcouple) < nstcmin)
3286 sprintf(warn_buf, "For proper integration of the %s thermostat, tau-t (%g) should be at least %d times larger than nsttcouple*dt (%g)",
3287 ETCOUPLTYPE(ir->etc),
3289 ir->nsttcouple*ir->delta_t);
3290 warning(wi, warn_buf);
3293 for (i = 0; (i < nr); i++)
3295 ir->opts.ref_t[i] = strtod(ptr2[i], NULL);
3296 if (ir->opts.ref_t[i] < 0)
3298 gmx_fatal(FARGS, "ref-t for group %d negative", i);
3301 /* set the lambda mc temperature to the md integrator temperature (which should be defined
3302 if we are in this conditional) if mc_temp is negative */
3303 if (ir->expandedvals->mc_temp < 0)
3305 ir->expandedvals->mc_temp = ir->opts.ref_t[0]; /*for now, set to the first reft */
3309 /* Simulated annealing for each group. There are nr groups */
3310 nSA = str_nelem(is->anneal, MAXPTR, ptr1);
3311 if (nSA == 1 && (ptr1[0][0] == 'n' || ptr1[0][0] == 'N'))
3315 if (nSA > 0 && nSA != nr)
3317 gmx_fatal(FARGS, "Not enough annealing values: %d (for %d groups)\n", nSA, nr);
3321 snew(ir->opts.annealing, nr);
3322 snew(ir->opts.anneal_npoints, nr);
3323 snew(ir->opts.anneal_time, nr);
3324 snew(ir->opts.anneal_temp, nr);
3325 for (i = 0; i < nr; i++)
3327 ir->opts.annealing[i] = eannNO;
3328 ir->opts.anneal_npoints[i] = 0;
3329 ir->opts.anneal_time[i] = NULL;
3330 ir->opts.anneal_temp[i] = NULL;
3335 for (i = 0; i < nr; i++)
3337 if (ptr1[i][0] == 'n' || ptr1[i][0] == 'N')
3339 ir->opts.annealing[i] = eannNO;
3341 else if (ptr1[i][0] == 's' || ptr1[i][0] == 'S')
3343 ir->opts.annealing[i] = eannSINGLE;
3346 else if (ptr1[i][0] == 'p' || ptr1[i][0] == 'P')
3348 ir->opts.annealing[i] = eannPERIODIC;
3354 /* Read the other fields too */
3355 nSA_points = str_nelem(is->anneal_npoints, MAXPTR, ptr1);
3356 if (nSA_points != nSA)
3358 gmx_fatal(FARGS, "Found %d annealing-npoints values for %d groups\n", nSA_points, nSA);
3360 for (k = 0, i = 0; i < nr; i++)
3362 ir->opts.anneal_npoints[i] = strtol(ptr1[i], NULL, 10);
3363 if (ir->opts.anneal_npoints[i] == 1)
3365 gmx_fatal(FARGS, "Please specify at least a start and an end point for annealing\n");
3367 snew(ir->opts.anneal_time[i], ir->opts.anneal_npoints[i]);
3368 snew(ir->opts.anneal_temp[i], ir->opts.anneal_npoints[i]);
3369 k += ir->opts.anneal_npoints[i];
3372 nSA_time = str_nelem(is->anneal_time, MAXPTR, ptr1);
3375 gmx_fatal(FARGS, "Found %d annealing-time values, wanter %d\n", nSA_time, k);
3377 nSA_temp = str_nelem(is->anneal_temp, MAXPTR, ptr2);
3380 gmx_fatal(FARGS, "Found %d annealing-temp values, wanted %d\n", nSA_temp, k);
3383 for (i = 0, k = 0; i < nr; i++)
3386 for (j = 0; j < ir->opts.anneal_npoints[i]; j++)
3388 ir->opts.anneal_time[i][j] = strtod(ptr1[k], NULL);
3389 ir->opts.anneal_temp[i][j] = strtod(ptr2[k], NULL);
3392 if (ir->opts.anneal_time[i][0] > (ir->init_t+GMX_REAL_EPS))
3394 gmx_fatal(FARGS, "First time point for annealing > init_t.\n");
3400 if (ir->opts.anneal_time[i][j] < ir->opts.anneal_time[i][j-1])
3402 gmx_fatal(FARGS, "Annealing timepoints out of order: t=%f comes after t=%f\n",
3403 ir->opts.anneal_time[i][j], ir->opts.anneal_time[i][j-1]);
3406 if (ir->opts.anneal_temp[i][j] < 0)
3408 gmx_fatal(FARGS, "Found negative temperature in annealing: %f\n", ir->opts.anneal_temp[i][j]);
3413 /* Print out some summary information, to make sure we got it right */
3414 for (i = 0, k = 0; i < nr; i++)
3416 if (ir->opts.annealing[i] != eannNO)
3418 j = groups->grps[egcTC].nm_ind[i];
3419 fprintf(stderr, "Simulated annealing for group %s: %s, %d timepoints\n",
3420 *(groups->grpname[j]), eann_names[ir->opts.annealing[i]],
3421 ir->opts.anneal_npoints[i]);
3422 fprintf(stderr, "Time (ps) Temperature (K)\n");
3423 /* All terms except the last one */
3424 for (j = 0; j < (ir->opts.anneal_npoints[i]-1); j++)
3426 fprintf(stderr, "%9.1f %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3429 /* Finally the last one */
3430 j = ir->opts.anneal_npoints[i]-1;
3431 if (ir->opts.annealing[i] == eannSINGLE)
3433 fprintf(stderr, "%9.1f- %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3437 fprintf(stderr, "%9.1f %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3438 if (fabs(ir->opts.anneal_temp[i][j]-ir->opts.anneal_temp[i][0]) > GMX_REAL_EPS)
3440 warning_note(wi, "There is a temperature jump when your annealing loops back.\n");
3449 if (ir->ePull != epullNO)
3451 make_pull_groups(ir->pull, is->pull_grp, grps, gnames);
3453 make_pull_coords(ir->pull);
3458 make_rotation_groups(ir->rot, is->rot_grp, grps, gnames);
3461 if (ir->eSwapCoords != eswapNO)
3463 make_swap_groups(ir->swap, swapgrp, splitgrp0, splitgrp1, solgrp, grps, gnames);
3466 /* Make indices for IMD session */
3469 make_IMD_group(ir->imd, is->imd_grp, grps, gnames);
3472 nacc = str_nelem(is->acc, MAXPTR, ptr1);
3473 nacg = str_nelem(is->accgrps, MAXPTR, ptr2);
3474 if (nacg*DIM != nacc)
3476 gmx_fatal(FARGS, "Invalid Acceleration input: %d groups and %d acc. values",
3479 do_numbering(natoms, groups, nacg, ptr2, grps, gnames, egcACC,
3480 restnm, egrptpALL_GENREST, bVerbose, wi);
3481 nr = groups->grps[egcACC].nr;
3482 snew(ir->opts.acc, nr);
3483 ir->opts.ngacc = nr;
3485 for (i = k = 0; (i < nacg); i++)
3487 for (j = 0; (j < DIM); j++, k++)
3489 ir->opts.acc[i][j] = strtod(ptr1[k], NULL);
3492 for (; (i < nr); i++)
3494 for (j = 0; (j < DIM); j++)
3496 ir->opts.acc[i][j] = 0;
3500 nfrdim = str_nelem(is->frdim, MAXPTR, ptr1);
3501 nfreeze = str_nelem(is->freeze, MAXPTR, ptr2);
3502 if (nfrdim != DIM*nfreeze)
3504 gmx_fatal(FARGS, "Invalid Freezing input: %d groups and %d freeze values",
3507 do_numbering(natoms, groups, nfreeze, ptr2, grps, gnames, egcFREEZE,
3508 restnm, egrptpALL_GENREST, bVerbose, wi);
3509 nr = groups->grps[egcFREEZE].nr;
3510 ir->opts.ngfrz = nr;
3511 snew(ir->opts.nFreeze, nr);
3512 for (i = k = 0; (i < nfreeze); i++)
3514 for (j = 0; (j < DIM); j++, k++)
3516 ir->opts.nFreeze[i][j] = (gmx_strncasecmp(ptr1[k], "Y", 1) == 0);
3517 if (!ir->opts.nFreeze[i][j])
3519 if (gmx_strncasecmp(ptr1[k], "N", 1) != 0)
3521 sprintf(warnbuf, "Please use Y(ES) or N(O) for freezedim only "
3522 "(not %s)", ptr1[k]);
3523 warning(wi, warn_buf);
3528 for (; (i < nr); i++)
3530 for (j = 0; (j < DIM); j++)
3532 ir->opts.nFreeze[i][j] = 0;
3536 nenergy = str_nelem(is->energy, MAXPTR, ptr1);
3537 do_numbering(natoms, groups, nenergy, ptr1, grps, gnames, egcENER,
3538 restnm, egrptpALL_GENREST, bVerbose, wi);
3539 add_wall_energrps(groups, ir->nwall, symtab);
3540 ir->opts.ngener = groups->grps[egcENER].nr;
3541 nvcm = str_nelem(is->vcm, MAXPTR, ptr1);
3543 do_numbering(natoms, groups, nvcm, ptr1, grps, gnames, egcVCM,
3544 restnm, nvcm == 0 ? egrptpALL_GENREST : egrptpPART, bVerbose, wi);
3547 warning(wi, "Some atoms are not part of any center of mass motion removal group.\n"
3548 "This may lead to artifacts.\n"
3549 "In most cases one should use one group for the whole system.");
3552 /* Now we have filled the freeze struct, so we can calculate NRDF */
3553 calc_nrdf(mtop, ir, gnames);
3559 /* Must check per group! */
3560 for (i = 0; (i < ir->opts.ngtc); i++)
3562 ntot += ir->opts.nrdf[i];
3564 if (ntot != (DIM*natoms))
3566 fac = sqrt(ntot/(DIM*natoms));
3569 fprintf(stderr, "Scaling velocities by a factor of %.3f to account for constraints\n"
3570 "and removal of center of mass motion\n", fac);
3572 for (i = 0; (i < natoms); i++)
3574 svmul(fac, v[i], v[i]);
3579 nuser = str_nelem(is->user1, MAXPTR, ptr1);
3580 do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcUser1,
3581 restnm, egrptpALL_GENREST, bVerbose, wi);
3582 nuser = str_nelem(is->user2, MAXPTR, ptr1);
3583 do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcUser2,
3584 restnm, egrptpALL_GENREST, bVerbose, wi);
3585 nuser = str_nelem(is->x_compressed_groups, MAXPTR, ptr1);
3586 do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcCompressedX,
3587 restnm, egrptpONE, bVerbose, wi);
3588 nofg = str_nelem(is->orirefitgrp, MAXPTR, ptr1);
3589 do_numbering(natoms, groups, nofg, ptr1, grps, gnames, egcORFIT,
3590 restnm, egrptpALL_GENREST, bVerbose, wi);
3592 /* QMMM input processing */
3593 nQMg = str_nelem(is->QMMM, MAXPTR, ptr1);
3594 nQMmethod = str_nelem(is->QMmethod, MAXPTR, ptr2);
3595 nQMbasis = str_nelem(is->QMbasis, MAXPTR, ptr3);
3596 if ((nQMmethod != nQMg) || (nQMbasis != nQMg))
3598 gmx_fatal(FARGS, "Invalid QMMM input: %d groups %d basissets"
3599 " and %d methods\n", nQMg, nQMbasis, nQMmethod);
3601 /* group rest, if any, is always MM! */
3602 do_numbering(natoms, groups, nQMg, ptr1, grps, gnames, egcQMMM,
3603 restnm, egrptpALL_GENREST, bVerbose, wi);
3604 nr = nQMg; /*atoms->grps[egcQMMM].nr;*/
3605 ir->opts.ngQM = nQMg;
3606 snew(ir->opts.QMmethod, nr);
3607 snew(ir->opts.QMbasis, nr);
3608 for (i = 0; i < nr; i++)
3610 /* input consists of strings: RHF CASSCF PM3 .. These need to be
3611 * converted to the corresponding enum in names.c
3613 ir->opts.QMmethod[i] = search_QMstring(ptr2[i], eQMmethodNR,
3615 ir->opts.QMbasis[i] = search_QMstring(ptr3[i], eQMbasisNR,
3619 nQMmult = str_nelem(is->QMmult, MAXPTR, ptr1);
3620 nQMcharge = str_nelem(is->QMcharge, MAXPTR, ptr2);
3621 nbSH = str_nelem(is->bSH, MAXPTR, ptr3);
3622 snew(ir->opts.QMmult, nr);
3623 snew(ir->opts.QMcharge, nr);
3624 snew(ir->opts.bSH, nr);
3626 for (i = 0; i < nr; i++)
3628 ir->opts.QMmult[i] = strtol(ptr1[i], NULL, 10);
3629 ir->opts.QMcharge[i] = strtol(ptr2[i], NULL, 10);
3630 ir->opts.bSH[i] = (gmx_strncasecmp(ptr3[i], "Y", 1) == 0);
3633 nCASelec = str_nelem(is->CASelectrons, MAXPTR, ptr1);
3634 nCASorb = str_nelem(is->CASorbitals, MAXPTR, ptr2);
3635 snew(ir->opts.CASelectrons, nr);
3636 snew(ir->opts.CASorbitals, nr);
3637 for (i = 0; i < nr; i++)
3639 ir->opts.CASelectrons[i] = strtol(ptr1[i], NULL, 10);
3640 ir->opts.CASorbitals[i] = strtol(ptr2[i], NULL, 10);
3642 /* special optimization options */
3644 nbOPT = str_nelem(is->bOPT, MAXPTR, ptr1);
3645 nbTS = str_nelem(is->bTS, MAXPTR, ptr2);
3646 snew(ir->opts.bOPT, nr);
3647 snew(ir->opts.bTS, nr);
3648 for (i = 0; i < nr; i++)
3650 ir->opts.bOPT[i] = (gmx_strncasecmp(ptr1[i], "Y", 1) == 0);
3651 ir->opts.bTS[i] = (gmx_strncasecmp(ptr2[i], "Y", 1) == 0);
3653 nSAon = str_nelem(is->SAon, MAXPTR, ptr1);
3654 nSAoff = str_nelem(is->SAoff, MAXPTR, ptr2);
3655 nSAsteps = str_nelem(is->SAsteps, MAXPTR, ptr3);
3656 snew(ir->opts.SAon, nr);
3657 snew(ir->opts.SAoff, nr);
3658 snew(ir->opts.SAsteps, nr);
3660 for (i = 0; i < nr; i++)
3662 ir->opts.SAon[i] = strtod(ptr1[i], NULL);
3663 ir->opts.SAoff[i] = strtod(ptr2[i], NULL);
3664 ir->opts.SAsteps[i] = strtol(ptr3[i], NULL, 10);
3666 /* end of QMMM input */
3670 for (i = 0; (i < egcNR); i++)
3672 fprintf(stderr, "%-16s has %d element(s):", gtypes[i], groups->grps[i].nr);
3673 for (j = 0; (j < groups->grps[i].nr); j++)
3675 fprintf(stderr, " %s", *(groups->grpname[groups->grps[i].nm_ind[j]]));
3677 fprintf(stderr, "\n");
3681 nr = groups->grps[egcENER].nr;
3682 snew(ir->opts.egp_flags, nr*nr);
3684 bExcl = do_egp_flag(ir, groups, "energygrp-excl", is->egpexcl, EGP_EXCL);
3685 if (bExcl && ir->cutoff_scheme == ecutsVERLET)
3687 warning_error(wi, "Energy group exclusions are not (yet) implemented for the Verlet scheme");
3689 if (bExcl && EEL_FULL(ir->coulombtype))
3691 warning(wi, "Can not exclude the lattice Coulomb energy between energy groups");
3694 bTable = do_egp_flag(ir, groups, "energygrp-table", is->egptable, EGP_TABLE);
3695 if (bTable && !(ir->vdwtype == evdwUSER) &&
3696 !(ir->coulombtype == eelUSER) && !(ir->coulombtype == eelPMEUSER) &&
3697 !(ir->coulombtype == eelPMEUSERSWITCH))
3699 gmx_fatal(FARGS, "Can only have energy group pair tables in combination with user tables for VdW and/or Coulomb");
3702 decode_cos(is->efield_x, &(ir->ex[XX]));
3703 decode_cos(is->efield_xt, &(ir->et[XX]));
3704 decode_cos(is->efield_y, &(ir->ex[YY]));
3705 decode_cos(is->efield_yt, &(ir->et[YY]));
3706 decode_cos(is->efield_z, &(ir->ex[ZZ]));
3707 decode_cos(is->efield_zt, &(ir->et[ZZ]));
3711 do_adress_index(ir->adress, groups, gnames, &(ir->opts), wi);
3714 for (i = 0; (i < grps->nr); i++)
3726 static void check_disre(gmx_mtop_t *mtop)
3728 gmx_ffparams_t *ffparams;
3729 t_functype *functype;
3731 int i, ndouble, ftype;
3732 int label, old_label;
3734 if (gmx_mtop_ftype_count(mtop, F_DISRES) > 0)
3736 ffparams = &mtop->ffparams;
3737 functype = ffparams->functype;
3738 ip = ffparams->iparams;
3741 for (i = 0; i < ffparams->ntypes; i++)
3743 ftype = functype[i];
3744 if (ftype == F_DISRES)
3746 label = ip[i].disres.label;
3747 if (label == old_label)
3749 fprintf(stderr, "Distance restraint index %d occurs twice\n", label);
3757 gmx_fatal(FARGS, "Found %d double distance restraint indices,\n"
3758 "probably the parameters for multiple pairs in one restraint "
3759 "are not identical\n", ndouble);
3764 static gmx_bool absolute_reference(t_inputrec *ir, gmx_mtop_t *sys,
3765 gmx_bool posres_only,
3769 gmx_mtop_ilistloop_t iloop;
3779 for (d = 0; d < DIM; d++)
3781 AbsRef[d] = (d < ndof_com(ir) ? 0 : 1);
3783 /* Check for freeze groups */
3784 for (g = 0; g < ir->opts.ngfrz; g++)
3786 for (d = 0; d < DIM; d++)
3788 if (ir->opts.nFreeze[g][d] != 0)
3796 /* Check for position restraints */
3797 iloop = gmx_mtop_ilistloop_init(sys);
3798 while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
3801 (AbsRef[XX] == 0 || AbsRef[YY] == 0 || AbsRef[ZZ] == 0))
3803 for (i = 0; i < ilist[F_POSRES].nr; i += 2)
3805 pr = &sys->ffparams.iparams[ilist[F_POSRES].iatoms[i]];
3806 for (d = 0; d < DIM; d++)
3808 if (pr->posres.fcA[d] != 0)
3814 for (i = 0; i < ilist[F_FBPOSRES].nr; i += 2)
3816 /* Check for flat-bottom posres */
3817 pr = &sys->ffparams.iparams[ilist[F_FBPOSRES].iatoms[i]];
3818 if (pr->fbposres.k != 0)
3820 switch (pr->fbposres.geom)
3822 case efbposresSPHERE:
3823 AbsRef[XX] = AbsRef[YY] = AbsRef[ZZ] = 1;
3825 case efbposresCYLINDER:
3826 AbsRef[XX] = AbsRef[YY] = 1;
3828 case efbposresX: /* d=XX */
3829 case efbposresY: /* d=YY */
3830 case efbposresZ: /* d=ZZ */
3831 d = pr->fbposres.geom - efbposresX;
3835 gmx_fatal(FARGS, " Invalid geometry for flat-bottom position restraint.\n"
3836 "Expected nr between 1 and %d. Found %d\n", efbposresNR-1,
3844 return (AbsRef[XX] != 0 && AbsRef[YY] != 0 && AbsRef[ZZ] != 0);
3848 check_combination_rule_differences(const gmx_mtop_t *mtop, int state,
3849 gmx_bool *bC6ParametersWorkWithGeometricRules,
3850 gmx_bool *bC6ParametersWorkWithLBRules,
3851 gmx_bool *bLBRulesPossible)
3853 int ntypes, tpi, tpj, thisLBdiff, thisgeomdiff;
3856 double geometricdiff, LBdiff;
3857 double c6i, c6j, c12i, c12j;
3858 double c6, c6_geometric, c6_LB;
3859 double sigmai, sigmaj, epsi, epsj;
3860 gmx_bool bCanDoLBRules, bCanDoGeometricRules;
3863 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
3864 * force-field floating point parameters.
3867 ptr = getenv("GMX_LJCOMB_TOL");
3872 sscanf(ptr, "%lf", &dbl);
3876 *bC6ParametersWorkWithLBRules = TRUE;
3877 *bC6ParametersWorkWithGeometricRules = TRUE;
3878 bCanDoLBRules = TRUE;
3879 bCanDoGeometricRules = TRUE;
3880 ntypes = mtop->ffparams.atnr;
3881 snew(typecount, ntypes);
3882 gmx_mtop_count_atomtypes(mtop, state, typecount);
3883 geometricdiff = LBdiff = 0.0;
3884 *bLBRulesPossible = TRUE;
3885 for (tpi = 0; tpi < ntypes; ++tpi)
3887 c6i = mtop->ffparams.iparams[(ntypes + 1) * tpi].lj.c6;
3888 c12i = mtop->ffparams.iparams[(ntypes + 1) * tpi].lj.c12;
3889 for (tpj = tpi; tpj < ntypes; ++tpj)
3891 c6j = mtop->ffparams.iparams[(ntypes + 1) * tpj].lj.c6;
3892 c12j = mtop->ffparams.iparams[(ntypes + 1) * tpj].lj.c12;
3893 c6 = mtop->ffparams.iparams[ntypes * tpi + tpj].lj.c6;
3894 c6_geometric = sqrt(c6i * c6j);
3895 if (!gmx_numzero(c6_geometric))
3897 if (!gmx_numzero(c12i) && !gmx_numzero(c12j))
3899 sigmai = pow(c12i / c6i, 1.0/6.0);
3900 sigmaj = pow(c12j / c6j, 1.0/6.0);
3901 epsi = c6i * c6i /(4.0 * c12i);
3902 epsj = c6j * c6j /(4.0 * c12j);
3903 c6_LB = 4.0 * pow(epsi * epsj, 1.0/2.0) * pow(0.5 * (sigmai + sigmaj), 6);
3907 *bLBRulesPossible = FALSE;
3908 c6_LB = c6_geometric;
3910 bCanDoLBRules = gmx_within_tol(c6_LB, c6, tol);
3913 if (FALSE == bCanDoLBRules)
3915 *bC6ParametersWorkWithLBRules = FALSE;
3918 bCanDoGeometricRules = gmx_within_tol(c6_geometric, c6, tol);
3920 if (FALSE == bCanDoGeometricRules)
3922 *bC6ParametersWorkWithGeometricRules = FALSE;
3930 check_combination_rules(const t_inputrec *ir, const gmx_mtop_t *mtop,
3934 gmx_bool bLBRulesPossible, bC6ParametersWorkWithGeometricRules, bC6ParametersWorkWithLBRules;
3936 check_combination_rule_differences(mtop, 0,
3937 &bC6ParametersWorkWithGeometricRules,
3938 &bC6ParametersWorkWithLBRules,
3940 if (ir->ljpme_combination_rule == eljpmeLB)
3942 if (FALSE == bC6ParametersWorkWithLBRules || FALSE == bLBRulesPossible)
3944 warning(wi, "You are using arithmetic-geometric combination rules "
3945 "in LJ-PME, but your non-bonded C6 parameters do not "
3946 "follow these rules.");
3951 if (FALSE == bC6ParametersWorkWithGeometricRules)
3953 if (ir->eDispCorr != edispcNO)
3955 warning_note(wi, "You are using geometric combination rules in "
3956 "LJ-PME, but your non-bonded C6 parameters do "
3957 "not follow these rules. "
3958 "This will introduce very small errors in the forces and energies in "
3959 "your simulations. Dispersion correction will correct total energy "
3960 "and/or pressure for isotropic systems, but not forces or surface tensions.");
3964 warning_note(wi, "You are using geometric combination rules in "
3965 "LJ-PME, but your non-bonded C6 parameters do "
3966 "not follow these rules. "
3967 "This will introduce very small errors in the forces and energies in "
3968 "your simulations. If your system is homogeneous, consider using dispersion correction "
3969 "for the total energy and pressure.");
3975 void triple_check(const char *mdparin, t_inputrec *ir, gmx_mtop_t *sys,
3978 char err_buf[STRLEN];
3979 int i, m, c, nmol, npct;
3980 gmx_bool bCharge, bAcc;
3981 real gdt_max, *mgrp, mt;
3983 gmx_mtop_atomloop_block_t aloopb;
3984 gmx_mtop_atomloop_all_t aloop;
3987 char warn_buf[STRLEN];
3989 set_warning_line(wi, mdparin, -1);
3991 if (ir->cutoff_scheme == ecutsVERLET &&
3992 ir->verletbuf_tol > 0 &&
3994 ((EI_MD(ir->eI) || EI_SD(ir->eI)) &&
3995 (ir->etc == etcVRESCALE || ir->etc == etcBERENDSEN)))
3997 /* Check if a too small Verlet buffer might potentially
3998 * cause more drift than the thermostat can couple off.
4000 /* Temperature error fraction for warning and suggestion */
4001 const real T_error_warn = 0.002;
4002 const real T_error_suggest = 0.001;
4003 /* For safety: 2 DOF per atom (typical with constraints) */
4004 const real nrdf_at = 2;
4005 real T, tau, max_T_error;
4010 for (i = 0; i < ir->opts.ngtc; i++)
4012 T = max(T, ir->opts.ref_t[i]);
4013 tau = max(tau, ir->opts.tau_t[i]);
4017 /* This is a worst case estimate of the temperature error,
4018 * assuming perfect buffer estimation and no cancelation
4019 * of errors. The factor 0.5 is because energy distributes
4020 * equally over Ekin and Epot.
4022 max_T_error = 0.5*tau*ir->verletbuf_tol/(nrdf_at*BOLTZ*T);
4023 if (max_T_error > T_error_warn)
4025 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.",
4026 ir->verletbuf_tol, T, tau,
4028 100*T_error_suggest,
4029 ir->verletbuf_tol*T_error_suggest/max_T_error);
4030 warning(wi, warn_buf);
4035 if (ETC_ANDERSEN(ir->etc))
4039 for (i = 0; i < ir->opts.ngtc; i++)
4041 sprintf(err_buf, "all tau_t must currently be equal using Andersen temperature control, violated for group %d", i);
4042 CHECK(ir->opts.tau_t[0] != ir->opts.tau_t[i]);
4043 sprintf(err_buf, "all tau_t must be postive using Andersen temperature control, tau_t[%d]=%10.6f",
4044 i, ir->opts.tau_t[i]);
4045 CHECK(ir->opts.tau_t[i] < 0);
4048 for (i = 0; i < ir->opts.ngtc; i++)
4050 int nsteps = (int)(ir->opts.tau_t[i]/ir->delta_t);
4051 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);
4052 CHECK((nsteps % ir->nstcomm) && (ir->etc == etcANDERSENMASSIVE));
4056 if (EI_DYNAMICS(ir->eI) && !EI_SD(ir->eI) && ir->eI != eiBD &&
4057 ir->comm_mode == ecmNO &&
4058 !(absolute_reference(ir, sys, FALSE, AbsRef) || ir->nsteps <= 10) &&
4059 !ETC_ANDERSEN(ir->etc))
4061 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");
4064 /* Check for pressure coupling with absolute position restraints */
4065 if (ir->epc != epcNO && ir->refcoord_scaling == erscNO)
4067 absolute_reference(ir, sys, TRUE, AbsRef);
4069 for (m = 0; m < DIM; m++)
4071 if (AbsRef[m] && norm2(ir->compress[m]) > 0)
4073 warning(wi, "You are using pressure coupling with absolute position restraints, this will give artifacts. Use the refcoord_scaling option.");
4081 aloopb = gmx_mtop_atomloop_block_init(sys);
4082 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
4084 if (atom->q != 0 || atom->qB != 0)
4092 if (EEL_FULL(ir->coulombtype))
4095 "You are using full electrostatics treatment %s for a system without charges.\n"
4096 "This costs a lot of performance for just processing zeros, consider using %s instead.\n",
4097 EELTYPE(ir->coulombtype), EELTYPE(eelCUT));
4098 warning(wi, err_buf);
4103 if (ir->coulombtype == eelCUT && ir->rcoulomb > 0 && !ir->implicit_solvent)
4106 "You are using a plain Coulomb cut-off, which might produce artifacts.\n"
4107 "You might want to consider using %s electrostatics.\n",
4109 warning_note(wi, err_buf);
4113 /* Check if combination rules used in LJ-PME are the same as in the force field */
4114 if (EVDW_PME(ir->vdwtype))
4116 check_combination_rules(ir, sys, wi);
4119 /* Generalized reaction field */
4120 if (ir->opts.ngtc == 0)
4122 sprintf(err_buf, "No temperature coupling while using coulombtype %s",
4124 CHECK(ir->coulombtype == eelGRF);
4128 sprintf(err_buf, "When using coulombtype = %s"
4129 " ref-t for temperature coupling should be > 0",
4131 CHECK((ir->coulombtype == eelGRF) && (ir->opts.ref_t[0] <= 0));
4134 if (ir->eI == eiSD2)
4136 sprintf(warn_buf, "The stochastic dynamics integrator %s is deprecated, since\n"
4137 "it is slower than integrator %s and is slightly less accurate\n"
4138 "with constraints. Use the %s integrator.",
4139 ei_names[ir->eI], ei_names[eiSD1], ei_names[eiSD1]);
4140 warning_note(wi, warn_buf);
4144 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4146 for (m = 0; (m < DIM); m++)
4148 if (fabs(ir->opts.acc[i][m]) > 1e-6)
4157 snew(mgrp, sys->groups.grps[egcACC].nr);
4158 aloop = gmx_mtop_atomloop_all_init(sys);
4159 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
4161 mgrp[ggrpnr(&sys->groups, egcACC, i)] += atom->m;
4164 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4166 for (m = 0; (m < DIM); m++)
4168 acc[m] += ir->opts.acc[i][m]*mgrp[i];
4172 for (m = 0; (m < DIM); m++)
4174 if (fabs(acc[m]) > 1e-6)
4176 const char *dim[DIM] = { "X", "Y", "Z" };
4178 "Net Acceleration in %s direction, will %s be corrected\n",
4179 dim[m], ir->nstcomm != 0 ? "" : "not");
4180 if (ir->nstcomm != 0 && m < ndof_com(ir))
4183 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4185 ir->opts.acc[i][m] -= acc[m];
4193 if (ir->efep != efepNO && ir->fepvals->sc_alpha != 0 &&
4194 !gmx_within_tol(sys->ffparams.reppow, 12.0, 10*GMX_DOUBLE_EPS))
4196 gmx_fatal(FARGS, "Soft-core interactions are only supported with VdW repulsion power 12");
4199 if (ir->ePull != epullNO)
4201 gmx_bool bPullAbsoluteRef;
4203 bPullAbsoluteRef = FALSE;
4204 for (i = 0; i < ir->pull->ncoord; i++)
4206 bPullAbsoluteRef = bPullAbsoluteRef ||
4207 ir->pull->coord[i].group[0] == 0 ||
4208 ir->pull->coord[i].group[1] == 0;
4210 if (bPullAbsoluteRef)
4212 absolute_reference(ir, sys, FALSE, AbsRef);
4213 for (m = 0; m < DIM; m++)
4215 if (ir->pull->dim[m] && !AbsRef[m])
4217 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.");
4223 if (ir->pull->eGeom == epullgDIRPBC)
4225 for (i = 0; i < 3; i++)
4227 for (m = 0; m <= i; m++)
4229 if ((ir->epc != epcNO && ir->compress[i][m] != 0) ||
4230 ir->deform[i][m] != 0)
4232 for (c = 0; c < ir->pull->ncoord; c++)
4234 if (ir->pull->coord[c].vec[m] != 0)
4236 gmx_fatal(FARGS, "Can not have dynamic box while using pull geometry '%s' (dim %c)", EPULLGEOM(ir->pull->eGeom), 'x'+m);
4248 void double_check(t_inputrec *ir, matrix box, gmx_bool bConstr, warninp_t wi)
4252 char warn_buf[STRLEN];
4255 ptr = check_box(ir->ePBC, box);
4258 warning_error(wi, ptr);
4261 if (bConstr && ir->eConstrAlg == econtSHAKE)
4263 if (ir->shake_tol <= 0.0)
4265 sprintf(warn_buf, "ERROR: shake-tol must be > 0 instead of %g\n",
4267 warning_error(wi, warn_buf);
4270 if (IR_TWINRANGE(*ir) && ir->nstlist > 1)
4272 sprintf(warn_buf, "With twin-range cut-off's and SHAKE the virial and the pressure are incorrect.");
4273 if (ir->epc == epcNO)
4275 warning(wi, warn_buf);
4279 warning_error(wi, warn_buf);
4284 if ( (ir->eConstrAlg == econtLINCS) && bConstr)
4286 /* If we have Lincs constraints: */
4287 if (ir->eI == eiMD && ir->etc == etcNO &&
4288 ir->eConstrAlg == econtLINCS && ir->nLincsIter == 1)
4290 sprintf(warn_buf, "For energy conservation with LINCS, lincs_iter should be 2 or larger.\n");
4291 warning_note(wi, warn_buf);
4294 if ((ir->eI == eiCG || ir->eI == eiLBFGS) && (ir->nProjOrder < 8))
4296 sprintf(warn_buf, "For accurate %s with LINCS constraints, lincs-order should be 8 or more.", ei_names[ir->eI]);
4297 warning_note(wi, warn_buf);
4299 if (ir->epc == epcMTTK)
4301 warning_error(wi, "MTTK not compatible with lincs -- use shake instead.");
4305 if (bConstr && ir->epc == epcMTTK)
4307 warning_note(wi, "MTTK with constraints is deprecated, and will be removed in GROMACS 5.1");
4310 if (ir->LincsWarnAngle > 90.0)
4312 sprintf(warn_buf, "lincs-warnangle can not be larger than 90 degrees, setting it to 90.\n");
4313 warning(wi, warn_buf);
4314 ir->LincsWarnAngle = 90.0;
4317 if (ir->ePBC != epbcNONE)
4319 if (ir->nstlist == 0)
4321 warning(wi, "With nstlist=0 atoms are only put into the box at step 0, therefore drifting atoms might cause the simulation to crash.");
4323 bTWIN = (ir->rlistlong > ir->rlist);
4324 if (ir->ns_type == ensGRID)
4326 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
4328 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",
4329 bTWIN ? (ir->rcoulomb == ir->rlistlong ? "rcoulomb" : "rvdw") : "rlist");
4330 warning_error(wi, warn_buf);
4335 min_size = min(box[XX][XX], min(box[YY][YY], box[ZZ][ZZ]));
4336 if (2*ir->rlistlong >= min_size)
4338 sprintf(warn_buf, "ERROR: One of the box lengths is smaller than twice the cut-off length. Increase the box size or decrease rlist.");
4339 warning_error(wi, warn_buf);
4342 fprintf(stderr, "Grid search might allow larger cut-off's than simple search with triclinic boxes.");
4349 void check_chargegroup_radii(const gmx_mtop_t *mtop, const t_inputrec *ir,
4353 real rvdw1, rvdw2, rcoul1, rcoul2;
4354 char warn_buf[STRLEN];
4356 calc_chargegroup_radii(mtop, x, &rvdw1, &rvdw2, &rcoul1, &rcoul2);
4360 printf("Largest charge group radii for Van der Waals: %5.3f, %5.3f nm\n",
4365 printf("Largest charge group radii for Coulomb: %5.3f, %5.3f nm\n",
4371 if (rvdw1 + rvdw2 > ir->rlist ||
4372 rcoul1 + rcoul2 > ir->rlist)
4375 "The sum of the two largest charge group radii (%f) "
4376 "is larger than rlist (%f)\n",
4377 max(rvdw1+rvdw2, rcoul1+rcoul2), ir->rlist);
4378 warning(wi, warn_buf);
4382 /* Here we do not use the zero at cut-off macro,
4383 * since user defined interactions might purposely
4384 * not be zero at the cut-off.
4386 if (ir_vdw_is_zero_at_cutoff(ir) &&
4387 rvdw1 + rvdw2 > ir->rlistlong - ir->rvdw)
4389 sprintf(warn_buf, "The sum of the two largest charge group "
4390 "radii (%f) is larger than %s (%f) - rvdw (%f).\n"
4391 "With exact cut-offs, better performance can be "
4392 "obtained with cutoff-scheme = %s, because it "
4393 "does not use charge groups at all.",
4395 ir->rlistlong > ir->rlist ? "rlistlong" : "rlist",
4396 ir->rlistlong, ir->rvdw,
4397 ecutscheme_names[ecutsVERLET]);
4400 warning(wi, warn_buf);
4404 warning_note(wi, warn_buf);
4407 if (ir_coulomb_is_zero_at_cutoff(ir) &&
4408 rcoul1 + rcoul2 > ir->rlistlong - ir->rcoulomb)
4410 sprintf(warn_buf, "The sum of the two largest charge group radii (%f) is larger than %s (%f) - rcoulomb (%f).\n"
4411 "With exact cut-offs, better performance can be obtained with cutoff-scheme = %s, because it does not use charge groups at all.",
4413 ir->rlistlong > ir->rlist ? "rlistlong" : "rlist",
4414 ir->rlistlong, ir->rcoulomb,
4415 ecutscheme_names[ecutsVERLET]);
4418 warning(wi, warn_buf);
4422 warning_note(wi, warn_buf);