Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
MatParams.cpp
Go to the documentation of this file.
1#include "MatParams.hpp"
2
5#include <iostream>
6
7namespace polyfem::assembler
8{
9 namespace
10 {
11 double convert_to_lambda(const bool is_volume, const double E, const double nu)
12 {
13 if (is_volume)
14 return (E * nu) / ((1.0 + nu) * (1.0 - 2.0 * nu));
15
16 return (nu * E) / (1.0 - nu * nu);
17 }
18
19 double convert_to_mu(const double E, const double nu)
20 {
21 return E / (2.0 * (1.0 + nu));
22 }
23 } // namespace
24
25 GenericMatParam::GenericMatParam(const std::string &param_name)
26 : param_name_(param_name)
27 {
28 }
29
30 void GenericMatParam::add_multimaterial(const int index, const json &params, const std::string &unit_type)
31 {
32 for (int i = param_.size(); i <= index; ++i)
33 {
34 param_.emplace_back();
35 }
36
37 if (params.count(param_name_))
38 {
39 param_[index].init(params[param_name_]);
40 param_[index].set_unit_type(unit_type);
41 }
42 }
43
44 double GenericMatParam::operator()(const RowVectorNd &p, double t, int index) const
45 {
46 const double x = p(0);
47 const double y = p(1);
48 const double z = p.size() == 3 ? p(2) : 0;
49
50 return (*this)(x, y, z, t, index);
51 }
52
53 double GenericMatParam::operator()(double x, double y, double z, double t, int index) const
54 {
55 assert(param_.size() == 1 || index < param_.size());
56
57 const auto &tmp_param = param_.size() == 1 ? param_[0] : param_[index];
58
59 return tmp_param(x, y, z, t, index);
60 }
61
62 GenericMatParams::GenericMatParams(const std::string &param_name)
63 : param_name_(param_name)
64 {
65 }
66
67 void GenericMatParams::add_multimaterial(const int index, const json &params, const std::string &unit_type)
68 {
69 if (!params.contains(param_name_))
70 return;
71
72 std::vector<json> params_array = utils::json_as_array(params[param_name_]);
73 assert(params_array.size() == params_.size() || params_.empty());
74
75 if (params_.empty())
76 for (int i = 0; i < params_array.size(); ++i)
77 params_.emplace_back(param_name_ + "_" + std::to_string(i));
78
79 for (int i = 0; i < params_.size(); ++i)
80 {
81 for (int j = params_.at(i).param_.size(); j <= index; ++j)
82 {
83 params_.at(i).param_.emplace_back();
84 }
85
86 params_.at(i).param_[index].init(params_array[i]);
87 params_.at(i).param_[index].set_unit_type(unit_type);
88 }
89 }
90
91 void ElasticityTensor::resize(const int size)
92 {
93 if (size == 2)
94 stiffness_tensor_.resize(3, 3);
95 else
96 stiffness_tensor_.resize(6, 6);
97
98 stiffness_tensor_.setZero();
99
100 size_ = size;
101 }
102
103 double ElasticityTensor::operator()(int i, int j) const
104 {
105 if (j < i)
106 {
107 std::swap(i, j);
108 }
109
110 assert(j >= i);
111 return stiffness_tensor_(i, j);
112 }
113
114 double &ElasticityTensor::operator()(int i, int j)
115 {
116 if (j < i)
117 {
118 std::swap(i, j);
119 }
120
121 assert(j >= i);
122 return stiffness_tensor_(i, j);
123 }
124
125 void ElasticityTensor::set_from_entries(const std::vector<double> &entries, const std::string &stress_units)
126 {
127 if (size_ == 2)
128 {
129 if (entries.size() == 4)
130 {
132 entries[0],
133 entries[1],
134 entries[2],
135 entries[3], stress_units);
136
137 return;
138 }
139
140 assert(entries.size() >= 6);
141
142 (*this)(0, 0) = entries[0];
143 (*this)(0, 1) = entries[1];
144 (*this)(0, 2) = entries[2];
145
146 (*this)(1, 1) = entries[3];
147 (*this)(1, 2) = entries[4];
148
149 (*this)(2, 2) = entries[5];
150 }
151 else
152 {
153 if (entries.size() == 5)
154 {
156 entries[0],
157 entries[1],
158 entries[2],
159 entries[3],
160 entries[4],
161 stress_units);
162
163 return;
164 }
165 else if (entries.size() == 9)
166 {
168 entries[0],
169 entries[1],
170 entries[2],
171 entries[3],
172 entries[4],
173 entries[5],
174 entries[6],
175 entries[7],
176 entries[8], stress_units);
177
178 return;
179 }
180 assert(entries.size() >= 21);
181
182 (*this)(0, 0) = entries[0];
183 (*this)(0, 1) = entries[1];
184 (*this)(0, 2) = entries[2];
185 (*this)(0, 3) = entries[3];
186 (*this)(0, 4) = entries[4];
187 (*this)(0, 5) = entries[5];
188
189 (*this)(1, 1) = entries[6];
190 (*this)(1, 2) = entries[7];
191 (*this)(1, 3) = entries[8];
192 (*this)(1, 4) = entries[9];
193 (*this)(1, 5) = entries[10];
194
195 (*this)(2, 2) = entries[11];
196 (*this)(2, 3) = entries[12];
197 (*this)(2, 4) = entries[13];
198 (*this)(2, 5) = entries[14];
199
200 (*this)(3, 3) = entries[15];
201 (*this)(3, 4) = entries[16];
202 (*this)(3, 5) = entries[17];
203
204 (*this)(4, 4) = entries[18];
205 (*this)(4, 5) = entries[19];
206
207 (*this)(5, 5) = entries[20];
208 }
209 }
210
211 void ElasticityTensor::set_from_lambda_mu(const double lambda, const double mu, const std::string &stress_units)
212 {
213 if (size_ == 2)
214 {
215 (*this)(0, 0) = 2 * mu + lambda;
216 (*this)(0, 1) = lambda;
217 (*this)(0, 2) = 0;
218
219 (*this)(1, 1) = 2 * mu + lambda;
220 (*this)(1, 2) = 0;
221
222 (*this)(2, 2) = mu;
223 }
224 else
225 {
226 (*this)(0, 0) = 2 * mu + lambda;
227 (*this)(0, 1) = lambda;
228 (*this)(0, 2) = lambda;
229 (*this)(0, 3) = 0;
230 (*this)(0, 4) = 0;
231 (*this)(0, 5) = 0;
232
233 (*this)(1, 1) = 2 * mu + lambda;
234 (*this)(1, 2) = lambda;
235 (*this)(1, 3) = 0;
236 (*this)(1, 4) = 0;
237 (*this)(1, 5) = 0;
238
239 (*this)(2, 2) = 2 * mu + lambda;
240 (*this)(2, 3) = 0;
241 (*this)(2, 4) = 0;
242 (*this)(2, 5) = 0;
243
244 (*this)(3, 3) = mu;
245 (*this)(3, 4) = 0;
246 (*this)(3, 5) = 0;
247
248 (*this)(4, 4) = mu;
249 (*this)(4, 5) = 0;
250
251 (*this)(5, 5) = mu;
252 }
253 }
254
255 void ElasticityTensor::set_from_young_poisson(const double young, const double nu, const std::string &stress_units)
256 {
257 if (size_ == 2)
258 {
259 stiffness_tensor_ << 1.0, nu, 0.0,
260 nu, 1.0, 0.0,
261 0.0, 0.0, (1.0 - nu) / 2.0;
262 stiffness_tensor_ *= young / (1.0 - nu * nu);
263 }
264 else
265 {
266 assert(size_ == 3);
267 const double v = nu;
268 stiffness_tensor_ << 1. - v, v, v, 0, 0, 0,
269 v, 1. - v, v, 0, 0, 0,
270 v, v, 1. - v, 0, 0, 0,
271 0, 0, 0, (1. - 2. * v) / 2., 0, 0,
272 0, 0, 0, 0, (1. - 2. * v) / 2., 0,
273 0, 0, 0, 0, 0, (1. - 2. * v) / 2.;
274 stiffness_tensor_ *= young / ((1. + v) * (1. - 2. * v));
275 }
276 }
277
279 double Ex, double Ey, double Ez,
280 double nuXY, double nuXZ, double nuYZ,
281 double muYZ, double muZX, double muXY, const std::string &stress_units)
282 {
283 assert(size_ == 3);
284
285 // from https://www.efunda.com/formulae/solid_mechanics/mat_mechanics/hooke_orthotropic.cfm
286 double nuYX = nuXY * Ey / Ex;
287 double nuZX = nuXZ * Ez / Ex;
288 double nuZY = nuYZ * Ez / Ey;
289
290 Eigen::MatrixXd compliance;
291 compliance.setZero(6, 6);
292 compliance << 1 / Ex, -nuYX / Ey, -nuZX / Ez, 0, 0, 0,
293 -nuXY / Ex, 1 / Ey, -nuZY / Ez, 0, 0, 0,
294 -nuXZ / Ex, -nuYZ / Ey, 1 / Ez, 0, 0, 0,
295 0, 0, 0, 1 / (2 * muYZ), 0, 0,
296 0, 0, 0, 0, 1 / (2 * muZX), 0,
297 0, 0, 0, 0, 0, 1 / (2 * muXY);
298 stiffness_tensor_ = compliance.inverse();
299 }
300
301 void ElasticityTensor::set_orthotropic(double Ex, double Ey, double nuXY, double muXY, const std::string &stress_units)
302 {
303 assert(size_ == 2);
304
305 // from https://www.efunda.com/formulae/solid_mechanics/mat_mechanics/hooke_orthotropic.cfm
306 double nuYX = nuXY * Ey / Ex;
307
308 Eigen::MatrixXd compliance;
309 compliance.setZero(3, 3);
310 compliance << 1.0 / Ex, -nuYX / Ey, 0.0,
311 -nuXY / Ex, 1.0 / Ey, 0.0,
312 0.0, 0.0, 1.0 / (2 * muXY);
313 stiffness_tensor_ = compliance.inverse();
314 }
315
317 double Et, double Ea,
318 double nu_t, double nu_a,
319 double Ga, const std::string &stress_units)
320 {
321 assert(size_ == 3);
322
323 // from https://osupdocs.forestry.oregonstate.edu/index.php/Transversely_Isotropic_Material
324 Eigen::MatrixXd compliance;
325 compliance.setZero(6, 6);
326 compliance << 1 / Et, -nu_t / Et, -nu_a / Ea, 0, 0, 0,
327 -nu_t / Et, 1 / Et, -nu_a / Ea, 0, 0, 0,
328 -nu_a / Ea, -nu_a / Ea, 1 / Ea, 0, 0, 0,
329 0, 0, 0, 1 / Ga, 0, 0,
330 0, 0, 0, 0, 1 / Ga, 0,
331 0, 0, 0, 0, 0, (2 * (1 + nu_t)) / Et;
332 stiffness_tensor_ = compliance.inverse();
333 }
334
335 template <int DIM>
336 double ElasticityTensor::compute_stress(const std::array<double, DIM> &strain, const int j) const
337 {
338 double res = 0;
339
340 for (int k = 0; k < DIM; ++k)
341 res += (*this)(j, k) * strain[k];
342
343 return res;
344 }
345
346 void ElasticityTensor::rotate_stiffness(const Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 6, 6> &rotation_mtx_voigt)
347 {
349 stiffness_tensor_ = rotation_mtx_voigt * stiffness_tensor_ * rotation_mtx_voigt.transpose();
350 }
351
356
358 {
359 lambda_or_E_.emplace_back();
360 lambda_or_E_.back().init(1.0);
361
362 mu_or_nu_.emplace_back();
363 mu_or_nu_.back().init(1.0);
364 size_ = -1;
365 is_lambda_mu_ = true;
366 }
367
368 void LameParameters::lambda_mu(double px, double py, double pz, double x, double y, double z, double t, int el_id, double &lambda, double &mu) const
369 {
370 assert(lambda_or_E_.size() == 1 || el_id < lambda_or_E_.size());
371 assert(mu_or_nu_.size() == 1 || el_id < mu_or_nu_.size());
372 assert(size_ == 2 || size_ == 3);
373
374 const auto &tmp1 = lambda_or_E_.size() == 1 ? lambda_or_E_[0] : lambda_or_E_[el_id];
375 const auto &tmp2 = mu_or_nu_.size() == 1 ? mu_or_nu_[0] : mu_or_nu_[el_id];
376
377 double llambda = tmp1(x, y, z, t, el_id);
378 double mmu = tmp2(x, y, z, t, el_id);
379
380 if (!is_lambda_mu_)
381 {
382 lambda = convert_to_lambda(size_ == 3, llambda, mmu);
383 mu = convert_to_mu(llambda, mmu);
384 }
385 else
386 {
387 lambda = llambda;
388 mu = mmu;
389 }
390
391 if (lambda_mat_.size() > el_id && mu_mat_.size() > el_id)
392 {
393 lambda = lambda_mat_(el_id);
394 mu = mu_mat_(el_id);
395 }
396
397 assert(!std::isnan(lambda));
398 assert(!std::isnan(mu));
399 assert(!std::isinf(lambda));
400 assert(!std::isinf(mu));
401 }
402
403 void LameParameters::add_multimaterial(const int index, const json &params, const bool is_volume, const std::string &stress_unit)
404 {
405 const int size = is_volume ? 3 : 2;
406 assert(size_ == -1 || size == size_);
407 size_ = size;
408
409 for (int i = lambda_or_E_.size(); i <= index; ++i)
410 {
411 lambda_or_E_.emplace_back();
412 mu_or_nu_.emplace_back();
413 }
414
415 if (params.count("young"))
416 {
417 set_e_nu(index, params["young"], params["nu"], stress_unit);
418 }
419 else if (params.count("E"))
420 {
421 set_e_nu(index, params["E"], params["nu"], stress_unit);
422 }
423 else if (params.count("lambda"))
424 {
425 lambda_or_E_[index].init(params["lambda"]);
426 mu_or_nu_[index].init(params["mu"]);
427
428 lambda_or_E_[index].set_unit_type(stress_unit);
429 mu_or_nu_[index].set_unit_type(stress_unit);
430 is_lambda_mu_ = true;
431 }
432 }
433
434 void LameParameters::set_e_nu(const int index, const json &E, const json &nu, const std::string &stress_unit)
435 {
436 // TODO: conversion is always called
437 is_lambda_mu_ = false;
438 lambda_or_E_[index].init(E);
439 mu_or_nu_[index].init(nu);
440
441 lambda_or_E_[index].set_unit_type(stress_unit);
442 // nu has no unit
443 mu_or_nu_[index].set_unit_type("");
444 }
445
447 {
448 rho_.emplace_back();
449 rho_.back().init(1.0);
450 }
451
452 double Density::operator()(double px, double py, double pz, double x, double y, double z, double t, int el_id) const
453 {
454 assert(rho_.size() == 1 || el_id < rho_.size());
455
456 const auto &tmp = rho_.size() == 1 ? rho_[0] : rho_[el_id];
457 const double res = tmp(x, y, z, t, el_id);
458 assert(!std::isnan(res));
459 assert(!std::isinf(res));
460 return res;
461 }
462
463 void Density::add_multimaterial(const int index, const json &params, const std::string &density_unit)
464 {
465 for (int i = rho_.size(); i <= index; ++i)
466 {
467 rho_.emplace_back();
468 }
469
470 if (params.count("rho"))
471 {
472 rho_[index].init(params["rho"]);
473 }
474 else if (params.count("density"))
475 {
476 rho_[index].init(params["density"]);
477 }
478
479 rho_[index].set_unit_type(density_unit);
480 }
481
485
486 void FiberDirection::resize(const int size)
487 {
488 assert(size == 2 || size == 3);
489 size_ = size;
490 if (!dir_.empty())
491 {
492 for (const auto &m : dir_)
493 {
494 assert(m.rows() == size && m.cols() == size);
495 }
496 }
497 }
498
499 Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 1, 3, 3> FiberDirection::operator()(double px, double py, double pz, double x, double y, double z, double t, int el_id) const
500 {
501 assert(dir_.size() == 1 || el_id < dir_.size());
502
503 const auto &tmp = dir_.size() == 1 ? dir_[0] : dir_[el_id];
504 Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 1, 3, 3> res;
505 res.resize(tmp.rows(), tmp.cols());
506 for (int i = 0; i < tmp.rows(); ++i)
507 {
508 for (int j = 0; j < tmp.cols(); ++j)
509 {
510 res(i, j) = tmp(i, j)(x, y, z, t, el_id);
511
512 assert(!std::isnan(res(i, j)));
513 assert(!std::isinf(res(i, j)));
514 }
515 }
516 return res;
517 }
518
519 Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 1, 6, 6> FiberDirection::stiffness_rotation_voigt(double px, double py, double pz, double x, double y, double z, double t, int el_id) const
520 {
521 // Rotate stiffness mtx in voigt notation according to:
522 // https://scicomp.stackexchange.com/questions/35600/4th-order-tensor-rotation-sources-to-refer
523
524 Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 1, 3, 3> rot = (*this)(px, py, pz, x, y, z, t, el_id);
525 Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 1, 6, 6> res;
526
527 int dim = rot.rows();
528
529 static const double sqrt2 = std::sqrt(2.0);
530
531 if (dim == 2)
532 {
533 res.resize(3, 3);
534 // Still need to compute for 2d
535 assert(false);
536 // res << rot(0, 0) * rot(0, 0), rot(0, 1) * rot(0, 1), 0,
537 // rot(1, 0) * rot(1, 0), rot(1, 1) * rot(1, 1), 0,
538 // 0, 0, 1;
539 }
540 else
541 {
542 assert(dim == 3);
543 res.resize(6, 6);
544 res << rot(0, 0) * rot(0, 0), rot(0, 1) * rot(0, 1), rot(0, 2) * rot(0, 2), sqrt2 * rot(0, 1) * rot(0, 2), sqrt2 * rot(0, 0) * rot(0, 2), sqrt2 * rot(0, 0) * rot(0, 1),
545 rot(1, 0) * rot(1, 0), rot(1, 1) * rot(1, 1), rot(1, 2) * rot(1, 2), sqrt2 * rot(1, 1) * rot(1, 2), sqrt2 * rot(1, 0) * rot(1, 2), sqrt2 * rot(1, 0) * rot(1, 1),
546 rot(2, 0) * rot(2, 0), rot(2, 1) * rot(2, 1), rot(2, 2) * rot(2, 2), sqrt2 * rot(2, 1) * rot(2, 2), sqrt2 * rot(2, 0) * rot(2, 2), sqrt2 * rot(2, 0) * rot(2, 1),
547 sqrt2 * rot(1, 0) * rot(2, 0), sqrt2 * rot(1, 1) * rot(2, 1), sqrt2 * rot(1, 2) * rot(2, 2), rot(1, 1) * rot(2, 2) + rot(1, 2) * rot(2, 1), rot(1, 0) * rot(2, 2) + rot(1, 2) * rot(2, 0), rot(1, 0) * rot(2, 1) + rot(1, 1) * rot(2, 0),
548 sqrt2 * rot(0, 0) * rot(2, 0), sqrt2 * rot(0, 1) * rot(2, 1), sqrt2 * rot(0, 2) * rot(2, 2), rot(0, 1) * rot(2, 2) + rot(0, 2) * rot(2, 1), rot(0, 0) * rot(2, 2) + rot(0, 2) * rot(2, 0), rot(0, 0) * rot(2, 1) + rot(0, 1) * rot(2, 0),
549 sqrt2 * rot(0, 0) * rot(1, 0), sqrt2 * rot(0, 1) * rot(1, 1), sqrt2 * rot(0, 2) * rot(1, 2), rot(0, 1) * rot(1, 2) + rot(0, 2) * rot(1, 1), rot(0, 0) * rot(1, 2) + rot(0, 2) * rot(1, 0), rot(0, 0) * rot(1, 1) + rot(0, 1) * rot(1, 0);
550 }
551
552 return res;
553 }
554
555 void FiberDirection::add_multimaterial(const int index, const json &dir, const std::string &unit)
556 {
557 for (int i = dir_.size(); i <= index; ++i)
558 {
559 dir_.emplace_back();
560 }
561
562 if (dir.size() == 3 || dir.size() == 2)
563 {
564 const int size = dir.size();
565 assert(size == size_);
566 dir_[index].resize(size, size);
567 for (int i = 0; i < size; ++i)
568 {
569 if (!dir[i].is_array() || dir[i].size() != size)
570 {
571 log_and_throw_error(fmt::format("Fiber must be {}x{}, row {} is {}", size, size, i, dir[i].dump()));
572 }
573 for (int j = 0; j < size; ++j)
574 {
575 dir_[index](i, j).init(dir[i][j]);
576 dir_[index](i, j).set_unit_type(unit);
577 }
578 }
579 has_rotation_ = true;
580 }
581 else if (dir.size() == 9 || dir.size() == 4)
582 {
583 const int size = dir.size() == 9 ? 3 : 2;
584 assert(size == size_);
585 dir_[index].resize(size, size);
586 for (int i = 0; i < size; ++i)
587 {
588 for (int j = 0; j < size; ++j)
589 {
590 dir_[index](i, j).init(dir[i * size + j]);
591 dir_[index](i, j).set_unit_type(unit);
592 }
593 }
594 has_rotation_ = true;
595 }
596 else if (dir.empty())
597 {
598 dir_[index].resize(size_, size_);
599 for (int i = 0; i < size_; ++i)
600 {
601 for (int j = 0; j < size_; ++j)
602 {
603 dir_[index](i, j).init(i == j ? 1.0 : 0.0);
604 dir_[index](i, j).set_unit_type(unit);
605 }
606 }
607 has_rotation_ = false;
608 }
609 else
610 {
611 log_and_throw_error("Fiber direction must be a 3x3 or 2x2 matrix");
612 }
613 }
614
615 // template instantiation
616 template double ElasticityTensor::compute_stress<3>(const std::array<double, 3> &strain, const int j) const;
617 template double ElasticityTensor::compute_stress<6>(const std::array<double, 6> &strain, const int j) const;
618
619} // namespace polyfem::assembler
std::vector< Eigen::Triplet< double > > entries
int y
int z
int x
std::vector< utils::ExpressionValue > rho_
virtual void add_multimaterial(const int index, const json &params, const std::string &density_unit)
virtual double operator()(double px, double py, double pz, double x, double y, double z, double t, int el_id) const
double operator()(int i, int j) const
void rotate_stiffness(const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 6, 6 > &rotation_mtx_voigt)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 6, 6 > stiffness_tensor_
Definition MatParams.hpp:70
void set_from_entries(const std::vector< double > &entries, const std::string &stress_unit)
void set_from_young_poisson(const double young, const double poisson, const std::string &stress_unit)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 6, 6 > reference_stiffness_tensor_
Definition MatParams.hpp:71
void set_orthotropic(double Ex, double Ey, double Ez, double nuXY, double nuXZ, double nuYZ, double muYZ, double muZX, double muXY, const std::string &stress_unit)
double compute_stress(const std::array< double, DIM > &strain, const int j) const
void set_transversely_isotropic(double Et, double Ea, double nu_t, double nu_a, double Ga, const std::string &stress_units)
void set_from_lambda_mu(const double lambda, const double mu, const std::string &stress_unit)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 1, 6, 6 > stiffness_rotation_voigt(double px, double py, double pz, double x, double y, double z, double t, int el_id) const
std::vector< Eigen::Matrix< utils::ExpressionValue, Eigen::Dynamic, Eigen::Dynamic, 1, 3, 3 > > dir_
void add_multimaterial(const int index, const json &params, const std::string &unit)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 1, 3, 3 > operator()(double px, double py, double pz, double x, double y, double z, double t, int el_id) const
void add_multimaterial(const int index, const json &params, const std::string &unit_type)
Definition MatParams.cpp:30
double operator()(const RowVectorNd &p, double t, int index) const
Definition MatParams.cpp:44
GenericMatParam(const std::string &param_name)
Definition MatParams.cpp:25
std::vector< utils::ExpressionValue > param_
Definition MatParams.hpp:21
void add_multimaterial(const int index, const json &params, const std::string &unit_type)
Definition MatParams.cpp:67
std::vector< GenericMatParam > params_
Definition MatParams.hpp:38
GenericMatParams(const std::string &param_name)
Definition MatParams.cpp:62
std::vector< utils::ExpressionValue > mu_or_nu_
void lambda_mu(double px, double py, double pz, double x, double y, double z, double t, int el_id, double &lambda, double &mu) const
void set_e_nu(const int index, const json &E, const json &nu, const std::string &stress_unit)
std::vector< utils::ExpressionValue > lambda_or_E_
void add_multimaterial(const int index, const json &params, const bool is_volume, const std::string &stress_unit)
Used for test only.
std::vector< T > json_as_array(const json &j)
Return the value of a json object as an array.
Definition JSONUtils.hpp:38
nlohmann::json json
Definition Common.hpp:9
double convert_to_mu(const double E, const double nu)
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13
double convert_to_lambda(const bool is_volume, const double E, const double nu)
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:73