45#include "brains/Velocitizer.hpp"
50#include "brains/Thermo.hpp"
51#include "flucq/FluctuatingChargeConstraints.hpp"
55#include "types/FluctuatingChargeAdapter.hpp"
56#include "utils/Constants.hpp"
57#include "utils/RandNumGen.hpp"
61 Velocitizer::Velocitizer(
SimInfo* info) :
62 info_(info), thermo_(info), globals_(info->getSimParams()),
63 randNumGen_(info->getRandomNumberGenerator()) {}
66 SimInfo::MoleculeIterator mi;
67 Molecule::IntegrableObjectIterator ioi;
74 for (sd = mol->beginIntegrableObject(ioi); sd != NULL;
75 sd = mol->nextIntegrableObject(ioi)) {
107 SimInfo::MoleculeIterator mi;
108 Molecule::IntegrableObjectIterator ioi;
112 std::normal_distribution<RealType> normalDistribution {0.0, 1.0};
114 kebar = Constants::kB * temperature * info_->
getNdfRaw() /
118 for (sd = mol->beginIntegrableObject(ioi); sd != NULL;
119 sd = mol->nextIntegrableObject(ioi)) {
122 av2 = 2.0 * kebar / sd->
getMass();
128 for (
int k = 0; k < 3; k++) {
129 v[k] = vbar * normalDistribution(*randNumGen_);
142 jbar = sqrt(2.0 * kebar * I(m, m));
143 j[m] = jbar * normalDistribution(*randNumGen_);
144 jbar = sqrt(2.0 * kebar * I(n, n));
145 j[n] = jbar * normalDistribution(*randNumGen_);
147 for (
int k = 0; k < 3; k++) {
148 jbar = sqrt(2.0 * kebar * I(k, k));
149 j[k] = jbar * normalDistribution(*randNumGen_);
171 SimInfo::MoleculeIterator mi;
172 Molecule::IntegrableObjectIterator ioi;
178 std::normal_distribution<RealType> normalDistribution {0.0, 1.0};
180 Globals* simParams = info_->getSimParams();
181 fqParams = simParams->getFluctuatingChargeParameters();
184 fqConstraints->setConstrainRegions(fqParams->getConstrainRegions());
188 ->getNumberOfFlucQConstraints();
189 int dfRaw = fqConstraints->getNumberOfFlucQAtoms();
190 int dfActual = dfRaw - nConstrain;
191 kebar = dfRaw * Constants::kb * temperature / (2 * dfActual);
195 for (sd = mol->beginIntegrableObject(ioi); sd != NULL;
196 sd = mol->nextIntegrableObject(ioi)) {
198 Atom* atom =
static_cast<Atom*
>(sd);
201 if (fqa.isFluctuatingCharge()) {
204 aw2 = 2.0 * kebar / atom->getChargeMass();
209 atom->
setFlucQVel(wbar * normalDistribution(*randNumGen_));
217 vector<Atom*> atomList;
220 for (
size_t i = 0; i < atomList.size(); ++i) {
221 Atom* atom = atomList[i];
224 if (fqa.isFluctuatingCharge()) {
226 aw2 = 2.0 * kebar / atom->getChargeMass();
230 atom->
setFlucQVel(wbar * normalDistribution(*randNumGen_));
236 fqConstraints->applyConstraintsOnChargeVelocities();
243 SimInfo::MoleculeIterator mi;
244 Molecule::IntegrableObjectIterator ioi;
252 for (sd = mol->beginIntegrableObject(ioi); sd != NULL;
253 sd = mol->nextIntegrableObject(ioi)) {
274 inertiaTensor = inertiaTensor.
inverse();
275 omega = inertiaTensor * angularMomentum;
277 SimInfo::MoleculeIterator mi;
278 Molecule::IntegrableObjectIterator ioi;
288 for (sd = mol->beginIntegrableObject(ioi); sd != NULL;
289 sd = mol->nextIntegrableObject(ioi)) {
290 tempComPos = sd->
getPos() - com;
AtomType * getAtomType()
Returns the AtomType of this Atom.
AtomType is what OpenMD looks to for unchanging data about an atom.
std::vector< Atom * > getAtoms()
Returns the atoms of this rigid body.
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
int getNdf()
Returns the number of degrees of freedom.
int getNdfRaw()
Returns the number of raw degrees of freedom.
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.
SquareMatrix3< Real > inverse() const
Sets the value of this matrix to the inverse of itself.
"Don't move, or you're dead! Stand up! Captain, we've got them!"
void setFlucQVel(RealType cvel)
Sets the current charge velocity of this stuntDouble.
Vector3d getVel()
Returns the current velocity of this stuntDouble.
int linearAxis()
Returns the linear axis of the rigidbody, atom and directional atom will always return -1.
RealType getMass()
Returns the mass of this stuntDouble.
virtual Mat3x3d getI()=0
Returns the inertia tensor of this stuntDouble.
bool isLinear()
Tests the if this stuntDouble is a linear rigidbody.
Vector3d getPos()
Returns the current position of this stuntDouble.
void setVel(const Vector3d &vel)
Sets the current velocity of this stuntDouble.
bool isRigidBody()
Tests if this stuntDouble is a rigid body.
Vector3d getJ()
Returns the current angular momentum of this stuntDouble (body -fixed).
bool isAtom()
Tests if this stuntDouble is an atom.
bool isDirectional()
Tests if this stuntDouble is a directional one.
void setJ(const Vector3d &angMom)
Sets the current angular momentum of this stuntDouble (body-fixed).
void getComAll(Vector3d &com, Vector3d &comVel)
Returns the center of the mass and Center of Mass velocity of the whole system.
void getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum)
Returns the inertia tensor and the total angular momentum for for the entire system.
Vector3d getComVel()
Returns the velocity of center of mass of the whole system.
void randomize(RealType ct)
Resamples velocities and angular momenta Resamples velocities and angular momenta from a Maxwell-Bolt...
void randomizeChargeVelocity(RealType ct)
Resamples charge velocities Resamples charge velocities from a Maxwell-Boltzmann distribution.
void removeComDrift()
Removes Center of Mass Drift Velocity Removes the center of mass drift velocity (required for accurat...
void scale(RealType lambda)
Scales velocities and angular momenta by a scaling factor Rescales velocity (and angular momenta) by ...
void removeAngularDrift()
Removes Center of Mass Angular momentum Removes the center of mass angular momentum (particularly use...
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
Vector3< Real > cross(const Vector3< Real > &v1, const Vector3< Real > &v2)
Returns the cross product of two Vectors.