source: src/Dynamics/ForceAnnealing.hpp@ b03d7d

ForceAnnealing_tocheck
Last change on this file since b03d7d was b03d7d, checked in by Frederik Heber <frederik.heber@…>, 8 years ago

tempcommit: Merge with 9190e6a

  • Property mode set to 100644
File size: 11.0 KB
Line 
1/*
2 * ForceAnnealing.hpp
3 *
4 * Created on: Aug 02, 2014
5 * Author: heber
6 */
7
8#ifndef FORCEANNEALING_HPP_
9#define FORCEANNEALING_HPP_
10
11// include config.h
12#ifdef HAVE_CONFIG_H
13#include <config.h>
14#endif
15
16#include "Atom/atom.hpp"
17#include "Atom/AtomSet.hpp"
18#include "CodePatterns/Assert.hpp"
19#include "CodePatterns/Info.hpp"
20#include "CodePatterns/Log.hpp"
21#include "CodePatterns/Verbose.hpp"
22#include "Descriptors/AtomIdDescriptor.hpp"
23#include "Dynamics/AtomicForceManipulator.hpp"
24#include "Fragmentation/ForceMatrix.hpp"
25#include "Graph/BoostGraphCreator.hpp"
26#include "Graph/BoostGraphHelpers.hpp"
27#include "Graph/BreadthFirstSearchGatherer.hpp"
28#include "Helpers/helpers.hpp"
29#include "Helpers/defs.hpp"
30#include "LinearAlgebra/Vector.hpp"
31#include "Thermostats/ThermoStatContainer.hpp"
32#include "Thermostats/Thermostat.hpp"
33#include "World.hpp"
34
35/** This class is the essential build block for performing structural optimization.
36 *
37 * Sadly, we have to use some static instances as so far values cannot be passed
38 * between actions. Hence, we need to store the current step and the adaptive-
39 * step width (we cannot perform a line search, as we have no control over the
40 * calculation of the forces).
41 *
42 * However, we do use the bond graph, i.e. if a single atom needs to be shifted
43 * to the left, then the whole molecule left of it is shifted, too. This is
44 * controlled by the \a max_distance parameter.
45 */
46template <class T>
47class ForceAnnealing : public AtomicForceManipulator<T>
48{
49public:
50 /** Constructor of class ForceAnnealing.
51 *
52 * \note We use a fixed delta t of 1.
53 *
54 * \param _atoms set of atoms to integrate
55 * \param _Deltat time step width in atomic units
56 * \param _IsAngstroem whether length units are in angstroem or bohr radii
57 * \param _maxSteps number of optimization steps to perform
58 * \param _max_distance up to this bond order is bond graph taken into account.
59 */
60 ForceAnnealing(
61 AtomSetMixin<T> &_atoms,
62 bool _IsAngstroem,
63 const size_t _maxSteps,
64 const int _max_distance) :
65 AtomicForceManipulator<T>(_atoms, 1., _IsAngstroem),
66 maxSteps(_maxSteps),
67 max_distance(_max_distance),
68 damping_factor(0.5)
69 {}
70 /** Destructor of class ForceAnnealing.
71 *
72 */
73 ~ForceAnnealing()
74 {}
75
76 /** Performs Gradient optimization.
77 *
78 * We assume that forces have just been calculated.
79 *
80 *
81 * \param NextStep current time step (i.e. \f$ t + \Delta t \f$ in the sense of the velocity verlet)
82 * \param offset offset in matrix file to the first force component
83 * \todo This is not yet checked if it is correctly working with DoConstrainedMD set >0.
84 */
85 void operator()(const int NextStep, const size_t offset)
86 {
87 // make sum of forces equal zero
88 AtomicForceManipulator<T>::correctForceMatrixForFixedCenterOfMass(offset,NextStep);
89
90 // are we in initial step? Then set static entities
91 if (currentStep == 0) {
92 currentDeltat = AtomicForceManipulator<T>::Deltat;
93 currentStep = 1;
94 LOG(2, "DEBUG: Initial step, setting values, current step is #" << currentStep);
95 } else {
96 ++currentStep;
97 LOG(2, "DEBUG: current step is #" << currentStep);
98 }
99
100 // get nodes on either side of selected bond via BFS discovery
101// std::vector<atomId_t> atomids;
102// for(typename AtomSetMixin<T>::iterator iter = AtomicForceManipulator<T>::atoms.begin();
103// iter != AtomicForceManipulator<T>::atoms.end(); ++iter) {
104// atomids.push_back((*iter)->getId());
105// }
106// ASSERT( atomids.size() == AtomicForceManipulator<T>::atoms.size(),
107// "ForceAnnealing() - could not gather all atomic ids?");
108 BoostGraphCreator BGcreator;
109 BGcreator.createFromRange(
110 AtomicForceManipulator<T>::atoms.begin(),
111 AtomicForceManipulator<T>::atoms.end(),
112 AtomicForceManipulator<T>::atoms.size(),
113 BreadthFirstSearchGatherer::AlwaysTruePredicate);
114 BreadthFirstSearchGatherer NodeGatherer(BGcreator);
115
116 std::map<atomId_t, Vector> GatheredUpdates; //!< gathers all updates which are applied at the end
117 Vector maxComponents(zeroVec);
118 for(typename AtomSetMixin<T>::iterator iter = AtomicForceManipulator<T>::atoms.begin();
119 iter != AtomicForceManipulator<T>::atoms.end(); ++iter) {
120 // atom's force vector gives steepest descent direction
121 const Vector oldPosition = (*iter)->getPositionAtStep(NextStep-2 >= 0 ? NextStep - 2 : 0);
122 const Vector currentPosition = (*iter)->getPosition();
123 const Vector oldGradient = (*iter)->getAtomicForceAtStep(NextStep-2 >= 0 ? NextStep - 2 : 0);
124 const Vector currentGradient = (*iter)->getAtomicForce();
125 LOG(4, "DEBUG: Force for atom " << **iter << " is " << currentGradient);
126
127 // we use Barzilai-Borwein update with position reversed to get descent
128 const Vector GradientDifference = (currentGradient - oldGradient);
129 const double stepwidth =
130 fabs((currentPosition - oldPosition).ScalarProduct(GradientDifference))/
131 GradientDifference.NormSquared();
132 Vector PositionUpdate = stepwidth * currentGradient;
133 if (fabs(stepwidth) < 1e-10) {
134 // dont' warn in first step, deltat usage normal
135 if (currentStep != 1)
136 ELOG(1, "INFO: Barzilai-Borwein stepwidth is zero, using deltat " << currentDeltat << " instead.");
137 PositionUpdate = currentDeltat * currentGradient;
138 }
139 LOG(3, "DEBUG: Update would be " << PositionUpdate);
140
141// // add update to central atom
142// const atomId_t atomid = (*iter)->getId();
143// if (GatheredUpdates.count(atomid)) {
144// GatheredUpdates[atomid] += PositionUpdate;
145// } else
146// GatheredUpdates.insert( std::make_pair(atomid, PositionUpdate) );
147
148 // We assume that a force is local, i.e. a bond is too short yet and hence
149 // the atom needs to be moved. However, all the adjacent (bound) atoms might
150 // already be at the perfect distance. If we just move the atom alone, we ruin
151 // all the other bonds. Hence, it would be sensible to move every atom found
152 // through the bond graph in the direction of the force as well by the same
153 // PositionUpdate. This is just what we are going to do.
154
155 /// get all nodes from bonds in the direction of the current force
156
157 // remove edges facing in the wrong direction
158 std::vector<bond::ptr> removed_bonds;
159 const BondList& ListOfBonds = (*iter)->getListOfBonds();
160 for(BondList::const_iterator bonditer = ListOfBonds.begin();
161 bonditer != ListOfBonds.end(); ++bonditer) {
162 const bond &current_bond = *(*bonditer);
163 LOG(2, "DEBUG: Looking at bond " << current_bond);
164 Vector BondVector = (*iter)->getPosition();
165 BondVector -= ((*iter)->getId() == current_bond.rightatom->getId())
166 ? current_bond.rightatom->getPosition() : current_bond.leftatom->getPosition();
167 BondVector.Normalize();
168 if (BondVector.ScalarProduct(currentGradient) < 0) {
169 removed_bonds.push_back(*bonditer);
170#ifndef NDEBUG
171 const bool status =
172#endif
173 BGcreator.removeEdge(current_bond.leftatom->getId(), current_bond.rightatom->getId());
174 ASSERT( status, "ForceAnnealing() - edge to found bond is not present?");
175 }
176 }
177 BoostGraphHelpers::Nodeset_t bondside_set = NodeGatherer((*iter)->getId(), max_distance);
178 const BreadthFirstSearchGatherer::distance_map_t &distance_map = NodeGatherer.getDistances();
179 std::sort(bondside_set.begin(), bondside_set.end());
180 // re-add those edges
181 for (std::vector<bond::ptr>::const_iterator bonditer = removed_bonds.begin();
182 bonditer != removed_bonds.end(); ++bonditer)
183 BGcreator.addEdge((*bonditer)->leftatom->getId(), (*bonditer)->rightatom->getId());
184
185 // apply PositionUpdate to all nodes in the bondside_set
186 for (BoostGraphHelpers::Nodeset_t::const_iterator setiter = bondside_set.begin();
187 setiter != bondside_set.end(); ++setiter) {
188 const BreadthFirstSearchGatherer::distance_map_t::const_iterator diter
189 = distance_map.find(*setiter);
190 ASSERT( diter != distance_map.end(),
191 "ForceAnnealing() - could not find distance to an atom.");
192 const double factor = pow(damping_factor, diter->second);
193 LOG(3, "DEBUG: Update for atom #" << *setiter << " will be "
194 << factor << "*" << PositionUpdate);
195 if (GatheredUpdates.count((*setiter))) {
196 GatheredUpdates[(*setiter)] += factor*PositionUpdate;
197 } else {
198 GatheredUpdates.insert(
199 std::make_pair(
200 (*setiter),
201 factor*PositionUpdate) );
202 }
203 }
204
205 // extract largest components for showing progress of annealing
206 for(size_t i=0;i<NDIM;++i)
207 if (currentGradient[i] > maxComponents[i])
208 maxComponents[i] = currentGradient[i];
209
210 // reset force vector for next step except on final one
211 if (currentStep != maxSteps)
212 (*iter)->setAtomicForce(zeroVec);
213 }
214 // apply the gathered updates
215 for (std::map<atomId_t, Vector>::const_iterator iter = GatheredUpdates.begin();
216 iter != GatheredUpdates.end(); ++iter) {
217 const atomId_t &atomid = iter->first;
218 const Vector &update = iter->second;
219 atom* const walker = World::getInstance().getAtom(AtomById(atomid));
220 ASSERT( walker != NULL,
221 "ForceAnnealing() - walker with id "+toString(atomid)+" has suddenly disappeared.");
222 walker->setPosition( walker->getPosition() + update );
223 walker->setAtomicVelocity(update);
224 LOG(3, "DEBUG: Applying update " << update << " to atom " << *walker);
225 }
226
227 LOG(1, "STATUS: Largest remaining force components at step #"
228 << currentStep << " are " << maxComponents);
229
230 // are we in final step? Remember to reset static entities
231 if (currentStep == maxSteps) {
232 LOG(2, "DEBUG: Final step, resetting values");
233 reset();
234 }
235 }
236
237 /** Reset function to unset static entities and artificial velocities.
238 *
239 */
240 void reset()
241 {
242 currentDeltat = 0.;
243 currentStep = 0;
244
245 // reset (artifical) velocities
246 for(typename AtomSetMixin<T>::iterator iter = AtomicForceManipulator<T>::atoms.begin();
247 iter != AtomicForceManipulator<T>::atoms.end(); ++iter)
248 (*iter)->setAtomicVelocity(zeroVec);
249 }
250
251private:
252 //!> contains the current step in relation to maxsteps
253 static size_t currentStep;
254 //!> contains the maximum number of steps, determines initial and final step with currentStep
255 size_t maxSteps;
256 static double currentDeltat;
257 //!> minimum deltat for internal while loop (adaptive step width)
258 static double MinimumDeltat;
259 //!> contains the maximum bond graph distance up to which shifts of a single atom are spread
260 const int max_distance;
261 //!> the shifted is dampened by this factor with the power of the bond graph distance to the shift causing atom
262 const double damping_factor;
263};
264
265template <class T>
266double ForceAnnealing<T>::currentDeltat = 0.;
267template <class T>
268size_t ForceAnnealing<T>::currentStep = 0;
269template <class T>
270double ForceAnnealing<T>::MinimumDeltat = 1e-8;
271
272#endif /* FORCEANNEALING_HPP_ */
Note: See TracBrowser for help on using the repository browser.