source: src/molecule_dynamics.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 ebe4d6, checked in by Frederik Heber <heber@…>, 15 years ago

BUGFIX - molecule::MinimiseConstrainedPotential() OldPotential was used outside of loop in output message.

Signed-off-by: Frederik Heber <heber@tabletINS.(none)>

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