24#ifndef DUMUX_MIMETICPRESSURE2P_HH
25#define DUMUX_MIMETICPRESSURE2P_HH
67 using SpatialParams =
typename GET_PROP_TYPE(TypeTag, SpatialParams);
68 using MaterialLaw =
typename SpatialParams::MaterialLaw;
70 using Indices =
typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
72 using FluidSystem =
typename GET_PROP_TYPE(TypeTag, FluidSystem);
73 using FluidState =
typename GET_PROP_TYPE(TypeTag, FluidState);
77 dim = GridView::dimension, dimWorld = GridView::dimensionworld
81 pw = Indices::pressureW,
82 pn = Indices::pressureNw,
83 pGlobal = Indices::pressureGlobal,
84 sw = Indices::saturationW,
85 sn = Indices::saturationNw,
86 vw = Indices::velocityW,
87 vn = Indices::velocityNw,
95 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
99 using Element =
typename GridView::Traits::template Codim<0>::Entity;
100 using Grid =
typename GridView::Grid;
102 using Geometry =
typename Element::Geometry;
103 using JacobianTransposed =
typename Geometry::JacobianTransposed ;
105 using LocalStiffness =
typename GET_PROP_TYPE(TypeTag, LocalStiffness);
106 using TraceType = Dune::BlockVector<Dune::FieldVector<Scalar, 1> >;
110 using SolutionTypes =
typename GET_PROP(TypeTag, SolutionTypes);
111 using ScalarSolutionType =
typename SolutionTypes::ScalarSolution;
116 using DimVector = Dune::FieldVector<Scalar, dim>;
119 void initializeMatrix();
122 void assemble(
bool first)
124 Scalar timeStep = problem_.timeManager().timeStepSize();
125 Scalar maxError = 0.0;
126 int size = problem_.gridView().size(0);
127 for (
int i = 0; i < size; i++)
131 switch (saturationType)
134 sat = problem_.variables().cellData(i).saturation(wPhaseIdx);
137 sat = problem_.variables().cellData(i).saturation(nPhaseIdx);
142 maxError = max(maxError, (sat - 1.0) / timeStep);
146 maxError = max(maxError, (-sat) / timeStep);
150 lstiff_.setErrorInfo(maxError, timeStep);
151 A_.
assemble(lstiff_, pressTrace_, f_);
179 const auto element = *problem_.gridView().template begin<0>();
180 FluidState fluidState;
181 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
182 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
183 fluidState.setTemperature(problem_.temperature(element));
184 fluidState.setSaturation(wPhaseIdx, 1.);
185 fluidState.setSaturation(nPhaseIdx, 0.);
193 pressTrace_.resize(problem_.gridView().size(1));
194 f_.resize(problem_.gridView().size(1));
195 lstiff_.initialize();
237 template<
class MultiWriter>
240 int size = problem_.gridView().size(0);
241 ScalarSolutionType *potential = writer.allocateManagedBuffer(size);
243 ScalarSolutionType *pressureSecond = 0;
244 ScalarSolutionType *potentialSecond = 0;
245 Dune::BlockVector < DimVector > *velocityWetting = 0;
246 Dune::BlockVector < DimVector > *velocityNonwetting = 0;
248 if (vtkOutputLevel_ > 0)
250 pressure = writer.allocateManagedBuffer(size);
251 pressureSecond = writer.allocateManagedBuffer(size);
252 potentialSecond = writer.allocateManagedBuffer(size);
253 velocityWetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
254 velocityNonwetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
258 for (
const auto& element : elements(problem_.gridView()))
260 int eIdxGlobal = problem_.variables().index(element);
261 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
263 if (pressureType == pw)
265 (*potential)[eIdxGlobal] = cellData.potential(wPhaseIdx);
268 if (pressureType == pn)
270 (*potential)[eIdxGlobal] = cellData.potential(nPhaseIdx);
273 if (vtkOutputLevel_ > 0)
276 if (pressureType == pw)
278 (*pressure)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
279 (*potentialSecond)[eIdxGlobal] = cellData.potential(nPhaseIdx);
280 (*pressureSecond)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
283 if (pressureType == pn)
285 (*pressure)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
286 (*potentialSecond)[eIdxGlobal] = cellData.potential(wPhaseIdx);
287 (*pressureSecond)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
290 const typename Element::Geometry& geometry = element.geometry();
293 using ReferenceElements = Dune::ReferenceElements<Scalar, dim>;
294 const auto refElement = ReferenceElements::general(geometry.type());
296 const int numberOfFaces=refElement.size(1);
297 std::vector<Scalar> fluxW(numberOfFaces,0);
298 std::vector<Scalar> fluxNw(numberOfFaces,0);
301 for (
const auto& intersection :
intersections(problem_.gridView(), element))
303 int isIndex = intersection.indexInInside();
305 fluxW[isIndex] += intersection.geometry().volume()
306 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(wPhaseIdx, isIndex));
307 fluxNw[isIndex] += intersection.geometry().volume()
308 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(nPhaseIdx, isIndex));
313 Dune::FieldVector<Scalar, dim> refVelocity;
315 if (refElement.type().isSimplex()) {
316 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
318 refVelocity[dimIdx] = -fluxW[dim - 1 - dimIdx];
319 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
321 refVelocity[dimIdx] += fluxW[fIdx]/(dim + 1);
326 else if (refElement.type().isCube()){
327 for (
int i = 0; i < dim; i++)
328 refVelocity[i] = 0.5 * (fluxW[2*i + 1] - fluxW[2*i]);
332 DUNE_THROW(Dune::NotImplemented,
"velocity output for prism/pyramid not implemented");
335 const DimVector& localPos = refElement.position(0, 0);
338 const JacobianTransposed jacobianT = geometry.jacobianTransposed(localPos);
341 DimVector elementVelocity(0);
342 jacobianT.umtv(refVelocity, elementVelocity);
343 elementVelocity /= geometry.integrationElement(localPos);
345 (*velocityWetting)[eIdxGlobal] = elementVelocity;
350 if (refElement.type().isSimplex()) {
351 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
353 refVelocity[dimIdx] = -fluxNw[dim - 1 - dimIdx];
354 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
356 refVelocity[dimIdx] += fluxNw[fIdx]/(dim + 1);
361 else if (refElement.type().isCube()){
362 for (
int i = 0; i < dim; i++)
363 refVelocity[i] = 0.5 * (fluxNw[2*i + 1] - fluxNw[2*i]);
367 DUNE_THROW(Dune::NotImplemented,
"velocity output for prism/pyramid not implemented");
372 jacobianT.umtv(refVelocity, elementVelocity);
373 elementVelocity /= geometry.integrationElement(localPos);
375 (*velocityNonwetting)[eIdxGlobal] = elementVelocity;
379 if (pressureType == pw)
381 writer.attachCellData(*potential,
"wetting potential");
384 if (pressureType == pn)
386 writer.attachCellData(*potential,
"nonwetting potential");
389 if (vtkOutputLevel_ > 0)
391 if (pressureType == pw)
393 writer.attachCellData(*
pressure,
"wetting pressure");
394 writer.attachCellData(*pressureSecond,
"nonwetting pressure");
395 writer.attachCellData(*potentialSecond,
"nonwetting potential");
398 if (pressureType == pn)
400 writer.attachCellData(*
pressure,
"nonwetting pressure");
401 writer.attachCellData(*pressureSecond,
"wetting pressure");
402 writer.attachCellData(*potentialSecond,
"wetting potential");
405 writer.attachCellData(*velocityWetting,
"wetting-velocity", dim);
406 writer.attachCellData(*velocityNonwetting,
"non-wetting-velocity", dim);
420 int numFaces = element.subEntities(1);
421 for (
int i=0; i < numFaces; i++)
423 int fIdxGlobal = A_.
faceMapper().subIndex(element, i, 1);
424 outstream << pressTrace_[fIdxGlobal][0];
436 int numFaces = element.subEntities(1);
437 for (
int i=0; i < numFaces; i++)
439 int fIdxGlobal = A_.
faceMapper().subIndex(element, i, 1);
440 instream >> pressTrace_[fIdxGlobal][0];
451 A_(problem.gridView()), lstiff_(problem_, false, problem_.gridView())
453 if (pressureType != pw && pressureType != pn)
455 DUNE_THROW(Dune::NotImplemented,
"Pressure type not supported!");
457 if (saturationType != sw && saturationType != sn)
459 DUNE_THROW(Dune::NotImplemented,
"Saturation type not supported!");
463 DUNE_THROW(Dune::NotImplemented,
"Compressibility not supported!");
466 density_[wPhaseIdx] = 0.0;
467 density_[nPhaseIdx] = 0.0;
468 viscosity_[wPhaseIdx] = 0.0;
469 viscosity_[nPhaseIdx] = 0.0;
471 vtkOutputLevel_ = getParam<int>(
"Vtk.OutputLevel");
476 TraceType pressTrace_;
478 OperatorAssembler A_;
481 Scalar density_[numPhases];
482 Scalar viscosity_[numPhases];
488template<
class TypeTag>
489void MimeticPressure2P<TypeTag>::solve()
491 using Solver =
typename GET_PROP_TYPE(TypeTag, LinearSolver);
493 auto verboseLevelSolver = getParam<int>(
"LinearSolver.Verbosity");
495 if (verboseLevelSolver)
496 std::cout <<
"MimeticPressure2P: solve for pressure" << std::endl;
498 auto solver = getSolver<Solver>(problem_);
499 solver.solve(*A_, pressTrace_, f_);
505template<
class TypeTag>
509 for (
const auto& element : elements(problem_.gridView()))
511 int eIdxGlobal = problem_.variables().index(element);
513 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
515 Scalar satW = cellData.saturation(wPhaseIdx);
518 Scalar mobilityW = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satW)
519 / viscosity_[wPhaseIdx];
520 Scalar mobilityNw = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satW)
521 / viscosity_[nPhaseIdx];
524 cellData.setMobility(wPhaseIdx, mobilityW);
525 cellData.setMobility(nPhaseIdx, mobilityNw);
528 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
529 cellData.setFracFlowFunc(nPhaseIdx, mobilityNw / (mobilityW + mobilityNw));
#define GET_PROP_VALUE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:282
#define GET_PROP(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:281
#define GET_PROP_TYPE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:283
Local stiffness matrix for the diffusion equation discretized by mimetic FD.
An assembler for the Jacobian matrix based on mimetic FD.
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Property tag PressureRHSVector
Type of the right hand side vector given to the linear solver.
Definition: sequential/pressureproperties.hh:60
Property tag EnableCompressibility
Returns whether compressibility is allowed.
Definition: porousmediumflow/2p/sequential/properties.hh:55
Property tag SaturationFormulation
The formulation of the saturation model.
Definition: porousmediumflow/2p/sequential/properties.hh:53
Property tag NumPhases
Number of phases in the system.
Definition: porousmediumflow/sequential/properties.hh:69
Property tag PressureCoefficientMatrix
Gives maximum number of intersections of an element and neighboring elements.
Definition: porousmediumflow/sequential/properties.hh:74
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
const FaceMapper & faceMapper()
Definition: croperator.hh:195
void assemble(LocalStiffness &loc, Vector &u, Vector &f)
Assembles global stiffness matrix.
Definition: croperator.hh:223
void initialize()
Initialize the CR operator assembler.
Definition: croperator.hh:112
Base class for local assemblers.
Definition: localstiffness.hh:60
Levelwise assembler.
Definition: operator.hh:45
void calculatePressure(LocalStiffness &loc, Vector &u, Problem &problem)
Definition: operator.hh:92
Mimetic method for the pressure equation.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:62
void updateVelocity()
Velocity update.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:213
void initialize(bool solveTwice=true)
Initializes the model.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:177
void addOutputVtkFields(MultiWriter &writer)
Write data file.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:238
void deserializeEntity(std::istream &instream, const Element &element)
General method for deserialization.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:434
void serializeEntity(std::ostream &outstream, const Element &element)
General method for serialization, output.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:418
void update()
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:220
void updateMaterialLaws()
Constitutive functions are initialized and stored in the variables object.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:506
MimeticPressure2P(Problem &problem)
Constructs a MimeticPressure2P object.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:449
Defines the basic properties required for a mimetic method.
Finite Volume Diffusion Model.