RASPA3 3.0.12
A molecular simulation code for computing adsorption and diffusion in nanoporous materials
Loading...
Searching...
No Matches
RunningEnergy Struct Referenceexport

Accumulates energy components during simulation. More...

Public Member Functions

 RunningEnergy ()
 Default constructor initializing all energy components to zero.
 
std::string printMC () const
 Generates a formatted string of Monte Carlo energy components.
 
std::string printMD () const
 Generates a formatted string of molecular dynamics energy components.
 
std::string printMCDiff (RunningEnergy &other) const
 Generates a formatted string comparing energy components between two Monte Carlo states.
 
std::string printMDDiff (RunningEnergy &other) const
 Generates a formatted string comparing energy components between two molecular dynamics states.
 
std::string printMC (const std::string &label) const
 Generates a labeled formatted string of Monte Carlo energy components.
 
std::string printMD (const std::string &label, double referenceEnergy) const
 Generates a labeled formatted string of molecular dynamics energy components, including conserved energy and drift.
 
nlohmann::json jsonMC () const
 Generates a JSON object of Monte Carlo energy components.
 
nlohmann::json jsonMD () const
 Generates a JSON object of molecular dynamics energy components.
 
double potentialEnergy () const
 Computes the total potential energy.
 
double CoulombEnergy () const
 Computes the total Coulombic energy.
 
double kineticEnergy () const
 Computes the total kinetic energy.
 
double ewaldFourier () const
 Computes the total Ewald Fourier energy.
 
double conservedEnergy () const
 Computes the conserved energy in molecular dynamics simulations.
 
double dudlambda (double lambda) const
 Computes the derivative of the potential energy with respect to lambda.
 
void zero ()
 Resets all energy components to zero.
 
RunningEnergyoperator+= (const RunningEnergy &b)
 
RunningEnergyoperator-= (const RunningEnergy &b)
 
RunningEnergy operator- () const
 

Public Attributes

std::uint64_t versionNumber {1}
 Version number for serialization.
 
double externalFieldVDW
 Energy from van der Waals interactions with external field.
 
double frameworkMoleculeVDW
 Energy from van der Waals interactions between framework and molecules.
 
double moleculeMoleculeVDW
 Energy from van der Waals interactions between molecules.
 
double externalFieldCharge
 Energy from Coulomb interactions with external field.
 
double frameworkMoleculeCharge
 Energy from Coulomb interactions between framework and molecules.
 
double moleculeMoleculeCharge
 Energy from Coulomb interactions between molecules.
 
double ewald_fourier
 Fourier component of Ewald sum for Coulomb interactions.
 
double ewald_self
 Self-interaction correction in Ewald sum.
 
double ewald_exclusion
 Exclusion term in Ewald sum for Coulomb interactions.
 
double bond
 
double ureyBradley
 
double bend
 
double inversionBend
 
double outOfPlaneBend
 
double torsion
 
double improperTorsion
 
double bondBond
 
double bondBend
 
double bondTorsion
 
double bendBend
 
double bendTorsion
 
double intraVDW
 Intramolecular van der Waals energy.
 
double intraCoul
 Intramolecular Coulomb energy.
 
double tail
 Tail correction energy for van der Waals interactions.
 
double polarization
 Energy contribution from polarization effects.
 
double dudlambdaVDW
 Derivative of van der Waals energy with respect to lambda.
 
double dudlambdaCharge
 Derivative of Coulomb energy with respect to lambda (real space).
 
double dudlambdaEwald
 Derivative of Coulomb energy with respect to lambda (Ewald sum).
 
double translationalKineticEnergy
 Translational kinetic energy.
 
double rotationalKineticEnergy
 Rotational kinetic energy.
 
double NoseHooverEnergy
 Energy associated with Nose-Hoover thermostat/barostat.
 

Friends

Archive< std::ofstream > & operator<< (Archive< std::ofstream > &archive, const RunningEnergy &c)
 
Archive< std::ifstream > & operator>> (Archive< std::ifstream > &archive, RunningEnergy &c)
 

Detailed Description

Accumulates energy components during simulation.

The RunningEnergy struct stores various components of energy calculated during a simulation run, including potential and kinetic energies, as well as contributions from various interactions such as van der Waals and Coulomb forces. It provides methods to calculate total energies, perform arithmetic operations, and output energy information.

Constructor & Destructor Documentation

◆ RunningEnergy()

RunningEnergy::RunningEnergy ( )
inline

Default constructor initializing all energy components to zero.

Initializes all the energy components and related terms to zero, preparing the RunningEnergy object for accumulation during simulation.

Member Function Documentation

◆ conservedEnergy()

double RunningEnergy::conservedEnergy ( ) const
inline

Computes the conserved energy in molecular dynamics simulations.

Sums up the total potential energy, kinetic energy, and energy associated with the Nose-Hoover thermostat/barostat to compute the conserved energy in molecular dynamics simulations.

Returns
The conserved energy.

◆ CoulombEnergy()

double RunningEnergy::CoulombEnergy ( ) const
inline

Computes the total Coulombic energy.

Sums up the Coulombic energy components, including framework-molecule Coulomb interactions, molecule-molecule Coulomb interactions, and Ewald sums.

Returns
The total Coulombic energy.

◆ dudlambda()

double RunningEnergy::dudlambda ( double  lambda) const
inline

Computes the derivative of the potential energy with respect to lambda.

Calculates the derivative of the potential energy with respect to the scaling parameter lambda, considering both van der Waals and Coulombic interactions.

Parameters
lambdaThe scaling parameter.
Returns
The derivative of the potential energy with respect to lambda.

◆ ewaldFourier()

double RunningEnergy::ewaldFourier ( ) const
inline

Computes the total Ewald Fourier energy.

Sums up the Fourier components of the Ewald sum, including the Fourier term, self-interaction correction, and exclusion term.

Returns
The total Ewald Fourier energy.

◆ jsonMC()

nlohmann::json RunningEnergy::jsonMC ( ) const

Generates a JSON object of Monte Carlo energy components.

Returns a JSON object containing the total potential energy and individual energy components relevant for Monte Carlo simulations.

Returns
A JSON object representing the energy components in Monte Carlo simulations.

◆ jsonMD()

nlohmann::json RunningEnergy::jsonMD ( ) const

Generates a JSON object of molecular dynamics energy components.

Returns a JSON object containing the total potential and kinetic energies, as well as individual energy components relevant for molecular dynamics simulations.

Returns
A JSON object representing the energy components in molecular dynamics simulations.

◆ kineticEnergy()

double RunningEnergy::kineticEnergy ( ) const
inline

Computes the total kinetic energy.

Sums up the translational and rotational kinetic energies.

Returns
The total kinetic energy.

◆ potentialEnergy()

double RunningEnergy::potentialEnergy ( ) const
inline

Computes the total potential energy.

Sums up all the potential energy components, including interactions with external fields, framework-molecule interactions, molecule-molecule interactions, Ewald sums, intramolecular energies, tail corrections, and polarization.

Returns
The total potential energy.

◆ printMC() [1/2]

std::string RunningEnergy::printMC ( ) const

Generates a formatted string of Monte Carlo energy components.

Returns a string containing the total potential energy and individual energy components relevant for Monte Carlo simulations, formatted for display.

Returns
A string representing the energy components in Monte Carlo simulations.

◆ printMC() [2/2]

std::string RunningEnergy::printMC ( const std::string &  label) const

Generates a labeled formatted string of Monte Carlo energy components.

Returns a string containing the total potential energy and individual energy components relevant for Monte Carlo simulations, with an additional label for identification.

Parameters
labelA string label to include in the output.
Returns
A labeled string representing the energy components in Monte Carlo simulations.

◆ printMCDiff()

std::string RunningEnergy::printMCDiff ( RunningEnergy other) const

Generates a formatted string comparing energy components between two Monte Carlo states.

Computes the difference between the current RunningEnergy and another instance, and returns a string that shows the comparison of energy components, highlighting any drifts or discrepancies.

Parameters
otherAnother RunningEnergy instance to compare with.
Returns
A string representing the differences in energy components between two Monte Carlo states.

◆ printMD() [1/2]

std::string RunningEnergy::printMD ( ) const

Generates a formatted string of molecular dynamics energy components.

Returns a string containing the total potential and kinetic energies, as well as individual energy components relevant for molecular dynamics simulations, formatted for display.

Returns
A string representing the energy components in molecular dynamics simulations.

◆ printMD() [2/2]

std::string RunningEnergy::printMD ( const std::string &  label,
double  referenceEnergy 
) const

Generates a labeled formatted string of molecular dynamics energy components, including conserved energy and drift.

Returns a string containing the total potential and kinetic energies, conserved energy, drift from a reference energy, and individual energy components relevant for molecular dynamics simulations, with an additional label for identification.

Parameters
labelA string label to include in the output.
referenceEnergyThe reference energy to compute the drift.
Returns
A labeled string representing the energy components in molecular dynamics simulations.

◆ printMDDiff()

std::string RunningEnergy::printMDDiff ( RunningEnergy other) const

Generates a formatted string comparing energy components between two molecular dynamics states.

Computes the difference between the current RunningEnergy and another instance, and returns a string that shows the comparison of energy components, highlighting any drifts or discrepancies in molecular dynamics simulations.

Parameters
otherAnother RunningEnergy instance to compare with.
Returns
A string representing the differences in energy components between two molecular dynamics states.

◆ zero()

void RunningEnergy::zero ( )
inline

Resets all energy components to zero.

Sets all the energy components and related terms to zero, effectively clearing the accumulated energies.


The documentation for this struct was generated from the following file: