Gromacs  2021-sycl
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members | Public Member Functions
ForeignLambdaTerms Class Reference

#include <gromacs/mdtypes/enerdata.h>

Description

Accumulates free-energy foreign lambda energies and dH/dlamba.

Public Member Functions

 ForeignLambdaTerms (int numLambdas)
 Constructor. More...
 
int numLambdas () const
 Returns the number of foreign lambda values.
 
double deltaH (int lambdaIndex) const
 Returns the H(lambdaIndex) - H(lambda_current)
 
gmx::ArrayRef< double > energies ()
 Returns a list of partial energies, the part which depends on lambda), current lambda in entry 0, foreign lambda i in entry 1+i. More...
 
gmx::ArrayRef< const double > energies () const
 Returns a list of partial energies, the part which depends on lambda), current lambda in entry 0, foreign lambda i in entry 1+i. More...
 
void accumulate (int listIndex, double energy, double dvdl)
 Adds an energy and dV/dl constribution to lambda list index listIndex. More...
 
void finalizePotentialContributions (gmx::ArrayRef< const double > dvdlLinear, gmx::ArrayRef< const real > lambda, const t_lambda &fepvals)
 Finalizes the potential (non-kinetic) terms. More...
 
void finalizeKineticContributions (gmx::ArrayRef< const real > energyTerms, double dhdlMass, gmx::ArrayRef< const real > lambda, const t_lambda &fepvals)
 Accumulates the kinetic and constraint free-energy contributions. More...
 
std::pair< std::vector< double >
, std::vector< double > > 
getTerms (const t_commrec *cr) const
 Returns a pair of lists of deltaH and dH/dlambda. More...
 
void zeroAllTerms ()
 Sets all terms to 0.
 

Constructor & Destructor Documentation

ForeignLambdaTerms::ForeignLambdaTerms ( int  numLambdas)

Constructor.

Parameters
[in]numLambdasThe number of foreign lambda values

Member Function Documentation

void ForeignLambdaTerms::accumulate ( int  listIndex,
double  energy,
double  dvdl 
)
inline

Adds an energy and dV/dl constribution to lambda list index listIndex.

This should only be used for terms with non-linear dependence on lambda The value passed as listIndex should be 0 for the current lambda and 1+i for foreign lambda index i.

gmx::ArrayRef<double> ForeignLambdaTerms::energies ( )
inline

Returns a list of partial energies, the part which depends on lambda), current lambda in entry 0, foreign lambda i in entry 1+i.

Note: the potential terms needs to be finalized before calling this method.

gmx::ArrayRef<const double> ForeignLambdaTerms::energies ( ) const
inline

Returns a list of partial energies, the part which depends on lambda), current lambda in entry 0, foreign lambda i in entry 1+i.

Note: the potential terms needs to be finalized before calling this method.

void ForeignLambdaTerms::finalizeKineticContributions ( gmx::ArrayRef< const real energyTerms,
double  dhdlMass,
gmx::ArrayRef< const real lambda,
const t_lambda &  fepvals 
)

Accumulates the kinetic and constraint free-energy contributions.

Parameters
[in]energyTermsList of energy terms, pass term in gmx_enerdata_t
[in]dhdlMassThe mass dependent contribution to dH/dlambda
[in]lambdaLambda values for the efptNR contribution types
[in]fepvalsFree-energy parameters
void ForeignLambdaTerms::finalizePotentialContributions ( gmx::ArrayRef< const double >  dvdlLinear,
gmx::ArrayRef< const real lambda,
const t_lambda &  fepvals 
)

Finalizes the potential (non-kinetic) terms.

Note: This can be called multiple times during the same force calculations without affecting the results.

Parameters
[in]dvdlLinearList of dV/dlambda contributions of size efptNR with depend linearly on lambda
[in]lambdaLambda values for the efptNR contribution types
[in]fepvalsFree-energy parameters
std::pair< std::vector< double >, std::vector< double > > ForeignLambdaTerms::getTerms ( const t_commrec *  cr) const

Returns a pair of lists of deltaH and dH/dlambda.

Both lists are of size numLambdas() and are indexed with the lambda index.

Note: should only be called after the object has been finalized by a call to accumulateLinearPotentialComponents() (is asserted).

Parameters
[in]crCommunication record, used to reduce the terms when !=nullptr

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