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