Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
/**
* \copyright
* Copyright (c) 2012-2016, OpenGeoSys Community (http://www.opengeosys.org)
* Distributed under a Modified BSD License.
* See accompanying file LICENSE.txt or
* http://www.opengeosys.org/project/license
*
*/
#ifndef PROCESSLIB_HYDROMECHANICS_LOCALDATAINITIALIZER_H_
#define PROCESSLIB_HYDROMECHANICS_LOCALDATAINITIALIZER_H_
#include <functional>
#include <memory>
#include <typeindex>
#include <typeinfo>
#include <type_traits>
#include <unordered_map>
#include "MeshLib/Elements/Elements.h"
#include "NumLib/DOF/LocalToGlobalIndexMap.h"
#include "NumLib/Fem/FiniteElement/LowerDimShapeTable.h"
#include "NumLib/Fem/Integration/GaussIntegrationPolicy.h"
#ifndef OGS_MAX_ELEMENT_DIM
static_assert(false, "The macro OGS_MAX_ELEMENT_DIM is undefined.");
#endif
#ifndef OGS_MAX_ELEMENT_ORDER
static_assert(false, "The macro OGS_MAX_ELEMENT_ORDER is undefined.");
#endif
// The following macros decide which element types will be compiled, i.e.
// which element types will be available for use in simulations.
#ifdef OGS_ENABLE_ELEMENT_SIMPLEX
#define ENABLED_ELEMENT_TYPE_SIMPLEX 1u
#else
#define ENABLED_ELEMENT_TYPE_SIMPLEX 0u
#endif
#ifdef OGS_ENABLE_ELEMENT_CUBOID
#define ENABLED_ELEMENT_TYPE_CUBOID 1u << 1
#else
#define ENABLED_ELEMENT_TYPE_CUBOID 0u
#endif
#ifdef OGS_ENABLE_ELEMENT_PRISM
#define ENABLED_ELEMENT_TYPE_PRISM 1u << 2
#else
#define ENABLED_ELEMENT_TYPE_PRISM 0u
#endif
#ifdef OGS_ENABLE_ELEMENT_PYRAMID
#define ENABLED_ELEMENT_TYPE_PYRAMID 1u << 3
#else
#define ENABLED_ELEMENT_TYPE_PYRAMID 0u
#endif
// Dependent element types.
// Faces of tets, pyramids and prisms are triangles
#define ENABLED_ELEMENT_TYPE_TRI \
((ENABLED_ELEMENT_TYPE_SIMPLEX) | (ENABLED_ELEMENT_TYPE_PYRAMID) | \
(ENABLED_ELEMENT_TYPE_PRISM))
// Faces of hexes, pyramids and prisms are quads
#define ENABLED_ELEMENT_TYPE_QUAD \
((ENABLED_ELEMENT_TYPE_CUBOID) | (ENABLED_ELEMENT_TYPE_PYRAMID) | \
(ENABLED_ELEMENT_TYPE_PRISM))
// All enabled element types
#define OGS_ENABLED_ELEMENTS \
((ENABLED_ELEMENT_TYPE_SIMPLEX) | (ENABLED_ELEMENT_TYPE_CUBOID) | \
(ENABLED_ELEMENT_TYPE_PYRAMID) | (ENABLED_ELEMENT_TYPE_PRISM))
// Include only what is needed (Well, the conditions are not sharp).
#if (OGS_ENABLED_ELEMENTS & ENABLED_ELEMENT_TYPE_SIMPLEX) != 0
#include "NumLib/Fem/ShapeFunction/ShapeTet10.h"
#include "NumLib/Fem/ShapeFunction/ShapeTet4.h"
#endif
#if (OGS_ENABLED_ELEMENTS & ENABLED_ELEMENT_TYPE_TRI) != 0
#include "NumLib/Fem/ShapeFunction/ShapeTri3.h"
#include "NumLib/Fem/ShapeFunction/ShapeTri6.h"
#endif
#if (OGS_ENABLED_ELEMENTS & ENABLED_ELEMENT_TYPE_CUBOID) != 0
#include "NumLib/Fem/ShapeFunction/ShapeHex20.h"
#include "NumLib/Fem/ShapeFunction/ShapeHex8.h"
#endif
#if (OGS_ENABLED_ELEMENTS & ENABLED_ELEMENT_TYPE_QUAD) != 0
#include "NumLib/Fem/ShapeFunction/ShapeQuad4.h"
#include "NumLib/Fem/ShapeFunction/ShapeQuad8.h"
#include "NumLib/Fem/ShapeFunction/ShapeQuad9.h"
#endif
#if (OGS_ENABLED_ELEMENTS & ENABLED_ELEMENT_TYPE_PRISM) != 0
#include "NumLib/Fem/ShapeFunction/ShapePrism15.h"
#include "NumLib/Fem/ShapeFunction/ShapePrism6.h"
#endif
#if (OGS_ENABLED_ELEMENTS & ENABLED_ELEMENT_TYPE_PYRAMID) != 0
#include "NumLib/Fem/ShapeFunction/ShapePyra13.h"
#include "NumLib/Fem/ShapeFunction/ShapePyra5.h"
#endif
namespace ProcessLib
{
namespace HydroMechanics
{
/// The LocalDataInitializer is a functor creating a local assembler data with
/// corresponding to the mesh element type shape functions and calling
/// initialization of the new local assembler data.
/// For example for MeshLib::Quad a local assembler data with template argument
/// NumLib::ShapeQuad4 is created.
///
/// \attention This is modified version of the ProcessLib::LocalDataInitializer
/// class which does not include line elements, allows only shapefunction of
/// order 2, and provides additional template argument DisplacementDim.
template <typename LocalAssemblerInterface,
template <typename, typename, typename, unsigned, int>
class LocalAssemblerData,
unsigned GlobalDim, int DisplacementDim, typename... ConstructorArgs>
class LocalDataInitializer final
{
public:
using LADataIntfPtr = std::unique_ptr<LocalAssemblerInterface>;
LocalDataInitializer(NumLib::LocalToGlobalIndexMap const& dof_table,
const unsigned shapefunction_order)
: _dof_table(dof_table)
{
if (shapefunction_order != 2)
OGS_FATAL(
"The given shape function order %d is not supported.\nOnly "
"shape functions of order 2 are supported.",
shapefunction_order);
// /// Quads and Hexahedra ///////////////////////////////////
#if (OGS_ENABLED_ELEMENTS & ENABLED_ELEMENT_TYPE_QUAD) != 0 && \
OGS_MAX_ELEMENT_DIM >= 2 && OGS_MAX_ELEMENT_ORDER >= 2
_builder[std::type_index(typeid(MeshLib::Quad8))] =
makeLocalAssemblerBuilder<NumLib::ShapeQuad8>();
_builder[std::type_index(typeid(MeshLib::Quad9))] =
makeLocalAssemblerBuilder<NumLib::ShapeQuad9>();
#endif
#if (OGS_ENABLED_ELEMENTS & ENABLED_ELEMENT_TYPE_CUBOID) != 0 && \
OGS_MAX_ELEMENT_DIM >= 3 && OGS_MAX_ELEMENT_ORDER >= 2
_builder[std::type_index(typeid(MeshLib::Hex20))] =
makeLocalAssemblerBuilder<NumLib::ShapeHex20>();
#endif
// /// Simplices ////////////////////////////////////////////////
#if (OGS_ENABLED_ELEMENTS & ENABLED_ELEMENT_TYPE_TRI) != 0 && \
OGS_MAX_ELEMENT_DIM >= 2 && OGS_MAX_ELEMENT_ORDER >= 2
_builder[std::type_index(typeid(MeshLib::Tri6))] =
makeLocalAssemblerBuilder<NumLib::ShapeTri6>();
#endif
#if (OGS_ENABLED_ELEMENTS & ENABLED_ELEMENT_TYPE_SIMPLEX) != 0 && \
OGS_MAX_ELEMENT_DIM >= 3 && OGS_MAX_ELEMENT_ORDER >= 2
_builder[std::type_index(typeid(MeshLib::Tet10))] =
makeLocalAssemblerBuilder<NumLib::ShapeTet10>();
#endif
// /// Prisms ////////////////////////////////////////////////////
#if (OGS_ENABLED_ELEMENTS & ENABLED_ELEMENT_TYPE_PRISM) != 0 && \
OGS_MAX_ELEMENT_DIM >= 3 && OGS_MAX_ELEMENT_ORDER >= 2
_builder[std::type_index(typeid(MeshLib::Prism15))] =
makeLocalAssemblerBuilder<NumLib::ShapePrism15>();
#endif
// /// Pyramids //////////////////////////////////////////////////
#if (OGS_ENABLED_ELEMENTS & ENABLED_ELEMENT_TYPE_PYRAMID) != 0 && \
OGS_MAX_ELEMENT_DIM >= 3 && OGS_MAX_ELEMENT_ORDER >= 2
_builder[std::type_index(typeid(MeshLib::Pyramid13))] =
makeLocalAssemblerBuilder<NumLib::ShapePyra13>();
#endif
}
/// Sets the provided \c data_ptr to the newly created local assembler data.
///
/// \attention
/// The index \c id is not necessarily the mesh item's id. Especially when
/// having multiple meshes it will differ from the latter.
void operator()(std::size_t const id,
MeshLib::Element const& mesh_item,
LADataIntfPtr& data_ptr,
ConstructorArgs&&... args) const
{
auto const type_idx = std::type_index(typeid(mesh_item));
auto const it = _builder.find(type_idx);
if (it != _builder.end())
{
auto const num_local_dof = _dof_table.getNumberOfElementDOF(id);
data_ptr = it->second(mesh_item, num_local_dof,
std::forward<ConstructorArgs>(args)...);
}
else
{
OGS_FATAL(
"You are trying to build a local assembler for an unknown mesh "
"element type (%s)."
" Maybe you have disabled this mesh element type in your build "
"configuration.",
type_idx.name());
}
}
private:
using LADataBuilder =
std::function<LADataIntfPtr(MeshLib::Element const& e,
std::size_t const local_matrix_size,
ConstructorArgs&&...)>;
template <typename ShapeFunctionDisplacement>
using IntegrationMethod = typename NumLib::GaussIntegrationPolicy<
typename ShapeFunctionDisplacement::MeshElement>::IntegrationMethod;
template <typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure>
using LAData =
LocalAssemblerData<ShapeFunctionDisplacement, ShapeFunctionPressure,
IntegrationMethod<ShapeFunctionDisplacement>,
GlobalDim, DisplacementDim>;
/// A helper forwarding to the correct version of makeLocalAssemblerBuilder
/// depending whether the global dimension is less than the shape function's
/// dimension or not.
template <typename ShapeFunctionDisplacement>
static LADataBuilder makeLocalAssemblerBuilder()
{
return makeLocalAssemblerBuilder<ShapeFunctionDisplacement>(
static_cast<std::integral_constant<
bool, (GlobalDim >= ShapeFunctionDisplacement::DIM)>*>(
nullptr));
}
/// Mapping of element types to local assembler constructors.
std::unordered_map<std::type_index, LADataBuilder> _builder;
NumLib::LocalToGlobalIndexMap const& _dof_table;
// local assembler builder implementations.
private:
/// Generates a function that creates a new LocalAssembler of type
/// LAData<ShapeFunctionDisplacement>. Only functions with shape function's
/// dimension less or equal to the global dimension are instantiated, e.g.
/// following combinations of shape functions and global dimensions: (Line2,
/// 1),
/// (Line2, 2), (Line2, 3), (Hex20, 3) but not (Hex20, 2) or (Hex20, 1).
template <typename ShapeFunctionDisplacement>
static LADataBuilder makeLocalAssemblerBuilder(std::true_type*)
{
// (Lower order elements = Order(ShapeFunctionDisplacement) - 1).
using ShapeFunctionPressure =
typename NumLib::LowerDim<ShapeFunctionDisplacement>::type;
return [](MeshLib::Element const& e,
std::size_t const local_matrix_size,
ConstructorArgs&&... args) {
return LADataIntfPtr{
new LAData<ShapeFunctionDisplacement, ShapeFunctionPressure>{
e, local_matrix_size,
std::forward<ConstructorArgs>(args)...}};
};
}
/// Returns nullptr for shape functions whose dimensions are less than the
/// global dimension.
template <typename ShapeFunctionDisplacement>
static LADataBuilder makeLocalAssemblerBuilder(std::false_type*)
{
return nullptr;
}
};
} // namespace HydroMechanics
} // namespace ProcessLib
#undef ENABLED_ELEMENT_TYPE_SIMPLEX
#undef ENABLED_ELEMENT_TYPE_CUBOID
#undef ENABLED_ELEMENT_TYPE_PYRAMID
#undef ENABLED_ELEMENT_TYPE_PRISM
#undef ENABLED_ELEMENT_TYPE_TRI
#undef ENABLED_ELEMENT_TYPE_QUAD
#undef OGS_ENABLED_ELEMENTS
#endif // PROCESSLIB_HYDROMECHANICS_LOCALDATAINITIALIZER_H_