3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
staggered.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 *****************************************************************************/
27#ifndef DUMUX_DISCRETIZATION_STAGGERD_HH
28#define DUMUX_DISCRETIZATION_STAGGERD_HH
29
32
36
39
47
48#include <dune/istl/multitypeblockvector.hh>
49#include <dune/istl/multitypeblockmatrix.hh>
51
52namespace Dumux {
53
54// forward declarations
55class CCElementBoundaryTypes;
56
57namespace Properties
58{
60// Create new type tags
61namespace TTag {
62struct StaggeredModel { using InheritsFrom = std::tuple<FiniteVolumeModel>; };
63} // end namespace TTag
64
66template<class TypeTag>
67struct GridFaceVariables<TypeTag, TTag::StaggeredModel>
68{
69private:
72 static constexpr auto enableCache = getPropValue<TypeTag, Properties::EnableGridFaceVariablesCache>();
73public:
75};
76
78template<class TypeTag>
79struct EnableGridFaceVariablesCache<TypeTag, TTag::StaggeredModel> { static constexpr bool value = true; };
80
82template<class TypeTag>
83struct GridFluxVariablesCache<TypeTag, TTag::StaggeredModel>
84{
85private:
89 static constexpr auto enableCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
90 static constexpr auto upwindSchemeOrder = getPropValue<TypeTag, Properties::UpwindSchemeOrder>();
91public:
93};
94
96template<class TypeTag>
97struct StaggeredFaceSolution<TypeTag, TTag::StaggeredModel>
98{
99private:
101public:
103};
104
106template<class TypeTag>
107struct GridVariables<TypeTag, TTag::StaggeredModel>
108{
109private:
114public:
116};
117
119template<class TypeTag>
120struct ElementBoundaryTypes<TypeTag, TTag::StaggeredModel> { using type = CCElementBoundaryTypes; };
121
123template<class TypeTag>
124struct BaseLocalResidual<TypeTag, TTag::StaggeredModel> { using type = StaggeredLocalResidual<TypeTag>; };
125
127template<class TypeTag>
128struct CellCenterPrimaryVariables<TypeTag, TTag::StaggeredModel>
129{
130 using type = Dune::FieldVector<GetPropType<TypeTag, Properties::Scalar>,
131 getPropValue<TypeTag, Properties::NumEqCellCenter>()>;
132};
133
135template<class TypeTag>
136struct FacePrimaryVariables<TypeTag, TTag::StaggeredModel>
137{
138 using type = Dune::FieldVector<GetPropType<TypeTag, Properties::Scalar>,
139 getPropValue<TypeTag, Properties::NumEqFace>()>;
140};
141
142DUNE_NO_DEPRECATED_BEGIN
144template<class TypeTag>
146DUNE_NO_DEPRECATED_END
147// TODO: bundle SolutionVector, JacobianMatrix
148// in LinearAlgebra traits
149
151template<class TypeTag>
152struct CellCenterSolutionVector<TypeTag, TTag::StaggeredModel>
153{ using type = Dune::BlockVector<GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>>; };
154
156template<class TypeTag>
157struct FaceSolutionVector<TypeTag, TTag::StaggeredModel>
158{ using type = Dune::BlockVector<GetPropType<TypeTag, Properties::FacePrimaryVariables>>; };
159
161template<class TypeTag>
162struct SolutionVector<TypeTag, TTag::StaggeredModel>
163{
164private:
166 using CellCenterSolutionVector = GetPropType<TypeTag, Properties::CellCenterSolutionVector>;
167public:
168 using type = Dune::MultiTypeBlockVector<FaceSolutionVector, CellCenterSolutionVector>;
169};
170
172template<class TypeTag>
173struct JacobianMatrix<TypeTag, TTag::StaggeredModel>
174{
175private:
177
178 static constexpr auto numEqCellCenter = getPropValue<TypeTag, Properties::NumEqCellCenter>();
179 static constexpr auto numEqFace = getPropValue<TypeTag, Properties::NumEqFace>();
180
181public:
182 // the sub-blocks
183 using MatrixLittleBlockCCToCC = typename Dune::FieldMatrix<Scalar, numEqCellCenter, numEqCellCenter>;
184 using MatrixLittleBlockCCToFace = typename Dune::FieldMatrix<Scalar, numEqCellCenter, numEqFace>;
185
186 using MatrixLittleBlockFaceToFace = typename Dune::FieldMatrix<Scalar, numEqFace, numEqFace>;
187 using MatrixLittleBlockFaceToCC = typename Dune::FieldMatrix<Scalar, numEqFace, numEqCellCenter>;
188
189 // the BCRS matrices of the subproblems as big blocks
192
195
196 // the row types
197 using RowFace = typename Dune::MultiTypeBlockVector<MatrixBlockFaceToFace, MatrixBlockFaceToCC>;
198 using RowCellCenter = typename Dune::MultiTypeBlockVector<MatrixBlockCCToFace, MatrixBlockCCToCC>;
199
200 // the jacobian matrix
202};
203
204} // namespace Properties
205
206namespace Detail {
207
208template<class Problem>
210{
211private:
212 using GG = std::decay_t<decltype(std::declval<Problem>().gridGeometry())>;
213 using Element = typename GG::GridView::template Codim<0>::Entity;
214 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
215public:
216 using GridGeometry = GG;
217 // BoundaryTypes is whatever the problem returns from boundaryTypes(element, scvf)
218 using BoundaryTypes = std::decay_t<decltype(std::declval<Problem>().boundaryTypes(std::declval<Element>(), std::declval<SubControlVolumeFace>()))>;
219};
220
221} // end namespace Detail
222
223} // namespace Dumux
224
225#endif
Calculates the element-wise residual for the staggered FV scheme.
Helpers for deprecation.
The available discretization methods in Dumux.
Declares properties required for finite-volume models models.
Classes related to flux variables caching.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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
Calculates the element-wise residual for the staggered FV scheme.
Definition: staggeredlocalresidual.hh:41
Class to specify the type of a boundary.
Definition: common/boundarytypes.hh:38
Definition: common/pdesolver.hh:37
The type of the base class of the local residual (specific to a discretization scheme)
Definition: common/properties.hh:65
Type of the global jacobian matrix.
Definition: common/properties.hh:67
Vector containing all primary variable vector of the grid.
Definition: common/properties.hh:69
Stores the boundary types of a single degree of freedom.
Definition: common/properties.hh:72
Stores the boundary types on an element.
Definition: common/properties.hh:98
The global vector of flux variable containers.
Definition: common/properties.hh:118
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:122
The solution vector type for cell-centered dofs.
Definition: common/properties.hh:213
The solution vector type for face dofs.
Definition: common/properties.hh:215
Global vector containing face-related data.
Definition: common/properties.hh:217
The primary variables container type for cell-centered dofs.
Definition: common/properties.hh:219
The primary variables container type for face dofs.
Definition: common/properties.hh:221
A vector containing the solution for a face (similar to ElementSolution)
Definition: common/properties.hh:233
Switch on/off caching of face variables.
Definition: common/properties.hh:235
Definition: matrix.hh:32
Definition: common/typetraits/problem.hh:35
Boundary types gathered on an element.
Definition: cellcentered/elementboundarytypes.hh:38
Definition: staggered.hh:62
std::tuple< FiniteVolumeModel > InheritsFrom
Definition: staggered.hh:62
Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, getPropValue< TypeTag, Properties::NumEqCellCenter >()> type
Definition: staggered.hh:131
Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, getPropValue< TypeTag, Properties::NumEqFace >()> type
Definition: staggered.hh:139
Dune::BlockVector< GetPropType< TypeTag, Properties::CellCenterPrimaryVariables > > type
Definition: staggered.hh:153
Dune::BlockVector< GetPropType< TypeTag, Properties::FacePrimaryVariables > > type
Definition: staggered.hh:158
Dune::MultiTypeBlockVector< FaceSolutionVector, CellCenterSolutionVector > type
Definition: staggered.hh:168
typename Dune::BCRSMatrix< MatrixLittleBlockCCToFace > MatrixBlockCCToFace
Definition: staggered.hh:191
typename Dune::MultiTypeBlockVector< MatrixBlockFaceToFace, MatrixBlockFaceToCC > RowFace
Definition: staggered.hh:197
typename Dune::BCRSMatrix< MatrixLittleBlockFaceToFace > MatrixBlockFaceToFace
Definition: staggered.hh:193
typename Dune::FieldMatrix< Scalar, numEqCellCenter, numEqFace > MatrixLittleBlockCCToFace
Definition: staggered.hh:184
typename Dune::MultiTypeBlockMatrix< RowFace, RowCellCenter > type
Definition: staggered.hh:201
typename Dune::BCRSMatrix< MatrixLittleBlockCCToCC > MatrixBlockCCToCC
Definition: staggered.hh:190
typename Dune::FieldMatrix< Scalar, numEqFace, numEqCellCenter > MatrixLittleBlockFaceToCC
Definition: staggered.hh:187
typename Dune::FieldMatrix< Scalar, numEqFace, numEqFace > MatrixLittleBlockFaceToFace
Definition: staggered.hh:186
typename Dune::MultiTypeBlockVector< MatrixBlockCCToFace, MatrixBlockCCToCC > RowCellCenter
Definition: staggered.hh:198
typename Dune::BCRSMatrix< MatrixLittleBlockFaceToCC > MatrixBlockFaceToCC
Definition: staggered.hh:194
typename Dune::FieldMatrix< Scalar, numEqCellCenter, numEqCellCenter > MatrixLittleBlockCCToCC
Definition: staggered.hh:183
std::decay_t< decltype(std::declval< Problem >().boundaryTypes(std::declval< Element >(), std::declval< SubControlVolumeFace >()))> BoundaryTypes
Definition: staggered.hh:218
The global face variables class for staggered models.
Definition: facesolution.hh:40
Face variables cache class for staggered models.
Definition: gridfacevariables.hh:59
Flux variables cache class for staggered models.
Definition: staggered/gridfluxvariablescache.hh:63
Class storing data associated to scvs and scvfs.
Definition: discretization/staggered/gridvariables.hh:195
Declares all properties used in Dumux.
Type traits for problem classes.
Boundary types gathered on an element.
Sub control volumes for cell-centered discretization schemes.