quadrature_build.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2014 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 // C++ includes
20 
21 
22 
23 // Local includes
29 #include "libmesh/string_to_enum.h"
30 
31 namespace libMesh
32 {
33 
34 
35 
36 //---------------------------------------------------------------
37 AutoPtr<QBase> QBase::build (const std::string &type,
38  const unsigned int _dim,
39  const Order _order)
40 {
41  return QBase::build (Utility::string_to_enum<QuadratureType> (type),
42  _dim,
43  _order);
44 }
45 
46 
47 
49  const unsigned int _dim,
50  const Order _order)
51 {
52  switch (_qt)
53  {
54 
55  case QCLOUGH:
56  {
57 #ifdef DEBUG
58  if (_order > TWENTYTHIRD)
59  {
60  libMesh::out << "WARNING: Clough quadrature implemented" << std::endl
61  << " up to TWENTYTHIRD order." << std::endl;
62  }
63 #endif
64 
65  AutoPtr<QBase> ap(new QClough(_dim, _order));
66  return ap;
67  }
68 
69  case QGAUSS:
70  {
71 
72 #ifdef DEBUG
73  if (_order > FORTYTHIRD)
74  {
75  libMesh::out << "WARNING: Gauss quadrature implemented" << std::endl
76  << " up to FORTYTHIRD order." << std::endl;
77  }
78 #endif
79 
80  AutoPtr<QBase> ap(new QGauss(_dim, _order));
81  return ap;
82  }
83 
84  case QJACOBI_1_0:
85  {
86 
87 #ifdef DEBUG
88  if (_order > TWENTYTHIRD)
89  {
90  libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl
91  << " up to TWENTYTHIRD order." << std::endl;
92  }
93 
94  if (_dim > 1)
95  {
96  libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl
97  << " in 1D only." << std::endl;
98  }
99 #endif
100 
101  AutoPtr<QBase> ap(new QJacobi(_dim, _order, 1, 0));
102  return ap;
103  }
104 
105  case QJACOBI_2_0:
106  {
107 
108 #ifdef DEBUG
109  if (_order > TWENTYTHIRD)
110  {
111  libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl
112  << " up to TWENTYTHIRD order." << std::endl;
113  }
114 
115  if (_dim > 1)
116  {
117  libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl
118  << " in 1D only." << std::endl;
119  }
120 #endif
121 
122  AutoPtr<QBase> ap(new QJacobi(_dim, _order, 2, 0));
123  return ap;
124  }
125 
126  case QSIMPSON:
127  {
128 
129 #ifdef DEBUG
130  if (_order > THIRD)
131  {
132  libMesh::out << "WARNING: Simpson rule provides only" << std::endl
133  << " THIRD order!" << std::endl;
134  }
135 #endif
136 
137  AutoPtr<QBase> ap(new QSimpson(_dim));
138  return ap;
139  }
140 
141  case QTRAP:
142  {
143 
144 #ifdef DEBUG
145  if (_order > FIRST)
146  {
147  libMesh::out << "WARNING: Trapezoidal rule provides only" << std::endl
148  << " FIRST order!" << std::endl;
149  }
150 #endif
151 
152  AutoPtr<QBase> ap(new QTrap(_dim));
153  return ap;
154  }
155 
156 
157  default:
158  {
159  libMesh::err << "ERROR: Bad qt=" << _qt << std::endl;
160  libmesh_error();
161  }
162  }
163 
164 
165  libmesh_error();
166  AutoPtr<QBase> ap(NULL);
167  return ap;
168 }
169 
170 } // namespace libMesh

Site Created By: libMesh Developers
Last modified: February 07 2014 16:57:06 UTC

Hosted By:
SourceForge.net Logo