-
Notifications
You must be signed in to change notification settings - Fork 0
Expand file tree
/
Copy pathFinite_Element_imp.h
More file actions
153 lines (130 loc) · 5.14 KB
/
Finite_Element_imp.h
File metadata and controls
153 lines (130 loc) · 5.14 KB
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
//
// Created by simonepanzeri on 26/11/2021.
//
#ifndef DEV_FDAPDE_FINITE_ELEMENT_IMP_H
#define DEV_FDAPDE_FINITE_ELEMENT_IMP_H
// Template auxiliary function forward declaration
// This function evaluate ^phi_i's at quadrature nodes
template<UInt NBASES, UInt mydim>
Eigen::Matrix<Real, NBASES, 1> reference_eval_point(const Point<mydim>& node);
// This function evaluate nabla ^phi_i's at quadrature nodes
template<UInt NBASES, UInt mydim>
Eigen::Matrix<Real, NBASES, mydim> reference_eval_der_point(const Point<mydim>& node);
//! FEData implementation
template <UInt ORDER, UInt mydim, UInt ndim>
FiniteElementData<ORDER, mydim, ndim>::FiniteElementData()
{
setPhi();
setPhiDer();
}
template <UInt ORDER, UInt mydim, UInt ndim>
FiniteElementData<ORDER, mydim, ndim>::~FiniteElementData()
{
}
template <UInt ORDER, UInt mydim, UInt ndim>
void FiniteElementData<ORDER, mydim, ndim>::updateElement(const Element<NBASES,mydim,ndim>& t)
{
t_ = t;
setElementPhiDer();
}
template <UInt ORDER, UInt mydim, UInt ndim>
void FiniteElementData<ORDER, mydim, ndim>::setPhi()
{
for(UInt iq=0; iq<Integrator::NNODES; ++iq)
referencePhi.row(iq) = reference_eval_point<NBASES,mydim>(Integrator::NODES[iq]).transpose();
}
template <UInt ORDER, UInt mydim, UInt ndim>
void FiniteElementData<ORDER, mydim, ndim>::setPhiDer()
{
// .template is needed! See Eigen documentation regarding
// the template and typename keywords in C++
for(UInt iq=0; iq<Integrator::NNODES; ++iq)
referencePhiDer.template block<mydim, NBASES>(0, NBASES*iq) = reference_eval_der_point<NBASES,mydim>(Integrator::NODES[iq]).transpose();
}
template <UInt ORDER, UInt mydim, UInt ndim>
void FiniteElementData<ORDER, mydim, ndim>::setElementPhiDer()
{
// we need J^(-T) nabla( phi)
for (UInt iq=0; iq < Integrator::NNODES; ++iq)
elementPhiDer.template block<ndim, NBASES>(0, NBASES*iq).noalias() = t_.getM_invJ().transpose() * referencePhiDer.template block<mydim, NBASES>(0, NBASES*iq);
}
// Templates for auxiliary functions
// This function evaluates the basis function on the reference element
// at the quadrature nodes
template<>
inline Eigen::Matrix<Real, 3, 1> reference_eval_point(const Point<2>& node){
return makeBaryCoord(node.eigenConstView());
}
template<>
inline Eigen::Matrix<Real, 4, 1> reference_eval_point(const Point<3>& node){
return makeBaryCoord(node.eigenConstView());
}
template<>
inline Eigen::Matrix<Real, 6, 1> reference_eval_point<6,2>(const Point<2>& node){
Eigen::Matrix<Real, 6, 1> phi;
phi << (1-node[0]-node[1]) * (1-2*node[0]-2*node[1]),
node[0] * (2*node[0]-1),
node[1] * (2*node[1]-1),
4*node[0] * node[1],
4*node[1] * (1-node[0]-node[1]),
4*node[0] * (1-node[0]-node[1]);
return phi;
}
template<>
inline Eigen::Matrix<Real, 10, 1> reference_eval_point<10,3>(const Point<3>& node){
Eigen::Matrix<Real, 10, 1> phi;
phi << (1-node[0]-node[1]-node[2]) * (1-2*node[0]-2*node[1]-2*node[2]),
node[0] * (2*node[0]-1),
node[1] * (2*node[1]-1),
node[2] * (2*node[2]-1),
4*node[0] * (1-node[0]-node[1]-node[2]),
4*node[1] * (1-node[0]-node[1]-node[2]),
4*node[2] * (1-node[0]-node[1]-node[2]),
4*node[0] * node[1],
4*node[1] * node[2],
4*node[2] * node[0];
return phi;
}
// This function evaluates the ndim-gradient of basis function on the reference element
// at the quadrature nodes
template<>
inline Eigen::Matrix<Real, 3,2> reference_eval_der_point(const Point<2>& node){
Eigen::Matrix<Real,3,2> B1;
B1.row(0).setConstant(-1);
B1.bottomRows(2).setIdentity();
return B1;
}
template<>
inline Eigen::Matrix<Real, 4,3> reference_eval_der_point(const Point<3>& node){
Eigen::Matrix<Real,4,3> B1;
B1.row(0).setConstant(-1);
B1.bottomRows(3).setIdentity();
return B1;
}
template<>
inline Eigen::Matrix<Real, 6,2> reference_eval_der_point<6,2>(const Point<2>& node){
Eigen::Matrix<Real,6,2> B2;
B2 << 1-4*(1-node[0]-node[1]), 1-4*(1-node[0]-node[1]),
4*node[0]-1, 0,
0, 4*node[1]-1,
4*node[1], 4*node[0],
-4*node[1], 4*(1-node[0]-2*node[1]),
4*(1-2*node[0]-node[1]), -4*node[0];
return B2;
}
template<>
inline Eigen::Matrix<Real, 10,3> reference_eval_der_point<10,3>(const Point<3>& node){
Eigen::Matrix<Real,10,3> B2;
B2 << 1-4*(1-node[0]-node[1]-node[2]), 1-4*(1-node[0]-node[1]-node[2]), 1-4*(1-node[0]-node[1]-node[2]),
4*node[0]-1, 0, 0,
0, 4*node[1]-1, 0,
0, 0, 4*node[2]-1,
4*(1-2*node[0]-node[1]-node[2]), -4*node[0], -4*node[0],
-4*node[1], 4*(1-node[0]-2*node[1]-node[2]), -4*node[1],
-4*node[2], -4*node[2], 4*(1-node[0]-node[1]-2*node[2]),
4*node[1], 4*node[0], 0,
0, 4*node[2], 4*node[1],
4*node[2], 0, 4*node[0];
return B2;
}
#endif //DEV_FDAPDE_FINITE_ELEMENT_IMP_H