3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
staggeredfvproblem.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_STAGGERD_FV_PROBLEM_HH
25#define DUMUX_STAGGERD_FV_PROBLEM_HH
26
27#include <dune/common/rangeutilities.hh>
28
31
32namespace Dumux {
33
43template<class TypeTag>
44class StaggeredFVProblem : public FVProblem<TypeTag>
45{
47 using Implementation = GetPropType<TypeTag, Properties::Problem>;
49 using Element = typename GridView::template Codim<0>::Entity;
50
52 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
53 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
54 using GridFaceVariables = typename GridVariables::GridFaceVariables;
55 using ElementFaceVariables = typename GridFaceVariables::LocalView;
56
59 using FVElementGeometry = typename GridGeometry::LocalView;
60 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
61 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
63
64 using CoordScalar = typename GridView::ctype;
65 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
66
67 static constexpr auto cellCenterIdx = GridGeometry::cellCenterIdx();
68 static constexpr auto faceIdx = GridGeometry::faceIdx();
69
70 static constexpr auto numEqCellCenter = getPropValue<TypeTag, Properties::NumEqCellCenter>();
71 static constexpr auto numEqFace = getPropValue<TypeTag, Properties::NumEqFace>();
72
73public:
79 StaggeredFVProblem(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
81 { }
82
91 bool isDirichletCell(const Element& element,
92 const FVElementGeometry& fvGeometry,
93 const SubControlVolume& scv,
94 int pvIdx) const
95 { return false; }
96
115 template<class ElementVolumeVariables, class ElementFaceVariables, class Entity>
116 NumEqVector source(const Element &element,
117 const FVElementGeometry& fvGeometry,
118 const ElementVolumeVariables& elemVolVars,
119 const ElementFaceVariables& elementFaceVars,
120 const Entity &e) const
121 {
122 // forward to solution independent, fully-implicit specific interface
123 return asImp_().sourceAtPos(e.center());
124 }
125
141 NumEqVector neumann(const Element& element,
142 const FVElementGeometry& fvGeometry,
143 const ElementVolumeVariables& elemVolVars,
144 const ElementFaceVariables& elemFaceVars,
145 const SubControlVolumeFace& scvf) const
146 {
147 // forward it to the interface with only the global position
148 return asImp_().neumannAtPos(scvf.ipGlobal());
149 }
150
157 template<class Entity>
158 PrimaryVariables initial(const Entity& entity) const
159 {
160 return asImp_().initialAtPos(entity.center());
161 }
162
167 template<class SolutionVector>
168 void applyInitialSolution(SolutionVector& sol) const
169 {
170 sol[cellCenterIdx].resize(this->gridGeometry().numCellCenterDofs());
171 sol[faceIdx].resize(this->gridGeometry().numFaceDofs());
172
173 for (const auto& element : elements(this->gridGeometry().gridView()))
174 {
175 auto fvGeometry = localView(this->gridGeometry());
176 fvGeometry.bindElement(element);
177
178 // loop over sub control volumes
179 for (auto&& scv : scvs(fvGeometry))
180 {
181 // let the problem do the dirty work of nailing down
182 // the initial solution.
183 auto initPriVars = asImp_().initial(scv);
184 asImp_().applyInitialCellCenterSolution(sol, scv, initPriVars);
185 }
186
187 // loop over faces
188 for(auto&& scvf : scvfs(fvGeometry))
189 {
190 auto initPriVars = asImp_().initial(scvf);
191 asImp_().applyInitialFaceSolution(sol, scvf, initPriVars);
192 }
193 }
194 }
195
196
198 template<class SolutionVector>
199 void applyInitialCellCenterSolution(SolutionVector& sol,
200 const SubControlVolume& scv,
201 const PrimaryVariables& initSol) const
202 {
203 // while the container within the actual solution vector holds numEqCellCenter
204 // elements, we need to specify an offset to get the correct entry of the initial solution
205 static constexpr auto offset = PrimaryVariables::dimension - numEqCellCenter;
206
207 for(int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
208 sol[cellCenterIdx][scv.dofIndex()][pvIdx] = initSol[pvIdx + offset];
209 }
210
212 template<class SolutionVector>
213 void applyInitialFaceSolution(SolutionVector& sol,
214 const SubControlVolumeFace& scvf,
215 const PrimaryVariables& initSol) const
216 {
217 for(int pvIdx = 0; pvIdx < numEqFace; ++pvIdx)
218 sol[faceIdx][scvf.dofIndex()][pvIdx] = initSol[pvIdx];
219 }
220
221protected:
223 Implementation &asImp_()
224 { return *static_cast<Implementation *>(this); }
225
227 const Implementation &asImp_() const
228 { return *static_cast<const Implementation *>(this); }
229};
230
231} // end namespace Dumux
232
233#endif
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
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:592
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: common/fvproblem.hh:588
Base class for all staggered finite-volume problems.
Definition: staggeredfvproblem.hh:45
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:141
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:91
const Implementation & asImp_() const
Returns the implementation of the problem (i.e. static polymorphism)
Definition: staggeredfvproblem.hh:227
StaggeredFVProblem(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
Constructor.
Definition: staggeredfvproblem.hh:79
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:158
Implementation & asImp_()
Returns the implementation of the problem (i.e. static polymorphism)
Definition: staggeredfvproblem.hh:223
void applyInitialFaceSolution(SolutionVector &sol, const SubControlVolumeFace &scvf, const PrimaryVariables &initSol) const
Applys the initial face solution.
Definition: staggeredfvproblem.hh:213
void applyInitialCellCenterSolution(SolutionVector &sol, const SubControlVolume &scv, const PrimaryVariables &initSol) const
Applys the initial cell center solution.
Definition: staggeredfvproblem.hh:199
void applyInitialSolution(SolutionVector &sol) const
Applies the initial solution for all degrees of freedom of the grid.
Definition: staggeredfvproblem.hh:168
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:116
Base class for all finite volume problems.
Declares all properties used in Dumux.