Intrepid2
Classes | Public Types | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Static Private Member Functions | Private Attributes | List of all members
Intrepid2::Data< DataScalar, DeviceType > Class Template Reference

Wrapper around a Kokkos::View that allows data that is constant or repeating in various logical dimensions to be stored just once, while providing a similar interface to that of View. More...

#include <Intrepid2_Data.hpp>

Classes

struct  bool_pack
 

Public Types

using value_type = DataScalar
 
using execution_space = typename DeviceType::execution_space
 
using reference_type = typename ScalarView< DataScalar, DeviceType >::reference_type
 
using const_reference_type = typename ScalarView< const DataScalar, DeviceType >::reference_type
 
template<bool... v>
using all_true = std::is_same< bool_pack< true, v...>, bool_pack< v..., true >>
 
template<class... IntArgs>
using valid_args = all_true< std::is_integral< IntArgs >{}...>
 

Public Member Functions

template<class UnaryOperator >
void applyOperator (UnaryOperator unaryOperator)
 applies the specified unary operator to each entry
 
template<class... IntArgs>
KOKKOS_INLINE_FUNCTION
reference_type 
getWritableEntryWithPassThroughOption (const bool &passThroughBlockDiagonalArgs, const IntArgs...intArgs) const
 Returns an l-value reference to the specified logical entry in the underlying view. Note that for variation types other than GENERAL, multiple valid argument sets will refer to the same memory location. Intended for Intrepid2 developers and expert users only. If passThroughBlockDiagonalArgs is TRUE, the corresponding arguments are interpreted as entries in the 1D packed matrix rather than as logical 2D matrix row and column.
 
template<class... IntArgs>
KOKKOS_INLINE_FUNCTION
reference_type 
getWritableEntry (const IntArgs...intArgs) const
 Returns an l-value reference to the specified logical entry in the underlying view. Note that for variation types other than GENERAL, multiple valid argument sets will refer to the same memory location. Intended for Intrepid2 developers and expert users only. If passThroughBlockDiagonalArgs is TRUE, the corresponding arguments are interpreted as entries in the 1D packed matrix rather than as logical 2D matrix row and column.
 
void allocateAndCopyFromDynRankView (ScalarView< DataScalar, DeviceType > data)
 allocate an underlying View that matches the provided DynRankView in dimensions, and copy. Called by constructors that accept a DynRankView as argument.
 
 Data (std::vector< DimensionInfo > dimInfoVector)
 Constructor in terms of DimensionInfo for each logical dimension; does not require a View to be specified. Will allocate a View of appropriate rank, zero-filled.
 
 variationType_ ({CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT})
 
 blockPlusDiagonalLastNonDiagonal_ (-1)
 
 rank_ (dimInfoVector.size())
 
 Data (const ScalarView< DataScalar, DeviceType > &data, int rank, Kokkos::Array< int, 7 > extents, Kokkos::Array< DataVariationType, 7 > variationType, const int blockPlusDiagonalLastNonDiagonal=-1)
 DynRankView constructor. Will copy to a View of appropriate rank.
 
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>
 Data (const Data< DataScalar, OtherDeviceType > &data)
 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>
 Data (const Data< DataScalar, OtherDeviceType > &data)
 copy-like constructor for differing execution spaces. This does a deep_copy of the underlying view.
 
 Data (ScalarView< DataScalar, DeviceType > data)
 copy constructor modeled after the copy-like constructor above. Not as efficient as the implicit copy constructor, so this should only be uncommented to emulate the copy-like constructor above in situations where host and device space are the same, and the copy-like constructor does not exist, or to provide a debugging breakpoint to assess when copies are being constructed. More...
 
 extents_ ({1, 1, 1, 1, 1, 1, 1})
 
 variationType_ ({CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT})
 
 blockPlusDiagonalLastNonDiagonal_ (blockPlusDiagonalLastNonDiagonal)
 
 rank_ (rank)
 
template<size_t rank, class... ViewProperties>
 Data (Kokkos::View< DataScalar *, DeviceType, ViewProperties...> data, Kokkos::Array< int, rank > extents, Kokkos::Array< DataVariationType, rank > variationType, const int blockPlusDiagonalLastNonDiagonal=-1)
 
 variationType_ ({CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT})
 
 blockPlusDiagonalLastNonDiagonal_ (blockPlusDiagonalLastNonDiagonal)
 
 rank_ (rank)
 
template<size_t rank, class... ViewProperties>
 Data (Kokkos::View< DataScalar **, DeviceType, ViewProperties...> data, Kokkos::Array< int, rank > extents, Kokkos::Array< DataVariationType, rank > variationType, const int blockPlusDiagonalLastNonDiagonal=-1)
 
 variationType_ ({CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT})
 
 blockPlusDiagonalLastNonDiagonal_ (blockPlusDiagonalLastNonDiagonal)
 
 rank_ (rank)
 
template<size_t rank, class... ViewProperties>
 Data (Kokkos::View< DataScalar ***, DeviceType, ViewProperties...> data, Kokkos::Array< int, rank > extents, Kokkos::Array< DataVariationType, rank > variationType, const int blockPlusDiagonalLastNonDiagonal=-1)
 
 variationType_ ({CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT})
 
 blockPlusDiagonalLastNonDiagonal_ (blockPlusDiagonalLastNonDiagonal)
 
 rank_ (rank)
 
template<size_t rank, class... ViewProperties>
 Data (Kokkos::View< DataScalar ****, DeviceType, ViewProperties...> data, Kokkos::Array< int, rank > extents, Kokkos::Array< DataVariationType, rank > variationType, const int blockPlusDiagonalLastNonDiagonal=-1)
 
 variationType_ ({CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT})
 
 blockPlusDiagonalLastNonDiagonal_ (blockPlusDiagonalLastNonDiagonal)
 
 rank_ (rank)
 
template<size_t rank, class... ViewProperties>
 Data (Kokkos::View< DataScalar *****, DeviceType, ViewProperties...> data, Kokkos::Array< int, rank > extents, Kokkos::Array< DataVariationType, rank > variationType, const int blockPlusDiagonalLastNonDiagonal=-1)
 
 variationType_ ({CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT})
 
 blockPlusDiagonalLastNonDiagonal_ (blockPlusDiagonalLastNonDiagonal)
 
 rank_ (rank)
 
template<size_t rank, class... ViewProperties>
 Data (Kokkos::View< DataScalar ******, DeviceType, ViewProperties...> data, Kokkos::Array< int, rank > extents, Kokkos::Array< DataVariationType, rank > variationType, const int blockPlusDiagonalLastNonDiagonal=-1)
 
 variationType_ ({CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT})
 
 blockPlusDiagonalLastNonDiagonal_ (blockPlusDiagonalLastNonDiagonal)
 
 rank_ (rank)
 
template<size_t rank, class... ViewProperties>
 Data (Kokkos::View< DataScalar *******, DeviceType, ViewProperties...> data, Kokkos::Array< int, rank > extents, Kokkos::Array< DataVariationType, rank > variationType, const int blockPlusDiagonalLastNonDiagonal=-1)
 
 variationType_ ({CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT})
 
 blockPlusDiagonalLastNonDiagonal_ (blockPlusDiagonalLastNonDiagonal)
 
 rank_ (rank)
 
template<class ViewScalar , class... ViewProperties>
 Data (const unsigned rank, Kokkos::View< ViewScalar, DeviceType, ViewProperties...> data, Kokkos::Array< int, 7 > extents, Kokkos::Array< DataVariationType, 7 > variationType, const int blockPlusDiagonalLastNonDiagonal=-1)
 constructor with run-time rank (requires full-length extents, variationType inputs; those beyond the rank will be ignored).
 
 variationType_ ({CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT})
 
 blockPlusDiagonalLastNonDiagonal_ (blockPlusDiagonalLastNonDiagonal)
 
 rank_ (rank)
 
template<size_t rank>
 Data (DataScalar constantValue, Kokkos::Array< int, rank > extents)
 constructor for everywhere-constant data
 
 variationType_ ({CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT})
 
 blockPlusDiagonalLastNonDiagonal_ (-1)
 
 rank_ (rank)
 
 Data ()
 default constructor (empty data)
 
 variationType_ ({CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT, CONSTANT})
 
 blockPlusDiagonalLastNonDiagonal_ (-1)
 
 rank_ (0)
 
KOKKOS_INLINE_FUNCTION const int & blockPlusDiagonalLastNonDiagonal () const
 For a Data object containing data with variation type BLOCK_PLUS_DIAGONAL, returns the row and column (which must match) of the last non-diagonal entry. For fully diagonal matrices, this is -1.
 
KOKKOS_INLINE_FUNCTION
Kokkos::Array< int, 7 > 
getExtents () const
 Returns an array containing the logical extents in each dimension.
 
KOKKOS_INLINE_FUNCTION
DimensionInfo 
getDimensionInfo (const int &dim) const
 Returns an object fully specifying the indicated dimension. This is used in determining appropriate sizing of Data objects that depend on other Data objects (e.g., when two matrix Data objects are multiplied together).
 
KOKKOS_INLINE_FUNCTION
DimensionInfo 
combinedDataDimensionInfo (const Data &otherData, const int &dim) const
 Returns (DataVariationType, data extent) in the specified dimension for a Data container that combines (through multiplication, say, or addition) this container with otherData.
 
template<int rank>
KOKKOS_INLINE_FUNCTION
enable_if_t< rank==1, const
Kokkos::View< typename
RankExpander< DataScalar, rank >
::value_type, DeviceType > & > 
getUnderlyingView () const
 Returns the underlying view. Throws an exception if the underlying view is not rank 1.
 
template<int rank>
KOKKOS_INLINE_FUNCTION
enable_if_t< rank==2, const
Kokkos::View< typename
RankExpander< DataScalar, rank >
::value_type, DeviceType > & > 
getUnderlyingView () const
 Returns the underlying view. Throws an exception if the underlying view is not rank 2.
 
template<int rank>
KOKKOS_INLINE_FUNCTION
enable_if_t< rank==3, const
Kokkos::View< typename
RankExpander< DataScalar, rank >
::value_type, DeviceType > & > 
getUnderlyingView () const
 Returns the underlying view. Throws an exception if the underlying view is not rank 3.
 
template<int rank>
KOKKOS_INLINE_FUNCTION
enable_if_t< rank==4, const
Kokkos::View< typename
RankExpander< DataScalar, rank >
::value_type, DeviceType > & > 
getUnderlyingView () const
 Returns the underlying view. Throws an exception if the underlying view is not rank 4.
 
template<int rank>
KOKKOS_INLINE_FUNCTION
enable_if_t< rank==5, const
Kokkos::View< typename
RankExpander< DataScalar, rank >
::value_type, DeviceType > & > 
getUnderlyingView () const
 Returns the underlying view. Throws an exception if the underlying view is not rank 5.
 
template<int rank>
KOKKOS_INLINE_FUNCTION
enable_if_t< rank==6, const
Kokkos::View< typename
RankExpander< DataScalar, rank >
::value_type, DeviceType > & > 
getUnderlyingView () const
 Returns the underlying view. Throws an exception if the underlying view is not rank 6.
 
template<int rank>
KOKKOS_INLINE_FUNCTION
enable_if_t< rank==7, const
Kokkos::View< typename
RankExpander< DataScalar, rank >
::value_type, DeviceType > & > 
getUnderlyingView () const
 Returns the underlying view. Throws an exception if the underlying view is not rank 7.
 
KOKKOS_INLINE_FUNCTION const
Kokkos::View< DataScalar
*, DeviceType > & 
getUnderlyingView1 () const
 returns the View that stores the unique data. For rank-1 underlying containers.
 
KOKKOS_INLINE_FUNCTION const
Kokkos::View< DataScalar
**, DeviceType > & 
getUnderlyingView2 () const
 returns the View that stores the unique data. For rank-2 underlying containers.
 
KOKKOS_INLINE_FUNCTION const
Kokkos::View< DataScalar
***, DeviceType > & 
getUnderlyingView3 () const
 returns the View that stores the unique data. For rank-3 underlying containers.
 
KOKKOS_INLINE_FUNCTION const
Kokkos::View< DataScalar
****, DeviceType > & 
getUnderlyingView4 () const
 returns the View that stores the unique data. For rank-4 underlying containers.
 
KOKKOS_INLINE_FUNCTION const
Kokkos::View< DataScalar
*****, DeviceType > & 
getUnderlyingView5 () const
 returns the View that stores the unique data. For rank-5 underlying containers.
 
KOKKOS_INLINE_FUNCTION const
Kokkos::View< DataScalar
******, DeviceType > & 
getUnderlyingView6 () const
 returns the View that stores the unique data. For rank-6 underlying containers.
 
KOKKOS_INLINE_FUNCTION const
Kokkos::View< DataScalar
*******, DeviceType > & 
getUnderlyingView7 () const
 returns the View that stores the unique data. For rank-7 underlying containers.
 
KOKKOS_INLINE_FUNCTION void setUnderlyingView1 (const Kokkos::View< DataScalar *, DeviceType > &view)
 sets the View that stores the unique data. For rank-1 underlying containers.
 
KOKKOS_INLINE_FUNCTION void setUnderlyingView2 (const Kokkos::View< DataScalar **, DeviceType > &view)
 sets the View that stores the unique data. For rank-2 underlying containers.
 
KOKKOS_INLINE_FUNCTION void setUnderlyingView3 (const Kokkos::View< DataScalar ***, DeviceType > &view)
 sets the View that stores the unique data. For rank-3 underlying containers.
 
KOKKOS_INLINE_FUNCTION void setUnderlyingView4 (const Kokkos::View< DataScalar ****, DeviceType > &view)
 sets the View that stores the unique data. For rank-4 underlying containers.
 
KOKKOS_INLINE_FUNCTION void setUnderlyingView5 (const Kokkos::View< DataScalar *****, DeviceType > &view)
 sets the View that stores the unique data. For rank-5 underlying containers.
 
KOKKOS_INLINE_FUNCTION void setUnderlyingView6 (const Kokkos::View< DataScalar ******, DeviceType > &view)
 sets the View that stores the unique data. For rank-6 underlying containers.
 
KOKKOS_INLINE_FUNCTION void setUnderlyingView7 (const Kokkos::View< DataScalar *******, DeviceType > &view)
 sets the View that stores the unique data. For rank-7 underlying containers.
 
template<int underlyingRank, class ViewScalar >
KOKKOS_INLINE_FUNCTION void setUnderlyingView (const Kokkos::View< ViewScalar, DeviceType > &view)
 
ScalarView< DataScalar,
DeviceType > 
getUnderlyingView () const
 Returns a DynRankView constructed atop the same underlying data as the fixed-rank Kokkos::View used internally.
 
KOKKOS_INLINE_FUNCTION ordinal_type getUnderlyingViewRank () const
 returns the rank of the View that stores the unique data
 
KOKKOS_INLINE_FUNCTION ordinal_type getUnderlyingViewSize () const
 returns the number of entries in the View that stores the unique data
 
ScalarView< DataScalar,
DeviceType > 
allocateDynRankViewMatchingUnderlying () const
 Returns a DynRankView that matches the underlying Kokkos::View object in value_type, layout, and dimension.
 
template<class... DimArgs>
ScalarView< DataScalar,
DeviceType > 
allocateDynRankViewMatchingUnderlying (DimArgs...dims) const
 Returns a DynRankView that matches the underlying Kokkos::View object value_type and layout, but with the specified dimensions.
 
void clear () const
 Copies 0.0 to the underlying View.
 
void copyDataFromDynRankViewMatchingUnderlying (const ScalarView< DataScalar, DeviceType > &dynRankView) const
 Copies from the provided DynRankView into the underlying Kokkos::View container storing the unique data.
 
KOKKOS_INLINE_FUNCTION int getDataExtent (const ordinal_type &d) const
 returns the true extent of the data corresponding to the logical dimension provided; if the data does not vary in that dimension, returns 1
 
KOKKOS_INLINE_FUNCTION int getVariationModulus (const int &d) const
 Variation modulus accessor. More...
 
KOKKOS_INLINE_FUNCTION const
Kokkos::Array
< DataVariationType, 7 > & 
getVariationTypes () const
 Returns an array with the variation types in each logical dimension.
 
template<class... IntArgs>
KOKKOS_INLINE_FUNCTION return_type getEntryWithPassThroughOption (const bool &passThroughBlockDiagonalArgs, const IntArgs &...intArgs) const
 Returns a (read-only) value corresponding to the specified logical data location. If passThroughBlockDiagonalArgs is TRUE, the corresponding arguments are interpreted as entries in the 1D packed matrix rather than as logical 2D matrix row and column.
 
template<class... IntArgs>
KOKKOS_INLINE_FUNCTION return_type getEntry (const IntArgs &...intArgs) const
 Returns a (read-only) value corresponding to the specified logical data location.
 
template<class... IntArgs>
KOKKOS_INLINE_FUNCTION
enable_if_t< valid_args
< IntArgs...>::value
&&(sizeof...(IntArgs)
<=7), return_type > 
operator() (const IntArgs &...intArgs) const
 Returns a value corresponding to the specified logical data location.
 
KOKKOS_INLINE_FUNCTION int extent_int (const int &r) const
 Returns the logical extent in the specified dimension.
 
template<typename iType >
KOKKOS_INLINE_FUNCTION
constexpr std::enable_if
< std::is_integral< iType >
::value, size_t >::type 
extent (const iType &r) const
 
KOKKOS_INLINE_FUNCTION bool isDiagonal () const
 returns true for containers that have two dimensions marked as BLOCK_PLUS_DIAGONAL for which the non-diagonal block is empty or size 1.
 
Data< DataScalar, DeviceType > allocateConstantData (const DataScalar &value)
 
template<int rank>
enable_if_t<(rank!=1)&&(rank!=7),
Kokkos::MDRangePolicy
< typename
DeviceType::execution_space,
Kokkos::Rank< rank > > > 
dataExtentRangePolicy ()
 returns an MDRangePolicy over the underlying data extents (but with the logical shape).
 
template<int rank>
enable_if_t< rank==7,
Kokkos::MDRangePolicy
< typename
DeviceType::execution_space,
Kokkos::Rank< 6 > > > 
dataExtentRangePolicy ()
 returns an MDRangePolicy over the first six underlying data extents (but with the logical shape).
 
template<int rank>
enable_if_t< rank==1,
Kokkos::RangePolicy< typename
DeviceType::execution_space > > 
dataExtentRangePolicy ()
 
Data shallowCopy (const int rank, const Kokkos::Array< int, 7 > &extents, const Kokkos::Array< DataVariationType, 7 > &variationTypes) const
 Creates a new Data object with the same underlying view, but with the specified logical rank, extents, and variation types.
 
template<class BinaryOperator >
void storeInPlaceCombination (const Data< DataScalar, DeviceType > &A, const Data< DataScalar, DeviceType > &B, BinaryOperator binaryOperator)
 Places the result of an in-place combination (e.g., entrywise sum) into this data container.
 
void storeInPlaceSum (const Data< DataScalar, DeviceType > &A, const Data< DataScalar, DeviceType > &B)
 stores the in-place (entrywise) sum, A .+ B, into this container.
 
void storeInPlaceProduct (const Data< DataScalar, DeviceType > &A, const Data< DataScalar, DeviceType > &B)
 stores the in-place (entrywise) product, A .* B, into this container.
 
void storeInPlaceDifference (const Data< DataScalar, DeviceType > &A, const Data< DataScalar, DeviceType > &B)
 stores the in-place (entrywise) difference, A .- B, into this container.
 
void storeInPlaceQuotient (const Data< DataScalar, DeviceType > &A, const Data< DataScalar, DeviceType > &B)
 stores the in-place (entrywise) quotient, A ./ B, into this container.
 
void storeMatVec (const Data< DataScalar, DeviceType > &matData, const Data< DataScalar, DeviceType > &vecData)
 Places the result of a matrix-vector multiply corresponding to the two provided containers into this Data container. This Data container should have been constructed by a call to allocateMatVecResult(), or should match such a container in underlying data extent and variation types.
 
void storeMatMat (const bool transposeA, const Data< DataScalar, DeviceType > &A_MatData, const bool transposeB, const Data< DataScalar, DeviceType > &B_MatData)
 
KOKKOS_INLINE_FUNCTION
constexpr bool 
isValid () const
 returns true for containers that have data; false for those that don't (namely, those that have been constructed by the default constructor).
 
KOKKOS_INLINE_FUNCTION unsigned rank () const
 Returns the logical rank of the Data container.
 
void setExtent (const ordinal_type &d, const ordinal_type &newExtent)
 sets the logical extent in the specified dimension. If needed, the underlying data container is resized. More...
 
KOKKOS_INLINE_FUNCTION bool underlyingMatchesLogical () const
 Returns true if the underlying container has exactly the same rank and extents as the logical container.
 

Static Public Member Functions

template<class ToContainer , class FromContainer >
static void copyContainer (ToContainer to, FromContainer from)
 Generic data copying method to allow construction of Data object from DynRankViews for which deep_copy() to the underlying view would be disallowed. This method made public to allow CUDA compilation (because it contains a Kokkos lambda).
 
static Data< DataScalar,
DeviceType > 
allocateInPlaceCombinationResult (const Data< DataScalar, DeviceType > &A, const Data< DataScalar, DeviceType > &B)
 
static Data< DataScalar,
DeviceType > 
allocateMatMatResult (const bool transposeA, const Data< DataScalar, DeviceType > &A_MatData, const bool transposeB, const Data< DataScalar, DeviceType > &B_MatData)
 
static Data< DataScalar,
DeviceType > 
allocateMatVecResult (const Data< DataScalar, DeviceType > &matData, const Data< DataScalar, DeviceType > &vecData)
 

Private Types

using return_type = const_reference_type
 

Private Member Functions

KOKKOS_INLINE_FUNCTION int getUnderlyingViewExtent (const int &dim) const
 Returns the extent of the underlying view in the specified dimension.
 
void setActiveDims ()
 class initialization method. Called by constructors.
 

Static Private Member Functions

static KOKKOS_INLINE_FUNCTION int blockPlusDiagonalNumNondiagonalEntries (const int &lastNondiagonal)
 Returns the number of non-diagonal entries based on the last non-diagonal. Only applicable for BLOCK_PLUS_DIAGONAL DataVariationType.
 
static KOKKOS_INLINE_FUNCTION int blockPlusDiagonalBlockEntryIndex (const int &lastNondiagonal, const int &numNondiagonalEntries, const int &i, const int &j)
 //! Returns flattened index of the specified (i,j) matrix entry, assuming that i,j ≤ lastNondiagonal. Only applicable for BLOCK_PLUS_DIAGONAL DataVariationType.
 
static KOKKOS_INLINE_FUNCTION int blockPlusDiagonalDiagonalEntryIndex (const int &lastNondiagonal, const int &numNondiagonalEntries, const int &i)
 Returns flattened index of the specified (i,i) matrix entry, assuming that i > lastNondiagonal. Only applicable for BLOCK_PLUS_DIAGONAL DataVariationType.
 

Private Attributes

ordinal_type dataRank_
 
Kokkos::View< DataScalar
*, DeviceType > 
data1_
 
Kokkos::View< DataScalar
**, DeviceType > 
data2_
 
Kokkos::View< DataScalar
***, DeviceType > 
data3_
 
Kokkos::View< DataScalar
****, DeviceType > 
data4_
 
Kokkos::View< DataScalar
*****, DeviceType > 
data5_
 
Kokkos::View< DataScalar
******, DeviceType > 
data6_
 
Kokkos::View< DataScalar
*******, DeviceType > 
data7_
 
Kokkos::Array< int, 7 > extents_
 
Kokkos::Array
< DataVariationType, 7 > 
variationType_
 
Kokkos::Array< int, 7 > variationModulus_
 
int blockPlusDiagonalLastNonDiagonal_ = -1
 
bool hasNontrivialModulusUNUSED_
 
bool underlyingMatchesLogical_
 
Kokkos::Array< ordinal_type, 7 > activeDims_
 
int numActiveDims_
 
ordinal_type rank_
 
ScalarView< DataScalar,
DeviceType > 
zeroView_
 

Detailed Description

template<class DataScalar, typename DeviceType>
class Intrepid2::Data< DataScalar, DeviceType >

Wrapper around a Kokkos::View that allows data that is constant or repeating in various logical dimensions to be stored just once, while providing a similar interface to that of View.

The Data class distinguishes between the logical extent and the data extent. For example, one could construct a data container corresponding to constant (cell, point) data with 100 cells and 25 points per cell as follows: auto cpData = Data(value, Kokkos::Array<int>{100,25}); The data extent of the container is 1 in every dimension, while the logical extent is 100 in the first dimension, and 25 in the second. Similarly, the logical rank of the container is 2, but the rank of the underlying View is 1.

There are four possible variation types in a logical dimension:

Definition at line 69 of file Intrepid2_Data.hpp.

Constructor & Destructor Documentation

template<class DataScalar, typename DeviceType>
Intrepid2::Data< DataScalar, DeviceType >::Data ( ScalarView< DataScalar, DeviceType >  data)
inline

copy constructor modeled after the copy-like constructor above. Not as efficient as the implicit copy constructor, so this should only be uncommented to emulate the copy-like constructor above in situations where host and device space are the same, and the copy-like constructor does not exist, or to provide a debugging breakpoint to assess when copies are being constructed.

constructor for fully varying data container, no expressed redundancy/repetition. Copies the data to a new Kokkos::View of matching rank.

Definition at line 680 of file Intrepid2_Data.hpp.

Member Function Documentation

template<class DataScalar, typename DeviceType>
Data<DataScalar,DeviceType> Intrepid2::Data< DataScalar, DeviceType >::allocateConstantData ( const DataScalar &  value)
inline

Constructs a container with extents matching this, with a single-value underlying View, CONSTANT in each dimension.

Parameters
value[in] - the constant value.
Returns
A container with the same logical shape as this, with a single-value underlying View.

Definition at line 1319 of file Intrepid2_Data.hpp.

Referenced by Intrepid2::CellTools< DeviceType >::setJacobianDetInv().

template<class DataScalar, typename DeviceType>
static Data<DataScalar,DeviceType> Intrepid2::Data< DataScalar, DeviceType >::allocateInPlaceCombinationResult ( const Data< DataScalar, DeviceType > &  A,
const Data< DataScalar, DeviceType > &  B 
)
inlinestatic

Constructs a container suitable for storing the result of an in-place combination of the two provided data containers. The two containers must have the same logical shape.

See Also
storeInPlaceCombination()
Parameters
A[in] - the first data container.
B[in] - the second data container. Must have the same logical shape as A.
Returns
A container with the same logical shape as A and B, with underlying View storage sufficient to store the result of A + B (or any other in-place combination).

Definition at line 1329 of file Intrepid2_Data.hpp.

Referenced by Intrepid2::IntegrationTools< DeviceType >::integrate(), Intrepid2::DataTools::multiplyByCPWeights(), and Intrepid2::TransformedBasisValues< Scalar, DeviceType >::multiplyByPointwiseWeights().

template<class DataScalar, typename DeviceType>
static Data<DataScalar,DeviceType> Intrepid2::Data< DataScalar, DeviceType >::allocateMatMatResult ( const bool  transposeA,
const Data< DataScalar, DeviceType > &  A_MatData,
const bool  transposeB,
const Data< DataScalar, DeviceType > &  B_MatData 
)
inlinestatic

Constructs a container suitable for storing the result of a matrix-vector multiply corresponding to the two provided containers.

See Also
storeMatMat()
Parameters
A_MatData[in] - logically (...,D1,D2)-dimensioned container, where D1,D2 correspond to matrix dimensions.
transposeA[in] - if true, A will be transposed prior to being multiplied by B (or B's transpose).
B_MatData[in] - logically (...,D3,D4)-dimensioned container, where D3,D4 correspond to matrix dimensions.
transposeB[in] - if true, B will be transposed prior to the multiplication by A (or A's transpose).

Definition at line 1349 of file Intrepid2_Data.hpp.

Referenced by Intrepid2::IntegrationTools< DeviceType >::integrate().

template<class DataScalar, typename DeviceType>
static Data<DataScalar,DeviceType> Intrepid2::Data< DataScalar, DeviceType >::allocateMatVecResult ( const Data< DataScalar, DeviceType > &  matData,
const Data< DataScalar, DeviceType > &  vecData 
)
inlinestatic

Constructs a container suitable for storing the result of a matrix-vector multiply corresponding to the two provided containers.

See Also
storeMatVec()

Definition at line 1531 of file Intrepid2_Data.hpp.

Referenced by Intrepid2::ProjectedGeometry< spaceDim, PointScalar, DeviceType >::projectOntoHGRADBasis().

template<class DataScalar, typename DeviceType>
KOKKOS_INLINE_FUNCTION int Intrepid2::Data< DataScalar, DeviceType >::getVariationModulus ( const int &  d) const
inline

Variation modulus accessor.

Parameters
[in]d- the logical dimension whose variation modulus is requested.
Returns
the variation modulus.

The variation modulus is defined as the number of unique entries in the specified dimension. This is defined as follows:

  • for CONSTANT variation, the variation modulus is 1
  • for MODULAR variation, the variation modulus is exactly the modulus by which the data repeats in the specified dimension
  • for GENERAL variation, the variation modulus is the extent in the specified dimension
  • for BLOCK_PLUS_DIAGONAL, the variation modulus in the first logical dimension of the matrix is the number of nonzeros in the matrix; in the second logical dimension the variation modulus is 1.

Definition at line 1232 of file Intrepid2_Data.hpp.

Referenced by Intrepid2::Data< Intrepid2::Orientation, DeviceType >::allocateMatMatResult(), and Intrepid2::Data< Intrepid2::Orientation, DeviceType >::allocateMatVecResult().

template<class DataScalar, typename DeviceType>
void Intrepid2::Data< DataScalar, DeviceType >::setExtent ( const ordinal_type &  d,
const ordinal_type &  newExtent 
)
inline

sets the logical extent in the specified dimension. If needed, the underlying data container is resized.

Parameters
[in]d- the logical dimension in which the extent is to be changed
[in]newExtent- the new extent
Note
Not supported for dimensions in which the variation type is BLOCK_PLUS_DIAGONAL.
If the variation type is MODULAR, the existing modulus must evenly divide the new extent; the underlying data structure will not be resized in this case.

Definition at line 1958 of file Intrepid2_Data.hpp.

template<class DataScalar, typename DeviceType>
void Intrepid2::Data< DataScalar, DeviceType >::storeMatMat ( const bool  transposeA,
const Data< DataScalar, DeviceType > &  A_MatData,
const bool  transposeB,
const Data< DataScalar, DeviceType > &  B_MatData 
)
inline

Places the result of a matrix-matrix multiply corresponding to the two provided containers into this Data container. This Data container should have been constructed by a call to allocateMatMatResult(), or should match such a container in underlying data extent and variation types.

See Also
allocateMatMat()
Parameters
A_MatData[in] - logically (...,D1,D2)-dimensioned container, where D1,D2 correspond to matrix dimensions.
transposeA[in] - if true, A will be transposed prior to being multiplied by B (or B's transpose).
B_MatData[in] - logically (...,D3,D4)-dimensioned container, where D3,D4 correspond to matrix dimensions.
transposeB[in] - if true, B will be transposed prior to the multiplication by A (or A's transpose).

Definition at line 1819 of file Intrepid2_Data.hpp.

Referenced by Intrepid2::IntegrationTools< DeviceType >::integrate().


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