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

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 c76b8b was 94453f1, checked in by Frederik Heber <heber@…>, 11 years ago

Added getCoordinator() to EmpiricalPotential interface.

  • implemented with all specific potentials.
  • Property mode set to 100644
File size: 7.1 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
55class Fragment;
56
57// static definitions
58const PairPotential_LennardJones::ParameterNames_t
59PairPotential_LennardJones::ParameterNames =
60 boost::assign::list_of<std::string>
61 ("epsilon")
62 ("sigma")
63 ;
64const std::string PairPotential_LennardJones::potential_token("lennardjones");
65Coordinator::ptr PairPotential_LennardJones::coordinator(new TwoBody_Length());
66
67void PairPotential_LennardJones::setDefaultParameters()
68{
69 params[epsilon] = 1e-5;
70 params[sigma] = 8.2;
71}
72
73PairPotential_LennardJones::PairPotential_LennardJones() :
74 EmpiricalPotential(),
75 params(parameters_t(MAXPARAMS, 0.))
76{
77 // have some decent defaults for parameter_derivative checking
78 setDefaultParameters();
79}
80
81PairPotential_LennardJones::PairPotential_LennardJones(
82 const ParticleTypes_t &_ParticleTypes
83 ) :
84 EmpiricalPotential(_ParticleTypes),
85 params(parameters_t(MAXPARAMS, 0.))
86{
87 // have some decent defaults for parameter_derivative checking
88 setDefaultParameters();
89}
90
91PairPotential_LennardJones::PairPotential_LennardJones(
92 const ParticleTypes_t &_ParticleTypes,
93 const double _epsilon,
94 const double _sigma
95 ) :
96 EmpiricalPotential(_ParticleTypes),
97 params(parameters_t(MAXPARAMS, 0.))
98{
99 params[epsilon] = _epsilon;
100 params[sigma] = _sigma;
101}
102
103void PairPotential_LennardJones::setParameters(const parameters_t &_params)
104{
105 const size_t paramsDim = _params.size();
106 ASSERT( paramsDim <= getParameterDimension(),
107 "PairPotential_LennardJones::setParameters() - we need not more than "
108 +toString(getParameterDimension())+" parameters.");
109 for(size_t i=0;i<paramsDim;++i)
110 params[i] = _params[i];
111
112#ifndef NDEBUG
113 parameters_t check_params(getParameters());
114 check_params.resize(paramsDim); // truncate to same size
115 ASSERT( check_params == _params,
116 "PairPotential_LennardJones::setParameters() - failed, mismatch in to be set "
117 +toString(_params)+" and set "+toString(check_params)+" params.");
118#endif
119}
120
121PairPotential_LennardJones::results_t
122PairPotential_LennardJones::operator()(
123 const arguments_t &arguments
124 ) const
125{
126 ASSERT( arguments.size() == 1,
127 "PairPotential_LennardJones::operator() - requires one argument.");
128 ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
129 arguments, getParticleTypes()),
130 "PairPotential_LennardJones::operator() - types don't match with ones in arguments.");
131 const double &r = arguments[0].distance;
132 const double temp = Helpers::pow(params[sigma]/r, 6);
133 const result_t result = 4.*params[epsilon] * (temp*temp - temp);
134 return std::vector<result_t>(1, result);
135}
136
137PairPotential_LennardJones::derivative_components_t
138PairPotential_LennardJones::derivative(
139 const arguments_t &arguments
140 ) const
141{
142 ASSERT( arguments.size() == 1,
143 "PairPotential_LennardJones::operator() - requires no argument.");
144 ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
145 arguments, getParticleTypes()),
146 "PairPotential_LennardJones::operator() - types don't match with ones in arguments.");
147 const double &r = arguments[0].distance;
148 const double sigma6 = Helpers::pow(params[sigma], 6);
149 const result_t result =
150 4.*params[epsilon] * (
151 sigma6*sigma6*(-12.) / Helpers::pow(r,13)
152 - sigma6*(-6.) /Helpers::pow(r,7)
153 );
154 derivative_components_t results(1, result);
155 return results;
156}
157
158PairPotential_LennardJones::results_t
159PairPotential_LennardJones::parameter_derivative(
160 const arguments_t &arguments,
161 const size_t index
162 ) const
163{
164 ASSERT( arguments.size() == 1,
165 "PairPotential_LennardJones::parameter_derivative() - requires no argument.");
166 ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
167 arguments, getParticleTypes()),
168 "PairPotential_LennardJones::operator() - types don't match with ones in arguments.");
169 const double &r = arguments[0].distance;
170 switch (index) {
171 case epsilon:
172 {
173 const double temp = Helpers::pow(params[sigma]/r, 6);
174 const result_t result = 4. * (temp*temp - temp);
175 return std::vector<result_t>(1, result);
176 break;
177 }
178 case sigma:
179 {
180 const double r6 = Helpers::pow(r, 6);
181 const result_t result =
182 4.*params[epsilon] * (
183 12. * Helpers::pow(params[sigma],11)/(r6*r6)
184 - 6. * Helpers::pow(params[sigma],5)/r6
185 );
186 return std::vector<result_t>(1, result);
187 break;
188 }
189 default:
190 break;
191 }
192 return std::vector<result_t>(1, 0.);
193}
194
195FunctionModel::extractor_t
196PairPotential_LennardJones::getSpecificExtractor() const
197{
198 Fragment::charges_t charges;
199 charges.resize(getParticleTypes().size());
200 std::transform(getParticleTypes().begin(), getParticleTypes().end(),
201 charges.begin(), boost::lambda::_1);
202 FunctionModel::extractor_t returnfunction =
203 boost::bind(&Extractors::gatherDistancesFromFragment,
204 boost::bind(&Fragment::getPositions, _1),
205 boost::bind(&Fragment::getCharges, _1),
206 charges,
207 _2);
208 return returnfunction;
209}
210
211FunctionModel::filter_t PairPotential_LennardJones::getSpecificFilter() const
212{
213 FunctionModel::filter_t returnfunction =
214 boost::bind(&Extractors::filterArgumentsByParticleTypes,
215 _1,
216 getParticleTypes());
217 return returnfunction;
218}
219
220void
221PairPotential_LennardJones::setParametersToRandomInitialValues(
222 const TrainingData &data)
223{
224 params[PairPotential_LennardJones::epsilon] = 1e-2*rand()/(double)RAND_MAX;
225 params[PairPotential_LennardJones::sigma] = (3.+10.*rand()/(double)RAND_MAX);// 0.5;
226}
227
Note: See TracBrowser for help on using the repository browser.