Gromacs  2024.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members | Public Types | Public Member Functions | Public Attributes
gmx::BiasParams Class Reference

#include <gromacs/applied_forces/awh/biasparams.h>

Description

Constant parameters for the bias.

Public Types

enum  DisableUpdateSkips { DisableUpdateSkips::no, DisableUpdateSkips::yes }
 Switch to turn off update skips, useful for testing. More...
 

Public Member Functions

bool skipUpdates () const
 Check if the parameters permit skipping updates. More...
 
const awh_iveccoverRadius () const
 Returns the radius that needs to be sampled around a point before it is considered covered.
 
bool isSampleCoordStep (int64_t step) const
 Returns whether we should sample the coordinate. More...
 
bool isUpdateFreeEnergyStep (int64_t step) const
 Returns whether we should update the free energy. More...
 
bool isUpdateTargetStep (int64_t step) const
 Returns whether we should update the target distribution. More...
 
bool isCheckCoveringStep (int64_t step) const
 Returns if to do checks for covering in the initial stage. More...
 
bool isCheckHistogramForAnomaliesStep (int64_t step) const
 Returns if to perform checks for anomalies in the histogram. More...
 
 BiasParams (const AwhParams &awhParams, const AwhBiasParams &awhBiasParams, ArrayRef< const DimParams > dimParams, double beta, double mdTimeStep, DisableUpdateSkips disableUpdateSkips, int numSharingSimulations, ArrayRef< const GridAxis > gridAxis, int biasIndex)
 Constructor. More...
 

Public Attributes

const double invBeta
 1/beta = kT in kJ/mol
 
const int numSamplesUpdateFreeEnergy_
 Number of samples per free energy update. More...
 
const AwhTargetType eTarget
 Type of target distribution. More...
 
const bool scaleTargetByMetric
 Scale the target distribution based on the friction metric?
 
const double targetMetricScalingLimit
 The upper limit for the scaling factor based on the friction metric. More...
 
const double freeEnergyCutoffInKT
 Free energy cut-off in kT for cut-off target distribution. More...
 
const double temperatureScaleFactor
 Temperature scaling factor for temperature scaled targed distributions. More...
 
const bool idealWeighthistUpdate
 Update reference weighthistogram using the target distribution? Otherwise use the realized distribution. More...
 
const int numSharedUpdate
 The number of (multi-)simulations sharing the bias update.
 
const double updateWeight
 The probability weight accumulated for each update. More...
 
const double localWeightScaling
 Scaling factor applied to a sample before adding it to the reference weight histogram (= 1, usually). More...
 
const double initialErrorInKT
 Estimated initial free energy error in kT. More...
 
const double initialHistogramSize
 Initial reference weight histogram size. More...
 
const bool convolveForce
 True if we convolve the force, false means use MC between umbrellas. More...
 
const int biasIndex_
 Index of the bias, used as a second random seed and for priting. More...
 

Member Enumeration Documentation

Switch to turn off update skips, useful for testing.

Enumerator
no 

Allow update skips (when supported by the method)

yes 

Disable update skips.

Constructor & Destructor Documentation

gmx::BiasParams::BiasParams ( const AwhParams awhParams,
const AwhBiasParams &  awhBiasParams,
ArrayRef< const DimParams dimParams,
double  beta,
double  mdTimeStep,
DisableUpdateSkips  disableUpdateSkips,
int  numSharingSimulations,
ArrayRef< const GridAxis gridAxis,
int  biasIndex 
)

Constructor.

The local Boltzmann target distibution is defined by 1) Adding the sampled weights instead of the target weights to the reference weight histogram. 2) Scaling the weights of these samples by the beta scaling factor. 3) Setting the target distribution equal the reference weight histogram. This requires the following special update settings: localWeightScaling = targetParam idealWeighthistUpdate = false Note: these variables could in principle be set to something else also for other target distribution types. However, localWeightScaling < 1 is in general expected to give lower efficiency and, except for local Boltzmann, idealWeightHistUpdate = false gives (in my experience) unstable, non-converging results.

Parameters
[in]awhParamsAWH parameters.
[in]awhBiasParamsBias parameters.
[in]dimParamsBias dimension parameters.
[in]beta1/(k_B T) in units of 1/(kJ/mol), should be > 0.
[in]mdTimeStepThe MD time step.
[in]numSharingSimulationsThe number of simulations to share the bias across.
[in]gridAxisThe grid axes.
[in]disableUpdateSkipsIf to disable update skips, useful for testing.
[in]biasIndexIndex of the bias.

Member Function Documentation

bool gmx::BiasParams::isCheckCoveringStep ( int64_t  step) const
inline

Returns if to do checks for covering in the initial stage.

To avoid overhead due to expensive checks, we do not check at every free energy update. However, if checks are performed too rarely the detection of coverings will be delayed, ultimately affecting free energy convergence.

Parameters
[in]stepTime step.
Returns
true at steps where checks should be performed.
Note
Only returns true at free energy update steps.
bool gmx::BiasParams::isCheckHistogramForAnomaliesStep ( int64_t  step) const
inline

Returns if to perform checks for anomalies in the histogram.

To avoid overhead due to expensive checks, we do not check at every free energy update. These checks are only used for warning the user and can be made as infrequently as neccessary without affecting the algorithm itself.

Parameters
[in]stepTime step.
Returns
true at steps where checks should be performed.
Note
Only returns true at free energy update steps.
Todo:
Currently this function just calls isCheckCoveringStep but the checks could be done less frequently.
bool gmx::BiasParams::isSampleCoordStep ( int64_t  step) const
inline

Returns whether we should sample the coordinate.

Parameters
[in]stepThe MD step number.
bool gmx::BiasParams::isUpdateFreeEnergyStep ( int64_t  step) const
inline

Returns whether we should update the free energy.

Parameters
[in]stepThe MD step number.
bool gmx::BiasParams::isUpdateTargetStep ( int64_t  step) const
inline

Returns whether we should update the target distribution.

Parameters
[in]stepThe MD step number.
bool gmx::BiasParams::skipUpdates ( ) const
inline

Check if the parameters permit skipping updates.

Generally, we can skip updates of points that are non-local at the time of the update if we for later times, when the points with skipped updates have become local, know exactly how to apply the previous updates. The free energy updates only depend on local sampling, but the histogram rescaling factors generally depend on the histogram size (all samples). If the histogram size is kept constant or the scaling factors are trivial, this is not a problem. However, if the histogram growth is scaled down by some factor the size at the time of the update needs to be known. It would be fairly simple to, for a deterministically growing histogram, backtrack and calculate this value, but currently we just disallow this case. This is not a restriction because it only affects the local Boltzmann target type for which every update is currently anyway global because the target is always updated globally.

Returns
true when we can skip updates.

Member Data Documentation

const int gmx::BiasParams::biasIndex_

Index of the bias, used as a second random seed and for priting.

const bool gmx::BiasParams::convolveForce

True if we convolve the force, false means use MC between umbrellas.

const AwhTargetType gmx::BiasParams::eTarget

Type of target distribution.

const double gmx::BiasParams::freeEnergyCutoffInKT

Free energy cut-off in kT for cut-off target distribution.

const bool gmx::BiasParams::idealWeighthistUpdate

Update reference weighthistogram using the target distribution? Otherwise use the realized distribution.

const double gmx::BiasParams::initialErrorInKT

Estimated initial free energy error in kT.

const double gmx::BiasParams::initialHistogramSize

Initial reference weight histogram size.

const double gmx::BiasParams::localWeightScaling

Scaling factor applied to a sample before adding it to the reference weight histogram (= 1, usually).

const int gmx::BiasParams::numSamplesUpdateFreeEnergy_

Number of samples per free energy update.

const double gmx::BiasParams::targetMetricScalingLimit

The upper limit for the scaling factor based on the friction metric.

The lower limit is the inverse.

const double gmx::BiasParams::temperatureScaleFactor

Temperature scaling factor for temperature scaled targed distributions.

const double gmx::BiasParams::updateWeight

The probability weight accumulated for each update.


The documentation for this class was generated from the following files: