45#include "FluctuatingChargeNVE.hpp"
48#include "utils/simError.h"
52 FluctuatingChargeNVE::FluctuatingChargeNVE(
SimInfo* info) :
54 snap(info->getSnapshotManager()->getCurrentSnapshot()), thermo(info) {}
56 void FluctuatingChargeNVE::initialize() {
57 FluctuatingChargePropagator::initialize();
59 if (info_->getSimParams()->haveDt()) {
60 dt_ = info_->getSimParams()->getDt();
63 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
64 "FluctuatingChargeNVE Error: dt is not set\n");
69 if (!info_->getSimParams()->getUseIntialExtendedSystemState()) {
70 snap->setElectronicThermostat(make_pair(0.0, 0.0));
75 void FluctuatingChargeNVE::moveA() {
76 if (!hasFlucQ_)
return;
78 SimInfo::MoleculeIterator i;
79 Molecule::FluctuatingChargeIterator j;
82 RealType cvel, cpos, cfrc, cmass;
86 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
87 atom = mol->nextFluctuatingCharge(j)) {
91 cmass = atom->getChargeMass();
94 cvel += dt2_ * cfrc / cmass;
104 void FluctuatingChargeNVE::moveB() {
105 if (!hasFlucQ_)
return;
106 SimInfo::MoleculeIterator i;
107 Molecule::FluctuatingChargeIterator j;
110 RealType cfrc, cvel, cmass;
114 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
115 atom = mol->nextFluctuatingCharge(j)) {
118 cmass = atom->getChargeMass();
121 cvel += dt2_ * cfrc / cmass;
126 void FluctuatingChargeNVE::VelocityStep(RealType dt) {
127 if (!hasFlucQ_)
return;
129 SimInfo::MoleculeIterator i;
130 Molecule::FluctuatingChargeIterator j;
133 RealType cvel, cfrc, cmass;
137 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
138 atom = mol->nextFluctuatingCharge(j)) {
141 cmass = atom->getChargeMass();
144 cvel += dt * cfrc / cmass;
150 void FluctuatingChargeNVE::PositionStep(RealType dt) {
151 if (!hasFlucQ_)
return;
153 SimInfo::MoleculeIterator i;
154 Molecule::FluctuatingChargeIterator j;
161 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
162 atom = mol->nextFluctuatingCharge(j)) {
abstract class for propagating fluctuating charge variables
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
Molecule * beginMolecule(MoleculeIterator &i)
Returns the first molecule in this SimInfo and intialize the iterator.
Molecule * nextMolecule(MoleculeIterator &i)
Returns the next avaliable Molecule based on the iterator.
void setFlucQVel(RealType cvel)
Sets the current charge velocity of this stuntDouble.
void setFlucQPos(RealType charge)
Sets the current fluctuating charge of this stuntDouble.
RealType getFlucQFrc()
Returns the current charge force of this stuntDouble.
RealType getFlucQPos()
Returns the current fluctuating charge of this stuntDouble.
RealType getFlucQVel()
Returns the current charge velocity of this stuntDouble.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.