Intrepid2
Public Types | Public Member Functions | Public Attributes | List of all members
Intrepid2::VectorData< Scalar, DeviceType > Class Template Reference

Reference-space field values for a basis, designed to support typical vector-valued bases. More...

#include <Intrepid2_VectorData.hpp>

Public Types

using VectorArray = Kokkos::Array< TensorData< Scalar, DeviceType >, Parameters::MaxVectorComponents >
 
using FamilyVectorArray = Kokkos::Array< VectorArray, Parameters::MaxTensorComponents >
 

Public Member Functions

void initialize ()
 Initialize members based on constructor parameters; all constructors should call this after populating numFamilies_, numComponents_, and vectorComponents_.
 
template<size_t numFamilies, size_t numComponents>
 VectorData (Kokkos::Array< Kokkos::Array< TensorData< Scalar, DeviceType >, numComponents >, numFamilies > vectorComponents)
 Standard constructor for the arbitrary case, accepting a fixed-length array argument. More...
 
 VectorData (const std::vector< std::vector< TensorData< Scalar, DeviceType > > > &vectorComponents)
 Standard constructor for the arbitrary case, accepting a variable-length std::vector argument. More...
 
template<size_t numComponents>
 VectorData (Kokkos::Array< TensorData< Scalar, DeviceType >, numComponents > vectorComponents, bool axialComponents)
 Simplified constructor for gradients of HGRAD, and values of HDIV and HCURL vector bases. More...
 
 VectorData (std::vector< TensorData< Scalar, DeviceType > > vectorComponents, bool axialComponents)
 Simplified constructor for gradients of HGRAD, and values of HDIV and HCURL vector bases. More...
 
template<typename OtherDeviceType , class = typename std::enable_if< std::is_same<typename DeviceType::memory_space, typename OtherDeviceType::memory_space>::value>::type, class = typename std::enable_if<!std::is_same<DeviceType,OtherDeviceType>::value>::type>
 VectorData (const VectorData< Scalar, OtherDeviceType > &vectorData)
 copy-like constructor for differing device type, but same memory space. This does a shallow copy of the underlying view.
 
template<typename OtherDeviceType , class = typename std::enable_if<!std::is_same<typename DeviceType::memory_space, typename OtherDeviceType::memory_space>::value>::type>
 VectorData (const VectorData< Scalar, OtherDeviceType > &vectorData)
 copy-like constructor for differing execution spaces. This does a deep copy of underlying views.
 
 VectorData (TensorData< Scalar, DeviceType > data)
 Simple 1-argument constructor for the case of trivial tensor product structure. The argument should have shape (F,P,D) where D has extent equal to the spatial dimension.
 
 VectorData (Data< Scalar, DeviceType > data)
 Simple 1-argument constructor for the case of trivial tensor product structure. The argument should have shape (F,P,D) where D has extent equal to the spatial dimension.
 
 numComponents_ (0)
 
KOKKOS_INLINE_FUNCTION bool axialComponents () const
 Returns true only if the families are so structured that the first family has nonzeros only in the x component, the second only in the y component, etc.
 
KOKKOS_INLINE_FUNCTION int numDimsForComponent (int &componentOrdinal) const
 Returns the number of dimensions corresponding to the specified component.
 
KOKKOS_INLINE_FUNCTION int numFields () const
 Returns the total number of fields; corresponds to the first dimension of this container.
 
KOKKOS_INLINE_FUNCTION int familyFieldOrdinalOffset (const int &familyOrdinal) const
 Returns the field ordinal offset for the specified family.
 
KOKKOS_INLINE_FUNCTION int numPoints () const
 Returns the number of points; corresponds to the second dimension of this container.
 
KOKKOS_INLINE_FUNCTION int spaceDim () const
 Returns the spatial dimension; corresponds to the third dimension of this container.
 
KOKKOS_INLINE_FUNCTION Scalar operator() (const int &fieldOrdinal, const int &pointOrdinal, const int &dim) const
 Accessor for the container, which has shape (F,P,D).
 
KOKKOS_INLINE_FUNCTION const
TensorData< Scalar, DeviceType > & 
getComponent (const int &componentOrdinal) const
 Single-argument component accessor for the axial-component or the single-family case; in this case, one argument suffices to uniquely identify the component. More...
 
KOKKOS_INLINE_FUNCTION const
TensorData< Scalar, DeviceType > & 
getComponent (const int &familyOrdinal, const int &componentOrdinal) const
 General component accessor. More...
 
KOKKOS_INLINE_FUNCTION int extent_int (const int &r) const
 Returns the extent in the specified dimension as an int.
 
KOKKOS_INLINE_FUNCTION unsigned rank () const
 Returns the rank of this container, which is 3.
 
KOKKOS_INLINE_FUNCTION int numComponents () const
 returns the number of components
 
KOKKOS_INLINE_FUNCTION int numFamilies () const
 returns the number of families
 
KOKKOS_INLINE_FUNCTION int familyForFieldOrdinal (const int &fieldOrdinal) const
 Returns the family ordinal corresponding to the indicated field ordinal.
 
KOKKOS_INLINE_FUNCTION int numFieldsInFamily (const unsigned &familyOrdinal) const
 returns the number of fields in the specified family
 
KOKKOS_INLINE_FUNCTION
constexpr bool 
isValid () const
 returns true for containers that have data; false for those that don't (e.g., those that have been constructed by the default constructor).
 

Public Attributes

FamilyVectorArray vectorComponents_
 
bool axialComponents_
 
int totalDimension_
 
Kokkos::Array< int,
Parameters::MaxVectorComponents
dimToComponent_
 
Kokkos::Array< int,
Parameters::MaxVectorComponents
dimToComponentDim_
 
Kokkos::Array< int,
Parameters::MaxVectorComponents
numDimsForComponent_
 
Kokkos::Array< int,
Parameters::MaxTensorComponents
familyFieldUpperBound_
 
unsigned numFamilies_
 
unsigned numComponents_
 
unsigned numPoints_
 
 true
 

Detailed Description

template<class Scalar, typename DeviceType>
class Intrepid2::VectorData< Scalar, DeviceType >

Reference-space field values for a basis, designed to support typical vector-valued bases.

VectorData is designed with typical HDIV/HCURL bases in mind; these often involve reference-space basis functions each of which is zero in every dimension but one. Moreover, on tensor product topologies each nonzero scalar function can be expressed as a product of functions defined on the topological components. Also supported: gradients of typical H^1 bases, which have a full vector, each component of which is a product of basis functions defined on lower-dimensional topologies.

Typically, HDIV and HCURL bases have multiple families corresponding to the nonzero structure of the vector components. VectorData therefore supports multiple families. (Gradients of H^1 are expressed as a single family.)

Definition at line 65 of file Intrepid2_VectorData.hpp.

Constructor & Destructor Documentation

template<class Scalar, typename DeviceType>
template<size_t numFamilies, size_t numComponents>
Intrepid2::VectorData< Scalar, DeviceType >::VectorData ( Kokkos::Array< Kokkos::Array< TensorData< Scalar, DeviceType >, numComponents >, numFamilies vectorComponents)
inline

Standard constructor for the arbitrary case, accepting a fixed-length array argument.

Parameters
[in]vectorComponents- an array of arrays, where the outer dimension corresponds to the family, and the inner to the number of components in each vector.

Outer dimension: number of families; inner dimension: number of components in each vector. Use empty/invalid TensorData objects to indicate zeroes. Each family, and each vector component dimension, must have at least one valid entry, and the number of points in each valid entry must agree with each other. The field count within components of a family must also agree; across families these may differ. Vector components are allowed to span multiple spatial dimensions, but when they do, any valid TensorData object in that vector position must agree in the dimension across families.

Definition at line 184 of file Intrepid2_VectorData.hpp.

Referenced by Intrepid2::VectorData< Scalar, ExecSpaceType >::VectorData().

template<class Scalar, typename DeviceType>
Intrepid2::VectorData< Scalar, DeviceType >::VectorData ( const std::vector< std::vector< TensorData< Scalar, DeviceType > > > &  vectorComponents)
inline

Standard constructor for the arbitrary case, accepting a variable-length std::vector argument.

Parameters
[in]vectorComponents- an array of arrays, where the outer dimension corresponds to the family, and the inner to the number of components in each vector.

Outer dimension: number of families; inner dimension: number of components in each vector. Use empty/invalid TensorData objects to indicate zeroes. Each family, and each vector component dimension, must have at least one valid entry, and the number of points in each valid entry must agree with each other. The field count within components of a family must also agree; across families these may differ. Vector components are allowed to span multiple spatial dimensions, but when they do, any valid TensorData object in that vector position must agree in the dimension across families.

Definition at line 207 of file Intrepid2_VectorData.hpp.

template<class Scalar, typename DeviceType>
template<size_t numComponents>
Intrepid2::VectorData< Scalar, DeviceType >::VectorData ( Kokkos::Array< TensorData< Scalar, DeviceType >, numComponents vectorComponents,
bool  axialComponents 
)
inline

Simplified constructor for gradients of HGRAD, and values of HDIV and HCURL vector bases.

Parameters
[in]vectorComponents- an array of components; see below for interpretation
[in]axialComponents- use false for a "full" vector, such as the gradient of a scalar HGRAD basis; use true for all-but-one zero entries in each vector, as for the values in a typical HDIV or HCURL basis

For the gradient use case, VectorData will have a single family. For the HDIV and HCURL use cases, will have the same number of families as there are components in the vectorComponents argument; each family will consist of vectors that have one entry filled, the others zeroes; this is what we mean by "axial components."

Definition at line 237 of file Intrepid2_VectorData.hpp.

template<class Scalar, typename DeviceType>
Intrepid2::VectorData< Scalar, DeviceType >::VectorData ( std::vector< TensorData< Scalar, DeviceType > >  vectorComponents,
bool  axialComponents 
)
inline

Simplified constructor for gradients of HGRAD, and values of HDIV and HCURL vector bases.

Parameters
[in]vectorComponents- an array of components; see below for interpretation
[in]axialComponents- use false for a "full" vector, such as the gradient of a scalar HGRAD basis; use true for all-but-one zero entries in each vector, as for the values in a typical HDIV or HCURL basis

For the gradient use case, VectorData will have a single family. For the HDIV and HCURL use cases, will have the same number of families as there are components in the vectorComponents argument; each family will consist of vectors that have one entry filled, the others zeroes; this is what we mean by "axial components."

Definition at line 267 of file Intrepid2_VectorData.hpp.

Member Function Documentation

template<class Scalar, typename DeviceType>
KOKKOS_INLINE_FUNCTION const TensorData<Scalar,DeviceType>& Intrepid2::VectorData< Scalar, DeviceType >::getComponent ( const int &  componentOrdinal) const
inline

Single-argument component accessor for the axial-component or the single-family case; in this case, one argument suffices to uniquely identify the component.

Parameters
[in]componentOrdinal- the vector component ordinal.

Definition at line 447 of file Intrepid2_VectorData.hpp.

Referenced by Intrepid2::BasisValues< Scalar, DeviceType >::basisValuesForFields(), and Intrepid2::VectorData< Scalar, ExecSpaceType >::VectorData().

template<class Scalar, typename DeviceType>
KOKKOS_INLINE_FUNCTION const TensorData<Scalar,DeviceType>& Intrepid2::VectorData< Scalar, DeviceType >::getComponent ( const int &  familyOrdinal,
const int &  componentOrdinal 
) const
inline

General component accessor.

Parameters
[in]familyOrdinal- the family ordinal for the requested component.
[in]componentOrdinal- the vector component ordinal.

Definition at line 471 of file Intrepid2_VectorData.hpp.

Member Data Documentation

template<class Scalar, typename DeviceType>
Intrepid2::VectorData< Scalar, DeviceType >::true
Initial value:
{}
:
numFamilies_(0)

Definition at line 340 of file Intrepid2_VectorData.hpp.


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