version 3.8
staggeredtraits.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-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_STAGGERED_MULTIDOMAIN_TRAITS_HH
15#define DUMUX_STAGGERED_MULTIDOMAIN_TRAITS_HH
16
17#include <type_traits>
18#include <tuple>
19#include <utility>
20
21#include <dune/common/fmatrix.hh>
22#include <dune/common/indices.hh>
23
24#include <dune/istl/bcrsmatrix.hh>
25#include <dune/istl/multitypeblockvector.hh>
26#include <dune/istl/multitypeblockmatrix.hh>
27
30
32
33#include "traits.hh"
34
35namespace Dumux {
36namespace Detail {
37namespace Staggered {
38
40template<template<std::size_t> class SubDomainTypeTag, std::size_t i>
42{ using type = GetPropType<SubDomainTypeTag<i>, Properties::GridGeometry>; };
43
44template<template<std::size_t> class SubDomainTypeTag>
45struct SubDomainFVGridGeometryImpl<SubDomainTypeTag, 0>
46{ using type = typename GetPropType<SubDomainTypeTag<0>, Properties::GridGeometry>::FaceFVGridGeometryType; };
47
48template<template<std::size_t> class SubDomainTypeTag>
49struct SubDomainFVGridGeometryImpl<SubDomainTypeTag, 1>
50{ using type = typename GetPropType<SubDomainTypeTag<0>, Properties::GridGeometry>::CellCenterFVGridGeometryType; };
52
54template<template<std::size_t> class SubDomainTypeTag, std::size_t i>
56{ using type = GetPropType<SubDomainTypeTag<i>, Properties::GridVariables>; };
57
58template<template<std::size_t> class SubDomainTypeTag>
59struct SubDomainGridVariablesImpl<SubDomainTypeTag, 0>
60{ using type = typename GetPropType<SubDomainTypeTag<0>, Properties::GridVariables>::FaceGridVariablesType; };
61
62template<template<std::size_t> class SubDomainTypeTag>
63struct SubDomainGridVariablesImpl<SubDomainTypeTag, 1>
64{ using type = typename GetPropType<SubDomainTypeTag<0>, Properties::GridVariables>::CellCenterGridVariablesType; };
66
68template<template<std::size_t> class SubDomainTypeTag, std::size_t i>
70{ using type = GetPropType<SubDomainTypeTag<i>, Properties::PrimaryVariables>; };
71
72template<template<std::size_t> class SubDomainTypeTag>
73struct SubDomainPrimaryVariablesImpl<SubDomainTypeTag, 0>
74{ using type = GetPropType<SubDomainTypeTag<0>, Properties::FacePrimaryVariables>; };
75
76template<template<std::size_t> class SubDomainTypeTag>
77struct SubDomainPrimaryVariablesImpl<SubDomainTypeTag, 1>
78{ using type = GetPropType<SubDomainTypeTag<0>, Properties::CellCenterPrimaryVariables>; };
80
81template<class Scalar, int numEq>
83{
84 private:
85 using MatrixBlock = typename Dune::FieldMatrix<Scalar, numEq, numEq>;
86 public:
88};
89
91template<template<std::size_t> class SubDomainTypeTag, std::size_t i>
93{ using type = GetPropType<SubDomainTypeTag<i>, Properties::JacobianMatrix>; };
94
95template<template<std::size_t> class SubDomainTypeTag>
96struct SubDomainJacobianMatrixImpl<SubDomainTypeTag, 0>
97{ using type = typename JacobianTypeImpl<GetPropType<SubDomainTypeTag<0>, Properties::Scalar>,
98 getPropValue<SubDomainTypeTag<0>, Properties::NumEqFace>()>::type; };
99
100template<template<std::size_t> class SubDomainTypeTag>
101struct SubDomainJacobianMatrixImpl<SubDomainTypeTag, 1>
102{ using type = typename JacobianTypeImpl<GetPropType<SubDomainTypeTag<1>, Properties::Scalar>,
103 getPropValue<SubDomainTypeTag<0>, Properties::NumEqCellCenter>()>::type; };
105
107template<template<std::size_t> class SubDomainTypeTag, std::size_t i>
109{ using type = GetPropType<SubDomainTypeTag<i>, Properties::SolutionVector>; };
110
111template<template<std::size_t> class SubDomainTypeTag>
112struct SubDomainSolutionVectorImpl<SubDomainTypeTag, 0>
113{ using type = GetPropType<SubDomainTypeTag<0>, Properties::FaceSolutionVector>; };
114
115template<template<std::size_t> class SubDomainTypeTag>
116struct SubDomainSolutionVectorImpl<SubDomainTypeTag, 1>
117{ using type = GetPropType<SubDomainTypeTag<0>, Properties::CellCenterSolutionVector>; };
119
120} // end namespace Staggered
121} // end namespace Detail
122
123/*
124 * \ingroup MultiDomain
125 * \ingroup StaggeredDiscretization
126 * \brief A traits class every multidomain model has to provide
127 * \tparam SubDomainTypeTags the TypeTags of the sub domain problems
128 * \note should export the types
129 * \code
130 * //! the type tag of the sub domain problem with id
131 * template<std::size_t id>
132 * using SubDomainTypeTag = ...
133 *
134 * //! the index to access sub domain matrices and vectors
135 * //! to use with multitype matrices and vectors
136 * template<std::size_t id>
137 * using DomainIdx = ...
138 *
139 * //! the scalar type
140 * using Scalar = ...
141 *
142 * //! the solution vector type
143 * using SolutionVector = ...
144 *
145 * //! the residual vector type
146 * using ResidualVector = ...
147 *
148 * //! the jacobian type
149 * using JacobianMatrix = ...
150 * \endcode
151 */
152template<typename... SubDomainTypeTags>
154{
156 static constexpr std::size_t numSubDomains = sizeof...(SubDomainTypeTags);
157
158private:
159
161 template<std::size_t id>
162 using SubDomainTypeTag = typename std::tuple_element_t<id, std::tuple<SubDomainTypeTags...>>;
163
165 using Indices = std::make_index_sequence<numSubDomains>;
166
168 template<std::size_t id>
169 using SubDomainScalar = GetPropType<SubDomainTypeTag<id>, Properties::Scalar>;
170
171 template<std::size_t id>
173
174 template<std::size_t id>
176
177 template<std::size_t id>
178 using SubDomainResidualVector = typename Detail::NativeDuneVectorType<SubDomainSolutionVector<id>>::type;
179
180public:
181
182 /*
183 * \brief sub domain types
184 */
185 //\{
186
187 template<std::size_t id>
189 {
190 using Index = Dune::index_constant<id>;
191 using TypeTag = SubDomainTypeTag<id>;
192 using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
193 using LocalResidual = GetPropType<SubDomainTypeTag<id>, Properties::LocalResidual>;
199 };
200
201 //\}
202
203 /*
204 * \brief multi domain types
205 */
206 //\{
207
210
213
216
219
220 //\}
221
222 /*
223 * \brief helper aliases to construct derived tuple types
224 */
225 //\{
226
228 template<template<std::size_t> class T>
230
232 template<template<std::size_t> class T>
234
236 template<template<std::size_t> class T>
238
239 //\}
240};
241
242} //end namespace Dumux
243
244#endif
Definition: matrix.hh:20
Defines all properties used in Dumux.
Helper to extract native Dune vector types from particular Dumux types.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Type traits to be used with matrix types.
Traits for multidomain problems.
Definition: adapt.hh:17
typename makeFromIndexedType< M, SubDomainDiagBlocks, Indices >::type type
Definition: multidomain/traits.hh:99
typename makeFromIndexedType< std::tuple, PtrType, Indices >::type type
Definition: multidomain/traits.hh:89
typename makeFromIndexedType< std::tuple, PtrType, Indices >::type type
Definition: multidomain/traits.hh:79
Definition: dunevectors.hh:56
typename NativeDuneVectorTypeImpl< V, Dune::Std::is_detected< Detail::DuneVectors::StateDetector, V >{} >::type type
Definition: dunevectors.hh:59
Definition: staggeredtraits.hh:83
typename Dune::BCRSMatrix< MatrixBlock > type
Definition: staggeredtraits.hh:87
typename GetPropType< SubDomainTypeTag< 0 >, Properties::GridGeometry >::FaceFVGridGeometryType type
Definition: staggeredtraits.hh:46
typename GetPropType< SubDomainTypeTag< 0 >, Properties::GridGeometry >::CellCenterFVGridGeometryType type
Definition: staggeredtraits.hh:50
GetPropType< SubDomainTypeTag< i >, Properties::GridGeometry > type
Definition: staggeredtraits.hh:42
typename GetPropType< SubDomainTypeTag< 0 >, Properties::GridVariables >::FaceGridVariablesType type
Definition: staggeredtraits.hh:60
typename GetPropType< SubDomainTypeTag< 0 >, Properties::GridVariables >::CellCenterGridVariablesType type
Definition: staggeredtraits.hh:64
GetPropType< SubDomainTypeTag< i >, Properties::GridVariables > type
Definition: staggeredtraits.hh:56
typename JacobianTypeImpl< GetPropType< SubDomainTypeTag< 0 >, Properties::Scalar >, getPropValue< SubDomainTypeTag< 0 >, Properties::NumEqFace >()>::type type
Definition: staggeredtraits.hh:98
typename JacobianTypeImpl< GetPropType< SubDomainTypeTag< 1 >, Properties::Scalar >, getPropValue< SubDomainTypeTag< 0 >, Properties::NumEqCellCenter >()>::type type
Definition: staggeredtraits.hh:103
GetPropType< SubDomainTypeTag< i >, Properties::JacobianMatrix > type
Definition: staggeredtraits.hh:93
GetPropType< SubDomainTypeTag< 0 >, Properties::FacePrimaryVariables > type
Definition: staggeredtraits.hh:74
GetPropType< SubDomainTypeTag< 0 >, Properties::CellCenterPrimaryVariables > type
Definition: staggeredtraits.hh:78
GetPropType< SubDomainTypeTag< i >, Properties::PrimaryVariables > type
Definition: staggeredtraits.hh:70
GetPropType< SubDomainTypeTag< 0 >, Properties::FaceSolutionVector > type
Definition: staggeredtraits.hh:113
GetPropType< SubDomainTypeTag< 0 >, Properties::CellCenterSolutionVector > type
Definition: staggeredtraits.hh:117
GetPropType< SubDomainTypeTag< i >, Properties::SolutionVector > type
Definition: staggeredtraits.hh:109
Definition: staggeredtraits.hh:189
SubDomainTypeTag< id > TypeTag
Definition: staggeredtraits.hh:191
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
Definition: staggeredtraits.hh:193
typename Detail::Staggered::SubDomainFVGridGeometryImpl< SubDomainTypeTag, id >::type GridGeometry
Definition: staggeredtraits.hh:194
GetPropType< SubDomainTypeTag< id >, Properties::Problem > Problem
Definition: staggeredtraits.hh:192
typename Detail::Staggered::SubDomainSolutionVectorImpl< SubDomainTypeTag, id >::type SolutionVector
Definition: staggeredtraits.hh:196
Dune::index_constant< id > Index
Definition: staggeredtraits.hh:190
typename Detail::Staggered::SubDomainGridVariablesImpl< SubDomainTypeTag, id >::type GridVariables
Definition: staggeredtraits.hh:195
typename Detail::NativeDuneVectorType< SolutionVector >::type ResidualVector
Definition: staggeredtraits.hh:197
typename Detail::Staggered::SubDomainPrimaryVariablesImpl< SubDomainTypeTag, id >::type PrimaryVariables
Definition: staggeredtraits.hh:198
Definition: staggeredtraits.hh:154
typename makeFromIndexedType< Dune::MultiTypeBlockVector, SubDomainResidualVector, Indices >::type ResidualVector
the residual vector type
Definition: staggeredtraits.hh:215
typename Detail::MultiDomainTupleSharedPtr< T, Indices >::type TupleOfSharedPtr
helper alias to create tuple<std::shared_ptr<...>> from indexed type
Definition: staggeredtraits.hh:233
typename Detail::MultiDomainMatrixType< SubDomainJacobianMatrix, Indices, Scalar >::type JacobianMatrix
the jacobian type
Definition: staggeredtraits.hh:218
typename makeFromIndexedType< std::tuple, T, Indices >::type Tuple
helper alias to create tuple<...> from indexed type
Definition: staggeredtraits.hh:229
typename makeFromIndexedType< Dune::MultiTypeBlockVector, SubDomainSolutionVector, Indices >::type SolutionVector
the solution vector type
Definition: staggeredtraits.hh:212
typename makeFromIndexedType< std::common_type_t, SubDomainScalar, Indices >::type Scalar
the scalar type
Definition: staggeredtraits.hh:209
typename Detail::MultiDomainTupleSharedPtrConst< T, Indices >::type TupleOfSharedPtrConst
helper alias to create tuple<std::shared_ptr<const ...>> from indexed type
Definition: staggeredtraits.hh:237
static constexpr std::size_t numSubDomains
the number of subdomains
Definition: staggeredtraits.hh:156
Definition: utility.hh:28