source: src/atom_trajectoryparticle.cpp@ ddc85b

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 ddc85b was ddc85b, checked in by Tillmann Crueger <crueger@…>, 14 years ago

Added a method that allows calculation of the total temperature of a set of atoms at a given time step

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