Simpatico  v1.10
List of all members | Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Static Protected Attributes
Simp::Species Class Reference

Detailed Description

A Species represents a set of chemically similar molecules.

A Species has:

The chemical structure of each molecule within a Species is defined by specifying the type of each atom and constructing a list of bond, angle and dihedral covalent groups. Information about covalent groups is stored internally in a set of SpeciesGroup<N> objects with N=2,3, and 4. This chemical structure is defined within the virtual protected readSpeciesParam(std::istream&) function, which is called by readParam.

Implementations of readSpeciesParam() must read whatever parameters are required to define the chemical structure of some class of molecular structures. The default implementation of Species::readSpeciesParam() is general enough to describe any molecular structure, and simply reads all of the required information about atom types and groups from file. Subclasses of Species generally hard-code some or all of the information required to describe a structure into the class definition, and require correspondingly less information to to be specified in the parameter file. For example, the parameter file format for a class that represents a flexible linear homopolymer might require only the number of monomers per chain, a single atom type index and a single bond type index. A subclass that represents a specific chemical structure (e.g., a water molecule) might not require any data from the parameter file.

A subclass of Species that is defined in the McMd namespace may also represent a "mutable" species. A mutable species is one in which each molecule can be in any of a finite number of different internal states, in which the chemical structure is slightly different in different states. A mutable species could, for example, be used implement a semi-grand ensemble for a polymer mixture, in which each polymer can be of either of two types. A subclass that represents a mutable species must have an associated McMd::SpeciesMutator object. See documentation for the Species::setMutatorPtr() function.

See also
parameter file format

Definition at line 71 of file simp/species/Species.h.

#include <Species.h>

Inheritance diagram for Simp::Species:
Util::ParamComposite Util::ParamComponent Util::Serializable Util::MpiFileIo Simp::Linear Simp::Point Simp::Ring McMd::HomopolymerSG McMd::LinearSG Simp::Diblock Simp::Homopolymer Simp::Multiblock Simp::HomoRing

Public Types

typedef SpeciesGroup< 2 > SpeciesBond
 A SpeciesBond has the local atom ids and a type id for one bond. More...
 
typedef FSArray< int, MaxBondPerAtomAtomBondIdArray
 An array of local integer bond ids for all bonds containing one atom. More...
 
typedef SpeciesGroup< 3 > SpeciesAngle
 A SpeciesAngle has local atom ids and a type id for one angle. More...
 
typedef FSArray< int, MaxAnglePerAtomAtomAngleIdArray
 An array of local angle ids for all angles containing one atom. More...
 
typedef SpeciesGroup< 4 > SpeciesDihedral
 A SpeciesDihedral has local atom ids and a type id for one dihedral. More...
 
typedef FSArray< int, MaxDihedralPerAtomAtomDihedralIdArray
 An array of local angle ids for all dihedrals containing one atom. More...
 
- 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...
 

Public Member Functions

 Species ()
 Constructor. More...
 
virtual ~Species ()
 Destructor. More...
 
Initialization
void setId (int id)
 Set integer id for this Species. More...
 
virtual void readParameters (std::istream &in)
 Read parameters and initialize structure for this species. More...
 
virtual void loadParameters (Serializable::IArchive &ar)
 Load internal state from an archive. More...
 
void readStructure (std::istream &in)
 Read structure from config/topology file format. More...
 
Output to File
virtual void save (Serializable::OArchive &ar)
 Save internal state to an archive. More...
 
void writeStructure (std::ostream &out, std::string indent=std::string())
 Write molecular structure in config/topology file format. More...
 
bool matchStructure (std::istream &in)
 Read structure, return true iff it matches existing structure. More...
 
Chemical Structure Accessors
int nAtom () const
 Get number of atoms per molecule for this Species. More...
 
int atomTypeId (int iAtom) const
 Get atom type index for a specific atom, by local atom index. More...
 
int nBond () const
 Get number of bonds per molecule for this Species. More...
 
const SpeciesBondspeciesBond (int iBond) const
 Get a specific SpeciesBond object, by local bond index. More...
 
const AtomBondIdArrayatomBondIds (int atomId) const
 Get array of ids for Bonds that contain one Atom. More...
 
int nAngle () const
 Get number of angles per molecule for this Species. More...
 
const SpeciesAnglespeciesAngle (int iAngle) const
 Get a specific SpeciesAngle object, by local angle index. More...
 
const AtomAngleIdArrayatomAngleIds (int atomId) const
 Get array of ids for angles that contain one Atom. More...
 
int nDihedral () const
 Get number of dihedrals per molecule for this Species. More...
 
const SpeciesDihedralspeciesDihedral (int iDihedral) const
 Get a specific SpeciesDihedral object, by local angle index. More...
 
const AtomDihedralIdArrayatomDihedralIds (int atomId) const
 Get array of ids for dihedrals that contain one Atom. More...
 
Interface for mutable species.
bool isMutable () const
 Is this a mutable Species? More...
 
McMd::SpeciesMutatormutator ()
 Return the species mutator object by reference. More...
 
Miscellaneous Accessors
int id () const
 Get integer id of this Species. More...
 
int capacity () const
 Maximum allowed number of molecules for this Species. More...
 
bool isValid () const
 Return true if Species is valid, or throw an Exception. 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...
 

Static Public Attributes

static const int MaxBondPerAtom = 4
 Maximum number of bonds that can be connected to one atom. More...
 
static const int MaxAnglePerAtom = 18
 Maximum number of angles groups that can contain one atom. More...
 
static const int MaxDihedralPerAtom = 72
 Maximum number of dihedral groups that can contain one atom. More...
 

Protected Member Functions

virtual void readSpeciesParam (std::istream &in)
 Define chemical structure for this Species. More...
 
virtual void loadSpeciesParam (Serializable::IArchive &ar)
 Define chemical structure for this Species. More...
 
void allocate ()
 Allocate chemical structure arrays. More...
 
void allocateAtoms ()
 Allocate and initialize array of atom type Ids. More...
 
void setAtomType (int atomId, int atomType)
 Set the type for one atom in a generic molecule of this Species. More...
 
void allocateBonds ()
 Allocate arrays associated with Bonds. More...
 
void makeBond (int bondId, int atomId1, int atomId2, int bondType)
 Add a bond to the chemical structure of a generic molecule. More...
 
void allocateAngles ()
 Allocate arrays associated with angles. More...
 
void makeAngle (int angleId, int atomId1, int atomId2, int atomId3, int angleType)
 Add an angle to the chemical structure of a generic molecule. More...
 
void allocateDihedrals ()
 Allocate arrays associated with dihedrals. More...
 
void makeDihedral (int dihedralId, int atomId1, int atomId2, int atomId3, int atomId4, int dihedralType)
 Add a dihedral to the chemical structure of a generic molecule. More...
 
void initializeAtomGroupIdArrays ()
 Initialize all atom groupId arrays (point from atoms to groups). More...
 
void setMutatorPtr (McMd::SpeciesMutator *mutatorPtr)
 Set pointer to associated McMd::SpeciesMutator for a mutable species. 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

int id_
 Integer index for this Species. More...
 
int moleculeCapacity_
 Number of molecules associated with the species. More...
 
int nAtom_
 Number of atoms per molecule. More...
 
DArray< int > atomTypeIds_
 Array of atom type Ids, indexed by local atom id. More...
 
int nBond_
 Number of bonds per molecule. More...
 
DArray< SpeciesBondspeciesBonds_
 Array of SpeciesBonds for all bonds, indexed by local bond id. More...
 
int nAngle_
 Number of angles per molecule. More...
 
DArray< SpeciesAnglespeciesAngles_
 Array of SpeciesAngles for all angles, indexed by local angle id. More...
 
int nDihedral_
 Number of dihedrals per molecule. More...
 
DArray< SpeciesDihedralspeciesDihedrals_
 Array of SpeciesDihedrals, indexed by local dihedral id. More...
 

Static Protected Attributes

static const int NullIndex = -1
 Null (unknown) value for any non-negative index. More...
 

Additional Inherited Members

- 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...
 

Member Typedef Documentation

A SpeciesBond has the local atom ids and a type id for one bond.

Definition at line 81 of file simp/species/Species.h.

An array of local integer bond ids for all bonds containing one atom.

Definition at line 84 of file simp/species/Species.h.

A SpeciesAngle has local atom ids and a type id for one angle.

Definition at line 92 of file simp/species/Species.h.

An array of local angle ids for all angles containing one atom.

Definition at line 95 of file simp/species/Species.h.

A SpeciesDihedral has local atom ids and a type id for one dihedral.

Definition at line 103 of file simp/species/Species.h.

An array of local angle ids for all dihedrals containing one atom.

Definition at line 106 of file simp/species/Species.h.

Constructor & Destructor Documentation

Tools::Species::Species ( )

Constructor.

Definition at line 22 of file simp/species/Species.cpp.

References Util::ParamComposite::setClassName().

Tools::Species::~Species ( )
virtual

Destructor.

Definition at line 45 of file simp/species/Species.cpp.

Member Function Documentation

void Tools::Species::setId ( int  id)

Set integer id for this Species.

Parameters
idinteger id

Definition at line 452 of file simp/species/Species.cpp.

References id(), id_, and UTIL_THROW.

Referenced by McMd::Simulation::loadParameters(), and McMd::Simulation::readParameters().

void Simp::Species::readParameters ( std::istream &  in)
virtual

Read parameters and initialize structure for this species.

This function reads the parameter moleculeCapacity (the maximum allowed number of molecules of this Species) and invokes the virtual Species::readSpeciesParam() function to define the chemical structure of the Species.

Parameters
ininput stream.

Reimplemented from Util::ParamComposite.

Reimplemented in McMd::LinearSG, and McMd::HomopolymerSG.

Definition at line 51 of file simp/species/Species.cpp.

References isValid(), moleculeCapacity_, and readSpeciesParam().

Referenced by McMd::HomopolymerSG::readParameters(), and McMd::LinearSG::readParameters().

void Simp::Species::loadParameters ( Serializable::IArchive ar)
virtual

Load internal state from an archive.

Parameters
arinput/loading archive

Reimplemented from Util::ParamComposite.

Definition at line 65 of file simp/species/Species.cpp.

References isValid(), loadSpeciesParam(), and moleculeCapacity_.

void Simp::Species::readStructure ( std::istream &  in)

Read structure from config/topology file format.

Parameters
infile from which to read structure.

Definition at line 160 of file simp/species/Species.cpp.

References allocateAngles(), allocateAtoms(), allocateBonds(), allocateDihedrals(), atomTypeIds_, initializeAtomGroupIdArrays(), Util::Label::isClear(), isValid(), nAngle_, nAtom_, nBond_, nDihedral_, speciesAngles_, speciesBonds_, speciesDihedrals_, and UTIL_CHECK.

void Simp::Species::save ( Serializable::OArchive ar)
virtual

Save internal state to an archive.

Parameters
aroutput/saving archive

Reimplemented from Util::ParamComposite.

Reimplemented in McMd::LinearSG, McMd::HomopolymerSG, Simp::Homopolymer, Simp::HomoRing, Simp::Diblock, Simp::Multiblock, and Simp::Point.

Definition at line 312 of file simp/species/Species.cpp.

References atomTypeIds_, moleculeCapacity_, nAngle_, nAtom_, nBond_, nDihedral_, speciesAngles_, speciesBonds_, and speciesDihedrals_.

void Simp::Species::writeStructure ( std::ostream &  out,
std::string  indent = std::string() 
)

Write molecular structure in config/topology file format.

Parameters
outfile to which to write structure.
indentindentation string (sequence of spaces)

Definition at line 340 of file simp/species/Species.cpp.

References atomTypeIds_, Util::ParamComponent::indent(), nAngle_, nAtom_, nBond_, nDihedral_, speciesAngles_, speciesBonds_, and speciesDihedrals_.

Referenced by McMd::SmpConfigIo::write().

bool Simp::Species::matchStructure ( std::istream &  in)

Read structure, return true iff it matches existing structure.

Throw exception if structure is not parseable. Return false if it is parseable but does not match.

Parameters
infile from which to read structure.

Definition at line 387 of file simp/species/Species.cpp.

References atomTypeIds_, nAngle(), nAngle_, nAtom(), nAtom_, nBond(), nBond_, nDihedral(), nDihedral_, speciesAngles_, speciesBonds_, and speciesDihedrals_.

Referenced by McMd::SmpConfigIo::read().

int Tools::Species::nAtom ( ) const
inline

Get number of atoms per molecule for this Species.

Definition at line 591 of file simp/species/Species.h.

References nAtom_.

Referenced by McMd::TrajectoryReader::addMolecules(), McMd::LinearGenerator::attemptPlaceMolecule(), Simp::Linear::buildLinear(), Simp::Ring::buildRing(), Tools::Species::initialize(), McMd::Simulation::isValid(), McMd::AtomMSD::loadParameters(), McMd::RigidDisplaceMove::loadParameters(), McMd::CfbLinearEndMove::loadParameters(), McMd::EndSwapMove::loadParameters(), McMd::RingTetraRebridgeMove::loadParameters(), McMd::RingRouseAutoCorr::loadParameters(), McMd::LinearRouseAutoCorr::loadParameters(), McMd::McMuExchange::loadParameters(), McMd::IntraPairAutoCorr::loadParameters(), McMd::RadiusGyration::loadParameters(), McMd::BlockRadiusGyration::loadParameters(), matchStructure(), Tools::Species::molecule(), McMd::SystemInterface::nAtom(), McMd::System::nAtom(), McMd::LinkLTPos::output(), McMd::SmpConfigIo::read(), McMd::DdMdConfigIo::read(), McMd::LammpsConfigIo::read(), McMd::RigidDisplaceMove::readParameters(), McMd::AtomMSD::readParameters(), McMd::RingTetraRebridgeMove::readParameters(), McMd::CfbLinearEndMove::readParameters(), McMd::G1MSD::readParameters(), McMd::CfbReptationMove::readParameters(), McMd::EndSwapMove::readParameters(), Tools::AtomMSD::readParameters(), McMd::ComMSD::readParameters(), McMd::RingRouseAutoCorr::readParameters(), McMd::McMuExchange::readParameters(), McMd::LinearRouseAutoCorr::readParameters(), McMd::EndtoEndXYZ::readParameters(), McMd::IntraPairAutoCorr::readParameters(), McMd::RadiusGyration::readParameters(), McMd::EndtoEnd::readParameters(), McMd::BlockRadiusGyration::readParameters(), McMd::LinkLTPos::sample(), McMd::Simulation::save(), Tools::ConfigReader::setAtomContexts(), McMd::HomopolymerSG::setMoleculeState(), McMd::LinearSG::setMoleculeState(), McMd::LinkLTPos::setup(), McMd::RingRouseAutoCorr::setup(), McMd::SmpConfigIo::write(), McMd::DdMdConfigIo::write(), and McMd::LammpsConfigIo::write().

int Simp::Species::atomTypeId ( int  iAtom) const
inline

Get atom type index for a specific atom, by local atom index.

Parameters
iAtomlocal index for an atom within a molecule.

Definition at line 597 of file simp/species/Species.h.

References atomTypeIds_.

Referenced by McMd::CfbLinear::addAtom(), McMd::Simulation::isValid(), McMd::McMuExchange::loadParameters(), McMd::SmpConfigIo::read(), McMd::CfbReptationMove::readParameters(), McMd::McMuExchange::readParameters(), and McMd::Simulation::save().

int Simp::Species::nBond ( ) const
inline
const Species::SpeciesBond & Simp::Species::speciesBond ( int  iBond) const
inline

Get a specific SpeciesBond object, by local bond index.

Parameters
iBondlocal index of a SpeciesBond within a molecule.

Definition at line 611 of file simp/species/Species.h.

References speciesBonds_.

Referenced by initializeAtomGroupIdArrays(), isValid(), McMd::Simulation::isValid(), McMd::CfbReptationMove::loadParameters(), McMd::CfbReptateMove::loadParameters(), McMd::CfbReptationMove::readParameters(), and McMd::Simulation::save().

const Species::AtomBondIdArray & Simp::Species::atomBondIds ( int  atomId) const
inline

Get array of ids for Bonds that contain one Atom.

Parameters
atomIdlocal index for atom of interest

Definition at line 618 of file simp/species/Species.h.

Referenced by McMd::Activate::deactivate(), McMd::getAtomBonds(), and McMd::Activate::reactivate().

int Simp::Species::nAngle ( ) const
inline
const Species::SpeciesAngle & Simp::Species::speciesAngle ( int  iAngle) const
inline

Get a specific SpeciesAngle object, by local angle index.

Parameters
iAnglelocal index of a SpeciesAngle within a molecule.

Definition at line 633 of file simp/species/Species.h.

References speciesAngles_.

Referenced by initializeAtomGroupIdArrays(), isValid(), McMd::Simulation::isValid(), and McMd::Simulation::save().

const Species::AtomAngleIdArray & Simp::Species::atomAngleIds ( int  atomId) const
inline

Get array of ids for angles that contain one Atom.

Parameters
atomIdlocal index for atom of interest

Definition at line 640 of file simp/species/Species.h.

Referenced by McMd::Activate::deactivate(), McMd::getAtomAngles(), and McMd::Activate::reactivate().

int Simp::Species::nDihedral ( ) const
inline
const Species::SpeciesDihedral & Simp::Species::speciesDihedral ( int  iDihedral) const
inline

Get a specific SpeciesDihedral object, by local angle index.

Parameters
iDihedrallocal index of a SpeciesDihedral within a molecule.

Definition at line 655 of file simp/species/Species.h.

References speciesDihedrals_.

Referenced by initializeAtomGroupIdArrays(), isValid(), McMd::Simulation::isValid(), and McMd::Simulation::save().

const Species::AtomDihedralIdArray & Simp::Species::atomDihedralIds ( int  atomId) const
inline

Get array of ids for dihedrals that contain one Atom.

Parameters
atomIdlocal index for atom of interest

Definition at line 663 of file simp/species/Species.h.

Referenced by McMd::Activate::deactivate(), McMd::getAtomDihedrals(), and McMd::Activate::reactivate().

bool Simp::Species::isMutable ( ) const
inline
McMd::SpeciesMutator & Simp::Species::mutator ( )
inline
int Tools::Species::id ( ) const
inline
int Tools::Species::capacity ( ) const
inline

Maximum allowed number of molecules for this Species.

Definition at line 585 of file simp/species/Species.h.

References moleculeCapacity_.

Referenced by McMd::TrajectoryReader::addMolecules(), McMd::Simulation::allocateMoleculeSet(), Tools::Species::initialize(), McMd::ClusterIdentifier::initialize(), McMd::Simulation::isValid(), McMd::SemiGrandDistribution::loadParameters(), McMd::AtomMSD::loadParameters(), McMd::ComMSD::loadParameters(), McMd::RingRouseAutoCorr::loadParameters(), McMd::IntraBondStressAutoCorr< SystemType >::loadParameters(), McMd::LinearRouseAutoCorr::loadParameters(), McMd::MdPairEnergyCoefficients::loadParameters(), McMd::IntraBondTensorAutoCorr< SystemType >::loadParameters(), McMd::McMuExchange::loadParameters(), McMd::IntraPairAutoCorr::loadParameters(), McMd::HomopolymerSG::loadSpeciesParam(), McMd::LinearSG::loadSpeciesParam(), McMd::SmpConfigIo::read(), McMd::DdMdConfigIo::read(), McMd::LammpsConfigIo::read(), McMd::LammpsDumpReader::readFrame(), McMd::DdMdTrajectoryReader::readFrame(), McMd::DCDTrajectoryReader::readFrame(), McMd::AtomMSD::readParameters(), McMd::CfbReptationMove::readParameters(), McMd::CfbReptateMove::readParameters(), Tools::AtomMSD::readParameters(), McMd::ComMSD::readParameters(), McMd::McMuExchange::readParameters(), McMd::RingRouseAutoCorr::readParameters(), McMd::IntraBondStressAutoCorr< SystemType >::readParameters(), McMd::LinearRouseAutoCorr::readParameters(), McMd::MdPairEnergyCoefficients::readParameters(), McMd::IntraBondTensorAutoCorr< SystemType >::readParameters(), McMd::IntraPairAutoCorr::readParameters(), McMd::HomopolymerSG::readSpeciesParam(), McMd::LinearSG::readSpeciesParam(), McMd::MdPairEnergyCoefficients::sample(), McMd::Simulation::save(), Tools::ConfigReader::setAtomContexts(), Tools::Species::size(), and McMd::MdPairEnergyCoefficients::~MdPairEnergyCoefficients().

bool Tools::Species::isValid ( ) const
void Simp::Species::readSpeciesParam ( std::istream &  in)
protectedvirtual

Define chemical structure for this Species.

This virtual function must define the structure of a molecule of this species, and read any data required to do so. The default implementation Species::readSpeciesParam reads nAtom_, nBond_, nAngle_, and the elements of the arrays atomTypeIds_ and speciesBonds_, speciesAngles_.

Re-implementations of this function by subclasses may hard-code some or all of the information contained in these variables, and may define a more compact file format for any parameters that are read from file.

Parameters
ininput parameter file stream.

Reimplemented in Simp::Homopolymer, Simp::Multiblock, McMd::LinearSG, McMd::HomopolymerSG, Simp::Diblock, Simp::HomoRing, and Simp::Point.

Definition at line 79 of file simp/species/Species.cpp.

References allocateAngles(), allocateAtoms(), allocateBonds(), allocateDihedrals(), atomTypeIds_, initializeAtomGroupIdArrays(), nAngle_, nAtom_, nBond_, nDihedral_, speciesAngles_, speciesBonds_, and speciesDihedrals_.

Referenced by readParameters().

void Simp::Species::loadSpeciesParam ( Serializable::IArchive ar)
protectedvirtual

Define chemical structure for this Species.

Analogous to readSpeciesParam, but reads data from a Serializable::IArchive rather than an input stream. This function must define the same parameter file format as readSpeciesParam.

Parameters
arinput parameter file stream.

Reimplemented in Simp::Homopolymer, Simp::Multiblock, McMd::LinearSG, McMd::HomopolymerSG, Simp::Diblock, Simp::HomoRing, and Simp::Point.

Definition at line 121 of file simp/species/Species.cpp.

References allocateAngles(), allocateAtoms(), allocateBonds(), allocateDihedrals(), atomTypeIds_, initializeAtomGroupIdArrays(), nAngle_, nAtom_, nBond_, nDihedral_, speciesAngles_, speciesBonds_, and speciesDihedrals_.

Referenced by loadParameters().

void Simp::Species::allocate ( )
protected

Allocate chemical structure arrays.

This function allocates the arrays that are used to define the chemical structure of a generic molecule, such as atomTypeIds_, speciesBonds_, atomBondIdArrays_, speciesAngles_, etc.

Precondition: nAtom_, nBond_, nAngles_, etc. must have been assigned final values on entry.

Definition at line 470 of file simp/species/Species.cpp.

References allocateAngles(), allocateAtoms(), allocateBonds(), and allocateDihedrals().

Referenced by Simp::Linear::buildLinear(), Simp::Ring::buildRing(), Simp::Point::loadSpeciesParam(), and Simp::Point::readSpeciesParam().

void Simp::Species::allocateAtoms ( )
protected

Allocate and initialize array of atom type Ids.

Precondition
nAtom_ must be assigned a positive value

Definition at line 487 of file simp/species/Species.cpp.

References Util::DArray< Data >::allocate(), atomTypeIds_, Util::DArray< Data >::isAllocated(), nAtom_, and UTIL_CHECK.

Referenced by allocate(), loadSpeciesParam(), readSpeciesParam(), and readStructure().

void Simp::Species::setAtomType ( int  atomId,
int  atomType 
)
protected

Set the type for one atom in a generic molecule of this Species.

Parameters
atomIdlocal atom id, in range 0,..., nAtom - 1;
atomTypeatom type index

Definition at line 503 of file simp/species/Species.cpp.

References atomTypeIds_.

Referenced by Simp::Linear::buildLinear(), and Simp::Ring::buildRing().

void Simp::Species::allocateBonds ( )
protected

Allocate arrays associated with Bonds.

Precondition
nAtom_ > 0 and nBond_ >= 0.

Definition at line 510 of file simp/species/Species.cpp.

References nAtom_, nBond_, speciesBonds_, and UTIL_CHECK.

Referenced by allocate(), loadSpeciesParam(), readSpeciesParam(), and readStructure().

void Simp::Species::makeBond ( int  bondId,
int  atomId1,
int  atomId2,
int  bondType 
)
protected

Add a bond to the chemical structure of a generic molecule.

This function creates and adds a SpeciesBond object, and also adds a reference to the list of bonds that are connected to each atom.

Parameters
bondIdlocal index of bond within a molecule
atomId1local index of 1st atom in bond
atomId2local index of 2nd atom in bond
bondTypebond type index

Definition at line 526 of file simp/species/Species.cpp.

References nAtom_, nBond_, and speciesBonds_.

Referenced by Simp::Linear::buildLinear(), and Simp::Ring::buildRing().

void Simp::Species::allocateAngles ( )
protected

Allocate arrays associated with angles.

Precondition
nAtom_ > 0 and nAngle_ >= 0.

Definition at line 553 of file simp/species/Species.cpp.

References nAngle_, nAtom_, speciesAngles_, and UTIL_CHECK.

Referenced by allocate(), loadSpeciesParam(), readSpeciesParam(), and readStructure().

void Simp::Species::makeAngle ( int  angleId,
int  atomId1,
int  atomId2,
int  atomId3,
int  angleType 
)
protected

Add an angle to the chemical structure of a generic molecule.

This function creates and adds a SpeciesAngle object, and also adds a reference to the list of angles that are connected to each atom.

Parameters
angleIdlocal index of angle within a molecule
atomId1local index of 1st atom in angle
atomId2local index of 2nd atom in angle
atomId3local index of 3rd atom in angle
angleTypeangle type index

Definition at line 569 of file simp/species/Species.cpp.

References nAngle_, nAtom_, and speciesAngles_.

Referenced by Simp::Linear::buildLinear(), and Simp::Ring::buildRing().

void Simp::Species::allocateDihedrals ( )
protected

Allocate arrays associated with dihedrals.

Definition at line 599 of file simp/species/Species.cpp.

References nAtom_, nDihedral_, speciesDihedrals_, and UTIL_CHECK.

Referenced by allocate(), loadSpeciesParam(), readSpeciesParam(), and readStructure().

void Simp::Species::makeDihedral ( int  dihedralId,
int  atomId1,
int  atomId2,
int  atomId3,
int  atomId4,
int  dihedralType 
)
protected

Add a dihedral to the chemical structure of a generic molecule.

This function creates and adds a SpeciesDihedral object, and also adds a reference to the list of dihedrals that are connected to each atom.

Parameters
dihedralIdlocal index of dihedral within a molecule
atomId1local index of 1st atom in dihedral
atomId2local index of 2nd atom in dihedral
atomId3local index of 3rd atom in dihedral
atomId4local index of 4th atom in dihedral
dihedralTypedihedral type index

Definition at line 616 of file simp/species/Species.cpp.

References nAtom_, nDihedral_, and speciesDihedrals_.

Referenced by Simp::Linear::buildLinear(), and Simp::Ring::buildRing().

void Simp::Species::initializeAtomGroupIdArrays ( )
protected

Initialize all atom groupId arrays (point from atoms to groups).

Precondition
Species::allocate() function must have been invoked.
All speciesGroup arrays (speciesBonds, etc.) must be initialized.

Definition at line 229 of file simp/species/Species.cpp.

References Simp::SpeciesGroup< NAtom >::atomId(), atomTypeIds_, Util::Array< Data >::capacity(), Util::DArray< Data >::isAllocated(), nAngle_, nAtom_, nBond_, nDihedral_, speciesAngle(), speciesAngles_, speciesBond(), speciesBonds_, speciesDihedral(), speciesDihedrals_, and UTIL_CHECK.

Referenced by loadSpeciesParam(), readSpeciesParam(), and readStructure().

void Simp::Species::setMutatorPtr ( McMd::SpeciesMutator mutatorPtr)
protected

Set pointer to associated McMd::SpeciesMutator for a mutable species.

A mutable subclass of Species must have an associated SpeciesMutator object. The constructor of each such subclass should pass a pointer to the SpeciesMutator to setMutatorPtr. After this function is called, Species::isMutator() will return true, the associated SpeciesMutator will be accessible via the Species::mutator() function. If a mutable Species subclass is derived from both Species and SpeciesMutator (recommended), the function setMutatorPtr(this) should be invoked within the subclass constructor to pass the address of the SpeciesMutator sub-object to the Species subobject.

Parameters
mutatorPtrpointer to an associated SpeciesMutator object

Definition at line 464 of file simp/species/Species.cpp.

Referenced by McMd::HomopolymerSG::HomopolymerSG(), and McMd::LinearSG::LinearSG().

Member Data Documentation

const int Simp::Species::MaxBondPerAtom = 4
static

Maximum number of bonds that can be connected to one atom.

Definition at line 78 of file simp/species/Species.h.

const int Simp::Species::MaxAnglePerAtom = 18
static

Maximum number of angles groups that can contain one atom.

Definition at line 89 of file simp/species/Species.h.

const int Simp::Species::MaxDihedralPerAtom = 72
static

Maximum number of dihedral groups that can contain one atom.

Definition at line 100 of file simp/species/Species.h.

const int Simp::Species::NullIndex = -1
staticprotected

Null (unknown) value for any non-negative index.

Definition at line 304 of file simp/species/Species.h.

Referenced by McMd::HomopolymerSG::calculateAtomTypeId(), McMd::LinearSG::calculateAtomTypeId(), and Simp::Diblock::Diblock().

int Simp::Species::id_
protected
int Simp::Species::moleculeCapacity_
protected
int Simp::Species::nAtom_
protected
DArray<int> Simp::Species::atomTypeIds_
protected

Array of atom type Ids, indexed by local atom id.

Element atomTypeIds_[id] is the atom type index for atom id of any molecule of this species, where 0 <= id < nAtom_.

Definition at line 328 of file simp/species/Species.h.

Referenced by allocateAtoms(), atomTypeId(), initializeAtomGroupIdArrays(), isValid(), Simp::Point::loadSpeciesParam(), loadSpeciesParam(), matchStructure(), Simp::Point::readSpeciesParam(), readSpeciesParam(), readStructure(), save(), setAtomType(), and writeStructure().

int Simp::Species::nBond_
protected
DArray<SpeciesBond> Simp::Species::speciesBonds_
protected

Array of SpeciesBonds for all bonds, indexed by local bond id.

Element speciesBonds_[id] is the SpeciesBond object for bond number id in any molecule of this Species, where 0 <= id < nBond_.

Definition at line 342 of file simp/species/Species.h.

Referenced by allocateBonds(), initializeAtomGroupIdArrays(), isValid(), loadSpeciesParam(), makeBond(), matchStructure(), readSpeciesParam(), readStructure(), save(), speciesBond(), and writeStructure().

int Simp::Species::nAngle_
protected
DArray<SpeciesAngle> Simp::Species::speciesAngles_
protected

Array of SpeciesAngles for all angles, indexed by local angle id.

Element speciesAngles_[id] is the SpeciesAngle object for angle number id in any molecule of this Species, where 0 <= id < nAngle_.

Definition at line 357 of file simp/species/Species.h.

Referenced by allocateAngles(), initializeAtomGroupIdArrays(), isValid(), loadSpeciesParam(), makeAngle(), matchStructure(), readSpeciesParam(), readStructure(), save(), speciesAngle(), and writeStructure().

int Simp::Species::nDihedral_
protected
DArray<SpeciesDihedral> Simp::Species::speciesDihedrals_
protected

Array of SpeciesDihedrals, indexed by local dihedral id.

Element speciesDihedrals_[id] is the SpeciesDihedral object for dihedral number id in any molecule of this Species, where the index satisfies 0 <= id < nAngle_.

Definition at line 373 of file simp/species/Species.h.

Referenced by allocateDihedrals(), initializeAtomGroupIdArrays(), isValid(), loadSpeciesParam(), makeDihedral(), matchStructure(), readSpeciesParam(), readStructure(), save(), speciesDihedral(), and writeStructure().


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