3.2-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 *****************************************************************************/
25#ifndef DUMUX_STAGGERD_FV_PROBLEM_HH
26#define DUMUX_STAGGERD_FV_PROBLEM_HH
27
28#include <dune/common/rangeutilities.hh>
29
32
33namespace Dumux {
34
45template<class TypeTag>
46class StaggeredFVProblem : public FVProblem<TypeTag>
47{
49 using Implementation = GetPropType<TypeTag, Properties::Problem>;
51 using Element = typename GridView::template Codim<0>::Entity;
52
54 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
55 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
56 using GridFaceVariables = typename GridVariables::GridFaceVariables;
57 using ElementFaceVariables = typename GridFaceVariables::LocalView;
58
61 using FVElementGeometry = typename GridGeometry::LocalView;
62 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
63 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
65
66 using CoordScalar = typename GridView::ctype;
67 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
68
69 static constexpr auto cellCenterIdx = GridGeometry::cellCenterIdx();
70 static constexpr auto faceIdx = GridGeometry::faceIdx();
71
72 static constexpr auto numEqCellCenter = getPropValue<TypeTag, Properties::NumEqCellCenter>();
73 static constexpr auto numEqFace = getPropValue<TypeTag, Properties::NumEqFace>();
74
75public:
81 StaggeredFVProblem(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
83 { }
84
93 bool isDirichletCell(const Element& element,
94 const FVElementGeometry& fvGeometry,
95 const SubControlVolume& scv,
96 int pvIdx) const
97 { return false; }
98
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
123 {
124 // forward to solution independent, fully-implicit specific interface
125 return asImp_().sourceAtPos(e.center());
126 }
127
143 NumEqVector neumann(const Element& element,
144 const FVElementGeometry& fvGeometry,
145 const ElementVolumeVariables& elemVolVars,
146 const ElementFaceVariables& elemFaceVars,
147 const SubControlVolumeFace& scvf) const
148 {
149 // forward it to the interface with only the global position
150 return asImp_().neumannAtPos(scvf.ipGlobal());
151 }
152
159 template<class Entity>
160 PrimaryVariables initial(const Entity& entity) const
161 {
162 return asImp_().initialAtPos(entity.center());
163 }
164
169 template<class SolutionVector>
170 void applyInitialSolution(SolutionVector& sol) const
171 {
172 sol[cellCenterIdx].resize(this->gridGeometry().numCellCenterDofs());
173 sol[faceIdx].resize(this->gridGeometry().numFaceDofs());
174
175 for (const auto& element : elements(this->gridGeometry().gridView()))
176 {
177 auto fvGeometry = localView(this->gridGeometry());
178 fvGeometry.bindElement(element);
179
180 // loop over sub control volumes
181 for (auto&& scv : scvs(fvGeometry))
182 {
183 // let the problem do the dirty work of nailing down
184 // the initial solution.
185 auto initPriVars = asImp_().initial(scv);
186 asImp_().applyInitialCellCenterSolution(sol, scv, initPriVars);
187 }
188
189 // loop over faces
190 for(auto&& scvf : scvfs(fvGeometry))
191 {
192 auto initPriVars = asImp_().initial(scvf);
193 asImp_().applyInitialFaceSolution(sol, scvf, initPriVars);
194 }
195 }
196 }
197
198
200 template<class SolutionVector>
201 void applyInitialCellCenterSolution(SolutionVector& sol,
202 const SubControlVolume& scv,
203 const PrimaryVariables& initSol) const
204 {
205 // while the container within the actual solution vector holds numEqCellCenter
206 // elements, we need to specify an offset to get the correct entry of the initial solution
207 static constexpr auto offset = PrimaryVariables::dimension - numEqCellCenter;
208
209 for(int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
210 sol[cellCenterIdx][scv.dofIndex()][pvIdx] = initSol[pvIdx + offset];
211 }
212
214 template<class SolutionVector>
215 void applyInitialFaceSolution(SolutionVector& sol,
216 const SubControlVolumeFace& scvf,
217 const PrimaryVariables& initSol) const
218 {
219 for(int pvIdx = 0; pvIdx < numEqFace; ++pvIdx)
220 sol[faceIdx][scvf.dofIndex()][pvIdx] = initSol[pvIdx];
221 }
222
223protected:
225 Implementation &asImp_()
226 { return *static_cast<Implementation *>(this); }
227
229 const Implementation &asImp_() const
230 { return *static_cast<const Implementation *>(this); }
231};
232
233} // end namespace Dumux
234
235#endif
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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: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 &paramGroup="")
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.