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