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;
65 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
71 pw = Indices::pressureW,
72 pn = Indices::pressureNw,
73 pressureType = getPropValue<TypeTag, Properties::PressureFormulation>(),
74 wPhaseIdx = Indices::wPhaseIdx,
75 nPhaseIdx = Indices::nPhaseIdx,
76 saturationIdx = Indices::saturationIdx,
77 satEqIdx = Indices::satEqIdx,
78 pressureEqIdx = Indices::pressureEqIdx
81 using FieldVector = Dune::FieldVector<Scalar, dimWorld>;
91 template<
class Vector>
94 Dune::FieldVector<Scalar, 2 * dim> velocityW(0);
95 Dune::FieldVector<Scalar, 2 * dim> velocityNw(0);
96 Dune::FieldVector<Scalar, 2 * dim> pressTrace(0);
97 Dune::FieldVector<Scalar, 2 * dim> gravPotTrace(0);
99 const auto firstElement = *this->
gridView_.template begin<0>();
100 FluidState fluidState;
101 fluidState.setPressure(wPhaseIdx, problem.referencePressure(firstElement));
102 fluidState.setPressure(nPhaseIdx, problem.referencePressure(firstElement));
103 fluidState.setTemperature(problem.temperature(firstElement));
104 fluidState.setSaturation(wPhaseIdx, 1.);
105 fluidState.setSaturation(nPhaseIdx, 0.);
111 for (
int i = 0; i < problem.gridView().size(0); i++)
113 problem.variables().cellData(i).fluxData().resetVelocity();
117 for (
const auto& element : elements(this->
gridView_))
119 int eIdxGlobal = problem.variables().index(element);
121 CellData& cellData = problem.variables().cellData(eIdxGlobal);
122 FieldVector globalPos = element.geometry().center();
125 for (
const auto& intersection : intersections(problem.gridView(), element))
127 int indexInInside = intersection.indexInInside();
129 int fIdxGlobal = this->
faceMapper_.subIndex(element, indexInInside, 1);
131 pressTrace[indexInInside] = u[fIdxGlobal];
133 gravPotTrace[indexInInside] = (problem.bBoxMax() - intersection.geometry().center()) * problem.gravity() * densityDiff;
136 switch (pressureType)
140 Scalar potW = loc.constructPressure(element, pressTrace);
141 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
142 Scalar potNw = potW + gravPot;
144 cellData.setPotential(wPhaseIdx, potW);
145 cellData.setPotential(nPhaseIdx, potNw);
147 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, wPhaseIdx);
149 cellData.setPressure(wPhaseIdx, potW - gravPot);
151 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, nPhaseIdx);
153 cellData.setPressure(nPhaseIdx, potNw - gravPot);
159 Scalar potNw = loc.constructPressure(element, pressTrace);
160 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
161 Scalar potW = potNw - gravPot;
163 cellData.setPotential(nPhaseIdx, potNw);
164 cellData.setPotential(wPhaseIdx, potW);
166 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, wPhaseIdx);
168 cellData.setPressure(wPhaseIdx, potW - gravPot);
170 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, nPhaseIdx);
172 cellData.setPressure(nPhaseIdx, potNw - gravPot);
181 switch (pressureType)
185 loc.constructVelocity(element, velocityW, pressTrace, cellData.potential(wPhaseIdx));
186 pressTrace += gravPotTrace;
187 loc.constructVelocity(element, velocityNw, pressTrace, cellData.potential(nPhaseIdx));
193 loc.constructVelocity(element, velocityW, pressTrace, cellData.potential(nPhaseIdx));
194 pressTrace -= gravPotTrace;
195 loc.constructVelocity(element, velocityNw, pressTrace, cellData.potential(wPhaseIdx));
201 for (
const auto& intersection : intersections(problem.gridView(), element))
203 int idxInInside = intersection.indexInInside();
205 cellData.fluxData().setUpwindPotential(wPhaseIdx, idxInInside, velocityW[idxInInside]);
206 cellData.fluxData().setUpwindPotential(nPhaseIdx, idxInInside, velocityNw[idxInInside]);
208 Scalar mobilityW = 0;
209 Scalar mobilityNw = 0;
211 if (intersection.neighbor())
213 int neighborIdx = problem.variables().index(intersection.outside());
215 CellData& cellDataNeighbor = problem.variables().cellData(neighborIdx);
218 (velocityW[idxInInside] >= 0.) ? cellData.mobility(wPhaseIdx) :
219 cellDataNeighbor.mobility(wPhaseIdx);
221 (velocityNw[idxInInside] >= 0.) ? cellData.mobility(nPhaseIdx) :
222 cellDataNeighbor.mobility(nPhaseIdx);
224 if (velocityW[idxInInside] >= 0.)
226 FieldVector velocity(intersection.centerUnitOuterNormal());
227 velocity *= mobilityW/(mobilityW+mobilityNw) * velocityW[idxInInside];
228 cellData.fluxData().addVelocity(wPhaseIdx, idxInInside, velocity);
229 cellDataNeighbor.fluxData().addVelocity(wPhaseIdx, intersection.indexInOutside(), velocity);
231 if (velocityNw[idxInInside] >= 0.)
233 FieldVector velocity(intersection.centerUnitOuterNormal());
234 velocity *= mobilityNw/(mobilityW+mobilityNw) * velocityNw[idxInInside];
235 cellData.fluxData().addVelocity(nPhaseIdx, idxInInside, velocity);
236 cellDataNeighbor.fluxData().addVelocity(nPhaseIdx, intersection.indexInOutside(), velocity);
239 cellData.fluxData().setVelocityMarker(idxInInside);
243 BoundaryTypes bctype;
244 problem.boundaryTypes(bctype, intersection);
245 if (bctype.isDirichlet(satEqIdx))
247 PrimaryVariables boundValues(0.0);
248 problem.dirichlet(boundValues, intersection);
250 if (velocityW[idxInInside] >= 0.)
252 mobilityW = cellData.mobility(wPhaseIdx);
256 mobilityW = MaterialLaw::krw(problem.spatialParams().materialLawParams(element),
257 boundValues[saturationIdx]) / viscosityW;
260 if (velocityNw[idxInInside] >= 0.)
262 mobilityNw = cellData.mobility(nPhaseIdx);
266 mobilityNw = MaterialLaw::krn(problem.spatialParams().materialLawParams(element),
267 boundValues[saturationIdx]) / viscosityNw;
272 mobilityW = cellData.mobility(wPhaseIdx);
273 mobilityNw = cellData.mobility(nPhaseIdx);
276 FieldVector velocity(intersection.centerUnitOuterNormal());
277 velocity *= mobilityW / (mobilityW + mobilityNw) * velocityW[idxInInside];
278 cellData.fluxData().setVelocity(wPhaseIdx, idxInInside, velocity);
281 velocity = intersection.centerUnitOuterNormal();
282 velocity *= mobilityNw / (mobilityW + mobilityNw) * velocityNw[idxInInside];
283 cellData.fluxData().setVelocity(nPhaseIdx, idxInInside, velocity);
284 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:45
MimeticOperatorAssemblerTwoP(const GridView &gridView)
Definition: operator.hh:85
void calculatePressure(LocalStiffness &loc, Vector &u, Problem &problem)
Definition: operator.hh:92
Specifies the properties for immiscible 2p diffusion/pressure models.
Defines the basic properties required for a mimetic method.