suanPan
tensorToolbox.h
Go to the documentation of this file.
1 /*******************************************************************************
2  * Copyright (C) 2017-2022 Theodore Chang
3  *
4  * This program is free software: you can redistribute it and/or modify
5  * it under the terms of the GNU General Public License as published by
6  * the Free Software Foundation, either version 3 of the License, or
7  * (at your option) any later version.
8  *
9  * This program is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with this program. If not, see <http://www.gnu.org/licenses/>.
16  ******************************************************************************/
17 
18 #ifndef TENSORTOOLBOX_H
19 #define TENSORTOOLBOX_H
20 
21 #include <Toolbox/utility.h>
22 
23 template<typename T> class Quaternion;
24 
25 namespace tensor {
26  mat isotropic_stiffness(double, double);
27  mat orthotropic_stiffness(const vec&, const vec&);
28 
32 
33  static const vec unit_tensor2{1., 1., 1., 0., 0., 0.};
34 
35  namespace stress {
36  // applies to 3D tensor only, either principal or not
37  double invariant1(const vec&);
38  // applies to 3D tensor only, either principal or not
39  double invariant2(const vec&);
40  // applies to 3D tensor only, either principal or not
41  double invariant3(const vec&);
42 
43  // compute lode angle based on input stress
44  double lode(vec);
45  // compute derivative of lode angle based on input stress
46  vec lode_der(vec);
47 
48  static const vec norm_weight{1., 1., 1., 2., 2., 2.};
49  } // namespace stress
50  namespace strain {
51  // applies to 3D tensor only, either principal or not
52  double invariant1(const vec&);
53  // applies to 3D tensor only, either principal or not
54  double invariant2(const vec&);
55  // applies to 3D tensor only, either principal or not
56  double invariant3(const vec&);
57 
58  // compute load angle based on input deviatoric strain
59  double lode(vec);
60 
61  static const vec norm_weight{1., 1., 1., .5, .5, .5};
62  } // namespace strain
63  double trace2(const vec&);
64  double trace3(const vec&);
65  double mean3(const vec&);
66  vec dev(const vec&);
67  vec dev(vec&&);
68 
69  mat dev(const mat&);
70  mat dev(mat&&);
71 
72  namespace strain {
73  mat to_green(mat&&);
74  mat to_green(const mat&);
75  mat to_tensor(const vec&);
76  vec to_voigt(const mat&);
77  double norm(const vec&);
78  double norm(vec&&);
79  double double_contraction(const vec&, const vec&);
80  double double_contraction(vec&&, vec&&);
81  } // namespace strain
82  namespace stress {
83  mat to_tensor(const vec&);
84  vec to_voigt(const mat&);
85  double norm(const vec&);
86  double norm(vec&&);
87  double double_contraction(const vec&, const vec&);
88  double double_contraction(vec&&, vec&&);
89  } // namespace stress
90 
91 } // namespace tensor
92 
93 namespace transform {
94  double atan2(const vec&);
97 
98  template<typename T> Mat<T> skew_symm(const Mat<T>& R) {
99  suanpan_debug([&] { if(R.n_elem != 3) throw invalid_argument("need 3 element vector"); });
100 
101  Mat<T> S(3, 3, fill::zeros);
102 
103  S(0, 1) = -(S(1, 0) = R(2));
104  S(2, 0) = -(S(0, 2) = R(1));
105  S(1, 2) = -(S(2, 1) = R(0));
106 
107  return S;
108  }
109 
110  template<typename T> Mat<T> skew_symm(const subview_col<T>& R) { return skew_symm(R.eval()); }
111 
112  template<typename T> Mat<T> rodrigues(const Mat<T>& R) { return arma::expmat(transform::skew_symm(R)); }
113 
114  template<typename T> Mat<T> rodrigues(const subview_col<T>& R) { return arma::expmat(transform::skew_symm(R)); }
115 
116  template<typename T> Quaternion<T> to_quaternion(const Mat<T>& R) {
117  if(3 == R.n_elem) {
118  const auto angle = arma::norm(R);
119 
120  if(suanpan::approx_equal(angle, 0.)) return {0., 0., 0., 0.};
121 
122  return {std::cos(.5 * angle), std::sin(.5 * angle) / angle * R};
123  }
124 
125  if(9 == R.n_elem) {
126  const auto tr_r = arma::trace(R);
127  const auto max_r = arma::max(R.diag());
128  if(tr_r >= max_r) {
129  const auto q0 = .5 * std::sqrt(tr_r + 1.);
130  const auto q1 = .25 / q0 * (R(2, 1) - R(1, 2));
131  const auto q2 = .25 / q0 * (R(0, 2) - R(2, 0));
132  const auto q3 = .25 / q0 * (R(1, 0) - R(0, 1));
133  return {q0, q1, q2, q3};
134  }
135 
136  for(auto I = 0; I < 3; ++I)
137  if(suanpan::approx_equal(R(I, I), max_r)) {
138  const auto J = (I + 1) % 3;
139  const auto K = (J + 1) % 3;
140  vec q(3, fill::none);
141  q(I) = std::sqrt(.5 * max_r + .25 * (1. - tr_r));
142  const double q0 = .25 / q(I) * (R(K, J) - R(J, K));
143  q(J) = .25 / q(I) * (R(J, I) + R(I, J));
144  q(K) = .25 / q(I) * (R(K, I) + R(I, K));
145 
146  return {q0, std::move(q)};
147  }
148  }
149 
150  throw invalid_argument("need either rotation vector or matrix");
151  }
152 
153  template<typename T> Col<T> to_pseudo(const Mat<T>& R) {
154  const Mat<T> S = arma::real(arma::logmat(R));
155 
156  return {S(2, 1), S(0, 2), S(1, 0)};
157  }
158 
159  namespace strain {
160  double angle(const vec&);
161  mat trans(double);
162  vec principal(const vec&);
163  vec rotate(const vec&, double);
164  } // namespace strain
165  namespace stress {
166  double angle(const vec&);
167  mat trans(double);
168  vec principal(const vec&);
169  vec rotate(const vec&, double);
170  } // namespace stress
171  namespace beam {
172  mat global_to_local(double, double, double);
173  mat global_to_local(const vec&, double);
174  } // namespace beam
175  namespace triangle {
176  vec to_area_coordinate(const vec&, const mat&);
177  }
178 
179 } // namespace transform
180 
181 namespace suanpan {
182  template<typename T> T ramp(const T in) { return in > T(0) ? in : T(0); }
183 } // namespace suanpan
184 
185 #endif
An Quaternion class.
Definition: Quaternion.hpp:34
Mat< T > stress(T X, T Y, unsigned S)
Definition: shapeFunction.h:475
Mat< T > strain(T X, T Y, T V, unsigned S)
Definition: shapeFunction.h:523
Col< T > beam(T int_pts, unsigned order, double length)
Definition: shapeFunction.h:151
T triangle(const Mat< T > &EC)
Definition: shapeFunction.h:107
Definition: MatrixModifier.hpp:36
T ramp(const T in)
Definition: tensorToolbox.h:182
std::enable_if_t<!std::numeric_limits< T >::is_integer, bool > approx_equal(T x, T y, int ulp=2)
Definition: utility.h:46
mat to_green(mat &&)
Definition: tensorToolbox.cpp:224
mat to_tensor(const vec &)
Definition: tensorToolbox.cpp:248
double invariant3(const vec &)
compute the third invariant of the given 3D strain tensor, could be either normal or deviatoric strai...
Definition: tensorToolbox.cpp:105
double norm(const vec &)
Definition: tensorToolbox.cpp:302
double invariant2(const vec &)
compute the second invariant of the given 3D strain tensor, could be either normal or deviatoric stra...
Definition: tensorToolbox.cpp:93
double invariant1(const vec &)
compute the first invariant of the given 3D strain tensor, could be either normal or deviatoric strai...
Definition: tensorToolbox.cpp:82
double double_contraction(const vec &, const vec &)
Definition: tensorToolbox.cpp:314
vec to_voigt(const mat &)
Definition: tensorToolbox.cpp:273
double lode(vec)
Definition: tensorToolbox.cpp:147
vec to_voigt(const mat &)
Definition: tensorToolbox.cpp:343
double invariant2(const vec &)
compute the second invariant of the given 3D stress tensor, could be either normal or deviatoric stre...
Definition: tensorToolbox.cpp:128
double double_contraction(const vec &, const vec &)
Definition: tensorToolbox.cpp:384
vec lode_der(vec)
Definition: tensorToolbox.cpp:165
double lode(vec)
Definition: tensorToolbox.cpp:156
double invariant3(const vec &)
compute the third invariant of the given 3D stress tensor, could be either normal or deviatoric stres...
Definition: tensorToolbox.cpp:140
double invariant1(const vec &)
compute the first invariant of the given 3D stress tensor, could be either normal or deviatoric stres...
Definition: tensorToolbox.cpp:117
double norm(const vec &)
Definition: tensorToolbox.cpp:372
mat to_tensor(const vec &)
Definition: tensorToolbox.cpp:318
double norm(vec &&)
Definition: tensorToolbox.cpp:378
Definition: tensorToolbox.h:25
double trace3(const vec &)
Only accepts 3D tensor!
Definition: tensorToolbox.cpp:189
mat unit_deviatoric_tensor4v2()
Definition: tensorToolbox.cpp:60
mat unit_deviatoric_tensor4()
Definition: tensorToolbox.cpp:50
double mean3(const vec &)
Definition: tensorToolbox.cpp:195
vec dev(const vec &)
Definition: tensorToolbox.cpp:197
mat unit_symmetric_tensor4()
Definition: tensorToolbox.cpp:68
mat isotropic_stiffness(double, double)
Definition: tensorToolbox.cpp:20
mat orthotropic_stiffness(const vec &, const vec &)
Definition: tensorToolbox.cpp:34
double trace2(const vec &)
Only accepts 2D tensor!
Definition: tensorToolbox.cpp:178
mat global_to_local(double, double, double)
Definition: tensorToolbox.cpp:541
vec principal(const vec &)
Definition: tensorToolbox.cpp:493
mat trans(double)
Definition: tensorToolbox.cpp:479
double angle(const vec &)
Definition: tensorToolbox.cpp:473
vec rotate(const vec &, double)
Definition: tensorToolbox.cpp:505
vec principal(const vec &)
Definition: tensorToolbox.cpp:527
mat trans(double)
Definition: tensorToolbox.cpp:513
double angle(const vec &)
Definition: tensorToolbox.cpp:511
vec rotate(const vec &, double)
Definition: tensorToolbox.cpp:539
vec to_area_coordinate(const vec &, const mat &)
Definition: tensorToolbox.cpp:458
Definition: tensorToolbox.h:93
Mat< T > rodrigues(const Mat< T > &R)
Definition: tensorToolbox.h:112
mat compute_jacobian_nominal_to_principal(const mat &)
Definition: tensorToolbox.cpp:390
Quaternion< T > to_quaternion(const Mat< T > &R)
Definition: tensorToolbox.h:116
Col< T > to_pseudo(const Mat< T > &R)
Definition: tensorToolbox.h:153
Mat< T > skew_symm(const Mat< T > &R)
Definition: tensorToolbox.h:98
double atan2(const vec &)
Definition: tensorToolbox.cpp:388
mat compute_jacobian_principal_to_nominal(const mat &)
Definition: tensorToolbox.cpp:424
void suanpan_debug(const char *M,...)
Definition: print.cpp:64