libMesh
fem_system.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #include "libmesh/dof_map.h"
21 #include "libmesh/elem.h"
22 #include "libmesh/equation_systems.h"
23 #include "libmesh/fe_base.h"
24 #include "libmesh/fem_context.h"
25 #include "libmesh/fem_system.h"
26 #include "libmesh/libmesh_logging.h"
27 #include "libmesh/mesh_base.h"
28 #include "libmesh/numeric_vector.h"
29 #include "libmesh/parallel_algebra.h"
30 #include "libmesh/parallel_ghost_sync.h"
31 #include "libmesh/quadrature.h"
32 #include "libmesh/sparse_matrix.h"
33 #include "libmesh/time_solver.h"
34 #include "libmesh/unsteady_solver.h" // For eulerian_residual
35 #include "libmesh/fe_interface.h"
36 
37 namespace {
38 using namespace libMesh;
39 
40 // give this guy some scope since there
41 // is underlying vector allocation upon
42 // creation/deletion
43 ConstElemRange elem_range;
44 
45 typedef Threads::spin_mutex femsystem_mutex;
46 femsystem_mutex assembly_mutex;
47 
48 void assemble_unconstrained_element_system(const FEMSystem & _sys,
49  const bool _get_jacobian,
50  const bool _constrain_heterogeneously,
51  FEMContext & _femcontext)
52 {
53  if (_sys.print_element_solutions)
54  {
55  std::streamsize old_precision = libMesh::out.precision();
57  if (_femcontext.has_elem())
58  libMesh::out << "U_elem " << _femcontext.get_elem().id();
59  else
60  libMesh::out << "U_scalar ";
61  libMesh::out << " = " << _femcontext.get_elem_solution() << std::endl;
62 
63  if (_sys.use_fixed_solution)
64  {
65  if (_femcontext.has_elem())
66  libMesh::out << "Ufixed_elem " << _femcontext.get_elem().id();
67  else
68  libMesh::out << "Ufixed_scalar ";
69  libMesh::out << " = " << _femcontext.get_elem_fixed_solution() << std::endl;
70  libMesh::out.precision(old_precision);
71  }
72  }
73 
74  // We need jacobians to do heterogeneous residual constraints
75  const bool need_jacobian =
76  (_get_jacobian || _constrain_heterogeneously);
77 
78  bool jacobian_computed =
79  _sys.time_solver->element_residual(need_jacobian, _femcontext);
80 
81  // Compute a numeric jacobian if we have to
82  if (need_jacobian && !jacobian_computed)
83  {
84  // Make sure we didn't compute a jacobian and lie about it
85  libmesh_assert_equal_to (_femcontext.get_elem_jacobian().l1_norm(), 0.0);
86  // Logging of numerical jacobians is done separately
87  _sys.numerical_elem_jacobian(_femcontext);
88  }
89 
90  // Compute a numeric jacobian if we're asked to verify the
91  // analytic jacobian we got
92  if (need_jacobian && jacobian_computed &&
93  _sys.verify_analytic_jacobians != 0.0)
94  {
95  DenseMatrix<Number> analytic_jacobian(_femcontext.get_elem_jacobian());
96 
97  _femcontext.get_elem_jacobian().zero();
98  // Logging of numerical jacobians is done separately
99  _sys.numerical_elem_jacobian(_femcontext);
100 
101  Real analytic_norm = analytic_jacobian.l1_norm();
102  Real numerical_norm = _femcontext.get_elem_jacobian().l1_norm();
103 
104  // If we can continue, we'll probably prefer the analytic jacobian
105  analytic_jacobian.swap(_femcontext.get_elem_jacobian());
106 
107  // The matrix "analytic_jacobian" will now hold the error matrix
108  analytic_jacobian.add(-1.0, _femcontext.get_elem_jacobian());
109  Real error_norm = analytic_jacobian.l1_norm();
110 
111  Real relative_error = error_norm /
112  std::max(analytic_norm, numerical_norm);
113 
114  if (relative_error > _sys.verify_analytic_jacobians)
115  {
116  libMesh::err << "Relative error " << relative_error
117  << " detected in analytic jacobian on element "
118  << _femcontext.get_elem().id() << '!' << std::endl;
119 
120  std::streamsize old_precision = libMesh::out.precision();
122  libMesh::out << "J_analytic " << _femcontext.get_elem().id() << " = "
123  << _femcontext.get_elem_jacobian() << std::endl;
124  analytic_jacobian.add(1.0, _femcontext.get_elem_jacobian());
125  libMesh::out << "J_numeric " << _femcontext.get_elem().id() << " = "
126  << analytic_jacobian << std::endl;
127 
128  libMesh::out.precision(old_precision);
129 
130  libmesh_error_msg("Relative error too large, exiting!");
131  }
132  }
133 
134  const unsigned char n_sides = _femcontext.get_elem().n_sides();
135  for (_femcontext.side = 0; _femcontext.side != n_sides;
136  ++_femcontext.side)
137  {
138  // Don't compute on non-boundary sides unless requested
139  if (!_sys.get_physics()->compute_internal_sides &&
140  _femcontext.get_elem().neighbor_ptr(_femcontext.side) != nullptr)
141  continue;
142 
143  // Any mesh movement has already been done (and restored,
144  // if the TimeSolver isn't broken), but
145  // reinitializing the side FE objects is still necessary
146  _femcontext.side_fe_reinit();
147 
148  DenseMatrix<Number> old_jacobian;
149  // If we're in DEBUG mode, we should always verify that the
150  // user's side_residual function doesn't alter our existing
151  // jacobian and then lie about it
152 #ifndef DEBUG
153  // Even if we're not in DEBUG mode, when we're verifying
154  // analytic jacobians we'll want to verify each side's
155  // jacobian contribution separately.
156  if (_sys.verify_analytic_jacobians != 0.0 && need_jacobian)
157 #endif // ifndef DEBUG
158  {
159  old_jacobian = _femcontext.get_elem_jacobian();
160  _femcontext.get_elem_jacobian().zero();
161  }
162 
163  jacobian_computed =
164  _sys.time_solver->side_residual(need_jacobian, _femcontext);
165 
166  // Compute a numeric jacobian if we have to
167  if (need_jacobian && !jacobian_computed)
168  {
169  // If we have already backed up old_jacobian,
170  // we can make sure side_residual didn't compute a
171  // jacobian and lie about it.
172  //
173  // If we haven't, then we need to, to let
174  // numerical_side_jacobian work.
175  if (old_jacobian.m())
176  libmesh_assert_equal_to (_femcontext.get_elem_jacobian().l1_norm(), 0.0);
177  else
178  {
179  old_jacobian = _femcontext.get_elem_jacobian();
180  _femcontext.get_elem_jacobian().zero();
181  }
182 
183  // Logging of numerical jacobians is done separately
184  _sys.numerical_side_jacobian(_femcontext);
185 
186  // Add back in element interior numerical Jacobian
187  _femcontext.get_elem_jacobian() += old_jacobian;
188  }
189 
190  // Compute a numeric jacobian if we're asked to verify the
191  // analytic jacobian we got
192  else if (need_jacobian && jacobian_computed &&
193  _sys.verify_analytic_jacobians != 0.0)
194  {
195  DenseMatrix<Number> analytic_jacobian(_femcontext.get_elem_jacobian());
196 
197  _femcontext.get_elem_jacobian().zero();
198  // Logging of numerical jacobians is done separately
199  _sys.numerical_side_jacobian(_femcontext);
200 
201  Real analytic_norm = analytic_jacobian.l1_norm();
202  Real numerical_norm = _femcontext.get_elem_jacobian().l1_norm();
203 
204  // If we can continue, we'll probably prefer the analytic jacobian
205  analytic_jacobian.swap(_femcontext.get_elem_jacobian());
206 
207  // The matrix "analytic_jacobian" will now hold the error matrix
208  analytic_jacobian.add(-1.0, _femcontext.get_elem_jacobian());
209  Real error_norm = analytic_jacobian.l1_norm();
210 
211  Real relative_error = error_norm /
212  std::max(analytic_norm, numerical_norm);
213 
214  if (relative_error > _sys.verify_analytic_jacobians)
215  {
216  libMesh::err << "Relative error " << relative_error
217  << " detected in analytic jacobian on element "
218  << _femcontext.get_elem().id()
219  << ", side "
220  << static_cast<unsigned int>(_femcontext.side) << '!' << std::endl;
221 
222  std::streamsize old_precision = libMesh::out.precision();
224  libMesh::out << "J_analytic " << _femcontext.get_elem().id() << " = "
225  << _femcontext.get_elem_jacobian() << std::endl;
226  analytic_jacobian.add(1.0, _femcontext.get_elem_jacobian());
227  libMesh::out << "J_numeric " << _femcontext.get_elem().id() << " = "
228  << analytic_jacobian << std::endl;
229  libMesh::out.precision(old_precision);
230 
231  libmesh_error_msg("Relative error too large, exiting!");
232  }
233  // Once we've verified a side, we'll want to add back the
234  // rest of the accumulated jacobian
235  _femcontext.get_elem_jacobian() += old_jacobian;
236  }
237 
238  // In DEBUG mode, we've set elem_jacobian == 0, and we
239  // may have yet to add the old jacobian back
240 #ifdef DEBUG
241  else
242  {
243  _femcontext.get_elem_jacobian() += old_jacobian;
244  }
245 #endif // ifdef DEBUG
246  }
247 }
248 
249 void add_element_system(const FEMSystem & _sys,
250  const bool _get_residual,
251  const bool _get_jacobian,
252  const bool _constrain_heterogeneously,
253  const bool _no_constraints,
254  FEMContext & _femcontext)
255 {
256 #ifdef LIBMESH_ENABLE_CONSTRAINTS
257  if (_get_residual && _sys.print_element_residuals)
258  {
259  std::streamsize old_precision = libMesh::out.precision();
261  if (_femcontext.has_elem())
262  libMesh::out << "Rraw_elem " << _femcontext.get_elem().id();
263  else
264  libMesh::out << "Rraw_scalar ";
265  libMesh::out << " = " << _femcontext.get_elem_residual() << std::endl;
266  libMesh::out.precision(old_precision);
267  }
268 
269  // We turn off the asymmetric constraint application;
270  // enforce_constraints_exactly() should be called in the solver
271  if (_get_residual && _get_jacobian)
272  {
273  if (_constrain_heterogeneously)
275  (_femcontext.get_elem_jacobian(),
276  _femcontext.get_elem_residual(),
277  _femcontext.get_dof_indices(), false);
278  else if (!_no_constraints)
280  (_femcontext.get_elem_jacobian(),
281  _femcontext.get_elem_residual(),
282  _femcontext.get_dof_indices(), false);
283  // Do nothing if (_no_constraints)
284  }
285  else if (_get_residual)
286  {
287  if (_constrain_heterogeneously)
289  (_femcontext.get_elem_jacobian(),
290  _femcontext.get_elem_residual(),
291  _femcontext.get_dof_indices(), false);
292  else if (!_no_constraints)
294  (_femcontext.get_elem_residual(), _femcontext.get_dof_indices(), false);
295  // Do nothing if (_no_constraints)
296  }
297  else if (_get_jacobian)
298  {
299  // Heterogeneous and homogeneous constraints are the same on the
300  // matrix
301  // Only get these contribs if we are applying some constraints
302  if (!_no_constraints)
304  _femcontext.get_dof_indices(),
305  false);
306  }
307 #else
308  libmesh_ignore(_constrain_heterogeneously, _no_constraints);
309 #endif // #ifdef LIBMESH_ENABLE_CONSTRAINTS
310 
311  if (_get_residual && _sys.print_element_residuals)
312  {
313  std::streamsize old_precision = libMesh::out.precision();
315  if (_femcontext.has_elem())
316  libMesh::out << "R_elem " << _femcontext.get_elem().id();
317  else
318  libMesh::out << "R_scalar ";
319  libMesh::out << " = " << _femcontext.get_elem_residual() << std::endl;
320  libMesh::out.precision(old_precision);
321  }
322 
323  if (_get_jacobian && _sys.print_element_jacobians)
324  {
325  std::streamsize old_precision = libMesh::out.precision();
327  if (_femcontext.has_elem())
328  libMesh::out << "J_elem " << _femcontext.get_elem().id();
329  else
330  libMesh::out << "J_scalar ";
331  libMesh::out << " = " << _femcontext.get_elem_jacobian() << std::endl;
332  libMesh::out.precision(old_precision);
333  }
334 
335  { // A lock is necessary around access to the global system
336  femsystem_mutex::scoped_lock lock(assembly_mutex);
337 
338  if (_get_jacobian)
339  _sys.matrix->add_matrix (_femcontext.get_elem_jacobian(),
340  _femcontext.get_dof_indices());
341  if (_get_residual)
342  _sys.rhs->add_vector (_femcontext.get_elem_residual(),
343  _femcontext.get_dof_indices());
344  } // Scope for assembly mutex
345 }
346 
347 
348 
349 class AssemblyContributions
350 {
351 public:
355  AssemblyContributions(FEMSystem & sys,
356  bool get_residual,
357  bool get_jacobian,
358  bool constrain_heterogeneously,
359  bool no_constraints) :
360  _sys(sys),
361  _get_residual(get_residual),
362  _get_jacobian(get_jacobian),
363  _constrain_heterogeneously(constrain_heterogeneously),
364  _no_constraints(no_constraints) {}
365 
369  void operator()(const ConstElemRange & range) const
370  {
371  std::unique_ptr<DiffContext> con = _sys.build_context();
372  FEMContext & _femcontext = cast_ref<FEMContext &>(*con);
373  _sys.init_context(_femcontext);
374 
375  for (const auto & elem : range)
376  {
377  _femcontext.pre_fe_reinit(_sys, elem);
378  _femcontext.elem_fe_reinit();
379 
380  assemble_unconstrained_element_system
381  (_sys, _get_jacobian, _constrain_heterogeneously, _femcontext);
382 
383  add_element_system
384  (_sys, _get_residual, _get_jacobian,
385  _constrain_heterogeneously, _no_constraints, _femcontext);
386  }
387  }
388 
389 private:
390 
391  FEMSystem & _sys;
392 
393  const bool _get_residual, _get_jacobian, _constrain_heterogeneously, _no_constraints;
394 };
395 
396 class PostprocessContributions
397 {
398 public:
402  explicit
403  PostprocessContributions(FEMSystem & sys) : _sys(sys) {}
404 
408  void operator()(const ConstElemRange & range) const
409  {
410  std::unique_ptr<DiffContext> con = _sys.build_context();
411  FEMContext & _femcontext = cast_ref<FEMContext &>(*con);
412  _sys.init_context(_femcontext);
413 
414  for (const auto & elem : range)
415  {
416  _femcontext.pre_fe_reinit(_sys, elem);
417 
418  // Optionally initialize all the interior FE objects on elem.
420  _femcontext.elem_fe_reinit();
421 
422  _sys.element_postprocess(_femcontext);
423 
424  const unsigned char n_sides = _femcontext.get_elem().n_sides();
425  for (_femcontext.side = 0; _femcontext.side != n_sides;
426  ++_femcontext.side)
427  {
428  // Don't compute on non-boundary sides unless requested
429  if (!_sys.postprocess_sides ||
431  _femcontext.get_elem().neighbor_ptr(_femcontext.side) != nullptr))
432  continue;
433 
434  // Optionally initialize all the FE objects on this side.
436  _femcontext.side_fe_reinit();
437 
438  _sys.side_postprocess(_femcontext);
439  }
440  }
441  }
442 
443 private:
444 
445  FEMSystem & _sys;
446 };
447 
448 class QoIContributions
449 {
450 public:
454  explicit
455  QoIContributions(FEMSystem & sys,
456  DifferentiableQoI & diff_qoi,
457  const QoISet & qoi_indices) :
458  qoi(sys.n_qois(), 0.), _sys(sys), _diff_qoi(diff_qoi),_qoi_indices(qoi_indices) {}
459 
463  QoIContributions(const QoIContributions & other,
464  Threads::split) :
465  qoi(other._sys.n_qois(), 0.), _sys(other._sys), _diff_qoi(other._diff_qoi) {}
466 
470  void operator()(const ConstElemRange & range)
471  {
472  std::unique_ptr<DiffContext> con = _sys.build_context();
473  FEMContext & _femcontext = cast_ref<FEMContext &>(*con);
474  _diff_qoi.init_context(_femcontext);
475 
476  bool have_some_heterogenous_qoi_bc = false;
477 #ifdef LIBMESH_ENABLE_CONSTRAINTS
478  std::vector<bool> have_heterogenous_qoi_bc(_sys.n_qois(), false);
479  for (auto q : IntRange<unsigned int>(0, _sys.n_qois()))
480  if (_qoi_indices.has_index(q) &&
482  {
483  have_heterogenous_qoi_bc[q] = true;
484  have_some_heterogenous_qoi_bc = true;
485  }
486 #endif
487 
488  if (have_some_heterogenous_qoi_bc)
489  _sys.init_context(_femcontext);
490 
491  for (const auto & elem : range)
492  {
493  _femcontext.pre_fe_reinit(_sys, elem);
494 
495  // We might have some heterogenous dofs here; let's see for
496  // certain
497  bool elem_has_some_heterogenous_qoi_bc = false;
498 
499 #ifdef LIBMESH_ENABLE_CONSTRAINTS
500  const unsigned int n_dofs =
501  cast_int<unsigned int>(_femcontext.get_dof_indices().size());
502 
503  std::vector<bool> elem_has_heterogenous_qoi_bc(_sys.n_qois(), false);
504  if (have_some_heterogenous_qoi_bc)
505  {
506  for (auto q : IntRange<unsigned int>(0, _sys.n_qois()))
507  {
508  if (have_heterogenous_qoi_bc[q])
509  {
510  for (auto d : IntRange<unsigned int>(0, n_dofs))
512  (q, _femcontext.get_dof_indices()[d]) != Number(0))
513  {
514  elem_has_some_heterogenous_qoi_bc = true;
515  elem_has_heterogenous_qoi_bc[q] = true;
516  break;
517  }
518  }
519  }
520  }
521 #endif
522 
523  if (_diff_qoi.assemble_qoi_elements ||
524  elem_has_some_heterogenous_qoi_bc)
525  _femcontext.elem_fe_reinit();
526 
527  if (_diff_qoi.assemble_qoi_elements)
528  _diff_qoi.element_qoi(_femcontext, _qoi_indices);
529 
530  // If we have some heterogenous dofs here, those are
531  // themselves part of a regularized flux QoI which the library
532  // handles integrating
533 #ifdef LIBMESH_ENABLE_CONSTRAINTS
534  if (elem_has_some_heterogenous_qoi_bc)
535  {
536  _sys.time_solver->element_residual(false, _femcontext);
537 
538  for (auto q : IntRange<unsigned int>(0, _sys.n_qois()))
539  {
540  if (elem_has_heterogenous_qoi_bc[q])
541  {
542  for (auto d : IntRange<unsigned int>(0, n_dofs))
543  this->qoi[q] -= _femcontext.get_elem_residual()(d) *
545 
546  }
547  }
548  }
549 #endif
550 
551  const unsigned char n_sides = _femcontext.get_elem().n_sides();
552  for (_femcontext.side = 0; _femcontext.side != n_sides;
553  ++_femcontext.side)
554  {
555  // Don't compute on non-boundary sides unless requested
556  if (!_diff_qoi.assemble_qoi_sides ||
557  (!_diff_qoi.assemble_qoi_internal_sides &&
558  _femcontext.get_elem().neighbor_ptr(_femcontext.side) != nullptr))
559  continue;
560 
561  _femcontext.side_fe_reinit();
562 
563  _diff_qoi.side_qoi(_femcontext, _qoi_indices);
564  }
565  }
566 
567  this->_diff_qoi.thread_join( this->qoi, _femcontext.get_qois(), _qoi_indices );
568  }
569 
570  void join (const QoIContributions & other)
571  {
572  libmesh_assert_equal_to (this->qoi.size(), other.qoi.size());
573  this->_diff_qoi.thread_join( this->qoi, other.qoi, _qoi_indices );
574  }
575 
576  std::vector<Number> qoi;
577 
578 private:
579 
580  FEMSystem & _sys;
581  DifferentiableQoI & _diff_qoi;
582 
583  const QoISet _qoi_indices;
584 };
585 
586 class QoIDerivativeContributions
587 {
588 public:
592  QoIDerivativeContributions(FEMSystem & sys,
593  const QoISet & qoi_indices,
594  DifferentiableQoI & qoi,
595  bool include_liftfunc,
596  bool apply_constraints) :
597  _sys(sys),
598  _qoi_indices(qoi_indices),
599  _qoi(qoi),
600  _include_liftfunc(include_liftfunc),
601  _apply_constraints(apply_constraints) {}
602 
606  void operator()(const ConstElemRange & range) const
607  {
608  std::unique_ptr<DiffContext> con = _sys.build_context();
609  FEMContext & _femcontext = cast_ref<FEMContext &>(*con);
610  _qoi.init_context(_femcontext);
611 
612  bool have_some_heterogenous_qoi_bc = false;
613 #ifdef LIBMESH_ENABLE_CONSTRAINTS
614  std::vector<bool> have_heterogenous_qoi_bc(_sys.n_qois(), false);
615  if (_include_liftfunc || _apply_constraints)
616  for (auto q : IntRange<unsigned int>(0, _sys.n_qois()))
617  if (_qoi_indices.has_index(q) &&
619  {
620  have_heterogenous_qoi_bc[q] = true;
621  have_some_heterogenous_qoi_bc = true;
622  }
623 #endif
624 
625  if (have_some_heterogenous_qoi_bc)
626  _sys.init_context(_femcontext);
627 
628  for (const auto & elem : range)
629  {
630  _femcontext.pre_fe_reinit(_sys, elem);
631 
632  // We might have some heterogenous dofs here; let's see for
633  // certain
634  bool elem_has_some_heterogenous_qoi_bc = false;
635 
636 #ifdef LIBMESH_ENABLE_CONSTRAINTS
637  const unsigned int n_dofs =
638  cast_int<unsigned int>(_femcontext.get_dof_indices().size());
639 
640  std::vector<bool> elem_has_heterogenous_qoi_bc(_sys.n_qois(), false);
641  if (have_some_heterogenous_qoi_bc)
642  {
643  for (auto q : IntRange<unsigned int>(0, _sys.n_qois()))
644  {
645  if (have_heterogenous_qoi_bc[q])
646  {
647  for (auto d : IntRange<unsigned int>(0, n_dofs))
649  (q, _femcontext.get_dof_indices()[d]) != Number(0))
650  {
651  elem_has_some_heterogenous_qoi_bc = true;
652  elem_has_heterogenous_qoi_bc[q] = true;
653  break;
654  }
655  }
656  }
657  }
658 #endif
659 
660  // If we're going to call a user integral, then we need FE
661  // information to call element_qoi.
662  // If we're going to evaluate lift-function-based components
663  // of a QoI, then we need FE information to assemble the
664  // element residual.
665  if (_qoi.assemble_qoi_elements ||
666  ((_include_liftfunc || _apply_constraints) &&
667  elem_has_some_heterogenous_qoi_bc))
668  _femcontext.elem_fe_reinit();
669 
670  if (_qoi.assemble_qoi_elements)
671  _qoi.element_qoi_derivative(_femcontext, _qoi_indices);
672 
673 #ifdef LIBMESH_ENABLE_CONSTRAINTS
674  // If we need to use heterogenous dofs here, we need the
675  // Jacobian either for the regularized flux QoI integration
676  // and/or for constraint application.
677  if ((_include_liftfunc || _apply_constraints) &&
678  elem_has_some_heterogenous_qoi_bc)
679  {
680  bool jacobian_computed = _sys.time_solver->element_residual(true, _femcontext);
681 
682  // If we're using numerical jacobians, above wont compute them
683  if (!jacobian_computed)
684  {
685  // Make sure we didn't compute a jacobian and lie about it
686  libmesh_assert_equal_to (_femcontext.get_elem_jacobian().l1_norm(), 0.0);
687  // Logging of numerical jacobians is done separately
688  _sys.numerical_elem_jacobian(_femcontext);
689  }
690  }
691 
692  // If we have some heterogenous dofs here, those are
693  // themselves part of a regularized flux QoI which the library
694  // may handle integrating
695  if (_include_liftfunc && elem_has_some_heterogenous_qoi_bc)
696  {
697  for (auto q : IntRange<unsigned int>(0, _sys.n_qois()))
698  {
699  if (elem_has_heterogenous_qoi_bc[q])
700  {
701  for (auto i : IntRange<unsigned int>(0, n_dofs))
702  {
703  Number liftfunc_val =
705 
706  if (liftfunc_val != Number(0))
707  {
708  for (auto j : IntRange<unsigned int>(0, n_dofs))
709  _femcontext.get_qoi_derivatives()[q](j) -=
710  _femcontext.get_elem_jacobian()(i,j) *
711  liftfunc_val;
712  }
713  }
714  }
715  }
716  }
717 #endif
718 
719 
720  const unsigned char n_sides = _femcontext.get_elem().n_sides();
721  for (_femcontext.side = 0; _femcontext.side != n_sides;
722  ++_femcontext.side)
723  {
724  // Don't compute on non-boundary sides unless requested
725  if (!_qoi.assemble_qoi_sides ||
726  (!_qoi.assemble_qoi_internal_sides &&
727  _femcontext.get_elem().neighbor_ptr(_femcontext.side) != nullptr))
728  continue;
729 
730  _femcontext.side_fe_reinit();
731 
732  _qoi.side_qoi_derivative(_femcontext, _qoi_indices);
733  }
734 
735  // We need some unmodified indices to use for constraining
736  // multiple vector
737  // FIXME - there should be a DofMap::constrain_element_vectors
738  // to do this more efficiently
739 #ifdef LIBMESH_ENABLE_CONSTRAINTS
740  std::vector<dof_id_type> original_dofs = _femcontext.get_dof_indices();
741 #endif
742 
743  { // A lock is necessary around access to the global system
744  femsystem_mutex::scoped_lock lock(assembly_mutex);
745 
746 #ifdef LIBMESH_ENABLE_CONSTRAINTS
747  // We'll need to see if any heterogenous constraints apply
748  // to the QoI dofs on this element *or* to any of the dofs
749  // they depend on, so let's get those dependencies
750  if (_apply_constraints)
751  _sys.get_dof_map().constrain_nothing(_femcontext.get_dof_indices());
752 #endif
753 
754  for (auto i : IntRange<unsigned int>(0, _sys.n_qois()))
755  if (_qoi_indices.has_index(i))
756  {
757 #ifdef LIBMESH_ENABLE_CONSTRAINTS
758  if (_apply_constraints)
759  {
760 #ifndef NDEBUG
761  bool has_heterogenous_constraint = false;
762  for (auto d : IntRange<unsigned int>(0, n_dofs))
764  (i, _femcontext.get_dof_indices()[d]) != Number(0))
765  {
766  has_heterogenous_constraint = true;
767  libmesh_assert(elem_has_heterogenous_qoi_bc[i]);
768  libmesh_assert(elem_has_some_heterogenous_qoi_bc);
769  break;
770  }
771 #else
772  bool has_heterogenous_constraint =
773  elem_has_heterogenous_qoi_bc[i];
774 #endif
775 
776  _femcontext.get_dof_indices() = original_dofs;
777 
778  if (has_heterogenous_constraint)
779  {
780  // Q_u gets used for *adjoint* solves, so we
781  // need K^T here.
782  DenseMatrix<Number> elem_jacobian_transpose;
783  _femcontext.get_elem_jacobian().get_transpose
784  (elem_jacobian_transpose);
785 
787  (elem_jacobian_transpose,
788  _femcontext.get_qoi_derivatives()[i],
789  _femcontext.get_dof_indices(), false, i);
790  }
791  else
792  {
794  (_femcontext.get_qoi_derivatives()[i],
795  _femcontext.get_dof_indices(), false);
796  }
797  }
798 #endif
799 
801  (_femcontext.get_qoi_derivatives()[i], _femcontext.get_dof_indices());
802  }
803  }
804  }
805  }
806 
807 private:
808 
809  FEMSystem & _sys;
810  const QoISet & _qoi_indices;
811  DifferentiableQoI & _qoi;
812  bool _include_liftfunc, _apply_constraints;
813 };
814 
815 
816 }
817 
818 
819 namespace libMesh
820 {
821 
822 
823 
824 
825 
827  const std::string & name_in,
828  const unsigned int number_in)
829  : Parent(es, name_in, number_in),
830  fe_reinit_during_postprocess(true),
831  numerical_jacobian_h(TOLERANCE),
832  verify_analytic_jacobians(0.0)
833 {
834 }
835 
836 
838 {
839 }
840 
841 
842 
844 {
845  // First initialize LinearImplicitSystem data
847 }
848 
849 
850 void FEMSystem::assembly (bool get_residual, bool get_jacobian,
851  bool apply_heterogeneous_constraints,
852  bool apply_no_constraints)
853 {
854  libmesh_assert(get_residual || get_jacobian);
855 
856  // Log residual and jacobian and combined performance separately
857 #ifdef LIBMESH_ENABLE_PERFORMANCE_LOGGING
858  const char * log_name;
859  if (get_residual && get_jacobian)
860  log_name = "assembly()";
861  else if (get_residual)
862  log_name = "assembly(get_residual)";
863  else
864  log_name = "assembly(get_jacobian)";
865 
866  LOG_SCOPE(log_name, "FEMSystem");
867 #endif
868 
869  const MeshBase & mesh = this->get_mesh();
870 
872  {
873  this->solution->close();
874 
875  std::streamsize old_precision = libMesh::out.precision();
877  libMesh::out << "|U| = "
878  << this->solution->l1_norm()
879  << std::endl;
880  libMesh::out.precision(old_precision);
881  }
882  if (print_solutions)
883  {
884  std::streamsize old_precision = libMesh::out.precision();
886  libMesh::out << "U = [" << *(this->solution)
887  << "];" << std::endl;
888  libMesh::out.precision(old_precision);
889  }
890 
891  // Is this definitely necessary? [RHS]
892  // Yes. [RHS 2012]
893  if (get_jacobian)
894  matrix->zero();
895  if (get_residual)
896  rhs->zero();
897 
898  // Stupid C++ lets you set *Real* verify_analytic_jacobians = true!
899  if (verify_analytic_jacobians > 0.5)
900  {
901  libMesh::err << "WARNING! verify_analytic_jacobians was set "
902  << "to absurdly large value of "
903  << verify_analytic_jacobians << std::endl;
904  libMesh::err << "Resetting to 1e-6!" << std::endl;
906  }
907 
908  // In time-dependent problems, the nonlinear function we're trying
909  // to solve at each timestep may depend on the particular solver
910  // we're using
912 
913  // Build the residual and jacobian contributions on every active
914  // mesh element on this processor
916  (elem_range.reset(mesh.active_local_elements_begin(),
918  AssemblyContributions(*this, get_residual, get_jacobian,
919  apply_heterogeneous_constraints,
920  apply_no_constraints));
921 
922  // Check and see if we have SCALAR variables
923  bool have_scalar = false;
924  for (auto i : IntRange<unsigned int>(0, this->n_variable_groups()))
925  {
926  if (this->variable_group(i).type().family == SCALAR)
927  {
928  have_scalar = true;
929  break;
930  }
931  }
932 
933  // SCALAR dofs are stored on the last processor, so we'll evaluate
934  // their equation terms there and only if we have a SCALAR variable
935  if (this->processor_id() == (this->n_processors()-1) && have_scalar)
936  {
937  std::unique_ptr<DiffContext> con = this->build_context();
938  FEMContext & _femcontext = cast_ref<FEMContext &>(*con);
939  this->init_context(_femcontext);
940  _femcontext.pre_fe_reinit(*this, nullptr);
941 
942  bool jacobian_computed =
943  this->time_solver->nonlocal_residual(get_jacobian, _femcontext);
944 
945  // Nonlocal residuals are likely to be length 0, in which case we
946  // don't need to do any more. And we shouldn't try to do any
947  // more; lots of DenseVector/DenseMatrix code assumes rank>0.
948  if (_femcontext.get_elem_residual().size())
949  {
950  // Compute a numeric jacobian if we have to
951  if (get_jacobian && !jacobian_computed)
952  {
953  // Make sure we didn't compute a jacobian and lie about it
954  libmesh_assert_equal_to (_femcontext.get_elem_jacobian().l1_norm(), 0.0);
955  // Logging of numerical jacobians is done separately
956  this->numerical_nonlocal_jacobian(_femcontext);
957  }
958 
959  // Compute a numeric jacobian if we're asked to verify the
960  // analytic jacobian we got
961  if (get_jacobian && jacobian_computed &&
962  this->verify_analytic_jacobians != 0.0)
963  {
964  DenseMatrix<Number> analytic_jacobian(_femcontext.get_elem_jacobian());
965 
966  _femcontext.get_elem_jacobian().zero();
967  // Logging of numerical jacobians is done separately
968  this->numerical_nonlocal_jacobian(_femcontext);
969 
970  Real analytic_norm = analytic_jacobian.l1_norm();
971  Real numerical_norm = _femcontext.get_elem_jacobian().l1_norm();
972 
973  // If we can continue, we'll probably prefer the analytic jacobian
974  analytic_jacobian.swap(_femcontext.get_elem_jacobian());
975 
976  // The matrix "analytic_jacobian" will now hold the error matrix
977  analytic_jacobian.add(-1.0, _femcontext.get_elem_jacobian());
978  Real error_norm = analytic_jacobian.l1_norm();
979 
980  Real relative_error = error_norm /
981  std::max(analytic_norm, numerical_norm);
982 
983  if (relative_error > this->verify_analytic_jacobians)
984  {
985  libMesh::err << "Relative error " << relative_error
986  << " detected in analytic jacobian on nonlocal dofs!"
987  << std::endl;
988 
989  std::streamsize old_precision = libMesh::out.precision();
991  libMesh::out << "J_analytic nonlocal = "
992  << _femcontext.get_elem_jacobian() << std::endl;
993  analytic_jacobian.add(1.0, _femcontext.get_elem_jacobian());
994  libMesh::out << "J_numeric nonlocal = "
995  << analytic_jacobian << std::endl;
996 
997  libMesh::out.precision(old_precision);
998 
999  libmesh_error_msg("Relative error too large, exiting!");
1000  }
1001  }
1002 
1003  add_element_system
1004  (*this, get_residual, get_jacobian,
1005  apply_heterogeneous_constraints, apply_no_constraints, _femcontext);
1006  }
1007  }
1008 
1009  if (get_residual && (print_residual_norms || print_residuals))
1010  this->rhs->close();
1011  if (get_residual && print_residual_norms)
1012  {
1013  std::streamsize old_precision = libMesh::out.precision();
1014  libMesh::out.precision(16);
1015  libMesh::out << "|F| = " << this->rhs->l1_norm() << std::endl;
1016  libMesh::out.precision(old_precision);
1017  }
1018  if (get_residual && print_residuals)
1019  {
1020  std::streamsize old_precision = libMesh::out.precision();
1021  libMesh::out.precision(16);
1022  libMesh::out << "F = [" << *(this->rhs) << "];" << std::endl;
1023  libMesh::out.precision(old_precision);
1024  }
1025 
1026  if (get_jacobian && (print_jacobian_norms || print_jacobians))
1027  this->matrix->close();
1028  if (get_jacobian && print_jacobian_norms)
1029  {
1030  std::streamsize old_precision = libMesh::out.precision();
1031  libMesh::out.precision(16);
1032  libMesh::out << "|J| = " << this->matrix->l1_norm() << std::endl;
1033  libMesh::out.precision(old_precision);
1034  }
1035  if (get_jacobian && print_jacobians)
1036  {
1037  std::streamsize old_precision = libMesh::out.precision();
1038  libMesh::out.precision(16);
1039  libMesh::out << "J = [" << *(this->matrix) << "];" << std::endl;
1040  libMesh::out.precision(old_precision);
1041  }
1042 }
1043 
1044 
1045 
1047 {
1048  // We are solving the primal problem
1049  Parent::solve();
1050 
1051  // On a moving mesh we want the mesh to reflect the new solution
1052  this->mesh_position_set();
1053 }
1054 
1055 
1056 
1058 {
1059  // If we don't need to move the mesh, we're done
1060  if (_mesh_sys != this)
1061  return;
1062 
1063  MeshBase & mesh = this->get_mesh();
1064 
1065  std::unique_ptr<DiffContext> con = this->build_context();
1066  FEMContext & _femcontext = cast_ref<FEMContext &>(*con);
1067  this->init_context(_femcontext);
1068 
1069  // Move every mesh element we can
1070  for (const auto & elem : mesh.active_local_element_ptr_range())
1071  {
1072  // We need the algebraic data
1073  _femcontext.pre_fe_reinit(*this, elem);
1074  // And when asserts are on, we also need the FE so
1075  // we can assert that the mesh data is of the right type.
1076 #ifndef NDEBUG
1077  _femcontext.elem_fe_reinit();
1078 #endif
1079 
1080  // This code won't handle moving subactive elements
1081  libmesh_assert(!_femcontext.get_elem().has_children());
1082 
1083  _femcontext.elem_position_set(1.);
1084  }
1085 
1086  // We've now got positions set on all local nodes (and some
1087  // semilocal nodes); let's request positions for non-local nodes
1088  // from their processors.
1089 
1090  SyncNodalPositions sync_object(mesh);
1092  (this->comm(), mesh.nodes_begin(), mesh.nodes_end(), sync_object);
1093 }
1094 
1095 
1096 
1098 {
1099  LOG_SCOPE("postprocess()", "FEMSystem");
1100 
1101  const MeshBase & mesh = this->get_mesh();
1102 
1103  this->update();
1104 
1105  // Get the time solver object associated with the system, and tell it that
1106  // we are not solving the adjoint problem
1107  this->get_time_solver().set_is_adjoint(false);
1108 
1109  // Loop over every active mesh element on this processor
1112  PostprocessContributions(*this));
1113 }
1114 
1115 
1116 
1117 void FEMSystem::assemble_qoi (const QoISet & qoi_indices)
1118 {
1119  LOG_SCOPE("assemble_qoi()", "FEMSystem");
1120 
1121  const MeshBase & mesh = this->get_mesh();
1122 
1123  this->update();
1124 
1125  const unsigned int Nq = this->n_qois();
1126 
1127  // the quantity of interest is assumed to be a sum of element and
1128  // side terms
1129  for (unsigned int i=0; i != Nq; ++i)
1130  if (qoi_indices.has_index(i))
1131  qoi[i] = 0;
1132 
1133  // Create a non-temporary qoi_contributions object, so we can query
1134  // its results after the reduction
1135  QoIContributions qoi_contributions(*this, *(this->diff_qoi), qoi_indices);
1136 
1137  // Loop over every active mesh element on this processor
1140  qoi_contributions);
1141 
1142  this->diff_qoi->parallel_op( this->comm(), this->qoi, qoi_contributions.qoi, qoi_indices );
1143 }
1144 
1145 
1146 
1147 void FEMSystem::assemble_qoi_derivative (const QoISet & qoi_indices,
1148  bool include_liftfunc,
1149  bool apply_constraints)
1150 {
1151  LOG_SCOPE("assemble_qoi_derivative()", "FEMSystem");
1152 
1153  const MeshBase & mesh = this->get_mesh();
1154 
1155  this->update();
1156 
1157  // The quantity of interest derivative assembly accumulates on
1158  // initially zero vectors
1159  for (auto i : IntRange<unsigned int>(0, this->n_qois()))
1160  if (qoi_indices.has_index(i))
1161  this->add_adjoint_rhs(i).zero();
1162 
1163  // Loop over every active mesh element on this processor
1166  QoIDerivativeContributions(*this, qoi_indices,
1167  *(this->diff_qoi),
1168  include_liftfunc,
1169  apply_constraints));
1170 
1171  for (auto i : IntRange<unsigned int>(0, this->n_qois()))
1172  if (qoi_indices.has_index(i))
1173  this->diff_qoi->finalize_derivative(this->get_adjoint_rhs(i),i);
1174 }
1175 
1176 
1177 
1178 void FEMSystem::numerical_jacobian (TimeSolverResPtr res,
1179  FEMContext & context) const
1180 {
1181  // Logging is done by numerical_elem_jacobian
1182  // or numerical_side_jacobian
1183 
1184  DenseVector<Number> original_residual(context.get_elem_residual());
1185  DenseVector<Number> backwards_residual(context.get_elem_residual());
1186  DenseMatrix<Number> numeric_jacobian(context.get_elem_jacobian());
1187 #ifdef DEBUG
1188  DenseMatrix<Number> old_jacobian(context.get_elem_jacobian());
1189 #endif
1190 
1191  Real numerical_point_h = 0.;
1192  if (_mesh_sys == this)
1193  numerical_point_h = numerical_jacobian_h * context.get_elem().hmin();
1194 
1195  const unsigned int n_dofs =
1196  cast_int<unsigned int>(context.get_dof_indices().size());
1197 
1198  for (auto v : IntRange<unsigned int>(0, context.n_vars()))
1199  {
1200  const Real my_h = this->numerical_jacobian_h_for_var(v);
1201 
1202  unsigned int j_offset = libMesh::invalid_uint;
1203 
1204  if (!context.get_dof_indices(v).empty())
1205  {
1206  for (auto i : IntRange<unsigned int>(0, n_dofs))
1207  if (context.get_dof_indices()[i] ==
1208  context.get_dof_indices(v)[0])
1209  j_offset = i;
1210 
1211  libmesh_assert_not_equal_to(j_offset, libMesh::invalid_uint);
1212  }
1213 
1214  for (auto j : IntRange<unsigned int>(0, context.get_dof_indices(v).size()))
1215  {
1216  const unsigned int total_j = j + j_offset;
1217 
1218  // Take the "minus" side of a central differenced first derivative
1219  Number original_solution = context.get_elem_solution(v)(j);
1220  context.get_elem_solution(v)(j) -= my_h;
1221 
1222  // Make sure to catch any moving mesh terms
1223  Real * coord = nullptr;
1224  if (_mesh_sys == this)
1225  {
1226  if (_mesh_x_var == v)
1227  coord = &(context.get_elem().point(j)(0));
1228  else if (_mesh_y_var == v)
1229  coord = &(context.get_elem().point(j)(1));
1230  else if (_mesh_z_var == v)
1231  coord = &(context.get_elem().point(j)(2));
1232  }
1233  if (coord)
1234  {
1235  // We have enough information to scale the perturbations
1236  // here appropriately
1237  context.get_elem_solution(v)(j) = original_solution - numerical_point_h;
1238  *coord = libmesh_real(context.get_elem_solution(v)(j));
1239  }
1240 
1241  context.get_elem_residual().zero();
1242  ((*time_solver).*(res))(false, context);
1243 #ifdef DEBUG
1244  libmesh_assert_equal_to (old_jacobian, context.get_elem_jacobian());
1245 #endif
1246  backwards_residual = context.get_elem_residual();
1247 
1248  // Take the "plus" side of a central differenced first derivative
1249  context.get_elem_solution(v)(j) = original_solution + my_h;
1250  if (coord)
1251  {
1252  context.get_elem_solution()(j) = original_solution + numerical_point_h;
1253  *coord = libmesh_real(context.get_elem_solution(v)(j));
1254  }
1255  context.get_elem_residual().zero();
1256  ((*time_solver).*(res))(false, context);
1257 #ifdef DEBUG
1258  libmesh_assert_equal_to (old_jacobian, context.get_elem_jacobian());
1259 #endif
1260 
1261  context.get_elem_solution(v)(j) = original_solution;
1262  if (coord)
1263  {
1264  *coord = libmesh_real(context.get_elem_solution(v)(j));
1265  for (auto i : IntRange<unsigned int>(0, n_dofs))
1266  {
1267  numeric_jacobian(i,total_j) =
1268  (context.get_elem_residual()(i) - backwards_residual(i)) /
1269  2. / numerical_point_h;
1270  }
1271  }
1272  else
1273  {
1274  for (auto i : IntRange<unsigned int>(0, n_dofs))
1275  {
1276  numeric_jacobian(i,total_j) =
1277  (context.get_elem_residual()(i) - backwards_residual(i)) /
1278  2. / my_h;
1279  }
1280  }
1281  }
1282  }
1283 
1284  context.get_elem_residual() = original_residual;
1285  context.get_elem_jacobian() = numeric_jacobian;
1286 }
1287 
1288 
1289 
1291 {
1292  LOG_SCOPE("numerical_elem_jacobian()", "FEMSystem");
1294 }
1295 
1296 
1297 
1299 {
1300  LOG_SCOPE("numerical_side_jacobian()", "FEMSystem");
1302 }
1303 
1304 
1305 
1307 {
1308  LOG_SCOPE("numerical_nonlocal_jacobian()", "FEMSystem");
1310 }
1311 
1312 
1313 
1314 std::unique_ptr<DiffContext> FEMSystem::build_context ()
1315 {
1316  FEMContext * fc = new FEMContext(*this);
1317 
1318  DifferentiablePhysics * phys = this->get_physics();
1319 
1320  libmesh_assert (phys);
1321 
1322  // If we are solving a moving mesh problem, tell that to the Context
1323  fc->set_mesh_system(phys->get_mesh_system());
1324  fc->set_mesh_x_var(phys->get_mesh_x_var());
1325  fc->set_mesh_y_var(phys->get_mesh_y_var());
1326  fc->set_mesh_z_var(phys->get_mesh_z_var());
1327 
1328  fc->set_deltat_pointer( &deltat );
1329 
1330  // If we are solving the adjoint problem, tell that to the Context
1331  fc->is_adjoint() = this->get_time_solver().is_adjoint();
1332 
1333  return std::unique_ptr<DiffContext>(fc);
1334 }
1335 
1336 
1337 
1339 {
1340  // Parent::init_context(c); // may be a good idea in derived classes
1341 
1342  // Although we do this in DiffSystem::build_context() and
1343  // FEMSystem::build_context() as well, we do it here just to be
1344  // extra sure that the deltat pointer gets set. Since the
1345  // intended behavior is for classes derived from FEMSystem to
1346  // call Parent::init_context() in their own init_context()
1347  // overloads, we can ensure that those classes get the correct
1348  // deltat pointers even if they have different build_context()
1349  // overloads.
1350  c.set_deltat_pointer ( &deltat );
1351 
1352  FEMContext & context = cast_ref<FEMContext &>(c);
1353 
1354  // Make sure we're prepared to do mass integration
1355  for (auto var : IntRange<unsigned int>(0, this->n_vars()))
1356  if (this->get_physics()->is_time_evolving(var))
1357  {
1358  // Request shape functions based on FEType
1359  switch( FEInterface::field_type( this->variable_type(var) ) )
1360  {
1361  case( TYPE_SCALAR ):
1362  {
1363  FEBase * elem_fe = nullptr;
1364  context.get_element_fe(var, elem_fe);
1365  elem_fe->get_JxW();
1366  elem_fe->get_phi();
1367  }
1368  break;
1369  case( TYPE_VECTOR ):
1370  {
1371  FEGenericBase<RealGradient> * elem_fe = nullptr;
1372  context.get_element_fe(var, elem_fe);
1373  elem_fe->get_JxW();
1374  elem_fe->get_phi();
1375  }
1376  break;
1377  default:
1378  libmesh_error_msg("Unrecognized field type!");
1379  }
1380  }
1381 }
1382 
1383 
1384 
1386 {
1387  // This function makes no sense unless we've already picked out some
1388  // variable(s) to reflect mesh position coordinates
1389  if (!_mesh_sys)
1390  libmesh_error_msg("_mesh_sys was nullptr!");
1391 
1392  // We currently assume mesh variables are in our own system
1393  if (_mesh_sys != this)
1394  libmesh_not_implemented();
1395 
1396  // Loop over every active mesh element on this processor
1397  const MeshBase & mesh = this->get_mesh();
1398 
1399  std::unique_ptr<DiffContext> con = this->build_context();
1400  FEMContext & _femcontext = cast_ref<FEMContext &>(*con);
1401  this->init_context(_femcontext);
1402 
1403  // Get the solution's mesh variables from every element
1404  for (const auto & elem : mesh.active_local_element_ptr_range())
1405  {
1406  _femcontext.pre_fe_reinit(*this, elem);
1407 
1408  _femcontext.elem_position_get();
1409 
1411  this->solution->insert(_femcontext.get_elem_solution(_mesh_x_var),
1412  _femcontext.get_dof_indices(_mesh_x_var) );
1414  this->solution->insert(_femcontext.get_elem_solution(_mesh_y_var),
1415  _femcontext.get_dof_indices(_mesh_y_var));
1417  this->solution->insert(_femcontext.get_elem_solution(_mesh_z_var),
1418  _femcontext.get_dof_indices(_mesh_z_var));
1419  }
1420 
1421  this->solution->close();
1422 
1423  // And make sure the current_local_solution is up to date too
1424  this->System::update();
1425 }
1426 
1427 } // namespace libMesh
libMesh::SparseMatrix::close
virtual void close()=0
Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across p...
libMesh::DifferentiableSystem::deltat
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time.
Definition: diff_system.h:260
libMesh::NumericVector::zero
virtual void zero()=0
Set all entries to zero.
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::DiffContext::get_elem_residual
const DenseVector< Number > & get_elem_residual() const
Const accessor for element residual.
Definition: diff_context.h:249
libMesh::DenseMatrix::get_transpose
void get_transpose(DenseMatrix< T > &dest) const
Put the tranposed matrix into dest.
Definition: dense_matrix_impl.h:612
libMesh::System::n_vars
unsigned int n_vars() const
Definition: system.h:2155
libMesh::FEMContext::side
unsigned char side
Current side for side_* to examine.
Definition: fem_context.h:1010
libMesh::FEMContext::elem_position_get
void elem_position_get()
Uses the geometry of elem to set the coordinate data specified by mesh_*_position configuration.
Definition: fem_context.C:1485
libMesh::FEMSystem::solve
virtual void solve() override
Invokes the solver associated with the system.
Definition: fem_system.C:1046
libMesh::invalid_uint
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value.
Definition: libmesh.h:249
libMesh::DifferentiablePhysics::_mesh_x_var
unsigned int _mesh_x_var
Variables from which to acquire moving mesh information.
Definition: diff_physics.h:548
libMesh::FEMSystem::init_context
virtual void init_context(DiffContext &) override
Definition: fem_system.C:1338
libMesh::DifferentiableSystem::print_jacobian_norms
bool print_jacobian_norms
Set print_jacobian_norms to true to print |J| whenever it is assembled.
Definition: diff_system.h:341
libMesh::DifferentiableSystem::print_element_jacobians
bool print_element_jacobians
Set print_element_jacobians to true to print each J_elem contribution.
Definition: diff_system.h:361
libMesh::ImplicitSystem
Manages consistently variables, degrees of freedom, coefficient vectors, and matrices for implicit sy...
Definition: implicit_system.h:57
libMesh::FEType::family
FEFamily family
The type of finite element.
Definition: fe_type.h:203
libMesh::ExplicitSystem::rhs
NumericVector< Number > * rhs
The system matrix.
Definition: explicit_system.h:114
libMesh::DifferentiableSystem::print_residual_norms
bool print_residual_norms
Set print_residual_norms to true to print |F| whenever it is assembled.
Definition: diff_system.h:331
libMesh::DiffContext::get_dof_indices
const std::vector< dof_id_type > & get_dof_indices() const
Accessor for element dof indices.
Definition: diff_context.h:367
libMesh::DofMap::constrain_element_vector
void constrain_element_vector(DenseVector< Number > &rhs, std::vector< dof_id_type > &dofs, bool asymmetric_constraint_rows=true) const
Constrains the element vector.
Definition: dof_map.h:2030
libMesh::FEMSystem::assemble_qoi_derivative
virtual void assemble_qoi_derivative(const QoISet &qoi_indices=QoISet(), bool include_liftfunc=true, bool apply_constraints=true) override
Runs a qoi derivative assembly loop over all elements, and if assemble_qoi_sides is true over all sid...
Definition: fem_system.C:1147
libMesh::DifferentiableSystem::print_element_solutions
bool print_element_solutions
Set print_element_solutions to true to print each U_elem input.
Definition: diff_system.h:351
libMesh::DifferentiablePhysics
This class provides a specific system class.
Definition: diff_physics.h:76
libMesh::DofMap::has_heterogenous_adjoint_constraints
bool has_heterogenous_adjoint_constraints(const unsigned int qoi_num) const
Definition: dof_map.h:1972
libMesh::FEMSystem::init_data
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
Definition: fem_system.C:843
libMesh::StoredRange::reset
StoredRange< iterator_type, object_type > & reset(const iterator_type &first, const iterator_type &last)
Resets the StoredRange to contain [first,last).
Definition: stored_range.h:222
libMesh::QoISet::has_index
bool has_index(std::size_t) const
Return whether or not this index is in the set to be calculated.
Definition: qoi_set.h:221
libMesh::FEMContext::side_fe_reinit
virtual void side_fe_reinit()
Reinitializes side FE objects on the current geometric element.
Definition: fem_context.C:1457
libMesh::libmesh_real
T libmesh_real(T a)
Definition: libmesh_common.h:166
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::MeshBase::active_local_elements_begin
virtual element_iterator active_local_elements_begin()=0
libMesh::TYPE_SCALAR
Definition: enum_fe_family.h:93
libMesh::DifferentiablePhysics::is_time_evolving
bool is_time_evolving(unsigned int var) const
Definition: diff_physics.h:278
libMesh::DifferentiableSystem::print_residuals
bool print_residuals
Set print_residuals to true to print F whenever it is assembled.
Definition: diff_system.h:336
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::DiffContext::n_vars
unsigned int n_vars() const
Number of variables in solution.
Definition: diff_context.h:99
libMesh::DenseVector::zero
virtual void zero() override
Set every element in the vector to 0.
Definition: dense_vector.h:379
libMesh::DenseMatrix::add
boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type add(const T2 factor, const DenseMatrix< T3 > &mat)
Adds factor times mat to this matrix.
Definition: dense_matrix.h:945
libMesh::FEMContext::pre_fe_reinit
virtual void pre_fe_reinit(const System &, const Elem *e)
Reinitializes local data vectors/matrices on the current geometric element.
Definition: fem_context.C:1642
libMesh::NumericVector::close
virtual void close()=0
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
libMesh::DiffContext::is_adjoint
bool is_adjoint() const
Accessor for querying whether we need to do a primal or adjoint solve.
Definition: diff_context.h:470
libMesh::FEGenericBase
This class forms the foundation from which generic finite elements may be derived.
Definition: exact_error_estimator.h:39
libMesh::FEAbstract::get_JxW
const std::vector< Real > & get_JxW() const
Definition: fe_abstract.h:244
libMesh::FEMContext::get_elem
const Elem & get_elem() const
Accessor for current Elem object.
Definition: fem_context.h:896
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::ParallelObject::comm
const Parallel::Communicator & comm() const
Definition: parallel_object.h:94
libMesh::DenseMatrix< Number >
libMesh::System::get_adjoint_rhs
NumericVector< Number > & get_adjoint_rhs(unsigned int i=0)
Definition: system.C:1019
libMesh::Threads::split
Dummy "splitting object" used to distinguish splitting constructors from copy constructors.
Definition: threads_none.h:63
libMesh::DifferentiableSystem::print_solution_norms
bool print_solution_norms
Set print_residual_norms to true to print |U| whenever it is used in an assembly() call.
Definition: diff_system.h:320
libMesh::Variable::type
const FEType & type() const
Definition: variable.h:119
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::DifferentiablePhysics::_mesh_sys
System * _mesh_sys
System from which to acquire moving mesh information.
Definition: diff_physics.h:543
libMesh::DofMap::has_heterogenous_adjoint_constraint
Number has_heterogenous_adjoint_constraint(const unsigned int qoi_num, const dof_id_type dof) const
Definition: dof_map.h:1986
libMesh::System::use_fixed_solution
bool use_fixed_solution
A boolean to be set to true by systems using elem_fixed_solution, for optional use by e....
Definition: system.h:1509
libMesh::DifferentiablePhysics::get_mesh_system
const System * get_mesh_system() const
Definition: diff_physics.h:625
libMesh::FEMSystem::FEMSystem
FEMSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
Definition: fem_system.C:826
libMesh::Elem::point
const Point & point(const unsigned int i) const
Definition: elem.h:1955
libMesh::DiffContext::get_qois
const std::vector< Number > & get_qois() const
Const accessor for QoI vector.
Definition: diff_context.h:319
libMesh::FEMSystem
This class provides a specific system class.
Definition: fem_system.h:53
libMesh::TimeSolver::element_residual
virtual bool element_residual(bool request_jacobian, DiffContext &)=0
This method uses the DifferentiablePhysics element_time_derivative(), element_constraint(),...
libMesh::DofMap::heterogenously_constrain_element_vector
void heterogenously_constrain_element_vector(const DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true, int qoi_index=-1) const
Constrains the element vector.
Definition: dof_map.h:2044
libMesh::DifferentiableSystem::time_solver
std::unique_ptr< TimeSolver > time_solver
A pointer to the solver object we're going to use.
Definition: diff_system.h:233
libMesh::System::n_qois
unsigned int n_qois() const
Number of currently active quantities of interest.
Definition: system.h:2328
libMesh::DenseMatrixBase::m
unsigned int m() const
Definition: dense_matrix_base.h:104
libMesh::FEMSystem::numerical_jacobian_h
Real numerical_jacobian_h
If calculating numeric jacobians is required, the FEMSystem will perturb each solution vector entry b...
Definition: fem_system.h:181
libMesh::System::add_adjoint_rhs
NumericVector< Number > & add_adjoint_rhs(unsigned int i=0)
Definition: system.C:1009
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::DifferentiablePhysics::_mesh_y_var
unsigned int _mesh_y_var
Definition: diff_physics.h:548
libMesh::IntRange
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
libMesh::FEMSystem::numerical_jacobian
void numerical_jacobian(TimeSolverResPtr res, FEMContext &context) const
Uses the results of multiple res calls to numerically differentiate the corresponding jacobian.
Definition: fem_system.C:1178
libMesh::NumericVector::add_vector
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i].
Definition: numeric_vector.C:363
libMesh::System::n_variable_groups
unsigned int n_variable_groups() const
Definition: system.h:2163
libMesh::TimeSolver::nonlocal_residual
virtual bool nonlocal_residual(bool request_jacobian, DiffContext &)=0
This method uses the DifferentiablePhysics nonlocal_time_derivative(), nonlocal_constraint(),...
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::Threads::parallel_for
void parallel_for(const Range &range, const Body &body)
Execute the provided function object in parallel on the specified range.
Definition: threads_none.h:73
libMesh::DiffContext::get_qoi_derivatives
const std::vector< DenseVector< Number > > & get_qoi_derivatives() const
Const accessor for QoI derivatives.
Definition: diff_context.h:331
libMesh::DifferentiableQoI
This class provides a specific system class.
Definition: diff_qoi.h:52
libMesh::System::qoi
std::vector< Number > qoi
Values of the quantities of interest.
Definition: system.h:1574
libMesh::ParallelObject::n_processors
processor_id_type n_processors() const
Definition: parallel_object.h:100
libMesh::System::variable_group
const VariableGroup & variable_group(unsigned int vg) const
Return a constant reference to VariableGroup vg.
Definition: system.h:2193
libMesh::FEMContext::has_elem
bool has_elem() const
Test for current Elem object.
Definition: fem_context.h:890
libMesh::QoISet
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
libMesh::FEMSystem::postprocess
virtual void postprocess() override
Runs a postprocessing loop over all elements, and if postprocess_sides is true over all sides.
Definition: fem_system.C:1097
libMesh::DofMap::constrain_nothing
void constrain_nothing(std::vector< dof_id_type > &dofs) const
Does not actually constrain anything, but modifies dofs in the same way as any of the constrain funct...
Definition: dof_map.h:2052
libMesh::FEMSystem::numerical_side_jacobian
void numerical_side_jacobian(FEMContext &context) const
Uses the results of multiple side_residual() calls to numerically differentiate the corresponding jac...
Definition: fem_system.C:1298
libMesh::FEMSystem::~FEMSystem
virtual ~FEMSystem()
Destructor.
Definition: fem_system.C:837
libMesh::DiffContext::set_deltat_pointer
void set_deltat_pointer(Real *dt)
Points the _deltat member of this class at a timestep value stored in the creating System,...
Definition: diff_context.C:103
libMesh::FEMSystem::mesh_position_get
void mesh_position_get()
Tells the FEMSystem to set the degree of freedom coefficients which should correspond to mesh nodal c...
Definition: fem_system.C:1385
libMesh::System::get_mesh
const MeshBase & get_mesh() const
Definition: system.h:2083
libMesh::ParallelObject::processor_id
processor_id_type processor_id() const
Definition: parallel_object.h:106
libMesh::libmesh_ignore
void libmesh_ignore(const Args &...)
Definition: libmesh_common.h:526
libMesh::DifferentiableSystem::element_postprocess
virtual void element_postprocess(DiffContext &)
Does any work that needs to be done on elem in a postprocessing loop.
Definition: diff_system.h:281
libMesh::DifferentiableSystem::print_jacobians
bool print_jacobians
Set print_jacobians to true to print J whenever it is assembled.
Definition: diff_system.h:346
libMesh::DifferentiableSystem::solve
virtual void solve() override
Invokes the solver associated with the system.
Definition: diff_system.C:152
libMesh::DifferentiableSystem::get_time_solver
TimeSolver & get_time_solver()
Definition: diff_system.h:416
libMesh::MeshBase::active_local_elements_end
virtual element_iterator active_local_elements_end()=0
libMesh::DifferentiablePhysics::get_mesh_x_var
unsigned int get_mesh_x_var() const
Definition: diff_physics.h:637
libMesh::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
libMesh::Threads::parallel_reduce
void parallel_reduce(const Range &range, Body &body)
Execute the provided reduction operation in parallel on the specified range.
Definition: threads_none.h:101
libMesh::TimeSolver::set_is_adjoint
void set_is_adjoint(bool _is_adjoint_value)
Accessor for setting whether we need to do a primal or adjoint solve.
Definition: time_solver.h:239
libMesh::DifferentiableQoI::finalize_derivative
virtual void finalize_derivative(NumericVector< Number > &derivatives, std::size_t qoi_index)
Method to finalize qoi derivatives which require more than just a simple sum of element contributions...
Definition: diff_qoi.C:53
libMesh::FEMContext::elem_position_set
void elem_position_set(Real theta)
Uses the coordinate data specified by mesh_*_position configuration to set the geometry of elem to th...
Definition: fem_context.h:1238
libMesh::DenseVector::size
virtual unsigned int size() const override
Definition: dense_vector.h:92
libMesh::NumericVector::l1_norm
virtual Real l1_norm() const =0
libMesh::SyncNodalPositions
Definition: parallel_ghost_sync.h:782
libMesh::FEMSystem::assembly
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false) override
Prepares matrix or rhs for matrix assembly.
Definition: fem_system.C:850
libMesh::DiffContext
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
libMesh::DifferentiablePhysics::get_mesh_z_var
unsigned int get_mesh_z_var() const
Definition: diff_physics.h:649
libMesh::DenseMatrix::zero
virtual void zero() override
Set every element in the matrix to 0.
Definition: dense_matrix.h:838
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::FEMSystem::build_context
virtual std::unique_ptr< DiffContext > build_context() override
Builds a FEMContext object with enough information to do evaluations on each element.
Definition: fem_system.C:1314
libMesh::MeshBase::nodes_end
virtual node_iterator nodes_end()=0
libMesh::DofMap::constrain_element_matrix_and_vector
void constrain_element_matrix_and_vector(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true) const
Constrains the element matrix and vector.
Definition: dof_map.h:2034
libMesh::DifferentiablePhysics::_mesh_z_var
unsigned int _mesh_z_var
Definition: diff_physics.h:548
libMesh::System::variable_type
const FEType & variable_type(const unsigned int i) const
Definition: system.h:2233
libMesh::FEMContext::set_mesh_system
virtual void set_mesh_system(System *sys)
Tells the FEMContext that system sys contains the isoparametric Lagrangian variables which correspond...
Definition: fem_context.h:830
libMesh::FEMSystem::mesh_position_set
void mesh_position_set()
Tells the FEMSystem to set the mesh nodal coordinates which should correspond to degree of freedom co...
Definition: fem_system.C:1057
libMesh::TimeSolver::is_adjoint
bool is_adjoint() const
Accessor for querying whether we need to do a primal or adjoint solve.
Definition: time_solver.h:232
libMesh::FEMContext::elem_fe_reinit
virtual void elem_fe_reinit(const std::vector< Point > *const pts=nullptr)
Reinitializes interior FE objects on the current geometric element.
Definition: fem_context.C:1436
libMesh::DofMap::constrain_element_matrix
void constrain_element_matrix(DenseMatrix< Number > &matrix, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true) const
Constrains the element matrix.
Definition: dof_map.h:2021
libMesh::Parallel::sync_dofobject_data_by_id
void sync_dofobject_data_by_id(const Communicator &comm, const Iterator &range_begin, const Iterator &range_end, SyncFunctor &sync)
Request data about a range of ghost dofobjects uniquely identified by their id.
Definition: parallel_ghost_sync.h:338
libMesh::DifferentiableSystem::diff_qoi
DifferentiableQoI * diff_qoi
Pointer to object to use for quantity of interest assembly evaluations.
Definition: diff_system.h:377
libMesh::DifferentiableQoI::parallel_op
virtual void parallel_op(const Parallel::Communicator &communicator, std::vector< Number > &sys_qoi, std::vector< Number > &local_qoi, const QoISet &qoi_indices)
Method to populate system qoi data structure with process-local qoi.
Definition: diff_qoi.C:41
libMesh::FEMContext::get_element_fe
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for interior finite element object for variable var for the largest dimension in the mesh.
Definition: fem_context.h:275
libMesh::System::solution
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1539
libMesh::SparseMatrix::add_matrix
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
Add the full matrix dm to the SparseMatrix.
libMesh::MeshBase::nodes_begin
virtual node_iterator nodes_begin()=0
Iterate over all the nodes in the Mesh.
libMesh::FEMContext::set_mesh_y_var
void set_mesh_y_var(unsigned int y_var)
Accessor for y-variable of moving mesh System.
Definition: fem_context.h:870
libMesh::FEMSystem::verify_analytic_jacobians
Real verify_analytic_jacobians
If verify_analytic_jacobian is equal to zero (as it is by default), no numeric jacobians will be calc...
Definition: fem_system.h:209
libMesh::DofMap::heterogenously_constrain_element_matrix_and_vector
void heterogenously_constrain_element_matrix_and_vector(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true, int qoi_index=-1) const
Constrains the element matrix and vector.
Definition: dof_map.h:2040
libMesh::Elem::hmin
virtual Real hmin() const
Definition: elem.C:359
libMesh::FEInterface::field_type
static FEFieldType field_type(const FEType &fe_type)
Definition: fe_interface.C:1683
libMesh::DofObject::id
dof_id_type id() const
Definition: dof_object.h:767
libMesh::TYPE_VECTOR
Definition: enum_fe_family.h:94
libMesh::DiffContext::get_elem_jacobian
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
Definition: diff_context.h:283
libMesh::FEMSystem::fe_reinit_during_postprocess
bool fe_reinit_during_postprocess
If fe_reinit_during_postprocess is true (it is true by default), FE objects will be reinit()ed with t...
Definition: fem_system.h:170
libMesh::DiffContext::get_elem_solution
const DenseVector< Number > & get_elem_solution() const
Accessor for element solution.
Definition: diff_context.h:111
libMesh::Elem::has_children
bool has_children() const
Definition: elem.h:2383
libMesh::StoredRange
The StoredRange class defines a contiguous, divisible set of objects.
Definition: stored_range.h:52
libMesh::BasicOStreamProxy::precision
std::streamsize precision() const
Get the associated write precision.
Definition: ostream_proxy.h:189
libMesh::DifferentiableSystem::get_physics
const DifferentiablePhysics * get_physics() const
Definition: diff_system.h:181
libMesh::DifferentiableSystem::init_data
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
Definition: diff_system.C:108
libMesh::DiffContext::get_elem_fixed_solution
const DenseVector< Number > & get_elem_fixed_solution() const
Accessor for element fixed solution.
Definition: diff_context.h:215
libMesh::FEMSystem::assemble_qoi
virtual void assemble_qoi(const QoISet &indices=QoISet()) override
Runs a qoi assembly loop over all elements, and if assemble_qoi_sides is true over all sides.
Definition: fem_system.C:1117
libMesh::DifferentiableSystem::postprocess_sides
bool postprocess_sides
If postprocess_sides is true (it is false by default), the postprocessing loop will loop over all sid...
Definition: diff_system.h:314
libMesh::DifferentiableSystem::print_solutions
bool print_solutions
Set print_solutions to true to print U whenever it is used in an assembly() call.
Definition: diff_system.h:326
libMesh::DifferentiableSystem::side_postprocess
virtual void side_postprocess(DiffContext &)
Does any work that needs to be done on side of elem in a postprocessing loop.
Definition: diff_system.h:287
libMesh::SparseMatrix::l1_norm
virtual Real l1_norm() const =0
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::System::n_dofs
dof_id_type n_dofs() const
Definition: system.C:150
libMesh::DenseMatrix::l1_norm
auto l1_norm() const -> decltype(std::abs(T(0)))
Definition: dense_matrix.h:1052
libMesh::err
OStreamProxy err
libMesh::SCALAR
Definition: enum_fe_family.h:58
libMesh::Elem::neighbor_ptr
const Elem * neighbor_ptr(unsigned int i) const
Definition: elem.h:2085
libMesh::FEGenericBase::get_phi
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:206
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::FEMSystem::numerical_jacobian_h_for_var
Real numerical_jacobian_h_for_var(unsigned int var_num) const
If numerical_jacobian_h_for_var(var_num) is changed from its default value (numerical_jacobian_h),...
Definition: fem_system.h:259
libMesh::FEMSystem::numerical_nonlocal_jacobian
void numerical_nonlocal_jacobian(FEMContext &context) const
Uses the results of multiple side_residual() calls to numerically differentiate the corresponding jac...
Definition: fem_system.C:1306
libMesh::Elem::n_sides
virtual unsigned int n_sides() const =0
libMesh::DifferentiablePhysics::get_mesh_y_var
unsigned int get_mesh_y_var() const
Definition: diff_physics.h:643
libMesh::TimeSolver::side_residual
virtual bool side_residual(bool request_jacobian, DiffContext &)=0
This method uses the DifferentiablePhysics side_time_derivative(), side_constraint(),...
libMesh::FEMSystem::numerical_elem_jacobian
void numerical_elem_jacobian(FEMContext &context) const
Uses the results of multiple element_residual() calls to numerically differentiate the corresponding ...
Definition: fem_system.C:1290
libMesh::out
OStreamProxy out
libMesh::DifferentiableSystem::print_element_residuals
bool print_element_residuals
Set print_element_residuals to true to print each R_elem contribution.
Definition: diff_system.h:356
libMesh::System::update
virtual void update()
Update the local values to reflect the solution on neighboring processors.
Definition: system.C:408
libMesh::SparseMatrix::zero
virtual void zero()=0
Set all entries to 0.
libMesh::DifferentiablePhysics::compute_internal_sides
bool compute_internal_sides
compute_internal_sides is false by default, indicating that side_* computations will only be done on ...
Definition: diff_physics.h:154
libMesh::Threads::spin_mutex
Spin mutex.
Definition: threads_none.h:127
libMesh::FEMContext::set_mesh_z_var
void set_mesh_z_var(unsigned int z_var)
Accessor for z-variable of moving mesh System.
Definition: fem_context.h:884
libMesh::DenseVector< Number >
libMesh::FEMContext::set_mesh_x_var
void set_mesh_x_var(unsigned int x_var)
Accessor for x-variable of moving mesh System.
Definition: fem_context.h:856
libMesh::FEMContext
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62