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