source: molecuilder/src/molecule_dynamics.cpp@ 17b3a5c

Last change on this file since 17b3a5c was 17b3a5c, checked in by Frederik Heber <heber@…>, 16 years ago

forward declarations used to untangle interdependet classes.

  • basically, everywhere in header files we removed '#include' lines were only pointer to the respective classes were used and the include line was moved to the implementation file.
  • as a sidenote, lots of funny errors happened because headers were included via a nesting over three other includes. Now, all should be declared directly as needed, as only very little include lines remain in header files.
  • Property mode set to 100644
File size: 38.8 KB
Line 
1/*
2 * molecule_dynamics.cpp
3 *
4 * Created on: Oct 5, 2009
5 * Author: heber
6 */
7
8#include "atom.hpp"
9#include "config.hpp"
10#include "element.hpp"
11#include "memoryallocator.hpp"
12#include "molecule.hpp"
13#include "parser.hpp"
14
15/************************************* Functions for class molecule *********************************/
16
17
18/** Evaluates the potential energy used for constrained molecular dynamics.
19 * \f$V_i^{con} = c^{bond} \cdot | r_{P(i)} - R_i | + sum_{i \neq j} C^{min} \cdot \frac{1}{C_{ij}} + C^{inj} \Bigl (1 - \theta \bigl (\prod_{i \neq j} (P(i) - P(j)) \bigr ) \Bigr )\f$
20 * where the first term points to the target in minimum distance, the second is a penalty for trajectories lying too close to each other (\f$C_{ij}\f$ is minimum distance between
21 * trajectories i and j) and the third term is a penalty for two atoms trying to each the same target point.
22 * Note that for the second term we have to solve the following linear system:
23 * \f$-c_1 \cdot n_1 + c_2 \cdot n_2 + C \cdot n_3 = - p_2 + p_1\f$, where \f$c_1\f$, \f$c_2\f$ and \f$C\f$ are constants,
24 * offset vector \f$p_1\f$ in direction \f$n_1\f$, offset vector \f$p_2\f$ in direction \f$n_2\f$,
25 * \f$n_3\f$ is the normal vector to both directions. \f$C\f$ would be the minimum distance between the two lines.
26 * \sa molecule::MinimiseConstrainedPotential(), molecule::VerletForceIntegration()
27 * \param *out output stream for debugging
28 * \param *PermutationMap gives target ptr for each atom, array of size molecule::AtomCount (this is "x" in \f$ V^{con}(x) \f$ )
29 * \param startstep start configuration (MDStep in molecule::trajectories)
30 * \param endstep end configuration (MDStep in molecule::trajectories)
31 * \param *constants constant in front of each term
32 * \param IsAngstroem whether coordinates are in angstroem (true) or bohrradius (false)
33 * \return potential energy
34 * \note This routine is scaling quadratically which is not optimal.
35 * \todo There's a bit double counting going on for the first time, bu nothing to worry really about.
36 */
37double molecule::ConstrainedPotential(ofstream *out, atom **PermutationMap, int startstep, int endstep, double *constants, bool IsAngstroem)
38{
39 double result = 0., tmp, Norm1, Norm2;
40 atom *Walker = NULL, *Runner = NULL, *Sprinter = NULL;
41 Vector trajectory1, trajectory2, normal, TestVector;
42 gsl_matrix *A = gsl_matrix_alloc(NDIM,NDIM);
43 gsl_vector *x = gsl_vector_alloc(NDIM);
44
45 // go through every atom
46 Walker = start;
47 while (Walker->next != end) {
48 Walker = Walker->next;
49 // first term: distance to target
50 Runner = PermutationMap[Walker->nr]; // find target point
51 tmp = (Walker->Trajectory.R.at(startstep).Distance(&Runner->Trajectory.R.at(endstep)));
52 tmp *= IsAngstroem ? 1. : 1./AtomicLengthToAngstroem;
53 result += constants[0] * tmp;
54 //*out << Verbose(4) << "Adding " << tmp*constants[0] << "." << endl;
55
56 // second term: sum of distances to other trajectories
57 Runner = start;
58 while (Runner->next != end) {
59 Runner = Runner->next;
60 if (Runner == Walker) // hence, we only go up to the Walker, not beyond (similar to i=0; i<j; i++)
61 break;
62 // determine normalized trajectories direction vector (n1, n2)
63 Sprinter = PermutationMap[Walker->nr]; // find first target point
64 trajectory1.CopyVector(&Sprinter->Trajectory.R.at(endstep));
65 trajectory1.SubtractVector(&Walker->Trajectory.R.at(startstep));
66 trajectory1.Normalize();
67 Norm1 = trajectory1.Norm();
68 Sprinter = PermutationMap[Runner->nr]; // find second target point
69 trajectory2.CopyVector(&Sprinter->Trajectory.R.at(endstep));
70 trajectory2.SubtractVector(&Runner->Trajectory.R.at(startstep));
71 trajectory2.Normalize();
72 Norm2 = trajectory1.Norm();
73 // check whether either is zero()
74 if ((Norm1 < MYEPSILON) && (Norm2 < MYEPSILON)) {
75 tmp = Walker->Trajectory.R.at(startstep).Distance(&Runner->Trajectory.R.at(startstep));
76 } else if (Norm1 < MYEPSILON) {
77 Sprinter = PermutationMap[Walker->nr]; // find first target point
78 trajectory1.CopyVector(&Sprinter->Trajectory.R.at(endstep)); // copy first offset
79 trajectory1.SubtractVector(&Runner->Trajectory.R.at(startstep)); // subtract second offset
80 trajectory2.Scale( trajectory1.ScalarProduct(&trajectory2) ); // trajectory2 is scaled to unity, hence we don't need to divide by anything
81 trajectory1.SubtractVector(&trajectory2); // project the part in norm direction away
82 tmp = trajectory1.Norm(); // remaining norm is distance
83 } else if (Norm2 < MYEPSILON) {
84 Sprinter = PermutationMap[Runner->nr]; // find second target point
85 trajectory2.CopyVector(&Sprinter->Trajectory.R.at(endstep)); // copy second offset
86 trajectory2.SubtractVector(&Walker->Trajectory.R.at(startstep)); // subtract first offset
87 trajectory1.Scale( trajectory2.ScalarProduct(&trajectory1) ); // trajectory1 is scaled to unity, hence we don't need to divide by anything
88 trajectory2.SubtractVector(&trajectory1); // project the part in norm direction away
89 tmp = trajectory2.Norm(); // remaining norm is distance
90 } else if ((fabs(trajectory1.ScalarProduct(&trajectory2)/Norm1/Norm2) - 1.) < MYEPSILON) { // check whether they're linear dependent
91// *out << Verbose(3) << "Both trajectories of " << *Walker << " and " << *Runner << " are linear dependent: ";
92// *out << trajectory1;
93// *out << " and ";
94// *out << trajectory2;
95 tmp = Walker->Trajectory.R.at(startstep).Distance(&Runner->Trajectory.R.at(startstep));
96// *out << " with distance " << tmp << "." << endl;
97 } else { // determine distance by finding minimum distance
98// *out << Verbose(3) << "Both trajectories of " << *Walker << " and " << *Runner << " are linear independent ";
99// *out << endl;
100// *out << "First Trajectory: ";
101// *out << trajectory1 << endl;
102// *out << "Second Trajectory: ";
103// *out << trajectory2 << endl;
104 // determine normal vector for both
105 normal.MakeNormalVector(&trajectory1, &trajectory2);
106 // print all vectors for debugging
107// *out << "Normal vector in between: ";
108// *out << normal << endl;
109 // setup matrix
110 for (int i=NDIM;i--;) {
111 gsl_matrix_set(A, 0, i, trajectory1.x[i]);
112 gsl_matrix_set(A, 1, i, trajectory2.x[i]);
113 gsl_matrix_set(A, 2, i, normal.x[i]);
114 gsl_vector_set(x,i, (Walker->Trajectory.R.at(startstep).x[i] - Runner->Trajectory.R.at(startstep).x[i]));
115 }
116 // solve the linear system by Householder transformations
117 gsl_linalg_HH_svx(A, x);
118 // distance from last component
119 tmp = gsl_vector_get(x,2);
120// *out << " with distance " << tmp << "." << endl;
121 // test whether we really have the intersection (by checking on c_1 and c_2)
122 TestVector.CopyVector(&Runner->Trajectory.R.at(startstep));
123 trajectory2.Scale(gsl_vector_get(x,1));
124 TestVector.AddVector(&trajectory2);
125 normal.Scale(gsl_vector_get(x,2));
126 TestVector.AddVector(&normal);
127 TestVector.SubtractVector(&Walker->Trajectory.R.at(startstep));
128 trajectory1.Scale(gsl_vector_get(x,0));
129 TestVector.SubtractVector(&trajectory1);
130 if (TestVector.Norm() < MYEPSILON) {
131// *out << Verbose(2) << "Test: ok.\tDistance of " << tmp << " is correct." << endl;
132 } else {
133// *out << Verbose(2) << "Test: failed.\tIntersection is off by ";
134// *out << TestVector;
135// *out << "." << endl;
136 }
137 }
138 // add up
139 tmp *= IsAngstroem ? 1. : 1./AtomicLengthToAngstroem;
140 if (fabs(tmp) > MYEPSILON) {
141 result += constants[1] * 1./tmp;
142 //*out << Verbose(4) << "Adding " << 1./tmp*constants[1] << "." << endl;
143 }
144 }
145
146 // third term: penalty for equal targets
147 Runner = start;
148 while (Runner->next != end) {
149 Runner = Runner->next;
150 if ((PermutationMap[Walker->nr] == PermutationMap[Runner->nr]) && (Walker->nr < Runner->nr)) {
151 Sprinter = PermutationMap[Walker->nr];
152// *out << *Walker << " and " << *Runner << " are heading to the same target at ";
153// *out << Sprinter->Trajectory.R.at(endstep);
154// *out << ", penalting." << endl;
155 result += constants[2];
156 //*out << Verbose(4) << "Adding " << constants[2] << "." << endl;
157 }
158 }
159 }
160
161 return result;
162};
163
164void PrintPermutationMap(ofstream *out, atom **PermutationMap, int Nr)
165{
166 stringstream zeile1, zeile2;
167 int *DoubleList = Malloc<int>(Nr, "PrintPermutationMap: *DoubleList");
168 int doubles = 0;
169 for (int i=0;i<Nr;i++)
170 DoubleList[i] = 0;
171 zeile1 << "PermutationMap: ";
172 zeile2 << " ";
173 for (int i=0;i<Nr;i++) {
174 DoubleList[PermutationMap[i]->nr]++;
175 zeile1 << i << " ";
176 zeile2 << PermutationMap[i]->nr << " ";
177 }
178 for (int i=0;i<Nr;i++)
179 if (DoubleList[i] > 1)
180 doubles++;
181 // *out << "Found " << doubles << " Doubles." << endl;
182 Free(&DoubleList);
183// *out << zeile1.str() << endl << zeile2.str() << endl;
184};
185
186/** Minimises the extra potential for constrained molecular dynamics and gives forces and the constrained potential energy.
187 * We do the following:
188 * -# Generate a distance list from all source to all target points
189 * -# Sort this per source point
190 * -# Take for each source point the target point with minimum distance, use this as initial permutation
191 * -# check whether molecule::ConstrainedPotential() is greater than injective penalty
192 * -# If so, we go through each source point, stepping down in the sorted target point distance list and re-checking potential.
193 * -# Next, we only apply transformations that keep the injectivity of the permutations list.
194 * -# Hence, for one source point we step down the ladder and seek the corresponding owner of this new target
195 * point and try to change it for one with lesser distance, or for the next one with greater distance, but only
196 * if this decreases the conditional potential.
197 * -# finished.
198 * -# Then, we calculate the forces by taking the spatial derivative, where we scale the potential to such a degree,
199 * that the total force is always pointing in direction of the constraint force (ensuring that we move in the
200 * right direction).
201 * -# Finally, we calculate the potential energy and return.
202 * \param *out output stream for debugging
203 * \param **PermutationMap on return: mapping between the atom label of the initial and the final configuration
204 * \param startstep current MD step giving initial position between which and \a endstep we perform the constrained MD (as further steps are always concatenated)
205 * \param endstep step giving final position in constrained MD
206 * \param IsAngstroem whether coordinates are in angstroem (true) or bohrradius (false)
207 * \sa molecule::VerletForceIntegration()
208 * \return potential energy (and allocated **PermutationMap (array of molecule::AtomCount ^2)
209 * \todo The constrained potential's constants are set to fixed values right now, but they should scale based on checks of the system in order
210 * to ensure they're properties (e.g. constants[2] always greater than the energy of the system).
211 * \bug this all is not O(N log N) but O(N^2)
212 */
213double molecule::MinimiseConstrainedPotential(ofstream *out, atom **&PermutationMap, int startstep, int endstep, bool IsAngstroem)
214{
215 double Potential, OldPotential, OlderPotential;
216 PermutationMap = Malloc<atom*>(AtomCount, "molecule::MinimiseConstrainedPotential: **PermutationMap");
217 DistanceMap **DistanceList = Malloc<DistanceMap*>(AtomCount, "molecule::MinimiseConstrainedPotential: **DistanceList");
218 DistanceMap::iterator *DistanceIterators = Malloc<DistanceMap::iterator>(AtomCount, "molecule::MinimiseConstrainedPotential: *DistanceIterators");
219 int *DoubleList = Malloc<int>(AtomCount, "molecule::MinimiseConstrainedPotential: *DoubleList");
220 DistanceMap::iterator *StepList = Malloc<DistanceMap::iterator>(AtomCount, "molecule::MinimiseConstrainedPotential: *StepList");
221 double constants[3];
222 int round;
223 atom *Walker = NULL, *Runner = NULL, *Sprinter = NULL;
224 DistanceMap::iterator Rider, Strider;
225
226 /// Minimise the potential
227 // set Lagrange multiplier constants
228 constants[0] = 10.;
229 constants[1] = 1.;
230 constants[2] = 1e+7; // just a huge penalty
231 // generate the distance list
232 *out << Verbose(1) << "Creating the distance list ... " << endl;
233 for (int i=AtomCount; i--;) {
234 DoubleList[i] = 0; // stores for how many atoms in startstep this atom is a target in endstep
235 DistanceList[i] = new DistanceMap; // is the distance sorted target list per atom
236 DistanceList[i]->clear();
237 }
238 *out << Verbose(1) << "Filling the distance list ... " << endl;
239 Walker = start;
240 while (Walker->next != end) {
241 Walker = Walker->next;
242 Runner = start;
243 while(Runner->next != end) {
244 Runner = Runner->next;
245 DistanceList[Walker->nr]->insert( DistancePair(Walker->Trajectory.R.at(startstep).Distance(&Runner->Trajectory.R.at(endstep)), Runner) );
246 }
247 }
248 // create the initial PermutationMap (source -> target)
249 Walker = start;
250 while (Walker->next != end) {
251 Walker = Walker->next;
252 StepList[Walker->nr] = DistanceList[Walker->nr]->begin(); // stores the step to the next iterator that could be a possible next target
253 PermutationMap[Walker->nr] = DistanceList[Walker->nr]->begin()->second; // always pick target with the smallest distance
254 DoubleList[DistanceList[Walker->nr]->begin()->second->nr]++; // increase this target's source count (>1? not injective)
255 DistanceIterators[Walker->nr] = DistanceList[Walker->nr]->begin(); // and remember which one we picked
256 *out << *Walker << " starts with distance " << DistanceList[Walker->nr]->begin()->first << "." << endl;
257 }
258 *out << Verbose(1) << "done." << endl;
259 // make the PermutationMap injective by checking whether we have a non-zero constants[2] term in it
260 *out << Verbose(1) << "Making the PermutationMap injective ... " << endl;
261 Walker = start;
262 DistanceMap::iterator NewBase;
263 OldPotential = fabs(ConstrainedPotential(out, PermutationMap, startstep, endstep, constants, IsAngstroem));
264 while ((OldPotential) > constants[2]) {
265 PrintPermutationMap(out, PermutationMap, AtomCount);
266 Walker = Walker->next;
267 if (Walker == end) // round-robin at the end
268 Walker = start->next;
269 if (DoubleList[DistanceIterators[Walker->nr]->second->nr] <= 1) // no need to make those injective that aren't
270 continue;
271 // now, try finding a new one
272 NewBase = DistanceIterators[Walker->nr]; // store old base
273 do {
274 NewBase++; // take next further distance in distance to targets list that's a target of no one
275 } while ((DoubleList[NewBase->second->nr] != 0) && (NewBase != DistanceList[Walker->nr]->end()));
276 if (NewBase != DistanceList[Walker->nr]->end()) {
277 PermutationMap[Walker->nr] = NewBase->second;
278 Potential = fabs(ConstrainedPotential(out, PermutationMap, startstep, endstep, constants, IsAngstroem));
279 if (Potential > OldPotential) { // undo
280 PermutationMap[Walker->nr] = DistanceIterators[Walker->nr]->second;
281 } else { // do
282 DoubleList[DistanceIterators[Walker->nr]->second->nr]--; // decrease the old entry in the doubles list
283 DoubleList[NewBase->second->nr]++; // increase the old entry in the doubles list
284 DistanceIterators[Walker->nr] = NewBase;
285 OldPotential = Potential;
286 *out << Verbose(3) << "Found a new permutation, new potential is " << OldPotential << "." << endl;
287 }
288 }
289 }
290 for (int i=AtomCount; i--;) // now each single entry in the DoubleList should be <=1
291 if (DoubleList[i] > 1) {
292 cerr << "Failed to create an injective PermutationMap!" << endl;
293 exit(1);
294 }
295 *out << Verbose(1) << "done." << endl;
296 Free(&DoubleList);
297 // argument minimise the constrained potential in this injective PermutationMap
298 *out << Verbose(1) << "Argument minimising the PermutationMap, at current potential " << OldPotential << " ... " << endl;
299 OldPotential = 1e+10;
300 round = 0;
301 do {
302 *out << "Starting round " << ++round << " ... " << endl;
303 OlderPotential = OldPotential;
304 do {
305 Walker = start;
306 while (Walker->next != end) { // pick one
307 Walker = Walker->next;
308 PrintPermutationMap(out, PermutationMap, AtomCount);
309 Sprinter = DistanceIterators[Walker->nr]->second; // store initial partner
310 Strider = DistanceIterators[Walker->nr]; //remember old iterator
311 DistanceIterators[Walker->nr] = StepList[Walker->nr];
312 if (DistanceIterators[Walker->nr] == DistanceList[Walker->nr]->end()) {// stop, before we run through the list and still on
313 DistanceIterators[Walker->nr] == DistanceList[Walker->nr]->begin();
314 break;
315 }
316 //*out << Verbose(2) << "Current Walker: " << *Walker << " with old/next candidate " << *Sprinter << "/" << *DistanceIterators[Walker->nr]->second << "." << endl;
317 // find source of the new target
318 Runner = start->next;
319 while(Runner != end) { // find the source whose toes we might be stepping on (Walker's new target should be in use by another already)
320 if (PermutationMap[Runner->nr] == DistanceIterators[Walker->nr]->second) {
321 //*out << Verbose(2) << "Found the corresponding owner " << *Runner << " to " << *PermutationMap[Runner->nr] << "." << endl;
322 break;
323 }
324 Runner = Runner->next;
325 }
326 if (Runner != end) { // we found the other source
327 // then look in its distance list for Sprinter
328 Rider = DistanceList[Runner->nr]->begin();
329 for (; Rider != DistanceList[Runner->nr]->end(); Rider++)
330 if (Rider->second == Sprinter)
331 break;
332 if (Rider != DistanceList[Runner->nr]->end()) { // if we have found one
333 //*out << Verbose(2) << "Current Other: " << *Runner << " with old/next candidate " << *PermutationMap[Runner->nr] << "/" << *Rider->second << "." << endl;
334 // exchange both
335 PermutationMap[Walker->nr] = DistanceIterators[Walker->nr]->second; // put next farther distance into PermutationMap
336 PermutationMap[Runner->nr] = Sprinter; // and hand the old target to its respective owner
337 PrintPermutationMap(out, PermutationMap, AtomCount);
338 // calculate the new potential
339 //*out << Verbose(2) << "Checking new potential ..." << endl;
340 Potential = ConstrainedPotential(out, PermutationMap, startstep, endstep, constants, IsAngstroem);
341 if (Potential > OldPotential) { // we made everything worse! Undo ...
342 //*out << Verbose(3) << "Nay, made the potential worse: " << Potential << " vs. " << OldPotential << "!" << endl;
343 //*out << Verbose(3) << "Setting " << *Runner << "'s source to " << *DistanceIterators[Runner->nr]->second << "." << endl;
344 // Undo for Runner (note, we haven't moved the iteration yet, we may use this)
345 PermutationMap[Runner->nr] = DistanceIterators[Runner->nr]->second;
346 // Undo for Walker
347 DistanceIterators[Walker->nr] = Strider; // take next farther distance target
348 //*out << Verbose(3) << "Setting " << *Walker << "'s source to " << *DistanceIterators[Walker->nr]->second << "." << endl;
349 PermutationMap[Walker->nr] = DistanceIterators[Walker->nr]->second;
350 } else {
351 DistanceIterators[Runner->nr] = Rider; // if successful also move the pointer in the iterator list
352 *out << Verbose(3) << "Found a better permutation, new potential is " << Potential << " vs." << OldPotential << "." << endl;
353 OldPotential = Potential;
354 }
355 if (Potential > constants[2]) {
356 cerr << "ERROR: The two-step permutation procedure did not maintain injectivity!" << endl;
357 exit(255);
358 }
359 //*out << endl;
360 } else {
361 cerr << "ERROR: " << *Runner << " was not the owner of " << *Sprinter << "!" << endl;
362 exit(255);
363 }
364 } else {
365 PermutationMap[Walker->nr] = DistanceIterators[Walker->nr]->second; // new target has no source!
366 }
367 StepList[Walker->nr]++; // take next farther distance target
368 }
369 } while (Walker->next != end);
370 } while ((OlderPotential - OldPotential) > 1e-3);
371 *out << Verbose(1) << "done." << endl;
372
373
374 /// free memory and return with evaluated potential
375 for (int i=AtomCount; i--;)
376 DistanceList[i]->clear();
377 Free(&DistanceList);
378 Free(&DistanceIterators);
379 return ConstrainedPotential(out, PermutationMap, startstep, endstep, constants, IsAngstroem);
380};
381
382/** Evaluates the (distance-related part) of the constrained potential for the constrained forces.
383 * \param *out output stream for debugging
384 * \param startstep current MD step giving initial position between which and \a endstep we perform the constrained MD (as further steps are always concatenated)
385 * \param endstep step giving final position in constrained MD
386 * \param **PermutationMap mapping between the atom label of the initial and the final configuration
387 * \param *Force ForceMatrix containing force vectors from the external energy functional minimisation.
388 * \todo the constant for the constrained potential distance part is hard-coded independently of the hard-coded value in MinimiseConstrainedPotential()
389 */
390void molecule::EvaluateConstrainedForces(ofstream *out, int startstep, int endstep, atom **PermutationMap, ForceMatrix *Force)
391{
392 double constant = 10.;
393 atom *Walker = NULL, *Sprinter = NULL;
394
395 /// evaluate forces (only the distance to target dependent part) with the final PermutationMap
396 *out << Verbose(1) << "Calculating forces and adding onto ForceMatrix ... " << endl;
397 Walker = start;
398 while (Walker->next != NULL) {
399 Walker = Walker->next;
400 Sprinter = PermutationMap[Walker->nr];
401 // set forces
402 for (int i=NDIM;i++;)
403 Force->Matrix[0][Walker->nr][5+i] += 2.*constant*sqrt(Walker->Trajectory.R.at(startstep).Distance(&Sprinter->Trajectory.R.at(endstep)));
404 }
405 *out << Verbose(1) << "done." << endl;
406};
407
408/** Performs a linear interpolation between two desired atomic configurations with a given number of steps.
409 * Note, step number is config::MaxOuterStep
410 * \param *out output stream for debugging
411 * \param startstep stating initial configuration in molecule::Trajectories
412 * \param endstep stating final configuration in molecule::Trajectories
413 * \param &config configuration structure
414 * \param MapByIdentity if true we just use the identity to map atoms in start config to end config, if not we find mapping by \sa MinimiseConstrainedPotential()
415 * \return true - success in writing step files, false - error writing files or only one step in molecule::Trajectories
416 */
417bool molecule::LinearInterpolationBetweenConfiguration(ofstream *out, int startstep, int endstep, const char *prefix, config &configuration, bool MapByIdentity)
418{
419 molecule *mol = NULL;
420 bool status = true;
421 int MaxSteps = configuration.MaxOuterStep;
422 MoleculeListClass *MoleculePerStep = new MoleculeListClass();
423 // Get the Permutation Map by MinimiseConstrainedPotential
424 atom **PermutationMap = NULL;
425 atom *Walker = NULL, *Sprinter = NULL;
426 if (!MapByIdentity)
427 MinimiseConstrainedPotential(out, PermutationMap, startstep, endstep, configuration.GetIsAngstroem());
428 else {
429 PermutationMap = Malloc<atom *>(AtomCount, "molecule::LinearInterpolationBetweenConfiguration: **PermutationMap");
430 Walker = start;
431 while (Walker->next != end) {
432 Walker = Walker->next;
433 PermutationMap[Walker->nr] = Walker; // always pick target with the smallest distance
434 }
435 }
436
437 // check whether we have sufficient space in Trajectories for each atom
438 Walker = start;
439 while (Walker->next != end) {
440 Walker = Walker->next;
441 if (Walker->Trajectory.R.size() <= (unsigned int)(MaxSteps)) {
442 //cout << "Increasing size for trajectory array of " << keyword << " to " << (MaxSteps+1) << "." << endl;
443 Walker->Trajectory.R.resize(MaxSteps+1);
444 Walker->Trajectory.U.resize(MaxSteps+1);
445 Walker->Trajectory.F.resize(MaxSteps+1);
446 }
447 }
448 // push endstep to last one
449 Walker = start;
450 while (Walker->next != end) { // remove the endstep (is now the very last one)
451 Walker = Walker->next;
452 for (int n=NDIM;n--;) {
453 Walker->Trajectory.R.at(MaxSteps).x[n] = Walker->Trajectory.R.at(endstep).x[n];
454 Walker->Trajectory.U.at(MaxSteps).x[n] = Walker->Trajectory.U.at(endstep).x[n];
455 Walker->Trajectory.F.at(MaxSteps).x[n] = Walker->Trajectory.F.at(endstep).x[n];
456 }
457 }
458 endstep = MaxSteps;
459
460 // go through all steps and add the molecular configuration to the list and to the Trajectories of \a this molecule
461 *out << Verbose(1) << "Filling intermediate " << MaxSteps << " steps with MDSteps of " << MDSteps << "." << endl;
462 for (int step = 0; step <= MaxSteps; step++) {
463 mol = new molecule(elemente);
464 MoleculePerStep->insert(mol);
465 Walker = start;
466 while (Walker->next != end) {
467 Walker = Walker->next;
468 // add to molecule list
469 Sprinter = mol->AddCopyAtom(Walker);
470 for (int n=NDIM;n--;) {
471 Sprinter->x.x[n] = Walker->Trajectory.R.at(startstep).x[n] + (PermutationMap[Walker->nr]->Trajectory.R.at(endstep).x[n] - Walker->Trajectory.R.at(startstep).x[n])*((double)step/(double)MaxSteps);
472 // add to Trajectories
473 //*out << Verbose(3) << step << ">=" << MDSteps-1 << endl;
474 if (step < MaxSteps) {
475 Walker->Trajectory.R.at(step).x[n] = Walker->Trajectory.R.at(startstep).x[n] + (PermutationMap[Walker->nr]->Trajectory.R.at(endstep).x[n] - Walker->Trajectory.R.at(startstep).x[n])*((double)step/(double)MaxSteps);
476 Walker->Trajectory.U.at(step).x[n] = 0.;
477 Walker->Trajectory.F.at(step).x[n] = 0.;
478 }
479 }
480 }
481 }
482 MDSteps = MaxSteps+1; // otherwise new Trajectories' points aren't stored on save&exit
483
484 // store the list to single step files
485 int *SortIndex = Malloc<int>(AtomCount, "molecule::LinearInterpolationBetweenConfiguration: *SortIndex");
486 for (int i=AtomCount; i--; )
487 SortIndex[i] = i;
488 status = MoleculePerStep->OutputConfigForListOfFragments(out, &configuration, SortIndex);
489
490 // free and return
491 Free(&PermutationMap);
492 delete(MoleculePerStep);
493 return status;
494};
495
496/** Parses nuclear forces from file and performs Verlet integration.
497 * Note that we assume the parsed forces to be in atomic units (hence, if coordinates are in angstroem, we
498 * have to transform them).
499 * This adds a new MD step to the config file.
500 * \param *out output stream for debugging
501 * \param *file filename
502 * \param config structure with config::Deltat, config::IsAngstroem, config::DoConstrained
503 * \param delta_t time step width in atomic units
504 * \param IsAngstroem whether coordinates are in angstroem (true) or bohrradius (false)
505 * \param DoConstrained whether we perform a constrained (>0, target step in molecule::trajectories) or unconstrained (0) molecular dynamics, \sa molecule::MinimiseConstrainedPotential()
506 * \return true - file found and parsed, false - file not found or imparsable
507 * \todo This is not yet checked if it is correctly working with DoConstrained set to true.
508 */
509bool molecule::VerletForceIntegration(ofstream *out, char *file, config &configuration)
510{
511 atom *walker = NULL;
512 ifstream input(file);
513 string token;
514 stringstream item;
515 double IonMass, Vector[NDIM], ConstrainedPotentialEnergy, ActualTemp;
516 ForceMatrix Force;
517
518 CountElements(); // make sure ElementsInMolecule is up to date
519
520 // check file
521 if (input == NULL) {
522 return false;
523 } else {
524 // parse file into ForceMatrix
525 if (!Force.ParseMatrix(file, 0,0,0)) {
526 cerr << "Could not parse Force Matrix file " << file << "." << endl;
527 return false;
528 }
529 if (Force.RowCounter[0] != AtomCount) {
530 cerr << "Mismatch between number of atoms in file " << Force.RowCounter[0] << " and in molecule " << AtomCount << "." << endl;
531 return false;
532 }
533 // correct Forces
534 for(int d=0;d<NDIM;d++)
535 Vector[d] = 0.;
536 for(int i=0;i<AtomCount;i++)
537 for(int d=0;d<NDIM;d++) {
538 Vector[d] += Force.Matrix[0][i][d+5];
539 }
540 for(int i=0;i<AtomCount;i++)
541 for(int d=0;d<NDIM;d++) {
542 Force.Matrix[0][i][d+5] -= Vector[d]/(double)AtomCount;
543 }
544 // solve a constrained potential if we are meant to
545 if (configuration.DoConstrainedMD) {
546 // calculate forces and potential
547 atom **PermutationMap = NULL;
548 ConstrainedPotentialEnergy = MinimiseConstrainedPotential(out, PermutationMap,configuration.DoConstrainedMD, 0, configuration.GetIsAngstroem());
549 EvaluateConstrainedForces(out, configuration.DoConstrainedMD, 0, PermutationMap, &Force);
550 Free(&PermutationMap);
551 }
552
553 // and perform Verlet integration for each atom with position, velocity and force vector
554 walker = start;
555 while (walker->next != end) { // go through every atom of this element
556 walker = walker->next;
557 //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
558 // check size of vectors
559 if (walker->Trajectory.R.size() <= (unsigned int)(MDSteps)) {
560 //out << "Increasing size for trajectory array of " << *walker << " to " << (size+10) << "." << endl;
561 walker->Trajectory.R.resize(MDSteps+10);
562 walker->Trajectory.U.resize(MDSteps+10);
563 walker->Trajectory.F.resize(MDSteps+10);
564 }
565
566 // Update R (and F)
567 for (int d=0; d<NDIM; d++) {
568 walker->Trajectory.F.at(MDSteps).x[d] = -Force.Matrix[0][walker->nr][d+5]*(configuration.GetIsAngstroem() ? AtomicLengthToAngstroem : 1.);
569 walker->Trajectory.R.at(MDSteps).x[d] = walker->Trajectory.R.at(MDSteps-1).x[d];
570 walker->Trajectory.R.at(MDSteps).x[d] += configuration.Deltat*(walker->Trajectory.U.at(MDSteps-1).x[d]); // s(t) = s(0) + v * deltat + 1/2 a * deltat^2
571 walker->Trajectory.R.at(MDSteps).x[d] += 0.5*configuration.Deltat*configuration.Deltat*(walker->Trajectory.F.at(MDSteps).x[d]/walker->type->mass); // F = m * a and s = 0.5 * F/m * t^2 = F * a * t
572 }
573 // Update U
574 for (int d=0; d<NDIM; d++) {
575 walker->Trajectory.U.at(MDSteps).x[d] = walker->Trajectory.U.at(MDSteps-1).x[d];
576 walker->Trajectory.U.at(MDSteps).x[d] += configuration.Deltat * (walker->Trajectory.F.at(MDSteps).x[d]+walker->Trajectory.F.at(MDSteps-1).x[d]/walker->type->mass); // v = F/m * t
577 }
578// out << "Integrated position&velocity of step " << (MDSteps) << ": (";
579// for (int d=0;d<NDIM;d++)
580// out << walker->Trajectory.R.at(MDSteps).x[d] << " "; // next step
581// out << ")\t(";
582// for (int d=0;d<NDIM;d++)
583// cout << walker->Trajectory.U.at(MDSteps).x[d] << " "; // next step
584// out << ")" << endl;
585 // next atom
586 }
587 }
588 // correct velocities (rather momenta) so that center of mass remains motionless
589 for(int d=0;d<NDIM;d++)
590 Vector[d] = 0.;
591 IonMass = 0.;
592 walker = start;
593 while (walker->next != end) { // go through every atom
594 walker = walker->next;
595 IonMass += walker->type->mass; // sum up total mass
596 for(int d=0;d<NDIM;d++) {
597 Vector[d] += walker->Trajectory.U.at(MDSteps).x[d]*walker->type->mass;
598 }
599 }
600 // correct velocities (rather momenta) so that center of mass remains motionless
601 for(int d=0;d<NDIM;d++)
602 Vector[d] /= IonMass;
603 ActualTemp = 0.;
604 walker = start;
605 while (walker->next != end) { // go through every atom of this element
606 walker = walker->next;
607 for(int d=0;d<NDIM;d++) {
608 walker->Trajectory.U.at(MDSteps).x[d] -= Vector[d];
609 ActualTemp += 0.5 * walker->type->mass * walker->Trajectory.U.at(MDSteps).x[d] * walker->Trajectory.U.at(MDSteps).x[d];
610 }
611 }
612 Thermostats(configuration, ActualTemp, Berendsen);
613 MDSteps++;
614
615
616 // exit
617 return true;
618};
619
620/** Implementation of various thermostats.
621 * All these thermostats apply an additional force which has the following forms:
622 * -# Woodcock
623 * \f$p_i \rightarrow \sqrt{\frac{T_0}{T}} \cdot p_i\f$
624 * -# Gaussian
625 * \f$ \frac{ \sum_i \frac{p_i}{m_i} \frac{\partial V}{\partial q_i}} {\sum_i \frac{p^2_i}{m_i}} \cdot p_i\f$
626 * -# Langevin
627 * \f$p_{i,n} \rightarrow \sqrt{1-\alpha^2} p_{i,0} + \alpha p_r\f$
628 * -# Berendsen
629 * \f$p_i \rightarrow \left [ 1+ \frac{\delta t}{\tau_T} \left ( \frac{T_0}{T} \right ) \right ]^{\frac{1}{2}} \cdot p_i\f$
630 * -# Nose-Hoover
631 * \f$\zeta p_i \f$ with \f$\frac{\partial \zeta}{\partial t} = \frac{1}{M_s} \left ( \sum^N_{i=1} \frac{p_i^2}{m_i} - g k_B T \right )\f$
632 * These Thermostats either simply rescale the velocities, thus this function should be called after ion velocities have been updated, and/or
633 * have a constraint force acting additionally on the ions. In the latter case, the ion speeds have to be modified
634 * belatedly and the constraint force set.
635 * \param *P Problem at hand
636 * \param i which of the thermostats to take: 0 - none, 1 - Woodcock, 2 - Gaussian, 3 - Langevin, 4 - Berendsen, 5 - Nose-Hoover
637 * \sa InitThermostat()
638 */
639void molecule::Thermostats(config &configuration, double ActualTemp, int Thermostat)
640{
641 double ekin = 0.;
642 double E = 0., G = 0.;
643 double delta_alpha = 0.;
644 double ScaleTempFactor;
645 double sigma;
646 double IonMass;
647 int d;
648 gsl_rng * r;
649 const gsl_rng_type * T;
650 double *U = NULL, *F = NULL, FConstraint[NDIM];
651 atom *walker = NULL;
652
653 // calculate scale configuration
654 ScaleTempFactor = configuration.TargetTemp/ActualTemp;
655
656 // differentating between the various thermostats
657 switch(Thermostat) {
658 case None:
659 cout << Verbose(2) << "Applying no thermostat..." << endl;
660 break;
661 case Woodcock:
662 if ((configuration.ScaleTempStep > 0) && ((MDSteps-1) % configuration.ScaleTempStep == 0)) {
663 cout << Verbose(2) << "Applying Woodcock thermostat..." << endl;
664 walker = start;
665 while (walker->next != end) { // go through every atom of this element
666 walker = walker->next;
667 IonMass = walker->type->mass;
668 U = walker->Trajectory.U.at(MDSteps).x;
669 if (walker->FixedIon == 0) // even FixedIon moves, only not by other's forces
670 for (d=0; d<NDIM; d++) {
671 U[d] *= sqrt(ScaleTempFactor);
672 ekin += 0.5*IonMass * U[d]*U[d];
673 }
674 }
675 }
676 break;
677 case Gaussian:
678 cout << Verbose(2) << "Applying Gaussian thermostat..." << endl;
679 walker = start;
680 while (walker->next != end) { // go through every atom of this element
681 walker = walker->next;
682 IonMass = walker->type->mass;
683 U = walker->Trajectory.U.at(MDSteps).x;
684 F = walker->Trajectory.F.at(MDSteps).x;
685 if (walker->FixedIon == 0) // even FixedIon moves, only not by other's forces
686 for (d=0; d<NDIM; d++) {
687 G += U[d] * F[d];
688 E += U[d]*U[d]*IonMass;
689 }
690 }
691 cout << Verbose(1) << "Gaussian Least Constraint constant is " << G/E << "." << endl;
692 walker = start;
693 while (walker->next != end) { // go through every atom of this element
694 walker = walker->next;
695 IonMass = walker->type->mass;
696 U = walker->Trajectory.U.at(MDSteps).x;
697 F = walker->Trajectory.F.at(MDSteps).x;
698 if (walker->FixedIon == 0) // even FixedIon moves, only not by other's forces
699 for (d=0; d<NDIM; d++) {
700 FConstraint[d] = (G/E) * (U[d]*IonMass);
701 U[d] += configuration.Deltat/IonMass * (FConstraint[d]);
702 ekin += IonMass * U[d]*U[d];
703 }
704 }
705 break;
706 case Langevin:
707 cout << Verbose(2) << "Applying Langevin thermostat..." << endl;
708 // init random number generator
709 gsl_rng_env_setup();
710 T = gsl_rng_default;
711 r = gsl_rng_alloc (T);
712 // Go through each ion
713 walker = start;
714 while (walker->next != end) { // go through every atom of this element
715 walker = walker->next;
716 IonMass = walker->type->mass;
717 sigma = sqrt(configuration.TargetTemp/IonMass); // sigma = (k_b T)/m (Hartree/atomicmass = atomiclength/atomictime)
718 U = walker->Trajectory.U.at(MDSteps).x;
719 F = walker->Trajectory.F.at(MDSteps).x;
720 if (walker->FixedIon == 0) { // even FixedIon moves, only not by other's forces
721 // throw a dice to determine whether it gets hit by a heat bath particle
722 if (((((rand()/(double)RAND_MAX))*configuration.TempFrequency) < 1.)) {
723 cout << Verbose(3) << "Particle " << *walker << " was hit (sigma " << sigma << "): " << sqrt(U[0]*U[0]+U[1]*U[1]+U[2]*U[2]) << " -> ";
724 // pick three random numbers from a Boltzmann distribution around the desired temperature T for each momenta axis
725 for (d=0; d<NDIM; d++) {
726 U[d] = gsl_ran_gaussian (r, sigma);
727 }
728 cout << sqrt(U[0]*U[0]+U[1]*U[1]+U[2]*U[2]) << endl;
729 }
730 for (d=0; d<NDIM; d++)
731 ekin += 0.5*IonMass * U[d]*U[d];
732 }
733 }
734 break;
735 case Berendsen:
736 cout << Verbose(2) << "Applying Berendsen-VanGunsteren thermostat..." << endl;
737 walker = start;
738 while (walker->next != end) { // go through every atom of this element
739 walker = walker->next;
740 IonMass = walker->type->mass;
741 U = walker->Trajectory.U.at(MDSteps).x;
742 F = walker->Trajectory.F.at(MDSteps).x;
743 if (walker->FixedIon == 0) { // even FixedIon moves, only not by other's forces
744 for (d=0; d<NDIM; d++) {
745 U[d] *= sqrt(1+(configuration.Deltat/configuration.TempFrequency)*(ScaleTempFactor-1));
746 ekin += 0.5*IonMass * U[d]*U[d];
747 }
748 }
749 }
750 break;
751 case NoseHoover:
752 cout << Verbose(2) << "Applying Nose-Hoover thermostat..." << endl;
753 // dynamically evolve alpha (the additional degree of freedom)
754 delta_alpha = 0.;
755 walker = start;
756 while (walker->next != end) { // go through every atom of this element
757 walker = walker->next;
758 IonMass = walker->type->mass;
759 U = walker->Trajectory.U.at(MDSteps).x;
760 if (walker->FixedIon == 0) { // even FixedIon moves, only not by other's forces
761 for (d=0; d<NDIM; d++) {
762 delta_alpha += U[d]*U[d]*IonMass;
763 }
764 }
765 }
766 delta_alpha = (delta_alpha - (3.*AtomCount+1.) * configuration.TargetTemp)/(configuration.HooverMass*Units2Electronmass);
767 configuration.alpha += delta_alpha*configuration.Deltat;
768 cout << Verbose(3) << "alpha = " << delta_alpha << " * " << configuration.Deltat << " = " << configuration.alpha << "." << endl;
769 // apply updated alpha as additional force
770 walker = start;
771 while (walker->next != end) { // go through every atom of this element
772 walker = walker->next;
773 IonMass = walker->type->mass;
774 U = walker->Trajectory.U.at(MDSteps).x;
775 if (walker->FixedIon == 0) { // even FixedIon moves, only not by other's forces
776 for (d=0; d<NDIM; d++) {
777 FConstraint[d] = - configuration.alpha * (U[d] * IonMass);
778 U[d] += configuration.Deltat/IonMass * (FConstraint[d]);
779 ekin += (0.5*IonMass) * U[d]*U[d];
780 }
781 }
782 }
783 break;
784 }
785 cout << Verbose(1) << "Kinetic energy is " << ekin << "." << endl;
786};
Note: See TracBrowser for help on using the repository browser.