Scine::Swoose  2.1.0
This is the SCINE module Swoose.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
SfamParameters.h
Go to the documentation of this file.
1 
8 #ifndef MOLECULARMECHANICS_SFAMPARAMETERS_H
9 #define MOLECULARMECHANICS_SFAMPARAMETERS_H
10 
11 #include "../MMParameters.h"
12 #include "../Parameters/AngleParameters.h"
13 #include "../Parameters/BondParameters.h"
14 #include "../Parameters/DihedralParameters.h"
15 #include "../Parameters/ImproperDihedralParameters.h"
16 #include "../Topology/AngleType.h"
17 #include "../Topology/BondType.h"
18 #include "../Topology/DihedralType.h"
19 #include "../Topology/ImproperDihedralType.h"
20 #include <Eigen/Core>
21 #include <map>
22 #include <vector>
23 
24 namespace Scine {
25 
26 namespace Utils {
27 class AtomCollection;
28 } // namespace Utils
29 
30 namespace MolecularMechanics {
31 class AtomTypesHolder;
37 class SfamParameters final : public MMParameters {
38  public:
40  bool sanityCheck(const AtomTypesHolder& atomTypes) const;
42  std::vector<double> getChargesForEachAtom(const AtomTypesHolder& atomTypes) const;
47  void prepareC6Matrix(const AtomTypesHolder& atomTypes);
49  float getC6(int indexOfAtomTypeA, int indexOfAtomTypeB) const;
51  float getC6(const std::string& a, const std::string& b) const;
53  std::vector<double> getNonCovalentParameters() const;
55  void setC6(int indexOfAtomTypeA, int indexOfAtomTypeB, float c6);
57  void setC6(const std::string& a, const std::string& b, float c6);
59  void resetC6Matrix();
61  void setNonCovalentParameters(std::vector<double> nonCovalentParameters);
66  std::vector<Dihedral> getMMDihedrals(std::string t1, std::string t2, std::string t3, std::string t4) const;
68  ImproperDihedral getMMImproperDihedral(std::string central, std::string t2, std::string t3, std::string t4) const;
73  const std::map<std::string, int>& getC6IndicesMap() const;
74 
75  // These functions add certain parameters of the MM model
76  void addCharge(std::string atomType, double charge);
77  void addDihedral(DihedralType dihedralType, DihedralParameters dihedralParameters);
78  void addImproperDihedral(ImproperDihedralType improperDihedralType, ImproperDihedralParameters improperDihedralParameters);
79 
80  int evaluateNumDistinctAtomTypes(const AtomTypesHolder& atomTypes) const;
81 
82  // These functions return all parameters of a certain type by reference
83  std::map<BondType, BondParameters>& getBonds();
84  std::map<AngleType, AngleParameters>& getAngles();
85  std::map<DihedralType, DihedralParameters>& getDihedrals();
86  std::map<ImproperDihedralType, ImproperDihedralParameters>& getImproperDihedrals();
87  std::map<std::string, double>& getCharges();
88 
89  // These functions return all parameters of a certain type by const. reference
90  const std::map<BondType, BondParameters>& getBonds() const;
91  const std::map<AngleType, AngleParameters>& getAngles() const;
92  const std::map<DihedralType, DihedralParameters>& getDihedrals() const;
93  const std::map<ImproperDihedralType, ImproperDihedralParameters>& getImproperDihedrals() const;
94  const std::map<std::string, double>& getCharges() const;
95 
96  private:
97  std::map<std::string, double> charges_;
98  Eigen::MatrixXf c6Matrix_; // store as floats for efficiency
99  std::map<std::string, int> c6IndicesMap_;
100  std::vector<double> nonCovalentParameters_;
101  std::map<DihedralType, DihedralParameters> dihedrals_;
102  std::map<ImproperDihedralType, ImproperDihedralParameters> improperDihedrals_;
103 };
104 
105 } // namespace MolecularMechanics
106 } // namespace Scine
107 
108 #endif // MOLECULARMECHANICS_SFAMPARAMETERS_H
void resetC6Matrix()
Resets the C6 matrix and the indices map.
Definition: SfamParameters.cpp:66
Definition: ImproperDihedralParameters.h:19
std::vector< double > getChargesForEachAtom(const AtomTypesHolder &atomTypes) const
Getter for the partial atomic charges for each atom.
Definition: SfamParameters.cpp:34
Class treating an improper dihedral interaction, based solely on the angle (in rad), i.e. in 1 dimension.
Definition: ImproperDihedral.h:19
void setNonCovalentParameters(std::vector< double > nonCovalentParameters)
Setter for the non-covalent parameters a1, s8, a2, beta and the atomic charges scaling factor...
Definition: SfamParameters.cpp:100
Unique descriptor a dihedral for given atom types. (useful for maps/unordered_maps) ...
Definition: DihedralType.h:24
Class containing the parameters for an MM dihedral.
Definition: DihedralParameters.h:19
std::vector< Dihedral > getMMDihedrals(std::string t1, std::string t2, std::string t3, std::string t4) const
Get MMDihedrals for the four atom types t1, t2, t3 and t4. (in principle there could be more than one...
Definition: SfamParameters.cpp:104
Class containing the MM atom types of the atoms in a molecular system.
Definition: AtomTypesHolder.h:21
ImproperDihedral getMMImproperDihedral(std::string central, std::string t2, std::string t3, std::string t4) const
Get MMImproperDihedral for the four atom types t1, t2, t3 and t4.
Definition: SfamParameters.cpp:124
void setC6(int indexOfAtomTypeA, int indexOfAtomTypeB, float c6)
Setter for the C6 dispersion coefficient with indices of atom types given.
Definition: SfamParameters.cpp:91
void prepareC6Matrix(const AtomTypesHolder &atomTypes)
Resizes the C6 dispersion coefficient matrix to size NxN (N is number of distinct atom types)...
Definition: SfamParameters.cpp:51
bool sanityCheck(const AtomTypesHolder &atomTypes) const
Checks whether the SFAM parameters are valid.
Definition: SfamParameters.cpp:17
Base class for the MM parameters classes, which holds and manages the bond and angle parameters as th...
Definition: MMParameters.h:30
std::vector< double > getNonCovalentParameters() const
Getter for the non-covalent parameters a1, s8, a2, beta and the atomic charges scaling factor...
Definition: SfamParameters.cpp:96
Unique descriptor an improper dihedral for given atom types. (useful for maps/unordered_maps) ...
Definition: ImproperDihedralType.h:24
Class containing the parameters for SFAM&#39;s MM model obtained after parsing a SFAM parameter file...
Definition: SfamParameters.h:37
const std::map< std::string, int > & getC6IndicesMap() const
Returns a const reference to the map containing the unique atom types of the system and their index i...
Definition: SfamParameters.cpp:185
float getC6(int indexOfAtomTypeA, int indexOfAtomTypeB) const
Getter for the C6 dispersion coefficient with indices of atom types given.
Definition: SfamParameters.cpp:79