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,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
48 #include "gromacs/math/veccompare.h"
49 #include "gromacs/math/vecdump.h"
50 #include "gromacs/mdtypes/awh_params.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/mdtypes/pull_params.h"
53 #include "gromacs/pbcutil/pbc.h"
54 #include "gromacs/utility/compare.h"
55 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/keyvaluetree.h"
58 #include "gromacs/utility/smalloc.h"
59 #include "gromacs/utility/snprintf.h"
60 #include "gromacs/utility/strconvert.h"
61 #include "gromacs/utility/stringutil.h"
62 #include "gromacs/utility/textwriter.h"
63 #include "gromacs/utility/txtdump.h"
65 //! Macro to select a bool name
66 #define EBOOL(e) gmx::boolToString(e)
68 /* The minimum number of integration steps required for reasonably accurate
69 * integration of first and second order coupling algorithms.
71 const int nstmin_berendsen_tcouple = 5;
72 const int nstmin_berendsen_pcouple = 10;
73 const int nstmin_harmonic = 20;
75 t_inputrec::t_inputrec()
77 // TODO When this memset is removed, remove the suppression of
78 // gcc -Wno-class-memaccess in a CMakeLists.txt file.
79 std::memset(this, 0, sizeof(*this)); // NOLINT(bugprone-undefined-memory-manipulation)
81 snew(expandedvals, 1);
85 t_inputrec::~t_inputrec()
90 static int nst_wanted(const t_inputrec* ir)
102 int ir_optimal_nstcalcenergy(const t_inputrec* ir)
104 return nst_wanted(ir);
107 int tcouple_min_integration_steps(int etc)
113 case etcNO: n = 0; break;
115 case etcYES: n = nstmin_berendsen_tcouple; break;
117 /* V-rescale supports instantaneous rescaling */
120 case etcNOSEHOOVER: n = nstmin_harmonic; break;
122 case etcANDERSENMASSIVE: n = 1; break;
123 default: gmx_incons("Unknown etc value");
129 int ir_optimal_nsttcouple(const t_inputrec* ir)
131 int nmin, nwanted, n;
135 nmin = tcouple_min_integration_steps(ir->etc);
137 nwanted = nst_wanted(ir);
140 if (ir->etc != etcNO)
142 for (g = 0; g < ir->opts.ngtc; g++)
144 if (ir->opts.tau_t[g] > 0)
146 tau_min = std::min(tau_min, ir->opts.tau_t[g]);
151 if (nmin == 0 || ir->delta_t * nwanted <= tau_min)
157 n = static_cast<int>(tau_min / (ir->delta_t * nmin) + 0.001);
162 while (nwanted % n != 0)
171 int pcouple_min_integration_steps(int epc)
177 case epcNO: n = 0; break;
179 case epcISOTROPIC: n = nstmin_berendsen_pcouple; break;
180 case epcPARRINELLORAHMAN:
181 case epcMTTK: n = nstmin_harmonic; break;
182 default: gmx_incons("Unknown epc value");
188 int ir_optimal_nstpcouple(const t_inputrec* ir)
190 int nmin, nwanted, n;
192 nmin = pcouple_min_integration_steps(ir->epc);
194 nwanted = nst_wanted(ir);
196 if (nmin == 0 || ir->delta_t * nwanted <= ir->tau_p)
202 n = static_cast<int>(ir->tau_p / (ir->delta_t * nmin) + 0.001);
207 while (nwanted % n != 0)
216 gmx_bool ir_coulomb_switched(const t_inputrec* ir)
218 return (ir->coulombtype == eelSWITCH || ir->coulombtype == eelSHIFT
219 || ir->coulombtype == eelPMESWITCH || ir->coulombtype == eelPMEUSERSWITCH
220 || ir->coulomb_modifier == eintmodPOTSWITCH || ir->coulomb_modifier == eintmodFORCESWITCH);
223 gmx_bool ir_coulomb_is_zero_at_cutoff(const t_inputrec* ir)
225 return (ir->cutoff_scheme == ecutsVERLET || ir_coulomb_switched(ir)
226 || ir->coulomb_modifier != eintmodNONE || ir->coulombtype == eelRF_ZERO);
229 gmx_bool ir_coulomb_might_be_zero_at_cutoff(const t_inputrec* ir)
231 return (ir_coulomb_is_zero_at_cutoff(ir) || ir->coulombtype == eelUSER || ir->coulombtype == eelPMEUSER);
234 gmx_bool ir_vdw_switched(const t_inputrec* ir)
236 return (ir->vdwtype == evdwSWITCH || ir->vdwtype == evdwSHIFT
237 || ir->vdw_modifier == eintmodPOTSWITCH || ir->vdw_modifier == eintmodFORCESWITCH);
240 gmx_bool ir_vdw_is_zero_at_cutoff(const t_inputrec* ir)
242 return (ir->cutoff_scheme == ecutsVERLET || ir_vdw_switched(ir) || ir->vdw_modifier != eintmodNONE);
245 gmx_bool ir_vdw_might_be_zero_at_cutoff(const t_inputrec* ir)
247 return (ir_vdw_is_zero_at_cutoff(ir) || ir->vdwtype == evdwUSER);
250 static void done_pull_group(t_pull_group* pgrp)
259 static void done_pull_params(pull_params_t* pull)
263 for (i = 0; i < pull->ngroup + 1; i++)
265 done_pull_group(pull->group);
272 static void done_lambdas(t_lambda* fep)
274 if (fep->n_lambda > 0)
276 for (int i = 0; i < efptNR; i++)
278 sfree(fep->all_lambda[i]);
281 sfree(fep->all_lambda);
284 void done_inputrec(t_inputrec* ir)
286 sfree(ir->opts.nrdf);
287 sfree(ir->opts.ref_t);
288 for (int i = 0; i < ir->opts.ngtc; i++)
290 sfree(ir->opts.anneal_time[i]);
291 sfree(ir->opts.anneal_temp[i]);
293 sfree(ir->opts.annealing);
294 sfree(ir->opts.anneal_npoints);
295 sfree(ir->opts.anneal_time);
296 sfree(ir->opts.anneal_temp);
297 sfree(ir->opts.tau_t);
299 sfree(ir->opts.nFreeze);
300 sfree(ir->opts.egp_flags);
301 done_lambdas(ir->fepvals);
303 sfree(ir->expandedvals);
304 sfree(ir->simtempvals);
308 done_pull_params(ir->pull);
314 static void pr_grp_opts(FILE* out, int indent, const char* title, const t_grpopts* opts, gmx_bool bMDPformat)
320 fprintf(out, "%s:\n", title);
323 pr_indent(out, indent);
324 fprintf(out, "nrdf%s", bMDPformat ? " = " : ":");
325 for (i = 0; (i < opts->ngtc); i++)
327 fprintf(out, " %10g", opts->nrdf[i]);
331 pr_indent(out, indent);
332 fprintf(out, "ref-t%s", bMDPformat ? " = " : ":");
333 for (i = 0; (i < opts->ngtc); i++)
335 fprintf(out, " %10g", opts->ref_t[i]);
339 pr_indent(out, indent);
340 fprintf(out, "tau-t%s", bMDPformat ? " = " : ":");
341 for (i = 0; (i < opts->ngtc); i++)
343 fprintf(out, " %10g", opts->tau_t[i]);
347 /* Pretty-print the simulated annealing info */
348 fprintf(out, "annealing%s", bMDPformat ? " = " : ":");
349 for (i = 0; (i < opts->ngtc); i++)
351 fprintf(out, " %10s", EANNEAL(opts->annealing[i]));
355 fprintf(out, "annealing-npoints%s", bMDPformat ? " = " : ":");
356 for (i = 0; (i < opts->ngtc); i++)
358 fprintf(out, " %10d", opts->anneal_npoints[i]);
362 for (i = 0; (i < opts->ngtc); i++)
364 if (opts->anneal_npoints[i] > 0)
366 fprintf(out, "annealing-time [%d]:\t", i);
367 for (j = 0; (j < opts->anneal_npoints[i]); j++)
369 fprintf(out, " %10.1f", opts->anneal_time[i][j]);
372 fprintf(out, "annealing-temp [%d]:\t", i);
373 for (j = 0; (j < opts->anneal_npoints[i]); j++)
375 fprintf(out, " %10.1f", opts->anneal_temp[i][j]);
381 pr_indent(out, indent);
382 fprintf(out, "acc:\t");
383 for (i = 0; (i < opts->ngacc); i++)
385 for (m = 0; (m < DIM); m++)
387 fprintf(out, " %10g", opts->acc[i][m]);
392 pr_indent(out, indent);
393 fprintf(out, "nfreeze:");
394 for (i = 0; (i < opts->ngfrz); i++)
396 for (m = 0; (m < DIM); m++)
398 fprintf(out, " %10s", opts->nFreeze[i][m] ? "Y" : "N");
404 for (i = 0; (i < opts->ngener); i++)
406 pr_indent(out, indent);
407 fprintf(out, "energygrp-flags[%3d]:", i);
408 for (m = 0; (m < opts->ngener); m++)
410 fprintf(out, " %d", opts->egp_flags[opts->ngener * i + m]);
418 static void pr_matrix(FILE* fp, int indent, const char* title, const rvec* m, gmx_bool bMDPformat)
422 fprintf(fp, "%-10s = %g %g %g %g %g %g\n", title, m[XX][XX], m[YY][YY], m[ZZ][ZZ],
423 m[XX][YY], m[XX][ZZ], m[YY][ZZ]);
427 pr_rvecs(fp, indent, title, m, DIM);
431 #define PS(t, s) pr_str(fp, indent, t, s)
432 #define PI(t, s) pr_int(fp, indent, t, s)
433 #define PSTEP(t, s) pr_int64(fp, indent, t, s)
434 #define PR(t, s) pr_real(fp, indent, t, s)
435 #define PD(t, s) pr_double(fp, indent, t, s)
437 static void pr_pull_group(FILE* fp, int indent, int g, const t_pull_group* pgrp)
439 pr_indent(fp, indent);
440 fprintf(fp, "pull-group %d:\n", g);
442 pr_ivec_block(fp, indent, "atom", pgrp->ind, pgrp->nat, TRUE);
443 pr_rvec(fp, indent, "weight", pgrp->weight, pgrp->nweight, TRUE);
444 PI("pbcatom", pgrp->pbcatom);
447 static void pr_pull_coord(FILE* fp, int indent, int c, const t_pull_coord* pcrd)
451 pr_indent(fp, indent);
452 fprintf(fp, "pull-coord %d:\n", c);
453 PS("type", EPULLTYPE(pcrd->eType));
454 if (pcrd->eType == epullEXTERNAL)
456 PS("potential-provider", pcrd->externalPotentialProvider);
458 PS("geometry", EPULLGEOM(pcrd->eGeom));
459 for (g = 0; g < pcrd->ngroup; g++)
463 sprintf(buf, "group[%d]", g);
464 PI(buf, pcrd->group[g]);
466 pr_ivec(fp, indent, "dim", pcrd->dim, DIM, TRUE);
467 pr_rvec(fp, indent, "origin", pcrd->origin, DIM, TRUE);
468 pr_rvec(fp, indent, "vec", pcrd->vec, DIM, TRUE);
469 PS("start", EBOOL(pcrd->bStart));
470 PR("init", pcrd->init);
471 PR("rate", pcrd->rate);
476 static void pr_simtempvals(FILE* fp, int indent, const t_simtemp* simtemp, int n_lambda)
478 PS("simulated-tempering-scaling", ESIMTEMP(simtemp->eSimTempScale));
479 PR("sim-temp-low", simtemp->simtemp_low);
480 PR("sim-temp-high", simtemp->simtemp_high);
481 pr_rvec(fp, indent, "simulated tempering temperatures", simtemp->temperatures, n_lambda, TRUE);
484 static void pr_expandedvals(FILE* fp, int indent, const t_expanded* expand, int n_lambda)
487 PI("nstexpanded", expand->nstexpanded);
488 PS("lmc-stats", elamstats_names[expand->elamstats]);
489 PS("lmc-move", elmcmove_names[expand->elmcmove]);
490 PS("lmc-weights-equil", elmceq_names[expand->elmceq]);
491 if (expand->elmceq == elmceqNUMATLAM)
493 PI("weight-equil-number-all-lambda", expand->equil_n_at_lam);
495 if (expand->elmceq == elmceqSAMPLES)
497 PI("weight-equil-number-samples", expand->equil_samples);
499 if (expand->elmceq == elmceqSTEPS)
501 PI("weight-equil-number-steps", expand->equil_steps);
503 if (expand->elmceq == elmceqWLDELTA)
505 PR("weight-equil-wl-delta", expand->equil_wl_delta);
507 if (expand->elmceq == elmceqRATIO)
509 PR("weight-equil-count-ratio", expand->equil_ratio);
511 PI("lmc-seed", expand->lmc_seed);
512 PR("mc-temperature", expand->mc_temp);
513 PI("lmc-repeats", expand->lmc_repeats);
514 PI("lmc-gibbsdelta", expand->gibbsdeltalam);
515 PI("lmc-forced-nstart", expand->lmc_forced_nstart);
516 PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
517 PI("nst-transition-matrix", expand->nstTij);
518 PI("mininum-var-min", expand->minvarmin); /*default is reasonable */
519 PI("weight-c-range", expand->c_range); /* default is just C=0 */
520 PR("wl-scale", expand->wl_scale);
521 PR("wl-ratio", expand->wl_ratio);
522 PR("init-wl-delta", expand->init_wl_delta);
523 PS("wl-oneovert", EBOOL(expand->bWLoneovert));
525 pr_indent(fp, indent);
526 pr_rvec(fp, indent, "init-lambda-weights", expand->init_lambda_weights, n_lambda, TRUE);
527 PS("init-weights", EBOOL(expand->bInit_weights));
530 static void pr_fepvals(FILE* fp, int indent, const t_lambda* fep, gmx_bool bMDPformat)
534 PR("init-lambda", fep->init_lambda);
535 PI("init-lambda-state", fep->init_fep_state);
536 PR("delta-lambda", fep->delta_lambda);
537 PI("nstdhdl", fep->nstdhdl);
541 PI("n-lambdas", fep->n_lambda);
543 if (fep->n_lambda > 0)
545 pr_indent(fp, indent);
546 fprintf(fp, "separate-dvdl%s\n", bMDPformat ? " = " : ":");
547 for (i = 0; i < efptNR; i++)
549 fprintf(fp, "%18s = ", efpt_names[i]);
550 if (fep->separate_dvdl[i])
552 fprintf(fp, " TRUE");
556 fprintf(fp, " FALSE");
560 fprintf(fp, "all-lambdas%s\n", bMDPformat ? " = " : ":");
561 for (i = 0; i < efptNR; i++)
563 fprintf(fp, "%18s = ", efpt_names[i]);
564 for (j = 0; j < fep->n_lambda; j++)
566 fprintf(fp, " %10g", fep->all_lambda[i][j]);
571 PI("calc-lambda-neighbors", fep->lambda_neighbors);
572 PS("dhdl-print-energy", edHdLPrintEnergy_names[fep->edHdLPrintEnergy]);
573 PR("sc-alpha", fep->sc_alpha);
574 PI("sc-power", fep->sc_power);
575 PR("sc-r-power", fep->sc_r_power);
576 PR("sc-sigma", fep->sc_sigma);
577 PR("sc-sigma-min", fep->sc_sigma_min);
578 PS("sc-coul", EBOOL(fep->bScCoul));
579 PI("dh-hist-size", fep->dh_hist_size);
580 PD("dh-hist-spacing", fep->dh_hist_spacing);
581 PS("separate-dhdl-file", SEPDHDLFILETYPE(fep->separate_dhdl_file));
582 PS("dhdl-derivatives", DHDLDERIVATIVESTYPE(fep->dhdl_derivatives));
585 static void pr_pull(FILE* fp, int indent, const pull_params_t* pull)
589 PR("pull-cylinder-r", pull->cylinder_r);
590 PR("pull-constr-tol", pull->constr_tol);
591 PS("pull-print-COM", EBOOL(pull->bPrintCOM));
592 PS("pull-print-ref-value", EBOOL(pull->bPrintRefValue));
593 PS("pull-print-components", EBOOL(pull->bPrintComp));
594 PI("pull-nstxout", pull->nstxout);
595 PI("pull-nstfout", pull->nstfout);
596 PS("pull-pbc-ref-prev-step-com", EBOOL(pull->bSetPbcRefToPrevStepCOM));
597 PS("pull-xout-average", EBOOL(pull->bXOutAverage));
598 PS("pull-fout-average", EBOOL(pull->bFOutAverage));
599 PI("pull-ngroups", pull->ngroup);
600 for (g = 0; g < pull->ngroup; g++)
602 pr_pull_group(fp, indent, g, &pull->group[g]);
604 PI("pull-ncoords", pull->ncoord);
605 for (g = 0; g < pull->ncoord; g++)
607 pr_pull_coord(fp, indent, g, &pull->coord[g]);
611 static void pr_awh_bias_dim(FILE* fp, int indent, gmx::AwhDimParams* awhDimParams, const char* prefix)
613 pr_indent(fp, indent);
615 fprintf(fp, "%s:\n", prefix);
616 PS("coord-provider", EAWHCOORDPROVIDER(awhDimParams->eCoordProvider));
617 PI("coord-index", awhDimParams->coordIndex + 1);
618 PR("start", awhDimParams->origin);
619 PR("end", awhDimParams->end);
620 PR("period", awhDimParams->period);
621 PR("force-constant", awhDimParams->forceConstant);
622 PR("diffusion", awhDimParams->diffusion);
623 PR("cover-diameter", awhDimParams->coverDiameter);
626 static void pr_awh_bias(FILE* fp, int indent, gmx::AwhBiasParams* awhBiasParams, const char* prefix)
630 sprintf(opt, "%s-error-init", prefix);
631 PR(opt, awhBiasParams->errorInitial);
632 sprintf(opt, "%s-growth", prefix);
633 PS(opt, EAWHGROWTH(awhBiasParams->eGrowth));
634 sprintf(opt, "%s-target", prefix);
635 PS(opt, EAWHTARGET(awhBiasParams->eTarget));
636 sprintf(opt, "%s-target-beta-scalng", prefix);
637 PR(opt, awhBiasParams->targetBetaScaling);
638 sprintf(opt, "%s-target-cutoff", prefix);
639 PR(opt, awhBiasParams->targetCutoff);
640 sprintf(opt, "%s-user-data", prefix);
641 PS(opt, EBOOL(awhBiasParams->bUserData));
642 sprintf(opt, "%s-share-group", prefix);
643 PI(opt, awhBiasParams->shareGroup);
644 sprintf(opt, "%s-equilibrate-histogram", prefix);
645 PS(opt, EBOOL(awhBiasParams->equilibrateHistogram));
646 sprintf(opt, "%s-ndim", prefix);
647 PI(opt, awhBiasParams->ndim);
649 for (int d = 0; d < awhBiasParams->ndim; d++)
651 char prefixdim[STRLEN];
652 sprintf(prefixdim, "%s-dim%d", prefix, d + 1);
653 pr_awh_bias_dim(fp, indent, &awhBiasParams->dimParams[d], prefixdim);
657 static void pr_awh(FILE* fp, int indent, gmx::AwhParams* awhParams)
659 PS("awh-potential", EAWHPOTENTIAL(awhParams->ePotential));
660 PI("awh-seed", awhParams->seed);
661 PI("awh-nstout", awhParams->nstOut);
662 PI("awh-nstsample", awhParams->nstSampleCoord);
663 PI("awh-nsamples-update", awhParams->numSamplesUpdateFreeEnergy);
664 PS("awh-share-bias-multisim", EBOOL(awhParams->shareBiasMultisim));
665 PI("awh-nbias", awhParams->numBias);
667 for (int k = 0; k < awhParams->numBias; k++)
669 auto prefix = gmx::formatString("awh%d", k + 1);
670 pr_awh_bias(fp, indent, &awhParams->awhBiasParams[k], prefix.c_str());
674 static void pr_rotgrp(FILE* fp, int indent, int g, const t_rotgrp* rotg)
676 pr_indent(fp, indent);
677 fprintf(fp, "rot-group %d:\n", g);
679 PS("rot-type", EROTGEOM(rotg->eType));
680 PS("rot-massw", EBOOL(rotg->bMassW));
681 pr_ivec_block(fp, indent, "atom", rotg->ind, rotg->nat, TRUE);
682 pr_rvecs(fp, indent, "x-ref", rotg->x_ref, rotg->nat);
683 pr_rvec(fp, indent, "rot-vec", rotg->inputVec, DIM, TRUE);
684 pr_rvec(fp, indent, "rot-pivot", rotg->pivot, DIM, TRUE);
685 PR("rot-rate", rotg->rate);
686 PR("rot-k", rotg->k);
687 PR("rot-slab-dist", rotg->slab_dist);
688 PR("rot-min-gauss", rotg->min_gaussian);
689 PR("rot-eps", rotg->eps);
690 PS("rot-fit-method", EROTFIT(rotg->eFittype));
691 PI("rot-potfit-nstep", rotg->PotAngle_nstep);
692 PR("rot-potfit-step", rotg->PotAngle_step);
695 static void pr_rot(FILE* fp, int indent, const t_rot* rot)
699 PI("rot-nstrout", rot->nstrout);
700 PI("rot-nstsout", rot->nstsout);
701 PI("rot-ngroups", rot->ngrp);
702 for (g = 0; g < rot->ngrp; g++)
704 pr_rotgrp(fp, indent, g, &rot->grp[g]);
709 static void pr_swap(FILE* fp, int indent, const t_swapcoords* swap)
713 /* Enums for better readability of the code */
721 PI("swap-frequency", swap->nstswap);
723 /* The split groups that define the compartments */
724 for (int j = 0; j < 2; j++)
726 snprintf(str, STRLEN, "massw_split%d", j);
727 PS(str, EBOOL(swap->massw_split[j]));
728 snprintf(str, STRLEN, "split atoms group %d", j);
729 pr_ivec_block(fp, indent, str, swap->grp[j].ind, swap->grp[j].nat, TRUE);
732 /* The solvent group */
733 snprintf(str, STRLEN, "solvent group %s", swap->grp[eGrpSolvent].molname);
734 pr_ivec_block(fp, indent, str, swap->grp[eGrpSolvent].ind, swap->grp[eGrpSolvent].nat, TRUE);
736 /* Now print the indices for all the ion groups: */
737 for (int ig = eSwapFixedGrpNR; ig < swap->ngrp; ig++)
739 snprintf(str, STRLEN, "ion group %s", swap->grp[ig].molname);
740 pr_ivec_block(fp, indent, str, swap->grp[ig].ind, swap->grp[ig].nat, TRUE);
743 PR("cyl0-r", swap->cyl0r);
744 PR("cyl0-up", swap->cyl0u);
745 PR("cyl0-down", swap->cyl0l);
746 PR("cyl1-r", swap->cyl1r);
747 PR("cyl1-up", swap->cyl1u);
748 PR("cyl1-down", swap->cyl1l);
749 PI("coupl-steps", swap->nAverage);
751 /* Print the requested ion counts for both compartments */
752 for (int ic = eCompA; ic <= eCompB; ic++)
754 for (int ig = eSwapFixedGrpNR; ig < swap->ngrp; ig++)
756 snprintf(str, STRLEN, "%s-in-%c", swap->grp[ig].molname, 'A' + ic);
757 PI(str, swap->grp[ig].nmolReq[ic]);
761 PR("threshold", swap->threshold);
762 PR("bulk-offsetA", swap->bulkOffset[eCompA]);
763 PR("bulk-offsetB", swap->bulkOffset[eCompB]);
767 static void pr_imd(FILE* fp, int indent, const t_IMD* imd)
769 PI("IMD-atoms", imd->nat);
770 pr_ivec_block(fp, indent, "atom", imd->ind, imd->nat, TRUE);
774 void pr_inputrec(FILE* fp, int indent, const char* title, const t_inputrec* ir, gmx_bool bMDPformat)
776 const char* infbuf = "inf";
778 if (available(fp, ir, indent, title))
782 indent = pr_title(fp, indent, title);
784 /* Try to make this list appear in the same order as the
785 * options are written in the default mdout.mdp, and with
786 * the same user-exposed names to facilitate debugging.
788 PS("integrator", EI(ir->eI));
789 PR("tinit", ir->init_t);
790 PR("dt", ir->delta_t);
791 PSTEP("nsteps", ir->nsteps);
792 PSTEP("init-step", ir->init_step);
793 PI("simulation-part", ir->simulation_part);
794 PS("comm-mode", ECOM(ir->comm_mode));
795 PI("nstcomm", ir->nstcomm);
797 /* Langevin dynamics */
798 PR("bd-fric", ir->bd_fric);
799 PSTEP("ld-seed", ir->ld_seed);
801 /* Energy minimization */
802 PR("emtol", ir->em_tol);
803 PR("emstep", ir->em_stepsize);
804 PI("niter", ir->niter);
805 PR("fcstep", ir->fc_stepsize);
806 PI("nstcgsteep", ir->nstcgsteep);
807 PI("nbfgscorr", ir->nbfgscorr);
809 /* Test particle insertion */
810 PR("rtpi", ir->rtpi);
813 PI("nstxout", ir->nstxout);
814 PI("nstvout", ir->nstvout);
815 PI("nstfout", ir->nstfout);
816 PI("nstlog", ir->nstlog);
817 PI("nstcalcenergy", ir->nstcalcenergy);
818 PI("nstenergy", ir->nstenergy);
819 PI("nstxout-compressed", ir->nstxout_compressed);
820 PR("compressed-x-precision", ir->x_compression_precision);
822 /* Neighborsearching parameters */
823 PS("cutoff-scheme", ECUTSCHEME(ir->cutoff_scheme));
824 PI("nstlist", ir->nstlist);
825 PS("pbc", c_pbcTypeNames[ir->pbcType].c_str());
826 PS("periodic-molecules", EBOOL(ir->bPeriodicMols));
827 PR("verlet-buffer-tolerance", ir->verletbuf_tol);
828 PR("rlist", ir->rlist);
830 /* Options for electrostatics and VdW */
831 PS("coulombtype", EELTYPE(ir->coulombtype));
832 PS("coulomb-modifier", INTMODIFIER(ir->coulomb_modifier));
833 PR("rcoulomb-switch", ir->rcoulomb_switch);
834 PR("rcoulomb", ir->rcoulomb);
835 if (ir->epsilon_r != 0)
837 PR("epsilon-r", ir->epsilon_r);
841 PS("epsilon-r", infbuf);
843 if (ir->epsilon_rf != 0)
845 PR("epsilon-rf", ir->epsilon_rf);
849 PS("epsilon-rf", infbuf);
851 PS("vdw-type", EVDWTYPE(ir->vdwtype));
852 PS("vdw-modifier", INTMODIFIER(ir->vdw_modifier));
853 PR("rvdw-switch", ir->rvdw_switch);
854 PR("rvdw", ir->rvdw);
855 PS("DispCorr", EDISPCORR(ir->eDispCorr));
856 PR("table-extension", ir->tabext);
858 PR("fourierspacing", ir->fourier_spacing);
859 PI("fourier-nx", ir->nkx);
860 PI("fourier-ny", ir->nky);
861 PI("fourier-nz", ir->nkz);
862 PI("pme-order", ir->pme_order);
863 PR("ewald-rtol", ir->ewald_rtol);
864 PR("ewald-rtol-lj", ir->ewald_rtol_lj);
865 PS("lj-pme-comb-rule", ELJPMECOMBNAMES(ir->ljpme_combination_rule));
866 PR("ewald-geometry", ir->ewald_geometry);
867 PR("epsilon-surface", ir->epsilon_surface);
869 /* Options for weak coupling algorithms */
870 PS("tcoupl", ETCOUPLTYPE(ir->etc));
871 PI("nsttcouple", ir->nsttcouple);
872 PI("nh-chain-length", ir->opts.nhchainlength);
873 PS("print-nose-hoover-chain-variables", EBOOL(ir->bPrintNHChains));
875 PS("pcoupl", EPCOUPLTYPE(ir->epc));
876 PS("pcoupltype", EPCOUPLTYPETYPE(ir->epct));
877 PI("nstpcouple", ir->nstpcouple);
878 PR("tau-p", ir->tau_p);
879 pr_matrix(fp, indent, "compressibility", ir->compress, bMDPformat);
880 pr_matrix(fp, indent, "ref-p", ir->ref_p, bMDPformat);
881 PS("refcoord-scaling", EREFSCALINGTYPE(ir->refcoord_scaling));
885 fprintf(fp, "posres-com = %g %g %g\n", ir->posres_com[XX], ir->posres_com[YY],
887 fprintf(fp, "posres-comB = %g %g %g\n", ir->posres_comB[XX], ir->posres_comB[YY],
888 ir->posres_comB[ZZ]);
892 pr_rvec(fp, indent, "posres-com", ir->posres_com, DIM, TRUE);
893 pr_rvec(fp, indent, "posres-comB", ir->posres_comB, DIM, TRUE);
897 PS("QMMM", EBOOL(ir->bQMMM));
898 fprintf(fp, "%s:\n", "qm-opts");
899 pr_int(fp, indent, "ngQM", ir->opts.ngQM);
901 /* CONSTRAINT OPTIONS */
902 PS("constraint-algorithm", ECONSTRTYPE(ir->eConstrAlg));
903 PS("continuation", EBOOL(ir->bContinuation));
905 PS("Shake-SOR", EBOOL(ir->bShakeSOR));
906 PR("shake-tol", ir->shake_tol);
907 PI("lincs-order", ir->nProjOrder);
908 PI("lincs-iter", ir->nLincsIter);
909 PR("lincs-warnangle", ir->LincsWarnAngle);
912 PI("nwall", ir->nwall);
913 PS("wall-type", EWALLTYPE(ir->wall_type));
914 PR("wall-r-linpot", ir->wall_r_linpot);
916 PI("wall-atomtype[0]", ir->wall_atomtype[0]);
917 PI("wall-atomtype[1]", ir->wall_atomtype[1]);
919 PR("wall-density[0]", ir->wall_density[0]);
920 PR("wall-density[1]", ir->wall_density[1]);
921 PR("wall-ewald-zfac", ir->wall_ewald_zfac);
924 PS("pull", EBOOL(ir->bPull));
927 pr_pull(fp, indent, ir->pull);
931 PS("awh", EBOOL(ir->bDoAwh));
934 pr_awh(fp, indent, ir->awhParams);
937 /* ENFORCED ROTATION */
938 PS("rotation", EBOOL(ir->bRot));
941 pr_rot(fp, indent, ir->rot);
945 PS("interactiveMD", EBOOL(ir->bIMD));
948 pr_imd(fp, indent, ir->imd);
951 /* NMR refinement stuff */
952 PS("disre", EDISRETYPE(ir->eDisre));
953 PS("disre-weighting", EDISREWEIGHTING(ir->eDisreWeighting));
954 PS("disre-mixed", EBOOL(ir->bDisreMixed));
955 PR("dr-fc", ir->dr_fc);
956 PR("dr-tau", ir->dr_tau);
957 PR("nstdisreout", ir->nstdisreout);
959 PR("orire-fc", ir->orires_fc);
960 PR("orire-tau", ir->orires_tau);
961 PR("nstorireout", ir->nstorireout);
963 /* FREE ENERGY VARIABLES */
964 PS("free-energy", EFEPTYPE(ir->efep));
965 if (ir->efep != efepNO || ir->bSimTemp)
967 pr_fepvals(fp, indent, ir->fepvals, bMDPformat);
971 pr_expandedvals(fp, indent, ir->expandedvals, ir->fepvals->n_lambda);
974 /* NON-equilibrium MD stuff */
975 PR("cos-acceleration", ir->cos_accel);
976 pr_matrix(fp, indent, "deform", ir->deform, bMDPformat);
978 /* SIMULATED TEMPERING */
979 PS("simulated-tempering", EBOOL(ir->bSimTemp));
982 pr_simtempvals(fp, indent, ir->simtempvals, ir->fepvals->n_lambda);
985 /* ION/WATER SWAPPING FOR COMPUTATIONAL ELECTROPHYSIOLOGY */
986 PS("swapcoords", ESWAPTYPE(ir->eSwapCoords));
987 if (ir->eSwapCoords != eswapNO)
989 pr_swap(fp, indent, ir->swap);
992 /* USER-DEFINED THINGIES */
993 PI("userint1", ir->userint1);
994 PI("userint2", ir->userint2);
995 PI("userint3", ir->userint3);
996 PI("userint4", ir->userint4);
997 PR("userreal1", ir->userreal1);
998 PR("userreal2", ir->userreal2);
999 PR("userreal3", ir->userreal3);
1000 PR("userreal4", ir->userreal4);
1004 gmx::TextWriter writer(fp);
1005 writer.wrapperSettings().setIndent(indent);
1006 gmx::dumpKeyValueTree(&writer, *ir->params);
1009 pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
1016 static void cmp_grpopts(FILE* fp, const t_grpopts* opt1, const t_grpopts* opt2, real ftol, real abstol)
1019 char buf1[256], buf2[256];
1021 cmp_int(fp, "inputrec->grpopts.ngtc", -1, opt1->ngtc, opt2->ngtc);
1022 cmp_int(fp, "inputrec->grpopts.ngacc", -1, opt1->ngacc, opt2->ngacc);
1023 cmp_int(fp, "inputrec->grpopts.ngfrz", -1, opt1->ngfrz, opt2->ngfrz);
1024 cmp_int(fp, "inputrec->grpopts.ngener", -1, opt1->ngener, opt2->ngener);
1025 for (i = 0; (i < std::min(opt1->ngtc, opt2->ngtc)); i++)
1027 cmp_real(fp, "inputrec->grpopts.nrdf", i, opt1->nrdf[i], opt2->nrdf[i], ftol, abstol);
1028 cmp_real(fp, "inputrec->grpopts.ref_t", i, opt1->ref_t[i], opt2->ref_t[i], ftol, abstol);
1029 cmp_real(fp, "inputrec->grpopts.tau_t", i, opt1->tau_t[i], opt2->tau_t[i], ftol, abstol);
1030 cmp_int(fp, "inputrec->grpopts.annealing", i, opt1->annealing[i], opt2->annealing[i]);
1031 cmp_int(fp, "inputrec->grpopts.anneal_npoints", i, opt1->anneal_npoints[i],
1032 opt2->anneal_npoints[i]);
1033 if (opt1->anneal_npoints[i] == opt2->anneal_npoints[i])
1035 sprintf(buf1, "inputrec->grpopts.anneal_time[%d]", i);
1036 sprintf(buf2, "inputrec->grpopts.anneal_temp[%d]", i);
1037 for (j = 0; j < opt1->anneal_npoints[i]; j++)
1039 cmp_real(fp, buf1, j, opt1->anneal_time[i][j], opt2->anneal_time[i][j], ftol, abstol);
1040 cmp_real(fp, buf2, j, opt1->anneal_temp[i][j], opt2->anneal_temp[i][j], ftol, abstol);
1044 if (opt1->ngener == opt2->ngener)
1046 for (i = 0; i < opt1->ngener; i++)
1048 for (j = i; j < opt1->ngener; j++)
1050 sprintf(buf1, "inputrec->grpopts.egp_flags[%d]", i);
1051 cmp_int(fp, buf1, j, opt1->egp_flags[opt1->ngener * i + j],
1052 opt2->egp_flags[opt1->ngener * i + j]);
1056 for (i = 0; (i < std::min(opt1->ngacc, opt2->ngacc)); i++)
1058 cmp_rvec(fp, "inputrec->grpopts.acc", i, opt1->acc[i], opt2->acc[i], ftol, abstol);
1060 for (i = 0; (i < std::min(opt1->ngfrz, opt2->ngfrz)); i++)
1062 cmp_ivec(fp, "inputrec->grpopts.nFreeze", i, opt1->nFreeze[i], opt2->nFreeze[i]);
1066 static void cmp_pull(FILE* fp)
1069 "WARNING: Both files use COM pulling, but comparing of the pull struct is not "
1070 "implemented (yet). The pull parameters could be the same or different.\n");
1073 static void cmp_awhDimParams(FILE* fp,
1074 const gmx::AwhDimParams* dimp1,
1075 const gmx::AwhDimParams* dimp2,
1080 /* Note that we have double index here, but the compare functions only
1081 * support one index, so here we only print the dim index and not the bias.
1083 cmp_int(fp, "inputrec.awhParams->bias?->dim->coord_index", dimIndex, dimp1->coordIndex,
1085 cmp_double(fp, "inputrec->awhParams->bias?->dim->period", dimIndex, dimp1->period,
1086 dimp2->period, ftol, abstol);
1087 cmp_double(fp, "inputrec->awhParams->bias?->dim->diffusion", dimIndex, dimp1->diffusion,
1088 dimp2->diffusion, ftol, abstol);
1089 cmp_double(fp, "inputrec->awhParams->bias?->dim->origin", dimIndex, dimp1->origin,
1090 dimp2->origin, ftol, abstol);
1091 cmp_double(fp, "inputrec->awhParams->bias?->dim->end", dimIndex, dimp1->end, dimp2->end, ftol, abstol);
1092 cmp_double(fp, "inputrec->awhParams->bias?->dim->coord_value_init", dimIndex,
1093 dimp1->coordValueInit, dimp2->coordValueInit, ftol, abstol);
1094 cmp_double(fp, "inputrec->awhParams->bias?->dim->coverDiameter", dimIndex, dimp1->coverDiameter,
1095 dimp2->coverDiameter, ftol, abstol);
1098 static void cmp_awhBiasParams(FILE* fp,
1099 const gmx::AwhBiasParams* bias1,
1100 const gmx::AwhBiasParams* bias2,
1105 cmp_int(fp, "inputrec->awhParams->ndim", biasIndex, bias1->ndim, bias2->ndim);
1106 cmp_int(fp, "inputrec->awhParams->biaseTarget", biasIndex, bias1->eTarget, bias2->eTarget);
1107 cmp_double(fp, "inputrec->awhParams->biastargetBetaScaling", biasIndex,
1108 bias1->targetBetaScaling, bias2->targetBetaScaling, ftol, abstol);
1109 cmp_double(fp, "inputrec->awhParams->biastargetCutoff", biasIndex, bias1->targetCutoff,
1110 bias2->targetCutoff, ftol, abstol);
1111 cmp_int(fp, "inputrec->awhParams->biaseGrowth", biasIndex, bias1->eGrowth, bias2->eGrowth);
1112 cmp_bool(fp, "inputrec->awhParams->biasbUserData", biasIndex, bias1->bUserData != 0,
1113 bias2->bUserData != 0);
1114 cmp_double(fp, "inputrec->awhParams->biaserror_initial", biasIndex, bias1->errorInitial,
1115 bias2->errorInitial, ftol, abstol);
1116 cmp_int(fp, "inputrec->awhParams->biasShareGroup", biasIndex, bias1->shareGroup, bias2->shareGroup);
1118 for (int dim = 0; dim < std::min(bias1->ndim, bias2->ndim); dim++)
1120 cmp_awhDimParams(fp, &bias1->dimParams[dim], &bias2->dimParams[dim], dim, ftol, abstol);
1124 static void cmp_awhParams(FILE* fp, const gmx::AwhParams* awh1, const gmx::AwhParams* awh2, real ftol, real abstol)
1126 cmp_int(fp, "inputrec->awhParams->nbias", -1, awh1->numBias, awh2->numBias);
1127 cmp_int64(fp, "inputrec->awhParams->seed", awh1->seed, awh2->seed);
1128 cmp_int(fp, "inputrec->awhParams->nstout", -1, awh1->nstOut, awh2->nstOut);
1129 cmp_int(fp, "inputrec->awhParams->nstsample_coord", -1, awh1->nstSampleCoord, awh2->nstSampleCoord);
1130 cmp_int(fp, "inputrec->awhParams->nsamples_update_free_energy", -1,
1131 awh1->numSamplesUpdateFreeEnergy, awh2->numSamplesUpdateFreeEnergy);
1132 cmp_int(fp, "inputrec->awhParams->ePotential", -1, awh1->ePotential, awh2->ePotential);
1133 cmp_bool(fp, "inputrec->awhParams->shareBiasMultisim", -1, awh1->shareBiasMultisim,
1134 awh2->shareBiasMultisim);
1136 if (awh1->numBias == awh2->numBias)
1138 for (int bias = 0; bias < awh1->numBias; bias++)
1140 cmp_awhBiasParams(fp, &awh1->awhBiasParams[bias], &awh2->awhBiasParams[bias], bias, ftol, abstol);
1145 static void cmp_simtempvals(FILE* fp,
1146 const t_simtemp* simtemp1,
1147 const t_simtemp* simtemp2,
1153 cmp_int(fp, "inputrec->simtempvals->eSimTempScale", -1, simtemp1->eSimTempScale, simtemp2->eSimTempScale);
1154 cmp_real(fp, "inputrec->simtempvals->simtemp_high", -1, simtemp1->simtemp_high,
1155 simtemp2->simtemp_high, ftol, abstol);
1156 cmp_real(fp, "inputrec->simtempvals->simtemp_low", -1, simtemp1->simtemp_low,
1157 simtemp2->simtemp_low, ftol, abstol);
1158 for (i = 0; i < n_lambda; i++)
1160 cmp_real(fp, "inputrec->simtempvals->temperatures", -1, simtemp1->temperatures[i],
1161 simtemp2->temperatures[i], ftol, abstol);
1165 static void cmp_expandedvals(FILE* fp,
1166 const t_expanded* expand1,
1167 const t_expanded* expand2,
1174 cmp_bool(fp, "inputrec->fepvals->bInit_weights", -1, expand1->bInit_weights, expand2->bInit_weights);
1175 cmp_bool(fp, "inputrec->fepvals->bWLoneovert", -1, expand1->bWLoneovert, expand2->bWLoneovert);
1177 for (i = 0; i < n_lambda; i++)
1179 cmp_real(fp, "inputrec->expandedvals->init_lambda_weights", -1,
1180 expand1->init_lambda_weights[i], expand2->init_lambda_weights[i], ftol, abstol);
1183 cmp_int(fp, "inputrec->expandedvals->lambda-stats", -1, expand1->elamstats, expand2->elamstats);
1184 cmp_int(fp, "inputrec->expandedvals->lambda-mc-move", -1, expand1->elmcmove, expand2->elmcmove);
1185 cmp_int(fp, "inputrec->expandedvals->lmc-repeats", -1, expand1->lmc_repeats, expand2->lmc_repeats);
1186 cmp_int(fp, "inputrec->expandedvals->lmc-gibbsdelta", -1, expand1->gibbsdeltalam, expand2->gibbsdeltalam);
1187 cmp_int(fp, "inputrec->expandedvals->lmc-forced-nstart", -1, expand1->lmc_forced_nstart,
1188 expand2->lmc_forced_nstart);
1189 cmp_int(fp, "inputrec->expandedvals->lambda-weights-equil", -1, expand1->elmceq, expand2->elmceq);
1190 cmp_int(fp, "inputrec->expandedvals->,weight-equil-number-all-lambda", -1,
1191 expand1->equil_n_at_lam, expand2->equil_n_at_lam);
1192 cmp_int(fp, "inputrec->expandedvals->weight-equil-number-samples", -1, expand1->equil_samples,
1193 expand2->equil_samples);
1194 cmp_int(fp, "inputrec->expandedvals->weight-equil-number-steps", -1, expand1->equil_steps,
1195 expand2->equil_steps);
1196 cmp_real(fp, "inputrec->expandedvals->weight-equil-wl-delta", -1, expand1->equil_wl_delta,
1197 expand2->equil_wl_delta, ftol, abstol);
1198 cmp_real(fp, "inputrec->expandedvals->weight-equil-count-ratio", -1, expand1->equil_ratio,
1199 expand2->equil_ratio, ftol, abstol);
1200 cmp_bool(fp, "inputrec->expandedvals->symmetrized-transition-matrix", -1,
1201 expand1->bSymmetrizedTMatrix, expand2->bSymmetrizedTMatrix);
1202 cmp_int(fp, "inputrec->expandedvals->nstTij", -1, expand1->nstTij, expand2->nstTij);
1203 cmp_int(fp, "inputrec->expandedvals->mininum-var-min", -1, expand1->minvarmin,
1204 expand2->minvarmin); /*default is reasonable */
1205 cmp_int(fp, "inputrec->expandedvals->weight-c-range", -1, expand1->c_range, expand2->c_range); /* default is just C=0 */
1206 cmp_real(fp, "inputrec->expandedvals->wl-scale", -1, expand1->wl_scale, expand2->wl_scale, ftol, abstol);
1207 cmp_real(fp, "inputrec->expandedvals->init-wl-delta", -1, expand1->init_wl_delta,
1208 expand2->init_wl_delta, ftol, abstol);
1209 cmp_real(fp, "inputrec->expandedvals->wl-ratio", -1, expand1->wl_ratio, expand2->wl_ratio, ftol, abstol);
1210 cmp_int(fp, "inputrec->expandedvals->nstexpanded", -1, expand1->nstexpanded, expand2->nstexpanded);
1211 cmp_int(fp, "inputrec->expandedvals->lmc-seed", -1, expand1->lmc_seed, expand2->lmc_seed);
1212 cmp_real(fp, "inputrec->expandedvals->mc-temperature", -1, expand1->mc_temp, expand2->mc_temp,
1216 static void cmp_fepvals(FILE* fp, const t_lambda* fep1, const t_lambda* fep2, real ftol, real abstol)
1219 cmp_int(fp, "inputrec->nstdhdl", -1, fep1->nstdhdl, fep2->nstdhdl);
1220 cmp_double(fp, "inputrec->fepvals->init_fep_state", -1, fep1->init_fep_state,
1221 fep2->init_fep_state, ftol, abstol);
1222 cmp_double(fp, "inputrec->fepvals->delta_lambda", -1, fep1->delta_lambda, fep2->delta_lambda,
1224 cmp_int(fp, "inputrec->fepvals->n_lambda", -1, fep1->n_lambda, fep2->n_lambda);
1225 for (i = 0; i < efptNR; i++)
1227 for (j = 0; j < std::min(fep1->n_lambda, fep2->n_lambda); j++)
1229 cmp_double(fp, "inputrec->fepvals->all_lambda", -1, fep1->all_lambda[i][j],
1230 fep2->all_lambda[i][j], ftol, abstol);
1233 cmp_int(fp, "inputrec->fepvals->lambda_neighbors", 1, fep1->lambda_neighbors, fep2->lambda_neighbors);
1234 cmp_real(fp, "inputrec->fepvals->sc_alpha", -1, fep1->sc_alpha, fep2->sc_alpha, ftol, abstol);
1235 cmp_int(fp, "inputrec->fepvals->sc_power", -1, fep1->sc_power, fep2->sc_power);
1236 cmp_real(fp, "inputrec->fepvals->sc_r_power", -1, fep1->sc_r_power, fep2->sc_r_power, ftol, abstol);
1237 cmp_real(fp, "inputrec->fepvals->sc_sigma", -1, fep1->sc_sigma, fep2->sc_sigma, ftol, abstol);
1238 cmp_int(fp, "inputrec->fepvals->edHdLPrintEnergy", -1, fep1->edHdLPrintEnergy, fep1->edHdLPrintEnergy);
1239 cmp_bool(fp, "inputrec->fepvals->bScCoul", -1, fep1->bScCoul, fep1->bScCoul);
1240 cmp_int(fp, "inputrec->separate_dhdl_file", -1, fep1->separate_dhdl_file, fep2->separate_dhdl_file);
1241 cmp_int(fp, "inputrec->dhdl_derivatives", -1, fep1->dhdl_derivatives, fep2->dhdl_derivatives);
1242 cmp_int(fp, "inputrec->dh_hist_size", -1, fep1->dh_hist_size, fep2->dh_hist_size);
1243 cmp_double(fp, "inputrec->dh_hist_spacing", -1, fep1->dh_hist_spacing, fep2->dh_hist_spacing,
1247 void cmp_inputrec(FILE* fp, const t_inputrec* ir1, const t_inputrec* ir2, real ftol, real abstol)
1249 fprintf(fp, "comparing inputrec\n");
1251 /* gcc 2.96 doesnt like these defines at all, but issues a huge list
1252 * of warnings. Maybe it will change in future versions, but for the
1253 * moment I've spelled them out instead. /EL 000820
1254 * #define CIB(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
1255 * #define CII(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
1256 * #define CIR(s) cmp_real(fp,"inputrec->"#s,0,ir1->##s,ir2->##s,ftol)
1258 cmp_int(fp, "inputrec->eI", -1, ir1->eI, ir2->eI);
1259 cmp_int64(fp, "inputrec->nsteps", ir1->nsteps, ir2->nsteps);
1260 cmp_int64(fp, "inputrec->init_step", ir1->init_step, ir2->init_step);
1261 cmp_int(fp, "inputrec->simulation_part", -1, ir1->simulation_part, ir2->simulation_part);
1262 cmp_int(fp, "inputrec->pbcType", -1, static_cast<int>(ir1->pbcType), static_cast<int>(ir2->pbcType));
1263 cmp_bool(fp, "inputrec->bPeriodicMols", -1, ir1->bPeriodicMols, ir2->bPeriodicMols);
1264 cmp_int(fp, "inputrec->cutoff_scheme", -1, ir1->cutoff_scheme, ir2->cutoff_scheme);
1265 cmp_int(fp, "inputrec->nstlist", -1, ir1->nstlist, ir2->nstlist);
1266 cmp_int(fp, "inputrec->nstcomm", -1, ir1->nstcomm, ir2->nstcomm);
1267 cmp_int(fp, "inputrec->comm_mode", -1, ir1->comm_mode, ir2->comm_mode);
1268 cmp_int(fp, "inputrec->nstlog", -1, ir1->nstlog, ir2->nstlog);
1269 cmp_int(fp, "inputrec->nstxout", -1, ir1->nstxout, ir2->nstxout);
1270 cmp_int(fp, "inputrec->nstvout", -1, ir1->nstvout, ir2->nstvout);
1271 cmp_int(fp, "inputrec->nstfout", -1, ir1->nstfout, ir2->nstfout);
1272 cmp_int(fp, "inputrec->nstcalcenergy", -1, ir1->nstcalcenergy, ir2->nstcalcenergy);
1273 cmp_int(fp, "inputrec->nstenergy", -1, ir1->nstenergy, ir2->nstenergy);
1274 cmp_int(fp, "inputrec->nstxout_compressed", -1, ir1->nstxout_compressed, ir2->nstxout_compressed);
1275 cmp_double(fp, "inputrec->init_t", -1, ir1->init_t, ir2->init_t, ftol, abstol);
1276 cmp_double(fp, "inputrec->delta_t", -1, ir1->delta_t, ir2->delta_t, ftol, abstol);
1277 cmp_real(fp, "inputrec->x_compression_precision", -1, ir1->x_compression_precision,
1278 ir2->x_compression_precision, ftol, abstol);
1279 cmp_real(fp, "inputrec->fourierspacing", -1, ir1->fourier_spacing, ir2->fourier_spacing, ftol, abstol);
1280 cmp_int(fp, "inputrec->nkx", -1, ir1->nkx, ir2->nkx);
1281 cmp_int(fp, "inputrec->nky", -1, ir1->nky, ir2->nky);
1282 cmp_int(fp, "inputrec->nkz", -1, ir1->nkz, ir2->nkz);
1283 cmp_int(fp, "inputrec->pme_order", -1, ir1->pme_order, ir2->pme_order);
1284 cmp_real(fp, "inputrec->ewald_rtol", -1, ir1->ewald_rtol, ir2->ewald_rtol, ftol, abstol);
1285 cmp_int(fp, "inputrec->ewald_geometry", -1, ir1->ewald_geometry, ir2->ewald_geometry);
1286 cmp_real(fp, "inputrec->epsilon_surface", -1, ir1->epsilon_surface, ir2->epsilon_surface, ftol, abstol);
1287 cmp_int(fp, "inputrec->bContinuation", -1, static_cast<int>(ir1->bContinuation),
1288 static_cast<int>(ir2->bContinuation));
1289 cmp_int(fp, "inputrec->bShakeSOR", -1, static_cast<int>(ir1->bShakeSOR),
1290 static_cast<int>(ir2->bShakeSOR));
1291 cmp_int(fp, "inputrec->etc", -1, ir1->etc, ir2->etc);
1292 cmp_int(fp, "inputrec->bPrintNHChains", -1, static_cast<int>(ir1->bPrintNHChains),
1293 static_cast<int>(ir2->bPrintNHChains));
1294 cmp_int(fp, "inputrec->epc", -1, ir1->epc, ir2->epc);
1295 cmp_int(fp, "inputrec->epct", -1, ir1->epct, ir2->epct);
1296 cmp_real(fp, "inputrec->tau_p", -1, ir1->tau_p, ir2->tau_p, ftol, abstol);
1297 cmp_rvec(fp, "inputrec->ref_p(x)", -1, ir1->ref_p[XX], ir2->ref_p[XX], ftol, abstol);
1298 cmp_rvec(fp, "inputrec->ref_p(y)", -1, ir1->ref_p[YY], ir2->ref_p[YY], ftol, abstol);
1299 cmp_rvec(fp, "inputrec->ref_p(z)", -1, ir1->ref_p[ZZ], ir2->ref_p[ZZ], ftol, abstol);
1300 cmp_rvec(fp, "inputrec->compress(x)", -1, ir1->compress[XX], ir2->compress[XX], ftol, abstol);
1301 cmp_rvec(fp, "inputrec->compress(y)", -1, ir1->compress[YY], ir2->compress[YY], ftol, abstol);
1302 cmp_rvec(fp, "inputrec->compress(z)", -1, ir1->compress[ZZ], ir2->compress[ZZ], ftol, abstol);
1303 cmp_int(fp, "refcoord_scaling", -1, ir1->refcoord_scaling, ir2->refcoord_scaling);
1304 cmp_rvec(fp, "inputrec->posres_com", -1, ir1->posres_com, ir2->posres_com, ftol, abstol);
1305 cmp_rvec(fp, "inputrec->posres_comB", -1, ir1->posres_comB, ir2->posres_comB, ftol, abstol);
1306 cmp_real(fp, "inputrec->verletbuf_tol", -1, ir1->verletbuf_tol, ir2->verletbuf_tol, ftol, abstol);
1307 cmp_real(fp, "inputrec->rlist", -1, ir1->rlist, ir2->rlist, ftol, abstol);
1308 cmp_real(fp, "inputrec->rtpi", -1, ir1->rtpi, ir2->rtpi, ftol, abstol);
1309 cmp_int(fp, "inputrec->coulombtype", -1, ir1->coulombtype, ir2->coulombtype);
1310 cmp_int(fp, "inputrec->coulomb_modifier", -1, ir1->coulomb_modifier, ir2->coulomb_modifier);
1311 cmp_real(fp, "inputrec->rcoulomb_switch", -1, ir1->rcoulomb_switch, ir2->rcoulomb_switch, ftol, abstol);
1312 cmp_real(fp, "inputrec->rcoulomb", -1, ir1->rcoulomb, ir2->rcoulomb, ftol, abstol);
1313 cmp_int(fp, "inputrec->vdwtype", -1, ir1->vdwtype, ir2->vdwtype);
1314 cmp_int(fp, "inputrec->vdw_modifier", -1, ir1->vdw_modifier, ir2->vdw_modifier);
1315 cmp_real(fp, "inputrec->rvdw_switch", -1, ir1->rvdw_switch, ir2->rvdw_switch, ftol, abstol);
1316 cmp_real(fp, "inputrec->rvdw", -1, ir1->rvdw, ir2->rvdw, ftol, abstol);
1317 cmp_real(fp, "inputrec->epsilon_r", -1, ir1->epsilon_r, ir2->epsilon_r, ftol, abstol);
1318 cmp_real(fp, "inputrec->epsilon_rf", -1, ir1->epsilon_rf, ir2->epsilon_rf, ftol, abstol);
1319 cmp_real(fp, "inputrec->tabext", -1, ir1->tabext, ir2->tabext, ftol, abstol);
1321 cmp_int(fp, "inputrec->eDispCorr", -1, ir1->eDispCorr, ir2->eDispCorr);
1322 cmp_real(fp, "inputrec->shake_tol", -1, ir1->shake_tol, ir2->shake_tol, ftol, abstol);
1323 cmp_int(fp, "inputrec->efep", -1, ir1->efep, ir2->efep);
1324 cmp_fepvals(fp, ir1->fepvals, ir2->fepvals, ftol, abstol);
1325 cmp_int(fp, "inputrec->bSimTemp", -1, static_cast<int>(ir1->bSimTemp), static_cast<int>(ir2->bSimTemp));
1326 if ((ir1->bSimTemp == ir2->bSimTemp) && (ir1->bSimTemp))
1328 cmp_simtempvals(fp, ir1->simtempvals, ir2->simtempvals,
1329 std::min(ir1->fepvals->n_lambda, ir2->fepvals->n_lambda), ftol, abstol);
1331 cmp_int(fp, "inputrec->bExpanded", -1, static_cast<int>(ir1->bExpanded),
1332 static_cast<int>(ir2->bExpanded));
1333 if ((ir1->bExpanded == ir2->bExpanded) && (ir1->bExpanded))
1335 cmp_expandedvals(fp, ir1->expandedvals, ir2->expandedvals,
1336 std::min(ir1->fepvals->n_lambda, ir2->fepvals->n_lambda), ftol, abstol);
1338 cmp_int(fp, "inputrec->nwall", -1, ir1->nwall, ir2->nwall);
1339 cmp_int(fp, "inputrec->wall_type", -1, ir1->wall_type, ir2->wall_type);
1340 cmp_int(fp, "inputrec->wall_atomtype[0]", -1, ir1->wall_atomtype[0], ir2->wall_atomtype[0]);
1341 cmp_int(fp, "inputrec->wall_atomtype[1]", -1, ir1->wall_atomtype[1], ir2->wall_atomtype[1]);
1342 cmp_real(fp, "inputrec->wall_density[0]", -1, ir1->wall_density[0], ir2->wall_density[0], ftol, abstol);
1343 cmp_real(fp, "inputrec->wall_density[1]", -1, ir1->wall_density[1], ir2->wall_density[1], ftol, abstol);
1344 cmp_real(fp, "inputrec->wall_ewald_zfac", -1, ir1->wall_ewald_zfac, ir2->wall_ewald_zfac, ftol, abstol);
1346 cmp_bool(fp, "inputrec->bPull", -1, ir1->bPull, ir2->bPull);
1347 if (ir1->bPull && ir2->bPull)
1352 cmp_bool(fp, "inputrec->bDoAwh", -1, ir1->bDoAwh, ir2->bDoAwh);
1353 if (ir1->bDoAwh && ir2->bDoAwh)
1355 cmp_awhParams(fp, ir1->awhParams, ir2->awhParams, ftol, abstol);
1358 cmp_int(fp, "inputrec->eDisre", -1, ir1->eDisre, ir2->eDisre);
1359 cmp_real(fp, "inputrec->dr_fc", -1, ir1->dr_fc, ir2->dr_fc, ftol, abstol);
1360 cmp_int(fp, "inputrec->eDisreWeighting", -1, ir1->eDisreWeighting, ir2->eDisreWeighting);
1361 cmp_int(fp, "inputrec->bDisreMixed", -1, static_cast<int>(ir1->bDisreMixed),
1362 static_cast<int>(ir2->bDisreMixed));
1363 cmp_int(fp, "inputrec->nstdisreout", -1, ir1->nstdisreout, ir2->nstdisreout);
1364 cmp_real(fp, "inputrec->dr_tau", -1, ir1->dr_tau, ir2->dr_tau, ftol, abstol);
1365 cmp_real(fp, "inputrec->orires_fc", -1, ir1->orires_fc, ir2->orires_fc, ftol, abstol);
1366 cmp_real(fp, "inputrec->orires_tau", -1, ir1->orires_tau, ir2->orires_tau, ftol, abstol);
1367 cmp_int(fp, "inputrec->nstorireout", -1, ir1->nstorireout, ir2->nstorireout);
1368 cmp_real(fp, "inputrec->em_stepsize", -1, ir1->em_stepsize, ir2->em_stepsize, ftol, abstol);
1369 cmp_real(fp, "inputrec->em_tol", -1, ir1->em_tol, ir2->em_tol, ftol, abstol);
1370 cmp_int(fp, "inputrec->niter", -1, ir1->niter, ir2->niter);
1371 cmp_real(fp, "inputrec->fc_stepsize", -1, ir1->fc_stepsize, ir2->fc_stepsize, ftol, abstol);
1372 cmp_int(fp, "inputrec->nstcgsteep", -1, ir1->nstcgsteep, ir2->nstcgsteep);
1373 cmp_int(fp, "inputrec->nbfgscorr", 0, ir1->nbfgscorr, ir2->nbfgscorr);
1374 cmp_int(fp, "inputrec->eConstrAlg", -1, ir1->eConstrAlg, ir2->eConstrAlg);
1375 cmp_int(fp, "inputrec->nProjOrder", -1, ir1->nProjOrder, ir2->nProjOrder);
1376 cmp_real(fp, "inputrec->LincsWarnAngle", -1, ir1->LincsWarnAngle, ir2->LincsWarnAngle, ftol, abstol);
1377 cmp_int(fp, "inputrec->nLincsIter", -1, ir1->nLincsIter, ir2->nLincsIter);
1378 cmp_real(fp, "inputrec->bd_fric", -1, ir1->bd_fric, ir2->bd_fric, ftol, abstol);
1379 cmp_int64(fp, "inputrec->ld_seed", ir1->ld_seed, ir2->ld_seed);
1380 cmp_real(fp, "inputrec->cos_accel", -1, ir1->cos_accel, ir2->cos_accel, ftol, abstol);
1381 cmp_rvec(fp, "inputrec->deform(a)", -1, ir1->deform[XX], ir2->deform[XX], ftol, abstol);
1382 cmp_rvec(fp, "inputrec->deform(b)", -1, ir1->deform[YY], ir2->deform[YY], ftol, abstol);
1383 cmp_rvec(fp, "inputrec->deform(c)", -1, ir1->deform[ZZ], ir2->deform[ZZ], ftol, abstol);
1386 cmp_int(fp, "inputrec->userint1", -1, ir1->userint1, ir2->userint1);
1387 cmp_int(fp, "inputrec->userint2", -1, ir1->userint2, ir2->userint2);
1388 cmp_int(fp, "inputrec->userint3", -1, ir1->userint3, ir2->userint3);
1389 cmp_int(fp, "inputrec->userint4", -1, ir1->userint4, ir2->userint4);
1390 cmp_real(fp, "inputrec->userreal1", -1, ir1->userreal1, ir2->userreal1, ftol, abstol);
1391 cmp_real(fp, "inputrec->userreal2", -1, ir1->userreal2, ir2->userreal2, ftol, abstol);
1392 cmp_real(fp, "inputrec->userreal3", -1, ir1->userreal3, ir2->userreal3, ftol, abstol);
1393 cmp_real(fp, "inputrec->userreal4", -1, ir1->userreal4, ir2->userreal4, ftol, abstol);
1394 cmp_grpopts(fp, &(ir1->opts), &(ir2->opts), ftol, abstol);
1395 gmx::TextWriter writer(fp);
1396 gmx::compareKeyValueTrees(&writer, *ir1->params, *ir2->params, ftol, abstol);
1399 void comp_pull_AB(FILE* fp, pull_params_t* pull, real ftol, real abstol)
1403 for (i = 0; i < pull->ncoord; i++)
1405 fprintf(fp, "comparing pull coord %d\n", i);
1406 cmp_real(fp, "pull-coord->k", -1, pull->coord[i].k, pull->coord[i].kB, ftol, abstol);
1410 gmx_bool inputrecDeform(const t_inputrec* ir)
1412 return (ir->deform[XX][XX] != 0 || ir->deform[YY][YY] != 0 || ir->deform[ZZ][ZZ] != 0
1413 || ir->deform[YY][XX] != 0 || ir->deform[ZZ][XX] != 0 || ir->deform[ZZ][YY] != 0);
1416 gmx_bool inputrecDynamicBox(const t_inputrec* ir)
1418 return (ir->epc != epcNO || ir->eI == eiTPI || inputrecDeform(ir));
1421 gmx_bool inputrecPreserveShape(const t_inputrec* ir)
1423 return (ir->epc != epcNO && ir->deform[XX][XX] == 0
1424 && (ir->epct == epctISOTROPIC || ir->epct == epctSEMIISOTROPIC));
1427 gmx_bool inputrecNeedMutot(const t_inputrec* ir)
1429 return ((ir->coulombtype == eelEWALD || EEL_PME(ir->coulombtype))
1430 && (ir->ewald_geometry == eewg3DC || ir->epsilon_surface != 0));
1433 gmx_bool inputrecExclForces(const t_inputrec* ir)
1435 return (EEL_FULL(ir->coulombtype) || (EEL_RF(ir->coulombtype)));
1438 gmx_bool inputrecNptTrotter(const t_inputrec* ir)
1440 return (((ir->eI == eiVV) || (ir->eI == eiVVAK)) && (ir->epc == epcMTTK) && (ir->etc == etcNOSEHOOVER));
1443 gmx_bool inputrecNvtTrotter(const t_inputrec* ir)
1445 return (((ir->eI == eiVV) || (ir->eI == eiVVAK)) && (ir->epc != epcMTTK) && (ir->etc == etcNOSEHOOVER));
1448 gmx_bool inputrecNphTrotter(const t_inputrec* ir)
1450 return (((ir->eI == eiVV) || (ir->eI == eiVVAK)) && (ir->epc == epcMTTK) && (ir->etc != etcNOSEHOOVER));
1453 bool inputrecPbcXY2Walls(const t_inputrec* ir)
1455 return (ir->pbcType == PbcType::XY && ir->nwall == 2);
1458 bool integratorHasConservedEnergyQuantity(const t_inputrec* ir)
1461 { // NOLINT bugprone-branch-clone
1462 // Energy minimization or stochastic integrator: no conservation
1465 else if (ir->etc == etcNO && ir->epc == epcNO)
1467 // The total energy is conserved, no additional conserved quanitity
1472 // Shear stress with Parrinello-Rahman is not supported (tedious)
1474 ((ir->epc == epcPARRINELLORAHMAN || ir->epc == epcMTTK)
1475 && (ir->ref_p[YY][XX] != 0 || ir->ref_p[ZZ][XX] != 0 || ir->ref_p[ZZ][YY] != 0));
1477 return !ETC_ANDERSEN(ir->etc) && !shearWithPR;
1481 bool integratorHasReferenceTemperature(const t_inputrec* ir)
1483 return ((ir->etc != etcNO) || EI_SD(ir->eI) || (ir->eI == eiBD) || EI_TPI(ir->eI));
1486 int inputrec2nboundeddim(const t_inputrec* ir)
1488 if (inputrecPbcXY2Walls(ir))
1494 return numPbcDimensions(ir->pbcType);
1498 int ndof_com(const t_inputrec* ir)
1502 switch (ir->pbcType)
1505 case PbcType::No: n = 3; break;
1506 case PbcType::XY: n = (ir->nwall == 0 ? 3 : 2); break;
1507 case PbcType::Screw: n = 1; break;
1508 default: gmx_incons("Unknown pbc in calc_nrdf");
1514 real maxReferenceTemperature(const t_inputrec& ir)
1516 if (EI_ENERGY_MINIMIZATION(ir.eI) || ir.eI == eiNM)
1521 if (EI_MD(ir.eI) && ir.etc == etcNO)
1526 /* SD and BD also use ref_t and tau_t for setting the reference temperature.
1527 * TPI can be treated as MD, since it needs an ensemble temperature.
1530 real maxTemperature = 0;
1531 for (int i = 0; i < ir.opts.ngtc; i++)
1533 if (ir.opts.tau_t[i] >= 0)
1535 maxTemperature = std::max(maxTemperature, ir.opts.ref_t[i]);
1539 return maxTemperature;
1542 bool haveEwaldSurfaceContribution(const t_inputrec& ir)
1544 return EEL_PME_EWALD(ir.coulombtype) && (ir.ewald_geometry == eewg3DC || ir.epsilon_surface != 0);
1547 bool haveFreeEnergyType(const t_inputrec& ir, const int fepType)
1549 for (int i = 0; i < ir.fepvals->n_lambda; i++)
1551 if (ir.fepvals->all_lambda[fepType][i] > 0)