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