source: src/Potentials/Specifics/PairPotential_LennardJones.cpp@ f4496d

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 f4496d was 0932c2, checked in by Frederik Heber <heber@…>, 11 years ago

FIX: RandomNumberEngine_Encapsulation::clone() did not set the seed obtained from its prototype.

  • this caused all engines to always start with default seed of 1, irrespective of the set parameter. The engine is always cloned from prototypes in the factory that can be manipulated to change their parameters. Hence, clone needs to set the seed (the engine's only parameter).
  • EmpiricalPotentials now use RandomNumberGenerator to obtain random starting values.
  • TESTFIX: fit-potential regression tests now use fixed seed in order to always take the same amount of time (some test runs in debug mode take very long because of ill-chosen random values).
  • Property mode set to 100644
File size: 7.3 KB
Line 
1/*
2 * Project: MoleCuilder
3 * Description: creates and alters molecular systems
4 * Copyright (C) 2013 Frederik Heber. All rights reserved.
5 * Please see the COPYING file or "Copyright notice" in builder.cpp for details.
6 *
7 *
8 * This file is part of MoleCuilder.
9 *
10 * MoleCuilder is free software: you can redistribute it and/or modify
11 * it under the terms of the GNU General Public License as published by
12 * the Free Software Foundation, either version 2 of the License, or
13 * (at your option) any later version.
14 *
15 * MoleCuilder is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU General Public License for more details.
19 *
20 * You should have received a copy of the GNU General Public License
21 * along with MoleCuilder. If not, see <http://www.gnu.org/licenses/>.
22 */
23
24/*
25 * PairPotential_LennardJones.cpp
26 *
27 * Created on: Jul 05, 2013
28 * Author: heber
29 */
30
31
32// include config.h
33#ifdef HAVE_CONFIG_H
34#include <config.h>
35#endif
36
37#include "CodePatterns/MemDebug.hpp"
38
39#include "PairPotential_LennardJones.hpp"
40
41#include <boost/assign/list_of.hpp> // for 'map_list_of()'
42#include <boost/bind.hpp>
43#include <boost/lambda/lambda.hpp>
44#include <cmath>
45#include <string>
46
47#include "CodePatterns/Assert.hpp"
48
49#include "FunctionApproximation/Extractors.hpp"
50#include "FunctionApproximation/TrainingData.hpp"
51#include "Potentials/helpers.hpp"
52#include "Potentials/InternalCoordinates/TwoBody_Length.hpp"
53#include "Potentials/ParticleTypeCheckers.hpp"
54#include "RandomNumbers/RandomNumberGeneratorFactory.hpp"
55#include "RandomNumbers/RandomNumberGenerator.hpp"
56
57class Fragment;
58
59// static definitions
60const PairPotential_LennardJones::ParameterNames_t
61PairPotential_LennardJones::ParameterNames =
62 boost::assign::list_of<std::string>
63 ("epsilon")
64 ("sigma")
65 ;
66const std::string PairPotential_LennardJones::potential_token("lennardjones");
67Coordinator::ptr PairPotential_LennardJones::coordinator(new TwoBody_Length());
68
69void PairPotential_LennardJones::setDefaultParameters()
70{
71 params[epsilon] = 1e-5;
72 params[sigma] = 8.2;
73}
74
75PairPotential_LennardJones::PairPotential_LennardJones() :
76 EmpiricalPotential(),
77 params(parameters_t(MAXPARAMS, 0.))
78{
79 // have some decent defaults for parameter_derivative checking
80 setDefaultParameters();
81}
82
83PairPotential_LennardJones::PairPotential_LennardJones(
84 const ParticleTypes_t &_ParticleTypes
85 ) :
86 EmpiricalPotential(_ParticleTypes),
87 params(parameters_t(MAXPARAMS, 0.))
88{
89 // have some decent defaults for parameter_derivative checking
90 setDefaultParameters();
91}
92
93PairPotential_LennardJones::PairPotential_LennardJones(
94 const ParticleTypes_t &_ParticleTypes,
95 const double _epsilon,
96 const double _sigma
97 ) :
98 EmpiricalPotential(_ParticleTypes),
99 params(parameters_t(MAXPARAMS, 0.))
100{
101 params[epsilon] = _epsilon;
102 params[sigma] = _sigma;
103}
104
105void PairPotential_LennardJones::setParameters(const parameters_t &_params)
106{
107 const size_t paramsDim = _params.size();
108 ASSERT( paramsDim <= getParameterDimension(),
109 "PairPotential_LennardJones::setParameters() - we need not more than "
110 +toString(getParameterDimension())+" parameters.");
111 for(size_t i=0;i<paramsDim;++i)
112 params[i] = _params[i];
113
114#ifndef NDEBUG
115 parameters_t check_params(getParameters());
116 check_params.resize(paramsDim); // truncate to same size
117 ASSERT( check_params == _params,
118 "PairPotential_LennardJones::setParameters() - failed, mismatch in to be set "
119 +toString(_params)+" and set "+toString(check_params)+" params.");
120#endif
121}
122
123PairPotential_LennardJones::results_t
124PairPotential_LennardJones::operator()(
125 const arguments_t &arguments
126 ) const
127{
128 ASSERT( arguments.size() == 1,
129 "PairPotential_LennardJones::operator() - requires one argument.");
130 ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
131 arguments, getParticleTypes()),
132 "PairPotential_LennardJones::operator() - types don't match with ones in arguments.");
133 const double &r = arguments[0].distance;
134 const double temp = Helpers::pow(params[sigma]/r, 6);
135 const result_t result = 4.*params[epsilon] * (temp*temp - temp);
136 return std::vector<result_t>(1, result);
137}
138
139PairPotential_LennardJones::derivative_components_t
140PairPotential_LennardJones::derivative(
141 const arguments_t &arguments
142 ) const
143{
144 ASSERT( arguments.size() == 1,
145 "PairPotential_LennardJones::operator() - requires no argument.");
146 ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
147 arguments, getParticleTypes()),
148 "PairPotential_LennardJones::operator() - types don't match with ones in arguments.");
149 const double &r = arguments[0].distance;
150 const double sigma6 = Helpers::pow(params[sigma], 6);
151 const result_t result =
152 4.*params[epsilon] * (
153 sigma6*sigma6*(-12.) / Helpers::pow(r,13)
154 - sigma6*(-6.) /Helpers::pow(r,7)
155 );
156 derivative_components_t results(1, result);
157 return results;
158}
159
160PairPotential_LennardJones::results_t
161PairPotential_LennardJones::parameter_derivative(
162 const arguments_t &arguments,
163 const size_t index
164 ) const
165{
166 ASSERT( arguments.size() == 1,
167 "PairPotential_LennardJones::parameter_derivative() - requires no argument.");
168 ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
169 arguments, getParticleTypes()),
170 "PairPotential_LennardJones::operator() - types don't match with ones in arguments.");
171 const double &r = arguments[0].distance;
172 switch (index) {
173 case epsilon:
174 {
175 const double temp = Helpers::pow(params[sigma]/r, 6);
176 const result_t result = 4. * (temp*temp - temp);
177 return std::vector<result_t>(1, result);
178 break;
179 }
180 case sigma:
181 {
182 const double r6 = Helpers::pow(r, 6);
183 const result_t result =
184 4.*params[epsilon] * (
185 12. * Helpers::pow(params[sigma],11)/(r6*r6)
186 - 6. * Helpers::pow(params[sigma],5)/r6
187 );
188 return std::vector<result_t>(1, result);
189 break;
190 }
191 default:
192 break;
193 }
194 return std::vector<result_t>(1, 0.);
195}
196
197FunctionModel::extractor_t
198PairPotential_LennardJones::getSpecificExtractor() const
199{
200 Fragment::charges_t charges;
201 charges.resize(getParticleTypes().size());
202 std::transform(getParticleTypes().begin(), getParticleTypes().end(),
203 charges.begin(), boost::lambda::_1);
204 FunctionModel::extractor_t returnfunction =
205 boost::bind(&Extractors::gatherDistancesFromFragment,
206 boost::bind(&Fragment::getPositions, _1),
207 boost::bind(&Fragment::getCharges, _1),
208 charges,
209 _2);
210 return returnfunction;
211}
212
213FunctionModel::filter_t PairPotential_LennardJones::getSpecificFilter() const
214{
215 FunctionModel::filter_t returnfunction =
216 boost::bind(&Extractors::filterArgumentsByParticleTypes,
217 _1,
218 getParticleTypes());
219 return returnfunction;
220}
221
222void
223PairPotential_LennardJones::setParametersToRandomInitialValues(
224 const TrainingData &data)
225{
226 RandomNumberGenerator &random = RandomNumberGeneratorFactory::getInstance().makeRandomNumberGenerator();
227 const double rng_min = random.min();
228 const double rng_max = random.max();
229 params[PairPotential_LennardJones::epsilon] = 1e-2*(random()/(rng_max-rng_min));
230 params[PairPotential_LennardJones::sigma] = (3.+10.*(random()/(rng_max-rng_min)));// 0.5;
231}
232
Note: See TracBrowser for help on using the repository browser.