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