source: src/Potentials/Specifics/PairPotential_Angle.cpp@ b15ae7

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 b15ae7 was 1e242a, checked in by Frederik Heber <heber@…>, 11 years ago

Removed energy_offset from PairPotential_Angle.

  • Property mode set to 100644
File size: 7.2 KB
Line 
1/*
2 * Project: MoleCuilder
3 * Description: creates and alters molecular systems
4 * Copyright (C) 2012 University of Bonn. 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_Angle.cpp
26 *
27 * Created on: Oct 11, 2012
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_Angle.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 <string>
45
46#include "CodePatterns/Assert.hpp"
47
48#include "FunctionApproximation/Extractors.hpp"
49#include "FunctionApproximation/TrainingData.hpp"
50#include "Potentials/helpers.hpp"
51#include "Potentials/ParticleTypeCheckers.hpp"
52
53class Fragment;
54
55// static definitions
56const PairPotential_Angle::ParameterNames_t
57PairPotential_Angle::ParameterNames =
58 boost::assign::list_of<std::string>
59 ("spring_constant")
60 ("equilibrium_distance")
61 ;
62const std::string PairPotential_Angle::potential_token("harmonic_angle");
63
64PairPotential_Angle::PairPotential_Angle(
65 const ParticleTypes_t &_ParticleTypes
66 ) :
67 EmpiricalPotential(_ParticleTypes),
68 params(parameters_t(MAXPARAMS, 0.))
69{
70 // have some decent defaults for parameter_derivative checking
71 params[spring_constant] = 1.;
72 params[equilibrium_distance] = 0.1;
73}
74
75PairPotential_Angle::PairPotential_Angle(
76 const ParticleTypes_t &_ParticleTypes,
77 const double _spring_constant,
78 const double _equilibrium_distance) :
79 EmpiricalPotential(_ParticleTypes),
80 params(parameters_t(MAXPARAMS, 0.))
81{
82 params[spring_constant] = _spring_constant;
83 params[equilibrium_distance] = _equilibrium_distance;
84}
85
86void PairPotential_Angle::setParameters(const parameters_t &_params)
87{
88 const size_t paramsDim = _params.size();
89 ASSERT( paramsDim <= getParameterDimension(),
90 "PairPotential_Angle::setParameters() - we need not more than "
91 +toString(getParameterDimension())+" parameters.");
92 for(size_t i=0;i<paramsDim;++i)
93 params[i] = _params[i];
94
95#ifndef NDEBUG
96 parameters_t check_params(getParameters());
97 check_params.resize(paramsDim); // truncate to same size
98 ASSERT( check_params == _params,
99 "PairPotential_Angle::setParameters() - failed, mismatch in to be set "
100 +toString(_params)+" and set "+toString(check_params)+" params.");
101#endif
102}
103
104PairPotential_Angle::result_t
105PairPotential_Angle::function_theta(
106 const double &r_ij,
107 const double &r_jk,
108 const double &r_ik
109 ) const
110{
111// Info info(__func__);
112 const double angle = Helpers::pow(r_ij,2) + Helpers::pow(r_jk,2) - Helpers::pow(r_ik,2);
113 const double divisor = 2.* r_ij * r_jk;
114
115// LOG(2, "DEBUG: cos(theta)= " << angle/divisor);
116 if (divisor == 0.)
117 return 0.;
118 else
119 return angle/divisor;
120}
121
122PairPotential_Angle::results_t
123PairPotential_Angle::operator()(
124 const arguments_t &arguments
125 ) const
126{
127 ASSERT( arguments.size() == 3,
128 "PairPotential_Angle::operator() - requires exactly three arguments.");
129 ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
130 arguments, getParticleTypes()),
131 "PairPotential_Angle::operator() - types don't match with ones in arguments.");
132 const argument_t &r_ij = arguments[0]; // 01
133 const argument_t &r_jk = arguments[2]; // 12
134 const argument_t &r_ik = arguments[1]; // 02
135 const result_t result =
136 params[spring_constant]
137 * Helpers::pow( function_theta(r_ij.distance, r_jk.distance, r_ik.distance) - params[equilibrium_distance], 2 );
138 return std::vector<result_t>(1, result);
139}
140
141PairPotential_Angle::derivative_components_t
142PairPotential_Angle::derivative(
143 const arguments_t &arguments
144 ) const
145{
146 ASSERT( arguments.size() == 3,
147 "PairPotential_Angle::operator() - requires exactly three arguments.");
148 ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
149 arguments, getParticleTypes()),
150 "PairPotential_Angle::operator() - types don't match with ones in arguments.");
151 derivative_components_t result;
152 const argument_t &r_ij = arguments[0]; //01
153 const argument_t &r_jk = arguments[2]; //12
154 const argument_t &r_ik = arguments[1]; //02
155 result.push_back( 2. * params[spring_constant] * ( function_theta(r_ij.distance, r_jk.distance, r_ik.distance) - params[equilibrium_distance]) );
156 ASSERT( result.size() == 1,
157 "PairPotential_Angle::operator() - we did not create exactly one component.");
158 return result;
159}
160
161PairPotential_Angle::results_t
162PairPotential_Angle::parameter_derivative(
163 const arguments_t &arguments,
164 const size_t index
165 ) const
166{
167 ASSERT( arguments.size() == 3,
168 "PairPotential_Angle::parameter_derivative() - requires exactly three arguments.");
169 ASSERT( ParticleTypeChecker::checkArgumentsAgainstParticleTypes(
170 arguments, getParticleTypes()),
171 "PairPotential_Angle::operator() - types don't match with ones in arguments.");
172 const argument_t &r_ij = arguments[0]; //01
173 const argument_t &r_jk = arguments[2]; //12
174 const argument_t &r_ik = arguments[1]; //02
175 switch (index) {
176 case spring_constant:
177 {
178 const result_t result =
179 Helpers::pow( function_theta(r_ij.distance, r_jk.distance, r_ik.distance) - params[equilibrium_distance], 2 );
180 return std::vector<result_t>(1, result);
181 break;
182 }
183 case equilibrium_distance:
184 {
185 const result_t result =
186 -2. * params[spring_constant]
187 * ( function_theta(r_ij.distance, r_jk.distance, r_ik.distance) - params[equilibrium_distance]);
188 return std::vector<result_t>(1, result);
189 break;
190 }
191 default:
192 ASSERT(0, "PairPotential_Angle::parameter_derivative() - derivative to unknown parameter desired.");
193 break;
194 }
195}
196
197FunctionModel::extractor_t
198PairPotential_Angle::getFragmentSpecificExtractor() const
199{
200 Fragment::charges_t charges;
201 charges.resize(getParticleTypes().size());
202 std::transform(getParticleTypes().begin(), getParticleTypes().end(),
203 charges.begin(), boost::lambda::_1);
204 FunctionModel::extractor_t returnfunction =
205 boost::bind(&Extractors::gatherDistancesFromFragment,
206 boost::bind(&Fragment::getPositions, _1),
207 boost::bind(&Fragment::getCharges, _1),
208 charges,
209 _2);
210 return returnfunction;
211}
212
213void
214PairPotential_Angle::setParametersToRandomInitialValues(
215 const TrainingData &data)
216{
217 params[PairPotential_Angle::spring_constant] = 1e+0*rand()/(double)RAND_MAX;// 0.2;
218 params[PairPotential_Angle::equilibrium_distance] = -0.3;//2e+0*rand()/(double)RAND_MAX - 1.;// 1.;
219}
220
Note: See TracBrowser for help on using the repository browser.