libMesh
quadrature_build.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 // Local includes
20 #include "libmesh/quadrature_clough.h"
21 #include "libmesh/quadrature_gauss.h"
22 #include "libmesh/quadrature_gm.h"
23 #include "libmesh/quadrature_grid.h"
24 #include "libmesh/quadrature_jacobi.h"
25 #include "libmesh/quadrature_monomial.h"
26 #include "libmesh/quadrature_simpson.h"
27 #include "libmesh/quadrature_trap.h"
28 #include "libmesh/quadrature_gauss_lobatto.h"
29 #include "libmesh/quadrature_conical.h"
30 #include "libmesh/quadrature_nodal.h"
31 #include "libmesh/string_to_enum.h"
32 #include "libmesh/enum_quadrature_type.h"
33 
34 // C++ Includes
35 #include <memory>
36 
37 namespace libMesh
38 {
39 
40 
41 
42 //---------------------------------------------------------------
43 std::unique_ptr<QBase> QBase::build (std::string_view type,
44  const unsigned int _dim,
45  const Order _order)
46 {
47  return QBase::build (Utility::string_to_enum<QuadratureType> (type),
48  _dim,
49  _order);
50 }
51 
52 
53 
54 std::unique_ptr<QBase> QBase::build(const QuadratureType _qt,
55  const unsigned int _dim,
56  const Order _order)
57 {
58  switch (_qt)
59  {
60 
61  case QCLOUGH:
62  {
63 #ifdef DEBUG
64  if (_order > TWENTYTHIRD)
65  {
66  libMesh::out << "WARNING: Clough quadrature implemented" << std::endl
67  << " up to TWENTYTHIRD order." << std::endl;
68  }
69 #endif
70 
71  return std::make_unique<QClough>(_dim, _order);
72  }
73 
74  case QGAUSS:
75  {
76 
77 #ifdef DEBUG
78  if (_order > FORTYTHIRD)
79  {
80  libMesh::out << "WARNING: Gauss quadrature implemented" << std::endl
81  << " up to FORTYTHIRD order." << std::endl;
82  }
83 #endif
84 
85  return std::make_unique<QGauss>(_dim, _order);
86  }
87 
88  case QJACOBI_1_0:
89  {
90 
91 #ifdef DEBUG
92  if (_order > FORTYTHIRD)
93  {
94  libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl
95  << " up to FORTYTHIRD order." << std::endl;
96  }
97 
98  if (_dim > 1)
99  {
100  libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl
101  << " in 1D only." << std::endl;
102  }
103 #endif
104 
105  return std::make_unique<QJacobi>(_dim, _order, 1, 0);
106  }
107 
108  case QJACOBI_2_0:
109  {
110 
111 #ifdef DEBUG
112  if (_order > FORTYTHIRD)
113  {
114  libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl
115  << " up to FORTYTHIRD order." << std::endl;
116  }
117 
118  if (_dim > 1)
119  {
120  libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl
121  << " in 1D only." << std::endl;
122  }
123 #endif
124 
125  return std::make_unique<QJacobi>(_dim, _order, 2, 0);
126  }
127 
128  case QSIMPSON:
129  {
130 
131 #ifdef DEBUG
132  if (_order > THIRD)
133  {
134  libMesh::out << "WARNING: Simpson rule provides only" << std::endl
135  << " THIRD order!" << std::endl;
136  }
137 #endif
138 
139  return std::make_unique<QSimpson>(_dim);
140  }
141 
142  case QTRAP:
143  {
144 
145 #ifdef DEBUG
146  if (_order > FIRST)
147  {
148  libMesh::out << "WARNING: Trapezoidal rule provides only" << std::endl
149  << " FIRST order!" << std::endl;
150  }
151 #endif
152 
153  return std::make_unique<QTrap>(_dim);
154  }
155 
156  case QGRID:
157  return std::make_unique<QGrid>(_dim, _order);
158 
159  case QGRUNDMANN_MOLLER:
160  return std::make_unique<QGrundmann_Moller>(_dim, _order);
161 
162  case QMONOMIAL:
163  return std::make_unique<QMonomial>(_dim, _order);
164 
165  case QGAUSS_LOBATTO:
166  return std::make_unique<QGaussLobatto>(_dim, _order);
167 
168  case QCONICAL:
169  return std::make_unique<QConical>(_dim, _order);
170 
171  case QNODAL:
172  return std::make_unique<QNodal>(_dim, _order);
173 
174  default:
175  libmesh_error_msg("ERROR: Bad qt=" << _qt);
176  }
177 }
178 
179 } // namespace libMesh
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
QuadratureType
Defines an enum for currently available quadrature rules.
unsigned int _dim
The spatial dimension of the quadrature rule.
Definition: quadrature.h:379
The libMesh namespace provides an interface to certain functionality in the library.
virtual QuadratureType type() const =0
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly. ...
Definition: quadrature.h:385
OStreamProxy out
static std::unique_ptr< QBase > build(std::string_view name, const unsigned int dim, const Order order=INVALID_ORDER)
Builds a specific quadrature rule based on the name string.