3.5-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
33
34namespace Dumux {
35
46template<class TypeTag>
48{
50 using Implementation = GetPropType<TypeTag, Properties::Problem>;
52 using Element = typename GridView::template Codim<0>::Entity;
53
55 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
56 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
57 using GridFaceVariables = typename GridVariables::GridFaceVariables;
58 using ElementFaceVariables = typename GridFaceVariables::LocalView;
59
62 using FVElementGeometry = typename GridGeometry::LocalView;
63 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
64 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
66
67 using CoordScalar = typename GridView::ctype;
68 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
69
70 static constexpr auto cellCenterIdx = GridGeometry::cellCenterIdx();
71 static constexpr auto faceIdx = GridGeometry::faceIdx();
72
73 static constexpr auto numEqCellCenter = getPropValue<TypeTag, Properties::NumEqCellCenter>();
74 static constexpr auto numEqFace = getPropValue<TypeTag, Properties::NumEqFace>();
75
76public:
82 StaggeredFVProblem(std::shared_ptr<const GridGeometry> gridGeometry,
83 const std::string& paramGroup = "")
85 { }
86
95 bool isDirichletCell(const Element& element,
96 const FVElementGeometry& fvGeometry,
97 const SubControlVolume& scv,
98 int pvIdx) const
99 { return false; }
100
119 template<class ElementVolumeVariables, class ElementFaceVariables, class Entity>
120 NumEqVector source(const Element &element,
121 const FVElementGeometry& fvGeometry,
122 const ElementVolumeVariables& elemVolVars,
123 const ElementFaceVariables& elementFaceVars,
124 const Entity &e) const
125 {
126 // forward to solution independent, fully-implicit specific interface
127 return this->asImp_().sourceAtPos(e.center());
128 }
129
145 NumEqVector neumann(const Element& element,
146 const FVElementGeometry& fvGeometry,
147 const ElementVolumeVariables& elemVolVars,
148 const ElementFaceVariables& elemFaceVars,
149 const SubControlVolumeFace& scvf) const
150 {
151 // forward it to the interface with only the global position
152 return this->asImp_().neumannAtPos(scvf.ipGlobal());
153 }
154
161 template<class Entity>
162 PrimaryVariables initial(const Entity& entity) const
163 {
164 return this->asImp_().initialAtPos(entity.center());
165 }
166
171 template<class SolutionVector>
172 void applyInitialSolution(SolutionVector& sol) const
173 {
174 sol[cellCenterIdx].resize(this->gridGeometry().numCellCenterDofs());
175 sol[faceIdx].resize(this->gridGeometry().numFaceDofs());
176
177 auto fvGeometry = localView(this->gridGeometry());
178 for (const auto& element : elements(this->gridGeometry().gridView()))
179 {
180 fvGeometry.bindElement(element);
181
182 // loop over sub control volumes
183 for (auto&& scv : scvs(fvGeometry))
184 {
185 // let the problem do the dirty work of nailing down
186 // the initial solution.
187 auto initPriVars = this->asImp_().initial(scv);
188 this->asImp_().applyInitialCellCenterSolution(sol, scv, initPriVars);
189 }
190
191 // loop over faces
192 for(auto&& scvf : scvfs(fvGeometry))
193 {
194 auto initPriVars = this->asImp_().initial(scvf);
195 this->asImp_().applyInitialFaceSolution(sol, scvf, initPriVars);
196 }
197 }
198 }
199
200
202 template<class SolutionVector>
203 void applyInitialCellCenterSolution(SolutionVector& sol,
204 const SubControlVolume& scv,
205 const PrimaryVariables& initSol) const
206 {
207 // while the container within the actual solution vector holds numEqCellCenter
208 // elements, we need to specify an offset to get the correct entry of the initial solution
209 static constexpr auto offset = PrimaryVariables::dimension - numEqCellCenter;
210
211 for(int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
212 sol[cellCenterIdx][scv.dofIndex()][pvIdx] = initSol[pvIdx + offset];
213 }
214
216 template<class SolutionVector>
217 void applyInitialFaceSolution(SolutionVector& sol,
218 const SubControlVolumeFace& scvf,
219 const PrimaryVariables& initSol) const
220 {
221 for(int pvIdx = 0; pvIdx < numEqFace; ++pvIdx)
222 sol[faceIdx][scvf.dofIndex()][pvIdx] = initSol[pvIdx];
223 }
224
225};
226
227} // end namespace Dumux
228
229#endif
A helper to deduce a vector with the same size as numbers of equations.
Base class for all finite volume problems that are parameterized.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:46
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:55
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
Implementation & asImp_()
Returns the implementation of the problem (i.e. static polymorphism)
Definition: common/fvproblem.hh:591
Base class for all finite-volume problems using spatial parameters.
Definition: fvproblemwithspatialparams.hh:41
Base class for all staggered finite-volume problems.
Definition: staggeredfvproblem.hh:48
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:145
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:95
StaggeredFVProblem(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
Constructor.
Definition: staggeredfvproblem.hh:82
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:162
void applyInitialFaceSolution(SolutionVector &sol, const SubControlVolumeFace &scvf, const PrimaryVariables &initSol) const
Applys the initial face solution.
Definition: staggeredfvproblem.hh:217
void applyInitialCellCenterSolution(SolutionVector &sol, const SubControlVolume &scv, const PrimaryVariables &initSol) const
Applys the initial cell center solution.
Definition: staggeredfvproblem.hh:203
void applyInitialSolution(SolutionVector &sol) const
Applies the initial solution for all degrees of freedom of the grid.
Definition: staggeredfvproblem.hh:172
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:120
Declares all properties used in Dumux.