source: src/Dynamics/ForceAnnealing.hpp@ 866dec

AutomationFragmentation_failures Candidate_v1.6.1 ChemicalSpaceEvaluator Enhanced_StructuralOptimization_continued Exclude_Hydrogens_annealWithBondGraph ForceAnnealing_with_BondGraph ForceAnnealing_with_BondGraph_contraction-expansion Gui_displays_atomic_force_velocity PythonUI_with_named_parameters StoppableMakroAction TremoloParser_IncreasedPrecision
Last change on this file since 866dec was 866dec, checked in by Frederik Heber <frederik.heber@…>, 7 years ago

ForceAnnealing no longer resets forces, AnalyseFragmentResults now sets forces.

  • before we added onto force vector in AnalyseFragmentResults.
  • if we do not reset forces in ForceAnnealing we can inspect the force components in the output file along the annealing run.
  • TESTFIX: data files for regression test AnalyseFragmentResults OBC and PBC changed because of modified way of resetting forces.
  • Property mode set to 100644
File size: 15.3 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 const double _Deltat,
63 bool _IsAngstroem,
64 const size_t _maxSteps,
65 const int _max_distance,
66 const double _damping_factor) :
67 AtomicForceManipulator<T>(_atoms, _Deltat, _IsAngstroem),
68 maxSteps(_maxSteps),
69 max_distance(_max_distance),
70 damping_factor(_damping_factor)
71 {}
72
73 /** Destructor of class ForceAnnealing.
74 *
75 */
76 ~ForceAnnealing()
77 {}
78
79 /** Performs Gradient optimization.
80 *
81 * We assume that forces have just been calculated.
82 *
83 *
84 * \param CurrentTimeStep current time step (i.e. \f$ t + \Delta t \f$ in the sense of the velocity verlet)
85 * \param offset offset in matrix file to the first force component
86 * \todo This is not yet checked if it is correctly working with DoConstrainedMD set >0.
87 */
88 void operator()(
89 const int _CurrentTimeStep,
90 const size_t _offset,
91 const bool _UseBondgraph)
92 {
93 // make sum of forces equal zero
94 AtomicForceManipulator<T>::correctForceMatrixForFixedCenterOfMass(_offset, _CurrentTimeStep);
95
96 // are we in initial step? Then set static entities
97 Vector maxComponents(zeroVec);
98 if (currentStep == 0) {
99 currentDeltat = AtomicForceManipulator<T>::Deltat;
100 currentStep = 1;
101 LOG(2, "DEBUG: Initial step, setting values, current step is #" << currentStep);
102
103 // always use atomic annealing on first step
104 anneal(_CurrentTimeStep, _offset, maxComponents);
105 } else {
106 ++currentStep;
107 LOG(2, "DEBUG: current step is #" << currentStep);
108
109 if (_UseBondgraph)
110 annealWithBondGraph(_CurrentTimeStep, _offset, maxComponents);
111 else
112 anneal(_CurrentTimeStep, _offset, maxComponents);
113 }
114
115 LOG(1, "STATUS: Largest remaining force components at step #"
116 << currentStep << " are " << maxComponents);
117
118 // are we in final step? Remember to reset static entities
119 if (currentStep == maxSteps) {
120 LOG(2, "DEBUG: Final step, resetting values");
121 reset();
122 }
123 }
124
125 /** Performs Gradient optimization on the atoms.
126 *
127 * We assume that forces have just been calculated.
128 *
129 * \param CurrentTimeStep current time step (i.e. \f$ t + \Delta t \f$ in the sense of the velocity verlet)
130 * \param offset offset in matrix file to the first force component
131 * \param maxComponents to be filled with maximum force component over all atoms
132 */
133 void anneal(
134 const int CurrentTimeStep,
135 const size_t offset,
136 Vector &maxComponents)
137 {
138 for(typename AtomSetMixin<T>::iterator iter = AtomicForceManipulator<T>::atoms.begin();
139 iter != AtomicForceManipulator<T>::atoms.end(); ++iter) {
140 // atom's force vector gives steepest descent direction
141 const Vector oldPosition = (*iter)->getPositionAtStep(CurrentTimeStep-2 >= 0 ? CurrentTimeStep - 2 : 0);
142 const Vector currentPosition = (*iter)->getPosition();
143 const Vector oldGradient = (*iter)->getAtomicForceAtStep(CurrentTimeStep-2 >= 0 ? CurrentTimeStep - 2 : 0);
144 const Vector currentGradient = (*iter)->getAtomicForce();
145 LOG(4, "DEBUG: oldPosition for atom " << **iter << " is " << oldPosition);
146 LOG(4, "DEBUG: currentPosition for atom " << **iter << " is " << currentPosition);
147 LOG(4, "DEBUG: oldGradient for atom " << **iter << " is " << oldGradient);
148 LOG(4, "DEBUG: currentGradient for atom " << **iter << " is " << currentGradient);
149// LOG(4, "DEBUG: Force for atom " << **iter << " is " << currentGradient);
150
151 // we use Barzilai-Borwein update with position reversed to get descent
152 const Vector PositionDifference = currentPosition - oldPosition;
153 const Vector GradientDifference = (currentGradient - oldGradient);
154 double stepwidth = 0.;
155 if (GradientDifference.Norm() > MYEPSILON)
156 stepwidth = fabs(PositionDifference.ScalarProduct(GradientDifference))/
157 GradientDifference.NormSquared();
158 if (fabs(stepwidth) < 1e-10) {
159 // dont' warn in first step, deltat usage normal
160 if (currentStep != 1)
161 ELOG(1, "INFO: Barzilai-Borwein stepwidth is zero, using deltat " << currentDeltat << " instead.");
162 stepwidth = currentDeltat;
163 }
164 Vector PositionUpdate = stepwidth * currentGradient;
165 LOG(3, "DEBUG: Update would be " << stepwidth << "*" << currentGradient << " = " << PositionUpdate);
166
167 // extract largest components for showing progress of annealing
168 for(size_t i=0;i<NDIM;++i)
169 if (currentGradient[i] > maxComponents[i])
170 maxComponents[i] = currentGradient[i];
171
172 // are we in initial step? Then don't check against velocity
173 if ((currentStep > 1) && (!(*iter)->getAtomicVelocity().IsZero()))
174 // update with currentDelta tells us how the current gradient relates to
175 // the last one: If it has become larger, reduce currentDelta
176 if ((PositionUpdate.ScalarProduct((*iter)->getAtomicVelocity()) < 0)
177 && (currentDeltat > MinimumDeltat)) {
178 currentDeltat = .5*currentDeltat;
179 LOG(2, "DEBUG: Upgrade in other direction: " << PositionUpdate.NormSquared()
180 << " > " << (*iter)->getAtomicVelocity().NormSquared()
181 << ", decreasing deltat: " << currentDeltat);
182 PositionUpdate = currentDeltat * currentGradient;
183 }
184 // finally set new values
185 (*iter)->setPosition(currentPosition + PositionUpdate);
186 (*iter)->setAtomicVelocity(PositionUpdate);
187 //std::cout << "Id of atom is " << (*iter)->getId() << std::endl;
188// (*iter)->VelocityVerletUpdateU((*iter)->getId(), CurrentTimeStep-1, Deltat, IsAngstroem);
189 }
190 }
191
192 /** Performs Gradient optimization on the bonds.
193 *
194 * We assume that forces have just been calculated. These forces are projected
195 * onto the bonds and these are annealed subsequently by moving atoms in the
196 * bond neighborhood on either side conjunctively.
197 *
198 *
199 * \param CurrentTimeStep current time step (i.e. \f$ t + \Delta t \f$ in the sense of the velocity verlet)
200 * \param offset offset in matrix file to the first force component
201 * \param maxComponents to be filled with maximum force component over all atoms
202 */
203 void annealWithBondGraph(
204 const int CurrentTimeStep,
205 const size_t offset,
206 Vector &maxComponents)
207 {
208 // get nodes on either side of selected bond via BFS discovery
209// std::vector<atomId_t> atomids;
210// for(typename AtomSetMixin<T>::iterator iter = AtomicForceManipulator<T>::atoms.begin();
211// iter != AtomicForceManipulator<T>::atoms.end(); ++iter) {
212// atomids.push_back((*iter)->getId());
213// }
214// ASSERT( atomids.size() == AtomicForceManipulator<T>::atoms.size(),
215// "ForceAnnealing() - could not gather all atomic ids?");
216 BoostGraphCreator BGcreator;
217 BGcreator.createFromRange(
218 AtomicForceManipulator<T>::atoms.begin(),
219 AtomicForceManipulator<T>::atoms.end(),
220 AtomicForceManipulator<T>::atoms.size(),
221 BreadthFirstSearchGatherer::AlwaysTruePredicate);
222 BreadthFirstSearchGatherer NodeGatherer(BGcreator);
223
224 std::map<atomId_t, Vector> GatheredUpdates; //!< gathers all updates which are applied at the end
225 for(typename AtomSetMixin<T>::iterator iter = AtomicForceManipulator<T>::atoms.begin();
226 iter != AtomicForceManipulator<T>::atoms.end(); ++iter) {
227 // atom's force vector gives steepest descent direction
228 const Vector oldPosition = (*iter)->getPositionAtStep(CurrentTimeStep-2 >= 0 ? CurrentTimeStep - 2 : 0);
229 const Vector currentPosition = (*iter)->getPosition();
230 const Vector oldGradient = (*iter)->getAtomicForceAtStep(CurrentTimeStep-2 >= 0 ? CurrentTimeStep - 2 : 0);
231 const Vector currentGradient = (*iter)->getAtomicForce();
232 LOG(4, "DEBUG: Force for atom " << **iter << " is " << currentGradient);
233
234 // we use Barzilai-Borwein update with position reversed to get descent
235 const Vector GradientDifference = (currentGradient - oldGradient);
236 const double stepwidth =
237 fabs((currentPosition - oldPosition).ScalarProduct(GradientDifference))/
238 GradientDifference.NormSquared();
239 Vector PositionUpdate = stepwidth * currentGradient;
240 if (fabs(stepwidth) < 1e-10) {
241 // dont' warn in first step, deltat usage normal
242 if (currentStep != 1)
243 ELOG(1, "INFO: Barzilai-Borwein stepwidth is zero, using deltat " << currentDeltat << " instead.");
244 PositionUpdate = currentDeltat * currentGradient;
245 }
246 LOG(3, "DEBUG: Update would be " << PositionUpdate);
247
248// // add update to central atom
249// const atomId_t atomid = (*iter)->getId();
250// if (GatheredUpdates.count(atomid)) {
251// GatheredUpdates[atomid] += PositionUpdate;
252// } else
253// GatheredUpdates.insert( std::make_pair(atomid, PositionUpdate) );
254
255 // We assume that a force is local, i.e. a bond is too short yet and hence
256 // the atom needs to be moved. However, all the adjacent (bound) atoms might
257 // already be at the perfect distance. If we just move the atom alone, we ruin
258 // all the other bonds. Hence, it would be sensible to move every atom found
259 // through the bond graph in the direction of the force as well by the same
260 // PositionUpdate. This is just what we are going to do.
261
262 /// get all nodes from bonds in the direction of the current force
263
264 // remove edges facing in the wrong direction
265 std::vector<bond::ptr> removed_bonds;
266 const BondList& ListOfBonds = (*iter)->getListOfBonds();
267 for(BondList::const_iterator bonditer = ListOfBonds.begin();
268 bonditer != ListOfBonds.end(); ++bonditer) {
269 const bond &current_bond = *(*bonditer);
270 LOG(2, "DEBUG: Looking at bond " << current_bond);
271 Vector BondVector = (*iter)->getPosition();
272 BondVector -= ((*iter)->getId() == current_bond.rightatom->getId())
273 ? current_bond.rightatom->getPosition() : current_bond.leftatom->getPosition();
274 BondVector.Normalize();
275 if (BondVector.ScalarProduct(currentGradient) < 0) {
276 removed_bonds.push_back(*bonditer);
277#ifndef NDEBUG
278 const bool status =
279#endif
280 BGcreator.removeEdge(current_bond.leftatom->getId(), current_bond.rightatom->getId());
281 ASSERT( status, "ForceAnnealing() - edge to found bond is not present?");
282 }
283 }
284 BoostGraphHelpers::Nodeset_t bondside_set = NodeGatherer((*iter)->getId(), max_distance);
285 const BreadthFirstSearchGatherer::distance_map_t &distance_map = NodeGatherer.getDistances();
286 std::sort(bondside_set.begin(), bondside_set.end());
287 // re-add those edges
288 for (std::vector<bond::ptr>::const_iterator bonditer = removed_bonds.begin();
289 bonditer != removed_bonds.end(); ++bonditer)
290 BGcreator.addEdge((*bonditer)->leftatom->getId(), (*bonditer)->rightatom->getId());
291
292 // apply PositionUpdate to all nodes in the bondside_set
293 for (BoostGraphHelpers::Nodeset_t::const_iterator setiter = bondside_set.begin();
294 setiter != bondside_set.end(); ++setiter) {
295 const BreadthFirstSearchGatherer::distance_map_t::const_iterator diter
296 = distance_map.find(*setiter);
297 ASSERT( diter != distance_map.end(),
298 "ForceAnnealing() - could not find distance to an atom.");
299 const double factor = pow(damping_factor, diter->second);
300 LOG(3, "DEBUG: Update for atom #" << *setiter << " will be "
301 << factor << "*" << PositionUpdate);
302 if (GatheredUpdates.count((*setiter))) {
303 GatheredUpdates[(*setiter)] += factor*PositionUpdate;
304 } else {
305 GatheredUpdates.insert(
306 std::make_pair(
307 (*setiter),
308 factor*PositionUpdate) );
309 }
310 }
311
312 // extract largest components for showing progress of annealing
313 for(size_t i=0;i<NDIM;++i)
314 if (currentGradient[i] > maxComponents[i])
315 maxComponents[i] = currentGradient[i];
316 }
317 // apply the gathered updates
318 for (std::map<atomId_t, Vector>::const_iterator iter = GatheredUpdates.begin();
319 iter != GatheredUpdates.end(); ++iter) {
320 const atomId_t &atomid = iter->first;
321 const Vector &update = iter->second;
322 atom* const walker = World::getInstance().getAtom(AtomById(atomid));
323 ASSERT( walker != NULL,
324 "ForceAnnealing() - walker with id "+toString(atomid)+" has suddenly disappeared.");
325 LOG(3, "DEBUG: Applying update " << update << " to atom #" << atomid
326 << ", namely " << *walker);
327 walker->setPosition( walker->getPosition() + update );
328 }
329 }
330
331 /** Reset function to unset static entities and artificial velocities.
332 *
333 */
334 void reset()
335 {
336 currentDeltat = 0.;
337 currentStep = 0;
338 }
339
340private:
341 //!> contains the current step in relation to maxsteps
342 static size_t currentStep;
343 //!> contains the maximum number of steps, determines initial and final step with currentStep
344 size_t maxSteps;
345 static double currentDeltat;
346 //!> minimum deltat for internal while loop (adaptive step width)
347 static double MinimumDeltat;
348 //!> contains the maximum bond graph distance up to which shifts of a single atom are spread
349 const int max_distance;
350 //!> the shifted is dampened by this factor with the power of the bond graph distance to the shift causing atom
351 const double damping_factor;
352};
353
354template <class T>
355double ForceAnnealing<T>::currentDeltat = 0.;
356template <class T>
357size_t ForceAnnealing<T>::currentStep = 0;
358template <class T>
359double ForceAnnealing<T>::MinimumDeltat = 1e-8;
360
361#endif /* FORCEANNEALING_HPP_ */
Note: See TracBrowser for help on using the repository browser.