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