https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members
SymmetricRankTwoTensorTempl< T > Class Template Reference

SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic material. More...

#include <SymmetricRankTwoTensor.h>

Public Types

enum  InitMethod { initNone, initIdentity }
 
enum  FillMethod { autodetect = 0, isotropic1 = 1, diagonal3 = 3, symmetric6 = 6 }
 To fill up the 6 entries in the 2nd-order tensor, fillFromInputVector is called with one of the following fill_methods. More...
 
typedef T value_type
 For generic programming. More...
 

Public Member Functions

 SymmetricRankTwoTensorTempl ()
 Default constructor; fills to zero. More...
 
 SymmetricRankTwoTensorTempl (const InitMethod)
 Select specific initialization pattern. More...
 
 SymmetricRankTwoTensorTempl (const std::vector< T > &input)
 Constructor that proxies the fillFromInputVector method. More...
 
 SymmetricRankTwoTensorTempl (const T &S11, const T &S22, const T &S33, const T &S23, const T &S13, const T &S12)
 Initialization list replacement constructors, 6 arguments. More...
 
 operator RankTwoTensorTempl< T > ()
 
 SymmetricRankTwoTensorTempl (const TensorValue< T > &a)
 Copy constructor from TensorValue<T> More...
 
 SymmetricRankTwoTensorTempl (const TypeTensor< T > &a)
 Copy constructor from TypeTensor<T> More...
 
template<typename T2 >
 SymmetricRankTwoTensorTempl (const SymmetricRankTwoTensorTempl< T2 > &a)
 Construct from other template. More...
 
void fillFromInputVector (const std::vector< T > &input, FillMethod fill_method=autodetect)
 fillFromInputVector takes 1, 3, or 6 inputs to fill in the symmmetric Rank-2 tensor. More...
 
void fillFromScalarVariable (const VariableValue &scalar_variable)
 fillFromScalarVariable takes FIRST/THIRD/SIXTH order scalar variable to fill in the Rank-2 tensor. More...
 
T & operator() (unsigned int i)
 Gets the raw value for the index specified. Takes index = 0,1,2,3,4,5. More...
 
operator() (unsigned int i) const
 Gets the raw value for the index specified. More...
 
operator() (unsigned int i, unsigned int j) const
 Gets the value for the index specified. More...
 
libMesh::VectorValue< T > row (const unsigned int n) const
 get the specified row of the tensor More...
 
libMesh::VectorValue< T > column (const unsigned int n) const
 get the specified column of the tensor More...
 
SymmetricRankTwoTensorTempl< T > square () const
 Returns the matrix squared. More...
 
tr () const
 Returns the trace. More...
 
void zero ()
 Set all components to zero. More...
 
void rotate (const TypeTensor< T > &R)
 rotates the tensor data given a rank two tensor rotation tensor _vals[i][j] = R_ij * R_jl * _vals[k][l] More...
 
SymmetricRankTwoTensorTempl< T > transpose () const
 Returns a matrix that is the transpose of the matrix this was called on. More...
 
template<typename T2 >
SymmetricRankTwoTensorTempl< T > & operator= (const SymmetricRankTwoTensorTempl< T2 > &a)
 sets _vals to a, and returns _vals More...
 
template<typename Scalar >
boostcopy::enable_if_c< libMesh::ScalarTraits< Scalar >::value, SymmetricRankTwoTensorTempl & >::type operator= (const Scalar &libmesh_dbg_var(p))
 Assignment-from-scalar operator. More...
 
SymmetricRankTwoTensorTempl< T > & operator+= (const SymmetricRankTwoTensorTempl< T > &a)
 adds a to _vals More...
 
template<typename T2 >
SymmetricRankTwoTensorTempl< typename libMesh::CompareTypes< T, T2 >::supertype > operator+ (const SymmetricRankTwoTensorTempl< T2 > &a) const
 returns _vals + a More...
 
SymmetricRankTwoTensorTempl< T > & operator-= (const SymmetricRankTwoTensorTempl< T > &a)
 sets _vals -= a and returns vals More...
 
template<typename T2 >
SymmetricRankTwoTensorTempl< typename libMesh::CompareTypes< T, T2 >::supertype > operator- (const SymmetricRankTwoTensorTempl< T2 > &a) const
 returns _vals - a More...
 
SymmetricRankTwoTensorTempl< T > operator- () const
 returns -_vals More...
 
SymmetricRankTwoTensorTempl< T > & operator*= (const T &a)
 performs _vals *= a More...
 
template<typename T2 >
auto operator* (const T2 &a) const -> typename std::enable_if< libMesh::ScalarTraits< T2 >::value, SymmetricRankTwoTensorTempl< decltype(T() *T2())>>::type
 returns _vals*a More...
 
SymmetricRankTwoTensorTempl< T > & operator/= (const T &a)
 performs _vals /= a More...
 
template<typename T2 >
auto operator/ (const T2 &a) const -> typename std::enable_if< libMesh::ScalarTraits< T2 >::value, SymmetricRankTwoTensorTempl< decltype(T()/T2())>>::type
 returns _vals/a More...
 
template<typename T2 >
TypeVector< typename libMesh::CompareTypes< T, T2 >::supertype > operator* (const TypeVector< T2 > &a) const
 Defines multiplication with a vector to get a vector. More...
 
template<typename T2 >
bool operator== (const SymmetricRankTwoTensorTempl< T2 > &a) const
 Defines logical equality with another SymmetricRankTwoTensorTempl<T2> More...
 
bool isSymmetric () const
 Test for symmetry. Surprisingly this is always true. More...
 
template<typename T2 >
bool operator!= (const SymmetricRankTwoTensorTempl< T2 > &a) const
 Defines logical inequality with another SymmetricRankTwoTensorTempl<T2> More...
 
SymmetricRankTwoTensorTempl< T > & operator= (const ColumnMajorMatrixTempl< T > &a)
 Sets _vals to the values in a ColumnMajorMatrix (must be 3x3) More...
 
doubleContraction (const SymmetricRankTwoTensorTempl< T > &a) const
 returns _vals_ij * a_ij (sum on i, j) More...
 
SymmetricRankFourTensorTempl< T > outerProduct (const SymmetricRankTwoTensorTempl< T > &a) const
 returns C_ijkl = a_ij * b_kl More...
 
SymmetricRankFourTensorTempl< T > positiveProjectionEigenDecomposition (std::vector< T > &, RankTwoTensorTempl< T > &) const
 return positive projection tensor of eigen-decomposition More...
 
SymmetricRankTwoTensorTempl< T > deviatoric () const
 returns A_ij - de_ij*tr(A)/3, where A are the _vals More...
 
trace () const
 returns the trace of the tensor, ie _vals[i][i] (sum i = 0, 1, 2) More...
 
SymmetricRankTwoTensorTempl< T > inverse () const
 retuns the inverse of the tensor More...
 
SymmetricRankTwoTensorTempl< T > dtrace () const
 Denote the _vals[i][j] by A_ij, then this returns d(trace)/dA_ij. More...
 
generalSecondInvariant () const
 Calculates the second invariant (I2) of a tensor. More...
 
secondInvariant () const
 Denote the _vals[i][j] by A_ij, then S_ij = A_ij - de_ij*tr(A)/3 Then this returns (S_ij + S_ji)*(S_ij + S_ji)/8 Note the explicit symmeterisation. More...
 
SymmetricRankTwoTensorTempl< T > dsecondInvariant () const
 Denote the _vals[i][j] by A_ij, then this returns d(secondInvariant)/dA_ij. More...
 
SymmetricRankFourTensorTempl< T > d2secondInvariant () const
 Denote the _vals[i][j] by A_ij, then this returns d^2(secondInvariant)/dA_ij/dA_kl. More...
 
sin3Lode (const T &r0, const T &r0_value) const
 Sin(3*Lode_angle) If secondInvariant() <= r0 then return r0_value This is to gaurd against precision-loss errors. More...
 
SymmetricRankTwoTensorTempl< T > dsin3Lode (const T &r0) const
 d(sin3Lode)/dA_ij If secondInvariant() <= r0 then return zero This is to gaurd against precision-loss errors. More...
 
thirdInvariant () const
 Denote the _vals[i][j] by A_ij, then S_ij = A_ij - de_ij*tr(A)/3 Then this returns det(S + S.transpose())/2 Note the explicit symmeterisation. More...
 
SymmetricRankTwoTensorTempl< T > dthirdInvariant () const
 Denote the _vals[i][j] by A_ij, then this returns d(thirdInvariant()/dA_ij. More...
 
SymmetricRankFourTensorTempl< T > d2thirdInvariant () const
 Denote the _vals[i][j] by A_ij, then this returns d^2(thirdInvariant)/dA_ij/dA_kl. More...
 
det () const
 
SymmetricRankTwoTensorTempl< T > ddet () const
 Denote the _vals[i][j] by A_ij, then this returns d(det)/dA_ij. More...
 
void print (std::ostream &stm=Moose::out) const
 Print the rank two tensor. More...
 
void printReal (std::ostream &stm=Moose::out) const
 Print the Real part of the ADReal rank two tensor. More...
 
void printADReal (unsigned int nDual, std::ostream &stm=Moose::out) const
 Print the Real part of the ADReal rank two tensor along with its first nDual dual numbers. More...
 
void addIa (const T &a)
 Add identity times a to _vals. More...
 
L2norm () const
 Sqrt(_vals[i][j]*_vals[i][j]) More...
 
void surfaceFillFromInputVector (const std::vector< T > &input)
 sets _vals[0][0], _vals[0][1], _vals[1][0], _vals[1][1] to input, and the remainder to zero More...
 
void symmetricEigenvalues (std::vector< T > &eigvals) const
 computes eigenvalues, assuming tens is symmetric, and places them in ascending order in eigvals More...
 
void symmetricEigenvaluesEigenvectors (std::vector< T > &eigvals, RankTwoTensorTempl< T > &eigvecs) const
 computes eigenvalues and eigenvectors, assuming tens is symmetric, and places them in ascending order in eigvals. More...
 
SymmetricRankTwoTensorTempl< T > initialContraction (const SymmetricRankFourTensorTempl< T > &b) const
 returns this_ij * b_ijkl More...
 
void setToIdentity ()
 set the tensor to the identity matrix More...
 

Static Public Member Functions

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. More...
 
static SymmetricRankTwoTensorTempl identity ()
 
static SymmetricRankTwoTensorTempl initializeSymmetric (const TypeVector< T > &v0, const TypeVector< T > &v1, const TypeVector< T > &v2)
 named constructor for initializing symmetrically More...
 
static MooseEnum fillMethodEnum ()
 Static method for use in validParams for getting the "fill_method". More...
 
static SymmetricRankTwoTensorTempl< T > timesTranspose (const RankTwoTensorTempl< T > &)
 return the matrix multiplied with its transpose A*A^T (guaranteed symmetric) More...
 
static SymmetricRankTwoTensorTempl< T > timesTranspose (const SymmetricRankTwoTensorTempl< T > &)
 
static SymmetricRankTwoTensorTempl< T > transposeTimes (const RankTwoTensorTempl< T > &)
 return the matrix multiplied with its transpose A^T*A (guaranteed symmetric) More...
 
static SymmetricRankTwoTensorTempl< T > transposeTimes (const SymmetricRankTwoTensorTempl< T > &)
 
static SymmetricRankTwoTensorTempl< T > plusTranspose (const RankTwoTensorTempl< T > &)
 return the matrix plus its transpose A-A^T (guaranteed symmetric) More...
 
static SymmetricRankTwoTensorTempl< T > plusTranspose (const SymmetricRankTwoTensorTempl< T > &)
 
static void initRandom (unsigned int)
 This function initializes random seed based on a user-defined number. More...
 
static SymmetricRankTwoTensorTempl< T > genRandomSymmTensor (T, T)
 This function generates a random symmetric rank two tensor. More...
 
static SymmetricRankTwoTensorTempl< T > selfOuterProduct (const TypeVector< T > &)
 SymmetricRankTwoTensorTempl<T> from outer product of a vector with itself. More...
 

Static Public Attributes

static constexpr unsigned int full_index [6][2] = {{0, 0}, {1, 1}, {2, 2}, {1, 2}, {0, 2}, {0, 1}}
 
static constexpr unsigned int reverse_index [3][3] = {{0, 5, 4}, {5, 1, 3}, {4, 3, 2}}
 
static constexpr unsigned int Ndim = Moose::dim
 tensor dimension and Mandel vector length More...
 
static constexpr unsigned int N = Ndim + (Ndim * (Ndim - 1)) / 2
 

Protected Member Functions

void syev (const char *calculation_type, std::vector< T > &eigvals, std::vector< T > &a) const
 Uses the petscblaslapack.h LAPACKsyev_ routine to find, for symmetric _vals: (1) the eigenvalues (if calculation_type == "N") (2) the eigenvalues and eigenvectors (if calculation_type == "V") More...
 

Private Member Functions

 SymmetricRankTwoTensorTempl (const T &S11, const T &S21, const T &S31, const T &S12, const T &S22, const T &S32, const T &S13, const T &S23, const T &S33)
 Initialization list replacement constructors, 9 arguments (for internal use only) More...
 

Static Private Member Functions

static SymmetricRankTwoTensorTempl fromRawComponents (const T &S11, const T &S22, const T &S33, const T &S23, const T &S13, const T &S12)
 

Private Attributes

std::array< T, N_vals
 

Static Private Attributes

static constexpr std::array< Real, NidentityCoords = {{1, 1, 1, 0, 0, 0}}
 

Friends

template<class T2 >
class SymmetricRankFourTensorTempl
 
std::ostream & operator<< (std::ostream &os, const SymmetricRankTwoTensorTempl< T > &t)
 
template<class T2 >
void dataStore (std::ostream &, SymmetricRankTwoTensorTempl< T2 > &, void *)
 
template<class T2 >
void dataLoad (std::istream &, SymmetricRankTwoTensorTempl< T2 > &, void *)
 

Detailed Description

template<typename T>
class SymmetricRankTwoTensorTempl< T >

SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic material.

It is designed to reduce the redundancies of the Complete tensor classes for regular mechanics problems and to enable use of the Mandel notation.

Definition at line 75 of file SymmetricRankTwoTensor.h.

Member Typedef Documentation

◆ value_type

template<typename T>
typedef T SymmetricRankTwoTensorTempl< T >::value_type

For generic programming.

Definition at line 79 of file SymmetricRankTwoTensor.h.

Member Enumeration Documentation

◆ FillMethod

template<typename T>
enum SymmetricRankTwoTensorTempl::FillMethod

To fill up the 6 entries in the 2nd-order tensor, fillFromInputVector is called with one of the following fill_methods.

See the fill*FromInputVector functions for more details

Enumerator
autodetect 
isotropic1 
diagonal3 
symmetric6 

Definition at line 113 of file SymmetricRankTwoTensor.h.

◆ InitMethod

template<typename T>
enum SymmetricRankTwoTensorTempl::InitMethod
Enumerator
initNone 
initIdentity 

Definition at line 96 of file SymmetricRankTwoTensor.h.

Constructor & Destructor Documentation

◆ SymmetricRankTwoTensorTempl() [1/8]

Default constructor; fills to zero.

Definition at line 56 of file SymmetricRankTwoTensorImplementation.h.

Referenced by SymmetricRankTwoTensorTempl< Real >::identity().

57 {
58  std::fill(_vals.begin(), _vals.end(), 0.0);
59 }

◆ SymmetricRankTwoTensorTempl() [2/8]

template<typename T >
SymmetricRankTwoTensorTempl< T >::SymmetricRankTwoTensorTempl ( const InitMethod  init)

Select specific initialization pattern.

Definition at line 62 of file SymmetricRankTwoTensorImplementation.h.

63 {
64  switch (init)
65  {
66  case initNone:
67  break;
68 
69  case initIdentity:
70  setToIdentity();
71  break;
72 
73  default:
74  mooseError("Unknown SymmetricRankTwoTensorTempl initialization pattern.");
75  }
76 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
void setToIdentity()
set the tensor to the identity matrix

◆ SymmetricRankTwoTensorTempl() [3/8]

template<typename T>
SymmetricRankTwoTensorTempl< T >::SymmetricRankTwoTensorTempl ( const std::vector< T > &  input)
inline

Constructor that proxies the fillFromInputVector method.

Definition at line 122 of file SymmetricRankTwoTensor.h.

122 { this->fillFromInputVector(input); };
void fillFromInputVector(const std::vector< T > &input, FillMethod fill_method=autodetect)
fillFromInputVector takes 1, 3, or 6 inputs to fill in the symmmetric Rank-2 tensor.

◆ SymmetricRankTwoTensorTempl() [4/8]

template<typename T>
SymmetricRankTwoTensorTempl< T >::SymmetricRankTwoTensorTempl ( const T &  S11,
const T &  S22,
const T &  S33,
const T &  S23,
const T &  S13,
const T &  S12 
)

Initialization list replacement constructors, 6 arguments.

Definition at line 79 of file SymmetricRankTwoTensorImplementation.h.

81 {
82  _vals[0] = S11;
83  _vals[1] = S22;
84  _vals[2] = S33;
85  _vals[3] = mandelFactor(3) * S23;
86  _vals[4] = mandelFactor(4) * S13;
87  _vals[5] = mandelFactor(5) * S12;
88 }
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...

◆ SymmetricRankTwoTensorTempl() [5/8]

template<typename T>
SymmetricRankTwoTensorTempl< T >::SymmetricRankTwoTensorTempl ( const T &  S11,
const T &  S21,
const T &  S31,
const T &  S12,
const T &  S22,
const T &  S32,
const T &  S13,
const T &  S23,
const T &  S33 
)
private

Initialization list replacement constructors, 9 arguments (for internal use only)

Definition at line 105 of file SymmetricRankTwoTensorImplementation.h.

114 {
115  _vals[0] = S11;
116  _vals[1] = S22;
117  _vals[2] = S33;
118  _vals[3] = (S23 + S32) / mandelFactor(3);
119  _vals[4] = (S13 + S31) / mandelFactor(4);
120  _vals[5] = (S12 + S21) / mandelFactor(5);
121 }
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...

◆ SymmetricRankTwoTensorTempl() [6/8]

template<typename T>
SymmetricRankTwoTensorTempl< T >::SymmetricRankTwoTensorTempl ( const TensorValue< T > &  a)
explicit

Copy constructor from TensorValue<T>

Definition at line 139 of file SymmetricRankTwoTensorImplementation.h.

140 {
141  _vals[0] = a(0, 0);
142  _vals[1] = a(1, 1);
143  _vals[2] = a(2, 2);
144  _vals[3] = (a(1, 2) + a(2, 1)) / mandelFactor(3);
145  _vals[4] = (a(0, 2) + a(2, 0)) / mandelFactor(4);
146  _vals[5] = (a(0, 1) + a(1, 0)) / mandelFactor(5);
147 }
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...

◆ SymmetricRankTwoTensorTempl() [7/8]

template<typename T>
SymmetricRankTwoTensorTempl< T >::SymmetricRankTwoTensorTempl ( const TypeTensor< T > &  a)
explicit

Copy constructor from TypeTensor<T>

Definition at line 150 of file SymmetricRankTwoTensorImplementation.h.

151 {
152  _vals[0] = a(0, 0);
153  _vals[1] = a(1, 1);
154  _vals[2] = a(2, 2);
155  _vals[3] = (a(1, 2) + a(2, 1)) / mandelFactor(3);
156  _vals[4] = (a(0, 2) + a(2, 0)) / mandelFactor(4);
157  _vals[5] = (a(0, 1) + a(1, 0)) / mandelFactor(5);
158 }
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...

◆ SymmetricRankTwoTensorTempl() [8/8]

template<typename T>
template<typename T2 >
SymmetricRankTwoTensorTempl< T >::SymmetricRankTwoTensorTempl ( const SymmetricRankTwoTensorTempl< T2 > &  a)
inline

Construct from other template.

Definition at line 156 of file SymmetricRankTwoTensor.h.

157  {
158  for (const auto i : make_range(N))
159  _vals[i] = a(i);
160  }
IntRange< T > make_range(T beg, T end)
static constexpr unsigned int N

Member Function Documentation

◆ addIa()

template<typename T>
void SymmetricRankTwoTensorTempl< T >::addIa ( const T &  a)

Add identity times a to _vals.

Definition at line 657 of file SymmetricRankTwoTensorImplementation.h.

Referenced by SymmetricRankTwoTensorTempl< Real >::deviatoric().

658 {
659  for (unsigned int i = 0; i < 3; ++i)
660  _vals[i] += a;
661 }

◆ column()

template<typename T>
libMesh::VectorValue<T> SymmetricRankTwoTensorTempl< T >::column ( const unsigned int  n) const
inline

get the specified column of the tensor

Definition at line 206 of file SymmetricRankTwoTensor.h.

206 { return row(n); }
libMesh::VectorValue< T > row(const unsigned int n) const
get the specified row of the tensor

◆ d2secondInvariant()

template<typename T >
SymmetricRankFourTensorTempl< T > SymmetricRankTwoTensorTempl< T >::d2secondInvariant ( ) const

Denote the _vals[i][j] by A_ij, then this returns d^2(secondInvariant)/dA_ij/dA_kl.

Definition at line 501 of file SymmetricRankTwoTensorImplementation.h.

502 {
504 
505  for (auto i : make_range(N))
506  for (auto j : make_range(N))
507  result(i, j) = 0.5 * (i == j) + 0.5 * (i == j) - (1.0 / 3.0) * (i < 3) * (j < 3);
508 
509  return result;
510 }
SymmetricRankFourTensorTempl is designed to handle an N-dimensional fourth order tensor with minor sy...
IntRange< T > make_range(T beg, T end)
static constexpr unsigned int N

◆ d2thirdInvariant()

template<typename T >
SymmetricRankFourTensorTempl< T > SymmetricRankTwoTensorTempl< T >::d2thirdInvariant ( ) const

Denote the _vals[i][j] by A_ij, then this returns d^2(thirdInvariant)/dA_ij/dA_kl.

Definition at line 571 of file SymmetricRankTwoTensorImplementation.h.

572 {
574 
576  for (auto i : make_range(N))
577  for (auto j : make_range(N))
578  d2(i, j) = Real(i < 3) * s(j) / 3.0 / (j < 3 ? 1 : MathUtils::sqrt2) +
579  Real(j < 3) * s(i) / 3.0 / (i < 3 ? 1 : MathUtils::sqrt2);
580 
581  d2(0, 1) += s(2);
582  d2(0, 2) += s(1);
583  d2(0, 3) -= s(3) / MathUtils::sqrt2;
584 
585  d2(1, 0) += s(2);
586  d2(1, 2) += s(0);
587  d2(1, 4) -= s(4) / MathUtils::sqrt2;
588 
589  d2(2, 0) += s(1);
590  d2(2, 1) += s(0);
591  d2(2, 5) -= s(5) / MathUtils::sqrt2;
592 
593  d2(3, 0) -= s(3) / MathUtils::sqrt2;
594  d2(3, 3) -= s(0) / 2.0;
595  d2(3, 4) += s(5) / 2.0 / MathUtils::sqrt2;
596  d2(3, 5) += s(4) / 2.0 / MathUtils::sqrt2;
597 
598  d2(4, 1) -= s(4) / MathUtils::sqrt2;
599  d2(4, 3) += s(5) / 2.0 / MathUtils::sqrt2;
600  d2(4, 4) -= s(1) / 2.0;
601  d2(4, 5) += s(3) / 2.0 / MathUtils::sqrt2;
602 
603  d2(5, 2) -= s(5) / MathUtils::sqrt2;
604  d2(5, 3) += s(4) / 2.0 / MathUtils::sqrt2;
605  d2(5, 4) += s(3) / 2.0 / MathUtils::sqrt2;
606  d2(5, 5) -= s(2) / 2.0;
607 
608  for (auto i : make_range(N))
609  for (auto j : make_range(N))
611 
612  return d2;
613 }
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...
SymmetricRankTwoTensorTempl< T > deviatoric() const
returns A_ij - de_ij*tr(A)/3, where A are the _vals
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static SymmetricRankTwoTensorTempl< T > plusTranspose(const RankTwoTensorTempl< T > &)
return the matrix plus its transpose A-A^T (guaranteed symmetric)
SymmetricRankFourTensorTempl is designed to handle an N-dimensional fourth order tensor with minor sy...
IntRange< T > make_range(T beg, T end)
static constexpr unsigned int N
static constexpr Real sqrt2
std::sqrt is not constexpr, so we add sqrt(2) as a constant (used in Mandel notation) ...
Definition: MathUtils.h:25

◆ ddet()

template<typename T >
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::ddet ( ) const

Denote the _vals[i][j] by A_ij, then this returns d(det)/dA_ij.

Definition at line 628 of file SymmetricRankTwoTensorImplementation.h.

629 {
631  _vals[2] * _vals[1] - _vals[3] * _vals[3] / 2.0,
632  _vals[0] * _vals[2] - _vals[4] * _vals[4] / 2.0,
633  _vals[0] * _vals[1] - _vals[5] * _vals[5] / 2.0,
634  _vals[4] * _vals[5] / 2.0 - _vals[0] * _vals[3] / MathUtils::sqrt2,
635  _vals[5] * _vals[3] / 2.0 - _vals[4] * _vals[1] / MathUtils::sqrt2,
636  _vals[4] * _vals[3] / 2.0 - _vals[5] * _vals[2] / MathUtils::sqrt2);
637 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...
static constexpr Real sqrt2
std::sqrt is not constexpr, so we add sqrt(2) as a constant (used in Mandel notation) ...
Definition: MathUtils.h:25

◆ det()

template<typename T >
T SymmetricRankTwoTensorTempl< T >::det ( ) const

Definition at line 617 of file SymmetricRankTwoTensorImplementation.h.

618 {
619  return _vals[0] * (_vals[2] * _vals[1] - _vals[3] * _vals[3] / 2.0) -
620  _vals[5] / MathUtils::sqrt2 *
621  (_vals[2] * _vals[5] / MathUtils::sqrt2 - _vals[3] * _vals[4] / 2.0) +
622  _vals[4] / MathUtils::sqrt2 *
623  (_vals[3] * _vals[5] / 2.0 - _vals[1] * _vals[4] / MathUtils::sqrt2);
624 }
static constexpr Real sqrt2
std::sqrt is not constexpr, so we add sqrt(2) as a constant (used in Mandel notation) ...
Definition: MathUtils.h:25

◆ deviatoric()

template<typename T >
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::deviatoric ( ) const

returns A_ij - de_ij*tr(A)/3, where A are the _vals

Definition at line 463 of file SymmetricRankTwoTensorImplementation.h.

464 {
466  deviatoric.addIa(-1.0 / 3.0 * this->tr());
467  return deviatoric;
468 }
T tr() const
Returns the trace.
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...
SymmetricRankTwoTensorTempl< T > deviatoric() const
returns A_ij - de_ij*tr(A)/3, where A are the _vals

◆ doubleContraction()

template<typename T>
T SymmetricRankTwoTensorTempl< T >::doubleContraction ( const SymmetricRankTwoTensorTempl< T > &  a) const

returns _vals_ij * a_ij (sum on i, j)

Definition at line 438 of file SymmetricRankTwoTensorImplementation.h.

439 {
440  T sum = 0;
441  for (unsigned int i = 0; i < N; ++i)
442  sum += _vals[i] * b._vals[i];
443  return sum;
444 }
static constexpr unsigned int N

◆ dsecondInvariant()

template<typename T >
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::dsecondInvariant ( ) const

Denote the _vals[i][j] by A_ij, then this returns d(secondInvariant)/dA_ij.

Definition at line 494 of file SymmetricRankTwoTensorImplementation.h.

495 {
497 }
SymmetricRankTwoTensorTempl< T > deviatoric() const
returns A_ij - de_ij*tr(A)/3, where A are the _vals
static SymmetricRankTwoTensorTempl< T > plusTranspose(const RankTwoTensorTempl< T > &)
return the matrix plus its transpose A-A^T (guaranteed symmetric)

◆ dsin3Lode()

template<typename T>
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::dsin3Lode ( const T &  r0) const

d(sin3Lode)/dA_ij If secondInvariant() <= r0 then return zero This is to gaurd against precision-loss errors.

Note that sin(3*Lode_angle) is not defined for secondInvariant() = 0

Definition at line 668 of file SymmetricRankTwoTensor.h.

669 {
671  {
672  T bar = secondInvariant();
673  if (bar <= r0)
675  else
676  return -1.5 * std::sqrt(3.0) *
677  (dthirdInvariant() / std::pow(bar, 1.5) -
678  1.5 * dsecondInvariant() * thirdInvariant() / std::pow(bar, 2.5));
679  }
680  else
681  mooseError("dsin3Lode is only available for ordered tensor component types");
682 }
SymmetricRankTwoTensorTempl< T > dsecondInvariant() const
Denote the _vals[i][j] by A_ij, then this returns d(secondInvariant)/dA_ij.
T thirdInvariant() const
Denote the _vals[i][j] by A_ij, then S_ij = A_ij - de_ij*tr(A)/3 Then this returns det(S + S...
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
Custom type trait that has a value of true for types that cam be use interchangeably with Real...
Definition: MooseUtils.h:985
T secondInvariant() const
Denote the _vals[i][j] by A_ij, then S_ij = A_ij - de_ij*tr(A)/3 Then this returns (S_ij + S_ji)*(S_i...
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...
SymmetricRankTwoTensorTempl< T > dthirdInvariant() const
Denote the _vals[i][j] by A_ij, then this returns d(thirdInvariant()/dA_ij.
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sqrt(_arg)) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tanh
MooseUnits pow(const MooseUnits &, int)
Definition: Units.C:537

◆ dthirdInvariant()

template<typename T >
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::dthirdInvariant ( ) const

Denote the _vals[i][j] by A_ij, then this returns d(thirdInvariant()/dA_ij.

Definition at line 557 of file SymmetricRankTwoTensorImplementation.h.

558 {
560  const T s3 = secondInvariant() / 3.0;
561  return SymmetricRankTwoTensorTempl<T>(s(1) * s(2) - s(3) * s(3) / 2.0 + s3,
562  s(0) * s(2) - s(4) * s(4) / 2.0 + s3,
563  s(0) * s(1) - s(5) * s(5) / 2.0 + s3,
564  s(5) * s(4) / 2.0 - s(0) * s(3) / MathUtils::sqrt2,
565  s(5) * s(3) / 2.0 - s(1) * s(4) / MathUtils::sqrt2,
566  s(4) * s(3) / 2.0 - s(5) * s(2) / MathUtils::sqrt2);
567 }
T secondInvariant() const
Denote the _vals[i][j] by A_ij, then S_ij = A_ij - de_ij*tr(A)/3 Then this returns (S_ij + S_ji)*(S_i...
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...
SymmetricRankTwoTensorTempl< T > deviatoric() const
returns A_ij - de_ij*tr(A)/3, where A are the _vals
static SymmetricRankTwoTensorTempl< T > plusTranspose(const RankTwoTensorTempl< T > &)
return the matrix plus its transpose A-A^T (guaranteed symmetric)
static constexpr Real sqrt2
std::sqrt is not constexpr, so we add sqrt(2) as a constant (used in Mandel notation) ...
Definition: MathUtils.h:25

◆ dtrace()

template<typename T >
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::dtrace ( ) const

Denote the _vals[i][j] by A_ij, then this returns d(trace)/dA_ij.

Definition at line 538 of file SymmetricRankTwoTensorImplementation.h.

539 {
540  return SymmetricRankTwoTensorTempl<T>(1, 1, 1, 0, 0, 0);
541 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...

◆ fillFromInputVector()

template<typename T>
void SymmetricRankTwoTensorTempl< T >::fillFromInputVector ( const std::vector< T > &  input,
FillMethod  fill_method = autodetect 
)

fillFromInputVector takes 1, 3, or 6 inputs to fill in the symmmetric Rank-2 tensor.

Definition at line 173 of file SymmetricRankTwoTensorImplementation.h.

Referenced by GenericConstantSymmetricRankTwoTensorTempl< is_ad >::GenericConstantSymmetricRankTwoTensorTempl(), and SymmetricRankTwoTensorTempl< Real >::SymmetricRankTwoTensorTempl().

175 {
176  if (fill_method != autodetect && fill_method != input.size())
177  mooseError("Expected an input vector size of ",
178  fill_method,
179  " to fill the SymmetricRankTwoTensorTempl");
180 
181  switch (input.size())
182  {
183  case 1:
184  _vals[0] = input[0];
185  _vals[1] = input[0];
186  _vals[2] = input[0];
187  _vals[3] = 0.0;
188  _vals[4] = 0.0;
189  _vals[5] = 0.0;
190  break;
191 
192  case 3:
193  _vals[0] = input[0];
194  _vals[1] = input[1];
195  _vals[2] = input[2];
196  _vals[3] = 0.0;
197  _vals[4] = 0.0;
198  _vals[5] = 0.0;
199  break;
200 
201  case 6:
202  _vals[0] = input[0];
203  _vals[1] = input[1];
204  _vals[2] = input[2];
205  _vals[3] = mandelFactor(3) * input[3];
206  _vals[4] = mandelFactor(4) * input[4];
207  _vals[5] = mandelFactor(5) * input[5];
208  break;
209 
210  default:
211  mooseError("Please check the number of entries in the input vector for building "
212  "a SymmetricRankTwoTensorTempl. It must be 1, 3, 6");
213  }
214 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
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...

◆ fillFromScalarVariable()

template<typename T >
void SymmetricRankTwoTensorTempl< T >::fillFromScalarVariable ( const VariableValue scalar_variable)

fillFromScalarVariable takes FIRST/THIRD/SIXTH order scalar variable to fill in the Rank-2 tensor.

Definition at line 218 of file SymmetricRankTwoTensorImplementation.h.

219 {
220  switch (scalar_variable.size())
221  {
222  case 1:
223  _vals[0] = scalar_variable[0];
224  _vals[1] = 0.0;
225  _vals[2] = 0.0;
226  _vals[3] = 0.0;
227  _vals[4] = 0.0;
228  _vals[5] = 0.0;
229  break;
230 
231  case 3:
232  _vals[0] = scalar_variable[0];
233  _vals[1] = scalar_variable[1];
234  _vals[2] = 0.0;
235  _vals[3] = 0.0;
236  _vals[4] = 0.0;
237  _vals[5] = mandelFactor(5) * scalar_variable[2];
238  break;
239 
240  case 6:
241  _vals[0] = scalar_variable[0];
242  _vals[1] = scalar_variable[1];
243  _vals[2] = scalar_variable[2];
244  _vals[3] = mandelFactor(3) * scalar_variable[3];
245  _vals[4] = mandelFactor(4) * scalar_variable[4];
246  _vals[5] = mandelFactor(5) * scalar_variable[5];
247  break;
248 
249  default:
250  mooseError("Only FIRST, THIRD, or SIXTH order scalar variable can be used to build "
251  "a SymmetricRankTwoTensorTempl.");
252  }
253 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
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...

◆ fillMethodEnum()

template<typename T >
MooseEnum SymmetricRankTwoTensorTempl< T >::fillMethodEnum ( )
static

Static method for use in validParams for getting the "fill_method".

Definition at line 50 of file SymmetricRankTwoTensorImplementation.h.

51 {
52  return MooseEnum("autodetect=0 isotropic1=1 diagonal3=3 symmetric6=6", "autodetect");
53 }
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33

◆ fromRawComponents()

template<typename T>
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::fromRawComponents ( const T &  S11,
const T &  S22,
const T &  S33,
const T &  S23,
const T &  S13,
const T &  S12 
)
staticprivate

Definition at line 125 of file SymmetricRankTwoTensorImplementation.h.

Referenced by SymmetricRankTwoTensorTempl< Real >::timesTranspose().

127 {
129  ret._vals[0] = S11;
130  ret._vals[1] = S22;
131  ret._vals[2] = S33;
132  ret._vals[3] = S23;
133  ret._vals[4] = S13;
134  ret._vals[5] = S12;
135  return ret;
136 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...

◆ generalSecondInvariant()

template<typename T >
T SymmetricRankTwoTensorTempl< T >::generalSecondInvariant ( ) const

Calculates the second invariant (I2) of a tensor.

Definition at line 472 of file SymmetricRankTwoTensorImplementation.h.

473 {
474  return _vals[0] * _vals[1] + _vals[0] * _vals[2] + _vals[1] * _vals[2] -
475  (_vals[3] * _vals[3] + _vals[4] * _vals[4] + _vals[5] * _vals[5]) / 2.0;
476 }

◆ genRandomSymmTensor()

template<typename T>
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::genRandomSymmTensor ( scale,
offset 
)
static

This function generates a random symmetric rank two tensor.

The first real scales the random number. The second real offsets the uniform random number

Definition at line 718 of file SymmetricRankTwoTensorImplementation.h.

719 {
720  auto r = [&]() { return (MooseRandom::rand() + offset) * scale; };
721  return SymmetricRankTwoTensorTempl<T>(r(), r(), r(), r(), r(), r());
722 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...
static Real rand()
This method returns the next random number (Real format) from the generator.
Definition: MooseRandom.h:50

◆ identity()

template<typename T>
static SymmetricRankTwoTensorTempl SymmetricRankTwoTensorTempl< T >::identity ( )
inlinestatic

Definition at line 163 of file SymmetricRankTwoTensor.h.

◆ initialContraction()

template<typename T>
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::initialContraction ( const SymmetricRankFourTensorTempl< T > &  b) const

returns this_ij * b_ijkl

Definition at line 734 of file SymmetricRankTwoTensorImplementation.h.

735 {
737  for (auto i : make_range(N))
738  for (auto j : make_range(N))
739  result(j) += (*this)(i)*b(i, j);
740  return result;
741 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...
IntRange< T > make_range(T beg, T end)
static constexpr unsigned int N

◆ initializeSymmetric()

template<typename T>
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::initializeSymmetric ( const TypeVector< T > &  v0,
const TypeVector< T > &  v1,
const TypeVector< T > &  v2 
)
static

named constructor for initializing symmetrically

Definition at line 163 of file SymmetricRankTwoTensorImplementation.h.

166 {
168  v0(0), v1(0), v2(0), v0(1), v1(1), v2(1), v0(2), v1(2), v2(2));
169 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...

◆ initRandom()

template<typename T >
void SymmetricRankTwoTensorTempl< T >::initRandom ( unsigned int  rand_seed)
static

This function initializes random seed based on a user-defined number.

Definition at line 711 of file SymmetricRankTwoTensorImplementation.h.

712 {
713  MooseRandom::seed(rand_seed);
714 }
static void seed(unsigned int seed)
The method seeds the random number generator.
Definition: MooseRandom.h:44

◆ inverse()

template<typename T >
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::inverse ( ) const

retuns the inverse of the tensor

Definition at line 521 of file SymmetricRankTwoTensorImplementation.h.

522 {
523  const auto d = det();
524  if (d == 0.0)
525  mooseException("Matrix not invertible");
527  _vals[2] * _vals[1] - _vals[3] * _vals[3] / 2.0,
528  _vals[2] * _vals[0] - _vals[4] * _vals[4] / 2.0,
529  _vals[0] * _vals[1] - _vals[5] * _vals[5] / 2.0,
530  _vals[5] * _vals[4] / 2.0 - _vals[0] * _vals[3] / MathUtils::sqrt2,
531  _vals[5] * _vals[3] / 2.0 - _vals[4] * _vals[1] / MathUtils::sqrt2,
532  _vals[4] * _vals[3] / 2.0 - _vals[2] * _vals[5] / MathUtils::sqrt2);
533  return inv * 1.0 / d;
534 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...
static constexpr Real sqrt2
std::sqrt is not constexpr, so we add sqrt(2) as a constant (used in Mandel notation) ...
Definition: MathUtils.h:25

◆ isSymmetric()

template<typename T>
bool SymmetricRankTwoTensorTempl< T >::isSymmetric ( ) const
inline

Test for symmetry. Surprisingly this is always true.

Definition at line 309 of file SymmetricRankTwoTensor.h.

309 { return true; }

◆ L2norm()

template<typename T >
T SymmetricRankTwoTensorTempl< T >::L2norm ( ) const

Sqrt(_vals[i][j]*_vals[i][j])

Definition at line 665 of file SymmetricRankTwoTensorImplementation.h.

666 {
667  T norm = _vals[0] * _vals[0];
668  for (unsigned int i = 1; i < N; ++i)
669  norm += _vals[i] * _vals[i];
670  return norm == 0.0 ? 0.0 : std::sqrt(norm);
671 }
auto norm(const T &a) -> decltype(std::abs(a))
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sqrt(_arg)) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tanh
static constexpr unsigned int N

◆ mandelFactor()

template<typename T>
static constexpr Real SymmetricRankTwoTensorTempl< T >::mandelFactor ( unsigned int  i)
inlinestatic

returns the 1 or sqrt(2) prefactor in the Mandel notation for the index i ranging from 0-5.

Definition at line 93 of file SymmetricRankTwoTensor.h.

Referenced by SymmetricRankFourTensorTempl< T >::mandelFactor(), and SymmetricRankTwoTensorTempl< Real >::operator()().

93 { return i < Ndim ? 1.0 : MathUtils::sqrt2; }
static constexpr unsigned int Ndim
tensor dimension and Mandel vector length
static constexpr Real sqrt2
std::sqrt is not constexpr, so we add sqrt(2) as a constant (used in Mandel notation) ...
Definition: MathUtils.h:25

◆ operator RankTwoTensorTempl< T >()

template<typename T >
SymmetricRankTwoTensorTempl< T >::operator RankTwoTensorTempl< T > ( )
explicit

Definition at line 91 of file SymmetricRankTwoTensorImplementation.h.

92 {
93  return RankTwoTensorTempl<T>(_vals[0],
94  _vals[5] / mandelFactor(5),
95  _vals[4] / mandelFactor(4),
96  _vals[5] / mandelFactor(5),
97  _vals[1],
98  _vals[3] / mandelFactor(3),
99  _vals[4] / mandelFactor(4),
100  _vals[3] / mandelFactor(3),
101  _vals[2]);
102 }
RankTwoTensorTempl is designed to handle the Stress or Strain Tensor for a fully anisotropic material...
Definition: RankTwoTensor.h:87
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...

◆ operator!=()

template<typename T >
template<typename T2 >
bool SymmetricRankTwoTensorTempl< T >::operator!= ( const SymmetricRankTwoTensorTempl< T2 > &  a) const

Defines logical inequality with another SymmetricRankTwoTensorTempl<T2>

Definition at line 573 of file SymmetricRankTwoTensor.h.

574 {
575  for (std::size_t i = 0; i < N; ++i)
576  if (_vals[i] != a._vals[i])
577  return true;
578  return false;
579 }
static constexpr unsigned int N

◆ operator()() [1/3]

template<typename T>
T& SymmetricRankTwoTensorTempl< T >::operator() ( unsigned int  i)
inline

Gets the raw value for the index specified. Takes index = 0,1,2,3,4,5.

Definition at line 185 of file SymmetricRankTwoTensor.h.

185 { return _vals[i]; }

◆ operator()() [2/3]

template<typename T>
T SymmetricRankTwoTensorTempl< T >::operator() ( unsigned int  i) const
inline

Gets the raw value for the index specified.

Takes index = 0,1,2 used for const

Definition at line 191 of file SymmetricRankTwoTensor.h.

191 { return _vals[i]; }

◆ operator()() [3/3]

template<typename T>
T SymmetricRankTwoTensorTempl< T >::operator() ( unsigned int  i,
unsigned int  j 
) const
inline

Gets the value for the index specified.

Takes index = 0,1,2

Definition at line 196 of file SymmetricRankTwoTensor.h.

197  {
198  const auto a = reverse_index[i][j];
199  return _vals[a] / mandelFactor(a);
200  }
static constexpr unsigned int reverse_index[3][3]
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...

◆ operator*() [1/2]

template<typename T>
template<typename T2 >
auto SymmetricRankTwoTensorTempl< T >::operator* ( const T2 &  a) const -> typename std::enable_if<libMesh::ScalarTraits<T2>::value, SymmetricRankTwoTensorTempl<decltype(T() * T2())>>::type

returns _vals*a

Definition at line 524 of file SymmetricRankTwoTensor.h.

527 {
529  for (const auto i : make_range(N))
530  result._vals[i] = _vals[i] * a;
531  return result;
532 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...
IntRange< T > make_range(T beg, T end)
static constexpr unsigned int N

◆ operator*() [2/2]

template<typename T>
template<typename T2 >
TypeVector< typename libMesh::CompareTypes< T, T2 >::supertype > SymmetricRankTwoTensorTempl< T >::operator* ( const TypeVector< T2 > &  a) const

Defines multiplication with a vector to get a vector.

Definition at line 551 of file SymmetricRankTwoTensor.h.

552 {
553  TypeVector<typename libMesh::CompareTypes<T, T2>::supertype> ret;
554  ret(0) = a(0) * _vals[0] + a(1) * _vals[5] + a(2) * _vals[4];
555  ret(1) = a(0) * _vals[5] + a(1) * _vals[1] + a(2) * _vals[3];
556  ret(2) = a(0) * _vals[4] + a(1) * _vals[3] + a(2) * _vals[2];
557 }

◆ operator*=()

template<typename T>
SymmetricRankTwoTensorTempl< T > & SymmetricRankTwoTensorTempl< T >::operator*= ( const T &  a)

performs _vals *= a

Definition at line 420 of file SymmetricRankTwoTensorImplementation.h.

421 {
422  for (std::size_t i = 0; i < N; ++i)
423  _vals[i] *= a;
424  return *this;
425 }
static constexpr unsigned int N

◆ operator+()

template<typename T >
template<typename T2 >
SymmetricRankTwoTensorTempl< typename libMesh::CompareTypes< T, T2 >::supertype > SymmetricRankTwoTensorTempl< T >::operator+ ( const SymmetricRankTwoTensorTempl< T2 > &  a) const

returns _vals + a

Definition at line 390 of file SymmetricRankTwoTensorImplementation.h.

391 {
393  for (std::size_t i = 0; i < N; ++i)
394  result(i) = _vals[i] + a(i);
395  return result;
396 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...
static constexpr unsigned int N

◆ operator+=()

template<typename T>
SymmetricRankTwoTensorTempl< T > & SymmetricRankTwoTensorTempl< T >::operator+= ( const SymmetricRankTwoTensorTempl< T > &  a)

adds a to _vals

Definition at line 371 of file SymmetricRankTwoTensorImplementation.h.

372 {
373  for (std::size_t i = 0; i < N; ++i)
374  _vals[i] += a._vals[i];
375  return *this;
376 }
static constexpr unsigned int N

◆ operator-() [1/2]

template<typename T >
template<typename T2 >
SymmetricRankTwoTensorTempl< typename libMesh::CompareTypes< T, T2 >::supertype > SymmetricRankTwoTensorTempl< T >::operator- ( const SymmetricRankTwoTensorTempl< T2 > &  a) const

returns _vals - a

Definition at line 401 of file SymmetricRankTwoTensorImplementation.h.

402 {
405  for (std::size_t i = 0; i < N; ++i)
406  result(i) = _vals[i] - a(i);
407  return result;
408 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...
static constexpr unsigned int N

◆ operator-() [2/2]

template<typename T >
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::operator- ( ) const

returns -_vals

Definition at line 413 of file SymmetricRankTwoTensorImplementation.h.

414 {
415  return (*this) * -1.0;
416 }

◆ operator-=()

template<typename T>
SymmetricRankTwoTensorTempl< T > & SymmetricRankTwoTensorTempl< T >::operator-= ( const SymmetricRankTwoTensorTempl< T > &  a)

sets _vals -= a and returns vals

Definition at line 380 of file SymmetricRankTwoTensorImplementation.h.

381 {
382  for (std::size_t i = 0; i < N; ++i)
383  _vals[i] -= a._vals[i];
384  return *this;
385 }
static constexpr unsigned int N

◆ operator/()

template<typename T>
template<typename T2 >
auto SymmetricRankTwoTensorTempl< T >::operator/ ( const T2 &  a) const -> typename std::enable_if<libMesh::ScalarTraits<T2>::value, SymmetricRankTwoTensorTempl<decltype(T() / T2())>>::type

returns _vals/a

Definition at line 537 of file SymmetricRankTwoTensor.h.

540 {
541  SymmetricRankTwoTensorTempl<decltype(T() / T2())> result;
542  for (const auto i : make_range(N))
543  result._vals[i] = _vals[i] / a;
544  return result;
545 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...
IntRange< T > make_range(T beg, T end)
static constexpr unsigned int N

◆ operator/=()

template<typename T>
SymmetricRankTwoTensorTempl< T > & SymmetricRankTwoTensorTempl< T >::operator/= ( const T &  a)

performs _vals /= a

Definition at line 429 of file SymmetricRankTwoTensorImplementation.h.

430 {
431  for (std::size_t i = 0; i < N; ++i)
432  _vals[i] /= a;
433  return *this;
434 }
static constexpr unsigned int N

◆ operator=() [1/3]

template<typename T >
template<typename T2 >
SymmetricRankTwoTensorTempl< T > & SymmetricRankTwoTensorTempl< T >::operator= ( const SymmetricRankTwoTensorTempl< T2 > &  a)

sets _vals to a, and returns _vals

Definition at line 687 of file SymmetricRankTwoTensor.h.

688 {
689  for (const auto i : make_range(N))
690  (*this)(i) = a(i);
691  return *this;
692 }
IntRange< T > make_range(T beg, T end)
static constexpr unsigned int N

◆ operator=() [2/3]

template<typename T>
template<typename Scalar >
boostcopy::enable_if_c<libMesh::ScalarTraits<Scalar>::value, SymmetricRankTwoTensorTempl &>::type SymmetricRankTwoTensorTempl< T >::operator= ( const Scalar &  libmesh_dbg_varp)
inline

Assignment-from-scalar operator.

Used only to zero out vectors.

Definition at line 255 of file SymmetricRankTwoTensor.h.

256  {
257  libmesh_assert_equal_to(p, Scalar(0));
258  this->zero();
259  return *this;
260  }
void zero()
Set all components to zero.

◆ operator=() [3/3]

template<typename T>
SymmetricRankTwoTensorTempl<T>& SymmetricRankTwoTensorTempl< T >::operator= ( const ColumnMajorMatrixTempl< T > &  a)

Sets _vals to the values in a ColumnMajorMatrix (must be 3x3)

◆ operator==()

template<typename T >
template<typename T2 >
bool SymmetricRankTwoTensorTempl< T >::operator== ( const SymmetricRankTwoTensorTempl< T2 > &  a) const

Defines logical equality with another SymmetricRankTwoTensorTempl<T2>

Definition at line 562 of file SymmetricRankTwoTensor.h.

563 {
564  for (std::size_t i = 0; i < N; ++i)
565  if (_vals[i] != a._vals[i])
566  return false;
567  return true;
568 }
static constexpr unsigned int N

◆ outerProduct()

template<typename T>
SymmetricRankFourTensorTempl< T > SymmetricRankTwoTensorTempl< T >::outerProduct ( const SymmetricRankTwoTensorTempl< T > &  a) const

returns C_ijkl = a_ij * b_kl

Definition at line 448 of file SymmetricRankTwoTensorImplementation.h.

449 {
451  unsigned int index = 0;
452  for (unsigned int i = 0; i < N; ++i)
453  {
454  const T & a = _vals[i];
455  for (unsigned int j = 0; j < N; ++j)
456  result._vals[index++] = a * b._vals[j];
457  }
458  return result;
459 }
std::array< T, N2 > _vals
The values of the rank-four tensor.
SymmetricRankFourTensorTempl is designed to handle an N-dimensional fourth order tensor with minor sy...
static constexpr unsigned int N

◆ plusTranspose() [1/2]

template<typename T>
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::plusTranspose ( const RankTwoTensorTempl< T > &  a)
static

return the matrix plus its transpose A-A^T (guaranteed symmetric)

Definition at line 340 of file SymmetricRankTwoTensorImplementation.h.

Referenced by SymmetricRankTwoTensorTempl< Real >::d2thirdInvariant(), SymmetricRankTwoTensorTempl< Real >::dsecondInvariant(), SymmetricRankTwoTensorTempl< Real >::dthirdInvariant(), and SymmetricRankTwoTensorTempl< Real >::thirdInvariant().

341 {
343  a(0, 0), a(0, 1), a(0, 2), a(1, 0), a(1, 1), a(1, 2), a(2, 0), a(2, 1), a(2, 2)) *
344  2.0;
345 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...

◆ plusTranspose() [2/2]

template<typename T>
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::plusTranspose ( const SymmetricRankTwoTensorTempl< T > &  a)
static

Definition at line 349 of file SymmetricRankTwoTensorImplementation.h.

350 {
351  return a * 2.0;
352 }

◆ positiveProjectionEigenDecomposition()

template<typename T>
SymmetricRankFourTensorTempl< T > SymmetricRankTwoTensorTempl< T >::positiveProjectionEigenDecomposition ( std::vector< T > &  eigval,
RankTwoTensorTempl< T > &  eigvec 
) const

return positive projection tensor of eigen-decomposition

Definition at line 583 of file SymmetricRankTwoTensor.h.

585 {
586  // The calculate of projection tensor follows
587  // C. Miehe and M. Lambrecht, Commun. Numer. Meth. Engng 2001; 17:337~353
588 
589  // Compute eigenvectors and eigenvalues of this tensor
591  {
592  this->symmetricEigenvaluesEigenvectors(eigval, eigvec);
593 
594  // Separate out positive and negative eigen values
595  std::array<T, N> epos;
596  std::array<T, N> d;
597  for (unsigned int i = 0; i < Ndim; ++i)
598  {
599  epos[i] = (std::abs(eigval[i]) + eigval[i]) / 2.0;
600  d[i] = 0 < eigval[i] ? 1.0 : 0.0;
601  }
602 
603  // projection tensor
605 
606  for (unsigned int a = 0; a < Ndim; ++a)
607  {
608  const auto Ma = SymmetricRankTwoTensorTempl<T>::selfOuterProduct(eigvec.column(a));
609  proj_pos += d[a] * Ma.outerProduct(Ma);
610  }
611 
612  for (const auto a : make_range(Ndim))
613  for (const auto b : make_range(a))
614  {
615  const auto Ma = SymmetricRankTwoTensorTempl<T>::selfOuterProduct(eigvec.column(a));
616  const auto Mb = SymmetricRankTwoTensorTempl<T>::selfOuterProduct(eigvec.column(b));
617 
619  for (const auto aa : make_range(N))
620  for (const auto bb : make_range(N))
621  {
622  const auto i = SymmetricRankFourTensorTempl<T>::full_index[aa][bb][0];
623  const auto j = SymmetricRankFourTensorTempl<T>::full_index[aa][bb][1];
624  const auto k = SymmetricRankFourTensorTempl<T>::full_index[aa][bb][2];
625  const auto l = SymmetricRankFourTensorTempl<T>::full_index[aa][bb][3];
626 
627  Gabba(aa, bb) = (Ma(i, k) * Mb(j, l) + Ma(i, l) * Mb(j, k) + Ma(j, l) * Mb(i, k) +
628  Ma(j, k) * Mb(i, l)) *
630  }
631 
632  T theta_ab;
633  if (!MooseUtils::relativeFuzzyEqual(eigval[a], eigval[b]))
634  theta_ab = 0.5 * (epos[a] - epos[b]) / (eigval[a] - eigval[b]);
635  else
636  theta_ab = 0.25 * (d[a] + d[b]);
637 
638  proj_pos += theta_ab * Gabba;
639  }
640  return proj_pos;
641  }
642  else
643  mooseError("positiveProjectionEigenDecomposition is only available for ordered tensor "
644  "component types");
645 }
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
Definition: EigenADReal.h:42
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
void symmetricEigenvaluesEigenvectors(std::vector< T > &eigvals, RankTwoTensorTempl< T > &eigvecs) const
computes eigenvalues and eigenvectors, assuming tens is symmetric, and places them in ascending order...
static constexpr unsigned int Ndim
tensor dimension and Mandel vector length
Custom type trait that has a value of true for types that cam be use interchangeably with Real...
Definition: MooseUtils.h:985
libMesh::VectorValue< T > column(const unsigned int i) const
Get the i-th column of the second order tensor.
static SymmetricRankTwoTensorTempl< T > selfOuterProduct(const TypeVector< T > &)
SymmetricRankTwoTensorTempl<T> from outer product of a vector with itself.
bool relativeFuzzyEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
Function to check whether two variables are equal within a relative tolerance.
Definition: MooseUtils.h:494
SymmetricRankFourTensorTempl is designed to handle an N-dimensional fourth order tensor with minor sy...
IntRange< T > make_range(T beg, T end)
static constexpr unsigned int N

◆ print()

template<typename T >
void SymmetricRankTwoTensorTempl< T >::print ( std::ostream &  stm = Moose::out) const

Print the rank two tensor.

Definition at line 641 of file SymmetricRankTwoTensorImplementation.h.

642 {
643  for (std::size_t i = 0; i < N; ++i)
644  stm << std::setw(15) << _vals[i] << std::endl;
645 }
static constexpr unsigned int N

◆ printADReal()

template<typename T>
void SymmetricRankTwoTensorTempl< T >::printADReal ( unsigned int  nDual,
std::ostream &  stm = Moose::out 
) const

Print the Real part of the ADReal rank two tensor along with its first nDual dual numbers.

◆ printReal()

template<typename T >
void SymmetricRankTwoTensorTempl< T >::printReal ( std::ostream &  stm = Moose::out) const

Print the Real part of the ADReal rank two tensor.

Definition at line 649 of file SymmetricRankTwoTensorImplementation.h.

650 {
651  for (std::size_t i = 0; i < N; ++i)
652  stm << std::setw(15) << MetaPhysicL::raw_value(_vals[i]) << std::endl;
653 }
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
static constexpr unsigned int N

◆ rotate()

template<typename T>
void SymmetricRankTwoTensorTempl< T >::rotate ( const TypeTensor< T > &  R)

rotates the tensor data given a rank two tensor rotation tensor _vals[i][j] = R_ij * R_jl * _vals[k][l]

Parameters
Rrotation matrix as a TypeTensor

Definition at line 257 of file SymmetricRankTwoTensorImplementation.h.

258 {
260  *this = M * (*this);
261 }
static SymmetricRankFourTensorTempl< T > rotationMatrix(const TypeTensor< T > &R)
Build a 6x6 rotation matrix MEHRABADI, MORTEZA M.

◆ row()

template<typename T >
libMesh::VectorValue< T > SymmetricRankTwoTensorTempl< T >::row ( const unsigned int  n) const

get the specified row of the tensor

multiply vector v with row n of this tensor

Definition at line 273 of file SymmetricRankTwoTensorImplementation.h.

Referenced by SymmetricRankTwoTensorTempl< Real >::column().

274 {
275  switch (n)
276  {
277  case 0:
279  _vals[0], _vals[5] / mandelFactor(5), _vals[4] / mandelFactor(4));
280  case 1:
282  _vals[5] / mandelFactor(5), _vals[1], _vals[3] / mandelFactor(3));
283  case 2:
285  _vals[4] / mandelFactor(4), _vals[3] / mandelFactor(3), _vals[2]);
286  default:
287  mooseError("Invalid row");
288  }
289 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
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...

◆ secondInvariant()

template<typename T >
T SymmetricRankTwoTensorTempl< T >::secondInvariant ( ) const

Denote the _vals[i][j] by A_ij, then S_ij = A_ij - de_ij*tr(A)/3 Then this returns (S_ij + S_ji)*(S_ij + S_ji)/8 Note the explicit symmeterisation.

Definition at line 480 of file SymmetricRankTwoTensorImplementation.h.

481 {
482  T result;
483  result = Utility::pow<2>(_vals[0] - _vals[1]) / 6.0;
484  result += Utility::pow<2>(_vals[0] - _vals[2]) / 6.0;
485  result += Utility::pow<2>(_vals[1] - _vals[2]) / 6.0;
486  result += Utility::pow<2>(_vals[5]) / 2.0;
487  result += Utility::pow<2>(_vals[4]) / 2.0;
488  result += Utility::pow<2>(_vals[3]) / 2.0;
489  return result;
490 }
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template pow< 2 >(tan(_arg))+1.0) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(sqrt

◆ selfOuterProduct()

template<typename T>
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::selfOuterProduct ( const TypeVector< T > &  v)
static

SymmetricRankTwoTensorTempl<T> from outer product of a vector with itself.

Definition at line 726 of file SymmetricRankTwoTensorImplementation.h.

Referenced by SymmetricRankTwoTensorTempl< Real >::positiveProjectionEigenDecomposition().

727 {
729  v(0) * v(0), v(1) * v(1), v(2) * v(2), v(1) * v(2), v(0) * v(2), v(0) * v(1));
730 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...

◆ setToIdentity()

template<typename T >
void SymmetricRankTwoTensorTempl< T >::setToIdentity ( )

set the tensor to the identity matrix

Definition at line 745 of file SymmetricRankTwoTensorImplementation.h.

746 {
747  std::copy(identityCoords.begin(), identityCoords.end(), _vals.begin());
748 }
static constexpr std::array< Real, N > identityCoords

◆ sin3Lode()

template<typename T>
T SymmetricRankTwoTensorTempl< T >::sin3Lode ( const T &  r0,
const T &  r0_value 
) const

Sin(3*Lode_angle) If secondInvariant() <= r0 then return r0_value This is to gaurd against precision-loss errors.

Note that sin(3*Lode_angle) is not defined for secondInvariant() = 0

Definition at line 649 of file SymmetricRankTwoTensor.h.

650 {
652  {
653  T bar = secondInvariant();
654  if (bar <= r0)
655  // in this case the Lode angle is not defined
656  return r0_value;
657  else
658  // the min and max here gaurd against precision-loss when bar is tiny but nonzero.
659  return std::max(std::min(thirdInvariant() * -1.5 * std::sqrt(3.0) / std::pow(bar, 1.5), 1.0),
660  -1.0);
661  }
662  else
663  mooseError("sin3Lode is only available for ordered tensor component types");
664 }
T thirdInvariant() const
Denote the _vals[i][j] by A_ij, then S_ij = A_ij - de_ij*tr(A)/3 Then this returns det(S + S...
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
auto max(const L &left, const R &right)
Custom type trait that has a value of true for types that cam be use interchangeably with Real...
Definition: MooseUtils.h:985
T secondInvariant() const
Denote the _vals[i][j] by A_ij, then S_ij = A_ij - de_ij*tr(A)/3 Then this returns (S_ij + S_ji)*(S_i...
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sqrt(_arg)) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tanh
auto min(const L &left, const R &right)
MooseUnits pow(const MooseUnits &, int)
Definition: Units.C:537

◆ square()

template<typename T >
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::square ( ) const

Returns the matrix squared.

Definition at line 356 of file SymmetricRankTwoTensorImplementation.h.

357 {
359 }
static SymmetricRankTwoTensorTempl< T > timesTranspose(const RankTwoTensorTempl< T > &)
return the matrix multiplied with its transpose A*A^T (guaranteed symmetric)

◆ surfaceFillFromInputVector()

template<typename T>
void SymmetricRankTwoTensorTempl< T >::surfaceFillFromInputVector ( const std::vector< T > &  input)

sets _vals[0][0], _vals[0][1], _vals[1][0], _vals[1][1] to input, and the remainder to zero

◆ syev()

template<typename T>
void SymmetricRankTwoTensorTempl< T >::syev ( const char *  calculation_type,
std::vector< T > &  eigvals,
std::vector< T > &  a 
) const
protected

Uses the petscblaslapack.h LAPACKsyev_ routine to find, for symmetric _vals: (1) the eigenvalues (if calculation_type == "N") (2) the eigenvalues and eigenvectors (if calculation_type == "V")

Parameters
calculation_typeIf "N" then calculation eigenvalues only
eigvalsEigenvalues are placed in this array, in ascending order
aEigenvectors are placed in this array if calculation_type == "V". See code in dsymmetricEigenvalues for extracting eigenvectors from the a output.

Definition at line 675 of file SymmetricRankTwoTensorImplementation.h.

676 {
677  mooseError("The syev method is only supported for Real valued tensors");
678 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323

◆ symmetricEigenvalues()

template<typename T>
void SymmetricRankTwoTensorTempl< T >::symmetricEigenvalues ( std::vector< T > &  eigvals) const

computes eigenvalues, assuming tens is symmetric, and places them in ascending order in eigvals

Definition at line 687 of file SymmetricRankTwoTensorImplementation.h.

688 {
691 }
void symmetricEigenvaluesEigenvectors(std::vector< T > &eigvals, RankTwoTensorTempl< T > &eigvecs) const
computes eigenvalues and eigenvectors, assuming tens is symmetric, and places them in ascending order...
RankTwoTensorTempl is designed to handle the Stress or Strain Tensor for a fully anisotropic material...
Definition: RankTwoTensor.h:87

◆ symmetricEigenvaluesEigenvectors()

template<typename T>
void SymmetricRankTwoTensorTempl< T >::symmetricEigenvaluesEigenvectors ( std::vector< T > &  eigvals,
RankTwoTensorTempl< T > &  eigvecs 
) const

computes eigenvalues and eigenvectors, assuming tens is symmetric, and places them in ascending order in eigvals.

eigvecs is a matrix with the first column being the first eigenvector, the second column being the second, etc.

Definition at line 698 of file SymmetricRankTwoTensorImplementation.h.

700 {
701  mooseError(
702  "symmetricEigenvaluesEigenvectors is only available for ordered tensor component types");
703 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323

◆ thirdInvariant()

template<typename T >
T SymmetricRankTwoTensorTempl< T >::thirdInvariant ( ) const

Denote the _vals[i][j] by A_ij, then S_ij = A_ij - de_ij*tr(A)/3 Then this returns det(S + S.transpose())/2 Note the explicit symmeterisation.

Definition at line 545 of file SymmetricRankTwoTensorImplementation.h.

546 {
548  return s(0) * (s(1) * s(2) - s(3) / MathUtils::sqrt2 * s(3) / MathUtils::sqrt2) -
549  s(5) / MathUtils::sqrt2 *
550  (s(5) / MathUtils::sqrt2 * s(2) - s(3) / MathUtils::sqrt2 * s(4) / MathUtils::sqrt2) +
551  s(4) / MathUtils::sqrt2 *
552  (s(5) / MathUtils::sqrt2 * s(3) / MathUtils::sqrt2 - s(1) * s(4) / MathUtils::sqrt2);
553 }
SymmetricRankTwoTensorTempl< T > deviatoric() const
returns A_ij - de_ij*tr(A)/3, where A are the _vals
static SymmetricRankTwoTensorTempl< T > plusTranspose(const RankTwoTensorTempl< T > &)
return the matrix plus its transpose A-A^T (guaranteed symmetric)
static constexpr Real sqrt2
std::sqrt is not constexpr, so we add sqrt(2) as a constant (used in Mandel notation) ...
Definition: MathUtils.h:25

◆ timesTranspose() [1/2]

template<typename T>
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::timesTranspose ( const RankTwoTensorTempl< T > &  a)
static

return the matrix multiplied with its transpose A*A^T (guaranteed symmetric)

Definition at line 293 of file SymmetricRankTwoTensorImplementation.h.

Referenced by SymmetricRankTwoTensorTempl< Real >::square().

294 {
295  return SymmetricRankTwoTensorTempl<T>(a(0, 0) * a(0, 0) + a(0, 1) * a(0, 1) + a(0, 2) * a(0, 2),
296  a(1, 0) * a(1, 0) + a(1, 1) * a(1, 1) + a(1, 2) * a(1, 2),
297  a(2, 0) * a(2, 0) + a(2, 1) * a(2, 1) + a(2, 2) * a(2, 2),
298  a(1, 0) * a(2, 0) + a(1, 1) * a(2, 1) + a(1, 2) * a(2, 2),
299  a(0, 0) * a(2, 0) + a(0, 1) * a(2, 1) + a(0, 2) * a(2, 2),
300  a(0, 0) * a(1, 0) + a(0, 1) * a(1, 1) + a(0, 2) * a(1, 2));
301 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...

◆ timesTranspose() [2/2]

template<typename T>
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::timesTranspose ( const SymmetricRankTwoTensorTempl< T > &  a)
static

Definition at line 305 of file SymmetricRankTwoTensorImplementation.h.

306 {
308  a._vals[0] * a._vals[0] + a._vals[4] * a._vals[4] / 2.0 + a._vals[5] * a._vals[5] / 2.0,
309  a._vals[1] * a._vals[1] + a._vals[3] * a._vals[3] / 2.0 + a._vals[5] * a._vals[5] / 2.0,
310  a._vals[2] * a._vals[2] + a._vals[3] * a._vals[3] / 2.0 + a._vals[4] * a._vals[4] / 2.0,
311  a._vals[1] * a._vals[3] + a._vals[2] * a._vals[3] +
312  a._vals[4] * a._vals[5] / MathUtils::sqrt2,
313  a._vals[0] * a._vals[4] + a._vals[2] * a._vals[4] +
314  a._vals[3] * a._vals[5] / MathUtils::sqrt2,
315  a._vals[0] * a._vals[5] + a._vals[1] * a._vals[5] +
316  a._vals[3] * a._vals[4] / MathUtils::sqrt2);
317 }
static SymmetricRankTwoTensorTempl fromRawComponents(const T &S11, const T &S22, const T &S33, const T &S23, const T &S13, const T &S12)
static constexpr Real sqrt2
std::sqrt is not constexpr, so we add sqrt(2) as a constant (used in Mandel notation) ...
Definition: MathUtils.h:25

◆ tr()

template<typename T>
T SymmetricRankTwoTensorTempl< T >::tr ( ) const
inline

Returns the trace.

Definition at line 227 of file SymmetricRankTwoTensor.h.

227 { return _vals[0] + _vals[1] + _vals[2]; }

◆ trace()

template<typename T >
T SymmetricRankTwoTensorTempl< T >::trace ( ) const

returns the trace of the tensor, ie _vals[i][i] (sum i = 0, 1, 2)

Definition at line 514 of file SymmetricRankTwoTensorImplementation.h.

515 {
516  return tr();
517 }
T tr() const
Returns the trace.

◆ transpose()

template<typename T >
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::transpose ( ) const

Returns a matrix that is the transpose of the matrix this was called on.

This is a non-operation.

Definition at line 265 of file SymmetricRankTwoTensorImplementation.h.

266 {
267  return *this;
268 }

◆ transposeTimes() [1/2]

template<typename T>
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::transposeTimes ( const RankTwoTensorTempl< T > &  a)
static

return the matrix multiplied with its transpose A^T*A (guaranteed symmetric)

Definition at line 321 of file SymmetricRankTwoTensorImplementation.h.

322 {
323  return SymmetricRankTwoTensorTempl<T>(a(0, 0) * a(0, 0) + a(1, 0) * a(1, 0) + a(2, 0) * a(2, 0),
324  a(0, 1) * a(0, 1) + a(1, 1) * a(1, 1) + a(2, 1) * a(2, 1),
325  a(0, 2) * a(0, 2) + a(1, 2) * a(1, 2) + a(2, 2) * a(2, 2),
326  a(0, 1) * a(0, 2) + a(1, 1) * a(1, 2) + a(2, 1) * a(2, 2),
327  a(0, 0) * a(0, 2) + a(1, 0) * a(1, 2) + a(2, 0) * a(2, 2),
328  a(0, 0) * a(0, 1) + a(1, 0) * a(1, 1) + a(2, 0) * a(2, 1));
329 }
SymmetricRankTwoTensorTempl is designed to handle the Stress or Strain Tensor for an anisotropic mate...

◆ transposeTimes() [2/2]

template<typename T>
SymmetricRankTwoTensorTempl< T > SymmetricRankTwoTensorTempl< T >::transposeTimes ( const SymmetricRankTwoTensorTempl< T > &  a)
static

Definition at line 333 of file SymmetricRankTwoTensorImplementation.h.

334 {
335  return timesTranspose(a);
336 }
static SymmetricRankTwoTensorTempl< T > timesTranspose(const RankTwoTensorTempl< T > &)
return the matrix multiplied with its transpose A*A^T (guaranteed symmetric)

◆ zero()

template<typename T >
void SymmetricRankTwoTensorTempl< T >::zero ( )

Set all components to zero.

Definition at line 363 of file SymmetricRankTwoTensorImplementation.h.

Referenced by SymmetricRankTwoTensorTempl< Real >::operator=().

364 {
365  for (std::size_t i = 0; i < N; ++i)
366  _vals[i] = 0.0;
367 }
static constexpr unsigned int N

Friends And Related Function Documentation

◆ dataLoad

template<typename T>
template<class T2 >
void dataLoad ( std::istream &  ,
SymmetricRankTwoTensorTempl< T2 > &  ,
void  
)
friend

◆ dataStore

template<typename T>
template<class T2 >
void dataStore ( std::ostream &  ,
SymmetricRankTwoTensorTempl< T2 > &  ,
void  
)
friend

◆ operator<<

template<typename T>
std::ostream& operator<< ( std::ostream &  os,
const SymmetricRankTwoTensorTempl< T > &  t 
)
friend

Definition at line 422 of file SymmetricRankTwoTensor.h.

423  {
424  t.print(os);
425  return os;
426  }
std::basic_ostream< charT, traits > * os
Definition: InfixIterator.h:33
void print(std::ostream &stm=Moose::out) const
Print the rank two tensor.

◆ SymmetricRankFourTensorTempl

template<typename T>
template<class T2 >
friend class SymmetricRankFourTensorTempl
friend

Definition at line 500 of file SymmetricRankTwoTensor.h.

Member Data Documentation

◆ _vals

template<typename T>
std::array<T, N> SymmetricRankTwoTensorTempl< T >::_vals
private

◆ full_index

template<typename T>
constexpr unsigned int SymmetricRankTwoTensorTempl< T >::full_index[6][2] = {{0, 0}, {1, 1}, {2, 2}, {1, 2}, {0, 2}, {0, 1}}
static

Definition at line 87 of file SymmetricRankTwoTensor.h.

◆ identityCoords

template<typename T>
constexpr std::array< Real, SymmetricRankTwoTensorTempl< T >::N > SymmetricRankTwoTensorTempl< T >::identityCoords = {{1, 1, 1, 0, 0, 0}}
staticprivate

Definition at line 489 of file SymmetricRankTwoTensor.h.

◆ N

template<typename T>
constexpr unsigned int SymmetricRankTwoTensorTempl< T >::N = Ndim + (Ndim * (Ndim - 1)) / 2
static

◆ Ndim

template<typename T>
constexpr unsigned int SymmetricRankTwoTensorTempl< T >::Ndim = Moose::dim
static

tensor dimension and Mandel vector length

Definition at line 82 of file SymmetricRankTwoTensor.h.

Referenced by SymmetricRankTwoTensorTempl< Real >::mandelFactor().

◆ reverse_index

template<typename T>
constexpr unsigned int SymmetricRankTwoTensorTempl< T >::reverse_index[3][3] = {{0, 5, 4}, {5, 1, 3}, {4, 3, 2}}
static

The documentation for this class was generated from the following files: