Intrepid2
Intrepid2_HGRAD_LINE_C2_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
14#ifndef __INTREPID2_HGRAD_LINE_C2_FEM_HPP__
15#define __INTREPID2_HGRAD_LINE_C2_FEM_HPP__
16
17#include "Intrepid2_Basis.hpp"
18
19namespace Intrepid2 {
20
45 namespace Impl {
46
51 public:
52 typedef struct Line<3> cell_topology_type;
56 template<EOperator opType>
57 struct Serial {
58 template<typename OutputViewType,
59 typename inputViewType>
60 KOKKOS_INLINE_FUNCTION
61 static void
62 getValues( OutputViewType output,
63 const inputViewType input );
64
65 };
66
67 template<typename DeviceType,
68 typename outputValueValueType, class ...outputValueProperties,
69 typename inputPointValueType, class ...inputPointProperties>
70 static void
71 getValues( const typename DeviceType::execution_space& space,
72 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
73 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
74 const EOperator operatorType);
75
79 template<typename outputValueViewType,
80 typename inputPointViewType,
81 EOperator opType>
82 struct Functor {
83 outputValueViewType _outputValues;
84 const inputPointViewType _inputPoints;
85
86 KOKKOS_INLINE_FUNCTION
87 Functor( outputValueViewType outputValues_,
88 inputPointViewType inputPoints_ )
89 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
90
91 KOKKOS_INLINE_FUNCTION
92 void operator()(const ordinal_type pt) const {
93 switch (opType) {
94 case OPERATOR_VALUE : {
95 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
96 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
97 Serial<opType>::getValues( output, input );
98 break;
99 }
100 case OPERATOR_GRAD :
101 case OPERATOR_MAX : {
102 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
103 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
104 Serial<opType>::getValues( output, input );
105 break;
106 }
107 default: {
108 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
109 opType != OPERATOR_GRAD &&
110 opType != OPERATOR_MAX,
111 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C2_FEM::Serial::getValues) operator is not supported");
112
113 }
114 }
115 }
116 };
117 };
118 }
119
120 template<typename DeviceType = void,
121 typename outputValueType = double,
122 typename pointValueType = double>
124 : public Basis<DeviceType,outputValueType,pointValueType> {
125 public:
127 using typename BasisBase::ExecutionSpace;
128
132
133 using typename BasisBase::OutputViewType;
134 using typename BasisBase::PointViewType ;
135 using typename BasisBase::ScalarViewType;
136
140
142
143 virtual
144 void
146 OutputViewType outputValues,
147 const PointViewType inputPoints,
148 const EOperator operatorType = OPERATOR_VALUE ) const override {
149#ifdef HAVE_INTREPID2_DEBUG
150 // Verify arguments
152 inputPoints,
153 operatorType,
154 this->getBaseCellTopology(),
155 this->getCardinality() );
156#endif
157 Impl::Basis_HGRAD_LINE_C2_FEM::
158 getValues<DeviceType>(space,
159 outputValues,
160 inputPoints,
161 operatorType);
162 }
163
164 virtual void
165 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
166 ordinal_type& perThreadSpaceSize,
167 const PointViewType inputPointsconst,
168 const EOperator operatorType = OPERATOR_VALUE) const override;
169
170 KOKKOS_INLINE_FUNCTION
171 virtual void
172 getValues(
173 OutputViewType outputValues,
174 const PointViewType inputPoints,
175 const EOperator operatorType,
176 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
177 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
178 const ordinal_type subcellDim = -1,
179 const ordinal_type subcellOrdinal = -1) const override;
180
181 virtual
182 void
183 getDofCoords( ScalarViewType dofCoords ) const override {
184#ifdef HAVE_INTREPID2_DEBUG
185 // Verify rank of output array.
186 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
187 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C2_FEM::getDofCoords) rank = 2 required for dofCoords array");
188 // Verify 0th dimension of output array.
189 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
190 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
191 // Verify 1st dimension of output array.
192 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
193 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
194#endif
195 Kokkos::deep_copy(dofCoords, this->dofCoords_);
196 }
197
198 virtual
199 void
200 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
201#ifdef HAVE_INTREPID2_DEBUG
202 // Verify rank of output array.
203 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
204 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
205 // Verify 0th dimension of output array.
206 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
207 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
208#endif
209 Kokkos::deep_copy(dofCoeffs, 1.0);
210 }
211
212 virtual
213 const char*
214 getName() const override {
215 return "Intrepid2_HGRAD_LINE_C2_FEM";
216 }
217
222
223 };
224
225}// namespace Intrepid2
226
228
229#endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HGRAD_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 HGRAD-conforming FEM basis....
Definition file for FEM basis functions of degree 1 for H(grad) functions on a Line.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line 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...
virtual void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
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 getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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...
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::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.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
See Intrepid2::Basis_HGRAD_LINE_C2_FEM.