45#include "optimization/PotentialEnergyObjectiveFunction.hpp"
53 PotentialEnergyObjectiveFunction::PotentialEnergyObjectiveFunction(
56 forceMan_(forceMan), thermo(info), hasFlucQ_(false) {
57 shake_ =
new Shake(info_);
59 if (info_->usesFluctuatingCharges()) {
60 if (info_->getNFluctuatingCharges() > 0) {
63 bool cr = info_->getSimParams()
64 ->getFluctuatingChargeParameters()
65 ->getConstrainRegions();
66 fqConstraints_->setConstrainRegions(cr);
74 shake_->constraintR();
75 forceMan_->calcForces();
76 if (hasFlucQ_) fqConstraints_->applyConstraints();
77 shake_->constraintF();
78 return thermo.getPotential();
84 shake_->constraintR();
85 forceMan_->calcForces();
86 if (hasFlucQ_) fqConstraints_->applyConstraints();
87 shake_->constraintF();
94 shake_->constraintR();
95 forceMan_->calcForces();
96 if (hasFlucQ_) fqConstraints_->applyConstraints();
97 shake_->constraintF();
99 return thermo.getPotential();
102 void PotentialEnergyObjectiveFunction::setCoor(
106 SimInfo::MoleculeIterator i;
107 Molecule::IntegrableObjectIterator j;
108 Molecule::AtomIterator ai;
117 index = displacements_[myrank_];
126 for (sd = mol->beginIntegrableObject(j); sd != NULL;
127 sd = mol->nextIntegrableObject(j)) {
128 position[0] = x[index++];
129 position[1] = x[index++];
130 position[2] = x[index++];
135 eulerAngle[0] = x[index++];
136 eulerAngle[1] = x[index++];
137 eulerAngle[2] = x[index++];
152 for (atom = mol->beginFluctuatingCharge(ai); atom != NULL;
153 atom = mol->nextFluctuatingCharge(ai)) {
160 void PotentialEnergyObjectiveFunction::getGrad(
162 SimInfo::MoleculeIterator i;
163 Molecule::IntegrableObjectIterator j;
164 Molecule::AtomIterator ai;
168 std::vector<RealType> myGrad;
172 index = displacements_[myrank_];
180 for (sd = mol->beginIntegrableObject(j); sd != NULL;
181 sd = mol->nextIntegrableObject(j)) {
184 for (
size_t k = 0; k < myGrad.size(); ++k) {
185 grad[index++] = myGrad[k];
193 for (atom = mol->beginFluctuatingCharge(ai); atom != NULL;
194 atom = mol->nextFluctuatingCharge(ai)) {
200 MPI_Allreduce(MPI_IN_PLACE, &grad[0], ndf_, MPI_REALTYPE, MPI_SUM,
207 MPI_Comm_size(MPI_COMM_WORLD, &nproc_);
208 MPI_Comm_rank(MPI_COMM_WORLD, &myrank_);
209 std::vector<int> onProc(nproc_, 0);
211 displacements_.clear();
212 displacements_.resize(nproc_, 0);
215 SimInfo::MoleculeIterator i;
216 Molecule::IntegrableObjectIterator j;
217 Molecule::AtomIterator ai;
229 for (sd = mol->beginIntegrableObject(j); sd != NULL;
230 sd = mol->nextIntegrableObject(j)) {
240 for (atom = mol->beginFluctuatingCharge(ai); atom != NULL;
241 atom = mol->nextFluctuatingCharge(ai)) {
248 MPI_Allgather(&ndf_, 1, MPI_INT, &onProc[0], 1, MPI_INT, MPI_COMM_WORLD);
251 for (
int iproc = 0; iproc < nproc_; iproc++) {
252 ndf_ += onProc[iproc];
255 displacements_[0] = 0;
256 for (
int iproc = 1; iproc < nproc_; iproc++) {
257 displacements_[iproc] = displacements_[iproc - 1] + onProc[iproc - 1];
265 index = displacements_[myrank_];
272 for (sd = mol->beginIntegrableObject(j); sd != NULL;
273 sd = mol->nextIntegrableObject(j)) {
275 xinit[index++] = pos[0];
276 xinit[index++] = pos[1];
277 xinit[index++] = pos[2];
281 xinit[index++] = eulerAngle[0];
282 xinit[index++] = eulerAngle[1];
283 xinit[index++] = eulerAngle[2];
291 for (atom = mol->beginFluctuatingCharge(ai); atom != NULL;
292 atom = mol->nextFluctuatingCharge(ai)) {
Dynamically-sized vector class.
void setZero()
zero out the vector
ForceManager is responsible for calculating both the short range (bonded) interactions and long range...
RealType value(const DynamicVector< RealType > &x)
method to overload to compute the objective function value in x
void gradient(DynamicVector< RealType > &grad, const DynamicVector< RealType > &x)
method to overload to compute grad_f, the first derivative of
RealType valueAndGradient(DynamicVector< RealType > &grad, const DynamicVector< RealType > &x)
method to overload to compute grad_f, the first derivative
void updateAtoms()
update the positions of atoms belong to this rigidbody
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.
"Don't move, or you're dead! Stand up! Captain, we've got them!"
Vector3d getEuler()
Returns the current euler angles of this stuntDouble.
void setFlucQPos(RealType charge)
Sets the current fluctuating charge of this stuntDouble.
Vector3d getPos()
Returns the current position of this stuntDouble.
RealType getFlucQFrc()
Returns the current charge force of this stuntDouble.
RealType getFlucQPos()
Returns the current fluctuating charge of this stuntDouble.
void setPos(const Vector3d &pos)
Sets the current position of this stuntDouble.
void setEuler(const Vector3d &euler)
Sets the current euler angles of this stuntDouble.
bool isRigidBody()
Tests if this stuntDouble is a rigid body.
bool isDirectional()
Tests if this stuntDouble is a directional one.
virtual std::vector< RealType > getGrad()=0
Returns the gradient of this stuntDouble.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.