Intrepid2
Intrepid2_NodalBasisFamily.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 
49 #ifndef Intrepid2_NodalBasisFamily_h
50 #define Intrepid2_NodalBasisFamily_h
51 
53 
56 
61 
66 
71 
76 
77 namespace Intrepid2 {
78  // the following defines a family of nodal basis functions, derived from a the standard high-order Intrepid2 bases on the line
79  // note that because the standard H(curl) basis uses a lower-order H(grad) basis in place of the H(vol) that is arguably more natural,
80  // the following will *not* match the standard nodal basis declared below. (Similarly for H(div).)
81 
89  template<typename ExecutionSpace=Kokkos::DefaultExecutionSpace,
90  typename OutputScalar = double,
91  typename PointScalar = double>
92  using DerivedNodalBasisFamily = DerivedBasisFamily< Basis_HGRAD_LINE_Cn_FEM<ExecutionSpace,OutputScalar,PointScalar>,
93  Basis_HVOL_LINE_Cn_FEM <ExecutionSpace,OutputScalar,PointScalar> >;
94 
102  template<typename ExecSpace=Kokkos::DefaultExecutionSpace,
103  typename OutputScalar = double,
104  typename PointScalar = double>
106  {
107  public:
108  using ExecutionSpace = ExecSpace;
109  using OutputValueType = OutputScalar;
110  using PointValueType = PointScalar;
111 
113  using BasisPtr = Teuchos::RCP<BasisType>;
114 
115  // line bases
118 
119  // quadrilateral bases
124 
125  // triangle bases
130 
131  // hexahedron bases
136 
137  // tetrahedron bases
142 
143  };
144 }
145 
146 #endif /* Intrepid2_NodalBasisFamily_h */
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Triangle cell...
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Header file for the Intrepid2::Basis_HCURL_HEX_In_FEM class.
Implementation of the default H(div)-compatible FEM basis on Quadrilateral cell.
Implementation of the default H(curl)-compatible FEM basis on Quadrilateral cell. ...
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Te...
Stateless class representing a family of basis functions, templated on H(vol) and H(grad) on the line...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Header file for the Intrepid2::Basis_HDIV_HEX_In_FEM class.
Implementation of the default H(curl)-compatible FEM basis on Hexahedron cell.
Header file for the Intrepid2::Basis_HVOL_HEX_Cn_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell...
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM class.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Header file for the Intrepid2::Basis_HCURL_TET_In_FEM class.
Header file for the Intrepid2::Basis_HVOL_TET_Cn_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM class.
Header file for the Intrepid2::Basis_HDIV_QUAD_In_FEM class.
Header file for the Intrepid2::Basis_HCURL_QUAD_In_FEM class.
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Tetrahedr...
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Tetrahedron ce...
Header file for the Intrepid2::Basis_HDIV_TRI_In_FEM class.
Header file for the Intrepid2::Basis_HVOL_QUAD_Cn_FEM class.
Implementation of the default HVOL-compatible FEM basis of degree n on Quadrilateral cell Implements ...
Header file for the Intrepid2::Basis_HGRAD_QUAD_Cn_FEM class.
A family of nodal basis functions representing the higher-order Lagrangian basis family that Intrepid...
Header file for the Intrepid2::Basis_HDIV_TET_In_FEM class.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Tr...
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Triangle ...
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell...
Implementation of the default H(div)-compatible FEM basis on Hexahedron cell.
Implementation of the default HVOL-compatible FEM basis of degree n on Hexahedron cell...
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Tetrahedron cell...
Implementation of the default H(grad)-compatible FEM basis of degree n on Quadrilateral cell Implemen...
Header file for the Intrepid2::Basis_HVOL_TRI_Cn_FEM class.
Header file for the Intrepid2::Basis_HCURL_TRI_In_FEM class.
Header file for the Intrepid2::Basis_HGRAD_HEX_Cn_FEM class.