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] Vardeman & Gezelter, in progress (2009). |
40 |
*/ |
41 |
|
42 |
#include <math.h> |
43 |
#include <iostream> |
44 |
|
45 |
#ifdef IS_MPI |
46 |
#include <mpi.h> |
47 |
#endif //is_mpi |
48 |
|
49 |
#include "brains/Thermo.hpp" |
50 |
#include "primitives/Molecule.hpp" |
51 |
#include "utils/simError.h" |
52 |
#include "utils/PhysicalConstants.hpp" |
53 |
|
54 |
namespace OpenMD { |
55 |
|
56 |
RealType Thermo::getKinetic() { |
57 |
SimInfo::MoleculeIterator miter; |
58 |
std::vector<StuntDouble*>::iterator iiter; |
59 |
Molecule* mol; |
60 |
StuntDouble* integrableObject; |
61 |
Vector3d vel; |
62 |
Vector3d angMom; |
63 |
Mat3x3d I; |
64 |
int i; |
65 |
int j; |
66 |
int k; |
67 |
RealType mass; |
68 |
RealType kinetic = 0.0; |
69 |
RealType kinetic_global = 0.0; |
70 |
|
71 |
for (mol = info_->beginMolecule(miter); mol != NULL; mol = info_->nextMolecule(miter)) { |
72 |
for (integrableObject = mol->beginIntegrableObject(iiter); integrableObject != NULL; |
73 |
integrableObject = mol->nextIntegrableObject(iiter)) { |
74 |
|
75 |
mass = integrableObject->getMass(); |
76 |
vel = integrableObject->getVel(); |
77 |
|
78 |
kinetic += mass * (vel[0]*vel[0] + vel[1]*vel[1] + vel[2]*vel[2]); |
79 |
|
80 |
if (integrableObject->isDirectional()) { |
81 |
angMom = integrableObject->getJ(); |
82 |
I = integrableObject->getI(); |
83 |
|
84 |
if (integrableObject->isLinear()) { |
85 |
i = integrableObject->linearAxis(); |
86 |
j = (i + 1) % 3; |
87 |
k = (i + 2) % 3; |
88 |
kinetic += angMom[j] * angMom[j] / I(j, j) + angMom[k] * angMom[k] / I(k, k); |
89 |
} else { |
90 |
kinetic += angMom[0]*angMom[0]/I(0, 0) + angMom[1]*angMom[1]/I(1, 1) |
91 |
+ angMom[2]*angMom[2]/I(2, 2); |
92 |
} |
93 |
} |
94 |
|
95 |
} |
96 |
} |
97 |
|
98 |
#ifdef IS_MPI |
99 |
|
100 |
MPI_Allreduce(&kinetic, &kinetic_global, 1, MPI_REALTYPE, MPI_SUM, |
101 |
MPI_COMM_WORLD); |
102 |
kinetic = kinetic_global; |
103 |
|
104 |
#endif //is_mpi |
105 |
|
106 |
kinetic = kinetic * 0.5 / PhysicalConstants::energyConvert; |
107 |
|
108 |
return kinetic; |
109 |
} |
110 |
|
111 |
RealType Thermo::getPotential() { |
112 |
RealType potential = 0.0; |
113 |
Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot(); |
114 |
RealType shortRangePot_local = curSnapshot->statData[Stats::SHORT_RANGE_POTENTIAL] ; |
115 |
|
116 |
// Get total potential for entire system from MPI. |
117 |
|
118 |
#ifdef IS_MPI |
119 |
|
120 |
MPI_Allreduce(&shortRangePot_local, &potential, 1, MPI_REALTYPE, MPI_SUM, |
121 |
MPI_COMM_WORLD); |
122 |
potential += curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL]; |
123 |
|
124 |
#else |
125 |
|
126 |
potential = shortRangePot_local + curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL]; |
127 |
|
128 |
#endif // is_mpi |
129 |
|
130 |
return potential; |
131 |
} |
132 |
|
133 |
RealType Thermo::getTotalE() { |
134 |
RealType total; |
135 |
|
136 |
total = this->getKinetic() + this->getPotential(); |
137 |
return total; |
138 |
} |
139 |
|
140 |
RealType Thermo::getTemperature() { |
141 |
|
142 |
RealType temperature = ( 2.0 * this->getKinetic() ) / (info_->getNdf()* PhysicalConstants::kb ); |
143 |
return temperature; |
144 |
} |
145 |
|
146 |
RealType Thermo::getVolume() { |
147 |
Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot(); |
148 |
return curSnapshot->getVolume(); |
149 |
} |
150 |
|
151 |
RealType Thermo::getPressure() { |
152 |
|
153 |
// Relies on the calculation of the full molecular pressure tensor |
154 |
|
155 |
|
156 |
Mat3x3d tensor; |
157 |
RealType pressure; |
158 |
|
159 |
tensor = getPressureTensor(); |
160 |
|
161 |
pressure = PhysicalConstants::pressureConvert * (tensor(0, 0) + tensor(1, 1) + tensor(2, 2)) / 3.0; |
162 |
|
163 |
return pressure; |
164 |
} |
165 |
|
166 |
RealType Thermo::getPressure(int direction) { |
167 |
|
168 |
// Relies on the calculation of the full molecular pressure tensor |
169 |
|
170 |
|
171 |
Mat3x3d tensor; |
172 |
RealType pressure; |
173 |
|
174 |
tensor = getPressureTensor(); |
175 |
|
176 |
pressure = PhysicalConstants::pressureConvert * tensor(direction, direction); |
177 |
|
178 |
return pressure; |
179 |
} |
180 |
|
181 |
Mat3x3d Thermo::getPressureTensor() { |
182 |
// returns pressure tensor in units amu*fs^-2*Ang^-1 |
183 |
// routine derived via viral theorem description in: |
184 |
// Paci, E. and Marchi, M. J.Phys.Chem. 1996, 100, 4314-4322 |
185 |
Mat3x3d pressureTensor; |
186 |
Mat3x3d p_local(0.0); |
187 |
Mat3x3d p_global(0.0); |
188 |
|
189 |
SimInfo::MoleculeIterator i; |
190 |
std::vector<StuntDouble*>::iterator j; |
191 |
Molecule* mol; |
192 |
StuntDouble* integrableObject; |
193 |
for (mol = info_->beginMolecule(i); mol != NULL; mol = info_->nextMolecule(i)) { |
194 |
for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL; |
195 |
integrableObject = mol->nextIntegrableObject(j)) { |
196 |
|
197 |
RealType mass = integrableObject->getMass(); |
198 |
Vector3d vcom = integrableObject->getVel(); |
199 |
p_local += mass * outProduct(vcom, vcom); |
200 |
} |
201 |
} |
202 |
|
203 |
#ifdef IS_MPI |
204 |
MPI_Allreduce(p_local.getArrayPointer(), p_global.getArrayPointer(), 9, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD); |
205 |
#else |
206 |
p_global = p_local; |
207 |
#endif // is_mpi |
208 |
|
209 |
RealType volume = this->getVolume(); |
210 |
Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot(); |
211 |
Mat3x3d tau = curSnapshot->statData.getTau(); |
212 |
|
213 |
pressureTensor = (p_global + PhysicalConstants::energyConvert* tau)/volume; |
214 |
|
215 |
return pressureTensor; |
216 |
} |
217 |
|
218 |
|
219 |
void Thermo::saveStat(){ |
220 |
Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot(); |
221 |
Stats& stat = currSnapshot->statData; |
222 |
|
223 |
stat[Stats::KINETIC_ENERGY] = getKinetic(); |
224 |
stat[Stats::POTENTIAL_ENERGY] = getPotential(); |
225 |
stat[Stats::TOTAL_ENERGY] = stat[Stats::KINETIC_ENERGY] + stat[Stats::POTENTIAL_ENERGY] ; |
226 |
stat[Stats::TEMPERATURE] = getTemperature(); |
227 |
stat[Stats::PRESSURE] = getPressure(); |
228 |
stat[Stats::VOLUME] = getVolume(); |
229 |
|
230 |
Mat3x3d tensor =getPressureTensor(); |
231 |
stat[Stats::PRESSURE_TENSOR_XX] = tensor(0, 0); |
232 |
stat[Stats::PRESSURE_TENSOR_XY] = tensor(0, 1); |
233 |
stat[Stats::PRESSURE_TENSOR_XZ] = tensor(0, 2); |
234 |
stat[Stats::PRESSURE_TENSOR_YX] = tensor(1, 0); |
235 |
stat[Stats::PRESSURE_TENSOR_YY] = tensor(1, 1); |
236 |
stat[Stats::PRESSURE_TENSOR_YZ] = tensor(1, 2); |
237 |
stat[Stats::PRESSURE_TENSOR_ZX] = tensor(2, 0); |
238 |
stat[Stats::PRESSURE_TENSOR_ZY] = tensor(2, 1); |
239 |
stat[Stats::PRESSURE_TENSOR_ZZ] = tensor(2, 2); |
240 |
Vector3d GKappa_t = getThermalHelfand(); |
241 |
stat[Stats::THERMAL_HELFANDMOMENT_X] = GKappa_t.x(); |
242 |
stat[Stats::THERMAL_HELFANDMOMENT_Y] = GKappa_t.y(); |
243 |
stat[Stats::THERMAL_HELFANDMOMENT_Z] = GKappa_t.z(); |
244 |
Vector3d HeatFlux_J = getHeatFlux(); |
245 |
stat[Stats::HEATFLUX_X] = HeatFlux_J.x(); |
246 |
stat[Stats::HEATFLUX_Y] = HeatFlux_J.y(); |
247 |
stat[Stats::HEATFLUX_Z] = HeatFlux_J.z(); |
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 |
|
317 |
Vector3d Thermo::getBoxDipole() { |
318 |
Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot(); |
319 |
SimInfo::MoleculeIterator miter; |
320 |
std::vector<Atom*>::iterator aiter; |
321 |
Molecule* mol; |
322 |
Atom* atom; |
323 |
RealType charge; |
324 |
RealType moment(0.0); |
325 |
Vector3d ri(0.0); |
326 |
Vector3d dipoleVector(0.0); |
327 |
Vector3d nPos(0.0); |
328 |
Vector3d pPos(0.0); |
329 |
RealType nChg(0.0); |
330 |
RealType pChg(0.0); |
331 |
int nCount = 0; |
332 |
int pCount = 0; |
333 |
|
334 |
RealType chargeToC = 1.60217733e-19; |
335 |
RealType angstromToM = 1.0e-10; 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 |
|
428 |
Vector3d Thermo::getThermalHelfand() { |
429 |
Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot(); |
430 |
SimInfo::MoleculeIterator miter; |
431 |
std::vector<Atom*>::iterator aiter; |
432 |
Molecule* mol; |
433 |
Atom* atom; |
434 |
RealType mass; |
435 |
Vector3d velocity; |
436 |
Vector3d x_a; |
437 |
RealType kinetic; |
438 |
RealType potential; |
439 |
RealType eatom; |
440 |
RealType AvgE_a_ = 0; |
441 |
Vector3d GKappa_t = V3Zero; |
442 |
Vector3d ThermalHelfandMoment; |
443 |
|
444 |
for (mol = info_->beginMolecule(miter); mol != NULL; |
445 |
mol = info_->nextMolecule(miter)) { |
446 |
|
447 |
for (atom = mol->beginAtom(aiter); atom != NULL; |
448 |
atom = mol->nextAtom(aiter)) { |
449 |
|
450 |
mass = atom->getMass(); |
451 |
velocity = atom->getVel(); |
452 |
kinetic = mass * (velocity[0]*velocity[0] + velocity[1]*velocity[1] + |
453 |
velocity[2]*velocity[2]) / PhysicalConstants::energyConvert; |
454 |
potential = atom->getParticlePot(); |
455 |
eatom += (kinetic + potential)/2.0; |
456 |
} |
457 |
} |
458 |
|
459 |
int natoms = info_->getNGlobalAtoms(); |
460 |
#ifdef IS_MPI |
461 |
|
462 |
MPI_Allreduce(&eatom, &AvgE_a_, 1, MPI_REALTYPE, MPI_SUM, |
463 |
MPI_COMM_WORLD); |
464 |
#else |
465 |
AvgE_a_ = eatom; |
466 |
#endif |
467 |
AvgE_a_ = AvgE_a_/RealType(natoms); |
468 |
|
469 |
for (mol = info_->beginMolecule(miter); mol != NULL; |
470 |
mol = info_->nextMolecule(miter)) { |
471 |
|
472 |
for (atom = mol->beginAtom(aiter); atom != NULL; |
473 |
atom = mol->nextAtom(aiter)) { |
474 |
|
475 |
/* We think that x_a is relative to the total box and should be a wrapped coordinate */ |
476 |
x_a = atom->getPos(); |
477 |
currSnapshot->wrapVector(x_a); |
478 |
potential = atom->getParticlePot(); |
479 |
velocity = atom->getVel(); |
480 |
kinetic = mass * (velocity[0]*velocity[0] + velocity[1]*velocity[1] + |
481 |
velocity[2]*velocity[2]) / PhysicalConstants::energyConvert; |
482 |
eatom += (kinetic + potential)/2.0; |
483 |
GKappa_t += x_a*(eatom-AvgE_a_); // UNITS KCAL/MOL |
484 |
} |
485 |
} |
486 |
#ifdef IS_MPI |
487 |
MPI_Allreduce(GKappa_t.getArrayPointer(), ThermalHelfandMoment.getArrayPointer(), 3, |
488 |
MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD); |
489 |
#else |
490 |
ThermalHelfandMoment = GKappa_t; |
491 |
#endif |
492 |
return ThermalHelfandMoment; |
493 |
|
494 |
} |
495 |
|
496 |
// Returns the Heat Flux Vector S for the system |
497 |
Vector3d Thermo::getHeatFlux(){ |
498 |
Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot(); |
499 |
SimInfo::MoleculeIterator miter; |
500 |
std::vector<Atom*>::iterator aiter; |
501 |
Molecule* mol; |
502 |
Atom* atom; |
503 |
RealType mass; |
504 |
Vector3d velocity; |
505 |
Vector3d x_a; |
506 |
RealType kinetic; |
507 |
RealType potential; |
508 |
RealType eatom; |
509 |
RealType AvgE_a_ = 0; |
510 |
// Conductive portion of the heat flux |
511 |
Vector3d heatFlux_Jc = V3Zero; |
512 |
Vector3d heatFlux_J; |
513 |
Vector3d heatFlux_Jc_local = V3Zero; |
514 |
|
515 |
/* Calculate conductive portion of the heat flux */ |
516 |
for (mol = info_->beginMolecule(miter); mol != NULL; |
517 |
mol = info_->nextMolecule(miter)) { |
518 |
|
519 |
for (atom = mol->beginAtom(aiter); atom != NULL; |
520 |
atom = mol->nextAtom(aiter)) { |
521 |
|
522 |
mass = atom->getMass(); |
523 |
velocity = atom->getVel(); |
524 |
kinetic = mass * (velocity[0]*velocity[0] + velocity[1]*velocity[1] + |
525 |
velocity[2]*velocity[2]); |
526 |
potential = atom->getParticlePot() * PhysicalConstants::energyConvert; // amu A^2/fs^2 |
527 |
// The potential may not be a 1/2 factor |
528 |
eatom = (kinetic + potential)/2.0; // amu A^2/fs^2 |
529 |
heatFlux_Jc_local[0] += eatom*velocity[0]; // amu A^3/fs^3 |
530 |
heatFlux_Jc_local[1] += eatom*velocity[1]; // amu A^3/fs^3 |
531 |
heatFlux_Jc_local[2] += eatom*velocity[2]; // amu A^3/fs^3 |
532 |
} |
533 |
} |
534 |
//std::cerr << "Heat flux heatFlux_Jc_local is: " << heatFlux_Jc_local << std::endl; |
535 |
/* The J_v vector is reduced in fortan so everyone has the global Jv. Jc is computed over |
536 |
* the local atoms and must be reduced among all processors. |
537 |
*/ |
538 |
#ifdef IS_MPI |
539 |
MPI_Allreduce(heatFlux_Jc_local.getArrayPointer(), heatFlux_Jc.getArrayPointer(), 3, |
540 |
MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD); |
541 |
#else |
542 |
heatFlux_Jc = heatFlux_Jc_local; |
543 |
#endif |
544 |
|
545 |
// (kcal/mol * A/fs) * conversion => (amu A^3)/fs^3 |
546 |
Vector3d heatFlux_Jv = currSnapshot->statData.getJv() * PhysicalConstants::energyConvert; |
547 |
|
548 |
//std::cerr << "Heat flux J_c is: " << heatFlux_Jc << std::endl; |
549 |
//std::cerr << "Heat flux J_v is: " << heatFlux_Jv << std::endl; |
550 |
|
551 |
heatFlux_J = heatFlux_Jv + heatFlux_Jc; // amu A^3/fs^3 + (amu A^3)/fs^3 |
552 |
// Correct for the fact the flux is 1/V (Jc + Jv) |
553 |
RealType volume = this->getVolume(); |
554 |
heatFlux_J = heatFlux_J/volume; // amu/fs^3 |
555 |
|
556 |
|
557 |
|
558 |
|
559 |
return heatFlux_J; |
560 |
|
561 |
|
562 |
} |
563 |
|
564 |
|
565 |
|
566 |
|
567 |
} //end namespace OpenMD |