Intrepid2
Intrepid2_TensorData.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Intrepid2 Package
4 //
5 // Copyright 2007 NTESS and the Intrepid2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
16 #ifndef Intrepid2_TensorData_h
17 #define Intrepid2_TensorData_h
18 
19 #include "Intrepid2_Data.hpp"
20 
21 #include "Intrepid2_ScalarView.hpp"
22 #include "Intrepid2_Types.hpp"
23 
24 namespace Intrepid2
25 {
29  template<class Scalar, typename DeviceType>
30  class TensorData {
31  public:
32  using value_type = Scalar;
33  using execution_space = typename DeviceType::execution_space;
34  protected:
35  Kokkos::Array< Data<Scalar,DeviceType>, Parameters::MaxTensorComponents> tensorComponents_;
36  Kokkos::Array<ordinal_type, 7> extents_;
37  Kokkos::Array<Kokkos::Array<ordinal_type, Parameters::MaxTensorComponents>, 7> entryModulus_;
38  ordinal_type rank_;
39  bool separateFirstComponent_ = false; // true supported only for rank 1 components; uses a rank-2 operator() in that case, where the first argument corresponds to the index in the first component, while the second argument corresponds to the tensorially-multiplied remaining components
40  ordinal_type numTensorComponents_ = 0;
41 
45  void initialize()
46  {
47  ordinal_type maxComponentRank = -1;
48  for (ordinal_type r=0; r<numTensorComponents_; r++)
49  {
50  const ordinal_type componentRank = tensorComponents_[r].rank();
51  maxComponentRank = (maxComponentRank > componentRank) ? maxComponentRank : componentRank;
52  }
53  INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE(separateFirstComponent_ && (maxComponentRank != 1), std::invalid_argument, "separateFirstComponent = true only supported if all components have rank 1");
54  ordinal_type rd_start; // where to begin the extents_ and entryModulus_ loops below
55  if ((maxComponentRank == 1) && separateFirstComponent_)
56  {
57  rank_ = 2;
58  rd_start = 1;
59  extents_[0] = tensorComponents_[0].extent_int(0);
60  entryModulus_[0][0] = extents_[0]; // should not be used
61  }
62  else
63  {
64  rank_ = maxComponentRank;
65  rd_start = 0;
66  }
67 
68  for (ordinal_type d=rd_start; d<7; d++)
69  {
70  extents_[d] = 1;
71  for (ordinal_type r=rd_start; r<numTensorComponents_; r++)
72  {
73  extents_[d] *= tensorComponents_[r].extent_int(d-rd_start);
74  }
75  ordinal_type entryModulus = extents_[d];
76  for (ordinal_type r=rd_start; r<numTensorComponents_; r++)
77  {
78  entryModulus /= tensorComponents_[r].extent_int(d-rd_start);
79  entryModulus_[d][r] = entryModulus;
80  }
81  }
82  }
83  public:
93  template<size_t numTensorComponents>
94  TensorData(Kokkos::Array< Data<Scalar,DeviceType>, numTensorComponents> tensorComponents, bool separateFirstComponent = false)
95  :
96  separateFirstComponent_(separateFirstComponent),
97  numTensorComponents_(numTensorComponents)
98  {
99  for (size_t r=0; r< numTensorComponents; r++)
100  {
101  tensorComponents_[r] = tensorComponents[r];
102  }
103 
104  initialize();
105  }
106 
116  TensorData(std::vector< Data<Scalar,DeviceType> > tensorComponents, bool separateFirstComponent = false)
117  :
118  separateFirstComponent_(separateFirstComponent),
119  numTensorComponents_(tensorComponents.size())
120  {
121  for (ordinal_type r=0; r<numTensorComponents_; r++)
122  {
123  tensorComponents_[r] = tensorComponents[r];
124  }
125 
126  initialize();
127  }
128 
137  TensorData(const TensorData &first, const TensorData &second, bool separateFirstComponent = false)
138  :
139  separateFirstComponent_(separateFirstComponent),
140  numTensorComponents_(first.numTensorComponents() + second.numTensorComponents())
141  {
142  ordinal_type r = 0;
143  for (ordinal_type r1=0; r1<first.numTensorComponents(); r1++, r++)
144  {
145  tensorComponents_[r] = first.getTensorComponent(r1);
146  }
147  for (ordinal_type r2=0; r2<second.numTensorComponents(); r2++, r++)
148  {
149  tensorComponents_[r] = second.getTensorComponent(r2);
150  }
151 
152  initialize();
153  }
154 
162  :
163  TensorData(Kokkos::Array< Data<Scalar,DeviceType>, 1>({tensorComponent}), false)
164  {}
165 
172  :
173  extents_({0,0,0,0,0,0,0}),
174  rank_(0)
175  {}
176 
184  TensorData(TensorData otherTensorData, std::vector<int> whichComps)
185  :
186  numTensorComponents_(whichComps.size())
187  {
188  int r = 0;
189  for (const auto & componentOrdinal : whichComps)
190  {
191  tensorComponents_[r++] = otherTensorData.getTensorComponent(componentOrdinal);
192  }
193 
194  initialize();
195  }
196 
198  template<typename OtherDeviceType, class = typename std::enable_if< std::is_same<typename DeviceType::memory_space, typename OtherDeviceType::memory_space>::value>::type,
199  class = typename std::enable_if<!std::is_same<DeviceType,OtherDeviceType>::value>::type>
200  TensorData(const TensorData<Scalar,OtherDeviceType> &tensorData)
201  {
202  if (tensorData.isValid())
203  {
204  numTensorComponents_ = tensorData.numTensorComponents();
205  for (ordinal_type r=0; r<numTensorComponents_; r++)
206  {
207  Data<Scalar,OtherDeviceType> otherTensorComponent = tensorData.getTensorComponent(r);
208  tensorComponents_[r] = Data<Scalar,DeviceType>(otherTensorComponent);
209  }
210  initialize();
211  }
212  else
213  {
214  extents_ = Kokkos::Array<ordinal_type,7>{0,0,0,0,0,0,0};
215  rank_ = 0;
216  }
217  }
218 
222  template<typename OtherDeviceType, class = typename std::enable_if<!std::is_same<typename DeviceType::memory_space, typename OtherDeviceType::memory_space>::value>::type>
224  {
225  if (tensorData.isValid())
226  {
227  numTensorComponents_ = tensorData.numTensorComponents();
228  for (ordinal_type r=0; r<numTensorComponents_; r++)
229  {
230  Data<Scalar,OtherDeviceType> otherTensorComponent = tensorData.getTensorComponent(r);
231  tensorComponents_[r] = Data<Scalar,DeviceType>(otherTensorComponent);
232  }
233  initialize();
234  }
235  else
236  {
237  extents_ = Kokkos::Array<ordinal_type,7>{0,0,0,0,0,0,0};
238  rank_ = 0;
239  }
240  }
241 
247  KOKKOS_INLINE_FUNCTION
248  const Data<Scalar,DeviceType> & getTensorComponent(const ordinal_type &r) const
249  {
250  return tensorComponents_[r];
251  }
252 
258  template <typename iType0>
259  KOKKOS_INLINE_FUNCTION typename std::enable_if<std::is_integral<iType0>::value, Scalar>::type
260  operator()(const iType0& tensorEntryIndex) const {
261 #ifdef HAVE_INTREPID2_DEBUG
262  INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE(rank_ != 1, std::invalid_argument, "This method is only valid for rank 1 containers.");
263 #endif
264  Scalar value = 1.0;
265  iType0 remainingEntryOrdinal = tensorEntryIndex;
266  for (ordinal_type r=0; r<numTensorComponents_; r++)
267  {
268  const ordinal_type componentEntryCount = tensorComponents_[r].extent_int(0);
269  const ordinal_type componentEntryOrdinal = remainingEntryOrdinal % componentEntryCount;
270  remainingEntryOrdinal /= componentEntryCount;
271 
272  value *= tensorComponents_[r](componentEntryOrdinal);
273  }
274 
275  return value;
276  }
277 
283  template <typename iType0, ordinal_type numTensorComponents>
284  KOKKOS_INLINE_FUNCTION typename std::enable_if<std::is_integral<iType0>::value, Scalar>::type
285  operator()(const Kokkos::Array<iType0,numTensorComponents>& entryComponents) const {
286 #ifdef HAVE_INTREPID2_DEBUG
287  INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE(rank_ != 1, std::invalid_argument, "This method is only valid for rank 1 containers.");
288  INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE(numTensorComponents_ != numTensorComponents, std::invalid_argument, "Tensorial component count mismatch");
289 #endif
290  Scalar value = 1.0;
291  for (ordinal_type r=0; r<numTensorComponents; r++)
292  {
293  value *= tensorComponents_[r](entryComponents[r]);
294  }
295  return value;
296  }
297 
308  template <typename iType0, typename iType1>
309  KOKKOS_INLINE_FUNCTION typename std::enable_if<
310  (std::is_integral<iType0>::value && std::is_integral<iType1>::value),
311  Scalar>::type
312  operator()(const iType0& tensorEntryIndex0, const iType1& tensorEntryIndex1) const {
313 #ifdef HAVE_INTREPID2_DEBUG
314  INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE(rank_ != 2, std::invalid_argument, "This method is only valid for rank 2 containers.");
315 #endif
316 
317  if (numTensorComponents_ == 1)
318  {
319  return tensorComponents_[0](tensorEntryIndex0,tensorEntryIndex1);
320  }
321 
322  if (!separateFirstComponent_)
323  {
324  Scalar value = 1.0;
325  iType0 remainingEntryOrdinal0 = tensorEntryIndex0;
326  iType1 remainingEntryOrdinal1 = tensorEntryIndex1;
327  for (ordinal_type r=0; r<numTensorComponents_; r++)
328  {
329  auto & component = tensorComponents_[r];
330  const ordinal_type componentEntryCount0 = component.extent_int(0);
331  const ordinal_type componentEntryCount1 = component.extent_int(1);
332  const iType0 componentEntryOrdinal0 = remainingEntryOrdinal0 % componentEntryCount0;
333  const iType1 componentEntryOrdinal1 = remainingEntryOrdinal1 % componentEntryCount1;
334  remainingEntryOrdinal0 /= componentEntryCount0;
335  remainingEntryOrdinal1 /= componentEntryCount1;
336 
337  const ordinal_type componentRank = component.rank();
338 
339  if (componentRank == 2)
340  {
341  value *= component(componentEntryOrdinal0,componentEntryOrdinal1);
342  }
343  else if (componentRank == 1)
344  {
345  value *= component(componentEntryOrdinal0);
346  }
347  else
348  {
349  INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE(true, std::logic_error, "unsupported component rank encountered");
350  }
351  }
352 
353  return value;
354  }
355  else
356  {
357  Scalar value = tensorComponents_[0](tensorEntryIndex0);
358  iType0 remainingEntryOrdinal = tensorEntryIndex1;
359  for (ordinal_type r=1; r<numTensorComponents_; r++)
360  {
361  const ordinal_type componentEntryCount = tensorComponents_[r].extent_int(0);
362  const ordinal_type componentEntryOrdinal = remainingEntryOrdinal % componentEntryCount;
363  remainingEntryOrdinal /= componentEntryCount;
364 
365  value *= tensorComponents_[r](componentEntryOrdinal);
366  }
367  return value;
368  }
369  }
370 
372  KOKKOS_INLINE_FUNCTION
373  ordinal_type getTensorComponentIndex(const ordinal_type &tensorComponent, const ordinal_type &dim, const ordinal_type &enumerationIndex) const
374  {
375  ordinal_type remainingEntryOrdinal = enumerationIndex;
376  for (ordinal_type r=0; r<tensorComponent; r++)
377  {
378  const auto & component = tensorComponents_[r];
379  const ordinal_type & componentEntryCount = component.extent_int(dim);
380 
381  remainingEntryOrdinal /= componentEntryCount;
382  }
383  return remainingEntryOrdinal % tensorComponents_[tensorComponent].extent_int(dim);
384  }
385 
395  template <typename iType0, typename iType1, typename iType2>
396  KOKKOS_INLINE_FUNCTION typename std::enable_if<
397  (std::is_integral<iType0>::value && std::is_integral<iType1>::value && std::is_integral<iType2>::value),
398  Scalar>::type
399  operator()(const iType0& tensorEntryIndex0, const iType1& tensorEntryIndex1, const iType2& tensorEntryIndex2) const
400  {
401 #ifdef HAVE_INTREPID2_DEBUG
402  INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE(rank_ != 3, std::invalid_argument, "This method is only valid for rank 3 containers.");
403  INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE(separateFirstComponent_, std::logic_error, "This method should never be called when separateFirstComponent is true");
404 #endif
405 
406  Scalar value = 1.0;
407  Kokkos::Array<ordinal_type,3> remainingEntryOrdinal {tensorEntryIndex0, tensorEntryIndex1, tensorEntryIndex2};
408  for (ordinal_type r=0; r<numTensorComponents_; r++)
409  {
410  auto & component = tensorComponents_[r];
411  const ordinal_type componentEntryCount0 = component.extent_int(0);
412  const ordinal_type componentEntryCount1 = component.extent_int(1);
413  const ordinal_type componentEntryCount2 = component.extent_int(2);
414  const ordinal_type componentEntryOrdinal0 = remainingEntryOrdinal[0] % componentEntryCount0;
415  const ordinal_type componentEntryOrdinal1 = remainingEntryOrdinal[1] % componentEntryCount1;
416  const ordinal_type componentEntryOrdinal2 = remainingEntryOrdinal[2] % componentEntryCount2;
417  remainingEntryOrdinal[0] /= componentEntryCount0;
418  remainingEntryOrdinal[1] /= componentEntryCount1;
419  remainingEntryOrdinal[2] /= componentEntryCount2;
420 
421  const ordinal_type componentRank = component.rank();
422 
423  if (componentRank == 3)
424  {
425  value *= component(componentEntryOrdinal0,componentEntryOrdinal1,componentEntryOrdinal2);
426  }
427  else if (componentRank == 2)
428  {
429  value *= component(componentEntryOrdinal0,componentEntryOrdinal1);
430  }
431  else if (componentRank == 1)
432  {
433  value *= component(componentEntryOrdinal0);
434  }
435  else
436  {
437  INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE(true, std::logic_error, "unsupported component rank encountered");
438  }
439  }
440  return value;
441  }
442 
450  template <typename iType0, typename iType1, ordinal_type numTensorComponents>
451  KOKKOS_INLINE_FUNCTION typename std::enable_if<
452  (std::is_integral<iType0>::value && std::is_integral<iType1>::value),
453  Scalar>::type
454  operator()(const Kokkos::Array<iType0,numTensorComponents>& entryComponents0, const Kokkos::Array<iType1,numTensorComponents>& entryComponents1) const {
455 #ifdef HAVE_INTREPID2_DEBUG
456  INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE(numTensorComponents_ != numTensorComponents, std::invalid_argument, "Tensorial component count mismatch");
457 #endif
458  Scalar value = 1.0;
459  for (ordinal_type r=0; r<numTensorComponents; r++)
460  {
461  auto & component = tensorComponents_[r];
462  const ordinal_type componentRank = component.rank();
463  if (componentRank == 2)
464  {
465  value *= component(entryComponents0[r],entryComponents1[r]);
466  }
467  else if (componentRank == 1)
468  {
469  value *= component(entryComponents0[r]);
470  }
471  else
472  {
473  INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE(true, std::logic_error, "unsupported component rank encountered");
474  }
475  }
476  return value;
477  }
478 
484  template <typename iType>
485  KOKKOS_INLINE_FUNCTION
486  typename std::enable_if<std::is_integral<iType>::value, ordinal_type>::type
487  extent_int(const iType& d) const {
488  return extents_[d];
489  }
490 
496  template <typename iType>
497  KOKKOS_INLINE_FUNCTION constexpr
498  typename std::enable_if<std::is_integral<iType>::value, size_t>::type
499  extent(const iType& d) const {
500  return extents_[d];
501  }
502 
504  KOKKOS_INLINE_FUNCTION constexpr bool isValid() const
505  {
506  return extents_[0] > 0;
507  }
508 
510  KOKKOS_INLINE_FUNCTION
511  ordinal_type rank() const
512  {
513  return rank_;
514  }
515 
517  KOKKOS_INLINE_FUNCTION
518  ordinal_type numTensorComponents() const
519  {
520  return numTensorComponents_;
521  }
522 
524  KOKKOS_INLINE_FUNCTION
526  {
527  return separateFirstComponent_;
528  }
529 
531  void setFirstComponentExtentInDimension0(const ordinal_type &newExtent)
532  {
533  INTREPID2_TEST_FOR_EXCEPTION(!separateFirstComponent_ && (numTensorComponents_ != 1), std::invalid_argument, "setFirstComponentExtent() is only allowed when separateFirstComponent_ is true, or there is only one component");
534  tensorComponents_[0].setExtent(0,newExtent);
535  extents_[0] = newExtent;
536  }
537  };
538 }
539 
540 #endif /* Intrepid2_TensorData_h */
TensorData(Kokkos::Array< Data< Scalar, DeviceType >, numTensorComponents > tensorComponents, bool separateFirstComponent=false)
Constructor with fixed-length Kokkos::Array argument.
#define INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE(test, x, msg)
Defines the Data class, a wrapper around a Kokkos::View that allows data that is constant or repeatin...
TensorData(const TensorData &first, const TensorData &second, bool separateFirstComponent=false)
Constructor to combine two other TensorData objects.
KOKKOS_INLINE_FUNCTION ordinal_type getTensorComponentIndex(const ordinal_type &tensorComponent, const ordinal_type &dim, const ordinal_type &enumerationIndex) const
return the index into the specified tensorial component in the dimension specified corresponding to t...
KOKKOS_INLINE_FUNCTION const Data< Scalar, DeviceType > & getTensorComponent(const ordinal_type &r) const
Returns the requested tensor component.
KOKKOS_INLINE_FUNCTION constexpr bool isValid() const
Returns true for containers that have data; false for those that don&#39;t (e.g., those that have been co...
void initialize()
Initialize members based on constructor parameters.
TensorData(const TensorData< Scalar, OtherDeviceType > &tensorData)
Copy-like constructor for differing execution spaces. This performs a deep copy of the underlying dat...
KOKKOS_INLINE_FUNCTION ordinal_type rank() const
Returns the rank of the container.
TensorData(TensorData otherTensorData, std::vector< int > whichComps)
Constructor that takes a subset of the tensorial components of another TensorData container...
TensorData(std::vector< Data< Scalar, DeviceType > > tensorComponents, bool separateFirstComponent=false)
Constructor with variable-length std::vector containing the components.
Contains definitions of custom data types in Intrepid2.
KOKKOS_INLINE_FUNCTION std::enable_if< std::is_integral< iType0 >::value, Scalar >::type operator()(const iType0 &tensorEntryIndex) const
Accessor for rank-1 objects.
TensorData()
Default constructor.
void setFirstComponentExtentInDimension0(const ordinal_type &newExtent)
Sets the extent of the first component. Only valid when either there is only one component, or when separateFirstComponent() returns true. The intended use case is when the 0 dimension in first component represents a cell index, and the container is resized to match a workset size that does not evenly divide the number of cells.
KOKKOS_INLINE_FUNCTION std::enable_if< std::is_integral< iType0 >::value, Scalar >::type operator()(const Kokkos::Array< iType0, numTensorComponents > &entryComponents) const
Accessor that accepts a fixed-length array with entries corresponding to component indices...
View-like interface to tensor data; tensor components are stored separately and multiplied together a...
KOKKOS_INLINE_FUNCTION std::enable_if< std::is_integral< iType >::value, ordinal_type >::type extent_int(const iType &d) const
Returns the logical extent in the requested dimension.
KOKKOS_INLINE_FUNCTION bool separateFirstComponent() const
Returns true if the first component is indexed separately; false if not.
KOKKOS_INLINE_FUNCTION ordinal_type numTensorComponents() const
Return the number of tensorial components.
TensorData(Data< Scalar, DeviceType > tensorComponent)
Simple constructor for the case of trivial tensor-product structure (single component) ...
KOKKOS_INLINE_FUNCTION constexpr std::enable_if< std::is_integral< iType >::value, size_t >::type extent(const iType &d) const
Returns the logical extent in the requested dimension.
static constexpr ordinal_type MaxTensorComponents
Maximum number of tensor/Cartesian products that can be taken: this allows hypercube basis in 7D to b...