9e39ad6edea032c6d4e8b4e87b9f97a6fb602bbc
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / readir.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "readir.h"
40
41 #include <ctype.h>
42 #include <limits.h>
43 #include <stdlib.h>
44
45 #include "gromacs/gmxpreprocess/toputil.h"
46 #include "gromacs/legacyheaders/chargegroup.h"
47 #include "gromacs/legacyheaders/inputrec.h"
48 #include "gromacs/legacyheaders/macros.h"
49 #include "gromacs/legacyheaders/names.h"
50 #include "gromacs/legacyheaders/network.h"
51 #include "gromacs/legacyheaders/readinp.h"
52 #include "gromacs/legacyheaders/typedefs.h"
53 #include "gromacs/legacyheaders/warninp.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/mdlib/calc_verletbuf.h"
57 #include "gromacs/pbcutil/pbc.h"
58 #include "gromacs/topology/block.h"
59 #include "gromacs/topology/index.h"
60 #include "gromacs/topology/mtop_util.h"
61 #include "gromacs/topology/symtab.h"
62 #include "gromacs/utility/cstringutil.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/smalloc.h"
65
66 #define MAXPTR 254
67 #define NOGID  255
68
69 /* Resource parameters
70  * Do not change any of these until you read the instruction
71  * in readinp.h. Some cpp's do not take spaces after the backslash
72  * (like the c-shell), which will give you a very weird compiler
73  * message.
74  */
75
76 typedef struct t_inputrec_strings
77 {
78     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], x_compressed_groups[STRLEN],
81          couple_moltype[STRLEN], orirefitgrp[STRLEN], egptable[STRLEN], egpexcl[STRLEN],
82          wall_atomtype[STRLEN], wall_density[STRLEN], deform[STRLEN], QMMM[STRLEN],
83          imd_grp[STRLEN];
84     char   fep_lambda[efptNR][STRLEN];
85     char   lambda_weights[STRLEN];
86     char **pull_grp;
87     char **rot_grp;
88     char   anneal[STRLEN], anneal_npoints[STRLEN],
89            anneal_time[STRLEN], anneal_temp[STRLEN];
90     char   QMmethod[STRLEN], QMbasis[STRLEN], QMcharge[STRLEN], QMmult[STRLEN],
91            bSH[STRLEN], CASorbitals[STRLEN], CASelectrons[STRLEN], SAon[STRLEN],
92            SAoff[STRLEN], SAsteps[STRLEN], bTS[STRLEN], bOPT[STRLEN];
93     char efield_x[STRLEN], efield_xt[STRLEN], efield_y[STRLEN],
94          efield_yt[STRLEN], efield_z[STRLEN], efield_zt[STRLEN];
95
96 } gmx_inputrec_strings;
97
98 static gmx_inputrec_strings *is = NULL;
99
100 void init_inputrec_strings()
101 {
102     if (is)
103     {
104         gmx_incons("Attempted to call init_inputrec_strings before calling done_inputrec_strings. Only one inputrec (i.e. .mdp file) can be parsed at a time.");
105     }
106     snew(is, 1);
107 }
108
109 void done_inputrec_strings()
110 {
111     sfree(is);
112     is = NULL;
113 }
114
115 static char swapgrp[STRLEN], splitgrp0[STRLEN], splitgrp1[STRLEN], solgrp[STRLEN];
116
117 enum {
118     egrptpALL,         /* All particles have to be a member of a group.     */
119     egrptpALL_GENREST, /* A rest group with name is generated for particles *
120                         * that are not part of any group.                   */
121     egrptpPART,        /* As egrptpALL_GENREST, but no name is generated    *
122                         * for the rest group.                               */
123     egrptpONE          /* Merge all selected groups into one group,         *
124                         * make a rest group for the remaining particles.    */
125 };
126
127 static const char *constraints[eshNR+1]    = {
128     "none", "h-bonds", "all-bonds", "h-angles", "all-angles", NULL
129 };
130
131 static const char *couple_lam[ecouplamNR+1]    = {
132     "vdw-q", "vdw", "q", "none", NULL
133 };
134
135 void init_ir(t_inputrec *ir, t_gromppopts *opts)
136 {
137     snew(opts->include, STRLEN);
138     snew(opts->define, STRLEN);
139     snew(ir->fepvals, 1);
140     snew(ir->expandedvals, 1);
141     snew(ir->simtempvals, 1);
142 }
143
144 static void GetSimTemps(int ntemps, t_simtemp *simtemp, double *temperature_lambdas)
145 {
146
147     int i;
148
149     for (i = 0; i < ntemps; i++)
150     {
151         /* simple linear scaling -- allows more control */
152         if (simtemp->eSimTempScale == esimtempLINEAR)
153         {
154             simtemp->temperatures[i] = simtemp->simtemp_low + (simtemp->simtemp_high-simtemp->simtemp_low)*temperature_lambdas[i];
155         }
156         else if (simtemp->eSimTempScale == esimtempGEOMETRIC)  /* should give roughly equal acceptance for constant heat capacity . . . */
157         {
158             simtemp->temperatures[i] = simtemp->simtemp_low * pow(simtemp->simtemp_high/simtemp->simtemp_low, (1.0*i)/(ntemps-1));
159         }
160         else if (simtemp->eSimTempScale == esimtempEXPONENTIAL)
161         {
162             simtemp->temperatures[i] = simtemp->simtemp_low + (simtemp->simtemp_high-simtemp->simtemp_low)*(gmx_expm1(temperature_lambdas[i])/gmx_expm1(1.0));
163         }
164         else
165         {
166             char errorstr[128];
167             sprintf(errorstr, "eSimTempScale=%d not defined", simtemp->eSimTempScale);
168             gmx_fatal(FARGS, errorstr);
169         }
170     }
171 }
172
173
174
175 static void _low_check(gmx_bool b, char *s, warninp_t wi)
176 {
177     if (b)
178     {
179         warning_error(wi, s);
180     }
181 }
182
183 static void check_nst(const char *desc_nst, int nst,
184                       const char *desc_p, int *p,
185                       warninp_t wi)
186 {
187     char buf[STRLEN];
188
189     if (*p > 0 && *p % nst != 0)
190     {
191         /* Round up to the next multiple of nst */
192         *p = ((*p)/nst + 1)*nst;
193         sprintf(buf, "%s should be a multiple of %s, changing %s to %d\n",
194                 desc_p, desc_nst, desc_p, *p);
195         warning(wi, buf);
196     }
197 }
198
199 static gmx_bool ir_NVE(const t_inputrec *ir)
200 {
201     return ((ir->eI == eiMD || EI_VV(ir->eI)) && ir->etc == etcNO);
202 }
203
204 static int lcd(int n1, int n2)
205 {
206     int d, i;
207
208     d = 1;
209     for (i = 2; (i <= n1 && i <= n2); i++)
210     {
211         if (n1 % i == 0 && n2 % i == 0)
212         {
213             d = i;
214         }
215     }
216
217     return d;
218 }
219
220 static void process_interaction_modifier(const t_inputrec *ir, int *eintmod)
221 {
222     if (*eintmod == eintmodPOTSHIFT_VERLET)
223     {
224         if (ir->cutoff_scheme == ecutsVERLET)
225         {
226             *eintmod = eintmodPOTSHIFT;
227         }
228         else
229         {
230             *eintmod = eintmodNONE;
231         }
232     }
233 }
234
235 void check_ir(const char *mdparin, t_inputrec *ir, t_gromppopts *opts,
236               warninp_t wi)
237 /* Check internal consistency.
238  * NOTE: index groups are not set here yet, don't check things
239  * like temperature coupling group options here, but in triple_check
240  */
241 {
242     /* Strange macro: first one fills the err_buf, and then one can check
243      * the condition, which will print the message and increase the error
244      * counter.
245      */
246 #define CHECK(b) _low_check(b, err_buf, wi)
247     char        err_buf[256], warn_buf[STRLEN];
248     int         i, j;
249     int         ns_type  = 0;
250     real        dt_coupl = 0;
251     real        dt_pcoupl;
252     int         nstcmin;
253     t_lambda   *fep    = ir->fepvals;
254     t_expanded *expand = ir->expandedvals;
255
256     set_warning_line(wi, mdparin, -1);
257
258     /* BASIC CUT-OFF STUFF */
259     if (ir->rcoulomb < 0)
260     {
261         warning_error(wi, "rcoulomb should be >= 0");
262     }
263     if (ir->rvdw < 0)
264     {
265         warning_error(wi, "rvdw should be >= 0");
266     }
267     if (ir->rlist < 0 &&
268         !(ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0))
269     {
270         warning_error(wi, "rlist should be >= 0");
271     }
272     sprintf(err_buf, "nstlist can not be smaller than 0. (If you were trying to use the heuristic neighbour-list update scheme for efficient buffering for improved energy conservation, please use the Verlet cut-off scheme instead.)");
273     CHECK(ir->nstlist < 0);
274
275     process_interaction_modifier(ir, &ir->coulomb_modifier);
276     process_interaction_modifier(ir, &ir->vdw_modifier);
277
278     if (ir->cutoff_scheme == ecutsGROUP)
279     {
280         warning_note(wi,
281                      "The group cutoff scheme is deprecated since GROMACS 5.0 and will be removed in a future "
282                      "release when all interaction forms are supported for the verlet scheme. The verlet "
283                      "scheme already scales better, and it is compatible with GPUs and other accelerators.");
284
285         /* BASIC CUT-OFF STUFF */
286         if (ir->rlist == 0 ||
287             !((ir_coulomb_might_be_zero_at_cutoff(ir) && ir->rcoulomb > ir->rlist) ||
288               (ir_vdw_might_be_zero_at_cutoff(ir)     && ir->rvdw     > ir->rlist)))
289         {
290             /* No switched potential and/or no twin-range:
291              * we can set the long-range cut-off to the maximum of the other cut-offs.
292              */
293             ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
294         }
295         else if (ir->rlistlong < 0)
296         {
297             ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
298             sprintf(warn_buf, "rlistlong was not set, setting it to %g (no buffer)",
299                     ir->rlistlong);
300             warning(wi, warn_buf);
301         }
302         if (ir->rlistlong == 0 && ir->ePBC != epbcNONE)
303         {
304             warning_error(wi, "Can not have an infinite cut-off with PBC");
305         }
306         if (ir->rlistlong > 0 && (ir->rlist == 0 || ir->rlistlong < ir->rlist))
307         {
308             warning_error(wi, "rlistlong can not be shorter than rlist");
309         }
310         if (IR_TWINRANGE(*ir) && ir->nstlist == 0)
311         {
312             warning_error(wi, "Can not have nstlist == 0 with twin-range interactions");
313         }
314     }
315
316     if (ir->rlistlong == ir->rlist)
317     {
318         ir->nstcalclr = 0;
319     }
320     else if (ir->rlistlong > ir->rlist && ir->nstcalclr == 0)
321     {
322         warning_error(wi, "With different cutoffs for electrostatics and VdW, nstcalclr must be -1 or a positive number");
323     }
324
325     if (ir->cutoff_scheme == ecutsVERLET)
326     {
327         real rc_max;
328
329         /* Normal Verlet type neighbor-list, currently only limited feature support */
330         if (inputrec2nboundeddim(ir) < 3)
331         {
332             warning_error(wi, "With Verlet lists only full pbc or pbc=xy with walls is supported");
333         }
334         if (ir->rcoulomb != ir->rvdw)
335         {
336             warning_error(wi, "With Verlet lists rcoulomb!=rvdw is not supported");
337         }
338         if (ir->vdwtype == evdwSHIFT || ir->vdwtype == evdwSWITCH)
339         {
340             if (ir->vdw_modifier == eintmodNONE ||
341                 ir->vdw_modifier == eintmodPOTSHIFT)
342             {
343                 ir->vdw_modifier = (ir->vdwtype == evdwSHIFT ? eintmodFORCESWITCH : eintmodPOTSWITCH);
344
345                 sprintf(warn_buf, "Replacing vdwtype=%s by the equivalent combination of vdwtype=%s and vdw_modifier=%s", evdw_names[ir->vdwtype], evdw_names[evdwCUT], eintmod_names[ir->vdw_modifier]);
346                 warning_note(wi, warn_buf);
347
348                 ir->vdwtype = evdwCUT;
349             }
350             else
351             {
352                 sprintf(warn_buf, "Unsupported combination of vdwtype=%s and vdw_modifier=%s", evdw_names[ir->vdwtype], eintmod_names[ir->vdw_modifier]);
353                 warning_error(wi, warn_buf);
354             }
355         }
356
357         if (!(ir->vdwtype == evdwCUT || ir->vdwtype == evdwPME))
358         {
359             warning_error(wi, "With Verlet lists only cut-off and PME LJ interactions are supported");
360         }
361         if (!(ir->coulombtype == eelCUT ||
362               (EEL_RF(ir->coulombtype) && ir->coulombtype != eelRF_NEC) ||
363               EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD))
364         {
365             warning_error(wi, "With Verlet lists only cut-off, reaction-field, PME and Ewald electrostatics are supported");
366         }
367         if (!(ir->coulomb_modifier == eintmodNONE ||
368               ir->coulomb_modifier == eintmodPOTSHIFT))
369         {
370             sprintf(warn_buf, "coulomb_modifier=%s is not supported with the Verlet cut-off scheme", eintmod_names[ir->coulomb_modifier]);
371             warning_error(wi, warn_buf);
372         }
373
374         if (ir->implicit_solvent != eisNO)
375         {
376             warning_error(wi, "Implicit solvent is not (yet) supported with the with Verlet lists.");
377         }
378
379         if (ir->nstlist <= 0)
380         {
381             warning_error(wi, "With Verlet lists nstlist should be larger than 0");
382         }
383
384         if (ir->nstlist < 10)
385         {
386             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.");
387         }
388
389         rc_max = max(ir->rvdw, ir->rcoulomb);
390
391         if (ir->verletbuf_tol <= 0)
392         {
393             if (ir->verletbuf_tol == 0)
394             {
395                 warning_error(wi, "Can not have Verlet buffer tolerance of exactly 0");
396             }
397
398             if (ir->rlist < rc_max)
399             {
400                 warning_error(wi, "With verlet lists rlist can not be smaller than rvdw or rcoulomb");
401             }
402
403             if (ir->rlist == rc_max && ir->nstlist > 1)
404             {
405                 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.");
406             }
407         }
408         else
409         {
410             if (ir->rlist > rc_max)
411             {
412                 warning_note(wi, "You have set rlist larger than the interaction cut-off, but you also have verlet-buffer-tolerance > 0. Will set rlist using verlet-buffer-tolerance.");
413             }
414
415             if (ir->nstlist == 1)
416             {
417                 /* No buffer required */
418                 ir->rlist = rc_max;
419             }
420             else
421             {
422                 if (EI_DYNAMICS(ir->eI))
423                 {
424                     if (inputrec2nboundeddim(ir) < 3)
425                     {
426                         warning_error(wi, "The box volume is required for calculating rlist from the energy drift with verlet-buffer-tolerance > 0. You are using at least one unbounded dimension, so no volume can be computed. Either use a finite box, or set rlist yourself together with verlet-buffer-tolerance = -1.");
427                     }
428                     /* Set rlist temporarily so we can continue processing */
429                     ir->rlist = rc_max;
430                 }
431                 else
432                 {
433                     /* Set the buffer to 5% of the cut-off */
434                     ir->rlist = (1.0 + verlet_buffer_ratio_nodynamics)*rc_max;
435                 }
436             }
437         }
438
439         /* No twin-range calculations with Verlet lists */
440         ir->rlistlong = ir->rlist;
441     }
442
443     if (ir->nstcalclr == -1)
444     {
445         /* if rlist=rlistlong, this will later be changed to nstcalclr=0 */
446         ir->nstcalclr = ir->nstlist;
447     }
448     else if (ir->nstcalclr > 0)
449     {
450         if (ir->nstlist > 0 && (ir->nstlist % ir->nstcalclr != 0))
451         {
452             warning_error(wi, "nstlist must be evenly divisible by nstcalclr. Use nstcalclr = -1 to automatically follow nstlist");
453         }
454     }
455     else if (ir->nstcalclr < -1)
456     {
457         warning_error(wi, "nstcalclr must be a positive number (divisor of nstcalclr), or -1 to follow nstlist.");
458     }
459
460     if (EEL_PME(ir->coulombtype) && ir->rcoulomb > ir->rlist && ir->nstcalclr > 1)
461     {
462         warning_error(wi, "When used with PME, the long-range component of twin-range interactions must be updated every step (nstcalclr)");
463     }
464
465     /* GENERAL INTEGRATOR STUFF */
466     if (!(ir->eI == eiMD || EI_VV(ir->eI)))
467     {
468         ir->etc = etcNO;
469     }
470     if (ir->eI == eiVVAK)
471     {
472         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]);
473         warning_note(wi, warn_buf);
474     }
475     if (!EI_DYNAMICS(ir->eI))
476     {
477         ir->epc = epcNO;
478     }
479     if (EI_DYNAMICS(ir->eI))
480     {
481         if (ir->nstcalcenergy < 0)
482         {
483             ir->nstcalcenergy = ir_optimal_nstcalcenergy(ir);
484             if (ir->nstenergy != 0 && ir->nstenergy < ir->nstcalcenergy)
485             {
486                 /* nstcalcenergy larger than nstener does not make sense.
487                  * We ideally want nstcalcenergy=nstener.
488                  */
489                 if (ir->nstlist > 0)
490                 {
491                     ir->nstcalcenergy = lcd(ir->nstenergy, ir->nstlist);
492                 }
493                 else
494                 {
495                     ir->nstcalcenergy = ir->nstenergy;
496                 }
497             }
498         }
499         else if ( (ir->nstenergy > 0 && ir->nstcalcenergy > ir->nstenergy) ||
500                   (ir->efep != efepNO && ir->fepvals->nstdhdl > 0 &&
501                    (ir->nstcalcenergy > ir->fepvals->nstdhdl) ) )
502
503         {
504             const char *nsten    = "nstenergy";
505             const char *nstdh    = "nstdhdl";
506             const char *min_name = nsten;
507             int         min_nst  = ir->nstenergy;
508
509             /* find the smallest of ( nstenergy, nstdhdl ) */
510             if (ir->efep != efepNO && ir->fepvals->nstdhdl > 0 &&
511                 (ir->nstenergy == 0 || ir->fepvals->nstdhdl < ir->nstenergy))
512             {
513                 min_nst  = ir->fepvals->nstdhdl;
514                 min_name = nstdh;
515             }
516             /* If the user sets nstenergy small, we should respect that */
517             sprintf(warn_buf,
518                     "Setting nstcalcenergy (%d) equal to %s (%d)",
519                     ir->nstcalcenergy, min_name, min_nst);
520             warning_note(wi, warn_buf);
521             ir->nstcalcenergy = min_nst;
522         }
523
524         if (ir->epc != epcNO)
525         {
526             if (ir->nstpcouple < 0)
527             {
528                 ir->nstpcouple = ir_optimal_nstpcouple(ir);
529             }
530         }
531         if (IR_TWINRANGE(*ir))
532         {
533             check_nst("nstlist", ir->nstlist,
534                       "nstcalcenergy", &ir->nstcalcenergy, wi);
535             if (ir->epc != epcNO)
536             {
537                 check_nst("nstlist", ir->nstlist,
538                           "nstpcouple", &ir->nstpcouple, wi);
539             }
540         }
541
542         if (ir->nstcalcenergy > 0)
543         {
544             if (ir->efep != efepNO)
545             {
546                 /* nstdhdl should be a multiple of nstcalcenergy */
547                 check_nst("nstcalcenergy", ir->nstcalcenergy,
548                           "nstdhdl", &ir->fepvals->nstdhdl, wi);
549                 /* nstexpanded should be a multiple of nstcalcenergy */
550                 check_nst("nstcalcenergy", ir->nstcalcenergy,
551                           "nstexpanded", &ir->expandedvals->nstexpanded, wi);
552             }
553             /* for storing exact averages nstenergy should be
554              * a multiple of nstcalcenergy
555              */
556             check_nst("nstcalcenergy", ir->nstcalcenergy,
557                       "nstenergy", &ir->nstenergy, wi);
558         }
559     }
560
561     if (ir->nsteps == 0 && !ir->bContinuation)
562     {
563         warning_note(wi, "For a correct single-point energy evaluation with nsteps = 0, use continuation = yes to avoid constraining the input coordinates.");
564     }
565
566     /* LD STUFF */
567     if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
568         ir->bContinuation && ir->ld_seed != -1)
569     {
570         warning_note(wi, "You are doing a continuation with SD or BD, make sure that ld_seed is different from the previous run (using ld_seed=-1 will ensure this)");
571     }
572
573     /* TPI STUFF */
574     if (EI_TPI(ir->eI))
575     {
576         sprintf(err_buf, "TPI only works with pbc = %s", epbc_names[epbcXYZ]);
577         CHECK(ir->ePBC != epbcXYZ);
578         sprintf(err_buf, "TPI only works with ns = %s", ens_names[ensGRID]);
579         CHECK(ir->ns_type != ensGRID);
580         sprintf(err_buf, "with TPI nstlist should be larger than zero");
581         CHECK(ir->nstlist <= 0);
582         sprintf(err_buf, "TPI does not work with full electrostatics other than PME");
583         CHECK(EEL_FULL(ir->coulombtype) && !EEL_PME(ir->coulombtype));
584         sprintf(err_buf, "TPI does not work (yet) with the Verlet cut-off scheme");
585         CHECK(ir->cutoff_scheme == ecutsVERLET);
586     }
587
588     /* SHAKE / LINCS */
589     if ( (opts->nshake > 0) && (opts->bMorse) )
590     {
591         sprintf(warn_buf,
592                 "Using morse bond-potentials while constraining bonds is useless");
593         warning(wi, warn_buf);
594     }
595
596     if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
597         ir->bContinuation && ir->ld_seed != -1)
598     {
599         warning_note(wi, "You are doing a continuation with SD or BD, make sure that ld_seed is different from the previous run (using ld_seed=-1 will ensure this)");
600     }
601     /* verify simulated tempering options */
602
603     if (ir->bSimTemp)
604     {
605         gmx_bool bAllTempZero = TRUE;
606         for (i = 0; i < fep->n_lambda; i++)
607         {
608             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]);
609             CHECK((fep->all_lambda[efptTEMPERATURE][i] < 0) || (fep->all_lambda[efptTEMPERATURE][i] > 1));
610             if (fep->all_lambda[efptTEMPERATURE][i] > 0)
611             {
612                 bAllTempZero = FALSE;
613             }
614         }
615         sprintf(err_buf, "if simulated tempering is on, temperature-lambdas may not be all zero");
616         CHECK(bAllTempZero == TRUE);
617
618         sprintf(err_buf, "Simulated tempering is currently only compatible with md-vv");
619         CHECK(ir->eI != eiVV);
620
621         /* check compatability of the temperature coupling with simulated tempering */
622
623         if (ir->etc == etcNOSEHOOVER)
624         {
625             sprintf(warn_buf, "Nose-Hoover based temperature control such as [%s] my not be entirelyconsistent with simulated tempering", etcoupl_names[ir->etc]);
626             warning_note(wi, warn_buf);
627         }
628
629         /* check that the temperatures make sense */
630
631         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);
632         CHECK(ir->simtempvals->simtemp_high <= ir->simtempvals->simtemp_low);
633
634         sprintf(err_buf, "Higher simulated tempering temperature (%g) must be >= zero", ir->simtempvals->simtemp_high);
635         CHECK(ir->simtempvals->simtemp_high <= 0);
636
637         sprintf(err_buf, "Lower simulated tempering temperature (%g) must be >= zero", ir->simtempvals->simtemp_low);
638         CHECK(ir->simtempvals->simtemp_low <= 0);
639     }
640
641     /* verify free energy options */
642
643     if (ir->efep != efepNO)
644     {
645         fep = ir->fepvals;
646         sprintf(err_buf, "The soft-core power is %d and can only be 1 or 2",
647                 fep->sc_power);
648         CHECK(fep->sc_alpha != 0 && fep->sc_power != 1 && fep->sc_power != 2);
649
650         sprintf(err_buf, "The soft-core sc-r-power is %d and can only be 6 or 48",
651                 (int)fep->sc_r_power);
652         CHECK(fep->sc_alpha != 0 && fep->sc_r_power != 6.0 && fep->sc_r_power != 48.0);
653
654         sprintf(err_buf, "Can't use postive delta-lambda (%g) if initial state/lambda does not start at zero", fep->delta_lambda);
655         CHECK(fep->delta_lambda > 0 && ((fep->init_fep_state > 0) ||  (fep->init_lambda > 0)));
656
657         sprintf(err_buf, "Can't use postive delta-lambda (%g) with expanded ensemble simulations", fep->delta_lambda);
658         CHECK(fep->delta_lambda > 0 && (ir->efep == efepEXPANDED));
659
660         sprintf(err_buf, "Can only use expanded ensemble with md-vv (for now)");
661         CHECK(!(EI_VV(ir->eI)) && (ir->efep == efepEXPANDED));
662
663         sprintf(err_buf, "Free-energy not implemented for Ewald");
664         CHECK(ir->coulombtype == eelEWALD);
665
666         /* check validty of lambda inputs */
667         if (fep->n_lambda == 0)
668         {
669             /* Clear output in case of no states:*/
670             sprintf(err_buf, "init-lambda-state set to %d: no lambda states are defined.", fep->init_fep_state);
671             CHECK((fep->init_fep_state >= 0) && (fep->n_lambda == 0));
672         }
673         else
674         {
675             sprintf(err_buf, "initial thermodynamic state %d does not exist, only goes to %d", fep->init_fep_state, fep->n_lambda-1);
676             CHECK((fep->init_fep_state >= fep->n_lambda));
677         }
678
679         sprintf(err_buf, "Lambda state must be set, either with init-lambda-state or with init-lambda");
680         CHECK((fep->init_fep_state < 0) && (fep->init_lambda < 0));
681
682         sprintf(err_buf, "init-lambda=%g while init-lambda-state=%d. Lambda state must be set either with init-lambda-state or with init-lambda, but not both",
683                 fep->init_lambda, fep->init_fep_state);
684         CHECK((fep->init_fep_state >= 0) && (fep->init_lambda >= 0));
685
686
687
688         if ((fep->init_lambda >= 0) && (fep->delta_lambda == 0))
689         {
690             int n_lambda_terms;
691             n_lambda_terms = 0;
692             for (i = 0; i < efptNR; i++)
693             {
694                 if (fep->separate_dvdl[i])
695                 {
696                     n_lambda_terms++;
697                 }
698             }
699             if (n_lambda_terms > 1)
700             {
701                 sprintf(warn_buf, "If lambda vector states (fep-lambdas, coul-lambdas etc.) are set, don't use init-lambda to set lambda state (except for slow growth). Use init-lambda-state instead.");
702                 warning(wi, warn_buf);
703             }
704
705             if (n_lambda_terms < 2 && fep->n_lambda > 0)
706             {
707                 warning_note(wi,
708                              "init-lambda is deprecated for setting lambda state (except for slow growth). Use init-lambda-state instead.");
709             }
710         }
711
712         for (j = 0; j < efptNR; j++)
713         {
714             for (i = 0; i < fep->n_lambda; i++)
715             {
716                 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]);
717                 CHECK((fep->all_lambda[j][i] < 0) || (fep->all_lambda[j][i] > 1));
718             }
719         }
720
721         if ((fep->sc_alpha > 0) && (!fep->bScCoul))
722         {
723             for (i = 0; i < fep->n_lambda; i++)
724             {
725                 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],
726                         fep->all_lambda[efptCOUL][i]);
727                 CHECK((fep->sc_alpha > 0) &&
728                       (((fep->all_lambda[efptCOUL][i] > 0.0) &&
729                         (fep->all_lambda[efptCOUL][i] < 1.0)) &&
730                        ((fep->all_lambda[efptVDW][i] > 0.0) &&
731                         (fep->all_lambda[efptVDW][i] < 1.0))));
732             }
733         }
734
735         if ((fep->bScCoul) && (EEL_PME(ir->coulombtype)))
736         {
737             real sigma, lambda, r_sc;
738
739             sigma  = 0.34;
740             /* Maximum estimate for A and B charges equal with lambda power 1 */
741             lambda = 0.5;
742             r_sc   = pow(lambda*fep->sc_alpha*pow(sigma/ir->rcoulomb, fep->sc_r_power) + 1.0, 1.0/fep->sc_r_power);
743             sprintf(warn_buf, "With PME there is a minor soft core effect present at the cut-off, proportional to (LJsigma/rcoulomb)^%g. This could have a minor effect on energy conservation, but usually other effects dominate. With a common sigma value of %g nm the fraction of the particle-particle potential at the cut-off at lambda=%g is around %.1e, while ewald-rtol is %.1e.",
744                     fep->sc_r_power,
745                     sigma, lambda, r_sc - 1.0, ir->ewald_rtol);
746             warning_note(wi, warn_buf);
747         }
748
749         /*  Free Energy Checks -- In an ideal world, slow growth and FEP would
750             be treated differently, but that's the next step */
751
752         for (i = 0; i < efptNR; i++)
753         {
754             for (j = 0; j < fep->n_lambda; j++)
755             {
756                 sprintf(err_buf, "%s[%d] must be between 0 and 1", efpt_names[i], j);
757                 CHECK((fep->all_lambda[i][j] < 0) || (fep->all_lambda[i][j] > 1));
758             }
759         }
760     }
761
762     if ((ir->bSimTemp) || (ir->efep == efepEXPANDED))
763     {
764         fep    = ir->fepvals;
765         expand = ir->expandedvals;
766
767         /* checking equilibration of weights inputs for validity */
768
769         sprintf(err_buf, "weight-equil-number-all-lambda (%d) is ignored if lmc-weights-equil is not equal to %s",
770                 expand->equil_n_at_lam, elmceq_names[elmceqNUMATLAM]);
771         CHECK((expand->equil_n_at_lam > 0) && (expand->elmceq != elmceqNUMATLAM));
772
773         sprintf(err_buf, "weight-equil-number-samples (%d) is ignored if lmc-weights-equil is not equal to %s",
774                 expand->equil_samples, elmceq_names[elmceqSAMPLES]);
775         CHECK((expand->equil_samples > 0) && (expand->elmceq != elmceqSAMPLES));
776
777         sprintf(err_buf, "weight-equil-number-steps (%d) is ignored if lmc-weights-equil is not equal to %s",
778                 expand->equil_steps, elmceq_names[elmceqSTEPS]);
779         CHECK((expand->equil_steps > 0) && (expand->elmceq != elmceqSTEPS));
780
781         sprintf(err_buf, "weight-equil-wl-delta (%d) is ignored if lmc-weights-equil is not equal to %s",
782                 expand->equil_samples, elmceq_names[elmceqWLDELTA]);
783         CHECK((expand->equil_wl_delta > 0) && (expand->elmceq != elmceqWLDELTA));
784
785         sprintf(err_buf, "weight-equil-count-ratio (%f) is ignored if lmc-weights-equil is not equal to %s",
786                 expand->equil_ratio, elmceq_names[elmceqRATIO]);
787         CHECK((expand->equil_ratio > 0) && (expand->elmceq != elmceqRATIO));
788
789         sprintf(err_buf, "weight-equil-number-all-lambda (%d) must be a positive integer if lmc-weights-equil=%s",
790                 expand->equil_n_at_lam, elmceq_names[elmceqNUMATLAM]);
791         CHECK((expand->equil_n_at_lam <= 0) && (expand->elmceq == elmceqNUMATLAM));
792
793         sprintf(err_buf, "weight-equil-number-samples (%d) must be a positive integer if lmc-weights-equil=%s",
794                 expand->equil_samples, elmceq_names[elmceqSAMPLES]);
795         CHECK((expand->equil_samples <= 0) && (expand->elmceq == elmceqSAMPLES));
796
797         sprintf(err_buf, "weight-equil-number-steps (%d) must be a positive integer if lmc-weights-equil=%s",
798                 expand->equil_steps, elmceq_names[elmceqSTEPS]);
799         CHECK((expand->equil_steps <= 0) && (expand->elmceq == elmceqSTEPS));
800
801         sprintf(err_buf, "weight-equil-wl-delta (%f) must be > 0 if lmc-weights-equil=%s",
802                 expand->equil_wl_delta, elmceq_names[elmceqWLDELTA]);
803         CHECK((expand->equil_wl_delta <= 0) && (expand->elmceq == elmceqWLDELTA));
804
805         sprintf(err_buf, "weight-equil-count-ratio (%f) must be > 0 if lmc-weights-equil=%s",
806                 expand->equil_ratio, elmceq_names[elmceqRATIO]);
807         CHECK((expand->equil_ratio <= 0) && (expand->elmceq == elmceqRATIO));
808
809         sprintf(err_buf, "lmc-weights-equil=%s only possible when lmc-stats = %s or lmc-stats %s",
810                 elmceq_names[elmceqWLDELTA], elamstats_names[elamstatsWL], elamstats_names[elamstatsWWL]);
811         CHECK((expand->elmceq == elmceqWLDELTA) && (!EWL(expand->elamstats)));
812
813         sprintf(err_buf, "lmc-repeats (%d) must be greater than 0", expand->lmc_repeats);
814         CHECK((expand->lmc_repeats <= 0));
815         sprintf(err_buf, "minimum-var-min (%d) must be greater than 0", expand->minvarmin);
816         CHECK((expand->minvarmin <= 0));
817         sprintf(err_buf, "weight-c-range (%d) must be greater or equal to 0", expand->c_range);
818         CHECK((expand->c_range < 0));
819         sprintf(err_buf, "init-lambda-state (%d) must be zero if lmc-forced-nstart (%d)> 0 and lmc-move != 'no'",
820                 fep->init_fep_state, expand->lmc_forced_nstart);
821         CHECK((fep->init_fep_state != 0) && (expand->lmc_forced_nstart > 0) && (expand->elmcmove != elmcmoveNO));
822         sprintf(err_buf, "lmc-forced-nstart (%d) must not be negative", expand->lmc_forced_nstart);
823         CHECK((expand->lmc_forced_nstart < 0));
824         sprintf(err_buf, "init-lambda-state (%d) must be in the interval [0,number of lambdas)", fep->init_fep_state);
825         CHECK((fep->init_fep_state < 0) || (fep->init_fep_state >= fep->n_lambda));
826
827         sprintf(err_buf, "init-wl-delta (%f) must be greater than or equal to 0", expand->init_wl_delta);
828         CHECK((expand->init_wl_delta < 0));
829         sprintf(err_buf, "wl-ratio (%f) must be between 0 and 1", expand->wl_ratio);
830         CHECK((expand->wl_ratio <= 0) || (expand->wl_ratio >= 1));
831         sprintf(err_buf, "wl-scale (%f) must be between 0 and 1", expand->wl_scale);
832         CHECK((expand->wl_scale <= 0) || (expand->wl_scale >= 1));
833
834         /* if there is no temperature control, we need to specify an MC temperature */
835         sprintf(err_buf, "If there is no temperature control, and lmc-mcmove!= 'no',mc_temperature must be set to a positive number");
836         if (expand->nstTij > 0)
837         {
838             sprintf(err_buf, "nst-transition-matrix (%d) must be an integer multiple of nstlog (%d)",
839                     expand->nstTij, ir->nstlog);
840             CHECK((mod(expand->nstTij, ir->nstlog) != 0));
841         }
842     }
843
844     /* PBC/WALLS */
845     sprintf(err_buf, "walls only work with pbc=%s", epbc_names[epbcXY]);
846     CHECK(ir->nwall && ir->ePBC != epbcXY);
847
848     /* VACUUM STUFF */
849     if (ir->ePBC != epbcXYZ && ir->nwall != 2)
850     {
851         if (ir->ePBC == epbcNONE)
852         {
853             if (ir->epc != epcNO)
854             {
855                 warning(wi, "Turning off pressure coupling for vacuum system");
856                 ir->epc = epcNO;
857             }
858         }
859         else
860         {
861             sprintf(err_buf, "Can not have pressure coupling with pbc=%s",
862                     epbc_names[ir->ePBC]);
863             CHECK(ir->epc != epcNO);
864         }
865         sprintf(err_buf, "Can not have Ewald with pbc=%s", epbc_names[ir->ePBC]);
866         CHECK(EEL_FULL(ir->coulombtype));
867
868         sprintf(err_buf, "Can not have dispersion correction with pbc=%s",
869                 epbc_names[ir->ePBC]);
870         CHECK(ir->eDispCorr != edispcNO);
871     }
872
873     if (ir->rlist == 0.0)
874     {
875         sprintf(err_buf, "can only have neighborlist cut-off zero (=infinite)\n"
876                 "with coulombtype = %s or coulombtype = %s\n"
877                 "without periodic boundary conditions (pbc = %s) and\n"
878                 "rcoulomb and rvdw set to zero",
879                 eel_names[eelCUT], eel_names[eelUSER], epbc_names[epbcNONE]);
880         CHECK(((ir->coulombtype != eelCUT) && (ir->coulombtype != eelUSER)) ||
881               (ir->ePBC     != epbcNONE) ||
882               (ir->rcoulomb != 0.0)      || (ir->rvdw != 0.0));
883
884         if (ir->nstlist > 0)
885         {
886             warning_note(wi, "Simulating without cut-offs can be (slightly) faster with nstlist=0, nstype=simple and only one MPI rank");
887         }
888     }
889
890     /* COMM STUFF */
891     if (ir->nstcomm == 0)
892     {
893         ir->comm_mode = ecmNO;
894     }
895     if (ir->comm_mode != ecmNO)
896     {
897         if (ir->nstcomm < 0)
898         {
899             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");
900             ir->nstcomm = abs(ir->nstcomm);
901         }
902
903         if (ir->nstcalcenergy > 0 && ir->nstcomm < ir->nstcalcenergy)
904         {
905             warning_note(wi, "nstcomm < nstcalcenergy defeats the purpose of nstcalcenergy, setting nstcomm to nstcalcenergy");
906             ir->nstcomm = ir->nstcalcenergy;
907         }
908
909         if (ir->comm_mode == ecmANGULAR)
910         {
911             sprintf(err_buf, "Can not remove the rotation around the center of mass with periodic molecules");
912             CHECK(ir->bPeriodicMols);
913             if (ir->ePBC != epbcNONE)
914             {
915                 warning(wi, "Removing the rotation around the center of mass in a periodic system, this can lead to artifacts. Only use this on a single (cluster of) molecules. This cluster should not cross periodic boundaries.");
916             }
917         }
918     }
919
920     if (EI_STATE_VELOCITY(ir->eI) && ir->ePBC == epbcNONE && ir->comm_mode != ecmANGULAR)
921     {
922         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.");
923     }
924
925     sprintf(err_buf, "Twin-range neighbour searching (NS) with simple NS"
926             " algorithm not implemented");
927     CHECK(((ir->rcoulomb > ir->rlist) || (ir->rvdw > ir->rlist))
928           && (ir->ns_type == ensSIMPLE));
929
930     /* TEMPERATURE COUPLING */
931     if (ir->etc == etcYES)
932     {
933         ir->etc = etcBERENDSEN;
934         warning_note(wi, "Old option for temperature coupling given: "
935                      "changing \"yes\" to \"Berendsen\"\n");
936     }
937
938     if ((ir->etc == etcNOSEHOOVER) || (ir->epc == epcMTTK))
939     {
940         if (ir->opts.nhchainlength < 1)
941         {
942             sprintf(warn_buf, "number of Nose-Hoover chains (currently %d) cannot be less than 1,reset to 1\n", ir->opts.nhchainlength);
943             ir->opts.nhchainlength = 1;
944             warning(wi, warn_buf);
945         }
946
947         if (ir->etc == etcNOSEHOOVER && !EI_VV(ir->eI) && ir->opts.nhchainlength > 1)
948         {
949             warning_note(wi, "leapfrog does not yet support Nose-Hoover chains, nhchainlength reset to 1");
950             ir->opts.nhchainlength = 1;
951         }
952     }
953     else
954     {
955         ir->opts.nhchainlength = 0;
956     }
957
958     if (ir->eI == eiVVAK)
959     {
960         sprintf(err_buf, "%s implemented primarily for validation, and requires nsttcouple = 1 and nstpcouple = 1.",
961                 ei_names[eiVVAK]);
962         CHECK((ir->nsttcouple != 1) || (ir->nstpcouple != 1));
963     }
964
965     if (ETC_ANDERSEN(ir->etc))
966     {
967         sprintf(err_buf, "%s temperature control not supported for integrator %s.", etcoupl_names[ir->etc], ei_names[ir->eI]);
968         CHECK(!(EI_VV(ir->eI)));
969
970         if (ir->nstcomm > 0 && (ir->etc == etcANDERSEN))
971         {
972             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]);
973             warning_note(wi, warn_buf);
974         }
975
976         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]);
977         CHECK(ir->nstcomm > 1 && (ir->etc == etcANDERSEN));
978     }
979
980     if (ir->etc == etcBERENDSEN)
981     {
982         sprintf(warn_buf, "The %s thermostat does not generate the correct kinetic energy distribution. You might want to consider using the %s thermostat.",
983                 ETCOUPLTYPE(ir->etc), ETCOUPLTYPE(etcVRESCALE));
984         warning_note(wi, warn_buf);
985     }
986
987     if ((ir->etc == etcNOSEHOOVER || ETC_ANDERSEN(ir->etc))
988         && ir->epc == epcBERENDSEN)
989     {
990         sprintf(warn_buf, "Using Berendsen pressure coupling invalidates the "
991                 "true ensemble for the thermostat");
992         warning(wi, warn_buf);
993     }
994
995     /* PRESSURE COUPLING */
996     if (ir->epc == epcISOTROPIC)
997     {
998         ir->epc = epcBERENDSEN;
999         warning_note(wi, "Old option for pressure coupling given: "
1000                      "changing \"Isotropic\" to \"Berendsen\"\n");
1001     }
1002
1003     if (ir->epc != epcNO)
1004     {
1005         dt_pcoupl = ir->nstpcouple*ir->delta_t;
1006
1007         sprintf(err_buf, "tau-p must be > 0 instead of %g\n", ir->tau_p);
1008         CHECK(ir->tau_p <= 0);
1009
1010         if (ir->tau_p/dt_pcoupl < pcouple_min_integration_steps(ir->epc) - 10*GMX_REAL_EPS)
1011         {
1012             sprintf(warn_buf, "For proper integration of the %s barostat, tau-p (%g) should be at least %d times larger than nstpcouple*dt (%g)",
1013                     EPCOUPLTYPE(ir->epc), ir->tau_p, pcouple_min_integration_steps(ir->epc), dt_pcoupl);
1014             warning(wi, warn_buf);
1015         }
1016
1017         sprintf(err_buf, "compressibility must be > 0 when using pressure"
1018                 " coupling %s\n", EPCOUPLTYPE(ir->epc));
1019         CHECK(ir->compress[XX][XX] < 0 || ir->compress[YY][YY] < 0 ||
1020               ir->compress[ZZ][ZZ] < 0 ||
1021               (trace(ir->compress) == 0 && ir->compress[YY][XX] <= 0 &&
1022                ir->compress[ZZ][XX] <= 0 && ir->compress[ZZ][YY] <= 0));
1023
1024         if (epcPARRINELLORAHMAN == ir->epc && opts->bGenVel)
1025         {
1026             sprintf(warn_buf,
1027                     "You are generating velocities so I am assuming you "
1028                     "are equilibrating a system. You are using "
1029                     "%s pressure coupling, but this can be "
1030                     "unstable for equilibration. If your system crashes, try "
1031                     "equilibrating first with Berendsen pressure coupling. If "
1032                     "you are not equilibrating the system, you can probably "
1033                     "ignore this warning.",
1034                     epcoupl_names[ir->epc]);
1035             warning(wi, warn_buf);
1036         }
1037     }
1038
1039     if (EI_VV(ir->eI))
1040     {
1041         if (ir->epc > epcNO)
1042         {
1043             if ((ir->epc != epcBERENDSEN) && (ir->epc != epcMTTK))
1044             {
1045                 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.");
1046             }
1047         }
1048     }
1049     else
1050     {
1051         if (ir->epc == epcMTTK)
1052         {
1053             warning_error(wi, "MTTK pressure coupling requires a Velocity-verlet integrator");
1054         }
1055     }
1056
1057     /* ELECTROSTATICS */
1058     /* More checks are in triple check (grompp.c) */
1059
1060     if (ir->coulombtype == eelSWITCH)
1061     {
1062         sprintf(warn_buf, "coulombtype = %s is only for testing purposes and can lead to serious "
1063                 "artifacts, advice: use coulombtype = %s",
1064                 eel_names[ir->coulombtype],
1065                 eel_names[eelRF_ZERO]);
1066         warning(wi, warn_buf);
1067     }
1068
1069     if (ir->epsilon_r != 1 && ir->implicit_solvent == eisGBSA)
1070     {
1071         sprintf(warn_buf, "epsilon-r = %g with GB implicit solvent, will use this value for inner dielectric", ir->epsilon_r);
1072         warning_note(wi, warn_buf);
1073     }
1074
1075     if (EEL_RF(ir->coulombtype) && ir->epsilon_rf == 1 && ir->epsilon_r != 1)
1076     {
1077         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);
1078         warning(wi, warn_buf);
1079         ir->epsilon_rf = ir->epsilon_r;
1080         ir->epsilon_r  = 1.0;
1081     }
1082
1083     if (getenv("GMX_DO_GALACTIC_DYNAMICS") == NULL)
1084     {
1085         sprintf(err_buf, "epsilon-r must be >= 0 instead of %g\n", ir->epsilon_r);
1086         CHECK(ir->epsilon_r < 0);
1087     }
1088
1089     if (EEL_RF(ir->coulombtype))
1090     {
1091         /* reaction field (at the cut-off) */
1092
1093         if (ir->coulombtype == eelRF_ZERO)
1094         {
1095             sprintf(warn_buf, "With coulombtype = %s, epsilon-rf must be 0, assuming you meant epsilon_rf=0",
1096                     eel_names[ir->coulombtype]);
1097             CHECK(ir->epsilon_rf != 0);
1098             ir->epsilon_rf = 0.0;
1099         }
1100
1101         sprintf(err_buf, "epsilon-rf must be >= epsilon-r");
1102         CHECK((ir->epsilon_rf < ir->epsilon_r && ir->epsilon_rf != 0) ||
1103               (ir->epsilon_r == 0));
1104         if (ir->epsilon_rf == ir->epsilon_r)
1105         {
1106             sprintf(warn_buf, "Using epsilon-rf = epsilon-r with %s does not make sense",
1107                     eel_names[ir->coulombtype]);
1108             warning(wi, warn_buf);
1109         }
1110     }
1111     /* Allow rlist>rcoulomb for tabulated long range stuff. This just
1112      * means the interaction is zero outside rcoulomb, but it helps to
1113      * provide accurate energy conservation.
1114      */
1115     if (ir_coulomb_might_be_zero_at_cutoff(ir))
1116     {
1117         if (ir_coulomb_switched(ir))
1118         {
1119             sprintf(err_buf,
1120                     "With coulombtype = %s rcoulomb_switch must be < rcoulomb. Or, better: Use the potential modifier options!",
1121                     eel_names[ir->coulombtype]);
1122             CHECK(ir->rcoulomb_switch >= ir->rcoulomb);
1123         }
1124     }
1125     else if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
1126     {
1127         if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE)
1128         {
1129             sprintf(err_buf, "With coulombtype = %s, rcoulomb should be >= rlist unless you use a potential modifier",
1130                     eel_names[ir->coulombtype]);
1131             CHECK(ir->rlist > ir->rcoulomb);
1132         }
1133     }
1134
1135     if (ir->coulombtype == eelSWITCH || ir->coulombtype == eelSHIFT)
1136     {
1137         sprintf(err_buf,
1138                 "Explicit switch/shift coulomb interactions cannot be used in combination with a secondary coulomb-modifier.");
1139         CHECK( ir->coulomb_modifier != eintmodNONE);
1140     }
1141     if (ir->vdwtype == evdwSWITCH || ir->vdwtype == evdwSHIFT)
1142     {
1143         sprintf(err_buf,
1144                 "Explicit switch/shift vdw interactions cannot be used in combination with a secondary vdw-modifier.");
1145         CHECK( ir->vdw_modifier != eintmodNONE);
1146     }
1147
1148     if (ir->coulombtype == eelSWITCH || ir->coulombtype == eelSHIFT ||
1149         ir->vdwtype == evdwSWITCH || ir->vdwtype == evdwSHIFT)
1150     {
1151         sprintf(warn_buf,
1152                 "The switch/shift interaction settings are just for compatibility; you will get better "
1153                 "performance from applying potential modifiers to your interactions!\n");
1154         warning_note(wi, warn_buf);
1155     }
1156
1157     if (ir->coulombtype == eelPMESWITCH || ir->coulomb_modifier == eintmodPOTSWITCH)
1158     {
1159         if (ir->rcoulomb_switch/ir->rcoulomb < 0.9499)
1160         {
1161             real percentage  = 100*(ir->rcoulomb-ir->rcoulomb_switch)/ir->rcoulomb;
1162             sprintf(warn_buf, "The switching range should be 5%% or less (currently %.2f%% using a switching range of %4f-%4f) for accurate electrostatic energies, energy conservation will be good regardless, since ewald_rtol = %g.",
1163                     percentage, ir->rcoulomb_switch, ir->rcoulomb, ir->ewald_rtol);
1164             warning(wi, warn_buf);
1165         }
1166     }
1167
1168     if (ir->vdwtype == evdwSWITCH || ir->vdw_modifier == eintmodPOTSWITCH)
1169     {
1170         if (ir->rvdw_switch == 0)
1171         {
1172             sprintf(warn_buf, "rvdw-switch is equal 0 even though you are using a switched Lennard-Jones potential.  This suggests it was not set in the mdp, which can lead to large energy errors.  In GROMACS, 0.05 to 0.1 nm is often a reasonable vdw switching range.");
1173             warning(wi, warn_buf);
1174         }
1175     }
1176
1177     if (EEL_FULL(ir->coulombtype))
1178     {
1179         if (ir->coulombtype == eelPMESWITCH || ir->coulombtype == eelPMEUSER ||
1180             ir->coulombtype == eelPMEUSERSWITCH)
1181         {
1182             sprintf(err_buf, "With coulombtype = %s, rcoulomb must be <= rlist",
1183                     eel_names[ir->coulombtype]);
1184             CHECK(ir->rcoulomb > ir->rlist);
1185         }
1186         else if (ir->cutoff_scheme == ecutsGROUP && ir->coulomb_modifier == eintmodNONE)
1187         {
1188             if (ir->coulombtype == eelPME || ir->coulombtype == eelP3M_AD)
1189             {
1190                 sprintf(err_buf,
1191                         "With coulombtype = %s (without modifier), rcoulomb must be equal to rlist,\n"
1192                         "or rlistlong if nstcalclr=1. For optimal energy conservation,consider using\n"
1193                         "a potential modifier.", eel_names[ir->coulombtype]);
1194                 if (ir->nstcalclr == 1)
1195                 {
1196                     CHECK(ir->rcoulomb != ir->rlist && ir->rcoulomb != ir->rlistlong);
1197                 }
1198                 else
1199                 {
1200                     CHECK(ir->rcoulomb != ir->rlist);
1201                 }
1202             }
1203         }
1204     }
1205
1206     if (EEL_PME(ir->coulombtype) || EVDW_PME(ir->vdwtype))
1207     {
1208         if (ir->pme_order < 3)
1209         {
1210             warning_error(wi, "pme-order can not be smaller than 3");
1211         }
1212     }
1213
1214     if (ir->nwall == 2 && EEL_FULL(ir->coulombtype))
1215     {
1216         if (ir->ewald_geometry == eewg3D)
1217         {
1218             sprintf(warn_buf, "With pbc=%s you should use ewald-geometry=%s",
1219                     epbc_names[ir->ePBC], eewg_names[eewg3DC]);
1220             warning(wi, warn_buf);
1221         }
1222         /* This check avoids extra pbc coding for exclusion corrections */
1223         sprintf(err_buf, "wall-ewald-zfac should be >= 2");
1224         CHECK(ir->wall_ewald_zfac < 2);
1225     }
1226     if ((ir->ewald_geometry == eewg3DC) && (ir->ePBC != epbcXY) &&
1227         EEL_FULL(ir->coulombtype))
1228     {
1229         sprintf(warn_buf, "With %s and ewald_geometry = %s you should use pbc = %s",
1230                 eel_names[ir->coulombtype], eewg_names[eewg3DC], epbc_names[epbcXY]);
1231         warning(wi, warn_buf);
1232     }
1233     if ((ir->epsilon_surface != 0) && EEL_FULL(ir->coulombtype))
1234     {
1235         if (ir->cutoff_scheme == ecutsVERLET)
1236         {
1237             sprintf(warn_buf, "Since molecules/charge groups are broken using the Verlet scheme, you can not use a dipole correction to the %s electrostatics.",
1238                     eel_names[ir->coulombtype]);
1239             warning(wi, warn_buf);
1240         }
1241         else
1242         {
1243             sprintf(warn_buf, "Dipole corrections to %s electrostatics only work if all charge groups that can cross PBC boundaries are dipoles. If this is not the case set epsilon_surface to 0",
1244                     eel_names[ir->coulombtype]);
1245             warning_note(wi, warn_buf);
1246         }
1247     }
1248
1249     if (ir_vdw_switched(ir))
1250     {
1251         sprintf(err_buf, "With switched vdw forces or potentials, rvdw-switch must be < rvdw");
1252         CHECK(ir->rvdw_switch >= ir->rvdw);
1253
1254         if (ir->rvdw_switch < 0.5*ir->rvdw)
1255         {
1256             sprintf(warn_buf, "You are applying a switch function to vdw forces or potentials from %g to %g nm, which is more than half the interaction range, whereas switch functions are intended to act only close to the cut-off.",
1257                     ir->rvdw_switch, ir->rvdw);
1258             warning_note(wi, warn_buf);
1259         }
1260     }
1261     else if (ir->vdwtype == evdwCUT || ir->vdwtype == evdwPME)
1262     {
1263         if (ir->cutoff_scheme == ecutsGROUP && ir->vdw_modifier == eintmodNONE)
1264         {
1265             sprintf(err_buf, "With vdwtype = %s, rvdw must be >= rlist unless you use a potential modifier", evdw_names[ir->vdwtype]);
1266             CHECK(ir->rlist > ir->rvdw);
1267         }
1268     }
1269
1270     if (ir->vdwtype == evdwPME)
1271     {
1272         if (!(ir->vdw_modifier == eintmodNONE || ir->vdw_modifier == eintmodPOTSHIFT))
1273         {
1274             sprintf(err_buf, "With vdwtype = %s, the only supported modifiers are %s a\
1275 nd %s",
1276                     evdw_names[ir->vdwtype],
1277                     eintmod_names[eintmodPOTSHIFT],
1278                     eintmod_names[eintmodNONE]);
1279         }
1280     }
1281
1282     if (ir->cutoff_scheme == ecutsGROUP)
1283     {
1284         if (((ir->coulomb_modifier != eintmodNONE && ir->rcoulomb == ir->rlist) ||
1285              (ir->vdw_modifier != eintmodNONE && ir->rvdw == ir->rlist)))
1286         {
1287             warning_note(wi, "With exact cut-offs, rlist should be "
1288                          "larger than rcoulomb and rvdw, so that there "
1289                          "is a buffer region for particle motion "
1290                          "between neighborsearch steps");
1291         }
1292
1293         if (ir_coulomb_is_zero_at_cutoff(ir) && ir->rlistlong <= ir->rcoulomb)
1294         {
1295             sprintf(warn_buf, "For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rcoulomb.",
1296                     IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
1297             warning_note(wi, warn_buf);
1298         }
1299         if (ir_vdw_switched(ir) && (ir->rlistlong <= ir->rvdw))
1300         {
1301             sprintf(warn_buf, "For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rvdw.",
1302                     IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
1303             warning_note(wi, warn_buf);
1304         }
1305     }
1306
1307     if (ir->vdwtype == evdwUSER && ir->eDispCorr != edispcNO)
1308     {
1309         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.");
1310     }
1311
1312     if (ir->eI == eiLBFGS && (ir->coulombtype == eelCUT || ir->vdwtype == evdwCUT)
1313         && ir->rvdw != 0)
1314     {
1315         warning(wi, "For efficient BFGS minimization, use switch/shift/pme instead of cut-off.");
1316     }
1317
1318     if (ir->eI == eiLBFGS && ir->nbfgscorr <= 0)
1319     {
1320         warning(wi, "Using L-BFGS with nbfgscorr<=0 just gets you steepest descent.");
1321     }
1322
1323     /* ENERGY CONSERVATION */
1324     if (ir_NVE(ir) && ir->cutoff_scheme == ecutsGROUP)
1325     {
1326         if (!ir_vdw_might_be_zero_at_cutoff(ir) && ir->rvdw > 0 && ir->vdw_modifier == eintmodNONE)
1327         {
1328             sprintf(warn_buf, "You are using a cut-off for VdW interactions with NVE, for good energy conservation use vdwtype = %s (possibly with DispCorr)",
1329                     evdw_names[evdwSHIFT]);
1330             warning_note(wi, warn_buf);
1331         }
1332         if (!ir_coulomb_might_be_zero_at_cutoff(ir) && ir->rcoulomb > 0)
1333         {
1334             sprintf(warn_buf, "You are using a cut-off for electrostatics with NVE, for good energy conservation use coulombtype = %s or %s",
1335                     eel_names[eelPMESWITCH], eel_names[eelRF_ZERO]);
1336             warning_note(wi, warn_buf);
1337         }
1338     }
1339
1340     if (EI_VV(ir->eI) && IR_TWINRANGE(*ir) && ir->nstlist > 1)
1341     {
1342         sprintf(warn_buf, "Twin-range multiple time stepping does not work with integrator %s.", ei_names[ir->eI]);
1343         warning_error(wi, warn_buf);
1344     }
1345
1346     /* IMPLICIT SOLVENT */
1347     if (ir->coulombtype == eelGB_NOTUSED)
1348     {
1349         sprintf(warn_buf, "Invalid option %s for coulombtype",
1350                 eel_names[ir->coulombtype]);
1351         warning_error(wi, warn_buf);
1352     }
1353
1354     if (ir->sa_algorithm == esaSTILL)
1355     {
1356         sprintf(err_buf, "Still SA algorithm not available yet, use %s or %s instead\n", esa_names[esaAPPROX], esa_names[esaNO]);
1357         CHECK(ir->sa_algorithm == esaSTILL);
1358     }
1359
1360     if (ir->implicit_solvent == eisGBSA)
1361     {
1362         sprintf(err_buf, "With GBSA implicit solvent, rgbradii must be equal to rlist.");
1363         CHECK(ir->rgbradii != ir->rlist);
1364
1365         if (ir->coulombtype != eelCUT)
1366         {
1367             sprintf(err_buf, "With GBSA, coulombtype must be equal to %s\n", eel_names[eelCUT]);
1368             CHECK(ir->coulombtype != eelCUT);
1369         }
1370         if (ir->vdwtype != evdwCUT)
1371         {
1372             sprintf(err_buf, "With GBSA, vdw-type must be equal to %s\n", evdw_names[evdwCUT]);
1373             CHECK(ir->vdwtype != evdwCUT);
1374         }
1375         if (ir->nstgbradii < 1)
1376         {
1377             sprintf(warn_buf, "Using GBSA with nstgbradii<1, setting nstgbradii=1");
1378             warning_note(wi, warn_buf);
1379             ir->nstgbradii = 1;
1380         }
1381         if (ir->sa_algorithm == esaNO)
1382         {
1383             sprintf(warn_buf, "No SA (non-polar) calculation requested together with GB. Are you sure this is what you want?\n");
1384             warning_note(wi, warn_buf);
1385         }
1386         if (ir->sa_surface_tension < 0 && ir->sa_algorithm != esaNO)
1387         {
1388             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");
1389             warning_note(wi, warn_buf);
1390
1391             if (ir->gb_algorithm == egbSTILL)
1392             {
1393                 ir->sa_surface_tension = 0.0049 * CAL2JOULE * 100;
1394             }
1395             else
1396             {
1397                 ir->sa_surface_tension = 0.0054 * CAL2JOULE * 100;
1398             }
1399         }
1400         if (ir->sa_surface_tension == 0 && ir->sa_algorithm != esaNO)
1401         {
1402             sprintf(err_buf, "Surface tension set to 0 while SA-calculation requested\n");
1403             CHECK(ir->sa_surface_tension == 0 && ir->sa_algorithm != esaNO);
1404         }
1405
1406     }
1407
1408     if (ir->bAdress)
1409     {
1410         if (ir->cutoff_scheme != ecutsGROUP)
1411         {
1412             warning_error(wi, "AdresS simulation supports only cutoff-scheme=group");
1413         }
1414         if (!EI_SD(ir->eI))
1415         {
1416             warning_error(wi, "AdresS simulation supports only stochastic dynamics");
1417         }
1418         if (ir->epc != epcNO)
1419         {
1420             warning_error(wi, "AdresS simulation does not support pressure coupling");
1421         }
1422         if (EEL_FULL(ir->coulombtype))
1423         {
1424             warning_error(wi, "AdresS simulation does not support long-range electrostatics");
1425         }
1426     }
1427 }
1428
1429 /* count the number of text elemets separated by whitespace in a string.
1430     str = the input string
1431     maxptr = the maximum number of allowed elements
1432     ptr = the output array of pointers to the first character of each element
1433     returns: the number of elements. */
1434 int str_nelem(const char *str, int maxptr, char *ptr[])
1435 {
1436     int   np = 0;
1437     char *copy0, *copy;
1438
1439     copy0 = gmx_strdup(str);
1440     copy  = copy0;
1441     ltrim(copy);
1442     while (*copy != '\0')
1443     {
1444         if (np >= maxptr)
1445         {
1446             gmx_fatal(FARGS, "Too many groups on line: '%s' (max is %d)",
1447                       str, maxptr);
1448         }
1449         if (ptr)
1450         {
1451             ptr[np] = copy;
1452         }
1453         np++;
1454         while ((*copy != '\0') && !isspace(*copy))
1455         {
1456             copy++;
1457         }
1458         if (*copy != '\0')
1459         {
1460             *copy = '\0';
1461             copy++;
1462         }
1463         ltrim(copy);
1464     }
1465     if (ptr == NULL)
1466     {
1467         sfree(copy0);
1468     }
1469
1470     return np;
1471 }
1472
1473 /* interpret a number of doubles from a string and put them in an array,
1474    after allocating space for them.
1475    str = the input string
1476    n = the (pre-allocated) number of doubles read
1477    r = the output array of doubles. */
1478 static void parse_n_real(char *str, int *n, real **r)
1479 {
1480     char *ptr[MAXPTR];
1481     int   i;
1482
1483     *n = str_nelem(str, MAXPTR, ptr);
1484
1485     snew(*r, *n);
1486     for (i = 0; i < *n; i++)
1487     {
1488         (*r)[i] = strtod(ptr[i], NULL);
1489     }
1490 }
1491
1492 static void do_fep_params(t_inputrec *ir, char fep_lambda[][STRLEN], char weights[STRLEN])
1493 {
1494
1495     int         i, j, max_n_lambda, nweights, nfep[efptNR];
1496     t_lambda   *fep    = ir->fepvals;
1497     t_expanded *expand = ir->expandedvals;
1498     real      **count_fep_lambdas;
1499     gmx_bool    bOneLambda = TRUE;
1500
1501     snew(count_fep_lambdas, efptNR);
1502
1503     /* FEP input processing */
1504     /* first, identify the number of lambda values for each type.
1505        All that are nonzero must have the same number */
1506
1507     for (i = 0; i < efptNR; i++)
1508     {
1509         parse_n_real(fep_lambda[i], &(nfep[i]), &(count_fep_lambdas[i]));
1510     }
1511
1512     /* now, determine the number of components.  All must be either zero, or equal. */
1513
1514     max_n_lambda = 0;
1515     for (i = 0; i < efptNR; i++)
1516     {
1517         if (nfep[i] > max_n_lambda)
1518         {
1519             max_n_lambda = nfep[i];  /* here's a nonzero one.  All of them
1520                                         must have the same number if its not zero.*/
1521             break;
1522         }
1523     }
1524
1525     for (i = 0; i < efptNR; i++)
1526     {
1527         if (nfep[i] == 0)
1528         {
1529             ir->fepvals->separate_dvdl[i] = FALSE;
1530         }
1531         else if (nfep[i] == max_n_lambda)
1532         {
1533             if (i != efptTEMPERATURE)  /* we treat this differently -- not really a reason to compute the derivative with
1534                                           respect to the temperature currently */
1535             {
1536                 ir->fepvals->separate_dvdl[i] = TRUE;
1537             }
1538         }
1539         else
1540         {
1541             gmx_fatal(FARGS, "Number of lambdas (%d) for FEP type %s not equal to number of other types (%d)",
1542                       nfep[i], efpt_names[i], max_n_lambda);
1543         }
1544     }
1545     /* we don't print out dhdl if the temperature is changing, since we can't correctly define dhdl in this case */
1546     ir->fepvals->separate_dvdl[efptTEMPERATURE] = FALSE;
1547
1548     /* the number of lambdas is the number we've read in, which is either zero
1549        or the same for all */
1550     fep->n_lambda = max_n_lambda;
1551
1552     /* allocate space for the array of lambda values */
1553     snew(fep->all_lambda, efptNR);
1554     /* if init_lambda is defined, we need to set lambda */
1555     if ((fep->init_lambda > 0) && (fep->n_lambda == 0))
1556     {
1557         ir->fepvals->separate_dvdl[efptFEP] = TRUE;
1558     }
1559     /* otherwise allocate the space for all of the lambdas, and transfer the data */
1560     for (i = 0; i < efptNR; i++)
1561     {
1562         snew(fep->all_lambda[i], fep->n_lambda);
1563         if (nfep[i] > 0)  /* if it's zero, then the count_fep_lambda arrays
1564                              are zero */
1565         {
1566             for (j = 0; j < fep->n_lambda; j++)
1567             {
1568                 fep->all_lambda[i][j] = (double)count_fep_lambdas[i][j];
1569             }
1570             sfree(count_fep_lambdas[i]);
1571         }
1572     }
1573     sfree(count_fep_lambdas);
1574
1575     /* "fep-vals" is either zero or the full number. If zero, we'll need to define fep-lambdas for internal
1576        bookkeeping -- for now, init_lambda */
1577
1578     if ((nfep[efptFEP] == 0) && (fep->init_lambda >= 0))
1579     {
1580         for (i = 0; i < fep->n_lambda; i++)
1581         {
1582             fep->all_lambda[efptFEP][i] = fep->init_lambda;
1583         }
1584     }
1585
1586     /* check to see if only a single component lambda is defined, and soft core is defined.
1587        In this case, turn on coulomb soft core */
1588
1589     if (max_n_lambda == 0)
1590     {
1591         bOneLambda = TRUE;
1592     }
1593     else
1594     {
1595         for (i = 0; i < efptNR; i++)
1596         {
1597             if ((nfep[i] != 0) && (i != efptFEP))
1598             {
1599                 bOneLambda = FALSE;
1600             }
1601         }
1602     }
1603     if ((bOneLambda) && (fep->sc_alpha > 0))
1604     {
1605         fep->bScCoul = TRUE;
1606     }
1607
1608     /* Fill in the others with the efptFEP if they are not explicitly
1609        specified (i.e. nfep[i] == 0).  This means if fep is not defined,
1610        they are all zero. */
1611
1612     for (i = 0; i < efptNR; i++)
1613     {
1614         if ((nfep[i] == 0) && (i != efptFEP))
1615         {
1616             for (j = 0; j < fep->n_lambda; j++)
1617             {
1618                 fep->all_lambda[i][j] = fep->all_lambda[efptFEP][j];
1619             }
1620         }
1621     }
1622
1623
1624     /* make it easier if sc_r_power = 48 by increasing it to the 4th power, to be in the right scale. */
1625     if (fep->sc_r_power == 48)
1626     {
1627         if (fep->sc_alpha > 0.1)
1628         {
1629             gmx_fatal(FARGS, "sc_alpha (%f) for sc_r_power = 48 should usually be between 0.001 and 0.004", fep->sc_alpha);
1630         }
1631     }
1632
1633     expand = ir->expandedvals;
1634     /* now read in the weights */
1635     parse_n_real(weights, &nweights, &(expand->init_lambda_weights));
1636     if (nweights == 0)
1637     {
1638         snew(expand->init_lambda_weights, fep->n_lambda); /* initialize to zero */
1639     }
1640     else if (nweights != fep->n_lambda)
1641     {
1642         gmx_fatal(FARGS, "Number of weights (%d) is not equal to number of lambda values (%d)",
1643                   nweights, fep->n_lambda);
1644     }
1645     if ((expand->nstexpanded < 0) && (ir->efep != efepNO))
1646     {
1647         expand->nstexpanded = fep->nstdhdl;
1648         /* if you don't specify nstexpanded when doing expanded ensemble free energy calcs, it is set to nstdhdl */
1649     }
1650     if ((expand->nstexpanded < 0) && ir->bSimTemp)
1651     {
1652         expand->nstexpanded = 2*(int)(ir->opts.tau_t[0]/ir->delta_t);
1653         /* if you don't specify nstexpanded when doing expanded ensemble simulated tempering, it is set to
1654            2*tau_t just to be careful so it's not to frequent  */
1655     }
1656 }
1657
1658
1659 static void do_simtemp_params(t_inputrec *ir)
1660 {
1661
1662     snew(ir->simtempvals->temperatures, ir->fepvals->n_lambda);
1663     GetSimTemps(ir->fepvals->n_lambda, ir->simtempvals, ir->fepvals->all_lambda[efptTEMPERATURE]);
1664
1665     return;
1666 }
1667
1668 static void do_wall_params(t_inputrec *ir,
1669                            char *wall_atomtype, char *wall_density,
1670                            t_gromppopts *opts)
1671 {
1672     int    nstr, i;
1673     char  *names[MAXPTR];
1674     double dbl;
1675
1676     opts->wall_atomtype[0] = NULL;
1677     opts->wall_atomtype[1] = NULL;
1678
1679     ir->wall_atomtype[0] = -1;
1680     ir->wall_atomtype[1] = -1;
1681     ir->wall_density[0]  = 0;
1682     ir->wall_density[1]  = 0;
1683
1684     if (ir->nwall > 0)
1685     {
1686         nstr = str_nelem(wall_atomtype, MAXPTR, names);
1687         if (nstr != ir->nwall)
1688         {
1689             gmx_fatal(FARGS, "Expected %d elements for wall_atomtype, found %d",
1690                       ir->nwall, nstr);
1691         }
1692         for (i = 0; i < ir->nwall; i++)
1693         {
1694             opts->wall_atomtype[i] = gmx_strdup(names[i]);
1695         }
1696
1697         if (ir->wall_type == ewt93 || ir->wall_type == ewt104)
1698         {
1699             nstr = str_nelem(wall_density, MAXPTR, names);
1700             if (nstr != ir->nwall)
1701             {
1702                 gmx_fatal(FARGS, "Expected %d elements for wall-density, found %d", ir->nwall, nstr);
1703             }
1704             for (i = 0; i < ir->nwall; i++)
1705             {
1706                 sscanf(names[i], "%lf", &dbl);
1707                 if (dbl <= 0)
1708                 {
1709                     gmx_fatal(FARGS, "wall-density[%d] = %f\n", i, dbl);
1710                 }
1711                 ir->wall_density[i] = dbl;
1712             }
1713         }
1714     }
1715 }
1716
1717 static void add_wall_energrps(gmx_groups_t *groups, int nwall, t_symtab *symtab)
1718 {
1719     int     i;
1720     t_grps *grps;
1721     char    str[STRLEN];
1722
1723     if (nwall > 0)
1724     {
1725         srenew(groups->grpname, groups->ngrpname+nwall);
1726         grps = &(groups->grps[egcENER]);
1727         srenew(grps->nm_ind, grps->nr+nwall);
1728         for (i = 0; i < nwall; i++)
1729         {
1730             sprintf(str, "wall%d", i);
1731             groups->grpname[groups->ngrpname] = put_symtab(symtab, str);
1732             grps->nm_ind[grps->nr++]          = groups->ngrpname++;
1733         }
1734     }
1735 }
1736
1737 void read_expandedparams(int *ninp_p, t_inpfile **inp_p,
1738                          t_expanded *expand, warninp_t wi)
1739 {
1740     int        ninp, nerror = 0;
1741     t_inpfile *inp;
1742
1743     ninp   = *ninp_p;
1744     inp    = *inp_p;
1745
1746     /* read expanded ensemble parameters */
1747     CCTYPE ("expanded ensemble variables");
1748     ITYPE ("nstexpanded", expand->nstexpanded, -1);
1749     EETYPE("lmc-stats", expand->elamstats, elamstats_names);
1750     EETYPE("lmc-move", expand->elmcmove, elmcmove_names);
1751     EETYPE("lmc-weights-equil", expand->elmceq, elmceq_names);
1752     ITYPE ("weight-equil-number-all-lambda", expand->equil_n_at_lam, -1);
1753     ITYPE ("weight-equil-number-samples", expand->equil_samples, -1);
1754     ITYPE ("weight-equil-number-steps", expand->equil_steps, -1);
1755     RTYPE ("weight-equil-wl-delta", expand->equil_wl_delta, -1);
1756     RTYPE ("weight-equil-count-ratio", expand->equil_ratio, -1);
1757     CCTYPE("Seed for Monte Carlo in lambda space");
1758     ITYPE ("lmc-seed", expand->lmc_seed, -1);
1759     RTYPE ("mc-temperature", expand->mc_temp, -1);
1760     ITYPE ("lmc-repeats", expand->lmc_repeats, 1);
1761     ITYPE ("lmc-gibbsdelta", expand->gibbsdeltalam, -1);
1762     ITYPE ("lmc-forced-nstart", expand->lmc_forced_nstart, 0);
1763     EETYPE("symmetrized-transition-matrix", expand->bSymmetrizedTMatrix, yesno_names);
1764     ITYPE("nst-transition-matrix", expand->nstTij, -1);
1765     ITYPE ("mininum-var-min", expand->minvarmin, 100); /*default is reasonable */
1766     ITYPE ("weight-c-range", expand->c_range, 0);      /* default is just C=0 */
1767     RTYPE ("wl-scale", expand->wl_scale, 0.8);
1768     RTYPE ("wl-ratio", expand->wl_ratio, 0.8);
1769     RTYPE ("init-wl-delta", expand->init_wl_delta, 1.0);
1770     EETYPE("wl-oneovert", expand->bWLoneovert, yesno_names);
1771
1772     *ninp_p   = ninp;
1773     *inp_p    = inp;
1774
1775     return;
1776 }
1777
1778 /*! \brief Return whether an end state with the given coupling-lambda
1779  * value describes fully-interacting VDW.
1780  *
1781  * \param[in]  couple_lambda_value  Enumeration ecouplam value describing the end state
1782  * \return                          Whether VDW is on (i.e. the user chose vdw or vdw-q in the .mdp file)
1783  */
1784 static gmx_bool couple_lambda_has_vdw_on(int couple_lambda_value)
1785 {
1786     return (couple_lambda_value == ecouplamVDW ||
1787             couple_lambda_value == ecouplamVDWQ);
1788 }
1789
1790 void get_ir(const char *mdparin, const char *mdparout,
1791             t_inputrec *ir, t_gromppopts *opts,
1792             warninp_t wi)
1793 {
1794     char       *dumstr[2];
1795     double      dumdub[2][6];
1796     t_inpfile  *inp;
1797     const char *tmp;
1798     int         i, j, m, ninp;
1799     char        warn_buf[STRLEN];
1800     t_lambda   *fep    = ir->fepvals;
1801     t_expanded *expand = ir->expandedvals;
1802
1803     init_inputrec_strings();
1804     inp = read_inpfile(mdparin, &ninp, wi);
1805
1806     snew(dumstr[0], STRLEN);
1807     snew(dumstr[1], STRLEN);
1808
1809     if (-1 == search_einp(ninp, inp, "cutoff-scheme"))
1810     {
1811         sprintf(warn_buf,
1812                 "%s did not specify a value for the .mdp option "
1813                 "\"cutoff-scheme\". Probably it was first intended for use "
1814                 "with GROMACS before 4.6. In 4.6, the Verlet scheme was "
1815                 "introduced, but the group scheme was still the default. "
1816                 "The default is now the Verlet scheme, so you will observe "
1817                 "different behaviour.", mdparin);
1818         warning_note(wi, warn_buf);
1819     }
1820
1821     /* ignore the following deprecated commands */
1822     REM_TYPE("title");
1823     REM_TYPE("cpp");
1824     REM_TYPE("domain-decomposition");
1825     REM_TYPE("andersen-seed");
1826     REM_TYPE("dihre");
1827     REM_TYPE("dihre-fc");
1828     REM_TYPE("dihre-tau");
1829     REM_TYPE("nstdihreout");
1830     REM_TYPE("nstcheckpoint");
1831     REM_TYPE("optimize-fft");
1832
1833     /* replace the following commands with the clearer new versions*/
1834     REPL_TYPE("unconstrained-start", "continuation");
1835     REPL_TYPE("foreign-lambda", "fep-lambdas");
1836     REPL_TYPE("verlet-buffer-drift", "verlet-buffer-tolerance");
1837     REPL_TYPE("nstxtcout", "nstxout-compressed");
1838     REPL_TYPE("xtc-grps", "compressed-x-grps");
1839     REPL_TYPE("xtc-precision", "compressed-x-precision");
1840
1841     CCTYPE ("VARIOUS PREPROCESSING OPTIONS");
1842     CTYPE ("Preprocessor information: use cpp syntax.");
1843     CTYPE ("e.g.: -I/home/joe/doe -I/home/mary/roe");
1844     STYPE ("include", opts->include,  NULL);
1845     CTYPE ("e.g.: -DPOSRES -DFLEXIBLE (note these variable names are case sensitive)");
1846     STYPE ("define",  opts->define,   NULL);
1847
1848     CCTYPE ("RUN CONTROL PARAMETERS");
1849     EETYPE("integrator",  ir->eI,         ei_names);
1850     CTYPE ("Start time and timestep in ps");
1851     RTYPE ("tinit",   ir->init_t, 0.0);
1852     RTYPE ("dt",      ir->delta_t,    0.001);
1853     STEPTYPE ("nsteps",   ir->nsteps,     0);
1854     CTYPE ("For exact run continuation or redoing part of a run");
1855     STEPTYPE ("init-step", ir->init_step,  0);
1856     CTYPE ("Part index is updated automatically on checkpointing (keeps files separate)");
1857     ITYPE ("simulation-part", ir->simulation_part, 1);
1858     CTYPE ("mode for center of mass motion removal");
1859     EETYPE("comm-mode",   ir->comm_mode,  ecm_names);
1860     CTYPE ("number of steps for center of mass motion removal");
1861     ITYPE ("nstcomm", ir->nstcomm,    100);
1862     CTYPE ("group(s) for center of mass motion removal");
1863     STYPE ("comm-grps",   is->vcm,            NULL);
1864
1865     CCTYPE ("LANGEVIN DYNAMICS OPTIONS");
1866     CTYPE ("Friction coefficient (amu/ps) and random seed");
1867     RTYPE ("bd-fric",     ir->bd_fric,    0.0);
1868     STEPTYPE ("ld-seed",  ir->ld_seed,    -1);
1869
1870     /* Em stuff */
1871     CCTYPE ("ENERGY MINIMIZATION OPTIONS");
1872     CTYPE ("Force tolerance and initial step-size");
1873     RTYPE ("emtol",       ir->em_tol,     10.0);
1874     RTYPE ("emstep",      ir->em_stepsize, 0.01);
1875     CTYPE ("Max number of iterations in relax-shells");
1876     ITYPE ("niter",       ir->niter,      20);
1877     CTYPE ("Step size (ps^2) for minimization of flexible constraints");
1878     RTYPE ("fcstep",      ir->fc_stepsize, 0);
1879     CTYPE ("Frequency of steepest descents steps when doing CG");
1880     ITYPE ("nstcgsteep",  ir->nstcgsteep, 1000);
1881     ITYPE ("nbfgscorr",   ir->nbfgscorr,  10);
1882
1883     CCTYPE ("TEST PARTICLE INSERTION OPTIONS");
1884     RTYPE ("rtpi",    ir->rtpi,   0.05);
1885
1886     /* Output options */
1887     CCTYPE ("OUTPUT CONTROL OPTIONS");
1888     CTYPE ("Output frequency for coords (x), velocities (v) and forces (f)");
1889     ITYPE ("nstxout", ir->nstxout,    0);
1890     ITYPE ("nstvout", ir->nstvout,    0);
1891     ITYPE ("nstfout", ir->nstfout,    0);
1892     CTYPE ("Output frequency for energies to log file and energy file");
1893     ITYPE ("nstlog",  ir->nstlog, 1000);
1894     ITYPE ("nstcalcenergy", ir->nstcalcenergy, 100);
1895     ITYPE ("nstenergy",   ir->nstenergy,  1000);
1896     CTYPE ("Output frequency and precision for .xtc file");
1897     ITYPE ("nstxout-compressed", ir->nstxout_compressed,  0);
1898     RTYPE ("compressed-x-precision", ir->x_compression_precision, 1000.0);
1899     CTYPE ("This selects the subset of atoms for the compressed");
1900     CTYPE ("trajectory file. You can select multiple groups. By");
1901     CTYPE ("default, all atoms will be written.");
1902     STYPE ("compressed-x-grps", is->x_compressed_groups, NULL);
1903     CTYPE ("Selection of energy groups");
1904     STYPE ("energygrps",  is->energy,         NULL);
1905
1906     /* Neighbor searching */
1907     CCTYPE ("NEIGHBORSEARCHING PARAMETERS");
1908     CTYPE ("cut-off scheme (Verlet: particle based cut-offs, group: using charge groups)");
1909     EETYPE("cutoff-scheme",     ir->cutoff_scheme,    ecutscheme_names);
1910     CTYPE ("nblist update frequency");
1911     ITYPE ("nstlist", ir->nstlist,    10);
1912     CTYPE ("ns algorithm (simple or grid)");
1913     EETYPE("ns-type",     ir->ns_type,    ens_names);
1914     CTYPE ("Periodic boundary conditions: xyz, no, xy");
1915     EETYPE("pbc",         ir->ePBC,       epbc_names);
1916     EETYPE("periodic-molecules", ir->bPeriodicMols, yesno_names);
1917     CTYPE ("Allowed energy error due to the Verlet buffer in kJ/mol/ps per atom,");
1918     CTYPE ("a value of -1 means: use rlist");
1919     RTYPE("verlet-buffer-tolerance", ir->verletbuf_tol,    0.005);
1920     CTYPE ("nblist cut-off");
1921     RTYPE ("rlist",   ir->rlist,  1.0);
1922     CTYPE ("long-range cut-off for switched potentials");
1923     RTYPE ("rlistlong",   ir->rlistlong,  -1);
1924     ITYPE ("nstcalclr",   ir->nstcalclr,  -1);
1925
1926     /* Electrostatics */
1927     CCTYPE ("OPTIONS FOR ELECTROSTATICS AND VDW");
1928     CTYPE ("Method for doing electrostatics");
1929     EETYPE("coulombtype", ir->coulombtype,    eel_names);
1930     EETYPE("coulomb-modifier",    ir->coulomb_modifier,    eintmod_names);
1931     CTYPE ("cut-off lengths");
1932     RTYPE ("rcoulomb-switch", ir->rcoulomb_switch,    0.0);
1933     RTYPE ("rcoulomb",    ir->rcoulomb,   1.0);
1934     CTYPE ("Relative dielectric constant for the medium and the reaction field");
1935     RTYPE ("epsilon-r",   ir->epsilon_r,  1.0);
1936     RTYPE ("epsilon-rf",  ir->epsilon_rf, 0.0);
1937     CTYPE ("Method for doing Van der Waals");
1938     EETYPE("vdw-type",    ir->vdwtype,    evdw_names);
1939     EETYPE("vdw-modifier",    ir->vdw_modifier,    eintmod_names);
1940     CTYPE ("cut-off lengths");
1941     RTYPE ("rvdw-switch", ir->rvdw_switch,    0.0);
1942     RTYPE ("rvdw",    ir->rvdw,   1.0);
1943     CTYPE ("Apply long range dispersion corrections for Energy and Pressure");
1944     EETYPE("DispCorr",    ir->eDispCorr,  edispc_names);
1945     CTYPE ("Extension of the potential lookup tables beyond the cut-off");
1946     RTYPE ("table-extension", ir->tabext, 1.0);
1947     CTYPE ("Separate tables between energy group pairs");
1948     STYPE ("energygrp-table", is->egptable,   NULL);
1949     CTYPE ("Spacing for the PME/PPPM FFT grid");
1950     RTYPE ("fourierspacing", ir->fourier_spacing, 0.12);
1951     CTYPE ("FFT grid size, when a value is 0 fourierspacing will be used");
1952     ITYPE ("fourier-nx",  ir->nkx,         0);
1953     ITYPE ("fourier-ny",  ir->nky,         0);
1954     ITYPE ("fourier-nz",  ir->nkz,         0);
1955     CTYPE ("EWALD/PME/PPPM parameters");
1956     ITYPE ("pme-order",   ir->pme_order,   4);
1957     RTYPE ("ewald-rtol",  ir->ewald_rtol, 0.00001);
1958     RTYPE ("ewald-rtol-lj", ir->ewald_rtol_lj, 0.001);
1959     EETYPE("lj-pme-comb-rule", ir->ljpme_combination_rule, eljpme_names);
1960     EETYPE("ewald-geometry", ir->ewald_geometry, eewg_names);
1961     RTYPE ("epsilon-surface", ir->epsilon_surface, 0.0);
1962
1963     CCTYPE("IMPLICIT SOLVENT ALGORITHM");
1964     EETYPE("implicit-solvent", ir->implicit_solvent, eis_names);
1965
1966     CCTYPE ("GENERALIZED BORN ELECTROSTATICS");
1967     CTYPE ("Algorithm for calculating Born radii");
1968     EETYPE("gb-algorithm", ir->gb_algorithm, egb_names);
1969     CTYPE ("Frequency of calculating the Born radii inside rlist");
1970     ITYPE ("nstgbradii", ir->nstgbradii, 1);
1971     CTYPE ("Cutoff for Born radii calculation; the contribution from atoms");
1972     CTYPE ("between rlist and rgbradii is updated every nstlist steps");
1973     RTYPE ("rgbradii",  ir->rgbradii, 1.0);
1974     CTYPE ("Dielectric coefficient of the implicit solvent");
1975     RTYPE ("gb-epsilon-solvent", ir->gb_epsilon_solvent, 80.0);
1976     CTYPE ("Salt concentration in M for Generalized Born models");
1977     RTYPE ("gb-saltconc",  ir->gb_saltconc, 0.0);
1978     CTYPE ("Scaling factors used in the OBC GB model. Default values are OBC(II)");
1979     RTYPE ("gb-obc-alpha", ir->gb_obc_alpha, 1.0);
1980     RTYPE ("gb-obc-beta", ir->gb_obc_beta, 0.8);
1981     RTYPE ("gb-obc-gamma", ir->gb_obc_gamma, 4.85);
1982     RTYPE ("gb-dielectric-offset", ir->gb_dielectric_offset, 0.009);
1983     EETYPE("sa-algorithm", ir->sa_algorithm, esa_names);
1984     CTYPE ("Surface tension (kJ/mol/nm^2) for the SA (nonpolar surface) part of GBSA");
1985     CTYPE ("The value -1 will set default value for Still/HCT/OBC GB-models.");
1986     RTYPE ("sa-surface-tension", ir->sa_surface_tension, -1);
1987
1988     /* Coupling stuff */
1989     CCTYPE ("OPTIONS FOR WEAK COUPLING ALGORITHMS");
1990     CTYPE ("Temperature coupling");
1991     EETYPE("tcoupl",  ir->etc,        etcoupl_names);
1992     ITYPE ("nsttcouple", ir->nsttcouple,  -1);
1993     ITYPE("nh-chain-length",     ir->opts.nhchainlength, 10);
1994     EETYPE("print-nose-hoover-chain-variables", ir->bPrintNHChains, yesno_names);
1995     CTYPE ("Groups to couple separately");
1996     STYPE ("tc-grps",     is->tcgrps,         NULL);
1997     CTYPE ("Time constant (ps) and reference temperature (K)");
1998     STYPE ("tau-t",   is->tau_t,      NULL);
1999     STYPE ("ref-t",   is->ref_t,      NULL);
2000     CTYPE ("pressure coupling");
2001     EETYPE("pcoupl",  ir->epc,        epcoupl_names);
2002     EETYPE("pcoupltype",  ir->epct,       epcoupltype_names);
2003     ITYPE ("nstpcouple", ir->nstpcouple,  -1);
2004     CTYPE ("Time constant (ps), compressibility (1/bar) and reference P (bar)");
2005     RTYPE ("tau-p",   ir->tau_p,  1.0);
2006     STYPE ("compressibility", dumstr[0],  NULL);
2007     STYPE ("ref-p",       dumstr[1],      NULL);
2008     CTYPE ("Scaling of reference coordinates, No, All or COM");
2009     EETYPE ("refcoord-scaling", ir->refcoord_scaling, erefscaling_names);
2010
2011     /* QMMM */
2012     CCTYPE ("OPTIONS FOR QMMM calculations");
2013     EETYPE("QMMM", ir->bQMMM, yesno_names);
2014     CTYPE ("Groups treated Quantum Mechanically");
2015     STYPE ("QMMM-grps",  is->QMMM,          NULL);
2016     CTYPE ("QM method");
2017     STYPE("QMmethod",     is->QMmethod, NULL);
2018     CTYPE ("QMMM scheme");
2019     EETYPE("QMMMscheme",  ir->QMMMscheme,    eQMMMscheme_names);
2020     CTYPE ("QM basisset");
2021     STYPE("QMbasis",      is->QMbasis, NULL);
2022     CTYPE ("QM charge");
2023     STYPE ("QMcharge",    is->QMcharge, NULL);
2024     CTYPE ("QM multiplicity");
2025     STYPE ("QMmult",      is->QMmult, NULL);
2026     CTYPE ("Surface Hopping");
2027     STYPE ("SH",          is->bSH, NULL);
2028     CTYPE ("CAS space options");
2029     STYPE ("CASorbitals",      is->CASorbitals,   NULL);
2030     STYPE ("CASelectrons",     is->CASelectrons,  NULL);
2031     STYPE ("SAon", is->SAon, NULL);
2032     STYPE ("SAoff", is->SAoff, NULL);
2033     STYPE ("SAsteps", is->SAsteps, NULL);
2034     CTYPE ("Scale factor for MM charges");
2035     RTYPE ("MMChargeScaleFactor", ir->scalefactor, 1.0);
2036     CTYPE ("Optimization of QM subsystem");
2037     STYPE ("bOPT",          is->bOPT, NULL);
2038     STYPE ("bTS",          is->bTS, NULL);
2039
2040     /* Simulated annealing */
2041     CCTYPE("SIMULATED ANNEALING");
2042     CTYPE ("Type of annealing for each temperature group (no/single/periodic)");
2043     STYPE ("annealing",   is->anneal,      NULL);
2044     CTYPE ("Number of time points to use for specifying annealing in each group");
2045     STYPE ("annealing-npoints", is->anneal_npoints, NULL);
2046     CTYPE ("List of times at the annealing points for each group");
2047     STYPE ("annealing-time",       is->anneal_time,       NULL);
2048     CTYPE ("Temp. at each annealing point, for each group.");
2049     STYPE ("annealing-temp",  is->anneal_temp,  NULL);
2050
2051     /* Startup run */
2052     CCTYPE ("GENERATE VELOCITIES FOR STARTUP RUN");
2053     EETYPE("gen-vel",     opts->bGenVel,  yesno_names);
2054     RTYPE ("gen-temp",    opts->tempi,    300.0);
2055     ITYPE ("gen-seed",    opts->seed,     -1);
2056
2057     /* Shake stuff */
2058     CCTYPE ("OPTIONS FOR BONDS");
2059     EETYPE("constraints", opts->nshake,   constraints);
2060     CTYPE ("Type of constraint algorithm");
2061     EETYPE("constraint-algorithm",  ir->eConstrAlg, econstr_names);
2062     CTYPE ("Do not constrain the start configuration");
2063     EETYPE("continuation", ir->bContinuation, yesno_names);
2064     CTYPE ("Use successive overrelaxation to reduce the number of shake iterations");
2065     EETYPE("Shake-SOR", ir->bShakeSOR, yesno_names);
2066     CTYPE ("Relative tolerance of shake");
2067     RTYPE ("shake-tol", ir->shake_tol, 0.0001);
2068     CTYPE ("Highest order in the expansion of the constraint coupling matrix");
2069     ITYPE ("lincs-order", ir->nProjOrder, 4);
2070     CTYPE ("Number of iterations in the final step of LINCS. 1 is fine for");
2071     CTYPE ("normal simulations, but use 2 to conserve energy in NVE runs.");
2072     CTYPE ("For energy minimization with constraints it should be 4 to 8.");
2073     ITYPE ("lincs-iter", ir->nLincsIter, 1);
2074     CTYPE ("Lincs will write a warning to the stderr if in one step a bond");
2075     CTYPE ("rotates over more degrees than");
2076     RTYPE ("lincs-warnangle", ir->LincsWarnAngle, 30.0);
2077     CTYPE ("Convert harmonic bonds to morse potentials");
2078     EETYPE("morse",       opts->bMorse, yesno_names);
2079
2080     /* Energy group exclusions */
2081     CCTYPE ("ENERGY GROUP EXCLUSIONS");
2082     CTYPE ("Pairs of energy groups for which all non-bonded interactions are excluded");
2083     STYPE ("energygrp-excl", is->egpexcl,     NULL);
2084
2085     /* Walls */
2086     CCTYPE ("WALLS");
2087     CTYPE ("Number of walls, type, atom types, densities and box-z scale factor for Ewald");
2088     ITYPE ("nwall", ir->nwall, 0);
2089     EETYPE("wall-type",     ir->wall_type,   ewt_names);
2090     RTYPE ("wall-r-linpot", ir->wall_r_linpot, -1);
2091     STYPE ("wall-atomtype", is->wall_atomtype, NULL);
2092     STYPE ("wall-density",  is->wall_density,  NULL);
2093     RTYPE ("wall-ewald-zfac", ir->wall_ewald_zfac, 3);
2094
2095     /* COM pulling */
2096     CCTYPE("COM PULLING");
2097     EETYPE("pull",          ir->bPull, yesno_names);
2098     if (ir->bPull)
2099     {
2100         snew(ir->pull, 1);
2101         is->pull_grp = read_pullparams(&ninp, &inp, ir->pull, wi);
2102     }
2103
2104     /* Enforced rotation */
2105     CCTYPE("ENFORCED ROTATION");
2106     CTYPE("Enforced rotation: No or Yes");
2107     EETYPE("rotation",       ir->bRot, yesno_names);
2108     if (ir->bRot)
2109     {
2110         snew(ir->rot, 1);
2111         is->rot_grp = read_rotparams(&ninp, &inp, ir->rot, wi);
2112     }
2113
2114     /* Interactive MD */
2115     ir->bIMD = FALSE;
2116     CCTYPE("Group to display and/or manipulate in interactive MD session");
2117     STYPE ("IMD-group", is->imd_grp, NULL);
2118     if (is->imd_grp[0] != '\0')
2119     {
2120         snew(ir->imd, 1);
2121         ir->bIMD = TRUE;
2122     }
2123
2124     /* Refinement */
2125     CCTYPE("NMR refinement stuff");
2126     CTYPE ("Distance restraints type: No, Simple or Ensemble");
2127     EETYPE("disre",       ir->eDisre,     edisre_names);
2128     CTYPE ("Force weighting of pairs in one distance restraint: Conservative or Equal");
2129     EETYPE("disre-weighting", ir->eDisreWeighting, edisreweighting_names);
2130     CTYPE ("Use sqrt of the time averaged times the instantaneous violation");
2131     EETYPE("disre-mixed", ir->bDisreMixed, yesno_names);
2132     RTYPE ("disre-fc",    ir->dr_fc,  1000.0);
2133     RTYPE ("disre-tau",   ir->dr_tau, 0.0);
2134     CTYPE ("Output frequency for pair distances to energy file");
2135     ITYPE ("nstdisreout", ir->nstdisreout, 100);
2136     CTYPE ("Orientation restraints: No or Yes");
2137     EETYPE("orire",       opts->bOrire,   yesno_names);
2138     CTYPE ("Orientation restraints force constant and tau for time averaging");
2139     RTYPE ("orire-fc",    ir->orires_fc,  0.0);
2140     RTYPE ("orire-tau",   ir->orires_tau, 0.0);
2141     STYPE ("orire-fitgrp", is->orirefitgrp,    NULL);
2142     CTYPE ("Output frequency for trace(SD) and S to energy file");
2143     ITYPE ("nstorireout", ir->nstorireout, 100);
2144
2145     /* free energy variables */
2146     CCTYPE ("Free energy variables");
2147     EETYPE("free-energy", ir->efep, efep_names);
2148     STYPE ("couple-moltype",  is->couple_moltype,  NULL);
2149     EETYPE("couple-lambda0", opts->couple_lam0, couple_lam);
2150     EETYPE("couple-lambda1", opts->couple_lam1, couple_lam);
2151     EETYPE("couple-intramol", opts->bCoupleIntra, yesno_names);
2152
2153     RTYPE ("init-lambda", fep->init_lambda, -1); /* start with -1 so
2154                                                     we can recognize if
2155                                                     it was not entered */
2156     ITYPE ("init-lambda-state", fep->init_fep_state, -1);
2157     RTYPE ("delta-lambda", fep->delta_lambda, 0.0);
2158     ITYPE ("nstdhdl", fep->nstdhdl, 50);
2159     STYPE ("fep-lambdas", is->fep_lambda[efptFEP], NULL);
2160     STYPE ("mass-lambdas", is->fep_lambda[efptMASS], NULL);
2161     STYPE ("coul-lambdas", is->fep_lambda[efptCOUL], NULL);
2162     STYPE ("vdw-lambdas", is->fep_lambda[efptVDW], NULL);
2163     STYPE ("bonded-lambdas", is->fep_lambda[efptBONDED], NULL);
2164     STYPE ("restraint-lambdas", is->fep_lambda[efptRESTRAINT], NULL);
2165     STYPE ("temperature-lambdas", is->fep_lambda[efptTEMPERATURE], NULL);
2166     ITYPE ("calc-lambda-neighbors", fep->lambda_neighbors, 1);
2167     STYPE ("init-lambda-weights", is->lambda_weights, NULL);
2168     EETYPE("dhdl-print-energy", fep->edHdLPrintEnergy, edHdLPrintEnergy_names);
2169     RTYPE ("sc-alpha", fep->sc_alpha, 0.0);
2170     ITYPE ("sc-power", fep->sc_power, 1);
2171     RTYPE ("sc-r-power", fep->sc_r_power, 6.0);
2172     RTYPE ("sc-sigma", fep->sc_sigma, 0.3);
2173     EETYPE("sc-coul", fep->bScCoul, yesno_names);
2174     ITYPE ("dh_hist_size", fep->dh_hist_size, 0);
2175     RTYPE ("dh_hist_spacing", fep->dh_hist_spacing, 0.1);
2176     EETYPE("separate-dhdl-file", fep->separate_dhdl_file,
2177            separate_dhdl_file_names);
2178     EETYPE("dhdl-derivatives", fep->dhdl_derivatives, dhdl_derivatives_names);
2179     ITYPE ("dh_hist_size", fep->dh_hist_size, 0);
2180     RTYPE ("dh_hist_spacing", fep->dh_hist_spacing, 0.1);
2181
2182     /* Non-equilibrium MD stuff */
2183     CCTYPE("Non-equilibrium MD stuff");
2184     STYPE ("acc-grps",    is->accgrps,        NULL);
2185     STYPE ("accelerate",  is->acc,            NULL);
2186     STYPE ("freezegrps",  is->freeze,         NULL);
2187     STYPE ("freezedim",   is->frdim,          NULL);
2188     RTYPE ("cos-acceleration", ir->cos_accel, 0);
2189     STYPE ("deform",      is->deform,         NULL);
2190
2191     /* simulated tempering variables */
2192     CCTYPE("simulated tempering variables");
2193     EETYPE("simulated-tempering", ir->bSimTemp, yesno_names);
2194     EETYPE("simulated-tempering-scaling", ir->simtempvals->eSimTempScale, esimtemp_names);
2195     RTYPE("sim-temp-low", ir->simtempvals->simtemp_low, 300.0);
2196     RTYPE("sim-temp-high", ir->simtempvals->simtemp_high, 300.0);
2197
2198     /* expanded ensemble variables */
2199     if (ir->efep == efepEXPANDED || ir->bSimTemp)
2200     {
2201         read_expandedparams(&ninp, &inp, expand, wi);
2202     }
2203
2204     /* Electric fields */
2205     CCTYPE("Electric fields");
2206     CTYPE ("Format is number of terms (int) and for all terms an amplitude (real)");
2207     CTYPE ("and a phase angle (real)");
2208     STYPE ("E-x",     is->efield_x,   NULL);
2209     CTYPE ("Time dependent (pulsed) electric field. Format is omega, time for pulse");
2210     CTYPE ("peak, and sigma (width) for pulse. Sigma = 0 removes pulse, leaving");
2211     CTYPE ("the field to be a cosine function.");
2212     STYPE ("E-xt",    is->efield_xt,  NULL);
2213     STYPE ("E-y",     is->efield_y,   NULL);
2214     STYPE ("E-yt",    is->efield_yt,  NULL);
2215     STYPE ("E-z",     is->efield_z,   NULL);
2216     STYPE ("E-zt",    is->efield_zt,  NULL);
2217
2218     CCTYPE("Ion/water position swapping for computational electrophysiology setups");
2219     CTYPE("Swap positions along direction: no, X, Y, Z");
2220     EETYPE("swapcoords", ir->eSwapCoords, eSwapTypes_names);
2221     if (ir->eSwapCoords != eswapNO)
2222     {
2223         snew(ir->swap, 1);
2224         CTYPE("Swap attempt frequency");
2225         ITYPE("swap-frequency", ir->swap->nstswap, 1);
2226         CTYPE("Two index groups that contain the compartment-partitioning atoms");
2227         STYPE("split-group0", splitgrp0, NULL);
2228         STYPE("split-group1", splitgrp1, NULL);
2229         CTYPE("Use center of mass of split groups (yes/no), otherwise center of geometry is used");
2230         EETYPE("massw-split0", ir->swap->massw_split[0], yesno_names);
2231         EETYPE("massw-split1", ir->swap->massw_split[1], yesno_names);
2232
2233         CTYPE("Group name of ions that can be exchanged with solvent molecules");
2234         STYPE("swap-group", swapgrp, NULL);
2235         CTYPE("Group name of solvent molecules");
2236         STYPE("solvent-group", solgrp, NULL);
2237
2238         CTYPE("Split cylinder: radius, upper and lower extension (nm) (this will define the channels)");
2239         CTYPE("Note that the split cylinder settings do not have an influence on the swapping protocol,");
2240         CTYPE("however, if correctly defined, the ion permeation events are counted per channel");
2241         RTYPE("cyl0-r", ir->swap->cyl0r, 2.0);
2242         RTYPE("cyl0-up", ir->swap->cyl0u, 1.0);
2243         RTYPE("cyl0-down", ir->swap->cyl0l, 1.0);
2244         RTYPE("cyl1-r", ir->swap->cyl1r, 2.0);
2245         RTYPE("cyl1-up", ir->swap->cyl1u, 1.0);
2246         RTYPE("cyl1-down", ir->swap->cyl1l, 1.0);
2247
2248         CTYPE("Average the number of ions per compartment over these many swap attempt steps");
2249         ITYPE("coupl-steps", ir->swap->nAverage, 10);
2250         CTYPE("Requested number of anions and cations for each of the two compartments");
2251         CTYPE("-1 means fix the numbers as found in time step 0");
2252         ITYPE("anionsA", ir->swap->nanions[0], -1);
2253         ITYPE("cationsA", ir->swap->ncations[0], -1);
2254         ITYPE("anionsB", ir->swap->nanions[1], -1);
2255         ITYPE("cationsB", ir->swap->ncations[1], -1);
2256         CTYPE("Start to swap ions if threshold difference to requested count is reached");
2257         RTYPE("threshold", ir->swap->threshold, 1.0);
2258     }
2259
2260     /* AdResS defined thingies */
2261     CCTYPE ("AdResS parameters");
2262     EETYPE("adress",       ir->bAdress, yesno_names);
2263     if (ir->bAdress)
2264     {
2265         snew(ir->adress, 1);
2266         read_adressparams(&ninp, &inp, ir->adress, wi);
2267     }
2268
2269     /* User defined thingies */
2270     CCTYPE ("User defined thingies");
2271     STYPE ("user1-grps",  is->user1,          NULL);
2272     STYPE ("user2-grps",  is->user2,          NULL);
2273     ITYPE ("userint1",    ir->userint1,   0);
2274     ITYPE ("userint2",    ir->userint2,   0);
2275     ITYPE ("userint3",    ir->userint3,   0);
2276     ITYPE ("userint4",    ir->userint4,   0);
2277     RTYPE ("userreal1",   ir->userreal1,  0);
2278     RTYPE ("userreal2",   ir->userreal2,  0);
2279     RTYPE ("userreal3",   ir->userreal3,  0);
2280     RTYPE ("userreal4",   ir->userreal4,  0);
2281 #undef CTYPE
2282
2283     write_inpfile(mdparout, ninp, inp, FALSE, wi);
2284     for (i = 0; (i < ninp); i++)
2285     {
2286         sfree(inp[i].name);
2287         sfree(inp[i].value);
2288     }
2289     sfree(inp);
2290
2291     /* Process options if necessary */
2292     for (m = 0; m < 2; m++)
2293     {
2294         for (i = 0; i < 2*DIM; i++)
2295         {
2296             dumdub[m][i] = 0.0;
2297         }
2298         if (ir->epc)
2299         {
2300             switch (ir->epct)
2301             {
2302                 case epctISOTROPIC:
2303                     if (sscanf(dumstr[m], "%lf", &(dumdub[m][XX])) != 1)
2304                     {
2305                         warning_error(wi, "Pressure coupling not enough values (I need 1)");
2306                     }
2307                     dumdub[m][YY] = dumdub[m][ZZ] = dumdub[m][XX];
2308                     break;
2309                 case epctSEMIISOTROPIC:
2310                 case epctSURFACETENSION:
2311                     if (sscanf(dumstr[m], "%lf%lf",
2312                                &(dumdub[m][XX]), &(dumdub[m][ZZ])) != 2)
2313                     {
2314                         warning_error(wi, "Pressure coupling not enough values (I need 2)");
2315                     }
2316                     dumdub[m][YY] = dumdub[m][XX];
2317                     break;
2318                 case epctANISOTROPIC:
2319                     if (sscanf(dumstr[m], "%lf%lf%lf%lf%lf%lf",
2320                                &(dumdub[m][XX]), &(dumdub[m][YY]), &(dumdub[m][ZZ]),
2321                                &(dumdub[m][3]), &(dumdub[m][4]), &(dumdub[m][5])) != 6)
2322                     {
2323                         warning_error(wi, "Pressure coupling not enough values (I need 6)");
2324                     }
2325                     break;
2326                 default:
2327                     gmx_fatal(FARGS, "Pressure coupling type %s not implemented yet",
2328                               epcoupltype_names[ir->epct]);
2329             }
2330         }
2331     }
2332     clear_mat(ir->ref_p);
2333     clear_mat(ir->compress);
2334     for (i = 0; i < DIM; i++)
2335     {
2336         ir->ref_p[i][i]    = dumdub[1][i];
2337         ir->compress[i][i] = dumdub[0][i];
2338     }
2339     if (ir->epct == epctANISOTROPIC)
2340     {
2341         ir->ref_p[XX][YY] = dumdub[1][3];
2342         ir->ref_p[XX][ZZ] = dumdub[1][4];
2343         ir->ref_p[YY][ZZ] = dumdub[1][5];
2344         if (ir->ref_p[XX][YY] != 0 && ir->ref_p[XX][ZZ] != 0 && ir->ref_p[YY][ZZ] != 0)
2345         {
2346             warning(wi, "All off-diagonal reference pressures are non-zero. Are you sure you want to apply a threefold shear stress?\n");
2347         }
2348         ir->compress[XX][YY] = dumdub[0][3];
2349         ir->compress[XX][ZZ] = dumdub[0][4];
2350         ir->compress[YY][ZZ] = dumdub[0][5];
2351         for (i = 0; i < DIM; i++)
2352         {
2353             for (m = 0; m < i; m++)
2354             {
2355                 ir->ref_p[i][m]    = ir->ref_p[m][i];
2356                 ir->compress[i][m] = ir->compress[m][i];
2357             }
2358         }
2359     }
2360
2361     if (ir->comm_mode == ecmNO)
2362     {
2363         ir->nstcomm = 0;
2364     }
2365
2366     opts->couple_moltype = NULL;
2367     if (strlen(is->couple_moltype) > 0)
2368     {
2369         if (ir->efep != efepNO)
2370         {
2371             opts->couple_moltype = gmx_strdup(is->couple_moltype);
2372             if (opts->couple_lam0 == opts->couple_lam1)
2373             {
2374                 warning(wi, "The lambda=0 and lambda=1 states for coupling are identical");
2375             }
2376             if (ir->eI == eiMD && (opts->couple_lam0 == ecouplamNONE ||
2377                                    opts->couple_lam1 == ecouplamNONE))
2378             {
2379                 warning(wi, "For proper sampling of the (nearly) decoupled state, stochastic dynamics should be used");
2380             }
2381         }
2382         else
2383         {
2384             warning_note(wi, "Free energy is turned off, so we will not decouple the molecule listed in your input.");
2385         }
2386     }
2387     /* FREE ENERGY AND EXPANDED ENSEMBLE OPTIONS */
2388     if (ir->efep != efepNO)
2389     {
2390         if (fep->delta_lambda > 0)
2391         {
2392             ir->efep = efepSLOWGROWTH;
2393         }
2394     }
2395
2396     if (fep->edHdLPrintEnergy == edHdLPrintEnergyYES)
2397     {
2398         fep->edHdLPrintEnergy = edHdLPrintEnergyTOTAL;
2399         warning_note(wi, "Old option for dhdl-print-energy given: "
2400                      "changing \"yes\" to \"total\"\n");
2401     }
2402
2403     if (ir->bSimTemp && (fep->edHdLPrintEnergy == edHdLPrintEnergyNO))
2404     {
2405         /* always print out the energy to dhdl if we are doing
2406            expanded ensemble, since we need the total energy for
2407            analysis if the temperature is changing. In some
2408            conditions one may only want the potential energy, so
2409            we will allow that if the appropriate mdp setting has
2410            been enabled. Otherwise, total it is:
2411          */
2412         fep->edHdLPrintEnergy = edHdLPrintEnergyTOTAL;
2413     }
2414
2415     if ((ir->efep != efepNO) || ir->bSimTemp)
2416     {
2417         ir->bExpanded = FALSE;
2418         if ((ir->efep == efepEXPANDED) || ir->bSimTemp)
2419         {
2420             ir->bExpanded = TRUE;
2421         }
2422         do_fep_params(ir, is->fep_lambda, is->lambda_weights);
2423         if (ir->bSimTemp) /* done after fep params */
2424         {
2425             do_simtemp_params(ir);
2426         }
2427
2428         /* Because sc-coul (=FALSE by default) only acts on the lambda state
2429          * setup and not on the old way of specifying the free-energy setup,
2430          * we should check for using soft-core when not needed, since that
2431          * can complicate the sampling significantly.
2432          * Note that we only check for the automated coupling setup.
2433          * If the (advanced) user does FEP through manual topology changes,
2434          * this check will not be triggered.
2435          */
2436         if (ir->efep != efepNO && ir->fepvals->n_lambda == 0 &&
2437             ir->fepvals->sc_alpha != 0 &&
2438             (couple_lambda_has_vdw_on(opts->couple_lam0) &&
2439              couple_lambda_has_vdw_on(opts->couple_lam1)))
2440         {
2441             warning(wi, "You are using soft-core interactions while the Van der Waals interactions are not decoupled (note that the sc-coul option is only active when using lambda states). Although this will not lead to errors, you will need much more sampling than without soft-core interactions. Consider using sc-alpha=0.");
2442         }
2443     }
2444     else
2445     {
2446         ir->fepvals->n_lambda = 0;
2447     }
2448
2449     /* WALL PARAMETERS */
2450
2451     do_wall_params(ir, is->wall_atomtype, is->wall_density, opts);
2452
2453     /* ORIENTATION RESTRAINT PARAMETERS */
2454
2455     if (opts->bOrire && str_nelem(is->orirefitgrp, MAXPTR, NULL) != 1)
2456     {
2457         warning_error(wi, "ERROR: Need one orientation restraint fit group\n");
2458     }
2459
2460     /* DEFORMATION PARAMETERS */
2461
2462     clear_mat(ir->deform);
2463     for (i = 0; i < 6; i++)
2464     {
2465         dumdub[0][i] = 0;
2466     }
2467     m = sscanf(is->deform, "%lf %lf %lf %lf %lf %lf",
2468                &(dumdub[0][0]), &(dumdub[0][1]), &(dumdub[0][2]),
2469                &(dumdub[0][3]), &(dumdub[0][4]), &(dumdub[0][5]));
2470     for (i = 0; i < 3; i++)
2471     {
2472         ir->deform[i][i] = dumdub[0][i];
2473     }
2474     ir->deform[YY][XX] = dumdub[0][3];
2475     ir->deform[ZZ][XX] = dumdub[0][4];
2476     ir->deform[ZZ][YY] = dumdub[0][5];
2477     if (ir->epc != epcNO)
2478     {
2479         for (i = 0; i < 3; i++)
2480         {
2481             for (j = 0; j <= i; j++)
2482             {
2483                 if (ir->deform[i][j] != 0 && ir->compress[i][j] != 0)
2484                 {
2485                     warning_error(wi, "A box element has deform set and compressibility > 0");
2486                 }
2487             }
2488         }
2489         for (i = 0; i < 3; i++)
2490         {
2491             for (j = 0; j < i; j++)
2492             {
2493                 if (ir->deform[i][j] != 0)
2494                 {
2495                     for (m = j; m < DIM; m++)
2496                     {
2497                         if (ir->compress[m][j] != 0)
2498                         {
2499                             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.");
2500                             warning(wi, warn_buf);
2501                         }
2502                     }
2503                 }
2504             }
2505         }
2506     }
2507
2508     /* Ion/water position swapping checks */
2509     if (ir->eSwapCoords != eswapNO)
2510     {
2511         if (ir->swap->nstswap < 1)
2512         {
2513             warning_error(wi, "swap_frequency must be 1 or larger when ion swapping is requested");
2514         }
2515         if (ir->swap->nAverage < 1)
2516         {
2517             warning_error(wi, "coupl_steps must be 1 or larger.\n");
2518         }
2519         if (ir->swap->threshold < 1.0)
2520         {
2521             warning_error(wi, "Ion count threshold must be at least 1.\n");
2522         }
2523     }
2524
2525     sfree(dumstr[0]);
2526     sfree(dumstr[1]);
2527 }
2528
2529 static int search_QMstring(const char *s, int ng, const char *gn[])
2530 {
2531     /* same as normal search_string, but this one searches QM strings */
2532     int i;
2533
2534     for (i = 0; (i < ng); i++)
2535     {
2536         if (gmx_strcasecmp(s, gn[i]) == 0)
2537         {
2538             return i;
2539         }
2540     }
2541
2542     gmx_fatal(FARGS, "this QM method or basisset (%s) is not implemented\n!", s);
2543
2544     return -1;
2545
2546 } /* search_QMstring */
2547
2548 /* We would like gn to be const as well, but C doesn't allow this */
2549 /* TODO this is utility functionality (search for the index of a
2550    string in a collection), so should be refactored and located more
2551    centrally. */
2552 int search_string(const char *s, int ng, char *gn[])
2553 {
2554     int i;
2555
2556     for (i = 0; (i < ng); i++)
2557     {
2558         if (gmx_strcasecmp(s, gn[i]) == 0)
2559         {
2560             return i;
2561         }
2562     }
2563
2564     gmx_fatal(FARGS,
2565               "Group %s referenced in the .mdp file was not found in the index file.\n"
2566               "Group names must match either [moleculetype] names or custom index group\n"
2567               "names, in which case you must supply an index file to the '-n' option\n"
2568               "of grompp.",
2569               s);
2570
2571     return -1;
2572 }
2573
2574 static gmx_bool do_numbering(int natoms, gmx_groups_t *groups, int ng, char *ptrs[],
2575                              t_blocka *block, char *gnames[],
2576                              int gtype, int restnm,
2577                              int grptp, gmx_bool bVerbose,
2578                              warninp_t wi)
2579 {
2580     unsigned short *cbuf;
2581     t_grps         *grps = &(groups->grps[gtype]);
2582     int             i, j, gid, aj, ognr, ntot = 0;
2583     const char     *title;
2584     gmx_bool        bRest;
2585     char            warn_buf[STRLEN];
2586
2587     if (debug)
2588     {
2589         fprintf(debug, "Starting numbering %d groups of type %d\n", ng, gtype);
2590     }
2591
2592     title = gtypes[gtype];
2593
2594     snew(cbuf, natoms);
2595     /* Mark all id's as not set */
2596     for (i = 0; (i < natoms); i++)
2597     {
2598         cbuf[i] = NOGID;
2599     }
2600
2601     snew(grps->nm_ind, ng+1); /* +1 for possible rest group */
2602     for (i = 0; (i < ng); i++)
2603     {
2604         /* Lookup the group name in the block structure */
2605         gid = search_string(ptrs[i], block->nr, gnames);
2606         if ((grptp != egrptpONE) || (i == 0))
2607         {
2608             grps->nm_ind[grps->nr++] = gid;
2609         }
2610         if (debug)
2611         {
2612             fprintf(debug, "Found gid %d for group %s\n", gid, ptrs[i]);
2613         }
2614
2615         /* Now go over the atoms in the group */
2616         for (j = block->index[gid]; (j < block->index[gid+1]); j++)
2617         {
2618
2619             aj = block->a[j];
2620
2621             /* Range checking */
2622             if ((aj < 0) || (aj >= natoms))
2623             {
2624                 gmx_fatal(FARGS, "Invalid atom number %d in indexfile", aj);
2625             }
2626             /* Lookup up the old group number */
2627             ognr = cbuf[aj];
2628             if (ognr != NOGID)
2629             {
2630                 gmx_fatal(FARGS, "Atom %d in multiple %s groups (%d and %d)",
2631                           aj+1, title, ognr+1, i+1);
2632             }
2633             else
2634             {
2635                 /* Store the group number in buffer */
2636                 if (grptp == egrptpONE)
2637                 {
2638                     cbuf[aj] = 0;
2639                 }
2640                 else
2641                 {
2642                     cbuf[aj] = i;
2643                 }
2644                 ntot++;
2645             }
2646         }
2647     }
2648
2649     /* Now check whether we have done all atoms */
2650     bRest = FALSE;
2651     if (ntot != natoms)
2652     {
2653         if (grptp == egrptpALL)
2654         {
2655             gmx_fatal(FARGS, "%d atoms are not part of any of the %s groups",
2656                       natoms-ntot, title);
2657         }
2658         else if (grptp == egrptpPART)
2659         {
2660             sprintf(warn_buf, "%d atoms are not part of any of the %s groups",
2661                     natoms-ntot, title);
2662             warning_note(wi, warn_buf);
2663         }
2664         /* Assign all atoms currently unassigned to a rest group */
2665         for (j = 0; (j < natoms); j++)
2666         {
2667             if (cbuf[j] == NOGID)
2668             {
2669                 cbuf[j] = grps->nr;
2670                 bRest   = TRUE;
2671             }
2672         }
2673         if (grptp != egrptpPART)
2674         {
2675             if (bVerbose)
2676             {
2677                 fprintf(stderr,
2678                         "Making dummy/rest group for %s containing %d elements\n",
2679                         title, natoms-ntot);
2680             }
2681             /* Add group name "rest" */
2682             grps->nm_ind[grps->nr] = restnm;
2683
2684             /* Assign the rest name to all atoms not currently assigned to a group */
2685             for (j = 0; (j < natoms); j++)
2686             {
2687                 if (cbuf[j] == NOGID)
2688                 {
2689                     cbuf[j] = grps->nr;
2690                 }
2691             }
2692             grps->nr++;
2693         }
2694     }
2695
2696     if (grps->nr == 1 && (ntot == 0 || ntot == natoms))
2697     {
2698         /* All atoms are part of one (or no) group, no index required */
2699         groups->ngrpnr[gtype] = 0;
2700         groups->grpnr[gtype]  = NULL;
2701     }
2702     else
2703     {
2704         groups->ngrpnr[gtype] = natoms;
2705         snew(groups->grpnr[gtype], natoms);
2706         for (j = 0; (j < natoms); j++)
2707         {
2708             groups->grpnr[gtype][j] = cbuf[j];
2709         }
2710     }
2711
2712     sfree(cbuf);
2713
2714     return (bRest && grptp == egrptpPART);
2715 }
2716
2717 static void calc_nrdf(gmx_mtop_t *mtop, t_inputrec *ir, char **gnames)
2718 {
2719     t_grpopts              *opts;
2720     gmx_groups_t           *groups;
2721     t_pull                 *pull;
2722     int                     natoms, ai, aj, i, j, d, g, imin, jmin;
2723     t_iatom                *ia;
2724     int                    *nrdf2, *na_vcm, na_tot;
2725     double                 *nrdf_tc, *nrdf_vcm, nrdf_uc, n_sub = 0;
2726     gmx_mtop_atomloop_all_t aloop;
2727     t_atom                 *atom;
2728     int                     mb, mol, ftype, as;
2729     gmx_molblock_t         *molb;
2730     gmx_moltype_t          *molt;
2731
2732     /* Calculate nrdf.
2733      * First calc 3xnr-atoms for each group
2734      * then subtract half a degree of freedom for each constraint
2735      *
2736      * Only atoms and nuclei contribute to the degrees of freedom...
2737      */
2738
2739     opts = &ir->opts;
2740
2741     groups = &mtop->groups;
2742     natoms = mtop->natoms;
2743
2744     /* Allocate one more for a possible rest group */
2745     /* We need to sum degrees of freedom into doubles,
2746      * since floats give too low nrdf's above 3 million atoms.
2747      */
2748     snew(nrdf_tc, groups->grps[egcTC].nr+1);
2749     snew(nrdf_vcm, groups->grps[egcVCM].nr+1);
2750     snew(na_vcm, groups->grps[egcVCM].nr+1);
2751
2752     for (i = 0; i < groups->grps[egcTC].nr; i++)
2753     {
2754         nrdf_tc[i] = 0;
2755     }
2756     for (i = 0; i < groups->grps[egcVCM].nr+1; i++)
2757     {
2758         nrdf_vcm[i] = 0;
2759     }
2760
2761     snew(nrdf2, natoms);
2762     aloop = gmx_mtop_atomloop_all_init(mtop);
2763     while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
2764     {
2765         nrdf2[i] = 0;
2766         if (atom->ptype == eptAtom || atom->ptype == eptNucleus)
2767         {
2768             g = ggrpnr(groups, egcFREEZE, i);
2769             /* Double count nrdf for particle i */
2770             for (d = 0; d < DIM; d++)
2771             {
2772                 if (opts->nFreeze[g][d] == 0)
2773                 {
2774                     nrdf2[i] += 2;
2775                 }
2776             }
2777             nrdf_tc [ggrpnr(groups, egcTC, i)]  += 0.5*nrdf2[i];
2778             nrdf_vcm[ggrpnr(groups, egcVCM, i)] += 0.5*nrdf2[i];
2779         }
2780     }
2781
2782     as = 0;
2783     for (mb = 0; mb < mtop->nmolblock; mb++)
2784     {
2785         molb = &mtop->molblock[mb];
2786         molt = &mtop->moltype[molb->type];
2787         atom = molt->atoms.atom;
2788         for (mol = 0; mol < molb->nmol; mol++)
2789         {
2790             for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
2791             {
2792                 ia = molt->ilist[ftype].iatoms;
2793                 for (i = 0; i < molt->ilist[ftype].nr; )
2794                 {
2795                     /* Subtract degrees of freedom for the constraints,
2796                      * if the particles still have degrees of freedom left.
2797                      * If one of the particles is a vsite or a shell, then all
2798                      * constraint motion will go there, but since they do not
2799                      * contribute to the constraints the degrees of freedom do not
2800                      * change.
2801                      */
2802                     ai = as + ia[1];
2803                     aj = as + ia[2];
2804                     if (((atom[ia[1]].ptype == eptNucleus) ||
2805                          (atom[ia[1]].ptype == eptAtom)) &&
2806                         ((atom[ia[2]].ptype == eptNucleus) ||
2807                          (atom[ia[2]].ptype == eptAtom)))
2808                     {
2809                         if (nrdf2[ai] > 0)
2810                         {
2811                             jmin = 1;
2812                         }
2813                         else
2814                         {
2815                             jmin = 2;
2816                         }
2817                         if (nrdf2[aj] > 0)
2818                         {
2819                             imin = 1;
2820                         }
2821                         else
2822                         {
2823                             imin = 2;
2824                         }
2825                         imin       = min(imin, nrdf2[ai]);
2826                         jmin       = min(jmin, nrdf2[aj]);
2827                         nrdf2[ai] -= imin;
2828                         nrdf2[aj] -= jmin;
2829                         nrdf_tc [ggrpnr(groups, egcTC, ai)]  -= 0.5*imin;
2830                         nrdf_tc [ggrpnr(groups, egcTC, aj)]  -= 0.5*jmin;
2831                         nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2832                         nrdf_vcm[ggrpnr(groups, egcVCM, aj)] -= 0.5*jmin;
2833                     }
2834                     ia += interaction_function[ftype].nratoms+1;
2835                     i  += interaction_function[ftype].nratoms+1;
2836                 }
2837             }
2838             ia = molt->ilist[F_SETTLE].iatoms;
2839             for (i = 0; i < molt->ilist[F_SETTLE].nr; )
2840             {
2841                 /* Subtract 1 dof from every atom in the SETTLE */
2842                 for (j = 0; j < 3; j++)
2843                 {
2844                     ai         = as + ia[1+j];
2845                     imin       = min(2, nrdf2[ai]);
2846                     nrdf2[ai] -= imin;
2847                     nrdf_tc [ggrpnr(groups, egcTC, ai)]  -= 0.5*imin;
2848                     nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2849                 }
2850                 ia += 4;
2851                 i  += 4;
2852             }
2853             as += molt->atoms.nr;
2854         }
2855     }
2856
2857     if (ir->bPull)
2858     {
2859         /* Correct nrdf for the COM constraints.
2860          * We correct using the TC and VCM group of the first atom
2861          * in the reference and pull group. If atoms in one pull group
2862          * belong to different TC or VCM groups it is anyhow difficult
2863          * to determine the optimal nrdf assignment.
2864          */
2865         pull = ir->pull;
2866
2867         for (i = 0; i < pull->ncoord; i++)
2868         {
2869             if (pull->coord[i].eType != epullCONSTRAINT)
2870             {
2871                 continue;
2872             }
2873
2874             imin = 1;
2875
2876             for (j = 0; j < 2; j++)
2877             {
2878                 const t_pull_group *pgrp;
2879
2880                 pgrp = &pull->group[pull->coord[i].group[j]];
2881
2882                 if (pgrp->nat > 0)
2883                 {
2884                     /* Subtract 1/2 dof from each group */
2885                     ai = pgrp->ind[0];
2886                     nrdf_tc [ggrpnr(groups, egcTC, ai)]  -= 0.5*imin;
2887                     nrdf_vcm[ggrpnr(groups, egcVCM, ai)] -= 0.5*imin;
2888                     if (nrdf_tc[ggrpnr(groups, egcTC, ai)] < 0)
2889                     {
2890                         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)]]);
2891                     }
2892                 }
2893                 else
2894                 {
2895                     /* We need to subtract the whole DOF from group j=1 */
2896                     imin += 1;
2897                 }
2898             }
2899         }
2900     }
2901
2902     if (ir->nstcomm != 0)
2903     {
2904         /* Subtract 3 from the number of degrees of freedom in each vcm group
2905          * when com translation is removed and 6 when rotation is removed
2906          * as well.
2907          */
2908         switch (ir->comm_mode)
2909         {
2910             case ecmLINEAR:
2911                 n_sub = ndof_com(ir);
2912                 break;
2913             case ecmANGULAR:
2914                 n_sub = 6;
2915                 break;
2916             default:
2917                 n_sub = 0;
2918                 gmx_incons("Checking comm_mode");
2919         }
2920
2921         for (i = 0; i < groups->grps[egcTC].nr; i++)
2922         {
2923             /* Count the number of atoms of TC group i for every VCM group */
2924             for (j = 0; j < groups->grps[egcVCM].nr+1; j++)
2925             {
2926                 na_vcm[j] = 0;
2927             }
2928             na_tot = 0;
2929             for (ai = 0; ai < natoms; ai++)
2930             {
2931                 if (ggrpnr(groups, egcTC, ai) == i)
2932                 {
2933                     na_vcm[ggrpnr(groups, egcVCM, ai)]++;
2934                     na_tot++;
2935                 }
2936             }
2937             /* Correct for VCM removal according to the fraction of each VCM
2938              * group present in this TC group.
2939              */
2940             nrdf_uc = nrdf_tc[i];
2941             if (debug)
2942             {
2943                 fprintf(debug, "T-group[%d] nrdf_uc = %g, n_sub = %g\n",
2944                         i, nrdf_uc, n_sub);
2945             }
2946             nrdf_tc[i] = 0;
2947             for (j = 0; j < groups->grps[egcVCM].nr+1; j++)
2948             {
2949                 if (nrdf_vcm[j] > n_sub)
2950                 {
2951                     nrdf_tc[i] += nrdf_uc*((double)na_vcm[j]/(double)na_tot)*
2952                         (nrdf_vcm[j] - n_sub)/nrdf_vcm[j];
2953                 }
2954                 if (debug)
2955                 {
2956                     fprintf(debug, "  nrdf_vcm[%d] = %g, nrdf = %g\n",
2957                             j, nrdf_vcm[j], nrdf_tc[i]);
2958                 }
2959             }
2960         }
2961     }
2962     for (i = 0; (i < groups->grps[egcTC].nr); i++)
2963     {
2964         opts->nrdf[i] = nrdf_tc[i];
2965         if (opts->nrdf[i] < 0)
2966         {
2967             opts->nrdf[i] = 0;
2968         }
2969         fprintf(stderr,
2970                 "Number of degrees of freedom in T-Coupling group %s is %.2f\n",
2971                 gnames[groups->grps[egcTC].nm_ind[i]], opts->nrdf[i]);
2972     }
2973
2974     sfree(nrdf2);
2975     sfree(nrdf_tc);
2976     sfree(nrdf_vcm);
2977     sfree(na_vcm);
2978 }
2979
2980 static void decode_cos(char *s, t_cosines *cosine)
2981 {
2982     char   *t;
2983     char    format[STRLEN], f1[STRLEN];
2984     double  a, phi;
2985     int     i;
2986
2987     t = gmx_strdup(s);
2988     trim(t);
2989
2990     cosine->n   = 0;
2991     cosine->a   = NULL;
2992     cosine->phi = NULL;
2993     if (strlen(t))
2994     {
2995         sscanf(t, "%d", &(cosine->n));
2996         if (cosine->n <= 0)
2997         {
2998             cosine->n = 0;
2999         }
3000         else
3001         {
3002             snew(cosine->a, cosine->n);
3003             snew(cosine->phi, cosine->n);
3004
3005             sprintf(format, "%%*d");
3006             for (i = 0; (i < cosine->n); i++)
3007             {
3008                 strcpy(f1, format);
3009                 strcat(f1, "%lf%lf");
3010                 if (sscanf(t, f1, &a, &phi) < 2)
3011                 {
3012                     gmx_fatal(FARGS, "Invalid input for electric field shift: '%s'", t);
3013                 }
3014                 cosine->a[i]   = a;
3015                 cosine->phi[i] = phi;
3016                 strcat(format, "%*lf%*lf");
3017             }
3018         }
3019     }
3020     sfree(t);
3021 }
3022
3023 static gmx_bool do_egp_flag(t_inputrec *ir, gmx_groups_t *groups,
3024                             const char *option, const char *val, int flag)
3025 {
3026     /* The maximum number of energy group pairs would be MAXPTR*(MAXPTR+1)/2.
3027      * But since this is much larger than STRLEN, such a line can not be parsed.
3028      * The real maximum is the number of names that fit in a string: STRLEN/2.
3029      */
3030 #define EGP_MAX (STRLEN/2)
3031     int      nelem, i, j, k, nr;
3032     char    *names[EGP_MAX];
3033     char  ***gnames;
3034     gmx_bool bSet;
3035
3036     gnames = groups->grpname;
3037
3038     nelem = str_nelem(val, EGP_MAX, names);
3039     if (nelem % 2 != 0)
3040     {
3041         gmx_fatal(FARGS, "The number of groups for %s is odd", option);
3042     }
3043     nr   = groups->grps[egcENER].nr;
3044     bSet = FALSE;
3045     for (i = 0; i < nelem/2; i++)
3046     {
3047         j = 0;
3048         while ((j < nr) &&
3049                gmx_strcasecmp(names[2*i], *(gnames[groups->grps[egcENER].nm_ind[j]])))
3050         {
3051             j++;
3052         }
3053         if (j == nr)
3054         {
3055             gmx_fatal(FARGS, "%s in %s is not an energy group\n",
3056                       names[2*i], option);
3057         }
3058         k = 0;
3059         while ((k < nr) &&
3060                gmx_strcasecmp(names[2*i+1], *(gnames[groups->grps[egcENER].nm_ind[k]])))
3061         {
3062             k++;
3063         }
3064         if (k == nr)
3065         {
3066             gmx_fatal(FARGS, "%s in %s is not an energy group\n",
3067                       names[2*i+1], option);
3068         }
3069         if ((j < nr) && (k < nr))
3070         {
3071             ir->opts.egp_flags[nr*j+k] |= flag;
3072             ir->opts.egp_flags[nr*k+j] |= flag;
3073             bSet = TRUE;
3074         }
3075     }
3076
3077     return bSet;
3078 }
3079
3080
3081 static void make_swap_groups(
3082         t_swapcoords *swap,
3083         char         *swapgname,
3084         char         *splitg0name,
3085         char         *splitg1name,
3086         char         *solgname,
3087         t_blocka     *grps,
3088         char        **gnames)
3089 {
3090     int   ig = -1, i = 0, j;
3091     char *splitg;
3092
3093
3094     /* Just a quick check here, more thorough checks are in mdrun */
3095     if (strcmp(splitg0name, splitg1name) == 0)
3096     {
3097         gmx_fatal(FARGS, "The split groups can not both be '%s'.", splitg0name);
3098     }
3099
3100     /* First get the swap group index atoms */
3101     ig        = search_string(swapgname, grps->nr, gnames);
3102     swap->nat = grps->index[ig+1] - grps->index[ig];
3103     if (swap->nat > 0)
3104     {
3105         fprintf(stderr, "Swap group '%s' contains %d atoms.\n", swapgname, swap->nat);
3106         snew(swap->ind, swap->nat);
3107         for (i = 0; i < swap->nat; i++)
3108         {
3109             swap->ind[i] = grps->a[grps->index[ig]+i];
3110         }
3111     }
3112     else
3113     {
3114         gmx_fatal(FARGS, "You defined an empty group of atoms for swapping.");
3115     }
3116
3117     /* Now do so for the split groups */
3118     for (j = 0; j < 2; j++)
3119     {
3120         if (j == 0)
3121         {
3122             splitg = splitg0name;
3123         }
3124         else
3125         {
3126             splitg = splitg1name;
3127         }
3128
3129         ig                 = search_string(splitg, grps->nr, gnames);
3130         swap->nat_split[j] = grps->index[ig+1] - grps->index[ig];
3131         if (swap->nat_split[j] > 0)
3132         {
3133             fprintf(stderr, "Split group %d '%s' contains %d atom%s.\n",
3134                     j, splitg, swap->nat_split[j], (swap->nat_split[j] > 1) ? "s" : "");
3135             snew(swap->ind_split[j], swap->nat_split[j]);
3136             for (i = 0; i < swap->nat_split[j]; i++)
3137             {
3138                 swap->ind_split[j][i] = grps->a[grps->index[ig]+i];
3139             }
3140         }
3141         else
3142         {
3143             gmx_fatal(FARGS, "Split group %d has to contain at least 1 atom!", j);
3144         }
3145     }
3146
3147     /* Now get the solvent group index atoms */
3148     ig            = search_string(solgname, grps->nr, gnames);
3149     swap->nat_sol = grps->index[ig+1] - grps->index[ig];
3150     if (swap->nat_sol > 0)
3151     {
3152         fprintf(stderr, "Solvent group '%s' contains %d atoms.\n", solgname, swap->nat_sol);
3153         snew(swap->ind_sol, swap->nat_sol);
3154         for (i = 0; i < swap->nat_sol; i++)
3155         {
3156             swap->ind_sol[i] = grps->a[grps->index[ig]+i];
3157         }
3158     }
3159     else
3160     {
3161         gmx_fatal(FARGS, "You defined an empty group of solvent. Cannot exchange ions.");
3162     }
3163 }
3164
3165
3166 void make_IMD_group(t_IMD *IMDgroup, char *IMDgname, t_blocka *grps, char **gnames)
3167 {
3168     int      ig, i;
3169
3170
3171     ig            = search_string(IMDgname, grps->nr, gnames);
3172     IMDgroup->nat = grps->index[ig+1] - grps->index[ig];
3173
3174     if (IMDgroup->nat > 0)
3175     {
3176         fprintf(stderr, "Group '%s' with %d atoms can be activated for interactive molecular dynamics (IMD).\n",
3177                 IMDgname, IMDgroup->nat);
3178         snew(IMDgroup->ind, IMDgroup->nat);
3179         for (i = 0; i < IMDgroup->nat; i++)
3180         {
3181             IMDgroup->ind[i] = grps->a[grps->index[ig]+i];
3182         }
3183     }
3184 }
3185
3186
3187 void do_index(const char* mdparin, const char *ndx,
3188               gmx_mtop_t *mtop,
3189               gmx_bool bVerbose,
3190               t_inputrec *ir, rvec *v,
3191               warninp_t wi)
3192 {
3193     t_blocka     *grps;
3194     gmx_groups_t *groups;
3195     int           natoms;
3196     t_symtab     *symtab;
3197     t_atoms       atoms_all;
3198     char          warnbuf[STRLEN], **gnames;
3199     int           nr, ntcg, ntau_t, nref_t, nacc, nofg, nSA, nSA_points, nSA_time, nSA_temp;
3200     real          tau_min;
3201     int           nstcmin;
3202     int           nacg, nfreeze, nfrdim, nenergy, nvcm, nuser;
3203     char         *ptr1[MAXPTR], *ptr2[MAXPTR], *ptr3[MAXPTR];
3204     int           i, j, k, restnm;
3205     real          SAtime;
3206     gmx_bool      bExcl, bTable, bSetTCpar, bAnneal, bRest;
3207     int           nQMmethod, nQMbasis, nQMcharge, nQMmult, nbSH, nCASorb, nCASelec,
3208                   nSAon, nSAoff, nSAsteps, nQMg, nbOPT, nbTS;
3209     char          warn_buf[STRLEN];
3210
3211     if (bVerbose)
3212     {
3213         fprintf(stderr, "processing index file...\n");
3214     }
3215     debug_gmx();
3216     if (ndx == NULL)
3217     {
3218         snew(grps, 1);
3219         snew(grps->index, 1);
3220         snew(gnames, 1);
3221         atoms_all = gmx_mtop_global_atoms(mtop);
3222         analyse(&atoms_all, grps, &gnames, FALSE, TRUE);
3223         free_t_atoms(&atoms_all, FALSE);
3224     }
3225     else
3226     {
3227         grps = init_index(ndx, &gnames);
3228     }
3229
3230     groups = &mtop->groups;
3231     natoms = mtop->natoms;
3232     symtab = &mtop->symtab;
3233
3234     snew(groups->grpname, grps->nr+1);
3235
3236     for (i = 0; (i < grps->nr); i++)
3237     {
3238         groups->grpname[i] = put_symtab(symtab, gnames[i]);
3239     }
3240     groups->grpname[i] = put_symtab(symtab, "rest");
3241     restnm             = i;
3242     srenew(gnames, grps->nr+1);
3243     gnames[restnm]   = *(groups->grpname[i]);
3244     groups->ngrpname = grps->nr+1;
3245
3246     set_warning_line(wi, mdparin, -1);
3247
3248     ntau_t = str_nelem(is->tau_t, MAXPTR, ptr1);
3249     nref_t = str_nelem(is->ref_t, MAXPTR, ptr2);
3250     ntcg   = str_nelem(is->tcgrps, MAXPTR, ptr3);
3251     if ((ntau_t != ntcg) || (nref_t != ntcg))
3252     {
3253         gmx_fatal(FARGS, "Invalid T coupling input: %d groups, %d ref-t values and "
3254                   "%d tau-t values", ntcg, nref_t, ntau_t);
3255     }
3256
3257     bSetTCpar = (ir->etc || EI_SD(ir->eI) || ir->eI == eiBD || EI_TPI(ir->eI));
3258     do_numbering(natoms, groups, ntcg, ptr3, grps, gnames, egcTC,
3259                  restnm, bSetTCpar ? egrptpALL : egrptpALL_GENREST, bVerbose, wi);
3260     nr            = groups->grps[egcTC].nr;
3261     ir->opts.ngtc = nr;
3262     snew(ir->opts.nrdf, nr);
3263     snew(ir->opts.tau_t, nr);
3264     snew(ir->opts.ref_t, nr);
3265     if (ir->eI == eiBD && ir->bd_fric == 0)
3266     {
3267         fprintf(stderr, "bd-fric=0, so tau-t will be used as the inverse friction constant(s)\n");
3268     }
3269
3270     if (bSetTCpar)
3271     {
3272         if (nr != nref_t)
3273         {
3274             gmx_fatal(FARGS, "Not enough ref-t and tau-t values!");
3275         }
3276
3277         tau_min = 1e20;
3278         for (i = 0; (i < nr); i++)
3279         {
3280             ir->opts.tau_t[i] = strtod(ptr1[i], NULL);
3281             if ((ir->eI == eiBD || ir->eI == eiSD2) && ir->opts.tau_t[i] <= 0)
3282             {
3283                 sprintf(warn_buf, "With integrator %s tau-t should be larger than 0", ei_names[ir->eI]);
3284                 warning_error(wi, warn_buf);
3285             }
3286
3287             if (ir->etc != etcVRESCALE && ir->opts.tau_t[i] == 0)
3288             {
3289                 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.");
3290             }
3291
3292             if (ir->opts.tau_t[i] >= 0)
3293             {
3294                 tau_min = min(tau_min, ir->opts.tau_t[i]);
3295             }
3296         }
3297         if (ir->etc != etcNO && ir->nsttcouple == -1)
3298         {
3299             ir->nsttcouple = ir_optimal_nsttcouple(ir);
3300         }
3301
3302         if (EI_VV(ir->eI))
3303         {
3304             if ((ir->etc == etcNOSEHOOVER) && (ir->epc == epcBERENDSEN))
3305             {
3306                 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");
3307             }
3308             if ((ir->epc == epcMTTK) && (ir->etc > etcNO))
3309             {
3310                 if (ir->nstpcouple != ir->nsttcouple)
3311                 {
3312                     int mincouple = min(ir->nstpcouple, ir->nsttcouple);
3313                     ir->nstpcouple = ir->nsttcouple = mincouple;
3314                     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);
3315                     warning_note(wi, warn_buf);
3316                 }
3317             }
3318         }
3319         /* velocity verlet with averaged kinetic energy KE = 0.5*(v(t+1/2) - v(t-1/2)) is implemented
3320            primarily for testing purposes, and does not work with temperature coupling other than 1 */
3321
3322         if (ETC_ANDERSEN(ir->etc))
3323         {
3324             if (ir->nsttcouple != 1)
3325             {
3326                 ir->nsttcouple = 1;
3327                 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");
3328                 warning_note(wi, warn_buf);
3329             }
3330         }
3331         nstcmin = tcouple_min_integration_steps(ir->etc);
3332         if (nstcmin > 1)
3333         {
3334             if (tau_min/(ir->delta_t*ir->nsttcouple) < nstcmin - 10*GMX_REAL_EPS)
3335             {
3336                 sprintf(warn_buf, "For proper integration of the %s thermostat, tau-t (%g) should be at least %d times larger than nsttcouple*dt (%g)",
3337                         ETCOUPLTYPE(ir->etc),
3338                         tau_min, nstcmin,
3339                         ir->nsttcouple*ir->delta_t);
3340                 warning(wi, warn_buf);
3341             }
3342         }
3343         for (i = 0; (i < nr); i++)
3344         {
3345             ir->opts.ref_t[i] = strtod(ptr2[i], NULL);
3346             if (ir->opts.ref_t[i] < 0)
3347             {
3348                 gmx_fatal(FARGS, "ref-t for group %d negative", i);
3349             }
3350         }
3351         /* set the lambda mc temperature to the md integrator temperature (which should be defined
3352            if we are in this conditional) if mc_temp is negative */
3353         if (ir->expandedvals->mc_temp < 0)
3354         {
3355             ir->expandedvals->mc_temp = ir->opts.ref_t[0]; /*for now, set to the first reft */
3356         }
3357     }
3358
3359     /* Simulated annealing for each group. There are nr groups */
3360     nSA = str_nelem(is->anneal, MAXPTR, ptr1);
3361     if (nSA == 1 && (ptr1[0][0] == 'n' || ptr1[0][0] == 'N'))
3362     {
3363         nSA = 0;
3364     }
3365     if (nSA > 0 && nSA != nr)
3366     {
3367         gmx_fatal(FARGS, "Not enough annealing values: %d (for %d groups)\n", nSA, nr);
3368     }
3369     else
3370     {
3371         snew(ir->opts.annealing, nr);
3372         snew(ir->opts.anneal_npoints, nr);
3373         snew(ir->opts.anneal_time, nr);
3374         snew(ir->opts.anneal_temp, nr);
3375         for (i = 0; i < nr; i++)
3376         {
3377             ir->opts.annealing[i]      = eannNO;
3378             ir->opts.anneal_npoints[i] = 0;
3379             ir->opts.anneal_time[i]    = NULL;
3380             ir->opts.anneal_temp[i]    = NULL;
3381         }
3382         if (nSA > 0)
3383         {
3384             bAnneal = FALSE;
3385             for (i = 0; i < nr; i++)
3386             {
3387                 if (ptr1[i][0] == 'n' || ptr1[i][0] == 'N')
3388                 {
3389                     ir->opts.annealing[i] = eannNO;
3390                 }
3391                 else if (ptr1[i][0] == 's' || ptr1[i][0] == 'S')
3392                 {
3393                     ir->opts.annealing[i] = eannSINGLE;
3394                     bAnneal               = TRUE;
3395                 }
3396                 else if (ptr1[i][0] == 'p' || ptr1[i][0] == 'P')
3397                 {
3398                     ir->opts.annealing[i] = eannPERIODIC;
3399                     bAnneal               = TRUE;
3400                 }
3401             }
3402             if (bAnneal)
3403             {
3404                 /* Read the other fields too */
3405                 nSA_points = str_nelem(is->anneal_npoints, MAXPTR, ptr1);
3406                 if (nSA_points != nSA)
3407                 {
3408                     gmx_fatal(FARGS, "Found %d annealing-npoints values for %d groups\n", nSA_points, nSA);
3409                 }
3410                 for (k = 0, i = 0; i < nr; i++)
3411                 {
3412                     ir->opts.anneal_npoints[i] = strtol(ptr1[i], NULL, 10);
3413                     if (ir->opts.anneal_npoints[i] == 1)
3414                     {
3415                         gmx_fatal(FARGS, "Please specify at least a start and an end point for annealing\n");
3416                     }
3417                     snew(ir->opts.anneal_time[i], ir->opts.anneal_npoints[i]);
3418                     snew(ir->opts.anneal_temp[i], ir->opts.anneal_npoints[i]);
3419                     k += ir->opts.anneal_npoints[i];
3420                 }
3421
3422                 nSA_time = str_nelem(is->anneal_time, MAXPTR, ptr1);
3423                 if (nSA_time != k)
3424                 {
3425                     gmx_fatal(FARGS, "Found %d annealing-time values, wanter %d\n", nSA_time, k);
3426                 }
3427                 nSA_temp = str_nelem(is->anneal_temp, MAXPTR, ptr2);
3428                 if (nSA_temp != k)
3429                 {
3430                     gmx_fatal(FARGS, "Found %d annealing-temp values, wanted %d\n", nSA_temp, k);
3431                 }
3432
3433                 for (i = 0, k = 0; i < nr; i++)
3434                 {
3435
3436                     for (j = 0; j < ir->opts.anneal_npoints[i]; j++)
3437                     {
3438                         ir->opts.anneal_time[i][j] = strtod(ptr1[k], NULL);
3439                         ir->opts.anneal_temp[i][j] = strtod(ptr2[k], NULL);
3440                         if (j == 0)
3441                         {
3442                             if (ir->opts.anneal_time[i][0] > (ir->init_t+GMX_REAL_EPS))
3443                             {
3444                                 gmx_fatal(FARGS, "First time point for annealing > init_t.\n");
3445                             }
3446                         }
3447                         else
3448                         {
3449                             /* j>0 */
3450                             if (ir->opts.anneal_time[i][j] < ir->opts.anneal_time[i][j-1])
3451                             {
3452                                 gmx_fatal(FARGS, "Annealing timepoints out of order: t=%f comes after t=%f\n",
3453                                           ir->opts.anneal_time[i][j], ir->opts.anneal_time[i][j-1]);
3454                             }
3455                         }
3456                         if (ir->opts.anneal_temp[i][j] < 0)
3457                         {
3458                             gmx_fatal(FARGS, "Found negative temperature in annealing: %f\n", ir->opts.anneal_temp[i][j]);
3459                         }
3460                         k++;
3461                     }
3462                 }
3463                 /* Print out some summary information, to make sure we got it right */
3464                 for (i = 0, k = 0; i < nr; i++)
3465                 {
3466                     if (ir->opts.annealing[i] != eannNO)
3467                     {
3468                         j = groups->grps[egcTC].nm_ind[i];
3469                         fprintf(stderr, "Simulated annealing for group %s: %s, %d timepoints\n",
3470                                 *(groups->grpname[j]), eann_names[ir->opts.annealing[i]],
3471                                 ir->opts.anneal_npoints[i]);
3472                         fprintf(stderr, "Time (ps)   Temperature (K)\n");
3473                         /* All terms except the last one */
3474                         for (j = 0; j < (ir->opts.anneal_npoints[i]-1); j++)
3475                         {
3476                             fprintf(stderr, "%9.1f      %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3477                         }
3478
3479                         /* Finally the last one */
3480                         j = ir->opts.anneal_npoints[i]-1;
3481                         if (ir->opts.annealing[i] == eannSINGLE)
3482                         {
3483                             fprintf(stderr, "%9.1f-     %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3484                         }
3485                         else
3486                         {
3487                             fprintf(stderr, "%9.1f      %5.1f\n", ir->opts.anneal_time[i][j], ir->opts.anneal_temp[i][j]);
3488                             if (fabs(ir->opts.anneal_temp[i][j]-ir->opts.anneal_temp[i][0]) > GMX_REAL_EPS)
3489                             {
3490                                 warning_note(wi, "There is a temperature jump when your annealing loops back.\n");
3491                             }
3492                         }
3493                     }
3494                 }
3495             }
3496         }
3497     }
3498
3499     if (ir->bPull)
3500     {
3501         make_pull_groups(ir->pull, is->pull_grp, grps, gnames);
3502
3503         make_pull_coords(ir->pull);
3504     }
3505
3506     if (ir->bRot)
3507     {
3508         make_rotation_groups(ir->rot, is->rot_grp, grps, gnames);
3509     }
3510
3511     if (ir->eSwapCoords != eswapNO)
3512     {
3513         make_swap_groups(ir->swap, swapgrp, splitgrp0, splitgrp1, solgrp, grps, gnames);
3514     }
3515
3516     /* Make indices for IMD session */
3517     if (ir->bIMD)
3518     {
3519         make_IMD_group(ir->imd, is->imd_grp, grps, gnames);
3520     }
3521
3522     nacc = str_nelem(is->acc, MAXPTR, ptr1);
3523     nacg = str_nelem(is->accgrps, MAXPTR, ptr2);
3524     if (nacg*DIM != nacc)
3525     {
3526         gmx_fatal(FARGS, "Invalid Acceleration input: %d groups and %d acc. values",
3527                   nacg, nacc);
3528     }
3529     do_numbering(natoms, groups, nacg, ptr2, grps, gnames, egcACC,
3530                  restnm, egrptpALL_GENREST, bVerbose, wi);
3531     nr = groups->grps[egcACC].nr;
3532     snew(ir->opts.acc, nr);
3533     ir->opts.ngacc = nr;
3534
3535     for (i = k = 0; (i < nacg); i++)
3536     {
3537         for (j = 0; (j < DIM); j++, k++)
3538         {
3539             ir->opts.acc[i][j] = strtod(ptr1[k], NULL);
3540         }
3541     }
3542     for (; (i < nr); i++)
3543     {
3544         for (j = 0; (j < DIM); j++)
3545         {
3546             ir->opts.acc[i][j] = 0;
3547         }
3548     }
3549
3550     nfrdim  = str_nelem(is->frdim, MAXPTR, ptr1);
3551     nfreeze = str_nelem(is->freeze, MAXPTR, ptr2);
3552     if (nfrdim != DIM*nfreeze)
3553     {
3554         gmx_fatal(FARGS, "Invalid Freezing input: %d groups and %d freeze values",
3555                   nfreeze, nfrdim);
3556     }
3557     do_numbering(natoms, groups, nfreeze, ptr2, grps, gnames, egcFREEZE,
3558                  restnm, egrptpALL_GENREST, bVerbose, wi);
3559     nr             = groups->grps[egcFREEZE].nr;
3560     ir->opts.ngfrz = nr;
3561     snew(ir->opts.nFreeze, nr);
3562     for (i = k = 0; (i < nfreeze); i++)
3563     {
3564         for (j = 0; (j < DIM); j++, k++)
3565         {
3566             ir->opts.nFreeze[i][j] = (gmx_strncasecmp(ptr1[k], "Y", 1) == 0);
3567             if (!ir->opts.nFreeze[i][j])
3568             {
3569                 if (gmx_strncasecmp(ptr1[k], "N", 1) != 0)
3570                 {
3571                     sprintf(warnbuf, "Please use Y(ES) or N(O) for freezedim only "
3572                             "(not %s)", ptr1[k]);
3573                     warning(wi, warn_buf);
3574                 }
3575             }
3576         }
3577     }
3578     for (; (i < nr); i++)
3579     {
3580         for (j = 0; (j < DIM); j++)
3581         {
3582             ir->opts.nFreeze[i][j] = 0;
3583         }
3584     }
3585
3586     nenergy = str_nelem(is->energy, MAXPTR, ptr1);
3587     do_numbering(natoms, groups, nenergy, ptr1, grps, gnames, egcENER,
3588                  restnm, egrptpALL_GENREST, bVerbose, wi);
3589     add_wall_energrps(groups, ir->nwall, symtab);
3590     ir->opts.ngener = groups->grps[egcENER].nr;
3591     nvcm            = str_nelem(is->vcm, MAXPTR, ptr1);
3592     bRest           =
3593         do_numbering(natoms, groups, nvcm, ptr1, grps, gnames, egcVCM,
3594                      restnm, nvcm == 0 ? egrptpALL_GENREST : egrptpPART, bVerbose, wi);
3595     if (bRest)
3596     {
3597         warning(wi, "Some atoms are not part of any center of mass motion removal group.\n"
3598                 "This may lead to artifacts.\n"
3599                 "In most cases one should use one group for the whole system.");
3600     }
3601
3602     /* Now we have filled the freeze struct, so we can calculate NRDF */
3603     calc_nrdf(mtop, ir, gnames);
3604
3605     if (v && NULL)
3606     {
3607         real fac, ntot = 0;
3608
3609         /* Must check per group! */
3610         for (i = 0; (i < ir->opts.ngtc); i++)
3611         {
3612             ntot += ir->opts.nrdf[i];
3613         }
3614         if (ntot != (DIM*natoms))
3615         {
3616             fac = sqrt(ntot/(DIM*natoms));
3617             if (bVerbose)
3618             {
3619                 fprintf(stderr, "Scaling velocities by a factor of %.3f to account for constraints\n"
3620                         "and removal of center of mass motion\n", fac);
3621             }
3622             for (i = 0; (i < natoms); i++)
3623             {
3624                 svmul(fac, v[i], v[i]);
3625             }
3626         }
3627     }
3628
3629     nuser = str_nelem(is->user1, MAXPTR, ptr1);
3630     do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcUser1,
3631                  restnm, egrptpALL_GENREST, bVerbose, wi);
3632     nuser = str_nelem(is->user2, MAXPTR, ptr1);
3633     do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcUser2,
3634                  restnm, egrptpALL_GENREST, bVerbose, wi);
3635     nuser = str_nelem(is->x_compressed_groups, MAXPTR, ptr1);
3636     do_numbering(natoms, groups, nuser, ptr1, grps, gnames, egcCompressedX,
3637                  restnm, egrptpONE, bVerbose, wi);
3638     nofg = str_nelem(is->orirefitgrp, MAXPTR, ptr1);
3639     do_numbering(natoms, groups, nofg, ptr1, grps, gnames, egcORFIT,
3640                  restnm, egrptpALL_GENREST, bVerbose, wi);
3641
3642     /* QMMM input processing */
3643     nQMg          = str_nelem(is->QMMM, MAXPTR, ptr1);
3644     nQMmethod     = str_nelem(is->QMmethod, MAXPTR, ptr2);
3645     nQMbasis      = str_nelem(is->QMbasis, MAXPTR, ptr3);
3646     if ((nQMmethod != nQMg) || (nQMbasis != nQMg))
3647     {
3648         gmx_fatal(FARGS, "Invalid QMMM input: %d groups %d basissets"
3649                   " and %d methods\n", nQMg, nQMbasis, nQMmethod);
3650     }
3651     /* group rest, if any, is always MM! */
3652     do_numbering(natoms, groups, nQMg, ptr1, grps, gnames, egcQMMM,
3653                  restnm, egrptpALL_GENREST, bVerbose, wi);
3654     nr            = nQMg; /*atoms->grps[egcQMMM].nr;*/
3655     ir->opts.ngQM = nQMg;
3656     snew(ir->opts.QMmethod, nr);
3657     snew(ir->opts.QMbasis, nr);
3658     for (i = 0; i < nr; i++)
3659     {
3660         /* input consists of strings: RHF CASSCF PM3 .. These need to be
3661          * converted to the corresponding enum in names.c
3662          */
3663         ir->opts.QMmethod[i] = search_QMstring(ptr2[i], eQMmethodNR,
3664                                                eQMmethod_names);
3665         ir->opts.QMbasis[i]  = search_QMstring(ptr3[i], eQMbasisNR,
3666                                                eQMbasis_names);
3667
3668     }
3669     nQMmult   = str_nelem(is->QMmult, MAXPTR, ptr1);
3670     nQMcharge = str_nelem(is->QMcharge, MAXPTR, ptr2);
3671     nbSH      = str_nelem(is->bSH, MAXPTR, ptr3);
3672     snew(ir->opts.QMmult, nr);
3673     snew(ir->opts.QMcharge, nr);
3674     snew(ir->opts.bSH, nr);
3675
3676     for (i = 0; i < nr; i++)
3677     {
3678         ir->opts.QMmult[i]   = strtol(ptr1[i], NULL, 10);
3679         ir->opts.QMcharge[i] = strtol(ptr2[i], NULL, 10);
3680         ir->opts.bSH[i]      = (gmx_strncasecmp(ptr3[i], "Y", 1) == 0);
3681     }
3682
3683     nCASelec  = str_nelem(is->CASelectrons, MAXPTR, ptr1);
3684     nCASorb   = str_nelem(is->CASorbitals, MAXPTR, ptr2);
3685     snew(ir->opts.CASelectrons, nr);
3686     snew(ir->opts.CASorbitals, nr);
3687     for (i = 0; i < nr; i++)
3688     {
3689         ir->opts.CASelectrons[i] = strtol(ptr1[i], NULL, 10);
3690         ir->opts.CASorbitals[i]  = strtol(ptr2[i], NULL, 10);
3691     }
3692     /* special optimization options */
3693
3694     nbOPT = str_nelem(is->bOPT, MAXPTR, ptr1);
3695     nbTS  = str_nelem(is->bTS, MAXPTR, ptr2);
3696     snew(ir->opts.bOPT, nr);
3697     snew(ir->opts.bTS, nr);
3698     for (i = 0; i < nr; i++)
3699     {
3700         ir->opts.bOPT[i] = (gmx_strncasecmp(ptr1[i], "Y", 1) == 0);
3701         ir->opts.bTS[i]  = (gmx_strncasecmp(ptr2[i], "Y", 1) == 0);
3702     }
3703     nSAon     = str_nelem(is->SAon, MAXPTR, ptr1);
3704     nSAoff    = str_nelem(is->SAoff, MAXPTR, ptr2);
3705     nSAsteps  = str_nelem(is->SAsteps, MAXPTR, ptr3);
3706     snew(ir->opts.SAon, nr);
3707     snew(ir->opts.SAoff, nr);
3708     snew(ir->opts.SAsteps, nr);
3709
3710     for (i = 0; i < nr; i++)
3711     {
3712         ir->opts.SAon[i]    = strtod(ptr1[i], NULL);
3713         ir->opts.SAoff[i]   = strtod(ptr2[i], NULL);
3714         ir->opts.SAsteps[i] = strtol(ptr3[i], NULL, 10);
3715     }
3716     /* end of QMMM input */
3717
3718     if (bVerbose)
3719     {
3720         for (i = 0; (i < egcNR); i++)
3721         {
3722             fprintf(stderr, "%-16s has %d element(s):", gtypes[i], groups->grps[i].nr);
3723             for (j = 0; (j < groups->grps[i].nr); j++)
3724             {
3725                 fprintf(stderr, " %s", *(groups->grpname[groups->grps[i].nm_ind[j]]));
3726             }
3727             fprintf(stderr, "\n");
3728         }
3729     }
3730
3731     nr = groups->grps[egcENER].nr;
3732     snew(ir->opts.egp_flags, nr*nr);
3733
3734     bExcl = do_egp_flag(ir, groups, "energygrp-excl", is->egpexcl, EGP_EXCL);
3735     if (bExcl && ir->cutoff_scheme == ecutsVERLET)
3736     {
3737         warning_error(wi, "Energy group exclusions are not (yet) implemented for the Verlet scheme");
3738     }
3739     if (bExcl && EEL_FULL(ir->coulombtype))
3740     {
3741         warning(wi, "Can not exclude the lattice Coulomb energy between energy groups");
3742     }
3743
3744     bTable = do_egp_flag(ir, groups, "energygrp-table", is->egptable, EGP_TABLE);
3745     if (bTable && !(ir->vdwtype == evdwUSER) &&
3746         !(ir->coulombtype == eelUSER) && !(ir->coulombtype == eelPMEUSER) &&
3747         !(ir->coulombtype == eelPMEUSERSWITCH))
3748     {
3749         gmx_fatal(FARGS, "Can only have energy group pair tables in combination with user tables for VdW and/or Coulomb");
3750     }
3751
3752     decode_cos(is->efield_x, &(ir->ex[XX]));
3753     decode_cos(is->efield_xt, &(ir->et[XX]));
3754     decode_cos(is->efield_y, &(ir->ex[YY]));
3755     decode_cos(is->efield_yt, &(ir->et[YY]));
3756     decode_cos(is->efield_z, &(ir->ex[ZZ]));
3757     decode_cos(is->efield_zt, &(ir->et[ZZ]));
3758
3759     if (ir->bAdress)
3760     {
3761         do_adress_index(ir->adress, groups, gnames, &(ir->opts), wi);
3762     }
3763
3764     for (i = 0; (i < grps->nr); i++)
3765     {
3766         sfree(gnames[i]);
3767     }
3768     sfree(gnames);
3769     done_blocka(grps);
3770     sfree(grps);
3771
3772 }
3773
3774
3775
3776 static void check_disre(gmx_mtop_t *mtop)
3777 {
3778     gmx_ffparams_t *ffparams;
3779     t_functype     *functype;
3780     t_iparams      *ip;
3781     int             i, ndouble, ftype;
3782     int             label, old_label;
3783
3784     if (gmx_mtop_ftype_count(mtop, F_DISRES) > 0)
3785     {
3786         ffparams  = &mtop->ffparams;
3787         functype  = ffparams->functype;
3788         ip        = ffparams->iparams;
3789         ndouble   = 0;
3790         old_label = -1;
3791         for (i = 0; i < ffparams->ntypes; i++)
3792         {
3793             ftype = functype[i];
3794             if (ftype == F_DISRES)
3795             {
3796                 label = ip[i].disres.label;
3797                 if (label == old_label)
3798                 {
3799                     fprintf(stderr, "Distance restraint index %d occurs twice\n", label);
3800                     ndouble++;
3801                 }
3802                 old_label = label;
3803             }
3804         }
3805         if (ndouble > 0)
3806         {
3807             gmx_fatal(FARGS, "Found %d double distance restraint indices,\n"
3808                       "probably the parameters for multiple pairs in one restraint "
3809                       "are not identical\n", ndouble);
3810         }
3811     }
3812 }
3813
3814 static gmx_bool absolute_reference(t_inputrec *ir, gmx_mtop_t *sys,
3815                                    gmx_bool posres_only,
3816                                    ivec AbsRef)
3817 {
3818     int                  d, g, i;
3819     gmx_mtop_ilistloop_t iloop;
3820     t_ilist             *ilist;
3821     int                  nmol;
3822     t_iparams           *pr;
3823
3824     clear_ivec(AbsRef);
3825
3826     if (!posres_only)
3827     {
3828         /* Check the COM */
3829         for (d = 0; d < DIM; d++)
3830         {
3831             AbsRef[d] = (d < ndof_com(ir) ? 0 : 1);
3832         }
3833         /* Check for freeze groups */
3834         for (g = 0; g < ir->opts.ngfrz; g++)
3835         {
3836             for (d = 0; d < DIM; d++)
3837             {
3838                 if (ir->opts.nFreeze[g][d] != 0)
3839                 {
3840                     AbsRef[d] = 1;
3841                 }
3842             }
3843         }
3844     }
3845
3846     /* Check for position restraints */
3847     iloop = gmx_mtop_ilistloop_init(sys);
3848     while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
3849     {
3850         if (nmol > 0 &&
3851             (AbsRef[XX] == 0 || AbsRef[YY] == 0 || AbsRef[ZZ] == 0))
3852         {
3853             for (i = 0; i < ilist[F_POSRES].nr; i += 2)
3854             {
3855                 pr = &sys->ffparams.iparams[ilist[F_POSRES].iatoms[i]];
3856                 for (d = 0; d < DIM; d++)
3857                 {
3858                     if (pr->posres.fcA[d] != 0)
3859                     {
3860                         AbsRef[d] = 1;
3861                     }
3862                 }
3863             }
3864             for (i = 0; i < ilist[F_FBPOSRES].nr; i += 2)
3865             {
3866                 /* Check for flat-bottom posres */
3867                 pr = &sys->ffparams.iparams[ilist[F_FBPOSRES].iatoms[i]];
3868                 if (pr->fbposres.k != 0)
3869                 {
3870                     switch (pr->fbposres.geom)
3871                     {
3872                         case efbposresSPHERE:
3873                             AbsRef[XX] = AbsRef[YY] = AbsRef[ZZ] = 1;
3874                             break;
3875                         case efbposresCYLINDERX:
3876                             AbsRef[YY] = AbsRef[ZZ] = 1;
3877                             break;
3878                         case efbposresCYLINDERY:
3879                             AbsRef[XX] = AbsRef[ZZ] = 1;
3880                             break;
3881                         case efbposresCYLINDER:
3882                         /* efbposres is a synonym for efbposresCYLINDERZ for backwards compatibility */
3883                         case efbposresCYLINDERZ:
3884                             AbsRef[XX] = AbsRef[YY] = 1;
3885                             break;
3886                         case efbposresX: /* d=XX */
3887                         case efbposresY: /* d=YY */
3888                         case efbposresZ: /* d=ZZ */
3889                             d         = pr->fbposres.geom - efbposresX;
3890                             AbsRef[d] = 1;
3891                             break;
3892                         default:
3893                             gmx_fatal(FARGS, " Invalid geometry for flat-bottom position restraint.\n"
3894                                       "Expected nr between 1 and %d. Found %d\n", efbposresNR-1,
3895                                       pr->fbposres.geom);
3896                     }
3897                 }
3898             }
3899         }
3900     }
3901
3902     return (AbsRef[XX] != 0 && AbsRef[YY] != 0 && AbsRef[ZZ] != 0);
3903 }
3904
3905 static void
3906 check_combination_rule_differences(const gmx_mtop_t *mtop, int state,
3907                                    gmx_bool *bC6ParametersWorkWithGeometricRules,
3908                                    gmx_bool *bC6ParametersWorkWithLBRules,
3909                                    gmx_bool *bLBRulesPossible)
3910 {
3911     int           ntypes, tpi, tpj, thisLBdiff, thisgeomdiff;
3912     int          *typecount;
3913     real          tol;
3914     double        geometricdiff, LBdiff;
3915     double        c6i, c6j, c12i, c12j;
3916     double        c6, c6_geometric, c6_LB;
3917     double        sigmai, sigmaj, epsi, epsj;
3918     gmx_bool      bCanDoLBRules, bCanDoGeometricRules;
3919     const char   *ptr;
3920
3921     /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
3922      * force-field floating point parameters.
3923      */
3924     tol = 1e-5;
3925     ptr = getenv("GMX_LJCOMB_TOL");
3926     if (ptr != NULL)
3927     {
3928         double dbl;
3929
3930         sscanf(ptr, "%lf", &dbl);
3931         tol = dbl;
3932     }
3933
3934     *bC6ParametersWorkWithLBRules         = TRUE;
3935     *bC6ParametersWorkWithGeometricRules  = TRUE;
3936     bCanDoLBRules                         = TRUE;
3937     bCanDoGeometricRules                  = TRUE;
3938     ntypes                                = mtop->ffparams.atnr;
3939     snew(typecount, ntypes);
3940     gmx_mtop_count_atomtypes(mtop, state, typecount);
3941     geometricdiff           = LBdiff = 0.0;
3942     *bLBRulesPossible       = TRUE;
3943     for (tpi = 0; tpi < ntypes; ++tpi)
3944     {
3945         c6i  = mtop->ffparams.iparams[(ntypes + 1) * tpi].lj.c6;
3946         c12i = mtop->ffparams.iparams[(ntypes + 1) * tpi].lj.c12;
3947         for (tpj = tpi; tpj < ntypes; ++tpj)
3948         {
3949             c6j          = mtop->ffparams.iparams[(ntypes + 1) * tpj].lj.c6;
3950             c12j         = mtop->ffparams.iparams[(ntypes + 1) * tpj].lj.c12;
3951             c6           = mtop->ffparams.iparams[ntypes * tpi + tpj].lj.c6;
3952             c6_geometric = sqrt(c6i * c6j);
3953             if (!gmx_numzero(c6_geometric))
3954             {
3955                 if (!gmx_numzero(c12i) && !gmx_numzero(c12j))
3956                 {
3957                     sigmai   = pow(c12i / c6i, 1.0/6.0);
3958                     sigmaj   = pow(c12j / c6j, 1.0/6.0);
3959                     epsi     = c6i * c6i /(4.0 * c12i);
3960                     epsj     = c6j * c6j /(4.0 * c12j);
3961                     c6_LB    = 4.0 * pow(epsi * epsj, 1.0/2.0) * pow(0.5 * (sigmai + sigmaj), 6);
3962                 }
3963                 else
3964                 {
3965                     *bLBRulesPossible = FALSE;
3966                     c6_LB             = c6_geometric;
3967                 }
3968                 bCanDoLBRules = gmx_within_tol(c6_LB, c6, tol);
3969             }
3970
3971             if (FALSE == bCanDoLBRules)
3972             {
3973                 *bC6ParametersWorkWithLBRules = FALSE;
3974             }
3975
3976             bCanDoGeometricRules = gmx_within_tol(c6_geometric, c6, tol);
3977
3978             if (FALSE == bCanDoGeometricRules)
3979             {
3980                 *bC6ParametersWorkWithGeometricRules = FALSE;
3981             }
3982         }
3983     }
3984     sfree(typecount);
3985 }
3986
3987 static void
3988 check_combination_rules(const t_inputrec *ir, const gmx_mtop_t *mtop,
3989                         warninp_t wi)
3990 {
3991     char     err_buf[256];
3992     gmx_bool bLBRulesPossible, bC6ParametersWorkWithGeometricRules, bC6ParametersWorkWithLBRules;
3993
3994     check_combination_rule_differences(mtop, 0,
3995                                        &bC6ParametersWorkWithGeometricRules,
3996                                        &bC6ParametersWorkWithLBRules,
3997                                        &bLBRulesPossible);
3998     if (ir->ljpme_combination_rule == eljpmeLB)
3999     {
4000         if (FALSE == bC6ParametersWorkWithLBRules || FALSE == bLBRulesPossible)
4001         {
4002             warning(wi, "You are using arithmetic-geometric combination rules "
4003                     "in LJ-PME, but your non-bonded C6 parameters do not "
4004                     "follow these rules.");
4005         }
4006     }
4007     else
4008     {
4009         if (FALSE == bC6ParametersWorkWithGeometricRules)
4010         {
4011             if (ir->eDispCorr != edispcNO)
4012             {
4013                 warning_note(wi, "You are using geometric combination rules in "
4014                              "LJ-PME, but your non-bonded C6 parameters do "
4015                              "not follow these rules. "
4016                              "This will introduce very small errors in the forces and energies in "
4017                              "your simulations. Dispersion correction will correct total energy "
4018                              "and/or pressure for isotropic systems, but not forces or surface tensions.");
4019             }
4020             else
4021             {
4022                 warning_note(wi, "You are using geometric combination rules in "
4023                              "LJ-PME, but your non-bonded C6 parameters do "
4024                              "not follow these rules. "
4025                              "This will introduce very small errors in the forces and energies in "
4026                              "your simulations. If your system is homogeneous, consider using dispersion correction "
4027                              "for the total energy and pressure.");
4028             }
4029         }
4030     }
4031 }
4032
4033 void triple_check(const char *mdparin, t_inputrec *ir, gmx_mtop_t *sys,
4034                   warninp_t wi)
4035 {
4036     char                      err_buf[STRLEN];
4037     int                       i, m, c, nmol, npct;
4038     gmx_bool                  bCharge, bAcc;
4039     real                      gdt_max, *mgrp, mt;
4040     rvec                      acc;
4041     gmx_mtop_atomloop_block_t aloopb;
4042     gmx_mtop_atomloop_all_t   aloop;
4043     t_atom                   *atom;
4044     ivec                      AbsRef;
4045     char                      warn_buf[STRLEN];
4046
4047     set_warning_line(wi, mdparin, -1);
4048
4049     if (ir->cutoff_scheme == ecutsVERLET &&
4050         ir->verletbuf_tol > 0 &&
4051         ir->nstlist > 1 &&
4052         ((EI_MD(ir->eI) || EI_SD(ir->eI)) &&
4053          (ir->etc == etcVRESCALE || ir->etc == etcBERENDSEN)))
4054     {
4055         /* Check if a too small Verlet buffer might potentially
4056          * cause more drift than the thermostat can couple off.
4057          */
4058         /* Temperature error fraction for warning and suggestion */
4059         const real T_error_warn    = 0.002;
4060         const real T_error_suggest = 0.001;
4061         /* For safety: 2 DOF per atom (typical with constraints) */
4062         const real nrdf_at         = 2;
4063         real       T, tau, max_T_error;
4064         int        i;
4065
4066         T   = 0;
4067         tau = 0;
4068         for (i = 0; i < ir->opts.ngtc; i++)
4069         {
4070             T   = max(T, ir->opts.ref_t[i]);
4071             tau = max(tau, ir->opts.tau_t[i]);
4072         }
4073         if (T > 0)
4074         {
4075             /* This is a worst case estimate of the temperature error,
4076              * assuming perfect buffer estimation and no cancelation
4077              * of errors. The factor 0.5 is because energy distributes
4078              * equally over Ekin and Epot.
4079              */
4080             max_T_error = 0.5*tau*ir->verletbuf_tol/(nrdf_at*BOLTZ*T);
4081             if (max_T_error > T_error_warn)
4082             {
4083                 sprintf(warn_buf, "With a verlet-buffer-tolerance of %g kJ/mol/ps, a reference temperature of %g and a tau_t of %g, your temperature might be off by up to %.1f%%. To ensure the error is below %.1f%%, decrease verlet-buffer-tolerance to %.0e or decrease tau_t.",
4084                         ir->verletbuf_tol, T, tau,
4085                         100*max_T_error,
4086                         100*T_error_suggest,
4087                         ir->verletbuf_tol*T_error_suggest/max_T_error);
4088                 warning(wi, warn_buf);
4089             }
4090         }
4091     }
4092
4093     if (ETC_ANDERSEN(ir->etc))
4094     {
4095         int i;
4096
4097         for (i = 0; i < ir->opts.ngtc; i++)
4098         {
4099             sprintf(err_buf, "all tau_t must currently be equal using Andersen temperature control, violated for group %d", i);
4100             CHECK(ir->opts.tau_t[0] != ir->opts.tau_t[i]);
4101             sprintf(err_buf, "all tau_t must be postive using Andersen temperature control, tau_t[%d]=%10.6f",
4102                     i, ir->opts.tau_t[i]);
4103             CHECK(ir->opts.tau_t[i] < 0);
4104         }
4105
4106         for (i = 0; i < ir->opts.ngtc; i++)
4107         {
4108             int nsteps = (int)(ir->opts.tau_t[i]/ir->delta_t);
4109             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);
4110             CHECK((nsteps % ir->nstcomm) && (ir->etc == etcANDERSENMASSIVE));
4111         }
4112     }
4113
4114     if (EI_DYNAMICS(ir->eI) && !EI_SD(ir->eI) && ir->eI != eiBD &&
4115         ir->comm_mode == ecmNO &&
4116         !(absolute_reference(ir, sys, FALSE, AbsRef) || ir->nsteps <= 10) &&
4117         !ETC_ANDERSEN(ir->etc))
4118     {
4119         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");
4120     }
4121
4122     /* Check for pressure coupling with absolute position restraints */
4123     if (ir->epc != epcNO && ir->refcoord_scaling == erscNO)
4124     {
4125         absolute_reference(ir, sys, TRUE, AbsRef);
4126         {
4127             for (m = 0; m < DIM; m++)
4128             {
4129                 if (AbsRef[m] && norm2(ir->compress[m]) > 0)
4130                 {
4131                     warning(wi, "You are using pressure coupling with absolute position restraints, this will give artifacts. Use the refcoord_scaling option.");
4132                     break;
4133                 }
4134             }
4135         }
4136     }
4137
4138     bCharge = FALSE;
4139     aloopb  = gmx_mtop_atomloop_block_init(sys);
4140     while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
4141     {
4142         if (atom->q != 0 || atom->qB != 0)
4143         {
4144             bCharge = TRUE;
4145         }
4146     }
4147
4148     if (!bCharge)
4149     {
4150         if (EEL_FULL(ir->coulombtype))
4151         {
4152             sprintf(err_buf,
4153                     "You are using full electrostatics treatment %s for a system without charges.\n"
4154                     "This costs a lot of performance for just processing zeros, consider using %s instead.\n",
4155                     EELTYPE(ir->coulombtype), EELTYPE(eelCUT));
4156             warning(wi, err_buf);
4157         }
4158     }
4159     else
4160     {
4161         if (ir->coulombtype == eelCUT && ir->rcoulomb > 0 && !ir->implicit_solvent)
4162         {
4163             sprintf(err_buf,
4164                     "You are using a plain Coulomb cut-off, which might produce artifacts.\n"
4165                     "You might want to consider using %s electrostatics.\n",
4166                     EELTYPE(eelPME));
4167             warning_note(wi, err_buf);
4168         }
4169     }
4170
4171     /* Check if combination rules used in LJ-PME are the same as in the force field */
4172     if (EVDW_PME(ir->vdwtype))
4173     {
4174         check_combination_rules(ir, sys, wi);
4175     }
4176
4177     /* Generalized reaction field */
4178     if (ir->opts.ngtc == 0)
4179     {
4180         sprintf(err_buf, "No temperature coupling while using coulombtype %s",
4181                 eel_names[eelGRF]);
4182         CHECK(ir->coulombtype == eelGRF);
4183     }
4184     else
4185     {
4186         sprintf(err_buf, "When using coulombtype = %s"
4187                 " ref-t for temperature coupling should be > 0",
4188                 eel_names[eelGRF]);
4189         CHECK((ir->coulombtype == eelGRF) && (ir->opts.ref_t[0] <= 0));
4190     }
4191
4192     if (ir->eI == eiSD2)
4193     {
4194         sprintf(warn_buf, "The stochastic dynamics integrator %s is deprecated, since\n"
4195                 "it is slower than integrator %s and is slightly less accurate\n"
4196                 "with constraints. Use the %s integrator.",
4197                 ei_names[ir->eI], ei_names[eiSD1], ei_names[eiSD1]);
4198         warning_note(wi, warn_buf);
4199     }
4200
4201     bAcc = FALSE;
4202     for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4203     {
4204         for (m = 0; (m < DIM); m++)
4205         {
4206             if (fabs(ir->opts.acc[i][m]) > 1e-6)
4207             {
4208                 bAcc = TRUE;
4209             }
4210         }
4211     }
4212     if (bAcc)
4213     {
4214         clear_rvec(acc);
4215         snew(mgrp, sys->groups.grps[egcACC].nr);
4216         aloop = gmx_mtop_atomloop_all_init(sys);
4217         while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
4218         {
4219             mgrp[ggrpnr(&sys->groups, egcACC, i)] += atom->m;
4220         }
4221         mt = 0.0;
4222         for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4223         {
4224             for (m = 0; (m < DIM); m++)
4225             {
4226                 acc[m] += ir->opts.acc[i][m]*mgrp[i];
4227             }
4228             mt += mgrp[i];
4229         }
4230         for (m = 0; (m < DIM); m++)
4231         {
4232             if (fabs(acc[m]) > 1e-6)
4233             {
4234                 const char *dim[DIM] = { "X", "Y", "Z" };
4235                 fprintf(stderr,
4236                         "Net Acceleration in %s direction, will %s be corrected\n",
4237                         dim[m], ir->nstcomm != 0 ? "" : "not");
4238                 if (ir->nstcomm != 0 && m < ndof_com(ir))
4239                 {
4240                     acc[m] /= mt;
4241                     for (i = 0; (i < sys->groups.grps[egcACC].nr); i++)
4242                     {
4243                         ir->opts.acc[i][m] -= acc[m];
4244                     }
4245                 }
4246             }
4247         }
4248         sfree(mgrp);
4249     }
4250
4251     if (ir->efep != efepNO && ir->fepvals->sc_alpha != 0 &&
4252         !gmx_within_tol(sys->ffparams.reppow, 12.0, 10*GMX_DOUBLE_EPS))
4253     {
4254         gmx_fatal(FARGS, "Soft-core interactions are only supported with VdW repulsion power 12");
4255     }
4256
4257     if (ir->bPull)
4258     {
4259         gmx_bool bWarned;
4260
4261         bWarned = FALSE;
4262         for (i = 0; i < ir->pull->ncoord && !bWarned; i++)
4263         {
4264             if (ir->pull->coord[i].group[0] == 0 ||
4265                 ir->pull->coord[i].group[1] == 0)
4266             {
4267                 absolute_reference(ir, sys, FALSE, AbsRef);
4268                 for (m = 0; m < DIM; m++)
4269                 {
4270                     if (ir->pull->coord[i].dim[m] && !AbsRef[m])
4271                     {
4272                         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.");
4273                         bWarned = TRUE;
4274                         break;
4275                     }
4276                 }
4277             }
4278         }
4279
4280         for (i = 0; i < 3; i++)
4281         {
4282             for (m = 0; m <= i; m++)
4283             {
4284                 if ((ir->epc != epcNO && ir->compress[i][m] != 0) ||
4285                     ir->deform[i][m] != 0)
4286                 {
4287                     for (c = 0; c < ir->pull->ncoord; c++)
4288                     {
4289                         if (ir->pull->coord[c].eGeom == epullgDIRPBC &&
4290                             ir->pull->coord[c].vec[m] != 0)
4291                         {
4292                             gmx_fatal(FARGS, "Can not have dynamic box while using pull geometry '%s' (dim %c)", EPULLGEOM(ir->pull->coord[c].eGeom), 'x'+m);
4293                         }
4294                     }
4295                 }
4296             }
4297         }
4298     }
4299
4300     check_disre(sys);
4301 }
4302
4303 void double_check(t_inputrec *ir, matrix box,
4304                   gmx_bool bHasNormalConstraints,
4305                   gmx_bool bHasAnyConstraints,
4306                   warninp_t wi)
4307 {
4308     real        min_size;
4309     gmx_bool    bTWIN;
4310     char        warn_buf[STRLEN];
4311     const char *ptr;
4312
4313     ptr = check_box(ir->ePBC, box);
4314     if (ptr)
4315     {
4316         warning_error(wi, ptr);
4317     }
4318
4319     if (bHasNormalConstraints && ir->eConstrAlg == econtSHAKE)
4320     {
4321         if (ir->shake_tol <= 0.0)
4322         {
4323             sprintf(warn_buf, "ERROR: shake-tol must be > 0 instead of %g\n",
4324                     ir->shake_tol);
4325             warning_error(wi, warn_buf);
4326         }
4327
4328         if (IR_TWINRANGE(*ir) && ir->nstlist > 1)
4329         {
4330             sprintf(warn_buf, "With twin-range cut-off's and SHAKE the virial and the pressure are incorrect.");
4331             if (ir->epc == epcNO)
4332             {
4333                 warning(wi, warn_buf);
4334             }
4335             else
4336             {
4337                 warning_error(wi, warn_buf);
4338             }
4339         }
4340     }
4341
4342     if ( (ir->eConstrAlg == econtLINCS) && bHasNormalConstraints)
4343     {
4344         /* If we have Lincs constraints: */
4345         if (ir->eI == eiMD && ir->etc == etcNO &&
4346             ir->eConstrAlg == econtLINCS && ir->nLincsIter == 1)
4347         {
4348             sprintf(warn_buf, "For energy conservation with LINCS, lincs_iter should be 2 or larger.\n");
4349             warning_note(wi, warn_buf);
4350         }
4351
4352         if ((ir->eI == eiCG || ir->eI == eiLBFGS) && (ir->nProjOrder < 8))
4353         {
4354             sprintf(warn_buf, "For accurate %s with LINCS constraints, lincs-order should be 8 or more.", ei_names[ir->eI]);
4355             warning_note(wi, warn_buf);
4356         }
4357         if (ir->epc == epcMTTK)
4358         {
4359             warning_error(wi, "MTTK not compatible with lincs -- use shake instead.");
4360         }
4361     }
4362
4363     if (bHasAnyConstraints && ir->epc == epcMTTK)
4364     {
4365         warning_error(wi, "Constraints are not implemented with MTTK pressure control.");
4366     }
4367
4368     if (ir->LincsWarnAngle > 90.0)
4369     {
4370         sprintf(warn_buf, "lincs-warnangle can not be larger than 90 degrees, setting it to 90.\n");
4371         warning(wi, warn_buf);
4372         ir->LincsWarnAngle = 90.0;
4373     }
4374
4375     if (ir->ePBC != epbcNONE)
4376     {
4377         if (ir->nstlist == 0)
4378         {
4379             warning(wi, "With nstlist=0 atoms are only put into the box at step 0, therefore drifting atoms might cause the simulation to crash.");
4380         }
4381         bTWIN = (ir->rlistlong > ir->rlist);
4382         if (ir->ns_type == ensGRID)
4383         {
4384             if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
4385             {
4386                 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",
4387                         bTWIN ? (ir->rcoulomb == ir->rlistlong ? "rcoulomb" : "rvdw") : "rlist");
4388                 warning_error(wi, warn_buf);
4389             }
4390         }
4391         else
4392         {
4393             min_size = min(box[XX][XX], min(box[YY][YY], box[ZZ][ZZ]));
4394             if (2*ir->rlistlong >= min_size)
4395             {
4396                 sprintf(warn_buf, "ERROR: One of the box lengths is smaller than twice the cut-off length. Increase the box size or decrease rlist.");
4397                 warning_error(wi, warn_buf);
4398                 if (TRICLINIC(box))
4399                 {
4400                     fprintf(stderr, "Grid search might allow larger cut-off's than simple search with triclinic boxes.");
4401                 }
4402             }
4403         }
4404     }
4405 }
4406
4407 void check_chargegroup_radii(const gmx_mtop_t *mtop, const t_inputrec *ir,
4408                              rvec *x,
4409                              warninp_t wi)
4410 {
4411     real rvdw1, rvdw2, rcoul1, rcoul2;
4412     char warn_buf[STRLEN];
4413
4414     calc_chargegroup_radii(mtop, x, &rvdw1, &rvdw2, &rcoul1, &rcoul2);
4415
4416     if (rvdw1 > 0)
4417     {
4418         printf("Largest charge group radii for Van der Waals: %5.3f, %5.3f nm\n",
4419                rvdw1, rvdw2);
4420     }
4421     if (rcoul1 > 0)
4422     {
4423         printf("Largest charge group radii for Coulomb:       %5.3f, %5.3f nm\n",
4424                rcoul1, rcoul2);
4425     }
4426
4427     if (ir->rlist > 0)
4428     {
4429         if (rvdw1  + rvdw2  > ir->rlist ||
4430             rcoul1 + rcoul2 > ir->rlist)
4431         {
4432             sprintf(warn_buf,
4433                     "The sum of the two largest charge group radii (%f) "
4434                     "is larger than rlist (%f)\n",
4435                     max(rvdw1+rvdw2, rcoul1+rcoul2), ir->rlist);
4436             warning(wi, warn_buf);
4437         }
4438         else
4439         {
4440             /* Here we do not use the zero at cut-off macro,
4441              * since user defined interactions might purposely
4442              * not be zero at the cut-off.
4443              */
4444             if (ir_vdw_is_zero_at_cutoff(ir) &&
4445                 rvdw1 + rvdw2 > ir->rlistlong - ir->rvdw)
4446             {
4447                 sprintf(warn_buf, "The sum of the two largest charge group "
4448                         "radii (%f) is larger than %s (%f) - rvdw (%f).\n"
4449                         "With exact cut-offs, better performance can be "
4450                         "obtained with cutoff-scheme = %s, because it "
4451                         "does not use charge groups at all.",
4452                         rvdw1+rvdw2,
4453                         ir->rlistlong > ir->rlist ? "rlistlong" : "rlist",
4454                         ir->rlistlong, ir->rvdw,
4455                         ecutscheme_names[ecutsVERLET]);
4456                 if (ir_NVE(ir))
4457                 {
4458                     warning(wi, warn_buf);
4459                 }
4460                 else
4461                 {
4462                     warning_note(wi, warn_buf);
4463                 }
4464             }
4465             if (ir_coulomb_is_zero_at_cutoff(ir) &&
4466                 rcoul1 + rcoul2 > ir->rlistlong - ir->rcoulomb)
4467             {
4468                 sprintf(warn_buf, "The sum of the two largest charge group radii (%f) is larger than %s (%f) - rcoulomb (%f).\n"
4469                         "With exact cut-offs, better performance can be obtained with cutoff-scheme = %s, because it does not use charge groups at all.",
4470                         rcoul1+rcoul2,
4471                         ir->rlistlong > ir->rlist ? "rlistlong" : "rlist",
4472                         ir->rlistlong, ir->rcoulomb,
4473                         ecutscheme_names[ecutsVERLET]);
4474                 if (ir_NVE(ir))
4475                 {
4476                     warning(wi, warn_buf);
4477                 }
4478                 else
4479                 {
4480                     warning_note(wi, warn_buf);
4481                 }
4482             }
4483         }
4484     }
4485 }