source: src/atom_trajectoryparticle.cpp@ fc1b24

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 fc1b24 was 717e0c, checked in by Frederik Heber <heber@…>, 15 years ago

Verbosity corrected for ERROR and WARNING

  • present ERROR and WARNING prefixes removed and placed by eLog() and respective Verbosity().
  • -v... is scanned for number of 'v's and verbosity is set accordingly
  • standard verbosity is now 0.

Signed-off-by: Frederik Heber <heber@…>

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