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