Intrepid2
Intrepid2_Basis.hpp
Go to the documentation of this file.
1// @HEADER
2// ************************************************************************
3//
4// Intrepid2 Package
5// Copyright (2007) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Kyungjoo Kim (kyukim@sandia.gov),
38// Mauro Perego (mperego@sandia.gov), or
39// Nate Roberts (nvrober@sandia.gov)
40//
41// ************************************************************************
42// @HEADER
43
50#ifndef __INTREPID2_BASIS_HPP__
51#define __INTREPID2_BASIS_HPP__
52
53#include "Intrepid2_ConfigDefs.hpp"
54#include "Intrepid2_Types.hpp"
55#include "Intrepid2_Utils.hpp"
56
60#include "Kokkos_Vector.hpp"
61#include "Shards_CellTopology.hpp"
62#include <Teuchos_RCPDecl.hpp>
63
64#include <vector>
65
66namespace Intrepid2 {
67
68template<typename DeviceType = void,
69 typename OutputType = double,
70 typename PointType = double>
71class Basis;
72
75template <typename DeviceType = void, typename OutputType = double, typename PointType = double>
76using BasisPtr = Teuchos::RCP<Basis<DeviceType,OutputType,PointType> >;
77
80template <typename OutputType = double, typename PointType = double>
82
119 template<typename Device,
120 typename outputValueType,
121 typename pointValueType>
122 class Basis {
123 public:
126 using DeviceType = Device;
127
130 using ExecutionSpace = typename DeviceType::execution_space;
131
132
135 using OutputValueType = outputValueType;
136
139 using PointValueType = pointValueType;
140
143 using OrdinalViewType = Kokkos::View<ordinal_type,DeviceType>;
144
147 using EBasisViewType = Kokkos::View<EBasis,DeviceType>;
148
151 using ECoordinatesViewType = Kokkos::View<ECoordinates,DeviceType>;
152
153 // ** tag interface
154 // - tag interface is not decorated with Kokkos inline so it should be allocated on hostspace
155
158 using OrdinalTypeArray1DHost = Kokkos::View<ordinal_type*,typename ExecutionSpace::array_layout,Kokkos::HostSpace>;
159
162 using OrdinalTypeArray2DHost = Kokkos::View<ordinal_type**,typename ExecutionSpace::array_layout,Kokkos::HostSpace>;
163
166 using OrdinalTypeArray3DHost = Kokkos::View<ordinal_type***,typename ExecutionSpace::array_layout,Kokkos::HostSpace>;
167
170 using OrdinalTypeArrayStride1DHost = Kokkos::View<ordinal_type*, Kokkos::LayoutStride, Kokkos::HostSpace>;
171
174 using OrdinalTypeArray1D = Kokkos::View<ordinal_type*,DeviceType>;
175
178 using OrdinalTypeArray2D = Kokkos::View<ordinal_type**,DeviceType>;
179
182 using OrdinalTypeArray3D = Kokkos::View<ordinal_type***,DeviceType>;
183
186 using OrdinalTypeArrayStride1D = Kokkos::View<ordinal_type*, Kokkos::LayoutStride, DeviceType>;
187
190 typedef typename ScalarTraits<pointValueType>::scalar_type scalarType;
191 protected:
192
195 ordinal_type basisCardinality_;
196
199 ordinal_type basisDegree_;
200
205 shards::CellTopology basisCellTopology_;
206
210
213 ECoordinates basisCoordinates_;
214
217 EFunctionSpace functionSpace_ = FUNCTION_SPACE_MAX;
218
221 //Kokkos::View<bool,DeviceType> basisTagsAreSet_;
222
235
248
260 template<typename OrdinalTypeView3D,
261 typename OrdinalTypeView2D,
262 typename OrdinalTypeView1D>
263 void setOrdinalTagData( OrdinalTypeView3D &tagToOrdinal,
264 OrdinalTypeView2D &ordinalToTag,
265 const OrdinalTypeView1D tags,
266 const ordinal_type basisCard,
267 const ordinal_type tagSize,
268 const ordinal_type posScDim,
269 const ordinal_type posScOrd,
270 const ordinal_type posDfOrd ) {
271 // Create ordinalToTag
272 ordinalToTag = OrdinalTypeView2D("ordinalToTag", basisCard, tagSize);
273
274 // Initialize with -1
275 Kokkos::deep_copy( ordinalToTag, -1 );
276
277 // Copy tags
278 for (ordinal_type i=0;i<basisCard;++i)
279 for (ordinal_type j=0;j<tagSize;++j)
280 ordinalToTag(i, j) = tags(i*tagSize + j);
281
282 // Find out dimension of tagToOrdinal
283 auto maxScDim = 0; // first dimension of tagToOrdinal
284 for (ordinal_type i=0;i<basisCard;++i)
285 if (maxScDim < tags(i*tagSize + posScDim))
286 maxScDim = tags(i*tagSize + posScDim);
287 ++maxScDim;
288
289 auto maxScOrd = 0; // second dimension of tagToOrdinal
290 for (ordinal_type i=0;i<basisCard;++i)
291 if (maxScOrd < tags(i*tagSize + posScOrd))
292 maxScOrd = tags(i*tagSize + posScOrd);
293 ++maxScOrd;
294
295 auto maxDfOrd = 0; // third dimension of tagToOrdinal
296 for (ordinal_type i=0;i<basisCard;++i)
297 if (maxDfOrd < tags(i*tagSize + posDfOrd))
298 maxDfOrd = tags(i*tagSize + posDfOrd);
299 ++maxDfOrd;
300
301 // Create tagToOrdinal
302 tagToOrdinal = OrdinalTypeView3D("tagToOrdinal", maxScDim, maxScOrd, maxDfOrd);
303
304 // Initialize with -1
305 Kokkos::deep_copy( tagToOrdinal, -1 );
306
307 // Overwrite elements of the array corresponding to tags with local DoF Id's, leave all other = -1
308 for (ordinal_type i=0;i<basisCard;++i)
309 tagToOrdinal(tags(i*tagSize), tags(i*tagSize+1), tags(i*tagSize+2)) = i;
310 }
311
312 // dof coords
315 Kokkos::DynRankView<scalarType,DeviceType> dofCoords_;
316
317 // dof coeffs
325 Kokkos::DynRankView<scalarType,DeviceType> dofCoeffs_;
326
335 public:
336
337 Basis() = default;
338 virtual~Basis() = default;
339
340 // receives input arguments
343 OutputValueType getDummyOutputValue() { return outputValueType(); }
344
347 PointValueType getDummyPointValue() { return pointValueType(); }
348
351 using OutputViewType = Kokkos::DynRankView<OutputValueType,Kokkos::LayoutStride,DeviceType>;
352
355 using PointViewType = Kokkos::DynRankView<PointValueType,Kokkos::LayoutStride,DeviceType>;
356
359 using ScalarViewType = Kokkos::DynRankView<scalarType,Kokkos::LayoutStride,DeviceType>;
360
365 Kokkos::DynRankView<OutputValueType,DeviceType> allocateOutputView( const int numPoints, const EOperator operatorType = OPERATOR_VALUE) const;
366
372 virtual BasisValues<OutputValueType,DeviceType> allocateBasisValues( TensorPoints<PointValueType,DeviceType> points, const EOperator operatorType = OPERATOR_VALUE) const
373 {
374 const bool operatorIsDk = (operatorType >= OPERATOR_D1) && (operatorType <= OPERATOR_D10);
375 const bool operatorSupported = (operatorType == OPERATOR_VALUE) || (operatorType == OPERATOR_GRAD) || (operatorType == OPERATOR_CURL) || (operatorType == OPERATOR_DIV) || operatorIsDk;
376 INTREPID2_TEST_FOR_EXCEPTION(!operatorSupported, std::invalid_argument, "operator is not supported by allocateBasisValues");
377
378// // this default implementation employs a trivial tensor-product structure; make sure that points also have a trivial tensor product structure:
379// INTREPID2_TEST_FOR_EXCEPTION(points.numTensorComponents() != 1, std::invalid_argument, "default implementation of allocateBasisValues() only supports a trivial tensor product structure (one tensor component)");
380
381 const int numPoints = points.extent_int(0);
382
383 using Scalar = OutputValueType;
384
385 auto dataView = allocateOutputView(numPoints, operatorType);
386 Data<Scalar,DeviceType> data(dataView);
387
388 bool useVectorData = (dataView.rank() == 3);
389
390 if (useVectorData)
391 {
392 VectorData<Scalar,DeviceType> vectorData(data);
393 return BasisValues<Scalar,DeviceType>(vectorData);
394 }
395 else
396 {
397 TensorData<Scalar,DeviceType> tensorData(data);
398 return BasisValues<Scalar,DeviceType>(tensorData);
399 }
400 }
401
420 virtual
421 void
422 getValues( OutputViewType /* outputValues */,
423 const PointViewType /* inputPoints */,
424 const EOperator /* operatorType */ = OPERATOR_VALUE ) const {
425 INTREPID2_TEST_FOR_EXCEPTION( true, std::logic_error,
426 ">>> ERROR (Basis::getValues): this method (FEM) is not supported or should be overridden accordingly by derived classes.");
427 }
428
440 virtual
441 void
444 const EOperator operatorType = OPERATOR_VALUE ) const {
445 // note the extra allocation/copy here (this is one reason, among several, to override this method):
446 auto rawExpandedPoints = inputPoints.allocateAndFillExpandedRawPointView();
447
448 OutputViewType rawOutputView;
450 if (outputValues.numTensorDataFamilies() > 0)
451 {
452 INTREPID2_TEST_FOR_EXCEPTION(outputValues.tensorData(0).numTensorComponents() != 1, std::invalid_argument, "default implementation of getValues() only supports outputValues with trivial tensor-product structure");
453 outputData = outputValues.tensorData().getTensorComponent(0);
454 }
455 else if (outputValues.vectorData().isValid())
456 {
457 INTREPID2_TEST_FOR_EXCEPTION(outputValues.vectorData().numComponents() != 1, std::invalid_argument, "default implementation of getValues() only supports outputValues with trivial tensor-product structure");
458 INTREPID2_TEST_FOR_EXCEPTION(outputValues.vectorData().getComponent(0).numTensorComponents() != 1, std::invalid_argument, "default implementation of getValues() only supports outputValues with trivial tensor-product structure");
459 outputData = outputValues.vectorData().getComponent(0).getTensorComponent(0);
460 }
461
462 this->getValues(outputData.getUnderlyingView(), rawExpandedPoints, operatorType);
463 }
464
484 virtual
485 void
486 getValues( OutputViewType /* outputValues */,
487 const PointViewType /* inputPoints */,
488 const PointViewType /* cellVertices */,
489 const EOperator /* operatorType */ = OPERATOR_VALUE ) const {
490 INTREPID2_TEST_FOR_EXCEPTION( true, std::logic_error,
491 ">>> ERROR (Basis::getValues): this method (FVM) is not supported or should be overridden accordingly by derived classes.");
492 }
493
494
498 virtual
499 void
500 getDofCoords( ScalarViewType /* dofCoords */ ) const {
501 INTREPID2_TEST_FOR_EXCEPTION( true, std::logic_error,
502 ">>> ERROR (Basis::getDofCoords): this method is not supported or should be overridden accordingly by derived classes.");
503 }
504
513 virtual
514 void
515 getDofCoeffs( ScalarViewType /* dofCoeffs */ ) const {
516 INTREPID2_TEST_FOR_EXCEPTION( true, std::logic_error,
517 ">>> ERROR (Basis::getDofCoeffs): this method is not supported or should be overridden accordingly by derived classes.");
518 }
519
528 {
529 INTREPID2_TEST_FOR_EXCEPTION( basisType_ != BASIS_FEM_HIERARCHICAL, std::logic_error,
530 ">>> ERROR (Basis::getFieldOrdinalsForDegree): this method is not supported for non-hierarchical bases.");
531 int degreeEntryLength = fieldOrdinalPolynomialDegree_.extent_int(1);
532 int requestedDegreeLength = degrees.extent_int(0);
533 INTREPID2_TEST_FOR_EXCEPTION(degreeEntryLength != requestedDegreeLength, std::invalid_argument, "length of degrees does not match the entries in fieldOrdinalPolynomialDegree_");
534 std::vector<int> fieldOrdinalsVector;
535 for (int basisOrdinal=0; basisOrdinal<fieldOrdinalPolynomialDegree_.extent_int(0); basisOrdinal++)
536 {
537 bool matches = true;
538 for (int d=0; d<degreeEntryLength; d++)
539 {
540 if (fieldOrdinalPolynomialDegree_(basisOrdinal,d) > degrees(d)) matches = false;
541 }
542 if (matches) fieldOrdinalsVector.push_back(basisOrdinal);
543 }
544 OrdinalTypeArray1DHost fieldOrdinals("fieldOrdinalsForDegree",fieldOrdinalsVector.size());
545 for (unsigned i=0; i<fieldOrdinalsVector.size(); i++)
546 {
547 fieldOrdinals(i) = fieldOrdinalsVector[i];
548 }
549 return fieldOrdinals;
550 }
551
563 std::vector<int> getFieldOrdinalsForDegree(std::vector<int> &degrees) const
564 {
565 INTREPID2_TEST_FOR_EXCEPTION( basisType_ != BASIS_FEM_HIERARCHICAL, std::logic_error,
566 ">>> ERROR (Basis::getFieldOrdinalsForDegree): this method is not supported for non-hierarchical bases.");
567 OrdinalTypeArray1DHost degreesView("degrees",degrees.size());
568 for (unsigned d=0; d<degrees.size(); d++)
569 {
570 degreesView(d) = degrees[d];
571 }
572 auto fieldOrdinalsView = getFieldOrdinalsForDegree(degreesView);
573 std::vector<int> fieldOrdinalsVector(fieldOrdinalsView.extent_int(0));
574 for (int i=0; i<fieldOrdinalsView.extent_int(0); i++)
575 {
576 fieldOrdinalsVector[i] = fieldOrdinalsView(i);
577 }
578 return fieldOrdinalsVector;
579 }
580
588 {
589 INTREPID2_TEST_FOR_EXCEPTION( basisType_ != BASIS_FEM_HIERARCHICAL, std::logic_error,
590 ">>> ERROR (Basis::getPolynomialDegreeOfField): this method is not supported for non-hierarchical bases.");
591 INTREPID2_TEST_FOR_EXCEPTION(fieldOrdinal < 0, std::invalid_argument, "field ordinal must be non-negative");
592 INTREPID2_TEST_FOR_EXCEPTION(fieldOrdinal >= fieldOrdinalPolynomialDegree_.extent_int(0), std::invalid_argument, "field ordinal out of bounds");
593
594 int polyDegreeLength = getPolynomialDegreeLength();
595 OrdinalTypeArray1DHost polyDegree("polynomial degree", polyDegreeLength);
596 for (int d=0; d<polyDegreeLength; d++)
597 {
598 polyDegree(d) = fieldOrdinalPolynomialDegree_(fieldOrdinal,d);
599 }
600 return polyDegree;
601 }
602
610 std::vector<int> getPolynomialDegreeOfFieldAsVector(int fieldOrdinal) const
611 {
612 INTREPID2_TEST_FOR_EXCEPTION( basisType_ != BASIS_FEM_HIERARCHICAL, std::logic_error,
613 ">>> ERROR (Basis::getPolynomialDegreeOfFieldAsVector): this method is not supported for non-hierarchical bases.");
614 auto polynomialDegreeView = getPolynomialDegreeOfField(fieldOrdinal);
615 std::vector<int> polynomialDegree(polynomialDegreeView.extent_int(0));
616
617 for (unsigned d=0; d<polynomialDegree.size(); d++)
618 {
619 polynomialDegree[d] = polynomialDegreeView(d);
620 }
621 return polynomialDegree;
622 }
623
627 {
628 INTREPID2_TEST_FOR_EXCEPTION( basisType_ != BASIS_FEM_HIERARCHICAL, std::logic_error,
629 ">>> ERROR (Basis::getPolynomialDegreeLength): this method is not supported for non-hierarchical bases.");
630 return fieldOrdinalPolynomialDegree_.extent_int(1);
631 }
632
637 virtual
638 const char*
639 getName() const {
640 return "Intrepid2_Basis";
641 }
642
645 virtual
646 bool
648 return false;
649 }
650
655 ordinal_type
657 return basisCardinality_;
658 }
659
660
665 ordinal_type
666 getDegree() const {
667 return basisDegree_;
668 }
669
674 EFunctionSpace
676 return functionSpace_;
677 }
678
684 shards::CellTopology
686 return basisCellTopology_;
687 }
688
689
694 EBasis
695 getBasisType() const {
696 return basisType_;
697 }
698
699
704 ECoordinates
706 return basisCoordinates_;
707 }
708
716 ordinal_type
717 getDofCount( const ordinal_type subcDim,
718 const ordinal_type subcOrd ) const {
719 if ( subcDim >= 0 && subcDim < static_cast<ordinal_type>(tagToOrdinal_.extent(0)) &&
720 subcOrd >= 0 && subcOrd < static_cast<ordinal_type>(tagToOrdinal_.extent(1)) )
721 {
722 int firstDofOrdinal = tagToOrdinal_(subcDim, subcOrd, 0); // will be -1 if no dofs for subcell
723 if (firstDofOrdinal == -1) return static_cast<ordinal_type>(0);
724 // otherwise, lookup its tag and return the dof count stored there
725 return static_cast<ordinal_type>(this->getDofTag(firstDofOrdinal)[3]);
726 }
727 else
728 {
729 // subcDim and/or subcOrd out of bounds -> no dofs associated with subcell
730 return static_cast<ordinal_type>(0);
731 }
732 }
733
742 ordinal_type
743 getDofOrdinal( const ordinal_type subcDim,
744 const ordinal_type subcOrd,
745 const ordinal_type subcDofOrd ) const {
746 // this should be abort and able to be called as a device function
747#ifdef HAVE_INTREPID2_DEBUG
748 INTREPID2_TEST_FOR_EXCEPTION( subcDim < 0 || subcDim >= static_cast<ordinal_type>(tagToOrdinal_.extent(0)), std::out_of_range,
749 ">>> ERROR (Basis::getDofOrdinal): subcDim is out of range");
750 INTREPID2_TEST_FOR_EXCEPTION( subcOrd < 0 || subcOrd >= static_cast<ordinal_type>(tagToOrdinal_.extent(1)), std::out_of_range,
751 ">>> ERROR (Basis::getDofOrdinal): subcOrd is out of range");
752 INTREPID2_TEST_FOR_EXCEPTION( subcDofOrd < 0 || subcDofOrd >= static_cast<ordinal_type>(tagToOrdinal_.extent(2)), std::out_of_range,
753 ">>> ERROR (Basis::getDofOrdinal): subcDofOrd is out of range");
754#endif
755 ordinal_type r_val = -1;
756 if ( subcDim < static_cast<ordinal_type>(tagToOrdinal_.extent(0)) &&
757 subcOrd < static_cast<ordinal_type>(tagToOrdinal_.extent(1)) &&
758 subcDofOrd < static_cast<ordinal_type>(tagToOrdinal_.extent(2)) )
759 r_val = tagToOrdinal_(subcDim, subcOrd, subcDofOrd);
760#ifdef HAVE_INTREPID2_DEBUG
761 INTREPID2_TEST_FOR_EXCEPTION( r_val == -1, std::runtime_error,
762 ">>> ERROR (Basis::getDofOrdinal): Invalid DoF tag is found.");
763#endif
764 return r_val;
765 }
766
770 return tagToOrdinal_;
771 }
772
773
785 getDofTag( const ordinal_type dofOrd ) const {
786#ifdef HAVE_INTREPID2_DEBUG
787 INTREPID2_TEST_FOR_EXCEPTION( dofOrd < 0 || dofOrd >= static_cast<ordinal_type>(ordinalToTag_.extent(0)), std::out_of_range,
788 ">>> ERROR (Basis::getDofTag): dofOrd is out of range");
789#endif
790 return Kokkos::subview(ordinalToTag_, dofOrd, Kokkos::ALL());
791 }
792
793
804 return ordinalToTag_;
805 }
806
822 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const {
823 INTREPID2_TEST_FOR_EXCEPTION( true, std::logic_error,
824 ">>> ERROR (Basis::getSubCellRefBasis): this method is not supported or should be overridden accordingly by derived classes.");
825 }
826
832 getHostBasis() const {
833 INTREPID2_TEST_FOR_EXCEPTION( true, std::logic_error,
834 ">>> ERROR (Basis::getHostBasis): this method is not supported or should be overridden accordingly by derived classes.");
835 }
836
837 }; // class Basis
838
839 //--------------------------------------------------------------------------------------------//
840 // //
841 // Helper functions of the Basis class //
842 // //
843 //--------------------------------------------------------------------------------------------//
844
845 //
846 // functions for orders, cardinality, etc.
847 //
848
849
861 KOKKOS_INLINE_FUNCTION
862 ordinal_type getFieldRank(const EFunctionSpace spaceType);
863
899 KOKKOS_INLINE_FUNCTION
900 ordinal_type getOperatorRank(const EFunctionSpace spaceType,
901 const EOperator operatorType,
902 const ordinal_type spaceDim);
903
909 KOKKOS_INLINE_FUNCTION
910 ordinal_type getOperatorOrder(const EOperator operatorType);
911
912 template<EOperator operatorType>
913 KOKKOS_INLINE_FUNCTION
914 constexpr ordinal_type getOperatorOrder();
915
939 template<ordinal_type spaceDim>
940 KOKKOS_INLINE_FUNCTION
941 ordinal_type getDkEnumeration(const ordinal_type xMult,
942 const ordinal_type yMult = -1,
943 const ordinal_type zMult = -1);
944
945
956 template<ordinal_type spaceDim>
957 KOKKOS_INLINE_FUNCTION
958 ordinal_type getPnEnumeration(const ordinal_type p,
959 const ordinal_type q = 0,
960 const ordinal_type r = 0);
961
962
963
982template<typename value_type>
983KOKKOS_INLINE_FUNCTION
984void getJacobyRecurrenceCoeffs (
985 value_type &an,
986 value_type &bn,
987 value_type &cn,
988 const ordinal_type alpha,
989 const ordinal_type beta ,
990 const ordinal_type n);
991
992
993
994
995
996 // /** \brief Returns multiplicities of dx, dy, and dz based on the enumeration of the partial
997 // derivative, its order and the space dimension. Inverse of the getDkEnumeration() method.
998
999 // \param partialMult [out] - array with the multiplicities f dx, dy and dz
1000 // \param derivativeEnum [in] - enumeration of the partial derivative
1001 // \param operatorType [in] - k-th partial derivative Dk
1002 // \param spaceDim [in] - space dimension
1003 // */
1004 // template<typename OrdinalArrayType>
1005 // KOKKOS_INLINE_FUNCTION
1006 // void getDkMultiplicities(OrdinalArrayType partialMult,
1007 // const ordinal_type derivativeEnum,
1008 // const EOperator operatorType,
1009 // const ordinal_type spaceDim);
1010
1029 KOKKOS_INLINE_FUNCTION
1030 ordinal_type getDkCardinality(const EOperator operatorType,
1031 const ordinal_type spaceDim);
1032
1033 template<EOperator operatorType, ordinal_type spaceDim>
1034 KOKKOS_INLINE_FUNCTION
1035 constexpr ordinal_type getDkCardinality();
1036
1037
1038
1048 template<ordinal_type spaceDim>
1049 KOKKOS_INLINE_FUNCTION
1050 ordinal_type getPnCardinality (ordinal_type n);
1051
1052 template<ordinal_type spaceDim, ordinal_type n>
1053 KOKKOS_INLINE_FUNCTION
1054 constexpr ordinal_type getPnCardinality ();
1055
1056
1057
1058 //
1059 // Argument check
1060 //
1061
1062
1073 template<typename outputValueViewType,
1074 typename inputPointViewType>
1075 void getValues_HGRAD_Args( const outputValueViewType outputValues,
1076 const inputPointViewType inputPoints,
1077 const EOperator operatorType,
1078 const shards::CellTopology cellTopo,
1079 const ordinal_type basisCard );
1080
1091 template<typename outputValueViewType,
1092 typename inputPointViewType>
1093 void getValues_HCURL_Args( const outputValueViewType outputValues,
1094 const inputPointViewType inputPoints,
1095 const EOperator operatorType,
1096 const shards::CellTopology cellTopo,
1097 const ordinal_type basisCard );
1098
1109 template<typename outputValueViewType,
1110 typename inputPointViewType>
1111 void getValues_HDIV_Args( const outputValueViewType outputValues,
1112 const inputPointViewType inputPoints,
1113 const EOperator operatorType,
1114 const shards::CellTopology cellTopo,
1115 const ordinal_type basisCard );
1116
1127 template<typename outputValueViewType,
1128 typename inputPointViewType>
1129 void getValues_HVOL_Args( const outputValueViewType outputValues,
1130 const inputPointViewType inputPoints,
1131 const EOperator operatorType,
1132 const shards::CellTopology cellTopo,
1133 const ordinal_type basisCard );
1134
1135}// namespace Intrepid2
1136
1137#include <Intrepid2_BasisDef.hpp>
1138
1139//--------------------------------------------------------------------------------------------//
1140// //
1141// D O C U M E N T A T I O N P A G E S //
1142// //
1143//--------------------------------------------------------------------------------------------//
1245#endif
Implementation file for the abstract base class Intrepid2::Basis.
Header file for the data-wrapper class Intrepid2::BasisValues.
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type),...
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Definition of cell topology information.
View-like interface to tensor points; point components are stored separately; the appropriate coordin...
Contains definitions of custom data types in Intrepid2.
Header function for Intrepid2::Util class and other utility functions.
The data containers in Intrepid2 that support sum factorization and other reduced-data optimizations ...
const VectorDataType & vectorData() const
VectorData accessor.
TensorDataType & tensorData()
TensorData accessor for single-family scalar data.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
OrdinalTypeArray1DHost getFieldOrdinalsForDegree(OrdinalTypeArray1DHost &degrees) const
For hierarchical bases, returns the field ordinals that have at most the specified degree in each dim...
ECoordinates basisCoordinates_
The coordinate system for which the basis is defined.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
const OrdinalTypeArrayStride1DHost getDofTag(const ordinal_type dofOrd) const
DoF ordinal to DoF tag lookup.
outputValueType OutputValueType
Output value type for basis; default is double.
ordinal_type getDofOrdinal(const ordinal_type subcDim, const ordinal_type subcOrd, const ordinal_type subcDofOrd) const
DoF tag to ordinal lookup.
ordinal_type getDofCount(const ordinal_type subcDim, const ordinal_type subcOrd) const
DoF count for specified subcell.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
const OrdinalTypeArray2DHost getAllDofTags() const
Retrieves all DoF tags.
EBasis basisType_
Type of the basis.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
Kokkos::View< ordinal_type *, Kokkos::LayoutStride, DeviceType > OrdinalTypeArrayStride1D
View type for 1d device array.
Kokkos::View< ordinal_type *, Kokkos::LayoutStride, Kokkos::HostSpace > OrdinalTypeArrayStride1DHost
View type for 1d host array.
ordinal_type getDegree() const
Returns the degree of the basis.
void setOrdinalTagData(OrdinalTypeView3D &tagToOrdinal, OrdinalTypeView2D &ordinalToTag, const OrdinalTypeView1D tags, const ordinal_type basisCard, const ordinal_type tagSize, const ordinal_type posScDim, const ordinal_type posScOrd, const ordinal_type posDfOrd)
Fills ordinalToTag_ and tagToOrdinal_ by basis-specific tag data.
virtual HostBasisPtr< OutputValueType, PointValueType > getHostBasis() const
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
ordinal_type getCardinality() const
Returns cardinality of the basis.
virtual BasisValues< OutputValueType, DeviceType > allocateBasisValues(TensorPoints< PointValueType, DeviceType > points, const EOperator operatorType=OPERATOR_VALUE) const
Allocate BasisValues container suitable for passing to the getValues() variant that takes a TensorPoi...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
OutputValueType getDummyOutputValue()
Dummy array to receive input arguments.
ECoordinates getCoordinateSystem() const
Returns the type of coordinate system for which the basis is defined.
OrdinalTypeArray2DHost ordinalToTag_
"true" if tagToOrdinal_ and ordinalToTag_ have been initialized
PointValueType getDummyPointValue()
Dummy array to receive input arguments.
int getPolynomialDegreeLength() const
For hierarchical bases, returns the number of entries required to specify the polynomial degree of a ...
virtual bool requireOrientation() const
True if orientation is required.
Kokkos::View< ordinal_type, DeviceType > OrdinalViewType
View type for ordinal.
const OrdinalTypeArray3DHost getAllDofOrdinal() const
DoF tag to ordinal data structure.
Device DeviceType
(Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_dev...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::View< ECoordinates, DeviceType > ECoordinatesViewType
View for coordinate system type.
virtual void getValues(BasisValues< OutputValueType, DeviceType > outputValues, const TensorPoints< PointValueType, DeviceType > inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell, using point and output value containers that allow pre...
pointValueType PointValueType
Point value type for basis; default is double.
Kokkos::View< ordinal_type **, DeviceType > OrdinalTypeArray2D
View type for 2d device array.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
OrdinalTypeArray1DHost getPolynomialDegreeOfField(int fieldOrdinal) const
For hierarchical bases, returns the polynomial degree (which may have multiple values in higher spati...
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
virtual void getDofCoords(ScalarViewType) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom.
OrdinalTypeArray3DHost tagToOrdinal_
DoF tag to ordinal lookup table.
EFunctionSpace getFunctionSpace() const
Returns the function space for the basis.
Kokkos::DynRankView< OutputValueType, DeviceType > allocateOutputView(const int numPoints, const EOperator operatorType=OPERATOR_VALUE) const
Allocate a View container suitable for passing to the getValues() variant that accepts Kokkos DynRank...
virtual void getValues(OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
virtual void getDofCoeffs(ScalarViewType) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
EBasis getBasisType() const
Returns the basis type.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Kokkos::View< ordinal_type ***, DeviceType > OrdinalTypeArray3D
View type for 3d device array.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...
Kokkos::View< ordinal_type *, DeviceType > OrdinalTypeArray1D
View type for 1d device array.
virtual void getValues(OutputViewType, const PointViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of an FVD basis evaluation on a physical cell.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
virtual const char * getName() const
Returns basis name.
std::vector< int > getFieldOrdinalsForDegree(std::vector< int > &degrees) const
For hierarchical bases, returns the field ordinals that have at most the specified degree in each dim...
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now....
virtual BasisPtr< DeviceType, OutputValueType, PointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const
returns the basis associated to a subCell.
EFunctionSpace functionSpace_
The function space in which the basis is defined.
std::vector< int > getPolynomialDegreeOfFieldAsVector(int fieldOrdinal) const
For hierarchical bases, returns the polynomial degree (which may have multiple values in higher spati...
Kokkos::View< EBasis, DeviceType > EBasisViewType
View for basis type.
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.
View-like interface to tensor data; tensor components are stored separately and multiplied together a...
KOKKOS_INLINE_FUNCTION const Data< Scalar, DeviceType > & getTensorComponent(const ordinal_type &r) const
Returns the requested tensor component.
KOKKOS_INLINE_FUNCTION ordinal_type numTensorComponents() const
Return the number of tensorial components.
View-like interface to tensor points; point components are stored separately; the appropriate coordin...
KOKKOS_INLINE_FUNCTION std::enable_if< std::is_integral< iType >::value, int >::type extent_int(const iType &r) const
Returns the logical extent in the requested dimension.
ScalarView< PointScalar, DeviceType > allocateAndFillExpandedRawPointView() const
This method is for compatibility with existing methods that take raw point views. Note that in genera...
Reference-space field values for a basis, designed to support typical vector-valued bases.
KOKKOS_INLINE_FUNCTION int numComponents() const
returns the number of components
KOKKOS_INLINE_FUNCTION const TensorData< Scalar, DeviceType > & getComponent(const int &componentOrdinal) const
Single-argument component accessor for the axial-component or the single-family case; in this case,...
KOKKOS_INLINE_FUNCTION constexpr bool isValid() const
returns true for containers that have data; false for those that don't (e.g., those that have been co...