source: src/molecule_dynamics.cpp@ 3839e5

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 3839e5 was aafd77, checked in by Tillmann Crueger <crueger@…>, 14 years ago

Removed all inclusions of GSL-Headers from molecule.hpp

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