|
Gromacs
2026.0-dev-20251111-8c1ac59
|
#include <gromacs/mdtypes/enerdata.h>
Accumulates free-energy foreign lambda energies and dH/dlamba.
Public Member Functions | |
| ForeignLambdaTerms (const gmx::EnumerationArray< FreeEnergyPerturbationCouplingType, std::vector< double >> *allLambdas) | |
| Constructor. More... | |
| int | numLambdas () const |
| Returns the number of foreign lambda values. | |
| gmx::ArrayRef< const double > | foreignLambdas (FreeEnergyPerturbationCouplingType couplingType) |
| Returns the foreign lambda values for the given perturbation type. | |
| 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, FreeEnergyPerturbationCouplingType couplingType, double energy, real dvdl) |
Adds an energy and dV/dl constribution to lambda list index listIndex. More... | |
| void | accumulate (int listIndex, double energy, const gmx::EnumerationArray< FreeEnergyPerturbationCouplingType, real > &dvdl) |
Adds an energy and dV/dl constribution to lambda list index listIndex. More... | |
| void | finalizePotentialContributions (const gmx::EnumerationArray< FreeEnergyPerturbationCouplingType, double > &dvdlLinear, gmx::ArrayRef< const real > lambda, const t_lambda &fepvals) |
| Finalizes the potential (non-kinetic) terms. More... | |
| void | finalizeKineticContributions (const gmx::EnumerationArray< InteractionFunction, 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 gmx::MpiComm &mpiComm) const |
| Returns a pair of lists of deltaH and dH/dlambda. More... | |
| void | zeroAllTerms () |
| Sets all terms to 0. | |
| ForeignLambdaTerms::ForeignLambdaTerms | ( | const gmx::EnumerationArray< FreeEnergyPerturbationCouplingType, std::vector< double >> * | allLambdas | ) |
Constructor.
| [in] | allLambdas | The list of lambda values for all lambda components, can be nullptr |
|
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.
|
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.
|
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.
|
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 | ( | const gmx::EnumerationArray< InteractionFunction, real > & | energyTerms, |
| double | dhdlMass, | ||
| gmx::ArrayRef< const real > | lambda, | ||
| const t_lambda & | fepvals | ||
| ) |
Accumulates the kinetic and constraint free-energy contributions.
| [in] | energyTerms | List of energy terms, pass term in gmx_enerdata_t |
| [in] | dhdlMass | The mass dependent contribution to dH/dlambda |
| [in] | lambda | Lambda values for the efptNR contribution types |
| [in] | fepvals | Free-energy parameters |
| void ForeignLambdaTerms::finalizePotentialContributions | ( | const gmx::EnumerationArray< FreeEnergyPerturbationCouplingType, 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.
| [in] | dvdlLinear | List of dV/dlambda contributions of size efptNR with depend linearly on lambda |
| [in] | lambda | Lambda values for the efptNR contribution types |
| [in] | fepvals | Free-energy parameters |
| std::pair< std::vector< double >, std::vector< double > > ForeignLambdaTerms::getTerms | ( | const gmx::MpiComm & | mpiComm | ) | 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).
| [in] | mpiComm | Communication object for my group |
1.8.5