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.
49 #include "gromacs/math/veccompare.h"
50 #include "gromacs/math/vecdump.h"
51 #include "gromacs/mdtypes/awh_params.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/mdtypes/multipletimestepping.h"
54 #include "gromacs/mdtypes/pull_params.h"
55 #include "gromacs/pbcutil/pbc.h"
56 #include "gromacs/utility/compare.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/keyvaluetree.h"
60 #include "gromacs/utility/smalloc.h"
61 #include "gromacs/utility/snprintf.h"
62 #include "gromacs/utility/strconvert.h"
63 #include "gromacs/utility/stringutil.h"
64 #include "gromacs/utility/textwriter.h"
65 #include "gromacs/utility/txtdump.h"
67 //! Macro to select a bool name
68 #define EBOOL(e) gmx::boolToString(e)
70 /* Default values for nstcalcenergy, used when the are no other restrictions. */
71 constexpr int c_defaultNstCalcEnergy
= 10;
73 /* The minimum number of integration steps required for reasonably accurate
74 * integration of first and second order coupling algorithms.
76 const int nstmin_berendsen_tcouple
= 5;
77 const int nstmin_berendsen_pcouple
= 10;
78 const int nstmin_harmonic
= 20;
80 /* Default values for T- and P- coupling intervals, used when the are no other
83 constexpr int c_defaultNstTCouple
= 10;
84 constexpr int c_defaultNstPCouple
= 10;
86 t_inputrec::t_inputrec()
88 // TODO When this memset is removed, remove the suppression of
89 // gcc -Wno-class-memaccess in a CMakeLists.txt file.
90 std::memset(this, 0, sizeof(*this)); // NOLINT(bugprone-undefined-memory-manipulation)
92 snew(expandedvals
, 1);
96 t_inputrec::~t_inputrec()
101 int ir_optimal_nstcalcenergy(const t_inputrec
* ir
)
111 nst
= c_defaultNstCalcEnergy
;
116 nst
= std::lcm(nst
, ir
->mtsLevels
.back().stepFactor
);
122 int tcouple_min_integration_steps(int etc
)
128 case etcNO
: n
= 0; break;
130 case etcYES
: n
= nstmin_berendsen_tcouple
; break;
132 /* V-rescale supports instantaneous rescaling */
135 case etcNOSEHOOVER
: n
= nstmin_harmonic
; break;
137 case etcANDERSENMASSIVE
: n
= 1; break;
138 default: gmx_incons("Unknown etc value");
144 int ir_optimal_nsttcouple(const t_inputrec
* ir
)
146 int nmin
, nwanted
, n
;
150 nmin
= tcouple_min_integration_steps(ir
->etc
);
152 nwanted
= c_defaultNstTCouple
;
155 if (ir
->etc
!= etcNO
)
157 for (g
= 0; g
< ir
->opts
.ngtc
; g
++)
159 if (ir
->opts
.tau_t
[g
] > 0)
161 tau_min
= std::min(tau_min
, ir
->opts
.tau_t
[g
]);
166 if (nmin
== 0 || ir
->delta_t
* nwanted
<= tau_min
)
172 n
= static_cast<int>(tau_min
/ (ir
->delta_t
* nmin
) + 0.001);
177 while (nwanted
% n
!= 0)
186 int pcouple_min_integration_steps(int epc
)
192 case epcNO
: n
= 0; break;
195 case epcISOTROPIC
: n
= nstmin_berendsen_pcouple
; break;
196 case epcPARRINELLORAHMAN
:
197 case epcMTTK
: n
= nstmin_harmonic
; break;
198 default: gmx_incons("Unknown epc value");
204 int ir_optimal_nstpcouple(const t_inputrec
* ir
)
206 const int minIntegrationSteps
= pcouple_min_integration_steps(ir
->epc
);
208 const int nwanted
= c_defaultNstPCouple
;
210 // With multiple time-stepping we can only compute the pressure at slowest steps
211 const int minNstPCouple
= (ir
->useMts
? ir
->mtsLevels
.back().stepFactor
: 1);
214 if (minIntegrationSteps
== 0 || ir
->delta_t
* nwanted
<= ir
->tau_p
)
220 n
= static_cast<int>(ir
->tau_p
/ (ir
->delta_t
* minIntegrationSteps
) + 0.001);
221 if (n
< minNstPCouple
)
225 // Without MTS we try to make nstpcouple a "nice" number
228 while (nwanted
% n
!= 0)
235 // With MTS, nstpcouple should be a multiple of the slowest MTS interval
238 n
= n
- (n
% minNstPCouple
);
244 gmx_bool
ir_coulomb_switched(const t_inputrec
* ir
)
246 return (ir
->coulombtype
== eelSWITCH
|| ir
->coulombtype
== eelSHIFT
247 || ir
->coulombtype
== eelPMESWITCH
|| ir
->coulombtype
== eelPMEUSERSWITCH
248 || ir
->coulomb_modifier
== eintmodPOTSWITCH
|| ir
->coulomb_modifier
== eintmodFORCESWITCH
);
251 gmx_bool
ir_coulomb_is_zero_at_cutoff(const t_inputrec
* ir
)
253 return (ir
->cutoff_scheme
== ecutsVERLET
|| ir_coulomb_switched(ir
)
254 || ir
->coulomb_modifier
!= eintmodNONE
|| ir
->coulombtype
== eelRF_ZERO
);
257 gmx_bool
ir_coulomb_might_be_zero_at_cutoff(const t_inputrec
* ir
)
259 return (ir_coulomb_is_zero_at_cutoff(ir
) || ir
->coulombtype
== eelUSER
|| ir
->coulombtype
== eelPMEUSER
);
262 gmx_bool
ir_vdw_switched(const t_inputrec
* ir
)
264 return (ir
->vdwtype
== evdwSWITCH
|| ir
->vdwtype
== evdwSHIFT
265 || ir
->vdw_modifier
== eintmodPOTSWITCH
|| ir
->vdw_modifier
== eintmodFORCESWITCH
);
268 gmx_bool
ir_vdw_is_zero_at_cutoff(const t_inputrec
* ir
)
270 return (ir
->cutoff_scheme
== ecutsVERLET
|| ir_vdw_switched(ir
) || ir
->vdw_modifier
!= eintmodNONE
);
273 gmx_bool
ir_vdw_might_be_zero_at_cutoff(const t_inputrec
* ir
)
275 return (ir_vdw_is_zero_at_cutoff(ir
) || ir
->vdwtype
== evdwUSER
);
278 static void done_pull_group(t_pull_group
* pgrp
)
287 static void done_pull_params(pull_params_t
* pull
)
291 for (i
= 0; i
< pull
->ngroup
+ 1; i
++)
293 done_pull_group(pull
->group
);
300 static void done_lambdas(t_lambda
* fep
)
302 if (fep
->n_lambda
> 0)
304 for (int i
= 0; i
< efptNR
; i
++)
306 sfree(fep
->all_lambda
[i
]);
309 sfree(fep
->all_lambda
);
312 static void done_t_rot(t_rot
* rot
)
318 if (rot
->grp
!= nullptr)
320 for (int i
= 0; i
< rot
->ngrp
; i
++)
322 sfree(rot
->grp
[i
].ind
);
323 sfree(rot
->grp
[i
].x_ref
);
330 void done_inputrec(t_inputrec
* ir
)
332 sfree(ir
->opts
.nrdf
);
333 sfree(ir
->opts
.ref_t
);
334 for (int i
= 0; i
< ir
->opts
.ngtc
; i
++)
336 sfree(ir
->opts
.anneal_time
[i
]);
337 sfree(ir
->opts
.anneal_temp
[i
]);
339 sfree(ir
->opts
.annealing
);
340 sfree(ir
->opts
.anneal_npoints
);
341 sfree(ir
->opts
.anneal_time
);
342 sfree(ir
->opts
.anneal_temp
);
343 sfree(ir
->opts
.tau_t
);
345 sfree(ir
->opts
.nFreeze
);
346 sfree(ir
->opts
.egp_flags
);
347 done_lambdas(ir
->fepvals
);
349 sfree(ir
->expandedvals
);
350 sfree(ir
->simtempvals
);
354 done_pull_params(ir
->pull
);
361 static void pr_grp_opts(FILE* out
, int indent
, const char* title
, const t_grpopts
* opts
, gmx_bool bMDPformat
)
367 fprintf(out
, "%s:\n", title
);
370 pr_indent(out
, indent
);
371 fprintf(out
, "nrdf%s", bMDPformat
? " = " : ":");
372 for (i
= 0; (i
< opts
->ngtc
); i
++)
374 fprintf(out
, " %10g", opts
->nrdf
[i
]);
378 pr_indent(out
, indent
);
379 fprintf(out
, "ref-t%s", bMDPformat
? " = " : ":");
380 for (i
= 0; (i
< opts
->ngtc
); i
++)
382 fprintf(out
, " %10g", opts
->ref_t
[i
]);
386 pr_indent(out
, indent
);
387 fprintf(out
, "tau-t%s", bMDPformat
? " = " : ":");
388 for (i
= 0; (i
< opts
->ngtc
); i
++)
390 fprintf(out
, " %10g", opts
->tau_t
[i
]);
394 /* Pretty-print the simulated annealing info */
395 fprintf(out
, "annealing%s", bMDPformat
? " = " : ":");
396 for (i
= 0; (i
< opts
->ngtc
); i
++)
398 fprintf(out
, " %10s", EANNEAL(opts
->annealing
[i
]));
402 fprintf(out
, "annealing-npoints%s", bMDPformat
? " = " : ":");
403 for (i
= 0; (i
< opts
->ngtc
); i
++)
405 fprintf(out
, " %10d", opts
->anneal_npoints
[i
]);
409 for (i
= 0; (i
< opts
->ngtc
); i
++)
411 if (opts
->anneal_npoints
[i
] > 0)
413 fprintf(out
, "annealing-time [%d]:\t", i
);
414 for (j
= 0; (j
< opts
->anneal_npoints
[i
]); j
++)
416 fprintf(out
, " %10.1f", opts
->anneal_time
[i
][j
]);
419 fprintf(out
, "annealing-temp [%d]:\t", i
);
420 for (j
= 0; (j
< opts
->anneal_npoints
[i
]); j
++)
422 fprintf(out
, " %10.1f", opts
->anneal_temp
[i
][j
]);
428 pr_indent(out
, indent
);
429 fprintf(out
, "acc:\t");
430 for (i
= 0; (i
< opts
->ngacc
); i
++)
432 for (m
= 0; (m
< DIM
); m
++)
434 fprintf(out
, " %10g", opts
->acc
[i
][m
]);
439 pr_indent(out
, indent
);
440 fprintf(out
, "nfreeze:");
441 for (i
= 0; (i
< opts
->ngfrz
); i
++)
443 for (m
= 0; (m
< DIM
); m
++)
445 fprintf(out
, " %10s", opts
->nFreeze
[i
][m
] ? "Y" : "N");
451 for (i
= 0; (i
< opts
->ngener
); i
++)
453 pr_indent(out
, indent
);
454 fprintf(out
, "energygrp-flags[%3d]:", i
);
455 for (m
= 0; (m
< opts
->ngener
); m
++)
457 fprintf(out
, " %d", opts
->egp_flags
[opts
->ngener
* i
+ m
]);
465 static void pr_matrix(FILE* fp
, int indent
, const char* title
, const rvec
* m
, gmx_bool bMDPformat
)
469 fprintf(fp
, "%-10s = %g %g %g %g %g %g\n", title
, m
[XX
][XX
], m
[YY
][YY
], m
[ZZ
][ZZ
],
470 m
[XX
][YY
], m
[XX
][ZZ
], m
[YY
][ZZ
]);
474 pr_rvecs(fp
, indent
, title
, m
, DIM
);
478 #define PS(t, s) pr_str(fp, indent, t, s)
479 #define PI(t, s) pr_int(fp, indent, t, s)
480 #define PSTEP(t, s) pr_int64(fp, indent, t, s)
481 #define PR(t, s) pr_real(fp, indent, t, s)
482 #define PD(t, s) pr_double(fp, indent, t, s)
484 static void pr_pull_group(FILE* fp
, int indent
, int g
, const t_pull_group
* pgrp
)
486 pr_indent(fp
, indent
);
487 fprintf(fp
, "pull-group %d:\n", g
);
489 pr_ivec_block(fp
, indent
, "atom", pgrp
->ind
, pgrp
->nat
, TRUE
);
490 pr_rvec(fp
, indent
, "weight", pgrp
->weight
, pgrp
->nweight
, TRUE
);
491 PI("pbcatom", pgrp
->pbcatom
);
494 static void pr_pull_coord(FILE* fp
, int indent
, int c
, const t_pull_coord
* pcrd
)
498 pr_indent(fp
, indent
);
499 fprintf(fp
, "pull-coord %d:\n", c
);
500 PS("type", EPULLTYPE(pcrd
->eType
));
501 if (pcrd
->eType
== epullEXTERNAL
)
503 PS("potential-provider", pcrd
->externalPotentialProvider
);
505 PS("geometry", EPULLGEOM(pcrd
->eGeom
));
506 for (g
= 0; g
< pcrd
->ngroup
; g
++)
510 sprintf(buf
, "group[%d]", g
);
511 PI(buf
, pcrd
->group
[g
]);
513 pr_ivec(fp
, indent
, "dim", pcrd
->dim
, DIM
, TRUE
);
514 pr_rvec(fp
, indent
, "origin", pcrd
->origin
, DIM
, TRUE
);
515 pr_rvec(fp
, indent
, "vec", pcrd
->vec
, DIM
, TRUE
);
516 PS("start", EBOOL(pcrd
->bStart
));
517 PR("init", pcrd
->init
);
518 PR("rate", pcrd
->rate
);
523 static void pr_simtempvals(FILE* fp
, int indent
, const t_simtemp
* simtemp
, int n_lambda
)
525 PS("simulated-tempering-scaling", ESIMTEMP(simtemp
->eSimTempScale
));
526 PR("sim-temp-low", simtemp
->simtemp_low
);
527 PR("sim-temp-high", simtemp
->simtemp_high
);
528 pr_rvec(fp
, indent
, "simulated tempering temperatures", simtemp
->temperatures
, n_lambda
, TRUE
);
531 static void pr_expandedvals(FILE* fp
, int indent
, const t_expanded
* expand
, int n_lambda
)
534 PI("nstexpanded", expand
->nstexpanded
);
535 PS("lmc-stats", elamstats_names
[expand
->elamstats
]);
536 PS("lmc-move", elmcmove_names
[expand
->elmcmove
]);
537 PS("lmc-weights-equil", elmceq_names
[expand
->elmceq
]);
538 if (expand
->elmceq
== elmceqNUMATLAM
)
540 PI("weight-equil-number-all-lambda", expand
->equil_n_at_lam
);
542 if (expand
->elmceq
== elmceqSAMPLES
)
544 PI("weight-equil-number-samples", expand
->equil_samples
);
546 if (expand
->elmceq
== elmceqSTEPS
)
548 PI("weight-equil-number-steps", expand
->equil_steps
);
550 if (expand
->elmceq
== elmceqWLDELTA
)
552 PR("weight-equil-wl-delta", expand
->equil_wl_delta
);
554 if (expand
->elmceq
== elmceqRATIO
)
556 PR("weight-equil-count-ratio", expand
->equil_ratio
);
558 PI("lmc-seed", expand
->lmc_seed
);
559 PR("mc-temperature", expand
->mc_temp
);
560 PI("lmc-repeats", expand
->lmc_repeats
);
561 PI("lmc-gibbsdelta", expand
->gibbsdeltalam
);
562 PI("lmc-forced-nstart", expand
->lmc_forced_nstart
);
563 PS("symmetrized-transition-matrix", EBOOL(expand
->bSymmetrizedTMatrix
));
564 PI("nst-transition-matrix", expand
->nstTij
);
565 PI("mininum-var-min", expand
->minvarmin
); /*default is reasonable */
566 PI("weight-c-range", expand
->c_range
); /* default is just C=0 */
567 PR("wl-scale", expand
->wl_scale
);
568 PR("wl-ratio", expand
->wl_ratio
);
569 PR("init-wl-delta", expand
->init_wl_delta
);
570 PS("wl-oneovert", EBOOL(expand
->bWLoneovert
));
572 pr_indent(fp
, indent
);
573 pr_rvec(fp
, indent
, "init-lambda-weights", expand
->init_lambda_weights
, n_lambda
, TRUE
);
574 PS("init-weights", EBOOL(expand
->bInit_weights
));
577 static void pr_fepvals(FILE* fp
, int indent
, const t_lambda
* fep
, gmx_bool bMDPformat
)
581 PR("init-lambda", fep
->init_lambda
);
582 PI("init-lambda-state", fep
->init_fep_state
);
583 PR("delta-lambda", fep
->delta_lambda
);
584 PI("nstdhdl", fep
->nstdhdl
);
588 PI("n-lambdas", fep
->n_lambda
);
590 if (fep
->n_lambda
> 0)
592 pr_indent(fp
, indent
);
593 fprintf(fp
, "separate-dvdl%s\n", bMDPformat
? " = " : ":");
594 for (i
= 0; i
< efptNR
; i
++)
596 fprintf(fp
, "%18s = ", efpt_names
[i
]);
597 if (fep
->separate_dvdl
[i
])
599 fprintf(fp
, " TRUE");
603 fprintf(fp
, " FALSE");
607 fprintf(fp
, "all-lambdas%s\n", bMDPformat
? " = " : ":");
608 for (i
= 0; i
< efptNR
; i
++)
610 fprintf(fp
, "%18s = ", efpt_names
[i
]);
611 for (j
= 0; j
< fep
->n_lambda
; j
++)
613 fprintf(fp
, " %10g", fep
->all_lambda
[i
][j
]);
618 PI("calc-lambda-neighbors", fep
->lambda_neighbors
);
619 PS("dhdl-print-energy", edHdLPrintEnergy_names
[fep
->edHdLPrintEnergy
]);
620 PR("sc-alpha", fep
->sc_alpha
);
621 PI("sc-power", fep
->sc_power
);
622 PR("sc-r-power", fep
->sc_r_power
);
623 PR("sc-sigma", fep
->sc_sigma
);
624 PR("sc-sigma-min", fep
->sc_sigma_min
);
625 PS("sc-coul", EBOOL(fep
->bScCoul
));
626 PI("dh-hist-size", fep
->dh_hist_size
);
627 PD("dh-hist-spacing", fep
->dh_hist_spacing
);
628 PS("separate-dhdl-file", SEPDHDLFILETYPE(fep
->separate_dhdl_file
));
629 PS("dhdl-derivatives", DHDLDERIVATIVESTYPE(fep
->dhdl_derivatives
));
632 static void pr_pull(FILE* fp
, int indent
, const pull_params_t
* pull
)
636 PR("pull-cylinder-r", pull
->cylinder_r
);
637 PR("pull-constr-tol", pull
->constr_tol
);
638 PS("pull-print-COM", EBOOL(pull
->bPrintCOM
));
639 PS("pull-print-ref-value", EBOOL(pull
->bPrintRefValue
));
640 PS("pull-print-components", EBOOL(pull
->bPrintComp
));
641 PI("pull-nstxout", pull
->nstxout
);
642 PI("pull-nstfout", pull
->nstfout
);
643 PS("pull-pbc-ref-prev-step-com", EBOOL(pull
->bSetPbcRefToPrevStepCOM
));
644 PS("pull-xout-average", EBOOL(pull
->bXOutAverage
));
645 PS("pull-fout-average", EBOOL(pull
->bFOutAverage
));
646 PI("pull-ngroups", pull
->ngroup
);
647 for (g
= 0; g
< pull
->ngroup
; g
++)
649 pr_pull_group(fp
, indent
, g
, &pull
->group
[g
]);
651 PI("pull-ncoords", pull
->ncoord
);
652 for (g
= 0; g
< pull
->ncoord
; g
++)
654 pr_pull_coord(fp
, indent
, g
, &pull
->coord
[g
]);
658 static void pr_awh_bias_dim(FILE* fp
, int indent
, gmx::AwhDimParams
* awhDimParams
, const char* prefix
)
660 pr_indent(fp
, 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("cover-diameter", awhDimParams
->coverDiameter
);
673 static void pr_awh_bias(FILE* fp
, int indent
, gmx::AwhBiasParams
* awhBiasParams
, const char* prefix
)
677 sprintf(opt
, "%s-error-init", prefix
);
678 PR(opt
, awhBiasParams
->errorInitial
);
679 sprintf(opt
, "%s-growth", prefix
);
680 PS(opt
, EAWHGROWTH(awhBiasParams
->eGrowth
));
681 sprintf(opt
, "%s-target", prefix
);
682 PS(opt
, EAWHTARGET(awhBiasParams
->eTarget
));
683 sprintf(opt
, "%s-target-beta-scalng", prefix
);
684 PR(opt
, awhBiasParams
->targetBetaScaling
);
685 sprintf(opt
, "%s-target-cutoff", prefix
);
686 PR(opt
, awhBiasParams
->targetCutoff
);
687 sprintf(opt
, "%s-user-data", prefix
);
688 PS(opt
, EBOOL(awhBiasParams
->bUserData
));
689 sprintf(opt
, "%s-share-group", prefix
);
690 PI(opt
, awhBiasParams
->shareGroup
);
691 sprintf(opt
, "%s-equilibrate-histogram", prefix
);
692 PS(opt
, EBOOL(awhBiasParams
->equilibrateHistogram
));
693 sprintf(opt
, "%s-ndim", prefix
);
694 PI(opt
, awhBiasParams
->ndim
);
696 for (int d
= 0; d
< awhBiasParams
->ndim
; d
++)
698 char prefixdim
[STRLEN
];
699 sprintf(prefixdim
, "%s-dim%d", prefix
, d
+ 1);
700 pr_awh_bias_dim(fp
, indent
, &awhBiasParams
->dimParams
[d
], prefixdim
);
704 static void pr_awh(FILE* fp
, int indent
, gmx::AwhParams
* awhParams
)
706 PS("awh-potential", EAWHPOTENTIAL(awhParams
->ePotential
));
707 PI("awh-seed", awhParams
->seed
);
708 PI("awh-nstout", awhParams
->nstOut
);
709 PI("awh-nstsample", awhParams
->nstSampleCoord
);
710 PI("awh-nsamples-update", awhParams
->numSamplesUpdateFreeEnergy
);
711 PS("awh-share-bias-multisim", EBOOL(awhParams
->shareBiasMultisim
));
712 PI("awh-nbias", awhParams
->numBias
);
714 for (int k
= 0; k
< awhParams
->numBias
; k
++)
716 auto prefix
= gmx::formatString("awh%d", k
+ 1);
717 pr_awh_bias(fp
, indent
, &awhParams
->awhBiasParams
[k
], prefix
.c_str());
721 static void pr_rotgrp(FILE* fp
, int indent
, int g
, const t_rotgrp
* rotg
)
723 pr_indent(fp
, indent
);
724 fprintf(fp
, "rot-group %d:\n", g
);
726 PS("rot-type", EROTGEOM(rotg
->eType
));
727 PS("rot-massw", EBOOL(rotg
->bMassW
));
728 pr_ivec_block(fp
, indent
, "atom", rotg
->ind
, rotg
->nat
, TRUE
);
729 pr_rvecs(fp
, indent
, "x-ref", rotg
->x_ref
, rotg
->nat
);
730 pr_rvec(fp
, indent
, "rot-vec", rotg
->inputVec
, DIM
, TRUE
);
731 pr_rvec(fp
, indent
, "rot-pivot", rotg
->pivot
, DIM
, TRUE
);
732 PR("rot-rate", rotg
->rate
);
733 PR("rot-k", rotg
->k
);
734 PR("rot-slab-dist", rotg
->slab_dist
);
735 PR("rot-min-gauss", rotg
->min_gaussian
);
736 PR("rot-eps", rotg
->eps
);
737 PS("rot-fit-method", EROTFIT(rotg
->eFittype
));
738 PI("rot-potfit-nstep", rotg
->PotAngle_nstep
);
739 PR("rot-potfit-step", rotg
->PotAngle_step
);
742 static void pr_rot(FILE* fp
, int indent
, const t_rot
* rot
)
746 PI("rot-nstrout", rot
->nstrout
);
747 PI("rot-nstsout", rot
->nstsout
);
748 PI("rot-ngroups", rot
->ngrp
);
749 for (g
= 0; g
< rot
->ngrp
; g
++)
751 pr_rotgrp(fp
, indent
, g
, &rot
->grp
[g
]);
756 static void pr_swap(FILE* fp
, int indent
, const t_swapcoords
* swap
)
760 /* Enums for better readability of the code */
768 PI("swap-frequency", swap
->nstswap
);
770 /* The split groups that define the compartments */
771 for (int j
= 0; j
< 2; j
++)
773 snprintf(str
, STRLEN
, "massw_split%d", j
);
774 PS(str
, EBOOL(swap
->massw_split
[j
]));
775 snprintf(str
, STRLEN
, "split atoms group %d", j
);
776 pr_ivec_block(fp
, indent
, str
, swap
->grp
[j
].ind
, swap
->grp
[j
].nat
, TRUE
);
779 /* The solvent group */
780 snprintf(str
, STRLEN
, "solvent group %s", swap
->grp
[eGrpSolvent
].molname
);
781 pr_ivec_block(fp
, indent
, str
, swap
->grp
[eGrpSolvent
].ind
, swap
->grp
[eGrpSolvent
].nat
, TRUE
);
783 /* Now print the indices for all the ion groups: */
784 for (int ig
= eSwapFixedGrpNR
; ig
< swap
->ngrp
; ig
++)
786 snprintf(str
, STRLEN
, "ion group %s", swap
->grp
[ig
].molname
);
787 pr_ivec_block(fp
, indent
, str
, swap
->grp
[ig
].ind
, swap
->grp
[ig
].nat
, TRUE
);
790 PR("cyl0-r", swap
->cyl0r
);
791 PR("cyl0-up", swap
->cyl0u
);
792 PR("cyl0-down", swap
->cyl0l
);
793 PR("cyl1-r", swap
->cyl1r
);
794 PR("cyl1-up", swap
->cyl1u
);
795 PR("cyl1-down", swap
->cyl1l
);
796 PI("coupl-steps", swap
->nAverage
);
798 /* Print the requested ion counts for both compartments */
799 for (int ic
= eCompA
; ic
<= eCompB
; ic
++)
801 for (int ig
= eSwapFixedGrpNR
; ig
< swap
->ngrp
; ig
++)
803 snprintf(str
, STRLEN
, "%s-in-%c", swap
->grp
[ig
].molname
, 'A' + ic
);
804 PI(str
, swap
->grp
[ig
].nmolReq
[ic
]);
808 PR("threshold", swap
->threshold
);
809 PR("bulk-offsetA", swap
->bulkOffset
[eCompA
]);
810 PR("bulk-offsetB", swap
->bulkOffset
[eCompB
]);
814 static void pr_imd(FILE* fp
, int indent
, const t_IMD
* imd
)
816 PI("IMD-atoms", imd
->nat
);
817 pr_ivec_block(fp
, indent
, "atom", imd
->ind
, imd
->nat
, TRUE
);
821 void pr_inputrec(FILE* fp
, int indent
, const char* title
, const t_inputrec
* ir
, gmx_bool bMDPformat
)
823 const char* infbuf
= "inf";
825 if (available(fp
, ir
, indent
, title
))
829 indent
= pr_title(fp
, indent
, title
);
831 /* Try to make this list appear in the same order as the
832 * options are written in the default mdout.mdp, and with
833 * the same user-exposed names to facilitate debugging.
835 PS("integrator", EI(ir
->eI
));
836 PR("tinit", ir
->init_t
);
837 PR("dt", ir
->delta_t
);
838 PSTEP("nsteps", ir
->nsteps
);
839 PSTEP("init-step", ir
->init_step
);
840 PI("simulation-part", ir
->simulation_part
);
841 PS("mts", EBOOL(ir
->useMts
));
844 for (int mtsIndex
= 1; mtsIndex
< static_cast<int>(ir
->mtsLevels
.size()); mtsIndex
++)
846 const auto& mtsLevel
= ir
->mtsLevels
[mtsIndex
];
847 const std::string forceKey
= gmx::formatString("mts-level%d-forces", mtsIndex
+ 1);
848 std::string forceGroups
;
849 for (int i
= 0; i
< static_cast<int>(gmx::MtsForceGroups::Count
); i
++)
851 if (mtsLevel
.forceGroups
[i
])
853 if (!forceGroups
.empty())
857 forceGroups
+= gmx::mtsForceGroupNames
[i
];
860 PS(forceKey
.c_str(), forceGroups
.c_str());
861 const std::string factorKey
= gmx::formatString("mts-level%d-factor", mtsIndex
+ 1);
862 PI(factorKey
.c_str(), mtsLevel
.stepFactor
);
865 PS("comm-mode", ECOM(ir
->comm_mode
));
866 PI("nstcomm", ir
->nstcomm
);
868 /* Langevin dynamics */
869 PR("bd-fric", ir
->bd_fric
);
870 PSTEP("ld-seed", ir
->ld_seed
);
872 /* Energy minimization */
873 PR("emtol", ir
->em_tol
);
874 PR("emstep", ir
->em_stepsize
);
875 PI("niter", ir
->niter
);
876 PR("fcstep", ir
->fc_stepsize
);
877 PI("nstcgsteep", ir
->nstcgsteep
);
878 PI("nbfgscorr", ir
->nbfgscorr
);
880 /* Test particle insertion */
881 PR("rtpi", ir
->rtpi
);
884 PI("nstxout", ir
->nstxout
);
885 PI("nstvout", ir
->nstvout
);
886 PI("nstfout", ir
->nstfout
);
887 PI("nstlog", ir
->nstlog
);
888 PI("nstcalcenergy", ir
->nstcalcenergy
);
889 PI("nstenergy", ir
->nstenergy
);
890 PI("nstxout-compressed", ir
->nstxout_compressed
);
891 PR("compressed-x-precision", ir
->x_compression_precision
);
893 /* Neighborsearching parameters */
894 PS("cutoff-scheme", ECUTSCHEME(ir
->cutoff_scheme
));
895 PI("nstlist", ir
->nstlist
);
896 PS("pbc", c_pbcTypeNames
[ir
->pbcType
].c_str());
897 PS("periodic-molecules", EBOOL(ir
->bPeriodicMols
));
898 PR("verlet-buffer-tolerance", ir
->verletbuf_tol
);
899 PR("rlist", ir
->rlist
);
901 /* Options for electrostatics and VdW */
902 PS("coulombtype", EELTYPE(ir
->coulombtype
));
903 PS("coulomb-modifier", INTMODIFIER(ir
->coulomb_modifier
));
904 PR("rcoulomb-switch", ir
->rcoulomb_switch
);
905 PR("rcoulomb", ir
->rcoulomb
);
906 if (ir
->epsilon_r
!= 0)
908 PR("epsilon-r", ir
->epsilon_r
);
912 PS("epsilon-r", infbuf
);
914 if (ir
->epsilon_rf
!= 0)
916 PR("epsilon-rf", ir
->epsilon_rf
);
920 PS("epsilon-rf", infbuf
);
922 PS("vdw-type", EVDWTYPE(ir
->vdwtype
));
923 PS("vdw-modifier", INTMODIFIER(ir
->vdw_modifier
));
924 PR("rvdw-switch", ir
->rvdw_switch
);
925 PR("rvdw", ir
->rvdw
);
926 PS("DispCorr", EDISPCORR(ir
->eDispCorr
));
927 PR("table-extension", ir
->tabext
);
929 PR("fourierspacing", ir
->fourier_spacing
);
930 PI("fourier-nx", ir
->nkx
);
931 PI("fourier-ny", ir
->nky
);
932 PI("fourier-nz", ir
->nkz
);
933 PI("pme-order", ir
->pme_order
);
934 PR("ewald-rtol", ir
->ewald_rtol
);
935 PR("ewald-rtol-lj", ir
->ewald_rtol_lj
);
936 PS("lj-pme-comb-rule", ELJPMECOMBNAMES(ir
->ljpme_combination_rule
));
937 PR("ewald-geometry", ir
->ewald_geometry
);
938 PR("epsilon-surface", ir
->epsilon_surface
);
940 /* Options for weak coupling algorithms */
941 PS("tcoupl", ETCOUPLTYPE(ir
->etc
));
942 PI("nsttcouple", ir
->nsttcouple
);
943 PI("nh-chain-length", ir
->opts
.nhchainlength
);
944 PS("print-nose-hoover-chain-variables", EBOOL(ir
->bPrintNHChains
));
946 PS("pcoupl", EPCOUPLTYPE(ir
->epc
));
947 PS("pcoupltype", EPCOUPLTYPETYPE(ir
->epct
));
948 PI("nstpcouple", ir
->nstpcouple
);
949 PR("tau-p", ir
->tau_p
);
950 pr_matrix(fp
, indent
, "compressibility", ir
->compress
, bMDPformat
);
951 pr_matrix(fp
, indent
, "ref-p", ir
->ref_p
, bMDPformat
);
952 PS("refcoord-scaling", EREFSCALINGTYPE(ir
->refcoord_scaling
));
956 fprintf(fp
, "posres-com = %g %g %g\n", ir
->posres_com
[XX
], ir
->posres_com
[YY
],
958 fprintf(fp
, "posres-comB = %g %g %g\n", ir
->posres_comB
[XX
], ir
->posres_comB
[YY
],
959 ir
->posres_comB
[ZZ
]);
963 pr_rvec(fp
, indent
, "posres-com", ir
->posres_com
, DIM
, TRUE
);
964 pr_rvec(fp
, indent
, "posres-comB", ir
->posres_comB
, DIM
, TRUE
);
968 PS("QMMM", EBOOL(ir
->bQMMM
));
969 fprintf(fp
, "%s:\n", "qm-opts");
970 pr_int(fp
, indent
, "ngQM", ir
->opts
.ngQM
);
972 /* CONSTRAINT OPTIONS */
973 PS("constraint-algorithm", ECONSTRTYPE(ir
->eConstrAlg
));
974 PS("continuation", EBOOL(ir
->bContinuation
));
976 PS("Shake-SOR", EBOOL(ir
->bShakeSOR
));
977 PR("shake-tol", ir
->shake_tol
);
978 PI("lincs-order", ir
->nProjOrder
);
979 PI("lincs-iter", ir
->nLincsIter
);
980 PR("lincs-warnangle", ir
->LincsWarnAngle
);
983 PI("nwall", ir
->nwall
);
984 PS("wall-type", EWALLTYPE(ir
->wall_type
));
985 PR("wall-r-linpot", ir
->wall_r_linpot
);
987 PI("wall-atomtype[0]", ir
->wall_atomtype
[0]);
988 PI("wall-atomtype[1]", ir
->wall_atomtype
[1]);
990 PR("wall-density[0]", ir
->wall_density
[0]);
991 PR("wall-density[1]", ir
->wall_density
[1]);
992 PR("wall-ewald-zfac", ir
->wall_ewald_zfac
);
995 PS("pull", EBOOL(ir
->bPull
));
998 pr_pull(fp
, indent
, ir
->pull
);
1002 PS("awh", EBOOL(ir
->bDoAwh
));
1005 pr_awh(fp
, indent
, ir
->awhParams
);
1008 /* ENFORCED ROTATION */
1009 PS("rotation", EBOOL(ir
->bRot
));
1012 pr_rot(fp
, indent
, ir
->rot
);
1015 /* INTERACTIVE MD */
1016 PS("interactiveMD", EBOOL(ir
->bIMD
));
1019 pr_imd(fp
, indent
, ir
->imd
);
1022 /* NMR refinement stuff */
1023 PS("disre", EDISRETYPE(ir
->eDisre
));
1024 PS("disre-weighting", EDISREWEIGHTING(ir
->eDisreWeighting
));
1025 PS("disre-mixed", EBOOL(ir
->bDisreMixed
));
1026 PR("dr-fc", ir
->dr_fc
);
1027 PR("dr-tau", ir
->dr_tau
);
1028 PR("nstdisreout", ir
->nstdisreout
);
1030 PR("orire-fc", ir
->orires_fc
);
1031 PR("orire-tau", ir
->orires_tau
);
1032 PR("nstorireout", ir
->nstorireout
);
1034 /* FREE ENERGY VARIABLES */
1035 PS("free-energy", EFEPTYPE(ir
->efep
));
1036 if (ir
->efep
!= efepNO
|| ir
->bSimTemp
)
1038 pr_fepvals(fp
, indent
, ir
->fepvals
, bMDPformat
);
1042 pr_expandedvals(fp
, indent
, ir
->expandedvals
, ir
->fepvals
->n_lambda
);
1045 /* NON-equilibrium MD stuff */
1046 PR("cos-acceleration", ir
->cos_accel
);
1047 pr_matrix(fp
, indent
, "deform", ir
->deform
, bMDPformat
);
1049 /* SIMULATED TEMPERING */
1050 PS("simulated-tempering", EBOOL(ir
->bSimTemp
));
1053 pr_simtempvals(fp
, indent
, ir
->simtempvals
, ir
->fepvals
->n_lambda
);
1056 /* ION/WATER SWAPPING FOR COMPUTATIONAL ELECTROPHYSIOLOGY */
1057 PS("swapcoords", ESWAPTYPE(ir
->eSwapCoords
));
1058 if (ir
->eSwapCoords
!= eswapNO
)
1060 pr_swap(fp
, indent
, ir
->swap
);
1063 /* USER-DEFINED THINGIES */
1064 PI("userint1", ir
->userint1
);
1065 PI("userint2", ir
->userint2
);
1066 PI("userint3", ir
->userint3
);
1067 PI("userint4", ir
->userint4
);
1068 PR("userreal1", ir
->userreal1
);
1069 PR("userreal2", ir
->userreal2
);
1070 PR("userreal3", ir
->userreal3
);
1071 PR("userreal4", ir
->userreal4
);
1075 gmx::TextWriter
writer(fp
);
1076 writer
.wrapperSettings().setIndent(indent
);
1077 gmx::dumpKeyValueTree(&writer
, *ir
->params
);
1080 pr_grp_opts(fp
, indent
, "grpopts", &(ir
->opts
), bMDPformat
);
1087 static void cmp_grpopts(FILE* fp
, const t_grpopts
* opt1
, const t_grpopts
* opt2
, real ftol
, real abstol
)
1090 char buf1
[256], buf2
[256];
1092 cmp_int(fp
, "inputrec->grpopts.ngtc", -1, opt1
->ngtc
, opt2
->ngtc
);
1093 cmp_int(fp
, "inputrec->grpopts.ngacc", -1, opt1
->ngacc
, opt2
->ngacc
);
1094 cmp_int(fp
, "inputrec->grpopts.ngfrz", -1, opt1
->ngfrz
, opt2
->ngfrz
);
1095 cmp_int(fp
, "inputrec->grpopts.ngener", -1, opt1
->ngener
, opt2
->ngener
);
1096 for (i
= 0; (i
< std::min(opt1
->ngtc
, opt2
->ngtc
)); i
++)
1098 cmp_real(fp
, "inputrec->grpopts.nrdf", i
, opt1
->nrdf
[i
], opt2
->nrdf
[i
], ftol
, abstol
);
1099 cmp_real(fp
, "inputrec->grpopts.ref_t", i
, opt1
->ref_t
[i
], opt2
->ref_t
[i
], ftol
, abstol
);
1100 cmp_real(fp
, "inputrec->grpopts.tau_t", i
, opt1
->tau_t
[i
], opt2
->tau_t
[i
], ftol
, abstol
);
1101 cmp_int(fp
, "inputrec->grpopts.annealing", i
, opt1
->annealing
[i
], opt2
->annealing
[i
]);
1102 cmp_int(fp
, "inputrec->grpopts.anneal_npoints", i
, opt1
->anneal_npoints
[i
],
1103 opt2
->anneal_npoints
[i
]);
1104 if (opt1
->anneal_npoints
[i
] == opt2
->anneal_npoints
[i
])
1106 sprintf(buf1
, "inputrec->grpopts.anneal_time[%d]", i
);
1107 sprintf(buf2
, "inputrec->grpopts.anneal_temp[%d]", i
);
1108 for (j
= 0; j
< opt1
->anneal_npoints
[i
]; j
++)
1110 cmp_real(fp
, buf1
, j
, opt1
->anneal_time
[i
][j
], opt2
->anneal_time
[i
][j
], ftol
, abstol
);
1111 cmp_real(fp
, buf2
, j
, opt1
->anneal_temp
[i
][j
], opt2
->anneal_temp
[i
][j
], ftol
, abstol
);
1115 if (opt1
->ngener
== opt2
->ngener
)
1117 for (i
= 0; i
< opt1
->ngener
; i
++)
1119 for (j
= i
; j
< opt1
->ngener
; j
++)
1121 sprintf(buf1
, "inputrec->grpopts.egp_flags[%d]", i
);
1122 cmp_int(fp
, buf1
, j
, opt1
->egp_flags
[opt1
->ngener
* i
+ j
],
1123 opt2
->egp_flags
[opt1
->ngener
* i
+ j
]);
1127 for (i
= 0; (i
< std::min(opt1
->ngacc
, opt2
->ngacc
)); i
++)
1129 cmp_rvec(fp
, "inputrec->grpopts.acc", i
, opt1
->acc
[i
], opt2
->acc
[i
], ftol
, abstol
);
1131 for (i
= 0; (i
< std::min(opt1
->ngfrz
, opt2
->ngfrz
)); i
++)
1133 cmp_ivec(fp
, "inputrec->grpopts.nFreeze", i
, opt1
->nFreeze
[i
], opt2
->nFreeze
[i
]);
1137 static void cmp_pull(FILE* fp
)
1140 "WARNING: Both files use COM pulling, but comparing of the pull struct is not "
1141 "implemented (yet). The pull parameters could be the same or different.\n");
1144 static void cmp_awhDimParams(FILE* fp
,
1145 const gmx::AwhDimParams
* dimp1
,
1146 const gmx::AwhDimParams
* dimp2
,
1151 /* Note that we have double index here, but the compare functions only
1152 * support one index, so here we only print the dim index and not the bias.
1154 cmp_int(fp
, "inputrec.awhParams->bias?->dim->coord_index", dimIndex
, dimp1
->coordIndex
,
1156 cmp_double(fp
, "inputrec->awhParams->bias?->dim->period", dimIndex
, dimp1
->period
,
1157 dimp2
->period
, ftol
, abstol
);
1158 cmp_double(fp
, "inputrec->awhParams->bias?->dim->diffusion", dimIndex
, dimp1
->diffusion
,
1159 dimp2
->diffusion
, ftol
, abstol
);
1160 cmp_double(fp
, "inputrec->awhParams->bias?->dim->origin", dimIndex
, dimp1
->origin
,
1161 dimp2
->origin
, ftol
, abstol
);
1162 cmp_double(fp
, "inputrec->awhParams->bias?->dim->end", dimIndex
, dimp1
->end
, dimp2
->end
, ftol
, abstol
);
1163 cmp_double(fp
, "inputrec->awhParams->bias?->dim->coord_value_init", dimIndex
,
1164 dimp1
->coordValueInit
, dimp2
->coordValueInit
, ftol
, abstol
);
1165 cmp_double(fp
, "inputrec->awhParams->bias?->dim->coverDiameter", dimIndex
, dimp1
->coverDiameter
,
1166 dimp2
->coverDiameter
, ftol
, abstol
);
1169 static void cmp_awhBiasParams(FILE* fp
,
1170 const gmx::AwhBiasParams
* bias1
,
1171 const gmx::AwhBiasParams
* bias2
,
1176 cmp_int(fp
, "inputrec->awhParams->ndim", biasIndex
, bias1
->ndim
, bias2
->ndim
);
1177 cmp_int(fp
, "inputrec->awhParams->biaseTarget", biasIndex
, bias1
->eTarget
, bias2
->eTarget
);
1178 cmp_double(fp
, "inputrec->awhParams->biastargetBetaScaling", biasIndex
,
1179 bias1
->targetBetaScaling
, bias2
->targetBetaScaling
, ftol
, abstol
);
1180 cmp_double(fp
, "inputrec->awhParams->biastargetCutoff", biasIndex
, bias1
->targetCutoff
,
1181 bias2
->targetCutoff
, ftol
, abstol
);
1182 cmp_int(fp
, "inputrec->awhParams->biaseGrowth", biasIndex
, bias1
->eGrowth
, bias2
->eGrowth
);
1183 cmp_bool(fp
, "inputrec->awhParams->biasbUserData", biasIndex
, bias1
->bUserData
!= 0,
1184 bias2
->bUserData
!= 0);
1185 cmp_double(fp
, "inputrec->awhParams->biaserror_initial", biasIndex
, bias1
->errorInitial
,
1186 bias2
->errorInitial
, ftol
, abstol
);
1187 cmp_int(fp
, "inputrec->awhParams->biasShareGroup", biasIndex
, bias1
->shareGroup
, bias2
->shareGroup
);
1189 for (int dim
= 0; dim
< std::min(bias1
->ndim
, bias2
->ndim
); dim
++)
1191 cmp_awhDimParams(fp
, &bias1
->dimParams
[dim
], &bias2
->dimParams
[dim
], dim
, ftol
, abstol
);
1195 static void cmp_awhParams(FILE* fp
, const gmx::AwhParams
* awh1
, const gmx::AwhParams
* awh2
, real ftol
, real abstol
)
1197 cmp_int(fp
, "inputrec->awhParams->nbias", -1, awh1
->numBias
, awh2
->numBias
);
1198 cmp_int64(fp
, "inputrec->awhParams->seed", awh1
->seed
, awh2
->seed
);
1199 cmp_int(fp
, "inputrec->awhParams->nstout", -1, awh1
->nstOut
, awh2
->nstOut
);
1200 cmp_int(fp
, "inputrec->awhParams->nstsample_coord", -1, awh1
->nstSampleCoord
, awh2
->nstSampleCoord
);
1201 cmp_int(fp
, "inputrec->awhParams->nsamples_update_free_energy", -1,
1202 awh1
->numSamplesUpdateFreeEnergy
, awh2
->numSamplesUpdateFreeEnergy
);
1203 cmp_int(fp
, "inputrec->awhParams->ePotential", -1, awh1
->ePotential
, awh2
->ePotential
);
1204 cmp_bool(fp
, "inputrec->awhParams->shareBiasMultisim", -1, awh1
->shareBiasMultisim
,
1205 awh2
->shareBiasMultisim
);
1207 if (awh1
->numBias
== awh2
->numBias
)
1209 for (int bias
= 0; bias
< awh1
->numBias
; bias
++)
1211 cmp_awhBiasParams(fp
, &awh1
->awhBiasParams
[bias
], &awh2
->awhBiasParams
[bias
], bias
, ftol
, abstol
);
1216 static void cmp_simtempvals(FILE* fp
,
1217 const t_simtemp
* simtemp1
,
1218 const t_simtemp
* simtemp2
,
1224 cmp_int(fp
, "inputrec->simtempvals->eSimTempScale", -1, simtemp1
->eSimTempScale
, simtemp2
->eSimTempScale
);
1225 cmp_real(fp
, "inputrec->simtempvals->simtemp_high", -1, simtemp1
->simtemp_high
,
1226 simtemp2
->simtemp_high
, ftol
, abstol
);
1227 cmp_real(fp
, "inputrec->simtempvals->simtemp_low", -1, simtemp1
->simtemp_low
,
1228 simtemp2
->simtemp_low
, ftol
, abstol
);
1229 for (i
= 0; i
< n_lambda
; i
++)
1231 cmp_real(fp
, "inputrec->simtempvals->temperatures", -1, simtemp1
->temperatures
[i
],
1232 simtemp2
->temperatures
[i
], ftol
, abstol
);
1236 static void cmp_expandedvals(FILE* fp
,
1237 const t_expanded
* expand1
,
1238 const t_expanded
* expand2
,
1245 cmp_bool(fp
, "inputrec->fepvals->bInit_weights", -1, expand1
->bInit_weights
, expand2
->bInit_weights
);
1246 cmp_bool(fp
, "inputrec->fepvals->bWLoneovert", -1, expand1
->bWLoneovert
, expand2
->bWLoneovert
);
1248 for (i
= 0; i
< n_lambda
; i
++)
1250 cmp_real(fp
, "inputrec->expandedvals->init_lambda_weights", -1,
1251 expand1
->init_lambda_weights
[i
], expand2
->init_lambda_weights
[i
], ftol
, abstol
);
1254 cmp_int(fp
, "inputrec->expandedvals->lambda-stats", -1, expand1
->elamstats
, expand2
->elamstats
);
1255 cmp_int(fp
, "inputrec->expandedvals->lambda-mc-move", -1, expand1
->elmcmove
, expand2
->elmcmove
);
1256 cmp_int(fp
, "inputrec->expandedvals->lmc-repeats", -1, expand1
->lmc_repeats
, expand2
->lmc_repeats
);
1257 cmp_int(fp
, "inputrec->expandedvals->lmc-gibbsdelta", -1, expand1
->gibbsdeltalam
, expand2
->gibbsdeltalam
);
1258 cmp_int(fp
, "inputrec->expandedvals->lmc-forced-nstart", -1, expand1
->lmc_forced_nstart
,
1259 expand2
->lmc_forced_nstart
);
1260 cmp_int(fp
, "inputrec->expandedvals->lambda-weights-equil", -1, expand1
->elmceq
, expand2
->elmceq
);
1261 cmp_int(fp
, "inputrec->expandedvals->,weight-equil-number-all-lambda", -1,
1262 expand1
->equil_n_at_lam
, expand2
->equil_n_at_lam
);
1263 cmp_int(fp
, "inputrec->expandedvals->weight-equil-number-samples", -1, expand1
->equil_samples
,
1264 expand2
->equil_samples
);
1265 cmp_int(fp
, "inputrec->expandedvals->weight-equil-number-steps", -1, expand1
->equil_steps
,
1266 expand2
->equil_steps
);
1267 cmp_real(fp
, "inputrec->expandedvals->weight-equil-wl-delta", -1, expand1
->equil_wl_delta
,
1268 expand2
->equil_wl_delta
, ftol
, abstol
);
1269 cmp_real(fp
, "inputrec->expandedvals->weight-equil-count-ratio", -1, expand1
->equil_ratio
,
1270 expand2
->equil_ratio
, ftol
, abstol
);
1271 cmp_bool(fp
, "inputrec->expandedvals->symmetrized-transition-matrix", -1,
1272 expand1
->bSymmetrizedTMatrix
, expand2
->bSymmetrizedTMatrix
);
1273 cmp_int(fp
, "inputrec->expandedvals->nstTij", -1, expand1
->nstTij
, expand2
->nstTij
);
1274 cmp_int(fp
, "inputrec->expandedvals->mininum-var-min", -1, expand1
->minvarmin
,
1275 expand2
->minvarmin
); /*default is reasonable */
1276 cmp_int(fp
, "inputrec->expandedvals->weight-c-range", -1, expand1
->c_range
, expand2
->c_range
); /* default is just C=0 */
1277 cmp_real(fp
, "inputrec->expandedvals->wl-scale", -1, expand1
->wl_scale
, expand2
->wl_scale
, ftol
, abstol
);
1278 cmp_real(fp
, "inputrec->expandedvals->init-wl-delta", -1, expand1
->init_wl_delta
,
1279 expand2
->init_wl_delta
, ftol
, abstol
);
1280 cmp_real(fp
, "inputrec->expandedvals->wl-ratio", -1, expand1
->wl_ratio
, expand2
->wl_ratio
, ftol
, abstol
);
1281 cmp_int(fp
, "inputrec->expandedvals->nstexpanded", -1, expand1
->nstexpanded
, expand2
->nstexpanded
);
1282 cmp_int(fp
, "inputrec->expandedvals->lmc-seed", -1, expand1
->lmc_seed
, expand2
->lmc_seed
);
1283 cmp_real(fp
, "inputrec->expandedvals->mc-temperature", -1, expand1
->mc_temp
, expand2
->mc_temp
,
1287 static void cmp_fepvals(FILE* fp
, const t_lambda
* fep1
, const t_lambda
* fep2
, real ftol
, real abstol
)
1290 cmp_int(fp
, "inputrec->nstdhdl", -1, fep1
->nstdhdl
, fep2
->nstdhdl
);
1291 cmp_double(fp
, "inputrec->fepvals->init_fep_state", -1, fep1
->init_fep_state
,
1292 fep2
->init_fep_state
, ftol
, abstol
);
1293 cmp_double(fp
, "inputrec->fepvals->delta_lambda", -1, fep1
->delta_lambda
, fep2
->delta_lambda
,
1295 cmp_int(fp
, "inputrec->fepvals->n_lambda", -1, fep1
->n_lambda
, fep2
->n_lambda
);
1296 for (i
= 0; i
< efptNR
; i
++)
1298 for (j
= 0; j
< std::min(fep1
->n_lambda
, fep2
->n_lambda
); j
++)
1300 cmp_double(fp
, "inputrec->fepvals->all_lambda", -1, fep1
->all_lambda
[i
][j
],
1301 fep2
->all_lambda
[i
][j
], ftol
, abstol
);
1304 cmp_int(fp
, "inputrec->fepvals->lambda_neighbors", 1, fep1
->lambda_neighbors
, fep2
->lambda_neighbors
);
1305 cmp_real(fp
, "inputrec->fepvals->sc_alpha", -1, fep1
->sc_alpha
, fep2
->sc_alpha
, ftol
, abstol
);
1306 cmp_int(fp
, "inputrec->fepvals->sc_power", -1, fep1
->sc_power
, fep2
->sc_power
);
1307 cmp_real(fp
, "inputrec->fepvals->sc_r_power", -1, fep1
->sc_r_power
, fep2
->sc_r_power
, ftol
, abstol
);
1308 cmp_real(fp
, "inputrec->fepvals->sc_sigma", -1, fep1
->sc_sigma
, fep2
->sc_sigma
, ftol
, abstol
);
1309 cmp_int(fp
, "inputrec->fepvals->edHdLPrintEnergy", -1, fep1
->edHdLPrintEnergy
, fep1
->edHdLPrintEnergy
);
1310 cmp_bool(fp
, "inputrec->fepvals->bScCoul", -1, fep1
->bScCoul
, fep1
->bScCoul
);
1311 cmp_int(fp
, "inputrec->separate_dhdl_file", -1, fep1
->separate_dhdl_file
, fep2
->separate_dhdl_file
);
1312 cmp_int(fp
, "inputrec->dhdl_derivatives", -1, fep1
->dhdl_derivatives
, fep2
->dhdl_derivatives
);
1313 cmp_int(fp
, "inputrec->dh_hist_size", -1, fep1
->dh_hist_size
, fep2
->dh_hist_size
);
1314 cmp_double(fp
, "inputrec->dh_hist_spacing", -1, fep1
->dh_hist_spacing
, fep2
->dh_hist_spacing
,
1318 void cmp_inputrec(FILE* fp
, const t_inputrec
* ir1
, const t_inputrec
* ir2
, real ftol
, real abstol
)
1320 fprintf(fp
, "comparing inputrec\n");
1322 /* gcc 2.96 doesnt like these defines at all, but issues a huge list
1323 * of warnings. Maybe it will change in future versions, but for the
1324 * moment I've spelled them out instead. /EL 000820
1325 * #define CIB(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
1326 * #define CII(s) cmp_int(fp,"inputrec->"#s,0,ir1->##s,ir2->##s)
1327 * #define CIR(s) cmp_real(fp,"inputrec->"#s,0,ir1->##s,ir2->##s,ftol)
1329 cmp_int(fp
, "inputrec->eI", -1, ir1
->eI
, ir2
->eI
);
1330 cmp_int64(fp
, "inputrec->nsteps", ir1
->nsteps
, ir2
->nsteps
);
1331 cmp_int64(fp
, "inputrec->init_step", ir1
->init_step
, ir2
->init_step
);
1332 cmp_int(fp
, "inputrec->simulation_part", -1, ir1
->simulation_part
, ir2
->simulation_part
);
1333 cmp_int(fp
, "inputrec->mts", -1, static_cast<int>(ir1
->useMts
), static_cast<int>(ir2
->useMts
));
1334 if (ir1
->useMts
&& ir2
->useMts
)
1336 cmp_int(fp
, "inputrec->mts-levels", -1, ir1
->mtsLevels
.size(), ir2
->mtsLevels
.size());
1337 cmp_int(fp
, "inputrec->mts-level2-forces", -1, ir1
->mtsLevels
[1].forceGroups
.to_ulong(),
1338 ir2
->mtsLevels
[1].forceGroups
.to_ulong());
1339 cmp_int(fp
, "inputrec->mts-level2-factor", -1, ir1
->mtsLevels
[1].stepFactor
,
1340 ir2
->mtsLevels
[1].stepFactor
);
1342 cmp_int(fp
, "inputrec->pbcType", -1, static_cast<int>(ir1
->pbcType
), static_cast<int>(ir2
->pbcType
));
1343 cmp_bool(fp
, "inputrec->bPeriodicMols", -1, ir1
->bPeriodicMols
, ir2
->bPeriodicMols
);
1344 cmp_int(fp
, "inputrec->cutoff_scheme", -1, ir1
->cutoff_scheme
, ir2
->cutoff_scheme
);
1345 cmp_int(fp
, "inputrec->nstlist", -1, ir1
->nstlist
, ir2
->nstlist
);
1346 cmp_int(fp
, "inputrec->nstcomm", -1, ir1
->nstcomm
, ir2
->nstcomm
);
1347 cmp_int(fp
, "inputrec->comm_mode", -1, ir1
->comm_mode
, ir2
->comm_mode
);
1348 cmp_int(fp
, "inputrec->nstlog", -1, ir1
->nstlog
, ir2
->nstlog
);
1349 cmp_int(fp
, "inputrec->nstxout", -1, ir1
->nstxout
, ir2
->nstxout
);
1350 cmp_int(fp
, "inputrec->nstvout", -1, ir1
->nstvout
, ir2
->nstvout
);
1351 cmp_int(fp
, "inputrec->nstfout", -1, ir1
->nstfout
, ir2
->nstfout
);
1352 cmp_int(fp
, "inputrec->nstcalcenergy", -1, ir1
->nstcalcenergy
, ir2
->nstcalcenergy
);
1353 cmp_int(fp
, "inputrec->nstenergy", -1, ir1
->nstenergy
, ir2
->nstenergy
);
1354 cmp_int(fp
, "inputrec->nstxout_compressed", -1, ir1
->nstxout_compressed
, ir2
->nstxout_compressed
);
1355 cmp_double(fp
, "inputrec->init_t", -1, ir1
->init_t
, ir2
->init_t
, ftol
, abstol
);
1356 cmp_double(fp
, "inputrec->delta_t", -1, ir1
->delta_t
, ir2
->delta_t
, ftol
, abstol
);
1357 cmp_real(fp
, "inputrec->x_compression_precision", -1, ir1
->x_compression_precision
,
1358 ir2
->x_compression_precision
, ftol
, abstol
);
1359 cmp_real(fp
, "inputrec->fourierspacing", -1, ir1
->fourier_spacing
, ir2
->fourier_spacing
, ftol
, abstol
);
1360 cmp_int(fp
, "inputrec->nkx", -1, ir1
->nkx
, ir2
->nkx
);
1361 cmp_int(fp
, "inputrec->nky", -1, ir1
->nky
, ir2
->nky
);
1362 cmp_int(fp
, "inputrec->nkz", -1, ir1
->nkz
, ir2
->nkz
);
1363 cmp_int(fp
, "inputrec->pme_order", -1, ir1
->pme_order
, ir2
->pme_order
);
1364 cmp_real(fp
, "inputrec->ewald_rtol", -1, ir1
->ewald_rtol
, ir2
->ewald_rtol
, ftol
, abstol
);
1365 cmp_int(fp
, "inputrec->ewald_geometry", -1, ir1
->ewald_geometry
, ir2
->ewald_geometry
);
1366 cmp_real(fp
, "inputrec->epsilon_surface", -1, ir1
->epsilon_surface
, ir2
->epsilon_surface
, ftol
, abstol
);
1367 cmp_int(fp
, "inputrec->bContinuation", -1, static_cast<int>(ir1
->bContinuation
),
1368 static_cast<int>(ir2
->bContinuation
));
1369 cmp_int(fp
, "inputrec->bShakeSOR", -1, static_cast<int>(ir1
->bShakeSOR
),
1370 static_cast<int>(ir2
->bShakeSOR
));
1371 cmp_int(fp
, "inputrec->etc", -1, ir1
->etc
, ir2
->etc
);
1372 cmp_int(fp
, "inputrec->bPrintNHChains", -1, static_cast<int>(ir1
->bPrintNHChains
),
1373 static_cast<int>(ir2
->bPrintNHChains
));
1374 cmp_int(fp
, "inputrec->epc", -1, ir1
->epc
, ir2
->epc
);
1375 cmp_int(fp
, "inputrec->epct", -1, ir1
->epct
, ir2
->epct
);
1376 cmp_real(fp
, "inputrec->tau_p", -1, ir1
->tau_p
, ir2
->tau_p
, ftol
, abstol
);
1377 cmp_rvec(fp
, "inputrec->ref_p(x)", -1, ir1
->ref_p
[XX
], ir2
->ref_p
[XX
], ftol
, abstol
);
1378 cmp_rvec(fp
, "inputrec->ref_p(y)", -1, ir1
->ref_p
[YY
], ir2
->ref_p
[YY
], ftol
, abstol
);
1379 cmp_rvec(fp
, "inputrec->ref_p(z)", -1, ir1
->ref_p
[ZZ
], ir2
->ref_p
[ZZ
], ftol
, abstol
);
1380 cmp_rvec(fp
, "inputrec->compress(x)", -1, ir1
->compress
[XX
], ir2
->compress
[XX
], ftol
, abstol
);
1381 cmp_rvec(fp
, "inputrec->compress(y)", -1, ir1
->compress
[YY
], ir2
->compress
[YY
], ftol
, abstol
);
1382 cmp_rvec(fp
, "inputrec->compress(z)", -1, ir1
->compress
[ZZ
], ir2
->compress
[ZZ
], ftol
, abstol
);
1383 cmp_int(fp
, "refcoord_scaling", -1, ir1
->refcoord_scaling
, ir2
->refcoord_scaling
);
1384 cmp_rvec(fp
, "inputrec->posres_com", -1, ir1
->posres_com
, ir2
->posres_com
, ftol
, abstol
);
1385 cmp_rvec(fp
, "inputrec->posres_comB", -1, ir1
->posres_comB
, ir2
->posres_comB
, ftol
, abstol
);
1386 cmp_real(fp
, "inputrec->verletbuf_tol", -1, ir1
->verletbuf_tol
, ir2
->verletbuf_tol
, ftol
, abstol
);
1387 cmp_real(fp
, "inputrec->rlist", -1, ir1
->rlist
, ir2
->rlist
, ftol
, abstol
);
1388 cmp_real(fp
, "inputrec->rtpi", -1, ir1
->rtpi
, ir2
->rtpi
, ftol
, abstol
);
1389 cmp_int(fp
, "inputrec->coulombtype", -1, ir1
->coulombtype
, ir2
->coulombtype
);
1390 cmp_int(fp
, "inputrec->coulomb_modifier", -1, ir1
->coulomb_modifier
, ir2
->coulomb_modifier
);
1391 cmp_real(fp
, "inputrec->rcoulomb_switch", -1, ir1
->rcoulomb_switch
, ir2
->rcoulomb_switch
, ftol
, abstol
);
1392 cmp_real(fp
, "inputrec->rcoulomb", -1, ir1
->rcoulomb
, ir2
->rcoulomb
, ftol
, abstol
);
1393 cmp_int(fp
, "inputrec->vdwtype", -1, ir1
->vdwtype
, ir2
->vdwtype
);
1394 cmp_int(fp
, "inputrec->vdw_modifier", -1, ir1
->vdw_modifier
, ir2
->vdw_modifier
);
1395 cmp_real(fp
, "inputrec->rvdw_switch", -1, ir1
->rvdw_switch
, ir2
->rvdw_switch
, ftol
, abstol
);
1396 cmp_real(fp
, "inputrec->rvdw", -1, ir1
->rvdw
, ir2
->rvdw
, ftol
, abstol
);
1397 cmp_real(fp
, "inputrec->epsilon_r", -1, ir1
->epsilon_r
, ir2
->epsilon_r
, ftol
, abstol
);
1398 cmp_real(fp
, "inputrec->epsilon_rf", -1, ir1
->epsilon_rf
, ir2
->epsilon_rf
, ftol
, abstol
);
1399 cmp_real(fp
, "inputrec->tabext", -1, ir1
->tabext
, ir2
->tabext
, ftol
, abstol
);
1401 cmp_int(fp
, "inputrec->eDispCorr", -1, ir1
->eDispCorr
, ir2
->eDispCorr
);
1402 cmp_real(fp
, "inputrec->shake_tol", -1, ir1
->shake_tol
, ir2
->shake_tol
, ftol
, abstol
);
1403 cmp_int(fp
, "inputrec->efep", -1, ir1
->efep
, ir2
->efep
);
1404 cmp_fepvals(fp
, ir1
->fepvals
, ir2
->fepvals
, ftol
, abstol
);
1405 cmp_int(fp
, "inputrec->bSimTemp", -1, static_cast<int>(ir1
->bSimTemp
), static_cast<int>(ir2
->bSimTemp
));
1406 if ((ir1
->bSimTemp
== ir2
->bSimTemp
) && (ir1
->bSimTemp
))
1408 cmp_simtempvals(fp
, ir1
->simtempvals
, ir2
->simtempvals
,
1409 std::min(ir1
->fepvals
->n_lambda
, ir2
->fepvals
->n_lambda
), ftol
, abstol
);
1411 cmp_int(fp
, "inputrec->bExpanded", -1, static_cast<int>(ir1
->bExpanded
),
1412 static_cast<int>(ir2
->bExpanded
));
1413 if ((ir1
->bExpanded
== ir2
->bExpanded
) && (ir1
->bExpanded
))
1415 cmp_expandedvals(fp
, ir1
->expandedvals
, ir2
->expandedvals
,
1416 std::min(ir1
->fepvals
->n_lambda
, ir2
->fepvals
->n_lambda
), ftol
, abstol
);
1418 cmp_int(fp
, "inputrec->nwall", -1, ir1
->nwall
, ir2
->nwall
);
1419 cmp_int(fp
, "inputrec->wall_type", -1, ir1
->wall_type
, ir2
->wall_type
);
1420 cmp_int(fp
, "inputrec->wall_atomtype[0]", -1, ir1
->wall_atomtype
[0], ir2
->wall_atomtype
[0]);
1421 cmp_int(fp
, "inputrec->wall_atomtype[1]", -1, ir1
->wall_atomtype
[1], ir2
->wall_atomtype
[1]);
1422 cmp_real(fp
, "inputrec->wall_density[0]", -1, ir1
->wall_density
[0], ir2
->wall_density
[0], ftol
, abstol
);
1423 cmp_real(fp
, "inputrec->wall_density[1]", -1, ir1
->wall_density
[1], ir2
->wall_density
[1], ftol
, abstol
);
1424 cmp_real(fp
, "inputrec->wall_ewald_zfac", -1, ir1
->wall_ewald_zfac
, ir2
->wall_ewald_zfac
, ftol
, abstol
);
1426 cmp_bool(fp
, "inputrec->bPull", -1, ir1
->bPull
, ir2
->bPull
);
1427 if (ir1
->bPull
&& ir2
->bPull
)
1432 cmp_bool(fp
, "inputrec->bDoAwh", -1, ir1
->bDoAwh
, ir2
->bDoAwh
);
1433 if (ir1
->bDoAwh
&& ir2
->bDoAwh
)
1435 cmp_awhParams(fp
, ir1
->awhParams
, ir2
->awhParams
, ftol
, abstol
);
1438 cmp_int(fp
, "inputrec->eDisre", -1, ir1
->eDisre
, ir2
->eDisre
);
1439 cmp_real(fp
, "inputrec->dr_fc", -1, ir1
->dr_fc
, ir2
->dr_fc
, ftol
, abstol
);
1440 cmp_int(fp
, "inputrec->eDisreWeighting", -1, ir1
->eDisreWeighting
, ir2
->eDisreWeighting
);
1441 cmp_int(fp
, "inputrec->bDisreMixed", -1, static_cast<int>(ir1
->bDisreMixed
),
1442 static_cast<int>(ir2
->bDisreMixed
));
1443 cmp_int(fp
, "inputrec->nstdisreout", -1, ir1
->nstdisreout
, ir2
->nstdisreout
);
1444 cmp_real(fp
, "inputrec->dr_tau", -1, ir1
->dr_tau
, ir2
->dr_tau
, ftol
, abstol
);
1445 cmp_real(fp
, "inputrec->orires_fc", -1, ir1
->orires_fc
, ir2
->orires_fc
, ftol
, abstol
);
1446 cmp_real(fp
, "inputrec->orires_tau", -1, ir1
->orires_tau
, ir2
->orires_tau
, ftol
, abstol
);
1447 cmp_int(fp
, "inputrec->nstorireout", -1, ir1
->nstorireout
, ir2
->nstorireout
);
1448 cmp_real(fp
, "inputrec->em_stepsize", -1, ir1
->em_stepsize
, ir2
->em_stepsize
, ftol
, abstol
);
1449 cmp_real(fp
, "inputrec->em_tol", -1, ir1
->em_tol
, ir2
->em_tol
, ftol
, abstol
);
1450 cmp_int(fp
, "inputrec->niter", -1, ir1
->niter
, ir2
->niter
);
1451 cmp_real(fp
, "inputrec->fc_stepsize", -1, ir1
->fc_stepsize
, ir2
->fc_stepsize
, ftol
, abstol
);
1452 cmp_int(fp
, "inputrec->nstcgsteep", -1, ir1
->nstcgsteep
, ir2
->nstcgsteep
);
1453 cmp_int(fp
, "inputrec->nbfgscorr", 0, ir1
->nbfgscorr
, ir2
->nbfgscorr
);
1454 cmp_int(fp
, "inputrec->eConstrAlg", -1, ir1
->eConstrAlg
, ir2
->eConstrAlg
);
1455 cmp_int(fp
, "inputrec->nProjOrder", -1, ir1
->nProjOrder
, ir2
->nProjOrder
);
1456 cmp_real(fp
, "inputrec->LincsWarnAngle", -1, ir1
->LincsWarnAngle
, ir2
->LincsWarnAngle
, ftol
, abstol
);
1457 cmp_int(fp
, "inputrec->nLincsIter", -1, ir1
->nLincsIter
, ir2
->nLincsIter
);
1458 cmp_real(fp
, "inputrec->bd_fric", -1, ir1
->bd_fric
, ir2
->bd_fric
, ftol
, abstol
);
1459 cmp_int64(fp
, "inputrec->ld_seed", ir1
->ld_seed
, ir2
->ld_seed
);
1460 cmp_real(fp
, "inputrec->cos_accel", -1, ir1
->cos_accel
, ir2
->cos_accel
, ftol
, abstol
);
1461 cmp_rvec(fp
, "inputrec->deform(a)", -1, ir1
->deform
[XX
], ir2
->deform
[XX
], ftol
, abstol
);
1462 cmp_rvec(fp
, "inputrec->deform(b)", -1, ir1
->deform
[YY
], ir2
->deform
[YY
], ftol
, abstol
);
1463 cmp_rvec(fp
, "inputrec->deform(c)", -1, ir1
->deform
[ZZ
], ir2
->deform
[ZZ
], ftol
, abstol
);
1466 cmp_int(fp
, "inputrec->userint1", -1, ir1
->userint1
, ir2
->userint1
);
1467 cmp_int(fp
, "inputrec->userint2", -1, ir1
->userint2
, ir2
->userint2
);
1468 cmp_int(fp
, "inputrec->userint3", -1, ir1
->userint3
, ir2
->userint3
);
1469 cmp_int(fp
, "inputrec->userint4", -1, ir1
->userint4
, ir2
->userint4
);
1470 cmp_real(fp
, "inputrec->userreal1", -1, ir1
->userreal1
, ir2
->userreal1
, ftol
, abstol
);
1471 cmp_real(fp
, "inputrec->userreal2", -1, ir1
->userreal2
, ir2
->userreal2
, ftol
, abstol
);
1472 cmp_real(fp
, "inputrec->userreal3", -1, ir1
->userreal3
, ir2
->userreal3
, ftol
, abstol
);
1473 cmp_real(fp
, "inputrec->userreal4", -1, ir1
->userreal4
, ir2
->userreal4
, ftol
, abstol
);
1474 cmp_grpopts(fp
, &(ir1
->opts
), &(ir2
->opts
), ftol
, abstol
);
1475 gmx::TextWriter
writer(fp
);
1476 gmx::compareKeyValueTrees(&writer
, *ir1
->params
, *ir2
->params
, ftol
, abstol
);
1479 void comp_pull_AB(FILE* fp
, pull_params_t
* pull
, real ftol
, real abstol
)
1483 for (i
= 0; i
< pull
->ncoord
; i
++)
1485 fprintf(fp
, "comparing pull coord %d\n", i
);
1486 cmp_real(fp
, "pull-coord->k", -1, pull
->coord
[i
].k
, pull
->coord
[i
].kB
, ftol
, abstol
);
1490 gmx_bool
inputrecDeform(const t_inputrec
* ir
)
1492 return (ir
->deform
[XX
][XX
] != 0 || ir
->deform
[YY
][YY
] != 0 || ir
->deform
[ZZ
][ZZ
] != 0
1493 || ir
->deform
[YY
][XX
] != 0 || ir
->deform
[ZZ
][XX
] != 0 || ir
->deform
[ZZ
][YY
] != 0);
1496 gmx_bool
inputrecDynamicBox(const t_inputrec
* ir
)
1498 return (ir
->epc
!= epcNO
|| ir
->eI
== eiTPI
|| inputrecDeform(ir
));
1501 gmx_bool
inputrecPreserveShape(const t_inputrec
* ir
)
1503 return (ir
->epc
!= epcNO
&& ir
->deform
[XX
][XX
] == 0
1504 && (ir
->epct
== epctISOTROPIC
|| ir
->epct
== epctSEMIISOTROPIC
));
1507 gmx_bool
inputrecNeedMutot(const t_inputrec
* ir
)
1509 return ((ir
->coulombtype
== eelEWALD
|| EEL_PME(ir
->coulombtype
))
1510 && (ir
->ewald_geometry
== eewg3DC
|| ir
->epsilon_surface
!= 0));
1513 gmx_bool
inputrecExclForces(const t_inputrec
* ir
)
1515 return (EEL_FULL(ir
->coulombtype
) || (EEL_RF(ir
->coulombtype
)));
1518 gmx_bool
inputrecNptTrotter(const t_inputrec
* ir
)
1520 return (((ir
->eI
== eiVV
) || (ir
->eI
== eiVVAK
)) && (ir
->epc
== epcMTTK
) && (ir
->etc
== etcNOSEHOOVER
));
1523 gmx_bool
inputrecNvtTrotter(const t_inputrec
* ir
)
1525 return (((ir
->eI
== eiVV
) || (ir
->eI
== eiVVAK
)) && (ir
->epc
!= epcMTTK
) && (ir
->etc
== etcNOSEHOOVER
));
1528 gmx_bool
inputrecNphTrotter(const t_inputrec
* ir
)
1530 return (((ir
->eI
== eiVV
) || (ir
->eI
== eiVVAK
)) && (ir
->epc
== epcMTTK
) && (ir
->etc
!= etcNOSEHOOVER
));
1533 bool inputrecPbcXY2Walls(const t_inputrec
* ir
)
1535 return (ir
->pbcType
== PbcType::XY
&& ir
->nwall
== 2);
1538 bool integratorHasConservedEnergyQuantity(const t_inputrec
* ir
)
1541 { // NOLINT bugprone-branch-clone
1542 // Energy minimization or stochastic integrator: no conservation
1545 else if (ir
->etc
== etcNO
&& ir
->epc
== epcNO
)
1547 // The total energy is conserved, no additional conserved quanitity
1552 // Shear stress with Parrinello-Rahman is not supported (tedious)
1554 ((ir
->epc
== epcPARRINELLORAHMAN
|| ir
->epc
== epcMTTK
)
1555 && (ir
->ref_p
[YY
][XX
] != 0 || ir
->ref_p
[ZZ
][XX
] != 0 || ir
->ref_p
[ZZ
][YY
] != 0));
1557 return !ETC_ANDERSEN(ir
->etc
) && !shearWithPR
;
1561 bool integratorHasReferenceTemperature(const t_inputrec
* ir
)
1563 return ((ir
->etc
!= etcNO
) || EI_SD(ir
->eI
) || (ir
->eI
== eiBD
) || EI_TPI(ir
->eI
));
1566 int inputrec2nboundeddim(const t_inputrec
* ir
)
1568 if (inputrecPbcXY2Walls(ir
))
1574 return numPbcDimensions(ir
->pbcType
);
1578 int ndof_com(const t_inputrec
* ir
)
1582 switch (ir
->pbcType
)
1585 case PbcType::No
: n
= 3; break;
1586 case PbcType::XY
: n
= (ir
->nwall
== 0 ? 3 : 2); break;
1587 case PbcType::Screw
: n
= 1; break;
1588 default: gmx_incons("Unknown pbc in calc_nrdf");
1594 real
maxReferenceTemperature(const t_inputrec
& ir
)
1596 if (EI_ENERGY_MINIMIZATION(ir
.eI
) || ir
.eI
== eiNM
)
1601 if (EI_MD(ir
.eI
) && ir
.etc
== etcNO
)
1606 /* SD and BD also use ref_t and tau_t for setting the reference temperature.
1607 * TPI can be treated as MD, since it needs an ensemble temperature.
1610 real maxTemperature
= 0;
1611 for (int i
= 0; i
< ir
.opts
.ngtc
; i
++)
1613 if (ir
.opts
.tau_t
[i
] >= 0)
1615 maxTemperature
= std::max(maxTemperature
, ir
.opts
.ref_t
[i
]);
1619 return maxTemperature
;
1622 bool haveEwaldSurfaceContribution(const t_inputrec
& ir
)
1624 return EEL_PME_EWALD(ir
.coulombtype
) && (ir
.ewald_geometry
== eewg3DC
|| ir
.epsilon_surface
!= 0);
1627 bool haveFreeEnergyType(const t_inputrec
& ir
, const int fepType
)
1629 for (int i
= 0; i
< ir
.fepvals
->n_lambda
; i
++)
1631 if (ir
.fepvals
->all_lambda
[fepType
][i
] > 0)