3.6-git
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:
87 using FluxVariablesCache = GetPropTypeOr<TypeTag,
89 >;
90 using FluxVariablesCacheFiller = GetPropTypeOr<TypeTag,
92 >;
93 static constexpr auto enableCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
94 static constexpr auto upwindSchemeOrder = getPropValue<TypeTag, Properties::UpwindSchemeOrder>();
95public:
97};
98
100template<class TypeTag>
101struct StaggeredFaceSolution<TypeTag, TTag::StaggeredModel>
102{
103private:
105public:
107};
108
110template<class TypeTag>
111struct GridVariables<TypeTag, TTag::StaggeredModel>
112{
113private:
118public:
120};
121
123template<class TypeTag>
124struct ElementBoundaryTypes<TypeTag, TTag::StaggeredModel> { using type = CCElementBoundaryTypes; };
125
127template<class TypeTag>
128struct BaseLocalResidual<TypeTag, TTag::StaggeredModel> { using type = StaggeredLocalResidual<TypeTag>; };
129
131template<class TypeTag>
132struct CellCenterPrimaryVariables<TypeTag, TTag::StaggeredModel>
133{
134 using type = Dune::FieldVector<GetPropType<TypeTag, Properties::Scalar>,
135 getPropValue<TypeTag, Properties::NumEqCellCenter>()>;
136};
137
139template<class TypeTag>
140struct FacePrimaryVariables<TypeTag, TTag::StaggeredModel>
141{
142 using type = Dune::FieldVector<GetPropType<TypeTag, Properties::Scalar>,
143 getPropValue<TypeTag, Properties::NumEqFace>()>;
144};
145
146// TODO: bundle SolutionVector, JacobianMatrix
147// in LinearAlgebra traits
148
150template<class TypeTag>
151struct CellCenterSolutionVector<TypeTag, TTag::StaggeredModel>
152{ using type = Dune::BlockVector<GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>>; };
153
155template<class TypeTag>
156struct FaceSolutionVector<TypeTag, TTag::StaggeredModel>
157{ using type = Dune::BlockVector<GetPropType<TypeTag, Properties::FacePrimaryVariables>>; };
158
160template<class TypeTag>
161struct SolutionVector<TypeTag, TTag::StaggeredModel>
162{
163private:
165 using CellCenterSolutionVector = GetPropType<TypeTag, Properties::CellCenterSolutionVector>;
166public:
168};
169
171template<class TypeTag>
172struct JacobianMatrix<TypeTag, TTag::StaggeredModel>
173{
174private:
176
177 static constexpr auto numEqCellCenter = getPropValue<TypeTag, Properties::NumEqCellCenter>();
178 static constexpr auto numEqFace = getPropValue<TypeTag, Properties::NumEqFace>();
179
180public:
181 // the sub-blocks
182 using MatrixLittleBlockCCToCC = typename Dune::FieldMatrix<Scalar, numEqCellCenter, numEqCellCenter>;
183 using MatrixLittleBlockCCToFace = typename Dune::FieldMatrix<Scalar, numEqCellCenter, numEqFace>;
184
185 using MatrixLittleBlockFaceToFace = typename Dune::FieldMatrix<Scalar, numEqFace, numEqFace>;
186 using MatrixLittleBlockFaceToCC = typename Dune::FieldMatrix<Scalar, numEqFace, numEqCellCenter>;
187
188 // the BCRS matrices of the subproblems as big blocks
191
194
195 // the row types
198
199 // the jacobian matrix
201};
202
203} // namespace Properties
204
205namespace Detail {
206
207template<class Problem>
208struct ProblemTraits<Problem, DiscretizationMethods::Staggered>
209{
210private:
211 using GG = std::decay_t<decltype(std::declval<Problem>().gridGeometry())>;
212 using Element = typename GG::GridView::template Codim<0>::Entity;
213 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
214public:
215 using GridGeometry = GG;
216 // BoundaryTypes is whatever the problem returns from boundaryTypes(element, scvf)
217 using BoundaryTypes = std::decay_t<decltype(std::declval<Problem>().boundaryTypes(std::declval<Element>(), std::declval<SubControlVolumeFace>()))>;
218};
219
220} // end namespace Detail
221
222} // namespace Dumux
223
224#endif
Calculates the element-wise residual for the staggered FV scheme.
Classes related to flux variables caching.
The available discretization methods in Dumux.
Declares properties required for finite-volume models models.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
typename GetPropOr< TypeTag, Property, T >::type GetPropTypeOr
get the type alias defined in the property or the type T if the property is undefined
Definition: propertysystem.hh:184
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:63
Type of the global jacobian matrix.
Definition: common/properties.hh:65
Vector containing all primary variable vector of the grid.
Definition: common/properties.hh:67
Stores the boundary types on an element.
Definition: common/properties.hh:97
Stores data associated with flux vars.
Definition: common/properties.hh:113
The engine behind the global flux cache (how to fill caches for the stencil)
Definition: common/properties.hh:115
The global vector of flux variable containers.
Definition: common/properties.hh:117
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:121
The solution vector type for cell-centered dofs.
Definition: common/properties.hh:218
The solution vector type for face dofs.
Definition: common/properties.hh:220
Global vector containing face-related data.
Definition: common/properties.hh:222
The primary variables container type for cell-centered dofs.
Definition: common/properties.hh:224
The primary variables container type for face dofs.
Definition: common/properties.hh:226
A vector containing the solution for a face (similar to ElementSolution)
Definition: common/properties.hh:238
Switch on/off caching of face variables.
Definition: common/properties.hh:240
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:135
Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, getPropValue< TypeTag, Properties::NumEqFace >()> type
Definition: staggered.hh:143
Dune::BlockVector< GetPropType< TypeTag, Properties::CellCenterPrimaryVariables > > type
Definition: staggered.hh:152
Dune::BlockVector< GetPropType< TypeTag, Properties::FacePrimaryVariables > > type
Definition: staggered.hh:157
typename Dune::BCRSMatrix< MatrixLittleBlockCCToFace > MatrixBlockCCToFace
Definition: staggered.hh:190
typename Dune::MultiTypeBlockVector< MatrixBlockFaceToFace, MatrixBlockFaceToCC > RowFace
Definition: staggered.hh:196
typename Dune::BCRSMatrix< MatrixLittleBlockFaceToFace > MatrixBlockFaceToFace
Definition: staggered.hh:192
typename Dune::FieldMatrix< Scalar, numEqCellCenter, numEqFace > MatrixLittleBlockCCToFace
Definition: staggered.hh:183
typename Dune::MultiTypeBlockMatrix< RowFace, RowCellCenter > type
Definition: staggered.hh:200
typename Dune::BCRSMatrix< MatrixLittleBlockCCToCC > MatrixBlockCCToCC
Definition: staggered.hh:189
typename Dune::FieldMatrix< Scalar, numEqFace, numEqCellCenter > MatrixLittleBlockFaceToCC
Definition: staggered.hh:186
typename Dune::FieldMatrix< Scalar, numEqFace, numEqFace > MatrixLittleBlockFaceToFace
Definition: staggered.hh:185
typename Dune::MultiTypeBlockVector< MatrixBlockCCToFace, MatrixBlockCCToCC > RowCellCenter
Definition: staggered.hh:197
typename Dune::BCRSMatrix< MatrixLittleBlockFaceToCC > MatrixBlockFaceToCC
Definition: staggered.hh:193
typename Dune::FieldMatrix< Scalar, numEqCellCenter, numEqCellCenter > MatrixLittleBlockCCToCC
Definition: staggered.hh:182
std::decay_t< decltype(std::declval< Problem >().boundaryTypes(std::declval< Element >(), std::declval< SubControlVolumeFace >()))> BoundaryTypes
Definition: staggered.hh:217
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
The empty filler class corresponding to EmptyCache.
Definition: fluxvariablescaching.hh:32
An empty flux variables cache.
Definition: fluxvariablescaching.hh:47
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.