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