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