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,2021, 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/enumerationhelpers.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/smalloc.h"
61 #include "gromacs/utility/stringutil.h"
63 #include "biasparams.h"
64 #include "biassharing.h"
69 const char* enumValueToString(AwhTargetType enumValue)
71 constexpr gmx::EnumerationArray<AwhTargetType, const char*> awhTargetTypeNames = {
72 "constant", "cutoff", "boltzmann", "local-boltzmann"
74 return awhTargetTypeNames[enumValue];
77 const char* enumValueToString(AwhHistogramGrowthType enumValue)
79 constexpr gmx::EnumerationArray<AwhHistogramGrowthType, const char*> awhHistogramGrowthTypeNames = {
80 "exp-linear", "linear"
82 return awhHistogramGrowthTypeNames[enumValue];
85 const char* enumValueToString(AwhPotentialType enumValue)
87 constexpr gmx::EnumerationArray<AwhPotentialType, const char*> awhPotentialTypeNames = {
88 "convolved", "umbrella"
90 return awhPotentialTypeNames[enumValue];
93 const char* enumValueToString(AwhCoordinateProviderType enumValue)
95 constexpr gmx::EnumerationArray<AwhCoordinateProviderType, const char*> awhCoordinateProviderTypeNames = {
98 return awhCoordinateProviderTypeNames[enumValue];
105 * Check multiple time-stepping consistency between AWH and pull and/or FEP
107 * \param[in] inputrec Inputput parameter struct.
108 * \param[in,out] wi Struct for bookeeping warnings.
110 void checkMtsConsistency(const t_inputrec& inputrec, warninp_t wi)
112 if (!inputrec.useMts)
117 GMX_RELEASE_ASSERT(inputrec.mtsLevels.size() == 2, "Only 2 MTS levels supported here");
119 bool usesPull = false;
120 bool usesFep = false;
121 for (int b = 0; b < inputrec.awhParams->numBias; b++)
123 const auto& biasParams = inputrec.awhParams->awhBiasParams[b];
124 for (int d = 0; d < biasParams.ndim; d++)
126 switch (biasParams.dimParams[d].eCoordProvider)
128 case AwhCoordinateProviderType::Pull: usesPull = true; break;
129 case AwhCoordinateProviderType::FreeEnergyLambda: usesFep = true; break;
130 default: GMX_RELEASE_ASSERT(false, "Unsupported coord provider");
134 const int awhMtsLevel = forceGroupMtsLevel(inputrec.mtsLevels, MtsForceGroups::Awh);
135 if (usesPull && forceGroupMtsLevel(inputrec.mtsLevels, MtsForceGroups::Pull) != awhMtsLevel)
138 "When AWH is applied to pull coordinates, pull and AWH should be computed at "
139 "the same MTS level");
141 if (usesFep && awhMtsLevel != ssize(inputrec.mtsLevels) - 1)
144 "When AWH is applied to the free-energy lambda with MTS, AWH should be "
145 "computed at the slow MTS level");
148 if (inputrec.awhParams->nstSampleCoord % inputrec.mtsLevels[awhMtsLevel].stepFactor != 0)
151 "With MTS applied to AWH, awh-nstsample should be a multiple of mts-factor");
156 * Check the parameters of an AWH bias pull dimension.
158 * \param[in] prefix Prefix for dimension parameters.
159 * \param[in,out] dimParams AWH dimensional parameters.
160 * \param[in] pull_params Pull parameters.
161 * \param[in,out] wi Struct for bookeeping warnings.
163 void checkPullDimParams(const std::string& prefix,
164 AwhDimParams* dimParams,
165 const pull_params_t& pull_params,
168 if (dimParams->coordIndex < 0)
171 "Failed to read a valid coordinate index for %s-coord-index. "
172 "Note that the pull coordinate indexing starts at 1.",
175 if (dimParams->coordIndex >= pull_params.ncoord)
178 "The given AWH coordinate index (%d) is larger than the number of pull "
180 dimParams->coordIndex + 1,
183 if (pull_params.coord[dimParams->coordIndex].rate != 0)
185 auto message = formatString(
186 "Setting pull-coord%d-rate (%g) is incompatible with AWH biasing this coordinate",
187 dimParams->coordIndex + 1,
188 pull_params.coord[dimParams->coordIndex].rate);
189 warning_error(wi, message);
192 if (gmx_within_tol(dimParams->end - dimParams->origin, 0, GMX_REAL_EPS))
194 auto message = formatString(
195 "The given interval length given by %s-start (%g) and %s-end (%g) is zero. "
196 "This will result in only one point along this axis in the coordinate value grid.",
201 warning(wi, message);
204 if (dimParams->forceConstant <= 0)
206 warning_error(wi, "The force AWH bias force constant should be > 0");
209 /* Grid params for each axis */
210 PullGroupGeometry eGeom = pull_params.coord[dimParams->coordIndex].eGeom;
212 /* Check that the requested interval is in allowed range */
213 if (eGeom == PullGroupGeometry::Distance)
215 if (dimParams->origin < 0 || dimParams->end < 0)
218 "%s-start (%g) or %s-end (%g) set to a negative value. With pull "
219 "geometry distance coordinate values are non-negative. "
220 "Perhaps you want to use geometry %s instead?",
225 enumValueToString(PullGroupGeometry::Direction));
228 else if (eGeom == PullGroupGeometry::Angle || eGeom == PullGroupGeometry::AngleAxis)
230 if (dimParams->origin < 0 || dimParams->end > 180)
233 "%s-start (%g) and %s-end (%g) are outside of the allowed range "
234 "0 to 180 deg for pull geometries %s and %s ",
239 enumValueToString(PullGroupGeometry::Angle),
240 enumValueToString(PullGroupGeometry::AngleAxis));
243 else if (eGeom == PullGroupGeometry::Dihedral)
245 if (dimParams->origin < -180 || dimParams->end > 180)
248 "%s-start (%g) and %s-end (%g) are outside of the allowed range "
249 "-180 to 180 deg for pull geometry %s. ",
254 enumValueToString(PullGroupGeometry::Dihedral));
260 * Check parameters of an AWH bias in a free energy lambda state dimension.
262 * \param[in] prefix Prefix for dimension parameters.
263 * \param[in,out] dimParams AWH dimensional parameters.
264 * \param[in] lambdaParams The free energy lambda related parameters.
265 * \param[in] efep This is the type of FEP calculation (efep enumerator).
266 * \param[in,out] wi Struct for bookeeping warnings.
268 void checkFepLambdaDimParams(const std::string& prefix,
269 const AwhDimParams* dimParams,
270 const t_lambda* lambdaParams,
271 const FreeEnergyPerturbationType efep,
279 "There must be free energy input if using AWH to steer the free energy lambda "
283 if (lambdaParams->lambda_neighbors != -1)
286 "When running AWH coupled to the free energy lambda state all lambda states "
287 "should be used as neighbors in order to get correct probabilities, i.e. "
288 "calc-lambda-neighbors (%d) must be %d.",
289 lambdaParams->lambda_neighbors,
293 if (efep == FreeEnergyPerturbationType::SlowGrowth || lambdaParams->delta_lambda != 0)
296 "AWH coupled to the free energy lambda state is not compatible with slow-growth "
297 "and delta-lambda must be 0.");
300 if (efep == FreeEnergyPerturbationType::Expanded)
303 "AWH is not treated like other expanded ensemble methods. Do not use expanded.");
306 if (dimParams->origin < 0)
308 opt = prefix + "-start";
310 "When running AWH coupled to the free energy lambda state the lower lambda state "
311 "for AWH, %s (%.0f), must be >= 0.",
315 if (dimParams->end >= lambdaParams->n_lambda)
317 opt = prefix + "-end";
319 "When running AWH coupled to the free energy lambda state the upper lambda state "
320 "for AWH, %s (%.0f), must be < n_lambda (%d).",
323 lambdaParams->n_lambda);
325 if (gmx_within_tol(dimParams->end - dimParams->origin, 0, GMX_REAL_EPS))
327 auto message = formatString(
328 "The given interval length given by %s-start (%g) and %s-end (%g) is zero. "
329 "This will result in only one lambda point along this free energy lambda state "
330 "axis in the coordinate value grid.",
335 warning(wi, message);
338 if (dimParams->forceConstant != 0)
342 "The force AWH bias force constant is not used with free energy lambda state as "
343 "coordinate provider.");
348 * Check that AWH FEP is not combined with incompatible decoupling.
350 * \param[in] mtop System topology.
351 * \param[in,out] wi Struct for bookeeping warnings.
353 void checkFepLambdaDimDecouplingConsistency(const gmx_mtop_t& mtop, warninp_t wi)
355 if (haveFepPerturbedMasses(mtop))
358 "Masses may not be perturbed when using the free energy lambda state as AWH "
359 "coordinate provider. If you are using fep-lambdas to specify lambda states "
360 "make sure that you also specify mass-lambdas without perturbation.");
362 if (havePerturbedConstraints(mtop))
365 "Constraints may not be perturbed when using the free energy lambda state as AWH "
366 "coordinate provider. If you are using fep-lambdas to specify lambda states "
367 "make sure that you also specify mass-lambdas without perturbation.");
372 * Read parameters of an AWH bias dimension.
374 * \param[in,out] inp Input file entries.
375 * \param[in] prefix Prefix for dimension parameters.
376 * \param[in,out] dimParams AWH dimensional parameters.
377 * \param[in,out] wi Struct for bookeeping warnings.
378 * \param[in] bComment True if comments should be printed.
380 void readDimParams(std::vector<t_inpfile>* inp,
381 const std::string& prefix,
382 AwhDimParams* dimParams,
389 printStringNoNewline(
391 "The provider of the reaction coordinate, "
392 "currently only 'pull' and 'fep-lambda' (free energy lambda state) is supported");
394 opt = prefix + "-coord-provider";
395 dimParams->eCoordProvider = getEnum<AwhCoordinateProviderType>(inp, opt.c_str(), wi);
399 printStringNoNewline(inp, "The coordinate index for this dimension");
401 opt = prefix + "-coord-index";
403 coordIndexInput = get_eint(inp, opt, 1, wi);
405 /* The pull coordinate indices start at 1 in the input file, at 0 internally */
406 dimParams->coordIndex = coordIndexInput - 1;
410 printStringNoNewline(inp, "Start and end values for each coordinate dimension");
412 opt = prefix + "-start";
413 dimParams->origin = get_ereal(inp, opt, 0., wi);
414 opt = prefix + "-end";
415 dimParams->end = get_ereal(inp, opt, 0., wi);
419 printStringNoNewline(
420 inp, "The force constant for this coordinate (kJ/mol/nm^2 or kJ/mol/rad^2)");
422 opt = prefix + "-force-constant";
423 dimParams->forceConstant = get_ereal(inp, opt, 0, wi);
427 printStringNoNewline(inp, "Estimated diffusion constant (nm^2/ps, rad^2/ps or ps^-1)");
429 opt = prefix + "-diffusion";
430 dimParams->diffusion = get_ereal(inp, opt, 0, wi);
434 printStringNoNewline(inp,
435 "Diameter that needs to be sampled around a point before it is "
436 "considered covered. In FEP dimensions the cover diameter is "
437 "specified in lambda states.");
439 opt = prefix + "-cover-diameter";
440 dimParams->coverDiameter = get_ereal(inp, opt, 0, wi);
444 * Check the parameters of an AWH bias dimension.
446 * \param[in] prefix Prefix for dimension parameters.
447 * \param[in,out] dimParams AWH dimensional parameters.
448 * \param[in] ir Input parameter struct.
449 * \param[in,out] wi Struct for bookeeping warnings.
451 void checkDimParams(const std::string& prefix, AwhDimParams* dimParams, const t_inputrec* ir, warninp_t wi)
453 if (dimParams->eCoordProvider == AwhCoordinateProviderType::Pull)
458 "AWH biasing along a pull dimension is only compatible with COM pulling "
461 checkPullDimParams(prefix, dimParams, *ir->pull, wi);
463 else if (dimParams->eCoordProvider == AwhCoordinateProviderType::FreeEnergyLambda)
465 if (ir->efep == FreeEnergyPerturbationType::No)
468 "AWH biasing along a free energy lambda state dimension is only compatible "
469 "with free energy turned on");
471 checkFepLambdaDimParams(prefix, dimParams, ir->fepvals.get(), ir->efep, wi);
476 "AWH biasing can only be applied to pull and free energy lambda state "
482 * Check consistency of input at the AWH bias level.
484 * \param[in] awhBiasParams AWH bias parameters.
485 * \param[in,out] wi Struct for bookkeeping warnings.
487 void checkInputConsistencyAwhBias(const AwhBiasParams& awhBiasParams, warninp_t wi)
489 /* Covering diameter and sharing warning. */
490 for (int d = 0; d < awhBiasParams.ndim; d++)
492 double coverDiameter = awhBiasParams.dimParams[d].coverDiameter;
493 if (awhBiasParams.shareGroup <= 0 && coverDiameter > 0)
496 "The covering diameter is only relevant to set for bias sharing simulations.");
502 * Read parameters of an AWH bias.
504 * \param[in,out] inp Input file entries.
505 * \param[in,out] awhBiasParams AWH dimensional parameters.
506 * \param[in] prefix Prefix for bias parameters.
507 * \param[in,out] wi Struct for bookeeping warnings.
508 * \param[in] bComment True if comments should be printed.
510 void readBiasParams(std::vector<t_inpfile>* inp,
511 AwhBiasParams* awhBiasParams,
512 const std::string& prefix,
518 printStringNoNewline(inp, "Estimated initial PMF error (kJ/mol)");
521 std::string opt = prefix + "-error-init";
522 awhBiasParams->errorInitial = get_ereal(inp, opt, 10, wi);
526 printStringNoNewline(inp,
527 "Growth rate of the reference histogram determining the bias update "
528 "size: exp-linear or linear");
530 opt = prefix + "-growth";
531 awhBiasParams->eGrowth = getEnum<AwhHistogramGrowthType>(inp, opt.c_str(), wi);
535 printStringNoNewline(inp,
536 "Start the simulation by equilibrating histogram towards the target "
537 "distribution: no or yes");
539 opt = prefix + "-equilibrate-histogram";
540 awhBiasParams->equilibrateHistogram = (getEnum<Boolean>(inp, opt.c_str(), wi) != Boolean::No);
544 printStringNoNewline(
545 inp, "Target distribution type: constant, cutoff, boltzmann or local-boltzmann");
547 opt = prefix + "-target";
548 awhBiasParams->eTarget = getEnum<AwhTargetType>(inp, opt.c_str(), wi);
552 printStringNoNewline(inp,
553 "Boltzmann beta scaling factor for target distribution types "
554 "'boltzmann' and 'boltzmann-local'");
556 opt = prefix + "-target-beta-scaling";
557 awhBiasParams->targetBetaScaling = get_ereal(inp, opt, 0, wi);
561 printStringNoNewline(inp, "Free energy cutoff value for target distribution type 'cutoff'");
563 opt = prefix + "-target-cutoff";
564 awhBiasParams->targetCutoff = get_ereal(inp, opt, 0, wi);
568 printStringNoNewline(inp, "Initialize PMF and target with user data: no or yes");
570 opt = prefix + "-user-data";
571 awhBiasParams->bUserData = getEnum<Boolean>(inp, opt.c_str(), wi) != Boolean::No;
575 printStringNoNewline(inp, "Group index to share the bias with, 0 means not shared");
577 opt = prefix + "-share-group";
578 awhBiasParams->shareGroup = get_eint(inp, opt, 0, wi);
582 printStringNoNewline(inp, "Dimensionality of the coordinate");
584 opt = prefix + "-ndim";
585 awhBiasParams->ndim = get_eint(inp, opt, 0, wi);
587 /* Check this before starting to read the AWH dimension parameters. */
588 if (awhBiasParams->ndim <= 0 || awhBiasParams->ndim > c_biasMaxNumDim)
590 gmx_fatal(FARGS, "%s (%d) needs to be > 0 and at most %d\n", opt.c_str(), awhBiasParams->ndim, c_biasMaxNumDim);
592 snew(awhBiasParams->dimParams, awhBiasParams->ndim);
593 for (int d = 0; d < awhBiasParams->ndim; d++)
595 bComment = bComment && d == 0;
596 std::string prefixdim = prefix + formatString("-dim%d", d + 1);
597 readDimParams(inp, prefixdim, &awhBiasParams->dimParams[d], wi, bComment);
602 * Check the parameters of an AWH bias.
604 * \param[in] awhBiasParams AWH dimensional parameters.
605 * \param[in] prefix Prefix for bias parameters.
606 * \param[in] ir Input parameter struct.
607 * \param[in,out] wi Struct for bookeeping warnings.
609 void checkBiasParams(const AwhBiasParams* awhBiasParams, const std::string& prefix, const t_inputrec* ir, warninp_t wi)
611 std::string opt = prefix + "-error-init";
612 if (awhBiasParams->errorInitial <= 0)
614 gmx_fatal(FARGS, "%s needs to be > 0.", opt.c_str());
617 opt = prefix + "-equilibrate-histogram";
618 if (awhBiasParams->equilibrateHistogram && awhBiasParams->eGrowth != AwhHistogramGrowthType::ExponentialLinear)
621 formatString("Option %s will only have an effect for histogram growth type '%s'.",
623 enumValueToString(AwhHistogramGrowthType::ExponentialLinear));
624 warning(wi, message);
627 if ((awhBiasParams->eTarget == AwhTargetType::LocalBoltzmann)
628 && (awhBiasParams->eGrowth == AwhHistogramGrowthType::ExponentialLinear))
630 auto message = formatString(
631 "Target type '%s' combined with histogram growth type '%s' is not "
632 "expected to give stable bias updates. You probably want to use growth type "
634 enumValueToString(AwhTargetType::LocalBoltzmann),
635 enumValueToString(AwhHistogramGrowthType::ExponentialLinear),
636 enumValueToString(AwhHistogramGrowthType::Linear));
637 warning(wi, message);
640 opt = prefix + "-target-beta-scaling";
641 switch (awhBiasParams->eTarget)
643 case AwhTargetType::Boltzmann:
644 case AwhTargetType::LocalBoltzmann:
645 if (awhBiasParams->targetBetaScaling < 0 || awhBiasParams->targetBetaScaling > 1)
648 "%s = %g is not useful for target type %s.",
650 awhBiasParams->targetBetaScaling,
651 enumValueToString(awhBiasParams->eTarget));
655 if (awhBiasParams->targetBetaScaling != 0)
659 "Value for %s (%g) set explicitly but will not be used for target type %s.",
661 awhBiasParams->targetBetaScaling,
662 enumValueToString(awhBiasParams->eTarget));
667 opt = prefix + "-target-cutoff";
668 switch (awhBiasParams->eTarget)
670 case AwhTargetType::Cutoff:
671 if (awhBiasParams->targetCutoff <= 0)
674 "%s = %g is not useful for target type %s.",
676 awhBiasParams->targetCutoff,
677 enumValueToString(awhBiasParams->eTarget));
681 if (awhBiasParams->targetCutoff != 0)
685 "Value for %s (%g) set explicitly but will not be used for target type %s.",
687 awhBiasParams->targetCutoff,
688 enumValueToString(awhBiasParams->eTarget));
693 opt = prefix + "-share-group";
694 if (awhBiasParams->shareGroup < 0)
696 warning_error(wi, "AWH bias share-group should be >= 0");
699 opt = prefix + "-ndim";
700 if (awhBiasParams->ndim <= 0 || awhBiasParams->ndim > c_biasMaxNumDim)
702 gmx_fatal(FARGS, "%s (%d) needs to be > 0 and at most %d\n", opt.c_str(), awhBiasParams->ndim, c_biasMaxNumDim);
704 if (awhBiasParams->ndim > 2)
707 "For awh-dim > 2 the estimate based on the diffusion and the initial error is "
708 "currently only a rough guideline."
709 " You should verify its usefulness for your system before production runs!");
711 for (int d = 0; d < awhBiasParams->ndim; d++)
713 std::string prefixdim = prefix + formatString("-dim%d", d + 1);
714 checkDimParams(prefixdim, &awhBiasParams->dimParams[d], ir, wi);
717 /* Check consistencies here that cannot be checked at read time at a lower level. */
718 checkInputConsistencyAwhBias(*awhBiasParams, wi);
722 * Check consistency of input at the AWH level.
724 * \param[in] awhParams AWH parameters.
725 * \param[in,out] wi Struct for bookkeeping warnings.
727 void checkInputConsistencyAwh(const AwhParams& awhParams, warninp_t wi)
729 /* Each pull coord can map to at most 1 AWH coord.
730 * Check that we have a shared bias when requesting multisim sharing.
732 bool haveSharedBias = false;
733 for (int k1 = 0; k1 < awhParams.numBias; k1++)
735 const AwhBiasParams& awhBiasParams1 = awhParams.awhBiasParams[k1];
737 if (awhBiasParams1.shareGroup > 0)
739 haveSharedBias = true;
742 /* k1 is the reference AWH, k2 is the AWH we compare with (can be equal to k1) */
743 for (int k2 = k1; k2 < awhParams.numBias; k2++)
745 for (int d1 = 0; d1 < awhBiasParams1.ndim; d1++)
747 if (awhBiasParams1.dimParams[d1].eCoordProvider == AwhCoordinateProviderType::FreeEnergyLambda)
751 const AwhBiasParams& awhBiasParams2 = awhParams.awhBiasParams[k2];
753 /* d1 is the reference dimension of the reference AWH. d2 is the dim index of the AWH to compare with. */
754 for (int d2 = 0; d2 < awhBiasParams2.ndim; d2++)
756 if (awhBiasParams2.dimParams[d2].eCoordProvider == AwhCoordinateProviderType::FreeEnergyLambda)
760 /* Give an error if (d1, k1) is different from (d2, k2) but the pull coordinate is the same */
761 if ((d1 != d2 || k1 != k2)
762 && (awhBiasParams1.dimParams[d1].coordIndex == awhBiasParams2.dimParams[d2].coordIndex))
764 char errormsg[STRLEN];
766 "One pull coordinate (%d) cannot be mapped to two separate AWH "
767 "dimensions (awh%d-dim%d and awh%d-dim%d). "
768 "If this is really what you want to do you will have to duplicate "
769 "this pull coordinate.",
770 awhBiasParams1.dimParams[d1].coordIndex + 1,
775 gmx_fatal(FARGS, "%s", errormsg);
782 if (awhParams.shareBiasMultisim && !haveSharedBias)
785 "Sharing of biases over multiple simulations is requested, but no bias is marked "
786 "as shared (share-group > 0)");
789 /* mdrun does not support this (yet), but will check again */
790 if (haveBiasSharingWithinSimulation(awhParams))
793 "You have shared biases within a single simulation, but mdrun does not support "
799 AwhParams* readAwhParams(std::vector<t_inpfile>* inp, warninp_t wi)
801 AwhParams* awhParams;
805 /* Parameters common for all biases */
807 printStringNoNewline(inp, "The way to apply the biasing potential: convolved or umbrella");
808 opt = "awh-potential";
809 awhParams->ePotential = getEnum<AwhPotentialType>(inp, opt.c_str(), wi);
811 printStringNoNewline(inp,
812 "The random seed used for sampling the umbrella center in the case of "
813 "umbrella type potential");
815 awhParams->seed = get_eint(inp, opt, -1, wi);
816 if (awhParams->seed == -1)
818 awhParams->seed = static_cast<int>(gmx::makeRandomSeed());
819 fprintf(stderr, "Setting the AWH bias MC random seed to %" PRId64 "\n", awhParams->seed);
822 printStringNoNewline(inp, "Data output interval in number of steps");
824 awhParams->nstOut = get_eint(inp, opt, 100000, wi);
826 printStringNoNewline(inp, "Coordinate sampling interval in number of steps");
827 opt = "awh-nstsample";
828 awhParams->nstSampleCoord = get_eint(inp, opt, 10, wi);
830 printStringNoNewline(inp, "Free energy and bias update interval in number of samples");
831 opt = "awh-nsamples-update";
832 awhParams->numSamplesUpdateFreeEnergy = get_eint(inp, opt, 10, wi);
834 printStringNoNewline(
835 inp, "When true, biases with share-group>0 are shared between multiple simulations");
836 opt = "awh-share-multisim";
837 awhParams->shareBiasMultisim = (getEnum<Boolean>(inp, opt.c_str(), wi) != Boolean::No);
839 printStringNoNewline(inp, "The number of independent AWH biases");
841 awhParams->numBias = get_eint(inp, opt, 1, wi);
842 /* Check this before starting to read the AWH biases. */
843 if (awhParams->numBias <= 0)
845 gmx_fatal(FARGS, "%s needs to be an integer > 0", opt.c_str());
848 /* Read the parameters specific to each AWH bias */
849 snew(awhParams->awhBiasParams, awhParams->numBias);
851 for (int k = 0; k < awhParams->numBias; k++)
853 bool bComment = (k == 0);
854 std::string prefixawh = formatString("awh%d", k + 1);
855 readBiasParams(inp, &awhParams->awhBiasParams[k], prefixawh, wi, bComment);
861 void checkAwhParams(const AwhParams* awhParams, const t_inputrec* ir, warninp_t wi)
865 checkMtsConsistency(*ir, wi);
868 if (awhParams->nstOut <= 0)
870 auto message = formatString("Not writing AWH output with AWH (%s = %d) does not make sense",
873 warning_error(wi, message);
875 /* This restriction can be removed by changing a flag of print_ebin() */
876 if (ir->nstenergy == 0 || awhParams->nstOut % ir->nstenergy != 0)
878 auto message = formatString(
879 "%s (%d) should be a multiple of nstenergy (%d)", opt.c_str(), awhParams->nstOut, ir->nstenergy);
880 warning_error(wi, message);
883 opt = "awh-nsamples-update";
884 if (awhParams->numSamplesUpdateFreeEnergy <= 0)
886 warning_error(wi, opt + " needs to be an integer > 0");
889 for (int k = 0; k < awhParams->numBias; k++)
891 std::string prefixawh = formatString("awh%d", k + 1);
892 checkBiasParams(&awhParams->awhBiasParams[k], prefixawh, ir, wi);
895 /* Do a final consistency check before returning */
896 checkInputConsistencyAwh(*awhParams, wi);
898 if (ir->init_step != 0)
900 warning_error(wi, "With AWH init-step should be 0");
905 * Gets the period of a pull coordinate.
907 * \param[in] pullCoordParams The parameters for the pull coordinate.
908 * \param[in] pbc The PBC setup
909 * \param[in] intervalLength The length of the AWH interval for this pull coordinate
910 * \returns the period (or 0 if not periodic).
912 static double get_pull_coord_period(const t_pull_coord& pullCoordParams, const t_pbc& pbc, const real intervalLength)
916 if (pullCoordParams.eGeom == PullGroupGeometry::Direction)
918 const real margin = 0.001;
919 // Make dims periodic when the interval covers > 95%
920 const real periodicFraction = 0.95;
922 // Check if the pull direction is along a box vector
923 for (int dim = 0; dim < pbc.ndim_ePBC; dim++)
925 const real boxLength = norm(pbc.box[dim]);
926 const real innerProduct = iprod(pullCoordParams.vec, pbc.box[dim]);
927 if (innerProduct >= (1 - margin) * boxLength && innerProduct <= (1 + margin) * boxLength)
929 if (intervalLength > (1 + margin) * boxLength)
932 "The AWH interval (%f nm) for a pull coordinate is larger than the "
938 if (intervalLength > periodicFraction * boxLength)
945 else if (pullCoordParams.eGeom == PullGroupGeometry::Dihedral)
947 /* The dihedral angle is periodic in -180 to 180 deg */
955 * Checks if the given interval is defined in the correct periodic interval.
957 * \param[in] origin Start value of interval.
958 * \param[in] end End value of interval.
959 * \param[in] period Period (or 0 if not periodic).
960 * \returns true if the end point values are in the correct periodic interval.
962 static bool intervalIsInPeriodicInterval(double origin, double end, double period)
964 return (period == 0) || (std::fabs(origin) <= 0.5 * period && std::fabs(end) <= 0.5 * period);
968 * Checks if a value is within an interval.
970 * \param[in] origin Start value of interval.
971 * \param[in] end End value of interval.
972 * \param[in] period Period (or 0 if not periodic).
973 * \param[in] value Value to check.
974 * \returns true if the value is within the interval.
976 static bool valueIsInInterval(double origin, double end, double period, double value)
984 /* The interval closes within the periodic interval */
985 bIn_interval = (value >= origin) && (value <= end);
989 /* The interval wraps around the periodic boundary */
990 bIn_interval = ((value >= origin) && (value <= 0.5 * period))
991 || ((value >= -0.5 * period) && (value <= end));
996 bIn_interval = (value >= origin) && (value <= end);
1003 * Check if the starting configuration is consistent with the given interval.
1005 * \param[in] awhParams AWH parameters.
1006 * \param[in,out] wi Struct for bookeeping warnings.
1008 static void checkInputConsistencyInterval(const AwhParams* awhParams, warninp_t wi)
1010 for (int k = 0; k < awhParams->numBias; k++)
1012 AwhBiasParams* awhBiasParams = &awhParams->awhBiasParams[k];
1013 for (int d = 0; d < awhBiasParams->ndim; d++)
1015 AwhDimParams* dimParams = &awhBiasParams->dimParams[d];
1016 int coordIndex = dimParams->coordIndex;
1017 double origin = dimParams->origin, end = dimParams->end, period = dimParams->period;
1018 double coordValueInit = dimParams->coordValueInit;
1020 if ((period == 0) && (origin > end))
1023 "For the non-periodic pull coordinates awh%d-dim%d-start (%f) cannot be "
1024 "larger than awh%d-dim%d-end (%f)",
1033 /* Currently we assume symmetric periodic intervals, meaning we use [-period/2, period/2] as the reference interval.
1034 Make sure the AWH interval is within this reference interval.
1036 Note: we could fairly simply allow using a more general interval (e.g. [x, x + period]) but it complicates
1037 things slightly and I don't see that there is a great need for it. It would also mean that the interval would
1038 depend on AWH input. Also, for dihedral angles you would always want the reference interval to be -180, +180,
1039 independent of AWH parameters.
1041 if (!intervalIsInPeriodicInterval(origin, end, period))
1044 "When using AWH with periodic pull coordinate geometries "
1045 "awh%d-dim%d-start (%.8g) and "
1046 "awh%d-dim%d-end (%.8g) should cover at most one period (%.8g) and take "
1047 "values in between "
1048 "minus half a period and plus half a period, i.e. in the interval [%.8g, "
1061 /* Warn if the pull initial coordinate value is not in the grid */
1062 if (!valueIsInInterval(origin, end, period, coordValueInit))
1064 auto message = formatString(
1065 "The initial coordinate value (%.8g) for pull coordinate index %d falls "
1067 "of the sampling nterval awh%d-dim%d-start (%.8g) to awh%d-dim%d-end "
1069 "This can lead to large initial forces pulling the coordinate towards the "
1070 "sampling interval.",
1079 warning(wi, message);
1086 * Sets AWH parameters, for one AWH pull dimension.
1088 * \param[in,out] dimParams AWH dimension parameters.
1089 * \param[in] biasIndex The index of the bias containing this AWH pull dimension.
1090 * \param[in] dimIndex The index of this AWH pull dimension.
1091 * \param[in] pull_params Pull parameters.
1092 * \param[in,out] pull_work Pull working struct to register AWH bias in.
1093 * \param[in] pbc A pbc information structure.
1094 * \param[in] compressibility Compressibility matrix for pressure coupling,
1095 * pass all 0 without pressure coupling.
1096 * \param[in,out] wi Struct for bookeeping warnings.
1099 static void setStateDependentAwhPullDimParams(AwhDimParams* dimParams,
1100 const int biasIndex,
1102 const pull_params_t* pull_params,
1105 const tensor& compressibility,
1108 const t_pull_coord& pullCoordParams = pull_params->coord[dimParams->coordIndex];
1110 if (pullCoordParams.eGeom == PullGroupGeometry::DirectionPBC)
1113 "AWH does not support pull geometry '%s'. "
1114 "If the maximum distance between the groups is always "
1115 "less than half the box size, "
1116 "you can use geometry '%s' instead.",
1117 enumValueToString(PullGroupGeometry::DirectionPBC),
1118 enumValueToString(PullGroupGeometry::Direction));
1121 dimParams->period = get_pull_coord_period(pullCoordParams, pbc, dimParams->end - dimParams->origin);
1122 // We would like to check for scaling, but we don't have the full inputrec available here
1123 if (dimParams->period > 0
1124 && !(pullCoordParams.eGeom == PullGroupGeometry::Angle
1125 || pullCoordParams.eGeom == PullGroupGeometry::Dihedral))
1127 bool coordIsScaled = false;
1128 for (int d2 = 0; d2 < DIM; d2++)
1130 if (pullCoordParams.vec[d2] != 0 && norm2(compressibility[d2]) != 0)
1132 coordIsScaled = true;
1137 std::string mesg = gmx::formatString(
1138 "AWH dimension %d of bias %d is periodic with pull geometry '%s', "
1139 "while you should are applying pressure scaling to the "
1140 "corresponding box vector, this is not supported.",
1143 enumValueToString(pullCoordParams.eGeom));
1144 warning(wi, mesg.c_str());
1148 /* The initial coordinate value, converted to external user units. */
1149 dimParams->coordValueInit = get_pull_coord_value(pull_work, dimParams->coordIndex, &pbc);
1151 dimParams->coordValueInit *= pull_conversion_factor_internal2userinput(&pullCoordParams);
1154 void setStateDependentAwhParams(AwhParams* awhParams,
1155 const pull_params_t& pull_params,
1159 const tensor& compressibility,
1160 const t_grpopts* inputrecGroupOptions,
1161 const real initLambda,
1162 const gmx_mtop_t& mtop,
1165 /* The temperature is not really state depenendent but is not known
1166 * when read_awhParams is called (in get ir).
1167 * It is known first after do_index has been called in grompp.cpp.
1169 if (inputrecGroupOptions->ref_t == nullptr || inputrecGroupOptions->ref_t[0] <= 0)
1171 gmx_fatal(FARGS, "AWH biasing is only supported for temperatures > 0");
1173 for (int i = 1; i < inputrecGroupOptions->ngtc; i++)
1175 if (inputrecGroupOptions->ref_t[i] != inputrecGroupOptions->ref_t[0])
1178 "AWH biasing is currently only supported for identical temperatures for all "
1179 "temperature coupling groups");
1184 set_pbc(&pbc, pbcType, box);
1186 for (int k = 0; k < awhParams->numBias; k++)
1188 AwhBiasParams* awhBiasParams = &awhParams->awhBiasParams[k];
1189 for (int d = 0; d < awhBiasParams->ndim; d++)
1191 AwhDimParams* dimParams = &awhBiasParams->dimParams[d];
1192 if (dimParams->eCoordProvider == AwhCoordinateProviderType::Pull)
1194 setStateDependentAwhPullDimParams(
1195 dimParams, k, d, &pull_params, pull_work, pbc, compressibility, wi);
1199 dimParams->coordValueInit = initLambda;
1200 checkFepLambdaDimDecouplingConsistency(mtop, wi);
1204 checkInputConsistencyInterval(awhParams, wi);
1206 /* Register AWH as external potential with pull (for AWH dimensions that use pull as
1207 * reaction coordinate provider) to check consistency. */
1208 Awh::registerAwhWithPull(*awhParams, pull_work);