25#ifndef DUMUX_FVPRESSURE1P_HH
26#define DUMUX_FVPRESSURE1P_HH
64 using Indices =
typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
66 using FluidSystem =
typename GET_PROP_TYPE(TypeTag, FluidSystem);
68 using BoundaryTypes =
typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
69 using SolutionTypes =
typename GET_PROP(TypeTag, SolutionTypes);
70 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
72 using ScalarSolutionType =
typename SolutionTypes::ScalarSolution;
76 dim = GridView::dimension, dimWorld = GridView::dimensionworld
81 pressEqIdx = Indices::pressureEqIdx
89 using Element =
typename GridView::Traits::template Codim<0>::Entity;
90 using Intersection =
typename GridView::Intersection;
92 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
93 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
94 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
99 void getSource(Dune::FieldVector<Scalar, 2>&,
const Element&,
const CellData&,
const bool);
102 void getStorage(Dune::FieldVector<Scalar, 2>& entry,
const Element& element,
103 const CellData& cellData,
const bool first)
109 void getFlux(Dune::FieldVector<Scalar, 2>&,
const Intersection&,
const CellData&,
const bool);
112 const Intersection&,
const CellData&,
const bool);
151 int size = problem_.gridView().size(0);
152 for (
int i = 0; i < size; i++)
154 CellData& cellData = problem_.variables().cellData(i);
156 cellData.fluxData().resetVelocity();
168 Scalar press = this->
pressure()[eIdxGlobal];
170 cellData.setPressure(press);
181 template<
class MultiWriter>
184 ScalarSolutionType *
pressure = writer.allocateManagedBuffer (
185 problem_.gridView().size(0));
189 writer.attachCellData(*
pressure,
"pressure");
204 auto element = *problem_.gridView().template begin<0> ();
205 Scalar
temperature = problem_.temperature(element);
206 Scalar referencePress = problem_.referencePressure(element);
214 const GravityVector& gravity_;
226template<
class TypeTag>
228 ,
const CellData& cellData,
const bool first)
231 Scalar volume = element.geometry().volume();
234 PrimaryVariables sourcePhase(0.0);
235 problem_.source(sourcePhase, element);
236 sourcePhase /= density_;
238 entry[rhs] = volume * sourcePhase;
248template<
class TypeTag>
250 ,
const CellData& cellData,
const bool first)
252 auto elementI = intersection.inside();
253 auto elementJ = intersection.outside();
256 const GlobalPosition& globalPosI = elementI.geometry().center();
257 const GlobalPosition& globalPosJ = elementJ.geometry().center();
260 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
263 Scalar faceArea = intersection.geometry().volume();
266 GlobalPosition distVec = globalPosJ - globalPosI;
269 Scalar dist = distVec.two_norm();
272 DimMatrix meanPermeability(0);
274 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(elementI),
275 problem_.spatialParams().intrinsicPermeability(elementJ));
283 entry[matrix] = ((
permeability * unitOuterNormal) / dist) * faceArea;
286 entry[rhs] = density_ * (
permeability * gravity_) * faceArea;
299template<
class TypeTag>
301const Intersection& intersection,
const CellData& cellData,
const bool first)
303 auto element = intersection.inside();
306 const GlobalPosition& globalPosI = element.geometry().center();
309 const GlobalPosition& globalPosJ = intersection.geometry().center();
312 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
315 Scalar faceArea = intersection.geometry().volume();
318 GlobalPosition distVec = globalPosJ - globalPosI;
321 Scalar dist = distVec.two_norm();
323 BoundaryTypes bcType;
324 problem_.boundaryTypes(bcType, intersection);
325 PrimaryVariables boundValues(0.0);
327 if (bcType.isDirichlet(pressEqIdx))
329 problem_.dirichlet(boundValues, intersection);
333 DimMatrix meanPermeability(0);
335 problem_.spatialParams().meanK(meanPermeability,
336 problem_.spatialParams().intrinsicPermeability(element));
344 Scalar pressBound = boundValues;
347 entry[matrix] = ((
permeability * unitOuterNormal) / dist) * faceArea;
348 entry[rhs] = entry[matrix] * pressBound;
351 entry[rhs] -= density_ * (
permeability * gravity_) * faceArea;
355 else if (bcType.isNeumann(pressEqIdx))
357 problem_.neumann(boundValues, intersection);
358 Scalar J = boundValues /= density_;
360 entry[rhs] = -(J * faceArea);
364 DUNE_THROW(Dune::NotImplemented,
"No valid boundary condition type defined for pressure equation!");
#define GET_PROP(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:281
#define GET_PROP_TYPE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:283
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string permeability() noexcept
I/O name of permeability.
Definition: name.hh:143
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Sequential OneP Model solving the equations for pressure and velocity separately.
Definition: 1p/sequential/diffusion/cellcentered/pressure.hh:57
void initialize(bool solveTwice=true)
Initializes the pressure model.
Definition: 1p/sequential/diffusion/cellcentered/pressure.hh:121
FVPressure1P(Problem &problem)
Constructs a FVPressure1P object.
Definition: 1p/sequential/diffusion/cellcentered/pressure.hh:199
void storePressureSolution(int eIdxGlobal, CellData &cellData)
Stores the pressure solution of a cell.
Definition: 1p/sequential/diffusion/cellcentered/pressure.hh:166
void getSource(Dune::FieldVector< Scalar, 2 > &, const Element &, const CellData &, const bool)
Function which calculates the source entry.
Definition: 1p/sequential/diffusion/cellcentered/pressure.hh:227
void storePressureSolution()
Globally stores the pressure solution.
Definition: 1p/sequential/diffusion/cellcentered/pressure.hh:149
void update()
Pressure update.
Definition: 1p/sequential/diffusion/cellcentered/pressure.hh:140
void addOutputVtkFields(MultiWriter &writer)
Adds pressure output to the output file.
Definition: 1p/sequential/diffusion/cellcentered/pressure.hh:182
void getFluxOnBoundary(Dune::FieldVector< Scalar, 2 > &, const Intersection &, const CellData &, const bool)
Function which calculates the flux entry at a boundary.
Definition: 1p/sequential/diffusion/cellcentered/pressure.hh:300
void getFlux(Dune::FieldVector< Scalar, 2 > &, const Intersection &, const CellData &, const bool)
Function which calculates the flux entry.
Definition: 1p/sequential/diffusion/cellcentered/pressure.hh:249
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:48
void assemble(bool first)
Function which assembles the system of equations to be solved.
Definition: sequential/cellcentered/pressure.hh:400
void getStorage(EntryType &entry, const Element &element, const CellData &cellData, const bool first)
Function which calculates the storage entry.
void initialize()
Initialize pressure model.
Definition: sequential/cellcentered/pressure.hh:212
PressureSolution & pressure()
Returns the vector containing the pressure solution.
Definition: sequential/cellcentered/pressure.hh:119
@ rhs
index for the right hand side entry
Definition: sequential/cellcentered/pressure.hh:87
@ matrix
index for the global matrix entry
Definition: sequential/cellcentered/pressure.hh:88
void solve()
Solves the global system of equations to get the spatial distribution of the pressure.
Definition: sequential/cellcentered/pressure.hh:526
void update()
Pressure update.
Definition: sequential/cellcentered/pressure.hh:227
Defines the properties required for the single phase sequential model.
Finite Volume Diffusion Model.