source: src/Potentials/Specifics/PairPotential_LennardJones.cpp@ 35d171

Action_Thermostats Add_AtomRandomPerturbation Add_FitFragmentPartialChargesAction Add_RotateAroundBondAction Add_SelectAtomByNameAction Added_ParseSaveFragmentResults Adding_Graph_to_ChangeBondActions Adding_MD_integration_tests Adding_StructOpt_integration_tests 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_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 FragmentMolecule_checks_bonddegrees GeometryObjects Gui_Fixes Gui_displays_atomic_force_velocity IndependentFragmentGrids IndependentFragmentGrids_IndividualZeroInstances IndependentFragmentGrids_IntegrationTest IndependentFragmentGrids_Sole_NN_Calculation JobMarket_RobustOnKillsSegFaults JobMarket_StableWorkerPool JobMarket_unresolvable_hostname_fix ODR_violation_mpqc_open PartialCharges_OrthogonalSummation PythonUI_with_named_parameters QtGui_reactivate_TimeChanged_changes Recreated_GuiChecks 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 35d171 was 7d320c, checked in by Frederik Heber <heber@…>, 10 years ago

MEMFIX: Static coordinator instances of specificv potentials need to be Memory::ignore()'d.

  • 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(Memory::ignore(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 list_of_arguments_t &listarguments
126 ) const
127{
128 result_t result = 0.;
129 for(list_of_arguments_t::const_iterator iter = listarguments.begin();
130 iter != listarguments.end(); ++iter) {
131 const arguments_t &arguments = *iter;
132 ASSERT( arguments.size() == 1,
133 "PairPotential_LennardJones::operator() - requires one argument.");
134 ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
135 arguments, getParticleTypes()),
136 "PairPotential_LennardJones::operator() - types don't match with ones in arguments.");
137 const double &r = arguments[0].distance;
138 const double temp = Helpers::pow(params[sigma]/r, 6);
139 result += 4.*params[epsilon] * (temp*temp - temp);
140 }
141 return results_t(1, result);
142}
143
144PairPotential_LennardJones::derivative_components_t
145PairPotential_LennardJones::derivative(
146 const list_of_arguments_t &listarguments
147 ) const
148{
149 const double sigma6 = Helpers::pow(params[sigma], 6);
150 result_t result = 0.;
151 for(list_of_arguments_t::const_iterator iter = listarguments.begin();
152 iter != listarguments.end(); ++iter) {
153 const arguments_t &arguments = *iter;
154 ASSERT( arguments.size() == 1,
155 "PairPotential_LennardJones::operator() - requires no argument.");
156 ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
157 arguments, getParticleTypes()),
158 "PairPotential_LennardJones::operator() - types don't match with ones in arguments.");
159 const double &r = arguments[0].distance;
160 result +=
161 4.*params[epsilon] * (
162 sigma6*sigma6*(-12.) / Helpers::pow(r,13)
163 - sigma6*(-6.) /Helpers::pow(r,7)
164 );
165 }
166 return derivative_components_t(1, result);
167}
168
169PairPotential_LennardJones::results_t
170PairPotential_LennardJones::parameter_derivative(
171 const list_of_arguments_t &listarguments,
172 const size_t index
173 ) const
174{
175 result_t result = 0.;
176 for(list_of_arguments_t::const_iterator iter = listarguments.begin();
177 iter != listarguments.end(); ++iter) {
178 const arguments_t &arguments = *iter;
179 ASSERT( arguments.size() == 1,
180 "PairPotential_LennardJones::parameter_derivative() - requires no argument.");
181 ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
182 arguments, getParticleTypes()),
183 "PairPotential_LennardJones::operator() - types don't match with ones in arguments.");
184 const double &r = arguments[0].distance;
185 switch (index) {
186 case epsilon:
187 {
188 const double temp = Helpers::pow(params[sigma]/r, 6);
189 result += 4. * (temp*temp - temp);
190 break;
191 }
192 case sigma:
193 {
194 const double r6 = Helpers::pow(r, 6);
195 result +=
196 4.*params[epsilon] * (
197 12. * Helpers::pow(params[sigma],11)/(r6*r6)
198 - 6. * Helpers::pow(params[sigma],5)/r6
199 );
200 break;
201 }
202 default:
203 break;
204 }
205 }
206 return results_t(1, result);
207}
208
209FunctionModel::filter_t PairPotential_LennardJones::getSpecificFilter() const
210{
211 FunctionModel::filter_t returnfunction =
212 boost::bind(&Extractors::filterArgumentsByParticleTypes,
213 _1,
214 getParticleTypes());
215 return returnfunction;
216}
217
218void
219PairPotential_LennardJones::setParametersToRandomInitialValues(
220 const TrainingData &data)
221{
222 RandomNumberGenerator &random = RandomNumberGeneratorFactory::getInstance().makeRandomNumberGenerator();
223 const double rng_min = random.min();
224 const double rng_max = random.max();
225 params[PairPotential_LennardJones::epsilon] = 1e-2*(random()/(rng_max-rng_min));
226 params[PairPotential_LennardJones::sigma] = (3.+10.*(random()/(rng_max-rng_min)));// 0.5;
227}
228
Note: See TracBrowser for help on using the repository browser.