OpenMD
2.6
Molecular Dynamics in the Open
|
Namespaces | |
Constants | |
Functions | |
template<class Cont , class Predict > | |
void | swap_if (Cont &b1, Cont &b2, Predict predict) |
std::ostream & | operator<< (std::ostream &o, HydrodynamicsModelFactory &factory) |
std::ostream & | operator<< (std::ostream &o, PairList &e) |
void | registerIntegrators () |
void | registerOptimizers () |
void | registerLattice () |
void | registerAll () |
int | getGlobalCountOfType (AtomType *atype) |
ostream & | operator<< (ostream &o, SimInfo &info) |
bool | pairCompare (const pair< RealType, int > &l, const pair< RealType, int > &r) |
std::ostream & | operator<< (std::ostream &o, IntegratorFactory &factory) |
NotEmptyConstraint | isNotEmpty () |
ZeroConstraint | isZero () |
ParamConstraintFacade< NonZeroConstraint > | isNonZero () |
PositiveConstraint | isPositive () |
NonPositiveConstraint | isNonPositive () |
NegativeConstraint | isNegative () |
NonNegativeConstraint | isNonNegative () |
EqualIgnoreCaseConstraint | isEqualIgnoreCase (std::string str) |
EvenConstraint | isEven () |
template<typename Cons1T , typename Cons2T > | |
AndParamConstraint< Cons1T, Cons2T > | operator && (const ParamConstraintFacade< Cons1T > &cons1, const ParamConstraintFacade< Cons2T > &cons2) |
template<typename Cons1T , typename Cons2T > | |
OrParamConstraint< Cons1T, Cons2T > | operator|| (const ParamConstraintFacade< Cons1T > &cons1, const ParamConstraintFacade< Cons2T > &cons2) |
template<typename ConsT > | |
NotParamConstraint< ConsT > | operator ! (const ParamConstraintFacade< ConsT > &cons) |
template<typename T > | |
LessThanConstraint< T > | isLessThan (T &v) |
template<typename T > | |
LessThanOrEqualToConstraint< T > | isLessThanOrEqualTo (T &v) |
template<typename T > | |
EqualConstraint< T > | isEqual (T &v) |
template<typename T > | |
GreaterThanConstraint< T > | isGreaterThan (T &v) |
template<typename T > | |
GreaterThanOrEqualTo< T > | isGreaterThanOrEqualTo (T &v) |
std::ostream & | operator<< (std::ostream &o, LatticeFactory &factory) |
template<class MatrixType > | |
void | CholeskyDecomposition (MatrixType &A, MatrixType &L) |
template<typename Real > | |
DynamicRectMatrix< Real > | operator - (const DynamicRectMatrix< Real > &m) |
template<typename Real > | |
DynamicRectMatrix< Real > | operator+ (const DynamicRectMatrix< Real > &m1, const DynamicRectMatrix< Real > &m2) |
template<typename Real > | |
DynamicRectMatrix< Real > | operator - (const DynamicRectMatrix< Real > &m1, const DynamicRectMatrix< Real > &m2) |
template<typename Real > | |
DynamicRectMatrix< Real > | operator * (const DynamicRectMatrix< Real > &m, Real s) |
template<typename Real > | |
DynamicRectMatrix< Real > | operator * (Real s, const DynamicRectMatrix< Real > &m) |
template<typename Real > | |
DynamicRectMatrix< Real > | operator * (const DynamicRectMatrix< Real > &m1, const DynamicRectMatrix< Real > &m2) |
template<typename Real > | |
DynamicVector< Real > | operator * (const DynamicRectMatrix< Real > &m, const DynamicVector< Real > &v) |
template<typename Real > | |
DynamicRectMatrix< Real > | operator/ (const DynamicRectMatrix< Real > &m, Real s) |
template<typename Real > | |
std::ostream & | operator<< (std::ostream &o, const DynamicRectMatrix< Real > &m) |
template<typename Real > | |
DynamicVector< Real > | operator - (const DynamicVector< Real > &v1) |
template<typename Real > | |
DynamicVector< Real > | operator+ (const DynamicVector< Real > &v1, const DynamicVector< Real > &v2) |
template<typename Real > | |
DynamicVector< Real > | operator - (const DynamicVector< Real > &v1, const DynamicVector< Real > &v2) |
template<typename Real > | |
DynamicVector< Real > | operator * (const DynamicVector< Real > &v1, Real s) |
template<typename Real > | |
DynamicVector< Real > | operator * (Real s, const DynamicVector< Real > &v1) |
template<typename Real > | |
DynamicVector< Real > | operator/ (const DynamicVector< Real > &v1, Real s) |
template<typename Real > | |
Real | dot (const DynamicVector< Real > &v1, const DynamicVector< Real > &v2) |
template<typename Real > | |
Real | distance (const DynamicVector< Real > &v1, const DynamicVector< Real > &v2) |
template<typename Real > | |
Real | distanceSquare (const DynamicVector< Real > &v1, const DynamicVector< Real > &v2) |
template<typename Real > | |
std::ostream & | operator<< (std::ostream &o, const DynamicVector< Real > &v) |
template<class MatrixType > | |
bool | invertMatrix (MatrixType &A, MatrixType &AI) |
template<class MatrixType > | |
bool | invertMatrix (MatrixType &A, MatrixType &AI, unsigned int size, int *tmp1Size, typename MatrixType::ElemPoinerType tmp2Size) |
template<class MatrixType > | |
int | LUFactorLinearSystem (MatrixType &A, int *index, int size, typename MatrixType::ElemPoinerType tmpSize) |
template<class MatrixType > | |
void | LUSolveLinearSystem (MatrixType &A, int *index, typename MatrixType::ElemPoinerType x, int size) |
std::ostream & | operator<< (std::ostream &os, const MTRand &mtrand) |
std::istream & | operator>> (std::istream &is, MTRand &mtrand) |
template<typename Real > | |
Real | fastpow (Real x, int N) |
template<typename Real > | |
Polynomial< Real > | operator * (const Polynomial< Real > &p1, const Polynomial< Real > &p2) |
template<typename Real > | |
Polynomial< Real > | operator * (const Polynomial< Real > &p, const Real v) |
template<typename Real > | |
Polynomial< Real > | operator * (const Real v, const Polynomial< Real > &p) |
template<typename Real > | |
Polynomial< Real > | operator+ (const Polynomial< Real > &p1, const Polynomial< Real > &p2) |
template<typename Real > | |
Polynomial< Real > | operator - (const Polynomial< Real > &p1, const Polynomial< Real > &p2) |
template<typename Real > | |
Polynomial< Real > * | getDerivative (const Polynomial< Real > &p1) |
template<typename Real > | |
bool | equal (const Polynomial< Real > &p1, const Polynomial< Real > &p2) |
template<typename Real , unsigned int Dim> | |
Quaternion< Real > | operator * (const Quaternion< Real > &q, Real s) |
template<typename Real , unsigned int Dim> | |
Quaternion< Real > | operator * (const Real &s, const Quaternion< Real > &q) |
template<typename Real > | |
Quaternion< Real > | operator * (const Quaternion< Real > &q1, const Quaternion< Real > &q2) |
template<typename Real > | |
Quaternion< Real > | operator/ (Quaternion< Real > &q1, Quaternion< Real > &q2) |
template<typename Real > | |
Quaternion< Real > | operator/ (const Real &s, Quaternion< Real > &q) |
template<class T > | |
bool | operator== (const Quaternion< T > &lhs, const Quaternion< T > &rhs) |
template<typename Real , unsigned int Row, unsigned int Col> | |
RectMatrix< Real, Row, Col > | operator - (const RectMatrix< Real, Row, Col > &m) |
template<typename Real , unsigned int Row, unsigned int Col> | |
RectMatrix< Real, Row, Col > | operator+ (const RectMatrix< Real, Row, Col > &m1, const RectMatrix< Real, Row, Col > &m2) |
template<typename Real , unsigned int Row, unsigned int Col> | |
RectMatrix< Real, Row, Col > | operator - (const RectMatrix< Real, Row, Col > &m1, const RectMatrix< Real, Row, Col > &m2) |
template<typename Real , unsigned int Row, unsigned int Col> | |
RectMatrix< Real, Row, Col > | operator * (const RectMatrix< Real, Row, Col > &m, Real s) |
template<typename Real , unsigned int Row, unsigned int Col> | |
RectMatrix< Real, Row, Col > | operator * (Real s, const RectMatrix< Real, Row, Col > &m) |
template<typename Real , unsigned int Row, unsigned int Col, unsigned int SameDim> | |
RectMatrix< Real, Row, Col > | operator * (const RectMatrix< Real, Row, SameDim > &m1, const RectMatrix< Real, SameDim, Col > &m2) |
template<typename Real , unsigned int Row, unsigned int Col> | |
Vector< Real, Row > | operator * (const RectMatrix< Real, Row, Col > &m, const Vector< Real, Col > &v) |
template<typename Real , unsigned int Row, unsigned int Col> | |
Vector< Real, Col > | operator * (const Vector< Real, Row > &v, const RectMatrix< Real, Row, Col > &m) |
template<typename Real , unsigned int Row, unsigned int Col> | |
RectMatrix< Real, Row, Col > | operator/ (const RectMatrix< Real, Row, Col > &m, Real s) |
template<typename Real , unsigned int Row, unsigned int Col> | |
Real | doubleDot (const RectMatrix< Real, Row, Col > &t1, const RectMatrix< Real, Row, Col > &t2) |
template<typename Real , unsigned int Row, unsigned int Col> | |
Vector< Real, Row > | mCross (const RectMatrix< Real, Row, Col > &t1, const RectMatrix< Real, Row, Col > &t2) |
template<typename Real , unsigned int Row, unsigned int Col> | |
std::ostream & | operator<< (std::ostream &o, const RectMatrix< Real, Row, Col > &m) |
template<typename Real > | |
SquareMatrix3< Real > | operator * (const SquareMatrix3< Real > &m1, const SquareMatrix3< Real > &m2) |
template<typename Real > | |
SquareMatrix3< Real > | outProduct (const Vector3< Real > &v1, const Vector3< Real > &v2) |
template<typename T > | |
bool | equal (T e1, T e2) |
template<> | |
bool | equal (RealType e1, RealType e2) |
template<typename Real , unsigned int Dim> | |
Vector< Real, Dim > | operator - (const Vector< Real, Dim > &v1) |
template<typename Real , unsigned int Dim> | |
Vector< Real, Dim > | operator+ (const Vector< Real, Dim > &v1, const Vector< Real, Dim > &v2) |
template<typename Real , unsigned int Dim> | |
Vector< Real, Dim > | operator - (const Vector< Real, Dim > &v1, const Vector< Real, Dim > &v2) |
template<typename Real , unsigned int Dim> | |
Vector< Real, Dim > | operator * (const Vector< Real, Dim > &v1, Real s) |
template<typename Real , unsigned int Dim> | |
Vector< Real, Dim > | operator * (Real s, const Vector< Real, Dim > &v1) |
template<typename Real , unsigned int Dim> | |
Vector< Real, Dim > | operator/ (const Vector< Real, Dim > &v1, Real s) |
template<typename Real , unsigned int Dim> | |
Real | dot (const Vector< Real, Dim > &v1, const Vector< Real, Dim > &v2) |
template<typename Real , unsigned int Dim> | |
Real | dot (const Vector< Real, Dim > &v1, const Vector< Real, Dim > &v2, const Vector< Real, Dim > &v3) |
template<typename Real , unsigned int Dim> | |
Real | distance (const Vector< Real, Dim > &v1, const Vector< Real, Dim > &v2) |
template<typename Real , unsigned int Dim> | |
Real | distanceSquare (const Vector< Real, Dim > &v1, const Vector< Real, Dim > &v2) |
template<typename Real , unsigned int Dim> | |
std::ostream & | operator<< (std::ostream &o, const Vector< Real, Dim > &v) |
template<typename Real > | |
Vector3< Real > | cross (const Vector3< Real > &v1, const Vector3< Real > &v2) |
template<typename Real > | |
Real | Vlinear (const Vector3< Real > &p, const Vector3< Real > &s) |
std::ostream & | operator<< (std::ostream &o, OptimizationFactory &factory) |
std::ostream & | operator<< (std::ostream &o, Molecule &mol) |
SelectionManager | operator| (const SelectionManager &sman1, const SelectionManager &sman2) |
SelectionManager | operator & (const SelectionManager &sman1, const SelectionManager &sman2) |
SelectionManager | operator^ (const SelectionManager &sman1, const SelectionManager &sman2) |
SelectionManager | operator- (const SelectionManager &sman1, const SelectionManager &sman2) |
SelectionSet | operator| (const SelectionSet &ss1, const SelectionSet &ss2) |
SelectionSet | operator & (const SelectionSet &ss1, const SelectionSet &ss2) |
SelectionSet | operator^ (const SelectionSet &ss1, const SelectionSet &ss2) |
SelectionSet | operator- (const SelectionSet &ss1, const SelectionSet &ss2) |
bool | operator== (const SelectionSet &ss1, const SelectionSet &ss2) |
std::ostream & | operator<< (std::ostream &os, const SelectionSet &ss) |
std::ostream & | operator<< (std::ostream &os, AmberImproperTorsionType &ott) |
std::ostream & | operator<< (std::ostream &os, HarmonicInversionType &hit) |
std::ostream & | operator<< (std::ostream &os, HarmonicTorsionType &htt) |
template<class ContainerType > | |
bool | hasDuplicateElement (const ContainerType &cont) |
std::ostream & | operator<< (std::ostream &os, OplsTorsionType &ott) |
std::ostream & | operator<< (std::ostream &os, PolynomialBendType &pbt) |
std::ostream & | operator<< (std::ostream &os, PolynomialBondType &pbt) |
std::ostream & | operator<< (std::ostream &os, TrappeTorsionType &ttt) |
template<typename InputIterator , typename OutputIterator , typename Predicate > | |
OutputIterator | copy_if (InputIterator first, InputIterator last, OutputIterator result, Predicate p) |
template<typename Container > | |
void | toLower (Container &cont, const std::locale &loc=std::locale()) |
template<typename Container > | |
Container | toLowerCopy (const Container &cont, const std::locale &loc=std::locale()) |
template<typename Container > | |
void | toUpper (Container &cont, const std::locale &loc=std::locale()) |
template<typename Container > | |
Container | toUpperCopy (const Container &cont, const std::locale &loc=std::locale()) |
template<typename O , typename I , typename C > | |
std::ostream & | operator<< (std::ostream &o, GenericFactory< O, I, C > &factory) |
template<class RandomAccessIterator > | |
bool | next_combination (std::vector< RandomAccessIterator > &iterContainer, RandomAccessIterator first, RandomAccessIterator last) |
STL next_permuation-like combination sequence generator. Given the first and last iterator of a sequence, next_combination iteratively generates all possible combinations. More... | |
OpenMDBitSet | operator| (const OpenMDBitSet &bs1, const OpenMDBitSet &bs2) |
OpenMDBitSet | operator & (const OpenMDBitSet &bs1, const OpenMDBitSet &bs2) |
OpenMDBitSet | operator^ (const OpenMDBitSet &bs1, const OpenMDBitSet &bs2) |
OpenMDBitSet | operator- (const OpenMDBitSet &bs1, const OpenMDBitSet &bs2) |
bool | operator== (const OpenMDBitSet &bs1, const OpenMDBitSet &bs2) |
std::ostream & | operator<< (std::ostream &os, const OpenMDBitSet &bs) |
CharClassification | isSpace (const std::locale &loc=std::locale()) |
CharClassification | isAlnum (const std::locale &loc=std::locale()) |
CharClassification | isAlpha (const std::locale &loc=std::locale()) |
CharClassification | isCntrl (const std::locale &loc=std::locale()) |
CharClassification | isDigit (const std::locale &loc=std::locale()) |
CharClassification | isGraph (const std::locale &loc=std::locale()) |
CharClassification | isLower (const std::locale &loc=std::locale()) |
CharClassification | isPrint (const std::locale &loc=std::locale()) |
CharClassification | isPunct (const std::locale &loc=std::locale()) |
CharClassification | isUpper (const std::locale &loc=std::locale()) |
CharClassification | isXDigit (const std::locale &loc=std::locale()) |
template<typename CharT > | |
FromRangeFunctor< CharT > | isFromRange (CharT from, CharT to) |
template<typename Pred1T , typename Pred2T > | |
PredAndFunctor< Pred1T, Pred2T > | operator && (const PredFacade< Pred1T > &Pred1, const PredFacade< Pred2T > &Pred2) |
template<typename Pred1T , typename Pred2T > | |
PredOrFunctor< Pred1T, Pred2T > | operator|| (const PredFacade< Pred1T > &Pred1, const PredFacade< Pred2T > &Pred2) |
template<typename PredT > | |
PredNotFunctor< PredT > | operator! (const PredFacade< PredT > &Pred) |
std::string | UpperCase (const std::string &S) |
std::string | LowerCase (const std::string &S) |
char * | trimSpaces (char *str) |
int | findBegin (std::istream &theStream, const char *startText) |
int | countTokens (char *line, char *delimiters) |
int | isEndLine (char *line) |
std::string | OpenMD_itoa (int value, unsigned int base) |
std::string | getPrefix (const std::string &str) |
std::string | getSuffix (const std::string &str) |
bool | isInteger (const std::string &str) |
bool | CaseInsensitiveEquals (const char ch1, const char ch2) |
size_t | CaseInsensitiveFind (const std::string &str1, const std::string &str2) |
unsigned long long | memparse (char *ptr, char **retptr) |
template<typename T > | |
std::string | toString (const T &v) |
template<typename T > | |
T | lexi_cast (const std::string &str) |
template<typename T > | |
bool | isType (const std::string &str) |
template<typename T > | |
std::string | to_string (T value) |
template<class ContainerType > | |
std::string | containerToString (const ContainerType &cont) |
void | trimLeft (std::string &str) |
void | trimRight (std::string &str) |
void | trim (std::string &str) |
std::string | trimLeftCopy (const std::string &input) |
std::string | trimRightCopy (const std::string &input) |
std::string | trimCopy (const std::string &input) |
template<typename P > | |
void | trimLeftIf (std::string &str, P pred) |
template<typename P > | |
void | trimRightIf (std::string &str, P pred) |
template<typename P > | |
void | trimIf (std::string &str, P pred) |
template<typename P > | |
std::string | trimLeftCopyIf (const std::string &input, P pred) |
template<typename P > | |
std::string | trimRightCopyIf (const std::string &input, P pred) |
template<typename P > | |
std::string | trimCopyIf (const std::string &input, P pred) |
template<class T1 , class T2 , class T3 > | |
tuple3< T1, T2, T3 > | make_tuple3 (T1 t1, T2 t2, T3 t3) |
template<class T1 , class T2 , class T3 , class T4 > | |
tuple4< T1, T2, T3, T4 > | make_tuple4 (T1 t1, T2 t2, T3 t3, T4 t4) |
template<class T1 , class T2 , class T3 > | |
bool | operator< (const tuple3< T1, T2, T3 > &t1, const tuple3< T1, T2, T3 > &t2) |
bool | operator< (const tuple3< int, int, std::vector< std::string > > &t1, const tuple3< int, int, std::vector< std::string > > &t2) |
template<class T1 , class T2 , class T3 , class T4 > | |
bool | operator< (const tuple4< T1, T2, T3, T4 > &t1, const tuple4< T1, T2, T3, T4 > &t2) |
bool | replaceWithWildCard (vector< vector< string >::iterator > &cont, vector< string > &sequence, vector< string > &result, const string &wildCard="X") |
iteratively replace the sequence with wild cards More... | |
RealType | roundMe (const RealType &x) |
Variables | |
const Mat3x3d | M3Zero (0.0) |
static const RealType | epsilon = 0.000001 |
const Vector3d | V3Zero (0.0, 0.0, 0.0) |
const Vector3d | V3X (1.0, 0.0, 0.0) |
const Vector3d | V3Y (0.0, 1.0, 0.0) |
const Vector3d | V3Z (0.0, 0.0, 1.0) |
static const int | ELECTROSTATIC_INTERACTION = (1 << 0) |
static const int | LJ_INTERACTION = (1 << 1) |
static const int | EAM_INTERACTION = (1 << 2) |
static const int | SC_INTERACTION = (1 << 3) |
static const int | STICKY_INTERACTION = (1 << 4) |
static const int | GB_INTERACTION = (1 << 5) |
static const int | MORSE_INTERACTION = (1 << 6) |
static const int | REPULSIVEPOWER_INTERACTION = (1 << 7) |
static const int | MAW_INTERACTION = (1 << 8) |
static const int | MIE_INTERACTION = (1 << 9) |
static const int | BUCKINGHAM_INTERACTION = (1 << 10) |
static const int | INVERSEPOWERSERIES_INTERACTION = (1 << 11) |
string const | DirectionalTypeID = "Directional" |
string const | EAMtypeID = "EAM" |
string const | FuncflTypeID = "FUNCFL" |
string const | ZhouTypeID = "ZHOU" |
string const | ZhouRoseTypeID = "ZHOUROSE" |
string const | FCtypeID = "Charge" |
string const | FQtypeID = "FlucQ" |
string const | GBtypeID = "GB" |
string const | LJtypeID = "LJ" |
string const | MultipoleTypeID = "Multipole" |
string const | PolTypeID = "Polarizable" |
string const | StickyTypeID = "Sticky" |
string const | SCtypeID = "SC" |
ElementsTable | etab |
const char * | progressSpinner_ = "|/-\\" |
typedef SimpleTypeData<bool> OpenMD::BoolGenericData |
BoolGenericData is a generic data type contains a bool variable
Definition at line 144 of file GenericData.hpp.
typedef DataStorageSnapshot::* OpenMD::DataStoragePointer |
Definition at line 331 of file Snapshot.hpp.
Definition at line 58 of file DirectionalAdapter.hpp.
typedef SimpleTypeData<RealType> OpenMD::DoubleGenericData |
DoubleGenericData is a generic data type contains a RealType variable
Definition at line 153 of file GenericData.hpp.
typedef Polynomial<RealType> OpenMD::DoublePolynomial |
Definition at line 669 of file Polynomial.hpp.
typedef VectorTypeData<RealType> OpenMD::DoubleVectorGenericData |
Definition at line 223 of file GenericData.hpp.
typedef SimpleTypeData<EAMParameters*> OpenMD::EAMData |
Definition at line 124 of file EAMAdapter.hpp.
Definition at line 57 of file FixedChargeAdapter.hpp.
typedef SimpleTypeData<float> OpenMD::FloatGenericData |
FloatGenericData is a generic data type contains a float variable
Definition at line 150 of file GenericData.hpp.
typedef VectorTypeData<float> OpenMD::FloatVectorGenericData |
Definition at line 217 of file GenericData.hpp.
Definition at line 69 of file FluctuatingChargeAdapter.hpp.
typedef SimpleTypeData<FuncflParameters*> OpenMD::FuncflData |
Definition at line 125 of file EAMAdapter.hpp.
Definition at line 62 of file GayBerneAdapter.hpp.
typedef SimpleTypeData<int> OpenMD::IntGenericData |
IntGenericData is a generic data type contains an integer variable
Definition at line 147 of file GenericData.hpp.
typedef tuple3<int, int, int> OpenMD::IntTuple3 |
typedef tuple4<int, int, int, int> OpenMD::IntTuple4 |
A generic data type contains a std::vector<int> variable.
A generic data type contains a std::vector<float> variable.
A generic data type contains a std::vector<RealType> variable.
Definition at line 211 of file GenericData.hpp.
Definition at line 59 of file LennardJonesAdapter.hpp.
typedef SquareMatrix3<RealType> OpenMD::Mat3x3d |
Definition at line 635 of file SquareMatrix3.hpp.
typedef SquareMatrix<RealType, 6> OpenMD::Mat6x6d |
Definition at line 452 of file SquareMatrix.hpp.
Definition at line 61 of file MultipoleAdapter.hpp.
Definition at line 57 of file PolarizableAdapter.hpp.
typedef Vector<RealType, N_INTERACTION_FAMILIES> OpenMD::potVec |
Definition at line 87 of file NonBondedInteraction.hpp.
typedef Quaternion<RealType> OpenMD::Quat4d |
Definition at line 651 of file Quaternion.hpp.
typedef SimpleTypeData<Restraint*> OpenMD::RestraintData |
Definition at line 174 of file Restraint.hpp.
typedef SquareMatrix3<RealType> OpenMD::RotMat3x3d |
Definition at line 636 of file SquareMatrix3.hpp.
Definition at line 61 of file SuttonChenAdapter.hpp.
Definition at line 108 of file ShapeAtomType.hpp.
typedef std::pair<int, int> OpenMD::SnapshotBlock |
Definition at line 50 of file BlockSnapshotManager.hpp.
Definition at line 64 of file StickyAdapter.hpp.
A generic data type contains a std::string variable
Definition at line 176 of file GenericData.hpp.
A generic data type contains a std::vector<string> variable.
Definition at line 253 of file GenericData.hpp.
typedef basic_teebuf<char> OpenMD::TeeBuf |
Definition at line 117 of file basic_teebuf.hpp.
typedef Vector3<RealType> OpenMD::Vector3d |
Definition at line 160 of file Vector3.hpp.
typedef Vector3<int> OpenMD::Vector3i |
Definition at line 158 of file Vector3.hpp.
typedef std::list<std::pair<BaseVisitor*, int> >::iterator OpenMD::VisitorIterator |
Definition at line 53 of file CompositeVisitor.hpp.
typedef SimpleTypeData<ZhouParameters*> OpenMD::ZhouData |
Definition at line 126 of file EAMAdapter.hpp.
Enumerator | |
---|---|
btTraditional | |
btModified | |
btUnknown |
Definition at line 68 of file BuckinghamInteractionType.hpp.
Enumerator | |
---|---|
Global | |
Row | |
Column |
Definition at line 61 of file Communicator.hpp.
enum OpenMD::CutoffMethod |
Enumerator | |
---|---|
HARD | |
SWITCHED | |
SHIFTED_POTENTIAL | |
SHIFTED_FORCE | |
TAYLOR_SHIFTED | |
EWALD_FULL |
Definition at line 48 of file Cutoffs.hpp.
enum OpenMD::EAMiType |
Enumerator | |
---|---|
eamitTabulated | |
eamitZhou |
Definition at line 63 of file EAMInteractionType.hpp.
enum OpenMD::EAMType |
Enumerator | |
---|---|
eamFuncfl | |
eamZhou2001 | |
eamZhou2004 | |
eamZhou2005 | |
eamZhou2005Oxygen | |
eamZhouRose | |
eamUnknown |
Definition at line 57 of file EAMAdapter.hpp.
Enumerator | |
---|---|
UNDAMPED | |
DAMPED |
Definition at line 85 of file Electrostatic.hpp.
Definition at line 73 of file Electrostatic.hpp.
The InteractionFamily enum.
This is used to sort different types of non-bonded interaction and to prevent multiple interactions in the same family from being applied to any given pair of atom types.
Definition at line 58 of file NonBondedInteraction.hpp.
enum OpenMD::InversionKey |
Enumerator | |
---|---|
itAngle | |
itCosAngle |
Definition at line 55 of file InversionType.hpp.
enum OpenMD::MorseType |
Enumerator | |
---|---|
mtShifted | |
mtRepulsive | |
mtUnknown |
Definition at line 63 of file MorseInteractionType.hpp.
Enumerator | |
---|---|
odhAverage | |
odhTotal | |
odhLastValue | |
odhMax | |
odhUnknownDataHandling |
Definition at line 58 of file StaticAnalyser.hpp.
Enumerator | |
---|---|
odtReal | |
odtVector3 | |
odtArray2d | |
odtUnknownDataType |
Definition at line 51 of file StaticAnalyser.hpp.
The SelectionType enum.
This is used to sort different types of selections by object type
Enumerator | |
---|---|
STUNTDOUBLE | StuntDoubles (Atoms & RigidBodies) |
BOND | Bonds |
BEND | Bends |
TORSION | Torsions |
INVERSION | Inversions |
MOLECULE | Molecules |
N_SELECTIONTYPES |
Definition at line 57 of file SelectionSet.hpp.
Enumerator | |
---|---|
cubic | |
fifth_order_poly |
Definition at line 51 of file SwitchingFunction.hpp.
bool OpenMD::CaseInsensitiveEquals | ( | const char | ch1, |
const char | ch2 | ||
) |
Definition at line 253 of file StringUtils.cpp.
Referenced by CaseInsensitiveFind().
size_t OpenMD::CaseInsensitiveFind | ( | const std::string & | str1, |
const std::string & | str2 | ||
) |
Definition at line 257 of file StringUtils.cpp.
References CaseInsensitiveEquals().
Referenced by OpenMD::SimCreator::createSim().
void OpenMD::CholeskyDecomposition | ( | MatrixType & | A, |
MatrixType & | L | ||
) |
Definition at line 52 of file CholeskyDecomposition.hpp.
References epsilon.
Referenced by OpenMD::HydroProp::complete(), OpenMD::BAOAB::getMomentData(), OpenMD::HydroProp::HydroProp(), and OpenMD::LangevinHullForceManager::postCalculation().
std::string OpenMD::containerToString | ( | const ContainerType & | cont | ) |
Definition at line 168 of file StringUtils.hpp.
Referenced by OpenMD::MoleculeStamp::checkBends(), OpenMD::MoleculeStamp::checkInversions(), OpenMD::MoleculeStamp::checkTorsions(), OpenMD::BondStamp::setMembers(), OpenMD::ConstraintStamp::setMembers(), OpenMD::TorsionStamp::setMembers(), OpenMD::BendStamp::setMembers(), OpenMD::AtomStamp::setOrientation(), and OpenMD::AtomStamp::setPosition().
OutputIterator OpenMD::copy_if | ( | InputIterator | first, |
InputIterator | last, | ||
OutputIterator | result, | ||
Predicate | p | ||
) |
copy_if is one of the missing functions in STL. copy_if copies elements from the range [first, last) to the range [result, result + (last-first)), except that any element for which pred is false is not copied.
Definition at line 53 of file Algorithm.hpp.
Referenced by OpenMD::BlockSnapshotManager::getActiveBlocks().
int OpenMD::countTokens | ( | char * | line, |
char * | delimiters | ||
) |
Counts the number of tokens on line which are delimited by the characters listed in delimiters
line | |
delimiters |
Definition at line 157 of file StringUtils.cpp.
|
inline |
Returns the cross product of two Vectors
v1 | first vector |
v2 | second vector |
Definition at line 135 of file Vector3.hpp.
References OpenMD::Vector3< Real >::x(), OpenMD::Vector3< Real >::y(), and OpenMD::Vector3< Real >::z().
Referenced by OpenMD::Quaternion< Real >::align(), OpenMD::MagneticField::applyPerturbation(), OpenMD::UniformField::applyPerturbation(), OpenMD::UniformGradient::applyPerturbation(), OpenMD::ForceMatrixDecomposition::buildNeighborList(), OpenMD::SHAPES::calcForce(), OpenMD::GhostBend::calcForce(), OpenMD::GhostTorsion::calcForce(), OpenMD::GB::calcForce(), OpenMD::Inversion::calcForce(), OpenMD::Torsion::calcForce(), OpenMD::MolecularRestraint::calcForce(), OpenMD::Electrostatic::calcForce(), OpenMD::SCDElem::calcSCD(), OpenMD::Electrostatic::calcSurfaceTerm(), OpenMD::RNEMD::collectData(), OpenMD::Triangle::computeHydrodynamicTensor(), OpenMD::Triangle::computeNormal(), OpenMD::MomAngMomCorrFunc::computeProperty2(), OpenMD::Triangle::computeUnitNormal(), OpenMD::SquareMatrix3< RealType >::diagonalize(), OpenMD::RNEMD::doVSS(), OpenMD::Field< Vector3d >::Field(), OpenMD::Quaternion< Real >::fromShortestArc(), OpenMD::Thermo::getAngularMomentum(), OpenMD::Thermo::getInertiaTensor(), OpenMD::Inversion::getValue(), OpenMD::Torsion::getValue(), OpenMD::Snapshot::getXYarea(), OpenMD::Snapshot::getXZarea(), OpenMD::Snapshot::getYZarea(), OpenMD::GofXyz::initializeHistogram(), OpenMD::LDForceManager::postCalculation(), OpenMD::Field< Vector3d >::processFrame(), OpenMD::RNEMDR::processFrame(), OpenMD::RNEMDRTheta::processFrame(), OpenMD::Electrostatic::ReciprocalSpaceSum(), OpenMD::Velocitizer::removeAngularDrift(), and OpenMD::LipidTransVisitor::update().
|
inline |
Returns the distance between two DynamicVectors
v1 | first vector |
v2 | second vector |
Definition at line 429 of file DynamicVector.hpp.
References OpenMD::DynamicVector< Real, Alloc >::length().
Referenced by OpenMD::SelectionCompiler::clauseWithin(), OpenMD::RNEMD::collectData(), OpenMD::GofRAngle2::collectHistogram(), OpenMD::GofR::collectHistogram(), OpenMD::GofZ::collectHistogram(), OpenMD::Kirkwood::collectHistogram(), OpenMD::TwoDGofR::collectHistogram(), OpenMD::GofRZ::collectHistogram(), OpenMD::GofRAngle::collectHistogram(), OpenMD::KirkwoodQuadrupoles::collectHistogram(), OpenMD::MoleculeCreator::createMolecule(), OpenMD::GCNSeq::doFrame(), OpenMD::DistanceFinder::find(), OpenMD::shapedLatticeRod::isInterior(), main(), OpenMD::MultipoleSum::process(), OpenMD::RhoR::process(), OpenMD::AngleR::process(), OpenMD::CoordinationNumber::process(), and OpenMD::SelectionEvaluator::withinInstruction().
|
inline |
Returns the distance between two Vectors
v1 | first vector |
v2 | second vector |
Definition at line 556 of file Vector.hpp.
References OpenMD::Vector< Real, Dim >::length().
|
inline |
Returns the squared distance between two DynamicVectors
v1 | first vector |
v2 | second vector |
Definition at line 441 of file DynamicVector.hpp.
References OpenMD::DynamicVector< Real, Alloc >::lengthSquare().
|
inline |
Returns the squared distance between two Vectors
v1 | first vector |
v2 | second vector |
Definition at line 568 of file Vector.hpp.
References OpenMD::Vector< Real, Dim >::lengthSquare().
|
inline |
Returns the dot product of two DynamicVectors
v1 | first vector |
v2 | second vector |
Definition at line 412 of file DynamicVector.hpp.
Referenced by OpenMD::Quaternion< Real >::align(), OpenMD::UniformField::applyPerturbation(), OpenMD::UniformGradient::applyPerturbation(), OpenMD::ForceMatrixDecomposition::buildNeighborList(), OpenMD::VCorrFunc::calcCorrVal(), OpenMD::ThetaCorrFunc::calcCorrVal(), OpenMD::DipoleCorrFunc::calcCorrVal(), OpenMD::MomAngMomCorrFunc::calcCorrVal(), OpenMD::LegendreCorrFunc::calcCorrVal(), OpenMD::LegendreCorrFuncZ::calcCorrVal(), OpenMD::FreqFlucCorrFunc::calcCorrVal(), OpenMD::DirectionalRCorrFunc::calcCorrVals(), OpenMD::COHZ::calcCorrVals(), OpenMD::SHAPES::calcForce(), OpenMD::GhostBend::calcForce(), OpenMD::GhostTorsion::calcForce(), OpenMD::GB::calcForce(), OpenMD::Inversion::calcForce(), OpenMD::Torsion::calcForce(), OpenMD::MolecularRestraint::calcForce(), OpenMD::Bend::calcForce(), OpenMD::Electrostatic::calcForce(), OpenMD::SCDElem::calcSCD(), OpenMD::Electrostatic::calcSurfaceTerm(), OpenMD::RMSD::calculate_rmsd(), OpenMD::GofRAngle2::collectHistogram(), OpenMD::Kirkwood::collectHistogram(), OpenMD::GofAngle2::collectHistogram(), OpenMD::HBondPersistence::computeFrame(), QuantLib::Problem::computeGradientNormValue(), OpenMD::VCorrFuncR::computeProperty1(), OpenMD::Rattle::constraintPairA(), OpenMD::Rattle::constraintPairB(), OpenMD::Shake::constraintPairF(), OpenMD::Shake::constraintPairR(), OpenMD::SystemDipoleCorrFunc::correlateFrames(), QuantLib::Problem::DotProduct(), OpenMD::RNEMD::doVSS(), OpenMD::GofRTheta::evaluateAngle(), OpenMD::GofROmega::evaluateAngle(), OpenMD::Field< Vector3d >::Field(), OpenMD::Quaternion< Real >::fromShortestArc(), OpenMD::RNEMDRTheta::getBins(), OpenMD::NanoLength::getLength(), OpenMD::LDForceManager::getMomentData(), OpenMD::BAOAB::getMomentData(), OpenMD::Electrostatic::getSitePotentials(), OpenMD::Inversion::getValue(), OpenMD::Torsion::getValue(), OpenMD::Bend::getValue(), OpenMD::UniformGradient::initialize(), OpenMD::HBondJump::isHBond(), OpenMD::DynamicVector< Real >::lengthSquare(), OpenMD::Vector< RealType, N_INTERACTION_FAMILIES >::lengthSquare(), OpenMD::ForceManager::longRangeInteractions(), OpenMD::LDForceManager::postCalculation(), OpenMD::FreqFlucCorrFunc::preCorrelate(), OpenMD::MultipoleSum::process(), OpenMD::P2OrderParameter::process(), OpenMD::NitrileFrequencyMap::process(), OpenMD::AngleR::process(), OpenMD::BondAngleDistribution::process(), OpenMD::pAngle::process(), OpenMD::HBondGeometric::process(), OpenMD::TetrahedralityHBMatrix::process(), OpenMD::TetrahedralityParamXYZ::process(), OpenMD::TetrahedralityParam::process(), OpenMD::TetrahedralityParamZ::process(), OpenMD::TetrahedralityParamDens::process(), OpenMD::Field< Vector3d >::processFrame(), OpenMD::RNEMDRTheta::processFrame(), OpenMD::Electrostatic::ReciprocalSpaceSum(), and OpenMD::ForceManager::selectedLongRangeInteractions().
|
inline |
Returns the dot product of two Vectors
v1 | first vector |
v2 | second vector |
Definition at line 515 of file Vector.hpp.
|
inline |
Returns the wide dot product of three Vectors. Compare with Rapaport's VWDot function.
v1 | first vector |
v2 | second vector |
v3 | third vector |
Definition at line 538 of file Vector.hpp.
|
inline |
Returns the tensor contraction (double dot product) of two rank 2 tensors (or Matrices)
t1 | first tensor |
t2 | second tensor |
Definition at line 580 of file RectMatrix.hpp.
References Row.
Referenced by OpenMD::UniformGradient::applyPerturbation(), OpenMD::KirkwoodQuadrupoles::collectHistogram(), and OpenMD::MultipoleSum::process().
|
inline |
Definition at line 62 of file Vector.hpp.
|
inline |
Definition at line 72 of file Vector.hpp.
References epsilon.
bool OpenMD::equal | ( | const Polynomial< Real > & | p1, |
const Polynomial< Real > & | p2 | ||
) |
Tests if two polynomial have the same exponents
p1 | the first polynomial |
p2 | the second polynomial |
Definition at line 649 of file Polynomial.hpp.
References OpenMD::Polynomial< Real >::begin(), OpenMD::Polynomial< Real >::end(), and OpenMD::Polynomial< Real >::size().
Referenced by _cmdline_parser_configfile(), OpenMD::DynamicVector< Real >::isNormalized(), OpenMD::Vector< RealType, N_INTERACTION_FAMILIES >::isNormalized(), OpenMD::DynamicVector< Real >::operator==(), OpenMD::Quaternion< Real >::operator==(), OpenMD::Vector< RealType, N_INTERACTION_FAMILIES >::operator==(), OpenMD::RectMatrix< Real, Dim, Dim >::operator==(), OpenMD::DynamicRectMatrix< RealType >::operator==(), and operator==().
Real OpenMD::fastpow | ( | Real | x, |
int | N | ||
) |
Definition at line 63 of file Polynomial.hpp.
Referenced by OpenMD::RNEMD::doNIVS(), OpenMD::Polynomial< RealType >::evaluate(), and OpenMD::Polynomial< RealType >::evaluateDerivative().
int OpenMD::findBegin | ( | std::istream & | theStream, |
const char * | startText | ||
) |
Finds the location of the string "begin <startText>" in an input stream.
theStream | |
startText |
Definition at line 107 of file StringUtils.cpp.
Referenced by OpenMD::ShapeAtomTypesSectionParser::parseShapeFile().
Polynomial<Real>* OpenMD::getDerivative | ( | const Polynomial< Real > & | p1 | ) |
Returns the first derivative of this polynomial.
Definition at line 625 of file Polynomial.hpp.
References OpenMD::Polynomial< Real >::begin(), OpenMD::Polynomial< Real >::end(), and OpenMD::Polynomial< Real >::setCoefficient().
int OpenMD::getGlobalCountOfType | ( | AtomType * | atype | ) |
Definition at line 803 of file SimInfo.cpp.
std::string OpenMD::getPrefix | ( | const std::string & | str | ) |
Definition at line 229 of file StringUtils.cpp.
Referenced by OpenMD::ActionCorrFunc::ActionCorrFunc(), OpenMD::AngleR::AngleR(), OpenMD::BondAngleDistribution::BondAngleDistribution(), OpenMD::BondCorrFunc::BondCorrFunc(), OpenMD::BondOrderParameter::BondOrderParameter(), OpenMD::BOPofR::BOPofR(), OpenMD::CenterOfMass::CenterOfMass(), OpenMD::ChargeField::ChargeField(), OpenMD::ChargeHistogram::ChargeHistogram(), OpenMD::COHZ::COHZ(), OpenMD::COMVel::COMVel(), OpenMD::ContactAngle1::ContactAngle1(), OpenMD::ContactAngle2::ContactAngle2(), OpenMD::CoordinationNumber::CoordinationNumber(), OpenMD::CurrentDensity::CurrentDensity(), OpenMD::DensityField::DensityField(), OpenMD::DensityHistogram::DensityHistogram(), OpenMD::DensityPlot::DensityPlot(), OpenMD::DipoleCorrFunc::DipoleCorrFunc(), OpenMD::DipoleField::DipoleField(), OpenMD::DirectionalRCorrFunc::DirectionalRCorrFunc(), OpenMD::Displacement::Displacement(), OpenMD::DisplacementZ::DisplacementZ(), OpenMD::EnergyCorrFunc::EnergyCorrFunc(), OpenMD::Equipartition::Equipartition(), OpenMD::Field< Vector3d >::Field(), OpenMD::ForceAutoCorrFunc::ForceAutoCorrFunc(), OpenMD::ForTorCorrFunc::ForTorCorrFunc(), OpenMD::FreqFlucCorrFunc::FreqFlucCorrFunc(), OpenMD::SimCreator::gatherParameters(), OpenMD::GCN::GCN(), OpenMD::GCNSeq::GCNSeq(), OpenMD::GofAngle2::GofAngle2(), OpenMD::GofR::GofR(), OpenMD::GofRAngle2::GofRAngle2(), OpenMD::GofROmega::GofROmega(), OpenMD::GofRTheta::GofRTheta(), OpenMD::GofRZ::GofRZ(), OpenMD::GofXyz::GofXyz(), OpenMD::GofZ::GofZ(), OpenMD::HBondGeometric::HBondGeometric(), OpenMD::HBondJump::HBondJump(), OpenMD::HBondJumpZ::HBondJumpZ(), OpenMD::HBondPersistence::HBondPersistence(), OpenMD::Hxy::Hxy(), OpenMD::Kirkwood::Kirkwood(), OpenMD::KirkwoodQuadrupoles::KirkwoodQuadrupoles(), OpenMD::LegendreCorrFunc::LegendreCorrFunc(), OpenMD::LegendreCorrFuncZ::LegendreCorrFuncZ(), main(), OpenMD::MomAngMomCorrFunc::MomAngMomCorrFunc(), OpenMD::MomentumCorrFunc::MomentumCorrFunc(), OpenMD::MomentumHistogram::MomentumHistogram(), OpenMD::MultipoleSum::MultipoleSum(), OpenMD::NanoLength::NanoLength(), OpenMD::NanoVolume::NanoVolume(), OpenMD::NitrileFrequencyMap::NitrileFrequencyMap(), OpenMD::ObjectCount::ObjectCount(), OpenMD::P2OrderParameter::P2OrderParameter(), OpenMD::pAngle::pAngle(), OpenMD::PipeDensity::PipeDensity(), OpenMD::PotDiff::PotDiff(), OpenMD::SurfaceDiffusion::process(), OpenMD::Rattle::Rattle(), OpenMD::RCorrFunc::RCorrFunc(), OpenMD::RCorrFuncR::RCorrFuncR(), OpenMD::RCorrFuncZ::RCorrFuncZ(), OpenMD::RestraintForceManager::RestraintForceManager(), OpenMD::RhoR::RhoR(), OpenMD::RhoZ::RhoZ(), OpenMD::RippleOP::RippleOP(), OpenMD::RNEMDR::RNEMDR(), OpenMD::RNEMDRTheta::RNEMDRTheta(), OpenMD::RNEMDZ::RNEMDZ(), OpenMD::SCDOrderParameter::SCDOrderParameter(), OpenMD::SCN::SCN(), OpenMD::SelectionCorrFunc::SelectionCorrFunc(), OpenMD::Shake::Shake(), OpenMD::SpatialStatistics::SpatialStatistics(), OpenMD::StressCorrFunc::StressCorrFunc(), OpenMD::SystemDipoleCorrFunc::SystemDipoleCorrFunc(), OpenMD::TetrahedralityHBMatrix::TetrahedralityHBMatrix(), OpenMD::TetrahedralityParam::TetrahedralityParam(), OpenMD::TetrahedralityParamDens::TetrahedralityParamDens(), OpenMD::TetrahedralityParamXYZ::TetrahedralityParamXYZ(), OpenMD::TetrahedralityParamZ::TetrahedralityParamZ(), OpenMD::ThetaCorrFunc::ThetaCorrFunc(), OpenMD::TorForCorrFunc::TorForCorrFunc(), OpenMD::TorqueAutoCorrFunc::TorqueAutoCorrFunc(), OpenMD::TwoDGofR::TwoDGofR(), OpenMD::VCorrFunc::VCorrFunc(), OpenMD::VCorrFuncR::VCorrFuncR(), OpenMD::VCorrFuncZ::VCorrFuncZ(), OpenMD::VelocityField::VelocityField(), OpenMD::VelocityZ::VelocityZ(), OpenMD::WCorrFunc::WCorrFunc(), OpenMD::ZConsReader::ZConsReader(), OpenMD::ZconstraintForceManager::ZconstraintForceManager(), and OpenMD::ZConsVisitor::ZConsVisitor().
std::string OpenMD::getSuffix | ( | const std::string & | str | ) |
Definition at line 233 of file StringUtils.cpp.
bool OpenMD::hasDuplicateElement | ( | const ContainerType & | cont | ) |
Definition at line 52 of file MoleculeStamp.cpp.
Referenced by OpenMD::MoleculeStamp::checkBends(), OpenMD::MoleculeStamp::checkInversions(), and OpenMD::MoleculeStamp::checkTorsions().
bool OpenMD::invertMatrix | ( | MatrixType & | A, |
MatrixType & | AI | ||
) |
Invert input square matrix A into matrix AI.
A | input square matrix |
AI | output square matrix |
Definition at line 94 of file LU.hpp.
Referenced by OpenMD::ApproximationModel::calcHydroPropsAtCD(), OpenMD::ApproximationModel::calcHydroPropsAtCR(), OpenMD::Sphere::getHydroProp(), OpenMD::Ellipsoid::getHydroProp(), OpenMD::SquareMatrix< RealType, 6 >::inverse(), and main().
bool OpenMD::invertMatrix | ( | MatrixType & | A, |
MatrixType & | AI, | ||
unsigned int | size, | ||
int * | tmp1Size, | ||
typename MatrixType::ElemPoinerType | tmp2Size | ||
) |
Invert input square matrix A into matrix AI (Thread safe versions).
A | input square matrix |
AI | output square matrix |
size | size of the matrix and temporary arrays |
tmp1Size | temporary array |
tmp2Size | temporary array |
Definition at line 137 of file LU.hpp.
References LUFactorLinearSystem(), and LUSolveLinearSystem().
|
inline |
Definition at line 107 of file Predicate.hpp.
|
inline |
Definition at line 111 of file Predicate.hpp.
|
inline |
Definition at line 115 of file Predicate.hpp.
|
inline |
Definition at line 119 of file Predicate.hpp.
int OpenMD::isEndLine | ( | char * | line | ) |
discovers whether or not the line contains the "end" token
line | The line to test |
Definition at line 177 of file StringUtils.cpp.
Referenced by OpenMD::ShapeAtomTypesSectionParser::parseShapeFile().
|
inline |
Definition at line 368 of file ParamConstraint.hpp.
EqualIgnoreCaseConstraint OpenMD::isEqualIgnoreCase | ( | std::string | str | ) |
Definition at line 74 of file ParamConstraint.cpp.
Referenced by OpenMD::MinimizerParameters::validate(), OpenMD::RestraintStamp::validate(), OpenMD::FluctuatingChargeParameters::validate(), OpenMD::RNEMDParameters::validate(), OpenMD::Globals::validate(), and OpenMD::ForceFieldOptions::validateOptions().
EvenConstraint OpenMD::isEven | ( | ) |
Definition at line 78 of file ParamConstraint.cpp.
|
inline |
Definition at line 149 of file Predicate.hpp.
|
inline |
Definition at line 123 of file Predicate.hpp.
|
inline |
Definition at line 373 of file ParamConstraint.hpp.
Referenced by OpenMD::MinimizerParameters::validate(), and OpenMD::FluctuatingChargeParameters::validate().
|
inline |
Definition at line 378 of file ParamConstraint.hpp.
Referenced by OpenMD::MinimizerParameters::validate().
bool OpenMD::isInteger | ( | const std::string & | str | ) |
Definition at line 237 of file StringUtils.cpp.
Referenced by OpenMD::SelectionCompiler::lookingAtInteger(), and OpenMD::ForceFieldOptions::setData().
|
inline |
Definition at line 358 of file ParamConstraint.hpp.
|
inline |
Definition at line 363 of file ParamConstraint.hpp.
Referenced by OpenMD::MinimizerParameters::validate(), and OpenMD::FluctuatingChargeParameters::validate().
|
inline |
Definition at line 127 of file Predicate.hpp.
NegativeConstraint OpenMD::isNegative | ( | ) |
Definition at line 66 of file ParamConstraint.cpp.
NonNegativeConstraint OpenMD::isNonNegative | ( | ) |
Definition at line 70 of file ParamConstraint.cpp.
Referenced by OpenMD::Component::validate(), OpenMD::ZConsStamp::validate(), OpenMD::RestraintStamp::validate(), OpenMD::FluctuatingChargeParameters::validate(), OpenMD::ConstraintStamp::validate(), OpenMD::BendStamp::validate(), OpenMD::TorsionStamp::validate(), and OpenMD::Globals::validate().
NonPositiveConstraint OpenMD::isNonPositive | ( | ) |
Definition at line 62 of file ParamConstraint.cpp.
ParamConstraintFacade< NonZeroConstraint > OpenMD::isNonZero | ( | ) |
Definition at line 54 of file ParamConstraint.cpp.
NotEmptyConstraint OpenMD::isNotEmpty | ( | ) |
Definition at line 46 of file ParamConstraint.cpp.
Referenced by OpenMD::FragmentStamp::validate(), OpenMD::Component::validate(), OpenMD::AtomStamp::validate(), and OpenMD::Globals::validate().
PositiveConstraint OpenMD::isPositive | ( | ) |
Definition at line 58 of file ParamConstraint.cpp.
Referenced by OpenMD::Component::validate(), OpenMD::MinimizerParameters::validate(), OpenMD::FluctuatingChargeParameters::validate(), OpenMD::RNEMDParameters::validate(), and OpenMD::Globals::validate().
|
inline |
Definition at line 131 of file Predicate.hpp.
|
inline |
Definition at line 135 of file Predicate.hpp.
|
inline |
Definition at line 103 of file Predicate.hpp.
Referenced by trim(), trimCopy(), trimLeft(), trimLeftCopy(), trimRight(), and trimRightCopy().
bool OpenMD::isType | ( | const std::string & | str | ) |
Definition at line 139 of file StringUtils.hpp.
|
inline |
Definition at line 139 of file Predicate.hpp.
|
inline |
Definition at line 143 of file Predicate.hpp.
ZeroConstraint OpenMD::isZero | ( | ) |
Definition at line 50 of file ParamConstraint.cpp.
T OpenMD::lexi_cast | ( | const std::string & | str | ) |
Definition at line 128 of file StringUtils.hpp.
std::string OpenMD::LowerCase | ( | const std::string & | S | ) |
int OpenMD::LUFactorLinearSystem | ( | MatrixType & | A, |
int * | index, | ||
int | size, | ||
typename MatrixType::ElemPoinerType | tmpSize | ||
) |
Factor linear equations Ax = b using LU decompostion A = LU where L is lower triangular matrix and U is upper triangular matrix.
A | input square matrix |
index | pivot indices |
size | size of the matrix and temporary arrays |
tmpSize | temporary array |
Definition at line 183 of file LU.hpp.
References epsilon.
Referenced by invertMatrix().
void OpenMD::LUSolveLinearSystem | ( | MatrixType & | A, |
int * | index, | ||
typename MatrixType::ElemPoinerType | x, | ||
int | size | ||
) |
Solve linear equations Ax = b using LU decompostion A = LU where L is lower triangular matrix and U is upper triangular matrix.
A | input square matrix |
index | pivot indices |
x | vector |
size | size of the matrix and temporary arrays |
Definition at line 277 of file LU.hpp.
Referenced by invertMatrix().
tuple3<T1,T2,T3> OpenMD::make_tuple3 | ( | T1 | t1, |
T2 | t2, | ||
T3 | t3 | ||
) |
Definition at line 61 of file Tuple.hpp.
Referenced by OpenMD::ForceField::getBendType(), OpenMD::ForceField::getInversionType(), OpenMD::ForceField::getTorsionType(), OpenMD::Icosahedron::Icosahedron(), and OpenMD::SCDElem::SCDElem().
tuple4<T1,T2,T3,T4> OpenMD::make_tuple4 | ( | T1 | t1, |
T2 | t2, | ||
T3 | t3, | ||
T4 | t4 | ||
) |
|
inline |
Returns the vector (cross) product of two matrices. This operation is defined in:
W. Smith, "Point Multipoles in the Ewald Summation (Revisited)," CCP5 Newsletter No 46., pp. 18-30.
Equation 21 defines:
where
and
are regarded as cyclic permuations of the matrix indices (i.e. for a 3x3 matrix, when
,
, and
).
t1 | first matrix |
t2 | second matrix |
Definition at line 616 of file RectMatrix.hpp.
References Row.
Referenced by OpenMD::UniformGradient::applyPerturbation(), and OpenMD::Electrostatic::calcForce().
unsigned long long OpenMD::memparse | ( | char * | ptr, |
char ** | retptr | ||
) |
memparse - parse a string with mem suffixes into a number
ptr | Where parse begins |
retptr | (output) Pointer to next char after parse completes |
Parses a string into a number. The number stored at ptr is potentially suffixed with K (for kilobytes, or 1024 bytes), M (for megabytes, or 1048576 bytes), or G (for gigabytes, or 1073741824). If the number is suffixed with K, M, or G, then the return value is the number multiplied by one kilobyte, one megabyte, or one gigabyte, respectively.
Definition at line 277 of file StringUtils.cpp.
Referenced by main().
bool OpenMD::next_combination | ( | std::vector< RandomAccessIterator > & | iterContainer, |
RandomAccessIterator | first, | ||
RandomAccessIterator | last | ||
) |
STL next_permuation-like combination sequence generator. Given the first and last iterator of a sequence, next_combination iteratively generates all possible combinations.
iterContainer | iterator container |
first | the first iterator |
last | the last iterator |
random access iterators . And all of the iteratos in iterContainer must be within the range [first, last)
Definition at line 87 of file next_combination.hpp.
Referenced by replaceWithWildCard().
std::string OpenMD::OpenMD_itoa | ( | int | value, |
unsigned int | base | ||
) |
Definition at line 198 of file StringUtils.cpp.
References MATPACK::sign().
Referenced by OpenMD::MoleculeCreator::createBend(), OpenMD::MoleculeCreator::createBond(), OpenMD::MoleculeCreator::createInversion(), OpenMD::MoleculeCreator::createRigidBody(), OpenMD::MoleculeCreator::createTorsion(), and OpenMD::SCDOrderParameter::SCDOrderParameter().
|
inline |
Definition at line 339 of file ParamConstraint.hpp.
OpenMDBitSet OpenMD::operator & | ( | const OpenMDBitSet & | bs1, |
const OpenMDBitSet & | bs2 | ||
) |
Definition at line 171 of file OpenMDBitSet.cpp.
References OpenMD::OpenMDBitSet::size().
SelectionSet OpenMD::operator & | ( | const SelectionSet & | ss1, |
const SelectionSet & | ss2 | ||
) |
Definition at line 196 of file SelectionSet.cpp.
References OpenMD::SelectionSet::bitsets_, and N_SELECTIONTYPES.
SelectionManager OpenMD::operator & | ( | const SelectionManager & | sman1, |
const SelectionManager & | sman2 | ||
) |
Definition at line 615 of file SelectionManager.cpp.
|
inline |
Definition at line 155 of file Predicate.hpp.
|
inline |
Definition at line 320 of file ParamConstraint.hpp.
DynamicVector<Real> OpenMD::operator * | ( | const DynamicVector< Real > & | v1, |
Real | s | ||
) |
Returns the vaule of scalar multiplication of this vector v1 (v1 * r).
v1 | the source vector |
s | the scalar value |
Definition at line 373 of file DynamicVector.hpp.
References OpenMD::DynamicVector< Real, Alloc >::mul().
DynamicVector<Real> OpenMD::operator * | ( | Real | s, |
const DynamicVector< Real > & | v1 | ||
) |
Returns the vaule of scalar multiplication of this vector v1 (v1 * r).
s | the scalar value |
v1 | the source vector |
Definition at line 386 of file DynamicVector.hpp.
References OpenMD::DynamicVector< Real, Alloc >::mul().
Vector<Real, Dim> OpenMD::operator * | ( | const Vector< Real, Dim > & | v1, |
Real | s | ||
) |
Returns the vaule of scalar multiplication of this vector v1 (v1 * r).
v1 | the source vector |
s | the scalar value |
Definition at line 476 of file Vector.hpp.
References OpenMD::Vector< Real, Dim >::mul().
|
inline |
Return the multiplication of scalra and matrix (m * s).
m | the matrix |
s | the scalar |
Definition at line 478 of file RectMatrix.hpp.
References OpenMD::RectMatrix< Real, Row, Col >::mul().
Vector<Real, Dim> OpenMD::operator * | ( | Real | s, |
const Vector< Real, Dim > & | v1 | ||
) |
Returns the vaule of scalar multiplication of this vector v1 (v1 * r).
s | the scalar value |
v1 | the source vector |
Definition at line 489 of file Vector.hpp.
References OpenMD::Vector< Real, Dim >::mul().
|
inline |
Return the multiplication of a scalra and a matrix (s * m).
s | the scalar |
m | the matrix |
Definition at line 493 of file RectMatrix.hpp.
References OpenMD::RectMatrix< Real, Row, Col >::mul().
|
inline |
Return the multiplication of two matrixes (m1 * m2).
m1 | the first matrix |
m2 | the second matrix |
Definition at line 508 of file RectMatrix.hpp.
References Row.
|
inline |
Return the multiplication of scalar and matrix (m * s).
m | the matrix |
s | the scalar |
Definition at line 511 of file DynamicRectMatrix.hpp.
References OpenMD::DynamicRectMatrix< Real >::getNCol(), OpenMD::DynamicRectMatrix< Real >::getNRow(), and OpenMD::DynamicRectMatrix< Real >::mul().
|
inline |
Return the multiplication of a scalar and a matrix (s * m).
s | the scalar |
m | the matrix |
Definition at line 526 of file DynamicRectMatrix.hpp.
References OpenMD::DynamicRectMatrix< Real >::getNCol(), OpenMD::DynamicRectMatrix< Real >::getNRow(), and OpenMD::DynamicRectMatrix< Real >::mul().
|
inline |
Returns the multiplication of a matrix and a vector (m * v).
m | the matrix |
v | the vector |
Definition at line 526 of file RectMatrix.hpp.
References Row.
|
inline |
Return the multiplication of two matrixes (m1 * m2).
m1 | the first matrix |
m2 | the second matrix |
Definition at line 541 of file DynamicRectMatrix.hpp.
References OpenMD::DynamicRectMatrix< Real >::getNCol(), and OpenMD::DynamicRectMatrix< Real >::getNRow().
|
inline |
Returns the multiplication of a vector transpose and a matrix (v^T * m).
v | the vector |
m | the matrix |
Definition at line 543 of file RectMatrix.hpp.
References Row.
Polynomial<Real> OpenMD::operator * | ( | const Polynomial< Real > & | p1, |
const Polynomial< Real > & | p2 | ||
) |
Generates and returns the product of two given Polynomials.
Definition at line 543 of file Polynomial.hpp.
References OpenMD::Polynomial< Real >::addCoefficient(), OpenMD::Polynomial< Real >::begin(), and OpenMD::Polynomial< Real >::end().
Polynomial<Real> OpenMD::operator * | ( | const Polynomial< Real > & | p, |
const Real | v | ||
) |
Definition at line 558 of file Polynomial.hpp.
References OpenMD::Polynomial< Real >::begin(), OpenMD::Polynomial< Real >::end(), and OpenMD::Polynomial< Real >::setCoefficient().
|
inline |
Return the multiplication of a matrix and a vector (m * v).
m | the matrix |
v | the vector |
Definition at line 562 of file DynamicRectMatrix.hpp.
References OpenMD::DynamicRectMatrix< Real >::getNCol(), and OpenMD::DynamicRectMatrix< Real >::getNRow().
Polynomial<Real> OpenMD::operator * | ( | const Real | v, |
const Polynomial< Real > & | p | ||
) |
Definition at line 570 of file Polynomial.hpp.
References OpenMD::Polynomial< Real >::begin(), OpenMD::Polynomial< Real >::end(), and OpenMD::Polynomial< Real >::setCoefficient().
Quaternion<Real> OpenMD::operator * | ( | const Quaternion< Real > & | q, |
Real | s | ||
) |
Returns the vaule of scalar multiplication of this quaterion q (q * s).
q | the source quaternion |
s | the scalar value |
Definition at line 587 of file Quaternion.hpp.
References OpenMD::Quaternion< Real >::mul().
Quaternion<Real> OpenMD::operator * | ( | const Real & | s, |
const Quaternion< Real > & | q | ||
) |
Returns the vaule of scalar multiplication of this quaterion q (q * s).
s | the scalar value |
q | the source quaternion |
Definition at line 600 of file Quaternion.hpp.
References OpenMD::Quaternion< Real >::mul().
|
inline |
Return the multiplication of two matrixes (m1 * m2).
m1 | the first matrix |
m2 | the second matrix |
Definition at line 610 of file SquareMatrix3.hpp.
|
inline |
Returns the multiplication of two quaternion
q1 | the first quaternion |
q2 | the second quaternion |
Definition at line 613 of file Quaternion.hpp.
|
inline |
unary minus
Definition at line 332 of file DynamicVector.hpp.
References OpenMD::DynamicVector< Real, Alloc >::negate().
DynamicVector<Real> OpenMD::operator - | ( | const DynamicVector< Real > & | v1, |
const DynamicVector< Real > & | v2 | ||
) |
Return the difference of two vectors (v1 - v2).
v1 | the first vector |
v2 | the second vector |
Definition at line 359 of file DynamicVector.hpp.
References OpenMD::DynamicVector< Real, Alloc >::sub().
|
inline |
Negate the value of every element of this matrix.
Definition at line 433 of file RectMatrix.hpp.
References OpenMD::RectMatrix< Real, Row, Col >::negate().
|
inline |
unary minus
Definition at line 436 of file Vector.hpp.
References OpenMD::Vector< Real, Dim >::negate().
|
inline |
Return the difference of two matrixes (m1 - m2).
m1 | the first matrix |
m2 | the second matrix |
Definition at line 463 of file RectMatrix.hpp.
References OpenMD::RectMatrix< Real, Row, Col >::sub().
Vector<Real, Dim> OpenMD::operator - | ( | const Vector< Real, Dim > & | v1, |
const Vector< Real, Dim > & | v2 | ||
) |
Return the difference of two vectors (v1 - v2).
v1 | the first vector |
v2 | the second vector |
Definition at line 463 of file Vector.hpp.
References OpenMD::Vector< Real, Dim >::sub().
|
inline |
Negate the value of every element of this matrix.
Definition at line 465 of file DynamicRectMatrix.hpp.
References OpenMD::DynamicRectMatrix< Real >::negate().
|
inline |
Return the difference of two matrixes (m1 - m2).
m1 | the first matrix |
m2 | the second matrix |
Definition at line 496 of file DynamicRectMatrix.hpp.
References OpenMD::DynamicRectMatrix< Real >::getNCol(), OpenMD::DynamicRectMatrix< Real >::getNRow(), and OpenMD::DynamicRectMatrix< Real >::sub().
Polynomial<Real> OpenMD::operator - | ( | const Polynomial< Real > & | p1, |
const Polynomial< Real > & | p2 | ||
) |
Generates and returns the difference of two given Polynomials.
p1 | the first polynomial |
p2 | the second polynomial |
Definition at line 607 of file Polynomial.hpp.
References OpenMD::Polynomial< Real >::addCoefficient(), OpenMD::Polynomial< Real >::begin(), and OpenMD::Polynomial< Real >::end().
|
inline |
Definition at line 178 of file Predicate.hpp.
|
inline |
Return the sum of two vectors (v1 - v2).
v1 | the first vector |
v2 | the second vector |
Definition at line 345 of file DynamicVector.hpp.
References OpenMD::DynamicVector< Real, Alloc >::add().
|
inline |
Return the sum of two matrixes (m1 + m2).
m1 | the first matrix |
m2 | the second matrix |
Definition at line 448 of file RectMatrix.hpp.
References OpenMD::RectMatrix< Real, Row, Col >::add().
|
inline |
Return the sum of two vectors (v1 - v2).
v1 | the first vector |
v2 | the second vector |
Definition at line 449 of file Vector.hpp.
References OpenMD::Vector< Real, Dim >::add().
|
inline |
Return the sum of two matrixes (m1 + m2).
m1 | the first matrix |
m2 | the second matrix |
Definition at line 480 of file DynamicRectMatrix.hpp.
References OpenMD::DynamicRectMatrix< Real >::add(), OpenMD::DynamicRectMatrix< Real >::getNCol(), and OpenMD::DynamicRectMatrix< Real >::getNRow().
Polynomial<Real> OpenMD::operator+ | ( | const Polynomial< Real > & | p1, |
const Polynomial< Real > & | p2 | ||
) |
Generates and returns the sum of two given Polynomials.
p1 | the first polynomial |
p2 | the second polynomial |
Definition at line 587 of file Polynomial.hpp.
References OpenMD::Polynomial< Real >::addCoefficient(), OpenMD::Polynomial< Real >::begin(), and OpenMD::Polynomial< Real >::end().
OpenMDBitSet OpenMD::operator- | ( | const OpenMDBitSet & | bs1, |
const OpenMDBitSet & | bs2 | ||
) |
Definition at line 187 of file OpenMDBitSet.cpp.
References OpenMD::OpenMDBitSet::size().
SelectionSet OpenMD::operator- | ( | const SelectionSet & | ss1, |
const SelectionSet & | ss2 | ||
) |
Definition at line 210 of file SelectionSet.cpp.
References OpenMD::SelectionSet::bitsets_, and N_SELECTIONTYPES.
SelectionManager OpenMD::operator- | ( | const SelectionManager & | sman1, |
const SelectionManager & | sman2 | ||
) |
Definition at line 629 of file SelectionManager.cpp.
DynamicVector<Real> OpenMD::operator/ | ( | const DynamicVector< Real > & | v1, |
Real | s | ||
) |
Returns the value of division of a vector by a scalar.
v1 | the source vector |
s | the scalar value |
Definition at line 399 of file DynamicVector.hpp.
References OpenMD::DynamicVector< Real, Alloc >::div().
Vector<Real, Dim> OpenMD::operator/ | ( | const Vector< Real, Dim > & | v1, |
Real | s | ||
) |
Returns the value of division of a vector by a scalar.
v1 | the source vector |
s | the scalar value |
Definition at line 502 of file Vector.hpp.
References OpenMD::Vector< Real, Dim >::div().
|
inline |
Return the scalar division of matrix (m / s).
m | the matrix |
s | the scalar |
Definition at line 560 of file RectMatrix.hpp.
References OpenMD::RectMatrix< Real, Row, Col >::div().
|
inline |
Return the scalar division of matrix (m / s).
m | the matrix |
s | the scalar |
Definition at line 582 of file DynamicRectMatrix.hpp.
References OpenMD::DynamicRectMatrix< Real >::div(), OpenMD::DynamicRectMatrix< Real >::getNCol(), and OpenMD::DynamicRectMatrix< Real >::getNRow().
|
inline |
Returns the division of two quaternion
q1 | divisor |
q2 | dividen |
Definition at line 626 of file Quaternion.hpp.
References OpenMD::Quaternion< Real >::inverse().
Quaternion<Real> OpenMD::operator/ | ( | const Real & | s, |
Quaternion< Real > & | q | ||
) |
Returns the value of the division of a scalar by a quaternion
s | scalar |
q | quaternion |
Definition at line 638 of file Quaternion.hpp.
References OpenMD::Quaternion< Real >::inverse().
|
inline |
Definition at line 88 of file Tuple.hpp.
References OpenMD::tuple3< T1, T2, T3 >::first, OpenMD::tuple3< T1, T2, T3 >::second, and OpenMD::tuple3< T1, T2, T3 >::third.
|
inline |
Definition at line 115 of file Tuple.hpp.
References OpenMD::tuple4< T1, T2, T3, T4 >::first, OpenMD::tuple4< T1, T2, T3, T4 >::fourth, OpenMD::tuple4< T1, T2, T3, T4 >::second, and OpenMD::tuple4< T1, T2, T3, T4 >::third.
std::ostream& OpenMD::operator<< | ( | std::ostream & | os, |
AmberImproperTorsionType & | ott | ||
) |
Definition at line 86 of file AmberImproperTorsionType.hpp.
References OpenMD::AmberImproperTorsionType::v2_.
std::ostream& OpenMD::operator<< | ( | std::ostream & | os, |
HarmonicInversionType & | hit | ||
) |
Definition at line 86 of file HarmonicInversionType.hpp.
References OpenMD::HarmonicInversionType::d0_, and OpenMD::HarmonicInversionType::phi0_.
std::ostream & OpenMD::operator<< | ( | std::ostream & | o, |
IntegratorFactory & | factory | ||
) |
write out all of the type identifiers to an output stream
Definition at line 89 of file IntegratorFactory.cpp.
References OpenMD::IntegratorFactory::getIdents().
std::ostream & OpenMD::operator<< | ( | std::ostream & | o, |
OptimizationFactory & | factory | ||
) |
write out all of the type identifiers to an output stream
Definition at line 89 of file OptimizationFactory.cpp.
References OpenMD::OptimizationFactory::getIdents().
std::ostream & OpenMD::operator<< | ( | std::ostream & | o, |
LatticeFactory & | factory | ||
) |
write out all of the type identifiers to an output stream
Definition at line 89 of file LatticeFactory.cpp.
References OpenMD::LatticeFactory::getIdents().
std::ostream & OpenMD::operator<< | ( | std::ostream & | o, |
HydrodynamicsModelFactory & | factory | ||
) |
write out all of the type identifiers to an output stream
Definition at line 90 of file HydrodynamicsModelFactory.cpp.
References OpenMD::HydrodynamicsModelFactory::getIdents().
std::ostream& OpenMD::operator<< | ( | std::ostream & | os, |
PolynomialBendType & | pbt | ||
) |
Definition at line 90 of file PolynomialBendType.hpp.
References OpenMD::Polynomial< Real >::begin(), OpenMD::Polynomial< Real >::end(), OpenMD::BendType::getTheta(), and OpenMD::PolynomialBendType::polynomial_.
std::ostream& OpenMD::operator<< | ( | std::ostream & | os, |
PolynomialBondType & | pbt | ||
) |
Definition at line 90 of file PolynomialBondType.hpp.
References OpenMD::Polynomial< Real >::begin(), OpenMD::Polynomial< Real >::end(), OpenMD::BondType::getEquilibriumBondLength(), and OpenMD::PolynomialBondType::polynomial_.
std::ostream& OpenMD::operator<< | ( | std::ostream & | os, |
HarmonicTorsionType & | htt | ||
) |
Definition at line 111 of file HarmonicTorsionType.hpp.
References OpenMD::HarmonicTorsionType::d0_, and OpenMD::HarmonicTorsionType::phi0_.
std::ostream& OpenMD::operator<< | ( | std::ostream & | os, |
OplsTorsionType & | ott | ||
) |
Definition at line 113 of file OplsTorsionType.hpp.
References OpenMD::OplsTorsionType::v1_, OpenMD::OplsTorsionType::v2_, and OpenMD::OplsTorsionType::v3_.
std::ostream& OpenMD::operator<< | ( | std::ostream & | os, |
TrappeTorsionType & | ttt | ||
) |
Definition at line 123 of file TrappeTorsionType.hpp.
References OpenMD::TrappeTorsionType::c0_, OpenMD::TrappeTorsionType::c1_, OpenMD::TrappeTorsionType::c2_, and OpenMD::TrappeTorsionType::c3_.
std::ostream& OpenMD::operator<< | ( | std::ostream & | o, |
PairList & | e | ||
) |
write out the exclusion list to an ostream
Definition at line 153 of file PairList.cpp.
References OpenMD::PairList::pairSet_.
std::ostream& OpenMD::operator<< | ( | std::ostream & | o, |
GenericFactory< O, I, C > & | factory | ||
) |
write out all of the type identifiers to an output stream
Definition at line 220 of file GenericFactory.hpp.
References OpenMD::GenericFactory< Object, IdentType, Creator >::getIdents().
std::ostream& OpenMD::operator<< | ( | std::ostream & | os, |
const OpenMDBitSet & | bs | ||
) |
Definition at line 233 of file OpenMDBitSet.cpp.
References OpenMD::OpenMDBitSet::bitset_.
std::ostream& OpenMD::operator<< | ( | std::ostream & | os, |
const SelectionSet & | ss | ||
) |
Definition at line 243 of file SelectionSet.cpp.
References OpenMD::SelectionSet::bitsets_, and N_SELECTIONTYPES.
std::ostream& OpenMD::operator<< | ( | std::ostream & | o, |
Molecule & | mol | ||
) |
Definition at line 413 of file Molecule.cpp.
References OpenMD::Molecule::getGlobalIndex(), OpenMD::Molecule::getNAtoms(), OpenMD::Molecule::getNBends(), OpenMD::Molecule::getNBonds(), OpenMD::Molecule::getNConstraintPairs(), OpenMD::Molecule::getNCutoffGroups(), OpenMD::Molecule::getNFluctuatingCharges(), OpenMD::Molecule::getNIntegrableObjects(), OpenMD::Molecule::getNInversions(), OpenMD::Molecule::getNRigidBodies(), and OpenMD::Molecule::getNTorsions().
|
inline |
Definition at line 425 of file MersenneTwister.hpp.
References OpenMD::MTRand::left, OpenMD::MTRand::N, and OpenMD::MTRand::state.
std::ostream& OpenMD::operator<< | ( | std::ostream & | o, |
const DynamicVector< Real > & | v | ||
) |
Write to an output stream
Definition at line 450 of file DynamicVector.hpp.
std::ostream& OpenMD::operator<< | ( | std::ostream & | o, |
const Vector< Real, Dim > & | v | ||
) |
Write to an output stream
Definition at line 577 of file Vector.hpp.
std::ostream& OpenMD::operator<< | ( | std::ostream & | o, |
const DynamicRectMatrix< Real > & | m | ||
) |
Write to an output stream
Definition at line 594 of file DynamicRectMatrix.hpp.
References OpenMD::DynamicRectMatrix< Real >::getNCol(), and OpenMD::DynamicRectMatrix< Real >::getNRow().
std::ostream& OpenMD::operator<< | ( | std::ostream & | o, |
const RectMatrix< Real, Row, Col > & | m | ||
) |
ostream& OpenMD::operator<< | ( | ostream & | o, |
SimInfo & | info | ||
) |
Definition at line 1060 of file SimInfo.cpp.
bool OpenMD::operator== | ( | const OpenMDBitSet & | bs1, |
const OpenMDBitSet & | bs2 | ||
) |
Definition at line 195 of file OpenMDBitSet.cpp.
References OpenMD::OpenMDBitSet::bitset_, equal(), and OpenMD::OpenMDBitSet::size().
bool OpenMD::operator== | ( | const SelectionSet & | ss1, |
const SelectionSet & | ss2 | ||
) |
Definition at line 217 of file SelectionSet.cpp.
References OpenMD::SelectionSet::bitsets_, and N_SELECTIONTYPES.
|
inline |
Definition at line 647 of file Quaternion.hpp.
References equal().
|
inline |
Definition at line 434 of file MersenneTwister.hpp.
References OpenMD::MTRand::left, OpenMD::MTRand::N, OpenMD::MTRand::pNext, and OpenMD::MTRand::state.
OpenMDBitSet OpenMD::operator^ | ( | const OpenMDBitSet & | bs1, |
const OpenMDBitSet & | bs2 | ||
) |
Definition at line 179 of file OpenMDBitSet.cpp.
References OpenMD::OpenMDBitSet::size().
SelectionSet OpenMD::operator^ | ( | const SelectionSet & | ss1, |
const SelectionSet & | ss2 | ||
) |
Definition at line 203 of file SelectionSet.cpp.
References OpenMD::SelectionSet::bitsets_, and N_SELECTIONTYPES.
SelectionManager OpenMD::operator^ | ( | const SelectionManager & | sman1, |
const SelectionManager & | sman2 | ||
) |
Definition at line 622 of file SelectionManager.cpp.
OpenMDBitSet OpenMD::operator| | ( | const OpenMDBitSet & | bs1, |
const OpenMDBitSet & | bs2 | ||
) |
Definition at line 163 of file OpenMDBitSet.cpp.
References OpenMD::OpenMDBitSet::size().
SelectionSet OpenMD::operator| | ( | const SelectionSet & | ss1, |
const SelectionSet & | ss2 | ||
) |
Definition at line 189 of file SelectionSet.cpp.
References OpenMD::SelectionSet::bitsets_, and N_SELECTIONTYPES.
SelectionManager OpenMD::operator| | ( | const SelectionManager & | sman1, |
const SelectionManager & | sman2 | ||
) |
Definition at line 608 of file SelectionManager.cpp.
|
inline |
Definition at line 167 of file Predicate.hpp.
|
inline |
Definition at line 329 of file ParamConstraint.hpp.
|
inline |
Definition at line 622 of file SquareMatrix3.hpp.
Referenced by OpenMD::ForceAutoCorrFunc::calcCorrVal(), OpenMD::TorqueAutoCorrFunc::calcCorrVal(), OpenMD::TorForCorrFunc::calcCorrVal(), OpenMD::ForTorCorrFunc::calcCorrVal(), OpenMD::MolecularRestraint::calcForce(), OpenMD::ApproximationModel::calcHydroPropsAtCD(), OpenMD::ApproximationModel::calcHydroPropsAtCR(), OpenMD::RigidBody::calcRefCoords(), OpenMD::RMSD::calculate_rmsd(), OpenMD::RNEMD::collectData(), OpenMD::ActionCorrFunc::correlateFrames(), OpenMD::StressCorrFunc::correlateFrames(), OpenMD::MomentumCorrFunc::correlateFrames(), OpenMD::RNEMD::doVSS(), OpenMD::NanoLength::getLength(), OpenMD::LDForceManager::getMomentData(), OpenMD::BAOAB::getMomentData(), OpenMD::Thermo::getPressureTensor(), OpenMD::Thermo::getSystemQuadrupole(), OpenMD::Triangle::hydro_tensor(), OpenMD::ForceManager::longRangeInteractions(), OpenMD::RMSD::optimal_superposition(), OpenMD::ForceAutoCorrFunc::postCorrelate(), OpenMD::TorqueAutoCorrFunc::postCorrelate(), OpenMD::TorForCorrFunc::postCorrelate(), OpenMD::ForTorCorrFunc::postCorrelate(), OpenMD::RippleOP::process(), OpenMD::P2OrderParameter::process(), OpenMD::RNEMDR::processFrame(), OpenMD::RNEMDRTheta::processFrame(), OpenMD::Electrostatic::ReciprocalSpaceSum(), and OpenMD::ForceManager::selectedLongRangeInteractions().
bool OpenMD::pairCompare | ( | const pair< RealType, int > & | l, |
const pair< RealType, int > & | r | ||
) |
Definition at line 53 of file Cuboctahedron.cpp.
Referenced by OpenMD::Cuboctahedron::getPoints().
void OpenMD::registerAll | ( | ) |
register force fields, integrators and optimizers
Definition at line 112 of file Register.cpp.
References registerIntegrators(), and registerOptimizers().
Referenced by main().
void OpenMD::registerIntegrators | ( | ) |
void OpenMD::registerLattice | ( | ) |
Register all lattice
Definition at line 106 of file Register.cpp.
Referenced by main(), and OpenMD::shapedLattice::shapedLattice().
void OpenMD::registerOptimizers | ( | ) |
bool OpenMD::replaceWithWildCard | ( | vector< vector< string >::iterator > & | cont, |
vector< string > & | sequence, | ||
vector< string > & | result, | ||
const string & | wildCard = "X" |
||
) |
iteratively replace the sequence with wild cards
cont | iterator container, if expect the whole series of combinations, pass an empty iterator container. The user should not modify this iterator container |
sequence | the whole sequence used to generate combination |
result | a possible combination sequence which is set on return |
wildCard | the wild card string. Its value is "X" by default |
Definition at line 48 of file Utility.cpp.
References next_combination().
Referenced by OpenMD::TypeContainer< BendType, 3 >::find().
|
inline |
Definition at line 52 of file Utility.hpp.
Referenced by OpenMD::ForceMatrixDecomposition::buildNeighborList(), OpenMD::Hxy::process(), OpenMD::Field< Vector3d >::processFrame(), and OpenMD::Snapshot::wrapVector().
void OpenMD::swap_if | ( | Cont & | b1, |
Cont & | b2, | ||
Predict | predict | ||
) |
Definition at line 65 of file CompositeShape.cpp.
Referenced by OpenMD::CompositeShape::getBoundingBox().
std::string OpenMD::to_string | ( | T | value | ) |
Definition at line 150 of file StringUtils.hpp.
Referenced by OpenMD::BendTypeParser::parseTypeAndPars(), OpenMD::BondTypeParser::parseTypeAndPars(), OpenMD::InversionTypeParser::parseTypeAndPars(), and OpenMD::TorsionTypeParser::parseTypeAndPars().
void OpenMD::toLower | ( | Container & | cont, |
const std::locale & | loc = std::locale() |
||
) |
Definition at line 74 of file CaseConversion.hpp.
Referenced by ParameterTraits< bool >::convert(), main(), OpenMD::ShapeAtomTypesSectionParser::parseShapeFile(), and toLowerCopy().
Container OpenMD::toLowerCopy | ( | const Container & | cont, |
const std::locale & | loc = std::locale() |
||
) |
Definition at line 79 of file CaseConversion.hpp.
References toLower().
std::string OpenMD::toString | ( | const T & | v | ) |
Convert a variable to a string
T | data type |
v | data to be converted |
Definition at line 119 of file StringUtils.hpp.
Referenced by OpenMD::SelectionEvaluator::invalidIndex(), and OpenMD::SelectionEvaluator::invalidIndexRange().
void OpenMD::toUpper | ( | Container & | cont, |
const std::locale & | loc = std::locale() |
||
) |
Definition at line 86 of file CaseConversion.hpp.
Referenced by OpenMD::GB::addType(), OpenMD::SimCreator::createSim(), OpenMD::Cuboctahedron::Cuboctahedron(), OpenMD::SC::getAlpha(), OpenMD::LJ::getSigma(), OpenMD::Electrostatic::initialize(), OpenMD::EAM::initialize(), OpenMD::EAMAdapter::makeZhou2001(), OpenMD::EAMAdapter::makeZhou2004(), OpenMD::EAMAdapter::makeZhou2005(), OpenMD::TorsionTypesSectionParser::parseLine(), OpenMD::EAMAtomTypesSectionParser::parseLine(), OpenMD::NonBondedInteractionsSectionParser::parseLine(), OpenMD::RNEMD::parseOutputFileFormat(), OpenMD::Stats::parseStatFileFormat(), OpenMD::ForceManager::setupCutoffs(), and toUpperCopy().
Container OpenMD::toUpperCopy | ( | const Container & | cont, |
const std::locale & | loc = std::locale() |
||
) |
Definition at line 91 of file CaseConversion.hpp.
References toUpper().
Referenced by OpenMD::ForceManager::initialize(), OpenMD::Integrator::Integrator(), main(), OpenMD::EqualIgnoreCaseConstraint::operator()(), OpenMD::ContainsConstraint::operator()(), OpenMD::RestraintForceManager::RestraintForceManager(), and OpenMD::ForceManager::setupCutoffs().
void OpenMD::trim | ( | std::string & | str | ) |
std::string OpenMD::trimCopy | ( | const std::string & | input | ) |
Remove all leading and trailing spaces from the input.
input | An input sequence |
Definition at line 66 of file Trim.cpp.
References isSpace(), and trimCopyIf().
std::string OpenMD::trimCopyIf | ( | const std::string & | input, |
P | pred | ||
) |
Remove all leading and trailing spaces from the input. The supplied predicate is used to determine which characters are considered spaces
input | An input sequence |
pred | The unary predicate identifying spaces |
Definition at line 149 of file Trim.hpp.
References trimIf().
Referenced by trimCopy().
void OpenMD::trimIf | ( | std::string & | str, |
P | pred | ||
) |
Remove all leading and trailing spaces in-place. The supplied predicate is used to determine which characters are considered spaces
str | An input sequence |
pred | The unary predicate identifying spaces |
Definition at line 108 of file Trim.hpp.
References trimLeftIf(), and trimRightIf().
Referenced by trim(), and trimCopyIf().
void OpenMD::trimLeft | ( | std::string & | str | ) |
Remove all leading spaces in-place.
str | An input sequence |
Definition at line 46 of file Trim.cpp.
References isSpace(), and trimLeftIf().
std::string OpenMD::trimLeftCopy | ( | const std::string & | input | ) |
Remove all leading spaces from the input.
input | An input sequence |
Definition at line 58 of file Trim.cpp.
References isSpace(), and trimLeftCopyIf().
Referenced by OpenMD::SimCreator::createSim(), OpenMD::SimplePreprocessor::doPreprocess(), OpenMD::SectionParser::parse(), and OpenMD::SectionParserManager::parse().
std::string OpenMD::trimLeftCopyIf | ( | const std::string & | input, |
P | pred | ||
) |
Remove all leading spaces from the input. The supplied predicate is used to determine which characters are considered spaces
input | An input sequence |
pred | The unary predicate identifying spaces |
Definition at line 121 of file Trim.hpp.
References trimLeftIf().
Referenced by trimLeftCopy().
void OpenMD::trimLeftIf | ( | std::string & | str, |
P | pred | ||
) |
Remove all leading spaces in-place. The supplied predicate is used to determine which characters are considered spaces
str | An input sequence |
pred | The unary predicate identifying spaces |
Definition at line 69 of file Trim.hpp.
Referenced by trimIf(), trimLeft(), and trimLeftCopyIf().
void OpenMD::trimRight | ( | std::string & | str | ) |
Remove all trailing spaces in-place.
str | An input sequence |
Definition at line 50 of file Trim.cpp.
References isSpace(), and trimRightIf().
std::string OpenMD::trimRightCopy | ( | const std::string & | input | ) |
Remove all trailing spaces from the input.
input | An input sequence |
Definition at line 62 of file Trim.cpp.
References isSpace(), and trimRightCopyIf().
std::string OpenMD::trimRightCopyIf | ( | const std::string & | input, |
P | pred | ||
) |
Remove all trailing spaces from the input. The supplied predicate is used to determine which characters are considered spaces
input | An input sequence |
pred | The unary predicate identifying spaces |
Definition at line 135 of file Trim.hpp.
References trimRightIf().
Referenced by trimRightCopy().
void OpenMD::trimRightIf | ( | std::string & | str, |
P | pred | ||
) |
Remove all trailing spaces in-place. The supplied predicate is used to determine which characters are considered spaces
str | An input sequence |
pred | The unary predicate identifying spaces |
Definition at line 88 of file Trim.hpp.
Referenced by trimIf(), trimRight(), and trimRightCopyIf().
char * OpenMD::trimSpaces | ( | char * | str | ) |
Removes left and right spaces from a string
str | String to trim |
Definition at line 79 of file StringUtils.cpp.
std::string OpenMD::UpperCase | ( | const std::string & | S | ) |
Converts a string to UPPER CASE
S |
Definition at line 59 of file StringUtils.cpp.
Referenced by OpenMD::Molecule::complete().
|
inline |
Returns the linear indexing for integer vectors. Compare to Rapaport's VLinear
p | first vector |
s | second vector |
Definition at line 154 of file Vector3.hpp.
References OpenMD::Vector3< Real >::x(), OpenMD::Vector3< Real >::y(), and OpenMD::Vector3< Real >::z().
Referenced by OpenMD::ForceMatrixDecomposition::buildNeighborList().
|
static |
Definition at line 84 of file NonBondedInteraction.hpp.
Referenced by OpenMD::Buckingham::getHash().
string const OpenMD::DirectionalTypeID = "Directional" |
Definition at line 53 of file DirectionalAdapter.hpp.
Referenced by OpenMD::DirectionalAdapter::getDirectionalParam(), OpenMD::DirectionalAdapter::isDirectional(), and OpenMD::DirectionalAdapter::makeDirectional().
|
static |
Definition at line 76 of file NonBondedInteraction.hpp.
Referenced by OpenMD::InteractionManager::doPair(), OpenMD::InteractionManager::doPreForce(), OpenMD::InteractionManager::doPrePair(), OpenMD::EAM::getHash(), and OpenMD::InteractionManager::initialize().
string const OpenMD::EAMtypeID = "EAM" |
Definition at line 52 of file EAMAdapter.hpp.
Referenced by OpenMD::EAMAdapter::getEAMParam(), OpenMD::EAMAdapter::isEAM(), OpenMD::EAMAdapter::makeFuncfl(), OpenMD::EAMAdapter::makeZhou2001(), OpenMD::EAMAdapter::makeZhou2004(), OpenMD::EAMAdapter::makeZhou2005(), OpenMD::EAMAdapter::makeZhou2005Oxygen(), and OpenMD::EAMAdapter::makeZhouRose().
|
static |
Boolean flags for the iHash_ and sHash_ data structures. These are used to greatly increase the speed of looking up the low-level interaction for any given pair:
Definition at line 74 of file NonBondedInteraction.hpp.
Referenced by OpenMD::InteractionManager::doPair(), OpenMD::InteractionManager::doSelfCorrection(), OpenMD::ElectrostaticInteraction::getHash(), and OpenMD::InteractionManager::initialize().
|
static |
Definition at line 59 of file Vector.hpp.
Referenced by OpenMD::LJ::addExplicitInteraction(), OpenMD::RepulsivePower::addExplicitInteraction(), OpenMD::Mie::addExplicitInteraction(), OpenMD::Buckingham::addExplicitInteraction(), OpenMD::SC::addExplicitInteraction(), OpenMD::GhostTorsion::calcForce(), OpenMD::LJ::calcForce(), OpenMD::RepulsivePower::calcForce(), OpenMD::Mie::calcForce(), OpenMD::Buckingham::calcForce(), OpenMD::Torsion::calcForce(), OpenMD::RigidBody::calcRefCoords(), CholeskyDecomposition(), OpenMD::DirectionalAtom::DirectionalAtom(), equal(), getSTOZeta(), OpenMD::Torsion::getValue(), OpenMD::LJ::initialize(), OpenMD::RepulsivePower::initialize(), OpenMD::Mie::initialize(), OpenMD::SquareMatrix3< RealType >::inverse(), OpenMD::SquareMatrix< RealType, 6 >::isDiagonal(), OpenMD::SquareMatrix< RealType, 6 >::isSymmetric(), OpenMD::SquareMatrix< RealType, 6 >::isUnitMatrix(), LUFactorLinearSystem(), main(), OpenMD::LennardJonesAdapter::makeLennardJones(), OpenMD::SuttonChenAdapter::makeSuttonChen(), QuantLib::LineSearchBasedMethod::minimize(), OpenMD::DumpReader::parseDumpLine(), OpenMD::RestReader::parseDumpLine(), OpenMD::LennardJonesAtomTypesSectionParser::parseLine(), OpenMD::BendTypeParser::parseLine(), OpenMD::BondTypeParser::parseLine(), OpenMD::SCAtomTypesSectionParser::parseLine(), OpenMD::NonBondedInteractionsSectionParser::parseLine(), OpenMD::Integrator::preStep(), sSTOCoulInt(), sSTOOvInt(), and OpenMD::SquareMatrix3< RealType >::toQuaternion().
ElementsTable OpenMD::etab |
Definition at line 55 of file ElementsTable.cpp.
Referenced by OpenMD::BAOAB::BAOAB(), OpenMD::BeadModel::createSingleBead(), OpenMD::ShapeBuilder::internalCreateShape(), and OpenMD::LDForceManager::LDForceManager().
string const OpenMD::FCtypeID = "Charge" |
Definition at line 52 of file FixedChargeAdapter.hpp.
Referenced by OpenMD::FixedChargeAdapter::getFixedChargeParam(), OpenMD::FixedChargeAdapter::isFixedCharge(), and OpenMD::FixedChargeAdapter::makeFixedCharge().
string const OpenMD::FQtypeID = "FlucQ" |
Definition at line 54 of file FluctuatingChargeAdapter.hpp.
Referenced by OpenMD::FluctuatingChargeAdapter::getFluctuatingChargeParam(), OpenMD::FluctuatingChargeAdapter::isFluctuatingCharge(), and OpenMD::FluctuatingChargeAdapter::makeFluctuatingCharge().
string const OpenMD::FuncflTypeID = "FUNCFL" |
Definition at line 53 of file EAMAdapter.hpp.
Referenced by OpenMD::EAMAdapter::getFuncflParam(), and OpenMD::EAMAdapter::makeFuncfl().
|
static |
Definition at line 79 of file NonBondedInteraction.hpp.
Referenced by OpenMD::InteractionManager::doPair(), OpenMD::GB::getHash(), and OpenMD::InteractionManager::initialize().
string const OpenMD::GBtypeID = "GB" |
Definition at line 52 of file GayBerneAdapter.hpp.
Referenced by OpenMD::GayBerneAdapter::getGayBerneParam(), OpenMD::GayBerneAdapter::isGayBerne(), and OpenMD::GayBerneAdapter::makeGayBerne().
|
static |
Definition at line 85 of file NonBondedInteraction.hpp.
Referenced by OpenMD::InteractionManager::doPair(), OpenMD::InversePowerSeries::getHash(), and OpenMD::InteractionManager::initialize().
|
static |
Definition at line 75 of file NonBondedInteraction.hpp.
Referenced by OpenMD::InteractionManager::doPair(), OpenMD::LJ::getHash(), and OpenMD::InteractionManager::initialize().
string const OpenMD::LJtypeID = "LJ" |
Definition at line 52 of file LennardJonesAdapter.hpp.
Referenced by OpenMD::LennardJonesAdapter::getLJParam(), OpenMD::LennardJonesAdapter::isLennardJones(), and OpenMD::LennardJonesAdapter::makeLennardJones().
const Mat3x3d OpenMD::M3Zero(0.0) |
|
static |
Definition at line 82 of file NonBondedInteraction.hpp.
Referenced by OpenMD::InteractionManager::doPair(), OpenMD::MAW::getHash(), and OpenMD::InteractionManager::initialize().
|
static |
Definition at line 83 of file NonBondedInteraction.hpp.
Referenced by OpenMD::InteractionManager::doPair(), OpenMD::Mie::getHash(), and OpenMD::InteractionManager::initialize().
|
static |
Definition at line 80 of file NonBondedInteraction.hpp.
Referenced by OpenMD::InteractionManager::doPair(), OpenMD::Morse::getHash(), and OpenMD::InteractionManager::initialize().
string const OpenMD::MultipoleTypeID = "Multipole" |
Definition at line 53 of file MultipoleAdapter.hpp.
Referenced by OpenMD::MultipoleAdapter::getMultipoleParam(), OpenMD::MultipoleAdapter::isMultipole(), and OpenMD::MultipoleAdapter::makeMultipole().
string const OpenMD::PolTypeID = "Polarizable" |
Definition at line 52 of file PolarizableAdapter.hpp.
Referenced by OpenMD::PolarizableAdapter::getPolarizableParam(), OpenMD::PolarizableAdapter::isPolarizable(), and OpenMD::PolarizableAdapter::makePolarizable().
const char* OpenMD::progressSpinner_ = "|/-\\" |
Definition at line 68 of file ProgressBar.cpp.
Referenced by OpenMD::ProgressBar::update().
|
static |
Definition at line 81 of file NonBondedInteraction.hpp.
Referenced by OpenMD::InteractionManager::doPair(), OpenMD::RepulsivePower::getHash(), and OpenMD::InteractionManager::initialize().
|
static |
Definition at line 77 of file NonBondedInteraction.hpp.
Referenced by OpenMD::InteractionManager::doPair(), OpenMD::InteractionManager::doPreForce(), OpenMD::InteractionManager::doPrePair(), OpenMD::SC::getHash(), and OpenMD::InteractionManager::initialize().
string const OpenMD::SCtypeID = "SC" |
Definition at line 52 of file SuttonChenAdapter.hpp.
Referenced by OpenMD::SuttonChenAdapter::getSuttonChenParam(), OpenMD::SuttonChenAdapter::isSuttonChen(), and OpenMD::SuttonChenAdapter::makeSuttonChen().
|
static |
Definition at line 78 of file NonBondedInteraction.hpp.
Referenced by OpenMD::InteractionManager::doPair(), OpenMD::Sticky::getHash(), and OpenMD::InteractionManager::initialize().
string const OpenMD::StickyTypeID = "Sticky" |
Definition at line 52 of file StickyAdapter.hpp.
Referenced by OpenMD::StickyAdapter::getStickyParam(), OpenMD::StickyAdapter::isSticky(), and OpenMD::StickyAdapter::makeSticky().
const Vector3d OpenMD::V3X(1.0, 0.0, 0.0) |
const Vector3d OpenMD::V3Y(0.0, 1.0, 0.0) |
const Vector3d OpenMD::V3Z(0.0, 0.0, 1.0) |
Referenced by OpenMD::MolecularRestraint::calcForce(), OpenMD::ObjectRestraint::calcForce(), OpenMD::GofRAngle2::collectHistogram(), OpenMD::GofAngle2::collectHistogram(), OpenMD::GofRTheta::evaluateAngle(), OpenMD::GofROmega::evaluateAngle(), OpenMD::Quaternion< Real >::getTwistSwingAxisAngle(), OpenMD::GofXyz::initializeHistogram(), OpenMD::MultipoleAtomTypesSectionParser::parseLine(), OpenMD::RippleOP::process(), OpenMD::P2OrderParameter::process(), OpenMD::AngleR::process(), OpenMD::pAngle::process(), OpenMD::Quaternion< Real >::RemoveTwist(), OpenMD::RNEMDRTheta::RNEMDRTheta(), OpenMD::LipidTransVisitor::update(), and OpenMD::DefaultAtomVisitor::visit().
const Vector3d OpenMD::V3Zero(0.0, 0.0, 0.0) |
Referenced by OpenMD::BAOAB::BAOAB(), OpenMD::SHAPES::calcForce(), OpenMD::GB::calcForce(), OpenMD::Sticky::calcForce(), OpenMD::RMSD::calculate_rmsd(), OpenMD::VectorAccumulator::clear(), OpenMD::Snapshot::clearDerivedProperties(), OpenMD::RNEMD::collectData(), OpenMD::ForceMatrixDecomposition::collectData(), OpenMD::ForceMatrixDecomposition::collectIntermediateData(), OpenMD::ConvexHull::computeHull(), OpenMD::AlphaHull::computeHull(), OpenMD::COMVel::doFrame(), OpenMD::CenterOfMass::doFrame(), OpenMD::ContactAngle1::doFrame(), OpenMD::RNEMD::doVSS(), OpenMD::RNEMD::doVSSCurrent(), OpenMD::ForceAutoCorrFunc::ForceAutoCorrFunc(), OpenMD::ForTorCorrFunc::ForTorCorrFunc(), OpenMD::StuntDouble::freeze(), OpenMD::LangevinHullForceManager::genTriangleForces(), OpenMD::Thermo::getHeatFlux(), OpenMD::Sphere::getHydroProp(), OpenMD::CompositeShape::getHydroProp(), OpenMD::Ellipsoid::getHydroProp(), OpenMD::Lattice::Lattice(), OpenMD::LDForceManager::LDForceManager(), main(), OpenMD::NitrileFrequencyMap::NitrileFrequencyMap(), OpenMD::RMSD::optimal_superposition(), OpenMD::MultipoleAtomTypesSectionParser::parseLine(), OpenMD::MoLocator::placeMol(), OpenMD::LangevinHullForceManager::postCalculation(), OpenMD::MultipoleSum::process(), OpenMD::RNEMDZ::processFrame(), OpenMD::RNEMDR::processFrame(), OpenMD::Electrostatic::ReciprocalSpaceSum(), OpenMD::RMSD::set_reference_structure(), OpenMD::ShellStatistics::ShellStatistics(), OpenMD::TorForCorrFunc::TorForCorrFunc(), OpenMD::TorqueAutoCorrFunc::TorqueAutoCorrFunc(), OpenMD::CutoffGroup::updateCOM(), OpenMD::ReplacementVisitor::visit(), OpenMD::DefaultAtomVisitor::visit(), OpenMD::RNEMD::writeOutputFile(), OpenMD::StuntDouble::zeroForcesAndTorques(), and OpenMD::ForceMatrixDecomposition::zeroWorkArrays().
string const OpenMD::ZhouRoseTypeID = "ZHOUROSE" |
Definition at line 55 of file EAMAdapter.hpp.
string const OpenMD::ZhouTypeID = "ZHOU" |
Definition at line 54 of file EAMAdapter.hpp.
Referenced by OpenMD::EAMAdapter::getZhouParam(), OpenMD::EAMAdapter::makeZhou2001(), OpenMD::EAMAdapter::makeZhou2004(), OpenMD::EAMAdapter::makeZhou2005(), OpenMD::EAMAdapter::makeZhou2005Oxygen(), and OpenMD::EAMAdapter::makeZhouRose().