libMesh
Public Types | Public Member Functions | Protected Attributes | List of all members
libMesh::TensorValue< T > Class Template Reference

This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space. More...

#include <exact_solution.h>

Inheritance diagram for libMesh::TensorValue< T >:
[legend]

Public Types

typedef T value_type
 Helper typedef for C++98 generic programming. More...
 
typedef std::tuple< unsigned int, unsigned intindex_type
 Helper typedef for generic index programming. More...
 

Public Member Functions

 TensorValue ()
 Empty constructor. More...
 
 TensorValue (const T &xx, const T &xy=0, const T &xz=0, const T &yx=0, const T &yy=0, const T &yz=0, const T &zx=0, const T &zy=0, const T &zz=0)
 Constructor-from-T. More...
 
template<typename Scalar >
 TensorValue (const Scalar &xx, const Scalar &xy=0, const Scalar &xz=0, const Scalar &yx=0, const Scalar &yy=0, const Scalar &yz=0, const Scalar &zx=0, const Scalar &zy=0, typename boostcopy::enable_if_c< ScalarTraits< Scalar >::value, const Scalar >::type &zz=0)
 Constructor-from-scalars. More...
 
template<typename T2 >
 TensorValue (const TypeVector< T2 > &vx)
 Constructor. More...
 
template<typename T2 >
 TensorValue (const TypeVector< T2 > &vx, const TypeVector< T2 > &vy)
 Constructor. More...
 
template<typename T2 >
 TensorValue (const TypeVector< T2 > &vx, const TypeVector< T2 > &vy, const TypeVector< T2 > &vz)
 Constructor. More...
 
template<typename T2 >
 TensorValue (const TensorValue< T2 > &p)
 Copy-constructor. More...
 
template<typename T2 >
 TensorValue (const TypeTensor< T2 > &p)
 Copy-constructor. More...
 
 TensorValue (const TypeTensor< Real > &p_re, const TypeTensor< Real > &p_im)
 Constructor that takes two TypeTensor<Real> representing the real and imaginary part as arguments. More...
 
template<typename Scalar >
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TensorValue & >::type operator= (const Scalar &libmesh_dbg_var(p))
 Assignment-from-scalar operator. More...
 
template<typename T2 >
void assign (const TypeTensor< T2 > &)
 Assign to this tensor without creating a temporary. More...
 
const T & operator() (const unsigned int i, const unsigned int j) const
 
T & operator() (const unsigned int i, const unsigned int j)
 
ConstTypeTensorColumn< T > slice (const unsigned int i) const
 
TypeTensorColumn< T > slice (const unsigned int i)
 
TypeVector< T > row (const unsigned int r) const
 
template<typename T2 >
TypeTensor< typename CompareTypes< T, T2 >::supertype > operator+ (const TypeTensor< T2 > &) const
 Add another tensor to this tensor. More...
 
template<typename T2 >
const TypeTensor< T > & operator+= (const TypeTensor< T2 > &)
 Add to this tensor. More...
 
template<typename T2 >
void add (const TypeTensor< T2 > &)
 Add to this tensor without creating a temporary. More...
 
template<typename T2 >
void add_scaled (const TypeTensor< T2 > &, const T &)
 Add a scaled tensor to this tensor without creating a temporary. More...
 
template<typename T2 >
TypeTensor< typename CompareTypes< T, T2 >::supertype > operator- (const TypeTensor< T2 > &) const
 Subtract a tensor from this tensor. More...
 
TypeTensor< T > operator- () const
 
template<typename T2 >
const TypeTensor< T > & operator-= (const TypeTensor< T2 > &)
 Subtract from this tensor. More...
 
template<typename T2 >
void subtract (const TypeTensor< T2 > &)
 Subtract from this tensor without creating a temporary. More...
 
template<typename T2 >
void subtract_scaled (const TypeTensor< T2 > &, const T &)
 Subtract a scaled value from this tensor without creating a temporary. More...
 
template<typename Scalar >
auto operator* (const Scalar &scalar) const -> typename boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeTensor< decltype(T() *scalar)>>::type
 Multiply this tensor by a scalar value. More...
 
template<typename T2 >
TypeTensor< typename CompareTypes< T, T2 >::supertype > operator* (const TypeTensor< T2 > &) const
 Multiply 2 tensors together, i.e. More...
 
template<typename T2 >
TypeVector< typename CompareTypes< T, T2 >::supertype > operator* (const TypeVector< T2 > &) const
 Right-multiply this tensor by a vector, i.e. More...
 
template<typename Scalar , typename boostcopy::enable_if_c< ScalarTraits< Scalar >::value, int >::type = 0>
const TypeTensor< T > & operator*= (const Scalar &factor)
 Multiply this tensor by a scalar value in place. More...
 
template<typename T2 >
const TypeTensor< T > & operator*= (const TypeTensor< T2 > &)
 Multiply this tensor by a tensor value in place. More...
 
template<typename Scalar >
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeTensor< typename CompareTypes< T, Scalar >::supertype > >::type operator/ (const Scalar &) const
 Divide each entry of this tensor by a scalar value. More...
 
const TypeTensor< T > & operator/= (const T &)
 Divide each entry of this tensor by a scalar value. More...
 
template<typename T2 >
CompareTypes< T, T2 >::supertype contract (const TypeTensor< T2 > &) const
 Multiply 2 tensors together to return a scalar, i.e. More...
 
template<typename T2 >
TypeVector< typename CompareTypes< T, T2 >::supertype > left_multiply (const TypeVector< T2 > &p) const
 Left-multiply this tensor by a vector, i.e. More...
 
TypeTensor< T > transpose () const
 
TypeTensor< T > inverse () const
 
void solve (const TypeVector< T > &b, TypeVector< T > &x) const
 Solve the 2x2 or 3x3 system of equations A*x = b for x by directly inverting A and multiplying it by b. More...
 
auto size () const -> decltype(std::norm(T()))
 
auto norm () const -> decltype(std::norm(T()))
 
auto size_sq () const -> decltype(std::norm(T()))
 
auto norm_sq () const -> decltype(std::norm(T()))
 
bool is_zero () const
 
det () const
 
tr () const
 
void zero ()
 Set all entries of the tensor to 0. More...
 
bool operator== (const TypeTensor< T > &rhs) const
 
bool operator< (const TypeTensor< T > &rhs) const
 
template<>
bool operator< (const TypeTensor< Real > &rhs) const
 
template<>
bool operator< (const TypeTensor< Complex > &rhs) const
 
bool operator> (const TypeTensor< T > &rhs) const
 
template<>
bool operator> (const TypeTensor< Real > &rhs) const
 
template<>
bool operator> (const TypeTensor< Complex > &rhs) const
 
void print (std::ostream &os=libMesh::out) const
 Formatted print to a stream which defaults to libMesh::out. More...
 
void write_unformatted (std::ostream &out, const bool newline=true) const
 Unformatted print to the stream out. More...
 

Protected Attributes

_coords [LIBMESH_DIM *LIBMESH_DIM]
 The coordinates of the TypeTensor. More...
 

Detailed Description

template<typename T>
class libMesh::TensorValue< T >

This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.

The typedef RealTensorValue always defines a real-valued tensor, and NumberTensorValue defines a real or complex-valued tensor depending on how the library was configured.

Author
Roy H. Stogner
Date
2004

Definition at line 53 of file exact_solution.h.

Member Typedef Documentation

◆ index_type

template<typename T>
typedef std::tuple<unsigned int, unsigned int> libMesh::TypeTensor< T >::index_type
inherited

Helper typedef for generic index programming.

Definition at line 144 of file type_tensor.h.

◆ value_type

template<typename T>
typedef T libMesh::TypeTensor< T >::value_type
inherited

Helper typedef for C++98 generic programming.

Definition at line 139 of file type_tensor.h.

Constructor & Destructor Documentation

◆ TensorValue() [1/9]

template<typename T >
libMesh::TensorValue< T >::TensorValue ( )
inline

Empty constructor.

Gives the tensor 0 in LIBMESH_DIM dimensional T space.

Definition at line 156 of file tensor_value.h.

156  :
157  TypeTensor<T> ()
158 {
159 }

◆ TensorValue() [2/9]

template<typename T >
libMesh::TensorValue< T >::TensorValue ( const T &  xx,
const T &  xy = 0,
const T &  xz = 0,
const T &  yx = 0,
const T &  yy = 0,
const T &  yz = 0,
const T &  zx = 0,
const T &  zy = 0,
const T &  zz = 0 
)
inlineexplicit

Constructor-from-T.

By default sets higher dimensional entries to 0.

Definition at line 165 of file tensor_value.h.

173  :
174  TypeTensor<T> (xx,xy,xz,yx,yy,yz,zx,zy,zz)
175 {
176 }

◆ TensorValue() [3/9]

template<typename T >
template<typename Scalar >
libMesh::TensorValue< T >::TensorValue ( const Scalar &  xx,
const Scalar &  xy = 0,
const Scalar &  xz = 0,
const Scalar &  yx = 0,
const Scalar &  yy = 0,
const Scalar &  yz = 0,
const Scalar &  zx = 0,
const Scalar &  zy = 0,
typename boostcopy::enable_if_c< ScalarTraits< Scalar >::value, const Scalar >::type &  zz = 0 
)
inlineexplicit

Constructor-from-scalars.

By default sets higher dimensional entries to 0.

Definition at line 182 of file tensor_value.h.

192  :
193  TypeTensor<T> (xx,xy,xz,yx,yy,yz,zx,zy,zz)
194 {
195 }

◆ TensorValue() [4/9]

template<typename T >
template<typename T2 >
libMesh::TensorValue< T >::TensorValue ( const TypeVector< T2 > &  vx)
inline

Constructor.

Takes 1 row vector for LIBMESH_DIM=1

Definition at line 212 of file tensor_value.h.

212  :
213  TypeTensor<T> (vx)
214 {
215 }

◆ TensorValue() [5/9]

template<typename T >
template<typename T2 >
libMesh::TensorValue< T >::TensorValue ( const TypeVector< T2 > &  vx,
const TypeVector< T2 > &  vy 
)
inline

Constructor.

Takes 2 row vectors for LIBMESH_DIM=2

Definition at line 222 of file tensor_value.h.

223  :
224  TypeTensor<T> (vx, vy)
225 {
226 }

◆ TensorValue() [6/9]

template<typename T >
template<typename T2 >
libMesh::TensorValue< T >::TensorValue ( const TypeVector< T2 > &  vx,
const TypeVector< T2 > &  vy,
const TypeVector< T2 > &  vz 
)
inline

Constructor.

Takes 3 row vectors for LIBMESH_DIM=3

Definition at line 233 of file tensor_value.h.

235  :
236  TypeTensor<T> (vx, vy, vz)
237 {
238 }

◆ TensorValue() [7/9]

template<typename T >
template<typename T2 >
libMesh::TensorValue< T >::TensorValue ( const TensorValue< T2 > &  p)
inline

Copy-constructor.

Definition at line 202 of file tensor_value.h.

202  :
203  TypeTensor<T> (p)
204 {
205 }

◆ TensorValue() [8/9]

template<typename T >
template<typename T2 >
libMesh::TensorValue< T >::TensorValue ( const TypeTensor< T2 > &  p)
inline

Copy-constructor.

Definition at line 245 of file tensor_value.h.

245  :
246  TypeTensor<T> (p)
247 {
248 }

◆ TensorValue() [9/9]

template<typename T >
libMesh::TensorValue< T >::TensorValue ( const TypeTensor< Real > &  p_re,
const TypeTensor< Real > &  p_im 
)
inline

Constructor that takes two TypeTensor<Real> representing the real and imaginary part as arguments.

Definition at line 254 of file tensor_value.h.

255  :
256  TypeTensor<T> (Complex (p_re(0,0), p_im(0,0)),
257  Complex (p_re(0,1), p_im(0,1)),
258  Complex (p_re(0,2), p_im(0,2)),
259  Complex (p_re(1,0), p_im(1,0)),
260  Complex (p_re(1,1), p_im(1,1)),
261  Complex (p_re(1,2), p_im(1,2)),
262  Complex (p_re(2,0), p_im(2,0)),
263  Complex (p_re(2,1), p_im(2,1)),
264  Complex (p_re(2,2), p_im(2,2)))
265 {
266 }

Member Function Documentation

◆ add()

template<typename T >
template<typename T2 >
void libMesh::TypeTensor< T >::add ( const TypeTensor< T2 > &  p)
inlineinherited

Add to this tensor without creating a temporary.

Definition at line 858 of file type_tensor.h.

859 {
860  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
861  _coords[i] += p._coords[i];
862 }

References libMesh::TypeTensor< T >::_coords.

◆ add_scaled()

template<typename T >
template<typename T2 >
void libMesh::TypeTensor< T >::add_scaled ( const TypeTensor< T2 > &  p,
const T &  factor 
)
inlineinherited

Add a scaled tensor to this tensor without creating a temporary.

Definition at line 869 of file type_tensor.h.

870 {
871  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
872  _coords[i] += factor*p._coords[i];
873 
874 }

References libMesh::TypeTensor< T >::_coords.

Referenced by libMesh::HPCoarsenTest::add_projection(), libMesh::System::calculate_norm(), libMesh::MeshFunction::hessian(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::System::point_hessian(), and libMesh::HPCoarsenTest::select_refinement().

◆ assign()

template<typename T >
template<typename T2 >
void libMesh::TypeTensor< T >::assign ( const TypeTensor< T2 > &  p)
inlineinherited

Assign to this tensor without creating a temporary.

Definition at line 729 of file type_tensor.h.

730 {
731  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
732  _coords[i] = p._coords[i];
733 }

References libMesh::TypeTensor< T >::_coords.

◆ contract()

template<typename T >
template<typename T2 >
CompareTypes< T, T2 >::supertype libMesh::TypeTensor< T >::contract ( const TypeTensor< T2 > &  t) const
inlineinherited

Multiply 2 tensors together to return a scalar, i.e.

Multiply 2 tensors together, i.e.

\( \sum_{ij} A_{ij} B_{ij} \) The tensors may contain different numeric types. Also known as the "double inner product" or "double dot product" of tensors.

Returns
The scalar-valued result, this tensor is unchanged.

sum Aij*Bij. The tensors may be of different types.

Definition at line 1286 of file type_tensor.h.

1287 {
1288  typename CompareTypes<T,T2>::supertype sum = 0.;
1289  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1290  sum += _coords[i]*t._coords[i];
1291  return sum;
1292 }

References libMesh::TypeTensor< T >::_coords.

Referenced by libMesh::HPCoarsenTest::add_projection(), libMesh::TensorTools::inner_product(), and libMesh::HPCoarsenTest::select_refinement().

◆ det()

template<typename T >
T libMesh::TypeTensor< T >::det ( ) const
inlineinherited
Returns
The determinant of the tensor.

Because these are 3x3 tensors at most, we don't do an LU decomposition like DenseMatrix does.

Definition at line 1328 of file type_tensor.h.

1329 {
1330 #if LIBMESH_DIM == 1
1331  return _coords[0];
1332 #endif
1333 
1334 #if LIBMESH_DIM == 2
1335  return (_coords[0] * _coords[3]
1336  - _coords[1] * _coords[2]);
1337 #endif
1338 
1339 #if LIBMESH_DIM == 3
1340  return (_coords[0] * _coords[4] * _coords[8]
1341  + _coords[1] * _coords[5] * _coords[6]
1342  + _coords[2] * _coords[3] * _coords[7]
1343  - _coords[0] * _coords[5] * _coords[7]
1344  - _coords[1] * _coords[3] * _coords[8]
1345  - _coords[2] * _coords[4] * _coords[6]);
1346 #endif
1347 }

Referenced by libMesh::Sphere::Sphere().

◆ inverse()

template<typename T >
TypeTensor< T > libMesh::TypeTensor< T >::inverse ( ) const
inlineinherited
Returns
The inverse of this tensor as an independent object.

Definition at line 1098 of file type_tensor.h.

1099 {
1100 #if LIBMESH_DIM == 1
1101  if (_coords[0] == static_cast<T>(0.))
1102  libmesh_convergence_failure();
1103  return TypeTensor(1. / _coords[0]);
1104 #endif
1105 
1106 #if LIBMESH_DIM == 2
1107  // Get convenient reference to this.
1108  const TypeTensor<T> & A = *this;
1109 
1110  // Use temporary variables, avoid multiple accesses
1111  T a = A(0,0), b = A(0,1),
1112  c = A(1,0), d = A(1,1);
1113 
1114  // Make sure det = ad - bc is not zero
1115  T my_det = a*d - b*c;
1116 
1117  if (my_det == static_cast<T>(0.))
1118  libmesh_convergence_failure();
1119 
1120  return TypeTensor(d/my_det, -b/my_det, -c/my_det, a/my_det);
1121 #endif
1122 
1123 #if LIBMESH_DIM == 3
1124  // Get convenient reference to this.
1125  const TypeTensor<T> & A = *this;
1126 
1127  T a11 = A(0,0), a12 = A(0,1), a13 = A(0,2),
1128  a21 = A(1,0), a22 = A(1,1), a23 = A(1,2),
1129  a31 = A(2,0), a32 = A(2,1), a33 = A(2,2);
1130 
1131  T my_det = a11*(a33*a22-a32*a23) - a21*(a33*a12-a32*a13) + a31*(a23*a12-a22*a13);
1132 
1133  if (my_det == static_cast<T>(0.))
1134  libmesh_convergence_failure();
1135 
1136  // Inline comment characters are for lining up columns.
1137  return TypeTensor( (a33*a22-a32*a23)/my_det, -(a33*a12-a32*a13)/my_det, (a23*a12-a22*a13)/my_det,
1138  -(a33*a21-a31*a23)/my_det, (a33*a11-a31*a13)/my_det, -(a23*a11-a21*a13)/my_det,
1139  (a32*a21-a31*a22)/my_det, -(a32*a11-a31*a12)/my_det, (a22*a11-a21*a12)/my_det);
1140 #endif
1141 }

References A.

Referenced by NonlinearNeoHookeCurrentConfig::calculate_stress(), and TypeTensorTest::testInverse().

◆ is_zero()

template<typename T >
bool libMesh::TypeTensor< T >::is_zero ( ) const
inlineinherited
Returns
True if all values in the tensor are zero

Definition at line 1318 of file type_tensor.h.

1319 {
1320  for (const auto & val : _coords)
1321  if (val != T(0))
1322  return false;
1323  return true;
1324 }

Referenced by TypeTensorTest::testIsZero().

◆ left_multiply()

template<typename T >
template<typename T2 >
TypeVector< typename CompareTypes< T, T2 >::supertype > libMesh::TypeTensor< T >::left_multiply ( const TypeVector< T2 > &  p) const
inlineinherited

Left-multiply this tensor by a vector, i.e.

matrix-vector product. The tensor and vector may contain different numeric types.

Returns
A copy of the result vector, this tensor is unchanged.

Definition at line 1229 of file type_tensor.h.

1230 {
1231  TypeVector<typename CompareTypes<T,T2>::supertype> returnval;
1232  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1233  for (unsigned int j=0; j<LIBMESH_DIM; j++)
1234  returnval(i) += p(j)*(*this)(j,i);
1235 
1236  return returnval;
1237 }

Referenced by libMesh::operator*().

◆ norm()

template<typename T >
auto libMesh::TypeTensor< T >::norm ( ) const -> decltype(std::norm(T()))
inlineinherited
Returns
The Frobenius norm of the tensor, i.e. the square-root of the sum of the elements squared.

Definition at line 1310 of file type_tensor.h.

1311 {
1312  return std::sqrt(this->norm_sq());
1313 }

References libMesh::TensorTools::norm_sq(), and std::sqrt().

Referenced by libMesh::System::calculate_norm().

◆ norm_sq()

template<typename T >
auto libMesh::TypeTensor< T >::norm_sq ( ) const -> decltype(std::norm(T()))
inlineinherited
Returns
The Frobenius norm of the tensor squared, i.e. sum of the element magnitudes squared.

Definition at line 1390 of file type_tensor.h.

1391 {
1392  Real sum = 0.;
1393  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1394  sum += TensorTools::norm_sq(_coords[i]);
1395  return sum;
1396 }

References libMesh::TensorTools::norm_sq(), and libMesh::Real.

Referenced by libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::System::calculate_norm(), libMesh::ExactErrorEstimator::find_squared_element_error(), and libMesh::HPCoarsenTest::select_refinement().

◆ operator()() [1/2]

template<typename T >
T & libMesh::TypeTensor< T >::operator() ( const unsigned int  i,
const unsigned int  j 
)
inlineinherited
Returns
A writable reference to the (i,j) entry of the tensor.

Definition at line 758 of file type_tensor.h.

760 {
761 #if LIBMESH_DIM < 3
762 
763  if (i >= LIBMESH_DIM || j >= LIBMESH_DIM)
764  libmesh_error_msg("ERROR: You are assigning to a tensor component that is out of range for the compiled LIBMESH_DIM!");
765 
766 #endif
767 
768  libmesh_assert_less (i, LIBMESH_DIM);
769  libmesh_assert_less (j, LIBMESH_DIM);
770 
771  return _coords[i*LIBMESH_DIM+j];
772 }

◆ operator()() [2/2]

template<typename T >
const T & libMesh::TypeTensor< T >::operator() ( const unsigned int  i,
const unsigned int  j 
) const
inlineinherited
Returns
A const reference to the (i,j) entry of the tensor.

Definition at line 739 of file type_tensor.h.

741 {
742  libmesh_assert_less (i, 3);
743  libmesh_assert_less (j, 3);
744 
745 #if LIBMESH_DIM < 3
746  const static T my_zero = 0;
747  if (i >= LIBMESH_DIM || j >= LIBMESH_DIM)
748  return my_zero;
749 #endif
750 
751  return _coords[i*LIBMESH_DIM+j];
752 }

◆ operator*() [1/3]

template<typename T >
template<typename Scalar >
auto libMesh::TypeTensor< T >::operator* ( const Scalar &  scalar) const -> typename boostcopy::enable_if_c< ScalarTraits<Scalar>::value, TypeTensor<decltype(T() * scalar)>>::type
inlineinherited

Multiply this tensor by a scalar value.

Returns
A copy of the result, this tensor is unchanged.

Definition at line 983 of file type_tensor.h.

986 {
987  typedef decltype((*this)(0, 0) * factor) TS;
988 
989 
990 #if LIBMESH_DIM == 1
991  return TypeTensor<TS>(_coords[0]*factor);
992 #endif
993 
994 #if LIBMESH_DIM == 2
995  return TypeTensor<TS>(_coords[0]*factor,
996  _coords[1]*factor,
997  _coords[2]*factor,
998  _coords[3]*factor);
999 #endif
1000 
1001 #if LIBMESH_DIM == 3
1002  return TypeTensor<TS>(_coords[0]*factor,
1003  _coords[1]*factor,
1004  _coords[2]*factor,
1005  _coords[3]*factor,
1006  _coords[4]*factor,
1007  _coords[5]*factor,
1008  _coords[6]*factor,
1009  _coords[7]*factor,
1010  _coords[8]*factor);
1011 #endif
1012 }

◆ operator*() [2/3]

template<typename T >
template<typename T2 >
TypeTensor< typename CompareTypes< T, T2 >::supertype > libMesh::TypeTensor< T >::operator* ( const TypeTensor< T2 > &  p) const
inlineinherited

Multiply 2 tensors together, i.e.

matrix-matrix product. The tensors may contain different numeric types.

Returns
A copy of the result, this tensor is unchanged.

Definition at line 1251 of file type_tensor.h.

1252 {
1253  TypeTensor<typename CompareTypes<T, T2>::supertype> returnval;
1254  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1255  for (unsigned int j=0; j<LIBMESH_DIM; j++)
1256  for (unsigned int k=0; k<LIBMESH_DIM; k++)
1257  returnval(i,j) += (*this)(i,k)*p(k,j);
1258 
1259  return returnval;
1260 }

◆ operator*() [3/3]

template<typename T >
template<typename T2 >
TypeVector< typename CompareTypes< T, T2 >::supertype > libMesh::TypeTensor< T >::operator* ( const TypeVector< T2 > &  p) const
inlineinherited

Right-multiply this tensor by a vector, i.e.

matrix-vector product. The tensor and vector may contain different numeric types.

Returns
A copy of the result vector, this tensor is unchanged.

Definition at line 1215 of file type_tensor.h.

1216 {
1217  TypeVector<typename CompareTypes<T,T2>::supertype> returnval;
1218  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1219  for (unsigned int j=0; j<LIBMESH_DIM; j++)
1220  returnval(i) += (*this)(i,j)*p(j);
1221 
1222  return returnval;
1223 }

◆ operator*=() [1/2]

template<typename T>
template<typename Scalar , typename boostcopy::enable_if_c< ScalarTraits< Scalar >::value, int >::type = 0>
const TypeTensor<T>& libMesh::TypeTensor< T >::operator*= ( const Scalar &  factor)
inlineinherited

Multiply this tensor by a scalar value in place.

Returns
A reference to *this.

Definition at line 282 of file type_tensor.h.

283  {
284  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
285  _coords[i] *= factor;
286 
287  return *this;
288  }

References libMesh::TypeTensor< T >::_coords.

◆ operator*=() [2/2]

template<typename T >
template<typename T2 >
const TypeTensor< T > & libMesh::TypeTensor< T >::operator*= ( const TypeTensor< T2 > &  p)
inlineinherited

Multiply this tensor by a tensor value in place.

Returns
A reference to *this

Definition at line 1265 of file type_tensor.h.

1266 {
1267  TypeTensor<T> temp;
1268  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1269  for (unsigned int j=0; j<LIBMESH_DIM; j++)
1270  for (unsigned int k=0; k<LIBMESH_DIM; k++)
1271  temp(i,j) += (*this)(i,k)*p(k,j);
1272 
1273  this->assign(temp);
1274  return *this;
1275 }

◆ operator+()

template<typename T >
template<typename T2 >
TypeTensor< typename CompareTypes< T, T2 >::supertype > libMesh::TypeTensor< T >::operator+ ( const TypeTensor< T2 > &  p) const
inlineinherited

Add another tensor to this tensor.

Returns
A copy of the result, this tensor is unchanged.

Definition at line 812 of file type_tensor.h.

813 {
814 
815 #if LIBMESH_DIM == 1
816  return TypeTensor<typename CompareTypes<T, T2>::supertype>(_coords[0] + p._coords[0]);
817 #endif
818 
819 #if LIBMESH_DIM == 2
820  return TypeTensor<typename CompareTypes<T, T2>::supertype>(_coords[0] + p._coords[0],
821  _coords[1] + p._coords[1],
822  0.,
823  _coords[2] + p._coords[2],
824  _coords[3] + p._coords[3]);
825 #endif
826 
827 #if LIBMESH_DIM == 3
828  return TypeTensor<typename CompareTypes<T, T2>::supertype>(_coords[0] + p._coords[0],
829  _coords[1] + p._coords[1],
830  _coords[2] + p._coords[2],
831  _coords[3] + p._coords[3],
832  _coords[4] + p._coords[4],
833  _coords[5] + p._coords[5],
834  _coords[6] + p._coords[6],
835  _coords[7] + p._coords[7],
836  _coords[8] + p._coords[8]);
837 #endif
838 
839 }

References libMesh::TypeTensor< T >::_coords.

◆ operator+=()

template<typename T >
template<typename T2 >
const TypeTensor< T > & libMesh::TypeTensor< T >::operator+= ( const TypeTensor< T2 > &  p)
inlineinherited

Add to this tensor.

Returns
A reference to *this.

Definition at line 846 of file type_tensor.h.

847 {
848  this->add (p);
849 
850  return *this;
851 }

◆ operator-() [1/2]

template<typename T >
TypeTensor< T > libMesh::TypeTensor< T >::operator- ( ) const
inlineinherited
Returns
The negative of this tensor in a separate copy.

Definition at line 949 of file type_tensor.h.

950 {
951 
952 #if LIBMESH_DIM == 1
953  return TypeTensor(-_coords[0]);
954 #endif
955 
956 #if LIBMESH_DIM == 2
957  return TypeTensor(-_coords[0],
958  -_coords[1],
959  -_coords[2],
960  -_coords[3]);
961 #endif
962 
963 #if LIBMESH_DIM == 3
964  return TypeTensor(-_coords[0],
965  -_coords[1],
966  -_coords[2],
967  -_coords[3],
968  -_coords[4],
969  -_coords[5],
970  -_coords[6],
971  -_coords[7],
972  -_coords[8]);
973 #endif
974 
975 }

◆ operator-() [2/2]

template<typename T >
template<typename T2 >
TypeTensor< typename CompareTypes< T, T2 >::supertype > libMesh::TypeTensor< T >::operator- ( const TypeTensor< T2 > &  p) const
inlineinherited

Subtract a tensor from this tensor.

Returns
A copy of the result, this tensor is unchanged.

Definition at line 882 of file type_tensor.h.

883 {
884 
885 #if LIBMESH_DIM == 1
886  return TypeTensor<typename CompareTypes<T, T2>::supertype>(_coords[0] - p._coords[0]);
887 #endif
888 
889 #if LIBMESH_DIM == 2
890  return TypeTensor<typename CompareTypes<T, T2>::supertype>(_coords[0] - p._coords[0],
891  _coords[1] - p._coords[1],
892  0.,
893  _coords[2] - p._coords[2],
894  _coords[3] - p._coords[3]);
895 #endif
896 
897 #if LIBMESH_DIM == 3
898  return TypeTensor<typename CompareTypes<T, T2>::supertype>(_coords[0] - p._coords[0],
899  _coords[1] - p._coords[1],
900  _coords[2] - p._coords[2],
901  _coords[3] - p._coords[3],
902  _coords[4] - p._coords[4],
903  _coords[5] - p._coords[5],
904  _coords[6] - p._coords[6],
905  _coords[7] - p._coords[7],
906  _coords[8] - p._coords[8]);
907 #endif
908 
909 }

References libMesh::TypeTensor< T >::_coords.

◆ operator-=()

template<typename T >
template<typename T2 >
const TypeTensor< T > & libMesh::TypeTensor< T >::operator-= ( const TypeTensor< T2 > &  p)
inlineinherited

Subtract from this tensor.

Returns
A reference to *this.

Definition at line 916 of file type_tensor.h.

917 {
918  this->subtract (p);
919 
920  return *this;
921 }

◆ operator/()

template<typename T >
template<typename Scalar >
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeTensor< typename CompareTypes< T, Scalar >::supertype > >::type libMesh::TypeTensor< T >::operator/ ( const Scalar &  factor) const
inlineinherited

Divide each entry of this tensor by a scalar value.

Returns
A copy of the result, this tensor is unchanged.

Definition at line 1033 of file type_tensor.h.

1034 {
1035  libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
1036 
1037  typedef typename CompareTypes<T, Scalar>::supertype TS;
1038 
1039 #if LIBMESH_DIM == 1
1040  return TypeTensor<TS>(_coords[0]/factor);
1041 #endif
1042 
1043 #if LIBMESH_DIM == 2
1044  return TypeTensor<TS>(_coords[0]/factor,
1045  _coords[1]/factor,
1046  _coords[2]/factor,
1047  _coords[3]/factor);
1048 #endif
1049 
1050 #if LIBMESH_DIM == 3
1051  return TypeTensor<TS>(_coords[0]/factor,
1052  _coords[1]/factor,
1053  _coords[2]/factor,
1054  _coords[3]/factor,
1055  _coords[4]/factor,
1056  _coords[5]/factor,
1057  _coords[6]/factor,
1058  _coords[7]/factor,
1059  _coords[8]/factor);
1060 #endif
1061 
1062 }

◆ operator/=()

template<typename T >
const TypeTensor< T > & libMesh::TypeTensor< T >::operator/= ( const T &  factor)
inlineinherited

Divide each entry of this tensor by a scalar value.

Returns
A reference to *this.

Definition at line 1198 of file type_tensor.h.

1199 {
1200  libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
1201 
1202  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1203  _coords[i] /= factor;
1204 
1205  return *this;
1206 }

◆ operator<() [1/3]

template<>
bool libMesh::TypeTensor< Complex >::operator< ( const TypeTensor< Complex > &  rhs) const
inherited

Definition at line 146 of file type_tensor.C.

147 {
148  for (unsigned int i=0; i<LIBMESH_DIM; i++)
149  for (unsigned int j=0; j<LIBMESH_DIM; j++)
150  {
151  if ((*this)(i,j).real() < rhs(i,j).real())
152  return true;
153  if ((*this)(i,j).real() > rhs(i,j).real())
154  return false;
155  if ((*this)(i,j).imag() < rhs(i,j).imag())
156  return true;
157  if ((*this)(i,j).imag() > rhs(i,j).imag())
158  return false;
159  }
160  return false;
161 }

References std::imag(), and std::real().

◆ operator<() [2/3]

template<>
bool libMesh::TypeTensor< Real >::operator< ( const TypeTensor< Real > &  rhs) const
inherited

Definition at line 113 of file type_tensor.C.

114 {
115  for (unsigned int i=0; i<LIBMESH_DIM; i++)
116  for (unsigned int j=0; j<LIBMESH_DIM; j++)
117  {
118  if ((*this)(i,j) < rhs(i,j))
119  return true;
120  if ((*this)(i,j) > rhs(i,j))
121  return false;
122  }
123  return false;
124 }

◆ operator<() [3/3]

template<typename T>
bool libMesh::TypeTensor< T >::operator< ( const TypeTensor< T > &  rhs) const
inherited
Returns
true if this tensor is "less" than another. Useful for sorting.

◆ operator=()

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

Assignment-from-scalar operator.

Used only to zero out tensors.

Definition at line 135 of file tensor_value.h.

136  { libmesh_assert_equal_to (p, Scalar(0)); this->zero(); return *this; }

References libMesh::TypeTensor< T >::zero().

◆ operator==()

template<typename T >
bool libMesh::TypeTensor< T >::operator== ( const TypeTensor< T > &  rhs) const
inlineinherited
Returns
true if two tensors are equal, false otherwise.

Definition at line 1402 of file type_tensor.h.

1403 {
1404 #if LIBMESH_DIM == 1
1405  return (std::abs(_coords[0] - rhs._coords[0])
1406  < TOLERANCE);
1407 #endif
1408 
1409 #if LIBMESH_DIM == 2
1410  return ((std::abs(_coords[0] - rhs._coords[0]) +
1411  std::abs(_coords[1] - rhs._coords[1]) +
1412  std::abs(_coords[2] - rhs._coords[2]) +
1413  std::abs(_coords[3] - rhs._coords[3]))
1414  < 4.*TOLERANCE);
1415 #endif
1416 
1417 #if LIBMESH_DIM == 3
1418  return ((std::abs(_coords[0] - rhs._coords[0]) +
1419  std::abs(_coords[1] - rhs._coords[1]) +
1420  std::abs(_coords[2] - rhs._coords[2]) +
1421  std::abs(_coords[3] - rhs._coords[3]) +
1422  std::abs(_coords[4] - rhs._coords[4]) +
1423  std::abs(_coords[5] - rhs._coords[5]) +
1424  std::abs(_coords[6] - rhs._coords[6]) +
1425  std::abs(_coords[7] - rhs._coords[7]) +
1426  std::abs(_coords[8] - rhs._coords[8]))
1427  < 9.*TOLERANCE);
1428 #endif
1429 
1430 }

References libMesh::TypeTensor< T >::_coords, std::abs(), and libMesh::TOLERANCE.

◆ operator>() [1/3]

template<>
bool libMesh::TypeTensor< Complex >::operator> ( const TypeTensor< Complex > &  rhs) const
inherited

Definition at line 166 of file type_tensor.C.

167 {
168  for (unsigned int i=0; i<LIBMESH_DIM; i++)
169  for (unsigned int j=0; j<LIBMESH_DIM; j++)
170  {
171  if ((*this)(i,j).real() > rhs(i,j).real())
172  return true;
173  if ((*this)(i,j).real() < rhs(i,j).real())
174  return false;
175  if ((*this)(i,j).imag() > rhs(i,j).imag())
176  return true;
177  if ((*this)(i,j).imag() < rhs(i,j).imag())
178  return false;
179  }
180  return false;
181 }

References std::imag(), and std::real().

◆ operator>() [2/3]

template<>
bool libMesh::TypeTensor< Real >::operator> ( const TypeTensor< Real > &  rhs) const
inherited

Definition at line 129 of file type_tensor.C.

130 {
131  for (unsigned int i=0; i<LIBMESH_DIM; i++)
132  for (unsigned int j=0; j<LIBMESH_DIM; j++)
133  {
134  if ((*this)(i,j) > rhs(i,j))
135  return true;
136  if ((*this)(i,j) < rhs(i,j))
137  return false;
138  }
139  return false;
140 }

◆ operator>() [3/3]

template<typename T>
bool libMesh::TypeTensor< T >::operator> ( const TypeTensor< T > &  rhs) const
inherited
Returns
true if this tensor is "greater" than another.

◆ print()

template<typename T >
void libMesh::TypeTensor< T >::print ( std::ostream &  os = libMesh::out) const
inherited

Formatted print to a stream which defaults to libMesh::out.

Definition at line 39 of file type_tensor.C.

40 {
41 #if LIBMESH_DIM == 1
42 
43  os << "x=" << (*this)(0,0) << std::endl;
44 
45 #endif
46 #if LIBMESH_DIM == 2
47 
48  os << "(xx,xy)=("
49  << std::setw(8) << (*this)(0,0) << ", "
50  << std::setw(8) << (*this)(0,1) << ")"
51  << std::endl;
52  os << "(yx,yy)=("
53  << std::setw(8) << (*this)(1,0) << ", "
54  << std::setw(8) << (*this)(1,1) << ")"
55  << std::endl;
56 
57 #endif
58 #if LIBMESH_DIM == 3
59 
60  os << "(xx,xy,xz)=("
61  << std::setw(8) << (*this)(0,0) << ", "
62  << std::setw(8) << (*this)(0,1) << ", "
63  << std::setw(8) << (*this)(0,2) << ")"
64  << std::endl;
65  os << "(yx,yy,yz)=("
66  << std::setw(8) << (*this)(1,0) << ", "
67  << std::setw(8) << (*this)(1,1) << ", "
68  << std::setw(8) << (*this)(1,2) << ")"
69  << std::endl;
70  os << "(zx,zy,zz)=("
71  << std::setw(8) << (*this)(2,0) << ", "
72  << std::setw(8) << (*this)(2,1) << ", "
73  << std::setw(8) << (*this)(2,2) << ")"
74  << std::endl;
75 #endif
76 }

◆ row()

template<typename T >
TypeVector< T > libMesh::TypeTensor< T >::row ( const unsigned int  r) const
inlineinherited
Returns
A copy of one row of the tensor as a TypeVector.

Definition at line 798 of file type_tensor.h.

799 {
800  TypeVector<T> return_vector;
801 
802  for (unsigned int j=0; j<LIBMESH_DIM; j++)
803  return_vector._coords[j] = _coords[r*LIBMESH_DIM + j];
804 
805  return return_vector;
806 }

References libMesh::TypeVector< T >::_coords.

◆ size()

template<typename T >
auto libMesh::TypeTensor< T >::size ( ) const -> decltype(std::norm(T()))
inlineinherited
Returns
The Frobenius norm of the tensor, i.e. the square-root of the sum of the elements squared.

Definition at line 1299 of file type_tensor.h.

1300 {
1301  libmesh_deprecated();
1302  return this->norm();
1303 }

References std::norm().

◆ size_sq()

template<typename T >
auto libMesh::TypeTensor< T >::size_sq ( ) const -> decltype(std::norm(T()))
inlineinherited
Returns
The Frobenius norm of the tensor squared, i.e. sum of the element magnitudes squared.

Definition at line 1379 of file type_tensor.h.

1380 {
1381  libmesh_deprecated();
1382  return this->norm_sq();
1383 }

References libMesh::TensorTools::norm_sq().

◆ slice() [1/2]

template<typename T >
TypeTensorColumn< T > libMesh::TypeTensor< T >::slice ( const unsigned int  i)
inlineinherited
Returns
A writable proxy for the \( i^{th} \) column of the tensor.

Definition at line 788 of file type_tensor.h.

789 {
790  libmesh_assert_less (i, LIBMESH_DIM);
791  return TypeTensorColumn<T>(*this, i);
792 }

◆ slice() [2/2]

template<typename T >
ConstTypeTensorColumn< T > libMesh::TypeTensor< T >::slice ( const unsigned int  i) const
inlineinherited
Returns
A proxy for the \( i^{th} \) column of the tensor.

Definition at line 778 of file type_tensor.h.

779 {
780  libmesh_assert_less (i, LIBMESH_DIM);
781  return ConstTypeTensorColumn<T>(*this, i);
782 }

◆ solve()

template<typename T >
void libMesh::TypeTensor< T >::solve ( const TypeVector< T > &  b,
TypeVector< T > &  x 
) const
inlineinherited

Solve the 2x2 or 3x3 system of equations A*x = b for x by directly inverting A and multiplying it by b.

Returns
The solution in the x vector.

Definition at line 1147 of file type_tensor.h.

1148 {
1149 #if LIBMESH_DIM == 1
1150  if (_coords[0] == static_cast<T>(0.))
1151  libmesh_convergence_failure();
1152  x(0) = b(0) / _coords[0];
1153 #endif
1154 
1155 #if LIBMESH_DIM == 2
1156  T my_det = _coords[0]*_coords[3] - _coords[1]*_coords[2];
1157 
1158  if (my_det == static_cast<T>(0.))
1159  libmesh_convergence_failure();
1160 
1161  T my_det_inv = 1./my_det;
1162 
1163  x(0) = my_det_inv*( _coords[3]*b(0) - _coords[1]*b(1));
1164  x(1) = my_det_inv*(-_coords[2]*b(0) + _coords[0]*b(1));
1165 #endif
1166 
1167 #if LIBMESH_DIM == 3
1168  T my_det =
1169  // a11*(a33 *a22 - a32 *a23)
1170  _coords[0]*(_coords[8]*_coords[4] - _coords[7]*_coords[5])
1171  // -a21*(a33 *a12 - a32 *a13)
1172  -_coords[3]*(_coords[8]*_coords[1] - _coords[7]*_coords[2]) +
1173  // +a31*(a23 *a12 - a22 *a13)
1174  +_coords[6]*(_coords[5]*_coords[1] - _coords[4]*_coords[2]);
1175 
1176  if (my_det == static_cast<T>(0.))
1177  libmesh_convergence_failure();
1178 
1179  T my_det_inv = 1./my_det;
1180  x(0) = my_det_inv*((_coords[8]*_coords[4] - _coords[7]*_coords[5])*b(0) -
1181  (_coords[8]*_coords[1] - _coords[7]*_coords[2])*b(1) +
1182  (_coords[5]*_coords[1] - _coords[4]*_coords[2])*b(2));
1183 
1184  x(1) = my_det_inv*((_coords[6]*_coords[5] - _coords[8]*_coords[3])*b(0) +
1185  (_coords[8]*_coords[0] - _coords[6]*_coords[2])*b(1) -
1186  (_coords[5]*_coords[0] - _coords[3]*_coords[2])*b(2));
1187 
1188  x(2) = my_det_inv*((_coords[7]*_coords[3] - _coords[6]*_coords[4])*b(0) -
1189  (_coords[7]*_coords[0] - _coords[6]*_coords[1])*b(1) +
1190  (_coords[4]*_coords[0] - _coords[3]*_coords[1])*b(2));
1191 #endif
1192 }

◆ subtract()

template<typename T >
template<typename T2 >
void libMesh::TypeTensor< T >::subtract ( const TypeTensor< T2 > &  p)
inlineinherited

Subtract from this tensor without creating a temporary.

Definition at line 928 of file type_tensor.h.

929 {
930  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
931  _coords[i] -= p._coords[i];
932 }

References libMesh::TypeTensor< T >::_coords.

◆ subtract_scaled()

template<typename T >
template<typename T2 >
void libMesh::TypeTensor< T >::subtract_scaled ( const TypeTensor< T2 > &  p,
const T &  factor 
)
inlineinherited

Subtract a scaled value from this tensor without creating a temporary.

Definition at line 939 of file type_tensor.h.

940 {
941  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
942  _coords[i] -= factor*p._coords[i];
943 }

References libMesh::TypeTensor< T >::_coords.

Referenced by libMesh::HPCoarsenTest::select_refinement().

◆ tr()

template<typename T >
T libMesh::TypeTensor< T >::tr ( ) const
inlineinherited
Returns
The trace of the tensor.

Definition at line 1351 of file type_tensor.h.

1352 {
1353 #if LIBMESH_DIM == 1
1354  return _coords[0];
1355 #endif
1356 
1357 #if LIBMESH_DIM == 2
1358  return _coords[0] + _coords[3];
1359 #endif
1360 
1361 #if LIBMESH_DIM == 3
1362  return _coords[0] + _coords[4] + _coords[8];
1363 #endif
1364 }

◆ transpose()

template<typename T >
TypeTensor< T > libMesh::TypeTensor< T >::transpose ( ) const
inlineinherited
Returns
The transpose of this tensor (with complex numbers not conjugated).

Definition at line 1068 of file type_tensor.h.

1069 {
1070 #if LIBMESH_DIM == 1
1071  return TypeTensor(_coords[0]);
1072 #endif
1073 
1074 #if LIBMESH_DIM == 2
1075  return TypeTensor(_coords[0],
1076  _coords[2],
1077  _coords[1],
1078  _coords[3]);
1079 #endif
1080 
1081 #if LIBMESH_DIM == 3
1082  return TypeTensor(_coords[0],
1083  _coords[3],
1084  _coords[6],
1085  _coords[1],
1086  _coords[4],
1087  _coords[7],
1088  _coords[2],
1089  _coords[5],
1090  _coords[8]);
1091 #endif
1092 }

Referenced by assemble_shell().

◆ write_unformatted()

template<typename T >
void libMesh::TypeTensor< T >::write_unformatted ( std::ostream &  out,
const bool  newline = true 
) const
inherited

Unformatted print to the stream out.

Simply prints the elements of the tensor separated by spaces and newlines.

Definition at line 83 of file type_tensor.C.

85 {
86  libmesh_assert (out_stream);
87 
88  out_stream << std::setiosflags(std::ios::showpoint)
89  << (*this)(0,0) << " "
90  << (*this)(0,1) << " "
91  << (*this)(0,2) << " ";
92  if (newline)
93  out_stream << '\n';
94 
95  out_stream << std::setiosflags(std::ios::showpoint)
96  << (*this)(1,0) << " "
97  << (*this)(1,1) << " "
98  << (*this)(1,2) << " ";
99  if (newline)
100  out_stream << '\n';
101 
102  out_stream << std::setiosflags(std::ios::showpoint)
103  << (*this)(2,0) << " "
104  << (*this)(2,1) << " "
105  << (*this)(2,2) << " ";
106  if (newline)
107  out_stream << '\n';
108 }

References libMesh::libmesh_assert().

◆ zero()

template<typename T >
void libMesh::TypeTensor< T >::zero ( )
inlineinherited

Set all entries of the tensor to 0.

Definition at line 1368 of file type_tensor.h.

1369 {
1370  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1371  _coords[i] = 0.;
1372 }

Referenced by NonlinearNeoHookeCurrentConfig::init_for_qp(), libMesh::TensorValue< T >::operator=(), and libMesh::TypeTensor< T >::operator=().

Member Data Documentation

◆ _coords

template<typename T>
T libMesh::TypeTensor< T >::_coords[LIBMESH_DIM *LIBMESH_DIM]
protectedinherited

The documentation for this class was generated from the following files:
std::sqrt
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::TypeTensor::TypeTensor
TypeTensor()
Empty constructor.
Definition: type_tensor.h:543
libMesh::TypeTensor::add
void add(const TypeTensor< T2 > &)
Add to this tensor without creating a temporary.
Definition: type_tensor.h:858
libMesh::libmesh_assert
libmesh_assert(ctx)
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::TypeTensor::norm_sq
auto norm_sq() const -> decltype(std::norm(T()))
Definition: type_tensor.h:1390
A
static PetscErrorCode Mat * A
Definition: petscdmlibmeshimpl.C:1026
libMesh::CompareTypes::supertype
void supertype
Definition: compare_types.h:100
libMesh::TypeTensor::assign
void assign(const TypeTensor< T2 > &)
Assign to this tensor without creating a temporary.
Definition: type_tensor.h:729
libMesh::Complex
std::complex< Real > Complex
Definition: libmesh_common.h:160
libMesh::TypeTensor::subtract
void subtract(const TypeTensor< T2 > &)
Subtract from this tensor without creating a temporary.
Definition: type_tensor.h:928
libMesh::TensorTools::norm_sq
T norm_sq(std::complex< T > a)
Definition: tensor_tools.h:85
libMesh::TypeTensor::_coords
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477
libMesh::TypeTensor::zero
void zero()
Set all entries of the tensor to 0.
Definition: type_tensor.h:1368
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::TypeTensor::norm
auto norm() const -> decltype(std::norm(T()))
Definition: type_tensor.h:1310
std::imag
boost::multiprecision::float128 imag(const boost::multiprecision::float128)
Definition: float128_shims.h:83
std::real
boost::multiprecision::float128 real(const boost::multiprecision::float128 in)
Definition: float128_shims.h:77