24#ifndef DUMUX_FVPRESSURE2P2C_MULTIPHYSICS_HH
25#define DUMUX_FVPRESSURE2P2C_MULTIPHYSICS_HH
27#include <dune/common/float_cmp.hh>
68template<
class TypeTag>
78 using MaterialLaw =
typename SpatialParams::MaterialLaw;
89 dim = GridView::dimension, dimWorld = GridView::dimensionworld
93 pw = Indices::pressureW
97 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
98 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
99 contiWEqIdx = Indices::contiWEqIdx, contiNEqIdx = Indices::contiNEqIdx
103 using Element =
typename GridView::Traits::template Codim<0>::Entity;
104 using Grid =
typename GridView::Grid;
105 using Intersection =
typename GridView::Intersection;
108 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
109 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
110 using PhaseVector = Dune::FieldVector<Scalar, getPropValue<TypeTag, Properties::NumPhases>()>;
114 using EntryType = Dune::FieldVector<Scalar, 2>;
124 void get1pSource(EntryType& sourceEntry,
const Element& elementI,
const CellData& cellDataI);
126 void get1pStorage(EntryType& storageEntry,
const Element& elementI, CellData& cellDataI);
128 void get1pFlux(EntryType& entries,
const Intersection& intersection,
129 const CellData& cellDataI);
132 const Intersection& intersection,
133 const CellData& cellDataI);
139 int size = this->
problem().gridView().size(0);
140 for (
int i = 0; i < size; i++)
142 CellData& cellData = this->
problem().variables().cellData(i);
143 cellData.subdomain() = 2;
160 int eIdxGlobal =
problem().variables().index(element);
161 CellData& cellData =
problem().variables().cellData(eIdxGlobal);
162 outstream <<
" "<< cellData.subdomain();
177 int eIdxGlobal =
problem().variables().index(element);
178 CellData& cellData =
problem().variables().cellData(eIdxGlobal);
180 instream >> subdomainIdx;
181 cellData.setSubdomainAndFluidStateType(subdomainIdx);
194 template<
class MultiWriter>
199 if(
problem().vtkOutputLevel()>=1)
201 int size =
problem().gridView().size(0);
203 Dune::BlockVector<Dune::FieldVector<int,1> >* subdomainPtr = writer.template allocateManagedBuffer<int, 1> (size);
204 for (
int i = 0; i < size; i++)
206 CellData& cellData =
problem().variables().cellData(i);
207 (*subdomainPtr)[i] = cellData.subdomain();
209 writer.attachCellData(*subdomainPtr,
"subdomain");
226 using DataHandle = VectorExchange<ElementMapper, Dune::BlockVector<Dune::FieldVector<int, 1> > >;
263template<
class TypeTag>
275 for (
const auto& element : elements(problem().gridView()))
278 int eIdxGlobalI = problem().variables().index(element);
281 if (element.partitionType() == Dune::InteriorEntity)
284 CellData& cellDataI = problem().variables().cellData(eIdxGlobalI);
286 Dune::FieldVector<Scalar, 2> entries(0.);
289 if(cellDataI.subdomain() != 2)
290 problem().pressureModel().get1pSource(entries,element, cellDataI);
292 problem().pressureModel().getSource(entries,element, cellDataI, first);
294 this->
f_[eIdxGlobalI] = entries[
rhs];
298 for (
const auto& intersection : intersections(problem().gridView(), element))
301 if (intersection.neighbor())
303 int eIdxGlobalJ = problem().variables().index(intersection.outside());
305 if (cellDataI.subdomain() != 2
306 or problem().variables().cellData(eIdxGlobalJ).subdomain() != 2)
307 get1pFlux(entries, intersection, cellDataI);
309 problem().pressureModel().getFlux(entries, intersection, cellDataI, first);
312 this->
f_[eIdxGlobalI] -= entries[
rhs];
314 this->
A_[eIdxGlobalI][eIdxGlobalI] += entries[
matrix];
316 this->
A_[eIdxGlobalI][eIdxGlobalJ] = -entries[
matrix];
323 if (cellDataI.subdomain() != 2)
324 problem().pressureModel().get1pFluxOnBoundary(entries, intersection, cellDataI);
326 problem().pressureModel().getFluxOnBoundary(entries, intersection, cellDataI, first);
329 this->
f_[eIdxGlobalI] += entries[
rhs];
331 this->
A_[eIdxGlobalI][eIdxGlobalI] += entries[
matrix];
337 if (cellDataI.subdomain() != 2)
338 problem().pressureModel().get1pStorage(entries, element, cellDataI);
340 problem().pressureModel().getStorage(entries, element, cellDataI, first);
342 this->
f_[eIdxGlobalI] += entries[
rhs];
344 this->
A_[eIdxGlobalI][eIdxGlobalI] += entries[
matrix];
349 this->
A_[eIdxGlobalI] = 0.0;
350 this->
A_[eIdxGlobalI][eIdxGlobalI] = 1.0;
351 this->
f_[eIdxGlobalI] = this->
pressure()[eIdxGlobalI];
370template<
class TypeTag>
372 const Element& elementI,
const CellData& cellDataI)
377 Scalar volume = elementI.geometry().volume();
378 int subdomainIdx = cellDataI.subdomain();
381 PrimaryVariables source(NAN);
382 problem().source(source, elementI);
383 source[1+subdomainIdx] /= cellDataI.density(subdomainIdx);
385 sourceEntry[1] = volume * source[1+subdomainIdx];
403template<
class TypeTag>
405 const Element& elementI,
410 int eIdxGlobalI = problem().variables().index(elementI);
411 int presentPhaseIdx = cellDataI.subdomain();
412 Scalar volume = elementI.geometry().volume();
415 Scalar timestep_ = problem().timeManager().timeStepSize();
420 Scalar& incp = this->
incp_;
423 PhaseVector p_(incp);
424 p_[nPhaseIdx] += cellDataI.pressure(nPhaseIdx);
425 p_[wPhaseIdx] += cellDataI.pressure(wPhaseIdx);
427 Scalar sumC = (cellDataI.massConcentration(wCompIdx) + cellDataI.massConcentration(nCompIdx));
428 Scalar Z0 = cellDataI.massConcentration(wCompIdx) / sumC;
435 Scalar v_ = 1. / pseudoFluidState.
density(presentPhaseIdx);
436 cellDataI.dv_dp() = (sumC * ( v_ - (1. /cellDataI.density(presentPhaseIdx)))) /incp;
438 if (cellDataI.dv_dp()>0)
441 Dune::dinfo <<
"dv_dp larger 0 at Idx " << eIdxGlobalI <<
" , try and invert secant"<< std::endl;
446 v_ = 1. / pseudoFluidState.
density(presentPhaseIdx);
447 cellDataI.dv_dp() = (sumC * ( v_ - (1. /cellDataI.density(presentPhaseIdx)))) /incp;
449 if (cellDataI.dv_dp()>0)
451 Dune::dinfo <<__FILE__<<
"dv_dp still larger 0 after inverting secant. regularize"<< std::endl;
452 cellDataI.dv_dp() *= -1;
457 Scalar compress_term = cellDataI.dv_dp() / timestep_;
459 storageEntry[
matrix] -= compress_term*volume;
460 storageEntry[
rhs] -= cellDataI.pressure(
pressureType) * compress_term * volume;
464 if (isnan(compress_term) ||isinf(compress_term))
465 DUNE_THROW(Dune::MathError,
"Compressibility term leads to NAN matrix entry at index " << eIdxGlobalI);
468 DUNE_THROW(Dune::NotImplemented,
"Compressibility is switched off???");
473 if( problem().timeManager().episodeWillBeFinished()
474 || problem().timeManager().willBeFinished())
476 problem().variables().cellData(eIdxGlobalI).errorCorrection() = 0.;
482 problem().variables().cellData(eIdxGlobalI).volumeError() /= timestep_;
484 Scalar erri = fabs(cellDataI.volumeError());
493 if ((erri*timestep_ > 5e-5) && (erri > x_lo * maxError) && (!problem().timeManager().willBeFinished()))
495 if (erri <= x_mi * maxError)
497 problem().variables().cellData(eIdxGlobalI).errorCorrection() =
498 fac* (1-x_mi*(lofac-1)/(x_lo-x_mi) + (lofac-1)/(x_lo-x_mi)*erri/maxError)
499 * cellDataI.volumeError() * volume;
502 problem().variables().cellData(eIdxGlobalI).errorCorrection() =
503 fac * (1 + x_mi - hifac*x_mi/(1-x_mi) + (hifac/(1-x_mi)-1)*erri/maxError)
504 * cellDataI.volumeError() * volume;
507 problem().variables().cellData(eIdxGlobalI).errorCorrection()=0 ;
526template<
class TypeTag>
528 const Intersection& intersection,
const CellData& cellDataI)
531 auto elementI = intersection.inside();
534 const GlobalPosition& globalPos = elementI.geometry().center();
537 DimMatrix permeabilityI(problem().spatialParams().intrinsicPermeability(elementI));
540 const GlobalPosition& unitOuterNormal = intersection.centerUnitOuterNormal();
543 Scalar faceArea = intersection.geometry().volume();
546 auto neighbor = intersection.outside();
547 int eIdxGlobalJ = problem().variables().index(neighbor);
548 CellData& cellDataJ = problem().variables().cellData(eIdxGlobalJ);
551 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
554 GlobalPosition distVec = globalPosNeighbor - globalPos;
557 Scalar dist = distVec.two_norm();
559 GlobalPosition unitDistVec(distVec);
562 DimMatrix permeabilityJ
563 = problem().spatialParams().intrinsicPermeability(neighbor);
566 DimMatrix meanPermeability(0);
569 Dune::FieldVector<Scalar, dim> permeability(0);
570 meanPermeability.mv(unitDistVec, permeability);
575 int phaseIdx = min(cellDataI.subdomain(), cellDataJ.subdomain());
577 Scalar rhoMean = 0.5 * (cellDataI.density(phaseIdx) + cellDataJ.density(phaseIdx));
580 Scalar potential = (cellDataI.pressure(phaseIdx) - cellDataJ.pressure(phaseIdx)) / dist;
582 potential += rhoMean * (unitDistVec *
gravity_);
588 lambda = cellDataI.mobility(phaseIdx);
589 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiWEqIdx,
false);
590 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiNEqIdx,
false);
592 else if (potential < 0.)
594 lambda = cellDataJ.mobility(phaseIdx);
595 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiWEqIdx,
true);
596 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiNEqIdx,
true);
600 lambda =
harmonicMean(cellDataI.mobility(phaseIdx) , cellDataJ.mobility(phaseIdx));
601 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiWEqIdx,
false);
602 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiNEqIdx,
false);
605 entries[0] = lambda * faceArea * fabs(permeability * unitOuterNormal) / (dist);
606 entries[1] = rhoMean * lambda;
607 entries[1] *= faceArea * fabs(permeability * unitOuterNormal) * (unitDistVec *
gravity_);
629template<
class TypeTag>
631 const Intersection& intersection,
const CellData& cellDataI)
635 auto elementI = intersection.inside();
636 const GlobalPosition& globalPos = elementI.geometry().center();
638 int phaseIdx = cellDataI.subdomain();
641 const GlobalPosition& unitOuterNormal = intersection.centerUnitOuterNormal();
643 Scalar faceArea = intersection.geometry().volume();
646 const GlobalPosition& globalPosFace = intersection.geometry().center();
649 GlobalPosition distVec(globalPosFace - globalPos);
650 Scalar dist = distVec.two_norm();
651 GlobalPosition unitDistVec(distVec);
655 BoundaryTypes bcType;
656 problem().boundaryTypes(bcType, intersection);
659 PhaseVector pressBC(0.);
662 if (bcType.isDirichlet(Indices::pressureEqIdx))
665 DimMatrix permeabilityI(problem().spatialParams().intrinsicPermeability(elementI));
668 int axis = intersection.indexInInside() / 2;
673 Scalar lambdaI = cellDataI.mobility(phaseIdx);
676 Dune::FieldVector<Scalar, dim> permeability(0);
677 permeabilityI.mv(unitDistVec, permeability);
680 FluidState BCfluidState;
683 PrimaryVariables primaryVariablesOnBoundary(NAN);
684 problem().dirichlet(primaryVariablesOnBoundary, intersection);
688 problem().transportModel().evalBoundary(globalPosFace,
694 Scalar densityBound =
695 FluidSystem::density(BCfluidState, phaseIdx);
696 Scalar viscosityBound =
697 FluidSystem::viscosity(BCfluidState, phaseIdx);
700 Scalar lambdaBound = 0.;
703 case Indices::satDependent:
705 lambdaBound = BCfluidState.saturation(phaseIdx)
709 case Indices::permDependent:
711 if (phaseIdx == wPhaseIdx)
712 lambdaBound = MaterialLaw::krw(
713 problem().spatialParams().materialLawParams(elementI), BCfluidState.saturation(wPhaseIdx))
716 lambdaBound = MaterialLaw::krn(
717 problem().spatialParams().materialLawParams(elementI), BCfluidState.saturation(wPhaseIdx))
722 Scalar rhoMean = 0.5 * (cellDataI.density(phaseIdx) + densityBound);
724 Scalar potential = 0;
727 potential = (cellDataI.pressure(phaseIdx) - pressBC[phaseIdx]) / dist;
729 potential += rhoMean * (unitDistVec *
gravity_);
733 if (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potential, 0.0, 1.0e-30))
735 lambda = 0.5*(lambdaI + lambdaBound);
737 else if (potential > 0.)
743 lambda = lambdaBound;
747 Scalar entry(0.), rightEntry(0.);
748 entry = lambda * (fabs(permeability * unitOuterNormal) / dist) * faceArea;
751 rightEntry = lambda * rhoMean * fabs(permeability * unitOuterNormal)
756 entries[1] += entry * primaryVariablesOnBoundary[Indices::pressureEqIdx];
757 entries[1] -= rightEntry * (
gravity_ * unitDistVec);
764 else if(bcType.isNeumann(Indices::pressureEqIdx))
766 PrimaryVariables J(NAN);
767 problem().neumann(J, intersection);
768 J[1+phaseIdx] /= cellDataI.density(phaseIdx);
770 entries[1] -= J[1+phaseIdx] * faceArea;
773 DUNE_THROW(Dune::NotImplemented,
"Boundary Condition neither Dirichlet nor Neumann!");
791template<
class TypeTag>
795 Scalar maxError = 0.;
798 if (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(problem().timeManager().time(), 0.0, 1.0e-30))
804 for (
const auto& element : elements(problem().gridView()))
807 int eIdxGlobal = problem().variables().index(element);
808 CellData& cellData = problem().variables().cellData(eIdxGlobal);
810 if(cellData.subdomain() == 2)
818 PrimaryVariables source(NAN);
819 problem().source(source, element);
821 if ((cellData.saturation(wPhaseIdx) > 0.0 && cellData.saturation(wPhaseIdx) < 1.0)
822 || Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(source.one_norm(), 0.0, 1.0e-30))
828 for (
const auto& intersection : intersections(problem().gridView(), element))
830 if (intersection.neighbor())
832 int eIdxGlobalJ = problem().variables().index(intersection.outside());
840 if(Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellData.saturation(wPhaseIdx), 0.0, 1.0e-30))
842 else if (Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellData.saturation(nPhaseIdx), 0.0, 1.0e-30))
858 problem().gridView().template communicate<DataHandle>(dataHandle,
859 Dune::InteriorBorder_All_Interface,
860 Dune::ForwardCommunication);
865 for (
const auto& element : elements(problem().gridView()))
867 int eIdxGlobal = problem().variables().index(element);
868 CellData& cellData = problem().variables().cellData(eIdxGlobal);
871 int oldSubdomainI = cellData.subdomain();
875 if(oldSubdomainI != 2
883 else if(oldSubdomainI != 2
895 maxError = max(maxError, fabs(cellData.volumeError()));
897 this->
maxError_ = maxError/problem().timeManager().timeStepSize();
901 if(problem().timeManager().willBeFinished() or problem().timeManager().episodeWillBeFinished())
902 Dune::dinfo <<
"Subdomain routines took " <<
timer_.elapsed() <<
" seconds" << std::endl;
917template<
class TypeTag>
921 GlobalPosition globalPos = elementI.geometry().center();
922 int eIdxGlobal = problem().variables().index(elementI);
925 int presentPhaseIdx = cellData.subdomain();
932 auto& pseudoFluidState = cellData.manipulateSimpleFluidState();
940 pc = MaterialLaw::pc(problem().spatialParams().materialLawParams(elementI),
941 ((presentPhaseIdx == wPhaseIdx) ? 1. : 0.));
944 pressure[wPhaseIdx] = this->pressure(eIdxGlobal);
945 pressure[nPhaseIdx] = this->pressure(eIdxGlobal)+pc;
949 pressure[wPhaseIdx] = this->pressure(eIdxGlobal)-pc;
950 pressure[nPhaseIdx] = this->pressure(eIdxGlobal);
954 Scalar sumConc = cellData.massConcentration(wCompIdx)
955 + cellData.massConcentration(nCompIdx);
956 Scalar Z0 = cellData.massConcentration(wCompIdx)/ sumConc;
962 assert(presentPhaseIdx == pseudoFluidState.presentPhaseIdx());
967 cellData.setViscosity(presentPhaseIdx, FluidSystem::viscosity(pseudoFluidState, presentPhaseIdx));
970 if(presentPhaseIdx == wPhaseIdx)
972 cellData.setMobility(wPhaseIdx,
973 MaterialLaw::krw(problem().spatialParams().materialLawParams(elementI), pseudoFluidState.saturation(wPhaseIdx))
974 / cellData.viscosity(wPhaseIdx));
975 cellData.setMobility(nPhaseIdx, 0.);
979 cellData.setMobility(nPhaseIdx,
980 MaterialLaw::krn(problem().spatialParams().materialLawParams(elementI), pseudoFluidState.saturation(wPhaseIdx))
981 / cellData.viscosity(nPhaseIdx));
982 cellData.setMobility(wPhaseIdx, 0.);
987 vol = sumConc / pseudoFluidState.density(presentPhaseIdx);
989 if (Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(problem().timeManager().timeStepSize(), 0.0, 1.0e-30))
990 cellData.volumeError() = (vol - problem().spatialParams().porosity(elementI));
Contains a class to exchange entries of a vector.
Determines the pressures and saturations of all fluid phases given the total mass of all components.
Finite volume 2p2c pressure model.
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition math.hh:68
void harmonicMeanMatrix(Dune::FieldMatrix< Scalar, m, n > &K, const Dune::FieldMatrix< Scalar, m, n > &Ki, const Dune::FieldMatrix< Scalar, m, n > &Kj)
Calculate the harmonic mean of a fixed-size matrix.
Definition math.hh:108
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:153
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition propertysystem.hh:149
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property (equivalent to old macro GET_PROP(...))
Definition propertysystem.hh:140
Flash calculation routines for compositional sequential models.
Definition compositionalflash.hh:51
static void concentrationFlash1p2c(FluidState1p2c &fluidState, const Scalar &Z0, const Dune::FieldVector< Scalar, numPhases > phasePressure, const int presentPhaseIdx, const Scalar &temperature)
The simplest possible update routine for 1p2c "flash" calculations.
Definition compositionalflash.hh:181
Container for compositional variables in a 1p2c situation.
Definition pseudo1p2c.hh:44
Scalar density(int phaseIdx) const
Set the density of a phase .
Definition pseudo1p2c.hh:109
Scalar temperature(int phaseIdx) const
Returns the temperature of the fluids .
Definition pseudo1p2c.hh:189
Scalar ErrorTermLowerBound_
Handling of error term: lower bound for error dampening.
Definition fvpressure.hh:187
Scalar ErrorTermFactor_
Handling of error term: relaxation factor.
Definition fvpressure.hh:186
Problem & problem_
Definition fvpressure.hh:182
bool regulateBoundaryPermeability
Enables regulation of permeability in the direction of a Dirichlet Boundary Condition.
Definition fvpressure.hh:184
Problem & problem()
Definition fvpressure.hh:138
Scalar minimalBoundaryPermeability
Minimal limit for the boundary permeability.
Definition fvpressure.hh:185
Scalar ErrorTermUpperBound_
Definition fvpressure.hh:188
FVPressure2P2C(Problem &problem)
Constructs a FVPressure2P2C object.
Definition fvpressure.hh:164
void updateMaterialLawsInElement(const Element &elementI, bool postTimeStep)
Updates secondary variables of one cell.
Definition fvpressure.hh:904
Scalar incp_
Increment for the volume derivative w.r.t pressure.
Definition fvpressurecompositional.hh:399
Scalar maxError_
Maximum volume error of all cells.
Definition fvpressurecompositional.hh:398
void addOutputVtkFields(MultiWriter &writer)
Write data files.
Definition fvpressurecompositional.hh:172
Dune::BlockVector< Dune::FieldVector< int, 1 > > nextSubdomain
vector holding next subdomain
Definition fvpressuremultiphysics.hh:230
void get1pFluxOnBoundary(EntryType &entries, const Intersection &intersection, const CellData &cellDataI)
The compositional single-phase flux in the multiphysics framework.
Definition fvpressuremultiphysics.hh:630
const GlobalPosition & gravity_
Definition fvpressuremultiphysics.hh:231
void update1pMaterialLawsInElement(const Element &elementI, CellData &cellData, bool postTimeStep)
updates secondary variables of one single phase cell
Definition fvpressuremultiphysics.hh:918
void assemble(bool first)
function which assembles the system of equations to be solved
Definition fvpressuremultiphysics.hh:264
void addOutputVtkFields(MultiWriter &writer)
Write data files.
Definition fvpressuremultiphysics.hh:195
FVPressure2P2CMultiPhysics(Problem &problem)
Constructs a FVPressure2P2CPC object.
Definition fvpressuremultiphysics.hh:219
@ rhs
index for the right hand side entry
Definition fvpressuremultiphysics.hh:245
@ matrix
index for the global matrix entry
Definition fvpressuremultiphysics.hh:246
VectorExchange< ElementMapper, Dune::BlockVector< Dune::FieldVector< int, 1 > > > DataHandle
Definition fvpressuremultiphysics.hh:226
void get1pStorage(EntryType &storageEntry, const Element &elementI, CellData &cellDataI)
Assembles the storage term for a 1p cell in a multiphysics framework.
Definition fvpressuremultiphysics.hh:404
void updateMaterialLaws(bool postTimeStep=false)
constitutive functions are updated once if new concentrations are calculated and stored in the variab...
Definition fvpressuremultiphysics.hh:792
void get1pSource(EntryType &sourceEntry, const Element &elementI, const CellData &cellDataI)
Assembles the source term.
Definition fvpressuremultiphysics.hh:371
Dune::Timer timer_
A timer for the time spent on the multiphysics framework.
Definition fvpressuremultiphysics.hh:234
void serializeEntity(std::ostream &outstream, const Element &element)
Function for serialization of the pressure field.
Definition fvpressuremultiphysics.hh:157
void deserializeEntity(std::istream &instream, const Element &element)
Function for deserialization of the pressure field.
Definition fvpressuremultiphysics.hh:173
static constexpr int pressureType
gives kind of pressure used ( , , )
Definition fvpressuremultiphysics.hh:233
void initialize(bool solveTwice=false)
Definition fvpressuremultiphysics.hh:136
typename SolutionTypes::ElementMapper ElementMapper
Definition fvpressuremultiphysics.hh:225
void get1pFlux(EntryType &entries, const Intersection &intersection, const CellData &cellDataI)
The compositional single-phase flux in the multiphysics framework.
Definition fvpressuremultiphysics.hh:527
void assemble(bool first)
Function which assembles the system of equations to be solved.
Definition sequential/cellcentered/pressure.hh:401
void initialize()
Initialize pressure model.
Definition sequential/cellcentered/pressure.hh:213
void deserializeEntity(std::istream &instream, const Element &element)
Function for deserialization of the pressure field.
Definition sequential/cellcentered/pressure.hh:264
PressureSolution & pressure()
Returns the vector containing the pressure solution.
Definition sequential/cellcentered/pressure.hh:120
void serializeEntity(std::ostream &outstream, const Element &element)
Function for serialization of the pressure field.
Definition sequential/cellcentered/pressure.hh:251
Matrix A_
Global stiffnes matrix (sparse matrix which is build by the initializeMatrix() function).
Definition sequential/cellcentered/pressure.hh:324
RHSVector f_
Right hand side vector.
Definition sequential/cellcentered/pressure.hh:325