Intrepid
Intrepid_HGRAD_LINE_Hermite_FEM.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid 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 Pavel Bochev (pbboche@sandia.gov)
38 // Denis Ridzal (dridzal@sandia.gov), or
39 // Kara Peterson (kjpeter@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
49 #ifndef INTREPID_HGRAD_LINE_HERMITE_FEM_HPP
50 #define INTREPID_HGRAD_LINE_HERMITE_FEM_HPP
51 
52 #include "Intrepid_Basis.hpp"
53 #include "Intrepid_PointTools.hpp"
54 #include "Intrepid_Polylib.hpp"
55 
56 #include "Teuchos_SerialDenseMatrix.hpp"
57 #include "Teuchos_SerialDenseSolver.hpp"
58 #include <ostream>
59 #include <cstdint>
60 
61 namespace Intrepid {
62 
92 template<class Scalar, class ArrayScalar>
93 class Basis_HGRAD_LINE_Hermite_FEM: public Basis<Scalar, ArrayScalar> {
94 
95  template<typename T> using RCP = Teuchos::RCP<T>;
96  using SerialDenseMatrix = Teuchos::SerialDenseMatrix<int,Scalar>;
97 
98  using TAGS = std::vector<int>;
99 
100 private:
101 
104 
106  mutable SerialDenseMatrix V_;
107 
108  mutable Teuchos::SerialDenseSolver<int,Scalar> solver_;
109 
110  mutable TAGS tags_;
111 
112  mutable bool isFactored_;
113 
116  void initializeTags();
117 
120  void recurrence( ArrayScalar &P, ArrayScalar &Px, const Scalar x ) const;
121 
125  void legendre_d( ArrayScalar &Pm, ArrayScalar &Pm1, const int m, const Scalar pt ) const;
126 
129  void setupVandermonde( bool factor=true );
130 
131 
132 public:
133 
138 
141  Basis_HGRAD_LINE_Hermite_FEM( const ArrayScalar &pts );
142 
146  Basis_HGRAD_LINE_Hermite_FEM( int numPoints , const EPointType &pointType );
147 
148 
159  void getValues(ArrayScalar & outputValues,
160  const ArrayScalar & inputPoints,
161  const EOperator operatorType) const;
162 
163 
166  void getValues(ArrayScalar & outputValues,
167  const ArrayScalar & inputPoints,
168  const ArrayScalar & cellVertices,
169  const EOperator operatorType = OPERATOR_VALUE) const;
170 
172  virtual void getDofCoords( ArrayScalar & DofCoords ) const;
173 
174 
175  void printTags( std::ostream &os );
176 
177 }; // class Basis_HGRAD_LINE_Hermite_FEM
178 
179 }// namespace Intrepid
180 
182 
183 #endif // INTREPID_HGRAD_LINE_HERMITE_FEM_HPP
virtual void getDofCoords(ArrayScalar &DofCoords) const
implements the dofcoords interface
Implements Hermite interpolant basis of degree n on the reference Line cell. The basis has cardinalit...
FieldContainer< Scalar > latticePts_
Holds the points defining the Hermite basis.
Header file for utility class to provide point tools, such as barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.
void legendre_d(ArrayScalar &Pm, ArrayScalar &Pm1, const int m, const Scalar pt) const
Evaluates and at a particular point .
void recurrence(ArrayScalar &P, ArrayScalar &Px, const Scalar x) const
Evaluates and at a particular point .
Definition file for Hermite FEM basis functions of degree 2n for H(grad) functions on a Line...
SerialDenseMatrix V_
Contains the values of the Legendre polynomials and their derivatives.
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Line cell.
void initializeTags()
Initializes tagToOrdinal_ and ordinalToTag_ lookup arrays.
Header file for the abstract base class Intrepid::Basis.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
void setupVandermonde(bool factor=true)
Form the Legendre/Derivative Vandermonde matrix at the given lattice points and have the linear solve...
Basis_HGRAD_LINE_Hermite_FEM()
Default Constructor assumes the two interpolation points are the cell vertices. Cubic Hermite Interpo...
Header file for a set of functions providing orthogonal polynomial polynomial calculus and interpolat...