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/utility/smalloc.h"
49 #include "gmx_fatal.h"
53 #include "gromacs/utility/cstringutil.h"
62 #include "mtop_util.h"
63 #include "chargegroup.h"
65 #include "calc_verletbuf.h"
70 /* Resource parameters
71 * Do not change any of these until you read the instruction
72 * in readinp.h. Some cpp's do not take spaces after the backslash
73 * (like the c-shell), which will give you a very weird compiler
77 typedef struct t_inputrec_strings
79 char tcgrps[STRLEN], tau_t[STRLEN], ref_t[STRLEN],
80 acc[STRLEN], accgrps[STRLEN], freeze[STRLEN], frdim[STRLEN],
81 energy[STRLEN], user1[STRLEN], user2[STRLEN], vcm[STRLEN], x_compressed_groups[STRLEN],
82 couple_moltype[STRLEN], orirefitgrp[STRLEN], egptable[STRLEN], egpexcl[STRLEN],
83 wall_atomtype[STRLEN], wall_density[STRLEN], deform[STRLEN], QMMM[STRLEN],
85 char fep_lambda[efptNR][STRLEN];
86 char lambda_weights[STRLEN];
89 char anneal[STRLEN], anneal_npoints[STRLEN],
90 anneal_time[STRLEN], anneal_temp[STRLEN];
91 char QMmethod[STRLEN], QMbasis[STRLEN], QMcharge[STRLEN], QMmult[STRLEN],
92 bSH[STRLEN], CASorbitals[STRLEN], CASelectrons[STRLEN], SAon[STRLEN],
93 SAoff[STRLEN], SAsteps[STRLEN], bTS[STRLEN], bOPT[STRLEN];
94 char efield_x[STRLEN], efield_xt[STRLEN], efield_y[STRLEN],
95 efield_yt[STRLEN], efield_z[STRLEN], efield_zt[STRLEN];
97 } gmx_inputrec_strings;
99 static gmx_inputrec_strings *is = NULL;
101 void init_inputrec_strings()
105 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.");
110 void done_inputrec_strings()
116 static char swapgrp[STRLEN], splitgrp0[STRLEN], splitgrp1[STRLEN], solgrp[STRLEN];
119 egrptpALL, /* All particles have to be a member of a group. */
120 egrptpALL_GENREST, /* A rest group with name is generated for particles *
121 * that are not part of any group. */
122 egrptpPART, /* As egrptpALL_GENREST, but no name is generated *
123 * for the rest group. */
124 egrptpONE /* Merge all selected groups into one group, *
125 * make a rest group for the remaining particles. */
128 static const char *constraints[eshNR+1] = {
129 "none", "h-bonds", "all-bonds", "h-angles", "all-angles", NULL
132 static const char *couple_lam[ecouplamNR+1] = {
133 "vdw-q", "vdw", "q", "none", NULL
136 void init_ir(t_inputrec *ir, t_gromppopts *opts)
138 snew(opts->include, STRLEN);
139 snew(opts->define, STRLEN);
140 snew(ir->fepvals, 1);
141 snew(ir->expandedvals, 1);
142 snew(ir->simtempvals, 1);
145 static void GetSimTemps(int ntemps, t_simtemp *simtemp, double *temperature_lambdas)
150 for (i = 0; i < ntemps; i++)
152 /* simple linear scaling -- allows more control */
153 if (simtemp->eSimTempScale == esimtempLINEAR)
155 simtemp->temperatures[i] = simtemp->simtemp_low + (simtemp->simtemp_high-simtemp->simtemp_low)*temperature_lambdas[i];
157 else if (simtemp->eSimTempScale == esimtempGEOMETRIC) /* should give roughly equal acceptance for constant heat capacity . . . */
159 simtemp->temperatures[i] = simtemp->simtemp_low * pow(simtemp->simtemp_high/simtemp->simtemp_low, (1.0*i)/(ntemps-1));
161 else if (simtemp->eSimTempScale == esimtempEXPONENTIAL)
163 simtemp->temperatures[i] = simtemp->simtemp_low + (simtemp->simtemp_high-simtemp->simtemp_low)*((exp(temperature_lambdas[i])-1)/(exp(1.0)-1));
168 sprintf(errorstr, "eSimTempScale=%d not defined", simtemp->eSimTempScale);
169 gmx_fatal(FARGS, errorstr);
176 static void _low_check(gmx_bool b, char *s, warninp_t wi)
180 warning_error(wi, s);
184 static void check_nst(const char *desc_nst, int nst,
185 const char *desc_p, int *p,
190 if (*p > 0 && *p % nst != 0)
192 /* Round up to the next multiple of nst */
193 *p = ((*p)/nst + 1)*nst;
194 sprintf(buf, "%s should be a multiple of %s, changing %s to %d\n",
195 desc_p, desc_nst, desc_p, *p);
200 static gmx_bool ir_NVE(const t_inputrec *ir)
202 return ((ir->eI == eiMD || EI_VV(ir->eI)) && ir->etc == etcNO);
205 static int lcd(int n1, int n2)
210 for (i = 2; (i <= n1 && i <= n2); i++)
212 if (n1 % i == 0 && n2 % i == 0)
221 static void process_interaction_modifier(const t_inputrec *ir, int *eintmod)
223 if (*eintmod == eintmodPOTSHIFT_VERLET)
225 if (ir->cutoff_scheme == ecutsVERLET)
227 *eintmod = eintmodPOTSHIFT;
231 *eintmod = eintmodNONE;
236 void check_ir(const char *mdparin, t_inputrec *ir, t_gromppopts *opts,
238 /* Check internal consistency.
239 * NOTE: index groups are not set here yet, don't check things
240 * like temperature coupling group options here, but in triple_check
243 /* Strange macro: first one fills the err_buf, and then one can check
244 * the condition, which will print the message and increase the error
247 #define CHECK(b) _low_check(b, err_buf, wi)
248 char err_buf[256], warn_buf[STRLEN];
254 t_lambda *fep = ir->fepvals;
255 t_expanded *expand = ir->expandedvals;
257 set_warning_line(wi, mdparin, -1);
259 /* BASIC CUT-OFF STUFF */
260 if (ir->rcoulomb < 0)
262 warning_error(wi, "rcoulomb should be >= 0");
266 warning_error(wi, "rvdw should be >= 0");
269 !(ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0))
271 warning_error(wi, "rlist should be >= 0");
274 process_interaction_modifier(ir, &ir->coulomb_modifier);
275 process_interaction_modifier(ir, &ir->vdw_modifier);
277 if (ir->cutoff_scheme == ecutsGROUP)
280 "The group cutoff scheme is deprecated in Gromacs 5.0 and will be removed in a future "
281 "release when all interaction forms are supported for the verlet scheme. The verlet "
282 "scheme already scales better, and it is compatible with GPUs and other accelerators.");
284 /* BASIC CUT-OFF STUFF */
285 if (ir->rlist == 0 ||
286 !((ir_coulomb_might_be_zero_at_cutoff(ir) && ir->rcoulomb > ir->rlist) ||
287 (ir_vdw_might_be_zero_at_cutoff(ir) && ir->rvdw > ir->rlist)))
289 /* No switched potential and/or no twin-range:
290 * we can set the long-range cut-off to the maximum of the other cut-offs.
292 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
294 else if (ir->rlistlong < 0)
296 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
297 sprintf(warn_buf, "rlistlong was not set, setting it to %g (no buffer)",
299 warning(wi, warn_buf);
301 if (ir->rlistlong == 0 && ir->ePBC != epbcNONE)
303 warning_error(wi, "Can not have an infinite cut-off with PBC");
305 if (ir->rlistlong > 0 && (ir->rlist == 0 || ir->rlistlong < ir->rlist))
307 warning_error(wi, "rlistlong can not be shorter than rlist");
309 if (IR_TWINRANGE(*ir) && ir->nstlist <= 0)
311 warning_error(wi, "Can not have nstlist<=0 with twin-range interactions");
315 if (ir->rlistlong == ir->rlist)
319 else if (ir->rlistlong > ir->rlist && ir->nstcalclr == 0)
321 warning_error(wi, "With different cutoffs for electrostatics and VdW, nstcalclr must be -1 or a positive number");
324 if (ir->cutoff_scheme == ecutsVERLET)
328 /* Normal Verlet type neighbor-list, currently only limited feature support */
329 if (inputrec2nboundeddim(ir) < 3)
331 warning_error(wi, "With Verlet lists only full pbc or pbc=xy with walls is supported");
333 if (ir->rcoulomb != ir->rvdw)
335 warning_error(wi, "With Verlet lists rcoulomb!=rvdw is not supported");
337 if (ir->vdwtype == evdwSHIFT || ir->vdwtype == evdwSWITCH)
339 if (ir->vdw_modifier == eintmodNONE ||
340 ir->vdw_modifier == eintmodPOTSHIFT)
342 ir->vdw_modifier = (ir->vdwtype == evdwSHIFT ? eintmodFORCESWITCH : eintmodPOTSWITCH);
344 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]);
345 warning_note(wi, warn_buf);
347 ir->vdwtype = evdwCUT;
351 sprintf(warn_buf, "Unsupported combination of vdwtype=%s and vdw_modifier=%s", evdw_names[ir->vdwtype], eintmod_names[ir->vdw_modifier]);
352 warning_error(wi, warn_buf);
356 if (!(ir->vdwtype == evdwCUT || ir->vdwtype == evdwPME))
358 warning_error(wi, "With Verlet lists only cut-off and PME LJ interactions are supported");
360 if (!(ir->coulombtype == eelCUT ||
361 (EEL_RF(ir->coulombtype) && ir->coulombtype != eelRF_NEC) ||
362 EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD))
364 warning_error(wi, "With Verlet lists only cut-off, reaction-field, PME and Ewald electrostatics are supported");
366 if (!(ir->coulomb_modifier == eintmodNONE ||
367 ir->coulomb_modifier == eintmodPOTSHIFT))
369 sprintf(warn_buf, "coulomb_modifier=%s is not supported with the Verlet cut-off scheme", eintmod_names[ir->coulomb_modifier]);
370 warning_error(wi, warn_buf);
373 if (ir->implicit_solvent != eisNO)
375 warning_error(wi, "Implicit solvent is not (yet) supported with the with Verlet lists.");
378 if (ir->nstlist <= 0)
380 warning_error(wi, "With Verlet lists nstlist should be larger than 0");
383 if (ir->nstlist < 10)
385 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.");
388 rc_max = max(ir->rvdw, ir->rcoulomb);
390 if (ir->verletbuf_tol <= 0)
392 if (ir->verletbuf_tol == 0)
394 warning_error(wi, "Can not have Verlet buffer tolerance of exactly 0");
397 if (ir->rlist < rc_max)
399 warning_error(wi, "With verlet lists rlist can not be smaller than rvdw or rcoulomb");
402 if (ir->rlist == rc_max && ir->nstlist > 1)
404 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.");
409 if (ir->rlist > rc_max)
411 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.");
414 if (ir->nstlist == 1)
416 /* No buffer required */
421 if (EI_DYNAMICS(ir->eI))
423 if (inputrec2nboundeddim(ir) < 3)
425 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.");
427 /* Set rlist temporarily so we can continue processing */
432 /* Set the buffer to 5% of the cut-off */
433 ir->rlist = (1.0 + verlet_buffer_ratio_nodynamics)*rc_max;
438 /* No twin-range calculations with Verlet lists */
439 ir->rlistlong = ir->rlist;
442 if (ir->nstcalclr == -1)
444 /* if rlist=rlistlong, this will later be changed to nstcalclr=0 */
445 ir->nstcalclr = ir->nstlist;
447 else if (ir->nstcalclr > 0)
449 if (ir->nstlist > 0 && (ir->nstlist % ir->nstcalclr != 0))
451 warning_error(wi, "nstlist must be evenly divisible by nstcalclr. Use nstcalclr = -1 to automatically follow nstlist");
454 else if (ir->nstcalclr < -1)
456 warning_error(wi, "nstcalclr must be a positive number (divisor of nstcalclr), or -1 to follow nstlist.");
459 if (EEL_PME(ir->coulombtype) && ir->rcoulomb > ir->rvdw && ir->nstcalclr > 1)
461 warning_error(wi, "When used with PME, the long-range component of twin-range interactions must be updated every step (nstcalclr)");
464 /* GENERAL INTEGRATOR STUFF */
465 if (!(ir->eI == eiMD || EI_VV(ir->eI)))
469 if (ir->eI == eiVVAK)
471 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]);
472 warning_note(wi, warn_buf);
474 if (!EI_DYNAMICS(ir->eI))
478 if (EI_DYNAMICS(ir->eI))
480 if (ir->nstcalcenergy < 0)
482 ir->nstcalcenergy = ir_optimal_nstcalcenergy(ir);
483 if (ir->nstenergy != 0 && ir->nstenergy < ir->nstcalcenergy)
485 /* nstcalcenergy larger than nstener does not make sense.
486 * We ideally want nstcalcenergy=nstener.
490 ir->nstcalcenergy = lcd(ir->nstenergy, ir->nstlist);
494 ir->nstcalcenergy = ir->nstenergy;
498 else if ( (ir->nstenergy > 0 && ir->nstcalcenergy > ir->nstenergy) ||
499 (ir->efep != efepNO && ir->fepvals->nstdhdl > 0 &&
500 (ir->nstcalcenergy > ir->fepvals->nstdhdl) ) )
503 const char *nsten = "nstenergy";
504 const char *nstdh = "nstdhdl";
505 const char *min_name = nsten;
506 int min_nst = ir->nstenergy;
508 /* find the smallest of ( nstenergy, nstdhdl ) */
509 if (ir->efep != efepNO && ir->fepvals->nstdhdl > 0 &&
510 (ir->nstenergy == 0 || ir->fepvals->nstdhdl < ir->nstenergy))
512 min_nst = ir->fepvals->nstdhdl;
515 /* If the user sets nstenergy small, we should respect that */
517 "Setting nstcalcenergy (%d) equal to %s (%d)",
518 ir->nstcalcenergy, min_name, min_nst);
519 warning_note(wi, warn_buf);
520 ir->nstcalcenergy = min_nst;
523 if (ir->epc != epcNO)
525 if (ir->nstpcouple < 0)
527 ir->nstpcouple = ir_optimal_nstpcouple(ir);
530 if (IR_TWINRANGE(*ir))
532 check_nst("nstlist", ir->nstlist,
533 "nstcalcenergy", &ir->nstcalcenergy, wi);
534 if (ir->epc != epcNO)
536 check_nst("nstlist", ir->nstlist,
537 "nstpcouple", &ir->nstpcouple, wi);
541 if (ir->nstcalcenergy > 0)
543 if (ir->efep != efepNO)
545 /* nstdhdl should be a multiple of nstcalcenergy */
546 check_nst("nstcalcenergy", ir->nstcalcenergy,
547 "nstdhdl", &ir->fepvals->nstdhdl, wi);
548 /* nstexpanded should be a multiple of nstcalcenergy */
549 check_nst("nstcalcenergy", ir->nstcalcenergy,
550 "nstexpanded", &ir->expandedvals->nstexpanded, wi);
552 /* for storing exact averages nstenergy should be
553 * a multiple of nstcalcenergy
555 check_nst("nstcalcenergy", ir->nstcalcenergy,
556 "nstenergy", &ir->nstenergy, wi);
560 if (ir->nsteps == 0 && !ir->bContinuation)
562 warning_note(wi, "For a correct single-point energy evaluation with nsteps = 0, use continuation = yes to avoid constraining the input coordinates.");
566 if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
567 ir->bContinuation && ir->ld_seed != -1)
569 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)");
575 sprintf(err_buf, "TPI only works with pbc = %s", epbc_names[epbcXYZ]);
576 CHECK(ir->ePBC != epbcXYZ);
577 sprintf(err_buf, "TPI only works with ns = %s", ens_names[ensGRID]);
578 CHECK(ir->ns_type != ensGRID);
579 sprintf(err_buf, "with TPI nstlist should be larger than zero");
580 CHECK(ir->nstlist <= 0);
581 sprintf(err_buf, "TPI does not work with full electrostatics other than PME");
582 CHECK(EEL_FULL(ir->coulombtype) && !EEL_PME(ir->coulombtype));
583 sprintf(err_buf, "TPI does not work (yet) with the Verlet cut-off scheme");
584 CHECK(ir->cutoff_scheme == ecutsVERLET);
588 if ( (opts->nshake > 0) && (opts->bMorse) )
591 "Using morse bond-potentials while constraining bonds is useless");
592 warning(wi, warn_buf);
595 if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
596 ir->bContinuation && ir->ld_seed != -1)
598 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)");
600 /* verify simulated tempering options */
604 gmx_bool bAllTempZero = TRUE;
605 for (i = 0; i < fep->n_lambda; i++)
607 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]);
608 CHECK((fep->all_lambda[efptTEMPERATURE][i] < 0) || (fep->all_lambda[efptTEMPERATURE][i] > 1));
609 if (fep->all_lambda[efptTEMPERATURE][i] > 0)
611 bAllTempZero = FALSE;
614 sprintf(err_buf, "if simulated tempering is on, temperature-lambdas may not be all zero");
615 CHECK(bAllTempZero == TRUE);
617 sprintf(err_buf, "Simulated tempering is currently only compatible with md-vv");
618 CHECK(ir->eI != eiVV);
620 /* check compatability of the temperature coupling with simulated tempering */
622 if (ir->etc == etcNOSEHOOVER)
624 sprintf(warn_buf, "Nose-Hoover based temperature control such as [%s] my not be entirelyconsistent with simulated tempering", etcoupl_names[ir->etc]);
625 warning_note(wi, warn_buf);
628 /* check that the temperatures make sense */
630 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);
631 CHECK(ir->simtempvals->simtemp_high <= ir->simtempvals->simtemp_low);
633 sprintf(err_buf, "Higher simulated tempering temperature (%g) must be >= zero", ir->simtempvals->simtemp_high);
634 CHECK(ir->simtempvals->simtemp_high <= 0);
636 sprintf(err_buf, "Lower simulated tempering temperature (%g) must be >= zero", ir->simtempvals->simtemp_low);
637 CHECK(ir->simtempvals->simtemp_low <= 0);
640 /* verify free energy options */
642 if (ir->efep != efepNO)
645 sprintf(err_buf, "The soft-core power is %d and can only be 1 or 2",
647 CHECK(fep->sc_alpha != 0 && fep->sc_power != 1 && fep->sc_power != 2);
649 sprintf(err_buf, "The soft-core sc-r-power is %d and can only be 6 or 48",
650 (int)fep->sc_r_power);
651 CHECK(fep->sc_alpha != 0 && fep->sc_r_power != 6.0 && fep->sc_r_power != 48.0);
653 sprintf(err_buf, "Can't use postive delta-lambda (%g) if initial state/lambda does not start at zero", fep->delta_lambda);
654 CHECK(fep->delta_lambda > 0 && ((fep->init_fep_state > 0) || (fep->init_lambda > 0)));
656 sprintf(err_buf, "Can't use postive delta-lambda (%g) with expanded ensemble simulations", fep->delta_lambda);
657 CHECK(fep->delta_lambda > 0 && (ir->efep == efepEXPANDED));
659 sprintf(err_buf, "Can only use expanded ensemble with md-vv for now; should be supported for other integrators in 5.0");
660 CHECK(!(EI_VV(ir->eI)) && (ir->efep == efepEXPANDED));
662 sprintf(err_buf, "Free-energy not implemented for Ewald");
663 CHECK(ir->coulombtype == eelEWALD);
665 /* check validty of lambda inputs */
666 if (fep->n_lambda == 0)
668 /* Clear output in case of no states:*/
669 sprintf(err_buf, "init-lambda-state set to %d: no lambda states are defined.", fep->init_fep_state);
670 CHECK((fep->init_fep_state >= 0) && (fep->n_lambda == 0));
674 sprintf(err_buf, "initial thermodynamic state %d does not exist, only goes to %d", fep->init_fep_state, fep->n_lambda-1);
675 CHECK((fep->init_fep_state >= fep->n_lambda));
678 sprintf(err_buf, "Lambda state must be set, either with init-lambda-state or with init-lambda");
679 CHECK((fep->init_fep_state < 0) && (fep->init_lambda < 0));
681 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",
682 fep->init_lambda, fep->init_fep_state);
683 CHECK((fep->init_fep_state >= 0) && (fep->init_lambda >= 0));
687 if ((fep->init_lambda >= 0) && (fep->delta_lambda == 0))
691 for (i = 0; i < efptNR; i++)
693 if (fep->separate_dvdl[i])
698 if (n_lambda_terms > 1)
700 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.");
701 warning(wi, warn_buf);
704 if (n_lambda_terms < 2 && fep->n_lambda > 0)
707 "init-lambda is deprecated for setting lambda state (except for slow growth). Use init-lambda-state instead.");
711 for (j = 0; j < efptNR; j++)
713 for (i = 0; i < fep->n_lambda; i++)
715 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]);
716 CHECK((fep->all_lambda[j][i] < 0) || (fep->all_lambda[j][i] > 1));
720 if ((fep->sc_alpha > 0) && (!fep->bScCoul))
722 for (i = 0; i < fep->n_lambda; i++)
724 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],
725 fep->all_lambda[efptCOUL][i]);
726 CHECK((fep->sc_alpha > 0) &&
727 (((fep->all_lambda[efptCOUL][i] > 0.0) &&
728 (fep->all_lambda[efptCOUL][i] < 1.0)) &&
729 ((fep->all_lambda[efptVDW][i] > 0.0) &&
730 (fep->all_lambda[efptVDW][i] < 1.0))));
734 if ((fep->bScCoul) && (EEL_PME(ir->coulombtype)))
736 real sigma, lambda, r_sc;
739 /* Maximum estimate for A and B charges equal with lambda power 1 */
741 r_sc = pow(lambda*fep->sc_alpha*pow(sigma/ir->rcoulomb, fep->sc_r_power) + 1.0, 1.0/fep->sc_r_power);
742 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.",
744 sigma, lambda, r_sc - 1.0, ir->ewald_rtol);
745 warning_note(wi, warn_buf);
748 /* Free Energy Checks -- In an ideal world, slow growth and FEP would
749 be treated differently, but that's the next step */
751 for (i = 0; i < efptNR; i++)
753 for (j = 0; j < fep->n_lambda; j++)
755 sprintf(err_buf, "%s[%d] must be between 0 and 1", efpt_names[i], j);
756 CHECK((fep->all_lambda[i][j] < 0) || (fep->all_lambda[i][j] > 1));
761 if ((ir->bSimTemp) || (ir->efep == efepEXPANDED))
764 expand = ir->expandedvals;
766 /* checking equilibration of weights inputs for validity */
768 sprintf(err_buf, "weight-equil-number-all-lambda (%d) is ignored if lmc-weights-equil is not equal to %s",
769 expand->equil_n_at_lam, elmceq_names[elmceqNUMATLAM]);
770 CHECK((expand->equil_n_at_lam > 0) && (expand->elmceq != elmceqNUMATLAM));
772 sprintf(err_buf, "weight-equil-number-samples (%d) is ignored if lmc-weights-equil is not equal to %s",
773 expand->equil_samples, elmceq_names[elmceqSAMPLES]);
774 CHECK((expand->equil_samples > 0) && (expand->elmceq != elmceqSAMPLES));
776 sprintf(err_buf, "weight-equil-number-steps (%d) is ignored if lmc-weights-equil is not equal to %s",
777 expand->equil_steps, elmceq_names[elmceqSTEPS]);
778 CHECK((expand->equil_steps > 0) && (expand->elmceq != elmceqSTEPS));
780 sprintf(err_buf, "weight-equil-wl-delta (%d) is ignored if lmc-weights-equil is not equal to %s",
781 expand->equil_samples, elmceq_names[elmceqWLDELTA]);
782 CHECK((expand->equil_wl_delta > 0) && (expand->elmceq != elmceqWLDELTA));
784 sprintf(err_buf, "weight-equil-count-ratio (%f) is ignored if lmc-weights-equil is not equal to %s",
785 expand->equil_ratio, elmceq_names[elmceqRATIO]);
786 CHECK((expand->equil_ratio > 0) && (expand->elmceq != elmceqRATIO));
788 sprintf(err_buf, "weight-equil-number-all-lambda (%d) must be a positive integer if lmc-weights-equil=%s",
789 expand->equil_n_at_lam, elmceq_names[elmceqNUMATLAM]);
790 CHECK((expand->equil_n_at_lam <= 0) && (expand->elmceq == elmceqNUMATLAM));
792 sprintf(err_buf, "weight-equil-number-samples (%d) must be a positive integer if lmc-weights-equil=%s",
793 expand->equil_samples, elmceq_names[elmceqSAMPLES]);
794 CHECK((expand->equil_samples <= 0) && (expand->elmceq == elmceqSAMPLES));
796 sprintf(err_buf, "weight-equil-number-steps (%d) must be a positive integer if lmc-weights-equil=%s",
797 expand->equil_steps, elmceq_names[elmceqSTEPS]);
798 CHECK((expand->equil_steps <= 0) && (expand->elmceq == elmceqSTEPS));
800 sprintf(err_buf, "weight-equil-wl-delta (%f) must be > 0 if lmc-weights-equil=%s",
801 expand->equil_wl_delta, elmceq_names[elmceqWLDELTA]);
802 CHECK((expand->equil_wl_delta <= 0) && (expand->elmceq == elmceqWLDELTA));
804 sprintf(err_buf, "weight-equil-count-ratio (%f) must be > 0 if lmc-weights-equil=%s",
805 expand->equil_ratio, elmceq_names[elmceqRATIO]);
806 CHECK((expand->equil_ratio <= 0) && (expand->elmceq == elmceqRATIO));
808 sprintf(err_buf, "lmc-weights-equil=%s only possible when lmc-stats = %s or lmc-stats %s",
809 elmceq_names[elmceqWLDELTA], elamstats_names[elamstatsWL], elamstats_names[elamstatsWWL]);
810 CHECK((expand->elmceq == elmceqWLDELTA) && (!EWL(expand->elamstats)));
812 sprintf(err_buf, "lmc-repeats (%d) must be greater than 0", expand->lmc_repeats);
813 CHECK((expand->lmc_repeats <= 0));
814 sprintf(err_buf, "minimum-var-min (%d) must be greater than 0", expand->minvarmin);
815 CHECK((expand->minvarmin <= 0));
816 sprintf(err_buf, "weight-c-range (%d) must be greater or equal to 0", expand->c_range);
817 CHECK((expand->c_range < 0));
818 sprintf(err_buf, "init-lambda-state (%d) must be zero if lmc-forced-nstart (%d)> 0 and lmc-move != 'no'",
819 fep->init_fep_state, expand->lmc_forced_nstart);
820 CHECK((fep->init_fep_state != 0) && (expand->lmc_forced_nstart > 0) && (expand->elmcmove != elmcmoveNO));
821 sprintf(err_buf, "lmc-forced-nstart (%d) must not be negative", expand->lmc_forced_nstart);
822 CHECK((expand->lmc_forced_nstart < 0));
823 sprintf(err_buf, "init-lambda-state (%d) must be in the interval [0,number of lambdas)", fep->init_fep_state);
824 CHECK((fep->init_fep_state < 0) || (fep->init_fep_state >= fep->n_lambda));
826 sprintf(err_buf, "init-wl-delta (%f) must be greater than or equal to 0", expand->init_wl_delta);
827 CHECK((expand->init_wl_delta < 0));
828 sprintf(err_buf, "wl-ratio (%f) must be between 0 and 1", expand->wl_ratio);
829 CHECK((expand->wl_ratio <= 0) || (expand->wl_ratio >= 1));
830 sprintf(err_buf, "wl-scale (%f) must be between 0 and 1", expand->wl_scale);
831 CHECK((expand->wl_scale <= 0) || (expand->wl_scale >= 1));
833 /* if there is no temperature control, we need to specify an MC temperature */
834 sprintf(err_buf, "If there is no temperature control, and lmc-mcmove!= 'no',mc_temperature must be set to a positive number");
835 if (expand->nstTij > 0)
837 sprintf(err_buf, "nst-transition-matrix (%d) must be an integer multiple of nstlog (%d)",
838 expand->nstTij, ir->nstlog);
839 CHECK((mod(expand->nstTij, ir->nstlog) != 0));
844 sprintf(err_buf, "walls only work with pbc=%s", epbc_names[epbcXY]);
845 CHECK(ir->nwall && ir->ePBC != epbcXY);
848 if (ir->ePBC != epbcXYZ && ir->nwall != 2)
850 if (ir->ePBC == epbcNONE)
852 if (ir->epc != epcNO)
854 warning(wi, "Turning off pressure coupling for vacuum system");
860 sprintf(err_buf, "Can not have pressure coupling with pbc=%s",
861 epbc_names[ir->ePBC]);
862 CHECK(ir->epc != epcNO);
864 sprintf(err_buf, "Can not have Ewald with pbc=%s", epbc_names[ir->ePBC]);
865 CHECK(EEL_FULL(ir->coulombtype));
867 sprintf(err_buf, "Can not have dispersion correction with pbc=%s",
868 epbc_names[ir->ePBC]);
869 CHECK(ir->eDispCorr != edispcNO);
872 if (ir->rlist == 0.0)
874 sprintf(err_buf, "can only have neighborlist cut-off zero (=infinite)\n"
875 "with coulombtype = %s or coulombtype = %s\n"
876 "without periodic boundary conditions (pbc = %s) and\n"
877 "rcoulomb and rvdw set to zero",
878 eel_names[eelCUT], eel_names[eelUSER], epbc_names[epbcNONE]);
879 CHECK(((ir->coulombtype != eelCUT) && (ir->coulombtype != eelUSER)) ||
880 (ir->ePBC != epbcNONE) ||
881 (ir->rcoulomb != 0.0) || (ir->rvdw != 0.0));
885 warning_error(wi, "Can not have heuristic neighborlist updates without cut-off");
889 warning_note(wi, "Simulating without cut-offs can be (slightly) faster with nstlist=0, nstype=simple and only one MPI rank");
894 if (ir->nstcomm == 0)
896 ir->comm_mode = ecmNO;
898 if (ir->comm_mode != ecmNO)
902 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");
903 ir->nstcomm = abs(ir->nstcomm);
906 if (ir->nstcalcenergy > 0 && ir->nstcomm < ir->nstcalcenergy)
908 warning_note(wi, "nstcomm < nstcalcenergy defeats the purpose of nstcalcenergy, setting nstcomm to nstcalcenergy");
909 ir->nstcomm = ir->nstcalcenergy;
912 if (ir->comm_mode == ecmANGULAR)
914 sprintf(err_buf, "Can not remove the rotation around the center of mass with periodic molecules");
915 CHECK(ir->bPeriodicMols);
916 if (ir->ePBC != epbcNONE)
918 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).");
923 if (EI_STATE_VELOCITY(ir->eI) && ir->ePBC == epbcNONE && ir->comm_mode != ecmANGULAR)
925 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.");
928 sprintf(err_buf, "Twin-range neighbour searching (NS) with simple NS"
929 " algorithm not implemented");
930 CHECK(((ir->rcoulomb > ir->rlist) || (ir->rvdw > ir->rlist))
931 && (ir->ns_type == ensSIMPLE));
933 /* TEMPERATURE COUPLING */
934 if (ir->etc == etcYES)
936 ir->etc = etcBERENDSEN;
937 warning_note(wi, "Old option for temperature coupling given: "
938 "changing \"yes\" to \"Berendsen\"\n");
941 if ((ir->etc == etcNOSEHOOVER) || (ir->epc == epcMTTK))
943 if (ir->opts.nhchainlength < 1)
945 sprintf(warn_buf, "number of Nose-Hoover chains (currently %d) cannot be less than 1,reset to 1\n", ir->opts.nhchainlength);
946 ir->opts.nhchainlength = 1;
947 warning(wi, warn_buf);
950 if (ir->etc == etcNOSEHOOVER && !EI_VV(ir->eI) && ir->opts.nhchainlength > 1)
952 warning_note(wi, "leapfrog does not yet support Nose-Hoover chains, nhchainlength reset to 1");
953 ir->opts.nhchainlength = 1;
958 ir->opts.nhchainlength = 0;
961 if (ir->eI == eiVVAK)
963 sprintf(err_buf, "%s implemented primarily for validation, and requires nsttcouple = 1 and nstpcouple = 1.",
965 CHECK((ir->nsttcouple != 1) || (ir->nstpcouple != 1));
968 if (ETC_ANDERSEN(ir->etc))
970 sprintf(err_buf, "%s temperature control not supported for integrator %s.", etcoupl_names[ir->etc], ei_names[ir->eI]);
971 CHECK(!(EI_VV(ir->eI)));
973 if (ir->nstcomm > 0 && (ir->etc == etcANDERSEN))
975 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]);
976 warning_note(wi, warn_buf);
979 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]);
980 CHECK(ir->nstcomm > 1 && (ir->etc == etcANDERSEN));
983 if (ir->etc == etcBERENDSEN)
985 sprintf(warn_buf, "The %s thermostat does not generate the correct kinetic energy distribution. You might want to consider using the %s thermostat.",
986 ETCOUPLTYPE(ir->etc), ETCOUPLTYPE(etcVRESCALE));
987 warning_note(wi, warn_buf);
990 if ((ir->etc == etcNOSEHOOVER || ETC_ANDERSEN(ir->etc))
991 && ir->epc == epcBERENDSEN)
993 sprintf(warn_buf, "Using Berendsen pressure coupling invalidates the "
994 "true ensemble for the thermostat");
995 warning(wi, warn_buf);
998 /* PRESSURE COUPLING */
999 if (ir->epc == epcISOTROPIC)
1001 ir->epc = epcBERENDSEN;
1002 warning_note(wi, "Old option for pressure coupling given: "
1003 "changing \"Isotropic\" to \"Berendsen\"\n");
1006 if (ir->epc != epcNO)
1008 dt_pcoupl = ir->nstpcouple*ir->delta_t;
1010 sprintf(err_buf, "tau-p must be > 0 instead of %g\n", ir->tau_p);
1011 CHECK(ir->tau_p <= 0);
1013 if (ir->tau_p/dt_pcoupl < pcouple_min_integration_steps(ir->epc) - 10*GMX_REAL_EPS)
1015 sprintf(warn_buf, "For proper integration of the %s barostat, tau-p (%g) should be at least %d times larger than nstpcouple*dt (%g)",
1016 EPCOUPLTYPE(ir->epc), ir->tau_p, pcouple_min_integration_steps(ir->epc), dt_pcoupl);
1017 warning(wi, warn_buf);
1020 sprintf(err_buf, "compressibility must be > 0 when using pressure"
1021 " coupling %s\n", EPCOUPLTYPE(ir->epc));
1022 CHECK(ir->compress[XX][XX] < 0 || ir->compress[YY][YY] < 0 ||
1023 ir->compress[ZZ][ZZ] < 0 ||
1024 (trace(ir->compress) == 0 && ir->compress[YY][XX] <= 0 &&
1025 ir->compress[ZZ][XX] <= 0 && ir->compress[ZZ][YY] <= 0));
1027 if (epcPARRINELLORAHMAN == ir->epc && opts->bGenVel)
1030 "You are generating velocities so I am assuming you "
1031 "are equilibrating a system. You are using "
1032 "%s pressure coupling, but this can be "
1033 "unstable for equilibration. If your system crashes, try "
1034 "equilibrating first with Berendsen pressure coupling. If "
1035 "you are not equilibrating the system, you can probably "
1036 "ignore this warning.",
1037 epcoupl_names[ir->epc]);
1038 warning(wi, warn_buf);
1044 if (ir->epc > epcNO)
1046 if ((ir->epc != epcBERENDSEN) && (ir->epc != epcMTTK))
1048 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.");
1054 if (ir->epc == epcMTTK)
1056 warning_error(wi, "MTTK pressure coupling requires a Velocity-verlet integrator");
1060 /* ELECTROSTATICS */
1061 /* More checks are in triple check (grompp.c) */
1063 if (ir->coulombtype == eelSWITCH)
1065 sprintf(warn_buf, "coulombtype = %s is only for testing purposes and can lead to serious "
1066 "artifacts, advice: use coulombtype = %s",
1067 eel_names[ir->coulombtype],
1068 eel_names[eelRF_ZERO]);
1069 warning(wi, warn_buf);
1072 if (ir->epsilon_r != 1 && ir->implicit_solvent == eisGBSA)
1074 sprintf(warn_buf, "epsilon-r = %g with GB implicit solvent, will use this value for inner dielectric", ir->epsilon_r);
1075 warning_note(wi, warn_buf);
1078 if (EEL_RF(ir->coulombtype) && ir->epsilon_rf == 1 && ir->epsilon_r != 1)
1080 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);
1081 warning(wi, warn_buf);
1082 ir->epsilon_rf = ir->epsilon_r;
1083 ir->epsilon_r = 1.0;
1086 if (getenv("GMX_DO_GALACTIC_DYNAMICS") == NULL)
1088 sprintf(err_buf, "epsilon-r must be >= 0 instead of %g\n", ir->epsilon_r);
1089 CHECK(ir->epsilon_r < 0);
1092 if (EEL_RF(ir->coulombtype))
1094 /* reaction field (at the cut-off) */
1096 if (ir->coulombtype == eelRF_ZERO)
1098 sprintf(warn_buf, "With coulombtype = %s, epsilon-rf must be 0, assuming you meant epsilon_rf=0",
1099 eel_names[ir->coulombtype]);
1100 CHECK(ir->epsilon_rf != 0);
1101 ir->epsilon_rf = 0.0;
1104 sprintf(err_buf, "epsilon-rf must be >= epsilon-r");
1105 CHECK((ir->epsilon_rf < ir->epsilon_r && ir->epsilon_rf != 0) ||
1106 (ir->epsilon_r == 0));
1107 if (ir->epsilon_rf == ir->epsilon_r)
1109 sprintf(warn_buf, "Using epsilon-rf = epsilon-r with %s does not make sense",
1110 eel_names[ir->coulombtype]);
1111 warning(wi, warn_buf);
1114 /* Allow rlist>rcoulomb for tabulated long range stuff. This just
1115 * means the interaction is zero outside rcoulomb, but it helps to
1116 * provide accurate energy conservation.
1118 if (ir_coulomb_might_be_zero_at_cutoff(ir))
1120 if (ir_coulomb_switched(ir))
1123 "With coulombtype = %s rcoulomb_switch must be < rcoulomb. Or, better: Use the potential modifier options!",
1124 eel_names[ir->coulombtype]);
1125 CHECK(ir->rcoulomb_switch >= ir->rcoulomb);
1128 else if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
1130 if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE)
1132 sprintf(err_buf, "With coulombtype = %s, rcoulomb should be >= rlist unless you use a potential modifier",
1133 eel_names[ir->coulombtype]);
1134 CHECK(ir->rlist > ir->rcoulomb);
1138 if (ir->coulombtype == eelSWITCH || ir->coulombtype == eelSHIFT)
1141 "Explicit switch/shift coulomb interactions cannot be used in combination with a secondary coulomb-modifier.");
1142 CHECK( ir->coulomb_modifier != eintmodNONE);
1144 if (ir->vdwtype == evdwSWITCH || ir->vdwtype == evdwSHIFT)
1147 "Explicit switch/shift vdw interactions cannot be used in combination with a secondary vdw-modifier.");
1148 CHECK( ir->vdw_modifier != eintmodNONE);
1151 if (ir->coulombtype == eelSWITCH || ir->coulombtype == eelSHIFT ||
1152 ir->vdwtype == evdwSWITCH || ir->vdwtype == evdwSHIFT)
1155 "The switch/shift interaction settings are just for compatibility; you will get better "
1156 "performance from applying potential modifiers to your interactions!\n");
1157 warning_note(wi, warn_buf);
1160 if (ir->coulombtype == eelPMESWITCH || ir->coulomb_modifier == eintmodPOTSWITCH)
1162 if (ir->rcoulomb_switch/ir->rcoulomb < 0.9499)
1164 real percentage = 100*(ir->rcoulomb-ir->rcoulomb_switch)/ir->rcoulomb;
1165 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.",
1166 percentage, ir->rcoulomb_switch, ir->rcoulomb, ir->ewald_rtol);
1167 warning(wi, warn_buf);
1171 if (ir->vdwtype == evdwSWITCH || ir->vdw_modifier == eintmodPOTSWITCH)
1173 if (ir->rvdw_switch == 0)
1175 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.");
1176 warning(wi, warn_buf);
1180 if (EEL_FULL(ir->coulombtype))
1182 if (ir->coulombtype == eelPMESWITCH || ir->coulombtype == eelPMEUSER ||
1183 ir->coulombtype == eelPMEUSERSWITCH)
1185 sprintf(err_buf, "With coulombtype = %s, rcoulomb must be <= rlist",
1186 eel_names[ir->coulombtype]);
1187 CHECK(ir->rcoulomb > ir->rlist);
1189 else if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE)
1191 if (ir->coulombtype == eelPME || ir->coulombtype == eelP3M_AD)
1194 "With coulombtype = %s (without modifier), rcoulomb must be equal to rlist,\n"
1195 "or rlistlong if nstcalclr=1. For optimal energy conservation,consider using\n"
1196 "a potential modifier.", eel_names[ir->coulombtype]);
1197 if (ir->nstcalclr == 1)
1199 CHECK(ir->rcoulomb != ir->rlist && ir->rcoulomb != ir->rlistlong);
1203 CHECK(ir->rcoulomb != ir->rlist);
1209 if (EEL_PME(ir->coulombtype) || EVDW_PME(ir->vdwtype))
1211 if (ir->pme_order < 3)
1213 warning_error(wi, "pme-order can not be smaller than 3");
1217 if (ir->nwall == 2 && EEL_FULL(ir->coulombtype))
1219 if (ir->ewald_geometry == eewg3D)
1221 sprintf(warn_buf, "With pbc=%s you should use ewald-geometry=%s",
1222 epbc_names[ir->ePBC], eewg_names[eewg3DC]);
1223 warning(wi, warn_buf);
1225 /* This check avoids extra pbc coding for exclusion corrections */
1226 sprintf(err_buf, "wall-ewald-zfac should be >= 2");
1227 CHECK(ir->wall_ewald_zfac < 2);
1230 if (ir_vdw_switched(ir))
1232 sprintf(err_buf, "With switched vdw forces or potentials, rvdw-switch must be < rvdw");
1233 CHECK(ir->rvdw_switch >= ir->rvdw);
1235 if (ir->rvdw_switch < 0.5*ir->rvdw)
1237 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.",
1238 ir->rvdw_switch, ir->rvdw);
1239 warning_note(wi, warn_buf);
1242 else if (ir->vdwtype == evdwCUT || ir->vdwtype == evdwPME)
1244 if (ir->cutoff_scheme == ecutsGROUP && ir->vdw_modifier == eintmodNONE)
1246 sprintf(err_buf, "With vdwtype = %s, rvdw must be >= rlist unless you use a potential modifier", evdw_names[ir->vdwtype]);
1247 CHECK(ir->rlist > ir->rvdw);
1251 if (ir->vdwtype == evdwPME)
1253 if (!(ir->vdw_modifier == eintmodNONE || ir->vdw_modifier == eintmodPOTSHIFT))
1255 sprintf(err_buf, "With vdwtype = %s, the only supported modifiers are %s a\
1257 evdw_names[ir->vdwtype],
1258 eintmod_names[eintmodPOTSHIFT],
1259 eintmod_names[eintmodNONE]);
1263 if (ir->cutoff_scheme == ecutsGROUP)
1265 if (((ir->coulomb_modifier != eintmodNONE && ir->rcoulomb == ir->rlist) ||
1266 (ir->vdw_modifier != eintmodNONE && ir->rvdw == ir->rlist)) &&
1269 warning_note(wi, "With exact cut-offs, rlist should be "
1270 "larger than rcoulomb and rvdw, so that there "
1271 "is a buffer region for particle motion "
1272 "between neighborsearch steps");
1275 if (ir_coulomb_is_zero_at_cutoff(ir) && ir->rlistlong <= ir->rcoulomb)
1277 sprintf(warn_buf, "For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rcoulomb.",
1278 IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
1279 warning_note(wi, warn_buf);
1281 if (ir_vdw_switched(ir) && (ir->rlistlong <= ir->rvdw))
1283 sprintf(warn_buf, "For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rvdw.",
1284 IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
1285 warning_note(wi, warn_buf);
1289 if (ir->vdwtype == evdwUSER && ir->eDispCorr != edispcNO)
1291 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.");
1294 if (ir->nstlist == -1)
1296 sprintf(err_buf, "With nstlist=-1 rvdw and rcoulomb should be smaller than rlist to account for diffusion and possibly charge-group radii");
1297 CHECK(ir->rvdw >= ir->rlist || ir->rcoulomb >= ir->rlist);
1299 sprintf(err_buf, "nstlist can not be smaller than -1");
1300 CHECK(ir->nstlist < -1);
1302 if (ir->eI == eiLBFGS && (ir->coulombtype == eelCUT || ir->vdwtype == evdwCUT)
1305 warning(wi, "For efficient BFGS minimization, use switch/shift/pme instead of cut-off.");
1308 if (ir->eI == eiLBFGS && ir->nbfgscorr <= 0)
1310 warning(wi, "Using L-BFGS with nbfgscorr<=0 just gets you steepest descent.");
1313 /* ENERGY CONSERVATION */
1314 if (ir_NVE(ir) && ir->cutoff_scheme == ecutsGROUP)
1316 if (!ir_vdw_might_be_zero_at_cutoff(ir) && ir->rvdw > 0 && ir->vdw_modifier == eintmodNONE)
1318 sprintf(warn_buf, "You are using a cut-off for VdW interactions with NVE, for good energy conservation use vdwtype = %s (possibly with DispCorr)",
1319 evdw_names[evdwSHIFT]);
1320 warning_note(wi, warn_buf);
1322 if (!ir_coulomb_might_be_zero_at_cutoff(ir) && ir->rcoulomb > 0)
1324 sprintf(warn_buf, "You are using a cut-off for electrostatics with NVE, for good energy conservation use coulombtype = %s or %s",
1325 eel_names[eelPMESWITCH], eel_names[eelRF_ZERO]);
1326 warning_note(wi, warn_buf);
1330 if (EI_VV(ir->eI) && IR_TWINRANGE(*ir) && ir->nstlist > 1)
1332 sprintf(warn_buf, "Twin-range multiple time stepping does not work with integrator %s.", ei_names[ir->eI]);
1333 warning_error(wi, warn_buf);
1336 /* IMPLICIT SOLVENT */
1337 if (ir->coulombtype == eelGB_NOTUSED)
1339 sprintf(warn_buf, "Invalid option %s for coulombtype",
1340 eel_names[ir->coulombtype]);
1341 warning_error(wi, warn_buf);
1344 if (ir->sa_algorithm == esaSTILL)
1346 sprintf(err_buf, "Still SA algorithm not available yet, use %s or %s instead\n", esa_names[esaAPPROX], esa_names[esaNO]);
1347 CHECK(ir->sa_algorithm == esaSTILL);
1350 if (ir->implicit_solvent == eisGBSA)
1352 sprintf(err_buf, "With GBSA implicit solvent, rgbradii must be equal to rlist.");
1353 CHECK(ir->rgbradii != ir->rlist);
1355 if (ir->coulombtype != eelCUT)
1357 sprintf(err_buf, "With GBSA, coulombtype must be equal to %s\n", eel_names[eelCUT]);
1358 CHECK(ir->coulombtype != eelCUT);
1360 if (ir->vdwtype != evdwCUT)
1362 sprintf(err_buf, "With GBSA, vdw-type must be equal to %s\n", evdw_names[evdwCUT]);
1363 CHECK(ir->vdwtype != evdwCUT);
1365 if (ir->nstgbradii < 1)
1367 sprintf(warn_buf, "Using GBSA with nstgbradii<1, setting nstgbradii=1");
1368 warning_note(wi, warn_buf);
1371 if (ir->sa_algorithm == esaNO)
1373 sprintf(warn_buf, "No SA (non-polar) calculation requested together with GB. Are you sure this is what you want?\n");
1374 warning_note(wi, warn_buf);
1376 if (ir->sa_surface_tension < 0 && ir->sa_algorithm != esaNO)
1378 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");
1379 warning_note(wi, warn_buf);
1381 if (ir->gb_algorithm == egbSTILL)
1383 ir->sa_surface_tension = 0.0049 * CAL2JOULE * 100;
1387 ir->sa_surface_tension = 0.0054 * CAL2JOULE * 100;
1390 if (ir->sa_surface_tension == 0 && ir->sa_algorithm != esaNO)
1392 sprintf(err_buf, "Surface tension set to 0 while SA-calculation requested\n");
1393 CHECK(ir->sa_surface_tension == 0 && ir->sa_algorithm != esaNO);
1400 if (ir->cutoff_scheme != ecutsGROUP)
1402 warning_error(wi, "AdresS simulation supports only cutoff-scheme=group");
1406 warning_error(wi, "AdresS simulation supports only stochastic dynamics");
1408 if (ir->epc != epcNO)
1410 warning_error(wi, "AdresS simulation does not support pressure coupling");
1412 if (EEL_FULL(ir->coulombtype))
1414 warning_error(wi, "AdresS simulation does not support long-range electrostatics");
1419 /* count the number of text elemets separated by whitespace in a string.
1420 str = the input string
1421 maxptr = the maximum number of allowed elements
1422 ptr = the output array of pointers to the first character of each element
1423 returns: the number of elements. */
1424 int str_nelem(const char *str, int maxptr, char *ptr[])
1429 copy0 = strdup(str);
1432 while (*copy != '\0')
1436 gmx_fatal(FARGS, "Too many groups on line: '%s' (max is %d)",
1444 while ((*copy != '\0') && !isspace(*copy))
1463 /* interpret a number of doubles from a string and put them in an array,
1464 after allocating space for them.
1465 str = the input string
1466 n = the (pre-allocated) number of doubles read
1467 r = the output array of doubles. */
1468 static void parse_n_real(char *str, int *n, real **r)
1473 *n = str_nelem(str, MAXPTR, ptr);
1476 for (i = 0; i < *n; i++)
1478 (*r)[i] = strtod(ptr[i], NULL);
1482 static void do_fep_params(t_inputrec *ir, char fep_lambda[][STRLEN], char weights[STRLEN])
1485 int i, j, max_n_lambda, nweights, nfep[efptNR];
1486 t_lambda *fep = ir->fepvals;
1487 t_expanded *expand = ir->expandedvals;
1488 real **count_fep_lambdas;
1489 gmx_bool bOneLambda = TRUE;
1491 snew(count_fep_lambdas, efptNR);
1493 /* FEP input processing */
1494 /* first, identify the number of lambda values for each type.
1495 All that are nonzero must have the same number */
1497 for (i = 0; i < efptNR; i++)
1499 parse_n_real(fep_lambda[i], &(nfep[i]), &(count_fep_lambdas[i]));
1502 /* now, determine the number of components. All must be either zero, or equal. */
1505 for (i = 0; i < efptNR; i++)
1507 if (nfep[i] > max_n_lambda)
1509 max_n_lambda = nfep[i]; /* here's a nonzero one. All of them
1510 must have the same number if its not zero.*/
1515 for (i = 0; i < efptNR; i++)
1519 ir->fepvals->separate_dvdl[i] = FALSE;
1521 else if (nfep[i] == max_n_lambda)
1523 if (i != efptTEMPERATURE) /* we treat this differently -- not really a reason to compute the derivative with
1524 respect to the temperature currently */
1526 ir->fepvals->separate_dvdl[i] = TRUE;
1531 gmx_fatal(FARGS, "Number of lambdas (%d) for FEP type %s not equal to number of other types (%d)",
1532 nfep[i], efpt_names[i], max_n_lambda);
1535 /* we don't print out dhdl if the temperature is changing, since we can't correctly define dhdl in this case */
1536 ir->fepvals->separate_dvdl[efptTEMPERATURE] = FALSE;
1538 /* the number of lambdas is the number we've read in, which is either zero
1539 or the same for all */
1540 fep->n_lambda = max_n_lambda;
1542 /* allocate space for the array of lambda values */
1543 snew(fep->all_lambda, efptNR);
1544 /* if init_lambda is defined, we need to set lambda */
1545 if ((fep->init_lambda > 0) && (fep->n_lambda == 0))
1547 ir->fepvals->separate_dvdl[efptFEP] = TRUE;
1549 /* otherwise allocate the space for all of the lambdas, and transfer the data */
1550 for (i = 0; i < efptNR; i++)
1552 snew(fep->all_lambda[i], fep->n_lambda);
1553 if (nfep[i] > 0) /* if it's zero, then the count_fep_lambda arrays
1556 for (j = 0; j < fep->n_lambda; j++)
1558 fep->all_lambda[i][j] = (double)count_fep_lambdas[i][j];
1560 sfree(count_fep_lambdas[i]);
1563 sfree(count_fep_lambdas);
1565 /* "fep-vals" is either zero or the full number. If zero, we'll need to define fep-lambdas for internal
1566 bookkeeping -- for now, init_lambda */
1568 if ((nfep[efptFEP] == 0) && (fep->init_lambda >= 0))
1570 for (i = 0; i < fep->n_lambda; i++)
1572 fep->all_lambda[efptFEP][i] = fep->init_lambda;
1576 /* check to see if only a single component lambda is defined, and soft core is defined.
1577 In this case, turn on coulomb soft core */
1579 if (max_n_lambda == 0)
1585 for (i = 0; i < efptNR; i++)
1587 if ((nfep[i] != 0) && (i != efptFEP))
1593 if ((bOneLambda) && (fep->sc_alpha > 0))
1595 fep->bScCoul = TRUE;
1598 /* Fill in the others with the efptFEP if they are not explicitly
1599 specified (i.e. nfep[i] == 0). This means if fep is not defined,
1600 they are all zero. */
1602 for (i = 0; i < efptNR; i++)
1604 if ((nfep[i] == 0) && (i != efptFEP))
1606 for (j = 0; j < fep->n_lambda; j++)
1608 fep->all_lambda[i][j] = fep->all_lambda[efptFEP][j];
1614 /* make it easier if sc_r_power = 48 by increasing it to the 4th power, to be in the right scale. */
1615 if (fep->sc_r_power == 48)
1617 if (fep->sc_alpha > 0.1)
1619 gmx_fatal(FARGS, "sc_alpha (%f) for sc_r_power = 48 should usually be between 0.001 and 0.004", fep->sc_alpha);
1623 expand = ir->expandedvals;
1624 /* now read in the weights */
1625 parse_n_real(weights, &nweights, &(expand->init_lambda_weights));
1628 snew(expand->init_lambda_weights, fep->n_lambda); /* initialize to zero */
1630 else if (nweights != fep->n_lambda)
1632 gmx_fatal(FARGS, "Number of weights (%d) is not equal to number of lambda values (%d)",
1633 nweights, fep->n_lambda);
1635 if ((expand->nstexpanded < 0) && (ir->efep != efepNO))
1637 expand->nstexpanded = fep->nstdhdl;
1638 /* if you don't specify nstexpanded when doing expanded ensemble free energy calcs, it is set to nstdhdl */
1640 if ((expand->nstexpanded < 0) && ir->bSimTemp)
1642 expand->nstexpanded = 2*(int)(ir->opts.tau_t[0]/ir->delta_t);
1643 /* if you don't specify nstexpanded when doing expanded ensemble simulated tempering, it is set to
1644 2*tau_t just to be careful so it's not to frequent */
1649 static void do_simtemp_params(t_inputrec *ir)
1652 snew(ir->simtempvals->temperatures, ir->fepvals->n_lambda);
1653 GetSimTemps(ir->fepvals->n_lambda, ir->simtempvals, ir->fepvals->all_lambda[efptTEMPERATURE]);
1658 static void do_wall_params(t_inputrec *ir,
1659 char *wall_atomtype, char *wall_density,
1663 char *names[MAXPTR];
1666 opts->wall_atomtype[0] = NULL;
1667 opts->wall_atomtype[1] = NULL;
1669 ir->wall_atomtype[0] = -1;
1670 ir->wall_atomtype[1] = -1;
1671 ir->wall_density[0] = 0;
1672 ir->wall_density[1] = 0;
1676 nstr = str_nelem(wall_atomtype, MAXPTR, names);
1677 if (nstr != ir->nwall)
1679 gmx_fatal(FARGS, "Expected %d elements for wall_atomtype, found %d",
1682 for (i = 0; i < ir->nwall; i++)
1684 opts->wall_atomtype[i] = strdup(names[i]);
1687 if (ir->wall_type == ewt93 || ir->wall_type == ewt104)
1689 nstr = str_nelem(wall_density, MAXPTR, names);
1690 if (nstr != ir->nwall)
1692 gmx_fatal(FARGS, "Expected %d elements for wall-density, found %d", ir->nwall, nstr);
1694 for (i = 0; i < ir->nwall; i++)
1696 sscanf(names[i], "%lf", &dbl);
1699 gmx_fatal(FARGS, "wall-density[%d] = %f\n", i, dbl);
1701 ir->wall_density[i] = dbl;
1707 static void add_wall_energrps(gmx_groups_t *groups, int nwall, t_symtab *symtab)
1715 srenew(groups->grpname, groups->ngrpname+nwall);
1716 grps = &(groups->grps[egcENER]);
1717 srenew(grps->nm_ind, grps->nr+nwall);
1718 for (i = 0; i < nwall; i++)
1720 sprintf(str, "wall%d", i);
1721 groups->grpname[groups->ngrpname] = put_symtab(symtab, str);
1722 grps->nm_ind[grps->nr++] = groups->ngrpname++;
1727 void read_expandedparams(int *ninp_p, t_inpfile **inp_p,
1728 t_expanded *expand, warninp_t wi)
1730 int ninp, nerror = 0;
1736 /* read expanded ensemble parameters */
1737 CCTYPE ("expanded ensemble variables");
1738 ITYPE ("nstexpanded", expand->nstexpanded, -1);
1739 EETYPE("lmc-stats", expand->elamstats, elamstats_names);
1740 EETYPE("lmc-move", expand->elmcmove, elmcmove_names);
1741 EETYPE("lmc-weights-equil", expand->elmceq, elmceq_names);
1742 ITYPE ("weight-equil-number-all-lambda", expand->equil_n_at_lam, -1);
1743 ITYPE ("weight-equil-number-samples", expand->equil_samples, -1);
1744 ITYPE ("weight-equil-number-steps", expand->equil_steps, -1);
1745 RTYPE ("weight-equil-wl-delta", expand->equil_wl_delta, -1);
1746 RTYPE ("weight-equil-count-ratio", expand->equil_ratio, -1);
1747 CCTYPE("Seed for Monte Carlo in lambda space");
1748 ITYPE ("lmc-seed", expand->lmc_seed, -1);
1749 RTYPE ("mc-temperature", expand->mc_temp, -1);
1750 ITYPE ("lmc-repeats", expand->lmc_repeats, 1);
1751 ITYPE ("lmc-gibbsdelta", expand->gibbsdeltalam, -1);
1752 ITYPE ("lmc-forced-nstart", expand->lmc_forced_nstart, 0);
1753 EETYPE("symmetrized-transition-matrix", expand->bSymmetrizedTMatrix, yesno_names);
1754 ITYPE("nst-transition-matrix", expand->nstTij, -1);
1755 ITYPE ("mininum-var-min", expand->minvarmin, 100); /*default is reasonable */
1756 ITYPE ("weight-c-range", expand->c_range, 0); /* default is just C=0 */
1757 RTYPE ("wl-scale", expand->wl_scale, 0.8);
1758 RTYPE ("wl-ratio", expand->wl_ratio, 0.8);
1759 RTYPE ("init-wl-delta", expand->init_wl_delta, 1.0);
1760 EETYPE("wl-oneovert", expand->bWLoneovert, yesno_names);
1768 void get_ir(const char *mdparin, const char *mdparout,
1769 t_inputrec *ir, t_gromppopts *opts,
1773 double dumdub[2][6];
1777 char warn_buf[STRLEN];
1778 t_lambda *fep = ir->fepvals;
1779 t_expanded *expand = ir->expandedvals;
1781 init_inputrec_strings();
1782 inp = read_inpfile(mdparin, &ninp, wi);
1784 snew(dumstr[0], STRLEN);
1785 snew(dumstr[1], STRLEN);
1787 if (-1 == search_einp(ninp, inp, "cutoff-scheme"))
1790 "%s did not specify a value for the .mdp option "
1791 "\"cutoff-scheme\". Probably it was first intended for use "
1792 "with GROMACS before 4.6. In 4.6, the Verlet scheme was "
1793 "introduced, but the group scheme was still the default. "
1794 "The default is now the Verlet scheme, so you will observe "
1795 "different behaviour.", mdparin);
1796 warning_note(wi, warn_buf);
1799 /* ignore the following deprecated commands */
1802 REM_TYPE("domain-decomposition");
1803 REM_TYPE("andersen-seed");
1805 REM_TYPE("dihre-fc");
1806 REM_TYPE("dihre-tau");
1807 REM_TYPE("nstdihreout");
1808 REM_TYPE("nstcheckpoint");
1809 REM_TYPE("optimize-fft");
1811 /* replace the following commands with the clearer new versions*/
1812 REPL_TYPE("unconstrained-start", "continuation");
1813 REPL_TYPE("foreign-lambda", "fep-lambdas");
1814 REPL_TYPE("verlet-buffer-drift", "verlet-buffer-tolerance");
1815 REPL_TYPE("nstxtcout", "nstxout-compressed");
1816 REPL_TYPE("xtc-grps", "compressed-x-grps");
1817 REPL_TYPE("xtc-precision", "compressed-x-precision");
1819 CCTYPE ("VARIOUS PREPROCESSING OPTIONS");
1820 CTYPE ("Preprocessor information: use cpp syntax.");
1821 CTYPE ("e.g.: -I/home/joe/doe -I/home/mary/roe");
1822 STYPE ("include", opts->include, NULL);
1823 CTYPE ("e.g.: -DPOSRES -DFLEXIBLE (note these variable names are case sensitive)");
1824 STYPE ("define", opts->define, NULL);
1826 CCTYPE ("RUN CONTROL PARAMETERS");
1827 EETYPE("integrator", ir->eI, ei_names);
1828 CTYPE ("Start time and timestep in ps");
1829 RTYPE ("tinit", ir->init_t, 0.0);
1830 RTYPE ("dt", ir->delta_t, 0.001);
1831 STEPTYPE ("nsteps", ir->nsteps, 0);
1832 CTYPE ("For exact run continuation or redoing part of a run");
1833 STEPTYPE ("init-step", ir->init_step, 0);
1834 CTYPE ("Part index is updated automatically on checkpointing (keeps files separate)");
1835 ITYPE ("simulation-part", ir->simulation_part, 1);
1836 CTYPE ("mode for center of mass motion removal");
1837 EETYPE("comm-mode", ir->comm_mode, ecm_names);
1838 CTYPE ("number of steps for center of mass motion removal");
1839 ITYPE ("nstcomm", ir->nstcomm, 100);
1840 CTYPE ("group(s) for center of mass motion removal");
1841 STYPE ("comm-grps", is->vcm, NULL);
1843 CCTYPE ("LANGEVIN DYNAMICS OPTIONS");
1844 CTYPE ("Friction coefficient (amu/ps) and random seed");
1845 RTYPE ("bd-fric", ir->bd_fric, 0.0);
1846 STEPTYPE ("ld-seed", ir->ld_seed, -1);
1849 CCTYPE ("ENERGY MINIMIZATION OPTIONS");
1850 CTYPE ("Force tolerance and initial step-size");
1851 RTYPE ("emtol", ir->em_tol, 10.0);
1852 RTYPE ("emstep", ir->em_stepsize, 0.01);
1853 CTYPE ("Max number of iterations in relax-shells");
1854 ITYPE ("niter", ir->niter, 20);
1855 CTYPE ("Step size (ps^2) for minimization of flexible constraints");
1856 RTYPE ("fcstep", ir->fc_stepsize, 0);
1857 CTYPE ("Frequency of steepest descents steps when doing CG");
1858 ITYPE ("nstcgsteep", ir->nstcgsteep, 1000);
1859 ITYPE ("nbfgscorr", ir->nbfgscorr, 10);
1861 CCTYPE ("TEST PARTICLE INSERTION OPTIONS");
1862 RTYPE ("rtpi", ir->rtpi, 0.05);
1864 /* Output options */
1865 CCTYPE ("OUTPUT CONTROL OPTIONS");
1866 CTYPE ("Output frequency for coords (x), velocities (v) and forces (f)");
1867 ITYPE ("nstxout", ir->nstxout, 0);
1868 ITYPE ("nstvout", ir->nstvout, 0);
1869 ITYPE ("nstfout", ir->nstfout, 0);
1870 CTYPE ("Output frequency for energies to log file and energy file");
1871 ITYPE ("nstlog", ir->nstlog, 1000);
1872 ITYPE ("nstcalcenergy", ir->nstcalcenergy, 100);
1873 ITYPE ("nstenergy", ir->nstenergy, 1000);
1874 CTYPE ("Output frequency and precision for .xtc file");
1875 ITYPE ("nstxout-compressed", ir->nstxout_compressed, 0);
1876 RTYPE ("compressed-x-precision", ir->x_compression_precision, 1000.0);
1877 CTYPE ("This selects the subset of atoms for the compressed");
1878 CTYPE ("trajectory file. You can select multiple groups. By");
1879 CTYPE ("default, all atoms will be written.");
1880 STYPE ("compressed-x-grps", is->x_compressed_groups, NULL);
1881 CTYPE ("Selection of energy groups");
1882 STYPE ("energygrps", is->energy, NULL);
1884 /* Neighbor searching */
1885 CCTYPE ("NEIGHBORSEARCHING PARAMETERS");
1886 CTYPE ("cut-off scheme (Verlet: particle based cut-offs, group: using charge groups)");
1887 EETYPE("cutoff-scheme", ir->cutoff_scheme, ecutscheme_names);
1888 CTYPE ("nblist update frequency");
1889 ITYPE ("nstlist", ir->nstlist, 10);
1890 CTYPE ("ns algorithm (simple or grid)");
1891 EETYPE("ns-type", ir->ns_type, ens_names);
1892 CTYPE ("Periodic boundary conditions: xyz, no, xy");
1893 EETYPE("pbc", ir->ePBC, epbc_names);
1894 EETYPE("periodic-molecules", ir->bPeriodicMols, yesno_names);
1895 CTYPE ("Allowed energy error due to the Verlet buffer in kJ/mol/ps per atom,");
1896 CTYPE ("a value of -1 means: use rlist");
1897 RTYPE("verlet-buffer-tolerance", ir->verletbuf_tol, 0.005);
1898 CTYPE ("nblist cut-off");
1899 RTYPE ("rlist", ir->rlist, 1.0);
1900 CTYPE ("long-range cut-off for switched potentials");
1901 RTYPE ("rlistlong", ir->rlistlong, -1);
1902 ITYPE ("nstcalclr", ir->nstcalclr, -1);
1904 /* Electrostatics */
1905 CCTYPE ("OPTIONS FOR ELECTROSTATICS AND VDW");
1906 CTYPE ("Method for doing electrostatics");
1907 EETYPE("coulombtype", ir->coulombtype, eel_names);
1908 EETYPE("coulomb-modifier", ir->coulomb_modifier, eintmod_names);
1909 CTYPE ("cut-off lengths");
1910 RTYPE ("rcoulomb-switch", ir->rcoulomb_switch, 0.0);
1911 RTYPE ("rcoulomb", ir->rcoulomb, 1.0);
1912 CTYPE ("Relative dielectric constant for the medium and the reaction field");
1913 RTYPE ("epsilon-r", ir->epsilon_r, 1.0);
1914 RTYPE ("epsilon-rf", ir->epsilon_rf, 0.0);
1915 CTYPE ("Method for doing Van der Waals");
1916 EETYPE("vdw-type", ir->vdwtype, evdw_names);
1917 EETYPE("vdw-modifier", ir->vdw_modifier, eintmod_names);
1918 CTYPE ("cut-off lengths");
1919 RTYPE ("rvdw-switch", ir->rvdw_switch, 0.0);
1920 RTYPE ("rvdw", ir->rvdw, 1.0);
1921 CTYPE ("Apply long range dispersion corrections for Energy and Pressure");
1922 EETYPE("DispCorr", ir->eDispCorr, edispc_names);
1923 CTYPE ("Extension of the potential lookup tables beyond the cut-off");
1924 RTYPE ("table-extension", ir->tabext, 1.0);
1925 CTYPE ("Separate tables between energy group pairs");
1926 STYPE ("energygrp-table", is->egptable, NULL);
1927 CTYPE ("Spacing for the PME/PPPM FFT grid");
1928 RTYPE ("fourierspacing", ir->fourier_spacing, 0.12);
1929 CTYPE ("FFT grid size, when a value is 0 fourierspacing will be used");
1930 ITYPE ("fourier-nx", ir->nkx, 0);
1931 ITYPE ("fourier-ny", ir->nky, 0);
1932 ITYPE ("fourier-nz", ir->nkz, 0);
1933 CTYPE ("EWALD/PME/PPPM parameters");
1934 ITYPE ("pme-order", ir->pme_order, 4);
1935 RTYPE ("ewald-rtol", ir->ewald_rtol, 0.00001);
1936 RTYPE ("ewald-rtol-lj", ir->ewald_rtol_lj, 0.001);
1937 EETYPE("lj-pme-comb-rule", ir->ljpme_combination_rule, eljpme_names);
1938 EETYPE("ewald-geometry", ir->ewald_geometry, eewg_names);
1939 RTYPE ("epsilon-surface", ir->epsilon_surface, 0.0);
1941 CCTYPE("IMPLICIT SOLVENT ALGORITHM");
1942 EETYPE("implicit-solvent", ir->implicit_solvent, eis_names);
1944 CCTYPE ("GENERALIZED BORN ELECTROSTATICS");
1945 CTYPE ("Algorithm for calculating Born radii");
1946 EETYPE("gb-algorithm", ir->gb_algorithm, egb_names);
1947 CTYPE ("Frequency of calculating the Born radii inside rlist");
1948 ITYPE ("nstgbradii", ir->nstgbradii, 1);
1949 CTYPE ("Cutoff for Born radii calculation; the contribution from atoms");
1950 CTYPE ("between rlist and rgbradii is updated every nstlist steps");
1951 RTYPE ("rgbradii", ir->rgbradii, 1.0);
1952 CTYPE ("Dielectric coefficient of the implicit solvent");
1953 RTYPE ("gb-epsilon-solvent", ir->gb_epsilon_solvent, 80.0);
1954 CTYPE ("Salt concentration in M for Generalized Born models");
1955 RTYPE ("gb-saltconc", ir->gb_saltconc, 0.0);
1956 CTYPE ("Scaling factors used in the OBC GB model. Default values are OBC(II)");
1957 RTYPE ("gb-obc-alpha", ir->gb_obc_alpha, 1.0);
1958 RTYPE ("gb-obc-beta", ir->gb_obc_beta, 0.8);
1959 RTYPE ("gb-obc-gamma", ir->gb_obc_gamma, 4.85);
1960 RTYPE ("gb-dielectric-offset", ir->gb_dielectric_offset, 0.009);
1961 EETYPE("sa-algorithm", ir->sa_algorithm, esa_names);
1962 CTYPE ("Surface tension (kJ/mol/nm^2) for the SA (nonpolar surface) part of GBSA");
1963 CTYPE ("The value -1 will set default value for Still/HCT/OBC GB-models.");
1964 RTYPE ("sa-surface-tension", ir->sa_surface_tension, -1);
1966 /* Coupling stuff */
1967 CCTYPE ("OPTIONS FOR WEAK COUPLING ALGORITHMS");
1968 CTYPE ("Temperature coupling");
1969 EETYPE("tcoupl", ir->etc, etcoupl_names);
1970 ITYPE ("nsttcouple", ir->nsttcouple, -1);
1971 ITYPE("nh-chain-length", ir->opts.nhchainlength, 10);
1972 EETYPE("print-nose-hoover-chain-variables", ir->bPrintNHChains, yesno_names);
1973 CTYPE ("Groups to couple separately");
1974 STYPE ("tc-grps", is->tcgrps, NULL);
1975 CTYPE ("Time constant (ps) and reference temperature (K)");
1976 STYPE ("tau-t", is->tau_t, NULL);
1977 STYPE ("ref-t", is->ref_t, NULL);
1978 CTYPE ("pressure coupling");
1979 EETYPE("pcoupl", ir->epc, epcoupl_names);
1980 EETYPE("pcoupltype", ir->epct, epcoupltype_names);
1981 ITYPE ("nstpcouple", ir->nstpcouple, -1);
1982 CTYPE ("Time constant (ps), compressibility (1/bar) and reference P (bar)");
1983 RTYPE ("tau-p", ir->tau_p, 1.0);
1984 STYPE ("compressibility", dumstr[0], NULL);
1985 STYPE ("ref-p", dumstr[1], NULL);
1986 CTYPE ("Scaling of reference coordinates, No, All or COM");
1987 EETYPE ("refcoord-scaling", ir->refcoord_scaling, erefscaling_names);
1990 CCTYPE ("OPTIONS FOR QMMM calculations");
1991 EETYPE("QMMM", ir->bQMMM, yesno_names);
1992 CTYPE ("Groups treated Quantum Mechanically");
1993 STYPE ("QMMM-grps", is->QMMM, NULL);
1994 CTYPE ("QM method");
1995 STYPE("QMmethod", is->QMmethod, NULL);
1996 CTYPE ("QMMM scheme");
1997 EETYPE("QMMMscheme", ir->QMMMscheme, eQMMMscheme_names);
1998 CTYPE ("QM basisset");
1999 STYPE("QMbasis", is->QMbasis, NULL);
2000 CTYPE ("QM charge");
2001 STYPE ("QMcharge", is->QMcharge, NULL);
2002 CTYPE ("QM multiplicity");
2003 STYPE ("QMmult", is->QMmult, NULL);
2004 CTYPE ("Surface Hopping");
2005 STYPE ("SH", is->bSH, NULL);
2006 CTYPE ("CAS space options");
2007 STYPE ("CASorbitals", is->CASorbitals, NULL);
2008 STYPE ("CASelectrons", is->CASelectrons, NULL);
2009 STYPE ("SAon", is->SAon, NULL);
2010 STYPE ("SAoff", is->SAoff, NULL);
2011 STYPE ("SAsteps", is->SAsteps, NULL);
2012 CTYPE ("Scale factor for MM charges");
2013 RTYPE ("MMChargeScaleFactor", ir->scalefactor, 1.0);
2014 CTYPE ("Optimization of QM subsystem");
2015 STYPE ("bOPT", is->bOPT, NULL);
2016 STYPE ("bTS", is->bTS, NULL);
2018 /* Simulated annealing */
2019 CCTYPE("SIMULATED ANNEALING");
2020 CTYPE ("Type of annealing for each temperature group (no/single/periodic)");
2021 STYPE ("annealing", is->anneal, NULL);
2022 CTYPE ("Number of time points to use for specifying annealing in each group");
2023 STYPE ("annealing-npoints", is->anneal_npoints, NULL);
2024 CTYPE ("List of times at the annealing points for each group");
2025 STYPE ("annealing-time", is->anneal_time, NULL);
2026 CTYPE ("Temp. at each annealing point, for each group.");
2027 STYPE ("annealing-temp", is->anneal_temp, NULL);
2030 CCTYPE ("GENERATE VELOCITIES FOR STARTUP RUN");
2031 EETYPE("gen-vel", opts->bGenVel, yesno_names);
2032 RTYPE ("gen-temp", opts->tempi, 300.0);
2033 ITYPE ("gen-seed", opts->seed, -1);
2036 CCTYPE ("OPTIONS FOR BONDS");
2037 EETYPE("constraints", opts->nshake, constraints);
2038 CTYPE ("Type of constraint algorithm");
2039 EETYPE("constraint-algorithm", ir->eConstrAlg, econstr_names);
2040 CTYPE ("Do not constrain the start configuration");
2041 EETYPE("continuation", ir->bContinuation, yesno_names);
2042 CTYPE ("Use successive overrelaxation to reduce the number of shake iterations");
2043 EETYPE("Shake-SOR", ir->bShakeSOR, yesno_names);
2044 CTYPE ("Relative tolerance of shake");
2045 RTYPE ("shake-tol", ir->shake_tol, 0.0001);
2046 CTYPE ("Highest order in the expansion of the constraint coupling matrix");
2047 ITYPE ("lincs-order", ir->nProjOrder, 4);
2048 CTYPE ("Number of iterations in the final step of LINCS. 1 is fine for");
2049 CTYPE ("normal simulations, but use 2 to conserve energy in NVE runs.");
2050 CTYPE ("For energy minimization with constraints it should be 4 to 8.");
2051 ITYPE ("lincs-iter", ir->nLincsIter, 1);
2052 CTYPE ("Lincs will write a warning to the stderr if in one step a bond");
2053 CTYPE ("rotates over more degrees than");
2054 RTYPE ("lincs-warnangle", ir->LincsWarnAngle, 30.0);
2055 CTYPE ("Convert harmonic bonds to morse potentials");
2056 EETYPE("morse", opts->bMorse, yesno_names);
2058 /* Energy group exclusions */
2059 CCTYPE ("ENERGY GROUP EXCLUSIONS");
2060 CTYPE ("Pairs of energy groups for which all non-bonded interactions are excluded");
2061 STYPE ("energygrp-excl", is->egpexcl, NULL);
2065 CTYPE ("Number of walls, type, atom types, densities and box-z scale factor for Ewald");
2066 ITYPE ("nwall", ir->nwall, 0);
2067 EETYPE("wall-type", ir->wall_type, ewt_names);
2068 RTYPE ("wall-r-linpot", ir->wall_r_linpot, -1);
2069 STYPE ("wall-atomtype", is->wall_atomtype, NULL);
2070 STYPE ("wall-density", is->wall_density, NULL);
2071 RTYPE ("wall-ewald-zfac", ir->wall_ewald_zfac, 3);
2074 CCTYPE("COM PULLING");
2075 CTYPE("Pull type: no, umbrella, constraint or constant-force");
2076 EETYPE("pull", ir->ePull, epull_names);
2077 if (ir->ePull != epullNO)
2080 is->pull_grp = read_pullparams(&ninp, &inp, ir->pull, &opts->pull_start, wi);
2083 /* Enforced rotation */
2084 CCTYPE("ENFORCED ROTATION");
2085 CTYPE("Enforced rotation: No or Yes");
2086 EETYPE("rotation", ir->bRot, yesno_names);
2090 is->rot_grp = read_rotparams(&ninp, &inp, ir->rot, wi);
2093 /* Interactive MD */
2095 CCTYPE("Group to display and/or manipulate in interactive MD session");
2096 STYPE ("IMD-group", is->imd_grp, NULL);
2097 if (is->imd_grp[0] != '\0')
2104 CCTYPE("NMR refinement stuff");
2105 CTYPE ("Distance restraints type: No, Simple or Ensemble");
2106 EETYPE("disre", ir->eDisre, edisre_names);
2107 CTYPE ("Force weighting of pairs in one distance restraint: Conservative or Equal");
2108 EETYPE("disre-weighting", ir->eDisreWeighting, edisreweighting_names);
2109 CTYPE ("Use sqrt of the time averaged times the instantaneous violation");
2110 EETYPE("disre-mixed", ir->bDisreMixed, yesno_names);
2111 RTYPE ("disre-fc", ir->dr_fc, 1000.0);
2112 RTYPE ("disre-tau", ir->dr_tau, 0.0);
2113 CTYPE ("Output frequency for pair distances to energy file");
2114 ITYPE ("nstdisreout", ir->nstdisreout, 100);
2115 CTYPE ("Orientation restraints: No or Yes");
2116 EETYPE("orire", opts->bOrire, yesno_names);
2117 CTYPE ("Orientation restraints force constant and tau for time averaging");
2118 RTYPE ("orire-fc", ir->orires_fc, 0.0);
2119 RTYPE ("orire-tau", ir->orires_tau, 0.0);
2120 STYPE ("orire-fitgrp", is->orirefitgrp, NULL);
2121 CTYPE ("Output frequency for trace(SD) and S to energy file");
2122 ITYPE ("nstorireout", ir->nstorireout, 100);
2124 /* free energy variables */
2125 CCTYPE ("Free energy variables");
2126 EETYPE("free-energy", ir->efep, efep_names);
2127 STYPE ("couple-moltype", is->couple_moltype, NULL);
2128 EETYPE("couple-lambda0", opts->couple_lam0, couple_lam);
2129 EETYPE("couple-lambda1", opts->couple_lam1, couple_lam);
2130 EETYPE("couple-intramol", opts->bCoupleIntra, yesno_names);
2132 RTYPE ("init-lambda", fep->init_lambda, -1); /* start with -1 so
2134 it was not entered */
2135 ITYPE ("init-lambda-state", fep->init_fep_state, -1);
2136 RTYPE ("delta-lambda", fep->delta_lambda, 0.0);
2137 ITYPE ("nstdhdl", fep->nstdhdl, 50);
2138 STYPE ("fep-lambdas", is->fep_lambda[efptFEP], NULL);
2139 STYPE ("mass-lambdas", is->fep_lambda[efptMASS], NULL);
2140 STYPE ("coul-lambdas", is->fep_lambda[efptCOUL], NULL);
2141 STYPE ("vdw-lambdas", is->fep_lambda[efptVDW], NULL);
2142 STYPE ("bonded-lambdas", is->fep_lambda[efptBONDED], NULL);
2143 STYPE ("restraint-lambdas", is->fep_lambda[efptRESTRAINT], NULL);
2144 STYPE ("temperature-lambdas", is->fep_lambda[efptTEMPERATURE], NULL);
2145 ITYPE ("calc-lambda-neighbors", fep->lambda_neighbors, 1);
2146 STYPE ("init-lambda-weights", is->lambda_weights, NULL);
2147 EETYPE("dhdl-print-energy", fep->edHdLPrintEnergy, edHdLPrintEnergy_names);
2148 RTYPE ("sc-alpha", fep->sc_alpha, 0.0);
2149 ITYPE ("sc-power", fep->sc_power, 1);
2150 RTYPE ("sc-r-power", fep->sc_r_power, 6.0);
2151 RTYPE ("sc-sigma", fep->sc_sigma, 0.3);
2152 EETYPE("sc-coul", fep->bScCoul, yesno_names);
2153 ITYPE ("dh_hist_size", fep->dh_hist_size, 0);
2154 RTYPE ("dh_hist_spacing", fep->dh_hist_spacing, 0.1);
2155 EETYPE("separate-dhdl-file", fep->separate_dhdl_file,
2156 separate_dhdl_file_names);
2157 EETYPE("dhdl-derivatives", fep->dhdl_derivatives, dhdl_derivatives_names);
2158 ITYPE ("dh_hist_size", fep->dh_hist_size, 0);
2159 RTYPE ("dh_hist_spacing", fep->dh_hist_spacing, 0.1);
2161 /* Non-equilibrium MD stuff */
2162 CCTYPE("Non-equilibrium MD stuff");
2163 STYPE ("acc-grps", is->accgrps, NULL);
2164 STYPE ("accelerate", is->acc, NULL);
2165 STYPE ("freezegrps", is->freeze, NULL);
2166 STYPE ("freezedim", is->frdim, NULL);
2167 RTYPE ("cos-acceleration", ir->cos_accel, 0);
2168 STYPE ("deform", is->deform, NULL);
2170 /* simulated tempering variables */
2171 CCTYPE("simulated tempering variables");
2172 EETYPE("simulated-tempering", ir->bSimTemp, yesno_names);
2173 EETYPE("simulated-tempering-scaling", ir->simtempvals->eSimTempScale, esimtemp_names);
2174 RTYPE("sim-temp-low", ir->simtempvals->simtemp_low, 300.0);
2175 RTYPE("sim-temp-high", ir->simtempvals->simtemp_high, 300.0);
2177 /* expanded ensemble variables */
2178 if (ir->efep == efepEXPANDED || ir->bSimTemp)
2180 read_expandedparams(&ninp, &inp, expand, wi);
2183 /* Electric fields */
2184 CCTYPE("Electric fields");
2185 CTYPE ("Format is number of terms (int) and for all terms an amplitude (real)");
2186 CTYPE ("and a phase angle (real)");
2187 STYPE ("E-x", is->efield_x, NULL);
2188 STYPE ("E-xt", is->efield_xt, NULL);
2189 STYPE ("E-y", is->efield_y, NULL);
2190 STYPE ("E-yt", is->efield_yt, NULL);
2191 STYPE ("E-z", is->efield_z, NULL);
2192 STYPE ("E-zt", is->efield_zt, NULL);
2194 CCTYPE("Ion/water position swapping for computational electrophysiology setups");
2195 CTYPE("Swap positions along direction: no, X, Y, Z");
2196 EETYPE("swapcoords", ir->eSwapCoords, eSwapTypes_names);
2197 if (ir->eSwapCoords != eswapNO)
2200 CTYPE("Swap attempt frequency");
2201 ITYPE("swap-frequency", ir->swap->nstswap, 1);
2202 CTYPE("Two index groups that contain the compartment-partitioning atoms");
2203 STYPE("split-group0", splitgrp0, NULL);
2204 STYPE("split-group1", splitgrp1, NULL);
2205 CTYPE("Use center of mass of split groups (yes/no), otherwise center of geometry is used");
2206 EETYPE("massw-split0", ir->swap->massw_split[0], yesno_names);
2207 EETYPE("massw-split1", ir->swap->massw_split[1], yesno_names);
2209 CTYPE("Group name of ions that can be exchanged with solvent molecules");
2210 STYPE("swap-group", swapgrp, NULL);
2211 CTYPE("Group name of solvent molecules");
2212 STYPE("solvent-group", solgrp, NULL);
2214 CTYPE("Split cylinder: radius, upper and lower extension (nm) (this will define the channels)");
2215 CTYPE("Note that the split cylinder settings do not have an influence on the swapping protocol,");
2216 CTYPE("however, if correctly defined, the ion permeation events are counted per channel");
2217 RTYPE("cyl0-r", ir->swap->cyl0r, 2.0);
2218 RTYPE("cyl0-up", ir->swap->cyl0u, 1.0);
2219 RTYPE("cyl0-down", ir->swap->cyl0l, 1.0);
2220 RTYPE("cyl1-r", ir->swap->cyl1r, 2.0);
2221 RTYPE("cyl1-up", ir->swap->cyl1u, 1.0);
2222 RTYPE("cyl1-down", ir->swap->cyl1l, 1.0);
2224 CTYPE("Average the number of ions per compartment over these many swap attempt steps");
2225 ITYPE("coupl-steps", ir->swap->nAverage, 10);
2226 CTYPE("Requested number of anions and cations for each of the two compartments");
2227 CTYPE("-1 means fix the numbers as found in time step 0");
2228 ITYPE("anionsA", ir->swap->nanions[0], -1);
2229 ITYPE("cationsA", ir->swap->ncations[0], -1);
2230 ITYPE("anionsB", ir->swap->nanions[1], -1);
2231 ITYPE("cationsB", ir->swap->ncations[1], -1);
2232 CTYPE("Start to swap ions if threshold difference to requested count is reached");
2233 RTYPE("threshold", ir->swap->threshold, 1.0);
2236 /* AdResS defined thingies */
2237 CCTYPE ("AdResS parameters");
2238 EETYPE("adress", ir->bAdress, yesno_names);
2241 snew(ir->adress, 1);
2242 read_adressparams(&ninp, &inp, ir->adress, wi);
2245 /* User defined thingies */
2246 CCTYPE ("User defined thingies");
2247 STYPE ("user1-grps", is->user1, NULL);
2248 STYPE ("user2-grps", is->user2, NULL);
2249 ITYPE ("userint1", ir->userint1, 0);
2250 ITYPE ("userint2", ir->userint2, 0);
2251 ITYPE ("userint3", ir->userint3, 0);
2252 ITYPE ("userint4", ir->userint4, 0);
2253 RTYPE ("userreal1", ir->userreal1, 0);
2254 RTYPE ("userreal2", ir->userreal2, 0);
2255 RTYPE ("userreal3", ir->userreal3, 0);
2256 RTYPE ("userreal4", ir->userreal4, 0);
2259 write_inpfile(mdparout, ninp, inp, FALSE, wi);
2260 for (i = 0; (i < ninp); i++)
2263 sfree(inp[i].value);
2267 /* Process options if necessary */
2268 for (m = 0; m < 2; m++)
2270 for (i = 0; i < 2*DIM; i++)
2279 if (sscanf(dumstr[m], "%lf", &(dumdub[m][XX])) != 1)
2281 warning_error(wi, "Pressure coupling not enough values (I need 1)");
2283 dumdub[m][YY] = dumdub[m][ZZ] = dumdub[m][XX];
2285 case epctSEMIISOTROPIC:
2286 case epctSURFACETENSION:
2287 if (sscanf(dumstr[m], "%lf%lf",
2288 &(dumdub[m][XX]), &(dumdub[m][ZZ])) != 2)
2290 warning_error(wi, "Pressure coupling not enough values (I need 2)");
2292 dumdub[m][YY] = dumdub[m][XX];
2294 case epctANISOTROPIC:
2295 if (sscanf(dumstr[m], "%lf%lf%lf%lf%lf%lf",
2296 &(dumdub[m][XX]), &(dumdub[m][YY]), &(dumdub[m][ZZ]),
2297 &(dumdub[m][3]), &(dumdub[m][4]), &(dumdub[m][5])) != 6)
2299 warning_error(wi, "Pressure coupling not enough values (I need 6)");
2303 gmx_fatal(FARGS, "Pressure coupling type %s not implemented yet",
2304 epcoupltype_names[ir->epct]);
2308 clear_mat(ir->ref_p);
2309 clear_mat(ir->compress);
2310 for (i = 0; i < DIM; i++)
2312 ir->ref_p[i][i] = dumdub[1][i];
2313 ir->compress[i][i] = dumdub[0][i];
2315 if (ir->epct == epctANISOTROPIC)
2317 ir->ref_p[XX][YY] = dumdub[1][3];
2318 ir->ref_p[XX][ZZ] = dumdub[1][4];
2319 ir->ref_p[YY][ZZ] = dumdub[1][5];
2320 if (ir->ref_p[XX][YY] != 0 && ir->ref_p[XX][ZZ] != 0 && ir->ref_p[YY][ZZ] != 0)
2322 warning(wi, "All off-diagonal reference pressures are non-zero. Are you sure you want to apply a threefold shear stress?\n");
2324 ir->compress[XX][YY] = dumdub[0][3];
2325 ir->compress[XX][ZZ] = dumdub[0][4];
2326 ir->compress[YY][ZZ] = dumdub[0][5];
2327 for (i = 0; i < DIM; i++)
2329 for (m = 0; m < i; m++)
2331 ir->ref_p[i][m] = ir->ref_p[m][i];
2332 ir->compress[i][m] = ir->compress[m][i];
2337 if (ir->comm_mode == ecmNO)
2342 opts->couple_moltype = NULL;
2343 if (strlen(is->couple_moltype) > 0)
2345 if (ir->efep != efepNO)
2347 opts->couple_moltype = strdup(is->couple_moltype);
2348 if (opts->couple_lam0 == opts->couple_lam1)
2350 warning(wi, "The lambda=0 and lambda=1 states for coupling are identical");
2352 if (ir->eI == eiMD && (opts->couple_lam0 == ecouplamNONE ||
2353 opts->couple_lam1 == ecouplamNONE))
2355 warning(wi, "For proper sampling of the (nearly) decoupled state, stochastic dynamics should be used");
2360 warning_note(wi, "Free energy is turned off, so we will not decouple the molecule listed in your input.");
2363 /* FREE ENERGY AND EXPANDED ENSEMBLE OPTIONS */
2364 if (ir->efep != efepNO)
2366 if (fep->delta_lambda > 0)
2368 ir->efep = efepSLOWGROWTH;
2372 if (fep->edHdLPrintEnergy == edHdLPrintEnergyYES)
2374 fep->edHdLPrintEnergy = edHdLPrintEnergyTOTAL;
2375 warning_note(wi, "Old option for dhdl-print-energy given: "
2376 "changing \"yes\" to \"total\"\n");
2379 if (ir->bSimTemp && (fep->edHdLPrintEnergy == edHdLPrintEnergyNO))
2381 /* always print out the energy to dhdl if we are doing
2382 expanded ensemble, since we need the total energy for
2383 analysis if the temperature is changing. In some
2384 conditions one may only want the potential energy, so
2385 we will allow that if the appropriate mdp setting has
2386 been enabled. Otherwise, total it is:
2388 fep->edHdLPrintEnergy = edHdLPrintEnergyTOTAL;
2391 if ((ir->efep != efepNO) || ir->bSimTemp)
2393 ir->bExpanded = FALSE;
2394 if ((ir->efep == efepEXPANDED) || ir->bSimTemp)
2396 ir->bExpanded = TRUE;
2398 do_fep_params(ir, is->fep_lambda, is->lambda_weights);
2399 if (ir->bSimTemp) /* done after fep params */
2401 do_simtemp_params(ir);
2406 ir->fepvals->n_lambda = 0;
2409 /* WALL PARAMETERS */
2411 do_wall_params(ir, is->wall_atomtype, is->wall_density, opts);
2413 /* ORIENTATION RESTRAINT PARAMETERS */
2415 if (opts->bOrire && str_nelem(is->orirefitgrp, MAXPTR, NULL) != 1)
2417 warning_error(wi, "ERROR: Need one orientation restraint fit group\n");
2420 /* DEFORMATION PARAMETERS */
2422 clear_mat(ir->deform);
2423 for (i = 0; i < 6; i++)
2427 m = sscanf(is->deform, "%lf %lf %lf %lf %lf %lf",
2428 &(dumdub[0][0]), &(dumdub[0][1]), &(dumdub[0][2]),
2429 &(dumdub[0][3]), &(dumdub[0][4]), &(dumdub[0][5]));
2430 for (i = 0; i < 3; i++)
2432 ir->deform[i][i] = dumdub[0][i];
2434 ir->deform[YY][XX] = dumdub[0][3];
2435 ir->deform[ZZ][XX] = dumdub[0][4];
2436 ir->deform[ZZ][YY] = dumdub[0][5];
2437 if (ir->epc != epcNO)
2439 for (i = 0; i < 3; i++)
2441 for (j = 0; j <= i; j++)
2443 if (ir->deform[i][j] != 0 && ir->compress[i][j] != 0)
2445 warning_error(wi, "A box element has deform set and compressibility > 0");
2449 for (i = 0; i < 3; i++)
2451 for (j = 0; j < i; j++)
2453 if (ir->deform[i][j] != 0)
2455 for (m = j; m < DIM; m++)
2457 if (ir->compress[m][j] != 0)
2459 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.");
2460 warning(wi, warn_buf);
2468 /* Ion/water position swapping checks */
2469 if (ir->eSwapCoords != eswapNO)
2471 if (ir->swap->nstswap < 1)
2473 warning_error(wi, "swap_frequency must be 1 or larger when ion swapping is requested");
2475 if (ir->swap->nAverage < 1)
2477 warning_error(wi, "coupl_steps must be 1 or larger.\n");
2479 if (ir->swap->threshold < 1.0)
2481 warning_error(wi, "Ion count threshold must be at least 1.\n");
2489 static int search_QMstring(const char *s, int ng, const char *gn[])
2491 /* same as normal search_string, but this one searches QM strings */
2494 for (i = 0; (i < ng); i++)
2496 if (gmx_strcasecmp(s, gn[i]) == 0)
2502 gmx_fatal(FARGS, "this QM method or basisset (%s) is not implemented\n!", s);
2506 } /* search_QMstring */
2508 /* We would like gn to be const as well, but C doesn't allow this */
2509 int search_string(const char *s, int ng, char *gn[])
2513 for (i = 0; (i < ng); i++)
2515 if (gmx_strcasecmp(s, gn[i]) == 0)
2522 "Group %s referenced in the .mdp file was not found in the index file.\n"
2523 "Group names must match either [moleculetype] names or custom index group\n"
2524 "names, in which case you must supply an index file to the '-n' option\n"
2531 static gmx_bool do_numbering(int natoms, gmx_groups_t *groups, int ng, char *ptrs[],
2532 t_blocka *block, char *gnames[],
2533 int gtype, int restnm,
2534 int grptp, gmx_bool bVerbose,
2537 unsigned short *cbuf;
2538 t_grps *grps = &(groups->grps[gtype]);
2539 int i, j, gid, aj, ognr, ntot = 0;
2542 char warn_buf[STRLEN];
2546 fprintf(debug, "Starting numbering %d groups of type %d\n", ng, gtype);
2549 title = gtypes[gtype];
2552 /* Mark all id's as not set */
2553 for (i = 0; (i < natoms); i++)
2558 snew(grps->nm_ind, ng+1); /* +1 for possible rest group */
2559 for (i = 0; (i < ng); i++)
2561 /* Lookup the group name in the block structure */
2562 gid = search_string(ptrs[i], block->nr, gnames);
2563 if ((grptp != egrptpONE) || (i == 0))
2565 grps->nm_ind[grps->nr++] = gid;
2569 fprintf(debug, "Found gid %d for group %s\n", gid, ptrs[i]);
2572 /* Now go over the atoms in the group */
2573 for (j = block->index[gid]; (j < block->index[gid+1]); j++)
2578 /* Range checking */
2579 if ((aj < 0) || (aj >= natoms))
2581 gmx_fatal(FARGS, "Invalid atom number %d in indexfile", aj);
2583 /* Lookup up the old group number */
2587 gmx_fatal(FARGS, "Atom %d in multiple %s groups (%d and %d)",
2588 aj+1, title, ognr+1, i+1);
2592 /* Store the group number in buffer */
2593 if (grptp == egrptpONE)
2606 /* Now check whether we have done all atoms */
2610 if (grptp == egrptpALL)
2612 gmx_fatal(FARGS, "%d atoms are not part of any of the %s groups",
2613 natoms-ntot, title);
2615 else if (grptp == egrptpPART)
2617 sprintf(warn_buf, "%d atoms are not part of any of the %s groups",
2618 natoms-ntot, title);
2619 warning_note(wi, warn_buf);
2621 /* Assign all atoms currently unassigned to a rest group */
2622 for (j = 0; (j < natoms); j++)
2624 if (cbuf[j] == NOGID)
2630 if (grptp != egrptpPART)
2635 "Making dummy/rest group for %s containing %d elements\n",
2636 title, natoms-ntot);
2638 /* Add group name "rest" */
2639 grps->nm_ind[grps->nr] = restnm;
2641 /* Assign the rest name to all atoms not currently assigned to a group */
2642 for (j = 0; (j < natoms); j++)
2644 if (cbuf[j] == NOGID)
2653 if (grps->nr == 1 && (ntot == 0 || ntot == natoms))
2655 /* All atoms are part of one (or no) group, no index required */
2656 groups->ngrpnr[gtype] = 0;
2657 groups->grpnr[gtype] = NULL;
2661 groups->ngrpnr[gtype] = natoms;
2662 snew(groups->grpnr[gtype], natoms);
2663 for (j = 0; (j < natoms); j++)
2665 groups->grpnr[gtype][j] = cbuf[j];
2671 return (bRest && grptp == egrptpPART);
2674 static void calc_nrdf(gmx_mtop_t *mtop, t_inputrec *ir, char **gnames)
2677 gmx_groups_t *groups;
2679 int natoms, ai, aj, i, j, d, g, imin, jmin;
2681 int *nrdf2, *na_vcm, na_tot;
2682 double *nrdf_tc, *nrdf_vcm, nrdf_uc, n_sub = 0;
2683 gmx_mtop_atomloop_all_t aloop;
2685 int mb, mol, ftype, as;
2686 gmx_molblock_t *molb;
2687 gmx_moltype_t *molt;
2690 * First calc 3xnr-atoms for each group
2691 * then subtract half a degree of freedom for each constraint
2693 * Only atoms and nuclei contribute to the degrees of freedom...
2698 groups = &mtop->groups;
2699 natoms = mtop->natoms;
2701 /* Allocate one more for a possible rest group */
2702 /* We need to sum degrees of freedom into doubles,
2703 * since floats give too low nrdf's above 3 million atoms.
2705 snew(nrdf_tc, groups->grps[egcTC].nr+1);
2706 snew(nrdf_vcm, groups->grps[egcVCM].nr+1);
2707 snew(na_vcm, groups->grps[egcVCM].nr+1);
2709 for (i = 0; i < groups->grps[egcTC].nr; i++)
2713 for (i = 0; i < groups->grps[egcVCM].nr+1; i++)
2718 snew(nrdf2, natoms);
2719 aloop = gmx_mtop_atomloop_all_init(mtop);
2720 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
2723 if (atom->ptype == eptAtom || atom->ptype == eptNucleus)
2725 g = ggrpnr(groups, egcFREEZE, i);
2726 /* Double count nrdf for particle i */
2727 for (d = 0; d < DIM; d++)
2729 if (opts->nFreeze[g][d] == 0)
2734 nrdf_tc [ggrpnr(groups, egcTC, i)] += 0.5*nrdf2[i];
2735 nrdf_vcm[ggrpnr(groups, egcVCM, i)] += 0.5*nrdf2[i];
2740 for (mb = 0; mb < mtop->nmolblock; mb++)
2742 molb = &mtop->molblock[mb];
2743 molt = &mtop->moltype[molb->type];
2744 atom = molt->atoms.atom;
2745 for (mol = 0; mol < molb->nmol; mol++)
2747 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
2749 ia = molt->ilist[ftype].iatoms;
2750 for (i = 0; i < molt->ilist[ftype].nr; )
2752 /* Subtract degrees of freedom for the constraints,
2753 * if the particles still have degrees of freedom left.
2754 * If one of the particles is a vsite or a shell, then all
2755 * constraint motion will go there, but since they do not
2756 * contribute to the constraints the degrees of freedom do not
2761 if (((atom[ia[1]].ptype == eptNucleus) ||
2762 (atom[ia[1]].ptype == eptAtom)) &&
2763 ((atom[ia[2]].ptype == eptNucleus) ||
2764 (atom[ia[2]].ptype == eptAtom)))
2782 imin = min(imin, nrdf2[ai]);
2783 jmin = min(jmin, nrdf2[aj]);
2786 nrdf_tc [ggrpnr(groups, egcTC, ai)] -= 0.5*imin;
2787 nrdf_tc [ggrpnr(groups, egcTC, aj)] -= 0.5*jmin;
2788 nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2789 nrdf_vcm[ggrpnr(groups, egcVCM, aj)] -= 0.5*jmin;
2791 ia += interaction_function[ftype].nratoms+1;
2792 i += interaction_function[ftype].nratoms+1;
2795 ia = molt->ilist[F_SETTLE].iatoms;
2796 for (i = 0; i < molt->ilist[F_SETTLE].nr; )
2798 /* Subtract 1 dof from every atom in the SETTLE */
2799 for (j = 0; j < 3; j++)
2802 imin = min(2, nrdf2[ai]);
2804 nrdf_tc [ggrpnr(groups, egcTC, ai)] -= 0.5*imin;
2805 nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2810 as += molt->atoms.nr;
2814 if (ir->ePull == epullCONSTRAINT)
2816 /* Correct nrdf for the COM constraints.
2817 * We correct using the TC and VCM group of the first atom
2818 * in the reference and pull group. If atoms in one pull group
2819 * belong to different TC or VCM groups it is anyhow difficult
2820 * to determine the optimal nrdf assignment.
2824 for (i = 0; i < pull->ncoord; i++)
2828 for (j = 0; j < 2; j++)
2830 const t_pull_group *pgrp;
2832 pgrp = &pull->group[pull->coord[i].group[j]];
2836 /* Subtract 1/2 dof from each group */
2838 nrdf_tc [ggrpnr(groups, egcTC, ai)] -= 0.5*imin;
2839 nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2840 if (nrdf_tc[ggrpnr(groups, egcTC, ai)] < 0)
2842 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)]]);
2847 /* We need to subtract the whole DOF from group j=1 */
2854 if (ir->nstcomm != 0)
2856 /* Subtract 3 from the number of degrees of freedom in each vcm group
2857 * when com translation is removed and 6 when rotation is removed
2860 switch (ir->comm_mode)
2863 n_sub = ndof_com(ir);
2870 gmx_incons("Checking comm_mode");
2873 for (i = 0; i < groups->grps[egcTC].nr; i++)
2875 /* Count the number of atoms of TC group i for every VCM group */
2876 for (j = 0; j < groups->grps[egcVCM].nr+1; j++)
2881 for (ai = 0; ai < natoms; ai++)
2883 if (ggrpnr(groups, egcTC, ai) == i)
2885 na_vcm[ggrpnr(groups, egcVCM, ai)]++;
2889 /* Correct for VCM removal according to the fraction of each VCM
2890 * group present in this TC group.
2892 nrdf_uc = nrdf_tc[i];
2895 fprintf(debug, "T-group[%d] nrdf_uc = %g, n_sub = %g\n",
2899 for (j = 0; j < groups->grps[egcVCM].nr+1; j++)
2901 if (nrdf_vcm[j] > n_sub)
2903 nrdf_tc[i] += nrdf_uc*((double)na_vcm[j]/(double)na_tot)*
2904 (nrdf_vcm[j] - n_sub)/nrdf_vcm[j];
2908 fprintf(debug, " nrdf_vcm[%d] = %g, nrdf = %g\n",
2909 j, nrdf_vcm[j], nrdf_tc[i]);
2914 for (i = 0; (i < groups->grps[egcTC].nr); i++)
2916 opts->nrdf[i] = nrdf_tc[i];
2917 if (opts->nrdf[i] < 0)
2922 "Number of degrees of freedom in T-Coupling group %s is %.2f\n",
2923 gnames[groups->grps[egcTC].nm_ind[i]], opts->nrdf[i]);
2932 static void decode_cos(char *s, t_cosines *cosine)
2935 char format[STRLEN], f1[STRLEN];
2947 sscanf(t, "%d", &(cosine->n));
2954 snew(cosine->a, cosine->n);
2955 snew(cosine->phi, cosine->n);
2957 sprintf(format, "%%*d");
2958 for (i = 0; (i < cosine->n); i++)
2961 strcat(f1, "%lf%lf");
2962 if (sscanf(t, f1, &a, &phi) < 2)
2964 gmx_fatal(FARGS, "Invalid input for electric field shift: '%s'", t);
2967 cosine->phi[i] = phi;
2968 strcat(format, "%*lf%*lf");
2975 static gmx_bool do_egp_flag(t_inputrec *ir, gmx_groups_t *groups,
2976 const char *option, const char *val, int flag)
2978 /* The maximum number of energy group pairs would be MAXPTR*(MAXPTR+1)/2.
2979 * But since this is much larger than STRLEN, such a line can not be parsed.
2980 * The real maximum is the number of names that fit in a string: STRLEN/2.
2982 #define EGP_MAX (STRLEN/2)
2983 int nelem, i, j, k, nr;
2984 char *names[EGP_MAX];
2988 gnames = groups->grpname;
2990 nelem = str_nelem(val, EGP_MAX, names);
2993 gmx_fatal(FARGS, "The number of groups for %s is odd", option);
2995 nr = groups->grps[egcENER].nr;
2997 for (i = 0; i < nelem/2; i++)
3001 gmx_strcasecmp(names[2*i], *(gnames[groups->grps[egcENER].nm_ind[j]])))
3007 gmx_fatal(FARGS, "%s in %s is not an energy group\n",
3008 names[2*i], option);
3012 gmx_strcasecmp(names[2*i+1], *(gnames[groups->grps[egcENER].nm_ind[k]])))
3018 gmx_fatal(FARGS, "%s in %s is not an energy group\n",
3019 names[2*i+1], option);
3021 if ((j < nr) && (k < nr))
3023 ir->opts.egp_flags[nr*j+k] |= flag;
3024 ir->opts.egp_flags[nr*k+j] |= flag;
3033 static void make_swap_groups(
3042 int ig = -1, i = 0, j;
3046 /* Just a quick check here, more thorough checks are in mdrun */
3047 if (strcmp(splitg0name, splitg1name) == 0)
3049 gmx_fatal(FARGS, "The split groups can not both be '%s'.", splitg0name);
3052 /* First get the swap group index atoms */
3053 ig = search_string(swapgname, grps->nr, gnames);
3054 swap->nat = grps->index[ig+1] - grps->index[ig];
3057 fprintf(stderr, "Swap group '%s' contains %d atoms.\n", swapgname, swap->nat);
3058 snew(swap->ind, swap->nat);
3059 for (i = 0; i < swap->nat; i++)
3061 swap->ind[i] = grps->a[grps->index[ig]+i];
3066 gmx_fatal(FARGS, "You defined an empty group of atoms for swapping.");
3069 /* Now do so for the split groups */
3070 for (j = 0; j < 2; j++)
3074 splitg = splitg0name;
3078 splitg = splitg1name;
3081 ig = search_string(splitg, grps->nr, gnames);
3082 swap->nat_split[j] = grps->index[ig+1] - grps->index[ig];
3083 if (swap->nat_split[j] > 0)
3085 fprintf(stderr, "Split group %d '%s' contains %d atom%s.\n",
3086 j, splitg, swap->nat_split[j], (swap->nat_split[j] > 1) ? "s" : "");
3087 snew(swap->ind_split[j], swap->nat_split[j]);
3088 for (i = 0; i < swap->nat_split[j]; i++)
3090 swap->ind_split[j][i] = grps->a[grps->index[ig]+i];
3095 gmx_fatal(FARGS, "Split group %d has to contain at least 1 atom!", j);
3099 /* Now get the solvent group index atoms */
3100 ig = search_string(solgname, grps->nr, gnames);
3101 swap->nat_sol = grps->index[ig+1] - grps->index[ig];
3102 if (swap->nat_sol > 0)
3104 fprintf(stderr, "Solvent group '%s' contains %d atoms.\n", solgname, swap->nat_sol);
3105 snew(swap->ind_sol, swap->nat_sol);
3106 for (i = 0; i < swap->nat_sol; i++)
3108 swap->ind_sol[i] = grps->a[grps->index[ig]+i];
3113 gmx_fatal(FARGS, "You defined an empty group of solvent. Cannot exchange ions.");
3118 void make_IMD_group(t_IMD *IMDgroup, char *IMDgname, t_blocka *grps, char **gnames)
3123 ig = search_string(IMDgname, grps->nr, gnames);
3124 IMDgroup->nat = grps->index[ig+1] - grps->index[ig];
3126 if (IMDgroup->nat > 0)
3128 fprintf(stderr, "Group '%s' with %d atoms can be activated for interactive molecular dynamics (IMD).\n",
3129 IMDgname, IMDgroup->nat);
3130 snew(IMDgroup->ind, IMDgroup->nat);
3131 for (i = 0; i < IMDgroup->nat; i++)
3133 IMDgroup->ind[i] = grps->a[grps->index[ig]+i];
3139 void do_index(const char* mdparin, const char *ndx,
3142 t_inputrec *ir, rvec *v,
3146 gmx_groups_t *groups;
3150 char warnbuf[STRLEN], **gnames;
3151 int nr, ntcg, ntau_t, nref_t, nacc, nofg, nSA, nSA_points, nSA_time, nSA_temp;
3154 int nacg, nfreeze, nfrdim, nenergy, nvcm, nuser;
3155 char *ptr1[MAXPTR], *ptr2[MAXPTR], *ptr3[MAXPTR];
3156 int i, j, k, restnm;
3158 gmx_bool bExcl, bTable, bSetTCpar, bAnneal, bRest;
3159 int nQMmethod, nQMbasis, nQMcharge, nQMmult, nbSH, nCASorb, nCASelec,
3160 nSAon, nSAoff, nSAsteps, nQMg, nbOPT, nbTS;
3161 char warn_buf[STRLEN];
3165 fprintf(stderr, "processing index file...\n");
3171 snew(grps->index, 1);
3173 atoms_all = gmx_mtop_global_atoms(mtop);
3174 analyse(&atoms_all, grps, &gnames, FALSE, TRUE);
3175 free_t_atoms(&atoms_all, FALSE);
3179 grps = init_index(ndx, &gnames);
3182 groups = &mtop->groups;
3183 natoms = mtop->natoms;
3184 symtab = &mtop->symtab;
3186 snew(groups->grpname, grps->nr+1);
3188 for (i = 0; (i < grps->nr); i++)
3190 groups->grpname[i] = put_symtab(symtab, gnames[i]);
3192 groups->grpname[i] = put_symtab(symtab, "rest");
3194 srenew(gnames, grps->nr+1);
3195 gnames[restnm] = *(groups->grpname[i]);
3196 groups->ngrpname = grps->nr+1;
3198 set_warning_line(wi, mdparin, -1);
3200 ntau_t = str_nelem(is->tau_t, MAXPTR, ptr1);
3201 nref_t = str_nelem(is->ref_t, MAXPTR, ptr2);
3202 ntcg = str_nelem(is->tcgrps, MAXPTR, ptr3);
3203 if ((ntau_t != ntcg) || (nref_t != ntcg))
3205 gmx_fatal(FARGS, "Invalid T coupling input: %d groups, %d ref-t values and "
3206 "%d tau-t values", ntcg, nref_t, ntau_t);
3209 bSetTCpar = (ir->etc || EI_SD(ir->eI) || ir->eI == eiBD || EI_TPI(ir->eI));
3210 do_numbering(natoms, groups, ntcg, ptr3, grps, gnames, egcTC,
3211 restnm, bSetTCpar ? egrptpALL : egrptpALL_GENREST, bVerbose, wi);
3212 nr = groups->grps[egcTC].nr;
3214 snew(ir->opts.nrdf, nr);
3215 snew(ir->opts.tau_t, nr);
3216 snew(ir->opts.ref_t, nr);
3217 if (ir->eI == eiBD && ir->bd_fric == 0)
3219 fprintf(stderr, "bd-fric=0, so tau-t will be used as the inverse friction constant(s)\n");
3226 gmx_fatal(FARGS, "Not enough ref-t and tau-t values!");
3230 for (i = 0; (i < nr); i++)
3232 ir->opts.tau_t[i] = strtod(ptr1[i], NULL);
3233 if ((ir->eI == eiBD || ir->eI == eiSD2) && ir->opts.tau_t[i] <= 0)
3235 sprintf(warn_buf, "With integrator %s tau-t should be larger than 0", ei_names[ir->eI]);
3236 warning_error(wi, warn_buf);
3239 if (ir->etc != etcVRESCALE && ir->opts.tau_t[i] == 0)
3241 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.");
3244 if (ir->opts.tau_t[i] >= 0)
3246 tau_min = min(tau_min, ir->opts.tau_t[i]);
3249 if (ir->etc != etcNO && ir->nsttcouple == -1)
3251 ir->nsttcouple = ir_optimal_nsttcouple(ir);
3256 if ((ir->etc == etcNOSEHOOVER) && (ir->epc == epcBERENDSEN))
3258 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");
3260 if ((ir->epc == epcMTTK) && (ir->etc > etcNO))
3262 if (ir->nstpcouple != ir->nsttcouple)
3264 int mincouple = min(ir->nstpcouple, ir->nsttcouple);
3265 ir->nstpcouple = ir->nsttcouple = mincouple;
3266 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);
3267 warning_note(wi, warn_buf);
3271 /* velocity verlet with averaged kinetic energy KE = 0.5*(v(t+1/2) - v(t-1/2)) is implemented
3272 primarily for testing purposes, and does not work with temperature coupling other than 1 */
3274 if (ETC_ANDERSEN(ir->etc))
3276 if (ir->nsttcouple != 1)
3279 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");
3280 warning_note(wi, warn_buf);
3283 nstcmin = tcouple_min_integration_steps(ir->etc);
3286 if (tau_min/(ir->delta_t*ir->nsttcouple) < nstcmin - 10*GMX_REAL_EPS)
3288 sprintf(warn_buf, "For proper integration of the %s thermostat, tau-t (%g) should be at least %d times larger than nsttcouple*dt (%g)",
3289 ETCOUPLTYPE(ir->etc),
3291 ir->nsttcouple*ir->delta_t);
3292 warning(wi, warn_buf);
3295 for (i = 0; (i < nr); i++)
3297 ir->opts.ref_t[i] = strtod(ptr2[i], NULL);
3298 if (ir->opts.ref_t[i] < 0)
3300 gmx_fatal(FARGS, "ref-t for group %d negative", i);
3303 /* set the lambda mc temperature to the md integrator temperature (which should be defined
3304 if we are in this conditional) if mc_temp is negative */
3305 if (ir->expandedvals->mc_temp < 0)
3307 ir->expandedvals->mc_temp = ir->opts.ref_t[0]; /*for now, set to the first reft */
3311 /* Simulated annealing for each group. There are nr groups */
3312 nSA = str_nelem(is->anneal, MAXPTR, ptr1);
3313 if (nSA == 1 && (ptr1[0][0] == 'n' || ptr1[0][0] == 'N'))
3317 if (nSA > 0 && nSA != nr)
3319 gmx_fatal(FARGS, "Not enough annealing values: %d (for %d groups)\n", nSA, nr);
3323 snew(ir->opts.annealing, nr);
3324 snew(ir->opts.anneal_npoints, nr);
3325 snew(ir->opts.anneal_time, nr);
3326 snew(ir->opts.anneal_temp, nr);
3327 for (i = 0; i < nr; i++)
3329 ir->opts.annealing[i] = eannNO;
3330 ir->opts.anneal_npoints[i] = 0;
3331 ir->opts.anneal_time[i] = NULL;
3332 ir->opts.anneal_temp[i] = NULL;
3337 for (i = 0; i < nr; i++)
3339 if (ptr1[i][0] == 'n' || ptr1[i][0] == 'N')
3341 ir->opts.annealing[i] = eannNO;
3343 else if (ptr1[i][0] == 's' || ptr1[i][0] == 'S')
3345 ir->opts.annealing[i] = eannSINGLE;
3348 else if (ptr1[i][0] == 'p' || ptr1[i][0] == 'P')
3350 ir->opts.annealing[i] = eannPERIODIC;
3356 /* Read the other fields too */
3357 nSA_points = str_nelem(is->anneal_npoints, MAXPTR, ptr1);
3358 if (nSA_points != nSA)
3360 gmx_fatal(FARGS, "Found %d annealing-npoints values for %d groups\n", nSA_points, nSA);
3362 for (k = 0, i = 0; i < nr; i++)
3364 ir->opts.anneal_npoints[i] = strtol(ptr1[i], NULL, 10);
3365 if (ir->opts.anneal_npoints[i] == 1)
3367 gmx_fatal(FARGS, "Please specify at least a start and an end point for annealing\n");
3369 snew(ir->opts.anneal_time[i], ir->opts.anneal_npoints[i]);
3370 snew(ir->opts.anneal_temp[i], ir->opts.anneal_npoints[i]);
3371 k += ir->opts.anneal_npoints[i];
3374 nSA_time = str_nelem(is->anneal_time, MAXPTR, ptr1);
3377 gmx_fatal(FARGS, "Found %d annealing-time values, wanter %d\n", nSA_time, k);
3379 nSA_temp = str_nelem(is->anneal_temp, MAXPTR, ptr2);
3382 gmx_fatal(FARGS, "Found %d annealing-temp values, wanted %d\n", nSA_temp, k);
3385 for (i = 0, k = 0; i < nr; i++)
3388 for (j = 0; j < ir->opts.anneal_npoints[i]; j++)
3390 ir->opts.anneal_time[i][j] = strtod(ptr1[k], NULL);
3391 ir->opts.anneal_temp[i][j] = strtod(ptr2[k], NULL);
3394 if (ir->opts.anneal_time[i][0] > (ir->init_t+GMX_REAL_EPS))
3396 gmx_fatal(FARGS, "First time point for annealing > init_t.\n");
3402 if (ir->opts.anneal_time[i][j] < ir->opts.anneal_time[i][j-1])
3404 gmx_fatal(FARGS, "Annealing timepoints out of order: t=%f comes after t=%f\n",
3405 ir->opts.anneal_time[i][j], ir->opts.anneal_time[i][j-1]);
3408 if (ir->opts.anneal_temp[i][j] < 0)
3410 gmx_fatal(FARGS, "Found negative temperature in annealing: %f\n", ir->opts.anneal_temp[i][j]);
3415 /* Print out some summary information, to make sure we got it right */
3416 for (i = 0, k = 0; i < nr; i++)
3418 if (ir->opts.annealing[i] != eannNO)
3420 j = groups->grps[egcTC].nm_ind[i];
3421 fprintf(stderr, "Simulated annealing for group %s: %s, %d timepoints\n",
3422 *(groups->grpname[j]), eann_names[ir->opts.annealing[i]],
3423 ir->opts.anneal_npoints[i]);
3424 fprintf(stderr, "Time (ps) Temperature (K)\n");
3425 /* All terms except the last one */
3426 for (j = 0; j < (ir->opts.anneal_npoints[i]-1); j++)
3428 fprintf(stderr, "%9.1f %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3431 /* Finally the last one */
3432 j = ir->opts.anneal_npoints[i]-1;
3433 if (ir->opts.annealing[i] == eannSINGLE)
3435 fprintf(stderr, "%9.1f- %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3439 fprintf(stderr, "%9.1f %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3440 if (fabs(ir->opts.anneal_temp[i][j]-ir->opts.anneal_temp[i][0]) > GMX_REAL_EPS)
3442 warning_note(wi, "There is a temperature jump when your annealing loops back.\n");
3451 if (ir->ePull != epullNO)
3453 make_pull_groups(ir->pull, is->pull_grp, grps, gnames);
3455 make_pull_coords(ir->pull);
3460 make_rotation_groups(ir->rot, is->rot_grp, grps, gnames);
3463 if (ir->eSwapCoords != eswapNO)
3465 make_swap_groups(ir->swap, swapgrp, splitgrp0, splitgrp1, solgrp, grps, gnames);
3468 /* Make indices for IMD session */
3471 make_IMD_group(ir->imd, is->imd_grp, grps, gnames);
3474 nacc = str_nelem(is->acc, MAXPTR, ptr1);
3475 nacg = str_nelem(is->accgrps, MAXPTR, ptr2);
3476 if (nacg*DIM != nacc)
3478 gmx_fatal(FARGS, "Invalid Acceleration input: %d groups and %d acc. values",
3481 do_numbering(natoms, groups, nacg, ptr2, grps, gnames, egcACC,
3482 restnm, egrptpALL_GENREST, bVerbose, wi);
3483 nr = groups->grps[egcACC].nr;
3484 snew(ir->opts.acc, nr);
3485 ir->opts.ngacc = nr;
3487 for (i = k = 0; (i < nacg); i++)
3489 for (j = 0; (j < DIM); j++, k++)
3491 ir->opts.acc[i][j] = strtod(ptr1[k], NULL);
3494 for (; (i < nr); i++)
3496 for (j = 0; (j < DIM); j++)
3498 ir->opts.acc[i][j] = 0;
3502 nfrdim = str_nelem(is->frdim, MAXPTR, ptr1);
3503 nfreeze = str_nelem(is->freeze, MAXPTR, ptr2);
3504 if (nfrdim != DIM*nfreeze)
3506 gmx_fatal(FARGS, "Invalid Freezing input: %d groups and %d freeze values",
3509 do_numbering(natoms, groups, nfreeze, ptr2, grps, gnames, egcFREEZE,
3510 restnm, egrptpALL_GENREST, bVerbose, wi);
3511 nr = groups->grps[egcFREEZE].nr;
3512 ir->opts.ngfrz = nr;
3513 snew(ir->opts.nFreeze, nr);
3514 for (i = k = 0; (i < nfreeze); i++)
3516 for (j = 0; (j < DIM); j++, k++)
3518 ir->opts.nFreeze[i][j] = (gmx_strncasecmp(ptr1[k], "Y", 1) == 0);
3519 if (!ir->opts.nFreeze[i][j])
3521 if (gmx_strncasecmp(ptr1[k], "N", 1) != 0)
3523 sprintf(warnbuf, "Please use Y(ES) or N(O) for freezedim only "
3524 "(not %s)", ptr1[k]);
3525 warning(wi, warn_buf);
3530 for (; (i < nr); i++)
3532 for (j = 0; (j < DIM); j++)
3534 ir->opts.nFreeze[i][j] = 0;
3538 nenergy = str_nelem(is->energy, MAXPTR, ptr1);
3539 do_numbering(natoms, groups, nenergy, ptr1, grps, gnames, egcENER,
3540 restnm, egrptpALL_GENREST, bVerbose, wi);
3541 add_wall_energrps(groups, ir->nwall, symtab);
3542 ir->opts.ngener = groups->grps[egcENER].nr;
3543 nvcm = str_nelem(is->vcm, MAXPTR, ptr1);
3545 do_numbering(natoms, groups, nvcm, ptr1, grps, gnames, egcVCM,
3546 restnm, nvcm == 0 ? egrptpALL_GENREST : egrptpPART, bVerbose, wi);
3549 warning(wi, "Some atoms are not part of any center of mass motion removal group.\n"
3550 "This may lead to artifacts.\n"
3551 "In most cases one should use one group for the whole system.");
3554 /* Now we have filled the freeze struct, so we can calculate NRDF */
3555 calc_nrdf(mtop, ir, gnames);
3561 /* Must check per group! */
3562 for (i = 0; (i < ir->opts.ngtc); i++)
3564 ntot += ir->opts.nrdf[i];
3566 if (ntot != (DIM*natoms))
3568 fac = sqrt(ntot/(DIM*natoms));
3571 fprintf(stderr, "Scaling velocities by a factor of %.3f to account for constraints\n"
3572 "and removal of center of mass motion\n", fac);
3574 for (i = 0; (i < natoms); i++)
3576 svmul(fac, v[i], v[i]);
3581 nuser = str_nelem(is->user1, MAXPTR, ptr1);
3582 do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcUser1,
3583 restnm, egrptpALL_GENREST, bVerbose, wi);
3584 nuser = str_nelem(is->user2, MAXPTR, ptr1);
3585 do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcUser2,
3586 restnm, egrptpALL_GENREST, bVerbose, wi);
3587 nuser = str_nelem(is->x_compressed_groups, MAXPTR, ptr1);
3588 do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcCompressedX,
3589 restnm, egrptpONE, bVerbose, wi);
3590 nofg = str_nelem(is->orirefitgrp, MAXPTR, ptr1);
3591 do_numbering(natoms, groups, nofg, ptr1, grps, gnames, egcORFIT,
3592 restnm, egrptpALL_GENREST, bVerbose, wi);
3594 /* QMMM input processing */
3595 nQMg = str_nelem(is->QMMM, MAXPTR, ptr1);
3596 nQMmethod = str_nelem(is->QMmethod, MAXPTR, ptr2);
3597 nQMbasis = str_nelem(is->QMbasis, MAXPTR, ptr3);
3598 if ((nQMmethod != nQMg) || (nQMbasis != nQMg))
3600 gmx_fatal(FARGS, "Invalid QMMM input: %d groups %d basissets"
3601 " and %d methods\n", nQMg, nQMbasis, nQMmethod);
3603 /* group rest, if any, is always MM! */
3604 do_numbering(natoms, groups, nQMg, ptr1, grps, gnames, egcQMMM,
3605 restnm, egrptpALL_GENREST, bVerbose, wi);
3606 nr = nQMg; /*atoms->grps[egcQMMM].nr;*/
3607 ir->opts.ngQM = nQMg;
3608 snew(ir->opts.QMmethod, nr);
3609 snew(ir->opts.QMbasis, nr);
3610 for (i = 0; i < nr; i++)
3612 /* input consists of strings: RHF CASSCF PM3 .. These need to be
3613 * converted to the corresponding enum in names.c
3615 ir->opts.QMmethod[i] = search_QMstring(ptr2[i], eQMmethodNR,
3617 ir->opts.QMbasis[i] = search_QMstring(ptr3[i], eQMbasisNR,
3621 nQMmult = str_nelem(is->QMmult, MAXPTR, ptr1);
3622 nQMcharge = str_nelem(is->QMcharge, MAXPTR, ptr2);
3623 nbSH = str_nelem(is->bSH, MAXPTR, ptr3);
3624 snew(ir->opts.QMmult, nr);
3625 snew(ir->opts.QMcharge, nr);
3626 snew(ir->opts.bSH, nr);
3628 for (i = 0; i < nr; i++)
3630 ir->opts.QMmult[i] = strtol(ptr1[i], NULL, 10);
3631 ir->opts.QMcharge[i] = strtol(ptr2[i], NULL, 10);
3632 ir->opts.bSH[i] = (gmx_strncasecmp(ptr3[i], "Y", 1) == 0);
3635 nCASelec = str_nelem(is->CASelectrons, MAXPTR, ptr1);
3636 nCASorb = str_nelem(is->CASorbitals, MAXPTR, ptr2);
3637 snew(ir->opts.CASelectrons, nr);
3638 snew(ir->opts.CASorbitals, nr);
3639 for (i = 0; i < nr; i++)
3641 ir->opts.CASelectrons[i] = strtol(ptr1[i], NULL, 10);
3642 ir->opts.CASorbitals[i] = strtol(ptr2[i], NULL, 10);
3644 /* special optimization options */
3646 nbOPT = str_nelem(is->bOPT, MAXPTR, ptr1);
3647 nbTS = str_nelem(is->bTS, MAXPTR, ptr2);
3648 snew(ir->opts.bOPT, nr);
3649 snew(ir->opts.bTS, nr);
3650 for (i = 0; i < nr; i++)
3652 ir->opts.bOPT[i] = (gmx_strncasecmp(ptr1[i], "Y", 1) == 0);
3653 ir->opts.bTS[i] = (gmx_strncasecmp(ptr2[i], "Y", 1) == 0);
3655 nSAon = str_nelem(is->SAon, MAXPTR, ptr1);
3656 nSAoff = str_nelem(is->SAoff, MAXPTR, ptr2);
3657 nSAsteps = str_nelem(is->SAsteps, MAXPTR, ptr3);
3658 snew(ir->opts.SAon, nr);
3659 snew(ir->opts.SAoff, nr);
3660 snew(ir->opts.SAsteps, nr);
3662 for (i = 0; i < nr; i++)
3664 ir->opts.SAon[i] = strtod(ptr1[i], NULL);
3665 ir->opts.SAoff[i] = strtod(ptr2[i], NULL);
3666 ir->opts.SAsteps[i] = strtol(ptr3[i], NULL, 10);
3668 /* end of QMMM input */
3672 for (i = 0; (i < egcNR); i++)
3674 fprintf(stderr, "%-16s has %d element(s):", gtypes[i], groups->grps[i].nr);
3675 for (j = 0; (j < groups->grps[i].nr); j++)
3677 fprintf(stderr, " %s", *(groups->grpname[groups->grps[i].nm_ind[j]]));
3679 fprintf(stderr, "\n");
3683 nr = groups->grps[egcENER].nr;
3684 snew(ir->opts.egp_flags, nr*nr);
3686 bExcl = do_egp_flag(ir, groups, "energygrp-excl", is->egpexcl, EGP_EXCL);
3687 if (bExcl && ir->cutoff_scheme == ecutsVERLET)
3689 warning_error(wi, "Energy group exclusions are not (yet) implemented for the Verlet scheme");
3691 if (bExcl && EEL_FULL(ir->coulombtype))
3693 warning(wi, "Can not exclude the lattice Coulomb energy between energy groups");
3696 bTable = do_egp_flag(ir, groups, "energygrp-table", is->egptable, EGP_TABLE);
3697 if (bTable && !(ir->vdwtype == evdwUSER) &&
3698 !(ir->coulombtype == eelUSER) && !(ir->coulombtype == eelPMEUSER) &&
3699 !(ir->coulombtype == eelPMEUSERSWITCH))
3701 gmx_fatal(FARGS, "Can only have energy group pair tables in combination with user tables for VdW and/or Coulomb");
3704 decode_cos(is->efield_x, &(ir->ex[XX]));
3705 decode_cos(is->efield_xt, &(ir->et[XX]));
3706 decode_cos(is->efield_y, &(ir->ex[YY]));
3707 decode_cos(is->efield_yt, &(ir->et[YY]));
3708 decode_cos(is->efield_z, &(ir->ex[ZZ]));
3709 decode_cos(is->efield_zt, &(ir->et[ZZ]));
3713 do_adress_index(ir->adress, groups, gnames, &(ir->opts), wi);
3716 for (i = 0; (i < grps->nr); i++)
3728 static void check_disre(gmx_mtop_t *mtop)
3730 gmx_ffparams_t *ffparams;
3731 t_functype *functype;
3733 int i, ndouble, ftype;
3734 int label, old_label;
3736 if (gmx_mtop_ftype_count(mtop, F_DISRES) > 0)
3738 ffparams = &mtop->ffparams;
3739 functype = ffparams->functype;
3740 ip = ffparams->iparams;
3743 for (i = 0; i < ffparams->ntypes; i++)
3745 ftype = functype[i];
3746 if (ftype == F_DISRES)
3748 label = ip[i].disres.label;
3749 if (label == old_label)
3751 fprintf(stderr, "Distance restraint index %d occurs twice\n", label);
3759 gmx_fatal(FARGS, "Found %d double distance restraint indices,\n"
3760 "probably the parameters for multiple pairs in one restraint "
3761 "are not identical\n", ndouble);
3766 static gmx_bool absolute_reference(t_inputrec *ir, gmx_mtop_t *sys,
3767 gmx_bool posres_only,
3771 gmx_mtop_ilistloop_t iloop;
3781 for (d = 0; d < DIM; d++)
3783 AbsRef[d] = (d < ndof_com(ir) ? 0 : 1);
3785 /* Check for freeze groups */
3786 for (g = 0; g < ir->opts.ngfrz; g++)
3788 for (d = 0; d < DIM; d++)
3790 if (ir->opts.nFreeze[g][d] != 0)
3798 /* Check for position restraints */
3799 iloop = gmx_mtop_ilistloop_init(sys);
3800 while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
3803 (AbsRef[XX] == 0 || AbsRef[YY] == 0 || AbsRef[ZZ] == 0))
3805 for (i = 0; i < ilist[F_POSRES].nr; i += 2)
3807 pr = &sys->ffparams.iparams[ilist[F_POSRES].iatoms[i]];
3808 for (d = 0; d < DIM; d++)
3810 if (pr->posres.fcA[d] != 0)
3816 for (i = 0; i < ilist[F_FBPOSRES].nr; i += 2)
3818 /* Check for flat-bottom posres */
3819 pr = &sys->ffparams.iparams[ilist[F_FBPOSRES].iatoms[i]];
3820 if (pr->fbposres.k != 0)
3822 switch (pr->fbposres.geom)
3824 case efbposresSPHERE:
3825 AbsRef[XX] = AbsRef[YY] = AbsRef[ZZ] = 1;
3827 case efbposresCYLINDER:
3828 AbsRef[XX] = AbsRef[YY] = 1;
3830 case efbposresX: /* d=XX */
3831 case efbposresY: /* d=YY */
3832 case efbposresZ: /* d=ZZ */
3833 d = pr->fbposres.geom - efbposresX;
3837 gmx_fatal(FARGS, " Invalid geometry for flat-bottom position restraint.\n"
3838 "Expected nr between 1 and %d. Found %d\n", efbposresNR-1,
3846 return (AbsRef[XX] != 0 && AbsRef[YY] != 0 && AbsRef[ZZ] != 0);
3850 check_combination_rule_differences(const gmx_mtop_t *mtop, int state,
3851 gmx_bool *bC6ParametersWorkWithGeometricRules,
3852 gmx_bool *bC6ParametersWorkWithLBRules,
3853 gmx_bool *bLBRulesPossible)
3855 int ntypes, tpi, tpj, thisLBdiff, thisgeomdiff;
3858 double geometricdiff, LBdiff;
3859 double c6i, c6j, c12i, c12j;
3860 double c6, c6_geometric, c6_LB;
3861 double sigmai, sigmaj, epsi, epsj;
3862 gmx_bool bCanDoLBRules, bCanDoGeometricRules;
3865 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
3866 * force-field floating point parameters.
3869 ptr = getenv("GMX_LJCOMB_TOL");
3874 sscanf(ptr, "%lf", &dbl);
3878 *bC6ParametersWorkWithLBRules = TRUE;
3879 *bC6ParametersWorkWithGeometricRules = TRUE;
3880 bCanDoLBRules = TRUE;
3881 bCanDoGeometricRules = TRUE;
3882 ntypes = mtop->ffparams.atnr;
3883 snew(typecount, ntypes);
3884 gmx_mtop_count_atomtypes(mtop, state, typecount);
3885 geometricdiff = LBdiff = 0.0;
3886 *bLBRulesPossible = TRUE;
3887 for (tpi = 0; tpi < ntypes; ++tpi)
3889 c6i = mtop->ffparams.iparams[(ntypes + 1) * tpi].lj.c6;
3890 c12i = mtop->ffparams.iparams[(ntypes + 1) * tpi].lj.c12;
3891 for (tpj = tpi; tpj < ntypes; ++tpj)
3893 c6j = mtop->ffparams.iparams[(ntypes + 1) * tpj].lj.c6;
3894 c12j = mtop->ffparams.iparams[(ntypes + 1) * tpj].lj.c12;
3895 c6 = mtop->ffparams.iparams[ntypes * tpi + tpj].lj.c6;
3896 c6_geometric = sqrt(c6i * c6j);
3897 if (!gmx_numzero(c6_geometric))
3899 if (!gmx_numzero(c12i) && !gmx_numzero(c12j))
3901 sigmai = pow(c12i / c6i, 1.0/6.0);
3902 sigmaj = pow(c12j / c6j, 1.0/6.0);
3903 epsi = c6i * c6i /(4.0 * c12i);
3904 epsj = c6j * c6j /(4.0 * c12j);
3905 c6_LB = 4.0 * pow(epsi * epsj, 1.0/2.0) * pow(0.5 * (sigmai + sigmaj), 6);
3909 *bLBRulesPossible = FALSE;
3910 c6_LB = c6_geometric;
3912 bCanDoLBRules = gmx_within_tol(c6_LB, c6, tol);
3915 if (FALSE == bCanDoLBRules)
3917 *bC6ParametersWorkWithLBRules = FALSE;
3920 bCanDoGeometricRules = gmx_within_tol(c6_geometric, c6, tol);
3922 if (FALSE == bCanDoGeometricRules)
3924 *bC6ParametersWorkWithGeometricRules = FALSE;
3932 check_combination_rules(const t_inputrec *ir, const gmx_mtop_t *mtop,
3936 gmx_bool bLBRulesPossible, bC6ParametersWorkWithGeometricRules, bC6ParametersWorkWithLBRules;
3938 check_combination_rule_differences(mtop, 0,
3939 &bC6ParametersWorkWithGeometricRules,
3940 &bC6ParametersWorkWithLBRules,
3942 if (ir->ljpme_combination_rule == eljpmeLB)
3944 if (FALSE == bC6ParametersWorkWithLBRules || FALSE == bLBRulesPossible)
3946 warning(wi, "You are using arithmetic-geometric combination rules "
3947 "in LJ-PME, but your non-bonded C6 parameters do not "
3948 "follow these rules.");
3953 if (FALSE == bC6ParametersWorkWithGeometricRules)
3955 if (ir->eDispCorr != edispcNO)
3957 warning_note(wi, "You are using geometric combination rules in "
3958 "LJ-PME, but your non-bonded C6 parameters do "
3959 "not follow these rules. "
3960 "This will introduce very small errors in the forces and energies in "
3961 "your simulations. Dispersion correction will correct total energy "
3962 "and/or pressure for isotropic systems, but not forces or surface tensions.");
3966 warning_note(wi, "You are using geometric combination rules in "
3967 "LJ-PME, but your non-bonded C6 parameters do "
3968 "not follow these rules. "
3969 "This will introduce very small errors in the forces and energies in "
3970 "your simulations. If your system is homogeneous, consider using dispersion correction "
3971 "for the total energy and pressure.");
3977 void triple_check(const char *mdparin, t_inputrec *ir, gmx_mtop_t *sys,
3980 char err_buf[STRLEN];
3981 int i, m, c, nmol, npct;
3982 gmx_bool bCharge, bAcc;
3983 real gdt_max, *mgrp, mt;
3985 gmx_mtop_atomloop_block_t aloopb;
3986 gmx_mtop_atomloop_all_t aloop;
3989 char warn_buf[STRLEN];
3991 set_warning_line(wi, mdparin, -1);
3993 if (ir->cutoff_scheme == ecutsVERLET &&
3994 ir->verletbuf_tol > 0 &&
3996 ((EI_MD(ir->eI) || EI_SD(ir->eI)) &&
3997 (ir->etc == etcVRESCALE || ir->etc == etcBERENDSEN)))
3999 /* Check if a too small Verlet buffer might potentially
4000 * cause more drift than the thermostat can couple off.
4002 /* Temperature error fraction for warning and suggestion */
4003 const real T_error_warn = 0.002;
4004 const real T_error_suggest = 0.001;
4005 /* For safety: 2 DOF per atom (typical with constraints) */
4006 const real nrdf_at = 2;
4007 real T, tau, max_T_error;
4012 for (i = 0; i < ir->opts.ngtc; i++)
4014 T = max(T, ir->opts.ref_t[i]);
4015 tau = max(tau, ir->opts.tau_t[i]);
4019 /* This is a worst case estimate of the temperature error,
4020 * assuming perfect buffer estimation and no cancelation
4021 * of errors. The factor 0.5 is because energy distributes
4022 * equally over Ekin and Epot.
4024 max_T_error = 0.5*tau*ir->verletbuf_tol/(nrdf_at*BOLTZ*T);
4025 if (max_T_error > T_error_warn)
4027 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.",
4028 ir->verletbuf_tol, T, tau,
4030 100*T_error_suggest,
4031 ir->verletbuf_tol*T_error_suggest/max_T_error);
4032 warning(wi, warn_buf);
4037 if (ETC_ANDERSEN(ir->etc))
4041 for (i = 0; i < ir->opts.ngtc; i++)
4043 sprintf(err_buf, "all tau_t must currently be equal using Andersen temperature control, violated for group %d", i);
4044 CHECK(ir->opts.tau_t[0] != ir->opts.tau_t[i]);
4045 sprintf(err_buf, "all tau_t must be postive using Andersen temperature control, tau_t[%d]=%10.6f",
4046 i, ir->opts.tau_t[i]);
4047 CHECK(ir->opts.tau_t[i] < 0);
4050 for (i = 0; i < ir->opts.ngtc; i++)
4052 int nsteps = (int)(ir->opts.tau_t[i]/ir->delta_t);
4053 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);
4054 CHECK((nsteps % ir->nstcomm) && (ir->etc == etcANDERSENMASSIVE));
4058 if (EI_DYNAMICS(ir->eI) && !EI_SD(ir->eI) && ir->eI != eiBD &&
4059 ir->comm_mode == ecmNO &&
4060 !(absolute_reference(ir, sys, FALSE, AbsRef) || ir->nsteps <= 10) &&
4061 !ETC_ANDERSEN(ir->etc))
4063 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");
4066 /* Check for pressure coupling with absolute position restraints */
4067 if (ir->epc != epcNO && ir->refcoord_scaling == erscNO)
4069 absolute_reference(ir, sys, TRUE, AbsRef);
4071 for (m = 0; m < DIM; m++)
4073 if (AbsRef[m] && norm2(ir->compress[m]) > 0)
4075 warning(wi, "You are using pressure coupling with absolute position restraints, this will give artifacts. Use the refcoord_scaling option.");
4083 aloopb = gmx_mtop_atomloop_block_init(sys);
4084 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
4086 if (atom->q != 0 || atom->qB != 0)
4094 if (EEL_FULL(ir->coulombtype))
4097 "You are using full electrostatics treatment %s for a system without charges.\n"
4098 "This costs a lot of performance for just processing zeros, consider using %s instead.\n",
4099 EELTYPE(ir->coulombtype), EELTYPE(eelCUT));
4100 warning(wi, err_buf);
4105 if (ir->coulombtype == eelCUT && ir->rcoulomb > 0 && !ir->implicit_solvent)
4108 "You are using a plain Coulomb cut-off, which might produce artifacts.\n"
4109 "You might want to consider using %s electrostatics.\n",
4111 warning_note(wi, err_buf);
4115 /* Check if combination rules used in LJ-PME are the same as in the force field */
4116 if (EVDW_PME(ir->vdwtype))
4118 check_combination_rules(ir, sys, wi);
4121 /* Generalized reaction field */
4122 if (ir->opts.ngtc == 0)
4124 sprintf(err_buf, "No temperature coupling while using coulombtype %s",
4126 CHECK(ir->coulombtype == eelGRF);
4130 sprintf(err_buf, "When using coulombtype = %s"
4131 " ref-t for temperature coupling should be > 0",
4133 CHECK((ir->coulombtype == eelGRF) && (ir->opts.ref_t[0] <= 0));
4136 if (ir->eI == eiSD2)
4138 sprintf(warn_buf, "The stochastic dynamics integrator %s is deprecated, since\n"
4139 "it is slower than integrator %s and is slightly less accurate\n"
4140 "with constraints. Use the %s integrator.",
4141 ei_names[ir->eI], ei_names[eiSD1], ei_names[eiSD1]);
4142 warning_note(wi, warn_buf);
4146 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4148 for (m = 0; (m < DIM); m++)
4150 if (fabs(ir->opts.acc[i][m]) > 1e-6)
4159 snew(mgrp, sys->groups.grps[egcACC].nr);
4160 aloop = gmx_mtop_atomloop_all_init(sys);
4161 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
4163 mgrp[ggrpnr(&sys->groups, egcACC, i)] += atom->m;
4166 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4168 for (m = 0; (m < DIM); m++)
4170 acc[m] += ir->opts.acc[i][m]*mgrp[i];
4174 for (m = 0; (m < DIM); m++)
4176 if (fabs(acc[m]) > 1e-6)
4178 const char *dim[DIM] = { "X", "Y", "Z" };
4180 "Net Acceleration in %s direction, will %s be corrected\n",
4181 dim[m], ir->nstcomm != 0 ? "" : "not");
4182 if (ir->nstcomm != 0 && m < ndof_com(ir))
4185 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4187 ir->opts.acc[i][m] -= acc[m];
4195 if (ir->efep != efepNO && ir->fepvals->sc_alpha != 0 &&
4196 !gmx_within_tol(sys->ffparams.reppow, 12.0, 10*GMX_DOUBLE_EPS))
4198 gmx_fatal(FARGS, "Soft-core interactions are only supported with VdW repulsion power 12");
4201 if (ir->ePull != epullNO)
4203 gmx_bool bPullAbsoluteRef;
4205 bPullAbsoluteRef = FALSE;
4206 for (i = 0; i < ir->pull->ncoord; i++)
4208 bPullAbsoluteRef = bPullAbsoluteRef ||
4209 ir->pull->coord[i].group[0] == 0 ||
4210 ir->pull->coord[i].group[1] == 0;
4212 if (bPullAbsoluteRef)
4214 absolute_reference(ir, sys, FALSE, AbsRef);
4215 for (m = 0; m < DIM; m++)
4217 if (ir->pull->dim[m] && !AbsRef[m])
4219 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.");
4225 if (ir->pull->eGeom == epullgDIRPBC)
4227 for (i = 0; i < 3; i++)
4229 for (m = 0; m <= i; m++)
4231 if ((ir->epc != epcNO && ir->compress[i][m] != 0) ||
4232 ir->deform[i][m] != 0)
4234 for (c = 0; c < ir->pull->ncoord; c++)
4236 if (ir->pull->coord[c].vec[m] != 0)
4238 gmx_fatal(FARGS, "Can not have dynamic box while using pull geometry '%s' (dim %c)", EPULLGEOM(ir->pull->eGeom), 'x'+m);
4250 void double_check(t_inputrec *ir, matrix box, gmx_bool bConstr, warninp_t wi)
4254 char warn_buf[STRLEN];
4257 ptr = check_box(ir->ePBC, box);
4260 warning_error(wi, ptr);
4263 if (bConstr && ir->eConstrAlg == econtSHAKE)
4265 if (ir->shake_tol <= 0.0)
4267 sprintf(warn_buf, "ERROR: shake-tol must be > 0 instead of %g\n",
4269 warning_error(wi, warn_buf);
4272 if (IR_TWINRANGE(*ir) && ir->nstlist > 1)
4274 sprintf(warn_buf, "With twin-range cut-off's and SHAKE the virial and the pressure are incorrect.");
4275 if (ir->epc == epcNO)
4277 warning(wi, warn_buf);
4281 warning_error(wi, warn_buf);
4286 if ( (ir->eConstrAlg == econtLINCS) && bConstr)
4288 /* If we have Lincs constraints: */
4289 if (ir->eI == eiMD && ir->etc == etcNO &&
4290 ir->eConstrAlg == econtLINCS && ir->nLincsIter == 1)
4292 sprintf(warn_buf, "For energy conservation with LINCS, lincs_iter should be 2 or larger.\n");
4293 warning_note(wi, warn_buf);
4296 if ((ir->eI == eiCG || ir->eI == eiLBFGS) && (ir->nProjOrder < 8))
4298 sprintf(warn_buf, "For accurate %s with LINCS constraints, lincs-order should be 8 or more.", ei_names[ir->eI]);
4299 warning_note(wi, warn_buf);
4301 if (ir->epc == epcMTTK)
4303 warning_error(wi, "MTTK not compatible with lincs -- use shake instead.");
4307 if (bConstr && ir->epc == epcMTTK)
4309 warning_note(wi, "MTTK with constraints is deprecated, and will be removed in GROMACS 5.1");
4312 if (ir->LincsWarnAngle > 90.0)
4314 sprintf(warn_buf, "lincs-warnangle can not be larger than 90 degrees, setting it to 90.\n");
4315 warning(wi, warn_buf);
4316 ir->LincsWarnAngle = 90.0;
4319 if (ir->ePBC != epbcNONE)
4321 if (ir->nstlist == 0)
4323 warning(wi, "With nstlist=0 atoms are only put into the box at step 0, therefore drifting atoms might cause the simulation to crash.");
4325 bTWIN = (ir->rlistlong > ir->rlist);
4326 if (ir->ns_type == ensGRID)
4328 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
4330 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",
4331 bTWIN ? (ir->rcoulomb == ir->rlistlong ? "rcoulomb" : "rvdw") : "rlist");
4332 warning_error(wi, warn_buf);
4337 min_size = min(box[XX][XX], min(box[YY][YY], box[ZZ][ZZ]));
4338 if (2*ir->rlistlong >= min_size)
4340 sprintf(warn_buf, "ERROR: One of the box lengths is smaller than twice the cut-off length. Increase the box size or decrease rlist.");
4341 warning_error(wi, warn_buf);
4344 fprintf(stderr, "Grid search might allow larger cut-off's than simple search with triclinic boxes.");
4351 void check_chargegroup_radii(const gmx_mtop_t *mtop, const t_inputrec *ir,
4355 real rvdw1, rvdw2, rcoul1, rcoul2;
4356 char warn_buf[STRLEN];
4358 calc_chargegroup_radii(mtop, x, &rvdw1, &rvdw2, &rcoul1, &rcoul2);
4362 printf("Largest charge group radii for Van der Waals: %5.3f, %5.3f nm\n",
4367 printf("Largest charge group radii for Coulomb: %5.3f, %5.3f nm\n",
4373 if (rvdw1 + rvdw2 > ir->rlist ||
4374 rcoul1 + rcoul2 > ir->rlist)
4377 "The sum of the two largest charge group radii (%f) "
4378 "is larger than rlist (%f)\n",
4379 max(rvdw1+rvdw2, rcoul1+rcoul2), ir->rlist);
4380 warning(wi, warn_buf);
4384 /* Here we do not use the zero at cut-off macro,
4385 * since user defined interactions might purposely
4386 * not be zero at the cut-off.
4388 if (ir_vdw_is_zero_at_cutoff(ir) &&
4389 rvdw1 + rvdw2 > ir->rlistlong - ir->rvdw)
4391 sprintf(warn_buf, "The sum of the two largest charge group "
4392 "radii (%f) is larger than %s (%f) - rvdw (%f).\n"
4393 "With exact cut-offs, better performance can be "
4394 "obtained with cutoff-scheme = %s, because it "
4395 "does not use charge groups at all.",
4397 ir->rlistlong > ir->rlist ? "rlistlong" : "rlist",
4398 ir->rlistlong, ir->rvdw,
4399 ecutscheme_names[ecutsVERLET]);
4402 warning(wi, warn_buf);
4406 warning_note(wi, warn_buf);
4409 if (ir_coulomb_is_zero_at_cutoff(ir) &&
4410 rcoul1 + rcoul2 > ir->rlistlong - ir->rcoulomb)
4412 sprintf(warn_buf, "The sum of the two largest charge group radii (%f) is larger than %s (%f) - rcoulomb (%f).\n"
4413 "With exact cut-offs, better performance can be obtained with cutoff-scheme = %s, because it does not use charge groups at all.",
4415 ir->rlistlong > ir->rlist ? "rlistlong" : "rlist",
4416 ir->rlistlong, ir->rcoulomb,
4417 ecutscheme_names[ecutsVERLET]);
4420 warning(wi, warn_buf);
4424 warning_note(wi, warn_buf);