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