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