ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/Thermo.cpp
Revision: 1760
Committed: Thu Jun 21 19:26:46 2012 UTC (12 years, 10 months ago) by gezelter
File size: 18499 byte(s)
Log Message:
Some bugfixes (CholeskyDecomposition), more work on fluctuating charges,
migrating stats stuff into frameData

File Contents

# User Rev Content
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 gezelter 1665 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40     * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 gezelter 246 */
42    
43 gezelter 2 #include <math.h>
44     #include <iostream>
45    
46     #ifdef IS_MPI
47     #include <mpi.h>
48     #endif //is_mpi
49    
50 tim 3 #include "brains/Thermo.hpp"
51 gezelter 246 #include "primitives/Molecule.hpp"
52 tim 3 #include "utils/simError.h"
53 gezelter 1390 #include "utils/PhysicalConstants.hpp"
54 gezelter 1710 #include "types/MultipoleAdapter.hpp"
55 gezelter 2
56 gezelter 1390 namespace OpenMD {
57 gezelter 2
58 tim 963 RealType Thermo::getKinetic() {
59 gezelter 246 SimInfo::MoleculeIterator miter;
60     std::vector<StuntDouble*>::iterator iiter;
61     Molecule* mol;
62     StuntDouble* integrableObject;
63     Vector3d vel;
64     Vector3d angMom;
65     Mat3x3d I;
66     int i;
67     int j;
68     int k;
69 chrisfen 998 RealType mass;
70 tim 963 RealType kinetic = 0.0;
71     RealType kinetic_global = 0.0;
72 gezelter 246
73     for (mol = info_->beginMolecule(miter); mol != NULL; mol = info_->nextMolecule(miter)) {
74 gezelter 507 for (integrableObject = mol->beginIntegrableObject(iiter); integrableObject != NULL;
75     integrableObject = mol->nextIntegrableObject(iiter)) {
76 gezelter 945
77 chrisfen 998 mass = integrableObject->getMass();
78     vel = integrableObject->getVel();
79 gezelter 945
80 gezelter 507 kinetic += mass * (vel[0]*vel[0] + vel[1]*vel[1] + vel[2]*vel[2]);
81 gezelter 945
82 gezelter 507 if (integrableObject->isDirectional()) {
83     angMom = integrableObject->getJ();
84     I = integrableObject->getI();
85 gezelter 2
86 gezelter 507 if (integrableObject->isLinear()) {
87     i = integrableObject->linearAxis();
88     j = (i + 1) % 3;
89     k = (i + 2) % 3;
90     kinetic += angMom[j] * angMom[j] / I(j, j) + angMom[k] * angMom[k] / I(k, k);
91     } else {
92     kinetic += angMom[0]*angMom[0]/I(0, 0) + angMom[1]*angMom[1]/I(1, 1)
93     + angMom[2]*angMom[2]/I(2, 2);
94     }
95     }
96 gezelter 246
97 gezelter 507 }
98 gezelter 246 }
99    
100     #ifdef IS_MPI
101 gezelter 2
102 tim 963 MPI_Allreduce(&kinetic, &kinetic_global, 1, MPI_REALTYPE, MPI_SUM,
103 gezelter 246 MPI_COMM_WORLD);
104     kinetic = kinetic_global;
105 gezelter 2
106 gezelter 246 #endif //is_mpi
107 gezelter 2
108 gezelter 1390 kinetic = kinetic * 0.5 / PhysicalConstants::energyConvert;
109 gezelter 2
110 gezelter 246 return kinetic;
111 gezelter 507 }
112 gezelter 2
113 tim 963 RealType Thermo::getPotential() {
114     RealType potential = 0.0;
115 gezelter 1760
116 gezelter 246 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
117 gezelter 1760 potential = curSnapshot->getShortRangePotential() + curSnapshot->getLongRangePotential();
118 gezelter 246 return potential;
119 gezelter 507 }
120 gezelter 2
121 tim 963 RealType Thermo::getTotalE() {
122     RealType total;
123 gezelter 2
124 gezelter 246 total = this->getKinetic() + this->getPotential();
125     return total;
126 gezelter 507 }
127 gezelter 2
128 tim 963 RealType Thermo::getTemperature() {
129 gezelter 246
130 gezelter 1390 RealType temperature = ( 2.0 * this->getKinetic() ) / (info_->getNdf()* PhysicalConstants::kb );
131 gezelter 246 return temperature;
132 gezelter 507 }
133 gezelter 2
134 gezelter 1715 RealType Thermo::getElectronicTemperature() {
135     SimInfo::MoleculeIterator miter;
136     std::vector<Atom*>::iterator iiter;
137     Molecule* mol;
138     Atom* atom;
139     RealType cvel;
140     RealType cmass;
141     RealType kinetic = 0.0;
142     RealType kinetic_global = 0.0;
143    
144     for (mol = info_->beginMolecule(miter); mol != NULL; mol = info_->nextMolecule(miter)) {
145     for (atom = mol->beginFluctuatingCharge(iiter); atom != NULL;
146     atom = mol->nextFluctuatingCharge(iiter)) {
147     cmass = atom->getChargeMass();
148     cvel = atom->getFlucQVel();
149    
150     kinetic += cmass * cvel * cvel;
151    
152     }
153     }
154    
155     #ifdef IS_MPI
156    
157     MPI_Allreduce(&kinetic, &kinetic_global, 1, MPI_REALTYPE, MPI_SUM,
158     MPI_COMM_WORLD);
159     kinetic = kinetic_global;
160    
161     #endif //is_mpi
162    
163 jmichalk 1733 kinetic = kinetic * 0.5;
164 gezelter 1715 return ( 2.0 * kinetic) / (info_->getNFluctuatingCharges()* PhysicalConstants::kb );
165     }
166    
167    
168    
169    
170 tim 963 RealType Thermo::getVolume() {
171 gezelter 246 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
172     return curSnapshot->getVolume();
173 gezelter 507 }
174 gezelter 2
175 tim 963 RealType Thermo::getPressure() {
176 gezelter 2
177 gezelter 246 // Relies on the calculation of the full molecular pressure tensor
178 gezelter 2
179    
180 gezelter 246 Mat3x3d tensor;
181 tim 963 RealType pressure;
182 gezelter 2
183 gezelter 246 tensor = getPressureTensor();
184 gezelter 2
185 gezelter 1390 pressure = PhysicalConstants::pressureConvert * (tensor(0, 0) + tensor(1, 1) + tensor(2, 2)) / 3.0;
186 gezelter 2
187 gezelter 246 return pressure;
188 gezelter 507 }
189 gezelter 2
190 tim 963 RealType Thermo::getPressure(int direction) {
191 tim 538
192     // Relies on the calculation of the full molecular pressure tensor
193    
194    
195     Mat3x3d tensor;
196 tim 963 RealType pressure;
197 tim 538
198     tensor = getPressureTensor();
199    
200 gezelter 1390 pressure = PhysicalConstants::pressureConvert * tensor(direction, direction);
201 tim 538
202     return pressure;
203     }
204    
205 gezelter 507 Mat3x3d Thermo::getPressureTensor() {
206 gezelter 246 // returns pressure tensor in units amu*fs^-2*Ang^-1
207     // routine derived via viral theorem description in:
208     // Paci, E. and Marchi, M. J.Phys.Chem. 1996, 100, 4314-4322
209     Mat3x3d pressureTensor;
210     Mat3x3d p_local(0.0);
211     Mat3x3d p_global(0.0);
212 gezelter 2
213 gezelter 246 SimInfo::MoleculeIterator i;
214     std::vector<StuntDouble*>::iterator j;
215     Molecule* mol;
216     StuntDouble* integrableObject;
217     for (mol = info_->beginMolecule(i); mol != NULL; mol = info_->nextMolecule(i)) {
218 gezelter 507 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
219     integrableObject = mol->nextIntegrableObject(j)) {
220 gezelter 2
221 tim 963 RealType mass = integrableObject->getMass();
222 gezelter 507 Vector3d vcom = integrableObject->getVel();
223     p_local += mass * outProduct(vcom, vcom);
224     }
225 gezelter 246 }
226 gezelter 2
227     #ifdef IS_MPI
228 tim 963 MPI_Allreduce(p_local.getArrayPointer(), p_global.getArrayPointer(), 9, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
229 gezelter 2 #else
230 gezelter 246 p_global = p_local;
231 gezelter 2 #endif // is_mpi
232    
233 tim 963 RealType volume = this->getVolume();
234 gezelter 246 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
235 gezelter 1723 Mat3x3d stressTensor = curSnapshot->getStressTensor();
236 gezelter 1126
237 gezelter 1723 pressureTensor = (p_global +
238     PhysicalConstants::energyConvert * stressTensor)/volume;
239 chrisfen 998
240 gezelter 246 return pressureTensor;
241 gezelter 507 }
242 gezelter 2
243 chrisfen 998
244 gezelter 507 void Thermo::saveStat(){
245 gezelter 246 Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
246     Stats& stat = currSnapshot->statData;
247 gezelter 2
248 gezelter 246 stat[Stats::KINETIC_ENERGY] = getKinetic();
249     stat[Stats::POTENTIAL_ENERGY] = getPotential();
250     stat[Stats::TOTAL_ENERGY] = stat[Stats::KINETIC_ENERGY] + stat[Stats::POTENTIAL_ENERGY] ;
251     stat[Stats::TEMPERATURE] = getTemperature();
252     stat[Stats::PRESSURE] = getPressure();
253     stat[Stats::VOLUME] = getVolume();
254 gezelter 2
255 tim 541 Mat3x3d tensor =getPressureTensor();
256 gezelter 1126 stat[Stats::PRESSURE_TENSOR_XX] = tensor(0, 0);
257     stat[Stats::PRESSURE_TENSOR_XY] = tensor(0, 1);
258     stat[Stats::PRESSURE_TENSOR_XZ] = tensor(0, 2);
259     stat[Stats::PRESSURE_TENSOR_YX] = tensor(1, 0);
260     stat[Stats::PRESSURE_TENSOR_YY] = tensor(1, 1);
261     stat[Stats::PRESSURE_TENSOR_YZ] = tensor(1, 2);
262     stat[Stats::PRESSURE_TENSOR_ZX] = tensor(2, 0);
263     stat[Stats::PRESSURE_TENSOR_ZY] = tensor(2, 1);
264     stat[Stats::PRESSURE_TENSOR_ZZ] = tensor(2, 2);
265 tim 541
266 gezelter 1503 // grab the simulation box dipole moment if specified
267     if (info_->getCalcBoxDipole()){
268     Vector3d totalDipole = getBoxDipole();
269     stat[Stats::BOX_DIPOLE_X] = totalDipole(0);
270     stat[Stats::BOX_DIPOLE_Y] = totalDipole(1);
271     stat[Stats::BOX_DIPOLE_Z] = totalDipole(2);
272     }
273 tim 541
274 gezelter 1291 Globals* simParams = info_->getSimParams();
275 gezelter 1723 // grab the heat flux if desired
276     if (simParams->havePrintHeatFlux()) {
277     if (simParams->getPrintHeatFlux()){
278     Vector3d heatFlux = getHeatFlux();
279     stat[Stats::HEATFLUX_X] = heatFlux(0);
280     stat[Stats::HEATFLUX_Y] = heatFlux(1);
281     stat[Stats::HEATFLUX_Z] = heatFlux(2);
282     }
283     }
284 gezelter 1291
285     if (simParams->haveTaggedAtomPair() &&
286     simParams->havePrintTaggedPairDistance()) {
287     if ( simParams->getPrintTaggedPairDistance()) {
288    
289     std::pair<int, int> tap = simParams->getTaggedAtomPair();
290     Vector3d pos1, pos2, rab;
291    
292     #ifdef IS_MPI
293 gezelter 1313 std::cerr << "tap = " << tap.first << " " << tap.second << std::endl;
294 gezelter 1291
295 chuckv 1292 int mol1 = info_->getGlobalMolMembership(tap.first);
296     int mol2 = info_->getGlobalMolMembership(tap.second);
297 gezelter 1313 std::cerr << "mols = " << mol1 << " " << mol2 << std::endl;
298    
299 gezelter 1291 int proc1 = info_->getMolToProc(mol1);
300     int proc2 = info_->getMolToProc(mol2);
301    
302 gezelter 1313 std::cerr << " procs = " << proc1 << " " <<proc2 <<std::endl;
303    
304 chuckv 1292 RealType data[3];
305 gezelter 1291 if (proc1 == worldRank) {
306     StuntDouble* sd1 = info_->getIOIndexToIntegrableObject(tap.first);
307 gezelter 1313 std::cerr << " on proc " << proc1 << ", sd1 has global index= " << sd1->getGlobalIndex() << std::endl;
308 gezelter 1291 pos1 = sd1->getPos();
309     data[0] = pos1.x();
310     data[1] = pos1.y();
311     data[2] = pos1.z();
312     MPI_Bcast(data, 3, MPI_REALTYPE, proc1, MPI_COMM_WORLD);
313     } else {
314     MPI_Bcast(data, 3, MPI_REALTYPE, proc1, MPI_COMM_WORLD);
315     pos1 = Vector3d(data);
316     }
317 chuckv 1292
318    
319 gezelter 1291 if (proc2 == worldRank) {
320     StuntDouble* sd2 = info_->getIOIndexToIntegrableObject(tap.second);
321 gezelter 1313 std::cerr << " on proc " << proc2 << ", sd2 has global index= " << sd2->getGlobalIndex() << std::endl;
322 gezelter 1291 pos2 = sd2->getPos();
323     data[0] = pos2.x();
324     data[1] = pos2.y();
325     data[2] = pos2.z();
326     MPI_Bcast(data, 3, MPI_REALTYPE, proc2, MPI_COMM_WORLD);
327     } else {
328     MPI_Bcast(data, 3, MPI_REALTYPE, proc2, MPI_COMM_WORLD);
329     pos2 = Vector3d(data);
330     }
331     #else
332     StuntDouble* at1 = info_->getIOIndexToIntegrableObject(tap.first);
333     StuntDouble* at2 = info_->getIOIndexToIntegrableObject(tap.second);
334     pos1 = at1->getPos();
335     pos2 = at2->getPos();
336     #endif
337     rab = pos2 - pos1;
338     currSnapshot->wrapVector(rab);
339     stat[Stats::TAGGED_PAIR_DISTANCE] = rab.length();
340     }
341     }
342    
343 gezelter 246 /**@todo need refactorying*/
344     //Conserved Quantity is set by integrator and time is set by setTime
345 gezelter 2
346 gezelter 507 }
347 gezelter 2
348 gezelter 1503
349     Vector3d Thermo::getBoxDipole() {
350     Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
351     SimInfo::MoleculeIterator miter;
352     std::vector<Atom*>::iterator aiter;
353     Molecule* mol;
354     Atom* atom;
355     RealType charge;
356     RealType moment(0.0);
357     Vector3d ri(0.0);
358     Vector3d dipoleVector(0.0);
359     Vector3d nPos(0.0);
360     Vector3d pPos(0.0);
361     RealType nChg(0.0);
362     RealType pChg(0.0);
363     int nCount = 0;
364     int pCount = 0;
365    
366     RealType chargeToC = 1.60217733e-19;
367     RealType angstromToM = 1.0e-10;
368     RealType debyeToCm = 3.33564095198e-30;
369    
370     for (mol = info_->beginMolecule(miter); mol != NULL;
371     mol = info_->nextMolecule(miter)) {
372    
373     for (atom = mol->beginAtom(aiter); atom != NULL;
374     atom = mol->nextAtom(aiter)) {
375    
376     if (atom->isCharge() ) {
377     charge = 0.0;
378     GenericData* data = atom->getAtomType()->getPropertyByName("Charge");
379     if (data != NULL) {
380    
381     charge = (dynamic_cast<DoubleGenericData*>(data))->getData();
382     charge *= chargeToC;
383    
384     ri = atom->getPos();
385     currSnapshot->wrapVector(ri);
386     ri *= angstromToM;
387    
388     if (charge < 0.0) {
389     nPos += ri;
390     nChg -= charge;
391     nCount++;
392     } else if (charge > 0.0) {
393     pPos += ri;
394     pChg += charge;
395     pCount++;
396     }
397     }
398     }
399    
400 gezelter 1710 MultipoleAdapter ma = MultipoleAdapter(atom->getAtomType());
401     if (ma.isDipole() ) {
402 gezelter 1503 Vector3d u_i = atom->getElectroFrame().getColumn(2);
403 gezelter 1710 moment = ma.getDipoleMoment();
404     moment *= debyeToCm;
405     dipoleVector += u_i * moment;
406 gezelter 1503 }
407     }
408     }
409    
410    
411     #ifdef IS_MPI
412     RealType pChg_global, nChg_global;
413     int pCount_global, nCount_global;
414     Vector3d pPos_global, nPos_global, dipVec_global;
415    
416     MPI_Allreduce(&pChg, &pChg_global, 1, MPI_REALTYPE, MPI_SUM,
417     MPI_COMM_WORLD);
418     pChg = pChg_global;
419     MPI_Allreduce(&nChg, &nChg_global, 1, MPI_REALTYPE, MPI_SUM,
420     MPI_COMM_WORLD);
421     nChg = nChg_global;
422     MPI_Allreduce(&pCount, &pCount_global, 1, MPI_INTEGER, MPI_SUM,
423     MPI_COMM_WORLD);
424     pCount = pCount_global;
425     MPI_Allreduce(&nCount, &nCount_global, 1, MPI_INTEGER, MPI_SUM,
426     MPI_COMM_WORLD);
427     nCount = nCount_global;
428     MPI_Allreduce(pPos.getArrayPointer(), pPos_global.getArrayPointer(), 3,
429     MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
430     pPos = pPos_global;
431     MPI_Allreduce(nPos.getArrayPointer(), nPos_global.getArrayPointer(), 3,
432     MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
433     nPos = nPos_global;
434     MPI_Allreduce(dipoleVector.getArrayPointer(),
435     dipVec_global.getArrayPointer(), 3,
436     MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
437     dipoleVector = dipVec_global;
438     #endif //is_mpi
439    
440     // first load the accumulated dipole moment (if dipoles were present)
441     Vector3d boxDipole = dipoleVector;
442     // now include the dipole moment due to charges
443     // use the lesser of the positive and negative charge totals
444     RealType chg_value = nChg <= pChg ? nChg : pChg;
445    
446     // find the average positions
447     if (pCount > 0 && nCount > 0 ) {
448     pPos /= pCount;
449     nPos /= nCount;
450     }
451    
452     // dipole is from the negative to the positive (physics notation)
453     boxDipole += (pPos - nPos) * chg_value;
454    
455     return boxDipole;
456     }
457 gezelter 1723
458     // Returns the Heat Flux Vector for the system
459     Vector3d Thermo::getHeatFlux(){
460     Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
461     SimInfo::MoleculeIterator miter;
462     std::vector<StuntDouble*>::iterator iiter;
463     Molecule* mol;
464     StuntDouble* integrableObject;
465     RigidBody::AtomIterator ai;
466     Atom* atom;
467     Vector3d vel;
468     Vector3d angMom;
469     Mat3x3d I;
470     int i;
471     int j;
472     int k;
473     RealType mass;
474    
475     Vector3d x_a;
476     RealType kinetic;
477     RealType potential;
478     RealType eatom;
479     RealType AvgE_a_ = 0;
480     // Convective portion of the heat flux
481     Vector3d heatFluxJc = V3Zero;
482    
483     /* Calculate convective portion of the heat flux */
484     for (mol = info_->beginMolecule(miter); mol != NULL;
485     mol = info_->nextMolecule(miter)) {
486    
487     for (integrableObject = mol->beginIntegrableObject(iiter);
488     integrableObject != NULL;
489     integrableObject = mol->nextIntegrableObject(iiter)) {
490    
491     mass = integrableObject->getMass();
492     vel = integrableObject->getVel();
493    
494     kinetic = mass * (vel[0]*vel[0] + vel[1]*vel[1] + vel[2]*vel[2]);
495    
496     if (integrableObject->isDirectional()) {
497     angMom = integrableObject->getJ();
498     I = integrableObject->getI();
499    
500     if (integrableObject->isLinear()) {
501     i = integrableObject->linearAxis();
502     j = (i + 1) % 3;
503     k = (i + 2) % 3;
504     kinetic += angMom[j] * angMom[j] / I(j, j) + angMom[k] * angMom[k] / I(k, k);
505     } else {
506     kinetic += angMom[0]*angMom[0]/I(0, 0) + angMom[1]*angMom[1]/I(1, 1)
507     + angMom[2]*angMom[2]/I(2, 2);
508     }
509     }
510    
511     potential = 0.0;
512    
513     if (integrableObject->isRigidBody()) {
514     RigidBody* rb = dynamic_cast<RigidBody*>(integrableObject);
515     for (atom = rb->beginAtom(ai); atom != NULL;
516     atom = rb->nextAtom(ai)) {
517     potential += atom->getParticlePot();
518     }
519     } else {
520     potential = integrableObject->getParticlePot();
521     cerr << "ppot = " << potential << "\n";
522     }
523    
524     potential *= PhysicalConstants::energyConvert; // amu A^2/fs^2
525     // The potential may not be a 1/2 factor
526     eatom = (kinetic + potential)/2.0; // amu A^2/fs^2
527     heatFluxJc[0] += eatom*vel[0]; // amu A^3/fs^3
528     heatFluxJc[1] += eatom*vel[1]; // amu A^3/fs^3
529     heatFluxJc[2] += eatom*vel[2]; // amu A^3/fs^3
530     }
531     }
532    
533     std::cerr << "Heat flux heatFluxJc is: " << heatFluxJc << std::endl;
534    
535     /* The J_v vector is reduced in fortan so everyone has the global
536     * Jv. Jc is computed over the local atoms and must be reduced
537     * among all processors.
538     */
539     #ifdef IS_MPI
540     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &heatFluxJc[0], 3, MPI::REALTYPE,
541     MPI::SUM);
542     #endif
543    
544     // (kcal/mol * A/fs) * conversion => (amu A^3)/fs^3
545    
546     Vector3d heatFluxJv = currSnapshot->getConductiveHeatFlux() *
547     PhysicalConstants::energyConvert;
548    
549     std::cerr << "Heat flux Jc is: " << heatFluxJc << std::endl;
550     std::cerr << "Heat flux Jv is: " << heatFluxJv << std::endl;
551    
552     // Correct for the fact the flux is 1/V (Jc + Jv)
553     return (heatFluxJv + heatFluxJc) / this->getVolume(); // amu / fs^3
554     }
555 gezelter 1390 } //end namespace OpenMD

Properties

Name Value
svn:keywords Author Id Revision Date