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>
51template<
class TypeTag>
63 dim = GridView::dimension
67 wPhaseIdx = Indices::wPhaseIdx,
68 nPhaseIdx = Indices::nPhaseIdx,
69 pressureIdx = Indices::pressureIdx,
70 saturationIdx = Indices::saturationIdx,
71 pressureEqIdx = Indices::pressureEqIdx,
72 satEqIdx = Indices::satEqIdx,
73 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
77 pw = Indices::pressureW,
78 pn = Indices::pressureNw,
79 pGlobal = Indices::pressureGlobal,
80 Sw = Indices::saturationW,
81 Sn = Indices::saturationNw,
82 vw = Indices::velocityW,
83 vn = Indices::velocityNw,
85 pressureType = getPropValue<TypeTag, Properties::PressureFormulation>(),
87 saturationType = getPropValue<TypeTag, Properties::SaturationFormulation>()
90 using Grid =
typename GridView::Grid;
91 using Element =
typename GridView::Traits::template Codim<0>::Entity;
95 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
115 const GridView& gridView,
bool procBoundaryAsDirichlet =
true) :
116 problem_(
problem), gridView_(gridView), maxError_(0), timeStep_(1)
118 ErrorTermFactor_ = getParam<Scalar>(
"Impet.ErrorTermFactor");
119 ErrorTermLowerBound_ = getParam<Scalar>(
"Impet.ErrorTermLowerBound");
120 ErrorTermUpperBound_ = getParam<Scalar>(
"Impet.ErrorTermUpperBound");
122 density_[wPhaseIdx] = 0.0;
123 density_[nPhaseIdx] = 0.0;
124 viscosity_[wPhaseIdx] =0.0;
125 viscosity_[nPhaseIdx] = 0.0;
130 const auto element = *problem_.gridView().template begin<0>();
131 FluidState fluidState;
132 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
133 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
134 fluidState.setTemperature(problem_.temperature(element));
135 fluidState.setSaturation(wPhaseIdx, 1.);
136 fluidState.setSaturation(nPhaseIdx, 0.);
142 int size = gridView_.size(0);
143 rhs_.resize(
size , 0.);
175 unsigned int numFaces = element.subEntities(1);
179 for (
unsigned int i = 0; i < numFaces; i++)
183 for (
unsigned int j = 0; j < numFaces; j++)
187 assembleV(element, k);
188 assembleBC(element, k);
196 void assemble(
const Element &cell,
const Dune::BlockVector<VBlockType>& localSolution,
int orderOfShapeFns = 1)
214 unsigned int numFaces = element.subEntities(1);
218 for (
unsigned int i = 0; i < numFaces; i++)
224 assembleBC(element, k);
228 template<
class Vector>
229 void completeRHS(
const Element& element, Dune::FieldVector<int, 2*dim>& local2Global, Vector& f)
231 int eIdxGlobal = problem_.variables().index(element);
232 unsigned int numFaces = element.subEntities(1);
234 Dune::FieldVector<Scalar, 2 * dim> F(0.);
238 for (
unsigned int i = 0; i < numFaces; i++)
240 if (!this->
bc(i).isDirichlet(pressureEqIdx))
241 f[local2Global[i]][0] += (dInv * F[i] * rhs_[eIdxGlobal]);
248 int eIdxGlobal = problem_.variables().index(element);
249 Scalar volume = element.geometry().volume();
251 PrimaryVariables sourceVec(0.0);
252 problem_.source(sourceVec, element);
253 Scalar qmean = volume * (sourceVec[wPhaseIdx]/density_[wPhaseIdx] + sourceVec[nPhaseIdx]/density_[nPhaseIdx]);
254 qmean += rhs_[eIdxGlobal];
256 Dune::FieldVector<Scalar, 2 * dim> F(0.);
260 return (dInv*(qmean + (F*pressTrace)));
265 Dune::FieldVector<Scalar,2*dim>& pressTrace, Scalar press)
267 int eIdxGlobal = problem_.variables().index(element);
269 Dune::FieldVector<Scalar, 2 * dim> faceVol(0);
270 for (
const auto& intersection : intersections(gridView_, element))
272 faceVol[intersection.indexInInside()] = intersection.geometry().volume();
276 for (
int i = 0; i < 2*dim; i++)
277 for (
int j = 0; j < 2*dim; j++)
278 vel[i] += W_[eIdxGlobal][i][j]*faceVol[j]*(press - pressTrace[j]);
283 Dune::FieldVector<Scalar,2*dim>& pressTrace, Scalar press)
285 int eIdxGlobal = problem_.variables().index(element);
287 Dune::FieldVector<Scalar, 2 * dim> faceVol(0);
288 for (
const auto& intersection : intersections(gridView_, element))
290 faceVol[intersection.indexInInside()] = intersection.geometry().volume();
294 for (
int j = 0; j < 2*dim; j++)
295 vel += W_[eIdxGlobal][fIdx][j]*faceVol[j]*(press - pressTrace[j]);
300 const Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim>& W,
301 Dune::FieldVector<Scalar, 2 * dim>& F, Scalar& dInv)
303 Dune::FieldVector<Scalar, 2 * dim> c(0);
304 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> Pi(0);
306 for (
const auto& intersection : intersections(gridView_, element))
309 int idx = intersection.indexInInside();
311 Scalar faceVol = intersection.geometry().volume();
320 Pi[idx][idx] = faceVol;
324 Dune::FieldVector<Scalar, 2 * dim> Wc(0);
326 dInv = 1.0 / (c * Wc);
334 Dune::FieldVector<Scalar, 2 * dim>& F,
335 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim>& Pi,
346 void assembleV(
const Element& element,
int k = 1);
348 void assembleBC(
const Element& element,
int k = 1);
351 Scalar evaluateErrorTerm(CellData& cellData)
356 switch (saturationType)
359 sat = cellData.saturation(wPhaseIdx);
362 sat = cellData.saturation(nPhaseIdx);
365 DUNE_THROW(Dune::NotImplemented,
"Only saturation formulation Sw and Sn are implemented!");
368 Scalar error = (sat > 1.0) ? sat - 1.0 : 0.0;
376 Scalar errorAbs = abs(error);
378 if ((errorAbs * timeStep_ > 1e-6) && (errorAbs > ErrorTermLowerBound_ * maxError_)
379 && (!problem_.timeManager().willBeFinished()))
381 return ErrorTermFactor_ * error;
389 Dune::DynamicVector<Scalar> rhs_;
390 std::vector<Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> > W_;
393 Scalar ErrorTermFactor_;
394 Scalar ErrorTermLowerBound_;
395 Scalar ErrorTermUpperBound_;
397 Scalar density_[numPhases];
398 Scalar viscosity_[numPhases];
402template<
class TypeTag>
403void MimeticTwoPLocalStiffness<TypeTag>::assembleV(
const Element& element,
int)
405 unsigned int numFaces = element.subEntities(1);
406 this->setcurrentsize(numFaces);
408 int eIdxGlobal = problem_.variables().index(element);
413 Dune::FieldVector<Scalar, 2 * dim> faceVol(0);
414 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> Pi(0.);
415 Dune::FieldVector<Scalar, 2 * dim> F(0.);
418 this->assembleElementMatrices(element, faceVol, F, Pi, dInv, qmean);
421 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> PiWPiT(W_[eIdxGlobal]);
422 PiWPiT.rightmultiply(Pi);
423 PiWPiT.leftmultiply(Pi);
426 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> FDinvFT(0);
427 for (
unsigned int i = 0; i < numFaces; i++)
428 for (
unsigned int j = 0; j < numFaces; j++)
429 FDinvFT[i][j] = dInv * F[i] * F[j];
432 for (
unsigned int i = 0; i < numFaces; i++)
433 for (
unsigned int j = 0; j < numFaces; j++)
434 this->A[i][j] = PiWPiT[i][j] - FDinvFT[i][j];
438 Scalar factor = dInv * qmean;
439 for (
unsigned int i = 0; i < numFaces; i++)
440 this->b[i] = F[i] * factor;
450template<
class TypeTag>
452 Dune::FieldVector<Scalar, 2 * dim>& faceVol,
453 Dune::FieldVector<Scalar, 2 * dim>& F,
454 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim>& Pi,
458 unsigned int numFaces = element.subEntities(1);
459 this->setcurrentsize(numFaces);
462 Dune::FieldVector<Scalar, dim> centerGlobal = element.geometry().center();
464 int eIdxGlobal = problem_.variables().index(element);
466 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
469 Scalar volume = element.geometry().volume();
472 Dune::FieldMatrix<Scalar, 2 * dim, dim> R(0), N(0);
477 Scalar gravPotFace[2*dim];
479 Scalar gravPot = (problem_.bBoxMax() - centerGlobal) * problem_.gravity() * (density_[nPhaseIdx] - density_[wPhaseIdx]);
481 for (
const auto& intersection : intersections(gridView_, element))
484 int fIdx = intersection.indexInInside();
486 Dune::FieldVector<Scalar, dim> faceGlobal = intersection.geometry().center();
487 faceVol[fIdx] = intersection.geometry().volume();
490 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
492 N[fIdx] = unitOuterNormal;
494 for (
int k = 0; k < dim; k++)
496 R[fIdx][k] = faceVol[fIdx] * (faceGlobal[k] - centerGlobal[k]);
498 gravPotFace[fIdx] = (problem_.bBoxMax() - faceGlobal) * problem_.gravity() * (density_[nPhaseIdx] - density_[wPhaseIdx]);
504 Scalar norm = R[0][0] * R[0][0];
506 for (
unsigned int i = 1; i < numFaces; i++)
507 norm += R[i][0] * R[i][0];
509 for (
unsigned int i = 0; i < numFaces; i++)
511 Scalar weight = R[0][1] * R[0][0];
512 for (
unsigned int i = 1; i < numFaces; i++)
513 weight += R[i][1] * R[i][0];
514 for (
unsigned int i = 0; i < numFaces; i++)
515 R[i][1] -= weight * R[i][0];
516 norm = R[0][1] * R[0][1];
517 for (
unsigned int i = 1; i < numFaces; i++)
518 norm += R[i][1] * R[i][1];
520 for (
unsigned int i = 0; i < numFaces; i++)
524 Scalar weight1 = R[0][2] * R[0][0];
525 Scalar weight2 = R[0][2] * R[0][1];
526 for (
unsigned int i = 1; i < numFaces; i++)
528 weight1 += R[i][2] * R[i][0];
529 weight2 += R[i][2] * R[i][1];
531 for (
unsigned int i = 0; i < numFaces; i++)
532 R[i][2] -= weight1 * R[i][0] + weight2 * R[i][1];
533 norm = R[0][2] * R[0][2];
534 for (
unsigned int i = 1; i < numFaces; i++)
535 norm += R[i][2] * R[i][2];
538 for (
unsigned int i = 0; i < numFaces; i++)
544 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> D(0);
545 for (
unsigned int s = 0; s < numFaces; s++)
547 Dune::FieldVector<Scalar, 2 * dim> es(0);
549 for (
unsigned int k = 0; k < numFaces; k++)
552 for (
unsigned int i = 0; i < dim; i++)
554 D[k][s] -= R[s][i] * R[k][i];
562 Dune::FieldMatrix<Scalar, dim, dim>
mobility(problem_.spatialParams().intrinsicPermeability(element));
563 mobility *= (cellData.mobility(wPhaseIdx) + cellData.mobility(nPhaseIdx));
566 for (
int i = 1; i < dim; i++)
568 D *= 2.0 * traceMob / volume;
572 Dune::FieldMatrix<Scalar, 2 * dim, dim> NK(N);
575 for (
unsigned int i = 0; i < numFaces; i++)
577 for (
unsigned int j = 0; j < numFaces; j++)
579 W_[eIdxGlobal][i][j] = NK[i][0] * N[j][0];
580 for (
unsigned int k = 1; k < dim; k++)
581 W_[eIdxGlobal][i][j] += NK[i][k] * N[j][k];
585 W_[eIdxGlobal] /= volume;
598 Dune::FieldVector<Scalar, 2 * dim> c(0);
599 for (
unsigned int i = 0; i < numFaces; i++)
604 for (
unsigned int i = 0; i < numFaces; i++)
605 Pi[i][i] = faceVol[i];
608 Dune::FieldVector<Scalar, 2 * dim> Wc(0);
609 W_[eIdxGlobal].umv(c, Wc);
610 dInv = 1.0 / (c * Wc);
618 for (
const auto& intersection : intersections(gridView_, element))
620 int idx = intersection.indexInInside();
625 for (
int j = 0; j < 2 * dim; j++)
626 flux += W_[eIdxGlobal][idx][j] * faceVol[j] * (gravPot - gravPotFace[j]);
629 if (intersection.neighbor())
631 int eIdxGlobalNeighbor = problem_.variables().index(intersection.outside());
634 switch (pressureType)
638 fracFlow = cellData.fracFlowFunc(nPhaseIdx);
643 fracFlow = -cellData.fracFlowFunc(wPhaseIdx);
647 DUNE_THROW(Dune::NotImplemented,
"Only pressure formulation pw and pn are implemented!");
651 rhs_[eIdxGlobal] -= (faceVol[idx] * fracFlow * flux);
652 rhs_[eIdxGlobalNeighbor] += (faceVol[idx] * fracFlow * flux);
655 else if (intersection.boundary())
657 BoundaryTypes bctype;
658 problem_.boundaryTypes(bctype, intersection);
660 if (bctype.isDirichlet(pressureEqIdx))
662 if (flux > 0. || !bctype.isDirichlet(satEqIdx))
664 switch (pressureType)
668 fracFlow = cellData.fracFlowFunc(nPhaseIdx);
673 fracFlow = -cellData.fracFlowFunc(wPhaseIdx);
677 DUNE_THROW(Dune::NotImplemented,
"Only pressure formulation pw and pn are implemented!");
680 else if (flux < 0. && bctype.isDirichlet(satEqIdx))
682 PrimaryVariables boundValues(0.0);
683 problem_.dirichlet(boundValues, intersection);
685 Scalar krw = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element),
686 boundValues[saturationIdx]);
687 Scalar krn = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element),
688 boundValues[saturationIdx]);
690 switch (pressureType)
694 fracFlow = krn / viscosity_[nPhaseIdx]
695 / (krw / viscosity_[wPhaseIdx] + krn / viscosity_[nPhaseIdx]);
700 fracFlow = -krw / viscosity_[wPhaseIdx]
701 / (krw / viscosity_[wPhaseIdx] + krn / viscosity_[nPhaseIdx]);
705 DUNE_THROW(Dune::NotImplemented,
"Only pressure formulation pw and pn are implemented!");
709 rhs_[eIdxGlobal] -= (faceVol[idx] * fracFlow * flux);
714 PrimaryVariables sourceVec(0.0);
715 problem_.source(sourceVec, element);
716 qmean = volume * (sourceVec[wPhaseIdx]/density_[wPhaseIdx] + sourceVec[nPhaseIdx]/density_[nPhaseIdx]);
718 qmean += evaluateErrorTerm(cellData) * volume;
722template<
class TypeTag>
726 for (
const auto& intersection : intersections(gridView_, element))
728 int faceIndex = intersection.indexInInside();
729 if (intersection.boundary())
731 problem_.boundaryTypes(this->bctype[faceIndex], intersection);
732 PrimaryVariables boundValues(0.0);
734 if (this->bctype[faceIndex].isNeumann(pressureEqIdx))
736 problem_.neumann(boundValues, intersection);
737 Scalar J = (boundValues[wPhaseIdx]/density_[wPhaseIdx] + boundValues[nPhaseIdx]/density_[nPhaseIdx]);
738 this->b[faceIndex] -= J * intersection.geometry().volume();
740 else if (this->bctype[faceIndex].isDirichlet(pressureEqIdx))
742 problem_.dirichlet(boundValues, intersection);
743 if (pressureType == pw)
744 this->b[faceIndex] = boundValues[pressureIdx] +
745 (problem_.bBoxMax() - intersection.geometry().center()) * problem_.gravity() * density_[wPhaseIdx];
746 else if (pressureType == pn)
747 this->b[faceIndex] = boundValues[pressureIdx] +
748 (problem_.bBoxMax() - intersection.geometry().center()) * problem_.gravity() * density_[nPhaseIdx];
750 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:53
void reset()
Definition: mimetic.hh:147
void assemble(const Element &element, int k=1)
Assembles local stiffness matrix for given element and order.
Definition: mimetic.hh:173
MimeticTwoPLocalStiffness(Problem &problem, bool levelBoundaryAsDirichlet, const GridView &gridView, bool procBoundaryAsDirichlet=true)
Constructor.
Definition: mimetic.hh:114
@ size
Definition: mimetic.hh:110
void constructVelocity(const Element &element, int fIdx, Scalar &vel, Dune::FieldVector< Scalar, 2 *dim > &pressTrace, Scalar press)
Definition: mimetic.hh:282
void assembleBoundaryCondition(const Element &element, int k=1)
Assembles only boundary conditions for given element.
Definition: mimetic.hh:212
void constructVelocity(const Element &element, Dune::FieldVector< Scalar, 2 *dim > &vel, Dune::FieldVector< Scalar, 2 *dim > &pressTrace, Scalar press)
Definition: mimetic.hh:264
const Problem & problem() const
Definition: mimetic.hh:340
void completeRHS(const Element &element, Dune::FieldVector< int, 2 *dim > &local2Global, Vector &f)
Definition: mimetic.hh:229
void initialize()
Definition: mimetic.hh:128
void assemble(const Element &cell, const Dune::BlockVector< VBlockType > &localSolution, int orderOfShapeFns=1)
Definition: mimetic.hh:196
void setErrorInfo(Scalar maxErr, Scalar dt)
Definition: mimetic.hh:154
Scalar constructPressure(const Element &element, Dune::FieldVector< Scalar, 2 *dim > &pressTrace)
Definition: mimetic.hh:246
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:451
@ m
Definition: mimetic.hh:106
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:299
Dune::FieldVector< Scalar, m > VBlockType
Definition: mimetic.hh:195
Specifies the properties for immiscible 2p diffusion/pressure models.