Intrepid2
Intrepid2_HDIV_TET_I1_FEM.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_HDIV_TET_I1_FEM_HPP__
17#define __INTREPID2_HDIV_TET_I1_FEM_HPP__
18
19#include "Intrepid2_Basis.hpp"
21
22namespace Intrepid2 {
23
70 namespace Impl {
71
76 public:
77 typedef struct Tetrahedron<4> cell_topology_type;
81 template<EOperator opType>
82 struct Serial {
83 template<typename OutputViewType,
84 typename inputViewType>
85 KOKKOS_INLINE_FUNCTION
86 static void
87 getValues( OutputViewType output,
88 const inputViewType input );
89
90 };
91
92 template<typename DeviceType,
93 typename outputValueValueType, class ...outputValueProperties,
94 typename inputPointValueType, class ...inputPointProperties>
95 static void
96 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
97 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
98 const EOperator operatorType);
99
103 template<typename outputValueViewType,
104 typename inputPointViewType,
105 EOperator opType>
106 struct Functor {
107 outputValueViewType _outputValues;
108 const inputPointViewType _inputPoints;
109
110 KOKKOS_INLINE_FUNCTION
111 Functor( outputValueViewType outputValues_,
112 inputPointViewType inputPoints_ )
113 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
114
115 KOKKOS_INLINE_FUNCTION
116 void operator()(const ordinal_type pt) const {
117 switch (opType) {
118 case OPERATOR_VALUE : {
119 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
120 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
121 Serial<opType>::getValues( output, input );
122 break;
123 }
124 case OPERATOR_DIV : {
125 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
126 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
127 Serial<opType>::getValues( output, input );
128 break;
129 }
130 default: {
131 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
132 opType != OPERATOR_DIV,
133 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::Serial::getValues) operator is not supported");
134 }
135 }
136 }
137 };
138
139 };
140 }
141
142 template<typename DeviceType = void,
143 typename outputValueType = double,
144 typename pointValueType = double>
145 class Basis_HDIV_TET_I1_FEM: public Basis<DeviceType,outputValueType,pointValueType> {
146 public:
148
152
153 using typename BasisBase::OutputViewType;
154 using typename BasisBase::PointViewType ;
155 using typename BasisBase::ScalarViewType;
156
160
162
163 virtual
164 void
165 getValues( OutputViewType outputValues,
166 const PointViewType inputPoints,
167 const EOperator operatorType = OPERATOR_VALUE ) const override {
168#ifdef HAVE_INTREPID2_DEBUG
169 // Verify arguments
171 inputPoints,
172 operatorType,
173 this->getBaseCellTopology(),
174 this->getCardinality() );
175#endif
176 Impl::Basis_HDIV_TET_I1_FEM::
177 getValues<DeviceType>( outputValues,
178 inputPoints,
179 operatorType );
180 }
181
182 virtual void
183 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
184 ordinal_type& perThreadSpaceSize,
185 const PointViewType inputPointsconst,
186 const EOperator operatorType = OPERATOR_VALUE) const override;
187
188 KOKKOS_INLINE_FUNCTION
189 virtual void
190 getValues(
191 OutputViewType outputValues,
192 const PointViewType inputPoints,
193 const EOperator operatorType,
194 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
195 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
196 const ordinal_type subcellDim = -1,
197 const ordinal_type subcellOrdinal = -1) const override;
198
199 virtual
200 void
201 getDofCoords( ScalarViewType dofCoords ) const override {
202#ifdef HAVE_INTREPID2_DEBUG
203 // Verify rank of output array.
204 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
205 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
206 // Verify 0th dimension of output array.
207 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
208 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
209 // Verify 1st dimension of output array.
210 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
211 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
212#endif
213 Kokkos::deep_copy(dofCoords, this->dofCoords_);
214 }
215
216 virtual
217 void
218 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
219#ifdef HAVE_INTREPID2_DEBUG
220 // Verify rank of output array.
221 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
222 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
223 // Verify 0th dimension of output array.
224 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
225 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
226 // Verify 1st dimension of output array.
227 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
228 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
229#endif
230 Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
231 }
232
233 virtual
234 const char*
235 getName() const override {
236 return "Intrepid2_HDIV_TET_I1_FEM";
237 }
238
239 virtual
240 bool
241 requireOrientation() const override {
242 return true;
243 }
244
255 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
256
257 if(subCellDim == 2) {
258 return Teuchos::rcp(new
259 Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Triangle<3> >()));
260 }
261 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
262 }
263
268 };
269
270}// namespace Intrepid2
271
273
274#endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HDIV_Args(const outputValueViewType outputValues, const inputPointViewType inputPoints, const EOperator operatorType, const shards::CellTopology cellTopo, const ordinal_type basisCard)
Runtime check of the arguments for the getValues method in an HDIV-conforming FEM basis....
Definition file for FEM basis functions of degree 1 for H(div) functions on TET cells.
Header file for the Intrepid2::Basis_HVOL_C0_FEM class.
Implementation of the default H(div)-compatible FEM basis of degree 1 on a Tetrahedron cell.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual const char * getName() const override
Returns basis name.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPointsconst, const EOperator operatorType=OPERATOR_VALUE) const override
Return the size of the scratch space, in bytes, needed for using the team-level implementation of get...
virtual bool requireOrientation() const override
True if orientation is required.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
BasisPtr< typename Kokkos::HostSpace::device_type, outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral,...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
virtual KOKKOS_INLINE_FUNCTION void getValues(OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
Team-level evaluation of basis functions on a reference cell.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host 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.
See Intrepid2::Basis_HDIV_TET_I1_FEM.