Intrepid2
Intrepid2_HGRAD_LINE_C1_FEM.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), or
38// Mauro Perego (mperego@sandia.gov)
39//
40// ************************************************************************
41// @HEADER
42
49#ifndef __INTREPID2_HGRAD_LINE_C1_FEM_HPP__
50#define __INTREPID2_HGRAD_LINE_C1_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
53
54namespace Intrepid2 {
55
78 namespace Impl {
79
84 public:
85 typedef struct Line<2> cell_topology_type;
89 template<EOperator opType>
90 struct Serial {
91 template<typename OutputViewType,
92 typename inputViewType>
93 KOKKOS_INLINE_FUNCTION
94 static void
95 getValues( OutputViewType output,
96 const inputViewType input );
97
98 };
99
100 template<typename DeviceType,
101 typename outputValueValueType, class ...outputValueProperties,
102 typename inputPointValueType, class ...inputPointProperties>
103 static void
104 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
105 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
106 const EOperator operatorType);
107
111 template<typename outputValueViewType,
112 typename inputPointViewType,
113 EOperator opType>
114 struct Functor {
115 outputValueViewType _outputValues;
116 const inputPointViewType _inputPoints;
117
118 KOKKOS_INLINE_FUNCTION
119 Functor( outputValueViewType outputValues_,
120 inputPointViewType inputPoints_ )
121 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
122
123 KOKKOS_INLINE_FUNCTION
124 void operator()(const ordinal_type pt) const {
125 switch (opType) {
126 case OPERATOR_VALUE : {
127 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
128 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
129 Serial<opType>::getValues( output, input );
130 break;
131 }
132 case OPERATOR_GRAD :
133 case OPERATOR_MAX : {
134 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
135 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
136 Serial<opType>::getValues( output, input );
137 break;
138 }
139 default: {
140 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
141 opType != OPERATOR_GRAD &&
142 opType != OPERATOR_MAX,
143 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::Serial::getValues) operator is not supported");
144
145 }
146 }
147 }
148 };
149 };
150 }
151
152 template<typename DeviceType = void,
153 typename outputValueType = double,
154 typename pointValueType = double>
156 : public Basis<DeviceType,outputValueType,pointValueType> {
157 public:
159
160 using OrdinalTypeArray1DHost = typename BasisBase::OrdinalTypeArray1DHost;
161 using OrdinalTypeArray2DHost = typename BasisBase::OrdinalTypeArray2DHost;
162 using OrdinalTypeArray3DHost = typename BasisBase::OrdinalTypeArray3DHost;
163
164 using OutputViewType = typename BasisBase::OutputViewType;
165 using PointViewType = typename BasisBase::PointViewType ;
166 using ScalarViewType = typename BasisBase::ScalarViewType;
167
171
173
174 virtual
175 void
176 getValues( OutputViewType outputValues,
177 const PointViewType inputPoints,
178 const EOperator operatorType = OPERATOR_VALUE ) const override {
179#ifdef HAVE_INTREPID2_DEBUG
180 // Verify arguments
182 inputPoints,
183 operatorType,
184 this->getBaseCellTopology(),
185 this->getCardinality() );
186#endif
187 Impl::Basis_HGRAD_LINE_C1_FEM::
188 getValues<DeviceType>( outputValues,
189 inputPoints,
190 operatorType );
191 }
192
193 virtual
194 void
195 getDofCoords( ScalarViewType dofCoords ) const override {
196#ifdef HAVE_INTREPID2_DEBUG
197 // Verify rank of output array.
198 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
199 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
200 // Verify 0th dimension of output array.
201 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
202 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
203 // Verify 1st dimension of output array.
204 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
205 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
206#endif
207 Kokkos::deep_copy(dofCoords, this->dofCoords_);
208 }
209
210 virtual
211 void
212 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
213#ifdef HAVE_INTREPID2_DEBUG
214 // Verify rank of output array.
215 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
216 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
217 // Verify 0th dimension of output array.
218 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
219 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
220#endif
221 Kokkos::deep_copy(dofCoeffs, 1.0);
222 }
223
224 virtual
225 const char*
226 getName() const override {
227 return "Intrepid2_HGRAD_LINE_C1_FEM";
228 }
229
231 getHostBasis() const override{
233 }
234
235 };
236
237}// namespace Intrepid2
238
240
241#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 1 on Line cell.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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...
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 const char * getName() const override
Returns basis name.
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.
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.
virtual void getValues(OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
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_HGRAD_LINE_C1_FEM.