Update advice about grompp -t
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / grompp.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 "grompp.h"
38
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42
43 #include <sys/types.h>
44 #include <math.h>
45 #include <string.h>
46 #include <errno.h>
47 #include <limits.h>
48 #include <assert.h>
49
50 #include "sysstuff.h"
51 #include "gromacs/utility/smalloc.h"
52 #include "macros.h"
53 #include "readir.h"
54 #include "toputil.h"
55 #include "topio.h"
56 #include "gromacs/fileio/confio.h"
57 #include "readir.h"
58 #include "symtab.h"
59 #include "names.h"
60 #include "grompp-impl.h"
61 #include "gromacs/random/random.h"
62 #include "gromacs/gmxpreprocess/gen_maxwell_velocities.h"
63 #include "vec.h"
64 #include "gromacs/fileio/futil.h"
65 #include "gromacs/commandline/pargs.h"
66 #include "splitter.h"
67 #include "gromacs/gmxpreprocess/sortwater.h"
68 #include "convparm.h"
69 #include "gmx_fatal.h"
70 #include "warninp.h"
71 #include "index.h"
72 #include "gromacs/fileio/gmxfio.h"
73 #include "gromacs/fileio/trnio.h"
74 #include "gromacs/fileio/tpxio.h"
75 #include "gromacs/fileio/trxio.h"
76 #include "vsite_parm.h"
77 #include "txtdump.h"
78 #include "calcgrid.h"
79 #include "add_par.h"
80 #include "gromacs/fileio/enxio.h"
81 #include "perf_est.h"
82 #include "compute_io.h"
83 #include "gpp_atomtype.h"
84 #include "mtop_util.h"
85 #include "genborn.h"
86 #include "calc_verletbuf.h"
87 #include "tomorse.h"
88 #include "gromacs/imd/imd.h"
89 #include "gromacs/utility/cstringutil.h"
90
91
92 static int rm_interactions(int ifunc, int nrmols, t_molinfo mols[])
93 {
94     int  i, n;
95
96     n = 0;
97     /* For all the molecule types */
98     for (i = 0; i < nrmols; i++)
99     {
100         n += mols[i].plist[ifunc].nr;
101         mols[i].plist[ifunc].nr = 0;
102     }
103     return n;
104 }
105
106 static int check_atom_names(const char *fn1, const char *fn2,
107                             gmx_mtop_t *mtop, t_atoms *at)
108 {
109     int      mb, m, i, j, nmismatch;
110     t_atoms *tat;
111 #define MAXMISMATCH 20
112
113     if (mtop->natoms != at->nr)
114     {
115         gmx_incons("comparing atom names");
116     }
117
118     nmismatch = 0;
119     i         = 0;
120     for (mb = 0; mb < mtop->nmolblock; mb++)
121     {
122         tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
123         for (m = 0; m < mtop->molblock[mb].nmol; m++)
124         {
125             for (j = 0; j < tat->nr; j++)
126             {
127                 if (strcmp( *(tat->atomname[j]), *(at->atomname[i]) ) != 0)
128                 {
129                     if (nmismatch < MAXMISMATCH)
130                     {
131                         fprintf(stderr,
132                                 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
133                                 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
134                     }
135                     else if (nmismatch == MAXMISMATCH)
136                     {
137                         fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
138                     }
139                     nmismatch++;
140                 }
141                 i++;
142             }
143         }
144     }
145
146     return nmismatch;
147 }
148
149 static void check_eg_vs_cg(gmx_mtop_t *mtop)
150 {
151     int            astart, mb, m, cg, j, firstj;
152     unsigned char  firsteg, eg;
153     gmx_moltype_t *molt;
154
155     /* Go through all the charge groups and make sure all their
156      * atoms are in the same energy group.
157      */
158
159     astart = 0;
160     for (mb = 0; mb < mtop->nmolblock; mb++)
161     {
162         molt = &mtop->moltype[mtop->molblock[mb].type];
163         for (m = 0; m < mtop->molblock[mb].nmol; m++)
164         {
165             for (cg = 0; cg < molt->cgs.nr; cg++)
166             {
167                 /* Get the energy group of the first atom in this charge group */
168                 firstj  = astart + molt->cgs.index[cg];
169                 firsteg = ggrpnr(&mtop->groups, egcENER, firstj);
170                 for (j = molt->cgs.index[cg]+1; j < molt->cgs.index[cg+1]; j++)
171                 {
172                     eg = ggrpnr(&mtop->groups, egcENER, astart+j);
173                     if (eg != firsteg)
174                     {
175                         gmx_fatal(FARGS, "atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
176                                   firstj+1, astart+j+1, cg+1, *molt->name);
177                     }
178                 }
179             }
180             astart += molt->atoms.nr;
181         }
182     }
183 }
184
185 static void check_cg_sizes(const char *topfn, t_block *cgs, warninp_t wi)
186 {
187     int  maxsize, cg;
188     char warn_buf[STRLEN];
189
190     maxsize = 0;
191     for (cg = 0; cg < cgs->nr; cg++)
192     {
193         maxsize = max(maxsize, cgs->index[cg+1]-cgs->index[cg]);
194     }
195
196     if (maxsize > MAX_CHARGEGROUP_SIZE)
197     {
198         gmx_fatal(FARGS, "The largest charge group contains %d atoms. The maximum is %d.", maxsize, MAX_CHARGEGROUP_SIZE);
199     }
200     else if (maxsize > 10)
201     {
202         set_warning_line(wi, topfn, -1);
203         sprintf(warn_buf,
204                 "The largest charge group contains %d atoms.\n"
205                 "Since atoms only see each other when the centers of geometry of the charge groups they belong to are within the cut-off distance, too large charge groups can lead to serious cut-off artifacts.\n"
206                 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
207                 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
208                 maxsize);
209         warning_note(wi, warn_buf);
210     }
211 }
212
213 static void check_bonds_timestep(gmx_mtop_t *mtop, double dt, warninp_t wi)
214 {
215     /* This check is not intended to ensure accurate integration,
216      * rather it is to signal mistakes in the mdp settings.
217      * A common mistake is to forget to turn on constraints
218      * for MD after energy minimization with flexible bonds.
219      * This check can also detect too large time steps for flexible water
220      * models, but such errors will often be masked by the constraints
221      * mdp options, which turns flexible water into water with bond constraints,
222      * but without an angle constraint. Unfortunately such incorrect use
223      * of water models can not easily be detected without checking
224      * for specific model names.
225      *
226      * The stability limit of leap-frog or velocity verlet is 4.44 steps
227      * per oscillational period.
228      * But accurate bonds distributions are lost far before that limit.
229      * To allow relatively common schemes (although not common with Gromacs)
230      * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
231      * we set the note limit to 10.
232      */
233     int            min_steps_warn = 5;
234     int            min_steps_note = 10;
235     t_iparams     *ip;
236     int            molt;
237     gmx_moltype_t *moltype, *w_moltype;
238     t_atom        *atom;
239     t_ilist       *ilist, *ilb, *ilc, *ils;
240     int            ftype;
241     int            i, a1, a2, w_a1, w_a2, j;
242     real           twopi2, limit2, fc, re, m1, m2, period2, w_period2;
243     gmx_bool       bFound, bWater, bWarn;
244     char           warn_buf[STRLEN];
245
246     ip = mtop->ffparams.iparams;
247
248     twopi2 = sqr(2*M_PI);
249
250     limit2 = sqr(min_steps_note*dt);
251
252     w_a1      = w_a2 = -1;
253     w_period2 = -1.0;
254
255     w_moltype = NULL;
256     for (molt = 0; molt < mtop->nmoltype; molt++)
257     {
258         moltype = &mtop->moltype[molt];
259         atom    = moltype->atoms.atom;
260         ilist   = moltype->ilist;
261         ilc     = &ilist[F_CONSTR];
262         ils     = &ilist[F_SETTLE];
263         for (ftype = 0; ftype < F_NRE; ftype++)
264         {
265             if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
266             {
267                 continue;
268             }
269
270             ilb = &ilist[ftype];
271             for (i = 0; i < ilb->nr; i += 3)
272             {
273                 fc = ip[ilb->iatoms[i]].harmonic.krA;
274                 re = ip[ilb->iatoms[i]].harmonic.rA;
275                 if (ftype == F_G96BONDS)
276                 {
277                     /* Convert squared sqaure fc to harmonic fc */
278                     fc = 2*fc*re;
279                 }
280                 a1 = ilb->iatoms[i+1];
281                 a2 = ilb->iatoms[i+2];
282                 m1 = atom[a1].m;
283                 m2 = atom[a2].m;
284                 if (fc > 0 && m1 > 0 && m2 > 0)
285                 {
286                     period2 = twopi2*m1*m2/((m1 + m2)*fc);
287                 }
288                 else
289                 {
290                     period2 = GMX_FLOAT_MAX;
291                 }
292                 if (debug)
293                 {
294                     fprintf(debug, "fc %g m1 %g m2 %g period %g\n",
295                             fc, m1, m2, sqrt(period2));
296                 }
297                 if (period2 < limit2)
298                 {
299                     bFound = FALSE;
300                     for (j = 0; j < ilc->nr; j += 3)
301                     {
302                         if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
303                             (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
304                         {
305                             bFound = TRUE;
306                         }
307                     }
308                     for (j = 0; j < ils->nr; j += 4)
309                     {
310                         if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
311                             (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
312                         {
313                             bFound = TRUE;
314                         }
315                     }
316                     if (!bFound &&
317                         (w_moltype == NULL || period2 < w_period2))
318                     {
319                         w_moltype = moltype;
320                         w_a1      = a1;
321                         w_a2      = a2;
322                         w_period2 = period2;
323                     }
324                 }
325             }
326         }
327     }
328
329     if (w_moltype != NULL)
330     {
331         bWarn = (w_period2 < sqr(min_steps_warn*dt));
332         /* A check that would recognize most water models */
333         bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
334                   w_moltype->atoms.nr <= 5);
335         sprintf(warn_buf, "The bond in molecule-type %s between atoms %d %s and %d %s has an estimated oscillational period of %.1e ps, which is less than %d times the time step of %.1e ps.\n"
336                 "%s",
337                 *w_moltype->name,
338                 w_a1+1, *w_moltype->atoms.atomname[w_a1],
339                 w_a2+1, *w_moltype->atoms.atomname[w_a2],
340                 sqrt(w_period2), bWarn ? min_steps_warn : min_steps_note, dt,
341                 bWater ?
342                 "Maybe you asked for fexible water." :
343                 "Maybe you forgot to change the constraints mdp option.");
344         if (bWarn)
345         {
346             warning(wi, warn_buf);
347         }
348         else
349         {
350             warning_note(wi, warn_buf);
351         }
352     }
353 }
354
355 static void check_vel(gmx_mtop_t *mtop, rvec v[])
356 {
357     gmx_mtop_atomloop_all_t aloop;
358     t_atom                 *atom;
359     int                     a;
360
361     aloop = gmx_mtop_atomloop_all_init(mtop);
362     while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
363     {
364         if (atom->ptype == eptShell ||
365             atom->ptype == eptBond  ||
366             atom->ptype == eptVSite)
367         {
368             clear_rvec(v[a]);
369         }
370     }
371 }
372
373 static void check_shells_inputrec(gmx_mtop_t *mtop,
374                                   t_inputrec *ir,
375                                   warninp_t   wi)
376 {
377     gmx_mtop_atomloop_all_t aloop;
378     t_atom                 *atom;
379     int                     a, nshells = 0;
380     char                    warn_buf[STRLEN];
381
382     aloop = gmx_mtop_atomloop_all_init(mtop);
383     while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
384     {
385         if (atom->ptype == eptShell ||
386             atom->ptype == eptBond)
387         {
388             nshells++;
389         }
390     }
391     if (IR_TWINRANGE(*ir) && (nshells > 0))
392     {
393         snprintf(warn_buf, STRLEN,
394                  "The combination of using shells and a twin-range cut-off is not supported");
395         warning_error(wi, warn_buf);
396     }
397     if ((nshells > 0) && (ir->nstcalcenergy != 1))
398     {
399         set_warning_line(wi, "unknown", -1);
400         snprintf(warn_buf, STRLEN,
401                  "There are %d shells, changing nstcalcenergy from %d to 1",
402                  nshells, ir->nstcalcenergy);
403         ir->nstcalcenergy = 1;
404         warning(wi, warn_buf);
405     }
406 }
407
408 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
409 {
410     int nint, mb;
411
412     nint = 0;
413     for (mb = 0; mb < mtop->nmolblock; mb++)
414     {
415         nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
416     }
417
418     return nint;
419 }
420
421 /* This routine reorders the molecule type array
422  * in the order of use in the molblocks,
423  * unused molecule types are deleted.
424  */
425 static void renumber_moltypes(gmx_mtop_t *sys,
426                               int *nmolinfo, t_molinfo **molinfo)
427 {
428     int       *order, norder, i;
429     int        mb, mi;
430     t_molinfo *minew;
431
432     snew(order, *nmolinfo);
433     norder = 0;
434     for (mb = 0; mb < sys->nmolblock; mb++)
435     {
436         for (i = 0; i < norder; i++)
437         {
438             if (order[i] == sys->molblock[mb].type)
439             {
440                 break;
441             }
442         }
443         if (i == norder)
444         {
445             /* This type did not occur yet, add it */
446             order[norder] = sys->molblock[mb].type;
447             /* Renumber the moltype in the topology */
448             norder++;
449         }
450         sys->molblock[mb].type = i;
451     }
452
453     /* We still need to reorder the molinfo structs */
454     snew(minew, norder);
455     for (mi = 0; mi < *nmolinfo; mi++)
456     {
457         for (i = 0; i < norder; i++)
458         {
459             if (order[i] == mi)
460             {
461                 break;
462             }
463         }
464         if (i == norder)
465         {
466             done_mi(&(*molinfo)[mi]);
467         }
468         else
469         {
470             minew[i] = (*molinfo)[mi];
471         }
472     }
473     sfree(*molinfo);
474
475     *nmolinfo = norder;
476     *molinfo  = minew;
477 }
478
479 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
480 {
481     int            m;
482     gmx_moltype_t *molt;
483
484     mtop->nmoltype = nmi;
485     snew(mtop->moltype, nmi);
486     for (m = 0; m < nmi; m++)
487     {
488         molt        = &mtop->moltype[m];
489         molt->name  = mi[m].name;
490         molt->atoms = mi[m].atoms;
491         /* ilists are copied later */
492         molt->cgs   = mi[m].cgs;
493         molt->excls = mi[m].excls;
494     }
495 }
496
497 static void
498 new_status(const char *topfile, const char *topppfile, const char *confin,
499            t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
500            gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
501            gpp_atomtype_t atype, gmx_mtop_t *sys,
502            int *nmi, t_molinfo **mi, t_params plist[],
503            int *comb, double *reppow, real *fudgeQQ,
504            gmx_bool bMorse,
505            warninp_t wi)
506 {
507     t_molinfo      *molinfo = NULL;
508     int             nmolblock;
509     gmx_molblock_t *molblock, *molbs;
510     t_atoms        *confat;
511     int             mb, i, nrmols, nmismatch;
512     char            buf[STRLEN];
513     gmx_bool        bGB = FALSE;
514     char            warn_buf[STRLEN];
515
516     init_mtop(sys);
517
518     /* Set gmx_boolean for GB */
519     if (ir->implicit_solvent)
520     {
521         bGB = TRUE;
522     }
523
524     /* TOPOLOGY processing */
525     sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
526                        plist, comb, reppow, fudgeQQ,
527                        atype, &nrmols, &molinfo, ir,
528                        &nmolblock, &molblock, bGB,
529                        wi);
530
531     sys->nmolblock = 0;
532     snew(sys->molblock, nmolblock);
533
534     sys->natoms = 0;
535     for (mb = 0; mb < nmolblock; mb++)
536     {
537         if (sys->nmolblock > 0 &&
538             molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
539         {
540             /* Merge consecutive blocks with the same molecule type */
541             sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
542             sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
543         }
544         else if (molblock[mb].nmol > 0)
545         {
546             /* Add a new molblock to the topology */
547             molbs             = &sys->molblock[sys->nmolblock];
548             *molbs            = molblock[mb];
549             molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
550             molbs->nposres_xA = 0;
551             molbs->nposres_xB = 0;
552             sys->natoms      += molbs->nmol*molbs->natoms_mol;
553             sys->nmolblock++;
554         }
555     }
556     if (sys->nmolblock == 0)
557     {
558         gmx_fatal(FARGS, "No molecules were defined in the system");
559     }
560
561     renumber_moltypes(sys, &nrmols, &molinfo);
562
563     if (bMorse)
564     {
565         convert_harmonics(nrmols, molinfo, atype);
566     }
567
568     if (ir->eDisre == edrNone)
569     {
570         i = rm_interactions(F_DISRES, nrmols, molinfo);
571         if (i > 0)
572         {
573             set_warning_line(wi, "unknown", -1);
574             sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
575             warning_note(wi, warn_buf);
576         }
577     }
578     if (opts->bOrire == FALSE)
579     {
580         i = rm_interactions(F_ORIRES, nrmols, molinfo);
581         if (i > 0)
582         {
583             set_warning_line(wi, "unknown", -1);
584             sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
585             warning_note(wi, warn_buf);
586         }
587     }
588
589     /* Copy structures from msys to sys */
590     molinfo2mtop(nrmols, molinfo, sys);
591
592     gmx_mtop_finalize(sys);
593
594     /* COORDINATE file processing */
595     if (bVerbose)
596     {
597         fprintf(stderr, "processing coordinates...\n");
598     }
599
600     get_stx_coordnum(confin, &state->natoms);
601     if (state->natoms != sys->natoms)
602     {
603         gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
604                   "             does not match topology (%s, %d)",
605                   confin, state->natoms, topfile, sys->natoms);
606     }
607     else
608     {
609         /* make space for coordinates and velocities */
610         char title[STRLEN];
611         snew(confat, 1);
612         init_t_atoms(confat, state->natoms, FALSE);
613         init_state(state, state->natoms, 0, 0, 0, 0);
614         read_stx_conf(confin, title, confat, state->x, state->v, NULL, state->box);
615         /* This call fixes the box shape for runs with pressure scaling */
616         set_box_rel(ir, state);
617
618         nmismatch = check_atom_names(topfile, confin, sys, confat);
619         free_t_atoms(confat, TRUE);
620         sfree(confat);
621
622         if (nmismatch)
623         {
624             sprintf(buf, "%d non-matching atom name%s\n"
625                     "atom names from %s will be used\n"
626                     "atom names from %s will be ignored\n",
627                     nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
628             warning(wi, buf);
629         }
630         if (bVerbose)
631         {
632             fprintf(stderr, "double-checking input for internal consistency...\n");
633         }
634         double_check(ir, state->box, nint_ftype(sys, molinfo, F_CONSTR), wi);
635     }
636
637     if (bGenVel)
638     {
639         real                   *mass;
640         gmx_mtop_atomloop_all_t aloop;
641         t_atom                 *atom;
642         unsigned int            useed;
643
644         snew(mass, state->natoms);
645         aloop = gmx_mtop_atomloop_all_init(sys);
646         while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
647         {
648             mass[i] = atom->m;
649         }
650
651         useed = opts->seed;
652         if (opts->seed == -1)
653         {
654             useed = (int)gmx_rng_make_seed();
655             fprintf(stderr, "Setting gen_seed to %u\n", useed);
656         }
657         maxwell_speed(opts->tempi, useed, sys, state->v);
658
659         stop_cm(stdout, state->natoms, mass, state->x, state->v);
660         sfree(mass);
661     }
662
663     *nmi = nrmols;
664     *mi  = molinfo;
665 }
666
667 static void copy_state(const char *slog, t_trxframe *fr,
668                        gmx_bool bReadVel, t_state *state,
669                        double *use_time)
670 {
671     int i;
672
673     if (fr->not_ok & FRAME_NOT_OK)
674     {
675         gmx_fatal(FARGS, "Can not start from an incomplete frame");
676     }
677     if (!fr->bX)
678     {
679         gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
680                   slog);
681     }
682
683     for (i = 0; i < state->natoms; i++)
684     {
685         copy_rvec(fr->x[i], state->x[i]);
686     }
687     if (bReadVel)
688     {
689         if (!fr->bV)
690         {
691             gmx_incons("Trajecory frame unexpectedly does not contain velocities");
692         }
693         for (i = 0; i < state->natoms; i++)
694         {
695             copy_rvec(fr->v[i], state->v[i]);
696         }
697     }
698     if (fr->bBox)
699     {
700         copy_mat(fr->box, state->box);
701     }
702
703     *use_time = fr->time;
704 }
705
706 static void cont_status(const char *slog, const char *ener,
707                         gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
708                         t_inputrec *ir, t_state *state,
709                         gmx_mtop_t *sys,
710                         const output_env_t oenv)
711 /* If fr_time == -1 read the last frame available which is complete */
712 {
713     gmx_bool     bReadVel;
714     t_trxframe   fr;
715     t_trxstatus *fp;
716     int          i;
717     double       use_time;
718
719     bReadVel = (bNeedVel && !bGenVel);
720
721     fprintf(stderr,
722             "Reading Coordinates%s and Box size from old trajectory\n",
723             bReadVel ? ", Velocities" : "");
724     if (fr_time == -1)
725     {
726         fprintf(stderr, "Will read whole trajectory\n");
727     }
728     else
729     {
730         fprintf(stderr, "Will read till time %g\n", fr_time);
731     }
732     if (!bReadVel)
733     {
734         if (bGenVel)
735         {
736             fprintf(stderr, "Velocities generated: "
737                     "ignoring velocities in input trajectory\n");
738         }
739         read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
740     }
741     else
742     {
743         read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
744
745         if (!fr.bV)
746         {
747             fprintf(stderr,
748                     "\n"
749                     "WARNING: Did not find a frame with velocities in file %s,\n"
750                     "         all velocities will be set to zero!\n\n", slog);
751             for (i = 0; i < sys->natoms; i++)
752             {
753                 clear_rvec(state->v[i]);
754             }
755             close_trj(fp);
756             /* Search for a frame without velocities */
757             bReadVel = FALSE;
758             read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
759         }
760     }
761
762     state->natoms = fr.natoms;
763
764     if (sys->natoms != state->natoms)
765     {
766         gmx_fatal(FARGS, "Number of atoms in Topology "
767                   "is not the same as in Trajectory");
768     }
769     copy_state(slog, &fr, bReadVel, state, &use_time);
770
771     /* Find the appropriate frame */
772     while ((fr_time == -1 || fr.time < fr_time) &&
773            read_next_frame(oenv, fp, &fr))
774     {
775         copy_state(slog, &fr, bReadVel, state, &use_time);
776     }
777
778     close_trj(fp);
779
780     /* Set the relative box lengths for preserving the box shape.
781      * Note that this call can lead to differences in the last bit
782      * with respect to using gmx convert-tpr to create a [TT].tpx[tt] file.
783      */
784     set_box_rel(ir, state);
785
786     fprintf(stderr, "Using frame at t = %g ps\n", use_time);
787     fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
788
789     if ((ir->epc != epcNO  || ir->etc == etcNOSEHOOVER) && ener)
790     {
791         get_enx_state(ener, use_time, &sys->groups, ir, state);
792         preserve_box_shape(ir, state->box_rel, state->boxv);
793     }
794 }
795
796 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
797                         char *fn,
798                         int rc_scaling, int ePBC,
799                         rvec com,
800                         warninp_t wi)
801 {
802     gmx_bool        bFirst = TRUE, *hadAtom;
803     rvec           *x, *v, *xp;
804     dvec            sum;
805     double          totmass;
806     t_atoms         dumat;
807     matrix          box, invbox;
808     int             natoms, npbcdim = 0;
809     char            warn_buf[STRLEN], title[STRLEN];
810     int             a, i, ai, j, k, mb, nat_molb;
811     gmx_molblock_t *molb;
812     t_params       *pr, *prfb;
813     t_atom         *atom;
814
815     get_stx_coordnum(fn, &natoms);
816     if (natoms != mtop->natoms)
817     {
818         sprintf(warn_buf, "The number of atoms in %s (%d) does not match the number of atoms in the topology (%d). Will assume that the first %d atoms in the topology and %s match.", fn, natoms, mtop->natoms, min(mtop->natoms, natoms), fn);
819         warning(wi, warn_buf);
820     }
821     snew(x, natoms);
822     snew(v, natoms);
823     init_t_atoms(&dumat, natoms, FALSE);
824     read_stx_conf(fn, title, &dumat, x, v, NULL, box);
825
826     npbcdim = ePBC2npbcdim(ePBC);
827     clear_rvec(com);
828     if (rc_scaling != erscNO)
829     {
830         copy_mat(box, invbox);
831         for (j = npbcdim; j < DIM; j++)
832         {
833             clear_rvec(invbox[j]);
834             invbox[j][j] = 1;
835         }
836         m_inv_ur0(invbox, invbox);
837     }
838
839     /* Copy the reference coordinates to mtop */
840     clear_dvec(sum);
841     totmass = 0;
842     a       = 0;
843     snew(hadAtom, natoms);
844     for (mb = 0; mb < mtop->nmolblock; mb++)
845     {
846         molb     = &mtop->molblock[mb];
847         nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
848         pr       = &(molinfo[molb->type].plist[F_POSRES]);
849         prfb     = &(molinfo[molb->type].plist[F_FBPOSRES]);
850         if (pr->nr > 0 || prfb->nr > 0)
851         {
852             atom = mtop->moltype[molb->type].atoms.atom;
853             for (i = 0; (i < pr->nr); i++)
854             {
855                 ai = pr->param[i].AI;
856                 if (ai >= natoms)
857                 {
858                     gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
859                               ai+1, *molinfo[molb->type].name, fn, natoms);
860                 }
861                 hadAtom[ai] = TRUE;
862                 if (rc_scaling == erscCOM)
863                 {
864                     /* Determine the center of mass of the posres reference coordinates */
865                     for (j = 0; j < npbcdim; j++)
866                     {
867                         sum[j] += atom[ai].m*x[a+ai][j];
868                     }
869                     totmass  += atom[ai].m;
870                 }
871             }
872             /* Same for flat-bottomed posres, but do not count an atom twice for COM */
873             for (i = 0; (i < prfb->nr); i++)
874             {
875                 ai = prfb->param[i].AI;
876                 if (ai >= natoms)
877                 {
878                     gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
879                               ai+1, *molinfo[molb->type].name, fn, natoms);
880                 }
881                 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE)
882                 {
883                     /* Determine the center of mass of the posres reference coordinates */
884                     for (j = 0; j < npbcdim; j++)
885                     {
886                         sum[j] += atom[ai].m*x[a+ai][j];
887                     }
888                     totmass  += atom[ai].m;
889                 }
890             }
891             if (!bTopB)
892             {
893                 molb->nposres_xA = nat_molb;
894                 snew(molb->posres_xA, molb->nposres_xA);
895                 for (i = 0; i < nat_molb; i++)
896                 {
897                     copy_rvec(x[a+i], molb->posres_xA[i]);
898                 }
899             }
900             else
901             {
902                 molb->nposres_xB = nat_molb;
903                 snew(molb->posres_xB, molb->nposres_xB);
904                 for (i = 0; i < nat_molb; i++)
905                 {
906                     copy_rvec(x[a+i], molb->posres_xB[i]);
907                 }
908             }
909         }
910         a += nat_molb;
911     }
912     if (rc_scaling == erscCOM)
913     {
914         if (totmass == 0)
915         {
916             gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
917         }
918         for (j = 0; j < npbcdim; j++)
919         {
920             com[j] = sum[j]/totmass;
921         }
922         fprintf(stderr, "The center of mass of the position restraint coord's is %6.3f %6.3f %6.3f\n", com[XX], com[YY], com[ZZ]);
923     }
924
925     if (rc_scaling != erscNO)
926     {
927         assert(npbcdim <= DIM);
928
929         for (mb = 0; mb < mtop->nmolblock; mb++)
930         {
931             molb     = &mtop->molblock[mb];
932             nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
933             if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
934             {
935                 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
936                 for (i = 0; i < nat_molb; i++)
937                 {
938                     for (j = 0; j < npbcdim; j++)
939                     {
940                         if (rc_scaling == erscALL)
941                         {
942                             /* Convert from Cartesian to crystal coordinates */
943                             xp[i][j] *= invbox[j][j];
944                             for (k = j+1; k < npbcdim; k++)
945                             {
946                                 xp[i][j] += invbox[k][j]*xp[i][k];
947                             }
948                         }
949                         else if (rc_scaling == erscCOM)
950                         {
951                             /* Subtract the center of mass */
952                             xp[i][j] -= com[j];
953                         }
954                     }
955                 }
956             }
957         }
958
959         if (rc_scaling == erscCOM)
960         {
961             /* Convert the COM from Cartesian to crystal coordinates */
962             for (j = 0; j < npbcdim; j++)
963             {
964                 com[j] *= invbox[j][j];
965                 for (k = j+1; k < npbcdim; k++)
966                 {
967                     com[j] += invbox[k][j]*com[k];
968                 }
969             }
970         }
971     }
972
973     free_t_atoms(&dumat, TRUE);
974     sfree(x);
975     sfree(v);
976     sfree(hadAtom);
977 }
978
979 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
980                        char *fnA, char *fnB,
981                        int rc_scaling, int ePBC,
982                        rvec com, rvec comB,
983                        warninp_t wi)
984 {
985     int i, j;
986
987     read_posres  (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
988     /* It is safer to simply read the b-state posres rather than trying
989      * to be smart and copy the positions.
990      */
991     read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
992 }
993
994 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
995                               t_inputrec *ir, warninp_t wi)
996 {
997     int  i;
998     char warn_buf[STRLEN];
999
1000     if (ir->nwall > 0)
1001     {
1002         fprintf(stderr, "Searching the wall atom type(s)\n");
1003     }
1004     for (i = 0; i < ir->nwall; i++)
1005     {
1006         ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
1007         if (ir->wall_atomtype[i] == NOTSET)
1008         {
1009             sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
1010             warning_error(wi, warn_buf);
1011         }
1012     }
1013 }
1014
1015 static int nrdf_internal(t_atoms *atoms)
1016 {
1017     int i, nmass, nrdf;
1018
1019     nmass = 0;
1020     for (i = 0; i < atoms->nr; i++)
1021     {
1022         /* Vsite ptype might not be set here yet, so also check the mass */
1023         if ((atoms->atom[i].ptype == eptAtom ||
1024              atoms->atom[i].ptype == eptNucleus)
1025             && atoms->atom[i].m > 0)
1026         {
1027             nmass++;
1028         }
1029     }
1030     switch (nmass)
1031     {
1032         case 0:  nrdf = 0; break;
1033         case 1:  nrdf = 0; break;
1034         case 2:  nrdf = 1; break;
1035         default: nrdf = nmass*3 - 6; break;
1036     }
1037
1038     return nrdf;
1039 }
1040
1041 void
1042 spline1d( double        dx,
1043           double *      y,
1044           int           n,
1045           double *      u,
1046           double *      y2 )
1047 {
1048     int    i;
1049     double p, q;
1050
1051     y2[0] = 0.0;
1052     u[0]  = 0.0;
1053
1054     for (i = 1; i < n-1; i++)
1055     {
1056         p     = 0.5*y2[i-1]+2.0;
1057         y2[i] = -0.5/p;
1058         q     = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1059         u[i]  = (3.0*q/dx-0.5*u[i-1])/p;
1060     }
1061
1062     y2[n-1] = 0.0;
1063
1064     for (i = n-2; i >= 0; i--)
1065     {
1066         y2[i] = y2[i]*y2[i+1]+u[i];
1067     }
1068 }
1069
1070
1071 void
1072 interpolate1d( double     xmin,
1073                double     dx,
1074                double *   ya,
1075                double *   y2a,
1076                double     x,
1077                double *   y,
1078                double *   y1)
1079 {
1080     int    ix;
1081     double a, b;
1082
1083     ix = (x-xmin)/dx;
1084
1085     a = (xmin+(ix+1)*dx-x)/dx;
1086     b = (x-xmin-ix*dx)/dx;
1087
1088     *y  = a*ya[ix]+b*ya[ix+1]+((a*a*a-a)*y2a[ix]+(b*b*b-b)*y2a[ix+1])*(dx*dx)/6.0;
1089     *y1 = (ya[ix+1]-ya[ix])/dx-(3.0*a*a-1.0)/6.0*dx*y2a[ix]+(3.0*b*b-1.0)/6.0*dx*y2a[ix+1];
1090 }
1091
1092
1093 void
1094 setup_cmap (int              grid_spacing,
1095             int              nc,
1096             real *           grid,
1097             gmx_cmap_t *     cmap_grid)
1098 {
1099     double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1100
1101     int     i, j, k, ii, jj, kk, idx;
1102     int     offset;
1103     double  dx, xmin, v, v1, v2, v12;
1104     double  phi, psi;
1105
1106     snew(tmp_u, 2*grid_spacing);
1107     snew(tmp_u2, 2*grid_spacing);
1108     snew(tmp_yy, 2*grid_spacing);
1109     snew(tmp_y1, 2*grid_spacing);
1110     snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1111     snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1112
1113     dx   = 360.0/grid_spacing;
1114     xmin = -180.0-dx*grid_spacing/2;
1115
1116     for (kk = 0; kk < nc; kk++)
1117     {
1118         /* Compute an offset depending on which cmap we are using
1119          * Offset will be the map number multiplied with the
1120          * grid_spacing * grid_spacing * 2
1121          */
1122         offset = kk * grid_spacing * grid_spacing * 2;
1123
1124         for (i = 0; i < 2*grid_spacing; i++)
1125         {
1126             ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1127
1128             for (j = 0; j < 2*grid_spacing; j++)
1129             {
1130                 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1131                 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1132             }
1133         }
1134
1135         for (i = 0; i < 2*grid_spacing; i++)
1136         {
1137             spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1138         }
1139
1140         for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1141         {
1142             ii  = i-grid_spacing/2;
1143             phi = ii*dx-180.0;
1144
1145             for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1146             {
1147                 jj  = j-grid_spacing/2;
1148                 psi = jj*dx-180.0;
1149
1150                 for (k = 0; k < 2*grid_spacing; k++)
1151                 {
1152                     interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1153                                   &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1154                 }
1155
1156                 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1157                 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1158                 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1159                 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1160
1161                 idx = ii*grid_spacing+jj;
1162                 cmap_grid->cmapdata[kk].cmap[idx*4]   = grid[offset+ii*grid_spacing+jj];
1163                 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1164                 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1165                 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1166             }
1167         }
1168     }
1169 }
1170
1171 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1172 {
1173     int i, k, nelem;
1174
1175     cmap_grid->ngrid        = ngrid;
1176     cmap_grid->grid_spacing = grid_spacing;
1177     nelem                   = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1178
1179     snew(cmap_grid->cmapdata, ngrid);
1180
1181     for (i = 0; i < cmap_grid->ngrid; i++)
1182     {
1183         snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1184     }
1185 }
1186
1187
1188 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1189 {
1190     int             count, count_mol, i, mb;
1191     gmx_molblock_t *molb;
1192     t_params       *plist;
1193     char            buf[STRLEN];
1194
1195     count = 0;
1196     for (mb = 0; mb < mtop->nmolblock; mb++)
1197     {
1198         count_mol = 0;
1199         molb      = &mtop->molblock[mb];
1200         plist     = mi[molb->type].plist;
1201
1202         for (i = 0; i < F_NRE; i++)
1203         {
1204             if (i == F_SETTLE)
1205             {
1206                 count_mol += 3*plist[i].nr;
1207             }
1208             else if (interaction_function[i].flags & IF_CONSTRAINT)
1209             {
1210                 count_mol += plist[i].nr;
1211             }
1212         }
1213
1214         if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1215         {
1216             sprintf(buf,
1217                     "Molecule type '%s' has %d constraints.\n"
1218                     "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1219                     *mi[molb->type].name, count_mol,
1220                     nrdf_internal(&mi[molb->type].atoms));
1221             warning(wi, buf);
1222         }
1223         count += molb->nmol*count_mol;
1224     }
1225
1226     return count;
1227 }
1228
1229 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1230 {
1231     int            i, nmiss, natoms, mt;
1232     real           q;
1233     const t_atoms *atoms;
1234
1235     nmiss = 0;
1236     for (mt = 0; mt < sys->nmoltype; mt++)
1237     {
1238         atoms  = &sys->moltype[mt].atoms;
1239         natoms = atoms->nr;
1240
1241         for (i = 0; i < natoms; i++)
1242         {
1243             q = atoms->atom[i].q;
1244             if ((get_atomtype_radius(atoms->atom[i].type, atype)    == 0  ||
1245                  get_atomtype_vol(atoms->atom[i].type, atype)       == 0  ||
1246                  get_atomtype_surftens(atoms->atom[i].type, atype)  == 0  ||
1247                  get_atomtype_gb_radius(atoms->atom[i].type, atype) == 0  ||
1248                  get_atomtype_S_hct(atoms->atom[i].type, atype)     == 0) &&
1249                 q != 0)
1250             {
1251                 fprintf(stderr, "\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1252                         get_atomtype_name(atoms->atom[i].type, atype), q);
1253                 nmiss++;
1254             }
1255         }
1256     }
1257
1258     if (nmiss > 0)
1259     {
1260         gmx_fatal(FARGS, "Can't do GB electrostatics; the implicit_genborn_params section of the forcefield has parameters with value zero for %d atomtypes that occur as charged atoms.", nmiss);
1261     }
1262 }
1263
1264
1265 static void check_gbsa_params(gpp_atomtype_t atype)
1266 {
1267     int  nmiss, i;
1268
1269     /* If we are doing GBSA, check that we got the parameters we need
1270      * This checking is to see if there are GBSA paratmeters for all
1271      * atoms in the force field. To go around this for testing purposes
1272      * comment out the nerror++ counter temporarily
1273      */
1274     nmiss = 0;
1275     for (i = 0; i < get_atomtype_ntypes(atype); i++)
1276     {
1277         if (get_atomtype_radius(i, atype)    < 0 ||
1278             get_atomtype_vol(i, atype)       < 0 ||
1279             get_atomtype_surftens(i, atype)  < 0 ||
1280             get_atomtype_gb_radius(i, atype) < 0 ||
1281             get_atomtype_S_hct(i, atype)     < 0)
1282         {
1283             fprintf(stderr, "\nGB parameter(s) missing or negative for atom type '%s'\n",
1284                     get_atomtype_name(i, atype));
1285             nmiss++;
1286         }
1287     }
1288
1289     if (nmiss > 0)
1290     {
1291         gmx_fatal(FARGS, "Can't do GB electrostatics; the implicit_genborn_params section of the forcefield is missing parameters for %d atomtypes or they might be negative.", nmiss);
1292     }
1293
1294 }
1295
1296 static real calc_temp(const gmx_mtop_t *mtop,
1297                       const t_inputrec *ir,
1298                       rvec             *v)
1299 {
1300     double                  sum_mv2;
1301     gmx_mtop_atomloop_all_t aloop;
1302     t_atom                 *atom;
1303     int                     a;
1304     int                     nrdf, g;
1305
1306     sum_mv2 = 0;
1307
1308     aloop = gmx_mtop_atomloop_all_init(mtop);
1309     while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
1310     {
1311         sum_mv2 += atom->m*norm2(v[a]);
1312     }
1313
1314     nrdf = 0;
1315     for (g = 0; g < ir->opts.ngtc; g++)
1316     {
1317         nrdf += ir->opts.nrdf[g];
1318     }
1319
1320     return sum_mv2/(nrdf*BOLTZ);
1321 }
1322
1323 static real get_max_reference_temp(const t_inputrec *ir,
1324                                    warninp_t         wi)
1325 {
1326     real     ref_t;
1327     int      i;
1328     gmx_bool bNoCoupl;
1329
1330     ref_t    = 0;
1331     bNoCoupl = FALSE;
1332     for (i = 0; i < ir->opts.ngtc; i++)
1333     {
1334         if (ir->opts.tau_t[i] < 0)
1335         {
1336             bNoCoupl = TRUE;
1337         }
1338         else
1339         {
1340             ref_t = max(ref_t, ir->opts.ref_t[i]);
1341         }
1342     }
1343
1344     if (bNoCoupl)
1345     {
1346         char buf[STRLEN];
1347
1348         sprintf(buf, "Some temperature coupling groups do not use temperature coupling. We will assume their temperature is not more than %.3f K. If their temperature is higher, the energy error and the Verlet buffer might be underestimated.",
1349                 ref_t);
1350         warning(wi, buf);
1351     }
1352
1353     return ref_t;
1354 }
1355
1356 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1357                               t_inputrec       *ir,
1358                               real              buffer_temp,
1359                               matrix            box,
1360                               warninp_t         wi)
1361 {
1362     int                    i;
1363     verletbuf_list_setup_t ls;
1364     real                   rlist_1x1;
1365     int                    n_nonlin_vsite;
1366     char                   warn_buf[STRLEN];
1367
1368     printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol, buffer_temp);
1369
1370     /* Calculate the buffer size for simple atom vs atoms list */
1371     ls.cluster_size_i = 1;
1372     ls.cluster_size_j = 1;
1373     calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1374                             &ls, &n_nonlin_vsite, &rlist_1x1);
1375
1376     /* Set the pair-list buffer size in ir */
1377     verletbuf_get_list_setup(FALSE, FALSE, &ls);
1378     calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1379                             &ls, &n_nonlin_vsite, &ir->rlist);
1380
1381     if (n_nonlin_vsite > 0)
1382     {
1383         sprintf(warn_buf, "There are %d non-linear virtual site constructions. Their contribution to the energy error is approximated. In most cases this does not affect the error significantly.", n_nonlin_vsite);
1384         warning_note(wi, warn_buf);
1385     }
1386
1387     printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1388            1, 1, rlist_1x1, rlist_1x1-max(ir->rvdw, ir->rcoulomb));
1389
1390     ir->rlistlong = ir->rlist;
1391     printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1392            ls.cluster_size_i, ls.cluster_size_j,
1393            ir->rlist, ir->rlist-max(ir->rvdw, ir->rcoulomb));
1394
1395     printf("Note that mdrun will redetermine rlist based on the actual pair-list setup\n");
1396
1397     if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
1398     {
1399         gmx_fatal(FARGS, "The pair-list cut-off (%g nm) is longer than half the shortest box vector or longer than the smallest box diagonal element (%g nm). Increase the box size or decrease nstlist or increase verlet-buffer-tolerance.", ir->rlistlong, sqrt(max_cutoff2(ir->ePBC, box)));
1400     }
1401 }
1402
1403 int gmx_grompp(int argc, char *argv[])
1404 {
1405     static const char *desc[] = {
1406         "[THISMODULE] (the gromacs preprocessor)",
1407         "reads a molecular topology file, checks the validity of the",
1408         "file, expands the topology from a molecular description to an atomic",
1409         "description. The topology file contains information about",
1410         "molecule types and the number of molecules, the preprocessor",
1411         "copies each molecule as needed. ",
1412         "There is no limitation on the number of molecule types. ",
1413         "Bonds and bond-angles can be converted into constraints, separately",
1414         "for hydrogens and heavy atoms.",
1415         "Then a coordinate file is read and velocities can be generated",
1416         "from a Maxwellian distribution if requested.",
1417         "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1418         "(eg. number of MD steps, time step, cut-off), and others such as",
1419         "NEMD parameters, which are corrected so that the net acceleration",
1420         "is zero.",
1421         "Eventually a binary file is produced that can serve as the sole input",
1422         "file for the MD program.[PAR]",
1423
1424         "[THISMODULE] uses the atom names from the topology file. The atom names",
1425         "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1426         "warnings when they do not match the atom names in the topology.",
1427         "Note that the atom names are irrelevant for the simulation as",
1428         "only the atom types are used for generating interaction parameters.[PAR]",
1429
1430         "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1431         "etc. The preprocessor supports the following keywords:[PAR]",
1432         "#ifdef VARIABLE[BR]",
1433         "#ifndef VARIABLE[BR]",
1434         "#else[BR]",
1435         "#endif[BR]",
1436         "#define VARIABLE[BR]",
1437         "#undef VARIABLE[BR]"
1438         "#include \"filename\"[BR]",
1439         "#include <filename>[PAR]",
1440         "The functioning of these statements in your topology may be modulated by",
1441         "using the following two flags in your [TT].mdp[tt] file:[PAR]",
1442         "[TT]define = -DVARIABLE1 -DVARIABLE2[BR]",
1443         "include = -I/home/john/doe[tt][BR]",
1444         "For further information a C-programming textbook may help you out.",
1445         "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1446         "topology file written out so that you can verify its contents.[PAR]",
1447
1448         /* cpp has been unnecessary for some time, hasn't it?
1449             "If your system does not have a C-preprocessor, you can still",
1450             "use [TT]grompp[tt], but you do not have access to the features ",
1451             "from the cpp. Command line options to the C-preprocessor can be given",
1452             "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1453          */
1454
1455         "When using position restraints a file with restraint coordinates",
1456         "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1457         "with respect to the conformation from the [TT]-c[tt] option.",
1458         "For free energy calculation the the coordinates for the B topology",
1459         "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1460         "those of the A topology.[PAR]",
1461
1462         "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1463         "The last frame with coordinates and velocities will be read,",
1464         "unless the [TT]-time[tt] option is used. Only if this information",
1465         "is absent will the coordinates in the [TT]-c[tt] file be used.",
1466         "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1467         "in your [TT].mdp[tt] file. An energy file can be supplied with",
1468         "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1469         "variables.[PAR]",
1470
1471         "[THISMODULE] can be used to restart simulations (preserving",
1472         "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1473         "However, for simply changing the number of run steps to extend",
1474         "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1475         "You then supply the old checkpoint file directly to [gmx-mdrun]",
1476         "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1477         "like output frequency, then supplying the checkpoint file to",
1478         "[THISMODULE] with [TT]-t[tt] along with a new [TT].mdp[tt] file",
1479         "with [TT]-f[tt] is the recommended procedure. Actually preserving",
1480         "the ensemble (if possible) still requires passing the checkpoint",
1481         "file to [gmx-mdrun] [TT]-cpi[tt].[PAR]",
1482
1483         "By default, all bonded interactions which have constant energy due to",
1484         "virtual site constructions will be removed. If this constant energy is",
1485         "not zero, this will result in a shift in the total energy. All bonded",
1486         "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1487         "all constraints for distances which will be constant anyway because",
1488         "of virtual site constructions will be removed. If any constraints remain",
1489         "which involve virtual sites, a fatal error will result.[PAR]"
1490
1491         "To verify your run input file, please take note of all warnings",
1492         "on the screen, and correct where necessary. Do also look at the contents",
1493         "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1494         "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1495         "with the [TT]-debug[tt] option which will give you more information",
1496         "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1497         "can see the contents of the run input file with the [gmx-dump]",
1498         "program. [gmx-check] can be used to compare the contents of two",
1499         "run input files.[PAR]"
1500
1501         "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1502         "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1503         "harmless, but usually they are not. The user is advised to carefully",
1504         "interpret the output messages before attempting to bypass them with",
1505         "this option."
1506     };
1507     t_gromppopts      *opts;
1508     gmx_mtop_t        *sys;
1509     int                nmi;
1510     t_molinfo         *mi;
1511     gpp_atomtype_t     atype;
1512     t_inputrec        *ir;
1513     int                natoms, nvsite, comb, mt;
1514     t_params          *plist;
1515     t_state           *state;
1516     matrix             box;
1517     real               max_spacing, fudgeQQ;
1518     double             reppow;
1519     char               fn[STRLEN], fnB[STRLEN];
1520     const char        *mdparin;
1521     int                ntype;
1522     gmx_bool           bNeedVel, bGenVel;
1523     gmx_bool           have_atomnumber;
1524     int                n12, n13, n14;
1525     t_params          *gb_plist = NULL;
1526     gmx_genborn_t     *born     = NULL;
1527     output_env_t       oenv;
1528     gmx_bool           bVerbose = FALSE;
1529     warninp_t          wi;
1530     char               warn_buf[STRLEN];
1531     unsigned int       useed;
1532     t_atoms            IMDatoms;   /* Atoms to be operated on interactively (IMD) */
1533
1534     t_filenm           fnm[] = {
1535         { efMDP, NULL,  NULL,        ffREAD  },
1536         { efMDP, "-po", "mdout",     ffWRITE },
1537         { efSTX, "-c",  NULL,        ffREAD  },
1538         { efSTX, "-r",  NULL,        ffOPTRD },
1539         { efSTX, "-rb", NULL,        ffOPTRD },
1540         { efNDX, NULL,  NULL,        ffOPTRD },
1541         { efTOP, NULL,  NULL,        ffREAD  },
1542         { efTOP, "-pp", "processed", ffOPTWR },
1543         { efTPX, "-o",  NULL,        ffWRITE },
1544         { efTRN, "-t",  NULL,        ffOPTRD },
1545         { efEDR, "-e",  NULL,        ffOPTRD },
1546         /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1547         { efGRO, "-imd", "imdgroup", ffOPTWR },
1548         { efTRN, "-ref", "rotref",   ffOPTRW }
1549     };
1550 #define NFILE asize(fnm)
1551
1552     /* Command line options */
1553     static gmx_bool bRenum   = TRUE;
1554     static gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1555     static int      i, maxwarn = 0;
1556     static real     fr_time = -1;
1557     t_pargs         pa[]    = {
1558         { "-v",       FALSE, etBOOL, {&bVerbose},
1559           "Be loud and noisy" },
1560         { "-time",    FALSE, etREAL, {&fr_time},
1561           "Take frame at or first after this time." },
1562         { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1563           "Remove constant bonded interactions with virtual sites" },
1564         { "-maxwarn", FALSE, etINT,  {&maxwarn},
1565           "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1566         { "-zero",    FALSE, etBOOL, {&bZero},
1567           "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1568         { "-renum",   FALSE, etBOOL, {&bRenum},
1569           "Renumber atomtypes and minimize number of atomtypes" }
1570     };
1571
1572     /* Initiate some variables */
1573     snew(ir, 1);
1574     snew(opts, 1);
1575     init_ir(ir, opts);
1576
1577     /* Parse the command line */
1578     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1579                            asize(desc), desc, 0, NULL, &oenv))
1580     {
1581         return 0;
1582     }
1583
1584     wi = init_warning(TRUE, maxwarn);
1585
1586     /* PARAMETER file processing */
1587     mdparin = opt2fn("-f", NFILE, fnm);
1588     set_warning_line(wi, mdparin, -1);
1589     get_ir(mdparin, opt2fn("-po", NFILE, fnm), ir, opts, wi);
1590
1591     if (bVerbose)
1592     {
1593         fprintf(stderr, "checking input for internal consistency...\n");
1594     }
1595     check_ir(mdparin, ir, opts, wi);
1596
1597     if (ir->ld_seed == -1)
1598     {
1599         ir->ld_seed = (gmx_int64_t)gmx_rng_make_seed();
1600         fprintf(stderr, "Setting the LD random seed to %"GMX_PRId64 "\n", ir->ld_seed);
1601     }
1602
1603     if (ir->expandedvals->lmc_seed == -1)
1604     {
1605         ir->expandedvals->lmc_seed = (int)gmx_rng_make_seed();
1606         fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1607     }
1608
1609     bNeedVel = EI_STATE_VELOCITY(ir->eI);
1610     bGenVel  = (bNeedVel && opts->bGenVel);
1611     if (bGenVel && ir->bContinuation)
1612     {
1613         sprintf(warn_buf,
1614                 "Generating velocities is inconsistent with attempting "
1615                 "to continue a previous run. Choose only one of "
1616                 "gen-vel = yes and continuation = yes.");
1617         warning_error(wi, warn_buf);
1618     }
1619
1620     snew(plist, F_NRE);
1621     init_plist(plist);
1622     snew(sys, 1);
1623     atype = init_atomtype();
1624     if (debug)
1625     {
1626         pr_symtab(debug, 0, "Just opened", &sys->symtab);
1627     }
1628
1629     strcpy(fn, ftp2fn(efTOP, NFILE, fnm));
1630     if (!gmx_fexist(fn))
1631     {
1632         gmx_fatal(FARGS, "%s does not exist", fn);
1633     }
1634     snew(state, 1);
1635     new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1636                opts, ir, bZero, bGenVel, bVerbose, state,
1637                atype, sys, &nmi, &mi, plist, &comb, &reppow, &fudgeQQ,
1638                opts->bMorse,
1639                wi);
1640
1641     if (debug)
1642     {
1643         pr_symtab(debug, 0, "After new_status", &sys->symtab);
1644     }
1645
1646     nvsite = 0;
1647     /* set parameters for virtual site construction (not for vsiten) */
1648     for (mt = 0; mt < sys->nmoltype; mt++)
1649     {
1650         nvsite +=
1651             set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1652     }
1653     /* now throw away all obsolete bonds, angles and dihedrals: */
1654     /* note: constraints are ALWAYS removed */
1655     if (nvsite)
1656     {
1657         for (mt = 0; mt < sys->nmoltype; mt++)
1658         {
1659             clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1660         }
1661     }
1662
1663     if (nvsite && ir->eI == eiNM)
1664     {
1665         gmx_fatal(FARGS, "Normal Mode analysis is not supported with virtual sites.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
1666     }
1667
1668     if (ir->cutoff_scheme == ecutsVERLET)
1669     {
1670         fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1671                 ecutscheme_names[ir->cutoff_scheme]);
1672
1673         /* Remove all charge groups */
1674         gmx_mtop_remove_chargegroups(sys);
1675     }
1676
1677     if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1678     {
1679         if (ir->eI == eiCG || ir->eI == eiLBFGS)
1680         {
1681             sprintf(warn_buf, "Can not do %s with %s, use %s",
1682                     EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1683             warning_error(wi, warn_buf);
1684         }
1685         if (ir->bPeriodicMols)
1686         {
1687             sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1688                     econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1689             warning_error(wi, warn_buf);
1690         }
1691     }
1692
1693     if (EI_SD (ir->eI) &&  ir->etc != etcNO)
1694     {
1695         warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1696     }
1697
1698     /* If we are doing QM/MM, check that we got the atom numbers */
1699     have_atomnumber = TRUE;
1700     for (i = 0; i < get_atomtype_ntypes(atype); i++)
1701     {
1702         have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1703     }
1704     if (!have_atomnumber && ir->bQMMM)
1705     {
1706         warning_error(wi,
1707                       "\n"
1708                       "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1709                       "field you are using does not contain atom numbers fields. This is an\n"
1710                       "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1711                       "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1712                       "an integer just before the mass column in ffXXXnb.itp.\n"
1713                       "NB: United atoms have the same atom numbers as normal ones.\n\n");
1714     }
1715
1716     if (ir->bAdress)
1717     {
1718         if ((ir->adress->const_wf > 1) || (ir->adress->const_wf < 0))
1719         {
1720             warning_error(wi, "AdResS contant weighting function should be between 0 and 1\n\n");
1721         }
1722         /** TODO check size of ex+hy width against box size */
1723     }
1724
1725     /* Check for errors in the input now, since they might cause problems
1726      * during processing further down.
1727      */
1728     check_warning_error(wi, FARGS);
1729
1730     if (opt2bSet("-r", NFILE, fnm))
1731     {
1732         sprintf(fn, "%s", opt2fn("-r", NFILE, fnm));
1733     }
1734     else
1735     {
1736         sprintf(fn, "%s", opt2fn("-c", NFILE, fnm));
1737     }
1738     if (opt2bSet("-rb", NFILE, fnm))
1739     {
1740         sprintf(fnB, "%s", opt2fn("-rb", NFILE, fnm));
1741     }
1742     else
1743     {
1744         strcpy(fnB, fn);
1745     }
1746
1747     if (nint_ftype(sys, mi, F_POSRES) > 0 || nint_ftype(sys, mi, F_FBPOSRES) > 0)
1748     {
1749         if (bVerbose)
1750         {
1751             fprintf(stderr, "Reading position restraint coords from %s", fn);
1752             if (strcmp(fn, fnB) == 0)
1753             {
1754                 fprintf(stderr, "\n");
1755             }
1756             else
1757             {
1758                 fprintf(stderr, " and %s\n", fnB);
1759             }
1760         }
1761         gen_posres(sys, mi, fn, fnB,
1762                    ir->refcoord_scaling, ir->ePBC,
1763                    ir->posres_com, ir->posres_comB,
1764                    wi);
1765     }
1766
1767     /* If we are using CMAP, setup the pre-interpolation grid */
1768     if (plist[F_CMAP].ncmap > 0)
1769     {
1770         init_cmap_grid(&sys->ffparams.cmap_grid, plist[F_CMAP].nc, plist[F_CMAP].grid_spacing);
1771         setup_cmap(plist[F_CMAP].grid_spacing, plist[F_CMAP].nc, plist[F_CMAP].cmap, &sys->ffparams.cmap_grid);
1772     }
1773
1774     set_wall_atomtype(atype, opts, ir, wi);
1775     if (bRenum)
1776     {
1777         renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1778         ntype = get_atomtype_ntypes(atype);
1779     }
1780
1781     if (ir->implicit_solvent != eisNO)
1782     {
1783         /* Now we have renumbered the atom types, we can check the GBSA params */
1784         check_gbsa_params(atype);
1785
1786         /* Check that all atoms that have charge and/or LJ-parameters also have
1787          * sensible GB-parameters
1788          */
1789         check_gbsa_params_charged(sys, atype);
1790     }
1791
1792     /* PELA: Copy the atomtype data to the topology atomtype list */
1793     copy_atomtype_atomtypes(atype, &(sys->atomtypes));
1794
1795     if (debug)
1796     {
1797         pr_symtab(debug, 0, "After renum_atype", &sys->symtab);
1798     }
1799
1800     if (bVerbose)
1801     {
1802         fprintf(stderr, "converting bonded parameters...\n");
1803     }
1804
1805     ntype = get_atomtype_ntypes(atype);
1806     convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1807
1808     if (debug)
1809     {
1810         pr_symtab(debug, 0, "After convert_params", &sys->symtab);
1811     }
1812
1813     /* set ptype to VSite for virtual sites */
1814     for (mt = 0; mt < sys->nmoltype; mt++)
1815     {
1816         set_vsites_ptype(FALSE, &sys->moltype[mt]);
1817     }
1818     if (debug)
1819     {
1820         pr_symtab(debug, 0, "After virtual sites", &sys->symtab);
1821     }
1822     /* Check velocity for virtual sites and shells */
1823     if (bGenVel)
1824     {
1825         check_vel(sys, state->v);
1826     }
1827
1828     /* check for shells and inpurecs */
1829     check_shells_inputrec(sys, ir, wi);
1830
1831     /* check masses */
1832     check_mol(sys, wi);
1833
1834     for (i = 0; i < sys->nmoltype; i++)
1835     {
1836         check_cg_sizes(ftp2fn(efTOP, NFILE, fnm), &sys->moltype[i].cgs, wi);
1837     }
1838
1839     if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1840     {
1841         check_bonds_timestep(sys, ir->delta_t, wi);
1842     }
1843
1844     if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
1845     {
1846         warning_note(wi, "Zero-step energy minimization will alter the coordinates before calculating the energy. If you just want the energy of a single point, try zero-step MD (with unconstrained_start = yes). To do multiple single-point energy evaluations of different configurations of the same topology, use mdrun -rerun.");
1847     }
1848
1849     check_warning_error(wi, FARGS);
1850
1851     if (bVerbose)
1852     {
1853         fprintf(stderr, "initialising group options...\n");
1854     }
1855     do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
1856              sys, bVerbose, ir,
1857              bGenVel ? state->v : NULL,
1858              wi);
1859
1860     if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0 &&
1861         ir->nstlist > 1)
1862     {
1863         if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
1864         {
1865             real buffer_temp;
1866
1867             if (EI_MD(ir->eI) && ir->etc == etcNO)
1868             {
1869                 if (bGenVel)
1870                 {
1871                     buffer_temp = opts->tempi;
1872                 }
1873                 else
1874                 {
1875                     buffer_temp = calc_temp(sys, ir, state->v);
1876                 }
1877                 if (buffer_temp > 0)
1878                 {
1879                     sprintf(warn_buf, "NVE simulation: will use the initial temperature of %.3f K for determining the Verlet buffer size", buffer_temp);
1880                     warning_note(wi, warn_buf);
1881                 }
1882                 else
1883                 {
1884                     sprintf(warn_buf, "NVE simulation with an initial temperature of zero: will use a Verlet buffer of %d%%. Check your energy drift!",
1885                             (int)(verlet_buffer_ratio_NVE_T0*100 + 0.5));
1886                     warning_note(wi, warn_buf);
1887                 }
1888             }
1889             else
1890             {
1891                 buffer_temp = get_max_reference_temp(ir, wi);
1892             }
1893
1894             if (EI_MD(ir->eI) && ir->etc == etcNO && buffer_temp == 0)
1895             {
1896                 /* NVE with initial T=0: we add a fixed ratio to rlist.
1897                  * Since we don't actually use verletbuf_tol, we set it to -1
1898                  * so it can't be misused later.
1899                  */
1900                 ir->rlist         *= 1.0 + verlet_buffer_ratio_NVE_T0;
1901                 ir->verletbuf_tol  = -1;
1902             }
1903             else
1904             {
1905                 /* We warn for NVE simulations with >1(.1)% drift tolerance */
1906                 const real drift_tol = 0.01;
1907                 real       ener_runtime;
1908
1909                 /* We use 2 DOF per atom = 2kT pot+kin energy, to be on
1910                  * the safe side with constraints (without constraints: 3 DOF).
1911                  */
1912                 ener_runtime = 2*BOLTZ*buffer_temp/(ir->nsteps*ir->delta_t);
1913
1914                 if (EI_MD(ir->eI) && ir->etc == etcNO && ir->nstlist > 1 &&
1915                     ir->nsteps > 0 &&
1916                     ir->verletbuf_tol > 1.1*drift_tol*ener_runtime)
1917                 {
1918                     sprintf(warn_buf, "You are using a Verlet buffer tolerance of %g kJ/mol/ps for an NVE simulation of length %g ps, which can give a final drift of %d%%. For conserving energy to %d%%, you might need to set verlet-buffer-tolerance to %.1e.",
1919                             ir->verletbuf_tol, ir->nsteps*ir->delta_t,
1920                             (int)(ir->verletbuf_tol/ener_runtime*100 + 0.5),
1921                             (int)(100*drift_tol + 0.5),
1922                             drift_tol*ener_runtime);
1923                     warning_note(wi, warn_buf);
1924                 }
1925
1926                 set_verlet_buffer(sys, ir, buffer_temp, state->box, wi);
1927             }
1928         }
1929     }
1930
1931     /* Init the temperature coupling state */
1932     init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
1933
1934     if (bVerbose)
1935     {
1936         fprintf(stderr, "Checking consistency between energy and charge groups...\n");
1937     }
1938     check_eg_vs_cg(sys);
1939
1940     if (debug)
1941     {
1942         pr_symtab(debug, 0, "After index", &sys->symtab);
1943     }
1944
1945     triple_check(mdparin, ir, sys, wi);
1946     close_symtab(&sys->symtab);
1947     if (debug)
1948     {
1949         pr_symtab(debug, 0, "After close", &sys->symtab);
1950     }
1951
1952     /* make exclusions between QM atoms */
1953     if (ir->bQMMM)
1954     {
1955         if (ir->QMMMscheme == eQMMMschemenormal && ir->ns_type == ensSIMPLE)
1956         {
1957             gmx_fatal(FARGS, "electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
1958         }
1959         else
1960         {
1961             generate_qmexcl(sys, ir, wi);
1962         }
1963     }
1964
1965     if (ftp2bSet(efTRN, NFILE, fnm))
1966     {
1967         if (bVerbose)
1968         {
1969             fprintf(stderr, "getting data from old trajectory ...\n");
1970         }
1971         cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
1972                     bNeedVel, bGenVel, fr_time, ir, state, sys, oenv);
1973     }
1974
1975     if (ir->ePBC == epbcXY && ir->nwall != 2)
1976     {
1977         clear_rvec(state->box[ZZ]);
1978     }
1979
1980     if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
1981     {
1982         set_warning_line(wi, mdparin, -1);
1983         check_chargegroup_radii(sys, ir, state->x, wi);
1984     }
1985
1986     if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
1987     {
1988         /* Calculate the optimal grid dimensions */
1989         copy_mat(state->box, box);
1990         if (ir->ePBC == epbcXY && ir->nwall == 2)
1991         {
1992             svmul(ir->wall_ewald_zfac, box[ZZ], box[ZZ]);
1993         }
1994         if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
1995         {
1996             /* Mark fourier_spacing as not used */
1997             ir->fourier_spacing = 0;
1998         }
1999         else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
2000         {
2001             set_warning_line(wi, mdparin, -1);
2002             warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
2003         }
2004         max_spacing = calc_grid(stdout, box, ir->fourier_spacing,
2005                                 &(ir->nkx), &(ir->nky), &(ir->nkz));
2006     }
2007
2008     /* MRS: eventually figure out better logic for initializing the fep
2009        values that makes declaring the lambda and declaring the state not
2010        potentially conflict if not handled correctly. */
2011     if (ir->efep != efepNO)
2012     {
2013         state->fep_state = ir->fepvals->init_fep_state;
2014         for (i = 0; i < efptNR; i++)
2015         {
2016             /* init_lambda trumps state definitions*/
2017             if (ir->fepvals->init_lambda >= 0)
2018             {
2019                 state->lambda[i] = ir->fepvals->init_lambda;
2020             }
2021             else
2022             {
2023                 if (ir->fepvals->all_lambda[i] == NULL)
2024                 {
2025                     gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
2026                 }
2027                 else
2028                 {
2029                     state->lambda[i] = ir->fepvals->all_lambda[i][state->fep_state];
2030                 }
2031             }
2032         }
2033     }
2034
2035     if (ir->ePull != epullNO)
2036     {
2037         set_pull_init(ir, sys, state->x, state->box, state->lambda[efptMASS], oenv, opts->pull_start);
2038     }
2039
2040     if (ir->bRot)
2041     {
2042         set_reference_positions(ir->rot, state->x, state->box,
2043                                 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
2044                                 wi);
2045     }
2046
2047     /*  reset_multinr(sys); */
2048
2049     if (EEL_PME(ir->coulombtype))
2050     {
2051         float ratio = pme_load_estimate(sys, ir, state->box);
2052         fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
2053         /* With free energy we might need to do PME both for the A and B state
2054          * charges. This will double the cost, but the optimal performance will
2055          * then probably be at a slightly larger cut-off and grid spacing.
2056          */
2057         if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
2058             (ir->efep != efepNO && ratio > 2.0/3.0))
2059         {
2060             warning_note(wi,
2061                          "The optimal PME mesh load for parallel simulations is below 0.5\n"
2062                          "and for highly parallel simulations between 0.25 and 0.33,\n"
2063                          "for higher performance, increase the cut-off and the PME grid spacing.\n");
2064             if (ir->efep != efepNO)
2065             {
2066                 warning_note(wi,
2067                              "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
2068             }
2069         }
2070     }
2071
2072     {
2073         char   warn_buf[STRLEN];
2074         double cio = compute_io(ir, sys->natoms, &sys->groups, F_NRE, 1);
2075         sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
2076         if (cio > 2000)
2077         {
2078             set_warning_line(wi, mdparin, -1);
2079             warning_note(wi, warn_buf);
2080         }
2081         else
2082         {
2083             printf("%s\n", warn_buf);
2084         }
2085     }
2086
2087     if (bVerbose)
2088     {
2089         fprintf(stderr, "writing run input file...\n");
2090     }
2091
2092     done_warning(wi, FARGS);
2093     write_tpx_state(ftp2fn(efTPX, NFILE, fnm), ir, state, sys);
2094
2095     /* Output IMD group, if bIMD is TRUE */
2096     write_IMDgroup_to_file(ir->bIMD, ir, state, sys, NFILE, fnm);
2097
2098     done_state(state);
2099     sfree(state);
2100     done_atomtype(atype);
2101     done_mtop(sys, TRUE);
2102     done_inputrec_strings();
2103
2104     return 0;
2105 }