3.2-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
31
35
38
46
47#include <dune/istl/multitypeblockvector.hh>
48#include <dune/istl/multitypeblockmatrix.hh>
49
50namespace Dumux {
51
52// forward declarations
53class CCElementBoundaryTypes;
54
55namespace Properties
56{
58// Create new type tags
59namespace TTag {
60struct StaggeredModel { using InheritsFrom = std::tuple<FiniteVolumeModel>; };
61} // end namespace TTag
62
64template<class TypeTag>
65struct GridFaceVariables<TypeTag, TTag::StaggeredModel>
66{
67private:
70 static constexpr auto enableCache = getPropValue<TypeTag, Properties::EnableGridFaceVariablesCache>();
71public:
73};
74
76template<class TypeTag>
77struct EnableGridFaceVariablesCache<TypeTag, TTag::StaggeredModel> { static constexpr bool value = true; };
78
80template<class TypeTag>
81struct GridFluxVariablesCache<TypeTag, TTag::StaggeredModel>
82{
83private:
87 static constexpr auto enableCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
88 static constexpr auto upwindSchemeOrder = getPropValue<TypeTag, Properties::UpwindSchemeOrder>();
89public:
91};
92
94template<class TypeTag>
95struct StaggeredFaceSolution<TypeTag, TTag::StaggeredModel>
96{
97private:
99public:
101};
102
104template<class TypeTag>
105struct GridVariables<TypeTag, TTag::StaggeredModel>
106{
107private:
112public:
114};
115
117template<class TypeTag>
118struct ElementBoundaryTypes<TypeTag, TTag::StaggeredModel> { using type = CCElementBoundaryTypes; };
119
121template<class TypeTag>
122struct BaseLocalResidual<TypeTag, TTag::StaggeredModel> { using type = StaggeredLocalResidual<TypeTag>; };
123
125template<class TypeTag>
126struct CellCenterPrimaryVariables<TypeTag, TTag::StaggeredModel>
127{
128 using type = Dune::FieldVector<GetPropType<TypeTag, Properties::Scalar>,
129 getPropValue<TypeTag, Properties::NumEqCellCenter>()>;
130};
131
133template<class TypeTag>
134struct FacePrimaryVariables<TypeTag, TTag::StaggeredModel>
135{
136 using type = Dune::FieldVector<GetPropType<TypeTag, Properties::Scalar>,
137 getPropValue<TypeTag, Properties::NumEqFace>()>;
138};
139
141template<class TypeTag>
143
144// TODO: bundle SolutionVector, JacobianMatrix
145// in LinearAlgebra traits
146
148template<class TypeTag>
149struct CellCenterSolutionVector<TypeTag, TTag::StaggeredModel>
150{ using type = Dune::BlockVector<GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>>; };
151
153template<class TypeTag>
154struct FaceSolutionVector<TypeTag, TTag::StaggeredModel>
155{ using type = Dune::BlockVector<GetPropType<TypeTag, Properties::FacePrimaryVariables>>; };
156
158template<class TypeTag>
159struct SolutionVector<TypeTag, TTag::StaggeredModel>
160{
161private:
163 using CellCenterSolutionVector = GetPropType<TypeTag, Properties::CellCenterSolutionVector>;
164public:
165 using type = Dune::MultiTypeBlockVector<FaceSolutionVector, CellCenterSolutionVector>;
166};
167
169template<class TypeTag>
170struct JacobianMatrix<TypeTag, TTag::StaggeredModel>
171{
172private:
174
175 static constexpr auto numEqCellCenter = getPropValue<TypeTag, Properties::NumEqCellCenter>();
176 static constexpr auto numEqFace = getPropValue<TypeTag, Properties::NumEqFace>();
177
178public:
179 // the sub-blocks
180 using MatrixLittleBlockCCToCC = typename Dune::FieldMatrix<Scalar, numEqCellCenter, numEqCellCenter>;
181 using MatrixLittleBlockCCToFace = typename Dune::FieldMatrix<Scalar, numEqCellCenter, numEqFace>;
182
183 using MatrixLittleBlockFaceToFace = typename Dune::FieldMatrix<Scalar, numEqFace, numEqFace>;
184 using MatrixLittleBlockFaceToCC = typename Dune::FieldMatrix<Scalar, numEqFace, numEqCellCenter>;
185
186 // the BCRS matrices of the subproblems as big blocks
189
192
193 // the row types
194 using RowFace = typename Dune::MultiTypeBlockVector<MatrixBlockFaceToFace, MatrixBlockFaceToCC>;
195 using RowCellCenter = typename Dune::MultiTypeBlockVector<MatrixBlockCCToFace, MatrixBlockCCToCC>;
196
197 // the jacobian matrix
199};
200
201} // namespace Properties
202} // namespace Dumux
203
204#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.
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Calculates the element-wise residual for the staggered FV scheme.
Definition: staggeredlocalresidual.hh:40
Class to specify the type of a boundary.
Definition: common/boundarytypes.hh:38
Definition: common/pdesolver.hh:37
The type of the base class of the local residual (specific to a discretization scheme)
Definition: common/properties.hh:78
Type of the global jacobian matrix.
Definition: common/properties.hh:80
Vector containing all primary variable vector of the grid.
Definition: common/properties.hh:82
Stores the boundary types of a single degree of freedom.
Definition: common/properties.hh:84
Stores the boundary types on an element.
Definition: common/properties.hh:110
The global vector of flux variable containers.
Definition: common/properties.hh:130
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:134
The solution vector type for cell-centered dofs.
Definition: common/properties.hh:225
The solution vector type for face dofs.
Definition: common/properties.hh:227
Global vector containing face-related data.
Definition: common/properties.hh:229
The primary variables container type for cell-centered dofs.
Definition: common/properties.hh:231
The primary variables container type for face dofs.
Definition: common/properties.hh:233
A vector containing the solution for a face (similar to ElementSolution)
Definition: common/properties.hh:245
Switch on/off caching of face variables.
Definition: common/properties.hh:247
Definition: matrix.hh:32
Boundary types gathered on an element.
Definition: discretization/cellcentered/elementboundarytypes.hh:38
Definition: staggered.hh:60
std::tuple< FiniteVolumeModel > InheritsFrom
Definition: staggered.hh:60
Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, getPropValue< TypeTag, Properties::NumEqCellCenter >()> type
Definition: staggered.hh:129
Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, getPropValue< TypeTag, Properties::NumEqFace >()> type
Definition: staggered.hh:137
Dune::BlockVector< GetPropType< TypeTag, Properties::CellCenterPrimaryVariables > > type
Definition: staggered.hh:150
Dune::BlockVector< GetPropType< TypeTag, Properties::FacePrimaryVariables > > type
Definition: staggered.hh:155
Dune::MultiTypeBlockVector< FaceSolutionVector, CellCenterSolutionVector > type
Definition: staggered.hh:165
typename Dune::BCRSMatrix< MatrixLittleBlockCCToFace > MatrixBlockCCToFace
Definition: staggered.hh:188
typename Dune::MultiTypeBlockVector< MatrixBlockFaceToFace, MatrixBlockFaceToCC > RowFace
Definition: staggered.hh:194
typename Dune::BCRSMatrix< MatrixLittleBlockFaceToFace > MatrixBlockFaceToFace
Definition: staggered.hh:190
typename Dune::FieldMatrix< Scalar, numEqCellCenter, numEqFace > MatrixLittleBlockCCToFace
Definition: staggered.hh:181
typename Dune::MultiTypeBlockMatrix< RowFace, RowCellCenter > type
Definition: staggered.hh:198
typename Dune::BCRSMatrix< MatrixLittleBlockCCToCC > MatrixBlockCCToCC
Definition: staggered.hh:187
typename Dune::FieldMatrix< Scalar, numEqFace, numEqCellCenter > MatrixLittleBlockFaceToCC
Definition: staggered.hh:184
typename Dune::FieldMatrix< Scalar, numEqFace, numEqFace > MatrixLittleBlockFaceToFace
Definition: staggered.hh:183
typename Dune::MultiTypeBlockVector< MatrixBlockCCToFace, MatrixBlockCCToCC > RowCellCenter
Definition: staggered.hh:195
typename Dune::BCRSMatrix< MatrixLittleBlockFaceToCC > MatrixBlockFaceToCC
Definition: staggered.hh:191
typename Dune::FieldMatrix< Scalar, numEqCellCenter, numEqCellCenter > MatrixLittleBlockCCToCC
Definition: staggered.hh:180
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: staggered/gridfluxvariablescache.hh:63
Class storing data associated to scvs and scvfs.
Definition: discretization/staggered/gridvariables.hh:195
Declares all properties used in Dumux.
Boundary types gathered on an element.
Sub control volumes for cell-centered discretization schemes.