STORMM Source Documentation
Loading...
Searching...
No Matches
stormm::structure::MeshForceField< T > Class Template Reference

A class to hold the rules by which the system underlying a mesh object interacts with its surroundings. The mesh object itself may have one probe or an array of them (for an array of potential surfaces), but the constants, rules, and modifying potentials in this object will determine how those probes interact with the underlying system. More...

#include <mesh_forcefield.h>

Public Member Functions

 MeshForceField (VdwCombiningRule mixing_protocol_in=VdwCombiningRule::LORENTZ_BERTHELOT, double coulomb_constant_in=amber_ancient_bioq, double clash_ratio_in=0.0, double clash_distance_in=0.0, int lj_type_count_in=0)
 The object is a convenient place to store constants and arrays of coefficients, nothing more. The constructor accepts these constants as well as sizing information for its arrays and allocates the appropriate amount of memory. Filling the arrays is a matter for the particular mesh object "creating" (perhaps more precise to say "populating") the force field.
 
 MeshForceField (VdwCombiningRule mixing_protocol_in, double clash_ratio_in, double clash_distance_in, const AtomGraph *ag)
 
 MeshForceField (VdwCombiningRule mixing_protocol_in, double clash_ratio_in, double clash_distance_in, const AtomGraph &ag)
 
 MeshForceField (VdwCombiningRule mixing_protocol_in, double coulomb_constant_in, double clash_ratio_in, double clash_distance_in, const std::vector< double > &probe_sigma, const std::vector< double > &probe_epsilon)
 
VdwCombiningRule getCombiningRule () const
 Get the Lennard-Jones mixing rule.
 
double getClashRatio () const
 Get the ratio of the van-der Waals (Lennard-Jones) sigma parameter for interacting pairs of particles at which softcore van-der Waals interactions take over.
 
double getClashDistance () const
 Get the absolute distance at which softcore electrostatic interactions take over.
 
const MeshFFKit< T > data (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get the abstract in whatever templated type (float or double) the object is cast in.
 
const MeshFFKit< double > referenceData (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get the refernce data in double precision, regardless of how the object is cast.
 
const MeshFFKit< void > templateFreeData (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get a template-free form of the object's abstract in its native type.
 
void setCombiningRule (VdwCombiningRule mixing_protocol_in)
 Set the combining rule that will be used to make the probe interact with the receptor on any mesh (with the exception of an electrostatic field). Checks on the validity of this setting must be made by the mesh that assembles the MeshForceField object.
 
void setClashDistance (double clash_distance_in)
 Set the absolute distance at which electrostatic charges will begin to be damped according to the softcore scheme (a harmonic potential with a relative maximum at r = -1.0 takes over).
 
void setClashRatio (double clash_ratio_in)
 Set the ratio of the van-der Waals (Lennard-Jones) sigma parameters at which softcore interactions take over (the interaction switches to a quartic potential with a relative maximum at r = -1.0).
 
void setElecSoftcoreParameter (double coef, int pos)
 Set a specific coefficient in the softcore electrostatic polynomial.
 
void setElecSoftcoreParameters (Interpolant stencil_kind)
 Compute the polynomial for an electrostatic softcore potential suitable for a certain interpolation strategy.
 
void setLJCoefficients (int index, double lja_in, double ljb_in)
 Set the Lennard-Jones pairwise coefficients to describe the interaction of one atom type in the mesh's underlying topology with its probe.
 
void setLJSoftcoreParameter (int index, double coef, int pos)
 Set a specific softcore polynomial coefficient.
 
void setLJSoftcoreParameters (Interpolant stencil_kind)
 Compute softcore polynomials for the loaded Lennard-Jones parameters, ratio of the sigma parameter at which the handoff is to occur, and a given interpolation strategy.
 
void reallocate (const int lj_type_count_in)
 Reallocate the object to hold a different number of Lennard-Jones types. The electrostatic softcore potential coefficients will be held in a temporary array during this process so that they may be preserved.
 
 MeshForceField (const MeshForceField< T > &original)
 The copy and move constructors as well as assignment operators must be written out to repair various pointers.
 
 MeshForceField (MeshForceField< T > &&original)
 
MeshForceFieldoperator= (const MeshForceField< T > &other)
 
MeshForceFieldoperator= (MeshForceField< T > &&other)
 
void setLJCoefficients (const AtomGraph *ag, const std::vector< double > &probe_radius, const std::vector< double > &probe_well_depth)
 Set the Lennard-Jones pairwise coefficients based on the mesh's underlying topology.
 
void setLJCoefficients (const AtomGraph &ag, const std::vector< double > &probe_radius, const std::vector< double > &probe_well_depth)
 
void setLJCoefficients (const AtomGraph *ag, double probe_radius, double probe_well_depth)
 
void setLJCoefficients (const AtomGraph &ag, double probe_radius, double probe_well_depth)
 

Detailed Description

template<typename T>
class stormm::structure::MeshForceField< T >

A class to hold the rules by which the system underlying a mesh object interacts with its surroundings. The mesh object itself may have one probe or an array of them (for an array of potential surfaces), but the constants, rules, and modifying potentials in this object will determine how those probes interact with the underlying system.

Constructor & Destructor Documentation

◆ MeshForceField()

template<typename T>
stormm::structure::MeshForceField< T >::MeshForceField ( const MeshForceField< T > & original)

The copy and move constructors as well as assignment operators must be written out to repair various pointers.

Parameters
originalThe original object to copy or move
otherAn existing object fulfilling the right hand side of an assignment statement

Member Function Documentation

◆ data()

template<typename T>
const MeshFFKit< T > stormm::structure::MeshForceField< T >::data ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get the abstract in whatever templated type (float or double) the object is cast in.

Parameters
tierObtain the data at the level of the CPU host or GPU device

◆ reallocate()

template<typename T>
void stormm::structure::MeshForceField< T >::reallocate ( const int lj_type_count_in)

Reallocate the object to hold a different number of Lennard-Jones types. The electrostatic softcore potential coefficients will be held in a temporary array during this process so that they may be preserved.

Parameters
lj_type_count_inThe new number of Lennard-Jones types to prepare for

◆ referenceData()

template<typename T>
const MeshFFKit< double > stormm::structure::MeshForceField< T >::referenceData ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get the refernce data in double precision, regardless of how the object is cast.

Parameters
tierObtain the data at the level of the CPU host or GPU device

◆ setClashDistance()

template<typename T>
void stormm::structure::MeshForceField< T >::setClashDistance ( double clash_distance_in)

Set the absolute distance at which electrostatic charges will begin to be damped according to the softcore scheme (a harmonic potential with a relative maximum at r = -1.0 takes over).

Parameters
clash_distance_inThe handoff range to set

◆ setClashRatio()

template<typename T>
void stormm::structure::MeshForceField< T >::setClashRatio ( double clash_ratio_in)

Set the ratio of the van-der Waals (Lennard-Jones) sigma parameters at which softcore interactions take over (the interaction switches to a quartic potential with a relative maximum at r = -1.0).

Parameters
clash_ratio_inThe handoff sigma ratio to set

◆ setCombiningRule()

template<typename T>
void stormm::structure::MeshForceField< T >::setCombiningRule ( VdwCombiningRule mixing_protocol_in)

Set the combining rule that will be used to make the probe interact with the receptor on any mesh (with the exception of an electrostatic field). Checks on the validity of this setting must be made by the mesh that assembles the MeshForceField object.

Parameters
mixing_protocol_inThe method to use

◆ setElecSoftcoreParameter()

template<typename T>
void stormm::structure::MeshForceField< T >::setElecSoftcoreParameter ( double coef,
int pos )

Set a specific coefficient in the softcore electrostatic polynomial.

Parameters
coefThe coefficient value to apply
posOrder of the coefficient, e.g. in Ax^3 + Bx^2 + Cx + D, D has order zero and C has order 1

◆ setElecSoftcoreParameters()

template<typename T>
void stormm::structure::MeshForceField< T >::setElecSoftcoreParameters ( Interpolant stencil_kind)

Compute the polynomial for an electrostatic softcore potential suitable for a certain interpolation strategy.

Parameters
stencil_kindPrepare for computing energies that will feed into this type of interpolation stencil. SMOOTHNESS interpolants require a fifth-order softcore function providing for three layers of continuous derivatives at the handoff, FUNCTION_VALUE interpolants, while more expensive to compute, require only third-order softcore polynomials with a single continuous derivative at the handoff.

◆ setLJCoefficients() [1/2]

template<typename T>
void stormm::structure::MeshForceField< T >::setLJCoefficients ( const AtomGraph * ag,
const std::vector< double > & probe_radius,
const std::vector< double > & probe_well_depth )

Set the Lennard-Jones pairwise coefficients based on the mesh's underlying topology.

Overloaded:

  • Provide the topology as a const pointer or const reference
  • Provide the a single value for the probe radius (sigma) and well depth (epsilon) to be combined with the topology's parameters according to the MeshForceField object's rule
Parameters
agThe mesh's underlying topology
probe_radiusThe spherical probe's radius (its sigma parameter if computing a Lennard-Jones non-bonded potential
probe_well_depthThe probe epsilon parameter

◆ setLJCoefficients() [2/2]

template<typename T>
void stormm::structure::MeshForceField< T >::setLJCoefficients ( int index,
double lja_in,
double ljb_in )

Set the Lennard-Jones pairwise coefficients to describe the interaction of one atom type in the mesh's underlying topology with its probe.

Parameters
indexIndex of the atom type in the underlying topology
lja_inThe Lennard-Jones A coefficient, as in U(LJ) = A/r^12 - B/r^6
ljb_inThe Lennard-Jones B coefficient

◆ setLJSoftcoreParameter()

template<typename T>
void stormm::structure::MeshForceField< T >::setLJSoftcoreParameter ( int index,
double coef,
int pos )

Set a specific softcore polynomial coefficient.

Parameters
indexLennard-Jones type index of a mesh atom type for which to apply the coefficient
coefThe coefficient value to apply
posOrder of the coefficient, e.g. in Ax^3 + Bx^2 + Cx + D, D has order zero and C has order 1

◆ setLJSoftcoreParameters()

template<typename T>
void stormm::structure::MeshForceField< T >::setLJSoftcoreParameters ( Interpolant stencil_kind)

Compute softcore polynomials for the loaded Lennard-Jones parameters, ratio of the sigma parameter at which the handoff is to occur, and a given interpolation strategy.

Parameters
stencil_kindPrepare for computing energies that will feed into this type of interpolation stencil.

◆ templateFreeData()

template<typename T>
const MeshFFKit< void > stormm::structure::MeshForceField< T >::templateFreeData ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get a template-free form of the object's abstract in its native type.

Parameters
tierObtain the data at the level of the CPU host or GPU device

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