24#ifndef DUMUX_MIMETICOPERATOR2PADAPTIVE_HH
25#define DUMUX_MIMETICOPERATOR2PADAPTIVE_HH
43template<
class TypeTag>
54 dim=GridView::dimension,
55 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
81 using FieldVector = Dune::FieldVector<Scalar, dimWorld>;
89 template<
class Vector>
92 Dune::DynamicVector<Scalar> velocityW(2*dim);
93 Dune::DynamicVector<Scalar> velocityNw(2*dim);
94 Dune::DynamicVector<Scalar> pressTraceW(2*dim);
95 Dune::DynamicVector<Scalar> pressTraceNw(2*dim);
97 const auto firstElement = *problem.gridView().template begin<0>();
98 FluidState fluidState;
99 fluidState.setPressure(wPhaseIdx, problem.referencePressure(firstElement));
100 fluidState.setPressure(nPhaseIdx, problem.referencePressure(firstElement));
101 fluidState.setTemperature(problem.temperature(firstElement));
102 fluidState.setSaturation(wPhaseIdx, 1.);
103 fluidState.setSaturation(nPhaseIdx, 0.);
109 for (
int i = 0; i < problem.gridView().size(0); i++)
111 problem.variables().cellData(i).fluxData().resetVelocity();
115 for (
const auto& element : elements(this->
gridView_))
117 int eIdxGlobal = problem.variables().index(element);
122 velocityW.resize(numFaces);
123 velocityNw.resize(numFaces);
124 pressTraceW.resize(numFaces);
125 pressTraceNw.resize(numFaces);
127 CellData& cellData = problem.variables().cellData(eIdxGlobal);
128 FieldVector globalPos = element.geometry().center();
130 int intersectionIdx = -1;
132 for (
const auto& intersection : intersections(problem.gridView(), element))
138 Scalar pcPotFace = (problem.bBoxMax() - intersection.geometry().center()) * problem.gravity() * densityDiff;
140 switch (pressureType)
144 pressTraceW[intersectionIdx] = u[fIdxGlobal];
145 pressTraceNw[intersectionIdx] = u[fIdxGlobal] + pcPotFace;
150 pressTraceNw[intersectionIdx] = u[fIdxGlobal];
151 pressTraceW[intersectionIdx] = u[fIdxGlobal] - pcPotFace;
158 switch (pressureType)
162 Scalar potW = loc.constructPressure(element, pressTraceW);
163 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
164 Scalar potNw = potW + gravPot;
166 cellData.setPotential(wPhaseIdx, potW);
167 cellData.setPotential(nPhaseIdx, potNw);
169 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, wPhaseIdx);
171 cellData.setPressure(wPhaseIdx, potW - gravPot);
173 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, nPhaseIdx);
175 cellData.setPressure(nPhaseIdx, potNw - gravPot);
181 Scalar potNw = loc.constructPressure(element, pressTraceNw);
182 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
183 Scalar potW = potNw - gravPot;
185 cellData.setPotential(nPhaseIdx, potNw);
186 cellData.setPotential(wPhaseIdx, potW);
188 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, wPhaseIdx);
190 cellData.setPressure(wPhaseIdx, potW - gravPot);
192 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, nPhaseIdx);
194 cellData.setPressure(nPhaseIdx, potNw - gravPot);
203 loc.constructVelocity(element, velocityW, pressTraceW, cellData.potential(wPhaseIdx));
204 loc.constructVelocity(element, velocityNw, pressTraceNw, cellData.potential(nPhaseIdx));
206 intersectionIdx = -1;
207 for (
const auto& intersection : intersections(problem.gridView(), element))
210 int idxInInside = intersection.indexInInside();
212 cellData.fluxData().addUpwindPotential(wPhaseIdx, idxInInside, velocityW[intersectionIdx]);
213 cellData.fluxData().addUpwindPotential(nPhaseIdx, idxInInside, velocityNw[intersectionIdx]);
215 Scalar mobilityW = 0;
216 Scalar mobilityNw = 0;
218 if (intersection.neighbor())
220 int neighborIdx = problem.variables().index(intersection.outside());
222 CellData& cellDataNeighbor = problem.variables().cellData(neighborIdx);
225 (velocityW[intersectionIdx] >= 0.) ? cellData.mobility(wPhaseIdx) :
226 cellDataNeighbor.mobility(wPhaseIdx);
228 (velocityNw[intersectionIdx] >= 0.) ? cellData.mobility(nPhaseIdx) :
229 cellDataNeighbor.mobility(nPhaseIdx);
231 if (velocityW[intersectionIdx] >= 0.)
233 FieldVector velocity(intersection.centerUnitOuterNormal());
234 velocity *= mobilityW/(mobilityW+mobilityNw) * velocityW[intersectionIdx];
235 cellData.fluxData().addVelocity(wPhaseIdx, idxInInside, velocity);
236 cellDataNeighbor.fluxData().addVelocity(wPhaseIdx, intersection.indexInOutside(), velocity);
238 if (velocityNw[intersectionIdx] >= 0.)
240 FieldVector velocity(intersection.centerUnitOuterNormal());
241 velocity *= mobilityNw/(mobilityW+mobilityNw) * velocityNw[intersectionIdx];
242 cellData.fluxData().addVelocity(nPhaseIdx, idxInInside, velocity);
243 cellDataNeighbor.fluxData().addVelocity(nPhaseIdx, intersection.indexInOutside(), velocity);
246 cellData.fluxData().setVelocityMarker(idxInInside);
250 BoundaryTypes bctype;
251 problem.boundaryTypes(bctype, intersection);
252 if (bctype.isDirichlet(satEqIdx))
254 PrimaryVariables boundValues(0.0);
255 problem.dirichlet(boundValues, intersection);
257 const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
258 if (velocityW[intersectionIdx] >= 0.)
260 mobilityW = cellData.mobility(wPhaseIdx);
264 mobilityW = fluidMatrixInteraction.krw(boundValues[saturationIdx]) / viscosityW;
267 if (velocityNw[intersectionIdx] >= 0.)
269 mobilityNw = cellData.mobility(nPhaseIdx);
273 mobilityNw = fluidMatrixInteraction.krn(boundValues[saturationIdx]) / viscosityNw;
278 mobilityW = cellData.mobility(wPhaseIdx);
279 mobilityNw = cellData.mobility(nPhaseIdx);
282 FieldVector velocity(intersection.centerUnitOuterNormal());
283 velocity *= mobilityW/(mobilityW+mobilityNw) * velocityW[intersectionIdx];
284 cellData.fluxData().addVelocity(wPhaseIdx, idxInInside, velocity);
287 velocity = intersection.centerUnitOuterNormal();
288 velocity *= mobilityNw/(mobilityW+mobilityNw) * velocityNw[intersectionIdx];
289 cellData.fluxData().addVelocity(nPhaseIdx, idxInInside, velocity);
290 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
int subIndex(int elemIdx, int fIdx)
Map interface fIdx'th interface of element index to array index.
Definition: intersectionmapper.hh:275
unsigned int size() const
Definition: intersectionmapper.hh:331
Extends CROperatorBase by a generic methods to assemble global stiffness matrix from local stiffness ...
Definition: croperatoradaptive.hh:69
const GridView gridView_
Definition: croperatoradaptive.hh:163
IntersectionMapper intersectionMapper_
Definition: croperatoradaptive.hh:165
Levelwise assembler.
Definition: operatoradaptive.hh:45
void calculatePressure(LocalStiffness &loc, Vector &u, Problem &problem)
Definition: operatoradaptive.hh:90
MimeticOperatorAssemblerTwoPAdaptive(const GridView &gridView)
Definition: operatoradaptive.hh:84
Specifies the properties for immiscible 2p diffusion/pressure models.
Defines the basic properties required for a mimetic method.