Clean up grompp memory leaks in rotation and pull structures.
[gromacs.git] / src / gromacs / mdtypes / inputrec.cpp
blobd712869ddef8533e431b4f204dc38c1b49bf10d3
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2010, The GROMACS development team.
6 * Copyright (c) 2012,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #include "gmxpre.h"
40 #include "inputrec.h"
42 #include <cstdio>
43 #include <cstdlib>
44 #include <cstring>
46 #include <algorithm>
48 #include "gromacs/math/veccompare.h"
49 #include "gromacs/math/vecdump.h"
50 #include "gromacs/mdtypes/awh_params.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/mdtypes/pull_params.h"
53 #include "gromacs/pbcutil/pbc.h"
54 #include "gromacs/utility/compare.h"
55 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/keyvaluetree.h"
58 #include "gromacs/utility/smalloc.h"
59 #include "gromacs/utility/snprintf.h"
60 #include "gromacs/utility/strconvert.h"
61 #include "gromacs/utility/stringutil.h"
62 #include "gromacs/utility/textwriter.h"
63 #include "gromacs/utility/txtdump.h"
65 //! Macro to select a bool name
66 #define EBOOL(e) gmx::boolToString(e)
68 /* The minimum number of integration steps required for reasonably accurate
69 * integration of first and second order coupling algorithms.
71 const int nstmin_berendsen_tcouple = 5;
72 const int nstmin_berendsen_pcouple = 10;
73 const int nstmin_harmonic = 20;
75 t_inputrec::t_inputrec()
77 // TODO When this memset is removed, remove the suppression of
78 // gcc -Wno-class-memaccess in a CMakeLists.txt file.
79 std::memset(this, 0, sizeof(*this)); // NOLINT(bugprone-undefined-memory-manipulation)
80 snew(fepvals, 1);
81 snew(expandedvals, 1);
82 snew(simtempvals, 1);
85 t_inputrec::~t_inputrec()
87 done_inputrec(this);
90 static int nst_wanted(const t_inputrec* ir)
92 if (ir->nstlist > 0)
94 return ir->nstlist;
96 else
98 return 10;
102 int ir_optimal_nstcalcenergy(const t_inputrec* ir)
104 return nst_wanted(ir);
107 int tcouple_min_integration_steps(int etc)
109 int n;
111 switch (etc)
113 case etcNO: n = 0; break;
114 case etcBERENDSEN:
115 case etcYES: n = nstmin_berendsen_tcouple; break;
116 case etcVRESCALE:
117 /* V-rescale supports instantaneous rescaling */
118 n = 0;
119 break;
120 case etcNOSEHOOVER: n = nstmin_harmonic; break;
121 case etcANDERSEN:
122 case etcANDERSENMASSIVE: n = 1; break;
123 default: gmx_incons("Unknown etc value");
126 return n;
129 int ir_optimal_nsttcouple(const t_inputrec* ir)
131 int nmin, nwanted, n;
132 real tau_min;
133 int g;
135 nmin = tcouple_min_integration_steps(ir->etc);
137 nwanted = nst_wanted(ir);
139 tau_min = 1e20;
140 if (ir->etc != etcNO)
142 for (g = 0; g < ir->opts.ngtc; g++)
144 if (ir->opts.tau_t[g] > 0)
146 tau_min = std::min(tau_min, ir->opts.tau_t[g]);
151 if (nmin == 0 || ir->delta_t * nwanted <= tau_min)
153 n = nwanted;
155 else
157 n = static_cast<int>(tau_min / (ir->delta_t * nmin) + 0.001);
158 if (n < 1)
160 n = 1;
162 while (nwanted % n != 0)
164 n--;
168 return n;
171 int pcouple_min_integration_steps(int epc)
173 int n;
175 switch (epc)
177 case epcNO: n = 0; break;
178 case etcBERENDSEN:
179 case epcISOTROPIC: n = nstmin_berendsen_pcouple; break;
180 case epcPARRINELLORAHMAN:
181 case epcMTTK: n = nstmin_harmonic; break;
182 default: gmx_incons("Unknown epc value");
185 return n;
188 int ir_optimal_nstpcouple(const t_inputrec* ir)
190 int nmin, nwanted, n;
192 nmin = pcouple_min_integration_steps(ir->epc);
194 nwanted = nst_wanted(ir);
196 if (nmin == 0 || ir->delta_t * nwanted <= ir->tau_p)
198 n = nwanted;
200 else
202 n = static_cast<int>(ir->tau_p / (ir->delta_t * nmin) + 0.001);
203 if (n < 1)
205 n = 1;
207 while (nwanted % n != 0)
209 n--;
213 return n;
216 gmx_bool ir_coulomb_switched(const t_inputrec* ir)
218 return (ir->coulombtype == eelSWITCH || ir->coulombtype == eelSHIFT
219 || ir->coulombtype == eelPMESWITCH || ir->coulombtype == eelPMEUSERSWITCH
220 || ir->coulomb_modifier == eintmodPOTSWITCH || ir->coulomb_modifier == eintmodFORCESWITCH);
223 gmx_bool ir_coulomb_is_zero_at_cutoff(const t_inputrec* ir)
225 return (ir->cutoff_scheme == ecutsVERLET || ir_coulomb_switched(ir)
226 || ir->coulomb_modifier != eintmodNONE || ir->coulombtype == eelRF_ZERO);
229 gmx_bool ir_coulomb_might_be_zero_at_cutoff(const t_inputrec* ir)
231 return (ir_coulomb_is_zero_at_cutoff(ir) || ir->coulombtype == eelUSER || ir->coulombtype == eelPMEUSER);
234 gmx_bool ir_vdw_switched(const t_inputrec* ir)
236 return (ir->vdwtype == evdwSWITCH || ir->vdwtype == evdwSHIFT
237 || ir->vdw_modifier == eintmodPOTSWITCH || ir->vdw_modifier == eintmodFORCESWITCH);
240 gmx_bool ir_vdw_is_zero_at_cutoff(const t_inputrec* ir)
242 return (ir->cutoff_scheme == ecutsVERLET || ir_vdw_switched(ir) || ir->vdw_modifier != eintmodNONE);
245 gmx_bool ir_vdw_might_be_zero_at_cutoff(const t_inputrec* ir)
247 return (ir_vdw_is_zero_at_cutoff(ir) || ir->vdwtype == evdwUSER);
250 static void done_pull_group(t_pull_group* pgrp)
252 if (pgrp->nat > 0)
254 sfree(pgrp->ind);
255 sfree(pgrp->weight);
259 static void done_pull_params(pull_params_t* pull)
261 int i;
263 for (i = 0; i < pull->ngroup + 1; i++)
265 done_pull_group(pull->group);
268 sfree(pull->group);
269 sfree(pull->coord);
272 static void done_lambdas(t_lambda* fep)
274 if (fep->n_lambda > 0)
276 for (int i = 0; i < efptNR; i++)
278 sfree(fep->all_lambda[i]);
281 sfree(fep->all_lambda);
284 static void done_t_rot(t_rot* rot)
286 if (rot == nullptr)
288 return;
290 if (rot->grp != nullptr)
292 for (int i = 0; i < rot->ngrp; i++)
294 sfree(rot->grp[i].ind);
295 sfree(rot->grp[i].x_ref);
297 sfree(rot->grp);
299 sfree(rot);
302 void done_inputrec(t_inputrec* ir)
304 sfree(ir->opts.nrdf);
305 sfree(ir->opts.ref_t);
306 for (int i = 0; i < ir->opts.ngtc; i++)
308 sfree(ir->opts.anneal_time[i]);
309 sfree(ir->opts.anneal_temp[i]);
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.egp_flags);
319 done_lambdas(ir->fepvals);
320 sfree(ir->fepvals);
321 sfree(ir->expandedvals);
322 sfree(ir->simtempvals);
324 if (ir->pull)
326 done_pull_params(ir->pull);
327 sfree(ir->pull);
329 done_t_rot(ir->rot);
330 delete ir->params;
333 static void pr_grp_opts(FILE* out, int indent, const char* title, const t_grpopts* opts, gmx_bool bMDPformat)
335 int i, m, j;
337 if (!bMDPformat)
339 fprintf(out, "%s:\n", title);
342 pr_indent(out, indent);
343 fprintf(out, "nrdf%s", bMDPformat ? " = " : ":");
344 for (i = 0; (i < opts->ngtc); i++)
346 fprintf(out, " %10g", opts->nrdf[i]);
348 fprintf(out, "\n");
350 pr_indent(out, indent);
351 fprintf(out, "ref-t%s", bMDPformat ? " = " : ":");
352 for (i = 0; (i < opts->ngtc); i++)
354 fprintf(out, " %10g", opts->ref_t[i]);
356 fprintf(out, "\n");
358 pr_indent(out, indent);
359 fprintf(out, "tau-t%s", bMDPformat ? " = " : ":");
360 for (i = 0; (i < opts->ngtc); i++)
362 fprintf(out, " %10g", opts->tau_t[i]);
364 fprintf(out, "\n");
366 /* Pretty-print the simulated annealing info */
367 fprintf(out, "annealing%s", bMDPformat ? " = " : ":");
368 for (i = 0; (i < opts->ngtc); i++)
370 fprintf(out, " %10s", EANNEAL(opts->annealing[i]));
372 fprintf(out, "\n");
374 fprintf(out, "annealing-npoints%s", bMDPformat ? " = " : ":");
375 for (i = 0; (i < opts->ngtc); i++)
377 fprintf(out, " %10d", opts->anneal_npoints[i]);
379 fprintf(out, "\n");
381 for (i = 0; (i < opts->ngtc); i++)
383 if (opts->anneal_npoints[i] > 0)
385 fprintf(out, "annealing-time [%d]:\t", i);
386 for (j = 0; (j < opts->anneal_npoints[i]); j++)
388 fprintf(out, " %10.1f", opts->anneal_time[i][j]);
390 fprintf(out, "\n");
391 fprintf(out, "annealing-temp [%d]:\t", i);
392 for (j = 0; (j < opts->anneal_npoints[i]); j++)
394 fprintf(out, " %10.1f", opts->anneal_temp[i][j]);
396 fprintf(out, "\n");
400 pr_indent(out, indent);
401 fprintf(out, "acc:\t");
402 for (i = 0; (i < opts->ngacc); i++)
404 for (m = 0; (m < DIM); m++)
406 fprintf(out, " %10g", opts->acc[i][m]);
409 fprintf(out, "\n");
411 pr_indent(out, indent);
412 fprintf(out, "nfreeze:");
413 for (i = 0; (i < opts->ngfrz); i++)
415 for (m = 0; (m < DIM); m++)
417 fprintf(out, " %10s", opts->nFreeze[i][m] ? "Y" : "N");
420 fprintf(out, "\n");
423 for (i = 0; (i < opts->ngener); i++)
425 pr_indent(out, indent);
426 fprintf(out, "energygrp-flags[%3d]:", i);
427 for (m = 0; (m < opts->ngener); m++)
429 fprintf(out, " %d", opts->egp_flags[opts->ngener * i + m]);
431 fprintf(out, "\n");
434 fflush(out);
437 static void pr_matrix(FILE* fp, int indent, const char* title, const rvec* m, gmx_bool bMDPformat)
439 if (bMDPformat)
441 fprintf(fp, "%-10s = %g %g %g %g %g %g\n", title, m[XX][XX], m[YY][YY], m[ZZ][ZZ],
442 m[XX][YY], m[XX][ZZ], m[YY][ZZ]);
444 else
446 pr_rvecs(fp, indent, title, m, DIM);
450 #define PS(t, s) pr_str(fp, indent, t, s)
451 #define PI(t, s) pr_int(fp, indent, t, s)
452 #define PSTEP(t, s) pr_int64(fp, indent, t, s)
453 #define PR(t, s) pr_real(fp, indent, t, s)
454 #define PD(t, s) pr_double(fp, indent, t, s)
456 static void pr_pull_group(FILE* fp, int indent, int g, const t_pull_group* pgrp)
458 pr_indent(fp, indent);
459 fprintf(fp, "pull-group %d:\n", g);
460 indent += 2;
461 pr_ivec_block(fp, indent, "atom", pgrp->ind, pgrp->nat, TRUE);
462 pr_rvec(fp, indent, "weight", pgrp->weight, pgrp->nweight, TRUE);
463 PI("pbcatom", pgrp->pbcatom);
466 static void pr_pull_coord(FILE* fp, int indent, int c, const t_pull_coord* pcrd)
468 int g;
470 pr_indent(fp, indent);
471 fprintf(fp, "pull-coord %d:\n", c);
472 PS("type", EPULLTYPE(pcrd->eType));
473 if (pcrd->eType == epullEXTERNAL)
475 PS("potential-provider", pcrd->externalPotentialProvider);
477 PS("geometry", EPULLGEOM(pcrd->eGeom));
478 for (g = 0; g < pcrd->ngroup; g++)
480 char buf[10];
482 sprintf(buf, "group[%d]", g);
483 PI(buf, pcrd->group[g]);
485 pr_ivec(fp, indent, "dim", pcrd->dim, DIM, TRUE);
486 pr_rvec(fp, indent, "origin", pcrd->origin, DIM, TRUE);
487 pr_rvec(fp, indent, "vec", pcrd->vec, DIM, TRUE);
488 PS("start", EBOOL(pcrd->bStart));
489 PR("init", pcrd->init);
490 PR("rate", pcrd->rate);
491 PR("k", pcrd->k);
492 PR("kB", pcrd->kB);
495 static void pr_simtempvals(FILE* fp, int indent, const t_simtemp* simtemp, int n_lambda)
497 PS("simulated-tempering-scaling", ESIMTEMP(simtemp->eSimTempScale));
498 PR("sim-temp-low", simtemp->simtemp_low);
499 PR("sim-temp-high", simtemp->simtemp_high);
500 pr_rvec(fp, indent, "simulated tempering temperatures", simtemp->temperatures, n_lambda, TRUE);
503 static void pr_expandedvals(FILE* fp, int indent, const t_expanded* expand, int n_lambda)
506 PI("nstexpanded", expand->nstexpanded);
507 PS("lmc-stats", elamstats_names[expand->elamstats]);
508 PS("lmc-move", elmcmove_names[expand->elmcmove]);
509 PS("lmc-weights-equil", elmceq_names[expand->elmceq]);
510 if (expand->elmceq == elmceqNUMATLAM)
512 PI("weight-equil-number-all-lambda", expand->equil_n_at_lam);
514 if (expand->elmceq == elmceqSAMPLES)
516 PI("weight-equil-number-samples", expand->equil_samples);
518 if (expand->elmceq == elmceqSTEPS)
520 PI("weight-equil-number-steps", expand->equil_steps);
522 if (expand->elmceq == elmceqWLDELTA)
524 PR("weight-equil-wl-delta", expand->equil_wl_delta);
526 if (expand->elmceq == elmceqRATIO)
528 PR("weight-equil-count-ratio", expand->equil_ratio);
530 PI("lmc-seed", expand->lmc_seed);
531 PR("mc-temperature", expand->mc_temp);
532 PI("lmc-repeats", expand->lmc_repeats);
533 PI("lmc-gibbsdelta", expand->gibbsdeltalam);
534 PI("lmc-forced-nstart", expand->lmc_forced_nstart);
535 PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
536 PI("nst-transition-matrix", expand->nstTij);
537 PI("mininum-var-min", expand->minvarmin); /*default is reasonable */
538 PI("weight-c-range", expand->c_range); /* default is just C=0 */
539 PR("wl-scale", expand->wl_scale);
540 PR("wl-ratio", expand->wl_ratio);
541 PR("init-wl-delta", expand->init_wl_delta);
542 PS("wl-oneovert", EBOOL(expand->bWLoneovert));
544 pr_indent(fp, indent);
545 pr_rvec(fp, indent, "init-lambda-weights", expand->init_lambda_weights, n_lambda, TRUE);
546 PS("init-weights", EBOOL(expand->bInit_weights));
549 static void pr_fepvals(FILE* fp, int indent, const t_lambda* fep, gmx_bool bMDPformat)
551 int i, j;
553 PR("init-lambda", fep->init_lambda);
554 PI("init-lambda-state", fep->init_fep_state);
555 PR("delta-lambda", fep->delta_lambda);
556 PI("nstdhdl", fep->nstdhdl);
558 if (!bMDPformat)
560 PI("n-lambdas", fep->n_lambda);
562 if (fep->n_lambda > 0)
564 pr_indent(fp, indent);
565 fprintf(fp, "separate-dvdl%s\n", bMDPformat ? " = " : ":");
566 for (i = 0; i < efptNR; i++)
568 fprintf(fp, "%18s = ", efpt_names[i]);
569 if (fep->separate_dvdl[i])
571 fprintf(fp, " TRUE");
573 else
575 fprintf(fp, " FALSE");
577 fprintf(fp, "\n");
579 fprintf(fp, "all-lambdas%s\n", bMDPformat ? " = " : ":");
580 for (i = 0; i < efptNR; i++)
582 fprintf(fp, "%18s = ", efpt_names[i]);
583 for (j = 0; j < fep->n_lambda; j++)
585 fprintf(fp, " %10g", fep->all_lambda[i][j]);
587 fprintf(fp, "\n");
590 PI("calc-lambda-neighbors", fep->lambda_neighbors);
591 PS("dhdl-print-energy", edHdLPrintEnergy_names[fep->edHdLPrintEnergy]);
592 PR("sc-alpha", fep->sc_alpha);
593 PI("sc-power", fep->sc_power);
594 PR("sc-r-power", fep->sc_r_power);
595 PR("sc-sigma", fep->sc_sigma);
596 PR("sc-sigma-min", fep->sc_sigma_min);
597 PS("sc-coul", EBOOL(fep->bScCoul));
598 PI("dh-hist-size", fep->dh_hist_size);
599 PD("dh-hist-spacing", fep->dh_hist_spacing);
600 PS("separate-dhdl-file", SEPDHDLFILETYPE(fep->separate_dhdl_file));
601 PS("dhdl-derivatives", DHDLDERIVATIVESTYPE(fep->dhdl_derivatives));
604 static void pr_pull(FILE* fp, int indent, const pull_params_t* pull)
606 int g;
608 PR("pull-cylinder-r", pull->cylinder_r);
609 PR("pull-constr-tol", pull->constr_tol);
610 PS("pull-print-COM", EBOOL(pull->bPrintCOM));
611 PS("pull-print-ref-value", EBOOL(pull->bPrintRefValue));
612 PS("pull-print-components", EBOOL(pull->bPrintComp));
613 PI("pull-nstxout", pull->nstxout);
614 PI("pull-nstfout", pull->nstfout);
615 PS("pull-pbc-ref-prev-step-com", EBOOL(pull->bSetPbcRefToPrevStepCOM));
616 PS("pull-xout-average", EBOOL(pull->bXOutAverage));
617 PS("pull-fout-average", EBOOL(pull->bFOutAverage));
618 PI("pull-ngroups", pull->ngroup);
619 for (g = 0; g < pull->ngroup; g++)
621 pr_pull_group(fp, indent, g, &pull->group[g]);
623 PI("pull-ncoords", pull->ncoord);
624 for (g = 0; g < pull->ncoord; g++)
626 pr_pull_coord(fp, indent, g, &pull->coord[g]);
630 static void pr_awh_bias_dim(FILE* fp, int indent, gmx::AwhDimParams* awhDimParams, const char* prefix)
632 pr_indent(fp, indent);
633 indent++;
634 fprintf(fp, "%s:\n", prefix);
635 PS("coord-provider", EAWHCOORDPROVIDER(awhDimParams->eCoordProvider));
636 PI("coord-index", awhDimParams->coordIndex + 1);
637 PR("start", awhDimParams->origin);
638 PR("end", awhDimParams->end);
639 PR("period", awhDimParams->period);
640 PR("force-constant", awhDimParams->forceConstant);
641 PR("diffusion", awhDimParams->diffusion);
642 PR("cover-diameter", awhDimParams->coverDiameter);
645 static void pr_awh_bias(FILE* fp, int indent, gmx::AwhBiasParams* awhBiasParams, const char* prefix)
647 char opt[STRLEN];
649 sprintf(opt, "%s-error-init", prefix);
650 PR(opt, awhBiasParams->errorInitial);
651 sprintf(opt, "%s-growth", prefix);
652 PS(opt, EAWHGROWTH(awhBiasParams->eGrowth));
653 sprintf(opt, "%s-target", prefix);
654 PS(opt, EAWHTARGET(awhBiasParams->eTarget));
655 sprintf(opt, "%s-target-beta-scalng", prefix);
656 PR(opt, awhBiasParams->targetBetaScaling);
657 sprintf(opt, "%s-target-cutoff", prefix);
658 PR(opt, awhBiasParams->targetCutoff);
659 sprintf(opt, "%s-user-data", prefix);
660 PS(opt, EBOOL(awhBiasParams->bUserData));
661 sprintf(opt, "%s-share-group", prefix);
662 PI(opt, awhBiasParams->shareGroup);
663 sprintf(opt, "%s-equilibrate-histogram", prefix);
664 PS(opt, EBOOL(awhBiasParams->equilibrateHistogram));
665 sprintf(opt, "%s-ndim", prefix);
666 PI(opt, awhBiasParams->ndim);
668 for (int d = 0; d < awhBiasParams->ndim; d++)
670 char prefixdim[STRLEN];
671 sprintf(prefixdim, "%s-dim%d", prefix, d + 1);
672 pr_awh_bias_dim(fp, indent, &awhBiasParams->dimParams[d], prefixdim);
676 static void pr_awh(FILE* fp, int indent, gmx::AwhParams* awhParams)
678 PS("awh-potential", EAWHPOTENTIAL(awhParams->ePotential));
679 PI("awh-seed", awhParams->seed);
680 PI("awh-nstout", awhParams->nstOut);
681 PI("awh-nstsample", awhParams->nstSampleCoord);
682 PI("awh-nsamples-update", awhParams->numSamplesUpdateFreeEnergy);
683 PS("awh-share-bias-multisim", EBOOL(awhParams->shareBiasMultisim));
684 PI("awh-nbias", awhParams->numBias);
686 for (int k = 0; k < awhParams->numBias; k++)
688 auto prefix = gmx::formatString("awh%d", k + 1);
689 pr_awh_bias(fp, indent, &awhParams->awhBiasParams[k], prefix.c_str());
693 static void pr_rotgrp(FILE* fp, int indent, int g, const t_rotgrp* rotg)
695 pr_indent(fp, indent);
696 fprintf(fp, "rot-group %d:\n", g);
697 indent += 2;
698 PS("rot-type", EROTGEOM(rotg->eType));
699 PS("rot-massw", EBOOL(rotg->bMassW));
700 pr_ivec_block(fp, indent, "atom", rotg->ind, rotg->nat, TRUE);
701 pr_rvecs(fp, indent, "x-ref", rotg->x_ref, rotg->nat);
702 pr_rvec(fp, indent, "rot-vec", rotg->inputVec, DIM, TRUE);
703 pr_rvec(fp, indent, "rot-pivot", rotg->pivot, DIM, TRUE);
704 PR("rot-rate", rotg->rate);
705 PR("rot-k", rotg->k);
706 PR("rot-slab-dist", rotg->slab_dist);
707 PR("rot-min-gauss", rotg->min_gaussian);
708 PR("rot-eps", rotg->eps);
709 PS("rot-fit-method", EROTFIT(rotg->eFittype));
710 PI("rot-potfit-nstep", rotg->PotAngle_nstep);
711 PR("rot-potfit-step", rotg->PotAngle_step);
714 static void pr_rot(FILE* fp, int indent, const t_rot* rot)
716 int g;
718 PI("rot-nstrout", rot->nstrout);
719 PI("rot-nstsout", rot->nstsout);
720 PI("rot-ngroups", rot->ngrp);
721 for (g = 0; g < rot->ngrp; g++)
723 pr_rotgrp(fp, indent, g, &rot->grp[g]);
728 static void pr_swap(FILE* fp, int indent, const t_swapcoords* swap)
730 char str[STRLEN];
732 /* Enums for better readability of the code */
733 enum
735 eCompA = 0,
736 eCompB
740 PI("swap-frequency", swap->nstswap);
742 /* The split groups that define the compartments */
743 for (int j = 0; j < 2; j++)
745 snprintf(str, STRLEN, "massw_split%d", j);
746 PS(str, EBOOL(swap->massw_split[j]));
747 snprintf(str, STRLEN, "split atoms group %d", j);
748 pr_ivec_block(fp, indent, str, swap->grp[j].ind, swap->grp[j].nat, TRUE);
751 /* The solvent group */
752 snprintf(str, STRLEN, "solvent group %s", swap->grp[eGrpSolvent].molname);
753 pr_ivec_block(fp, indent, str, swap->grp[eGrpSolvent].ind, swap->grp[eGrpSolvent].nat, TRUE);
755 /* Now print the indices for all the ion groups: */
756 for (int ig = eSwapFixedGrpNR; ig < swap->ngrp; ig++)
758 snprintf(str, STRLEN, "ion group %s", swap->grp[ig].molname);
759 pr_ivec_block(fp, indent, str, swap->grp[ig].ind, swap->grp[ig].nat, TRUE);
762 PR("cyl0-r", swap->cyl0r);
763 PR("cyl0-up", swap->cyl0u);
764 PR("cyl0-down", swap->cyl0l);
765 PR("cyl1-r", swap->cyl1r);
766 PR("cyl1-up", swap->cyl1u);
767 PR("cyl1-down", swap->cyl1l);
768 PI("coupl-steps", swap->nAverage);
770 /* Print the requested ion counts for both compartments */
771 for (int ic = eCompA; ic <= eCompB; ic++)
773 for (int ig = eSwapFixedGrpNR; ig < swap->ngrp; ig++)
775 snprintf(str, STRLEN, "%s-in-%c", swap->grp[ig].molname, 'A' + ic);
776 PI(str, swap->grp[ig].nmolReq[ic]);
780 PR("threshold", swap->threshold);
781 PR("bulk-offsetA", swap->bulkOffset[eCompA]);
782 PR("bulk-offsetB", swap->bulkOffset[eCompB]);
786 static void pr_imd(FILE* fp, int indent, const t_IMD* imd)
788 PI("IMD-atoms", imd->nat);
789 pr_ivec_block(fp, indent, "atom", imd->ind, imd->nat, TRUE);
793 void pr_inputrec(FILE* fp, int indent, const char* title, const t_inputrec* ir, gmx_bool bMDPformat)
795 const char* infbuf = "inf";
797 if (available(fp, ir, indent, title))
799 if (!bMDPformat)
801 indent = pr_title(fp, indent, title);
803 /* Try to make this list appear in the same order as the
804 * options are written in the default mdout.mdp, and with
805 * the same user-exposed names to facilitate debugging.
807 PS("integrator", EI(ir->eI));
808 PR("tinit", ir->init_t);
809 PR("dt", ir->delta_t);
810 PSTEP("nsteps", ir->nsteps);
811 PSTEP("init-step", ir->init_step);
812 PI("simulation-part", ir->simulation_part);
813 PS("comm-mode", ECOM(ir->comm_mode));
814 PI("nstcomm", ir->nstcomm);
816 /* Langevin dynamics */
817 PR("bd-fric", ir->bd_fric);
818 PSTEP("ld-seed", ir->ld_seed);
820 /* Energy minimization */
821 PR("emtol", ir->em_tol);
822 PR("emstep", ir->em_stepsize);
823 PI("niter", ir->niter);
824 PR("fcstep", ir->fc_stepsize);
825 PI("nstcgsteep", ir->nstcgsteep);
826 PI("nbfgscorr", ir->nbfgscorr);
828 /* Test particle insertion */
829 PR("rtpi", ir->rtpi);
831 /* Output control */
832 PI("nstxout", ir->nstxout);
833 PI("nstvout", ir->nstvout);
834 PI("nstfout", ir->nstfout);
835 PI("nstlog", ir->nstlog);
836 PI("nstcalcenergy", ir->nstcalcenergy);
837 PI("nstenergy", ir->nstenergy);
838 PI("nstxout-compressed", ir->nstxout_compressed);
839 PR("compressed-x-precision", ir->x_compression_precision);
841 /* Neighborsearching parameters */
842 PS("cutoff-scheme", ECUTSCHEME(ir->cutoff_scheme));
843 PI("nstlist", ir->nstlist);
844 PS("pbc", c_pbcTypeNames[ir->pbcType].c_str());
845 PS("periodic-molecules", EBOOL(ir->bPeriodicMols));
846 PR("verlet-buffer-tolerance", ir->verletbuf_tol);
847 PR("rlist", ir->rlist);
849 /* Options for electrostatics and VdW */
850 PS("coulombtype", EELTYPE(ir->coulombtype));
851 PS("coulomb-modifier", INTMODIFIER(ir->coulomb_modifier));
852 PR("rcoulomb-switch", ir->rcoulomb_switch);
853 PR("rcoulomb", ir->rcoulomb);
854 if (ir->epsilon_r != 0)
856 PR("epsilon-r", ir->epsilon_r);
858 else
860 PS("epsilon-r", infbuf);
862 if (ir->epsilon_rf != 0)
864 PR("epsilon-rf", ir->epsilon_rf);
866 else
868 PS("epsilon-rf", infbuf);
870 PS("vdw-type", EVDWTYPE(ir->vdwtype));
871 PS("vdw-modifier", INTMODIFIER(ir->vdw_modifier));
872 PR("rvdw-switch", ir->rvdw_switch);
873 PR("rvdw", ir->rvdw);
874 PS("DispCorr", EDISPCORR(ir->eDispCorr));
875 PR("table-extension", ir->tabext);
877 PR("fourierspacing", ir->fourier_spacing);
878 PI("fourier-nx", ir->nkx);
879 PI("fourier-ny", ir->nky);
880 PI("fourier-nz", ir->nkz);
881 PI("pme-order", ir->pme_order);
882 PR("ewald-rtol", ir->ewald_rtol);
883 PR("ewald-rtol-lj", ir->ewald_rtol_lj);
884 PS("lj-pme-comb-rule", ELJPMECOMBNAMES(ir->ljpme_combination_rule));
885 PR("ewald-geometry", ir->ewald_geometry);
886 PR("epsilon-surface", ir->epsilon_surface);
888 /* Options for weak coupling algorithms */
889 PS("tcoupl", ETCOUPLTYPE(ir->etc));
890 PI("nsttcouple", ir->nsttcouple);
891 PI("nh-chain-length", ir->opts.nhchainlength);
892 PS("print-nose-hoover-chain-variables", EBOOL(ir->bPrintNHChains));
894 PS("pcoupl", EPCOUPLTYPE(ir->epc));
895 PS("pcoupltype", EPCOUPLTYPETYPE(ir->epct));
896 PI("nstpcouple", ir->nstpcouple);
897 PR("tau-p", ir->tau_p);
898 pr_matrix(fp, indent, "compressibility", ir->compress, bMDPformat);
899 pr_matrix(fp, indent, "ref-p", ir->ref_p, bMDPformat);
900 PS("refcoord-scaling", EREFSCALINGTYPE(ir->refcoord_scaling));
902 if (bMDPformat)
904 fprintf(fp, "posres-com = %g %g %g\n", ir->posres_com[XX], ir->posres_com[YY],
905 ir->posres_com[ZZ]);
906 fprintf(fp, "posres-comB = %g %g %g\n", ir->posres_comB[XX], ir->posres_comB[YY],
907 ir->posres_comB[ZZ]);
909 else
911 pr_rvec(fp, indent, "posres-com", ir->posres_com, DIM, TRUE);
912 pr_rvec(fp, indent, "posres-comB", ir->posres_comB, DIM, TRUE);
915 /* QMMM */
916 PS("QMMM", EBOOL(ir->bQMMM));
917 fprintf(fp, "%s:\n", "qm-opts");
918 pr_int(fp, indent, "ngQM", ir->opts.ngQM);
920 /* CONSTRAINT OPTIONS */
921 PS("constraint-algorithm", ECONSTRTYPE(ir->eConstrAlg));
922 PS("continuation", EBOOL(ir->bContinuation));
924 PS("Shake-SOR", EBOOL(ir->bShakeSOR));
925 PR("shake-tol", ir->shake_tol);
926 PI("lincs-order", ir->nProjOrder);
927 PI("lincs-iter", ir->nLincsIter);
928 PR("lincs-warnangle", ir->LincsWarnAngle);
930 /* Walls */
931 PI("nwall", ir->nwall);
932 PS("wall-type", EWALLTYPE(ir->wall_type));
933 PR("wall-r-linpot", ir->wall_r_linpot);
934 /* wall-atomtype */
935 PI("wall-atomtype[0]", ir->wall_atomtype[0]);
936 PI("wall-atomtype[1]", ir->wall_atomtype[1]);
937 /* wall-density */
938 PR("wall-density[0]", ir->wall_density[0]);
939 PR("wall-density[1]", ir->wall_density[1]);
940 PR("wall-ewald-zfac", ir->wall_ewald_zfac);
942 /* COM PULLING */
943 PS("pull", EBOOL(ir->bPull));
944 if (ir->bPull)
946 pr_pull(fp, indent, ir->pull);
949 /* AWH BIASING */
950 PS("awh", EBOOL(ir->bDoAwh));
951 if (ir->bDoAwh)
953 pr_awh(fp, indent, ir->awhParams);
956 /* ENFORCED ROTATION */
957 PS("rotation", EBOOL(ir->bRot));
958 if (ir->bRot)
960 pr_rot(fp, indent, ir->rot);
963 /* INTERACTIVE MD */
964 PS("interactiveMD", EBOOL(ir->bIMD));
965 if (ir->bIMD)
967 pr_imd(fp, indent, ir->imd);
970 /* NMR refinement stuff */
971 PS("disre", EDISRETYPE(ir->eDisre));
972 PS("disre-weighting", EDISREWEIGHTING(ir->eDisreWeighting));
973 PS("disre-mixed", EBOOL(ir->bDisreMixed));
974 PR("dr-fc", ir->dr_fc);
975 PR("dr-tau", ir->dr_tau);
976 PR("nstdisreout", ir->nstdisreout);
978 PR("orire-fc", ir->orires_fc);
979 PR("orire-tau", ir->orires_tau);
980 PR("nstorireout", ir->nstorireout);
982 /* FREE ENERGY VARIABLES */
983 PS("free-energy", EFEPTYPE(ir->efep));
984 if (ir->efep != efepNO || ir->bSimTemp)
986 pr_fepvals(fp, indent, ir->fepvals, bMDPformat);
988 if (ir->bExpanded)
990 pr_expandedvals(fp, indent, ir->expandedvals, ir->fepvals->n_lambda);
993 /* NON-equilibrium MD stuff */
994 PR("cos-acceleration", ir->cos_accel);
995 pr_matrix(fp, indent, "deform", ir->deform, bMDPformat);
997 /* SIMULATED TEMPERING */
998 PS("simulated-tempering", EBOOL(ir->bSimTemp));
999 if (ir->bSimTemp)
1001 pr_simtempvals(fp, indent, ir->simtempvals, ir->fepvals->n_lambda);
1004 /* ION/WATER SWAPPING FOR COMPUTATIONAL ELECTROPHYSIOLOGY */
1005 PS("swapcoords", ESWAPTYPE(ir->eSwapCoords));
1006 if (ir->eSwapCoords != eswapNO)
1008 pr_swap(fp, indent, ir->swap);
1011 /* USER-DEFINED THINGIES */
1012 PI("userint1", ir->userint1);
1013 PI("userint2", ir->userint2);
1014 PI("userint3", ir->userint3);
1015 PI("userint4", ir->userint4);
1016 PR("userreal1", ir->userreal1);
1017 PR("userreal2", ir->userreal2);
1018 PR("userreal3", ir->userreal3);
1019 PR("userreal4", ir->userreal4);
1021 if (!bMDPformat)
1023 gmx::TextWriter writer(fp);
1024 writer.wrapperSettings().setIndent(indent);
1025 gmx::dumpKeyValueTree(&writer, *ir->params);
1028 pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
1031 #undef PS
1032 #undef PR
1033 #undef PI
1035 static void cmp_grpopts(FILE* fp, const t_grpopts* opt1, const t_grpopts* opt2, real ftol, real abstol)
1037 int i, j;
1038 char buf1[256], buf2[256];
1040 cmp_int(fp, "inputrec->grpopts.ngtc", -1, opt1->ngtc, opt2->ngtc);
1041 cmp_int(fp, "inputrec->grpopts.ngacc", -1, opt1->ngacc, opt2->ngacc);
1042 cmp_int(fp, "inputrec->grpopts.ngfrz", -1, opt1->ngfrz, opt2->ngfrz);
1043 cmp_int(fp, "inputrec->grpopts.ngener", -1, opt1->ngener, opt2->ngener);
1044 for (i = 0; (i < std::min(opt1->ngtc, opt2->ngtc)); i++)
1046 cmp_real(fp, "inputrec->grpopts.nrdf", i, opt1->nrdf[i], opt2->nrdf[i], ftol, abstol);
1047 cmp_real(fp, "inputrec->grpopts.ref_t", i, opt1->ref_t[i], opt2->ref_t[i], ftol, abstol);
1048 cmp_real(fp, "inputrec->grpopts.tau_t", i, opt1->tau_t[i], opt2->tau_t[i], ftol, abstol);
1049 cmp_int(fp, "inputrec->grpopts.annealing", i, opt1->annealing[i], opt2->annealing[i]);
1050 cmp_int(fp, "inputrec->grpopts.anneal_npoints", i, opt1->anneal_npoints[i],
1051 opt2->anneal_npoints[i]);
1052 if (opt1->anneal_npoints[i] == opt2->anneal_npoints[i])
1054 sprintf(buf1, "inputrec->grpopts.anneal_time[%d]", i);
1055 sprintf(buf2, "inputrec->grpopts.anneal_temp[%d]", i);
1056 for (j = 0; j < opt1->anneal_npoints[i]; j++)
1058 cmp_real(fp, buf1, j, opt1->anneal_time[i][j], opt2->anneal_time[i][j], ftol, abstol);
1059 cmp_real(fp, buf2, j, opt1->anneal_temp[i][j], opt2->anneal_temp[i][j], ftol, abstol);
1063 if (opt1->ngener == opt2->ngener)
1065 for (i = 0; i < opt1->ngener; i++)
1067 for (j = i; j < opt1->ngener; j++)
1069 sprintf(buf1, "inputrec->grpopts.egp_flags[%d]", i);
1070 cmp_int(fp, buf1, j, opt1->egp_flags[opt1->ngener * i + j],
1071 opt2->egp_flags[opt1->ngener * i + j]);
1075 for (i = 0; (i < std::min(opt1->ngacc, opt2->ngacc)); i++)
1077 cmp_rvec(fp, "inputrec->grpopts.acc", i, opt1->acc[i], opt2->acc[i], ftol, abstol);
1079 for (i = 0; (i < std::min(opt1->ngfrz, opt2->ngfrz)); i++)
1081 cmp_ivec(fp, "inputrec->grpopts.nFreeze", i, opt1->nFreeze[i], opt2->nFreeze[i]);
1085 static void cmp_pull(FILE* fp)
1087 fprintf(fp,
1088 "WARNING: Both files use COM pulling, but comparing of the pull struct is not "
1089 "implemented (yet). The pull parameters could be the same or different.\n");
1092 static void cmp_awhDimParams(FILE* fp,
1093 const gmx::AwhDimParams* dimp1,
1094 const gmx::AwhDimParams* dimp2,
1095 int dimIndex,
1096 real ftol,
1097 real abstol)
1099 /* Note that we have double index here, but the compare functions only
1100 * support one index, so here we only print the dim index and not the bias.
1102 cmp_int(fp, "inputrec.awhParams->bias?->dim->coord_index", dimIndex, dimp1->coordIndex,
1103 dimp2->coordIndex);
1104 cmp_double(fp, "inputrec->awhParams->bias?->dim->period", dimIndex, dimp1->period,
1105 dimp2->period, ftol, abstol);
1106 cmp_double(fp, "inputrec->awhParams->bias?->dim->diffusion", dimIndex, dimp1->diffusion,
1107 dimp2->diffusion, ftol, abstol);
1108 cmp_double(fp, "inputrec->awhParams->bias?->dim->origin", dimIndex, dimp1->origin,
1109 dimp2->origin, ftol, abstol);
1110 cmp_double(fp, "inputrec->awhParams->bias?->dim->end", dimIndex, dimp1->end, dimp2->end, ftol, abstol);
1111 cmp_double(fp, "inputrec->awhParams->bias?->dim->coord_value_init", dimIndex,
1112 dimp1->coordValueInit, dimp2->coordValueInit, ftol, abstol);
1113 cmp_double(fp, "inputrec->awhParams->bias?->dim->coverDiameter", dimIndex, dimp1->coverDiameter,
1114 dimp2->coverDiameter, ftol, abstol);
1117 static void cmp_awhBiasParams(FILE* fp,
1118 const gmx::AwhBiasParams* bias1,
1119 const gmx::AwhBiasParams* bias2,
1120 int biasIndex,
1121 real ftol,
1122 real abstol)
1124 cmp_int(fp, "inputrec->awhParams->ndim", biasIndex, bias1->ndim, bias2->ndim);
1125 cmp_int(fp, "inputrec->awhParams->biaseTarget", biasIndex, bias1->eTarget, bias2->eTarget);
1126 cmp_double(fp, "inputrec->awhParams->biastargetBetaScaling", biasIndex,
1127 bias1->targetBetaScaling, bias2->targetBetaScaling, ftol, abstol);
1128 cmp_double(fp, "inputrec->awhParams->biastargetCutoff", biasIndex, bias1->targetCutoff,
1129 bias2->targetCutoff, ftol, abstol);
1130 cmp_int(fp, "inputrec->awhParams->biaseGrowth", biasIndex, bias1->eGrowth, bias2->eGrowth);
1131 cmp_bool(fp, "inputrec->awhParams->biasbUserData", biasIndex, bias1->bUserData != 0,
1132 bias2->bUserData != 0);
1133 cmp_double(fp, "inputrec->awhParams->biaserror_initial", biasIndex, bias1->errorInitial,
1134 bias2->errorInitial, ftol, abstol);
1135 cmp_int(fp, "inputrec->awhParams->biasShareGroup", biasIndex, bias1->shareGroup, bias2->shareGroup);
1137 for (int dim = 0; dim < std::min(bias1->ndim, bias2->ndim); dim++)
1139 cmp_awhDimParams(fp, &bias1->dimParams[dim], &bias2->dimParams[dim], dim, ftol, abstol);
1143 static void cmp_awhParams(FILE* fp, const gmx::AwhParams* awh1, const gmx::AwhParams* awh2, real ftol, real abstol)
1145 cmp_int(fp, "inputrec->awhParams->nbias", -1, awh1->numBias, awh2->numBias);
1146 cmp_int64(fp, "inputrec->awhParams->seed", awh1->seed, awh2->seed);
1147 cmp_int(fp, "inputrec->awhParams->nstout", -1, awh1->nstOut, awh2->nstOut);
1148 cmp_int(fp, "inputrec->awhParams->nstsample_coord", -1, awh1->nstSampleCoord, awh2->nstSampleCoord);
1149 cmp_int(fp, "inputrec->awhParams->nsamples_update_free_energy", -1,
1150 awh1->numSamplesUpdateFreeEnergy, awh2->numSamplesUpdateFreeEnergy);
1151 cmp_int(fp, "inputrec->awhParams->ePotential", -1, awh1->ePotential, awh2->ePotential);
1152 cmp_bool(fp, "inputrec->awhParams->shareBiasMultisim", -1, awh1->shareBiasMultisim,
1153 awh2->shareBiasMultisim);
1155 if (awh1->numBias == awh2->numBias)
1157 for (int bias = 0; bias < awh1->numBias; bias++)
1159 cmp_awhBiasParams(fp, &awh1->awhBiasParams[bias], &awh2->awhBiasParams[bias], bias, ftol, abstol);
1164 static void cmp_simtempvals(FILE* fp,
1165 const t_simtemp* simtemp1,
1166 const t_simtemp* simtemp2,
1167 int n_lambda,
1168 real ftol,
1169 real abstol)
1171 int i;
1172 cmp_int(fp, "inputrec->simtempvals->eSimTempScale", -1, simtemp1->eSimTempScale, simtemp2->eSimTempScale);
1173 cmp_real(fp, "inputrec->simtempvals->simtemp_high", -1, simtemp1->simtemp_high,
1174 simtemp2->simtemp_high, ftol, abstol);
1175 cmp_real(fp, "inputrec->simtempvals->simtemp_low", -1, simtemp1->simtemp_low,
1176 simtemp2->simtemp_low, ftol, abstol);
1177 for (i = 0; i < n_lambda; i++)
1179 cmp_real(fp, "inputrec->simtempvals->temperatures", -1, simtemp1->temperatures[i],
1180 simtemp2->temperatures[i], ftol, abstol);
1184 static void cmp_expandedvals(FILE* fp,
1185 const t_expanded* expand1,
1186 const t_expanded* expand2,
1187 int n_lambda,
1188 real ftol,
1189 real abstol)
1191 int i;
1193 cmp_bool(fp, "inputrec->fepvals->bInit_weights", -1, expand1->bInit_weights, expand2->bInit_weights);
1194 cmp_bool(fp, "inputrec->fepvals->bWLoneovert", -1, expand1->bWLoneovert, expand2->bWLoneovert);
1196 for (i = 0; i < n_lambda; i++)
1198 cmp_real(fp, "inputrec->expandedvals->init_lambda_weights", -1,
1199 expand1->init_lambda_weights[i], expand2->init_lambda_weights[i], ftol, abstol);
1202 cmp_int(fp, "inputrec->expandedvals->lambda-stats", -1, expand1->elamstats, expand2->elamstats);
1203 cmp_int(fp, "inputrec->expandedvals->lambda-mc-move", -1, expand1->elmcmove, expand2->elmcmove);
1204 cmp_int(fp, "inputrec->expandedvals->lmc-repeats", -1, expand1->lmc_repeats, expand2->lmc_repeats);
1205 cmp_int(fp, "inputrec->expandedvals->lmc-gibbsdelta", -1, expand1->gibbsdeltalam, expand2->gibbsdeltalam);
1206 cmp_int(fp, "inputrec->expandedvals->lmc-forced-nstart", -1, expand1->lmc_forced_nstart,
1207 expand2->lmc_forced_nstart);
1208 cmp_int(fp, "inputrec->expandedvals->lambda-weights-equil", -1, expand1->elmceq, expand2->elmceq);
1209 cmp_int(fp, "inputrec->expandedvals->,weight-equil-number-all-lambda", -1,
1210 expand1->equil_n_at_lam, expand2->equil_n_at_lam);
1211 cmp_int(fp, "inputrec->expandedvals->weight-equil-number-samples", -1, expand1->equil_samples,
1212 expand2->equil_samples);
1213 cmp_int(fp, "inputrec->expandedvals->weight-equil-number-steps", -1, expand1->equil_steps,
1214 expand2->equil_steps);
1215 cmp_real(fp, "inputrec->expandedvals->weight-equil-wl-delta", -1, expand1->equil_wl_delta,
1216 expand2->equil_wl_delta, ftol, abstol);
1217 cmp_real(fp, "inputrec->expandedvals->weight-equil-count-ratio", -1, expand1->equil_ratio,
1218 expand2->equil_ratio, ftol, abstol);
1219 cmp_bool(fp, "inputrec->expandedvals->symmetrized-transition-matrix", -1,
1220 expand1->bSymmetrizedTMatrix, expand2->bSymmetrizedTMatrix);
1221 cmp_int(fp, "inputrec->expandedvals->nstTij", -1, expand1->nstTij, expand2->nstTij);
1222 cmp_int(fp, "inputrec->expandedvals->mininum-var-min", -1, expand1->minvarmin,
1223 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,
1227 expand2->init_wl_delta, ftol, abstol);
1228 cmp_real(fp, "inputrec->expandedvals->wl-ratio", -1, expand1->wl_ratio, expand2->wl_ratio, ftol, abstol);
1229 cmp_int(fp, "inputrec->expandedvals->nstexpanded", -1, expand1->nstexpanded, expand2->nstexpanded);
1230 cmp_int(fp, "inputrec->expandedvals->lmc-seed", -1, expand1->lmc_seed, expand2->lmc_seed);
1231 cmp_real(fp, "inputrec->expandedvals->mc-temperature", -1, expand1->mc_temp, expand2->mc_temp,
1232 ftol, abstol);
1235 static void cmp_fepvals(FILE* fp, const t_lambda* fep1, const t_lambda* fep2, real ftol, real abstol)
1237 int i, j;
1238 cmp_int(fp, "inputrec->nstdhdl", -1, fep1->nstdhdl, fep2->nstdhdl);
1239 cmp_double(fp, "inputrec->fepvals->init_fep_state", -1, fep1->init_fep_state,
1240 fep2->init_fep_state, ftol, abstol);
1241 cmp_double(fp, "inputrec->fepvals->delta_lambda", -1, fep1->delta_lambda, fep2->delta_lambda,
1242 ftol, abstol);
1243 cmp_int(fp, "inputrec->fepvals->n_lambda", -1, fep1->n_lambda, fep2->n_lambda);
1244 for (i = 0; i < efptNR; i++)
1246 for (j = 0; j < std::min(fep1->n_lambda, fep2->n_lambda); j++)
1248 cmp_double(fp, "inputrec->fepvals->all_lambda", -1, fep1->all_lambda[i][j],
1249 fep2->all_lambda[i][j], ftol, abstol);
1252 cmp_int(fp, "inputrec->fepvals->lambda_neighbors", 1, fep1->lambda_neighbors, fep2->lambda_neighbors);
1253 cmp_real(fp, "inputrec->fepvals->sc_alpha", -1, fep1->sc_alpha, fep2->sc_alpha, ftol, abstol);
1254 cmp_int(fp, "inputrec->fepvals->sc_power", -1, fep1->sc_power, fep2->sc_power);
1255 cmp_real(fp, "inputrec->fepvals->sc_r_power", -1, fep1->sc_r_power, fep2->sc_r_power, ftol, abstol);
1256 cmp_real(fp, "inputrec->fepvals->sc_sigma", -1, fep1->sc_sigma, fep2->sc_sigma, ftol, abstol);
1257 cmp_int(fp, "inputrec->fepvals->edHdLPrintEnergy", -1, fep1->edHdLPrintEnergy, fep1->edHdLPrintEnergy);
1258 cmp_bool(fp, "inputrec->fepvals->bScCoul", -1, fep1->bScCoul, fep1->bScCoul);
1259 cmp_int(fp, "inputrec->separate_dhdl_file", -1, fep1->separate_dhdl_file, fep2->separate_dhdl_file);
1260 cmp_int(fp, "inputrec->dhdl_derivatives", -1, fep1->dhdl_derivatives, fep2->dhdl_derivatives);
1261 cmp_int(fp, "inputrec->dh_hist_size", -1, fep1->dh_hist_size, fep2->dh_hist_size);
1262 cmp_double(fp, "inputrec->dh_hist_spacing", -1, fep1->dh_hist_spacing, fep2->dh_hist_spacing,
1263 ftol, abstol);
1266 void cmp_inputrec(FILE* fp, const t_inputrec* ir1, const t_inputrec* ir2, real ftol, real abstol)
1268 fprintf(fp, "comparing inputrec\n");
1270 /* gcc 2.96 doesnt like these defines at all, but issues a huge list
1271 * of warnings. Maybe it will change in future versions, but for the
1272 * moment I've spelled them out instead. /EL 000820
1273 * #define CIB(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
1274 * #define CII(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
1275 * #define CIR(s) cmp_real(fp,"inputrec->"#s,0,ir1->##s,ir2->##s,ftol)
1277 cmp_int(fp, "inputrec->eI", -1, ir1->eI, ir2->eI);
1278 cmp_int64(fp, "inputrec->nsteps", ir1->nsteps, ir2->nsteps);
1279 cmp_int64(fp, "inputrec->init_step", ir1->init_step, ir2->init_step);
1280 cmp_int(fp, "inputrec->simulation_part", -1, ir1->simulation_part, ir2->simulation_part);
1281 cmp_int(fp, "inputrec->pbcType", -1, static_cast<int>(ir1->pbcType), static_cast<int>(ir2->pbcType));
1282 cmp_bool(fp, "inputrec->bPeriodicMols", -1, ir1->bPeriodicMols, ir2->bPeriodicMols);
1283 cmp_int(fp, "inputrec->cutoff_scheme", -1, ir1->cutoff_scheme, ir2->cutoff_scheme);
1284 cmp_int(fp, "inputrec->nstlist", -1, ir1->nstlist, ir2->nstlist);
1285 cmp_int(fp, "inputrec->nstcomm", -1, ir1->nstcomm, ir2->nstcomm);
1286 cmp_int(fp, "inputrec->comm_mode", -1, ir1->comm_mode, ir2->comm_mode);
1287 cmp_int(fp, "inputrec->nstlog", -1, ir1->nstlog, ir2->nstlog);
1288 cmp_int(fp, "inputrec->nstxout", -1, ir1->nstxout, ir2->nstxout);
1289 cmp_int(fp, "inputrec->nstvout", -1, ir1->nstvout, ir2->nstvout);
1290 cmp_int(fp, "inputrec->nstfout", -1, ir1->nstfout, ir2->nstfout);
1291 cmp_int(fp, "inputrec->nstcalcenergy", -1, ir1->nstcalcenergy, ir2->nstcalcenergy);
1292 cmp_int(fp, "inputrec->nstenergy", -1, ir1->nstenergy, ir2->nstenergy);
1293 cmp_int(fp, "inputrec->nstxout_compressed", -1, ir1->nstxout_compressed, ir2->nstxout_compressed);
1294 cmp_double(fp, "inputrec->init_t", -1, ir1->init_t, ir2->init_t, ftol, abstol);
1295 cmp_double(fp, "inputrec->delta_t", -1, ir1->delta_t, ir2->delta_t, ftol, abstol);
1296 cmp_real(fp, "inputrec->x_compression_precision", -1, ir1->x_compression_precision,
1297 ir2->x_compression_precision, ftol, abstol);
1298 cmp_real(fp, "inputrec->fourierspacing", -1, ir1->fourier_spacing, ir2->fourier_spacing, ftol, abstol);
1299 cmp_int(fp, "inputrec->nkx", -1, ir1->nkx, ir2->nkx);
1300 cmp_int(fp, "inputrec->nky", -1, ir1->nky, ir2->nky);
1301 cmp_int(fp, "inputrec->nkz", -1, ir1->nkz, ir2->nkz);
1302 cmp_int(fp, "inputrec->pme_order", -1, ir1->pme_order, ir2->pme_order);
1303 cmp_real(fp, "inputrec->ewald_rtol", -1, ir1->ewald_rtol, ir2->ewald_rtol, ftol, abstol);
1304 cmp_int(fp, "inputrec->ewald_geometry", -1, ir1->ewald_geometry, ir2->ewald_geometry);
1305 cmp_real(fp, "inputrec->epsilon_surface", -1, ir1->epsilon_surface, ir2->epsilon_surface, ftol, abstol);
1306 cmp_int(fp, "inputrec->bContinuation", -1, static_cast<int>(ir1->bContinuation),
1307 static_cast<int>(ir2->bContinuation));
1308 cmp_int(fp, "inputrec->bShakeSOR", -1, static_cast<int>(ir1->bShakeSOR),
1309 static_cast<int>(ir2->bShakeSOR));
1310 cmp_int(fp, "inputrec->etc", -1, ir1->etc, ir2->etc);
1311 cmp_int(fp, "inputrec->bPrintNHChains", -1, static_cast<int>(ir1->bPrintNHChains),
1312 static_cast<int>(ir2->bPrintNHChains));
1313 cmp_int(fp, "inputrec->epc", -1, ir1->epc, ir2->epc);
1314 cmp_int(fp, "inputrec->epct", -1, ir1->epct, ir2->epct);
1315 cmp_real(fp, "inputrec->tau_p", -1, ir1->tau_p, ir2->tau_p, ftol, abstol);
1316 cmp_rvec(fp, "inputrec->ref_p(x)", -1, ir1->ref_p[XX], ir2->ref_p[XX], ftol, abstol);
1317 cmp_rvec(fp, "inputrec->ref_p(y)", -1, ir1->ref_p[YY], ir2->ref_p[YY], ftol, abstol);
1318 cmp_rvec(fp, "inputrec->ref_p(z)", -1, ir1->ref_p[ZZ], ir2->ref_p[ZZ], ftol, abstol);
1319 cmp_rvec(fp, "inputrec->compress(x)", -1, ir1->compress[XX], ir2->compress[XX], ftol, abstol);
1320 cmp_rvec(fp, "inputrec->compress(y)", -1, ir1->compress[YY], ir2->compress[YY], ftol, abstol);
1321 cmp_rvec(fp, "inputrec->compress(z)", -1, ir1->compress[ZZ], ir2->compress[ZZ], ftol, abstol);
1322 cmp_int(fp, "refcoord_scaling", -1, ir1->refcoord_scaling, ir2->refcoord_scaling);
1323 cmp_rvec(fp, "inputrec->posres_com", -1, ir1->posres_com, ir2->posres_com, ftol, abstol);
1324 cmp_rvec(fp, "inputrec->posres_comB", -1, ir1->posres_comB, ir2->posres_comB, ftol, abstol);
1325 cmp_real(fp, "inputrec->verletbuf_tol", -1, ir1->verletbuf_tol, ir2->verletbuf_tol, ftol, abstol);
1326 cmp_real(fp, "inputrec->rlist", -1, ir1->rlist, ir2->rlist, ftol, abstol);
1327 cmp_real(fp, "inputrec->rtpi", -1, ir1->rtpi, ir2->rtpi, ftol, abstol);
1328 cmp_int(fp, "inputrec->coulombtype", -1, ir1->coulombtype, ir2->coulombtype);
1329 cmp_int(fp, "inputrec->coulomb_modifier", -1, ir1->coulomb_modifier, ir2->coulomb_modifier);
1330 cmp_real(fp, "inputrec->rcoulomb_switch", -1, ir1->rcoulomb_switch, ir2->rcoulomb_switch, ftol, abstol);
1331 cmp_real(fp, "inputrec->rcoulomb", -1, ir1->rcoulomb, ir2->rcoulomb, ftol, abstol);
1332 cmp_int(fp, "inputrec->vdwtype", -1, ir1->vdwtype, ir2->vdwtype);
1333 cmp_int(fp, "inputrec->vdw_modifier", -1, ir1->vdw_modifier, ir2->vdw_modifier);
1334 cmp_real(fp, "inputrec->rvdw_switch", -1, ir1->rvdw_switch, ir2->rvdw_switch, ftol, abstol);
1335 cmp_real(fp, "inputrec->rvdw", -1, ir1->rvdw, ir2->rvdw, ftol, abstol);
1336 cmp_real(fp, "inputrec->epsilon_r", -1, ir1->epsilon_r, ir2->epsilon_r, ftol, abstol);
1337 cmp_real(fp, "inputrec->epsilon_rf", -1, ir1->epsilon_rf, ir2->epsilon_rf, ftol, abstol);
1338 cmp_real(fp, "inputrec->tabext", -1, ir1->tabext, ir2->tabext, ftol, abstol);
1340 cmp_int(fp, "inputrec->eDispCorr", -1, ir1->eDispCorr, ir2->eDispCorr);
1341 cmp_real(fp, "inputrec->shake_tol", -1, ir1->shake_tol, ir2->shake_tol, ftol, abstol);
1342 cmp_int(fp, "inputrec->efep", -1, ir1->efep, ir2->efep);
1343 cmp_fepvals(fp, ir1->fepvals, ir2->fepvals, ftol, abstol);
1344 cmp_int(fp, "inputrec->bSimTemp", -1, static_cast<int>(ir1->bSimTemp), static_cast<int>(ir2->bSimTemp));
1345 if ((ir1->bSimTemp == ir2->bSimTemp) && (ir1->bSimTemp))
1347 cmp_simtempvals(fp, ir1->simtempvals, ir2->simtempvals,
1348 std::min(ir1->fepvals->n_lambda, ir2->fepvals->n_lambda), ftol, abstol);
1350 cmp_int(fp, "inputrec->bExpanded", -1, static_cast<int>(ir1->bExpanded),
1351 static_cast<int>(ir2->bExpanded));
1352 if ((ir1->bExpanded == ir2->bExpanded) && (ir1->bExpanded))
1354 cmp_expandedvals(fp, ir1->expandedvals, ir2->expandedvals,
1355 std::min(ir1->fepvals->n_lambda, ir2->fepvals->n_lambda), ftol, abstol);
1357 cmp_int(fp, "inputrec->nwall", -1, ir1->nwall, ir2->nwall);
1358 cmp_int(fp, "inputrec->wall_type", -1, ir1->wall_type, ir2->wall_type);
1359 cmp_int(fp, "inputrec->wall_atomtype[0]", -1, ir1->wall_atomtype[0], ir2->wall_atomtype[0]);
1360 cmp_int(fp, "inputrec->wall_atomtype[1]", -1, ir1->wall_atomtype[1], ir2->wall_atomtype[1]);
1361 cmp_real(fp, "inputrec->wall_density[0]", -1, ir1->wall_density[0], ir2->wall_density[0], ftol, abstol);
1362 cmp_real(fp, "inputrec->wall_density[1]", -1, ir1->wall_density[1], ir2->wall_density[1], ftol, abstol);
1363 cmp_real(fp, "inputrec->wall_ewald_zfac", -1, ir1->wall_ewald_zfac, ir2->wall_ewald_zfac, ftol, abstol);
1365 cmp_bool(fp, "inputrec->bPull", -1, ir1->bPull, ir2->bPull);
1366 if (ir1->bPull && ir2->bPull)
1368 cmp_pull(fp);
1371 cmp_bool(fp, "inputrec->bDoAwh", -1, ir1->bDoAwh, ir2->bDoAwh);
1372 if (ir1->bDoAwh && ir2->bDoAwh)
1374 cmp_awhParams(fp, ir1->awhParams, ir2->awhParams, ftol, abstol);
1377 cmp_int(fp, "inputrec->eDisre", -1, ir1->eDisre, ir2->eDisre);
1378 cmp_real(fp, "inputrec->dr_fc", -1, ir1->dr_fc, ir2->dr_fc, ftol, abstol);
1379 cmp_int(fp, "inputrec->eDisreWeighting", -1, ir1->eDisreWeighting, ir2->eDisreWeighting);
1380 cmp_int(fp, "inputrec->bDisreMixed", -1, static_cast<int>(ir1->bDisreMixed),
1381 static_cast<int>(ir2->bDisreMixed));
1382 cmp_int(fp, "inputrec->nstdisreout", -1, ir1->nstdisreout, ir2->nstdisreout);
1383 cmp_real(fp, "inputrec->dr_tau", -1, ir1->dr_tau, ir2->dr_tau, ftol, abstol);
1384 cmp_real(fp, "inputrec->orires_fc", -1, ir1->orires_fc, ir2->orires_fc, ftol, abstol);
1385 cmp_real(fp, "inputrec->orires_tau", -1, ir1->orires_tau, ir2->orires_tau, ftol, abstol);
1386 cmp_int(fp, "inputrec->nstorireout", -1, ir1->nstorireout, ir2->nstorireout);
1387 cmp_real(fp, "inputrec->em_stepsize", -1, ir1->em_stepsize, ir2->em_stepsize, ftol, abstol);
1388 cmp_real(fp, "inputrec->em_tol", -1, ir1->em_tol, ir2->em_tol, ftol, abstol);
1389 cmp_int(fp, "inputrec->niter", -1, ir1->niter, ir2->niter);
1390 cmp_real(fp, "inputrec->fc_stepsize", -1, ir1->fc_stepsize, ir2->fc_stepsize, ftol, abstol);
1391 cmp_int(fp, "inputrec->nstcgsteep", -1, ir1->nstcgsteep, ir2->nstcgsteep);
1392 cmp_int(fp, "inputrec->nbfgscorr", 0, ir1->nbfgscorr, ir2->nbfgscorr);
1393 cmp_int(fp, "inputrec->eConstrAlg", -1, ir1->eConstrAlg, ir2->eConstrAlg);
1394 cmp_int(fp, "inputrec->nProjOrder", -1, ir1->nProjOrder, ir2->nProjOrder);
1395 cmp_real(fp, "inputrec->LincsWarnAngle", -1, ir1->LincsWarnAngle, ir2->LincsWarnAngle, ftol, abstol);
1396 cmp_int(fp, "inputrec->nLincsIter", -1, ir1->nLincsIter, ir2->nLincsIter);
1397 cmp_real(fp, "inputrec->bd_fric", -1, ir1->bd_fric, ir2->bd_fric, ftol, abstol);
1398 cmp_int64(fp, "inputrec->ld_seed", ir1->ld_seed, ir2->ld_seed);
1399 cmp_real(fp, "inputrec->cos_accel", -1, ir1->cos_accel, ir2->cos_accel, ftol, abstol);
1400 cmp_rvec(fp, "inputrec->deform(a)", -1, ir1->deform[XX], ir2->deform[XX], ftol, abstol);
1401 cmp_rvec(fp, "inputrec->deform(b)", -1, ir1->deform[YY], ir2->deform[YY], ftol, abstol);
1402 cmp_rvec(fp, "inputrec->deform(c)", -1, ir1->deform[ZZ], ir2->deform[ZZ], ftol, abstol);
1405 cmp_int(fp, "inputrec->userint1", -1, ir1->userint1, ir2->userint1);
1406 cmp_int(fp, "inputrec->userint2", -1, ir1->userint2, ir2->userint2);
1407 cmp_int(fp, "inputrec->userint3", -1, ir1->userint3, ir2->userint3);
1408 cmp_int(fp, "inputrec->userint4", -1, ir1->userint4, ir2->userint4);
1409 cmp_real(fp, "inputrec->userreal1", -1, ir1->userreal1, ir2->userreal1, ftol, abstol);
1410 cmp_real(fp, "inputrec->userreal2", -1, ir1->userreal2, ir2->userreal2, ftol, abstol);
1411 cmp_real(fp, "inputrec->userreal3", -1, ir1->userreal3, ir2->userreal3, ftol, abstol);
1412 cmp_real(fp, "inputrec->userreal4", -1, ir1->userreal4, ir2->userreal4, ftol, abstol);
1413 cmp_grpopts(fp, &(ir1->opts), &(ir2->opts), ftol, abstol);
1414 gmx::TextWriter writer(fp);
1415 gmx::compareKeyValueTrees(&writer, *ir1->params, *ir2->params, ftol, abstol);
1418 void comp_pull_AB(FILE* fp, pull_params_t* pull, real ftol, real abstol)
1420 int i;
1422 for (i = 0; i < pull->ncoord; i++)
1424 fprintf(fp, "comparing pull coord %d\n", i);
1425 cmp_real(fp, "pull-coord->k", -1, pull->coord[i].k, pull->coord[i].kB, ftol, abstol);
1429 gmx_bool inputrecDeform(const t_inputrec* ir)
1431 return (ir->deform[XX][XX] != 0 || ir->deform[YY][YY] != 0 || ir->deform[ZZ][ZZ] != 0
1432 || ir->deform[YY][XX] != 0 || ir->deform[ZZ][XX] != 0 || ir->deform[ZZ][YY] != 0);
1435 gmx_bool inputrecDynamicBox(const t_inputrec* ir)
1437 return (ir->epc != epcNO || ir->eI == eiTPI || inputrecDeform(ir));
1440 gmx_bool inputrecPreserveShape(const t_inputrec* ir)
1442 return (ir->epc != epcNO && ir->deform[XX][XX] == 0
1443 && (ir->epct == epctISOTROPIC || ir->epct == epctSEMIISOTROPIC));
1446 gmx_bool inputrecNeedMutot(const t_inputrec* ir)
1448 return ((ir->coulombtype == eelEWALD || EEL_PME(ir->coulombtype))
1449 && (ir->ewald_geometry == eewg3DC || ir->epsilon_surface != 0));
1452 gmx_bool inputrecExclForces(const t_inputrec* ir)
1454 return (EEL_FULL(ir->coulombtype) || (EEL_RF(ir->coulombtype)));
1457 gmx_bool inputrecNptTrotter(const t_inputrec* ir)
1459 return (((ir->eI == eiVV) || (ir->eI == eiVVAK)) && (ir->epc == epcMTTK) && (ir->etc == etcNOSEHOOVER));
1462 gmx_bool inputrecNvtTrotter(const t_inputrec* ir)
1464 return (((ir->eI == eiVV) || (ir->eI == eiVVAK)) && (ir->epc != epcMTTK) && (ir->etc == etcNOSEHOOVER));
1467 gmx_bool inputrecNphTrotter(const t_inputrec* ir)
1469 return (((ir->eI == eiVV) || (ir->eI == eiVVAK)) && (ir->epc == epcMTTK) && (ir->etc != etcNOSEHOOVER));
1472 bool inputrecPbcXY2Walls(const t_inputrec* ir)
1474 return (ir->pbcType == PbcType::XY && ir->nwall == 2);
1477 bool integratorHasConservedEnergyQuantity(const t_inputrec* ir)
1479 if (!EI_MD(ir->eI))
1480 { // NOLINT bugprone-branch-clone
1481 // Energy minimization or stochastic integrator: no conservation
1482 return false;
1484 else if (ir->etc == etcNO && ir->epc == epcNO)
1486 // The total energy is conserved, no additional conserved quanitity
1487 return false;
1489 else
1491 // Shear stress with Parrinello-Rahman is not supported (tedious)
1492 bool shearWithPR =
1493 ((ir->epc == epcPARRINELLORAHMAN || ir->epc == epcMTTK)
1494 && (ir->ref_p[YY][XX] != 0 || ir->ref_p[ZZ][XX] != 0 || ir->ref_p[ZZ][YY] != 0));
1496 return !ETC_ANDERSEN(ir->etc) && !shearWithPR;
1500 bool integratorHasReferenceTemperature(const t_inputrec* ir)
1502 return ((ir->etc != etcNO) || EI_SD(ir->eI) || (ir->eI == eiBD) || EI_TPI(ir->eI));
1505 int inputrec2nboundeddim(const t_inputrec* ir)
1507 if (inputrecPbcXY2Walls(ir))
1509 return 3;
1511 else
1513 return numPbcDimensions(ir->pbcType);
1517 int ndof_com(const t_inputrec* ir)
1519 int n = 0;
1521 switch (ir->pbcType)
1523 case PbcType::Xyz:
1524 case PbcType::No: n = 3; break;
1525 case PbcType::XY: n = (ir->nwall == 0 ? 3 : 2); break;
1526 case PbcType::Screw: n = 1; break;
1527 default: gmx_incons("Unknown pbc in calc_nrdf");
1530 return n;
1533 real maxReferenceTemperature(const t_inputrec& ir)
1535 if (EI_ENERGY_MINIMIZATION(ir.eI) || ir.eI == eiNM)
1537 return 0;
1540 if (EI_MD(ir.eI) && ir.etc == etcNO)
1542 return -1;
1545 /* SD and BD also use ref_t and tau_t for setting the reference temperature.
1546 * TPI can be treated as MD, since it needs an ensemble temperature.
1549 real maxTemperature = 0;
1550 for (int i = 0; i < ir.opts.ngtc; i++)
1552 if (ir.opts.tau_t[i] >= 0)
1554 maxTemperature = std::max(maxTemperature, ir.opts.ref_t[i]);
1558 return maxTemperature;
1561 bool haveEwaldSurfaceContribution(const t_inputrec& ir)
1563 return EEL_PME_EWALD(ir.coulombtype) && (ir.ewald_geometry == eewg3DC || ir.epsilon_surface != 0);
1566 bool haveFreeEnergyType(const t_inputrec& ir, const int fepType)
1568 for (int i = 0; i < ir.fepvals->n_lambda; i++)
1570 if (ir.fepvals->all_lambda[fepType][i] > 0)
1572 return true;
1575 return false;