Scine::Swoose  2.1.0
This is the SCINE module Swoose.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
UpdateFunctionManager.h
Go to the documentation of this file.
1 
8 #ifndef MMPARAMETRIZATION_UPDATEFUNCTIONMANAGER_H
9 #define MMPARAMETRIZATION_UPDATEFUNCTIONMANAGER_H
10 
13 #include <memory>
14 
15 namespace Scine {
16 
17 namespace Utils {
18 class Settings;
19 } // namespace Utils
20 
21 namespace MMParametrization {
22 struct ParametrizationData;
23 
28 enum class ParameterToOptimize { Bond = (1 << 0), Angle = (1 << 1), Dihedral = (1 << 2), ImproperDihedral = (1 << 3) };
32 constexpr inline ParameterToOptimize operator|(ParameterToOptimize p1, ParameterToOptimize p2) {
33  using utype = std::underlying_type<ParameterToOptimize>::type;
34  return static_cast<ParameterToOptimize>(static_cast<utype>(p1) | static_cast<utype>(p2));
35 }
36 
38  public:
39  int hessianCounter;
40  UpdateFunctionManager(ParametrizationData& data, std::shared_ptr<Utils::Settings> settings,
41  ParameterToOptimize typeOfParameter);
45  void updateErrors(const Eigen::VectorXd& parameters, Eigen::VectorXd& errors) override;
52  int getNumberOfDataPoints(const Eigen::VectorXd& parameters) const override;
56  MolecularMechanics::SfamParameters vectorToParameters(const Eigen::VectorXd& parameters);
60  Eigen::VectorXd parametersToVector(const MolecularMechanics::SfamParameters& parameters);
64  void setInitialParameters(Eigen::VectorXd initialParameters);
65 
66  private:
67  // The data used within all MM parametrization classes
68  ParametrizationData& data_;
69  // The settings
70  std::shared_ptr<Utils::Settings> settings_;
71  // Returns a list of atom indices for which the Hessian contributions need to be calculated.
72  std::vector<int> getAtomsToConsiderForHessian() const;
73  // The type of parameter that is going to be optimized here
74  ParameterToOptimize typeOfParameter_;
75  // The MM calculator
76  std::unique_ptr<MolecularMechanics::SfamMolecularMechanicsCalculator> mmCalculator_;
77  // Determines whether the given parameter type is a subset of the member typeOfParameter_
78  bool typeOfParameterShouldBeOptimized(const ParameterToOptimize& typeToCheck) const {
79  auto combined = typeToCheck | typeOfParameter_;
80  return combined == typeOfParameter_;
81  }
82  // The initial parameters as an Eigen::VectorXd, this is saved so that the constraints can be set
83  Eigen::VectorXd initialParameters_;
84  // Whether to add constraints in the parameter optimization
85  bool addConstraints_;
86  // Holds the minimum and maximum values of the parameters as a scaling factor of the initial parameter.
87  Eigen::VectorXd minimumParameters_;
88  Eigen::VectorXd maximumParameters_;
89 };
90 
91 } // namespace MMParametrization
92 } // namespace Scine
93 
94 #endif // MMPARAMETRIZATION_UPDATEFUNCTIONMANAGER_H
void updateErrors(const Eigen::VectorXd &parameters, Eigen::VectorXd &errors) override
Update the errors vector for the least squares optimization for a given parameter.
Definition: UpdateFunctionManager.cpp:38
int getNumberOfDataPoints(const Eigen::VectorXd &parameters) const override
This function returns the number of data points present in the least squares optimization, i.e., the number of Hessian elements for the parameter to be optimized and optionally 2 data points to implement constraints on the parameter.
Definition: UpdateFunctionManager.cpp:216
MolecularMechanics::SfamParameters vectorToParameters(const Eigen::VectorXd &parameters)
Function that converts the variables of the optimization to the parameters object.
Definition: UpdateFunctionManager.cpp:100
void setInitialParameters(Eigen::VectorXd initialParameters)
Set the initial variables of the optimization (variables are MM parameters).
Definition: UpdateFunctionManager.cpp:228
Eigen::VectorXd parametersToVector(const MolecularMechanics::SfamParameters &parameters)
Function that converts the parameters object to the variables of the optimization.
Definition: UpdateFunctionManager.cpp:167
Definition: UpdateFunctionManager.h:37
This struct holds all objects used inside the MM parametrization algorithm.
Definition: ParametrizationData.h:29
Class containing the parameters for SFAM&#39;s MM model obtained after parsing a SFAM parameter file...
Definition: SfamParameters.h:37