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;
67 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
73 pw = Indices::pressureW,
74 pn = Indices::pressureNw,
75 pressureType = getPropValue<TypeTag, Properties::PressureFormulation>(),
76 wPhaseIdx = Indices::wPhaseIdx,
77 nPhaseIdx = Indices::nPhaseIdx,
78 saturationIdx = Indices::saturationIdx,
79 satEqIdx = Indices::satEqIdx
82 using FieldVector = Dune::FieldVector<Scalar, dimWorld>;
90 template<
class Vector>
93 Dune::DynamicVector<Scalar> velocityW(2*dim);
94 Dune::DynamicVector<Scalar> velocityNw(2*dim);
95 Dune::DynamicVector<Scalar> pressTraceW(2*dim);
96 Dune::DynamicVector<Scalar> pressTraceNw(2*dim);
98 const auto firstElement = *problem.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);
123 velocityW.resize(numFaces);
124 velocityNw.resize(numFaces);
125 pressTraceW.resize(numFaces);
126 pressTraceNw.resize(numFaces);
128 CellData& cellData = problem.variables().cellData(eIdxGlobal);
129 FieldVector globalPos = element.geometry().center();
131 int intersectionIdx = -1;
133 for (
const auto& intersection : intersections(problem.gridView(), element))
139 Scalar pcPotFace = (problem.bBoxMax() - intersection.geometry().center()) * problem.gravity() * densityDiff;
141 switch (pressureType)
145 pressTraceW[intersectionIdx] = u[fIdxGlobal];
146 pressTraceNw[intersectionIdx] = u[fIdxGlobal] + pcPotFace;
151 pressTraceNw[intersectionIdx] = u[fIdxGlobal];
152 pressTraceW[intersectionIdx] = u[fIdxGlobal] - pcPotFace;
159 switch (pressureType)
163 Scalar potW = loc.constructPressure(element, pressTraceW);
164 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
165 Scalar potNw = potW + gravPot;
167 cellData.setPotential(wPhaseIdx, potW);
168 cellData.setPotential(nPhaseIdx, potNw);
170 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, wPhaseIdx);
172 cellData.setPressure(wPhaseIdx, potW - gravPot);
174 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, nPhaseIdx);
176 cellData.setPressure(nPhaseIdx, potNw - gravPot);
182 Scalar potNw = loc.constructPressure(element, pressTraceNw);
183 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
184 Scalar potW = potNw - gravPot;
186 cellData.setPotential(nPhaseIdx, potNw);
187 cellData.setPotential(wPhaseIdx, potW);
189 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, wPhaseIdx);
191 cellData.setPressure(wPhaseIdx, potW - gravPot);
193 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() *
FluidSystem::density(fluidState, nPhaseIdx);
195 cellData.setPressure(nPhaseIdx, potNw - gravPot);
204 loc.constructVelocity(element, velocityW, pressTraceW, cellData.potential(wPhaseIdx));
205 loc.constructVelocity(element, velocityNw, pressTraceNw, cellData.potential(nPhaseIdx));
207 intersectionIdx = -1;
208 for (
const auto& intersection : intersections(problem.gridView(), element))
211 int idxInInside = intersection.indexInInside();
213 cellData.fluxData().addUpwindPotential(wPhaseIdx, idxInInside, velocityW[intersectionIdx]);
214 cellData.fluxData().addUpwindPotential(nPhaseIdx, idxInInside, velocityNw[intersectionIdx]);
216 Scalar mobilityW = 0;
217 Scalar mobilityNw = 0;
219 if (intersection.neighbor())
221 int neighborIdx = problem.variables().index(intersection.outside());
223 CellData& cellDataNeighbor = problem.variables().cellData(neighborIdx);
226 (velocityW[intersectionIdx] >= 0.) ? cellData.mobility(wPhaseIdx) :
227 cellDataNeighbor.mobility(wPhaseIdx);
229 (velocityNw[intersectionIdx] >= 0.) ? cellData.mobility(nPhaseIdx) :
230 cellDataNeighbor.mobility(nPhaseIdx);
232 if (velocityW[intersectionIdx] >= 0.)
234 FieldVector velocity(intersection.centerUnitOuterNormal());
235 velocity *= mobilityW/(mobilityW+mobilityNw) * velocityW[intersectionIdx];
236 cellData.fluxData().addVelocity(wPhaseIdx, idxInInside, velocity);
237 cellDataNeighbor.fluxData().addVelocity(wPhaseIdx, intersection.indexInOutside(), velocity);
239 if (velocityNw[intersectionIdx] >= 0.)
241 FieldVector velocity(intersection.centerUnitOuterNormal());
242 velocity *= mobilityNw/(mobilityW+mobilityNw) * velocityNw[intersectionIdx];
243 cellData.fluxData().addVelocity(nPhaseIdx, idxInInside, velocity);
244 cellDataNeighbor.fluxData().addVelocity(nPhaseIdx, intersection.indexInOutside(), velocity);
247 cellData.fluxData().setVelocityMarker(idxInInside);
251 BoundaryTypes bctype;
252 problem.boundaryTypes(bctype, intersection);
253 if (bctype.isDirichlet(satEqIdx))
255 PrimaryVariables boundValues(0.0);
256 problem.dirichlet(boundValues, intersection);
258 if (velocityW[intersectionIdx] >= 0.)
260 mobilityW = cellData.mobility(wPhaseIdx);
264 mobilityW = MaterialLaw::krw(problem.spatialParams().materialLawParams(element),
265 boundValues[saturationIdx]) / viscosityW;
268 if (velocityNw[intersectionIdx] >= 0.)
270 mobilityNw = cellData.mobility(nPhaseIdx);
274 mobilityNw = MaterialLaw::krn(problem.spatialParams().materialLawParams(element),
275 boundValues[saturationIdx]) / viscosityNw;
280 mobilityW = cellData.mobility(wPhaseIdx);
281 mobilityNw = cellData.mobility(nPhaseIdx);
284 FieldVector velocity(intersection.centerUnitOuterNormal());
285 velocity *= mobilityW/(mobilityW+mobilityNw) * velocityW[intersectionIdx];
286 cellData.fluxData().addVelocity(wPhaseIdx, idxInInside, velocity);
289 velocity = intersection.centerUnitOuterNormal();
290 velocity *= mobilityNw/(mobilityW+mobilityNw) * velocityNw[intersectionIdx];
291 cellData.fluxData().addVelocity(nPhaseIdx, idxInInside, velocity);
292 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
int subIndex(int elemIdx, int fIdx)
Map interface fIdx'th interface of element index to array index.
Definition: intersectionmapper.hh:236
unsigned int size() const
Definition: intersectionmapper.hh:292
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:91
MimeticOperatorAssemblerTwoPAdaptive(const GridView &gridView)
Definition: operatoradaptive.hh:85
Specifies the properties for immiscible 2p diffusion/pressure models.
Defines the basic properties required for a mimetic method.