45#include "parallel/ForceDecomposition.hpp"
56 ForceDecomposition::ForceDecomposition(
SimInfo* info,
59 interactionMan_(iMan), needVelocities_(false) {
60 sman_ = info_->getSnapshotManager();
61 atomStorageLayout_ = sman_->getAtomStorageLayout();
62 rigidBodyStorageLayout_ = sman_->getRigidBodyStorageLayout();
63 cutoffGroupStorageLayout_ = sman_->getCutoffGroupStorageLayout();
64 ff_ = info_->getForceField();
66 usePeriodicBoundaryConditions_ =
67 info->getSimParams()->getUsePeriodicBoundaryConditions();
69 Globals* simParams_ = info_->getSimParams();
70 if (simParams_->havePrintHeatFlux()) {
71 if (simParams_->getPrintHeatFlux()) { needVelocities_ =
true; }
74 if (simParams_->haveSkinThickness()) {
75 skinThickness_ = simParams_->getSkinThickness();
78 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
79 "ForceDecomposition: No value was set for the skinThickness.\n"
80 "\tOpenMD will use a default value of %f Angstroms\n"
81 "\tfor this simulation\n",
83 painCave.severity = OPENMD_INFO;
91 cellOffsets_.push_back(
Vector3i(0, 0, 0));
92 cellOffsets_.push_back(
Vector3i(1, 0, 0));
93 cellOffsets_.push_back(
Vector3i(1, 1, 0));
94 cellOffsets_.push_back(
Vector3i(0, 1, 0));
95 cellOffsets_.push_back(
Vector3i(-1, 1, 0));
96 cellOffsets_.push_back(
Vector3i(0, 0, 1));
97 cellOffsets_.push_back(
Vector3i(1, 0, 1));
98 cellOffsets_.push_back(
Vector3i(1, 1, 1));
99 cellOffsets_.push_back(
Vector3i(0, 1, 1));
100 cellOffsets_.push_back(
Vector3i(-1, 1, 1));
101 cellOffsets_.push_back(
Vector3i(-1, 0, 1));
102 cellOffsets_.push_back(
Vector3i(-1, -1, 1));
103 cellOffsets_.push_back(
Vector3i(0, -1, 1));
104 cellOffsets_.push_back(
Vector3i(1, -1, 1));
107 void ForceDecomposition::setCutoffRadius(RealType rcut) {
110 rListSq_ = rList_ * rList_;
113 void ForceDecomposition::fillPreForceData(
SelfData& sdat,
int atom) {
114 sdat.
atid = idents[atom];
116 sdat.
selePot = selectedSelfPot;
118 if (atomStorageLayout_ & DataStorage::dslDensity) {
122 if (atomStorageLayout_ & DataStorage::dslParticlePot) {
127 void ForceDecomposition::fillSelfData(
SelfData& sdat,
int atom) {
128 sdat.
atid = idents[atom];
130 sdat.
selePot = selectedSelfPot;
132 if (atomStorageLayout_ & DataStorage::dslSkippedCharge) {
136 if (atomStorageLayout_ & DataStorage::dslParticlePot) {
140 if (atomStorageLayout_ & DataStorage::dslFlucQPosition) {
144 if (atomStorageLayout_ & DataStorage::dslFlucQForce) {
149 void ForceDecomposition::unpackPreForceData(
SelfData& sdat,
int atom) {
151 selectedSelfPot = sdat.
selePot;
153 if (atomStorageLayout_ & DataStorage::dslFunctional) {
157 if (atomStorageLayout_ & DataStorage::dslFunctionalDerivative) {
162 (atomStorageLayout_ & DataStorage::dslParticlePot)) {
167 void ForceDecomposition::unpackSelfData(
SelfData& sdat,
int atom) {
169 selectedSelfPot = sdat.
selePot;
171 if (atomStorageLayout_ & DataStorage::dslFlucQForce) {
176 bool ForceDecomposition::checkNeighborList(vector<Vector3d> savedPositions) {
178 std::size_t nGroups = snap_->cgData.position.size();
181 DataStorage::dslVelocity);
186 if (savedPositions.size() != nGroups)
return true;
188 RealType dispmax = 0.0;
191 for (std::size_t i = 0; i < nGroups; i++) {
192 disp = snap_->cgData.position[i] - savedPositions[i];
197 MPI_Allreduce(MPI_IN_PLACE, &dispmax, 1, MPI_REALTYPE, MPI_MAX,
201 return (dispmax > st2) ? true :
false;
204 void ForceDecomposition::addToHeatFlux(
Vector3d hf) {
205 Vector3d chf = snap_->getConductiveHeatFlux();
207 snap_->setConductiveHeatFlux(chf);
209 void ForceDecomposition::setHeatFlux(
Vector3d hf) {
210 snap_->setConductiveHeatFlux(hf);
vector< RealType > functional
electron density
vector< RealType > flucQFrc
fluctuating charge velocities
vector< RealType > particlePot
torque array
vector< RealType > functionalDerivative
density functional
vector< RealType > flucQPos
charge skipped during normal pairwise calculation
vector< RealType > density
particle potential arrray
vector< RealType > skippedCharge
local electric field
void setStorageLayout(int layout)
Sets the storage layout
RealType skinThickness_
Verlet neighbor list skin thickness.
InteractionManager is responsible for keeping track of the non-bonded interactions (C++)
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
Real lengthSquare()
Returns the squared length of this vector.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
potVec selePot
potential energy of the selected site
RealType flucQfrc
fluctuating charge derivative
RealType frho
value of density functional for atom
RealType dfrhodrho
derivative of density functional for atom
potVec selfPot
total potential (including embedding energy)
RealType particlePot
contribution to potential from this particle
RealType rho
electron density
RealType flucQ
current value of atom's fluctuating charge
RealType skippedCharge
charge skipped in pairwise interaction loop
bool doParticlePot
should we bother with the particle pot?
int atid
atomType ident for the atom