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