24#ifndef DUMUX_MIMETICPRESSURE2PADAPTIVE_HH
25#define DUMUX_MIMETICPRESSURE2PADAPTIVE_HH
27#include <dune/common/exceptions.hh>
78 dim = GridView::dimension, dimWorld = GridView::dimensionworld
82 pw = Indices::pressureW,
83 pn = Indices::pressureNw,
84 pGlobal = Indices::pressureGlobal,
85 Sw = Indices::saturationW,
86 Sn = Indices::saturationNw,
87 vw = Indices::velocityW,
88 vn = Indices::velocityNw,
90 pressureType = getPropValue<TypeTag, Properties::PressureFormulation>(),
92 saturationType = getPropValue<TypeTag, Properties::SaturationFormulation>(),
96 wPhaseIdx = Indices::wPhaseIdx,
97 nPhaseIdx = Indices::nPhaseIdx,
98 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
101 using Element =
typename GridView::Traits::template Codim<0>::Entity;
102 using Grid =
typename GridView::Grid;
104 using Geometry =
typename Element::Geometry;
105 using JacobianTransposed =
typename Geometry::JacobianTransposed ;
107 using DimVector = Dune::FieldVector<Scalar, dim>;
110 using TraceType = Dune::BlockVector<Dune::FieldVector<Scalar, 1> >;
115 using ScalarSolutionType =
typename SolutionTypes::ScalarSolution;
121 void initializeMatrix();
124 void assemble(
bool first)
126 Scalar timeStep = problem_.timeManager().timeStepSize();
127 Scalar maxError = 0.0;
128 int size = problem_.gridView().size(0);
129 for (
int i = 0; i < size; i++)
133 switch (saturationType)
136 sat = problem_.variables().cellData(i).saturation(wPhaseIdx);
139 sat = problem_.variables().cellData(i).saturation(nPhaseIdx);
142 DUNE_THROW(Dune::NotImplemented,
"Only saturation formulation Sw and Sn are implemented!");
146 maxError = max(maxError, (sat - 1.0) / timeStep);
150 maxError = max(maxError, (-sat) / timeStep);
154 lstiff_.setErrorInfo(maxError, timeStep);
155 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 lstiff_.initialize();
233 if (problem_.gridAdapt().wasAdapted())
254 template<
class MultiWriter>
257 int size = problem_.gridView().size(0);
258 ScalarSolutionType *potential = writer.allocateManagedBuffer(size);
260 ScalarSolutionType *pressureSecond = 0;
261 ScalarSolutionType *potentialSecond = 0;
262 Dune::BlockVector < DimVector > *velocityWetting = 0;
263 Dune::BlockVector < DimVector > *velocityNonwetting = 0;
265 if (vtkOutputLevel_ > 0)
267 pressure = writer.allocateManagedBuffer(size);
268 pressureSecond = writer.allocateManagedBuffer(size);
269 potentialSecond = writer.allocateManagedBuffer(size);
270 velocityWetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
271 velocityNonwetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
275 for (
const auto& element : elements(problem_.gridView()))
277 int eIdxGlobal = problem_.variables().index(element);
278 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
280 if (pressureType == pw)
282 (*potential)[eIdxGlobal] = cellData.potential(wPhaseIdx);
285 if (pressureType == pn)
287 (*potential)[eIdxGlobal] = cellData.potential(nPhaseIdx);
290 if (vtkOutputLevel_ > 0)
293 if (pressureType == pw)
295 (*pressure)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
296 (*potentialSecond)[eIdxGlobal] = cellData.potential(nPhaseIdx);
297 (*pressureSecond)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
300 if (pressureType == pn)
302 (*pressure)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
303 (*potentialSecond)[eIdxGlobal] = cellData.potential(wPhaseIdx);
304 (*pressureSecond)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
307 const typename Element::Geometry& geometry = element.geometry();
310 const auto refElement = referenceElement(geometry);
312 const int numberOfFaces=refElement.size(1);
313 std::vector<Scalar> fluxW(numberOfFaces,0);
314 std::vector<Scalar> fluxNw(numberOfFaces,0);
317 for (
const auto& intersection : intersections(problem_.gridView(), element))
319 int isIndex = intersection.indexInInside();
321 fluxW[isIndex] += intersection.geometry().volume()
322 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(wPhaseIdx, isIndex));
323 fluxNw[isIndex] += intersection.geometry().volume()
324 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(nPhaseIdx, isIndex));
329 Dune::FieldVector<Scalar, dim> refVelocity;
331 if (refElement.type().isSimplex()) {
332 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
334 refVelocity[dimIdx] = -fluxW[dim - 1 - dimIdx];
335 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
337 refVelocity[dimIdx] += fluxW[fIdx]/(dim + 1);
342 else if (refElement.type().isCube()){
343 for (
int i = 0; i < dim; i++)
344 refVelocity[i] = 0.5 * (fluxW[2*i + 1] - fluxW[2*i]);
348 DUNE_THROW(Dune::NotImplemented,
"velocity output for prism/pyramid not implemented");
351 const DimVector& localPos = refElement.position(0, 0);
354 const JacobianTransposed jacobianT = geometry.jacobianTransposed(localPos);
357 DimVector elementVelocity(0);
358 jacobianT.umtv(refVelocity, elementVelocity);
359 elementVelocity /= geometry.integrationElement(localPos);
361 (*velocityWetting)[eIdxGlobal] = elementVelocity;
366 if (refElement.type().isSimplex()) {
367 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
369 refVelocity[dimIdx] = -fluxNw[dim - 1 - dimIdx];
370 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
372 refVelocity[dimIdx] += fluxNw[fIdx]/(dim + 1);
377 else if (refElement.type().isCube()){
378 for (
int i = 0; i < dim; i++)
379 refVelocity[i] = 0.5 * (fluxNw[2*i + 1] - fluxNw[2*i]);
383 DUNE_THROW(Dune::NotImplemented,
"velocity output for prism/pyramid not implemented");
388 jacobianT.umtv(refVelocity, elementVelocity);
389 elementVelocity /= geometry.integrationElement(localPos);
391 (*velocityNonwetting)[eIdxGlobal] = elementVelocity;
395 if (pressureType == pw)
397 writer.attachCellData(*potential,
"wetting potential");
400 if (pressureType == pn)
402 writer.attachCellData(*potential,
"nonwetting potential");
405 if (vtkOutputLevel_ > 0)
407 if (pressureType == pw)
409 writer.attachCellData(*
pressure,
"wetting pressure");
410 writer.attachCellData(*pressureSecond,
"nonwetting pressure");
411 writer.attachCellData(*potentialSecond,
"nonwetting potential");
414 if (pressureType == pn)
416 writer.attachCellData(*
pressure,
"nonwetting pressure");
417 writer.attachCellData(*pressureSecond,
"wetting pressure");
418 writer.attachCellData(*potentialSecond,
"wetting potential");
421 writer.attachCellData(*velocityWetting,
"wetting-velocity", dim);
422 writer.attachCellData(*velocityNonwetting,
"nonwetting-velocity", dim);
436 int numFaces = element.subEntities(1);
437 for (
int i=0; i < numFaces; i++)
440 outstream << pressTrace_[isIdxGlobal][0];
452 int numFaces = element.subEntities(1);
453 for (
int i=0; i < numFaces; i++)
456 instream >> pressTrace_[isIdxGlobal][0];
467 A_(problem.gridView()), lstiff_(problem_, false, problem_.gridView(), A_.intersectionMapper())
469 if (pressureType != pw && pressureType != pn)
471 DUNE_THROW(Dune::NotImplemented,
"Pressure type not supported!");
473 if (saturationType != Sw && saturationType != Sn)
475 DUNE_THROW(Dune::NotImplemented,
"Saturation type not supported!");
477 if (getPropValue<TypeTag, Properties::EnableCompressibility>())
479 DUNE_THROW(Dune::NotImplemented,
"Compressibility not supported!");
482 density_[wPhaseIdx] = 0.0;
483 density_[nPhaseIdx] = 0.0;
484 viscosity_[wPhaseIdx] = 0.0;
485 viscosity_[nPhaseIdx] = 0.0;
487 vtkOutputLevel_ = getParam<int>(
"Vtk.OutputLevel");
492 TraceType pressTrace_;
494 OperatorAssembler A_;
497 Scalar density_[numPhases];
498 Scalar viscosity_[numPhases];
504template<
class TypeTag>
505void MimeticPressure2PAdaptive<TypeTag>::solve()
507 using Solver = GetPropType<TypeTag, Properties::LinearSolver>;
509 int verboseLevelSolver = getParam<int>(
"LinearSolver.Verbosity", 0);
511 if (verboseLevelSolver)
512 std::cout <<
"MimeticPressure2PAdaptive: solve for pressure" << std::endl;
514 auto solver = getSolver<Solver>(problem_);
515 solver.solve(*A_, pressTrace_, f_);
524template<
class TypeTag>
528 for (
const auto& element : elements(problem_.gridView()))
530 int eIdxGlobal = problem_.variables().index(element);
532 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
534 const Scalar satW = cellData.saturation(wPhaseIdx);
536 const auto fluidMatrixInteraction = problem_.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
539 const Scalar mobilityW = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
540 const Scalar mobilityNw = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
543 cellData.setMobility(wPhaseIdx, mobilityW);
544 cellData.setMobility(nPhaseIdx, mobilityNw);
547 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
548 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
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 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
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
const IntersectionMapper & intersectionMapper()
Definition: croperatoradaptive.hh:130
void assemble(LocalStiffness &loc, Vector &u, Vector &f)
Assembles global stiffness matrix.
Definition: croperatoradaptive.hh:232
void adapt()
Definition: croperatoradaptive.hh:111
Base class for local assemblers.
Definition: localstiffness.hh:60
Levelwise assembler.
Definition: operatoradaptive.hh:45
void calculatePressure(LocalStiffness &loc, Vector &u, Problem &problem)
Definition: operatoradaptive.hh:90
Mimetic method for the pressure equation.
Definition: mimetic/pressureadaptive.hh:64
void adapt()
Definition: mimetic/pressureadaptive.hh:209
void updateVelocity()
Velocity update.
Definition: mimetic/pressureadaptive.hh:224
void initialize(bool solveTwice=true)
Initializes the model.
Definition: mimetic/pressureadaptive.hh:182
MimeticPressure2PAdaptive(Problem &problem)
Constructs a MimeticPressure2PAdaptive object.
Definition: mimetic/pressureadaptive.hh:465
void update()
updates the model
Definition: mimetic/pressureadaptive.hh:231
void updateMaterialLaws()
Constitutive functions are initialized and stored in the variables object.
Definition: mimetic/pressureadaptive.hh:525
void serializeEntity(std::ostream &outstream, const Element &element)
General method for serialization, output.
Definition: mimetic/pressureadaptive.hh:434
void addOutputVtkFields(MultiWriter &writer)
Write data file.
Definition: mimetic/pressureadaptive.hh:255
void deserializeEntity(std::istream &instream, const Element &element)
General method for deserialization.
Definition: mimetic/pressureadaptive.hh:450
Defines the basic properties required for a mimetic method.
Finite Volume Diffusion Model.