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->nstlist <= 0)
375 warning_error(wi, "With Verlet lists nstlist should be larger than 0");
378 if (ir->nstlist < 10)
380 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.");
383 rc_max = max(ir->rvdw, ir->rcoulomb);
385 if (ir->verletbuf_tol <= 0)
387 if (ir->verletbuf_tol == 0)
389 warning_error(wi, "Can not have Verlet buffer tolerance of exactly 0");
392 if (ir->rlist < rc_max)
394 warning_error(wi, "With verlet lists rlist can not be smaller than rvdw or rcoulomb");
397 if (ir->rlist == rc_max && ir->nstlist > 1)
399 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.");
404 if (ir->rlist > rc_max)
406 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.");
409 if (ir->nstlist == 1)
411 /* No buffer required */
416 if (EI_DYNAMICS(ir->eI))
418 if (inputrec2nboundeddim(ir) < 3)
420 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.");
422 /* Set rlist temporarily so we can continue processing */
427 /* Set the buffer to 5% of the cut-off */
428 ir->rlist = (1.0 + verlet_buffer_ratio_nodynamics)*rc_max;
433 /* No twin-range calculations with Verlet lists */
434 ir->rlistlong = ir->rlist;
437 if (ir->nstcalclr == -1)
439 /* if rlist=rlistlong, this will later be changed to nstcalclr=0 */
440 ir->nstcalclr = ir->nstlist;
442 else if (ir->nstcalclr > 0)
444 if (ir->nstlist > 0 && (ir->nstlist % ir->nstcalclr != 0))
446 warning_error(wi, "nstlist must be evenly divisible by nstcalclr. Use nstcalclr = -1 to automatically follow nstlist");
449 else if (ir->nstcalclr < -1)
451 warning_error(wi, "nstcalclr must be a positive number (divisor of nstcalclr), or -1 to follow nstlist.");
454 if (EEL_PME(ir->coulombtype) && ir->rcoulomb > ir->rvdw && ir->nstcalclr > 1)
456 warning_error(wi, "When used with PME, the long-range component of twin-range interactions must be updated every step (nstcalclr)");
459 /* GENERAL INTEGRATOR STUFF */
460 if (!(ir->eI == eiMD || EI_VV(ir->eI)))
464 if (ir->eI == eiVVAK)
466 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]);
467 warning_note(wi, warn_buf);
469 if (!EI_DYNAMICS(ir->eI))
473 if (EI_DYNAMICS(ir->eI))
475 if (ir->nstcalcenergy < 0)
477 ir->nstcalcenergy = ir_optimal_nstcalcenergy(ir);
478 if (ir->nstenergy != 0 && ir->nstenergy < ir->nstcalcenergy)
480 /* nstcalcenergy larger than nstener does not make sense.
481 * We ideally want nstcalcenergy=nstener.
485 ir->nstcalcenergy = lcd(ir->nstenergy, ir->nstlist);
489 ir->nstcalcenergy = ir->nstenergy;
493 else if ( (ir->nstenergy > 0 && ir->nstcalcenergy > ir->nstenergy) ||
494 (ir->efep != efepNO && ir->fepvals->nstdhdl > 0 &&
495 (ir->nstcalcenergy > ir->fepvals->nstdhdl) ) )
498 const char *nsten = "nstenergy";
499 const char *nstdh = "nstdhdl";
500 const char *min_name = nsten;
501 int min_nst = ir->nstenergy;
503 /* find the smallest of ( nstenergy, nstdhdl ) */
504 if (ir->efep != efepNO && ir->fepvals->nstdhdl > 0 &&
505 (ir->nstenergy == 0 || ir->fepvals->nstdhdl < ir->nstenergy))
507 min_nst = ir->fepvals->nstdhdl;
510 /* If the user sets nstenergy small, we should respect that */
512 "Setting nstcalcenergy (%d) equal to %s (%d)",
513 ir->nstcalcenergy, min_name, min_nst);
514 warning_note(wi, warn_buf);
515 ir->nstcalcenergy = min_nst;
518 if (ir->epc != epcNO)
520 if (ir->nstpcouple < 0)
522 ir->nstpcouple = ir_optimal_nstpcouple(ir);
525 if (IR_TWINRANGE(*ir))
527 check_nst("nstlist", ir->nstlist,
528 "nstcalcenergy", &ir->nstcalcenergy, wi);
529 if (ir->epc != epcNO)
531 check_nst("nstlist", ir->nstlist,
532 "nstpcouple", &ir->nstpcouple, wi);
536 if (ir->nstcalcenergy > 0)
538 if (ir->efep != efepNO)
540 /* nstdhdl should be a multiple of nstcalcenergy */
541 check_nst("nstcalcenergy", ir->nstcalcenergy,
542 "nstdhdl", &ir->fepvals->nstdhdl, wi);
543 /* nstexpanded should be a multiple of nstcalcenergy */
544 check_nst("nstcalcenergy", ir->nstcalcenergy,
545 "nstexpanded", &ir->expandedvals->nstexpanded, wi);
547 /* for storing exact averages nstenergy should be
548 * a multiple of nstcalcenergy
550 check_nst("nstcalcenergy", ir->nstcalcenergy,
551 "nstenergy", &ir->nstenergy, wi);
556 if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
557 ir->bContinuation && ir->ld_seed != -1)
559 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)");
565 sprintf(err_buf, "TPI only works with pbc = %s", epbc_names[epbcXYZ]);
566 CHECK(ir->ePBC != epbcXYZ);
567 sprintf(err_buf, "TPI only works with ns = %s", ens_names[ensGRID]);
568 CHECK(ir->ns_type != ensGRID);
569 sprintf(err_buf, "with TPI nstlist should be larger than zero");
570 CHECK(ir->nstlist <= 0);
571 sprintf(err_buf, "TPI does not work with full electrostatics other than PME");
572 CHECK(EEL_FULL(ir->coulombtype) && !EEL_PME(ir->coulombtype));
576 if ( (opts->nshake > 0) && (opts->bMorse) )
579 "Using morse bond-potentials while constraining bonds is useless");
580 warning(wi, warn_buf);
583 if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
584 ir->bContinuation && ir->ld_seed != -1)
586 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)");
588 /* verify simulated tempering options */
592 gmx_bool bAllTempZero = TRUE;
593 for (i = 0; i < fep->n_lambda; i++)
595 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]);
596 CHECK((fep->all_lambda[efptTEMPERATURE][i] < 0) || (fep->all_lambda[efptTEMPERATURE][i] > 1));
597 if (fep->all_lambda[efptTEMPERATURE][i] > 0)
599 bAllTempZero = FALSE;
602 sprintf(err_buf, "if simulated tempering is on, temperature-lambdas may not be all zero");
603 CHECK(bAllTempZero == TRUE);
605 sprintf(err_buf, "Simulated tempering is currently only compatible with md-vv");
606 CHECK(ir->eI != eiVV);
608 /* check compatability of the temperature coupling with simulated tempering */
610 if (ir->etc == etcNOSEHOOVER)
612 sprintf(warn_buf, "Nose-Hoover based temperature control such as [%s] my not be entirelyconsistent with simulated tempering", etcoupl_names[ir->etc]);
613 warning_note(wi, warn_buf);
616 /* check that the temperatures make sense */
618 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);
619 CHECK(ir->simtempvals->simtemp_high <= ir->simtempvals->simtemp_low);
621 sprintf(err_buf, "Higher simulated tempering temperature (%g) must be >= zero", ir->simtempvals->simtemp_high);
622 CHECK(ir->simtempvals->simtemp_high <= 0);
624 sprintf(err_buf, "Lower simulated tempering temperature (%g) must be >= zero", ir->simtempvals->simtemp_low);
625 CHECK(ir->simtempvals->simtemp_low <= 0);
628 /* verify free energy options */
630 if (ir->efep != efepNO)
633 sprintf(err_buf, "The soft-core power is %d and can only be 1 or 2",
635 CHECK(fep->sc_alpha != 0 && fep->sc_power != 1 && fep->sc_power != 2);
637 sprintf(err_buf, "The soft-core sc-r-power is %d and can only be 6 or 48",
638 (int)fep->sc_r_power);
639 CHECK(fep->sc_alpha != 0 && fep->sc_r_power != 6.0 && fep->sc_r_power != 48.0);
641 sprintf(err_buf, "Can't use postive delta-lambda (%g) if initial state/lambda does not start at zero", fep->delta_lambda);
642 CHECK(fep->delta_lambda > 0 && ((fep->init_fep_state > 0) || (fep->init_lambda > 0)));
644 sprintf(err_buf, "Can't use postive delta-lambda (%g) with expanded ensemble simulations", fep->delta_lambda);
645 CHECK(fep->delta_lambda > 0 && (ir->efep == efepEXPANDED));
647 sprintf(err_buf, "Can only use expanded ensemble with md-vv for now; should be supported for other integrators in 5.0");
648 CHECK(!(EI_VV(ir->eI)) && (ir->efep == efepEXPANDED));
650 sprintf(err_buf, "Free-energy not implemented for Ewald");
651 CHECK(ir->coulombtype == eelEWALD);
653 /* check validty of lambda inputs */
654 if (fep->n_lambda == 0)
656 /* Clear output in case of no states:*/
657 sprintf(err_buf, "init-lambda-state set to %d: no lambda states are defined.", fep->init_fep_state);
658 CHECK((fep->init_fep_state >= 0) && (fep->n_lambda == 0));
662 sprintf(err_buf, "initial thermodynamic state %d does not exist, only goes to %d", fep->init_fep_state, fep->n_lambda-1);
663 CHECK((fep->init_fep_state >= fep->n_lambda));
666 sprintf(err_buf, "Lambda state must be set, either with init-lambda-state or with init-lambda");
667 CHECK((fep->init_fep_state < 0) && (fep->init_lambda < 0));
669 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",
670 fep->init_lambda, fep->init_fep_state);
671 CHECK((fep->init_fep_state >= 0) && (fep->init_lambda >= 0));
675 if ((fep->init_lambda >= 0) && (fep->delta_lambda == 0))
679 for (i = 0; i < efptNR; i++)
681 if (fep->separate_dvdl[i])
686 if (n_lambda_terms > 1)
688 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.");
689 warning(wi, warn_buf);
692 if (n_lambda_terms < 2 && fep->n_lambda > 0)
695 "init-lambda is deprecated for setting lambda state (except for slow growth). Use init-lambda-state instead.");
699 for (j = 0; j < efptNR; j++)
701 for (i = 0; i < fep->n_lambda; i++)
703 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]);
704 CHECK((fep->all_lambda[j][i] < 0) || (fep->all_lambda[j][i] > 1));
708 if ((fep->sc_alpha > 0) && (!fep->bScCoul))
710 for (i = 0; i < fep->n_lambda; i++)
712 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],
713 fep->all_lambda[efptCOUL][i]);
714 CHECK((fep->sc_alpha > 0) &&
715 (((fep->all_lambda[efptCOUL][i] > 0.0) &&
716 (fep->all_lambda[efptCOUL][i] < 1.0)) &&
717 ((fep->all_lambda[efptVDW][i] > 0.0) &&
718 (fep->all_lambda[efptVDW][i] < 1.0))));
722 if ((fep->bScCoul) && (EEL_PME(ir->coulombtype)))
724 real sigma, lambda, r_sc;
727 /* Maximum estimate for A and B charges equal with lambda power 1 */
729 r_sc = pow(lambda*fep->sc_alpha*pow(sigma/ir->rcoulomb, fep->sc_r_power) + 1.0, 1.0/fep->sc_r_power);
730 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.",
732 sigma, lambda, r_sc - 1.0, ir->ewald_rtol);
733 warning_note(wi, warn_buf);
736 /* Free Energy Checks -- In an ideal world, slow growth and FEP would
737 be treated differently, but that's the next step */
739 for (i = 0; i < efptNR; i++)
741 for (j = 0; j < fep->n_lambda; j++)
743 sprintf(err_buf, "%s[%d] must be between 0 and 1", efpt_names[i], j);
744 CHECK((fep->all_lambda[i][j] < 0) || (fep->all_lambda[i][j] > 1));
749 if ((ir->bSimTemp) || (ir->efep == efepEXPANDED))
752 expand = ir->expandedvals;
754 /* checking equilibration of weights inputs for validity */
756 sprintf(err_buf, "weight-equil-number-all-lambda (%d) is ignored if lmc-weights-equil is not equal to %s",
757 expand->equil_n_at_lam, elmceq_names[elmceqNUMATLAM]);
758 CHECK((expand->equil_n_at_lam > 0) && (expand->elmceq != elmceqNUMATLAM));
760 sprintf(err_buf, "weight-equil-number-samples (%d) is ignored if lmc-weights-equil is not equal to %s",
761 expand->equil_samples, elmceq_names[elmceqSAMPLES]);
762 CHECK((expand->equil_samples > 0) && (expand->elmceq != elmceqSAMPLES));
764 sprintf(err_buf, "weight-equil-number-steps (%d) is ignored if lmc-weights-equil is not equal to %s",
765 expand->equil_steps, elmceq_names[elmceqSTEPS]);
766 CHECK((expand->equil_steps > 0) && (expand->elmceq != elmceqSTEPS));
768 sprintf(err_buf, "weight-equil-wl-delta (%d) is ignored if lmc-weights-equil is not equal to %s",
769 expand->equil_samples, elmceq_names[elmceqWLDELTA]);
770 CHECK((expand->equil_wl_delta > 0) && (expand->elmceq != elmceqWLDELTA));
772 sprintf(err_buf, "weight-equil-count-ratio (%f) is ignored if lmc-weights-equil is not equal to %s",
773 expand->equil_ratio, elmceq_names[elmceqRATIO]);
774 CHECK((expand->equil_ratio > 0) && (expand->elmceq != elmceqRATIO));
776 sprintf(err_buf, "weight-equil-number-all-lambda (%d) must be a positive integer if lmc-weights-equil=%s",
777 expand->equil_n_at_lam, elmceq_names[elmceqNUMATLAM]);
778 CHECK((expand->equil_n_at_lam <= 0) && (expand->elmceq == elmceqNUMATLAM));
780 sprintf(err_buf, "weight-equil-number-samples (%d) must be a positive integer if lmc-weights-equil=%s",
781 expand->equil_samples, elmceq_names[elmceqSAMPLES]);
782 CHECK((expand->equil_samples <= 0) && (expand->elmceq == elmceqSAMPLES));
784 sprintf(err_buf, "weight-equil-number-steps (%d) must be a positive integer if lmc-weights-equil=%s",
785 expand->equil_steps, elmceq_names[elmceqSTEPS]);
786 CHECK((expand->equil_steps <= 0) && (expand->elmceq == elmceqSTEPS));
788 sprintf(err_buf, "weight-equil-wl-delta (%f) must be > 0 if lmc-weights-equil=%s",
789 expand->equil_wl_delta, elmceq_names[elmceqWLDELTA]);
790 CHECK((expand->equil_wl_delta <= 0) && (expand->elmceq == elmceqWLDELTA));
792 sprintf(err_buf, "weight-equil-count-ratio (%f) must be > 0 if lmc-weights-equil=%s",
793 expand->equil_ratio, elmceq_names[elmceqRATIO]);
794 CHECK((expand->equil_ratio <= 0) && (expand->elmceq == elmceqRATIO));
796 sprintf(err_buf, "lmc-weights-equil=%s only possible when lmc-stats = %s or lmc-stats %s",
797 elmceq_names[elmceqWLDELTA], elamstats_names[elamstatsWL], elamstats_names[elamstatsWWL]);
798 CHECK((expand->elmceq == elmceqWLDELTA) && (!EWL(expand->elamstats)));
800 sprintf(err_buf, "lmc-repeats (%d) must be greater than 0", expand->lmc_repeats);
801 CHECK((expand->lmc_repeats <= 0));
802 sprintf(err_buf, "minimum-var-min (%d) must be greater than 0", expand->minvarmin);
803 CHECK((expand->minvarmin <= 0));
804 sprintf(err_buf, "weight-c-range (%d) must be greater or equal to 0", expand->c_range);
805 CHECK((expand->c_range < 0));
806 sprintf(err_buf, "init-lambda-state (%d) must be zero if lmc-forced-nstart (%d)> 0 and lmc-move != 'no'",
807 fep->init_fep_state, expand->lmc_forced_nstart);
808 CHECK((fep->init_fep_state != 0) && (expand->lmc_forced_nstart > 0) && (expand->elmcmove != elmcmoveNO));
809 sprintf(err_buf, "lmc-forced-nstart (%d) must not be negative", expand->lmc_forced_nstart);
810 CHECK((expand->lmc_forced_nstart < 0));
811 sprintf(err_buf, "init-lambda-state (%d) must be in the interval [0,number of lambdas)", fep->init_fep_state);
812 CHECK((fep->init_fep_state < 0) || (fep->init_fep_state >= fep->n_lambda));
814 sprintf(err_buf, "init-wl-delta (%f) must be greater than or equal to 0", expand->init_wl_delta);
815 CHECK((expand->init_wl_delta < 0));
816 sprintf(err_buf, "wl-ratio (%f) must be between 0 and 1", expand->wl_ratio);
817 CHECK((expand->wl_ratio <= 0) || (expand->wl_ratio >= 1));
818 sprintf(err_buf, "wl-scale (%f) must be between 0 and 1", expand->wl_scale);
819 CHECK((expand->wl_scale <= 0) || (expand->wl_scale >= 1));
821 /* if there is no temperature control, we need to specify an MC temperature */
822 sprintf(err_buf, "If there is no temperature control, and lmc-mcmove!= 'no',mc_temperature must be set to a positive number");
823 if (expand->nstTij > 0)
825 sprintf(err_buf, "nst-transition-matrix (%d) must be an integer multiple of nstlog (%d)",
826 expand->nstTij, ir->nstlog);
827 CHECK((mod(expand->nstTij, ir->nstlog) != 0));
832 sprintf(err_buf, "walls only work with pbc=%s", epbc_names[epbcXY]);
833 CHECK(ir->nwall && ir->ePBC != epbcXY);
836 if (ir->ePBC != epbcXYZ && ir->nwall != 2)
838 if (ir->ePBC == epbcNONE)
840 if (ir->epc != epcNO)
842 warning(wi, "Turning off pressure coupling for vacuum system");
848 sprintf(err_buf, "Can not have pressure coupling with pbc=%s",
849 epbc_names[ir->ePBC]);
850 CHECK(ir->epc != epcNO);
852 sprintf(err_buf, "Can not have Ewald with pbc=%s", epbc_names[ir->ePBC]);
853 CHECK(EEL_FULL(ir->coulombtype));
855 sprintf(err_buf, "Can not have dispersion correction with pbc=%s",
856 epbc_names[ir->ePBC]);
857 CHECK(ir->eDispCorr != edispcNO);
860 if (ir->rlist == 0.0)
862 sprintf(err_buf, "can only have neighborlist cut-off zero (=infinite)\n"
863 "with coulombtype = %s or coulombtype = %s\n"
864 "without periodic boundary conditions (pbc = %s) and\n"
865 "rcoulomb and rvdw set to zero",
866 eel_names[eelCUT], eel_names[eelUSER], epbc_names[epbcNONE]);
867 CHECK(((ir->coulombtype != eelCUT) && (ir->coulombtype != eelUSER)) ||
868 (ir->ePBC != epbcNONE) ||
869 (ir->rcoulomb != 0.0) || (ir->rvdw != 0.0));
873 warning_error(wi, "Can not have heuristic neighborlist updates without cut-off");
877 warning_note(wi, "Simulating without cut-offs can be (slightly) faster with nstlist=0, nstype=simple and only one MPI rank");
882 if (ir->nstcomm == 0)
884 ir->comm_mode = ecmNO;
886 if (ir->comm_mode != ecmNO)
890 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");
891 ir->nstcomm = abs(ir->nstcomm);
894 if (ir->nstcalcenergy > 0 && ir->nstcomm < ir->nstcalcenergy)
896 warning_note(wi, "nstcomm < nstcalcenergy defeats the purpose of nstcalcenergy, setting nstcomm to nstcalcenergy");
897 ir->nstcomm = ir->nstcalcenergy;
900 if (ir->comm_mode == ecmANGULAR)
902 sprintf(err_buf, "Can not remove the rotation around the center of mass with periodic molecules");
903 CHECK(ir->bPeriodicMols);
904 if (ir->ePBC != epbcNONE)
906 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).");
911 if (EI_STATE_VELOCITY(ir->eI) && ir->ePBC == epbcNONE && ir->comm_mode != ecmANGULAR)
913 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.");
916 sprintf(err_buf, "Twin-range neighbour searching (NS) with simple NS"
917 " algorithm not implemented");
918 CHECK(((ir->rcoulomb > ir->rlist) || (ir->rvdw > ir->rlist))
919 && (ir->ns_type == ensSIMPLE));
921 /* TEMPERATURE COUPLING */
922 if (ir->etc == etcYES)
924 ir->etc = etcBERENDSEN;
925 warning_note(wi, "Old option for temperature coupling given: "
926 "changing \"yes\" to \"Berendsen\"\n");
929 if ((ir->etc == etcNOSEHOOVER) || (ir->epc == epcMTTK))
931 if (ir->opts.nhchainlength < 1)
933 sprintf(warn_buf, "number of Nose-Hoover chains (currently %d) cannot be less than 1,reset to 1\n", ir->opts.nhchainlength);
934 ir->opts.nhchainlength = 1;
935 warning(wi, warn_buf);
938 if (ir->etc == etcNOSEHOOVER && !EI_VV(ir->eI) && ir->opts.nhchainlength > 1)
940 warning_note(wi, "leapfrog does not yet support Nose-Hoover chains, nhchainlength reset to 1");
941 ir->opts.nhchainlength = 1;
946 ir->opts.nhchainlength = 0;
949 if (ir->eI == eiVVAK)
951 sprintf(err_buf, "%s implemented primarily for validation, and requires nsttcouple = 1 and nstpcouple = 1.",
953 CHECK((ir->nsttcouple != 1) || (ir->nstpcouple != 1));
956 if (ETC_ANDERSEN(ir->etc))
958 sprintf(err_buf, "%s temperature control not supported for integrator %s.", etcoupl_names[ir->etc], ei_names[ir->eI]);
959 CHECK(!(EI_VV(ir->eI)));
961 if (ir->nstcomm > 0 && (ir->etc == etcANDERSEN))
963 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]);
964 warning_note(wi, warn_buf);
967 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]);
968 CHECK(ir->nstcomm > 1 && (ir->etc == etcANDERSEN));
971 if (ir->etc == etcBERENDSEN)
973 sprintf(warn_buf, "The %s thermostat does not generate the correct kinetic energy distribution. You might want to consider using the %s thermostat.",
974 ETCOUPLTYPE(ir->etc), ETCOUPLTYPE(etcVRESCALE));
975 warning_note(wi, warn_buf);
978 if ((ir->etc == etcNOSEHOOVER || ETC_ANDERSEN(ir->etc))
979 && ir->epc == epcBERENDSEN)
981 sprintf(warn_buf, "Using Berendsen pressure coupling invalidates the "
982 "true ensemble for the thermostat");
983 warning(wi, warn_buf);
986 /* PRESSURE COUPLING */
987 if (ir->epc == epcISOTROPIC)
989 ir->epc = epcBERENDSEN;
990 warning_note(wi, "Old option for pressure coupling given: "
991 "changing \"Isotropic\" to \"Berendsen\"\n");
994 if (ir->epc != epcNO)
996 dt_pcoupl = ir->nstpcouple*ir->delta_t;
998 sprintf(err_buf, "tau-p must be > 0 instead of %g\n", ir->tau_p);
999 CHECK(ir->tau_p <= 0);
1001 if (ir->tau_p/dt_pcoupl < pcouple_min_integration_steps(ir->epc))
1003 sprintf(warn_buf, "For proper integration of the %s barostat, tau-p (%g) should be at least %d times larger than nstpcouple*dt (%g)",
1004 EPCOUPLTYPE(ir->epc), ir->tau_p, pcouple_min_integration_steps(ir->epc), dt_pcoupl);
1005 warning(wi, warn_buf);
1008 sprintf(err_buf, "compressibility must be > 0 when using pressure"
1009 " coupling %s\n", EPCOUPLTYPE(ir->epc));
1010 CHECK(ir->compress[XX][XX] < 0 || ir->compress[YY][YY] < 0 ||
1011 ir->compress[ZZ][ZZ] < 0 ||
1012 (trace(ir->compress) == 0 && ir->compress[YY][XX] <= 0 &&
1013 ir->compress[ZZ][XX] <= 0 && ir->compress[ZZ][YY] <= 0));
1015 if (epcPARRINELLORAHMAN == ir->epc && opts->bGenVel)
1018 "You are generating velocities so I am assuming you "
1019 "are equilibrating a system. You are using "
1020 "%s pressure coupling, but this can be "
1021 "unstable for equilibration. If your system crashes, try "
1022 "equilibrating first with Berendsen pressure coupling. If "
1023 "you are not equilibrating the system, you can probably "
1024 "ignore this warning.",
1025 epcoupl_names[ir->epc]);
1026 warning(wi, warn_buf);
1032 if (ir->epc > epcNO)
1034 if ((ir->epc != epcBERENDSEN) && (ir->epc != epcMTTK))
1036 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.");
1042 if (ir->epc == epcMTTK)
1044 warning_error(wi, "MTTK pressure coupling requires a Velocity-verlet integrator");
1048 /* ELECTROSTATICS */
1049 /* More checks are in triple check (grompp.c) */
1051 if (ir->coulombtype == eelSWITCH)
1053 sprintf(warn_buf, "coulombtype = %s is only for testing purposes and can lead to serious "
1054 "artifacts, advice: use coulombtype = %s",
1055 eel_names[ir->coulombtype],
1056 eel_names[eelRF_ZERO]);
1057 warning(wi, warn_buf);
1060 if (ir->epsilon_r != 1 && ir->implicit_solvent == eisGBSA)
1062 sprintf(warn_buf, "epsilon-r = %g with GB implicit solvent, will use this value for inner dielectric", ir->epsilon_r);
1063 warning_note(wi, warn_buf);
1066 if (EEL_RF(ir->coulombtype) && ir->epsilon_rf == 1 && ir->epsilon_r != 1)
1068 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);
1069 warning(wi, warn_buf);
1070 ir->epsilon_rf = ir->epsilon_r;
1071 ir->epsilon_r = 1.0;
1074 if (getenv("GMX_DO_GALACTIC_DYNAMICS") == NULL)
1076 sprintf(err_buf, "epsilon-r must be >= 0 instead of %g\n", ir->epsilon_r);
1077 CHECK(ir->epsilon_r < 0);
1080 if (EEL_RF(ir->coulombtype))
1082 /* reaction field (at the cut-off) */
1084 if (ir->coulombtype == eelRF_ZERO)
1086 sprintf(warn_buf, "With coulombtype = %s, epsilon-rf must be 0, assuming you meant epsilon_rf=0",
1087 eel_names[ir->coulombtype]);
1088 CHECK(ir->epsilon_rf != 0);
1089 ir->epsilon_rf = 0.0;
1092 sprintf(err_buf, "epsilon-rf must be >= epsilon-r");
1093 CHECK((ir->epsilon_rf < ir->epsilon_r && ir->epsilon_rf != 0) ||
1094 (ir->epsilon_r == 0));
1095 if (ir->epsilon_rf == ir->epsilon_r)
1097 sprintf(warn_buf, "Using epsilon-rf = epsilon-r with %s does not make sense",
1098 eel_names[ir->coulombtype]);
1099 warning(wi, warn_buf);
1102 /* Allow rlist>rcoulomb for tabulated long range stuff. This just
1103 * means the interaction is zero outside rcoulomb, but it helps to
1104 * provide accurate energy conservation.
1106 if (ir_coulomb_might_be_zero_at_cutoff(ir))
1108 if (ir_coulomb_switched(ir))
1111 "With coulombtype = %s rcoulomb_switch must be < rcoulomb. Or, better: Use the potential modifier options!",
1112 eel_names[ir->coulombtype]);
1113 CHECK(ir->rcoulomb_switch >= ir->rcoulomb);
1116 else if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
1118 if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE)
1120 sprintf(err_buf, "With coulombtype = %s, rcoulomb should be >= rlist unless you use a potential modifier",
1121 eel_names[ir->coulombtype]);
1122 CHECK(ir->rlist > ir->rcoulomb);
1126 if (ir->coulombtype == eelSWITCH || ir->coulombtype == eelSHIFT ||
1127 ir->vdwtype == evdwSWITCH || ir->vdwtype == evdwSHIFT)
1130 "The switch/shift interaction settings are just for compatibility; you will get better "
1131 "performance from applying potential modifiers to your interactions!\n");
1132 warning_note(wi, warn_buf);
1135 if (ir->coulombtype == eelPMESWITCH || ir->coulomb_modifier == eintmodPOTSWITCH)
1137 if (ir->rcoulomb_switch/ir->rcoulomb < 0.9499)
1139 real percentage = 100*(ir->rcoulomb-ir->rcoulomb_switch)/ir->rcoulomb;
1140 sprintf(warn_buf, "The switching range for should be 5%% or less (currently %.2f%% using a switching range of %4f-%4f) for accurate electrostatic energies, energy conservation will be good regardless, since ewald_rtol = %g.",
1141 percentage, ir->rcoulomb_switch, ir->rcoulomb, ir->ewald_rtol);
1142 warning(wi, warn_buf);
1146 if (ir->vdwtype == evdwSWITCH || ir->vdw_modifier == eintmodPOTSWITCH)
1148 if (ir->rvdw_switch == 0)
1150 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.");
1151 warning(wi, warn_buf);
1155 if (EEL_FULL(ir->coulombtype))
1157 if (ir->coulombtype == eelPMESWITCH || ir->coulombtype == eelPMEUSER ||
1158 ir->coulombtype == eelPMEUSERSWITCH)
1160 sprintf(err_buf, "With coulombtype = %s, rcoulomb must be <= rlist",
1161 eel_names[ir->coulombtype]);
1162 CHECK(ir->rcoulomb > ir->rlist);
1164 else if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE)
1166 if (ir->coulombtype == eelPME || ir->coulombtype == eelP3M_AD)
1169 "With coulombtype = %s (without modifier), rcoulomb must be equal to rlist,\n"
1170 "or rlistlong if nstcalclr=1. For optimal energy conservation,consider using\n"
1171 "a potential modifier.", eel_names[ir->coulombtype]);
1172 if (ir->nstcalclr == 1)
1174 CHECK(ir->rcoulomb != ir->rlist && ir->rcoulomb != ir->rlistlong);
1178 CHECK(ir->rcoulomb != ir->rlist);
1184 if (EEL_PME(ir->coulombtype) || EVDW_PME(ir->vdwtype))
1186 if (ir->pme_order < 3)
1188 warning_error(wi, "pme-order can not be smaller than 3");
1192 if (ir->nwall == 2 && EEL_FULL(ir->coulombtype))
1194 if (ir->ewald_geometry == eewg3D)
1196 sprintf(warn_buf, "With pbc=%s you should use ewald-geometry=%s",
1197 epbc_names[ir->ePBC], eewg_names[eewg3DC]);
1198 warning(wi, warn_buf);
1200 /* This check avoids extra pbc coding for exclusion corrections */
1201 sprintf(err_buf, "wall-ewald-zfac should be >= 2");
1202 CHECK(ir->wall_ewald_zfac < 2);
1205 if (ir_vdw_switched(ir))
1207 sprintf(err_buf, "With switched vdw forces or potentials, rvdw-switch must be < rvdw");
1208 CHECK(ir->rvdw_switch >= ir->rvdw);
1210 if (ir->rvdw_switch < 0.5*ir->rvdw)
1212 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.",
1213 ir->rvdw_switch, ir->rvdw);
1214 warning_note(wi, warn_buf);
1217 else if (ir->vdwtype == evdwCUT || ir->vdwtype == evdwPME)
1219 if (ir->cutoff_scheme == ecutsGROUP && ir->vdw_modifier == eintmodNONE)
1221 sprintf(err_buf, "With vdwtype = %s, rvdw must be >= rlist unless you use a potential modifier", evdw_names[ir->vdwtype]);
1222 CHECK(ir->rlist > ir->rvdw);
1226 if (ir->vdwtype == evdwPME)
1228 if (!(ir->vdw_modifier == eintmodNONE || ir->vdw_modifier == eintmodPOTSHIFT))
1230 sprintf(err_buf, "With vdwtype = %s, the only supported modifiers are %s a\
1232 evdw_names[ir->vdwtype],
1233 eintmod_names[eintmodPOTSHIFT],
1234 eintmod_names[eintmodNONE]);
1238 if (ir->cutoff_scheme == ecutsGROUP)
1240 if (((ir->coulomb_modifier != eintmodNONE && ir->rcoulomb == ir->rlist) ||
1241 (ir->vdw_modifier != eintmodNONE && ir->rvdw == ir->rlist)) &&
1244 warning_note(wi, "With exact cut-offs, rlist should be "
1245 "larger than rcoulomb and rvdw, so that there "
1246 "is a buffer region for particle motion "
1247 "between neighborsearch steps");
1250 if (ir_coulomb_is_zero_at_cutoff(ir) && ir->rlistlong <= ir->rcoulomb)
1252 sprintf(warn_buf, "For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rcoulomb.",
1253 IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
1254 warning_note(wi, warn_buf);
1256 if (ir_vdw_switched(ir) && (ir->rlistlong <= ir->rvdw))
1258 sprintf(warn_buf, "For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rvdw.",
1259 IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
1260 warning_note(wi, warn_buf);
1264 if (ir->vdwtype == evdwUSER && ir->eDispCorr != edispcNO)
1266 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.");
1269 if (ir->nstlist == -1)
1271 sprintf(err_buf, "With nstlist=-1 rvdw and rcoulomb should be smaller than rlist to account for diffusion and possibly charge-group radii");
1272 CHECK(ir->rvdw >= ir->rlist || ir->rcoulomb >= ir->rlist);
1274 sprintf(err_buf, "nstlist can not be smaller than -1");
1275 CHECK(ir->nstlist < -1);
1277 if (ir->eI == eiLBFGS && (ir->coulombtype == eelCUT || ir->vdwtype == evdwCUT)
1280 warning(wi, "For efficient BFGS minimization, use switch/shift/pme instead of cut-off.");
1283 if (ir->eI == eiLBFGS && ir->nbfgscorr <= 0)
1285 warning(wi, "Using L-BFGS with nbfgscorr<=0 just gets you steepest descent.");
1288 /* ENERGY CONSERVATION */
1289 if (ir_NVE(ir) && ir->cutoff_scheme == ecutsGROUP)
1291 if (!ir_vdw_might_be_zero_at_cutoff(ir) && ir->rvdw > 0 && ir->vdw_modifier == eintmodNONE)
1293 sprintf(warn_buf, "You are using a cut-off for VdW interactions with NVE, for good energy conservation use vdwtype = %s (possibly with DispCorr)",
1294 evdw_names[evdwSHIFT]);
1295 warning_note(wi, warn_buf);
1297 if (!ir_coulomb_might_be_zero_at_cutoff(ir) && ir->rcoulomb > 0)
1299 sprintf(warn_buf, "You are using a cut-off for electrostatics with NVE, for good energy conservation use coulombtype = %s or %s",
1300 eel_names[eelPMESWITCH], eel_names[eelRF_ZERO]);
1301 warning_note(wi, warn_buf);
1305 if (EI_VV(ir->eI) && IR_TWINRANGE(*ir) && ir->nstlist > 1)
1307 sprintf(warn_buf, "Twin-range multiple time stepping does not work with integrator %s.", ei_names[ir->eI]);
1308 warning_error(wi, warn_buf);
1311 /* IMPLICIT SOLVENT */
1312 if (ir->coulombtype == eelGB_NOTUSED)
1314 ir->coulombtype = eelCUT;
1315 ir->implicit_solvent = eisGBSA;
1316 fprintf(stderr, "Note: Old option for generalized born electrostatics given:\n"
1317 "Changing coulombtype from \"generalized-born\" to \"cut-off\" and instead\n"
1318 "setting implicit-solvent value to \"GBSA\" in input section.\n");
1321 if (ir->sa_algorithm == esaSTILL)
1323 sprintf(err_buf, "Still SA algorithm not available yet, use %s or %s instead\n", esa_names[esaAPPROX], esa_names[esaNO]);
1324 CHECK(ir->sa_algorithm == esaSTILL);
1327 if (ir->implicit_solvent == eisGBSA)
1329 sprintf(err_buf, "With GBSA implicit solvent, rgbradii must be equal to rlist.");
1330 CHECK(ir->rgbradii != ir->rlist);
1332 if (ir->coulombtype != eelCUT)
1334 sprintf(err_buf, "With GBSA, coulombtype must be equal to %s\n", eel_names[eelCUT]);
1335 CHECK(ir->coulombtype != eelCUT);
1337 if (ir->vdwtype != evdwCUT)
1339 sprintf(err_buf, "With GBSA, vdw-type must be equal to %s\n", evdw_names[evdwCUT]);
1340 CHECK(ir->vdwtype != evdwCUT);
1342 if (ir->nstgbradii < 1)
1344 sprintf(warn_buf, "Using GBSA with nstgbradii<1, setting nstgbradii=1");
1345 warning_note(wi, warn_buf);
1348 if (ir->sa_algorithm == esaNO)
1350 sprintf(warn_buf, "No SA (non-polar) calculation requested together with GB. Are you sure this is what you want?\n");
1351 warning_note(wi, warn_buf);
1353 if (ir->sa_surface_tension < 0 && ir->sa_algorithm != esaNO)
1355 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");
1356 warning_note(wi, warn_buf);
1358 if (ir->gb_algorithm == egbSTILL)
1360 ir->sa_surface_tension = 0.0049 * CAL2JOULE * 100;
1364 ir->sa_surface_tension = 0.0054 * CAL2JOULE * 100;
1367 if (ir->sa_surface_tension == 0 && ir->sa_algorithm != esaNO)
1369 sprintf(err_buf, "Surface tension set to 0 while SA-calculation requested\n");
1370 CHECK(ir->sa_surface_tension == 0 && ir->sa_algorithm != esaNO);
1377 if (ir->cutoff_scheme != ecutsGROUP)
1379 warning_error(wi, "AdresS simulation supports only cutoff-scheme=group");
1383 warning_error(wi, "AdresS simulation supports only stochastic dynamics");
1385 if (ir->epc != epcNO)
1387 warning_error(wi, "AdresS simulation does not support pressure coupling");
1389 if (EEL_FULL(ir->coulombtype))
1391 warning_error(wi, "AdresS simulation does not support long-range electrostatics");
1396 /* count the number of text elemets separated by whitespace in a string.
1397 str = the input string
1398 maxptr = the maximum number of allowed elements
1399 ptr = the output array of pointers to the first character of each element
1400 returns: the number of elements. */
1401 int str_nelem(const char *str, int maxptr, char *ptr[])
1406 copy0 = strdup(str);
1409 while (*copy != '\0')
1413 gmx_fatal(FARGS, "Too many groups on line: '%s' (max is %d)",
1421 while ((*copy != '\0') && !isspace(*copy))
1440 /* interpret a number of doubles from a string and put them in an array,
1441 after allocating space for them.
1442 str = the input string
1443 n = the (pre-allocated) number of doubles read
1444 r = the output array of doubles. */
1445 static void parse_n_real(char *str, int *n, real **r)
1450 *n = str_nelem(str, MAXPTR, ptr);
1453 for (i = 0; i < *n; i++)
1455 (*r)[i] = strtod(ptr[i], NULL);
1459 static void do_fep_params(t_inputrec *ir, char fep_lambda[][STRLEN], char weights[STRLEN])
1462 int i, j, max_n_lambda, nweights, nfep[efptNR];
1463 t_lambda *fep = ir->fepvals;
1464 t_expanded *expand = ir->expandedvals;
1465 real **count_fep_lambdas;
1466 gmx_bool bOneLambda = TRUE;
1468 snew(count_fep_lambdas, efptNR);
1470 /* FEP input processing */
1471 /* first, identify the number of lambda values for each type.
1472 All that are nonzero must have the same number */
1474 for (i = 0; i < efptNR; i++)
1476 parse_n_real(fep_lambda[i], &(nfep[i]), &(count_fep_lambdas[i]));
1479 /* now, determine the number of components. All must be either zero, or equal. */
1482 for (i = 0; i < efptNR; i++)
1484 if (nfep[i] > max_n_lambda)
1486 max_n_lambda = nfep[i]; /* here's a nonzero one. All of them
1487 must have the same number if its not zero.*/
1492 for (i = 0; i < efptNR; i++)
1496 ir->fepvals->separate_dvdl[i] = FALSE;
1498 else if (nfep[i] == max_n_lambda)
1500 if (i != efptTEMPERATURE) /* we treat this differently -- not really a reason to compute the derivative with
1501 respect to the temperature currently */
1503 ir->fepvals->separate_dvdl[i] = TRUE;
1508 gmx_fatal(FARGS, "Number of lambdas (%d) for FEP type %s not equal to number of other types (%d)",
1509 nfep[i], efpt_names[i], max_n_lambda);
1512 /* we don't print out dhdl if the temperature is changing, since we can't correctly define dhdl in this case */
1513 ir->fepvals->separate_dvdl[efptTEMPERATURE] = FALSE;
1515 /* the number of lambdas is the number we've read in, which is either zero
1516 or the same for all */
1517 fep->n_lambda = max_n_lambda;
1519 /* allocate space for the array of lambda values */
1520 snew(fep->all_lambda, efptNR);
1521 /* if init_lambda is defined, we need to set lambda */
1522 if ((fep->init_lambda > 0) && (fep->n_lambda == 0))
1524 ir->fepvals->separate_dvdl[efptFEP] = TRUE;
1526 /* otherwise allocate the space for all of the lambdas, and transfer the data */
1527 for (i = 0; i < efptNR; i++)
1529 snew(fep->all_lambda[i], fep->n_lambda);
1530 if (nfep[i] > 0) /* if it's zero, then the count_fep_lambda arrays
1533 for (j = 0; j < fep->n_lambda; j++)
1535 fep->all_lambda[i][j] = (double)count_fep_lambdas[i][j];
1537 sfree(count_fep_lambdas[i]);
1540 sfree(count_fep_lambdas);
1542 /* "fep-vals" is either zero or the full number. If zero, we'll need to define fep-lambdas for internal
1543 bookkeeping -- for now, init_lambda */
1545 if ((nfep[efptFEP] == 0) && (fep->init_lambda >= 0))
1547 for (i = 0; i < fep->n_lambda; i++)
1549 fep->all_lambda[efptFEP][i] = fep->init_lambda;
1553 /* check to see if only a single component lambda is defined, and soft core is defined.
1554 In this case, turn on coulomb soft core */
1556 if (max_n_lambda == 0)
1562 for (i = 0; i < efptNR; i++)
1564 if ((nfep[i] != 0) && (i != efptFEP))
1570 if ((bOneLambda) && (fep->sc_alpha > 0))
1572 fep->bScCoul = TRUE;
1575 /* Fill in the others with the efptFEP if they are not explicitly
1576 specified (i.e. nfep[i] == 0). This means if fep is not defined,
1577 they are all zero. */
1579 for (i = 0; i < efptNR; i++)
1581 if ((nfep[i] == 0) && (i != efptFEP))
1583 for (j = 0; j < fep->n_lambda; j++)
1585 fep->all_lambda[i][j] = fep->all_lambda[efptFEP][j];
1591 /* make it easier if sc_r_power = 48 by increasing it to the 4th power, to be in the right scale. */
1592 if (fep->sc_r_power == 48)
1594 if (fep->sc_alpha > 0.1)
1596 gmx_fatal(FARGS, "sc_alpha (%f) for sc_r_power = 48 should usually be between 0.001 and 0.004", fep->sc_alpha);
1600 expand = ir->expandedvals;
1601 /* now read in the weights */
1602 parse_n_real(weights, &nweights, &(expand->init_lambda_weights));
1605 snew(expand->init_lambda_weights, fep->n_lambda); /* initialize to zero */
1607 else if (nweights != fep->n_lambda)
1609 gmx_fatal(FARGS, "Number of weights (%d) is not equal to number of lambda values (%d)",
1610 nweights, fep->n_lambda);
1612 if ((expand->nstexpanded < 0) && (ir->efep != efepNO))
1614 expand->nstexpanded = fep->nstdhdl;
1615 /* if you don't specify nstexpanded when doing expanded ensemble free energy calcs, it is set to nstdhdl */
1617 if ((expand->nstexpanded < 0) && ir->bSimTemp)
1619 expand->nstexpanded = 2*(int)(ir->opts.tau_t[0]/ir->delta_t);
1620 /* if you don't specify nstexpanded when doing expanded ensemble simulated tempering, it is set to
1621 2*tau_t just to be careful so it's not to frequent */
1626 static void do_simtemp_params(t_inputrec *ir)
1629 snew(ir->simtempvals->temperatures, ir->fepvals->n_lambda);
1630 GetSimTemps(ir->fepvals->n_lambda, ir->simtempvals, ir->fepvals->all_lambda[efptTEMPERATURE]);
1635 static void do_wall_params(t_inputrec *ir,
1636 char *wall_atomtype, char *wall_density,
1640 char *names[MAXPTR];
1643 opts->wall_atomtype[0] = NULL;
1644 opts->wall_atomtype[1] = NULL;
1646 ir->wall_atomtype[0] = -1;
1647 ir->wall_atomtype[1] = -1;
1648 ir->wall_density[0] = 0;
1649 ir->wall_density[1] = 0;
1653 nstr = str_nelem(wall_atomtype, MAXPTR, names);
1654 if (nstr != ir->nwall)
1656 gmx_fatal(FARGS, "Expected %d elements for wall_atomtype, found %d",
1659 for (i = 0; i < ir->nwall; i++)
1661 opts->wall_atomtype[i] = strdup(names[i]);
1664 if (ir->wall_type == ewt93 || ir->wall_type == ewt104)
1666 nstr = str_nelem(wall_density, MAXPTR, names);
1667 if (nstr != ir->nwall)
1669 gmx_fatal(FARGS, "Expected %d elements for wall-density, found %d", ir->nwall, nstr);
1671 for (i = 0; i < ir->nwall; i++)
1673 sscanf(names[i], "%lf", &dbl);
1676 gmx_fatal(FARGS, "wall-density[%d] = %f\n", i, dbl);
1678 ir->wall_density[i] = dbl;
1684 static void add_wall_energrps(gmx_groups_t *groups, int nwall, t_symtab *symtab)
1692 srenew(groups->grpname, groups->ngrpname+nwall);
1693 grps = &(groups->grps[egcENER]);
1694 srenew(grps->nm_ind, grps->nr+nwall);
1695 for (i = 0; i < nwall; i++)
1697 sprintf(str, "wall%d", i);
1698 groups->grpname[groups->ngrpname] = put_symtab(symtab, str);
1699 grps->nm_ind[grps->nr++] = groups->ngrpname++;
1704 void read_expandedparams(int *ninp_p, t_inpfile **inp_p,
1705 t_expanded *expand, warninp_t wi)
1707 int ninp, nerror = 0;
1713 /* read expanded ensemble parameters */
1714 CCTYPE ("expanded ensemble variables");
1715 ITYPE ("nstexpanded", expand->nstexpanded, -1);
1716 EETYPE("lmc-stats", expand->elamstats, elamstats_names);
1717 EETYPE("lmc-move", expand->elmcmove, elmcmove_names);
1718 EETYPE("lmc-weights-equil", expand->elmceq, elmceq_names);
1719 ITYPE ("weight-equil-number-all-lambda", expand->equil_n_at_lam, -1);
1720 ITYPE ("weight-equil-number-samples", expand->equil_samples, -1);
1721 ITYPE ("weight-equil-number-steps", expand->equil_steps, -1);
1722 RTYPE ("weight-equil-wl-delta", expand->equil_wl_delta, -1);
1723 RTYPE ("weight-equil-count-ratio", expand->equil_ratio, -1);
1724 CCTYPE("Seed for Monte Carlo in lambda space");
1725 ITYPE ("lmc-seed", expand->lmc_seed, -1);
1726 RTYPE ("mc-temperature", expand->mc_temp, -1);
1727 ITYPE ("lmc-repeats", expand->lmc_repeats, 1);
1728 ITYPE ("lmc-gibbsdelta", expand->gibbsdeltalam, -1);
1729 ITYPE ("lmc-forced-nstart", expand->lmc_forced_nstart, 0);
1730 EETYPE("symmetrized-transition-matrix", expand->bSymmetrizedTMatrix, yesno_names);
1731 ITYPE("nst-transition-matrix", expand->nstTij, -1);
1732 ITYPE ("mininum-var-min", expand->minvarmin, 100); /*default is reasonable */
1733 ITYPE ("weight-c-range", expand->c_range, 0); /* default is just C=0 */
1734 RTYPE ("wl-scale", expand->wl_scale, 0.8);
1735 RTYPE ("wl-ratio", expand->wl_ratio, 0.8);
1736 RTYPE ("init-wl-delta", expand->init_wl_delta, 1.0);
1737 EETYPE("wl-oneovert", expand->bWLoneovert, yesno_names);
1745 void get_ir(const char *mdparin, const char *mdparout,
1746 t_inputrec *ir, t_gromppopts *opts,
1750 double dumdub[2][6];
1754 char warn_buf[STRLEN];
1755 t_lambda *fep = ir->fepvals;
1756 t_expanded *expand = ir->expandedvals;
1758 init_inputrec_strings();
1759 inp = read_inpfile(mdparin, &ninp, wi);
1761 snew(dumstr[0], STRLEN);
1762 snew(dumstr[1], STRLEN);
1764 if (-1 == search_einp(ninp, inp, "cutoff-scheme"))
1767 "%s did not specify a value for the .mdp option "
1768 "\"cutoff-scheme\". Probably it was first intended for use "
1769 "with GROMACS before 4.6. In 4.6, the Verlet scheme was "
1770 "introduced, but the group scheme was still the default. "
1771 "The default is now the Verlet scheme, so you will observe "
1772 "different behaviour.", mdparin);
1773 warning_note(wi, warn_buf);
1776 /* remove the following deprecated commands */
1779 REM_TYPE("domain-decomposition");
1780 REM_TYPE("andersen-seed");
1782 REM_TYPE("dihre-fc");
1783 REM_TYPE("dihre-tau");
1784 REM_TYPE("nstdihreout");
1785 REM_TYPE("nstcheckpoint");
1787 /* replace the following commands with the clearer new versions*/
1788 REPL_TYPE("unconstrained-start", "continuation");
1789 REPL_TYPE("foreign-lambda", "fep-lambdas");
1790 REPL_TYPE("verlet-buffer-drift", "verlet-buffer-tolerance");
1791 REPL_TYPE("nstxtcout", "nstxout-compressed");
1792 REPL_TYPE("xtc-grps", "compressed-x-grps");
1793 REPL_TYPE("xtc-precision", "compressed-x-precision");
1795 CCTYPE ("VARIOUS PREPROCESSING OPTIONS");
1796 CTYPE ("Preprocessor information: use cpp syntax.");
1797 CTYPE ("e.g.: -I/home/joe/doe -I/home/mary/roe");
1798 STYPE ("include", opts->include, NULL);
1799 CTYPE ("e.g.: -DPOSRES -DFLEXIBLE (note these variable names are case sensitive)");
1800 STYPE ("define", opts->define, NULL);
1802 CCTYPE ("RUN CONTROL PARAMETERS");
1803 EETYPE("integrator", ir->eI, ei_names);
1804 CTYPE ("Start time and timestep in ps");
1805 RTYPE ("tinit", ir->init_t, 0.0);
1806 RTYPE ("dt", ir->delta_t, 0.001);
1807 STEPTYPE ("nsteps", ir->nsteps, 0);
1808 CTYPE ("For exact run continuation or redoing part of a run");
1809 STEPTYPE ("init-step", ir->init_step, 0);
1810 CTYPE ("Part index is updated automatically on checkpointing (keeps files separate)");
1811 ITYPE ("simulation-part", ir->simulation_part, 1);
1812 CTYPE ("mode for center of mass motion removal");
1813 EETYPE("comm-mode", ir->comm_mode, ecm_names);
1814 CTYPE ("number of steps for center of mass motion removal");
1815 ITYPE ("nstcomm", ir->nstcomm, 100);
1816 CTYPE ("group(s) for center of mass motion removal");
1817 STYPE ("comm-grps", is->vcm, NULL);
1819 CCTYPE ("LANGEVIN DYNAMICS OPTIONS");
1820 CTYPE ("Friction coefficient (amu/ps) and random seed");
1821 RTYPE ("bd-fric", ir->bd_fric, 0.0);
1822 STEPTYPE ("ld-seed", ir->ld_seed, -1);
1825 CCTYPE ("ENERGY MINIMIZATION OPTIONS");
1826 CTYPE ("Force tolerance and initial step-size");
1827 RTYPE ("emtol", ir->em_tol, 10.0);
1828 RTYPE ("emstep", ir->em_stepsize, 0.01);
1829 CTYPE ("Max number of iterations in relax-shells");
1830 ITYPE ("niter", ir->niter, 20);
1831 CTYPE ("Step size (ps^2) for minimization of flexible constraints");
1832 RTYPE ("fcstep", ir->fc_stepsize, 0);
1833 CTYPE ("Frequency of steepest descents steps when doing CG");
1834 ITYPE ("nstcgsteep", ir->nstcgsteep, 1000);
1835 ITYPE ("nbfgscorr", ir->nbfgscorr, 10);
1837 CCTYPE ("TEST PARTICLE INSERTION OPTIONS");
1838 RTYPE ("rtpi", ir->rtpi, 0.05);
1840 /* Output options */
1841 CCTYPE ("OUTPUT CONTROL OPTIONS");
1842 CTYPE ("Output frequency for coords (x), velocities (v) and forces (f)");
1843 ITYPE ("nstxout", ir->nstxout, 0);
1844 ITYPE ("nstvout", ir->nstvout, 0);
1845 ITYPE ("nstfout", ir->nstfout, 0);
1846 ir->nstcheckpoint = 1000;
1847 CTYPE ("Output frequency for energies to log file and energy file");
1848 ITYPE ("nstlog", ir->nstlog, 1000);
1849 ITYPE ("nstcalcenergy", ir->nstcalcenergy, 100);
1850 ITYPE ("nstenergy", ir->nstenergy, 1000);
1851 CTYPE ("Output frequency and precision for .xtc file");
1852 ITYPE ("nstxout-compressed", ir->nstxout_compressed, 0);
1853 RTYPE ("compressed-x-precision", ir->x_compression_precision, 1000.0);
1854 CTYPE ("This selects the subset of atoms for the compressed");
1855 CTYPE ("trajectory file. You can select multiple groups. By");
1856 CTYPE ("default, all atoms will be written.");
1857 STYPE ("compressed-x-grps", is->x_compressed_groups, NULL);
1858 CTYPE ("Selection of energy groups");
1859 STYPE ("energygrps", is->energy, NULL);
1861 /* Neighbor searching */
1862 CCTYPE ("NEIGHBORSEARCHING PARAMETERS");
1863 CTYPE ("cut-off scheme (Verlet: particle based cut-offs, group: using charge groups)");
1864 EETYPE("cutoff-scheme", ir->cutoff_scheme, ecutscheme_names);
1865 CTYPE ("nblist update frequency");
1866 ITYPE ("nstlist", ir->nstlist, 10);
1867 CTYPE ("ns algorithm (simple or grid)");
1868 EETYPE("ns-type", ir->ns_type, ens_names);
1869 /* set ndelta to the optimal value of 2 */
1871 CTYPE ("Periodic boundary conditions: xyz, no, xy");
1872 EETYPE("pbc", ir->ePBC, epbc_names);
1873 EETYPE("periodic-molecules", ir->bPeriodicMols, yesno_names);
1874 CTYPE ("Allowed energy error due to the Verlet buffer in kJ/mol/ps per atom,");
1875 CTYPE ("a value of -1 means: use rlist");
1876 RTYPE("verlet-buffer-tolerance", ir->verletbuf_tol, 0.005);
1877 CTYPE ("nblist cut-off");
1878 RTYPE ("rlist", ir->rlist, 1.0);
1879 CTYPE ("long-range cut-off for switched potentials");
1880 RTYPE ("rlistlong", ir->rlistlong, -1);
1881 ITYPE ("nstcalclr", ir->nstcalclr, -1);
1883 /* Electrostatics */
1884 CCTYPE ("OPTIONS FOR ELECTROSTATICS AND VDW");
1885 CTYPE ("Method for doing electrostatics");
1886 EETYPE("coulombtype", ir->coulombtype, eel_names);
1887 EETYPE("coulomb-modifier", ir->coulomb_modifier, eintmod_names);
1888 CTYPE ("cut-off lengths");
1889 RTYPE ("rcoulomb-switch", ir->rcoulomb_switch, 0.0);
1890 RTYPE ("rcoulomb", ir->rcoulomb, 1.0);
1891 CTYPE ("Relative dielectric constant for the medium and the reaction field");
1892 RTYPE ("epsilon-r", ir->epsilon_r, 1.0);
1893 RTYPE ("epsilon-rf", ir->epsilon_rf, 0.0);
1894 CTYPE ("Method for doing Van der Waals");
1895 EETYPE("vdw-type", ir->vdwtype, evdw_names);
1896 EETYPE("vdw-modifier", ir->vdw_modifier, eintmod_names);
1897 CTYPE ("cut-off lengths");
1898 RTYPE ("rvdw-switch", ir->rvdw_switch, 0.0);
1899 RTYPE ("rvdw", ir->rvdw, 1.0);
1900 CTYPE ("Apply long range dispersion corrections for Energy and Pressure");
1901 EETYPE("DispCorr", ir->eDispCorr, edispc_names);
1902 CTYPE ("Extension of the potential lookup tables beyond the cut-off");
1903 RTYPE ("table-extension", ir->tabext, 1.0);
1904 CTYPE ("Separate tables between energy group pairs");
1905 STYPE ("energygrp-table", is->egptable, NULL);
1906 CTYPE ("Spacing for the PME/PPPM FFT grid");
1907 RTYPE ("fourierspacing", ir->fourier_spacing, 0.12);
1908 CTYPE ("FFT grid size, when a value is 0 fourierspacing will be used");
1909 ITYPE ("fourier-nx", ir->nkx, 0);
1910 ITYPE ("fourier-ny", ir->nky, 0);
1911 ITYPE ("fourier-nz", ir->nkz, 0);
1912 CTYPE ("EWALD/PME/PPPM parameters");
1913 ITYPE ("pme-order", ir->pme_order, 4);
1914 RTYPE ("ewald-rtol", ir->ewald_rtol, 0.00001);
1915 RTYPE ("ewald-rtol-lj", ir->ewald_rtol_lj, 0.001);
1916 EETYPE("lj-pme-comb-rule", ir->ljpme_combination_rule, eljpme_names);
1917 EETYPE("ewald-geometry", ir->ewald_geometry, eewg_names);
1918 RTYPE ("epsilon-surface", ir->epsilon_surface, 0.0);
1919 EETYPE("optimize-fft", ir->bOptFFT, yesno_names);
1921 CCTYPE("IMPLICIT SOLVENT ALGORITHM");
1922 EETYPE("implicit-solvent", ir->implicit_solvent, eis_names);
1924 CCTYPE ("GENERALIZED BORN ELECTROSTATICS");
1925 CTYPE ("Algorithm for calculating Born radii");
1926 EETYPE("gb-algorithm", ir->gb_algorithm, egb_names);
1927 CTYPE ("Frequency of calculating the Born radii inside rlist");
1928 ITYPE ("nstgbradii", ir->nstgbradii, 1);
1929 CTYPE ("Cutoff for Born radii calculation; the contribution from atoms");
1930 CTYPE ("between rlist and rgbradii is updated every nstlist steps");
1931 RTYPE ("rgbradii", ir->rgbradii, 1.0);
1932 CTYPE ("Dielectric coefficient of the implicit solvent");
1933 RTYPE ("gb-epsilon-solvent", ir->gb_epsilon_solvent, 80.0);
1934 CTYPE ("Salt concentration in M for Generalized Born models");
1935 RTYPE ("gb-saltconc", ir->gb_saltconc, 0.0);
1936 CTYPE ("Scaling factors used in the OBC GB model. Default values are OBC(II)");
1937 RTYPE ("gb-obc-alpha", ir->gb_obc_alpha, 1.0);
1938 RTYPE ("gb-obc-beta", ir->gb_obc_beta, 0.8);
1939 RTYPE ("gb-obc-gamma", ir->gb_obc_gamma, 4.85);
1940 RTYPE ("gb-dielectric-offset", ir->gb_dielectric_offset, 0.009);
1941 EETYPE("sa-algorithm", ir->sa_algorithm, esa_names);
1942 CTYPE ("Surface tension (kJ/mol/nm^2) for the SA (nonpolar surface) part of GBSA");
1943 CTYPE ("The value -1 will set default value for Still/HCT/OBC GB-models.");
1944 RTYPE ("sa-surface-tension", ir->sa_surface_tension, -1);
1946 /* Coupling stuff */
1947 CCTYPE ("OPTIONS FOR WEAK COUPLING ALGORITHMS");
1948 CTYPE ("Temperature coupling");
1949 EETYPE("tcoupl", ir->etc, etcoupl_names);
1950 ITYPE ("nsttcouple", ir->nsttcouple, -1);
1951 ITYPE("nh-chain-length", ir->opts.nhchainlength, 10);
1952 EETYPE("print-nose-hoover-chain-variables", ir->bPrintNHChains, yesno_names);
1953 CTYPE ("Groups to couple separately");
1954 STYPE ("tc-grps", is->tcgrps, NULL);
1955 CTYPE ("Time constant (ps) and reference temperature (K)");
1956 STYPE ("tau-t", is->tau_t, NULL);
1957 STYPE ("ref-t", is->ref_t, NULL);
1958 CTYPE ("pressure coupling");
1959 EETYPE("pcoupl", ir->epc, epcoupl_names);
1960 EETYPE("pcoupltype", ir->epct, epcoupltype_names);
1961 ITYPE ("nstpcouple", ir->nstpcouple, -1);
1962 CTYPE ("Time constant (ps), compressibility (1/bar) and reference P (bar)");
1963 RTYPE ("tau-p", ir->tau_p, 1.0);
1964 STYPE ("compressibility", dumstr[0], NULL);
1965 STYPE ("ref-p", dumstr[1], NULL);
1966 CTYPE ("Scaling of reference coordinates, No, All or COM");
1967 EETYPE ("refcoord-scaling", ir->refcoord_scaling, erefscaling_names);
1970 CCTYPE ("OPTIONS FOR QMMM calculations");
1971 EETYPE("QMMM", ir->bQMMM, yesno_names);
1972 CTYPE ("Groups treated Quantum Mechanically");
1973 STYPE ("QMMM-grps", is->QMMM, NULL);
1974 CTYPE ("QM method");
1975 STYPE("QMmethod", is->QMmethod, NULL);
1976 CTYPE ("QMMM scheme");
1977 EETYPE("QMMMscheme", ir->QMMMscheme, eQMMMscheme_names);
1978 CTYPE ("QM basisset");
1979 STYPE("QMbasis", is->QMbasis, NULL);
1980 CTYPE ("QM charge");
1981 STYPE ("QMcharge", is->QMcharge, NULL);
1982 CTYPE ("QM multiplicity");
1983 STYPE ("QMmult", is->QMmult, NULL);
1984 CTYPE ("Surface Hopping");
1985 STYPE ("SH", is->bSH, NULL);
1986 CTYPE ("CAS space options");
1987 STYPE ("CASorbitals", is->CASorbitals, NULL);
1988 STYPE ("CASelectrons", is->CASelectrons, NULL);
1989 STYPE ("SAon", is->SAon, NULL);
1990 STYPE ("SAoff", is->SAoff, NULL);
1991 STYPE ("SAsteps", is->SAsteps, NULL);
1992 CTYPE ("Scale factor for MM charges");
1993 RTYPE ("MMChargeScaleFactor", ir->scalefactor, 1.0);
1994 CTYPE ("Optimization of QM subsystem");
1995 STYPE ("bOPT", is->bOPT, NULL);
1996 STYPE ("bTS", is->bTS, NULL);
1998 /* Simulated annealing */
1999 CCTYPE("SIMULATED ANNEALING");
2000 CTYPE ("Type of annealing for each temperature group (no/single/periodic)");
2001 STYPE ("annealing", is->anneal, NULL);
2002 CTYPE ("Number of time points to use for specifying annealing in each group");
2003 STYPE ("annealing-npoints", is->anneal_npoints, NULL);
2004 CTYPE ("List of times at the annealing points for each group");
2005 STYPE ("annealing-time", is->anneal_time, NULL);
2006 CTYPE ("Temp. at each annealing point, for each group.");
2007 STYPE ("annealing-temp", is->anneal_temp, NULL);
2010 CCTYPE ("GENERATE VELOCITIES FOR STARTUP RUN");
2011 EETYPE("gen-vel", opts->bGenVel, yesno_names);
2012 RTYPE ("gen-temp", opts->tempi, 300.0);
2013 ITYPE ("gen-seed", opts->seed, -1);
2016 CCTYPE ("OPTIONS FOR BONDS");
2017 EETYPE("constraints", opts->nshake, constraints);
2018 CTYPE ("Type of constraint algorithm");
2019 EETYPE("constraint-algorithm", ir->eConstrAlg, econstr_names);
2020 CTYPE ("Do not constrain the start configuration");
2021 EETYPE("continuation", ir->bContinuation, yesno_names);
2022 CTYPE ("Use successive overrelaxation to reduce the number of shake iterations");
2023 EETYPE("Shake-SOR", ir->bShakeSOR, yesno_names);
2024 CTYPE ("Relative tolerance of shake");
2025 RTYPE ("shake-tol", ir->shake_tol, 0.0001);
2026 CTYPE ("Highest order in the expansion of the constraint coupling matrix");
2027 ITYPE ("lincs-order", ir->nProjOrder, 4);
2028 CTYPE ("Number of iterations in the final step of LINCS. 1 is fine for");
2029 CTYPE ("normal simulations, but use 2 to conserve energy in NVE runs.");
2030 CTYPE ("For energy minimization with constraints it should be 4 to 8.");
2031 ITYPE ("lincs-iter", ir->nLincsIter, 1);
2032 CTYPE ("Lincs will write a warning to the stderr if in one step a bond");
2033 CTYPE ("rotates over more degrees than");
2034 RTYPE ("lincs-warnangle", ir->LincsWarnAngle, 30.0);
2035 CTYPE ("Convert harmonic bonds to morse potentials");
2036 EETYPE("morse", opts->bMorse, yesno_names);
2038 /* Energy group exclusions */
2039 CCTYPE ("ENERGY GROUP EXCLUSIONS");
2040 CTYPE ("Pairs of energy groups for which all non-bonded interactions are excluded");
2041 STYPE ("energygrp-excl", is->egpexcl, NULL);
2045 CTYPE ("Number of walls, type, atom types, densities and box-z scale factor for Ewald");
2046 ITYPE ("nwall", ir->nwall, 0);
2047 EETYPE("wall-type", ir->wall_type, ewt_names);
2048 RTYPE ("wall-r-linpot", ir->wall_r_linpot, -1);
2049 STYPE ("wall-atomtype", is->wall_atomtype, NULL);
2050 STYPE ("wall-density", is->wall_density, NULL);
2051 RTYPE ("wall-ewald-zfac", ir->wall_ewald_zfac, 3);
2054 CCTYPE("COM PULLING");
2055 CTYPE("Pull type: no, umbrella, constraint or constant-force");
2056 EETYPE("pull", ir->ePull, epull_names);
2057 if (ir->ePull != epullNO)
2060 is->pull_grp = read_pullparams(&ninp, &inp, ir->pull, &opts->pull_start, wi);
2063 /* Enforced rotation */
2064 CCTYPE("ENFORCED ROTATION");
2065 CTYPE("Enforced rotation: No or Yes");
2066 EETYPE("rotation", ir->bRot, yesno_names);
2070 is->rot_grp = read_rotparams(&ninp, &inp, ir->rot, wi);
2073 /* Interactive MD */
2075 CCTYPE("Group to display and/or manipulate in interactive MD session");
2076 STYPE ("IMD-group", is->imd_grp, NULL);
2077 if (is->imd_grp[0] != '\0')
2084 CCTYPE("NMR refinement stuff");
2085 CTYPE ("Distance restraints type: No, Simple or Ensemble");
2086 EETYPE("disre", ir->eDisre, edisre_names);
2087 CTYPE ("Force weighting of pairs in one distance restraint: Conservative or Equal");
2088 EETYPE("disre-weighting", ir->eDisreWeighting, edisreweighting_names);
2089 CTYPE ("Use sqrt of the time averaged times the instantaneous violation");
2090 EETYPE("disre-mixed", ir->bDisreMixed, yesno_names);
2091 RTYPE ("disre-fc", ir->dr_fc, 1000.0);
2092 RTYPE ("disre-tau", ir->dr_tau, 0.0);
2093 CTYPE ("Output frequency for pair distances to energy file");
2094 ITYPE ("nstdisreout", ir->nstdisreout, 100);
2095 CTYPE ("Orientation restraints: No or Yes");
2096 EETYPE("orire", opts->bOrire, yesno_names);
2097 CTYPE ("Orientation restraints force constant and tau for time averaging");
2098 RTYPE ("orire-fc", ir->orires_fc, 0.0);
2099 RTYPE ("orire-tau", ir->orires_tau, 0.0);
2100 STYPE ("orire-fitgrp", is->orirefitgrp, NULL);
2101 CTYPE ("Output frequency for trace(SD) and S to energy file");
2102 ITYPE ("nstorireout", ir->nstorireout, 100);
2104 /* free energy variables */
2105 CCTYPE ("Free energy variables");
2106 EETYPE("free-energy", ir->efep, efep_names);
2107 STYPE ("couple-moltype", is->couple_moltype, NULL);
2108 EETYPE("couple-lambda0", opts->couple_lam0, couple_lam);
2109 EETYPE("couple-lambda1", opts->couple_lam1, couple_lam);
2110 EETYPE("couple-intramol", opts->bCoupleIntra, yesno_names);
2112 RTYPE ("init-lambda", fep->init_lambda, -1); /* start with -1 so
2114 it was not entered */
2115 ITYPE ("init-lambda-state", fep->init_fep_state, -1);
2116 RTYPE ("delta-lambda", fep->delta_lambda, 0.0);
2117 ITYPE ("nstdhdl", fep->nstdhdl, 50);
2118 STYPE ("fep-lambdas", is->fep_lambda[efptFEP], NULL);
2119 STYPE ("mass-lambdas", is->fep_lambda[efptMASS], NULL);
2120 STYPE ("coul-lambdas", is->fep_lambda[efptCOUL], NULL);
2121 STYPE ("vdw-lambdas", is->fep_lambda[efptVDW], NULL);
2122 STYPE ("bonded-lambdas", is->fep_lambda[efptBONDED], NULL);
2123 STYPE ("restraint-lambdas", is->fep_lambda[efptRESTRAINT], NULL);
2124 STYPE ("temperature-lambdas", is->fep_lambda[efptTEMPERATURE], NULL);
2125 ITYPE ("calc-lambda-neighbors", fep->lambda_neighbors, 1);
2126 STYPE ("init-lambda-weights", is->lambda_weights, NULL);
2127 EETYPE("dhdl-print-energy", fep->bPrintEnergy, yesno_names);
2128 RTYPE ("sc-alpha", fep->sc_alpha, 0.0);
2129 ITYPE ("sc-power", fep->sc_power, 1);
2130 RTYPE ("sc-r-power", fep->sc_r_power, 6.0);
2131 RTYPE ("sc-sigma", fep->sc_sigma, 0.3);
2132 EETYPE("sc-coul", fep->bScCoul, yesno_names);
2133 ITYPE ("dh_hist_size", fep->dh_hist_size, 0);
2134 RTYPE ("dh_hist_spacing", fep->dh_hist_spacing, 0.1);
2135 EETYPE("separate-dhdl-file", fep->separate_dhdl_file,
2136 separate_dhdl_file_names);
2137 EETYPE("dhdl-derivatives", fep->dhdl_derivatives, dhdl_derivatives_names);
2138 ITYPE ("dh_hist_size", fep->dh_hist_size, 0);
2139 RTYPE ("dh_hist_spacing", fep->dh_hist_spacing, 0.1);
2141 /* Non-equilibrium MD stuff */
2142 CCTYPE("Non-equilibrium MD stuff");
2143 STYPE ("acc-grps", is->accgrps, NULL);
2144 STYPE ("accelerate", is->acc, NULL);
2145 STYPE ("freezegrps", is->freeze, NULL);
2146 STYPE ("freezedim", is->frdim, NULL);
2147 RTYPE ("cos-acceleration", ir->cos_accel, 0);
2148 STYPE ("deform", is->deform, NULL);
2150 /* simulated tempering variables */
2151 CCTYPE("simulated tempering variables");
2152 EETYPE("simulated-tempering", ir->bSimTemp, yesno_names);
2153 EETYPE("simulated-tempering-scaling", ir->simtempvals->eSimTempScale, esimtemp_names);
2154 RTYPE("sim-temp-low", ir->simtempvals->simtemp_low, 300.0);
2155 RTYPE("sim-temp-high", ir->simtempvals->simtemp_high, 300.0);
2157 /* expanded ensemble variables */
2158 if (ir->efep == efepEXPANDED || ir->bSimTemp)
2160 read_expandedparams(&ninp, &inp, expand, wi);
2163 /* Electric fields */
2164 CCTYPE("Electric fields");
2165 CTYPE ("Format is number of terms (int) and for all terms an amplitude (real)");
2166 CTYPE ("and a phase angle (real)");
2167 STYPE ("E-x", is->efield_x, NULL);
2168 STYPE ("E-xt", is->efield_xt, NULL);
2169 STYPE ("E-y", is->efield_y, NULL);
2170 STYPE ("E-yt", is->efield_yt, NULL);
2171 STYPE ("E-z", is->efield_z, NULL);
2172 STYPE ("E-zt", is->efield_zt, NULL);
2174 CCTYPE("Ion/water position swapping for computational electrophysiology setups");
2175 CTYPE("Swap positions along direction: no, X, Y, Z");
2176 EETYPE("swapcoords", ir->eSwapCoords, eSwapTypes_names);
2177 if (ir->eSwapCoords != eswapNO)
2180 CTYPE("Swap attempt frequency");
2181 ITYPE("swap-frequency", ir->swap->nstswap, 1);
2182 CTYPE("Two index groups that contain the compartment-partitioning atoms");
2183 STYPE("split-group0", splitgrp0, NULL);
2184 STYPE("split-group1", splitgrp1, NULL);
2185 CTYPE("Use center of mass of split groups (yes/no), otherwise center of geometry is used");
2186 EETYPE("massw-split0", ir->swap->massw_split[0], yesno_names);
2187 EETYPE("massw-split1", ir->swap->massw_split[1], yesno_names);
2189 CTYPE("Group name of ions that can be exchanged with solvent molecules");
2190 STYPE("swap-group", swapgrp, NULL);
2191 CTYPE("Group name of solvent molecules");
2192 STYPE("solvent-group", solgrp, NULL);
2194 CTYPE("Split cylinder: radius, upper and lower extension (nm) (this will define the channels)");
2195 CTYPE("Note that the split cylinder settings do not have an influence on the swapping protocol,");
2196 CTYPE("however, if correctly defined, the ion permeation events are counted per channel");
2197 RTYPE("cyl0-r", ir->swap->cyl0r, 2.0);
2198 RTYPE("cyl0-up", ir->swap->cyl0u, 1.0);
2199 RTYPE("cyl0-down", ir->swap->cyl0l, 1.0);
2200 RTYPE("cyl1-r", ir->swap->cyl1r, 2.0);
2201 RTYPE("cyl1-up", ir->swap->cyl1u, 1.0);
2202 RTYPE("cyl1-down", ir->swap->cyl1l, 1.0);
2204 CTYPE("Average the number of ions per compartment over these many swap attempt steps");
2205 ITYPE("coupl-steps", ir->swap->nAverage, 10);
2206 CTYPE("Requested number of anions and cations for each of the two compartments");
2207 CTYPE("-1 means fix the numbers as found in time step 0");
2208 ITYPE("anionsA", ir->swap->nanions[0], -1);
2209 ITYPE("cationsA", ir->swap->ncations[0], -1);
2210 ITYPE("anionsB", ir->swap->nanions[1], -1);
2211 ITYPE("cationsB", ir->swap->ncations[1], -1);
2212 CTYPE("Start to swap ions if threshold difference to requested count is reached");
2213 RTYPE("threshold", ir->swap->threshold, 1.0);
2216 /* AdResS defined thingies */
2217 CCTYPE ("AdResS parameters");
2218 EETYPE("adress", ir->bAdress, yesno_names);
2221 snew(ir->adress, 1);
2222 read_adressparams(&ninp, &inp, ir->adress, wi);
2225 /* User defined thingies */
2226 CCTYPE ("User defined thingies");
2227 STYPE ("user1-grps", is->user1, NULL);
2228 STYPE ("user2-grps", is->user2, NULL);
2229 ITYPE ("userint1", ir->userint1, 0);
2230 ITYPE ("userint2", ir->userint2, 0);
2231 ITYPE ("userint3", ir->userint3, 0);
2232 ITYPE ("userint4", ir->userint4, 0);
2233 RTYPE ("userreal1", ir->userreal1, 0);
2234 RTYPE ("userreal2", ir->userreal2, 0);
2235 RTYPE ("userreal3", ir->userreal3, 0);
2236 RTYPE ("userreal4", ir->userreal4, 0);
2239 write_inpfile(mdparout, ninp, inp, FALSE, wi);
2240 for (i = 0; (i < ninp); i++)
2243 sfree(inp[i].value);
2247 /* Process options if necessary */
2248 for (m = 0; m < 2; m++)
2250 for (i = 0; i < 2*DIM; i++)
2259 if (sscanf(dumstr[m], "%lf", &(dumdub[m][XX])) != 1)
2261 warning_error(wi, "Pressure coupling not enough values (I need 1)");
2263 dumdub[m][YY] = dumdub[m][ZZ] = dumdub[m][XX];
2265 case epctSEMIISOTROPIC:
2266 case epctSURFACETENSION:
2267 if (sscanf(dumstr[m], "%lf%lf",
2268 &(dumdub[m][XX]), &(dumdub[m][ZZ])) != 2)
2270 warning_error(wi, "Pressure coupling not enough values (I need 2)");
2272 dumdub[m][YY] = dumdub[m][XX];
2274 case epctANISOTROPIC:
2275 if (sscanf(dumstr[m], "%lf%lf%lf%lf%lf%lf",
2276 &(dumdub[m][XX]), &(dumdub[m][YY]), &(dumdub[m][ZZ]),
2277 &(dumdub[m][3]), &(dumdub[m][4]), &(dumdub[m][5])) != 6)
2279 warning_error(wi, "Pressure coupling not enough values (I need 6)");
2283 gmx_fatal(FARGS, "Pressure coupling type %s not implemented yet",
2284 epcoupltype_names[ir->epct]);
2288 clear_mat(ir->ref_p);
2289 clear_mat(ir->compress);
2290 for (i = 0; i < DIM; i++)
2292 ir->ref_p[i][i] = dumdub[1][i];
2293 ir->compress[i][i] = dumdub[0][i];
2295 if (ir->epct == epctANISOTROPIC)
2297 ir->ref_p[XX][YY] = dumdub[1][3];
2298 ir->ref_p[XX][ZZ] = dumdub[1][4];
2299 ir->ref_p[YY][ZZ] = dumdub[1][5];
2300 if (ir->ref_p[XX][YY] != 0 && ir->ref_p[XX][ZZ] != 0 && ir->ref_p[YY][ZZ] != 0)
2302 warning(wi, "All off-diagonal reference pressures are non-zero. Are you sure you want to apply a threefold shear stress?\n");
2304 ir->compress[XX][YY] = dumdub[0][3];
2305 ir->compress[XX][ZZ] = dumdub[0][4];
2306 ir->compress[YY][ZZ] = dumdub[0][5];
2307 for (i = 0; i < DIM; i++)
2309 for (m = 0; m < i; m++)
2311 ir->ref_p[i][m] = ir->ref_p[m][i];
2312 ir->compress[i][m] = ir->compress[m][i];
2317 if (ir->comm_mode == ecmNO)
2322 opts->couple_moltype = NULL;
2323 if (strlen(is->couple_moltype) > 0)
2325 if (ir->efep != efepNO)
2327 opts->couple_moltype = strdup(is->couple_moltype);
2328 if (opts->couple_lam0 == opts->couple_lam1)
2330 warning(wi, "The lambda=0 and lambda=1 states for coupling are identical");
2332 if (ir->eI == eiMD && (opts->couple_lam0 == ecouplamNONE ||
2333 opts->couple_lam1 == ecouplamNONE))
2335 warning(wi, "For proper sampling of the (nearly) decoupled state, stochastic dynamics should be used");
2340 warning(wi, "Can not couple a molecule with free_energy = no");
2343 /* FREE ENERGY AND EXPANDED ENSEMBLE OPTIONS */
2344 if (ir->efep != efepNO)
2346 if (fep->delta_lambda > 0)
2348 ir->efep = efepSLOWGROWTH;
2354 fep->bPrintEnergy = TRUE;
2355 /* always print out the energy to dhdl if we are doing expanded ensemble, since we need the total energy
2356 if the temperature is changing. */
2359 if ((ir->efep != efepNO) || ir->bSimTemp)
2361 ir->bExpanded = FALSE;
2362 if ((ir->efep == efepEXPANDED) || ir->bSimTemp)
2364 ir->bExpanded = TRUE;
2366 do_fep_params(ir, is->fep_lambda, is->lambda_weights);
2367 if (ir->bSimTemp) /* done after fep params */
2369 do_simtemp_params(ir);
2374 ir->fepvals->n_lambda = 0;
2377 /* WALL PARAMETERS */
2379 do_wall_params(ir, is->wall_atomtype, is->wall_density, opts);
2381 /* ORIENTATION RESTRAINT PARAMETERS */
2383 if (opts->bOrire && str_nelem(is->orirefitgrp, MAXPTR, NULL) != 1)
2385 warning_error(wi, "ERROR: Need one orientation restraint fit group\n");
2388 /* DEFORMATION PARAMETERS */
2390 clear_mat(ir->deform);
2391 for (i = 0; i < 6; i++)
2395 m = sscanf(is->deform, "%lf %lf %lf %lf %lf %lf",
2396 &(dumdub[0][0]), &(dumdub[0][1]), &(dumdub[0][2]),
2397 &(dumdub[0][3]), &(dumdub[0][4]), &(dumdub[0][5]));
2398 for (i = 0; i < 3; i++)
2400 ir->deform[i][i] = dumdub[0][i];
2402 ir->deform[YY][XX] = dumdub[0][3];
2403 ir->deform[ZZ][XX] = dumdub[0][4];
2404 ir->deform[ZZ][YY] = dumdub[0][5];
2405 if (ir->epc != epcNO)
2407 for (i = 0; i < 3; i++)
2409 for (j = 0; j <= i; j++)
2411 if (ir->deform[i][j] != 0 && ir->compress[i][j] != 0)
2413 warning_error(wi, "A box element has deform set and compressibility > 0");
2417 for (i = 0; i < 3; i++)
2419 for (j = 0; j < i; j++)
2421 if (ir->deform[i][j] != 0)
2423 for (m = j; m < DIM; m++)
2425 if (ir->compress[m][j] != 0)
2427 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.");
2428 warning(wi, warn_buf);
2436 /* Ion/water position swapping checks */
2437 if (ir->eSwapCoords != eswapNO)
2439 if (ir->swap->nstswap < 1)
2441 warning_error(wi, "swap_frequency must be 1 or larger when ion swapping is requested");
2443 if (ir->swap->nAverage < 1)
2445 warning_error(wi, "coupl_steps must be 1 or larger.\n");
2447 if (ir->swap->threshold < 1.0)
2449 warning_error(wi, "Ion count threshold must be at least 1.\n");
2457 static int search_QMstring(const char *s, int ng, const char *gn[])
2459 /* same as normal search_string, but this one searches QM strings */
2462 for (i = 0; (i < ng); i++)
2464 if (gmx_strcasecmp(s, gn[i]) == 0)
2470 gmx_fatal(FARGS, "this QM method or basisset (%s) is not implemented\n!", s);
2474 } /* search_QMstring */
2476 /* We would like gn to be const as well, but C doesn't allow this */
2477 int search_string(const char *s, int ng, char *gn[])
2481 for (i = 0; (i < ng); i++)
2483 if (gmx_strcasecmp(s, gn[i]) == 0)
2490 "Group %s referenced in the .mdp file was not found in the index file.\n"
2491 "Group names must match either [moleculetype] names or custom index group\n"
2492 "names, in which case you must supply an index file to the '-n' option\n"
2499 static gmx_bool do_numbering(int natoms, gmx_groups_t *groups, int ng, char *ptrs[],
2500 t_blocka *block, char *gnames[],
2501 int gtype, int restnm,
2502 int grptp, gmx_bool bVerbose,
2505 unsigned short *cbuf;
2506 t_grps *grps = &(groups->grps[gtype]);
2507 int i, j, gid, aj, ognr, ntot = 0;
2510 char warn_buf[STRLEN];
2514 fprintf(debug, "Starting numbering %d groups of type %d\n", ng, gtype);
2517 title = gtypes[gtype];
2520 /* Mark all id's as not set */
2521 for (i = 0; (i < natoms); i++)
2526 snew(grps->nm_ind, ng+1); /* +1 for possible rest group */
2527 for (i = 0; (i < ng); i++)
2529 /* Lookup the group name in the block structure */
2530 gid = search_string(ptrs[i], block->nr, gnames);
2531 if ((grptp != egrptpONE) || (i == 0))
2533 grps->nm_ind[grps->nr++] = gid;
2537 fprintf(debug, "Found gid %d for group %s\n", gid, ptrs[i]);
2540 /* Now go over the atoms in the group */
2541 for (j = block->index[gid]; (j < block->index[gid+1]); j++)
2546 /* Range checking */
2547 if ((aj < 0) || (aj >= natoms))
2549 gmx_fatal(FARGS, "Invalid atom number %d in indexfile", aj);
2551 /* Lookup up the old group number */
2555 gmx_fatal(FARGS, "Atom %d in multiple %s groups (%d and %d)",
2556 aj+1, title, ognr+1, i+1);
2560 /* Store the group number in buffer */
2561 if (grptp == egrptpONE)
2574 /* Now check whether we have done all atoms */
2578 if (grptp == egrptpALL)
2580 gmx_fatal(FARGS, "%d atoms are not part of any of the %s groups",
2581 natoms-ntot, title);
2583 else if (grptp == egrptpPART)
2585 sprintf(warn_buf, "%d atoms are not part of any of the %s groups",
2586 natoms-ntot, title);
2587 warning_note(wi, warn_buf);
2589 /* Assign all atoms currently unassigned to a rest group */
2590 for (j = 0; (j < natoms); j++)
2592 if (cbuf[j] == NOGID)
2598 if (grptp != egrptpPART)
2603 "Making dummy/rest group for %s containing %d elements\n",
2604 title, natoms-ntot);
2606 /* Add group name "rest" */
2607 grps->nm_ind[grps->nr] = restnm;
2609 /* Assign the rest name to all atoms not currently assigned to a group */
2610 for (j = 0; (j < natoms); j++)
2612 if (cbuf[j] == NOGID)
2621 if (grps->nr == 1 && (ntot == 0 || ntot == natoms))
2623 /* All atoms are part of one (or no) group, no index required */
2624 groups->ngrpnr[gtype] = 0;
2625 groups->grpnr[gtype] = NULL;
2629 groups->ngrpnr[gtype] = natoms;
2630 snew(groups->grpnr[gtype], natoms);
2631 for (j = 0; (j < natoms); j++)
2633 groups->grpnr[gtype][j] = cbuf[j];
2639 return (bRest && grptp == egrptpPART);
2642 static void calc_nrdf(gmx_mtop_t *mtop, t_inputrec *ir, char **gnames)
2645 gmx_groups_t *groups;
2647 int natoms, ai, aj, i, j, d, g, imin, jmin;
2649 int *nrdf2, *na_vcm, na_tot;
2650 double *nrdf_tc, *nrdf_vcm, nrdf_uc, n_sub = 0;
2651 gmx_mtop_atomloop_all_t aloop;
2653 int mb, mol, ftype, as;
2654 gmx_molblock_t *molb;
2655 gmx_moltype_t *molt;
2658 * First calc 3xnr-atoms for each group
2659 * then subtract half a degree of freedom for each constraint
2661 * Only atoms and nuclei contribute to the degrees of freedom...
2666 groups = &mtop->groups;
2667 natoms = mtop->natoms;
2669 /* Allocate one more for a possible rest group */
2670 /* We need to sum degrees of freedom into doubles,
2671 * since floats give too low nrdf's above 3 million atoms.
2673 snew(nrdf_tc, groups->grps[egcTC].nr+1);
2674 snew(nrdf_vcm, groups->grps[egcVCM].nr+1);
2675 snew(na_vcm, groups->grps[egcVCM].nr+1);
2677 for (i = 0; i < groups->grps[egcTC].nr; i++)
2681 for (i = 0; i < groups->grps[egcVCM].nr+1; i++)
2686 snew(nrdf2, natoms);
2687 aloop = gmx_mtop_atomloop_all_init(mtop);
2688 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
2691 if (atom->ptype == eptAtom || atom->ptype == eptNucleus)
2693 g = ggrpnr(groups, egcFREEZE, i);
2694 /* Double count nrdf for particle i */
2695 for (d = 0; d < DIM; d++)
2697 if (opts->nFreeze[g][d] == 0)
2702 nrdf_tc [ggrpnr(groups, egcTC, i)] += 0.5*nrdf2[i];
2703 nrdf_vcm[ggrpnr(groups, egcVCM, i)] += 0.5*nrdf2[i];
2708 for (mb = 0; mb < mtop->nmolblock; mb++)
2710 molb = &mtop->molblock[mb];
2711 molt = &mtop->moltype[molb->type];
2712 atom = molt->atoms.atom;
2713 for (mol = 0; mol < molb->nmol; mol++)
2715 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
2717 ia = molt->ilist[ftype].iatoms;
2718 for (i = 0; i < molt->ilist[ftype].nr; )
2720 /* Subtract degrees of freedom for the constraints,
2721 * if the particles still have degrees of freedom left.
2722 * If one of the particles is a vsite or a shell, then all
2723 * constraint motion will go there, but since they do not
2724 * contribute to the constraints the degrees of freedom do not
2729 if (((atom[ia[1]].ptype == eptNucleus) ||
2730 (atom[ia[1]].ptype == eptAtom)) &&
2731 ((atom[ia[2]].ptype == eptNucleus) ||
2732 (atom[ia[2]].ptype == eptAtom)))
2750 imin = min(imin, nrdf2[ai]);
2751 jmin = min(jmin, nrdf2[aj]);
2754 nrdf_tc [ggrpnr(groups, egcTC, ai)] -= 0.5*imin;
2755 nrdf_tc [ggrpnr(groups, egcTC, aj)] -= 0.5*jmin;
2756 nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2757 nrdf_vcm[ggrpnr(groups, egcVCM, aj)] -= 0.5*jmin;
2759 ia += interaction_function[ftype].nratoms+1;
2760 i += interaction_function[ftype].nratoms+1;
2763 ia = molt->ilist[F_SETTLE].iatoms;
2764 for (i = 0; i < molt->ilist[F_SETTLE].nr; )
2766 /* Subtract 1 dof from every atom in the SETTLE */
2767 for (j = 0; j < 3; j++)
2770 imin = min(2, nrdf2[ai]);
2772 nrdf_tc [ggrpnr(groups, egcTC, ai)] -= 0.5*imin;
2773 nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2778 as += molt->atoms.nr;
2782 if (ir->ePull == epullCONSTRAINT)
2784 /* Correct nrdf for the COM constraints.
2785 * We correct using the TC and VCM group of the first atom
2786 * in the reference and pull group. If atoms in one pull group
2787 * belong to different TC or VCM groups it is anyhow difficult
2788 * to determine the optimal nrdf assignment.
2792 for (i = 0; i < pull->ncoord; i++)
2796 for (j = 0; j < 2; j++)
2798 const t_pull_group *pgrp;
2800 pgrp = &pull->group[pull->coord[i].group[j]];
2804 /* Subtract 1/2 dof from each group */
2806 nrdf_tc [ggrpnr(groups, egcTC, ai)] -= 0.5*imin;
2807 nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2808 if (nrdf_tc[ggrpnr(groups, egcTC, ai)] < 0)
2810 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)]]);
2815 /* We need to subtract the whole DOF from group j=1 */
2822 if (ir->nstcomm != 0)
2824 /* Subtract 3 from the number of degrees of freedom in each vcm group
2825 * when com translation is removed and 6 when rotation is removed
2828 switch (ir->comm_mode)
2831 n_sub = ndof_com(ir);
2838 gmx_incons("Checking comm_mode");
2841 for (i = 0; i < groups->grps[egcTC].nr; i++)
2843 /* Count the number of atoms of TC group i for every VCM group */
2844 for (j = 0; j < groups->grps[egcVCM].nr+1; j++)
2849 for (ai = 0; ai < natoms; ai++)
2851 if (ggrpnr(groups, egcTC, ai) == i)
2853 na_vcm[ggrpnr(groups, egcVCM, ai)]++;
2857 /* Correct for VCM removal according to the fraction of each VCM
2858 * group present in this TC group.
2860 nrdf_uc = nrdf_tc[i];
2863 fprintf(debug, "T-group[%d] nrdf_uc = %g, n_sub = %g\n",
2867 for (j = 0; j < groups->grps[egcVCM].nr+1; j++)
2869 if (nrdf_vcm[j] > n_sub)
2871 nrdf_tc[i] += nrdf_uc*((double)na_vcm[j]/(double)na_tot)*
2872 (nrdf_vcm[j] - n_sub)/nrdf_vcm[j];
2876 fprintf(debug, " nrdf_vcm[%d] = %g, nrdf = %g\n",
2877 j, nrdf_vcm[j], nrdf_tc[i]);
2882 for (i = 0; (i < groups->grps[egcTC].nr); i++)
2884 opts->nrdf[i] = nrdf_tc[i];
2885 if (opts->nrdf[i] < 0)
2890 "Number of degrees of freedom in T-Coupling group %s is %.2f\n",
2891 gnames[groups->grps[egcTC].nm_ind[i]], opts->nrdf[i]);
2900 static void decode_cos(char *s, t_cosines *cosine)
2903 char format[STRLEN], f1[STRLEN];
2915 sscanf(t, "%d", &(cosine->n));
2922 snew(cosine->a, cosine->n);
2923 snew(cosine->phi, cosine->n);
2925 sprintf(format, "%%*d");
2926 for (i = 0; (i < cosine->n); i++)
2929 strcat(f1, "%lf%lf");
2930 if (sscanf(t, f1, &a, &phi) < 2)
2932 gmx_fatal(FARGS, "Invalid input for electric field shift: '%s'", t);
2935 cosine->phi[i] = phi;
2936 strcat(format, "%*lf%*lf");
2943 static gmx_bool do_egp_flag(t_inputrec *ir, gmx_groups_t *groups,
2944 const char *option, const char *val, int flag)
2946 /* The maximum number of energy group pairs would be MAXPTR*(MAXPTR+1)/2.
2947 * But since this is much larger than STRLEN, such a line can not be parsed.
2948 * The real maximum is the number of names that fit in a string: STRLEN/2.
2950 #define EGP_MAX (STRLEN/2)
2951 int nelem, i, j, k, nr;
2952 char *names[EGP_MAX];
2956 gnames = groups->grpname;
2958 nelem = str_nelem(val, EGP_MAX, names);
2961 gmx_fatal(FARGS, "The number of groups for %s is odd", option);
2963 nr = groups->grps[egcENER].nr;
2965 for (i = 0; i < nelem/2; i++)
2969 gmx_strcasecmp(names[2*i], *(gnames[groups->grps[egcENER].nm_ind[j]])))
2975 gmx_fatal(FARGS, "%s in %s is not an energy group\n",
2976 names[2*i], option);
2980 gmx_strcasecmp(names[2*i+1], *(gnames[groups->grps[egcENER].nm_ind[k]])))
2986 gmx_fatal(FARGS, "%s in %s is not an energy group\n",
2987 names[2*i+1], option);
2989 if ((j < nr) && (k < nr))
2991 ir->opts.egp_flags[nr*j+k] |= flag;
2992 ir->opts.egp_flags[nr*k+j] |= flag;
3001 static void make_swap_groups(
3010 int ig = -1, i = 0, j;
3014 /* Just a quick check here, more thorough checks are in mdrun */
3015 if (strcmp(splitg0name, splitg1name) == 0)
3017 gmx_fatal(FARGS, "The split groups can not both be '%s'.", splitg0name);
3020 /* First get the swap group index atoms */
3021 ig = search_string(swapgname, grps->nr, gnames);
3022 swap->nat = grps->index[ig+1] - grps->index[ig];
3025 fprintf(stderr, "Swap group '%s' contains %d atoms.\n", swapgname, swap->nat);
3026 snew(swap->ind, swap->nat);
3027 for (i = 0; i < swap->nat; i++)
3029 swap->ind[i] = grps->a[grps->index[ig]+i];
3034 gmx_fatal(FARGS, "You defined an empty group of atoms for swapping.");
3037 /* Now do so for the split groups */
3038 for (j = 0; j < 2; j++)
3042 splitg = splitg0name;
3046 splitg = splitg1name;
3049 ig = search_string(splitg, grps->nr, gnames);
3050 swap->nat_split[j] = grps->index[ig+1] - grps->index[ig];
3051 if (swap->nat_split[j] > 0)
3053 fprintf(stderr, "Split group %d '%s' contains %d atom%s.\n",
3054 j, splitg, swap->nat_split[j], (swap->nat_split[j] > 1) ? "s" : "");
3055 snew(swap->ind_split[j], swap->nat_split[j]);
3056 for (i = 0; i < swap->nat_split[j]; i++)
3058 swap->ind_split[j][i] = grps->a[grps->index[ig]+i];
3063 gmx_fatal(FARGS, "Split group %d has to contain at least 1 atom!", j);
3067 /* Now get the solvent group index atoms */
3068 ig = search_string(solgname, grps->nr, gnames);
3069 swap->nat_sol = grps->index[ig+1] - grps->index[ig];
3070 if (swap->nat_sol > 0)
3072 fprintf(stderr, "Solvent group '%s' contains %d atoms.\n", solgname, swap->nat_sol);
3073 snew(swap->ind_sol, swap->nat_sol);
3074 for (i = 0; i < swap->nat_sol; i++)
3076 swap->ind_sol[i] = grps->a[grps->index[ig]+i];
3081 gmx_fatal(FARGS, "You defined an empty group of solvent. Cannot exchange ions.");
3086 void make_IMD_group(t_IMD *IMDgroup, char *IMDgname, t_blocka *grps, char **gnames)
3091 ig = search_string(IMDgname, grps->nr, gnames);
3092 IMDgroup->nat = grps->index[ig+1] - grps->index[ig];
3094 if (IMDgroup->nat > 0)
3096 fprintf(stderr, "Group '%s' with %d atoms can be activated for interactive molecular dynamics (IMD).\n",
3097 IMDgname, IMDgroup->nat);
3098 snew(IMDgroup->ind, IMDgroup->nat);
3099 for (i = 0; i < IMDgroup->nat; i++)
3101 IMDgroup->ind[i] = grps->a[grps->index[ig]+i];
3107 void do_index(const char* mdparin, const char *ndx,
3110 t_inputrec *ir, rvec *v,
3114 gmx_groups_t *groups;
3118 char warnbuf[STRLEN], **gnames;
3119 int nr, ntcg, ntau_t, nref_t, nacc, nofg, nSA, nSA_points, nSA_time, nSA_temp;
3122 int nacg, nfreeze, nfrdim, nenergy, nvcm, nuser;
3123 char *ptr1[MAXPTR], *ptr2[MAXPTR], *ptr3[MAXPTR];
3124 int i, j, k, restnm;
3126 gmx_bool bExcl, bTable, bSetTCpar, bAnneal, bRest;
3127 int nQMmethod, nQMbasis, nQMcharge, nQMmult, nbSH, nCASorb, nCASelec,
3128 nSAon, nSAoff, nSAsteps, nQMg, nbOPT, nbTS;
3129 char warn_buf[STRLEN];
3133 fprintf(stderr, "processing index file...\n");
3139 snew(grps->index, 1);
3141 atoms_all = gmx_mtop_global_atoms(mtop);
3142 analyse(&atoms_all, grps, &gnames, FALSE, TRUE);
3143 free_t_atoms(&atoms_all, FALSE);
3147 grps = init_index(ndx, &gnames);
3150 groups = &mtop->groups;
3151 natoms = mtop->natoms;
3152 symtab = &mtop->symtab;
3154 snew(groups->grpname, grps->nr+1);
3156 for (i = 0; (i < grps->nr); i++)
3158 groups->grpname[i] = put_symtab(symtab, gnames[i]);
3160 groups->grpname[i] = put_symtab(symtab, "rest");
3162 srenew(gnames, grps->nr+1);
3163 gnames[restnm] = *(groups->grpname[i]);
3164 groups->ngrpname = grps->nr+1;
3166 set_warning_line(wi, mdparin, -1);
3168 ntau_t = str_nelem(is->tau_t, MAXPTR, ptr1);
3169 nref_t = str_nelem(is->ref_t, MAXPTR, ptr2);
3170 ntcg = str_nelem(is->tcgrps, MAXPTR, ptr3);
3171 if ((ntau_t != ntcg) || (nref_t != ntcg))
3173 gmx_fatal(FARGS, "Invalid T coupling input: %d groups, %d ref-t values and "
3174 "%d tau-t values", ntcg, nref_t, ntau_t);
3177 bSetTCpar = (ir->etc || EI_SD(ir->eI) || ir->eI == eiBD || EI_TPI(ir->eI));
3178 do_numbering(natoms, groups, ntcg, ptr3, grps, gnames, egcTC,
3179 restnm, bSetTCpar ? egrptpALL : egrptpALL_GENREST, bVerbose, wi);
3180 nr = groups->grps[egcTC].nr;
3182 snew(ir->opts.nrdf, nr);
3183 snew(ir->opts.tau_t, nr);
3184 snew(ir->opts.ref_t, nr);
3185 if (ir->eI == eiBD && ir->bd_fric == 0)
3187 fprintf(stderr, "bd-fric=0, so tau-t will be used as the inverse friction constant(s)\n");
3194 gmx_fatal(FARGS, "Not enough ref-t and tau-t values!");
3198 for (i = 0; (i < nr); i++)
3200 ir->opts.tau_t[i] = strtod(ptr1[i], NULL);
3201 if ((ir->eI == eiBD || ir->eI == eiSD2) && ir->opts.tau_t[i] <= 0)
3203 sprintf(warn_buf, "With integrator %s tau-t should be larger than 0", ei_names[ir->eI]);
3204 warning_error(wi, warn_buf);
3207 if (ir->etc != etcVRESCALE && ir->opts.tau_t[i] == 0)
3209 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.");
3212 if (ir->opts.tau_t[i] >= 0)
3214 tau_min = min(tau_min, ir->opts.tau_t[i]);
3217 if (ir->etc != etcNO && ir->nsttcouple == -1)
3219 ir->nsttcouple = ir_optimal_nsttcouple(ir);
3224 if ((ir->etc == etcNOSEHOOVER) && (ir->epc == epcBERENDSEN))
3226 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");
3228 if ((ir->epc == epcMTTK) && (ir->etc > etcNO))
3230 if (ir->nstpcouple != ir->nsttcouple)
3232 int mincouple = min(ir->nstpcouple, ir->nsttcouple);
3233 ir->nstpcouple = ir->nsttcouple = mincouple;
3234 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);
3235 warning_note(wi, warn_buf);
3239 /* velocity verlet with averaged kinetic energy KE = 0.5*(v(t+1/2) - v(t-1/2)) is implemented
3240 primarily for testing purposes, and does not work with temperature coupling other than 1 */
3242 if (ETC_ANDERSEN(ir->etc))
3244 if (ir->nsttcouple != 1)
3247 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");
3248 warning_note(wi, warn_buf);
3251 nstcmin = tcouple_min_integration_steps(ir->etc);
3254 if (tau_min/(ir->delta_t*ir->nsttcouple) < nstcmin)
3256 sprintf(warn_buf, "For proper integration of the %s thermostat, tau-t (%g) should be at least %d times larger than nsttcouple*dt (%g)",
3257 ETCOUPLTYPE(ir->etc),
3259 ir->nsttcouple*ir->delta_t);
3260 warning(wi, warn_buf);
3263 for (i = 0; (i < nr); i++)
3265 ir->opts.ref_t[i] = strtod(ptr2[i], NULL);
3266 if (ir->opts.ref_t[i] < 0)
3268 gmx_fatal(FARGS, "ref-t for group %d negative", i);
3271 /* set the lambda mc temperature to the md integrator temperature (which should be defined
3272 if we are in this conditional) if mc_temp is negative */
3273 if (ir->expandedvals->mc_temp < 0)
3275 ir->expandedvals->mc_temp = ir->opts.ref_t[0]; /*for now, set to the first reft */
3279 /* Simulated annealing for each group. There are nr groups */
3280 nSA = str_nelem(is->anneal, MAXPTR, ptr1);
3281 if (nSA == 1 && (ptr1[0][0] == 'n' || ptr1[0][0] == 'N'))
3285 if (nSA > 0 && nSA != nr)
3287 gmx_fatal(FARGS, "Not enough annealing values: %d (for %d groups)\n", nSA, nr);
3291 snew(ir->opts.annealing, nr);
3292 snew(ir->opts.anneal_npoints, nr);
3293 snew(ir->opts.anneal_time, nr);
3294 snew(ir->opts.anneal_temp, nr);
3295 for (i = 0; i < nr; i++)
3297 ir->opts.annealing[i] = eannNO;
3298 ir->opts.anneal_npoints[i] = 0;
3299 ir->opts.anneal_time[i] = NULL;
3300 ir->opts.anneal_temp[i] = NULL;
3305 for (i = 0; i < nr; i++)
3307 if (ptr1[i][0] == 'n' || ptr1[i][0] == 'N')
3309 ir->opts.annealing[i] = eannNO;
3311 else if (ptr1[i][0] == 's' || ptr1[i][0] == 'S')
3313 ir->opts.annealing[i] = eannSINGLE;
3316 else if (ptr1[i][0] == 'p' || ptr1[i][0] == 'P')
3318 ir->opts.annealing[i] = eannPERIODIC;
3324 /* Read the other fields too */
3325 nSA_points = str_nelem(is->anneal_npoints, MAXPTR, ptr1);
3326 if (nSA_points != nSA)
3328 gmx_fatal(FARGS, "Found %d annealing-npoints values for %d groups\n", nSA_points, nSA);
3330 for (k = 0, i = 0; i < nr; i++)
3332 ir->opts.anneal_npoints[i] = strtol(ptr1[i], NULL, 10);
3333 if (ir->opts.anneal_npoints[i] == 1)
3335 gmx_fatal(FARGS, "Please specify at least a start and an end point for annealing\n");
3337 snew(ir->opts.anneal_time[i], ir->opts.anneal_npoints[i]);
3338 snew(ir->opts.anneal_temp[i], ir->opts.anneal_npoints[i]);
3339 k += ir->opts.anneal_npoints[i];
3342 nSA_time = str_nelem(is->anneal_time, MAXPTR, ptr1);
3345 gmx_fatal(FARGS, "Found %d annealing-time values, wanter %d\n", nSA_time, k);
3347 nSA_temp = str_nelem(is->anneal_temp, MAXPTR, ptr2);
3350 gmx_fatal(FARGS, "Found %d annealing-temp values, wanted %d\n", nSA_temp, k);
3353 for (i = 0, k = 0; i < nr; i++)
3356 for (j = 0; j < ir->opts.anneal_npoints[i]; j++)
3358 ir->opts.anneal_time[i][j] = strtod(ptr1[k], NULL);
3359 ir->opts.anneal_temp[i][j] = strtod(ptr2[k], NULL);
3362 if (ir->opts.anneal_time[i][0] > (ir->init_t+GMX_REAL_EPS))
3364 gmx_fatal(FARGS, "First time point for annealing > init_t.\n");
3370 if (ir->opts.anneal_time[i][j] < ir->opts.anneal_time[i][j-1])
3372 gmx_fatal(FARGS, "Annealing timepoints out of order: t=%f comes after t=%f\n",
3373 ir->opts.anneal_time[i][j], ir->opts.anneal_time[i][j-1]);
3376 if (ir->opts.anneal_temp[i][j] < 0)
3378 gmx_fatal(FARGS, "Found negative temperature in annealing: %f\n", ir->opts.anneal_temp[i][j]);
3383 /* Print out some summary information, to make sure we got it right */
3384 for (i = 0, k = 0; i < nr; i++)
3386 if (ir->opts.annealing[i] != eannNO)
3388 j = groups->grps[egcTC].nm_ind[i];
3389 fprintf(stderr, "Simulated annealing for group %s: %s, %d timepoints\n",
3390 *(groups->grpname[j]), eann_names[ir->opts.annealing[i]],
3391 ir->opts.anneal_npoints[i]);
3392 fprintf(stderr, "Time (ps) Temperature (K)\n");
3393 /* All terms except the last one */
3394 for (j = 0; j < (ir->opts.anneal_npoints[i]-1); j++)
3396 fprintf(stderr, "%9.1f %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3399 /* Finally the last one */
3400 j = ir->opts.anneal_npoints[i]-1;
3401 if (ir->opts.annealing[i] == eannSINGLE)
3403 fprintf(stderr, "%9.1f- %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3407 fprintf(stderr, "%9.1f %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3408 if (fabs(ir->opts.anneal_temp[i][j]-ir->opts.anneal_temp[i][0]) > GMX_REAL_EPS)
3410 warning_note(wi, "There is a temperature jump when your annealing loops back.\n");
3419 if (ir->ePull != epullNO)
3421 make_pull_groups(ir->pull, is->pull_grp, grps, gnames);
3423 make_pull_coords(ir->pull);
3428 make_rotation_groups(ir->rot, is->rot_grp, grps, gnames);
3431 if (ir->eSwapCoords != eswapNO)
3433 make_swap_groups(ir->swap, swapgrp, splitgrp0, splitgrp1, solgrp, grps, gnames);
3436 /* Make indices for IMD session */
3439 make_IMD_group(ir->imd, is->imd_grp, grps, gnames);
3442 nacc = str_nelem(is->acc, MAXPTR, ptr1);
3443 nacg = str_nelem(is->accgrps, MAXPTR, ptr2);
3444 if (nacg*DIM != nacc)
3446 gmx_fatal(FARGS, "Invalid Acceleration input: %d groups and %d acc. values",
3449 do_numbering(natoms, groups, nacg, ptr2, grps, gnames, egcACC,
3450 restnm, egrptpALL_GENREST, bVerbose, wi);
3451 nr = groups->grps[egcACC].nr;
3452 snew(ir->opts.acc, nr);
3453 ir->opts.ngacc = nr;
3455 for (i = k = 0; (i < nacg); i++)
3457 for (j = 0; (j < DIM); j++, k++)
3459 ir->opts.acc[i][j] = strtod(ptr1[k], NULL);
3462 for (; (i < nr); i++)
3464 for (j = 0; (j < DIM); j++)
3466 ir->opts.acc[i][j] = 0;
3470 nfrdim = str_nelem(is->frdim, MAXPTR, ptr1);
3471 nfreeze = str_nelem(is->freeze, MAXPTR, ptr2);
3472 if (nfrdim != DIM*nfreeze)
3474 gmx_fatal(FARGS, "Invalid Freezing input: %d groups and %d freeze values",
3477 do_numbering(natoms, groups, nfreeze, ptr2, grps, gnames, egcFREEZE,
3478 restnm, egrptpALL_GENREST, bVerbose, wi);
3479 nr = groups->grps[egcFREEZE].nr;
3480 ir->opts.ngfrz = nr;
3481 snew(ir->opts.nFreeze, nr);
3482 for (i = k = 0; (i < nfreeze); i++)
3484 for (j = 0; (j < DIM); j++, k++)
3486 ir->opts.nFreeze[i][j] = (gmx_strncasecmp(ptr1[k], "Y", 1) == 0);
3487 if (!ir->opts.nFreeze[i][j])
3489 if (gmx_strncasecmp(ptr1[k], "N", 1) != 0)
3491 sprintf(warnbuf, "Please use Y(ES) or N(O) for freezedim only "
3492 "(not %s)", ptr1[k]);
3493 warning(wi, warn_buf);
3498 for (; (i < nr); i++)
3500 for (j = 0; (j < DIM); j++)
3502 ir->opts.nFreeze[i][j] = 0;
3506 nenergy = str_nelem(is->energy, MAXPTR, ptr1);
3507 do_numbering(natoms, groups, nenergy, ptr1, grps, gnames, egcENER,
3508 restnm, egrptpALL_GENREST, bVerbose, wi);
3509 add_wall_energrps(groups, ir->nwall, symtab);
3510 ir->opts.ngener = groups->grps[egcENER].nr;
3511 nvcm = str_nelem(is->vcm, MAXPTR, ptr1);
3513 do_numbering(natoms, groups, nvcm, ptr1, grps, gnames, egcVCM,
3514 restnm, nvcm == 0 ? egrptpALL_GENREST : egrptpPART, bVerbose, wi);
3517 warning(wi, "Some atoms are not part of any center of mass motion removal group.\n"
3518 "This may lead to artifacts.\n"
3519 "In most cases one should use one group for the whole system.");
3522 /* Now we have filled the freeze struct, so we can calculate NRDF */
3523 calc_nrdf(mtop, ir, gnames);
3529 /* Must check per group! */
3530 for (i = 0; (i < ir->opts.ngtc); i++)
3532 ntot += ir->opts.nrdf[i];
3534 if (ntot != (DIM*natoms))
3536 fac = sqrt(ntot/(DIM*natoms));
3539 fprintf(stderr, "Scaling velocities by a factor of %.3f to account for constraints\n"
3540 "and removal of center of mass motion\n", fac);
3542 for (i = 0; (i < natoms); i++)
3544 svmul(fac, v[i], v[i]);
3549 nuser = str_nelem(is->user1, MAXPTR, ptr1);
3550 do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcUser1,
3551 restnm, egrptpALL_GENREST, bVerbose, wi);
3552 nuser = str_nelem(is->user2, MAXPTR, ptr1);
3553 do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcUser2,
3554 restnm, egrptpALL_GENREST, bVerbose, wi);
3555 nuser = str_nelem(is->x_compressed_groups, MAXPTR, ptr1);
3556 do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcCompressedX,
3557 restnm, egrptpONE, bVerbose, wi);
3558 nofg = str_nelem(is->orirefitgrp, MAXPTR, ptr1);
3559 do_numbering(natoms, groups, nofg, ptr1, grps, gnames, egcORFIT,
3560 restnm, egrptpALL_GENREST, bVerbose, wi);
3562 /* QMMM input processing */
3563 nQMg = str_nelem(is->QMMM, MAXPTR, ptr1);
3564 nQMmethod = str_nelem(is->QMmethod, MAXPTR, ptr2);
3565 nQMbasis = str_nelem(is->QMbasis, MAXPTR, ptr3);
3566 if ((nQMmethod != nQMg) || (nQMbasis != nQMg))
3568 gmx_fatal(FARGS, "Invalid QMMM input: %d groups %d basissets"
3569 " and %d methods\n", nQMg, nQMbasis, nQMmethod);
3571 /* group rest, if any, is always MM! */
3572 do_numbering(natoms, groups, nQMg, ptr1, grps, gnames, egcQMMM,
3573 restnm, egrptpALL_GENREST, bVerbose, wi);
3574 nr = nQMg; /*atoms->grps[egcQMMM].nr;*/
3575 ir->opts.ngQM = nQMg;
3576 snew(ir->opts.QMmethod, nr);
3577 snew(ir->opts.QMbasis, nr);
3578 for (i = 0; i < nr; i++)
3580 /* input consists of strings: RHF CASSCF PM3 .. These need to be
3581 * converted to the corresponding enum in names.c
3583 ir->opts.QMmethod[i] = search_QMstring(ptr2[i], eQMmethodNR,
3585 ir->opts.QMbasis[i] = search_QMstring(ptr3[i], eQMbasisNR,
3589 nQMmult = str_nelem(is->QMmult, MAXPTR, ptr1);
3590 nQMcharge = str_nelem(is->QMcharge, MAXPTR, ptr2);
3591 nbSH = str_nelem(is->bSH, MAXPTR, ptr3);
3592 snew(ir->opts.QMmult, nr);
3593 snew(ir->opts.QMcharge, nr);
3594 snew(ir->opts.bSH, nr);
3596 for (i = 0; i < nr; i++)
3598 ir->opts.QMmult[i] = strtol(ptr1[i], NULL, 10);
3599 ir->opts.QMcharge[i] = strtol(ptr2[i], NULL, 10);
3600 ir->opts.bSH[i] = (gmx_strncasecmp(ptr3[i], "Y", 1) == 0);
3603 nCASelec = str_nelem(is->CASelectrons, MAXPTR, ptr1);
3604 nCASorb = str_nelem(is->CASorbitals, MAXPTR, ptr2);
3605 snew(ir->opts.CASelectrons, nr);
3606 snew(ir->opts.CASorbitals, nr);
3607 for (i = 0; i < nr; i++)
3609 ir->opts.CASelectrons[i] = strtol(ptr1[i], NULL, 10);
3610 ir->opts.CASorbitals[i] = strtol(ptr2[i], NULL, 10);
3612 /* special optimization options */
3614 nbOPT = str_nelem(is->bOPT, MAXPTR, ptr1);
3615 nbTS = str_nelem(is->bTS, MAXPTR, ptr2);
3616 snew(ir->opts.bOPT, nr);
3617 snew(ir->opts.bTS, nr);
3618 for (i = 0; i < nr; i++)
3620 ir->opts.bOPT[i] = (gmx_strncasecmp(ptr1[i], "Y", 1) == 0);
3621 ir->opts.bTS[i] = (gmx_strncasecmp(ptr2[i], "Y", 1) == 0);
3623 nSAon = str_nelem(is->SAon, MAXPTR, ptr1);
3624 nSAoff = str_nelem(is->SAoff, MAXPTR, ptr2);
3625 nSAsteps = str_nelem(is->SAsteps, MAXPTR, ptr3);
3626 snew(ir->opts.SAon, nr);
3627 snew(ir->opts.SAoff, nr);
3628 snew(ir->opts.SAsteps, nr);
3630 for (i = 0; i < nr; i++)
3632 ir->opts.SAon[i] = strtod(ptr1[i], NULL);
3633 ir->opts.SAoff[i] = strtod(ptr2[i], NULL);
3634 ir->opts.SAsteps[i] = strtol(ptr3[i], NULL, 10);
3636 /* end of QMMM input */
3640 for (i = 0; (i < egcNR); i++)
3642 fprintf(stderr, "%-16s has %d element(s):", gtypes[i], groups->grps[i].nr);
3643 for (j = 0; (j < groups->grps[i].nr); j++)
3645 fprintf(stderr, " %s", *(groups->grpname[groups->grps[i].nm_ind[j]]));
3647 fprintf(stderr, "\n");
3651 nr = groups->grps[egcENER].nr;
3652 snew(ir->opts.egp_flags, nr*nr);
3654 bExcl = do_egp_flag(ir, groups, "energygrp-excl", is->egpexcl, EGP_EXCL);
3655 if (bExcl && ir->cutoff_scheme == ecutsVERLET)
3657 warning_error(wi, "Energy group exclusions are not (yet) implemented for the Verlet scheme");
3659 if (bExcl && EEL_FULL(ir->coulombtype))
3661 warning(wi, "Can not exclude the lattice Coulomb energy between energy groups");
3664 bTable = do_egp_flag(ir, groups, "energygrp-table", is->egptable, EGP_TABLE);
3665 if (bTable && !(ir->vdwtype == evdwUSER) &&
3666 !(ir->coulombtype == eelUSER) && !(ir->coulombtype == eelPMEUSER) &&
3667 !(ir->coulombtype == eelPMEUSERSWITCH))
3669 gmx_fatal(FARGS, "Can only have energy group pair tables in combination with user tables for VdW and/or Coulomb");
3672 decode_cos(is->efield_x, &(ir->ex[XX]));
3673 decode_cos(is->efield_xt, &(ir->et[XX]));
3674 decode_cos(is->efield_y, &(ir->ex[YY]));
3675 decode_cos(is->efield_yt, &(ir->et[YY]));
3676 decode_cos(is->efield_z, &(ir->ex[ZZ]));
3677 decode_cos(is->efield_zt, &(ir->et[ZZ]));
3681 do_adress_index(ir->adress, groups, gnames, &(ir->opts), wi);
3684 for (i = 0; (i < grps->nr); i++)
3696 static void check_disre(gmx_mtop_t *mtop)
3698 gmx_ffparams_t *ffparams;
3699 t_functype *functype;
3701 int i, ndouble, ftype;
3702 int label, old_label;
3704 if (gmx_mtop_ftype_count(mtop, F_DISRES) > 0)
3706 ffparams = &mtop->ffparams;
3707 functype = ffparams->functype;
3708 ip = ffparams->iparams;
3711 for (i = 0; i < ffparams->ntypes; i++)
3713 ftype = functype[i];
3714 if (ftype == F_DISRES)
3716 label = ip[i].disres.label;
3717 if (label == old_label)
3719 fprintf(stderr, "Distance restraint index %d occurs twice\n", label);
3727 gmx_fatal(FARGS, "Found %d double distance restraint indices,\n"
3728 "probably the parameters for multiple pairs in one restraint "
3729 "are not identical\n", ndouble);
3734 static gmx_bool absolute_reference(t_inputrec *ir, gmx_mtop_t *sys,
3735 gmx_bool posres_only,
3739 gmx_mtop_ilistloop_t iloop;
3749 for (d = 0; d < DIM; d++)
3751 AbsRef[d] = (d < ndof_com(ir) ? 0 : 1);
3753 /* Check for freeze groups */
3754 for (g = 0; g < ir->opts.ngfrz; g++)
3756 for (d = 0; d < DIM; d++)
3758 if (ir->opts.nFreeze[g][d] != 0)
3766 /* Check for position restraints */
3767 iloop = gmx_mtop_ilistloop_init(sys);
3768 while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
3771 (AbsRef[XX] == 0 || AbsRef[YY] == 0 || AbsRef[ZZ] == 0))
3773 for (i = 0; i < ilist[F_POSRES].nr; i += 2)
3775 pr = &sys->ffparams.iparams[ilist[F_POSRES].iatoms[i]];
3776 for (d = 0; d < DIM; d++)
3778 if (pr->posres.fcA[d] != 0)
3784 for (i = 0; i < ilist[F_FBPOSRES].nr; i += 2)
3786 /* Check for flat-bottom posres */
3787 pr = &sys->ffparams.iparams[ilist[F_FBPOSRES].iatoms[i]];
3788 if (pr->fbposres.k != 0)
3790 switch (pr->fbposres.geom)
3792 case efbposresSPHERE:
3793 AbsRef[XX] = AbsRef[YY] = AbsRef[ZZ] = 1;
3795 case efbposresCYLINDER:
3796 AbsRef[XX] = AbsRef[YY] = 1;
3798 case efbposresX: /* d=XX */
3799 case efbposresY: /* d=YY */
3800 case efbposresZ: /* d=ZZ */
3801 d = pr->fbposres.geom - efbposresX;
3805 gmx_fatal(FARGS, " Invalid geometry for flat-bottom position restraint.\n"
3806 "Expected nr between 1 and %d. Found %d\n", efbposresNR-1,
3814 return (AbsRef[XX] != 0 && AbsRef[YY] != 0 && AbsRef[ZZ] != 0);
3818 check_combination_rule_differences(const gmx_mtop_t *mtop, int state,
3819 gmx_bool *bC6ParametersWorkWithGeometricRules,
3820 gmx_bool *bC6ParametersWorkWithLBRules,
3821 gmx_bool *bLBRulesPossible)
3823 int ntypes, tpi, tpj, thisLBdiff, thisgeomdiff;
3826 double geometricdiff, LBdiff;
3827 double c6i, c6j, c12i, c12j;
3828 double c6, c6_geometric, c6_LB;
3829 double sigmai, sigmaj, epsi, epsj;
3830 gmx_bool bCanDoLBRules, bCanDoGeometricRules;
3833 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
3834 * force-field floating point parameters.
3837 ptr = getenv("GMX_LJCOMB_TOL");
3842 sscanf(ptr, "%lf", &dbl);
3846 *bC6ParametersWorkWithLBRules = TRUE;
3847 *bC6ParametersWorkWithGeometricRules = TRUE;
3848 bCanDoLBRules = TRUE;
3849 bCanDoGeometricRules = TRUE;
3850 ntypes = mtop->ffparams.atnr;
3851 snew(typecount, ntypes);
3852 gmx_mtop_count_atomtypes(mtop, state, typecount);
3853 geometricdiff = LBdiff = 0.0;
3854 *bLBRulesPossible = TRUE;
3855 for (tpi = 0; tpi < ntypes; ++tpi)
3857 c6i = mtop->ffparams.iparams[(ntypes + 1) * tpi].lj.c6;
3858 c12i = mtop->ffparams.iparams[(ntypes + 1) * tpi].lj.c12;
3859 for (tpj = tpi; tpj < ntypes; ++tpj)
3861 c6j = mtop->ffparams.iparams[(ntypes + 1) * tpj].lj.c6;
3862 c12j = mtop->ffparams.iparams[(ntypes + 1) * tpj].lj.c12;
3863 c6 = mtop->ffparams.iparams[ntypes * tpi + tpj].lj.c6;
3864 c6_geometric = sqrt(c6i * c6j);
3865 if (!gmx_numzero(c6_geometric))
3867 if (!gmx_numzero(c12i) && !gmx_numzero(c12j))
3869 sigmai = pow(c12i / c6i, 1.0/6.0);
3870 sigmaj = pow(c12j / c6j, 1.0/6.0);
3871 epsi = c6i * c6i /(4.0 * c12i);
3872 epsj = c6j * c6j /(4.0 * c12j);
3873 c6_LB = 4.0 * pow(epsi * epsj, 1.0/2.0) * pow(0.5 * (sigmai + sigmaj), 6);
3877 *bLBRulesPossible = FALSE;
3878 c6_LB = c6_geometric;
3880 bCanDoLBRules = gmx_within_tol(c6_LB, c6, tol);
3883 if (FALSE == bCanDoLBRules)
3885 *bC6ParametersWorkWithLBRules = FALSE;
3888 bCanDoGeometricRules = gmx_within_tol(c6_geometric, c6, tol);
3890 if (FALSE == bCanDoGeometricRules)
3892 *bC6ParametersWorkWithGeometricRules = FALSE;
3900 check_combination_rules(const t_inputrec *ir, const gmx_mtop_t *mtop,
3904 gmx_bool bLBRulesPossible, bC6ParametersWorkWithGeometricRules, bC6ParametersWorkWithLBRules;
3906 check_combination_rule_differences(mtop, 0,
3907 &bC6ParametersWorkWithGeometricRules,
3908 &bC6ParametersWorkWithLBRules,
3910 if (ir->ljpme_combination_rule == eljpmeLB)
3912 if (FALSE == bC6ParametersWorkWithLBRules || FALSE == bLBRulesPossible)
3914 warning(wi, "You are using arithmetic-geometric combination rules "
3915 "in LJ-PME, but your non-bonded C6 parameters do not "
3916 "follow these rules.");
3921 if (FALSE == bC6ParametersWorkWithGeometricRules)
3923 if (ir->eDispCorr != edispcNO)
3925 warning_note(wi, "You are using geometric combination rules in "
3926 "LJ-PME, but your non-bonded C6 parameters do "
3927 "not follow these rules. "
3928 "This will introduce very small errors in the forces and energies in "
3929 "your simulations. Dispersion correction will correct total energy "
3930 "and/or pressure for isotropic systems, but not forces or surface tensions.");
3934 warning_note(wi, "You are using geometric combination rules in "
3935 "LJ-PME, but your non-bonded C6 parameters do "
3936 "not follow these rules. "
3937 "This will introduce very small errors in the forces and energies in "
3938 "your simulations. If your system is homogeneous, consider using dispersion correction "
3939 "for the total energy and pressure.");
3945 void triple_check(const char *mdparin, t_inputrec *ir, gmx_mtop_t *sys,
3948 char err_buf[STRLEN];
3949 int i, m, c, nmol, npct;
3950 gmx_bool bCharge, bAcc;
3951 real gdt_max, *mgrp, mt;
3953 gmx_mtop_atomloop_block_t aloopb;
3954 gmx_mtop_atomloop_all_t aloop;
3957 char warn_buf[STRLEN];
3959 set_warning_line(wi, mdparin, -1);
3961 if (ir->cutoff_scheme == ecutsVERLET &&
3962 ir->verletbuf_tol > 0 &&
3964 ((EI_MD(ir->eI) || EI_SD(ir->eI)) &&
3965 (ir->etc == etcVRESCALE || ir->etc == etcBERENDSEN)))
3967 /* Check if a too small Verlet buffer might potentially
3968 * cause more drift than the thermostat can couple off.
3970 /* Temperature error fraction for warning and suggestion */
3971 const real T_error_warn = 0.002;
3972 const real T_error_suggest = 0.001;
3973 /* For safety: 2 DOF per atom (typical with constraints) */
3974 const real nrdf_at = 2;
3975 real T, tau, max_T_error;
3980 for (i = 0; i < ir->opts.ngtc; i++)
3982 T = max(T, ir->opts.ref_t[i]);
3983 tau = max(tau, ir->opts.tau_t[i]);
3987 /* This is a worst case estimate of the temperature error,
3988 * assuming perfect buffer estimation and no cancelation
3989 * of errors. The factor 0.5 is because energy distributes
3990 * equally over Ekin and Epot.
3992 max_T_error = 0.5*tau*ir->verletbuf_tol/(nrdf_at*BOLTZ*T);
3993 if (max_T_error > T_error_warn)
3995 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.",
3996 ir->verletbuf_tol, T, tau,
3998 100*T_error_suggest,
3999 ir->verletbuf_tol*T_error_suggest/max_T_error);
4000 warning(wi, warn_buf);
4005 if (ETC_ANDERSEN(ir->etc))
4009 for (i = 0; i < ir->opts.ngtc; i++)
4011 sprintf(err_buf, "all tau_t must currently be equal using Andersen temperature control, violated for group %d", i);
4012 CHECK(ir->opts.tau_t[0] != ir->opts.tau_t[i]);
4013 sprintf(err_buf, "all tau_t must be postive using Andersen temperature control, tau_t[%d]=%10.6f",
4014 i, ir->opts.tau_t[i]);
4015 CHECK(ir->opts.tau_t[i] < 0);
4018 for (i = 0; i < ir->opts.ngtc; i++)
4020 int nsteps = (int)(ir->opts.tau_t[i]/ir->delta_t);
4021 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);
4022 CHECK((nsteps % ir->nstcomm) && (ir->etc == etcANDERSENMASSIVE));
4026 if (EI_DYNAMICS(ir->eI) && !EI_SD(ir->eI) && ir->eI != eiBD &&
4027 ir->comm_mode == ecmNO &&
4028 !(absolute_reference(ir, sys, FALSE, AbsRef) || ir->nsteps <= 10) &&
4029 !ETC_ANDERSEN(ir->etc))
4031 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");
4034 /* Check for pressure coupling with absolute position restraints */
4035 if (ir->epc != epcNO && ir->refcoord_scaling == erscNO)
4037 absolute_reference(ir, sys, TRUE, AbsRef);
4039 for (m = 0; m < DIM; m++)
4041 if (AbsRef[m] && norm2(ir->compress[m]) > 0)
4043 warning(wi, "You are using pressure coupling with absolute position restraints, this will give artifacts. Use the refcoord_scaling option.");
4051 aloopb = gmx_mtop_atomloop_block_init(sys);
4052 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
4054 if (atom->q != 0 || atom->qB != 0)
4062 if (EEL_FULL(ir->coulombtype))
4065 "You are using full electrostatics treatment %s for a system without charges.\n"
4066 "This costs a lot of performance for just processing zeros, consider using %s instead.\n",
4067 EELTYPE(ir->coulombtype), EELTYPE(eelCUT));
4068 warning(wi, err_buf);
4073 if (ir->coulombtype == eelCUT && ir->rcoulomb > 0 && !ir->implicit_solvent)
4076 "You are using a plain Coulomb cut-off, which might produce artifacts.\n"
4077 "You might want to consider using %s electrostatics.\n",
4079 warning_note(wi, err_buf);
4083 /* Check if combination rules used in LJ-PME are the same as in the force field */
4084 if (EVDW_PME(ir->vdwtype))
4086 check_combination_rules(ir, sys, wi);
4089 /* Generalized reaction field */
4090 if (ir->opts.ngtc == 0)
4092 sprintf(err_buf, "No temperature coupling while using coulombtype %s",
4094 CHECK(ir->coulombtype == eelGRF);
4098 sprintf(err_buf, "When using coulombtype = %s"
4099 " ref-t for temperature coupling should be > 0",
4101 CHECK((ir->coulombtype == eelGRF) && (ir->opts.ref_t[0] <= 0));
4104 if (ir->eI == eiSD1 &&
4105 (gmx_mtop_ftype_count(sys, F_CONSTR) > 0 ||
4106 gmx_mtop_ftype_count(sys, F_SETTLE) > 0))
4108 sprintf(warn_buf, "With constraints integrator %s is less accurate, consider using %s instead", ei_names[ir->eI], ei_names[eiSD2]);
4109 warning_note(wi, warn_buf);
4113 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4115 for (m = 0; (m < DIM); m++)
4117 if (fabs(ir->opts.acc[i][m]) > 1e-6)
4126 snew(mgrp, sys->groups.grps[egcACC].nr);
4127 aloop = gmx_mtop_atomloop_all_init(sys);
4128 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
4130 mgrp[ggrpnr(&sys->groups, egcACC, i)] += atom->m;
4133 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4135 for (m = 0; (m < DIM); m++)
4137 acc[m] += ir->opts.acc[i][m]*mgrp[i];
4141 for (m = 0; (m < DIM); m++)
4143 if (fabs(acc[m]) > 1e-6)
4145 const char *dim[DIM] = { "X", "Y", "Z" };
4147 "Net Acceleration in %s direction, will %s be corrected\n",
4148 dim[m], ir->nstcomm != 0 ? "" : "not");
4149 if (ir->nstcomm != 0 && m < ndof_com(ir))
4152 for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4154 ir->opts.acc[i][m] -= acc[m];
4162 if (ir->efep != efepNO && ir->fepvals->sc_alpha != 0 &&
4163 !gmx_within_tol(sys->ffparams.reppow, 12.0, 10*GMX_DOUBLE_EPS))
4165 gmx_fatal(FARGS, "Soft-core interactions are only supported with VdW repulsion power 12");
4168 if (ir->ePull != epullNO)
4170 gmx_bool bPullAbsoluteRef;
4172 bPullAbsoluteRef = FALSE;
4173 for (i = 0; i < ir->pull->ncoord; i++)
4175 bPullAbsoluteRef = bPullAbsoluteRef ||
4176 ir->pull->coord[i].group[0] == 0 ||
4177 ir->pull->coord[i].group[1] == 0;
4179 if (bPullAbsoluteRef)
4181 absolute_reference(ir, sys, FALSE, AbsRef);
4182 for (m = 0; m < DIM; m++)
4184 if (ir->pull->dim[m] && !AbsRef[m])
4186 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.");
4192 if (ir->pull->eGeom == epullgDIRPBC)
4194 for (i = 0; i < 3; i++)
4196 for (m = 0; m <= i; m++)
4198 if ((ir->epc != epcNO && ir->compress[i][m] != 0) ||
4199 ir->deform[i][m] != 0)
4201 for (c = 0; c < ir->pull->ncoord; c++)
4203 if (ir->pull->coord[c].vec[m] != 0)
4205 gmx_fatal(FARGS, "Can not have dynamic box while using pull geometry '%s' (dim %c)", EPULLGEOM(ir->pull->eGeom), 'x'+m);
4217 void double_check(t_inputrec *ir, matrix box, gmx_bool bConstr, warninp_t wi)
4221 char warn_buf[STRLEN];
4224 ptr = check_box(ir->ePBC, box);
4227 warning_error(wi, ptr);
4230 if (bConstr && ir->eConstrAlg == econtSHAKE)
4232 if (ir->shake_tol <= 0.0)
4234 sprintf(warn_buf, "ERROR: shake-tol must be > 0 instead of %g\n",
4236 warning_error(wi, warn_buf);
4239 if (IR_TWINRANGE(*ir) && ir->nstlist > 1)
4241 sprintf(warn_buf, "With twin-range cut-off's and SHAKE the virial and the pressure are incorrect.");
4242 if (ir->epc == epcNO)
4244 warning(wi, warn_buf);
4248 warning_error(wi, warn_buf);
4253 if ( (ir->eConstrAlg == econtLINCS) && bConstr)
4255 /* If we have Lincs constraints: */
4256 if (ir->eI == eiMD && ir->etc == etcNO &&
4257 ir->eConstrAlg == econtLINCS && ir->nLincsIter == 1)
4259 sprintf(warn_buf, "For energy conservation with LINCS, lincs_iter should be 2 or larger.\n");
4260 warning_note(wi, warn_buf);
4263 if ((ir->eI == eiCG || ir->eI == eiLBFGS) && (ir->nProjOrder < 8))
4265 sprintf(warn_buf, "For accurate %s with LINCS constraints, lincs-order should be 8 or more.", ei_names[ir->eI]);
4266 warning_note(wi, warn_buf);
4268 if (ir->epc == epcMTTK)
4270 warning_error(wi, "MTTK not compatible with lincs -- use shake instead.");
4274 if (bConstr && ir->epc == epcMTTK)
4276 warning_note(wi, "MTTK with constraints is deprecated, and will be removed in GROMACS 5.1");
4279 if (ir->LincsWarnAngle > 90.0)
4281 sprintf(warn_buf, "lincs-warnangle can not be larger than 90 degrees, setting it to 90.\n");
4282 warning(wi, warn_buf);
4283 ir->LincsWarnAngle = 90.0;
4286 if (ir->ePBC != epbcNONE)
4288 if (ir->nstlist == 0)
4290 warning(wi, "With nstlist=0 atoms are only put into the box at step 0, therefore drifting atoms might cause the simulation to crash.");
4292 bTWIN = (ir->rlistlong > ir->rlist);
4293 if (ir->ns_type == ensGRID)
4295 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
4297 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",
4298 bTWIN ? (ir->rcoulomb == ir->rlistlong ? "rcoulomb" : "rvdw") : "rlist");
4299 warning_error(wi, warn_buf);
4304 min_size = min(box[XX][XX], min(box[YY][YY], box[ZZ][ZZ]));
4305 if (2*ir->rlistlong >= min_size)
4307 sprintf(warn_buf, "ERROR: One of the box lengths is smaller than twice the cut-off length. Increase the box size or decrease rlist.");
4308 warning_error(wi, warn_buf);
4311 fprintf(stderr, "Grid search might allow larger cut-off's than simple search with triclinic boxes.");
4318 void check_chargegroup_radii(const gmx_mtop_t *mtop, const t_inputrec *ir,
4322 real rvdw1, rvdw2, rcoul1, rcoul2;
4323 char warn_buf[STRLEN];
4325 calc_chargegroup_radii(mtop, x, &rvdw1, &rvdw2, &rcoul1, &rcoul2);
4329 printf("Largest charge group radii for Van der Waals: %5.3f, %5.3f nm\n",
4334 printf("Largest charge group radii for Coulomb: %5.3f, %5.3f nm\n",
4340 if (rvdw1 + rvdw2 > ir->rlist ||
4341 rcoul1 + rcoul2 > ir->rlist)
4344 "The sum of the two largest charge group radii (%f) "
4345 "is larger than rlist (%f)\n",
4346 max(rvdw1+rvdw2, rcoul1+rcoul2), ir->rlist);
4347 warning(wi, warn_buf);
4351 /* Here we do not use the zero at cut-off macro,
4352 * since user defined interactions might purposely
4353 * not be zero at the cut-off.
4355 if (ir_vdw_is_zero_at_cutoff(ir) &&
4356 rvdw1 + rvdw2 > ir->rlistlong - ir->rvdw)
4358 sprintf(warn_buf, "The sum of the two largest charge group "
4359 "radii (%f) is larger than %s (%f) - rvdw (%f).\n"
4360 "With exact cut-offs, better performance can be "
4361 "obtained with cutoff-scheme = %s, because it "
4362 "does not use charge groups at all.",
4364 ir->rlistlong > ir->rlist ? "rlistlong" : "rlist",
4365 ir->rlistlong, ir->rvdw,
4366 ecutscheme_names[ecutsVERLET]);
4369 warning(wi, warn_buf);
4373 warning_note(wi, warn_buf);
4376 if (ir_coulomb_is_zero_at_cutoff(ir) &&
4377 rcoul1 + rcoul2 > ir->rlistlong - ir->rcoulomb)
4379 sprintf(warn_buf, "The sum of the two largest charge group radii (%f) is larger than %s (%f) - rcoulomb (%f).\n"
4380 "With exact cut-offs, better performance can be obtained with cutoff-scheme = %s, because it does not use charge groups at all.",
4382 ir->rlistlong > ir->rlist ? "rlistlong" : "rlist",
4383 ir->rlistlong, ir->rcoulomb,
4384 ecutscheme_names[ecutsVERLET]);
4387 warning(wi, warn_buf);
4391 warning_note(wi, warn_buf);