Intrepid
Classes | Static Public Member Functions | List of all members
Intrepid::RealSpaceTools< Scalar > Class Template Reference

Implementation of basic linear algebra functionality in Euclidean space. More...

#include <Intrepid_RealSpaceTools.hpp>

Classes

struct  detTempSpec
 

Static Public Member Functions

static void absval (Scalar *absArray, const Scalar *inArray, const int size)
 Computes absolute value of contiguous input data inArray of size size. More...
 
static void absval (Scalar *inoutArray, const int size)
 Computes absolute value of contiguous data inoutAbsArray of size size in place. More...
 
template<class ArrayAbs , class ArrayIn >
static void absval (ArrayAbs &absArray, const ArrayIn &inArray)
 Computes absolute value of an array. More...
 
template<class ArrayInOut >
static void absval (ArrayInOut &inoutAbsArray)
 Computes, in place, absolute value of an array. More...
 
static Scalar vectorNorm (const Scalar *inVec, const size_t dim, const ENorm normType)
 Computes norm (1, 2, infinity) of the vector inVec of size dim. More...
 
template<class ArrayIn >
static Scalar vectorNorm (const ArrayIn &inVec, const ENorm normType)
 Computes norm (1, 2, infinity) of a single vector stored in an array of rank 1. More...
 
template<class ArrayNorm , class ArrayIn >
static void vectorNorm (ArrayNorm &normArray, const ArrayIn &inVecs, const ENorm normType)
 Computes norms (1, 2, infinity) of vectors stored in a array of total rank 2 (array of vectors), indexed by (i0, D), or 3 (array of arrays of vectors), indexed by (i0, i1, D). More...
 
static void transpose (Scalar *transposeMat, const Scalar *inMat, const size_t dim)
 Computes transpose of the square matrix inMat of size dim by dim. More...
 
template<class ArrayTranspose , class ArrayIn >
static void transpose (ArrayTranspose &transposeMats, const ArrayIn &inMats)
 Computes transposes of square matrices stored in an array of total rank 2 (single matrix), indexed by (D, D), 3 (array of matrices), indexed by (i0, D, D), or 4 (array of arrays of matrices), indexed by (i0, i1, D, D). More...
 
static void inverse (Scalar *inverseMat, const Scalar *inMat, const size_t dim)
 Computes inverse of the square matrix inMat of size dim by dim. More...
 
template<class ArrayInverse , class ArrayIn >
static void inverse (ArrayInverse &inverseMats, const ArrayIn &inMats)
 Computes inverses of nonsingular matrices stored in an array of total rank 2 (single matrix), indexed by (D, D), 3 (array of matrices), indexed by (i0, D, D), or 4 (array of arrays of matrices), indexed by (i0, i1, D, D). More...
 
static Scalar det (const Scalar *inMat, const size_t dim)
 
template<class ArrayIn >
static Scalar det (const ArrayIn &inMat)
 Computes determinant of a single square matrix stored in an array of rank 2. More...
 
template<class ArrayDet , class ArrayIn >
static void det (ArrayDet &detArray, const ArrayIn &inMats)
 Computes determinants of matrices stored in an array of total rank 3 (array of matrices), indexed by (i0, D, D), or 4 (array of arrays of matrices), indexed by (i0, i1, D, D). More...
 
static void add (Scalar *sumArray, const Scalar *inArray1, const Scalar *inArray2, const int size)
 Adds contiguous data inArray1 and inArray2 of size size:
sumArray = inArray1 + inArray2. More...
 
static void add (Scalar *inoutSumArray, const Scalar *inArray, const int size)
 Adds, in place, contiguous data inArray into inoutSumArray of size size:
inoutSumArray = inoutSumArray + inArray. More...
 
template<class ArraySum , class ArrayIn1 , class ArrayIn2 >
static void add (ArraySum &sumArray, const ArrayIn1 &inArray1, const ArrayIn2 &inArray2)
 Adds arrays inArray1 and inArray2:
sumArray = inArray1 + inArray2. More...
 
template<class ArraySum , class ArrayIn >
static void add (ArraySum &inoutSumArray, const ArrayIn &inArray)
 Adds, in place, inArray into inoutSumArray:
inoutSumArray = inoutSumArray + inArray. More...
 
static void subtract (Scalar *diffArray, const Scalar *inArray1, const Scalar *inArray2, const int size)
 Subtracts contiguous data inArray2 from inArray1 of size size:
diffArray = inArray1 - inArray2. More...
 
static void subtract (Scalar *inoutDiffArray, const Scalar *inArray, const int size)
 Subtracts, in place, contiguous data inArray from inoutDiffArray of size size:
inoutDiffArray = inoutDiffArray - inArray. More...
 
template<class ArrayDiff , class ArrayIn1 , class ArrayIn2 >
static void subtract (ArrayDiff &diffArray, const ArrayIn1 &inArray1, const ArrayIn2 &inArray2)
 Subtracts inArray2 from inArray1:
diffArray = inArray1 - inArray2. More...
 
template<class ArrayDiff , class ArrayIn >
static void subtract (ArrayDiff &inoutDiffArray, const ArrayIn &inArray)
 Subtracts, in place, inArray from inoutDiffArray:
inoutDiffArray = inoutDiffArray - inArray. More...
 
template<class ArrayDiff , class ArrayIn >
static void subtractTemp (ArrayDiff &inoutDiffArray, const ArrayIn &inArray)
 
static void scale (Scalar *scaledArray, const Scalar *inArray, const int size, const Scalar scalar)
 Multiplies contiguous data inArray of size size by a scalar (componentwise):
scaledArray = scalar * inArray. More...
 
static void scale (Scalar *inoutScaledArray, const int size, const Scalar scalar)
 Multiplies, in place, contiguous data inoutScaledArray of size size by a scalar (componentwise):
inoutScaledArray = scalar * inoutScaledArray. More...
 
template<class ArrayScaled , class ArrayIn >
static void scale (ArrayScaled &scaledArray, const ArrayIn &inArray, const Scalar scalar)
 Multiplies array inArray by the scalar scalar (componentwise):
scaledArray = scalar * inArray. More...
 
template<class ArrayScaled >
static void scale (ArrayScaled &inoutScaledArray, const Scalar scalar)
 Multiplies, in place, array inoutScaledArray by the scalar scalar (componentwise):
inoutScaledArray = scalar * inoutScaledArray. More...
 
static Scalar dot (const Scalar *inArray1, const Scalar *inArray2, const int size)
 Computes dot product of contiguous data inArray1 and inArray2 of size size. More...
 
template<class ArrayVec1 , class ArrayVec2 >
static Scalar dot (const ArrayVec1 &inVec1, const ArrayVec2 &inVec2)
 Computes dot product of two vectors stored in arrays of rank 1. More...
 
template<class ArrayDot , class ArrayVec1 , class ArrayVec2 >
static void dot (ArrayDot &dotArray, const ArrayVec1 &inVecs1, const ArrayVec2 &inVecs2)
 Computes dot product of vectors stored in an array of total rank 2 (array of vectors), indexed by (i0, D), or 3 (array of arrays of vectors), indexed by (i0, i1, D). More...
 
static void matvec (Scalar *matVec, const Scalar *inMat, const Scalar *inVec, const size_t dim)
 Matrix-vector left multiply using contiguous data:
matVec = inMat * inVec. More...
 
template<class ArrayMatVec , class ArrayMat , class ArrayVec >
static void matvec (ArrayMatVec &matVecs, const ArrayMat &inMats, const ArrayVec &inVecs)
 Matrix-vector left multiply using multidimensional arrays:
matVec = inMat * inVec. More...
 
template<class ArrayVecProd , class ArrayIn1 , class ArrayIn2 >
static void vecprod (ArrayVecProd &vecProd, const ArrayIn1 &inLeft, const ArrayIn2 &inRight)
 Vector product using multidimensional arrays:
vecProd = inVecLeft x inVecRight More...
 

Detailed Description

template<class Scalar>
class Intrepid::RealSpaceTools< Scalar >

Implementation of basic linear algebra functionality in Euclidean space.

Definition at line 65 of file Intrepid_RealSpaceTools.hpp.

Member Function Documentation

template<class Scalar >
void Intrepid::RealSpaceTools< Scalar >::absval ( Scalar *  absArray,
const Scalar *  inArray,
const int  size 
)
static

Computes absolute value of contiguous input data inArray of size size.

Parameters
absArray[out] - output data
inArray[in] - input data
size[in] - size

Definition at line 55 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
void Intrepid::RealSpaceTools< Scalar >::absval ( Scalar *  inoutArray,
const int  size 
)
static

Computes absolute value of contiguous data inoutAbsArray of size size in place.

Parameters
inoutAbsArray[in/out] - input/output data
size[in] - size

Definition at line 64 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArrayAbs , class ArrayIn >
void Intrepid::RealSpaceTools< Scalar >::absval ( ArrayAbs &  absArray,
const ArrayIn &  inArray 
)
static

Computes absolute value of an array.

Parameters
outArray[out] - output array
inArray[in] - input array
Note
Requirements (checked at runtime, in debug mode):
  • rank(absArray) == rank(inArray)
  • dimensions(absArray) == dimensions(inArray)

Definition at line 74 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArrayInOut >
void Intrepid::RealSpaceTools< Scalar >::absval ( ArrayInOut &  inoutAbsArray)
static

Computes, in place, absolute value of an array.

Parameters
inoutAbsArray[in/out] - input/output array

Definition at line 130 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
void Intrepid::RealSpaceTools< Scalar >::add ( Scalar *  sumArray,
const Scalar *  inArray1,
const Scalar *  inArray2,
const int  size 
)
static

Adds contiguous data inArray1 and inArray2 of size size:
sumArray = inArray1 + inArray2.

Parameters
sumArray[out] - sum
inArray1[in] - first summand
inArray2[in] - second summand
size[in] - size of input/output data

Definition at line 1506 of file Intrepid_RealSpaceToolsDef.hpp.

Referenced by Intrepid::CellTools< Scalar >::mapToReferenceFrameInitGuess().

template<class Scalar >
void Intrepid::RealSpaceTools< Scalar >::add ( Scalar *  inoutSumArray,
const Scalar *  inArray,
const int  size 
)
static

Adds, in place, contiguous data inArray into inoutSumArray of size size:
inoutSumArray = inoutSumArray + inArray.

Parameters
inoutSumArray[in/out] - sum / first summand
inArray[in] - second summand
size[in] - size of input/output data

Definition at line 1515 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArraySum , class ArrayIn1 , class ArrayIn2 >
void Intrepid::RealSpaceTools< Scalar >::add ( ArraySum &  sumArray,
const ArrayIn1 &  inArray1,
const ArrayIn2 &  inArray2 
)
static

Adds arrays inArray1 and inArray2:
sumArray = inArray1 + inArray2.

Parameters
sumArray[out] - sum
inArray1[in] - first summand
inArray2[in] - second summand
Note
Requirements (checked at runtime, in debug mode):
  • rank(sumArray) == rank(inArray1) == rank(inArray2)
  • dimensions(sumArray) == dimensions(inArray1) == dimensions(inArray2)

Definition at line 1525 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArraySum , class ArrayIn >
void Intrepid::RealSpaceTools< Scalar >::add ( ArraySum &  inoutSumArray,
const ArrayIn &  inArray 
)
static

Adds, in place, inArray into inoutSumArray:
inoutSumArray = inoutSumArray + inArray.

Parameters
inoutSumArray[in/out] - sum/first summand
inArray[in] - second summand
Note
Requirements (checked at runtime, in debug mode):
  • rank(inoutSumArray) == rank(inArray)
  • dimensions(inoutSumArray) == dimensions(inArray)

Definition at line 1584 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArrayIn >
Scalar Intrepid::RealSpaceTools< Scalar >::det ( const ArrayIn &  inMat)
static

Computes determinant of a single square matrix stored in an array of rank 2.

Parameters
inMat[in] - array representing a single matrix, indexed by (D, D)
Note
Requirements (checked at runtime, in debug mode):
  • rank(inMats) == 2
  • matrix dimension is limited to 1, 2, and 3

Definition at line 1240 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArrayDet , class ArrayIn >
void Intrepid::RealSpaceTools< Scalar >::det ( ArrayDet &  detArray,
const ArrayIn &  inMats 
)
static

Computes determinants of matrices stored in an array of total rank 3 (array of matrices), indexed by (i0, D, D), or 4 (array of arrays of matrices), indexed by (i0, i1, D, D).

Parameters
detArray[out] - array of determinants indexed by (i0) or (i0, i1)
inMats[in] - array of matrices indexed by (i0, D, D) or (i0, i1, D, D)
Note
Requirements (checked at runtime, in debug mode):
  • rank(detArray) == rank(inMats) - 2
  • rank(inMats) == 3 or 4
  • dimensions i0, i1 of detArray and inMats must agree
  • matrix dimensions are limited to 1, 2, and 3

Definition at line 1321 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
Scalar Intrepid::RealSpaceTools< Scalar >::dot ( const Scalar *  inArray1,
const Scalar *  inArray2,
const int  size 
)
static

Computes dot product of contiguous data inArray1 and inArray2 of size size.

Parameters
inArray1[in] - first array
inArray2[in] - second array
size[in] - size of input arrays

Definition at line 1847 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArrayVec1 , class ArrayVec2 >
Scalar Intrepid::RealSpaceTools< Scalar >::dot ( const ArrayVec1 &  inVec1,
const ArrayVec2 &  inVec2 
)
static

Computes dot product of two vectors stored in arrays of rank 1.

Parameters
inVec1[in] - first vector
inVec2[in] - second vector
Note
Requirements (checked at runtime, in debug mode):
  • rank(inVec1) == rank(inVec2) == 1
  • inVec1 and inVec2 have same dimension

Definition at line 1859 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArrayDot , class ArrayVec1 , class ArrayVec2 >
void Intrepid::RealSpaceTools< Scalar >::dot ( ArrayDot &  dotArray,
const ArrayVec1 &  inVecs1,
const ArrayVec2 &  inVecs2 
)
static

Computes dot product of vectors stored in an array of total rank 2 (array of vectors), indexed by (i0, D), or 3 (array of arrays of vectors), indexed by (i0, i1, D).

Parameters
dotArray[out] - dot product array indexed by (i0) or (i0, i1)
inVecs1[in] - first array of vectors indexed by (i0, D) or (i0, i1, D)
inVecs2[in] - second array of vectors indexed by (i0, D) or (i0, i1, D)
Note
Requirements (checked at runtime, in debug mode):
  • rank(dotArray) == rank(inVecs1) - 1 == rank(inVecs2) - 1
  • rank(inVecs1) == 2 or 3
  • dimensions i0, i1 of dotArray and inVecs1 / inVecs2 must agree

Definition at line 1882 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
void Intrepid::RealSpaceTools< Scalar >::inverse ( Scalar *  inverseMat,
const Scalar *  inMat,
const size_t  dim 
)
static

Computes inverse of the square matrix inMat of size dim by dim.

Parameters
inverseMat[out] - matrix inverse
inMat[in] - matrix
dim[in] - matrix dimension

Definition at line 653 of file Intrepid_RealSpaceToolsDef.hpp.

Referenced by Intrepid::CellTools< Scalar >::setJacobianInv().

template<class Scalar >
template<class ArrayInverse , class ArrayIn >
void Intrepid::RealSpaceTools< Scalar >::inverse ( ArrayInverse &  inverseMats,
const ArrayIn &  inMats 
)
static

Computes inverses of nonsingular matrices stored in an array of total rank 2 (single matrix), indexed by (D, D), 3 (array of matrices), indexed by (i0, D, D), or 4 (array of arrays of matrices), indexed by (i0, i1, D, D).

Parameters
inverseMats[out] - array of inverses indexed by (D, D), (i0, D, D) or (i0, i1, D, D)
inMats[in] - array of matrices indexed by (D, D), (i0, D, D) or (i0, i1, D, D)
Note
Requirements (checked at runtime, in debug mode):
  • rank(inverseMats) == rank(inMats)
  • rank(inMats) == 3 or 4
  • dimensions(inverseMats) == dimensions(inMats)
  • matrices must be square
  • matrix dimensions are limited to 1, 2, and 3

Definition at line 757 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
void Intrepid::RealSpaceTools< Scalar >::matvec ( Scalar *  matVec,
const Scalar *  inMat,
const Scalar *  inVec,
const size_t  dim 
)
static

Matrix-vector left multiply using contiguous data:
matVec = inMat * inVec.

   A single "column" vector <b><var>inVec</var></b> of size <b><var>dim</var></b> is
   multiplied on the left by a square matrix <b><var>inMat</var></b> of size
   <b><var>dim</var></b> by <b><var>dim</var></b>.
Parameters
matVec[out] - matrix-vector product
inMat[in] - the matrix argument
inVec[in] - the vector argument
dim[in] - matrix/vector dimension

Definition at line 1954 of file Intrepid_RealSpaceToolsDef.hpp.

Referenced by Intrepid::CellTools< Scalar >::mapToReferenceFrameInitGuess().

template<class Scalar >
template<class ArrayMatVec , class ArrayMat , class ArrayVec >
void Intrepid::RealSpaceTools< Scalar >::matvec ( ArrayMatVec &  matVecs,
const ArrayMat &  inMats,
const ArrayVec &  inVecs 
)
static

Matrix-vector left multiply using multidimensional arrays:
matVec = inMat * inVec.

   An array (rank 1, 2 or 3) of "column" vectors, indexed by
   (D), (i0, D) or (i0, i1, D), is multiplied on the left by an
   array (rank 2, 3 or 4) of square matrices, indexed by (D, D),
   (i0, D, D) or (i0, i1, D, D).
Parameters
matVec[out] - matrix-vector product indexed by (D), (i0, D) or (i0, i1, D)
inMat[in] - the matrix argument indexed by (D, D), (i0, D, D) or (i0, i1, D, D)
inVec[in] - the vector argument indexed by (D), (i0, D) or (i0, i1, D)
Note
Requirements (checked at runtime, in debug mode):
  • rank(matVec) == rank(inVec) == rank(inMat) - 1
  • dimensions(matVec) == dimensions(inVec)
  • matrix and vector dimensions D, i0 and i1 must agree
  • matrices are square

Definition at line 1968 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
void Intrepid::RealSpaceTools< Scalar >::scale ( Scalar *  scaledArray,
const Scalar *  inArray,
const int  size,
const Scalar  scalar 
)
static

Multiplies contiguous data inArray of size size by a scalar (componentwise):
scaledArray = scalar * inArray.

Parameters
scaledArray[out] - scaled array
inArray[in] - input array
size[in] - size of the input array
scalar[in] - multiplier

Definition at line 1762 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
void Intrepid::RealSpaceTools< Scalar >::scale ( Scalar *  inoutScaledArray,
const int  size,
const Scalar  scalar 
)
static

Multiplies, in place, contiguous data inoutScaledArray of size size by a scalar (componentwise):
inoutScaledArray = scalar * inoutScaledArray.

Parameters
inoutScaledArray[in/out] - input/scaled array
size[in] - size of array
scalar[in] - multiplier

Definition at line 1771 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArrayScaled , class ArrayIn >
void Intrepid::RealSpaceTools< Scalar >::scale ( ArrayScaled &  scaledArray,
const ArrayIn &  inArray,
const Scalar  scalar 
)
static

Multiplies array inArray by the scalar scalar (componentwise):
scaledArray = scalar * inArray.

Parameters
scaledArray[out] - scaled array
inArray[in] - input array
scalar[in] - multiplier
Note
Requirements (checked at runtime, in debug mode):
  • rank(scaledArray) == rank(inArray)
  • dimensions(scaledArray) == dimensions(inArray)

Definition at line 1781 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArrayScaled >
void Intrepid::RealSpaceTools< Scalar >::scale ( ArrayScaled &  inoutScaledArray,
const Scalar  scalar 
)
static

Multiplies, in place, array inoutScaledArray by the scalar scalar (componentwise):
inoutScaledArray = scalar * inoutScaledArray.

Parameters
inoutScaledArray[in/out] - input/output array
scalar[in] - multiplier

Definition at line 1835 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
void Intrepid::RealSpaceTools< Scalar >::subtract ( Scalar *  diffArray,
const Scalar *  inArray1,
const Scalar *  inArray2,
const int  size 
)
static

Subtracts contiguous data inArray2 from inArray1 of size size:
diffArray = inArray1 - inArray2.

Parameters
diffArray[out] - difference
inArray1[in] - minuend
inArray2[in] - subtrahend
size[in] - size of input/output data

Definition at line 1638 of file Intrepid_RealSpaceToolsDef.hpp.

Referenced by Intrepid::CellTools< Scalar >::mapToReferenceFrameInitGuess().

template<class Scalar >
void Intrepid::RealSpaceTools< Scalar >::subtract ( Scalar *  inoutDiffArray,
const Scalar *  inArray,
const int  size 
)
static

Subtracts, in place, contiguous data inArray from inoutDiffArray of size size:
inoutDiffArray = inoutDiffArray - inArray.

Parameters
inoutDiffArray[in/out] - difference/minuend
inArray[in] - subtrahend
size[in] - size of input/output data

Definition at line 1647 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArrayDiff , class ArrayIn1 , class ArrayIn2 >
void Intrepid::RealSpaceTools< Scalar >::subtract ( ArrayDiff &  diffArray,
const ArrayIn1 &  inArray1,
const ArrayIn2 &  inArray2 
)
static

Subtracts inArray2 from inArray1:
diffArray = inArray1 - inArray2.

Parameters
diffArray[out] - difference
inArray1[in] - minuend
inArray2[in] - subtrahend
Note
Requirements (checked at runtime, in debug mode):
  • rank(sumArray) == rank(inArray1) == rank(inArray2)
  • dimensions(sumArray) == dimensions(inArray1) == dimensions(inArray2)

Definition at line 1657 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArrayDiff , class ArrayIn >
void Intrepid::RealSpaceTools< Scalar >::subtract ( ArrayDiff &  inoutDiffArray,
const ArrayIn &  inArray 
)
static

Subtracts, in place, inArray from inoutDiffArray:
inoutDiffArray = inoutDiffArray - inArray.

Parameters
inoutDiffArray[in/out] - difference/minuend
inArray[in] - subtrahend
Note
Requirements (checked at runtime, in debug mode):
  • rank(inoutDiffArray) == rank(inArray)
  • dimensions(inoutDiffArray) == dimensions(inArray)

Definition at line 1710 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
void Intrepid::RealSpaceTools< Scalar >::transpose ( Scalar *  transposeMat,
const Scalar *  inMat,
const size_t  dim 
)
static

Computes transpose of the square matrix inMat of size dim by dim.

Parameters
transposeMat[out] - matrix transpose
inMat[in] - matrix
dim[in] - matrix dimension

Definition at line 576 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArrayTranspose , class ArrayIn >
void Intrepid::RealSpaceTools< Scalar >::transpose ( ArrayTranspose &  transposeMats,
const ArrayIn &  inMats 
)
static

Computes transposes of square matrices stored in an array of total rank 2 (single matrix), indexed by (D, D), 3 (array of matrices), indexed by (i0, D, D), or 4 (array of arrays of matrices), indexed by (i0, i1, D, D).

Parameters
transposeMats[out] - array of transposes indexed by (D, D), (i0, D, D) or (i0, i1, D, D)
inMats[in] - array of matrices indexed by (D, D), (i0, D, D) or (i0, i1, D, D)
Note
Requirements (checked at runtime, in debug mode):
  • rank(transposeMats) == rank(inMats)
  • rank(inMats) == 3 or 4
  • dimensions(transposeMats) == dimensions(inMats)
  • matrices must be square

Definition at line 590 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArrayVecProd , class ArrayIn1 , class ArrayIn2 >
void Intrepid::RealSpaceTools< Scalar >::vecprod ( ArrayVecProd &  vecProd,
const ArrayIn1 &  inLeft,
const ArrayIn2 &  inRight 
)
static

Vector product using multidimensional arrays:
vecProd = inVecLeft x inVecRight

     Vector multiplication of two "column" vectors stored in arrays (rank 1, 2, or 3)
     indexed by (D), (i0, D) or (i0, i1, D). 
Parameters
vecProd[in] - vector product indexed by (D), (i0, D) or (i0, i1, D)
inLeft[in] - left vector argument indexed by (D), (i0, D) or (i0, i1, D)
inRight[in] - right vector argument indexed by (D), (i0, D) or (i0, i1, D)
Todo:
Need to decide on how to handle vecprod in 2D: is the result a vector, i.e., there's dimension D or a scalar?

Definition at line 2039 of file Intrepid_RealSpaceToolsDef.hpp.

Referenced by Intrepid::CellTools< Scalar >::getPhysicalFaceNormals(), and Intrepid::CellTools< Scalar >::getReferenceFaceNormal().

template<class Scalar >
Scalar Intrepid::RealSpaceTools< Scalar >::vectorNorm ( const Scalar *  inVec,
const size_t  dim,
const ENorm  normType 
)
static

Computes norm (1, 2, infinity) of the vector inVec of size dim.

Parameters
inVec[in] - vector
dim[in] - vector dimension
normType[in] - norm type

Definition at line 139 of file Intrepid_RealSpaceToolsDef.hpp.

Referenced by Intrepid::FunctionSpaceTools::computeEdgeMeasure(), Intrepid::FunctionSpaceTools::computeFaceMeasure(), and Intrepid::CellTools< Scalar >::mapToReferenceFrameInitGuess().

template<class Scalar >
template<class ArrayIn >
Scalar Intrepid::RealSpaceTools< Scalar >::vectorNorm ( const ArrayIn &  inVec,
const ENorm  normType 
)
static

Computes norm (1, 2, infinity) of a single vector stored in an array of rank 1.

Parameters
inVec[in] - array representing a single vector
normType[in] - norm type
Note
Requirements (checked at runtime, in debug mode):
  • rank(inVec) == 1

Definition at line 170 of file Intrepid_RealSpaceToolsDef.hpp.

template<class Scalar >
template<class ArrayNorm , class ArrayIn >
void Intrepid::RealSpaceTools< Scalar >::vectorNorm ( ArrayNorm &  normArray,
const ArrayIn &  inVecs,
const ENorm  normType 
)
static

Computes norms (1, 2, infinity) of vectors stored in a array of total rank 2 (array of vectors), indexed by (i0, D), or 3 (array of arrays of vectors), indexed by (i0, i1, D).

Parameters
normArray[out] - norm array indexed by (i0) or (i0, i1)
inVecs[in] - array of vectors indexed by (i0, D) or (i0, i1, D)
normType[in] - norm type
Note
Requirements (checked at runtime, in debug mode):
  • rank(normArray) == rank(inVecs) - 1
  • rank(inVecs) == 2 or 3
  • dimensions i0, i1 of normArray and inVecs must agree

Definition at line 338 of file Intrepid_RealSpaceToolsDef.hpp.


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