source: src/Potentials/Specifics/PairPotential_Angle.hpp@ 1e565c

Action_Thermostats Add_AtomRandomPerturbation Add_FitFragmentPartialChargesAction Add_RotateAroundBondAction Add_SelectAtomByNameAction Added_ParseSaveFragmentResults AddingActions_SaveParseParticleParameters Adding_Graph_to_ChangeBondActions Adding_MD_integration_tests Adding_ParticleName_to_Atom Adding_StructOpt_integration_tests AtomFragments Automaking_mpqc_open AutomationFragmentation_failures Candidate_v1.5.4 Candidate_v1.6.0 Candidate_v1.6.1 ChangeBugEmailaddress ChangingTestPorts ChemicalSpaceEvaluator CombiningParticlePotentialParsing Combining_Subpackages Debian_Package_split Debian_package_split_molecuildergui_only Disabling_MemDebug Docu_Python_wait EmpiricalPotential_contain_HomologyGraph EmpiricalPotential_contain_HomologyGraph_documentation Enable_parallel_make_install Enhance_userguide Enhanced_StructuralOptimization Enhanced_StructuralOptimization_continued Example_ManyWaysToTranslateAtom Exclude_Hydrogens_annealWithBondGraph FitPartialCharges_GlobalError Fix_BoundInBox_CenterInBox_MoleculeActions Fix_ChargeSampling_PBC Fix_ChronosMutex Fix_FitPartialCharges Fix_FitPotential_needs_atomicnumbers Fix_ForceAnnealing Fix_IndependentFragmentGrids Fix_ParseParticles Fix_ParseParticles_split_forward_backward_Actions Fix_PopActions Fix_QtFragmentList_sorted_selection Fix_Restrictedkeyset_FragmentMolecule Fix_StatusMsg Fix_StepWorldTime_single_argument Fix_Verbose_Codepatterns Fix_fitting_potentials Fixes ForceAnnealing_goodresults ForceAnnealing_oldresults ForceAnnealing_tocheck ForceAnnealing_with_BondGraph ForceAnnealing_with_BondGraph_continued ForceAnnealing_with_BondGraph_continued_betteresults ForceAnnealing_with_BondGraph_contraction-expansion FragmentAction_writes_AtomFragments FragmentMolecule_checks_bonddegrees GeometryObjects Gui_Fixes Gui_displays_atomic_force_velocity ImplicitCharges IndependentFragmentGrids IndependentFragmentGrids_IndividualZeroInstances IndependentFragmentGrids_IntegrationTest IndependentFragmentGrids_Sole_NN_Calculation JobMarket_RobustOnKillsSegFaults JobMarket_StableWorkerPool JobMarket_unresolvable_hostname_fix MoreRobust_FragmentAutomation ODR_violation_mpqc_open PartialCharges_OrthogonalSummation PdbParser_setsAtomName PythonUI_with_named_parameters QtGui_reactivate_TimeChanged_changes Recreated_GuiChecks Rewrite_FitPartialCharges RotateToPrincipalAxisSystem_UndoRedo SaturateAtoms_findBestMatching SaturateAtoms_singleDegree StoppableMakroAction Subpackage_CodePatterns Subpackage_JobMarket Subpackage_LinearAlgebra Subpackage_levmar Subpackage_mpqc_open Subpackage_vmg Switchable_LogView ThirdParty_MPQC_rebuilt_buildsystem TrajectoryDependenant_MaxOrder TremoloParser_IncreasedPrecision TremoloParser_MultipleTimesteps TremoloParser_setsAtomName Ubuntu_1604_changes stable
Last change on this file since 1e565c was 2bde72, checked in by Frederik Heber <heber@…>, 12 years ago

FIX: Removed CodePattern/Assert include from various PairPotential headers.

  • Property mode set to 100644
File size: 5.1 KB
RevLine 
[a63187]1/*
2 * PairPotential_Angle.hpp
3 *
4 * Created on: Oct 11, 2012
5 * Author: heber
6 */
7
8#ifndef PAIRPOTENTIAL_ANGLE_HPP_
9#define PAIRPOTENTIAL_ANGLE_HPP_
10
11
12// include config.h
13#ifdef HAVE_CONFIG_H
14#include <config.h>
15#endif
16
17#include <limits>
18
19#include "Potentials/EmpiricalPotential.hpp"
20
[d52819]21class TrainingData;
22
[a63187]23/** This is the implementation of a harmonic angle potential.
24 *
25 * This evaluates \f$ k \cdot (\theta -\theta_0)^2 \f$.
26 *
27 */
[ed2551]28class PairPotential_Angle :
[fdd23a]29 public EmpiricalPotential
[a63187]30{
31 //!> grant unit test access to internal parts
32 friend class PairPotential_AngleTest;
33 // some repeated typedefs to avoid ambiguities
34 typedef FunctionModel::arguments_t arguments_t;
35 typedef FunctionModel::result_t result_t;
36 typedef FunctionModel::results_t results_t;
37 typedef EmpiricalPotential::derivative_components_t derivative_components_t;
38 typedef FunctionModel::parameters_t parameters_t;
39public:
[ed2551]40 PairPotential_Angle(const ParticleTypes_t &_ParticleTypes);
[a63187]41 PairPotential_Angle(
[ed2551]42 const ParticleTypes_t &_ParticleTypes,
[a63187]43 const double _spring_constant,
44 const double _equilibrium_distance,
45 const double _energy_offset);
46 virtual ~PairPotential_Angle() {}
47
48 /** Setter for parameters as required by FunctionModel interface.
49 *
50 * \param _params given set of parameters
51 */
[086070]52 void setParameters(const parameters_t &_params);
[a63187]53
54 /** Getter for parameters as required by FunctionModel interface.
55 *
56 * \return set of parameters
57 */
58 parameters_t getParameters() const
59 {
60 return params;
61 }
62
[d52819]63 /** Sets the parameter randomly within the sensible range of each parameter.
64 *
65 * \param data container with training data for guesstimating range
66 */
67 void setParametersToRandomInitialValues(const TrainingData &data);
68
[a63187]69 /** Getter for the number of parameters of this model function.
70 *
71 * \return number of parameters
72 */
73 size_t getParameterDimension() const
74 {
75 return 3;
76 }
77
78 /** Evaluates the harmonic potential function for the given arguments.
79 *
80 * @param arguments single distance
81 * @return value of the potential function
82 */
83 results_t operator()(const arguments_t &arguments) const;
84
85 /** Evaluates the derivative of the potential function.
86 *
87 * @param arguments single distance
88 * @return vector with derivative with respect to the input degrees of freedom
89 */
90 derivative_components_t derivative(const arguments_t &arguments) const;
91
92 /** Evaluates the derivative of the function with the given \a arguments
93 * with respect to a specific parameter indicated by \a index.
94 *
95 * \param arguments set of arguments as input variables to the function
96 * \param index derivative of which parameter
97 * \return result vector containing the derivative with respect to the given
98 * input
99 */
100 results_t parameter_derivative(const arguments_t &arguments, const size_t index) const;
101
[6efcae]102 /** Return the token name of this specific potential.
[67cd3a]103 *
[6efcae]104 * \return token name of the potential
[67cd3a]105 */
[ed2551]106 const std::string& getToken() const
107 { return potential_token; }
108
109 /** Returns a vector of parameter names.
110 *
111 * This is required from the specific implementation
112 *
113 * \return vector of strings containing parameter names
114 */
115 const ParameterNames_t& getParameterNames() const
116 { return ParameterNames; }
[67cd3a]117
[a63187]118 /** States whether lower and upper boundaries should be used to constraint
119 * the parameter search for this function model.
120 *
121 * \return true - constraints should be used, false - else
122 */
123 bool isBoxConstraint() const {
124 return true;
125 }
126
127 /** Returns a vector which are the lower boundaries for each parameter_t
128 * of this FunctionModel.
129 *
130 * \return vector of parameter_t resembling lowest allowed values
131 */
132 parameters_t getLowerBoxConstraints() const {
133 parameters_t lowerbounds(getParameterDimension(), -std::numeric_limits<double>::max());
[9340ee]134 lowerbounds[equilibrium_distance] = -1.;
[a63187]135 return lowerbounds;
136 }
137
138 /** Returns a vector which are the upper boundaries for each parameter_t
139 * of this FunctionModel.
140 *
141 * \return vector of parameter_t resembling highest allowed values
142 */
143 parameters_t getUpperBoxConstraints() const {
[9340ee]144 parameters_t upperbounds(getParameterDimension(), std::numeric_limits<double>::max());
145 upperbounds[equilibrium_distance] = 1.;
146 return upperbounds;
[a63187]147 }
148
[7b019a]149 /** Returns a bound function to be used with TrainingData, extracting distances
150 * from a Fragment.
151 *
152 * \return bound function extracting distances from a fragment
153 */
[da2d5c]154 FunctionModel::extractor_t getFragmentSpecificExtractor() const;
[7b019a]155
[a63187]156 enum parameter_enum_t {
157 spring_constant=0,
158 equilibrium_distance=1,
159 energy_offset=2,
160 MAXPARAMS
161 };
162
163private:
164 result_t
165 function_theta(
166 const double &r_ij,
167 const double &r_ik,
168 const double &r_jk
169 ) const;
170private:
171 //!> parameter vector with parameters as in enum parameter_enum_t
172 parameters_t params;
[ed2551]173
174 //!> static definitions of the parameter name for this potential
175 static const ParameterNames_t ParameterNames;
176
177 //!> static token of this potential type
178 static const std::string potential_token;
[a63187]179};
180
181#endif /* PAIRPOTENTIAL_ANGLE_HPP_ */
Note: See TracBrowser for help on using the repository browser.