| [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] | 53 | AtomInfo::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] | 65 | AtomInfo::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] | 75 | AtomInfo::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 | */ | 
|---|
|  | 87 | AtomInfo::~AtomInfo() | 
|---|
|  | 88 | { | 
|---|
|  | 89 | }; | 
|---|
|  | 90 |  | 
|---|
| [8cc22f] | 91 | void 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 |  | 
|---|
| [8c6b68] | 103 | void AtomInfo::eraseInTrajctory( | 
|---|
|  | 104 | VectorTrajectory_t &_trajectory, | 
|---|
|  | 105 | const unsigned int _firststep, const unsigned int _laststep) | 
|---|
|  | 106 | { | 
|---|
|  | 107 | const VectorTrajectory_t::iterator firstiter = _trajectory.lower_bound(_firststep); | 
|---|
|  | 108 | const VectorTrajectory_t::iterator lastiter = _trajectory.upper_bound(_laststep); | 
|---|
|  | 109 | _trajectory.erase(firstiter, lastiter); | 
|---|
|  | 110 | } | 
|---|
|  | 111 |  | 
|---|
|  | 112 | void AtomInfo::removeTrajectorySteps(const unsigned int _firststep, const unsigned int _laststep) | 
|---|
| [7e51e1] | 113 | { | 
|---|
| [fe0ddf] | 114 | NOTIFY(TrajectoryChanged); | 
|---|
| [8c6b68] | 115 | eraseInTrajctory(AtomicPosition, _firststep, _laststep); | 
|---|
|  | 116 | eraseInTrajctory(AtomicVelocity, _firststep, _laststep); | 
|---|
|  | 117 | eraseInTrajctory(AtomicForce, _firststep, _laststep); | 
|---|
|  | 118 | LOG(5,"AtomInfo::removeTrajectorySteps() called, size is (" | 
|---|
| [7e51e1] | 119 | << AtomicPosition.size() << "," | 
|---|
|  | 120 | << AtomicVelocity.size() << "," | 
|---|
|  | 121 | << AtomicForce.size() << ")"); | 
|---|
|  | 122 | } | 
|---|
|  | 123 |  | 
|---|
| [d74077] | 124 | const element *AtomInfo::getType() const | 
|---|
|  | 125 | { | 
|---|
| [35a25a] | 126 | const element *elem = World::getInstance().getPeriode()->FindElement(AtomicElement); | 
|---|
|  | 127 | return elem; | 
|---|
| [d74077] | 128 | } | 
|---|
|  | 129 |  | 
|---|
| [08a0f52] | 130 | const element &AtomInfo::getElement() const | 
|---|
|  | 131 | { | 
|---|
|  | 132 | const element &elem = *World::getInstance().getPeriode()->FindElement(AtomicElement); | 
|---|
|  | 133 | return elem; | 
|---|
|  | 134 | } | 
|---|
|  | 135 |  | 
|---|
|  | 136 | atomicNumber_t AtomInfo::getElementNo() const | 
|---|
|  | 137 | { | 
|---|
|  | 138 | return AtomicElement; | 
|---|
|  | 139 | } | 
|---|
|  | 140 |  | 
|---|
| [843590] | 141 | const std::string &AtomInfo::getParticleName() const | 
|---|
|  | 142 | { | 
|---|
|  | 143 | return particlename; | 
|---|
|  | 144 | } | 
|---|
|  | 145 |  | 
|---|
|  | 146 | void AtomInfo::setParticleName(const std::string & _name) | 
|---|
|  | 147 | { | 
|---|
|  | 148 | particlename = _name; | 
|---|
|  | 149 | } | 
|---|
|  | 150 |  | 
|---|
| [d74077] | 151 | const double& AtomInfo::operator[](size_t i) const | 
|---|
|  | 152 | { | 
|---|
| [8cc22f] | 153 | return atStep(i, WorldTime::getTime()); | 
|---|
| [d74077] | 154 | } | 
|---|
|  | 155 |  | 
|---|
|  | 156 | const double& AtomInfo::at(size_t i) const | 
|---|
|  | 157 | { | 
|---|
| [8cc22f] | 158 | return atStep(i, WorldTime::getTime()); | 
|---|
| [d74077] | 159 | } | 
|---|
|  | 160 |  | 
|---|
| [6b020f] | 161 | const double& AtomInfo::atStep(size_t i, unsigned int _step) const | 
|---|
| [056e70] | 162 | { | 
|---|
| [8cc22f] | 163 | ASSERT(!AtomicPosition.empty(), | 
|---|
|  | 164 | "AtomInfo::operator[]() - AtomicPosition is empty."); | 
|---|
|  | 165 | VectorTrajectory_t::const_iterator iter = | 
|---|
|  | 166 | AtomicPosition.lower_bound(_step); | 
|---|
|  | 167 | return iter->second[i]; | 
|---|
| [056e70] | 168 | } | 
|---|
|  | 169 |  | 
|---|
| [d74077] | 170 | void AtomInfo::set(size_t i, const double value) | 
|---|
|  | 171 | { | 
|---|
| [7188b1] | 172 | OBSERVE; | 
|---|
|  | 173 | NOTIFY(AtomObservable::PositionChanged); | 
|---|
| [8cc22f] | 174 | VectorTrajectory_t::iterator iter = AtomicPosition.find(WorldTime::getTime()); | 
|---|
|  | 175 | if (iter !=  AtomicPosition.end()) { | 
|---|
|  | 176 | iter->second[i] = value; | 
|---|
|  | 177 | } else { | 
|---|
|  | 178 | Vector newPos; | 
|---|
|  | 179 | newPos[i] = value; | 
|---|
|  | 180 | #ifndef NDEBUG | 
|---|
|  | 181 | std::pair<VectorTrajectory_t::iterator, bool> inserter = | 
|---|
|  | 182 | #endif | 
|---|
|  | 183 | AtomicPosition.insert( std::make_pair(WorldTime::getTime(), newPos) ); | 
|---|
|  | 184 | ASSERT( inserter.second, | 
|---|
|  | 185 | "AtomInfo::set() - time step "+toString(WorldTime::getTime()) | 
|---|
|  | 186 | +" present after all?"); | 
|---|
|  | 187 | } | 
|---|
|  | 188 | } | 
|---|
|  | 189 |  | 
|---|
| [95b64f] | 190 | void AtomInfo::setAtStep(size_t i, unsigned int _step, const double value) | 
|---|
|  | 191 | { | 
|---|
|  | 192 | OBSERVE; | 
|---|
|  | 193 | NOTIFY(AtomObservable::PositionChanged); | 
|---|
|  | 194 | VectorTrajectory_t::iterator iter = AtomicPosition.find(_step); | 
|---|
|  | 195 | if (iter !=  AtomicPosition.end()) { | 
|---|
|  | 196 | iter->second[i] = value; | 
|---|
|  | 197 | } else { | 
|---|
|  | 198 | Vector newPos; | 
|---|
|  | 199 | newPos[i] = value; | 
|---|
|  | 200 | #ifndef NDEBUG | 
|---|
|  | 201 | std::pair<VectorTrajectory_t::iterator, bool> inserter = | 
|---|
|  | 202 | #endif | 
|---|
|  | 203 | AtomicPosition.insert( std::make_pair(_step, newPos) ); | 
|---|
|  | 204 | ASSERT( inserter.second, | 
|---|
|  | 205 | "AtomInfo::setAtStep() - time step "+toString(_step) | 
|---|
|  | 206 | +" present after all?"); | 
|---|
|  | 207 | } | 
|---|
|  | 208 | } | 
|---|
|  | 209 |  | 
|---|
| [8cc22f] | 210 | /** Helps to determine whether the current step really exists or getPosition() has just | 
|---|
|  | 211 | * delivered the one closest to it in the past. | 
|---|
|  | 212 | * | 
|---|
|  | 213 | * \param _step step to check | 
|---|
|  | 214 | * \param true - step exists, false - step does not exist, getPosition() delivers closest | 
|---|
|  | 215 | */ | 
|---|
|  | 216 | bool AtomInfo::isStepPresent(const unsigned int _step) const | 
|---|
|  | 217 | { | 
|---|
|  | 218 | VectorTrajectory_t::const_iterator iter = | 
|---|
|  | 219 | AtomicPosition.find(_step); | 
|---|
|  | 220 | return iter != AtomicPosition.end(); | 
|---|
| [d74077] | 221 | } | 
|---|
|  | 222 |  | 
|---|
|  | 223 | const Vector& AtomInfo::getPosition() const | 
|---|
|  | 224 | { | 
|---|
| [8cc22f] | 225 | return getPositionAtStep(WorldTime::getTime()); | 
|---|
| [f16a4b] | 226 | } | 
|---|
|  | 227 |  | 
|---|
| [6b020f] | 228 | const Vector& AtomInfo::getPositionAtStep(const unsigned int _step) const | 
|---|
| [6625c3] | 229 | { | 
|---|
| [8cc22f] | 230 | ASSERT(!AtomicPosition.empty(), | 
|---|
|  | 231 | "AtomInfo::operator[]() - AtomicPosition is empty."); | 
|---|
|  | 232 | VectorTrajectory_t::const_iterator iter = | 
|---|
|  | 233 | AtomicPosition.lower_bound(_step); | 
|---|
|  | 234 | return iter->second; | 
|---|
| [6625c3] | 235 | } | 
|---|
|  | 236 |  | 
|---|
| [7188b1] | 237 | void AtomInfo::setType(const element* _type) | 
|---|
|  | 238 | { | 
|---|
| [fe0ddf] | 239 | OBSERVE; | 
|---|
|  | 240 | NOTIFY(AtomObservable::ElementChanged); | 
|---|
|  | 241 | AtomicElement = _type->getAtomicNumber(); | 
|---|
| [f16a4b] | 242 | } | 
|---|
|  | 243 |  | 
|---|
| [7188b1] | 244 | void AtomInfo::setType(const int Z) | 
|---|
|  | 245 | { | 
|---|
| [ead4e6] | 246 | const element *elem = World::getInstance().getPeriode()->FindElement(Z); | 
|---|
| [8cc22f] | 247 | setType(elem); | 
|---|
| [f16a4b] | 248 | } | 
|---|
| [d74077] | 249 |  | 
|---|
| [bce72c] | 250 | const Vector& AtomInfo::getAtomicVelocity() const | 
|---|
|  | 251 | { | 
|---|
| [8cc22f] | 252 | return getAtomicVelocityAtStep(WorldTime::getTime()); | 
|---|
| [bce72c] | 253 | } | 
|---|
|  | 254 |  | 
|---|
| [6b020f] | 255 | const Vector& AtomInfo::getAtomicVelocityAtStep(const unsigned int _step) const | 
|---|
| [6625c3] | 256 | { | 
|---|
| [8cc22f] | 257 | ASSERT(!AtomicVelocity.empty(), | 
|---|
|  | 258 | "AtomInfo::operator[]() - AtomicVelocity is empty."); | 
|---|
|  | 259 | VectorTrajectory_t::const_iterator iter = | 
|---|
|  | 260 | AtomicVelocity.lower_bound(_step); | 
|---|
|  | 261 | // special, we only interpolate between present time steps not into the future | 
|---|
|  | 262 | if (_step > AtomicVelocity.begin()->first) | 
|---|
|  | 263 | return zeroVec; | 
|---|
|  | 264 | else | 
|---|
|  | 265 | return iter->second; | 
|---|
| [6625c3] | 266 | } | 
|---|
|  | 267 |  | 
|---|
| [bce72c] | 268 | void AtomInfo::setAtomicVelocity(const Vector &_newvelocity) | 
|---|
|  | 269 | { | 
|---|
| [8cc22f] | 270 | setAtomicVelocityAtStep(WorldTime::getTime(), _newvelocity); | 
|---|
| [bce72c] | 271 | } | 
|---|
|  | 272 |  | 
|---|
| [6b020f] | 273 | void AtomInfo::setAtomicVelocityAtStep(const unsigned int _step, const Vector &_newvelocity) | 
|---|
| [6625c3] | 274 | { | 
|---|
| [7188b1] | 275 | OBSERVE; | 
|---|
| [8cc22f] | 276 | VectorTrajectory_t::iterator iter = AtomicVelocity.find(_step); | 
|---|
|  | 277 | if (iter !=  AtomicVelocity.end()) { | 
|---|
|  | 278 | iter->second = _newvelocity; | 
|---|
|  | 279 | } else { | 
|---|
|  | 280 | #ifndef NDEBUG | 
|---|
|  | 281 | std::pair<VectorTrajectory_t::iterator, bool> inserter = | 
|---|
|  | 282 | #endif | 
|---|
|  | 283 | AtomicVelocity.insert( std::make_pair(_step, _newvelocity) ); | 
|---|
|  | 284 | ASSERT( inserter.second, | 
|---|
|  | 285 | "AtomInfo::set() - time step "+toString(_step) | 
|---|
|  | 286 | +" present after all?"); | 
|---|
|  | 287 | } | 
|---|
| [7188b1] | 288 | if (WorldTime::getTime() == _step) | 
|---|
|  | 289 | NOTIFY(AtomObservable::VelocityChanged); | 
|---|
| [6625c3] | 290 | } | 
|---|
|  | 291 |  | 
|---|
| [bce72c] | 292 | const Vector& AtomInfo::getAtomicForce() const | 
|---|
|  | 293 | { | 
|---|
| [8cc22f] | 294 | return getAtomicForceAtStep(WorldTime::getTime()); | 
|---|
| [bce72c] | 295 | } | 
|---|
|  | 296 |  | 
|---|
| [6b020f] | 297 | const Vector& AtomInfo::getAtomicForceAtStep(const unsigned int _step) const | 
|---|
| [6625c3] | 298 | { | 
|---|
| [8cc22f] | 299 | ASSERT(!AtomicForce.empty(), | 
|---|
|  | 300 | "AtomInfo::operator[]() - AtomicForce is empty."); | 
|---|
|  | 301 | VectorTrajectory_t::const_iterator iter = | 
|---|
|  | 302 | AtomicForce.lower_bound(_step); | 
|---|
|  | 303 | // special, we only interpolate between present time steps not into the future | 
|---|
|  | 304 | if (_step > AtomicForce.begin()->first) | 
|---|
|  | 305 | return zeroVec; | 
|---|
|  | 306 | else | 
|---|
|  | 307 | return iter->second; | 
|---|
| [6625c3] | 308 | } | 
|---|
|  | 309 |  | 
|---|
| [bce72c] | 310 | void AtomInfo::setAtomicForce(const Vector &_newforce) | 
|---|
|  | 311 | { | 
|---|
| [8cc22f] | 312 | setAtomicForceAtStep(WorldTime::getTime(), _newforce); | 
|---|
| [bce72c] | 313 | } | 
|---|
|  | 314 |  | 
|---|
| [6b020f] | 315 | void AtomInfo::setAtomicForceAtStep(const unsigned int _step, const Vector &_newforce) | 
|---|
| [6625c3] | 316 | { | 
|---|
| [7188b1] | 317 | OBSERVE; | 
|---|
| [8cc22f] | 318 | VectorTrajectory_t::iterator iter = AtomicForce.find(_step); | 
|---|
|  | 319 | if (iter !=  AtomicForce.end()) { | 
|---|
|  | 320 | iter->second = _newforce; | 
|---|
|  | 321 | } else { | 
|---|
|  | 322 | #ifndef NDEBUG | 
|---|
|  | 323 | std::pair<VectorTrajectory_t::iterator, bool> inserter = | 
|---|
|  | 324 | #endif | 
|---|
|  | 325 | AtomicForce.insert( std::make_pair(_step, _newforce) ); | 
|---|
|  | 326 | ASSERT( inserter.second, | 
|---|
|  | 327 | "AtomInfo::set() - time step "+toString(_step) | 
|---|
|  | 328 | +" present after all?"); | 
|---|
|  | 329 | } | 
|---|
| [7188b1] | 330 | if (WorldTime::getTime() == _step) | 
|---|
| [bcb593] | 331 | NOTIFY(AtomObservable::ForceChanged); | 
|---|
| [6625c3] | 332 | } | 
|---|
|  | 333 |  | 
|---|
|  | 334 | bool AtomInfo::getFixedIon() const | 
|---|
|  | 335 | { | 
|---|
|  | 336 | return FixedIon; | 
|---|
|  | 337 | } | 
|---|
|  | 338 |  | 
|---|
|  | 339 | void AtomInfo::setFixedIon(const bool _fixedion) | 
|---|
|  | 340 | { | 
|---|
| [7188b1] | 341 | OBSERVE; | 
|---|
|  | 342 | NOTIFY(AtomObservable::PropertyChanged); | 
|---|
| [6625c3] | 343 | FixedIon = _fixedion; | 
|---|
|  | 344 | } | 
|---|
|  | 345 |  | 
|---|
| [d74077] | 346 | void AtomInfo::setPosition(const Vector& _vector) | 
|---|
|  | 347 | { | 
|---|
| [8cc22f] | 348 | setPositionAtStep(WorldTime::getTime(), _vector); | 
|---|
| [d74077] | 349 | } | 
|---|
|  | 350 |  | 
|---|
| [6b020f] | 351 | void AtomInfo::setPositionAtStep(unsigned int _step, const Vector& _vector) | 
|---|
| [6625c3] | 352 | { | 
|---|
| [7188b1] | 353 | OBSERVE; | 
|---|
| [8cc22f] | 354 | VectorTrajectory_t::iterator iter = AtomicPosition.find(_step); | 
|---|
|  | 355 | if (iter !=  AtomicPosition.end()) { | 
|---|
|  | 356 | iter->second = _vector; | 
|---|
|  | 357 | } else { | 
|---|
|  | 358 | #ifndef NDEBUG | 
|---|
|  | 359 | std::pair<VectorTrajectory_t::iterator, bool> inserter = | 
|---|
|  | 360 | #endif | 
|---|
|  | 361 | AtomicPosition.insert( std::make_pair(_step, _vector) ); | 
|---|
|  | 362 | ASSERT( inserter.second, | 
|---|
|  | 363 | "AtomInfo::set() - time step "+toString(_step) | 
|---|
|  | 364 | +" present after all?"); | 
|---|
|  | 365 | } | 
|---|
| [7188b1] | 366 | if (WorldTime::getTime() == _step) | 
|---|
|  | 367 | NOTIFY(AtomObservable::PositionChanged); | 
|---|
| [6625c3] | 368 | } | 
|---|
|  | 369 |  | 
|---|
| [d74077] | 370 | const VectorInterface& AtomInfo::operator+=(const Vector& b) | 
|---|
|  | 371 | { | 
|---|
| [8cc22f] | 372 | setPosition(getPosition()+b); | 
|---|
| [d74077] | 373 | return *this; | 
|---|
|  | 374 | } | 
|---|
|  | 375 |  | 
|---|
|  | 376 | const VectorInterface& AtomInfo::operator-=(const Vector& b) | 
|---|
|  | 377 | { | 
|---|
| [8cc22f] | 378 | setPosition(getPosition()-b); | 
|---|
| [d74077] | 379 | return *this; | 
|---|
|  | 380 | } | 
|---|
|  | 381 |  | 
|---|
|  | 382 | Vector const AtomInfo::operator+(const Vector& b) const | 
|---|
|  | 383 | { | 
|---|
| [8cc22f] | 384 | Vector a(getPosition()); | 
|---|
| [d74077] | 385 | a += b; | 
|---|
|  | 386 | return a; | 
|---|
|  | 387 | } | 
|---|
|  | 388 |  | 
|---|
|  | 389 | Vector const AtomInfo::operator-(const Vector& b) const | 
|---|
|  | 390 | { | 
|---|
| [8cc22f] | 391 | Vector a(getPosition()); | 
|---|
| [d74077] | 392 | a -= b; | 
|---|
|  | 393 | return a; | 
|---|
|  | 394 | } | 
|---|
|  | 395 |  | 
|---|
|  | 396 | double AtomInfo::distance(const Vector &point) const | 
|---|
|  | 397 | { | 
|---|
| [8cc22f] | 398 | return getPosition().distance(point); | 
|---|
| [d74077] | 399 | } | 
|---|
|  | 400 |  | 
|---|
|  | 401 | double AtomInfo::DistanceSquared(const Vector &y) const | 
|---|
|  | 402 | { | 
|---|
| [8cc22f] | 403 | return getPosition().DistanceSquared(y); | 
|---|
| [d74077] | 404 | } | 
|---|
|  | 405 |  | 
|---|
|  | 406 | double AtomInfo::distance(const VectorInterface &_atom) const | 
|---|
|  | 407 | { | 
|---|
| [8cc22f] | 408 | return _atom.distance(getPosition()); | 
|---|
| [d74077] | 409 | } | 
|---|
|  | 410 |  | 
|---|
|  | 411 | double AtomInfo::DistanceSquared(const VectorInterface &_atom) const | 
|---|
|  | 412 | { | 
|---|
| [8cc22f] | 413 | return _atom.DistanceSquared(getPosition()); | 
|---|
| [d74077] | 414 | } | 
|---|
|  | 415 |  | 
|---|
|  | 416 | VectorInterface &AtomInfo::operator=(const Vector& _vector) | 
|---|
|  | 417 | { | 
|---|
| [8cc22f] | 418 | setPosition(_vector); | 
|---|
| [d74077] | 419 | return *this; | 
|---|
|  | 420 | } | 
|---|
|  | 421 |  | 
|---|
|  | 422 | void AtomInfo::ScaleAll(const double *factor) | 
|---|
|  | 423 | { | 
|---|
| [8cc22f] | 424 | Vector temp(getPosition()); | 
|---|
|  | 425 | temp.ScaleAll(factor); | 
|---|
|  | 426 | setPosition(temp); | 
|---|
| [d74077] | 427 | } | 
|---|
|  | 428 |  | 
|---|
|  | 429 | void AtomInfo::ScaleAll(const Vector &factor) | 
|---|
|  | 430 | { | 
|---|
| [8cc22f] | 431 | Vector temp(getPosition()); | 
|---|
|  | 432 | temp.ScaleAll(factor); | 
|---|
|  | 433 | setPosition(temp); | 
|---|
| [d74077] | 434 | } | 
|---|
|  | 435 |  | 
|---|
|  | 436 | void AtomInfo::Scale(const double factor) | 
|---|
|  | 437 | { | 
|---|
| [8cc22f] | 438 | Vector temp(getPosition()); | 
|---|
|  | 439 | temp.Scale(factor); | 
|---|
|  | 440 | setPosition(temp); | 
|---|
| [d74077] | 441 | } | 
|---|
|  | 442 |  | 
|---|
|  | 443 | void AtomInfo::Zero() | 
|---|
|  | 444 | { | 
|---|
| [8cc22f] | 445 | setPosition(zeroVec); | 
|---|
| [d74077] | 446 | } | 
|---|
|  | 447 |  | 
|---|
|  | 448 | void AtomInfo::One(const double one) | 
|---|
|  | 449 | { | 
|---|
| [8cc22f] | 450 | setPosition(Vector(one,one,one)); | 
|---|
| [d74077] | 451 | } | 
|---|
|  | 452 |  | 
|---|
|  | 453 | void AtomInfo::LinearCombinationOfVectors(const Vector &x1, const Vector &x2, const Vector &x3, const double * const factors) | 
|---|
|  | 454 | { | 
|---|
| [8cc22f] | 455 | Vector newPos; | 
|---|
|  | 456 | newPos.LinearCombinationOfVectors(x1,x2,x3,factors); | 
|---|
|  | 457 | setPosition(newPos); | 
|---|
| [d74077] | 458 | } | 
|---|
|  | 459 |  | 
|---|
| [6625c3] | 460 | /** | 
|---|
|  | 461 | *  returns the kinetic energy of this atom at a given time step | 
|---|
|  | 462 | */ | 
|---|
| [7188b1] | 463 | double AtomInfo::getKineticEnergy(const unsigned int _step) const | 
|---|
|  | 464 | { | 
|---|
| [8cc22f] | 465 | return getMass() * getAtomicVelocityAtStep(_step).NormSquared(); | 
|---|
| [6625c3] | 466 | } | 
|---|
|  | 467 |  | 
|---|
| [7188b1] | 468 | Vector AtomInfo::getMomentum(const unsigned int _step) const | 
|---|
|  | 469 | { | 
|---|
| [8cc22f] | 470 | return getMass() * getAtomicVelocityAtStep(_step); | 
|---|
| [6625c3] | 471 | } | 
|---|
|  | 472 |  | 
|---|
| [8cc22f] | 473 | /** Decrease the trajectory if given \a MaxSteps is smaller. | 
|---|
|  | 474 | * Does nothing if \a MaxSteps is larger than current size. | 
|---|
|  | 475 | * | 
|---|
| [6625c3] | 476 | * \param MaxSteps | 
|---|
|  | 477 | */ | 
|---|
|  | 478 | void AtomInfo::ResizeTrajectory(size_t MaxSteps) | 
|---|
|  | 479 | { | 
|---|
| [8cc22f] | 480 | // mind the reverse ordering due to std::greater, latest time steps are at beginning | 
|---|
|  | 481 | VectorTrajectory_t::iterator positer = AtomicPosition.lower_bound(MaxSteps); | 
|---|
|  | 482 | if (positer != AtomicPosition.begin()) { | 
|---|
|  | 483 | if (positer->first == MaxSteps) | 
|---|
|  | 484 | --positer; | 
|---|
|  | 485 | AtomicPosition.erase(AtomicPosition.begin(), positer); | 
|---|
|  | 486 | } | 
|---|
|  | 487 | VectorTrajectory_t::iterator veliter = AtomicVelocity.lower_bound(MaxSteps); | 
|---|
|  | 488 | if (veliter != AtomicVelocity.begin()) { | 
|---|
|  | 489 | if (veliter->first == MaxSteps) | 
|---|
|  | 490 | --veliter; | 
|---|
|  | 491 | AtomicVelocity.erase(AtomicVelocity.begin(), veliter); | 
|---|
|  | 492 | } | 
|---|
|  | 493 | VectorTrajectory_t::iterator forceiter = AtomicForce.lower_bound(MaxSteps); | 
|---|
|  | 494 | if (forceiter != AtomicForce.begin()) { | 
|---|
|  | 495 | if (forceiter->first == MaxSteps) | 
|---|
|  | 496 | --forceiter; | 
|---|
|  | 497 | AtomicForce.erase(AtomicForce.begin(), forceiter); | 
|---|
|  | 498 | } | 
|---|
| [e2373df] | 499 | } | 
|---|
| [6625c3] | 500 |  | 
|---|
|  | 501 | size_t AtomInfo::getTrajectorySize() const | 
|---|
|  | 502 | { | 
|---|
| [8cc22f] | 503 | // mind greater comp for map here: first element is latest in time steps! | 
|---|
|  | 504 | return AtomicPosition.begin()->first+1; | 
|---|
| [6625c3] | 505 | } | 
|---|
|  | 506 |  | 
|---|
| [35a25a] | 507 | double AtomInfo::getMass() const | 
|---|
|  | 508 | { | 
|---|
|  | 509 | return getType()->getMass(); | 
|---|
| [6625c3] | 510 | } | 
|---|
|  | 511 |  | 
|---|
| [8cc22f] | 512 | /** Helper function to either insert or assign, depending on the element being | 
|---|
|  | 513 | * present already. | 
|---|
|  | 514 | * | 
|---|
|  | 515 | * \param _trajectory vector of Vectors to assign | 
|---|
|  | 516 | * \param dest step to insert/assign to | 
|---|
|  | 517 | * \param _newvalue new Vector value | 
|---|
|  | 518 | */ | 
|---|
|  | 519 | void assignTrajectoryElement( | 
|---|
|  | 520 | std::map<unsigned int, Vector, std::greater<unsigned int> > &_trajectory, | 
|---|
|  | 521 | const unsigned int dest, | 
|---|
|  | 522 | const Vector &_newvalue) | 
|---|
|  | 523 | { | 
|---|
|  | 524 | std::pair<std::map<unsigned int, Vector, std::greater<unsigned int> >::iterator, bool> inserter = | 
|---|
|  | 525 | _trajectory.insert( std::make_pair(dest, _newvalue) ); | 
|---|
|  | 526 | if (!inserter.second) | 
|---|
|  | 527 | inserter.first->second = _newvalue; | 
|---|
|  | 528 | } | 
|---|
|  | 529 |  | 
|---|
| [6625c3] | 530 | /** Copies a given trajectory step \a src onto another \a dest | 
|---|
|  | 531 | * \param dest index of destination step | 
|---|
|  | 532 | * \param src index of source step | 
|---|
|  | 533 | */ | 
|---|
| [6b020f] | 534 | void AtomInfo::CopyStepOnStep(const unsigned int dest, const unsigned int src) | 
|---|
| [6625c3] | 535 | { | 
|---|
|  | 536 | if (dest == src)  // self assignment check | 
|---|
|  | 537 | return; | 
|---|
|  | 538 |  | 
|---|
| [7188b1] | 539 | if (WorldTime::getTime() == dest){ | 
|---|
|  | 540 | NOTIFY(AtomObservable::PositionChanged); | 
|---|
|  | 541 | NOTIFY(AtomObservable::VelocityChanged); | 
|---|
|  | 542 | NOTIFY(AtomObservable::ForceChanged); | 
|---|
|  | 543 | } | 
|---|
|  | 544 |  | 
|---|
| [8cc22f] | 545 | VectorTrajectory_t::iterator positer = AtomicPosition.find(src); | 
|---|
|  | 546 | ASSERT( positer != AtomicPosition.end(), | 
|---|
|  | 547 | "AtomInfo::CopyStepOnStep() - step " | 
|---|
|  | 548 | +toString(src)+" to copy from not present in AtomicPosition."); | 
|---|
| [f946b2] | 549 | assignTrajectoryElement(AtomicPosition, dest, positer->second); | 
|---|
| [8cc22f] | 550 | VectorTrajectory_t::iterator veliter = AtomicVelocity.find(src); | 
|---|
| [f946b2] | 551 | if (veliter != AtomicVelocity.end()) | 
|---|
|  | 552 | assignTrajectoryElement(AtomicVelocity, dest, veliter->second); | 
|---|
| [8cc22f] | 553 | VectorTrajectory_t::iterator forceiter = AtomicForce.find(src); | 
|---|
| [f946b2] | 554 | if (forceiter != AtomicForce.end()) | 
|---|
|  | 555 | assignTrajectoryElement(AtomicForce, dest, forceiter->second); | 
|---|
| [6625c3] | 556 | }; | 
|---|
|  | 557 |  | 
|---|
| [bcb593] | 558 | /** Performs a velocity verlet update of the position at \a NextStep from \a LastStep information only. | 
|---|
|  | 559 | * | 
|---|
|  | 560 | * We calculate \f$x(t + \delta t) = x(t) + v(t)* \delta t + .5 * \delta t * \delta t * F(t)/m \f$. | 
|---|
|  | 561 | * | 
|---|
|  | 562 | * | 
|---|
| [6625c3] | 563 | * \param NextStep index of sequential step to set | 
|---|
| [435065] | 564 | * \param Deltat time step width | 
|---|
|  | 565 | * \param IsAngstroem whether the force's underlying unit of length is angstroem or bohr radii | 
|---|
| [6625c3] | 566 | */ | 
|---|
| [bcb593] | 567 | void AtomInfo::VelocityVerletUpdateX(int nr, const unsigned int NextStep, double Deltat, bool IsAngstroem) | 
|---|
| [6625c3] | 568 | { | 
|---|
| [4882d5] | 569 | const unsigned int LastStep = NextStep == 0 ? 0 : NextStep-1; | 
|---|
|  | 570 |  | 
|---|
| [435065] | 571 | LOG(2, "INFO: Particle that currently " << *this); | 
|---|
| [bcb593] | 572 | LOG(2, "INFO: Integrating position with mass=" << getMass() << " and Deltat=" | 
|---|
| [435065] | 573 | << Deltat << " at NextStep=" << NextStep); | 
|---|
| [056e70] | 574 |  | 
|---|
|  | 575 | // update position | 
|---|
| [4882d5] | 576 | { | 
|---|
|  | 577 | Vector tempVector = getPositionAtStep(LastStep); | 
|---|
|  | 578 | LOG(4, "INFO: initial position from last step " << setprecision(4) << tempVector); | 
|---|
|  | 579 | tempVector += Deltat*(getAtomicVelocityAtStep(LastStep));     // s(t) = s(0) + v * deltat + 1/2 a * deltat^2 | 
|---|
|  | 580 | LOG(4, "INFO: position with velocity " << getAtomicVelocityAtStep(LastStep) << " from last step " << tempVector); | 
|---|
|  | 581 | tempVector += .5*Deltat*Deltat*(getAtomicForceAtStep(LastStep))*(1./getMass());     // F = m * a and s = | 
|---|
| [bcb593] | 582 | LOG(4, "INFO: position with force " << getAtomicForceAtStep(LastStep) << " from last step " << tempVector); | 
|---|
| [4882d5] | 583 | setPositionAtStep(NextStep, tempVector); | 
|---|
|  | 584 | LOG(3, "INFO: Position at step " << NextStep << " set to " << tempVector); | 
|---|
|  | 585 | } | 
|---|
| [bcb593] | 586 | }; | 
|---|
|  | 587 |  | 
|---|
|  | 588 | /** Performs a velocity verlet update of the velocity at \a NextStep. | 
|---|
|  | 589 | * | 
|---|
|  | 590 | * \note forces at NextStep should have been calculated based on position at NextStep prior | 
|---|
|  | 591 | * to calling this function. | 
|---|
|  | 592 | * | 
|---|
|  | 593 | * We calculate \f$v(t) = v(t - \delta t) + \delta _t * .5 * (F(t - \delta t) + F(t))/m \f$. | 
|---|
|  | 594 | * | 
|---|
|  | 595 | * Parameters are according to those in configuration class. | 
|---|
|  | 596 | * \param NextStep index of sequential step to set | 
|---|
|  | 597 | * \param Deltat time step width | 
|---|
|  | 598 | * \param IsAngstroem whether the force's underlying unit of length is angstroem or bohr radii | 
|---|
|  | 599 | */ | 
|---|
|  | 600 | void AtomInfo::VelocityVerletUpdateU(int nr, const unsigned int NextStep, double Deltat, bool IsAngstroem) | 
|---|
|  | 601 | { | 
|---|
|  | 602 | const unsigned int LastStep = NextStep == 0 ? 0 : NextStep-1; | 
|---|
|  | 603 |  | 
|---|
|  | 604 | LOG(2, "INFO: Particle that currently " << *this); | 
|---|
|  | 605 | LOG(2, "INFO: Integrating velocity with mass=" << getMass() << " and Deltat=" | 
|---|
|  | 606 | << Deltat << " at NextStep=" << NextStep); | 
|---|
| [056e70] | 607 |  | 
|---|
| [6625c3] | 608 | // Update U | 
|---|
| [4882d5] | 609 | { | 
|---|
|  | 610 | Vector tempVector = getAtomicVelocityAtStep(LastStep); | 
|---|
|  | 611 | LOG(4, "INFO: initial velocity from last step " << tempVector); | 
|---|
|  | 612 | tempVector += Deltat * .5*(getAtomicForceAtStep(LastStep)+getAtomicForceAtStep(NextStep))*(1./getMass()); // v = F/m * t | 
|---|
|  | 613 | LOG(4, "INFO: Velocity with force from last " << getAtomicForceAtStep(LastStep) | 
|---|
|  | 614 | << " and present " << getAtomicForceAtStep(NextStep) << " step " << tempVector); | 
|---|
|  | 615 | setAtomicVelocityAtStep(NextStep, tempVector); | 
|---|
|  | 616 | LOG(3, "INFO: Velocity at step " << NextStep << " set to " << tempVector); | 
|---|
|  | 617 | } | 
|---|
| [6625c3] | 618 | }; | 
|---|
|  | 619 |  | 
|---|
| [d74077] | 620 | std::ostream & AtomInfo::operator << (std::ostream &ost) const | 
|---|
|  | 621 | { | 
|---|
|  | 622 | return (ost << getPosition()); | 
|---|
|  | 623 | } | 
|---|
|  | 624 |  | 
|---|
|  | 625 | std::ostream & operator << (std::ostream &ost, const AtomInfo &a) | 
|---|
|  | 626 | { | 
|---|
| [b1a5d9] | 627 | const size_t terminalstep = a.getTrajectorySize()-1; | 
|---|
| [e34254] | 628 | if (terminalstep) { | 
|---|
|  | 629 | ost << "starts at " | 
|---|
|  | 630 | << a.getPositionAtStep(0) << " and ends at " | 
|---|
|  | 631 | << a.getPositionAtStep(terminalstep) | 
|---|
|  | 632 | << " at time step " << terminalstep; | 
|---|
|  | 633 | } else { | 
|---|
|  | 634 | ost << "is at " | 
|---|
|  | 635 | << a.getPositionAtStep(0) << " with a single time step only"; | 
|---|
|  | 636 | } | 
|---|
| [d74077] | 637 | return ost; | 
|---|
|  | 638 | } | 
|---|
|  | 639 |  | 
|---|