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