-
Notifications
You must be signed in to change notification settings - Fork 2
/
SVKElasticForcefield.cpp
419 lines (334 loc) · 19.2 KB
/
SVKElasticForcefield.cpp
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
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
#include "SVKElasticForcefield.h"
#include <sofa/core/ObjectFactory.h>
#include <sofa/core/visual/VisualParams.h>
SVKElasticForcefield::SVKElasticForcefield()
: d_youngModulus(initData(&d_youngModulus,
Real(1000), "youngModulus",
"Young's modulus of the material",
true /*displayed_in_GUI*/, false /*read_only_in_GUI*/))
, d_poissonRatio(initData(&d_poissonRatio,
Real(0.3), "poissonRatio",
"Poisson's ratio of the material",
true /*displayed_in_GUI*/, false /*read_only_in_GUI*/))
, d_topology_container(initLink(
"topology", "Topology that contains the elements on which this force will be computed."))
{
}
void SVKElasticForcefield::init() {
using Mat33 = Eigen::Matrix<double, 3, 3>;
ForceField::init();
if (!this->mstate.get() || !d_topology_container.get()) {
msg_error() << "Both a mechanical object and a topology container are required";
}
auto * state = this->mstate.get();
auto * topology = d_topology_container.get();
// Convert SOFA position vector to an Eigen matrix (nx3 for n nodes)
const auto sofa_x0 = state->readRestPositions();
Eigen::Map<const Eigen::Matrix<Real, Eigen::Dynamic, 3, Eigen::RowMajor>> x0 (sofa_x0.ref().data()->data(), state->getSize(), 3);
// Convert the node index vector from SOFA to an Eigen matrix (nxm for n elements of m nodes each)
Eigen::Map<const Eigen::Matrix<sofa::Index, Eigen::Dynamic, Element::NumberOfNodes, Eigen::RowMajor>> node_indices (
topology->getTetras().data()->data(), topology->getNbTetrahedra(), Element::NumberOfNodes
);
// Gather the integration points for each tetrahedron
const auto number_of_elements = topology->getNbTetrahedra();
p_gauss_nodes.resize(number_of_elements);
for (Eigen::Index element_id = 0; element_id < number_of_elements; ++element_id) {
// Position vector of each of the element nodes
Eigen::Matrix<double, Element::NumberOfNodes, 3, Eigen::RowMajor> node_positions;
for (Eigen::Index node_id = 0; node_id < Element::NumberOfNodes; ++node_id) {
node_positions.row(node_id) = x0.row(node_indices(element_id, node_id));
}
// Initialize each Gauss nodes
for (std::size_t gauss_node_id = 0; gauss_node_id < Element::NumberOfGaussNode; ++ gauss_node_id) {
const auto & gauss_node = Element::gauss_nodes()[gauss_node_id];
const auto dN_du = Element::dN(gauss_node.position[0], gauss_node.position[1], gauss_node.position[2]);
const Mat33 J = node_positions.transpose() * dN_du;
const Mat33 Jinv = J.inverse();
const auto detJ = J.determinant();
const Eigen::Matrix<double, Element::NumberOfNodes, 3> dN_dx = (Jinv.transpose() * dN_du.transpose()).transpose();
p_gauss_nodes[element_id][gauss_node_id].weight = gauss_node.weight;
p_gauss_nodes[element_id][gauss_node_id].jacobian_determinant = detJ;
p_gauss_nodes[element_id][gauss_node_id].dN_dx = dN_dx;
}
}
}
double SVKElasticForcefield::getPotentialEnergy(const sofa::core::MechanicalParams *,
const Data<sofa::type::vector<Coord>> & d_x) const {
using Mat33 = Eigen::Matrix<double, 3, 3>;
if (!this->mstate.get() || !d_topology_container.get()) {
return 0;
}
auto * state = this->mstate.get();
auto * topology = d_topology_container.get();
const auto poisson_ratio = d_poissonRatio.getValue();
const auto young_modulus = d_youngModulus.getValue();
const auto mu = young_modulus / (2.0 * (1.0 + poisson_ratio));
const auto lm = young_modulus * poisson_ratio / ((1.0 + poisson_ratio) * (1.0 - 2.0 * poisson_ratio));
static const auto Id = Mat33::Identity();
// Convert SOFA input position vector to an Eigen matrix (nx3 for n nodes)
auto sofa_x = sofa::helper::getReadAccessor(d_x);
Eigen::Map<const Eigen::Matrix<Real, Eigen::Dynamic, 3, Eigen::RowMajor>> x (sofa_x.ref().data()->data(), state->getSize(), 3);
// Convert the node index vector from SOFA to an Eigen matrix (nxm for n elements of m nodes each)
Eigen::Map<const Eigen::Matrix<sofa::Index, Eigen::Dynamic, Element::NumberOfNodes, Eigen::RowMajor>> node_indices (
topology->getTetras().data()->data(), topology->getNbTetrahedra(), Element::NumberOfNodes
);
double Psi = 0.;
const auto nb_elements = topology->getNbTetrahedra();
for (Eigen::Index element_id = 0; element_id < nb_elements; ++element_id) {
// Position vector of each of the element nodes
Eigen::Matrix<double, Element::NumberOfNodes, 3, Eigen::RowMajor> node_positions;
for (Eigen::Index node_id = 0; node_id < Element::NumberOfNodes; ++node_id) {
node_positions.row(node_id) = x.row(node_indices(element_id, node_id));
}
for (const GaussNode & gauss_node : p_gauss_nodes[element_id]) {
// Jacobian of the gauss node's transformation mapping from the elementary space to the world space
const auto & detJ = gauss_node.jacobian_determinant;
// Derivatives of the shape functions at the gauss node with respect to global coordinates x,y and z
const auto & dN_dx = gauss_node.dN_dx;
// Gauss quadrature node weight
const auto & w = gauss_node.weight;
// Deformation tensor at gauss node
const Mat33 F = node_positions.transpose()*dN_dx;
// Green-Lagrange strain tensor at gauss node
const Mat33 E = 1/2. * (F.transpose() * F - Id);
const double trE = E.trace();
const double trEE = (E*E).trace();
// Add the potential energy at gauss node
Psi += (detJ * w) * lm/2.*(trE*trE) + mu*trEE;
}
}
return Psi;
}
void SVKElasticForcefield::addForce(const sofa::core::MechanicalParams */*mparams*/,
Data<sofa::type::vector<Deriv>> &d_f,
const Data<sofa::type::vector<Coord>> &d_x,
const Data<sofa::type::vector<Deriv>> &/*d_v*/) {
using Mat33 = Eigen::Matrix<double, 3, 3>;
using Vec3 = Eigen::Matrix<double, 3, 1>;
if (!this->mstate.get() || !d_topology_container.get()) {
return;
}
auto * state = this->mstate.get();
auto * topology = d_topology_container.get();
const auto poisson_ratio = d_poissonRatio.getValue();
const auto young_modulus = d_youngModulus.getValue();
const auto mu = young_modulus / (2.0 * (1.0 + poisson_ratio));
const auto lm = young_modulus * poisson_ratio / ((1.0 + poisson_ratio) * (1.0 - 2.0 * poisson_ratio));
static const auto Id = Mat33::Identity();
// Convert SOFA input position vector to an Eigen matrix (nx3 for n nodes)
auto sofa_x = sofa::helper::getReadAccessor(d_x);
Eigen::Map<const Eigen::Matrix<Real, Eigen::Dynamic, 3, Eigen::RowMajor>> x (sofa_x.ref().data()->data(), state->getSize(), 3);
// Convert SOFA output residual vector to an Eigen matrix (nx3 for n nodes)
auto sofa_f = sofa::helper::getWriteAccessor(d_f);
Eigen::Map<Eigen::Matrix<Real, Eigen::Dynamic, 3, Eigen::RowMajor>> R (&(sofa_f[0][0]), state->getSize(), 3);
// Convert the node index vector from SOFA to an Eigen matrix (nxm for n elements of m nodes each)
Eigen::Map<const Eigen::Matrix<sofa::Index, Eigen::Dynamic, Element::NumberOfNodes, Eigen::RowMajor>> node_indices (
topology->getTetras().data()->data(), topology->getNbTetrahedra(), Element::NumberOfNodes
);
// Assemble the residual vector
const auto nb_elements = topology->getNbTetrahedra();
for (Eigen::Index element_id = 0; element_id < nb_elements; ++element_id) {
// Position vector of each of the element nodes
Eigen::Matrix<double, Element::NumberOfNodes, 3, Eigen::RowMajor> node_positions;
for (Eigen::Index node_id = 0; node_id < Element::NumberOfNodes; ++node_id) {
node_positions.row(node_id) = x.row(node_indices(element_id, node_id));
}
for (const GaussNode & gauss_node : p_gauss_nodes[element_id]) {
// Jacobian of the gauss node's transformation mapping from the elementary space to the world space
const auto & detJ = gauss_node.jacobian_determinant;
// Derivatives of the shape functions at the gauss node with respect to global coordinates x,y and z
const auto & dN_dx = gauss_node.dN_dx;
// Gauss quadrature node weight
const auto & w = gauss_node.weight;
// Deformation tensor at gauss node
const Mat33 F = node_positions.transpose()*dN_dx;
// Green-Lagrange strain tensor at gauss node
const Mat33 E = 1/2. * (F.transpose() * F - Id);
// Second Piola-Kirchhoff stress tensor at gauss node
const Mat33 S = lm*E.trace()*Id + 2*mu*E;
// Elastic forces (residual) w.r.t the gauss node applied on each nodes
for (Eigen::Index i = 0; i < Element::NumberOfNodes; ++i) {
const auto dx = dN_dx.row(i).transpose();
const Vec3 f_ = (detJ * w) * F*S*dx;
R.row(node_indices(element_id, i)).noalias() -= f_.transpose();
}
}
}
}
void SVKElasticForcefield::addKToMatrix(sofa::defaulttype::BaseMatrix * matrix,
double kFact,
unsigned int & offset) {
using Mat33 = Eigen::Matrix<double, 3, 3>;
using Vec3 = Eigen::Matrix<double, 3, 1>;
if (!this->mstate.get() || !d_topology_container.get()) {
return;
}
auto * state = this->mstate.get();
auto * topology = d_topology_container.get();
const auto poisson_ratio = d_poissonRatio.getValue();
const auto young_modulus = d_youngModulus.getValue();
const auto mu = young_modulus / (2.0 * (1.0 + poisson_ratio));
const auto lm = young_modulus * poisson_ratio / ((1.0 + poisson_ratio) * (1.0 - 2.0 * poisson_ratio));
static const auto Id = Mat33::Identity();
// Convert SOFA input position vector to an Eigen matrix (nx3 for n nodes)
auto sofa_x = state->readPositions();
Eigen::Map<const Eigen::Matrix<Real, Eigen::Dynamic, 3, Eigen::RowMajor>> x (sofa_x.ref().data()->data(), state->getSize(), 3);
// Convert the node index vector from SOFA to an Eigen matrix (nxm for n elements of m nodes each)
Eigen::Map<const Eigen::Matrix<sofa::Index, Eigen::Dynamic, Element::NumberOfNodes, Eigen::RowMajor>> node_indices (
topology->getTetras().data()->data(), topology->getNbTetrahedra(), Element::NumberOfNodes
);
// Assemble the stiffness matrix
const auto nb_elements = topology->getNbTetrahedra();
for (Eigen::Index element_id = 0; element_id < nb_elements; ++element_id) {
// Position vector of each of the element nodes
Eigen::Matrix<double, Element::NumberOfNodes, 3, Eigen::RowMajor> node_positions;
for (Eigen::Index node_id = 0; node_id < Element::NumberOfNodes; ++node_id) {
node_positions.row(node_id) = x.row(node_indices(element_id, node_id));
}
for (const GaussNode & gauss_node : p_gauss_nodes[element_id]) {
// Jacobian of the gauss node's transformation mapping from the elementary space to the world space
const auto & detJ = gauss_node.jacobian_determinant;
// Derivatives of the shape functions at the gauss node with respect to global coordinates x,y and z
const auto & dN_dx = gauss_node.dN_dx;
// Gauss quadrature node weight
const auto & w = gauss_node.weight;
// Deformation tensor at gauss node
// (Note: this is usually saved up for each Gauss node at the residual assembly step to improve a bit the computational cost)
const Mat33 F = node_positions.transpose()*dN_dx;
// Green-Lagrange strain tensor at gauss node
const Mat33 E = 1/2. * (F.transpose() * F - Id);
// Second Piola-Kirchhoff stress tensor at gauss node
const Mat33 S = lm*E.trace()*Id + 2*mu*E;
Eigen::Matrix<double, 6, 6> D;
D <<
lm + 2*mu, lm, lm, 0, 0, 0,
lm, lm + 2*mu, lm, 0, 0, 0,
lm, lm, lm + 2*mu, 0, 0, 0,
0, 0, 0, mu, 0, 0,
0, 0, 0, 0, mu, 0,
0, 0, 0, 0, 0, mu;
// Symmetric elemental tangent stiffness matrix
for (Eigen::Index i = 0; i < Element::NumberOfNodes; ++i) {
const Vec3 dxi = dN_dx.row(i).transpose();
const auto I = static_cast<int>(offset + node_indices(element_id, i) * 3);
Eigen::Matrix<double, 6,3> Bi;
Bi <<
F(0,0)*dxi[0], F(1,0)*dxi[0], F(2,0)*dxi[0],
F(0,1)*dxi[1], F(1,1)*dxi[1], F(2,1)*dxi[1],
F(0,2)*dxi[2], F(1,2)*dxi[2], F(2,2)*dxi[2],
F(0,0)*dxi[1] + F(0,1)*dxi[0], F(1,0)*dxi[1] + F(1,1)*dxi[0], F(2,0)*dxi[1] + F(2,1)*dxi[0],
F(0,1)*dxi[2] + F(0,2)*dxi[1], F(1,1)*dxi[2] + F(1,2)*dxi[1], F(2,1)*dxi[2] + F(2,2)*dxi[1],
F(0,0)*dxi[2] + F(0,2)*dxi[0], F(1,0)*dxi[2] + F(1,2)*dxi[0], F(2,0)*dxi[2] + F(2,2)*dxi[0];
// The 3x3 sub-matrix Kii is symmetric
Mat33 Kii = (dxi.dot(S*dxi)*Id + Bi.transpose()*D*Bi) * (detJ * w) * -kFact;
for (int m = 0; m < 3; ++m) {
matrix->add(I + m, I + m, Kii(m, m));
for (int n = m+1; n < 3; ++n) {
matrix->add(I + m, I + n, Kii(m, n));
matrix->add(I + n, I + m, Kii(m, n));
}
}
// We now loop only on the upper triangular part of the
// element stiffness matrix Ke since it is symmetric
for (Eigen::Index j = i+1; j < Element::NumberOfNodes; ++j) {
const Vec3 dxj = dN_dx.row(j).transpose();
const auto J = static_cast<int>(offset + node_indices(element_id,j) * 3);
Eigen::Matrix<double, 6,3> Bj;
Bj <<
F(0,0)*dxj[0], F(1,0)*dxj[0], F(2,0)*dxj[0],
F(0,1)*dxj[1], F(1,1)*dxj[1], F(2,1)*dxj[1],
F(0,2)*dxj[2], F(1,2)*dxj[2], F(2,2)*dxj[2],
F(0,0)*dxj[1] + F(0,1)*dxj[0], F(1,0)*dxj[1] + F(1,1)*dxj[0], F(2,0)*dxj[1] + F(2,1)*dxj[0],
F(0,1)*dxj[2] + F(0,2)*dxj[1], F(1,1)*dxj[2] + F(1,2)*dxj[1], F(2,1)*dxj[2] + F(2,2)*dxj[1],
F(0,0)*dxj[2] + F(0,2)*dxj[0], F(1,0)*dxj[2] + F(1,2)*dxj[0], F(2,0)*dxj[2] + F(2,2)*dxj[0];
// The 3x3 sub-matrix Kij is NOT symmetric, we store its full part
const Mat33 Kij = (dxi.dot(S*dxj)*Id + Bi.transpose()*D*Bj) * (detJ * w) * -kFact;
for (int m = 0; m < 3; ++m) {
for (int n = 0; n < 3; ++n) {
matrix->add(I + m, J + n, Kij(m, n));
matrix->add(J + n, I + m, Kij(m, n));
}
}
}
}
}
}
}
void SVKElasticForcefield::addDForce(const sofa::core::MechanicalParams * /*mparams*/,
SVKElasticForcefield::Data<sofa::type::vector<sofa::type::Vec3>> & /*d_df*/,
const SVKElasticForcefield::Data<sofa::type::vector<sofa::type::Vec3>> & /*d_dx*/) {
// Here you would compute df = K*dx
}
void SVKElasticForcefield::draw(const sofa::core::visual::VisualParams *vparams) {
if (!this->mstate.get() || !d_topology_container.get()) {
return;
}
if (!vparams->displayFlags().getShowForceFields()) {
return;
}
auto * state = this->mstate.get();
auto * topology = d_topology_container.get();
vparams->drawTool()->disableLighting();
const auto x = state->readPositions();
std::vector< sofa::type::Vec<3, double> > points[4];
const auto number_of_elements = topology->getNbTetrahedra();
for (sofa::core::topology::Topology::TetrahedronID i = 0 ; i<number_of_elements;++i) {
const auto t=topology->getTetra(i);
const auto & a = t[0];
const auto & b = t[1];
const auto & c = t[2];
const auto & d = t[3];
Coord center = (x[a]+x[b]+x[c]+x[d])*0.125;
Coord pa = (x[a]+center)*(Real)0.666667;
Coord pb = (x[b]+center)*(Real)0.666667;
Coord pc = (x[c]+center)*(Real)0.666667;
Coord pd = (x[d]+center)*(Real)0.666667;
points[0].push_back(pa);
points[0].push_back(pb);
points[0].push_back(pc);
points[1].push_back(pb);
points[1].push_back(pc);
points[1].push_back(pd);
points[2].push_back(pc);
points[2].push_back(pd);
points[2].push_back(pa);
points[3].push_back(pd);
points[3].push_back(pa);
points[3].push_back(pb);
}
sofa::type::RGBAColor face_colors[4] = {
{1.0, 0.0, 0.0, 1.0},
{1.0, 0.0, 0.5, 1.0},
{1.0, 1.0, 0.0, 1.0},
{1.0, 0.5, 1.0, 1.0}
};
vparams->drawTool()->drawTriangles(points[0], face_colors[0]);
vparams->drawTool()->drawTriangles(points[1], face_colors[1]);
vparams->drawTool()->drawTriangles(points[2], face_colors[2]);
vparams->drawTool()->drawTriangles(points[3], face_colors[3]);
if (vparams->displayFlags().getShowWireFrame())
vparams->drawTool()->setPolygonMode(0,false);
vparams->drawTool()->restoreLastState();
}
void SVKElasticForcefield::computeBBox(const sofa::core::ExecParams * /*params*/, bool onlyVisible) {
using namespace sofa::core::objectmodel;
if (!onlyVisible) return;
if (!this->mstate) return;
sofa::helper::ReadAccessor<Data < VecCoord>>
x = this->mstate->read(sofa::core::VecCoordId::position());
static const Real max_real = std::numeric_limits<Real>::max();
static const Real min_real = std::numeric_limits<Real>::lowest();
Real maxBBox[3] = {min_real, min_real, min_real};
Real minBBox[3] = {max_real, max_real, max_real};
for (size_t i = 0; i < x.size(); i++) {
for (int c = 0; c < 3; c++) {
if (x[i][c] > maxBBox[c]) maxBBox[c] = static_cast<Real>(x[i][c]);
else if (x[i][c] < minBBox[c]) minBBox[c] = static_cast<Real>(x[i][c]);
}
}
this->f_bbox.setValue(sofa::type::TBoundingBox<Real>(minBBox, maxBBox));
}
using sofa::core::RegisterObject;
[[maybe_unused]]
static int _c_ = RegisterObject("Simple implementation of a Saint-Venant-Kirchhoff force field for tetrahedral meshes.")
.add<SVKElasticForcefield>();