51#include "brains/ForceModifier.hpp"
52#include "nonbonded/NonBondedInteraction.hpp"
54#include "types/FixedChargeAdapter.hpp"
55#include "types/FluctuatingChargeAdapter.hpp"
56#include "types/MultipoleAdapter.hpp"
57#include "utils/Constants.hpp"
61 UniformField::UniformField(
SimInfo* info) :
62 ForceModifier {info}, initialized {false}, doUniformField {false},
63 doParticlePot {false} {
64 simParams = info_->getSimParams();
67 void UniformField::initialize() {
68 std::vector<RealType> ef;
70 if (simParams->haveElectricField()) {
71 doUniformField =
true;
72 ef = simParams->getElectricField();
74 if (simParams->haveUniformField()) {
75 doUniformField =
true;
76 ef = simParams->getUniformField();
79 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
80 "UniformField: Incorrect number of parameters specified.\n"
81 "\tthere should be 3 parameters, but %lu were specified.\n",
91 if (storageLayout_ & DataStorage::dslParticlePot) doParticlePot =
true;
95 void UniformField::modifyForces() {
96 if (!initialized) initialize();
98 SimInfo::MoleculeIterator i;
99 Molecule::AtomIterator j;
103 potVec longRangePotential(0.0);
115 if (doUniformField) {
121 for (atom = mol->beginAtom(j); atom != NULL; atom = mol->nextAtom(j)) {
136 if (fca.isFixedCharge()) {
142 if (fqa.isFluctuatingCharge()) {
149 f = EF * C * Constants::chargeFieldConvert;
159 D = atom->
getDipole() * Constants::dipoleFieldConvert;
173 MPI_Allreduce(MPI_IN_PLACE, &fPot, 1, MPI_REALTYPE, MPI_SUM,
178 longRangePotential = snap->getLongRangePotentials();
180 snap->setLongRangePotentials(longRangePotential);
AtomType * getAtomType()
Returns the AtomType of this Atom.
AtomType is what OpenMD looks to for unchanging data about an atom.
Abstract class for external ForceModifier classes.
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.
SnapshotManager * getSnapshotManager()
Returns the snapshot manager.
The Snapshot class is a repository storing dynamic data during a Simulation.
void wrapVector(Vector3d &v)
Wrapping the vector according to periodic boundary condition.
Snapshot * getCurrentSnapshot()
Returns the pointer of current snapshot.
void addFlucQFrc(RealType cfrc)
Adds charge force into the current charge force of this stuntDouble.
void addElectricField(const Vector3d &eField)
Adds electric field into the current electric field of this stuntDouble.
Vector3d getPos()
Returns the current position of this stuntDouble.
RealType getFlucQPos()
Returns the current fluctuating charge of this stuntDouble.
Vector3d getDipole()
Returns the current dipole vector of this stuntDouble.
void addTrq(const Vector3d &trq)
Adds torque into the current torque of this stuntDouble.
void addParticlePot(const RealType &particlePot)
Adds particlePot into the current particlePot of this stuntDouble.
void addFrc(const Vector3d &frc)
Adds force into the current force of this stuntDouble.
Real & z()
Returns reference of the third element of Vector3.
Real & x()
Returns reference of the first element of Vector3.
Real & y()
Returns reference of the second element of Vector3.
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.
Real dot(const DynamicVector< Real > &v1, const DynamicVector< Real > &v2)
Returns the dot product of two DynamicVectors.
@ ELECTROSTATIC_FAMILY
Coulombic and point-multipole interactions.