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