Merge branch 'master' of git://git.gromacs.org/gromacs
[gromacs/adressmacs.git] / src / kernel / grompp.c
blob410553cdadc210c73d2be58498f4032911d139bd
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.03
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include <sys/types.h>
41 #include <math.h>
42 #include <string.h>
43 #include <errno.h>
44 #include <limits.h>
46 #include "sysstuff.h"
47 #include "smalloc.h"
48 #include "macros.h"
49 #include "string2.h"
50 #include "readir.h"
51 #include "toputil.h"
52 #include "topio.h"
53 #include "confio.h"
54 #include "copyrite.h"
55 #include "readir.h"
56 #include "symtab.h"
57 #include "names.h"
58 #include "grompp.h"
59 #include "random.h"
60 #include "vec.h"
61 #include "futil.h"
62 #include "statutil.h"
63 #include "splitter.h"
64 #include "sortwater.h"
65 #include "convparm.h"
66 #include "gmx_fatal.h"
67 #include "warninp.h"
68 #include "index.h"
69 #include "gmxfio.h"
70 #include "trnio.h"
71 #include "tpxio.h"
72 #include "vsite_parm.h"
73 #include "txtdump.h"
74 #include "calcgrid.h"
75 #include "add_par.h"
76 #include "enxio.h"
77 #include "perf_est.h"
78 #include "compute_io.h"
79 #include "gpp_atomtype.h"
80 #include "gpp_tomorse.h"
81 #include "mtop_util.h"
82 #include "genborn.h"
84 static int rm_interactions(int ifunc,int nrmols,t_molinfo mols[])
86 int i,n;
88 n=0;
89 /* For all the molecule types */
90 for(i=0; i<nrmols; i++) {
91 n += mols[i].plist[ifunc].nr;
92 mols[i].plist[ifunc].nr=0;
94 return n;
97 static int check_atom_names(const char *fn1, const char *fn2,
98 gmx_mtop_t *mtop, t_atoms *at)
100 int mb,m,i,j,nmismatch;
101 t_atoms *tat;
102 #define MAXMISMATCH 20
104 if (mtop->natoms != at->nr)
105 gmx_incons("comparing atom names");
107 nmismatch=0;
108 i = 0;
109 for(mb=0; mb<mtop->nmolblock; mb++) {
110 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
111 for(m=0; m<mtop->molblock[mb].nmol; m++) {
112 for(j=0; j < tat->nr; j++) {
113 if (strcmp( *(tat->atomname[j]) , *(at->atomname[i]) ) != 0) {
114 if (nmismatch < MAXMISMATCH) {
115 fprintf(stderr,
116 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
117 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
118 } else if (nmismatch == MAXMISMATCH) {
119 fprintf(stderr,"(more than %d non-matching atom names)\n",MAXMISMATCH);
121 nmismatch++;
123 i++;
128 return nmismatch;
131 static void check_eg_vs_cg(gmx_mtop_t *mtop)
133 int astart,mb,m,cg,j,firstj;
134 unsigned char firsteg,eg;
135 gmx_moltype_t *molt;
137 /* Go through all the charge groups and make sure all their
138 * atoms are in the same energy group.
141 astart = 0;
142 for(mb=0; mb<mtop->nmolblock; mb++) {
143 molt = &mtop->moltype[mtop->molblock[mb].type];
144 for(m=0; m<mtop->molblock[mb].nmol; m++) {
145 for(cg=0; cg<molt->cgs.nr;cg++) {
146 /* Get the energy group of the first atom in this charge group */
147 firstj = astart + molt->cgs.index[cg];
148 firsteg = ggrpnr(&mtop->groups,egcENER,firstj);
149 for(j=molt->cgs.index[cg]+1;j<molt->cgs.index[cg+1];j++) {
150 eg = ggrpnr(&mtop->groups,egcENER,astart+j);
151 if(eg != firsteg) {
152 gmx_fatal(FARGS,"atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
153 firstj+1,astart+j+1,cg+1,*molt->name);
157 astart += molt->atoms.nr;
162 static void check_cg_sizes(const char *topfn,t_block *cgs,warninp_t wi)
164 int maxsize,cg;
165 char warn_buf[STRLEN];
167 maxsize = 0;
168 for(cg=0; cg<cgs->nr; cg++)
170 maxsize = max(maxsize,cgs->index[cg+1]-cgs->index[cg]);
173 if (maxsize > MAX_CHARGEGROUP_SIZE)
175 gmx_fatal(FARGS,"The largest charge group contains %d atoms. The maximum is %d.",maxsize,MAX_CHARGEGROUP_SIZE);
177 else if (maxsize > 10)
179 set_warning_line(wi,topfn,-1);
180 sprintf(warn_buf,
181 "The largest charge group contains %d atoms.\n"
182 "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"
183 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
184 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
185 maxsize);
186 warning_note(wi,warn_buf);
190 static void check_bonds_timestep(gmx_mtop_t *mtop,double dt,warninp_t wi)
192 /* This check is not intended to ensure accurate integration,
193 * rather it is to signal mistakes in the mdp settings.
194 * A common mistake is to forget to turn on constraints
195 * for MD after energy minimization with flexible bonds.
196 * This check can also detect too large time steps for flexible water
197 * models, but such errors will often be masked by the constraints
198 * mdp options, which turns flexible water into water with bond constraints,
199 * but without an angle constraint. Unfortunately such incorrect use
200 * of water models can not easily be detected without checking
201 * for specific model names.
203 * The stability limit of leap-frog or velocity verlet is 4.44 steps
204 * per oscillational period.
205 * But accurate bonds distributions are lost far before that limit.
206 * To allow relatively common schemes (although not common with Gromacs)
207 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
208 * we set the note limit to 10.
210 int min_steps_warn=5;
211 int min_steps_note=10;
212 t_iparams *ip;
213 int molt;
214 gmx_moltype_t *moltype,*w_moltype;
215 t_atom *atom;
216 t_ilist *ilist,*ilb,*ilc,*ils;
217 int ftype;
218 int i,a1,a2,w_a1,w_a2,j;
219 real twopi2,limit2,fc,re,m1,m2,period2,w_period2;
220 gmx_bool bFound,bWater,bWarn;
221 char warn_buf[STRLEN];
223 ip = mtop->ffparams.iparams;
225 twopi2 = sqr(2*M_PI);
227 limit2 = sqr(min_steps_note*dt);
229 w_a1 = w_a2 = -1;
230 w_period2 = -1.0;
232 w_moltype = NULL;
233 for(molt=0; molt<mtop->nmoltype; molt++)
235 moltype = &mtop->moltype[molt];
236 atom = moltype->atoms.atom;
237 ilist = moltype->ilist;
238 ilc = &ilist[F_CONSTR];
239 ils = &ilist[F_SETTLE];
240 for(ftype=0; ftype<F_NRE; ftype++)
242 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
244 continue;
247 ilb = &ilist[ftype];
248 for(i=0; i<ilb->nr; i+=3)
250 fc = ip[ilb->iatoms[i]].harmonic.krA;
251 re = ip[ilb->iatoms[i]].harmonic.rA;
252 if (ftype == F_G96BONDS)
254 /* Convert squared sqaure fc to harmonic fc */
255 fc = 2*fc*re;
257 a1 = ilb->iatoms[i+1];
258 a2 = ilb->iatoms[i+2];
259 m1 = atom[a1].m;
260 m2 = atom[a2].m;
261 if (fc > 0 && m1 > 0 && m2 > 0)
263 period2 = twopi2*m1*m2/((m1 + m2)*fc);
265 else
267 period2 = GMX_FLOAT_MAX;
269 if (debug)
271 fprintf(debug,"fc %g m1 %g m2 %g period %g\n",
272 fc,m1,m2,sqrt(period2));
274 if (period2 < limit2)
276 bFound = FALSE;
277 for(j=0; j<ilc->nr; j+=3)
279 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
280 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
282 bFound = TRUE;
285 for(j=0; j<ils->nr; j+=2)
287 if ((a1 >= ils->iatoms[j+1] && a1 < ils->iatoms[j+1]+3) &&
288 (a2 >= ils->iatoms[j+1] && a2 < ils->iatoms[j+1]+3))
290 bFound = TRUE;
293 if (!bFound &&
294 (w_moltype == NULL || period2 < w_period2))
296 w_moltype = moltype;
297 w_a1 = a1;
298 w_a2 = a2;
299 w_period2 = period2;
306 if (w_moltype != NULL)
308 bWarn = (w_period2 < sqr(min_steps_warn*dt));
309 /* A check that would recognize most water models */
310 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
311 w_moltype->atoms.nr <= 5);
312 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"
313 "%s",
314 *w_moltype->name,
315 w_a1+1,*w_moltype->atoms.atomname[w_a1],
316 w_a2+1,*w_moltype->atoms.atomname[w_a2],
317 sqrt(w_period2),bWarn ? min_steps_warn : min_steps_note,dt,
318 bWater ?
319 "Maybe you asked for fexible water." :
320 "Maybe you forgot to change the constraints mdp option.");
321 if (bWarn)
323 warning(wi,warn_buf);
325 else
327 warning_note(wi,warn_buf);
332 static void check_vel(gmx_mtop_t *mtop,rvec v[])
334 gmx_mtop_atomloop_all_t aloop;
335 t_atom *atom;
336 int a;
338 aloop = gmx_mtop_atomloop_all_init(mtop);
339 while (gmx_mtop_atomloop_all_next(aloop,&a,&atom)) {
340 if (atom->ptype == eptShell ||
341 atom->ptype == eptBond ||
342 atom->ptype == eptVSite) {
343 clear_rvec(v[a]);
348 static gmx_bool nint_ftype(gmx_mtop_t *mtop,t_molinfo *mi,int ftype)
350 int nint,mb;
352 nint = 0;
353 for(mb=0; mb<mtop->nmolblock; mb++) {
354 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
357 return nint;
360 /* This routine reorders the molecule type array
361 * in the order of use in the molblocks,
362 * unused molecule types are deleted.
364 static void renumber_moltypes(gmx_mtop_t *sys,
365 int *nmolinfo,t_molinfo **molinfo)
367 int *order,norder,i;
368 int mb,mi;
369 t_molinfo *minew;
371 snew(order,*nmolinfo);
372 norder = 0;
373 for(mb=0; mb<sys->nmolblock; mb++) {
374 for(i=0; i<norder; i++) {
375 if (order[i] == sys->molblock[mb].type) {
376 break;
379 if (i == norder) {
380 /* This type did not occur yet, add it */
381 order[norder] = sys->molblock[mb].type;
382 /* Renumber the moltype in the topology */
383 norder++;
385 sys->molblock[mb].type = i;
388 /* We still need to reorder the molinfo structs */
389 snew(minew,norder);
390 for(mi=0; mi<*nmolinfo; mi++) {
391 for(i=0; i<norder; i++) {
392 if (order[i] == mi) {
393 break;
396 if (i == norder) {
397 done_mi(&(*molinfo)[mi]);
398 } else {
399 minew[i] = (*molinfo)[mi];
402 sfree(*molinfo);
404 *nmolinfo = norder;
405 *molinfo = minew;
408 static void molinfo2mtop(int nmi,t_molinfo *mi,gmx_mtop_t *mtop)
410 int m;
411 gmx_moltype_t *molt;
413 mtop->nmoltype = nmi;
414 snew(mtop->moltype,nmi);
415 for(m=0; m<nmi; m++) {
416 molt = &mtop->moltype[m];
417 molt->name = mi[m].name;
418 molt->atoms = mi[m].atoms;
419 /* ilists are copied later */
420 molt->cgs = mi[m].cgs;
421 molt->excls = mi[m].excls;
425 static void
426 new_status(const char *topfile,const char *topppfile,const char *confin,
427 t_gromppopts *opts,t_inputrec *ir,gmx_bool bZero,
428 gmx_bool bGenVel,gmx_bool bVerbose,t_state *state,
429 gpp_atomtype_t atype,gmx_mtop_t *sys,
430 int *nmi,t_molinfo **mi,t_params plist[],
431 int *comb,double *reppow,real *fudgeQQ,
432 gmx_bool bMorse,
433 warninp_t wi)
435 t_molinfo *molinfo=NULL;
436 int nmolblock;
437 gmx_molblock_t *molblock,*molbs;
438 t_atoms *confat;
439 int mb,i,nrmols,nmismatch;
440 char buf[STRLEN];
441 gmx_bool bGB=FALSE;
442 char warn_buf[STRLEN];
444 init_mtop(sys);
446 /* Set gmx_boolean for GB */
447 if(ir->implicit_solvent)
448 bGB=TRUE;
450 /* TOPOLOGY processing */
451 sys->name = do_top(bVerbose,topfile,topppfile,opts,bZero,&(sys->symtab),
452 plist,comb,reppow,fudgeQQ,
453 atype,&nrmols,&molinfo,ir,
454 &nmolblock,&molblock,bGB,
455 wi);
457 sys->nmolblock = 0;
458 snew(sys->molblock,nmolblock);
460 sys->natoms = 0;
461 for(mb=0; mb<nmolblock; mb++) {
462 if (sys->nmolblock > 0 &&
463 molblock[mb].type == sys->molblock[sys->nmolblock-1].type) {
464 /* Merge consecutive blocks with the same molecule type */
465 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
466 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
467 } else if (molblock[mb].nmol > 0) {
468 /* Add a new molblock to the topology */
469 molbs = &sys->molblock[sys->nmolblock];
470 *molbs = molblock[mb];
471 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
472 molbs->nposres_xA = 0;
473 molbs->nposres_xB = 0;
474 sys->natoms += molbs->nmol*molbs->natoms_mol;
475 sys->nmolblock++;
478 if (sys->nmolblock == 0) {
479 gmx_fatal(FARGS,"No molecules were defined in the system");
482 renumber_moltypes(sys,&nrmols,&molinfo);
484 if (bMorse)
485 convert_harmonics(nrmols,molinfo,atype);
487 if (ir->eDisre == edrNone) {
488 i = rm_interactions(F_DISRES,nrmols,molinfo);
489 if (i > 0) {
490 set_warning_line(wi,"unknown",-1);
491 sprintf(warn_buf,"disre = no, removed %d distance restraints",i);
492 warning_note(wi,warn_buf);
495 if (opts->bOrire == FALSE) {
496 i = rm_interactions(F_ORIRES,nrmols,molinfo);
497 if (i > 0) {
498 set_warning_line(wi,"unknown",-1);
499 sprintf(warn_buf,"orire = no, removed %d orientation restraints",i);
500 warning_note(wi,warn_buf);
503 if (opts->bDihre == FALSE) {
504 i = rm_interactions(F_DIHRES,nrmols,molinfo);
505 if (i > 0) {
506 set_warning_line(wi,"unknown",-1);
507 sprintf(warn_buf,"dihre = no, removed %d dihedral restraints",i);
508 warning_note(wi,warn_buf);
512 /* Copy structures from msys to sys */
513 molinfo2mtop(nrmols,molinfo,sys);
515 gmx_mtop_finalize(sys);
517 /* COORDINATE file processing */
518 if (bVerbose)
519 fprintf(stderr,"processing coordinates...\n");
521 get_stx_coordnum(confin,&state->natoms);
522 if (state->natoms != sys->natoms)
523 gmx_fatal(FARGS,"number of coordinates in coordinate file (%s, %d)\n"
524 " does not match topology (%s, %d)",
525 confin,state->natoms,topfile,sys->natoms);
526 else {
527 /* make space for coordinates and velocities */
528 char title[STRLEN];
529 snew(confat,1);
530 init_t_atoms(confat,state->natoms,FALSE);
531 init_state(state,state->natoms,0,0,0);
532 read_stx_conf(confin,title,confat,state->x,state->v,NULL,state->box);
533 /* This call fixes the box shape for runs with pressure scaling */
534 set_box_rel(ir,state);
536 nmismatch = check_atom_names(topfile, confin, sys, confat);
537 free_t_atoms(confat,TRUE);
538 sfree(confat);
540 if (nmismatch) {
541 sprintf(buf,"%d non-matching atom name%s\n"
542 "atom names from %s will be used\n"
543 "atom names from %s will be ignored\n",
544 nmismatch,(nmismatch == 1) ? "" : "s",topfile,confin);
545 warning(wi,buf);
547 if (bVerbose)
548 fprintf(stderr,"double-checking input for internal consistency...\n");
549 double_check(ir,state->box,nint_ftype(sys,molinfo,F_CONSTR),wi);
552 if (bGenVel) {
553 real *mass;
554 gmx_mtop_atomloop_all_t aloop;
555 t_atom *atom;
557 snew(mass,state->natoms);
558 aloop = gmx_mtop_atomloop_all_init(sys);
559 while (gmx_mtop_atomloop_all_next(aloop,&i,&atom)) {
560 mass[i] = atom->m;
563 if (opts->seed == -1) {
564 opts->seed = make_seed();
565 fprintf(stderr,"Setting gen_seed to %d\n",opts->seed);
567 maxwell_speed(opts->tempi,opts->seed,sys,state->v);
569 stop_cm(stdout,state->natoms,mass,state->x,state->v);
570 sfree(mass);
573 *nmi = nrmols;
574 *mi = molinfo;
577 static void copy_state(const char *slog,t_trxframe *fr,
578 gmx_bool bReadVel,t_state *state,
579 double *use_time)
581 int i;
583 if (fr->not_ok & FRAME_NOT_OK)
585 gmx_fatal(FARGS,"Can not start from an incomplete frame");
587 if (!fr->bX)
589 gmx_fatal(FARGS,"Did not find a frame with coordinates in file %s",
590 slog);
593 for(i=0; i<state->natoms; i++)
595 copy_rvec(fr->x[i],state->x[i]);
597 if (bReadVel)
599 if (!fr->bV)
601 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
603 for(i=0; i<state->natoms; i++)
605 copy_rvec(fr->v[i],state->v[i]);
608 if (fr->bBox)
610 copy_mat(fr->box,state->box);
613 *use_time = fr->time;
616 static void cont_status(const char *slog,const char *ener,
617 gmx_bool bNeedVel,gmx_bool bGenVel, real fr_time,
618 t_inputrec *ir,t_state *state,
619 gmx_mtop_t *sys,
620 const output_env_t oenv)
621 /* If fr_time == -1 read the last frame available which is complete */
623 gmx_bool bReadVel;
624 t_trxframe fr;
625 t_trxstatus *fp;
626 int i;
627 double use_time;
629 bReadVel = (bNeedVel && !bGenVel);
631 fprintf(stderr,
632 "Reading Coordinates%s and Box size from old trajectory\n",
633 bReadVel ? ", Velocities" : "");
634 if (fr_time == -1)
636 fprintf(stderr,"Will read whole trajectory\n");
638 else
640 fprintf(stderr,"Will read till time %g\n",fr_time);
642 if (!bReadVel)
644 if (bGenVel)
646 fprintf(stderr,"Velocities generated: "
647 "ignoring velocities in input trajectory\n");
649 read_first_frame(oenv,&fp,slog,&fr,TRX_NEED_X);
651 else
653 read_first_frame(oenv,&fp,slog,&fr,TRX_NEED_X | TRX_NEED_V);
655 if (!fr.bV)
657 fprintf(stderr,
658 "\n"
659 "WARNING: Did not find a frame with velocities in file %s,\n"
660 " all velocities will be set to zero!\n\n",slog);
661 for(i=0; i<sys->natoms; i++)
663 clear_rvec(state->v[i]);
665 close_trj(fp);
666 /* Search for a frame without velocities */
667 bReadVel = FALSE;
668 read_first_frame(oenv,&fp,slog,&fr,TRX_NEED_X);
672 state->natoms = fr.natoms;
674 if (sys->natoms != state->natoms)
676 gmx_fatal(FARGS,"Number of atoms in Topology "
677 "is not the same as in Trajectory");
679 copy_state(slog,&fr,bReadVel,state,&use_time);
681 /* Find the appropriate frame */
682 while ((fr_time == -1 || fr.time < fr_time) &&
683 read_next_frame(oenv,fp,&fr))
685 copy_state(slog,&fr,bReadVel,state,&use_time);
688 close_trj(fp);
690 /* Set the relative box lengths for preserving the box shape.
691 * Note that this call can lead to differences in the last bit
692 * with respect to using tpbconv to create a tpx file.
694 set_box_rel(ir,state);
696 fprintf(stderr,"Using frame at t = %g ps\n",use_time);
697 fprintf(stderr,"Starting time for run is %g ps\n",ir->init_t);
699 if ((ir->epc != epcNO || ir->etc ==etcNOSEHOOVER) && ener)
701 get_enx_state(ener,use_time,&sys->groups,ir,state);
702 preserve_box_shape(ir,state->box_rel,state->boxv);
706 static void read_posres(gmx_mtop_t *mtop,t_molinfo *molinfo,gmx_bool bTopB,
707 char *fn,
708 int rc_scaling, int ePBC,
709 rvec com,
710 warninp_t wi)
712 gmx_bool bFirst = TRUE;
713 rvec *x,*v,*xp;
714 dvec sum;
715 double totmass;
716 t_atoms dumat;
717 matrix box,invbox;
718 int natoms,npbcdim=0;
719 char warn_buf[STRLEN],title[STRLEN];
720 int a,i,ai,j,k,mb,nat_molb;
721 gmx_molblock_t *molb;
722 t_params *pr;
723 t_atom *atom;
725 get_stx_coordnum(fn,&natoms);
726 if (natoms != mtop->natoms) {
727 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);
728 warning(wi,warn_buf);
730 snew(x,natoms);
731 snew(v,natoms);
732 init_t_atoms(&dumat,natoms,FALSE);
733 read_stx_conf(fn,title,&dumat,x,v,NULL,box);
735 npbcdim = ePBC2npbcdim(ePBC);
736 clear_rvec(com);
737 if (rc_scaling != erscNO) {
738 copy_mat(box,invbox);
739 for(j=npbcdim; j<DIM; j++) {
740 clear_rvec(invbox[j]);
741 invbox[j][j] = 1;
743 m_inv_ur0(invbox,invbox);
746 /* Copy the reference coordinates to mtop */
747 clear_dvec(sum);
748 totmass = 0;
749 a = 0;
750 for(mb=0; mb<mtop->nmolblock; mb++) {
751 molb = &mtop->molblock[mb];
752 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
753 pr = &(molinfo[molb->type].plist[F_POSRES]);
754 if (pr->nr > 0) {
755 atom = mtop->moltype[molb->type].atoms.atom;
756 for(i=0; (i<pr->nr); i++) {
757 ai=pr->param[i].AI;
758 if (ai >= natoms) {
759 gmx_fatal(FARGS,"Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
760 ai+1,*molinfo[molb->type].name,fn,natoms);
762 if (rc_scaling == erscCOM) {
763 /* Determine the center of mass of the posres reference coordinates */
764 for(j=0; j<npbcdim; j++) {
765 sum[j] += atom[ai].m*x[a+ai][j];
767 totmass += atom[ai].m;
770 if (!bTopB) {
771 molb->nposres_xA = nat_molb;
772 snew(molb->posres_xA,molb->nposres_xA);
773 for(i=0; i<nat_molb; i++) {
774 copy_rvec(x[a+i],molb->posres_xA[i]);
776 } else {
777 molb->nposres_xB = nat_molb;
778 snew(molb->posres_xB,molb->nposres_xB);
779 for(i=0; i<nat_molb; i++) {
780 copy_rvec(x[a+i],molb->posres_xB[i]);
784 a += nat_molb;
786 if (rc_scaling == erscCOM) {
787 if (totmass == 0)
788 gmx_fatal(FARGS,"The total mass of the position restraint atoms is 0");
789 for(j=0; j<npbcdim; j++)
790 com[j] = sum[j]/totmass;
791 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]);
794 if (rc_scaling != erscNO) {
795 for(mb=0; mb<mtop->nmolblock; mb++) {
796 molb = &mtop->molblock[mb];
797 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
798 if (molb->nposres_xA > 0 || molb->nposres_xB > 0) {
799 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
800 for(i=0; i<nat_molb; i++) {
801 for(j=0; j<npbcdim; j++) {
802 if (rc_scaling == erscALL) {
803 /* Convert from Cartesian to crystal coordinates */
804 xp[i][j] *= invbox[j][j];
805 for(k=j+1; k<npbcdim; k++) {
806 xp[i][j] += invbox[k][j]*xp[i][k];
808 } else if (rc_scaling == erscCOM) {
809 /* Subtract the center of mass */
810 xp[i][j] -= com[j];
817 if (rc_scaling == erscCOM) {
818 /* Convert the COM from Cartesian to crystal coordinates */
819 for(j=0; j<npbcdim; j++) {
820 com[j] *= invbox[j][j];
821 for(k=j+1; k<npbcdim; k++) {
822 com[j] += invbox[k][j]*com[k];
828 free_t_atoms(&dumat,TRUE);
829 sfree(x);
830 sfree(v);
833 static void gen_posres(gmx_mtop_t *mtop,t_molinfo *mi,
834 char *fnA, char *fnB,
835 int rc_scaling, int ePBC,
836 rvec com, rvec comB,
837 warninp_t wi)
839 int i,j;
841 read_posres (mtop,mi,FALSE,fnA,rc_scaling,ePBC,com,wi);
842 if (strcmp(fnA,fnB) != 0) {
843 read_posres(mtop,mi,TRUE ,fnB,rc_scaling,ePBC,comB,wi);
847 static void set_wall_atomtype(gpp_atomtype_t at,t_gromppopts *opts,
848 t_inputrec *ir)
850 int i;
852 if (ir->nwall > 0)
853 fprintf(stderr,"Searching the wall atom type(s)\n");
854 for(i=0; i<ir->nwall; i++)
855 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i],at);
858 static int nrdf_internal(t_atoms *atoms)
860 int i,nmass,nrdf;
862 nmass = 0;
863 for(i=0; i<atoms->nr; i++) {
864 /* Vsite ptype might not be set here yet, so also check the mass */
865 if ((atoms->atom[i].ptype == eptAtom ||
866 atoms->atom[i].ptype == eptNucleus)
867 && atoms->atom[i].m > 0) {
868 nmass++;
871 switch (nmass) {
872 case 0: nrdf = 0; break;
873 case 1: nrdf = 0; break;
874 case 2: nrdf = 1; break;
875 default: nrdf = nmass*3 - 6; break;
878 return nrdf;
881 void
882 spline1d( double dx,
883 double * y,
884 int n,
885 double * u,
886 double * y2 )
888 int i;
889 double p,q;
891 y2[0] = 0.0;
892 u[0] = 0.0;
894 for(i=1;i<n-1;i++)
896 p = 0.5*y2[i-1]+2.0;
897 y2[i] = -0.5/p;
898 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
899 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
902 y2[n-1] = 0.0;
904 for(i=n-2;i>=0;i--)
906 y2[i] = y2[i]*y2[i+1]+u[i];
911 void
912 interpolate1d( double xmin,
913 double dx,
914 double * ya,
915 double * y2a,
916 double x,
917 double * y,
918 double * y1)
920 int ix;
921 double a,b;
923 ix = (x-xmin)/dx;
925 a = (xmin+(ix+1)*dx-x)/dx;
926 b = (x-xmin-ix*dx)/dx;
928 *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;
929 *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];
933 void
934 setup_cmap (int grid_spacing,
935 int nc,
936 real * grid ,
937 gmx_cmap_t * cmap_grid)
939 double *tmp_u,*tmp_u2,*tmp_yy,*tmp_y1,*tmp_t2,*tmp_grid;
941 int i,j,k,ii,jj,kk,idx;
942 int offset;
943 double dx,xmin,v,v1,v2,v12;
944 double phi,psi;
946 snew(tmp_u,2*grid_spacing);
947 snew(tmp_u2,2*grid_spacing);
948 snew(tmp_yy,2*grid_spacing);
949 snew(tmp_y1,2*grid_spacing);
950 snew(tmp_t2,2*grid_spacing*2*grid_spacing);
951 snew(tmp_grid,2*grid_spacing*2*grid_spacing);
953 dx = 360.0/grid_spacing;
954 xmin = -180.0-dx*grid_spacing/2;
956 for(kk=0;kk<nc;kk++)
958 /* Compute an offset depending on which cmap we are using
959 * Offset will be the map number multiplied with the grid_spacing * grid_spacing * 2
961 offset = kk * grid_spacing * grid_spacing * 2;
963 for(i=0;i<2*grid_spacing;i++)
965 ii=(i+grid_spacing-grid_spacing/2)%grid_spacing;
967 for(j=0;j<2*grid_spacing;j++)
969 jj=(j+grid_spacing-grid_spacing/2)%grid_spacing;
970 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
974 for(i=0;i<2*grid_spacing;i++)
976 spline1d(dx,&(tmp_grid[2*grid_spacing*i]),2*grid_spacing,tmp_u,&(tmp_t2[2*grid_spacing*i]));
979 for(i=grid_spacing/2;i<grid_spacing+grid_spacing/2;i++)
981 ii = i-grid_spacing/2;
982 phi = ii*dx-180.0;
984 for(j=grid_spacing/2;j<grid_spacing+grid_spacing/2;j++)
986 jj = j-grid_spacing/2;
987 psi = jj*dx-180.0;
989 for(k=0;k<2*grid_spacing;k++)
991 interpolate1d(xmin,dx,&(tmp_grid[2*grid_spacing*k]),
992 &(tmp_t2[2*grid_spacing*k]),psi,&tmp_yy[k],&tmp_y1[k]);
995 spline1d(dx,tmp_yy,2*grid_spacing,tmp_u,tmp_u2);
996 interpolate1d(xmin,dx,tmp_yy,tmp_u2,phi,&v,&v1);
997 spline1d(dx,tmp_y1,2*grid_spacing,tmp_u,tmp_u2);
998 interpolate1d(xmin,dx,tmp_y1,tmp_u2,phi,&v2,&v12);
1000 idx = ii*grid_spacing+jj;
1001 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1002 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1003 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1004 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1010 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1012 int i,k,nelem;
1014 cmap_grid->ngrid = ngrid;
1015 cmap_grid->grid_spacing = grid_spacing;
1016 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1018 snew(cmap_grid->cmapdata,ngrid);
1020 for(i=0;i<cmap_grid->ngrid;i++)
1022 snew(cmap_grid->cmapdata[i].cmap,4*nelem);
1027 static int count_constraints(gmx_mtop_t *mtop,t_molinfo *mi,warninp_t wi)
1029 int count,count_mol,i,mb;
1030 gmx_molblock_t *molb;
1031 t_params *plist;
1032 char buf[STRLEN];
1034 count = 0;
1035 for(mb=0; mb<mtop->nmolblock; mb++) {
1036 count_mol = 0;
1037 molb = &mtop->molblock[mb];
1038 plist = mi[molb->type].plist;
1040 for(i=0; i<F_NRE; i++) {
1041 if (i == F_SETTLE)
1042 count_mol += 3*plist[i].nr;
1043 else if (interaction_function[i].flags & IF_CONSTRAINT)
1044 count_mol += plist[i].nr;
1047 if (count_mol > nrdf_internal(&mi[molb->type].atoms)) {
1048 sprintf(buf,
1049 "Molecule type '%s' has %d constraints.\n"
1050 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1051 *mi[molb->type].name,count_mol,
1052 nrdf_internal(&mi[molb->type].atoms));
1053 warning(wi,buf);
1055 count += molb->nmol*count_mol;
1058 return count;
1061 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1063 int i,nmiss,natoms,mt;
1064 real q;
1065 const t_atoms *atoms;
1067 nmiss = 0;
1068 for(mt=0;mt<sys->nmoltype;mt++)
1070 atoms = &sys->moltype[mt].atoms;
1071 natoms = atoms->nr;
1073 for(i=0;i<natoms;i++)
1075 q = atoms->atom[i].q;
1076 if ((get_atomtype_radius(atoms->atom[i].type,atype) == 0 ||
1077 get_atomtype_vol(atoms->atom[i].type,atype) == 0 ||
1078 get_atomtype_surftens(atoms->atom[i].type,atype) == 0 ||
1079 get_atomtype_gb_radius(atoms->atom[i].type,atype) == 0 ||
1080 get_atomtype_S_hct(atoms->atom[i].type,atype) == 0) &&
1081 q != 0)
1083 fprintf(stderr,"\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1084 get_atomtype_name(atoms->atom[i].type,atype),q);
1085 nmiss++;
1090 if (nmiss > 0)
1092 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);
1097 static void check_gbsa_params(t_inputrec *ir,gpp_atomtype_t atype)
1099 int nmiss,i;
1101 /* If we are doing GBSA, check that we got the parameters we need
1102 * This checking is to see if there are GBSA paratmeters for all
1103 * atoms in the force field. To go around this for testing purposes
1104 * comment out the nerror++ counter temporarily
1106 nmiss = 0;
1107 for(i=0;i<get_atomtype_ntypes(atype);i++)
1109 if (get_atomtype_radius(i,atype) < 0 ||
1110 get_atomtype_vol(i,atype) < 0 ||
1111 get_atomtype_surftens(i,atype) < 0 ||
1112 get_atomtype_gb_radius(i,atype) < 0 ||
1113 get_atomtype_S_hct(i,atype) < 0)
1115 fprintf(stderr,"\nGB parameter(s) missing or negative for atom type '%s'\n",
1116 get_atomtype_name(i,atype));
1117 nmiss++;
1121 if (nmiss > 0)
1123 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);
1128 int main (int argc, char *argv[])
1130 static const char *desc[] = {
1131 "The gromacs preprocessor",
1132 "reads a molecular topology file, checks the validity of the",
1133 "file, expands the topology from a molecular description to an atomic",
1134 "description. The topology file contains information about",
1135 "molecule types and the number of molecules, the preprocessor",
1136 "copies each molecule as needed. ",
1137 "There is no limitation on the number of molecule types. ",
1138 "Bonds and bond-angles can be converted into constraints, separately",
1139 "for hydrogens and heavy atoms.",
1140 "Then a coordinate file is read and velocities can be generated",
1141 "from a Maxwellian distribution if requested.",
1142 "grompp also reads parameters for the mdrun ",
1143 "(eg. number of MD steps, time step, cut-off), and others such as",
1144 "NEMD parameters, which are corrected so that the net acceleration",
1145 "is zero.",
1146 "Eventually a binary file is produced that can serve as the sole input",
1147 "file for the MD program.[PAR]",
1149 "grompp uses the atom names from the topology file. The atom names",
1150 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1151 "warnings when they do not match the atom names in the topology.",
1152 "Note that the atom names are irrelevant for the simulation as",
1153 "only the atom types are used for generating interaction parameters.[PAR]",
1155 "grompp uses a built-in preprocessor to resolve includes, macros ",
1156 "etcetera. The preprocessor supports the following keywords:[BR]",
1157 "#ifdef VARIABLE[BR]",
1158 "#ifndef VARIABLE[BR]",
1159 "#else[BR]",
1160 "#endif[BR]",
1161 "#define VARIABLE[BR]",
1162 "#undef VARIABLE[BR]"
1163 "#include \"filename\"[BR]",
1164 "#include <filename>[BR]",
1165 "The functioning of these statements in your topology may be modulated by",
1166 "using the following two flags in your [TT]mdp[tt] file:[BR]",
1167 "define = -DVARIABLE1 -DVARIABLE2[BR]",
1168 "include = -I/home/john/doe[BR]",
1169 "For further information a C-programming textbook may help you out.",
1170 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1171 "topology file written out so that you can verify its contents.[PAR]",
1173 "If your system does not have a c-preprocessor, you can still",
1174 "use grompp, but you do not have access to the features ",
1175 "from the cpp. Command line options to the c-preprocessor can be given",
1176 "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1178 "When using position restraints a file with restraint coordinates",
1179 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1180 "with respect to the conformation from the [TT]-c[tt] option.",
1181 "For free energy calculation the the coordinates for the B topology",
1182 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1183 "those of the A topology.[PAR]",
1185 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1186 "The last frame with coordinates and velocities will be read,",
1187 "unless the [TT]-time[tt] option is used.",
1188 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1189 "in your [TT].mdp[tt] file. An energy file can be supplied with",
1190 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1191 "variables. Note that for continuation it is better and easier to supply",
1192 "a checkpoint file directly to mdrun, since that always contains",
1193 "the complete state of the system and you don't need to generate",
1194 "a new run input file. Note that if you only want to change the number",
1195 "of run steps tpbconv is more convenient than grompp.[PAR]",
1197 "By default all bonded interactions which have constant energy due to",
1198 "virtual site constructions will be removed. If this constant energy is",
1199 "not zero, this will result in a shift in the total energy. All bonded",
1200 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1201 "all constraints for distances which will be constant anyway because",
1202 "of virtual site constructions will be removed. If any constraints remain",
1203 "which involve virtual sites, a fatal error will result.[PAR]"
1205 "To verify your run input file, please make notice of all warnings",
1206 "on the screen, and correct where necessary. Do also look at the contents",
1207 "of the [TT]mdout.mdp[tt] file, this contains comment lines, as well as",
1208 "the input that [TT]grompp[tt] has read. If in doubt you can start grompp",
1209 "with the [TT]-debug[tt] option which will give you more information",
1210 "in a file called grompp.log (along with real debug info). Finally, you",
1211 "can see the contents of the run input file with the [TT]gmxdump[tt]",
1212 "program."
1214 t_gromppopts *opts;
1215 gmx_mtop_t *sys;
1216 int nmi;
1217 t_molinfo *mi;
1218 gpp_atomtype_t atype;
1219 t_inputrec *ir;
1220 int natoms,nvsite,comb,mt;
1221 t_params *plist;
1222 t_state state;
1223 matrix box;
1224 real max_spacing,fudgeQQ;
1225 double reppow;
1226 char fn[STRLEN],fnB[STRLEN];
1227 const char *mdparin;
1228 int ntype;
1229 gmx_bool bNeedVel,bGenVel;
1230 gmx_bool have_atomnumber;
1231 int n12,n13,n14;
1232 t_params *gb_plist = NULL;
1233 gmx_genborn_t *born = NULL;
1234 output_env_t oenv;
1235 gmx_bool bVerbose = FALSE;
1236 warninp_t wi;
1237 char warn_buf[STRLEN];
1239 t_filenm fnm[] = {
1240 { efMDP, NULL, NULL, ffREAD },
1241 { efMDP, "-po", "mdout", ffWRITE },
1242 { efSTX, "-c", NULL, ffREAD },
1243 { efSTX, "-r", NULL, ffOPTRD },
1244 { efSTX, "-rb", NULL, ffOPTRD },
1245 { efNDX, NULL, NULL, ffOPTRD },
1246 { efTOP, NULL, NULL, ffREAD },
1247 { efTOP, "-pp", "processed", ffOPTWR },
1248 { efTPX, "-o", NULL, ffWRITE },
1249 { efTRN, "-t", NULL, ffOPTRD },
1250 { efEDR, "-e", NULL, ffOPTRD },
1251 { efTRN, "-ref","rotref", ffOPTRW }
1253 #define NFILE asize(fnm)
1255 /* Command line options */
1256 static gmx_bool bRenum=TRUE;
1257 static gmx_bool bRmVSBds=TRUE,bZero=FALSE;
1258 static int i,maxwarn=0;
1259 static real fr_time=-1;
1260 t_pargs pa[] = {
1261 { "-v", FALSE, etBOOL,{&bVerbose},
1262 "Be loud and noisy" },
1263 { "-time", FALSE, etREAL, {&fr_time},
1264 "Take frame at or first after this time." },
1265 { "-rmvsbds",FALSE, etBOOL, {&bRmVSBds},
1266 "Remove constant bonded interactions with virtual sites" },
1267 { "-maxwarn", FALSE, etINT, {&maxwarn},
1268 "Number of allowed warnings during input processing" },
1269 { "-zero", FALSE, etBOOL, {&bZero},
1270 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1271 { "-renum", FALSE, etBOOL, {&bRenum},
1272 "Renumber atomtypes and minimize number of atomtypes" }
1275 CopyRight(stdout,argv[0]);
1277 /* Initiate some variables */
1278 snew(ir,1);
1279 snew(opts,1);
1280 init_ir(ir,opts);
1282 /* Parse the command line */
1283 parse_common_args(&argc,argv,0,NFILE,fnm,asize(pa),pa,
1284 asize(desc),desc,0,NULL,&oenv);
1286 wi = init_warning(TRUE,maxwarn);
1288 /* PARAMETER file processing */
1289 mdparin = opt2fn("-f",NFILE,fnm);
1290 set_warning_line(wi,mdparin,-1);
1291 get_ir(mdparin,opt2fn("-po",NFILE,fnm),ir,opts,wi);
1293 if (bVerbose)
1294 fprintf(stderr,"checking input for internal consistency...\n");
1295 check_ir(mdparin,ir,opts,wi);
1297 if (ir->ld_seed == -1) {
1298 ir->ld_seed = make_seed();
1299 fprintf(stderr,"Setting the LD random seed to %d\n",ir->ld_seed);
1302 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1303 bGenVel = (bNeedVel && opts->bGenVel);
1305 snew(plist,F_NRE);
1306 init_plist(plist);
1307 snew(sys,1);
1308 atype = init_atomtype();
1309 if (debug)
1310 pr_symtab(debug,0,"Just opened",&sys->symtab);
1312 strcpy(fn,ftp2fn(efTOP,NFILE,fnm));
1313 if (!gmx_fexist(fn))
1314 gmx_fatal(FARGS,"%s does not exist",fn);
1315 new_status(fn,opt2fn_null("-pp",NFILE,fnm),opt2fn("-c",NFILE,fnm),
1316 opts,ir,bZero,bGenVel,bVerbose,&state,
1317 atype,sys,&nmi,&mi,plist,&comb,&reppow,&fudgeQQ,
1318 opts->bMorse,
1319 wi);
1321 if (debug)
1322 pr_symtab(debug,0,"After new_status",&sys->symtab);
1324 if (count_constraints(sys,mi,wi) && (ir->eConstrAlg == econtSHAKE)) {
1325 if (ir->eI == eiCG || ir->eI == eiLBFGS) {
1326 sprintf(warn_buf,"Can not do %s with %s, use %s",
1327 EI(ir->eI),econstr_names[econtSHAKE],econstr_names[econtLINCS]);
1328 warning_error(wi,warn_buf);
1330 if (ir->bPeriodicMols) {
1331 sprintf(warn_buf,"Can not do periodic molecules with %s, use %s",
1332 econstr_names[econtSHAKE],econstr_names[econtLINCS]);
1333 warning_error(wi,warn_buf);
1337 /* If we are doing QM/MM, check that we got the atom numbers */
1338 have_atomnumber = TRUE;
1339 for (i=0; i<get_atomtype_ntypes(atype); i++) {
1340 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i,atype) >= 0);
1342 if (!have_atomnumber && ir->bQMMM)
1344 warning_error(wi,
1345 "\n"
1346 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1347 "field you are using does not contain atom numbers fields. This is an\n"
1348 "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1349 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1350 "an integer just before the mass column in ffXXXnb.itp.\n"
1351 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1354 if (ir->adress_type!=eAdressOff) {
1355 if ((ir->adress_const_wf>1) || (ir->adress_const_wf<0)) {
1356 warning_error(wi,"AdResS contant weighting function should be between 0 and 1\n\n");
1358 /** \TODO check size of ex+hy width against box size */
1361 /* Check for errors in the input now, since they might cause problems
1362 * during processing further down.
1364 check_warning_error(wi,FARGS);
1366 if (opt2bSet("-r",NFILE,fnm))
1367 sprintf(fn,"%s",opt2fn("-r",NFILE,fnm));
1368 else
1369 sprintf(fn,"%s",opt2fn("-c",NFILE,fnm));
1370 if (opt2bSet("-rb",NFILE,fnm))
1371 sprintf(fnB,"%s",opt2fn("-rb",NFILE,fnm));
1372 else
1373 strcpy(fnB,fn);
1375 if (nint_ftype(sys,mi,F_POSRES) > 0)
1377 if (bVerbose)
1379 fprintf(stderr,"Reading position restraint coords from %s",fn);
1380 if (strcmp(fn,fnB) == 0)
1382 fprintf(stderr,"\n");
1384 else
1386 fprintf(stderr," and %s\n",fnB);
1387 if (ir->efep != efepNO && ir->n_flambda > 0)
1389 warning_error(wi,"Can not change the position restraint reference coordinates with lambda togther with foreign lambda calculation.");
1393 gen_posres(sys,mi,fn,fnB,
1394 ir->refcoord_scaling,ir->ePBC,
1395 ir->posres_com,ir->posres_comB,
1396 wi);
1399 nvsite = 0;
1400 /* set parameters for virtual site construction (not for vsiten) */
1401 for(mt=0; mt<sys->nmoltype; mt++) {
1402 nvsite +=
1403 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1405 /* now throw away all obsolete bonds, angles and dihedrals: */
1406 /* note: constraints are ALWAYS removed */
1407 if (nvsite) {
1408 for(mt=0; mt<sys->nmoltype; mt++) {
1409 clean_vsite_bondeds(mi[mt].plist,sys->moltype[mt].atoms.nr,bRmVSBds);
1413 /* If we are using CMAP, setup the pre-interpolation grid */
1414 if(plist->ncmap>0)
1416 init_cmap_grid(&sys->ffparams.cmap_grid, plist->nc, plist->grid_spacing);
1417 setup_cmap(plist->grid_spacing, plist->nc, plist->cmap,&sys->ffparams.cmap_grid);
1420 set_wall_atomtype(atype,opts,ir);
1421 if (bRenum) {
1422 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1423 ntype = get_atomtype_ntypes(atype);
1426 if (ir->implicit_solvent != eisNO)
1428 /* Now we have renumbered the atom types, we can check the GBSA params */
1429 check_gbsa_params(ir,atype);
1431 /* Check that all atoms that have charge and/or LJ-parameters also have
1432 * sensible GB-parameters
1434 check_gbsa_params_charged(sys,atype);
1437 /* PELA: Copy the atomtype data to the topology atomtype list */
1438 copy_atomtype_atomtypes(atype,&(sys->atomtypes));
1440 if (debug)
1441 pr_symtab(debug,0,"After renum_atype",&sys->symtab);
1443 if (bVerbose)
1444 fprintf(stderr,"converting bonded parameters...\n");
1446 ntype = get_atomtype_ntypes(atype);
1447 convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1449 if (debug)
1450 pr_symtab(debug,0,"After convert_params",&sys->symtab);
1452 /* set ptype to VSite for virtual sites */
1453 for(mt=0; mt<sys->nmoltype; mt++) {
1454 set_vsites_ptype(FALSE,&sys->moltype[mt]);
1456 if (debug) {
1457 pr_symtab(debug,0,"After virtual sites",&sys->symtab);
1459 /* Check velocity for virtual sites and shells */
1460 if (bGenVel) {
1461 check_vel(sys,state.v);
1464 /* check masses */
1465 check_mol(sys,wi);
1467 for(i=0; i<sys->nmoltype; i++) {
1468 check_cg_sizes(ftp2fn(efTOP,NFILE,fnm),&sys->moltype[i].cgs,wi);
1471 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1473 check_bonds_timestep(sys,ir->delta_t,wi);
1476 check_warning_error(wi,FARGS);
1478 if (bVerbose)
1479 fprintf(stderr,"initialising group options...\n");
1480 do_index(mdparin,ftp2fn_null(efNDX,NFILE,fnm),
1481 sys,bVerbose,ir,
1482 bGenVel ? state.v : NULL,
1483 wi);
1485 /* Init the temperature coupling state */
1486 init_gtc_state(&state,ir->opts.ngtc,0,ir->opts.nhchainlength);
1488 if (bVerbose)
1489 fprintf(stderr,"Checking consistency between energy and charge groups...\n");
1490 check_eg_vs_cg(sys);
1492 if (debug)
1493 pr_symtab(debug,0,"After index",&sys->symtab);
1494 triple_check(mdparin,ir,sys,wi);
1495 close_symtab(&sys->symtab);
1496 if (debug)
1497 pr_symtab(debug,0,"After close",&sys->symtab);
1499 /* make exclusions between QM atoms */
1500 if (ir->bQMMM) {
1501 generate_qmexcl(sys,ir);
1504 if (ftp2bSet(efTRN,NFILE,fnm)) {
1505 if (bVerbose)
1506 fprintf(stderr,"getting data from old trajectory ...\n");
1507 cont_status(ftp2fn(efTRN,NFILE,fnm),ftp2fn_null(efEDR,NFILE,fnm),
1508 bNeedVel,bGenVel,fr_time,ir,&state,sys,oenv);
1511 if (ir->ePBC==epbcXY && ir->nwall!=2)
1513 clear_rvec(state.box[ZZ]);
1516 if (ir->rlist > 0)
1518 set_warning_line(wi,mdparin,-1);
1519 check_chargegroup_radii(sys,ir,state.x,wi);
1522 if (EEL_FULL(ir->coulombtype)) {
1523 /* Calculate the optimal grid dimensions */
1524 copy_mat(state.box,box);
1525 if (ir->ePBC==epbcXY && ir->nwall==2)
1526 svmul(ir->wall_ewald_zfac,box[ZZ],box[ZZ]);
1527 max_spacing = calc_grid(stdout,box,opts->fourierspacing,
1528 &(ir->nkx),&(ir->nky),&(ir->nkz));
1529 if ((ir->coulombtype == eelPPPM) && (max_spacing > 0.1)) {
1530 set_warning_line(wi,mdparin,-1);
1531 warning_note(wi,"Grid spacing larger then 0.1 while using PPPM.");
1535 if (ir->ePull != epullNO)
1536 set_pull_init(ir,sys,state.x,state.box,oenv,opts->pull_start);
1538 if (ir->bRot)
1540 set_reference_positions(ir->rot,sys,state.x,state.box,
1541 opt2fn("-ref",NFILE,fnm),opt2bSet("-ref",NFILE,fnm),
1542 wi);
1545 /* reset_multinr(sys); */
1547 if (EEL_PME(ir->coulombtype)) {
1548 float ratio = pme_load_estimate(sys,ir,state.box);
1549 fprintf(stderr,"Estimate for the relative computational load of the PME mesh part: %.2f\n",ratio);
1550 /* With free energy we might need to do PME both for the A and B state
1551 * charges. This will double the cost, but the optimal performance will
1552 * then probably be at a slightly larger cut-off and grid spacing.
1554 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
1555 (ir->efep != efepNO && ratio > 2.0/3.0)) {
1556 warning_note(wi,
1557 "The optimal PME mesh load for parallel simulations is below 0.5\n"
1558 "and for highly parallel simulations between 0.25 and 0.33,\n"
1559 "for higher performance, increase the cut-off and the PME grid spacing");
1564 char warn_buf[STRLEN];
1565 double cio = compute_io(ir,sys->natoms,&sys->groups,F_NRE,1);
1566 sprintf(warn_buf,"This run will generate roughly %.0f Mb of data",cio);
1567 if (cio > 2000) {
1568 set_warning_line(wi,mdparin,-1);
1569 warning_note(wi,warn_buf);
1570 } else {
1571 printf("%s\n",warn_buf);
1575 if (bVerbose)
1576 fprintf(stderr,"writing run input file...\n");
1578 done_warning(wi,FARGS);
1580 state.lambda = ir->init_lambda;
1581 write_tpx_state(ftp2fn(efTPX,NFILE,fnm),ir,&state,sys);
1583 thanx(stderr);
1585 return 0;