source: src/Potentials/Specifics/PairPotential_LennardJones.hpp@ ff90e3

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 ff90e3 was 1413f4, checked in by Frederik Heber <heber@…>, 12 years ago

Added LennardJones pair potential.

  • this is for later fitting between order-1 fragments.
  • also added regression test for LJ fitting.
  • Property mode set to 100644
File size: 5.6 KB
Line 
1/*
2 * PairPotential_LennardJones.hpp
3 *
4 * Created on: Jul 05, 2013
5 * Author: heber
6 */
7
8#ifndef PAIRPOTENTIAL_LENNARDJONES_HPP_
9#define PAIRPOTENTIAL_LENNARDJONES_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
21class PotentialFactory;
22class TrainingData;
23
24/** This is the implementation of a constant potential to adapt to any offset
25 * in the potential energy.
26 *
27 * This evaluates \f$ k \f$.
28 *
29 */
30class PairPotential_LennardJones :
31 public EmpiricalPotential
32{
33 //!> grant unit test access to internal parts
34 friend class PairPotential_LennardJonesTest;
35 //!> grant PotentialFactory access to default cstor
36 friend class PotentialFactory;
37 // some repeated typedefs to avoid ambiguities
38 typedef FunctionModel::arguments_t arguments_t;
39 typedef FunctionModel::result_t result_t;
40 typedef FunctionModel::results_t results_t;
41 typedef EmpiricalPotential::derivative_components_t derivative_components_t;
42 typedef FunctionModel::parameters_t parameters_t;
43private:
44 /** Private default constructor.
45 *
46 * This prevents creation of potential without set ParticleTypes_t.
47 *
48 * \note PotentialFactory may use this default cstor
49 *
50 */
51 PairPotential_LennardJones();
52
53public:
54 PairPotential_LennardJones(const ParticleTypes_t &_ParticleTypes);
55 PairPotential_LennardJones(
56 const ParticleTypes_t &_ParticleTypes,
57 const double _epsilon,
58 const double _sigma
59 );
60 virtual ~PairPotential_LennardJones() {}
61
62 /** Setter for parameters as required by FunctionModel interface.
63 *
64 * \param _params given set of parameters
65 */
66 void setParameters(const parameters_t &_params);
67
68 /** Getter for parameters as required by FunctionModel interface.
69 *
70 * \return set of parameters
71 */
72 parameters_t getParameters() const
73 {
74 return params;
75 }
76
77 /** Sets the parameter randomly within the sensible range of each parameter.
78 *
79 * \param data container with training data for guesstimating range
80 */
81 void setParametersToRandomInitialValues(const TrainingData &data);
82
83 /** Getter for the number of parameters of this model function.
84 *
85 * \return number of parameters
86 */
87 size_t getParameterDimension() const
88 {
89 return 2;
90 }
91
92 /** Evaluates the harmonic potential function for the given arguments.
93 *
94 * @param arguments single distance
95 * @return value of the potential function
96 */
97 results_t operator()(const arguments_t &arguments) const;
98
99 /** Evaluates the derivative of the potential function.
100 *
101 * @param arguments single distance
102 * @return vector with derivative with respect to the input degrees of freedom
103 */
104 derivative_components_t derivative(const arguments_t &arguments) const;
105
106 /** Evaluates the derivative of the function with the given \a arguments
107 * with respect to a specific parameter indicated by \a index.
108 *
109 * \param arguments set of arguments as input variables to the function
110 * \param index derivative of which parameter
111 * \return result vector containing the derivative with respect to the given
112 * input
113 */
114 results_t parameter_derivative(const arguments_t &arguments, const size_t index) const;
115
116 /** Return the token name of this specific potential.
117 *
118 * \return token name of the potential
119 */
120 const std::string& getToken() const
121 { return potential_token; }
122
123 /** Returns a vector of parameter names.
124 *
125 * This is required from the specific implementation
126 *
127 * \return vector of strings containing parameter names
128 */
129 const ParameterNames_t& getParameterNames() const
130 { return ParameterNames; }
131
132 /** States whether lower and upper boundaries should be used to constraint
133 * the parameter search for this function model.
134 *
135 * \return true - constraints should be used, false - else
136 */
137 bool isBoxConstraint() const {
138 return true;
139 }
140
141 /** Returns a vector which are the lower boundaries for each parameter_t
142 * of this FunctionModel.
143 *
144 * \return vector of parameter_t resembling lowest allowed values
145 */
146 parameters_t getLowerBoxConstraints() const {
147 parameters_t lowerbound(getParameterDimension(), 0.);
148 return lowerbound;
149 }
150
151 /** Returns a vector which are the upper boundaries for each parameter_t
152 * of this FunctionModel.
153 *
154 * \return vector of parameter_t resembling highest allowed values
155 */
156 parameters_t getUpperBoxConstraints() const {
157 return parameters_t(getParameterDimension(), std::numeric_limits<double>::max());
158 }
159
160 /** Returns a bound function to be used with TrainingData, extracting distances
161 * from a Fragment.
162 *
163 * \return bound function extracting distances from a fragment
164 */
165 FunctionModel::extractor_t getSpecificExtractor() const;
166
167 /** Returns a bound function to be used with TrainingData, extracting distances
168 * from a Fragment.
169 *
170 * \return bound function extracting distances from a fragment
171 */
172 FunctionModel::filter_t getSpecificFilter() const;
173
174 /** Returns the number of arguments the underlying function requires.
175 *
176 * \return number of arguments of the function
177 */
178 size_t getSpecificArgumentCount() const
179 { return 0; }
180
181 enum parameter_enum_t {
182 epsilon=0,
183 sigma=1,
184 MAXPARAMS
185 };
186
187private:
188 /** Sets some sensible default parameter values.
189 *
190 */
191 void setDefaultParameters();
192
193private:
194 //!> parameter vector with parameters as in enum parameter_enum_t
195 parameters_t params;
196
197 //!> static definitions of the parameter name for this potential
198 static const ParameterNames_t ParameterNames;
199
200 //!> static token of this potential type
201 static const std::string potential_token;
202};
203
204#endif /* PAIRPOTENTIAL_LENNARDJONES_HPP_ */
Note: See TracBrowser for help on using the repository browser.