24#ifndef DUMUX_MIMETICPRESSURE2P_HH
25#define DUMUX_MIMETICPRESSURE2P_HH
27#include <dune/common/exceptions.hh>
80 dim = GridView::dimension, dimWorld = GridView::dimensionworld
84 pw = Indices::pressureW,
85 pn = Indices::pressureNw,
86 pGlobal = Indices::pressureGlobal,
87 sw = Indices::saturationW,
88 sn = Indices::saturationNw,
89 vw = Indices::velocityW,
90 vn = Indices::velocityNw,
92 pressureType = getPropValue<TypeTag, Properties::PressureFormulation>(),
94 saturationType = getPropValue<TypeTag, Properties::SaturationFormulation>(),
98 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
99 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
102 using Element =
typename GridView::Traits::template Codim<0>::Entity;
103 using Grid =
typename GridView::Grid;
105 using Geometry =
typename Element::Geometry;
106 using JacobianTransposed =
typename Geometry::JacobianTransposed ;
109 using TraceType = Dune::BlockVector<Dune::FieldVector<Scalar, 1> >;
114 using ScalarSolutionType =
typename SolutionTypes::ScalarSolution;
119 using DimVector = Dune::FieldVector<Scalar, dim>;
122 void initializeMatrix();
125 void assemble(
bool first)
127 Scalar timeStep = problem_.timeManager().timeStepSize();
128 Scalar maxError = 0.0;
129 int size = problem_.gridView().size(0);
130 for (
int i = 0; i < size; i++)
134 switch (saturationType)
137 sat = problem_.variables().cellData(i).saturation(wPhaseIdx);
140 sat = problem_.variables().cellData(i).saturation(nPhaseIdx);
143 DUNE_THROW(Dune::NotImplemented,
"Only saturation formulation sw and sn are implemented!");
147 maxError = max(maxError, (sat - 1.0) / timeStep);
151 maxError = max(maxError, (-sat) / timeStep);
155 lstiff_.setErrorInfo(maxError, timeStep);
156 A_.
assemble(lstiff_, pressTrace_, f_);
184 const auto element = *problem_.gridView().template begin<0>();
185 FluidState fluidState;
186 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
187 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
188 fluidState.setTemperature(problem_.temperature(element));
189 fluidState.setSaturation(wPhaseIdx, 1.);
190 fluidState.setSaturation(nPhaseIdx, 0.);
198 pressTrace_.resize(problem_.gridView().size(1));
199 f_.resize(problem_.gridView().size(1));
200 lstiff_.initialize();
242 template<
class MultiWriter>
245 int size = problem_.gridView().size(0);
246 ScalarSolutionType *potential = writer.allocateManagedBuffer(size);
248 ScalarSolutionType *pressureSecond = 0;
249 ScalarSolutionType *potentialSecond = 0;
250 Dune::BlockVector < DimVector > *velocityWetting = 0;
251 Dune::BlockVector < DimVector > *velocityNonwetting = 0;
253 if (vtkOutputLevel_ > 0)
255 pressure = writer.allocateManagedBuffer(size);
256 pressureSecond = writer.allocateManagedBuffer(size);
257 potentialSecond = writer.allocateManagedBuffer(size);
258 velocityWetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
259 velocityNonwetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
263 for (
const auto& element : elements(problem_.gridView()))
265 int eIdxGlobal = problem_.variables().index(element);
266 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
268 if (pressureType == pw)
270 (*potential)[eIdxGlobal] = cellData.potential(wPhaseIdx);
273 if (pressureType == pn)
275 (*potential)[eIdxGlobal] = cellData.potential(nPhaseIdx);
278 if (vtkOutputLevel_ > 0)
281 if (pressureType == pw)
283 (*pressure)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
284 (*potentialSecond)[eIdxGlobal] = cellData.potential(nPhaseIdx);
285 (*pressureSecond)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
288 if (pressureType == pn)
290 (*pressure)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
291 (*potentialSecond)[eIdxGlobal] = cellData.potential(wPhaseIdx);
292 (*pressureSecond)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
295 const typename Element::Geometry& geometry = element.geometry();
298 const auto refElement = referenceElement(geometry);
300 const int numberOfFaces=refElement.size(1);
301 std::vector<Scalar> fluxW(numberOfFaces,0);
302 std::vector<Scalar> fluxNw(numberOfFaces,0);
305 for (
const auto& intersection : intersections(problem_.gridView(), element))
307 int isIndex = intersection.indexInInside();
309 fluxW[isIndex] += intersection.geometry().volume()
310 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(wPhaseIdx, isIndex));
311 fluxNw[isIndex] += intersection.geometry().volume()
312 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(nPhaseIdx, isIndex));
317 Dune::FieldVector<Scalar, dim> refVelocity;
319 if (refElement.type().isSimplex()) {
320 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
322 refVelocity[dimIdx] = -fluxW[dim - 1 - dimIdx];
323 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
325 refVelocity[dimIdx] += fluxW[fIdx]/(dim + 1);
330 else if (refElement.type().isCube()){
331 for (
int i = 0; i < dim; i++)
332 refVelocity[i] = 0.5 * (fluxW[2*i + 1] - fluxW[2*i]);
336 DUNE_THROW(Dune::NotImplemented,
"velocity output for prism/pyramid not implemented");
339 const DimVector& localPos = refElement.position(0, 0);
342 const JacobianTransposed jacobianT = geometry.jacobianTransposed(localPos);
345 DimVector elementVelocity(0);
346 jacobianT.umtv(refVelocity, elementVelocity);
347 elementVelocity /= geometry.integrationElement(localPos);
349 (*velocityWetting)[eIdxGlobal] = elementVelocity;
354 if (refElement.type().isSimplex()) {
355 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
357 refVelocity[dimIdx] = -fluxNw[dim - 1 - dimIdx];
358 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
360 refVelocity[dimIdx] += fluxNw[fIdx]/(dim + 1);
365 else if (refElement.type().isCube()){
366 for (
int i = 0; i < dim; i++)
367 refVelocity[i] = 0.5 * (fluxNw[2*i + 1] - fluxNw[2*i]);
371 DUNE_THROW(Dune::NotImplemented,
"velocity output for prism/pyramid not implemented");
376 jacobianT.umtv(refVelocity, elementVelocity);
377 elementVelocity /= geometry.integrationElement(localPos);
379 (*velocityNonwetting)[eIdxGlobal] = elementVelocity;
383 if (pressureType == pw)
385 writer.attachCellData(*potential,
"wetting potential");
388 if (pressureType == pn)
390 writer.attachCellData(*potential,
"nonwetting potential");
393 if (vtkOutputLevel_ > 0)
395 if (pressureType == pw)
397 writer.attachCellData(*
pressure,
"wetting pressure");
398 writer.attachCellData(*pressureSecond,
"nonwetting pressure");
399 writer.attachCellData(*potentialSecond,
"nonwetting potential");
402 if (pressureType == pn)
404 writer.attachCellData(*
pressure,
"nonwetting pressure");
405 writer.attachCellData(*pressureSecond,
"wetting pressure");
406 writer.attachCellData(*potentialSecond,
"wetting potential");
409 writer.attachCellData(*velocityWetting,
"wetting-velocity", dim);
410 writer.attachCellData(*velocityNonwetting,
"nonwetting-velocity", dim);
424 int numFaces = element.subEntities(1);
425 for (
int i=0; i < numFaces; i++)
427 int fIdxGlobal = A_.
faceMapper().subIndex(element, i, 1);
428 outstream << pressTrace_[fIdxGlobal][0];
440 int numFaces = element.subEntities(1);
441 for (
int i=0; i < numFaces; i++)
443 int fIdxGlobal = A_.
faceMapper().subIndex(element, i, 1);
444 instream >> pressTrace_[fIdxGlobal][0];
455 A_(problem.gridView()), lstiff_(problem_, false, problem_.gridView())
457 if (pressureType != pw && pressureType != pn)
459 DUNE_THROW(Dune::NotImplemented,
"Pressure type not supported!");
461 if (saturationType != sw && saturationType != sn)
463 DUNE_THROW(Dune::NotImplemented,
"Saturation type not supported!");
465 if (getPropValue<TypeTag, Properties::EnableCompressibility>())
467 DUNE_THROW(Dune::NotImplemented,
"Compressibility not supported!");
470 density_[wPhaseIdx] = 0.0;
471 density_[nPhaseIdx] = 0.0;
472 viscosity_[wPhaseIdx] = 0.0;
473 viscosity_[nPhaseIdx] = 0.0;
475 vtkOutputLevel_ = getParam<int>(
"Vtk.OutputLevel");
480 TraceType pressTrace_;
482 OperatorAssembler A_;
485 Scalar density_[numPhases];
486 Scalar viscosity_[numPhases];
492template<
class TypeTag>
493void MimeticPressure2P<TypeTag>::solve()
495 using Solver = GetPropType<TypeTag, Properties::LinearSolver>;
497 auto verboseLevelSolver = getParam<int>(
"LinearSolver.Verbosity", 0);
499 if (verboseLevelSolver)
500 std::cout <<
"MimeticPressure2P: solve for pressure" << std::endl;
502 auto solver = getSolver<Solver>(problem_);
503 solver.solve(*A_, pressTrace_, f_);
509template<
class TypeTag>
513 for (
const auto& element : elements(problem_.gridView()))
515 int eIdxGlobal = problem_.variables().index(element);
517 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
519 const Scalar satW = cellData.saturation(wPhaseIdx);
524 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
527 const Scalar mobilityW = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
528 const Scalar mobilityNw = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
531 cellData.setMobility(wPhaseIdx, mobilityW);
532 cellData.setMobility(nPhaseIdx, mobilityNw);
535 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
536 cellData.setFracFlowFunc(nPhaseIdx, mobilityNw / (mobilityW + mobilityNw));
Local stiffness matrix for the diffusion equation discretized by mimetic FD.
An assembler for the Jacobian matrix based on mimetic FD.
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 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:47
void calculatePressure(LocalStiffness &loc, Vector &u, Problem &problem)
Definition: operator.hh:93
Mimetic method for the pressure equation.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:66
void updateVelocity()
Velocity update.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:218
void initialize(bool solveTwice=true)
Initializes the model.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:182
void addOutputVtkFields(MultiWriter &writer)
Write data file.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:243
void deserializeEntity(std::istream &instream, const Element &element)
General method for deserialization.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:438
void serializeEntity(std::ostream &outstream, const Element &element)
General method for serialization, output.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:422
void update()
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:225
void updateMaterialLaws()
Constitutive functions are initialized and stored in the variables object.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:510
MimeticPressure2P(Problem &problem)
Constructs a MimeticPressure2P object.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:453
Defines the basic properties required for a mimetic method.
Finite Volume Diffusion Model.