24#ifndef DUMUX_MIMETIC2P_HH
25#define DUMUX_MIMETIC2P_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>
43#include <dune/common/dynvector.hh>
53template<
class TypeTag>
65 dim = GridView::dimension
69 wPhaseIdx = Indices::wPhaseIdx,
70 nPhaseIdx = Indices::nPhaseIdx,
71 pressureIdx = Indices::pressureIdx,
72 saturationIdx = Indices::saturationIdx,
73 pressureEqIdx = Indices::pressureEqIdx,
74 satEqIdx = Indices::satEqIdx,
75 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
79 pw = Indices::pressureW,
80 pn = Indices::pressureNw,
81 pGlobal = Indices::pressureGlobal,
82 Sw = Indices::saturationW,
83 Sn = Indices::saturationNw,
84 vw = Indices::velocityW,
85 vn = Indices::velocityNw,
87 pressureType = getPropValue<TypeTag, Properties::PressureFormulation>(),
89 saturationType = getPropValue<TypeTag, Properties::SaturationFormulation>()
92 using Grid =
typename GridView::Grid;
93 using Element =
typename GridView::Traits::template Codim<0>::Entity;
97 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
116 const GridView& gridView,
bool procBoundaryAsDirichlet =
true) :
117 problem_(
problem), gridView_(gridView), maxError_(0), timeStep_(1)
119 ErrorTermFactor_ = getParam<Scalar>(
"Impet.ErrorTermFactor");
120 ErrorTermLowerBound_ = getParam<Scalar>(
"Impet.ErrorTermLowerBound");
121 ErrorTermUpperBound_ = getParam<Scalar>(
"Impet.ErrorTermUpperBound");
123 density_[wPhaseIdx] = 0.0;
124 density_[nPhaseIdx] = 0.0;
125 viscosity_[wPhaseIdx] =0.0;
126 viscosity_[nPhaseIdx] = 0.0;
131 const auto element = *problem_.gridView().template begin<0>();
132 FluidState fluidState;
133 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
134 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
135 fluidState.setTemperature(problem_.temperature(element));
136 fluidState.setSaturation(wPhaseIdx, 1.);
137 fluidState.setSaturation(nPhaseIdx, 0.);
143 int size = gridView_.size(0);
144 rhs_.resize(
size , 0.);
176 unsigned int numFaces = element.subEntities(1);
180 for (
unsigned int i = 0; i < numFaces; i++)
184 for (
unsigned int j = 0; j < numFaces; j++)
188 assembleV(element, k);
189 assembleBC(element, k);
197 void assemble(
const Element &cell,
const Dune::BlockVector<VBlockType>& localSolution,
int orderOfShapeFns = 1)
215 unsigned int numFaces = element.subEntities(1);
219 for (
unsigned int i = 0; i < numFaces; i++)
225 assembleBC(element, k);
229 template<
class Vector>
230 void completeRHS(
const Element& element, Dune::FieldVector<int, 2*dim>& local2Global, Vector& f)
232 int eIdxGlobal = problem_.variables().index(element);
233 unsigned int numFaces = element.subEntities(1);
235 Dune::FieldVector<Scalar, 2 * dim> F(0.);
239 for (
unsigned int i = 0; i < numFaces; i++)
241 if (!this->
bc(i).isDirichlet(pressureEqIdx))
242 f[local2Global[i]][0] += (dInv * F[i] * rhs_[eIdxGlobal]);
249 int eIdxGlobal = problem_.variables().index(element);
250 Scalar volume = element.geometry().volume();
252 PrimaryVariables sourceVec(0.0);
253 problem_.source(sourceVec, element);
254 Scalar qmean = volume * (sourceVec[wPhaseIdx]/density_[wPhaseIdx] + sourceVec[nPhaseIdx]/density_[nPhaseIdx]);
255 qmean += rhs_[eIdxGlobal];
257 Dune::FieldVector<Scalar, 2 * dim> F(0.);
261 return (dInv*(qmean + (F*pressTrace)));
266 Dune::FieldVector<Scalar,2*dim>& pressTrace, Scalar press)
268 int eIdxGlobal = problem_.variables().index(element);
270 Dune::FieldVector<Scalar, 2 * dim> faceVol(0);
271 for (
const auto& intersection : intersections(gridView_, element))
273 faceVol[intersection.indexInInside()] = intersection.geometry().volume();
277 for (
int i = 0; i < 2*dim; i++)
278 for (
int j = 0; j < 2*dim; j++)
279 vel[i] += W_[eIdxGlobal][i][j]*faceVol[j]*(press - pressTrace[j]);
284 Dune::FieldVector<Scalar,2*dim>& pressTrace, Scalar press)
286 int eIdxGlobal = problem_.variables().index(element);
288 Dune::FieldVector<Scalar, 2 * dim> faceVol(0);
289 for (
const auto& intersection : intersections(gridView_, element))
291 faceVol[intersection.indexInInside()] = intersection.geometry().volume();
295 for (
int j = 0; j < 2*dim; j++)
296 vel += W_[eIdxGlobal][fIdx][j]*faceVol[j]*(press - pressTrace[j]);
301 const Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim>& W,
302 Dune::FieldVector<Scalar, 2 * dim>& F, Scalar& dInv)
304 Dune::FieldVector<Scalar, 2 * dim> c(0);
305 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> Pi(0);
307 for (
const auto& intersection : intersections(gridView_, element))
310 int idx = intersection.indexInInside();
312 Scalar faceVol = intersection.geometry().volume();
321 Pi[idx][idx] = faceVol;
325 Dune::FieldVector<Scalar, 2 * dim> Wc(0);
327 dInv = 1.0 / (c * Wc);
335 Dune::FieldVector<Scalar, 2 * dim>& F,
336 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim>& Pi,
347 void assembleV(
const Element& element,
int k = 1);
349 void assembleBC(
const Element& element,
int k = 1);
352 Scalar evaluateErrorTerm(CellData& cellData)
357 switch (saturationType)
360 sat = cellData.saturation(wPhaseIdx);
363 sat = cellData.saturation(nPhaseIdx);
366 DUNE_THROW(Dune::NotImplemented,
"Only saturation formulation Sw and Sn are implemented!");
369 Scalar error = (sat > 1.0) ? sat - 1.0 : 0.0;
377 Scalar errorAbs = abs(error);
379 if ((errorAbs * timeStep_ > 1e-6) && (errorAbs > ErrorTermLowerBound_ * maxError_)
380 && (!problem_.timeManager().willBeFinished()))
382 return ErrorTermFactor_ * error;
390 Dune::DynamicVector<Scalar> rhs_;
391 std::vector<Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> > W_;
394 Scalar ErrorTermFactor_;
395 Scalar ErrorTermLowerBound_;
396 Scalar ErrorTermUpperBound_;
398 Scalar density_[numPhases];
399 Scalar viscosity_[numPhases];
403template<
class TypeTag>
404void MimeticTwoPLocalStiffness<TypeTag>::assembleV(
const Element& element,
int)
406 unsigned int numFaces = element.subEntities(1);
407 this->setcurrentsize(numFaces);
409 int eIdxGlobal = problem_.variables().index(element);
414 Dune::FieldVector<Scalar, 2 * dim> faceVol(0);
415 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> Pi(0.);
416 Dune::FieldVector<Scalar, 2 * dim> F(0.);
419 this->assembleElementMatrices(element, faceVol, F, Pi, dInv, qmean);
422 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> PiWPiT(W_[eIdxGlobal]);
423 PiWPiT.rightmultiply(Pi);
424 PiWPiT.leftmultiply(Pi);
427 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> FDinvFT(0);
428 for (
unsigned int i = 0; i < numFaces; i++)
429 for (
unsigned int j = 0; j < numFaces; j++)
430 FDinvFT[i][j] = dInv * F[i] * F[j];
433 for (
unsigned int i = 0; i < numFaces; i++)
434 for (
unsigned int j = 0; j < numFaces; j++)
435 this->A[i][j] = PiWPiT[i][j] - FDinvFT[i][j];
439 Scalar factor = dInv * qmean;
440 for (
unsigned int i = 0; i < numFaces; i++)
441 this->b[i] = F[i] * factor;
451template<
class TypeTag>
453 Dune::FieldVector<Scalar, 2 * dim>& faceVol,
454 Dune::FieldVector<Scalar, 2 * dim>& F,
455 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim>& Pi,
459 unsigned int numFaces = element.subEntities(1);
460 this->setcurrentsize(numFaces);
463 Dune::FieldVector<Scalar, dim> centerGlobal = element.geometry().center();
465 int eIdxGlobal = problem_.variables().index(element);
467 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
470 Scalar volume = element.geometry().volume();
473 Dune::FieldMatrix<Scalar, 2 * dim, dim> R(0), N(0);
478 Scalar gravPotFace[2*dim];
480 Scalar gravPot = (problem_.bBoxMax() - centerGlobal) * problem_.gravity() * (density_[nPhaseIdx] - density_[wPhaseIdx]);
482 for (
const auto& intersection : intersections(gridView_, element))
485 int fIdx = intersection.indexInInside();
487 Dune::FieldVector<Scalar, dim> faceGlobal = intersection.geometry().center();
488 faceVol[fIdx] = intersection.geometry().volume();
491 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
493 N[fIdx] = unitOuterNormal;
495 for (
int k = 0; k < dim; k++)
497 R[fIdx][k] = faceVol[fIdx] * (faceGlobal[k] - centerGlobal[k]);
499 gravPotFace[fIdx] = (problem_.bBoxMax() - faceGlobal) * problem_.gravity() * (density_[nPhaseIdx] - density_[wPhaseIdx]);
505 Scalar norm = R[0][0] * R[0][0];
507 for (
unsigned int i = 1; i < numFaces; i++)
508 norm += R[i][0] * R[i][0];
510 for (
unsigned int i = 0; i < numFaces; i++)
512 Scalar weight = R[0][1] * R[0][0];
513 for (
unsigned int i = 1; i < numFaces; i++)
514 weight += R[i][1] * R[i][0];
515 for (
unsigned int i = 0; i < numFaces; i++)
516 R[i][1] -= weight * R[i][0];
517 norm = R[0][1] * R[0][1];
518 for (
unsigned int i = 1; i < numFaces; i++)
519 norm += R[i][1] * R[i][1];
521 for (
unsigned int i = 0; i < numFaces; i++)
525 Scalar weight1 = R[0][2] * R[0][0];
526 Scalar weight2 = R[0][2] * R[0][1];
527 for (
unsigned int i = 1; i < numFaces; i++)
529 weight1 += R[i][2] * R[i][0];
530 weight2 += R[i][2] * R[i][1];
532 for (
unsigned int i = 0; i < numFaces; i++)
533 R[i][2] -= weight1 * R[i][0] + weight2 * R[i][1];
534 norm = R[0][2] * R[0][2];
535 for (
unsigned int i = 1; i < numFaces; i++)
536 norm += R[i][2] * R[i][2];
539 for (
unsigned int i = 0; i < numFaces; i++)
545 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> D(0);
546 for (
unsigned int s = 0; s < numFaces; s++)
548 Dune::FieldVector<Scalar, 2 * dim> es(0);
550 for (
unsigned int k = 0; k < numFaces; k++)
553 for (
unsigned int i = 0; i < dim; i++)
555 D[k][s] -= R[s][i] * R[k][i];
563 Dune::FieldMatrix<Scalar, dim, dim>
mobility(problem_.spatialParams().intrinsicPermeability(element));
564 mobility *= (cellData.mobility(wPhaseIdx) + cellData.mobility(nPhaseIdx));
567 for (
int i = 1; i < dim; i++)
569 D *= 2.0 * traceMob / volume;
573 Dune::FieldMatrix<Scalar, 2 * dim, dim> NK(N);
576 for (
unsigned int i = 0; i < numFaces; i++)
578 for (
unsigned int j = 0; j < numFaces; j++)
580 W_[eIdxGlobal][i][j] = NK[i][0] * N[j][0];
581 for (
unsigned int k = 1; k < dim; k++)
582 W_[eIdxGlobal][i][j] += NK[i][k] * N[j][k];
586 W_[eIdxGlobal] /= volume;
599 Dune::FieldVector<Scalar, 2 * dim> c(0);
600 for (
unsigned int i = 0; i < numFaces; i++)
605 for (
unsigned int i = 0; i < numFaces; i++)
606 Pi[i][i] = faceVol[i];
609 Dune::FieldVector<Scalar, 2 * dim> Wc(0);
610 W_[eIdxGlobal].umv(c, Wc);
611 dInv = 1.0 / (c * Wc);
619 for (
const auto& intersection : intersections(gridView_, element))
621 int idx = intersection.indexInInside();
626 for (
int j = 0; j < 2 * dim; j++)
627 flux += W_[eIdxGlobal][idx][j] * faceVol[j] * (gravPot - gravPotFace[j]);
630 if (intersection.neighbor())
632 int eIdxGlobalNeighbor = problem_.variables().index(intersection.outside());
635 switch (pressureType)
639 fracFlow = cellData.fracFlowFunc(nPhaseIdx);
644 fracFlow = -cellData.fracFlowFunc(wPhaseIdx);
648 DUNE_THROW(Dune::NotImplemented,
"Only pressure formulation pw and pn are implemented!");
652 rhs_[eIdxGlobal] -= (faceVol[idx] * fracFlow * flux);
653 rhs_[eIdxGlobalNeighbor] += (faceVol[idx] * fracFlow * flux);
656 else if (intersection.boundary())
658 BoundaryTypes bctype;
659 problem_.boundaryTypes(bctype, intersection);
661 if (bctype.isDirichlet(pressureEqIdx))
663 if (flux > 0. || !bctype.isDirichlet(satEqIdx))
665 switch (pressureType)
669 fracFlow = cellData.fracFlowFunc(nPhaseIdx);
674 fracFlow = -cellData.fracFlowFunc(wPhaseIdx);
678 DUNE_THROW(Dune::NotImplemented,
"Only pressure formulation pw and pn are implemented!");
681 else if (flux < 0. && bctype.isDirichlet(satEqIdx))
683 PrimaryVariables boundValues(0.0);
684 problem_.dirichlet(boundValues, intersection);
689 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
691 const Scalar krw = fluidMatrixInteraction.krw(boundValues[saturationIdx]);
692 const Scalar krn = fluidMatrixInteraction.krn(boundValues[saturationIdx]);
694 switch (pressureType)
698 fracFlow = krn / viscosity_[nPhaseIdx]
699 / (krw / viscosity_[wPhaseIdx] + krn / viscosity_[nPhaseIdx]);
704 fracFlow = -krw / viscosity_[wPhaseIdx]
705 / (krw / viscosity_[wPhaseIdx] + krn / viscosity_[nPhaseIdx]);
709 DUNE_THROW(Dune::NotImplemented,
"Only pressure formulation pw and pn are implemented!");
713 rhs_[eIdxGlobal] -= (faceVol[idx] * fracFlow * flux);
718 PrimaryVariables sourceVec(0.0);
719 problem_.source(sourceVec, element);
720 qmean = volume * (sourceVec[wPhaseIdx]/density_[wPhaseIdx] + sourceVec[nPhaseIdx]/density_[nPhaseIdx]);
722 qmean += evaluateErrorTerm(cellData) * volume;
726template<
class TypeTag>
730 for (
const auto& intersection : intersections(gridView_, element))
732 int faceIndex = intersection.indexInInside();
733 if (intersection.boundary())
735 problem_.boundaryTypes(this->bctype[faceIndex], intersection);
736 PrimaryVariables boundValues(0.0);
738 if (this->bctype[faceIndex].isNeumann(pressureEqIdx))
740 problem_.neumann(boundValues, intersection);
741 Scalar J = (boundValues[wPhaseIdx]/density_[wPhaseIdx] + boundValues[nPhaseIdx]/density_[nPhaseIdx]);
742 this->b[faceIndex] -= J * intersection.geometry().volume();
744 else if (this->bctype[faceIndex].isDirichlet(pressureEqIdx))
746 problem_.dirichlet(boundValues, intersection);
747 if (pressureType == pw)
748 this->b[faceIndex] = boundValues[pressureIdx] +
749 (problem_.bBoxMax() - intersection.geometry().center()) * problem_.gravity() * density_[wPhaseIdx];
750 else if (pressureType == pn)
751 this->b[faceIndex] = boundValues[pressureIdx] +
752 (problem_.bBoxMax() - intersection.geometry().center()) * problem_.gravity() * density_[nPhaseIdx];
754 this->b[faceIndex] = boundValues[pressureIdx];
Definition of boundary condition types, extend if necessary.
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
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
const BoundaryTypes & bc(int i) const
Accesses boundary condition for each dof.
Definition: localstiffness.hh:223
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: mimetic.hh:55
void reset()
Definition: mimetic.hh:148
void assemble(const Element &element, int k=1)
Assembles local stiffness matrix for given element and order.
Definition: mimetic.hh:174
MimeticTwoPLocalStiffness(Problem &problem, bool levelBoundaryAsDirichlet, const GridView &gridView, bool procBoundaryAsDirichlet=true)
Constructor.
Definition: mimetic.hh:115
void constructVelocity(const Element &element, int fIdx, Scalar &vel, Dune::FieldVector< Scalar, 2 *dim > &pressTrace, Scalar press)
Definition: mimetic.hh:283
void assembleBoundaryCondition(const Element &element, int k=1)
Assembles only boundary conditions for given element.
Definition: mimetic.hh:213
void constructVelocity(const Element &element, Dune::FieldVector< Scalar, 2 *dim > &vel, Dune::FieldVector< Scalar, 2 *dim > &pressTrace, Scalar press)
Definition: mimetic.hh:265
const Problem & problem() const
Definition: mimetic.hh:341
@ size
Definition: mimetic.hh:111
void completeRHS(const Element &element, Dune::FieldVector< int, 2 *dim > &local2Global, Vector &f)
Definition: mimetic.hh:230
void initialize()
Definition: mimetic.hh:129
void assemble(const Element &cell, const Dune::BlockVector< VBlockType > &localSolution, int orderOfShapeFns=1)
Definition: mimetic.hh:197
void setErrorInfo(Scalar maxErr, Scalar dt)
Definition: mimetic.hh:155
Scalar constructPressure(const Element &element, Dune::FieldVector< Scalar, 2 *dim > &pressTrace)
Definition: mimetic.hh:247
void assembleElementMatrices(const Element &element, Dune::FieldVector< Scalar, 2 *dim > &faceVol, Dune::FieldVector< Scalar, 2 *dim > &F, Dune::FieldMatrix< Scalar, 2 *dim, 2 *dim > &Pi, Scalar &dInv, Scalar &qmean)
Definition: mimetic.hh:452
void computeReconstructionMatrices(const Element &element, const Dune::FieldMatrix< Scalar, 2 *dim, 2 *dim > &W, Dune::FieldVector< Scalar, 2 *dim > &F, Scalar &dInv)
Definition: mimetic.hh:300
@ m
Definition: mimetic.hh:107
Dune::FieldVector< Scalar, m > VBlockType
Definition: mimetic.hh:196
Specifies the properties for immiscible 2p diffusion/pressure models.