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"
60 UniformGradient::UniformGradient(
SimInfo* info) :
61 ForceModifier {info}, initialized {false}, doUniformGradient {false},
62 doParticlePot {false} {
63 simParams = info_->getSimParams();
66 void UniformGradient::initialize() {
71 if (simParams->haveUniformGradientDirection1()) {
72 std::vector<RealType> d1 = simParams->getUniformGradientDirection1();
74 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
75 "uniformGradientDirection1: Incorrect number of parameters\n"
76 "\tspecified. There should be 3 parameters, but %lu were\n"
90 if (simParams->haveUniformGradientDirection2()) {
91 std::vector<RealType> d2 = simParams->getUniformGradientDirection2();
93 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
94 "uniformGradientDirection2: Incorrect number of parameters\n"
95 "\tspecified. There should be 3 parameters, but %lu were\n"
109 if (simParams->haveUniformGradientStrength()) {
110 g_ = simParams->getUniformGradientStrength();
114 if (haveA && haveB && haveG) {
115 doUniformGradient =
true;
118 Grad_(0, 0) = 2.0 * (a_.
x() * b_.
x() - cpsi_ / 3.0);
119 Grad_(0, 1) = a_.
x() * b_.
y() + a_.
y() * b_.
x();
120 Grad_(0, 2) = a_.
x() * b_.
z() + a_.
z() * b_.
x();
121 Grad_(1, 0) = Grad_(0, 1);
122 Grad_(1, 1) = 2.0 * (a_.
y() * b_.
y() - cpsi_ / 3.0);
123 Grad_(1, 2) = a_.
y() * b_.
z() + a_.
z() * b_.
y();
124 Grad_(2, 0) = Grad_(0, 2);
125 Grad_(2, 1) = Grad_(1, 2);
126 Grad_(2, 2) = 2.0 * (a_.
z() * b_.
z() - cpsi_ / 3.0);
132 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
133 "UniformGradient: uniformGradientDirection1 not specified.\n");
134 painCave.isFatal = 1;
138 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
139 "UniformGradient: uniformGradientDirection2 not specified.\n");
140 painCave.isFatal = 1;
144 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
145 "UniformGradient: uniformGradientStrength not specified.\n");
146 painCave.isFatal = 1;
152 if (storageLayout_ & DataStorage::dslParticlePot) doParticlePot =
true;
156 void UniformGradient::modifyForces() {
157 if (!initialized) initialize();
159 SimInfo::MoleculeIterator i;
160 Molecule::AtomIterator j;
164 potVec longRangePotential(0.0);
180 if (doUniformGradient) {
186 for (atom = mol->beginAtom(j); atom != NULL; atom = mol->nextAtom(j)) {
200 if (fca.isFixedCharge()) {
206 if (fqa.isFluctuatingCharge()) {
213 f = EF * C * Constants::chargeFieldConvert;
223 D = atom->
getDipole() * Constants::dipoleFieldConvert;
236 if (ma.isQuadrupole()) {
239 t = 2.0 *
mCross(Q, Grad_);
250 MPI_Allreduce(MPI_IN_PLACE, &fPot, 1, MPI_REALTYPE, MPI_SUM,
255 longRangePotential = snap->getLongRangePotentials();
257 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.
Mat3x3d getQuadrupole()
Returns the current quadrupole tensor 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.
void normalize()
Normalizes this vector in place.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
Real doubleDot(const RectMatrix< Real, Row, Col > &t1, const RectMatrix< Real, Row, Col > &t2)
Returns the tensor contraction (double dot product) of two rank 2 tensors (or Matrices)
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.
Vector< Real, Row > mCross(const RectMatrix< Real, Row, Col > &t1, const RectMatrix< Real, Row, Col > &t2)
Returns the vector (cross) product of two matrices.
@ ELECTROSTATIC_FAMILY
Coulombic and point-multipole interactions.