Intrepid2
Intrepid2_HierarchicalBasis_HDIV_TRI.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 
51 #ifndef Intrepid2_HierarchicalBasis_HDIV_TRI_h
52 #define Intrepid2_HierarchicalBasis_HDIV_TRI_h
53 
54 #include <Kokkos_DynRankView.hpp>
55 
56 #include <Intrepid2_config.h>
57 
58 #include "Intrepid2_Basis.hpp"
61 #include "Intrepid2_Utils.hpp"
62 
63 namespace Intrepid2
64 {
65 
74  template<typename DeviceType,
75  typename OutputScalar = double,
76  typename PointScalar = double,
77  bool useCGBasis = true> // if useCGBasis is false, all basis functions will be associated with the interior
79  : public HierarchicalBasis_HCURL_TRI<DeviceType,OutputScalar,PointScalar,useCGBasis>
80  {
81  public:
83 
85 
86  using typename BasisBase::OrdinalTypeArray1DHost;
87  using typename BasisBase::OrdinalTypeArray2DHost;
88 
89  using typename BasisBase::OutputViewType;
90  using typename BasisBase::PointViewType;
91  using typename BasisBase::ScalarViewType;
92 
93  using typename BasisBase::ExecutionSpace;
94 
95  public:
100  HierarchicalBasis_HDIV_TRI(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
101  :
102  CurlBasis(polyOrder,pointType)
103  {
104  this->functionSpace_ = FUNCTION_SPACE_HDIV; // CurlBasis will set the other Basis member variables correctly…
105  }
106 
111  const char* getName() const override {
112  return "Intrepid2_HierarchicalBasis_HDIV_TRI";
113  }
114 
117  virtual bool requireOrientation() const override {
118  return (this->getDegree() > 2);
119  }
120 
121  // since the getValues() below only overrides the FEM variant, we specify that
122  // we use the base class's getValues(), which implements the FVD variant by throwing an exception.
123  // (It's an error to use the FVD variant on this basis.)
124  using BasisBase::getValues;
125 
144  virtual void getValues( OutputViewType outputValues, const PointViewType inputPoints,
145  const EOperator operatorType = OPERATOR_VALUE ) const override
146  {
147  // H(div) functions are a rotation of the corresponding H(curl) functions
148  if (operatorType == OPERATOR_DIV)
149  {
150  this->CurlBasis::getValues(outputValues, inputPoints, OPERATOR_CURL);
151  }
152  else if (operatorType == OPERATOR_VALUE)
153  {
154  this->CurlBasis::getValues(outputValues, inputPoints, OPERATOR_VALUE);
155 
156  const ordinal_type numFields = outputValues.extent_int(0);
157  const ordinal_type numPoints = outputValues.extent_int(1);
158 
159  auto policy = Kokkos::MDRangePolicy<ExecutionSpace,Kokkos::Rank<2>>(Kokkos::Array<int,2>{0,0},Kokkos::Array<int,2>{numFields,numPoints});
160 
161  // multiply by [[0 1] [-1 0]]
162  Kokkos::parallel_for("apply rotation to H(curl) values", policy,
163  KOKKOS_LAMBDA (const int &fieldOrdinal, const int &pointOrdinal) {
164  const auto tempValue = outputValues(fieldOrdinal,pointOrdinal,0);
165  outputValues(fieldOrdinal,pointOrdinal,0) = outputValues(fieldOrdinal,pointOrdinal,1);
166  outputValues(fieldOrdinal,pointOrdinal,1) = -tempValue;
167  });
168  }
169  }
170 
176  getHostBasis() const override {
177  using HostDeviceType = typename Kokkos::HostSpace::device_type;
179  return Teuchos::rcp( new HostBasisType(this->CurlBasis::polyOrder_) );
180  }
181  };
182 } // end namespace Intrepid2
183 
184 #endif /* Intrepid2_HierarchicalBasis_HDIV_TRI_h */
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
const char * getName() const override
Returns basis name.
EFunctionSpace functionSpace_
The function space in which the basis is defined.
virtual bool requireOrientation() const override
True if orientation is required.
Free functions, callable from device code, that implement various polynomials useful in basis definit...
H(curl) basis on the triangle using a construction involving Legendre and integrated Jacobi polynomia...
Header function for Intrepid2::Util class and other utility functions.
For mathematical details of the construction, see:
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
HierarchicalBasis_HDIV_TRI(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
virtual void getValues(OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
EPointType
Enumeration of types of point distributions in Intrepid.
ordinal_type getDegree() const
Returns the degree of the basis.
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.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Header file for the abstract base class Intrepid2::Basis.
For mathematical details of the construction, see:
virtual BasisPtr< typename Kokkos::HostSpace::device_type, OutputScalar, PointScalar > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...