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