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 |
|
---|
54 | class Fragment;
|
---|
55 |
|
---|
56 | // static definitions
|
---|
57 | const PairPotential_LennardJones::ParameterNames_t
|
---|
58 | PairPotential_LennardJones::ParameterNames =
|
---|
59 | boost::assign::list_of<std::string>
|
---|
60 | ("epsilon")
|
---|
61 | ("sigma")
|
---|
62 | ;
|
---|
63 | const std::string PairPotential_LennardJones::potential_token("lennardjones");
|
---|
64 |
|
---|
65 | void PairPotential_LennardJones::setDefaultParameters()
|
---|
66 | {
|
---|
67 | params[epsilon] = 1e-5;
|
---|
68 | params[sigma] = 8.2;
|
---|
69 | }
|
---|
70 |
|
---|
71 | PairPotential_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 |
|
---|
79 | PairPotential_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 |
|
---|
89 | PairPotential_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 |
|
---|
101 | void 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 |
|
---|
119 | PairPotential_LennardJones::results_t
|
---|
120 | PairPotential_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 |
|
---|
135 | PairPotential_LennardJones::derivative_components_t
|
---|
136 | PairPotential_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 |
|
---|
156 | PairPotential_LennardJones::results_t
|
---|
157 | PairPotential_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 |
|
---|
193 | FunctionModel::extractor_t
|
---|
194 | PairPotential_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 |
|
---|
209 | FunctionModel::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 |
|
---|
218 | void
|
---|
219 | PairPotential_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 |
|
---|