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-2004, The GROMACS development team.
6 * Copyright (c) 2013,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.
40 #include "read_params.h"
42 #include "gromacs/applied_forces/awh/awh.h"
43 #include "gromacs/fileio/readinp.h"
44 #include "gromacs/fileio/warninp.h"
45 #include "gromacs/math/units.h"
46 #include "gromacs/math/utilities.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/mdtypes/awh_params.h"
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/mdtypes/md_enums.h"
51 #include "gromacs/mdtypes/multipletimestepping.h"
52 #include "gromacs/mdtypes/pull_params.h"
53 #include "gromacs/pbcutil/pbc.h"
54 #include "gromacs/pulling/pull.h"
55 #include "gromacs/random/seed.h"
56 #include "gromacs/topology/mtop_util.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/smalloc.h"
60 #include "gromacs/utility/stringutil.h"
62 #include "biasparams.h"
63 #include "biassharing.h"
68 const char* eawhtarget_names[eawhtargetNR + 1] = { "constant", "cutoff", "boltzmann",
69 "local-boltzmann", nullptr };
71 const char* eawhgrowth_names[eawhgrowthNR + 1] = { "exp-linear", "linear", nullptr };
73 const char* eawhpotential_names[eawhpotentialNR + 1] = { "convolved", "umbrella", nullptr };
75 const char* eawhcoordprovider_names[eawhcoordproviderNR + 1] = { "pull", "fep-lambda", nullptr };
81 * Check multiple time-stepping consistency between AWH and pull and/or FEP
83 * \param[in] inputrec Inputput parameter struct.
84 * \param[in,out] wi Struct for bookeeping warnings.
86 void checkMtsConsistency(const t_inputrec& inputrec, warninp_t wi)
93 GMX_RELEASE_ASSERT(inputrec.mtsLevels.size() == 2, "Only 2 MTS levels supported here");
95 bool usesPull = false;
97 for (int b = 0; b < inputrec.awhParams->numBias; b++)
99 const auto& biasParams = inputrec.awhParams->awhBiasParams[b];
100 for (int d = 0; d < biasParams.ndim; d++)
102 switch (biasParams.dimParams[d].eCoordProvider)
104 case eawhcoordproviderPULL: usesPull = true; break;
105 case eawhcoordproviderFREE_ENERGY_LAMBDA: usesFep = true; break;
106 default: GMX_RELEASE_ASSERT(false, "Unsupported coord provider");
110 const int awhMtsLevel = forceGroupMtsLevel(inputrec.mtsLevels, MtsForceGroups::Awh);
111 if (usesPull && forceGroupMtsLevel(inputrec.mtsLevels, MtsForceGroups::Pull) != awhMtsLevel)
114 "When AWH is applied to pull coordinates, pull and AWH should be computed at "
115 "the same MTS level");
117 if (usesFep && awhMtsLevel != ssize(inputrec.mtsLevels) - 1)
120 "When AWH is applied to the free-energy lambda with MTS, AWH should be "
121 "computed at the slow MTS level");
124 if (inputrec.awhParams->nstSampleCoord % inputrec.mtsLevels[awhMtsLevel].stepFactor != 0)
127 "With MTS applied to AWH, awh-nstsample should be a multiple of mts-factor");
132 * Check the parameters of an AWH bias pull dimension.
134 * \param[in] prefix Prefix for dimension parameters.
135 * \param[in,out] dimParams AWH dimensional parameters.
136 * \param[in] pull_params Pull parameters.
137 * \param[in,out] wi Struct for bookeeping warnings.
139 void checkPullDimParams(const std::string& prefix,
140 AwhDimParams* dimParams,
141 const pull_params_t& pull_params,
144 if (dimParams->coordIndex < 0)
147 "Failed to read a valid coordinate index for %s-coord-index. "
148 "Note that the pull coordinate indexing starts at 1.",
151 if (dimParams->coordIndex >= pull_params.ncoord)
154 "The given AWH coordinate index (%d) is larger than the number of pull "
156 dimParams->coordIndex + 1, pull_params.ncoord);
158 if (pull_params.coord[dimParams->coordIndex].rate != 0)
160 auto message = formatString(
161 "Setting pull-coord%d-rate (%g) is incompatible with AWH biasing this coordinate",
162 dimParams->coordIndex + 1, pull_params.coord[dimParams->coordIndex].rate);
163 warning_error(wi, message);
166 if (gmx_within_tol(dimParams->end - dimParams->origin, 0, GMX_REAL_EPS))
168 auto message = formatString(
169 "The given interval length given by %s-start (%g) and %s-end (%g) is zero. "
170 "This will result in only one point along this axis in the coordinate value grid.",
171 prefix.c_str(), dimParams->origin, prefix.c_str(), dimParams->end);
172 warning(wi, message);
175 if (dimParams->forceConstant <= 0)
177 warning_error(wi, "The force AWH bias force constant should be > 0");
180 /* Grid params for each axis */
181 int eGeom = pull_params.coord[dimParams->coordIndex].eGeom;
183 /* Check that the requested interval is in allowed range */
184 if (eGeom == epullgDIST)
186 if (dimParams->origin < 0 || dimParams->end < 0)
189 "%s-start (%g) or %s-end (%g) set to a negative value. With pull "
190 "geometry distance coordinate values are non-negative. "
191 "Perhaps you want to use geometry %s instead?",
192 prefix.c_str(), dimParams->origin, prefix.c_str(), dimParams->end,
193 EPULLGEOM(epullgDIR));
196 else if (eGeom == epullgANGLE || eGeom == epullgANGLEAXIS)
198 if (dimParams->origin < 0 || dimParams->end > 180)
201 "%s-start (%g) and %s-end (%g) are outside of the allowed range "
202 "0 to 180 deg for pull geometries %s and %s ",
203 prefix.c_str(), dimParams->origin, prefix.c_str(), dimParams->end,
204 EPULLGEOM(epullgANGLE), EPULLGEOM(epullgANGLEAXIS));
207 else if (eGeom == epullgDIHEDRAL)
209 if (dimParams->origin < -180 || dimParams->end > 180)
212 "%s-start (%g) and %s-end (%g) are outside of the allowed range "
213 "-180 to 180 deg for pull geometry %s. ",
214 prefix.c_str(), dimParams->origin, prefix.c_str(), dimParams->end,
215 EPULLGEOM(epullgDIHEDRAL));
221 * Check parameters of an AWH bias in a free energy lambda state dimension.
223 * \param[in] prefix Prefix for dimension parameters.
224 * \param[in,out] dimParams AWH dimensional parameters.
225 * \param[in] lambdaParams The free energy lambda related parameters.
226 * \param[in] efep This is the type of FEP calculation (efep enumerator).
227 * \param[in,out] wi Struct for bookeeping warnings.
229 void checkFepLambdaDimParams(const std::string& prefix,
230 const AwhDimParams* dimParams,
231 const t_lambda* lambdaParams,
240 "There must be free energy input if using AWH to steer the free energy lambda "
244 if (lambdaParams->lambda_neighbors != -1)
247 "When running AWH coupled to the free energy lambda state all lambda states "
248 "should be used as neighbors in order to get correct probabilities, i.e. "
249 "calc-lambda-neighbors (%d) must be %d.",
250 lambdaParams->lambda_neighbors, -1);
253 if (efep == efepSLOWGROWTH || lambdaParams->delta_lambda != 0)
256 "AWH coupled to the free energy lambda state is not compatible with slow-growth "
257 "and delta-lambda must be 0.");
260 if (efep == efepEXPANDED)
263 "AWH is not treated like other expanded ensemble methods. Do not use expanded.");
266 if (dimParams->origin < 0)
268 opt = prefix + "-start";
270 "When running AWH coupled to the free energy lambda state the lower lambda state "
271 "for AWH, %s (%.0f), must be >= 0.",
272 opt.c_str(), dimParams->origin);
274 if (dimParams->end >= lambdaParams->n_lambda)
276 opt = prefix + "-end";
278 "When running AWH coupled to the free energy lambda state the upper lambda state "
279 "for AWH, %s (%.0f), must be < n_lambda (%d).",
280 opt.c_str(), dimParams->origin, lambdaParams->n_lambda);
282 if (gmx_within_tol(dimParams->end - dimParams->origin, 0, GMX_REAL_EPS))
284 auto message = formatString(
285 "The given interval length given by %s-start (%g) and %s-end (%g) is zero. "
286 "This will result in only one lambda point along this free energy lambda state "
287 "axis in the coordinate value grid.",
288 prefix.c_str(), dimParams->origin, prefix.c_str(), dimParams->end);
289 warning(wi, message);
292 if (dimParams->forceConstant != 0)
296 "The force AWH bias force constant is not used with free energy lambda state as "
297 "coordinate provider.");
302 * Check that AWH FEP is not combined with incompatible decoupling.
304 * \param[in] mtop System topology.
305 * \param[in,out] wi Struct for bookeeping warnings.
307 void checkFepLambdaDimDecouplingConsistency(const gmx_mtop_t& mtop, warninp_t wi)
309 if (haveFepPerturbedMasses(mtop))
312 "Masses may not be perturbed when using the free energy lambda state as AWH "
313 "coordinate provider. If you are using fep-lambdas to specify lambda states "
314 "make sure that you also specify mass-lambdas without perturbation.");
316 if (havePerturbedConstraints(mtop))
319 "Constraints may not be perturbed when using the free energy lambda state as AWH "
320 "coordinate provider. If you are using fep-lambdas to specify lambda states "
321 "make sure that you also specify mass-lambdas without perturbation.");
326 * Read parameters of an AWH bias dimension.
328 * \param[in,out] inp Input file entries.
329 * \param[in] prefix Prefix for dimension parameters.
330 * \param[in,out] dimParams AWH dimensional parameters.
331 * \param[in,out] wi Struct for bookeeping warnings.
332 * \param[in] bComment True if comments should be printed.
334 void readDimParams(std::vector<t_inpfile>* inp,
335 const std::string& prefix,
336 AwhDimParams* dimParams,
343 printStringNoNewline(
345 "The provider of the reaction coordinate, "
346 "currently only 'pull' and 'fep-lambda' (free energy lambda state) is supported");
348 opt = prefix + "-coord-provider";
349 dimParams->eCoordProvider = get_eeenum(inp, opt, eawhcoordprovider_names, wi);
353 printStringNoNewline(inp, "The coordinate index for this dimension");
355 opt = prefix + "-coord-index";
357 coordIndexInput = get_eint(inp, opt, 1, wi);
359 /* The pull coordinate indices start at 1 in the input file, at 0 internally */
360 dimParams->coordIndex = coordIndexInput - 1;
364 printStringNoNewline(inp, "Start and end values for each coordinate dimension");
366 opt = prefix + "-start";
367 dimParams->origin = get_ereal(inp, opt, 0., wi);
368 opt = prefix + "-end";
369 dimParams->end = get_ereal(inp, opt, 0., wi);
373 printStringNoNewline(
374 inp, "The force constant for this coordinate (kJ/mol/nm^2 or kJ/mol/rad^2)");
376 opt = prefix + "-force-constant";
377 dimParams->forceConstant = get_ereal(inp, opt, 0, wi);
381 printStringNoNewline(inp, "Estimated diffusion constant (nm^2/ps, rad^2/ps or ps^-1)");
383 opt = prefix + "-diffusion";
384 dimParams->diffusion = get_ereal(inp, opt, 0, wi);
388 printStringNoNewline(inp,
389 "Diameter that needs to be sampled around a point before it is "
390 "considered covered. In FEP dimensions the cover diameter is "
391 "specified in lambda states.");
393 opt = prefix + "-cover-diameter";
394 dimParams->coverDiameter = get_ereal(inp, opt, 0, wi);
398 * Check the parameters of an AWH bias dimension.
400 * \param[in] prefix Prefix for dimension parameters.
401 * \param[in,out] dimParams AWH dimensional parameters.
402 * \param[in] ir Input parameter struct.
403 * \param[in,out] wi Struct for bookeeping warnings.
405 void checkDimParams(const std::string& prefix, AwhDimParams* dimParams, const t_inputrec* ir, warninp_t wi)
407 if (dimParams->eCoordProvider == eawhcoordproviderPULL)
412 "AWH biasing along a pull dimension is only compatible with COM pulling "
415 checkPullDimParams(prefix, dimParams, *ir->pull, wi);
417 else if (dimParams->eCoordProvider == eawhcoordproviderFREE_ENERGY_LAMBDA)
419 if (ir->efep == efepNO)
422 "AWH biasing along a free energy lambda state dimension is only compatible "
423 "with free energy turned on");
425 checkFepLambdaDimParams(prefix, dimParams, ir->fepvals, ir->efep, wi);
430 "AWH biasing can only be applied to pull and free energy lambda state "
436 * Check consistency of input at the AWH bias level.
438 * \param[in] awhBiasParams AWH bias parameters.
439 * \param[in,out] wi Struct for bookkeeping warnings.
441 void checkInputConsistencyAwhBias(const AwhBiasParams& awhBiasParams, warninp_t wi)
443 /* Covering diameter and sharing warning. */
444 for (int d = 0; d < awhBiasParams.ndim; d++)
446 double coverDiameter = awhBiasParams.dimParams[d].coverDiameter;
447 if (awhBiasParams.shareGroup <= 0 && coverDiameter > 0)
450 "The covering diameter is only relevant to set for bias sharing simulations.");
456 * Read parameters of an AWH bias.
458 * \param[in,out] inp Input file entries.
459 * \param[in,out] awhBiasParams AWH dimensional parameters.
460 * \param[in] prefix Prefix for bias parameters.
461 * \param[in,out] wi Struct for bookeeping warnings.
462 * \param[in] bComment True if comments should be printed.
464 void readBiasParams(std::vector<t_inpfile>* inp,
465 AwhBiasParams* awhBiasParams,
466 const std::string& prefix,
472 printStringNoNewline(inp, "Estimated initial PMF error (kJ/mol)");
475 std::string opt = prefix + "-error-init";
476 awhBiasParams->errorInitial = get_ereal(inp, opt, 10, wi);
480 printStringNoNewline(inp,
481 "Growth rate of the reference histogram determining the bias update "
482 "size: exp-linear or linear");
484 opt = prefix + "-growth";
485 awhBiasParams->eGrowth = get_eeenum(inp, opt, eawhgrowth_names, wi);
489 printStringNoNewline(inp,
490 "Start the simulation by equilibrating histogram towards the target "
491 "distribution: no or yes");
493 opt = prefix + "-equilibrate-histogram";
494 awhBiasParams->equilibrateHistogram = (get_eeenum(inp, opt, yesno_names, wi) != 0);
498 printStringNoNewline(
499 inp, "Target distribution type: constant, cutoff, boltzmann or local-boltzmann");
501 opt = prefix + "-target";
502 awhBiasParams->eTarget = get_eeenum(inp, opt, eawhtarget_names, wi);
506 printStringNoNewline(inp,
507 "Boltzmann beta scaling factor for target distribution types "
508 "'boltzmann' and 'boltzmann-local'");
510 opt = prefix + "-target-beta-scaling";
511 awhBiasParams->targetBetaScaling = get_ereal(inp, opt, 0, wi);
515 printStringNoNewline(inp, "Free energy cutoff value for target distribution type 'cutoff'");
517 opt = prefix + "-target-cutoff";
518 awhBiasParams->targetCutoff = get_ereal(inp, opt, 0, wi);
522 printStringNoNewline(inp, "Initialize PMF and target with user data: no or yes");
524 opt = prefix + "-user-data";
525 awhBiasParams->bUserData = get_eeenum(inp, opt, yesno_names, wi);
529 printStringNoNewline(inp, "Group index to share the bias with, 0 means not shared");
531 opt = prefix + "-share-group";
532 awhBiasParams->shareGroup = get_eint(inp, opt, 0, wi);
536 printStringNoNewline(inp, "Dimensionality of the coordinate");
538 opt = prefix + "-ndim";
539 awhBiasParams->ndim = get_eint(inp, opt, 0, wi);
541 /* Check this before starting to read the AWH dimension parameters. */
542 if (awhBiasParams->ndim <= 0 || awhBiasParams->ndim > c_biasMaxNumDim)
544 gmx_fatal(FARGS, "%s (%d) needs to be > 0 and at most %d\n", opt.c_str(),
545 awhBiasParams->ndim, c_biasMaxNumDim);
547 snew(awhBiasParams->dimParams, awhBiasParams->ndim);
548 for (int d = 0; d < awhBiasParams->ndim; d++)
550 bComment = bComment && d == 0;
551 std::string prefixdim = prefix + formatString("-dim%d", d + 1);
552 readDimParams(inp, prefixdim, &awhBiasParams->dimParams[d], wi, bComment);
557 * Check the parameters of an AWH bias.
559 * \param[in] awhBiasParams AWH dimensional parameters.
560 * \param[in] prefix Prefix for bias parameters.
561 * \param[in] ir Input parameter struct.
562 * \param[in,out] wi Struct for bookeeping warnings.
564 void checkBiasParams(const AwhBiasParams* awhBiasParams, const std::string& prefix, const t_inputrec* ir, warninp_t wi)
566 std::string opt = prefix + "-error-init";
567 if (awhBiasParams->errorInitial <= 0)
569 gmx_fatal(FARGS, "%s needs to be > 0.", opt.c_str());
572 opt = prefix + "-equilibrate-histogram";
573 if (awhBiasParams->equilibrateHistogram && awhBiasParams->eGrowth != eawhgrowthEXP_LINEAR)
576 formatString("Option %s will only have an effect for histogram growth type '%s'.",
577 opt.c_str(), EAWHGROWTH(eawhgrowthEXP_LINEAR));
578 warning(wi, message);
581 if ((awhBiasParams->eTarget == eawhtargetLOCALBOLTZMANN)
582 && (awhBiasParams->eGrowth == eawhgrowthEXP_LINEAR))
584 auto message = formatString(
585 "Target type '%s' combined with histogram growth type '%s' is not "
586 "expected to give stable bias updates. You probably want to use growth type "
588 EAWHTARGET(eawhtargetLOCALBOLTZMANN), EAWHGROWTH(eawhgrowthEXP_LINEAR),
589 EAWHGROWTH(eawhgrowthLINEAR));
590 warning(wi, message);
593 opt = prefix + "-target-beta-scaling";
594 switch (awhBiasParams->eTarget)
596 case eawhtargetBOLTZMANN:
597 case eawhtargetLOCALBOLTZMANN:
598 if (awhBiasParams->targetBetaScaling < 0 || awhBiasParams->targetBetaScaling > 1)
600 gmx_fatal(FARGS, "%s = %g is not useful for target type %s.", opt.c_str(),
601 awhBiasParams->targetBetaScaling, EAWHTARGET(awhBiasParams->eTarget));
605 if (awhBiasParams->targetBetaScaling != 0)
609 "Value for %s (%g) set explicitly but will not be used for target type %s.",
610 opt.c_str(), awhBiasParams->targetBetaScaling,
611 EAWHTARGET(awhBiasParams->eTarget));
616 opt = prefix + "-target-cutoff";
617 switch (awhBiasParams->eTarget)
619 case eawhtargetCUTOFF:
620 if (awhBiasParams->targetCutoff <= 0)
622 gmx_fatal(FARGS, "%s = %g is not useful for target type %s.", opt.c_str(),
623 awhBiasParams->targetCutoff, EAWHTARGET(awhBiasParams->eTarget));
627 if (awhBiasParams->targetCutoff != 0)
631 "Value for %s (%g) set explicitly but will not be used for target type %s.",
632 opt.c_str(), awhBiasParams->targetCutoff, EAWHTARGET(awhBiasParams->eTarget));
637 opt = prefix + "-share-group";
638 if (awhBiasParams->shareGroup < 0)
640 warning_error(wi, "AWH bias share-group should be >= 0");
643 opt = prefix + "-ndim";
644 if (awhBiasParams->ndim <= 0 || awhBiasParams->ndim > c_biasMaxNumDim)
646 gmx_fatal(FARGS, "%s (%d) needs to be > 0 and at most %d\n", opt.c_str(),
647 awhBiasParams->ndim, c_biasMaxNumDim);
649 if (awhBiasParams->ndim > 2)
652 "For awh-dim > 2 the estimate based on the diffusion and the initial error is "
653 "currently only a rough guideline."
654 " You should verify its usefulness for your system before production runs!");
656 for (int d = 0; d < awhBiasParams->ndim; d++)
658 std::string prefixdim = prefix + formatString("-dim%d", d + 1);
659 checkDimParams(prefixdim, &awhBiasParams->dimParams[d], ir, wi);
662 /* Check consistencies here that cannot be checked at read time at a lower level. */
663 checkInputConsistencyAwhBias(*awhBiasParams, wi);
667 * Check consistency of input at the AWH level.
669 * \param[in] awhParams AWH parameters.
670 * \param[in,out] wi Struct for bookkeeping warnings.
672 void checkInputConsistencyAwh(const AwhParams& awhParams, warninp_t wi)
674 /* Each pull coord can map to at most 1 AWH coord.
675 * Check that we have a shared bias when requesting multisim sharing.
677 bool haveSharedBias = false;
678 for (int k1 = 0; k1 < awhParams.numBias; k1++)
680 const AwhBiasParams& awhBiasParams1 = awhParams.awhBiasParams[k1];
682 if (awhBiasParams1.shareGroup > 0)
684 haveSharedBias = true;
687 /* k1 is the reference AWH, k2 is the AWH we compare with (can be equal to k1) */
688 for (int k2 = k1; k2 < awhParams.numBias; k2++)
690 for (int d1 = 0; d1 < awhBiasParams1.ndim; d1++)
692 if (awhBiasParams1.dimParams[d1].eCoordProvider == eawhcoordproviderFREE_ENERGY_LAMBDA)
696 const AwhBiasParams& awhBiasParams2 = awhParams.awhBiasParams[k2];
698 /* d1 is the reference dimension of the reference AWH. d2 is the dim index of the AWH to compare with. */
699 for (int d2 = 0; d2 < awhBiasParams2.ndim; d2++)
701 if (awhBiasParams2.dimParams[d2].eCoordProvider == eawhcoordproviderFREE_ENERGY_LAMBDA)
705 /* Give an error if (d1, k1) is different from (d2, k2) but the pull coordinate is the same */
706 if ((d1 != d2 || k1 != k2)
707 && (awhBiasParams1.dimParams[d1].coordIndex == awhBiasParams2.dimParams[d2].coordIndex))
709 char errormsg[STRLEN];
711 "One pull coordinate (%d) cannot be mapped to two separate AWH "
712 "dimensions (awh%d-dim%d and awh%d-dim%d). "
713 "If this is really what you want to do you will have to duplicate "
714 "this pull coordinate.",
715 awhBiasParams1.dimParams[d1].coordIndex + 1, k1 + 1, d1 + 1, k2 + 1,
717 gmx_fatal(FARGS, "%s", errormsg);
724 if (awhParams.shareBiasMultisim && !haveSharedBias)
727 "Sharing of biases over multiple simulations is requested, but no bias is marked "
728 "as shared (share-group > 0)");
731 /* mdrun does not support this (yet), but will check again */
732 if (haveBiasSharingWithinSimulation(awhParams))
735 "You have shared biases within a single simulation, but mdrun does not support "
741 AwhParams* readAwhParams(std::vector<t_inpfile>* inp, warninp_t wi)
743 AwhParams* awhParams;
747 /* Parameters common for all biases */
749 printStringNoNewline(inp, "The way to apply the biasing potential: convolved or umbrella");
750 opt = "awh-potential";
751 awhParams->ePotential = get_eeenum(inp, opt, eawhpotential_names, wi);
753 printStringNoNewline(inp,
754 "The random seed used for sampling the umbrella center in the case of "
755 "umbrella type potential");
757 awhParams->seed = get_eint(inp, opt, -1, wi);
758 if (awhParams->seed == -1)
760 awhParams->seed = static_cast<int>(gmx::makeRandomSeed());
761 fprintf(stderr, "Setting the AWH bias MC random seed to %" PRId64 "\n", awhParams->seed);
764 printStringNoNewline(inp, "Data output interval in number of steps");
766 awhParams->nstOut = get_eint(inp, opt, 100000, wi);
768 printStringNoNewline(inp, "Coordinate sampling interval in number of steps");
769 opt = "awh-nstsample";
770 awhParams->nstSampleCoord = get_eint(inp, opt, 10, wi);
772 printStringNoNewline(inp, "Free energy and bias update interval in number of samples");
773 opt = "awh-nsamples-update";
774 awhParams->numSamplesUpdateFreeEnergy = get_eint(inp, opt, 10, wi);
776 printStringNoNewline(
777 inp, "When true, biases with share-group>0 are shared between multiple simulations");
778 opt = "awh-share-multisim";
779 awhParams->shareBiasMultisim = (get_eeenum(inp, opt, yesno_names, wi) != 0);
781 printStringNoNewline(inp, "The number of independent AWH biases");
783 awhParams->numBias = get_eint(inp, opt, 1, wi);
784 /* Check this before starting to read the AWH biases. */
785 if (awhParams->numBias <= 0)
787 gmx_fatal(FARGS, "%s needs to be an integer > 0", opt.c_str());
790 /* Read the parameters specific to each AWH bias */
791 snew(awhParams->awhBiasParams, awhParams->numBias);
793 for (int k = 0; k < awhParams->numBias; k++)
795 bool bComment = (k == 0);
796 std::string prefixawh = formatString("awh%d", k + 1);
797 readBiasParams(inp, &awhParams->awhBiasParams[k], prefixawh, wi, bComment);
803 void checkAwhParams(const AwhParams* awhParams, const t_inputrec* ir, warninp_t wi)
807 checkMtsConsistency(*ir, wi);
810 if (awhParams->nstOut <= 0)
812 auto message = formatString("Not writing AWH output with AWH (%s = %d) does not make sense",
813 opt.c_str(), awhParams->nstOut);
814 warning_error(wi, message);
816 /* This restriction can be removed by changing a flag of print_ebin() */
817 if (ir->nstenergy == 0 || awhParams->nstOut % ir->nstenergy != 0)
819 auto message = formatString("%s (%d) should be a multiple of nstenergy (%d)", opt.c_str(),
820 awhParams->nstOut, ir->nstenergy);
821 warning_error(wi, message);
824 opt = "awh-nsamples-update";
825 if (awhParams->numSamplesUpdateFreeEnergy <= 0)
827 warning_error(wi, opt + " needs to be an integer > 0");
830 for (int k = 0; k < awhParams->numBias; k++)
832 std::string prefixawh = formatString("awh%d", k + 1);
833 checkBiasParams(&awhParams->awhBiasParams[k], prefixawh, ir, wi);
836 /* Do a final consistency check before returning */
837 checkInputConsistencyAwh(*awhParams, wi);
839 if (ir->init_step != 0)
841 warning_error(wi, "With AWH init-step should be 0");
846 * Gets the period of a pull coordinate.
848 * \param[in] pullCoordParams The parameters for the pull coordinate.
849 * \param[in] pbc The PBC setup
850 * \param[in] intervalLength The length of the AWH interval for this pull coordinate
851 * \returns the period (or 0 if not periodic).
853 static double get_pull_coord_period(const t_pull_coord& pullCoordParams, const t_pbc& pbc, const real intervalLength)
857 if (pullCoordParams.eGeom == epullgDIR)
859 const real margin = 0.001;
860 // Make dims periodic when the interval covers > 95%
861 const real periodicFraction = 0.95;
863 // Check if the pull direction is along a box vector
864 for (int dim = 0; dim < pbc.ndim_ePBC; dim++)
866 const real boxLength = norm(pbc.box[dim]);
867 const real innerProduct = iprod(pullCoordParams.vec, pbc.box[dim]);
868 if (innerProduct >= (1 - margin) * boxLength && innerProduct <= (1 + margin) * boxLength)
870 if (intervalLength > (1 + margin) * boxLength)
873 "The AWH interval (%f nm) for a pull coordinate is larger than the "
875 intervalLength, boxLength);
878 if (intervalLength > periodicFraction * boxLength)
885 else if (pullCoordParams.eGeom == epullgDIHEDRAL)
887 /* The dihedral angle is periodic in -180 to 180 deg */
895 * Checks if the given interval is defined in the correct periodic interval.
897 * \param[in] origin Start value of interval.
898 * \param[in] end End value of interval.
899 * \param[in] period Period (or 0 if not periodic).
900 * \returns true if the end point values are in the correct periodic interval.
902 static bool intervalIsInPeriodicInterval(double origin, double end, double period)
904 return (period == 0) || (std::fabs(origin) <= 0.5 * period && std::fabs(end) <= 0.5 * period);
908 * Checks if a value is within an interval.
910 * \param[in] origin Start value of interval.
911 * \param[in] end End value of interval.
912 * \param[in] period Period (or 0 if not periodic).
913 * \param[in] value Value to check.
914 * \returns true if the value is within the interval.
916 static bool valueIsInInterval(double origin, double end, double period, double value)
924 /* The interval closes within the periodic interval */
925 bIn_interval = (value >= origin) && (value <= end);
929 /* The interval wraps around the periodic boundary */
930 bIn_interval = ((value >= origin) && (value <= 0.5 * period))
931 || ((value >= -0.5 * period) && (value <= end));
936 bIn_interval = (value >= origin) && (value <= end);
943 * Check if the starting configuration is consistent with the given interval.
945 * \param[in] awhParams AWH parameters.
946 * \param[in,out] wi Struct for bookeeping warnings.
948 static void checkInputConsistencyInterval(const AwhParams* awhParams, warninp_t wi)
950 for (int k = 0; k < awhParams->numBias; k++)
952 AwhBiasParams* awhBiasParams = &awhParams->awhBiasParams[k];
953 for (int d = 0; d < awhBiasParams->ndim; d++)
955 AwhDimParams* dimParams = &awhBiasParams->dimParams[d];
956 int coordIndex = dimParams->coordIndex;
957 double origin = dimParams->origin, end = dimParams->end, period = dimParams->period;
958 double coordValueInit = dimParams->coordValueInit;
960 if ((period == 0) && (origin > end))
963 "For the non-periodic pull coordinates awh%d-dim%d-start (%f) cannot be "
964 "larger than awh%d-dim%d-end (%f)",
965 k + 1, d + 1, origin, k + 1, d + 1, end);
968 /* Currently we assume symmetric periodic intervals, meaning we use [-period/2, period/2] as the reference interval.
969 Make sure the AWH interval is within this reference interval.
971 Note: we could fairly simply allow using a more general interval (e.g. [x, x + period]) but it complicates
972 things slightly and I don't see that there is a great need for it. It would also mean that the interval would
973 depend on AWH input. Also, for dihedral angles you would always want the reference interval to be -180, +180,
974 independent of AWH parameters.
976 if (!intervalIsInPeriodicInterval(origin, end, period))
979 "When using AWH with periodic pull coordinate geometries "
980 "awh%d-dim%d-start (%.8g) and "
981 "awh%d-dim%d-end (%.8g) should cover at most one period (%.8g) and take "
983 "minus half a period and plus half a period, i.e. in the interval [%.8g, "
985 k + 1, d + 1, origin, k + 1, d + 1, end, period, -0.5 * period, 0.5 * period);
988 /* Warn if the pull initial coordinate value is not in the grid */
989 if (!valueIsInInterval(origin, end, period, coordValueInit))
991 auto message = formatString(
992 "The initial coordinate value (%.8g) for pull coordinate index %d falls "
994 "of the sampling nterval awh%d-dim%d-start (%.8g) to awh%d-dim%d-end "
996 "This can lead to large initial forces pulling the coordinate towards the "
997 "sampling interval.",
998 coordValueInit, coordIndex + 1, k + 1, d + 1, origin, k + 1, d + 1, end);
999 warning(wi, message);
1006 * Sets AWH parameters, for one AWH pull dimension.
1008 * \param[in,out] dimParams AWH dimension parameters.
1009 * \param[in] biasIndex The index of the bias containing this AWH pull dimension.
1010 * \param[in] dimIndex The index of this AWH pull dimension.
1011 * \param[in] pull_params Pull parameters.
1012 * \param[in,out] pull_work Pull working struct to register AWH bias in.
1013 * \param[in] pbc A pbc information structure.
1014 * \param[in] compressibility Compressibility matrix for pressure coupling,
1015 * pass all 0 without pressure coupling.
1016 * \param[in,out] wi Struct for bookeeping warnings.
1019 static void setStateDependentAwhPullDimParams(AwhDimParams* dimParams,
1020 const int biasIndex,
1022 const pull_params_t* pull_params,
1025 const tensor& compressibility,
1028 const t_pull_coord& pullCoordParams = pull_params->coord[dimParams->coordIndex];
1030 if (pullCoordParams.eGeom == epullgDIRPBC)
1033 "AWH does not support pull geometry '%s'. "
1034 "If the maximum distance between the groups is always "
1035 "less than half the box size, "
1036 "you can use geometry '%s' instead.",
1037 EPULLGEOM(epullgDIRPBC), EPULLGEOM(epullgDIR));
1040 dimParams->period = get_pull_coord_period(pullCoordParams, pbc, dimParams->end - dimParams->origin);
1041 // We would like to check for scaling, but we don't have the full inputrec available here
1042 if (dimParams->period > 0
1043 && !(pullCoordParams.eGeom == epullgANGLE || pullCoordParams.eGeom == epullgDIHEDRAL))
1045 bool coordIsScaled = false;
1046 for (int d2 = 0; d2 < DIM; d2++)
1048 if (pullCoordParams.vec[d2] != 0 && norm2(compressibility[d2]) != 0)
1050 coordIsScaled = true;
1055 std::string mesg = gmx::formatString(
1056 "AWH dimension %d of bias %d is periodic with pull geometry '%s', "
1057 "while you should are applying pressure scaling to the "
1058 "corresponding box vector, this is not supported.",
1059 dimIndex + 1, biasIndex + 1, EPULLGEOM(pullCoordParams.eGeom));
1060 warning(wi, mesg.c_str());
1064 /* The initial coordinate value, converted to external user units. */
1065 dimParams->coordValueInit = get_pull_coord_value(pull_work, dimParams->coordIndex, &pbc);
1067 dimParams->coordValueInit *= pull_conversion_factor_internal2userinput(&pullCoordParams);
1070 void setStateDependentAwhParams(AwhParams* awhParams,
1071 const pull_params_t& pull_params,
1075 const tensor& compressibility,
1076 const t_grpopts* inputrecGroupOptions,
1077 const real initLambda,
1078 const gmx_mtop_t& mtop,
1081 /* The temperature is not really state depenendent but is not known
1082 * when read_awhParams is called (in get ir).
1083 * It is known first after do_index has been called in grompp.cpp.
1085 if (inputrecGroupOptions->ref_t == nullptr || inputrecGroupOptions->ref_t[0] <= 0)
1087 gmx_fatal(FARGS, "AWH biasing is only supported for temperatures > 0");
1089 for (int i = 1; i < inputrecGroupOptions->ngtc; i++)
1091 if (inputrecGroupOptions->ref_t[i] != inputrecGroupOptions->ref_t[0])
1094 "AWH biasing is currently only supported for identical temperatures for all "
1095 "temperature coupling groups");
1100 set_pbc(&pbc, pbcType, box);
1102 for (int k = 0; k < awhParams->numBias; k++)
1104 AwhBiasParams* awhBiasParams = &awhParams->awhBiasParams[k];
1105 for (int d = 0; d < awhBiasParams->ndim; d++)
1107 AwhDimParams* dimParams = &awhBiasParams->dimParams[d];
1108 if (dimParams->eCoordProvider == eawhcoordproviderPULL)
1110 setStateDependentAwhPullDimParams(dimParams, k, d, &pull_params, pull_work, pbc,
1111 compressibility, wi);
1115 dimParams->coordValueInit = initLambda;
1116 checkFepLambdaDimDecouplingConsistency(mtop, wi);
1120 checkInputConsistencyInterval(awhParams, wi);
1122 /* Register AWH as external potential with pull (for AWH dimensions that use pull as
1123 * reaction coordinate provider) to check consistency. */
1124 Awh::registerAwhWithPull(*awhParams, pull_work);