C++ API Reference

Normally, users of MELD will not use the C++ API. This documentation is primarily for developers.

All inputs to the MELD plugin are assumed to be in units of:

  • kilojoules per mole

  • nanometers

  • degrees

class CalcMeldForceKernel : public KernelImpl
#include <MeldKernels.h>

Public Functions

inline CalcMeldForceKernel(std::string name, const OpenMM::Platform &platform)
virtual void initialize(const OpenMM::System &system, const MeldForce &force) = 0

Initialize the kernel.

Parameters
  • system – the System this kernel will be applied to

  • force – the MeldForce this kernel will be used for

virtual double execute(OpenMM::ContextImpl &context, bool includeForces, bool includeEnergy) = 0

Execute the kernel to calculate the forces and/or energy.

Parameters
  • context – the context in which to execute this kernel

  • includeForces – true if forces should be calculated

  • includeEnergy – true if the energy should be calculated

Returns

the potential energy due to the force

virtual void copyParametersToContext(OpenMM::ContextImpl &context, const MeldForce &force) = 0

Public Static Functions

static inline std::string Name()
class CollectionInfo

Public Functions

inline CollectionInfo()
inline CollectionInfo(std::vector<int> group_indices, int n_active)

Public Members

std::vector<int> group_indices
int n_active
class DistanceRestraintInfo

Public Functions

inline DistanceRestraintInfo()
inline DistanceRestraintInfo(int particle1, int particle2, float r1, float r2, float r3, float r4, float force_constant, int global_index)

Public Members

int particle1
int particle2
float r1
float r2
float r3
float r4
float force_constant
int global_index
class DistProfileRestraintInfo

Public Functions

inline DistProfileRestraintInfo()
inline DistProfileRestraintInfo(int atom1, int atom2, float rMin, float rMax, int nBins, std::vector<double> a0, std::vector<double> a1, std::vector<double> a2, std::vector<double> a3, float scaleFactor, int globalIndex)

Public Members

int atom1
int atom2
int nBins
float rMin
float rMax
float scaleFactor
std::vector<double> a0
std::vector<double> a1
std::vector<double> a2
std::vector<double> a3
int globalIndex
class GMMRestraintInfo

Public Functions

inline GMMRestraintInfo()
inline GMMRestraintInfo(int nPairs, int nComponents, int globalIndex, float scale, std::vector<int> atomIndices, std::vector<double> weights, std::vector<double> means, std::vector<double> precisionOnDiagonal, std::vector<double> precisionOffDiagonal)

Public Members

int nPairs
int nComponents
int globalIndex
float scale
std::vector<int> atomIndices
std::vector<double> weights
std::vector<double> means
std::vector<double> precisionOnDiagonal
std::vector<double> precisionOffDiagonal
class GridPotentialInfo

Public Functions

inline GridPotentialInfo()
inline GridPotentialInfo(std::vector<double> potential, double originx, double originy, double originz, double gridx, double gridy, double gridz, int nx, int ny, int nz, int density_index)

Public Members

std::vector<double> potential
double originx
double originy
double originz
double gridx
double gridy
double gridz
int nx
int ny
int nz
int density_index
class GridPotentialRestraintInfo

Public Functions

inline GridPotentialRestraintInfo()
inline GridPotentialRestraintInfo(std::vector<int> particle, std::vector<double> mu, std::vector<double> gridpos_x, std::vector<double> gridpos_y, std::vector<double> gridpos_z, int globalIndex)

Public Members

std::vector<int> particle
std::vector<double> mu
std::vector<double> gridpos_x
std::vector<double> gridpos_y
std::vector<double> gridpos_z
int globalIndex
class GroupInfo

Public Functions

inline GroupInfo()
inline GroupInfo(std::vector<int> restraint_indices, int n_active)

Public Members

std::vector<int> restraint_indices
int n_active
class HyperbolicDistanceRestraintInfo

Public Functions

inline HyperbolicDistanceRestraintInfo()
inline HyperbolicDistanceRestraintInfo(int particle1, int particle2, float r1, float r2, float r3, float r4, float force_constant, float asymptote, int global_index)

Public Members

int particle1
int particle2
float r1
float r2
float r3
float r4
float force_constant
float asymptote
int global_index
class MeldForce : public Force
#include <MeldForce.h>

This is the MELD Force.

Public Functions

MeldForce()

Constructors.

MeldForce(int numRDCAlignments, float rdcScaleFactor)
void updateParametersInContext(OpenMM::Context &context)

Update the per-restraint parameters in a Context to match those stored in this Force object.

This method provides an efficient method to update certain parameters in an existing Context without needing to reinitialize it. Simply call modifyDistanceRestaint(), modifyTorsionRestraint(), modifyDistProfileRestraint(), or modifyTorsProfileRestraint() to modify the parameters of a restraint, then call updateParametersInContext() to copy them over to the Context.

This method has several limitations. The only information it updates is the values of per-restraint parameters. All other aspects of the Force (such as the energy function) are unaffected and can only be changed by reinitializing the Context. The set of particles involved in a restraint cannot be changed, nor can new restraints be added.

int getNumRDCAlignments() const
Returns

The number of RDC alignment tensors in this force.

float getRDCScaleFactor() const
Returns

The RDC scale factor.

bool containsParticle(int particle) const
Returns

A bool indicating if particle is involved in MELD force

int getNumRDCRestraints() const
Returns

The number of RDC restraints.

int getNumDistRestraints() const
Returns

The number of distance restraints.

int getNumHyperbolicDistRestraints() const
Returns

The number of hyperbolic distance restraints.

int getNumTorsionRestraints() const
Returns

The number of torsion restraints.

int getNumDistProfileRestraints() const
Returns

The number of distance profile restraints.

int getNumDistProfileRestParams() const
Returns

The number of distance profile restraint parameters.

int getNumTorsProfileRestraints() const
Returns

The number of torsion profile restraints.

int getNumTorsProfileRestParams() const
Returns

The number of torsion profile restraint parameters.

int getNumGMMRestraints() const
Returns

The number of GMM restraints

int getNumGridPotentials() const
Returns

The number of grid potentials

int getNumGridPotentialRestraints() const
Returns

The number of grid potential restraints

int getNumTotalRestraints() const
Returns

The total number of distance and torsion restraints.

int getNumGroups() const
Returns

The number of restraint groups.

int getNumCollections() const
Returns

The number of collections of restraint groups.

void getRDCRestraintParameters(int index, int &particle1, int &particle2, int &alignment, float &kappa, float &obs, float &tol, float &quad_cut, float &force_constant, int &globalIndex) const

Get the parameters of an RDC restraint.

Parameters
  • index – the index to retrieve

  • particle1 – the first particle in the RDC

  • particle2 – the second particle in the RDC

  • alignment – the index of the alignment tensor to use

  • kappa – the kappa parameter in Hz nm^3

  • obs – the observed dipolar coupling in Hz

  • tol – the tolerance in Hz

  • quad_cut – number of Hz for transition from quadratic to linear

  • force_constant – the force constant in kJ mol^-1 Hz^-2

  • globalIndex – the global index of the restraint

void getDistanceRestraintParams(int index, int &atom1, int &atom2, float &r1, float &r2, float &r3, float &r4, float &forceConstant, int &globalIndex) const

Get the parameters for a distance restraint.

See addDistanceRestraint() for more details about the parameters.

Parameters
  • index – the index of the restraint

  • atom1 – the first atom

  • atom2 – the second atom

  • r1 – the upper bound of region 1

  • r2 – the upper bound of region 2

  • r3 – the upper bound of region 3

  • r4 – the upper bound of region 4

  • forceConstant – the force constant

  • globalIndex – the global index of the restraint

void getHyperbolicDistanceRestraintParams(int index, int &atom1, int &atom2, float &r1, float &r2, float &r3, float &r4, float &forceConstant, float &asymptote, int &globalIndex) const

Get the parameters for a hyperbolic distance restraint.

See addHyperbolicDistanceRestraint() for more details about the parameters.

Parameters
  • index – the index of the restraint

  • atom1 – the first atom

  • atom2 – the second atom

  • r1 – the upper bound of region 1

  • r2 – the upper bound of region 2

  • r3 – the upper bound of region 3

  • r4 – the upper bound of region 4

  • forceConstant – the force constant for region 1

  • asymptote – the asymptotic energy in region 4

  • globalIndex – the global index of the restraint

void getTorsionRestraintParams(int index, int &atom1, int &atom2, int &atom3, int &atom4, float &phi, float &deltaPhi, float &forceConstant, int &globalIndex) const

Get the parameters for a torsion restraint.

See addTorsionRestraint() for more details about the parameters.

Parameters
  • index – the index of the restraint

  • atom1 – the first atom

  • atom2 – the second atom

  • atom3 – the third atom

  • atom4 – the fourth atom

  • phi – the equilibrium torsion (degrees)

  • deltaPhi – the deltaPhi parameter (degrees)

  • forceConstant – the force constant

  • globalIndex – the global index of the restraint

void getDistProfileRestraintParams(int index, int &atom1, int &atom2, float &rMin, float &rMax, int &nBins, std::vector<double> &a0, std::vector<double> &a1, std::vector<double> &a2, std::vector<double> &a3, float &scaleFactor, int &globalIndex) const

Get the parameters for a distance profile restraint.

See addDistProfileRestraint() for more details about the parameters.

Parameters
  • index – the index of the restraint

  • atom1 – the first atom

  • atom2 – the second atom

  • rMin – the lower bound of the restraint

  • rMax – the upper bound of the restraint

  • nBins – the number of bins

  • aN – the Nth spline parameter where N is in (0,1,2,3)

  • scaleFactor – the scale factor

  • globalIndex – the global index of the restraint

void getTorsProfileRestraintParams(int index, int &atom1, int &atom2, int &atom3, int &atom4, int &atom5, int &atom6, int &atom7, int &atom8, int &nBins, std::vector<double> &a0, std::vector<double> &a1, std::vector<double> &a2, std::vector<double> &a3, std::vector<double> &a4, std::vector<double> &a5, std::vector<double> &a6, std::vector<double> &a7, std::vector<double> &a8, std::vector<double> &a9, std::vector<double> &a10, std::vector<double> &a11, std::vector<double> &a12, std::vector<double> &a13, std::vector<double> &a14, std::vector<double> &a15, float &scaleFactor, int &globalIndex) const

Get the parameters for a torsion profile restraint.

Parameters
  • index – the index of the restraint

  • atom1 – the first atom

  • atom2 – the second atom

  • atom3 – the third atom

  • atom4 – the fourth atom

  • atom5 – the fifth atom

  • atom6 – the sixth atom

  • atom7 – the seventh atom

  • atom8 – the eighth atom

  • nBins – the number of bins

  • aN – the Nth spline parameter where N is in (0,1,…,14,15)

  • scaleFactor – the scale factor

  • globalIndex – the global index of the restraint

void getGMMRestraintParams(int index, int &nPairs, int &nComponents, float &scale, std::vector<int> &atomIndices, std::vector<double> &weights, std::vector<double> &means, std::vector<double> &precisionOnDiagonal, std::vector<double> &precisionOffDiagonal, int &globalIndex) const

Get the parameters for a GMM restraint.

See addGMMRestraint() for more details about the parameters.

Parameters
  • index – the index of the restraint

  • nPairs – the number of atom pairs

  • nComponents – the number of GMM components

  • scale – the overall scaling applied to the forces and energies

  • atomIndices – the vector of atom indices

  • weights – the vector of weights

  • means – the vector of means in nm

  • precisionOnDiagonal – the diagonals of the precision matrix in nm^(-2)

  • precisionOffDiagonal – the off-diagonals of the precision matrix in nm^(-2)

  • globalIndex – the global index of the restraint

void getGroupParams(int index, std::vector<int> &indices, int &numActive) const

Get the parameters for a group of restraints.

Parameters
  • index – the index of the group

  • indices – the indices of the restraints in the group

  • numActive – the number of active restraints in the group

void getCollectionParams(int index, std::vector<int> &indices, int &numActive) const

Get the parameters for a collection of restraint groups.

Parameters
  • index – the index of the collection

  • indices – the indices of the groups in the collection

  • numActive – the number of active groups in the collection

void getGridPotentialParams(int index, std::vector<double> &potential, float &originx, float &originy, float &originz, float &gridx, float &gridy, float &gridz, int &nx, int &ny, int &nz) const
void getGridPotentialRestraintParams(int index, std::vector<int> &atom, std::vector<double> &mu, std::vector<double> &gridpos_x, std::vector<double> &gridpos_y, std::vector<double> &gridpos_z, int &globalIndex) const
int addRDCRestraint(int particle1, int particle2, int alignment, float kappa, float obs, float tol, float quad_cut, float force_constant)
Parameters
  • particle1 – the first particle in the RDC

  • particle2 – the second particle in the RDC

  • alignment – the index of the alignment tensor to use

  • kappa – the kappa parameter in Hz nm^3

  • obs – the observed dipolar coupling in Hz

  • tol – the tolerance in Hz

  • quad_cut – number of Hz for transition from quadratic to linear

  • force_constant – the force constant in kJ mol^-1 Hz^-2

Returns

the index of the restraint that was created

void modifyRDCRestraint(int index, int particle1, int particle2, int alignment, float kappa, float obs, float tol, float quad_cut, float force_constant)
Parameters
  • index – the index of the restraint to modify

  • particle1 – the first particle in the RDC

  • particle2 – the second particle in the RDC

  • alignment – the index of the alignment tensor to use

  • kappa – the kappa parameter in Hz nm^3

  • obs – the observed dipolar coupling in Hz

  • tol – the tolerance in Hz

  • quad_cut – number of Hz for transition from quadratic to linear

  • force_constant – the force constant in kJ mol^-1 Hz^-2

int addDistanceRestraint(int particle1, int particle2, float r1, float r2, float r3, float r4, float force_constant)

Create a new distance restraint.

There are five regions:

I: r < r1

II: r1 < r < r2

III: r2 < r < r3

IV: r3 < r < r4

V: r4 < r

The energy is linear in regions I and V, quadratic in II and IV, and zero in III.

There is special handling of the restraints if particle1=-1. In this case, the energy and force will be zero. However, when selecting active restraints within a group, the energy will be treated as MAX_FLOAT. This behavior is to facilitate support for peak mapping, which needs a way for a restraint to be effectively ignored.

Parameters
  • particle1 – the first atom

  • particle2 – the second atom

  • r1 – the upper bound of region 1

  • r2 – the upper bound of region 2

  • r3 – the upper bound of region 3

  • r4 – the upper bound of region 4

  • forceConstant – the force constant

Returns

the index of the restraint that was created

void modifyDistanceRestraint(int index, int particle1, int particle2, float r1, float r2, float r3, float r4, float force_constant)

Modify an existing distance restraint.

See addDistanceRestraint() for more details about the parameters.

Parameters
  • index – the index of the restraint

  • particle1 – the first atom

  • particle2 – the second atom

  • r1 – the upper bound of region 1

  • r2 – the upper bound of region 2

  • r3 – the upper bound of region 3

  • r4 – the upper bound of region 4

  • forceConstant – the force constant

int addHyperbolicDistanceRestraint(int particle1, int particle2, float r1, float r2, float r3, float r4, float force_constant, float asymptote)

Create a new hyperbolic distance restraint.

There are five regions:

I: r < r1

II: r1 < r < r2

III: r2 < r < r3

IV: r3 < r < r4

V: r4 < r

The energy is linear in region I, quadratic in II and IV, and zero in III.

The energy is hyperbolic in region V, with an asymptotic value set by the parameter asymptote. The energy will be 1/3 of the asymptotic value at r=r4. The distance between r3 and r4 controls the steepness of the potential.

Parameters
  • particle1 – the first atom

  • particle2 – the second atom

  • r1 – the upper bound of region 1

  • r2 – the upper bound of region 2

  • r3 – the upper bound of region 3

  • r4 – the upper bound of region 4

  • forceConstant – the force constant in regions I and II

  • asymptote – the asymptotic value in region V, also controls the steepness in region IV.

Returns

the index of the restraint that was created

void modifyHyperbolicDistanceRestraint(int index, int particle1, int particle2, float r1, float r2, float r3, float r4, float force_constant, float asymptote)

Modify an existing hyperbolic distance restraint.

See addHyperbolicDistanceRestraint() for more details about the parameters.

Parameters
  • index – the index of the restraint

  • particle1 – the first atom

  • particle2 – the second atom

  • r1 – the upper bound of region 1

  • r2 – the upper bound of region 2

  • r3 – the upper bound of region 3

  • r4 – the upper bound of region 4

  • forceConstant – the force constant

  • asymptote – the asymptotic value

int addTorsionRestraint(int atom1, int atom2, int atom3, int atom4, float phi, float deltaPhi, float forceConstant)

Create a new torsion restraint.

If (x - phi) < -deltaPhi: E = 1/2 * forceConstant * (x - phi + deltaPhi)^2

Else if (x - phi) > deltaPhi: E = 1/2 * forceConstant * (x - phi - deltaPhi)^2

Else: E = 0

Parameters
  • atom1 – the first atom

  • atom2 – the second atom

  • atom3 – the third atom

  • atom4 – the fourth atom

  • phi – the equilibrium torsion (degrees)

  • deltaPhi – the deltaPhi parameter (degrees)

  • forceConstant – the force constant

Returns

the index of the restraint that was created

void modifyTorsionRestraint(int index, int atom1, int atom2, int atom3, int atom4, float phi, float deltaPhi, float forceConstant)

Modify an existing torsion restraint.

See addTorsionRestraint() for more details about the parameters.

Parameters
  • index – the index of the restraint

  • atom1 – the first atom

  • atom2 – the second atom

  • atom3 – the third atom

  • atom4 – the fourth atom

  • phi – the equilibrium torsion (degrees)

  • deltaPhi – the deltaPhi parameter (degrees)

  • forceConstant – the force constant

int addDistProfileRestraint(int atom1, int atom2, float rMin, float rMax, int nBins, std::vector<double> a0, std::vector<double> a1, std::vector<double> a2, std::vector<double> a3, float scaleFactor)

Create a new distance profile restraint.

bin = floor( (r - rMin) / (rMax - rMin) * nBins) )

binWidth = (rMax - rMin) / nBins

t = (r - bin * binWidth + rMin) / binWidth;

E = scaleFactor * (a0 + a1 * t + a2 * t^2 + a3 * t^3)

Parameters
  • atom1 – the first atom

  • atom2 – the second atom

  • rMin – the lower bound of the restraint

  • rMax – the upper bound of the restraint

  • nBins – the number of bins

  • aN – the Nth spline parameter where N is in (0,1,2,3)

  • scaleFactor – the scale factor

Returns

the index of the restraint that was created

void modifyDistProfileRestraint(int index, int atom1, int atom2, float rMin, float rMax, int nBins, std::vector<double> a0, std::vector<double> a1, std::vector<double> a2, std::vector<double> a3, float scaleFactor)

Modify an existing distance profile restraint.

See addDistProfileRestraint() for more details about the parameters.

Parameters
  • index – the index of the restraint

  • atom1 – the first atom

  • atom2 – the second atom

  • rMin – the lower bound of the restraint

  • rMax – the upper bound of the restraint

  • nBins – the number of bins

  • aN – the Nth spline parameter where N is in (0,1,2,3)

  • scaleFactor – the scale factor

int addTorsProfileRestraint(int atom1, int atom2, int atom3, int atom4, int atom5, int atom6, int atom7, int atom8, int nBins, std::vector<double> a0, std::vector<double> a1, std::vector<double> a2, std::vector<double> a3, std::vector<double> a4, std::vector<double> a5, std::vector<double> a6, std::vector<double> a7, std::vector<double> a8, std::vector<double> a9, std::vector<double> a10, std::vector<double> a11, std::vector<double> a12, std::vector<double> a13, std::vector<double> a14, std::vector<double> a15, float scaleFactor)

Create a new torsion profile restraint.

Parameters
  • atom1 – the first atom

  • atom2 – the second atom

  • atom3 – the third atom

  • atom4 – the fourth atom

  • atom5 – the fifth atom

  • atom6 – the sixth atom

  • atom7 – the seventh atom

  • atom8 – the eighth atom

  • nBins – the number of bins

  • aN – the Nth spline parameter where N is in (0,1,…,14,15)

  • scaleFactor – the scale factor

Returns

the index of the restraint that was created

void modifyTorsProfileRestraint(int index, int atom1, int atom2, int atom3, int atom4, int atom5, int atom6, int atom7, int atom8, int nBins, std::vector<double> a0, std::vector<double> a1, std::vector<double> a2, std::vector<double> a3, std::vector<double> a4, std::vector<double> a5, std::vector<double> a6, std::vector<double> a7, std::vector<double> a8, std::vector<double> a9, std::vector<double> a10, std::vector<double> a11, std::vector<double> a12, std::vector<double> a13, std::vector<double> a14, std::vector<double> a15, float scaleFactor)

Modify an existing torsion profile restraint.

Parameters
  • index – the index of the restraint

  • atom1 – the first atom

  • atom2 – the second atom

  • atom3 – the third atom

  • atom4 – the fourth atom

  • atom5 – the fifth atom

  • atom6 – the sixth atom

  • atom7 – the seventh atom

  • atom8 – the eighth atom

  • nBins – the number of bins

  • aN – the Nth spline parameter where N is in (0,1,…,14,15)

  • scaleFactor – the scale factor

int addGMMRestraint(int nPairs, int nComponents, float scale, std::vector<int> atomIndices, std::vector<double> weights, std::vector<double> means, std::vector<double> precisionOnDiagonal, std::vector<double> precisionOffDiagonal)

Create a new GMM restraint.

This is a Gaussian Mixture Model restraint involving a series of distances.

The energy has the form:

E = w1 N1 exp(-0.5 (r-u1)^T P1 (r-u1)) + w2 N2 exp(-0.5 (r-u2)^T P2 (r-u2)) + …

where: w1, w2, … are the nComponents weights, N1, N2, … are automatically calculated normalization factors r is the vector of distances for the atom pairs in atomIndices u1, u2, … are the mean vectors for each component P1, P2, … are the precision (inverse covariance) matrices for each component

The memory layout is as follows: atomIndices -> [pair1_atom_1, pair1_atom_2, pair2_atom_1, pair2_atom_2, …] weights -> [wa, wb, …] means -> [ma1, ma2, …, mb1, mb2, …] precisionOnDiagonal -> [Pa11, Pa22, …, Pb11, Pb22, …] precisionOffDiagonal -> [Pa12, Pa13, …, Pa23, …, Pb12, Pb13, …, Pb23, …]

where a, b, … are the different components and 1, 2, … are the different distances.

atomIndices.size() must be 2 * nPairs weights.size() must be nComponents means.size() must be nPairs * nComponents precisionOnDiagonal.size() must be nPairs * nComponents precisionOffDiagonal.size() must be nPairs * nComponents * (nComponents - 1) / 2

Parameters
  • nPairs – the number of atom pairs

  • nComponents – the number of GMM components

  • scale – the overall scaling applied to the forces and energies

  • atomIndices – the vector of atom indices

  • weights – the vector of weights

  • means – the vector of means in nm

  • precisionOnDiagonal – the diagonals of the precision matrix in nm^(-2)

  • precisionOffDiagonal – the off-diagonals of the precision matrix in nm^(-2)

Returns

the index of the restraint that was created

void modifyGMMRestraint(int index, int nPairs, int nComponents, float scale, std::vector<int> atomIndices, std::vector<double> weights, std::vector<double> means, std::vector<double> precisionOnDiagonal, std::vector<double> precisionOffDiagonal)

Modify an existing GMM restraint.

See addGMMRestraint() for more details about the parameters.

Parameters
  • index – the index of the restraint

  • nPairs – the number of atom pairs

  • nComponents – the number of GMM components

  • scale – the overall scaling applied to the forces and energies

  • atomIndices – the vector of atom indices

  • weights – the vector of weights

  • means – the vector of means in nm

  • precisionOnDiagonal – the diagonals of the precision matrix in nm^(-2)

  • precisionOffDiagonal – the off-diagonals of the precision matrix in nm^(-2)

int addGroup(std::vector<int> restraint_indices, int n_active)

Create a new group of restraints.

Parameters
  • restraint_indices – the indices of the restraints in the group

  • n_active – the number of active restraints in the group

Returns

the index of the group that was created

void modifyGroupNumActive(int index, int n_active)

Modify the number of active restraints in a group.

Parameters
  • index – the index of the group to modify

  • n_active – the new number of active restraints

int addCollection(std::vector<int> group_indices, int n_active)

Create a new collection of restraint groups.

Parameters
  • group_indices – the indices of the groups in the collection

  • n_active – the number of active groups in the collection

Returns

the index of the collection that was created

void modifyCollectionNumActive(int index, int n_active)

Modify the number of active groups in a collection.

Parameters
  • index – the index of the collection to modify

  • n_active – the new number of active groups

bool usesPeriodicBoundaryConditions() const override
void addGridPotential(std::vector<double> potential, float originx, float originy, float originz, float gridx, float gridy, float gridz, int nx, int ny, int nz, int density_index)

Add a new grid potential to the system.

void modifyGridPotential(int index, std::vector<double> potential, float originx, float originy, float originz, float gridx, float gridy, float gridz, int nx, int ny, int nz)

Modify a grid potential.

int addGridPotentialRestraint(std::vector<int> particle, std::vector<double> mu, std::vector<double> gridpos_x, std::vector<double> gridpos_y, std::vector<double> gridpos_z)

Add a grid potential restraint.

void modifyGridPotentialRestraint(int index, std::vector<int> particle, std::vector<double> mu, std::vector<double> gridpos_x, std::vector<double> gridpos_y, std::vector<double> gridpos_z)
std::vector<std::pair<int, int>> getBondedParticles() const

Protected Functions

OpenMM::ForceImpl *createImpl() const

Private Functions

void updateMeldParticleSet()

Private Members

int n_restraints
int n_rdc_alignments
float rdcScaleFactor
std::vector<RDCRestraintInfo> rdcRestraints
std::vector<DistanceRestraintInfo> distanceRestraints
std::vector<HyperbolicDistanceRestraintInfo> hyperbolicDistanceRestraints
std::vector<TorsionRestraintInfo> torsions
std::vector<DistProfileRestraintInfo> distProfileRestraints
std::vector<TorsProfileRestraintInfo> torsProfileRestraints
std::vector<GMMRestraintInfo> gmmRestraints
std::vector<GridPotentialInfo> gridPotentials
std::vector<GridPotentialRestraintInfo> gridPotentialRestraints
std::vector<GroupInfo> groups
std::vector<CollectionInfo> collections
std::set<int> meldParticleSet
bool isDirty
class MeldForceImpl : public ForceImpl
#include <MeldForceImpl.h>

Public Functions

MeldForceImpl(const MeldForce &owner)
~MeldForceImpl()
void initialize(OpenMM::ContextImpl &context)
inline const MeldForce &getOwner() const
inline void updateContextState(OpenMM::ContextImpl &context)
double calcForcesAndEnergy(OpenMM::ContextImpl &context, bool includeForces, bool includeEnergy, int groups)
inline std::map<std::string, double> getDefaultParameters()
std::vector<std::string> getKernelNames()
void updateParametersInContext(OpenMM::ContextImpl &context)
inline std::vector<std::pair<int, int>> getBondedParticles() const

Private Members

const MeldForce &owner
OpenMM::Kernel kernel
class MeldForceProxy : public SerializationProxy
#include <MeldForceProxy.h>

This is a proxy for serializing MeldForce objects.

Public Functions

MeldForceProxy()
void serialize(const void *object, SerializationNode &node) const
void *deserialize(const SerializationNode &node) const
class RDCRestraintInfo

Public Functions

inline RDCRestraintInfo()
inline RDCRestraintInfo(int particle1, int particle2, int alignment, float kappa, float obs, float tol, float quad_cut, float force_constant, int global_index)

Public Members

int particle1
int particle2
int alignment
float kappa
float obs
float tol
float quad_cut
float force_constant
int global_index
class TorsionRestraintInfo

Public Functions

inline TorsionRestraintInfo()
inline TorsionRestraintInfo(int atom1, int atom2, int atom3, int atom4, float phi, float deltaPhi, float forceConstant, int globalIndex)

Public Members

int atom1
int atom2
int atom3
int atom4
float phi
float deltaPhi
float forceConstant
int globalIndex
class TorsProfileRestraintInfo

Public Functions

inline TorsProfileRestraintInfo()
inline TorsProfileRestraintInfo(int atom1, int atom2, int atom3, int atom4, int atom5, int atom6, int atom7, int atom8, int nBins, std::vector<double> a0, std::vector<double> a1, std::vector<double> a2, std::vector<double> a3, std::vector<double> a4, std::vector<double> a5, std::vector<double> a6, std::vector<double> a7, std::vector<double> a8, std::vector<double> a9, std::vector<double> a10, std::vector<double> a11, std::vector<double> a12, std::vector<double> a13, std::vector<double> a14, std::vector<double> a15, float scaleFactor, int globalIndex)

Public Members

int atom1
int atom2
int atom3
int atom4
int atom5
int atom6
int atom7
int atom8
int nBins
int globalIndex
float scaleFactor
std::vector<double> a0
std::vector<double> a1
std::vector<double> a2
std::vector<double> a3
std::vector<double> a4
std::vector<double> a5
std::vector<double> a6
std::vector<double> a7
std::vector<double> a8
std::vector<double> a9
std::vector<double> a10
std::vector<double> a11
std::vector<double> a12
std::vector<double> a13
std::vector<double> a14
std::vector<double> a15
namespace MeldPlugin
namespace OpenMM
file MeldForce.h
#include “openmm/Force.h”
#include “openmm/Vec3.h”
#include <map>
#include <vector>
#include <set>
file MeldForceImpl.h
#include “MeldForce.h
#include “openmm/internal/ForceImpl.h”
#include “openmm/Kernel.h”
#include <utility>
#include <set>
#include <vector>
#include <string>
file MeldForceProxy.h
#include “openmm/serialization/SerializationProxy.h”
file MeldKernels.h
#include “MeldForce.h
#include “openmm/KernelImpl.h”
#include “openmm/System.h”
#include “openmm/Platform.h”
#include <set>
#include <string>
#include <vector>
file windowsExportMeld.h

Defines

OPENMM_EXPORT_MELD
dir include
dir include
dir internal
dir openmmapi
dir serialization