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