Split up compare.*
[gromacs.git] / src / gromacs / mdtypes / inputrec.cpp
blobaa4b2502526beb5d775af8a9866c530e96af9db1
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2010, The GROMACS development team.
6 * Copyright (c) 2012,2014,2015,2016, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "gmxpre.h"
39 #include "inputrec.h"
41 #include <cstdio>
42 #include <cstring>
44 #include <algorithm>
46 #include "gromacs/math/veccompare.h"
47 #include "gromacs/math/vecdump.h"
48 #include "gromacs/mdtypes/md_enums.h"
49 #include "gromacs/mdtypes/pull-params.h"
50 #include "gromacs/pbcutil/pbc.h"
51 #include "gromacs/utility/compare.h"
52 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/smalloc.h"
55 #include "gromacs/utility/snprintf.h"
56 #include "gromacs/utility/stringutil.h"
57 #include "gromacs/utility/txtdump.h"
59 //! Macro to select a bool name
60 #define EBOOL(e) gmx::boolToString(e)
62 /* The minimum number of integration steps required for reasonably accurate
63 * integration of first and second order coupling algorithms.
65 const int nstmin_berendsen_tcouple = 5;
66 const int nstmin_berendsen_pcouple = 10;
67 const int nstmin_harmonic = 20;
69 static int nst_wanted(const t_inputrec *ir)
71 if (ir->nstlist > 0)
73 return ir->nstlist;
75 else
77 return 10;
81 int ir_optimal_nstcalcenergy(const t_inputrec *ir)
83 return nst_wanted(ir);
86 int tcouple_min_integration_steps(int etc)
88 int n;
90 switch (etc)
92 case etcNO:
93 n = 0;
94 break;
95 case etcBERENDSEN:
96 case etcYES:
97 n = nstmin_berendsen_tcouple;
98 break;
99 case etcVRESCALE:
100 /* V-rescale supports instantaneous rescaling */
101 n = 0;
102 break;
103 case etcNOSEHOOVER:
104 n = nstmin_harmonic;
105 break;
106 case etcANDERSEN:
107 case etcANDERSENMASSIVE:
108 n = 1;
109 break;
110 default:
111 gmx_incons("Unknown etc value");
112 n = 0;
115 return n;
118 int ir_optimal_nsttcouple(const t_inputrec *ir)
120 int nmin, nwanted, n;
121 real tau_min;
122 int g;
124 nmin = tcouple_min_integration_steps(ir->etc);
126 nwanted = nst_wanted(ir);
128 tau_min = 1e20;
129 if (ir->etc != etcNO)
131 for (g = 0; g < ir->opts.ngtc; g++)
133 if (ir->opts.tau_t[g] > 0)
135 tau_min = std::min(tau_min, ir->opts.tau_t[g]);
140 if (nmin == 0 || ir->delta_t*nwanted <= tau_min)
142 n = nwanted;
144 else
146 n = (int)(tau_min/(ir->delta_t*nmin) + 0.001);
147 if (n < 1)
149 n = 1;
151 while (nwanted % n != 0)
153 n--;
157 return n;
160 int pcouple_min_integration_steps(int epc)
162 int n;
164 switch (epc)
166 case epcNO:
167 n = 0;
168 break;
169 case etcBERENDSEN:
170 case epcISOTROPIC:
171 n = nstmin_berendsen_pcouple;
172 break;
173 case epcPARRINELLORAHMAN:
174 case epcMTTK:
175 n = nstmin_harmonic;
176 break;
177 default:
178 gmx_incons("Unknown epc value");
179 n = 0;
182 return n;
185 int ir_optimal_nstpcouple(const t_inputrec *ir)
187 int nmin, nwanted, n;
189 nmin = pcouple_min_integration_steps(ir->epc);
191 nwanted = nst_wanted(ir);
193 if (nmin == 0 || ir->delta_t*nwanted <= ir->tau_p)
195 n = nwanted;
197 else
199 n = static_cast<int>(ir->tau_p/(ir->delta_t*nmin) + 0.001);
200 if (n < 1)
202 n = 1;
204 while (nwanted % n != 0)
206 n--;
210 return n;
213 gmx_bool ir_coulomb_switched(const t_inputrec *ir)
215 return (ir->coulombtype == eelSWITCH ||
216 ir->coulombtype == eelSHIFT ||
217 ir->coulombtype == eelENCADSHIFT ||
218 ir->coulombtype == eelPMESWITCH ||
219 ir->coulombtype == eelPMEUSERSWITCH ||
220 ir->coulomb_modifier == eintmodPOTSWITCH ||
221 ir->coulomb_modifier == eintmodFORCESWITCH);
224 gmx_bool ir_coulomb_is_zero_at_cutoff(const t_inputrec *ir)
226 return (ir->cutoff_scheme == ecutsVERLET ||
227 ir_coulomb_switched(ir) || ir->coulomb_modifier != eintmodNONE ||
228 ir->coulombtype == eelRF_ZERO);
231 gmx_bool ir_coulomb_might_be_zero_at_cutoff(const t_inputrec *ir)
233 return (ir_coulomb_is_zero_at_cutoff(ir) || ir->coulombtype == eelUSER || ir->coulombtype == eelPMEUSER);
236 gmx_bool ir_vdw_switched(const t_inputrec *ir)
238 return (ir->vdwtype == evdwSWITCH ||
239 ir->vdwtype == evdwSHIFT ||
240 ir->vdwtype == evdwENCADSHIFT ||
241 ir->vdw_modifier == eintmodPOTSWITCH ||
242 ir->vdw_modifier == eintmodFORCESWITCH);
245 gmx_bool ir_vdw_is_zero_at_cutoff(const t_inputrec *ir)
247 return (ir->cutoff_scheme == ecutsVERLET ||
248 ir_vdw_switched(ir) || ir->vdw_modifier != eintmodNONE);
251 gmx_bool ir_vdw_might_be_zero_at_cutoff(const t_inputrec *ir)
253 return (ir_vdw_is_zero_at_cutoff(ir) || ir->vdwtype == evdwUSER);
256 static void done_pull_group(t_pull_group *pgrp)
258 if (pgrp->nat > 0)
260 sfree(pgrp->ind);
261 sfree(pgrp->weight);
265 static void done_pull_params(pull_params_t *pull)
267 int i;
269 for (i = 0; i < pull->ngroup+1; i++)
271 done_pull_group(pull->group);
274 sfree(pull->group);
275 sfree(pull->coord);
278 void done_inputrec(t_inputrec *ir)
280 int m;
282 for (m = 0; (m < DIM); m++)
284 if (ir->ex[m].a)
286 sfree(ir->ex[m].a);
288 if (ir->ex[m].phi)
290 sfree(ir->ex[m].phi);
292 if (ir->et[m].a)
294 sfree(ir->et[m].a);
296 if (ir->et[m].phi)
298 sfree(ir->et[m].phi);
302 sfree(ir->opts.nrdf);
303 sfree(ir->opts.ref_t);
304 sfree(ir->opts.annealing);
305 sfree(ir->opts.anneal_npoints);
306 sfree(ir->opts.anneal_time);
307 sfree(ir->opts.anneal_temp);
308 sfree(ir->opts.tau_t);
309 sfree(ir->opts.acc);
310 sfree(ir->opts.nFreeze);
311 sfree(ir->opts.QMmethod);
312 sfree(ir->opts.QMbasis);
313 sfree(ir->opts.QMcharge);
314 sfree(ir->opts.QMmult);
315 sfree(ir->opts.bSH);
316 sfree(ir->opts.CASorbitals);
317 sfree(ir->opts.CASelectrons);
318 sfree(ir->opts.SAon);
319 sfree(ir->opts.SAoff);
320 sfree(ir->opts.SAsteps);
321 sfree(ir->opts.bOPT);
322 sfree(ir->opts.bTS);
324 if (ir->pull)
326 done_pull_params(ir->pull);
327 sfree(ir->pull);
331 static void pr_qm_opts(FILE *fp, int indent, const char *title, const t_grpopts *opts)
333 fprintf(fp, "%s:\n", title);
335 pr_int(fp, indent, "ngQM", opts->ngQM);
336 if (opts->ngQM > 0)
338 pr_ivec(fp, indent, "QMmethod", opts->QMmethod, opts->ngQM, FALSE);
339 pr_ivec(fp, indent, "QMbasis", opts->QMbasis, opts->ngQM, FALSE);
340 pr_ivec(fp, indent, "QMcharge", opts->QMcharge, opts->ngQM, FALSE);
341 pr_ivec(fp, indent, "QMmult", opts->QMmult, opts->ngQM, FALSE);
342 pr_bvec(fp, indent, "SH", opts->bSH, opts->ngQM, FALSE);
343 pr_ivec(fp, indent, "CASorbitals", opts->CASorbitals, opts->ngQM, FALSE);
344 pr_ivec(fp, indent, "CASelectrons", opts->CASelectrons, opts->ngQM, FALSE);
345 pr_rvec(fp, indent, "SAon", opts->SAon, opts->ngQM, FALSE);
346 pr_rvec(fp, indent, "SAoff", opts->SAoff, opts->ngQM, FALSE);
347 pr_ivec(fp, indent, "SAsteps", opts->SAsteps, opts->ngQM, FALSE);
348 pr_bvec(fp, indent, "bOPT", opts->bOPT, opts->ngQM, FALSE);
349 pr_bvec(fp, indent, "bTS", opts->bTS, opts->ngQM, FALSE);
353 static void pr_grp_opts(FILE *out, int indent, const char *title, const t_grpopts *opts,
354 gmx_bool bMDPformat)
356 int i, m, j;
358 if (!bMDPformat)
360 fprintf(out, "%s:\n", title);
363 pr_indent(out, indent);
364 fprintf(out, "nrdf%s", bMDPformat ? " = " : ":");
365 for (i = 0; (i < opts->ngtc); i++)
367 fprintf(out, " %10g", opts->nrdf[i]);
369 fprintf(out, "\n");
371 pr_indent(out, indent);
372 fprintf(out, "ref-t%s", bMDPformat ? " = " : ":");
373 for (i = 0; (i < opts->ngtc); i++)
375 fprintf(out, " %10g", opts->ref_t[i]);
377 fprintf(out, "\n");
379 pr_indent(out, indent);
380 fprintf(out, "tau-t%s", bMDPformat ? " = " : ":");
381 for (i = 0; (i < opts->ngtc); i++)
383 fprintf(out, " %10g", opts->tau_t[i]);
385 fprintf(out, "\n");
387 /* Pretty-print the simulated annealing info */
388 fprintf(out, "annealing%s", bMDPformat ? " = " : ":");
389 for (i = 0; (i < opts->ngtc); i++)
391 fprintf(out, " %10s", EANNEAL(opts->annealing[i]));
393 fprintf(out, "\n");
395 fprintf(out, "annealing-npoints%s", bMDPformat ? " = " : ":");
396 for (i = 0; (i < opts->ngtc); i++)
398 fprintf(out, " %10d", opts->anneal_npoints[i]);
400 fprintf(out, "\n");
402 for (i = 0; (i < opts->ngtc); i++)
404 if (opts->anneal_npoints[i] > 0)
406 fprintf(out, "annealing-time [%d]:\t", i);
407 for (j = 0; (j < opts->anneal_npoints[i]); j++)
409 fprintf(out, " %10.1f", opts->anneal_time[i][j]);
411 fprintf(out, "\n");
412 fprintf(out, "annealing-temp [%d]:\t", i);
413 for (j = 0; (j < opts->anneal_npoints[i]); j++)
415 fprintf(out, " %10.1f", opts->anneal_temp[i][j]);
417 fprintf(out, "\n");
421 pr_indent(out, indent);
422 fprintf(out, "acc:\t");
423 for (i = 0; (i < opts->ngacc); i++)
425 for (m = 0; (m < DIM); m++)
427 fprintf(out, " %10g", opts->acc[i][m]);
430 fprintf(out, "\n");
432 pr_indent(out, indent);
433 fprintf(out, "nfreeze:");
434 for (i = 0; (i < opts->ngfrz); i++)
436 for (m = 0; (m < DIM); m++)
438 fprintf(out, " %10s", opts->nFreeze[i][m] ? "Y" : "N");
441 fprintf(out, "\n");
444 for (i = 0; (i < opts->ngener); i++)
446 pr_indent(out, indent);
447 fprintf(out, "energygrp-flags[%3d]:", i);
448 for (m = 0; (m < opts->ngener); m++)
450 fprintf(out, " %d", opts->egp_flags[opts->ngener*i+m]);
452 fprintf(out, "\n");
455 fflush(out);
458 static void pr_matrix(FILE *fp, int indent, const char *title, const rvec *m,
459 gmx_bool bMDPformat)
461 if (bMDPformat)
463 fprintf(fp, "%-10s = %g %g %g %g %g %g\n", title,
464 m[XX][XX], m[YY][YY], m[ZZ][ZZ], m[XX][YY], m[XX][ZZ], m[YY][ZZ]);
466 else
468 pr_rvecs(fp, indent, title, m, DIM);
472 static void pr_cosine(FILE *fp, int indent, const char *title, const t_cosines *cos,
473 gmx_bool bMDPformat)
475 int j;
477 if (bMDPformat)
479 fprintf(fp, "%s = %d\n", title, cos->n);
481 else
483 indent = pr_title(fp, indent, title);
484 pr_indent(fp, indent);
485 fprintf(fp, "n = %d\n", cos->n);
486 if (cos->n > 0)
488 pr_indent(fp, indent+2);
489 fprintf(fp, "a =");
490 for (j = 0; (j < cos->n); j++)
492 fprintf(fp, " %e", cos->a[j]);
494 fprintf(fp, "\n");
495 pr_indent(fp, indent+2);
496 fprintf(fp, "phi =");
497 for (j = 0; (j < cos->n); j++)
499 fprintf(fp, " %e", cos->phi[j]);
501 fprintf(fp, "\n");
506 #define PS(t, s) pr_str(fp, indent, t, s)
507 #define PI(t, s) pr_int(fp, indent, t, s)
508 #define PSTEP(t, s) pr_int64(fp, indent, t, s)
509 #define PR(t, s) pr_real(fp, indent, t, s)
510 #define PD(t, s) pr_double(fp, indent, t, s)
512 static void pr_pull_group(FILE *fp, int indent, int g, const t_pull_group *pgrp)
514 pr_indent(fp, indent);
515 fprintf(fp, "pull-group %d:\n", g);
516 indent += 2;
517 pr_ivec_block(fp, indent, "atom", pgrp->ind, pgrp->nat, TRUE);
518 pr_rvec(fp, indent, "weight", pgrp->weight, pgrp->nweight, TRUE);
519 PI("pbcatom", pgrp->pbcatom);
522 static void pr_pull_coord(FILE *fp, int indent, int c, const t_pull_coord *pcrd)
524 int g;
526 pr_indent(fp, indent);
527 fprintf(fp, "pull-coord %d:\n", c);
528 PS("type", EPULLTYPE(pcrd->eType));
529 if (pcrd->eType == epullEXTERNAL)
531 PS("potential-provider", pcrd->externalPotentialProvider);
533 PS("geometry", EPULLGEOM(pcrd->eGeom));
534 for (g = 0; g < pcrd->ngroup; g++)
536 char buf[10];
538 sprintf(buf, "group[%d]", g);
539 PI(buf, pcrd->group[g]);
541 pr_ivec(fp, indent, "dim", pcrd->dim, DIM, TRUE);
542 pr_rvec(fp, indent, "origin", pcrd->origin, DIM, TRUE);
543 pr_rvec(fp, indent, "vec", pcrd->vec, DIM, TRUE);
544 PS("start", EBOOL(pcrd->bStart));
545 PR("init", pcrd->init);
546 PR("rate", pcrd->rate);
547 PR("k", pcrd->k);
548 PR("kB", pcrd->kB);
551 static void pr_simtempvals(FILE *fp, int indent, const t_simtemp *simtemp, int n_lambda)
553 PS("simulated-tempering-scaling", ESIMTEMP(simtemp->eSimTempScale));
554 PR("sim-temp-low", simtemp->simtemp_low);
555 PR("sim-temp-high", simtemp->simtemp_high);
556 pr_rvec(fp, indent, "simulated tempering temperatures", simtemp->temperatures, n_lambda, TRUE);
559 static void pr_expandedvals(FILE *fp, int indent, const t_expanded *expand, int n_lambda)
562 PI("nstexpanded", expand->nstexpanded);
563 PS("lmc-stats", elamstats_names[expand->elamstats]);
564 PS("lmc-move", elmcmove_names[expand->elmcmove]);
565 PS("lmc-weights-equil", elmceq_names[expand->elmceq]);
566 if (expand->elmceq == elmceqNUMATLAM)
568 PI("weight-equil-number-all-lambda", expand->equil_n_at_lam);
570 if (expand->elmceq == elmceqSAMPLES)
572 PI("weight-equil-number-samples", expand->equil_samples);
574 if (expand->elmceq == elmceqSTEPS)
576 PI("weight-equil-number-steps", expand->equil_steps);
578 if (expand->elmceq == elmceqWLDELTA)
580 PR("weight-equil-wl-delta", expand->equil_wl_delta);
582 if (expand->elmceq == elmceqRATIO)
584 PR("weight-equil-count-ratio", expand->equil_ratio);
586 PI("lmc-seed", expand->lmc_seed);
587 PR("mc-temperature", expand->mc_temp);
588 PI("lmc-repeats", expand->lmc_repeats);
589 PI("lmc-gibbsdelta", expand->gibbsdeltalam);
590 PI("lmc-forced-nstart", expand->lmc_forced_nstart);
591 PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
592 PI("nst-transition-matrix", expand->nstTij);
593 PI("mininum-var-min", expand->minvarmin); /*default is reasonable */
594 PI("weight-c-range", expand->c_range); /* default is just C=0 */
595 PR("wl-scale", expand->wl_scale);
596 PR("wl-ratio", expand->wl_ratio);
597 PR("init-wl-delta", expand->init_wl_delta);
598 PS("wl-oneovert", EBOOL(expand->bWLoneovert));
600 pr_indent(fp, indent);
601 pr_rvec(fp, indent, "init-lambda-weights", expand->init_lambda_weights, n_lambda, TRUE);
602 PS("init-weights", EBOOL(expand->bInit_weights));
605 static void pr_fepvals(FILE *fp, int indent, const t_lambda *fep, gmx_bool bMDPformat)
607 int i, j;
609 PR("init-lambda", fep->init_lambda);
610 PI("init-lambda-state", fep->init_fep_state);
611 PR("delta-lambda", fep->delta_lambda);
612 PI("nstdhdl", fep->nstdhdl);
614 if (!bMDPformat)
616 PI("n-lambdas", fep->n_lambda);
618 if (fep->n_lambda > 0)
620 pr_indent(fp, indent);
621 fprintf(fp, "separate-dvdl%s\n", bMDPformat ? " = " : ":");
622 for (i = 0; i < efptNR; i++)
624 fprintf(fp, "%18s = ", efpt_names[i]);
625 if (fep->separate_dvdl[i])
627 fprintf(fp, " TRUE");
629 else
631 fprintf(fp, " FALSE");
633 fprintf(fp, "\n");
635 fprintf(fp, "all-lambdas%s\n", bMDPformat ? " = " : ":");
636 for (i = 0; i < efptNR; i++)
638 fprintf(fp, "%18s = ", efpt_names[i]);
639 for (j = 0; j < fep->n_lambda; j++)
641 fprintf(fp, " %10g", fep->all_lambda[i][j]);
643 fprintf(fp, "\n");
646 PI("calc-lambda-neighbors", fep->lambda_neighbors);
647 PS("dhdl-print-energy", edHdLPrintEnergy_names[fep->edHdLPrintEnergy]);
648 PR("sc-alpha", fep->sc_alpha);
649 PI("sc-power", fep->sc_power);
650 PR("sc-r-power", fep->sc_r_power);
651 PR("sc-sigma", fep->sc_sigma);
652 PR("sc-sigma-min", fep->sc_sigma_min);
653 PS("sc-coul", EBOOL(fep->bScCoul));
654 PI("dh-hist-size", fep->dh_hist_size);
655 PD("dh-hist-spacing", fep->dh_hist_spacing);
656 PS("separate-dhdl-file", SEPDHDLFILETYPE(fep->separate_dhdl_file));
657 PS("dhdl-derivatives", DHDLDERIVATIVESTYPE(fep->dhdl_derivatives));
660 static void pr_pull(FILE *fp, int indent, const pull_params_t *pull)
662 int g;
664 PR("pull-cylinder-r", pull->cylinder_r);
665 PR("pull-constr-tol", pull->constr_tol);
666 PS("pull-print-COM", EBOOL(pull->bPrintCOM));
667 PS("pull-print-ref-value", EBOOL(pull->bPrintRefValue));
668 PS("pull-print-components", EBOOL(pull->bPrintComp));
669 PI("pull-nstxout", pull->nstxout);
670 PI("pull-nstfout", pull->nstfout);
671 PI("pull-ngroups", pull->ngroup);
672 for (g = 0; g < pull->ngroup; g++)
674 pr_pull_group(fp, indent, g, &pull->group[g]);
676 PI("pull-ncoords", pull->ncoord);
677 for (g = 0; g < pull->ncoord; g++)
679 pr_pull_coord(fp, indent, g, &pull->coord[g]);
683 static void pr_rotgrp(FILE *fp, int indent, int g, const t_rotgrp *rotg)
685 pr_indent(fp, indent);
686 fprintf(fp, "rot-group %d:\n", g);
687 indent += 2;
688 PS("rot-type", EROTGEOM(rotg->eType));
689 PS("rot-massw", EBOOL(rotg->bMassW));
690 pr_ivec_block(fp, indent, "atom", rotg->ind, rotg->nat, TRUE);
691 pr_rvecs(fp, indent, "x-ref", rotg->x_ref, rotg->nat);
692 pr_rvec(fp, indent, "rot-vec", rotg->vec, DIM, TRUE);
693 pr_rvec(fp, indent, "rot-pivot", rotg->pivot, DIM, TRUE);
694 PR("rot-rate", rotg->rate);
695 PR("rot-k", rotg->k);
696 PR("rot-slab-dist", rotg->slab_dist);
697 PR("rot-min-gauss", rotg->min_gaussian);
698 PR("rot-eps", rotg->eps);
699 PS("rot-fit-method", EROTFIT(rotg->eFittype));
700 PI("rot-potfit-nstep", rotg->PotAngle_nstep);
701 PR("rot-potfit-step", rotg->PotAngle_step);
704 static void pr_rot(FILE *fp, int indent, const t_rot *rot)
706 int g;
708 PI("rot-nstrout", rot->nstrout);
709 PI("rot-nstsout", rot->nstsout);
710 PI("rot-ngroups", rot->ngrp);
711 for (g = 0; g < rot->ngrp; g++)
713 pr_rotgrp(fp, indent, g, &rot->grp[g]);
718 static void pr_swap(FILE *fp, int indent, const t_swapcoords *swap)
720 char str[STRLEN];
722 /* Enums for better readability of the code */
723 enum {
724 eCompA = 0, eCompB
728 PI("swap-frequency", swap->nstswap);
730 /* The split groups that define the compartments */
731 for (int j = 0; j < 2; j++)
733 snprintf(str, STRLEN, "massw_split%d", j);
734 PS(str, EBOOL(swap->massw_split[j]));
735 snprintf(str, STRLEN, "split atoms group %d", j);
736 pr_ivec_block(fp, indent, str, swap->grp[j].ind, swap->grp[j].nat, TRUE);
739 /* The solvent group */
740 snprintf(str, STRLEN, "solvent group %s", swap->grp[eGrpSolvent].molname);
741 pr_ivec_block(fp, indent, str, swap->grp[eGrpSolvent].ind, swap->grp[eGrpSolvent].nat, TRUE);
743 /* Now print the indices for all the ion groups: */
744 for (int ig = eSwapFixedGrpNR; ig < swap->ngrp; ig++)
746 snprintf(str, STRLEN, "ion group %s", swap->grp[ig].molname);
747 pr_ivec_block(fp, indent, str, swap->grp[ig].ind, swap->grp[ig].nat, TRUE);
750 PR("cyl0-r", swap->cyl0r);
751 PR("cyl0-up", swap->cyl0u);
752 PR("cyl0-down", swap->cyl0l);
753 PR("cyl1-r", swap->cyl1r);
754 PR("cyl1-up", swap->cyl1u);
755 PR("cyl1-down", swap->cyl1l);
756 PI("coupl-steps", swap->nAverage);
758 /* Print the requested ion counts for both compartments */
759 for (int ic = eCompA; ic <= eCompB; ic++)
761 for (int ig = eSwapFixedGrpNR; ig < swap->ngrp; ig++)
763 snprintf(str, STRLEN, "%s-in-%c", swap->grp[ig].molname, 'A'+ic);
764 PI(str, swap->grp[ig].nmolReq[ic]);
768 PR("threshold", swap->threshold);
769 PR("bulk-offsetA", swap->bulkOffset[eCompA]);
770 PR("bulk-offsetB", swap->bulkOffset[eCompB]);
774 static void pr_imd(FILE *fp, int indent, const t_IMD *imd)
776 PI("IMD-atoms", imd->nat);
777 pr_ivec_block(fp, indent, "atom", imd->ind, imd->nat, TRUE);
781 void pr_inputrec(FILE *fp, int indent, const char *title, const t_inputrec *ir,
782 gmx_bool bMDPformat)
784 const char *infbuf = "inf";
786 if (available(fp, ir, indent, title))
788 if (!bMDPformat)
790 indent = pr_title(fp, indent, title);
792 /* Try to make this list appear in the same order as the
793 * options are written in the default mdout.mdp, and with
794 * the same user-exposed names to facilitate debugging.
796 PS("integrator", EI(ir->eI));
797 PR("tinit", ir->init_t);
798 PR("dt", ir->delta_t);
799 PSTEP("nsteps", ir->nsteps);
800 PSTEP("init-step", ir->init_step);
801 PI("simulation-part", ir->simulation_part);
802 PS("comm-mode", ECOM(ir->comm_mode));
803 PI("nstcomm", ir->nstcomm);
805 /* Langevin dynamics */
806 PR("bd-fric", ir->bd_fric);
807 PSTEP("ld-seed", ir->ld_seed);
809 /* Energy minimization */
810 PR("emtol", ir->em_tol);
811 PR("emstep", ir->em_stepsize);
812 PI("niter", ir->niter);
813 PR("fcstep", ir->fc_stepsize);
814 PI("nstcgsteep", ir->nstcgsteep);
815 PI("nbfgscorr", ir->nbfgscorr);
817 /* Test particle insertion */
818 PR("rtpi", ir->rtpi);
820 /* Output control */
821 PI("nstxout", ir->nstxout);
822 PI("nstvout", ir->nstvout);
823 PI("nstfout", ir->nstfout);
824 PI("nstlog", ir->nstlog);
825 PI("nstcalcenergy", ir->nstcalcenergy);
826 PI("nstenergy", ir->nstenergy);
827 PI("nstxout-compressed", ir->nstxout_compressed);
828 PR("compressed-x-precision", ir->x_compression_precision);
830 /* Neighborsearching parameters */
831 PS("cutoff-scheme", ECUTSCHEME(ir->cutoff_scheme));
832 PI("nstlist", ir->nstlist);
833 PS("ns-type", ENS(ir->ns_type));
834 PS("pbc", epbc_names[ir->ePBC]);
835 PS("periodic-molecules", EBOOL(ir->bPeriodicMols));
836 PR("verlet-buffer-tolerance", ir->verletbuf_tol);
837 PR("rlist", ir->rlist);
839 /* Options for electrostatics and VdW */
840 PS("coulombtype", EELTYPE(ir->coulombtype));
841 PS("coulomb-modifier", INTMODIFIER(ir->coulomb_modifier));
842 PR("rcoulomb-switch", ir->rcoulomb_switch);
843 PR("rcoulomb", ir->rcoulomb);
844 if (ir->epsilon_r != 0)
846 PR("epsilon-r", ir->epsilon_r);
848 else
850 PS("epsilon-r", infbuf);
852 if (ir->epsilon_rf != 0)
854 PR("epsilon-rf", ir->epsilon_rf);
856 else
858 PS("epsilon-rf", infbuf);
860 PS("vdw-type", EVDWTYPE(ir->vdwtype));
861 PS("vdw-modifier", INTMODIFIER(ir->vdw_modifier));
862 PR("rvdw-switch", ir->rvdw_switch);
863 PR("rvdw", ir->rvdw);
864 PS("DispCorr", EDISPCORR(ir->eDispCorr));
865 PR("table-extension", ir->tabext);
867 PR("fourierspacing", ir->fourier_spacing);
868 PI("fourier-nx", ir->nkx);
869 PI("fourier-ny", ir->nky);
870 PI("fourier-nz", ir->nkz);
871 PI("pme-order", ir->pme_order);
872 PR("ewald-rtol", ir->ewald_rtol);
873 PR("ewald-rtol-lj", ir->ewald_rtol_lj);
874 PS("lj-pme-comb-rule", ELJPMECOMBNAMES(ir->ljpme_combination_rule));
875 PR("ewald-geometry", ir->ewald_geometry);
876 PR("epsilon-surface", ir->epsilon_surface);
878 /* Implicit solvent */
879 PS("implicit-solvent", EIMPLICITSOL(ir->implicit_solvent));
881 /* Generalized born electrostatics */
882 PS("gb-algorithm", EGBALGORITHM(ir->gb_algorithm));
883 PI("nstgbradii", ir->nstgbradii);
884 PR("rgbradii", ir->rgbradii);
885 PR("gb-epsilon-solvent", ir->gb_epsilon_solvent);
886 PR("gb-saltconc", ir->gb_saltconc);
887 PR("gb-obc-alpha", ir->gb_obc_alpha);
888 PR("gb-obc-beta", ir->gb_obc_beta);
889 PR("gb-obc-gamma", ir->gb_obc_gamma);
890 PR("gb-dielectric-offset", ir->gb_dielectric_offset);
891 PS("sa-algorithm", ESAALGORITHM(ir->sa_algorithm));
892 PR("sa-surface-tension", ir->sa_surface_tension);
894 /* Options for weak coupling algorithms */
895 PS("tcoupl", ETCOUPLTYPE(ir->etc));
896 PI("nsttcouple", ir->nsttcouple);
897 PI("nh-chain-length", ir->opts.nhchainlength);
898 PS("print-nose-hoover-chain-variables", EBOOL(ir->bPrintNHChains));
900 PS("pcoupl", EPCOUPLTYPE(ir->epc));
901 PS("pcoupltype", EPCOUPLTYPETYPE(ir->epct));
902 PI("nstpcouple", ir->nstpcouple);
903 PR("tau-p", ir->tau_p);
904 pr_matrix(fp, indent, "compressibility", ir->compress, bMDPformat);
905 pr_matrix(fp, indent, "ref-p", ir->ref_p, bMDPformat);
906 PS("refcoord-scaling", EREFSCALINGTYPE(ir->refcoord_scaling));
908 if (bMDPformat)
910 fprintf(fp, "posres-com = %g %g %g\n", ir->posres_com[XX],
911 ir->posres_com[YY], ir->posres_com[ZZ]);
912 fprintf(fp, "posres-comB = %g %g %g\n", ir->posres_comB[XX],
913 ir->posres_comB[YY], ir->posres_comB[ZZ]);
915 else
917 pr_rvec(fp, indent, "posres-com", ir->posres_com, DIM, TRUE);
918 pr_rvec(fp, indent, "posres-comB", ir->posres_comB, DIM, TRUE);
921 /* QMMM */
922 PS("QMMM", EBOOL(ir->bQMMM));
923 PI("QMconstraints", ir->QMconstraints);
924 PI("QMMMscheme", ir->QMMMscheme);
925 PR("MMChargeScaleFactor", ir->scalefactor);
926 pr_qm_opts(fp, indent, "qm-opts", &(ir->opts));
928 /* CONSTRAINT OPTIONS */
929 PS("constraint-algorithm", ECONSTRTYPE(ir->eConstrAlg));
930 PS("continuation", EBOOL(ir->bContinuation));
932 PS("Shake-SOR", EBOOL(ir->bShakeSOR));
933 PR("shake-tol", ir->shake_tol);
934 PI("lincs-order", ir->nProjOrder);
935 PI("lincs-iter", ir->nLincsIter);
936 PR("lincs-warnangle", ir->LincsWarnAngle);
938 /* Walls */
939 PI("nwall", ir->nwall);
940 PS("wall-type", EWALLTYPE(ir->wall_type));
941 PR("wall-r-linpot", ir->wall_r_linpot);
942 /* wall-atomtype */
943 PI("wall-atomtype[0]", ir->wall_atomtype[0]);
944 PI("wall-atomtype[1]", ir->wall_atomtype[1]);
945 /* wall-density */
946 PR("wall-density[0]", ir->wall_density[0]);
947 PR("wall-density[1]", ir->wall_density[1]);
948 PR("wall-ewald-zfac", ir->wall_ewald_zfac);
950 /* COM PULLING */
951 PS("pull", EBOOL(ir->bPull));
952 if (ir->bPull)
954 pr_pull(fp, indent, ir->pull);
957 /* ENFORCED ROTATION */
958 PS("rotation", EBOOL(ir->bRot));
959 if (ir->bRot)
961 pr_rot(fp, indent, ir->rot);
964 /* INTERACTIVE MD */
965 PS("interactiveMD", EBOOL(ir->bIMD));
966 if (ir->bIMD)
968 pr_imd(fp, indent, ir->imd);
971 /* NMR refinement stuff */
972 PS("disre", EDISRETYPE(ir->eDisre));
973 PS("disre-weighting", EDISREWEIGHTING(ir->eDisreWeighting));
974 PS("disre-mixed", EBOOL(ir->bDisreMixed));
975 PR("dr-fc", ir->dr_fc);
976 PR("dr-tau", ir->dr_tau);
977 PR("nstdisreout", ir->nstdisreout);
979 PR("orire-fc", ir->orires_fc);
980 PR("orire-tau", ir->orires_tau);
981 PR("nstorireout", ir->nstorireout);
983 /* FREE ENERGY VARIABLES */
984 PS("free-energy", EFEPTYPE(ir->efep));
985 if (ir->efep != efepNO || ir->bSimTemp)
987 pr_fepvals(fp, indent, ir->fepvals, bMDPformat);
989 if (ir->bExpanded)
991 pr_expandedvals(fp, indent, ir->expandedvals, ir->fepvals->n_lambda);
994 /* NON-equilibrium MD stuff */
995 PR("cos-acceleration", ir->cos_accel);
996 pr_matrix(fp, indent, "deform", ir->deform, bMDPformat);
998 /* SIMULATED TEMPERING */
999 PS("simulated-tempering", EBOOL(ir->bSimTemp));
1000 if (ir->bSimTemp)
1002 pr_simtempvals(fp, indent, ir->simtempvals, ir->fepvals->n_lambda);
1005 /* ELECTRIC FIELDS */
1006 pr_cosine(fp, indent, "E-x", &(ir->ex[XX]), bMDPformat);
1007 pr_cosine(fp, indent, "E-xt", &(ir->et[XX]), bMDPformat);
1008 pr_cosine(fp, indent, "E-y", &(ir->ex[YY]), bMDPformat);
1009 pr_cosine(fp, indent, "E-yt", &(ir->et[YY]), bMDPformat);
1010 pr_cosine(fp, indent, "E-z", &(ir->ex[ZZ]), bMDPformat);
1011 pr_cosine(fp, indent, "E-zt", &(ir->et[ZZ]), bMDPformat);
1013 /* ION/WATER SWAPPING FOR COMPUTATIONAL ELECTROPHYSIOLOGY */
1014 PS("swapcoords", ESWAPTYPE(ir->eSwapCoords));
1015 if (ir->eSwapCoords != eswapNO)
1017 pr_swap(fp, indent, ir->swap);
1020 /* USER-DEFINED THINGIES */
1021 PI("userint1", ir->userint1);
1022 PI("userint2", ir->userint2);
1023 PI("userint3", ir->userint3);
1024 PI("userint4", ir->userint4);
1025 PR("userreal1", ir->userreal1);
1026 PR("userreal2", ir->userreal2);
1027 PR("userreal3", ir->userreal3);
1028 PR("userreal4", ir->userreal4);
1030 pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
1033 #undef PS
1034 #undef PR
1035 #undef PI
1037 static void cmp_grpopts(FILE *fp, const t_grpopts *opt1, const t_grpopts *opt2, real ftol, real abstol)
1039 int i, j;
1040 char buf1[256], buf2[256];
1042 cmp_int(fp, "inputrec->grpopts.ngtc", -1, opt1->ngtc, opt2->ngtc);
1043 cmp_int(fp, "inputrec->grpopts.ngacc", -1, opt1->ngacc, opt2->ngacc);
1044 cmp_int(fp, "inputrec->grpopts.ngfrz", -1, opt1->ngfrz, opt2->ngfrz);
1045 cmp_int(fp, "inputrec->grpopts.ngener", -1, opt1->ngener, opt2->ngener);
1046 for (i = 0; (i < std::min(opt1->ngtc, opt2->ngtc)); i++)
1048 cmp_real(fp, "inputrec->grpopts.nrdf", i, opt1->nrdf[i], opt2->nrdf[i], ftol, abstol);
1049 cmp_real(fp, "inputrec->grpopts.ref_t", i, opt1->ref_t[i], opt2->ref_t[i], ftol, abstol);
1050 cmp_real(fp, "inputrec->grpopts.tau_t", i, opt1->tau_t[i], opt2->tau_t[i], ftol, abstol);
1051 cmp_int(fp, "inputrec->grpopts.annealing", i, opt1->annealing[i], opt2->annealing[i]);
1052 cmp_int(fp, "inputrec->grpopts.anneal_npoints", i,
1053 opt1->anneal_npoints[i], opt2->anneal_npoints[i]);
1054 if (opt1->anneal_npoints[i] == opt2->anneal_npoints[i])
1056 sprintf(buf1, "inputrec->grpopts.anneal_time[%d]", i);
1057 sprintf(buf2, "inputrec->grpopts.anneal_temp[%d]", i);
1058 for (j = 0; j < opt1->anneal_npoints[i]; j++)
1060 cmp_real(fp, buf1, j, opt1->anneal_time[i][j], opt2->anneal_time[i][j], ftol, abstol);
1061 cmp_real(fp, buf2, j, opt1->anneal_temp[i][j], opt2->anneal_temp[i][j], ftol, abstol);
1065 if (opt1->ngener == opt2->ngener)
1067 for (i = 0; i < opt1->ngener; i++)
1069 for (j = i; j < opt1->ngener; j++)
1071 sprintf(buf1, "inputrec->grpopts.egp_flags[%d]", i);
1072 cmp_int(fp, buf1, j,
1073 opt1->egp_flags[opt1->ngener*i+j],
1074 opt2->egp_flags[opt1->ngener*i+j]);
1078 for (i = 0; (i < std::min(opt1->ngacc, opt2->ngacc)); i++)
1080 cmp_rvec(fp, "inputrec->grpopts.acc", i, opt1->acc[i], opt2->acc[i], ftol, abstol);
1082 for (i = 0; (i < std::min(opt1->ngfrz, opt2->ngfrz)); i++)
1084 cmp_ivec(fp, "inputrec->grpopts.nFreeze", i, opt1->nFreeze[i], opt2->nFreeze[i]);
1089 static void cmp_cosines(FILE *fp, const char *s, const t_cosines c1[DIM], const t_cosines c2[DIM], real ftol, real abstol)
1091 int i, m;
1092 char buf[256];
1094 for (m = 0; (m < DIM); m++)
1096 sprintf(buf, "inputrec->%s[%d]", s, m);
1097 cmp_int(fp, buf, 0, c1->n, c2->n);
1098 for (i = 0; (i < std::min(c1->n, c2->n)); i++)
1100 cmp_real(fp, buf, i, c1->a[i], c2->a[i], ftol, abstol);
1101 cmp_real(fp, buf, i, c1->phi[i], c2->phi[i], ftol, abstol);
1105 static void cmp_pull(FILE *fp)
1107 fprintf(fp, "WARNING: Both files use COM pulling, but comparing of the pull struct is not implemented (yet). The pull parameters could be the same or different.\n");
1110 static void cmp_simtempvals(FILE *fp, const t_simtemp *simtemp1, const t_simtemp *simtemp2, int n_lambda, real ftol, real abstol)
1112 int i;
1113 cmp_int(fp, "inputrec->simtempvals->eSimTempScale", -1, simtemp1->eSimTempScale, simtemp2->eSimTempScale);
1114 cmp_real(fp, "inputrec->simtempvals->simtemp_high", -1, simtemp1->simtemp_high, simtemp2->simtemp_high, ftol, abstol);
1115 cmp_real(fp, "inputrec->simtempvals->simtemp_low", -1, simtemp1->simtemp_low, simtemp2->simtemp_low, ftol, abstol);
1116 for (i = 0; i < n_lambda; i++)
1118 cmp_real(fp, "inputrec->simtempvals->temperatures", -1, simtemp1->temperatures[i], simtemp2->temperatures[i], ftol, abstol);
1122 static void cmp_expandedvals(FILE *fp, const t_expanded *expand1, const t_expanded *expand2, int n_lambda, real ftol, real abstol)
1124 int i;
1126 cmp_bool(fp, "inputrec->fepvals->bInit_weights", -1, expand1->bInit_weights, expand2->bInit_weights);
1127 cmp_bool(fp, "inputrec->fepvals->bWLoneovert", -1, expand1->bWLoneovert, expand2->bWLoneovert);
1129 for (i = 0; i < n_lambda; i++)
1131 cmp_real(fp, "inputrec->expandedvals->init_lambda_weights", -1,
1132 expand1->init_lambda_weights[i], expand2->init_lambda_weights[i], ftol, abstol);
1135 cmp_int(fp, "inputrec->expandedvals->lambda-stats", -1, expand1->elamstats, expand2->elamstats);
1136 cmp_int(fp, "inputrec->expandedvals->lambda-mc-move", -1, expand1->elmcmove, expand2->elmcmove);
1137 cmp_int(fp, "inputrec->expandedvals->lmc-repeats", -1, expand1->lmc_repeats, expand2->lmc_repeats);
1138 cmp_int(fp, "inputrec->expandedvals->lmc-gibbsdelta", -1, expand1->gibbsdeltalam, expand2->gibbsdeltalam);
1139 cmp_int(fp, "inputrec->expandedvals->lmc-forced-nstart", -1, expand1->lmc_forced_nstart, expand2->lmc_forced_nstart);
1140 cmp_int(fp, "inputrec->expandedvals->lambda-weights-equil", -1, expand1->elmceq, expand2->elmceq);
1141 cmp_int(fp, "inputrec->expandedvals->,weight-equil-number-all-lambda", -1, expand1->equil_n_at_lam, expand2->equil_n_at_lam);
1142 cmp_int(fp, "inputrec->expandedvals->weight-equil-number-samples", -1, expand1->equil_samples, expand2->equil_samples);
1143 cmp_int(fp, "inputrec->expandedvals->weight-equil-number-steps", -1, expand1->equil_steps, expand2->equil_steps);
1144 cmp_real(fp, "inputrec->expandedvals->weight-equil-wl-delta", -1, expand1->equil_wl_delta, expand2->equil_wl_delta, ftol, abstol);
1145 cmp_real(fp, "inputrec->expandedvals->weight-equil-count-ratio", -1, expand1->equil_ratio, expand2->equil_ratio, ftol, abstol);
1146 cmp_bool(fp, "inputrec->expandedvals->symmetrized-transition-matrix", -1, expand1->bSymmetrizedTMatrix, expand2->bSymmetrizedTMatrix);
1147 cmp_int(fp, "inputrec->expandedvals->nstTij", -1, expand1->nstTij, expand2->nstTij);
1148 cmp_int(fp, "inputrec->expandedvals->mininum-var-min", -1, expand1->minvarmin, expand2->minvarmin); /*default is reasonable */
1149 cmp_int(fp, "inputrec->expandedvals->weight-c-range", -1, expand1->c_range, expand2->c_range); /* default is just C=0 */
1150 cmp_real(fp, "inputrec->expandedvals->wl-scale", -1, expand1->wl_scale, expand2->wl_scale, ftol, abstol);
1151 cmp_real(fp, "inputrec->expandedvals->init-wl-delta", -1, expand1->init_wl_delta, expand2->init_wl_delta, ftol, abstol);
1152 cmp_real(fp, "inputrec->expandedvals->wl-ratio", -1, expand1->wl_ratio, expand2->wl_ratio, ftol, abstol);
1153 cmp_int(fp, "inputrec->expandedvals->nstexpanded", -1, expand1->nstexpanded, expand2->nstexpanded);
1154 cmp_int(fp, "inputrec->expandedvals->lmc-seed", -1, expand1->lmc_seed, expand2->lmc_seed);
1155 cmp_real(fp, "inputrec->expandedvals->mc-temperature", -1, expand1->mc_temp, expand2->mc_temp, ftol, abstol);
1158 static void cmp_fepvals(FILE *fp, const t_lambda *fep1, const t_lambda *fep2, real ftol, real abstol)
1160 int i, j;
1161 cmp_int(fp, "inputrec->nstdhdl", -1, fep1->nstdhdl, fep2->nstdhdl);
1162 cmp_double(fp, "inputrec->fepvals->init_fep_state", -1, fep1->init_fep_state, fep2->init_fep_state, ftol, abstol);
1163 cmp_double(fp, "inputrec->fepvals->delta_lambda", -1, fep1->delta_lambda, fep2->delta_lambda, ftol, abstol);
1164 cmp_int(fp, "inputrec->fepvals->n_lambda", -1, fep1->n_lambda, fep2->n_lambda);
1165 for (i = 0; i < efptNR; i++)
1167 for (j = 0; j < std::min(fep1->n_lambda, fep2->n_lambda); j++)
1169 cmp_double(fp, "inputrec->fepvals->all_lambda", -1, fep1->all_lambda[i][j], fep2->all_lambda[i][j], ftol, abstol);
1172 cmp_int(fp, "inputrec->fepvals->lambda_neighbors", 1, fep1->lambda_neighbors,
1173 fep2->lambda_neighbors);
1174 cmp_real(fp, "inputrec->fepvals->sc_alpha", -1, fep1->sc_alpha, fep2->sc_alpha, ftol, abstol);
1175 cmp_int(fp, "inputrec->fepvals->sc_power", -1, fep1->sc_power, fep2->sc_power);
1176 cmp_real(fp, "inputrec->fepvals->sc_r_power", -1, fep1->sc_r_power, fep2->sc_r_power, ftol, abstol);
1177 cmp_real(fp, "inputrec->fepvals->sc_sigma", -1, fep1->sc_sigma, fep2->sc_sigma, ftol, abstol);
1178 cmp_int(fp, "inputrec->fepvals->edHdLPrintEnergy", -1, fep1->edHdLPrintEnergy, fep1->edHdLPrintEnergy);
1179 cmp_bool(fp, "inputrec->fepvals->bScCoul", -1, fep1->bScCoul, fep1->bScCoul);
1180 cmp_int(fp, "inputrec->separate_dhdl_file", -1, fep1->separate_dhdl_file, fep2->separate_dhdl_file);
1181 cmp_int(fp, "inputrec->dhdl_derivatives", -1, fep1->dhdl_derivatives, fep2->dhdl_derivatives);
1182 cmp_int(fp, "inputrec->dh_hist_size", -1, fep1->dh_hist_size, fep2->dh_hist_size);
1183 cmp_double(fp, "inputrec->dh_hist_spacing", -1, fep1->dh_hist_spacing, fep2->dh_hist_spacing, ftol, abstol);
1186 void cmp_inputrec(FILE *fp, const t_inputrec *ir1, const t_inputrec *ir2, real ftol, real abstol)
1188 fprintf(fp, "comparing inputrec\n");
1190 /* gcc 2.96 doesnt like these defines at all, but issues a huge list
1191 * of warnings. Maybe it will change in future versions, but for the
1192 * moment I've spelled them out instead. /EL 000820
1193 * #define CIB(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
1194 * #define CII(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
1195 * #define CIR(s) cmp_real(fp,"inputrec->"#s,0,ir1->##s,ir2->##s,ftol)
1197 cmp_int(fp, "inputrec->eI", -1, ir1->eI, ir2->eI);
1198 cmp_int64(fp, "inputrec->nsteps", ir1->nsteps, ir2->nsteps);
1199 cmp_int64(fp, "inputrec->init_step", ir1->init_step, ir2->init_step);
1200 cmp_int(fp, "inputrec->simulation_part", -1, ir1->simulation_part, ir2->simulation_part);
1201 cmp_int(fp, "inputrec->ePBC", -1, ir1->ePBC, ir2->ePBC);
1202 cmp_int(fp, "inputrec->bPeriodicMols", -1, ir1->bPeriodicMols, ir2->bPeriodicMols);
1203 cmp_int(fp, "inputrec->cutoff_scheme", -1, ir1->cutoff_scheme, ir2->cutoff_scheme);
1204 cmp_int(fp, "inputrec->ns_type", -1, ir1->ns_type, ir2->ns_type);
1205 cmp_int(fp, "inputrec->nstlist", -1, ir1->nstlist, ir2->nstlist);
1206 cmp_int(fp, "inputrec->nstcomm", -1, ir1->nstcomm, ir2->nstcomm);
1207 cmp_int(fp, "inputrec->comm_mode", -1, ir1->comm_mode, ir2->comm_mode);
1208 cmp_int(fp, "inputrec->nstlog", -1, ir1->nstlog, ir2->nstlog);
1209 cmp_int(fp, "inputrec->nstxout", -1, ir1->nstxout, ir2->nstxout);
1210 cmp_int(fp, "inputrec->nstvout", -1, ir1->nstvout, ir2->nstvout);
1211 cmp_int(fp, "inputrec->nstfout", -1, ir1->nstfout, ir2->nstfout);
1212 cmp_int(fp, "inputrec->nstcalcenergy", -1, ir1->nstcalcenergy, ir2->nstcalcenergy);
1213 cmp_int(fp, "inputrec->nstenergy", -1, ir1->nstenergy, ir2->nstenergy);
1214 cmp_int(fp, "inputrec->nstxout_compressed", -1, ir1->nstxout_compressed, ir2->nstxout_compressed);
1215 cmp_double(fp, "inputrec->init_t", -1, ir1->init_t, ir2->init_t, ftol, abstol);
1216 cmp_double(fp, "inputrec->delta_t", -1, ir1->delta_t, ir2->delta_t, ftol, abstol);
1217 cmp_real(fp, "inputrec->x_compression_precision", -1, ir1->x_compression_precision, ir2->x_compression_precision, ftol, abstol);
1218 cmp_real(fp, "inputrec->fourierspacing", -1, ir1->fourier_spacing, ir2->fourier_spacing, ftol, abstol);
1219 cmp_int(fp, "inputrec->nkx", -1, ir1->nkx, ir2->nkx);
1220 cmp_int(fp, "inputrec->nky", -1, ir1->nky, ir2->nky);
1221 cmp_int(fp, "inputrec->nkz", -1, ir1->nkz, ir2->nkz);
1222 cmp_int(fp, "inputrec->pme_order", -1, ir1->pme_order, ir2->pme_order);
1223 cmp_real(fp, "inputrec->ewald_rtol", -1, ir1->ewald_rtol, ir2->ewald_rtol, ftol, abstol);
1224 cmp_int(fp, "inputrec->ewald_geometry", -1, ir1->ewald_geometry, ir2->ewald_geometry);
1225 cmp_real(fp, "inputrec->epsilon_surface", -1, ir1->epsilon_surface, ir2->epsilon_surface, ftol, abstol);
1226 cmp_int(fp, "inputrec->bContinuation", -1, ir1->bContinuation, ir2->bContinuation);
1227 cmp_int(fp, "inputrec->bShakeSOR", -1, ir1->bShakeSOR, ir2->bShakeSOR);
1228 cmp_int(fp, "inputrec->etc", -1, ir1->etc, ir2->etc);
1229 cmp_int(fp, "inputrec->bPrintNHChains", -1, ir1->bPrintNHChains, ir2->bPrintNHChains);
1230 cmp_int(fp, "inputrec->epc", -1, ir1->epc, ir2->epc);
1231 cmp_int(fp, "inputrec->epct", -1, ir1->epct, ir2->epct);
1232 cmp_real(fp, "inputrec->tau_p", -1, ir1->tau_p, ir2->tau_p, ftol, abstol);
1233 cmp_rvec(fp, "inputrec->ref_p(x)", -1, ir1->ref_p[XX], ir2->ref_p[XX], ftol, abstol);
1234 cmp_rvec(fp, "inputrec->ref_p(y)", -1, ir1->ref_p[YY], ir2->ref_p[YY], ftol, abstol);
1235 cmp_rvec(fp, "inputrec->ref_p(z)", -1, ir1->ref_p[ZZ], ir2->ref_p[ZZ], ftol, abstol);
1236 cmp_rvec(fp, "inputrec->compress(x)", -1, ir1->compress[XX], ir2->compress[XX], ftol, abstol);
1237 cmp_rvec(fp, "inputrec->compress(y)", -1, ir1->compress[YY], ir2->compress[YY], ftol, abstol);
1238 cmp_rvec(fp, "inputrec->compress(z)", -1, ir1->compress[ZZ], ir2->compress[ZZ], ftol, abstol);
1239 cmp_int(fp, "refcoord_scaling", -1, ir1->refcoord_scaling, ir2->refcoord_scaling);
1240 cmp_rvec(fp, "inputrec->posres_com", -1, ir1->posres_com, ir2->posres_com, ftol, abstol);
1241 cmp_rvec(fp, "inputrec->posres_comB", -1, ir1->posres_comB, ir2->posres_comB, ftol, abstol);
1242 cmp_real(fp, "inputrec->verletbuf_tol", -1, ir1->verletbuf_tol, ir2->verletbuf_tol, ftol, abstol);
1243 cmp_real(fp, "inputrec->rlist", -1, ir1->rlist, ir2->rlist, ftol, abstol);
1244 cmp_real(fp, "inputrec->rtpi", -1, ir1->rtpi, ir2->rtpi, ftol, abstol);
1245 cmp_int(fp, "inputrec->coulombtype", -1, ir1->coulombtype, ir2->coulombtype);
1246 cmp_int(fp, "inputrec->coulomb_modifier", -1, ir1->coulomb_modifier, ir2->coulomb_modifier);
1247 cmp_real(fp, "inputrec->rcoulomb_switch", -1, ir1->rcoulomb_switch, ir2->rcoulomb_switch, ftol, abstol);
1248 cmp_real(fp, "inputrec->rcoulomb", -1, ir1->rcoulomb, ir2->rcoulomb, ftol, abstol);
1249 cmp_int(fp, "inputrec->vdwtype", -1, ir1->vdwtype, ir2->vdwtype);
1250 cmp_int(fp, "inputrec->vdw_modifier", -1, ir1->vdw_modifier, ir2->vdw_modifier); cmp_real(fp, "inputrec->rvdw_switch", -1, ir1->rvdw_switch, ir2->rvdw_switch, ftol, abstol);
1251 cmp_real(fp, "inputrec->rvdw", -1, ir1->rvdw, ir2->rvdw, ftol, abstol);
1252 cmp_real(fp, "inputrec->epsilon_r", -1, ir1->epsilon_r, ir2->epsilon_r, ftol, abstol);
1253 cmp_real(fp, "inputrec->epsilon_rf", -1, ir1->epsilon_rf, ir2->epsilon_rf, ftol, abstol);
1254 cmp_real(fp, "inputrec->tabext", -1, ir1->tabext, ir2->tabext, ftol, abstol);
1255 cmp_int(fp, "inputrec->implicit_solvent", -1, ir1->implicit_solvent, ir2->implicit_solvent);
1256 cmp_int(fp, "inputrec->gb_algorithm", -1, ir1->gb_algorithm, ir2->gb_algorithm);
1257 cmp_int(fp, "inputrec->nstgbradii", -1, ir1->nstgbradii, ir2->nstgbradii);
1258 cmp_real(fp, "inputrec->rgbradii", -1, ir1->rgbradii, ir2->rgbradii, ftol, abstol);
1259 cmp_real(fp, "inputrec->gb_saltconc", -1, ir1->gb_saltconc, ir2->gb_saltconc, ftol, abstol);
1260 cmp_real(fp, "inputrec->gb_epsilon_solvent", -1, ir1->gb_epsilon_solvent, ir2->gb_epsilon_solvent, ftol, abstol);
1261 cmp_real(fp, "inputrec->gb_obc_alpha", -1, ir1->gb_obc_alpha, ir2->gb_obc_alpha, ftol, abstol);
1262 cmp_real(fp, "inputrec->gb_obc_beta", -1, ir1->gb_obc_beta, ir2->gb_obc_beta, ftol, abstol);
1263 cmp_real(fp, "inputrec->gb_obc_gamma", -1, ir1->gb_obc_gamma, ir2->gb_obc_gamma, ftol, abstol);
1264 cmp_real(fp, "inputrec->gb_dielectric_offset", -1, ir1->gb_dielectric_offset, ir2->gb_dielectric_offset, ftol, abstol);
1265 cmp_int(fp, "inputrec->sa_algorithm", -1, ir1->sa_algorithm, ir2->sa_algorithm);
1266 cmp_real(fp, "inputrec->sa_surface_tension", -1, ir1->sa_surface_tension, ir2->sa_surface_tension, ftol, abstol);
1268 cmp_int(fp, "inputrec->eDispCorr", -1, ir1->eDispCorr, ir2->eDispCorr);
1269 cmp_real(fp, "inputrec->shake_tol", -1, ir1->shake_tol, ir2->shake_tol, ftol, abstol);
1270 cmp_int(fp, "inputrec->efep", -1, ir1->efep, ir2->efep);
1271 cmp_fepvals(fp, ir1->fepvals, ir2->fepvals, ftol, abstol);
1272 cmp_int(fp, "inputrec->bSimTemp", -1, ir1->bSimTemp, ir2->bSimTemp);
1273 if ((ir1->bSimTemp == ir2->bSimTemp) && (ir1->bSimTemp))
1275 cmp_simtempvals(fp, ir1->simtempvals, ir2->simtempvals, std::min(ir1->fepvals->n_lambda, ir2->fepvals->n_lambda), ftol, abstol);
1277 cmp_int(fp, "inputrec->bExpanded", -1, ir1->bExpanded, ir2->bExpanded);
1278 if ((ir1->bExpanded == ir2->bExpanded) && (ir1->bExpanded))
1280 cmp_expandedvals(fp, ir1->expandedvals, ir2->expandedvals, std::min(ir1->fepvals->n_lambda, ir2->fepvals->n_lambda), ftol, abstol);
1282 cmp_int(fp, "inputrec->nwall", -1, ir1->nwall, ir2->nwall);
1283 cmp_int(fp, "inputrec->wall_type", -1, ir1->wall_type, ir2->wall_type);
1284 cmp_int(fp, "inputrec->wall_atomtype[0]", -1, ir1->wall_atomtype[0], ir2->wall_atomtype[0]);
1285 cmp_int(fp, "inputrec->wall_atomtype[1]", -1, ir1->wall_atomtype[1], ir2->wall_atomtype[1]);
1286 cmp_real(fp, "inputrec->wall_density[0]", -1, ir1->wall_density[0], ir2->wall_density[0], ftol, abstol);
1287 cmp_real(fp, "inputrec->wall_density[1]", -1, ir1->wall_density[1], ir2->wall_density[1], ftol, abstol);
1288 cmp_real(fp, "inputrec->wall_ewald_zfac", -1, ir1->wall_ewald_zfac, ir2->wall_ewald_zfac, ftol, abstol);
1290 cmp_bool(fp, "inputrec->bPull", -1, ir1->bPull, ir2->bPull);
1291 if (ir1->bPull && ir2->bPull)
1293 cmp_pull(fp);
1296 cmp_int(fp, "inputrec->eDisre", -1, ir1->eDisre, ir2->eDisre);
1297 cmp_real(fp, "inputrec->dr_fc", -1, ir1->dr_fc, ir2->dr_fc, ftol, abstol);
1298 cmp_int(fp, "inputrec->eDisreWeighting", -1, ir1->eDisreWeighting, ir2->eDisreWeighting);
1299 cmp_int(fp, "inputrec->bDisreMixed", -1, ir1->bDisreMixed, ir2->bDisreMixed);
1300 cmp_int(fp, "inputrec->nstdisreout", -1, ir1->nstdisreout, ir2->nstdisreout);
1301 cmp_real(fp, "inputrec->dr_tau", -1, ir1->dr_tau, ir2->dr_tau, ftol, abstol);
1302 cmp_real(fp, "inputrec->orires_fc", -1, ir1->orires_fc, ir2->orires_fc, ftol, abstol);
1303 cmp_real(fp, "inputrec->orires_tau", -1, ir1->orires_tau, ir2->orires_tau, ftol, abstol);
1304 cmp_int(fp, "inputrec->nstorireout", -1, ir1->nstorireout, ir2->nstorireout);
1305 cmp_real(fp, "inputrec->em_stepsize", -1, ir1->em_stepsize, ir2->em_stepsize, ftol, abstol);
1306 cmp_real(fp, "inputrec->em_tol", -1, ir1->em_tol, ir2->em_tol, ftol, abstol);
1307 cmp_int(fp, "inputrec->niter", -1, ir1->niter, ir2->niter);
1308 cmp_real(fp, "inputrec->fc_stepsize", -1, ir1->fc_stepsize, ir2->fc_stepsize, ftol, abstol);
1309 cmp_int(fp, "inputrec->nstcgsteep", -1, ir1->nstcgsteep, ir2->nstcgsteep);
1310 cmp_int(fp, "inputrec->nbfgscorr", 0, ir1->nbfgscorr, ir2->nbfgscorr);
1311 cmp_int(fp, "inputrec->eConstrAlg", -1, ir1->eConstrAlg, ir2->eConstrAlg);
1312 cmp_int(fp, "inputrec->nProjOrder", -1, ir1->nProjOrder, ir2->nProjOrder);
1313 cmp_real(fp, "inputrec->LincsWarnAngle", -1, ir1->LincsWarnAngle, ir2->LincsWarnAngle, ftol, abstol);
1314 cmp_int(fp, "inputrec->nLincsIter", -1, ir1->nLincsIter, ir2->nLincsIter);
1315 cmp_real(fp, "inputrec->bd_fric", -1, ir1->bd_fric, ir2->bd_fric, ftol, abstol);
1316 cmp_int64(fp, "inputrec->ld_seed", ir1->ld_seed, ir2->ld_seed);
1317 cmp_real(fp, "inputrec->cos_accel", -1, ir1->cos_accel, ir2->cos_accel, ftol, abstol);
1318 cmp_rvec(fp, "inputrec->deform(a)", -1, ir1->deform[XX], ir2->deform[XX], ftol, abstol);
1319 cmp_rvec(fp, "inputrec->deform(b)", -1, ir1->deform[YY], ir2->deform[YY], ftol, abstol);
1320 cmp_rvec(fp, "inputrec->deform(c)", -1, ir1->deform[ZZ], ir2->deform[ZZ], ftol, abstol);
1323 cmp_int(fp, "inputrec->userint1", -1, ir1->userint1, ir2->userint1);
1324 cmp_int(fp, "inputrec->userint2", -1, ir1->userint2, ir2->userint2);
1325 cmp_int(fp, "inputrec->userint3", -1, ir1->userint3, ir2->userint3);
1326 cmp_int(fp, "inputrec->userint4", -1, ir1->userint4, ir2->userint4);
1327 cmp_real(fp, "inputrec->userreal1", -1, ir1->userreal1, ir2->userreal1, ftol, abstol);
1328 cmp_real(fp, "inputrec->userreal2", -1, ir1->userreal2, ir2->userreal2, ftol, abstol);
1329 cmp_real(fp, "inputrec->userreal3", -1, ir1->userreal3, ir2->userreal3, ftol, abstol);
1330 cmp_real(fp, "inputrec->userreal4", -1, ir1->userreal4, ir2->userreal4, ftol, abstol);
1331 cmp_grpopts(fp, &(ir1->opts), &(ir2->opts), ftol, abstol);
1332 cmp_cosines(fp, "ex", ir1->ex, ir2->ex, ftol, abstol);
1333 cmp_cosines(fp, "et", ir1->et, ir2->et, ftol, abstol);
1336 void comp_pull_AB(FILE *fp, pull_params_t *pull, real ftol, real abstol)
1338 int i;
1340 for (i = 0; i < pull->ncoord; i++)
1342 fprintf(fp, "comparing pull coord %d\n", i);
1343 cmp_real(fp, "pull-coord->k", -1, pull->coord[i].k, pull->coord[i].kB, ftol, abstol);
1347 gmx_bool inputrecDeform(const t_inputrec *ir)
1349 return (ir->deform[XX][XX] != 0 || ir->deform[YY][YY] != 0 || ir->deform[ZZ][ZZ] != 0 ||
1350 ir->deform[YY][XX] != 0 || ir->deform[ZZ][XX] != 0 || ir->deform[ZZ][YY] != 0);
1353 gmx_bool inputrecDynamicBox(const t_inputrec *ir)
1355 return (ir->epc != epcNO || ir->eI == eiTPI || inputrecDeform(ir));
1358 gmx_bool inputrecPreserveShape(const t_inputrec *ir)
1360 return (ir->epc != epcNO && ir->deform[XX][XX] == 0 &&
1361 (ir->epct == epctISOTROPIC || ir->epct == epctSEMIISOTROPIC));
1364 gmx_bool inputrecNeedMutot(const t_inputrec *ir)
1366 return ((ir->coulombtype == eelEWALD || EEL_PME(ir->coulombtype)) &&
1367 (ir->ewald_geometry == eewg3DC || ir->epsilon_surface != 0));
1370 gmx_bool inputrecElecField(const t_inputrec *ir)
1372 return (ir->ex[XX].n > 0 || ir->ex[YY].n > 0 || ir->ex[ZZ].n > 0);
1375 gmx_bool inputrecExclForces(const t_inputrec *ir)
1377 return (EEL_FULL(ir->coulombtype) || (EEL_RF(ir->coulombtype)) ||
1378 ir->implicit_solvent != eisNO);
1381 gmx_bool inputrecNptTrotter(const t_inputrec *ir)
1383 return ( ( (ir->eI == eiVV) || (ir->eI == eiVVAK) ) &&
1384 (ir->epc == epcMTTK) && (ir->etc == etcNOSEHOOVER) );
1387 gmx_bool inputrecNvtTrotter(const t_inputrec *ir)
1389 return ( ( (ir->eI == eiVV) || (ir->eI == eiVVAK) ) &&
1390 (ir->epc != epcMTTK) && (ir->etc == etcNOSEHOOVER) );
1393 gmx_bool inputrecNphTrotter(const t_inputrec *ir)
1395 return ( ( (ir->eI == eiVV) || (ir->eI == eiVVAK) ) &&
1396 (ir->epc == epcMTTK) && (ir->etc != etcNOSEHOOVER) );