https://mooseframework.inl.gov
FactorizedRankTwoTensor.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 "RankTwoTensor.h"
13 #include "SymmetricRankTwoTensor.h"
14 
15 // forward declarations
16 template <typename>
18 
34 template <typename T>
36 {
37 public:
39  typedef typename T::value_type value_type;
40 
43 
46 
49  {
50  if (!A.isSymmetric())
51  mooseError("The tensor is not symmetric.");
52  A.symmetricEigenvaluesEigenvectors(_eigvals, _eigvecs);
53  }
54 
55  // Construct from the factorization
56  // Assume that regardless of the type of T, the eigenvectors are always stored as as a full second
57  // order tensor, e.g. the eigenvector matrix of a symmetric second order tensor isn't symmetric in
58  // general. Hence we don't take advantage of orthonormal and/or unitary second order tensors.
59  FactorizedRankTwoTensorTempl(const std::vector<typename T::value_type> & eigvals,
62  {
63  }
64 
65  // @{ getters
66  template <typename T2 = T>
67  T2 get() const
68  {
69  return static_cast<T2>(assemble());
70  }
71  template <typename T2 = T>
72  T2 get()
73  {
74  return static_cast<T2>(assemble());
75  }
76  const std::vector<typename T::value_type> & eigvals() const { return _eigvals; }
77  std::vector<typename T::value_type> eigvals() { return _eigvals; }
80  // @}
81 
82  void print(std::ostream & stm = Moose::out) const;
83 
84  // Returns _A rotated by R.
87 
90 
91  // @{ Assignment operators
94  // @}
95 
97  FactorizedRankTwoTensorTempl<T> & operator*=(const typename T::value_type & a);
98 
100  template <typename T2>
101  FactorizedRankTwoTensorTempl<T> operator*(const T2 & a) const;
102 
104  FactorizedRankTwoTensorTempl<T> & operator/=(const typename T::value_type & a);
105 
107  template <typename T2>
108  FactorizedRankTwoTensorTempl<T> operator/(const T2 & a) const;
109 
111  bool operator==(const T & A) const;
112 
114  bool operator==(const FactorizedRankTwoTensorTempl<T> & A) const;
115 
118 
120  void addIa(const typename T::value_type & a);
121 
123  typename T::value_type trace() const;
124 
126  typename T::value_type det() const;
127 
128 private:
129  // Assemble the tensor from the factorization
131  {
133  return _eigvals[0] * R2T::selfOuterProduct(_eigvecs.column(0)) +
134  _eigvals[1] * R2T::selfOuterProduct(_eigvecs.column(1)) +
135  _eigvals[2] * R2T::selfOuterProduct(_eigvecs.column(2));
136  }
137 
138  // The eigen values of _A;
139  std::vector<typename T::value_type> _eigvals;
140 
141  // The eigen vectors of _A;
143 };
144 
145 namespace MathUtils
146 {
147 #define FactorizedRankTwoTensorOperatorMapBody(operator) \
148  { \
149  using std::log, std::exp, std::sqrt, std::cbrt, std::pow; \
150  std::vector<typename T::value_type> op_eigvals; \
151  for (const auto & eigval : A.eigvals()) \
152  op_eigvals.push_back(operator); \
153  return FactorizedRankTwoTensorTempl<T>(op_eigvals, A.eigvecs()); \
154  }
155 
156 #define FactorizedRankTwoTensorOperatorMapUnary(operatorname, operator) \
157  template <typename T> \
158  FactorizedRankTwoTensorTempl<T> operatorname(const FactorizedRankTwoTensorTempl<T> & A) \
159  { \
160  FactorizedRankTwoTensorOperatorMapBody(operator) \
161  } \
162  static_assert(true, "")
163 
164 #define FactorizedRankTwoTensorOperatorMapBinary(operatorname, operator) \
165  template <typename T, typename T2> \
166  FactorizedRankTwoTensorTempl<T> operatorname(const FactorizedRankTwoTensorTempl<T> & A, \
167  const T2 & arg) \
168  { \
169  if constexpr (libMesh::ScalarTraits<T2>::value) \
170  { \
171  FactorizedRankTwoTensorOperatorMapBody(operator) \
172  } \
173  } \
174  static_assert(true, "")
175 
176 #define FactorizedRankTwoTensorOperatorMapDerivativeBody(operator, derivative) \
177  { \
178  using std::log, std::exp, std::sqrt, std::cbrt, std::pow; \
179  std::vector<typename T::value_type> op_eigvals, op_derivs; \
180  for (const auto & eigval : A.eigvals()) \
181  { \
182  op_eigvals.push_back(operator); \
183  op_derivs.push_back(derivative); \
184  } \
185  \
186  RankFourTensorTempl<typename T::value_type> P, Gab, Gba; \
187  typedef RankTwoTensorTempl<typename T::value_type> R2T; \
188  \
189  for (auto a : make_range(3)) \
190  { \
191  const auto Ma = R2T::selfOuterProduct(A.eigvecs().column(a)); \
192  P += op_derivs[a] * Ma.outerProduct(Ma); \
193  } \
194  \
195  for (auto a : make_range(3)) \
196  for (auto b : make_range(a)) \
197  { \
198  const auto Ma = R2T::selfOuterProduct(A.eigvecs().column(a)); \
199  const auto Mb = R2T::selfOuterProduct(A.eigvecs().column(b)); \
200  \
201  usingTensorIndices(i_, j_, k_, l_); \
202  Gab = Ma.template times<i_, k_, j_, l_>(Mb) + Ma.template times<i_, l_, j_, k_>(Mb); \
203  Gba = Mb.template times<i_, k_, j_, l_>(Ma) + Mb.template times<i_, l_, j_, k_>(Ma); \
204  \
205  Real theta_ab; \
206  if (!MooseUtils::absoluteFuzzyEqual(A.eigvals()[a], A.eigvals()[b])) \
207  theta_ab = 0.5 * (op_eigvals[a] - op_eigvals[b]) / (A.eigvals()[a] - A.eigvals()[b]); \
208  else \
209  theta_ab = 0.25 * (op_derivs[a] + op_derivs[b]); \
210  \
211  P += theta_ab * (Gab + Gba); \
212  } \
213  return P; \
214  }
215 
216 #define FactorizedRankTwoTensorOperatorMapDerivativeUnary(derivativename, operator, derivative) \
217  template <typename T> \
218  RankFourTensorTempl<typename T::value_type> derivativename( \
219  const FactorizedRankTwoTensorTempl<T> & A) \
220  { \
221  FactorizedRankTwoTensorOperatorMapDerivativeBody(operator, derivative) \
222  } \
223  static_assert(true, "")
224 
225 #define FactorizedRankTwoTensorOperatorMapDerivativeBinary(derivativename, operator, derivative) \
226  template <typename T, typename T2> \
227  RankFourTensorTempl<typename T::value_type> derivativename( \
228  const FactorizedRankTwoTensorTempl<T> & A, const T2 & arg) \
229  { \
230  if constexpr (libMesh::ScalarTraits<T2>::value) \
231  { \
232  FactorizedRankTwoTensorOperatorMapDerivativeBody(operator, derivative) \
233  } \
234  } \
235  static_assert(true, "")
236 
237 // TODO: While the macros are here, in the future we could instantiate other operator maps like
238 // trignometry functions.
239 // @{
241 FactorizedRankTwoTensorOperatorMapDerivativeUnary(dlog, log(eigval), 1 / eigval);
242 
245 
247 FactorizedRankTwoTensorOperatorMapDerivativeUnary(dsqrt, sqrt(eigval), pow(eigval, -1. / 2.) / 2.);
248 
249 FactorizedRankTwoTensorOperatorMapUnary(cbrt, cbrt(eigval));
250 FactorizedRankTwoTensorOperatorMapDerivativeUnary(dcbrt, cbrt(eigval), pow(eigval, -2. / 3.) / 3.);
251 
254  pow(eigval, arg),
255  arg * pow(eigval, arg - 1));
256 // @}
257 } // end namespace MathUtils
258 
259 template <typename T>
260 template <typename T2>
263 {
264  if constexpr (libMesh::ScalarTraits<T2>::value)
265  {
267  for (auto & eigval : A._eigvals)
268  eigval *= a;
269  return A;
270  }
271 }
272 
273 template <typename T>
274 template <typename T2>
277 {
278  if constexpr (libMesh::ScalarTraits<T2>::value)
279  {
281  for (auto & eigval : A._eigvals)
282  eigval /= a;
283  return A;
284  }
285 }
286 
T::value_type trace() const
trace of _A
std::vector< typename T::value_type > eigvals()
T::value_type value_type
For generic programming.
FactorizedRankTwoTensorTempl< T > & operator/=(const typename T::value_type &a)
performs _A /= a in place, also updates eigen values
FactorizedRankTwoTensorTempl< T > & operator*=(const typename T::value_type &a)
performs _A *= a in place, also updates eigen values
auto exp(const T &)
FactorizedRankTwoTensorTempl< T > rotated(const RankTwoTensorTempl< typename T::value_type > &R) const
T::value_type det() const
determinant of _A
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
const RankTwoTensorTempl< typename T::value_type > & eigvecs() const
FactorizedRankTwoTensorOperatorMapUnary(log, log(eigval))
RankTwoTensorTempl< typename T::value_type > _eigvecs
std::vector< typename T::value_type > _eigvals
FactorizedRankTwoTensorOperatorMapDerivativeUnary(dlog, log(eigval), 1/eigval)
FactorizedRankTwoTensorTempl< T > operator*(const T2 &a) const
returns _A * a, also updates eigen values
FactorizedRankTwoTensorOperatorMapDerivativeBinary(dpow, pow(eigval, arg), arg *pow(eigval, arg - 1))
FactorizedRankTwoTensorTempl< ADSymmetricRankTwoTensor > ADFactorizedSymmetricRankTwoTensor
FactorizedRankTwoTensorTempl(const std::vector< typename T::value_type > &eigvals, const RankTwoTensorTempl< typename T::value_type > &eigvecs)
FactorizedRankTwoTensorTempl< T > transpose() const
Returns the transpose of _A.
FactorizedRankTwoTensorOperatorMapBinary(pow, pow(eigval, arg))
FactorizedRankTwoTensorTempl< SymmetricRankTwoTensor > FactorizedSymmetricRankTwoTensor
FactorizedRankTwoTensorTempl is designed to perform the spectral decomposition of an underlying symme...
RankTwoTensorTempl< typename T::value_type > eigvecs()
libMesh::VectorValue< T > column(const unsigned int i) const
Get the i-th column of the second order tensor.
const std::vector< typename T::value_type > & eigvals() const
auto log(const T &)
FactorizedRankTwoTensorTempl< RankTwoTensor > FactorizedRankTwoTensor
FactorizedRankTwoTensorTempl()=delete
No default constructor.
void print(std::ostream &stm=Moose::out) const
void addIa(const typename T::value_type &a)
add identity times a to _A
FactorizedRankTwoTensorTempl(const T &A)
Constructor if the factorization isn&#39;t known a priori.
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sqrt(_arg)) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tanh
FactorizedRankTwoTensorTempl< T > & operator=(const FactorizedRankTwoTensorTempl< T > &A)
bool operator==(const T &A) const
Defines logical equality with another second order tensor.
RankTwoTensorTempl< typename T::value_type > assemble() const
T pow(T x, int e)
Definition: MathUtils.h:91
FactorizedRankTwoTensorTempl< T > inverse() const
inverse of _A
FactorizedRankTwoTensorTempl< T > operator/(const T2 &a) const
returns _A / a, also updates eigen values
FactorizedRankTwoTensorTempl< ADRankTwoTensor > ADFactorizedRankTwoTensor