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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
50 #include "gmx_fatal.h"
63 #include "mtop_util.h"
64 #include "chargegroup.h"
69 #define MAXLAMBDAS 1024
71 /* Resource parameters
72 * Do not change any of these until you read the instruction
73 * in readinp.h. Some cpp's do not take spaces after the backslash
74 * (like the c-shell), which will give you a very weird compiler
78 static char tcgrps[STRLEN],tau_t[STRLEN],ref_t[STRLEN],
79 acc[STRLEN],accgrps[STRLEN],freeze[STRLEN],frdim[STRLEN],
80 energy[STRLEN],user1[STRLEN],user2[STRLEN],vcm[STRLEN],xtc_grps[STRLEN],
81 couple_moltype[STRLEN],orirefitgrp[STRLEN],egptable[STRLEN],egpexcl[STRLEN],
82 wall_atomtype[STRLEN],wall_density[STRLEN],deform[STRLEN],QMMM[STRLEN];
83 static char fep_lambda[efptNR][STRLEN];
84 static char lambda_weights[STRLEN];
85 static char **pull_grp;
86 static char **rot_grp;
87 static char anneal[STRLEN],anneal_npoints[STRLEN],
88 anneal_time[STRLEN],anneal_temp[STRLEN];
89 static char QMmethod[STRLEN],QMbasis[STRLEN],QMcharge[STRLEN],QMmult[STRLEN],
90 bSH[STRLEN],CASorbitals[STRLEN], CASelectrons[STRLEN],SAon[STRLEN],
91 SAoff[STRLEN],SAsteps[STRLEN],bTS[STRLEN],bOPT[STRLEN];
92 static char efield_x[STRLEN],efield_xt[STRLEN],efield_y[STRLEN],
93 efield_yt[STRLEN],efield_z[STRLEN],efield_zt[STRLEN];
96 egrptpALL, /* All particles have to be a member of a group. */
97 egrptpALL_GENREST, /* A rest group with name is generated for particles *
98 * that are not part of any group. */
99 egrptpPART, /* As egrptpALL_GENREST, but no name is generated *
100 * for the rest group. */
101 egrptpONE /* Merge all selected groups into one group, *
102 * make a rest group for the remaining particles. */
106 void init_ir(t_inputrec *ir, t_gromppopts *opts)
108 snew(opts->include,STRLEN);
109 snew(opts->define,STRLEN);
111 snew(ir->expandedvals,1);
112 snew(ir->simtempvals,1);
115 static void GetSimTemps(int ntemps, t_simtemp *simtemp, double *temperature_lambdas)
120 for (i=0;i<ntemps;i++)
122 /* simple linear scaling -- allows more control */
123 if (simtemp->eSimTempScale == esimtempLINEAR)
125 simtemp->temperatures[i] = simtemp->simtemp_low + (simtemp->simtemp_high-simtemp->simtemp_low)*temperature_lambdas[i];
127 else if (simtemp->eSimTempScale == esimtempGEOMETRIC) /* should give roughly equal acceptance for constant heat capacity . . . */
129 simtemp->temperatures[i] = simtemp->simtemp_low * pow(simtemp->simtemp_high/simtemp->simtemp_low,(1.0*i)/(ntemps-1));
131 else if (simtemp->eSimTempScale == esimtempEXPONENTIAL)
133 simtemp->temperatures[i] = simtemp->simtemp_low + (simtemp->simtemp_high-simtemp->simtemp_low)*((exp(temperature_lambdas[i])-1)/(exp(1.0)-1));
138 sprintf(errorstr,"eSimTempScale=%d not defined",simtemp->eSimTempScale);
139 gmx_fatal(FARGS,errorstr);
146 static void _low_check(gmx_bool b,char *s,warninp_t wi)
154 static void check_nst(const char *desc_nst,int nst,
155 const char *desc_p,int *p,
160 if (*p > 0 && *p % nst != 0)
162 /* Round up to the next multiple of nst */
163 *p = ((*p)/nst + 1)*nst;
164 sprintf(buf,"%s should be a multiple of %s, changing %s to %d\n",
165 desc_p,desc_nst,desc_p,*p);
170 static gmx_bool ir_NVE(const t_inputrec *ir)
172 return ((ir->eI == eiMD || EI_VV(ir->eI)) && ir->etc == etcNO);
175 static int lcd(int n1,int n2)
180 for(i=2; (i<=n1 && i<=n2); i++)
182 if (n1 % i == 0 && n2 % i == 0)
191 static void process_interaction_modifier(const t_inputrec *ir,int *eintmod)
193 if (*eintmod == eintmodPOTSHIFT_VERLET)
195 if (ir->cutoff_scheme == ecutsVERLET)
197 *eintmod = eintmodPOTSHIFT;
201 *eintmod = eintmodNONE;
206 void check_ir(const char *mdparin,t_inputrec *ir, t_gromppopts *opts,
208 /* Check internal consistency */
210 /* Strange macro: first one fills the err_buf, and then one can check
211 * the condition, which will print the message and increase the error
214 #define CHECK(b) _low_check(b,err_buf,wi)
215 char err_buf[256],warn_buf[STRLEN];
221 t_lambda *fep = ir->fepvals;
222 t_expanded *expand = ir->expandedvals;
224 set_warning_line(wi,mdparin,-1);
226 /* BASIC CUT-OFF STUFF */
227 if (ir->rcoulomb < 0)
229 warning_error(wi,"rcoulomb should be >= 0");
233 warning_error(wi,"rvdw should be >= 0");
236 !(ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_drift > 0))
238 warning_error(wi,"rlist should be >= 0");
241 process_interaction_modifier(ir,&ir->coulomb_modifier);
242 process_interaction_modifier(ir,&ir->vdw_modifier);
244 if (ir->cutoff_scheme == ecutsGROUP)
246 /* BASIC CUT-OFF STUFF */
247 if (ir->rlist == 0 ||
248 !((EEL_MIGHT_BE_ZERO_AT_CUTOFF(ir->coulombtype) && ir->rcoulomb > ir->rlist) ||
249 (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(ir->vdwtype) && ir->rvdw > ir->rlist))) {
250 /* No switched potential and/or no twin-range:
251 * we can set the long-range cut-off to the maximum of the other cut-offs.
253 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
255 else if (ir->rlistlong < 0)
257 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
258 sprintf(warn_buf,"rlistlong was not set, setting it to %g (no buffer)",
260 warning(wi,warn_buf);
262 if (ir->rlistlong == 0 && ir->ePBC != epbcNONE)
264 warning_error(wi,"Can not have an infinite cut-off with PBC");
266 if (ir->rlistlong > 0 && (ir->rlist == 0 || ir->rlistlong < ir->rlist))
268 warning_error(wi,"rlistlong can not be shorter than rlist");
270 if (IR_TWINRANGE(*ir) && ir->nstlist <= 0)
272 warning_error(wi,"Can not have nstlist<=0 with twin-range interactions");
276 if(ir->rlistlong == ir->rlist)
280 else if(ir->rlistlong>ir->rlist && ir->nstcalclr==0)
282 warning_error(wi,"With different cutoffs for electrostatics and VdW, nstcalclr must be -1 or a positive number");
285 if (ir->cutoff_scheme == ecutsVERLET)
289 /* Normal Verlet type neighbor-list, currently only limited feature support */
290 if (inputrec2nboundeddim(ir) < 3)
292 warning_error(wi,"With Verlet lists only full pbc or pbc=xy with walls is supported");
294 if (ir->rcoulomb != ir->rvdw)
296 warning_error(wi,"With Verlet lists rcoulomb!=rvdw is not supported");
298 if (ir->vdwtype != evdwCUT)
300 warning_error(wi,"With Verlet lists only cut-off LJ interactions are supported");
302 if (!(ir->coulombtype == eelCUT ||
303 (EEL_RF(ir->coulombtype) && ir->coulombtype != eelRF_NEC) ||
304 EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD))
306 warning_error(wi,"With Verlet lists only cut-off, reaction-field, PME and Ewald electrostatics are supported");
309 if (ir->nstlist <= 0)
311 warning_error(wi,"With Verlet lists nstlist should be larger than 0");
314 if (ir->nstlist < 10)
316 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.");
319 rc_max = max(ir->rvdw,ir->rcoulomb);
321 if (ir->verletbuf_drift <= 0)
323 if (ir->verletbuf_drift == 0)
325 warning_error(wi,"Can not have an energy drift of exactly 0");
328 if (ir->rlist < rc_max)
330 warning_error(wi,"With verlet lists rlist can not be smaller than rvdw or rcoulomb");
333 if (ir->rlist == rc_max && ir->nstlist > 1)
335 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.");
340 if (ir->rlist > rc_max)
342 warning_note(wi,"You have set rlist larger than the interaction cut-off, but you also have verlet-buffer-drift > 0. Will set rlist using verlet-buffer-drift.");
345 if (ir->nstlist == 1)
347 /* No buffer required */
352 if (EI_DYNAMICS(ir->eI))
354 if (EI_MD(ir->eI) && ir->etc == etcNO)
356 warning_error(wi,"Temperature coupling is required for calculating rlist using the energy drift with verlet-buffer-drift > 0. Either use temperature coupling or set rlist yourself together with verlet-buffer-drift = -1.");
359 if (inputrec2nboundeddim(ir) < 3)
361 warning_error(wi,"The box volume is required for calculating rlist from the energy drift with verlet-buffer-drift > 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-drift = -1.");
363 /* Set rlist temporarily so we can continue processing */
368 /* Set the buffer to 5% of the cut-off */
369 ir->rlist = 1.05*rc_max;
374 /* No twin-range calculations with Verlet lists */
375 ir->rlistlong = ir->rlist;
378 if(ir->nstcalclr==-1)
380 /* if rlist=rlistlong, this will later be changed to nstcalclr=0 */
381 ir->nstcalclr = ir->nstlist;
383 else if(ir->nstcalclr>0)
385 if(ir->nstlist>0 && (ir->nstlist % ir->nstcalclr != 0))
387 warning_error(wi,"nstlist must be evenly divisible by nstcalclr. Use nstcalclr = -1 to automatically follow nstlist");
390 else if(ir->nstcalclr<-1)
392 warning_error(wi,"nstcalclr must be a positive number (divisor of nstcalclr), or -1 to follow nstlist.");
395 if(EEL_PME(ir->coulombtype) && ir->rcoulomb > ir->rvdw && ir->nstcalclr>1)
397 warning_error(wi,"When used with PME, the long-range component of twin-range interactions must be updated every step (nstcalclr)");
400 /* GENERAL INTEGRATOR STUFF */
401 if (!(ir->eI == eiMD || EI_VV(ir->eI)))
405 if (ir->eI == eiVVAK) {
406 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]);
407 warning_note(wi,warn_buf);
409 if (!EI_DYNAMICS(ir->eI))
413 if (EI_DYNAMICS(ir->eI))
415 if (ir->nstcalcenergy < 0)
417 ir->nstcalcenergy = ir_optimal_nstcalcenergy(ir);
418 if (ir->nstenergy != 0 && ir->nstenergy < ir->nstcalcenergy)
420 /* nstcalcenergy larger than nstener does not make sense.
421 * We ideally want nstcalcenergy=nstener.
425 ir->nstcalcenergy = lcd(ir->nstenergy,ir->nstlist);
429 ir->nstcalcenergy = ir->nstenergy;
433 else if (ir->nstenergy > 0 && ir->nstcalcenergy > ir->nstenergy)
435 /* If the user sets nstenergy small, we should respect that */
436 sprintf(warn_buf,"Setting nstcalcenergy (%d) equal to nstenergy (%d)",ir->nstcalcenergy,ir->nstenergy);
437 ir->nstcalcenergy = ir->nstenergy;
440 if (ir->epc != epcNO)
442 if (ir->nstpcouple < 0)
444 ir->nstpcouple = ir_optimal_nstpcouple(ir);
447 if (IR_TWINRANGE(*ir))
449 check_nst("nstlist",ir->nstlist,
450 "nstcalcenergy",&ir->nstcalcenergy,wi);
451 if (ir->epc != epcNO)
453 check_nst("nstlist",ir->nstlist,
454 "nstpcouple",&ir->nstpcouple,wi);
458 if (ir->nstcalcenergy > 1)
460 /* for storing exact averages nstenergy should be
461 * a multiple of nstcalcenergy
463 check_nst("nstcalcenergy",ir->nstcalcenergy,
464 "nstenergy",&ir->nstenergy,wi);
465 if (ir->efep != efepNO)
467 /* nstdhdl should be a multiple of nstcalcenergy */
468 check_nst("nstcalcenergy",ir->nstcalcenergy,
469 "nstdhdl",&ir->fepvals->nstdhdl,wi);
470 /* nstexpanded should be a multiple of nstcalcenergy */
471 check_nst("nstcalcenergy",ir->nstcalcenergy,
472 "nstdhdl",&ir->expandedvals->nstexpanded,wi);
478 if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
479 ir->bContinuation && ir->ld_seed != -1) {
480 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)");
484 if (EI_TPI(ir->eI)) {
485 sprintf(err_buf,"TPI only works with pbc = %s",epbc_names[epbcXYZ]);
486 CHECK(ir->ePBC != epbcXYZ);
487 sprintf(err_buf,"TPI only works with ns = %s",ens_names[ensGRID]);
488 CHECK(ir->ns_type != ensGRID);
489 sprintf(err_buf,"with TPI nstlist should be larger than zero");
490 CHECK(ir->nstlist <= 0);
491 sprintf(err_buf,"TPI does not work with full electrostatics other than PME");
492 CHECK(EEL_FULL(ir->coulombtype) && !EEL_PME(ir->coulombtype));
496 if ( (opts->nshake > 0) && (opts->bMorse) ) {
498 "Using morse bond-potentials while constraining bonds is useless");
499 warning(wi,warn_buf);
502 if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
503 ir->bContinuation && ir->ld_seed != -1) {
504 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)");
506 /* verify simulated tempering options */
509 gmx_bool bAllTempZero = TRUE;
510 for (i=0;i<fep->n_lambda;i++)
512 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]);
513 CHECK((fep->all_lambda[efptTEMPERATURE][i] < 0) || (fep->all_lambda[efptTEMPERATURE][i] > 1));
514 if (fep->all_lambda[efptTEMPERATURE][i] > 0)
516 bAllTempZero = FALSE;
519 sprintf(err_buf,"if simulated tempering is on, temperature-lambdas may not be all zero");
520 CHECK(bAllTempZero==TRUE);
522 sprintf(err_buf,"Simulated tempering is currently only compatible with md-vv");
523 CHECK(ir->eI != eiVV);
525 /* check compatability of the temperature coupling with simulated tempering */
527 if (ir->etc == etcNOSEHOOVER) {
528 sprintf(warn_buf,"Nose-Hoover based temperature control such as [%s] my not be entirelyconsistent with simulated tempering",etcoupl_names[ir->etc]);
529 warning_note(wi,warn_buf);
532 /* check that the temperatures make sense */
534 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);
535 CHECK(ir->simtempvals->simtemp_high <= ir->simtempvals->simtemp_low);
537 sprintf(err_buf,"Higher simulated tempering temperature (%g) must be >= zero",ir->simtempvals->simtemp_high);
538 CHECK(ir->simtempvals->simtemp_high <= 0);
540 sprintf(err_buf,"Lower simulated tempering temperature (%g) must be >= zero",ir->simtempvals->simtemp_low);
541 CHECK(ir->simtempvals->simtemp_low <= 0);
544 /* verify free energy options */
546 if (ir->efep != efepNO) {
548 sprintf(err_buf,"The soft-core power is %d and can only be 1 or 2",
550 CHECK(fep->sc_alpha!=0 && fep->sc_power!=1 && fep->sc_power!=2);
552 sprintf(err_buf,"The soft-core sc-r-power is %d and can only be 6 or 48",
553 (int)fep->sc_r_power);
554 CHECK(fep->sc_alpha!=0 && fep->sc_r_power!=6.0 && fep->sc_r_power!=48.0);
556 /* check validity of options */
557 if (fep->n_lambda > 0 && ir->rlist < max(ir->rvdw,ir->rcoulomb))
560 "For foreign lambda free energy differences it is assumed that the soft-core interactions have no effect beyond the neighborlist cut-off");
561 warning(wi,warn_buf);
564 sprintf(err_buf,"Can't use postive delta-lambda (%g) if initial state/lambda does not start at zero",fep->delta_lambda);
565 CHECK(fep->delta_lambda > 0 && ((fep->init_fep_state !=0) || (fep->init_lambda !=0)));
567 sprintf(err_buf,"Can't use postive delta-lambda (%g) with expanded ensemble simulations",fep->delta_lambda);
568 CHECK(fep->delta_lambda > 0 && (ir->efep == efepEXPANDED));
570 sprintf(err_buf,"Free-energy not implemented for Ewald");
571 CHECK(ir->coulombtype==eelEWALD);
573 /* check validty of lambda inputs */
574 sprintf(err_buf,"initial thermodynamic state %d does not exist, only goes to %d",fep->init_fep_state,fep->n_lambda);
575 CHECK((fep->init_fep_state > fep->n_lambda));
577 for (j=0;j<efptNR;j++)
579 for (i=0;i<fep->n_lambda;i++)
581 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]);
582 CHECK((fep->all_lambda[j][i] < 0) || (fep->all_lambda[j][i] > 1));
586 if ((fep->sc_alpha>0) && (!fep->bScCoul))
588 for (i=0;i<fep->n_lambda;i++)
590 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],
591 fep->all_lambda[efptCOUL][i]);
592 CHECK((fep->sc_alpha>0) &&
593 (((fep->all_lambda[efptCOUL][i] > 0.0) &&
594 (fep->all_lambda[efptCOUL][i] < 1.0)) &&
595 ((fep->all_lambda[efptVDW][i] > 0.0) &&
596 (fep->all_lambda[efptVDW][i] < 1.0))));
600 if ((fep->bScCoul) && (EEL_PME(ir->coulombtype)))
602 sprintf(warn_buf,"With coulomb soft core, the reciprocal space calculation will not necessarily cancel. It may be necessary to decrease the reciprocal space energy, and increase the cutoff radius to get sufficiently close matches to energies with free energy turned off.");
603 warning(wi, warn_buf);
606 /* Free Energy Checks -- In an ideal world, slow growth and FEP would
607 be treated differently, but that's the next step */
609 for (i=0;i<efptNR;i++) {
610 for (j=0;j<fep->n_lambda;j++) {
611 sprintf(err_buf,"%s[%d] must be between 0 and 1",efpt_names[i],j);
612 CHECK((fep->all_lambda[i][j] < 0) || (fep->all_lambda[i][j] > 1));
617 if ((ir->bSimTemp) || (ir->efep == efepEXPANDED)) {
619 expand = ir->expandedvals;
621 /* checking equilibration of weights inputs for validity */
623 sprintf(err_buf,"weight-equil-number-all-lambda (%d) is ignored if lmc-weights-equil is not equal to %s",
624 expand->equil_n_at_lam,elmceq_names[elmceqNUMATLAM]);
625 CHECK((expand->equil_n_at_lam>0) && (expand->elmceq!=elmceqNUMATLAM));
627 sprintf(err_buf,"weight-equil-number-samples (%d) is ignored if lmc-weights-equil is not equal to %s",
628 expand->equil_samples,elmceq_names[elmceqSAMPLES]);
629 CHECK((expand->equil_samples>0) && (expand->elmceq!=elmceqSAMPLES));
631 sprintf(err_buf,"weight-equil-number-steps (%d) is ignored if lmc-weights-equil is not equal to %s",
632 expand->equil_steps,elmceq_names[elmceqSTEPS]);
633 CHECK((expand->equil_steps>0) && (expand->elmceq!=elmceqSTEPS));
635 sprintf(err_buf,"weight-equil-wl-delta (%d) is ignored if lmc-weights-equil is not equal to %s",
636 expand->equil_samples,elmceq_names[elmceqWLDELTA]);
637 CHECK((expand->equil_wl_delta>0) && (expand->elmceq!=elmceqWLDELTA));
639 sprintf(err_buf,"weight-equil-count-ratio (%f) is ignored if lmc-weights-equil is not equal to %s",
640 expand->equil_ratio,elmceq_names[elmceqRATIO]);
641 CHECK((expand->equil_ratio>0) && (expand->elmceq!=elmceqRATIO));
643 sprintf(err_buf,"weight-equil-number-all-lambda (%d) must be a positive integer if lmc-weights-equil=%s",
644 expand->equil_n_at_lam,elmceq_names[elmceqNUMATLAM]);
645 CHECK((expand->equil_n_at_lam<=0) && (expand->elmceq==elmceqNUMATLAM));
647 sprintf(err_buf,"weight-equil-number-samples (%d) must be a positive integer if lmc-weights-equil=%s",
648 expand->equil_samples,elmceq_names[elmceqSAMPLES]);
649 CHECK((expand->equil_samples<=0) && (expand->elmceq==elmceqSAMPLES));
651 sprintf(err_buf,"weight-equil-number-steps (%d) must be a positive integer if lmc-weights-equil=%s",
652 expand->equil_steps,elmceq_names[elmceqSTEPS]);
653 CHECK((expand->equil_steps<=0) && (expand->elmceq==elmceqSTEPS));
655 sprintf(err_buf,"weight-equil-wl-delta (%f) must be > 0 if lmc-weights-equil=%s",
656 expand->equil_wl_delta,elmceq_names[elmceqWLDELTA]);
657 CHECK((expand->equil_wl_delta<=0) && (expand->elmceq==elmceqWLDELTA));
659 sprintf(err_buf,"weight-equil-count-ratio (%f) must be > 0 if lmc-weights-equil=%s",
660 expand->equil_ratio,elmceq_names[elmceqRATIO]);
661 CHECK((expand->equil_ratio<=0) && (expand->elmceq==elmceqRATIO));
663 sprintf(err_buf,"lmc-weights-equil=%s only possible when lmc-stats = %s or lmc-stats %s",
664 elmceq_names[elmceqWLDELTA],elamstats_names[elamstatsWL],elamstats_names[elamstatsWWL]);
665 CHECK((expand->elmceq==elmceqWLDELTA) && (!EWL(expand->elamstats)));
667 sprintf(err_buf,"lmc-repeats (%d) must be greater than 0",expand->lmc_repeats);
668 CHECK((expand->lmc_repeats <= 0));
669 sprintf(err_buf,"minimum-var-min (%d) must be greater than 0",expand->minvarmin);
670 CHECK((expand->minvarmin <= 0));
671 sprintf(err_buf,"weight-c-range (%d) must be greater or equal to 0",expand->c_range);
672 CHECK((expand->c_range < 0));
673 sprintf(err_buf,"init-lambda-state (%d) must be zero if lmc-forced-nstart (%d)> 0 and lmc-move != 'no'",
674 fep->init_fep_state, expand->lmc_forced_nstart);
675 CHECK((fep->init_fep_state!=0) && (expand->lmc_forced_nstart>0) && (expand->elmcmove!=elmcmoveNO));
676 sprintf(err_buf,"lmc-forced-nstart (%d) must not be negative",expand->lmc_forced_nstart);
677 CHECK((expand->lmc_forced_nstart < 0));
678 sprintf(err_buf,"init-lambda-state (%d) must be in the interval [0,number of lambdas)",fep->init_fep_state);
679 CHECK((fep->init_fep_state < 0) || (fep->init_fep_state >= fep->n_lambda));
681 sprintf(err_buf,"init-wl-delta (%f) must be greater than or equal to 0",expand->init_wl_delta);
682 CHECK((expand->init_wl_delta < 0));
683 sprintf(err_buf,"wl-ratio (%f) must be between 0 and 1",expand->wl_ratio);
684 CHECK((expand->wl_ratio <= 0) || (expand->wl_ratio >= 1));
685 sprintf(err_buf,"wl-scale (%f) must be between 0 and 1",expand->wl_scale);
686 CHECK((expand->wl_scale <= 0) || (expand->wl_scale >= 1));
688 /* if there is no temperature control, we need to specify an MC temperature */
689 sprintf(err_buf,"If there is no temperature control, and lmc-mcmove!= 'no',mc_temperature must be set to a positive number");
690 if (expand->nstTij > 0)
692 sprintf(err_buf,"nst-transition-matrix (%d) must be an integer multiple of nstlog (%d)",
693 expand->nstTij,ir->nstlog);
694 CHECK((mod(expand->nstTij,ir->nstlog)!=0));
699 sprintf(err_buf,"walls only work with pbc=%s",epbc_names[epbcXY]);
700 CHECK(ir->nwall && ir->ePBC!=epbcXY);
703 if (ir->ePBC != epbcXYZ && ir->nwall != 2) {
704 if (ir->ePBC == epbcNONE) {
705 if (ir->epc != epcNO) {
706 warning(wi,"Turning off pressure coupling for vacuum system");
710 sprintf(err_buf,"Can not have pressure coupling with pbc=%s",
711 epbc_names[ir->ePBC]);
712 CHECK(ir->epc != epcNO);
714 sprintf(err_buf,"Can not have Ewald with pbc=%s",epbc_names[ir->ePBC]);
715 CHECK(EEL_FULL(ir->coulombtype));
717 sprintf(err_buf,"Can not have dispersion correction with pbc=%s",
718 epbc_names[ir->ePBC]);
719 CHECK(ir->eDispCorr != edispcNO);
722 if (ir->rlist == 0.0) {
723 sprintf(err_buf,"can only have neighborlist cut-off zero (=infinite)\n"
724 "with coulombtype = %s or coulombtype = %s\n"
725 "without periodic boundary conditions (pbc = %s) and\n"
726 "rcoulomb and rvdw set to zero",
727 eel_names[eelCUT],eel_names[eelUSER],epbc_names[epbcNONE]);
728 CHECK(((ir->coulombtype != eelCUT) && (ir->coulombtype != eelUSER)) ||
729 (ir->ePBC != epbcNONE) ||
730 (ir->rcoulomb != 0.0) || (ir->rvdw != 0.0));
732 if (ir->nstlist < 0) {
733 warning_error(wi,"Can not have heuristic neighborlist updates without cut-off");
735 if (ir->nstlist > 0) {
736 warning_note(wi,"Simulating without cut-offs is usually (slightly) faster with nstlist=0, nstype=simple and particle decomposition");
741 if (ir->nstcomm == 0) {
742 ir->comm_mode = ecmNO;
744 if (ir->comm_mode != ecmNO) {
745 if (ir->nstcomm < 0) {
746 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");
747 ir->nstcomm = abs(ir->nstcomm);
750 if (ir->nstcalcenergy > 0 && ir->nstcomm < ir->nstcalcenergy) {
751 warning_note(wi,"nstcomm < nstcalcenergy defeats the purpose of nstcalcenergy, setting nstcomm to nstcalcenergy");
752 ir->nstcomm = ir->nstcalcenergy;
755 if (ir->comm_mode == ecmANGULAR) {
756 sprintf(err_buf,"Can not remove the rotation around the center of mass with periodic molecules");
757 CHECK(ir->bPeriodicMols);
758 if (ir->ePBC != epbcNONE)
759 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).");
763 if (EI_STATE_VELOCITY(ir->eI) && ir->ePBC == epbcNONE && ir->comm_mode != ecmANGULAR) {
764 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.");
767 sprintf(err_buf,"Twin-range neighbour searching (NS) with simple NS"
768 " algorithm not implemented");
769 CHECK(((ir->rcoulomb > ir->rlist) || (ir->rvdw > ir->rlist))
770 && (ir->ns_type == ensSIMPLE));
772 /* TEMPERATURE COUPLING */
773 if (ir->etc == etcYES)
775 ir->etc = etcBERENDSEN;
776 warning_note(wi,"Old option for temperature coupling given: "
777 "changing \"yes\" to \"Berendsen\"\n");
780 if ((ir->etc == etcNOSEHOOVER) || (ir->epc == epcMTTK))
782 if (ir->opts.nhchainlength < 1)
784 sprintf(warn_buf,"number of Nose-Hoover chains (currently %d) cannot be less than 1,reset to 1\n",ir->opts.nhchainlength);
785 ir->opts.nhchainlength =1;
786 warning(wi,warn_buf);
789 if (ir->etc==etcNOSEHOOVER && !EI_VV(ir->eI) && ir->opts.nhchainlength > 1)
791 warning_note(wi,"leapfrog does not yet support Nose-Hoover chains, nhchainlength reset to 1");
792 ir->opts.nhchainlength = 1;
797 ir->opts.nhchainlength = 0;
800 if (ir->eI == eiVVAK) {
801 sprintf(err_buf,"%s implemented primarily for validation, and requires nsttcouple = 1 and nstpcouple = 1.",
803 CHECK((ir->nsttcouple != 1) || (ir->nstpcouple != 1));
806 if (ETC_ANDERSEN(ir->etc))
808 sprintf(err_buf,"%s temperature control not supported for integrator %s.",etcoupl_names[ir->etc],ei_names[ir->eI]);
809 CHECK(!(EI_VV(ir->eI)));
811 for (i=0;i<ir->opts.ngtc;i++)
813 sprintf(err_buf,"all tau_t must currently be equal using Andersen temperature control, violated for group %d",i);
814 CHECK(ir->opts.tau_t[0] != ir->opts.tau_t[i]);
815 sprintf(err_buf,"all tau_t must be postive using Andersen temperature control, tau_t[%d]=%10.6f",
816 i,ir->opts.tau_t[i]);
817 CHECK(ir->opts.tau_t[i]<0);
819 if (ir->nstcomm > 0 && (ir->etc == etcANDERSEN)) {
820 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]);
821 warning_note(wi,warn_buf);
824 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]);
825 CHECK(ir->nstcomm > 1 && (ir->etc == etcANDERSEN));
827 for (i=0;i<ir->opts.ngtc;i++)
829 int nsteps = (int)(ir->opts.tau_t[i]/ir->delta_t);
830 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);
831 CHECK((nsteps % ir->nstcomm) && (ir->etc == etcANDERSENMASSIVE));
834 if (ir->etc == etcBERENDSEN)
836 sprintf(warn_buf,"The %s thermostat does not generate the correct kinetic energy distribution. You might want to consider using the %s thermostat.",
837 ETCOUPLTYPE(ir->etc),ETCOUPLTYPE(etcVRESCALE));
838 warning_note(wi,warn_buf);
841 if ((ir->etc==etcNOSEHOOVER || ETC_ANDERSEN(ir->etc))
842 && ir->epc==epcBERENDSEN)
844 sprintf(warn_buf,"Using Berendsen pressure coupling invalidates the "
845 "true ensemble for the thermostat");
846 warning(wi,warn_buf);
849 /* PRESSURE COUPLING */
850 if (ir->epc == epcISOTROPIC)
852 ir->epc = epcBERENDSEN;
853 warning_note(wi,"Old option for pressure coupling given: "
854 "changing \"Isotropic\" to \"Berendsen\"\n");
857 if (ir->epc != epcNO)
859 dt_pcoupl = ir->nstpcouple*ir->delta_t;
861 sprintf(err_buf,"tau-p must be > 0 instead of %g\n",ir->tau_p);
862 CHECK(ir->tau_p <= 0);
864 if (ir->tau_p/dt_pcoupl < pcouple_min_integration_steps(ir->epc))
866 sprintf(warn_buf,"For proper integration of the %s barostat, tau-p (%g) should be at least %d times larger than nstpcouple*dt (%g)",
867 EPCOUPLTYPE(ir->epc),ir->tau_p,pcouple_min_integration_steps(ir->epc),dt_pcoupl);
868 warning(wi,warn_buf);
871 sprintf(err_buf,"compressibility must be > 0 when using pressure"
872 " coupling %s\n",EPCOUPLTYPE(ir->epc));
873 CHECK(ir->compress[XX][XX] < 0 || ir->compress[YY][YY] < 0 ||
874 ir->compress[ZZ][ZZ] < 0 ||
875 (trace(ir->compress) == 0 && ir->compress[YY][XX] <= 0 &&
876 ir->compress[ZZ][XX] <= 0 && ir->compress[ZZ][YY] <= 0));
878 if (epcPARRINELLORAHMAN == ir->epc && opts->bGenVel)
881 "You are generating velocities so I am assuming you "
882 "are equilibrating a system. You are using "
883 "%s pressure coupling, but this can be "
884 "unstable for equilibration. If your system crashes, try "
885 "equilibrating first with Berendsen pressure coupling. If "
886 "you are not equilibrating the system, you can probably "
887 "ignore this warning.",
888 epcoupl_names[ir->epc]);
889 warning(wi,warn_buf);
897 if ((ir->epc!=epcBERENDSEN) && (ir->epc!=epcMTTK))
899 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.");
905 /* More checks are in triple check (grompp.c) */
907 if (ir->coulombtype == eelSWITCH) {
908 sprintf(warn_buf,"coulombtype = %s is only for testing purposes and can lead to serious "
909 "artifacts, advice: use coulombtype = %s",
910 eel_names[ir->coulombtype],
911 eel_names[eelRF_ZERO]);
912 warning(wi,warn_buf);
915 if (ir->epsilon_r!=1 && ir->implicit_solvent==eisGBSA) {
916 sprintf(warn_buf,"epsilon-r = %g with GB implicit solvent, will use this value for inner dielectric",ir->epsilon_r);
917 warning_note(wi,warn_buf);
920 if (EEL_RF(ir->coulombtype) && ir->epsilon_rf==1 && ir->epsilon_r!=1) {
921 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);
922 warning(wi,warn_buf);
923 ir->epsilon_rf = ir->epsilon_r;
927 if (getenv("GALACTIC_DYNAMICS") == NULL) {
928 sprintf(err_buf,"epsilon-r must be >= 0 instead of %g\n",ir->epsilon_r);
929 CHECK(ir->epsilon_r < 0);
932 if (EEL_RF(ir->coulombtype)) {
933 /* reaction field (at the cut-off) */
935 if (ir->coulombtype == eelRF_ZERO) {
936 sprintf(warn_buf,"With coulombtype = %s, epsilon-rf must be 0, assuming you meant epsilon_rf=0",
937 eel_names[ir->coulombtype]);
938 CHECK(ir->epsilon_rf != 0);
939 ir->epsilon_rf = 0.0;
942 sprintf(err_buf,"epsilon-rf must be >= epsilon-r");
943 CHECK((ir->epsilon_rf < ir->epsilon_r && ir->epsilon_rf != 0) ||
944 (ir->epsilon_r == 0));
945 if (ir->epsilon_rf == ir->epsilon_r) {
946 sprintf(warn_buf,"Using epsilon-rf = epsilon-r with %s does not make sense",
947 eel_names[ir->coulombtype]);
948 warning(wi,warn_buf);
951 /* Allow rlist>rcoulomb for tabulated long range stuff. This just
952 * means the interaction is zero outside rcoulomb, but it helps to
953 * provide accurate energy conservation.
955 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(ir->coulombtype)) {
956 if (EEL_SWITCHED(ir->coulombtype)) {
958 "With coulombtype = %s rcoulomb_switch must be < rcoulomb. Or, better: Use the potential modifier options!",
959 eel_names[ir->coulombtype]);
960 CHECK(ir->rcoulomb_switch >= ir->rcoulomb);
962 } else if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype)) {
963 if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE) {
964 sprintf(err_buf,"With coulombtype = %s, rcoulomb should be >= rlist unless you use a potential modifier",
965 eel_names[ir->coulombtype]);
966 CHECK(ir->rlist > ir->rcoulomb);
970 if(ir->coulombtype==eelSWITCH || ir->coulombtype==eelSHIFT ||
971 ir->vdwtype==evdwSWITCH || ir->vdwtype==evdwSHIFT)
974 "The switch/shift interaction settings are just for compatibility; you will get better"
975 "performance from applying potential modifiers to your interactions!\n");
976 warning_note(wi,warn_buf);
979 if (EEL_FULL(ir->coulombtype))
981 if (ir->coulombtype==eelPMESWITCH || ir->coulombtype==eelPMEUSER ||
982 ir->coulombtype==eelPMEUSERSWITCH)
984 sprintf(err_buf,"With coulombtype = %s, rcoulomb must be <= rlist",
985 eel_names[ir->coulombtype]);
986 CHECK(ir->rcoulomb > ir->rlist);
988 else if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE)
990 if (ir->coulombtype == eelPME || ir->coulombtype == eelP3M_AD)
993 "With coulombtype = %s (without modifier), rcoulomb must be equal to rlist,\n"
994 "or rlistlong if nstcalclr=1. For optimal energy conservation,consider using\n"
995 "a potential modifier.",eel_names[ir->coulombtype]);
998 CHECK(ir->rcoulomb != ir->rlist && ir->rcoulomb != ir->rlistlong);
1002 CHECK(ir->rcoulomb != ir->rlist);
1008 if (EEL_PME(ir->coulombtype)) {
1009 if (ir->pme_order < 3) {
1010 warning_error(wi,"pme-order can not be smaller than 3");
1014 if (ir->nwall==2 && EEL_FULL(ir->coulombtype)) {
1015 if (ir->ewald_geometry == eewg3D) {
1016 sprintf(warn_buf,"With pbc=%s you should use ewald-geometry=%s",
1017 epbc_names[ir->ePBC],eewg_names[eewg3DC]);
1018 warning(wi,warn_buf);
1020 /* This check avoids extra pbc coding for exclusion corrections */
1021 sprintf(err_buf,"wall-ewald-zfac should be >= 2");
1022 CHECK(ir->wall_ewald_zfac < 2);
1025 if (EVDW_SWITCHED(ir->vdwtype)) {
1026 sprintf(err_buf,"With vdwtype = %s rvdw-switch must be < rvdw. Or, better - use a potential modifier.",
1027 evdw_names[ir->vdwtype]);
1028 CHECK(ir->rvdw_switch >= ir->rvdw);
1029 } else if (ir->vdwtype == evdwCUT) {
1030 if (ir->cutoff_scheme == ecutsGROUP && ir->vdw_modifier == eintmodNONE) {
1031 sprintf(err_buf,"With vdwtype = %s, rvdw must be >= rlist unless you use a potential modifier",evdw_names[ir->vdwtype]);
1032 CHECK(ir->rlist > ir->rvdw);
1035 if (ir->cutoff_scheme == ecutsGROUP)
1037 if (EEL_IS_ZERO_AT_CUTOFF(ir->coulombtype)
1038 && (ir->rlistlong <= ir->rcoulomb))
1040 sprintf(warn_buf,"For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rcoulomb.",
1041 IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
1042 warning_note(wi,warn_buf);
1044 if (EVDW_SWITCHED(ir->vdwtype) && (ir->rlistlong <= ir->rvdw))
1046 sprintf(warn_buf,"For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rvdw.",
1047 IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
1048 warning_note(wi,warn_buf);
1052 if (ir->vdwtype == evdwUSER && ir->eDispCorr != edispcNO) {
1053 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.");
1056 if (ir->nstlist == -1) {
1057 sprintf(err_buf,"With nstlist=-1 rvdw and rcoulomb should be smaller than rlist to account for diffusion and possibly charge-group radii");
1058 CHECK(ir->rvdw >= ir->rlist || ir->rcoulomb >= ir->rlist);
1060 sprintf(err_buf,"nstlist can not be smaller than -1");
1061 CHECK(ir->nstlist < -1);
1063 if (ir->eI == eiLBFGS && (ir->coulombtype==eelCUT || ir->vdwtype==evdwCUT)
1065 warning(wi,"For efficient BFGS minimization, use switch/shift/pme instead of cut-off.");
1068 if (ir->eI == eiLBFGS && ir->nbfgscorr <= 0) {
1069 warning(wi,"Using L-BFGS with nbfgscorr<=0 just gets you steepest descent.");
1072 /* ENERGY CONSERVATION */
1073 if (ir_NVE(ir) && ir->cutoff_scheme == ecutsGROUP)
1075 if (!EVDW_MIGHT_BE_ZERO_AT_CUTOFF(ir->vdwtype) && ir->rvdw > 0 && ir->vdw_modifier == eintmodNONE)
1077 sprintf(warn_buf,"You are using a cut-off for VdW interactions with NVE, for good energy conservation use vdwtype = %s (possibly with DispCorr)",
1078 evdw_names[evdwSHIFT]);
1079 warning_note(wi,warn_buf);
1081 if (!EEL_MIGHT_BE_ZERO_AT_CUTOFF(ir->coulombtype) && ir->rcoulomb > 0 && ir->coulomb_modifier == eintmodNONE)
1083 sprintf(warn_buf,"You are using a cut-off for electrostatics with NVE, for good energy conservation use coulombtype = %s or %s",
1084 eel_names[eelPMESWITCH],eel_names[eelRF_ZERO]);
1085 warning_note(wi,warn_buf);
1089 /* IMPLICIT SOLVENT */
1090 if(ir->coulombtype==eelGB_NOTUSED)
1092 ir->coulombtype=eelCUT;
1093 ir->implicit_solvent=eisGBSA;
1094 fprintf(stderr,"Note: Old option for generalized born electrostatics given:\n"
1095 "Changing coulombtype from \"generalized-born\" to \"cut-off\" and instead\n"
1096 "setting implicit-solvent value to \"GBSA\" in input section.\n");
1099 if(ir->sa_algorithm==esaSTILL)
1101 sprintf(err_buf,"Still SA algorithm not available yet, use %s or %s instead\n",esa_names[esaAPPROX],esa_names[esaNO]);
1102 CHECK(ir->sa_algorithm == esaSTILL);
1105 if(ir->implicit_solvent==eisGBSA)
1107 sprintf(err_buf,"With GBSA implicit solvent, rgbradii must be equal to rlist.");
1108 CHECK(ir->rgbradii != ir->rlist);
1110 if(ir->coulombtype!=eelCUT)
1112 sprintf(err_buf,"With GBSA, coulombtype must be equal to %s\n",eel_names[eelCUT]);
1113 CHECK(ir->coulombtype!=eelCUT);
1115 if(ir->vdwtype!=evdwCUT)
1117 sprintf(err_buf,"With GBSA, vdw-type must be equal to %s\n",evdw_names[evdwCUT]);
1118 CHECK(ir->vdwtype!=evdwCUT);
1120 if(ir->nstgbradii<1)
1122 sprintf(warn_buf,"Using GBSA with nstgbradii<1, setting nstgbradii=1");
1123 warning_note(wi,warn_buf);
1126 if(ir->sa_algorithm==esaNO)
1128 sprintf(warn_buf,"No SA (non-polar) calculation requested together with GB. Are you sure this is what you want?\n");
1129 warning_note(wi,warn_buf);
1131 if(ir->sa_surface_tension<0 && ir->sa_algorithm!=esaNO)
1133 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");
1134 warning_note(wi,warn_buf);
1136 if(ir->gb_algorithm==egbSTILL)
1138 ir->sa_surface_tension = 0.0049 * CAL2JOULE * 100;
1142 ir->sa_surface_tension = 0.0054 * CAL2JOULE * 100;
1145 if(ir->sa_surface_tension==0 && ir->sa_algorithm!=esaNO)
1147 sprintf(err_buf, "Surface tension set to 0 while SA-calculation requested\n");
1148 CHECK(ir->sa_surface_tension==0 && ir->sa_algorithm!=esaNO);
1155 if (ir->cutoff_scheme != ecutsGROUP)
1157 warning_error(wi,"AdresS simulation supports only cutoff-scheme=group");
1161 warning_error(wi,"AdresS simulation supports only stochastic dynamics");
1163 if (ir->epc != epcNO)
1165 warning_error(wi,"AdresS simulation does not support pressure coupling");
1167 if (EEL_FULL(ir->coulombtype))
1169 warning_error(wi,"AdresS simulation does not support long-range electrostatics");
1174 /* count the number of text elemets separated by whitespace in a string.
1175 str = the input string
1176 maxptr = the maximum number of allowed elements
1177 ptr = the output array of pointers to the first character of each element
1178 returns: the number of elements. */
1179 int str_nelem(const char *str,int maxptr,char *ptr[])
1187 while (*copy != '\0') {
1189 gmx_fatal(FARGS,"Too many groups on line: '%s' (max is %d)",
1194 while ((*copy != '\0') && !isspace(*copy))
1196 if (*copy != '\0') {
1208 /* interpret a number of doubles from a string and put them in an array,
1209 after allocating space for them.
1210 str = the input string
1211 n = the (pre-allocated) number of doubles read
1212 r = the output array of doubles. */
1213 static void parse_n_real(char *str,int *n,real **r)
1218 *n = str_nelem(str,MAXPTR,ptr);
1221 for(i=0; i<*n; i++) {
1222 (*r)[i] = strtod(ptr[i],NULL);
1226 static void do_fep_params(t_inputrec *ir, char fep_lambda[][STRLEN],char weights[STRLEN]) {
1228 int i,j,max_n_lambda,nweights,nfep[efptNR];
1229 t_lambda *fep = ir->fepvals;
1230 t_expanded *expand = ir->expandedvals;
1231 real **count_fep_lambdas;
1232 gmx_bool bOneLambda = TRUE;
1234 snew(count_fep_lambdas,efptNR);
1236 /* FEP input processing */
1237 /* first, identify the number of lambda values for each type.
1238 All that are nonzero must have the same number */
1240 for (i=0;i<efptNR;i++)
1242 parse_n_real(fep_lambda[i],&(nfep[i]),&(count_fep_lambdas[i]));
1245 /* now, determine the number of components. All must be either zero, or equal. */
1248 for (i=0;i<efptNR;i++)
1250 if (nfep[i] > max_n_lambda) {
1251 max_n_lambda = nfep[i]; /* here's a nonzero one. All of them
1252 must have the same number if its not zero.*/
1257 for (i=0;i<efptNR;i++)
1261 ir->fepvals->separate_dvdl[i] = FALSE;
1263 else if (nfep[i] == max_n_lambda)
1265 if (i!=efptTEMPERATURE) /* we treat this differently -- not really a reason to compute the derivative with
1266 respect to the temperature currently */
1268 ir->fepvals->separate_dvdl[i] = TRUE;
1273 gmx_fatal(FARGS,"Number of lambdas (%d) for FEP type %s not equal to number of other types (%d)",
1274 nfep[i],efpt_names[i],max_n_lambda);
1277 /* we don't print out dhdl if the temperature is changing, since we can't correctly define dhdl in this case */
1278 ir->fepvals->separate_dvdl[efptTEMPERATURE] = FALSE;
1280 /* the number of lambdas is the number we've read in, which is either zero
1281 or the same for all */
1282 fep->n_lambda = max_n_lambda;
1284 /* allocate space for the array of lambda values */
1285 snew(fep->all_lambda,efptNR);
1286 /* if init_lambda is defined, we need to set lambda */
1287 if ((fep->init_lambda > 0) && (fep->n_lambda == 0))
1289 ir->fepvals->separate_dvdl[efptFEP] = TRUE;
1291 /* otherwise allocate the space for all of the lambdas, and transfer the data */
1292 for (i=0;i<efptNR;i++)
1294 snew(fep->all_lambda[i],fep->n_lambda);
1295 if (nfep[i] > 0) /* if it's zero, then the count_fep_lambda arrays
1298 for (j=0;j<fep->n_lambda;j++)
1300 fep->all_lambda[i][j] = (double)count_fep_lambdas[i][j];
1302 sfree(count_fep_lambdas[i]);
1305 sfree(count_fep_lambdas);
1307 /* "fep-vals" is either zero or the full number. If zero, we'll need to define fep-lambdas for internal
1308 bookkeeping -- for now, init_lambda */
1310 if ((nfep[efptFEP] == 0) && (fep->init_lambda >= 0) && (fep->init_lambda <= 1))
1312 for (i=0;i<fep->n_lambda;i++)
1314 fep->all_lambda[efptFEP][i] = fep->init_lambda;
1318 /* check to see if only a single component lambda is defined, and soft core is defined.
1319 In this case, turn on coulomb soft core */
1321 if (max_n_lambda == 0)
1327 for (i=0;i<efptNR;i++)
1329 if ((nfep[i] != 0) && (i!=efptFEP))
1335 if ((bOneLambda) && (fep->sc_alpha > 0))
1337 fep->bScCoul = TRUE;
1340 /* Fill in the others with the efptFEP if they are not explicitly
1341 specified (i.e. nfep[i] == 0). This means if fep is not defined,
1342 they are all zero. */
1344 for (i=0;i<efptNR;i++)
1346 if ((nfep[i] == 0) && (i!=efptFEP))
1348 for (j=0;j<fep->n_lambda;j++)
1350 fep->all_lambda[i][j] = fep->all_lambda[efptFEP][j];
1356 /* make it easier if sc_r_power = 48 by increasing it to the 4th power, to be in the right scale. */
1357 if (fep->sc_r_power == 48)
1359 if (fep->sc_alpha > 0.1)
1361 gmx_fatal(FARGS,"sc_alpha (%f) for sc_r_power = 48 should usually be between 0.001 and 0.004", fep->sc_alpha);
1365 expand = ir->expandedvals;
1366 /* now read in the weights */
1367 parse_n_real(weights,&nweights,&(expand->init_lambda_weights));
1370 expand->bInit_weights = FALSE;
1371 snew(expand->init_lambda_weights,fep->n_lambda); /* initialize to zero */
1373 else if (nweights != fep->n_lambda)
1375 gmx_fatal(FARGS,"Number of weights (%d) is not equal to number of lambda values (%d)",
1376 nweights,fep->n_lambda);
1380 expand->bInit_weights = TRUE;
1382 if ((expand->nstexpanded < 0) && (ir->efep != efepNO)) {
1383 expand->nstexpanded = fep->nstdhdl;
1384 /* if you don't specify nstexpanded when doing expanded ensemble free energy calcs, it is set to nstdhdl */
1386 if ((expand->nstexpanded < 0) && ir->bSimTemp) {
1387 expand->nstexpanded = 2*(int)(ir->opts.tau_t[0]/ir->delta_t);
1388 /* if you don't specify nstexpanded when doing expanded ensemble simulated tempering, it is set to
1389 2*tau_t just to be careful so it's not to frequent */
1394 static void do_simtemp_params(t_inputrec *ir) {
1396 snew(ir->simtempvals->temperatures,ir->fepvals->n_lambda);
1397 GetSimTemps(ir->fepvals->n_lambda,ir->simtempvals,ir->fepvals->all_lambda[efptTEMPERATURE]);
1402 static void do_wall_params(t_inputrec *ir,
1403 char *wall_atomtype, char *wall_density,
1407 char *names[MAXPTR];
1410 opts->wall_atomtype[0] = NULL;
1411 opts->wall_atomtype[1] = NULL;
1413 ir->wall_atomtype[0] = -1;
1414 ir->wall_atomtype[1] = -1;
1415 ir->wall_density[0] = 0;
1416 ir->wall_density[1] = 0;
1420 nstr = str_nelem(wall_atomtype,MAXPTR,names);
1421 if (nstr != ir->nwall)
1423 gmx_fatal(FARGS,"Expected %d elements for wall_atomtype, found %d",
1426 for(i=0; i<ir->nwall; i++)
1428 opts->wall_atomtype[i] = strdup(names[i]);
1431 if (ir->wall_type == ewt93 || ir->wall_type == ewt104) {
1432 nstr = str_nelem(wall_density,MAXPTR,names);
1433 if (nstr != ir->nwall)
1435 gmx_fatal(FARGS,"Expected %d elements for wall-density, found %d",ir->nwall,nstr);
1437 for(i=0; i<ir->nwall; i++)
1439 sscanf(names[i],"%lf",&dbl);
1442 gmx_fatal(FARGS,"wall-density[%d] = %f\n",i,dbl);
1444 ir->wall_density[i] = dbl;
1450 static void add_wall_energrps(gmx_groups_t *groups,int nwall,t_symtab *symtab)
1457 srenew(groups->grpname,groups->ngrpname+nwall);
1458 grps = &(groups->grps[egcENER]);
1459 srenew(grps->nm_ind,grps->nr+nwall);
1460 for(i=0; i<nwall; i++) {
1461 sprintf(str,"wall%d",i);
1462 groups->grpname[groups->ngrpname] = put_symtab(symtab,str);
1463 grps->nm_ind[grps->nr++] = groups->ngrpname++;
1468 void read_expandedparams(int *ninp_p,t_inpfile **inp_p,
1469 t_expanded *expand,warninp_t wi)
1477 /* read expanded ensemble parameters */
1478 CCTYPE ("expanded ensemble variables");
1479 ITYPE ("nstexpanded",expand->nstexpanded,-1);
1480 EETYPE("lmc-stats", expand->elamstats, elamstats_names);
1481 EETYPE("lmc-move", expand->elmcmove, elmcmove_names);
1482 EETYPE("lmc-weights-equil",expand->elmceq,elmceq_names);
1483 ITYPE ("weight-equil-number-all-lambda",expand->equil_n_at_lam,-1);
1484 ITYPE ("weight-equil-number-samples",expand->equil_samples,-1);
1485 ITYPE ("weight-equil-number-steps",expand->equil_steps,-1);
1486 RTYPE ("weight-equil-wl-delta",expand->equil_wl_delta,-1);
1487 RTYPE ("weight-equil-count-ratio",expand->equil_ratio,-1);
1488 CCTYPE("Seed for Monte Carlo in lambda space");
1489 ITYPE ("lmc-seed",expand->lmc_seed,-1);
1490 RTYPE ("mc-temperature",expand->mc_temp,-1);
1491 ITYPE ("lmc-repeats",expand->lmc_repeats,1);
1492 ITYPE ("lmc-gibbsdelta",expand->gibbsdeltalam,-1);
1493 ITYPE ("lmc-forced-nstart",expand->lmc_forced_nstart,0);
1494 EETYPE("symmetrized-transition-matrix", expand->bSymmetrizedTMatrix, yesno_names);
1495 ITYPE("nst-transition-matrix", expand->nstTij, -1);
1496 ITYPE ("mininum-var-min",expand->minvarmin, 100); /*default is reasonable */
1497 ITYPE ("weight-c-range",expand->c_range, 0); /* default is just C=0 */
1498 RTYPE ("wl-scale",expand->wl_scale,0.8);
1499 RTYPE ("wl-ratio",expand->wl_ratio,0.8);
1500 RTYPE ("init-wl-delta",expand->init_wl_delta,1.0);
1501 EETYPE("wl-oneovert",expand->bWLoneovert,yesno_names);
1509 void get_ir(const char *mdparin,const char *mdparout,
1510 t_inputrec *ir,t_gromppopts *opts,
1514 double dumdub[2][6];
1518 char warn_buf[STRLEN];
1519 t_lambda *fep = ir->fepvals;
1520 t_expanded *expand = ir->expandedvals;
1522 inp = read_inpfile(mdparin, &ninp, NULL, wi);
1524 snew(dumstr[0],STRLEN);
1525 snew(dumstr[1],STRLEN);
1527 /* remove the following deprecated commands */
1530 REM_TYPE("domain-decomposition");
1531 REM_TYPE("andersen-seed");
1533 REM_TYPE("dihre-fc");
1534 REM_TYPE("dihre-tau");
1535 REM_TYPE("nstdihreout");
1536 REM_TYPE("nstcheckpoint");
1538 /* replace the following commands with the clearer new versions*/
1539 REPL_TYPE("unconstrained-start","continuation");
1540 REPL_TYPE("foreign-lambda","fep-lambdas");
1542 CCTYPE ("VARIOUS PREPROCESSING OPTIONS");
1543 CTYPE ("Preprocessor information: use cpp syntax.");
1544 CTYPE ("e.g.: -I/home/joe/doe -I/home/mary/roe");
1545 STYPE ("include", opts->include, NULL);
1546 CTYPE ("e.g.: -DPOSRES -DFLEXIBLE (note these variable names are case sensitive)");
1547 STYPE ("define", opts->define, NULL);
1549 CCTYPE ("RUN CONTROL PARAMETERS");
1550 EETYPE("integrator", ir->eI, ei_names);
1551 CTYPE ("Start time and timestep in ps");
1552 RTYPE ("tinit", ir->init_t, 0.0);
1553 RTYPE ("dt", ir->delta_t, 0.001);
1554 STEPTYPE ("nsteps", ir->nsteps, 0);
1555 CTYPE ("For exact run continuation or redoing part of a run");
1556 STEPTYPE ("init-step",ir->init_step, 0);
1557 CTYPE ("Part index is updated automatically on checkpointing (keeps files separate)");
1558 ITYPE ("simulation-part", ir->simulation_part, 1);
1559 CTYPE ("mode for center of mass motion removal");
1560 EETYPE("comm-mode", ir->comm_mode, ecm_names);
1561 CTYPE ("number of steps for center of mass motion removal");
1562 ITYPE ("nstcomm", ir->nstcomm, 100);
1563 CTYPE ("group(s) for center of mass motion removal");
1564 STYPE ("comm-grps", vcm, NULL);
1566 CCTYPE ("LANGEVIN DYNAMICS OPTIONS");
1567 CTYPE ("Friction coefficient (amu/ps) and random seed");
1568 RTYPE ("bd-fric", ir->bd_fric, 0.0);
1569 ITYPE ("ld-seed", ir->ld_seed, 1993);
1572 CCTYPE ("ENERGY MINIMIZATION OPTIONS");
1573 CTYPE ("Force tolerance and initial step-size");
1574 RTYPE ("emtol", ir->em_tol, 10.0);
1575 RTYPE ("emstep", ir->em_stepsize,0.01);
1576 CTYPE ("Max number of iterations in relax-shells");
1577 ITYPE ("niter", ir->niter, 20);
1578 CTYPE ("Step size (ps^2) for minimization of flexible constraints");
1579 RTYPE ("fcstep", ir->fc_stepsize, 0);
1580 CTYPE ("Frequency of steepest descents steps when doing CG");
1581 ITYPE ("nstcgsteep", ir->nstcgsteep, 1000);
1582 ITYPE ("nbfgscorr", ir->nbfgscorr, 10);
1584 CCTYPE ("TEST PARTICLE INSERTION OPTIONS");
1585 RTYPE ("rtpi", ir->rtpi, 0.05);
1587 /* Output options */
1588 CCTYPE ("OUTPUT CONTROL OPTIONS");
1589 CTYPE ("Output frequency for coords (x), velocities (v) and forces (f)");
1590 ITYPE ("nstxout", ir->nstxout, 0);
1591 ITYPE ("nstvout", ir->nstvout, 0);
1592 ITYPE ("nstfout", ir->nstfout, 0);
1593 ir->nstcheckpoint = 1000;
1594 CTYPE ("Output frequency for energies to log file and energy file");
1595 ITYPE ("nstlog", ir->nstlog, 1000);
1596 ITYPE ("nstcalcenergy",ir->nstcalcenergy, 100);
1597 ITYPE ("nstenergy", ir->nstenergy, 1000);
1598 CTYPE ("Output frequency and precision for .xtc file");
1599 ITYPE ("nstxtcout", ir->nstxtcout, 0);
1600 RTYPE ("xtc-precision",ir->xtcprec, 1000.0);
1601 CTYPE ("This selects the subset of atoms for the .xtc file. You can");
1602 CTYPE ("select multiple groups. By default all atoms will be written.");
1603 STYPE ("xtc-grps", xtc_grps, NULL);
1604 CTYPE ("Selection of energy groups");
1605 STYPE ("energygrps", energy, NULL);
1607 /* Neighbor searching */
1608 CCTYPE ("NEIGHBORSEARCHING PARAMETERS");
1609 CTYPE ("cut-off scheme (group: using charge groups, Verlet: particle based cut-offs)");
1610 EETYPE("cutoff-scheme", ir->cutoff_scheme, ecutscheme_names);
1611 CTYPE ("nblist update frequency");
1612 ITYPE ("nstlist", ir->nstlist, 10);
1613 CTYPE ("ns algorithm (simple or grid)");
1614 EETYPE("ns-type", ir->ns_type, ens_names);
1615 /* set ndelta to the optimal value of 2 */
1617 CTYPE ("Periodic boundary conditions: xyz, no, xy");
1618 EETYPE("pbc", ir->ePBC, epbc_names);
1619 EETYPE("periodic-molecules", ir->bPeriodicMols, yesno_names);
1620 CTYPE ("Allowed energy drift due to the Verlet buffer in kJ/mol/ps per atom,");
1621 CTYPE ("a value of -1 means: use rlist");
1622 RTYPE("verlet-buffer-drift", ir->verletbuf_drift, 0.005);
1623 CTYPE ("nblist cut-off");
1624 RTYPE ("rlist", ir->rlist, -1);
1625 CTYPE ("long-range cut-off for switched potentials");
1626 RTYPE ("rlistlong", ir->rlistlong, -1);
1627 ITYPE ("nstcalclr", ir->nstcalclr, -1);
1629 /* Electrostatics */
1630 CCTYPE ("OPTIONS FOR ELECTROSTATICS AND VDW");
1631 CTYPE ("Method for doing electrostatics");
1632 EETYPE("coulombtype", ir->coulombtype, eel_names);
1633 EETYPE("coulomb-modifier", ir->coulomb_modifier, eintmod_names);
1634 CTYPE ("cut-off lengths");
1635 RTYPE ("rcoulomb-switch", ir->rcoulomb_switch, 0.0);
1636 RTYPE ("rcoulomb", ir->rcoulomb, -1);
1637 CTYPE ("Relative dielectric constant for the medium and the reaction field");
1638 RTYPE ("epsilon-r", ir->epsilon_r, 1.0);
1639 RTYPE ("epsilon-rf", ir->epsilon_rf, 0.0);
1640 CTYPE ("Method for doing Van der Waals");
1641 EETYPE("vdw-type", ir->vdwtype, evdw_names);
1642 EETYPE("vdw-modifier", ir->vdw_modifier, eintmod_names);
1643 CTYPE ("cut-off lengths");
1644 RTYPE ("rvdw-switch", ir->rvdw_switch, 0.0);
1645 RTYPE ("rvdw", ir->rvdw, -1);
1646 CTYPE ("Apply long range dispersion corrections for Energy and Pressure");
1647 EETYPE("DispCorr", ir->eDispCorr, edispc_names);
1648 CTYPE ("Extension of the potential lookup tables beyond the cut-off");
1649 RTYPE ("table-extension", ir->tabext, 1.0);
1650 CTYPE ("Seperate tables between energy group pairs");
1651 STYPE ("energygrp-table", egptable, NULL);
1652 CTYPE ("Spacing for the PME/PPPM FFT grid");
1653 RTYPE ("fourierspacing", ir->fourier_spacing,0.12);
1654 CTYPE ("FFT grid size, when a value is 0 fourierspacing will be used");
1655 ITYPE ("fourier-nx", ir->nkx, 0);
1656 ITYPE ("fourier-ny", ir->nky, 0);
1657 ITYPE ("fourier-nz", ir->nkz, 0);
1658 CTYPE ("EWALD/PME/PPPM parameters");
1659 ITYPE ("pme-order", ir->pme_order, 4);
1660 RTYPE ("ewald-rtol", ir->ewald_rtol, 0.00001);
1661 EETYPE("ewald-geometry", ir->ewald_geometry, eewg_names);
1662 RTYPE ("epsilon-surface", ir->epsilon_surface, 0.0);
1663 EETYPE("optimize-fft",ir->bOptFFT, yesno_names);
1665 CCTYPE("IMPLICIT SOLVENT ALGORITHM");
1666 EETYPE("implicit-solvent", ir->implicit_solvent, eis_names);
1668 CCTYPE ("GENERALIZED BORN ELECTROSTATICS");
1669 CTYPE ("Algorithm for calculating Born radii");
1670 EETYPE("gb-algorithm", ir->gb_algorithm, egb_names);
1671 CTYPE ("Frequency of calculating the Born radii inside rlist");
1672 ITYPE ("nstgbradii", ir->nstgbradii, 1);
1673 CTYPE ("Cutoff for Born radii calculation; the contribution from atoms");
1674 CTYPE ("between rlist and rgbradii is updated every nstlist steps");
1675 RTYPE ("rgbradii", ir->rgbradii, 1.0);
1676 CTYPE ("Dielectric coefficient of the implicit solvent");
1677 RTYPE ("gb-epsilon-solvent",ir->gb_epsilon_solvent, 80.0);
1678 CTYPE ("Salt concentration in M for Generalized Born models");
1679 RTYPE ("gb-saltconc", ir->gb_saltconc, 0.0);
1680 CTYPE ("Scaling factors used in the OBC GB model. Default values are OBC(II)");
1681 RTYPE ("gb-obc-alpha", ir->gb_obc_alpha, 1.0);
1682 RTYPE ("gb-obc-beta", ir->gb_obc_beta, 0.8);
1683 RTYPE ("gb-obc-gamma", ir->gb_obc_gamma, 4.85);
1684 RTYPE ("gb-dielectric-offset", ir->gb_dielectric_offset, 0.009);
1685 EETYPE("sa-algorithm", ir->sa_algorithm, esa_names);
1686 CTYPE ("Surface tension (kJ/mol/nm^2) for the SA (nonpolar surface) part of GBSA");
1687 CTYPE ("The value -1 will set default value for Still/HCT/OBC GB-models.");
1688 RTYPE ("sa-surface-tension", ir->sa_surface_tension, -1);
1690 /* Coupling stuff */
1691 CCTYPE ("OPTIONS FOR WEAK COUPLING ALGORITHMS");
1692 CTYPE ("Temperature coupling");
1693 EETYPE("tcoupl", ir->etc, etcoupl_names);
1694 ITYPE ("nsttcouple", ir->nsttcouple, -1);
1695 ITYPE("nh-chain-length", ir->opts.nhchainlength, NHCHAINLENGTH);
1696 EETYPE("print-nose-hoover-chain-variables", ir->bPrintNHChains, yesno_names);
1697 CTYPE ("Groups to couple separately");
1698 STYPE ("tc-grps", tcgrps, NULL);
1699 CTYPE ("Time constant (ps) and reference temperature (K)");
1700 STYPE ("tau-t", tau_t, NULL);
1701 STYPE ("ref-t", ref_t, NULL);
1702 CTYPE ("pressure coupling");
1703 EETYPE("pcoupl", ir->epc, epcoupl_names);
1704 EETYPE("pcoupltype", ir->epct, epcoupltype_names);
1705 ITYPE ("nstpcouple", ir->nstpcouple, -1);
1706 CTYPE ("Time constant (ps), compressibility (1/bar) and reference P (bar)");
1707 RTYPE ("tau-p", ir->tau_p, 1.0);
1708 STYPE ("compressibility", dumstr[0], NULL);
1709 STYPE ("ref-p", dumstr[1], NULL);
1710 CTYPE ("Scaling of reference coordinates, No, All or COM");
1711 EETYPE ("refcoord-scaling",ir->refcoord_scaling,erefscaling_names);
1714 CCTYPE ("OPTIONS FOR QMMM calculations");
1715 EETYPE("QMMM", ir->bQMMM, yesno_names);
1716 CTYPE ("Groups treated Quantum Mechanically");
1717 STYPE ("QMMM-grps", QMMM, NULL);
1718 CTYPE ("QM method");
1719 STYPE("QMmethod", QMmethod, NULL);
1720 CTYPE ("QMMM scheme");
1721 EETYPE("QMMMscheme", ir->QMMMscheme, eQMMMscheme_names);
1722 CTYPE ("QM basisset");
1723 STYPE("QMbasis", QMbasis, NULL);
1724 CTYPE ("QM charge");
1725 STYPE ("QMcharge", QMcharge,NULL);
1726 CTYPE ("QM multiplicity");
1727 STYPE ("QMmult", QMmult,NULL);
1728 CTYPE ("Surface Hopping");
1729 STYPE ("SH", bSH, NULL);
1730 CTYPE ("CAS space options");
1731 STYPE ("CASorbitals", CASorbitals, NULL);
1732 STYPE ("CASelectrons", CASelectrons, NULL);
1733 STYPE ("SAon", SAon, NULL);
1734 STYPE ("SAoff",SAoff,NULL);
1735 STYPE ("SAsteps", SAsteps, NULL);
1736 CTYPE ("Scale factor for MM charges");
1737 RTYPE ("MMChargeScaleFactor", ir->scalefactor, 1.0);
1738 CTYPE ("Optimization of QM subsystem");
1739 STYPE ("bOPT", bOPT, NULL);
1740 STYPE ("bTS", bTS, NULL);
1742 /* Simulated annealing */
1743 CCTYPE("SIMULATED ANNEALING");
1744 CTYPE ("Type of annealing for each temperature group (no/single/periodic)");
1745 STYPE ("annealing", anneal, NULL);
1746 CTYPE ("Number of time points to use for specifying annealing in each group");
1747 STYPE ("annealing-npoints", anneal_npoints, NULL);
1748 CTYPE ("List of times at the annealing points for each group");
1749 STYPE ("annealing-time", anneal_time, NULL);
1750 CTYPE ("Temp. at each annealing point, for each group.");
1751 STYPE ("annealing-temp", anneal_temp, NULL);
1754 CCTYPE ("GENERATE VELOCITIES FOR STARTUP RUN");
1755 EETYPE("gen-vel", opts->bGenVel, yesno_names);
1756 RTYPE ("gen-temp", opts->tempi, 300.0);
1757 ITYPE ("gen-seed", opts->seed, 173529);
1760 CCTYPE ("OPTIONS FOR BONDS");
1761 EETYPE("constraints", opts->nshake, constraints);
1762 CTYPE ("Type of constraint algorithm");
1763 EETYPE("constraint-algorithm", ir->eConstrAlg, econstr_names);
1764 CTYPE ("Do not constrain the start configuration");
1765 EETYPE("continuation", ir->bContinuation, yesno_names);
1766 CTYPE ("Use successive overrelaxation to reduce the number of shake iterations");
1767 EETYPE("Shake-SOR", ir->bShakeSOR, yesno_names);
1768 CTYPE ("Relative tolerance of shake");
1769 RTYPE ("shake-tol", ir->shake_tol, 0.0001);
1770 CTYPE ("Highest order in the expansion of the constraint coupling matrix");
1771 ITYPE ("lincs-order", ir->nProjOrder, 4);
1772 CTYPE ("Number of iterations in the final step of LINCS. 1 is fine for");
1773 CTYPE ("normal simulations, but use 2 to conserve energy in NVE runs.");
1774 CTYPE ("For energy minimization with constraints it should be 4 to 8.");
1775 ITYPE ("lincs-iter", ir->nLincsIter, 1);
1776 CTYPE ("Lincs will write a warning to the stderr if in one step a bond");
1777 CTYPE ("rotates over more degrees than");
1778 RTYPE ("lincs-warnangle", ir->LincsWarnAngle, 30.0);
1779 CTYPE ("Convert harmonic bonds to morse potentials");
1780 EETYPE("morse", opts->bMorse,yesno_names);
1782 /* Energy group exclusions */
1783 CCTYPE ("ENERGY GROUP EXCLUSIONS");
1784 CTYPE ("Pairs of energy groups for which all non-bonded interactions are excluded");
1785 STYPE ("energygrp-excl", egpexcl, NULL);
1789 CTYPE ("Number of walls, type, atom types, densities and box-z scale factor for Ewald");
1790 ITYPE ("nwall", ir->nwall, 0);
1791 EETYPE("wall-type", ir->wall_type, ewt_names);
1792 RTYPE ("wall-r-linpot", ir->wall_r_linpot, -1);
1793 STYPE ("wall-atomtype", wall_atomtype, NULL);
1794 STYPE ("wall-density", wall_density, NULL);
1795 RTYPE ("wall-ewald-zfac", ir->wall_ewald_zfac, 3);
1798 CCTYPE("COM PULLING");
1799 CTYPE("Pull type: no, umbrella, constraint or constant-force");
1800 EETYPE("pull", ir->ePull, epull_names);
1801 if (ir->ePull != epullNO) {
1803 pull_grp = read_pullparams(&ninp,&inp,ir->pull,&opts->pull_start,wi);
1806 /* Enforced rotation */
1807 CCTYPE("ENFORCED ROTATION");
1808 CTYPE("Enforced rotation: No or Yes");
1809 EETYPE("rotation", ir->bRot, yesno_names);
1812 rot_grp = read_rotparams(&ninp,&inp,ir->rot,wi);
1816 CCTYPE("NMR refinement stuff");
1817 CTYPE ("Distance restraints type: No, Simple or Ensemble");
1818 EETYPE("disre", ir->eDisre, edisre_names);
1819 CTYPE ("Force weighting of pairs in one distance restraint: Conservative or Equal");
1820 EETYPE("disre-weighting", ir->eDisreWeighting, edisreweighting_names);
1821 CTYPE ("Use sqrt of the time averaged times the instantaneous violation");
1822 EETYPE("disre-mixed", ir->bDisreMixed, yesno_names);
1823 RTYPE ("disre-fc", ir->dr_fc, 1000.0);
1824 RTYPE ("disre-tau", ir->dr_tau, 0.0);
1825 CTYPE ("Output frequency for pair distances to energy file");
1826 ITYPE ("nstdisreout", ir->nstdisreout, 100);
1827 CTYPE ("Orientation restraints: No or Yes");
1828 EETYPE("orire", opts->bOrire, yesno_names);
1829 CTYPE ("Orientation restraints force constant and tau for time averaging");
1830 RTYPE ("orire-fc", ir->orires_fc, 0.0);
1831 RTYPE ("orire-tau", ir->orires_tau, 0.0);
1832 STYPE ("orire-fitgrp",orirefitgrp, NULL);
1833 CTYPE ("Output frequency for trace(SD) and S to energy file");
1834 ITYPE ("nstorireout", ir->nstorireout, 100);
1836 /* free energy variables */
1837 CCTYPE ("Free energy variables");
1838 EETYPE("free-energy", ir->efep, efep_names);
1839 STYPE ("couple-moltype", couple_moltype, NULL);
1840 EETYPE("couple-lambda0", opts->couple_lam0, couple_lam);
1841 EETYPE("couple-lambda1", opts->couple_lam1, couple_lam);
1842 EETYPE("couple-intramol", opts->bCoupleIntra, yesno_names);
1844 RTYPE ("init-lambda", fep->init_lambda,-1); /* start with -1 so
1846 it was not entered */
1847 ITYPE ("init-lambda-state", fep->init_fep_state,0);
1848 RTYPE ("delta-lambda",fep->delta_lambda,0.0);
1849 ITYPE ("nstdhdl",fep->nstdhdl, 100);
1850 STYPE ("fep-lambdas", fep_lambda[efptFEP], NULL);
1851 STYPE ("mass-lambdas", fep_lambda[efptMASS], NULL);
1852 STYPE ("coul-lambdas", fep_lambda[efptCOUL], NULL);
1853 STYPE ("vdw-lambdas", fep_lambda[efptVDW], NULL);
1854 STYPE ("bonded-lambdas", fep_lambda[efptBONDED], NULL);
1855 STYPE ("restraint-lambdas", fep_lambda[efptRESTRAINT], NULL);
1856 STYPE ("temperature-lambdas", fep_lambda[efptTEMPERATURE], NULL);
1857 STYPE ("init-lambda-weights",lambda_weights,NULL);
1858 EETYPE("dhdl-print-energy", fep->bPrintEnergy, yesno_names);
1859 RTYPE ("sc-alpha",fep->sc_alpha,0.0);
1860 ITYPE ("sc-power",fep->sc_power,1);
1861 RTYPE ("sc-r-power",fep->sc_r_power,6.0);
1862 RTYPE ("sc-sigma",fep->sc_sigma,0.3);
1863 EETYPE("sc-coul",fep->bScCoul,yesno_names);
1864 ITYPE ("dh_hist_size", fep->dh_hist_size, 0);
1865 RTYPE ("dh_hist_spacing", fep->dh_hist_spacing, 0.1);
1866 EETYPE("separate-dhdl-file", fep->separate_dhdl_file,
1867 separate_dhdl_file_names);
1868 EETYPE("dhdl-derivatives", fep->dhdl_derivatives, dhdl_derivatives_names);
1869 ITYPE ("dh_hist_size", fep->dh_hist_size, 0);
1870 RTYPE ("dh_hist_spacing", fep->dh_hist_spacing, 0.1);
1872 /* Non-equilibrium MD stuff */
1873 CCTYPE("Non-equilibrium MD stuff");
1874 STYPE ("acc-grps", accgrps, NULL);
1875 STYPE ("accelerate", acc, NULL);
1876 STYPE ("freezegrps", freeze, NULL);
1877 STYPE ("freezedim", frdim, NULL);
1878 RTYPE ("cos-acceleration", ir->cos_accel, 0);
1879 STYPE ("deform", deform, NULL);
1881 /* simulated tempering variables */
1882 CCTYPE("simulated tempering variables");
1883 EETYPE("simulated-tempering",ir->bSimTemp,yesno_names);
1884 EETYPE("simulated-tempering-scaling",ir->simtempvals->eSimTempScale,esimtemp_names);
1885 RTYPE("sim-temp-low",ir->simtempvals->simtemp_low,300.0);
1886 RTYPE("sim-temp-high",ir->simtempvals->simtemp_high,300.0);
1888 /* expanded ensemble variables */
1889 if (ir->efep==efepEXPANDED || ir->bSimTemp)
1891 read_expandedparams(&ninp,&inp,expand,wi);
1894 /* Electric fields */
1895 CCTYPE("Electric fields");
1896 CTYPE ("Format is number of terms (int) and for all terms an amplitude (real)");
1897 CTYPE ("and a phase angle (real)");
1898 STYPE ("E-x", efield_x, NULL);
1899 STYPE ("E-xt", efield_xt, NULL);
1900 STYPE ("E-y", efield_y, NULL);
1901 STYPE ("E-yt", efield_yt, NULL);
1902 STYPE ("E-z", efield_z, NULL);
1903 STYPE ("E-zt", efield_zt, NULL);
1905 /* AdResS defined thingies */
1906 CCTYPE ("AdResS parameters");
1907 EETYPE("adress", ir->bAdress, yesno_names);
1910 read_adressparams(&ninp,&inp,ir->adress,wi);
1913 /* User defined thingies */
1914 CCTYPE ("User defined thingies");
1915 STYPE ("user1-grps", user1, NULL);
1916 STYPE ("user2-grps", user2, NULL);
1917 ITYPE ("userint1", ir->userint1, 0);
1918 ITYPE ("userint2", ir->userint2, 0);
1919 ITYPE ("userint3", ir->userint3, 0);
1920 ITYPE ("userint4", ir->userint4, 0);
1921 RTYPE ("userreal1", ir->userreal1, 0);
1922 RTYPE ("userreal2", ir->userreal2, 0);
1923 RTYPE ("userreal3", ir->userreal3, 0);
1924 RTYPE ("userreal4", ir->userreal4, 0);
1927 write_inpfile(mdparout,ninp,inp,FALSE,wi);
1928 for (i=0; (i<ninp); i++) {
1930 sfree(inp[i].value);
1934 /* Process options if necessary */
1935 for(m=0; m<2; m++) {
1936 for(i=0; i<2*DIM; i++)
1941 if (sscanf(dumstr[m],"%lf",&(dumdub[m][XX]))!=1) {
1942 warning_error(wi,"Pressure coupling not enough values (I need 1)");
1944 dumdub[m][YY]=dumdub[m][ZZ]=dumdub[m][XX];
1946 case epctSEMIISOTROPIC:
1947 case epctSURFACETENSION:
1948 if (sscanf(dumstr[m],"%lf%lf",
1949 &(dumdub[m][XX]),&(dumdub[m][ZZ]))!=2) {
1950 warning_error(wi,"Pressure coupling not enough values (I need 2)");
1952 dumdub[m][YY]=dumdub[m][XX];
1954 case epctANISOTROPIC:
1955 if (sscanf(dumstr[m],"%lf%lf%lf%lf%lf%lf",
1956 &(dumdub[m][XX]),&(dumdub[m][YY]),&(dumdub[m][ZZ]),
1957 &(dumdub[m][3]),&(dumdub[m][4]),&(dumdub[m][5]))!=6) {
1958 warning_error(wi,"Pressure coupling not enough values (I need 6)");
1962 gmx_fatal(FARGS,"Pressure coupling type %s not implemented yet",
1963 epcoupltype_names[ir->epct]);
1967 clear_mat(ir->ref_p);
1968 clear_mat(ir->compress);
1969 for(i=0; i<DIM; i++) {
1970 ir->ref_p[i][i] = dumdub[1][i];
1971 ir->compress[i][i] = dumdub[0][i];
1973 if (ir->epct == epctANISOTROPIC) {
1974 ir->ref_p[XX][YY] = dumdub[1][3];
1975 ir->ref_p[XX][ZZ] = dumdub[1][4];
1976 ir->ref_p[YY][ZZ] = dumdub[1][5];
1977 if (ir->ref_p[XX][YY]!=0 && ir->ref_p[XX][ZZ]!=0 && ir->ref_p[YY][ZZ]!=0) {
1978 warning(wi,"All off-diagonal reference pressures are non-zero. Are you sure you want to apply a threefold shear stress?\n");
1980 ir->compress[XX][YY] = dumdub[0][3];
1981 ir->compress[XX][ZZ] = dumdub[0][4];
1982 ir->compress[YY][ZZ] = dumdub[0][5];
1983 for(i=0; i<DIM; i++) {
1984 for(m=0; m<i; m++) {
1985 ir->ref_p[i][m] = ir->ref_p[m][i];
1986 ir->compress[i][m] = ir->compress[m][i];
1991 if (ir->comm_mode == ecmNO)
1994 opts->couple_moltype = NULL;
1995 if (strlen(couple_moltype) > 0)
1997 if (ir->efep != efepNO)
1999 opts->couple_moltype = strdup(couple_moltype);
2000 if (opts->couple_lam0 == opts->couple_lam1)
2002 warning(wi,"The lambda=0 and lambda=1 states for coupling are identical");
2004 if (ir->eI == eiMD && (opts->couple_lam0 == ecouplamNONE ||
2005 opts->couple_lam1 == ecouplamNONE))
2007 warning(wi,"For proper sampling of the (nearly) decoupled state, stochastic dynamics should be used");
2012 warning(wi,"Can not couple a molecule with free_energy = no");
2015 /* FREE ENERGY AND EXPANDED ENSEMBLE OPTIONS */
2016 if (ir->efep != efepNO) {
2017 if (fep->delta_lambda > 0) {
2018 ir->efep = efepSLOWGROWTH;
2023 fep->bPrintEnergy = TRUE;
2024 /* always print out the energy to dhdl if we are doing expanded ensemble, since we need the total energy
2025 if the temperature is changing. */
2028 if ((ir->efep != efepNO) || ir->bSimTemp)
2030 ir->bExpanded = FALSE;
2031 if ((ir->efep == efepEXPANDED) || ir->bSimTemp)
2033 ir->bExpanded = TRUE;
2035 do_fep_params(ir,fep_lambda,lambda_weights);
2036 if (ir->bSimTemp) { /* done after fep params */
2037 do_simtemp_params(ir);
2042 ir->fepvals->n_lambda = 0;
2045 /* WALL PARAMETERS */
2047 do_wall_params(ir,wall_atomtype,wall_density,opts);
2049 /* ORIENTATION RESTRAINT PARAMETERS */
2051 if (opts->bOrire && str_nelem(orirefitgrp,MAXPTR,NULL)!=1) {
2052 warning_error(wi,"ERROR: Need one orientation restraint fit group\n");
2055 /* DEFORMATION PARAMETERS */
2057 clear_mat(ir->deform);
2062 m = sscanf(deform,"%lf %lf %lf %lf %lf %lf",
2063 &(dumdub[0][0]),&(dumdub[0][1]),&(dumdub[0][2]),
2064 &(dumdub[0][3]),&(dumdub[0][4]),&(dumdub[0][5]));
2067 ir->deform[i][i] = dumdub[0][i];
2069 ir->deform[YY][XX] = dumdub[0][3];
2070 ir->deform[ZZ][XX] = dumdub[0][4];
2071 ir->deform[ZZ][YY] = dumdub[0][5];
2072 if (ir->epc != epcNO) {
2075 if (ir->deform[i][j]!=0 && ir->compress[i][j]!=0) {
2076 warning_error(wi,"A box element has deform set and compressibility > 0");
2080 if (ir->deform[i][j]!=0) {
2081 for(m=j; m<DIM; m++)
2082 if (ir->compress[m][j]!=0) {
2083 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.");
2084 warning(wi,warn_buf);
2093 static int search_QMstring(char *s,int ng,const char *gn[])
2095 /* same as normal search_string, but this one searches QM strings */
2098 for(i=0; (i<ng); i++)
2099 if (gmx_strcasecmp(s,gn[i]) == 0)
2102 gmx_fatal(FARGS,"this QM method or basisset (%s) is not implemented\n!",s);
2106 } /* search_QMstring */
2109 int search_string(char *s,int ng,char *gn[])
2113 for(i=0; (i<ng); i++)
2115 if (gmx_strcasecmp(s,gn[i]) == 0)
2122 "Group %s referenced in the .mdp file was not found in the index file.\n"
2123 "Group names must match either [moleculetype] names or custom index group\n"
2124 "names, in which case you must supply an index file to the '-n' option\n"
2131 static gmx_bool do_numbering(int natoms,gmx_groups_t *groups,int ng,char *ptrs[],
2132 t_blocka *block,char *gnames[],
2133 int gtype,int restnm,
2134 int grptp,gmx_bool bVerbose,
2137 unsigned short *cbuf;
2138 t_grps *grps=&(groups->grps[gtype]);
2139 int i,j,gid,aj,ognr,ntot=0;
2142 char warn_buf[STRLEN];
2146 fprintf(debug,"Starting numbering %d groups of type %d\n",ng,gtype);
2149 title = gtypes[gtype];
2152 /* Mark all id's as not set */
2153 for(i=0; (i<natoms); i++)
2158 snew(grps->nm_ind,ng+1); /* +1 for possible rest group */
2159 for(i=0; (i<ng); i++)
2161 /* Lookup the group name in the block structure */
2162 gid = search_string(ptrs[i],block->nr,gnames);
2163 if ((grptp != egrptpONE) || (i == 0))
2165 grps->nm_ind[grps->nr++]=gid;
2169 fprintf(debug,"Found gid %d for group %s\n",gid,ptrs[i]);
2172 /* Now go over the atoms in the group */
2173 for(j=block->index[gid]; (j<block->index[gid+1]); j++)
2178 /* Range checking */
2179 if ((aj < 0) || (aj >= natoms))
2181 gmx_fatal(FARGS,"Invalid atom number %d in indexfile",aj);
2183 /* Lookup up the old group number */
2187 gmx_fatal(FARGS,"Atom %d in multiple %s groups (%d and %d)",
2188 aj+1,title,ognr+1,i+1);
2192 /* Store the group number in buffer */
2193 if (grptp == egrptpONE)
2206 /* Now check whether we have done all atoms */
2210 if (grptp == egrptpALL)
2212 gmx_fatal(FARGS,"%d atoms are not part of any of the %s groups",
2215 else if (grptp == egrptpPART)
2217 sprintf(warn_buf,"%d atoms are not part of any of the %s groups",
2219 warning_note(wi,warn_buf);
2221 /* Assign all atoms currently unassigned to a rest group */
2222 for(j=0; (j<natoms); j++)
2224 if (cbuf[j] == NOGID)
2230 if (grptp != egrptpPART)
2235 "Making dummy/rest group for %s containing %d elements\n",
2238 /* Add group name "rest" */
2239 grps->nm_ind[grps->nr] = restnm;
2241 /* Assign the rest name to all atoms not currently assigned to a group */
2242 for(j=0; (j<natoms); j++)
2244 if (cbuf[j] == NOGID)
2253 if (grps->nr == 1 && (ntot == 0 || ntot == natoms))
2255 /* All atoms are part of one (or no) group, no index required */
2256 groups->ngrpnr[gtype] = 0;
2257 groups->grpnr[gtype] = NULL;
2261 groups->ngrpnr[gtype] = natoms;
2262 snew(groups->grpnr[gtype],natoms);
2263 for(j=0; (j<natoms); j++)
2265 groups->grpnr[gtype][j] = cbuf[j];
2271 return (bRest && grptp == egrptpPART);
2274 static void calc_nrdf(gmx_mtop_t *mtop,t_inputrec *ir,char **gnames)
2277 gmx_groups_t *groups;
2279 int natoms,ai,aj,i,j,d,g,imin,jmin,nc;
2281 int *nrdf2,*na_vcm,na_tot;
2282 double *nrdf_tc,*nrdf_vcm,nrdf_uc,n_sub=0;
2283 gmx_mtop_atomloop_all_t aloop;
2285 int mb,mol,ftype,as;
2286 gmx_molblock_t *molb;
2287 gmx_moltype_t *molt;
2290 * First calc 3xnr-atoms for each group
2291 * then subtract half a degree of freedom for each constraint
2293 * Only atoms and nuclei contribute to the degrees of freedom...
2298 groups = &mtop->groups;
2299 natoms = mtop->natoms;
2301 /* Allocate one more for a possible rest group */
2302 /* We need to sum degrees of freedom into doubles,
2303 * since floats give too low nrdf's above 3 million atoms.
2305 snew(nrdf_tc,groups->grps[egcTC].nr+1);
2306 snew(nrdf_vcm,groups->grps[egcVCM].nr+1);
2307 snew(na_vcm,groups->grps[egcVCM].nr+1);
2309 for(i=0; i<groups->grps[egcTC].nr; i++)
2311 for(i=0; i<groups->grps[egcVCM].nr+1; i++)
2315 aloop = gmx_mtop_atomloop_all_init(mtop);
2316 while (gmx_mtop_atomloop_all_next(aloop,&i,&atom)) {
2318 if (atom->ptype == eptAtom || atom->ptype == eptNucleus) {
2319 g = ggrpnr(groups,egcFREEZE,i);
2320 /* Double count nrdf for particle i */
2321 for(d=0; d<DIM; d++) {
2322 if (opts->nFreeze[g][d] == 0) {
2326 nrdf_tc [ggrpnr(groups,egcTC ,i)] += 0.5*nrdf2[i];
2327 nrdf_vcm[ggrpnr(groups,egcVCM,i)] += 0.5*nrdf2[i];
2332 for(mb=0; mb<mtop->nmolblock; mb++) {
2333 molb = &mtop->molblock[mb];
2334 molt = &mtop->moltype[molb->type];
2335 atom = molt->atoms.atom;
2336 for(mol=0; mol<molb->nmol; mol++) {
2337 for (ftype=F_CONSTR; ftype<=F_CONSTRNC; ftype++) {
2338 ia = molt->ilist[ftype].iatoms;
2339 for(i=0; i<molt->ilist[ftype].nr; ) {
2340 /* Subtract degrees of freedom for the constraints,
2341 * if the particles still have degrees of freedom left.
2342 * If one of the particles is a vsite or a shell, then all
2343 * constraint motion will go there, but since they do not
2344 * contribute to the constraints the degrees of freedom do not
2349 if (((atom[ia[1]].ptype == eptNucleus) ||
2350 (atom[ia[1]].ptype == eptAtom)) &&
2351 ((atom[ia[2]].ptype == eptNucleus) ||
2352 (atom[ia[2]].ptype == eptAtom))) {
2361 imin = min(imin,nrdf2[ai]);
2362 jmin = min(jmin,nrdf2[aj]);
2365 nrdf_tc [ggrpnr(groups,egcTC ,ai)] -= 0.5*imin;
2366 nrdf_tc [ggrpnr(groups,egcTC ,aj)] -= 0.5*jmin;
2367 nrdf_vcm[ggrpnr(groups,egcVCM,ai)] -= 0.5*imin;
2368 nrdf_vcm[ggrpnr(groups,egcVCM,aj)] -= 0.5*jmin;
2370 ia += interaction_function[ftype].nratoms+1;
2371 i += interaction_function[ftype].nratoms+1;
2374 ia = molt->ilist[F_SETTLE].iatoms;
2375 for(i=0; i<molt->ilist[F_SETTLE].nr; ) {
2376 /* Subtract 1 dof from every atom in the SETTLE */
2377 for(j=0; j<3; j++) {
2379 imin = min(2,nrdf2[ai]);
2381 nrdf_tc [ggrpnr(groups,egcTC ,ai)] -= 0.5*imin;
2382 nrdf_vcm[ggrpnr(groups,egcVCM,ai)] -= 0.5*imin;
2387 as += molt->atoms.nr;
2391 if (ir->ePull == epullCONSTRAINT) {
2392 /* Correct nrdf for the COM constraints.
2393 * We correct using the TC and VCM group of the first atom
2394 * in the reference and pull group. If atoms in one pull group
2395 * belong to different TC or VCM groups it is anyhow difficult
2396 * to determine the optimal nrdf assignment.
2399 if (pull->eGeom == epullgPOS) {
2401 for(i=0; i<DIM; i++) {
2408 for(i=0; i<pull->ngrp; i++) {
2410 if (pull->grp[0].nat > 0) {
2411 /* Subtract 1/2 dof from the reference group */
2412 ai = pull->grp[0].ind[0];
2413 if (nrdf_tc[ggrpnr(groups,egcTC,ai)] > 1) {
2414 nrdf_tc [ggrpnr(groups,egcTC ,ai)] -= 0.5;
2415 nrdf_vcm[ggrpnr(groups,egcVCM,ai)] -= 0.5;
2419 /* Subtract 1/2 dof from the pulled group */
2420 ai = pull->grp[1+i].ind[0];
2421 nrdf_tc [ggrpnr(groups,egcTC ,ai)] -= 0.5*imin;
2422 nrdf_vcm[ggrpnr(groups,egcVCM,ai)] -= 0.5*imin;
2423 if (nrdf_tc[ggrpnr(groups,egcTC,ai)] < 0)
2424 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)]]);
2428 if (ir->nstcomm != 0) {
2429 /* Subtract 3 from the number of degrees of freedom in each vcm group
2430 * when com translation is removed and 6 when rotation is removed
2433 switch (ir->comm_mode) {
2435 n_sub = ndof_com(ir);
2442 gmx_incons("Checking comm_mode");
2445 for(i=0; i<groups->grps[egcTC].nr; i++) {
2446 /* Count the number of atoms of TC group i for every VCM group */
2447 for(j=0; j<groups->grps[egcVCM].nr+1; j++)
2450 for(ai=0; ai<natoms; ai++)
2451 if (ggrpnr(groups,egcTC,ai) == i) {
2452 na_vcm[ggrpnr(groups,egcVCM,ai)]++;
2455 /* Correct for VCM removal according to the fraction of each VCM
2456 * group present in this TC group.
2458 nrdf_uc = nrdf_tc[i];
2460 fprintf(debug,"T-group[%d] nrdf_uc = %g, n_sub = %g\n",
2464 for(j=0; j<groups->grps[egcVCM].nr+1; j++) {
2465 if (nrdf_vcm[j] > n_sub) {
2466 nrdf_tc[i] += nrdf_uc*((double)na_vcm[j]/(double)na_tot)*
2467 (nrdf_vcm[j] - n_sub)/nrdf_vcm[j];
2470 fprintf(debug," nrdf_vcm[%d] = %g, nrdf = %g\n",
2471 j,nrdf_vcm[j],nrdf_tc[i]);
2476 for(i=0; (i<groups->grps[egcTC].nr); i++) {
2477 opts->nrdf[i] = nrdf_tc[i];
2478 if (opts->nrdf[i] < 0)
2481 "Number of degrees of freedom in T-Coupling group %s is %.2f\n",
2482 gnames[groups->grps[egcTC].nm_ind[i]],opts->nrdf[i]);
2491 static void decode_cos(char *s,t_cosines *cosine,gmx_bool bTime)
2494 char format[STRLEN],f1[STRLEN];
2505 sscanf(t,"%d",&(cosine->n));
2506 if (cosine->n <= 0) {
2509 snew(cosine->a,cosine->n);
2510 snew(cosine->phi,cosine->n);
2512 sprintf(format,"%%*d");
2513 for(i=0; (i<cosine->n); i++) {
2515 strcat(f1,"%lf%lf");
2516 if (sscanf(t,f1,&a,&phi) < 2)
2517 gmx_fatal(FARGS,"Invalid input for electric field shift: '%s'",t);
2520 strcat(format,"%*lf%*lf");
2527 static gmx_bool do_egp_flag(t_inputrec *ir,gmx_groups_t *groups,
2528 const char *option,const char *val,int flag)
2530 /* The maximum number of energy group pairs would be MAXPTR*(MAXPTR+1)/2.
2531 * But since this is much larger than STRLEN, such a line can not be parsed.
2532 * The real maximum is the number of names that fit in a string: STRLEN/2.
2534 #define EGP_MAX (STRLEN/2)
2536 char *names[EGP_MAX];
2540 gnames = groups->grpname;
2542 nelem = str_nelem(val,EGP_MAX,names);
2544 gmx_fatal(FARGS,"The number of groups for %s is odd",option);
2545 nr = groups->grps[egcENER].nr;
2547 for(i=0; i<nelem/2; i++) {
2550 gmx_strcasecmp(names[2*i],*(gnames[groups->grps[egcENER].nm_ind[j]])))
2553 gmx_fatal(FARGS,"%s in %s is not an energy group\n",
2557 gmx_strcasecmp(names[2*i+1],*(gnames[groups->grps[egcENER].nm_ind[k]])))
2560 gmx_fatal(FARGS,"%s in %s is not an energy group\n",
2561 names[2*i+1],option);
2562 if ((j < nr) && (k < nr)) {
2563 ir->opts.egp_flags[nr*j+k] |= flag;
2564 ir->opts.egp_flags[nr*k+j] |= flag;
2572 void do_index(const char* mdparin, const char *ndx,
2575 t_inputrec *ir,rvec *v,
2579 gmx_groups_t *groups;
2583 char warnbuf[STRLEN],**gnames;
2584 int nr,ntcg,ntau_t,nref_t,nacc,nofg,nSA,nSA_points,nSA_time,nSA_temp;
2587 int nacg,nfreeze,nfrdim,nenergy,nvcm,nuser;
2588 char *ptr1[MAXPTR],*ptr2[MAXPTR],*ptr3[MAXPTR];
2591 gmx_bool bExcl,bTable,bSetTCpar,bAnneal,bRest;
2592 int nQMmethod,nQMbasis,nQMcharge,nQMmult,nbSH,nCASorb,nCASelec,
2593 nSAon,nSAoff,nSAsteps,nQMg,nbOPT,nbTS;
2594 char warn_buf[STRLEN];
2597 fprintf(stderr,"processing index file...\n");
2601 snew(grps->index,1);
2603 atoms_all = gmx_mtop_global_atoms(mtop);
2604 analyse(&atoms_all,grps,&gnames,FALSE,TRUE);
2605 free_t_atoms(&atoms_all,FALSE);
2607 grps = init_index(ndx,&gnames);
2610 groups = &mtop->groups;
2611 natoms = mtop->natoms;
2612 symtab = &mtop->symtab;
2614 snew(groups->grpname,grps->nr+1);
2616 for(i=0; (i<grps->nr); i++) {
2617 groups->grpname[i] = put_symtab(symtab,gnames[i]);
2619 groups->grpname[i] = put_symtab(symtab,"rest");
2621 srenew(gnames,grps->nr+1);
2622 gnames[restnm] = *(groups->grpname[i]);
2623 groups->ngrpname = grps->nr+1;
2625 set_warning_line(wi,mdparin,-1);
2627 ntau_t = str_nelem(tau_t,MAXPTR,ptr1);
2628 nref_t = str_nelem(ref_t,MAXPTR,ptr2);
2629 ntcg = str_nelem(tcgrps,MAXPTR,ptr3);
2630 if ((ntau_t != ntcg) || (nref_t != ntcg)) {
2631 gmx_fatal(FARGS,"Invalid T coupling input: %d groups, %d ref-t values and "
2632 "%d tau-t values",ntcg,nref_t,ntau_t);
2635 bSetTCpar = (ir->etc || EI_SD(ir->eI) || ir->eI==eiBD || EI_TPI(ir->eI));
2636 do_numbering(natoms,groups,ntcg,ptr3,grps,gnames,egcTC,
2637 restnm,bSetTCpar ? egrptpALL : egrptpALL_GENREST,bVerbose,wi);
2638 nr = groups->grps[egcTC].nr;
2640 snew(ir->opts.nrdf,nr);
2641 snew(ir->opts.tau_t,nr);
2642 snew(ir->opts.ref_t,nr);
2643 if (ir->eI==eiBD && ir->bd_fric==0) {
2644 fprintf(stderr,"bd-fric=0, so tau-t will be used as the inverse friction constant(s)\n");
2651 gmx_fatal(FARGS,"Not enough ref-t and tau-t values!");
2655 for(i=0; (i<nr); i++)
2657 ir->opts.tau_t[i] = strtod(ptr1[i],NULL);
2658 if ((ir->eI == eiBD || ir->eI == eiSD2) && ir->opts.tau_t[i] <= 0)
2660 sprintf(warn_buf,"With integrator %s tau-t should be larger than 0",ei_names[ir->eI]);
2661 warning_error(wi,warn_buf);
2664 if (ir->etc != etcVRESCALE && ir->opts.tau_t[i] == 0)
2666 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.");
2669 if (ir->opts.tau_t[i] >= 0)
2671 tau_min = min(tau_min,ir->opts.tau_t[i]);
2674 if (ir->etc != etcNO && ir->nsttcouple == -1)
2676 ir->nsttcouple = ir_optimal_nsttcouple(ir);
2681 if ((ir->etc==etcNOSEHOOVER) && (ir->epc==epcBERENDSEN)) {
2682 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");
2684 if ((ir->epc==epcMTTK) && (ir->etc>etcNO))
2687 mincouple = ir->nsttcouple;
2688 if (ir->nstpcouple < mincouple)
2690 mincouple = ir->nstpcouple;
2692 ir->nstpcouple = mincouple;
2693 ir->nsttcouple = mincouple;
2694 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);
2695 warning_note(wi,warn_buf);
2698 /* velocity verlet with averaged kinetic energy KE = 0.5*(v(t+1/2) - v(t-1/2)) is implemented
2699 primarily for testing purposes, and does not work with temperature coupling other than 1 */
2701 if (ETC_ANDERSEN(ir->etc)) {
2702 if (ir->nsttcouple != 1) {
2704 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");
2705 warning_note(wi,warn_buf);
2708 nstcmin = tcouple_min_integration_steps(ir->etc);
2711 if (tau_min/(ir->delta_t*ir->nsttcouple) < nstcmin)
2713 sprintf(warn_buf,"For proper integration of the %s thermostat, tau-t (%g) should be at least %d times larger than nsttcouple*dt (%g)",
2714 ETCOUPLTYPE(ir->etc),
2716 ir->nsttcouple*ir->delta_t);
2717 warning(wi,warn_buf);
2720 for(i=0; (i<nr); i++)
2722 ir->opts.ref_t[i] = strtod(ptr2[i],NULL);
2723 if (ir->opts.ref_t[i] < 0)
2725 gmx_fatal(FARGS,"ref-t for group %d negative",i);
2728 /* set the lambda mc temperature to the md integrator temperature (which should be defined
2729 if we are in this conditional) if mc_temp is negative */
2730 if (ir->expandedvals->mc_temp < 0)
2732 ir->expandedvals->mc_temp = ir->opts.ref_t[0]; /*for now, set to the first reft */
2736 /* Simulated annealing for each group. There are nr groups */
2737 nSA = str_nelem(anneal,MAXPTR,ptr1);
2738 if (nSA == 1 && (ptr1[0][0]=='n' || ptr1[0][0]=='N'))
2740 if(nSA>0 && nSA != nr)
2741 gmx_fatal(FARGS,"Not enough annealing values: %d (for %d groups)\n",nSA,nr);
2743 snew(ir->opts.annealing,nr);
2744 snew(ir->opts.anneal_npoints,nr);
2745 snew(ir->opts.anneal_time,nr);
2746 snew(ir->opts.anneal_temp,nr);
2748 ir->opts.annealing[i]=eannNO;
2749 ir->opts.anneal_npoints[i]=0;
2750 ir->opts.anneal_time[i]=NULL;
2751 ir->opts.anneal_temp[i]=NULL;
2756 if(ptr1[i][0]=='n' || ptr1[i][0]=='N') {
2757 ir->opts.annealing[i]=eannNO;
2758 } else if(ptr1[i][0]=='s'|| ptr1[i][0]=='S') {
2759 ir->opts.annealing[i]=eannSINGLE;
2761 } else if(ptr1[i][0]=='p'|| ptr1[i][0]=='P') {
2762 ir->opts.annealing[i]=eannPERIODIC;
2767 /* Read the other fields too */
2768 nSA_points = str_nelem(anneal_npoints,MAXPTR,ptr1);
2770 gmx_fatal(FARGS,"Found %d annealing-npoints values for %d groups\n",nSA_points,nSA);
2771 for(k=0,i=0;i<nr;i++) {
2772 ir->opts.anneal_npoints[i]=strtol(ptr1[i],NULL,10);
2773 if(ir->opts.anneal_npoints[i]==1)
2774 gmx_fatal(FARGS,"Please specify at least a start and an end point for annealing\n");
2775 snew(ir->opts.anneal_time[i],ir->opts.anneal_npoints[i]);
2776 snew(ir->opts.anneal_temp[i],ir->opts.anneal_npoints[i]);
2777 k += ir->opts.anneal_npoints[i];
2780 nSA_time = str_nelem(anneal_time,MAXPTR,ptr1);
2782 gmx_fatal(FARGS,"Found %d annealing-time values, wanter %d\n",nSA_time,k);
2783 nSA_temp = str_nelem(anneal_temp,MAXPTR,ptr2);
2785 gmx_fatal(FARGS,"Found %d annealing-temp values, wanted %d\n",nSA_temp,k);
2787 for(i=0,k=0;i<nr;i++) {
2789 for(j=0;j<ir->opts.anneal_npoints[i];j++) {
2790 ir->opts.anneal_time[i][j]=strtod(ptr1[k],NULL);
2791 ir->opts.anneal_temp[i][j]=strtod(ptr2[k],NULL);
2793 if(ir->opts.anneal_time[i][0] > (ir->init_t+GMX_REAL_EPS))
2794 gmx_fatal(FARGS,"First time point for annealing > init_t.\n");
2797 if(ir->opts.anneal_time[i][j]<ir->opts.anneal_time[i][j-1])
2798 gmx_fatal(FARGS,"Annealing timepoints out of order: t=%f comes after t=%f\n",
2799 ir->opts.anneal_time[i][j],ir->opts.anneal_time[i][j-1]);
2801 if(ir->opts.anneal_temp[i][j]<0)
2802 gmx_fatal(FARGS,"Found negative temperature in annealing: %f\n",ir->opts.anneal_temp[i][j]);
2806 /* Print out some summary information, to make sure we got it right */
2807 for(i=0,k=0;i<nr;i++) {
2808 if(ir->opts.annealing[i]!=eannNO) {
2809 j = groups->grps[egcTC].nm_ind[i];
2810 fprintf(stderr,"Simulated annealing for group %s: %s, %d timepoints\n",
2811 *(groups->grpname[j]),eann_names[ir->opts.annealing[i]],
2812 ir->opts.anneal_npoints[i]);
2813 fprintf(stderr,"Time (ps) Temperature (K)\n");
2814 /* All terms except the last one */
2815 for(j=0;j<(ir->opts.anneal_npoints[i]-1);j++)
2816 fprintf(stderr,"%9.1f %5.1f\n",ir->opts.anneal_time[i][j],ir->opts.anneal_temp[i][j]);
2818 /* Finally the last one */
2819 j = ir->opts.anneal_npoints[i]-1;
2820 if(ir->opts.annealing[i]==eannSINGLE)
2821 fprintf(stderr,"%9.1f- %5.1f\n",ir->opts.anneal_time[i][j],ir->opts.anneal_temp[i][j]);
2823 fprintf(stderr,"%9.1f %5.1f\n",ir->opts.anneal_time[i][j],ir->opts.anneal_temp[i][j]);
2824 if(fabs(ir->opts.anneal_temp[i][j]-ir->opts.anneal_temp[i][0])>GMX_REAL_EPS)
2825 warning_note(wi,"There is a temperature jump when your annealing loops back.\n");
2833 if (ir->ePull != epullNO) {
2834 make_pull_groups(ir->pull,pull_grp,grps,gnames);
2838 make_rotation_groups(ir->rot,rot_grp,grps,gnames);
2841 nacc = str_nelem(acc,MAXPTR,ptr1);
2842 nacg = str_nelem(accgrps,MAXPTR,ptr2);
2843 if (nacg*DIM != nacc)
2844 gmx_fatal(FARGS,"Invalid Acceleration input: %d groups and %d acc. values",
2846 do_numbering(natoms,groups,nacg,ptr2,grps,gnames,egcACC,
2847 restnm,egrptpALL_GENREST,bVerbose,wi);
2848 nr = groups->grps[egcACC].nr;
2849 snew(ir->opts.acc,nr);
2852 for(i=k=0; (i<nacg); i++)
2853 for(j=0; (j<DIM); j++,k++)
2854 ir->opts.acc[i][j]=strtod(ptr1[k],NULL);
2856 for(j=0; (j<DIM); j++)
2857 ir->opts.acc[i][j]=0;
2859 nfrdim = str_nelem(frdim,MAXPTR,ptr1);
2860 nfreeze = str_nelem(freeze,MAXPTR,ptr2);
2861 if (nfrdim != DIM*nfreeze)
2862 gmx_fatal(FARGS,"Invalid Freezing input: %d groups and %d freeze values",
2864 do_numbering(natoms,groups,nfreeze,ptr2,grps,gnames,egcFREEZE,
2865 restnm,egrptpALL_GENREST,bVerbose,wi);
2866 nr = groups->grps[egcFREEZE].nr;
2868 snew(ir->opts.nFreeze,nr);
2869 for(i=k=0; (i<nfreeze); i++)
2870 for(j=0; (j<DIM); j++,k++) {
2871 ir->opts.nFreeze[i][j]=(gmx_strncasecmp(ptr1[k],"Y",1)==0);
2872 if (!ir->opts.nFreeze[i][j]) {
2873 if (gmx_strncasecmp(ptr1[k],"N",1) != 0) {
2874 sprintf(warnbuf,"Please use Y(ES) or N(O) for freezedim only "
2875 "(not %s)", ptr1[k]);
2876 warning(wi,warn_buf);
2881 for(j=0; (j<DIM); j++)
2882 ir->opts.nFreeze[i][j]=0;
2884 nenergy=str_nelem(energy,MAXPTR,ptr1);
2885 do_numbering(natoms,groups,nenergy,ptr1,grps,gnames,egcENER,
2886 restnm,egrptpALL_GENREST,bVerbose,wi);
2887 add_wall_energrps(groups,ir->nwall,symtab);
2888 ir->opts.ngener = groups->grps[egcENER].nr;
2889 nvcm=str_nelem(vcm,MAXPTR,ptr1);
2891 do_numbering(natoms,groups,nvcm,ptr1,grps,gnames,egcVCM,
2892 restnm,nvcm==0 ? egrptpALL_GENREST : egrptpPART,bVerbose,wi);
2894 warning(wi,"Some atoms are not part of any center of mass motion removal group.\n"
2895 "This may lead to artifacts.\n"
2896 "In most cases one should use one group for the whole system.");
2899 /* Now we have filled the freeze struct, so we can calculate NRDF */
2900 calc_nrdf(mtop,ir,gnames);
2905 /* Must check per group! */
2906 for(i=0; (i<ir->opts.ngtc); i++)
2907 ntot += ir->opts.nrdf[i];
2908 if (ntot != (DIM*natoms)) {
2909 fac = sqrt(ntot/(DIM*natoms));
2911 fprintf(stderr,"Scaling velocities by a factor of %.3f to account for constraints\n"
2912 "and removal of center of mass motion\n",fac);
2913 for(i=0; (i<natoms); i++)
2914 svmul(fac,v[i],v[i]);
2918 nuser=str_nelem(user1,MAXPTR,ptr1);
2919 do_numbering(natoms,groups,nuser,ptr1,grps,gnames,egcUser1,
2920 restnm,egrptpALL_GENREST,bVerbose,wi);
2921 nuser=str_nelem(user2,MAXPTR,ptr1);
2922 do_numbering(natoms,groups,nuser,ptr1,grps,gnames,egcUser2,
2923 restnm,egrptpALL_GENREST,bVerbose,wi);
2924 nuser=str_nelem(xtc_grps,MAXPTR,ptr1);
2925 do_numbering(natoms,groups,nuser,ptr1,grps,gnames,egcXTC,
2926 restnm,egrptpONE,bVerbose,wi);
2927 nofg = str_nelem(orirefitgrp,MAXPTR,ptr1);
2928 do_numbering(natoms,groups,nofg,ptr1,grps,gnames,egcORFIT,
2929 restnm,egrptpALL_GENREST,bVerbose,wi);
2931 /* QMMM input processing */
2932 nQMg = str_nelem(QMMM,MAXPTR,ptr1);
2933 nQMmethod = str_nelem(QMmethod,MAXPTR,ptr2);
2934 nQMbasis = str_nelem(QMbasis,MAXPTR,ptr3);
2935 if((nQMmethod != nQMg)||(nQMbasis != nQMg)){
2936 gmx_fatal(FARGS,"Invalid QMMM input: %d groups %d basissets"
2937 " and %d methods\n",nQMg,nQMbasis,nQMmethod);
2939 /* group rest, if any, is always MM! */
2940 do_numbering(natoms,groups,nQMg,ptr1,grps,gnames,egcQMMM,
2941 restnm,egrptpALL_GENREST,bVerbose,wi);
2942 nr = nQMg; /*atoms->grps[egcQMMM].nr;*/
2943 ir->opts.ngQM = nQMg;
2944 snew(ir->opts.QMmethod,nr);
2945 snew(ir->opts.QMbasis,nr);
2947 /* input consists of strings: RHF CASSCF PM3 .. These need to be
2948 * converted to the corresponding enum in names.c
2950 ir->opts.QMmethod[i] = search_QMstring(ptr2[i],eQMmethodNR,
2952 ir->opts.QMbasis[i] = search_QMstring(ptr3[i],eQMbasisNR,
2956 nQMmult = str_nelem(QMmult,MAXPTR,ptr1);
2957 nQMcharge = str_nelem(QMcharge,MAXPTR,ptr2);
2958 nbSH = str_nelem(bSH,MAXPTR,ptr3);
2959 snew(ir->opts.QMmult,nr);
2960 snew(ir->opts.QMcharge,nr);
2961 snew(ir->opts.bSH,nr);
2964 ir->opts.QMmult[i] = strtol(ptr1[i],NULL,10);
2965 ir->opts.QMcharge[i] = strtol(ptr2[i],NULL,10);
2966 ir->opts.bSH[i] = (gmx_strncasecmp(ptr3[i],"Y",1)==0);
2969 nCASelec = str_nelem(CASelectrons,MAXPTR,ptr1);
2970 nCASorb = str_nelem(CASorbitals,MAXPTR,ptr2);
2971 snew(ir->opts.CASelectrons,nr);
2972 snew(ir->opts.CASorbitals,nr);
2974 ir->opts.CASelectrons[i]= strtol(ptr1[i],NULL,10);
2975 ir->opts.CASorbitals[i] = strtol(ptr2[i],NULL,10);
2977 /* special optimization options */
2979 nbOPT = str_nelem(bOPT,MAXPTR,ptr1);
2980 nbTS = str_nelem(bTS,MAXPTR,ptr2);
2981 snew(ir->opts.bOPT,nr);
2982 snew(ir->opts.bTS,nr);
2984 ir->opts.bOPT[i] = (gmx_strncasecmp(ptr1[i],"Y",1)==0);
2985 ir->opts.bTS[i] = (gmx_strncasecmp(ptr2[i],"Y",1)==0);
2987 nSAon = str_nelem(SAon,MAXPTR,ptr1);
2988 nSAoff = str_nelem(SAoff,MAXPTR,ptr2);
2989 nSAsteps = str_nelem(SAsteps,MAXPTR,ptr3);
2990 snew(ir->opts.SAon,nr);
2991 snew(ir->opts.SAoff,nr);
2992 snew(ir->opts.SAsteps,nr);
2995 ir->opts.SAon[i] = strtod(ptr1[i],NULL);
2996 ir->opts.SAoff[i] = strtod(ptr2[i],NULL);
2997 ir->opts.SAsteps[i] = strtol(ptr3[i],NULL,10);
2999 /* end of QMMM input */
3002 for(i=0; (i<egcNR); i++) {
3003 fprintf(stderr,"%-16s has %d element(s):",gtypes[i],groups->grps[i].nr);
3004 for(j=0; (j<groups->grps[i].nr); j++)
3005 fprintf(stderr," %s",*(groups->grpname[groups->grps[i].nm_ind[j]]));
3006 fprintf(stderr,"\n");
3009 nr = groups->grps[egcENER].nr;
3010 snew(ir->opts.egp_flags,nr*nr);
3012 bExcl = do_egp_flag(ir,groups,"energygrp-excl",egpexcl,EGP_EXCL);
3013 if (bExcl && ir->cutoff_scheme == ecutsVERLET)
3015 warning_error(wi,"Energy group exclusions are not (yet) implemented for the Verlet scheme");
3017 if (bExcl && EEL_FULL(ir->coulombtype))
3018 warning(wi,"Can not exclude the lattice Coulomb energy between energy groups");
3020 bTable = do_egp_flag(ir,groups,"energygrp-table",egptable,EGP_TABLE);
3021 if (bTable && !(ir->vdwtype == evdwUSER) &&
3022 !(ir->coulombtype == eelUSER) && !(ir->coulombtype == eelPMEUSER) &&
3023 !(ir->coulombtype == eelPMEUSERSWITCH))
3024 gmx_fatal(FARGS,"Can only have energy group pair tables in combination with user tables for VdW and/or Coulomb");
3026 decode_cos(efield_x,&(ir->ex[XX]),FALSE);
3027 decode_cos(efield_xt,&(ir->et[XX]),TRUE);
3028 decode_cos(efield_y,&(ir->ex[YY]),FALSE);
3029 decode_cos(efield_yt,&(ir->et[YY]),TRUE);
3030 decode_cos(efield_z,&(ir->ex[ZZ]),FALSE);
3031 decode_cos(efield_zt,&(ir->et[ZZ]),TRUE);
3034 do_adress_index(ir->adress,groups,gnames,&(ir->opts),wi);
3036 for(i=0; (i<grps->nr); i++)
3046 static void check_disre(gmx_mtop_t *mtop)
3048 gmx_ffparams_t *ffparams;
3049 t_functype *functype;
3051 int i,ndouble,ftype;
3052 int label,old_label;
3054 if (gmx_mtop_ftype_count(mtop,F_DISRES) > 0) {
3055 ffparams = &mtop->ffparams;
3056 functype = ffparams->functype;
3057 ip = ffparams->iparams;
3060 for(i=0; i<ffparams->ntypes; i++) {
3061 ftype = functype[i];
3062 if (ftype == F_DISRES) {
3063 label = ip[i].disres.label;
3064 if (label == old_label) {
3065 fprintf(stderr,"Distance restraint index %d occurs twice\n",label);
3072 gmx_fatal(FARGS,"Found %d double distance restraint indices,\n"
3073 "probably the parameters for multiple pairs in one restraint "
3074 "are not identical\n",ndouble);
3078 static gmx_bool absolute_reference(t_inputrec *ir,gmx_mtop_t *sys,
3079 gmx_bool posres_only,
3083 gmx_mtop_ilistloop_t iloop;
3093 for(d=0; d<DIM; d++)
3095 AbsRef[d] = (d < ndof_com(ir) ? 0 : 1);
3097 /* Check for freeze groups */
3098 for(g=0; g<ir->opts.ngfrz; g++)
3100 for(d=0; d<DIM; d++)
3102 if (ir->opts.nFreeze[g][d] != 0)
3110 /* Check for position restraints */
3111 iloop = gmx_mtop_ilistloop_init(sys);
3112 while (gmx_mtop_ilistloop_next(iloop,&ilist,&nmol))
3115 (AbsRef[XX] == 0 || AbsRef[YY] == 0 || AbsRef[ZZ] == 0))
3117 for(i=0; i<ilist[F_POSRES].nr; i+=2)
3119 pr = &sys->ffparams.iparams[ilist[F_POSRES].iatoms[i]];
3120 for(d=0; d<DIM; d++)
3122 if (pr->posres.fcA[d] != 0)
3131 return (AbsRef[XX] != 0 && AbsRef[YY] != 0 && AbsRef[ZZ] != 0);
3134 void triple_check(const char *mdparin,t_inputrec *ir,gmx_mtop_t *sys,
3138 int i,m,g,nmol,npct;
3139 gmx_bool bCharge,bAcc;
3140 real gdt_max,*mgrp,mt;
3142 gmx_mtop_atomloop_block_t aloopb;
3143 gmx_mtop_atomloop_all_t aloop;
3146 char warn_buf[STRLEN];
3148 set_warning_line(wi,mdparin,-1);
3150 if (EI_DYNAMICS(ir->eI) && !EI_SD(ir->eI) && ir->eI != eiBD &&
3151 ir->comm_mode == ecmNO &&
3152 !(absolute_reference(ir,sys,FALSE,AbsRef) || ir->nsteps <= 10)) {
3153 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");
3156 /* Check for pressure coupling with absolute position restraints */
3157 if (ir->epc != epcNO && ir->refcoord_scaling == erscNO)
3159 absolute_reference(ir,sys,TRUE,AbsRef);
3161 for(m=0; m<DIM; m++)
3163 if (AbsRef[m] && norm2(ir->compress[m]) > 0)
3165 warning(wi,"You are using pressure coupling with absolute position restraints, this will give artifacts. Use the refcoord_scaling option.");
3173 aloopb = gmx_mtop_atomloop_block_init(sys);
3174 while (gmx_mtop_atomloop_block_next(aloopb,&atom,&nmol)) {
3175 if (atom->q != 0 || atom->qB != 0) {
3181 if (EEL_FULL(ir->coulombtype)) {
3183 "You are using full electrostatics treatment %s for a system without charges.\n"
3184 "This costs a lot of performance for just processing zeros, consider using %s instead.\n",
3185 EELTYPE(ir->coulombtype),EELTYPE(eelCUT));
3186 warning(wi,err_buf);
3189 if (ir->coulombtype == eelCUT && ir->rcoulomb > 0 && !ir->implicit_solvent) {
3191 "You are using a plain Coulomb cut-off, which might produce artifacts.\n"
3192 "You might want to consider using %s electrostatics.\n",
3194 warning_note(wi,err_buf);
3198 /* Generalized reaction field */
3199 if (ir->opts.ngtc == 0) {
3200 sprintf(err_buf,"No temperature coupling while using coulombtype %s",
3202 CHECK(ir->coulombtype == eelGRF);
3205 sprintf(err_buf,"When using coulombtype = %s"
3206 " ref-t for temperature coupling should be > 0",
3208 CHECK((ir->coulombtype == eelGRF) && (ir->opts.ref_t[0] <= 0));
3211 if (ir->eI == eiSD1 &&
3212 (gmx_mtop_ftype_count(sys,F_CONSTR) > 0 ||
3213 gmx_mtop_ftype_count(sys,F_SETTLE) > 0))
3215 sprintf(warn_buf,"With constraints integrator %s is less accurate, consider using %s instead",ei_names[ir->eI],ei_names[eiSD2]);
3216 warning_note(wi,warn_buf);
3220 for(i=0; (i<sys->groups.grps[egcACC].nr); i++) {
3221 for(m=0; (m<DIM); m++) {
3222 if (fabs(ir->opts.acc[i][m]) > 1e-6) {
3229 snew(mgrp,sys->groups.grps[egcACC].nr);
3230 aloop = gmx_mtop_atomloop_all_init(sys);
3231 while (gmx_mtop_atomloop_all_next(aloop,&i,&atom)) {
3232 mgrp[ggrpnr(&sys->groups,egcACC,i)] += atom->m;
3235 for(i=0; (i<sys->groups.grps[egcACC].nr); i++) {
3236 for(m=0; (m<DIM); m++)
3237 acc[m] += ir->opts.acc[i][m]*mgrp[i];
3240 for(m=0; (m<DIM); m++) {
3241 if (fabs(acc[m]) > 1e-6) {
3242 const char *dim[DIM] = { "X", "Y", "Z" };
3244 "Net Acceleration in %s direction, will %s be corrected\n",
3245 dim[m],ir->nstcomm != 0 ? "" : "not");
3246 if (ir->nstcomm != 0 && m < ndof_com(ir)) {
3248 for (i=0; (i<sys->groups.grps[egcACC].nr); i++)
3249 ir->opts.acc[i][m] -= acc[m];
3256 if (ir->efep != efepNO && ir->fepvals->sc_alpha != 0 &&
3257 !gmx_within_tol(sys->ffparams.reppow,12.0,10*GMX_DOUBLE_EPS)) {
3258 gmx_fatal(FARGS,"Soft-core interactions are only supported with VdW repulsion power 12");
3261 if (ir->ePull != epullNO) {
3262 if (ir->pull->grp[0].nat == 0) {
3263 absolute_reference(ir,sys,FALSE,AbsRef);
3264 for(m=0; m<DIM; m++) {
3265 if (ir->pull->dim[m] && !AbsRef[m]) {
3266 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.");
3272 if (ir->pull->eGeom == epullgDIRPBC) {
3273 for(i=0; i<3; i++) {
3274 for(m=0; m<=i; m++) {
3275 if ((ir->epc != epcNO && ir->compress[i][m] != 0) ||
3276 ir->deform[i][m] != 0) {
3277 for(g=1; g<ir->pull->ngrp; g++) {
3278 if (ir->pull->grp[g].vec[m] != 0) {
3279 gmx_fatal(FARGS,"Can not have dynamic box while using pull geometry '%s' (dim %c)",EPULLGEOM(ir->pull->eGeom),'x'+m);
3291 void double_check(t_inputrec *ir,matrix box,gmx_bool bConstr,warninp_t wi)
3295 char warn_buf[STRLEN];
3298 ptr = check_box(ir->ePBC,box);
3300 warning_error(wi,ptr);
3303 if (bConstr && ir->eConstrAlg == econtSHAKE) {
3304 if (ir->shake_tol <= 0.0) {
3305 sprintf(warn_buf,"ERROR: shake-tol must be > 0 instead of %g\n",
3307 warning_error(wi,warn_buf);
3310 if (IR_TWINRANGE(*ir) && ir->nstlist > 1) {
3311 sprintf(warn_buf,"With twin-range cut-off's and SHAKE the virial and the pressure are incorrect.");
3312 if (ir->epc == epcNO) {
3313 warning(wi,warn_buf);
3315 warning_error(wi,warn_buf);
3320 if( (ir->eConstrAlg == econtLINCS) && bConstr) {
3321 /* If we have Lincs constraints: */
3322 if(ir->eI==eiMD && ir->etc==etcNO &&
3323 ir->eConstrAlg==econtLINCS && ir->nLincsIter==1) {
3324 sprintf(warn_buf,"For energy conservation with LINCS, lincs_iter should be 2 or larger.\n");
3325 warning_note(wi,warn_buf);
3328 if ((ir->eI == eiCG || ir->eI == eiLBFGS) && (ir->nProjOrder<8)) {
3329 sprintf(warn_buf,"For accurate %s with LINCS constraints, lincs-order should be 8 or more.",ei_names[ir->eI]);
3330 warning_note(wi,warn_buf);
3332 if (ir->epc==epcMTTK) {
3333 warning_error(wi,"MTTK not compatible with lincs -- use shake instead.");
3337 if (ir->LincsWarnAngle > 90.0) {
3338 sprintf(warn_buf,"lincs-warnangle can not be larger than 90 degrees, setting it to 90.\n");
3339 warning(wi,warn_buf);
3340 ir->LincsWarnAngle = 90.0;
3343 if (ir->ePBC != epbcNONE) {
3344 if (ir->nstlist == 0) {
3345 warning(wi,"With nstlist=0 atoms are only put into the box at step 0, therefore drifting atoms might cause the simulation to crash.");
3347 bTWIN = (ir->rlistlong > ir->rlist);
3348 if (ir->ns_type == ensGRID) {
3349 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC,box)) {
3350 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",
3351 bTWIN ? (ir->rcoulomb==ir->rlistlong ? "rcoulomb" : "rvdw"):"rlist");
3352 warning_error(wi,warn_buf);
3355 min_size = min(box[XX][XX],min(box[YY][YY],box[ZZ][ZZ]));
3356 if (2*ir->rlistlong >= min_size) {
3357 sprintf(warn_buf,"ERROR: One of the box lengths is smaller than twice the cut-off length. Increase the box size or decrease rlist.");
3358 warning_error(wi,warn_buf);
3360 fprintf(stderr,"Grid search might allow larger cut-off's than simple search with triclinic boxes.");
3366 void check_chargegroup_radii(const gmx_mtop_t *mtop,const t_inputrec *ir,
3370 real rvdw1,rvdw2,rcoul1,rcoul2;
3371 char warn_buf[STRLEN];
3373 calc_chargegroup_radii(mtop,x,&rvdw1,&rvdw2,&rcoul1,&rcoul2);
3377 printf("Largest charge group radii for Van der Waals: %5.3f, %5.3f nm\n",
3382 printf("Largest charge group radii for Coulomb: %5.3f, %5.3f nm\n",
3388 if (rvdw1 + rvdw2 > ir->rlist ||
3389 rcoul1 + rcoul2 > ir->rlist)
3391 sprintf(warn_buf,"The sum of the two largest charge group radii (%f) is larger than rlist (%f)\n",max(rvdw1+rvdw2,rcoul1+rcoul2),ir->rlist);
3392 warning(wi,warn_buf);
3396 /* Here we do not use the zero at cut-off macro,
3397 * since user defined interactions might purposely
3398 * not be zero at the cut-off.
3400 if (EVDW_IS_ZERO_AT_CUTOFF(ir->vdwtype) &&
3401 rvdw1 + rvdw2 > ir->rlist - ir->rvdw)
3403 sprintf(warn_buf,"The sum of the two largest charge group radii (%f) is larger than rlist (%f) - rvdw (%f)\n",
3405 ir->rlist,ir->rvdw);
3408 warning(wi,warn_buf);
3412 warning_note(wi,warn_buf);
3415 if (EEL_IS_ZERO_AT_CUTOFF(ir->coulombtype) &&
3416 rcoul1 + rcoul2 > ir->rlistlong - ir->rcoulomb)
3418 sprintf(warn_buf,"The sum of the two largest charge group radii (%f) is larger than %s (%f) - rcoulomb (%f)\n",
3420 ir->rlistlong > ir->rlist ? "rlistlong" : "rlist",
3421 ir->rlistlong,ir->rcoulomb);
3424 warning(wi,warn_buf);
3428 warning_note(wi,warn_buf);