Intrepid2
Public Types | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
Intrepid2::TensorPoints< PointScalar, DeviceType > Class Template Reference

View-like interface to tensor points; point components are stored separately; the appropriate coordinate is determined from the composite point index and requested dimension at access time. More...

#include <Intrepid2_TensorPoints.hpp>

Public Types

using value_type = PointScalar
 

Public Member Functions

template<size_t numTensorComponents>
 TensorPoints (Kokkos::Array< ScalarView< PointScalar, DeviceType >, numTensorComponents > pointTensorComponents)
 Constructor with fixed-length Kokkos::Array argument.
 
 TensorPoints (TensorPoints otherPointsContainer, std::vector< int > whichDims)
 Constructor that takes a subset of the tensorial components of another points container.
 
 TensorPoints (std::vector< ScalarView< PointScalar, DeviceType > > pointTensorComponents)
 Constructor with variable-length std::vector argument.
 
 TensorPoints (ScalarView< PointScalar, DeviceType > points)
 Constructor for point set with trivial tensor structure.
 
template<class OtherPointsContainer >
void copyPointsContainer (ScalarView< PointScalar, DeviceType > toPoints, OtherPointsContainer fromPoints)
 Copy from one points container, which may be an arbitrary functor, to a DynRankView.
 
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>
 TensorPoints (const TensorPoints< PointScalar, OtherDeviceType > &tensorPoints)
 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>
 TensorPoints (const TensorPoints< PointScalar, OtherDeviceType > &tensorPoints)
 copy-like constructor for differing memory spaces. This does a deep_copy of the underlying view.
 
 TensorPoints ()
 Default constructor. TensorPoints::isValid() will return false.
 
ordinal_type componentPointCount (const ordinal_type &tensorComponentOrdinal) const
 Returns the number of points in the indicated component.
 
template<typename iType0 , typename iType1 >
KOKKOS_INLINE_FUNCTION std::enable_if<(std::is_integral< iType0 >::value &&std::is_integral< iType1 >::value), reference_type >::type operator() (const iType0 &tensorPointIndex, const iType1 &dim) const
 Accessor that accepts a composite point index.
 
template<typename iType0 , typename iType1 , size_t numTensorComponents>
KOKKOS_INLINE_FUNCTION std::enable_if<(std::is_integral< iType0 >::value &&std::is_integral< iType1 >::value), reference_type >::type operator() (const Kokkos::Array< iType0, numTensorComponents > &pointOrdinalComponents, const iType1 &dim) const
 Accessor that accepts a a fixed-length array with entries corresponding to component indices.
 
template<typename iType >
KOKKOS_INLINE_FUNCTION std::enable_if< std::is_integral< iType >::value, int >::type extent_int (const iType &r) const
 Returns the logical extent in the requested dimension.
 
template<typename iType >
KOKKOS_INLINE_FUNCTION constexpr std::enable_if< std::is_integral< iType >::value, size_t >::type extent (const iType &r) const
 Returns the logical extent in the requested dimension.
 
ScalarView< PointScalar, DeviceType > allocateAndFillExpandedRawPointView () const
 This method is for compatibility with existing methods that take raw point views. Note that in general it is probably better for performance to make the existing methods accept a TensorPoints object, since that can dramatically reduce memory footprint, and avoids an allocation here.
 
KOKKOS_INLINE_FUNCTION ScalarView< PointScalar, DeviceType > getTensorComponent (const ordinal_type &r) const
 Returns the requested tensor component.
 
KOKKOS_INLINE_FUNCTION 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).
 
KOKKOS_INLINE_FUNCTION ordinal_type numTensorComponents () const
 Returns the number of tensorial components.
 
KOKKOS_INLINE_FUNCTION constexpr ordinal_type rank () const
 Return the rank of the container, which is 2.
 

Protected Types

using reference_type = typename ScalarView< PointScalar, DeviceType >::reference_type
 

Protected Member Functions

void TEST_VALID_POINT_COMPONENTS ()
 
void initialize ()
 Initialize members based on constructor parameters.
 

Protected Attributes

Kokkos::Array< ScalarView< PointScalar, DeviceType >, Parameters::MaxTensorComponentspointTensorComponents_
 
ordinal_type numTensorComponents_
 
ordinal_type totalPointCount_
 
ordinal_type totalDimension_
 
Kokkos::View< ordinal_type *, DeviceType > dimToComponent_
 
Kokkos::View< ordinal_type *, DeviceType > dimToComponentDim_
 
Kokkos::Array< ordinal_type, Parameters::MaxTensorComponentspointModulus_
 
Kokkos::Array< ordinal_type, Parameters::MaxTensorComponentspointDivisor_
 
bool isValid_
 

Detailed Description

template<class PointScalar, typename DeviceType>
class Intrepid2::TensorPoints< PointScalar, DeviceType >

View-like interface to tensor points; point components are stored separately; the appropriate coordinate is determined from the composite point index and requested dimension at access time.

Definition at line 26 of file Intrepid2_TensorPoints.hpp.

Member Typedef Documentation

◆ reference_type

template<class PointScalar , typename DeviceType >
using Intrepid2::TensorPoints< PointScalar, DeviceType >::reference_type = typename ScalarView<PointScalar,DeviceType>::reference_type
protected

Definition at line 40 of file Intrepid2_TensorPoints.hpp.

◆ value_type

template<class PointScalar , typename DeviceType >
using Intrepid2::TensorPoints< PointScalar, DeviceType >::value_type = PointScalar

Definition at line 28 of file Intrepid2_TensorPoints.hpp.

Constructor & Destructor Documentation

◆ TensorPoints() [1/7]

template<class PointScalar , typename DeviceType >
template<size_t numTensorComponents>
Intrepid2::TensorPoints< PointScalar, DeviceType >::TensorPoints ( Kokkos::Array< ScalarView< PointScalar, DeviceType >, numTensorComponents pointTensorComponents)
inline

Constructor with fixed-length Kokkos::Array argument.

Parameters
[in]pointTensorComponents- the components representing the points.

TensorPoints has shape (P,D), where P is the product of the first dimensions of the component points, and D is the sum of the second dimensions.

Definition at line 107 of file Intrepid2_TensorPoints.hpp.

References Intrepid2::TensorPoints< PointScalar, DeviceType >::initialize(), and Intrepid2::TensorPoints< PointScalar, DeviceType >::numTensorComponents().

◆ TensorPoints() [2/7]

template<class PointScalar , typename DeviceType >
Intrepid2::TensorPoints< PointScalar, DeviceType >::TensorPoints ( TensorPoints< PointScalar, DeviceType >  otherPointsContainer,
std::vector< int >  whichDims 
)
inline

Constructor that takes a subset of the tensorial components of another points container.

Parameters
[in]otherPointsContainer- the original points container
[in]whichDims- the tensorial component indices to take from the other container.
Note
this does not copy the points.

Definition at line 127 of file Intrepid2_TensorPoints.hpp.

References Intrepid2::TensorPoints< PointScalar, DeviceType >::getTensorComponent(), and Intrepid2::TensorPoints< PointScalar, DeviceType >::initialize().

◆ TensorPoints() [3/7]

template<class PointScalar , typename DeviceType >
Intrepid2::TensorPoints< PointScalar, DeviceType >::TensorPoints ( std::vector< ScalarView< PointScalar, DeviceType > >  pointTensorComponents)
inline

Constructor with variable-length std::vector argument.

Parameters
[in]pointTensorComponents- the components representing the points.

TensorPoints has shape (P,D), where P is the product of the first dimensions of the component points, and D is the sum of the second dimensions.

Definition at line 147 of file Intrepid2_TensorPoints.hpp.

References Intrepid2::TensorPoints< PointScalar, DeviceType >::initialize().

◆ TensorPoints() [4/7]

template<class PointScalar , typename DeviceType >
Intrepid2::TensorPoints< PointScalar, DeviceType >::TensorPoints ( ScalarView< PointScalar, DeviceType >  points)
inline

Constructor for point set with trivial tensor structure.

Parameters
[in]points- the points, with shape (P,D).

TensorPoints has the same shape (P,D) as the input points.

Definition at line 166 of file Intrepid2_TensorPoints.hpp.

References Intrepid2::TensorPoints< PointScalar, DeviceType >::initialize().

◆ TensorPoints() [5/7]

template<class PointScalar , typename DeviceType >
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>
Intrepid2::TensorPoints< PointScalar, DeviceType >::TensorPoints ( const TensorPoints< PointScalar, OtherDeviceType > &  tensorPoints)
inline

copy-like constructor for differing device type, but same memory space. This does a shallow copy of the underlying view.

Definition at line 196 of file Intrepid2_TensorPoints.hpp.

References Intrepid2::TensorPoints< PointScalar, DeviceType >::getTensorComponent().

◆ TensorPoints() [6/7]

template<class PointScalar , typename DeviceType >
template<typename OtherDeviceType , class = typename std::enable_if<!std::is_same<typename DeviceType::memory_space, typename OtherDeviceType::memory_space>::value>::type>
Intrepid2::TensorPoints< PointScalar, DeviceType >::TensorPoints ( const TensorPoints< PointScalar, OtherDeviceType > &  tensorPoints)
inline

◆ TensorPoints() [7/7]

template<class PointScalar , typename DeviceType >
Intrepid2::TensorPoints< PointScalar, DeviceType >::TensorPoints ( )
inline

Default constructor. TensorPoints::isValid() will return false.

Definition at line 237 of file Intrepid2_TensorPoints.hpp.

Member Function Documentation

◆ allocateAndFillExpandedRawPointView()

template<class PointScalar , typename DeviceType >
ScalarView< PointScalar, DeviceType > Intrepid2::TensorPoints< PointScalar, DeviceType >::allocateAndFillExpandedRawPointView ( ) const
inline

This method is for compatibility with existing methods that take raw point views. Note that in general it is probably better for performance to make the existing methods accept a TensorPoints object, since that can dramatically reduce memory footprint, and avoids an allocation here.

Definition at line 327 of file Intrepid2_TensorPoints.hpp.

References Intrepid2::TensorPoints< PointScalar, DeviceType >::extent_int().

◆ componentPointCount()

template<class PointScalar , typename DeviceType >
ordinal_type Intrepid2::TensorPoints< PointScalar, DeviceType >::componentPointCount ( const ordinal_type &  tensorComponentOrdinal) const
inline

Returns the number of points in the indicated component.

Definition at line 245 of file Intrepid2_TensorPoints.hpp.

Referenced by Intrepid2::Basis_TensorBasis< BasisBaseClass >::allocateBasisValues().

◆ copyPointsContainer()

template<class PointScalar , typename DeviceType >
template<class OtherPointsContainer >
void Intrepid2::TensorPoints< PointScalar, DeviceType >::copyPointsContainer ( ScalarView< PointScalar, DeviceType >  toPoints,
OtherPointsContainer  fromPoints 
)
inline

Copy from one points container, which may be an arbitrary functor, to a DynRankView.

Parameters
[in]toPoints- the container to copy to.
[in]fromPoints- the container to copy from.

Definition at line 181 of file Intrepid2_TensorPoints.hpp.

Referenced by Intrepid2::TensorPoints< PointScalar, DeviceType >::TensorPoints().

◆ extent()

template<class PointScalar , typename DeviceType >
template<typename iType >
KOKKOS_INLINE_FUNCTION constexpr std::enable_if< std::is_integral< iType >::value, size_t >::type Intrepid2::TensorPoints< PointScalar, DeviceType >::extent ( const iType &  r) const
inlineconstexpr

Returns the logical extent in the requested dimension.

Parameters
[in]r- the dimension
Returns
the logical extent in the requested dimension. (This will be 1 for r >= 2, as this is a rank-2 container.)

Definition at line 319 of file Intrepid2_TensorPoints.hpp.

◆ extent_int()

template<class PointScalar , typename DeviceType >
template<typename iType >
KOKKOS_INLINE_FUNCTION std::enable_if< std::is_integral< iType >::value, int >::type Intrepid2::TensorPoints< PointScalar, DeviceType >::extent_int ( const iType &  r) const
inline

◆ getTensorComponent()

template<class PointScalar , typename DeviceType >
KOKKOS_INLINE_FUNCTION ScalarView< PointScalar, DeviceType > Intrepid2::TensorPoints< PointScalar, DeviceType >::getTensorComponent ( const ordinal_type &  r) const
inline

◆ initialize()

template<class PointScalar , typename DeviceType >
void Intrepid2::TensorPoints< PointScalar, DeviceType >::initialize ( )
inlineprotected

◆ isValid()

template<class PointScalar , typename DeviceType >
KOKKOS_INLINE_FUNCTION bool Intrepid2::TensorPoints< PointScalar, DeviceType >::isValid ( ) const
inline

Returns true for containers that have data; false for those that don't (e.g., those that have been constructed by the default constructor).

Definition at line 355 of file Intrepid2_TensorPoints.hpp.

◆ numTensorComponents()

template<class PointScalar , typename DeviceType >
KOKKOS_INLINE_FUNCTION ordinal_type Intrepid2::TensorPoints< PointScalar, DeviceType >::numTensorComponents ( ) const
inline

◆ operator()() [1/2]

template<class PointScalar , typename DeviceType >
template<typename iType0 , typename iType1 >
KOKKOS_INLINE_FUNCTION std::enable_if<(std::is_integral< iType0 >::value &&std::is_integral< iType1 >::value), reference_type >::type Intrepid2::TensorPoints< PointScalar, DeviceType >::operator() ( const iType0 &  tensorPointIndex,
const iType1 &  dim 
) const
inline

Accessor that accepts a composite point index.

Parameters
[in]tensorPointIndex- the composite point index.
[in]dim- the coordinate dimension. From tensorPointIndex, we can determine a point index for each component. From dim, we can determine which component contains the coordinate of interest, and the coordinate dimension in that component.
Returns
The component coordinate.

Definition at line 262 of file Intrepid2_TensorPoints.hpp.

◆ operator()() [2/2]

template<class PointScalar , typename DeviceType >
template<typename iType0 , typename iType1 , size_t numTensorComponents>
KOKKOS_INLINE_FUNCTION std::enable_if<(std::is_integral< iType0 >::value &&std::is_integral< iType1 >::value), reference_type >::type Intrepid2::TensorPoints< PointScalar, DeviceType >::operator() ( const Kokkos::Array< iType0, numTensorComponents > &  pointOrdinalComponents,
const iType1 &  dim 
) const
inline

Accessor that accepts a a fixed-length array with entries corresponding to component indices.

Parameters
[in]pointOrdinalComponents- the component point indices.
[in]dim- the coordinate dimension. From dim, we can determine which component contains the coordinate of interest, and the coordinate dimension in that component. We use the indicated point index for that component specified in pointOrdinalComponents to return the appropriate coordinate.
Returns
The component coordinate.

Definition at line 281 of file Intrepid2_TensorPoints.hpp.

◆ rank()

template<class PointScalar , typename DeviceType >
KOKKOS_INLINE_FUNCTION constexpr ordinal_type Intrepid2::TensorPoints< PointScalar, DeviceType >::rank ( ) const
inlineconstexpr

Return the rank of the container, which is 2.

Definition at line 369 of file Intrepid2_TensorPoints.hpp.

◆ TEST_VALID_POINT_COMPONENTS()

template<class PointScalar , typename DeviceType >
void Intrepid2::TensorPoints< PointScalar, DeviceType >::TEST_VALID_POINT_COMPONENTS ( )
inlineprotected

Definition at line 42 of file Intrepid2_TensorPoints.hpp.

Member Data Documentation

◆ dimToComponent_

template<class PointScalar , typename DeviceType >
Kokkos::View<ordinal_type*, DeviceType> Intrepid2::TensorPoints< PointScalar, DeviceType >::dimToComponent_
protected

Definition at line 34 of file Intrepid2_TensorPoints.hpp.

◆ dimToComponentDim_

template<class PointScalar , typename DeviceType >
Kokkos::View<ordinal_type*, DeviceType> Intrepid2::TensorPoints< PointScalar, DeviceType >::dimToComponentDim_
protected

Definition at line 35 of file Intrepid2_TensorPoints.hpp.

◆ isValid_

template<class PointScalar , typename DeviceType >
bool Intrepid2::TensorPoints< PointScalar, DeviceType >::isValid_
protected

Definition at line 39 of file Intrepid2_TensorPoints.hpp.

◆ numTensorComponents_

template<class PointScalar , typename DeviceType >
ordinal_type Intrepid2::TensorPoints< PointScalar, DeviceType >::numTensorComponents_
protected

Definition at line 31 of file Intrepid2_TensorPoints.hpp.

◆ pointDivisor_

template<class PointScalar , typename DeviceType >
Kokkos::Array<ordinal_type, Parameters::MaxTensorComponents> Intrepid2::TensorPoints< PointScalar, DeviceType >::pointDivisor_
protected

Definition at line 37 of file Intrepid2_TensorPoints.hpp.

◆ pointModulus_

template<class PointScalar , typename DeviceType >
Kokkos::Array<ordinal_type, Parameters::MaxTensorComponents> Intrepid2::TensorPoints< PointScalar, DeviceType >::pointModulus_
protected

Definition at line 36 of file Intrepid2_TensorPoints.hpp.

◆ pointTensorComponents_

template<class PointScalar , typename DeviceType >
Kokkos::Array< ScalarView<PointScalar,DeviceType>, Parameters::MaxTensorComponents> Intrepid2::TensorPoints< PointScalar, DeviceType >::pointTensorComponents_
protected

Definition at line 30 of file Intrepid2_TensorPoints.hpp.

◆ totalDimension_

template<class PointScalar , typename DeviceType >
ordinal_type Intrepid2::TensorPoints< PointScalar, DeviceType >::totalDimension_
protected

Definition at line 33 of file Intrepid2_TensorPoints.hpp.

◆ totalPointCount_

template<class PointScalar , typename DeviceType >
ordinal_type Intrepid2::TensorPoints< PointScalar, DeviceType >::totalPointCount_
protected

Definition at line 32 of file Intrepid2_TensorPoints.hpp.


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