ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/Thermo.cpp
Revision: 1709
Committed: Tue May 15 13:04:08 2012 UTC (12 years, 11 months ago) by gezelter
File size: 14402 byte(s)
Log Message:
Moving silly stuff out of Stats and into Snapshot.  Most of it should go
into a not-yet-implemented FrameData class.

File Contents

# Content
1 /*
2 * 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 * 1. Redistributions of source code must retain the above copyright
10 * notice, this list of conditions and the following disclaimer.
11 *
12 * 2. Redistributions in binary form must reproduce the above copyright
13 * 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 *
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] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40 * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 */
42
43 #include <math.h>
44 #include <iostream>
45
46 #ifdef IS_MPI
47 #include <mpi.h>
48 #endif //is_mpi
49
50 #include "brains/Thermo.hpp"
51 #include "primitives/Molecule.hpp"
52 #include "utils/simError.h"
53 #include "utils/PhysicalConstants.hpp"
54
55 namespace OpenMD {
56
57 RealType Thermo::getKinetic() {
58 SimInfo::MoleculeIterator miter;
59 std::vector<StuntDouble*>::iterator iiter;
60 Molecule* mol;
61 StuntDouble* integrableObject;
62 Vector3d vel;
63 Vector3d angMom;
64 Mat3x3d I;
65 int i;
66 int j;
67 int k;
68 RealType mass;
69 RealType kinetic = 0.0;
70 RealType kinetic_global = 0.0;
71
72 for (mol = info_->beginMolecule(miter); mol != NULL; mol = info_->nextMolecule(miter)) {
73 for (integrableObject = mol->beginIntegrableObject(iiter); integrableObject != NULL;
74 integrableObject = mol->nextIntegrableObject(iiter)) {
75
76 mass = integrableObject->getMass();
77 vel = integrableObject->getVel();
78
79 kinetic += mass * (vel[0]*vel[0] + vel[1]*vel[1] + vel[2]*vel[2]);
80
81 if (integrableObject->isDirectional()) {
82 angMom = integrableObject->getJ();
83 I = integrableObject->getI();
84
85 if (integrableObject->isLinear()) {
86 i = integrableObject->linearAxis();
87 j = (i + 1) % 3;
88 k = (i + 2) % 3;
89 kinetic += angMom[j] * angMom[j] / I(j, j) + angMom[k] * angMom[k] / I(k, k);
90 } else {
91 kinetic += angMom[0]*angMom[0]/I(0, 0) + angMom[1]*angMom[1]/I(1, 1)
92 + angMom[2]*angMom[2]/I(2, 2);
93 }
94 }
95
96 }
97 }
98
99 #ifdef IS_MPI
100
101 MPI_Allreduce(&kinetic, &kinetic_global, 1, MPI_REALTYPE, MPI_SUM,
102 MPI_COMM_WORLD);
103 kinetic = kinetic_global;
104
105 #endif //is_mpi
106
107 kinetic = kinetic * 0.5 / PhysicalConstants::energyConvert;
108
109 return kinetic;
110 }
111
112 RealType Thermo::getPotential() {
113 RealType potential = 0.0;
114 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
115 RealType shortRangePot_local = curSnapshot->statData[Stats::SHORT_RANGE_POTENTIAL] ;
116
117 // Get total potential for entire system from MPI.
118
119 #ifdef IS_MPI
120
121 MPI_Allreduce(&shortRangePot_local, &potential, 1, MPI_REALTYPE, MPI_SUM,
122 MPI_COMM_WORLD);
123 potential += curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL];
124
125 #else
126
127 potential = shortRangePot_local + curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL];
128
129 #endif // is_mpi
130
131 return potential;
132 }
133
134 RealType Thermo::getTotalE() {
135 RealType total;
136
137 total = this->getKinetic() + this->getPotential();
138 return total;
139 }
140
141 RealType Thermo::getTemperature() {
142
143 RealType temperature = ( 2.0 * this->getKinetic() ) / (info_->getNdf()* PhysicalConstants::kb );
144 return temperature;
145 }
146
147 RealType Thermo::getVolume() {
148 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
149 return curSnapshot->getVolume();
150 }
151
152 RealType Thermo::getPressure() {
153
154 // Relies on the calculation of the full molecular pressure tensor
155
156
157 Mat3x3d tensor;
158 RealType pressure;
159
160 tensor = getPressureTensor();
161
162 pressure = PhysicalConstants::pressureConvert * (tensor(0, 0) + tensor(1, 1) + tensor(2, 2)) / 3.0;
163
164 return pressure;
165 }
166
167 RealType Thermo::getPressure(int direction) {
168
169 // Relies on the calculation of the full molecular pressure tensor
170
171
172 Mat3x3d tensor;
173 RealType pressure;
174
175 tensor = getPressureTensor();
176
177 pressure = PhysicalConstants::pressureConvert * tensor(direction, direction);
178
179 return pressure;
180 }
181
182 Mat3x3d Thermo::getPressureTensor() {
183 // returns pressure tensor in units amu*fs^-2*Ang^-1
184 // routine derived via viral theorem description in:
185 // Paci, E. and Marchi, M. J.Phys.Chem. 1996, 100, 4314-4322
186 Mat3x3d pressureTensor;
187 Mat3x3d p_local(0.0);
188 Mat3x3d p_global(0.0);
189
190 SimInfo::MoleculeIterator i;
191 std::vector<StuntDouble*>::iterator j;
192 Molecule* mol;
193 StuntDouble* integrableObject;
194 for (mol = info_->beginMolecule(i); mol != NULL; mol = info_->nextMolecule(i)) {
195 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
196 integrableObject = mol->nextIntegrableObject(j)) {
197
198 RealType mass = integrableObject->getMass();
199 Vector3d vcom = integrableObject->getVel();
200 p_local += mass * outProduct(vcom, vcom);
201 }
202 }
203
204 #ifdef IS_MPI
205 MPI_Allreduce(p_local.getArrayPointer(), p_global.getArrayPointer(), 9, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
206 #else
207 p_global = p_local;
208 #endif // is_mpi
209
210 RealType volume = this->getVolume();
211 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
212 Mat3x3d tau = curSnapshot->getTau();
213
214 pressureTensor = (p_global + PhysicalConstants::energyConvert* tau)/volume;
215
216 return pressureTensor;
217 }
218
219
220 void Thermo::saveStat(){
221 Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
222 Stats& stat = currSnapshot->statData;
223
224 stat[Stats::KINETIC_ENERGY] = getKinetic();
225 stat[Stats::POTENTIAL_ENERGY] = getPotential();
226 stat[Stats::TOTAL_ENERGY] = stat[Stats::KINETIC_ENERGY] + stat[Stats::POTENTIAL_ENERGY] ;
227 stat[Stats::TEMPERATURE] = getTemperature();
228 stat[Stats::PRESSURE] = getPressure();
229 stat[Stats::VOLUME] = getVolume();
230
231 Mat3x3d tensor =getPressureTensor();
232 stat[Stats::PRESSURE_TENSOR_XX] = tensor(0, 0);
233 stat[Stats::PRESSURE_TENSOR_XY] = tensor(0, 1);
234 stat[Stats::PRESSURE_TENSOR_XZ] = tensor(0, 2);
235 stat[Stats::PRESSURE_TENSOR_YX] = tensor(1, 0);
236 stat[Stats::PRESSURE_TENSOR_YY] = tensor(1, 1);
237 stat[Stats::PRESSURE_TENSOR_YZ] = tensor(1, 2);
238 stat[Stats::PRESSURE_TENSOR_ZX] = tensor(2, 0);
239 stat[Stats::PRESSURE_TENSOR_ZY] = tensor(2, 1);
240 stat[Stats::PRESSURE_TENSOR_ZZ] = tensor(2, 2);
241
242 // grab the simulation box dipole moment if specified
243 if (info_->getCalcBoxDipole()){
244 Vector3d totalDipole = getBoxDipole();
245 stat[Stats::BOX_DIPOLE_X] = totalDipole(0);
246 stat[Stats::BOX_DIPOLE_Y] = totalDipole(1);
247 stat[Stats::BOX_DIPOLE_Z] = totalDipole(2);
248 }
249
250 Globals* simParams = info_->getSimParams();
251
252 if (simParams->haveTaggedAtomPair() &&
253 simParams->havePrintTaggedPairDistance()) {
254 if ( simParams->getPrintTaggedPairDistance()) {
255
256 std::pair<int, int> tap = simParams->getTaggedAtomPair();
257 Vector3d pos1, pos2, rab;
258
259 #ifdef IS_MPI
260 std::cerr << "tap = " << tap.first << " " << tap.second << std::endl;
261
262 int mol1 = info_->getGlobalMolMembership(tap.first);
263 int mol2 = info_->getGlobalMolMembership(tap.second);
264 std::cerr << "mols = " << mol1 << " " << mol2 << std::endl;
265
266 int proc1 = info_->getMolToProc(mol1);
267 int proc2 = info_->getMolToProc(mol2);
268
269 std::cerr << " procs = " << proc1 << " " <<proc2 <<std::endl;
270
271 RealType data[3];
272 if (proc1 == worldRank) {
273 StuntDouble* sd1 = info_->getIOIndexToIntegrableObject(tap.first);
274 std::cerr << " on proc " << proc1 << ", sd1 has global index= " << sd1->getGlobalIndex() << std::endl;
275 pos1 = sd1->getPos();
276 data[0] = pos1.x();
277 data[1] = pos1.y();
278 data[2] = pos1.z();
279 MPI_Bcast(data, 3, MPI_REALTYPE, proc1, MPI_COMM_WORLD);
280 } else {
281 MPI_Bcast(data, 3, MPI_REALTYPE, proc1, MPI_COMM_WORLD);
282 pos1 = Vector3d(data);
283 }
284
285
286 if (proc2 == worldRank) {
287 StuntDouble* sd2 = info_->getIOIndexToIntegrableObject(tap.second);
288 std::cerr << " on proc " << proc2 << ", sd2 has global index= " << sd2->getGlobalIndex() << std::endl;
289 pos2 = sd2->getPos();
290 data[0] = pos2.x();
291 data[1] = pos2.y();
292 data[2] = pos2.z();
293 MPI_Bcast(data, 3, MPI_REALTYPE, proc2, MPI_COMM_WORLD);
294 } else {
295 MPI_Bcast(data, 3, MPI_REALTYPE, proc2, MPI_COMM_WORLD);
296 pos2 = Vector3d(data);
297 }
298 #else
299 StuntDouble* at1 = info_->getIOIndexToIntegrableObject(tap.first);
300 StuntDouble* at2 = info_->getIOIndexToIntegrableObject(tap.second);
301 pos1 = at1->getPos();
302 pos2 = at2->getPos();
303 #endif
304 rab = pos2 - pos1;
305 currSnapshot->wrapVector(rab);
306 stat[Stats::TAGGED_PAIR_DISTANCE] = rab.length();
307 }
308 }
309
310 /**@todo need refactorying*/
311 //Conserved Quantity is set by integrator and time is set by setTime
312
313 }
314
315
316 Vector3d Thermo::getBoxDipole() {
317 Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
318 SimInfo::MoleculeIterator miter;
319 std::vector<Atom*>::iterator aiter;
320 Molecule* mol;
321 Atom* atom;
322 RealType charge;
323 RealType moment(0.0);
324 Vector3d ri(0.0);
325 Vector3d dipoleVector(0.0);
326 Vector3d nPos(0.0);
327 Vector3d pPos(0.0);
328 RealType nChg(0.0);
329 RealType pChg(0.0);
330 int nCount = 0;
331 int pCount = 0;
332
333 RealType chargeToC = 1.60217733e-19;
334 RealType angstromToM = 1.0e-10;
335 RealType debyeToCm = 3.33564095198e-30;
336
337 for (mol = info_->beginMolecule(miter); mol != NULL;
338 mol = info_->nextMolecule(miter)) {
339
340 for (atom = mol->beginAtom(aiter); atom != NULL;
341 atom = mol->nextAtom(aiter)) {
342
343 if (atom->isCharge() ) {
344 charge = 0.0;
345 GenericData* data = atom->getAtomType()->getPropertyByName("Charge");
346 if (data != NULL) {
347
348 charge = (dynamic_cast<DoubleGenericData*>(data))->getData();
349 charge *= chargeToC;
350
351 ri = atom->getPos();
352 currSnapshot->wrapVector(ri);
353 ri *= angstromToM;
354
355 if (charge < 0.0) {
356 nPos += ri;
357 nChg -= charge;
358 nCount++;
359 } else if (charge > 0.0) {
360 pPos += ri;
361 pChg += charge;
362 pCount++;
363 }
364 }
365 }
366
367 if (atom->isDipole() ) {
368 Vector3d u_i = atom->getElectroFrame().getColumn(2);
369 GenericData* data = dynamic_cast<DirectionalAtomType*>(atom->getAtomType())->getPropertyByName("Dipole");
370 if (data != NULL) {
371 moment = (dynamic_cast<DoubleGenericData*>(data))->getData();
372
373 moment *= debyeToCm;
374 dipoleVector += u_i * moment;
375 }
376 }
377 }
378 }
379
380
381 #ifdef IS_MPI
382 RealType pChg_global, nChg_global;
383 int pCount_global, nCount_global;
384 Vector3d pPos_global, nPos_global, dipVec_global;
385
386 MPI_Allreduce(&pChg, &pChg_global, 1, MPI_REALTYPE, MPI_SUM,
387 MPI_COMM_WORLD);
388 pChg = pChg_global;
389 MPI_Allreduce(&nChg, &nChg_global, 1, MPI_REALTYPE, MPI_SUM,
390 MPI_COMM_WORLD);
391 nChg = nChg_global;
392 MPI_Allreduce(&pCount, &pCount_global, 1, MPI_INTEGER, MPI_SUM,
393 MPI_COMM_WORLD);
394 pCount = pCount_global;
395 MPI_Allreduce(&nCount, &nCount_global, 1, MPI_INTEGER, MPI_SUM,
396 MPI_COMM_WORLD);
397 nCount = nCount_global;
398 MPI_Allreduce(pPos.getArrayPointer(), pPos_global.getArrayPointer(), 3,
399 MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
400 pPos = pPos_global;
401 MPI_Allreduce(nPos.getArrayPointer(), nPos_global.getArrayPointer(), 3,
402 MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
403 nPos = nPos_global;
404 MPI_Allreduce(dipoleVector.getArrayPointer(),
405 dipVec_global.getArrayPointer(), 3,
406 MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
407 dipoleVector = dipVec_global;
408 #endif //is_mpi
409
410 // first load the accumulated dipole moment (if dipoles were present)
411 Vector3d boxDipole = dipoleVector;
412 // now include the dipole moment due to charges
413 // use the lesser of the positive and negative charge totals
414 RealType chg_value = nChg <= pChg ? nChg : pChg;
415
416 // find the average positions
417 if (pCount > 0 && nCount > 0 ) {
418 pPos /= pCount;
419 nPos /= nCount;
420 }
421
422 // dipole is from the negative to the positive (physics notation)
423 boxDipole += (pPos - nPos) * chg_value;
424
425 return boxDipole;
426 }
427 } //end namespace OpenMD

Properties

Name Value
svn:keywords Author Id Revision Date