OpenMD 3.0
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
Thermo.cpp
1/*
2 * Copyright (c) 2004-present, The University of Notre Dame. All rights
3 * reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * 1. Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice,
12 * this list of conditions and the following disclaimer in the documentation
13 * and/or other materials provided with the distribution.
14 *
15 * 3. Neither the name of the copyright holder nor the names of its
16 * contributors may be used to endorse or promote products derived from
17 * this software without specific prior written permission.
18 *
19 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
20 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
21 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
22 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
23 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
24 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
25 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
26 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
27 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
28 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
29 * POSSIBILITY OF SUCH DAMAGE.
30 *
31 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
32 * research, please cite the appropriate papers when you publish your
33 * work. Good starting points are:
34 *
35 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
36 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
37 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
38 * [4] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
39 * [5] Kuang & Gezelter, Mol. Phys., 110, 691-701 (2012).
40 * [6] Lamichhane, Gezelter & Newman, J. Chem. Phys. 141, 134109 (2014).
41 * [7] Lamichhane, Newman & Gezelter, J. Chem. Phys. 141, 134110 (2014).
42 * [8] Bhattarai, Newman & Gezelter, Phys. Rev. B 99, 094106 (2019).
43 */
44
45#include "brains/Thermo.hpp"
46
47#include <cmath>
48#include <iostream>
49
50#ifdef IS_MPI
51#include <mpi.h>
52#endif
53
55#include "types/FixedChargeAdapter.hpp"
56#include "types/FluctuatingChargeAdapter.hpp"
57#include "types/MultipoleAdapter.hpp"
58#include "utils/Constants.hpp"
59#include "utils/simError.h"
60#ifdef HAVE_QHULL
61#include "math/AlphaHull.hpp"
62#include "math/ConvexHull.hpp"
63#endif
64
65using namespace std;
66namespace OpenMD {
67
68 RealType Thermo::getTranslationalKinetic() {
70
71 if (!snap->hasTranslationalKineticEnergy) {
72 SimInfo::MoleculeIterator miter;
73 vector<StuntDouble*>::iterator iiter;
74 Molecule* mol;
75 StuntDouble* sd;
76 Vector3d vel;
77 RealType mass;
78 RealType kinetic(0.0);
79
80 for (mol = info_->beginMolecule(miter); mol != NULL;
81 mol = info_->nextMolecule(miter)) {
82 for (sd = mol->beginIntegrableObject(iiter); sd != NULL;
83 sd = mol->nextIntegrableObject(iiter)) {
84 mass = sd->getMass();
85 vel = sd->getVel();
86
87 kinetic +=
88 mass * (vel[0] * vel[0] + vel[1] * vel[1] + vel[2] * vel[2]);
89 }
90 }
91
92#ifdef IS_MPI
93 MPI_Allreduce(MPI_IN_PLACE, &kinetic, 1, MPI_REALTYPE, MPI_SUM,
94 MPI_COMM_WORLD);
95#endif
96
97 kinetic = kinetic * 0.5 / Constants::energyConvert;
98
99 snap->setTranslationalKineticEnergy(kinetic);
100 }
101 return snap->getTranslationalKineticEnergy();
102 }
103
104 RealType Thermo::getRotationalKinetic() {
106
107 if (!snap->hasRotationalKineticEnergy) {
108 SimInfo::MoleculeIterator miter;
109 vector<StuntDouble*>::iterator iiter;
110 Molecule* mol;
111 StuntDouble* sd;
112 Vector3d angMom;
113 Mat3x3d I;
114 int i, j, k;
115 RealType kinetic(0.0);
116
117 for (mol = info_->beginMolecule(miter); mol != NULL;
118 mol = info_->nextMolecule(miter)) {
119 for (sd = mol->beginIntegrableObject(iiter); sd != NULL;
120 sd = mol->nextIntegrableObject(iiter)) {
121 if (sd->isDirectional()) {
122 angMom = sd->getJ();
123 I = sd->getI();
124
125 if (sd->isLinear()) {
126 i = sd->linearAxis();
127 j = (i + 1) % 3;
128 k = (i + 2) % 3;
129 kinetic += angMom[j] * angMom[j] / I(j, j) +
130 angMom[k] * angMom[k] / I(k, k);
131 } else {
132 kinetic += angMom[0] * angMom[0] / I(0, 0) +
133 angMom[1] * angMom[1] / I(1, 1) +
134 angMom[2] * angMom[2] / I(2, 2);
135 }
136 }
137 }
138 }
139
140#ifdef IS_MPI
141 MPI_Allreduce(MPI_IN_PLACE, &kinetic, 1, MPI_REALTYPE, MPI_SUM,
142 MPI_COMM_WORLD);
143#endif
144
145 kinetic = kinetic * 0.5 / Constants::energyConvert;
146
147 snap->setRotationalKineticEnergy(kinetic);
148 }
149 return snap->getRotationalKineticEnergy();
150 }
151
152 RealType Thermo::getKinetic() {
154
155 if (!snap->hasKineticEnergy) {
156 RealType ke = getTranslationalKinetic() + getRotationalKinetic() +
157 getElectronicKinetic();
158
159 snap->setKineticEnergy(ke);
160 }
161 return snap->getKineticEnergy();
162 }
163
164 RealType Thermo::getPotential() {
165 // ForceManager computes the potential and stores it in the
166 // Snapshot. All we have to do is report it.
167
169 return snap->getPotentialEnergy();
170 }
171
172 potVec Thermo::getSelectionPotentials() {
173 // ForceManager computes the selection potentials and stores them
174 // in the Snapshot. All we have to do is report them.
175
177 return snap->getSelectionPotentials();
178 }
179
180 RealType Thermo::getTotalEnergy() {
182
183 if (!snap->hasTotalEnergy) {
184 snap->setTotalEnergy(this->getKinetic() + this->getPotential());
185 }
186
187 return snap->getTotalEnergy();
188 }
189
190 /*
191 * Returns only the nuclear portion of the temperature - see
192 * getElectronicTemperature for the electronic portion */
193 RealType Thermo::getTemperature() {
195
196 if (!snap->hasTemperature) {
197 RealType nuclearKE =
198 this->getTranslationalKinetic() + this->getRotationalKinetic();
199
200 RealType temperature =
201 (2.0 * nuclearKE) / (info_->getNdf() * Constants::kb);
202
203 snap->setTemperature(temperature);
204 }
205
206 return snap->getTemperature();
207 }
208
209 RealType Thermo::getElectronicKinetic() {
211
212 if (!snap->hasElectronicKineticEnergy) {
213 SimInfo::MoleculeIterator miter;
214 vector<Atom*>::iterator iiter;
215 Molecule* mol;
216 Atom* atom;
217 RealType cvel;
218 RealType cmass;
219 RealType kinetic(0.0);
220
221 for (mol = info_->beginMolecule(miter); mol != NULL;
222 mol = info_->nextMolecule(miter)) {
223 for (atom = mol->beginFluctuatingCharge(iiter); atom != NULL;
224 atom = mol->nextFluctuatingCharge(iiter)) {
225 cmass = atom->getChargeMass();
226 cvel = atom->getFlucQVel();
227 kinetic += cmass * cvel * cvel;
228 }
229 }
230
231#ifdef IS_MPI
232 MPI_Allreduce(MPI_IN_PLACE, &kinetic, 1, MPI_REALTYPE, MPI_SUM,
233 MPI_COMM_WORLD);
234#endif
235
236 kinetic *= 0.5;
237 snap->setElectronicKineticEnergy(kinetic);
238 }
239
240 return snap->getElectronicKineticEnergy();
241 }
242
243 RealType Thermo::getElectronicTemperature() {
245
246 if (!snap->hasElectronicTemperature) {
247 RealType eTemp = (2.0 * this->getElectronicKinetic()) /
248 (info_->getNFluctuatingCharges() * Constants::kb);
249
250 snap->setElectronicTemperature(eTemp);
251 }
252
253 return snap->getElectronicTemperature();
254 }
255
256 RealType Thermo::getNetCharge() {
258
259 if (!snap->hasNetCharge) {
260 SimInfo::MoleculeIterator miter;
261 vector<Atom*>::iterator aiter;
262 Molecule* mol;
263 Atom* atom;
264 RealType charge;
265 RealType netCharge(0.0);
266
267 for (mol = info_->beginMolecule(miter); mol != NULL;
268 mol = info_->nextMolecule(miter)) {
269 for (atom = mol->beginAtom(aiter); atom != NULL;
270 atom = mol->nextAtom(aiter)) {
271 charge = 0.0;
272
273 FixedChargeAdapter fca = FixedChargeAdapter(atom->getAtomType());
274 if (fca.isFixedCharge()) { charge = fca.getCharge(); }
275
277 FluctuatingChargeAdapter(atom->getAtomType());
278 if (fqa.isFluctuatingCharge()) { charge += atom->getFlucQPos(); }
279
280 netCharge += charge;
281 }
282 }
283
284#ifdef IS_MPI
285 MPI_Allreduce(MPI_IN_PLACE, &netCharge, 1, MPI_REALTYPE, MPI_SUM,
286 MPI_COMM_WORLD);
287#endif
288
289 snap->setNetCharge(netCharge);
290 }
291
292 return snap->getNetCharge();
293 }
294
295 RealType Thermo::getChargeMomentum() {
297
298 if (!snap->hasChargeMomentum) {
299 SimInfo::MoleculeIterator miter;
300 vector<Atom*>::iterator iiter;
301 Molecule* mol;
302 Atom* atom;
303 RealType cvel;
304 RealType cmass;
305 RealType momentum(0.0);
306
307 for (mol = info_->beginMolecule(miter); mol != NULL;
308 mol = info_->nextMolecule(miter)) {
309 for (atom = mol->beginFluctuatingCharge(iiter); atom != NULL;
310 atom = mol->nextFluctuatingCharge(iiter)) {
311 cmass = atom->getChargeMass();
312 cvel = atom->getFlucQVel();
313
314 momentum += cmass * cvel;
315 }
316 }
317
318#ifdef IS_MPI
319 MPI_Allreduce(MPI_IN_PLACE, &momentum, 1, MPI_REALTYPE, MPI_SUM,
320 MPI_COMM_WORLD);
321#endif
322
323 snap->setChargeMomentum(momentum);
324 }
325
326 return snap->getChargeMomentum();
327 }
328
329 std::vector<Vector3d> Thermo::getCurrentDensity() {
331 AtomTypeSet simTypes = info_->getSimulatedAtomTypes();
332
333 SimInfo::MoleculeIterator miter;
334 std::vector<Atom*>::iterator iiter;
335 std::vector<RigidBody*>::iterator ri;
336 Molecule* mol;
337 RigidBody* rb;
338 Atom* atom;
339 AtomType* atype;
340 AtomTypeSet::iterator at;
341 Vector3d Jc(0.0);
342 std::vector<Vector3d> typeJc(simTypes.size(), V3Zero);
343
344 for (mol = info_->beginMolecule(miter); mol != NULL;
345 mol = info_->nextMolecule(miter)) {
346 // change the velocities of atoms which belong to the rigidbodies
347 for (rb = mol->beginRigidBody(ri); rb != NULL;
348 rb = mol->nextRigidBody(ri)) {
349 rb->updateAtomVel();
350 }
351
352 for (atom = mol->beginAtom(iiter); atom != NULL;
353 atom = mol->nextAtom(iiter)) {
354 Vector3d v = atom->getVel();
355 RealType q = 0.0;
356 int typeIndex(-1);
357
358 atype = atom->getAtomType();
360 if (fca.isFixedCharge()) q = fca.getCharge();
362 if (fqa.isFluctuatingCharge()) q += atom->getFlucQPos();
363
364 typeIndex = -1;
365 at = std::find(simTypes.begin(), simTypes.end(), atype);
366 if (at != simTypes.end()) {
367 typeIndex = std::distance(simTypes.begin(), at);
368 }
369
370 if (typeIndex != -1) { typeJc[typeIndex] += q * v; }
371 Jc += q * v;
372 }
373 }
374
375#ifdef IS_MPI
376 MPI_Allreduce(MPI_IN_PLACE, &Jc[0], 3, MPI_REALTYPE, MPI_SUM,
377 MPI_COMM_WORLD);
378 for (unsigned int j = 0; j < simTypes.size(); j++) {
379 MPI_Allreduce(MPI_IN_PLACE, &typeJc[j][0], 3, MPI_REALTYPE, MPI_SUM,
380 MPI_COMM_WORLD);
381 }
382#endif
383
384 RealType vol = snap->getVolume();
385
386 Jc /= (vol * Constants::currentDensityConvert);
387 for (unsigned int j = 0; j < simTypes.size(); j++) {
388 typeJc[j] /= (vol * Constants::currentDensityConvert);
389 }
390
391 std::vector<Vector3d> result;
392 result.clear();
393
394 result.push_back(Jc);
395 for (unsigned int j = 0; j < simTypes.size(); j++)
396 result.push_back(typeJc[j]);
397
398 return result;
399 }
400
401 RealType Thermo::getVolume() {
403 return snap->getVolume();
404 }
405
406 RealType Thermo::getPressure() {
408
409 if (!snap->hasPressure) {
410 // Relies on the calculation of the full molecular pressure tensor
411
412 Mat3x3d tensor;
413 RealType pressure;
414
415 tensor = getPressureTensor();
416
417 pressure = Constants::pressureConvert *
418 (tensor(0, 0) + tensor(1, 1) + tensor(2, 2)) / 3.0;
419
420 snap->setPressure(pressure);
421 }
422
423 return snap->getPressure();
424 }
425
426 RealType Thermo::getPressure(Snapshot* snap) {
427 if (!snap->hasPressure) {
428 // Relies on the calculation of the full molecular pressure tensor
429 Mat3x3d tensor;
430 RealType pressure;
431
432 tensor = getPressureTensor(snap);
433
434 pressure = Constants::pressureConvert *
435 (tensor(0, 0) + tensor(1, 1) + tensor(2, 2)) / 3.0;
436
437 snap->setPressure(pressure);
438 }
439
440 return snap->getPressure();
441 }
442
444 // returns pressure tensor in units amu*fs^-2*Ang^-1
445 // routine derived via viral theorem description in:
446 // Paci, E. and Marchi, M. J.Phys.Chem. 1996, 100, 4314-4322
448
449 if (!snap->hasPressureTensor) {
450 Mat3x3d pressureTensor;
451 Mat3x3d p_tens(0.0);
452 RealType mass;
453 Vector3d vcom;
454
455 SimInfo::MoleculeIterator i;
456 vector<StuntDouble*>::iterator j;
457 Molecule* mol;
458 StuntDouble* sd;
459 for (mol = info_->beginMolecule(i); mol != NULL;
460 mol = info_->nextMolecule(i)) {
461 for (sd = mol->beginIntegrableObject(j); sd != NULL;
462 sd = mol->nextIntegrableObject(j)) {
463 mass = sd->getMass();
464 vcom = sd->getVel();
465 p_tens += mass * outProduct(vcom, vcom);
466 }
467 }
468
469#ifdef IS_MPI
470 MPI_Allreduce(MPI_IN_PLACE, p_tens.getArrayPointer(), 9, MPI_REALTYPE,
471 MPI_SUM, MPI_COMM_WORLD);
472#endif
473
474 RealType volume = this->getVolume();
475 Mat3x3d virialTensor = snap->getVirialTensor();
476
477 pressureTensor =
478 (p_tens + Constants::energyConvert * virialTensor) / volume;
479
480 snap->setPressureTensor(pressureTensor);
481 }
482 return snap->getPressureTensor();
483 }
484
486 // returns pressure tensor in units amu*fs^-2*Ang^-1
487 // routine derived via viral theorem description in:
488 // Paci, E. and Marchi, M. J.Phys.Chem. 1996, 100, 4314-4322
489 if (!snap->hasPressureTensor) {
490 Mat3x3d pressureTensor;
491 Mat3x3d p_tens(0.0);
492 RealType mass;
493 Vector3d vcom;
494
495 SimInfo::MoleculeIterator i;
496 vector<StuntDouble*>::iterator j;
497 Molecule* mol;
498 StuntDouble* sd;
499 for (mol = info_->beginMolecule(i); mol != NULL;
500 mol = info_->nextMolecule(i)) {
501 for (sd = mol->beginIntegrableObject(j); sd != NULL;
502 sd = mol->nextIntegrableObject(j)) {
503 mass = sd->getMass();
504 vcom = sd->getVel(snap);
505 p_tens += mass * outProduct(vcom, vcom);
506 }
507 }
508
509#ifdef IS_MPI
510 MPI_Allreduce(MPI_IN_PLACE, p_tens.getArrayPointer(), 9, MPI_REALTYPE,
511 MPI_SUM, MPI_COMM_WORLD);
512#endif
513
514 RealType volume = this->getVolume();
515 Mat3x3d virialTensor = snap->getVirialTensor();
516
517 pressureTensor =
518 (p_tens + Constants::energyConvert * virialTensor) / volume;
519
520 snap->setPressureTensor(pressureTensor);
521 }
522 return snap->getPressureTensor();
523 }
524
527
528 if (!snap->hasSystemDipole) {
529 SimInfo::MoleculeIterator miter;
530 vector<Atom*>::iterator aiter;
531 Molecule* mol;
532 Atom* atom;
533 RealType charge;
534 Vector3d ri(0.0);
535 Vector3d dipoleVector(0.0);
536 Vector3d nPos(0.0);
537 Vector3d pPos(0.0);
538 RealType nChg(0.0);
539 RealType pChg(0.0);
540 int nCount = 0;
541 int pCount = 0;
542
543 RealType chargeToC = 1.60217733e-19;
544 RealType angstromToM = 1.0e-10;
545 RealType debyeToCm = 3.33564095198e-30;
546
547 for (mol = info_->beginMolecule(miter); mol != NULL;
548 mol = info_->nextMolecule(miter)) {
549 for (atom = mol->beginAtom(aiter); atom != NULL;
550 atom = mol->nextAtom(aiter)) {
551 charge = 0.0;
552
553 FixedChargeAdapter fca = FixedChargeAdapter(atom->getAtomType());
554 if (fca.isFixedCharge()) { charge = fca.getCharge(); }
555
557 FluctuatingChargeAdapter(atom->getAtomType());
558 if (fqa.isFluctuatingCharge()) { charge += atom->getFlucQPos(); }
559
560 charge *= chargeToC;
561
562 ri = atom->getPos();
563 snap->wrapVector(ri);
564 ri *= angstromToM;
565
566 if (charge < 0.0) {
567 nPos += ri;
568 nChg -= charge;
569 nCount++;
570 } else if (charge > 0.0) {
571 pPos += ri;
572 pChg += charge;
573 pCount++;
574 }
575
576 if (atom->isDipole()) {
577 dipoleVector += atom->getDipole() * debyeToCm;
578 }
579 }
580 }
581
582#ifdef IS_MPI
583 MPI_Allreduce(MPI_IN_PLACE, &pChg, 1, MPI_REALTYPE, MPI_SUM,
584 MPI_COMM_WORLD);
585 MPI_Allreduce(MPI_IN_PLACE, &nChg, 1, MPI_REALTYPE, MPI_SUM,
586 MPI_COMM_WORLD);
587
588 MPI_Allreduce(MPI_IN_PLACE, &pCount, 1, MPI_INTEGER, MPI_SUM,
589 MPI_COMM_WORLD);
590 MPI_Allreduce(MPI_IN_PLACE, &nCount, 1, MPI_INTEGER, MPI_SUM,
591 MPI_COMM_WORLD);
592
593 MPI_Allreduce(MPI_IN_PLACE, pPos.getArrayPointer(), 3, MPI_REALTYPE,
594 MPI_SUM, MPI_COMM_WORLD);
595 MPI_Allreduce(MPI_IN_PLACE, nPos.getArrayPointer(), 3, MPI_REALTYPE,
596 MPI_SUM, MPI_COMM_WORLD);
597
598 MPI_Allreduce(MPI_IN_PLACE, dipoleVector.getArrayPointer(), 3,
599 MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
600#endif
601
602 // first load the accumulated dipole moment (if dipoles were present)
603 Vector3d boxDipole = dipoleVector;
604 // now include the dipole moment due to charges
605 // use the lesser of the positive and negative charge totals
606 RealType chg_value = nChg <= pChg ? nChg : pChg;
607
608 // find the average positions
609 if (pCount > 0 && nCount > 0) {
610 pPos /= pCount;
611 nPos /= nCount;
612 }
613
614 // dipole is from the negative to the positive (physics notation)
615 boxDipole += (pPos - nPos) * chg_value;
616 snap->setSystemDipole(boxDipole);
617 }
618
619 return snap->getSystemDipole();
620 }
621
624
625 if (!snap->hasSystemQuadrupole) {
626 SimInfo::MoleculeIterator miter;
627 vector<Atom*>::iterator aiter;
628 Molecule* mol;
629 Atom* atom;
630 RealType charge;
631 Vector3d ri(0.0);
632 Vector3d dipole(0.0);
633 Mat3x3d qpole(0.0);
634
635 RealType chargeToC = 1.60217733e-19;
636 RealType angstromToM = 1.0e-10;
637 RealType debyeToCm = 3.33564095198e-30;
638
639 for (mol = info_->beginMolecule(miter); mol != NULL;
640 mol = info_->nextMolecule(miter)) {
641 for (atom = mol->beginAtom(aiter); atom != NULL;
642 atom = mol->nextAtom(aiter)) {
643 ri = atom->getPos();
644 snap->wrapVector(ri);
645 ri *= angstromToM;
646
647 charge = 0.0;
648
649 FixedChargeAdapter fca = FixedChargeAdapter(atom->getAtomType());
650 if (fca.isFixedCharge()) { charge = fca.getCharge(); }
651
653 FluctuatingChargeAdapter(atom->getAtomType());
654 if (fqa.isFluctuatingCharge()) { charge += atom->getFlucQPos(); }
655
656 charge *= chargeToC;
657
658 qpole += 0.5 * charge * outProduct(ri, ri);
659
660 MultipoleAdapter ma = MultipoleAdapter(atom->getAtomType());
661
662 if (ma.isDipole()) {
663 dipole = atom->getDipole() * debyeToCm;
664 qpole += 0.5 * outProduct(dipole, ri);
665 qpole += 0.5 * outProduct(ri, dipole);
666 }
667
668 if (ma.isQuadrupole()) {
669 qpole += atom->getQuadrupole() * debyeToCm * angstromToM;
670 }
671 }
672 }
673
674#ifdef IS_MPI
675 MPI_Allreduce(MPI_IN_PLACE, qpole.getArrayPointer(), 9, MPI_REALTYPE,
676 MPI_SUM, MPI_COMM_WORLD);
677#endif
678
679 snap->setSystemQuadrupole(qpole);
680 }
681
682 return snap->getSystemQuadrupole();
683 }
684
685 // Returns the Heat Flux Vector for the system
686 Vector3d Thermo::getHeatFlux() {
687 Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
688 SimInfo::MoleculeIterator miter;
689 vector<StuntDouble*>::iterator iiter;
690 Molecule* mol;
691 StuntDouble* sd;
692 RigidBody::AtomIterator ai;
693 Atom* atom;
694 Vector3d vel;
695 Vector3d angMom;
696 Mat3x3d I;
697 int i;
698 int j;
699 int k;
700 RealType mass;
701
702 Vector3d x_a;
703 RealType kinetic;
704 RealType potential;
705 RealType eatom;
706 // Convective portion of the heat flux
707 Vector3d heatFluxJc = V3Zero;
708
709 /* Calculate convective portion of the heat flux */
710 for (mol = info_->beginMolecule(miter); mol != NULL;
711 mol = info_->nextMolecule(miter)) {
712 for (sd = mol->beginIntegrableObject(iiter); sd != NULL;
713 sd = mol->nextIntegrableObject(iiter)) {
714 mass = sd->getMass();
715 vel = sd->getVel();
716
717 kinetic = mass * (vel[0] * vel[0] + vel[1] * vel[1] + vel[2] * vel[2]);
718
719 if (sd->isDirectional()) {
720 angMom = sd->getJ();
721 I = sd->getI();
722
723 if (sd->isLinear()) {
724 i = sd->linearAxis();
725 j = (i + 1) % 3;
726 k = (i + 2) % 3;
727 kinetic += angMom[j] * angMom[j] / I(j, j) +
728 angMom[k] * angMom[k] / I(k, k);
729 } else {
730 kinetic += angMom[0] * angMom[0] / I(0, 0) +
731 angMom[1] * angMom[1] / I(1, 1) +
732 angMom[2] * angMom[2] / I(2, 2);
733 }
734 }
735
736 potential = 0.0;
737
738 if (sd->isRigidBody()) {
739 RigidBody* rb = dynamic_cast<RigidBody*>(sd);
740 for (atom = rb->beginAtom(ai); atom != NULL;
741 atom = rb->nextAtom(ai)) {
742 potential += atom->getParticlePot();
743 }
744 } else {
745 potential = sd->getParticlePot();
746 }
747
748 potential *= Constants::energyConvert; // amu A^2/fs^2
749 // The potential may not be a 1/2 factor
750 eatom = (kinetic + potential) / 2.0; // amu A^2/fs^2
751 heatFluxJc[0] += eatom * vel[0]; // amu A^3/fs^3
752 heatFluxJc[1] += eatom * vel[1]; // amu A^3/fs^3
753 heatFluxJc[2] += eatom * vel[2]; // amu A^3/fs^3
754 }
755 }
756
757 /* The J_v vector is reduced in the forceManager so everyone has
758 * the global Jv. Jc is computed over the local atoms and must be
759 * reduced among all processors.
760 */
761#ifdef IS_MPI
762 MPI_Allreduce(MPI_IN_PLACE, &heatFluxJc[0], 3, MPI_REALTYPE, MPI_SUM,
763 MPI_COMM_WORLD);
764#endif
765
766 // (kcal/mol * A/fs) * conversion => (amu A^3)/fs^3
767
768 Vector3d heatFluxJv =
769 currSnapshot->getConductiveHeatFlux() * Constants::energyConvert;
770
771 // Correct for the fact the flux is 1/V (Jc + Jv)
772 return (heatFluxJv + heatFluxJc) / this->getVolume(); // amu / fs^3
773 }
774
777
778 if (!snap->hasCOMvel) {
779 SimInfo::MoleculeIterator i;
780 Molecule* mol;
781
782 Vector3d comVel(0.0);
783 RealType totalMass(0.0);
784
785 for (mol = info_->beginMolecule(i); mol != NULL;
786 mol = info_->nextMolecule(i)) {
787 RealType mass = mol->getMass();
788 totalMass += mass;
789 comVel += mass * mol->getComVel();
790 }
791
792#ifdef IS_MPI
793 MPI_Allreduce(MPI_IN_PLACE, &totalMass, 1, MPI_REALTYPE, MPI_SUM,
794 MPI_COMM_WORLD);
795 MPI_Allreduce(MPI_IN_PLACE, comVel.getArrayPointer(), 3, MPI_REALTYPE,
796 MPI_SUM, MPI_COMM_WORLD);
797#endif
798
799 comVel /= totalMass;
800 snap->setCOMvel(comVel);
801 }
802 return snap->getCOMvel();
803 }
804
807
808 if (!snap->hasCOM) {
809 SimInfo::MoleculeIterator i;
810 Molecule* mol;
811
812 Vector3d com(0.0);
813 RealType totalMass(0.0);
814
815 for (mol = info_->beginMolecule(i); mol != NULL;
816 mol = info_->nextMolecule(i)) {
817 RealType mass = mol->getMass();
818 totalMass += mass;
819 com += mass * mol->getCom();
820 }
821
822#ifdef IS_MPI
823 MPI_Allreduce(MPI_IN_PLACE, &totalMass, 1, MPI_REALTYPE, MPI_SUM,
824 MPI_COMM_WORLD);
825 MPI_Allreduce(MPI_IN_PLACE, com.getArrayPointer(), 3, MPI_REALTYPE,
826 MPI_SUM, MPI_COMM_WORLD);
827#endif
828
829 com /= totalMass;
830 snap->setCOM(com);
831 }
832 return snap->getCOM();
833 }
834
835 /**
836 * Returns center of mass and center of mass velocity in one
837 * function call.
838 */
839 void Thermo::getComAll(Vector3d& com, Vector3d& comVel) {
841
842 if (!(snap->hasCOM && snap->hasCOMvel)) {
843 SimInfo::MoleculeIterator i;
844 Molecule* mol;
845
846 RealType totalMass(0.0);
847
848 com = 0.0;
849 comVel = 0.0;
850
851 for (mol = info_->beginMolecule(i); mol != NULL;
852 mol = info_->nextMolecule(i)) {
853 RealType mass = mol->getMass();
854 totalMass += mass;
855 com += mass * mol->getCom();
856 comVel += mass * mol->getComVel();
857 }
858
859#ifdef IS_MPI
860 MPI_Allreduce(MPI_IN_PLACE, &totalMass, 1, MPI_REALTYPE, MPI_SUM,
861 MPI_COMM_WORLD);
862 MPI_Allreduce(MPI_IN_PLACE, com.getArrayPointer(), 3, MPI_REALTYPE,
863 MPI_SUM, MPI_COMM_WORLD);
864 MPI_Allreduce(MPI_IN_PLACE, comVel.getArrayPointer(), 3, MPI_REALTYPE,
865 MPI_SUM, MPI_COMM_WORLD);
866#endif
867
868 com /= totalMass;
869 comVel /= totalMass;
870 snap->setCOM(com);
871 snap->setCOMvel(comVel);
872 }
873 com = snap->getCOM();
874 comVel = snap->getCOMvel();
875 return;
876 }
877
878 /**
879 * \brief Return inertia tensor for entire system and angular momentum
880 * Vector.
881 *
882 *
883 *
884 * [ Ixx -Ixy -Ixz ]
885 * I =| -Iyx Iyy -Iyz |
886 * [ -Izx -Iyz Izz ]
887 */
889 Vector3d& angularMomentum) {
891 Molecule::IntegrableObjectIterator j;
892 StuntDouble* sd;
893
894 if (!(snap->hasInertiaTensor && snap->hasCOMw)) {
895 RealType xx = 0.0;
896 RealType yy = 0.0;
897 RealType zz = 0.0;
898 RealType xy = 0.0;
899 RealType xz = 0.0;
900 RealType yz = 0.0;
901 Vector3d com(0.0);
902 Vector3d comVel(0.0);
903
904 getComAll(com, comVel);
905
906 SimInfo::MoleculeIterator i;
907 Molecule* mol;
908
909 Vector3d r(0.0);
910 Vector3d v(0.0);
911
912 RealType m = 0.0;
913
914 for (mol = info_->beginMolecule(i); mol != NULL;
915 mol = info_->nextMolecule(i)) {
916 for (sd = mol->beginIntegrableObject(j); sd != NULL;
917 sd = mol->nextIntegrableObject(j)) {
918 r = sd->getPos() - com;
919 v = sd->getVel() - comVel;
920
921 m = sd->getMass();
922
923 // Compute moment of intertia coefficients.
924 xx += r[0] * r[0] * m;
925 yy += r[1] * r[1] * m;
926 zz += r[2] * r[2] * m;
927
928 // compute products of intertia
929 xy += r[0] * r[1] * m;
930 xz += r[0] * r[2] * m;
931 yz += r[1] * r[2] * m;
932
933 angularMomentum += cross(r, v) * m;
934 }
935 }
936
937 inertiaTensor(0, 0) = yy + zz;
938 inertiaTensor(0, 1) = -xy;
939 inertiaTensor(0, 2) = -xz;
940 inertiaTensor(1, 0) = -xy;
941 inertiaTensor(1, 1) = xx + zz;
942 inertiaTensor(1, 2) = -yz;
943 inertiaTensor(2, 0) = -xz;
944 inertiaTensor(2, 1) = -yz;
945 inertiaTensor(2, 2) = xx + yy;
946
947#ifdef IS_MPI
948 MPI_Allreduce(MPI_IN_PLACE, inertiaTensor.getArrayPointer(), 9,
949 MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
950 MPI_Allreduce(MPI_IN_PLACE, angularMomentum.getArrayPointer(), 3,
951 MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
952#endif
953
954 snap->setCOMw(angularMomentum);
955 snap->setInertiaTensor(inertiaTensor);
956 }
957
958 angularMomentum = snap->getCOMw();
959 inertiaTensor = snap->getInertiaTensor();
960
961 return;
962 }
963
966
967 if (!(snap->hasBoundingBox)) {
968 SimInfo::MoleculeIterator i;
969 Molecule::RigidBodyIterator ri;
970 Molecule::AtomIterator ai;
971 Molecule* mol;
972 RigidBody* rb;
973 Atom* atom;
974 Vector3d pos, bMax, bMin;
975 int index = 0;
976
977 for (mol = info_->beginMolecule(i); mol != NULL;
978 mol = info_->nextMolecule(i)) {
979 // change the positions of atoms which belong to the rigidbodies
980 for (rb = mol->beginRigidBody(ri); rb != NULL;
981 rb = mol->nextRigidBody(ri)) {
982 rb->updateAtoms();
983 }
984
985 for (atom = mol->beginAtom(ai); atom != NULL;
986 atom = mol->nextAtom(ai)) {
987 pos = atom->getPos();
988
989 if (index == 0) {
990 bMax = pos;
991 bMin = pos;
992 } else {
993 for (int i = 0; i < 3; i++) {
994 bMax[i] = max(bMax[i], pos[i]);
995 bMin[i] = min(bMin[i], pos[i]);
996 }
997 }
998 index++;
999 }
1000 }
1001
1002#ifdef IS_MPI
1003 MPI_Allreduce(MPI_IN_PLACE, &bMax[0], 3, MPI_REALTYPE, MPI_MAX,
1004 MPI_COMM_WORLD);
1005
1006 MPI_Allreduce(MPI_IN_PLACE, &bMin[0], 3, MPI_REALTYPE, MPI_MIN,
1007 MPI_COMM_WORLD);
1008#endif
1009 Mat3x3d bBox = Mat3x3d(0.0);
1010 for (int i = 0; i < 3; i++) {
1011 bBox(i, i) = bMax[i] - bMin[i];
1012 }
1013 snap->setBoundingBox(bBox);
1014 }
1015
1016 return snap->getBoundingBox();
1017 }
1018
1019 // Returns the angular momentum of the system
1022
1023 if (!snap->hasCOMw) {
1024 Vector3d com(0.0);
1025 Vector3d comVel(0.0);
1026 Vector3d angularMomentum(0.0);
1027
1028 getComAll(com, comVel);
1029
1030 SimInfo::MoleculeIterator i;
1031 Molecule* mol;
1032
1033 Vector3d thisr(0.0);
1034 Vector3d thisp(0.0);
1035
1036 RealType thisMass;
1037
1038 for (mol = info_->beginMolecule(i); mol != NULL;
1039 mol = info_->nextMolecule(i)) {
1040 thisMass = mol->getMass();
1041 thisr = mol->getCom() - com;
1042 thisp = (mol->getComVel() - comVel) * thisMass;
1043
1044 angularMomentum += cross(thisr, thisp);
1045 }
1046
1047#ifdef IS_MPI
1048 MPI_Allreduce(MPI_IN_PLACE, angularMomentum.getArrayPointer(), 3,
1049 MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
1050#endif
1051
1052 snap->setCOMw(angularMomentum);
1053 }
1054
1055 return snap->getCOMw();
1056 }
1057
1058 /**
1059 * Returns the Volume of the system based on a ellipsoid with
1060 * semi-axes based on the radius of gyration V=4/3*Pi*R_1*R_2*R_3
1061 * where R_i are related to the principle inertia moments
1062 * R_i = sqrt(C*I_i/N), this reduces to
1063 * V = 4/3*Pi*(C/N)^3/2*sqrt(det(I)).
1064 * See S.E. Baltazar et. al. Comp. Mat. Sci. 37 (2006) 526-536.
1065 */
1068
1069 if (!snap->hasGyrationalVolume) {
1070 Mat3x3d intTensor;
1071 RealType det;
1072 Vector3d dummyAngMom;
1073 RealType sysconstants;
1074 RealType geomCnst;
1075 RealType volume;
1076
1077 geomCnst = 3.0 / 2.0;
1078 /* Get the inertial tensor and angular momentum for free*/
1079 getInertiaTensor(intTensor, dummyAngMom);
1080
1081 det = intTensor.determinant();
1082 sysconstants =
1083 geomCnst / (RealType)(info_->getNGlobalIntegrableObjects());
1084 volume =
1085 4.0 / 3.0 * Constants::PI * pow(sysconstants, geomCnst) * sqrt(det);
1086
1087 snap->setGyrationalVolume(volume);
1088 }
1089 return snap->getGyrationalVolume();
1090 }
1091
1092 void Thermo::getGyrationalVolume(RealType& volume, RealType& detI) {
1094
1095 if (!(snap->hasInertiaTensor && snap->hasGyrationalVolume)) {
1096 Mat3x3d intTensor;
1097 Vector3d dummyAngMom;
1098 RealType sysconstants;
1099 RealType geomCnst;
1100
1101 geomCnst = 3.0 / 2.0;
1102 /* Get the inertia tensor and angular momentum for free*/
1103 this->getInertiaTensor(intTensor, dummyAngMom);
1104
1105 detI = intTensor.determinant();
1106 sysconstants =
1107 geomCnst / (RealType)(info_->getNGlobalIntegrableObjects());
1108 volume =
1109 4.0 / 3.0 * Constants::PI * pow(sysconstants, geomCnst) * sqrt(detI);
1110 snap->setGyrationalVolume(volume);
1111 } else {
1112 volume = snap->getGyrationalVolume();
1113 detI = snap->getInertiaTensor().determinant();
1114 }
1115 return;
1116 }
1117
1118 RealType Thermo::getTaggedAtomPairDistance() {
1119 Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
1120 Globals* simParams = info_->getSimParams();
1121
1122 if (simParams->haveTaggedAtomPair() &&
1123 simParams->havePrintTaggedPairDistance()) {
1124 if (simParams->getPrintTaggedPairDistance()) {
1125 pair<int, int> tap = simParams->getTaggedAtomPair();
1126 Vector3d pos1, pos2, rab;
1127
1128#ifdef IS_MPI
1129 int mol1 = info_->getGlobalMolMembership(tap.first);
1130 int mol2 = info_->getGlobalMolMembership(tap.second);
1131
1132 int proc1 = info_->getMolToProc(mol1);
1133 int proc2 = info_->getMolToProc(mol2);
1134
1135 RealType data[3];
1136 if (proc1 == worldRank) {
1137 StuntDouble* sd1 = info_->getIOIndexToIntegrableObject(tap.first);
1138 pos1 = sd1->getPos();
1139 data[0] = pos1.x();
1140 data[1] = pos1.y();
1141 data[2] = pos1.z();
1142 MPI_Bcast(data, 3, MPI_REALTYPE, proc1, MPI_COMM_WORLD);
1143 } else {
1144 MPI_Bcast(data, 3, MPI_REALTYPE, proc1, MPI_COMM_WORLD);
1145 pos1 = Vector3d(data);
1146 }
1147
1148 if (proc2 == worldRank) {
1149 StuntDouble* sd2 = info_->getIOIndexToIntegrableObject(tap.second);
1150 pos2 = sd2->getPos();
1151 data[0] = pos2.x();
1152 data[1] = pos2.y();
1153 data[2] = pos2.z();
1154 MPI_Bcast(data, 3, MPI_REALTYPE, proc2, MPI_COMM_WORLD);
1155 } else {
1156 MPI_Bcast(data, 3, MPI_REALTYPE, proc2, MPI_COMM_WORLD);
1157 pos2 = Vector3d(data);
1158 }
1159#else
1160 StuntDouble* at1 = info_->getIOIndexToIntegrableObject(tap.first);
1161 StuntDouble* at2 = info_->getIOIndexToIntegrableObject(tap.second);
1162 pos1 = at1->getPos();
1163 pos2 = at2->getPos();
1164#endif
1165 rab = pos2 - pos1;
1166 currSnapshot->wrapVector(rab);
1167 return rab.length();
1168 }
1169 return 0.0;
1170 }
1171 return 0.0;
1172 }
1173
1174 RealType Thermo::getHullVolume() {
1175#ifdef HAVE_QHULL
1177 if (!snap->hasHullVolume) {
1178 Hull* surfaceMesh_;
1179
1180 Globals* simParams = info_->getSimParams();
1181 const std::string ht = simParams->getHULL_Method();
1182
1183 if (ht == "Convex") {
1184 surfaceMesh_ = new ConvexHull();
1185 } else if (ht == "AlphaShape") {
1186 surfaceMesh_ = new AlphaHull(simParams->getAlpha());
1187 } else {
1188 return 0.0;
1189 }
1190
1191 // Build a vector of stunt doubles to determine if they are
1192 // surface atoms
1193 std::vector<StuntDouble*> localSites_;
1194 Molecule* mol;
1195 StuntDouble* sd;
1196 SimInfo::MoleculeIterator i;
1197 Molecule::IntegrableObjectIterator j;
1198
1199 for (mol = info_->beginMolecule(i); mol != NULL;
1200 mol = info_->nextMolecule(i)) {
1201 for (sd = mol->beginIntegrableObject(j); sd != NULL;
1202 sd = mol->nextIntegrableObject(j)) {
1203 localSites_.push_back(sd);
1204 }
1205 }
1206
1207 // Compute surface Mesh
1208 surfaceMesh_->computeHull(localSites_);
1209 snap->setHullVolume(surfaceMesh_->getVolume());
1210
1211 delete surfaceMesh_;
1212 }
1213
1214 return snap->getHullVolume();
1215#else
1216 return 0.0;
1217#endif
1218 }
1219
1220} // namespace OpenMD
Compute alpha complex or alpha shape.
Definition AlphaHull.hpp:79
AtomType is what OpenMD looks to for unchanging data about an atom.
Definition AtomType.hpp:66
RealType getMass()
get total mass of this molecule
Definition Molecule.cpp:266
Vector3d getComVel()
Returns the velocity of center of mass of this molecule.
Definition Molecule.cpp:351
Vector3d getCom()
Returns the current center of mass position of this molecule.
Definition Molecule.cpp:298
Real * getArrayPointer()
Returns the pointer of internal array.
int getNdf()
Returns the number of degrees of freedom.
Definition SimInfo.hpp:220
int getNGlobalIntegrableObjects()
Returns the total number of integrable objects (total number of rigid bodies plus the total number of...
Definition SimInfo.hpp:139
int getMolToProc(int globalIndex)
Finds the processor where a molecule resides.
Definition SimInfo.hpp:667
Molecule * beginMolecule(MoleculeIterator &i)
Returns the first molecule in this SimInfo and intialize the iterator.
Definition SimInfo.cpp:240
int getNFluctuatingCharges()
Returns the total number of fluctuating charges that are present.
Definition SimInfo.hpp:217
Molecule * nextMolecule(MoleculeIterator &i)
Returns the next avaliable Molecule based on the iterator.
Definition SimInfo.cpp:245
SnapshotManager * getSnapshotManager()
Returns the snapshot manager.
Definition SimInfo.hpp:248
StuntDouble * getIOIndexToIntegrableObject(int index)
return an integral objects by its global index.
Definition SimInfo.cpp:1027
AtomTypeSet getSimulatedAtomTypes()
Returns the set of atom types present in this simulation.
Definition SimInfo.cpp:712
The Snapshot class is a repository storing dynamic data during a Simulation.
Definition Snapshot.hpp:147
Mat3x3d getBoundingBox()
Returns the Bounding Box.
Definition Snapshot.cpp:281
void setBoundingBox(const Mat3x3d &m)
Sets the Bounding Box.
Definition Snapshot.cpp:284
void wrapVector(Vector3d &v)
Wrapping the vector according to periodic boundary condition.
Definition Snapshot.cpp:337
Snapshot * getCurrentSnapshot()
Returns the pointer of current snapshot.
Real determinant() const
Returns the determinant of this matrix.
"Don't move, or you're dead! Stand up! Captain, we've got them!"
Vector3d getVel()
Returns the current velocity of this stuntDouble.
int linearAxis()
Returns the linear axis of the rigidbody, atom and directional atom will always return -1.
RealType getMass()
Returns the mass of this stuntDouble.
virtual Mat3x3d getI()=0
Returns the inertia tensor of this stuntDouble.
bool isLinear()
Tests the if this stuntDouble is a linear rigidbody.
Vector3d getPos()
Returns the current position of this stuntDouble.
RealType getParticlePot()
Returns the current particlePot of this stuntDouble.
bool isRigidBody()
Tests if this stuntDouble is a rigid body.
Vector3d getJ()
Returns the current angular momentum of this stuntDouble (body -fixed).
bool isDirectional()
Tests if this stuntDouble is a directional one.
Mat3x3d getPressureTensor()
gives the pressure tensor in amu*fs^-2*Ang^-1
Definition Thermo.cpp:443
Vector3d getSystemDipole()
accumulate and return the simulation box dipole moment in C*m
Definition Thermo.cpp:525
void getComAll(Vector3d &com, Vector3d &comVel)
Returns the center of the mass and Center of Mass velocity of the whole system.
Definition Thermo.cpp:839
void getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum)
Returns the inertia tensor and the total angular momentum for for the entire system.
Definition Thermo.cpp:888
RealType getGyrationalVolume()
Returns volume of system as estimated by an ellipsoid defined by the radii of gyration.
Definition Thermo.cpp:1066
Vector3d getCom()
Returns the center of the mass of the whole system.
Definition Thermo.cpp:805
Mat3x3d getSystemQuadrupole()
accumulate and return the simulation box dipole moment in debye Angstroms
Definition Thermo.cpp:622
Vector3d getAngularMomentum()
Returns system angular momentum.
Definition Thermo.cpp:1020
Vector3d getComVel()
Returns the velocity of center of mass of the whole system.
Definition Thermo.cpp:775
Mat3x3d getBoundingBox()
Returns the Axis-aligned bounding box for the current system.
Definition Thermo.cpp:964
Real length()
Returns the length of this vector.
Definition Vector.hpp:393
Real * getArrayPointer()
Returns the pointer of internal array.
Definition Vector.hpp:171
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
Vector3< Real > cross(const Vector3< Real > &v1, const Vector3< Real > &v2)
Returns the cross product of two Vectors.
Definition Vector3.hpp:136