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