G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsNewtonCotesRule.h
Go to the documentation of this file.
1 
14 #pragma once
15 
16 #include <gsAssembler/gsQuadRule.h>
17 
18 namespace gismo
19 {
20 
26 template<class T>
27 class gsNewtonCotesRule GISMO_FINAL : public gsQuadRule<T>
28 {
29 public:
30 
31  typedef memory::unique_ptr<gsNewtonCotesRule> uPtr;
32 
35 
38  const unsigned digits = 0 )
39  {
40  gsNewtonCotesRule::setNodes(numNodes, digits);
41  }
42 
44  static uPtr make(gsVector<index_t> const & numNodes,
45  const unsigned digits = 0 )
46  { return uPtr( new gsNewtonCotesRule(numNodes,digits) ); }
47 
49  gsNewtonCotesRule(index_t numNodes, const unsigned digits = 0 )
50  {
51  this->setNodes(numNodes, digits);
52  }
53 
56  gsNewtonCotesRule(const gsBasis<T> & basis, const T quA, const index_t quB, short_t fixDir = -1);
57  //const unsigned digits = std::numeric_limits<T>::digits10 );
58 
62  gsNewtonCotesRule(const gsBasis<T> & basis, const gsOptionList & options, short_t fixDir = -1);
63  //const unsigned digits = std::numeric_limits<T>::digits10 );
64 
65  ~gsNewtonCotesRule() { }
66 
67 public:
68  // see gsQuadRule.h for documentation
69  void setNodes( gsVector<index_t> const & numNodes,
70  unsigned digits = 0 );
71 
72  using gsQuadRule<T>::setNodes; // unhide base
73 
74 private:
75 
76  void init(const gsBasis<T> & basis, const T quA,
77  const index_t quB, short_t fixDir);
78 
84  static void computeReference(index_t n, gsVector<T> & x, gsVector<T> & w);
85 
86 }; // class gsNewtonCotesRule
87 
88 
89 } // namespace gismo
90 
91 
92 #ifndef GISMO_BUILD_LIB
93 #include GISMO_HPP_HEADER(gsNewtonCotesRule.hpp)
94 #endif
Provides a base class for a quadrature rule.
Class representing a reference quadrature rule.
Definition: gsQuadRule.h:28
void setNodes(gsVector< index_t > const &numNodes, unsigned digits=0)
Initialize quadrature rule with numNodes number of quadrature points per integration variable...
Definition: gsNewtonCotesRule.hpp:92
#define short_t
Definition: gsConfig.h:35
#define index_t
Definition: gsConfig.h:32
static void computeReference(index_t n, gsVector< T > &x, gsVector< T > &w)
Computes the Newton-Cotes quadrature rule with n nodes in the interval [-1,1].
Definition: gsNewtonCotesRule.hpp:111
Class that represents the (tensor) Newton-Cotes quadrature rule.
Definition: gsNewtonCotesRule.h:27
gsNewtonCotesRule(gsVector< index_t > const &numNodes, const unsigned digits=0)
Initialize a tensor-product Newton-Cotes quadrature rule with numNodes (direction-wise) ...
Definition: gsNewtonCotesRule.h:37
gsNewtonCotesRule(index_t numNodes, const unsigned digits=0)
Initialize a 1D Newton-Cotes quadrature rule with numNodes.
Definition: gsNewtonCotesRule.h:49
gsNewtonCotesRule()
Default empty constructor.
Definition: gsNewtonCotesRule.h:34
Class which holds a list of parameters/options, and provides easy access to them. ...
Definition: gsOptionList.h:32
index_t numNodes() const
Number of nodes in the currently kept rule.
Definition: gsQuadRule.h:106
A basis represents a family of scalar basis functions defined over a common parameter domain...
Definition: gsBasis.h:78
static uPtr make(gsVector< index_t > const &numNodes, const unsigned digits=0)
Make function returning a smart pointer.
Definition: gsNewtonCotesRule.h:44