version 3.11-dev
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
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// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
15#ifndef DUMUX_DISCRETIZATION_STAGGERD_HH
16#define DUMUX_DISCRETIZATION_STAGGERD_HH
17
18#include <concepts>
19#include <type_traits>
20
23
28
31
39
40#include <dune/istl/multitypeblockvector.hh>
41#include <dune/istl/multitypeblockmatrix.hh>
42
43namespace Dumux {
44
45// forward declarations
46class CCElementBoundaryTypes;
47
48namespace Properties
49{
51// Create new type tags
52namespace TTag {
53struct StaggeredModel { using InheritsFrom = std::tuple<FiniteVolumeModel>; };
54} // end namespace TTag
55
57template<class TypeTag>
58struct GridFaceVariables<TypeTag, TTag::StaggeredModel>
59{
60private:
63 static constexpr auto enableCache = getPropValue<TypeTag, Properties::EnableGridFaceVariablesCache>();
64public:
66};
67
69template<class TypeTag>
70struct EnableGridFaceVariablesCache<TypeTag, TTag::StaggeredModel> { static constexpr bool value = true; };
71
73template<class TypeTag>
74struct GridFluxVariablesCache<TypeTag, TTag::StaggeredModel>
75{
76private:
79 using FluxVariablesCache = GetPropTypeOr<TypeTag,
80 Properties::FluxVariablesCache, FluxVariablesCaching::EmptyCache<Scalar>
81 >;
82 using FluxVariablesCacheFiller = GetPropTypeOr<TypeTag,
83 Properties::FluxVariablesCacheFiller, FluxVariablesCaching::EmptyCacheFiller
84 >;
85 static constexpr auto enableCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
86 static constexpr auto upwindSchemeOrder = getPropValue<TypeTag, Properties::UpwindSchemeOrder>();
87public:
89};
90
92template<class TypeTag>
93struct StaggeredFaceSolution<TypeTag, TTag::StaggeredModel>
94{
95private:
97public:
99};
100
102template<class TypeTag>
103struct GridVariables<TypeTag, TTag::StaggeredModel>
104{
105private:
110public:
112};
113
115template<class TypeTag>
116struct ElementBoundaryTypes<TypeTag, TTag::StaggeredModel> { using type = CCElementBoundaryTypes; };
117
119template<class TypeTag>
120struct CellCenterPrimaryVariables<TypeTag, TTag::StaggeredModel>
121{
122 using type = Dune::FieldVector<GetPropType<TypeTag, Properties::Scalar>,
123 getPropValue<TypeTag, Properties::NumEqCellCenter>()>;
124};
125
127template<class TypeTag>
128struct FacePrimaryVariables<TypeTag, TTag::StaggeredModel>
129{
130 using type = Dune::FieldVector<GetPropType<TypeTag, Properties::Scalar>,
131 getPropValue<TypeTag, Properties::NumEqFace>()>;
132};
133
134// TODO: bundle SolutionVector, JacobianMatrix
135// in LinearAlgebra traits
136
138template<class TypeTag>
139struct CellCenterSolutionVector<TypeTag, TTag::StaggeredModel>
140{ using type = Dune::BlockVector<GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>>; };
141
143template<class TypeTag>
144struct FaceSolutionVector<TypeTag, TTag::StaggeredModel>
145{ using type = Dune::BlockVector<GetPropType<TypeTag, Properties::FacePrimaryVariables>>; };
146
148template<class TypeTag>
149struct SolutionVector<TypeTag, TTag::StaggeredModel>
150{
151private:
153 using CellCenterSolutionVector = GetPropType<TypeTag, Properties::CellCenterSolutionVector>;
154public:
156};
157
159template<class TypeTag>
160struct JacobianMatrix<TypeTag, TTag::StaggeredModel>
161{
162private:
164
165 static constexpr auto numEqCellCenter = getPropValue<TypeTag, Properties::NumEqCellCenter>();
166 static constexpr auto numEqFace = getPropValue<TypeTag, Properties::NumEqFace>();
167
168public:
169 // the sub-blocks
170 using MatrixLittleBlockCCToCC = typename Dune::FieldMatrix<Scalar, numEqCellCenter, numEqCellCenter>;
171 using MatrixLittleBlockCCToFace = typename Dune::FieldMatrix<Scalar, numEqCellCenter, numEqFace>;
172
173 using MatrixLittleBlockFaceToFace = typename Dune::FieldMatrix<Scalar, numEqFace, numEqFace>;
174 using MatrixLittleBlockFaceToCC = typename Dune::FieldMatrix<Scalar, numEqFace, numEqCellCenter>;
175
176 // the BCRS matrices of the subproblems as big blocks
179
182
183 // the row types
186
187 // the jacobian matrix
189};
190
191} // namespace Properties
192
193namespace Detail {
194
195template<class Problem>
196struct ProblemTraits<Problem, DiscretizationMethods::Staggered>
197{
198private:
199 using GG = std::decay_t<decltype(std::declval<Problem>().gridGeometry())>;
200 using Element = typename GG::GridView::template Codim<0>::Entity;
201 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
202public:
203 using GridGeometry = GG;
204 // BoundaryTypes is whatever the problem returns from boundaryTypes(element, scvf)
205 using BoundaryTypes = std::decay_t<decltype(std::declval<Problem>().boundaryTypes(std::declval<Element>(), std::declval<SubControlVolumeFace>()))>;
206};
207
208template<class TypeTag>
209concept StaggeredModel = std::is_same_v<
210 typename GetPropType<TypeTag, Properties::GridGeometry>::DiscretizationMethod,
211 DiscretizationMethods::Staggered
212>;
213
214template<StaggeredModel TypeTag>
216{
218};
219
220} // end namespace Detail
221
222} // namespace Dumux
223
224#endif
Boundary types gathered on an element.
Boundary types gathered on an element.
Definition: cellcentered/elementboundarytypes.hh:26
The element-wise residual for control-volume finite element schemes.
Definition: cvfelocalresidual.hh:60
The global face variables class for staggered models.
Definition: facesolution.hh:28
Face variables cache class for staggered models.
Definition: gridfacevariables.hh:47
Flux variables cache class for staggered models.
Definition: discretization/staggered/gridfluxvariablescache.hh:51
Class storing data associated to scvs and scvfs.
Definition: discretization/staggered/gridvariables.hh:183
Calculates the element-wise residual for the staggered FV scheme.
Definition: staggeredlocalresidual.hh:29
Definition: matrix.hh:20
Definition: common/pdesolver.hh:26
Definition: variablesbackend.hh:34
Defines all properties used in Dumux.
Type traits for problem classes.
Definition: staggered.hh:209
The default local operator than can be specialized for each discretization scheme.
Sub control volumes for cell-centered discretization schemes.
Classes related to flux variables caching.
Declares properties required for finite-volume models models.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
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:303
The available discretization methods in Dumux.
Definition: adapt.hh:17
typename Detail::DiscretizationDefaultLocalOperator< TypeTag >::type DiscretizationDefaultLocalOperator
Definition: defaultlocaloperator.hh:27
Calculates the element-wise residual for the staggered FV scheme.
std::decay_t< decltype(std::declval< Problem >().boundaryTypes(std::declval< Element >(), std::declval< SubControlVolumeFace >()))> BoundaryTypes
Definition: staggered.hh:205
Definition: common/typetraits/problem.hh:23
The empty filler class corresponding to EmptyCache.
Definition: fluxvariablescaching.hh:20
An empty flux variables cache.
Definition: fluxvariablescaching.hh:35
Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, getPropValue< TypeTag, Properties::NumEqCellCenter >()> type
Definition: staggered.hh:123
Dune::BlockVector< GetPropType< TypeTag, Properties::CellCenterPrimaryVariables > > type
Definition: staggered.hh:140
Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, getPropValue< TypeTag, Properties::NumEqFace >()> type
Definition: staggered.hh:131
Dune::BlockVector< GetPropType< TypeTag, Properties::FacePrimaryVariables > > type
Definition: staggered.hh:145
typename Dune::BCRSMatrix< MatrixLittleBlockCCToFace > MatrixBlockCCToFace
Definition: staggered.hh:178
typename Dune::MultiTypeBlockVector< MatrixBlockFaceToFace, MatrixBlockFaceToCC > RowFace
Definition: staggered.hh:184
typename Dune::BCRSMatrix< MatrixLittleBlockFaceToFace > MatrixBlockFaceToFace
Definition: staggered.hh:180
typename Dune::FieldMatrix< Scalar, numEqCellCenter, numEqFace > MatrixLittleBlockCCToFace
Definition: staggered.hh:171
typename Dune::MultiTypeBlockMatrix< RowFace, RowCellCenter > type
Definition: staggered.hh:188
typename Dune::BCRSMatrix< MatrixLittleBlockCCToCC > MatrixBlockCCToCC
Definition: staggered.hh:177
typename Dune::FieldMatrix< Scalar, numEqFace, numEqCellCenter > MatrixLittleBlockFaceToCC
Definition: staggered.hh:174
typename Dune::FieldMatrix< Scalar, numEqFace, numEqFace > MatrixLittleBlockFaceToFace
Definition: staggered.hh:173
typename Dune::MultiTypeBlockVector< MatrixBlockCCToFace, MatrixBlockCCToCC > RowCellCenter
Definition: staggered.hh:185
typename Dune::BCRSMatrix< MatrixLittleBlockFaceToCC > MatrixBlockFaceToCC
Definition: staggered.hh:181
typename Dune::FieldMatrix< Scalar, numEqCellCenter, numEqCellCenter > MatrixLittleBlockCCToCC
Definition: staggered.hh:170
Definition: staggered.hh:53
std::tuple< FiniteVolumeModel > InheritsFrom
Definition: staggered.hh:53