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>
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;
119 const GridView& gridView,
const IntersectionMapper& intersectionMapper,
bool procBoundaryAsDirichlet =
true) :
120 problem_(
problem), gridView_(gridView), intersectionMapper_(intersectionMapper), maxError_(0), timeStep_(1)
122 int size = gridView_.size(0);
123 rhs_.resize(
size , 0.);
125 ErrorTermFactor_ = getParam<Scalar>(
"Impet.ErrorTermFactor");
126 ErrorTermLowerBound_ = getParam<Scalar>(
"Impet.ErrorTermLowerBound");
127 ErrorTermUpperBound_ = getParam<Scalar>(
"Impet.ErrorTermUpperBound");
129 density_[wPhaseIdx] = 0.0;
130 density_[nPhaseIdx] = 0.0;
131 viscosity_[wPhaseIdx] =0.0;
132 viscosity_[nPhaseIdx] = 0.0;
137 const auto element = *problem_.gridView().template begin<0>();
138 FluidState fluidState;
139 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
140 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
141 fluidState.setTemperature(problem_.temperature(element));
142 fluidState.setSaturation(wPhaseIdx, 1.);
143 fluidState.setSaturation(nPhaseIdx, 0.);
154 int size = gridView_.size(0);
155 rhs_.resize(
size , 0.);
174 return W_[eIdxGlobal].N();
191 int numFaces = intersectionMapper_.
size(element);
196 for (
int i = 0; i < numFaces; i++)
200 for (
int j = 0; j < numFaces; j++)
204 assembleV(element, numFaces, k);
205 assembleBC(element, k);
213 void assemble(
const Element &cell,
const Dune::BlockVector<VBlockType>& localSolution,
int orderOfShapeFns = 1)
230 int numFaces = intersectionMapper_.
size(element);
233 for (
int i = 0; i < numFaces; i++)
239 assembleBC(element, k);
242 template<
class Vector>
243 void completeRHS(
const Element& element, std::vector<int>& local2Global, Vector& f)
245 int eIdxGlobal = problem_.variables().index(element);
249 Dune::DynamicVector<Scalar> F(numFaces);
253 for (
int i = 0; i < numFaces; i++)
256 f[local2Global[i]][0] += (dInv * F[i] * rhs_[eIdxGlobal]);
262 int eIdxGlobal = problem_.variables().index(element);
266 Scalar volume = element.geometry().volume();
268 PrimaryVariables sourceVec(0.0);
269 problem_.source(sourceVec, element);
270 Scalar qmean = volume * (sourceVec[wPhaseIdx]/density_[wPhaseIdx] + sourceVec[nPhaseIdx]/density_[nPhaseIdx]);
271 qmean += rhs_[eIdxGlobal];
273 Dune::DynamicVector<Scalar> F(numFaces, 0.);
277 return (dInv*(qmean + (F*pressTrace)));
281 Dune::DynamicVector<Scalar>& pressTrace, Scalar press)
283 int eIdxGlobal = problem_.variables().index(element);
287 Dune::DynamicVector<Scalar> faceVol(numFaces);
289 Dune::FieldVector<Scalar, 2*dim> faceVolumeReal(0.0);
292 for (
const auto& intersection : intersections(gridView_, element))
294 faceVol[idx] = intersection.geometry().volume();
295 int indexInInside = intersection.indexInInside();
296 faceVolumeReal[indexInInside] += faceVol[idx];
302 for (
int i = 0; i < numFaces; i++)
303 for (
int j = 0; j < numFaces; j++)
304 vel[i] += W_[eIdxGlobal][i][j]*faceVol[j]*(press - pressTrace[j]);
306 for (
int i = 0; i < numFaces; i++)
308 vel[i] *= faceVol[i]/faceVolumeReal[intersectionMapper_.
maplocal(eIdxGlobal, i)];
314 Dune::DynamicVector<Scalar>& F, Scalar& dInv)
316 int eIdxGlobal = problem_.variables().index(element);
320 Dune::DynamicVector<Scalar> c(numFaces);
321 Dune::DynamicMatrix<Scalar> Pi(numFaces, numFaces);
324 for (
const auto& intersection : intersections(gridView_, element))
326 Scalar faceVol = intersection.geometry().volume();
335 Pi[idx][idx] = faceVol;
340 Dune::DynamicVector<Scalar> Wc(numFaces);
342 dInv = 1.0 / (c * Wc);
350 Dune::DynamicVector<Scalar>& faceVol,
351 Dune::DynamicVector<Scalar>& F,
352 Dune::DynamicMatrix<Scalar>& Pi,
363 void assembleV(
const Element& element,
int numFaces,
int k = 1);
365 void assembleBC(
const Element& element,
int k = 1);
367 Scalar evaluateErrorTerm(CellData& cellData)
372 switch (saturationType)
375 sat = cellData.saturation(wPhaseIdx);
378 sat = cellData.saturation(nPhaseIdx);
381 DUNE_THROW(Dune::NotImplemented,
"Only saturation formulation Sw and Sn are implemented!");
384 Scalar error = (sat > 1.0) ? sat - 1.0 : 0.0;
392 Scalar errorAbs = abs(error);
394 if ((errorAbs * timeStep_ > 1e-6) && (errorAbs > ErrorTermLowerBound_ * maxError_)
395 && (!problem_.timeManager().willBeFinished()))
397 return ErrorTermFactor_ * error;
404 const GridView gridView_;
405 const IntersectionMapper& intersectionMapper_;
406 Dune::DynamicVector<Scalar> rhs_;
407 std::vector<Dune::DynamicMatrix<Scalar> > W_;
410 Scalar ErrorTermFactor_;
411 Scalar ErrorTermLowerBound_;
412 Scalar ErrorTermUpperBound_;
414 Scalar density_[numPhases];
415 Scalar viscosity_[numPhases];
419template<
class TypeTag>
420void MimeticTwoPLocalStiffnessAdaptive<TypeTag>::assembleV(
const Element& element,
int numFaces,
int k)
422 int eIdxGlobal = problem_.variables().index(element);
427 W_[eIdxGlobal].resize(numFaces, numFaces);
428 Dune::DynamicVector<Scalar> faceVol(numFaces);
429 Dune::DynamicMatrix<Scalar> Pi(numFaces, numFaces);
430 Dune::DynamicVector<Scalar> F(numFaces);
433 this->assembleElementMatrices(element, faceVol, F, Pi, dInv, qmean);
436 Dune::DynamicMatrix<Scalar> PiWPiT(W_[eIdxGlobal]);
437 PiWPiT.rightmultiply(Pi);
438 PiWPiT.leftmultiply(Pi);
441 Dune::DynamicMatrix<Scalar> FDinvFT(numFaces, numFaces);
442 for (
int i = 0; i < numFaces; i++)
443 for (
int j = 0; j < numFaces; j++)
444 FDinvFT[i][j] = dInv * F[i] * F[j];
447 for (
int i = 0; i < numFaces; i++)
448 for (
int j = 0; j < numFaces; j++)
449 this->A[i][j] = PiWPiT[i][j] - FDinvFT[i][j];
453 Scalar factor = dInv * qmean;
454 for (
int i = 0; i < numFaces; i++)
455 this->b[i] = F[i] * factor;
465template<
class TypeTag>
467 Dune::DynamicVector<Scalar>& faceVol,
468 Dune::DynamicVector<Scalar>& F,
469 Dune::DynamicMatrix<Scalar>& Pi,
474 const Dune::FieldVector<Scalar, dim>& centerGlobal = element.geometry().center();
476 int eIdxGlobal = problem_.variables().index(element);
478 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
480 int numFaces = intersectionMapper_.size(eIdxGlobal);
483 Scalar volume = element.geometry().volume();
486 Dune::DynamicMatrix<Scalar> R(numFaces, dim, 0.);
487 Dune::DynamicMatrix<Scalar> N(numFaces, dim, 0.);
492 std::vector<Scalar> pcPotFace(numFaces);
493 Scalar pcPot = (problem_.bBoxMax() - element.geometry().center()) * problem_.gravity()
494 * (density_[nPhaseIdx] - density_[wPhaseIdx]);
497 for (
const auto& intersection : intersections(gridView_, element))
501 Dune::FieldVector<Scalar, dim> faceGlobal(intersection.geometry().center());
502 faceVol[idx] = intersection.geometry().volume();
505 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
508 for (
int k = 0; k < dim; k++)
511 R[idx][k] = faceVol[idx] * (faceGlobal[k] - centerGlobal[k]);
512 N[idx][k] = unitOuterNormal[k];
515 pcPotFace[idx] = (problem_.bBoxMax() - faceGlobal) * problem_.gravity()
516 * (density_[nPhaseIdx] - density_[wPhaseIdx]);
524 Scalar norm = R[0][0] * R[0][0];
526 for (
int i = 1; i < numFaces; i++)
527 norm += R[i][0] * R[i][0];
529 for (
int i = 0; i < numFaces; i++)
531 Scalar weight = R[0][1] * R[0][0];
532 for (
int i = 1; i < numFaces; i++)
533 weight += R[i][1] * R[i][0];
534 for (
int i = 0; i < numFaces; i++)
535 R[i][1] -= weight * R[i][0];
536 norm = R[0][1] * R[0][1];
537 for (
int i = 1; i < numFaces; i++)
538 norm += R[i][1] * R[i][1];
540 for (
int i = 0; i < numFaces; i++)
544 Scalar weight1 = R[0][2] * R[0][0];
545 Scalar weight2 = R[0][2] * R[0][1];
546 for (
int i = 1; i < numFaces; i++)
548 weight1 += R[i][2] * R[i][0];
549 weight2 += R[i][2] * R[i][1];
551 for (
int i = 0; i < numFaces; i++)
552 R[i][2] -= weight1 * R[i][0] + weight2 * R[i][1];
553 norm = R[0][2] * R[0][2];
554 for (
int i = 1; i < numFaces; i++)
555 norm += R[i][2] * R[i][2];
558 for (
int i = 0; i < numFaces; i++)
564 Dune::DynamicMatrix<Scalar> D(numFaces, numFaces, 0.);
565 for (
int s = 0; s < numFaces; s++)
567 Dune::DynamicVector<Scalar> es(numFaces, 0.);
569 for (
int k = 0; k < numFaces; k++)
572 for (
int i = 0; i < dim; i++)
574 D[k][s] -= R[s][i] * R[k][i];
581 Dune::FieldMatrix<Scalar, dim, dim>
mobility(problem_.spatialParams().intrinsicPermeability(element));
582 mobility *= (cellData.mobility(wPhaseIdx) + cellData.mobility(nPhaseIdx));
585 for (
int i = 1; i < dim; i++)
587 D *= 2.0 * traceMob / volume;
591 Dune::DynamicMatrix<Scalar> NK(N);
596 for (
int i = 0; i < numFaces; i++)
598 for (
int j = 0; j < numFaces; j++)
600 W_[eIdxGlobal][i][j] = NK[i][0] * N[j][0];
601 for (
int k = 1; k < dim; k++)
602 W_[eIdxGlobal][i][j] += NK[i][k] * N[j][k];
606 W_[eIdxGlobal] /= volume;
619 Dune::DynamicVector<Scalar> c(numFaces);
620 for (
int i = 0; i < numFaces; i++)
625 for (
int i = 0; i < numFaces; i++)
626 Pi[i][i] = faceVol[i];
629 Dune::DynamicVector<Scalar> Wc(numFaces);
630 W_[eIdxGlobal].umv(c, Wc);
631 dInv = 1.0 / (c * Wc);
640 for (
const auto& intersection : intersections(gridView_, element))
645 for (
int j = 0; j < numFaces; j++)
646 flux += W_[eIdxGlobal][idx][j] * faceVol[j] * (pcPot - pcPotFace[j]);
648 if (intersection.neighbor())
651 int eIdxGlobalNeighbor = problem_.variables().index(intersection.outside());
654 switch (pressureType)
658 fracFlow = cellData.fracFlowFunc(nPhaseIdx);
663 fracFlow = -cellData.fracFlowFunc(wPhaseIdx);
667 DUNE_THROW(Dune::NotImplemented,
"Only pressure formulation pw and pn are implemented!");
670 rhs_[eIdxGlobal] -= (faceVol[idx] * fracFlow * flux);
671 rhs_[eIdxGlobalNeighbor] += (faceVol[idx] * fracFlow * flux);
674 else if (intersection.boundary())
676 BoundaryTypes bctype;
677 problem_.boundaryTypes(bctype, intersection);
679 if (bctype.isDirichlet(pressureEqIdx))
681 if (flux > 0. || !bctype.isDirichlet(satEqIdx))
683 switch (pressureType)
687 fracFlow = cellData.fracFlowFunc(nPhaseIdx);
692 fracFlow = -cellData.fracFlowFunc(wPhaseIdx);
696 DUNE_THROW(Dune::NotImplemented,
"Only pressure formulation pw and pn are implemented!");
699 else if (flux < 0. && bctype.isDirichlet(satEqIdx))
701 PrimaryVariables boundValues(0.0);
702 problem_.dirichlet(boundValues, intersection);
704 Scalar krw = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element),
705 boundValues[saturationIdx]);
706 Scalar krn = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element),
707 boundValues[saturationIdx]);
709 switch (pressureType)
713 fracFlow = krn / viscosity_[nPhaseIdx]
714 / (krw / viscosity_[wPhaseIdx] + krn / viscosity_[nPhaseIdx]);
719 fracFlow = -krw / viscosity_[wPhaseIdx]
720 / (krw / viscosity_[wPhaseIdx] + krn / viscosity_[nPhaseIdx]);
724 DUNE_THROW(Dune::NotImplemented,
"Only pressure formulation pw and pn are implemented!");
728 rhs_[eIdxGlobal] -= (faceVol[idx] * fracFlow * flux);
734 PrimaryVariables sourceVec(0.0);
735 problem_.source(sourceVec, element);
736 qmean = volume * (sourceVec[wPhaseIdx]/density_[wPhaseIdx] + sourceVec[nPhaseIdx]/density_[nPhaseIdx]);
738 qmean += evaluateErrorTerm(cellData) * volume;
742template<
class TypeTag>
746 unsigned int faceIndex = 0;
747 for (
const auto& intersection : intersections(gridView_, element))
749 if (intersection.neighbor())
753 else if (intersection.boundary())
755 problem_.boundaryTypes(this->bctype[faceIndex], intersection);
756 PrimaryVariables boundValues(0.0);
758 if (this->bctype[faceIndex].isNeumann(pressureEqIdx))
760 problem_.neumann(boundValues, intersection);
761 Scalar J = (boundValues[wPhaseIdx]/density_[wPhaseIdx] + boundValues[nPhaseIdx]/density_[nPhaseIdx]);
762 this->b[faceIndex] -= J * intersection.geometry().volume();
764 else if (this->bctype[faceIndex].isDirichlet(pressureEqIdx))
766 problem_.dirichlet(boundValues, intersection);
767 if (pressureType == pw)
768 this->b[faceIndex] = boundValues[pressureIdx] + (problem_.bBoxMax()
769 - intersection.geometry().center()) * problem_.gravity() * density_[wPhaseIdx];
770 else if (pressureType == pn)
771 this->b[faceIndex] = boundValues[pressureIdx] + (problem_.bBoxMax()
772 - intersection.geometry().center()) * problem_.gravity() * density_[nPhaseIdx];
774 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:55
Scalar constructPressure(const Element &element, Dune::DynamicVector< Scalar > &pressTrace)
Definition: mimeticadaptive.hh:260
void assembleBoundaryCondition(const Element &element, int k=1)
Assembles only boundary conditions for given element.
Definition: mimeticadaptive.hh:228
void initialize()
Definition: mimeticadaptive.hh:135
void computeReconstructionMatrices(const Element &element, const Dune::DynamicMatrix< Scalar > &W, Dune::DynamicVector< Scalar > &F, Scalar &dInv)
Definition: mimeticadaptive.hh:313
void reset()
Definition: mimeticadaptive.hh:159
void setErrorInfo(Scalar maxErr, Scalar dt)
Definition: mimeticadaptive.hh:166
void assemble(const Element &cell, const Dune::BlockVector< VBlockType > &localSolution, int orderOfShapeFns=1)
Definition: mimeticadaptive.hh:213
@ m
Definition: mimeticadaptive.hh:110
void adapt()
Definition: mimeticadaptive.hh:152
void completeRHS(const Element &element, std::vector< int > &local2Global, Vector &f)
Definition: mimeticadaptive.hh:243
void assemble(const Element &element, int k=1)
Assembles local stiffness matrix for given element and order.
Definition: mimeticadaptive.hh:189
void constructVelocity(const Element &element, Dune::DynamicVector< Scalar > &vel, Dune::DynamicVector< Scalar > &pressTrace, Scalar press)
Definition: mimeticadaptive.hh:280
const Problem & problem() const
Definition: mimeticadaptive.hh:357
int numberOfFaces(int eIdxGlobal)
Definition: mimeticadaptive.hh:172
MimeticTwoPLocalStiffnessAdaptive(Problem &problem, bool levelBoundaryAsDirichlet, const GridView &gridView, const IntersectionMapper &intersectionMapper, bool procBoundaryAsDirichlet=true)
Constructor.
Definition: mimeticadaptive.hh:118
Dune::FieldVector< Scalar, m > VBlockType
Definition: mimeticadaptive.hh:212
@ size
Definition: mimeticadaptive.hh:114
void assembleElementMatrices(const Element &element, Dune::DynamicVector< Scalar > &faceVol, Dune::DynamicVector< Scalar > &F, Dune::DynamicMatrix< Scalar > &Pi, Scalar &dInv, Scalar &qmean)
Definition: mimeticadaptive.hh:466
Specifies the properties for immiscible 2p diffusion/pressure models.