ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/Thermo.cpp
(Generate patch)

Comparing trunk/OOPSE/libmdtools/Thermo.cpp (file contents):
Revision 445 by gezelter, Thu Apr 3 19:58:24 2003 UTC vs.
Revision 458 by gezelter, Fri Apr 4 19:47:19 2003 UTC

# Line 4 | Line 4 | using namespace std;
4  
5   #ifdef IS_MPI
6   #include <mpi.h>
7 #include <mpi++.h>
7   #endif //is_mpi
8  
9   #include "Thermo.hpp"
# Line 73 | Line 72 | double Thermo::getKinetic(){
72      }
73    }
74   #ifdef IS_MPI
75 <  MPI::COMM_WORLD.Allreduce(&kinetic,&kinetic_global,1,MPI_DOUBLE,MPI_SUM);
75 >  MPI_Allreduce(&kinetic,&kinetic_global,1,MPI_DOUBLE,
76 >                MPI_SUM, MPI_COMM_WORLD);
77    kinetic = kinetic_global;
78   #endif //is_mpi
79  
# Line 100 | Line 100 | double Thermo::getPotential(){
100      potential_local += molecules[el].getPotential();
101    }
102  
103 #ifdef IS_MPI
104  /*
105  std::cerr << "node " << worldRank << ": before LONG RANGE pot = " << entry_plug->lrPot
106            << "; pot_local = " << potential_local
107            << "; pot = " << potential << "\n";
108  */
109 #endif
110
103    // Get total potential for entire system from MPI.
104   #ifdef IS_MPI
105 <  MPI::COMM_WORLD.Allreduce(&potential_local,&potential,1,MPI_DOUBLE,MPI_SUM);
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
# Line 136 | Line 129 | double Thermo::getTemperature(){
129  
130    const double kb = 1.9872179E-3; // boltzman's constant in kcal/(mol K)
131    double temperature;
139  int ndf_local, ndf;
132    
133 <  ndf_local = 3 * entry_plug->n_atoms + 3 * entry_plug->n_oriented
142 <    - entry_plug->n_constraints;
143 <
144 < #ifdef IS_MPI
145 <  MPI::COMM_WORLD.Allreduce(&ndf_local,&ndf,1,MPI_INT,MPI_SUM);
146 < #else
147 <  ndf = ndf_local;
148 < #endif
149 <
150 <  ndf = ndf - 3;
151 <  
152 <  temperature = ( 2.0 * this->getKinetic() ) / ( ndf * kb );
133 >  temperature = ( 2.0 * this->getKinetic() ) / ((double)entry_plug->ndf * kb );
134    return temperature;
135   }
136  
# Line 172 | Line 153 | void Thermo::velocitize() {
153    const double kb = 8.31451e-7; // kb in amu, angstroms, fs, etc.
154    double av2;
155    double kebar;
175  int ndf, ndf_local; // number of degrees of freedom
176  int ndfRaw, ndfRaw_local; // the raw number of degrees of freedom
156    int n_atoms;
157    Atom** atoms;
158    DirectionalAtom* dAtom;
# Line 187 | Line 166 | void Thermo::velocitize() {
166    n_oriented    = entry_plug->n_oriented;
167    n_constraints = entry_plug->n_constraints;
168    
169 <  // Raw degrees of freedom that we have to set
170 <  ndfRaw_local = 3 * entry_plug->n_atoms + 3 * entry_plug->n_oriented;
192 <
193 <  // Degrees of freedom that can contain kinetic energy
194 <  ndf_local = 3 * entry_plug->n_atoms + 3 * entry_plug->n_oriented
195 <    - entry_plug->n_constraints;
169 >  kebar = kb * temperature * (double)entry_plug->ndf /
170 >    ( 2.0 * (double)entry_plug->ndfRaw );
171    
197 #ifdef IS_MPI
198  MPI::COMM_WORLD.Allreduce(&ndf_local,&ndf,1,MPI_INT,MPI_SUM);
199  MPI::COMM_WORLD.Allreduce(&ndfRaw_local,&ndfRaw,1,MPI_INT,MPI_SUM);
200 #else
201  ndfRaw = ndfRaw_local;
202  ndf = ndf_local;
203 #endif
204  ndf = ndf - 3;
205
206  kebar = kb * temperature * (double)ndf / ( 2.0 * (double)ndfRaw );
207  
172    for(vr = 0; vr < n_atoms; vr++){
173      
174      // uses equipartition theory to solve for vbar in angstrom/fs
# Line 260 | Line 224 | void Thermo::velocitize() {
224  
225          vbar = sqrt( 2.0 * kebar * dAtom->getIyy() );
226          jy = vbar * gaussStream->getGaussian();
227 <
227 >        
228          vbar = sqrt( 2.0 * kebar * dAtom->getIzz() );
229          jz = vbar * gaussStream->getGaussian();
230          
# Line 300 | Line 264 | void Thermo::getCOMVel(double vdrift[3]){
264    }
265  
266   #ifdef IS_MPI
267 <  MPI::COMM_WORLD.Allreduce(&mtot_local,&mtot,1,MPI_DOUBLE,MPI_SUM);
268 <  MPI::COMM_WORLD.Allreduce(vdrift_local,vdrift,3,MPI_DOUBLE,MPI_SUM);
267 >  MPI_Allreduce(&mtot_local,&mtot,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
268 >  MPI_Allreduce(vdrift_local,vdrift,3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
269   #else
270    mtot = mtot_local;
271    for(vd = 0; vd < 3; vd++) {

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines