|
Intrepid2
|
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. | |
| 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. | |
| template<size_t rank, class ... DynRankViewProperties> | |
| Data (const Kokkos::DynRankView< DataScalar, DeviceType, DynRankViewProperties... > &data, Kokkos::Array< int, rank > extents, Kokkos::Array< DataVariationType, rank > variationType, const int blockPlusDiagonalLastNonDiagonal=-1) | |
| Constructor that accepts a DynRankView as an argument. The data belonging to the DynRankView will be copied into a new View of matching dimensions. | |
| 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) | |
| 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) | |
| 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) | |
| 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) | |
| 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) | |
| 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) | |
| 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) | |
| 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). | |
| template<size_t rank> | |
| Data (DataScalar constantValue, Kokkos::Array< int, rank > extents) | |
| constructor for everywhere-constant data | |
| Data () | |
| default constructor (empty data) | |
| 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. | |
| 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. | |
| 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. | |
| void | storeDotProduct (const Data< DataScalar, DeviceType > &A, const Data< DataScalar, DeviceType > &B) |
| Places the result of a contraction along the final dimension of A and B into this data container. | |
| 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, const bool transposeMatrix=false) |
| 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. | |
| 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 > | allocateContractionResult (const Data< DataScalar, DeviceType > &A, const Data< DataScalar, DeviceType > &B, const int &numContractionDims) |
| static Data< DataScalar, DeviceType > | allocateDotProductResult (const Data< DataScalar, DeviceType > &A, const Data< DataScalar, DeviceType > &B) |
| static Data< DataScalar, DeviceType > | allocateMatVecResult (const Data< DataScalar, DeviceType > &matData, const Data< DataScalar, DeviceType > &vecData, const bool transposeMatrix=false) |
Private Types | |
| using | View1 = Kokkos::View< DataScalar *, DeviceType > |
| using | View2 = Kokkos::View< DataScalar **, DeviceType > |
| using | View3 = Kokkos::View< DataScalar ***, DeviceType > |
| using | View4 = Kokkos::View< DataScalar ****, DeviceType > |
| using | View5 = Kokkos::View< DataScalar *****, DeviceType > |
| using | View6 = Kokkos::View< DataScalar ******, DeviceType > |
| using | View7 = Kokkos::View< DataScalar *******, DeviceType > |
| using | return_type = const_reference_type |
Private Member Functions | |
| template<class T > | |
| KOKKOS_INLINE_FUNCTION const T & | get_fixed_view (const std::variant< View1, View2, View3, View4, View5, View6, View7 > &v) const |
| 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_ |
| std::variant< View1, View2, View3, View4, View5, View6, View7 > | underlyingView_ |
| Kokkos::Array< int, 7 > | extents_ |
| Kokkos::Array< DataVariationType, 7 > | variationType_ |
| Kokkos::Array< int, 7 > | variationModulus_ |
| int | blockPlusDiagonalLastNonDiagonal_ = -1 |
| bool | underlyingMatchesLogical_ |
| Kokkos::Array< ordinal_type, 7 > | activeDims_ |
| int | numActiveDims_ |
| ordinal_type | rank_ |
| ScalarView< DataScalar, DeviceType > | zeroView_ |
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 79 of file Intrepid2_Data.hpp.
| using Intrepid2::Data< DataScalar, DeviceType >::all_true = std::is_same<bool_pack<true, v...>, bool_pack<v..., true> > |
Definition at line 1192 of file Intrepid2_Data.hpp.
| using Intrepid2::Data< DataScalar, DeviceType >::const_reference_type = typename ScalarView<const DataScalar,DeviceType>::reference_type |
Definition at line 85 of file Intrepid2_Data.hpp.
| using Intrepid2::Data< DataScalar, DeviceType >::execution_space = typename DeviceType::execution_space |
Definition at line 82 of file Intrepid2_Data.hpp.
| using Intrepid2::Data< DataScalar, DeviceType >::reference_type = typename ScalarView<DataScalar,DeviceType>::reference_type |
Definition at line 84 of file Intrepid2_Data.hpp.
|
private |
Definition at line 108 of file Intrepid2_Data.hpp.
| using Intrepid2::Data< DataScalar, DeviceType >::valid_args = all_true<std::is_integral<IntArgs>{}...> |
Definition at line 1195 of file Intrepid2_Data.hpp.
| using Intrepid2::Data< DataScalar, DeviceType >::value_type = DataScalar |
Definition at line 81 of file Intrepid2_Data.hpp.
|
private |
Definition at line 88 of file Intrepid2_Data.hpp.
|
private |
Definition at line 89 of file Intrepid2_Data.hpp.
|
private |
Definition at line 90 of file Intrepid2_Data.hpp.
|
private |
Definition at line 91 of file Intrepid2_Data.hpp.
|
private |
Definition at line 92 of file Intrepid2_Data.hpp.
|
private |
Definition at line 93 of file Intrepid2_Data.hpp.
|
private |
Definition at line 94 of file Intrepid2_Data.hpp.
|
inline |
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.
Definition at line 528 of file Intrepid2_Data.hpp.
|
inline |
DynRankView constructor. Will copy to a View of appropriate rank.
Definition at line 579 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::allocateAndCopyFromDynRankView(), and Intrepid2::Data< DataScalar, DeviceType >::setActiveDims().
|
inline |
copy-like constructor for differing device type, but same memory space. This does a shallow copy of the underlying view.
Definition at line 590 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView().
|
inline |
copy-like constructor for differing execution spaces. This does a deep_copy of the underlying view.
Definition at line 615 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::copyContainer(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView1(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView2(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView3(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView4(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView5(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView6(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView7(), and Intrepid2::Data< DataScalar, DeviceType >::setActiveDims().
|
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 696 of file Intrepid2_Data.hpp.
|
inline |
Constructor that accepts a DynRankView as an argument. The data belonging to the DynRankView will be copied into a new View of matching dimensions.
Definition at line 706 of file Intrepid2_Data.hpp.
|
inline |
Definition at line 721 of file Intrepid2_Data.hpp.
|
inline |
Definition at line 735 of file Intrepid2_Data.hpp.
|
inline |
Definition at line 749 of file Intrepid2_Data.hpp.
|
inline |
Definition at line 763 of file Intrepid2_Data.hpp.
|
inline |
Definition at line 777 of file Intrepid2_Data.hpp.
|
inline |
Definition at line 791 of file Intrepid2_Data.hpp.
|
inline |
Definition at line 805 of file Intrepid2_Data.hpp.
|
inline |
constructor with run-time rank (requires full-length extents, variationType inputs; those beyond the rank will be ignored).
Definition at line 820 of file Intrepid2_Data.hpp.
|
inline |
constructor for everywhere-constant data
Definition at line 835 of file Intrepid2_Data.hpp.
|
inline |
default constructor (empty data)
Definition at line 849 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::shallowCopy().
|
inline |
allocate an underlying View that matches the provided DynRankView in dimensions, and copy. Called by constructors that accept a DynRankView as argument.
Definition at line 499 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::copyContainer().
Referenced by Intrepid2::Data< DataScalar, DeviceType >::Data().
|
inline |
Constructs a container with extents matching this, with a single-value underlying View, CONSTANT in each dimension.
| value | [in] - the constant value. |
Definition at line 1248 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::getExtents().
Referenced by Intrepid2::CellTools< DeviceType >::setJacobianDetInv().
|
inlinestatic |
Constructs a container suitable for storing the result of a contraction over the final dimensions of the two provided containers. The two containers must have the same logical shape.
| A | [in] - the first data container. |
| B | [in] - the second data container. Must have the same logical shape as A. |
| numContractionDims | [in] - the number of dimensions over which the contraction should take place. |
Definition at line 1458 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::combinedDataDimensionInfo(), Intrepid2::Data< DataScalar, DeviceType >::extent_int(), INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::Data< DataScalar, DeviceType >::rank().
Referenced by Intrepid2::Data< DataScalar, DeviceType >::allocateDotProductResult().
|
inlinestatic |
Constructs a container suitable for storing the result of a contraction over the final dimension of the two provided containers. The two containers must have the same logical shape.
| A | [in] - the first data container. |
| B | [in] - the second data container. Must have the same logical shape as A. |
Definition at line 1478 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::allocateContractionResult().
|
inline |
Returns a DynRankView that matches the underlying Kokkos::View object in value_type, layout, and dimension.
Definition at line 1067 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::extent_int(), and Intrepid2::getMatchingViewWithLabel().
Referenced by Intrepid2::Data< DataScalar, DeviceType >::allocateMatVecResult().
|
inline |
Returns a DynRankView that matches the underlying Kokkos::View object value_type and layout, but with the specified dimensions.
Definition at line 1084 of file Intrepid2_Data.hpp.
References Intrepid2::getMatchingViewWithLabel().
|
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.
| A | [in] - the first data container. |
| B | [in] - the second data container. Must have the same logical shape as A. |
Definition at line 1258 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::combinedDataDimensionInfo(), Intrepid2::Data< DataScalar, DeviceType >::extent_int(), INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::Data< DataScalar, DeviceType >::rank().
Referenced by Intrepid2::IntegrationTools< DeviceType >::integrate(), Intrepid2::DataTools::multiplyByCPWeights(), and Intrepid2::TransformedBasisValues< Scalar, DeviceType >::multiplyByPointwiseWeights().
|
inlinestatic |
Constructs a container suitable for storing the result of a matrix-vector multiply corresponding to the two provided containers.
| 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 1278 of file Intrepid2_Data.hpp.
References Intrepid2::BLOCK_PLUS_DIAGONAL, Intrepid2::Data< DataScalar, DeviceType >::blockPlusDiagonalLastNonDiagonal(), Intrepid2::CONSTANT, Intrepid2::Data< DataScalar, DeviceType >::extent_int(), Intrepid2::GENERAL, Intrepid2::getMatchingViewWithLabel(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView(), Intrepid2::Data< DataScalar, DeviceType >::getVariationModulus(), Intrepid2::Data< DataScalar, DeviceType >::getVariationTypes(), INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, Intrepid2::MODULAR, and Intrepid2::Data< DataScalar, DeviceType >::rank().
Referenced by Intrepid2::DataTools::transposeMatrix().
|
inlinestatic |
Constructs a container suitable for storing the result of a matrix-vector multiply corresponding to the two provided containers.
Definition at line 1485 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::allocateDynRankViewMatchingUnderlying(), Intrepid2::BLOCK_PLUS_DIAGONAL, Intrepid2::CONSTANT, Intrepid2::Data< DataScalar, DeviceType >::extent_int(), Intrepid2::GENERAL, Intrepid2::Data< DataScalar, DeviceType >::getVariationModulus(), Intrepid2::Data< DataScalar, DeviceType >::getVariationTypes(), INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, Intrepid2::MODULAR, and Intrepid2::Data< DataScalar, DeviceType >::rank().
Referenced by Intrepid2::ProjectedGeometry< spaceDim, PointScalar, DeviceType >::projectOntoHGRADBasis().
|
inline |
applies the specified unary operator to each entry
Definition at line 259 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::getDataExtent().
|
inlinestaticprivate |
//! Returns flattened index of the specified (i,j) matrix entry, assuming that i,j ≤ lastNondiagonal. Only applicable for BLOCK_PLUS_DIAGONAL DataVariationType.
Definition at line 121 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::getWritableEntryWithPassThroughOption().
|
inlinestaticprivate |
Returns flattened index of the specified (i,i) matrix entry, assuming that i > lastNondiagonal. Only applicable for BLOCK_PLUS_DIAGONAL DataVariationType.
Definition at line 128 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::getWritableEntryWithPassThroughOption().
|
inline |
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.
Definition at line 858 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::allocateMatMatResult(), Intrepid2::CellTools< DeviceType >::setJacobianDet(), and Intrepid2::CellTools< DeviceType >::setJacobianInv().
|
inlinestaticprivate |
Returns the number of non-diagonal entries based on the last non-diagonal. Only applicable for BLOCK_PLUS_DIAGONAL DataVariationType.
Definition at line 114 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::getWritableEntryWithPassThroughOption(), and Intrepid2::Data< DataScalar, DeviceType >::setActiveDims().
|
inline |
Copies 0.0 to the underlying View.
Definition at line 1100 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::IntegrationTools< DeviceType >::integrate().
|
inline |
Returns (DataVariationType, data extent) in the specified dimension for a Data container that combines (through multiplication, say, or addition) this container with otherData.
Definition at line 890 of file Intrepid2_Data.hpp.
References Intrepid2::combinedDimensionInfo(), and Intrepid2::Data< DataScalar, DeviceType >::getDimensionInfo().
Referenced by Intrepid2::Data< DataScalar, DeviceType >::allocateContractionResult(), and Intrepid2::Data< DataScalar, DeviceType >::allocateInPlaceCombinationResult().
|
inlinestatic |
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).
Definition at line 484 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::allocateAndCopyFromDynRankView(), Intrepid2::Data< DataScalar, DeviceType >::copyDataFromDynRankViewMatchingUnderlying(), and Intrepid2::Data< DataScalar, DeviceType >::Data().
|
inline |
Copies from the provided DynRankView into the underlying Kokkos::View container storing the unique data.
Definition at line 1116 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::copyContainer().
|
inline |
returns an MDRangePolicy over the underlying data extents (but with the logical shape).
Definition at line 1619 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::getDataExtent(), and Intrepid2::Data< DataScalar, DeviceType >::rank().
|
inline |
returns an MDRangePolicy over the first six underlying data extents (but with the logical shape).
Definition at line 1637 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::getDataExtent().
|
inline |
Definition at line 1655 of file Intrepid2_Data.hpp.
|
inlineconstexpr |
Definition at line 1224 of file Intrepid2_Data.hpp.
|
inline |
Returns the logical extent in the specified dimension.
Definition at line 1216 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::allocateCellMeasure(), Intrepid2::Data< DataScalar, DeviceType >::allocateContractionResult(), Intrepid2::Data< DataScalar, DeviceType >::allocateDynRankViewMatchingUnderlying(), Intrepid2::Data< DataScalar, DeviceType >::allocateInPlaceCombinationResult(), Intrepid2::Data< DataScalar, DeviceType >::allocateMatMatResult(), Intrepid2::Data< DataScalar, DeviceType >::allocateMatVecResult(), Intrepid2::Impl::F_Integrate< Scalar, DeviceType, integralViewRank >::approximateFlopCountPerCell(), Intrepid2::Impl::F_IntegratePointValueCache< Scalar, DeviceType, integralViewRank >::approximateFlopCountPerCell(), Intrepid2::Data< DataScalar, DeviceType >::getDimensionInfo(), Intrepid2::IntegrationTools< DeviceType >::integrate(), Intrepid2::DataTools::multiplyByCPWeights(), Intrepid2::DataTools::multiplyByCPWeights(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::operator()(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::operator()(), Intrepid2::Impl::F_Integrate< Scalar, DeviceType, integralViewRank >::runSpecialized3(), Intrepid2::Impl::F_IntegratePointValueCache< Scalar, DeviceType, integralViewRank >::runSpecialized3(), Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::setJacobianDataPrivate(), Intrepid2::CellTools< DeviceType >::setJacobianDet(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::spaceDim(), Intrepid2::Data< DataScalar, DeviceType >::storeDotProduct(), Intrepid2::Data< DataScalar, DeviceType >::storeMatMat(), Intrepid2::Data< DataScalar, DeviceType >::storeMatVec(), Intrepid2::Impl::F_IntegratePointValueCache< Scalar, DeviceType, integralViewRank >::team_shmem_size(), Intrepid2::Impl::F_Integrate< Scalar, DeviceType, integralViewRank >::team_shmem_size(), and Intrepid2::TensorArgumentIterator::TensorArgumentIterator().
|
inlineprivate |
Definition at line 134 of file Intrepid2_Data.hpp.
|
inline |
returns the true extent of the data corresponding to the logical dimension provided; if the data does not vary in that dimension, returns 1
Definition at line 1133 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingViewExtent().
Referenced by Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::allocateCellMeasure(), Intrepid2::Data< DataScalar, DeviceType >::applyOperator(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::cellDataExtent(), Intrepid2::Data< DataScalar, DeviceType >::dataExtentRangePolicy(), Intrepid2::Data< DataScalar, DeviceType >::getDimensionInfo(), Intrepid2::IntegrationTools< DeviceType >::integrate(), Intrepid2::Data< DataScalar, DeviceType >::storeDotProduct(), Intrepid2::DataCombiner< DataScalar, DeviceType, BinaryOperator >::storeInPlaceCombination(), Intrepid2::Data< DataScalar, DeviceType >::storeMatMat(), and Intrepid2::Data< DataScalar, DeviceType >::storeMatVec().
|
inline |
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).
Definition at line 872 of file Intrepid2_Data.hpp.
References Intrepid2::BLOCK_PLUS_DIAGONAL, Intrepid2::Data< DataScalar, DeviceType >::extent_int(), and Intrepid2::Data< DataScalar, DeviceType >::getDataExtent().
Referenced by Intrepid2::IntegrationTools< DeviceType >::allocateIntegralData(), and Intrepid2::Data< DataScalar, DeviceType >::combinedDataDimensionInfo().
|
inline |
Returns a (read-only) value corresponding to the specified logical data location.
Definition at line 1184 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::getEntryWithPassThroughOption().
Referenced by Intrepid2::Data< DataScalar, DeviceType >::operator()().
|
inline |
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.
Definition at line 1176 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::getWritableEntryWithPassThroughOption().
Referenced by Intrepid2::Data< DataScalar, DeviceType >::getEntry().
|
inline |
Returns an array containing the logical extents in each dimension.
Definition at line 865 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::allocateCellMeasure(), Intrepid2::Data< DataScalar, DeviceType >::allocateConstantData(), Intrepid2::CellTools< DeviceType >::allocateJacobianDet(), Intrepid2::CellTools< DeviceType >::allocateJacobianInv(), Intrepid2::IntegrationTools< DeviceType >::integrate(), Intrepid2::DataTools::multiplyByCPWeights(), Intrepid2::DataTools::multiplyByCPWeights(), and Intrepid2::DataTools::transposeMatrix().
|
inline |
Returns the underlying view. Throws an exception if the underlying view is not rank 1.
Definition at line 902 of file Intrepid2_Data.hpp.
References INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::Data< DataScalar, DeviceType >::rank().
Referenced by Intrepid2::Data< DataScalar, DeviceType >::allocateMatMatResult(), Intrepid2::Data< DataScalar, DeviceType >::Data(), Intrepid2::Cubature< DeviceType, pointValueType, weightValueType >::getCubature(), Intrepid2::Basis_TensorBasis< BasisBaseClass >::getValues(), Intrepid2::IntegrationTools< DeviceType >::integrate(), and Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::setJacobianDataPrivate().
|
inline |
Returns the underlying view. Throws an exception if the underlying view is not rank 2.
Definition at line 914 of file Intrepid2_Data.hpp.
References INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::Data< DataScalar, DeviceType >::rank().
|
inline |
Returns the underlying view. Throws an exception if the underlying view is not rank 3.
Definition at line 926 of file Intrepid2_Data.hpp.
References INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::Data< DataScalar, DeviceType >::rank().
|
inline |
Returns the underlying view. Throws an exception if the underlying view is not rank 4.
Definition at line 938 of file Intrepid2_Data.hpp.
References INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::Data< DataScalar, DeviceType >::rank().
|
inline |
Returns the underlying view. Throws an exception if the underlying view is not rank 5.
Definition at line 950 of file Intrepid2_Data.hpp.
References INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::Data< DataScalar, DeviceType >::rank().
|
inline |
Returns the underlying view. Throws an exception if the underlying view is not rank 6.
Definition at line 962 of file Intrepid2_Data.hpp.
References INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::Data< DataScalar, DeviceType >::rank().
|
inline |
Returns the underlying view. Throws an exception if the underlying view is not rank 7.
Definition at line 974 of file Intrepid2_Data.hpp.
References INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::Data< DataScalar, DeviceType >::rank().
|
inline |
Returns a DynRankView constructed atop the same underlying data as the fixed-rank Kokkos::View used internally.
Definition at line 1032 of file Intrepid2_Data.hpp.
|
inline |
returns the View that stores the unique data. For rank-1 underlying containers.
Definition at line 984 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::CellTools< DeviceType >::allocateJacobianDet(), Intrepid2::CellTools< DeviceType >::allocateJacobianInv(), Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::computeCellMeasure(), Intrepid2::Data< DataScalar, DeviceType >::Data(), Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::setJacobianDataPrivate(), Intrepid2::CellTools< DeviceType >::setJacobianDet(), and Intrepid2::CellTools< DeviceType >::setJacobianInv().
|
inline |
returns the View that stores the unique data. For rank-2 underlying containers.
Definition at line 991 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::CellTools< DeviceType >::allocateJacobianInv(), Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::computeCellMeasure(), Intrepid2::Data< DataScalar, DeviceType >::Data(), Intrepid2::IntegrationTools< DeviceType >::integrate(), Intrepid2::Impl::F_Integrate< Scalar, DeviceType, integralViewRank >::runSpecialized3(), Intrepid2::CellTools< DeviceType >::setJacobianDet(), and Intrepid2::CellTools< DeviceType >::setJacobianInv().
|
inline |
returns the View that stores the unique data. For rank-3 underlying containers.
Definition at line 998 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::CellTools< DeviceType >::allocateJacobianDet(), Intrepid2::CellTools< DeviceType >::allocateJacobianInv(), Intrepid2::Data< DataScalar, DeviceType >::Data(), Intrepid2::IntegrationTools< DeviceType >::integrate(), Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::setJacobianDataPrivate(), Intrepid2::CellTools< DeviceType >::setJacobianDet(), and Intrepid2::CellTools< DeviceType >::setJacobianInv().
|
inline |
returns the View that stores the unique data. For rank-4 underlying containers.
Definition at line 1005 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::CellTools< DeviceType >::allocateJacobianDet(), Intrepid2::CellTools< DeviceType >::allocateJacobianInv(), Intrepid2::Data< DataScalar, DeviceType >::Data(), Intrepid2::Impl::F_Integrate< Scalar, DeviceType, integralViewRank >::runSpecialized3(), Intrepid2::CellTools< DeviceType >::setJacobianDet(), Intrepid2::CellTools< DeviceType >::setJacobianInv(), and Intrepid2::Data< DataScalar, DeviceType >::storeMatMat().
|
inline |
returns the View that stores the unique data. For rank-5 underlying containers.
Definition at line 1012 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::Data().
|
inline |
returns the View that stores the unique data. For rank-6 underlying containers.
Definition at line 1019 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::Data().
|
inline |
returns the View that stores the unique data. For rank-7 underlying containers.
Definition at line 1026 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::Data().
|
inlineprivate |
Returns the extent of the underlying view in the specified dimension.
Definition at line 140 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::getDataExtent(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingViewSize(), Intrepid2::Data< DataScalar, DeviceType >::setActiveDims(), and Intrepid2::Data< DataScalar, DeviceType >::setExtent().
|
inline |
returns the rank of the View that stores the unique data
Definition at line 1049 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::CellTools< DeviceType >::allocateJacobianInv(), Intrepid2::IntegrationTools< DeviceType >::integrate(), Intrepid2::CellTools< DeviceType >::setJacobianDet(), Intrepid2::CellTools< DeviceType >::setJacobianInv(), and Intrepid2::DataCombiner< DataScalar, DeviceType, BinaryOperator >::storeInPlaceCombination().
|
inline |
returns the number of entries in the View that stores the unique data
Definition at line 1056 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingViewExtent().
Referenced by Intrepid2::IntegrationTools< DeviceType >::integrate(), and Intrepid2::DataCombiner< DataScalar, DeviceType, BinaryOperator >::storeInPlaceCombination().
|
inline |
Variation modulus accessor.
| [in] | d | - the logical dimension whose variation modulus is requested. |
The variation modulus is defined as the number of unique entries in the specified dimension. This is defined as follows:
Definition at line 1161 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::allocateMatMatResult(), and Intrepid2::Data< DataScalar, DeviceType >::allocateMatVecResult().
|
inline |
Returns an array with the variation types in each logical dimension.
Definition at line 1168 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::allocateCellMeasure(), Intrepid2::CellTools< DeviceType >::allocateJacobianDet(), Intrepid2::CellTools< DeviceType >::allocateJacobianInv(), Intrepid2::Data< DataScalar, DeviceType >::allocateMatMatResult(), Intrepid2::Data< DataScalar, DeviceType >::allocateMatVecResult(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::cellVariationType(), Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::computeCellMeasure(), Intrepid2::IntegrationTools< DeviceType >::integrate(), Intrepid2::DataTools::multiplyByCPWeights(), Intrepid2::DataTools::multiplyByCPWeights(), Intrepid2::CellTools< DeviceType >::setJacobianDet(), Intrepid2::CellTools< DeviceType >::setJacobianInv(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::TransformedBasisValues(), and Intrepid2::DataTools::transposeMatrix().
|
inline |
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.
Definition at line 477 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::getWritableEntryWithPassThroughOption().
Referenced by Intrepid2::Data< DataScalar, DeviceType >::storeDotProduct(), Intrepid2::Data< DataScalar, DeviceType >::storeMatMat(), and Intrepid2::Data< DataScalar, DeviceType >::storeMatVec().
|
inline |
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.
Definition at line 364 of file Intrepid2_Data.hpp.
References Intrepid2::BLOCK_PLUS_DIAGONAL, Intrepid2::Data< DataScalar, DeviceType >::blockPlusDiagonalBlockEntryIndex(), Intrepid2::Data< DataScalar, DeviceType >::blockPlusDiagonalDiagonalEntryIndex(), Intrepid2::Data< DataScalar, DeviceType >::blockPlusDiagonalNumNondiagonalEntries(), Intrepid2::CONSTANT, Intrepid2::GENERAL, INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::MODULAR.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::getEntryWithPassThroughOption(), and Intrepid2::Data< DataScalar, DeviceType >::getWritableEntry().
|
inline |
returns true for containers that have two dimensions marked as BLOCK_PLUS_DIAGONAL for which the non-diagonal block is empty or size 1.
Definition at line 1229 of file Intrepid2_Data.hpp.
References Intrepid2::BLOCK_PLUS_DIAGONAL, and INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE.
|
inlineconstexpr |
returns true for containers that have data; false for those that don't (namely, those that have been constructed by the default constructor).
Definition at line 1968 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::IntegrationTools< DeviceType >::allocateIntegralData(), Intrepid2::IntegrationTools< DeviceType >::integrate(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::multiplyByPointwiseWeights(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::operator()(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::operator()(), Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::setJacobianDataPrivate(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::spaceDim(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::transformWeight(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::transformWeight(), and Intrepid2::TransformedBasisValues< Scalar, DeviceType >::transformWeight().
|
inline |
Returns a value corresponding to the specified logical data location.
Definition at line 1210 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::getEntry().
|
inline |
Returns the logical rank of the Data container.
Definition at line 1975 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::allocateContractionResult(), Intrepid2::Data< DataScalar, DeviceType >::allocateInPlaceCombinationResult(), Intrepid2::Data< DataScalar, DeviceType >::allocateMatMatResult(), Intrepid2::Data< DataScalar, DeviceType >::allocateMatVecResult(), Intrepid2::Data< DataScalar, DeviceType >::dataExtentRangePolicy(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView(), Intrepid2::IntegrationTools< DeviceType >::integrate(), Intrepid2::DataTools::multiplyByCPWeights(), Intrepid2::DataTools::multiplyByCPWeights(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::multiplyByPointwiseWeights(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::operator()(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::operator()(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::rank(), Intrepid2::Impl::F_Integrate< Scalar, DeviceType, integralViewRank >::runSpecialized3(), Intrepid2::Impl::F_IntegratePointValueCache< Scalar, DeviceType, integralViewRank >::runSpecialized3(), Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::setJacobianDataPrivate(), Intrepid2::Data< DataScalar, DeviceType >::shallowCopy(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::spaceDim(), Intrepid2::Data< DataScalar, DeviceType >::storeDotProduct(), Intrepid2::Data< DataScalar, DeviceType >::storeMatMat(), Intrepid2::Data< DataScalar, DeviceType >::storeMatVec(), Intrepid2::TransformedBasisValues< Scalar, DeviceType >::TransformedBasisValues(), and Intrepid2::DataTools::transposeMatrix().
|
inlineprivate |
class initialization method. Called by constructors.
Definition at line 156 of file Intrepid2_Data.hpp.
References Intrepid2::BLOCK_PLUS_DIAGONAL, Intrepid2::Data< DataScalar, DeviceType >::blockPlusDiagonalNumNondiagonalEntries(), Intrepid2::GENERAL, Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingViewExtent(), and Intrepid2::MODULAR.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::Data(), and Intrepid2::Data< DataScalar, DeviceType >::Data().
|
inline |
sets the logical extent in the specified dimension. If needed, the underlying data container is resized.
| [in] | d | - the logical dimension in which the extent is to be changed |
| [in] | newExtent | - the new extent |
Definition at line 1986 of file Intrepid2_Data.hpp.
References Intrepid2::BLOCK_PLUS_DIAGONAL, Intrepid2::GENERAL, Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingViewExtent(), and Intrepid2::MODULAR.
|
inline |
Creates a new Data object with the same underlying view, but with the specified logical rank, extents, and variation types.
Definition at line 1663 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::Data(), and Intrepid2::Data< DataScalar, DeviceType >::rank().
Referenced by Intrepid2::IntegrationTools< DeviceType >::integrate(), Intrepid2::DataTools::multiplyByCPWeights(), and Intrepid2::DataTools::multiplyByCPWeights().
|
inline |
Places the result of a contraction along the final dimension of A and B into this data container.
Definition at line 1680 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::extent_int(), Intrepid2::Data< DataScalar, DeviceType >::getDataExtent(), Intrepid2::Data< DataScalar, DeviceType >::getWritableEntry(), INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::Data< DataScalar, DeviceType >::rank().
| void Intrepid2::Data< DataScalar, DeviceType >::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.
Definition at line 31 of file Intrepid2_DataDef.hpp.
Referenced by Intrepid2::Data< DataScalar, DeviceType >::storeInPlaceDifference(), Intrepid2::Data< DataScalar, DeviceType >::storeInPlaceProduct(), Intrepid2::Data< DataScalar, DeviceType >::storeInPlaceQuotient(), and Intrepid2::Data< DataScalar, DeviceType >::storeInPlaceSum().
|
inline |
stores the in-place (entrywise) difference, A .- B, into this container.
Definition at line 1756 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::storeInPlaceCombination().
|
inline |
stores the in-place (entrywise) product, A .* B, into this container.
Definition at line 1749 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::storeInPlaceCombination().
Referenced by Intrepid2::IntegrationTools< DeviceType >::integrate(), and Intrepid2::DataTools::multiplyByCPWeights().
|
inline |
stores the in-place (entrywise) quotient, A ./ B, into this container.
Definition at line 1763 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::storeInPlaceCombination().
Referenced by Intrepid2::CellTools< DeviceType >::setJacobianDetInv().
|
inline |
stores the in-place (entrywise) sum, A .+ B, into this container.
Definition at line 1742 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::storeInPlaceCombination().
|
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.
| 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 1847 of file Intrepid2_Data.hpp.
References Intrepid2::BLOCK_PLUS_DIAGONAL, Intrepid2::Data< DataScalar, DeviceType >::extent_int(), Intrepid2::Data< DataScalar, DeviceType >::getDataExtent(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView4(), Intrepid2::Data< DataScalar, DeviceType >::getWritableEntry(), INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::Data< DataScalar, DeviceType >::rank().
Referenced by Intrepid2::IntegrationTools< DeviceType >::integrate().
|
inline |
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.
Definition at line 1770 of file Intrepid2_Data.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::extent_int(), Intrepid2::Data< DataScalar, DeviceType >::getDataExtent(), Intrepid2::Data< DataScalar, DeviceType >::getWritableEntry(), INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::Data< DataScalar, DeviceType >::rank().
|
inline |
Returns true if the underlying container has exactly the same rank and extents as the logical container.
Definition at line 2039 of file Intrepid2_Data.hpp.
Referenced by Intrepid2::Impl::F_Integrate< Scalar, DeviceType, integralViewRank >::runSpecialized3(), and Intrepid2::DataCombiner< DataScalar, DeviceType, BinaryOperator >::storeInPlaceCombination().
|
private |
Definition at line 102 of file Intrepid2_Data.hpp.
|
private |
Definition at line 99 of file Intrepid2_Data.hpp.
|
private |
Definition at line 87 of file Intrepid2_Data.hpp.
|
private |
Definition at line 96 of file Intrepid2_Data.hpp.
|
private |
Definition at line 103 of file Intrepid2_Data.hpp.
|
private |
Definition at line 105 of file Intrepid2_Data.hpp.
|
private |
Definition at line 101 of file Intrepid2_Data.hpp.
|
private |
Definition at line 95 of file Intrepid2_Data.hpp.
|
private |
Definition at line 98 of file Intrepid2_Data.hpp.
|
private |
Definition at line 97 of file Intrepid2_Data.hpp.
|
private |
Definition at line 110 of file Intrepid2_Data.hpp.