Loading [MathJax]/extensions/tex2jax.js
OGS
Toggle main menu visibility
Main Page
Related Pages
Namespaces
Namespace List
Namespace Members
All
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Functions
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Variables
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
Typedefs
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
v
w
z
Enumerations
b
c
d
e
f
g
i
l
m
n
o
p
s
t
u
v
Enumerator
a
b
c
d
e
f
g
h
l
m
n
p
r
s
t
v
y
Classes
Class List
Class Index
Class Hierarchy
Files
File List
File Members
All
a
b
c
d
e
f
g
i
k
m
n
o
p
r
s
t
v
w
x
Functions
a
c
d
e
f
g
i
m
o
p
r
s
t
v
w
Variables
Typedefs
Enumerations
Macros
b
c
g
m
n
o
p
r
s
t
ExtrapolatableElementCollection.h
Go to the documentation of this file.
1
11
#pragma once
12
13
#include <functional>
14
#include <vector>
15
16
#include "
ExtrapolatableElement.h
"
17
#include "
MathLib/LinAlg/GlobalMatrixVectorTypes.h
"
18
19
namespace
NumLib
20
{
21
class
LocalToGlobalIndexMap;
22
29
class
ExtrapolatableElementCollection
30
{
31
public
:
34
virtual
Eigen::Map<const Eigen::RowVectorXd>
getShapeMatrix
(
35
std::size_t
const
id
,
unsigned
const
integration_point)
const
= 0;
36
63
virtual
std::vector<double>
const
&
getIntegrationPointValues
(
64
std::size_t
const
id
,
const
double
t,
65
std::vector<GlobalVector*>
const
& x,
66
std::vector<NumLib::LocalToGlobalIndexMap const*>
const
& dof_table,
67
std::vector<double>& cache)
const
= 0;
68
70
virtual
std::size_t
size
()
const
= 0;
71
72
virtual
~ExtrapolatableElementCollection
() =
default
;
73
};
29
class
ExtrapolatableElementCollection
{
…
};
74
75
template
<
typename
LocalAssemblerCollection>
76
class
ExtrapolatableLocalAssemblerCollection
77
:
public
ExtrapolatableElementCollection
78
{
79
public
:
81
using
LocalAssembler
=
82
typename
std::decay_t<decltype(*std::declval<LocalAssemblerCollection>()
83
[
static_cast<
std::size_t
>
(0)])>;
84
85
static_assert
(std::is_base_of<ExtrapolatableElement, LocalAssembler>::value,
86
"Local assemblers used for extrapolation must be derived "
87
"from ExtrapolatableElement."
);
88
97
using
IntegrationPointValuesMethod
=
98
std::function<std::vector<double>
const
&(
99
LocalAssembler
const
& loc_asm,
const
double
t,
100
std::vector<GlobalVector*>
const
& x,
101
std::vector<NumLib::LocalToGlobalIndexMap const*>
const
& dof_table,
102
std::vector<double>& cache)>;
103
111
ExtrapolatableLocalAssemblerCollection
(
112
LocalAssemblerCollection
const
& local_assemblers,
113
IntegrationPointValuesMethod
const
& integration_point_values_method)
114
:
_local_assemblers
(local_assemblers),
115
_integration_point_values_method
{integration_point_values_method}
116
{
117
}
111
ExtrapolatableLocalAssemblerCollection
( {
…
}
118
119
Eigen::Map<const Eigen::RowVectorXd>
getShapeMatrix
(
120
std::size_t
const
id
,
unsigned
const
integration_point)
const override
121
{
122
ExtrapolatableElement
const
& loc_asm = *
_local_assemblers
[id];
123
return
loc_asm.
getShapeMatrix
(integration_point);
124
}
119
Eigen::Map<const Eigen::RowVectorXd>
getShapeMatrix
( {
…
}
125
126
std::vector<double>
const
&
getIntegrationPointValues
(
127
std::size_t
const
id
,
const
double
t,
128
std::vector<GlobalVector*>
const
& x,
129
std::vector<NumLib::LocalToGlobalIndexMap const*>
const
& dof_table,
130
std::vector<double>& cache)
const override
131
{
132
auto
const
& loc_asm = *
_local_assemblers
[id];
133
return
_integration_point_values_method
(loc_asm, t, x, dof_table,
134
cache);
135
}
126
std::vector<double>
const
&
getIntegrationPointValues
( {
…
}
136
137
std::size_t
size
()
const override
{
return
_local_assemblers
.size(); }
138
139
private
:
140
LocalAssemblerCollection
const
&
_local_assemblers
;
141
IntegrationPointValuesMethod
const
_integration_point_values_method
;
142
};
76
class
ExtrapolatableLocalAssemblerCollection
{
…
};
143
146
template
<
typename
LocalAssemblerCollection,
147
typename
IntegrationPointValuesMethod>
148
ExtrapolatableLocalAssemblerCollection<LocalAssemblerCollection>
149
makeExtrapolatable
(LocalAssemblerCollection
const
& local_assemblers,
150
IntegrationPointValuesMethod integration_point_values_method)
151
{
152
return
ExtrapolatableLocalAssemblerCollection<LocalAssemblerCollection>
{
153
local_assemblers, integration_point_values_method};
154
}
149
makeExtrapolatable
(LocalAssemblerCollection
const
& local_assemblers, {
…
}
155
}
// namespace NumLib
ExtrapolatableElement.h
GlobalMatrixVectorTypes.h
NumLib::ExtrapolatableElementCollection
Definition
ExtrapolatableElementCollection.h:30
NumLib::ExtrapolatableElementCollection::getShapeMatrix
virtual Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(std::size_t const id, unsigned const integration_point) const =0
NumLib::ExtrapolatableElementCollection::size
virtual std::size_t size() const =0
Returns the number of elements whose properties shall be extrapolated.
NumLib::ExtrapolatableElementCollection::~ExtrapolatableElementCollection
virtual ~ExtrapolatableElementCollection()=default
NumLib::ExtrapolatableElementCollection::getIntegrationPointValues
virtual std::vector< double > const & getIntegrationPointValues(std::size_t const id, const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
NumLib::ExtrapolatableElement
Definition
ExtrapolatableElement.h:24
NumLib::ExtrapolatableElement::getShapeMatrix
virtual Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const =0
Provides the shape matrix at the given integration point.
NumLib::ExtrapolatableLocalAssemblerCollection
Definition
ExtrapolatableElementCollection.h:78
NumLib::ExtrapolatableLocalAssemblerCollection::IntegrationPointValuesMethod
std::function< std::vector< double > const &( LocalAssembler const &loc_asm, const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache)> IntegrationPointValuesMethod
Definition
ExtrapolatableElementCollection.h:97
NumLib::ExtrapolatableLocalAssemblerCollection::getIntegrationPointValues
std::vector< double > const & getIntegrationPointValues(std::size_t const id, const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
Definition
ExtrapolatableElementCollection.h:126
NumLib::ExtrapolatableLocalAssemblerCollection::LocalAssembler
typename std::decay_t< decltype(*std::declval< LocalAssemblerCollection >()[static_cast< std::size_t >(0)])> LocalAssembler
LocalAssemblerCollection contains many LocalAssembler's.
Definition
ExtrapolatableElementCollection.h:81
NumLib::ExtrapolatableLocalAssemblerCollection::_local_assemblers
LocalAssemblerCollection const & _local_assemblers
Definition
ExtrapolatableElementCollection.h:140
NumLib::ExtrapolatableLocalAssemblerCollection::getShapeMatrix
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(std::size_t const id, unsigned const integration_point) const override
Definition
ExtrapolatableElementCollection.h:119
NumLib::ExtrapolatableLocalAssemblerCollection::size
std::size_t size() const override
Returns the number of elements whose properties shall be extrapolated.
Definition
ExtrapolatableElementCollection.h:137
NumLib::ExtrapolatableLocalAssemblerCollection::_integration_point_values_method
IntegrationPointValuesMethod const _integration_point_values_method
Definition
ExtrapolatableElementCollection.h:141
NumLib::ExtrapolatableLocalAssemblerCollection::ExtrapolatableLocalAssemblerCollection
ExtrapolatableLocalAssemblerCollection(LocalAssemblerCollection const &local_assemblers, IntegrationPointValuesMethod const &integration_point_values_method)
Definition
ExtrapolatableElementCollection.h:111
NumLib
Definition
ProjectData.h:46
NumLib::makeExtrapolatable
ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection > makeExtrapolatable(LocalAssemblerCollection const &local_assemblers, IntegrationPointValuesMethod integration_point_values_method)
Definition
ExtrapolatableElementCollection.h:149
NumLib
Extrapolation
ExtrapolatableElementCollection.h
Generated by
1.12.0