-
Notifications
You must be signed in to change notification settings - Fork 1.2k
/
linear_simplex_element.cc
145 lines (136 loc) · 7.11 KB
/
linear_simplex_element.cc
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
#include "drake/multibody/fem/linear_simplex_element.h"
#include "drake/common/default_scalars.h"
namespace drake {
namespace multibody {
namespace fem {
namespace internal {
template <typename T, int natural_dimension_at_compile_time,
int spatial_dimension_at_compile_time,
int num_sample_locations_at_compile_time>
typename LinearSimplexElement<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>::
template ArrayNumSamples<typename IsoparametricElement<
LinearSimplexElement<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>,
LinearSimplexElementTraits<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>>::
GradientInSpatialCoordinates>
LinearSimplexElement<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>::
CalcGradientInSpatialCoordinates(
const Eigen::Ref<const Eigen::Matrix<T, spatial_dimension,
num_nodes>>& xa) const {
return this->DefaultCalcGradientInSpatialCoordinates(xa);
}
template <typename T, int natural_dimension_at_compile_time,
int spatial_dimension_at_compile_time,
int num_sample_locations_at_compile_time>
typename LinearSimplexElement<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>::
template ArrayNumSamples<typename IsoparametricElement<
LinearSimplexElement<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>,
LinearSimplexElementTraits<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>>::
JacobianMatrix>
LinearSimplexElement<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>::
CalcJacobian(const Eigen::Ref<
const Eigen::Matrix<T, spatial_dimension, num_nodes>>& xa)
const {
return this->DefaultCalcJacobian(xa);
}
template <typename T, int natural_dimension_at_compile_time,
int spatial_dimension_at_compile_time,
int num_sample_locations_at_compile_time>
typename LinearSimplexElement<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>::
template ArrayNumSamples<typename IsoparametricElement<
LinearSimplexElement<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>,
LinearSimplexElementTraits<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>>::
PseudoinverseJacobianMatrix>
LinearSimplexElement<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>::
CalcJacobianPseudoinverse(
const ArrayNumSamples<JacobianMatrix>& jacobian) const {
return this->DefaultCalcJacobianPseudoinverse(jacobian);
}
template <typename T, int natural_dimension_at_compile_time,
int spatial_dimension_at_compile_time,
int num_sample_locations_at_compile_time>
typename LinearSimplexElement<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>::
template ArrayNumSamples<typename LinearSimplexElement<
T, natural_dimension_at_compile_time, spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>::VectorNumNodes>
LinearSimplexElement<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>::
CalcShapeFunctions(
const ArrayNumSamples<Vector<double, natural_dimension>>& locations)
const {
ArrayNumSamples<VectorNumNodes> S;
for (int q = 0; q < num_sample_locations; ++q) {
VectorNumNodes Sq;
/* Sₐ = ξₐ₋₁ for a = 1, ..., num_nodes - 1. */
for (int a = 1; a < num_nodes; ++a) {
Sq(a) = locations[q](a - 1);
}
/* S₀ = 1−ξ₀ − ... − ξₙ₋₂. */
Sq(0) = 1 - Sq.template tail<num_nodes - 1>().sum();
S[q] = std::move(Sq);
}
return S;
}
template <typename T, int natural_dimension_at_compile_time,
int spatial_dimension_at_compile_time,
int num_sample_locations_at_compile_time>
const typename LinearSimplexElement<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>::
template ArrayNumSamples<typename LinearSimplexElement<
T, natural_dimension_at_compile_time, spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>::GradientInParentCoordinates>
LinearSimplexElement<T, natural_dimension_at_compile_time,
spatial_dimension_at_compile_time,
num_sample_locations_at_compile_time>::
CalcGradientInParentCoordinates() {
GradientInParentCoordinates dSdxi_q;
dSdxi_q.template topRows<1>() = -1.0 * Vector<T, natural_dimension>::Ones();
dSdxi_q.template bottomRows<natural_dimension>() =
Eigen::Matrix<T, natural_dimension, natural_dimension>::Identity();
ArrayNumSamples<GradientInParentCoordinates> dSdxi;
dSdxi.fill(dSdxi_q);
return dSdxi;
}
template class LinearSimplexElement<double, 2, 2, 2>;
template class LinearSimplexElement<double, 2, 2, 4>;
template class LinearSimplexElement<double, 2, 3, 4>;
template class LinearSimplexElement<double, 3, 3, 1>;
template class LinearSimplexElement<double, 3, 3, 2>;
template class LinearSimplexElement<double, 3, 3, 5>;
template class LinearSimplexElement<AutoDiffXd, 2, 2, 2>;
template class LinearSimplexElement<AutoDiffXd, 2, 2, 4>;
template class LinearSimplexElement<AutoDiffXd, 2, 3, 4>;
template class LinearSimplexElement<AutoDiffXd, 3, 3, 1>;
template class LinearSimplexElement<AutoDiffXd, 3, 3, 2>;
template class LinearSimplexElement<AutoDiffXd, 3, 3, 5>;
} // namespace internal
} // namespace fem
} // namespace multibody
} // namespace drake