MixedHybridCompositionalTrangenstein1985 Class Reference

#include <mixedhybridcompositionaltrangenstein1985.h>

Inheritance diagram for MixedHybridCompositionalTrangenstein1985:
Inheritance graph
[legend]
Collaboration diagram for MixedHybridCompositionalTrangenstein1985:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 MixedHybridCompositionalTrangenstein1985 (OrthoMesh &mesh, Function3D &fPInit, Function3D &fPrescribedVelocities, Function3D &fPrescribedPression, const VecDouble &K, const VecDouble &vPor, Function1D &fMobT, Function1D &fRelMobilities, double grav, double dt, FlashCompositional &flash, LinearSolver &solver, int debugLevel, bool bCheckNoBC=false)
virtual void iterate (TransportBase &trans)
virtual void printOutput ()

Protected Member Functions

void assemblySystem (SparseMatrix< double > &A, VecDouble &B, double dt, const VecDouble &vP, const VecDouble &vK, const VecDouble &vPor, FlashCompositional &flash, GeneralFunctionInterface &fMobT, GeneralFunctionInterface &fRelMobilities, double grav)
void getNormalVdtFromPressure (VecDouble &Vq, OrthoMesh &mesh, HybridFaceBC &bc, const VecDouble &K, const VecDouble &vP, const VecDouble &vMobT)
void getTotalMassFluxFromPressure (VecDouble &Vq, OrthoMesh &mesh, const HybridFaceBC &bc, const VecDouble &K, const VecDouble &vP, const Matrix &m_vC)

Detailed Description

Definition at line 8 of file mixedhybridcompositionaltrangenstein1985.h.


Constructor & Destructor Documentation

MixedHybridCompositionalTrangenstein1985::MixedHybridCompositionalTrangenstein1985 ( OrthoMesh mesh,
Function3D fPInit,
Function3D fPrescribedVelocities,
Function3D fPrescribedPression,
const VecDouble K,
const VecDouble vPor,
Function1D fMobT,
Function1D fRelMobilities,
double  grav,
double  dt,
FlashCompositional flash,
LinearSolver solver,
int  debugLevel,
bool  bCheckNoBC = false 
)

Definition at line 6 of file mixedhybridcompositionaltrangenstein1985.cpp.

00007  :MixedHybridCompositionalFull(mesh,fPInit,fPrescribedVelocities, fPrescribedPression, K, vPor,fMobT,fRelMobilities,grav,  dt,flash,solver,debugLevel, bCheckNoBC)
00008 {
00009 
00010 }


Member Function Documentation

void MixedHybridCompositionalTrangenstein1985::assemblySystem ( SparseMatrix< double > &  A,
VecDouble B,
double  dt,
const VecDouble vP,
const VecDouble vK,
const VecDouble vPor,
FlashCompositional flash,
GeneralFunctionInterface fMobT,
GeneralFunctionInterface fRelMobilities,
double  grav 
) [protected]

Reimplemented from MixedHybridCompositionalFull.

Definition at line 32 of file mixedhybridcompositionaltrangenstein1985.cpp.

00033 {
00034   assert((flash.numPhases() -1) ==  fRelMobilities.getImageDim());
00035   assert(fMobT.getDomainDim() == flash.numPhases()-1);
00036   assert(fMobT.getImageDim() == 1);
00037   assert(fRelMobilities.getDomainDim() == flash.numPhases()-1);
00038 
00039   static VecDouble totalComponents(flash.numComponents());
00040   static VecDouble dv_dm(flash.numComponents());
00041   static VecDouble phasesVol(flash.numPhases());
00042   static VecDouble phasesDensities(flash.numPhases());
00043   static VecDouble phasesViscosities(flash.numPhases());
00044   
00045   double maxStE;
00046   //Allocate the vector of saturations. Since the sum
00047   //must be one the saturation of the last phase is
00048   //a function of the others, so we just have numPhases()-1
00049   //Independent saturation variables
00050   static VecDouble vS(flash.numPhases()-1);
00051 
00052   
00053                            
00054 
00055   //ALLOCATE THE flash data
00056   static FlashData data(flash.numPhases(),flash.numComponents());
00057 
00058 
00059   //Get the cell iterator and the number of cells.
00060   OrthoMesh::Cell_It cell = m_mesh.begin_cell();
00061   Index nCells = m_mesh.numCells();
00062   //Write Zeros in the system
00063   A=0.0;
00064   B=0.0;
00065 
00066   //For each cell do
00067   double BetaG_dt;
00068   double St;
00069   for (Index cellIndex=0;cellIndex < nCells; cellIndex++,cell++)
00070   {
00071     //Compute/Get the flash data
00072     flash.flash(cellIndex,data);
00073     double P=vP(cellIndex);
00074 
00075     flash.getPhasesVolume(P,data,phasesVol);     
00076     St=phasesVol.l1_norm();
00077 
00078 
00079     if (fabs(1-St) > maxStE) maxStE = fabs(1-St); 
00080     
00081     BetaG_dt=vPor(cellIndex)*flash.getFluidCompressibility(P,data)/(dt);
00082     if (BetaG_dt < -0.0)
00083     {
00084       printf("Error BetaG_dt = %g\n",BetaG_dt);
00085       data.print();
00086       abort();
00087     }
00088     A.add(cellIndex,cellIndex,BetaG_dt);
00089     B(cellIndex)+=BetaG_dt*P; // +vPor(cellIndex)*(St-1.0)/dt;
00090 
00091 
00092                               
00093     
00094     
00095     //Get viscosities
00096     flash.getPhasesViscosities(P,data,phasesViscosities);
00097     fMobT.setParameters(phasesViscosities);
00098     fRelMobilities.setParameters(phasesViscosities);
00099 
00100 
00101 
00102     
00103     //Get vector of saturations, remember that the last saturation
00104     //is a function of the others because the sum of all saturations
00105     //must be one. So we dont consider the last saturation as a variable
00106     //for the mobilities
00107     for (unsigned k=0;k<vS.size();k++)
00108       vS(k) = phasesVol(k)/St;
00109     m_vMobT(cellIndex)=fMobT(vS);
00110     
00111     //Get total conservation components
00112     data.getTotalComponentsMoles(totalComponents);
00113 
00114     //Get Volume derivatives in terms of the total moles
00115     flash.getTotalVolumeDerivatives(P,data,dv_dm);
00116 
00117     m_vDvDm(cellIndex,0)=dv_dm(0);
00118     m_vDvDm(cellIndex,1)=dv_dm(1);
00119     assert(data.getMoles(1,0) == 0.0);
00120     m_vC(cellIndex,0)=_div(data.getMoles(0,0),phasesVol(0))*fRelMobilities(vS)*fMobT(vS);
00121     m_vC(cellIndex,1)=( _div(data.getMoles(0,1), phasesVol(0)) * fRelMobilities(vS) +
00122                         _div(data.getMoles(1,1), phasesVol(1))*(1-fRelMobilities(vS)))*fMobT(vS);
00123     //data.print();
00124     //printf("%g %g %g %g\n",phasesVol(0),phasesVol(1),fRelMobilities.value(vS),fMobT.value(vS));
00125     //printf("%g %g\n", m_vC(cellIndex,0), m_vC(cellIndex,1));
00126            
00127 
00128   }
00129 
00130 
00131   //Now with all the information gathered from the flash in each cell
00132   //build the system
00133 
00134   OrthoMesh::Face_It face = m_mesh.begin_face();
00135   unsigned nFaces = m_mesh.numFaces();
00136   
00137 
00138 
00139 
00140   //For each face
00141   for (unsigned faceIndex=0;faceIndex < nFaces; faceIndex++,face++)
00142   {
00143     unsigned c1,c2;
00144     face->getAdjCellIndices(c1,c2);
00145     
00146     //Make sure the face is not at boundary
00147     if (!face->at_boundary())
00148     {
00149       double inv_dx2 = m_vInvH2(face->getNormalOrientation());
00150 
00151 
00152       double K1 = vK(c1)*m_vC(c1,0);
00153       double K2 = vK(c2)*m_vC(c2,0);
00154       double Ceff=(K1+K2)/2.0;
00155       A.add(c1,c1, +m_vDvDm(c1,0)*Ceff * inv_dx2);
00156       A.add(c1,c2, -m_vDvDm(c1,0)*Ceff * inv_dx2);
00157       A.add(c2,c2, +m_vDvDm(c2,0)*Ceff * inv_dx2);
00158       A.add(c2,c1, -m_vDvDm(c2,0)*Ceff * inv_dx2);
00159 
00160 
00161       K1 = vK(c1)*m_vC(c1,1);
00162       K2 = vK(c2)*m_vC(c2,1);
00163       Ceff=(K1+K2)/2.0;
00164       A.add(c1,c1, +m_vDvDm(c1,1)*Ceff * inv_dx2);
00165       A.add(c1,c2, -m_vDvDm(c1,1)*Ceff * inv_dx2);
00166       A.add(c2,c2, +m_vDvDm(c2,1)*Ceff * inv_dx2);
00167       A.add(c2,c1, -m_vDvDm(c2,1)*Ceff * inv_dx2);
00168     }
00169   }
00170   //Now process the boundary conditions
00171     
00172   HybridFaceBC::iterator bcEnd = m_bc.end();
00173   for(HybridFaceBC::iterator bcIt = m_bc.begin();bcIt!=bcEnd;bcIt++)
00174   {
00175     OrthoMesh::Face_It face = bcIt->face;
00176     assert(face->at_boundary());
00177       
00178     unsigned c1,c2;
00179     face->getAdjCellIndices(c1,c2);
00180 
00181     if (bcIt->bcType == Dirichlet)
00182     {
00183       double inv_dx2 =  m_vInvH2(face->getNormalOrientation());
00184       double c = (c1 != OrthoMesh::invalidIndex()) ? c1 : c2;
00185       double Coeff = +2*vK(c)*(m_vC(c,0)*m_vDvDm(c,0) + m_vC(c,1)*m_vDvDm(c,1))*inv_dx2;
00186         
00187       A.add(c,c, Coeff);
00188       B(c) += Coeff*bcIt->value;
00189 
00190     }
00191     else //Neumann condition
00192     {
00193       if (bcIt->value != 0.0)
00194       {
00195         throw new Exception("The fully compositional model only accept null newmann condition");
00196                                  
00197       }
00198       if (c1 != OrthoMesh::invalidIndex())
00199         B(c1) += -bcIt->value;
00200       else
00201         B(c2) += +bcIt->value;
00202     }
00203   }
00204 
00205   HDF5OrthoWriter::getHDF5OrthoWriter().setVariable("St_Error",maxStE);
00206   printf("StError= %g\n",maxStE);
00207   maxStE=0;
00208   
00209 }

void MixedHybridCompositionalTrangenstein1985::getNormalVdtFromPressure ( VecDouble Vq,
OrthoMesh mesh,
HybridFaceBC bc,
const VecDouble K,
const VecDouble vP,
const VecDouble vMobT 
) [protected]

Reimplemented from MixedHybridCompositionalFull.

Definition at line 213 of file mixedhybridcompositionaltrangenstein1985.cpp.

00214 {
00215   assert(Vq.size() == mesh.numFaces());
00216   assert(vP.size() == mesh.numCells());
00217   assert(vMobT.size() == mesh.numCells());
00218   assert(K.size() == mesh.numCells());
00219   
00220   
00221   //Get the face iterator and the number of faces.
00222   OrthoMesh::Face_It face = mesh.begin_face();
00223   OrthoMesh::Face_It endf = mesh.end_face();
00224 
00225     
00226   //For each face
00227   for (;face!=endf;face++)
00228   {
00229     if (!face->at_boundary())
00230     {
00231       unsigned c1,c2;
00232       face->getAdjCellIndices(c1,c2);
00233       double K1 = K(c1)*vMobT(c1);
00234       double K2 = K(c2)*vMobT(c2);
00235       double Keff = (K1+K2)/2.0*face->areaPerCellVol();
00236       Vq(face->index()) = - Keff * (vP(c2) - vP(c1));
00237       
00238     }
00239   }
00240 
00241   //Now for boundary conditions
00242     
00243   HybridFaceBC::iterator bcEnd = bc.end();
00244   for(HybridFaceBC::iterator bcIt = bc.begin();bcIt!=bcEnd;bcIt++)
00245   {
00246     OrthoMesh::Face_It face = bcIt->face;
00247     assert(face->at_boundary());
00248       
00249     unsigned c1,c2;
00250     face->getAdjCellIndices(c1,c2);
00251 
00252     if (bcIt->bcType == Dirichlet)
00253     {
00254       if (c1 != OrthoMesh::invalidIndex())
00255       {
00256         double Keff = 2*K(c1)*vMobT(c1)*face->areaPerCellVol();
00257         Vq(face->index()) = - Keff*(bcIt->value - vP(c1)); 
00258         //      printf("Pressions R: %d %g %g = %g\n",c1,bcIt->value, m_sol(c1),Vq(bcIt->faceIndex));
00259       }
00260       else
00261       {
00262         double Keff = 2*K(c2)*vMobT(c2)*face->areaPerCellVol();
00263         Vq(face->index()) = Keff*(bcIt->value - vP(c2));
00264         //printf("Pressions L: %d %g %g = %g\n",c2,bcIt->value, m_sol(c2),Vq(bcIt->faceIndex));
00265       }
00266     }
00267     else
00268     {
00269       if (c1 != OrthoMesh::invalidIndex())
00270         Vq(face->index()) = bcIt->value/face->area();
00271       else
00272         Vq(face->index()) = bcIt->value/face->area();
00273     }
00274   }
00275 
00276 }

void MixedHybridCompositionalTrangenstein1985::getTotalMassFluxFromPressure ( VecDouble Vq,
OrthoMesh mesh,
const HybridFaceBC bc,
const VecDouble K,
const VecDouble vP,
const Matrix m_vC 
) [protected]

Reimplemented from MixedHybridCompositionalFull.

Definition at line 282 of file mixedhybridcompositionaltrangenstein1985.cpp.

00283 {
00284   assert(Vq.size() == mesh.numFaces());
00285   assert(vP.size() == mesh.numCells());
00286   assert(K.size() == mesh.numCells());
00287   
00288   
00289   //Get the face iterator and the number of faces.
00290   OrthoMesh::Face_It face = mesh.begin_face();
00291   OrthoMesh::Face_It endf = mesh.end_face();
00292 
00293     
00294   //For each face
00295   for (;face!=endf;face++)
00296   {
00297     if (!face->at_boundary())
00298     {
00299       unsigned c1,c2;
00300       face->getAdjCellIndices(c1,c2);
00301       double K1 = K(c1)*m_vC(c1,0) + K(c1)*m_vC(c1,1);
00302       double K2 = K(c2)*m_vC(c2,0) + K(c2)*m_vC(c2,1);
00303       double Keff = (K1+K2)/2.0*face->areaPerCellVol();
00304       Vq(face->index()) = - Keff * (vP(c2) - vP(c1));
00305       
00306     }
00307   }
00308 
00309   //Now for boundary conditions
00310     
00311   HybridFaceBC::const_iterator bcEnd = bc.end();
00312   for(HybridFaceBC::const_iterator bcIt = bc.begin();bcIt!=bcEnd;bcIt++)
00313   {
00314     OrthoMesh::Face_It face = bcIt->face;
00315     assert(face->at_boundary());
00316       
00317     unsigned c1,c2;
00318     face->getAdjCellIndices(c1,c2);
00319 
00320     if (bcIt->bcType == Dirichlet)
00321     {
00322       if (c1 != OrthoMesh::invalidIndex())
00323       {
00324         double Keff = 2*K(c1)*(m_vC(c1,0) + m_vC(c1,1))*face->areaPerCellVol();
00325         Vq(face->index()) = - Keff*(bcIt->value - vP(c1)); 
00326         //      printf("Pressions R: %d %g %g = %g\n",c1,bcIt->value, m_sol(c1),Vq(bcIt->faceIndex));
00327       }
00328       else
00329       {
00330         double Keff = 2*K(c2)*(m_vC(c2,0) + m_vC(c2,1))*face->areaPerCellVol();
00331         Vq(face->index()) = Keff*(bcIt->value - vP(c2));
00332         //printf("Pressions L: %d %g %g = %g\n",c2,bcIt->value, m_sol(c2),Vq(bcIt->faceIndex));
00333       }
00334     }
00335     else
00336     {
00337       if (c1 != OrthoMesh::invalidIndex())
00338         Vq(face->index()) = bcIt->value/face->area();
00339       else
00340         Vq(face->index()) = bcIt->value/face->area();
00341     }
00342   }
00343 
00344   
00345 }

void MixedHybridCompositionalTrangenstein1985::iterate ( TransportBase trans  )  [virtual]

Reimplemented from MixedHybridCompositionalFull.

Definition at line 12 of file mixedhybridcompositionaltrangenstein1985.cpp.

00013 {
00014   this->assemblySystem(m_A,m_B,m_dt,m_sol,m_vK,m_vPor,m_flash,m_fMobT,m_fFracFlux,m_grav);
00015 
00016   VecDouble PAnt=m_sol;
00017 
00018   m_solver.solve(m_A,m_sol,m_B);
00019 
00020   Point3D p(5,50,50);
00021   OrthoMesh::Cell_It cell = m_mesh.getCellAt(p);
00022   //testSolution(cell,m_dt, v1,m_vPor,m_vK,m_sol,PAnt,m_vMobT,v2); 
00023 
00024   //Get the fluxes
00025   this->getNormalVdtFromPressure(m_Vq,m_mesh,m_bc,m_vK,m_sol,m_vMobT);
00026 
00027 }

void MixedHybridCompositionalTrangenstein1985::printOutput (  )  [virtual]

Reimplemented from MixedHybridCompositionalFull.

Definition at line 351 of file mixedhybridcompositionaltrangenstein1985.cpp.

00352 {
00353   MixedHybridCompositional::printOutput();
00354 
00355   if (_debugLevel > 1)
00356   {
00357     Matrix VelM(m_mesh.numVertices(),3);
00358     VecDouble Vf(m_mesh.numFaces());
00359     VecDouble Vv(m_mesh.numVertices());
00360     VecDouble Vc(m_mesh.numCells());
00361     HDF5OrthoWriter &hdf5 = HDF5OrthoWriter::getHDF5OrthoWriter(); 
00362     MixedHybridCompositionalTrangenstein1985::getTotalMassFluxFromPressure(Vf,m_mesh,getFaceBC(),m_vK,m_sol,m_vC);
00363     m_mesh.projectVelocitiesInFacesToVertices(VelM,Vf);
00364     NumericMethods::matrixGetColumn(Vv,VelM,0);
00365     hdf5.writeScalarField(Vv,"VelM");
00366 
00367   }
00368   
00369 
00370 
00371   
00372 }


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
Generated on Sun Apr 8 23:13:22 2012 for CO2INJECTION by  doxygen 1.6.3