CompDiffusiveStep Class Reference

#include <compdiffusivestep.h>

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

List of all members.

Classes

struct  FlashBC

Public Member Functions

 CompDiffusiveStep (OrthoMesh &mesh, Function3D &fPrescribedPressure, Function3D &fTransportBC, const VecDouble &vK, const VecDouble &vPor, Function1D &fMobPc, Function1D &DfPc, FlashCompositional &flash, LinearSolver &solver)
virtual ~CompDiffusiveStep ()
void setSaqBC (VecIncAccess< FlashBC > &flashbc, ScalarBC &swbc, const Function3D &fPressureBC, const Function3D &fTransBC, FlashCompositional &flash)
virtual void iterate (double dt)
virtual void printOutput ()

Protected Attributes

VecDouble m_vSaq
VecDouble m_cV
VecDouble m_vV
VecDouble m_vSaqPrev
OrthoMeshm_mesh
ScalarBC m_swBC
SparsityPattern m_sparsePattern
SparseMatrix< double > m_A
VecDouble m_sol
VecDouble m_B
VecDouble m_vCoeff
VecDouble m_vCoeff2
LinearSolverm_solver
FlashCompositionalm_flash
ScalarBC m_bc
const VecDoublem_vK
const VecDoublem_vPor
Function1Dm_fMobPc
Function1Dm_DfPc

Private Attributes

VecIncAccess< FlashBCm_flashBC

Detailed Description

Definition at line 12 of file compdiffusivestep.h.


Constructor & Destructor Documentation

CompDiffusiveStep::CompDiffusiveStep ( OrthoMesh mesh,
Function3D fPrescribedPressure,
Function3D fTransportBC,
const VecDouble vK,
const VecDouble vPor,
Function1D fMobPc,
Function1D DfPc,
FlashCompositional flash,
LinearSolver solver 
)

Definition at line 7 of file compdiffusivestep.cpp.

00008   :m_mesh(mesh),m_solver(solver),m_flash(flash),m_vK(K),m_vPor(vPor),m_fMobPc(fMobPc),m_DfPc(DfPc)
00009  {
00010    m_vCoeff.reinit(mesh.numCells());
00011    m_vCoeff2.reinit(mesh.numCells());
00012    m_sol.reinit(mesh.numCells());
00013    m_B.reinit(mesh.numCells());
00014 
00015    m_vSaq.reinit(mesh.numCells());
00016    m_vSaqPrev.reinit(mesh.numCells());
00017    
00018    //setSaqBC(m_flashBC,m_swBC,fPrescribedPressure,fTransportBC,flash);
00019    
00020    buildSparsityPattern(mesh,m_sparsePattern);
00021    m_A.reinit(m_sparsePattern);
00022    m_cV.reinit(mesh.numCells());
00023    m_vV.reinit(mesh.numVertices());
00024 
00025  }

CompDiffusiveStep::~CompDiffusiveStep (  )  [virtual]

Definition at line 186 of file compdiffusivestep.cpp.

00187 {
00188 
00189 }


Member Function Documentation

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

Implements DiffusiveStep.

Definition at line 70 of file compdiffusivestep.cpp.

00071 {
00072   const VecDouble &vP=m_flash.getDynamicModule().getPressureAtCells();
00073   TransportBase &trans = m_flash.getTransportModule();
00074   assert(m_flash.numPhases()==2);
00075   
00076   static VecDouble phasesVol(m_flash.numPhases());
00077 
00078 
00079 
00080 
00081   //Allocate the flash data
00082   static FlashData data(m_flash.numPhases(),m_flash.numComponents());
00083 
00084 
00085   Index nCells = m_mesh.numCells();
00086   double cellVol=m_mesh.begin_cell()->volume();  
00087 
00088   //Write Zeros in the system
00089   m_A=0.0;
00090   m_B=0.0;
00091 
00092   for (Index cellIndex=0;cellIndex < nCells; cellIndex++)
00093   {
00094     //Compute/Get the flash data
00095     m_flash.flash(cellIndex,data);
00096     double P=vP(cellIndex);
00097     m_flash.getPhasesVolume(P,data,phasesVol);     
00098 
00099 
00100     double Saq=phasesVol(0)/(phasesVol(1)+phasesVol(0));
00101     m_vSaqPrev(cellIndex)=Saq;
00102     double por=m_vPor(cellIndex);
00103     //double molarVol = Saq/data.getPhaseTotalMoles(0);    
00104     double molarVol;
00105     if (fabs(Saq) == 0.0)
00106       molarVol=1.2733e-4;
00107     else
00108       molarVol = Saq/data.getPhaseTotalMoles(0);
00109     
00110     m_A.add(cellIndex,cellIndex,molarVol*por*cellVol/dt);
00111     m_B(cellIndex)+=molarVol*por*Saq/dt*cellVol;
00112 
00113     m_vCoeff(cellIndex)=-m_vK(cellIndex)*m_DfPc(Saq)*m_fMobPc(Saq);
00114     m_vCoeff2(cellIndex)=m_vCoeff(cellIndex)*molarVol;
00115   }
00116 
00117   //Now solve the system
00118   MixedHybridBase::assemblySystem(m_A,m_B,m_mesh,m_swBC,m_vCoeff2);
00119 
00120 
00121   m_solver.solve(m_A,m_vSaq,m_B);
00122   unsigned n_components = m_flash.numComponents();
00123 
00124   ArrayOfVecDouble &vSol = trans.getSolutionAtCells();
00125 
00126   //m_sol contains the saturation
00127 
00128 
00129   for (unsigned cmp=0;cmp<n_components;cmp++)
00130   {
00131     OrthoMesh::Face_It face = m_mesh.begin_face();
00132     OrthoMesh::Face_It end = m_mesh.end_face();
00133     Index c1,c2;
00134     FlashData data1(m_flash.numPhases(),m_flash.numComponents());
00135     FlashData data2(m_flash.numPhases(),m_flash.numComponents());
00136     VecIncAccess<FlashBC>::iterator it=m_flashBC.begin();
00137     for (;face!=end;face++)
00138     {
00139       if (!face->at_boundary())
00140       {
00141         face->getAdjCellIndices(c1,c2);
00142         double Saq1=m_vSaqPrev(c1);
00143         double Saq2=m_vSaqPrev(c2);
00144         m_flash.flash(c1,data1);
00145         m_flash.flash(c2,data2);
00146         double Kflux1 = m_vCoeff(c1)*(_div(data1.getMoles(0,cmp),Saq1) -  _div(data1.getMoles(1,cmp),(1-Saq1)));
00147         double Kflux2 = m_vCoeff(c2)*( _div(data2.getMoles(0,cmp),Saq2) - _div(data2.getMoles(1,cmp),(1-Saq2)));
00148         double flux =  - 2*_div(Kflux1*Kflux2, (Kflux1+Kflux2))*(m_vSaq(c2)-m_vSaq(c1))*face->areaPerCellVol();
00149         //double flux =  -(Kflux1+Kflux2)/2.0*(m_vSaq(c2)-m_vSaq(c1))*face->areaPerCellVol();
00150 
00151         //if (cmp == 0)
00152         //  printf("%d) %g %g %g\n",face->index(),flux,Kflux1,Kflux2);
00153         vSol(c1,cmp)-=flux/m_vPor(c1)*dt*face->areaPerCellVol();
00154         vSol(c2,cmp)+=flux/m_vPor(c2)*dt*face->areaPerCellVol();
00155       }
00156       else
00157       {
00158         Index faceIndex = face->index();
00159         const FlashBC *flash_bc = it.inc_search(faceIndex);
00160         if (flash_bc)
00161         {
00162           abort();
00163           face->getAdjCellIndices(c1,c2);
00164           double Saq = flash_bc->Saq;
00165           if (c1 != OrthoMesh::INVALID_INDEX)
00166           {  
00167             m_flash.flash(c1,data1);
00168             double Kflux = m_vCoeff(c1)*(_div(data1.getMoles(0,cmp),m_vSaqPrev(c1)) - _div(data1.getMoles(1,cmp),1-m_vSaqPrev(c1)));
00169             double flux= -2*Kflux*(Saq-m_vSaq(c1))*face->areaPerCellVol();
00170             vSol(c1,cmp)-=flux/m_vPor(c1)*dt*face->areaPerCellVol();
00171           }
00172           else
00173           {
00174             m_flash.flash(c2,data2);
00175             double Kflux = m_vCoeff(c2)*(_div(data2.getMoles(0,cmp),m_vSaqPrev(c2)) - _div(data2.getMoles(1,cmp),1-m_vSaqPrev(c2)));
00176             double flux= -2*Kflux*(m_vSaq(c2) - Saq)*face->areaPerCellVol();
00177             vSol(c2,cmp)+=flux/m_vPor(c2)*dt*face->areaPerCellVol();
00178           }
00179         }
00180       }
00181     }
00182   }
00183 }

virtual void CompDiffusiveStep::printOutput (  )  [inline, virtual]

Implements DiffusiveStep.

Definition at line 47 of file compdiffusivestep.h.

00047 {}

void CompDiffusiveStep::setSaqBC ( VecIncAccess< FlashBC > &  flashbc,
ScalarBC swbc,
const Function3D fPressureBC,
const Function3D fTransBC,
FlashCompositional flash 
)

Definition at line 27 of file compdiffusivestep.cpp.

00028 {
00029   assert(flash.numPhases() == 2);
00030   assert(fTransBC.getImageDim() == flash.numComponents());
00031   unsigned ncmps = flash.numComponents();
00032   
00033   VecDouble cmps(flash.numComponents());
00034   VecDouble phasesVol(flash.numPhases());
00035   FlashData flashData(flash.numPhases(),flash.numComponents());
00036   flashData.allocateOwnMemory();
00037 
00038   OrthoMesh::Face_It face = m_mesh.begin_face();
00039   OrthoMesh::Face_It endf = m_mesh.end_face();
00040   for (;face!=endf;face++)
00041   {
00042     if (face->at_boundary())
00043     {
00044       Point3D p = face->barycenter();
00045       if (fTransBC.isInDomain(p))
00046       {
00047         assert(fPressureBC.isInDomain(p));
00048         for (unsigned i=0;i<ncmps;i++)
00049         {
00050           cmps(i)=fTransBC(p,i);
00051         }
00052         double P = fPressureBC(p);
00053         flash.flash(P,cmps,flashData);
00054         flash.getPhasesVolume(P,flashData,phasesVol);
00055         double Saq = phasesVol(0)/(phasesVol(1)+phasesVol(0));
00056         FlashBC  _flashBC; //(flash.numPhases(),flash.numComponents());;
00057         _ScalarBC _sBC;
00058         _flashBC.Saq=Saq;
00059         _sBC.value=Saq;
00060         _sBC.face=face;
00061         swbc.push_back(_sBC);
00062         flashbc.push_back(face->index(),_flashBC);
00063       }
00064     }
00065   }
00066 }


Member Data Documentation

SparseMatrix<double> CompDiffusiveStep::m_A [protected]

Definition at line 32 of file compdiffusivestep.h.

Right hand side of the system

Definition at line 34 of file compdiffusivestep.h.

Definition at line 38 of file compdiffusivestep.h.

Definition at line 27 of file compdiffusivestep.h.

Definition at line 40 of file compdiffusivestep.h.

Definition at line 37 of file compdiffusivestep.h.

Definition at line 23 of file compdiffusivestep.h.

Definition at line 40 of file compdiffusivestep.h.

Definition at line 29 of file compdiffusivestep.h.

Definition at line 33 of file compdiffusivestep.h.

Definition at line 36 of file compdiffusivestep.h.

SparsityPattern CompDiffusiveStep::m_sparsePattern [protected]

Definition at line 31 of file compdiffusivestep.h.

Definition at line 30 of file compdiffusivestep.h.

Definition at line 35 of file compdiffusivestep.h.

Definition at line 35 of file compdiffusivestep.h.

const VecDouble& CompDiffusiveStep::m_vK [protected]

Definition at line 39 of file compdiffusivestep.h.

const VecDouble & CompDiffusiveStep::m_vPor [protected]

Definition at line 39 of file compdiffusivestep.h.

Definition at line 27 of file compdiffusivestep.h.

Definition at line 28 of file compdiffusivestep.h.

Definition at line 27 of file compdiffusivestep.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:59 2012 for CO2INJECTION by  doxygen 1.6.3