BiphasicDiff Class Reference

#include <biphasicdiff.h>

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

List of all members.

Public Member Functions

 BiphasicDiff (OrthoMesh &mesh, const VecDouble &vK, const VecDouble &vPor, Function1D &fMobPc, Function1D &DfPc, LinearSolver &solver)
virtual ~BiphasicDiff ()
virtual void iterate (double dt)
virtual void printOutput ()

Private Attributes

OrthoMeshm_mesh
SparsityPattern m_sparsePattern
SparseMatrix< double > m_A
VecDouble m_B
VecDouble m_vCoeff
LinearSolverm_solver
const VecDoublem_vK
const VecDoublem_vPor
Function1Dm_fMobPc
Function1Dm_DfPc
HybridFaceBC m_bc

Detailed Description

A Diffusive Step used for the classic biphasic uncompressible flow Note, the problem solves the parabolic problem with mixed hybrid method. This is the flow created by the gradient of pc, we assume tha no flux crossed the boundary of the domain due to the cappillary pressure gradient.

So the mixed hybrid method has a null flux boundary condition

Definition at line 17 of file biphasicdiff.h.


Constructor & Destructor Documentation

BiphasicDiff::BiphasicDiff ( OrthoMesh mesh,
const VecDouble vK,
const VecDouble vPor,
Function1D fMobPc,
Function1D DfPc,
LinearSolver solver 
)

Definition at line 5 of file biphasicdiff.cpp.

00006    :m_mesh(mesh),m_solver(solver),m_vK(vK),m_vPor(vPor),m_fMobPc(fMobPc),m_DfPc(DfPc)
00007  {
00008    m_B.reinit(mesh.numCells());
00009    buildSparsityPattern(mesh,m_sparsePattern);
00010    m_A.reinit(m_sparsePattern);
00011    m_vCoeff.reinit(mesh.numCells());
00012  }

virtual BiphasicDiff::~BiphasicDiff (  )  [inline, virtual]

Definition at line 35 of file biphasicdiff.h.

00035 {}


Member Function Documentation

void BiphasicDiff::iterate ( double  dt  )  [virtual]

Implements DiffusiveStep.

Definition at line 15 of file biphasicdiff.cpp.

00016 {
00017   ArrayOfVecDouble &transSol = getTransport().getSolutionAtCells();
00018   assert(transSol.vecs_size() == 1);
00019   VecDoubleRef vS(&(transSol(0,0)),transSol.size());
00020 
00021   
00022 
00023   Index nCells = m_mesh.numCells();
00024   double cellVol=m_mesh.begin_cell()->volume();  
00025 
00026   //Write Zeros in the system
00027   m_A=0.0;
00028   m_B=0.0;
00029 
00030   for (Index cellIndex=0;cellIndex < nCells; cellIndex++)
00031   {
00032 
00033 
00034     double Saq=vS(cellIndex);;
00035     double por=m_vPor(cellIndex);
00036     
00037     m_A.add(cellIndex,cellIndex,por*cellVol/dt);
00038     m_B(cellIndex)+=por*Saq/dt*cellVol;
00039 
00040 
00041     m_vCoeff(cellIndex)=-m_vK(cellIndex)*m_DfPc(Saq)*m_fMobPc(Saq);
00042     printf("Flux= %g %g %g\n",m_vK(cellIndex),m_DfPc(Saq),m_fMobPc(Saq));
00043 
00044   }
00045 
00046   //Now solve the system
00047   MixedHybridBase::assemblySystem(m_A,m_B,m_mesh,m_bc,m_vCoeff);
00048 
00049 
00050   m_solver.solve(m_A,vS,m_B);
00051 
00052 
00053   
00054 }

void BiphasicDiff::printOutput (  )  [virtual]

Implements DiffusiveStep.

Definition at line 57 of file biphasicdiff.cpp.

00058 {
00059 
00060   ArrayOfVecDouble &transSol = getTransport().getSolutionAtCells();
00061   VecDoubleRef vS(&(transSol(0,0)),transSol.size());
00062   VecDouble Vq(m_mesh.numFaces());
00063   VecDouble Vv(m_mesh.numVertices());
00064   MixedHybridBase::getNormalVelocitiesFromPressure(Vq,m_mesh,m_bc,m_vCoeff,vS);
00065   HDF5OrthoWriter &hdf5 = HDF5OrthoWriter::getHDF5OrthoWriter();
00066   Matrix M(m_mesh.numVertices(),3);
00067   m_mesh.projectVelocitiesInFacesToVertices(M,Vq);
00068   M.copyColumn(Vv,2);
00069   hdf5.writeScalarField(Vv,"DiffVelZ");
00070 }


Member Data Documentation

SparseMatrix<double> BiphasicDiff::m_A [private]

Definition at line 22 of file biphasicdiff.h.

Right hand side of the system

Definition at line 23 of file biphasicdiff.h.

A vector of boundary conditions since the boundaries use neumman condition for the flux this vector is empty

Definition at line 28 of file biphasicdiff.h.

The mobility and derivative of Pc

Definition at line 27 of file biphasicdiff.h.

Definition at line 27 of file biphasicdiff.h.

Definition at line 20 of file biphasicdiff.h.

Definition at line 25 of file biphasicdiff.h.

SparsityPattern BiphasicDiff::m_sparsePattern [private]

Definition at line 21 of file biphasicdiff.h.

Definition at line 24 of file biphasicdiff.h.

const VecDouble& BiphasicDiff::m_vK [private]

Definition at line 26 of file biphasicdiff.h.

const VecDouble & BiphasicDiff::m_vPor [private]

Store porosity and permeability

Definition at line 26 of file biphasicdiff.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:12:57 2012 for CO2INJECTION by  doxygen 1.6.3