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 
54 namespace 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( const typename DeviceType::execution_space& space,
105  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
106  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
107  const EOperator operatorType);
108 
112  template<typename outputValueViewType,
113  typename inputPointViewType,
114  EOperator opType>
115  struct Functor {
116  outputValueViewType _outputValues;
117  const inputPointViewType _inputPoints;
118 
119  KOKKOS_INLINE_FUNCTION
120  Functor( outputValueViewType outputValues_,
121  inputPointViewType inputPoints_ )
122  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
123 
124  KOKKOS_INLINE_FUNCTION
125  void operator()(const ordinal_type pt) const {
126  switch (opType) {
127  case OPERATOR_VALUE : {
128  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
129  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
130  Serial<opType>::getValues( output, input );
131  break;
132  }
133  case OPERATOR_GRAD :
134  case OPERATOR_MAX : {
135  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
136  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
137  Serial<opType>::getValues( output, input );
138  break;
139  }
140  default: {
141  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
142  opType != OPERATOR_GRAD &&
143  opType != OPERATOR_MAX,
144  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::Serial::getValues) operator is not supported");
145 
146  }
147  }
148  }
149  };
150  };
151  }
152 
153  template<typename DeviceType = void,
154  typename outputValueType = double,
155  typename pointValueType = double>
157  : public Basis<DeviceType,outputValueType,pointValueType> {
158  public:
160  using typename BasisBase::ExecutionSpace;
161 
162  using typename BasisBase::OrdinalTypeArray1DHost;
163  using typename BasisBase::OrdinalTypeArray2DHost;
164  using typename BasisBase::OrdinalTypeArray3DHost;
165 
166  using typename BasisBase::OutputViewType;
167  using typename BasisBase::PointViewType ;
168  using typename BasisBase::ScalarViewType;
169 
173 
174  using BasisBase::getValues;
175 
176  virtual
177  void
178  getValues( const ExecutionSpace& space,
179  OutputViewType outputValues,
180  const PointViewType inputPoints,
181  const EOperator operatorType = OPERATOR_VALUE ) const override {
182 #ifdef HAVE_INTREPID2_DEBUG
183  // Verify arguments
184  Intrepid2::getValues_HGRAD_Args(outputValues,
185  inputPoints,
186  operatorType,
187  this->getBaseCellTopology(),
188  this->getCardinality() );
189 #endif
190  Impl::Basis_HGRAD_LINE_C1_FEM::
191  getValues<DeviceType>(space,
192  outputValues,
193  inputPoints,
194  operatorType );
195  }
196 
197  virtual
198  void
199  getDofCoords( ScalarViewType dofCoords ) const override {
200 #ifdef HAVE_INTREPID2_DEBUG
201  // Verify rank of output array.
202  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
203  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
204  // Verify 0th dimension of output array.
205  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
206  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
207  // Verify 1st dimension of output array.
208  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
209  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
210 #endif
211  Kokkos::deep_copy(dofCoords, this->dofCoords_);
212  }
213 
214  virtual
215  void
216  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
217 #ifdef HAVE_INTREPID2_DEBUG
218  // Verify rank of output array.
219  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
220  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
221  // Verify 0th dimension of output array.
222  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
223  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
224 #endif
225  Kokkos::deep_copy(dofCoeffs, 1.0);
226  }
227 
228  virtual
229  const char*
230  getName() const override {
231  return "Intrepid2_HGRAD_LINE_C1_FEM";
232  }
233 
234  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
235  getHostBasis() const override{
237  }
238 
239  };
240 
241 }// namespace Intrepid2
242 
244 
245 #endif
See Intrepid2::Basis_HGRAD_LINE_C1_FEM.
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...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Definition file for FEM basis functions of degree 1 for H(grad) functions on a Line.
ordinal_type getCardinality() const
Returns cardinality of the basis.
virtual void getValues(const ExecutionSpace &, OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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.
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...
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
virtual const char * getName() const override
Returns basis name.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Line cell.
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.