24#ifndef DUMUX_MIMETICOPERATOR2P_HH
25#define DUMUX_MIMETICOPERATOR2P_HH
43template<
class TypeTag>
53 dim = GridView::dimension, dimWorld = GridView::dimensionworld,
57 using Element =
typename GridView::template Codim<0>::Entity;
64 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
70 pw = Indices::pressureW,
71 pn = Indices::pressureNw,
72 pressureType = getPropValue<TypeTag, Properties::PressureFormulation>(),
73 wPhaseIdx = Indices::wPhaseIdx,
74 nPhaseIdx = Indices::nPhaseIdx,
75 saturationIdx = Indices::saturationIdx,
76 satEqIdx = Indices::satEqIdx,
77 pressureEqIdx = Indices::pressureEqIdx
80 using FieldVector = Dune::FieldVector<Scalar, dimWorld>;
90 template<
class Vector>
93 Dune::FieldVector<Scalar, 2 * dim> velocityW(0);
94 Dune::FieldVector<Scalar, 2 * dim> velocityNw(0);
95 Dune::FieldVector<Scalar, 2 * dim> pressTrace(0);
96 Dune::FieldVector<Scalar, 2 * dim> gravPotTrace(0);
98 const auto firstElement = *this->
gridView_.template begin<0>();
99 FluidState fluidState;
100 fluidState.setPressure(wPhaseIdx, problem.referencePressure(firstElement));
101 fluidState.setPressure(nPhaseIdx, problem.referencePressure(firstElement));
102 fluidState.setTemperature(problem.temperature(firstElement));
103 fluidState.setSaturation(wPhaseIdx, 1.);
104 fluidState.setSaturation(nPhaseIdx, 0.);
110 for (
int i = 0; i < problem.gridView().size(0); i++)
112 problem.variables().cellData(i).fluxData().resetVelocity();
116 for (
const auto& element : elements(this->
gridView_))
118 int eIdxGlobal = problem.variables().index(element);
120 CellData& cellData = problem.variables().cellData(eIdxGlobal);
121 FieldVector globalPos = element.geometry().center();
124 for (
const auto& intersection : intersections(problem.gridView(), element))
126 int indexInInside = intersection.indexInInside();
128 int fIdxGlobal = this->
faceMapper_.subIndex(element, indexInInside, 1);
130 pressTrace[indexInInside] = u[fIdxGlobal];
132 gravPotTrace[indexInInside] = (problem.bBoxMax() - intersection.geometry().center()) * problem.gravity() * densityDiff;
135 switch (pressureType)
139 Scalar potW = loc.constructPressure(element, pressTrace);
140 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
141 Scalar potNw = potW + gravPot;
143 cellData.setPotential(wPhaseIdx, potW);
144 cellData.setPotential(nPhaseIdx, potNw);
146 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, wPhaseIdx);
148 cellData.setPressure(wPhaseIdx, potW - gravPot);
150 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, nPhaseIdx);
152 cellData.setPressure(nPhaseIdx, potNw - gravPot);
158 Scalar potNw = loc.constructPressure(element, pressTrace);
159 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
160 Scalar potW = potNw - gravPot;
162 cellData.setPotential(nPhaseIdx, potNw);
163 cellData.setPotential(wPhaseIdx, potW);
165 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, wPhaseIdx);
167 cellData.setPressure(wPhaseIdx, potW - gravPot);
169 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, nPhaseIdx);
171 cellData.setPressure(nPhaseIdx, potNw - gravPot);
180 switch (pressureType)
184 loc.constructVelocity(element, velocityW, pressTrace, cellData.potential(wPhaseIdx));
185 pressTrace += gravPotTrace;
186 loc.constructVelocity(element, velocityNw, pressTrace, cellData.potential(nPhaseIdx));
192 loc.constructVelocity(element, velocityW, pressTrace, cellData.potential(nPhaseIdx));
193 pressTrace -= gravPotTrace;
194 loc.constructVelocity(element, velocityNw, pressTrace, cellData.potential(wPhaseIdx));
200 for (
const auto& intersection : intersections(problem.gridView(), element))
202 int idxInInside = intersection.indexInInside();
204 cellData.fluxData().setUpwindPotential(wPhaseIdx, idxInInside, velocityW[idxInInside]);
205 cellData.fluxData().setUpwindPotential(nPhaseIdx, idxInInside, velocityNw[idxInInside]);
207 Scalar mobilityW = 0;
208 Scalar mobilityNw = 0;
210 if (intersection.neighbor())
212 int neighborIdx = problem.variables().index(intersection.outside());
214 CellData& cellDataNeighbor = problem.variables().cellData(neighborIdx);
217 (velocityW[idxInInside] >= 0.) ? cellData.mobility(wPhaseIdx) :
218 cellDataNeighbor.mobility(wPhaseIdx);
220 (velocityNw[idxInInside] >= 0.) ? cellData.mobility(nPhaseIdx) :
221 cellDataNeighbor.mobility(nPhaseIdx);
223 if (velocityW[idxInInside] >= 0.)
225 FieldVector velocity(intersection.centerUnitOuterNormal());
226 velocity *= mobilityW/(mobilityW+mobilityNw) * velocityW[idxInInside];
227 cellData.fluxData().addVelocity(wPhaseIdx, idxInInside, velocity);
228 cellDataNeighbor.fluxData().addVelocity(wPhaseIdx, intersection.indexInOutside(), velocity);
230 if (velocityNw[idxInInside] >= 0.)
232 FieldVector velocity(intersection.centerUnitOuterNormal());
233 velocity *= mobilityNw/(mobilityW+mobilityNw) * velocityNw[idxInInside];
234 cellData.fluxData().addVelocity(nPhaseIdx, idxInInside, velocity);
235 cellDataNeighbor.fluxData().addVelocity(nPhaseIdx, intersection.indexInOutside(), velocity);
238 cellData.fluxData().setVelocityMarker(idxInInside);
242 BoundaryTypes bctype;
243 problem.boundaryTypes(bctype, intersection);
244 if (bctype.isDirichlet(satEqIdx))
246 PrimaryVariables boundValues(0.0);
247 problem.dirichlet(boundValues, intersection);
249 const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
250 if (velocityW[idxInInside] >= 0.)
252 mobilityW = cellData.mobility(wPhaseIdx);
256 mobilityW = fluidMatrixInteraction.krw(boundValues[saturationIdx]) / viscosityW;
259 if (velocityNw[idxInInside] >= 0.)
261 mobilityNw = cellData.mobility(nPhaseIdx);
265 mobilityNw = fluidMatrixInteraction.krn(boundValues[saturationIdx]) / viscosityNw;
270 mobilityW = cellData.mobility(wPhaseIdx);
271 mobilityNw = cellData.mobility(nPhaseIdx);
274 FieldVector velocity(intersection.centerUnitOuterNormal());
275 velocity *= mobilityW / (mobilityW + mobilityNw) * velocityW[idxInInside];
276 cellData.fluxData().setVelocity(wPhaseIdx, idxInInside, velocity);
279 velocity = intersection.centerUnitOuterNormal();
280 velocity *= mobilityNw / (mobilityW + mobilityNw) * velocityNw[idxInInside];
281 cellData.fluxData().setVelocity(nPhaseIdx, idxInInside, velocity);
282 cellData.fluxData().setVelocityMarker(idxInInside);
Defines a class for Crozieux-Raviart piecewise linear finite element functions.
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property
Definition: propertysystem.hh:141
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Extends CROperatorBase by a generic methods to assemble global stiffness matrix from local stiffness ...
Definition: croperator.hh:67
const GridView gridView_
Definition: croperator.hh:331
FaceMapper faceMapper_
Definition: croperator.hh:332
Levelwise assembler.
Definition: operator.hh:45
MimeticOperatorAssemblerTwoP(const GridView &gridView)
Definition: operator.hh:84
void calculatePressure(LocalStiffness &loc, Vector &u, Problem &problem)
Definition: operator.hh:91
Specifies the properties for immiscible 2p diffusion/pressure models.
Defines the basic properties required for a mimetic method.