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.
43 #include <sys/types.h>
56 #include "gromacs/fileio/confio.h"
60 #include "grompp-impl.h"
63 #include "gromacs/fileio/futil.h"
64 #include "gromacs/commandline/pargs.h"
66 #include "sortwater.h"
68 #include "gmx_fatal.h"
71 #include "gromacs/fileio/gmxfio.h"
72 #include "gromacs/fileio/trnio.h"
73 #include "gromacs/fileio/tpxio.h"
74 #include "gromacs/fileio/trxio.h"
75 #include "vsite_parm.h"
79 #include "gromacs/fileio/enxio.h"
81 #include "compute_io.h"
82 #include "gpp_atomtype.h"
83 #include "mtop_util.h"
85 #include "calc_verletbuf.h"
89 static int rm_interactions(int ifunc, int nrmols, t_molinfo mols[])
94 /* For all the molecule types */
95 for (i = 0; i < nrmols; i++)
97 n += mols[i].plist[ifunc].nr;
98 mols[i].plist[ifunc].nr = 0;
103 static int check_atom_names(const char *fn1, const char *fn2,
104 gmx_mtop_t *mtop, t_atoms *at)
106 int mb, m, i, j, nmismatch;
108 #define MAXMISMATCH 20
110 if (mtop->natoms != at->nr)
112 gmx_incons("comparing atom names");
117 for (mb = 0; mb < mtop->nmolblock; mb++)
119 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
120 for (m = 0; m < mtop->molblock[mb].nmol; m++)
122 for (j = 0; j < tat->nr; j++)
124 if (strcmp( *(tat->atomname[j]), *(at->atomname[i]) ) != 0)
126 if (nmismatch < MAXMISMATCH)
129 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
130 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
132 else if (nmismatch == MAXMISMATCH)
134 fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
146 static void check_eg_vs_cg(gmx_mtop_t *mtop)
148 int astart, mb, m, cg, j, firstj;
149 unsigned char firsteg, eg;
152 /* Go through all the charge groups and make sure all their
153 * atoms are in the same energy group.
157 for (mb = 0; mb < mtop->nmolblock; mb++)
159 molt = &mtop->moltype[mtop->molblock[mb].type];
160 for (m = 0; m < mtop->molblock[mb].nmol; m++)
162 for (cg = 0; cg < molt->cgs.nr; cg++)
164 /* Get the energy group of the first atom in this charge group */
165 firstj = astart + molt->cgs.index[cg];
166 firsteg = ggrpnr(&mtop->groups, egcENER, firstj);
167 for (j = molt->cgs.index[cg]+1; j < molt->cgs.index[cg+1]; j++)
169 eg = ggrpnr(&mtop->groups, egcENER, astart+j);
172 gmx_fatal(FARGS, "atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
173 firstj+1, astart+j+1, cg+1, *molt->name);
177 astart += molt->atoms.nr;
182 static void check_cg_sizes(const char *topfn, t_block *cgs, warninp_t wi)
185 char warn_buf[STRLEN];
188 for (cg = 0; cg < cgs->nr; cg++)
190 maxsize = max(maxsize, cgs->index[cg+1]-cgs->index[cg]);
193 if (maxsize > MAX_CHARGEGROUP_SIZE)
195 gmx_fatal(FARGS, "The largest charge group contains %d atoms. The maximum is %d.", maxsize, MAX_CHARGEGROUP_SIZE);
197 else if (maxsize > 10)
199 set_warning_line(wi, topfn, -1);
201 "The largest charge group contains %d atoms.\n"
202 "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"
203 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
204 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
206 warning_note(wi, warn_buf);
210 static void check_bonds_timestep(gmx_mtop_t *mtop, double dt, warninp_t wi)
212 /* This check is not intended to ensure accurate integration,
213 * rather it is to signal mistakes in the mdp settings.
214 * A common mistake is to forget to turn on constraints
215 * for MD after energy minimization with flexible bonds.
216 * This check can also detect too large time steps for flexible water
217 * models, but such errors will often be masked by the constraints
218 * mdp options, which turns flexible water into water with bond constraints,
219 * but without an angle constraint. Unfortunately such incorrect use
220 * of water models can not easily be detected without checking
221 * for specific model names.
223 * The stability limit of leap-frog or velocity verlet is 4.44 steps
224 * per oscillational period.
225 * But accurate bonds distributions are lost far before that limit.
226 * To allow relatively common schemes (although not common with Gromacs)
227 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
228 * we set the note limit to 10.
230 int min_steps_warn = 5;
231 int min_steps_note = 10;
234 gmx_moltype_t *moltype, *w_moltype;
236 t_ilist *ilist, *ilb, *ilc, *ils;
238 int i, a1, a2, w_a1, w_a2, j;
239 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
240 gmx_bool bFound, bWater, bWarn;
241 char warn_buf[STRLEN];
243 ip = mtop->ffparams.iparams;
245 twopi2 = sqr(2*M_PI);
247 limit2 = sqr(min_steps_note*dt);
253 for (molt = 0; molt < mtop->nmoltype; molt++)
255 moltype = &mtop->moltype[molt];
256 atom = moltype->atoms.atom;
257 ilist = moltype->ilist;
258 ilc = &ilist[F_CONSTR];
259 ils = &ilist[F_SETTLE];
260 for (ftype = 0; ftype < F_NRE; ftype++)
262 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
268 for (i = 0; i < ilb->nr; i += 3)
270 fc = ip[ilb->iatoms[i]].harmonic.krA;
271 re = ip[ilb->iatoms[i]].harmonic.rA;
272 if (ftype == F_G96BONDS)
274 /* Convert squared sqaure fc to harmonic fc */
277 a1 = ilb->iatoms[i+1];
278 a2 = ilb->iatoms[i+2];
281 if (fc > 0 && m1 > 0 && m2 > 0)
283 period2 = twopi2*m1*m2/((m1 + m2)*fc);
287 period2 = GMX_FLOAT_MAX;
291 fprintf(debug, "fc %g m1 %g m2 %g period %g\n",
292 fc, m1, m2, sqrt(period2));
294 if (period2 < limit2)
297 for (j = 0; j < ilc->nr; j += 3)
299 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
300 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
305 for (j = 0; j < ils->nr; j += 4)
307 if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
308 (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
314 (w_moltype == NULL || period2 < w_period2))
326 if (w_moltype != NULL)
328 bWarn = (w_period2 < sqr(min_steps_warn*dt));
329 /* A check that would recognize most water models */
330 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
331 w_moltype->atoms.nr <= 5);
332 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 w_a1+1, *w_moltype->atoms.atomname[w_a1],
336 w_a2+1, *w_moltype->atoms.atomname[w_a2],
337 sqrt(w_period2), bWarn ? min_steps_warn : min_steps_note, dt,
339 "Maybe you asked for fexible water." :
340 "Maybe you forgot to change the constraints mdp option.");
343 warning(wi, warn_buf);
347 warning_note(wi, warn_buf);
352 static void check_vel(gmx_mtop_t *mtop, rvec v[])
354 gmx_mtop_atomloop_all_t aloop;
358 aloop = gmx_mtop_atomloop_all_init(mtop);
359 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
361 if (atom->ptype == eptShell ||
362 atom->ptype == eptBond ||
363 atom->ptype == eptVSite)
370 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
375 for (mb = 0; mb < mtop->nmolblock; mb++)
377 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
383 /* This routine reorders the molecule type array
384 * in the order of use in the molblocks,
385 * unused molecule types are deleted.
387 static void renumber_moltypes(gmx_mtop_t *sys,
388 int *nmolinfo, t_molinfo **molinfo)
390 int *order, norder, i;
394 snew(order, *nmolinfo);
396 for (mb = 0; mb < sys->nmolblock; mb++)
398 for (i = 0; i < norder; i++)
400 if (order[i] == sys->molblock[mb].type)
407 /* This type did not occur yet, add it */
408 order[norder] = sys->molblock[mb].type;
409 /* Renumber the moltype in the topology */
412 sys->molblock[mb].type = i;
415 /* We still need to reorder the molinfo structs */
417 for (mi = 0; mi < *nmolinfo; mi++)
419 for (i = 0; i < norder; i++)
428 done_mi(&(*molinfo)[mi]);
432 minew[i] = (*molinfo)[mi];
441 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
446 mtop->nmoltype = nmi;
447 snew(mtop->moltype, nmi);
448 for (m = 0; m < nmi; m++)
450 molt = &mtop->moltype[m];
451 molt->name = mi[m].name;
452 molt->atoms = mi[m].atoms;
453 /* ilists are copied later */
454 molt->cgs = mi[m].cgs;
455 molt->excls = mi[m].excls;
460 new_status(const char *topfile, const char *topppfile, const char *confin,
461 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
462 gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
463 gpp_atomtype_t atype, gmx_mtop_t *sys,
464 int *nmi, t_molinfo **mi, t_params plist[],
465 int *comb, double *reppow, real *fudgeQQ,
469 t_molinfo *molinfo = NULL;
471 gmx_molblock_t *molblock, *molbs;
473 int mb, i, nrmols, nmismatch;
475 gmx_bool bGB = FALSE;
476 char warn_buf[STRLEN];
480 /* Set gmx_boolean for GB */
481 if (ir->implicit_solvent)
486 /* TOPOLOGY processing */
487 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
488 plist, comb, reppow, fudgeQQ,
489 atype, &nrmols, &molinfo, ir,
490 &nmolblock, &molblock, bGB,
494 snew(sys->molblock, nmolblock);
497 for (mb = 0; mb < nmolblock; mb++)
499 if (sys->nmolblock > 0 &&
500 molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
502 /* Merge consecutive blocks with the same molecule type */
503 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
504 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
506 else if (molblock[mb].nmol > 0)
508 /* Add a new molblock to the topology */
509 molbs = &sys->molblock[sys->nmolblock];
510 *molbs = molblock[mb];
511 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
512 molbs->nposres_xA = 0;
513 molbs->nposres_xB = 0;
514 sys->natoms += molbs->nmol*molbs->natoms_mol;
518 if (sys->nmolblock == 0)
520 gmx_fatal(FARGS, "No molecules were defined in the system");
523 renumber_moltypes(sys, &nrmols, &molinfo);
527 convert_harmonics(nrmols, molinfo, atype);
530 if (ir->eDisre == edrNone)
532 i = rm_interactions(F_DISRES, nrmols, molinfo);
535 set_warning_line(wi, "unknown", -1);
536 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
537 warning_note(wi, warn_buf);
540 if (opts->bOrire == FALSE)
542 i = rm_interactions(F_ORIRES, nrmols, molinfo);
545 set_warning_line(wi, "unknown", -1);
546 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
547 warning_note(wi, warn_buf);
551 /* Copy structures from msys to sys */
552 molinfo2mtop(nrmols, molinfo, sys);
554 gmx_mtop_finalize(sys);
556 /* COORDINATE file processing */
559 fprintf(stderr, "processing coordinates...\n");
562 get_stx_coordnum(confin, &state->natoms);
563 if (state->natoms != sys->natoms)
565 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
566 " does not match topology (%s, %d)",
567 confin, state->natoms, topfile, sys->natoms);
571 /* make space for coordinates and velocities */
574 init_t_atoms(confat, state->natoms, FALSE);
575 init_state(state, state->natoms, 0, 0, 0, 0);
576 read_stx_conf(confin, title, confat, state->x, state->v, NULL, state->box);
577 /* This call fixes the box shape for runs with pressure scaling */
578 set_box_rel(ir, state);
580 nmismatch = check_atom_names(topfile, confin, sys, confat);
581 free_t_atoms(confat, TRUE);
586 sprintf(buf, "%d non-matching atom name%s\n"
587 "atom names from %s will be used\n"
588 "atom names from %s will be ignored\n",
589 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
594 fprintf(stderr, "double-checking input for internal consistency...\n");
596 double_check(ir, state->box, nint_ftype(sys, molinfo, F_CONSTR), wi);
602 gmx_mtop_atomloop_all_t aloop;
605 snew(mass, state->natoms);
606 aloop = gmx_mtop_atomloop_all_init(sys);
607 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
612 if (opts->seed == -1)
614 opts->seed = make_seed();
615 fprintf(stderr, "Setting gen_seed to %d\n", opts->seed);
617 maxwell_speed(opts->tempi, opts->seed, sys, state->v);
619 stop_cm(stdout, state->natoms, mass, state->x, state->v);
627 static void copy_state(const char *slog, t_trxframe *fr,
628 gmx_bool bReadVel, t_state *state,
633 if (fr->not_ok & FRAME_NOT_OK)
635 gmx_fatal(FARGS, "Can not start from an incomplete frame");
639 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
643 for (i = 0; i < state->natoms; i++)
645 copy_rvec(fr->x[i], state->x[i]);
651 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
653 for (i = 0; i < state->natoms; i++)
655 copy_rvec(fr->v[i], state->v[i]);
660 copy_mat(fr->box, state->box);
663 *use_time = fr->time;
666 static void cont_status(const char *slog, const char *ener,
667 gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
668 t_inputrec *ir, t_state *state,
670 const output_env_t oenv)
671 /* If fr_time == -1 read the last frame available which is complete */
679 bReadVel = (bNeedVel && !bGenVel);
682 "Reading Coordinates%s and Box size from old trajectory\n",
683 bReadVel ? ", Velocities" : "");
686 fprintf(stderr, "Will read whole trajectory\n");
690 fprintf(stderr, "Will read till time %g\n", fr_time);
696 fprintf(stderr, "Velocities generated: "
697 "ignoring velocities in input trajectory\n");
699 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
703 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
709 "WARNING: Did not find a frame with velocities in file %s,\n"
710 " all velocities will be set to zero!\n\n", slog);
711 for (i = 0; i < sys->natoms; i++)
713 clear_rvec(state->v[i]);
716 /* Search for a frame without velocities */
718 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
722 state->natoms = fr.natoms;
724 if (sys->natoms != state->natoms)
726 gmx_fatal(FARGS, "Number of atoms in Topology "
727 "is not the same as in Trajectory");
729 copy_state(slog, &fr, bReadVel, state, &use_time);
731 /* Find the appropriate frame */
732 while ((fr_time == -1 || fr.time < fr_time) &&
733 read_next_frame(oenv, fp, &fr))
735 copy_state(slog, &fr, bReadVel, state, &use_time);
740 /* Set the relative box lengths for preserving the box shape.
741 * Note that this call can lead to differences in the last bit
742 * with respect to using gmx convert-tpr to create a [TT].tpx[tt] file.
744 set_box_rel(ir, state);
746 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
747 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
749 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
751 get_enx_state(ener, use_time, &sys->groups, ir, state);
752 preserve_box_shape(ir, state->box_rel, state->boxv);
756 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
758 int rc_scaling, int ePBC,
762 gmx_bool bFirst = TRUE, *hadAtom;
768 int natoms, npbcdim = 0;
769 char warn_buf[STRLEN], title[STRLEN];
770 int a, i, ai, j, k, mb, nat_molb;
771 gmx_molblock_t *molb;
775 get_stx_coordnum(fn, &natoms);
776 if (natoms != mtop->natoms)
778 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);
779 warning(wi, warn_buf);
783 init_t_atoms(&dumat, natoms, FALSE);
784 read_stx_conf(fn, title, &dumat, x, v, NULL, box);
786 npbcdim = ePBC2npbcdim(ePBC);
788 if (rc_scaling != erscNO)
790 copy_mat(box, invbox);
791 for (j = npbcdim; j < DIM; j++)
793 clear_rvec(invbox[j]);
796 m_inv_ur0(invbox, invbox);
799 /* Copy the reference coordinates to mtop */
803 snew(hadAtom, natoms);
804 for (mb = 0; mb < mtop->nmolblock; mb++)
806 molb = &mtop->molblock[mb];
807 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
808 pr = &(molinfo[molb->type].plist[F_POSRES]);
809 prfb = &(molinfo[molb->type].plist[F_FBPOSRES]);
810 if (pr->nr > 0 || prfb->nr > 0)
812 atom = mtop->moltype[molb->type].atoms.atom;
813 for (i = 0; (i < pr->nr); i++)
815 ai = pr->param[i].AI;
818 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
819 ai+1, *molinfo[molb->type].name, fn, natoms);
822 if (rc_scaling == erscCOM)
824 /* Determine the center of mass of the posres reference coordinates */
825 for (j = 0; j < npbcdim; j++)
827 sum[j] += atom[ai].m*x[a+ai][j];
829 totmass += atom[ai].m;
832 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
833 for (i = 0; (i < prfb->nr); i++)
835 ai = prfb->param[i].AI;
838 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
839 ai+1, *molinfo[molb->type].name, fn, natoms);
841 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE)
843 /* Determine the center of mass of the posres reference coordinates */
844 for (j = 0; j < npbcdim; j++)
846 sum[j] += atom[ai].m*x[a+ai][j];
848 totmass += atom[ai].m;
853 molb->nposres_xA = nat_molb;
854 snew(molb->posres_xA, molb->nposres_xA);
855 for (i = 0; i < nat_molb; i++)
857 copy_rvec(x[a+i], molb->posres_xA[i]);
862 molb->nposres_xB = nat_molb;
863 snew(molb->posres_xB, molb->nposres_xB);
864 for (i = 0; i < nat_molb; i++)
866 copy_rvec(x[a+i], molb->posres_xB[i]);
872 if (rc_scaling == erscCOM)
876 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
878 for (j = 0; j < npbcdim; j++)
880 com[j] = sum[j]/totmass;
882 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]);
885 if (rc_scaling != erscNO)
887 for (mb = 0; mb < mtop->nmolblock; mb++)
889 molb = &mtop->molblock[mb];
890 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
891 if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
893 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
894 for (i = 0; i < nat_molb; i++)
896 for (j = 0; j < npbcdim; j++)
898 if (rc_scaling == erscALL)
900 /* Convert from Cartesian to crystal coordinates */
901 xp[i][j] *= invbox[j][j];
902 for (k = j+1; k < npbcdim; k++)
904 xp[i][j] += invbox[k][j]*xp[i][k];
907 else if (rc_scaling == erscCOM)
909 /* Subtract the center of mass */
917 if (rc_scaling == erscCOM)
919 /* Convert the COM from Cartesian to crystal coordinates */
920 for (j = 0; j < npbcdim; j++)
922 com[j] *= invbox[j][j];
923 for (k = j+1; k < npbcdim; k++)
925 com[j] += invbox[k][j]*com[k];
931 free_t_atoms(&dumat, TRUE);
937 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
938 char *fnA, char *fnB,
939 int rc_scaling, int ePBC,
945 read_posres (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
946 if (strcmp(fnA, fnB) != 0)
948 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
952 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
953 t_inputrec *ir, warninp_t wi)
956 char warn_buf[STRLEN];
960 fprintf(stderr, "Searching the wall atom type(s)\n");
962 for (i = 0; i < ir->nwall; i++)
964 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
965 if (ir->wall_atomtype[i] == NOTSET)
967 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
968 warning_error(wi, warn_buf);
973 static int nrdf_internal(t_atoms *atoms)
978 for (i = 0; i < atoms->nr; i++)
980 /* Vsite ptype might not be set here yet, so also check the mass */
981 if ((atoms->atom[i].ptype == eptAtom ||
982 atoms->atom[i].ptype == eptNucleus)
983 && atoms->atom[i].m > 0)
990 case 0: nrdf = 0; break;
991 case 1: nrdf = 0; break;
992 case 2: nrdf = 1; break;
993 default: nrdf = nmass*3 - 6; break;
1000 spline1d( double dx,
1012 for (i = 1; i < n-1; i++)
1014 p = 0.5*y2[i-1]+2.0;
1016 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1017 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1022 for (i = n-2; i >= 0; i--)
1024 y2[i] = y2[i]*y2[i+1]+u[i];
1030 interpolate1d( double xmin,
1043 a = (xmin+(ix+1)*dx-x)/dx;
1044 b = (x-xmin-ix*dx)/dx;
1046 *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;
1047 *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];
1052 setup_cmap (int grid_spacing,
1055 gmx_cmap_t * cmap_grid)
1057 double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1059 int i, j, k, ii, jj, kk, idx;
1061 double dx, xmin, v, v1, v2, v12;
1064 snew(tmp_u, 2*grid_spacing);
1065 snew(tmp_u2, 2*grid_spacing);
1066 snew(tmp_yy, 2*grid_spacing);
1067 snew(tmp_y1, 2*grid_spacing);
1068 snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1069 snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1071 dx = 360.0/grid_spacing;
1072 xmin = -180.0-dx*grid_spacing/2;
1074 for (kk = 0; kk < nc; kk++)
1076 /* Compute an offset depending on which cmap we are using
1077 * Offset will be the map number multiplied with the
1078 * grid_spacing * grid_spacing * 2
1080 offset = kk * grid_spacing * grid_spacing * 2;
1082 for (i = 0; i < 2*grid_spacing; i++)
1084 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1086 for (j = 0; j < 2*grid_spacing; j++)
1088 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1089 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1093 for (i = 0; i < 2*grid_spacing; i++)
1095 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1098 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1100 ii = i-grid_spacing/2;
1103 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1105 jj = j-grid_spacing/2;
1108 for (k = 0; k < 2*grid_spacing; k++)
1110 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1111 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1114 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1115 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1116 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1117 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1119 idx = ii*grid_spacing+jj;
1120 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1121 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1122 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1123 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1129 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1133 cmap_grid->ngrid = ngrid;
1134 cmap_grid->grid_spacing = grid_spacing;
1135 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1137 snew(cmap_grid->cmapdata, ngrid);
1139 for (i = 0; i < cmap_grid->ngrid; i++)
1141 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1146 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1148 int count, count_mol, i, mb;
1149 gmx_molblock_t *molb;
1154 for (mb = 0; mb < mtop->nmolblock; mb++)
1157 molb = &mtop->molblock[mb];
1158 plist = mi[molb->type].plist;
1160 for (i = 0; i < F_NRE; i++)
1164 count_mol += 3*plist[i].nr;
1166 else if (interaction_function[i].flags & IF_CONSTRAINT)
1168 count_mol += plist[i].nr;
1172 if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1175 "Molecule type '%s' has %d constraints.\n"
1176 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1177 *mi[molb->type].name, count_mol,
1178 nrdf_internal(&mi[molb->type].atoms));
1181 count += molb->nmol*count_mol;
1187 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1189 int i, nmiss, natoms, mt;
1191 const t_atoms *atoms;
1194 for (mt = 0; mt < sys->nmoltype; mt++)
1196 atoms = &sys->moltype[mt].atoms;
1199 for (i = 0; i < natoms; i++)
1201 q = atoms->atom[i].q;
1202 if ((get_atomtype_radius(atoms->atom[i].type, atype) == 0 ||
1203 get_atomtype_vol(atoms->atom[i].type, atype) == 0 ||
1204 get_atomtype_surftens(atoms->atom[i].type, atype) == 0 ||
1205 get_atomtype_gb_radius(atoms->atom[i].type, atype) == 0 ||
1206 get_atomtype_S_hct(atoms->atom[i].type, atype) == 0) &&
1209 fprintf(stderr, "\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1210 get_atomtype_name(atoms->atom[i].type, atype), q);
1218 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);
1223 static void check_gbsa_params(gpp_atomtype_t atype)
1227 /* If we are doing GBSA, check that we got the parameters we need
1228 * This checking is to see if there are GBSA paratmeters for all
1229 * atoms in the force field. To go around this for testing purposes
1230 * comment out the nerror++ counter temporarily
1233 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1235 if (get_atomtype_radius(i, atype) < 0 ||
1236 get_atomtype_vol(i, atype) < 0 ||
1237 get_atomtype_surftens(i, atype) < 0 ||
1238 get_atomtype_gb_radius(i, atype) < 0 ||
1239 get_atomtype_S_hct(i, atype) < 0)
1241 fprintf(stderr, "\nGB parameter(s) missing or negative for atom type '%s'\n",
1242 get_atomtype_name(i, atype));
1249 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);
1254 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1261 verletbuf_list_setup_t ls;
1264 char warn_buf[STRLEN];
1267 for (i = 0; i < ir->opts.ngtc; i++)
1269 if (ir->opts.ref_t[i] < 0)
1271 warning(wi, "Some atom groups do not use temperature coupling. This cannot be accounted for in the energy error estimation for the Verlet buffer size. The energy error and the Verlet buffer might be underestimated.");
1275 ref_T = max(ref_T, ir->opts.ref_t[i]);
1279 printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol, ref_T);
1281 for (i = 0; i < ir->opts.ngtc; i++)
1283 if (ir->opts.ref_t[i] >= 0 && ir->opts.ref_t[i] != ref_T)
1285 sprintf(warn_buf, "ref_T for group of %.1f DOFs is %g K, which is smaller than the maximum of %g K used for the buffer size calculation. The buffer size might be on the conservative (large) side.",
1286 ir->opts.nrdf[i], ir->opts.ref_t[i], ref_T);
1287 warning_note(wi, warn_buf);
1291 /* Calculate the buffer size for simple atom vs atoms list */
1292 ls.cluster_size_i = 1;
1293 ls.cluster_size_j = 1;
1294 calc_verlet_buffer_size(mtop, det(box), ir,
1295 &ls, &n_nonlin_vsite, &rlist_1x1);
1297 /* Set the pair-list buffer size in ir */
1298 verletbuf_get_list_setup(FALSE, &ls);
1299 calc_verlet_buffer_size(mtop, det(box), ir,
1300 &ls, &n_nonlin_vsite, &ir->rlist);
1302 if (n_nonlin_vsite > 0)
1304 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);
1305 warning_note(wi, warn_buf);
1308 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1309 1, 1, rlist_1x1, rlist_1x1-max(ir->rvdw, ir->rcoulomb));
1311 ir->rlistlong = ir->rlist;
1312 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1313 ls.cluster_size_i, ls.cluster_size_j,
1314 ir->rlist, ir->rlist-max(ir->rvdw, ir->rcoulomb));
1316 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
1318 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)));
1322 int gmx_grompp(int argc, char *argv[])
1324 static const char *desc[] = {
1325 "[THISMODULE] (the gromacs preprocessor)",
1326 "reads a molecular topology file, checks the validity of the",
1327 "file, expands the topology from a molecular description to an atomic",
1328 "description. The topology file contains information about",
1329 "molecule types and the number of molecules, the preprocessor",
1330 "copies each molecule as needed. ",
1331 "There is no limitation on the number of molecule types. ",
1332 "Bonds and bond-angles can be converted into constraints, separately",
1333 "for hydrogens and heavy atoms.",
1334 "Then a coordinate file is read and velocities can be generated",
1335 "from a Maxwellian distribution if requested.",
1336 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1337 "(eg. number of MD steps, time step, cut-off), and others such as",
1338 "NEMD parameters, which are corrected so that the net acceleration",
1340 "Eventually a binary file is produced that can serve as the sole input",
1341 "file for the MD program.[PAR]",
1343 "[THISMODULE] uses the atom names from the topology file. The atom names",
1344 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1345 "warnings when they do not match the atom names in the topology.",
1346 "Note that the atom names are irrelevant for the simulation as",
1347 "only the atom types are used for generating interaction parameters.[PAR]",
1349 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1350 "etc. The preprocessor supports the following keywords:[PAR]",
1351 "#ifdef VARIABLE[BR]",
1352 "#ifndef VARIABLE[BR]",
1355 "#define VARIABLE[BR]",
1356 "#undef VARIABLE[BR]"
1357 "#include \"filename\"[BR]",
1358 "#include <filename>[PAR]",
1359 "The functioning of these statements in your topology may be modulated by",
1360 "using the following two flags in your [TT].mdp[tt] file:[PAR]",
1361 "[TT]define = -DVARIABLE1 -DVARIABLE2[BR]",
1362 "include = -I/home/john/doe[tt][BR]",
1363 "For further information a C-programming textbook may help you out.",
1364 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1365 "topology file written out so that you can verify its contents.[PAR]",
1367 /* cpp has been unnecessary for some time, hasn't it?
1368 "If your system does not have a C-preprocessor, you can still",
1369 "use [TT]grompp[tt], but you do not have access to the features ",
1370 "from the cpp. Command line options to the C-preprocessor can be given",
1371 "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1374 "When using position restraints a file with restraint coordinates",
1375 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1376 "with respect to the conformation from the [TT]-c[tt] option.",
1377 "For free energy calculation the the coordinates for the B topology",
1378 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1379 "those of the A topology.[PAR]",
1381 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1382 "The last frame with coordinates and velocities will be read,",
1383 "unless the [TT]-time[tt] option is used. Only if this information",
1384 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1385 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1386 "in your [TT].mdp[tt] file. An energy file can be supplied with",
1387 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1390 "[THISMODULE] can be used to restart simulations (preserving",
1391 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1392 "However, for simply changing the number of run steps to extend",
1393 "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1394 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1395 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1396 "like output frequency, then supplying the checkpoint file to",
1397 "[THISMODULE] with [TT]-t[tt] along with a new [TT].mdp[tt] file",
1398 "with [TT]-f[tt] is the recommended procedure.[PAR]",
1400 "By default, all bonded interactions which have constant energy due to",
1401 "virtual site constructions will be removed. If this constant energy is",
1402 "not zero, this will result in a shift in the total energy. All bonded",
1403 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1404 "all constraints for distances which will be constant anyway because",
1405 "of virtual site constructions will be removed. If any constraints remain",
1406 "which involve virtual sites, a fatal error will result.[PAR]"
1408 "To verify your run input file, please take note of all warnings",
1409 "on the screen, and correct where necessary. Do also look at the contents",
1410 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1411 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1412 "with the [TT]-debug[tt] option which will give you more information",
1413 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1414 "can see the contents of the run input file with the [gmx-dump]",
1415 "program. [gmx-check] can be used to compare the contents of two",
1416 "run input files.[PAR]"
1418 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1419 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1420 "harmless, but usually they are not. The user is advised to carefully",
1421 "interpret the output messages before attempting to bypass them with",
1428 gpp_atomtype_t atype;
1430 int natoms, nvsite, comb, mt;
1434 real max_spacing, fudgeQQ;
1436 char fn[STRLEN], fnB[STRLEN];
1437 const char *mdparin;
1439 gmx_bool bNeedVel, bGenVel;
1440 gmx_bool have_atomnumber;
1442 t_params *gb_plist = NULL;
1443 gmx_genborn_t *born = NULL;
1445 gmx_bool bVerbose = FALSE;
1447 char warn_buf[STRLEN];
1450 { efMDP, NULL, NULL, ffREAD },
1451 { efMDP, "-po", "mdout", ffWRITE },
1452 { efSTX, "-c", NULL, ffREAD },
1453 { efSTX, "-r", NULL, ffOPTRD },
1454 { efSTX, "-rb", NULL, ffOPTRD },
1455 { efNDX, NULL, NULL, ffOPTRD },
1456 { efTOP, NULL, NULL, ffREAD },
1457 { efTOP, "-pp", "processed", ffOPTWR },
1458 { efTPX, "-o", NULL, ffWRITE },
1459 { efTRN, "-t", NULL, ffOPTRD },
1460 { efEDR, "-e", NULL, ffOPTRD },
1461 { efTRN, "-ref", "rotref", ffOPTRW }
1463 #define NFILE asize(fnm)
1465 /* Command line options */
1466 static gmx_bool bRenum = TRUE;
1467 static gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1468 static int i, maxwarn = 0;
1469 static real fr_time = -1;
1471 { "-v", FALSE, etBOOL, {&bVerbose},
1472 "Be loud and noisy" },
1473 { "-time", FALSE, etREAL, {&fr_time},
1474 "Take frame at or first after this time." },
1475 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1476 "Remove constant bonded interactions with virtual sites" },
1477 { "-maxwarn", FALSE, etINT, {&maxwarn},
1478 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1479 { "-zero", FALSE, etBOOL, {&bZero},
1480 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1481 { "-renum", FALSE, etBOOL, {&bRenum},
1482 "Renumber atomtypes and minimize number of atomtypes" }
1485 /* Initiate some variables */
1490 /* Parse the command line */
1491 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1492 asize(desc), desc, 0, NULL, &oenv))
1497 wi = init_warning(TRUE, maxwarn);
1499 /* PARAMETER file processing */
1500 mdparin = opt2fn("-f", NFILE, fnm);
1501 set_warning_line(wi, mdparin, -1);
1502 get_ir(mdparin, opt2fn("-po", NFILE, fnm), ir, opts, wi);
1506 fprintf(stderr, "checking input for internal consistency...\n");
1508 check_ir(mdparin, ir, opts, wi);
1510 if (ir->ld_seed == -1)
1512 ir->ld_seed = make_seed();
1513 fprintf(stderr, "Setting the LD random seed to %d\n", ir->ld_seed);
1516 if (ir->expandedvals->lmc_seed == -1)
1518 ir->expandedvals->lmc_seed = make_seed();
1519 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1522 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1523 bGenVel = (bNeedVel && opts->bGenVel);
1524 if (bGenVel && ir->bContinuation)
1527 "Generating velocities is inconsistent with attempting "
1528 "to continue a previous run. Choose only one of "
1529 "gen-vel = yes and continuation = yes.");
1530 warning_error(wi, warn_buf);
1536 atype = init_atomtype();
1539 pr_symtab(debug, 0, "Just opened", &sys->symtab);
1542 strcpy(fn, ftp2fn(efTOP, NFILE, fnm));
1543 if (!gmx_fexist(fn))
1545 gmx_fatal(FARGS, "%s does not exist", fn);
1547 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1548 opts, ir, bZero, bGenVel, bVerbose, &state,
1549 atype, sys, &nmi, &mi, plist, &comb, &reppow, &fudgeQQ,
1555 pr_symtab(debug, 0, "After new_status", &sys->symtab);
1559 /* set parameters for virtual site construction (not for vsiten) */
1560 for (mt = 0; mt < sys->nmoltype; mt++)
1563 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1565 /* now throw away all obsolete bonds, angles and dihedrals: */
1566 /* note: constraints are ALWAYS removed */
1569 for (mt = 0; mt < sys->nmoltype; mt++)
1571 clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1575 if (ir->cutoff_scheme == ecutsVERLET)
1577 fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1578 ecutscheme_names[ir->cutoff_scheme]);
1580 /* Remove all charge groups */
1581 gmx_mtop_remove_chargegroups(sys);
1583 if (EVDW_PME(ir->vdwtype))
1585 gmx_fatal(FARGS, "LJ-PME not implemented together with verlet-scheme!");
1589 if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1591 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1593 sprintf(warn_buf, "Can not do %s with %s, use %s",
1594 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1595 warning_error(wi, warn_buf);
1597 if (ir->bPeriodicMols)
1599 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1600 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1601 warning_error(wi, warn_buf);
1605 if (EI_SD (ir->eI) && ir->etc != etcNO)
1607 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1610 /* If we are doing QM/MM, check that we got the atom numbers */
1611 have_atomnumber = TRUE;
1612 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1614 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1616 if (!have_atomnumber && ir->bQMMM)
1620 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1621 "field you are using does not contain atom numbers fields. This is an\n"
1622 "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1623 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1624 "an integer just before the mass column in ffXXXnb.itp.\n"
1625 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1630 if ((ir->adress->const_wf > 1) || (ir->adress->const_wf < 0))
1632 warning_error(wi, "AdResS contant weighting function should be between 0 and 1\n\n");
1634 /** \TODO check size of ex+hy width against box size */
1637 /* Check for errors in the input now, since they might cause problems
1638 * during processing further down.
1640 check_warning_error(wi, FARGS);
1642 if (opt2bSet("-r", NFILE, fnm))
1644 sprintf(fn, "%s", opt2fn("-r", NFILE, fnm));
1648 sprintf(fn, "%s", opt2fn("-c", NFILE, fnm));
1650 if (opt2bSet("-rb", NFILE, fnm))
1652 sprintf(fnB, "%s", opt2fn("-rb", NFILE, fnm));
1659 if (nint_ftype(sys, mi, F_POSRES) > 0 || nint_ftype(sys, mi, F_FBPOSRES) > 0)
1663 fprintf(stderr, "Reading position restraint coords from %s", fn);
1664 if (strcmp(fn, fnB) == 0)
1666 fprintf(stderr, "\n");
1670 fprintf(stderr, " and %s\n", fnB);
1673 gen_posres(sys, mi, fn, fnB,
1674 ir->refcoord_scaling, ir->ePBC,
1675 ir->posres_com, ir->posres_comB,
1679 /* If we are using CMAP, setup the pre-interpolation grid */
1680 if (plist->ncmap > 0)
1682 init_cmap_grid(&sys->ffparams.cmap_grid, plist->nc, plist->grid_spacing);
1683 setup_cmap(plist->grid_spacing, plist->nc, plist->cmap, &sys->ffparams.cmap_grid);
1686 set_wall_atomtype(atype, opts, ir, wi);
1689 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1690 ntype = get_atomtype_ntypes(atype);
1693 if (ir->implicit_solvent != eisNO)
1695 /* Now we have renumbered the atom types, we can check the GBSA params */
1696 check_gbsa_params(atype);
1698 /* Check that all atoms that have charge and/or LJ-parameters also have
1699 * sensible GB-parameters
1701 check_gbsa_params_charged(sys, atype);
1704 /* PELA: Copy the atomtype data to the topology atomtype list */
1705 copy_atomtype_atomtypes(atype, &(sys->atomtypes));
1709 pr_symtab(debug, 0, "After renum_atype", &sys->symtab);
1714 fprintf(stderr, "converting bonded parameters...\n");
1717 ntype = get_atomtype_ntypes(atype);
1718 convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1722 pr_symtab(debug, 0, "After convert_params", &sys->symtab);
1725 /* set ptype to VSite for virtual sites */
1726 for (mt = 0; mt < sys->nmoltype; mt++)
1728 set_vsites_ptype(FALSE, &sys->moltype[mt]);
1732 pr_symtab(debug, 0, "After virtual sites", &sys->symtab);
1734 /* Check velocity for virtual sites and shells */
1737 check_vel(sys, state.v);
1743 for (i = 0; i < sys->nmoltype; i++)
1745 check_cg_sizes(ftp2fn(efTOP, NFILE, fnm), &sys->moltype[i].cgs, wi);
1748 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1750 check_bonds_timestep(sys, ir->delta_t, wi);
1753 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
1755 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.");
1758 check_warning_error(wi, FARGS);
1762 fprintf(stderr, "initialising group options...\n");
1764 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
1766 bGenVel ? state.v : NULL,
1769 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0 &&
1772 if (EI_DYNAMICS(ir->eI) &&
1773 !(EI_MD(ir->eI) && ir->etc == etcNO) &&
1774 inputrec2nboundeddim(ir) == 3)
1776 set_verlet_buffer(sys, ir, state.box, wi);
1780 /* Init the temperature coupling state */
1781 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
1785 fprintf(stderr, "Checking consistency between energy and charge groups...\n");
1787 check_eg_vs_cg(sys);
1791 pr_symtab(debug, 0, "After index", &sys->symtab);
1794 triple_check(mdparin, ir, sys, wi);
1795 close_symtab(&sys->symtab);
1798 pr_symtab(debug, 0, "After close", &sys->symtab);
1801 /* make exclusions between QM atoms */
1804 if (ir->QMMMscheme == eQMMMschemenormal && ir->ns_type == ensSIMPLE)
1806 gmx_fatal(FARGS, "electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
1810 generate_qmexcl(sys, ir, wi);
1814 if (ftp2bSet(efTRN, NFILE, fnm))
1818 fprintf(stderr, "getting data from old trajectory ...\n");
1820 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
1821 bNeedVel, bGenVel, fr_time, ir, &state, sys, oenv);
1824 if (ir->ePBC == epbcXY && ir->nwall != 2)
1826 clear_rvec(state.box[ZZ]);
1829 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
1831 set_warning_line(wi, mdparin, -1);
1832 check_chargegroup_radii(sys, ir, state.x, wi);
1835 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
1837 /* Calculate the optimal grid dimensions */
1838 copy_mat(state.box, box);
1839 if (ir->ePBC == epbcXY && ir->nwall == 2)
1841 svmul(ir->wall_ewald_zfac, box[ZZ], box[ZZ]);
1843 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
1845 /* Mark fourier_spacing as not used */
1846 ir->fourier_spacing = 0;
1848 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
1850 set_warning_line(wi, mdparin, -1);
1851 warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
1853 max_spacing = calc_grid(stdout, box, ir->fourier_spacing,
1854 &(ir->nkx), &(ir->nky), &(ir->nkz));
1857 /* MRS: eventually figure out better logic for initializing the fep
1858 values that makes declaring the lambda and declaring the state not
1859 potentially conflict if not handled correctly. */
1860 if (ir->efep != efepNO)
1862 if (EVDW_PME(ir->vdwtype))
1864 gmx_fatal(FARGS, "LJ-PME not implemented together with free energy calculations!");
1867 state.fep_state = ir->fepvals->init_fep_state;
1868 for (i = 0; i < efptNR; i++)
1870 /* init_lambda trumps state definitions*/
1871 if (ir->fepvals->init_lambda >= 0)
1873 state.lambda[i] = ir->fepvals->init_lambda;
1877 if (ir->fepvals->all_lambda[i] == NULL)
1879 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
1883 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
1889 if (ir->ePull != epullNO)
1891 set_pull_init(ir, sys, state.x, state.box, state.lambda[efptMASS], oenv, opts->pull_start);
1896 set_reference_positions(ir->rot, state.x, state.box,
1897 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
1901 /* reset_multinr(sys); */
1903 if (EEL_PME(ir->coulombtype))
1905 float ratio = pme_load_estimate(sys, ir, state.box);
1906 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
1907 /* With free energy we might need to do PME both for the A and B state
1908 * charges. This will double the cost, but the optimal performance will
1909 * then probably be at a slightly larger cut-off and grid spacing.
1911 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
1912 (ir->efep != efepNO && ratio > 2.0/3.0))
1915 "The optimal PME mesh load for parallel simulations is below 0.5\n"
1916 "and for highly parallel simulations between 0.25 and 0.33,\n"
1917 "for higher performance, increase the cut-off and the PME grid spacing.\n");
1918 if (ir->efep != efepNO)
1921 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
1927 char warn_buf[STRLEN];
1928 double cio = compute_io(ir, sys->natoms, &sys->groups, F_NRE, 1);
1929 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
1932 set_warning_line(wi, mdparin, -1);
1933 warning_note(wi, warn_buf);
1937 printf("%s\n", warn_buf);
1943 fprintf(stderr, "writing run input file...\n");
1946 done_warning(wi, FARGS);
1947 write_tpx_state(ftp2fn(efTPX, NFILE, fnm), ir, &state, sys);
1949 done_atomtype(atype);
1950 done_mtop(sys, TRUE);
1951 done_inputrec_strings();