45#include "FluctuatingChargeDamped.hpp"
48#include "utils/Constants.hpp"
49#include "utils/simError.h"
53 FluctuatingChargeDamped::FluctuatingChargeDamped(
SimInfo* info) :
55 snap(info->getSnapshotManager()->getCurrentSnapshot()) {}
57 void FluctuatingChargeDamped::initialize() {
58 FluctuatingChargePropagator::initialize();
60 if (info_->getSimParams()->haveDt()) {
61 dt_ = info_->getSimParams()->getDt();
64 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
65 "FluctuatingChargeDamped Error: dt is not set\n");
70 if (!fqParams_->haveDragCoefficient()) {
71 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
72 "If you use the FluctuatingChargeDamped\n"
73 "\tpropagator, you must set flucQ dragCoefficient .\n");
75 painCave.severity = OPENMD_ERROR;
79 drag_ = fqParams_->getDragCoefficient();
84 void FluctuatingChargeDamped::moveA() {
85 if (!hasFlucQ_)
return;
87 SimInfo::MoleculeIterator i;
88 Molecule::FluctuatingChargeIterator j;
91 RealType cvel, cpos, cfrc, cmass;
95 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
96 atom = mol->nextFluctuatingCharge(j)) {
100 cmass = atom->getChargeMass();
103 cvel += dt2_ * cfrc / cmass;
113 void FluctuatingChargeDamped::applyConstraints() {
114 if (!hasFlucQ_)
return;
116 SimInfo::MoleculeIterator i;
117 Molecule::FluctuatingChargeIterator j;
120 RealType cvel, cfrc, cmass, frictionForce;
121 RealType velStep, oldFF;
125 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
126 atom = mol->nextFluctuatingCharge(j)) {
140 cmass = atom->getChargeMass();
141 velStep = cvel + dt2_ * cfrc / cmass;
147 for (
int k = 0; k < maxIterNum_; k++) {
148 oldFF = frictionForce;
149 frictionForce = -drag_ * velStep;
152 velStep = cvel + dt2_ * (cfrc + frictionForce) / cmass;
156 if (fabs(frictionForce - oldFF) <= forceTolerance_)
162 fqConstraints_->applyConstraints();
165 void FluctuatingChargeDamped::moveB() {
166 if (!hasFlucQ_)
return;
167 SimInfo::MoleculeIterator i;
168 Molecule::FluctuatingChargeIterator j;
171 RealType cfrc, cvel, cmass;
175 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
176 atom = mol->nextFluctuatingCharge(j)) {
179 cmass = atom->getChargeMass();
182 cvel += (dt2_ * cfrc) / cmass;
189 void FluctuatingChargeDamped::updateSizes() {}
191 RealType FluctuatingChargeDamped::calcConservedQuantity() {
return 0.0; }
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 addFlucQFrc(RealType cfrc)
Adds charge force into the current charge force 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.