[1413f4] | 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 |
|
---|
[9eb71b3] | 37 | //#include "CodePatterns/MemDebug.hpp"
|
---|
[1413f4] | 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"
|
---|
[94453f1] | 52 | #include "Potentials/InternalCoordinates/TwoBody_Length.hpp"
|
---|
[1413f4] | 53 | #include "Potentials/ParticleTypeCheckers.hpp"
|
---|
[0932c2] | 54 | #include "RandomNumbers/RandomNumberGeneratorFactory.hpp"
|
---|
| 55 | #include "RandomNumbers/RandomNumberGenerator.hpp"
|
---|
[1413f4] | 56 |
|
---|
| 57 | class Fragment;
|
---|
| 58 |
|
---|
| 59 | // static definitions
|
---|
| 60 | const PairPotential_LennardJones::ParameterNames_t
|
---|
| 61 | PairPotential_LennardJones::ParameterNames =
|
---|
| 62 | boost::assign::list_of<std::string>
|
---|
| 63 | ("epsilon")
|
---|
| 64 | ("sigma")
|
---|
| 65 | ;
|
---|
| 66 | const std::string PairPotential_LennardJones::potential_token("lennardjones");
|
---|
[9eb71b3] | 67 | Coordinator::ptr PairPotential_LennardJones::coordinator( /* Memory::ignore( */ new TwoBody_Length() /* ) */ );
|
---|
[1413f4] | 68 |
|
---|
[9c793c] | 69 | static BindingModel generateBindingModel(const EmpiricalPotential::ParticleTypes_t &_ParticleTypes)
|
---|
[d5ca1a] | 70 | {
|
---|
| 71 | // fill nodes
|
---|
[9c793c] | 72 | BindingModel::vector_nodes_t nodes;
|
---|
[d5ca1a] | 73 | {
|
---|
| 74 | ASSERT( _ParticleTypes.size() == (size_t)2,
|
---|
| 75 | "generateBindingModel() - PairPotential_LennardJones needs two types.");
|
---|
[9c793c] | 76 | nodes.push_back( FragmentNode(_ParticleTypes[0], 1) );
|
---|
| 77 | nodes.push_back( FragmentNode(_ParticleTypes[1], 1) );
|
---|
[d5ca1a] | 78 | }
|
---|
| 79 |
|
---|
| 80 | // there are no edges
|
---|
| 81 | HomologyGraph::edges_t edges;
|
---|
| 82 |
|
---|
[9c793c] | 83 | return BindingModel(nodes, edges);
|
---|
[d5ca1a] | 84 | }
|
---|
| 85 |
|
---|
[1413f4] | 86 | void PairPotential_LennardJones::setDefaultParameters()
|
---|
| 87 | {
|
---|
| 88 | params[epsilon] = 1e-5;
|
---|
| 89 | params[sigma] = 8.2;
|
---|
| 90 | }
|
---|
| 91 |
|
---|
| 92 | PairPotential_LennardJones::PairPotential_LennardJones() :
|
---|
| 93 | EmpiricalPotential(),
|
---|
[d5ca1a] | 94 | params(parameters_t(MAXPARAMS, 0.)),
|
---|
[9c793c] | 95 | bindingmodel(BindingModel())
|
---|
[1413f4] | 96 | {
|
---|
| 97 | // have some decent defaults for parameter_derivative checking
|
---|
| 98 | setDefaultParameters();
|
---|
| 99 | }
|
---|
| 100 |
|
---|
| 101 | PairPotential_LennardJones::PairPotential_LennardJones(
|
---|
| 102 | const ParticleTypes_t &_ParticleTypes
|
---|
| 103 | ) :
|
---|
| 104 | EmpiricalPotential(_ParticleTypes),
|
---|
[d5ca1a] | 105 | params(parameters_t(MAXPARAMS, 0.)),
|
---|
| 106 | bindingmodel(generateBindingModel(_ParticleTypes))
|
---|
[1413f4] | 107 | {
|
---|
| 108 | // have some decent defaults for parameter_derivative checking
|
---|
| 109 | setDefaultParameters();
|
---|
| 110 | }
|
---|
| 111 |
|
---|
| 112 | PairPotential_LennardJones::PairPotential_LennardJones(
|
---|
| 113 | const ParticleTypes_t &_ParticleTypes,
|
---|
| 114 | const double _epsilon,
|
---|
| 115 | const double _sigma
|
---|
| 116 | ) :
|
---|
| 117 | EmpiricalPotential(_ParticleTypes),
|
---|
[d5ca1a] | 118 | params(parameters_t(MAXPARAMS, 0.)),
|
---|
| 119 | bindingmodel(generateBindingModel(_ParticleTypes))
|
---|
[1413f4] | 120 | {
|
---|
| 121 | params[epsilon] = _epsilon;
|
---|
| 122 | params[sigma] = _sigma;
|
---|
| 123 | }
|
---|
| 124 |
|
---|
| 125 | void PairPotential_LennardJones::setParameters(const parameters_t &_params)
|
---|
| 126 | {
|
---|
| 127 | const size_t paramsDim = _params.size();
|
---|
| 128 | ASSERT( paramsDim <= getParameterDimension(),
|
---|
| 129 | "PairPotential_LennardJones::setParameters() - we need not more than "
|
---|
| 130 | +toString(getParameterDimension())+" parameters.");
|
---|
| 131 | for(size_t i=0;i<paramsDim;++i)
|
---|
| 132 | params[i] = _params[i];
|
---|
| 133 |
|
---|
| 134 | #ifndef NDEBUG
|
---|
| 135 | parameters_t check_params(getParameters());
|
---|
| 136 | check_params.resize(paramsDim); // truncate to same size
|
---|
| 137 | ASSERT( check_params == _params,
|
---|
| 138 | "PairPotential_LennardJones::setParameters() - failed, mismatch in to be set "
|
---|
| 139 | +toString(_params)+" and set "+toString(check_params)+" params.");
|
---|
| 140 | #endif
|
---|
| 141 | }
|
---|
| 142 |
|
---|
| 143 | PairPotential_LennardJones::results_t
|
---|
| 144 | PairPotential_LennardJones::operator()(
|
---|
[e1fe7e] | 145 | const list_of_arguments_t &listarguments
|
---|
[1413f4] | 146 | ) const
|
---|
| 147 | {
|
---|
[e1fe7e] | 148 | result_t result = 0.;
|
---|
| 149 | for(list_of_arguments_t::const_iterator iter = listarguments.begin();
|
---|
| 150 | iter != listarguments.end(); ++iter) {
|
---|
| 151 | const arguments_t &arguments = *iter;
|
---|
| 152 | ASSERT( arguments.size() == 1,
|
---|
| 153 | "PairPotential_LennardJones::operator() - requires one argument.");
|
---|
| 154 | ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
|
---|
| 155 | arguments, getParticleTypes()),
|
---|
| 156 | "PairPotential_LennardJones::operator() - types don't match with ones in arguments.");
|
---|
| 157 | const double &r = arguments[0].distance;
|
---|
| 158 | const double temp = Helpers::pow(params[sigma]/r, 6);
|
---|
| 159 | result += 4.*params[epsilon] * (temp*temp - temp);
|
---|
| 160 | }
|
---|
| 161 | return results_t(1, result);
|
---|
[1413f4] | 162 | }
|
---|
| 163 |
|
---|
| 164 | PairPotential_LennardJones::derivative_components_t
|
---|
| 165 | PairPotential_LennardJones::derivative(
|
---|
[e1fe7e] | 166 | const list_of_arguments_t &listarguments
|
---|
[1413f4] | 167 | ) const
|
---|
| 168 | {
|
---|
| 169 | const double sigma6 = Helpers::pow(params[sigma], 6);
|
---|
[e1fe7e] | 170 | result_t result = 0.;
|
---|
| 171 | for(list_of_arguments_t::const_iterator iter = listarguments.begin();
|
---|
| 172 | iter != listarguments.end(); ++iter) {
|
---|
| 173 | const arguments_t &arguments = *iter;
|
---|
| 174 | ASSERT( arguments.size() == 1,
|
---|
| 175 | "PairPotential_LennardJones::operator() - requires no argument.");
|
---|
| 176 | ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
|
---|
| 177 | arguments, getParticleTypes()),
|
---|
| 178 | "PairPotential_LennardJones::operator() - types don't match with ones in arguments.");
|
---|
| 179 | const double &r = arguments[0].distance;
|
---|
| 180 | result +=
|
---|
| 181 | 4.*params[epsilon] * (
|
---|
| 182 | sigma6*sigma6*(-12.) / Helpers::pow(r,13)
|
---|
| 183 | - sigma6*(-6.) /Helpers::pow(r,7)
|
---|
| 184 | );
|
---|
| 185 | }
|
---|
| 186 | return derivative_components_t(1, result);
|
---|
[1413f4] | 187 | }
|
---|
| 188 |
|
---|
| 189 | PairPotential_LennardJones::results_t
|
---|
| 190 | PairPotential_LennardJones::parameter_derivative(
|
---|
[e1fe7e] | 191 | const list_of_arguments_t &listarguments,
|
---|
[1413f4] | 192 | const size_t index
|
---|
| 193 | ) const
|
---|
| 194 | {
|
---|
[e1fe7e] | 195 | result_t result = 0.;
|
---|
| 196 | for(list_of_arguments_t::const_iterator iter = listarguments.begin();
|
---|
| 197 | iter != listarguments.end(); ++iter) {
|
---|
| 198 | const arguments_t &arguments = *iter;
|
---|
| 199 | ASSERT( arguments.size() == 1,
|
---|
| 200 | "PairPotential_LennardJones::parameter_derivative() - requires no argument.");
|
---|
| 201 | ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
|
---|
| 202 | arguments, getParticleTypes()),
|
---|
| 203 | "PairPotential_LennardJones::operator() - types don't match with ones in arguments.");
|
---|
| 204 | const double &r = arguments[0].distance;
|
---|
| 205 | switch (index) {
|
---|
| 206 | case epsilon:
|
---|
| 207 | {
|
---|
| 208 | const double temp = Helpers::pow(params[sigma]/r, 6);
|
---|
| 209 | result += 4. * (temp*temp - temp);
|
---|
| 210 | break;
|
---|
| 211 | }
|
---|
| 212 | case sigma:
|
---|
| 213 | {
|
---|
| 214 | const double r6 = Helpers::pow(r, 6);
|
---|
| 215 | result +=
|
---|
| 216 | 4.*params[epsilon] * (
|
---|
| 217 | 12. * Helpers::pow(params[sigma],11)/(r6*r6)
|
---|
| 218 | - 6. * Helpers::pow(params[sigma],5)/r6
|
---|
| 219 | );
|
---|
| 220 | break;
|
---|
| 221 | }
|
---|
| 222 | default:
|
---|
| 223 | break;
|
---|
[1413f4] | 224 | }
|
---|
| 225 | }
|
---|
[e1fe7e] | 226 | return results_t(1, result);
|
---|
[1413f4] | 227 | }
|
---|
| 228 |
|
---|
| 229 | FunctionModel::filter_t PairPotential_LennardJones::getSpecificFilter() const
|
---|
| 230 | {
|
---|
| 231 | FunctionModel::filter_t returnfunction =
|
---|
| 232 | boost::bind(&Extractors::filterArgumentsByParticleTypes,
|
---|
[e60558] | 233 | _2, _1,
|
---|
[67044a] | 234 | boost::cref(getParticleTypes()), boost::cref(getBindingModel()));
|
---|
[1413f4] | 235 | return returnfunction;
|
---|
| 236 | }
|
---|
| 237 |
|
---|
| 238 | void
|
---|
| 239 | PairPotential_LennardJones::setParametersToRandomInitialValues(
|
---|
| 240 | const TrainingData &data)
|
---|
| 241 | {
|
---|
[0932c2] | 242 | RandomNumberGenerator &random = RandomNumberGeneratorFactory::getInstance().makeRandomNumberGenerator();
|
---|
| 243 | const double rng_min = random.min();
|
---|
| 244 | const double rng_max = random.max();
|
---|
| 245 | params[PairPotential_LennardJones::epsilon] = 1e-2*(random()/(rng_max-rng_min));
|
---|
| 246 | params[PairPotential_LennardJones::sigma] = (3.+10.*(random()/(rng_max-rng_min)));// 0.5;
|
---|
[1413f4] | 247 | }
|
---|
| 248 |
|
---|