1 |
gezelter |
507 |
/* |
2 |
gezelter |
246 |
* Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
|
|
* |
4 |
|
|
* The University of Notre Dame grants you ("Licensee") a |
5 |
|
|
* non-exclusive, royalty free, license to use, modify and |
6 |
|
|
* redistribute this software in source and binary code form, provided |
7 |
|
|
* that the following conditions are met: |
8 |
|
|
* |
9 |
gezelter |
1390 |
* 1. Redistributions of source code must retain the above copyright |
10 |
gezelter |
246 |
* notice, this list of conditions and the following disclaimer. |
11 |
|
|
* |
12 |
gezelter |
1390 |
* 2. Redistributions in binary form must reproduce the above copyright |
13 |
gezelter |
246 |
* notice, this list of conditions and the following disclaimer in the |
14 |
|
|
* documentation and/or other materials provided with the |
15 |
|
|
* distribution. |
16 |
|
|
* |
17 |
|
|
* This software is provided "AS IS," without a warranty of any |
18 |
|
|
* kind. All express or implied conditions, representations and |
19 |
|
|
* warranties, including any implied warranty of merchantability, |
20 |
|
|
* fitness for a particular purpose or non-infringement, are hereby |
21 |
|
|
* excluded. The University of Notre Dame and its licensors shall not |
22 |
|
|
* be liable for any damages suffered by licensee as a result of |
23 |
|
|
* using, modifying or distributing the software or its |
24 |
|
|
* derivatives. In no event will the University of Notre Dame or its |
25 |
|
|
* licensors be liable for any lost revenue, profit or data, or for |
26 |
|
|
* direct, indirect, special, consequential, incidental or punitive |
27 |
|
|
* damages, however caused and regardless of the theory of liability, |
28 |
|
|
* arising out of the use of or inability to use software, even if the |
29 |
|
|
* University of Notre Dame has been advised of the possibility of |
30 |
|
|
* such damages. |
31 |
gezelter |
1390 |
* |
32 |
|
|
* SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your |
33 |
|
|
* research, please cite the appropriate papers when you publish your |
34 |
|
|
* work. Good starting points are: |
35 |
|
|
* |
36 |
|
|
* [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). |
37 |
|
|
* [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). |
38 |
|
|
* [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008). |
39 |
|
|
* [4] Vardeman & Gezelter, in progress (2009). |
40 |
gezelter |
246 |
*/ |
41 |
|
|
|
42 |
|
|
/** |
43 |
|
|
* @file SimInfo.cpp |
44 |
|
|
* @author tlin |
45 |
|
|
* @date 11/02/2004 |
46 |
|
|
* @version 1.0 |
47 |
|
|
*/ |
48 |
gezelter |
2 |
|
49 |
gezelter |
246 |
#include <algorithm> |
50 |
|
|
#include <set> |
51 |
tim |
749 |
#include <map> |
52 |
gezelter |
2 |
|
53 |
tim |
3 |
#include "brains/SimInfo.hpp" |
54 |
gezelter |
246 |
#include "math/Vector3.hpp" |
55 |
|
|
#include "primitives/Molecule.hpp" |
56 |
tim |
1024 |
#include "primitives/StuntDouble.hpp" |
57 |
gezelter |
586 |
#include "UseTheForce/fCutoffPolicy.h" |
58 |
chrisfen |
726 |
#include "UseTheForce/DarkSide/fSwitchingFunctionType.h" |
59 |
gezelter |
246 |
#include "UseTheForce/doForces_interface.h" |
60 |
chuckv |
1095 |
#include "UseTheForce/DarkSide/neighborLists_interface.h" |
61 |
chrisfen |
726 |
#include "UseTheForce/DarkSide/switcheroo_interface.h" |
62 |
gezelter |
246 |
#include "utils/MemoryUtils.hpp" |
63 |
tim |
3 |
#include "utils/simError.h" |
64 |
tim |
316 |
#include "selection/SelectionManager.hpp" |
65 |
chuckv |
834 |
#include "io/ForceFieldOptions.hpp" |
66 |
|
|
#include "UseTheForce/ForceField.hpp" |
67 |
gezelter |
2 |
|
68 |
chuckv |
1095 |
|
69 |
gezelter |
246 |
#ifdef IS_MPI |
70 |
|
|
#include "UseTheForce/mpiComponentPlan.h" |
71 |
|
|
#include "UseTheForce/DarkSide/simParallel_interface.h" |
72 |
|
|
#endif |
73 |
gezelter |
2 |
|
74 |
gezelter |
1390 |
namespace OpenMD { |
75 |
tim |
749 |
std::set<int> getRigidSet(int index, std::map<int, std::set<int> >& container) { |
76 |
|
|
std::map<int, std::set<int> >::iterator i = container.find(index); |
77 |
|
|
std::set<int> result; |
78 |
|
|
if (i != container.end()) { |
79 |
|
|
result = i->second; |
80 |
|
|
} |
81 |
gezelter |
2 |
|
82 |
tim |
749 |
return result; |
83 |
|
|
} |
84 |
|
|
|
85 |
tim |
770 |
SimInfo::SimInfo(ForceField* ff, Globals* simParams) : |
86 |
|
|
forceField_(ff), simParams_(simParams), |
87 |
gezelter |
945 |
ndf_(0), fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0), |
88 |
gezelter |
507 |
nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0), |
89 |
|
|
nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0), |
90 |
gezelter |
1277 |
nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nInversions_(0), |
91 |
|
|
nRigidBodies_(0), nIntegrableObjects_(0), nCutoffGroups_(0), |
92 |
|
|
nConstraints_(0), sman_(NULL), fortranInitialized_(false), |
93 |
|
|
calcBoxDipole_(false), useAtomicVirial_(true) { |
94 |
gezelter |
2 |
|
95 |
gezelter |
1277 |
|
96 |
gezelter |
507 |
MoleculeStamp* molStamp; |
97 |
|
|
int nMolWithSameStamp; |
98 |
|
|
int nCutoffAtoms = 0; // number of atoms belong to cutoff groups |
99 |
chrisfen |
645 |
int nGroups = 0; //total cutoff groups defined in meta-data file |
100 |
gezelter |
507 |
CutoffGroupStamp* cgStamp; |
101 |
|
|
RigidBodyStamp* rbStamp; |
102 |
|
|
int nRigidAtoms = 0; |
103 |
gezelter |
1277 |
|
104 |
tim |
770 |
std::vector<Component*> components = simParams->getComponents(); |
105 |
|
|
|
106 |
|
|
for (std::vector<Component*>::iterator i = components.begin(); i !=components.end(); ++i) { |
107 |
|
|
molStamp = (*i)->getMoleculeStamp(); |
108 |
|
|
nMolWithSameStamp = (*i)->getNMol(); |
109 |
gezelter |
246 |
|
110 |
|
|
addMoleculeStamp(molStamp, nMolWithSameStamp); |
111 |
gezelter |
2 |
|
112 |
gezelter |
246 |
//calculate atoms in molecules |
113 |
|
|
nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp; |
114 |
gezelter |
2 |
|
115 |
gezelter |
246 |
//calculate atoms in cutoff groups |
116 |
|
|
int nAtomsInGroups = 0; |
117 |
|
|
int nCutoffGroupsInStamp = molStamp->getNCutoffGroups(); |
118 |
|
|
|
119 |
|
|
for (int j=0; j < nCutoffGroupsInStamp; j++) { |
120 |
tim |
770 |
cgStamp = molStamp->getCutoffGroupStamp(j); |
121 |
gezelter |
507 |
nAtomsInGroups += cgStamp->getNMembers(); |
122 |
gezelter |
246 |
} |
123 |
gezelter |
2 |
|
124 |
gezelter |
246 |
nGroups += nCutoffGroupsInStamp * nMolWithSameStamp; |
125 |
chrisfen |
645 |
|
126 |
gezelter |
246 |
nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp; |
127 |
gezelter |
2 |
|
128 |
gezelter |
246 |
//calculate atoms in rigid bodies |
129 |
|
|
int nAtomsInRigidBodies = 0; |
130 |
tim |
274 |
int nRigidBodiesInStamp = molStamp->getNRigidBodies(); |
131 |
gezelter |
246 |
|
132 |
|
|
for (int j=0; j < nRigidBodiesInStamp; j++) { |
133 |
tim |
770 |
rbStamp = molStamp->getRigidBodyStamp(j); |
134 |
gezelter |
507 |
nAtomsInRigidBodies += rbStamp->getNMembers(); |
135 |
gezelter |
246 |
} |
136 |
gezelter |
2 |
|
137 |
gezelter |
246 |
nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp; |
138 |
|
|
nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp; |
139 |
|
|
|
140 |
gezelter |
507 |
} |
141 |
chrisfen |
143 |
|
142 |
chrisfen |
645 |
//every free atom (atom does not belong to cutoff groups) is a cutoff |
143 |
|
|
//group therefore the total number of cutoff groups in the system is |
144 |
|
|
//equal to the total number of atoms minus number of atoms belong to |
145 |
|
|
//cutoff group defined in meta-data file plus the number of cutoff |
146 |
|
|
//groups defined in meta-data file |
147 |
gezelter |
507 |
nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups; |
148 |
gezelter |
2 |
|
149 |
chrisfen |
645 |
//every free atom (atom does not belong to rigid bodies) is an |
150 |
|
|
//integrable object therefore the total number of integrable objects |
151 |
|
|
//in the system is equal to the total number of atoms minus number of |
152 |
|
|
//atoms belong to rigid body defined in meta-data file plus the number |
153 |
|
|
//of rigid bodies defined in meta-data file |
154 |
|
|
nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms |
155 |
|
|
+ nGlobalRigidBodies_; |
156 |
|
|
|
157 |
gezelter |
507 |
nGlobalMols_ = molStampIds_.size(); |
158 |
|
|
molToProcMap_.resize(nGlobalMols_); |
159 |
|
|
} |
160 |
gezelter |
2 |
|
161 |
gezelter |
507 |
SimInfo::~SimInfo() { |
162 |
tim |
398 |
std::map<int, Molecule*>::iterator i; |
163 |
|
|
for (i = molecules_.begin(); i != molecules_.end(); ++i) { |
164 |
gezelter |
507 |
delete i->second; |
165 |
tim |
398 |
} |
166 |
|
|
molecules_.clear(); |
167 |
tim |
490 |
|
168 |
gezelter |
246 |
delete sman_; |
169 |
|
|
delete simParams_; |
170 |
|
|
delete forceField_; |
171 |
gezelter |
507 |
} |
172 |
gezelter |
2 |
|
173 |
gezelter |
507 |
int SimInfo::getNGlobalConstraints() { |
174 |
gezelter |
246 |
int nGlobalConstraints; |
175 |
|
|
#ifdef IS_MPI |
176 |
|
|
MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM, |
177 |
|
|
MPI_COMM_WORLD); |
178 |
|
|
#else |
179 |
|
|
nGlobalConstraints = nConstraints_; |
180 |
|
|
#endif |
181 |
|
|
return nGlobalConstraints; |
182 |
gezelter |
507 |
} |
183 |
gezelter |
2 |
|
184 |
gezelter |
507 |
bool SimInfo::addMolecule(Molecule* mol) { |
185 |
gezelter |
246 |
MoleculeIterator i; |
186 |
gezelter |
2 |
|
187 |
gezelter |
246 |
i = molecules_.find(mol->getGlobalIndex()); |
188 |
|
|
if (i == molecules_.end() ) { |
189 |
gezelter |
2 |
|
190 |
gezelter |
507 |
molecules_.insert(std::make_pair(mol->getGlobalIndex(), mol)); |
191 |
gezelter |
246 |
|
192 |
gezelter |
507 |
nAtoms_ += mol->getNAtoms(); |
193 |
|
|
nBonds_ += mol->getNBonds(); |
194 |
|
|
nBends_ += mol->getNBends(); |
195 |
|
|
nTorsions_ += mol->getNTorsions(); |
196 |
gezelter |
1277 |
nInversions_ += mol->getNInversions(); |
197 |
gezelter |
507 |
nRigidBodies_ += mol->getNRigidBodies(); |
198 |
|
|
nIntegrableObjects_ += mol->getNIntegrableObjects(); |
199 |
|
|
nCutoffGroups_ += mol->getNCutoffGroups(); |
200 |
|
|
nConstraints_ += mol->getNConstraintPairs(); |
201 |
gezelter |
2 |
|
202 |
gezelter |
1287 |
addInteractionPairs(mol); |
203 |
|
|
|
204 |
gezelter |
507 |
return true; |
205 |
gezelter |
246 |
} else { |
206 |
gezelter |
507 |
return false; |
207 |
gezelter |
246 |
} |
208 |
gezelter |
507 |
} |
209 |
gezelter |
2 |
|
210 |
gezelter |
507 |
bool SimInfo::removeMolecule(Molecule* mol) { |
211 |
gezelter |
246 |
MoleculeIterator i; |
212 |
|
|
i = molecules_.find(mol->getGlobalIndex()); |
213 |
gezelter |
2 |
|
214 |
gezelter |
246 |
if (i != molecules_.end() ) { |
215 |
gezelter |
2 |
|
216 |
gezelter |
507 |
assert(mol == i->second); |
217 |
gezelter |
246 |
|
218 |
gezelter |
507 |
nAtoms_ -= mol->getNAtoms(); |
219 |
|
|
nBonds_ -= mol->getNBonds(); |
220 |
|
|
nBends_ -= mol->getNBends(); |
221 |
|
|
nTorsions_ -= mol->getNTorsions(); |
222 |
gezelter |
1277 |
nInversions_ -= mol->getNInversions(); |
223 |
gezelter |
507 |
nRigidBodies_ -= mol->getNRigidBodies(); |
224 |
|
|
nIntegrableObjects_ -= mol->getNIntegrableObjects(); |
225 |
|
|
nCutoffGroups_ -= mol->getNCutoffGroups(); |
226 |
|
|
nConstraints_ -= mol->getNConstraintPairs(); |
227 |
gezelter |
2 |
|
228 |
gezelter |
1287 |
removeInteractionPairs(mol); |
229 |
gezelter |
507 |
molecules_.erase(mol->getGlobalIndex()); |
230 |
gezelter |
2 |
|
231 |
gezelter |
507 |
delete mol; |
232 |
gezelter |
246 |
|
233 |
gezelter |
507 |
return true; |
234 |
gezelter |
246 |
} else { |
235 |
gezelter |
507 |
return false; |
236 |
gezelter |
246 |
} |
237 |
|
|
|
238 |
|
|
|
239 |
gezelter |
507 |
} |
240 |
gezelter |
246 |
|
241 |
|
|
|
242 |
gezelter |
507 |
Molecule* SimInfo::beginMolecule(MoleculeIterator& i) { |
243 |
gezelter |
246 |
i = molecules_.begin(); |
244 |
|
|
return i == molecules_.end() ? NULL : i->second; |
245 |
gezelter |
507 |
} |
246 |
gezelter |
246 |
|
247 |
gezelter |
507 |
Molecule* SimInfo::nextMolecule(MoleculeIterator& i) { |
248 |
gezelter |
246 |
++i; |
249 |
|
|
return i == molecules_.end() ? NULL : i->second; |
250 |
gezelter |
507 |
} |
251 |
gezelter |
2 |
|
252 |
|
|
|
253 |
gezelter |
507 |
void SimInfo::calcNdf() { |
254 |
gezelter |
246 |
int ndf_local; |
255 |
|
|
MoleculeIterator i; |
256 |
|
|
std::vector<StuntDouble*>::iterator j; |
257 |
|
|
Molecule* mol; |
258 |
|
|
StuntDouble* integrableObject; |
259 |
gezelter |
2 |
|
260 |
gezelter |
246 |
ndf_local = 0; |
261 |
|
|
|
262 |
|
|
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
263 |
gezelter |
507 |
for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL; |
264 |
|
|
integrableObject = mol->nextIntegrableObject(j)) { |
265 |
gezelter |
2 |
|
266 |
gezelter |
507 |
ndf_local += 3; |
267 |
gezelter |
2 |
|
268 |
gezelter |
507 |
if (integrableObject->isDirectional()) { |
269 |
|
|
if (integrableObject->isLinear()) { |
270 |
|
|
ndf_local += 2; |
271 |
|
|
} else { |
272 |
|
|
ndf_local += 3; |
273 |
|
|
} |
274 |
|
|
} |
275 |
gezelter |
246 |
|
276 |
tim |
770 |
} |
277 |
|
|
} |
278 |
gezelter |
246 |
|
279 |
|
|
// n_constraints is local, so subtract them on each processor |
280 |
|
|
ndf_local -= nConstraints_; |
281 |
|
|
|
282 |
|
|
#ifdef IS_MPI |
283 |
|
|
MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD); |
284 |
|
|
#else |
285 |
|
|
ndf_ = ndf_local; |
286 |
|
|
#endif |
287 |
|
|
|
288 |
|
|
// nZconstraints_ is global, as are the 3 COM translations for the |
289 |
|
|
// entire system: |
290 |
|
|
ndf_ = ndf_ - 3 - nZconstraint_; |
291 |
|
|
|
292 |
gezelter |
507 |
} |
293 |
gezelter |
2 |
|
294 |
gezelter |
945 |
int SimInfo::getFdf() { |
295 |
|
|
#ifdef IS_MPI |
296 |
|
|
MPI_Allreduce(&fdf_local,&fdf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD); |
297 |
|
|
#else |
298 |
|
|
fdf_ = fdf_local; |
299 |
|
|
#endif |
300 |
|
|
return fdf_; |
301 |
|
|
} |
302 |
|
|
|
303 |
gezelter |
507 |
void SimInfo::calcNdfRaw() { |
304 |
gezelter |
246 |
int ndfRaw_local; |
305 |
gezelter |
2 |
|
306 |
gezelter |
246 |
MoleculeIterator i; |
307 |
|
|
std::vector<StuntDouble*>::iterator j; |
308 |
|
|
Molecule* mol; |
309 |
|
|
StuntDouble* integrableObject; |
310 |
|
|
|
311 |
|
|
// Raw degrees of freedom that we have to set |
312 |
|
|
ndfRaw_local = 0; |
313 |
|
|
|
314 |
|
|
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
315 |
gezelter |
507 |
for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL; |
316 |
|
|
integrableObject = mol->nextIntegrableObject(j)) { |
317 |
gezelter |
246 |
|
318 |
gezelter |
507 |
ndfRaw_local += 3; |
319 |
gezelter |
246 |
|
320 |
gezelter |
507 |
if (integrableObject->isDirectional()) { |
321 |
|
|
if (integrableObject->isLinear()) { |
322 |
|
|
ndfRaw_local += 2; |
323 |
|
|
} else { |
324 |
|
|
ndfRaw_local += 3; |
325 |
|
|
} |
326 |
|
|
} |
327 |
gezelter |
246 |
|
328 |
gezelter |
507 |
} |
329 |
gezelter |
246 |
} |
330 |
|
|
|
331 |
|
|
#ifdef IS_MPI |
332 |
|
|
MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD); |
333 |
|
|
#else |
334 |
|
|
ndfRaw_ = ndfRaw_local; |
335 |
|
|
#endif |
336 |
gezelter |
507 |
} |
337 |
gezelter |
2 |
|
338 |
gezelter |
507 |
void SimInfo::calcNdfTrans() { |
339 |
gezelter |
246 |
int ndfTrans_local; |
340 |
gezelter |
2 |
|
341 |
gezelter |
246 |
ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_; |
342 |
gezelter |
2 |
|
343 |
|
|
|
344 |
gezelter |
246 |
#ifdef IS_MPI |
345 |
|
|
MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD); |
346 |
|
|
#else |
347 |
|
|
ndfTrans_ = ndfTrans_local; |
348 |
|
|
#endif |
349 |
gezelter |
2 |
|
350 |
gezelter |
246 |
ndfTrans_ = ndfTrans_ - 3 - nZconstraint_; |
351 |
|
|
|
352 |
gezelter |
507 |
} |
353 |
gezelter |
2 |
|
354 |
gezelter |
1287 |
void SimInfo::addInteractionPairs(Molecule* mol) { |
355 |
|
|
ForceFieldOptions& options_ = forceField_->getForceFieldOptions(); |
356 |
gezelter |
246 |
std::vector<Bond*>::iterator bondIter; |
357 |
|
|
std::vector<Bend*>::iterator bendIter; |
358 |
|
|
std::vector<Torsion*>::iterator torsionIter; |
359 |
gezelter |
1277 |
std::vector<Inversion*>::iterator inversionIter; |
360 |
gezelter |
246 |
Bond* bond; |
361 |
|
|
Bend* bend; |
362 |
|
|
Torsion* torsion; |
363 |
gezelter |
1277 |
Inversion* inversion; |
364 |
gezelter |
246 |
int a; |
365 |
|
|
int b; |
366 |
|
|
int c; |
367 |
|
|
int d; |
368 |
tim |
749 |
|
369 |
gezelter |
1287 |
// atomGroups can be used to add special interaction maps between |
370 |
|
|
// groups of atoms that are in two separate rigid bodies. |
371 |
|
|
// However, most site-site interactions between two rigid bodies |
372 |
|
|
// are probably not special, just the ones between the physically |
373 |
|
|
// bonded atoms. Interactions *within* a single rigid body should |
374 |
|
|
// always be excluded. These are done at the bottom of this |
375 |
|
|
// function. |
376 |
|
|
|
377 |
tim |
749 |
std::map<int, std::set<int> > atomGroups; |
378 |
|
|
Molecule::RigidBodyIterator rbIter; |
379 |
|
|
RigidBody* rb; |
380 |
|
|
Molecule::IntegrableObjectIterator ii; |
381 |
|
|
StuntDouble* integrableObject; |
382 |
gezelter |
246 |
|
383 |
gezelter |
1287 |
for (integrableObject = mol->beginIntegrableObject(ii); |
384 |
|
|
integrableObject != NULL; |
385 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
386 |
|
|
|
387 |
tim |
749 |
if (integrableObject->isRigidBody()) { |
388 |
gezelter |
1287 |
rb = static_cast<RigidBody*>(integrableObject); |
389 |
|
|
std::vector<Atom*> atoms = rb->getAtoms(); |
390 |
|
|
std::set<int> rigidAtoms; |
391 |
|
|
for (int i = 0; i < static_cast<int>(atoms.size()); ++i) { |
392 |
|
|
rigidAtoms.insert(atoms[i]->getGlobalIndex()); |
393 |
|
|
} |
394 |
|
|
for (int i = 0; i < static_cast<int>(atoms.size()); ++i) { |
395 |
|
|
atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms)); |
396 |
|
|
} |
397 |
tim |
749 |
} else { |
398 |
|
|
std::set<int> oneAtomSet; |
399 |
|
|
oneAtomSet.insert(integrableObject->getGlobalIndex()); |
400 |
|
|
atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet)); |
401 |
|
|
} |
402 |
|
|
} |
403 |
gezelter |
1287 |
|
404 |
|
|
for (bond= mol->beginBond(bondIter); bond != NULL; |
405 |
|
|
bond = mol->nextBond(bondIter)) { |
406 |
tim |
749 |
|
407 |
gezelter |
1287 |
a = bond->getAtomA()->getGlobalIndex(); |
408 |
|
|
b = bond->getAtomB()->getGlobalIndex(); |
409 |
tim |
749 |
|
410 |
gezelter |
1287 |
if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { |
411 |
|
|
oneTwoInteractions_.addPair(a, b); |
412 |
|
|
} else { |
413 |
|
|
excludedInteractions_.addPair(a, b); |
414 |
|
|
} |
415 |
gezelter |
246 |
} |
416 |
gezelter |
2 |
|
417 |
gezelter |
1287 |
for (bend= mol->beginBend(bendIter); bend != NULL; |
418 |
|
|
bend = mol->nextBend(bendIter)) { |
419 |
|
|
|
420 |
gezelter |
507 |
a = bend->getAtomA()->getGlobalIndex(); |
421 |
|
|
b = bend->getAtomB()->getGlobalIndex(); |
422 |
|
|
c = bend->getAtomC()->getGlobalIndex(); |
423 |
gezelter |
1287 |
|
424 |
|
|
if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { |
425 |
|
|
oneTwoInteractions_.addPair(a, b); |
426 |
|
|
oneTwoInteractions_.addPair(b, c); |
427 |
|
|
} else { |
428 |
|
|
excludedInteractions_.addPair(a, b); |
429 |
|
|
excludedInteractions_.addPair(b, c); |
430 |
|
|
} |
431 |
gezelter |
2 |
|
432 |
gezelter |
1287 |
if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) { |
433 |
|
|
oneThreeInteractions_.addPair(a, c); |
434 |
|
|
} else { |
435 |
|
|
excludedInteractions_.addPair(a, c); |
436 |
|
|
} |
437 |
gezelter |
246 |
} |
438 |
gezelter |
2 |
|
439 |
gezelter |
1287 |
for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; |
440 |
|
|
torsion = mol->nextTorsion(torsionIter)) { |
441 |
|
|
|
442 |
gezelter |
507 |
a = torsion->getAtomA()->getGlobalIndex(); |
443 |
|
|
b = torsion->getAtomB()->getGlobalIndex(); |
444 |
|
|
c = torsion->getAtomC()->getGlobalIndex(); |
445 |
gezelter |
1287 |
d = torsion->getAtomD()->getGlobalIndex(); |
446 |
cli2 |
1290 |
|
447 |
gezelter |
1287 |
if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { |
448 |
|
|
oneTwoInteractions_.addPair(a, b); |
449 |
|
|
oneTwoInteractions_.addPair(b, c); |
450 |
|
|
oneTwoInteractions_.addPair(c, d); |
451 |
|
|
} else { |
452 |
|
|
excludedInteractions_.addPair(a, b); |
453 |
|
|
excludedInteractions_.addPair(b, c); |
454 |
|
|
excludedInteractions_.addPair(c, d); |
455 |
|
|
} |
456 |
gezelter |
2 |
|
457 |
gezelter |
1287 |
if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) { |
458 |
|
|
oneThreeInteractions_.addPair(a, c); |
459 |
|
|
oneThreeInteractions_.addPair(b, d); |
460 |
|
|
} else { |
461 |
|
|
excludedInteractions_.addPair(a, c); |
462 |
|
|
excludedInteractions_.addPair(b, d); |
463 |
|
|
} |
464 |
tim |
749 |
|
465 |
gezelter |
1287 |
if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) { |
466 |
|
|
oneFourInteractions_.addPair(a, d); |
467 |
|
|
} else { |
468 |
|
|
excludedInteractions_.addPair(a, d); |
469 |
|
|
} |
470 |
gezelter |
2 |
} |
471 |
|
|
|
472 |
gezelter |
1277 |
for (inversion= mol->beginInversion(inversionIter); inversion != NULL; |
473 |
|
|
inversion = mol->nextInversion(inversionIter)) { |
474 |
gezelter |
1287 |
|
475 |
gezelter |
1277 |
a = inversion->getAtomA()->getGlobalIndex(); |
476 |
|
|
b = inversion->getAtomB()->getGlobalIndex(); |
477 |
|
|
c = inversion->getAtomC()->getGlobalIndex(); |
478 |
|
|
d = inversion->getAtomD()->getGlobalIndex(); |
479 |
|
|
|
480 |
gezelter |
1287 |
if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { |
481 |
|
|
oneTwoInteractions_.addPair(a, b); |
482 |
|
|
oneTwoInteractions_.addPair(a, c); |
483 |
|
|
oneTwoInteractions_.addPair(a, d); |
484 |
|
|
} else { |
485 |
|
|
excludedInteractions_.addPair(a, b); |
486 |
|
|
excludedInteractions_.addPair(a, c); |
487 |
|
|
excludedInteractions_.addPair(a, d); |
488 |
|
|
} |
489 |
gezelter |
1277 |
|
490 |
gezelter |
1287 |
if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) { |
491 |
|
|
oneThreeInteractions_.addPair(b, c); |
492 |
|
|
oneThreeInteractions_.addPair(b, d); |
493 |
|
|
oneThreeInteractions_.addPair(c, d); |
494 |
|
|
} else { |
495 |
|
|
excludedInteractions_.addPair(b, c); |
496 |
|
|
excludedInteractions_.addPair(b, d); |
497 |
|
|
excludedInteractions_.addPair(c, d); |
498 |
|
|
} |
499 |
gezelter |
1277 |
} |
500 |
|
|
|
501 |
gezelter |
1287 |
for (rb = mol->beginRigidBody(rbIter); rb != NULL; |
502 |
|
|
rb = mol->nextRigidBody(rbIter)) { |
503 |
gezelter |
507 |
std::vector<Atom*> atoms = rb->getAtoms(); |
504 |
gezelter |
1287 |
for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) { |
505 |
|
|
for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) { |
506 |
gezelter |
507 |
a = atoms[i]->getGlobalIndex(); |
507 |
|
|
b = atoms[j]->getGlobalIndex(); |
508 |
gezelter |
1287 |
excludedInteractions_.addPair(a, b); |
509 |
gezelter |
507 |
} |
510 |
|
|
} |
511 |
tim |
430 |
} |
512 |
|
|
|
513 |
gezelter |
507 |
} |
514 |
gezelter |
246 |
|
515 |
gezelter |
1287 |
void SimInfo::removeInteractionPairs(Molecule* mol) { |
516 |
|
|
ForceFieldOptions& options_ = forceField_->getForceFieldOptions(); |
517 |
gezelter |
246 |
std::vector<Bond*>::iterator bondIter; |
518 |
|
|
std::vector<Bend*>::iterator bendIter; |
519 |
|
|
std::vector<Torsion*>::iterator torsionIter; |
520 |
gezelter |
1277 |
std::vector<Inversion*>::iterator inversionIter; |
521 |
gezelter |
246 |
Bond* bond; |
522 |
|
|
Bend* bend; |
523 |
|
|
Torsion* torsion; |
524 |
gezelter |
1277 |
Inversion* inversion; |
525 |
gezelter |
246 |
int a; |
526 |
|
|
int b; |
527 |
|
|
int c; |
528 |
|
|
int d; |
529 |
tim |
749 |
|
530 |
|
|
std::map<int, std::set<int> > atomGroups; |
531 |
|
|
Molecule::RigidBodyIterator rbIter; |
532 |
|
|
RigidBody* rb; |
533 |
|
|
Molecule::IntegrableObjectIterator ii; |
534 |
|
|
StuntDouble* integrableObject; |
535 |
gezelter |
246 |
|
536 |
gezelter |
1287 |
for (integrableObject = mol->beginIntegrableObject(ii); |
537 |
|
|
integrableObject != NULL; |
538 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
539 |
|
|
|
540 |
tim |
749 |
if (integrableObject->isRigidBody()) { |
541 |
gezelter |
1287 |
rb = static_cast<RigidBody*>(integrableObject); |
542 |
|
|
std::vector<Atom*> atoms = rb->getAtoms(); |
543 |
|
|
std::set<int> rigidAtoms; |
544 |
|
|
for (int i = 0; i < static_cast<int>(atoms.size()); ++i) { |
545 |
|
|
rigidAtoms.insert(atoms[i]->getGlobalIndex()); |
546 |
|
|
} |
547 |
|
|
for (int i = 0; i < static_cast<int>(atoms.size()); ++i) { |
548 |
|
|
atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms)); |
549 |
|
|
} |
550 |
tim |
749 |
} else { |
551 |
|
|
std::set<int> oneAtomSet; |
552 |
|
|
oneAtomSet.insert(integrableObject->getGlobalIndex()); |
553 |
|
|
atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet)); |
554 |
|
|
} |
555 |
|
|
} |
556 |
|
|
|
557 |
gezelter |
1287 |
for (bond= mol->beginBond(bondIter); bond != NULL; |
558 |
|
|
bond = mol->nextBond(bondIter)) { |
559 |
|
|
|
560 |
|
|
a = bond->getAtomA()->getGlobalIndex(); |
561 |
|
|
b = bond->getAtomB()->getGlobalIndex(); |
562 |
tim |
749 |
|
563 |
gezelter |
1287 |
if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { |
564 |
|
|
oneTwoInteractions_.removePair(a, b); |
565 |
|
|
} else { |
566 |
|
|
excludedInteractions_.removePair(a, b); |
567 |
|
|
} |
568 |
gezelter |
2 |
} |
569 |
gezelter |
246 |
|
570 |
gezelter |
1287 |
for (bend= mol->beginBend(bendIter); bend != NULL; |
571 |
|
|
bend = mol->nextBend(bendIter)) { |
572 |
|
|
|
573 |
gezelter |
507 |
a = bend->getAtomA()->getGlobalIndex(); |
574 |
|
|
b = bend->getAtomB()->getGlobalIndex(); |
575 |
|
|
c = bend->getAtomC()->getGlobalIndex(); |
576 |
gezelter |
1287 |
|
577 |
|
|
if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { |
578 |
|
|
oneTwoInteractions_.removePair(a, b); |
579 |
|
|
oneTwoInteractions_.removePair(b, c); |
580 |
|
|
} else { |
581 |
|
|
excludedInteractions_.removePair(a, b); |
582 |
|
|
excludedInteractions_.removePair(b, c); |
583 |
|
|
} |
584 |
gezelter |
246 |
|
585 |
gezelter |
1287 |
if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) { |
586 |
|
|
oneThreeInteractions_.removePair(a, c); |
587 |
|
|
} else { |
588 |
|
|
excludedInteractions_.removePair(a, c); |
589 |
|
|
} |
590 |
gezelter |
2 |
} |
591 |
gezelter |
246 |
|
592 |
gezelter |
1287 |
for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; |
593 |
|
|
torsion = mol->nextTorsion(torsionIter)) { |
594 |
|
|
|
595 |
gezelter |
507 |
a = torsion->getAtomA()->getGlobalIndex(); |
596 |
|
|
b = torsion->getAtomB()->getGlobalIndex(); |
597 |
|
|
c = torsion->getAtomC()->getGlobalIndex(); |
598 |
gezelter |
1287 |
d = torsion->getAtomD()->getGlobalIndex(); |
599 |
|
|
|
600 |
|
|
if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { |
601 |
|
|
oneTwoInteractions_.removePair(a, b); |
602 |
|
|
oneTwoInteractions_.removePair(b, c); |
603 |
|
|
oneTwoInteractions_.removePair(c, d); |
604 |
|
|
} else { |
605 |
|
|
excludedInteractions_.removePair(a, b); |
606 |
|
|
excludedInteractions_.removePair(b, c); |
607 |
|
|
excludedInteractions_.removePair(c, d); |
608 |
|
|
} |
609 |
gezelter |
246 |
|
610 |
gezelter |
1287 |
if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) { |
611 |
|
|
oneThreeInteractions_.removePair(a, c); |
612 |
|
|
oneThreeInteractions_.removePair(b, d); |
613 |
|
|
} else { |
614 |
|
|
excludedInteractions_.removePair(a, c); |
615 |
|
|
excludedInteractions_.removePair(b, d); |
616 |
|
|
} |
617 |
tim |
749 |
|
618 |
gezelter |
1287 |
if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) { |
619 |
|
|
oneFourInteractions_.removePair(a, d); |
620 |
|
|
} else { |
621 |
|
|
excludedInteractions_.removePair(a, d); |
622 |
|
|
} |
623 |
|
|
} |
624 |
tim |
749 |
|
625 |
gezelter |
1287 |
for (inversion= mol->beginInversion(inversionIter); inversion != NULL; |
626 |
|
|
inversion = mol->nextInversion(inversionIter)) { |
627 |
tim |
749 |
|
628 |
gezelter |
1277 |
a = inversion->getAtomA()->getGlobalIndex(); |
629 |
|
|
b = inversion->getAtomB()->getGlobalIndex(); |
630 |
|
|
c = inversion->getAtomC()->getGlobalIndex(); |
631 |
|
|
d = inversion->getAtomD()->getGlobalIndex(); |
632 |
|
|
|
633 |
gezelter |
1287 |
if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { |
634 |
|
|
oneTwoInteractions_.removePair(a, b); |
635 |
|
|
oneTwoInteractions_.removePair(a, c); |
636 |
|
|
oneTwoInteractions_.removePair(a, d); |
637 |
|
|
} else { |
638 |
|
|
excludedInteractions_.removePair(a, b); |
639 |
|
|
excludedInteractions_.removePair(a, c); |
640 |
|
|
excludedInteractions_.removePair(a, d); |
641 |
|
|
} |
642 |
gezelter |
1277 |
|
643 |
gezelter |
1287 |
if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) { |
644 |
|
|
oneThreeInteractions_.removePair(b, c); |
645 |
|
|
oneThreeInteractions_.removePair(b, d); |
646 |
|
|
oneThreeInteractions_.removePair(c, d); |
647 |
|
|
} else { |
648 |
|
|
excludedInteractions_.removePair(b, c); |
649 |
|
|
excludedInteractions_.removePair(b, d); |
650 |
|
|
excludedInteractions_.removePair(c, d); |
651 |
|
|
} |
652 |
gezelter |
1277 |
} |
653 |
|
|
|
654 |
gezelter |
1287 |
for (rb = mol->beginRigidBody(rbIter); rb != NULL; |
655 |
|
|
rb = mol->nextRigidBody(rbIter)) { |
656 |
gezelter |
507 |
std::vector<Atom*> atoms = rb->getAtoms(); |
657 |
gezelter |
1287 |
for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) { |
658 |
|
|
for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) { |
659 |
gezelter |
507 |
a = atoms[i]->getGlobalIndex(); |
660 |
|
|
b = atoms[j]->getGlobalIndex(); |
661 |
gezelter |
1287 |
excludedInteractions_.removePair(a, b); |
662 |
gezelter |
507 |
} |
663 |
|
|
} |
664 |
tim |
430 |
} |
665 |
gezelter |
1287 |
|
666 |
gezelter |
507 |
} |
667 |
gezelter |
1287 |
|
668 |
|
|
|
669 |
gezelter |
507 |
void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) { |
670 |
gezelter |
246 |
int curStampId; |
671 |
gezelter |
1287 |
|
672 |
gezelter |
246 |
//index from 0 |
673 |
|
|
curStampId = moleculeStamps_.size(); |
674 |
gezelter |
2 |
|
675 |
gezelter |
246 |
moleculeStamps_.push_back(molStamp); |
676 |
|
|
molStampIds_.insert(molStampIds_.end(), nmol, curStampId); |
677 |
gezelter |
507 |
} |
678 |
gezelter |
2 |
|
679 |
gezelter |
507 |
void SimInfo::update() { |
680 |
gezelter |
2 |
|
681 |
gezelter |
246 |
setupSimType(); |
682 |
gezelter |
2 |
|
683 |
gezelter |
246 |
#ifdef IS_MPI |
684 |
|
|
setupFortranParallel(); |
685 |
|
|
#endif |
686 |
gezelter |
2 |
|
687 |
gezelter |
246 |
setupFortranSim(); |
688 |
gezelter |
2 |
|
689 |
gezelter |
246 |
//setup fortran force field |
690 |
|
|
/** @deprecate */ |
691 |
|
|
int isError = 0; |
692 |
chrisfen |
598 |
|
693 |
chrisfen |
1045 |
setupCutoff(); |
694 |
|
|
|
695 |
chrisfen |
603 |
setupElectrostaticSummationMethod( isError ); |
696 |
chrisfen |
726 |
setupSwitchingFunction(); |
697 |
chrisfen |
998 |
setupAccumulateBoxDipole(); |
698 |
chrisfen |
598 |
|
699 |
gezelter |
246 |
if(isError){ |
700 |
gezelter |
507 |
sprintf( painCave.errMsg, |
701 |
|
|
"ForceField error: There was an error initializing the forceField in fortran.\n" ); |
702 |
|
|
painCave.isFatal = 1; |
703 |
|
|
simError(); |
704 |
gezelter |
246 |
} |
705 |
gezelter |
2 |
|
706 |
gezelter |
246 |
calcNdf(); |
707 |
|
|
calcNdfRaw(); |
708 |
|
|
calcNdfTrans(); |
709 |
|
|
|
710 |
|
|
fortranInitialized_ = true; |
711 |
gezelter |
507 |
} |
712 |
gezelter |
2 |
|
713 |
gezelter |
507 |
std::set<AtomType*> SimInfo::getUniqueAtomTypes() { |
714 |
gezelter |
246 |
SimInfo::MoleculeIterator mi; |
715 |
|
|
Molecule* mol; |
716 |
|
|
Molecule::AtomIterator ai; |
717 |
|
|
Atom* atom; |
718 |
|
|
std::set<AtomType*> atomTypes; |
719 |
gezelter |
2 |
|
720 |
gezelter |
246 |
for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { |
721 |
gezelter |
2 |
|
722 |
gezelter |
507 |
for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { |
723 |
|
|
atomTypes.insert(atom->getAtomType()); |
724 |
|
|
} |
725 |
gezelter |
246 |
|
726 |
|
|
} |
727 |
gezelter |
2 |
|
728 |
gezelter |
246 |
return atomTypes; |
729 |
gezelter |
507 |
} |
730 |
gezelter |
2 |
|
731 |
gezelter |
507 |
void SimInfo::setupSimType() { |
732 |
gezelter |
246 |
std::set<AtomType*>::iterator i; |
733 |
|
|
std::set<AtomType*> atomTypes; |
734 |
|
|
atomTypes = getUniqueAtomTypes(); |
735 |
gezelter |
2 |
|
736 |
gezelter |
246 |
int useLennardJones = 0; |
737 |
|
|
int useElectrostatic = 0; |
738 |
|
|
int useEAM = 0; |
739 |
chuckv |
734 |
int useSC = 0; |
740 |
gezelter |
246 |
int useCharge = 0; |
741 |
|
|
int useDirectional = 0; |
742 |
|
|
int useDipole = 0; |
743 |
|
|
int useGayBerne = 0; |
744 |
|
|
int useSticky = 0; |
745 |
chrisfen |
523 |
int useStickyPower = 0; |
746 |
gezelter |
246 |
int useShape = 0; |
747 |
|
|
int useFLARB = 0; //it is not in AtomType yet |
748 |
|
|
int useDirectionalAtom = 0; |
749 |
|
|
int useElectrostatics = 0; |
750 |
|
|
//usePBC and useRF are from simParams |
751 |
tim |
665 |
int usePBC = simParams_->getUsePeriodicBoundaryConditions(); |
752 |
chrisfen |
611 |
int useRF; |
753 |
chrisfen |
720 |
int useSF; |
754 |
chrisfen |
998 |
int useSP; |
755 |
|
|
int useBoxDipole; |
756 |
gezelter |
1126 |
|
757 |
tim |
665 |
std::string myMethod; |
758 |
gezelter |
2 |
|
759 |
chrisfen |
611 |
// set the useRF logical |
760 |
tim |
665 |
useRF = 0; |
761 |
chrisfen |
720 |
useSF = 0; |
762 |
gezelter |
1078 |
useSP = 0; |
763 |
gezelter |
1313 |
useBoxDipole = 0; |
764 |
chrisfen |
691 |
|
765 |
|
|
|
766 |
tim |
665 |
if (simParams_->haveElectrostaticSummationMethod()) { |
767 |
chrisfen |
691 |
std::string myMethod = simParams_->getElectrostaticSummationMethod(); |
768 |
|
|
toUpper(myMethod); |
769 |
chrisfen |
998 |
if (myMethod == "REACTION_FIELD"){ |
770 |
gezelter |
1078 |
useRF = 1; |
771 |
chrisfen |
998 |
} else if (myMethod == "SHIFTED_FORCE"){ |
772 |
|
|
useSF = 1; |
773 |
|
|
} else if (myMethod == "SHIFTED_POTENTIAL"){ |
774 |
|
|
useSP = 1; |
775 |
chrisfen |
691 |
} |
776 |
tim |
665 |
} |
777 |
chrisfen |
998 |
|
778 |
|
|
if (simParams_->haveAccumulateBoxDipole()) |
779 |
|
|
if (simParams_->getAccumulateBoxDipole()) |
780 |
|
|
useBoxDipole = 1; |
781 |
chrisfen |
611 |
|
782 |
gezelter |
1126 |
useAtomicVirial_ = simParams_->getUseAtomicVirial(); |
783 |
|
|
|
784 |
gezelter |
246 |
//loop over all of the atom types |
785 |
|
|
for (i = atomTypes.begin(); i != atomTypes.end(); ++i) { |
786 |
gezelter |
507 |
useLennardJones |= (*i)->isLennardJones(); |
787 |
|
|
useElectrostatic |= (*i)->isElectrostatic(); |
788 |
|
|
useEAM |= (*i)->isEAM(); |
789 |
chuckv |
734 |
useSC |= (*i)->isSC(); |
790 |
gezelter |
507 |
useCharge |= (*i)->isCharge(); |
791 |
|
|
useDirectional |= (*i)->isDirectional(); |
792 |
|
|
useDipole |= (*i)->isDipole(); |
793 |
|
|
useGayBerne |= (*i)->isGayBerne(); |
794 |
|
|
useSticky |= (*i)->isSticky(); |
795 |
chrisfen |
523 |
useStickyPower |= (*i)->isStickyPower(); |
796 |
gezelter |
507 |
useShape |= (*i)->isShape(); |
797 |
gezelter |
246 |
} |
798 |
gezelter |
2 |
|
799 |
chrisfen |
523 |
if (useSticky || useStickyPower || useDipole || useGayBerne || useShape) { |
800 |
gezelter |
507 |
useDirectionalAtom = 1; |
801 |
gezelter |
246 |
} |
802 |
gezelter |
2 |
|
803 |
gezelter |
246 |
if (useCharge || useDipole) { |
804 |
gezelter |
507 |
useElectrostatics = 1; |
805 |
gezelter |
246 |
} |
806 |
gezelter |
2 |
|
807 |
gezelter |
246 |
#ifdef IS_MPI |
808 |
|
|
int temp; |
809 |
gezelter |
2 |
|
810 |
gezelter |
246 |
temp = usePBC; |
811 |
|
|
MPI_Allreduce(&temp, &usePBC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
812 |
gezelter |
2 |
|
813 |
gezelter |
246 |
temp = useDirectionalAtom; |
814 |
|
|
MPI_Allreduce(&temp, &useDirectionalAtom, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
815 |
gezelter |
2 |
|
816 |
gezelter |
246 |
temp = useLennardJones; |
817 |
|
|
MPI_Allreduce(&temp, &useLennardJones, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
818 |
gezelter |
2 |
|
819 |
gezelter |
246 |
temp = useElectrostatics; |
820 |
|
|
MPI_Allreduce(&temp, &useElectrostatics, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
821 |
gezelter |
2 |
|
822 |
gezelter |
246 |
temp = useCharge; |
823 |
|
|
MPI_Allreduce(&temp, &useCharge, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
824 |
gezelter |
2 |
|
825 |
gezelter |
246 |
temp = useDipole; |
826 |
|
|
MPI_Allreduce(&temp, &useDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
827 |
gezelter |
2 |
|
828 |
gezelter |
246 |
temp = useSticky; |
829 |
|
|
MPI_Allreduce(&temp, &useSticky, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
830 |
gezelter |
2 |
|
831 |
chrisfen |
523 |
temp = useStickyPower; |
832 |
|
|
MPI_Allreduce(&temp, &useStickyPower, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
833 |
|
|
|
834 |
gezelter |
246 |
temp = useGayBerne; |
835 |
|
|
MPI_Allreduce(&temp, &useGayBerne, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
836 |
gezelter |
2 |
|
837 |
gezelter |
246 |
temp = useEAM; |
838 |
|
|
MPI_Allreduce(&temp, &useEAM, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
839 |
gezelter |
2 |
|
840 |
chuckv |
734 |
temp = useSC; |
841 |
|
|
MPI_Allreduce(&temp, &useSC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
842 |
|
|
|
843 |
gezelter |
246 |
temp = useShape; |
844 |
|
|
MPI_Allreduce(&temp, &useShape, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
845 |
|
|
|
846 |
|
|
temp = useFLARB; |
847 |
|
|
MPI_Allreduce(&temp, &useFLARB, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
848 |
|
|
|
849 |
chrisfen |
611 |
temp = useRF; |
850 |
|
|
MPI_Allreduce(&temp, &useRF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
851 |
|
|
|
852 |
chrisfen |
720 |
temp = useSF; |
853 |
chrisfen |
998 |
MPI_Allreduce(&temp, &useSF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
854 |
chrisfen |
705 |
|
855 |
chrisfen |
998 |
temp = useSP; |
856 |
|
|
MPI_Allreduce(&temp, &useSP, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
857 |
|
|
|
858 |
|
|
temp = useBoxDipole; |
859 |
|
|
MPI_Allreduce(&temp, &useBoxDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
860 |
|
|
|
861 |
gezelter |
1126 |
temp = useAtomicVirial_; |
862 |
|
|
MPI_Allreduce(&temp, &useAtomicVirial_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
863 |
|
|
|
864 |
gezelter |
2 |
#endif |
865 |
|
|
|
866 |
gezelter |
246 |
fInfo_.SIM_uses_PBC = usePBC; |
867 |
|
|
fInfo_.SIM_uses_DirectionalAtoms = useDirectionalAtom; |
868 |
|
|
fInfo_.SIM_uses_LennardJones = useLennardJones; |
869 |
|
|
fInfo_.SIM_uses_Electrostatics = useElectrostatics; |
870 |
|
|
fInfo_.SIM_uses_Charges = useCharge; |
871 |
|
|
fInfo_.SIM_uses_Dipoles = useDipole; |
872 |
|
|
fInfo_.SIM_uses_Sticky = useSticky; |
873 |
chrisfen |
523 |
fInfo_.SIM_uses_StickyPower = useStickyPower; |
874 |
gezelter |
246 |
fInfo_.SIM_uses_GayBerne = useGayBerne; |
875 |
|
|
fInfo_.SIM_uses_EAM = useEAM; |
876 |
chuckv |
734 |
fInfo_.SIM_uses_SC = useSC; |
877 |
gezelter |
246 |
fInfo_.SIM_uses_Shapes = useShape; |
878 |
|
|
fInfo_.SIM_uses_FLARB = useFLARB; |
879 |
chrisfen |
611 |
fInfo_.SIM_uses_RF = useRF; |
880 |
chrisfen |
720 |
fInfo_.SIM_uses_SF = useSF; |
881 |
chrisfen |
998 |
fInfo_.SIM_uses_SP = useSP; |
882 |
|
|
fInfo_.SIM_uses_BoxDipole = useBoxDipole; |
883 |
gezelter |
1126 |
fInfo_.SIM_uses_AtomicVirial = useAtomicVirial_; |
884 |
gezelter |
507 |
} |
885 |
gezelter |
2 |
|
886 |
gezelter |
507 |
void SimInfo::setupFortranSim() { |
887 |
gezelter |
246 |
int isError; |
888 |
gezelter |
1287 |
int nExclude, nOneTwo, nOneThree, nOneFour; |
889 |
gezelter |
246 |
std::vector<int> fortranGlobalGroupMembership; |
890 |
|
|
|
891 |
|
|
isError = 0; |
892 |
gezelter |
2 |
|
893 |
gezelter |
246 |
//globalGroupMembership_ is filled by SimCreator |
894 |
|
|
for (int i = 0; i < nGlobalAtoms_; i++) { |
895 |
gezelter |
507 |
fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1); |
896 |
gezelter |
246 |
} |
897 |
gezelter |
2 |
|
898 |
gezelter |
246 |
//calculate mass ratio of cutoff group |
899 |
tim |
963 |
std::vector<RealType> mfact; |
900 |
gezelter |
246 |
SimInfo::MoleculeIterator mi; |
901 |
|
|
Molecule* mol; |
902 |
|
|
Molecule::CutoffGroupIterator ci; |
903 |
|
|
CutoffGroup* cg; |
904 |
|
|
Molecule::AtomIterator ai; |
905 |
|
|
Atom* atom; |
906 |
tim |
963 |
RealType totalMass; |
907 |
gezelter |
246 |
|
908 |
|
|
//to avoid memory reallocation, reserve enough space for mfact |
909 |
|
|
mfact.reserve(getNCutoffGroups()); |
910 |
gezelter |
2 |
|
911 |
gezelter |
246 |
for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { |
912 |
gezelter |
507 |
for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) { |
913 |
gezelter |
2 |
|
914 |
gezelter |
507 |
totalMass = cg->getMass(); |
915 |
|
|
for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) { |
916 |
chrisfen |
645 |
// Check for massless groups - set mfact to 1 if true |
917 |
|
|
if (totalMass != 0) |
918 |
|
|
mfact.push_back(atom->getMass()/totalMass); |
919 |
|
|
else |
920 |
|
|
mfact.push_back( 1.0 ); |
921 |
gezelter |
507 |
} |
922 |
|
|
} |
923 |
gezelter |
246 |
} |
924 |
gezelter |
2 |
|
925 |
gezelter |
246 |
//fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!) |
926 |
|
|
std::vector<int> identArray; |
927 |
gezelter |
2 |
|
928 |
gezelter |
246 |
//to avoid memory reallocation, reserve enough space identArray |
929 |
|
|
identArray.reserve(getNAtoms()); |
930 |
|
|
|
931 |
|
|
for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { |
932 |
gezelter |
507 |
for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { |
933 |
|
|
identArray.push_back(atom->getIdent()); |
934 |
|
|
} |
935 |
gezelter |
246 |
} |
936 |
gezelter |
2 |
|
937 |
gezelter |
246 |
//fill molMembershipArray |
938 |
|
|
//molMembershipArray is filled by SimCreator |
939 |
|
|
std::vector<int> molMembershipArray(nGlobalAtoms_); |
940 |
|
|
for (int i = 0; i < nGlobalAtoms_; i++) { |
941 |
gezelter |
507 |
molMembershipArray[i] = globalMolMembership_[i] + 1; |
942 |
gezelter |
246 |
} |
943 |
|
|
|
944 |
|
|
//setup fortran simulation |
945 |
gezelter |
1287 |
|
946 |
|
|
nExclude = excludedInteractions_.getSize(); |
947 |
|
|
nOneTwo = oneTwoInteractions_.getSize(); |
948 |
|
|
nOneThree = oneThreeInteractions_.getSize(); |
949 |
|
|
nOneFour = oneFourInteractions_.getSize(); |
950 |
|
|
|
951 |
|
|
int* excludeList = excludedInteractions_.getPairList(); |
952 |
|
|
int* oneTwoList = oneTwoInteractions_.getPairList(); |
953 |
|
|
int* oneThreeList = oneThreeInteractions_.getPairList(); |
954 |
|
|
int* oneFourList = oneFourInteractions_.getPairList(); |
955 |
|
|
|
956 |
gezelter |
1241 |
setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0], |
957 |
gezelter |
1287 |
&nExclude, excludeList, |
958 |
|
|
&nOneTwo, oneTwoList, |
959 |
|
|
&nOneThree, oneThreeList, |
960 |
|
|
&nOneFour, oneFourList, |
961 |
gezelter |
1241 |
&molMembershipArray[0], &mfact[0], &nCutoffGroups_, |
962 |
|
|
&fortranGlobalGroupMembership[0], &isError); |
963 |
|
|
|
964 |
gezelter |
246 |
if( isError ){ |
965 |
gezelter |
1241 |
|
966 |
gezelter |
507 |
sprintf( painCave.errMsg, |
967 |
|
|
"There was an error setting the simulation information in fortran.\n" ); |
968 |
|
|
painCave.isFatal = 1; |
969 |
gezelter |
1390 |
painCave.severity = OPENMD_ERROR; |
970 |
gezelter |
507 |
simError(); |
971 |
gezelter |
246 |
} |
972 |
gezelter |
1241 |
|
973 |
|
|
|
974 |
gezelter |
246 |
sprintf( checkPointMsg, |
975 |
gezelter |
507 |
"succesfully sent the simulation information to fortran.\n"); |
976 |
gezelter |
1241 |
|
977 |
|
|
errorCheckPoint(); |
978 |
|
|
|
979 |
chuckv |
1095 |
// Setup number of neighbors in neighbor list if present |
980 |
|
|
if (simParams_->haveNeighborListNeighbors()) { |
981 |
chuckv |
1121 |
int nlistNeighbors = simParams_->getNeighborListNeighbors(); |
982 |
|
|
setNeighbors(&nlistNeighbors); |
983 |
chuckv |
1095 |
} |
984 |
|
|
|
985 |
|
|
|
986 |
gezelter |
507 |
} |
987 |
gezelter |
2 |
|
988 |
|
|
|
989 |
gezelter |
507 |
void SimInfo::setupFortranParallel() { |
990 |
gezelter |
1241 |
#ifdef IS_MPI |
991 |
gezelter |
246 |
//SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex |
992 |
|
|
std::vector<int> localToGlobalAtomIndex(getNAtoms(), 0); |
993 |
|
|
std::vector<int> localToGlobalCutoffGroupIndex; |
994 |
|
|
SimInfo::MoleculeIterator mi; |
995 |
|
|
Molecule::AtomIterator ai; |
996 |
|
|
Molecule::CutoffGroupIterator ci; |
997 |
|
|
Molecule* mol; |
998 |
|
|
Atom* atom; |
999 |
|
|
CutoffGroup* cg; |
1000 |
|
|
mpiSimData parallelData; |
1001 |
|
|
int isError; |
1002 |
gezelter |
2 |
|
1003 |
gezelter |
246 |
for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { |
1004 |
gezelter |
2 |
|
1005 |
gezelter |
507 |
//local index(index in DataStorge) of atom is important |
1006 |
|
|
for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { |
1007 |
|
|
localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1; |
1008 |
|
|
} |
1009 |
gezelter |
2 |
|
1010 |
gezelter |
507 |
//local index of cutoff group is trivial, it only depends on the order of travesing |
1011 |
|
|
for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) { |
1012 |
|
|
localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1); |
1013 |
|
|
} |
1014 |
gezelter |
246 |
|
1015 |
|
|
} |
1016 |
gezelter |
2 |
|
1017 |
gezelter |
246 |
//fill up mpiSimData struct |
1018 |
|
|
parallelData.nMolGlobal = getNGlobalMolecules(); |
1019 |
|
|
parallelData.nMolLocal = getNMolecules(); |
1020 |
|
|
parallelData.nAtomsGlobal = getNGlobalAtoms(); |
1021 |
|
|
parallelData.nAtomsLocal = getNAtoms(); |
1022 |
|
|
parallelData.nGroupsGlobal = getNGlobalCutoffGroups(); |
1023 |
|
|
parallelData.nGroupsLocal = getNCutoffGroups(); |
1024 |
|
|
parallelData.myNode = worldRank; |
1025 |
|
|
MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors)); |
1026 |
gezelter |
2 |
|
1027 |
gezelter |
246 |
//pass mpiSimData struct and index arrays to fortran |
1028 |
|
|
setFsimParallel(¶llelData, &(parallelData.nAtomsLocal), |
1029 |
|
|
&localToGlobalAtomIndex[0], &(parallelData.nGroupsLocal), |
1030 |
|
|
&localToGlobalCutoffGroupIndex[0], &isError); |
1031 |
gezelter |
2 |
|
1032 |
gezelter |
246 |
if (isError) { |
1033 |
gezelter |
507 |
sprintf(painCave.errMsg, |
1034 |
|
|
"mpiRefresh errror: fortran didn't like something we gave it.\n"); |
1035 |
|
|
painCave.isFatal = 1; |
1036 |
|
|
simError(); |
1037 |
gezelter |
246 |
} |
1038 |
gezelter |
2 |
|
1039 |
gezelter |
246 |
sprintf(checkPointMsg, " mpiRefresh successful.\n"); |
1040 |
gezelter |
1241 |
errorCheckPoint(); |
1041 |
gezelter |
2 |
|
1042 |
gezelter |
1241 |
#endif |
1043 |
gezelter |
507 |
} |
1044 |
chrisfen |
143 |
|
1045 |
gezelter |
764 |
void SimInfo::setupCutoff() { |
1046 |
gezelter |
2 |
|
1047 |
chuckv |
834 |
ForceFieldOptions& forceFieldOptions_ = forceField_->getForceFieldOptions(); |
1048 |
|
|
|
1049 |
gezelter |
764 |
// Check the cutoff policy |
1050 |
chuckv |
834 |
int cp = TRADITIONAL_CUTOFF_POLICY; // Set to traditional by default |
1051 |
|
|
|
1052 |
chrisfen |
1129 |
// Set LJ shifting bools to false |
1053 |
gezelter |
1386 |
ljsp_ = 0; |
1054 |
|
|
ljsf_ = 0; |
1055 |
chrisfen |
1129 |
|
1056 |
chuckv |
834 |
std::string myPolicy; |
1057 |
|
|
if (forceFieldOptions_.haveCutoffPolicy()){ |
1058 |
|
|
myPolicy = forceFieldOptions_.getCutoffPolicy(); |
1059 |
|
|
}else if (simParams_->haveCutoffPolicy()) { |
1060 |
|
|
myPolicy = simParams_->getCutoffPolicy(); |
1061 |
|
|
} |
1062 |
|
|
|
1063 |
|
|
if (!myPolicy.empty()){ |
1064 |
tim |
665 |
toUpper(myPolicy); |
1065 |
gezelter |
586 |
if (myPolicy == "MIX") { |
1066 |
|
|
cp = MIX_CUTOFF_POLICY; |
1067 |
|
|
} else { |
1068 |
|
|
if (myPolicy == "MAX") { |
1069 |
|
|
cp = MAX_CUTOFF_POLICY; |
1070 |
|
|
} else { |
1071 |
|
|
if (myPolicy == "TRADITIONAL") { |
1072 |
|
|
cp = TRADITIONAL_CUTOFF_POLICY; |
1073 |
|
|
} else { |
1074 |
|
|
// throw error |
1075 |
|
|
sprintf( painCave.errMsg, |
1076 |
|
|
"SimInfo error: Unknown cutoffPolicy. (Input file specified %s .)\n\tcutoffPolicy must be one of: \"Mix\", \"Max\", or \"Traditional\".", myPolicy.c_str() ); |
1077 |
|
|
painCave.isFatal = 1; |
1078 |
|
|
simError(); |
1079 |
|
|
} |
1080 |
|
|
} |
1081 |
|
|
} |
1082 |
gezelter |
764 |
} |
1083 |
|
|
notifyFortranCutoffPolicy(&cp); |
1084 |
chuckv |
629 |
|
1085 |
gezelter |
764 |
// Check the Skin Thickness for neighborlists |
1086 |
tim |
963 |
RealType skin; |
1087 |
gezelter |
764 |
if (simParams_->haveSkinThickness()) { |
1088 |
|
|
skin = simParams_->getSkinThickness(); |
1089 |
|
|
notifyFortranSkinThickness(&skin); |
1090 |
|
|
} |
1091 |
|
|
|
1092 |
|
|
// Check if the cutoff was set explicitly: |
1093 |
|
|
if (simParams_->haveCutoffRadius()) { |
1094 |
|
|
rcut_ = simParams_->getCutoffRadius(); |
1095 |
|
|
if (simParams_->haveSwitchingRadius()) { |
1096 |
|
|
rsw_ = simParams_->getSwitchingRadius(); |
1097 |
|
|
} else { |
1098 |
chrisfen |
878 |
if (fInfo_.SIM_uses_Charges | |
1099 |
|
|
fInfo_.SIM_uses_Dipoles | |
1100 |
|
|
fInfo_.SIM_uses_RF) { |
1101 |
|
|
|
1102 |
|
|
rsw_ = 0.85 * rcut_; |
1103 |
|
|
sprintf(painCave.errMsg, |
1104 |
|
|
"SimCreator Warning: No value was set for the switchingRadius.\n" |
1105 |
gezelter |
1390 |
"\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n" |
1106 |
chrisfen |
878 |
"\tswitchingRadius = %f. for this simulation\n", rsw_); |
1107 |
|
|
painCave.isFatal = 0; |
1108 |
|
|
simError(); |
1109 |
|
|
} else { |
1110 |
|
|
rsw_ = rcut_; |
1111 |
|
|
sprintf(painCave.errMsg, |
1112 |
|
|
"SimCreator Warning: No value was set for the switchingRadius.\n" |
1113 |
gezelter |
1390 |
"\tOpenMD will use the same value as the cutoffRadius.\n" |
1114 |
chrisfen |
878 |
"\tswitchingRadius = %f. for this simulation\n", rsw_); |
1115 |
|
|
painCave.isFatal = 0; |
1116 |
|
|
simError(); |
1117 |
|
|
} |
1118 |
chrisfen |
879 |
} |
1119 |
chrisfen |
1129 |
|
1120 |
|
|
if (simParams_->haveElectrostaticSummationMethod()) { |
1121 |
|
|
std::string myMethod = simParams_->getElectrostaticSummationMethod(); |
1122 |
|
|
toUpper(myMethod); |
1123 |
|
|
|
1124 |
|
|
if (myMethod == "SHIFTED_POTENTIAL") { |
1125 |
gezelter |
1386 |
ljsp_ = 1; |
1126 |
chrisfen |
1129 |
} else if (myMethod == "SHIFTED_FORCE") { |
1127 |
gezelter |
1386 |
ljsf_ = 1; |
1128 |
chrisfen |
1129 |
} |
1129 |
|
|
} |
1130 |
gezelter |
1386 |
|
1131 |
chrisfen |
1129 |
notifyFortranCutoffs(&rcut_, &rsw_, &ljsp_, &ljsf_); |
1132 |
chrisfen |
879 |
|
1133 |
gezelter |
764 |
} else { |
1134 |
|
|
|
1135 |
|
|
// For electrostatic atoms, we'll assume a large safe value: |
1136 |
|
|
if (fInfo_.SIM_uses_Charges | fInfo_.SIM_uses_Dipoles | fInfo_.SIM_uses_RF) { |
1137 |
|
|
sprintf(painCave.errMsg, |
1138 |
|
|
"SimCreator Warning: No value was set for the cutoffRadius.\n" |
1139 |
gezelter |
1390 |
"\tOpenMD will use a default value of 15.0 angstroms" |
1140 |
gezelter |
764 |
"\tfor the cutoffRadius.\n"); |
1141 |
|
|
painCave.isFatal = 0; |
1142 |
|
|
simError(); |
1143 |
|
|
rcut_ = 15.0; |
1144 |
|
|
|
1145 |
|
|
if (simParams_->haveElectrostaticSummationMethod()) { |
1146 |
|
|
std::string myMethod = simParams_->getElectrostaticSummationMethod(); |
1147 |
|
|
toUpper(myMethod); |
1148 |
gezelter |
1503 |
|
1149 |
|
|
// For the time being, we're tethering the LJ shifted behavior to the |
1150 |
|
|
// electrostaticSummationMethod keyword options |
1151 |
chrisfen |
1129 |
if (myMethod == "SHIFTED_POTENTIAL") { |
1152 |
gezelter |
1386 |
ljsp_ = 1; |
1153 |
chrisfen |
1129 |
} else if (myMethod == "SHIFTED_FORCE") { |
1154 |
gezelter |
1386 |
ljsf_ = 1; |
1155 |
chrisfen |
1129 |
} |
1156 |
|
|
if (myMethod == "SHIFTED_POTENTIAL" || myMethod == "SHIFTED_FORCE") { |
1157 |
gezelter |
764 |
if (simParams_->haveSwitchingRadius()){ |
1158 |
|
|
sprintf(painCave.errMsg, |
1159 |
|
|
"SimInfo Warning: A value was set for the switchingRadius\n" |
1160 |
|
|
"\teven though the electrostaticSummationMethod was\n" |
1161 |
|
|
"\tset to %s\n", myMethod.c_str()); |
1162 |
|
|
painCave.isFatal = 1; |
1163 |
|
|
simError(); |
1164 |
|
|
} |
1165 |
|
|
} |
1166 |
|
|
} |
1167 |
|
|
|
1168 |
|
|
if (simParams_->haveSwitchingRadius()){ |
1169 |
|
|
rsw_ = simParams_->getSwitchingRadius(); |
1170 |
|
|
} else { |
1171 |
|
|
sprintf(painCave.errMsg, |
1172 |
|
|
"SimCreator Warning: No value was set for switchingRadius.\n" |
1173 |
gezelter |
1390 |
"\tOpenMD will use a default value of\n" |
1174 |
gezelter |
764 |
"\t0.85 * cutoffRadius for the switchingRadius\n"); |
1175 |
|
|
painCave.isFatal = 0; |
1176 |
|
|
simError(); |
1177 |
|
|
rsw_ = 0.85 * rcut_; |
1178 |
|
|
} |
1179 |
chrisfen |
1129 |
|
1180 |
gezelter |
1503 |
Electrostatic::setElectrostaticCutoffRadius(rcut_, rsw_); |
1181 |
chrisfen |
1129 |
notifyFortranCutoffs(&rcut_, &rsw_, &ljsp_, &ljsf_); |
1182 |
|
|
|
1183 |
gezelter |
764 |
} else { |
1184 |
|
|
// We didn't set rcut explicitly, and we don't have electrostatic atoms, so |
1185 |
|
|
// We'll punt and let fortran figure out the cutoffs later. |
1186 |
|
|
|
1187 |
|
|
notifyFortranYouAreOnYourOwn(); |
1188 |
chuckv |
629 |
|
1189 |
gezelter |
764 |
} |
1190 |
chuckv |
629 |
} |
1191 |
gezelter |
507 |
} |
1192 |
gezelter |
2 |
|
1193 |
chrisfen |
603 |
void SimInfo::setupElectrostaticSummationMethod( int isError ) { |
1194 |
chrisfen |
598 |
|
1195 |
|
|
int errorOut; |
1196 |
gezelter |
1503 |
ElectrostaticSummationMethod esm = NONE; |
1197 |
|
|
ElectrostaticScreeningMethod sm = UNDAMPED; |
1198 |
tim |
963 |
RealType alphaVal; |
1199 |
|
|
RealType dielectric; |
1200 |
chrisfen |
1045 |
|
1201 |
chrisfen |
598 |
errorOut = isError; |
1202 |
|
|
|
1203 |
chrisfen |
603 |
if (simParams_->haveElectrostaticSummationMethod()) { |
1204 |
chrisfen |
604 |
std::string myMethod = simParams_->getElectrostaticSummationMethod(); |
1205 |
tim |
665 |
toUpper(myMethod); |
1206 |
chrisfen |
603 |
if (myMethod == "NONE") { |
1207 |
|
|
esm = NONE; |
1208 |
chrisfen |
598 |
} else { |
1209 |
chrisfen |
709 |
if (myMethod == "SWITCHING_FUNCTION") { |
1210 |
|
|
esm = SWITCHING_FUNCTION; |
1211 |
chrisfen |
598 |
} else { |
1212 |
chrisfen |
709 |
if (myMethod == "SHIFTED_POTENTIAL") { |
1213 |
|
|
esm = SHIFTED_POTENTIAL; |
1214 |
|
|
} else { |
1215 |
|
|
if (myMethod == "SHIFTED_FORCE") { |
1216 |
|
|
esm = SHIFTED_FORCE; |
1217 |
chrisfen |
598 |
} else { |
1218 |
chrisfen |
1050 |
if (myMethod == "REACTION_FIELD") { |
1219 |
chrisfen |
709 |
esm = REACTION_FIELD; |
1220 |
chrisfen |
1050 |
dielectric = simParams_->getDielectric(); |
1221 |
|
|
if (!simParams_->haveDielectric()) { |
1222 |
|
|
// throw warning |
1223 |
|
|
sprintf( painCave.errMsg, |
1224 |
|
|
"SimInfo warning: dielectric was not specified in the input file\n\tfor the reaction field correction method.\n" |
1225 |
|
|
"\tA default value of %f will be used for the dielectric.\n", dielectric); |
1226 |
|
|
painCave.isFatal = 0; |
1227 |
|
|
simError(); |
1228 |
|
|
} |
1229 |
chrisfen |
709 |
} else { |
1230 |
|
|
// throw error |
1231 |
|
|
sprintf( painCave.errMsg, |
1232 |
gezelter |
764 |
"SimInfo error: Unknown electrostaticSummationMethod.\n" |
1233 |
|
|
"\t(Input file specified %s .)\n" |
1234 |
|
|
"\telectrostaticSummationMethod must be one of: \"none\",\n" |
1235 |
|
|
"\t\"shifted_potential\", \"shifted_force\", or \n" |
1236 |
|
|
"\t\"reaction_field\".\n", myMethod.c_str() ); |
1237 |
chrisfen |
709 |
painCave.isFatal = 1; |
1238 |
|
|
simError(); |
1239 |
|
|
} |
1240 |
|
|
} |
1241 |
|
|
} |
1242 |
chrisfen |
598 |
} |
1243 |
|
|
} |
1244 |
|
|
} |
1245 |
chrisfen |
709 |
|
1246 |
chrisfen |
716 |
if (simParams_->haveElectrostaticScreeningMethod()) { |
1247 |
|
|
std::string myScreen = simParams_->getElectrostaticScreeningMethod(); |
1248 |
chrisfen |
709 |
toUpper(myScreen); |
1249 |
|
|
if (myScreen == "UNDAMPED") { |
1250 |
|
|
sm = UNDAMPED; |
1251 |
|
|
} else { |
1252 |
|
|
if (myScreen == "DAMPED") { |
1253 |
|
|
sm = DAMPED; |
1254 |
|
|
if (!simParams_->haveDampingAlpha()) { |
1255 |
chrisfen |
1045 |
// first set a cutoff dependent alpha value |
1256 |
|
|
// we assume alpha depends linearly with rcut from 0 to 20.5 ang |
1257 |
|
|
alphaVal = 0.5125 - rcut_* 0.025; |
1258 |
|
|
// for values rcut > 20.5, alpha is zero |
1259 |
|
|
if (alphaVal < 0) alphaVal = 0; |
1260 |
|
|
|
1261 |
|
|
// throw warning |
1262 |
chrisfen |
709 |
sprintf( painCave.errMsg, |
1263 |
gezelter |
764 |
"SimInfo warning: dampingAlpha was not specified in the input file.\n" |
1264 |
chrisfen |
1045 |
"\tA default value of %f (1/ang) will be used for the cutoff of\n\t%f (ang).\n", alphaVal, rcut_); |
1265 |
chrisfen |
709 |
painCave.isFatal = 0; |
1266 |
|
|
simError(); |
1267 |
chrisfen |
1089 |
} else { |
1268 |
|
|
alphaVal = simParams_->getDampingAlpha(); |
1269 |
chrisfen |
709 |
} |
1270 |
chrisfen |
1089 |
|
1271 |
chrisfen |
716 |
} else { |
1272 |
|
|
// throw error |
1273 |
|
|
sprintf( painCave.errMsg, |
1274 |
gezelter |
764 |
"SimInfo error: Unknown electrostaticScreeningMethod.\n" |
1275 |
|
|
"\t(Input file specified %s .)\n" |
1276 |
|
|
"\telectrostaticScreeningMethod must be one of: \"undamped\"\n" |
1277 |
|
|
"or \"damped\".\n", myScreen.c_str() ); |
1278 |
chrisfen |
716 |
painCave.isFatal = 1; |
1279 |
|
|
simError(); |
1280 |
chrisfen |
709 |
} |
1281 |
|
|
} |
1282 |
|
|
} |
1283 |
chrisfen |
716 |
|
1284 |
gezelter |
1503 |
|
1285 |
|
|
Electrostatic::setElectrostaticSummationMethod( esm ); |
1286 |
|
|
Electrostatic::setElectrostaticScreeningMethod( sm ); |
1287 |
|
|
Electrostatic::setDampingAlpha( alphaVal ); |
1288 |
|
|
Electrostatic::setReactionFieldDielectric( dielectric ); |
1289 |
gezelter |
764 |
initFortranFF( &errorOut ); |
1290 |
chrisfen |
598 |
} |
1291 |
|
|
|
1292 |
chrisfen |
726 |
void SimInfo::setupSwitchingFunction() { |
1293 |
|
|
int ft = CUBIC; |
1294 |
|
|
|
1295 |
|
|
if (simParams_->haveSwitchingFunctionType()) { |
1296 |
|
|
std::string funcType = simParams_->getSwitchingFunctionType(); |
1297 |
|
|
toUpper(funcType); |
1298 |
|
|
if (funcType == "CUBIC") { |
1299 |
|
|
ft = CUBIC; |
1300 |
|
|
} else { |
1301 |
|
|
if (funcType == "FIFTH_ORDER_POLYNOMIAL") { |
1302 |
|
|
ft = FIFTH_ORDER_POLY; |
1303 |
|
|
} else { |
1304 |
|
|
// throw error |
1305 |
|
|
sprintf( painCave.errMsg, |
1306 |
|
|
"SimInfo error: Unknown switchingFunctionType. (Input file specified %s .)\n\tswitchingFunctionType must be one of: \"cubic\" or \"fifth_order_polynomial\".", funcType.c_str() ); |
1307 |
|
|
painCave.isFatal = 1; |
1308 |
|
|
simError(); |
1309 |
|
|
} |
1310 |
|
|
} |
1311 |
|
|
} |
1312 |
|
|
|
1313 |
|
|
// send switching function notification to switcheroo |
1314 |
|
|
setFunctionType(&ft); |
1315 |
|
|
|
1316 |
|
|
} |
1317 |
|
|
|
1318 |
chrisfen |
998 |
void SimInfo::setupAccumulateBoxDipole() { |
1319 |
|
|
|
1320 |
|
|
// we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true |
1321 |
|
|
if ( simParams_->haveAccumulateBoxDipole() ) |
1322 |
|
|
if ( simParams_->getAccumulateBoxDipole() ) { |
1323 |
|
|
setAccumulateBoxDipole(); |
1324 |
|
|
calcBoxDipole_ = true; |
1325 |
|
|
} |
1326 |
|
|
|
1327 |
|
|
} |
1328 |
|
|
|
1329 |
gezelter |
507 |
void SimInfo::addProperty(GenericData* genData) { |
1330 |
gezelter |
246 |
properties_.addProperty(genData); |
1331 |
gezelter |
507 |
} |
1332 |
gezelter |
2 |
|
1333 |
gezelter |
507 |
void SimInfo::removeProperty(const std::string& propName) { |
1334 |
gezelter |
246 |
properties_.removeProperty(propName); |
1335 |
gezelter |
507 |
} |
1336 |
gezelter |
2 |
|
1337 |
gezelter |
507 |
void SimInfo::clearProperties() { |
1338 |
gezelter |
246 |
properties_.clearProperties(); |
1339 |
gezelter |
507 |
} |
1340 |
gezelter |
2 |
|
1341 |
gezelter |
507 |
std::vector<std::string> SimInfo::getPropertyNames() { |
1342 |
gezelter |
246 |
return properties_.getPropertyNames(); |
1343 |
gezelter |
507 |
} |
1344 |
gezelter |
246 |
|
1345 |
gezelter |
507 |
std::vector<GenericData*> SimInfo::getProperties() { |
1346 |
gezelter |
246 |
return properties_.getProperties(); |
1347 |
gezelter |
507 |
} |
1348 |
gezelter |
2 |
|
1349 |
gezelter |
507 |
GenericData* SimInfo::getPropertyByName(const std::string& propName) { |
1350 |
gezelter |
246 |
return properties_.getPropertyByName(propName); |
1351 |
gezelter |
507 |
} |
1352 |
gezelter |
2 |
|
1353 |
gezelter |
507 |
void SimInfo::setSnapshotManager(SnapshotManager* sman) { |
1354 |
tim |
432 |
if (sman_ == sman) { |
1355 |
gezelter |
507 |
return; |
1356 |
tim |
432 |
} |
1357 |
|
|
delete sman_; |
1358 |
gezelter |
246 |
sman_ = sman; |
1359 |
gezelter |
2 |
|
1360 |
gezelter |
246 |
Molecule* mol; |
1361 |
|
|
RigidBody* rb; |
1362 |
|
|
Atom* atom; |
1363 |
|
|
SimInfo::MoleculeIterator mi; |
1364 |
|
|
Molecule::RigidBodyIterator rbIter; |
1365 |
|
|
Molecule::AtomIterator atomIter;; |
1366 |
|
|
|
1367 |
|
|
for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { |
1368 |
|
|
|
1369 |
gezelter |
507 |
for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) { |
1370 |
|
|
atom->setSnapshotManager(sman_); |
1371 |
|
|
} |
1372 |
gezelter |
246 |
|
1373 |
gezelter |
507 |
for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) { |
1374 |
|
|
rb->setSnapshotManager(sman_); |
1375 |
|
|
} |
1376 |
gezelter |
246 |
} |
1377 |
gezelter |
2 |
|
1378 |
gezelter |
507 |
} |
1379 |
gezelter |
2 |
|
1380 |
gezelter |
507 |
Vector3d SimInfo::getComVel(){ |
1381 |
gezelter |
246 |
SimInfo::MoleculeIterator i; |
1382 |
|
|
Molecule* mol; |
1383 |
gezelter |
2 |
|
1384 |
gezelter |
246 |
Vector3d comVel(0.0); |
1385 |
tim |
963 |
RealType totalMass = 0.0; |
1386 |
gezelter |
2 |
|
1387 |
gezelter |
246 |
|
1388 |
|
|
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
1389 |
tim |
963 |
RealType mass = mol->getMass(); |
1390 |
gezelter |
507 |
totalMass += mass; |
1391 |
|
|
comVel += mass * mol->getComVel(); |
1392 |
gezelter |
246 |
} |
1393 |
gezelter |
2 |
|
1394 |
gezelter |
246 |
#ifdef IS_MPI |
1395 |
tim |
963 |
RealType tmpMass = totalMass; |
1396 |
gezelter |
246 |
Vector3d tmpComVel(comVel); |
1397 |
tim |
963 |
MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); |
1398 |
|
|
MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); |
1399 |
gezelter |
246 |
#endif |
1400 |
|
|
|
1401 |
|
|
comVel /= totalMass; |
1402 |
|
|
|
1403 |
|
|
return comVel; |
1404 |
gezelter |
507 |
} |
1405 |
gezelter |
2 |
|
1406 |
gezelter |
507 |
Vector3d SimInfo::getCom(){ |
1407 |
gezelter |
246 |
SimInfo::MoleculeIterator i; |
1408 |
|
|
Molecule* mol; |
1409 |
gezelter |
2 |
|
1410 |
gezelter |
246 |
Vector3d com(0.0); |
1411 |
tim |
963 |
RealType totalMass = 0.0; |
1412 |
gezelter |
246 |
|
1413 |
|
|
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
1414 |
tim |
963 |
RealType mass = mol->getMass(); |
1415 |
gezelter |
507 |
totalMass += mass; |
1416 |
|
|
com += mass * mol->getCom(); |
1417 |
gezelter |
246 |
} |
1418 |
gezelter |
2 |
|
1419 |
|
|
#ifdef IS_MPI |
1420 |
tim |
963 |
RealType tmpMass = totalMass; |
1421 |
gezelter |
246 |
Vector3d tmpCom(com); |
1422 |
tim |
963 |
MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); |
1423 |
|
|
MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); |
1424 |
gezelter |
2 |
#endif |
1425 |
|
|
|
1426 |
gezelter |
246 |
com /= totalMass; |
1427 |
gezelter |
2 |
|
1428 |
gezelter |
246 |
return com; |
1429 |
gezelter |
2 |
|
1430 |
gezelter |
507 |
} |
1431 |
gezelter |
246 |
|
1432 |
gezelter |
507 |
std::ostream& operator <<(std::ostream& o, SimInfo& info) { |
1433 |
gezelter |
246 |
|
1434 |
|
|
return o; |
1435 |
gezelter |
507 |
} |
1436 |
chuckv |
555 |
|
1437 |
|
|
|
1438 |
|
|
/* |
1439 |
|
|
Returns center of mass and center of mass velocity in one function call. |
1440 |
|
|
*/ |
1441 |
|
|
|
1442 |
|
|
void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){ |
1443 |
|
|
SimInfo::MoleculeIterator i; |
1444 |
|
|
Molecule* mol; |
1445 |
|
|
|
1446 |
|
|
|
1447 |
tim |
963 |
RealType totalMass = 0.0; |
1448 |
chuckv |
555 |
|
1449 |
gezelter |
246 |
|
1450 |
chuckv |
555 |
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
1451 |
tim |
963 |
RealType mass = mol->getMass(); |
1452 |
chuckv |
555 |
totalMass += mass; |
1453 |
|
|
com += mass * mol->getCom(); |
1454 |
|
|
comVel += mass * mol->getComVel(); |
1455 |
|
|
} |
1456 |
|
|
|
1457 |
|
|
#ifdef IS_MPI |
1458 |
tim |
963 |
RealType tmpMass = totalMass; |
1459 |
chuckv |
555 |
Vector3d tmpCom(com); |
1460 |
|
|
Vector3d tmpComVel(comVel); |
1461 |
tim |
963 |
MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); |
1462 |
|
|
MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); |
1463 |
|
|
MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); |
1464 |
chuckv |
555 |
#endif |
1465 |
|
|
|
1466 |
|
|
com /= totalMass; |
1467 |
|
|
comVel /= totalMass; |
1468 |
|
|
} |
1469 |
|
|
|
1470 |
|
|
/* |
1471 |
|
|
Return intertia tensor for entire system and angular momentum Vector. |
1472 |
chuckv |
557 |
|
1473 |
|
|
|
1474 |
|
|
[ Ixx -Ixy -Ixz ] |
1475 |
|
|
J =| -Iyx Iyy -Iyz | |
1476 |
|
|
[ -Izx -Iyz Izz ] |
1477 |
chuckv |
555 |
*/ |
1478 |
|
|
|
1479 |
|
|
void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){ |
1480 |
|
|
|
1481 |
|
|
|
1482 |
tim |
963 |
RealType xx = 0.0; |
1483 |
|
|
RealType yy = 0.0; |
1484 |
|
|
RealType zz = 0.0; |
1485 |
|
|
RealType xy = 0.0; |
1486 |
|
|
RealType xz = 0.0; |
1487 |
|
|
RealType yz = 0.0; |
1488 |
chuckv |
555 |
Vector3d com(0.0); |
1489 |
|
|
Vector3d comVel(0.0); |
1490 |
|
|
|
1491 |
|
|
getComAll(com, comVel); |
1492 |
|
|
|
1493 |
|
|
SimInfo::MoleculeIterator i; |
1494 |
|
|
Molecule* mol; |
1495 |
|
|
|
1496 |
|
|
Vector3d thisq(0.0); |
1497 |
|
|
Vector3d thisv(0.0); |
1498 |
|
|
|
1499 |
tim |
963 |
RealType thisMass = 0.0; |
1500 |
chuckv |
555 |
|
1501 |
|
|
|
1502 |
|
|
|
1503 |
|
|
|
1504 |
|
|
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
1505 |
|
|
|
1506 |
|
|
thisq = mol->getCom()-com; |
1507 |
|
|
thisv = mol->getComVel()-comVel; |
1508 |
|
|
thisMass = mol->getMass(); |
1509 |
|
|
// Compute moment of intertia coefficients. |
1510 |
|
|
xx += thisq[0]*thisq[0]*thisMass; |
1511 |
|
|
yy += thisq[1]*thisq[1]*thisMass; |
1512 |
|
|
zz += thisq[2]*thisq[2]*thisMass; |
1513 |
|
|
|
1514 |
|
|
// compute products of intertia |
1515 |
|
|
xy += thisq[0]*thisq[1]*thisMass; |
1516 |
|
|
xz += thisq[0]*thisq[2]*thisMass; |
1517 |
|
|
yz += thisq[1]*thisq[2]*thisMass; |
1518 |
|
|
|
1519 |
|
|
angularMomentum += cross( thisq, thisv ) * thisMass; |
1520 |
|
|
|
1521 |
|
|
} |
1522 |
|
|
|
1523 |
|
|
|
1524 |
|
|
inertiaTensor(0,0) = yy + zz; |
1525 |
|
|
inertiaTensor(0,1) = -xy; |
1526 |
|
|
inertiaTensor(0,2) = -xz; |
1527 |
|
|
inertiaTensor(1,0) = -xy; |
1528 |
chuckv |
557 |
inertiaTensor(1,1) = xx + zz; |
1529 |
chuckv |
555 |
inertiaTensor(1,2) = -yz; |
1530 |
|
|
inertiaTensor(2,0) = -xz; |
1531 |
|
|
inertiaTensor(2,1) = -yz; |
1532 |
|
|
inertiaTensor(2,2) = xx + yy; |
1533 |
|
|
|
1534 |
|
|
#ifdef IS_MPI |
1535 |
|
|
Mat3x3d tmpI(inertiaTensor); |
1536 |
|
|
Vector3d tmpAngMom; |
1537 |
tim |
963 |
MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); |
1538 |
|
|
MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); |
1539 |
chuckv |
555 |
#endif |
1540 |
|
|
|
1541 |
|
|
return; |
1542 |
|
|
} |
1543 |
|
|
|
1544 |
|
|
//Returns the angular momentum of the system |
1545 |
|
|
Vector3d SimInfo::getAngularMomentum(){ |
1546 |
|
|
|
1547 |
|
|
Vector3d com(0.0); |
1548 |
|
|
Vector3d comVel(0.0); |
1549 |
|
|
Vector3d angularMomentum(0.0); |
1550 |
|
|
|
1551 |
|
|
getComAll(com,comVel); |
1552 |
|
|
|
1553 |
|
|
SimInfo::MoleculeIterator i; |
1554 |
|
|
Molecule* mol; |
1555 |
|
|
|
1556 |
chuckv |
557 |
Vector3d thisr(0.0); |
1557 |
|
|
Vector3d thisp(0.0); |
1558 |
chuckv |
555 |
|
1559 |
tim |
963 |
RealType thisMass; |
1560 |
chuckv |
555 |
|
1561 |
|
|
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
1562 |
chuckv |
557 |
thisMass = mol->getMass(); |
1563 |
|
|
thisr = mol->getCom()-com; |
1564 |
|
|
thisp = (mol->getComVel()-comVel)*thisMass; |
1565 |
chuckv |
555 |
|
1566 |
chuckv |
557 |
angularMomentum += cross( thisr, thisp ); |
1567 |
|
|
|
1568 |
chuckv |
555 |
} |
1569 |
|
|
|
1570 |
|
|
#ifdef IS_MPI |
1571 |
|
|
Vector3d tmpAngMom; |
1572 |
tim |
963 |
MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); |
1573 |
chuckv |
555 |
#endif |
1574 |
|
|
|
1575 |
|
|
return angularMomentum; |
1576 |
|
|
} |
1577 |
|
|
|
1578 |
tim |
1024 |
StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) { |
1579 |
|
|
return IOIndexToIntegrableObject.at(index); |
1580 |
|
|
} |
1581 |
|
|
|
1582 |
|
|
void SimInfo::setIOIndexToIntegrableObject(const std::vector<StuntDouble*>& v) { |
1583 |
|
|
IOIndexToIntegrableObject= v; |
1584 |
|
|
} |
1585 |
|
|
|
1586 |
chuckv |
1103 |
/* Returns the Volume of the simulation based on a ellipsoid with semi-axes |
1587 |
|
|
based on the radius of gyration V=4/3*Pi*R_1*R_2*R_3 |
1588 |
|
|
where R_i are related to the principle inertia moments R_i = sqrt(C*I_i/N), this reduces to |
1589 |
|
|
V = 4/3*Pi*(C/N)^3/2*sqrt(det(I)). See S.E. Baltazar et. al. Comp. Mat. Sci. 37 (2006) 526-536. |
1590 |
|
|
*/ |
1591 |
|
|
void SimInfo::getGyrationalVolume(RealType &volume){ |
1592 |
|
|
Mat3x3d intTensor; |
1593 |
|
|
RealType det; |
1594 |
|
|
Vector3d dummyAngMom; |
1595 |
|
|
RealType sysconstants; |
1596 |
|
|
RealType geomCnst; |
1597 |
|
|
|
1598 |
|
|
geomCnst = 3.0/2.0; |
1599 |
|
|
/* Get the inertial tensor and angular momentum for free*/ |
1600 |
|
|
getInertiaTensor(intTensor,dummyAngMom); |
1601 |
|
|
|
1602 |
|
|
det = intTensor.determinant(); |
1603 |
|
|
sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_; |
1604 |
|
|
volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(det); |
1605 |
|
|
return; |
1606 |
|
|
} |
1607 |
|
|
|
1608 |
|
|
void SimInfo::getGyrationalVolume(RealType &volume, RealType &detI){ |
1609 |
|
|
Mat3x3d intTensor; |
1610 |
|
|
Vector3d dummyAngMom; |
1611 |
|
|
RealType sysconstants; |
1612 |
|
|
RealType geomCnst; |
1613 |
|
|
|
1614 |
|
|
geomCnst = 3.0/2.0; |
1615 |
|
|
/* Get the inertial tensor and angular momentum for free*/ |
1616 |
|
|
getInertiaTensor(intTensor,dummyAngMom); |
1617 |
|
|
|
1618 |
|
|
detI = intTensor.determinant(); |
1619 |
|
|
sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_; |
1620 |
|
|
volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(detI); |
1621 |
|
|
return; |
1622 |
|
|
} |
1623 |
tim |
1024 |
/* |
1624 |
|
|
void SimInfo::setStuntDoubleFromGlobalIndex(std::vector<StuntDouble*> v) { |
1625 |
|
|
assert( v.size() == nAtoms_ + nRigidBodies_); |
1626 |
|
|
sdByGlobalIndex_ = v; |
1627 |
|
|
} |
1628 |
|
|
|
1629 |
|
|
StuntDouble* SimInfo::getStuntDoubleFromGlobalIndex(int index) { |
1630 |
|
|
//assert(index < nAtoms_ + nRigidBodies_); |
1631 |
|
|
return sdByGlobalIndex_.at(index); |
1632 |
|
|
} |
1633 |
|
|
*/ |
1634 |
gezelter |
1390 |
}//end namespace OpenMD |
1635 |
gezelter |
246 |
|