ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-3.0/src/integrators/NPTxyz.cpp
Revision: 1883
Committed: Mon Dec 13 22:30:27 2004 UTC (19 years, 7 months ago) by tim
File size: 3740 byte(s)
Log Message:
MPI version is built

File Contents

# Content
1 #include "brains/SimInfo.hpp"
2 #include "brains/Thermo.hpp"
3 #include "integrators/IntegratorCreator.hpp"
4 #include "integrators/NPTxyz.hpp"
5 #include "primitives/Molecule.hpp"
6 #include "utils/OOPSEConstant.hpp"
7 #include "utils/simError.h"
8
9 // Basic non-isotropic thermostating and barostating via the Melchionna
10 // modification of the Hoover algorithm:
11 //
12 // Melchionna, S., Ciccotti, G., and Holian, B. L., 1993,
13 // Molec. Phys., 78, 533.
14 //
15 // and
16 //
17 // Hoover, W. G., 1986, Phys. Rev. A, 34, 2499.
18
19 namespace oopse {
20
21
22 double NPTxyz::calcConservedQuantity(){
23
24 // We need NkBT a lot, so just set it here: This is the RAW number
25 // of integrableObjects, so no subtraction or addition of constraints or
26 // orientational degrees of freedom:
27 NkBT = info_->getNGlobalIntegrableObjects()*OOPSEConstant::kB *targetTemp;
28
29 // fkBT is used because the thermostat operates on more degrees of freedom
30 // than the barostat (when there are particles with orientational degrees
31 // of freedom).
32 fkBT = info_->getNdf()*OOPSEConstant::kB *targetTemp;
33
34 double conservedQuantity;
35 double totalEnergy;
36 double thermostat_kinetic;
37 double thermostat_potential;
38 double barostat_kinetic;
39 double barostat_potential;
40 double trEta;
41
42 totalEnergy = thermo.getTotalE();
43
44 thermostat_kinetic = fkBT * tt2 * chi * chi /(2.0 * OOPSEConstant::energyConvert);
45
46 thermostat_potential = fkBT* integralOfChidt / OOPSEConstant::energyConvert;
47
48 SquareMatrix<double, 3> tmp = eta.transpose() * eta;
49 trEta = tmp.trace();
50
51 barostat_kinetic = NkBT * tb2 * trEta /(2.0 * OOPSEConstant::energyConvert);
52
53 barostat_potential = (targetPressure * thermo.getVolume() / OOPSEConstant::pressureConvert) /OOPSEConstant::energyConvert;
54
55 conservedQuantity = totalEnergy + thermostat_kinetic + thermostat_potential +
56 barostat_kinetic + barostat_potential;
57
58
59 return conservedQuantity;
60
61 }
62
63
64 void NPTxyz::scaleSimBox(){
65
66 int i,j,k;
67 Mat3x3d scaleMat;
68 double eta2ij, scaleFactor;
69 double bigScale, smallScale, offDiagMax;
70 Mat3x3d hm;
71 Mat3x3d hmnew;
72
73
74
75 // Scale the box after all the positions have been moved:
76
77 // Use a taylor expansion for eta products: Hmat = Hmat . exp(dt * etaMat)
78 // Hmat = Hmat . ( Ident + dt * etaMat + dt^2 * etaMat*etaMat / 2)
79
80 bigScale = 1.0;
81 smallScale = 1.0;
82 offDiagMax = 0.0;
83
84 for(i=0; i<3; i++){
85 for(j=0; j<3; j++){
86 scaleMat(i, j) = 0.0;
87 if(i==j) {
88 scaleMat(i, j) = 1.0;
89 }
90 }
91 }
92
93 for(i=0;i<3;i++){
94
95 // calculate the scaleFactors
96
97 scaleFactor = exp(dt*eta(i, i));
98
99 scaleMat(i, i) = scaleFactor;
100
101 if (scaleMat(i, i) > bigScale) {
102 bigScale = scaleMat(i, i);
103 }
104
105 if (scaleMat(i, i) < smallScale) {
106 smallScale = scaleMat(i, i);
107 }
108 }
109
110 if ((bigScale > 1.1) || (smallScale < 0.9)) {
111 sprintf( painCave.errMsg,
112 "NPTxyz error: Attempting a Box scaling of more than 10 percent.\n"
113 " Check your tauBarostat, as it is probably too small!\n\n"
114 " scaleMat = [%lf\t%lf\t%lf]\n"
115 " [%lf\t%lf\t%lf]\n"
116 " [%lf\t%lf\t%lf]\n",
117 scaleMat(0, 0),scaleMat(0, 1),scaleMat(0, 2),
118 scaleMat(1, 0),scaleMat(1, 1),scaleMat(1, 2),
119 scaleMat(2, 0),scaleMat(2, 1),scaleMat(2, 2));
120 painCave.isFatal = 1;
121 simError();
122 } else {
123
124 Mat3x3d hmat = currentSnapshot_->getHmat();
125 hmat = hmat *scaleMat;
126 currentSnapshot_->setHmat(hmat);
127 }
128 }
129
130 void NPTxyz::loadEta() {
131 eta= currentSnapshot_->getEta();
132 }
133
134 }