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

This class defines a tensor in LIBMESH_DIM dimensional space of type T. More...

#include <tensor_tools.h>

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

Public Types

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

Public Member Functions

 TypeTensor ()
 Empty constructor. More...
 
template<typename T2 >
 TypeTensor (const TypeTensor< T2 > &p)
 Copy-constructor. More...
 
 ~TypeTensor ()
 Destructor. More...
 
template<typename T2 >
void assign (const TypeTensor< T2 > &)
 Assign to this tensor without creating a temporary. More...
 
template<typename Scalar >
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeTensor & >::type operator= (const Scalar &libmesh_dbg_var(p))
 Assignment-from-scalar operator. 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...
 
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...
 
TypeTensor< T > operator- () const
 
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 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 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 >
TypeTensor< typename CompareTypes< T, T2 >::supertype > operator* (const TypeTensor< T2 > &) const
 Multiply 2 tensors together, i.e. More...
 
template<typename T2 >
const TypeTensor< T > & operator*= (const TypeTensor< T2 > &)
 Multiply this tensor by a tensor value in place. 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 > operator* (const TypeVector< T2 > &) const
 Right-multiply this tensor by a vector, 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
 
bool operator> (const TypeTensor< T > &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 Member Functions

 TypeTensor (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 >
 TypeTensor (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-Scalar. More...
 
template<typename T2 >
 TypeTensor (const TypeVector< T2 > &vx)
 Constructor. More...
 
template<typename T2 >
 TypeTensor (const TypeVector< T2 > &vx, const TypeVector< T2 > &vy)
 
template<typename T2 >
 TypeTensor (const TypeVector< T2 > &vx, const TypeVector< T2 > &vy, const TypeVector< T2 > &vz)
 

Protected Attributes

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

Friends

template<typename T2 >
class TypeTensor
 
std::ostream & operator<< (std::ostream &os, const TypeTensor< T > &t)
 Formatted print as above but supports the syntax: More...
 

Detailed Description

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

This class defines a tensor in LIBMESH_DIM dimensional space of type T.

T may either be Real or Complex.

Author
Roy Stogner
Date
2004

Definition at line 47 of file tensor_tools.h.

Member Typedef Documentation

◆ index_type

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

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

Helper typedef for C++98 generic programming.

Definition at line 139 of file type_tensor.h.

Constructor & Destructor Documentation

◆ TypeTensor() [1/7]

template<typename T >
libMesh::TypeTensor< T >::TypeTensor ( )

Empty constructor.

Gives the tensor 0 in LIBMESH_DIM dimensions.

Definition at line 543 of file type_tensor.h.

544 {
545  _coords[0] = 0;
546 
547 #if LIBMESH_DIM > 1
548  _coords[1] = 0;
549  _coords[2] = 0;
550  _coords[3] = 0;
551 #endif
552 
553 #if LIBMESH_DIM > 2
554  _coords[4] = 0;
555  _coords[5] = 0;
556  _coords[6] = 0;
557  _coords[7] = 0;
558  _coords[8] = 0;
559 #endif
560 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ TypeTensor() [2/7]

template<typename T >
libMesh::TypeTensor< T >::TypeTensor ( 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 
)
explicitprotected

Constructor-from-T.

By default sets higher dimensional entries to 0. This is a poor constructor for 2D tensors - if the default arguments are to be overridden it requires that the "xz = 0." etc. arguments also be given explicitly.

Definition at line 566 of file type_tensor.h.

575 {
576  _coords[0] = xx;
577 
578 #if LIBMESH_DIM == 2
579  _coords[1] = xy;
580  _coords[2] = yx;
581  _coords[3] = yy;
582 #elif LIBMESH_DIM == 1
583  libmesh_assert_equal_to (xy, 0);
584  libmesh_assert_equal_to (yx, 0);
585  libmesh_assert_equal_to (yy, 0);
586 #endif
587 
588 #if LIBMESH_DIM == 3
589  _coords[1] = xy;
590  _coords[2] = xz;
591  _coords[3] = yx;
592  _coords[4] = yy;
593  _coords[5] = yz;
594  _coords[6] = zx;
595  _coords[7] = zy;
596  _coords[8] = zz;
597 #else
598  libmesh_assert_equal_to (xz, 0);
599  libmesh_assert_equal_to (yz, 0);
600  libmesh_assert_equal_to (zx, 0);
601  libmesh_assert_equal_to (zy, 0);
602  libmesh_assert_equal_to (zz, 0);
603 #endif
604 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ TypeTensor() [3/7]

template<typename T >
template<typename Scalar >
libMesh::TypeTensor< T >::TypeTensor ( 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 
)
explicitprotected

Constructor-from-Scalar.

Definition at line 610 of file type_tensor.h.

621 {
622  _coords[0] = xx;
623 
624 #if LIBMESH_DIM == 2
625  _coords[1] = xy;
626  _coords[2] = yx;
627  _coords[3] = yy;
628 #endif
629 
630 #if LIBMESH_DIM == 3
631  _coords[1] = xy;
632  _coords[2] = xz;
633  _coords[3] = yx;
634  _coords[4] = yy;
635  _coords[5] = yz;
636  _coords[6] = zx;
637  _coords[7] = zy;
638  _coords[8] = zz;
639 #endif
640 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ TypeTensor() [4/7]

template<typename T >
template<typename T2 >
libMesh::TypeTensor< T >::TypeTensor ( const TypeVector< T2 > &  vx)
protected

Constructor.

Assigns each vector to a different row of the tensor. We're in LIBMESH_DIM space dimensions and so LIBMESH_DIM many vectors are needed.

Definition at line 657 of file type_tensor.h.

658 {
659  libmesh_assert_equal_to (LIBMESH_DIM, 1);
660  _coords[0] = vx(0);
661 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ TypeTensor() [5/7]

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

Definition at line 665 of file type_tensor.h.

667 {
668  libmesh_assert_equal_to (LIBMESH_DIM, 2);
669  _coords[0] = vx(0);
670  _coords[1] = vx(1);
671  _coords[2] = vy(0);
672  _coords[3] = vy(1);
673 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ TypeTensor() [6/7]

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

Definition at line 677 of file type_tensor.h.

680 {
681  libmesh_assert_equal_to (LIBMESH_DIM, 3);
682  _coords[0] = vx(0);
683  _coords[1] = vx(1);
684  _coords[2] = vx(2);
685  _coords[3] = vy(0);
686  _coords[4] = vy(1);
687  _coords[5] = vy(2);
688  _coords[6] = vz(0);
689  _coords[7] = vz(1);
690  _coords[8] = vz(2);
691 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ TypeTensor() [7/7]

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

Copy-constructor.

Definition at line 647 of file type_tensor.h.

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

648 {
649  // copy the nodes from vector p to me
650  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
651  _coords[i] = p._coords[i];
652 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ ~TypeTensor()

template<typename T >
libMesh::TypeTensor< T >::~TypeTensor ( )

Destructor.

Definition at line 698 of file type_tensor.h.

699 {
700 }

Member Function Documentation

◆ add()

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

Add to this tensor without creating a temporary.

Definition at line 836 of file type_tensor.h.

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

837 {
838  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
839  _coords[i] += p._coords[i];
840 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ add_scaled()

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

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

Definition at line 847 of file type_tensor.h.

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

848 {
849  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
850  _coords[i] += factor*p._coords[i];
851 
852 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ assign()

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

Assign to this tensor without creating a temporary.

Definition at line 707 of file type_tensor.h.

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

708 {
709  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
710  _coords[i] = p._coords[i];
711 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ contract()

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

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 1264 of file type_tensor.h.

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

Referenced by libMesh::TensorTools::inner_product().

1265 {
1266  typename CompareTypes<T,T2>::supertype sum = 0.;
1267  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1268  sum += _coords[i]*t._coords[i];
1269  return sum;
1270 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ det()

template<typename T >
T libMesh::TypeTensor< T >::det ( ) const
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 1306 of file type_tensor.h.

1307 {
1308 #if LIBMESH_DIM == 1
1309  return _coords[0];
1310 #endif
1311 
1312 #if LIBMESH_DIM == 2
1313  return (_coords[0] * _coords[3]
1314  - _coords[1] * _coords[2]);
1315 #endif
1316 
1317 #if LIBMESH_DIM == 3
1318  return (_coords[0] * _coords[4] * _coords[8]
1319  + _coords[1] * _coords[5] * _coords[6]
1320  + _coords[2] * _coords[3] * _coords[7]
1321  - _coords[0] * _coords[5] * _coords[7]
1322  - _coords[1] * _coords[3] * _coords[8]
1323  - _coords[2] * _coords[4] * _coords[6]);
1324 #endif
1325 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ inverse()

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

Definition at line 1076 of file type_tensor.h.

1077 {
1078 #if LIBMESH_DIM == 1
1079  if (_coords[0] == static_cast<T>(0.))
1080  libmesh_convergence_failure();
1081  return TypeTensor(1. / _coords[0]);
1082 #endif
1083 
1084 #if LIBMESH_DIM == 2
1085  // Get convenient reference to this.
1086  const TypeTensor<T> & A = *this;
1087 
1088  // Use temporary variables, avoid multiple accesses
1089  T a = A(0,0), b = A(0,1),
1090  c = A(1,0), d = A(1,1);
1091 
1092  // Make sure det = ad - bc is not zero
1093  T my_det = a*d - b*c;
1094 
1095  if (my_det == static_cast<T>(0.))
1096  libmesh_convergence_failure();
1097 
1098  return TypeTensor(d/my_det, -b/my_det, -c/my_det, a/my_det);
1099 #endif
1100 
1101 #if LIBMESH_DIM == 3
1102  // Get convenient reference to this.
1103  const TypeTensor<T> & A = *this;
1104 
1105  T a11 = A(0,0), a12 = A(0,1), a13 = A(0,2),
1106  a21 = A(1,0), a22 = A(1,1), a23 = A(1,2),
1107  a31 = A(2,0), a32 = A(2,1), a33 = A(2,2);
1108 
1109  T my_det = a11*(a33*a22-a32*a23) - a21*(a33*a12-a32*a13) + a31*(a23*a12-a22*a13);
1110 
1111  if (my_det == static_cast<T>(0.))
1112  libmesh_convergence_failure();
1113 
1114  // Inline comment characters are for lining up columns.
1115  return TypeTensor( (a33*a22-a32*a23)/my_det, -(a33*a12-a32*a13)/my_det, (a23*a12-a22*a13)/my_det,
1116  -(a33*a21-a31*a23)/my_det, (a33*a11-a31*a13)/my_det, -(a23*a11-a21*a13)/my_det,
1117  (a32*a21-a31*a22)/my_det, -(a32*a11-a31*a12)/my_det, (a22*a11-a21*a12)/my_det);
1118 #endif
1119 }
TypeTensor()
Empty constructor.
Definition: type_tensor.h:543
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ is_zero()

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

Definition at line 1296 of file type_tensor.h.

1297 {
1298  for (const auto & val : _coords)
1299  if (val != T(0))
1300  return false;
1301  return true;
1302 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ left_multiply()

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

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 1207 of file type_tensor.h.

Referenced by libMesh::operator*().

1208 {
1209  TypeVector<typename CompareTypes<T,T2>::supertype> returnval;
1210  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1211  for (unsigned int j=0; j<LIBMESH_DIM; j++)
1212  returnval(i) += p(j)*(*this)(j,i);
1213 
1214  return returnval;
1215 }

◆ norm()

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

Definition at line 1288 of file type_tensor.h.

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

1289 {
1290  return std::sqrt(this->norm_sq());
1291 }
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
auto norm_sq() const -> decltype(std::norm(T()))
Definition: type_tensor.h:1368

◆ norm_sq()

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

Definition at line 1368 of file type_tensor.h.

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

1369 {
1370  Real sum = 0.;
1371  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1372  sum += TensorTools::norm_sq(_coords[i]);
1373  return sum;
1374 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477
T norm_sq(std::complex< T > a)
Definition: tensor_tools.h:85
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ operator()() [1/2]

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

Definition at line 717 of file type_tensor.h.

719 {
720  libmesh_assert_less (i, 3);
721  libmesh_assert_less (j, 3);
722 
723 #if LIBMESH_DIM < 3
724  const static T my_zero = 0;
725  if (i >= LIBMESH_DIM || j >= LIBMESH_DIM)
726  return my_zero;
727 #endif
728 
729  return _coords[i*LIBMESH_DIM+j];
730 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ operator()() [2/2]

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

Definition at line 736 of file type_tensor.h.

738 {
739 #if LIBMESH_DIM < 3
740 
741  if (i >= LIBMESH_DIM || j >= LIBMESH_DIM)
742  libmesh_error_msg("ERROR: You are assigning to a tensor component that is out of range for the compiled LIBMESH_DIM!");
743 
744 #endif
745 
746  libmesh_assert_less (i, LIBMESH_DIM);
747  libmesh_assert_less (j, LIBMESH_DIM);
748 
749  return _coords[i*LIBMESH_DIM+j];
750 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ 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

Multiply this tensor by a scalar value.

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

Definition at line 961 of file type_tensor.h.

964 {
965  typedef decltype((*this)(0, 0) * factor) TS;
966 
967 
968 #if LIBMESH_DIM == 1
969  return TypeTensor<TS>(_coords[0]*factor);
970 #endif
971 
972 #if LIBMESH_DIM == 2
973  return TypeTensor<TS>(_coords[0]*factor,
974  _coords[1]*factor,
975  _coords[2]*factor,
976  _coords[3]*factor);
977 #endif
978 
979 #if LIBMESH_DIM == 3
980  return TypeTensor<TS>(_coords[0]*factor,
981  _coords[1]*factor,
982  _coords[2]*factor,
983  _coords[3]*factor,
984  _coords[4]*factor,
985  _coords[5]*factor,
986  _coords[6]*factor,
987  _coords[7]*factor,
988  _coords[8]*factor);
989 #endif
990 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ operator*() [2/3]

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

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 1229 of file type_tensor.h.

1230 {
1231  TypeTensor<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  for (unsigned int k=0; k<LIBMESH_DIM; k++)
1235  returnval(i,j) += (*this)(i,k)*p(k,j);
1236 
1237  return returnval;
1238 }

◆ operator*() [3/3]

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

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 1193 of file type_tensor.h.

1194 {
1195  TypeVector<typename CompareTypes<T,T2>::supertype> returnval;
1196  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1197  for (unsigned int j=0; j<LIBMESH_DIM; j++)
1198  returnval(i) += (*this)(i,j)*p(j);
1199 
1200  return returnval;
1201 }

◆ 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)

Multiply this tensor by a scalar value in place.

Returns
A reference to *this.

Definition at line 282 of file type_tensor.h.

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

283  {
284  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
285  _coords[i] *= factor;
286 
287  return *this;
288  }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ operator*=() [2/2]

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

Multiply this tensor by a tensor value in place.

Returns
A reference to *this

Definition at line 1243 of file type_tensor.h.

1244 {
1245  TypeTensor<T> temp;
1246  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1247  for (unsigned int j=0; j<LIBMESH_DIM; j++)
1248  for (unsigned int k=0; k<LIBMESH_DIM; k++)
1249  temp(i,j) += (*this)(i,k)*p(k,j);
1250 
1251  this->assign(temp);
1252  return *this;
1253 }
void assign(const TypeTensor< T2 > &)
Assign to this tensor without creating a temporary.
Definition: type_tensor.h:707

◆ operator+()

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

Add another tensor to this tensor.

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

Definition at line 790 of file type_tensor.h.

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

791 {
792 
793 #if LIBMESH_DIM == 1
794  return TypeTensor<typename CompareTypes<T, T2>::supertype>(_coords[0] + p._coords[0]);
795 #endif
796 
797 #if LIBMESH_DIM == 2
798  return TypeTensor<typename CompareTypes<T, T2>::supertype>(_coords[0] + p._coords[0],
799  _coords[1] + p._coords[1],
800  0.,
801  _coords[2] + p._coords[2],
802  _coords[3] + p._coords[3]);
803 #endif
804 
805 #if LIBMESH_DIM == 3
806  return TypeTensor<typename CompareTypes<T, T2>::supertype>(_coords[0] + p._coords[0],
807  _coords[1] + p._coords[1],
808  _coords[2] + p._coords[2],
809  _coords[3] + p._coords[3],
810  _coords[4] + p._coords[4],
811  _coords[5] + p._coords[5],
812  _coords[6] + p._coords[6],
813  _coords[7] + p._coords[7],
814  _coords[8] + p._coords[8]);
815 #endif
816 
817 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ operator+=()

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

Add to this tensor.

Returns
A reference to *this.

Definition at line 824 of file type_tensor.h.

825 {
826  this->add (p);
827 
828  return *this;
829 }
void add(const TypeTensor< T2 > &)
Add to this tensor without creating a temporary.
Definition: type_tensor.h:836

◆ operator-() [1/2]

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

Subtract a tensor from this tensor.

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

Definition at line 860 of file type_tensor.h.

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

861 {
862 
863 #if LIBMESH_DIM == 1
864  return TypeTensor<typename CompareTypes<T, T2>::supertype>(_coords[0] - p._coords[0]);
865 #endif
866 
867 #if LIBMESH_DIM == 2
868  return TypeTensor<typename CompareTypes<T, T2>::supertype>(_coords[0] - p._coords[0],
869  _coords[1] - p._coords[1],
870  0.,
871  _coords[2] - p._coords[2],
872  _coords[3] - p._coords[3]);
873 #endif
874 
875 #if LIBMESH_DIM == 3
876  return TypeTensor<typename CompareTypes<T, T2>::supertype>(_coords[0] - p._coords[0],
877  _coords[1] - p._coords[1],
878  _coords[2] - p._coords[2],
879  _coords[3] - p._coords[3],
880  _coords[4] - p._coords[4],
881  _coords[5] - p._coords[5],
882  _coords[6] - p._coords[6],
883  _coords[7] - p._coords[7],
884  _coords[8] - p._coords[8]);
885 #endif
886 
887 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ operator-() [2/2]

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

Definition at line 927 of file type_tensor.h.

928 {
929 
930 #if LIBMESH_DIM == 1
931  return TypeTensor(-_coords[0]);
932 #endif
933 
934 #if LIBMESH_DIM == 2
935  return TypeTensor(-_coords[0],
936  -_coords[1],
937  -_coords[2],
938  -_coords[3]);
939 #endif
940 
941 #if LIBMESH_DIM == 3
942  return TypeTensor(-_coords[0],
943  -_coords[1],
944  -_coords[2],
945  -_coords[3],
946  -_coords[4],
947  -_coords[5],
948  -_coords[6],
949  -_coords[7],
950  -_coords[8]);
951 #endif
952 
953 }
TypeTensor()
Empty constructor.
Definition: type_tensor.h:543
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ operator-=()

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

Subtract from this tensor.

Returns
A reference to *this.

Definition at line 894 of file type_tensor.h.

895 {
896  this->subtract (p);
897 
898  return *this;
899 }
void subtract(const TypeTensor< T2 > &)
Subtract from this tensor without creating a temporary.
Definition: type_tensor.h:906

◆ 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

Divide each entry of this tensor by a scalar value.

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

Definition at line 1011 of file type_tensor.h.

1012 {
1013  libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
1014 
1015  typedef typename CompareTypes<T, Scalar>::supertype TS;
1016 
1017 #if LIBMESH_DIM == 1
1018  return TypeTensor<TS>(_coords[0]/factor);
1019 #endif
1020 
1021 #if LIBMESH_DIM == 2
1022  return TypeTensor<TS>(_coords[0]/factor,
1023  _coords[1]/factor,
1024  _coords[2]/factor,
1025  _coords[3]/factor);
1026 #endif
1027 
1028 #if LIBMESH_DIM == 3
1029  return TypeTensor<TS>(_coords[0]/factor,
1030  _coords[1]/factor,
1031  _coords[2]/factor,
1032  _coords[3]/factor,
1033  _coords[4]/factor,
1034  _coords[5]/factor,
1035  _coords[6]/factor,
1036  _coords[7]/factor,
1037  _coords[8]/factor);
1038 #endif
1039 
1040 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ operator/=()

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

Divide each entry of this tensor by a scalar value.

Returns
A reference to *this.

Definition at line 1176 of file type_tensor.h.

1177 {
1178  libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
1179 
1180  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1181  _coords[i] /= factor;
1182 
1183  return *this;
1184 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ operator<()

template<typename T>
bool libMesh::TypeTensor< T >::operator< ( const TypeTensor< T > &  rhs) const
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, TypeTensor &>::type libMesh::TypeTensor< T >::operator= ( const Scalar &  libmesh_dbg_varp)

Assignment-from-scalar operator.

Used only to zero out vectors.

Returns
A reference to *this.

Definition at line 172 of file type_tensor.h.

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

173  { libmesh_assert_equal_to (p, Scalar(0)); this->zero(); return *this; }
void zero()
Set all entries of the tensor to 0.
Definition: type_tensor.h:1346

◆ operator==()

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

Definition at line 1380 of file type_tensor.h.

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

1381 {
1382 #if LIBMESH_DIM == 1
1383  return (std::abs(_coords[0] - rhs._coords[0])
1384  < TOLERANCE);
1385 #endif
1386 
1387 #if LIBMESH_DIM == 2
1388  return ((std::abs(_coords[0] - rhs._coords[0]) +
1389  std::abs(_coords[1] - rhs._coords[1]) +
1390  std::abs(_coords[2] - rhs._coords[2]) +
1391  std::abs(_coords[3] - rhs._coords[3]))
1392  < 4.*TOLERANCE);
1393 #endif
1394 
1395 #if LIBMESH_DIM == 3
1396  return ((std::abs(_coords[0] - rhs._coords[0]) +
1397  std::abs(_coords[1] - rhs._coords[1]) +
1398  std::abs(_coords[2] - rhs._coords[2]) +
1399  std::abs(_coords[3] - rhs._coords[3]) +
1400  std::abs(_coords[4] - rhs._coords[4]) +
1401  std::abs(_coords[5] - rhs._coords[5]) +
1402  std::abs(_coords[6] - rhs._coords[6]) +
1403  std::abs(_coords[7] - rhs._coords[7]) +
1404  std::abs(_coords[8] - rhs._coords[8]))
1405  < 9.*TOLERANCE);
1406 #endif
1407 
1408 }
double abs(double a)
static const Real TOLERANCE
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ operator>()

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

◆ print()

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

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

◆ row()

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

Definition at line 776 of file type_tensor.h.

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

777 {
778  TypeVector<T> return_vector;
779 
780  for (unsigned int j=0; j<LIBMESH_DIM; j++)
781  return_vector._coords[j] = _coords[r*LIBMESH_DIM + j];
782 
783  return return_vector;
784 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ size()

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

Definition at line 1277 of file type_tensor.h.

References std::norm().

1278 {
1279  libmesh_deprecated();
1280  return this->norm();
1281 }
auto norm() const -> decltype(std::norm(T()))
Definition: type_tensor.h:1288

◆ size_sq()

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

Definition at line 1357 of file type_tensor.h.

References libMesh::TensorTools::norm_sq().

1358 {
1359  libmesh_deprecated();
1360  return this->norm_sq();
1361 }
auto norm_sq() const -> decltype(std::norm(T()))
Definition: type_tensor.h:1368

◆ slice() [1/2]

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

Definition at line 756 of file type_tensor.h.

757 {
758  libmesh_assert_less (i, LIBMESH_DIM);
759  return ConstTypeTensorColumn<T>(*this, i);
760 }

◆ slice() [2/2]

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

Definition at line 766 of file type_tensor.h.

767 {
768  libmesh_assert_less (i, LIBMESH_DIM);
769  return TypeTensorColumn<T>(*this, i);
770 }

◆ solve()

template<typename T >
void libMesh::TypeTensor< T >::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.

Returns
The solution in the x vector.

Definition at line 1125 of file type_tensor.h.

1126 {
1127 #if LIBMESH_DIM == 1
1128  if (_coords[0] == static_cast<T>(0.))
1129  libmesh_convergence_failure();
1130  x(0) = b(0) / _coords[0];
1131 #endif
1132 
1133 #if LIBMESH_DIM == 2
1134  T my_det = _coords[0]*_coords[3] - _coords[1]*_coords[2];
1135 
1136  if (my_det == static_cast<T>(0.))
1137  libmesh_convergence_failure();
1138 
1139  T my_det_inv = 1./my_det;
1140 
1141  x(0) = my_det_inv*( _coords[3]*b(0) - _coords[1]*b(1));
1142  x(1) = my_det_inv*(-_coords[2]*b(0) + _coords[0]*b(1));
1143 #endif
1144 
1145 #if LIBMESH_DIM == 3
1146  T my_det =
1147  // a11*(a33 *a22 - a32 *a23)
1148  _coords[0]*(_coords[8]*_coords[4] - _coords[7]*_coords[5])
1149  // -a21*(a33 *a12 - a32 *a13)
1150  -_coords[3]*(_coords[8]*_coords[1] - _coords[7]*_coords[2]) +
1151  // +a31*(a23 *a12 - a22 *a13)
1152  +_coords[6]*(_coords[5]*_coords[1] - _coords[4]*_coords[2]);
1153 
1154  if (my_det == static_cast<T>(0.))
1155  libmesh_convergence_failure();
1156 
1157  T my_det_inv = 1./my_det;
1158  x(0) = my_det_inv*((_coords[8]*_coords[4] - _coords[7]*_coords[5])*b(0) -
1159  (_coords[8]*_coords[1] - _coords[7]*_coords[2])*b(1) +
1160  (_coords[5]*_coords[1] - _coords[4]*_coords[2])*b(2));
1161 
1162  x(1) = my_det_inv*((_coords[6]*_coords[5] - _coords[8]*_coords[3])*b(0) +
1163  (_coords[8]*_coords[0] - _coords[6]*_coords[2])*b(1) -
1164  (_coords[5]*_coords[0] - _coords[3]*_coords[2])*b(2));
1165 
1166  x(2) = my_det_inv*((_coords[7]*_coords[3] - _coords[6]*_coords[4])*b(0) -
1167  (_coords[7]*_coords[0] - _coords[6]*_coords[1])*b(1) +
1168  (_coords[4]*_coords[0] - _coords[3]*_coords[1])*b(2));
1169 #endif
1170 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ subtract()

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

Subtract from this tensor without creating a temporary.

Definition at line 906 of file type_tensor.h.

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

907 {
908  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
909  _coords[i] -= p._coords[i];
910 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ subtract_scaled()

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

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

Definition at line 917 of file type_tensor.h.

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

918 {
919  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
920  _coords[i] -= factor*p._coords[i];
921 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ tr()

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

Definition at line 1329 of file type_tensor.h.

1330 {
1331 #if LIBMESH_DIM == 1
1332  return _coords[0];
1333 #endif
1334 
1335 #if LIBMESH_DIM == 2
1336  return _coords[0] + _coords[3];
1337 #endif
1338 
1339 #if LIBMESH_DIM == 3
1340  return _coords[0] + _coords[4] + _coords[8];
1341 #endif
1342 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ transpose()

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

Definition at line 1046 of file type_tensor.h.

1047 {
1048 #if LIBMESH_DIM == 1
1049  return TypeTensor(_coords[0]);
1050 #endif
1051 
1052 #if LIBMESH_DIM == 2
1053  return TypeTensor(_coords[0],
1054  _coords[2],
1055  _coords[1],
1056  _coords[3]);
1057 #endif
1058 
1059 #if LIBMESH_DIM == 3
1060  return TypeTensor(_coords[0],
1061  _coords[3],
1062  _coords[6],
1063  _coords[1],
1064  _coords[4],
1065  _coords[7],
1066  _coords[2],
1067  _coords[5],
1068  _coords[8]);
1069 #endif
1070 }
TypeTensor()
Empty constructor.
Definition: type_tensor.h:543
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

◆ write_unformatted()

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

Unformatted print to the stream out.

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

◆ zero()

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

Set all entries of the tensor to 0.

Definition at line 1346 of file type_tensor.h.

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

1347 {
1348  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1349  _coords[i] = 0.;
1350 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:477

Friends And Related Function Documentation

◆ operator<<

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

Formatted print as above but supports the syntax:

std::cout << t << std::endl;

Definition at line 460 of file type_tensor.h.

461  {
462  t.print(os);
463  return os;
464  }

◆ TypeTensor

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

Definition at line 74 of file type_tensor.h.

Member Data Documentation

◆ _coords

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

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