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,2021, 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.
40 #include "gromacs/utility/enumerationhelpers.h"
51 #include "gromacs/applied_forces/awh/read_params.h"
52 #include "gromacs/math/veccompare.h"
53 #include "gromacs/math/vecdump.h"
54 #include "gromacs/mdlib/vcm.h"
55 #include "gromacs/mdtypes/awh_params.h"
56 #include "gromacs/mdtypes/md_enums.h"
57 #include "gromacs/mdtypes/multipletimestepping.h"
58 #include "gromacs/mdtypes/pull_params.h"
59 #include "gromacs/pbcutil/pbc.h"
60 #include "gromacs/utility/compare.h"
61 #include "gromacs/utility/cstringutil.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/keyvaluetree.h"
64 #include "gromacs/utility/smalloc.h"
65 #include "gromacs/utility/snprintf.h"
66 #include "gromacs/utility/strconvert.h"
67 #include "gromacs/utility/stringutil.h"
68 #include "gromacs/utility/textwriter.h"
69 #include "gromacs/utility/txtdump.h"
71 //! Macro to select a bool name
72 #define EBOOL(e) gmx::boolToString(e)
74 /* Default values for nstcalcenergy, used when the are no other restrictions. */
75 constexpr int c_defaultNstCalcEnergy = 10;
77 /* The minimum number of integration steps required for reasonably accurate
78 * integration of first and second order coupling algorithms.
80 const int nstmin_berendsen_tcouple = 5;
81 const int nstmin_berendsen_pcouple = 10;
82 const int nstmin_harmonic = 20;
84 /* Default values for T- and P- coupling intervals, used when the are no other
87 constexpr int c_defaultNstTCouple = 10;
88 constexpr int c_defaultNstPCouple = 10;
90 t_inputrec::t_inputrec()
92 // TODO When this memset is removed, remove the suppression of
93 // gcc -Wno-class-memaccess in a CMakeLists.txt file.
94 std::memset(this, 0, sizeof(*this)); // NOLINT(bugprone-undefined-memory-manipulation)
95 fepvals = std::make_unique<t_lambda>();
96 expandedvals = std::make_unique<t_expanded>();
97 simtempvals = std::make_unique<t_simtemp>();
100 t_inputrec::~t_inputrec()
105 int ir_optimal_nstcalcenergy(const t_inputrec* ir)
115 nst = c_defaultNstCalcEnergy;
120 nst = std::lcm(nst, ir->mtsLevels.back().stepFactor);
126 int tcouple_min_integration_steps(TemperatureCoupling etc)
132 case TemperatureCoupling::No: n = 0; break;
133 case TemperatureCoupling::Berendsen:
134 case TemperatureCoupling::Yes: n = nstmin_berendsen_tcouple; break;
135 case TemperatureCoupling::VRescale:
136 /* V-rescale supports instantaneous rescaling */
139 case TemperatureCoupling::NoseHoover: n = nstmin_harmonic; break;
140 case TemperatureCoupling::Andersen:
141 case TemperatureCoupling::AndersenMassive: n = 1; break;
142 default: gmx_incons("Unknown etc value");
148 int ir_optimal_nsttcouple(const t_inputrec* ir)
150 int nmin, nwanted, n;
154 nmin = tcouple_min_integration_steps(ir->etc);
156 nwanted = c_defaultNstTCouple;
159 if (ir->etc != TemperatureCoupling::No)
161 for (g = 0; g < ir->opts.ngtc; g++)
163 if (ir->opts.tau_t[g] > 0)
165 tau_min = std::min(tau_min, ir->opts.tau_t[g]);
170 if (nmin == 0 || ir->delta_t * nwanted <= tau_min)
176 n = static_cast<int>(tau_min / (ir->delta_t * nmin) + 0.001);
181 while (nwanted % n != 0)
190 int pcouple_min_integration_steps(PressureCoupling epc)
196 case PressureCoupling::No: n = 0; break;
197 case PressureCoupling::Berendsen:
198 case PressureCoupling::CRescale:
199 case PressureCoupling::Isotropic: n = nstmin_berendsen_pcouple; break;
200 case PressureCoupling::ParrinelloRahman:
201 case PressureCoupling::Mttk: n = nstmin_harmonic; break;
202 default: gmx_incons("Unknown epc value");
208 int ir_optimal_nstpcouple(const t_inputrec* ir)
210 const int minIntegrationSteps = pcouple_min_integration_steps(ir->epc);
212 const int nwanted = c_defaultNstPCouple;
214 // With multiple time-stepping we can only compute the pressure at slowest steps
215 const int minNstPCouple = (ir->useMts ? ir->mtsLevels.back().stepFactor : 1);
218 if (minIntegrationSteps == 0 || ir->delta_t * nwanted <= ir->tau_p)
224 n = static_cast<int>(ir->tau_p / (ir->delta_t * minIntegrationSteps) + 0.001);
225 if (n < minNstPCouple)
229 // Without MTS we try to make nstpcouple a "nice" number
232 while (nwanted % n != 0)
239 // With MTS, nstpcouple should be a multiple of the slowest MTS interval
242 n = n - (n % minNstPCouple);
248 gmx_bool ir_coulomb_switched(const t_inputrec* ir)
250 return (ir->coulombtype == CoulombInteractionType::Switch
251 || ir->coulombtype == CoulombInteractionType::Shift
252 || ir->coulombtype == CoulombInteractionType::PmeSwitch
253 || ir->coulombtype == CoulombInteractionType::PmeUserSwitch
254 || ir->coulomb_modifier == InteractionModifiers::PotSwitch
255 || ir->coulomb_modifier == InteractionModifiers::ForceSwitch);
258 gmx_bool ir_coulomb_is_zero_at_cutoff(const t_inputrec* ir)
260 return (ir->cutoff_scheme == CutoffScheme::Verlet || ir_coulomb_switched(ir)
261 || ir->coulomb_modifier != InteractionModifiers::None
262 || ir->coulombtype == CoulombInteractionType::RFZero);
265 gmx_bool ir_coulomb_might_be_zero_at_cutoff(const t_inputrec* ir)
267 return (ir_coulomb_is_zero_at_cutoff(ir) || ir->coulombtype == CoulombInteractionType::User
268 || ir->coulombtype == CoulombInteractionType::PmeUser);
271 gmx_bool ir_vdw_switched(const t_inputrec* ir)
273 return (ir->vdwtype == VanDerWaalsType::Switch || ir->vdwtype == VanDerWaalsType::Shift
274 || ir->vdw_modifier == InteractionModifiers::PotSwitch
275 || ir->vdw_modifier == InteractionModifiers::ForceSwitch);
278 gmx_bool ir_vdw_is_zero_at_cutoff(const t_inputrec* ir)
280 return (ir->cutoff_scheme == CutoffScheme::Verlet || ir_vdw_switched(ir)
281 || ir->vdw_modifier != InteractionModifiers::None);
284 gmx_bool ir_vdw_might_be_zero_at_cutoff(const t_inputrec* ir)
286 return (ir_vdw_is_zero_at_cutoff(ir) || ir->vdwtype == VanDerWaalsType::User);
289 static void done_t_rot(t_rot* rot)
295 if (rot->grp != nullptr)
297 for (int i = 0; i < rot->ngrp; i++)
299 sfree(rot->grp[i].ind);
300 sfree(rot->grp[i].x_ref);
307 static void done_t_swapCoords(t_swapcoords* swapCoords)
309 if (swapCoords == nullptr)
313 for (int i = 0; i < swapCoords->ngrp; i++)
315 sfree(swapCoords->grp[i].ind);
316 sfree(swapCoords->grp[i].molname);
318 sfree(swapCoords->grp);
322 void done_inputrec(t_inputrec* ir)
324 sfree(ir->opts.nrdf);
325 sfree(ir->opts.ref_t);
326 for (int i = 0; i < ir->opts.ngtc; i++)
328 sfree(ir->opts.anneal_time[i]);
329 sfree(ir->opts.anneal_temp[i]);
331 sfree(ir->opts.annealing);
332 sfree(ir->opts.anneal_npoints);
333 sfree(ir->opts.anneal_time);
334 sfree(ir->opts.anneal_temp);
335 sfree(ir->opts.tau_t);
336 sfree(ir->opts.nFreeze);
337 sfree(ir->opts.egp_flags);
339 done_t_swapCoords(ir->swap);
344 static void pr_grp_opts(FILE* out, int indent, const char* title, const t_grpopts* opts, gmx_bool bMDPformat)
350 fprintf(out, "%s:\n", title);
353 pr_indent(out, indent);
354 fprintf(out, "nrdf%s", bMDPformat ? " = " : ":");
355 for (i = 0; (i < opts->ngtc); i++)
357 fprintf(out, " %10g", opts->nrdf[i]);
361 pr_indent(out, indent);
362 fprintf(out, "ref-t%s", bMDPformat ? " = " : ":");
363 for (i = 0; (i < opts->ngtc); i++)
365 fprintf(out, " %10g", opts->ref_t[i]);
369 pr_indent(out, indent);
370 fprintf(out, "tau-t%s", bMDPformat ? " = " : ":");
371 for (i = 0; (i < opts->ngtc); i++)
373 fprintf(out, " %10g", opts->tau_t[i]);
377 /* Pretty-print the simulated annealing info */
378 fprintf(out, "annealing%s", bMDPformat ? " = " : ":");
379 for (i = 0; (i < opts->ngtc); i++)
381 fprintf(out, " %10s", enumValueToString(opts->annealing[i]));
385 fprintf(out, "annealing-npoints%s", bMDPformat ? " = " : ":");
386 for (i = 0; (i < opts->ngtc); i++)
388 fprintf(out, " %10d", opts->anneal_npoints[i]);
392 for (i = 0; (i < opts->ngtc); i++)
394 if (opts->anneal_npoints[i] > 0)
396 fprintf(out, "annealing-time [%d]:\t", i);
397 for (j = 0; (j < opts->anneal_npoints[i]); j++)
399 fprintf(out, " %10.1f", opts->anneal_time[i][j]);
402 fprintf(out, "annealing-temp [%d]:\t", i);
403 for (j = 0; (j < opts->anneal_npoints[i]); j++)
405 fprintf(out, " %10.1f", opts->anneal_temp[i][j]);
411 pr_indent(out, indent);
412 fprintf(out, "nfreeze:");
413 for (i = 0; (i < opts->ngfrz); i++)
415 for (m = 0; (m < DIM); m++)
417 fprintf(out, " %10s", opts->nFreeze[i][m] ? "Y" : "N");
423 for (i = 0; (i < opts->ngener); i++)
425 pr_indent(out, indent);
426 fprintf(out, "energygrp-flags[%3d]:", i);
427 for (m = 0; (m < opts->ngener); m++)
429 fprintf(out, " %d", opts->egp_flags[opts->ngener * i + m]);
437 static void pr_matrix(FILE* fp, int indent, const char* title, const rvec* m, gmx_bool bMDPformat)
442 "%-10s = %g %g %g %g %g %g\n",
453 pr_rvecs(fp, indent, title, m, DIM);
457 #define PS(t, s) pr_str(fp, indent, t, s)
458 #define PI(t, s) pr_int(fp, indent, t, s)
459 #define PSTEP(t, s) pr_int64(fp, indent, t, s)
460 #define PR(t, s) pr_real(fp, indent, t, s)
461 #define PD(t, s) pr_double(fp, indent, t, s)
463 static void pr_pull_group(FILE* fp, int indent, int g, const t_pull_group* pgrp)
465 pr_indent(fp, indent);
466 fprintf(fp, "pull-group %d:\n", g);
468 pr_ivec_block(fp, indent, "atom", pgrp->ind.data(), pgrp->ind.size(), TRUE);
469 pr_rvec(fp, indent, "weight", pgrp->weight.data(), pgrp->weight.size(), TRUE);
470 PI("pbcatom", pgrp->pbcatom);
473 static void pr_pull_coord(FILE* fp, int indent, int c, const t_pull_coord* pcrd)
477 pr_indent(fp, indent);
478 fprintf(fp, "pull-coord %d:\n", c);
479 PS("type", enumValueToString(pcrd->eType));
480 if (pcrd->eType == PullingAlgorithm::External)
482 PS("potential-provider", pcrd->externalPotentialProvider.c_str());
484 PS("geometry", enumValueToString(pcrd->eGeom));
485 for (g = 0; g < pcrd->ngroup; g++)
489 sprintf(buf, "group[%d]", g);
490 PI(buf, pcrd->group[g]);
492 pr_ivec(fp, indent, "dim", pcrd->dim, DIM, TRUE);
493 pr_rvec(fp, indent, "origin", pcrd->origin, DIM, TRUE);
494 pr_rvec(fp, indent, "vec", pcrd->vec, DIM, TRUE);
495 PS("start", EBOOL(pcrd->bStart));
496 PR("init", pcrd->init);
497 PR("rate", pcrd->rate);
502 static void pr_simtempvals(FILE* fp, int indent, const t_simtemp* simtemp, int n_lambda)
504 PS("simulated-tempering-scaling", enumValueToString(simtemp->eSimTempScale));
505 PR("sim-temp-low", simtemp->simtemp_low);
506 PR("sim-temp-high", simtemp->simtemp_high);
507 pr_rvec(fp, indent, "simulated tempering temperatures", simtemp->temperatures.data(), n_lambda, TRUE);
510 static void pr_expandedvals(FILE* fp, int indent, const t_expanded* expand, int n_lambda)
513 PI("nstexpanded", expand->nstexpanded);
514 PS("lmc-stats", enumValueToString(expand->elamstats));
515 PS("lmc-move", enumValueToString(expand->elmcmove));
516 PS("lmc-weights-equil", enumValueToString(expand->elmceq));
517 if (expand->elmceq == LambdaWeightWillReachEquilibrium::NumAtLambda)
519 PI("weight-equil-number-all-lambda", expand->equil_n_at_lam);
521 if (expand->elmceq == LambdaWeightWillReachEquilibrium::Samples)
523 PI("weight-equil-number-samples", expand->equil_samples);
525 if (expand->elmceq == LambdaWeightWillReachEquilibrium::Steps)
527 PI("weight-equil-number-steps", expand->equil_steps);
529 if (expand->elmceq == LambdaWeightWillReachEquilibrium::WLDelta)
531 PR("weight-equil-wl-delta", expand->equil_wl_delta);
533 if (expand->elmceq == LambdaWeightWillReachEquilibrium::Ratio)
535 PR("weight-equil-count-ratio", expand->equil_ratio);
537 PI("lmc-seed", expand->lmc_seed);
538 PR("mc-temperature", expand->mc_temp);
539 PI("lmc-repeats", expand->lmc_repeats);
540 PI("lmc-gibbsdelta", expand->gibbsdeltalam);
541 PI("lmc-forced-nstart", expand->lmc_forced_nstart);
542 PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
543 PI("nst-transition-matrix", expand->nstTij);
544 PI("mininum-var-min", expand->minvarmin); /*default is reasonable */
545 PI("weight-c-range", expand->c_range); /* default is just C=0 */
546 PR("wl-scale", expand->wl_scale);
547 PR("wl-ratio", expand->wl_ratio);
548 PR("init-wl-delta", expand->init_wl_delta);
549 PS("wl-oneovert", EBOOL(expand->bWLoneovert));
551 pr_indent(fp, indent);
552 pr_rvec(fp, indent, "init-lambda-weights", expand->init_lambda_weights.data(), n_lambda, TRUE);
553 PS("init-weights", EBOOL(expand->bInit_weights));
556 static void pr_fepvals(FILE* fp, int indent, const t_lambda* fep, gmx_bool bMDPformat)
560 PR("init-lambda", fep->init_lambda);
561 PI("init-lambda-state", fep->init_fep_state);
562 PR("delta-lambda", fep->delta_lambda);
563 PI("nstdhdl", fep->nstdhdl);
567 PI("n-lambdas", fep->n_lambda);
569 if (fep->n_lambda > 0)
571 pr_indent(fp, indent);
572 fprintf(fp, "separate-dvdl%s\n", bMDPformat ? " = " : ":");
573 for (auto i : gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, bool>::keys())
575 fprintf(fp, "%18s = ", enumValueToString(i));
576 if (fep->separate_dvdl[i])
578 fprintf(fp, " TRUE");
582 fprintf(fp, " FALSE");
586 fprintf(fp, "all-lambdas%s\n", bMDPformat ? " = " : ":");
587 for (auto key : gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, bool>::keys())
589 fprintf(fp, "%18s = ", enumValueToString(key));
590 int i = static_cast<int>(key);
591 for (j = 0; j < fep->n_lambda; j++)
593 fprintf(fp, " %10g", fep->all_lambda[i][j]);
598 PI("calc-lambda-neighbors", fep->lambda_neighbors);
599 PS("dhdl-print-energy", enumValueToString(fep->edHdLPrintEnergy));
600 PR("sc-alpha", fep->sc_alpha);
601 PI("sc-power", fep->sc_power);
602 PR("sc-r-power", fep->sc_r_power);
603 PR("sc-sigma", fep->sc_sigma);
604 PR("sc-sigma-min", fep->sc_sigma_min);
605 PS("sc-coul", EBOOL(fep->bScCoul));
606 PI("dh-hist-size", fep->dh_hist_size);
607 PD("dh-hist-spacing", fep->dh_hist_spacing);
608 PS("separate-dhdl-file", enumValueToString(fep->separate_dhdl_file));
609 PS("dhdl-derivatives", enumValueToString(fep->dhdl_derivatives));
612 static void pr_pull(FILE* fp, int indent, const pull_params_t& pull)
616 PR("pull-cylinder-r", pull.cylinder_r);
617 PR("pull-constr-tol", pull.constr_tol);
618 PS("pull-print-COM", EBOOL(pull.bPrintCOM));
619 PS("pull-print-ref-value", EBOOL(pull.bPrintRefValue));
620 PS("pull-print-components", EBOOL(pull.bPrintComp));
621 PI("pull-nstxout", pull.nstxout);
622 PI("pull-nstfout", pull.nstfout);
623 PS("pull-pbc-ref-prev-step-com", EBOOL(pull.bSetPbcRefToPrevStepCOM));
624 PS("pull-xout-average", EBOOL(pull.bXOutAverage));
625 PS("pull-fout-average", EBOOL(pull.bFOutAverage));
626 PI("pull-ngroups", pull.ngroup);
627 for (g = 0; g < pull.ngroup; g++)
629 pr_pull_group(fp, indent, g, &pull.group[g]);
631 PI("pull-ncoords", pull.ncoord);
632 for (g = 0; g < pull.ncoord; g++)
634 pr_pull_coord(fp, indent, g, &pull.coord[g]);
638 static void pr_awh_bias_dim(FILE* fp, int indent, const gmx::AwhDimParams& awhDimParams, const char* prefix)
640 pr_indent(fp, indent);
642 fprintf(fp, "%s:\n", prefix);
643 PS("coord-provider", enumValueToString(awhDimParams.coordinateProvider()));
644 PI("coord-index", awhDimParams.coordinateIndex() + 1);
645 PR("start", awhDimParams.origin());
646 PR("end", awhDimParams.end());
647 PR("period", awhDimParams.period());
648 PR("force-constant", awhDimParams.forceConstant());
649 PR("diffusion", awhDimParams.diffusion());
650 PR("cover-diameter", awhDimParams.coverDiameter());
653 static void pr_awh_bias(FILE* fp, int indent, const gmx::AwhBiasParams& awhBiasParams, const char* prefix)
657 sprintf(opt, "%s-error-init", prefix);
658 PR(opt, awhBiasParams.initialErrorEstimate());
659 sprintf(opt, "%s-growth", prefix);
660 PS(opt, enumValueToString(awhBiasParams.growthType()));
661 sprintf(opt, "%s-target", prefix);
662 PS(opt, enumValueToString(awhBiasParams.targetDistribution()));
663 sprintf(opt, "%s-target-beta-scalng", prefix);
664 PR(opt, awhBiasParams.targetBetaScaling());
665 sprintf(opt, "%s-target-cutoff", prefix);
666 PR(opt, awhBiasParams.targetCutoff());
667 sprintf(opt, "%s-user-data", prefix);
668 PS(opt, EBOOL(awhBiasParams.userPMFEstimate()));
669 sprintf(opt, "%s-share-group", prefix);
670 PI(opt, awhBiasParams.shareGroup());
671 sprintf(opt, "%s-equilibrate-histogram", prefix);
672 PS(opt, EBOOL(awhBiasParams.equilibrateHistogram()));
673 sprintf(opt, "%s-ndim", prefix);
674 PI(opt, awhBiasParams.ndim());
677 for (const auto& dimParam : awhBiasParams.dimParams())
679 char prefixdim[STRLEN];
680 sprintf(prefixdim, "%s-dim%d", prefix, d + 1);
681 pr_awh_bias_dim(fp, indent, dimParam, prefixdim);
686 static void pr_awh(FILE* fp, int indent, gmx::AwhParams* awhParams)
688 PS("awh-potential", enumValueToString(awhParams->potential()));
689 PI("awh-seed", awhParams->seed());
690 PI("awh-nstout", awhParams->nstout());
691 PI("awh-nstsample", awhParams->nstSampleCoord());
692 PI("awh-nsamples-update", awhParams->numSamplesUpdateFreeEnergy());
693 PS("awh-share-bias-multisim", EBOOL(awhParams->shareBiasMultisim()));
694 PI("awh-nbias", awhParams->numBias());
697 for (const auto& awhBiasParam : awhParams->awhBiasParams())
699 auto prefix = gmx::formatString("awh%d", k + 1);
700 pr_awh_bias(fp, indent, awhBiasParam, prefix.c_str());
705 static void pr_rotgrp(FILE* fp, int indent, int g, const t_rotgrp* rotg)
707 pr_indent(fp, indent);
708 fprintf(fp, "rot-group %d:\n", g);
710 PS("rot-type", enumValueToString(rotg->eType));
711 PS("rot-massw", EBOOL(rotg->bMassW));
712 pr_ivec_block(fp, indent, "atom", rotg->ind, rotg->nat, TRUE);
713 pr_rvecs(fp, indent, "x-ref", rotg->x_ref, rotg->nat);
714 pr_rvec(fp, indent, "rot-vec", rotg->inputVec, DIM, TRUE);
715 pr_rvec(fp, indent, "rot-pivot", rotg->pivot, DIM, TRUE);
716 PR("rot-rate", rotg->rate);
717 PR("rot-k", rotg->k);
718 PR("rot-slab-dist", rotg->slab_dist);
719 PR("rot-min-gauss", rotg->min_gaussian);
720 PR("rot-eps", rotg->eps);
721 PS("rot-fit-method", enumValueToString(rotg->eFittype));
722 PI("rot-potfit-nstep", rotg->PotAngle_nstep);
723 PR("rot-potfit-step", rotg->PotAngle_step);
726 static void pr_rot(FILE* fp, int indent, const t_rot* rot)
730 PI("rot-nstrout", rot->nstrout);
731 PI("rot-nstsout", rot->nstsout);
732 PI("rot-ngroups", rot->ngrp);
733 for (g = 0; g < rot->ngrp; g++)
735 pr_rotgrp(fp, indent, g, &rot->grp[g]);
740 static void pr_swap(FILE* fp, int indent, const t_swapcoords* swap)
744 /* Enums for better readability of the code */
752 PI("swap-frequency", swap->nstswap);
754 /* The split groups that define the compartments */
755 for (int j = 0; j < 2; j++)
757 snprintf(str, STRLEN, "massw_split%d", j);
758 PS(str, EBOOL(swap->massw_split[j]));
759 snprintf(str, STRLEN, "split atoms group %d", j);
760 pr_ivec_block(fp, indent, str, swap->grp[j].ind, swap->grp[j].nat, TRUE);
763 /* The solvent group */
767 swap->grp[static_cast<int>(SwapGroupSplittingType::Solvent)].molname);
771 swap->grp[static_cast<int>(SwapGroupSplittingType::Solvent)].ind,
772 swap->grp[static_cast<int>(SwapGroupSplittingType::Solvent)].nat,
775 /* Now print the indices for all the ion groups: */
776 for (int ig = static_cast<int>(SwapGroupSplittingType::Count); ig < swap->ngrp; ig++)
778 snprintf(str, STRLEN, "ion group %s", swap->grp[ig].molname);
779 pr_ivec_block(fp, indent, str, swap->grp[ig].ind, swap->grp[ig].nat, TRUE);
782 PR("cyl0-r", swap->cyl0r);
783 PR("cyl0-up", swap->cyl0u);
784 PR("cyl0-down", swap->cyl0l);
785 PR("cyl1-r", swap->cyl1r);
786 PR("cyl1-up", swap->cyl1u);
787 PR("cyl1-down", swap->cyl1l);
788 PI("coupl-steps", swap->nAverage);
790 /* Print the requested ion counts for both compartments */
791 for (int ic = eCompA; ic <= eCompB; ic++)
793 for (int ig = static_cast<int>(SwapGroupSplittingType::Count); ig < swap->ngrp; ig++)
795 snprintf(str, STRLEN, "%s-in-%c", swap->grp[ig].molname, 'A' + ic);
796 PI(str, swap->grp[ig].nmolReq[ic]);
800 PR("threshold", swap->threshold);
801 PR("bulk-offsetA", swap->bulkOffset[eCompA]);
802 PR("bulk-offsetB", swap->bulkOffset[eCompB]);
806 static void pr_imd(FILE* fp, int indent, const t_IMD* imd)
808 PI("IMD-atoms", imd->nat);
809 pr_ivec_block(fp, indent, "atom", imd->ind, imd->nat, TRUE);
813 void pr_inputrec(FILE* fp, int indent, const char* title, const t_inputrec* ir, gmx_bool bMDPformat)
815 const char* infbuf = "inf";
817 if (available(fp, ir, indent, title))
821 indent = pr_title(fp, indent, title);
823 /* Try to make this list appear in the same order as the
824 * options are written in the default mdout.mdp, and with
825 * the same user-exposed names to facilitate debugging.
827 PS("integrator", enumValueToString(ir->eI));
828 PR("tinit", ir->init_t);
829 PR("dt", ir->delta_t);
830 PSTEP("nsteps", ir->nsteps);
831 PSTEP("init-step", ir->init_step);
832 PI("simulation-part", ir->simulation_part);
833 PS("mts", EBOOL(ir->useMts));
836 for (int mtsIndex = 1; mtsIndex < static_cast<int>(ir->mtsLevels.size()); mtsIndex++)
838 const auto& mtsLevel = ir->mtsLevels[mtsIndex];
839 const std::string forceKey = gmx::formatString("mts-level%d-forces", mtsIndex + 1);
840 std::string forceGroups;
841 for (int i = 0; i < static_cast<int>(gmx::MtsForceGroups::Count); i++)
843 if (mtsLevel.forceGroups[i])
845 if (!forceGroups.empty())
849 forceGroups += gmx::mtsForceGroupNames[i];
852 PS(forceKey.c_str(), forceGroups.c_str());
853 const std::string factorKey = gmx::formatString("mts-level%d-factor", mtsIndex + 1);
854 PI(factorKey.c_str(), mtsLevel.stepFactor);
857 PS("comm-mode", enumValueToString(ir->comm_mode));
858 PI("nstcomm", ir->nstcomm);
860 /* Langevin dynamics */
861 PR("bd-fric", ir->bd_fric);
862 PSTEP("ld-seed", ir->ld_seed);
864 /* Energy minimization */
865 PR("emtol", ir->em_tol);
866 PR("emstep", ir->em_stepsize);
867 PI("niter", ir->niter);
868 PR("fcstep", ir->fc_stepsize);
869 PI("nstcgsteep", ir->nstcgsteep);
870 PI("nbfgscorr", ir->nbfgscorr);
872 /* Test particle insertion */
873 PR("rtpi", ir->rtpi);
876 PI("nstxout", ir->nstxout);
877 PI("nstvout", ir->nstvout);
878 PI("nstfout", ir->nstfout);
879 PI("nstlog", ir->nstlog);
880 PI("nstcalcenergy", ir->nstcalcenergy);
881 PI("nstenergy", ir->nstenergy);
882 PI("nstxout-compressed", ir->nstxout_compressed);
883 PR("compressed-x-precision", ir->x_compression_precision);
885 /* Neighborsearching parameters */
886 PS("cutoff-scheme", enumValueToString(ir->cutoff_scheme));
887 PI("nstlist", ir->nstlist);
888 PS("pbc", c_pbcTypeNames[ir->pbcType].c_str());
889 PS("periodic-molecules", EBOOL(ir->bPeriodicMols));
890 PR("verlet-buffer-tolerance", ir->verletbuf_tol);
891 PR("rlist", ir->rlist);
893 /* Options for electrostatics and VdW */
894 PS("coulombtype", enumValueToString(ir->coulombtype));
895 PS("coulomb-modifier", enumValueToString(ir->coulomb_modifier));
896 PR("rcoulomb-switch", ir->rcoulomb_switch);
897 PR("rcoulomb", ir->rcoulomb);
898 if (ir->epsilon_r != 0)
900 PR("epsilon-r", ir->epsilon_r);
904 PS("epsilon-r", infbuf);
906 if (ir->epsilon_rf != 0)
908 PR("epsilon-rf", ir->epsilon_rf);
912 PS("epsilon-rf", infbuf);
914 PS("vdw-type", enumValueToString(ir->vdwtype));
915 PS("vdw-modifier", enumValueToString(ir->vdw_modifier));
916 PR("rvdw-switch", ir->rvdw_switch);
917 PR("rvdw", ir->rvdw);
918 PS("DispCorr", enumValueToString(ir->eDispCorr));
919 PR("table-extension", ir->tabext);
921 PR("fourierspacing", ir->fourier_spacing);
922 PI("fourier-nx", ir->nkx);
923 PI("fourier-ny", ir->nky);
924 PI("fourier-nz", ir->nkz);
925 PI("pme-order", ir->pme_order);
926 PR("ewald-rtol", ir->ewald_rtol);
927 PR("ewald-rtol-lj", ir->ewald_rtol_lj);
928 PS("lj-pme-comb-rule", enumValueToString(ir->ljpme_combination_rule));
929 PS("ewald-geometry", enumValueToString(ir->ewald_geometry));
930 PR("epsilon-surface", ir->epsilon_surface);
932 /* Options for weak coupling algorithms */
933 PS("tcoupl", enumValueToString(ir->etc));
934 PI("nsttcouple", ir->nsttcouple);
935 PI("nh-chain-length", ir->opts.nhchainlength);
936 PS("print-nose-hoover-chain-variables", EBOOL(ir->bPrintNHChains));
938 PS("pcoupl", enumValueToString(ir->epc));
939 PS("pcoupltype", enumValueToString(ir->epct));
940 PI("nstpcouple", ir->nstpcouple);
941 PR("tau-p", ir->tau_p);
942 pr_matrix(fp, indent, "compressibility", ir->compress, bMDPformat);
943 pr_matrix(fp, indent, "ref-p", ir->ref_p, bMDPformat);
944 PS("refcoord-scaling", enumValueToString(ir->refcoord_scaling));
948 fprintf(fp, "posres-com = %g %g %g\n", ir->posres_com[XX], ir->posres_com[YY], ir->posres_com[ZZ]);
949 fprintf(fp, "posres-comB = %g %g %g\n", ir->posres_comB[XX], ir->posres_comB[YY], ir->posres_comB[ZZ]);
953 pr_rvec(fp, indent, "posres-com", ir->posres_com, DIM, TRUE);
954 pr_rvec(fp, indent, "posres-comB", ir->posres_comB, DIM, TRUE);
958 PS("QMMM", EBOOL(ir->bQMMM));
959 fprintf(fp, "%s:\n", "qm-opts");
960 pr_int(fp, indent, "ngQM", ir->opts.ngQM);
962 /* CONSTRAINT OPTIONS */
963 PS("constraint-algorithm", enumValueToString(ir->eConstrAlg));
964 PS("continuation", EBOOL(ir->bContinuation));
966 PS("Shake-SOR", EBOOL(ir->bShakeSOR));
967 PR("shake-tol", ir->shake_tol);
968 PI("lincs-order", ir->nProjOrder);
969 PI("lincs-iter", ir->nLincsIter);
970 PR("lincs-warnangle", ir->LincsWarnAngle);
973 PI("nwall", ir->nwall);
974 PS("wall-type", enumValueToString(ir->wall_type));
975 PR("wall-r-linpot", ir->wall_r_linpot);
977 PI("wall-atomtype[0]", ir->wall_atomtype[0]);
978 PI("wall-atomtype[1]", ir->wall_atomtype[1]);
980 PR("wall-density[0]", ir->wall_density[0]);
981 PR("wall-density[1]", ir->wall_density[1]);
982 PR("wall-ewald-zfac", ir->wall_ewald_zfac);
985 PS("pull", EBOOL(ir->bPull));
988 pr_pull(fp, indent, *ir->pull);
992 PS("awh", EBOOL(ir->bDoAwh));
995 pr_awh(fp, indent, ir->awhParams.get());
998 /* ENFORCED ROTATION */
999 PS("rotation", EBOOL(ir->bRot));
1002 pr_rot(fp, indent, ir->rot);
1005 /* INTERACTIVE MD */
1006 PS("interactiveMD", EBOOL(ir->bIMD));
1009 pr_imd(fp, indent, ir->imd);
1012 /* NMR refinement stuff */
1013 PS("disre", enumValueToString(ir->eDisre));
1014 PS("disre-weighting", enumValueToString(ir->eDisreWeighting));
1015 PS("disre-mixed", EBOOL(ir->bDisreMixed));
1016 PR("dr-fc", ir->dr_fc);
1017 PR("dr-tau", ir->dr_tau);
1018 PR("nstdisreout", ir->nstdisreout);
1020 PR("orire-fc", ir->orires_fc);
1021 PR("orire-tau", ir->orires_tau);
1022 PR("nstorireout", ir->nstorireout);
1024 /* FREE ENERGY VARIABLES */
1025 PS("free-energy", enumValueToString(ir->efep));
1026 if (ir->efep != FreeEnergyPerturbationType::No || ir->bSimTemp)
1028 pr_fepvals(fp, indent, ir->fepvals.get(), bMDPformat);
1032 pr_expandedvals(fp, indent, ir->expandedvals.get(), ir->fepvals->n_lambda);
1035 /* NON-equilibrium MD stuff */
1036 PR("cos-acceleration", ir->cos_accel);
1037 pr_matrix(fp, indent, "deform", ir->deform, bMDPformat);
1039 /* SIMULATED TEMPERING */
1040 PS("simulated-tempering", EBOOL(ir->bSimTemp));
1043 pr_simtempvals(fp, indent, ir->simtempvals.get(), ir->fepvals->n_lambda);
1046 /* ION/WATER SWAPPING FOR COMPUTATIONAL ELECTROPHYSIOLOGY */
1047 PS("swapcoords", enumValueToString(ir->eSwapCoords));
1048 if (ir->eSwapCoords != SwapType::No)
1050 pr_swap(fp, indent, ir->swap);
1053 /* USER-DEFINED THINGIES */
1054 PI("userint1", ir->userint1);
1055 PI("userint2", ir->userint2);
1056 PI("userint3", ir->userint3);
1057 PI("userint4", ir->userint4);
1058 PR("userreal1", ir->userreal1);
1059 PR("userreal2", ir->userreal2);
1060 PR("userreal3", ir->userreal3);
1061 PR("userreal4", ir->userreal4);
1065 gmx::TextWriter writer(fp);
1066 writer.wrapperSettings().setIndent(indent);
1067 gmx::dumpKeyValueTree(&writer, *ir->params);
1070 pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
1077 static void cmp_grpopts(FILE* fp, const t_grpopts* opt1, const t_grpopts* opt2, real ftol, real abstol)
1080 char buf1[256], buf2[256];
1082 cmp_int(fp, "inputrec->grpopts.ngtc", -1, opt1->ngtc, opt2->ngtc);
1083 cmp_int(fp, "inputrec->grpopts.ngfrz", -1, opt1->ngfrz, opt2->ngfrz);
1084 cmp_int(fp, "inputrec->grpopts.ngener", -1, opt1->ngener, opt2->ngener);
1085 for (i = 0; (i < std::min(opt1->ngtc, opt2->ngtc)); i++)
1087 cmp_real(fp, "inputrec->grpopts.nrdf", i, opt1->nrdf[i], opt2->nrdf[i], ftol, abstol);
1088 cmp_real(fp, "inputrec->grpopts.ref_t", i, opt1->ref_t[i], opt2->ref_t[i], ftol, abstol);
1089 cmp_real(fp, "inputrec->grpopts.tau_t", i, opt1->tau_t[i], opt2->tau_t[i], ftol, abstol);
1090 cmpEnum(fp, "inputrec->grpopts.annealing", opt1->annealing[i], opt2->annealing[i]);
1091 cmp_int(fp, "inputrec->grpopts.anneal_npoints", i, opt1->anneal_npoints[i], opt2->anneal_npoints[i]);
1092 if (opt1->anneal_npoints[i] == opt2->anneal_npoints[i])
1094 sprintf(buf1, "inputrec->grpopts.anneal_time[%d]", i);
1095 sprintf(buf2, "inputrec->grpopts.anneal_temp[%d]", i);
1096 for (j = 0; j < opt1->anneal_npoints[i]; j++)
1098 cmp_real(fp, buf1, j, opt1->anneal_time[i][j], opt2->anneal_time[i][j], ftol, abstol);
1099 cmp_real(fp, buf2, j, opt1->anneal_temp[i][j], opt2->anneal_temp[i][j], ftol, abstol);
1103 if (opt1->ngener == opt2->ngener)
1105 for (i = 0; i < opt1->ngener; i++)
1107 for (j = i; j < opt1->ngener; j++)
1109 sprintf(buf1, "inputrec->grpopts.egp_flags[%d]", i);
1110 cmp_int(fp, buf1, j, opt1->egp_flags[opt1->ngener * i + j], opt2->egp_flags[opt1->ngener * i + j]);
1114 for (i = 0; (i < std::min(opt1->ngfrz, opt2->ngfrz)); i++)
1116 cmp_ivec(fp, "inputrec->grpopts.nFreeze", i, opt1->nFreeze[i], opt2->nFreeze[i]);
1120 static void cmp_pull(FILE* fp)
1123 "WARNING: Both files use COM pulling, but comparing of the pull struct is not "
1124 "implemented (yet). The pull parameters could be the same or different.\n");
1127 static void cmp_awhDimParams(FILE* fp,
1128 const gmx::AwhDimParams& dimp1,
1129 const gmx::AwhDimParams& dimp2,
1134 /* Note that we have double index here, but the compare functions only
1135 * support one index, so here we only print the dim index and not the bias.
1138 "inputrec.awhParams->bias?->dim->coord_index",
1140 dimp1.coordinateIndex(),
1141 dimp2.coordinateIndex());
1142 cmp_double(fp, "inputrec->awhParams->bias?->dim->period", dimIndex, dimp1.period(), dimp2.period(), ftol, abstol);
1144 "inputrec->awhParams->bias?->dim->diffusion",
1150 cmp_double(fp, "inputrec->awhParams->bias?->dim->origin", dimIndex, dimp1.origin(), dimp2.origin(), ftol, abstol);
1151 cmp_double(fp, "inputrec->awhParams->bias?->dim->end", dimIndex, dimp1.end(), dimp2.end(), ftol, abstol);
1153 "inputrec->awhParams->bias?->dim->coord_value_init",
1155 dimp1.initialCoordinate(),
1156 dimp2.initialCoordinate(),
1160 "inputrec->awhParams->bias?->dim->coverDiameter",
1162 dimp1.coverDiameter(),
1163 dimp2.coverDiameter(),
1168 static void cmp_awhBiasParams(FILE* fp,
1169 const gmx::AwhBiasParams& bias1,
1170 const gmx::AwhBiasParams& bias2,
1175 cmp_int(fp, "inputrec->awhParams->ndim", biasIndex, bias1.ndim(), bias2.ndim());
1176 cmpEnum<gmx::AwhTargetType>(
1177 fp, "inputrec->awhParams->biaseTarget", bias1.targetDistribution(), bias2.targetDistribution());
1179 "inputrec->awhParams->biastargetBetaScaling",
1181 bias1.targetBetaScaling(),
1182 bias2.targetBetaScaling(),
1186 "inputrec->awhParams->biastargetCutoff",
1188 bias1.targetCutoff(),
1189 bias2.targetCutoff(),
1192 cmpEnum<gmx::AwhHistogramGrowthType>(
1193 fp, "inputrec->awhParams->biaseGrowth", bias1.growthType(), bias2.growthType());
1194 cmp_bool(fp, "inputrec->awhParams->biasbUserData", biasIndex, bias1.userPMFEstimate(), bias2.userPMFEstimate());
1196 "inputrec->awhParams->biaserror_initial",
1198 bias1.initialErrorEstimate(),
1199 bias2.initialErrorEstimate(),
1202 cmp_int(fp, "inputrec->awhParams->biasShareGroup", biasIndex, bias1.shareGroup(), bias2.shareGroup());
1204 const auto dimParams1 = bias1.dimParams();
1205 const auto dimParams2 = bias2.dimParams();
1206 for (int dim = 0; dim < std::min(bias1.ndim(), bias2.ndim()); dim++)
1208 cmp_awhDimParams(fp, dimParams1[dim], dimParams2[dim], dim, ftol, abstol);
1212 static void cmp_awhParams(FILE* fp, const gmx::AwhParams& awh1, const gmx::AwhParams& awh2, real ftol, real abstol)
1214 cmp_int(fp, "inputrec->awhParams->nbias", -1, awh1.numBias(), awh2.numBias());
1215 cmp_int64(fp, "inputrec->awhParams->seed", awh1.seed(), awh2.seed());
1216 cmp_int(fp, "inputrec->awhParams->nstout", -1, awh1.nstout(), awh2.nstout());
1217 cmp_int(fp, "inputrec->awhParams->nstsample_coord", -1, awh1.nstSampleCoord(), awh2.nstSampleCoord());
1219 "inputrec->awhParams->nsamples_update_free_energy",
1221 awh1.numSamplesUpdateFreeEnergy(),
1222 awh2.numSamplesUpdateFreeEnergy());
1223 cmpEnum<gmx::AwhPotentialType>(
1224 fp, "inputrec->awhParams->ePotential", awh1.potential(), awh2.potential());
1225 cmp_bool(fp, "inputrec->awhParams->shareBiasMultisim", -1, awh1.shareBiasMultisim(), awh2.shareBiasMultisim());
1227 if (awh1.numBias() == awh2.numBias())
1229 const auto awhBiasParams1 = awh1.awhBiasParams();
1230 const auto awhBiasParams2 = awh2.awhBiasParams();
1231 for (int bias = 0; bias < awh1.numBias(); bias++)
1233 cmp_awhBiasParams(fp, awhBiasParams1[bias], awhBiasParams2[bias], bias, ftol, abstol);
1238 static void cmp_simtempvals(FILE* fp,
1239 const t_simtemp* simtemp1,
1240 const t_simtemp* simtemp2,
1246 cmpEnum(fp, "inputrec->simtempvals->eSimTempScale", simtemp1->eSimTempScale, simtemp2->eSimTempScale);
1247 cmp_real(fp, "inputrec->simtempvals->simtemp_high", -1, simtemp1->simtemp_high, simtemp2->simtemp_high, ftol, abstol);
1248 cmp_real(fp, "inputrec->simtempvals->simtemp_low", -1, simtemp1->simtemp_low, simtemp2->simtemp_low, ftol, abstol);
1249 for (i = 0; i < n_lambda; i++)
1252 "inputrec->simtempvals->temperatures",
1254 simtemp1->temperatures[i],
1255 simtemp2->temperatures[i],
1261 static void cmp_expandedvals(FILE* fp,
1262 const t_expanded* expand1,
1263 const t_expanded* expand2,
1270 cmp_bool(fp, "inputrec->fepvals->bInit_weights", -1, expand1->bInit_weights, expand2->bInit_weights);
1271 cmp_bool(fp, "inputrec->fepvals->bWLoneovert", -1, expand1->bWLoneovert, expand2->bWLoneovert);
1273 for (i = 0; i < n_lambda; i++)
1276 "inputrec->expandedvals->init_lambda_weights",
1278 expand1->init_lambda_weights[i],
1279 expand2->init_lambda_weights[i],
1284 cmpEnum(fp, "inputrec->expandedvals->lambda-stats", expand1->elamstats, expand2->elamstats);
1285 cmpEnum(fp, "inputrec->expandedvals->lambda-mc-move", expand1->elmcmove, expand2->elmcmove);
1286 cmp_int(fp, "inputrec->expandedvals->lmc-repeats", -1, expand1->lmc_repeats, expand2->lmc_repeats);
1287 cmp_int(fp, "inputrec->expandedvals->lmc-gibbsdelta", -1, expand1->gibbsdeltalam, expand2->gibbsdeltalam);
1288 cmp_int(fp, "inputrec->expandedvals->lmc-forced-nstart", -1, expand1->lmc_forced_nstart, expand2->lmc_forced_nstart);
1289 cmpEnum(fp, "inputrec->expandedvals->lambda-weights-equil", expand1->elmceq, expand2->elmceq);
1291 "inputrec->expandedvals->,weight-equil-number-all-lambda",
1293 expand1->equil_n_at_lam,
1294 expand2->equil_n_at_lam);
1295 cmp_int(fp, "inputrec->expandedvals->weight-equil-number-samples", -1, expand1->equil_samples, expand2->equil_samples);
1296 cmp_int(fp, "inputrec->expandedvals->weight-equil-number-steps", -1, expand1->equil_steps, expand2->equil_steps);
1298 "inputrec->expandedvals->weight-equil-wl-delta",
1300 expand1->equil_wl_delta,
1301 expand2->equil_wl_delta,
1305 "inputrec->expandedvals->weight-equil-count-ratio",
1307 expand1->equil_ratio,
1308 expand2->equil_ratio,
1312 "inputrec->expandedvals->symmetrized-transition-matrix",
1314 expand1->bSymmetrizedTMatrix,
1315 expand2->bSymmetrizedTMatrix);
1316 cmp_int(fp, "inputrec->expandedvals->nstTij", -1, expand1->nstTij, expand2->nstTij);
1317 cmp_int(fp, "inputrec->expandedvals->mininum-var-min", -1, expand1->minvarmin, expand2->minvarmin); /*default is reasonable */
1318 cmp_int(fp, "inputrec->expandedvals->weight-c-range", -1, expand1->c_range, expand2->c_range); /* default is just C=0 */
1319 cmp_real(fp, "inputrec->expandedvals->wl-scale", -1, expand1->wl_scale, expand2->wl_scale, ftol, abstol);
1320 cmp_real(fp, "inputrec->expandedvals->init-wl-delta", -1, expand1->init_wl_delta, expand2->init_wl_delta, ftol, abstol);
1321 cmp_real(fp, "inputrec->expandedvals->wl-ratio", -1, expand1->wl_ratio, expand2->wl_ratio, ftol, abstol);
1322 cmp_int(fp, "inputrec->expandedvals->nstexpanded", -1, expand1->nstexpanded, expand2->nstexpanded);
1323 cmp_int(fp, "inputrec->expandedvals->lmc-seed", -1, expand1->lmc_seed, expand2->lmc_seed);
1324 cmp_real(fp, "inputrec->expandedvals->mc-temperature", -1, expand1->mc_temp, expand2->mc_temp, ftol, abstol);
1327 static void cmp_fepvals(FILE* fp, const t_lambda* fep1, const t_lambda* fep2, real ftol, real abstol)
1330 cmp_int(fp, "inputrec->nstdhdl", -1, fep1->nstdhdl, fep2->nstdhdl);
1331 cmp_double(fp, "inputrec->fepvals->init_fep_state", -1, fep1->init_fep_state, fep2->init_fep_state, ftol, abstol);
1332 cmp_double(fp, "inputrec->fepvals->delta_lambda", -1, fep1->delta_lambda, fep2->delta_lambda, ftol, abstol);
1333 cmp_int(fp, "inputrec->fepvals->n_lambda", -1, fep1->n_lambda, fep2->n_lambda);
1334 for (i = 0; i < static_cast<int>(FreeEnergyPerturbationCouplingType::Count); i++)
1336 for (j = 0; j < std::min(fep1->n_lambda, fep2->n_lambda); j++)
1339 "inputrec->fepvals->all_lambda",
1341 fep1->all_lambda[i][j],
1342 fep2->all_lambda[i][j],
1347 cmp_int(fp, "inputrec->fepvals->lambda_neighbors", 1, fep1->lambda_neighbors, fep2->lambda_neighbors);
1348 cmp_real(fp, "inputrec->fepvals->sc_alpha", -1, fep1->sc_alpha, fep2->sc_alpha, ftol, abstol);
1349 cmp_int(fp, "inputrec->fepvals->sc_power", -1, fep1->sc_power, fep2->sc_power);
1350 cmp_real(fp, "inputrec->fepvals->sc_r_power", -1, fep1->sc_r_power, fep2->sc_r_power, ftol, abstol);
1351 cmp_real(fp, "inputrec->fepvals->sc_sigma", -1, fep1->sc_sigma, fep2->sc_sigma, ftol, abstol);
1352 cmpEnum(fp, "inputrec->fepvals->edHdLPrintEnergy", fep1->edHdLPrintEnergy, fep1->edHdLPrintEnergy);
1353 cmp_bool(fp, "inputrec->fepvals->bScCoul", -1, fep1->bScCoul, fep1->bScCoul);
1354 cmpEnum(fp, "inputrec->separate_dhdl_file", fep1->separate_dhdl_file, fep2->separate_dhdl_file);
1355 cmpEnum(fp, "inputrec->dhdl_derivatives", fep1->dhdl_derivatives, fep2->dhdl_derivatives);
1356 cmp_int(fp, "inputrec->dh_hist_size", -1, fep1->dh_hist_size, fep2->dh_hist_size);
1357 cmp_double(fp, "inputrec->dh_hist_spacing", -1, fep1->dh_hist_spacing, fep2->dh_hist_spacing, ftol, abstol);
1360 void cmp_inputrec(FILE* fp, const t_inputrec* ir1, const t_inputrec* ir2, real ftol, real abstol)
1362 fprintf(fp, "comparing inputrec\n");
1364 /* gcc 2.96 doesnt like these defines at all, but issues a huge list
1365 * of warnings. Maybe it will change in future versions, but for the
1366 * moment I've spelled them out instead. /EL 000820
1367 * #define CIB(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
1368 * #define CII(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
1369 * #define CIR(s) cmp_real(fp,"inputrec->"#s,0,ir1->##s,ir2->##s,ftol)
1371 cmpEnum(fp, "inputrec->eI", ir1->eI, ir2->eI);
1372 cmp_int64(fp, "inputrec->nsteps", ir1->nsteps, ir2->nsteps);
1373 cmp_int64(fp, "inputrec->init_step", ir1->init_step, ir2->init_step);
1374 cmp_int(fp, "inputrec->simulation_part", -1, ir1->simulation_part, ir2->simulation_part);
1375 cmp_int(fp, "inputrec->mts", -1, static_cast<int>(ir1->useMts), static_cast<int>(ir2->useMts));
1376 if (ir1->useMts && ir2->useMts)
1378 cmp_int(fp, "inputrec->mts-levels", -1, ir1->mtsLevels.size(), ir2->mtsLevels.size());
1380 "inputrec->mts-level2-forces",
1382 ir1->mtsLevels[1].forceGroups.to_ulong(),
1383 ir2->mtsLevels[1].forceGroups.to_ulong());
1385 "inputrec->mts-level2-factor",
1387 ir1->mtsLevels[1].stepFactor,
1388 ir2->mtsLevels[1].stepFactor);
1390 cmp_int(fp, "inputrec->pbcType", -1, static_cast<int>(ir1->pbcType), static_cast<int>(ir2->pbcType));
1391 cmp_bool(fp, "inputrec->bPeriodicMols", -1, ir1->bPeriodicMols, ir2->bPeriodicMols);
1392 cmpEnum(fp, "inputrec->cutoff_scheme", ir1->cutoff_scheme, ir2->cutoff_scheme);
1393 cmp_int(fp, "inputrec->nstlist", -1, ir1->nstlist, ir2->nstlist);
1394 cmp_int(fp, "inputrec->nstcomm", -1, ir1->nstcomm, ir2->nstcomm);
1395 cmpEnum(fp, "inputrec->comm_mode", ir1->comm_mode, ir2->comm_mode);
1396 cmp_int(fp, "inputrec->nstlog", -1, ir1->nstlog, ir2->nstlog);
1397 cmp_int(fp, "inputrec->nstxout", -1, ir1->nstxout, ir2->nstxout);
1398 cmp_int(fp, "inputrec->nstvout", -1, ir1->nstvout, ir2->nstvout);
1399 cmp_int(fp, "inputrec->nstfout", -1, ir1->nstfout, ir2->nstfout);
1400 cmp_int(fp, "inputrec->nstcalcenergy", -1, ir1->nstcalcenergy, ir2->nstcalcenergy);
1401 cmp_int(fp, "inputrec->nstenergy", -1, ir1->nstenergy, ir2->nstenergy);
1402 cmp_int(fp, "inputrec->nstxout_compressed", -1, ir1->nstxout_compressed, ir2->nstxout_compressed);
1403 cmp_double(fp, "inputrec->init_t", -1, ir1->init_t, ir2->init_t, ftol, abstol);
1404 cmp_double(fp, "inputrec->delta_t", -1, ir1->delta_t, ir2->delta_t, ftol, abstol);
1406 "inputrec->x_compression_precision",
1408 ir1->x_compression_precision,
1409 ir2->x_compression_precision,
1412 cmp_real(fp, "inputrec->fourierspacing", -1, ir1->fourier_spacing, ir2->fourier_spacing, ftol, abstol);
1413 cmp_int(fp, "inputrec->nkx", -1, ir1->nkx, ir2->nkx);
1414 cmp_int(fp, "inputrec->nky", -1, ir1->nky, ir2->nky);
1415 cmp_int(fp, "inputrec->nkz", -1, ir1->nkz, ir2->nkz);
1416 cmp_int(fp, "inputrec->pme_order", -1, ir1->pme_order, ir2->pme_order);
1417 cmp_real(fp, "inputrec->ewald_rtol", -1, ir1->ewald_rtol, ir2->ewald_rtol, ftol, abstol);
1418 cmpEnum(fp, "inputrec->ewald_geometry", ir1->ewald_geometry, ir2->ewald_geometry);
1419 cmp_real(fp, "inputrec->epsilon_surface", -1, ir1->epsilon_surface, ir2->epsilon_surface, ftol, abstol);
1421 "inputrec->bContinuation",
1423 static_cast<int>(ir1->bContinuation),
1424 static_cast<int>(ir2->bContinuation));
1425 cmp_int(fp, "inputrec->bShakeSOR", -1, static_cast<int>(ir1->bShakeSOR), static_cast<int>(ir2->bShakeSOR));
1426 cmpEnum(fp, "inputrec->etc", ir1->etc, ir2->etc);
1428 "inputrec->bPrintNHChains",
1430 static_cast<int>(ir1->bPrintNHChains),
1431 static_cast<int>(ir2->bPrintNHChains));
1432 cmpEnum(fp, "inputrec->epc", ir1->epc, ir2->epc);
1433 cmpEnum(fp, "inputrec->epct", ir1->epct, ir2->epct);
1434 cmp_real(fp, "inputrec->tau_p", -1, ir1->tau_p, ir2->tau_p, ftol, abstol);
1435 cmp_rvec(fp, "inputrec->ref_p(x)", -1, ir1->ref_p[XX], ir2->ref_p[XX], ftol, abstol);
1436 cmp_rvec(fp, "inputrec->ref_p(y)", -1, ir1->ref_p[YY], ir2->ref_p[YY], ftol, abstol);
1437 cmp_rvec(fp, "inputrec->ref_p(z)", -1, ir1->ref_p[ZZ], ir2->ref_p[ZZ], ftol, abstol);
1438 cmp_rvec(fp, "inputrec->compress(x)", -1, ir1->compress[XX], ir2->compress[XX], ftol, abstol);
1439 cmp_rvec(fp, "inputrec->compress(y)", -1, ir1->compress[YY], ir2->compress[YY], ftol, abstol);
1440 cmp_rvec(fp, "inputrec->compress(z)", -1, ir1->compress[ZZ], ir2->compress[ZZ], ftol, abstol);
1441 cmpEnum(fp, "refcoord_scaling", ir1->refcoord_scaling, ir2->refcoord_scaling);
1442 cmp_rvec(fp, "inputrec->posres_com", -1, ir1->posres_com, ir2->posres_com, ftol, abstol);
1443 cmp_rvec(fp, "inputrec->posres_comB", -1, ir1->posres_comB, ir2->posres_comB, ftol, abstol);
1444 cmp_real(fp, "inputrec->verletbuf_tol", -1, ir1->verletbuf_tol, ir2->verletbuf_tol, ftol, abstol);
1445 cmp_real(fp, "inputrec->rlist", -1, ir1->rlist, ir2->rlist, ftol, abstol);
1446 cmp_real(fp, "inputrec->rtpi", -1, ir1->rtpi, ir2->rtpi, ftol, abstol);
1447 cmpEnum(fp, "inputrec->coulombtype", ir1->coulombtype, ir2->coulombtype);
1448 cmpEnum(fp, "inputrec->coulomb_modifier", ir1->coulomb_modifier, ir2->coulomb_modifier);
1449 cmp_real(fp, "inputrec->rcoulomb_switch", -1, ir1->rcoulomb_switch, ir2->rcoulomb_switch, ftol, abstol);
1450 cmp_real(fp, "inputrec->rcoulomb", -1, ir1->rcoulomb, ir2->rcoulomb, ftol, abstol);
1451 cmpEnum(fp, "inputrec->vdwtype", ir1->vdwtype, ir2->vdwtype);
1452 cmpEnum(fp, "inputrec->vdw_modifier", ir1->vdw_modifier, ir2->vdw_modifier);
1453 cmp_real(fp, "inputrec->rvdw_switch", -1, ir1->rvdw_switch, ir2->rvdw_switch, ftol, abstol);
1454 cmp_real(fp, "inputrec->rvdw", -1, ir1->rvdw, ir2->rvdw, ftol, abstol);
1455 cmp_real(fp, "inputrec->epsilon_r", -1, ir1->epsilon_r, ir2->epsilon_r, ftol, abstol);
1456 cmp_real(fp, "inputrec->epsilon_rf", -1, ir1->epsilon_rf, ir2->epsilon_rf, ftol, abstol);
1457 cmp_real(fp, "inputrec->tabext", -1, ir1->tabext, ir2->tabext, ftol, abstol);
1459 cmpEnum(fp, "inputrec->eDispCorr", ir1->eDispCorr, ir2->eDispCorr);
1460 cmp_real(fp, "inputrec->shake_tol", -1, ir1->shake_tol, ir2->shake_tol, ftol, abstol);
1461 cmpEnum(fp, "inputrec->efep", ir1->efep, ir2->efep);
1462 cmp_fepvals(fp, ir1->fepvals.get(), ir2->fepvals.get(), ftol, abstol);
1463 cmp_int(fp, "inputrec->bSimTemp", -1, static_cast<int>(ir1->bSimTemp), static_cast<int>(ir2->bSimTemp));
1464 if ((ir1->bSimTemp == ir2->bSimTemp) && (ir1->bSimTemp))
1467 ir1->simtempvals.get(),
1468 ir2->simtempvals.get(),
1469 std::min(ir1->fepvals->n_lambda, ir2->fepvals->n_lambda),
1473 cmp_int(fp, "inputrec->bExpanded", -1, static_cast<int>(ir1->bExpanded), static_cast<int>(ir2->bExpanded));
1474 if ((ir1->bExpanded == ir2->bExpanded) && (ir1->bExpanded))
1476 cmp_expandedvals(fp,
1477 ir1->expandedvals.get(),
1478 ir2->expandedvals.get(),
1479 std::min(ir1->fepvals->n_lambda, ir2->fepvals->n_lambda),
1483 cmp_int(fp, "inputrec->nwall", -1, ir1->nwall, ir2->nwall);
1484 cmpEnum(fp, "inputrec->wall_type", ir1->wall_type, ir2->wall_type);
1485 cmp_int(fp, "inputrec->wall_atomtype[0]", -1, ir1->wall_atomtype[0], ir2->wall_atomtype[0]);
1486 cmp_int(fp, "inputrec->wall_atomtype[1]", -1, ir1->wall_atomtype[1], ir2->wall_atomtype[1]);
1487 cmp_real(fp, "inputrec->wall_density[0]", -1, ir1->wall_density[0], ir2->wall_density[0], ftol, abstol);
1488 cmp_real(fp, "inputrec->wall_density[1]", -1, ir1->wall_density[1], ir2->wall_density[1], ftol, abstol);
1489 cmp_real(fp, "inputrec->wall_ewald_zfac", -1, ir1->wall_ewald_zfac, ir2->wall_ewald_zfac, ftol, abstol);
1491 cmp_bool(fp, "inputrec->bPull", -1, ir1->bPull, ir2->bPull);
1492 if (ir1->bPull && ir2->bPull)
1497 cmp_bool(fp, "inputrec->bDoAwh", -1, ir1->bDoAwh, ir2->bDoAwh);
1498 if (ir1->bDoAwh && ir2->bDoAwh)
1500 cmp_awhParams(fp, *ir1->awhParams, *ir2->awhParams, ftol, abstol);
1503 cmpEnum(fp, "inputrec->eDisre", ir1->eDisre, ir2->eDisre);
1504 cmp_real(fp, "inputrec->dr_fc", -1, ir1->dr_fc, ir2->dr_fc, ftol, abstol);
1505 cmpEnum(fp, "inputrec->eDisreWeighting", ir1->eDisreWeighting, ir2->eDisreWeighting);
1506 cmp_int(fp, "inputrec->bDisreMixed", -1, static_cast<int>(ir1->bDisreMixed), static_cast<int>(ir2->bDisreMixed));
1507 cmp_int(fp, "inputrec->nstdisreout", -1, ir1->nstdisreout, ir2->nstdisreout);
1508 cmp_real(fp, "inputrec->dr_tau", -1, ir1->dr_tau, ir2->dr_tau, ftol, abstol);
1509 cmp_real(fp, "inputrec->orires_fc", -1, ir1->orires_fc, ir2->orires_fc, ftol, abstol);
1510 cmp_real(fp, "inputrec->orires_tau", -1, ir1->orires_tau, ir2->orires_tau, ftol, abstol);
1511 cmp_int(fp, "inputrec->nstorireout", -1, ir1->nstorireout, ir2->nstorireout);
1512 cmp_real(fp, "inputrec->em_stepsize", -1, ir1->em_stepsize, ir2->em_stepsize, ftol, abstol);
1513 cmp_real(fp, "inputrec->em_tol", -1, ir1->em_tol, ir2->em_tol, ftol, abstol);
1514 cmp_int(fp, "inputrec->niter", -1, ir1->niter, ir2->niter);
1515 cmp_real(fp, "inputrec->fc_stepsize", -1, ir1->fc_stepsize, ir2->fc_stepsize, ftol, abstol);
1516 cmp_int(fp, "inputrec->nstcgsteep", -1, ir1->nstcgsteep, ir2->nstcgsteep);
1517 cmp_int(fp, "inputrec->nbfgscorr", 0, ir1->nbfgscorr, ir2->nbfgscorr);
1518 cmpEnum(fp, "inputrec->eConstrAlg", ir1->eConstrAlg, ir2->eConstrAlg);
1519 cmp_int(fp, "inputrec->nProjOrder", -1, ir1->nProjOrder, ir2->nProjOrder);
1520 cmp_real(fp, "inputrec->LincsWarnAngle", -1, ir1->LincsWarnAngle, ir2->LincsWarnAngle, ftol, abstol);
1521 cmp_int(fp, "inputrec->nLincsIter", -1, ir1->nLincsIter, ir2->nLincsIter);
1522 cmp_real(fp, "inputrec->bd_fric", -1, ir1->bd_fric, ir2->bd_fric, ftol, abstol);
1523 cmp_int64(fp, "inputrec->ld_seed", ir1->ld_seed, ir2->ld_seed);
1524 cmp_real(fp, "inputrec->cos_accel", -1, ir1->cos_accel, ir2->cos_accel, ftol, abstol);
1525 cmp_rvec(fp, "inputrec->deform(a)", -1, ir1->deform[XX], ir2->deform[XX], ftol, abstol);
1526 cmp_rvec(fp, "inputrec->deform(b)", -1, ir1->deform[YY], ir2->deform[YY], ftol, abstol);
1527 cmp_rvec(fp, "inputrec->deform(c)", -1, ir1->deform[ZZ], ir2->deform[ZZ], ftol, abstol);
1530 cmp_int(fp, "inputrec->userint1", -1, ir1->userint1, ir2->userint1);
1531 cmp_int(fp, "inputrec->userint2", -1, ir1->userint2, ir2->userint2);
1532 cmp_int(fp, "inputrec->userint3", -1, ir1->userint3, ir2->userint3);
1533 cmp_int(fp, "inputrec->userint4", -1, ir1->userint4, ir2->userint4);
1534 cmp_real(fp, "inputrec->userreal1", -1, ir1->userreal1, ir2->userreal1, ftol, abstol);
1535 cmp_real(fp, "inputrec->userreal2", -1, ir1->userreal2, ir2->userreal2, ftol, abstol);
1536 cmp_real(fp, "inputrec->userreal3", -1, ir1->userreal3, ir2->userreal3, ftol, abstol);
1537 cmp_real(fp, "inputrec->userreal4", -1, ir1->userreal4, ir2->userreal4, ftol, abstol);
1538 cmp_grpopts(fp, &(ir1->opts), &(ir2->opts), ftol, abstol);
1539 gmx::TextWriter writer(fp);
1540 gmx::compareKeyValueTrees(&writer, *ir1->params, *ir2->params, ftol, abstol);
1543 void comp_pull_AB(FILE* fp, const pull_params_t& pull, real ftol, real abstol)
1545 for (int i = 0; i < pull.ncoord; i++)
1547 fprintf(fp, "comparing pull coord %d\n", i);
1548 cmp_real(fp, "pull-coord->k", -1, pull.coord[i].k, pull.coord[i].kB, ftol, abstol);
1552 gmx_bool inputrecDeform(const t_inputrec* ir)
1554 return (ir->deform[XX][XX] != 0 || ir->deform[YY][YY] != 0 || ir->deform[ZZ][ZZ] != 0
1555 || ir->deform[YY][XX] != 0 || ir->deform[ZZ][XX] != 0 || ir->deform[ZZ][YY] != 0);
1558 gmx_bool inputrecDynamicBox(const t_inputrec* ir)
1560 return (ir->epc != PressureCoupling::No || ir->eI == IntegrationAlgorithm::TPI || inputrecDeform(ir));
1563 gmx_bool inputrecPreserveShape(const t_inputrec* ir)
1565 return (ir->epc != PressureCoupling::No && ir->deform[XX][XX] == 0
1566 && (ir->epct == PressureCouplingType::Isotropic
1567 || ir->epct == PressureCouplingType::SemiIsotropic));
1570 gmx_bool inputrecNeedMutot(const t_inputrec* ir)
1572 return ((ir->coulombtype == CoulombInteractionType::Ewald || EEL_PME(ir->coulombtype))
1573 && (ir->ewald_geometry == EwaldGeometry::ThreeDC || ir->epsilon_surface != 0));
1576 gmx_bool inputrecExclForces(const t_inputrec* ir)
1578 return (EEL_FULL(ir->coulombtype) || (EEL_RF(ir->coulombtype)));
1581 gmx_bool inputrecNptTrotter(const t_inputrec* ir)
1583 return (((ir->eI == IntegrationAlgorithm::VV) || (ir->eI == IntegrationAlgorithm::VVAK))
1584 && (ir->epc == PressureCoupling::Mttk) && (ir->etc == TemperatureCoupling::NoseHoover));
1587 gmx_bool inputrecNvtTrotter(const t_inputrec* ir)
1589 return (((ir->eI == IntegrationAlgorithm::VV) || (ir->eI == IntegrationAlgorithm::VVAK))
1590 && (ir->epc != PressureCoupling::Mttk) && (ir->etc == TemperatureCoupling::NoseHoover));
1593 gmx_bool inputrecNphTrotter(const t_inputrec* ir)
1595 return (((ir->eI == IntegrationAlgorithm::VV) || (ir->eI == IntegrationAlgorithm::VVAK))
1596 && (ir->epc == PressureCoupling::Mttk) && (ir->etc != TemperatureCoupling::NoseHoover));
1599 bool inputrecPbcXY2Walls(const t_inputrec* ir)
1601 return (ir->pbcType == PbcType::XY && ir->nwall == 2);
1604 bool inputrecFrozenAtoms(const t_inputrec* ir)
1606 return ((ir->opts.nFreeze != nullptr)
1607 && (ir->opts.ngfrz > 1 || ir->opts.nFreeze[0][XX] != 0 || ir->opts.nFreeze[0][YY] != 0
1608 || ir->opts.nFreeze[0][ZZ] != 0));
1611 bool integratorHasConservedEnergyQuantity(const t_inputrec* ir)
1614 { // NOLINT bugprone-branch-clone
1615 // Energy minimization or stochastic integrator: no conservation
1618 else if (ir->etc == TemperatureCoupling::No && ir->epc == PressureCoupling::No)
1620 // The total energy is conserved, no additional conserved quanitity
1625 // Shear stress with Parrinello-Rahman is not supported (tedious)
1627 ((ir->epc == PressureCoupling::ParrinelloRahman || ir->epc == PressureCoupling::Mttk)
1628 && (ir->ref_p[YY][XX] != 0 || ir->ref_p[ZZ][XX] != 0 || ir->ref_p[ZZ][YY] != 0));
1630 return !ETC_ANDERSEN(ir->etc) && !shearWithPR;
1634 bool integratorHasReferenceTemperature(const t_inputrec* ir)
1636 return ((ir->etc != TemperatureCoupling::No) || EI_SD(ir->eI)
1637 || (ir->eI == IntegrationAlgorithm::BD) || EI_TPI(ir->eI));
1640 int inputrec2nboundeddim(const t_inputrec* ir)
1642 if (inputrecPbcXY2Walls(ir))
1648 return numPbcDimensions(ir->pbcType);
1652 int ndof_com(const t_inputrec* ir)
1656 switch (ir->pbcType)
1659 case PbcType::No: n = 3; break;
1660 case PbcType::XY: n = (ir->nwall == 0 ? 3 : 2); break;
1661 case PbcType::Screw: n = 1; break;
1662 default: gmx_incons("Unknown pbc in calc_nrdf");
1668 real maxReferenceTemperature(const t_inputrec& ir)
1670 if (EI_ENERGY_MINIMIZATION(ir.eI) || ir.eI == IntegrationAlgorithm::NM)
1675 if (EI_MD(ir.eI) && ir.etc == TemperatureCoupling::No)
1680 /* SD and BD also use ref_t and tau_t for setting the reference temperature.
1681 * TPI can be treated as MD, since it needs an ensemble temperature.
1684 real maxTemperature = 0;
1685 for (int i = 0; i < ir.opts.ngtc; i++)
1687 if (ir.opts.tau_t[i] >= 0)
1689 maxTemperature = std::max(maxTemperature, ir.opts.ref_t[i]);
1693 return maxTemperature;
1696 bool haveEwaldSurfaceContribution(const t_inputrec& ir)
1698 return EEL_PME_EWALD(ir.coulombtype)
1699 && (ir.ewald_geometry == EwaldGeometry::ThreeDC || ir.epsilon_surface != 0);
1702 bool haveFreeEnergyType(const t_inputrec& ir, const int fepType)
1704 for (int i = 0; i < ir.fepvals->n_lambda; i++)
1706 if (ir.fepvals->all_lambda[fepType][i] > 0)