DUNE-FEM (unstable)

orthogonal.hh
1#ifndef DUNE_FEM_HPDG_SPACE_DISCONTINUOUSGALERKIN_ORTHOGONAL_HH
2#define DUNE_FEM_HPDG_SPACE_DISCONTINUOUSGALERKIN_ORTHOGONAL_HH
3
4#include <dune/grid/common/gridenums.hh>
5
6#include <dune/fem/operator/matrix/istlmatrixadapter.hh>
7
8#include <dune/fem/space/common/capabilities.hh>
9#include <dune/fem/space/common/commoperations.hh>
10#include <dune/fem/space/common/defaultcommhandler.hh>
11#include <dune/fem/space/common/localrestrictprolong.hh>
12#include <dune/fem/space/discontinuousgalerkin/localrestrictprolong.hh>
13
14#include <dune/fem/space/shapefunctionset/selectcaching.hh>
15#include <dune/fem/space/basisfunctionset/hpdg/orthogonal.hh>
16
17#include "blockmapper.hh"
18#include "space.hh"
19
20namespace Dune
21{
22
23 namespace Fem
24 {
25
26 namespace hpDG
27 {
28
29 // Internal forward declaration
30 // ----------------------------
31
32 template< class FunctionSpace, class GridPart, int order, class Storage = Fem::CachingStorage >
33 class OrthogonalDiscontinuousGalerkinSpace;
34
35
36
37#ifndef DOXYGEN
38
39 // OrthogonalDiscontinuousGalerkinSpaceTraits
40 // ------------------------------------------
41
42 template< class FunctionSpace, class GridPart, int order, class Storage >
43 struct OrthogonalDiscontinuousGalerkinSpaceTraits
44 {
45 using DiscreteFunctionSpaceType = hpDG::OrthogonalDiscontinuousGalerkinSpace< FunctionSpace, GridPart, order, Storage >;
46
47 using FunctionSpaceType = FunctionSpace;
48
49 using GridPartType = GridPart;
50
51 using BasisFunctionSetsType = hpDG::OrthogonalBasisFunctionSets< FunctionSpaceType, GridPartType, order, Storage >;
52 using BasisFunctionSetType = typename BasisFunctionSetsType::BasisFunctionSetType;
53
54 static const int codimension = BasisFunctionSetType::EntityType::codimension;
55
56 using BlockMapperType = hpDG::DiscontinuousGalerkinBlockMapper< GridPartType, BasisFunctionSetsType >;
57 static const int localBlockSize = BasisFunctionSetsType::localBlockSize;
58 static_assert( localBlockSize == FunctionSpace::dimRange, " dimRange prob ");
59
60 typedef Hybrid::IndexRange< int, localBlockSize > LocalBlockIndices;
61
62 template< class DiscreteFunction, class Operation = Dune::Fem::DFCommunicationOperation::Copy >
63 struct CommDataHandle
64 {
65 using OperationType = Operation;
67 };
68 };
69
70#endif // #ifndef DOXYGEN
71
72
73
74 // OrthogonalDiscontinuousGalerkinSpace
75 // ------------------------------------
76
86 template< class FunctionSpace, class GridPart, int order, class Storage >
88 : public hpDG::DiscontinuousGalerkinSpace< OrthogonalDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, Storage > >
89 {
91
92 public:
93
94 static const int polynomialOrder = order ;
95
96 using GridPartType = typename BaseType::GridPartType;
97 using EntityType = typename BaseType::EntityType;
98 using BasisFunctionSetsType = typename BaseType::BasisFunctionSetsType;
99
100 explicit OrthogonalDiscontinuousGalerkinSpace ( GridPartType &gridPart,
103 : BaseType( gridPart, BasisFunctionSetsType{}, order, interface, direction )
104 {}
105
106 template <class Function,
107 std::enable_if_t<
108 std::is_arithmetic<
109 decltype(Function(std::declval<const EntityType>()))>::value,int> i=0>
113 : BaseType( gridPart, BasisFunctionSetsType{}, order, function, interface, direction )
114 {}
115 };
116
117 } // namespace hpDG
118
120 template <class FunctionSpaceImp,
121 class GridPartImp,
122 int polOrd,
123 class Storage,
124 class VolumeQuadratureImp>
126 hpDG::OrthogonalDiscontinuousGalerkinSpace< FunctionSpaceImp, GridPartImp, polOrd, Storage >, VolumeQuadratureImp >
128 hpDG::OrthogonalDiscontinuousGalerkinSpace< FunctionSpaceImp, GridPartImp, polOrd, Storage >, VolumeQuadratureImp >
129 {
132 public:
133 using BaseType::BaseType;
134 };
135
136
137#ifndef DOXYGEN
138
139 // DefaultLocalRestrictProlong
140 // ---------------------------
141
142 template< class FunctionSpace, class GridPart, int order, class Storage >
143 class DefaultLocalRestrictProlong< hpDG::OrthogonalDiscontinuousGalerkinSpace< FunctionSpace, GridPart, order, Storage > >
144 : public DiscontinuousGalerkinLocalRestrictProlong< hpDG::OrthogonalDiscontinuousGalerkinSpace< FunctionSpace, GridPart, order, Storage >, false >
145 {
146 using BaseType = DiscontinuousGalerkinLocalRestrictProlong< hpDG::OrthogonalDiscontinuousGalerkinSpace< FunctionSpace, GridPart, order, Storage >, false >;
147
148 public:
149 explicit DefaultLocalRestrictProlong ( const typename BaseType::DiscreteFunctionSpaceType &space )
150 : BaseType( space )
151 {}
152 };
153
154
155#endif // #ifndef DOXYGEN
156
157
158 namespace Capabilities
159 {
161 // hpDG::OrthogonalDiscontinuousGalerkinSpace
163
164 template< class FunctionSpace, class GridPart, int polOrder, class Storage >
165 struct hasStaticPolynomialOrder< hpDG::OrthogonalDiscontinuousGalerkinSpace< FunctionSpace, GridPart, polOrder, Storage > >
166 {
167 static const bool v = true;
168 static const int order = polOrder;
169 };
170
171 template< class FunctionSpace, class GridPart, int polOrder, class Storage >
172 struct isLocalized< hpDG::OrthogonalDiscontinuousGalerkinSpace< FunctionSpace, GridPart, polOrder, Storage > >
173 {
174 static const bool v = true;
175 };
176
177 template< class FunctionSpace, class GridPart, int polOrder, class Storage >
178 struct isAdaptive< hpDG::OrthogonalDiscontinuousGalerkinSpace< FunctionSpace, GridPart, polOrder, Storage > >
179 {
180 static const bool v = true;
181 };
182
183 template< class FunctionSpace, class GridPart, int polOrder, class Storage >
184 struct viewThreadSafe< hpDG::OrthogonalDiscontinuousGalerkinSpace< FunctionSpace, GridPart, polOrder, Storage > >
185 {
186 static const bool v = true;
187 };
188
189 template< class FunctionSpace, class GridPart, int polOrder, class Storage >
190 struct isHierarchic< hpDG::OrthogonalDiscontinuousGalerkinSpace< FunctionSpace, GridPart, polOrder, Storage > >
191 {
192 static const bool v = true;
193 };
194
195 } // namespace Capabilities
196
197
198 } // namespace Fem
199
200} // namespace Dune
201
202#endif // #ifndef DUNE_FEM_HPDG_SPACE_DISCONTINUOUSGALERKIN_ORTHOGONAL_HH
Default communication handler for discrete functions.
Definition: defaultcommhandler.hh:38
@ dimRange
dimension of range vector space
Definition: functionspaceinterface.hh:48
A vector valued function space.
Definition: functionspace.hh:60
Abstract class representing a function.
Definition: function.hh:50
DG Local Mass Matrix for arbitrary spaces.
Definition: localmassmatrix.hh:987
Local Mass Matrix inversion implementation, select the correct method in your implementation.
Definition: localmassmatrix.hh:37
Local Mass Matrix for arbitrary spaces.
Definition: localmassmatrix.hh:967
Generic implementation of a -adaptive discontinuous finite element space.
Definition: space.hh:46
typename Traits::BasisFunctionSetsType BasisFunctionSetsType
basis function sets type
Definition: space.hh:56
Implementation of an -adaptive discrete function space using orthogonal polynomials.
Definition: orthogonal.hh:89
CommunicationDirection
Define a type for communication direction parameter.
Definition: gridenums.hh:170
InterfaceType
Parameter to be used for the communication functions.
Definition: gridenums.hh:86
@ ForwardCommunication
communicate as given in InterfaceType
Definition: gridenums.hh:171
@ InteriorBorder_All_Interface
send interior and border, receive all entities
Definition: gridenums.hh:88
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)