MixedHybridCompositional Class Reference

#include <mixedhybridcompositional.h>

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

List of all members.

Public Member Functions

 MixedHybridCompositional (OrthoMesh &mesh, Function3D &fPInit, Function3D &fPrescribedVelocities, Function3D &fPrescribedPression, const VecDouble &K, const VecDouble &vPor, GeneralFunctionInterface &fMobT, GeneralFunctionInterface &fRelMobilities, double grav, double dt, FlashCompositional &flash, LinearSolver &solver, int debugLevel, bool bCheckNoBC=false)
 ~MixedHybridCompositional ()
virtual void getVelocitiesAtFaces (Matrix &vel)
virtual void printOutput ()
virtual const VecDoublegetPressureAtCells ()
virtual void getNormalVelocityAtFaces (VecDouble &vFNC)
void getVelocitiesAtCells (Matrix &vel)
virtual void setDt (double dt)
const HybridFaceBCgetFaceBC ()
void testSolution ()

Protected Attributes

OrthoMeshm_mesh
GeneralFunctionInterfacem_fMobT
GeneralFunctionInterfacem_fFracFlux
const VecDoublem_vPor
const VecDoublem_vK
double m_grav
double m_dt
SparsityPattern m_sparsePattern
SparseMatrix< double > m_A
VecDouble m_B
VecDouble m_sol
VecDouble m_solAnt
VecDouble m_Vq
LinearSolverm_solver
FlashCompositionalm_flash
int _debugLevel
VecDouble m_vInvH2
double m_inv_hz
HybridFaceBC m_bc

Detailed Description

Definition at line 11 of file mixedhybridcompositional.h.


Constructor & Destructor Documentation

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

Definition at line 7 of file mixedhybridcompositional.cpp.

00008   :m_mesh(mesh),
00009    m_fMobT(fMobT),
00010    m_fFracFlux(fFracFlux),
00011    m_vPor(vPor),
00012    m_vK(K),
00013    m_grav(grav),
00014    m_dt(dt),
00015    m_solver(solver),
00016    m_flash(flash),
00017    _debugLevel(debugLevel)
00018 {
00019   
00020   
00021   this->setBoundaryCondition(m_bc,m_mesh,fPrescribedVelocities,fPrescribedPression,bCheckNoBC);
00022   this->buildSparsityPattern(mesh,m_sparsePattern);
00023   m_A.reinit(m_sparsePattern);
00024   m_B.reinit(m_mesh.numCells());
00025   m_sol.reinit(m_mesh.numCells(),NAN);
00026   m_Vq.reinit(m_mesh.numFaces());
00027   //Set the initial values for the pressure;
00028   mesh.setCentralValuesFromFunction(fInitP,m_sol);
00029 
00030 
00031   m_vInvH2.reinit(3);
00032   m_vInvH2(OrthoMesh::NORMAL_X)=1.0/(pow(mesh.get_dx(),2));
00033   m_vInvH2(OrthoMesh::NORMAL_Y)=1.0/(pow(mesh.get_dy(),2));
00034   m_vInvH2(OrthoMesh::NORMAL_Z)=1.0/(pow(mesh.get_dz(),2));
00035   m_inv_hz=1.0/mesh.get_dz();
00036 
00037   
00038 }

MixedHybridCompositional::~MixedHybridCompositional (  ) 

Definition at line 40 of file mixedhybridcompositional.cpp.

00041 {
00042 }


Member Function Documentation

const HybridFaceBC& MixedHybridCompositional::getFaceBC (  )  [inline]

Definition at line 71 of file mixedhybridcompositional.h.

00071 {return m_bc;}

void MixedHybridCompositional::getNormalVelocityAtFaces ( VecDouble vFNC  )  [virtual]

Reimplemented from DynamicBase.

Definition at line 121 of file mixedhybridcompositional.cpp.

00122 {
00123   vFNC=m_Vq;
00124 }

const VecDouble & MixedHybridCompositional::getPressureAtCells (  )  [virtual]

Reimplemented from DynamicBase.

Definition at line 127 of file mixedhybridcompositional.cpp.

00128 {
00129   return m_sol;
00130 }

void MixedHybridCompositional::getVelocitiesAtCells ( Matrix vel  ) 

Definition at line 70 of file mixedhybridcompositional.cpp.

00071 {
00072   m_mesh.projectVelocitiesInFacesToCells(vel,m_Vq);
00073 }

void MixedHybridCompositional::getVelocitiesAtFaces ( Matrix vel  )  [virtual]

Reimplemented from DynamicBase.

Definition at line 54 of file mixedhybridcompositional.cpp.

00055 {
00056   assert(vel.m() == m_mesh.numFaces() && vel.n() == 3);
00057   
00058   vel = 0.0;
00059   //Get the face iterator and the number of faces.
00060   OrthoMesh::Face_It face = m_mesh.begin_face();
00061   OrthoMesh::Face_It endf = m_mesh.end_face();
00062     
00063   //For each face
00064   for (;face!=endf;face++)
00065   {
00066     vel(face->index(),face->getNormalOrientation()) = m_Vq(face->index());
00067   }
00068 }

void MixedHybridCompositional::printOutput (  )  [virtual]

Implements DynamicBase.

Reimplemented in MixedHybridCompositionalFull, MixedHybridCompositionalSimple, and MixedHybridCompositionalTrangenstein1985.

Definition at line 77 of file mixedhybridcompositional.cpp.

00078 {
00079   if (_debugLevel > 0)
00080   {
00081     VecDouble v(m_mesh.numCells());
00082     VecDouble Vv(m_mesh.numVertices());
00083     HDF5OrthoWriter &hdf5 = HDF5OrthoWriter::getHDF5OrthoWriter(); 
00084 
00085     m_mesh.projectCentralValuesAtVertices(getPressureAtCells(),Vv);
00086     hdf5.writeScalarField(Vv,"P");
00087 
00088 
00089 
00090   
00091     if (_debugLevel > 1)
00092     {
00093       Matrix Mvel(m_mesh.numCells(),3);
00094       FEValuesExtractors::Vector velocities(0);
00095       this->getVelocitiesAtCells(Mvel);
00096       if (_debugLevel > 2)
00097       {
00098         //printf("\nVelocities at cells' barycenters\n\n");
00099         //Mvel.print_formatted(std::cout);
00100       }
00101       VecDouble v(m_mesh.numCells());
00102       for (unsigned i=0;i<v.size();i++)
00103       {
00104         v(i) = Mvel(i,0); 
00105       }
00106       VecDouble minVel(3);
00107       VecDouble maxVel(3);
00108       NumericMethods::getMinMaxValueByColumn(Mvel,minVel,maxVel);
00109       printf("Min Vel <%g,%g,%g> - Max Vel <%g,%g,%g>\n",minVel(0),minVel(1),minVel(2),maxVel(0),maxVel(1),maxVel(2));
00110   
00111       m_mesh.projectCentralValuesAtVertices(v,Vv);
00112       hdf5.writeScalarField(Vv,"Vel");
00113     }
00114   }
00115 
00116 }

virtual void MixedHybridCompositional::setDt ( double  dt  )  [inline, virtual]

Reimplemented from CompressibleDynamic.

Definition at line 70 of file mixedhybridcompositional.h.

00070 {m_dt=dt;} 

void MixedHybridCompositional::testSolution (  ) 

Member Data Documentation

Definition at line 44 of file mixedhybridcompositional.h.

Definition at line 54 of file mixedhybridcompositional.h.

double MixedHybridCompositional::m_dt [protected]

Time step

Definition at line 24 of file mixedhybridcompositional.h.

Function containing the fractional flux function

Definition at line 20 of file mixedhybridcompositional.h.

Definition at line 42 of file mixedhybridcompositional.h.

The total mobility function

Definition at line 19 of file mixedhybridcompositional.h.

Gravity acceleration m/s^2

Definition at line 23 of file mixedhybridcompositional.h.

The mesh

Definition at line 16 of file mixedhybridcompositional.h.

Store the permeability

Definition at line 22 of file mixedhybridcompositional.h.

Store the porosity

Definition at line 21 of file mixedhybridcompositional.h.


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:19 2012 for CO2INJECTION by  doxygen 1.6.3