libMesh
fem_context.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2024 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 
20 #include "libmesh/fem_context.h"
21 
22 #include "libmesh/boundary_info.h"
23 #include "libmesh/diff_system.h"
24 #include "libmesh/dof_map.h"
25 #include "libmesh/elem.h"
26 #include "libmesh/fe_base.h"
27 #include "libmesh/fe_interface.h"
28 #include "libmesh/libmesh_logging.h"
29 #include "libmesh/mesh_base.h"
30 #include "libmesh/numeric_vector.h"
31 #include "libmesh/quadrature.h"
32 #include "libmesh/system.h"
33 #include "libmesh/time_solver.h"
34 #include "libmesh/unsteady_solver.h" // For euler_residual
35 
36 namespace libMesh
37 {
38 
40  const std::vector<unsigned int> * active_vars,
41  bool allocate_local_matrices)
42  : FEMContext(sys, sys.extra_quadrature_order, active_vars,
43  allocate_local_matrices)
44 {
45  init_internal_data(sys);
46 }
47 
49  int extra_quadrature_order,
50  const std::vector<unsigned int> * active_vars,
51  bool allocate_local_matrices)
52  : DiffContext(sys, allocate_local_matrices),
53  side(0), edge(0),
54  _mesh_sys(nullptr),
55  _mesh_x_var(0),
56  _mesh_y_var(0),
57  _mesh_z_var(0),
58  _atype(CURRENT),
59  _custom_solution(nullptr),
60  _boundary_info(sys.get_mesh().get_boundary_info()),
61  _elem(nullptr),
62  _dim(cast_int<unsigned char>(sys.get_mesh().mesh_dimension())),
63  _elem_dim(0), /* This will be reset in set_elem(). */
64  _elem_dims(sys.get_mesh().elem_dimensions()),
65  _element_qrule(4),
66  _side_qrule(4),
67  _extra_quadrature_order(extra_quadrature_order)
68 {
69  if (active_vars)
70  {
71  libmesh_assert(!active_vars->empty());
72  auto vars_copy =
73  std::make_unique<std::vector<unsigned int>>(*active_vars);
74 
75  // We want to do quick binary_search later
76  std::sort(vars_copy->begin(), vars_copy->end());
77 
78  _active_vars = std::move(vars_copy);
79  }
80 
81  init_internal_data(sys);
82 }
83 
84 
85 
87 {
88  const System & sys = this->get_system();
89  FEType hardest_fe_type =
91  (*_active_vars)[0] : 0);
92 
93  auto check_var = [&hardest_fe_type, &sys](unsigned int v)
94  {
95  FEType fe_type = sys.variable_type(v);
96 
97  // Make sure we find a non-SCALAR FE family, even in the case
98  // where the first variable(s) weren't
99  if (hardest_fe_type.family == SCALAR)
100  {
101  hardest_fe_type.family = fe_type.family;
102  hardest_fe_type.order = fe_type.order;
103  }
104 
105  // FIXME - we don't yet handle mixed finite elements from
106  // different families which require different quadrature rules
107  // libmesh_assert_equal_to (fe_type.family, hardest_fe_type.family);
108 
109  // We need to detect SCALAR's so we can prepare FE objects for
110  // them, and so we don't mistake high order scalars as a reason
111  // to crank up the quadrature order on other types.
112  if (fe_type.family != SCALAR && fe_type.order > hardest_fe_type.order)
113  hardest_fe_type = fe_type;
114  };
115 
116  if (_active_vars)
117  for (auto v : *_active_vars)
118  check_var(v);
119  else
120  for (auto v : make_range(sys.n_vars()))
121  check_var(v);
122 
123  return hardest_fe_type;
124 }
125 
126 
128 {
129  const System & sys = this->get_system();
130 
131  auto attach_rules = [this, &sys](unsigned int v)
132  {
133  for (const auto & dim : _elem_dims)
134  {
135  FEType fe_type = sys.variable_type(v);
136 
137  _element_fe[dim][fe_type]->attach_quadrature_rule(_element_qrule[dim].get());
138  if (dim)
139  _side_fe[dim][fe_type]->attach_quadrature_rule(_side_qrule[dim].get());
140  if (dim == 3)
141  _edge_fe[fe_type]->attach_quadrature_rule(_edge_qrule.get());
142  };
143  };
144 
145  if (_active_vars)
146  for (auto v : *_active_vars)
147  attach_rules(v);
148  else
149  for (auto v : make_range(sys.n_vars()))
150  attach_rules(v);
151 }
152 
153 
154 
155 void FEMContext::use_default_quadrature_rules(int extra_quadrature_order)
156 {
157  _extra_quadrature_order = extra_quadrature_order;
158 
159  FEType hardest_fe_type = this->find_hardest_fe_type();
160 
161  for (const auto & dim : _elem_dims)
162  {
163  // Create an adequate quadrature rule
166  if (dim)
167  _side_qrule[dim] =
169  if (dim == 3)
171  }
172 
173  this->attach_quadrature_rules();
174 }
175 
176 
177 void FEMContext::use_unweighted_quadrature_rules(int extra_quadrature_order)
178 {
179  _extra_quadrature_order = extra_quadrature_order;
180 
181  FEType hardest_fe_type = this->find_hardest_fe_type();
182 
183  for (const auto & dim : _elem_dims)
184  {
185  // Create an adequate quadrature rule
188  _side_qrule[dim] =
190  if (dim == 3)
192  }
193 
194  this->attach_quadrature_rules();
195 }
196 
197 
199 {
200  // Reserve space for the FEAbstract and QBase objects for each
201  // element dimension possibility (0,1,2,3)
202 
203  // Note: we would simply resize() all four of these vectors, but
204  // some compilers (ICC 19, MSVC) generate a diagnostic about copying
205  // a std::unique_ptr in this case, so the following two lines are a
206  // workaround.
207  _element_fe = std::vector<std::map<FEType, std::unique_ptr<FEAbstract>>>(4);
208  _side_fe = std::vector<std::map<FEType, std::unique_ptr<FEAbstract>>>(4);
209  _element_fe_var.resize(4);
210  _side_fe_var.resize(4);
211 
212  // We need to know which of our variables has the hardest
213  // shape functions to numerically integrate.
214 
215  unsigned int nv = sys.n_vars();
216  libmesh_assert (nv);
217 
218  bool have_scalar = false;
219 
220  if (_active_vars)
221  {
222  for (auto v : *_active_vars)
223  if (sys.variable_type(v).family == SCALAR)
224  {
225  have_scalar = true;
226  break;
227  }
228  }
229  else
230  {
231  for (auto v : make_range(sys.n_vars()))
232  if (sys.variable_type(v).family == SCALAR)
233  {
234  have_scalar = true;
235  break;
236  }
237  }
238 
239  if (have_scalar)
240  // SCALAR FEs have dimension 0 by assumption
241  _elem_dims.insert(0);
242 
243  auto build_var_fe = [this, &sys](unsigned int dim,
244  unsigned int i)
245  {
246  FEType fe_type = sys.variable_type(i);
247  const auto & dof_map = sys.get_dof_map();
248  const bool add_p_level = dof_map.should_p_refine(dof_map.var_group_from_var_number(i));
249 
250  auto & element_fe = _element_fe[dim][fe_type];
251  auto & side_fe = _side_fe[dim][fe_type];
252  if (!element_fe)
253  {
254  element_fe = FEAbstract::build(dim, fe_type);
255  element_fe->add_p_level_in_reinit(add_p_level);
256  side_fe = FEAbstract::build(dim, fe_type);
257  side_fe->add_p_level_in_reinit(add_p_level);
258 
259  if (dim == 3)
260  {
261  auto & edge_fe = _edge_fe[fe_type];
262  edge_fe = FEAbstract::build(dim, fe_type);
263  edge_fe->add_p_level_in_reinit(add_p_level);
264  }
265  }
266 
267  _element_fe_var[dim][i] = element_fe.get();
268  _side_fe_var[dim][i] = side_fe.get();
269  if ((dim) == 3)
270  _edge_fe_var[i] = _edge_fe[fe_type].get();
271  };
272 
273  for (const auto & dim : _elem_dims)
274  {
275  // Create finite element objects
276  _element_fe_var[dim].resize(nv);
277  _side_fe_var[dim].resize(nv);
278  if (dim == 3)
279  _edge_fe_var.resize(nv);
280 
281  if (_active_vars)
282  for (auto v : *_active_vars)
283  build_var_fe(dim, v);
284  else
285  for (auto v : make_range(nv))
286  build_var_fe(dim, v);
287  }
288 
290 }
291 
293 {
294 }
295 
296 
297 
299 {
300  return _boundary_info.has_boundary_id(&(this->get_elem()), side, id);
301 }
302 
303 
304 
305 void FEMContext::side_boundary_ids(std::vector<boundary_id_type> & vec_to_fill) const
306 {
307  _boundary_info.boundary_ids(&(this->get_elem()), side, vec_to_fill);
308 }
309 
310 
311 
312 template<typename OutputType,
314  FEMContext::diff_subsolution_getter subsolution_getter>
315 void FEMContext::some_value(unsigned int var, unsigned int qp, OutputType & u) const
316 {
317  // Get local-to-global dof index lookup
318  const unsigned int n_dofs = cast_int<unsigned int>
319  (this->get_dof_indices(var).size());
320 
321  // Get current local coefficients
322  const DenseSubVector<Number> & coef = (this->*subsolution_getter)(var);
323 
324  // Get finite element object
325  typename FENeeded<OutputType>::value_base * fe = nullptr;
326  (this->*fe_getter)( var, fe, this->get_elem_dim() );
327 
328  // Get shape function values at quadrature point
329  const std::vector<std::vector
330  <typename FENeeded<OutputType>::value_shape>> & phi = fe->get_phi();
331 
332  // Accumulate solution value
333  u = 0.;
334 
335  for (unsigned int l=0; l != n_dofs; l++)
336  u += phi[l][qp] * coef(l);
337 }
338 
339 
340 
341 template<typename OutputType,
343  FEMContext::diff_subsolution_getter subsolution_getter>
344 void FEMContext::some_gradient(unsigned int var, unsigned int qp, OutputType & du) const
345 {
346  // Get local-to-global dof index lookup
347  const unsigned int n_dofs = cast_int<unsigned int>
348  (this->get_dof_indices(var).size());
349 
350  // Get current local coefficients
351  const DenseSubVector<Number> & coef = (this->*subsolution_getter)(var);
352 
353  // Get finite element object
354  typename FENeeded<OutputType>::grad_base * fe = nullptr;
355  (this->*fe_getter)( var, fe, this->get_elem_dim() );
356 
357  // Get shape function values at quadrature point
358  const std::vector<std::vector
360  & dphi = fe->get_dphi();
361 
362  // Accumulate solution derivatives
363  du = 0;
364 
365  for (unsigned int l=0; l != n_dofs; l++)
366  du.add_scaled(dphi[l][qp], coef(l));
367 
368  return;
369 }
370 
371 
372 
373 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
374 template<typename OutputType,
376  FEMContext::diff_subsolution_getter subsolution_getter>
377 void FEMContext::some_hessian(unsigned int var, unsigned int qp, OutputType & d2u) const
378 {
379  // Get local-to-global dof index lookup
380  const unsigned int n_dofs = cast_int<unsigned int>
381  (this->get_dof_indices(var).size());
382 
383  // Get current local coefficients
384  const DenseSubVector<Number> & coef = (this->*subsolution_getter)(var);
385 
386  // Get finite element object
387  typename FENeeded<OutputType>::hess_base * fe = nullptr;
388  (this->*fe_getter)( var, fe, this->get_elem_dim() );
389 
390  // Get shape function values at quadrature point
391  const std::vector<std::vector
393  & d2phi = fe->get_d2phi();
394 
395  // Accumulate solution second derivatives
396  d2u = 0.0;
397 
398  for (unsigned int l=0; l != n_dofs; l++)
399  d2u.add_scaled(d2phi[l][qp], coef(l));
400 
401  return;
402 }
403 #endif
404 
405 
406 
407 Number FEMContext::interior_value(unsigned int var, unsigned int qp) const
408 {
409  Number u;
410 
411  this->interior_value( var, qp, u );
412 
413  return u;
414 }
415 
416 template<typename OutputType>
417 void FEMContext::interior_value(unsigned int var, unsigned int qp,
418  OutputType & u) const
419 {
420  this->some_value<OutputType,
421  &FEMContext::get_element_fe<typename TensorTools::MakeReal<OutputType>::type>,
422  &DiffContext::get_elem_solution>(var, qp, u);
423 }
424 
425 
426 template<typename OutputType>
427 void FEMContext::interior_values (unsigned int var,
428  const NumericVector<Number> & _system_vector,
429  std::vector<OutputType> & u_vals) const
430 {
431  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
432 
433  // Get local-to-global dof index lookup
434  const unsigned int n_dofs = cast_int<unsigned int>
435  (this->get_dof_indices(var).size());
436 
437  // Get current local coefficients
438  const DenseSubVector<Number> & coef = get_localized_subvector(_system_vector, var);
439 
440  // Get the finite element object
441  FEGenericBase<OutputShape> * fe = nullptr;
442  this->get_element_fe<OutputShape>( var, fe, this->get_elem_dim() );
443 
444  // Get shape function values at quadrature point
445  const std::vector<std::vector<OutputShape>> & phi = fe->get_phi();
446 
447  // Loop over all the q_points on this element
448  for (auto qp : index_range(u_vals))
449  {
450  OutputType & u = u_vals[qp];
451 
452  // Compute the value at this q_point
453  u = 0.;
454 
455  for (unsigned int l=0; l != n_dofs; l++)
456  u += phi[l][qp] * coef(l);
457  }
458 
459  return;
460 }
461 
463  unsigned int qp) const
464 {
465  Gradient du;
466 
467  this->interior_gradient( var, qp, du );
468 
469  return du;
470 }
471 
472 
473 
474 template<typename OutputType>
475 void FEMContext::interior_gradient(unsigned int var,
476  unsigned int qp,
477  OutputType & du) const
478 {
479  this->some_gradient<OutputType,
482  <OutputType>::type>::type>,
483  &DiffContext::get_elem_solution>(var, qp, du);
484 }
485 
486 
487 
488 template<typename OutputType>
489 void FEMContext::interior_gradients(unsigned int var,
490  const NumericVector<Number> & _system_vector,
491  std::vector<OutputType> & du_vals) const
492 {
493  typedef typename TensorTools::MakeReal
495  OutputShape;
496 
497  // Get local-to-global dof index lookup
498  const unsigned int n_dofs = cast_int<unsigned int>
499  (this->get_dof_indices(var).size());
500 
501  // Get current local coefficients
502  const DenseSubVector<Number> & coef = get_localized_subvector(_system_vector, var);
503 
504  // Get finite element object
505  FEGenericBase<OutputShape> * fe = nullptr;
506  this->get_element_fe<OutputShape>( var, fe, this->get_elem_dim() );
507 
508  // Get shape function values at quadrature point
509  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient>> & dphi = fe->get_dphi();
510 
511  // Loop over all the q_points in this finite element
512  for (auto qp : index_range(du_vals))
513  {
514  OutputType & du = du_vals[qp];
515 
516  // Compute the gradient at this q_point
517  du = 0;
518 
519  for (unsigned int l=0; l != n_dofs; l++)
520  du.add_scaled(dphi[l][qp], coef(l));
521  }
522 
523  return;
524 }
525 
526 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
527 Tensor FEMContext::interior_hessian(unsigned int var, unsigned int qp) const
528 {
529  Tensor d2u;
530 
531  this->interior_hessian( var, qp, d2u );
532 
533  return d2u;
534 }
535 
536 template<typename OutputType>
537 void FEMContext::interior_hessian(unsigned int var, unsigned int qp,
538  OutputType & d2u) const
539 {
540  this->some_hessian<OutputType,
542  <typename TensorTools::MakeReal
545  <OutputType>::type>::type>::type>,
546  &DiffContext::get_elem_solution>(var, qp, d2u);
547 }
548 
549 
550 template<typename OutputType>
551 void FEMContext::interior_hessians(unsigned int var,
552  const NumericVector<Number> & _system_vector,
553  std::vector<OutputType> & d2u_vals) const
554 {
555  typedef typename TensorTools::DecrementRank<OutputType>::type Rank1Decrement;
556  typedef typename TensorTools::DecrementRank<Rank1Decrement>::type Rank2Decrement;
557  typedef typename TensorTools::MakeReal<Rank2Decrement>::type OutputShape;
558 
559  // Get local-to-global dof index lookup
560  const unsigned int n_dofs = cast_int<unsigned int>
561  (this->get_dof_indices(var).size());
562 
563  // Get current local coefficients
564  const DenseSubVector<Number> & coef = get_localized_subvector(_system_vector, var);
565 
566  // Get finite element object
567  FEGenericBase<OutputShape> * fe = nullptr;
568  this->get_element_fe<OutputShape>( var, fe, this->get_elem_dim() );
569 
570  // Get shape function values at quadrature point
571  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor>> & d2phi = fe->get_d2phi();
572 
573  // Loop over all the q_points in this finite element
574  for (auto qp : index_range(d2u_vals))
575  {
576  OutputType & d2u = d2u_vals[qp];
577 
578  // Compute the gradient at this q_point
579  d2u = 0;
580 
581  for (unsigned int l=0; l != n_dofs; l++)
582  d2u.add_scaled(d2phi[l][qp], coef(l));
583  }
584 
585  return;
586 }
587 
588 
589 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
590 
591 
592 template<typename OutputType>
593 void FEMContext::interior_curl(unsigned int var, unsigned int qp,
594  OutputType & curl_u) const
595 {
596  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
597 
598  // Get local-to-global dof index lookup
599  const unsigned int n_dofs = cast_int<unsigned int>
600  (this->get_dof_indices(var).size());
601 
602  // Get current local coefficients
603  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
604  const DenseSubVector<Number> & coef = this->get_elem_solution(var);
605 
606  // Get finite element object
607  FEGenericBase<OutputShape> * fe = nullptr;
608  this->get_element_fe<OutputShape>( var, fe, this->get_elem_dim() );
609 
610  // Get shape function values at quadrature point
611  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputShape>> & curl_phi = fe->get_curl_phi();
612 
613  // Accumulate solution curl
614  curl_u = 0.;
615 
616  for (unsigned int l=0; l != n_dofs; l++)
617  curl_u.add_scaled(curl_phi[l][qp], coef(l));
618 
619  return;
620 }
621 
622 
623 template<typename OutputType>
624 void FEMContext::interior_div(unsigned int var, unsigned int qp,
625  OutputType & div_u) const
626 {
627  typedef typename
629  <typename TensorTools::MakeReal<OutputType>::type>::type OutputShape;
630 
631  // Get local-to-global dof index lookup
632  const unsigned int n_dofs = cast_int<unsigned int>
633  (this->get_dof_indices(var).size());
634 
635  // Get current local coefficients
636  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
637  const DenseSubVector<Number> & coef = this->get_elem_solution(var);
638 
639  // Get finite element object
640  FEGenericBase<OutputShape> * fe = nullptr;
641  this->get_element_fe<OutputShape>( var, fe, this->get_elem_dim() );
642 
643  // Get shape function values at quadrature point
644  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputDivergence>> & div_phi = fe->get_div_phi();
645 
646  // Accumulate solution curl
647  div_u = 0.;
648 
649  for (unsigned int l=0; l != n_dofs; l++)
650  div_u += div_phi[l][qp] * coef(l);
651 
652  return;
653 }
654 
655 
657  unsigned int qp) const
658 {
659  Number u = 0.;
660 
661  this->side_value( var, qp, u );
662 
663  return u;
664 }
665 
666 
667 template<typename OutputType>
668 void FEMContext::side_value(unsigned int var,
669  unsigned int qp,
670  OutputType & u) const
671 {
672  this->some_value<OutputType,
673  &FEMContext::get_side_fe<typename TensorTools::MakeReal<OutputType>::type>,
674  &DiffContext::get_elem_solution>(var, qp, u);
675 }
676 
677 
678 template<typename OutputType>
679 void FEMContext::side_values(unsigned int var,
680  const NumericVector<Number> & _system_vector,
681  std::vector<OutputType> & u_vals) const
682 {
683  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
684 
685  // Get local-to-global dof index lookup
686  const unsigned int n_dofs = cast_int<unsigned int>
687  (this->get_dof_indices(var).size());
688 
689  // Get current local coefficients
690  const DenseSubVector<Number> & coef = get_localized_subvector(_system_vector, var);
691 
692  // Get the finite element object
693  FEGenericBase<OutputShape> * the_side_fe = nullptr;
694  this->get_side_fe<OutputShape>( var, the_side_fe, this->get_elem_dim() );
695 
696  // Get shape function values at quadrature point
697  const std::vector<std::vector<OutputShape>> & phi = the_side_fe->get_phi();
698 
699  // Loop over all the q_points on this element
700  for (auto qp : index_range(u_vals))
701  {
702  OutputType & u = u_vals[qp];
703 
704  // Compute the value at this q_point
705  u = 0.;
706 
707  for (unsigned int l=0; l != n_dofs; l++)
708  u += phi[l][qp] * coef(l);
709  }
710 
711  return;
712 }
713 
714 Gradient FEMContext::side_gradient(unsigned int var, unsigned int qp) const
715 {
716  Gradient du;
717 
718  this->side_gradient( var, qp, du );
719 
720  return du;
721 }
722 
723 
724 template<typename OutputType>
725 void FEMContext::side_gradient(unsigned int var, unsigned int qp,
726  OutputType & du) const
727 {
728  typedef typename TensorTools::MakeReal
730  OutputShape;
731 
732  // Get local-to-global dof index lookup
733  const unsigned int n_dofs = cast_int<unsigned int>
734  (this->get_dof_indices(var).size());
735 
736  // Get current local coefficients
737  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
738  const DenseSubVector<Number> & coef = this->get_elem_solution(var);
739 
740  // Get finite element object
741  FEGenericBase<OutputShape> * the_side_fe = nullptr;
742  this->get_side_fe<OutputShape>( var, the_side_fe, this->get_elem_dim() );
743 
744  // Get shape function values at quadrature point
745  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient>> & dphi = the_side_fe->get_dphi();
746 
747  // Accumulate solution derivatives
748  du = 0.;
749 
750  for (unsigned int l=0; l != n_dofs; l++)
751  du.add_scaled(dphi[l][qp], coef(l));
752 
753  return;
754 }
755 
756 
757 
758 template<typename OutputType>
759 void FEMContext::side_gradients(unsigned int var,
760  const NumericVector<Number> & _system_vector,
761  std::vector<OutputType> & du_vals) const
762 {
763  typedef typename TensorTools::MakeReal
765  OutputShape;
766 
767  // Get local-to-global dof index lookup
768  const unsigned int n_dofs = cast_int<unsigned int>
769  (this->get_dof_indices(var).size());
770 
771  // Get current local coefficients
772  const DenseSubVector<Number> & coef = get_localized_subvector(_system_vector, var);
773 
774  // Get finite element object
775  FEGenericBase<OutputShape> * the_side_fe = nullptr;
776  this->get_side_fe<OutputShape>( var, the_side_fe, this->get_elem_dim() );
777 
778  // Get shape function values at quadrature point
779  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient>> & dphi = the_side_fe->get_dphi();
780 
781  // Loop over all the q_points in this finite element
782  for (auto qp : index_range(du_vals))
783  {
784  OutputType & du = du_vals[qp];
785 
786  du = 0;
787 
788  // Compute the gradient at this q_point
789  for (unsigned int l=0; l != n_dofs; l++)
790  du.add_scaled(dphi[l][qp], coef(l));
791  }
792 
793  return;
794 }
795 
796 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
798  unsigned int qp) const
799 {
800  Tensor d2u;
801 
802  this->side_hessian( var, qp, d2u );
803 
804  return d2u;
805 }
806 
807 
808 
809 template<typename OutputType>
810 void FEMContext::side_hessian(unsigned int var,
811  unsigned int qp,
812  OutputType & d2u) const
813 {
814  this->some_hessian<OutputType,
816  <typename TensorTools::MakeReal
819  <OutputType>::type>::type>::type>,
820  &DiffContext::get_elem_solution>(var, qp, d2u);
821 }
822 
823 
824 
825 template<typename OutputType>
826 void FEMContext::side_hessians(unsigned int var,
827  const NumericVector<Number> & _system_vector,
828  std::vector<OutputType> & d2u_vals) const
829 {
830  typedef typename TensorTools::DecrementRank<OutputType>::type Rank1Decrement;
831  typedef typename TensorTools::DecrementRank<Rank1Decrement>::type Rank2Decrement;
832  typedef typename TensorTools::MakeReal<Rank2Decrement>::type OutputShape;
833 
834  // Get local-to-global dof index lookup
835  const unsigned int n_dofs = cast_int<unsigned int>
836  (this->get_dof_indices(var).size());
837 
838  // Get current local coefficients
839  const DenseSubVector<Number> & coef = get_localized_subvector(_system_vector, var);
840 
841  // Get finite element object
842  FEGenericBase<OutputShape> * the_side_fe = nullptr;
843  this->get_side_fe<OutputShape>( var, the_side_fe, this->get_elem_dim() );
844 
845  // Get shape function values at quadrature point
846  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor>> & d2phi = the_side_fe->get_d2phi();
847 
848  // Loop over all the q_points in this finite element
849  for (auto qp : index_range(d2u_vals))
850  {
851  OutputType & d2u = d2u_vals[qp];
852 
853  // Compute the gradient at this q_point
854  d2u = 0;
855 
856  for (unsigned int l=0; l != n_dofs; l++)
857  d2u.add_scaled(d2phi[l][qp], coef(l));
858  }
859 
860  return;
861 }
862 
863 
864 
865 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
866 
867 
868 
869 Number FEMContext::point_value(unsigned int var, const Point & p) const
870 {
871  Number u = 0.;
872 
873  this->point_value( var, p, u );
874 
875  return u;
876 }
877 
878 template<typename OutputType>
879 void FEMContext::point_value(unsigned int var,
880  const Point & p,
881  OutputType & u,
882  const Real tolerance) const
883 {
884  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
885 
886  // Get local-to-global dof index lookup
887  const unsigned int n_dofs = cast_int<unsigned int>
888  (this->get_dof_indices(var).size());
889 
890  // Get current local coefficients
891  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
892  const DenseSubVector<Number> & coef = this->get_elem_solution(var);
893 
894  // Get finite element object
895  FEGenericBase<OutputShape> * fe = nullptr;
896  this->get_element_fe<OutputShape>( var, fe, this->get_elem_dim() );
897 
898  // Build a FE for calculating u(p)
899  FEGenericBase<OutputShape> * fe_new =
900  this->build_new_fe( fe, p, tolerance, 0 );
901 
902  // Get the values of the shape function derivatives
903  const std::vector<std::vector<OutputShape>> & phi = fe_new->get_phi();
904 
905  u = 0.;
906 
907  for (unsigned int l=0; l != n_dofs; l++)
908  u += phi[l][0] * coef(l);
909 
910  return;
911 }
912 
913 
914 
915 Gradient FEMContext::point_gradient(unsigned int var, const Point & p) const
916 {
917  Gradient grad_u;
918 
919  this->point_gradient( var, p, grad_u );
920 
921  return grad_u;
922 }
923 
924 
925 
926 template<typename OutputType>
927 void FEMContext::point_gradient(unsigned int var,
928  const Point & p,
929  OutputType & grad_u,
930  const Real tolerance) const
931 {
932  typedef typename TensorTools::MakeReal
934  OutputShape;
935 
936  // Get local-to-global dof index lookup
937  const unsigned int n_dofs = cast_int<unsigned int>
938  (this->get_dof_indices(var).size());
939 
940  // Get current local coefficients
941  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
942  const DenseSubVector<Number> & coef = this->get_elem_solution(var);
943 
944  // Get finite element object
945  FEGenericBase<OutputShape> * fe = nullptr;
946  this->get_element_fe<OutputShape>( var, fe, this->get_elem_dim() );
947 
948  // Build a FE for calculating u(p)
949  FEGenericBase<OutputShape> * fe_new =
950  this->build_new_fe( fe, p, tolerance, 1 );
951 
952  // Get the values of the shape function derivatives
953  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient>> & dphi = fe_new->get_dphi();
954 
955  grad_u = 0.0;
956 
957  for (unsigned int l=0; l != n_dofs; l++)
958  grad_u.add_scaled(dphi[l][0], coef(l));
959 
960  return;
961 }
962 
963 
964 
965 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
966 
967 Tensor FEMContext::point_hessian(unsigned int var, const Point & p) const
968 {
969  Tensor hess_u;
970 
971  this->point_hessian( var, p, hess_u );
972 
973  return hess_u;
974 }
975 
976 
977 template<typename OutputType>
978 void FEMContext::point_hessian(unsigned int var,
979  const Point & p,
980  OutputType & hess_u,
981  const Real tolerance) const
982 {
983  typedef typename TensorTools::DecrementRank<OutputType>::type Rank1Decrement;
984  typedef typename TensorTools::DecrementRank<Rank1Decrement>::type Rank2Decrement;
985  typedef typename TensorTools::MakeReal<Rank2Decrement>::type OutputShape;
986 
987  // Get local-to-global dof index lookup
988  const unsigned int n_dofs = cast_int<unsigned int>
989  (this->get_dof_indices(var).size());
990 
991  // Get current local coefficients
992  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
993  const DenseSubVector<Number> & coef = this->get_elem_solution(var);
994 
995  // Get finite element object
996  FEGenericBase<OutputShape> * fe = nullptr;
997  this->get_element_fe<OutputShape>( var, fe, this->get_elem_dim() );
998 
999  // Build a FE for calculating u(p)
1000  FEGenericBase<OutputShape> * fe_new =
1001  this->build_new_fe( fe, p, tolerance, 2 );
1002 
1003  // Get the values of the shape function derivatives
1004  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor>> & d2phi = fe_new->get_d2phi();
1005 
1006  hess_u = 0.0;
1007 
1008  for (unsigned int l=0; l != n_dofs; l++)
1009  hess_u.add_scaled(d2phi[l][0], coef(l));
1010 
1011  return;
1012 }
1013 
1014 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
1015 
1016 
1017 template<typename OutputType>
1018 void FEMContext::point_curl(unsigned int var,
1019  const Point & p,
1020  OutputType & curl_u,
1021  const Real tolerance) const
1022 {
1023  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
1024 
1025  // Get local-to-global dof index lookup
1026  const unsigned int n_dofs = cast_int<unsigned int>
1027  (this->get_dof_indices(var).size());
1028 
1029  // Get current local coefficients
1030  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
1031  const DenseSubVector<Number> & coef = this->get_elem_solution(var);
1032 
1033  // Get finite element object
1034  FEGenericBase<OutputShape> * fe = nullptr;
1035  this->get_element_fe<OutputShape>( var, fe, this->get_elem_dim() );
1036 
1037  // Build a FE for calculating u(p)
1038  FEGenericBase<OutputShape> * fe_new =
1039  this->build_new_fe( fe, p, tolerance, 3 );
1040 
1041  // Get the values of the shape function derivatives
1042  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputShape>> & curl_phi = fe_new->get_curl_phi();
1043 
1044  curl_u = 0.0;
1045 
1046  for (unsigned int l=0; l != n_dofs; l++)
1047  curl_u.add_scaled(curl_phi[l][0], coef(l));
1048 
1049  return;
1050 }
1051 
1052 
1053 
1054 Number FEMContext::fixed_interior_value(unsigned int var, unsigned int qp) const
1055 {
1056  Number u = 0.;
1057 
1058  this->fixed_interior_value( var, qp, u );
1059 
1060  return u;
1061 }
1062 
1063 
1064 
1065 template<typename OutputType>
1066 void FEMContext::fixed_interior_value(unsigned int var, unsigned int qp,
1067  OutputType & u) const
1068 {
1069  this->some_value<OutputType,
1073 }
1074 
1075 
1076 
1077 Gradient FEMContext::fixed_interior_gradient(unsigned int var, unsigned int qp) const
1078 {
1079  Gradient du;
1080 
1081  this->fixed_interior_gradient( var, qp, du );
1082 
1083  return du;
1084 }
1085 
1086 
1087 template<typename OutputType>
1088 void FEMContext::fixed_interior_gradient(unsigned int var, unsigned int qp,
1089  OutputType & du) const
1090 {
1091  this->some_gradient
1092  <OutputType,
1094  <typename TensorTools::MakeReal
1095  <typename TensorTools::DecrementRank
1096  <OutputType>::type>::type>,
1098  (var, qp, du);
1099 }
1100 
1101 
1102 
1103 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1104 Tensor FEMContext::fixed_interior_hessian(unsigned int var, unsigned int qp) const
1105 {
1106  Tensor d2u;
1107 
1108  this->fixed_interior_hessian( var, qp, d2u );
1109 
1110  return d2u;
1111 }
1112 
1113 
1114 template<typename OutputType>
1115 void FEMContext::fixed_interior_hessian(unsigned int var, unsigned int qp,
1116  OutputType & d2u) const
1117 {
1118  this->some_hessian<OutputType,
1120  <typename TensorTools::MakeReal
1121  <typename TensorTools::DecrementRank
1122  <typename TensorTools::DecrementRank
1123  <OutputType>::type>::type>::type>,
1124  &DiffContext::get_elem_fixed_solution>(var, qp, d2u);
1125 }
1126 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1127 
1128 
1129 
1130 Number FEMContext::fixed_side_value(unsigned int var, unsigned int qp) const
1131 {
1132  Number u = 0.;
1133 
1134  this->fixed_side_value( var, qp, u );
1135 
1136  return u;
1137 }
1138 
1139 
1140 template<typename OutputType>
1141 void FEMContext::fixed_side_value(unsigned int var, unsigned int qp,
1142  OutputType & u) const
1143 {
1144  this->some_value
1145  <OutputType,
1149  (var, qp, u);
1150 }
1151 
1152 
1153 
1154 Gradient FEMContext::fixed_side_gradient(unsigned int var, unsigned int qp) const
1155 {
1156  Gradient du;
1157 
1158  this->fixed_side_gradient( var, qp, du );
1159 
1160  return du;
1161 }
1162 
1163 
1164 template<typename OutputType>
1165 void FEMContext::fixed_side_gradient(unsigned int var, unsigned int qp,
1166  OutputType & du) const
1167 {
1168  this->some_gradient<OutputType,
1170  <typename TensorTools::MakeReal
1171  <typename TensorTools::DecrementRank
1172  <OutputType>::type>::type>,
1173  &DiffContext::get_elem_fixed_solution>(var, qp, du);
1174 }
1175 
1176 
1177 
1178 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1179 Tensor FEMContext::fixed_side_hessian(unsigned int var, unsigned int qp) const
1180 {
1181  Tensor d2u;
1182 
1183  this->fixed_side_hessian( var, qp, d2u );
1184 
1185  return d2u;
1186 }
1187 
1188 template<typename OutputType>
1189 void FEMContext::fixed_side_hessian(unsigned int var, unsigned int qp,
1190  OutputType & d2u) const
1191 {
1192  this->some_hessian<OutputType,
1194  <typename TensorTools::MakeReal
1195  <typename TensorTools::DecrementRank
1196  <typename TensorTools::DecrementRank
1197  <OutputType>::type>::type>::type>,
1198  &DiffContext::get_elem_fixed_solution>(var, qp, d2u);
1199 }
1200 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1201 
1202 
1203 
1204 Number FEMContext::fixed_point_value(unsigned int var, const Point & p) const
1205 {
1206  Number u = 0.;
1207 
1208  this->fixed_point_value( var, p, u );
1209 
1210  return u;
1211 }
1212 
1213 template<typename OutputType>
1214 void FEMContext::fixed_point_value(unsigned int var,
1215  const Point & p,
1216  OutputType & u,
1217  const Real tolerance) const
1218 {
1219  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
1220 
1221  // Get local-to-global dof index lookup
1222  const unsigned int n_dofs = cast_int<unsigned int>
1223  (this->get_dof_indices(var).size());
1224 
1225  // Get current local coefficients
1226  libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
1227  const DenseSubVector<Number> & coef = this->get_elem_fixed_solution(var);
1228 
1229  // Get finite element object
1230  FEGenericBase<OutputShape> * fe = nullptr;
1231  this->get_element_fe<OutputShape>( var, fe, this->get_elem_dim() );
1232 
1233  // Build a FE for calculating u(p)
1234  FEGenericBase<OutputShape> * fe_new =
1235  this->build_new_fe( fe, p, tolerance, 0 );
1236 
1237  // Get the values of the shape function derivatives
1238  const std::vector<std::vector<OutputShape>> & phi = fe_new->get_phi();
1239 
1240  u = 0.;
1241 
1242  for (unsigned int l=0; l != n_dofs; l++)
1243  u += phi[l][0] * coef(l);
1244 
1245  return;
1246 }
1247 
1248 
1249 
1250 Gradient FEMContext::fixed_point_gradient(unsigned int var, const Point & p) const
1251 {
1252  Gradient grad_u;
1253 
1254  this->fixed_point_gradient( var, p, grad_u );
1255 
1256  return grad_u;
1257 }
1258 
1259 
1260 
1261 template<typename OutputType>
1262 void FEMContext::fixed_point_gradient(unsigned int var,
1263  const Point & p,
1264  OutputType & grad_u,
1265  const Real tolerance) const
1266 {
1267  typedef typename TensorTools::MakeReal
1269  OutputShape;
1270 
1271  // Get local-to-global dof index lookup
1272  const unsigned int n_dofs = cast_int<unsigned int>
1273  (this->get_dof_indices(var).size());
1274 
1275  // Get current local coefficients
1276  libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
1277  const DenseSubVector<Number> & coef = this->get_elem_fixed_solution(var);
1278 
1279  // Get finite element object
1280  FEGenericBase<OutputShape> * fe = nullptr;
1281  this->get_element_fe<OutputShape>( var, fe, this->get_elem_dim() );
1282 
1283  // Build a FE for calculating u(p)
1284  FEGenericBase<OutputShape> * fe_new =
1285  this->build_new_fe( fe, p, tolerance, 1 );
1286 
1287  // Get the values of the shape function derivatives
1288  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient>> & dphi = fe_new->get_dphi();
1289 
1290  grad_u = 0.0;
1291 
1292  for (unsigned int l=0; l != n_dofs; l++)
1293  grad_u.add_scaled(dphi[l][0], coef(l));
1294 
1295  return;
1296 }
1297 
1298 
1299 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1300 
1301 Tensor FEMContext::fixed_point_hessian(unsigned int var, const Point & p) const
1302 {
1303  Tensor hess_u;
1304 
1305  this->fixed_point_hessian( var, p, hess_u );
1306 
1307  return hess_u;
1308 }
1309 
1310 
1311 
1312 template<typename OutputType>
1313 void FEMContext::fixed_point_hessian(unsigned int var,
1314  const Point & p,
1315  OutputType & hess_u,
1316  const Real tolerance) const
1317 {
1318  typedef typename TensorTools::DecrementRank<OutputType>::type Rank1Decrement;
1319  typedef typename TensorTools::DecrementRank<Rank1Decrement>::type Rank2Decrement;
1320  typedef typename TensorTools::MakeReal<Rank2Decrement>::type OutputShape;
1321 
1322  // Get local-to-global dof index lookup
1323  const unsigned int n_dofs = cast_int<unsigned int>
1324  (this->get_dof_indices(var).size());
1325 
1326  // Get current local coefficients
1327  libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
1328  const DenseSubVector<Number> & coef = this->get_elem_fixed_solution(var);
1329 
1330  // Get finite element object
1331  FEGenericBase<OutputShape> * fe = nullptr;
1332  this->get_element_fe<OutputShape>( var, fe, this->get_elem_dim() );
1333 
1334  // Build a FE for calculating u(p)
1335  FEGenericBase<OutputShape> * fe_new =
1336  this->build_new_fe( fe, p, tolerance, 2 );
1337 
1338  // Get the values of the shape function derivatives
1339  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor>> & d2phi = fe_new->get_d2phi();
1340 
1341  hess_u = 0.0;
1342 
1343  for (unsigned int l=0; l != n_dofs; l++)
1344  hess_u.add_scaled(d2phi[l][0], coef(l));
1345 
1346  return;
1347 }
1348 
1349 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
1350 
1351 
1352 
1353 template<typename OutputType>
1354 void FEMContext::interior_rate(unsigned int var, unsigned int qp,
1355  OutputType & u) const
1356 {
1357  this->some_value<OutputType,
1361 }
1362 
1363 template<typename OutputType>
1364 void FEMContext::interior_rate_gradient(unsigned int var, unsigned int qp,
1365  OutputType & dudot) const
1366 {
1367  this->some_gradient<OutputType,
1369  <typename TensorTools::DecrementRank
1370  <OutputType>::type>::type>,
1371  &DiffContext::get_elem_solution_rate>(var, qp, dudot);
1372 }
1373 
1374 template<typename OutputType>
1375 void FEMContext::side_rate(unsigned int var, unsigned int qp,
1376  OutputType & u) const
1377 {
1378  this->some_value<OutputType,
1382 }
1383 
1384 template<typename OutputType>
1385 void FEMContext::interior_accel(unsigned int var, unsigned int qp,
1386  OutputType & u) const
1387 {
1388  this->some_value<OutputType,
1392 }
1393 
1394 
1395 
1396 template<typename OutputType>
1397 void FEMContext::side_accel(unsigned int var, unsigned int qp,
1398  OutputType & u) const
1399 {
1400  this->some_value<OutputType,
1404 }
1405 
1406 
1407 
1409 {
1410  // Update the "time" variable of this context object
1411  this->_update_time_from_system(theta);
1412 
1413  // Handle a moving element if necessary.
1414  if (_mesh_sys)
1415  {
1416  // We assume that the ``default'' state
1417  // of the mesh is its final, theta=1.0
1418  // position, so we don't bother with
1419  // mesh motion in that case.
1420  if (theta != 1.0)
1421  {
1422  // FIXME - ALE is not threadsafe yet!
1423  libmesh_assert_equal_to (libMesh::n_threads(), 1);
1424 
1425  elem_position_set(theta);
1426  }
1427  elem_fe_reinit();
1428  }
1429 }
1430 
1431 
1433 {
1434  // Update the "time" variable of this context object
1435  this->_update_time_from_system(theta);
1436 
1437  // Handle a moving element if necessary
1438  if (_mesh_sys)
1439  {
1440  // FIXME - not threadsafe yet!
1441  elem_position_set(theta);
1442  side_fe_reinit();
1443  }
1444 }
1445 
1446 
1448 {
1449  // Update the "time" variable of this context object
1450  this->_update_time_from_system(theta);
1451 
1452  // Handle a moving element if necessary
1453  if (_mesh_sys)
1454  {
1455  // FIXME - not threadsafe yet!
1456  elem_position_set(theta);
1457  edge_fe_reinit();
1458  }
1459 }
1460 
1461 
1463 {
1464  // Update the "time" variable of this context object
1465  this->_update_time_from_system(theta);
1466 
1467  // We can reuse the Elem FE safely here.
1468  elem_fe_reinit();
1469 }
1470 
1471 
1472 void FEMContext::elem_fe_reinit(const std::vector<Point> * const pts)
1473 {
1474  // Initialize all the interior FE objects on elem.
1475  // Logging of FE::reinit is done in the FE functions
1476  // We only reinit the FE objects for the current element
1477  // dimension
1478  const unsigned char dim = this->get_elem_dim();
1479 
1480  libmesh_assert( !_element_fe[dim].empty() );
1481 
1482  for (const auto & pr : _element_fe[dim])
1483  {
1484  if (this->has_elem())
1485  pr.second->reinit(&(this->get_elem()), pts);
1486  else
1487  // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
1488  pr.second->reinit(nullptr);
1489  }
1490 }
1491 
1492 
1494 {
1495  // Initialize all the side FE objects on elem/side.
1496  // Logging of FE::reinit is done in the FE functions
1497  // We only reinit the FE objects for the current element
1498  // dimension
1499  const unsigned char dim = this->get_elem_dim();
1500 
1501  libmesh_assert( !_side_fe[dim].empty() );
1502 
1503  for (auto & pr : _side_fe[dim])
1504  pr.second->reinit(&(this->get_elem()), this->get_side());
1505 }
1506 
1507 
1508 
1510 {
1511  libmesh_assert_equal_to (this->get_elem_dim(), 3);
1512 
1513  // Initialize all the interior FE objects on elem/edge.
1514  // Logging of FE::reinit is done in the FE functions
1515  for (auto & pr : _edge_fe)
1516  pr.second->edge_reinit(&(this->get_elem()), this->get_edge());
1517 }
1518 
1519 
1520 
1522 {
1523  // This is too expensive to call unless we've been asked to move the mesh
1525 
1526  // This will probably break with threading when two contexts are
1527  // operating on elements which share a node
1528  libmesh_assert_equal_to (libMesh::n_threads(), 1);
1529 
1530  // If the coordinate data is in our own system, it's already
1531  // been set up for us
1532  // if (_mesh_sys == this->number())
1533  // {
1534  unsigned int n_nodes = this->get_elem().n_nodes();
1535 
1536 #ifndef NDEBUG
1537  const unsigned char dim = this->get_elem_dim();
1538 
1539  // For simplicity we demand that mesh coordinates be stored
1540  // in a format that allows a direct copy
1542  (this->get_element_fe(this->get_mesh_x_var(), dim)->get_fe_type().family
1543  == FEMap::map_fe_type(this->get_elem()) &&
1544  this->get_element_fe(this->get_mesh_x_var(), dim)->get_fe_type().order.get_order()
1545  == this->get_elem().default_order()));
1547  (this->get_element_fe(this->get_mesh_y_var(), dim)->get_fe_type().family
1548  == FEMap::map_fe_type(this->get_elem()) &&
1549  this->get_element_fe(this->get_mesh_y_var(), dim)->get_fe_type().order.get_order()
1550  == this->get_elem().default_order()));
1552  (this->get_element_fe(this->get_mesh_z_var(), dim)->get_fe_type().family
1553  == FEMap::map_fe_type(this->get_elem()) &&
1554  this->get_element_fe(this->get_mesh_z_var(), dim)->get_fe_type().order.get_order()
1555  == this->get_elem().default_order()));
1556 #endif
1557 
1558  // Get degree of freedom coefficients from point coordinates
1559  if (this->get_mesh_x_var() != libMesh::invalid_uint)
1560  for (unsigned int i=0; i != n_nodes; ++i)
1561  (this->get_elem_solution(this->get_mesh_x_var()))(i) = this->get_elem().point(i)(0);
1562 
1563  if (this->get_mesh_y_var() != libMesh::invalid_uint)
1564  for (unsigned int i=0; i != n_nodes; ++i)
1565  (this->get_elem_solution(this->get_mesh_y_var()))(i) = this->get_elem().point(i)(1);
1566 
1567  if (this->get_mesh_z_var() != libMesh::invalid_uint)
1568  for (unsigned int i=0; i != n_nodes; ++i)
1569  (this->get_elem_solution(this->get_mesh_z_var()))(i) = this->get_elem().point(i)(2);
1570  // }
1571  // FIXME - If the coordinate data is not in our own system, someone
1572  // had better get around to implementing that... - RHS
1573  // else
1574  // {
1575  // libmesh_not_implemented();
1576  // }
1577 }
1578 
1579 
1580 
1582 {
1583  for (auto & m : _element_fe)
1584  for (auto & pr : m)
1585  pr.second->get_fe_map().set_jacobian_tolerance(tol);
1586 
1587  for (auto & m : _side_fe)
1588  for (auto & pr : m)
1589  pr.second->get_fe_map().set_jacobian_tolerance(tol);
1590 
1591  for (auto & pr : _edge_fe)
1592  pr.second->get_fe_map().set_jacobian_tolerance(tol);
1593 }
1594 
1595 
1596 
1597 // We can ignore the theta argument in the current use of this
1598 // function, because elem_subsolutions will already have been set to
1599 // the theta value.
1600 //
1601 // To enable loose mesh movement coupling things will need to change.
1603 {
1604  // This is too expensive to call unless we've been asked to move the mesh
1606 
1607  // This will probably break with threading when two contexts are
1608  // operating on elements which share a node
1609  libmesh_assert_equal_to (libMesh::n_threads(), 1);
1610 
1611  // If the coordinate data is in our own system, it's already
1612  // been set up for us, and we can ignore our input parameter theta
1613  // if (_mesh_sys == this->number())
1614  // {
1615  unsigned int n_nodes = this->get_elem().n_nodes();
1616 
1617 #ifndef NDEBUG
1618  const unsigned char dim = this->get_elem_dim();
1619 
1620  // For simplicity we demand that mesh coordinates be stored
1621  // in a format that allows a direct copy
1623  (this->get_element_fe(this->get_mesh_x_var(), dim)->get_fe_type().family
1624  == FEMap::map_fe_type(this->get_elem()) &&
1625  this->get_elem_solution(this->get_mesh_x_var()).size() == n_nodes));
1627  (this->get_element_fe(this->get_mesh_y_var(), dim)->get_fe_type().family
1628  == FEMap::map_fe_type(this->get_elem()) &&
1629  this->get_elem_solution(this->get_mesh_y_var()).size() == n_nodes));
1631  (this->get_element_fe(this->get_mesh_z_var(), dim)->get_fe_type().family
1632  == FEMap::map_fe_type(this->get_elem()) &&
1633  this->get_elem_solution(this->get_mesh_z_var()).size() == n_nodes));
1634 #endif
1635 
1636  // Set the new point coordinates
1637  if (this->get_mesh_x_var() != libMesh::invalid_uint)
1638  for (unsigned int i=0; i != n_nodes; ++i)
1639  const_cast<Elem &>(this->get_elem()).point(i)(0) =
1640  libmesh_real(this->get_elem_solution(this->get_mesh_x_var())(i));
1641 
1642  if (this->get_mesh_y_var() != libMesh::invalid_uint)
1643  for (unsigned int i=0; i != n_nodes; ++i)
1644  const_cast<Elem &>(this->get_elem()).point(i)(1) =
1645  libmesh_real(this->get_elem_solution(this->get_mesh_y_var())(i));
1646 
1647  if (this->get_mesh_z_var() != libMesh::invalid_uint)
1648  for (unsigned int i=0; i != n_nodes; ++i)
1649  const_cast<Elem &>(this->get_elem()).point(i)(2) =
1650  libmesh_real(this->get_elem_solution(this->get_mesh_z_var())(i));
1651  // }
1652  // FIXME - If the coordinate data is not in our own system, someone
1653  // had better get around to implementing that... - RHS
1654  // else
1655  // {
1656  // libmesh_not_implemented();
1657  // }
1658 }
1659 
1660 
1661 
1662 
1663 
1664 /*
1665  void FEMContext::reinit(const FEMSystem & sys, Elem * e)
1666  {
1667  // Initialize our elem pointer, algebraic objects
1668  this->pre_fe_reinit(e);
1669 
1670  // Moving the mesh may be necessary
1671  // Reinitializing the FE objects is definitely necessary
1672  this->elem_reinit(1.);
1673  }
1674 */
1675 
1676 
1677 
1678 void FEMContext::pre_fe_reinit(const System & sys, const Elem * e)
1679 {
1680  this->set_elem(e);
1681 
1682  if (algebraic_type() == CURRENT ||
1684  {
1685  // Initialize the per-element data for elem.
1686  if (this->has_elem())
1687  sys.get_dof_map().dof_indices (&(this->get_elem()), this->get_dof_indices());
1688  else
1689  // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
1690  sys.get_dof_map().dof_indices
1691  (static_cast<Elem*>(nullptr), this->get_dof_indices());
1692  }
1693 #ifdef LIBMESH_ENABLE_AMR
1694  else if (algebraic_type() == OLD ||
1696  {
1697  // Initialize the per-element data for elem.
1698  if (this->has_elem())
1699  sys.get_dof_map().old_dof_indices (&(this->get_elem()), this->get_dof_indices());
1700  else
1701  // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
1703  (static_cast<Elem*>(nullptr), this->get_dof_indices());
1704  }
1705 #endif // LIBMESH_ENABLE_AMR
1706 
1707  const unsigned int n_dofs = cast_int<unsigned int>
1708  (this->get_dof_indices().size());
1709  const unsigned int n_qoi = sys.n_qois();
1710 
1711  if (this->algebraic_type() != NONE &&
1712  this->algebraic_type() != DOFS_ONLY &&
1713  this->algebraic_type() != OLD_DOFS_ONLY)
1714  {
1715  // This also resizes elem_solution
1716  if (_custom_solution == nullptr)
1717  sys.current_local_solution->get(this->get_dof_indices(), this->get_elem_solution().get_values());
1718  else
1719  _custom_solution->get(this->get_dof_indices(), this->get_elem_solution().get_values());
1720 
1721  if (sys.use_fixed_solution)
1722  this->get_elem_fixed_solution().resize(n_dofs);
1723 
1724  // Only make space for these if we're using DiffSystem
1725  // This is assuming *only* DiffSystem is using elem_solution_rate/accel
1726  const DifferentiableSystem * diff_system = dynamic_cast<const DifferentiableSystem *>(&sys);
1727  if (diff_system)
1728  {
1729  // Now, we only need these if the solver is unsteady
1730  if (!diff_system->get_time_solver().is_steady())
1731  {
1732  this->get_elem_solution_rate().resize(n_dofs);
1733 
1734  // We only need accel space if the TimeSolver is second order
1735  const UnsteadySolver & time_solver = cast_ref<const UnsteadySolver &>(diff_system->get_time_solver());
1736 
1737  if (time_solver.time_order() >= 2 || !diff_system->get_second_order_vars().empty())
1738  this->get_elem_solution_accel().resize(n_dofs);
1739  }
1740  }
1741 
1742  if (algebraic_type() != OLD)
1743  {
1744  // These resize calls also zero out the residual and jacobian
1745  this->get_elem_residual().resize(n_dofs);
1746  if (this->_have_local_matrices)
1747  this->get_elem_jacobian().resize(n_dofs, n_dofs);
1748 
1749  this->get_qoi_derivatives().resize(n_qoi);
1750  this->_elem_qoi_subderivatives.resize(n_qoi);
1751  for (std::size_t q=0; q != n_qoi; ++q)
1752  (this->get_qoi_derivatives())[q].resize(n_dofs);
1753  }
1754  }
1755 
1756  // Initialize the per-variable data for elem.
1757  {
1758  unsigned int sub_dofs = 0;
1759  for (auto i : make_range(sys.n_vars()))
1760  {
1761  if (algebraic_type() == CURRENT ||
1763  {
1764  if (this->has_elem())
1765  sys.get_dof_map().dof_indices (&(this->get_elem()), this->get_dof_indices(i), i);
1766  else
1767  // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
1768  sys.get_dof_map().dof_indices
1769  (static_cast<Elem*>(nullptr), this->get_dof_indices(i), i);
1770  }
1771 #ifdef LIBMESH_ENABLE_AMR
1772  else if (algebraic_type() == OLD ||
1774  {
1775  if (this->has_elem())
1776  sys.get_dof_map().old_dof_indices (&(this->get_elem()), this->get_dof_indices(i), i);
1777  else
1778  // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
1780  (static_cast<Elem*>(nullptr), this->get_dof_indices(i), i);
1781  }
1782 #endif // LIBMESH_ENABLE_AMR
1783 
1784  if (this->algebraic_type() != NONE &&
1785  this->algebraic_type() != DOFS_ONLY &&
1786  this->algebraic_type() != OLD_DOFS_ONLY)
1787  {
1788  const unsigned int n_dofs_var = cast_int<unsigned int>
1789  (this->get_dof_indices(i).size());
1790 
1791 
1792  if (!_active_vars ||
1793  std::binary_search(_active_vars->begin(),
1794  _active_vars->end(), i))
1795  {
1796  this->get_elem_solution(i).reposition
1797  (sub_dofs, n_dofs_var);
1798 
1799  // Only make space for these if we're using DiffSystem
1800  // This is assuming *only* DiffSystem is using elem_solution_rate/accel
1801  const DifferentiableSystem * diff_system = dynamic_cast<const DifferentiableSystem *>(&sys);
1802  if (diff_system)
1803  {
1804  // Now, we only need these if the solver is unsteady
1805  if (!diff_system->get_time_solver().is_steady())
1806  {
1807  this->get_elem_solution_rate(i).reposition
1808  (sub_dofs, n_dofs_var);
1809 
1810  // We only need accel space if the TimeSolver is second order
1811  const UnsteadySolver & time_solver = cast_ref<const UnsteadySolver &>(diff_system->get_time_solver());
1812 
1813  if (time_solver.time_order() >= 2 || !diff_system->get_second_order_vars().empty())
1814  this->get_elem_solution_accel(i).reposition
1815  (sub_dofs, n_dofs_var);
1816  }
1817  }
1818 
1819  if (sys.use_fixed_solution)
1820  this->get_elem_fixed_solution(i).reposition
1821  (sub_dofs, n_dofs_var);
1822 
1823  if (algebraic_type() != OLD)
1824  {
1825  this->get_elem_residual(i).reposition
1826  (sub_dofs, n_dofs_var);
1827 
1828  for (std::size_t q=0; q != n_qoi; ++q)
1829  this->get_qoi_derivatives(q,i).reposition
1830  (sub_dofs, n_dofs_var);
1831 
1832  if (this->_have_local_matrices)
1833  {
1834  for (unsigned int j=0; j != i; ++j)
1835  {
1836  const unsigned int n_dofs_var_j =
1837  cast_int<unsigned int>
1838  (this->get_dof_indices(j).size());
1839 
1840  this->get_elem_jacobian(i,j).reposition
1841  (sub_dofs, this->get_elem_residual(j).i_off(),
1842  n_dofs_var, n_dofs_var_j);
1843  this->get_elem_jacobian(j,i).reposition
1844  (this->get_elem_residual(j).i_off(), sub_dofs,
1845  n_dofs_var_j, n_dofs_var);
1846  }
1847  this->get_elem_jacobian(i,i).reposition
1848  (sub_dofs, sub_dofs,
1849  n_dofs_var,
1850  n_dofs_var);
1851  }
1852  }
1853  }
1854 
1855  sub_dofs += n_dofs_var;
1856  }
1857  }
1858 
1859  if (this->algebraic_type() != NONE &&
1860  this->algebraic_type() != DOFS_ONLY &&
1861  this->algebraic_type() != OLD &&
1862  this->algebraic_type() != OLD_DOFS_ONLY)
1863  libmesh_assert_equal_to (sub_dofs, n_dofs);
1864  }
1865 
1866  // Now do the localization for the user requested vectors
1867  if (this->algebraic_type() != NONE &&
1868  this->algebraic_type() != DOFS_ONLY &&
1869  this->algebraic_type() != OLD_DOFS_ONLY)
1870  {
1871  DiffContext::localized_vectors_iterator localized_vec_it = this->_localized_vectors.begin();
1872  const DiffContext::localized_vectors_iterator localized_vec_end = this->_localized_vectors.end();
1873 
1874  for (; localized_vec_it != localized_vec_end; ++localized_vec_it)
1875  {
1876  const NumericVector<Number> & current_localized_vector = *localized_vec_it->first;
1877  DenseVector<Number> & target_vector = localized_vec_it->second.first;
1878 
1879  current_localized_vector.get(this->get_dof_indices(), target_vector.get_values());
1880 
1881  // Initialize the per-variable data for elem.
1882  unsigned int sub_dofs = 0;
1883  auto init_localized_var_data = [this, localized_vec_it, &sub_dofs](unsigned int i)
1884  {
1885  const unsigned int n_dofs_var = cast_int<unsigned int>
1886  (this->get_dof_indices(i).size());
1887 
1888  // This is redundant with earlier initialization, isn't it? - RHS
1889  // sys.get_dof_map().dof_indices (&(this->get_elem()), this->get_dof_indices(i), i);
1890 
1891  localized_vec_it->second.second[i].reposition
1892  (sub_dofs, n_dofs_var);
1893 
1894  sub_dofs += n_dofs_var;
1895  };
1896 
1897  if (_active_vars)
1898  for (auto v : *_active_vars)
1899  init_localized_var_data(v);
1900  else
1901  for (auto v : make_range(sys.n_vars()))
1902  init_localized_var_data(v);
1903 
1904  libmesh_assert_equal_to (sub_dofs, n_dofs);
1905  }
1906  }
1907 }
1908 
1909 void FEMContext::set_elem( const Elem * e )
1910 {
1911  this->_elem = e;
1912 
1913  // If e is nullptr, we assume it's SCALAR and set _elem_dim to 0.
1914  this->_elem_dim =
1915  cast_int<unsigned char>(this->_elem ? this->_elem->dim() : 0);
1916 }
1917 
1919 {
1920  // Update the "time" variable based on the value of theta. For this
1921  // to work, we need to know the value of deltat, a pointer to which is now
1922  // stored by our parent DiffContext class. Note: get_deltat_value() will
1923  // assert in debug mode if the requested pointer is nullptr.
1924  const Real deltat = this->get_deltat_value();
1925 
1926  this->set_time(theta*(this->get_system_time() + deltat) + (1.-theta)*this->get_system_time());
1927 }
1928 
1929 
1930 
1931 template<>
1933 FEMContext::cached_fe( const unsigned int elem_dim,
1934  const FEType fe_type,
1935  const int get_derivative_level ) const
1936 {
1937 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1938  const bool fe_needs_inf =
1939  this->has_elem() && this->get_elem().infinite();
1940 #endif
1941 
1942  if (!_real_fe ||
1943  elem_dim != _real_fe->get_dim() ||
1944  fe_type != _real_fe->get_fe_type() ||
1945  get_derivative_level != _real_fe_derivative_level)
1946  {
1947  _real_fe_derivative_level = get_derivative_level;
1948 
1949  _real_fe =
1950 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1951  fe_needs_inf ?
1952  FEGenericBase<Real>::build_InfFE(elem_dim, fe_type) :
1953 #endif
1954  FEGenericBase<Real>::build(elem_dim, fe_type);
1955  }
1956 
1957 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1958  else if (fe_needs_inf && !_real_fe_is_inf)
1959  _real_fe =
1960  FEGenericBase<Real>::build_InfFE(elem_dim, fe_type);
1961  else if (!fe_needs_inf && _real_fe_is_inf)
1962  _real_fe =
1963  FEGenericBase<Real>::build(elem_dim, fe_type);
1964 
1965  _real_fe_is_inf =
1966  (this->has_elem() && this->get_elem().infinite());
1967 #endif
1968 
1969  return _real_fe.get();
1970 }
1971 
1972 
1973 template<>
1975 FEMContext::cached_fe( const unsigned int elem_dim,
1976  const FEType fe_type,
1977  const int get_derivative_level ) const
1978 {
1979 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1980  const bool fe_needs_inf =
1981  this->has_elem() && this->get_elem().infinite();
1982 #endif
1983 
1984  if (!_real_grad_fe ||
1985  elem_dim != _real_grad_fe->get_dim() ||
1986  fe_type != _real_grad_fe->get_fe_type() ||
1987  get_derivative_level != _real_grad_fe_derivative_level)
1988  {
1989  _real_grad_fe_derivative_level = get_derivative_level;
1990 
1991  _real_grad_fe =
1992 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1993  fe_needs_inf ?
1994  FEGenericBase<RealGradient>::build_InfFE(elem_dim, fe_type) :
1995 #endif
1996  FEGenericBase<RealGradient>::build(elem_dim, fe_type);
1997  }
1998 
1999 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2000  else if (fe_needs_inf && !_real_grad_fe_is_inf)
2001  _real_grad_fe =
2002  FEGenericBase<RealGradient>::build_InfFE(elem_dim, fe_type);
2003  else if (!fe_needs_inf && _real_grad_fe_is_inf)
2004  _real_grad_fe =
2005  FEGenericBase<RealGradient>::build(elem_dim, fe_type);
2006 
2008  (this->has_elem() && this->get_elem().infinite());
2009 #endif
2010 
2011  return _real_grad_fe.get();
2012 }
2013 
2014 
2015 
2016 template<typename OutputShape>
2019  const Point & p,
2020  const Real tolerance,
2021  const int get_derivative_level) const
2022 {
2023  FEType fe_type = fe->get_fe_type();
2024 
2025  // If we don't have an Elem to evaluate on, then the only functions
2026  // we can sensibly evaluate are the scalar dofs which are the same
2027  // everywhere.
2028  libmesh_assert(this->has_elem() || fe_type.family == SCALAR);
2029 
2030 #ifdef LIBMESH_ENABLE_AMR
2031  const bool add_p_level = fe->add_p_level_in_reinit();
2032  if ((algebraic_type() == OLD) &&
2033  this->has_elem())
2034  {
2035  if (this->get_elem().p_refinement_flag() == Elem::JUST_REFINED)
2036  fe_type.order = static_cast<Order>(fe_type.order - add_p_level);
2037  else if (this->get_elem().p_refinement_flag() == Elem::JUST_COARSENED)
2038  fe_type.order = static_cast<Order>(fe_type.order + add_p_level);
2039  }
2040 #endif // LIBMESH_ENABLE_AMR
2041 
2042  const unsigned int elem_dim = this->has_elem() ? this->get_elem().dim() : 0;
2043 
2044  FEGenericBase<OutputShape>* fe_new =
2045  cached_fe<OutputShape>(elem_dim, fe_type, get_derivative_level);
2046 #ifdef LIBMESH_ENABLE_AMR
2047  fe_new->add_p_level_in_reinit(add_p_level);
2048 #endif // LIBMESH_ENABLE_AMR
2049 
2050  // Map the physical co-ordinates to the master co-ordinates using the inverse_map from fe_interface.h
2051  // Build a vector of point co-ordinates to send to reinit
2052  Point master_point = this->has_elem() ?
2053  FEMap::inverse_map (elem_dim, &this->get_elem(), p, tolerance) :
2054  Point(0);
2055 
2056  std::vector<Point> coor(1, master_point);
2057 
2058  switch (get_derivative_level)
2059  {
2060  case -1:
2061  fe_new->get_phi();
2062  fe_new->get_dphi();
2063 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
2064  fe_new->get_d2phi();
2065 #endif
2066  fe_new->get_curl_phi();
2067  break;
2068  case 0:
2069  fe_new->get_phi();
2070  break;
2071  case 1:
2072  fe_new->get_dphi();
2073  break;
2074  case 2:
2075 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
2076  fe_new->get_d2phi();
2077 #else
2078  // here a different configuration is required.
2079  libmesh_not_implemented();
2080 #endif
2081  break;
2082  case 3:
2083  fe_new->get_curl_phi();
2084  break;
2085  default:
2086  libmesh_error();
2087  }
2088 
2089  // Reinitialize the element and compute the shape function values at coor
2090  if (this->has_elem())
2091  fe_new->reinit (&this->get_elem(), &coor);
2092  else
2093  // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
2094  fe_new->reinit (nullptr, &coor);
2095 
2096  return fe_new;
2097 }
2098 
2099 
2100 
2101 
2102 
2103 // Instantiate member function templates
2104 template LIBMESH_EXPORT void FEMContext::interior_value<Number>(unsigned int, unsigned int, Number &) const;
2105 template LIBMESH_EXPORT void FEMContext::interior_values<Number>(unsigned int, const NumericVector<Number> &,
2106  std::vector<Number> &) const;
2107 template LIBMESH_EXPORT void FEMContext::interior_value<Gradient>(unsigned int, unsigned int, Gradient &) const;
2108 template LIBMESH_EXPORT void FEMContext::interior_values<Gradient>(unsigned int, const NumericVector<Number> &,
2109  std::vector<Gradient> &) const;
2110 
2111 template LIBMESH_EXPORT void FEMContext::interior_gradient<Gradient>(unsigned int, unsigned int, Gradient &) const;
2112 template LIBMESH_EXPORT void FEMContext::interior_gradients<Gradient>(unsigned int, const NumericVector<Number> &,
2113  std::vector<Gradient> &) const;
2114 template LIBMESH_EXPORT void FEMContext::interior_gradient<Tensor>(unsigned int, unsigned int, Tensor &) const;
2115 template LIBMESH_EXPORT void FEMContext::interior_gradients<Tensor>(unsigned int, const NumericVector<Number> &,
2116  std::vector<Tensor> &) const;
2117 
2118 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
2119 template LIBMESH_EXPORT void FEMContext::interior_hessian<Tensor>(unsigned int, unsigned int, Tensor &) const;
2120 template LIBMESH_EXPORT void FEMContext::interior_hessians<Tensor>(unsigned int, const NumericVector<Number> &,
2121  std::vector<Tensor> &) const;
2122 //FIXME: Not everything is implemented yet for second derivatives of RealGradients
2123 //template LIBMESH_EXPORT void FEMContext::interior_hessian<??>(unsigned int, unsigned int, ??&) const;
2124 //template LIBMESH_EXPORT void FEMContext::interior_hessians<??>(unsigned int, const NumericVector<Number> &,
2125 // std::vector<??> &) const;
2126 #endif
2127 
2128 template LIBMESH_EXPORT void FEMContext::interior_curl<Gradient>(unsigned int, unsigned int, Gradient &) const;
2129 
2130 template LIBMESH_EXPORT void FEMContext::interior_div<Number>(unsigned int, unsigned int, Number &) const;
2131 
2132 template LIBMESH_EXPORT void FEMContext::side_value<Number>(unsigned int, unsigned int, Number &) const;
2133 template LIBMESH_EXPORT void FEMContext::side_value<Gradient>(unsigned int, unsigned int, Gradient &) const;
2134 template LIBMESH_EXPORT void FEMContext::side_values<Number>(unsigned int, const NumericVector<Number> &,
2135  std::vector<Number> &) const;
2136 template LIBMESH_EXPORT void FEMContext::side_values<Gradient>(unsigned int, const NumericVector<Number> &,
2137  std::vector<Gradient> &) const;
2138 
2139 template LIBMESH_EXPORT void FEMContext::side_gradient<Gradient>(unsigned int, unsigned int, Gradient &) const;
2140 template LIBMESH_EXPORT void FEMContext::side_gradients<Gradient>(unsigned int, const NumericVector<Number> &,
2141  std::vector<Gradient> &) const;
2142 template LIBMESH_EXPORT void FEMContext::side_gradient<Tensor>(unsigned int, unsigned int, Tensor &) const;
2143 template LIBMESH_EXPORT void FEMContext::side_gradients<Tensor>(unsigned int, const NumericVector<Number> &,
2144  std::vector<Tensor> &) const;
2145 
2146 
2147 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
2148 template LIBMESH_EXPORT void FEMContext::side_hessian<Tensor>(unsigned int, unsigned int, Tensor &) const;
2149 template LIBMESH_EXPORT void FEMContext::side_hessians<Tensor>(unsigned int, const NumericVector<Number> &,
2150  std::vector<Tensor> &) const;
2151 //FIXME: Not everything is implemented yet for second derivatives of RealGradients
2152 //template LIBMESH_EXPORT void FEMContext::side_hessian<??>(unsigned int, unsigned int,
2153 // ??&) const;
2154 //template LIBMESH_EXPORT void FEMContext::side_hessians<??>(unsigned int, const NumericVector<Number> &,
2155 // std::vector<??> &) const;
2156 #endif
2157 
2158 template LIBMESH_EXPORT void FEMContext::point_value<Number>(unsigned int, const Point &, Number &, const Real) const;
2159 template LIBMESH_EXPORT void FEMContext::point_value<Gradient>(unsigned int, const Point &, Gradient &, const Real) const;
2160 
2161 template LIBMESH_EXPORT void FEMContext::point_gradient<Gradient>(unsigned int, const Point &, Gradient &, const Real) const;
2162 template LIBMESH_EXPORT void FEMContext::point_gradient<Tensor>(unsigned int, const Point &, Tensor &, const Real) const;
2163 
2164 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
2165 template LIBMESH_EXPORT void FEMContext::point_hessian<Tensor>(unsigned int, const Point &, Tensor &, const Real) const;
2166 //FIXME: Not everything is implemented yet for second derivatives of RealGradients
2167 //template LIBMESH_EXPORT void FEMContext::point_hessian<??>(unsigned int, const Point &, ??&) const;
2168 #endif
2169 
2170 template LIBMESH_EXPORT void FEMContext::point_curl<Gradient>(unsigned int, const Point &, Gradient &, const Real) const;
2171 
2172 template LIBMESH_EXPORT void FEMContext::fixed_interior_value<Number>(unsigned int, unsigned int, Number &) const;
2173 template LIBMESH_EXPORT void FEMContext::fixed_interior_value<Gradient>(unsigned int, unsigned int, Gradient &) const;
2174 
2175 template LIBMESH_EXPORT void FEMContext::fixed_interior_gradient<Gradient>(unsigned int, unsigned int, Gradient &) const;
2176 template LIBMESH_EXPORT void FEMContext::fixed_interior_gradient<Tensor>(unsigned int, unsigned int, Tensor &) const;
2177 
2178 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
2179 template LIBMESH_EXPORT void FEMContext::fixed_interior_hessian<Tensor>(unsigned int, unsigned int, Tensor &) const;
2180 //FIXME: Not everything is implemented yet for second derivatives of RealGradients
2181 //template LIBMESH_EXPORT void FEMContext::fixed_interior_hessian<??>(unsigned int, unsigned int, ??&) const;
2182 #endif
2183 
2184 template LIBMESH_EXPORT void FEMContext::fixed_side_value<Number>(unsigned int, unsigned int, Number &) const;
2185 template LIBMESH_EXPORT void FEMContext::fixed_side_value<Gradient>(unsigned int, unsigned int, Gradient &) const;
2186 
2187 template LIBMESH_EXPORT void FEMContext::fixed_side_gradient<Gradient>(unsigned int, unsigned int, Gradient &) const;
2188 template LIBMESH_EXPORT void FEMContext::fixed_side_gradient<Tensor>(unsigned int, unsigned int, Tensor &) const;
2189 
2190 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
2191 template LIBMESH_EXPORT void FEMContext::fixed_side_hessian<Tensor>(unsigned int, unsigned int, Tensor &) const;
2192 //FIXME: Not everything is implemented yet for second derivatives of RealGradients
2193 //template LIBMESH_EXPORT void FEMContext::fixed_side_hessian<??>(unsigned int, unsigned int, ??&) const;
2194 #endif
2195 
2196 template LIBMESH_EXPORT void FEMContext::fixed_point_value<Number>(unsigned int, const Point &, Number &, const Real) const;
2197 template LIBMESH_EXPORT void FEMContext::fixed_point_value<Gradient>(unsigned int, const Point &, Gradient &, const Real) const;
2198 
2199 template LIBMESH_EXPORT void FEMContext::fixed_point_gradient<Gradient>(unsigned int, const Point &, Gradient &, const Real) const;
2200 template LIBMESH_EXPORT void FEMContext::fixed_point_gradient<Tensor>(unsigned int, const Point &, Tensor &, const Real) const;
2201 
2202 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
2203 template LIBMESH_EXPORT void FEMContext::fixed_point_hessian<Tensor>(unsigned int, const Point &, Tensor &, const Real) const;
2204 //FIXME: Not everything is implemented yet for second derivatives of RealGradients
2205 //template LIBMESH_EXPORT void FEMContext::fixed_point_hessian<??>(unsigned int, const Point &, ??&) const;
2206 #endif
2207 
2208 template LIBMESH_EXPORT void FEMContext::interior_rate<Number>(unsigned int, unsigned int, Number &) const;
2209 template LIBMESH_EXPORT void FEMContext::interior_rate<Gradient>(unsigned int, unsigned int, Gradient &) const;
2210 
2211 template LIBMESH_EXPORT void FEMContext::interior_rate_gradient<Gradient>(unsigned int, unsigned int, Gradient &) const;
2212 template LIBMESH_EXPORT void FEMContext::interior_rate_gradient<Tensor>(unsigned int, unsigned int, Tensor &) const;
2213 
2214 template LIBMESH_EXPORT void FEMContext::side_rate<Number>(unsigned int, unsigned int, Number &) const;
2215 template LIBMESH_EXPORT void FEMContext::side_rate<Gradient>(unsigned int, unsigned int, Gradient &) const;
2216 
2217 template LIBMESH_EXPORT void FEMContext::interior_accel<Number>(unsigned int, unsigned int, Number &) const;
2218 template LIBMESH_EXPORT void FEMContext::interior_accel<Gradient>(unsigned int, unsigned int, Gradient &) const;
2219 
2220 template LIBMESH_EXPORT void FEMContext::side_accel<Number>(unsigned int, unsigned int, Number &) const;
2221 template LIBMESH_EXPORT void FEMContext::side_accel<Gradient>(unsigned int, unsigned int, Gradient &) const;
2222 
2223 template LIBMESH_EXPORT FEGenericBase<Real> *
2225  const Point &,
2226  const Real,
2227  const int) const;
2228 
2229 template LIBMESH_EXPORT FEGenericBase<RealGradient> *
2231  const Point &,
2232  const Real,
2233  const int) const;
2234 
2235 } // namespace libMesh
T libmesh_real(T a)
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:182
void interior_curl(unsigned int var, unsigned int qp, OutputType &curl_u) const
Definition: fem_context.C:593
FEFamily family
The type of finite element.
Definition: fe_type.h:207
unsigned char get_edge() const
Accessor for current edge of Elem object.
Definition: fem_context.h:928
virtual void nonlocal_reinit(Real theta) override
Gives derived classes the opportunity to reinitialize data needed for nonlocal calculations at a new ...
Definition: fem_context.C:1462
FEGenericBase< OutputShape > * cached_fe(const unsigned int elem_dim, const FEType fe_type, const int get_derivative_level) const
virtual_for_inffe const std::vector< std::vector< OutputDivergence > > & get_div_phi() const
Definition: fe_base.h:261
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
Definition: diff_context.h:274
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
Number fixed_side_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1130
Number side_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:656
Tensor fixed_point_hessian(unsigned int var, const Point &p) const
Definition: fem_context.C:1301
Number interior_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:407
unsigned int n_threads()
Definition: libmesh_base.h:96
void elem_position_set(Real theta)
Uses the coordinate data specified by mesh_*_position configuration to set the geometry of elem to th...
Definition: fem_context.h:1270
std::vector< DenseSubVector< Number > > _elem_fixed_subsolutions
Definition: diff_context.h:605
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
void get_side_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for edge/face (2D/3D) finite element object for variable var for the largest dimension in th...
Definition: fem_context.h:317
virtual void pre_fe_reinit(const System &, const Elem *e)
Reinitializes local data vectors/matrices on the current geometric element.
Definition: fem_context.C:1678
bool has_boundary_id(const Node *const node, const boundary_id_type id) const
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
unsigned int get_mesh_x_var() const
Accessor for x-variable of moving mesh System.
Definition: fem_context.h:860
std::unique_ptr< const std::vector< unsigned int > > _active_vars
Variables on which to enable calculations, or nullptr if all variables in the System are to be enable...
Definition: fem_context.h:1054
const NumericVector< Number > * _custom_solution
Data with which to do algebra reinitialization.
Definition: fem_context.h:1074
void set_elem(const Elem *e)
Helper function to promote accessor usage.
Definition: fem_context.C:1909
std::unique_ptr< FEGenericBase< RealGradient > > _real_grad_fe
Definition: fem_context.h:1077
virtual void get(const std::vector< numeric_index_type > &index, T *values) const
Access multiple components at once.
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
const DenseVector< Number > & get_elem_fixed_solution() const
Accessor for element fixed solution.
Definition: diff_context.h:210
const DenseVector< Number > & get_elem_solution_rate() const
Accessor for element solution rate of change w.r.t.
Definition: diff_context.h:144
unsigned char get_elem_dim() const
Definition: fem_context.h:944
const Elem & get_elem() const
Accessor for current Elem object.
Definition: fem_context.h:908
unsigned int dim
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
FEMContext(const System &sys, const std::vector< unsigned int > *active_vars=nullptr, bool allocate_local_matrices=true)
Constructor.
Definition: fem_context.C:39
unsigned int n_qois() const
Number of currently active quantities of interest.
Definition: system.h:2516
void side_hessians(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &d2u_vals) const
Fills a vector of hessians of the _system_vector at the all the quadrature points on the current elem...
Definition: fem_context.C:826
void interior_accel(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1385
TensorTools::MakeReal< OutputType >::type value_shape
Definition: fem_context.h:1117
const DenseSubVector< Number > &(DiffContext::* diff_subsolution_getter)(unsigned int) const
Helper typedef to simplify refactoring.
Definition: fem_context.h:1103
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:374
RefinementState p_refinement_flag() const
Definition: elem.h:3063
void interior_gradients(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &interior_gradients_vector) const
Fills a vector with the gradient of the solution variable var at all the quadrature points in the cur...
Definition: fem_context.C:489
std::vector< T > & get_values()
Definition: dense_vector.h:268
virtual void elem_edge_reinit(Real theta) override
Resets the current time in the context.
Definition: fem_context.C:1447
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
OrderWrapper order
The approximation order of the element.
Definition: fe_type.h:201
int _extra_quadrature_order
The extra quadrature order for this context.
Definition: fem_context.h:1238
void boundary_ids(const Node *node, std::vector< boundary_id_type > &vec_to_fill) const
Fills a user-provided std::vector with the boundary ids associated with Node node.
void should_p_refine(unsigned int g, bool p_refine)
Describe whether the given variable group should be p-refined.
Definition: dof_map.h:2322
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
std::vector< std::vector< DenseSubVector< Number > > > _elem_qoi_subderivatives
Definition: diff_context.h:627
virtual bool is_steady() const =0
Is this effectively a steady-state solver?
void elem_position_get()
Uses the geometry of elem to set the coordinate data specified by mesh_*_position configuration...
Definition: fem_context.C:1521
The libMesh namespace provides an interface to certain functionality in the library.
void interior_values(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &interior_values_vector) const
Fills a vector of values of the _system_vector at the all the quadrature points in the current elemen...
Definition: fem_context.C:427
bool has_side_boundary_id(boundary_id_type id) const
Reports if the boundary id is found on the current side.
Definition: fem_context.C:298
void use_default_quadrature_rules(int extra_quadrature_order=0)
Use quadrature rules designed to over-integrate a mass matrix, plus extra_quadrature_order.
Definition: fem_context.C:155
virtual ~FEMContext()
Destructor.
Definition: fem_context.C:292
void side_boundary_ids(std::vector< boundary_id_type > &vec_to_fill) const
As above, but fills in the std::set provided by the user.
Definition: fem_context.C:305
unsigned char get_side() const
Accessor for current side of Elem object.
Definition: fem_context.h:922
Number point_value(unsigned int var, const Point &p) const
Definition: fem_context.C:869
std::vector< std::map< FEType, std::unique_ptr< FEAbstract > > > _element_fe
Finite element objects for each variable&#39;s interior, sides and edges.
Definition: fem_context.h:1168
virtual void elem_fe_reinit(const std::vector< Point > *const pts=nullptr)
Reinitializes interior FE objects on the current geometric element.
Definition: fem_context.C:1472
const bool _have_local_matrices
Whether we have local matrices allocated/initialized.
Definition: diff_context.h:576
std::unique_ptr< FEGenericBase< Real > > _real_fe
Definition: fem_context.h:1076
Number fixed_interior_value(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1054
Tnew cast_int(Told oldvar)
void interior_rate_gradient(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1364
Defines a dense subvector for use in finite element computations.
FEType get_fe_type() const
Definition: fe_abstract.h:499
This class provides a specific system class.
Definition: diff_system.h:54
Tensor fixed_side_hessian(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1179
std::unique_ptr< QBase > default_quadrature_rule(const unsigned int dim, const int extraorder=0) const
Definition: fe_type.C:34
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Definition: fe_base.h:230
std::set< unsigned char > _elem_dims
Cached dimensions of elements in the mesh, plus dimension 0 if SCALAR variables are in use...
Definition: fem_context.h:1208
Gradient interior_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:462
void point_curl(unsigned int var, const Point &p, OutputType &curl_u, const Real tolerance=TOLERANCE) const
Definition: fem_context.C:1018
const dof_id_type n_nodes
Definition: tecplot_io.C:67
virtual unsigned int time_order() const =0
const std::set< unsigned int > & get_second_order_vars() const
Definition: diff_physics.h:519
const DenseVector< Number > & get_elem_solution() const
Accessor for element solution.
Definition: diff_context.h:112
int8_t boundary_id_type
Definition: id_types.h:51
unsigned char _elem_dim
Cached dimension of this->_elem.
Definition: fem_context.h:1202
void side_gradients(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &side_gradients_vector) const
Fills a vector with the gradient of the solution variable var at all the quadrature points on the cur...
Definition: fem_context.C:759
Tensor side_hessian(unsigned int var, unsigned int qp) const
Definition: fem_context.C:797
virtual void side_fe_reinit()
Reinitializes side FE objects on the current geometric element.
Definition: fem_context.C:1493
std::map< FEType, std::unique_ptr< FEAbstract > > _edge_fe
Definition: fem_context.h:1170
const System & get_system() const
Accessor for associated system.
Definition: diff_context.h:106
void _update_time_from_system(Real theta)
Update the time in the context object for the given value of theta, based on the values of "time" and...
Definition: fem_context.C:1918
virtual unsigned int n_nodes() const =0
void init_internal_data(const System &sys)
Helper function used in constructors to set up internal data.
Definition: fem_context.C:198
Gradient fixed_point_gradient(unsigned int var, const Point &p) const
Definition: fem_context.C:1250
bool use_fixed_solution
A boolean to be set to true by systems using elem_fixed_solution, for optional use by e...
Definition: system.h:1543
void side_accel(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1397
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
void set_jacobian_tolerance(Real tol)
Calls set_jacobian_tolerance() on all the FE objects controlled by this class.
Definition: fem_context.C:1581
System * _mesh_sys
System from which to acquire moving mesh information.
Definition: fem_context.h:1059
void interior_hessians(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &d2u_vals) const
Fills a vector of hessians of the _system_vector at the all the quadrature points in the current elem...
Definition: fem_context.C:551
void some_value(unsigned int var, unsigned int qp, OutputType &u) const
Helper function to reduce some code duplication in the *interior_value methods.
Definition: fem_context.C:315
void some_hessian(unsigned int var, unsigned int qp, OutputType &u) const
Helper function to reduce some code duplication in the *interior_hessian methods. ...
Definition: fem_context.C:377
FEType find_hardest_fe_type()
Helper function for creating quadrature rules.
Definition: fem_context.C:86
void add_p_level_in_reinit(bool value)
Indicate whether to add p-refinement levels in init/reinit methods.
Definition: fe_abstract.h:610
Tensor interior_hessian(unsigned int var, unsigned int qp) const
Definition: fem_context.C:527
const std::vector< unsigned int > * active_vars() const
Return a pointer to the vector of active variables being computed for, or a null pointer if all varia...
Definition: fem_context.h:1034
Tensor point_hessian(unsigned int var, const Point &p) const
Definition: fem_context.C:967
Gradient fixed_interior_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1077
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
libmesh_assert(ctx)
void use_unweighted_quadrature_rules(int extra_quadrature_order=0)
Use quadrature rules designed to exactly integrate unweighted undistorted basis functions, plus extra_quadrature_order.
Definition: fem_context.C:177
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr)=0
This is at the core of this class.
void side_values(unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &side_values_vector) const
Fills a vector of values of the _system_vector at the all the quadrature points on the current elemen...
Definition: fem_context.C:679
const std::vector< dof_id_type > & get_dof_indices() const
Accessor for element dof indices.
Definition: diff_context.h:363
const BoundaryInfo & _boundary_info
Saved reference to BoundaryInfo on the mesh for this System.
Definition: fem_context.h:1187
void set_time(Real time_in)
Set the time for which the current nonlinear_solution is defined.
Definition: diff_context.h:423
std::vector< std::vector< FEAbstract * > > _element_fe_var
Pointers to the same finite element objects, but indexed by variable number.
Definition: fem_context.h:1179
const std::vector< std::vector< OutputTensor > > & get_d2phi() const
Definition: fe_base.h:319
std::map< const NumericVector< Number > *, std::pair< DenseVector< Number >, std::vector< DenseSubVector< Number > > > >::iterator localized_vectors_iterator
Typedef for the localized_vectors iterator.
Definition: diff_context.h:540
std::map< const NumericVector< Number > *, std::pair< DenseVector< Number >, std::vector< DenseSubVector< Number > > > > _localized_vectors
Contains pointers to vectors the user has asked to be localized, keyed with pairs of element localize...
Definition: diff_context.h:571
int _real_grad_fe_derivative_level
Definition: fem_context.h:1079
virtual void elem_side_reinit(Real theta) override
Resets the current time in the context.
Definition: fem_context.C:1432
DenseSubVector< Number > & get_localized_subvector(const NumericVector< Number > &localized_vector, unsigned int var)
Return a reference to DenseSubVector localization of localized_vector at variable var contained in th...
Definition: diff_context.C:154
Real get_system_time() const
Accessor for the time variable stored in the system class.
Definition: diff_context.h:411
Number fixed_point_value(unsigned int var, const Point &p) const
Definition: fem_context.C:1204
virtual void edge_fe_reinit()
Reinitializes edge FE objects on the current geometric element.
Definition: fem_context.C:1509
std::vector< DenseSubVector< Number > > _elem_subsolutions
Definition: diff_context.h:583
std::vector< FEAbstract * > _edge_fe_var
Definition: fem_context.h:1181
const DenseVector< Number > & get_elem_residual() const
Const accessor for element residual.
Definition: diff_context.h:242
This is a generic class that defines a solver to handle time integration of DifferentiableSystems.
std::vector< std::unique_ptr< QBase > > _side_qrule
Quadrature rules for element sides The FEM context will try to find a quadrature rule that correctly ...
Definition: fem_context.h:1224
std::vector< std::map< FEType, std::unique_ptr< FEAbstract > > > _side_fe
Definition: fem_context.h:1169
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< std::unique_ptr< QBase > > _element_qrule
Quadrature rule for element interior.
Definition: fem_context.h:1216
std::unique_ptr< QBase > unweighted_quadrature_rule(const unsigned int dim, const int extraorder=0) const
Definition: fe_type.C:53
unsigned char side
Current side for side_* to examine.
Definition: fem_context.h:1013
virtual unsigned short dim() const =0
const DenseVector< Number > & get_elem_solution_accel() const
Accessor for element solution accel of change w.r.t.
Definition: diff_context.h:177
static std::unique_ptr< FEAbstract > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
Definition: fe_abstract.C:77
void resize(const unsigned int new_m, const unsigned int new_n)
Resizes the matrix to the specified size and calls zero().
Definition: dense_matrix.h:895
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
std::unique_ptr< QBase > _edge_qrule
Quadrature rules for element edges.
Definition: fem_context.h:1233
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
Definition: system.h:1585
Helper nested class for C++03-compatible "template typedef".
Definition: fem_context.h:1110
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for interior finite element object for variable var for the largest dimension in the mesh...
Definition: fem_context.h:277
static std::unique_ptr< FEGenericBase > build_InfFE(const unsigned int dim, const FEType &type)
Builds a specific infinite element type.
Gradient fixed_side_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1154
unsigned int get_mesh_y_var() const
Accessor for y-variable of moving mesh System.
Definition: fem_context.h:874
AlgebraicType algebraic_type() const
Definition: fem_context.h:992
FEGenericBase< OutputShape > * build_new_fe(const FEGenericBase< OutputShape > *fe, const Point &p, const Real tolerance=TOLERANCE, const int get_derivative_level=-1) const
Helper function to reduce some code duplication in the *_point_* methods.
Definition: fem_context.C:2018
void interior_div(unsigned int var, unsigned int qp, OutputType &div_u) const
Definition: fem_context.C:624
virtual bool infinite() const =0
void some_gradient(unsigned int var, unsigned int qp, OutputType &u) const
Helper function to reduce some code duplication in the *interior_gradient methods.
Definition: fem_context.C:344
std::vector< std::vector< FEAbstract * > > _side_fe_var
Definition: fem_context.h:1180
unsigned int n_vars() const
Definition: system.h:2349
Gradient point_gradient(unsigned int var, const Point &p) const
Definition: fem_context.C:915
const Elem * _elem
Current element for element_* to examine.
Definition: fem_context.h:1192
void attach_quadrature_rules()
Helper function for attaching quadrature rules.
Definition: fem_context.C:127
void _do_elem_position_set(Real theta)
Uses the coordinate data specified by mesh_*_position configuration to set the geometry of elem to th...
Definition: fem_context.C:1602
virtual Order default_order() const =0
const std::vector< DenseVector< Number > > & get_qoi_derivatives() const
Const accessor for QoI derivatives.
Definition: diff_context.h:329
bool has_elem() const
Test for current Elem object.
Definition: fem_context.h:902
Gradient side_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:714
const DofMap & get_dof_map() const
Definition: system.h:2293
void interior_rate(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1354
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
const Point & point(const unsigned int i) const
Definition: elem.h:2277
void ErrorVector unsigned int
Definition: adjoints_ex3.C:360
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:111
void side_rate(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1375
This class forms the foundation from which generic finite elements may be derived.
virtual void elem_reinit(Real theta) override
Resets the current time in the context.
Definition: fem_context.C:1408
unsigned int get_mesh_z_var() const
Accessor for z-variable of moving mesh System.
Definition: fem_context.h:888
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
void old_dof_indices(const Elem &elem, unsigned int n, std::vector< dof_id_type > &di, const unsigned int vn) const
Appends to the vector di the old global degree of freedom indices for elem.node_ref(n), for one variable vn.
Definition: dof_map.C:2283
virtual_for_inffe const std::vector< std::vector< OutputShape > > & get_curl_phi() const
Definition: fe_base.h:252
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:207
TimeSolver & get_time_solver()
Definition: diff_system.h:454
static FEFamily map_fe_type(const Elem &elem)
Definition: fe_map.C:45
Tensor fixed_interior_hessian(unsigned int var, unsigned int qp) const
Definition: fem_context.C:1104