libMesh
inf_fe_radial_test.C
Go to the documentation of this file.
1 // libmesh includes
2 #include "libmesh/libmesh.h"
3 #include "libmesh/replicated_mesh.h"
4 #include "libmesh/mesh_generation.h"
5 #include "libmesh/inf_elem_builder.h"
6 #include "libmesh/fe_type.h"
7 #include "libmesh/quadrature_gauss.h"
8 #include "libmesh/elem.h"
9 #include "libmesh/inf_fe.h"
10 
11 // unit test includes
12 #include "test_comm.h"
13 
14 // C++ includes
15 #include <map>
16 #include <vector>
17 
18 #include "libmesh_cppunit.h"
19 
20 using namespace libMesh;
21 
32 class InfFERadialTest : public CppUnit::TestCase
33 {
34 public:
35  CPPUNIT_TEST_SUITE( InfFERadialTest );
36  CPPUNIT_TEST( testDifferentOrders );
37  CPPUNIT_TEST_SUITE_END();
38 
39 public:
40  struct TabulatedVal
41  {
42  TabulatedVal(unsigned int i_in, unsigned int qp_in, Real val_in)
43  : i(i_in), qp(qp_in), val(val_in) {}
44  unsigned int i;
45  unsigned int qp;
47  };
48 
50  {
51  TabulatedGrad(unsigned int i_in, unsigned int qp_in, Point grad_in)
52  : i(i_in), qp(qp_in), grad(grad_in) {}
53  unsigned int i;
54  unsigned int qp;
56  };
57 
58  typedef std::vector<TabulatedVal> TabulatedVals;
59  typedef std::vector<TabulatedGrad> TabulatedGrads;
60  typedef std::pair<Order, FEFamily> KeyType;
61  std::map<KeyType, TabulatedVals> tabulated_vals;
62  std::map<KeyType, TabulatedGrads> tabulated_grads;
63 
65  {
66  testSingleOrder(FIRST, LEGENDRE);
67  testSingleOrder(SECOND, LEGENDRE);
68  testSingleOrder(THIRD, LEGENDRE);
69  testSingleOrder(FOURTH, LEGENDRE);
70 
71  testSingleOrder(FIRST, JACOBI_20_00);
72  testSingleOrder(SECOND, JACOBI_20_00);
73  testSingleOrder(THIRD, JACOBI_20_00);
74  testSingleOrder(FOURTH, JACOBI_20_00);
75 
76  testSingleOrder(FIRST, JACOBI_30_00);
77  testSingleOrder(SECOND, JACOBI_30_00);
78  testSingleOrder(THIRD, JACOBI_30_00);
79  testSingleOrder(FOURTH, JACOBI_30_00);
80  }
81 
82  void testSingleOrder (Order radial_order,
83  FEFamily radial_family)
84  {
85  // Avoid warnings when not compiled with infinite element support.
86  libmesh_ignore(radial_order, radial_family);
87 
88 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
89  // std::cout << "Called testSingleOrder with radial_order = "
90  // << radial_order
91  // << std::endl;
92 
95  (mesh,
96  /*nx=*/1, /*ny=*/1, /*nz=*/1,
97  /*xmin=*/-1., /*xmax=*/1.,
98  /*ymin=*/-1., /*ymax=*/1.,
99  /*zmin=*/-1., /*zmax=*/1.,
100  HEX8);
101 
102  // Add infinite elements to the mesh. The optional verbose flag only
103  // prints extra information in non-optimized mode.
104  InfElemBuilder builder(mesh);
105  builder.build_inf_elem(/*verbose=true*/);
106 
107  // Get pointer to the last infinite Elem. I originally intended to
108  // get a pointer to the *first* infinite Elem but then I computed
109  // all the reference values on the last Elem, so that's the one we
110  // are using now...
111  const Elem * infinite_elem = mesh.elem_ptr(mesh.n_elem() - 1);
112  if (!infinite_elem || !infinite_elem->infinite())
113  libmesh_error_msg("Error setting Elem pointer.");
114 
115  // We will construct FEs, etc. of the same dimension as the mesh elements.
116  auto dim = mesh.mesh_dimension();
117 
118  FEType fe_type(/*Order*/FIRST,
119  /*FEFamily*/LAGRANGE,
120  radial_order,
121  radial_family,
122  /*inf_map*/CARTESIAN);
123 
124  // Construct FE, quadrature rule, etc.
125  std::unique_ptr<FEBase> inf_fe (FEBase::build_InfFE(dim, fe_type));
126  QGauss qrule (dim, fe_type.default_quadrature_order());
127  inf_fe->attach_quadrature_rule(&qrule);
128  const std::vector<std::vector<Real>> & phi = inf_fe->get_phi();
129  const std::vector<std::vector<RealGradient>> & dphi = inf_fe->get_dphi();
130 
131  // Reinit on infinite elem
132  inf_fe->reinit(infinite_elem);
133 
134  // Check phi values against reference values. The number of
135  // quadrature points and shape functions seem to both increase by
136  // 4 as the radial_order is increased.
137  //
138  // radial_order | n_qp | n_sf
139  // --------------------------
140  // FIRST | 16 | 8
141  // SECOND | 20 | 12
142  // THIRD | 24 | 16
143  // FOURTH | 28 | 20
144  // FIFTH | 32 | 24
145 
146  // Debugging print values so we can tabulate them
147  // auto n_qp = inf_fe->n_quadrature_points();
148  // auto n_sf = inf_fe->n_shape_functions();
149  // for (unsigned int i=0; i<n_sf; ++i)
150  // for (unsigned qp=0; qp<n_qp; ++qp)
151  // {
152  // libMesh::out << "{" << i << "," << qp << "," << phi[i][qp] << "}," << std::endl;
153  // }
154  // for (unsigned int i=0; i<n_sf; ++i)
155  // for (unsigned qp=0; qp<n_qp; ++qp)
156  // {
157  // libMesh::out << "{" << i << "," << qp << ",Point("
158  // << dphi[i][qp](0)
159  // << ","
160  // << dphi[i][qp](1)
161  // << ","
162  // << dphi[i][qp](2)
163  // << ")},"
164  // << std::endl;
165  // }
166 
167  // Get reference vals for this Order/Family combination.
168  const TabulatedVals & val_table =
169  tabulated_vals[std::make_pair(radial_order, radial_family)];
170  const TabulatedGrads & grad_table =
171  tabulated_grads[std::make_pair(radial_order, radial_family)];
172 
173  // Test whether the computed values match the tabulated values to
174  // the specified accuracy.
175  for (const auto & t : val_table)
176  LIBMESH_ASSERT_FP_EQUAL(t.val, phi[t.i][t.qp], 1.e-4);
177  for (const auto & t : grad_table)
178  LIBMESH_ASSERT_FP_EQUAL(0., (dphi[t.i][t.qp] - t.grad).norm_sq(), 1.e-4);
179 
180  // Make sure there are actually reference values
181  if (val_table.empty())
182  libmesh_error_msg("No tabulated values found!");
183  if (grad_table.empty())
184  libmesh_error_msg("No tabulated gradients found!");
185 
186 #endif // LIBMESH_ENABLE_INFINITE_ELEMENTS
187  }
188 
189  void setUp()
190  {
192  // Arbitrarily selected LEGENDRE reference values.
194  tabulated_vals[std::make_pair(FIRST, LEGENDRE)] =
195  {
196  {0,9,0.0550016}, {1,5,0.41674}, {2,8,0.0147376}, {3,7,0.111665},
197  {4,6,0.0737011}, {5,1,0.0803773}, {6,11,0.275056}, {7,15,0.021537}
198  };
199 
200  tabulated_vals[std::make_pair(SECOND, LEGENDRE)] =
201  {
202  {0,3,0.0425633}, {1,6,0.0343526}, {2,2,0.158848}, {3,7,0.128206},
203  {4,12,0.220829}, {5,5,-0.136549}, {6,0,0.0149032}, {7,10,-0.0334936},
204  {8,16,0.00399329}, {9,18,-0.00209733}, {10,2,0.0556194}, {11,17,-0.000561977}
205  };
206 
207  tabulated_vals[std::make_pair(THIRD, LEGENDRE)] =
208  {
209  {12,9,0.136657}, {13,6,0.175034}, {14,20,-0.0011016}, {15,15,0.0428935}
210  };
211 
212  tabulated_vals[std::make_pair(FOURTH, LEGENDRE)] =
213  {
214  {16,3,0.00826618}, {17,10,-0.547813}, {18,14,0.311004}, {19,27,-0.00192085}
215  };
216 
218  // Arbitrarily selected LEGENDRE reference gradients.
220  tabulated_grads[std::make_pair(FIRST, LEGENDRE)] =
221  {
222  {0,9,Point(-0.0429458, -0.0115073, 0.)},
223  {1,5,Point(0.177013, -0.177013, -0.483609)},
224  {2,8,Point(0.0115073, 0.0115073, 0.00842393)},
225  {3,7,Point(-0.177013, 0.0474305, 0.)},
226  {4,6,Point(-0.031305, -0.116832, 0.10025)},
227  {5,1,Point(0.0474191, -0.0474191, 0.872917)},
228  {6,11,Point(0.0575466, 0.0575466, -0.11251)},
229  {7,15,Point(-0.00353805, 0.000948017, 0.000111572)}
230  };
231 
232  tabulated_grads[std::make_pair(SECOND, LEGENDRE)] =
233  {
234  {0,3,Point(-0.0959817, -0.0959817, 0.0702635)},
235  {1,6,Point(0.0625228, -0.0625228, 0.0457699)},
236  {2,2,Point(0.358209, 0.0959817, -2.77556e-17)},
237  {3,7,Point(-0.233338, 0.0625228, -6.93889e-17)},
238  {4,12,Point(-0.0323071, -0.0323071, -0.0729771)},
239  {5,5,Point(0.248523, 0.0665915, -0.245097)},
240  {6,0,Point(0.0336072, -0.00900502, 0.288589)},
241  {7,10,Point(-0.0396234, 0.0396234, -0.0290064)},
242  {8,16,Point(0.000443217, 0.000443217, 0.000333678)},
243  {9,18,Point(-0.000232783, -6.23741e-05, 9.35433e-05)},
244  {10,2,Point(-0.0336072, 0.0336072, 0.985214)},
245  {11,17,Point(6.23741e-05, -6.23741e-05, -2.05961e-05)},
246  };
247 
248  tabulated_grads[std::make_pair(THIRD, LEGENDRE)] =
249  {
250  {12,9,Point(0.0536552, 0.200244, -0.0849541)},
251  {13,6,Point(-0.0921697, 0.0921697, 0.461056)},
252  {14,20,Point(2.35811e-05, -8.80059e-05, 3.58959e-05)},
253  {15,15,Point(-0.0386352, 0.0103523, -0.0197323)},
254  };
255 
256  tabulated_grads[std::make_pair(FOURTH, LEGENDRE)] =
257  {
258  {16,3,Point(-0.0190603, 0.0051072, 0.308529)},
259  {17,10,Point(0.244125, -0.244125, 0.140907)},
260  {18,14,Point(-0.0985844, 0.0985844, -0.502591)},
261  {19,27,Point(0.000115647, -3.09874e-05, 4.30775e-05)}
262  };
263 
265  // Arbitrarily selected JACOBI_20_00 reference values.
267  tabulated_vals[std::make_pair(FIRST, JACOBI_20_00)] =
268  {
269  // These values are the same as Legendre
270  {0,9,0.0550016}, {1,5,0.41674}, {2,8,0.0147376}, {3,7,0.111665},
271  // These values are different from Legendre
272  {4,6,0.147402}, {5,1,0.160755}, {6,11,0.550112}, {7,15,0.043074}
273  };
274 
275  tabulated_vals[std::make_pair(SECOND, JACOBI_20_00)] =
276  {
277  // These values are the same as Legendre
278  {0,3,0.0425633}, {1,6,0.0343526}, {2,2,0.158848}, {3,7,0.128206},
279  // These values are different from Legendre
280  {4,12,0.441658}, {5,5,-0.193445}, {6,0,0.0298063}, {7,10,-0.0279114},
281  {8,16,0.00798659}, {9,18,0.0320146}, {10,2,0.111239}, {11,17,0.00857829}
282  };
283 
284  tabulated_vals[std::make_pair(THIRD, JACOBI_20_00)] =
285  {
286  {12,9,0.0837876}, {13,6,0.350068}, {14,20,0.0244336}, {15,15,0.0181532}
287  };
288 
289  tabulated_vals[std::make_pair(FOURTH, JACOBI_20_00)] =
290  {
291  {16,3,0.0165324}, {17,10,-0.720085}, {18,14,0.0777511}, {19,27,0.0453842}
292  };
293 
295  // Arbitrarily selected JACOBI_20_00 reference gradients.
297  tabulated_grads[std::make_pair(FIRST, JACOBI_20_00)] =
298  {
299  // These values are the same as Legendre
300  {0,9,Point(-0.0429458, -0.0115073, 0.)},
301  {1,5,Point(0.177013, -0.177013, -0.483609)},
302  {2,8,Point(0.0115073, 0.0115073, 0.00842393)},
303  {3,7,Point(-0.177013, 0.0474305, 0.)},
304  // These values are different from Legendre
305  {4,6,Point(-0.0626101, -0.233664, 0.2005)},
306  {5,1,Point(0.0948382, -0.0948382, 1.74583)},
307  {6,11,Point(0.115093, 0.115093, -0.22502)},
308  {7,15,Point(-0.0070761, 0.00189603, 0.000223144)}
309  };
310 
311  tabulated_grads[std::make_pair(SECOND, JACOBI_20_00)] =
312  {
313  // These values are the same as Legendre
314  {0,3,Point(-0.0959817, -0.0959817, 0.0702635)},
315  {1,6,Point(0.0625228, -0.0625228, 0.0457699)},
316  {2,2,Point(0.358209, 0.0959817, -2.77556e-17)},
317  {3,7,Point(-0.233338, 0.0625228, -6.93889e-17)},
318  // These values are different from Legendre
319  {4,12,Point(-0.0646142, -0.0646142, -0.145954)},
320  {5,5,Point(0.352076, 0.0943384, -0.233431)},
321  {6,0,Point(0.0672144, -0.01801, 0.577179)},
322  {7,10,Point(-0.0330195, 0.0330195, 0.00373942)},
323  {8,16,Point(0.000886435, 0.000886435, 0.000667355)},
324  {9,18,Point(0.00355332, 0.000952108, 0.000319882)},
325  {10,2,Point(-0.0672144, 0.0672144, 1.97043)},
326  {11,17,Point(-0.000952108, 0.000952108, 0.000782704)}
327  };
328 
329  tabulated_grads[std::make_pair(THIRD, JACOBI_20_00)] =
330  {
331  {12,9,Point(0.0328972,0.122774,-0.222472)},
332  {13,6,Point(-0.184339,0.184339,0.922112)},
333  {14,20,Point(-0.000523034,0.00195199,0.000121819)},
334  {15,15,Point(-0.016351,0.00438123,0.0404817)}
335  };
336 
337  tabulated_grads[std::make_pair(FOURTH, JACOBI_20_00)] =
338  {
339  {16,3,Point(-0.0381206,0.0102144,0.617059)},
340  {17,10,Point(0.320895,-0.320895,0.641729)},
341  {18,14,Point(-0.0246461,0.0246461,-0.300588)},
342  {19,27,Point(-0.0027324,0.000732145,0.000328402)}
343  };
344 
346  // Arbitrarily selected JACOBI_30_00 reference values.
348  tabulated_vals[std::make_pair(FIRST, JACOBI_30_00)] =
349  {
350  // These values are the same as Legendre
351  {0,9,0.0550016}, {1,5,0.41674}, {2,8,0.0147376}, {3,7,0.111665},
352  // These values are different from Legendre
353  {4,6,0.184253}, {5,1,0.200943}, {6,11,0.68764}, {7,15,0.0538426}
354  };
355 
356  tabulated_vals[std::make_pair(SECOND, JACOBI_30_00)] =
357  {
358  // These values are the same as Legendre
359  {0,3,0.0425633}, {1,6,0.0343526}, {2,2,0.158848}, {3,7,0.128206},
360  // These values are different from Legendre
361  {4,12,0.552072}, {5,5,-0.211652}, {6,0,0.0372579}, {7,10,-0.0167468},
362  {8,16,0.00998323}, {9,18,0.0597236}, {10,2,0.139049}, {11,17,0.0160029}
363  };
364 
365  tabulated_vals[std::make_pair(THIRD, JACOBI_30_00)] =
366  {
367  {12,9,0.0469849}, {13,6,0.437585}, {14,20,0.0450821}, {15,15,0.0469849}
368  };
369 
370  tabulated_vals[std::make_pair(FOURTH, JACOBI_30_00)] =
371  {
372  {16,3,0.0206655}, {17,10,-0.748341}, {18,14,0}, {19,27,0.115995}
373  };
374 
376  // Arbitrarily selected JACOBI_30_00 reference gradients.
378  tabulated_grads[std::make_pair(FIRST, JACOBI_30_00)] =
379  {
380  // These values are the same as Legendre
381  {0,9,Point(-0.0429458, -0.0115073, 0.)},
382  {1,5,Point(0.177013, -0.177013, -0.483609)},
383  {2,8,Point(0.0115073, 0.0115073, 0.00842393)},
384  {3,7,Point(-0.177013, 0.0474305, 0.)},
385  // These values are different from Legendre
386  {4,6,Point(-0.0782626,-0.29208,0.250625)},
387  {5,1,Point(0.118548,-0.118548,2.18229)},
388  {6,11,Point(0.143866,0.143866,-0.281275)},
389  {7,15,Point(-0.00884512,0.00237004,0.00027893)}
390  };
391 
392  tabulated_grads[std::make_pair(SECOND, JACOBI_30_00)] =
393  {
394  // These values are the same as Legendre
395  {0,3,Point(-0.0959817, -0.0959817, 0.0702635)},
396  {1,6,Point(0.0625228, -0.0625228, 0.0457699)},
397  {2,2,Point(0.358209, 0.0959817, -2.77556e-17)},
398  {3,7,Point(-0.233338, 0.0625228, -6.93889e-17)},
399  // These values are different from Legendre
400  {4,12,Point(-0.0807678,-0.0807678,-0.182443)},
401  {5,5,Point(0.385213,0.103218,-0.175079)},
402  {6,0,Point(0.0840179,-0.0225125,0.721474)},
403  {7,10,Point(-0.0198117,0.0198117,0.0357373)},
404  {8,16,Point(0.00110804,0.00110804,0.000834194)},
405  {9,18,Point(0.00662875,0.00177617,0.000482244)},
406  {10,2,Point(-0.0840179,0.0840179,2.46303)},
407  {11,17,Point(-0.00177617,0.00177617,0.00142946)},
408  };
409 
410  tabulated_grads[std::make_pair(THIRD, JACOBI_30_00)] =
411  {
412  {12,9,Point(0.0184475,0.0688471,-0.254748)},
413  {13,6,Point(-0.230424,0.230424,1.15264)},
414  {14,20,Point(-0.000965042,0.00360159,0.000183379)},
415  {15,15,Point(-0.0423204,0.0113397,0.12514)},
416  };
417 
418  tabulated_grads[std::make_pair(FOURTH, JACOBI_30_00)] =
419  {
420  {16,3,Point(-0.0476508,0.012768,0.771323)},
421  {17,10,Point(0.333487,-0.333487,1.01421)},
422  {18,14,Point(0,0,0)},
423  {19,27,Point(-0.00698362,0.00187126,0.000667664)},
424  };
425  }
426 
427  void tearDown() {}
428 
429 };
430 
libMesh::InfElemBuilder::build_inf_elem
const Point build_inf_elem(const bool be_verbose=false)
Build infinite elements atop a volume-based mesh, determine origin automatically.
Definition: inf_elem_builder.C:41
InfFERadialTest::TabulatedVal::qp
unsigned int qp
Definition: inf_fe_radial_test.C:45
InfFERadialTest::TabulatedVal::TabulatedVal
TabulatedVal(unsigned int i_in, unsigned int qp_in, Real val_in)
Definition: inf_fe_radial_test.C:42
InfFERadialTest
This class is for unit testing the "radial" basis function aspects of the InfFE.
Definition: inf_fe_radial_test.C:32
libMesh::QGauss
This class implements specific orders of Gauss quadrature.
Definition: quadrature_gauss.h:39
libMesh::HEX8
Definition: enum_elem_type.h:47
InfFERadialTest::testSingleOrder
void testSingleOrder(Order radial_order, FEFamily radial_family)
Definition: inf_fe_radial_test.C:82
libMesh::MeshBase::n_elem
virtual dof_id_type n_elem() const =0
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::JACOBI_20_00
Definition: enum_fe_family.h:49
CPPUNIT_TEST_SUITE_REGISTRATION
CPPUNIT_TEST_SUITE_REGISTRATION(InfFERadialTest)
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::Elem::infinite
virtual bool infinite() const =0
libMesh::MeshTools::Generation::build_cube
void build_cube(UnstructuredMesh &mesh, const unsigned int nx=0, const unsigned int ny=0, const unsigned int nz=0, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const Real zmin=0., const Real zmax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
Builds a (elements) cube.
Definition: mesh_generation.C:298
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
InfFERadialTest::tabulated_grads
std::map< KeyType, TabulatedGrads > tabulated_grads
Definition: inf_fe_radial_test.C:62
libMesh::SECOND
Definition: enum_order.h:43
InfFERadialTest::testDifferentOrders
void testDifferentOrders()
Definition: inf_fe_radial_test.C:64
InfFERadialTest::TabulatedVals
std::vector< TabulatedVal > TabulatedVals
Definition: inf_fe_radial_test.C:58
libMesh::MeshBase::elem_ptr
virtual const Elem * elem_ptr(const dof_id_type i) const =0
libMesh::FEGenericBase::build_InfFE
static std::unique_ptr< FEGenericBase > build_InfFE(const unsigned int dim, const FEType &type)
Builds a specific infinite element type.
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::CARTESIAN
Definition: enum_inf_map_type.h:35
libMesh::ReplicatedMesh
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
Definition: replicated_mesh.h:47
InfFERadialTest::KeyType
std::pair< Order, FEFamily > KeyType
Definition: inf_fe_radial_test.C:60
libMesh::JACOBI_30_00
Definition: enum_fe_family.h:50
InfFERadialTest::TabulatedGrads
std::vector< TabulatedGrad > TabulatedGrads
Definition: inf_fe_radial_test.C:59
libMesh::libmesh_ignore
void libmesh_ignore(const Args &...)
Definition: libmesh_common.h:526
TestCommWorld
libMesh::Parallel::Communicator * TestCommWorld
Definition: driver.C:111
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::InfElemBuilder
This class is used to build infinite elements on top of an existing mesh.
Definition: inf_elem_builder.h:53
InfFERadialTest::TabulatedVal::val
Real val
Definition: inf_fe_radial_test.C:46
libMesh::LEGENDRE
Definition: enum_fe_family.h:51
InfFERadialTest::TabulatedGrad::qp
unsigned int qp
Definition: inf_fe_radial_test.C:54
InfFERadialTest::TabulatedVal
Definition: inf_fe_radial_test.C:40
libMesh::FOURTH
Definition: enum_order.h:45
InfFERadialTest::TabulatedGrad::grad
Point grad
Definition: inf_fe_radial_test.C:55
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
InfFERadialTest::tabulated_vals
std::map< KeyType, TabulatedVals > tabulated_vals
Definition: inf_fe_radial_test.C:61
libMesh::FEFamily
FEFamily
Definition: enum_fe_family.h:34
libmesh_cppunit.h
InfFERadialTest::tearDown
void tearDown()
Definition: inf_fe_radial_test.C:427
InfFERadialTest::setUp
void setUp()
Definition: inf_fe_radial_test.C:189
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
InfFERadialTest::TabulatedGrad
Definition: inf_fe_radial_test.C:49
InfFERadialTest::TabulatedVal::i
unsigned int i
Definition: inf_fe_radial_test.C:44
libMesh::THIRD
Definition: enum_order.h:44
libMesh::LAGRANGE
Definition: enum_fe_family.h:36
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
test_comm.h
libMesh::FEType::default_quadrature_order
Order default_quadrature_order() const
Definition: fe_type.h:353
libMesh::FIRST
Definition: enum_order.h:42
InfFERadialTest::TabulatedGrad::i
unsigned int i
Definition: inf_fe_radial_test.C:53
InfFERadialTest::TabulatedGrad::TabulatedGrad
TabulatedGrad(unsigned int i_in, unsigned int qp_in, Point grad_in)
Definition: inf_fe_radial_test.C:51