| 1 |
#include <cmath> |
| 2 |
#include <iostream> |
| 3 |
using namespace std; |
| 4 |
|
| 5 |
#ifdef IS_MPI |
| 6 |
#include <mpi.h> |
| 7 |
#endif //is_mpi |
| 8 |
|
| 9 |
#include "Thermo.hpp" |
| 10 |
#include "SRI.hpp" |
| 11 |
#include "Integrator.hpp" |
| 12 |
#include "simError.h" |
| 13 |
|
| 14 |
#ifdef IS_MPI |
| 15 |
#define __C |
| 16 |
#include "mpiSimulation.hpp" |
| 17 |
#endif // is_mpi |
| 18 |
|
| 19 |
|
| 20 |
#define BASE_SEED 123456789 |
| 21 |
|
| 22 |
Thermo::Thermo( SimInfo* the_entry_plug ) { |
| 23 |
entry_plug = the_entry_plug; |
| 24 |
int baseSeed = BASE_SEED; |
| 25 |
|
| 26 |
gaussStream = new gaussianSPRNG( baseSeed ); |
| 27 |
} |
| 28 |
|
| 29 |
Thermo::~Thermo(){ |
| 30 |
delete gaussStream; |
| 31 |
} |
| 32 |
|
| 33 |
double Thermo::getKinetic(){ |
| 34 |
|
| 35 |
const double e_convert = 4.184E-4; // convert kcal/mol -> (amu A^2)/fs^2 |
| 36 |
double vx2, vy2, vz2; |
| 37 |
double kinetic, v_sqr; |
| 38 |
int kl; |
| 39 |
double jx2, jy2, jz2; // the square of the angular momentums |
| 40 |
|
| 41 |
DirectionalAtom *dAtom; |
| 42 |
|
| 43 |
int n_atoms; |
| 44 |
double kinetic_global; |
| 45 |
Atom** atoms; |
| 46 |
|
| 47 |
|
| 48 |
n_atoms = entry_plug->n_atoms; |
| 49 |
atoms = entry_plug->atoms; |
| 50 |
|
| 51 |
kinetic = 0.0; |
| 52 |
kinetic_global = 0.0; |
| 53 |
for( kl=0; kl < n_atoms; kl++ ){ |
| 54 |
|
| 55 |
vx2 = atoms[kl]->get_vx() * atoms[kl]->get_vx(); |
| 56 |
vy2 = atoms[kl]->get_vy() * atoms[kl]->get_vy(); |
| 57 |
vz2 = atoms[kl]->get_vz() * atoms[kl]->get_vz(); |
| 58 |
|
| 59 |
v_sqr = vx2 + vy2 + vz2; |
| 60 |
kinetic += atoms[kl]->getMass() * v_sqr; |
| 61 |
|
| 62 |
if( atoms[kl]->isDirectional() ){ |
| 63 |
|
| 64 |
dAtom = (DirectionalAtom *)atoms[kl]; |
| 65 |
|
| 66 |
jx2 = dAtom->getJx() * dAtom->getJx(); |
| 67 |
jy2 = dAtom->getJy() * dAtom->getJy(); |
| 68 |
jz2 = dAtom->getJz() * dAtom->getJz(); |
| 69 |
|
| 70 |
kinetic += (jx2 / dAtom->getIxx()) + (jy2 / dAtom->getIyy()) |
| 71 |
+ (jz2 / dAtom->getIzz()); |
| 72 |
} |
| 73 |
} |
| 74 |
#ifdef IS_MPI |
| 75 |
MPI_Allreduce(&kinetic,&kinetic_global,1,MPI_DOUBLE, |
| 76 |
MPI_SUM, MPI_COMM_WORLD); |
| 77 |
kinetic = kinetic_global; |
| 78 |
#endif //is_mpi |
| 79 |
|
| 80 |
kinetic = kinetic * 0.5 / e_convert; |
| 81 |
|
| 82 |
return kinetic; |
| 83 |
} |
| 84 |
|
| 85 |
double Thermo::getPotential(){ |
| 86 |
|
| 87 |
double potential_local; |
| 88 |
double potential; |
| 89 |
int el, nSRI; |
| 90 |
Molecule* molecules; |
| 91 |
|
| 92 |
molecules = entry_plug->molecules; |
| 93 |
nSRI = entry_plug->n_SRI; |
| 94 |
|
| 95 |
potential_local = 0.0; |
| 96 |
potential = 0.0; |
| 97 |
potential_local += entry_plug->lrPot; |
| 98 |
|
| 99 |
for( el=0; el<entry_plug->n_mol; el++ ){ |
| 100 |
potential_local += molecules[el].getPotential(); |
| 101 |
} |
| 102 |
|
| 103 |
// Get total potential for entire system from MPI. |
| 104 |
#ifdef IS_MPI |
| 105 |
MPI_Allreduce(&potential_local,&potential,1,MPI_DOUBLE, |
| 106 |
MPI_SUM, MPI_COMM_WORLD); |
| 107 |
#else |
| 108 |
potential = potential_local; |
| 109 |
#endif // is_mpi |
| 110 |
|
| 111 |
#ifdef IS_MPI |
| 112 |
/* |
| 113 |
std::cerr << "node " << worldRank << ": after pot = " << potential << "\n"; |
| 114 |
*/ |
| 115 |
#endif |
| 116 |
|
| 117 |
return potential; |
| 118 |
} |
| 119 |
|
| 120 |
double Thermo::getTotalE(){ |
| 121 |
|
| 122 |
double total; |
| 123 |
|
| 124 |
total = this->getKinetic() + this->getPotential(); |
| 125 |
return total; |
| 126 |
} |
| 127 |
|
| 128 |
double Thermo::getTemperature(){ |
| 129 |
|
| 130 |
const double kb = 1.9872179E-3; // boltzman's constant in kcal/(mol K) |
| 131 |
double temperature; |
| 132 |
|
| 133 |
temperature = ( 2.0 * this->getKinetic() ) / ((double)entry_plug->ndf * kb ); |
| 134 |
return temperature; |
| 135 |
} |
| 136 |
|
| 137 |
double Thermo::getPressure() { |
| 138 |
// returns the pressure in units of atm |
| 139 |
// Relies on the calculation of the full molecular pressure tensor |
| 140 |
|
| 141 |
const double p_convert = 1.63882576e8; |
| 142 |
double press[9]; |
| 143 |
double pressure; |
| 144 |
|
| 145 |
this->getPressureTensor(press); |
| 146 |
|
| 147 |
pressure = p_convert * (press[0] + press[4] + press[8]) / 3.0; |
| 148 |
|
| 149 |
return pressure; |
| 150 |
} |
| 151 |
|
| 152 |
|
| 153 |
void Thermo::getPressureTensor(double press[9]){ |
| 154 |
// returns pressure tensor in units amu*fs^-2*Ang^-1 |
| 155 |
// routine derived via viral theorem description in: |
| 156 |
// Paci, E. and Marchi, M. J.Phys.Chem. 1996, 100, 4314-4322 |
| 157 |
|
| 158 |
const double e_convert = 4.184e-4; |
| 159 |
|
| 160 |
double molmass, volume; |
| 161 |
double vcom[3]; |
| 162 |
double p_local[9], p_global[9]; |
| 163 |
double theBox[3]; |
| 164 |
double* tau; |
| 165 |
int i, nMols; |
| 166 |
Molecule* molecules; |
| 167 |
|
| 168 |
nMols = entry_plug->n_mol; |
| 169 |
molecules = entry_plug->molecules; |
| 170 |
tau = entry_plug->tau; |
| 171 |
|
| 172 |
// use velocities of molecular centers of mass and molecular masses: |
| 173 |
for (i=0; i < 9; i++) { |
| 174 |
p_local[i] = 0.0; |
| 175 |
p_global[i] = 0.0; |
| 176 |
} |
| 177 |
|
| 178 |
for (i=0; i < nMols; i++) { |
| 179 |
molmass = molecules[i].getCOMvel(vcom); |
| 180 |
|
| 181 |
p_local[0] += molmass * (vcom[0] * vcom[0]); |
| 182 |
p_local[1] += molmass * (vcom[0] * vcom[1]); |
| 183 |
p_local[2] += molmass * (vcom[0] * vcom[2]); |
| 184 |
p_local[3] += molmass * (vcom[1] * vcom[0]); |
| 185 |
p_local[4] += molmass * (vcom[1] * vcom[1]); |
| 186 |
p_local[5] += molmass * (vcom[1] * vcom[2]); |
| 187 |
p_local[6] += molmass * (vcom[2] * vcom[0]); |
| 188 |
p_local[7] += molmass * (vcom[2] * vcom[1]); |
| 189 |
p_local[8] += molmass * (vcom[2] * vcom[2]); |
| 190 |
} |
| 191 |
|
| 192 |
// Get total for entire system from MPI. |
| 193 |
|
| 194 |
#ifdef IS_MPI |
| 195 |
MPI_Allreduce(p_local,p_global,9,MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD); |
| 196 |
#else |
| 197 |
for (i=0; i<9; i++) { |
| 198 |
p_global[i] = p_local[i]; |
| 199 |
} |
| 200 |
#endif // is_mpi |
| 201 |
|
| 202 |
entry_plug->getBox(theBox); |
| 203 |
|
| 204 |
volume = theBox[0] * theBox[1] * theBox[2]; |
| 205 |
|
| 206 |
for(i=0; i<9; i++) { |
| 207 |
press[i] = (p_global[i] - tau[i]*e_convert) / volume; |
| 208 |
} |
| 209 |
} |
| 210 |
|
| 211 |
void Thermo::velocitize() { |
| 212 |
|
| 213 |
double x,y; |
| 214 |
double vx, vy, vz; |
| 215 |
double jx, jy, jz; |
| 216 |
int i, vr, vd; // velocity randomizer loop counters |
| 217 |
double vdrift[3]; |
| 218 |
double vbar; |
| 219 |
const double kb = 8.31451e-7; // kb in amu, angstroms, fs, etc. |
| 220 |
double av2; |
| 221 |
double kebar; |
| 222 |
int n_atoms; |
| 223 |
Atom** atoms; |
| 224 |
DirectionalAtom* dAtom; |
| 225 |
double temperature; |
| 226 |
int n_oriented; |
| 227 |
int n_constraints; |
| 228 |
|
| 229 |
atoms = entry_plug->atoms; |
| 230 |
n_atoms = entry_plug->n_atoms; |
| 231 |
temperature = entry_plug->target_temp; |
| 232 |
n_oriented = entry_plug->n_oriented; |
| 233 |
n_constraints = entry_plug->n_constraints; |
| 234 |
|
| 235 |
kebar = kb * temperature * (double)entry_plug->ndf / |
| 236 |
( 2.0 * (double)entry_plug->ndfRaw ); |
| 237 |
|
| 238 |
for(vr = 0; vr < n_atoms; vr++){ |
| 239 |
|
| 240 |
// uses equipartition theory to solve for vbar in angstrom/fs |
| 241 |
|
| 242 |
av2 = 2.0 * kebar / atoms[vr]->getMass(); |
| 243 |
vbar = sqrt( av2 ); |
| 244 |
|
| 245 |
// vbar = sqrt( 8.31451e-7 * temperature / atoms[vr]->getMass() ); |
| 246 |
|
| 247 |
// picks random velocities from a gaussian distribution |
| 248 |
// centered on vbar |
| 249 |
|
| 250 |
vx = vbar * gaussStream->getGaussian(); |
| 251 |
vy = vbar * gaussStream->getGaussian(); |
| 252 |
vz = vbar * gaussStream->getGaussian(); |
| 253 |
|
| 254 |
atoms[vr]->set_vx( vx ); |
| 255 |
atoms[vr]->set_vy( vy ); |
| 256 |
atoms[vr]->set_vz( vz ); |
| 257 |
} |
| 258 |
|
| 259 |
// Get the Center of Mass drift velocity. |
| 260 |
|
| 261 |
getCOMVel(vdrift); |
| 262 |
|
| 263 |
// Corrects for the center of mass drift. |
| 264 |
// sums all the momentum and divides by total mass. |
| 265 |
|
| 266 |
for(vd = 0; vd < n_atoms; vd++){ |
| 267 |
|
| 268 |
vx = atoms[vd]->get_vx(); |
| 269 |
vy = atoms[vd]->get_vy(); |
| 270 |
vz = atoms[vd]->get_vz(); |
| 271 |
|
| 272 |
vx -= vdrift[0]; |
| 273 |
vy -= vdrift[1]; |
| 274 |
vz -= vdrift[2]; |
| 275 |
|
| 276 |
atoms[vd]->set_vx(vx); |
| 277 |
atoms[vd]->set_vy(vy); |
| 278 |
atoms[vd]->set_vz(vz); |
| 279 |
} |
| 280 |
if( n_oriented ){ |
| 281 |
|
| 282 |
for( i=0; i<n_atoms; i++ ){ |
| 283 |
|
| 284 |
if( atoms[i]->isDirectional() ){ |
| 285 |
|
| 286 |
dAtom = (DirectionalAtom *)atoms[i]; |
| 287 |
|
| 288 |
vbar = sqrt( 2.0 * kebar * dAtom->getIxx() ); |
| 289 |
jx = vbar * gaussStream->getGaussian(); |
| 290 |
|
| 291 |
vbar = sqrt( 2.0 * kebar * dAtom->getIyy() ); |
| 292 |
jy = vbar * gaussStream->getGaussian(); |
| 293 |
|
| 294 |
vbar = sqrt( 2.0 * kebar * dAtom->getIzz() ); |
| 295 |
jz = vbar * gaussStream->getGaussian(); |
| 296 |
|
| 297 |
dAtom->setJx( jx ); |
| 298 |
dAtom->setJy( jy ); |
| 299 |
dAtom->setJz( jz ); |
| 300 |
} |
| 301 |
} |
| 302 |
} |
| 303 |
} |
| 304 |
|
| 305 |
void Thermo::getCOMVel(double vdrift[3]){ |
| 306 |
|
| 307 |
double mtot, mtot_local; |
| 308 |
double vdrift_local[3]; |
| 309 |
int vd, n_atoms; |
| 310 |
Atom** atoms; |
| 311 |
|
| 312 |
// We are very careless here with the distinction between n_atoms and n_local |
| 313 |
// We should really fix this before someone pokes an eye out. |
| 314 |
|
| 315 |
n_atoms = entry_plug->n_atoms; |
| 316 |
atoms = entry_plug->atoms; |
| 317 |
|
| 318 |
mtot_local = 0.0; |
| 319 |
vdrift_local[0] = 0.0; |
| 320 |
vdrift_local[1] = 0.0; |
| 321 |
vdrift_local[2] = 0.0; |
| 322 |
|
| 323 |
for(vd = 0; vd < n_atoms; vd++){ |
| 324 |
|
| 325 |
vdrift_local[0] += atoms[vd]->get_vx() * atoms[vd]->getMass(); |
| 326 |
vdrift_local[1] += atoms[vd]->get_vy() * atoms[vd]->getMass(); |
| 327 |
vdrift_local[2] += atoms[vd]->get_vz() * atoms[vd]->getMass(); |
| 328 |
|
| 329 |
mtot_local += atoms[vd]->getMass(); |
| 330 |
} |
| 331 |
|
| 332 |
#ifdef IS_MPI |
| 333 |
MPI_Allreduce(&mtot_local,&mtot,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
| 334 |
MPI_Allreduce(vdrift_local,vdrift,3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
| 335 |
#else |
| 336 |
mtot = mtot_local; |
| 337 |
for(vd = 0; vd < 3; vd++) { |
| 338 |
vdrift[vd] = vdrift_local[vd]; |
| 339 |
} |
| 340 |
#endif |
| 341 |
|
| 342 |
for (vd = 0; vd < 3; vd++) { |
| 343 |
vdrift[vd] = vdrift[vd] / mtot; |
| 344 |
} |
| 345 |
|
| 346 |
} |
| 347 |
|