ab292ec5231c3c768cae3b427614cdc12dec386e
[alexxy/gromacs.git] / src / gromacs / mdtypes / inputrec.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #include "gmxpre.h"
39
40 #include "gromacs/utility/enumerationhelpers.h"
41 #include "inputrec.h"
42
43 #include <cstdio>
44 #include <cstdlib>
45 #include <cstring>
46
47 #include <algorithm>
48 #include <memory>
49 #include <numeric>
50
51 #include "gromacs/applied_forces/awh/read_params.h"
52 #include "gromacs/math/veccompare.h"
53 #include "gromacs/math/vecdump.h"
54 #include "gromacs/mdtypes/awh_params.h"
55 #include "gromacs/mdtypes/md_enums.h"
56 #include "gromacs/mdtypes/multipletimestepping.h"
57 #include "gromacs/mdtypes/pull_params.h"
58 #include "gromacs/pbcutil/pbc.h"
59 #include "gromacs/utility/compare.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/keyvaluetree.h"
63 #include "gromacs/utility/smalloc.h"
64 #include "gromacs/utility/snprintf.h"
65 #include "gromacs/utility/strconvert.h"
66 #include "gromacs/utility/stringutil.h"
67 #include "gromacs/utility/textwriter.h"
68 #include "gromacs/utility/txtdump.h"
69
70 //! Macro to select a bool name
71 #define EBOOL(e) gmx::boolToString(e)
72
73 /* Default values for nstcalcenergy, used when the are no other restrictions. */
74 constexpr int c_defaultNstCalcEnergy = 10;
75
76 /* The minimum number of integration steps required for reasonably accurate
77  * integration of first and second order coupling algorithms.
78  */
79 const int nstmin_berendsen_tcouple = 5;
80 const int nstmin_berendsen_pcouple = 10;
81 const int nstmin_harmonic          = 20;
82
83 /* Default values for T- and P- coupling intervals, used when the are no other
84  * restrictions.
85  */
86 constexpr int c_defaultNstTCouple = 10;
87 constexpr int c_defaultNstPCouple = 10;
88
89 t_inputrec::t_inputrec()
90 {
91     // TODO When this memset is removed, remove the suppression of
92     // gcc -Wno-class-memaccess in a CMakeLists.txt file.
93     std::memset(this, 0, sizeof(*this)); // NOLINT(bugprone-undefined-memory-manipulation)
94     fepvals      = std::make_unique<t_lambda>();
95     expandedvals = std::make_unique<t_expanded>();
96     simtempvals  = std::make_unique<t_simtemp>();
97 }
98
99 t_inputrec::~t_inputrec()
100 {
101     done_inputrec(this);
102 }
103
104 int ir_optimal_nstcalcenergy(const t_inputrec* ir)
105 {
106     int nst;
107
108     if (ir->nstlist > 0)
109     {
110         nst = ir->nstlist;
111     }
112     else
113     {
114         nst = c_defaultNstCalcEnergy;
115     }
116
117     if (ir->useMts)
118     {
119         nst = std::lcm(nst, ir->mtsLevels.back().stepFactor);
120     }
121
122     return nst;
123 }
124
125 int tcouple_min_integration_steps(TemperatureCoupling etc)
126 {
127     int n;
128
129     switch (etc)
130     {
131         case TemperatureCoupling::No: n = 0; break;
132         case TemperatureCoupling::Berendsen:
133         case TemperatureCoupling::Yes: n = nstmin_berendsen_tcouple; break;
134         case TemperatureCoupling::VRescale:
135             /* V-rescale supports instantaneous rescaling */
136             n = 0;
137             break;
138         case TemperatureCoupling::NoseHoover: n = nstmin_harmonic; break;
139         case TemperatureCoupling::Andersen:
140         case TemperatureCoupling::AndersenMassive: n = 1; break;
141         default: gmx_incons("Unknown etc value");
142     }
143
144     return n;
145 }
146
147 int ir_optimal_nsttcouple(const t_inputrec* ir)
148 {
149     int  nmin, nwanted, n;
150     real tau_min;
151     int  g;
152
153     nmin = tcouple_min_integration_steps(ir->etc);
154
155     nwanted = c_defaultNstTCouple;
156
157     tau_min = 1e20;
158     if (ir->etc != TemperatureCoupling::No)
159     {
160         for (g = 0; g < ir->opts.ngtc; g++)
161         {
162             if (ir->opts.tau_t[g] > 0)
163             {
164                 tau_min = std::min(tau_min, ir->opts.tau_t[g]);
165             }
166         }
167     }
168
169     if (nmin == 0 || ir->delta_t * nwanted <= tau_min)
170     {
171         n = nwanted;
172     }
173     else
174     {
175         n = static_cast<int>(tau_min / (ir->delta_t * nmin) + 0.001);
176         if (n < 1)
177         {
178             n = 1;
179         }
180         while (nwanted % n != 0)
181         {
182             n--;
183         }
184     }
185
186     return n;
187 }
188
189 int pcouple_min_integration_steps(PressureCoupling epc)
190 {
191     int n;
192
193     switch (epc)
194     {
195         case PressureCoupling::No: n = 0; break;
196         case PressureCoupling::Berendsen:
197         case PressureCoupling::CRescale:
198         case PressureCoupling::Isotropic: n = nstmin_berendsen_pcouple; break;
199         case PressureCoupling::ParrinelloRahman:
200         case PressureCoupling::Mttk: n = nstmin_harmonic; break;
201         default: gmx_incons("Unknown epc value");
202     }
203
204     return n;
205 }
206
207 int ir_optimal_nstpcouple(const t_inputrec* ir)
208 {
209     const int minIntegrationSteps = pcouple_min_integration_steps(ir->epc);
210
211     const int nwanted = c_defaultNstPCouple;
212
213     // With multiple time-stepping we can only compute the pressure at slowest steps
214     const int minNstPCouple = (ir->useMts ? ir->mtsLevels.back().stepFactor : 1);
215
216     int n;
217     if (minIntegrationSteps == 0 || ir->delta_t * nwanted <= ir->tau_p)
218     {
219         n = nwanted;
220     }
221     else
222     {
223         n = static_cast<int>(ir->tau_p / (ir->delta_t * minIntegrationSteps) + 0.001);
224         if (n < minNstPCouple)
225         {
226             n = minNstPCouple;
227         }
228         // Without MTS we try to make nstpcouple a "nice" number
229         if (!ir->useMts)
230         {
231             while (nwanted % n != 0)
232             {
233                 n--;
234             }
235         }
236     }
237
238     // With MTS, nstpcouple should be a multiple of the slowest MTS interval
239     if (ir->useMts)
240     {
241         n = n - (n % minNstPCouple);
242     }
243
244     return n;
245 }
246
247 gmx_bool ir_coulomb_switched(const t_inputrec* ir)
248 {
249     return (ir->coulombtype == CoulombInteractionType::Switch
250             || ir->coulombtype == CoulombInteractionType::Shift
251             || ir->coulombtype == CoulombInteractionType::PmeSwitch
252             || ir->coulombtype == CoulombInteractionType::PmeUserSwitch
253             || ir->coulomb_modifier == InteractionModifiers::PotSwitch
254             || ir->coulomb_modifier == InteractionModifiers::ForceSwitch);
255 }
256
257 gmx_bool ir_coulomb_is_zero_at_cutoff(const t_inputrec* ir)
258 {
259     return (ir->cutoff_scheme == CutoffScheme::Verlet || ir_coulomb_switched(ir)
260             || ir->coulomb_modifier != InteractionModifiers::None
261             || ir->coulombtype == CoulombInteractionType::RFZero);
262 }
263
264 gmx_bool ir_coulomb_might_be_zero_at_cutoff(const t_inputrec* ir)
265 {
266     return (ir_coulomb_is_zero_at_cutoff(ir) || ir->coulombtype == CoulombInteractionType::User
267             || ir->coulombtype == CoulombInteractionType::PmeUser);
268 }
269
270 gmx_bool ir_vdw_switched(const t_inputrec* ir)
271 {
272     return (ir->vdwtype == VanDerWaalsType::Switch || ir->vdwtype == VanDerWaalsType::Shift
273             || ir->vdw_modifier == InteractionModifiers::PotSwitch
274             || ir->vdw_modifier == InteractionModifiers::ForceSwitch);
275 }
276
277 gmx_bool ir_vdw_is_zero_at_cutoff(const t_inputrec* ir)
278 {
279     return (ir->cutoff_scheme == CutoffScheme::Verlet || ir_vdw_switched(ir)
280             || ir->vdw_modifier != InteractionModifiers::None);
281 }
282
283 gmx_bool ir_vdw_might_be_zero_at_cutoff(const t_inputrec* ir)
284 {
285     return (ir_vdw_is_zero_at_cutoff(ir) || ir->vdwtype == VanDerWaalsType::User);
286 }
287
288 static void done_t_rot(t_rot* rot)
289 {
290     if (rot == nullptr)
291     {
292         return;
293     }
294     if (rot->grp != nullptr)
295     {
296         for (int i = 0; i < rot->ngrp; i++)
297         {
298             sfree(rot->grp[i].ind);
299             sfree(rot->grp[i].x_ref);
300         }
301         sfree(rot->grp);
302     }
303     sfree(rot);
304 }
305
306 void done_inputrec(t_inputrec* ir)
307 {
308     sfree(ir->opts.nrdf);
309     sfree(ir->opts.ref_t);
310     for (int i = 0; i < ir->opts.ngtc; i++)
311     {
312         sfree(ir->opts.anneal_time[i]);
313         sfree(ir->opts.anneal_temp[i]);
314     }
315     sfree(ir->opts.annealing);
316     sfree(ir->opts.anneal_npoints);
317     sfree(ir->opts.anneal_time);
318     sfree(ir->opts.anneal_temp);
319     sfree(ir->opts.tau_t);
320     sfree(ir->opts.nFreeze);
321     sfree(ir->opts.egp_flags);
322
323     done_t_rot(ir->rot);
324     delete ir->params;
325 }
326
327 static void pr_grp_opts(FILE* out, int indent, const char* title, const t_grpopts* opts, gmx_bool bMDPformat)
328 {
329     int i, m, j;
330
331     if (!bMDPformat)
332     {
333         fprintf(out, "%s:\n", title);
334     }
335
336     pr_indent(out, indent);
337     fprintf(out, "nrdf%s", bMDPformat ? " = " : ":");
338     for (i = 0; (i < opts->ngtc); i++)
339     {
340         fprintf(out, "  %10g", opts->nrdf[i]);
341     }
342     fprintf(out, "\n");
343
344     pr_indent(out, indent);
345     fprintf(out, "ref-t%s", bMDPformat ? " = " : ":");
346     for (i = 0; (i < opts->ngtc); i++)
347     {
348         fprintf(out, "  %10g", opts->ref_t[i]);
349     }
350     fprintf(out, "\n");
351
352     pr_indent(out, indent);
353     fprintf(out, "tau-t%s", bMDPformat ? " = " : ":");
354     for (i = 0; (i < opts->ngtc); i++)
355     {
356         fprintf(out, "  %10g", opts->tau_t[i]);
357     }
358     fprintf(out, "\n");
359
360     /* Pretty-print the simulated annealing info */
361     fprintf(out, "annealing%s", bMDPformat ? " = " : ":");
362     for (i = 0; (i < opts->ngtc); i++)
363     {
364         fprintf(out, "  %10s", enumValueToString(opts->annealing[i]));
365     }
366     fprintf(out, "\n");
367
368     fprintf(out, "annealing-npoints%s", bMDPformat ? " = " : ":");
369     for (i = 0; (i < opts->ngtc); i++)
370     {
371         fprintf(out, "  %10d", opts->anneal_npoints[i]);
372     }
373     fprintf(out, "\n");
374
375     for (i = 0; (i < opts->ngtc); i++)
376     {
377         if (opts->anneal_npoints[i] > 0)
378         {
379             fprintf(out, "annealing-time [%d]:\t", i);
380             for (j = 0; (j < opts->anneal_npoints[i]); j++)
381             {
382                 fprintf(out, "  %10.1f", opts->anneal_time[i][j]);
383             }
384             fprintf(out, "\n");
385             fprintf(out, "annealing-temp [%d]:\t", i);
386             for (j = 0; (j < opts->anneal_npoints[i]); j++)
387             {
388                 fprintf(out, "  %10.1f", opts->anneal_temp[i][j]);
389             }
390             fprintf(out, "\n");
391         }
392     }
393
394     pr_indent(out, indent);
395     fprintf(out, "nfreeze:");
396     for (i = 0; (i < opts->ngfrz); i++)
397     {
398         for (m = 0; (m < DIM); m++)
399         {
400             fprintf(out, "  %10s", opts->nFreeze[i][m] ? "Y" : "N");
401         }
402     }
403     fprintf(out, "\n");
404
405
406     for (i = 0; (i < opts->ngener); i++)
407     {
408         pr_indent(out, indent);
409         fprintf(out, "energygrp-flags[%3d]:", i);
410         for (m = 0; (m < opts->ngener); m++)
411         {
412             fprintf(out, " %d", opts->egp_flags[opts->ngener * i + m]);
413         }
414         fprintf(out, "\n");
415     }
416
417     fflush(out);
418 }
419
420 static void pr_matrix(FILE* fp, int indent, const char* title, const rvec* m, gmx_bool bMDPformat)
421 {
422     if (bMDPformat)
423     {
424         fprintf(fp,
425                 "%-10s    = %g %g %g %g %g %g\n",
426                 title,
427                 m[XX][XX],
428                 m[YY][YY],
429                 m[ZZ][ZZ],
430                 m[XX][YY],
431                 m[XX][ZZ],
432                 m[YY][ZZ]);
433     }
434     else
435     {
436         pr_rvecs(fp, indent, title, m, DIM);
437     }
438 }
439
440 #define PS(t, s) pr_str(fp, indent, t, s)
441 #define PI(t, s) pr_int(fp, indent, t, s)
442 #define PSTEP(t, s) pr_int64(fp, indent, t, s)
443 #define PR(t, s) pr_real(fp, indent, t, s)
444 #define PD(t, s) pr_double(fp, indent, t, s)
445
446 static void pr_pull_group(FILE* fp, int indent, int g, const t_pull_group* pgrp)
447 {
448     pr_indent(fp, indent);
449     fprintf(fp, "pull-group %d:\n", g);
450     indent += 2;
451     pr_ivec_block(fp, indent, "atom", pgrp->ind.data(), pgrp->ind.size(), TRUE);
452     pr_rvec(fp, indent, "weight", pgrp->weight.data(), pgrp->weight.size(), TRUE);
453     PI("pbcatom", pgrp->pbcatom);
454 }
455
456 static void pr_pull_coord(FILE* fp, int indent, int c, const t_pull_coord* pcrd)
457 {
458     int g;
459
460     pr_indent(fp, indent);
461     fprintf(fp, "pull-coord %d:\n", c);
462     PS("type", enumValueToString(pcrd->eType));
463     if (pcrd->eType == PullingAlgorithm::External)
464     {
465         PS("potential-provider", pcrd->externalPotentialProvider.c_str());
466     }
467     PS("geometry", enumValueToString(pcrd->eGeom));
468     for (g = 0; g < pcrd->ngroup; g++)
469     {
470         char buf[STRLEN];
471
472         sprintf(buf, "group[%d]", g);
473         PI(buf, pcrd->group[g]);
474     }
475     pr_ivec(fp, indent, "dim", pcrd->dim, DIM, TRUE);
476     pr_rvec(fp, indent, "origin", pcrd->origin, DIM, TRUE);
477     pr_rvec(fp, indent, "vec", pcrd->vec, DIM, TRUE);
478     PS("start", EBOOL(pcrd->bStart));
479     PR("init", pcrd->init);
480     PR("rate", pcrd->rate);
481     PR("k", pcrd->k);
482     PR("kB", pcrd->kB);
483 }
484
485 static void pr_simtempvals(FILE* fp, int indent, const t_simtemp* simtemp, int n_lambda)
486 {
487     PS("simulated-tempering-scaling", enumValueToString(simtemp->eSimTempScale));
488     PR("sim-temp-low", simtemp->simtemp_low);
489     PR("sim-temp-high", simtemp->simtemp_high);
490     pr_rvec(fp, indent, "simulated tempering temperatures", simtemp->temperatures, n_lambda, TRUE);
491 }
492
493 static void pr_expandedvals(FILE* fp, int indent, const t_expanded* expand, int n_lambda)
494 {
495
496     PI("nstexpanded", expand->nstexpanded);
497     PS("lmc-stats", enumValueToString(expand->elamstats));
498     PS("lmc-move", enumValueToString(expand->elmcmove));
499     PS("lmc-weights-equil", enumValueToString(expand->elmceq));
500     if (expand->elmceq == LambdaWeightWillReachEquilibrium::NumAtLambda)
501     {
502         PI("weight-equil-number-all-lambda", expand->equil_n_at_lam);
503     }
504     if (expand->elmceq == LambdaWeightWillReachEquilibrium::Samples)
505     {
506         PI("weight-equil-number-samples", expand->equil_samples);
507     }
508     if (expand->elmceq == LambdaWeightWillReachEquilibrium::Steps)
509     {
510         PI("weight-equil-number-steps", expand->equil_steps);
511     }
512     if (expand->elmceq == LambdaWeightWillReachEquilibrium::WLDelta)
513     {
514         PR("weight-equil-wl-delta", expand->equil_wl_delta);
515     }
516     if (expand->elmceq == LambdaWeightWillReachEquilibrium::Ratio)
517     {
518         PR("weight-equil-count-ratio", expand->equil_ratio);
519     }
520     PI("lmc-seed", expand->lmc_seed);
521     PR("mc-temperature", expand->mc_temp);
522     PI("lmc-repeats", expand->lmc_repeats);
523     PI("lmc-gibbsdelta", expand->gibbsdeltalam);
524     PI("lmc-forced-nstart", expand->lmc_forced_nstart);
525     PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
526     PI("nst-transition-matrix", expand->nstTij);
527     PI("mininum-var-min", expand->minvarmin); /*default is reasonable */
528     PI("weight-c-range", expand->c_range);    /* default is just C=0 */
529     PR("wl-scale", expand->wl_scale);
530     PR("wl-ratio", expand->wl_ratio);
531     PR("init-wl-delta", expand->init_wl_delta);
532     PS("wl-oneovert", EBOOL(expand->bWLoneovert));
533
534     pr_indent(fp, indent);
535     pr_rvec(fp, indent, "init-lambda-weights", expand->init_lambda_weights.data(), n_lambda, TRUE);
536     PS("init-weights", EBOOL(expand->bInit_weights));
537 }
538
539 static void pr_fepvals(FILE* fp, int indent, const t_lambda* fep, gmx_bool bMDPformat)
540 {
541     int j;
542
543     PR("init-lambda", fep->init_lambda);
544     PI("init-lambda-state", fep->init_fep_state);
545     PR("delta-lambda", fep->delta_lambda);
546     PI("nstdhdl", fep->nstdhdl);
547
548     if (!bMDPformat)
549     {
550         PI("n-lambdas", fep->n_lambda);
551     }
552     if (fep->n_lambda > 0)
553     {
554         pr_indent(fp, indent);
555         fprintf(fp, "separate-dvdl%s\n", bMDPformat ? " = " : ":");
556         for (auto i : gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, bool>::keys())
557         {
558             fprintf(fp, "%18s = ", enumValueToString(i));
559             if (fep->separate_dvdl[i])
560             {
561                 fprintf(fp, "  TRUE");
562             }
563             else
564             {
565                 fprintf(fp, "  FALSE");
566             }
567             fprintf(fp, "\n");
568         }
569         fprintf(fp, "all-lambdas%s\n", bMDPformat ? " = " : ":");
570         for (auto key : gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, bool>::keys())
571         {
572             fprintf(fp, "%18s = ", enumValueToString(key));
573             int i = static_cast<int>(key);
574             for (j = 0; j < fep->n_lambda; j++)
575             {
576                 fprintf(fp, "  %10g", fep->all_lambda[i][j]);
577             }
578             fprintf(fp, "\n");
579         }
580     }
581     PI("calc-lambda-neighbors", fep->lambda_neighbors);
582     PS("dhdl-print-energy", enumValueToString(fep->edHdLPrintEnergy));
583     PR("sc-alpha", fep->sc_alpha);
584     PI("sc-power", fep->sc_power);
585     PR("sc-r-power", fep->sc_r_power);
586     PR("sc-sigma", fep->sc_sigma);
587     PR("sc-sigma-min", fep->sc_sigma_min);
588     PS("sc-coul", EBOOL(fep->bScCoul));
589     PI("dh-hist-size", fep->dh_hist_size);
590     PD("dh-hist-spacing", fep->dh_hist_spacing);
591     PS("separate-dhdl-file", enumValueToString(fep->separate_dhdl_file));
592     PS("dhdl-derivatives", enumValueToString(fep->dhdl_derivatives));
593 };
594
595 static void pr_pull(FILE* fp, int indent, const pull_params_t& pull)
596 {
597     int g;
598
599     PR("pull-cylinder-r", pull.cylinder_r);
600     PR("pull-constr-tol", pull.constr_tol);
601     PS("pull-print-COM", EBOOL(pull.bPrintCOM));
602     PS("pull-print-ref-value", EBOOL(pull.bPrintRefValue));
603     PS("pull-print-components", EBOOL(pull.bPrintComp));
604     PI("pull-nstxout", pull.nstxout);
605     PI("pull-nstfout", pull.nstfout);
606     PS("pull-pbc-ref-prev-step-com", EBOOL(pull.bSetPbcRefToPrevStepCOM));
607     PS("pull-xout-average", EBOOL(pull.bXOutAverage));
608     PS("pull-fout-average", EBOOL(pull.bFOutAverage));
609     PI("pull-ngroups", pull.ngroup);
610     for (g = 0; g < pull.ngroup; g++)
611     {
612         pr_pull_group(fp, indent, g, &pull.group[g]);
613     }
614     PI("pull-ncoords", pull.ncoord);
615     for (g = 0; g < pull.ncoord; g++)
616     {
617         pr_pull_coord(fp, indent, g, &pull.coord[g]);
618     }
619 }
620
621 static void pr_awh_bias_dim(FILE* fp, int indent, const gmx::AwhDimParams& awhDimParams, const char* prefix)
622 {
623     pr_indent(fp, indent);
624     indent++;
625     fprintf(fp, "%s:\n", prefix);
626     PS("coord-provider", enumValueToString(awhDimParams.coordinateProvider()));
627     PI("coord-index", awhDimParams.coordinateIndex() + 1);
628     PR("start", awhDimParams.origin());
629     PR("end", awhDimParams.end());
630     PR("period", awhDimParams.period());
631     PR("force-constant", awhDimParams.forceConstant());
632     PR("diffusion", awhDimParams.diffusion());
633     PR("cover-diameter", awhDimParams.coverDiameter());
634 }
635
636 static void pr_awh_bias(FILE* fp, int indent, const gmx::AwhBiasParams& awhBiasParams, const char* prefix)
637 {
638     char opt[STRLEN];
639
640     sprintf(opt, "%s-error-init", prefix);
641     PR(opt, awhBiasParams.initialErrorEstimate());
642     sprintf(opt, "%s-growth", prefix);
643     PS(opt, enumValueToString(awhBiasParams.growthType()));
644     sprintf(opt, "%s-target", prefix);
645     PS(opt, enumValueToString(awhBiasParams.targetDistribution()));
646     sprintf(opt, "%s-target-beta-scalng", prefix);
647     PR(opt, awhBiasParams.targetBetaScaling());
648     sprintf(opt, "%s-target-cutoff", prefix);
649     PR(opt, awhBiasParams.targetCutoff());
650     sprintf(opt, "%s-user-data", prefix);
651     PS(opt, EBOOL(awhBiasParams.userPMFEstimate()));
652     sprintf(opt, "%s-share-group", prefix);
653     PI(opt, awhBiasParams.shareGroup());
654     sprintf(opt, "%s-equilibrate-histogram", prefix);
655     PS(opt, EBOOL(awhBiasParams.equilibrateHistogram()));
656     sprintf(opt, "%s-ndim", prefix);
657     PI(opt, awhBiasParams.ndim());
658
659     int d = 0;
660     for (const auto& dimParam : awhBiasParams.dimParams())
661     {
662         char prefixdim[STRLEN];
663         sprintf(prefixdim, "%s-dim%d", prefix, d + 1);
664         pr_awh_bias_dim(fp, indent, dimParam, prefixdim);
665         d++;
666     }
667 }
668
669 static void pr_awh(FILE* fp, int indent, gmx::AwhParams* awhParams)
670 {
671     PS("awh-potential", enumValueToString(awhParams->potential()));
672     PI("awh-seed", awhParams->seed());
673     PI("awh-nstout", awhParams->nstout());
674     PI("awh-nstsample", awhParams->nstSampleCoord());
675     PI("awh-nsamples-update", awhParams->numSamplesUpdateFreeEnergy());
676     PS("awh-share-bias-multisim", EBOOL(awhParams->shareBiasMultisim()));
677     PI("awh-nbias", awhParams->numBias());
678
679     int k = 0;
680     for (const auto& awhBiasParam : awhParams->awhBiasParams())
681     {
682         auto prefix = gmx::formatString("awh%d", k + 1);
683         pr_awh_bias(fp, indent, awhBiasParam, prefix.c_str());
684         k++;
685     }
686 }
687
688 static void pr_rotgrp(FILE* fp, int indent, int g, const t_rotgrp* rotg)
689 {
690     pr_indent(fp, indent);
691     fprintf(fp, "rot-group %d:\n", g);
692     indent += 2;
693     PS("rot-type", enumValueToString(rotg->eType));
694     PS("rot-massw", EBOOL(rotg->bMassW));
695     pr_ivec_block(fp, indent, "atom", rotg->ind, rotg->nat, TRUE);
696     pr_rvecs(fp, indent, "x-ref", rotg->x_ref, rotg->nat);
697     pr_rvec(fp, indent, "rot-vec", rotg->inputVec, DIM, TRUE);
698     pr_rvec(fp, indent, "rot-pivot", rotg->pivot, DIM, TRUE);
699     PR("rot-rate", rotg->rate);
700     PR("rot-k", rotg->k);
701     PR("rot-slab-dist", rotg->slab_dist);
702     PR("rot-min-gauss", rotg->min_gaussian);
703     PR("rot-eps", rotg->eps);
704     PS("rot-fit-method", enumValueToString(rotg->eFittype));
705     PI("rot-potfit-nstep", rotg->PotAngle_nstep);
706     PR("rot-potfit-step", rotg->PotAngle_step);
707 }
708
709 static void pr_rot(FILE* fp, int indent, const t_rot* rot)
710 {
711     int g;
712
713     PI("rot-nstrout", rot->nstrout);
714     PI("rot-nstsout", rot->nstsout);
715     PI("rot-ngroups", rot->ngrp);
716     for (g = 0; g < rot->ngrp; g++)
717     {
718         pr_rotgrp(fp, indent, g, &rot->grp[g]);
719     }
720 }
721
722
723 static void pr_swap(FILE* fp, int indent, const t_swapcoords* swap)
724 {
725     char str[STRLEN];
726
727     /* Enums for better readability of the code */
728     enum
729     {
730         eCompA = 0,
731         eCompB
732     };
733
734
735     PI("swap-frequency", swap->nstswap);
736
737     /* The split groups that define the compartments */
738     for (int j = 0; j < 2; j++)
739     {
740         snprintf(str, STRLEN, "massw_split%d", j);
741         PS(str, EBOOL(swap->massw_split[j]));
742         snprintf(str, STRLEN, "split atoms group %d", j);
743         pr_ivec_block(fp, indent, str, swap->grp[j].ind, swap->grp[j].nat, TRUE);
744     }
745
746     /* The solvent group */
747     snprintf(str,
748              STRLEN,
749              "solvent group %s",
750              swap->grp[static_cast<int>(SwapGroupSplittingType::Solvent)].molname);
751     pr_ivec_block(fp,
752                   indent,
753                   str,
754                   swap->grp[static_cast<int>(SwapGroupSplittingType::Solvent)].ind,
755                   swap->grp[static_cast<int>(SwapGroupSplittingType::Solvent)].nat,
756                   TRUE);
757
758     /* Now print the indices for all the ion groups: */
759     for (int ig = static_cast<int>(SwapGroupSplittingType::Count); ig < swap->ngrp; ig++)
760     {
761         snprintf(str, STRLEN, "ion group %s", swap->grp[ig].molname);
762         pr_ivec_block(fp, indent, str, swap->grp[ig].ind, swap->grp[ig].nat, TRUE);
763     }
764
765     PR("cyl0-r", swap->cyl0r);
766     PR("cyl0-up", swap->cyl0u);
767     PR("cyl0-down", swap->cyl0l);
768     PR("cyl1-r", swap->cyl1r);
769     PR("cyl1-up", swap->cyl1u);
770     PR("cyl1-down", swap->cyl1l);
771     PI("coupl-steps", swap->nAverage);
772
773     /* Print the requested ion counts for both compartments */
774     for (int ic = eCompA; ic <= eCompB; ic++)
775     {
776         for (int ig = static_cast<int>(SwapGroupSplittingType::Count); ig < swap->ngrp; ig++)
777         {
778             snprintf(str, STRLEN, "%s-in-%c", swap->grp[ig].molname, 'A' + ic);
779             PI(str, swap->grp[ig].nmolReq[ic]);
780         }
781     }
782
783     PR("threshold", swap->threshold);
784     PR("bulk-offsetA", swap->bulkOffset[eCompA]);
785     PR("bulk-offsetB", swap->bulkOffset[eCompB]);
786 }
787
788
789 static void pr_imd(FILE* fp, int indent, const t_IMD* imd)
790 {
791     PI("IMD-atoms", imd->nat);
792     pr_ivec_block(fp, indent, "atom", imd->ind, imd->nat, TRUE);
793 }
794
795
796 void pr_inputrec(FILE* fp, int indent, const char* title, const t_inputrec* ir, gmx_bool bMDPformat)
797 {
798     const char* infbuf = "inf";
799
800     if (available(fp, ir, indent, title))
801     {
802         if (!bMDPformat)
803         {
804             indent = pr_title(fp, indent, title);
805         }
806         /* Try to make this list appear in the same order as the
807          * options are written in the default mdout.mdp, and with
808          * the same user-exposed names to facilitate debugging.
809          */
810         PS("integrator", enumValueToString(ir->eI));
811         PR("tinit", ir->init_t);
812         PR("dt", ir->delta_t);
813         PSTEP("nsteps", ir->nsteps);
814         PSTEP("init-step", ir->init_step);
815         PI("simulation-part", ir->simulation_part);
816         PS("mts", EBOOL(ir->useMts));
817         if (ir->useMts)
818         {
819             for (int mtsIndex = 1; mtsIndex < static_cast<int>(ir->mtsLevels.size()); mtsIndex++)
820             {
821                 const auto&       mtsLevel = ir->mtsLevels[mtsIndex];
822                 const std::string forceKey = gmx::formatString("mts-level%d-forces", mtsIndex + 1);
823                 std::string       forceGroups;
824                 for (int i = 0; i < static_cast<int>(gmx::MtsForceGroups::Count); i++)
825                 {
826                     if (mtsLevel.forceGroups[i])
827                     {
828                         if (!forceGroups.empty())
829                         {
830                             forceGroups += " ";
831                         }
832                         forceGroups += gmx::mtsForceGroupNames[i];
833                     }
834                 }
835                 PS(forceKey.c_str(), forceGroups.c_str());
836                 const std::string factorKey = gmx::formatString("mts-level%d-factor", mtsIndex + 1);
837                 PI(factorKey.c_str(), mtsLevel.stepFactor);
838             }
839         }
840         PS("comm-mode", enumValueToString(ir->comm_mode));
841         PI("nstcomm", ir->nstcomm);
842
843         /* Langevin dynamics */
844         PR("bd-fric", ir->bd_fric);
845         PSTEP("ld-seed", ir->ld_seed);
846
847         /* Energy minimization */
848         PR("emtol", ir->em_tol);
849         PR("emstep", ir->em_stepsize);
850         PI("niter", ir->niter);
851         PR("fcstep", ir->fc_stepsize);
852         PI("nstcgsteep", ir->nstcgsteep);
853         PI("nbfgscorr", ir->nbfgscorr);
854
855         /* Test particle insertion */
856         PR("rtpi", ir->rtpi);
857
858         /* Output control */
859         PI("nstxout", ir->nstxout);
860         PI("nstvout", ir->nstvout);
861         PI("nstfout", ir->nstfout);
862         PI("nstlog", ir->nstlog);
863         PI("nstcalcenergy", ir->nstcalcenergy);
864         PI("nstenergy", ir->nstenergy);
865         PI("nstxout-compressed", ir->nstxout_compressed);
866         PR("compressed-x-precision", ir->x_compression_precision);
867
868         /* Neighborsearching parameters */
869         PS("cutoff-scheme", enumValueToString(ir->cutoff_scheme));
870         PI("nstlist", ir->nstlist);
871         PS("pbc", c_pbcTypeNames[ir->pbcType].c_str());
872         PS("periodic-molecules", EBOOL(ir->bPeriodicMols));
873         PR("verlet-buffer-tolerance", ir->verletbuf_tol);
874         PR("rlist", ir->rlist);
875
876         /* Options for electrostatics and VdW */
877         PS("coulombtype", enumValueToString(ir->coulombtype));
878         PS("coulomb-modifier", enumValueToString(ir->coulomb_modifier));
879         PR("rcoulomb-switch", ir->rcoulomb_switch);
880         PR("rcoulomb", ir->rcoulomb);
881         if (ir->epsilon_r != 0)
882         {
883             PR("epsilon-r", ir->epsilon_r);
884         }
885         else
886         {
887             PS("epsilon-r", infbuf);
888         }
889         if (ir->epsilon_rf != 0)
890         {
891             PR("epsilon-rf", ir->epsilon_rf);
892         }
893         else
894         {
895             PS("epsilon-rf", infbuf);
896         }
897         PS("vdw-type", enumValueToString(ir->vdwtype));
898         PS("vdw-modifier", enumValueToString(ir->vdw_modifier));
899         PR("rvdw-switch", ir->rvdw_switch);
900         PR("rvdw", ir->rvdw);
901         PS("DispCorr", enumValueToString(ir->eDispCorr));
902         PR("table-extension", ir->tabext);
903
904         PR("fourierspacing", ir->fourier_spacing);
905         PI("fourier-nx", ir->nkx);
906         PI("fourier-ny", ir->nky);
907         PI("fourier-nz", ir->nkz);
908         PI("pme-order", ir->pme_order);
909         PR("ewald-rtol", ir->ewald_rtol);
910         PR("ewald-rtol-lj", ir->ewald_rtol_lj);
911         PS("lj-pme-comb-rule", enumValueToString(ir->ljpme_combination_rule));
912         PS("ewald-geometry", enumValueToString(ir->ewald_geometry));
913         PR("epsilon-surface", ir->epsilon_surface);
914
915         /* Options for weak coupling algorithms */
916         PS("tcoupl", enumValueToString(ir->etc));
917         PI("nsttcouple", ir->nsttcouple);
918         PI("nh-chain-length", ir->opts.nhchainlength);
919         PS("print-nose-hoover-chain-variables", EBOOL(ir->bPrintNHChains));
920
921         PS("pcoupl", enumValueToString(ir->epc));
922         PS("pcoupltype", enumValueToString(ir->epct));
923         PI("nstpcouple", ir->nstpcouple);
924         PR("tau-p", ir->tau_p);
925         pr_matrix(fp, indent, "compressibility", ir->compress, bMDPformat);
926         pr_matrix(fp, indent, "ref-p", ir->ref_p, bMDPformat);
927         PS("refcoord-scaling", enumValueToString(ir->refcoord_scaling));
928
929         if (bMDPformat)
930         {
931             fprintf(fp, "posres-com  = %g %g %g\n", ir->posres_com[XX], ir->posres_com[YY], ir->posres_com[ZZ]);
932             fprintf(fp, "posres-comB = %g %g %g\n", ir->posres_comB[XX], ir->posres_comB[YY], ir->posres_comB[ZZ]);
933         }
934         else
935         {
936             pr_rvec(fp, indent, "posres-com", ir->posres_com, DIM, TRUE);
937             pr_rvec(fp, indent, "posres-comB", ir->posres_comB, DIM, TRUE);
938         }
939
940         /* QMMM */
941         PS("QMMM", EBOOL(ir->bQMMM));
942         fprintf(fp, "%s:\n", "qm-opts");
943         pr_int(fp, indent, "ngQM", ir->opts.ngQM);
944
945         /* CONSTRAINT OPTIONS */
946         PS("constraint-algorithm", enumValueToString(ir->eConstrAlg));
947         PS("continuation", EBOOL(ir->bContinuation));
948
949         PS("Shake-SOR", EBOOL(ir->bShakeSOR));
950         PR("shake-tol", ir->shake_tol);
951         PI("lincs-order", ir->nProjOrder);
952         PI("lincs-iter", ir->nLincsIter);
953         PR("lincs-warnangle", ir->LincsWarnAngle);
954
955         /* Walls */
956         PI("nwall", ir->nwall);
957         PS("wall-type", enumValueToString(ir->wall_type));
958         PR("wall-r-linpot", ir->wall_r_linpot);
959         /* wall-atomtype */
960         PI("wall-atomtype[0]", ir->wall_atomtype[0]);
961         PI("wall-atomtype[1]", ir->wall_atomtype[1]);
962         /* wall-density */
963         PR("wall-density[0]", ir->wall_density[0]);
964         PR("wall-density[1]", ir->wall_density[1]);
965         PR("wall-ewald-zfac", ir->wall_ewald_zfac);
966
967         /* COM PULLING */
968         PS("pull", EBOOL(ir->bPull));
969         if (ir->bPull)
970         {
971             pr_pull(fp, indent, *ir->pull);
972         }
973
974         /* AWH BIASING */
975         PS("awh", EBOOL(ir->bDoAwh));
976         if (ir->bDoAwh)
977         {
978             pr_awh(fp, indent, ir->awhParams.get());
979         }
980
981         /* ENFORCED ROTATION */
982         PS("rotation", EBOOL(ir->bRot));
983         if (ir->bRot)
984         {
985             pr_rot(fp, indent, ir->rot);
986         }
987
988         /* INTERACTIVE MD */
989         PS("interactiveMD", EBOOL(ir->bIMD));
990         if (ir->bIMD)
991         {
992             pr_imd(fp, indent, ir->imd);
993         }
994
995         /* NMR refinement stuff */
996         PS("disre", enumValueToString(ir->eDisre));
997         PS("disre-weighting", enumValueToString(ir->eDisreWeighting));
998         PS("disre-mixed", EBOOL(ir->bDisreMixed));
999         PR("dr-fc", ir->dr_fc);
1000         PR("dr-tau", ir->dr_tau);
1001         PR("nstdisreout", ir->nstdisreout);
1002
1003         PR("orire-fc", ir->orires_fc);
1004         PR("orire-tau", ir->orires_tau);
1005         PR("nstorireout", ir->nstorireout);
1006
1007         /* FREE ENERGY VARIABLES */
1008         PS("free-energy", enumValueToString(ir->efep));
1009         if (ir->efep != FreeEnergyPerturbationType::No || ir->bSimTemp)
1010         {
1011             pr_fepvals(fp, indent, ir->fepvals.get(), bMDPformat);
1012         }
1013         if (ir->bExpanded)
1014         {
1015             pr_expandedvals(fp, indent, ir->expandedvals.get(), ir->fepvals->n_lambda);
1016         }
1017
1018         /* NON-equilibrium MD stuff */
1019         PR("cos-acceleration", ir->cos_accel);
1020         pr_matrix(fp, indent, "deform", ir->deform, bMDPformat);
1021
1022         /* SIMULATED TEMPERING */
1023         PS("simulated-tempering", EBOOL(ir->bSimTemp));
1024         if (ir->bSimTemp)
1025         {
1026             pr_simtempvals(fp, indent, ir->simtempvals.get(), ir->fepvals->n_lambda);
1027         }
1028
1029         /* ION/WATER SWAPPING FOR COMPUTATIONAL ELECTROPHYSIOLOGY */
1030         PS("swapcoords", enumValueToString(ir->eSwapCoords));
1031         if (ir->eSwapCoords != SwapType::No)
1032         {
1033             pr_swap(fp, indent, ir->swap);
1034         }
1035
1036         /* USER-DEFINED THINGIES */
1037         PI("userint1", ir->userint1);
1038         PI("userint2", ir->userint2);
1039         PI("userint3", ir->userint3);
1040         PI("userint4", ir->userint4);
1041         PR("userreal1", ir->userreal1);
1042         PR("userreal2", ir->userreal2);
1043         PR("userreal3", ir->userreal3);
1044         PR("userreal4", ir->userreal4);
1045
1046         if (!bMDPformat)
1047         {
1048             gmx::TextWriter writer(fp);
1049             writer.wrapperSettings().setIndent(indent);
1050             gmx::dumpKeyValueTree(&writer, *ir->params);
1051         }
1052
1053         pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
1054     }
1055 }
1056 #undef PS
1057 #undef PR
1058 #undef PI
1059
1060 static void cmp_grpopts(FILE* fp, const t_grpopts* opt1, const t_grpopts* opt2, real ftol, real abstol)
1061 {
1062     int  i, j;
1063     char buf1[256], buf2[256];
1064
1065     cmp_int(fp, "inputrec->grpopts.ngtc", -1, opt1->ngtc, opt2->ngtc);
1066     cmp_int(fp, "inputrec->grpopts.ngfrz", -1, opt1->ngfrz, opt2->ngfrz);
1067     cmp_int(fp, "inputrec->grpopts.ngener", -1, opt1->ngener, opt2->ngener);
1068     for (i = 0; (i < std::min(opt1->ngtc, opt2->ngtc)); i++)
1069     {
1070         cmp_real(fp, "inputrec->grpopts.nrdf", i, opt1->nrdf[i], opt2->nrdf[i], ftol, abstol);
1071         cmp_real(fp, "inputrec->grpopts.ref_t", i, opt1->ref_t[i], opt2->ref_t[i], ftol, abstol);
1072         cmp_real(fp, "inputrec->grpopts.tau_t", i, opt1->tau_t[i], opt2->tau_t[i], ftol, abstol);
1073         cmpEnum(fp, "inputrec->grpopts.annealing", opt1->annealing[i], opt2->annealing[i]);
1074         cmp_int(fp, "inputrec->grpopts.anneal_npoints", i, opt1->anneal_npoints[i], opt2->anneal_npoints[i]);
1075         if (opt1->anneal_npoints[i] == opt2->anneal_npoints[i])
1076         {
1077             sprintf(buf1, "inputrec->grpopts.anneal_time[%d]", i);
1078             sprintf(buf2, "inputrec->grpopts.anneal_temp[%d]", i);
1079             for (j = 0; j < opt1->anneal_npoints[i]; j++)
1080             {
1081                 cmp_real(fp, buf1, j, opt1->anneal_time[i][j], opt2->anneal_time[i][j], ftol, abstol);
1082                 cmp_real(fp, buf2, j, opt1->anneal_temp[i][j], opt2->anneal_temp[i][j], ftol, abstol);
1083             }
1084         }
1085     }
1086     if (opt1->ngener == opt2->ngener)
1087     {
1088         for (i = 0; i < opt1->ngener; i++)
1089         {
1090             for (j = i; j < opt1->ngener; j++)
1091             {
1092                 sprintf(buf1, "inputrec->grpopts.egp_flags[%d]", i);
1093                 cmp_int(fp, buf1, j, opt1->egp_flags[opt1->ngener * i + j], opt2->egp_flags[opt1->ngener * i + j]);
1094             }
1095         }
1096     }
1097     for (i = 0; (i < std::min(opt1->ngfrz, opt2->ngfrz)); i++)
1098     {
1099         cmp_ivec(fp, "inputrec->grpopts.nFreeze", i, opt1->nFreeze[i], opt2->nFreeze[i]);
1100     }
1101 }
1102
1103 static void cmp_pull(FILE* fp)
1104 {
1105     fprintf(fp,
1106             "WARNING: Both files use COM pulling, but comparing of the pull struct is not "
1107             "implemented (yet). The pull parameters could be the same or different.\n");
1108 }
1109
1110 static void cmp_awhDimParams(FILE*                    fp,
1111                              const gmx::AwhDimParams& dimp1,
1112                              const gmx::AwhDimParams& dimp2,
1113                              int                      dimIndex,
1114                              real                     ftol,
1115                              real                     abstol)
1116 {
1117     /* Note that we have double index here, but the compare functions only
1118      * support one index, so here we only print the dim index and not the bias.
1119      */
1120     cmp_int(fp,
1121             "inputrec.awhParams->bias?->dim->coord_index",
1122             dimIndex,
1123             dimp1.coordinateIndex(),
1124             dimp2.coordinateIndex());
1125     cmp_double(fp, "inputrec->awhParams->bias?->dim->period", dimIndex, dimp1.period(), dimp2.period(), ftol, abstol);
1126     cmp_double(fp,
1127                "inputrec->awhParams->bias?->dim->diffusion",
1128                dimIndex,
1129                dimp1.diffusion(),
1130                dimp2.diffusion(),
1131                ftol,
1132                abstol);
1133     cmp_double(fp, "inputrec->awhParams->bias?->dim->origin", dimIndex, dimp1.origin(), dimp2.origin(), ftol, abstol);
1134     cmp_double(fp, "inputrec->awhParams->bias?->dim->end", dimIndex, dimp1.end(), dimp2.end(), ftol, abstol);
1135     cmp_double(fp,
1136                "inputrec->awhParams->bias?->dim->coord_value_init",
1137                dimIndex,
1138                dimp1.initialCoordinate(),
1139                dimp2.initialCoordinate(),
1140                ftol,
1141                abstol);
1142     cmp_double(fp,
1143                "inputrec->awhParams->bias?->dim->coverDiameter",
1144                dimIndex,
1145                dimp1.coverDiameter(),
1146                dimp2.coverDiameter(),
1147                ftol,
1148                abstol);
1149 }
1150
1151 static void cmp_awhBiasParams(FILE*                     fp,
1152                               const gmx::AwhBiasParams& bias1,
1153                               const gmx::AwhBiasParams& bias2,
1154                               int                       biasIndex,
1155                               real                      ftol,
1156                               real                      abstol)
1157 {
1158     cmp_int(fp, "inputrec->awhParams->ndim", biasIndex, bias1.ndim(), bias2.ndim());
1159     cmpEnum<gmx::AwhTargetType>(
1160             fp, "inputrec->awhParams->biaseTarget", bias1.targetDistribution(), bias2.targetDistribution());
1161     cmp_double(fp,
1162                "inputrec->awhParams->biastargetBetaScaling",
1163                biasIndex,
1164                bias1.targetBetaScaling(),
1165                bias2.targetBetaScaling(),
1166                ftol,
1167                abstol);
1168     cmp_double(fp,
1169                "inputrec->awhParams->biastargetCutoff",
1170                biasIndex,
1171                bias1.targetCutoff(),
1172                bias2.targetCutoff(),
1173                ftol,
1174                abstol);
1175     cmpEnum<gmx::AwhHistogramGrowthType>(
1176             fp, "inputrec->awhParams->biaseGrowth", bias1.growthType(), bias2.growthType());
1177     cmp_bool(fp, "inputrec->awhParams->biasbUserData", biasIndex, bias1.userPMFEstimate(), bias2.userPMFEstimate());
1178     cmp_double(fp,
1179                "inputrec->awhParams->biaserror_initial",
1180                biasIndex,
1181                bias1.initialErrorEstimate(),
1182                bias2.initialErrorEstimate(),
1183                ftol,
1184                abstol);
1185     cmp_int(fp, "inputrec->awhParams->biasShareGroup", biasIndex, bias1.shareGroup(), bias2.shareGroup());
1186
1187     const auto dimParams1 = bias1.dimParams();
1188     const auto dimParams2 = bias2.dimParams();
1189     for (int dim = 0; dim < std::min(bias1.ndim(), bias2.ndim()); dim++)
1190     {
1191         cmp_awhDimParams(fp, dimParams1[dim], dimParams2[dim], dim, ftol, abstol);
1192     }
1193 }
1194
1195 static void cmp_awhParams(FILE* fp, const gmx::AwhParams& awh1, const gmx::AwhParams& awh2, real ftol, real abstol)
1196 {
1197     cmp_int(fp, "inputrec->awhParams->nbias", -1, awh1.numBias(), awh2.numBias());
1198     cmp_int64(fp, "inputrec->awhParams->seed", awh1.seed(), awh2.seed());
1199     cmp_int(fp, "inputrec->awhParams->nstout", -1, awh1.nstout(), awh2.nstout());
1200     cmp_int(fp, "inputrec->awhParams->nstsample_coord", -1, awh1.nstSampleCoord(), awh2.nstSampleCoord());
1201     cmp_int(fp,
1202             "inputrec->awhParams->nsamples_update_free_energy",
1203             -1,
1204             awh1.numSamplesUpdateFreeEnergy(),
1205             awh2.numSamplesUpdateFreeEnergy());
1206     cmpEnum<gmx::AwhPotentialType>(
1207             fp, "inputrec->awhParams->ePotential", awh1.potential(), awh2.potential());
1208     cmp_bool(fp, "inputrec->awhParams->shareBiasMultisim", -1, awh1.shareBiasMultisim(), awh2.shareBiasMultisim());
1209
1210     if (awh1.numBias() == awh2.numBias())
1211     {
1212         const auto awhBiasParams1 = awh1.awhBiasParams();
1213         const auto awhBiasParams2 = awh2.awhBiasParams();
1214         for (int bias = 0; bias < awh1.numBias(); bias++)
1215         {
1216             cmp_awhBiasParams(fp, awhBiasParams1[bias], awhBiasParams2[bias], bias, ftol, abstol);
1217         }
1218     }
1219 }
1220
1221 static void cmp_simtempvals(FILE*            fp,
1222                             const t_simtemp* simtemp1,
1223                             const t_simtemp* simtemp2,
1224                             int              n_lambda,
1225                             real             ftol,
1226                             real             abstol)
1227 {
1228     int i;
1229     cmpEnum(fp, "inputrec->simtempvals->eSimTempScale", simtemp1->eSimTempScale, simtemp2->eSimTempScale);
1230     cmp_real(fp, "inputrec->simtempvals->simtemp_high", -1, simtemp1->simtemp_high, simtemp2->simtemp_high, ftol, abstol);
1231     cmp_real(fp, "inputrec->simtempvals->simtemp_low", -1, simtemp1->simtemp_low, simtemp2->simtemp_low, ftol, abstol);
1232     for (i = 0; i < n_lambda; i++)
1233     {
1234         cmp_real(fp,
1235                  "inputrec->simtempvals->temperatures",
1236                  -1,
1237                  simtemp1->temperatures[i],
1238                  simtemp2->temperatures[i],
1239                  ftol,
1240                  abstol);
1241     }
1242 }
1243
1244 static void cmp_expandedvals(FILE*             fp,
1245                              const t_expanded* expand1,
1246                              const t_expanded* expand2,
1247                              int               n_lambda,
1248                              real              ftol,
1249                              real              abstol)
1250 {
1251     int i;
1252
1253     cmp_bool(fp, "inputrec->fepvals->bInit_weights", -1, expand1->bInit_weights, expand2->bInit_weights);
1254     cmp_bool(fp, "inputrec->fepvals->bWLoneovert", -1, expand1->bWLoneovert, expand2->bWLoneovert);
1255
1256     for (i = 0; i < n_lambda; i++)
1257     {
1258         cmp_real(fp,
1259                  "inputrec->expandedvals->init_lambda_weights",
1260                  -1,
1261                  expand1->init_lambda_weights[i],
1262                  expand2->init_lambda_weights[i],
1263                  ftol,
1264                  abstol);
1265     }
1266
1267     cmpEnum(fp, "inputrec->expandedvals->lambda-stats", expand1->elamstats, expand2->elamstats);
1268     cmpEnum(fp, "inputrec->expandedvals->lambda-mc-move", expand1->elmcmove, expand2->elmcmove);
1269     cmp_int(fp, "inputrec->expandedvals->lmc-repeats", -1, expand1->lmc_repeats, expand2->lmc_repeats);
1270     cmp_int(fp, "inputrec->expandedvals->lmc-gibbsdelta", -1, expand1->gibbsdeltalam, expand2->gibbsdeltalam);
1271     cmp_int(fp, "inputrec->expandedvals->lmc-forced-nstart", -1, expand1->lmc_forced_nstart, expand2->lmc_forced_nstart);
1272     cmpEnum(fp, "inputrec->expandedvals->lambda-weights-equil", expand1->elmceq, expand2->elmceq);
1273     cmp_int(fp,
1274             "inputrec->expandedvals->,weight-equil-number-all-lambda",
1275             -1,
1276             expand1->equil_n_at_lam,
1277             expand2->equil_n_at_lam);
1278     cmp_int(fp, "inputrec->expandedvals->weight-equil-number-samples", -1, expand1->equil_samples, expand2->equil_samples);
1279     cmp_int(fp, "inputrec->expandedvals->weight-equil-number-steps", -1, expand1->equil_steps, expand2->equil_steps);
1280     cmp_real(fp,
1281              "inputrec->expandedvals->weight-equil-wl-delta",
1282              -1,
1283              expand1->equil_wl_delta,
1284              expand2->equil_wl_delta,
1285              ftol,
1286              abstol);
1287     cmp_real(fp,
1288              "inputrec->expandedvals->weight-equil-count-ratio",
1289              -1,
1290              expand1->equil_ratio,
1291              expand2->equil_ratio,
1292              ftol,
1293              abstol);
1294     cmp_bool(fp,
1295              "inputrec->expandedvals->symmetrized-transition-matrix",
1296              -1,
1297              expand1->bSymmetrizedTMatrix,
1298              expand2->bSymmetrizedTMatrix);
1299     cmp_int(fp, "inputrec->expandedvals->nstTij", -1, expand1->nstTij, expand2->nstTij);
1300     cmp_int(fp, "inputrec->expandedvals->mininum-var-min", -1, expand1->minvarmin, expand2->minvarmin); /*default is reasonable */
1301     cmp_int(fp, "inputrec->expandedvals->weight-c-range", -1, expand1->c_range, expand2->c_range); /* default is just C=0 */
1302     cmp_real(fp, "inputrec->expandedvals->wl-scale", -1, expand1->wl_scale, expand2->wl_scale, ftol, abstol);
1303     cmp_real(fp, "inputrec->expandedvals->init-wl-delta", -1, expand1->init_wl_delta, expand2->init_wl_delta, ftol, abstol);
1304     cmp_real(fp, "inputrec->expandedvals->wl-ratio", -1, expand1->wl_ratio, expand2->wl_ratio, ftol, abstol);
1305     cmp_int(fp, "inputrec->expandedvals->nstexpanded", -1, expand1->nstexpanded, expand2->nstexpanded);
1306     cmp_int(fp, "inputrec->expandedvals->lmc-seed", -1, expand1->lmc_seed, expand2->lmc_seed);
1307     cmp_real(fp, "inputrec->expandedvals->mc-temperature", -1, expand1->mc_temp, expand2->mc_temp, ftol, abstol);
1308 }
1309
1310 static void cmp_fepvals(FILE* fp, const t_lambda* fep1, const t_lambda* fep2, real ftol, real abstol)
1311 {
1312     int i, j;
1313     cmp_int(fp, "inputrec->nstdhdl", -1, fep1->nstdhdl, fep2->nstdhdl);
1314     cmp_double(fp, "inputrec->fepvals->init_fep_state", -1, fep1->init_fep_state, fep2->init_fep_state, ftol, abstol);
1315     cmp_double(fp, "inputrec->fepvals->delta_lambda", -1, fep1->delta_lambda, fep2->delta_lambda, ftol, abstol);
1316     cmp_int(fp, "inputrec->fepvals->n_lambda", -1, fep1->n_lambda, fep2->n_lambda);
1317     for (i = 0; i < static_cast<int>(FreeEnergyPerturbationCouplingType::Count); i++)
1318     {
1319         for (j = 0; j < std::min(fep1->n_lambda, fep2->n_lambda); j++)
1320         {
1321             cmp_double(fp,
1322                        "inputrec->fepvals->all_lambda",
1323                        -1,
1324                        fep1->all_lambda[i][j],
1325                        fep2->all_lambda[i][j],
1326                        ftol,
1327                        abstol);
1328         }
1329     }
1330     cmp_int(fp, "inputrec->fepvals->lambda_neighbors", 1, fep1->lambda_neighbors, fep2->lambda_neighbors);
1331     cmp_real(fp, "inputrec->fepvals->sc_alpha", -1, fep1->sc_alpha, fep2->sc_alpha, ftol, abstol);
1332     cmp_int(fp, "inputrec->fepvals->sc_power", -1, fep1->sc_power, fep2->sc_power);
1333     cmp_real(fp, "inputrec->fepvals->sc_r_power", -1, fep1->sc_r_power, fep2->sc_r_power, ftol, abstol);
1334     cmp_real(fp, "inputrec->fepvals->sc_sigma", -1, fep1->sc_sigma, fep2->sc_sigma, ftol, abstol);
1335     cmpEnum(fp, "inputrec->fepvals->edHdLPrintEnergy", fep1->edHdLPrintEnergy, fep1->edHdLPrintEnergy);
1336     cmp_bool(fp, "inputrec->fepvals->bScCoul", -1, fep1->bScCoul, fep1->bScCoul);
1337     cmpEnum(fp, "inputrec->separate_dhdl_file", fep1->separate_dhdl_file, fep2->separate_dhdl_file);
1338     cmpEnum(fp, "inputrec->dhdl_derivatives", fep1->dhdl_derivatives, fep2->dhdl_derivatives);
1339     cmp_int(fp, "inputrec->dh_hist_size", -1, fep1->dh_hist_size, fep2->dh_hist_size);
1340     cmp_double(fp, "inputrec->dh_hist_spacing", -1, fep1->dh_hist_spacing, fep2->dh_hist_spacing, ftol, abstol);
1341 }
1342
1343 void cmp_inputrec(FILE* fp, const t_inputrec* ir1, const t_inputrec* ir2, real ftol, real abstol)
1344 {
1345     fprintf(fp, "comparing inputrec\n");
1346
1347     /* gcc 2.96 doesnt like these defines at all, but issues a huge list
1348      * of warnings. Maybe it will change in future versions, but for the
1349      * moment I've spelled them out instead. /EL 000820
1350      * #define CIB(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
1351      * #define CII(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
1352      * #define CIR(s) cmp_real(fp,"inputrec->"#s,0,ir1->##s,ir2->##s,ftol)
1353      */
1354     cmpEnum(fp, "inputrec->eI", ir1->eI, ir2->eI);
1355     cmp_int64(fp, "inputrec->nsteps", ir1->nsteps, ir2->nsteps);
1356     cmp_int64(fp, "inputrec->init_step", ir1->init_step, ir2->init_step);
1357     cmp_int(fp, "inputrec->simulation_part", -1, ir1->simulation_part, ir2->simulation_part);
1358     cmp_int(fp, "inputrec->mts", -1, static_cast<int>(ir1->useMts), static_cast<int>(ir2->useMts));
1359     if (ir1->useMts && ir2->useMts)
1360     {
1361         cmp_int(fp, "inputrec->mts-levels", -1, ir1->mtsLevels.size(), ir2->mtsLevels.size());
1362         cmp_int(fp,
1363                 "inputrec->mts-level2-forces",
1364                 -1,
1365                 ir1->mtsLevels[1].forceGroups.to_ulong(),
1366                 ir2->mtsLevels[1].forceGroups.to_ulong());
1367         cmp_int(fp,
1368                 "inputrec->mts-level2-factor",
1369                 -1,
1370                 ir1->mtsLevels[1].stepFactor,
1371                 ir2->mtsLevels[1].stepFactor);
1372     }
1373     cmp_int(fp, "inputrec->pbcType", -1, static_cast<int>(ir1->pbcType), static_cast<int>(ir2->pbcType));
1374     cmp_bool(fp, "inputrec->bPeriodicMols", -1, ir1->bPeriodicMols, ir2->bPeriodicMols);
1375     cmpEnum(fp, "inputrec->cutoff_scheme", ir1->cutoff_scheme, ir2->cutoff_scheme);
1376     cmp_int(fp, "inputrec->nstlist", -1, ir1->nstlist, ir2->nstlist);
1377     cmp_int(fp, "inputrec->nstcomm", -1, ir1->nstcomm, ir2->nstcomm);
1378     cmpEnum(fp, "inputrec->comm_mode", ir1->comm_mode, ir2->comm_mode);
1379     cmp_int(fp, "inputrec->nstlog", -1, ir1->nstlog, ir2->nstlog);
1380     cmp_int(fp, "inputrec->nstxout", -1, ir1->nstxout, ir2->nstxout);
1381     cmp_int(fp, "inputrec->nstvout", -1, ir1->nstvout, ir2->nstvout);
1382     cmp_int(fp, "inputrec->nstfout", -1, ir1->nstfout, ir2->nstfout);
1383     cmp_int(fp, "inputrec->nstcalcenergy", -1, ir1->nstcalcenergy, ir2->nstcalcenergy);
1384     cmp_int(fp, "inputrec->nstenergy", -1, ir1->nstenergy, ir2->nstenergy);
1385     cmp_int(fp, "inputrec->nstxout_compressed", -1, ir1->nstxout_compressed, ir2->nstxout_compressed);
1386     cmp_double(fp, "inputrec->init_t", -1, ir1->init_t, ir2->init_t, ftol, abstol);
1387     cmp_double(fp, "inputrec->delta_t", -1, ir1->delta_t, ir2->delta_t, ftol, abstol);
1388     cmp_real(fp,
1389              "inputrec->x_compression_precision",
1390              -1,
1391              ir1->x_compression_precision,
1392              ir2->x_compression_precision,
1393              ftol,
1394              abstol);
1395     cmp_real(fp, "inputrec->fourierspacing", -1, ir1->fourier_spacing, ir2->fourier_spacing, ftol, abstol);
1396     cmp_int(fp, "inputrec->nkx", -1, ir1->nkx, ir2->nkx);
1397     cmp_int(fp, "inputrec->nky", -1, ir1->nky, ir2->nky);
1398     cmp_int(fp, "inputrec->nkz", -1, ir1->nkz, ir2->nkz);
1399     cmp_int(fp, "inputrec->pme_order", -1, ir1->pme_order, ir2->pme_order);
1400     cmp_real(fp, "inputrec->ewald_rtol", -1, ir1->ewald_rtol, ir2->ewald_rtol, ftol, abstol);
1401     cmpEnum(fp, "inputrec->ewald_geometry", ir1->ewald_geometry, ir2->ewald_geometry);
1402     cmp_real(fp, "inputrec->epsilon_surface", -1, ir1->epsilon_surface, ir2->epsilon_surface, ftol, abstol);
1403     cmp_int(fp,
1404             "inputrec->bContinuation",
1405             -1,
1406             static_cast<int>(ir1->bContinuation),
1407             static_cast<int>(ir2->bContinuation));
1408     cmp_int(fp, "inputrec->bShakeSOR", -1, static_cast<int>(ir1->bShakeSOR), static_cast<int>(ir2->bShakeSOR));
1409     cmpEnum(fp, "inputrec->etc", ir1->etc, ir2->etc);
1410     cmp_int(fp,
1411             "inputrec->bPrintNHChains",
1412             -1,
1413             static_cast<int>(ir1->bPrintNHChains),
1414             static_cast<int>(ir2->bPrintNHChains));
1415     cmpEnum(fp, "inputrec->epc", ir1->epc, ir2->epc);
1416     cmpEnum(fp, "inputrec->epct", ir1->epct, ir2->epct);
1417     cmp_real(fp, "inputrec->tau_p", -1, ir1->tau_p, ir2->tau_p, ftol, abstol);
1418     cmp_rvec(fp, "inputrec->ref_p(x)", -1, ir1->ref_p[XX], ir2->ref_p[XX], ftol, abstol);
1419     cmp_rvec(fp, "inputrec->ref_p(y)", -1, ir1->ref_p[YY], ir2->ref_p[YY], ftol, abstol);
1420     cmp_rvec(fp, "inputrec->ref_p(z)", -1, ir1->ref_p[ZZ], ir2->ref_p[ZZ], ftol, abstol);
1421     cmp_rvec(fp, "inputrec->compress(x)", -1, ir1->compress[XX], ir2->compress[XX], ftol, abstol);
1422     cmp_rvec(fp, "inputrec->compress(y)", -1, ir1->compress[YY], ir2->compress[YY], ftol, abstol);
1423     cmp_rvec(fp, "inputrec->compress(z)", -1, ir1->compress[ZZ], ir2->compress[ZZ], ftol, abstol);
1424     cmpEnum(fp, "refcoord_scaling", ir1->refcoord_scaling, ir2->refcoord_scaling);
1425     cmp_rvec(fp, "inputrec->posres_com", -1, ir1->posres_com, ir2->posres_com, ftol, abstol);
1426     cmp_rvec(fp, "inputrec->posres_comB", -1, ir1->posres_comB, ir2->posres_comB, ftol, abstol);
1427     cmp_real(fp, "inputrec->verletbuf_tol", -1, ir1->verletbuf_tol, ir2->verletbuf_tol, ftol, abstol);
1428     cmp_real(fp, "inputrec->rlist", -1, ir1->rlist, ir2->rlist, ftol, abstol);
1429     cmp_real(fp, "inputrec->rtpi", -1, ir1->rtpi, ir2->rtpi, ftol, abstol);
1430     cmpEnum(fp, "inputrec->coulombtype", ir1->coulombtype, ir2->coulombtype);
1431     cmpEnum(fp, "inputrec->coulomb_modifier", ir1->coulomb_modifier, ir2->coulomb_modifier);
1432     cmp_real(fp, "inputrec->rcoulomb_switch", -1, ir1->rcoulomb_switch, ir2->rcoulomb_switch, ftol, abstol);
1433     cmp_real(fp, "inputrec->rcoulomb", -1, ir1->rcoulomb, ir2->rcoulomb, ftol, abstol);
1434     cmpEnum(fp, "inputrec->vdwtype", ir1->vdwtype, ir2->vdwtype);
1435     cmpEnum(fp, "inputrec->vdw_modifier", ir1->vdw_modifier, ir2->vdw_modifier);
1436     cmp_real(fp, "inputrec->rvdw_switch", -1, ir1->rvdw_switch, ir2->rvdw_switch, ftol, abstol);
1437     cmp_real(fp, "inputrec->rvdw", -1, ir1->rvdw, ir2->rvdw, ftol, abstol);
1438     cmp_real(fp, "inputrec->epsilon_r", -1, ir1->epsilon_r, ir2->epsilon_r, ftol, abstol);
1439     cmp_real(fp, "inputrec->epsilon_rf", -1, ir1->epsilon_rf, ir2->epsilon_rf, ftol, abstol);
1440     cmp_real(fp, "inputrec->tabext", -1, ir1->tabext, ir2->tabext, ftol, abstol);
1441
1442     cmpEnum(fp, "inputrec->eDispCorr", ir1->eDispCorr, ir2->eDispCorr);
1443     cmp_real(fp, "inputrec->shake_tol", -1, ir1->shake_tol, ir2->shake_tol, ftol, abstol);
1444     cmpEnum(fp, "inputrec->efep", ir1->efep, ir2->efep);
1445     cmp_fepvals(fp, ir1->fepvals.get(), ir2->fepvals.get(), ftol, abstol);
1446     cmp_int(fp, "inputrec->bSimTemp", -1, static_cast<int>(ir1->bSimTemp), static_cast<int>(ir2->bSimTemp));
1447     if ((ir1->bSimTemp == ir2->bSimTemp) && (ir1->bSimTemp))
1448     {
1449         cmp_simtempvals(fp,
1450                         ir1->simtempvals.get(),
1451                         ir2->simtempvals.get(),
1452                         std::min(ir1->fepvals->n_lambda, ir2->fepvals->n_lambda),
1453                         ftol,
1454                         abstol);
1455     }
1456     cmp_int(fp, "inputrec->bExpanded", -1, static_cast<int>(ir1->bExpanded), static_cast<int>(ir2->bExpanded));
1457     if ((ir1->bExpanded == ir2->bExpanded) && (ir1->bExpanded))
1458     {
1459         cmp_expandedvals(fp,
1460                          ir1->expandedvals.get(),
1461                          ir2->expandedvals.get(),
1462                          std::min(ir1->fepvals->n_lambda, ir2->fepvals->n_lambda),
1463                          ftol,
1464                          abstol);
1465     }
1466     cmp_int(fp, "inputrec->nwall", -1, ir1->nwall, ir2->nwall);
1467     cmpEnum(fp, "inputrec->wall_type", ir1->wall_type, ir2->wall_type);
1468     cmp_int(fp, "inputrec->wall_atomtype[0]", -1, ir1->wall_atomtype[0], ir2->wall_atomtype[0]);
1469     cmp_int(fp, "inputrec->wall_atomtype[1]", -1, ir1->wall_atomtype[1], ir2->wall_atomtype[1]);
1470     cmp_real(fp, "inputrec->wall_density[0]", -1, ir1->wall_density[0], ir2->wall_density[0], ftol, abstol);
1471     cmp_real(fp, "inputrec->wall_density[1]", -1, ir1->wall_density[1], ir2->wall_density[1], ftol, abstol);
1472     cmp_real(fp, "inputrec->wall_ewald_zfac", -1, ir1->wall_ewald_zfac, ir2->wall_ewald_zfac, ftol, abstol);
1473
1474     cmp_bool(fp, "inputrec->bPull", -1, ir1->bPull, ir2->bPull);
1475     if (ir1->bPull && ir2->bPull)
1476     {
1477         cmp_pull(fp);
1478     }
1479
1480     cmp_bool(fp, "inputrec->bDoAwh", -1, ir1->bDoAwh, ir2->bDoAwh);
1481     if (ir1->bDoAwh && ir2->bDoAwh)
1482     {
1483         cmp_awhParams(fp, *ir1->awhParams, *ir2->awhParams, ftol, abstol);
1484     }
1485
1486     cmpEnum(fp, "inputrec->eDisre", ir1->eDisre, ir2->eDisre);
1487     cmp_real(fp, "inputrec->dr_fc", -1, ir1->dr_fc, ir2->dr_fc, ftol, abstol);
1488     cmpEnum(fp, "inputrec->eDisreWeighting", ir1->eDisreWeighting, ir2->eDisreWeighting);
1489     cmp_int(fp, "inputrec->bDisreMixed", -1, static_cast<int>(ir1->bDisreMixed), static_cast<int>(ir2->bDisreMixed));
1490     cmp_int(fp, "inputrec->nstdisreout", -1, ir1->nstdisreout, ir2->nstdisreout);
1491     cmp_real(fp, "inputrec->dr_tau", -1, ir1->dr_tau, ir2->dr_tau, ftol, abstol);
1492     cmp_real(fp, "inputrec->orires_fc", -1, ir1->orires_fc, ir2->orires_fc, ftol, abstol);
1493     cmp_real(fp, "inputrec->orires_tau", -1, ir1->orires_tau, ir2->orires_tau, ftol, abstol);
1494     cmp_int(fp, "inputrec->nstorireout", -1, ir1->nstorireout, ir2->nstorireout);
1495     cmp_real(fp, "inputrec->em_stepsize", -1, ir1->em_stepsize, ir2->em_stepsize, ftol, abstol);
1496     cmp_real(fp, "inputrec->em_tol", -1, ir1->em_tol, ir2->em_tol, ftol, abstol);
1497     cmp_int(fp, "inputrec->niter", -1, ir1->niter, ir2->niter);
1498     cmp_real(fp, "inputrec->fc_stepsize", -1, ir1->fc_stepsize, ir2->fc_stepsize, ftol, abstol);
1499     cmp_int(fp, "inputrec->nstcgsteep", -1, ir1->nstcgsteep, ir2->nstcgsteep);
1500     cmp_int(fp, "inputrec->nbfgscorr", 0, ir1->nbfgscorr, ir2->nbfgscorr);
1501     cmpEnum(fp, "inputrec->eConstrAlg", ir1->eConstrAlg, ir2->eConstrAlg);
1502     cmp_int(fp, "inputrec->nProjOrder", -1, ir1->nProjOrder, ir2->nProjOrder);
1503     cmp_real(fp, "inputrec->LincsWarnAngle", -1, ir1->LincsWarnAngle, ir2->LincsWarnAngle, ftol, abstol);
1504     cmp_int(fp, "inputrec->nLincsIter", -1, ir1->nLincsIter, ir2->nLincsIter);
1505     cmp_real(fp, "inputrec->bd_fric", -1, ir1->bd_fric, ir2->bd_fric, ftol, abstol);
1506     cmp_int64(fp, "inputrec->ld_seed", ir1->ld_seed, ir2->ld_seed);
1507     cmp_real(fp, "inputrec->cos_accel", -1, ir1->cos_accel, ir2->cos_accel, ftol, abstol);
1508     cmp_rvec(fp, "inputrec->deform(a)", -1, ir1->deform[XX], ir2->deform[XX], ftol, abstol);
1509     cmp_rvec(fp, "inputrec->deform(b)", -1, ir1->deform[YY], ir2->deform[YY], ftol, abstol);
1510     cmp_rvec(fp, "inputrec->deform(c)", -1, ir1->deform[ZZ], ir2->deform[ZZ], ftol, abstol);
1511
1512
1513     cmp_int(fp, "inputrec->userint1", -1, ir1->userint1, ir2->userint1);
1514     cmp_int(fp, "inputrec->userint2", -1, ir1->userint2, ir2->userint2);
1515     cmp_int(fp, "inputrec->userint3", -1, ir1->userint3, ir2->userint3);
1516     cmp_int(fp, "inputrec->userint4", -1, ir1->userint4, ir2->userint4);
1517     cmp_real(fp, "inputrec->userreal1", -1, ir1->userreal1, ir2->userreal1, ftol, abstol);
1518     cmp_real(fp, "inputrec->userreal2", -1, ir1->userreal2, ir2->userreal2, ftol, abstol);
1519     cmp_real(fp, "inputrec->userreal3", -1, ir1->userreal3, ir2->userreal3, ftol, abstol);
1520     cmp_real(fp, "inputrec->userreal4", -1, ir1->userreal4, ir2->userreal4, ftol, abstol);
1521     cmp_grpopts(fp, &(ir1->opts), &(ir2->opts), ftol, abstol);
1522     gmx::TextWriter writer(fp);
1523     gmx::compareKeyValueTrees(&writer, *ir1->params, *ir2->params, ftol, abstol);
1524 }
1525
1526 void comp_pull_AB(FILE* fp, const pull_params_t& pull, real ftol, real abstol)
1527 {
1528     for (int i = 0; i < pull.ncoord; i++)
1529     {
1530         fprintf(fp, "comparing pull coord %d\n", i);
1531         cmp_real(fp, "pull-coord->k", -1, pull.coord[i].k, pull.coord[i].kB, ftol, abstol);
1532     }
1533 }
1534
1535 gmx_bool inputrecDeform(const t_inputrec* ir)
1536 {
1537     return (ir->deform[XX][XX] != 0 || ir->deform[YY][YY] != 0 || ir->deform[ZZ][ZZ] != 0
1538             || ir->deform[YY][XX] != 0 || ir->deform[ZZ][XX] != 0 || ir->deform[ZZ][YY] != 0);
1539 }
1540
1541 gmx_bool inputrecDynamicBox(const t_inputrec* ir)
1542 {
1543     return (ir->epc != PressureCoupling::No || ir->eI == IntegrationAlgorithm::TPI || inputrecDeform(ir));
1544 }
1545
1546 gmx_bool inputrecPreserveShape(const t_inputrec* ir)
1547 {
1548     return (ir->epc != PressureCoupling::No && ir->deform[XX][XX] == 0
1549             && (ir->epct == PressureCouplingType::Isotropic
1550                 || ir->epct == PressureCouplingType::SemiIsotropic));
1551 }
1552
1553 gmx_bool inputrecNeedMutot(const t_inputrec* ir)
1554 {
1555     return ((ir->coulombtype == CoulombInteractionType::Ewald || EEL_PME(ir->coulombtype))
1556             && (ir->ewald_geometry == EwaldGeometry::ThreeDC || ir->epsilon_surface != 0));
1557 }
1558
1559 gmx_bool inputrecExclForces(const t_inputrec* ir)
1560 {
1561     return (EEL_FULL(ir->coulombtype) || (EEL_RF(ir->coulombtype)));
1562 }
1563
1564 gmx_bool inputrecNptTrotter(const t_inputrec* ir)
1565 {
1566     return (((ir->eI == IntegrationAlgorithm::VV) || (ir->eI == IntegrationAlgorithm::VVAK))
1567             && (ir->epc == PressureCoupling::Mttk) && (ir->etc == TemperatureCoupling::NoseHoover));
1568 }
1569
1570 gmx_bool inputrecNvtTrotter(const t_inputrec* ir)
1571 {
1572     return (((ir->eI == IntegrationAlgorithm::VV) || (ir->eI == IntegrationAlgorithm::VVAK))
1573             && (ir->epc != PressureCoupling::Mttk) && (ir->etc == TemperatureCoupling::NoseHoover));
1574 }
1575
1576 gmx_bool inputrecNphTrotter(const t_inputrec* ir)
1577 {
1578     return (((ir->eI == IntegrationAlgorithm::VV) || (ir->eI == IntegrationAlgorithm::VVAK))
1579             && (ir->epc == PressureCoupling::Mttk) && (ir->etc != TemperatureCoupling::NoseHoover));
1580 }
1581
1582 bool inputrecPbcXY2Walls(const t_inputrec* ir)
1583 {
1584     return (ir->pbcType == PbcType::XY && ir->nwall == 2);
1585 }
1586
1587 bool integratorHasConservedEnergyQuantity(const t_inputrec* ir)
1588 {
1589     if (!EI_MD(ir->eI))
1590     { // NOLINT bugprone-branch-clone
1591         // Energy minimization or stochastic integrator: no conservation
1592         return false;
1593     }
1594     else if (ir->etc == TemperatureCoupling::No && ir->epc == PressureCoupling::No)
1595     {
1596         // The total energy is conserved, no additional conserved quanitity
1597         return false;
1598     }
1599     else
1600     {
1601         // Shear stress with Parrinello-Rahman is not supported (tedious)
1602         bool shearWithPR =
1603                 ((ir->epc == PressureCoupling::ParrinelloRahman || ir->epc == PressureCoupling::Mttk)
1604                  && (ir->ref_p[YY][XX] != 0 || ir->ref_p[ZZ][XX] != 0 || ir->ref_p[ZZ][YY] != 0));
1605
1606         return !ETC_ANDERSEN(ir->etc) && !shearWithPR;
1607     }
1608 }
1609
1610 bool integratorHasReferenceTemperature(const t_inputrec* ir)
1611 {
1612     return ((ir->etc != TemperatureCoupling::No) || EI_SD(ir->eI)
1613             || (ir->eI == IntegrationAlgorithm::BD) || EI_TPI(ir->eI));
1614 }
1615
1616 int inputrec2nboundeddim(const t_inputrec* ir)
1617 {
1618     if (inputrecPbcXY2Walls(ir))
1619     {
1620         return 3;
1621     }
1622     else
1623     {
1624         return numPbcDimensions(ir->pbcType);
1625     }
1626 }
1627
1628 int ndof_com(const t_inputrec* ir)
1629 {
1630     int n = 0;
1631
1632     switch (ir->pbcType)
1633     {
1634         case PbcType::Xyz:
1635         case PbcType::No: n = 3; break;
1636         case PbcType::XY: n = (ir->nwall == 0 ? 3 : 2); break;
1637         case PbcType::Screw: n = 1; break;
1638         default: gmx_incons("Unknown pbc in calc_nrdf");
1639     }
1640
1641     return n;
1642 }
1643
1644 real maxReferenceTemperature(const t_inputrec& ir)
1645 {
1646     if (EI_ENERGY_MINIMIZATION(ir.eI) || ir.eI == IntegrationAlgorithm::NM)
1647     {
1648         return 0;
1649     }
1650
1651     if (EI_MD(ir.eI) && ir.etc == TemperatureCoupling::No)
1652     {
1653         return -1;
1654     }
1655
1656     /* SD and BD also use ref_t and tau_t for setting the reference temperature.
1657      * TPI can be treated as MD, since it needs an ensemble temperature.
1658      */
1659
1660     real maxTemperature = 0;
1661     for (int i = 0; i < ir.opts.ngtc; i++)
1662     {
1663         if (ir.opts.tau_t[i] >= 0)
1664         {
1665             maxTemperature = std::max(maxTemperature, ir.opts.ref_t[i]);
1666         }
1667     }
1668
1669     return maxTemperature;
1670 }
1671
1672 bool haveEwaldSurfaceContribution(const t_inputrec& ir)
1673 {
1674     return EEL_PME_EWALD(ir.coulombtype)
1675            && (ir.ewald_geometry == EwaldGeometry::ThreeDC || ir.epsilon_surface != 0);
1676 }
1677
1678 bool haveFreeEnergyType(const t_inputrec& ir, const int fepType)
1679 {
1680     for (int i = 0; i < ir.fepvals->n_lambda; i++)
1681     {
1682         if (ir.fepvals->all_lambda[fepType][i] > 0)
1683         {
1684             return true;
1685         }
1686     }
1687     return false;
1688 }