24#ifndef DUMUX_MIMETIC2PADAPTIVE_HH
25#define DUMUX_MIMETIC2PADAPTIVE_HH
34#include<dune/common/exceptions.hh>
35#include<dune/grid/common/grid.hh>
36#include<dune/geometry/type.hh>
37#include<dune/geometry/quadraturerules.hh>
44#include <dune/common/dynvector.hh>
45#include <dune/common/dynmatrix.hh>
55template<
class TypeTag>
67 dim = GridView::dimension
71 wPhaseIdx = Indices::wPhaseIdx,
72 nPhaseIdx = Indices::nPhaseIdx,
73 pressureIdx = Indices::pressureIdx,
74 saturationIdx = Indices::saturationIdx,
75 pressureEqIdx = Indices::pressureEqIdx,
76 satEqIdx = Indices::satEqIdx,
77 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
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,
89 pressureType = getPropValue<TypeTag, Properties::PressureFormulation>(),
91 saturationType = getPropValue<TypeTag, Properties::SaturationFormulation>()
94 using Grid =
typename GridView::Grid;
95 using Element =
typename GridView::Traits::template Codim<0>::Entity;
99 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
120 const GridView& gridView,
const IntersectionMapper& intersectionMapper,
bool procBoundaryAsDirichlet =
true) :
121 problem_(
problem), gridView_(gridView), intersectionMapper_(intersectionMapper), maxError_(0), timeStep_(1)
123 int size = gridView_.size(0);
124 rhs_.resize(
size , 0.);
126 ErrorTermFactor_ = getParam<Scalar>(
"Impet.ErrorTermFactor");
127 ErrorTermLowerBound_ = getParam<Scalar>(
"Impet.ErrorTermLowerBound");
128 ErrorTermUpperBound_ = getParam<Scalar>(
"Impet.ErrorTermUpperBound");
130 density_[wPhaseIdx] = 0.0;
131 density_[nPhaseIdx] = 0.0;
132 viscosity_[wPhaseIdx] =0.0;
133 viscosity_[nPhaseIdx] = 0.0;
138 const auto element = *problem_.gridView().template begin<0>();
139 FluidState fluidState;
140 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
141 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
142 fluidState.setTemperature(problem_.temperature(element));
143 fluidState.setSaturation(wPhaseIdx, 1.);
144 fluidState.setSaturation(nPhaseIdx, 0.);
155 int size = gridView_.size(0);
156 rhs_.resize(
size , 0.);
175 return W_[eIdxGlobal].N();
192 int numFaces = intersectionMapper_.
size(element);
197 for (
int i = 0; i < numFaces; i++)
201 for (
int j = 0; j < numFaces; j++)
205 assembleV(element, numFaces, k);
206 assembleBC(element, k);
214 void assemble(
const Element &cell,
const Dune::BlockVector<VBlockType>& localSolution,
int orderOfShapeFns = 1)
231 int numFaces = intersectionMapper_.
size(element);
234 for (
int i = 0; i < numFaces; i++)
240 assembleBC(element, k);
243 template<
class Vector>
244 void completeRHS(
const Element& element, std::vector<int>& local2Global, Vector& f)
246 int eIdxGlobal = problem_.variables().index(element);
250 Dune::DynamicVector<Scalar> F(numFaces);
254 for (
int i = 0; i < numFaces; i++)
257 f[local2Global[i]][0] += (dInv * F[i] * rhs_[eIdxGlobal]);
263 int eIdxGlobal = problem_.variables().index(element);
267 Scalar volume = element.geometry().volume();
269 PrimaryVariables sourceVec(0.0);
270 problem_.source(sourceVec, element);
271 Scalar qmean = volume * (sourceVec[wPhaseIdx]/density_[wPhaseIdx] + sourceVec[nPhaseIdx]/density_[nPhaseIdx]);
272 qmean += rhs_[eIdxGlobal];
274 Dune::DynamicVector<Scalar> F(numFaces, 0.);
278 return (dInv*(qmean + (F*pressTrace)));
282 Dune::DynamicVector<Scalar>& pressTrace, Scalar press)
284 int eIdxGlobal = problem_.variables().index(element);
288 Dune::DynamicVector<Scalar> faceVol(numFaces);
290 Dune::FieldVector<Scalar, 2*dim> faceVolumeReal(0.0);
293 for (
const auto& intersection : intersections(gridView_, element))
295 faceVol[idx] = intersection.geometry().volume();
296 int indexInInside = intersection.indexInInside();
297 faceVolumeReal[indexInInside] += faceVol[idx];
303 for (
int i = 0; i < numFaces; i++)
304 for (
int j = 0; j < numFaces; j++)
305 vel[i] += W_[eIdxGlobal][i][j]*faceVol[j]*(press - pressTrace[j]);
307 for (
int i = 0; i < numFaces; i++)
309 vel[i] *= faceVol[i]/faceVolumeReal[intersectionMapper_.
maplocal(eIdxGlobal, i)];
315 Dune::DynamicVector<Scalar>& F, Scalar& dInv)
317 int eIdxGlobal = problem_.variables().index(element);
321 Dune::DynamicVector<Scalar> c(numFaces);
322 Dune::DynamicMatrix<Scalar> Pi(numFaces, numFaces);
325 for (
const auto& intersection : intersections(gridView_, element))
327 Scalar faceVol = intersection.geometry().volume();
336 Pi[idx][idx] = faceVol;
341 Dune::DynamicVector<Scalar> Wc(numFaces);
343 dInv = 1.0 / (c * Wc);
351 Dune::DynamicVector<Scalar>& faceVol,
352 Dune::DynamicVector<Scalar>& F,
353 Dune::DynamicMatrix<Scalar>& Pi,
364 void assembleV(
const Element& element,
int numFaces,
int k = 1);
366 void assembleBC(
const Element& element,
int k = 1);
368 Scalar evaluateErrorTerm(CellData& cellData)
373 switch (saturationType)
376 sat = cellData.saturation(wPhaseIdx);
379 sat = cellData.saturation(nPhaseIdx);
382 DUNE_THROW(Dune::NotImplemented,
"Only saturation formulation Sw and Sn are implemented!");
385 Scalar error = (sat > 1.0) ? sat - 1.0 : 0.0;
393 Scalar errorAbs = abs(error);
395 if ((errorAbs * timeStep_ > 1e-6) && (errorAbs > ErrorTermLowerBound_ * maxError_)
396 && (!problem_.timeManager().willBeFinished()))
398 return ErrorTermFactor_ * error;
405 const GridView gridView_;
406 const IntersectionMapper& intersectionMapper_;
407 Dune::DynamicVector<Scalar> rhs_;
408 std::vector<Dune::DynamicMatrix<Scalar> > W_;
411 Scalar ErrorTermFactor_;
412 Scalar ErrorTermLowerBound_;
413 Scalar ErrorTermUpperBound_;
415 Scalar density_[numPhases];
416 Scalar viscosity_[numPhases];
420template<
class TypeTag>
421void MimeticTwoPLocalStiffnessAdaptive<TypeTag>::assembleV(
const Element& element,
int numFaces,
int k)
423 int eIdxGlobal = problem_.variables().index(element);
428 W_[eIdxGlobal].resize(numFaces, numFaces);
429 Dune::DynamicVector<Scalar> faceVol(numFaces);
430 Dune::DynamicMatrix<Scalar> Pi(numFaces, numFaces);
431 Dune::DynamicVector<Scalar> F(numFaces);
434 this->assembleElementMatrices(element, faceVol, F, Pi, dInv, qmean);
437 Dune::DynamicMatrix<Scalar> PiWPiT(W_[eIdxGlobal]);
438 PiWPiT.rightmultiply(Pi);
439 PiWPiT.leftmultiply(Pi);
442 Dune::DynamicMatrix<Scalar> FDinvFT(numFaces, numFaces);
443 for (
int i = 0; i < numFaces; i++)
444 for (
int j = 0; j < numFaces; j++)
445 FDinvFT[i][j] = dInv * F[i] * F[j];
448 for (
int i = 0; i < numFaces; i++)
449 for (
int j = 0; j < numFaces; j++)
450 this->A[i][j] = PiWPiT[i][j] - FDinvFT[i][j];
454 Scalar factor = dInv * qmean;
455 for (
int i = 0; i < numFaces; i++)
456 this->b[i] = F[i] * factor;
466template<
class TypeTag>
468 Dune::DynamicVector<Scalar>& faceVol,
469 Dune::DynamicVector<Scalar>& F,
470 Dune::DynamicMatrix<Scalar>& Pi,
475 const Dune::FieldVector<Scalar, dim>& centerGlobal = element.geometry().center();
477 int eIdxGlobal = problem_.variables().index(element);
479 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
481 int numFaces = intersectionMapper_.size(eIdxGlobal);
484 Scalar volume = element.geometry().volume();
487 Dune::DynamicMatrix<Scalar> R(numFaces, dim, 0.);
488 Dune::DynamicMatrix<Scalar> N(numFaces, dim, 0.);
493 std::vector<Scalar> pcPotFace(numFaces);
494 Scalar pcPot = (problem_.bBoxMax() - element.geometry().center()) * problem_.gravity()
495 * (density_[nPhaseIdx] - density_[wPhaseIdx]);
498 for (
const auto& intersection : intersections(gridView_, element))
502 Dune::FieldVector<Scalar, dim> faceGlobal(intersection.geometry().center());
503 faceVol[idx] = intersection.geometry().volume();
506 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
509 for (
int k = 0; k < dim; k++)
512 R[idx][k] = faceVol[idx] * (faceGlobal[k] - centerGlobal[k]);
513 N[idx][k] = unitOuterNormal[k];
516 pcPotFace[idx] = (problem_.bBoxMax() - faceGlobal) * problem_.gravity()
517 * (density_[nPhaseIdx] - density_[wPhaseIdx]);
525 Scalar norm = R[0][0] * R[0][0];
527 for (
int i = 1; i < numFaces; i++)
528 norm += R[i][0] * R[i][0];
530 for (
int i = 0; i < numFaces; i++)
532 Scalar weight = R[0][1] * R[0][0];
533 for (
int i = 1; i < numFaces; i++)
534 weight += R[i][1] * R[i][0];
535 for (
int i = 0; i < numFaces; i++)
536 R[i][1] -= weight * R[i][0];
537 norm = R[0][1] * R[0][1];
538 for (
int i = 1; i < numFaces; i++)
539 norm += R[i][1] * R[i][1];
541 for (
int i = 0; i < numFaces; i++)
545 Scalar weight1 = R[0][2] * R[0][0];
546 Scalar weight2 = R[0][2] * R[0][1];
547 for (
int i = 1; i < numFaces; i++)
549 weight1 += R[i][2] * R[i][0];
550 weight2 += R[i][2] * R[i][1];
552 for (
int i = 0; i < numFaces; i++)
553 R[i][2] -= weight1 * R[i][0] + weight2 * R[i][1];
554 norm = R[0][2] * R[0][2];
555 for (
int i = 1; i < numFaces; i++)
556 norm += R[i][2] * R[i][2];
559 for (
int i = 0; i < numFaces; i++)
565 Dune::DynamicMatrix<Scalar> D(numFaces, numFaces, 0.);
566 for (
int s = 0; s < numFaces; s++)
568 Dune::DynamicVector<Scalar> es(numFaces, 0.);
570 for (
int k = 0; k < numFaces; k++)
573 for (
int i = 0; i < dim; i++)
575 D[k][s] -= R[s][i] * R[k][i];
582 Dune::FieldMatrix<Scalar, dim, dim>
mobility(problem_.spatialParams().intrinsicPermeability(element));
583 mobility *= (cellData.mobility(wPhaseIdx) + cellData.mobility(nPhaseIdx));
586 for (
int i = 1; i < dim; i++)
588 D *= 2.0 * traceMob / volume;
592 Dune::DynamicMatrix<Scalar> NK(N);
597 for (
int i = 0; i < numFaces; i++)
599 for (
int j = 0; j < numFaces; j++)
601 W_[eIdxGlobal][i][j] = NK[i][0] * N[j][0];
602 for (
int k = 1; k < dim; k++)
603 W_[eIdxGlobal][i][j] += NK[i][k] * N[j][k];
607 W_[eIdxGlobal] /= volume;
620 Dune::DynamicVector<Scalar> c(numFaces);
621 for (
int i = 0; i < numFaces; i++)
626 for (
int i = 0; i < numFaces; i++)
627 Pi[i][i] = faceVol[i];
630 Dune::DynamicVector<Scalar> Wc(numFaces);
631 W_[eIdxGlobal].umv(c, Wc);
632 dInv = 1.0 / (c * Wc);
641 for (
const auto& intersection : intersections(gridView_, element))
646 for (
int j = 0; j < numFaces; j++)
647 flux += W_[eIdxGlobal][idx][j] * faceVol[j] * (pcPot - pcPotFace[j]);
649 if (intersection.neighbor())
652 int eIdxGlobalNeighbor = problem_.variables().index(intersection.outside());
655 switch (pressureType)
659 fracFlow = cellData.fracFlowFunc(nPhaseIdx);
664 fracFlow = -cellData.fracFlowFunc(wPhaseIdx);
668 DUNE_THROW(Dune::NotImplemented,
"Only pressure formulation pw and pn are implemented!");
671 rhs_[eIdxGlobal] -= (faceVol[idx] * fracFlow * flux);
672 rhs_[eIdxGlobalNeighbor] += (faceVol[idx] * fracFlow * flux);
675 else if (intersection.boundary())
677 BoundaryTypes bctype;
678 problem_.boundaryTypes(bctype, intersection);
680 if (bctype.isDirichlet(pressureEqIdx))
682 if (flux > 0. || !bctype.isDirichlet(satEqIdx))
684 switch (pressureType)
688 fracFlow = cellData.fracFlowFunc(nPhaseIdx);
693 fracFlow = -cellData.fracFlowFunc(wPhaseIdx);
697 DUNE_THROW(Dune::NotImplemented,
"Only pressure formulation pw and pn are implemented!");
700 else if (flux < 0. && bctype.isDirichlet(satEqIdx))
702 PrimaryVariables boundValues(0.0);
703 problem_.dirichlet(boundValues, intersection);
708 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
710 const Scalar krw = fluidMatrixInteraction.krw(boundValues[saturationIdx]);
711 const Scalar krn = fluidMatrixInteraction.krn(boundValues[saturationIdx]);
713 switch (pressureType)
717 fracFlow = krn / viscosity_[nPhaseIdx]
718 / (krw / viscosity_[wPhaseIdx] + krn / viscosity_[nPhaseIdx]);
723 fracFlow = -krw / viscosity_[wPhaseIdx]
724 / (krw / viscosity_[wPhaseIdx] + krn / viscosity_[nPhaseIdx]);
728 DUNE_THROW(Dune::NotImplemented,
"Only pressure formulation pw and pn are implemented!");
732 rhs_[eIdxGlobal] -= (faceVol[idx] * fracFlow * flux);
738 PrimaryVariables sourceVec(0.0);
739 problem_.source(sourceVec, element);
740 qmean = volume * (sourceVec[wPhaseIdx]/density_[wPhaseIdx] + sourceVec[nPhaseIdx]/density_[nPhaseIdx]);
742 qmean += evaluateErrorTerm(cellData) * volume;
746template<
class TypeTag>
750 unsigned int faceIndex = 0;
751 for (
const auto& intersection : intersections(gridView_, element))
753 if (intersection.neighbor())
757 else if (intersection.boundary())
759 problem_.boundaryTypes(this->bctype[faceIndex], intersection);
760 PrimaryVariables boundValues(0.0);
762 if (this->bctype[faceIndex].isNeumann(pressureEqIdx))
764 problem_.neumann(boundValues, intersection);
765 Scalar J = (boundValues[wPhaseIdx]/density_[wPhaseIdx] + boundValues[nPhaseIdx]/density_[nPhaseIdx]);
766 this->b[faceIndex] -= J * intersection.geometry().volume();
768 else if (this->bctype[faceIndex].isDirichlet(pressureEqIdx))
770 problem_.dirichlet(boundValues, intersection);
771 if (pressureType == pw)
772 this->b[faceIndex] = boundValues[pressureIdx] + (problem_.bBoxMax()
773 - intersection.geometry().center()) * problem_.gravity() * density_[wPhaseIdx];
774 else if (pressureType == pn)
775 this->b[faceIndex] = boundValues[pressureIdx] + (problem_.bBoxMax()
776 - intersection.geometry().center()) * problem_.gravity() * density_[nPhaseIdx];
778 this->b[faceIndex] = boundValues[pressureIdx];
Definition of boundary condition types, extend if necessary.
defines intersection mappers.
Base class for assembling local stiffness matrices.
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 mobility(int phaseIdx) noexcept
I/O name of mobility for multiphase systems.
Definition: name.hh:101
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
defines an intersection mapper for mapping of global DOFs assigned to faces which also works for adap...
Definition: intersectionmapper.hh:185
unsigned int size() const
Definition: intersectionmapper.hh:292
int maplocal(int elemIdx, int fIdx)
Definition: intersectionmapper.hh:270
Base class for local assemblers.
Definition: localstiffness.hh:60
void setcurrentsize(int s)
Sets the current size of the local stiffness matrix.
Definition: localstiffness.hh:233
Dune::Matrix< MBlockType > A
Definition: localstiffness.hh:250
std::vector< VBlockType > b
Definition: localstiffness.hh:251
std::vector< BoundaryTypes > bctype
Definition: localstiffness.hh:252
compute local stiffness matrix for conforming finite elements for the full 2-phase pressure equation
Definition: mimeticadaptive.hh:57
Scalar constructPressure(const Element &element, Dune::DynamicVector< Scalar > &pressTrace)
Definition: mimeticadaptive.hh:261
void assembleBoundaryCondition(const Element &element, int k=1)
Assembles only boundary conditions for given element.
Definition: mimeticadaptive.hh:229
void initialize()
Definition: mimeticadaptive.hh:136
void computeReconstructionMatrices(const Element &element, const Dune::DynamicMatrix< Scalar > &W, Dune::DynamicVector< Scalar > &F, Scalar &dInv)
Definition: mimeticadaptive.hh:314
void reset()
Definition: mimeticadaptive.hh:160
void setErrorInfo(Scalar maxErr, Scalar dt)
Definition: mimeticadaptive.hh:167
void assemble(const Element &cell, const Dune::BlockVector< VBlockType > &localSolution, int orderOfShapeFns=1)
Definition: mimeticadaptive.hh:214
void adapt()
Definition: mimeticadaptive.hh:153
@ m
Definition: mimeticadaptive.hh:111
void completeRHS(const Element &element, std::vector< int > &local2Global, Vector &f)
Definition: mimeticadaptive.hh:244
void assemble(const Element &element, int k=1)
Assembles local stiffness matrix for given element and order.
Definition: mimeticadaptive.hh:190
void constructVelocity(const Element &element, Dune::DynamicVector< Scalar > &vel, Dune::DynamicVector< Scalar > &pressTrace, Scalar press)
Definition: mimeticadaptive.hh:281
const Problem & problem() const
Definition: mimeticadaptive.hh:358
@ size
Definition: mimeticadaptive.hh:115
int numberOfFaces(int eIdxGlobal)
Definition: mimeticadaptive.hh:173
MimeticTwoPLocalStiffnessAdaptive(Problem &problem, bool levelBoundaryAsDirichlet, const GridView &gridView, const IntersectionMapper &intersectionMapper, bool procBoundaryAsDirichlet=true)
Constructor.
Definition: mimeticadaptive.hh:119
Dune::FieldVector< Scalar, m > VBlockType
Definition: mimeticadaptive.hh:213
void assembleElementMatrices(const Element &element, Dune::DynamicVector< Scalar > &faceVol, Dune::DynamicVector< Scalar > &F, Dune::DynamicMatrix< Scalar > &Pi, Scalar &dInv, Scalar &qmean)
Definition: mimeticadaptive.hh:467
Specifies the properties for immiscible 2p diffusion/pressure models.