Intrepid
Static Public Member Functions | Static Private Member Functions | List of all members
Intrepid::TensorProductSpaceTools Class Reference

Defines expert-level interfaces for the evaluation, differentiation and integration of finite element-functions defined by tensor products of one-dimensional spaces. These are useful in implementing spectral element methods. More...

#include <Intrepid_TensorProductSpaceTools.hpp>

Static Public Member Functions

template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
static void evaluate (ArrayTypeOut &vals, const ArrayTypeCoeffs &coeffs, const Array< RCP< ArrayTypeBasis > > &bases)
 Computes point values of a set of polynomials expressed in a tensor product basis at output points. The array coeffs is assumed to have dimensions (C,F1,F2), where F1 runs over the number of different polynomials per cell and F2 runs over the coefficients run over a tensor product basis (lowest space dimension runs fastest). The Teuchos::Array of (pointers to) Arrays bases have the one-dimensional bases tabulated at the one-dimensional points. The output array is (C,F1,P). More...
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
static void evaluateCollocated (ArrayTypeOut &vals, const ArrayTypeCoeffs &coeffs, const Array< RCP< ArrayTypeBasis > > &bases)
 Computes point values of a set of array-valued polynomials expressed in a tensor product basis at output points. The array coeffs is assumed to have dimensions (C,F1,F2), where F1 runs over the number of different polynomials per cell and F2 runs over the coefficients run over a tensor product basis (lowest space dimension runs fastest). The Teuchos::Array of (pointers to) Arrays bases have the one-dimensional bases tabulated at the one-dimensional points. The output array is (C,F1,P,D). This method assumes that the nodes for the basis coincide with the evaluation points, which leads to a big simplification. More...
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
static void evaluateGradient (ArrayTypeOut &vals, const ArrayTypeCoeffs &coeffs, const Array< RCP< ArrayTypeBasis > > &bases, const Array< RCP< ArrayTypeBasis > > &Dbases)
 Given a polynomial expressed in a tensor product basis, evaluates the gradient at a tensor product of points. The array coeffs is assumed to have dimensions (C,F1,F2), where F1 runs over the number of different polynomials per cell and F2 runs over the coefficients run over a tensor product basis (lowest space dimension runs fastest). The Teuchos::Array of (pointers to) Arrays bases and Dbases have the one-dimensional bases and their derivatives, respectively, tabulated at the one-dimensional points. The output array is (C,F1,P). More...
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
static void evaluateGradientCollocated (ArrayTypeOut &vals, const ArrayTypeCoeffs &coeffs, const Array< RCP< ArrayTypeBasis > > &bases, const Array< RCP< ArrayTypeBasis > > &Dbases)
 Given a polynomial expressed in a tensor product basis, evaluates the gradient at a tensor product of points. The array coeffs is assumed to have dimensions (C,F1,F2), where F1 runs over the number of different polynomials per cell and F2 runs over the coefficients run over a tensor product basis (lowest space dimension runs fastest). The Teuchos::Array of (pointers to) Arrays bases and Dbases have the one-dimensional bases and their derivatives, respectively, tabulated at the one-dimensional points. The output array is (C,F1,P). This method assumes that the basis nodes and the output points coincide, which leads to a considerable simplification. More...
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
static void moments (ArrayTypeOut &vals, const ArrayTypeData &data, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeWeights > > &wts)
 Computes the moments of a set of data integrated against a basis tabulated at points. More...
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
static void momentsCollocated (ArrayTypeOut &vals, const ArrayTypeData &data, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeWeights > > &wts)
 Computes the moments of a set of data integrated against a basis tabulated at points, assuming that the basis nodes and integration points coincide. More...
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
static void momentsGrad (ArrayTypeOut &vals, const ArrayTypeData &data, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeBasis > > &basisDVals, const Array< RCP< ArrayTypeWeights > > &wts)
 Computes the moments of a collection of F1 data integrated against a list of functions tabulated at points. F1 runs over the input data, F2 runs over the members of the basis. More...
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
static void momentsGradCollocated (ArrayTypeOut &vals, const ArrayTypeData &data, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeBasis > > &basisDVals, const Array< RCP< ArrayTypeWeights > > &wts)
 Computes the moments of a collection of F1 data integrated against a list of functions tabulated at points. F1 runs over the input data, F2 runs over the members of the basis. This assumes the basis nodes and integration points coincide. More...
 

Static Private Member Functions

template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
static void evaluate2D (ArrayTypeOut &vals, const ArrayTypeCoeffs &coeffs, const Array< RCP< ArrayTypeBasis > > &basisVals)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
static void evaluate3D (ArrayTypeOut &vals, const ArrayTypeCoeffs &coeffs, const Array< RCP< ArrayTypeBasis > > &basisDVals)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
static void evaluateCollocated2D (ArrayTypeOut &vals, const ArrayTypeCoeffs &coeffs, const Array< RCP< ArrayTypeBasis > > &basisVals)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
static void evaluateCollocated3D (ArrayTypeOut &vals, const ArrayTypeCoeffs &coeffs, const Array< RCP< ArrayTypeBasis > > &basisDVals)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
static void evaluateGradient2D (ArrayTypeOut &vals, const ArrayTypeCoeffs &coeffs, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeBasis > > &basisDVals)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
static void evaluateGradient3D (ArrayTypeOut &vals, const ArrayTypeCoeffs &coeffs, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeBasis > > &basisDVals)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
static void evaluateGradientCollocated2D (ArrayTypeOut &vals, const ArrayTypeCoeffs &coeffs, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeBasis > > &basisDVals)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
static void evaluateGradientCollocated3D (ArrayTypeOut &vals, const ArrayTypeCoeffs &coeffs, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeBasis > > &basisDVals)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
static void moments2D (ArrayTypeOut &vals, const ArrayTypeData &data, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeWeights > > &wts)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
static void moments3D (ArrayTypeOut &vals, const ArrayTypeData &data, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeWeights > > &wts)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
static void momentsCollocated2D (ArrayTypeOut &vals, const ArrayTypeData &data, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeWeights > > &wts)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
static void momentsCollocated3D (ArrayTypeOut &vals, const ArrayTypeData &data, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeWeights > > &wts)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
static void momentsGradCollocated2D (ArrayTypeOut &vals, const ArrayTypeData &data, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeBasis > > &basisDVals, const Array< RCP< ArrayTypeWeights > > &wts)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
static void momentsGradCollocated3D (ArrayTypeOut &vals, const ArrayTypeData &data, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeBasis > > &basisDVals, const Array< RCP< ArrayTypeWeights > > &wts)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
static void momentsGrad2D (ArrayTypeOut &vals, const ArrayTypeData &data, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeBasis > > &basisDVals, const Array< RCP< ArrayTypeWeights > > &wts)
 
template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
static void momentsGrad3D (ArrayTypeOut &vals, const ArrayTypeData &data, const Array< RCP< ArrayTypeBasis > > &basisVals, const Array< RCP< ArrayTypeBasis > > &basisDVals, const Array< RCP< ArrayTypeWeights > > &wts)
 

Detailed Description

Defines expert-level interfaces for the evaluation, differentiation and integration of finite element-functions defined by tensor products of one-dimensional spaces. These are useful in implementing spectral element methods.

Definition at line 55 of file Intrepid_TensorProductSpaceTools.hpp.

Member Function Documentation

template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
void Intrepid::TensorProductSpaceTools::evaluate ( ArrayTypeOut &  vals,
const ArrayTypeCoeffs &  coeffs,
const Array< RCP< ArrayTypeBasis > > &  bases 
)
static

Computes point values of a set of polynomials expressed in a tensor product basis at output points. The array coeffs is assumed to have dimensions (C,F1,F2), where F1 runs over the number of different polynomials per cell and F2 runs over the coefficients run over a tensor product basis (lowest space dimension runs fastest). The Teuchos::Array of (pointers to) Arrays bases have the one-dimensional bases tabulated at the one-dimensional points. The output array is (C,F1,P).

Parameters
vals[out] - output point values of the discrete function
coeffs[in] - coefficients of the input function
bases[in] - one-dimensional bases tabulated at points

Definition at line 53 of file Intrepid_TensorProductSpaceToolsDef.hpp.

template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
void Intrepid::TensorProductSpaceTools::evaluateCollocated ( ArrayTypeOut &  vals,
const ArrayTypeCoeffs &  coeffs,
const Array< RCP< ArrayTypeBasis > > &  bases 
)
static

Computes point values of a set of array-valued polynomials expressed in a tensor product basis at output points. The array coeffs is assumed to have dimensions (C,F1,F2), where F1 runs over the number of different polynomials per cell and F2 runs over the coefficients run over a tensor product basis (lowest space dimension runs fastest). The Teuchos::Array of (pointers to) Arrays bases have the one-dimensional bases tabulated at the one-dimensional points. The output array is (C,F1,P,D). This method assumes that the nodes for the basis coincide with the evaluation points, which leads to a big simplification.

Parameters
vals[out] - output point values of the discrete function
coeffs[in] - coefficients of the input function
bases[in] - one-dimensional bases tabulated at pointsComputes point values of a set of polynomials expressed in a tensor product basis at output points. The array coeffs is assumed to have dimensions (C,F1,F2), where F1 runs over the number of different polynomials per cell and F2 runs over the coefficients run over a tensor product basis (lowest space dimension runs fastest). The Teuchos::Array of (pointers to) Arrays bases have the one-dimensional bases tabulated at the one-dimensional points. The output array is (C,F1,P). This method assumes that the nodes for the basis coincide with the evaluation points, which leads to a big simplification.
vals[out] - output point values of the discrete function
coeffs[in] - coefficients of the input function
bases[in] - one-dimensional bases tabulated at points

Definition at line 108 of file Intrepid_TensorProductSpaceToolsDef.hpp.

template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
void Intrepid::TensorProductSpaceTools::evaluateGradient ( ArrayTypeOut &  vals,
const ArrayTypeCoeffs &  coeffs,
const Array< RCP< ArrayTypeBasis > > &  bases,
const Array< RCP< ArrayTypeBasis > > &  Dbases 
)
static

Given a polynomial expressed in a tensor product basis, evaluates the gradient at a tensor product of points. The array coeffs is assumed to have dimensions (C,F1,F2), where F1 runs over the number of different polynomials per cell and F2 runs over the coefficients run over a tensor product basis (lowest space dimension runs fastest). The Teuchos::Array of (pointers to) Arrays bases and Dbases have the one-dimensional bases and their derivatives, respectively, tabulated at the one-dimensional points. The output array is (C,F1,P).

Parameters
vals[out] - output point values of the discrete function
coeffs[in] - coefficients of the input function
bases[in] - one-dimensional bases tabulated at points
Dbases[in] - one-dimensional bases differentiated at points

Definition at line 224 of file Intrepid_TensorProductSpaceToolsDef.hpp.

template<class Scalar , class ArrayTypeOut , class ArrayTypeCoeffs , class ArrayTypeBasis >
void Intrepid::TensorProductSpaceTools::evaluateGradientCollocated ( ArrayTypeOut &  vals,
const ArrayTypeCoeffs &  coeffs,
const Array< RCP< ArrayTypeBasis > > &  bases,
const Array< RCP< ArrayTypeBasis > > &  Dbases 
)
static

Given a polynomial expressed in a tensor product basis, evaluates the gradient at a tensor product of points. The array coeffs is assumed to have dimensions (C,F1,F2), where F1 runs over the number of different polynomials per cell and F2 runs over the coefficients run over a tensor product basis (lowest space dimension runs fastest). The Teuchos::Array of (pointers to) Arrays bases and Dbases have the one-dimensional bases and their derivatives, respectively, tabulated at the one-dimensional points. The output array is (C,F1,P). This method assumes that the basis nodes and the output points coincide, which leads to a considerable simplification.

Parameters
vals[out] - output point values of the discrete function
coeffs[in] - coefficients of the input function
bases[in] - one-dimensional bases tabulated at points
Dbases[in] - one-dimensional bases differentiated at points

Definition at line 294 of file Intrepid_TensorProductSpaceToolsDef.hpp.

template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
void Intrepid::TensorProductSpaceTools::moments ( ArrayTypeOut &  vals,
const ArrayTypeData &  data,
const Array< RCP< ArrayTypeBasis > > &  basisVals,
const Array< RCP< ArrayTypeWeights > > &  wts 
)
static

Computes the moments of a set of data integrated against a basis tabulated at points.

Parameters
vals[out] - (C,F1,F2) output moments of the data against the basis functions
data[in] - (C,F1,P) data tabulated at the tensor product of points
basisvals[in] - one-dimensional bases tabulated at points
wts[in] - array of one-dimensional quadrature weights

Definition at line 364 of file Intrepid_TensorProductSpaceToolsDef.hpp.

template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
void Intrepid::TensorProductSpaceTools::momentsCollocated ( ArrayTypeOut &  vals,
const ArrayTypeData &  data,
const Array< RCP< ArrayTypeBasis > > &  basisVals,
const Array< RCP< ArrayTypeWeights > > &  wts 
)
static

Computes the moments of a set of data integrated against a basis tabulated at points, assuming that the basis nodes and integration points coincide.

Parameters
vals[out] - (C,F1,F2) output moments of the data against the basis functions
data[in] - (C,F1,P) data tabulated at the tensor product of points
basisvals[in] - one-dimensional bases tabulated at points
wts[in] - array of one-dimensional quadrature weights

Definition at line 429 of file Intrepid_TensorProductSpaceToolsDef.hpp.

template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
void Intrepid::TensorProductSpaceTools::momentsGrad ( ArrayTypeOut &  vals,
const ArrayTypeData &  data,
const Array< RCP< ArrayTypeBasis > > &  basisVals,
const Array< RCP< ArrayTypeBasis > > &  basisDVals,
const Array< RCP< ArrayTypeWeights > > &  wts 
)
static

Computes the moments of a collection of F1 data integrated against a list of functions tabulated at points. F1 runs over the input data, F2 runs over the members of the basis.

Parameters
vals[out] - (C,F1,F2) output moments of the data against the basis functions
data[in] - (C,F1,P,D) data tabulated at the tensor product of points
basisvals[in] - one-dimensional bases tabulated at points
basisDvals[in] - one-dimensional bases differentated at points
wts[in] - one-dimensional quadrature weights

Definition at line 494 of file Intrepid_TensorProductSpaceToolsDef.hpp.

template<class Scalar , class ArrayTypeOut , class ArrayTypeData , class ArrayTypeBasis , class ArrayTypeWeights >
void Intrepid::TensorProductSpaceTools::momentsGradCollocated ( ArrayTypeOut &  vals,
const ArrayTypeData &  data,
const Array< RCP< ArrayTypeBasis > > &  basisVals,
const Array< RCP< ArrayTypeBasis > > &  basisDVals,
const Array< RCP< ArrayTypeWeights > > &  wts 
)
static

Computes the moments of a collection of F1 data integrated against a list of functions tabulated at points. F1 runs over the input data, F2 runs over the members of the basis. This assumes the basis nodes and integration points coincide.

Parameters
vals[out] - (C,F1,F2) output moments of the data against the basis functions
data[in] - (C,F1,P,D) data tabulated at the tensor product of points
basisvals[in] - one-dimensional bases tabulated at points
basisDvals[in] - one-dimensional bases differentated at points
wts[in] - one-dimensional quadrature weights

Definition at line 563 of file Intrepid_TensorProductSpaceToolsDef.hpp.


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