24#ifndef DUMUX_MIMETICPRESSURE2PADAPTIVE_HH
25#define DUMUX_MIMETICPRESSURE2PADAPTIVE_HH
27#include <dune/common/exceptions.hh>
70 using MaterialLaw =
typename SpatialParams::MaterialLaw;
79 dim = GridView::dimension, dimWorld = GridView::dimensionworld
83 pw = Indices::pressureW,
84 pn = Indices::pressureNw,
85 pGlobal = Indices::pressureGlobal,
86 Sw = Indices::saturationW,
87 Sn = Indices::saturationNw,
88 vw = Indices::velocityW,
89 vn = Indices::velocityNw,
91 pressureType = getPropValue<TypeTag, Properties::PressureFormulation>(),
93 saturationType = getPropValue<TypeTag, Properties::SaturationFormulation>(),
97 wPhaseIdx = Indices::wPhaseIdx,
98 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 ;
108 using DimVector = Dune::FieldVector<Scalar, dim>;
111 using TraceType = Dune::BlockVector<Dune::FieldVector<Scalar, 1> >;
116 using ScalarSolutionType =
typename SolutionTypes::ScalarSolution;
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_);
185 const auto element = *problem_.gridView().template begin<0>();
186 FluidState fluidState;
187 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
188 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
189 fluidState.setTemperature(problem_.temperature(element));
190 fluidState.setSaturation(wPhaseIdx, 1.);
191 fluidState.setSaturation(nPhaseIdx, 0.);
199 lstiff_.initialize();
234 if (problem_.gridAdapt().wasAdapted())
255 template<
class MultiWriter>
258 int size = problem_.gridView().size(0);
259 ScalarSolutionType *potential = writer.allocateManagedBuffer(size);
261 ScalarSolutionType *pressureSecond = 0;
262 ScalarSolutionType *potentialSecond = 0;
263 Dune::BlockVector < DimVector > *velocityWetting = 0;
264 Dune::BlockVector < DimVector > *velocityNonwetting = 0;
266 if (vtkOutputLevel_ > 0)
268 pressure = writer.allocateManagedBuffer(size);
269 pressureSecond = writer.allocateManagedBuffer(size);
270 potentialSecond = writer.allocateManagedBuffer(size);
271 velocityWetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
272 velocityNonwetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
276 for (
const auto& element : elements(problem_.gridView()))
278 int eIdxGlobal = problem_.variables().index(element);
279 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
281 if (pressureType == pw)
283 (*potential)[eIdxGlobal] = cellData.potential(wPhaseIdx);
286 if (pressureType == pn)
288 (*potential)[eIdxGlobal] = cellData.potential(nPhaseIdx);
291 if (vtkOutputLevel_ > 0)
294 if (pressureType == pw)
296 (*pressure)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
297 (*potentialSecond)[eIdxGlobal] = cellData.potential(nPhaseIdx);
298 (*pressureSecond)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
301 if (pressureType == pn)
303 (*pressure)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
304 (*potentialSecond)[eIdxGlobal] = cellData.potential(wPhaseIdx);
305 (*pressureSecond)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
308 const typename Element::Geometry& geometry = element.geometry();
311 using ReferenceElements = Dune::ReferenceElements<Scalar, dim>;
312 const auto refElement = ReferenceElements::general(geometry.type());
314 const int numberOfFaces=refElement.size(1);
315 std::vector<Scalar> fluxW(numberOfFaces,0);
316 std::vector<Scalar> fluxNw(numberOfFaces,0);
319 for (
const auto& intersection : intersections(problem_.gridView(), element))
321 int isIndex = intersection.indexInInside();
323 fluxW[isIndex] += intersection.geometry().volume()
324 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(wPhaseIdx, isIndex));
325 fluxNw[isIndex] += intersection.geometry().volume()
326 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(nPhaseIdx, isIndex));
331 Dune::FieldVector<Scalar, dim> refVelocity;
333 if (refElement.type().isSimplex()) {
334 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
336 refVelocity[dimIdx] = -fluxW[dim - 1 - dimIdx];
337 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
339 refVelocity[dimIdx] += fluxW[fIdx]/(dim + 1);
344 else if (refElement.type().isCube()){
345 for (
int i = 0; i < dim; i++)
346 refVelocity[i] = 0.5 * (fluxW[2*i + 1] - fluxW[2*i]);
350 DUNE_THROW(Dune::NotImplemented,
"velocity output for prism/pyramid not implemented");
353 const DimVector& localPos = refElement.position(0, 0);
356 const JacobianTransposed jacobianT = geometry.jacobianTransposed(localPos);
359 DimVector elementVelocity(0);
360 jacobianT.umtv(refVelocity, elementVelocity);
361 elementVelocity /= geometry.integrationElement(localPos);
363 (*velocityWetting)[eIdxGlobal] = elementVelocity;
368 if (refElement.type().isSimplex()) {
369 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
371 refVelocity[dimIdx] = -fluxNw[dim - 1 - dimIdx];
372 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
374 refVelocity[dimIdx] += fluxNw[fIdx]/(dim + 1);
379 else if (refElement.type().isCube()){
380 for (
int i = 0; i < dim; i++)
381 refVelocity[i] = 0.5 * (fluxNw[2*i + 1] - fluxNw[2*i]);
385 DUNE_THROW(Dune::NotImplemented,
"velocity output for prism/pyramid not implemented");
390 jacobianT.umtv(refVelocity, elementVelocity);
391 elementVelocity /= geometry.integrationElement(localPos);
393 (*velocityNonwetting)[eIdxGlobal] = elementVelocity;
397 if (pressureType == pw)
399 writer.attachCellData(*potential,
"wetting potential");
402 if (pressureType == pn)
404 writer.attachCellData(*potential,
"nonwetting potential");
407 if (vtkOutputLevel_ > 0)
409 if (pressureType == pw)
411 writer.attachCellData(*
pressure,
"wetting pressure");
412 writer.attachCellData(*pressureSecond,
"nonwetting pressure");
413 writer.attachCellData(*potentialSecond,
"nonwetting potential");
416 if (pressureType == pn)
418 writer.attachCellData(*
pressure,
"nonwetting pressure");
419 writer.attachCellData(*pressureSecond,
"wetting pressure");
420 writer.attachCellData(*potentialSecond,
"wetting potential");
423 writer.attachCellData(*velocityWetting,
"wetting-velocity", dim);
424 writer.attachCellData(*velocityNonwetting,
"non-wetting-velocity", dim);
438 int numFaces = element.subEntities(1);
439 for (
int i=0; i < numFaces; i++)
442 outstream << pressTrace_[isIdxGlobal][0];
454 int numFaces = element.subEntities(1);
455 for (
int i=0; i < numFaces; i++)
458 instream >> pressTrace_[isIdxGlobal][0];
469 A_(problem.gridView()), lstiff_(problem_, false, problem_.gridView(), A_.intersectionMapper())
471 if (pressureType != pw && pressureType != pn)
473 DUNE_THROW(Dune::NotImplemented,
"Pressure type not supported!");
475 if (saturationType != Sw && saturationType != Sn)
477 DUNE_THROW(Dune::NotImplemented,
"Saturation type not supported!");
479 if (getPropValue<TypeTag, Properties::EnableCompressibility>())
481 DUNE_THROW(Dune::NotImplemented,
"Compressibility not supported!");
484 density_[wPhaseIdx] = 0.0;
485 density_[nPhaseIdx] = 0.0;
486 viscosity_[wPhaseIdx] = 0.0;
487 viscosity_[nPhaseIdx] = 0.0;
489 vtkOutputLevel_ = getParam<int>(
"Vtk.OutputLevel");
494 TraceType pressTrace_;
496 OperatorAssembler A_;
499 Scalar density_[numPhases];
500 Scalar viscosity_[numPhases];
506template<
class TypeTag>
507void MimeticPressure2PAdaptive<TypeTag>::solve()
509 using Solver = GetPropType<TypeTag, Properties::LinearSolver>;
511 int verboseLevelSolver = getParam<int>(
"LinearSolver.Verbosity", 0);
513 if (verboseLevelSolver)
514 std::cout <<
"MimeticPressure2PAdaptive: solve for pressure" << std::endl;
516 auto solver = getSolver<Solver>(problem_);
517 solver.solve(*A_, pressTrace_, f_);
526template<
class TypeTag>
530 for (
const auto& element : elements(problem_.gridView()))
532 int eIdxGlobal = problem_.variables().index(element);
534 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
536 Scalar satW = cellData.saturation(wPhaseIdx);
539 Scalar mobilityW = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satW)
540 / viscosity_[wPhaseIdx];
541 Scalar mobilityNw = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satW)
542 / viscosity_[nPhaseIdx];
545 cellData.setMobility(wPhaseIdx, mobilityW);
546 cellData.setMobility(nPhaseIdx, mobilityNw);
549 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
550 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
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
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:91
Mimetic method for the pressure equation.
Definition: mimetic/pressureadaptive.hh:64
void adapt()
Definition: mimetic/pressureadaptive.hh:210
void updateVelocity()
Velocity update.
Definition: mimetic/pressureadaptive.hh:225
void initialize(bool solveTwice=true)
Initializes the model.
Definition: mimetic/pressureadaptive.hh:183
MimeticPressure2PAdaptive(Problem &problem)
Constructs a MimeticPressure2PAdaptive object.
Definition: mimetic/pressureadaptive.hh:467
void update()
updates the model
Definition: mimetic/pressureadaptive.hh:232
void updateMaterialLaws()
Constitutive functions are initialized and stored in the variables object.
Definition: mimetic/pressureadaptive.hh:527
void serializeEntity(std::ostream &outstream, const Element &element)
General method for serialization, output.
Definition: mimetic/pressureadaptive.hh:436
void addOutputVtkFields(MultiWriter &writer)
Write data file.
Definition: mimetic/pressureadaptive.hh:256
void deserializeEntity(std::istream &instream, const Element &element)
General method for deserialization.
Definition: mimetic/pressureadaptive.hh:452
Defines the basic properties required for a mimetic method.
Finite Volume Diffusion Model.