source: src/Atom/atom_atominfo.cpp@ 05a2c5

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

FIX: Fixing the use of the trajectories that were changed to maps.

  • as the atom's trajectories are now stored inside a map, i.e. not every step necessarily needs to be present, there will still some places in the code that were meant for the old vector storing.
  • Atom::UpdateStep() updates the current step, not the future step.
  • AtomInfo::append/removeTrajectoryStep() now no longer ASSERT that the step is not the current one. Due to the trajectory in map concept (with current it simply sought by lower_bound()) this is no longer necessary.
  • HydrogenPool no longer adds all trajectory steps till current time.
  • World has commented out code for also copying the current atom positions in case a new time step is added, but this is not used as so far we do not recognize whether it is a new time step ... and also it is not really necessary.
  • Property mode set to 100644
File size: 17.4 KB
RevLine 
[bcf653]1/*
2 * Project: MoleCuilder
3 * Description: creates and alters molecular systems
[0aa122]4 * Copyright (C) 2010-2012 University of Bonn. All rights reserved.
[8cc22f]5 * Copyright (C) 2014 Frederik Heber. All rights reserved.
[94d5ac6]6 *
7 *
8 * This file is part of MoleCuilder.
9 *
10 * MoleCuilder is free software: you can redistribute it and/or modify
11 * it under the terms of the GNU General Public License as published by
12 * the Free Software Foundation, either version 2 of the License, or
13 * (at your option) any later version.
14 *
15 * MoleCuilder is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU General Public License for more details.
19 *
20 * You should have received a copy of the GNU General Public License
21 * along with MoleCuilder. If not, see <http://www.gnu.org/licenses/>.
[bcf653]22 */
23
[6b919f8]24/*
25 * atom_atominfo.cpp
26 *
27 * Created on: Oct 19, 2009
28 * Author: heber
29 */
30
[bf3817]31// include config.h
32#ifdef HAVE_CONFIG_H
33#include <config.h>
34#endif
35
[9eb71b3]36//#include "CodePatterns/MemDebug.hpp"
[112b09]37
[6625c3]38#include "CodePatterns/Verbose.hpp"
[08a0f52]39
40#include "atom_atominfo.hpp"
41#include "CodePatterns/Log.hpp"
[6625c3]42#include "config.hpp"
[3bdb6d]43#include "Element/element.hpp"
44#include "Element/periodentafel.hpp"
[a9b86d]45#include "Fragmentation/ForceMatrix.hpp"
[f16a4b]46#include "World.hpp"
[1b558c]47#include "WorldTime.hpp"
[6b919f8]48
[435065]49#include <iomanip>
50
[6b919f8]51/** Constructor of class AtomInfo.
52 */
[97b825]53AtomInfo::AtomInfo() :
[606c24]54 AtomicElement(1),
[2034f3]55 FixedIon(false),
56 charge(0.)
[54b42e]57{
[8cc22f]58 AtomicPosition.insert( std::make_pair(0, zeroVec) );
59 AtomicVelocity.insert( std::make_pair(0, zeroVec) );
60 AtomicForce.insert( std::make_pair(0, zeroVec) );
61}
[d74077]62
63/** Copy constructor of class AtomInfo.
64 */
[97b825]65AtomInfo::AtomInfo(const AtomInfo &_atom) :
66 AtomicPosition(_atom.AtomicPosition),
[8cc22f]67 AtomicVelocity(_atom.AtomicVelocity),
68 AtomicForce(_atom.AtomicForce),
[6625c3]69 AtomicElement(_atom.AtomicElement),
[2034f3]70 FixedIon(_atom.FixedIon),
71 charge(_atom.charge)
[6625c3]72{
[8cc22f]73}
[d74077]74
[97b825]75AtomInfo::AtomInfo(const VectorInterface &_v) :
[606c24]76 AtomicElement(1),
[2034f3]77 FixedIon(false),
78 charge(0.)
[54b42e]79{
[8cc22f]80 AtomicPosition.insert( std::make_pair(0, _v.getPosition()) );
81 AtomicVelocity.insert( std::make_pair(0, zeroVec) );
82 AtomicForce.insert( std::make_pair(0, zeroVec) );
[54b42e]83};
[6b919f8]84
85/** Destructor of class AtomInfo.
86 */
87AtomInfo::~AtomInfo()
88{
89};
90
[8cc22f]91void AtomInfo::AppendTrajectoryStep(const unsigned int _step)
[e2373df]92{
[fe0ddf]93 NOTIFY(TrajectoryChanged);
[8cc22f]94 AtomicPosition.insert( std::make_pair(_step, zeroVec) );
95 AtomicVelocity.insert( std::make_pair(_step, zeroVec) );
96 AtomicForce.insert( std::make_pair(_step, zeroVec) );
[fb9eba]97 LOG(5,"AtomInfo::AppendTrajectoryStep() called, size is ("
98 << AtomicPosition.size() << ","
99 << AtomicVelocity.size() << ","
100 << AtomicForce.size() << ")");
[e2373df]101}
102
[8cc22f]103void AtomInfo::removeTrajectoryStep(const unsigned int _step)
[7e51e1]104{
[fe0ddf]105 NOTIFY(TrajectoryChanged);
[8cc22f]106 AtomicPosition.erase(_step);
107 AtomicVelocity.erase(_step);
108 AtomicForce.erase(_step);
[7e51e1]109 LOG(5,"AtomInfo::removeTrajectoryStep() called, size is ("
110 << AtomicPosition.size() << ","
111 << AtomicVelocity.size() << ","
112 << AtomicForce.size() << ")");
113}
114
[d74077]115const element *AtomInfo::getType() const
116{
[35a25a]117 const element *elem = World::getInstance().getPeriode()->FindElement(AtomicElement);
118 return elem;
[d74077]119}
120
[08a0f52]121const element &AtomInfo::getElement() const
122{
123 const element &elem = *World::getInstance().getPeriode()->FindElement(AtomicElement);
124 return elem;
125}
126
127atomicNumber_t AtomInfo::getElementNo() const
128{
129 return AtomicElement;
130}
131
[843590]132const std::string &AtomInfo::getParticleName() const
133{
134 return particlename;
135}
136
137void AtomInfo::setParticleName(const std::string & _name)
138{
139 particlename = _name;
140}
141
[d74077]142const double& AtomInfo::operator[](size_t i) const
143{
[8cc22f]144 return atStep(i, WorldTime::getTime());
[d74077]145}
146
147const double& AtomInfo::at(size_t i) const
148{
[8cc22f]149 return atStep(i, WorldTime::getTime());
[d74077]150}
151
[6b020f]152const double& AtomInfo::atStep(size_t i, unsigned int _step) const
[056e70]153{
[8cc22f]154 ASSERT(!AtomicPosition.empty(),
155 "AtomInfo::operator[]() - AtomicPosition is empty.");
156 VectorTrajectory_t::const_iterator iter =
157 AtomicPosition.lower_bound(_step);
158 return iter->second[i];
[056e70]159}
160
[d74077]161void AtomInfo::set(size_t i, const double value)
162{
[7188b1]163 OBSERVE;
164 NOTIFY(AtomObservable::PositionChanged);
[8cc22f]165 VectorTrajectory_t::iterator iter = AtomicPosition.find(WorldTime::getTime());
166 if (iter != AtomicPosition.end()) {
167 iter->second[i] = value;
168 } else {
169 Vector newPos;
170 newPos[i] = value;
171#ifndef NDEBUG
172 std::pair<VectorTrajectory_t::iterator, bool> inserter =
173#endif
174 AtomicPosition.insert( std::make_pair(WorldTime::getTime(), newPos) );
175 ASSERT( inserter.second,
176 "AtomInfo::set() - time step "+toString(WorldTime::getTime())
177 +" present after all?");
178 }
179}
180
181/** Helps to determine whether the current step really exists or getPosition() has just
182 * delivered the one closest to it in the past.
183 *
184 * \param _step step to check
185 * \param true - step exists, false - step does not exist, getPosition() delivers closest
186 */
187bool AtomInfo::isStepPresent(const unsigned int _step) const
188{
189 VectorTrajectory_t::const_iterator iter =
190 AtomicPosition.find(_step);
191 return iter != AtomicPosition.end();
[d74077]192}
193
194const Vector& AtomInfo::getPosition() const
195{
[8cc22f]196 return getPositionAtStep(WorldTime::getTime());
[f16a4b]197}
198
[6b020f]199const Vector& AtomInfo::getPositionAtStep(const unsigned int _step) const
[6625c3]200{
[8cc22f]201 ASSERT(!AtomicPosition.empty(),
202 "AtomInfo::operator[]() - AtomicPosition is empty.");
203 VectorTrajectory_t::const_iterator iter =
204 AtomicPosition.lower_bound(_step);
205 return iter->second;
[6625c3]206}
207
[7188b1]208void AtomInfo::setType(const element* _type)
209{
[fe0ddf]210 OBSERVE;
211 NOTIFY(AtomObservable::ElementChanged);
212 AtomicElement = _type->getAtomicNumber();
[f16a4b]213}
214
[7188b1]215void AtomInfo::setType(const int Z)
216{
[ead4e6]217 const element *elem = World::getInstance().getPeriode()->FindElement(Z);
[8cc22f]218 setType(elem);
[f16a4b]219}
[d74077]220
[bce72c]221const Vector& AtomInfo::getAtomicVelocity() const
222{
[8cc22f]223 return getAtomicVelocityAtStep(WorldTime::getTime());
[bce72c]224}
225
[6b020f]226const Vector& AtomInfo::getAtomicVelocityAtStep(const unsigned int _step) const
[6625c3]227{
[8cc22f]228 ASSERT(!AtomicVelocity.empty(),
229 "AtomInfo::operator[]() - AtomicVelocity is empty.");
230 VectorTrajectory_t::const_iterator iter =
231 AtomicVelocity.lower_bound(_step);
232 // special, we only interpolate between present time steps not into the future
233 if (_step > AtomicVelocity.begin()->first)
234 return zeroVec;
235 else
236 return iter->second;
[6625c3]237}
238
[bce72c]239void AtomInfo::setAtomicVelocity(const Vector &_newvelocity)
240{
[8cc22f]241 setAtomicVelocityAtStep(WorldTime::getTime(), _newvelocity);
[bce72c]242}
243
[6b020f]244void AtomInfo::setAtomicVelocityAtStep(const unsigned int _step, const Vector &_newvelocity)
[6625c3]245{
[7188b1]246 OBSERVE;
[8cc22f]247 VectorTrajectory_t::iterator iter = AtomicVelocity.find(_step);
248 if (iter != AtomicVelocity.end()) {
249 iter->second = _newvelocity;
250 } else {
251#ifndef NDEBUG
252 std::pair<VectorTrajectory_t::iterator, bool> inserter =
253#endif
254 AtomicVelocity.insert( std::make_pair(_step, _newvelocity) );
255 ASSERT( inserter.second,
256 "AtomInfo::set() - time step "+toString(_step)
257 +" present after all?");
258 }
[7188b1]259 if (WorldTime::getTime() == _step)
260 NOTIFY(AtomObservable::VelocityChanged);
[6625c3]261}
262
[bce72c]263const Vector& AtomInfo::getAtomicForce() const
264{
[8cc22f]265 return getAtomicForceAtStep(WorldTime::getTime());
[bce72c]266}
267
[6b020f]268const Vector& AtomInfo::getAtomicForceAtStep(const unsigned int _step) const
[6625c3]269{
[8cc22f]270 ASSERT(!AtomicForce.empty(),
271 "AtomInfo::operator[]() - AtomicForce is empty.");
272 VectorTrajectory_t::const_iterator iter =
273 AtomicForce.lower_bound(_step);
274 // special, we only interpolate between present time steps not into the future
275 if (_step > AtomicForce.begin()->first)
276 return zeroVec;
277 else
278 return iter->second;
[6625c3]279}
280
[bce72c]281void AtomInfo::setAtomicForce(const Vector &_newforce)
282{
[8cc22f]283 setAtomicForceAtStep(WorldTime::getTime(), _newforce);
[bce72c]284}
285
[6b020f]286void AtomInfo::setAtomicForceAtStep(const unsigned int _step, const Vector &_newforce)
[6625c3]287{
[7188b1]288 OBSERVE;
[8cc22f]289 VectorTrajectory_t::iterator iter = AtomicForce.find(_step);
290 if (iter != AtomicForce.end()) {
291 iter->second = _newforce;
292 } else {
293#ifndef NDEBUG
294 std::pair<VectorTrajectory_t::iterator, bool> inserter =
295#endif
296 AtomicForce.insert( std::make_pair(_step, _newforce) );
297 ASSERT( inserter.second,
298 "AtomInfo::set() - time step "+toString(_step)
299 +" present after all?");
300 }
[7188b1]301 if (WorldTime::getTime() == _step)
[bcb593]302 NOTIFY(AtomObservable::ForceChanged);
[6625c3]303}
304
305bool AtomInfo::getFixedIon() const
306{
307 return FixedIon;
308}
309
310void AtomInfo::setFixedIon(const bool _fixedion)
311{
[7188b1]312 OBSERVE;
313 NOTIFY(AtomObservable::PropertyChanged);
[6625c3]314 FixedIon = _fixedion;
315}
316
[d74077]317void AtomInfo::setPosition(const Vector& _vector)
318{
[8cc22f]319 setPositionAtStep(WorldTime::getTime(), _vector);
[d74077]320}
321
[6b020f]322void AtomInfo::setPositionAtStep(unsigned int _step, const Vector& _vector)
[6625c3]323{
[7188b1]324 OBSERVE;
[8cc22f]325 VectorTrajectory_t::iterator iter = AtomicPosition.find(_step);
326 if (iter != AtomicPosition.end()) {
327 iter->second = _vector;
328 } else {
329#ifndef NDEBUG
330 std::pair<VectorTrajectory_t::iterator, bool> inserter =
331#endif
332 AtomicPosition.insert( std::make_pair(_step, _vector) );
333 ASSERT( inserter.second,
334 "AtomInfo::set() - time step "+toString(_step)
335 +" present after all?");
336 }
[7188b1]337 if (WorldTime::getTime() == _step)
338 NOTIFY(AtomObservable::PositionChanged);
[6625c3]339}
340
[d74077]341const VectorInterface& AtomInfo::operator+=(const Vector& b)
342{
[8cc22f]343 setPosition(getPosition()+b);
[d74077]344 return *this;
345}
346
347const VectorInterface& AtomInfo::operator-=(const Vector& b)
348{
[8cc22f]349 setPosition(getPosition()-b);
[d74077]350 return *this;
351}
352
353Vector const AtomInfo::operator+(const Vector& b) const
354{
[8cc22f]355 Vector a(getPosition());
[d74077]356 a += b;
357 return a;
358}
359
360Vector const AtomInfo::operator-(const Vector& b) const
361{
[8cc22f]362 Vector a(getPosition());
[d74077]363 a -= b;
364 return a;
365}
366
367double AtomInfo::distance(const Vector &point) const
368{
[8cc22f]369 return getPosition().distance(point);
[d74077]370}
371
372double AtomInfo::DistanceSquared(const Vector &y) const
373{
[8cc22f]374 return getPosition().DistanceSquared(y);
[d74077]375}
376
377double AtomInfo::distance(const VectorInterface &_atom) const
378{
[8cc22f]379 return _atom.distance(getPosition());
[d74077]380}
381
382double AtomInfo::DistanceSquared(const VectorInterface &_atom) const
383{
[8cc22f]384 return _atom.DistanceSquared(getPosition());
[d74077]385}
386
387VectorInterface &AtomInfo::operator=(const Vector& _vector)
388{
[8cc22f]389 setPosition(_vector);
[d74077]390 return *this;
391}
392
393void AtomInfo::ScaleAll(const double *factor)
394{
[8cc22f]395 Vector temp(getPosition());
396 temp.ScaleAll(factor);
397 setPosition(temp);
[d74077]398}
399
400void AtomInfo::ScaleAll(const Vector &factor)
401{
[8cc22f]402 Vector temp(getPosition());
403 temp.ScaleAll(factor);
404 setPosition(temp);
[d74077]405}
406
407void AtomInfo::Scale(const double factor)
408{
[8cc22f]409 Vector temp(getPosition());
410 temp.Scale(factor);
411 setPosition(temp);
[d74077]412}
413
414void AtomInfo::Zero()
415{
[8cc22f]416 setPosition(zeroVec);
[d74077]417}
418
419void AtomInfo::One(const double one)
420{
[8cc22f]421 setPosition(Vector(one,one,one));
[d74077]422}
423
424void AtomInfo::LinearCombinationOfVectors(const Vector &x1, const Vector &x2, const Vector &x3, const double * const factors)
425{
[8cc22f]426 Vector newPos;
427 newPos.LinearCombinationOfVectors(x1,x2,x3,factors);
428 setPosition(newPos);
[d74077]429}
430
[6625c3]431/**
432 * returns the kinetic energy of this atom at a given time step
433 */
[7188b1]434double AtomInfo::getKineticEnergy(const unsigned int _step) const
435{
[8cc22f]436 return getMass() * getAtomicVelocityAtStep(_step).NormSquared();
[6625c3]437}
438
[7188b1]439Vector AtomInfo::getMomentum(const unsigned int _step) const
440{
[8cc22f]441 return getMass() * getAtomicVelocityAtStep(_step);
[6625c3]442}
443
[8cc22f]444/** Decrease the trajectory if given \a MaxSteps is smaller.
445 * Does nothing if \a MaxSteps is larger than current size.
446 *
[6625c3]447 * \param MaxSteps
448 */
449void AtomInfo::ResizeTrajectory(size_t MaxSteps)
450{
[8cc22f]451 // mind the reverse ordering due to std::greater, latest time steps are at beginning
452 VectorTrajectory_t::iterator positer = AtomicPosition.lower_bound(MaxSteps);
453 if (positer != AtomicPosition.begin()) {
454 if (positer->first == MaxSteps)
455 --positer;
456 AtomicPosition.erase(AtomicPosition.begin(), positer);
457 }
458 VectorTrajectory_t::iterator veliter = AtomicVelocity.lower_bound(MaxSteps);
459 if (veliter != AtomicVelocity.begin()) {
460 if (veliter->first == MaxSteps)
461 --veliter;
462 AtomicVelocity.erase(AtomicVelocity.begin(), veliter);
463 }
464 VectorTrajectory_t::iterator forceiter = AtomicForce.lower_bound(MaxSteps);
465 if (forceiter != AtomicForce.begin()) {
466 if (forceiter->first == MaxSteps)
467 --forceiter;
468 AtomicForce.erase(AtomicForce.begin(), forceiter);
469 }
[e2373df]470}
[6625c3]471
472size_t AtomInfo::getTrajectorySize() const
473{
[8cc22f]474 // mind greater comp for map here: first element is latest in time steps!
475 return AtomicPosition.begin()->first+1;
[6625c3]476}
477
[35a25a]478double AtomInfo::getMass() const
479{
480 return getType()->getMass();
[6625c3]481}
482
[8cc22f]483/** Helper function to either insert or assign, depending on the element being
484 * present already.
485 *
486 * \param _trajectory vector of Vectors to assign
487 * \param dest step to insert/assign to
488 * \param _newvalue new Vector value
489 */
490void assignTrajectoryElement(
491 std::map<unsigned int, Vector, std::greater<unsigned int> > &_trajectory,
492 const unsigned int dest,
493 const Vector &_newvalue)
494{
495 std::pair<std::map<unsigned int, Vector, std::greater<unsigned int> >::iterator, bool> inserter =
496 _trajectory.insert( std::make_pair(dest, _newvalue) );
497 if (!inserter.second)
498 inserter.first->second = _newvalue;
499}
500
[6625c3]501/** Copies a given trajectory step \a src onto another \a dest
502 * \param dest index of destination step
503 * \param src index of source step
504 */
[6b020f]505void AtomInfo::CopyStepOnStep(const unsigned int dest, const unsigned int src)
[6625c3]506{
507 if (dest == src) // self assignment check
508 return;
509
[7188b1]510 if (WorldTime::getTime() == dest){
511 NOTIFY(AtomObservable::PositionChanged);
512 NOTIFY(AtomObservable::VelocityChanged);
513 NOTIFY(AtomObservable::ForceChanged);
514 }
515
[8cc22f]516 VectorTrajectory_t::iterator positer = AtomicPosition.find(src);
517 ASSERT( positer != AtomicPosition.end(),
518 "AtomInfo::CopyStepOnStep() - step "
519 +toString(src)+" to copy from not present in AtomicPosition.");
[f946b2]520 assignTrajectoryElement(AtomicPosition, dest, positer->second);
[8cc22f]521 VectorTrajectory_t::iterator veliter = AtomicVelocity.find(src);
[f946b2]522 if (veliter != AtomicVelocity.end())
523 assignTrajectoryElement(AtomicVelocity, dest, veliter->second);
[8cc22f]524 VectorTrajectory_t::iterator forceiter = AtomicForce.find(src);
[f946b2]525 if (forceiter != AtomicForce.end())
526 assignTrajectoryElement(AtomicForce, dest, forceiter->second);
[6625c3]527};
528
[bcb593]529/** Performs a velocity verlet update of the position at \a NextStep from \a LastStep information only.
530 *
531 * We calculate \f$x(t + \delta t) = x(t) + v(t)* \delta t + .5 * \delta t * \delta t * F(t)/m \f$.
532 *
533 *
[6625c3]534 * \param NextStep index of sequential step to set
[435065]535 * \param Deltat time step width
536 * \param IsAngstroem whether the force's underlying unit of length is angstroem or bohr radii
[6625c3]537 */
[bcb593]538void AtomInfo::VelocityVerletUpdateX(int nr, const unsigned int NextStep, double Deltat, bool IsAngstroem)
[6625c3]539{
[4882d5]540 const unsigned int LastStep = NextStep == 0 ? 0 : NextStep-1;
541
[435065]542 LOG(2, "INFO: Particle that currently " << *this);
[bcb593]543 LOG(2, "INFO: Integrating position with mass=" << getMass() << " and Deltat="
[435065]544 << Deltat << " at NextStep=" << NextStep);
[056e70]545
546 // update position
[4882d5]547 {
548 Vector tempVector = getPositionAtStep(LastStep);
549 LOG(4, "INFO: initial position from last step " << setprecision(4) << tempVector);
550 tempVector += Deltat*(getAtomicVelocityAtStep(LastStep)); // s(t) = s(0) + v * deltat + 1/2 a * deltat^2
551 LOG(4, "INFO: position with velocity " << getAtomicVelocityAtStep(LastStep) << " from last step " << tempVector);
552 tempVector += .5*Deltat*Deltat*(getAtomicForceAtStep(LastStep))*(1./getMass()); // F = m * a and s =
[bcb593]553 LOG(4, "INFO: position with force " << getAtomicForceAtStep(LastStep) << " from last step " << tempVector);
[4882d5]554 setPositionAtStep(NextStep, tempVector);
555 LOG(3, "INFO: Position at step " << NextStep << " set to " << tempVector);
556 }
[bcb593]557};
558
559/** Performs a velocity verlet update of the velocity at \a NextStep.
560 *
561 * \note forces at NextStep should have been calculated based on position at NextStep prior
562 * to calling this function.
563 *
564 * We calculate \f$v(t) = v(t - \delta t) + \delta _t * .5 * (F(t - \delta t) + F(t))/m \f$.
565 *
566 * Parameters are according to those in configuration class.
567 * \param NextStep index of sequential step to set
568 * \param Deltat time step width
569 * \param IsAngstroem whether the force's underlying unit of length is angstroem or bohr radii
570 */
571void AtomInfo::VelocityVerletUpdateU(int nr, const unsigned int NextStep, double Deltat, bool IsAngstroem)
572{
573 const unsigned int LastStep = NextStep == 0 ? 0 : NextStep-1;
574
575 LOG(2, "INFO: Particle that currently " << *this);
576 LOG(2, "INFO: Integrating velocity with mass=" << getMass() << " and Deltat="
577 << Deltat << " at NextStep=" << NextStep);
[056e70]578
[6625c3]579 // Update U
[4882d5]580 {
581 Vector tempVector = getAtomicVelocityAtStep(LastStep);
582 LOG(4, "INFO: initial velocity from last step " << tempVector);
583 tempVector += Deltat * .5*(getAtomicForceAtStep(LastStep)+getAtomicForceAtStep(NextStep))*(1./getMass()); // v = F/m * t
584 LOG(4, "INFO: Velocity with force from last " << getAtomicForceAtStep(LastStep)
585 << " and present " << getAtomicForceAtStep(NextStep) << " step " << tempVector);
586 setAtomicVelocityAtStep(NextStep, tempVector);
587 LOG(3, "INFO: Velocity at step " << NextStep << " set to " << tempVector);
588 }
[6625c3]589};
590
[d74077]591std::ostream & AtomInfo::operator << (std::ostream &ost) const
592{
593 return (ost << getPosition());
594}
595
596std::ostream & operator << (std::ostream &ost, const AtomInfo &a)
597{
[b1a5d9]598 const size_t terminalstep = a.getTrajectorySize()-1;
[e34254]599 if (terminalstep) {
600 ost << "starts at "
601 << a.getPositionAtStep(0) << " and ends at "
602 << a.getPositionAtStep(terminalstep)
603 << " at time step " << terminalstep;
604 } else {
605 ost << "is at "
606 << a.getPositionAtStep(0) << " with a single time step only";
607 }
[d74077]608 return ost;
609}
610
Note: See TracBrowser for help on using the repository browser.