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