Merge branch 'master' of git://git.gromacs.org/gromacs
[gromacs/adressmacs.git] / src / kernel / readir.c
blob946f46e12c7e2bd904dd5e8892674e8b70d6bcd6
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.0
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 <ctype.h>
41 #include <stdlib.h>
42 #include <limits.h>
43 #include "sysstuff.h"
44 #include "smalloc.h"
45 #include "typedefs.h"
46 #include "physics.h"
47 #include "names.h"
48 #include "gmx_fatal.h"
49 #include "macros.h"
50 #include "index.h"
51 #include "symtab.h"
52 #include "string2.h"
53 #include "readinp.h"
54 #include "warninp.h"
55 #include "readir.h"
56 #include "toputil.h"
57 #include "index.h"
58 #include "network.h"
59 #include "vec.h"
60 #include "pbc.h"
61 #include "mtop_util.h"
62 #include "chargegroup.h"
64 #define MAXPTR 254
65 #define NOGID 255
67 /* Resource parameters
68 * Do not change any of these until you read the instruction
69 * in readinp.h. Some cpp's do not take spaces after the backslash
70 * (like the c-shell), which will give you a very weird compiler
71 * message.
74 static char tcgrps[STRLEN],tau_t[STRLEN],ref_t[STRLEN],
75 acc[STRLEN],accgrps[STRLEN],freeze[STRLEN],frdim[STRLEN],
76 energy[STRLEN],user1[STRLEN],user2[STRLEN],vcm[STRLEN],xtc_grps[STRLEN],
77 couple_moltype[STRLEN],orirefitgrp[STRLEN],egptable[STRLEN],egpexcl[STRLEN],
78 wall_atomtype[STRLEN],wall_density[STRLEN],deform[STRLEN],QMMM[STRLEN],adress_refs[STRLEN],
79 adress_tf_grp_names[STRLEN];
80 static char foreign_lambda[STRLEN];
81 static char **pull_grp;
82 static char anneal[STRLEN],anneal_npoints[STRLEN],
83 anneal_time[STRLEN],anneal_temp[STRLEN];
84 static char QMmethod[STRLEN],QMbasis[STRLEN],QMcharge[STRLEN],QMmult[STRLEN],
85 bSH[STRLEN],CASorbitals[STRLEN], CASelectrons[STRLEN],SAon[STRLEN],
86 SAoff[STRLEN],SAsteps[STRLEN],bTS[STRLEN],bOPT[STRLEN];
87 static char efield_x[STRLEN],efield_xt[STRLEN],efield_y[STRLEN],
88 efield_yt[STRLEN],efield_z[STRLEN],efield_zt[STRLEN];
90 enum { egrptpALL, egrptpALL_GENREST, egrptpPART, egrptpONE };
93 /* Minimum number of time steps required for accurate coupling integration
94 * of first and second order thermo- and barostats:
96 int nstcmin1 = 10;
97 int nstcmin2 = 20;
100 void init_ir(t_inputrec *ir, t_gromppopts *opts)
102 snew(opts->include,STRLEN);
103 snew(opts->define,STRLEN);
106 static void _low_check(bool b,char *s,warninp_t wi)
108 if (b)
110 warning_error(wi,s);
114 static void check_nst(const char *desc_nst,int nst,
115 const char *desc_p,int *p,
116 warninp_t wi)
118 char buf[STRLEN];
120 if (*p > 0 && *p % nst != 0)
122 /* Round up to the next multiple of nst */
123 *p = ((*p)/nst + 1)*nst;
124 sprintf(buf,"%s should be a multiple of %s, changing %s to %d\n",
125 desc_p,desc_nst,desc_p,*p);
126 warning(wi,buf);
130 static bool ir_NVE(const t_inputrec *ir)
132 return ((ir->eI == eiMD || EI_VV(ir->eI)) && ir->etc == etcNO);
135 void check_ir(const char *mdparin,t_inputrec *ir, t_gromppopts *opts,
136 warninp_t wi)
137 /* Check internal consistency */
139 /* Strange macro: first one fills the err_buf, and then one can check
140 * the condition, which will print the message and increase the error
141 * counter.
143 #define CHECK(b) _low_check(b,err_buf,wi)
144 char err_buf[256],warn_buf[STRLEN];
145 int ns_type=0;
146 real dt_coupl=0;
147 int nstcmin;
149 set_warning_line(wi,mdparin,-1);
151 /* BASIC CUT-OFF STUFF */
152 if (ir->rlist == 0 ||
153 !((EEL_MIGHT_BE_ZERO_AT_CUTOFF(ir->coulombtype) && ir->rcoulomb > ir->rlist) ||
154 (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(ir->vdwtype) && ir->rvdw > ir->rlist))) {
155 /* No switched potential and/or no twin-range:
156 * we can set the long-range cut-off to the maximum of the other cut-offs.
158 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
159 } else if (ir->rlistlong < 0) {
160 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
161 sprintf(warn_buf,"rlistlong was not set, setting it to %g (no buffer)",
162 ir->rlistlong);
163 warning(wi,warn_buf);
165 if (ir->rlistlong == 0 && ir->ePBC != epbcNONE) {
166 warning_error(wi,"Can not have an infinite cut-off with PBC");
168 if (ir->rlistlong > 0 && (ir->rlist == 0 || ir->rlistlong < ir->rlist)) {
169 warning_error(wi,"rlistlong can not be shorter than rlist");
171 if (IR_TWINRANGE(*ir) && ir->nstlist <= 0) {
172 warning_error(wi,"Can not have nstlist<=0 with twin-range interactions");
175 /* GENERAL INTEGRATOR STUFF */
176 if (!(ir->eI == eiMD || EI_VV(ir->eI)))
178 ir->etc = etcNO;
180 if (!EI_DYNAMICS(ir->eI))
182 ir->epc = epcNO;
184 if (EI_DYNAMICS(ir->eI))
186 if (ir->nstcalcenergy < 0)
188 if (EI_VV(ir->eI))
190 /* VV coupling algorithms currently only support 1 */
191 ir->nstcalcenergy = 1;
193 else
195 if (ir->nstlist > 0)
197 ir->nstcalcenergy = ir->nstlist;
199 else
201 ir->nstcalcenergy = 10;
205 if (ir->etc != etcNO || ir->epc != epcNO)
207 if (ir->nstcalcenergy == 0)
209 gmx_fatal(FARGS,"Can not have nstcalcenergy=0 with global T/P-coupling");
211 if (EI_VV(ir->eI))
213 sprintf(err_buf,"T- and P-coupling with VV integrators currently only supports nstcalcenergy=1");
214 CHECK(ir->nstcalcenergy > 1);
217 if (IR_TWINRANGE(*ir))
219 check_nst("nstlist",ir->nstlist,
220 "nstcalcenergy",&ir->nstcalcenergy,wi);
222 dt_coupl = ir->nstcalcenergy*ir->delta_t;
224 if (ir->nstcalcenergy > 1)
226 /* for storing exact averages nstenergy should be
227 * a multiple of nstcalcenergy
229 check_nst("nstcalcenergy",ir->nstcalcenergy,
230 "nstenergy",&ir->nstenergy,wi);
231 if (ir->efep != efepNO)
233 /* nstdhdl should be a multiple of nstcalcenergy */
234 check_nst("nstcalcenergy",ir->nstcalcenergy,
235 "nstdhdl",&ir->nstdhdl,wi);
240 /* LD STUFF */
241 if ((EI_SD(ir->eI) || ir->eI == eiBD) &&
242 ir->bContinuation && ir->ld_seed != -1) {
243 warning_note(wi,"You are doing a continuation with SD or BD, make sure that ld_seed is different from the previous run (using ld_seed=-1 will ensure this)");
246 /* TPI STUFF */
247 if (EI_TPI(ir->eI)) {
248 sprintf(err_buf,"TPI only works with pbc = %s",epbc_names[epbcXYZ]);
249 CHECK(ir->ePBC != epbcXYZ);
250 sprintf(err_buf,"TPI only works with ns = %s",ens_names[ensGRID]);
251 CHECK(ir->ns_type != ensGRID);
252 sprintf(err_buf,"with TPI nstlist should be larger than zero");
253 CHECK(ir->nstlist <= 0);
254 sprintf(err_buf,"TPI does not work with full electrostatics other than PME");
255 CHECK(EEL_FULL(ir->coulombtype) && !EEL_PME(ir->coulombtype));
258 /* SHAKE / LINCS */
259 if ( (opts->nshake > 0) && (opts->bMorse) ) {
260 sprintf(warn_buf,
261 "Using morse bond-potentials while constraining bonds is useless");
262 warning(wi,warn_buf);
265 sprintf(err_buf,"shake_tol must be > 0 instead of %g while using shake",
266 ir->shake_tol);
267 CHECK(((ir->shake_tol <= 0.0) && (opts->nshake>0) &&
268 (ir->eConstrAlg == econtSHAKE)));
270 /* PBC/WALLS */
271 sprintf(err_buf,"walls only work with pbc=%s",epbc_names[epbcXY]);
272 CHECK(ir->nwall && ir->ePBC!=epbcXY);
274 /* VACUUM STUFF */
275 if (ir->ePBC != epbcXYZ && ir->nwall != 2) {
276 if (ir->ePBC == epbcNONE) {
277 if (ir->epc != epcNO) {
278 warning(wi,"Turning off pressure coupling for vacuum system");
279 ir->epc = epcNO;
281 } else {
282 sprintf(err_buf,"Can not have pressure coupling with pbc=%s",
283 epbc_names[ir->ePBC]);
284 CHECK(ir->epc != epcNO);
286 sprintf(err_buf,"Can not have Ewald with pbc=%s",epbc_names[ir->ePBC]);
287 CHECK(EEL_FULL(ir->coulombtype));
289 sprintf(err_buf,"Can not have dispersion correction with pbc=%s",
290 epbc_names[ir->ePBC]);
291 CHECK(ir->eDispCorr != edispcNO);
294 if (ir->rlist == 0.0) {
295 sprintf(err_buf,"can only have neighborlist cut-off zero (=infinite)\n"
296 "with coulombtype = %s or coulombtype = %s\n"
297 "without periodic boundary conditions (pbc = %s) and\n"
298 "rcoulomb and rvdw set to zero",
299 eel_names[eelCUT],eel_names[eelUSER],epbc_names[epbcNONE]);
300 CHECK(((ir->coulombtype != eelCUT) && (ir->coulombtype != eelUSER)) ||
301 (ir->ePBC != epbcNONE) ||
302 (ir->rcoulomb != 0.0) || (ir->rvdw != 0.0));
304 if (ir->nstlist < 0) {
305 warning_error(wi,"Can not have heuristic neighborlist updates without cut-off");
307 if (ir->nstlist > 0) {
308 warning_note(wi,"Simulating without cut-offs is usually (slightly) faster with nstlist=0, nstype=simple and particle decomposition");
312 /* COMM STUFF */
313 if (ir->nstcomm == 0) {
314 ir->comm_mode = ecmNO;
316 if (ir->comm_mode != ecmNO) {
317 if (ir->nstcomm < 0) {
318 warning(wi,"If you want to remove the rotation around the center of mass, you should set comm_mode = Angular instead of setting nstcomm < 0. nstcomm is modified to its absolute value");
319 ir->nstcomm = abs(ir->nstcomm);
322 if (ir->nstcalcenergy > 0 && ir->nstcomm < ir->nstcalcenergy) {
323 warning_note(wi,"nstcomm < nstcalcenergy defeats the purpose of nstcalcenergy, setting nstcomm to nstcalcenergy");
324 ir->nstcomm = ir->nstcalcenergy;
327 if (ir->comm_mode == ecmANGULAR) {
328 sprintf(err_buf,"Can not remove the rotation around the center of mass with periodic molecules");
329 CHECK(ir->bPeriodicMols);
330 if (ir->ePBC != epbcNONE)
331 warning(wi,"Removing the rotation around the center of mass in a periodic system (this is not a problem when you have only one molecule).");
335 if (EI_STATE_VELOCITY(ir->eI) && ir->ePBC == epbcNONE && ir->comm_mode != ecmANGULAR) {
336 warning_note(wi,"Tumbling and or flying ice-cubes: We are not removing rotation around center of mass in a non-periodic system. You should probably set comm_mode = ANGULAR.");
339 sprintf(err_buf,"Free-energy not implemented for Ewald and PPPM");
340 CHECK((ir->coulombtype==eelEWALD || ir->coulombtype==eelPPPM)
341 && (ir->efep!=efepNO));
343 sprintf(err_buf,"Twin-range neighbour searching (NS) with simple NS"
344 " algorithm not implemented");
345 CHECK(((ir->rcoulomb > ir->rlist) || (ir->rvdw > ir->rlist))
346 && (ir->ns_type == ensSIMPLE));
348 /* TEMPERATURE COUPLING */
349 if (ir->etc == etcYES)
351 ir->etc = etcBERENDSEN;
352 warning_note(wi,"Old option for temperature coupling given: "
353 "changing \"yes\" to \"Berendsen\"\n");
356 if (ir->etc == etcNOSEHOOVER)
358 if (ir->opts.nhchainlength < 1)
360 sprintf(warn_buf,"number of Nose-Hoover chains (currently %d) cannot be less than 1,reset to 1\n",ir->opts.nhchainlength);
361 ir->opts.nhchainlength =1;
362 warning(wi,warn_buf);
365 if (ir->etc==etcNOSEHOOVER && !EI_VV(ir->eI) && ir->opts.nhchainlength > 1)
367 warning_note(wi,"leapfrog does not yet support Nose-Hoover chains, nhchainlength reset to 1");
368 ir->opts.nhchainlength = 1;
371 else
373 ir->opts.nhchainlength = 0;
376 if (ir->etc == etcBERENDSEN)
378 sprintf(warn_buf,"The %s thermostat does not generate the correct kinetic energy distribution. You might want to consider using the %s thermostat.",
379 ETCOUPLTYPE(ir->etc),ETCOUPLTYPE(etcVRESCALE));
380 warning_note(wi,warn_buf);
383 if ((ir->etc==etcNOSEHOOVER || ir->etc==etcANDERSEN || ir->etc==etcANDERSENINTERVAL)
384 && ir->epc==epcBERENDSEN)
386 sprintf(warn_buf,"Using Berendsen pressure coupling invalidates the "
387 "true ensemble for the thermostat");
388 warning(wi,warn_buf);
391 /* PRESSURE COUPLING */
392 if (ir->epc == epcISOTROPIC)
394 ir->epc = epcBERENDSEN;
395 warning_note(wi,"Old option for pressure coupling given: "
396 "changing \"Isotropic\" to \"Berendsen\"\n");
399 if (ir->epc != epcNO)
401 sprintf(err_buf,"tau_p must be > 0 instead of %g\n",ir->tau_p);
402 CHECK(ir->tau_p <= 0);
404 nstcmin = (ir->epc == epcBERENDSEN ? nstcmin1 : nstcmin2);
405 if (ir->tau_p < nstcmin*dt_coupl)
407 sprintf(warn_buf,"For proper integration of the %s barostat, tau_p (%g) should be at least %d times larger than nstcalcenergy*dt (%g)",
408 EPCOUPLTYPE(ir->epc),ir->tau_p,nstcmin,dt_coupl);
409 warning(wi,warn_buf);
412 sprintf(err_buf,"compressibility must be > 0 when using pressure"
413 " coupling %s\n",EPCOUPLTYPE(ir->epc));
414 CHECK(ir->compress[XX][XX] < 0 || ir->compress[YY][YY] < 0 ||
415 ir->compress[ZZ][ZZ] < 0 ||
416 (trace(ir->compress) == 0 && ir->compress[YY][XX] <= 0 &&
417 ir->compress[ZZ][XX] <= 0 && ir->compress[ZZ][YY] <= 0));
419 sprintf(err_buf,"pressure coupling with PPPM not implemented, use PME");
420 CHECK(ir->coulombtype == eelPPPM);
423 else if (ir->coulombtype == eelPPPM)
425 sprintf(warn_buf,"The pressure with PPPM is incorrect, if you need the pressure use PME");
426 warning(wi,warn_buf);
429 if (EI_VV(ir->eI))
431 if (ir->epc > epcNO)
433 if (ir->epc!=epcMTTK)
435 warning_error(wi,"NPT only defined for vv using Martyna-Tuckerman-Tobias-Klein equations");
440 /* ELECTROSTATICS */
441 /* More checks are in triple check (grompp.c) */
442 if (ir->coulombtype == eelSWITCH) {
443 sprintf(warn_buf,"coulombtype = %s is only for testing purposes and can lead to serious artifacts, advice: use coulombtype = %s",
444 eel_names[ir->coulombtype],
445 eel_names[eelRF_ZERO]);
446 warning(wi,warn_buf);
449 if (ir->epsilon_r!=1 && ir->implicit_solvent==eisGBSA) {
450 sprintf(warn_buf,"epsilon_r = %g with GB implicit solvent, will use this value for inner dielectric",ir->epsilon_r);
451 warning_note(wi,warn_buf);
454 if (EEL_RF(ir->coulombtype) && ir->epsilon_rf==1 && ir->epsilon_r!=1) {
455 sprintf(warn_buf,"epsilon_r = %g and epsilon_rf = 1 with reaction field, assuming old format and exchanging epsilon_r and epsilon_rf",ir->epsilon_r);
456 warning(wi,warn_buf);
457 ir->epsilon_rf = ir->epsilon_r;
458 ir->epsilon_r = 1.0;
461 if (getenv("GALACTIC_DYNAMICS") == NULL) {
462 sprintf(err_buf,"epsilon_r must be >= 0 instead of %g\n",ir->epsilon_r);
463 CHECK(ir->epsilon_r < 0);
466 if (EEL_RF(ir->coulombtype)) {
467 /* reaction field (at the cut-off) */
469 if (ir->coulombtype == eelRF_ZERO) {
470 sprintf(err_buf,"With coulombtype = %s, epsilon_rf must be 0",
471 eel_names[ir->coulombtype]);
472 CHECK(ir->epsilon_rf != 0);
475 sprintf(err_buf,"epsilon_rf must be >= epsilon_r");
476 CHECK((ir->epsilon_rf < ir->epsilon_r && ir->epsilon_rf != 0) ||
477 (ir->epsilon_r == 0));
478 if (ir->epsilon_rf == ir->epsilon_r) {
479 sprintf(warn_buf,"Using epsilon_rf = epsilon_r with %s does not make sense",
480 eel_names[ir->coulombtype]);
481 warning(wi,warn_buf);
484 /* Allow rlist>rcoulomb for tabulated long range stuff. This just
485 * means the interaction is zero outside rcoulomb, but it helps to
486 * provide accurate energy conservation.
488 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(ir->coulombtype)) {
489 if (EEL_SWITCHED(ir->coulombtype)) {
490 sprintf(err_buf,
491 "With coulombtype = %s rcoulomb_switch must be < rcoulomb",
492 eel_names[ir->coulombtype]);
493 CHECK(ir->rcoulomb_switch >= ir->rcoulomb);
495 } else if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype)) {
496 sprintf(err_buf,"With coulombtype = %s, rcoulomb must be >= rlist",
497 eel_names[ir->coulombtype]);
498 CHECK(ir->rlist > ir->rcoulomb);
501 if (EEL_FULL(ir->coulombtype)) {
502 if (ir->coulombtype==eelPMESWITCH || ir->coulombtype==eelPMEUSER) {
503 sprintf(err_buf,"With coulombtype = %s, rcoulomb must be <= rlist",
504 eel_names[ir->coulombtype]);
505 CHECK(ir->rcoulomb > ir->rlist);
506 } else {
507 if (ir->coulombtype == eelPME) {
508 sprintf(err_buf,
509 "With coulombtype = %s, rcoulomb must be equal to rlist\n"
510 "If you want optimal energy conservation or exact integration use %s",
511 eel_names[ir->coulombtype],eel_names[eelPMESWITCH]);
512 } else {
513 sprintf(err_buf,
514 "With coulombtype = %s, rcoulomb must be equal to rlist",
515 eel_names[ir->coulombtype]);
517 CHECK(ir->rcoulomb != ir->rlist);
521 if (EEL_PME(ir->coulombtype)) {
522 if (ir->pme_order < 3) {
523 warning_error(wi,"pme_order can not be smaller than 3");
527 if (ir->nwall==2 && EEL_FULL(ir->coulombtype)) {
528 if (ir->ewald_geometry == eewg3D) {
529 sprintf(warn_buf,"With pbc=%s you should use ewald_geometry=%s",
530 epbc_names[ir->ePBC],eewg_names[eewg3DC]);
531 warning(wi,warn_buf);
533 /* This check avoids extra pbc coding for exclusion corrections */
534 sprintf(err_buf,"wall_ewald_zfac should be >= 2");
535 CHECK(ir->wall_ewald_zfac < 2);
538 if (EVDW_SWITCHED(ir->vdwtype)) {
539 sprintf(err_buf,"With vdwtype = %s rvdw_switch must be < rvdw",
540 evdw_names[ir->vdwtype]);
541 CHECK(ir->rvdw_switch >= ir->rvdw);
542 } else if (ir->vdwtype == evdwCUT) {
543 sprintf(err_buf,"With vdwtype = %s, rvdw must be >= rlist",evdw_names[ir->vdwtype]);
544 CHECK(ir->rlist > ir->rvdw);
546 if (EEL_IS_ZERO_AT_CUTOFF(ir->coulombtype)
547 && (ir->rlistlong <= ir->rcoulomb)) {
548 sprintf(warn_buf,"For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rcoulomb.",
549 IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
550 warning_note(wi,warn_buf);
552 if (EVDW_SWITCHED(ir->vdwtype) && (ir->rlistlong <= ir->rvdw)) {
553 sprintf(warn_buf,"For energy conservation with switch/shift potentials, %s should be 0.1 to 0.3 nm larger than rvdw.",
554 IR_TWINRANGE(*ir) ? "rlistlong" : "rlist");
555 warning_note(wi,warn_buf);
558 if (ir->vdwtype == evdwUSER && ir->eDispCorr != edispcNO) {
559 warning_note(wi,"You have selected user tables with dispersion correction, the dispersion will be corrected to -C6/r^6 beyond rvdw_switch (the tabulated interaction between rvdw_switch and rvdw will not be double counted). Make sure that you really want dispersion correction to -C6/r^6.");
562 if (ir->nstlist == -1) {
563 sprintf(err_buf,
564 "nstlist=-1 only works with switched or shifted potentials,\n"
565 "suggestion: use vdw-type=%s and coulomb-type=%s",
566 evdw_names[evdwSHIFT],eel_names[eelPMESWITCH]);
567 CHECK(!(EEL_MIGHT_BE_ZERO_AT_CUTOFF(ir->coulombtype) &&
568 EVDW_MIGHT_BE_ZERO_AT_CUTOFF(ir->vdwtype)));
570 sprintf(err_buf,"With nstlist=-1 rvdw and rcoulomb should be smaller than rlist to account for diffusion and possibly charge-group radii");
571 CHECK(ir->rvdw >= ir->rlist || ir->rcoulomb >= ir->rlist);
573 sprintf(err_buf,"nstlist can not be smaller than -1");
574 CHECK(ir->nstlist < -1);
576 if (ir->eI == eiLBFGS && (ir->coulombtype==eelCUT || ir->vdwtype==evdwCUT)
577 && ir->rvdw != 0) {
578 warning(wi,"For efficient BFGS minimization, use switch/shift/pme instead of cut-off.");
581 if (ir->eI == eiLBFGS && ir->nbfgscorr <= 0) {
582 warning(wi,"Using L-BFGS with nbfgscorr<=0 just gets you steepest descent.");
585 /* FREE ENERGY */
586 if (ir->efep != efepNO) {
587 sprintf(err_buf,"The soft-core power is %d and can only be 1 or 2",
588 ir->sc_power);
589 CHECK(ir->sc_alpha!=0 && ir->sc_power!=1 && ir->sc_power!=2);
592 /* ENERGY CONSERVATION */
593 if (ir_NVE(ir))
595 if (!EVDW_MIGHT_BE_ZERO_AT_CUTOFF(ir->vdwtype) && ir->rvdw > 0)
597 sprintf(warn_buf,"You are using a cut-off for VdW interactions with NVE, for good energy conservation use vdwtype = %s (possibly with DispCorr)",
598 evdw_names[evdwSHIFT]);
599 warning_note(wi,warn_buf);
601 if (!EEL_MIGHT_BE_ZERO_AT_CUTOFF(ir->coulombtype) && ir->rcoulomb > 0)
603 sprintf(warn_buf,"You are using a cut-off for electrostatics with NVE, for good energy conservation use coulombtype = %s or %s",
604 eel_names[eelPMESWITCH],eel_names[eelRF_ZERO]);
605 warning_note(wi,warn_buf);
609 /* IMPLICIT SOLVENT */
610 if(ir->coulombtype==eelGB_NOTUSED)
612 ir->coulombtype=eelCUT;
613 ir->implicit_solvent=eisGBSA;
614 fprintf(stderr,"Note: Old option for generalized born electrostatics given:\n"
615 "Changing coulombtype from \"generalized-born\" to \"cut-off\" and instead\n"
616 "setting implicit_solvent value to \"GBSA\" in input section.\n");
619 if(ir->implicit_solvent==eisGBSA)
621 sprintf(err_buf,"With GBSA implicit solvent, rgbradii must be equal to rlist.");
622 CHECK(ir->rgbradii != ir->rlist);
624 if(ir->coulombtype!=eelCUT)
626 sprintf(err_buf,"With GBSA, coulombtype must be equal to %s\n",eel_names[eelCUT]);
627 CHECK(ir->coulombtype!=eelCUT);
629 if(ir->vdwtype!=evdwCUT)
631 sprintf(err_buf,"With GBSA, vdw-type must be equal to %s\n",evdw_names[evdwCUT]);
632 CHECK(ir->vdwtype!=evdwCUT);
635 if(ir->nstgbradii<1)
637 sprintf(warn_buf,"Using GBSA with nstgbradii<1, setting nstgbradii=1");
638 warning_note(wi,warn_buf);
639 ir->nstgbradii=1;
644 static int str_nelem(const char *str,int maxptr,char *ptr[])
646 int np=0;
647 char *copy0,*copy;
649 copy0=strdup(str);
650 copy=copy0;
651 ltrim(copy);
652 while (*copy != '\0') {
653 if (np >= maxptr)
654 gmx_fatal(FARGS,"Too many groups on line: '%s' (max is %d)",
655 str,maxptr);
656 if (ptr)
657 ptr[np]=copy;
658 np++;
659 while ((*copy != '\0') && !isspace(*copy))
660 copy++;
661 if (*copy != '\0') {
662 *copy='\0';
663 copy++;
665 ltrim(copy);
667 if (ptr == NULL)
668 sfree(copy0);
670 return np;
673 static void parse_n_double(char *str,int *n,double **r)
675 char *ptr[MAXPTR];
676 int i;
678 *n = str_nelem(str,MAXPTR,ptr);
680 snew(*r,*n);
681 for(i=0; i<*n; i++) {
682 (*r)[i] = strtod(ptr[i],NULL);
686 static void do_wall_params(t_inputrec *ir,
687 char *wall_atomtype, char *wall_density,
688 t_gromppopts *opts)
690 int nstr,i;
691 char *names[MAXPTR];
692 double dbl;
694 opts->wall_atomtype[0] = NULL;
695 opts->wall_atomtype[1] = NULL;
697 ir->wall_atomtype[0] = -1;
698 ir->wall_atomtype[1] = -1;
699 ir->wall_density[0] = 0;
700 ir->wall_density[1] = 0;
702 if (ir->nwall > 0) {
703 nstr = str_nelem(wall_atomtype,MAXPTR,names);
704 if (nstr != ir->nwall)
705 gmx_fatal(FARGS,"Expected %d elements for wall_atomtype, found %d",
706 ir->nwall,nstr);
707 for(i=0; i<ir->nwall; i++)
708 opts->wall_atomtype[i] = strdup(names[i]);
710 if (ir->wall_type != ewtTABLE) {
711 nstr = str_nelem(wall_density,MAXPTR,names);
712 if (nstr != ir->nwall)
713 gmx_fatal(FARGS,"Expected %d elements for wall_density, found %d",
714 ir->nwall,nstr);
715 for(i=0; i<ir->nwall; i++) {
716 sscanf(names[i],"%lf",&dbl);
717 if (dbl <= 0)
718 gmx_fatal(FARGS,"wall_density[%d] = %f\n",i,dbl);
719 ir->wall_density[i] = dbl;
725 static void add_wall_energrps(gmx_groups_t *groups,int nwall,t_symtab *symtab)
727 int i;
728 t_grps *grps;
729 char str[STRLEN];
731 if (nwall > 0) {
732 srenew(groups->grpname,groups->ngrpname+nwall);
733 grps = &(groups->grps[egcENER]);
734 srenew(grps->nm_ind,grps->nr+nwall);
735 for(i=0; i<nwall; i++) {
736 sprintf(str,"wall%d",i);
737 groups->grpname[groups->ngrpname] = put_symtab(symtab,str);
738 grps->nm_ind[grps->nr++] = groups->ngrpname++;
743 void get_ir(const char *mdparin,const char *mdparout,
744 t_inputrec *ir,t_gromppopts *opts,
745 warninp_t wi)
747 char *dumstr[2];
748 double dumdub[2][6];
749 t_inpfile *inp;
750 const char *tmp;
751 int i,j,m,ninp;
752 char warn_buf[STRLEN];
754 inp = read_inpfile(mdparin, &ninp, NULL, wi);
756 snew(dumstr[0],STRLEN);
757 snew(dumstr[1],STRLEN);
759 REM_TYPE("title");
760 REM_TYPE("cpp");
761 REM_TYPE("domain-decomposition");
762 REPL_TYPE("unconstrained-start","continuation");
763 REM_TYPE("dihre-tau");
764 REM_TYPE("nstdihreout");
765 REM_TYPE("nstcheckpoint");
767 CCTYPE ("VARIOUS PREPROCESSING OPTIONS");
768 CTYPE ("Preprocessor information: use cpp syntax.");
769 CTYPE ("e.g.: -I/home/joe/doe -I/home/mary/roe");
770 STYPE ("include", opts->include, NULL);
771 CTYPE ("e.g.: -DI_Want_Cookies -DMe_Too");
772 STYPE ("define", opts->define, NULL);
774 CCTYPE ("RUN CONTROL PARAMETERS");
775 EETYPE("integrator", ir->eI, ei_names);
776 CTYPE ("Start time and timestep in ps");
777 RTYPE ("tinit", ir->init_t, 0.0);
778 RTYPE ("dt", ir->delta_t, 0.001);
779 STEPTYPE ("nsteps", ir->nsteps, 0);
780 CTYPE ("For exact run continuation or redoing part of a run");
781 STEPTYPE ("init_step",ir->init_step, 0);
782 CTYPE ("Part index is updated automatically on checkpointing (keeps files separate)");
783 ITYPE ("simulation_part", ir->simulation_part, 1);
784 CTYPE ("mode for center of mass motion removal");
785 CTYPE ("energy calculation and T/P-coupling frequency");
786 ITYPE ("nstcalcenergy",ir->nstcalcenergy, -1);
787 EETYPE("comm-mode", ir->comm_mode, ecm_names);
788 CTYPE ("number of steps for center of mass motion removal");
789 ITYPE ("nstcomm", ir->nstcomm, 10);
790 CTYPE ("group(s) for center of mass motion removal");
791 STYPE ("comm-grps", vcm, NULL);
793 CCTYPE ("LANGEVIN DYNAMICS OPTIONS");
794 CTYPE ("Friction coefficient (amu/ps) and random seed");
795 RTYPE ("bd-fric", ir->bd_fric, 0.0);
796 ITYPE ("ld-seed", ir->ld_seed, 1993);
798 /* Em stuff */
799 CCTYPE ("ENERGY MINIMIZATION OPTIONS");
800 CTYPE ("Force tolerance and initial step-size");
801 RTYPE ("emtol", ir->em_tol, 10.0);
802 RTYPE ("emstep", ir->em_stepsize,0.01);
803 CTYPE ("Max number of iterations in relax_shells");
804 ITYPE ("niter", ir->niter, 20);
805 CTYPE ("Step size (ps^2) for minimization of flexible constraints");
806 RTYPE ("fcstep", ir->fc_stepsize, 0);
807 CTYPE ("Frequency of steepest descents steps when doing CG");
808 ITYPE ("nstcgsteep", ir->nstcgsteep, 1000);
809 ITYPE ("nbfgscorr", ir->nbfgscorr, 10);
811 CCTYPE ("TEST PARTICLE INSERTION OPTIONS");
812 RTYPE ("rtpi", ir->rtpi, 0.05);
814 /* Output options */
815 CCTYPE ("OUTPUT CONTROL OPTIONS");
816 CTYPE ("Output frequency for coords (x), velocities (v) and forces (f)");
817 ITYPE ("nstxout", ir->nstxout, 100);
818 ITYPE ("nstvout", ir->nstvout, 100);
819 ITYPE ("nstfout", ir->nstfout, 0);
820 ir->nstcheckpoint = 1000;
821 CTYPE ("Output frequency for energies to log file and energy file");
822 ITYPE ("nstlog", ir->nstlog, 100);
823 ITYPE ("nstenergy", ir->nstenergy, 100);
824 CTYPE ("Output frequency and precision for xtc file");
825 ITYPE ("nstxtcout", ir->nstxtcout, 0);
826 RTYPE ("xtc-precision",ir->xtcprec, 1000.0);
827 CTYPE ("This selects the subset of atoms for the xtc file. You can");
828 CTYPE ("select multiple groups. By default all atoms will be written.");
829 STYPE ("xtc-grps", xtc_grps, NULL);
830 CTYPE ("Selection of energy groups");
831 STYPE ("energygrps", energy, NULL);
833 /* Neighbor searching */
834 CCTYPE ("NEIGHBORSEARCHING PARAMETERS");
835 CTYPE ("nblist update frequency");
836 ITYPE ("nstlist", ir->nstlist, 10);
837 CTYPE ("ns algorithm (simple or grid)");
838 EETYPE("ns-type", ir->ns_type, ens_names);
839 /* set ndelta to the optimal value of 2 */
840 ir->ndelta = 2;
841 CTYPE ("Periodic boundary conditions: xyz, no, xy");
842 EETYPE("pbc", ir->ePBC, epbc_names);
843 EETYPE("periodic_molecules", ir->bPeriodicMols, yesno_names);
844 CTYPE ("nblist cut-off");
845 RTYPE ("rlist", ir->rlist, 1.0);
846 CTYPE ("long-range cut-off for switched potentials");
847 RTYPE ("rlistlong", ir->rlistlong, -1);
849 /* Electrostatics */
850 CCTYPE ("OPTIONS FOR ELECTROSTATICS AND VDW");
851 CTYPE ("Method for doing electrostatics");
852 EETYPE("coulombtype", ir->coulombtype, eel_names);
853 CTYPE ("cut-off lengths");
854 RTYPE ("rcoulomb-switch", ir->rcoulomb_switch, 0.0);
855 RTYPE ("rcoulomb", ir->rcoulomb, 1.0);
856 CTYPE ("Relative dielectric constant for the medium and the reaction field");
857 RTYPE ("epsilon_r", ir->epsilon_r, 1.0);
858 RTYPE ("epsilon_rf", ir->epsilon_rf, 1.0);
859 CTYPE ("Method for doing Van der Waals");
860 EETYPE("vdw-type", ir->vdwtype, evdw_names);
861 CTYPE ("cut-off lengths");
862 RTYPE ("rvdw-switch", ir->rvdw_switch, 0.0);
863 RTYPE ("rvdw", ir->rvdw, 1.0);
864 CTYPE ("Apply long range dispersion corrections for Energy and Pressure");
865 EETYPE("DispCorr", ir->eDispCorr, edispc_names);
866 CTYPE ("Extension of the potential lookup tables beyond the cut-off");
867 RTYPE ("table-extension", ir->tabext, 1.0);
868 CTYPE ("Seperate tables between energy group pairs");
869 STYPE ("energygrp_table", egptable, NULL);
870 CTYPE ("Spacing for the PME/PPPM FFT grid");
871 RTYPE ("fourierspacing", opts->fourierspacing,0.12);
872 CTYPE ("FFT grid size, when a value is 0 fourierspacing will be used");
873 ITYPE ("fourier_nx", ir->nkx, 0);
874 ITYPE ("fourier_ny", ir->nky, 0);
875 ITYPE ("fourier_nz", ir->nkz, 0);
876 CTYPE ("EWALD/PME/PPPM parameters");
877 ITYPE ("pme_order", ir->pme_order, 4);
878 RTYPE ("ewald_rtol", ir->ewald_rtol, 0.00001);
879 EETYPE("ewald_geometry", ir->ewald_geometry, eewg_names);
880 RTYPE ("epsilon_surface", ir->epsilon_surface, 0.0);
881 EETYPE("optimize_fft",ir->bOptFFT, yesno_names);
883 CCTYPE("IMPLICIT SOLVENT ALGORITHM");
884 EETYPE("implicit_solvent", ir->implicit_solvent, eis_names);
886 CCTYPE ("GENERALIZED BORN ELECTROSTATICS");
887 CTYPE ("Algorithm for calculating Born radii");
888 EETYPE("gb_algorithm", ir->gb_algorithm, egb_names);
889 CTYPE ("Frequency of calculating the Born radii inside rlist");
890 ITYPE ("nstgbradii", ir->nstgbradii, 1);
891 CTYPE ("Cutoff for Born radii calculation; the contribution from atoms");
892 CTYPE ("between rlist and rgbradii is updated every nstlist steps");
893 RTYPE ("rgbradii", ir->rgbradii, 1.0);
894 CTYPE ("Dielectric coefficient of the implicit solvent");
895 RTYPE ("gb_epsilon_solvent",ir->gb_epsilon_solvent, 80.0);
896 CTYPE ("Salt concentration in M for Generalized Born models");
897 RTYPE ("gb_saltconc", ir->gb_saltconc, 0.0);
898 CTYPE ("Scaling factors used in the OBC GB model. Default values are OBC(II)");
899 RTYPE ("gb_obc_alpha", ir->gb_obc_alpha, 1.0);
900 RTYPE ("gb_obc_beta", ir->gb_obc_beta, 0.8);
901 RTYPE ("gb_obc_gamma", ir->gb_obc_gamma, 4.85);
902 RTYPE ("gb_dielectric_offset", ir->gb_dielectric_offset, 0.009);
903 EETYPE("sa_algorithm", ir->sa_algorithm, esa_names);
904 CTYPE ("Surface tension (kJ/mol/nm^2) for the SA (nonpolar surface) part of GBSA");
905 CTYPE ("The default value (2.092) corresponds to 0.005 kcal/mol/Angstrom^2.");
906 RTYPE ("sa_surface_tension", ir->sa_surface_tension, 2.092);
908 /* Coupling stuff */
909 CCTYPE ("OPTIONS FOR WEAK COUPLING ALGORITHMS");
910 CTYPE ("Temperature coupling");
911 EETYPE("tcoupl", ir->etc, etcoupl_names);
912 ITYPE("nh-chain-length", ir->opts.nhchainlength, NHCHAINLENGTH);
913 CTYPE ("Groups to couple separately");
914 STYPE ("tc-grps", tcgrps, NULL);
915 CTYPE ("Time constant (ps) and reference temperature (K)");
916 STYPE ("tau-t", tau_t, NULL);
917 STYPE ("ref-t", ref_t, NULL);
918 CTYPE ("Pressure coupling");
919 EETYPE("Pcoupl", ir->epc, epcoupl_names);
920 EETYPE("Pcoupltype", ir->epct, epcoupltype_names);
921 CTYPE ("Time constant (ps), compressibility (1/bar) and reference P (bar)");
922 RTYPE ("tau-p", ir->tau_p, 1.0);
923 STYPE ("compressibility", dumstr[0], NULL);
924 STYPE ("ref-p", dumstr[1], NULL);
925 CTYPE ("Scaling of reference coordinates, No, All or COM");
926 EETYPE ("refcoord_scaling",ir->refcoord_scaling,erefscaling_names);
928 CTYPE ("Random seed for Andersen thermostat");
929 ITYPE ("andersen_seed", ir->andersen_seed, 815131);
931 /* QMMM */
932 CCTYPE ("OPTIONS FOR QMMM calculations");
933 EETYPE("QMMM", ir->bQMMM, yesno_names);
934 CTYPE ("Groups treated Quantum Mechanically");
935 STYPE ("QMMM-grps", QMMM, NULL);
936 CTYPE ("QM method");
937 STYPE("QMmethod", QMmethod, NULL);
938 CTYPE ("QMMM scheme");
939 EETYPE("QMMMscheme", ir->QMMMscheme, eQMMMscheme_names);
940 CTYPE ("QM basisset");
941 STYPE("QMbasis", QMbasis, NULL);
942 CTYPE ("QM charge");
943 STYPE ("QMcharge", QMcharge,NULL);
944 CTYPE ("QM multiplicity");
945 STYPE ("QMmult", QMmult,NULL);
946 CTYPE ("Surface Hopping");
947 STYPE ("SH", bSH, NULL);
948 CTYPE ("CAS space options");
949 STYPE ("CASorbitals", CASorbitals, NULL);
950 STYPE ("CASelectrons", CASelectrons, NULL);
951 STYPE ("SAon", SAon, NULL);
952 STYPE ("SAoff",SAoff,NULL);
953 STYPE ("SAsteps", SAsteps, NULL);
954 CTYPE ("Scale factor for MM charges");
955 RTYPE ("MMChargeScaleFactor", ir->scalefactor, 1.0);
956 CTYPE ("Optimization of QM subsystem");
957 STYPE ("bOPT", bOPT, NULL);
958 STYPE ("bTS", bTS, NULL);
960 /* Simulated annealing */
961 CCTYPE("SIMULATED ANNEALING");
962 CTYPE ("Type of annealing for each temperature group (no/single/periodic)");
963 STYPE ("annealing", anneal, NULL);
964 CTYPE ("Number of time points to use for specifying annealing in each group");
965 STYPE ("annealing_npoints", anneal_npoints, NULL);
966 CTYPE ("List of times at the annealing points for each group");
967 STYPE ("annealing_time", anneal_time, NULL);
968 CTYPE ("Temp. at each annealing point, for each group.");
969 STYPE ("annealing_temp", anneal_temp, NULL);
971 /* Startup run */
972 CCTYPE ("GENERATE VELOCITIES FOR STARTUP RUN");
973 EETYPE("gen-vel", opts->bGenVel, yesno_names);
974 RTYPE ("gen-temp", opts->tempi, 300.0);
975 ITYPE ("gen-seed", opts->seed, 173529);
977 /* Shake stuff */
978 CCTYPE ("OPTIONS FOR BONDS");
979 EETYPE("constraints", opts->nshake, constraints);
980 CTYPE ("Type of constraint algorithm");
981 EETYPE("constraint-algorithm", ir->eConstrAlg, econstr_names);
982 CTYPE ("Do not constrain the start configuration");
983 EETYPE("continuation", ir->bContinuation, yesno_names);
984 CTYPE ("Use successive overrelaxation to reduce the number of shake iterations");
985 EETYPE("Shake-SOR", ir->bShakeSOR, yesno_names);
986 CTYPE ("Relative tolerance of shake");
987 RTYPE ("shake-tol", ir->shake_tol, 0.0001);
988 CTYPE ("Highest order in the expansion of the constraint coupling matrix");
989 ITYPE ("lincs-order", ir->nProjOrder, 4);
990 CTYPE ("Number of iterations in the final step of LINCS. 1 is fine for");
991 CTYPE ("normal simulations, but use 2 to conserve energy in NVE runs.");
992 CTYPE ("For energy minimization with constraints it should be 4 to 8.");
993 ITYPE ("lincs-iter", ir->nLincsIter, 1);
994 CTYPE ("Lincs will write a warning to the stderr if in one step a bond");
995 CTYPE ("rotates over more degrees than");
996 RTYPE ("lincs-warnangle", ir->LincsWarnAngle, 30.0);
997 CTYPE ("Convert harmonic bonds to morse potentials");
998 EETYPE("morse", opts->bMorse,yesno_names);
1000 /* Energy group exclusions */
1001 CCTYPE ("ENERGY GROUP EXCLUSIONS");
1002 CTYPE ("Pairs of energy groups for which all non-bonded interactions are excluded");
1003 STYPE ("energygrp_excl", egpexcl, NULL);
1005 /* Walls */
1006 CCTYPE ("WALLS");
1007 CTYPE ("Number of walls, type, atom types, densities and box-z scale factor for Ewald");
1008 ITYPE ("nwall", ir->nwall, 0);
1009 EETYPE("wall_type", ir->wall_type, ewt_names);
1010 RTYPE ("wall_r_linpot", ir->wall_r_linpot, -1);
1011 STYPE ("wall_atomtype", wall_atomtype, NULL);
1012 STYPE ("wall_density", wall_density, NULL);
1013 RTYPE ("wall_ewald_zfac", ir->wall_ewald_zfac, 3);
1015 /* COM pulling */
1016 CCTYPE("COM PULLING");
1017 CTYPE("Pull type: no, umbrella, constraint or constant_force");
1018 EETYPE("pull", ir->ePull, epull_names);
1019 if (ir->ePull != epullNO) {
1020 snew(ir->pull,1);
1021 pull_grp = read_pullparams(&ninp,&inp,ir->pull,&opts->pull_start,wi);
1024 /* Refinement */
1025 CCTYPE("NMR refinement stuff");
1026 CTYPE ("Distance restraints type: No, Simple or Ensemble");
1027 EETYPE("disre", ir->eDisre, edisre_names);
1028 CTYPE ("Force weighting of pairs in one distance restraint: Conservative or Equal");
1029 EETYPE("disre-weighting", ir->eDisreWeighting, edisreweighting_names);
1030 CTYPE ("Use sqrt of the time averaged times the instantaneous violation");
1031 EETYPE("disre-mixed", ir->bDisreMixed, yesno_names);
1032 RTYPE ("disre-fc", ir->dr_fc, 1000.0);
1033 RTYPE ("disre-tau", ir->dr_tau, 0.0);
1034 CTYPE ("Output frequency for pair distances to energy file");
1035 ITYPE ("nstdisreout", ir->nstdisreout, 100);
1036 CTYPE ("Orientation restraints: No or Yes");
1037 EETYPE("orire", opts->bOrire, yesno_names);
1038 CTYPE ("Orientation restraints force constant and tau for time averaging");
1039 RTYPE ("orire-fc", ir->orires_fc, 0.0);
1040 RTYPE ("orire-tau", ir->orires_tau, 0.0);
1041 STYPE ("orire-fitgrp",orirefitgrp, NULL);
1042 CTYPE ("Output frequency for trace(SD) and S to energy file");
1043 ITYPE ("nstorireout", ir->nstorireout, 100);
1044 CTYPE ("Dihedral angle restraints: No or Yes");
1045 EETYPE("dihre", opts->bDihre, yesno_names);
1046 RTYPE ("dihre-fc", ir->dihre_fc, 1000.0);
1048 /* Free energy stuff */
1049 CCTYPE ("Free energy control stuff");
1050 EETYPE("free-energy", ir->efep, efep_names);
1051 RTYPE ("init-lambda", ir->init_lambda,0.0);
1052 RTYPE ("delta-lambda",ir->delta_lambda,0.0);
1053 STYPE ("foreign_lambda", foreign_lambda, NULL);
1054 RTYPE ("sc-alpha",ir->sc_alpha,0.0);
1055 ITYPE ("sc-power",ir->sc_power,0);
1056 RTYPE ("sc-sigma",ir->sc_sigma,0.3);
1057 ITYPE ("nstdhdl", ir->nstdhdl, 10);
1058 STYPE ("couple-moltype", couple_moltype, NULL);
1059 EETYPE("couple-lambda0", opts->couple_lam0, couple_lam);
1060 EETYPE("couple-lambda1", opts->couple_lam1, couple_lam);
1061 EETYPE("couple-intramol", opts->bCoupleIntra, yesno_names);
1063 /* Non-equilibrium MD stuff */
1064 CCTYPE("Non-equilibrium MD stuff");
1065 STYPE ("acc-grps", accgrps, NULL);
1066 STYPE ("accelerate", acc, NULL);
1067 STYPE ("freezegrps", freeze, NULL);
1068 STYPE ("freezedim", frdim, NULL);
1069 RTYPE ("cos-acceleration", ir->cos_accel, 0);
1070 STYPE ("deform", deform, NULL);
1072 /* Electric fields */
1073 CCTYPE("Electric fields");
1074 CTYPE ("Format is number of terms (int) and for all terms an amplitude (real)");
1075 CTYPE ("and a phase angle (real)");
1076 STYPE ("E-x", efield_x, NULL);
1077 STYPE ("E-xt", efield_xt, NULL);
1078 STYPE ("E-y", efield_y, NULL);
1079 STYPE ("E-yt", efield_yt, NULL);
1080 STYPE ("E-z", efield_z, NULL);
1081 STYPE ("E-zt", efield_zt, NULL);
1083 /* AdResS defined thingies */
1084 CCTYPE ("AdResS parameters");
1085 EETYPE("adress_type", ir->adress_type, eAdresstype_names);
1086 EETYPE("adress_new_wf", ir->badress_new_wf, yesno_names);
1087 EETYPE("adress_chempot_dx", ir->badress_chempot_dx, yesno_names);
1088 EETYPE("adress_tf_full_box", ir->badress_tf_full_box, yesno_names);
1089 RTYPE ("adress_const_wf", ir->adress_const_wf, 1);
1090 RTYPE ("adress_ex_width", ir->adress_ex_width, 0);
1091 RTYPE ("adress_hy_width", ir->adress_hy_width, 0);
1092 EETYPE("adress_interface_correction",ir->adress_icor, eAdressICtype_names);
1093 EETYPE("adress_exvdw", ir->adress_ivdw, evdw_names);
1094 EETYPE("adress_site", ir->adress_site, eAdressSITEtype_names);
1095 STYPE ("adress_reference_coords", adress_refs, NULL);
1096 STYPE ("adress_tf_grp_names", adress_tf_grp_names, NULL);
1098 /* User defined thingies */
1099 CCTYPE ("User defined thingies");
1100 STYPE ("user1-grps", user1, NULL);
1101 STYPE ("user2-grps", user2, NULL);
1102 ITYPE ("userint1", ir->userint1, 0);
1103 ITYPE ("userint2", ir->userint2, 0);
1104 ITYPE ("userint3", ir->userint3, 0);
1105 ITYPE ("userint4", ir->userint4, 0);
1106 RTYPE ("userreal1", ir->userreal1, 0);
1107 RTYPE ("userreal2", ir->userreal2, 0);
1108 RTYPE ("userreal3", ir->userreal3, 0);
1109 RTYPE ("userreal4", ir->userreal4, 0);
1110 #undef CTYPE
1112 write_inpfile(mdparout,ninp,inp,FALSE,wi);
1113 for (i=0; (i<ninp); i++) {
1114 sfree(inp[i].name);
1115 sfree(inp[i].value);
1117 sfree(inp);
1119 /* Process options if necessary */
1120 for(m=0; m<2; m++) {
1121 for(i=0; i<2*DIM; i++)
1122 dumdub[m][i]=0.0;
1123 if(ir->epc) {
1124 switch (ir->epct) {
1125 case epctISOTROPIC:
1126 if (sscanf(dumstr[m],"%lf",&(dumdub[m][XX]))!=1) {
1127 warning_error(wi,"Pressure coupling not enough values (I need 1)");
1129 dumdub[m][YY]=dumdub[m][ZZ]=dumdub[m][XX];
1130 break;
1131 case epctSEMIISOTROPIC:
1132 case epctSURFACETENSION:
1133 if (sscanf(dumstr[m],"%lf%lf",
1134 &(dumdub[m][XX]),&(dumdub[m][ZZ]))!=2) {
1135 warning_error(wi,"Pressure coupling not enough values (I need 2)");
1137 dumdub[m][YY]=dumdub[m][XX];
1138 break;
1139 case epctANISOTROPIC:
1140 if (sscanf(dumstr[m],"%lf%lf%lf%lf%lf%lf",
1141 &(dumdub[m][XX]),&(dumdub[m][YY]),&(dumdub[m][ZZ]),
1142 &(dumdub[m][3]),&(dumdub[m][4]),&(dumdub[m][5]))!=6) {
1143 warning_error(wi,"Pressure coupling not enough values (I need 6)");
1145 break;
1146 default:
1147 gmx_fatal(FARGS,"Pressure coupling type %s not implemented yet",
1148 epcoupltype_names[ir->epct]);
1152 clear_mat(ir->ref_p);
1153 clear_mat(ir->compress);
1154 for(i=0; i<DIM; i++) {
1155 ir->ref_p[i][i] = dumdub[1][i];
1156 ir->compress[i][i] = dumdub[0][i];
1158 if (ir->epct == epctANISOTROPIC) {
1159 ir->ref_p[XX][YY] = dumdub[1][3];
1160 ir->ref_p[XX][ZZ] = dumdub[1][4];
1161 ir->ref_p[YY][ZZ] = dumdub[1][5];
1162 if (ir->ref_p[XX][YY]!=0 && ir->ref_p[XX][ZZ]!=0 && ir->ref_p[YY][ZZ]!=0) {
1163 warning(wi,"All off-diagonal reference pressures are non-zero. Are you sure you want to apply a threefold shear stress?\n");
1165 ir->compress[XX][YY] = dumdub[0][3];
1166 ir->compress[XX][ZZ] = dumdub[0][4];
1167 ir->compress[YY][ZZ] = dumdub[0][5];
1168 for(i=0; i<DIM; i++) {
1169 for(m=0; m<i; m++) {
1170 ir->ref_p[i][m] = ir->ref_p[m][i];
1171 ir->compress[i][m] = ir->compress[m][i];
1176 if (ir->comm_mode == ecmNO)
1177 ir->nstcomm = 0;
1179 opts->couple_moltype = NULL;
1180 if (strlen(couple_moltype) > 0) {
1181 if (ir->efep != efepNO) {
1182 opts->couple_moltype = strdup(couple_moltype);
1183 if (opts->couple_lam0 == opts->couple_lam1)
1184 warning(wi,"The lambda=0 and lambda=1 states for coupling are identical");
1185 if (ir->eI == eiMD && (opts->couple_lam0 == ecouplamNONE ||
1186 opts->couple_lam1 == ecouplamNONE)) {
1187 warning(wi,"For proper sampling of the (nearly) decoupled state, stochastic dynamics should be used");
1189 } else {
1190 warning(wi,"Can not couple a molecule with free_energy = no");
1194 do_wall_params(ir,wall_atomtype,wall_density,opts);
1196 if (opts->bOrire && str_nelem(orirefitgrp,MAXPTR,NULL)!=1) {
1197 warning_error(wi,"ERROR: Need one orientation restraint fit group\n");
1200 clear_mat(ir->deform);
1201 for(i=0; i<6; i++)
1202 dumdub[0][i] = 0;
1203 m = sscanf(deform,"%lf %lf %lf %lf %lf %lf",
1204 &(dumdub[0][0]),&(dumdub[0][1]),&(dumdub[0][2]),
1205 &(dumdub[0][3]),&(dumdub[0][4]),&(dumdub[0][5]));
1206 for(i=0; i<3; i++)
1207 ir->deform[i][i] = dumdub[0][i];
1208 ir->deform[YY][XX] = dumdub[0][3];
1209 ir->deform[ZZ][XX] = dumdub[0][4];
1210 ir->deform[ZZ][YY] = dumdub[0][5];
1211 if (ir->epc != epcNO) {
1212 for(i=0; i<3; i++)
1213 for(j=0; j<=i; j++)
1214 if (ir->deform[i][j]!=0 && ir->compress[i][j]!=0) {
1215 warning_error(wi,"A box element has deform set and compressibility > 0");
1217 for(i=0; i<3; i++)
1218 for(j=0; j<i; j++)
1219 if (ir->deform[i][j]!=0) {
1220 for(m=j; m<DIM; m++)
1221 if (ir->compress[m][j]!=0) {
1222 sprintf(warn_buf,"An off-diagonal box element has deform set while compressibility > 0 for the same component of another box vector, this might lead to spurious periodicity effects.");
1223 warning(wi,warn_buf);
1228 if (ir->efep != efepNO) {
1229 parse_n_double(foreign_lambda,&ir->n_flambda,&ir->flambda);
1230 if (ir->n_flambda > 0 && ir->rlist < max(ir->rvdw,ir->rcoulomb)) {
1231 warning_note(wi,"For foreign lambda free energy differences it is assumed that the soft-core interactions have no effect beyond the neighborlist cut-off");
1233 } else {
1234 ir->n_flambda = 0;
1237 sfree(dumstr[0]);
1238 sfree(dumstr[1]);
1241 static int search_QMstring(char *s,int ng,const char *gn[])
1243 /* same as normal search_string, but this one searches QM strings */
1244 int i;
1246 for(i=0; (i<ng); i++)
1247 if (strcasecmp(s,gn[i]) == 0)
1248 return i;
1250 gmx_fatal(FARGS,"this QM method or basisset (%s) is not implemented\n!",s);
1252 return -1;
1254 } /* search_QMstring */
1257 int search_string(char *s,int ng,char *gn[])
1259 int i;
1261 for(i=0; (i<ng); i++)
1262 if (strcasecmp(s,gn[i]) == 0)
1263 return i;
1265 gmx_fatal(FARGS,"Group %s not found in indexfile.\nMaybe you have non-default goups in your mdp file, while not using the '-n' option of grompp.\nIn that case use the '-n' option.\n",s);
1267 return -1;
1270 static bool do_numbering(int natoms,gmx_groups_t *groups,int ng,char *ptrs[],
1271 t_blocka *block,char *gnames[],
1272 int gtype,int restnm,
1273 int grptp,bool bVerbose,
1274 warninp_t wi)
1276 unsigned short *cbuf;
1277 t_grps *grps=&(groups->grps[gtype]);
1278 int i,j,gid,aj,ognr,ntot=0;
1279 const char *title;
1280 bool bRest;
1281 char warn_buf[STRLEN];
1283 if (debug)
1284 fprintf(debug,"Starting numbering %d groups of type %d\n",ng,gtype);
1286 title = gtypes[gtype];
1288 snew(cbuf,natoms);
1289 for(i=0; (i<natoms); i++)
1290 cbuf[i]=NOGID;
1292 snew(grps->nm_ind,ng+1); /* +1 for possible rest group */
1293 for(i=0; (i<ng); i++) {
1294 /* Lookup the group name in the block structure */
1295 gid = search_string(ptrs[i],block->nr,gnames);
1296 if ((grptp != egrptpONE) || (i == 0))
1297 grps->nm_ind[grps->nr++]=gid;
1298 if (debug)
1299 fprintf(debug,"Found gid %d for group %s\n",gid,ptrs[i]);
1301 /* Now go over the atoms in the group */
1302 for(j=block->index[gid]; (j<block->index[gid+1]); j++) {
1303 aj=block->a[j];
1305 /* Range checking */
1306 if ((aj < 0) || (aj >= natoms))
1307 gmx_fatal(FARGS,"Invalid atom number %d in indexfile",aj);
1309 /* Lookup up the old group number */
1310 ognr = cbuf[aj];
1311 if (ognr != NOGID)
1312 gmx_fatal(FARGS,"Atom %d in multiple %s groups (%d and %d)",
1313 aj+1,title,ognr+1,i+1);
1314 else {
1315 /* Store the group number in buffer */
1316 if (grptp == egrptpONE)
1317 cbuf[aj] = 0;
1318 else
1319 cbuf[aj] = i;
1320 ntot++;
1325 /* Now check whether we have done all atoms */
1326 bRest = FALSE;
1327 if (ntot != natoms) {
1328 if (grptp == egrptpALL) {
1329 gmx_fatal(FARGS,"%d atoms are not part of any of the %s groups",
1330 natoms-ntot,title);
1331 } else if (grptp == egrptpPART) {
1332 sprintf(warn_buf,"%d atoms are not part of any of the %s groups",
1333 natoms-ntot,title);
1334 warning_note(wi,warn_buf);
1336 /* Assign all atoms currently unassigned to a rest group */
1337 for(j=0; (j<natoms); j++) {
1338 if (cbuf[j] == NOGID) {
1339 cbuf[j] = grps->nr;
1340 bRest = TRUE;
1343 if (grptp != egrptpPART) {
1344 if (bVerbose)
1345 fprintf(stderr,
1346 "Making dummy/rest group for %s containing %d elements\n",
1347 title,natoms-ntot);
1348 /* Add group name "rest" */
1349 grps->nm_ind[grps->nr] = restnm;
1351 /* Assign the rest name to all atoms not currently assigned to a group */
1352 for(j=0; (j<natoms); j++) {
1353 if (cbuf[j] == NOGID)
1354 cbuf[j] = grps->nr;
1356 grps->nr++;
1360 if (grps->nr == 1) {
1361 groups->ngrpnr[gtype] = 0;
1362 groups->grpnr[gtype] = NULL;
1363 } else {
1364 groups->ngrpnr[gtype] = natoms;
1365 snew(groups->grpnr[gtype],natoms);
1366 for(j=0; (j<natoms); j++) {
1367 groups->grpnr[gtype][j] = cbuf[j];
1371 sfree(cbuf);
1373 return (bRest && grptp == egrptpPART);
1376 static void calc_nrdf(gmx_mtop_t *mtop,t_inputrec *ir,char **gnames)
1378 t_grpopts *opts;
1379 gmx_groups_t *groups;
1380 t_pull *pull;
1381 int natoms,ai,aj,i,j,d,g,imin,jmin,nc;
1382 t_iatom *ia;
1383 int *nrdf2,*na_vcm,na_tot;
1384 double *nrdf_tc,*nrdf_vcm,nrdf_uc,n_sub=0;
1385 gmx_mtop_atomloop_all_t aloop;
1386 t_atom *atom;
1387 int mb,mol,ftype,as;
1388 gmx_molblock_t *molb;
1389 gmx_moltype_t *molt;
1391 /* Calculate nrdf.
1392 * First calc 3xnr-atoms for each group
1393 * then subtract half a degree of freedom for each constraint
1395 * Only atoms and nuclei contribute to the degrees of freedom...
1398 opts = &ir->opts;
1400 groups = &mtop->groups;
1401 natoms = mtop->natoms;
1403 /* Allocate one more for a possible rest group */
1404 /* We need to sum degrees of freedom into doubles,
1405 * since floats give too low nrdf's above 3 million atoms.
1407 snew(nrdf_tc,groups->grps[egcTC].nr+1);
1408 snew(nrdf_vcm,groups->grps[egcVCM].nr+1);
1409 snew(na_vcm,groups->grps[egcVCM].nr+1);
1411 for(i=0; i<groups->grps[egcTC].nr; i++)
1412 nrdf_tc[i] = 0;
1413 for(i=0; i<groups->grps[egcVCM].nr+1; i++)
1414 nrdf_vcm[i] = 0;
1416 snew(nrdf2,natoms);
1417 aloop = gmx_mtop_atomloop_all_init(mtop);
1418 while (gmx_mtop_atomloop_all_next(aloop,&i,&atom)) {
1419 nrdf2[i] = 0;
1420 if (atom->ptype == eptAtom || atom->ptype == eptNucleus) {
1421 g = ggrpnr(groups,egcFREEZE,i);
1422 /* Double count nrdf for particle i */
1423 for(d=0; d<DIM; d++) {
1424 if (opts->nFreeze[g][d] == 0) {
1425 nrdf2[i] += 2;
1428 nrdf_tc [ggrpnr(groups,egcTC ,i)] += 0.5*nrdf2[i];
1429 nrdf_vcm[ggrpnr(groups,egcVCM,i)] += 0.5*nrdf2[i];
1433 as = 0;
1434 for(mb=0; mb<mtop->nmolblock; mb++) {
1435 molb = &mtop->molblock[mb];
1436 molt = &mtop->moltype[molb->type];
1437 atom = molt->atoms.atom;
1438 for(mol=0; mol<molb->nmol; mol++) {
1439 for (ftype=F_CONSTR; ftype<=F_CONSTRNC; ftype++) {
1440 ia = molt->ilist[ftype].iatoms;
1441 for(i=0; i<molt->ilist[ftype].nr; ) {
1442 /* Subtract degrees of freedom for the constraints,
1443 * if the particles still have degrees of freedom left.
1444 * If one of the particles is a vsite or a shell, then all
1445 * constraint motion will go there, but since they do not
1446 * contribute to the constraints the degrees of freedom do not
1447 * change.
1449 ai = as + ia[1];
1450 aj = as + ia[2];
1451 if (((atom[ia[1]].ptype == eptNucleus) ||
1452 (atom[ia[1]].ptype == eptAtom)) &&
1453 ((atom[ia[2]].ptype == eptNucleus) ||
1454 (atom[ia[2]].ptype == eptAtom))) {
1455 if (nrdf2[ai] > 0)
1456 jmin = 1;
1457 else
1458 jmin = 2;
1459 if (nrdf2[aj] > 0)
1460 imin = 1;
1461 else
1462 imin = 2;
1463 imin = min(imin,nrdf2[ai]);
1464 jmin = min(jmin,nrdf2[aj]);
1465 nrdf2[ai] -= imin;
1466 nrdf2[aj] -= jmin;
1467 nrdf_tc [ggrpnr(groups,egcTC ,ai)] -= 0.5*imin;
1468 nrdf_tc [ggrpnr(groups,egcTC ,aj)] -= 0.5*jmin;
1469 nrdf_vcm[ggrpnr(groups,egcVCM,ai)] -= 0.5*imin;
1470 nrdf_vcm[ggrpnr(groups,egcVCM,aj)] -= 0.5*jmin;
1472 ia += interaction_function[ftype].nratoms+1;
1473 i += interaction_function[ftype].nratoms+1;
1476 ia = molt->ilist[F_SETTLE].iatoms;
1477 for(i=0; i<molt->ilist[F_SETTLE].nr; ) {
1478 /* Subtract 1 dof from every atom in the SETTLE */
1479 for(ai=as+ia[1]; ai<as+ia[1]+3; ai++) {
1480 imin = min(2,nrdf2[ai]);
1481 nrdf2[ai] -= imin;
1482 nrdf_tc [ggrpnr(groups,egcTC ,ai)] -= 0.5*imin;
1483 nrdf_vcm[ggrpnr(groups,egcVCM,ai)] -= 0.5*imin;
1485 ia += 2;
1486 i += 2;
1488 as += molt->atoms.nr;
1492 if (ir->ePull == epullCONSTRAINT) {
1493 /* Correct nrdf for the COM constraints.
1494 * We correct using the TC and VCM group of the first atom
1495 * in the reference and pull group. If atoms in one pull group
1496 * belong to different TC or VCM groups it is anyhow difficult
1497 * to determine the optimal nrdf assignment.
1499 pull = ir->pull;
1500 if (pull->eGeom == epullgPOS) {
1501 nc = 0;
1502 for(i=0; i<DIM; i++) {
1503 if (pull->dim[i])
1504 nc++;
1506 } else {
1507 nc = 1;
1509 for(i=0; i<pull->ngrp; i++) {
1510 imin = 2*nc;
1511 if (pull->grp[0].nat > 0) {
1512 /* Subtract 1/2 dof from the reference group */
1513 ai = pull->grp[0].ind[0];
1514 if (nrdf_tc[ggrpnr(groups,egcTC,ai)] > 1) {
1515 nrdf_tc [ggrpnr(groups,egcTC ,ai)] -= 0.5;
1516 nrdf_vcm[ggrpnr(groups,egcVCM,ai)] -= 0.5;
1517 imin--;
1520 /* Subtract 1/2 dof from the pulled group */
1521 ai = pull->grp[1+i].ind[0];
1522 nrdf_tc [ggrpnr(groups,egcTC ,ai)] -= 0.5*imin;
1523 nrdf_vcm[ggrpnr(groups,egcVCM,ai)] -= 0.5*imin;
1524 if (nrdf_tc[ggrpnr(groups,egcTC,ai)] < 0)
1525 gmx_fatal(FARGS,"Center of mass pulling constraints caused the number of degrees of freedom for temperature coupling group %s to be negative",gnames[groups->grps[egcTC].nm_ind[ggrpnr(groups,egcTC,ai)]]);
1529 if (ir->nstcomm != 0) {
1530 /* Subtract 3 from the number of degrees of freedom in each vcm group
1531 * when com translation is removed and 6 when rotation is removed
1532 * as well.
1534 switch (ir->comm_mode) {
1535 case ecmLINEAR:
1536 n_sub = ndof_com(ir);
1537 break;
1538 case ecmANGULAR:
1539 n_sub = 6;
1540 break;
1541 default:
1542 n_sub = 0;
1543 gmx_incons("Checking comm_mode");
1546 for(i=0; i<groups->grps[egcTC].nr; i++) {
1547 /* Count the number of atoms of TC group i for every VCM group */
1548 for(j=0; j<groups->grps[egcVCM].nr+1; j++)
1549 na_vcm[j] = 0;
1550 na_tot = 0;
1551 for(ai=0; ai<natoms; ai++)
1552 if (ggrpnr(groups,egcTC,ai) == i) {
1553 na_vcm[ggrpnr(groups,egcVCM,ai)]++;
1554 na_tot++;
1556 /* Correct for VCM removal according to the fraction of each VCM
1557 * group present in this TC group.
1559 nrdf_uc = nrdf_tc[i];
1560 if (debug) {
1561 fprintf(debug,"T-group[%d] nrdf_uc = %g, n_sub = %g\n",
1562 i,nrdf_uc,n_sub);
1564 nrdf_tc[i] = 0;
1565 for(j=0; j<groups->grps[egcVCM].nr+1; j++) {
1566 if (nrdf_vcm[j] > n_sub) {
1567 nrdf_tc[i] += nrdf_uc*((double)na_vcm[j]/(double)na_tot)*
1568 (nrdf_vcm[j] - n_sub)/nrdf_vcm[j];
1570 if (debug) {
1571 fprintf(debug," nrdf_vcm[%d] = %g, nrdf = %g\n",
1572 j,nrdf_vcm[j],nrdf_tc[i]);
1577 for(i=0; (i<groups->grps[egcTC].nr); i++) {
1578 opts->nrdf[i] = nrdf_tc[i];
1579 if (opts->nrdf[i] < 0)
1580 opts->nrdf[i] = 0;
1581 fprintf(stderr,
1582 "Number of degrees of freedom in T-Coupling group %s is %.2f\n",
1583 gnames[groups->grps[egcTC].nm_ind[i]],opts->nrdf[i]);
1586 sfree(nrdf2);
1587 sfree(nrdf_tc);
1588 sfree(nrdf_vcm);
1589 sfree(na_vcm);
1592 static void decode_cos(char *s,t_cosines *cosine,bool bTime)
1594 char *t;
1595 char format[STRLEN],f1[STRLEN];
1596 double a,phi;
1597 int i;
1599 t=strdup(s);
1600 trim(t);
1602 cosine->n=0;
1603 cosine->a=NULL;
1604 cosine->phi=NULL;
1605 if (strlen(t)) {
1606 sscanf(t,"%d",&(cosine->n));
1607 if (cosine->n <= 0) {
1608 cosine->n=0;
1609 } else {
1610 snew(cosine->a,cosine->n);
1611 snew(cosine->phi,cosine->n);
1613 sprintf(format,"%%*d");
1614 for(i=0; (i<cosine->n); i++) {
1615 strcpy(f1,format);
1616 strcat(f1,"%lf%lf");
1617 if (sscanf(t,f1,&a,&phi) < 2)
1618 gmx_fatal(FARGS,"Invalid input for electric field shift: '%s'",t);
1619 cosine->a[i]=a;
1620 cosine->phi[i]=phi;
1621 strcat(format,"%*lf%*lf");
1625 sfree(t);
1628 static bool do_egp_flag(t_inputrec *ir,gmx_groups_t *groups,
1629 const char *option,const char *val,int flag)
1631 /* The maximum number of energy group pairs would be MAXPTR*(MAXPTR+1)/2.
1632 * But since this is much larger than STRLEN, such a line can not be parsed.
1633 * The real maximum is the number of names that fit in a string: STRLEN/2.
1635 #define EGP_MAX (STRLEN/2)
1636 int nelem,i,j,k,nr;
1637 char *names[EGP_MAX];
1638 char ***gnames;
1639 bool bSet;
1641 gnames = groups->grpname;
1643 nelem = str_nelem(val,EGP_MAX,names);
1644 if (nelem % 2 != 0)
1645 gmx_fatal(FARGS,"The number of groups for %s is odd",option);
1646 nr = groups->grps[egcENER].nr;
1647 bSet = FALSE;
1648 for(i=0; i<nelem/2; i++) {
1649 j = 0;
1650 while ((j < nr) &&
1651 strcasecmp(names[2*i],*(gnames[groups->grps[egcENER].nm_ind[j]])))
1652 j++;
1653 if (j == nr)
1654 gmx_fatal(FARGS,"%s in %s is not an energy group\n",
1655 names[2*i],option);
1656 k = 0;
1657 while ((k < nr) &&
1658 strcasecmp(names[2*i+1],*(gnames[groups->grps[egcENER].nm_ind[k]])))
1659 k++;
1660 if (k==nr)
1661 gmx_fatal(FARGS,"%s in %s is not an energy group\n",
1662 names[2*i+1],option);
1663 if ((j < nr) && (k < nr)) {
1664 ir->opts.egp_flags[nr*j+k] |= flag;
1665 ir->opts.egp_flags[nr*k+j] |= flag;
1666 bSet = TRUE;
1670 return bSet;
1673 void do_index(const char* mdparin, const char *ndx,
1674 gmx_mtop_t *mtop,
1675 bool bVerbose,
1676 t_inputrec *ir,rvec *v,
1677 warninp_t wi)
1679 t_blocka *grps;
1680 gmx_groups_t *groups;
1681 int natoms;
1682 t_symtab *symtab;
1683 t_atoms atoms_all;
1684 char warnbuf[STRLEN],**gnames;
1685 int nr,ntcg,ntau_t,nref_t,nacc,nofg,nSA,nSA_points,nSA_time,nSA_temp;
1686 int nstcmin;
1687 int nacg,nfreeze,nfrdim,nenergy,nvcm,nuser;
1688 int nadress_refs, nadress_tf_grp_names;
1689 char *ptr1[MAXPTR],*ptr2[MAXPTR],*ptr3[MAXPTR];
1690 int i,j,k,restnm;
1691 real SAtime;
1692 bool bExcl,bTable,bSetTCpar,bAnneal,bRest;
1693 int nQMmethod,nQMbasis,nQMcharge,nQMmult,nbSH,nCASorb,nCASelec,
1694 nSAon,nSAoff,nSAsteps,nQMg,nbOPT,nbTS;
1695 char warn_buf[STRLEN];
1697 if (bVerbose)
1698 fprintf(stderr,"processing index file...\n");
1699 debug_gmx();
1700 if (ndx == NULL) {
1701 snew(grps,1);
1702 snew(grps->index,1);
1703 snew(gnames,1);
1704 atoms_all = gmx_mtop_global_atoms(mtop);
1705 analyse(&atoms_all,grps,&gnames,FALSE,TRUE);
1706 free_t_atoms(&atoms_all,FALSE);
1707 } else {
1708 grps = init_index(ndx,&gnames);
1711 groups = &mtop->groups;
1712 natoms = mtop->natoms;
1713 symtab = &mtop->symtab;
1715 snew(groups->grpname,grps->nr+1);
1717 for(i=0; (i<grps->nr); i++) {
1718 groups->grpname[i] = put_symtab(symtab,gnames[i]);
1720 groups->grpname[i] = put_symtab(symtab,"rest");
1721 restnm=i;
1722 srenew(gnames,grps->nr+1);
1723 gnames[restnm] = *(groups->grpname[i]);
1724 groups->ngrpname = grps->nr+1;
1726 set_warning_line(wi,mdparin,-1);
1728 ntau_t = str_nelem(tau_t,MAXPTR,ptr1);
1729 nref_t = str_nelem(ref_t,MAXPTR,ptr2);
1730 ntcg = str_nelem(tcgrps,MAXPTR,ptr3);
1731 if ((ntau_t != ntcg) || (nref_t != ntcg)) {
1732 gmx_fatal(FARGS,"Invalid T coupling input: %d groups, %d ref_t values and "
1733 "%d tau_t values",ntcg,nref_t,ntau_t);
1736 bSetTCpar = (ir->etc || EI_SD(ir->eI) || ir->eI==eiBD || EI_TPI(ir->eI));
1737 do_numbering(natoms,groups,ntcg,ptr3,grps,gnames,egcTC,
1738 restnm,bSetTCpar ? egrptpALL : egrptpALL_GENREST,bVerbose,wi);
1739 nr = groups->grps[egcTC].nr;
1740 ir->opts.ngtc = nr;
1741 snew(ir->opts.nrdf,nr);
1742 snew(ir->opts.tau_t,nr);
1743 snew(ir->opts.ref_t,nr);
1744 if (ir->eI==eiBD && ir->bd_fric==0) {
1745 fprintf(stderr,"bd_fric=0, so tau_t will be used as the inverse friction constant(s)\n");
1748 if (bSetTCpar)
1750 if (nr != nref_t)
1752 gmx_fatal(FARGS,"Not enough ref_t and tau_t values!");
1754 nstcmin = (ir->etc == etcBERENDSEN ? nstcmin1 : nstcmin2);
1756 for(i=0; (i<nr); i++)
1758 ir->opts.tau_t[i] = strtod(ptr1[i],NULL);
1759 if (ir->opts.tau_t[i] < 0)
1761 gmx_fatal(FARGS,"tau_t for group %d negative",i);
1763 /* We check the relative magnitude of the coupling time tau_t.
1764 * V-rescale works correctly, even for tau_t=0.
1766 if ((ir->etc == etcBERENDSEN || ir->etc == etcNOSEHOOVER) &&
1767 ir->opts.tau_t[i] != 0 &&
1768 ir->opts.tau_t[i] < nstcmin*ir->nstcalcenergy*ir->delta_t)
1770 sprintf(warn_buf,"For proper integration of the %s thermostat, tau_t (%g) should be at least %d times larger than nstcalcenergy*dt (%g)",
1771 ETCOUPLTYPE(ir->etc),
1772 ir->opts.tau_t[i],nstcmin,
1773 ir->nstcalcenergy*ir->delta_t);
1774 warning(wi,warn_buf);
1777 for(i=0; (i<nr); i++)
1779 ir->opts.ref_t[i] = strtod(ptr2[i],NULL);
1780 if (ir->opts.ref_t[i] < 0)
1782 gmx_fatal(FARGS,"ref_t for group %d negative",i);
1787 /* Simulated annealing for each group. There are nr groups */
1788 nSA = str_nelem(anneal,MAXPTR,ptr1);
1789 if (nSA == 1 && (ptr1[0][0]=='n' || ptr1[0][0]=='N'))
1790 nSA = 0;
1791 if(nSA>0 && nSA != nr)
1792 gmx_fatal(FARGS,"Not enough annealing values: %d (for %d groups)\n",nSA,nr);
1793 else {
1794 snew(ir->opts.annealing,nr);
1795 snew(ir->opts.anneal_npoints,nr);
1796 snew(ir->opts.anneal_time,nr);
1797 snew(ir->opts.anneal_temp,nr);
1798 for(i=0;i<nr;i++) {
1799 ir->opts.annealing[i]=eannNO;
1800 ir->opts.anneal_npoints[i]=0;
1801 ir->opts.anneal_time[i]=NULL;
1802 ir->opts.anneal_temp[i]=NULL;
1804 if (nSA > 0) {
1805 bAnneal=FALSE;
1806 for(i=0;i<nr;i++) {
1807 if(ptr1[i][0]=='n' || ptr1[i][0]=='N') {
1808 ir->opts.annealing[i]=eannNO;
1809 } else if(ptr1[i][0]=='s'|| ptr1[i][0]=='S') {
1810 ir->opts.annealing[i]=eannSINGLE;
1811 bAnneal=TRUE;
1812 } else if(ptr1[i][0]=='p'|| ptr1[i][0]=='P') {
1813 ir->opts.annealing[i]=eannPERIODIC;
1814 bAnneal=TRUE;
1817 if(bAnneal) {
1818 /* Read the other fields too */
1819 nSA_points = str_nelem(anneal_npoints,MAXPTR,ptr1);
1820 if(nSA_points!=nSA)
1821 gmx_fatal(FARGS,"Found %d annealing_npoints values for %d groups\n",nSA_points,nSA);
1822 for(k=0,i=0;i<nr;i++) {
1823 ir->opts.anneal_npoints[i]=strtol(ptr1[i],NULL,10);
1824 if(ir->opts.anneal_npoints[i]==1)
1825 gmx_fatal(FARGS,"Please specify at least a start and an end point for annealing\n");
1826 snew(ir->opts.anneal_time[i],ir->opts.anneal_npoints[i]);
1827 snew(ir->opts.anneal_temp[i],ir->opts.anneal_npoints[i]);
1828 k += ir->opts.anneal_npoints[i];
1831 nSA_time = str_nelem(anneal_time,MAXPTR,ptr1);
1832 if(nSA_time!=k)
1833 gmx_fatal(FARGS,"Found %d annealing_time values, wanter %d\n",nSA_time,k);
1834 nSA_temp = str_nelem(anneal_temp,MAXPTR,ptr2);
1835 if(nSA_temp!=k)
1836 gmx_fatal(FARGS,"Found %d annealing_temp values, wanted %d\n",nSA_temp,k);
1838 for(i=0,k=0;i<nr;i++) {
1840 for(j=0;j<ir->opts.anneal_npoints[i];j++) {
1841 ir->opts.anneal_time[i][j]=strtod(ptr1[k],NULL);
1842 ir->opts.anneal_temp[i][j]=strtod(ptr2[k],NULL);
1843 if(j==0) {
1844 if(ir->opts.anneal_time[i][0] > (ir->init_t+GMX_REAL_EPS))
1845 gmx_fatal(FARGS,"First time point for annealing > init_t.\n");
1846 } else {
1847 /* j>0 */
1848 if(ir->opts.anneal_time[i][j]<ir->opts.anneal_time[i][j-1])
1849 gmx_fatal(FARGS,"Annealing timepoints out of order: t=%f comes after t=%f\n",
1850 ir->opts.anneal_time[i][j],ir->opts.anneal_time[i][j-1]);
1852 if(ir->opts.anneal_temp[i][j]<0)
1853 gmx_fatal(FARGS,"Found negative temperature in annealing: %f\n",ir->opts.anneal_temp[i][j]);
1854 k++;
1857 /* Print out some summary information, to make sure we got it right */
1858 for(i=0,k=0;i<nr;i++) {
1859 if(ir->opts.annealing[i]!=eannNO) {
1860 j = groups->grps[egcTC].nm_ind[i];
1861 fprintf(stderr,"Simulated annealing for group %s: %s, %d timepoints\n",
1862 *(groups->grpname[j]),eann_names[ir->opts.annealing[i]],
1863 ir->opts.anneal_npoints[i]);
1864 fprintf(stderr,"Time (ps) Temperature (K)\n");
1865 /* All terms except the last one */
1866 for(j=0;j<(ir->opts.anneal_npoints[i]-1);j++)
1867 fprintf(stderr,"%9.1f %5.1f\n",ir->opts.anneal_time[i][j],ir->opts.anneal_temp[i][j]);
1869 /* Finally the last one */
1870 j = ir->opts.anneal_npoints[i]-1;
1871 if(ir->opts.annealing[i]==eannSINGLE)
1872 fprintf(stderr,"%9.1f- %5.1f\n",ir->opts.anneal_time[i][j],ir->opts.anneal_temp[i][j]);
1873 else {
1874 fprintf(stderr,"%9.1f %5.1f\n",ir->opts.anneal_time[i][j],ir->opts.anneal_temp[i][j]);
1875 if(fabs(ir->opts.anneal_temp[i][j]-ir->opts.anneal_temp[i][0])>GMX_REAL_EPS)
1876 warning_note(wi,"There is a temperature jump when your annealing loops back.\n");
1884 if (ir->ePull != epullNO) {
1885 make_pull_groups(ir->pull,pull_grp,grps,gnames);
1888 nacc = str_nelem(acc,MAXPTR,ptr1);
1889 nacg = str_nelem(accgrps,MAXPTR,ptr2);
1890 if (nacg*DIM != nacc)
1891 gmx_fatal(FARGS,"Invalid Acceleration input: %d groups and %d acc. values",
1892 nacg,nacc);
1893 do_numbering(natoms,groups,nacg,ptr2,grps,gnames,egcACC,
1894 restnm,egrptpALL_GENREST,bVerbose,wi);
1895 nr = groups->grps[egcACC].nr;
1896 snew(ir->opts.acc,nr);
1897 ir->opts.ngacc=nr;
1899 for(i=k=0; (i<nacg); i++)
1900 for(j=0; (j<DIM); j++,k++)
1901 ir->opts.acc[i][j]=strtod(ptr1[k],NULL);
1902 for( ;(i<nr); i++)
1903 for(j=0; (j<DIM); j++)
1904 ir->opts.acc[i][j]=0;
1906 nfrdim = str_nelem(frdim,MAXPTR,ptr1);
1907 nfreeze = str_nelem(freeze,MAXPTR,ptr2);
1908 if (nfrdim != DIM*nfreeze)
1909 gmx_fatal(FARGS,"Invalid Freezing input: %d groups and %d freeze values",
1910 nfreeze,nfrdim);
1911 do_numbering(natoms,groups,nfreeze,ptr2,grps,gnames,egcFREEZE,
1912 restnm,egrptpALL_GENREST,bVerbose,wi);
1913 nr = groups->grps[egcFREEZE].nr;
1914 ir->opts.ngfrz=nr;
1915 snew(ir->opts.nFreeze,nr);
1916 for(i=k=0; (i<nfreeze); i++)
1917 for(j=0; (j<DIM); j++,k++) {
1918 ir->opts.nFreeze[i][j]=(strncasecmp(ptr1[k],"Y",1)==0);
1919 if (!ir->opts.nFreeze[i][j]) {
1920 if (strncasecmp(ptr1[k],"N",1) != 0) {
1921 sprintf(warnbuf,"Please use Y(ES) or N(O) for freezedim only "
1922 "(not %s)", ptr1[k]);
1923 warning(wi,warn_buf);
1927 for( ; (i<nr); i++)
1928 for(j=0; (j<DIM); j++)
1929 ir->opts.nFreeze[i][j]=0;
1931 nenergy=str_nelem(energy,MAXPTR,ptr1);
1932 do_numbering(natoms,groups,nenergy,ptr1,grps,gnames,egcENER,
1933 restnm,egrptpALL_GENREST,bVerbose,wi);
1934 add_wall_energrps(groups,ir->nwall,symtab);
1935 ir->opts.ngener = groups->grps[egcENER].nr;
1936 nvcm=str_nelem(vcm,MAXPTR,ptr1);
1937 bRest =
1938 do_numbering(natoms,groups,nvcm,ptr1,grps,gnames,egcVCM,
1939 restnm,nvcm==0 ? egrptpALL_GENREST : egrptpPART,bVerbose,wi);
1940 if (bRest) {
1941 warning(wi,"Some atoms are not part of any center of mass motion removal group.\n"
1942 "This may lead to artifacts.\n"
1943 "In most cases one should use one group for the whole system.");
1946 /* Now we have filled the freeze struct, so we can calculate NRDF */
1947 calc_nrdf(mtop,ir,gnames);
1949 if (v && NULL) {
1950 real fac,ntot=0;
1952 /* Must check per group! */
1953 for(i=0; (i<ir->opts.ngtc); i++)
1954 ntot += ir->opts.nrdf[i];
1955 if (ntot != (DIM*natoms)) {
1956 fac = sqrt(ntot/(DIM*natoms));
1957 if (bVerbose)
1958 fprintf(stderr,"Scaling velocities by a factor of %.3f to account for constraints\n"
1959 "and removal of center of mass motion\n",fac);
1960 for(i=0; (i<natoms); i++)
1961 svmul(fac,v[i],v[i]);
1965 nuser=str_nelem(user1,MAXPTR,ptr1);
1966 do_numbering(natoms,groups,nuser,ptr1,grps,gnames,egcUser1,
1967 restnm,egrptpALL_GENREST,bVerbose,wi);
1968 nuser=str_nelem(user2,MAXPTR,ptr1);
1969 do_numbering(natoms,groups,nuser,ptr1,grps,gnames,egcUser2,
1970 restnm,egrptpALL_GENREST,bVerbose,wi);
1971 nuser=str_nelem(xtc_grps,MAXPTR,ptr1);
1972 do_numbering(natoms,groups,nuser,ptr1,grps,gnames,egcXTC,
1973 restnm,egrptpONE,bVerbose,wi);
1974 nofg = str_nelem(orirefitgrp,MAXPTR,ptr1);
1975 do_numbering(natoms,groups,nofg,ptr1,grps,gnames,egcORFIT,
1976 restnm,egrptpALL_GENREST,bVerbose,wi);
1978 /* QMMM input processing */
1979 nQMg = str_nelem(QMMM,MAXPTR,ptr1);
1980 nQMmethod = str_nelem(QMmethod,MAXPTR,ptr2);
1981 nQMbasis = str_nelem(QMbasis,MAXPTR,ptr3);
1982 if((nQMmethod != nQMg)||(nQMbasis != nQMg)){
1983 gmx_fatal(FARGS,"Invalid QMMM input: %d groups %d basissets"
1984 " and %d methods\n",nQMg,nQMbasis,nQMmethod);
1986 /* group rest, if any, is always MM! */
1987 do_numbering(natoms,groups,nQMg,ptr1,grps,gnames,egcQMMM,
1988 restnm,egrptpALL_GENREST,bVerbose,wi);
1989 nr = nQMg; /*atoms->grps[egcQMMM].nr;*/
1990 ir->opts.ngQM = nQMg;
1991 snew(ir->opts.QMmethod,nr);
1992 snew(ir->opts.QMbasis,nr);
1993 for(i=0;i<nr;i++){
1994 /* input consists of strings: RHF CASSCF PM3 .. These need to be
1995 * converted to the corresponding enum in names.c
1997 ir->opts.QMmethod[i] = search_QMstring(ptr2[i],eQMmethodNR,
1998 eQMmethod_names);
1999 ir->opts.QMbasis[i] = search_QMstring(ptr3[i],eQMbasisNR,
2000 eQMbasis_names);
2003 nQMmult = str_nelem(QMmult,MAXPTR,ptr1);
2004 nQMcharge = str_nelem(QMcharge,MAXPTR,ptr2);
2005 nbSH = str_nelem(bSH,MAXPTR,ptr3);
2006 snew(ir->opts.QMmult,nr);
2007 snew(ir->opts.QMcharge,nr);
2008 snew(ir->opts.bSH,nr);
2010 for(i=0;i<nr;i++){
2011 ir->opts.QMmult[i] = strtol(ptr1[i],NULL,10);
2012 ir->opts.QMcharge[i] = strtol(ptr2[i],NULL,10);
2013 ir->opts.bSH[i] = (strncasecmp(ptr3[i],"Y",1)==0);
2016 nCASelec = str_nelem(CASelectrons,MAXPTR,ptr1);
2017 nCASorb = str_nelem(CASorbitals,MAXPTR,ptr2);
2018 snew(ir->opts.CASelectrons,nr);
2019 snew(ir->opts.CASorbitals,nr);
2020 for(i=0;i<nr;i++){
2021 ir->opts.CASelectrons[i]= strtol(ptr1[i],NULL,10);
2022 ir->opts.CASorbitals[i] = strtol(ptr2[i],NULL,10);
2024 /* special optimization options */
2026 nbOPT = str_nelem(bOPT,MAXPTR,ptr1);
2027 nbTS = str_nelem(bTS,MAXPTR,ptr2);
2028 snew(ir->opts.bOPT,nr);
2029 snew(ir->opts.bTS,nr);
2030 for(i=0;i<nr;i++){
2031 ir->opts.bOPT[i] = (strncasecmp(ptr1[i],"Y",1)==0);
2032 ir->opts.bTS[i] = (strncasecmp(ptr2[i],"Y",1)==0);
2034 nSAon = str_nelem(SAon,MAXPTR,ptr1);
2035 nSAoff = str_nelem(SAoff,MAXPTR,ptr2);
2036 nSAsteps = str_nelem(SAsteps,MAXPTR,ptr3);
2037 snew(ir->opts.SAon,nr);
2038 snew(ir->opts.SAoff,nr);
2039 snew(ir->opts.SAsteps,nr);
2041 for(i=0;i<nr;i++){
2042 ir->opts.SAon[i] = strtod(ptr1[i],NULL);
2043 ir->opts.SAoff[i] = strtod(ptr2[i],NULL);
2044 ir->opts.SAsteps[i] = strtol(ptr3[i],NULL,10);
2046 /* end of QMMM input */
2048 /* AdResS reference input */
2049 nadress_refs = str_nelem(adress_refs,MAXPTR,ptr1);
2051 for(i=0; (i<nadress_refs); i++)
2052 ir->adress_refs[i]=strtod(ptr1[i],NULL);
2053 for( ;(i<DIM); i++)
2054 ir->adress_refs[i]=0;
2056 /* End AdResS input */
2058 if (bVerbose)
2059 for(i=0; (i<egcNR); i++) {
2060 fprintf(stderr,"%-16s has %d element(s):",gtypes[i],groups->grps[i].nr);
2061 for(j=0; (j<groups->grps[i].nr); j++)
2062 fprintf(stderr," %s",*(groups->grpname[groups->grps[i].nm_ind[j]]));
2063 fprintf(stderr,"\n");
2066 nr = groups->grps[egcENER].nr;
2067 snew(ir->opts.egp_flags,nr*nr);
2069 bExcl = do_egp_flag(ir,groups,"energygrp_excl",egpexcl,EGP_EXCL);
2070 if (bExcl && EEL_FULL(ir->coulombtype))
2071 warning(wi,"Can not exclude the lattice Coulomb energy between energy groups");
2073 bTable = do_egp_flag(ir,groups,"energygrp_table",egptable,EGP_TABLE);
2074 if (bTable && !(ir->vdwtype == evdwUSER) &&
2075 !(ir->coulombtype == eelUSER) &&!(ir->coulombtype == eelPMEUSER))
2076 gmx_fatal(FARGS,"Can only have energy group pair tables in combination with user tables for VdW and/or Coulomb");
2078 decode_cos(efield_x,&(ir->ex[XX]),FALSE);
2079 decode_cos(efield_xt,&(ir->et[XX]),TRUE);
2080 decode_cos(efield_y,&(ir->ex[YY]),FALSE);
2081 decode_cos(efield_yt,&(ir->et[YY]),TRUE);
2082 decode_cos(efield_z,&(ir->ex[ZZ]),FALSE);
2083 decode_cos(efield_zt,&(ir->et[ZZ]),TRUE);
2085 /* AdResS multiple tf tables input */
2086 nadress_tf_grp_names = str_nelem(adress_tf_grp_names,MAXPTR,ptr1);
2087 ir->n_adress_tf_grps = nadress_tf_grp_names;
2088 snew(ir->adress_tf_table_index, nadress_tf_grp_names);
2090 nr = groups->grps[egcENER].nr;
2092 if (nadress_tf_grp_names > 0){
2093 for (i=0; i <nadress_tf_grp_names; i++){
2094 //search for the group name mathching the tf group name
2095 k = 0;
2096 while ((k < nr) &&
2097 strcasecmp(ptr1[i],(char*)(gnames[groups->grps[egcENER].nm_ind[k]])))
2098 k++;
2099 if (k==nr) gmx_fatal(FARGS,"Adress tf energy group %s not found\n",ptr1[i]);
2101 ir->adress_tf_table_index[i] = k;
2102 if (debug) fprintf(debug,"found tf group %s id %d \n",ptr1[i], k);
2106 /* end AdResS multiple tf tables input */
2108 for(i=0; (i<grps->nr); i++)
2109 sfree(gnames[i]);
2110 sfree(gnames);
2111 done_blocka(grps);
2112 sfree(grps);
2118 static void check_disre(gmx_mtop_t *mtop)
2120 gmx_ffparams_t *ffparams;
2121 t_functype *functype;
2122 t_iparams *ip;
2123 int i,ndouble,ftype;
2124 int label,old_label;
2126 if (gmx_mtop_ftype_count(mtop,F_DISRES) > 0) {
2127 ffparams = &mtop->ffparams;
2128 functype = ffparams->functype;
2129 ip = ffparams->iparams;
2130 ndouble = 0;
2131 old_label = -1;
2132 for(i=0; i<ffparams->ntypes; i++) {
2133 ftype = functype[i];
2134 if (ftype == F_DISRES) {
2135 label = ip[i].disres.label;
2136 if (label == old_label) {
2137 fprintf(stderr,"Distance restraint index %d occurs twice\n",label);
2138 ndouble++;
2140 old_label = label;
2143 if (ndouble>0)
2144 gmx_fatal(FARGS,"Found %d double distance restraint indices,\n"
2145 "probably the parameters for multiple pairs in one restraint "
2146 "are not identical\n",ndouble);
2150 static bool absolute_reference(t_inputrec *ir,gmx_mtop_t *sys,ivec AbsRef)
2152 int d,g,i;
2153 gmx_mtop_ilistloop_t iloop;
2154 t_ilist *ilist;
2155 int nmol;
2156 t_iparams *pr;
2158 /* Check the COM */
2159 for(d=0; d<DIM; d++) {
2160 AbsRef[d] = (d < ndof_com(ir) ? 0 : 1);
2162 /* Check for freeze groups */
2163 for(g=0; g<ir->opts.ngfrz; g++) {
2164 for(d=0; d<DIM; d++) {
2165 if (ir->opts.nFreeze[g][d] != 0) {
2166 AbsRef[d] = 1;
2170 /* Check for position restraints */
2171 iloop = gmx_mtop_ilistloop_init(sys);
2172 while (gmx_mtop_ilistloop_next(iloop,&ilist,&nmol)) {
2173 if (nmol > 0) {
2174 for(i=0; i<ilist[F_POSRES].nr; i+=2) {
2175 pr = &sys->ffparams.iparams[ilist[F_POSRES].iatoms[i]];
2176 for(d=0; d<DIM; d++) {
2177 if (pr->posres.fcA[d] != 0) {
2178 AbsRef[d] = 1;
2185 return (AbsRef[XX] != 0 && AbsRef[YY] != 0 && AbsRef[ZZ] != 0);
2188 void triple_check(const char *mdparin,t_inputrec *ir,gmx_mtop_t *sys,
2189 warninp_t wi)
2191 char err_buf[256];
2192 int i,m,g,nmol,npct;
2193 bool bCharge,bAcc;
2194 real gdt_max,*mgrp,mt;
2195 rvec acc;
2196 gmx_mtop_atomloop_block_t aloopb;
2197 gmx_mtop_atomloop_all_t aloop;
2198 t_atom *atom;
2199 ivec AbsRef;
2200 char warn_buf[STRLEN];
2202 set_warning_line(wi,mdparin,-1);
2204 if (EI_DYNAMICS(ir->eI) && !EI_SD(ir->eI) && ir->eI != eiBD &&
2205 ir->comm_mode == ecmNO &&
2206 !(absolute_reference(ir,sys,AbsRef) || ir->nsteps <= 10)) {
2207 warning(wi,"You are not using center of mass motion removal (mdp option comm-mode), numerical rounding errors can lead to build up of kinetic energy of the center of mass");
2210 bCharge = FALSE;
2211 aloopb = gmx_mtop_atomloop_block_init(sys);
2212 while (gmx_mtop_atomloop_block_next(aloopb,&atom,&nmol)) {
2213 if (atom->q != 0 || atom->qB != 0) {
2214 bCharge = TRUE;
2218 if (!bCharge) {
2219 if (EEL_FULL(ir->coulombtype)) {
2220 sprintf(err_buf,
2221 "You are using full electrostatics treatment %s for a system without charges.\n"
2222 "This costs a lot of performance for just processing zeros, consider using %s instead.\n",
2223 EELTYPE(ir->coulombtype),EELTYPE(eelCUT));
2224 warning(wi,err_buf);
2226 } else {
2227 if (ir->coulombtype == eelCUT && ir->rcoulomb > 0) {
2228 sprintf(err_buf,
2229 "You are using a plain Coulomb cut-off, which might produce artifacts.\n"
2230 "You might want to consider using %s electrostatics.\n",
2231 EELTYPE(eelPME));
2232 warning_note(wi,err_buf);
2236 /* Generalized reaction field */
2237 if (ir->opts.ngtc == 0) {
2238 sprintf(err_buf,"No temperature coupling while using coulombtype %s",
2239 eel_names[eelGRF]);
2240 CHECK(ir->coulombtype == eelGRF);
2242 else {
2243 sprintf(err_buf,"When using coulombtype = %s"
2244 " ref_t for temperature coupling should be > 0",
2245 eel_names[eelGRF]);
2246 CHECK((ir->coulombtype == eelGRF) && (ir->opts.ref_t[0] <= 0));
2249 if (ir->eI == eiSD1) {
2250 gdt_max = 0;
2251 for(i=0; (i<ir->opts.ngtc); i++)
2252 gdt_max = max(gdt_max,ir->delta_t/ir->opts.tau_t[i]);
2253 if (0.5*gdt_max > 0.0015) {
2254 sprintf(warn_buf,"The relative error with integrator %s is 0.5*delta_t/tau_t = %g, you might want to switch to integrator %s\n",
2255 ei_names[ir->eI],0.5*gdt_max,ei_names[eiSD2]);
2256 warning_note(wi,warn_buf);
2260 bAcc = FALSE;
2261 for(i=0; (i<sys->groups.grps[egcACC].nr); i++) {
2262 for(m=0; (m<DIM); m++) {
2263 if (fabs(ir->opts.acc[i][m]) > 1e-6) {
2264 bAcc = TRUE;
2268 if (bAcc) {
2269 clear_rvec(acc);
2270 snew(mgrp,sys->groups.grps[egcACC].nr);
2271 aloop = gmx_mtop_atomloop_all_init(sys);
2272 while (gmx_mtop_atomloop_all_next(aloop,&i,&atom)) {
2273 mgrp[ggrpnr(&sys->groups,egcACC,i)] += atom->m;
2275 mt = 0.0;
2276 for(i=0; (i<sys->groups.grps[egcACC].nr); i++) {
2277 for(m=0; (m<DIM); m++)
2278 acc[m] += ir->opts.acc[i][m]*mgrp[i];
2279 mt += mgrp[i];
2281 for(m=0; (m<DIM); m++) {
2282 if (fabs(acc[m]) > 1e-6) {
2283 const char *dim[DIM] = { "X", "Y", "Z" };
2284 fprintf(stderr,
2285 "Net Acceleration in %s direction, will %s be corrected\n",
2286 dim[m],ir->nstcomm != 0 ? "" : "not");
2287 if (ir->nstcomm != 0 && m < ndof_com(ir)) {
2288 acc[m] /= mt;
2289 for (i=0; (i<sys->groups.grps[egcACC].nr); i++)
2290 ir->opts.acc[i][m] -= acc[m];
2294 sfree(mgrp);
2297 if (ir->efep != efepNO && ir->sc_alpha != 0 &&
2298 !gmx_within_tol(sys->ffparams.reppow,12.0,10*GMX_DOUBLE_EPS)) {
2299 gmx_fatal(FARGS,"Soft-core interactions are only supported with VdW repulsion power 12");
2302 if (ir->ePull != epullNO) {
2303 if (ir->pull->grp[0].nat == 0) {
2304 absolute_reference(ir,sys,AbsRef);
2305 for(m=0; m<DIM; m++) {
2306 if (ir->pull->dim[m] && !AbsRef[m]) {
2307 warning(wi,"You are using an absolute reference for pulling, but the rest of the system does not have an absolute reference. This will lead to artifacts.");
2308 break;
2313 if (ir->pull->eGeom == epullgDIRPBC) {
2314 for(i=0; i<3; i++) {
2315 for(m=0; m<=i; m++) {
2316 if ((ir->epc != epcNO && ir->compress[i][m] != 0) ||
2317 ir->deform[i][m] != 0) {
2318 for(g=1; g<ir->pull->ngrp; g++) {
2319 if (ir->pull->grp[g].vec[m] != 0) {
2320 gmx_fatal(FARGS,"Can not have dynamic box while using pull geometry '%s' (dim %c)",EPULLGEOM(ir->pull->eGeom),'x'+m);
2329 check_disre(sys);
2332 void double_check(t_inputrec *ir,matrix box,bool bConstr,warninp_t wi)
2334 real min_size;
2335 bool bTWIN;
2336 char warn_buf[STRLEN];
2337 const char *ptr;
2339 ptr = check_box(ir->ePBC,box);
2340 if (ptr) {
2341 warning_error(wi,ptr);
2344 if (bConstr && ir->eConstrAlg == econtSHAKE) {
2345 if (ir->shake_tol <= 0.0) {
2346 sprintf(warn_buf,"ERROR: shake_tol must be > 0 instead of %g\n",
2347 ir->shake_tol);
2348 warning_error(wi,warn_buf);
2351 if (IR_TWINRANGE(*ir) && ir->nstlist > 1) {
2352 sprintf(warn_buf,"With twin-range cut-off's and SHAKE the virial and the pressure are incorrect.");
2353 if (ir->epc == epcNO) {
2354 warning(wi,warn_buf);
2355 } else {
2356 warning_error(wi,warn_buf);
2361 if( (ir->eConstrAlg == econtLINCS) && bConstr) {
2362 /* If we have Lincs constraints: */
2363 if(ir->eI==eiMD && ir->etc==etcNO &&
2364 ir->eConstrAlg==econtLINCS && ir->nLincsIter==1) {
2365 sprintf(warn_buf,"For energy conservation with LINCS, lincs_iter should be 2 or larger.\n");
2366 warning_note(wi,warn_buf);
2369 if ((ir->eI == eiCG || ir->eI == eiLBFGS) && (ir->nProjOrder<8)) {
2370 sprintf(warn_buf,"For accurate %s with LINCS constraints, lincs_order should be 8 or more.",ei_names[ir->eI]);
2371 warning_note(wi,warn_buf);
2373 if (ir->epc==epcMTTK) {
2374 warning_error(wi,"MTTK not compatible with lincs -- use shake instead.");
2378 if (ir->LincsWarnAngle > 90.0) {
2379 sprintf(warn_buf,"lincs-warnangle can not be larger than 90 degrees, setting it to 90.\n");
2380 warning(wi,warn_buf);
2381 ir->LincsWarnAngle = 90.0;
2384 if (ir->ePBC != epbcNONE) {
2385 if (ir->nstlist == 0) {
2386 warning(wi,"With nstlist=0 atoms are only put into the box at step 0, therefore drifting atoms might cause the simulation to crash.");
2388 bTWIN = (ir->rlistlong > ir->rlist);
2389 if (ir->ns_type == ensGRID) {
2390 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC,box)) {
2391 sprintf(warn_buf,"ERROR: The cut-off length is longer than half the shortest box vector or longer than the smallest box diagonal element. Increase the box size or decrease %s.\n",
2392 bTWIN ? (ir->rcoulomb==ir->rlistlong ? "rcoulomb" : "rvdw"):"rlist");
2393 warning_error(wi,warn_buf);
2395 } else {
2396 min_size = min(box[XX][XX],min(box[YY][YY],box[ZZ][ZZ]));
2397 if (2*ir->rlistlong >= min_size) {
2398 sprintf(warn_buf,"ERROR: One of the box lengths is smaller than twice the cut-off length. Increase the box size or decrease rlist.");
2399 warning_error(wi,warn_buf);
2400 if (TRICLINIC(box))
2401 fprintf(stderr,"Grid search might allow larger cut-off's than simple search with triclinic boxes.");
2407 void check_chargegroup_radii(const gmx_mtop_t *mtop,const t_inputrec *ir,
2408 rvec *x,
2409 warninp_t wi)
2411 real rvdw1,rvdw2,rcoul1,rcoul2;
2412 char warn_buf[STRLEN];
2414 calc_chargegroup_radii(mtop,x,&rvdw1,&rvdw2,&rcoul1,&rcoul2);
2416 if (rvdw1 > 0)
2418 printf("Largest charge group radii for Van der Waals: %5.3f, %5.3f nm\n",
2419 rvdw1,rvdw2);
2421 if (rcoul1 > 0)
2423 printf("Largest charge group radii for Coulomb: %5.3f, %5.3f nm\n",
2424 rcoul1,rcoul2);
2427 if (ir->rlist > 0)
2429 if (rvdw1 + rvdw2 > ir->rlist ||
2430 rcoul1 + rcoul2 > ir->rlist)
2432 sprintf(warn_buf,"The sum of the two largest charge group radii (%f) is larger than rlist (%f)\n",max(rvdw1+rvdw2,rcoul1+rcoul2),ir->rlist);
2433 warning(wi,warn_buf);
2435 else
2437 /* Here we do not use the zero at cut-off macro,
2438 * since user defined interactions might purposely
2439 * not be zero at the cut-off.
2441 if (EVDW_IS_ZERO_AT_CUTOFF(ir->vdwtype) &&
2442 rvdw1 + rvdw2 > ir->rlistlong - ir->rvdw)
2444 sprintf(warn_buf,"The sum of the two largest charge group radii (%f) is larger than %s (%f) - rvdw (%f)\n",
2445 rvdw1+rvdw2,
2446 ir->rlistlong > ir->rlist ? "rlistlong" : "rlist",
2447 ir->rlist,ir->rvdw);
2448 if (ir_NVE(ir))
2450 warning(wi,warn_buf);
2452 else
2454 warning_note(wi,warn_buf);
2457 if (EEL_IS_ZERO_AT_CUTOFF(ir->coulombtype) &&
2458 rcoul1 + rcoul2 > ir->rlistlong - ir->rcoulomb)
2460 sprintf(warn_buf,"The sum of the two largest charge group radii (%f) is larger than %s (%f) - rcoulomb (%f)\n",
2461 rcoul1+rcoul2,
2462 ir->rlistlong > ir->rlist ? "rlistlong" : "rlist",
2463 ir->rlistlong,ir->rcoulomb);
2464 if (ir_NVE(ir))
2466 warning(wi,warn_buf);
2468 else
2470 warning_note(wi,warn_buf);