https://mooseframework.inl.gov
SymmetricRankFourTensor.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "Moose.h"
13 #include "ADRankTwoTensorForward.h"
17 #include "MooseUtils.h"
18 #include "MathUtils.h"
19 
20 #include "libmesh/libmesh.h"
21 #include "libmesh/tuple_of.h"
22 #include "libmesh/int_range.h"
23 
24 #include "metaphysicl/raw_type.h"
25 
26 #include <petscsys.h>
27 
28 // Eigen includes
29 #include <Eigen/Core>
30 #include <Eigen/Dense>
31 
32 #include <array>
33 
34 using libMesh::Real;
35 using libMesh::tuple_of;
36 namespace libMesh
37 {
38 template <typename>
39 class TensorValue;
40 template <typename>
41 class TypeTensor;
42 template <typename>
43 class VectorValue;
44 }
45 
46 // Forward declarations
47 class MooseEnum;
48 
49 namespace boostcopy = libMesh::boostcopy;
50 
51 namespace MathUtils
52 {
57 template <>
59 template <>
61 }
62 
68 template <typename T>
70 {
71 public:
73  static constexpr unsigned int Ndim = LIBMESH_DIM;
74  static constexpr unsigned int N = Ndim + Ndim * (Ndim - 1) / 2;
75  static constexpr unsigned int N2 = N * N;
77 
78  // Full tensor indices in the Mandel/Voigt representation
79  static constexpr unsigned int full_index[6][6][4] = {
80  {{0, 0, 0, 0}, {0, 0, 1, 1}, {0, 0, 2, 2}, {0, 0, 1, 2}, {0, 0, 0, 2}, {0, 0, 0, 1}},
81  {{1, 1, 0, 0}, {1, 1, 1, 1}, {1, 1, 2, 2}, {1, 1, 1, 2}, {1, 1, 0, 2}, {1, 1, 0, 1}},
82  {{2, 2, 0, 0}, {2, 2, 1, 1}, {2, 2, 2, 2}, {2, 2, 1, 2}, {2, 2, 0, 2}, {2, 2, 0, 1}},
83  {{1, 2, 0, 0}, {1, 2, 1, 1}, {1, 2, 2, 2}, {1, 2, 1, 2}, {1, 2, 0, 2}, {1, 2, 0, 1}},
84  {{0, 2, 0, 0}, {0, 2, 1, 1}, {0, 2, 2, 2}, {0, 2, 1, 2}, {0, 2, 0, 2}, {0, 2, 0, 1}},
85  {{0, 1, 0, 0}, {0, 1, 1, 1}, {0, 1, 2, 2}, {0, 1, 1, 2}, {0, 1, 0, 2}, {0, 1, 0, 1}}};
86 
88  static constexpr Real mandelFactor(unsigned int i, unsigned int j)
89  {
92  }
93 
96  {
100  };
101 
108  {
109  symmetric9, // fillSymmetric9FromInputVector
110  symmetric21, // fillSymmetric21FromInputVector
111  symmetric_isotropic, // fillSymmetricIsotropicFromInputVector
112  symmetric_isotropic_E_nu, // fillSymmetricIsotropicEandNu
113  axisymmetric_rz, // fillAxisymmetricRZFromInputVector
114  principal, // fillPrincipalFromInputVector
115  orthotropic // fillGeneralOrthotropicFromInputVector
116  };
117 
118  template <template <typename> class Tensor, typename Scalar>
120  {
121  static const bool value = false;
122  };
123  template <typename Scalar>
125  {
127  };
128  template <typename Scalar>
129  struct TwoTensorMultTraits<TensorValue, Scalar>
130  {
132  };
133  template <typename Scalar>
134  struct TwoTensorMultTraits<TypeTensor, Scalar>
135  {
137  };
138 
141 
144 
146  SymmetricRankFourTensorTempl(const std::vector<T> &, FillMethod);
147 
150 
154  template <typename T2>
156 
159 
161  explicit operator RankFourTensorTempl<T>();
162 
163  // Named constructors
165  {
167  }
169  {
171  };
172 
174  inline T & operator()(unsigned int i, unsigned int j) { return _vals[i * N + j]; }
175 
180  inline const T & operator()(unsigned int i, unsigned int j) const { return _vals[i * N + j]; }
181 
183  void zero();
184 
186  void print(std::ostream & stm = Moose::out) const;
187 
189  void printReal(std::ostream & stm = Moose::out) const;
190 
191  friend std::ostream & operator<<(std::ostream & os, const SymmetricRankFourTensorTempl<T> & t)
192  {
193  t.print(os);
194  return os;
195  }
196 
199 
205  template <typename Scalar>
208  operator=(const Scalar & libmesh_dbg_var(p))
209  {
210  libmesh_assert_equal_to(p, Scalar(0));
211  this->zero();
212  return *this;
213  }
214 
216  template <typename T2>
217  auto operator*(const SymmetricRankTwoTensorTempl<T2> & b) const
219 
221  template <typename T2>
222  auto operator*(const T2 & a) const ->
223  typename std::enable_if<libMesh::ScalarTraits<T2>::value,
225 
228 
230  template <typename T2>
231  auto operator/(const T2 & a) const ->
232  typename std::enable_if<libMesh::ScalarTraits<T2>::value,
233  SymmetricRankFourTensorTempl<decltype(T() / T2())>>::type;
234 
237 
240 
242  template <typename T2>
243  auto operator+(const SymmetricRankFourTensorTempl<T2> & a) const
245 
248 
250  template <typename T2>
251  auto operator-(const SymmetricRankFourTensorTempl<T2> & a) const
252  -> SymmetricRankFourTensorTempl<decltype(T() - T2())>;
253 
256 
258  template <typename T2>
259  auto operator*(const SymmetricRankFourTensorTempl<T2> & a) const
261 
263  T L2norm() const;
264 
270 
277  static SymmetricRankFourTensorTempl<T> rotationMatrix(const TypeTensor<T> & R);
278 
283  void rotate(const TypeTensor<T> & R);
284 
291 
296 
298  static MooseEnum fillMethodEnum();
299 
306  void fillFromInputVector(const std::vector<T> & input, FillMethod fill_method);
307 
309  void fillSymmetricIsotropic(const T & i0, const T & i1);
310  void fillSymmetricIsotropicEandNu(const T & E, const T & nu);
312 
323  template <typename T2>
324  void fillSymmetric9FromInputVector(const T2 & input);
325 
334  template <typename T2>
335  void fillSymmetric21FromInputVector(const T2 & input);
336 
338  T sum3x3() const;
339 
342 
344  bool isSymmetric() const;
345 
347  bool isIsotropic() const;
348 
349 protected:
351  std::array<T, N2> _vals;
352 
360  void fillSymmetricIsotropicFromInputVector(const std::vector<T> & input);
361 
371  void fillSymmetricIsotropicEandNuFromInputVector(const std::vector<T> & input);
372 
380  void fillAxisymmetricRZFromInputVector(const std::vector<T> & input);
381 
396  void fillPrincipalFromInputVector(const std::vector<T> & input);
397 
404  void fillGeneralOrthotropicFromInputVector(const std::vector<T> & input);
405 
406  template <class T2>
407  friend void dataStore(std::ostream &, SymmetricRankFourTensorTempl<T2> &, void *);
408 
409  template <class T2>
410  friend void dataLoad(std::istream &, SymmetricRankFourTensorTempl<T2> &, void *);
411 
412  template <typename T2>
414  template <typename T2>
416  template <typename T2>
417  friend class RankThreeTensorTempl;
418 };
419 
420 namespace MetaPhysicL
421 {
422 template <typename T>
424 {
426 
428  {
429  value_type ret;
430  for (unsigned int i = 0; i < SymmetricRankFourTensorTempl<T>::N; ++i)
431  for (unsigned int j = 0; j < SymmetricRankFourTensorTempl<T>::N; ++j)
432  ret(i, j) = raw_value(in(i, j));
433 
434  return ret;
435  }
436 };
437 }
438 
439 template <typename T1, typename T2>
440 inline auto
441 operator*(const T1 & a, const SymmetricRankFourTensorTempl<T2> & b) ->
442  typename std::enable_if<libMesh::ScalarTraits<T1>::value,
444 {
445  return b * a;
446 }
447 
448 template <typename T>
449 template <typename T2>
452 {
453  for (const auto i : make_range(N2))
454  _vals[i] = copy._vals[i];
455 }
456 
457 template <typename T>
458 template <typename T2>
459 auto
461  typename std::enable_if<libMesh::ScalarTraits<T2>::value,
463 {
464  typedef decltype(T() * T2()) ValueType;
466 
467  for (const auto i : make_range(N2))
468  result._vals[i] = _vals[i] * b;
469 
470  return result;
471 }
472 
473 template <typename T>
474 template <typename T2>
475 auto
478 {
479  typedef decltype(T() * T2()) ValueType;
481 
482  std::size_t index = 0;
483  for (const auto i : make_range(N))
484  {
485  ValueType tmp = 0.0;
486  for (const auto j : make_range(N))
487  tmp += _vals[index++] * b._vals[j];
488  result._vals[i] = tmp;
489  }
490 
491  return result;
492 }
493 
494 template <typename T>
495 template <typename T2>
496 auto
498  typename std::enable_if<libMesh::ScalarTraits<T2>::value,
499  SymmetricRankFourTensorTempl<decltype(T() / T2())>>::type
500 {
501  SymmetricRankFourTensorTempl<decltype(T() / T2())> result;
502  for (const auto i : make_range(N2))
503  result._vals[i] = _vals[i] / b;
504  return result;
505 }
506 
507 template <typename T>
508 template <typename T2>
509 void
511 {
512  mooseAssert(LIBMESH_DIM == 3, "This method assumes LIBMESH_DIM == 3");
513  mooseAssert(input.size() == 9,
514  "To use fillSymmetric9FromInputVector, your input must have size 9.");
515  zero();
516 
517  _vals[0] = input[0]; // C1111
518  _vals[7] = input[3]; // C2222
519  _vals[14] = input[5]; // C3333
520 
521  _vals[1] = _vals[6] = input[1]; // C1122
522  _vals[2] = _vals[12] = input[2]; // C1133
523  _vals[8] = _vals[13] = input[4]; // C2233
524 
525  static constexpr std::size_t C2323 = 21;
526  static constexpr std::size_t C1313 = 28;
527  static constexpr std::size_t C1212 = 35;
528 
529  _vals[C2323] = 2.0 * input[6];
530  _vals[C1313] = 2.0 * input[7];
531  _vals[C1212] = 2.0 * input[8];
532 }
533 
534 template <typename T>
535 template <typename T2>
536 void
538 {
539  // C1111 C1122 C1133 C1123 C1113 C1112
540  // C2222 C2233 C2223 C2213 C2212
541  // C3333 C3323 C3313 C3312
542  // C2323 C2313 C2312
543  // C1313 C1312
544  // C1212
545 
546  mooseAssert(LIBMESH_DIM == 3, "This method assumes LIBMESH_DIM == 3");
547  mooseAssert(input.size() == 21,
548  "To use fillSymmetric21FromInputVector, your input must have size 21.");
549  std::size_t index = 0;
550  for (const auto i : make_range(N))
551  for (const auto j : make_range(i, N))
552  {
553  _vals[i + N * j] = mandelFactor(i, j) * input[index];
554  _vals[j + N * i] = mandelFactor(j, i) * input[index];
555  index++;
556  }
557 }
558 
559 template <typename T>
562 {
563  SymmetricRankFourTensorTempl<T> result(initNone);
564 
565  if constexpr (SymmetricRankFourTensorTempl<T>::N2 * sizeof(T) > EIGEN_STACK_ALLOCATION_LIMIT)
566  {
567  Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic> mat(N, N);
568  for (const auto i : make_range(N))
569  for (const auto j : make_range(N))
570  mat(i, j) = (*this)(i, j);
571  mat = mat.inverse();
572  for (const auto i : make_range(N))
573  for (const auto j : make_range(N))
574  result(i, j) = mat(i, j);
575  }
576  else
577  {
578  const Eigen::Map<const Eigen::Matrix<T, N, N, Eigen::RowMajor>> mat(&_vals[0]);
579  Eigen::Map<Eigen::Matrix<T, N, N, Eigen::RowMajor>> res(&result._vals[0]);
580  res = mat.inverse();
581  }
582 
583  return result;
584 }
void fillAxisymmetricRZFromInputVector(const std::vector< T > &input)
fillAxisymmetricRZFromInputVector takes 5 inputs to fill the axisymmetric Rank-4 tensor with the appr...
RankFourTensorTempl is designed to handle any N-dimensional fourth order tensor, C.
std::array< T, N2 > _vals
The values of the rank-four tensor.
libMesh::VectorValue< T > sum3x1() const
Calculates the vector a[i] = sum over j Ciijj for i and j varying from 0 to 2.
void fillSymmetricIsotropicEandNu(const T &E, const T &nu)
const T & operator()(unsigned int i, unsigned int j) const
Gets the value for the indices specified.
static MooseEnum fillMethodEnum()
Static method for use in validParams for getting the "fill_method".
T & operator()(unsigned int i, unsigned int j)
Gets the value for the indices specified. Takes indices ranging from 0-5 for i and j...
bool isIsotropic() const
checks if the tensor is isotropic
auto operator+(const SymmetricRankFourTensorTempl< T2 > &a) const -> SymmetricRankFourTensorTempl< decltype(T()+T2())>
C_ijkl + a_ijkl.
FillMethod
To fill up the 36 entries in the 4th-order tensor, fillFromInputVector is called with one of the foll...
void fillSymmetricIsotropic(const T &i0, const T &i1)
Vector-less fill API functions. See docs of the corresponding ...FromInputVector methods.
SymmetricRankFourTensorTempl< T > & operator-=(const SymmetricRankFourTensorTempl< T > &a)
C_ijkl -= a_ijkl.
void mooseSetToZero< SymmetricRankFourTensor >(SymmetricRankFourTensor &v)
Helper function template specialization to set an object to zero.
auto operator*(const SymmetricRankTwoTensorTempl< T2 > &b) const -> SymmetricRankTwoTensorTempl< decltype(T() *T2())>
C_ijkl*a_kl.
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
void fillSymmetric21FromInputVector(const T2 &input)
fillSymmetric21FromInputVector takes 21 inputs to fill in the Rank-4 tensor with the appropriate crys...
static constexpr Real mandelFactor(unsigned int i, unsigned int j)
returns the 1, sqrt(2), or 2 prefactor in the Mandel notation for the indices i,j ranging from 0-5...
We need to instantiate the following CompareTypes to tell the compiler that ADReal is a subtype of Ch...
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
std::basic_ostream< charT, traits > * os
Definition: InfixIterator.h:33
const Number zero
static constexpr unsigned int N
T sum3x3() const
Calculates the sum of Ciijj for i and j varying from 0 to 2.
static SymmetricRankFourTensorTempl< T > rotationMatrix(const TypeTensor< T > &R)
Build a 6x6 rotation matrix MEHRABADI, MORTEZA M.
RankThreeTensor is designed to handle any N-dimensional third order tensor, r.
static value_type value(const SymmetricRankFourTensorTempl< T > &in)
auto operator/(const T2 &a) const -> typename std::enable_if< libMesh::ScalarTraits< T2 >::value, SymmetricRankFourTensorTempl< decltype(T()/T2())>>::type
C_ijkl/a.
typename tuple_n< Index, T >::template type<> tuple_of
SymmetricRankFourTensorTempl< T > & operator+=(const SymmetricRankFourTensorTempl< T > &a)
C_ijkl += a_ijkl for all i, j, k, l.
SymmetricRankFourTensorTempl< typename RawType< T >::value_type > value_type
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...
void fillSymmetricIsotropicEandNuFromInputVector(const std::vector< T > &input)
fillSymmetricIsotropicEandNuFromInputVector is a variation of the fillSymmetricIsotropicFromInputVect...
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
void rotate(const TypeTensor< T > &R)
Rotate the tensor using C_ijkl = R_im R_jn R_ko R_lp C_mnop.
SymmetricRankFourTensorTempl< T > & operator*=(const T &a)
C_ijkl *= a.
static SymmetricRankFourTensorTempl< T > identitySymmetricFour()
void fillPrincipalFromInputVector(const std::vector< T > &input)
fillPrincipalFromInputVector takes 9 inputs to fill a Rank-4 tensor C1111 = input0 C1122 = input1 C11...
void printReal(std::ostream &stm=Moose::out) const
Print the values of the rank four tensor.
void fillSymmetricIsotropicFromInputVector(const std::vector< T > &input)
fillSymmetricIsotropicFromInputVector takes 2 inputs to fill the the symmetric Rank-4 tensor with the...
friend void dataStore(std::ostream &, SymmetricRankFourTensorTempl< T2 > &, void *)
void mooseSetToZero< ADSymmetricRankFourTensor >(ADSymmetricRankFourTensor &v)
SymmetricRankFourTensorTempl< T > transposeMajor() const
Transpose the tensor by swapping the first pair with the second pair of indices This amounts to a reg...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
SymmetricRankFourTensorTempl< T > & operator=(const SymmetricRankFourTensorTempl< T > &a)=default
copies values from a into this tensor
void print(std::ostream &stm=Moose::out) const
Print the rank four tensor.
void fillFromInputVector(const std::vector< T > &input, FillMethod fill_method)
fillFromInputVector takes some number of inputs to fill the Rank-4 tensor.
NumberTensorValue Tensor
SymmetricRankFourTensorTempl< T > operator-() const
-C_ijkl
friend void dataLoad(std::istream &, SymmetricRankFourTensorTempl< T2 > &, void *)
void fillGeneralOrthotropicFromInputVector(const std::vector< T > &input)
fillGeneralOrhotropicFromInputVector takes 10 inputs to fill the Rank-4 tensor It defines a general o...
static constexpr unsigned int full_index[6][6][4]
SymmetricRankFourTensorTempl is designed to handle an N-dimensional fourth order tensor with minor sy...
IntRange< T > make_range(T beg, T end)
boostcopy::enable_if_c< libMesh::ScalarTraits< Scalar >::value, SymmetricRankFourTensorTempl & >::type operator=(const Scalar &libmesh_dbg_var(p))
Assignment-from-scalar operator.
SymmetricRankFourTensorTempl< T > invSymm() const
This returns A_ijkl such that C_ijkl*A_klmn = 0.5*(de_im de_jn + de_in de_jm) This routine assumes th...
static SymmetricRankFourTensorTempl< T > identity()
bool isSymmetric() const
checks if the tensor is symmetric
static constexpr Real mandelFactor(unsigned int i)
returns the 1 or sqrt(2) prefactor in the Mandel notation for the index i ranging from 0-5...
static constexpr unsigned int Ndim
tensor dimension, Mandel matrix dimension, and Mandel matrix size
SymmetricRankFourTensorTempl< T > & operator/=(const T &a)
C_ijkl /= a for all i, j, k, l.
SymmetricRankFourTensorTempl()
Default constructor; fills to zero.
void fillSymmetric9FromInputVector(const T2 &input)
fillSymmetric9FromInputVector takes 9 inputs to fill in the Rank-4 tensor with the appropriate crysta...
auto operator*(const T1 &a, const SymmetricRankFourTensorTempl< T2 > &b) -> typename std::enable_if< libMesh::ScalarTraits< T1 >::value, SymmetricRankFourTensorTempl< decltype(T1() *T2())>>::type
SymmetricRankFourTensorTempl< T > transposeIj() const
Transpose the tensor by swapping the first two indices - a no-op.
static constexpr unsigned int N2