24#ifndef DUMUX_MIMETICOPERATOR2P_HH
25#define DUMUX_MIMETICOPERATOR2P_HH
45template<
class TypeTag>
55 dim = GridView::dimension, dimWorld = GridView::dimensionworld,
59 using Element =
typename GridView::template Codim<0>::Entity;
66 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
72 pw = Indices::pressureW,
73 pn = Indices::pressureNw,
74 pressureType = getPropValue<TypeTag, Properties::PressureFormulation>(),
75 wPhaseIdx = Indices::wPhaseIdx,
76 nPhaseIdx = Indices::nPhaseIdx,
77 saturationIdx = Indices::saturationIdx,
78 satEqIdx = Indices::satEqIdx,
79 pressureEqIdx = Indices::pressureEqIdx
82 using FieldVector = Dune::FieldVector<Scalar, dimWorld>;
92 template<
class Vector>
95 Dune::FieldVector<Scalar, 2 * dim> velocityW(0);
96 Dune::FieldVector<Scalar, 2 * dim> velocityNw(0);
97 Dune::FieldVector<Scalar, 2 * dim> pressTrace(0);
98 Dune::FieldVector<Scalar, 2 * dim> gravPotTrace(0);
100 const auto firstElement = *this->
gridView_.template begin<0>();
101 FluidState fluidState;
102 fluidState.setPressure(wPhaseIdx, problem.referencePressure(firstElement));
103 fluidState.setPressure(nPhaseIdx, problem.referencePressure(firstElement));
104 fluidState.setTemperature(problem.temperature(firstElement));
105 fluidState.setSaturation(wPhaseIdx, 1.);
106 fluidState.setSaturation(nPhaseIdx, 0.);
112 for (
int i = 0; i < problem.gridView().size(0); i++)
114 problem.variables().cellData(i).fluxData().resetVelocity();
118 for (
const auto& element : elements(this->
gridView_))
120 int eIdxGlobal = problem.variables().index(element);
122 CellData& cellData = problem.variables().cellData(eIdxGlobal);
123 FieldVector globalPos = element.geometry().center();
126 for (
const auto& intersection : intersections(problem.gridView(), element))
128 int indexInInside = intersection.indexInInside();
130 int fIdxGlobal = this->
faceMapper_.subIndex(element, indexInInside, 1);
132 pressTrace[indexInInside] = u[fIdxGlobal];
134 gravPotTrace[indexInInside] = (problem.bBoxMax() - intersection.geometry().center()) * problem.gravity() * densityDiff;
137 switch (pressureType)
141 Scalar potW = loc.constructPressure(element, pressTrace);
142 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
143 Scalar potNw = potW + gravPot;
145 cellData.setPotential(wPhaseIdx, potW);
146 cellData.setPotential(nPhaseIdx, potNw);
148 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, wPhaseIdx);
150 cellData.setPressure(wPhaseIdx, potW - gravPot);
152 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, nPhaseIdx);
154 cellData.setPressure(nPhaseIdx, potNw - gravPot);
160 Scalar potNw = loc.constructPressure(element, pressTrace);
161 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
162 Scalar potW = potNw - gravPot;
164 cellData.setPotential(nPhaseIdx, potNw);
165 cellData.setPotential(wPhaseIdx, potW);
167 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, wPhaseIdx);
169 cellData.setPressure(wPhaseIdx, potW - gravPot);
171 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, nPhaseIdx);
173 cellData.setPressure(nPhaseIdx, potNw - gravPot);
182 switch (pressureType)
186 loc.constructVelocity(element, velocityW, pressTrace, cellData.potential(wPhaseIdx));
187 pressTrace += gravPotTrace;
188 loc.constructVelocity(element, velocityNw, pressTrace, cellData.potential(nPhaseIdx));
194 loc.constructVelocity(element, velocityW, pressTrace, cellData.potential(nPhaseIdx));
195 pressTrace -= gravPotTrace;
196 loc.constructVelocity(element, velocityNw, pressTrace, cellData.potential(wPhaseIdx));
202 for (
const auto& intersection : intersections(problem.gridView(), element))
204 int idxInInside = intersection.indexInInside();
206 cellData.fluxData().setUpwindPotential(wPhaseIdx, idxInInside, velocityW[idxInInside]);
207 cellData.fluxData().setUpwindPotential(nPhaseIdx, idxInInside, velocityNw[idxInInside]);
209 Scalar mobilityW = 0;
210 Scalar mobilityNw = 0;
212 if (intersection.neighbor())
214 int neighborIdx = problem.variables().index(intersection.outside());
216 CellData& cellDataNeighbor = problem.variables().cellData(neighborIdx);
219 (velocityW[idxInInside] >= 0.) ? cellData.mobility(wPhaseIdx) :
220 cellDataNeighbor.mobility(wPhaseIdx);
222 (velocityNw[idxInInside] >= 0.) ? cellData.mobility(nPhaseIdx) :
223 cellDataNeighbor.mobility(nPhaseIdx);
225 if (velocityW[idxInInside] >= 0.)
227 FieldVector velocity(intersection.centerUnitOuterNormal());
228 velocity *= mobilityW/(mobilityW+mobilityNw) * velocityW[idxInInside];
229 cellData.fluxData().addVelocity(wPhaseIdx, idxInInside, velocity);
230 cellDataNeighbor.fluxData().addVelocity(wPhaseIdx, intersection.indexInOutside(), velocity);
232 if (velocityNw[idxInInside] >= 0.)
234 FieldVector velocity(intersection.centerUnitOuterNormal());
235 velocity *= mobilityNw/(mobilityW+mobilityNw) * velocityNw[idxInInside];
236 cellData.fluxData().addVelocity(nPhaseIdx, idxInInside, velocity);
237 cellDataNeighbor.fluxData().addVelocity(nPhaseIdx, intersection.indexInOutside(), velocity);
240 cellData.fluxData().setVelocityMarker(idxInInside);
244 BoundaryTypes bctype;
245 problem.boundaryTypes(bctype, intersection);
246 if (bctype.isDirichlet(satEqIdx))
248 PrimaryVariables boundValues(0.0);
249 problem.dirichlet(boundValues, intersection);
254 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem.spatialParams(), element);
256 if (velocityW[idxInInside] >= 0.)
258 mobilityW = cellData.mobility(wPhaseIdx);
262 mobilityW = fluidMatrixInteraction.krw(boundValues[saturationIdx]) / viscosityW;
265 if (velocityNw[idxInInside] >= 0.)
267 mobilityNw = cellData.mobility(nPhaseIdx);
271 mobilityNw = fluidMatrixInteraction.krn(boundValues[saturationIdx]) / viscosityNw;
276 mobilityW = cellData.mobility(wPhaseIdx);
277 mobilityNw = cellData.mobility(nPhaseIdx);
280 FieldVector velocity(intersection.centerUnitOuterNormal());
281 velocity *= mobilityW / (mobilityW + mobilityNw) * velocityW[idxInInside];
282 cellData.fluxData().setVelocity(wPhaseIdx, idxInInside, velocity);
285 velocity = intersection.centerUnitOuterNormal();
286 velocity *= mobilityNw / (mobilityW + mobilityNw) * velocityNw[idxInInside];
287 cellData.fluxData().setVelocity(nPhaseIdx, idxInInside, velocity);
288 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 (equivalent to old macro GET_PROP(...))
Definition: propertysystem.hh:140
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
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:66
const GridView gridView_
Definition: croperator.hh:327
FaceMapper faceMapper_
Definition: croperator.hh:328
Levelwise assembler.
Definition: operator.hh:47
MimeticOperatorAssemblerTwoP(const GridView &gridView)
Definition: operator.hh:86
void calculatePressure(LocalStiffness &loc, Vector &u, Problem &problem)
Definition: operator.hh:93
Specifies the properties for immiscible 2p diffusion/pressure models.
Defines the basic properties required for a mimetic method.