version 3.8
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-FileCopyrightInfo: 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
20
24
27
35
36#include <dune/istl/multitypeblockvector.hh>
37#include <dune/istl/multitypeblockmatrix.hh>
38
39namespace Dumux {
40
41// forward declarations
42class CCElementBoundaryTypes;
43
44namespace Properties
45{
47// Create new type tags
48namespace TTag {
49struct StaggeredModel { using InheritsFrom = std::tuple<FiniteVolumeModel>; };
50} // end namespace TTag
51
53template<class TypeTag>
54struct GridFaceVariables<TypeTag, TTag::StaggeredModel>
55{
56private:
59 static constexpr auto enableCache = getPropValue<TypeTag, Properties::EnableGridFaceVariablesCache>();
60public:
62};
63
65template<class TypeTag>
66struct EnableGridFaceVariablesCache<TypeTag, TTag::StaggeredModel> { static constexpr bool value = true; };
67
69template<class TypeTag>
70struct GridFluxVariablesCache<TypeTag, TTag::StaggeredModel>
71{
72private:
75 using FluxVariablesCache = GetPropTypeOr<TypeTag,
76 Properties::FluxVariablesCache, FluxVariablesCaching::EmptyCache<Scalar>
77 >;
78 using FluxVariablesCacheFiller = GetPropTypeOr<TypeTag,
79 Properties::FluxVariablesCacheFiller, FluxVariablesCaching::EmptyCacheFiller
80 >;
81 static constexpr auto enableCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
82 static constexpr auto upwindSchemeOrder = getPropValue<TypeTag, Properties::UpwindSchemeOrder>();
83public:
85};
86
88template<class TypeTag>
89struct StaggeredFaceSolution<TypeTag, TTag::StaggeredModel>
90{
91private:
93public:
95};
96
98template<class TypeTag>
99struct GridVariables<TypeTag, TTag::StaggeredModel>
100{
101private:
106public:
108};
109
111template<class TypeTag>
112struct ElementBoundaryTypes<TypeTag, TTag::StaggeredModel> { using type = CCElementBoundaryTypes; };
113
115template<class TypeTag>
116struct BaseLocalResidual<TypeTag, TTag::StaggeredModel> { using type = StaggeredLocalResidual<TypeTag>; };
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
208} // end namespace Detail
209
210} // namespace Dumux
211
212#endif
Boundary types gathered on an element.
Boundary types gathered on an element.
Definition: cellcentered/elementboundarytypes.hh:26
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:31
Defines all properties used in Dumux.
Type traits for problem classes.
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
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:49
std::tuple< FiniteVolumeModel > InheritsFrom
Definition: staggered.hh:49