libMesh
fe_test.h
Go to the documentation of this file.
1 #ifndef FE_TEST_H
2 #define FE_TEST_H
3 
4 #include "test_comm.h"
5 
6 #include <libmesh/dof_map.h>
7 #include <libmesh/elem.h>
8 #include <libmesh/equation_systems.h>
9 #include <libmesh/fe_base.h>
10 #include <libmesh/fe_interface.h>
11 #include <libmesh/function_base.h>
12 #include <libmesh/mesh.h>
13 #include <libmesh/mesh_generation.h>
14 #include <libmesh/mesh_modification.h>
15 #include <libmesh/numeric_vector.h>
16 #include <libmesh/system.h>
17 #include <libmesh/quadrature_gauss.h>
18 
19 #include <vector>
20 
21 #include "libmesh_cppunit.h"
22 
23 #define FETEST \
24  CPPUNIT_TEST( testFEInterface ); \
25  CPPUNIT_TEST( testU ); \
26  CPPUNIT_TEST( testPartitionOfUnity ); \
27  CPPUNIT_TEST( testGradU ); \
28  CPPUNIT_TEST( testGradUComp ); \
29  CPPUNIT_TEST( testHessU ); \
30  CPPUNIT_TEST( testHessUComp ); \
31  CPPUNIT_TEST( testDualDoesntScreamAndDie ); \
32  CPPUNIT_TEST( testCustomReinit );
33 
34 using namespace libMesh;
35 
36 
37 class SkewFunc : public FunctionBase<Real>
38 {
39  std::unique_ptr<FunctionBase<Real>> clone () const override
40  { return std::make_unique<SkewFunc>(); }
41 
42  Real operator() (const Point &,
43  const Real = 0.) override
44  { libmesh_not_implemented(); } // scalar-only API
45 
46  // Skew in x based on y, y based on z
47  void operator() (const Point & p,
48  const Real,
49  DenseVector<Real> & output)
50  {
51  output.resize(3);
52  output(0) = p(0);
53 #if LIBMESH_DIM > 0
54  output(0) += .1*p(1);
55  output(1) = p(1);
56 #endif
57 #if LIBMESH_DIM > 1
58  output(1) += .2*p(2);
59  output(2) = p(2);
60 #endif
61  }
62 };
63 
64 
65 inline
67  const Parameters&,
68  const std::string&,
69  const std::string&)
70 {
71  const Real & x = p(0);
72  const Real & y = (LIBMESH_DIM > 1) ? p(1) : 0;
73  const Real & z = (LIBMESH_DIM > 2) ? p(2) : 0;
74 
75  return x + 0.25*y + 0.0625*z;
76 }
77 
78 inline
80  const Parameters&,
81  const std::string&,
82  const std::string&)
83 {
84  Gradient grad = 1;
85  if (LIBMESH_DIM > 1)
86  grad(1) = 0.25;
87  if (LIBMESH_DIM > 2)
88  grad(2) = 0.0625;
89 
90  return grad;
91 }
92 
93 
94 inline
96  const Parameters&,
97  const std::string&,
98  const std::string&)
99 {
100  const Real & x = p(0);
101  const Real & y = (LIBMESH_DIM > 1) ? p(1) : 0;
102  const Real & z = (LIBMESH_DIM > 2) ? p(2) : 0;
103 
104  return x*x + 0.5*y*y + 0.25*z*z + 0.125*x*y + 0.0625*x*z + 0.03125*y*z;
105 }
106 
107 inline
109  const Parameters&,
110  const std::string&,
111  const std::string&)
112 {
113  const Real & x = p(0);
114  const Real & y = (LIBMESH_DIM > 1) ? p(1) : 0;
115  const Real & z = (LIBMESH_DIM > 2) ? p(2) : 0;
116 
117  Gradient grad = 2*x + 0.125*y + 0.0625*z;
118  if (LIBMESH_DIM > 1)
119  grad(1) = y + 0.125*x + 0.03125*z;
120  if (LIBMESH_DIM > 2)
121  grad(2) = 0.5*z + 0.0625*x + 0.03125*y;
122 
123  return grad;
124 }
125 
126 
127 inline
129  const Parameters&,
130  const std::string&,
131  const std::string&)
132 {
133  const Real & x = p(0);
134  const Real & y = (LIBMESH_DIM > 1) ? p(1) : 0;
135  const Real & z = (LIBMESH_DIM > 2) ? p(2) : 0;
136 
137  return x*(1-x)*(1-x) + x*x*(1-y) + x*(1-y)*(1-z) + y*(1-y)*z + z*(1-z)*(1-z);
138 }
139 
140 inline
142  const Parameters&,
143  const std::string&,
144  const std::string&)
145 {
146  const Real & x = p(0);
147  const Real & y = (LIBMESH_DIM > 1) ? p(1) : 0;
148  const Real & z = (LIBMESH_DIM > 2) ? p(2) : 0;
149 
150  Gradient grad = 3*x*x-4*x+1 + 2*x*(1-y) + (1-y)*(1-z);
151  if (LIBMESH_DIM > 1)
152  grad(1) = -x*x - x*(1-z) + (1-2*y)*z;
153  if (LIBMESH_DIM > 2)
154  grad(2) = -x*(1-y) + y*(1-y) + 3*z*z-4*z+1;
155 
156  return grad;
157 }
158 
159 
160 inline
162  const Parameters&,
163  const std::string&,
164  const std::string&)
165 {
166  const Real & x = p(0);
167  const Real & y = (LIBMESH_DIM > 1) ? p(1) : 0;
168  const Real & z = (LIBMESH_DIM > 2) ? p(2) : 0;
169 
170  return x*x*(1-x)*(1-x) + x*x*z*(1-y) + x*(1-x)*(1-y)*(1-z) + (1-x)*y*(1-y)*z + z*z*(1-z)*(1-z);
171 }
172 
173 inline
175  const Parameters&,
176  const std::string&,
177  const std::string&)
178 {
179  const Real & x = p(0);
180  const Real & y = (LIBMESH_DIM > 1) ? p(1) : 0;
181  const Real & z = (LIBMESH_DIM > 2) ? p(2) : 0;
182 
183  Gradient grad = 4*x*x*x-6*x*x+2*x + 2*x*z*(1-y) + (1-2*x)*(1-y)*(1-z) - y*(1-y)*z;
184  if (LIBMESH_DIM > 1)
185  grad(1) = -x*x*z - x*(1-x)*(1-z) + (1-x)*(1-2*y)*z;
186  if (LIBMESH_DIM > 2)
187  grad(2) = x*x*(1-y) - x*(1-x)*(1-y) + (1-x)*y*(1-y) + 4*z*z*z-6*z*z+2*z;
188 
189  return grad;
190 }
191 
192 
193 
194 // Higher order rational bases need uniform weights to exactly
195 // represent linears; we can easily try out other functions on
196 // tensor product elements.
197 static const Real rational_w = 0.75;
198 
199 inline
201  const Parameters&,
202  const std::string&,
203  const std::string&)
204 {
205  const Real & x = p(0);
206  const Real & y = (LIBMESH_DIM > 1) ? p(1) : 0;
207  const Real & z = (LIBMESH_DIM > 2) ? p(2) : 0;
208 
209  const Real denom = ((1-x)*(1-x)+x*x+2*rational_w*x*(1-x))*
210  ((1-y)*(1-y)+y*y+2*rational_w*y*(1-y))*
211  ((1-z)*(1-z)+z*z+2*rational_w*z*(1-z));
212 
213  return (x + 0.25*y + 0.0625*z)/denom;
214 }
215 
216 inline
218  const Parameters&,
219  const std::string&,
220  const std::string&)
221 {
222  const Real & x = p(0);
223  const Real & y = (LIBMESH_DIM > 1) ? p(1) : 0;
224  const Real & z = (LIBMESH_DIM > 2) ? p(2) : 0;
225 
226  const Real xpoly = (1-x)*(1-x)+x*x+2*rational_w*x*(1-x);
227  const Real xderiv = -2*(1-x)+2*x+2*rational_w*(1-2*x);
228  const Real ypoly = (1-y)*(1-y)+y*y+2*rational_w*y*(1-y);
229  const Real yderiv = -2*(1-y)+2*y+2*rational_w*(1-2*y);
230  const Real zpoly = (1-z)*(1-z)+z*z+2*rational_w*z*(1-z);
231  const Real zderiv = -2*(1-z)+2*z+2*rational_w*(1-2*z);
232 
233  const Real denom = xpoly * ypoly * zpoly;
234 
235  const Real numer = (x + 0.25*y + 0.0625*z);
236 
237  Gradient grad_n = 1, grad_d = xderiv * ypoly * zpoly;
238  if (LIBMESH_DIM > 1)
239  {
240  grad_n(1) = 0.25;
241  grad_d(1) = xpoly * yderiv * zpoly;
242  }
243  if (LIBMESH_DIM > 2)
244  {
245  grad_n(2) = 0.0625;
246  grad_d(2) = xpoly * ypoly * zderiv;
247  }
248 
249  Gradient grad = (grad_n - numer * grad_d / denom) / denom;
250 
251  return grad;
252 }
253 
254 
255 #define FE_CAN_TEST_CUBIC \
256  (((family != LAGRANGE && family != L2_LAGRANGE) || \
257  (elem_type != TRI7 && elem_type != TET14 && \
258  elem_type != PRISM20 && elem_type != PRISM21 && \
259  elem_type != PYRAMID18)) && order > 2)
260 
261 
262 template <Order order, FEFamily family, ElemType elem_type, unsigned int build_nx>
263 class FETestBase : public CppUnit::TestCase {
264 
265 protected:
266  std::string libmesh_suite_name;
267 
268  unsigned int _dim, _nx, _ny, _nz;
270  std::vector<dof_id_type> _dof_indices;
272  std::unique_ptr<Mesh> _mesh;
273  std::unique_ptr<EquationSystems> _es;
274  std::unique_ptr<FEBase> _fe;
275  std::unique_ptr<QGauss> _qrule;
276 
277  Real _value_tol, _grad_tol, _hess_tol;
278 
279 
281  {
282  Parameters dummy;
283 
284  Gradient true_grad;
285  RealGradient returnval;
286 
287  if (family == RATIONAL_BERNSTEIN && order > 1)
288  true_grad = rational_test_grad(p, dummy, "", "");
289  else if (order > 3)
290  true_grad = fe_quartic_test_grad(p, dummy, "", "");
291  else if (FE_CAN_TEST_CUBIC)
292  true_grad = fe_cubic_test_grad(p, dummy, "", "");
293  else if (order > 1)
294  {
295  const Real & x = p(0);
296  const Real & y = (LIBMESH_DIM > 1) ? p(1) : 0;
297  const Real & z = (LIBMESH_DIM > 2) ? p(2) : 0;
298 
299  true_grad = Gradient(2*x+0.125*y+0.0625*z,
300  y+0.125*x+0.03125*z,
301  0.5*z+0.0625*x+0.03125*y);
302  }
303  else
304  true_grad = Gradient(1.0, 0.25, 0.0625);
305 
306  for (unsigned int d=0; d != LIBMESH_DIM; ++d)
307  {
308  CPPUNIT_ASSERT(true_grad(d) ==
309  Number(libmesh_real(true_grad(d))));
310 
311  returnval(d) = libmesh_real(true_grad(d));
312  }
313 
314  return returnval;
315  }
316 
317 
319  {
320  const Real & x = p(0);
321  const Real & y = LIBMESH_DIM > 1 ? p(1) : 0;
322  const Real & z = LIBMESH_DIM > 2 ? p(2) : 0;
323 
324  if (order > 3)
325  return RealTensor
326  { 12*x*x-12*x+2+2*z*(1-y)-2*(1-y)*(1-z), -2*x*z-(1-2*x)*(1-z)-(1-2*y)*z, 2*x*(1-y)-(1-2*x)*(1-y)-y*(1-y),
327  -2*x*z-(1-2*x)*(1-z)-(1-2*y)*z, -2*(1-x)*z, -x*x+x*(1-x)+(1-x)*(1-2*y),
328  2*x*(1-y)-(1-2*x)*(1-y)-y*(1-y), -x*x+x*(1-x)+(1-x)*(1-2*y), 12*z*z-12*z+2 };
329  else if (FE_CAN_TEST_CUBIC)
330  return RealTensor
331  { 6*x-4+2*(1-y), -2*x+z-1, y-1,
332  -2*x+z-1, -2*z, x+1-2*y,
333  y-1, x+1-2*y, 6*z-4 };
334  else if (order > 1)
335  return RealTensor
336  { 2, 0.125, 0.0625,
337  0.125, 1, 0.03125,
338  0.0625, 0.03125, 0.5 };
339 
340  return RealTensor
341  { 0, 0, 0,
342  0, 0, 0,
343  0, 0, 0 };
344  }
345 
346 public:
347  void setUp()
348  {
349  _mesh = std::make_unique<Mesh>(*TestCommWorld);
350  const std::unique_ptr<Elem> test_elem = Elem::build(elem_type);
351  _dim = test_elem->dim();
352  const unsigned int build_ny = (_dim > 1) * build_nx;
353  const unsigned int build_nz = (_dim > 2) * build_nx;
354 
355  unsigned char weight_index = 0;
356 
357  if (family == RATIONAL_BERNSTEIN)
358  {
359  // Make sure we can handle non-zero weight indices
360  _mesh->add_node_integer("buffer integer");
361 
362  // Default weight to 1.0 so we don't get NaN from contains_point
363  // checks with a default GhostPointNeighbors ghosting
364  const Real default_weight = 1.0;
365  weight_index = cast_int<unsigned char>
366  (_mesh->add_node_datum<Real>("rational_weight", true,
367  &default_weight));
368  libmesh_assert_not_equal_to(weight_index, 0);
369 
370  // Set mapping data but not type, since here we're testing
371  // rational basis functions in the FE space but testing then
372  // with Lagrange bases for the mapping space.
373  _mesh->set_default_mapping_data(weight_index);
374  }
375 
377  build_nx, build_ny, build_nz,
378  0., 1., 0., 1., 0., 1.,
379  elem_type);
380 
381  // For debugging purposes it can be helpful to only consider one
382  // element even when we're using an element type that requires
383  // more than one element to fill out a square or cube.
384 #if 0
385  for (dof_id_type i = 0; i != _mesh->max_elem_id(); ++i)
386  {
387  Elem * elem = _mesh->query_elem_ptr(i);
388  if (elem && elem->id())
389  _mesh->delete_elem(elem);
390  }
391  _mesh->prepare_for_use();
392  CPPUNIT_ASSERT_EQUAL(_mesh->n_elem(), dof_id_type(1));
393 #endif
394 
395  // Permute our elements randomly and rotate and skew our mesh so
396  // we test all sorts of orientations ... except with Hermite
397  // elements, which are only designed to support meshes with a
398  // single orientation shared by all elements. We're also not
399  // rotating and/or skewing the rational elements, since our test
400  // solution was designed for a specific weighted mesh.
401  if (family != HERMITE &&
402  family != RATIONAL_BERNSTEIN)
403  {
405 
406  // Not yet testing manifolds embedded in higher-D space
407  if (_dim > 1)
409  8*(_dim>2), 16*(_dim>2));
410 
411  SkewFunc skew_func;
412  MeshTools::Modification::redistribute(*_mesh, skew_func);
413  }
414 
415  // Set rational weights so we can exactly match our test solution
416  if (family == RATIONAL_BERNSTEIN)
417  {
418  for (auto elem : _mesh->active_element_ptr_range())
419  {
420  const unsigned int nv = elem->n_vertices();
421  const unsigned int nn = elem->n_nodes();
422  // We want interiors in lower dimensional elements treated
423  // like edges or faces as appropriate.
424  const unsigned int n_edges =
425  (elem->type() == EDGE3) ? 1 : elem->n_edges();
426  const unsigned int n_faces =
427  (elem->type() == QUAD9) ? 1 : elem->n_faces();
428  const unsigned int nve = std::min(nv + n_edges, nn);
429  const unsigned int nvef = std::min(nve + n_faces, nn);
430 
431  for (unsigned int i = 0; i != nv; ++i)
432  elem->node_ref(i).set_extra_datum<Real>(weight_index, 1.);
433  for (unsigned int i = nv; i != nve; ++i)
434  elem->node_ref(i).set_extra_datum<Real>(weight_index, rational_w);
435  const Real w2 = rational_w * rational_w;
436  for (unsigned int i = nve; i != nvef; ++i)
437  elem->node_ref(i).set_extra_datum<Real>(weight_index, w2);
438  const Real w3 = rational_w * w2;
439  for (unsigned int i = nvef; i != nn; ++i)
440  elem->node_ref(i).set_extra_datum<Real>(weight_index, w3);
441  }
442  }
443 
444  _es = std::make_unique<EquationSystems>(*_mesh);
445  _sys = &(_es->add_system<System> ("SimpleSystem"));
446  _sys->add_variable("u", order, family);
447  _es->init();
448 
449  if (family == RATIONAL_BERNSTEIN && order > 1)
450  {
451  _sys->project_solution(rational_test, rational_test_grad, _es->parameters);
452  }
453  else if (order > 3)
454  {
455  _sys->project_solution(fe_quartic_test, fe_quartic_test_grad, _es->parameters);
456  }
457  // Lagrange "cubic" on Tet7 only supports a bubble function, not
458  // all of P^3
459  else if (FE_CAN_TEST_CUBIC)
460  {
461  _sys->project_solution(fe_cubic_test, fe_cubic_test_grad, _es->parameters);
462  }
463  else if (order > 1)
464  {
465  _sys->project_solution(quadratic_test, quadratic_test_grad, _es->parameters);
466  }
467  else
468  {
469  _sys->project_solution(linear_test, linear_test_grad, _es->parameters);
470  }
471 
472  FEType fe_type = _sys->variable_type(0);
473  _fe = FEBase::build(_dim, fe_type);
474 
475  // Create quadrature rule for use in computing dual shape coefficients
476  _qrule = std::make_unique<QGauss>(_dim, fe_type.default_quadrature_order());
477  _fe->attach_quadrature_rule(_qrule.get());
478 
479  auto rng = _mesh->active_local_element_ptr_range();
480  this->_elem = rng.begin() == rng.end() ? nullptr : *(rng.begin());
481 
482  _sys->get_dof_map().dof_indices(this->_elem, _dof_indices);
483 
484  _nx = 10;
485  _ny = (_dim > 1) ? _nx : 1;
486  _nz = (_dim > 2) ? _nx : 1;
487 
488  this->_value_tol = TOLERANCE * sqrt(TOLERANCE);
489 
490  // We see 6.5*tol*sqrt(tol) errors on cubic Hermites with the fe_cubic
491  // hermite test function
492  // On Tri7 we see 10*tol*sqrt(tol) errors, even!
493  this->_grad_tol = 12 * TOLERANCE * sqrt(TOLERANCE);
494 
495  this->_hess_tol = sqrt(TOLERANCE); // FIXME: we see some ~1e-5 errors?!?
496 
497  // Prerequest everything we'll want to calculate later.
498  _fe->get_phi();
499  _fe->get_dphi();
500  _fe->get_dphidx();
501 #if LIBMESH_DIM > 1
502  _fe->get_dphidy();
503 #endif
504 #if LIBMESH_DIM > 2
505  _fe->get_dphidz();
506 #endif
507 
508 #if LIBMESH_ENABLE_SECOND_DERIVATIVES
509 
510  // Szabab elements don't have second derivatives yet
511  if (family == SZABAB)
512  return;
513 
514  _fe->get_d2phi();
515  _fe->get_d2phidx2();
516 #if LIBMESH_DIM > 1
517  _fe->get_d2phidxdy();
518  _fe->get_d2phidy2();
519 #endif
520 #if LIBMESH_DIM > 2
521  _fe->get_d2phidxdz();
522  _fe->get_d2phidydz();
523  _fe->get_d2phidz2();
524 #endif
525 
526 #endif
527  }
528 
529  void tearDown() {}
530 };
531 
532 
533 
534 template <Order order, FEFamily family, ElemType elem_type>
535 class FETest : public FETestBase<order, family, elem_type, 1> {
536 
537 public:
538 
539  template <typename Functor>
540  void testLoop(Functor f)
541  {
542  // Handle the "more processors than elements" case
543  if (!this->_elem)
544  return;
545 
546  // These tests require exceptions to be enabled because a
547  // TypeTensor::solve() call down in Elem::contains_point()
548  // actually throws a non-fatal exception for a certain Point which
549  // is not in the Elem. When exceptions are not enabled, this test
550  // simply aborts.
551 #ifdef LIBMESH_ENABLE_EXCEPTIONS
552  for (unsigned int i=0; i != this->_nx; ++i)
553  for (unsigned int j=0; j != this->_ny; ++j)
554  for (unsigned int k=0; k != this->_nz; ++k)
555  {
556  Point p = Real(i)/this->_nx;
557  if (j > 0)
558  p(1) = Real(j)/this->_ny;
559  if (k > 0)
560  p(2) = Real(k)/this->_nz;
561  if (!this->_elem->contains_point(p))
562  continue;
563 
564  // If at a singular node, cannot use FEMap::map
565  if (this->_elem->local_singular_node(p) != invalid_uint)
566  continue;
567 
568  std::vector<Point> master_points
569  (1, FEMap::inverse_map(this->_dim, this->_elem, p));
570 
571  // Reinit at point to test against analytic solution
572  this->_fe->reinit(this->_elem, &master_points);
573 
574  f(p);
575  }
576 #endif // LIBMESH_ENABLE_EXCEPTIONS
577  }
578 
580  {
581  if (!this->_elem)
582  return;
583 
584  this->_fe->reinit(this->_elem);
585 
586  bool satisfies_partition_of_unity = true;
587  for (const auto qp : make_range(this->_qrule->n_points()))
588  {
589  Real phi_sum = 0;
590  for (std::size_t d = 0; d != this->_dof_indices.size(); ++d)
591  phi_sum += this->_fe->get_phi()[d][qp];
592  if (phi_sum < (1 - TOLERANCE) || phi_sum > (1 + TOLERANCE))
593  {
594  satisfies_partition_of_unity = false;
595  break;
596  }
597  }
598 
599  switch (this->_fe->get_family())
600  {
601  case MONOMIAL:
602  {
603  switch (this->_fe->get_order())
604  {
605  case CONSTANT:
606  CPPUNIT_ASSERT(satisfies_partition_of_unity);
607  break;
608 
609  default:
610  CPPUNIT_ASSERT(!satisfies_partition_of_unity);
611  break;
612  }
613  break;
614  }
615  case SZABAB:
616  case HIERARCHIC:
617  case L2_HIERARCHIC:
618  {
619  switch (this->_fe->get_order())
620  {
621  case FIRST:
622  CPPUNIT_ASSERT(satisfies_partition_of_unity);
623  break;
624 
625  default:
626  CPPUNIT_ASSERT(!satisfies_partition_of_unity);
627  break;
628  }
629  break;
630  }
631 
632  case XYZ:
633  case CLOUGH:
634  case HERMITE:
635  {
636  CPPUNIT_ASSERT(!satisfies_partition_of_unity);
637  break;
638  }
639 
640  case LAGRANGE:
641  case L2_LAGRANGE:
642  case BERNSTEIN:
643  case RATIONAL_BERNSTEIN:
644  {
645  CPPUNIT_ASSERT(satisfies_partition_of_unity);
646  break;
647  }
648 
649  default:
650  CPPUNIT_FAIL("Uncovered FEFamily");
651  }
652  }
653 
654 
656  {
657  LOG_UNIT_TEST;
658 
659  // Handle the "more processors than elements" case
660  if (!this->_elem)
661  return;
662 
663  this->_fe->reinit(this->_elem);
664 
665  const FEType fe_type = this->_sys->variable_type(0);
666 
667  CPPUNIT_ASSERT_EQUAL(
668  FEInterface::n_shape_functions(fe_type, this->_elem),
669  this->_fe->n_shape_functions());
670 
671  CPPUNIT_ASSERT_EQUAL(
673  this->_fe->get_continuity());
674 
675  CPPUNIT_ASSERT_EQUAL(
677  this->_fe->is_hierarchic());
678  }
679 
680  void testU()
681  {
682  LOG_UNIT_TEST;
683 
684  auto f = [this](Point p)
685  {
686  Parameters dummy;
687 
688  Number u = 0;
689  for (std::size_t d = 0; d != this->_dof_indices.size(); ++d)
690  u += this->_fe->get_phi()[d][0] * (*this->_sys->current_local_solution)(this->_dof_indices[d]);
691 
692  Number true_u;
693 
694  if (family == RATIONAL_BERNSTEIN && order > 1)
695  true_u = rational_test(p, dummy, "", "");
696  else if (order > 3)
697  true_u = fe_quartic_test(p, dummy, "", "");
698  else if (FE_CAN_TEST_CUBIC)
699  true_u = fe_cubic_test(p, dummy, "", "");
700  else if (order > 1)
701  true_u = p(0)*p(0) + 0.5*p(1)*p(1) + 0.25*p(2)*p(2) +
702  0.125*p(0)*p(1) + 0.0625*p(0)*p(2) + 0.03125*p(1)*p(2);
703  else
704  true_u = p(0) + 0.25*p(1) + 0.0625*p(2);
705 
706  LIBMESH_ASSERT_FP_EQUAL
707  (libmesh_real(true_u), libmesh_real(u), this->_value_tol);
708  };
709 
710  testLoop(f);
711  }
712 
714  {
715  LOG_UNIT_TEST;
716 
717  // Handle the "more processors than elements" case
718  if (!this->_elem)
719  return;
720 
721  // Request dual calculations
722  this->_fe->get_dual_phi();
723 
724  // reinit using the default quadrature rule in order to calculate the dual coefficients
725  this->_fe->reinit(this->_elem);
726  }
727 
728 
729  void testGradU()
730  {
731  LOG_UNIT_TEST;
732 
733  auto f = [this](Point p)
734  {
735  Parameters dummy;
736 
737  Gradient grad_u = 0;
738  for (std::size_t d = 0; d != this->_dof_indices.size(); ++d)
739  grad_u += this->_fe->get_dphi()[d][0] * (*this->_sys->current_local_solution)(this->_dof_indices[d]);
740 
741  RealGradient true_grad = this->true_gradient(p);
742 
743  LIBMESH_ASSERT_FP_EQUAL(libmesh_real(grad_u(0)),
744  true_grad(0), this->_grad_tol);
745  if (this->_dim > 1)
746  LIBMESH_ASSERT_FP_EQUAL(libmesh_real(grad_u(1)),
747  true_grad(1), this->_grad_tol);
748  if (this->_dim > 2)
749  LIBMESH_ASSERT_FP_EQUAL(libmesh_real(grad_u(2)),
750  true_grad(2), this->_grad_tol);
751  };
752 
753  testLoop(f);
754  }
755 
757  {
758  LOG_UNIT_TEST;
759 
760  auto f = [this](Point p)
761  {
762  Parameters dummy;
763 
764  Number grad_u_x = 0, grad_u_y = 0, grad_u_z = 0;
765  for (std::size_t d = 0; d != this->_dof_indices.size(); ++d)
766  {
767  grad_u_x += this->_fe->get_dphidx()[d][0] * (*this->_sys->current_local_solution)(this->_dof_indices[d]);
768 #if LIBMESH_DIM > 1
769  grad_u_y += this->_fe->get_dphidy()[d][0] * (*this->_sys->current_local_solution)(this->_dof_indices[d]);
770 #endif
771 #if LIBMESH_DIM > 2
772  grad_u_z += this->_fe->get_dphidz()[d][0] * (*this->_sys->current_local_solution)(this->_dof_indices[d]);
773 #endif
774  }
775 
776  RealGradient true_grad = this->true_gradient(p);
777 
778  LIBMESH_ASSERT_FP_EQUAL(libmesh_real(grad_u_x),
779  true_grad(0), this->_grad_tol);
780  if (this->_dim > 1)
781  LIBMESH_ASSERT_FP_EQUAL(libmesh_real(grad_u_y),
782  true_grad(1), this->_grad_tol);
783  if (this->_dim > 2)
784  LIBMESH_ASSERT_FP_EQUAL(libmesh_real(grad_u_z),
785  true_grad(2), this->_grad_tol);
786  };
787 
788  testLoop(f);
789  }
790 
791 
792  void testHessU()
793  {
794  LOG_UNIT_TEST;
795 
796  // Szabab elements don't have second derivatives yet
797  if (family == SZABAB)
798  return;
799 
800 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
801  auto f = [this](Point p)
802  {
803  Tensor hess_u;
804  for (std::size_t d = 0; d != this->_dof_indices.size(); ++d)
805  hess_u += this->_fe->get_d2phi()[d][0] * (*this->_sys->current_local_solution)(this->_dof_indices[d]);
806 
807  // TODO: Yeah we'll test the ugly expressions later.
808  if (family == RATIONAL_BERNSTEIN && order > 1)
809  return;
810 
811  RealTensor true_hess = this->true_hessian(p);
812 
813  LIBMESH_ASSERT_FP_EQUAL(true_hess(0,0), libmesh_real(hess_u(0,0)),
814  this->_hess_tol);
815  if (this->_dim > 1)
816  {
817  LIBMESH_ASSERT_FP_EQUAL(libmesh_real(hess_u(0,1)), libmesh_real(hess_u(1,0)),
818  this->_hess_tol);
819  LIBMESH_ASSERT_FP_EQUAL(true_hess(0,1), libmesh_real(hess_u(0,1)),
820  this->_hess_tol);
821  LIBMESH_ASSERT_FP_EQUAL(true_hess(1,1), libmesh_real(hess_u(1,1)),
822  this->_hess_tol);
823  }
824  if (this->_dim > 2)
825  {
826  LIBMESH_ASSERT_FP_EQUAL(libmesh_real(hess_u(0,2)), libmesh_real(hess_u(2,0)),
827  this->_hess_tol);
828  LIBMESH_ASSERT_FP_EQUAL(libmesh_real(hess_u(1,2)), libmesh_real(hess_u(2,1)),
829  this->_hess_tol);
830  LIBMESH_ASSERT_FP_EQUAL(true_hess(0,2), libmesh_real(hess_u(0,2)),
831  this->_hess_tol);
832  LIBMESH_ASSERT_FP_EQUAL(true_hess(1,2), libmesh_real(hess_u(1,2)),
833  this->_hess_tol);
834  LIBMESH_ASSERT_FP_EQUAL(true_hess(2,2), libmesh_real(hess_u(2,2)),
835  this->_hess_tol);
836  }
837  };
838 
839  testLoop(f);
840 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
841  }
842 
844  {
845 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
846  LOG_UNIT_TEST;
847 
848  // Szabab elements don't have second derivatives yet
849  if (family == SZABAB)
850  return;
851 
852  auto f = [this](Point p)
853  {
854  Number hess_u_xx = 0, hess_u_xy = 0, hess_u_yy = 0,
855  hess_u_xz = 0, hess_u_yz = 0, hess_u_zz = 0;
856  for (std::size_t d = 0; d != this->_dof_indices.size(); ++d)
857  {
858  hess_u_xx += this->_fe->get_d2phidx2()[d][0] * (*this->_sys->current_local_solution)(this->_dof_indices[d]);
859 #if LIBMESH_DIM > 1
860  hess_u_xy += this->_fe->get_d2phidxdy()[d][0] * (*this->_sys->current_local_solution)(this->_dof_indices[d]);
861  hess_u_yy += this->_fe->get_d2phidy2()[d][0] * (*this->_sys->current_local_solution)(this->_dof_indices[d]);
862 #endif
863 #if LIBMESH_DIM > 2
864  hess_u_xz += this->_fe->get_d2phidxdz()[d][0] * (*this->_sys->current_local_solution)(this->_dof_indices[d]);
865  hess_u_yz += this->_fe->get_d2phidydz()[d][0] * (*this->_sys->current_local_solution)(this->_dof_indices[d]);
866  hess_u_zz += this->_fe->get_d2phidz2()[d][0] * (*this->_sys->current_local_solution)(this->_dof_indices[d]);
867 #endif
868  }
869 
870  // TODO: Yeah we'll test the ugly expressions later.
871  if (family == RATIONAL_BERNSTEIN && order > 1)
872  return;
873 
874  RealTensor true_hess = this->true_hessian(p);
875 
876  LIBMESH_ASSERT_FP_EQUAL(true_hess(0,0), libmesh_real(hess_u_xx),
877  this->_hess_tol);
878  if (this->_dim > 1)
879  {
880  LIBMESH_ASSERT_FP_EQUAL(true_hess(0,1), libmesh_real(hess_u_xy),
881  this->_hess_tol);
882  LIBMESH_ASSERT_FP_EQUAL(true_hess(1,1), libmesh_real(hess_u_yy),
883  this->_hess_tol);
884  }
885  if (this->_dim > 2)
886  {
887  LIBMESH_ASSERT_FP_EQUAL(true_hess(0,2), libmesh_real(hess_u_xz),
888  this->_hess_tol);
889  LIBMESH_ASSERT_FP_EQUAL(true_hess(1,2), libmesh_real(hess_u_yz),
890  this->_hess_tol);
891  LIBMESH_ASSERT_FP_EQUAL(true_hess(2,2), libmesh_real(hess_u_zz),
892  this->_hess_tol);
893  }
894  };
895 
896  testLoop(f);
897 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
898  }
899 
901  {
902  LOG_UNIT_TEST;
903 
904  std::vector<Point> q_points;
905  std::vector<Real> weights;
906  q_points.resize(3); weights.resize(3);
907  q_points[0](0) = 0.0; q_points[0](1) = 0.0; weights[0] = Real(1)/6;
908  q_points[1](0) = 1.0; q_points[1](1) = 0.0; weights[1] = weights[0];
909  q_points[2](0) = 0.0; q_points[2](1) = 1.0; weights[2] = weights[0];
910 
911  FEType fe_type = this->_sys->variable_type(0);
912  std::unique_ptr<FEBase> fe (FEBase::build(this->_dim, fe_type));
913  const int extraorder = 3;
914  std::unique_ptr<QBase> qrule (fe_type.default_quadrature_rule (this->_dim, extraorder));
915  fe->attach_quadrature_rule (qrule.get());
916 
917  const std::vector<Point> & q_pos = fe->get_xyz();
918 
919  for (const auto & elem : this->_mesh->active_local_element_ptr_range()) {
920  fe->reinit (elem, &q_points, &weights);
921  CPPUNIT_ASSERT_EQUAL(q_points.size(), std::size_t(3));
922  CPPUNIT_ASSERT_EQUAL(q_pos.size(), std::size_t(3)); // 6? bug?
923  }
924  }
925 
926 };
927 
928 
929 #define INSTANTIATE_FETEST(order, family, elemtype) \
930  class FETest_##order##_##family##_##elemtype : public FETest<order, family, elemtype> { \
931  public: \
932  FETest_##order##_##family##_##elemtype() : \
933  FETest<order,family,elemtype>() { \
934  if (unitlog->summarized_logs_enabled()) \
935  this->libmesh_suite_name = "FETest"; \
936  else \
937  this->libmesh_suite_name = "FETest_" #order "_" #family "_" #elemtype; \
938  } \
939  CPPUNIT_TEST_SUITE( FETest_##order##_##family##_##elemtype ); \
940  FETEST \
941  CPPUNIT_TEST_SUITE_END(); \
942  }; \
943  \
944  CPPUNIT_TEST_SUITE_REGISTRATION( FETest_##order##_##family##_##elemtype );
945 
946 #endif
T libmesh_real(T a)
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:182
std::unique_ptr< FunctionBase< Real > > clone() const override
Definition: fe_test.h:39
unsigned int _nz
Definition: fe_test.h:268
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
Definition: libmesh.h:286
void testDualDoesntScreamAndDie()
Definition: fe_test.h:713
This class provides the ability to map between arbitrary, user-defined strings and several data types...
Definition: parameters.h:67
void testLoop(Functor f)
Definition: fe_test.h:540
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
Definition: dof_map.C:1992
static constexpr Real TOLERANCE
std::unique_ptr< QGauss > _qrule
Definition: fe_test.h:275
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true, const bool extra_checks=true)
Definition: fe_map.C:1626
Gradient quadratic_test_grad(const Point &p, const Parameters &, const std::string &, const std::string &)
Definition: fe_test.h:108
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:374
Elem * _elem
Definition: fe_test.h:269
void testGradUComp()
Definition: fe_test.h:756
void testFEInterface()
Definition: fe_test.h:655
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
Gradient rational_test_grad(const Point &p, const Parameters &, const std::string &, const std::string &)
Definition: fe_test.h:217
Order default_quadrature_order() const
Definition: fe_type.h:357
System * _sys
Definition: fe_test.h:271
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
Gradient fe_cubic_test_grad(const Point &p, const Parameters &, const std::string &, const std::string &)
Definition: fe_test.h:141
ADRealEigenVector< T, D, asd > sqrt(const ADRealEigenVector< T, D, asd > &)
Definition: type_vector.h:53
The libMesh namespace provides an interface to certain functionality in the library.
Gradient fe_quartic_test_grad(const Point &p, const Parameters &, const std::string &, const std::string &)
Definition: fe_test.h:174
Number fe_cubic_test(const Point &p, const Parameters &, const std::string &, const std::string &)
Definition: fe_test.h:128
static const Real rational_w
Definition: fe_test.h:197
void permute_elements(MeshBase &mesh)
Randomly permute the nodal ordering of each element (without twisting the element mapping)...
std::unique_ptr< EquationSystems > _es
Definition: fe_test.h:273
void tearDown()
Definition: fe_test.h:529
void project_solution(FunctionBase< Number > *f, FunctionBase< Gradient > *g=nullptr) const
Projects arbitrary functions onto the current solution.
std::unique_ptr< QBase > default_quadrature_rule(const unsigned int dim, const int extraorder=0) const
Definition: fe_type.C:34
Number linear_test(const Point &p, const Parameters &, const std::string &, const std::string &)
Definition: fe_test.h:66
void testHessUComp()
Definition: fe_test.h:843
Number rational_test(const Point &p, const Parameters &, const std::string &, const std::string &)
Definition: fe_test.h:200
dof_id_type id() const
Definition: dof_object.h:823
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:515
static RealTensor true_hessian(Point p)
Definition: fe_test.h:318
void testCustomReinit()
Definition: fe_test.h:900
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:273
NumberVectorValue Gradient
static bool is_hierarchic(const FEType &fe_type)
Returns whether or not the input FEType&#39;s higher-order shape functions are always hierarchic...
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
unsigned int add_variable(std::string_view var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
Definition: system.C:1305
void testHessU()
Definition: fe_test.h:792
Number quadratic_test(const Point &p, const Parameters &, const std::string &, const std::string &)
Definition: fe_test.h:95
std::string libmesh_suite_name
Definition: fe_test.h:266
Real _value_tol
Definition: fe_test.h:277
RealTensorValue rotate(MeshBase &mesh, const Real phi, const Real theta=0., const Real psi=0.)
Rotates the mesh in the xy plane.
std::unique_ptr< Mesh > _mesh
Definition: fe_test.h:272
static FEContinuity get_continuity(const FEType &fe_type)
Returns the input FEType&#39;s FEContinuity based on the underlying FEFamily and potentially the Order...
const FEType & variable_type(const unsigned int i) const
Definition: system.h:2427
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::vector< dof_id_type > _dof_indices
Definition: fe_test.h:270
static RealGradient true_gradient(Point p)
Definition: fe_test.h:280
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:134
void redistribute(MeshBase &mesh, const FunctionBase< Real > &mapfunc)
Deterministically perturb the nodal locations.
std::unique_ptr< FEBase > _fe
Definition: fe_test.h:274
Defines a dense vector for use in Finite Element-type computations.
Number fe_quartic_test(const Point &p, const Parameters &, const std::string &, const std::string &)
Definition: fe_test.h:161
Gradient linear_test_grad(const Point &, const Parameters &, const std::string &, const std::string &)
Definition: fe_test.h:79
Base class for functors that can be evaluated at a point and (optionally) time.
const DofMap & get_dof_map() const
Definition: system.h:2293
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
void testPartitionOfUnity()
Definition: fe_test.h:579
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.
void setUp()
Definition: fe_test.h:347
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
uint8_t dof_id_type
Definition: id_types.h:67
void testGradU()
Definition: fe_test.h:729
void testU()
Definition: fe_test.h:680