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()))
 
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 538 of file type_tensor.h.

539 {
540  _coords[0] = 0;
541 
542 #if LIBMESH_DIM > 1
543  _coords[1] = 0;
544  _coords[2] = 0;
545  _coords[3] = 0;
546 #endif
547 
548 #if LIBMESH_DIM > 2
549  _coords[4] = 0;
550  _coords[5] = 0;
551  _coords[6] = 0;
552  _coords[7] = 0;
553  _coords[8] = 0;
554 #endif
555 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

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

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

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

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

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

653 {
654  libmesh_assert_equal_to (LIBMESH_DIM, 1);
655  _coords[0] = vx(0);
656 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

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

662 {
663  libmesh_assert_equal_to (LIBMESH_DIM, 2);
664  _coords[0] = vx(0);
665  _coords[1] = vx(1);
666  _coords[2] = vy(0);
667  _coords[3] = vy(1);
668 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

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

675 {
676  libmesh_assert_equal_to (LIBMESH_DIM, 3);
677  _coords[0] = vx(0);
678  _coords[1] = vx(1);
679  _coords[2] = vx(2);
680  _coords[3] = vy(0);
681  _coords[4] = vy(1);
682  _coords[5] = vy(2);
683  _coords[6] = vz(0);
684  _coords[7] = vz(1);
685  _coords[8] = vz(2);
686 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

◆ TypeTensor() [7/7]

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

Copy-constructor.

Definition at line 642 of file type_tensor.h.

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

643 {
644  // copy the nodes from vector p to me
645  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
646  _coords[i] = p._coords[i];
647 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

◆ ~TypeTensor()

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

Destructor.

Definition at line 693 of file type_tensor.h.

694 {
695 }

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

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

832 {
833  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
834  _coords[i] += p._coords[i];
835 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

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

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

843 {
844  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
845  _coords[i] += factor*p._coords[i];
846 
847 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

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

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

703 {
704  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
705  _coords[i] = p._coords[i];
706 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

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

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

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

1260 {
1261  typename CompareTypes<T,T2>::supertype sum = 0.;
1262  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1263  sum += _coords[i]*t._coords[i];
1264  return sum;
1265 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

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

1293 {
1294 #if LIBMESH_DIM == 1
1295  return _coords[0];
1296 #endif
1297 
1298 #if LIBMESH_DIM == 2
1299  return (_coords[0] * _coords[3]
1300  - _coords[1] * _coords[2]);
1301 #endif
1302 
1303 #if LIBMESH_DIM == 3
1304  return (_coords[0] * _coords[4] * _coords[8]
1305  + _coords[1] * _coords[5] * _coords[6]
1306  + _coords[2] * _coords[3] * _coords[7]
1307  - _coords[0] * _coords[5] * _coords[7]
1308  - _coords[1] * _coords[3] * _coords[8]
1309  - _coords[2] * _coords[4] * _coords[6]);
1310 #endif
1311 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

◆ inverse()

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

Definition at line 1071 of file type_tensor.h.

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

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

Referenced by libMesh::operator*().

1203 {
1204  TypeVector<typename CompareTypes<T,T2>::supertype> returnval;
1205  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1206  for (unsigned int j=0; j<LIBMESH_DIM; j++)
1207  returnval(i) += p(i)*(*this)(i,j);
1208 
1209  return returnval;
1210 }

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

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

1284 {
1285  return std::sqrt(this->norm_sq());
1286 }
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
auto norm_sq() const -> decltype(std::norm(T()))
Definition: type_tensor.h:1354

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

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

1355 {
1356  Real sum = 0.;
1357  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1358  sum += TensorTools::norm_sq(_coords[i]);
1359  return sum;
1360 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472
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 712 of file type_tensor.h.

714 {
715  libmesh_assert_less (i, 3);
716  libmesh_assert_less (j, 3);
717 
718 #if LIBMESH_DIM < 3
719  const static T my_zero = 0;
720  if (i >= LIBMESH_DIM || j >= LIBMESH_DIM)
721  return my_zero;
722 #endif
723 
724  return _coords[i*LIBMESH_DIM+j];
725 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

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

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

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

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

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

1225 {
1226  TypeTensor<typename CompareTypes<T, T2>::supertype> returnval;
1227  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1228  for (unsigned int j=0; j<LIBMESH_DIM; j++)
1229  for (unsigned int k=0; k<LIBMESH_DIM; k++)
1230  returnval(i,j) += (*this)(i,k)*p(k,j);
1231 
1232  return returnval;
1233 }

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

1189 {
1190  TypeVector<typename CompareTypes<T,T2>::supertype> returnval;
1191  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1192  for (unsigned int j=0; j<LIBMESH_DIM; j++)
1193  returnval(i) += (*this)(i,j)*p(j);
1194 
1195  return returnval;
1196 }

◆ 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:472

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

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

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

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

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

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

820 {
821  this->add (p);
822 
823  return *this;
824 }
void add(const TypeTensor< T2 > &)
Add to this tensor without creating a temporary.
Definition: type_tensor.h:831

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

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

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

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

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

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

890 {
891  this->subtract (p);
892 
893  return *this;
894 }
void subtract(const TypeTensor< T2 > &)
Subtract from this tensor without creating a temporary.
Definition: type_tensor.h:901

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

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

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

1172 {
1173  libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
1174 
1175  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1176  _coords[i] /= factor;
1177 
1178  return *this;
1179 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

◆ 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:1332

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

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

1367 {
1368 #if LIBMESH_DIM == 1
1369  return (std::abs(_coords[0] - rhs._coords[0])
1370  < TOLERANCE);
1371 #endif
1372 
1373 #if LIBMESH_DIM == 2
1374  return ((std::abs(_coords[0] - rhs._coords[0]) +
1375  std::abs(_coords[1] - rhs._coords[1]) +
1376  std::abs(_coords[2] - rhs._coords[2]) +
1377  std::abs(_coords[3] - rhs._coords[3]))
1378  < 4.*TOLERANCE);
1379 #endif
1380 
1381 #if LIBMESH_DIM == 3
1382  return ((std::abs(_coords[0] - rhs._coords[0]) +
1383  std::abs(_coords[1] - rhs._coords[1]) +
1384  std::abs(_coords[2] - rhs._coords[2]) +
1385  std::abs(_coords[3] - rhs._coords[3]) +
1386  std::abs(_coords[4] - rhs._coords[4]) +
1387  std::abs(_coords[5] - rhs._coords[5]) +
1388  std::abs(_coords[6] - rhs._coords[6]) +
1389  std::abs(_coords[7] - rhs._coords[7]) +
1390  std::abs(_coords[8] - rhs._coords[8]))
1391  < 9.*TOLERANCE);
1392 #endif
1393 
1394 }
double abs(double a)
static const Real TOLERANCE
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

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

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

772 {
773  TypeVector<T> return_vector;
774 
775  for (unsigned int j=0; j<LIBMESH_DIM; j++)
776  return_vector._coords[j] = _coords[r*LIBMESH_DIM + j];
777 
778  return return_vector;
779 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

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

References std::norm().

1273 {
1274  libmesh_deprecated();
1275  return this->norm();
1276 }
auto norm() const -> decltype(std::norm(T()))
Definition: type_tensor.h:1283

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

References libMesh::TensorTools::norm_sq().

1344 {
1345  libmesh_deprecated();
1346  return this->norm_sq();
1347 }
auto norm_sq() const -> decltype(std::norm(T()))
Definition: type_tensor.h:1354

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

752 {
753  libmesh_assert_less (i, LIBMESH_DIM);
754  return ConstTypeTensorColumn<T>(*this, i);
755 }

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

762 {
763  libmesh_assert_less (i, LIBMESH_DIM);
764  return TypeTensorColumn<T>(*this, i);
765 }

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

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

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

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

902 {
903  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
904  _coords[i] -= p._coords[i];
905 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

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

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

913 {
914  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
915  _coords[i] -= factor*p._coords[i];
916 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

◆ tr()

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

Definition at line 1315 of file type_tensor.h.

1316 {
1317 #if LIBMESH_DIM == 1
1318  return _coords[0];
1319 #endif
1320 
1321 #if LIBMESH_DIM == 2
1322  return _coords[0] + _coords[3];
1323 #endif
1324 
1325 #if LIBMESH_DIM == 3
1326  return _coords[0] + _coords[4] + _coords[8];
1327 #endif
1328 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

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

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

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

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

1333 {
1334  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1335  _coords[i] = 0.;
1336 }
T _coords[LIBMESH_DIM *LIBMESH_DIM]
The coordinates of the TypeTensor.
Definition: type_tensor.h:472

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

456  {
457  t.print(os);
458  return os;
459  }

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