source: src/atom_trajectoryparticle.cpp@ 81a9bc

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 81a9bc was bf3817, checked in by Frederik Heber <heber@…>, 14 years ago

Added ifdef HAVE_CONFIG and config.h include to each and every cpp file.

  • is now topmost in front of MemDebug.hpp (and any other).
  • Property mode set to 100644
File size: 10.4 KB
Line 
1/*
2 * atom_trajectoryparticle.cpp
3 *
4 * Created on: Oct 19, 2009
5 * Author: heber
6 */
7
8// include config.h
9#ifdef HAVE_CONFIG_H
10#include <config.h>
11#endif
12
13#include "Helpers/MemDebug.hpp"
14
15#include "atom.hpp"
16#include "atom_trajectoryparticle.hpp"
17#include "config.hpp"
18#include "element.hpp"
19#include "Helpers/Info.hpp"
20#include "Helpers/Log.hpp"
21#include "parser.hpp"
22#include "ThermoStatContainer.hpp"
23#include "Helpers/Verbose.hpp"
24
25/** Constructor of class TrajectoryParticle.
26 */
27TrajectoryParticle::TrajectoryParticle()
28{
29};
30
31/** Destructor of class TrajectoryParticle.
32 */
33TrajectoryParticle::~TrajectoryParticle()
34{
35};
36
37
38/** Adds kinetic energy of this atom to given temperature value.
39 * \param *temperature add on this value
40 * \param step given step of trajectory to add
41 */
42void TrajectoryParticle::AddKineticToTemperature(double *temperature, int step) const
43{
44 for (int i=NDIM;i--;)
45 *temperature += getType()->mass * Trajectory.U.at(step)[i]* Trajectory.U.at(step)[i];
46};
47
48/** Evaluates some constraint potential if atom moves from \a startstep at once to \endstep in trajectory.
49 * \param startstep trajectory begins at
50 * \param endstep trajectory ends at
51 * \param **PermutationMap if atom switches places with some other atom, there is no translation but a permutaton noted here (not in the trajectories of ea
52 * \param *Force Force matrix to store result in
53 */
54void TrajectoryParticle::EvaluateConstrainedForce(int startstep, int endstep, atom **PermutationMap, ForceMatrix *Force) const
55{
56 double constant = 10.;
57 TrajectoryParticle *Sprinter = PermutationMap[nr];
58 // set forces
59 for (int i=NDIM;i++;)
60 Force->Matrix[0][nr][5+i] += 2.*constant*sqrt(Trajectory.R.at(startstep).distance(Sprinter->Trajectory.R.at(endstep)));
61};
62
63/** Correct velocity against the summed \a CoGVelocity for \a step.
64 * \param *ActualTemp sum up actual temperature meanwhile
65 * \param Step MD step in atom::Tracjetory
66 * \param *CoGVelocity remnant velocity (i.e. vector sum of all atom velocities)
67 */
68void TrajectoryParticle::CorrectVelocity(double *ActualTemp, int Step, Vector *CoGVelocity)
69{
70 for(int d=0;d<NDIM;d++) {
71 Trajectory.U.at(Step)[d] -= CoGVelocity->at(d);
72 *ActualTemp += 0.5 * getType()->mass * Trajectory.U.at(Step)[d] * Trajectory.U.at(Step)[d];
73 }
74};
75
76/** Extends the trajectory STL vector to the new size.
77 * Does nothing if \a MaxSteps is smaller than current size.
78 * \param MaxSteps
79 */
80void TrajectoryParticle::ResizeTrajectory(int MaxSteps)
81{
82 Info FunctionInfo(__func__);
83 if (Trajectory.R.size() <= (unsigned int)(MaxSteps)) {
84 DoLog(0) && (Log() << Verbose(0) << "Increasing size for trajectory array of " << nr << " from " << Trajectory.R.size() << " to " << (MaxSteps+1) << "." << endl);
85 Trajectory.R.resize(MaxSteps+1);
86 Trajectory.U.resize(MaxSteps+1);
87 Trajectory.F.resize(MaxSteps+1);
88 }
89};
90
91/** Copies a given trajectory step \a src onto another \a dest
92 * \param dest index of destination step
93 * \param src index of source step
94 */
95void TrajectoryParticle::CopyStepOnStep(int dest, int src)
96{
97 if (dest == src) // self assignment check
98 return;
99
100 for (int n=NDIM;n--;) {
101 Trajectory.R.at(dest)[n] = Trajectory.R.at(src)[n];
102 Trajectory.U.at(dest)[n] = Trajectory.U.at(src)[n];
103 Trajectory.F.at(dest)[n] = Trajectory.F.at(src)[n];
104 }
105};
106
107/** Performs a velocity verlet update of the trajectory.
108 * Parameters are according to those in configuration class.
109 * \param NextStep index of sequential step to set
110 * \param *configuration pointer to configuration with parameters
111 * \param *Force matrix with forces
112 */
113void TrajectoryParticle::VelocityVerletUpdate(int NextStep, config *configuration, ForceMatrix *Force, const size_t offset)
114{
115 //a = configuration.Deltat*0.5/walker->type->mass; // (F+F_old)/2m = a and thus: v = (F+F_old)/2m * t = (F + F_old) * a
116 for (int d=0; d<NDIM; d++) {
117 Trajectory.F.at(NextStep)[d] = -Force->Matrix[0][nr][d+offset]*(configuration->GetIsAngstroem() ? AtomicLengthToAngstroem : 1.);
118 Trajectory.R.at(NextStep)[d] = Trajectory.R.at(NextStep-1)[d];
119 Trajectory.R.at(NextStep)[d] += configuration->Deltat*(Trajectory.U.at(NextStep-1)[d]); // s(t) = s(0) + v * deltat + 1/2 a * deltat^2
120 Trajectory.R.at(NextStep)[d] += 0.5*configuration->Deltat*configuration->Deltat*(Trajectory.F.at(NextStep)[d]/getType()->mass); // F = m * a and s =
121 }
122 // Update U
123 for (int d=0; d<NDIM; d++) {
124 Trajectory.U.at(NextStep)[d] = Trajectory.U.at(NextStep-1)[d];
125 Trajectory.U.at(NextStep)[d] += configuration->Deltat * (Trajectory.F.at(NextStep)[d]+Trajectory.F.at(NextStep-1)[d]/getType()->mass); // v = F/m * t
126 }
127 // Update R (and F)
128// out << "Integrated position&velocity of step " << (NextStep) << ": (";
129// for (int d=0;d<NDIM;d++)
130// out << Trajectory.R.at(NextStep).x[d] << " "; // next step
131// out << ")\t(";
132// for (int d=0;d<NDIM;d++)
133// Log() << Verbose(0) << Trajectory.U.at(NextStep).x[d] << " "; // next step
134// out << ")" << endl;
135};
136
137/** Sums up mass and kinetics.
138 * \param Step step to sum for
139 * \param *TotalMass pointer to total mass sum
140 * \param *TotalVelocity pointer to tota velocity sum
141 */
142void TrajectoryParticle::SumUpKineticEnergy( int Step, double *TotalMass, Vector *TotalVelocity ) const
143{
144 *TotalMass += getType()->mass; // sum up total mass
145 for(int d=0;d<NDIM;d++) {
146 TotalVelocity->at(d) += Trajectory.U.at(Step)[d]*getType()->mass;
147 }
148};
149
150/** Scales velocity of atom according to Woodcock thermostat.
151 * \param ScaleTempFactor factor to scale the velocities with (i.e. sqrt of energy scale factor)
152 * \param Step MD step to scale
153 * \param *ekin sum of kinetic energy
154 */
155void TrajectoryParticle::Thermostat_Woodcock(double ScaleTempFactor, int Step, double *ekin)
156{
157 Vector &U = Trajectory.U.at(Step);
158 if (FixedIon == 0) // even FixedIon moves, only not by other's forces
159 for (int d=0; d<NDIM; d++) {
160 U[d] *= ScaleTempFactor;
161 *ekin += 0.5*getType()->mass * U[d]*U[d];
162 }
163};
164
165/** Scales velocity of atom according to Gaussian thermostat.
166 * \param Step MD step to scale
167 * \param *G
168 * \param *E
169 */
170void TrajectoryParticle::Thermostat_Gaussian_init(int Step, double *G, double *E)
171{
172 Vector &U = Trajectory.U.at(Step);
173 Vector &F = Trajectory.F.at(Step);
174 if (FixedIon == 0) // even FixedIon moves, only not by other's forces
175 for (int d=0; d<NDIM; d++) {
176 *G += U[d] * F[d];
177 *E += U[d]*U[d]*getType()->mass;
178 }
179};
180
181/** Determines scale factors according to Gaussian thermostat.
182 * \param Step MD step to scale
183 * \param GE G over E ratio
184 * \param *ekin sum of kinetic energy
185 * \param *configuration configuration class with TempFrequency and TargetTemp
186 */
187void TrajectoryParticle::Thermostat_Gaussian_least_constraint(int Step, double G_over_E, double *ekin, config *configuration)
188{
189 Vector &U = Trajectory.U.at(Step);
190 if (FixedIon == 0) // even FixedIon moves, only not by other's forces
191 for (int d=0; d<NDIM; d++) {
192 U[d] += configuration->Deltat/getType()->mass * ( (G_over_E) * (U[d]*getType()->mass) );
193 *ekin += getType()->mass * U[d]*U[d];
194 }
195};
196
197/** Scales velocity of atom according to Langevin thermostat.
198 * \param Step MD step to scale
199 * \param *r random number generator
200 * \param *ekin sum of kinetic energy
201 * \param *configuration configuration class with TempFrequency and TargetTemp
202 */
203void TrajectoryParticle::Thermostat_Langevin(int Step, gsl_rng * r, double *ekin, config *configuration)
204{
205 double sigma = sqrt(configuration->Thermostats->TargetTemp/getType()->mass); // sigma = (k_b T)/m (Hartree/atomicmass = atomiclength/atomictime)
206 Vector &U = Trajectory.U.at(Step);
207 if (FixedIon == 0) { // even FixedIon moves, only not by other's forces
208 // throw a dice to determine whether it gets hit by a heat bath particle
209 if (((((rand()/(double)RAND_MAX))*configuration->Thermostats->TempFrequency) < 1.)) {
210 DoLog(3) && (Log() << Verbose(3) << "Particle " << *this << " was hit (sigma " << sigma << "): " << sqrt(U[0]*U[0]+U[1]*U[1]+U[2]*U[2]) << " -> ");
211 // pick three random numbers from a Boltzmann distribution around the desired temperature T for each momenta axis
212 for (int d=0; d<NDIM; d++) {
213 U[d] = gsl_ran_gaussian (r, sigma);
214 }
215 DoLog(2) && (Log() << Verbose(2) << sqrt(U[0]*U[0]+U[1]*U[1]+U[2]*U[2]) << endl);
216 }
217 for (int d=0; d<NDIM; d++)
218 *ekin += 0.5*getType()->mass * U[d]*U[d];
219 }
220};
221
222/** Scales velocity of atom according to Berendsen thermostat.
223 * \param Step MD step to scale
224 * \param ScaleTempFactor factor to scale energy (not velocity!) with
225 * \param *ekin sum of kinetic energy
226 * \param *configuration configuration class with TempFrequency and Deltat
227 */
228void TrajectoryParticle::Thermostat_Berendsen(int Step, double ScaleTempFactor, double *ekin, config *configuration)
229{
230 Vector &U = Trajectory.U.at(Step);
231 if (FixedIon == 0) { // even FixedIon moves, only not by other's forces
232 for (int d=0; d<NDIM; d++) {
233 U[d] *= sqrt(1+(configuration->Deltat/configuration->Thermostats->TempFrequency)*(ScaleTempFactor-1));
234 *ekin += 0.5*getType()->mass * U[d]*U[d];
235 }
236 }
237};
238
239/** Initializes current run of NoseHoover thermostat.
240 * \param Step MD step to scale
241 * \param *delta_alpha additional sum of kinetic energy on return
242 */
243void TrajectoryParticle::Thermostat_NoseHoover_init(int Step, double *delta_alpha)
244{
245 Vector &U = Trajectory.U.at(Step);
246 if (FixedIon == 0) { // even FixedIon moves, only not by other's forces
247 for (int d=0; d<NDIM; d++) {
248 *delta_alpha += U[d]*U[d]*getType()->mass;
249 }
250 }
251};
252
253/** Initializes current run of NoseHoover thermostat.
254 * \param Step MD step to scale
255 * \param *ekin sum of kinetic energy
256 * \param *configuration configuration class with TempFrequency and Deltat
257 */
258void TrajectoryParticle::Thermostat_NoseHoover_scale(int Step, double *ekin, config *configuration)
259{
260 Vector &U = Trajectory.U.at(Step);
261 if (FixedIon == 0) { // even FixedIon moves, only not by other's forces
262 for (int d=0; d<NDIM; d++) {
263 U[d] += configuration->Deltat/getType()->mass * (configuration->Thermostats->alpha * (U[d] * getType()->mass));
264 *ekin += (0.5*getType()->mass) * U[d]*U[d];
265 }
266 }
267};
268
269
270std::ostream & TrajectoryParticle::operator << (std::ostream &ost) const
271{
272 ParticleInfo::operator<<(ost);
273 ost << "," << getPosition();
274 return ost;
275}
276
277std::ostream & operator << (std::ostream &ost, const TrajectoryParticle &a)
278{
279 a.ParticleInfo::operator<<(ost);
280 ost << "," << a.getPosition();
281 return ost;
282}
283
Note: See TracBrowser for help on using the repository browser.