25#ifndef DUMUX_STAGGERD_FV_PROBLEM_HH
26#define DUMUX_STAGGERD_FV_PROBLEM_HH
28#include <dune/common/rangeutilities.hh>
45template<
class TypeTag>
51 using Element =
typename GridView::template Codim<0>::Entity;
54 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
55 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
56 using GridFaceVariables =
typename GridVariables::GridFaceVariables;
57 using ElementFaceVariables =
typename GridFaceVariables::LocalView;
61 using FVElementGeometry =
typename GridGeometry::LocalView;
62 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
63 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
66 using CoordScalar =
typename GridView::ctype;
67 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
69 static constexpr auto cellCenterIdx = GridGeometry::cellCenterIdx();
70 static constexpr auto faceIdx = GridGeometry::faceIdx();
72 static constexpr auto numEqCellCenter = getPropValue<TypeTag, Properties::NumEqCellCenter>();
73 static constexpr auto numEqFace = getPropValue<TypeTag, Properties::NumEqFace>();
94 const FVElementGeometry& fvGeometry,
95 const SubControlVolume& scv,
117 template<
class ElementVolumeVariables,
class ElementFaceVariables,
class Entity>
118 NumEqVector
source(
const Element &element,
119 const FVElementGeometry& fvGeometry,
120 const ElementVolumeVariables& elemVolVars,
121 const ElementFaceVariables& elementFaceVars,
122 const Entity &e)
const
125 return asImp_().sourceAtPos(e.center());
144 const FVElementGeometry& fvGeometry,
145 const ElementVolumeVariables& elemVolVars,
146 const ElementFaceVariables& elemFaceVars,
147 const SubControlVolumeFace& scvf)
const
150 return asImp_().neumannAtPos(scvf.ipGlobal());
159 template<
class Entity>
160 PrimaryVariables
initial(
const Entity& entity)
const
162 return asImp_().initialAtPos(entity.center());
169 template<
class SolutionVector>
172 sol[cellCenterIdx].resize(this->
gridGeometry().numCellCenterDofs());
173 sol[faceIdx].resize(this->
gridGeometry().numFaceDofs());
175 for (
const auto& element : elements(this->
gridGeometry().gridView()))
178 fvGeometry.bindElement(element);
181 for (
auto&& scv : scvs(fvGeometry))
185 auto initPriVars =
asImp_().initial(scv);
186 asImp_().applyInitialCellCenterSolution(sol, scv, initPriVars);
190 for(
auto&& scvf : scvfs(fvGeometry))
192 auto initPriVars =
asImp_().initial(scvf);
193 asImp_().applyInitialFaceSolution(sol, scvf, initPriVars);
200 template<
class SolutionVector>
202 const SubControlVolume& scv,
203 const PrimaryVariables& initSol)
const
207 static constexpr auto offset = PrimaryVariables::dimension - numEqCellCenter;
209 for(
int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
210 sol[cellCenterIdx][scv.dofIndex()][pvIdx] = initSol[pvIdx + offset];
214 template<
class SolutionVector>
216 const SubControlVolumeFace& scvf,
217 const PrimaryVariables& initSol)
const
219 for(
int pvIdx = 0; pvIdx < numEqFace; ++pvIdx)
220 sol[faceIdx][scvf.dofIndex()][pvIdx] = initSol[pvIdx];
226 {
return *
static_cast<Implementation *
>(
this); }
230 {
return *
static_cast<const Implementation *
>(
this); }
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:50
const std::string & paramGroup() const
The parameter group in which to retrieve runtime parameters.
Definition: common/fvproblem.hh:586
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: common/fvproblem.hh:582
Base class for all staggered finite-volume problems.
Definition: staggeredfvproblem.hh:47
NumEqVector neumann(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf) const
Evaluate the boundary conditions for a neumann boundary segment.
Definition: staggeredfvproblem.hh:143
bool isDirichletCell(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolume &scv, int pvIdx) const
Returns whether a fixed Dirichlet value shall be used at a given cell.
Definition: staggeredfvproblem.hh:93
const Implementation & asImp_() const
Returns the implementation of the problem (i.e. static polymorphism)
Definition: staggeredfvproblem.hh:229
StaggeredFVProblem(std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="")
Constructor.
Definition: staggeredfvproblem.hh:81
PrimaryVariables initial(const Entity &entity) const
Evaluate the initial value for an element (for cell-centered primary variables) or face (for velociti...
Definition: staggeredfvproblem.hh:160
Implementation & asImp_()
Returns the implementation of the problem (i.e. static polymorphism)
Definition: staggeredfvproblem.hh:225
void applyInitialFaceSolution(SolutionVector &sol, const SubControlVolumeFace &scvf, const PrimaryVariables &initSol) const
Applys the initial face solution.
Definition: staggeredfvproblem.hh:215
void applyInitialCellCenterSolution(SolutionVector &sol, const SubControlVolume &scv, const PrimaryVariables &initSol) const
Applys the initial cell center solution.
Definition: staggeredfvproblem.hh:201
void applyInitialSolution(SolutionVector &sol) const
Applies the initial solution for all degrees of freedom of the grid.
Definition: staggeredfvproblem.hh:170
NumEqVector source(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elementFaceVars, const Entity &e) const
Evaluate the source term for all phases within a given sub-control-volume (-face).
Definition: staggeredfvproblem.hh:118
Base class for all finite volume problems.
Declares all properties used in Dumux.