Move symtab.* to topology/
[gromacs.git] / src / gromacs / gmxpreprocess / grompp.c
blob88a24a4243bbfdf444de8046d8c9d69e25e4172f
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "grompp.h"
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
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>
50 #include "macros.h"
51 #include "readir.h"
52 #include "toputil.h"
53 #include "topio.h"
54 #include "gromacs/fileio/confio.h"
55 #include "readir.h"
56 #include "names.h"
57 #include "grompp-impl.h"
58 #include "gromacs/random/random.h"
59 #include "gromacs/gmxpreprocess/gen_maxwell_velocities.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/commandline/pargs.h"
63 #include "splitter.h"
64 #include "gromacs/gmxpreprocess/sortwater.h"
65 #include "convparm.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "warninp.h"
68 #include "index.h"
69 #include "gromacs/fileio/gmxfio.h"
70 #include "gromacs/fileio/trnio.h"
71 #include "gromacs/fileio/tpxio.h"
72 #include "gromacs/fileio/trxio.h"
73 #include "vsite_parm.h"
74 #include "txtdump.h"
75 #include "calcgrid.h"
76 #include "add_par.h"
77 #include "gromacs/fileio/enxio.h"
78 #include "perf_est.h"
79 #include "compute_io.h"
80 #include "gpp_atomtype.h"
81 #include "mtop_util.h"
82 #include "genborn.h"
83 #include "calc_verletbuf.h"
84 #include "tomorse.h"
85 #include "gromacs/imd/imd.h"
87 #include "gromacs/pbcutil/pbc.h"
88 #include "gromacs/topology/symtab.h"
89 #include "gromacs/utility/smalloc.h"
91 static int rm_interactions(int ifunc, int nrmols, t_molinfo mols[])
93 int i, n;
95 n = 0;
96 /* For all the molecule types */
97 for (i = 0; i < nrmols; i++)
99 n += mols[i].plist[ifunc].nr;
100 mols[i].plist[ifunc].nr = 0;
102 return n;
105 static int check_atom_names(const char *fn1, const char *fn2,
106 gmx_mtop_t *mtop, t_atoms *at)
108 int mb, m, i, j, nmismatch;
109 t_atoms *tat;
110 #define MAXMISMATCH 20
112 if (mtop->natoms != at->nr)
114 gmx_incons("comparing atom names");
117 nmismatch = 0;
118 i = 0;
119 for (mb = 0; mb < mtop->nmolblock; mb++)
121 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
122 for (m = 0; m < mtop->molblock[mb].nmol; m++)
124 for (j = 0; j < tat->nr; j++)
126 if (strcmp( *(tat->atomname[j]), *(at->atomname[i]) ) != 0)
128 if (nmismatch < MAXMISMATCH)
130 fprintf(stderr,
131 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
132 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
134 else if (nmismatch == MAXMISMATCH)
136 fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
138 nmismatch++;
140 i++;
145 return nmismatch;
148 static void check_eg_vs_cg(gmx_mtop_t *mtop)
150 int astart, mb, m, cg, j, firstj;
151 unsigned char firsteg, eg;
152 gmx_moltype_t *molt;
154 /* Go through all the charge groups and make sure all their
155 * atoms are in the same energy group.
158 astart = 0;
159 for (mb = 0; mb < mtop->nmolblock; mb++)
161 molt = &mtop->moltype[mtop->molblock[mb].type];
162 for (m = 0; m < mtop->molblock[mb].nmol; m++)
164 for (cg = 0; cg < molt->cgs.nr; cg++)
166 /* Get the energy group of the first atom in this charge group */
167 firstj = astart + molt->cgs.index[cg];
168 firsteg = ggrpnr(&mtop->groups, egcENER, firstj);
169 for (j = molt->cgs.index[cg]+1; j < molt->cgs.index[cg+1]; j++)
171 eg = ggrpnr(&mtop->groups, egcENER, astart+j);
172 if (eg != firsteg)
174 gmx_fatal(FARGS, "atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
175 firstj+1, astart+j+1, cg+1, *molt->name);
179 astart += molt->atoms.nr;
184 static void check_cg_sizes(const char *topfn, t_block *cgs, warninp_t wi)
186 int maxsize, cg;
187 char warn_buf[STRLEN];
189 maxsize = 0;
190 for (cg = 0; cg < cgs->nr; cg++)
192 maxsize = max(maxsize, cgs->index[cg+1]-cgs->index[cg]);
195 if (maxsize > MAX_CHARGEGROUP_SIZE)
197 gmx_fatal(FARGS, "The largest charge group contains %d atoms. The maximum is %d.", maxsize, MAX_CHARGEGROUP_SIZE);
199 else if (maxsize > 10)
201 set_warning_line(wi, topfn, -1);
202 sprintf(warn_buf,
203 "The largest charge group contains %d atoms.\n"
204 "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"
205 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
206 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
207 maxsize);
208 warning_note(wi, warn_buf);
212 static void check_bonds_timestep(gmx_mtop_t *mtop, double dt, warninp_t wi)
214 /* This check is not intended to ensure accurate integration,
215 * rather it is to signal mistakes in the mdp settings.
216 * A common mistake is to forget to turn on constraints
217 * for MD after energy minimization with flexible bonds.
218 * This check can also detect too large time steps for flexible water
219 * models, but such errors will often be masked by the constraints
220 * mdp options, which turns flexible water into water with bond constraints,
221 * but without an angle constraint. Unfortunately such incorrect use
222 * of water models can not easily be detected without checking
223 * for specific model names.
225 * The stability limit of leap-frog or velocity verlet is 4.44 steps
226 * per oscillational period.
227 * But accurate bonds distributions are lost far before that limit.
228 * To allow relatively common schemes (although not common with Gromacs)
229 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
230 * we set the note limit to 10.
232 int min_steps_warn = 5;
233 int min_steps_note = 10;
234 t_iparams *ip;
235 int molt;
236 gmx_moltype_t *moltype, *w_moltype;
237 t_atom *atom;
238 t_ilist *ilist, *ilb, *ilc, *ils;
239 int ftype;
240 int i, a1, a2, w_a1, w_a2, j;
241 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
242 gmx_bool bFound, bWater, bWarn;
243 char warn_buf[STRLEN];
245 ip = mtop->ffparams.iparams;
247 twopi2 = sqr(2*M_PI);
249 limit2 = sqr(min_steps_note*dt);
251 w_a1 = w_a2 = -1;
252 w_period2 = -1.0;
254 w_moltype = NULL;
255 for (molt = 0; molt < mtop->nmoltype; molt++)
257 moltype = &mtop->moltype[molt];
258 atom = moltype->atoms.atom;
259 ilist = moltype->ilist;
260 ilc = &ilist[F_CONSTR];
261 ils = &ilist[F_SETTLE];
262 for (ftype = 0; ftype < F_NRE; ftype++)
264 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
266 continue;
269 ilb = &ilist[ftype];
270 for (i = 0; i < ilb->nr; i += 3)
272 fc = ip[ilb->iatoms[i]].harmonic.krA;
273 re = ip[ilb->iatoms[i]].harmonic.rA;
274 if (ftype == F_G96BONDS)
276 /* Convert squared sqaure fc to harmonic fc */
277 fc = 2*fc*re;
279 a1 = ilb->iatoms[i+1];
280 a2 = ilb->iatoms[i+2];
281 m1 = atom[a1].m;
282 m2 = atom[a2].m;
283 if (fc > 0 && m1 > 0 && m2 > 0)
285 period2 = twopi2*m1*m2/((m1 + m2)*fc);
287 else
289 period2 = GMX_FLOAT_MAX;
291 if (debug)
293 fprintf(debug, "fc %g m1 %g m2 %g period %g\n",
294 fc, m1, m2, sqrt(period2));
296 if (period2 < limit2)
298 bFound = FALSE;
299 for (j = 0; j < ilc->nr; j += 3)
301 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
302 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
304 bFound = TRUE;
307 for (j = 0; j < ils->nr; j += 4)
309 if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
310 (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
312 bFound = TRUE;
315 if (!bFound &&
316 (w_moltype == NULL || period2 < w_period2))
318 w_moltype = moltype;
319 w_a1 = a1;
320 w_a2 = a2;
321 w_period2 = period2;
328 if (w_moltype != NULL)
330 bWarn = (w_period2 < sqr(min_steps_warn*dt));
331 /* A check that would recognize most water models */
332 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
333 w_moltype->atoms.nr <= 5);
334 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"
335 "%s",
336 *w_moltype->name,
337 w_a1+1, *w_moltype->atoms.atomname[w_a1],
338 w_a2+1, *w_moltype->atoms.atomname[w_a2],
339 sqrt(w_period2), bWarn ? min_steps_warn : min_steps_note, dt,
340 bWater ?
341 "Maybe you asked for fexible water." :
342 "Maybe you forgot to change the constraints mdp option.");
343 if (bWarn)
345 warning(wi, warn_buf);
347 else
349 warning_note(wi, warn_buf);
354 static void check_vel(gmx_mtop_t *mtop, rvec v[])
356 gmx_mtop_atomloop_all_t aloop;
357 t_atom *atom;
358 int a;
360 aloop = gmx_mtop_atomloop_all_init(mtop);
361 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
363 if (atom->ptype == eptShell ||
364 atom->ptype == eptBond ||
365 atom->ptype == eptVSite)
367 clear_rvec(v[a]);
372 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
374 int nint, mb;
376 nint = 0;
377 for (mb = 0; mb < mtop->nmolblock; mb++)
379 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
382 return nint;
385 /* This routine reorders the molecule type array
386 * in the order of use in the molblocks,
387 * unused molecule types are deleted.
389 static void renumber_moltypes(gmx_mtop_t *sys,
390 int *nmolinfo, t_molinfo **molinfo)
392 int *order, norder, i;
393 int mb, mi;
394 t_molinfo *minew;
396 snew(order, *nmolinfo);
397 norder = 0;
398 for (mb = 0; mb < sys->nmolblock; mb++)
400 for (i = 0; i < norder; i++)
402 if (order[i] == sys->molblock[mb].type)
404 break;
407 if (i == norder)
409 /* This type did not occur yet, add it */
410 order[norder] = sys->molblock[mb].type;
411 /* Renumber the moltype in the topology */
412 norder++;
414 sys->molblock[mb].type = i;
417 /* We still need to reorder the molinfo structs */
418 snew(minew, norder);
419 for (mi = 0; mi < *nmolinfo; mi++)
421 for (i = 0; i < norder; i++)
423 if (order[i] == mi)
425 break;
428 if (i == norder)
430 done_mi(&(*molinfo)[mi]);
432 else
434 minew[i] = (*molinfo)[mi];
437 sfree(*molinfo);
439 *nmolinfo = norder;
440 *molinfo = minew;
443 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
445 int m;
446 gmx_moltype_t *molt;
448 mtop->nmoltype = nmi;
449 snew(mtop->moltype, nmi);
450 for (m = 0; m < nmi; m++)
452 molt = &mtop->moltype[m];
453 molt->name = mi[m].name;
454 molt->atoms = mi[m].atoms;
455 /* ilists are copied later */
456 molt->cgs = mi[m].cgs;
457 molt->excls = mi[m].excls;
461 static void
462 new_status(const char *topfile, const char *topppfile, const char *confin,
463 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
464 gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
465 gpp_atomtype_t atype, gmx_mtop_t *sys,
466 int *nmi, t_molinfo **mi, t_params plist[],
467 int *comb, double *reppow, real *fudgeQQ,
468 gmx_bool bMorse,
469 warninp_t wi)
471 t_molinfo *molinfo = NULL;
472 int nmolblock;
473 gmx_molblock_t *molblock, *molbs;
474 t_atoms *confat;
475 int mb, i, nrmols, nmismatch;
476 char buf[STRLEN];
477 gmx_bool bGB = FALSE;
478 char warn_buf[STRLEN];
480 init_mtop(sys);
482 /* Set gmx_boolean for GB */
483 if (ir->implicit_solvent)
485 bGB = TRUE;
488 /* TOPOLOGY processing */
489 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
490 plist, comb, reppow, fudgeQQ,
491 atype, &nrmols, &molinfo, ir,
492 &nmolblock, &molblock, bGB,
493 wi);
495 sys->nmolblock = 0;
496 snew(sys->molblock, nmolblock);
498 sys->natoms = 0;
499 for (mb = 0; mb < nmolblock; mb++)
501 if (sys->nmolblock > 0 &&
502 molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
504 /* Merge consecutive blocks with the same molecule type */
505 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
506 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
508 else if (molblock[mb].nmol > 0)
510 /* Add a new molblock to the topology */
511 molbs = &sys->molblock[sys->nmolblock];
512 *molbs = molblock[mb];
513 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
514 molbs->nposres_xA = 0;
515 molbs->nposres_xB = 0;
516 sys->natoms += molbs->nmol*molbs->natoms_mol;
517 sys->nmolblock++;
520 if (sys->nmolblock == 0)
522 gmx_fatal(FARGS, "No molecules were defined in the system");
525 renumber_moltypes(sys, &nrmols, &molinfo);
527 if (bMorse)
529 convert_harmonics(nrmols, molinfo, atype);
532 if (ir->eDisre == edrNone)
534 i = rm_interactions(F_DISRES, nrmols, molinfo);
535 if (i > 0)
537 set_warning_line(wi, "unknown", -1);
538 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
539 warning_note(wi, warn_buf);
542 if (opts->bOrire == FALSE)
544 i = rm_interactions(F_ORIRES, nrmols, molinfo);
545 if (i > 0)
547 set_warning_line(wi, "unknown", -1);
548 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
549 warning_note(wi, warn_buf);
553 /* Copy structures from msys to sys */
554 molinfo2mtop(nrmols, molinfo, sys);
556 gmx_mtop_finalize(sys);
558 /* COORDINATE file processing */
559 if (bVerbose)
561 fprintf(stderr, "processing coordinates...\n");
564 get_stx_coordnum(confin, &state->natoms);
565 if (state->natoms != sys->natoms)
567 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
568 " does not match topology (%s, %d)",
569 confin, state->natoms, topfile, sys->natoms);
571 else
573 /* make space for coordinates and velocities */
574 char title[STRLEN];
575 snew(confat, 1);
576 init_t_atoms(confat, state->natoms, FALSE);
577 init_state(state, state->natoms, 0, 0, 0, 0);
578 read_stx_conf(confin, title, confat, state->x, state->v, NULL, state->box);
579 /* This call fixes the box shape for runs with pressure scaling */
580 set_box_rel(ir, state);
582 nmismatch = check_atom_names(topfile, confin, sys, confat);
583 free_t_atoms(confat, TRUE);
584 sfree(confat);
586 if (nmismatch)
588 sprintf(buf, "%d non-matching atom name%s\n"
589 "atom names from %s will be used\n"
590 "atom names from %s will be ignored\n",
591 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
592 warning(wi, buf);
594 if (bVerbose)
596 fprintf(stderr, "double-checking input for internal consistency...\n");
598 double_check(ir, state->box, nint_ftype(sys, molinfo, F_CONSTR), wi);
601 if (bGenVel)
603 real *mass;
604 gmx_mtop_atomloop_all_t aloop;
605 t_atom *atom;
606 unsigned int useed;
608 snew(mass, state->natoms);
609 aloop = gmx_mtop_atomloop_all_init(sys);
610 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
612 mass[i] = atom->m;
615 useed = opts->seed;
616 if (opts->seed == -1)
618 useed = (int)gmx_rng_make_seed();
619 fprintf(stderr, "Setting gen_seed to %u\n", useed);
621 maxwell_speed(opts->tempi, useed, sys, state->v);
623 stop_cm(stdout, state->natoms, mass, state->x, state->v);
624 sfree(mass);
627 *nmi = nrmols;
628 *mi = molinfo;
631 static void copy_state(const char *slog, t_trxframe *fr,
632 gmx_bool bReadVel, t_state *state,
633 double *use_time)
635 int i;
637 if (fr->not_ok & FRAME_NOT_OK)
639 gmx_fatal(FARGS, "Can not start from an incomplete frame");
641 if (!fr->bX)
643 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
644 slog);
647 for (i = 0; i < state->natoms; i++)
649 copy_rvec(fr->x[i], state->x[i]);
651 if (bReadVel)
653 if (!fr->bV)
655 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
657 for (i = 0; i < state->natoms; i++)
659 copy_rvec(fr->v[i], state->v[i]);
662 if (fr->bBox)
664 copy_mat(fr->box, state->box);
667 *use_time = fr->time;
670 static void cont_status(const char *slog, const char *ener,
671 gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
672 t_inputrec *ir, t_state *state,
673 gmx_mtop_t *sys,
674 const output_env_t oenv)
675 /* If fr_time == -1 read the last frame available which is complete */
677 gmx_bool bReadVel;
678 t_trxframe fr;
679 t_trxstatus *fp;
680 int i;
681 double use_time;
683 bReadVel = (bNeedVel && !bGenVel);
685 fprintf(stderr,
686 "Reading Coordinates%s and Box size from old trajectory\n",
687 bReadVel ? ", Velocities" : "");
688 if (fr_time == -1)
690 fprintf(stderr, "Will read whole trajectory\n");
692 else
694 fprintf(stderr, "Will read till time %g\n", fr_time);
696 if (!bReadVel)
698 if (bGenVel)
700 fprintf(stderr, "Velocities generated: "
701 "ignoring velocities in input trajectory\n");
703 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
705 else
707 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
709 if (!fr.bV)
711 fprintf(stderr,
712 "\n"
713 "WARNING: Did not find a frame with velocities in file %s,\n"
714 " all velocities will be set to zero!\n\n", slog);
715 for (i = 0; i < sys->natoms; i++)
717 clear_rvec(state->v[i]);
719 close_trj(fp);
720 /* Search for a frame without velocities */
721 bReadVel = FALSE;
722 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
726 state->natoms = fr.natoms;
728 if (sys->natoms != state->natoms)
730 gmx_fatal(FARGS, "Number of atoms in Topology "
731 "is not the same as in Trajectory");
733 copy_state(slog, &fr, bReadVel, state, &use_time);
735 /* Find the appropriate frame */
736 while ((fr_time == -1 || fr.time < fr_time) &&
737 read_next_frame(oenv, fp, &fr))
739 copy_state(slog, &fr, bReadVel, state, &use_time);
742 close_trj(fp);
744 /* Set the relative box lengths for preserving the box shape.
745 * Note that this call can lead to differences in the last bit
746 * with respect to using gmx convert-tpr to create a [TT].tpx[tt] file.
748 set_box_rel(ir, state);
750 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
751 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
753 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
755 get_enx_state(ener, use_time, &sys->groups, ir, state);
756 preserve_box_shape(ir, state->box_rel, state->boxv);
760 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
761 char *fn,
762 int rc_scaling, int ePBC,
763 rvec com,
764 warninp_t wi)
766 gmx_bool bFirst = TRUE, *hadAtom;
767 rvec *x, *v, *xp;
768 dvec sum;
769 double totmass;
770 t_atoms dumat;
771 matrix box, invbox;
772 int natoms, npbcdim = 0;
773 char warn_buf[STRLEN], title[STRLEN];
774 int a, i, ai, j, k, mb, nat_molb;
775 gmx_molblock_t *molb;
776 t_params *pr, *prfb;
777 t_atom *atom;
779 get_stx_coordnum(fn, &natoms);
780 if (natoms != mtop->natoms)
782 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);
783 warning(wi, warn_buf);
785 snew(x, natoms);
786 snew(v, natoms);
787 init_t_atoms(&dumat, natoms, FALSE);
788 read_stx_conf(fn, title, &dumat, x, v, NULL, box);
790 npbcdim = ePBC2npbcdim(ePBC);
791 clear_rvec(com);
792 if (rc_scaling != erscNO)
794 copy_mat(box, invbox);
795 for (j = npbcdim; j < DIM; j++)
797 clear_rvec(invbox[j]);
798 invbox[j][j] = 1;
800 m_inv_ur0(invbox, invbox);
803 /* Copy the reference coordinates to mtop */
804 clear_dvec(sum);
805 totmass = 0;
806 a = 0;
807 snew(hadAtom, natoms);
808 for (mb = 0; mb < mtop->nmolblock; mb++)
810 molb = &mtop->molblock[mb];
811 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
812 pr = &(molinfo[molb->type].plist[F_POSRES]);
813 prfb = &(molinfo[molb->type].plist[F_FBPOSRES]);
814 if (pr->nr > 0 || prfb->nr > 0)
816 atom = mtop->moltype[molb->type].atoms.atom;
817 for (i = 0; (i < pr->nr); i++)
819 ai = pr->param[i].AI;
820 if (ai >= natoms)
822 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
823 ai+1, *molinfo[molb->type].name, fn, natoms);
825 hadAtom[ai] = TRUE;
826 if (rc_scaling == erscCOM)
828 /* Determine the center of mass of the posres reference coordinates */
829 for (j = 0; j < npbcdim; j++)
831 sum[j] += atom[ai].m*x[a+ai][j];
833 totmass += atom[ai].m;
836 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
837 for (i = 0; (i < prfb->nr); i++)
839 ai = prfb->param[i].AI;
840 if (ai >= natoms)
842 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
843 ai+1, *molinfo[molb->type].name, fn, natoms);
845 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE)
847 /* Determine the center of mass of the posres reference coordinates */
848 for (j = 0; j < npbcdim; j++)
850 sum[j] += atom[ai].m*x[a+ai][j];
852 totmass += atom[ai].m;
855 if (!bTopB)
857 molb->nposres_xA = nat_molb;
858 snew(molb->posres_xA, molb->nposres_xA);
859 for (i = 0; i < nat_molb; i++)
861 copy_rvec(x[a+i], molb->posres_xA[i]);
864 else
866 molb->nposres_xB = nat_molb;
867 snew(molb->posres_xB, molb->nposres_xB);
868 for (i = 0; i < nat_molb; i++)
870 copy_rvec(x[a+i], molb->posres_xB[i]);
874 a += nat_molb;
876 if (rc_scaling == erscCOM)
878 if (totmass == 0)
880 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
882 for (j = 0; j < npbcdim; j++)
884 com[j] = sum[j]/totmass;
886 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]);
889 if (rc_scaling != erscNO)
891 assert(npbcdim <= DIM);
893 for (mb = 0; mb < mtop->nmolblock; mb++)
895 molb = &mtop->molblock[mb];
896 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
897 if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
899 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
900 for (i = 0; i < nat_molb; i++)
902 for (j = 0; j < npbcdim; j++)
904 if (rc_scaling == erscALL)
906 /* Convert from Cartesian to crystal coordinates */
907 xp[i][j] *= invbox[j][j];
908 for (k = j+1; k < npbcdim; k++)
910 xp[i][j] += invbox[k][j]*xp[i][k];
913 else if (rc_scaling == erscCOM)
915 /* Subtract the center of mass */
916 xp[i][j] -= com[j];
923 if (rc_scaling == erscCOM)
925 /* Convert the COM from Cartesian to crystal coordinates */
926 for (j = 0; j < npbcdim; j++)
928 com[j] *= invbox[j][j];
929 for (k = j+1; k < npbcdim; k++)
931 com[j] += invbox[k][j]*com[k];
937 free_t_atoms(&dumat, TRUE);
938 sfree(x);
939 sfree(v);
940 sfree(hadAtom);
943 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
944 char *fnA, char *fnB,
945 int rc_scaling, int ePBC,
946 rvec com, rvec comB,
947 warninp_t wi)
949 int i, j;
951 read_posres (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
952 if (strcmp(fnA, fnB) != 0)
954 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
958 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
959 t_inputrec *ir, warninp_t wi)
961 int i;
962 char warn_buf[STRLEN];
964 if (ir->nwall > 0)
966 fprintf(stderr, "Searching the wall atom type(s)\n");
968 for (i = 0; i < ir->nwall; i++)
970 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
971 if (ir->wall_atomtype[i] == NOTSET)
973 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
974 warning_error(wi, warn_buf);
979 static int nrdf_internal(t_atoms *atoms)
981 int i, nmass, nrdf;
983 nmass = 0;
984 for (i = 0; i < atoms->nr; i++)
986 /* Vsite ptype might not be set here yet, so also check the mass */
987 if ((atoms->atom[i].ptype == eptAtom ||
988 atoms->atom[i].ptype == eptNucleus)
989 && atoms->atom[i].m > 0)
991 nmass++;
994 switch (nmass)
996 case 0: nrdf = 0; break;
997 case 1: nrdf = 0; break;
998 case 2: nrdf = 1; break;
999 default: nrdf = nmass*3 - 6; break;
1002 return nrdf;
1005 void
1006 spline1d( double dx,
1007 double * y,
1008 int n,
1009 double * u,
1010 double * y2 )
1012 int i;
1013 double p, q;
1015 y2[0] = 0.0;
1016 u[0] = 0.0;
1018 for (i = 1; i < n-1; i++)
1020 p = 0.5*y2[i-1]+2.0;
1021 y2[i] = -0.5/p;
1022 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1023 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1026 y2[n-1] = 0.0;
1028 for (i = n-2; i >= 0; i--)
1030 y2[i] = y2[i]*y2[i+1]+u[i];
1035 void
1036 interpolate1d( double xmin,
1037 double dx,
1038 double * ya,
1039 double * y2a,
1040 double x,
1041 double * y,
1042 double * y1)
1044 int ix;
1045 double a, b;
1047 ix = (x-xmin)/dx;
1049 a = (xmin+(ix+1)*dx-x)/dx;
1050 b = (x-xmin-ix*dx)/dx;
1052 *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;
1053 *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];
1057 void
1058 setup_cmap (int grid_spacing,
1059 int nc,
1060 real * grid,
1061 gmx_cmap_t * cmap_grid)
1063 double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1065 int i, j, k, ii, jj, kk, idx;
1066 int offset;
1067 double dx, xmin, v, v1, v2, v12;
1068 double phi, psi;
1070 snew(tmp_u, 2*grid_spacing);
1071 snew(tmp_u2, 2*grid_spacing);
1072 snew(tmp_yy, 2*grid_spacing);
1073 snew(tmp_y1, 2*grid_spacing);
1074 snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1075 snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1077 dx = 360.0/grid_spacing;
1078 xmin = -180.0-dx*grid_spacing/2;
1080 for (kk = 0; kk < nc; kk++)
1082 /* Compute an offset depending on which cmap we are using
1083 * Offset will be the map number multiplied with the
1084 * grid_spacing * grid_spacing * 2
1086 offset = kk * grid_spacing * grid_spacing * 2;
1088 for (i = 0; i < 2*grid_spacing; i++)
1090 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1092 for (j = 0; j < 2*grid_spacing; j++)
1094 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1095 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1099 for (i = 0; i < 2*grid_spacing; i++)
1101 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1104 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1106 ii = i-grid_spacing/2;
1107 phi = ii*dx-180.0;
1109 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1111 jj = j-grid_spacing/2;
1112 psi = jj*dx-180.0;
1114 for (k = 0; k < 2*grid_spacing; k++)
1116 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1117 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1120 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1121 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1122 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1123 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1125 idx = ii*grid_spacing+jj;
1126 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1127 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1128 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1129 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1135 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1137 int i, k, nelem;
1139 cmap_grid->ngrid = ngrid;
1140 cmap_grid->grid_spacing = grid_spacing;
1141 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1143 snew(cmap_grid->cmapdata, ngrid);
1145 for (i = 0; i < cmap_grid->ngrid; i++)
1147 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1152 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1154 int count, count_mol, i, mb;
1155 gmx_molblock_t *molb;
1156 t_params *plist;
1157 char buf[STRLEN];
1159 count = 0;
1160 for (mb = 0; mb < mtop->nmolblock; mb++)
1162 count_mol = 0;
1163 molb = &mtop->molblock[mb];
1164 plist = mi[molb->type].plist;
1166 for (i = 0; i < F_NRE; i++)
1168 if (i == F_SETTLE)
1170 count_mol += 3*plist[i].nr;
1172 else if (interaction_function[i].flags & IF_CONSTRAINT)
1174 count_mol += plist[i].nr;
1178 if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1180 sprintf(buf,
1181 "Molecule type '%s' has %d constraints.\n"
1182 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1183 *mi[molb->type].name, count_mol,
1184 nrdf_internal(&mi[molb->type].atoms));
1185 warning(wi, buf);
1187 count += molb->nmol*count_mol;
1190 return count;
1193 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1195 int i, nmiss, natoms, mt;
1196 real q;
1197 const t_atoms *atoms;
1199 nmiss = 0;
1200 for (mt = 0; mt < sys->nmoltype; mt++)
1202 atoms = &sys->moltype[mt].atoms;
1203 natoms = atoms->nr;
1205 for (i = 0; i < natoms; i++)
1207 q = atoms->atom[i].q;
1208 if ((get_atomtype_radius(atoms->atom[i].type, atype) == 0 ||
1209 get_atomtype_vol(atoms->atom[i].type, atype) == 0 ||
1210 get_atomtype_surftens(atoms->atom[i].type, atype) == 0 ||
1211 get_atomtype_gb_radius(atoms->atom[i].type, atype) == 0 ||
1212 get_atomtype_S_hct(atoms->atom[i].type, atype) == 0) &&
1213 q != 0)
1215 fprintf(stderr, "\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1216 get_atomtype_name(atoms->atom[i].type, atype), q);
1217 nmiss++;
1222 if (nmiss > 0)
1224 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);
1229 static void check_gbsa_params(gpp_atomtype_t atype)
1231 int nmiss, i;
1233 /* If we are doing GBSA, check that we got the parameters we need
1234 * This checking is to see if there are GBSA paratmeters for all
1235 * atoms in the force field. To go around this for testing purposes
1236 * comment out the nerror++ counter temporarily
1238 nmiss = 0;
1239 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1241 if (get_atomtype_radius(i, atype) < 0 ||
1242 get_atomtype_vol(i, atype) < 0 ||
1243 get_atomtype_surftens(i, atype) < 0 ||
1244 get_atomtype_gb_radius(i, atype) < 0 ||
1245 get_atomtype_S_hct(i, atype) < 0)
1247 fprintf(stderr, "\nGB parameter(s) missing or negative for atom type '%s'\n",
1248 get_atomtype_name(i, atype));
1249 nmiss++;
1253 if (nmiss > 0)
1255 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);
1260 static real calc_temp(const gmx_mtop_t *mtop,
1261 const t_inputrec *ir,
1262 rvec *v)
1264 double sum_mv2;
1265 gmx_mtop_atomloop_all_t aloop;
1266 t_atom *atom;
1267 int a;
1268 int nrdf, g;
1270 sum_mv2 = 0;
1272 aloop = gmx_mtop_atomloop_all_init(mtop);
1273 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
1275 sum_mv2 += atom->m*norm2(v[a]);
1278 nrdf = 0;
1279 for (g = 0; g < ir->opts.ngtc; g++)
1281 nrdf += ir->opts.nrdf[g];
1284 return sum_mv2/(nrdf*BOLTZ);
1287 static real get_max_reference_temp(const t_inputrec *ir,
1288 warninp_t wi)
1290 real ref_t;
1291 int i;
1292 gmx_bool bNoCoupl;
1294 ref_t = 0;
1295 bNoCoupl = FALSE;
1296 for (i = 0; i < ir->opts.ngtc; i++)
1298 if (ir->opts.tau_t[i] < 0)
1300 bNoCoupl = TRUE;
1302 else
1304 ref_t = max(ref_t, ir->opts.ref_t[i]);
1308 if (bNoCoupl)
1310 char buf[STRLEN];
1312 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.",
1313 ref_t);
1314 warning(wi, buf);
1317 return ref_t;
1320 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1321 t_inputrec *ir,
1322 real buffer_temp,
1323 matrix box,
1324 warninp_t wi)
1326 int i;
1327 verletbuf_list_setup_t ls;
1328 real rlist_1x1;
1329 int n_nonlin_vsite;
1330 char warn_buf[STRLEN];
1332 printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol, buffer_temp);
1334 /* Calculate the buffer size for simple atom vs atoms list */
1335 ls.cluster_size_i = 1;
1336 ls.cluster_size_j = 1;
1337 calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1338 &ls, &n_nonlin_vsite, &rlist_1x1);
1340 /* Set the pair-list buffer size in ir */
1341 verletbuf_get_list_setup(FALSE, &ls);
1342 calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1343 &ls, &n_nonlin_vsite, &ir->rlist);
1345 if (n_nonlin_vsite > 0)
1347 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);
1348 warning_note(wi, warn_buf);
1351 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1352 1, 1, rlist_1x1, rlist_1x1-max(ir->rvdw, ir->rcoulomb));
1354 ir->rlistlong = ir->rlist;
1355 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1356 ls.cluster_size_i, ls.cluster_size_j,
1357 ir->rlist, ir->rlist-max(ir->rvdw, ir->rcoulomb));
1359 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
1361 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)));
1365 int gmx_grompp(int argc, char *argv[])
1367 static const char *desc[] = {
1368 "[THISMODULE] (the gromacs preprocessor)",
1369 "reads a molecular topology file, checks the validity of the",
1370 "file, expands the topology from a molecular description to an atomic",
1371 "description. The topology file contains information about",
1372 "molecule types and the number of molecules, the preprocessor",
1373 "copies each molecule as needed. ",
1374 "There is no limitation on the number of molecule types. ",
1375 "Bonds and bond-angles can be converted into constraints, separately",
1376 "for hydrogens and heavy atoms.",
1377 "Then a coordinate file is read and velocities can be generated",
1378 "from a Maxwellian distribution if requested.",
1379 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1380 "(eg. number of MD steps, time step, cut-off), and others such as",
1381 "NEMD parameters, which are corrected so that the net acceleration",
1382 "is zero.",
1383 "Eventually a binary file is produced that can serve as the sole input",
1384 "file for the MD program.[PAR]",
1386 "[THISMODULE] uses the atom names from the topology file. The atom names",
1387 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1388 "warnings when they do not match the atom names in the topology.",
1389 "Note that the atom names are irrelevant for the simulation as",
1390 "only the atom types are used for generating interaction parameters.[PAR]",
1392 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1393 "etc. The preprocessor supports the following keywords:[PAR]",
1394 "#ifdef VARIABLE[BR]",
1395 "#ifndef VARIABLE[BR]",
1396 "#else[BR]",
1397 "#endif[BR]",
1398 "#define VARIABLE[BR]",
1399 "#undef VARIABLE[BR]"
1400 "#include \"filename\"[BR]",
1401 "#include <filename>[PAR]",
1402 "The functioning of these statements in your topology may be modulated by",
1403 "using the following two flags in your [TT].mdp[tt] file:[PAR]",
1404 "[TT]define = -DVARIABLE1 -DVARIABLE2[BR]",
1405 "include = -I/home/john/doe[tt][BR]",
1406 "For further information a C-programming textbook may help you out.",
1407 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1408 "topology file written out so that you can verify its contents.[PAR]",
1410 /* cpp has been unnecessary for some time, hasn't it?
1411 "If your system does not have a C-preprocessor, you can still",
1412 "use [TT]grompp[tt], but you do not have access to the features ",
1413 "from the cpp. Command line options to the C-preprocessor can be given",
1414 "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1417 "When using position restraints a file with restraint coordinates",
1418 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1419 "with respect to the conformation from the [TT]-c[tt] option.",
1420 "For free energy calculation the the coordinates for the B topology",
1421 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1422 "those of the A topology.[PAR]",
1424 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1425 "The last frame with coordinates and velocities will be read,",
1426 "unless the [TT]-time[tt] option is used. Only if this information",
1427 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1428 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1429 "in your [TT].mdp[tt] file. An energy file can be supplied with",
1430 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1431 "variables.[PAR]",
1433 "[THISMODULE] can be used to restart simulations (preserving",
1434 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1435 "However, for simply changing the number of run steps to extend",
1436 "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1437 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1438 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1439 "like output frequency, then supplying the checkpoint file to",
1440 "[THISMODULE] with [TT]-t[tt] along with a new [TT].mdp[tt] file",
1441 "with [TT]-f[tt] is the recommended procedure.[PAR]",
1443 "By default, all bonded interactions which have constant energy due to",
1444 "virtual site constructions will be removed. If this constant energy is",
1445 "not zero, this will result in a shift in the total energy. All bonded",
1446 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1447 "all constraints for distances which will be constant anyway because",
1448 "of virtual site constructions will be removed. If any constraints remain",
1449 "which involve virtual sites, a fatal error will result.[PAR]"
1451 "To verify your run input file, please take note of all warnings",
1452 "on the screen, and correct where necessary. Do also look at the contents",
1453 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1454 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1455 "with the [TT]-debug[tt] option which will give you more information",
1456 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1457 "can see the contents of the run input file with the [gmx-dump]",
1458 "program. [gmx-check] can be used to compare the contents of two",
1459 "run input files.[PAR]"
1461 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1462 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1463 "harmless, but usually they are not. The user is advised to carefully",
1464 "interpret the output messages before attempting to bypass them with",
1465 "this option."
1467 t_gromppopts *opts;
1468 gmx_mtop_t *sys;
1469 int nmi;
1470 t_molinfo *mi;
1471 gpp_atomtype_t atype;
1472 t_inputrec *ir;
1473 int natoms, nvsite, comb, mt;
1474 t_params *plist;
1475 t_state state;
1476 matrix box;
1477 real max_spacing, fudgeQQ;
1478 double reppow;
1479 char fn[STRLEN], fnB[STRLEN];
1480 const char *mdparin;
1481 int ntype;
1482 gmx_bool bNeedVel, bGenVel;
1483 gmx_bool have_atomnumber;
1484 int n12, n13, n14;
1485 t_params *gb_plist = NULL;
1486 gmx_genborn_t *born = NULL;
1487 output_env_t oenv;
1488 gmx_bool bVerbose = FALSE;
1489 warninp_t wi;
1490 char warn_buf[STRLEN];
1491 unsigned int useed;
1492 t_atoms IMDatoms; /* Atoms to be operated on interactively (IMD) */
1494 t_filenm fnm[] = {
1495 { efMDP, NULL, NULL, ffREAD },
1496 { efMDP, "-po", "mdout", ffWRITE },
1497 { efSTX, "-c", NULL, ffREAD },
1498 { efSTX, "-r", NULL, ffOPTRD },
1499 { efSTX, "-rb", NULL, ffOPTRD },
1500 { efNDX, NULL, NULL, ffOPTRD },
1501 { efTOP, NULL, NULL, ffREAD },
1502 { efTOP, "-pp", "processed", ffOPTWR },
1503 { efTPX, "-o", NULL, ffWRITE },
1504 { efTRN, "-t", NULL, ffOPTRD },
1505 { efEDR, "-e", NULL, ffOPTRD },
1506 /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1507 { efGRO, "-imd", "imdgroup", ffOPTWR },
1508 { efTRN, "-ref", "rotref", ffOPTRW }
1510 #define NFILE asize(fnm)
1512 /* Command line options */
1513 static gmx_bool bRenum = TRUE;
1514 static gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1515 static int i, maxwarn = 0;
1516 static real fr_time = -1;
1517 t_pargs pa[] = {
1518 { "-v", FALSE, etBOOL, {&bVerbose},
1519 "Be loud and noisy" },
1520 { "-time", FALSE, etREAL, {&fr_time},
1521 "Take frame at or first after this time." },
1522 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1523 "Remove constant bonded interactions with virtual sites" },
1524 { "-maxwarn", FALSE, etINT, {&maxwarn},
1525 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1526 { "-zero", FALSE, etBOOL, {&bZero},
1527 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1528 { "-renum", FALSE, etBOOL, {&bRenum},
1529 "Renumber atomtypes and minimize number of atomtypes" }
1532 /* Initiate some variables */
1533 snew(ir, 1);
1534 snew(opts, 1);
1535 init_ir(ir, opts);
1537 /* Parse the command line */
1538 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1539 asize(desc), desc, 0, NULL, &oenv))
1541 return 0;
1544 wi = init_warning(TRUE, maxwarn);
1546 /* PARAMETER file processing */
1547 mdparin = opt2fn("-f", NFILE, fnm);
1548 set_warning_line(wi, mdparin, -1);
1549 get_ir(mdparin, opt2fn("-po", NFILE, fnm), ir, opts, wi);
1551 if (bVerbose)
1553 fprintf(stderr, "checking input for internal consistency...\n");
1555 check_ir(mdparin, ir, opts, wi);
1557 if (ir->ld_seed == -1)
1559 ir->ld_seed = (gmx_int64_t)gmx_rng_make_seed();
1560 fprintf(stderr, "Setting the LD random seed to %"GMX_PRId64 "\n", ir->ld_seed);
1563 if (ir->expandedvals->lmc_seed == -1)
1565 ir->expandedvals->lmc_seed = (int)gmx_rng_make_seed();
1566 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1569 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1570 bGenVel = (bNeedVel && opts->bGenVel);
1571 if (bGenVel && ir->bContinuation)
1573 sprintf(warn_buf,
1574 "Generating velocities is inconsistent with attempting "
1575 "to continue a previous run. Choose only one of "
1576 "gen-vel = yes and continuation = yes.");
1577 warning_error(wi, warn_buf);
1580 snew(plist, F_NRE);
1581 init_plist(plist);
1582 snew(sys, 1);
1583 atype = init_atomtype();
1584 if (debug)
1586 pr_symtab(debug, 0, "Just opened", &sys->symtab);
1589 strcpy(fn, ftp2fn(efTOP, NFILE, fnm));
1590 if (!gmx_fexist(fn))
1592 gmx_fatal(FARGS, "%s does not exist", fn);
1594 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1595 opts, ir, bZero, bGenVel, bVerbose, &state,
1596 atype, sys, &nmi, &mi, plist, &comb, &reppow, &fudgeQQ,
1597 opts->bMorse,
1598 wi);
1600 if (debug)
1602 pr_symtab(debug, 0, "After new_status", &sys->symtab);
1605 nvsite = 0;
1606 /* set parameters for virtual site construction (not for vsiten) */
1607 for (mt = 0; mt < sys->nmoltype; mt++)
1609 nvsite +=
1610 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1612 /* now throw away all obsolete bonds, angles and dihedrals: */
1613 /* note: constraints are ALWAYS removed */
1614 if (nvsite)
1616 for (mt = 0; mt < sys->nmoltype; mt++)
1618 clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1622 if (ir->cutoff_scheme == ecutsVERLET)
1624 fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1625 ecutscheme_names[ir->cutoff_scheme]);
1627 /* Remove all charge groups */
1628 gmx_mtop_remove_chargegroups(sys);
1631 if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1633 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1635 sprintf(warn_buf, "Can not do %s with %s, use %s",
1636 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1637 warning_error(wi, warn_buf);
1639 if (ir->bPeriodicMols)
1641 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1642 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1643 warning_error(wi, warn_buf);
1647 if (EI_SD (ir->eI) && ir->etc != etcNO)
1649 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1652 /* If we are doing QM/MM, check that we got the atom numbers */
1653 have_atomnumber = TRUE;
1654 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1656 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1658 if (!have_atomnumber && ir->bQMMM)
1660 warning_error(wi,
1661 "\n"
1662 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1663 "field you are using does not contain atom numbers fields. This is an\n"
1664 "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1665 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1666 "an integer just before the mass column in ffXXXnb.itp.\n"
1667 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1670 if (ir->bAdress)
1672 if ((ir->adress->const_wf > 1) || (ir->adress->const_wf < 0))
1674 warning_error(wi, "AdResS contant weighting function should be between 0 and 1\n\n");
1676 /** TODO check size of ex+hy width against box size */
1679 /* Check for errors in the input now, since they might cause problems
1680 * during processing further down.
1682 check_warning_error(wi, FARGS);
1684 if (opt2bSet("-r", NFILE, fnm))
1686 sprintf(fn, "%s", opt2fn("-r", NFILE, fnm));
1688 else
1690 sprintf(fn, "%s", opt2fn("-c", NFILE, fnm));
1692 if (opt2bSet("-rb", NFILE, fnm))
1694 sprintf(fnB, "%s", opt2fn("-rb", NFILE, fnm));
1696 else
1698 strcpy(fnB, fn);
1701 if (nint_ftype(sys, mi, F_POSRES) > 0 || nint_ftype(sys, mi, F_FBPOSRES) > 0)
1703 if (bVerbose)
1705 fprintf(stderr, "Reading position restraint coords from %s", fn);
1706 if (strcmp(fn, fnB) == 0)
1708 fprintf(stderr, "\n");
1710 else
1712 fprintf(stderr, " and %s\n", fnB);
1715 gen_posres(sys, mi, fn, fnB,
1716 ir->refcoord_scaling, ir->ePBC,
1717 ir->posres_com, ir->posres_comB,
1718 wi);
1721 /* If we are using CMAP, setup the pre-interpolation grid */
1722 if (plist->ncmap > 0)
1724 init_cmap_grid(&sys->ffparams.cmap_grid, plist->nc, plist->grid_spacing);
1725 setup_cmap(plist->grid_spacing, plist->nc, plist->cmap, &sys->ffparams.cmap_grid);
1728 set_wall_atomtype(atype, opts, ir, wi);
1729 if (bRenum)
1731 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1732 ntype = get_atomtype_ntypes(atype);
1735 if (ir->implicit_solvent != eisNO)
1737 /* Now we have renumbered the atom types, we can check the GBSA params */
1738 check_gbsa_params(atype);
1740 /* Check that all atoms that have charge and/or LJ-parameters also have
1741 * sensible GB-parameters
1743 check_gbsa_params_charged(sys, atype);
1746 /* PELA: Copy the atomtype data to the topology atomtype list */
1747 copy_atomtype_atomtypes(atype, &(sys->atomtypes));
1749 if (debug)
1751 pr_symtab(debug, 0, "After renum_atype", &sys->symtab);
1754 if (bVerbose)
1756 fprintf(stderr, "converting bonded parameters...\n");
1759 ntype = get_atomtype_ntypes(atype);
1760 convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1762 if (debug)
1764 pr_symtab(debug, 0, "After convert_params", &sys->symtab);
1767 /* set ptype to VSite for virtual sites */
1768 for (mt = 0; mt < sys->nmoltype; mt++)
1770 set_vsites_ptype(FALSE, &sys->moltype[mt]);
1772 if (debug)
1774 pr_symtab(debug, 0, "After virtual sites", &sys->symtab);
1776 /* Check velocity for virtual sites and shells */
1777 if (bGenVel)
1779 check_vel(sys, state.v);
1782 /* check masses */
1783 check_mol(sys, wi);
1785 for (i = 0; i < sys->nmoltype; i++)
1787 check_cg_sizes(ftp2fn(efTOP, NFILE, fnm), &sys->moltype[i].cgs, wi);
1790 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1792 check_bonds_timestep(sys, ir->delta_t, wi);
1795 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
1797 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.");
1800 check_warning_error(wi, FARGS);
1802 if (bVerbose)
1804 fprintf(stderr, "initialising group options...\n");
1806 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
1807 sys, bVerbose, ir,
1808 bGenVel ? state.v : NULL,
1809 wi);
1811 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0 &&
1812 ir->nstlist > 1)
1814 if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
1816 real buffer_temp;
1818 if (EI_MD(ir->eI) && ir->etc == etcNO)
1820 if (bGenVel)
1822 buffer_temp = opts->tempi;
1824 else
1826 buffer_temp = calc_temp(sys, ir, state.v);
1828 if (buffer_temp > 0)
1830 sprintf(warn_buf, "NVE simulation: will use the initial temperature of %.3f K for determining the Verlet buffer size", buffer_temp);
1831 warning_note(wi, warn_buf);
1833 else
1835 sprintf(warn_buf, "NVE simulation with an initial temperature of zero: will use a Verlet buffer of %d%%. Check your energy drift!",
1836 (int)(verlet_buffer_ratio_NVE_T0*100 + 0.5));
1837 warning_note(wi, warn_buf);
1840 else
1842 buffer_temp = get_max_reference_temp(ir, wi);
1845 if (EI_MD(ir->eI) && ir->etc == etcNO && buffer_temp == 0)
1847 /* NVE with initial T=0: we add a fixed ratio to rlist.
1848 * Since we don't actually use verletbuf_tol, we set it to -1
1849 * so it can't be misused later.
1851 ir->rlist *= 1.0 + verlet_buffer_ratio_NVE_T0;
1852 ir->verletbuf_tol = -1;
1854 else
1856 /* We warn for NVE simulations with >1(.1)% drift tolerance */
1857 const real drift_tol = 0.01;
1858 real ener_runtime;
1860 /* We use 2 DOF per atom = 2kT pot+kin energy, to be on
1861 * the safe side with constraints (without constraints: 3 DOF).
1863 ener_runtime = 2*BOLTZ*buffer_temp/(ir->nsteps*ir->delta_t);
1865 if (EI_MD(ir->eI) && ir->etc == etcNO && ir->nstlist > 1 &&
1866 ir->nsteps > 0 &&
1867 ir->verletbuf_tol > 1.1*drift_tol*ener_runtime)
1869 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.",
1870 ir->verletbuf_tol, ir->nsteps*ir->delta_t,
1871 (int)(ir->verletbuf_tol/ener_runtime*100 + 0.5),
1872 (int)(100*drift_tol + 0.5),
1873 drift_tol*ener_runtime);
1874 warning_note(wi, warn_buf);
1877 set_verlet_buffer(sys, ir, buffer_temp, state.box, wi);
1882 /* Init the temperature coupling state */
1883 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
1885 if (bVerbose)
1887 fprintf(stderr, "Checking consistency between energy and charge groups...\n");
1889 check_eg_vs_cg(sys);
1891 if (debug)
1893 pr_symtab(debug, 0, "After index", &sys->symtab);
1896 triple_check(mdparin, ir, sys, wi);
1897 close_symtab(&sys->symtab);
1898 if (debug)
1900 pr_symtab(debug, 0, "After close", &sys->symtab);
1903 /* make exclusions between QM atoms */
1904 if (ir->bQMMM)
1906 if (ir->QMMMscheme == eQMMMschemenormal && ir->ns_type == ensSIMPLE)
1908 gmx_fatal(FARGS, "electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
1910 else
1912 generate_qmexcl(sys, ir, wi);
1916 if (ftp2bSet(efTRN, NFILE, fnm))
1918 if (bVerbose)
1920 fprintf(stderr, "getting data from old trajectory ...\n");
1922 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
1923 bNeedVel, bGenVel, fr_time, ir, &state, sys, oenv);
1926 if (ir->ePBC == epbcXY && ir->nwall != 2)
1928 clear_rvec(state.box[ZZ]);
1931 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
1933 set_warning_line(wi, mdparin, -1);
1934 check_chargegroup_radii(sys, ir, state.x, wi);
1937 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
1939 /* Calculate the optimal grid dimensions */
1940 copy_mat(state.box, box);
1941 if (ir->ePBC == epbcXY && ir->nwall == 2)
1943 svmul(ir->wall_ewald_zfac, box[ZZ], box[ZZ]);
1945 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
1947 /* Mark fourier_spacing as not used */
1948 ir->fourier_spacing = 0;
1950 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
1952 set_warning_line(wi, mdparin, -1);
1953 warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
1955 max_spacing = calc_grid(stdout, box, ir->fourier_spacing,
1956 &(ir->nkx), &(ir->nky), &(ir->nkz));
1959 /* MRS: eventually figure out better logic for initializing the fep
1960 values that makes declaring the lambda and declaring the state not
1961 potentially conflict if not handled correctly. */
1962 if (ir->efep != efepNO)
1964 state.fep_state = ir->fepvals->init_fep_state;
1965 for (i = 0; i < efptNR; i++)
1967 /* init_lambda trumps state definitions*/
1968 if (ir->fepvals->init_lambda >= 0)
1970 state.lambda[i] = ir->fepvals->init_lambda;
1972 else
1974 if (ir->fepvals->all_lambda[i] == NULL)
1976 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
1978 else
1980 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
1986 if (ir->ePull != epullNO)
1988 set_pull_init(ir, sys, state.x, state.box, state.lambda[efptMASS], oenv, opts->pull_start);
1991 if (ir->bRot)
1993 set_reference_positions(ir->rot, state.x, state.box,
1994 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
1995 wi);
1998 /* reset_multinr(sys); */
2000 if (EEL_PME(ir->coulombtype))
2002 float ratio = pme_load_estimate(sys, ir, state.box);
2003 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
2004 /* With free energy we might need to do PME both for the A and B state
2005 * charges. This will double the cost, but the optimal performance will
2006 * then probably be at a slightly larger cut-off and grid spacing.
2008 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
2009 (ir->efep != efepNO && ratio > 2.0/3.0))
2011 warning_note(wi,
2012 "The optimal PME mesh load for parallel simulations is below 0.5\n"
2013 "and for highly parallel simulations between 0.25 and 0.33,\n"
2014 "for higher performance, increase the cut-off and the PME grid spacing.\n");
2015 if (ir->efep != efepNO)
2017 warning_note(wi,
2018 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
2024 char warn_buf[STRLEN];
2025 double cio = compute_io(ir, sys->natoms, &sys->groups, F_NRE, 1);
2026 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
2027 if (cio > 2000)
2029 set_warning_line(wi, mdparin, -1);
2030 warning_note(wi, warn_buf);
2032 else
2034 printf("%s\n", warn_buf);
2038 if (bVerbose)
2040 fprintf(stderr, "writing run input file...\n");
2043 done_warning(wi, FARGS);
2044 write_tpx_state(ftp2fn(efTPX, NFILE, fnm), ir, &state, sys);
2046 /* Output IMD group, if bIMD is TRUE */
2047 write_IMDgroup_to_file(ir->bIMD, ir, &state, sys, NFILE, fnm);
2049 done_atomtype(atype);
2050 done_mtop(sys, TRUE);
2051 done_inputrec_strings();
2053 return 0;