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