Simpatico  v1.10
List of all members | Public Member Functions | Protected Member Functions | Protected Attributes
DdMd::PairPotential Class Referenceabstract

Detailed Description

Abstract base class for computing nonbonded pair forces and energies.

A PairPotential has a private CellList and PairList which it can use to calculate nonbonded pair forces. Concrete subclasses of PairPotential are constructed as instances of the PairPotentialImpl class template.

Definition at line 44 of file ddMd/potentials/pair/PairPotential.h.

#include <PairPotential.h>

Inheritance diagram for DdMd::PairPotential:
DdMd::Potential Util::ParamComposite Util::ParamComponent Util::Serializable Util::MpiFileIo DdMd::PairPotentialImpl< Interaction >

Public Member Functions

 PairPotential (Simulation &simulation)
 Constructor. More...
 
 PairPotential ()
 Default constructor. More...
 
virtual ~PairPotential ()
 Destructor. More...
 
double skin () const
 Get value of the pair list skin. More...
 
double cutoff () const
 Get value of the pair list cutoff (maxPairCutoff + skin). More...
 
int methodId () const
 Return integer id for algorithm (0=PAIR, 1=CELL, 2=NSQ) More...
 
Initialization and Serialization
void associate (Domain &domain, Boundary &boundary, AtomStorage &storage)
 Associate with related objects. More...
 
virtual void setNAtomType (int nAtomType)=0
 Set the maximum number of atom types. More...
 
void initialize (const Boundary &maxBoundary, double skin, int pairCapacity)
 Set parameters and allocate memory. More...
 
void readParameters (std::istream &in)
 Initialize, by reading parameters and allocating memory for PairList. More...
 
virtual void loadParameters (Serializable::IArchive &ar)
 Load parameters for PairList from archive, and allocate memory. More...
 
virtual void save (Serializable::OArchive &ar)
 Save internal state to an output archive. More...
 
void setMethodId (int methodId)
 Set id to specify algorithm for energy, force calculations. More...
 
Interaction interface
virtual double pairEnergy (double rsq, int iAtomType, int jAtomType) const =0
 Return energy for a single pair. More...
 
virtual double pairForceOverR (double rsq, int iAtomType, int jAtomType) const =0
 Return force / separation for a single pair. More...
 
virtual void set (std::string name, int i, int j, double value)=0
 Modify a parameter, identified by a string. More...
 
virtual double get (std::string name, int i, int j) const =0
 Get a parameter value, identified by a string. More...
 
virtual double maxPairCutoff () const =0
 Return maximum cutoff. More...
 
virtual std::string interactionClassName () const =0
 Return pair interaction class name (e.g., "LJPair"). More...
 
Pair and Cell Lists.
void buildCellList ()
 Build a cell list. More...
 
void buildPairList ()
 Build the Verlet Pair list. More...
 
virtual void computePairEnergies (MPI::Intracomm &communicator)=0
 Compute pair energies on all processors. More...
 
DMatrix< double > pairEnergies () const
 Return total pair energies, from all processors. More...
 
void unsetPairEnergies ()
 Mark pair energy as unknown (nullify). More...
 
void computeNPair (MPI::Intracomm &communicator)
 Compute twice the number of pairs within the force cutoff. More...
 
int nPair () const
 Return twice the number of pairs within the specified force cutoff. More...
 
CellListcellList ()
 Get the CellList by reference. More...
 
PairListpairList ()
 Get the PairList by reference. More...
 
- Public Member Functions inherited from DdMd::Potential
 Potential ()
 Constructor. More...
 
virtual ~Potential ()
 Destructor. More...
 
void setReverseUpdateFlag (bool reverseUpdateFlag)
 Set flag to identify if reverse communication is enabled. More...
 
bool reverseUpdateFlag () const
 Get flag to identify if reverse communication is enabled. More...
 
virtual void computeForces ()=0
 Add force contributions to all atomic forces. More...
 
virtual void computeEnergy (MPI::Intracomm &communicator)=0
 Compute potential energy on all processors. More...
 
double energy () const
 Return the total potential, from all processors. More...
 
void unsetEnergy ()
 Mark the energy as unknown (nullify). More...
 
bool isEnergySet () const
 Is the energy set (known)? More...
 
virtual void computeStress (MPI::Intracomm &communicator)
 Compute stress on all processors. More...
 
virtual void computeForcesAndStress (MPI::Intracomm &communicator)
 Compute forces and stress for all processors. More...
 
Tensor stress () const
 Return the stress tensor. More...
 
double pressure () const
 Return the pressure. More...
 
void unsetStress ()
 Mark the stress as unknown (nullify). More...
 
bool isStressSet () const
 Is the stress set (known)? More...
 
virtual bool isValid (MPI::Intracomm &communicator) const
 Is the potential in a valid internal state? More...
 
- Public Member Functions inherited from Util::ParamComposite
 ParamComposite ()
 Constructor. More...
 
 ParamComposite (const ParamComposite &other)
 Copy constructor. More...
 
 ParamComposite (int capacity)
 Constructor. More...
 
virtual ~ParamComposite ()
 Virtual destructor. More...
 
void resetParam ()
 Resets ParamComposite to its empty state. More...
 
virtual void readParam (std::istream &in)
 Read the parameter file block. More...
 
virtual void readParamOptional (std::istream &in)
 Read optional parameter file block. More...
 
virtual void writeParam (std::ostream &out)
 Write all parameters to an output stream. More...
 
virtual void load (Serializable::IArchive &ar)
 Load all parameters from an input archive. More...
 
virtual void loadOptional (Serializable::IArchive &ar)
 Load an optional ParamComposite. More...
 
void saveOptional (Serializable::OArchive &ar)
 Saves isActive flag, and then calls save() iff isActive is true. More...
 
void readParamComposite (std::istream &in, ParamComposite &child, bool next=true)
 Add and read a required child ParamComposite. More...
 
void readParamCompositeOptional (std::istream &in, ParamComposite &child, bool next=true)
 Add and attempt to read an optional child ParamComposite. More...
 
template<typename Type >
ScalarParam< Type > & read (std::istream &in, const char *label, Type &value)
 Add and read a new required ScalarParam < Type > object. More...
 
template<typename Type >
ScalarParam< Type > & readOptional (std::istream &in, const char *label, Type &value)
 Add and read a new optional ScalarParam < Type > object. More...
 
template<typename Type >
CArrayParam< Type > & readCArray (std::istream &in, const char *label, Type *value, int n)
 Add and read a required C array parameter. More...
 
template<typename Type >
CArrayParam< Type > & readOptionalCArray (std::istream &in, const char *label, Type *value, int n)
 Add and read an optional C array parameter. More...
 
template<typename Type >
DArrayParam< Type > & readDArray (std::istream &in, const char *label, DArray< Type > &array, int n)
 Add and read a required DArray < Type > parameter. More...
 
template<typename Type >
DArrayParam< Type > & readOptionalDArray (std::istream &in, const char *label, DArray< Type > &array, int n)
 Add and read an optional DArray < Type > parameter. More...
 
template<typename Type , int N>
FArrayParam< Type, N > & readFArray (std::istream &in, const char *label, FArray< Type, N > &array)
 Add and read a required FArray < Type, N > array parameter. More...
 
template<typename Type , int N>
FArrayParam< Type, N > & readOptionalFArray (std::istream &in, const char *label, FArray< Type, N > &array)
 Add and read an optional FArray < Type, N > array parameter. More...
 
template<typename Type >
CArray2DParam< Type > & readCArray2D (std::istream &in, const char *label, Type *value, int m, int n, int np)
 Add and read a required CArray2DParam < Type > 2D C-array. More...
 
template<typename Type >
CArray2DParam< Type > & readOptionalCArray2D (std::istream &in, const char *label, Type *value, int m, int n, int np)
 Add and read an optional CArray2DParam < Type > 2D C-array parameter. More...
 
template<typename Type >
DMatrixParam< Type > & readDMatrix (std::istream &in, const char *label, DMatrix< Type > &matrix, int m, int n)
 Add and read a required DMatrix < Type > matrix parameter. More...
 
template<typename Type >
DMatrixParam< Type > & readOptionalDMatrix (std::istream &in, const char *label, DMatrix< Type > &matrix, int m, int n)
 Add and read an optional DMatrix < Type > matrix parameter. More...
 
template<typename Type >
DSymmMatrixParam< Type > & readDSymmMatrix (std::istream &in, const char *label, DMatrix< Type > &matrix, int n)
 Add and read a required symmetrix DMatrix. More...
 
template<typename Type >
DSymmMatrixParam< Type > & readOptionalDSymmMatrix (std::istream &in, const char *label, DMatrix< Type > &matrix, int n)
 Add and read an optional DMatrix matrix parameter. More...
 
BeginreadBegin (std::istream &in, const char *label, bool isRequired=true)
 Add and read a class label and opening bracket. More...
 
EndreadEnd (std::istream &in)
 Add and read the closing bracket. More...
 
BlankreadBlank (std::istream &in)
 Add and read a new Blank object, representing a blank line. More...
 
void loadParamComposite (Serializable::IArchive &ar, ParamComposite &child, bool next=true)
 Add and load a required child ParamComposite. More...
 
void loadParamCompositeOptional (Serializable::IArchive &ar, ParamComposite &child, bool next=true)
 Add and load an optional child ParamComposite if isActive. More...
 
template<typename Type >
ScalarParam< Type > & loadParameter (Serializable::IArchive &ar, const char *label, Type &value, bool isRequired)
 Add and load a new ScalarParam < Type > object. More...
 
template<typename Type >
ScalarParam< Type > & loadParameter (Serializable::IArchive &ar, const char *label, Type &value)
 Add and load new required ScalarParam < Type > object. More...
 
template<typename Type >
CArrayParam< Type > & loadCArray (Serializable::IArchive &ar, const char *label, Type *value, int n, bool isRequired)
 Add a C array parameter and load its elements. More...
 
template<typename Type >
CArrayParam< Type > & loadCArray (Serializable::IArchive &ar, const char *label, Type *value, int n)
 Add and load a required CArrayParam< Type > array parameter. More...
 
template<typename Type >
DArrayParam< Type > & loadDArray (Serializable::IArchive &ar, const char *label, DArray< Type > &array, int n, bool isRequired)
 Add an load a DArray < Type > array parameter. More...
 
template<typename Type >
DArrayParam< Type > & loadDArray (Serializable::IArchive &ar, const char *label, DArray< Type > &array, int n)
 Add and load a required DArray< Type > array parameter. More...
 
template<typename Type , int N>
FArrayParam< Type, N > & loadFArray (Serializable::IArchive &ar, const char *label, FArray< Type, N > &array, bool isRequired)
 Add and load an FArray < Type, N > fixed-size array parameter. More...
 
template<typename Type , int N>
FArrayParam< Type, N > & loadFArray (Serializable::IArchive &ar, const char *label, FArray< Type, N > &array)
 Add and load a required FArray < Type > array parameter. More...
 
template<typename Type >
CArray2DParam< Type > & loadCArray2D (Serializable::IArchive &ar, const char *label, Type *value, int m, int n, int np, bool isRequired)
 Add and load a CArray2DParam < Type > C 2D array parameter. More...
 
template<typename Type >
CArray2DParam< Type > & loadCArray2D (Serializable::IArchive &ar, const char *label, Type *value, int m, int n, int np)
 Add and load a required < Type > matrix parameter. More...
 
template<typename Type >
DMatrixParam< Type > & loadDMatrix (Serializable::IArchive &ar, const char *label, DMatrix< Type > &matrix, int m, int n, bool isRequired)
 Add and load a DMatrixParam < Type > matrix parameter. More...
 
template<typename Type >
DMatrixParam< Type > & loadDMatrix (Serializable::IArchive &ar, const char *label, DMatrix< Type > &matrix, int m, int n)
 Add and load a required DMatrixParam < Type > matrix parameter. More...
 
template<typename Type >
DSymmMatrixParam< Type > & loadDSymmMatrix (Serializable::IArchive &ar, const char *label, DMatrix< Type > &matrix, int n, bool isRequired)
 Add and load a symmetric DSymmMatrixParam < Type > matrix parameter. More...
 
template<typename Type >
DSymmMatrixParam< Type > & loadDSymmMatrix (Serializable::IArchive &ar, const char *label, DMatrix< Type > &matrix, int n)
 Add and load a required DSymmMatrixParam < Type > matrix parameter. More...
 
void addParamComposite (ParamComposite &child, bool next=true)
 Add a child ParamComposite object to the format array. More...
 
BeginaddBegin (const char *label)
 Add a Begin object representing a class name and bracket. More...
 
EndaddEnd ()
 Add a closing bracket. More...
 
BlankaddBlank ()
 Create and add a new Blank object, representing a blank line. More...
 
std::string className () const
 Get class name string. More...
 
bool isRequired () const
 Is this ParamComposite required in the input file? More...
 
bool isActive () const
 Is this parameter active? More...
 
- Public Member Functions inherited from Util::ParamComponent
virtual ~ParamComponent ()
 Destructor. More...
 
void setIndent (const ParamComponent &parent, bool next=true)
 Set indent level. More...
 
std::string indent () const
 Return indent string for this object (string of spaces). More...
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 Serialize this ParamComponent as a string. More...
 
- Public Member Functions inherited from Util::Serializable
virtual ~Serializable ()
 Destructor. More...
 
- Public Member Functions inherited from Util::MpiFileIo
 MpiFileIo ()
 Constructor. More...
 
 MpiFileIo (const MpiFileIo &other)
 Copy constructor. More...
 
bool isIoProcessor () const
 Can this processor do file I/O ? More...
 
void setIoCommunicator (MPI::Intracomm &communicator)
 Set the communicator. More...
 
void clearCommunicator ()
 Clear (nullify) the communicator. More...
 
bool hasIoCommunicator () const
 Does this object have an associated MPI communicator? More...
 
MPI::Intracomm & ioCommunicator () const
 Get the MPI communicator by reference. More...
 

Protected Member Functions

Boundaryboundary ()
 Get the PairList by const reference. More...
 
Domaindomain ()
 Get the PairList by const reference. More...
 
AtomStoragestorage ()
 Get the AtomStorage by reference. More...
 
void setPairEnergies (DMatrix< double > pairEnergies)
 Set values for pair energies. More...
 
- Protected Member Functions inherited from DdMd::Potential
void setEnergy (double energy)
 Set a value for the total energy. More...
 
void setStress (const Tensor &stress)
 Set a value for the total stress. More...
 
void incrementPairStress (const Vector &f, const Vector &dr, Tensor &stress) const
 Add a pair contribution to the stress tensor. More...
 
void reduceEnergy (double localEnergy, MPI::Intracomm &communicator)
 Add local energies from all processors, set energy on master. More...
 
void reduceStress (Tensor &localStress, MPI::Intracomm &communicator)
 Add local stresses from all processors, set total on master. More...
 
- Protected Member Functions inherited from Util::ParamComposite
void setClassName (const char *className)
 Set class name string. More...
 
void setIsRequired (bool isRequired)
 Set or unset the isActive flag. More...
 
void setIsActive (bool isActive)
 Set or unset the isActive flag. More...
 
void setParent (ParamComponent &param, bool next=true)
 Set this to the parent of a child component. More...
 
void addComponent (ParamComponent &param, bool isLeaf=true)
 Add a new ParamComponent object to the format array. More...
 
template<typename Type >
ScalarParam< Type > & add (std::istream &in, const char *label, Type &value, bool isRequired=true)
 Add a new required ScalarParam < Type > object. More...
 
template<typename Type >
CArrayParam< Type > & addCArray (std::istream &in, const char *label, Type *value, int n, bool isRequired=true)
 Add (but do not read) a required C array parameter. More...
 
template<typename Type >
DArrayParam< Type > & addDArray (std::istream &in, const char *label, DArray< Type > &array, int n, bool isRequired=true)
 Add (but do not read) a DArray < Type > parameter. More...
 
template<typename Type , int N>
FArrayParam< Type, N > & addFArray (std::istream &in, const char *label, FArray< Type, N > &array, bool isRequired=true)
 Add (but do not read) a FArray < Type, N > array parameter. More...
 
template<typename Type >
CArray2DParam< Type > & addCArray2D (std::istream &in, const char *label, Type *value, int m, int n, int np, bool isRequired=true)
 Add (but do not read) a CArray2DParam < Type > 2D C-array. More...
 
template<typename Type >
DMatrixParam< Type > & addDMatrix (std::istream &in, const char *label, DMatrix< Type > &matrix, int m, int n, bool isRequired=true)
 Add and read a required DMatrix < Type > matrix parameter. More...
 
- Protected Member Functions inherited from Util::ParamComponent
 ParamComponent ()
 Constructor. More...
 
 ParamComponent (const ParamComponent &other)
 Copy constructor. More...
 

Protected Attributes

CellList cellList_
 CellList to construct PairList or calculate nonbonded pair forces. More...
 
PairList pairList_
 Verlet pair list, to calculate nonbonded pair forces. More...
 
Boundary maxBoundary_
 Boundary used to allocate space for the cell list. More...
 
double skin_
 Difference between pairlist cutoff and pair potential cutoff. More...
 
double cutoff_
 Minimum cell size = pair potential cutoff + skin. More...
 
int nCellCut_
 Approximate number of cells per cutoff distance in each direction. More...
 
int pairCapacity_
 Maximum number of nonbonded pairs in pair list. More...
 

Additional Inherited Members

- Public Types inherited from Util::Serializable
typedef BinaryFileOArchive OArchive
 Type of output archive used by save method. More...
 
typedef BinaryFileIArchive IArchive
 Type of input archive used by load method. More...
 
- Static Public Member Functions inherited from Util::ParamComponent
static void initStatic ()
 Initialize static echo member to false. More...
 
static void setEcho (bool echo=true)
 Enable or disable echoing for all subclasses of ParamComponent. More...
 
static bool echo ()
 Get echo parameter. More...
 

Constructor & Destructor Documentation

DdMd::PairPotential::PairPotential ( Simulation simulation)

Constructor.

This is the constructor used during a simulation, and is passed the parent Simulation object as a parameter.

Parameters
simulationparent Simulation object.

Definition at line 43 of file PairPotential.cpp.

References Util::ParamComposite::setClassName().

DdMd::PairPotential::PairPotential ( )

Default constructor.

This constructor is provided to simplify unit testing.

Definition at line 28 of file PairPotential.cpp.

References Util::ParamComposite::setClassName().

DdMd::PairPotential::~PairPotential ( )
virtual

Destructor.

Definition at line 69 of file PairPotential.cpp.

Member Function Documentation

void DdMd::PairPotential::associate ( Domain domain,
Boundary boundary,
AtomStorage storage 
)

Associate with related objects.

Required iff object instantiated with default constructor. If required, it must be called before initialize, readParameters, or loadParameters.

Parameters
domainassociated Domain object.
boundaryassociated Boundary object.
storageassociated AtomStorage object.

Definition at line 58 of file PairPotential.cpp.

References boundary(), domain(), and storage().

virtual void DdMd::PairPotential::setNAtomType ( int  nAtomType)
pure virtual

Set the maximum number of atom types.

Required iff object instantiated with default constructor. If required, it must be called before initialize, readParameters, or loadParameters.

Parameters
nAtomTypemaximum number of atom types (max typeId + 1)

Implemented in DdMd::PairPotentialImpl< Interaction >.

void DdMd::PairPotential::initialize ( const Boundary maxBoundary,
double  skin,
int  pairCapacity 
)

Set parameters and allocate memory.

This method initializes the object by assigning parameter values and allocating memory. It provides a programmatic alternative to readParameters.

Parameters
maxBoundarylargest expected Boundary (used for allocation)
skinpair list skin length
pairCapacitymaximum number of pairs per processor

Definition at line 75 of file PairPotential.cpp.

References cutoff_, maxBoundary_, maxPairCutoff(), pairCapacity_, skin(), and skin_.

void DdMd::PairPotential::readParameters ( std::istream &  in)
virtual

Initialize, by reading parameters and allocating memory for PairList.

Parameters
ininput parameter stream.

Reimplemented from Util::ParamComposite.

Reimplemented in DdMd::PairPotentialImpl< Interaction >.

Definition at line 89 of file PairPotential.cpp.

References cutoff_, maxBoundary_, maxPairCutoff(), nCellCut_, pairCapacity_, and skin_.

Referenced by DdMd::PairPotentialImpl< Interaction >::readParameters().

void DdMd::PairPotential::loadParameters ( Serializable::IArchive ar)
virtual

Load parameters for PairList from archive, and allocate memory.

Parameters
arinput/loading archive

Reimplemented from Util::ParamComposite.

Reimplemented in DdMd::PairPotentialImpl< Interaction >.

Definition at line 103 of file PairPotential.cpp.

References cutoff_, maxBoundary_, nCellCut_, pairCapacity_, and skin_.

Referenced by DdMd::PairPotentialImpl< Interaction >::loadParameters().

void DdMd::PairPotential::save ( Serializable::OArchive ar)
virtual
void DdMd::PairPotential::setMethodId ( int  methodId)
inline

Set id to specify algorithm for energy, force calculations.

The default algorithm, methodId=0, uses a pair list that is constructed using a cell list.

Parameters
methodIdalgorithm id: 0=pair list, 1=cell list, 2=N^2 loop.

Definition at line 396 of file ddMd/potentials/pair/PairPotential.h.

virtual double DdMd::PairPotential::pairEnergy ( double  rsq,
int  iAtomType,
int  jAtomType 
) const
pure virtual

Return energy for a single pair.

Parameters
rsqsquare distance between atoms in pair
iAtomTypeatom type index of 1st atom
jAtomTypeatom type index of 2nd atom
Returns
energy of pair

Implemented in DdMd::PairPotentialImpl< Interaction >.

virtual double DdMd::PairPotential::pairForceOverR ( double  rsq,
int  iAtomType,
int  jAtomType 
) const
pure virtual

Return force / separation for a single pair.

Parameters
rsqsquare distance between atoms in pair
iAtomTypeatom type index of 1st atom
jAtomTypeatom type index of 2nd atom
Returns
repulsive force (< 0 if attractive) over distance

Implemented in DdMd::PairPotentialImpl< Interaction >.

virtual void DdMd::PairPotential::set ( std::string  name,
int  i,
int  j,
double  value 
)
pure virtual

Modify a parameter, identified by a string.

Parameters
nameparameter name
itype index of first atom
jtype index of first atom
valuenew value of parameter

Implemented in DdMd::PairPotentialImpl< Interaction >.

virtual double DdMd::PairPotential::get ( std::string  name,
int  i,
int  j 
) const
pure virtual

Get a parameter value, identified by a string.

Parameters
nameparameter name
itype index of first atom
jtype index of first atom

Implemented in DdMd::PairPotentialImpl< Interaction >.

virtual double DdMd::PairPotential::maxPairCutoff ( ) const
pure virtual

Return maximum cutoff.

Implemented in DdMd::PairPotentialImpl< Interaction >.

Referenced by computeNPair(), initialize(), and readParameters().

virtual std::string DdMd::PairPotential::interactionClassName ( ) const
pure virtual

Return pair interaction class name (e.g., "LJPair").

Implemented in DdMd::PairPotentialImpl< Interaction >.

void DdMd::PairPotential::buildCellList ( )

Build a cell list.

This method adds all local and ghost atoms to a cell list that is used in buildPairList() to build a Verlet pair list. On entry to this method, atomic positions must be generalized. This method should be followed by AtomStorage::transformGenToCart(Boundary& ) to transform all positions back to Cartesian coordinates.

Definition at line 159 of file PairPotential.cpp.

References DdMd::AtomStorage::begin(), DdMd::CellList::build(), cellList_, DdMd::CellList::clear(), cutoff_, Util::Dimension, domain(), DdMd::Domain::domainBound(), DdMd::AtomStorage::isCartesian(), DdMd::CellList::isValid(), Simp::OrthorhombicBoundary::length(), DdMd::CellList::makeGrid(), DdMd::AtomStorage::nAtom(), DdMd::CellList::nAtom(), nCellCut_, DdMd::AtomStorage::nGhost(), Util::PArrayIterator< Data >::notEnd(), DdMd::CellList::nReject(), DdMd::CellList::placeAtom(), storage(), and UTIL_THROW.

Referenced by DdMd::TwoStepIntegrator::run(), and DdMd::Integrator::setupAtoms().

void DdMd::PairPotential::buildPairList ( )

Build the Verlet Pair list.

Preconditions:

Definition at line 206 of file PairPotential.cpp.

References DdMd::PairList::build(), cellList_, pairList_, DdMd::Potential::reverseUpdateFlag(), storage(), and UTIL_THROW.

Referenced by DdMd::TwoStepIntegrator::run(), and DdMd::Integrator::setupAtoms().

virtual void DdMd::PairPotential::computePairEnergies ( MPI::Intracomm &  communicator)
pure virtual

Compute pair energies on all processors.

This method must be called on all processors. The resulting total is stored on the master processor, and may be retrieved by calling the pairEnergies() function on the master processor.

Implemented in DdMd::PairPotentialImpl< Interaction >.

Referenced by DdMd::PairEnergyAnalyzer::compute().

DMatrix< double > DdMd::PairPotential::pairEnergies ( ) const

Return total pair energies, from all processors.

This method should only be called on the master (rank 0) processor, after a previous call to computePairEnergies.

Definition at line 217 of file PairPotential.cpp.

Referenced by DdMd::PairEnergyAnalyzer::value().

void DdMd::PairPotential::unsetPairEnergies ( )

Mark pair energy as unknown (nullify).

Definition at line 229 of file PairPotential.cpp.

Referenced by DdMd::PairPotentialImpl< Interaction >::computePairEnergies().

void DdMd::PairPotential::computeNPair ( MPI::Intracomm &  communicator)

Compute twice the number of pairs within the force cutoff.

This method must be called on all processors.

The method for distributing pairs among processors is the same as is used to distribute the energy, and depends on the value of reverseUpdateFlag().

Definition at line 236 of file PairPotential.cpp.

References maxPairCutoff(), and methodId().

int DdMd::PairPotential::nPair ( ) const
CellList & DdMd::PairPotential::cellList ( )
inline

Get the CellList by reference.

Definition at line 375 of file ddMd/potentials/pair/PairPotential.h.

PairList & DdMd::PairPotential::pairList ( )
inline

Get the PairList by reference.

Definition at line 378 of file ddMd/potentials/pair/PairPotential.h.

Referenced by DdMd::Integrator::clear(), and DdMd::Integrator::outputStatistics().

double DdMd::PairPotential::skin ( ) const
inline

Get value of the pair list skin.

Definition at line 381 of file ddMd/potentials/pair/PairPotential.h.

Referenced by initialize().

double DdMd::PairPotential::cutoff ( ) const
inline

Get value of the pair list cutoff (maxPairCutoff + skin).

Definition at line 384 of file ddMd/potentials/pair/PairPotential.h.

int DdMd::PairPotential::methodId ( ) const
inline

Return integer id for algorithm (0=PAIR, 1=CELL, 2=NSQ)

Definition at line 399 of file ddMd/potentials/pair/PairPotential.h.

Referenced by DdMd::PairPotentialImpl< Interaction >::computeEnergy(), DdMd::PairPotentialImpl< Interaction >::computeForces(), and computeNPair().

Boundary & DdMd::PairPotential::boundary ( )
inlineprotected
Domain & DdMd::PairPotential::domain ( )
inlineprotected

Get the PairList by const reference.

Definition at line 390 of file ddMd/potentials/pair/PairPotential.h.

Referenced by associate(), buildCellList(), and save().

AtomStorage & DdMd::PairPotential::storage ( )
inlineprotected
void DdMd::PairPotential::setPairEnergies ( DMatrix< double >  pairEnergies)
protected

Set values for pair energies.

Definition at line 223 of file PairPotential.cpp.

Referenced by DdMd::PairPotentialImpl< Interaction >::computePairEnergies().

Member Data Documentation

CellList DdMd::PairPotential::cellList_
protected

CellList to construct PairList or calculate nonbonded pair forces.

Definition at line 303 of file ddMd/potentials/pair/PairPotential.h.

Referenced by buildCellList(), buildPairList(), DdMd::PairPotentialImpl< Interaction >::computeEnergy(), nPair(), and save().

PairList DdMd::PairPotential::pairList_
protected
Boundary DdMd::PairPotential::maxBoundary_
protected

Boundary used to allocate space for the cell list.

Definition at line 309 of file ddMd/potentials/pair/PairPotential.h.

Referenced by initialize(), loadParameters(), readParameters(), and save().

double DdMd::PairPotential::skin_
protected

Difference between pairlist cutoff and pair potential cutoff.

Definition at line 312 of file ddMd/potentials/pair/PairPotential.h.

Referenced by initialize(), loadParameters(), readParameters(), and save().

double DdMd::PairPotential::cutoff_
protected

Minimum cell size = pair potential cutoff + skin.

Definition at line 315 of file ddMd/potentials/pair/PairPotential.h.

Referenced by buildCellList(), initialize(), loadParameters(), readParameters(), and save().

int DdMd::PairPotential::nCellCut_
protected

Approximate number of cells per cutoff distance in each direction.

Definition at line 318 of file ddMd/potentials/pair/PairPotential.h.

Referenced by buildCellList(), loadParameters(), readParameters(), and save().

int DdMd::PairPotential::pairCapacity_
protected

Maximum number of nonbonded pairs in pair list.

Definition at line 321 of file ddMd/potentials/pair/PairPotential.h.

Referenced by initialize(), loadParameters(), readParameters(), and save().


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