3.4
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>
50
51namespace Dumux {
52
53// forward declarations
54class CCElementBoundaryTypes;
55
56namespace Properties
57{
59// Create new type tags
60namespace TTag {
61struct StaggeredModel { using InheritsFrom = std::tuple<FiniteVolumeModel>; };
62} // end namespace TTag
63
65template<class TypeTag>
66struct GridFaceVariables<TypeTag, TTag::StaggeredModel>
67{
68private:
71 static constexpr auto enableCache = getPropValue<TypeTag, Properties::EnableGridFaceVariablesCache>();
72public:
74};
75
77template<class TypeTag>
78struct EnableGridFaceVariablesCache<TypeTag, TTag::StaggeredModel> { static constexpr bool value = true; };
79
81template<class TypeTag>
82struct GridFluxVariablesCache<TypeTag, TTag::StaggeredModel>
83{
84private:
88 static constexpr auto enableCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
89 static constexpr auto upwindSchemeOrder = getPropValue<TypeTag, Properties::UpwindSchemeOrder>();
90public:
92};
93
95template<class TypeTag>
96struct StaggeredFaceSolution<TypeTag, TTag::StaggeredModel>
97{
98private:
100public:
102};
103
105template<class TypeTag>
106struct GridVariables<TypeTag, TTag::StaggeredModel>
107{
108private:
113public:
115};
116
118template<class TypeTag>
119struct ElementBoundaryTypes<TypeTag, TTag::StaggeredModel> { using type = CCElementBoundaryTypes; };
120
122template<class TypeTag>
123struct BaseLocalResidual<TypeTag, TTag::StaggeredModel> { using type = StaggeredLocalResidual<TypeTag>; };
124
126template<class TypeTag>
127struct CellCenterPrimaryVariables<TypeTag, TTag::StaggeredModel>
128{
129 using type = Dune::FieldVector<GetPropType<TypeTag, Properties::Scalar>,
130 getPropValue<TypeTag, Properties::NumEqCellCenter>()>;
131};
132
134template<class TypeTag>
135struct FacePrimaryVariables<TypeTag, TTag::StaggeredModel>
136{
137 using type = Dune::FieldVector<GetPropType<TypeTag, Properties::Scalar>,
138 getPropValue<TypeTag, Properties::NumEqFace>()>;
139};
140
141// TODO: bundle SolutionVector, JacobianMatrix
142// in LinearAlgebra traits
143
145template<class TypeTag>
146struct CellCenterSolutionVector<TypeTag, TTag::StaggeredModel>
147{ using type = Dune::BlockVector<GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>>; };
148
150template<class TypeTag>
151struct FaceSolutionVector<TypeTag, TTag::StaggeredModel>
152{ using type = Dune::BlockVector<GetPropType<TypeTag, Properties::FacePrimaryVariables>>; };
153
155template<class TypeTag>
156struct SolutionVector<TypeTag, TTag::StaggeredModel>
157{
158private:
160 using CellCenterSolutionVector = GetPropType<TypeTag, Properties::CellCenterSolutionVector>;
161public:
163};
164
166template<class TypeTag>
167struct JacobianMatrix<TypeTag, TTag::StaggeredModel>
168{
169private:
171
172 static constexpr auto numEqCellCenter = getPropValue<TypeTag, Properties::NumEqCellCenter>();
173 static constexpr auto numEqFace = getPropValue<TypeTag, Properties::NumEqFace>();
174
175public:
176 // the sub-blocks
177 using MatrixLittleBlockCCToCC = typename Dune::FieldMatrix<Scalar, numEqCellCenter, numEqCellCenter>;
178 using MatrixLittleBlockCCToFace = typename Dune::FieldMatrix<Scalar, numEqCellCenter, numEqFace>;
179
180 using MatrixLittleBlockFaceToFace = typename Dune::FieldMatrix<Scalar, numEqFace, numEqFace>;
181 using MatrixLittleBlockFaceToCC = typename Dune::FieldMatrix<Scalar, numEqFace, numEqCellCenter>;
182
183 // the BCRS matrices of the subproblems as big blocks
186
189
190 // the row types
193
194 // the jacobian matrix
196};
197
198} // namespace Properties
199
200namespace Detail {
201
202template<class Problem>
204{
205private:
206 using GG = std::decay_t<decltype(std::declval<Problem>().gridGeometry())>;
207 using Element = typename GG::GridView::template Codim<0>::Entity;
208 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
209public:
210 using GridGeometry = GG;
211 // BoundaryTypes is whatever the problem returns from boundaryTypes(element, scvf)
212 using BoundaryTypes = std::decay_t<decltype(std::declval<Problem>().boundaryTypes(std::declval<Element>(), std::declval<SubControlVolumeFace>()))>;
213};
214
215} // end namespace Detail
216
217} // namespace Dumux
218
219#endif
Calculates the element-wise residual for the staggered FV scheme.
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
Definition: propertysystem.hh:150
Calculates the element-wise residual for the staggered FV scheme.
Definition: staggeredlocalresidual.hh:41
Definition: common/pdesolver.hh:38
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 on an element.
Definition: common/properties.hh:95
The global vector of flux variable containers.
Definition: common/properties.hh:115
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:119
The solution vector type for cell-centered dofs.
Definition: common/properties.hh:210
The solution vector type for face dofs.
Definition: common/properties.hh:212
Global vector containing face-related data.
Definition: common/properties.hh:214
The primary variables container type for cell-centered dofs.
Definition: common/properties.hh:216
The primary variables container type for face dofs.
Definition: common/properties.hh:218
A vector containing the solution for a face (similar to ElementSolution)
Definition: common/properties.hh:230
Switch on/off caching of face variables.
Definition: common/properties.hh:232
Definition: matrix.hh:32
Definition: common/typetraits/problem.hh:35
Definition: variablesbackend.hh:43
Boundary types gathered on an element.
Definition: cellcentered/elementboundarytypes.hh:38
Definition: staggered.hh:61
std::tuple< FiniteVolumeModel > InheritsFrom
Definition: staggered.hh:61
Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, getPropValue< TypeTag, Properties::NumEqCellCenter >()> type
Definition: staggered.hh:130
Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, getPropValue< TypeTag, Properties::NumEqFace >()> type
Definition: staggered.hh:138
Dune::BlockVector< GetPropType< TypeTag, Properties::CellCenterPrimaryVariables > > type
Definition: staggered.hh:147
Dune::BlockVector< GetPropType< TypeTag, Properties::FacePrimaryVariables > > type
Definition: staggered.hh:152
typename Dune::BCRSMatrix< MatrixLittleBlockCCToFace > MatrixBlockCCToFace
Definition: staggered.hh:185
typename Dune::MultiTypeBlockVector< MatrixBlockFaceToFace, MatrixBlockFaceToCC > RowFace
Definition: staggered.hh:191
typename Dune::BCRSMatrix< MatrixLittleBlockFaceToFace > MatrixBlockFaceToFace
Definition: staggered.hh:187
typename Dune::FieldMatrix< Scalar, numEqCellCenter, numEqFace > MatrixLittleBlockCCToFace
Definition: staggered.hh:178
typename Dune::MultiTypeBlockMatrix< RowFace, RowCellCenter > type
Definition: staggered.hh:195
typename Dune::BCRSMatrix< MatrixLittleBlockCCToCC > MatrixBlockCCToCC
Definition: staggered.hh:184
typename Dune::FieldMatrix< Scalar, numEqFace, numEqCellCenter > MatrixLittleBlockFaceToCC
Definition: staggered.hh:181
typename Dune::FieldMatrix< Scalar, numEqFace, numEqFace > MatrixLittleBlockFaceToFace
Definition: staggered.hh:180
typename Dune::MultiTypeBlockVector< MatrixBlockCCToFace, MatrixBlockCCToCC > RowCellCenter
Definition: staggered.hh:192
typename Dune::BCRSMatrix< MatrixLittleBlockFaceToCC > MatrixBlockFaceToCC
Definition: staggered.hh:188
typename Dune::FieldMatrix< Scalar, numEqCellCenter, numEqCellCenter > MatrixLittleBlockCCToCC
Definition: staggered.hh:177
std::decay_t< decltype(std::declval< Problem >().boundaryTypes(std::declval< Element >(), std::declval< SubControlVolumeFace >()))> BoundaryTypes
Definition: staggered.hh:212
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: discretization/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.