source: src/molecule_dynamics.cpp@ f2bb0f

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 f2bb0f was 9879f6, checked in by Frederik Heber <heber@…>, 15 years ago

Huge Refactoring due to class molecule now being an STL container.

  • molecule::start and molecule::end were dropped. Hence, the usual construct Walker = start while (Walker->next != end) {

Walker = walker->next
...

}
was changed to
for (molecule::iterator iter = begin(); iter != end(); ++iter) {

...

}
and (*iter) used instead of Walker.

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