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