24#ifndef DUMUX_FVPRESSURE2P2C_MULTIPHYSICS_HH
25#define DUMUX_FVPRESSURE2P2C_MULTIPHYSICS_HH
27#include <dune/common/float_cmp.hh>
68template<
class TypeTag>
88 dim = GridView::dimension, dimWorld = GridView::dimensionworld
92 pw = Indices::pressureW
96 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
97 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
98 contiWEqIdx = Indices::contiWEqIdx, contiNEqIdx = Indices::contiNEqIdx
102 using Element =
typename GridView::Traits::template Codim<0>::Entity;
103 using Grid =
typename GridView::Grid;
104 using Intersection =
typename GridView::Intersection;
107 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
108 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
109 using PhaseVector = Dune::FieldVector<Scalar, getPropValue<TypeTag, Properties::NumPhases>()>;
113 using EntryType = Dune::FieldVector<Scalar, 2>;
123 void get1pSource(EntryType& sourceEntry,
const Element& elementI,
const CellData& cellDataI);
125 void get1pStorage(EntryType& storageEntry,
const Element& elementI, CellData& cellDataI);
127 void get1pFlux(EntryType& entries,
const Intersection& intersection,
128 const CellData& cellDataI);
131 const Intersection& intersection,
132 const CellData& cellDataI);
138 int size = this->
problem().gridView().size(0);
139 for (
int i = 0; i < size; i++)
141 CellData& cellData = this->
problem().variables().cellData(i);
142 cellData.subdomain() = 2;
159 int eIdxGlobal =
problem().variables().index(element);
160 CellData& cellData =
problem().variables().cellData(eIdxGlobal);
161 outstream <<
" "<< cellData.subdomain();
176 int eIdxGlobal =
problem().variables().index(element);
177 CellData& cellData =
problem().variables().cellData(eIdxGlobal);
179 instream >> subdomainIdx;
180 cellData.setSubdomainAndFluidStateType(subdomainIdx);
193 template<
class MultiWriter>
198 if(
problem().vtkOutputLevel()>=1)
200 int size =
problem().gridView().size(0);
202 Dune::BlockVector<Dune::FieldVector<int,1> >* subdomainPtr = writer.template allocateManagedBuffer<int, 1> (size);
203 for (
int i = 0; i < size; i++)
205 CellData& cellData =
problem().variables().cellData(i);
206 (*subdomainPtr)[i] = cellData.subdomain();
208 writer.attachCellData(*subdomainPtr,
"subdomain");
262template<
class TypeTag>
274 for (
const auto& element : elements(problem().gridView()))
277 int eIdxGlobalI = problem().variables().index(element);
280 if (element.partitionType() == Dune::InteriorEntity)
283 CellData& cellDataI = problem().variables().cellData(eIdxGlobalI);
285 Dune::FieldVector<Scalar, 2> entries(0.);
288 if(cellDataI.subdomain() != 2)
289 problem().pressureModel().get1pSource(entries,element, cellDataI);
291 problem().pressureModel().getSource(entries,element, cellDataI, first);
293 this->
f_[eIdxGlobalI] = entries[
rhs];
297 for (
const auto& intersection : intersections(problem().gridView(), element))
300 if (intersection.neighbor())
302 int eIdxGlobalJ = problem().variables().index(intersection.outside());
304 if (cellDataI.subdomain() != 2
305 or problem().variables().cellData(eIdxGlobalJ).subdomain() != 2)
306 get1pFlux(entries, intersection, cellDataI);
308 problem().pressureModel().getFlux(entries, intersection, cellDataI, first);
311 this->
f_[eIdxGlobalI] -= entries[
rhs];
313 this->
A_[eIdxGlobalI][eIdxGlobalI] += entries[
matrix];
315 this->
A_[eIdxGlobalI][eIdxGlobalJ] = -entries[
matrix];
322 if (cellDataI.subdomain() != 2)
323 problem().pressureModel().get1pFluxOnBoundary(entries, intersection, cellDataI);
325 problem().pressureModel().getFluxOnBoundary(entries, intersection, cellDataI, first);
328 this->
f_[eIdxGlobalI] += entries[
rhs];
330 this->
A_[eIdxGlobalI][eIdxGlobalI] += entries[
matrix];
336 if (cellDataI.subdomain() != 2)
337 problem().pressureModel().get1pStorage(entries, element, cellDataI);
339 problem().pressureModel().getStorage(entries, element, cellDataI, first);
341 this->
f_[eIdxGlobalI] += entries[
rhs];
343 this->
A_[eIdxGlobalI][eIdxGlobalI] += entries[
matrix];
348 this->
A_[eIdxGlobalI] = 0.0;
349 this->
A_[eIdxGlobalI][eIdxGlobalI] = 1.0;
350 this->
f_[eIdxGlobalI] = this->
pressure()[eIdxGlobalI];
369template<
class TypeTag>
371 const Element& elementI,
const CellData& cellDataI)
376 Scalar volume = elementI.geometry().volume();
377 int subdomainIdx = cellDataI.subdomain();
380 PrimaryVariables source(NAN);
381 problem().source(source, elementI);
382 source[1+subdomainIdx] /= cellDataI.density(subdomainIdx);
384 sourceEntry[1] = volume * source[1+subdomainIdx];
402template<
class TypeTag>
404 const Element& elementI,
409 int eIdxGlobalI = problem().variables().index(elementI);
410 int presentPhaseIdx = cellDataI.subdomain();
411 Scalar volume = elementI.geometry().volume();
414 Scalar timestep_ = problem().timeManager().timeStepSize();
419 Scalar& incp = this->
incp_;
422 PhaseVector p_(incp);
423 p_[nPhaseIdx] += cellDataI.pressure(nPhaseIdx);
424 p_[wPhaseIdx] += cellDataI.pressure(wPhaseIdx);
426 Scalar sumC = (cellDataI.massConcentration(wCompIdx) + cellDataI.massConcentration(nCompIdx));
427 Scalar Z0 = cellDataI.massConcentration(wCompIdx) / sumC;
434 Scalar v_ = 1. / pseudoFluidState.
density(presentPhaseIdx);
435 cellDataI.dv_dp() = (sumC * ( v_ - (1. /cellDataI.density(presentPhaseIdx)))) /incp;
437 if (cellDataI.dv_dp()>0)
440 Dune::dinfo <<
"dv_dp larger 0 at Idx " << eIdxGlobalI <<
" , try and invert secant"<< std::endl;
445 v_ = 1. / pseudoFluidState.
density(presentPhaseIdx);
446 cellDataI.dv_dp() = (sumC * ( v_ - (1. /cellDataI.density(presentPhaseIdx)))) /incp;
448 if (cellDataI.dv_dp()>0)
450 Dune::dinfo <<__FILE__<<
"dv_dp still larger 0 after inverting secant. regularize"<< std::endl;
451 cellDataI.dv_dp() *= -1;
456 Scalar compress_term = cellDataI.dv_dp() / timestep_;
458 storageEntry[
matrix] -= compress_term*volume;
459 storageEntry[
rhs] -= cellDataI.pressure(
pressureType) * compress_term * volume;
463 if (isnan(compress_term) ||isinf(compress_term))
464 DUNE_THROW(Dune::MathError,
"Compressibility term leads to NAN matrix entry at index " << eIdxGlobalI);
467 DUNE_THROW(Dune::NotImplemented,
"Compressibility is switched off???");
472 if( problem().timeManager().episodeWillBeFinished()
473 || problem().timeManager().willBeFinished())
475 problem().variables().cellData(eIdxGlobalI).errorCorrection() = 0.;
481 problem().variables().cellData(eIdxGlobalI).volumeError() /= timestep_;
483 Scalar erri = fabs(cellDataI.volumeError());
492 if ((erri*timestep_ > 5e-5) && (erri > x_lo * maxError) && (!problem().timeManager().willBeFinished()))
494 if (erri <= x_mi * maxError)
496 problem().variables().cellData(eIdxGlobalI).errorCorrection() =
497 fac* (1-x_mi*(lofac-1)/(x_lo-x_mi) + (lofac-1)/(x_lo-x_mi)*erri/maxError)
498 * cellDataI.volumeError() * volume;
501 problem().variables().cellData(eIdxGlobalI).errorCorrection() =
502 fac * (1 + x_mi - hifac*x_mi/(1-x_mi) + (hifac/(1-x_mi)-1)*erri/maxError)
503 * cellDataI.volumeError() * volume;
506 problem().variables().cellData(eIdxGlobalI).errorCorrection()=0 ;
525template<
class TypeTag>
527 const Intersection& intersection,
const CellData& cellDataI)
530 auto elementI = intersection.inside();
533 const GlobalPosition& globalPos = elementI.geometry().center();
536 DimMatrix permeabilityI(problem().spatialParams().intrinsicPermeability(elementI));
539 const GlobalPosition& unitOuterNormal = intersection.centerUnitOuterNormal();
542 Scalar faceArea = intersection.geometry().volume();
545 auto neighbor = intersection.outside();
546 int eIdxGlobalJ = problem().variables().index(neighbor);
547 CellData& cellDataJ = problem().variables().cellData(eIdxGlobalJ);
550 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
553 GlobalPosition distVec = globalPosNeighbor - globalPos;
556 Scalar dist = distVec.two_norm();
558 GlobalPosition unitDistVec(distVec);
561 DimMatrix permeabilityJ
562 = problem().spatialParams().intrinsicPermeability(neighbor);
565 DimMatrix meanPermeability(0);
568 Dune::FieldVector<Scalar, dim> permeability(0);
569 meanPermeability.mv(unitDistVec, permeability);
574 int phaseIdx = min(cellDataI.subdomain(), cellDataJ.subdomain());
576 Scalar rhoMean = 0.5 * (cellDataI.density(phaseIdx) + cellDataJ.density(phaseIdx));
579 Scalar potential = (cellDataI.pressure(phaseIdx) - cellDataJ.pressure(phaseIdx)) / dist;
581 potential += rhoMean * (unitDistVec *
gravity_);
587 lambda = cellDataI.mobility(phaseIdx);
588 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiWEqIdx,
false);
589 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiNEqIdx,
false);
591 else if (potential < 0.)
593 lambda = cellDataJ.mobility(phaseIdx);
594 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiWEqIdx,
true);
595 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiNEqIdx,
true);
599 lambda =
harmonicMean(cellDataI.mobility(phaseIdx) , cellDataJ.mobility(phaseIdx));
600 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiWEqIdx,
false);
601 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiNEqIdx,
false);
604 entries[0] = lambda * faceArea * fabs(permeability * unitOuterNormal) / (dist);
605 entries[1] = rhoMean * lambda;
606 entries[1] *= faceArea * fabs(permeability * unitOuterNormal) * (unitDistVec *
gravity_);
628template<
class TypeTag>
630 const Intersection& intersection,
const CellData& cellDataI)
634 auto elementI = intersection.inside();
635 const GlobalPosition& globalPos = elementI.geometry().center();
637 int phaseIdx = cellDataI.subdomain();
640 const GlobalPosition& unitOuterNormal = intersection.centerUnitOuterNormal();
642 Scalar faceArea = intersection.geometry().volume();
645 const GlobalPosition& globalPosFace = intersection.geometry().center();
648 GlobalPosition distVec(globalPosFace - globalPos);
649 Scalar dist = distVec.two_norm();
650 GlobalPosition unitDistVec(distVec);
654 BoundaryTypes bcType;
655 problem().boundaryTypes(bcType, intersection);
658 PhaseVector pressBC(0.);
661 if (bcType.isDirichlet(Indices::pressureEqIdx))
664 DimMatrix permeabilityI(problem().spatialParams().intrinsicPermeability(elementI));
667 int axis = intersection.indexInInside() / 2;
672 Scalar lambdaI = cellDataI.mobility(phaseIdx);
675 Dune::FieldVector<Scalar, dim> permeability(0);
676 permeabilityI.mv(unitDistVec, permeability);
679 FluidState BCfluidState;
682 PrimaryVariables primaryVariablesOnBoundary(NAN);
683 problem().dirichlet(primaryVariablesOnBoundary, intersection);
687 problem().transportModel().evalBoundary(globalPosFace,
693 Scalar densityBound =
694 FluidSystem::density(BCfluidState, phaseIdx);
695 Scalar viscosityBound =
696 FluidSystem::viscosity(BCfluidState, phaseIdx);
699 Scalar lambdaBound = 0.;
702 case Indices::satDependent:
704 lambdaBound = BCfluidState.saturation(phaseIdx) / viscosityBound;
707 case Indices::permDependent:
709 const auto fluidMatrixInteraction = problem().spatialParams().fluidMatrixInteractionAtPos(elementI.geometry().center());
710 if (phaseIdx == wPhaseIdx)
711 lambdaBound = fluidMatrixInteraction.krw(BCfluidState.saturation(wPhaseIdx)) / viscosityBound;
713 lambdaBound = fluidMatrixInteraction.krn(BCfluidState.saturation(wPhaseIdx)) / viscosityBound;
717 Scalar rhoMean = 0.5 * (cellDataI.density(phaseIdx) + densityBound);
719 Scalar potential = 0;
722 potential = (cellDataI.pressure(phaseIdx) - pressBC[phaseIdx]) / dist;
724 potential += rhoMean * (unitDistVec *
gravity_);
728 if (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potential, 0.0, 1.0e-30))
730 lambda = 0.5*(lambdaI + lambdaBound);
732 else if (potential > 0.)
738 lambda = lambdaBound;
742 Scalar entry(0.), rightEntry(0.);
743 entry = lambda * (fabs(permeability * unitOuterNormal) / dist) * faceArea;
746 rightEntry = lambda * rhoMean * fabs(permeability * unitOuterNormal)
751 entries[1] += entry * primaryVariablesOnBoundary[Indices::pressureEqIdx];
752 entries[1] -= rightEntry * (
gravity_ * unitDistVec);
759 else if(bcType.isNeumann(Indices::pressureEqIdx))
761 PrimaryVariables J(NAN);
762 problem().neumann(J, intersection);
763 J[1+phaseIdx] /= cellDataI.density(phaseIdx);
765 entries[1] -= J[1+phaseIdx] * faceArea;
768 DUNE_THROW(Dune::NotImplemented,
"Boundary Condition neither Dirichlet nor Neumann!");
786template<
class TypeTag>
790 Scalar maxError = 0.;
793 if (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(problem().timeManager().time(), 0.0, 1.0e-30))
799 for (
const auto& element : elements(problem().gridView()))
802 int eIdxGlobal = problem().variables().index(element);
803 CellData& cellData = problem().variables().cellData(eIdxGlobal);
805 if(cellData.subdomain() == 2)
813 PrimaryVariables source(NAN);
814 problem().source(source, element);
816 if ((cellData.saturation(wPhaseIdx) > 0.0 && cellData.saturation(wPhaseIdx) < 1.0)
817 || Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(source.one_norm(), 0.0, 1.0e-30))
823 for (
const auto& intersection : intersections(problem().gridView(), element))
825 if (intersection.neighbor())
827 int eIdxGlobalJ = problem().variables().index(intersection.outside());
835 if(Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellData.saturation(wPhaseIdx), 0.0, 1.0e-30))
837 else if (Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellData.saturation(nPhaseIdx), 0.0, 1.0e-30))
853 problem().gridView().template communicate<DataHandle>(dataHandle,
854 Dune::InteriorBorder_All_Interface,
855 Dune::ForwardCommunication);
860 for (
const auto& element : elements(problem().gridView()))
862 int eIdxGlobal = problem().variables().index(element);
863 CellData& cellData = problem().variables().cellData(eIdxGlobal);
866 int oldSubdomainI = cellData.subdomain();
870 if(oldSubdomainI != 2
878 else if(oldSubdomainI != 2
890 maxError = max(maxError, fabs(cellData.volumeError()));
892 this->
maxError_ = maxError/problem().timeManager().timeStepSize();
896 if(problem().timeManager().willBeFinished() or problem().timeManager().episodeWillBeFinished())
897 Dune::dinfo <<
"Subdomain routines took " <<
timer_.elapsed() <<
" seconds" << std::endl;
912template<
class TypeTag>
916 GlobalPosition globalPos = elementI.geometry().center();
917 int eIdxGlobal = problem().variables().index(elementI);
920 int presentPhaseIdx = cellData.subdomain();
927 auto& pseudoFluidState = cellData.manipulateSimpleFluidState();
929 const auto fluidMatrixInteraction = problem().spatialParams().fluidMatrixInteractionAtPos(elementI.geometry().center());
937 pc = fluidMatrixInteraction.pc(((presentPhaseIdx == wPhaseIdx) ? 1. : 0.));
940 pressure[wPhaseIdx] = this->pressure(eIdxGlobal);
941 pressure[nPhaseIdx] = this->pressure(eIdxGlobal)+pc;
945 pressure[wPhaseIdx] = this->pressure(eIdxGlobal)-pc;
946 pressure[nPhaseIdx] = this->pressure(eIdxGlobal);
950 Scalar sumConc = cellData.massConcentration(wCompIdx)
951 + cellData.massConcentration(nCompIdx);
952 Scalar Z0 = cellData.massConcentration(wCompIdx)/ sumConc;
958 assert(presentPhaseIdx == pseudoFluidState.presentPhaseIdx());
963 cellData.setViscosity(presentPhaseIdx, FluidSystem::viscosity(pseudoFluidState, presentPhaseIdx));
966 if(presentPhaseIdx == wPhaseIdx)
968 cellData.setMobility(wPhaseIdx,
969 fluidMatrixInteraction.krw(pseudoFluidState.saturation(wPhaseIdx)) / cellData.viscosity(wPhaseIdx));
970 cellData.setMobility(nPhaseIdx, 0.);
974 cellData.setMobility(nPhaseIdx,
975 fluidMatrixInteraction.krn(pseudoFluidState.saturation(wPhaseIdx)) / cellData.viscosity(nPhaseIdx));
976 cellData.setMobility(wPhaseIdx, 0.);
981 vol = sumConc / pseudoFluidState.density(presentPhaseIdx);
983 if (Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(problem().timeManager().timeStepSize(), 0.0, 1.0e-30))
984 cellData.volumeError() = (vol - problem().spatialParams().porosity(elementI));
Determines the pressures and saturations of all fluid phases given the total mass of all components.
Contains a class to exchange entries of a vector.
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:69
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:109
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:154
VectorCommDataHandle< Mapper, Vector, codim, Detail::SetEqual, DataType > VectorCommDataHandleEqual
Definition vectorcommdatahandle.hh:126
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition propertysystem.hh:150
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property
Definition propertysystem.hh:141
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:179
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:186
Scalar ErrorTermFactor_
Handling of error term: relaxation factor.
Definition fvpressure.hh:185
Problem & problem_
Definition fvpressure.hh:181
bool regulateBoundaryPermeability
Enables regulation of permeability in the direction of a Dirichlet Boundary Condition.
Definition fvpressure.hh:183
Problem & problem()
Definition fvpressure.hh:137
Scalar minimalBoundaryPermeability
Minimal limit for the boundary permeability.
Definition fvpressure.hh:184
Scalar ErrorTermUpperBound_
Definition fvpressure.hh:187
FVPressure2P2C(Problem &problem)
Constructs a FVPressure2P2C object.
Definition fvpressure.hh:163
void updateMaterialLawsInElement(const Element &elementI, bool postTimeStep)
Updates secondary variables of one cell.
Definition fvpressure.hh:903
Scalar incp_
Increment for the volume derivative w.r.t pressure.
Definition fvpressurecompositional.hh:396
Scalar maxError_
Maximum volume error of all cells.
Definition fvpressurecompositional.hh:395
void addOutputVtkFields(MultiWriter &writer)
Write data files.
Definition fvpressurecompositional.hh:169
Dune::BlockVector< Dune::FieldVector< int, 1 > > nextSubdomain
vector holding next subdomain
Definition fvpressuremultiphysics.hh:229
void get1pFluxOnBoundary(EntryType &entries, const Intersection &intersection, const CellData &cellDataI)
The compositional single-phase flux in the multiphysics framework.
Definition fvpressuremultiphysics.hh:629
const GlobalPosition & gravity_
Definition fvpressuremultiphysics.hh:230
void update1pMaterialLawsInElement(const Element &elementI, CellData &cellData, bool postTimeStep)
updates secondary variables of one single phase cell
Definition fvpressuremultiphysics.hh:913
void assemble(bool first)
function which assembles the system of equations to be solved
Definition fvpressuremultiphysics.hh:263
void addOutputVtkFields(MultiWriter &writer)
Write data files.
Definition fvpressuremultiphysics.hh:194
VectorCommDataHandleEqual< ElementMapper, Dune::BlockVector< Dune::FieldVector< int, 1 > >, 0 > DataHandle
Definition fvpressuremultiphysics.hh:225
FVPressure2P2CMultiPhysics(Problem &problem)
Constructs a FVPressure2P2CPC object.
Definition fvpressuremultiphysics.hh:218
@ rhs
index for the right hand side entry
Definition fvpressuremultiphysics.hh:244
@ matrix
index for the global matrix entry
Definition fvpressuremultiphysics.hh:245
void get1pStorage(EntryType &storageEntry, const Element &elementI, CellData &cellDataI)
Assembles the storage term for a 1p cell in a multiphysics framework.
Definition fvpressuremultiphysics.hh:403
void updateMaterialLaws(bool postTimeStep=false)
constitutive functions are updated once if new concentrations are calculated and stored in the variab...
Definition fvpressuremultiphysics.hh:787
void get1pSource(EntryType &sourceEntry, const Element &elementI, const CellData &cellDataI)
Assembles the source term.
Definition fvpressuremultiphysics.hh:370
Dune::Timer timer_
A timer for the time spent on the multiphysics framework.
Definition fvpressuremultiphysics.hh:233
void serializeEntity(std::ostream &outstream, const Element &element)
Function for serialization of the pressure field.
Definition fvpressuremultiphysics.hh:156
void deserializeEntity(std::istream &instream, const Element &element)
Function for deserialization of the pressure field.
Definition fvpressuremultiphysics.hh:172
static constexpr int pressureType
gives kind of pressure used ( , , )
Definition fvpressuremultiphysics.hh:232
void initialize(bool solveTwice=false)
Definition fvpressuremultiphysics.hh:135
typename SolutionTypes::ElementMapper ElementMapper
Definition fvpressuremultiphysics.hh:224
void get1pFlux(EntryType &entries, const Intersection &intersection, const CellData &cellDataI)
The compositional single-phase flux in the multiphysics framework.
Definition fvpressuremultiphysics.hh:526
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