libMesh
system.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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 // Local includes
21 #include "libmesh/dof_map.h"
22 #include "libmesh/equation_systems.h"
23 #include "libmesh/int_range.h"
24 #include "libmesh/libmesh_logging.h"
25 #include "libmesh/mesh_base.h"
26 #include "libmesh/numeric_vector.h"
27 #include "libmesh/parameter_vector.h"
28 #include "libmesh/point.h" // For point_value
29 #include "libmesh/point_locator_base.h" // For point_value
30 #include "libmesh/qoi_set.h"
31 #include "libmesh/enum_to_string.h"
32 #include "libmesh/sparse_matrix.h"
33 #include "libmesh/system.h"
34 #include "libmesh/system_norm.h"
35 #include "libmesh/utility.h"
36 #include "libmesh/elem.h"
37 #include "libmesh/fe_type.h"
38 #include "libmesh/fe_interface.h"
39 #include "libmesh/fe_compute_data.h"
40 #include "libmesh/static_condensation.h"
41 
42 // includes for calculate_norm, point_*
43 #include "libmesh/fe_base.h"
44 #include "libmesh/fe_interface.h"
45 #include "libmesh/parallel.h"
46 #include "libmesh/parallel_algebra.h"
47 #include "libmesh/quadrature.h"
48 #include "libmesh/tensor_value.h"
49 #include "libmesh/vector_value.h"
50 #include "libmesh/tensor_tools.h"
51 #include "libmesh/enum_norm_type.h"
52 #include "libmesh/enum_fe_family.h"
53 
54 // C++ includes
55 #include <sstream> // for std::ostringstream
56 
57 namespace libMesh
58 {
59 
60 
61 // ------------------------------------------------------------
62 // System implementation
64  const std::string & name_in,
65  const unsigned int number_in) :
66 
67  ParallelObject (es),
68  assemble_before_solve (true),
69  use_fixed_solution (false),
70  extra_quadrature_order (0),
71  solution (NumericVector<Number>::build(this->comm())),
72  current_local_solution (NumericVector<Number>::build(this->comm())),
73  time (0.),
74  qoi (0),
75  qoi_error_estimates (0),
76  _init_system_function (nullptr),
77  _init_system_object (nullptr),
78  _assemble_system_function (nullptr),
79  _assemble_system_object (nullptr),
80  _constrain_system_function (nullptr),
81  _constrain_system_object (nullptr),
82  _qoi_evaluate_function (nullptr),
83  _qoi_evaluate_object (nullptr),
84  _qoi_evaluate_derivative_function (nullptr),
85  _qoi_evaluate_derivative_object (nullptr),
86  _dof_map (std::make_unique<DofMap>(number_in, es.get_mesh())),
87  _equation_systems (es),
88  _mesh (es.get_mesh()),
89  _sys_name (name_in),
90  _sys_number (number_in),
91  _active (true),
92  _matrices_initialized (false),
93  _solution_projection (true),
94  _basic_system_only (false),
95  _is_initialized (false),
96  _additional_data_written (false),
97  adjoint_already_solved (false),
98  _hide_output (false),
99  project_with_constraints (true),
100  _prefer_hash_table_matrix_assembly(false),
101  _require_sparsity_pattern (false),
102  _prefix_with_name (false)
103 {
104  if (libMesh::on_command_line("--solver-system-names"))
105  this->prefix_with_name(true);
106  if (libMesh::on_command_line("--" + name_in + "-static-condensation"))
108 }
109 
110 
111 
113 {
114  libmesh_exceptionless_assert (!libMesh::closed());
115 }
116 
117 
118 
120 {
121  return _dof_map->n_dofs();
122 }
123 
124 
125 
127 {
128 #ifdef LIBMESH_ENABLE_CONSTRAINTS
129 
130  return _dof_map->n_constrained_dofs();
131 
132 #else
133 
134  return 0;
135 
136 #endif
137 }
138 
139 
140 
142 {
143 #ifdef LIBMESH_ENABLE_CONSTRAINTS
144 
145  return _dof_map->n_local_constrained_dofs();
146 
147 #else
148 
149  return 0;
150 
151 #endif
152 }
153 
154 
155 
157 {
158  return _dof_map->n_local_dofs();
159 }
160 
161 
162 
163 Number System::current_solution (const dof_id_type global_dof_number) const
164 {
165  // Check the sizes
166  libmesh_assert_less (global_dof_number, _dof_map->n_dofs());
167  libmesh_assert_less (global_dof_number, current_local_solution->size());
168 
169  return (*current_local_solution)(global_dof_number);
170 }
171 
172 
173 
175 {
176  _dof_map->clear ();
177  solution->clear ();
178  current_local_solution->clear ();
179 
180  // clear any user-added vectors
181  _vectors.clear();
182  _vector_projections.clear();
183  _vector_is_adjoint.clear();
184  _is_initialized = false;
185 
186  // clear any user-added matrices
187  _matrices.clear();
188  _matrices_initialized = false;
189 }
190 
191 
192 
194 {
195  // Calling init() twice on the same system currently works evil
196  // magic, whether done directly or via EquationSystems::read()
197  libmesh_assert(!this->is_initialized());
198 
199  this->reinit_mesh();
200 }
201 
202 
203 
205 {
206  parallel_object_only();
207 
208  MeshBase & mesh = this->get_mesh();
209 
210  // Distribute the degrees of freedom on the mesh
211  auto total_dofs = _dof_map->distribute_dofs (mesh);
212 
213  // Throw an error if the total number of DOFs is not capable of
214  // being indexed by our solution vector.
215  auto max_allowed_id = solution->max_allowed_id();
216  libmesh_error_msg_if(total_dofs > max_allowed_id,
217  "Cannot allocate a NumericVector with " << total_dofs << " degrees of freedom. "
218  "The vector can only index up to " << max_allowed_id << " entries.");
219 
220  // Recreate any user or internal constraints
221  this->reinit_constraints();
222 
223  // Even if there weren't any constraint changes,
224  // reinit_constraints() did prepare_send_list() for us.
225 
226  // Now finally after dof distribution and construction of any
227  // possible constraints, we may init any static condensation
228  // data
229  _dof_map->reinit_static_condensation();
230 
231  // Resize the solution conformal to the current mesh
232  solution->init (this->n_dofs(), this->n_local_dofs(), false, PARALLEL);
233 
234  // Resize the current_local_solution for the current mesh
235 #ifdef LIBMESH_ENABLE_GHOSTED
236  current_local_solution->init (this->n_dofs(), this->n_local_dofs(),
237  _dof_map->get_send_list(), /*fast=*/false,
238  GHOSTED);
239 #else
240  current_local_solution->init (this->n_dofs(), false, SERIAL);
241 #endif
242 
243  // from now on, adding additional vectors or variables can't be done
244  // without immediately initializing them
245  _is_initialized = true;
246 
247  // initialize & zero other vectors, if necessary
248  for (auto & [vec_name, vec] : _vectors)
249  {
250  libmesh_ignore(vec_name); // spurious warning from old gcc
251  const ParallelType type = vec->type();
252 
253  if (type == GHOSTED)
254  {
255 #ifdef LIBMESH_ENABLE_GHOSTED
256  vec->init (this->n_dofs(), this->n_local_dofs(),
257  _dof_map->get_send_list(), /*fast=*/false,
258  GHOSTED);
259 #else
260  libmesh_error_msg("Cannot initialize ghosted vectors when they are not enabled.");
261 #endif
262  }
263  else if (type == SERIAL)
264  {
265  vec->init (this->n_dofs(), false, type);
266  }
267  else
268  {
269  libmesh_assert_equal_to(type, PARALLEL);
270  vec->init (this->n_dofs(), this->n_local_dofs(), false, type);
271  }
272  }
273 
274  // Add matrices
275  this->add_matrices();
276 
277  // Clear any existing matrices
278  for (auto & pr : _matrices)
279  pr.second->clear();
280 
281  // Initialize the matrices for the system
282  if (!_basic_system_only)
283  this->init_matrices();
284 }
285 
287 {
288  parallel_object_only();
289 
290  // First initialize any required data:
291  // either only the basic System data
292  if (_basic_system_only)
294  // or all the derived class' data too
295  else
296  this->init_data();
297 
298  // If no variables have been added to this system
299  // don't do anything
300  if (!this->n_vars())
301  return;
302 
303  // Then call the user-provided initialization function
304  this->user_initialization();
305 
306 }
307 
309 {
310  parallel_object_only();
311 
312  // No matrices to init
313  if (_matrices.empty())
314  {
315  // any future matrices to be added will need their own
316  // initialization
317  _matrices_initialized = true;
318 
319  return;
320  }
321 
322  // Check for quick return in case the first matrix
323  // (and by extension all the matrices) has already
324  // been initialized
325  if (_matrices.begin()->second->initialized())
326  {
328  return;
329  }
330 
331  _matrices_initialized = true;
332 
333  // Tell the matrices about the dof map, and vice versa
334  for (auto & pr : _matrices)
335  {
336  SparseMatrix<Number> & m = *(pr.second);
338 
339  // We want to allow repeated init() on systems, but we don't
340  // want to attach the same matrix to the DofMap twice
341  if (!this->get_dof_map().is_attached(m))
342  this->get_dof_map().attach_matrix(m);
343 
344  // If the user has already explicitly requested that this matrix use a hash table, then we
345  // always honor that
346  const bool use_hash =
347  pr.second->use_hash_table() ||
348  (this->_prefer_hash_table_matrix_assembly && pr.second->supports_hash_table());
349  pr.second->use_hash_table(use_hash);
350  // Make this call after we've determined whether the matrix is using a hash table
351  if (pr.second->require_sparsity_pattern())
352  this->_require_sparsity_pattern = true;
353  }
354 
355  // Compute the sparsity pattern for the current
356  // mesh and DOF distribution. This also updates
357  // additional matrices, \p DofMap now knows them
358  if (this->_require_sparsity_pattern)
359  this->get_dof_map().compute_sparsity(this->get_mesh());
360 
361  // Initialize matrices and set to zero
362  for (auto & [name, mat] : _matrices)
363  {
364  mat->init(_matrix_types[name]);
365  mat->zero();
366  }
367 }
368 
369 
370 
372 {
373  parallel_object_only();
374 
375 #ifdef LIBMESH_ENABLE_AMR
376  // Restrict the _vectors on the coarsened cells
377  for (auto & [vec_name, vec] : _vectors)
378  {
379  NumericVector<Number> * v = vec.get();
380 
381  if (_vector_projections[vec_name])
382  {
383  this->project_vector (*v, this->vector_is_adjoint(vec_name));
384  }
385  else
386  {
387  const ParallelType type = vec->type();
388 
389  if (type == GHOSTED)
390  {
391 #ifdef LIBMESH_ENABLE_GHOSTED
392  vec->init (this->n_dofs(), this->n_local_dofs(),
393  _dof_map->get_send_list(), /*fast=*/false,
394  GHOSTED);
395 #else
396  libmesh_error_msg("Cannot initialize ghosted vectors when they are not enabled.");
397 #endif
398  }
399  else
400  vec->init (this->n_dofs(), this->n_local_dofs(), false, type);
401  }
402  }
403 
404  const std::vector<dof_id_type> & send_list = _dof_map->get_send_list ();
405 
406  // Restrict the solution on the coarsened cells
408  this->project_vector (*solution);
409  // Or at least make sure the solution vector is the correct size
410  else
411  solution->init (this->n_dofs(), this->n_local_dofs(), true, PARALLEL);
412 
413 #ifdef LIBMESH_ENABLE_GHOSTED
414  current_local_solution->init(this->n_dofs(),
415  this->n_local_dofs(), send_list,
416  false, GHOSTED);
417 #else
418  current_local_solution->init(this->n_dofs());
419 #endif
420 
422  solution->localize (*current_local_solution, send_list);
423 
424 #endif // LIBMESH_ENABLE_AMR
425 }
426 
427 
428 
430 {
431 #ifdef LIBMESH_ENABLE_AMR
432  // Currently project_vector handles both restriction and prolongation
433  this->restrict_vectors();
434 #endif
435 }
436 
437 
438 
440 {
441  parallel_object_only();
442 
443  // project_vector handles vector initialization now
444  libmesh_assert_equal_to (solution->size(), current_local_solution->size());
445 
446  // Make sure our static condensation dof map is up-to-date before we init any
447  // static condensation matrices
449 
450  if (!_matrices.empty() && !_basic_system_only)
451  {
452  // Clear the matrices
453  for (auto & pr : _matrices)
454  {
455  pr.second->clear();
456  pr.second->attach_dof_map(this->get_dof_map());
457  }
458 
459  if (this->_require_sparsity_pattern)
460  {
461  // Clear the sparsity pattern
462  this->get_dof_map().clear_sparsity();
463 
464  // Compute the sparsity pattern for the current
465  // mesh and DOF distribution. This also updates
466  // additional matrices, \p DofMap now knows them
467  this->get_dof_map().compute_sparsity (this->get_mesh());
468  }
469 
470  // Initialize matrices and set to zero
471  for (auto & pr : _matrices)
472  {
473  pr.second->init();
474  pr.second->zero();
475  }
476  }
477 }
478 
479 
481 {
482  parallel_object_only();
483 
484 #ifdef LIBMESH_ENABLE_CONSTRAINTS
486  user_constrain();
488  if (libMesh::on_command_line ("--print-constraints"))
490 #endif
492 }
493 
494 
496 {
497  parallel_object_only();
498 
499  libmesh_assert(solution->closed());
500 
501  const std::vector<dof_id_type> & send_list = _dof_map->get_send_list ();
502 
503  // Check sizes
504  libmesh_assert_equal_to (current_local_solution->size(), solution->size());
505  // More processors than elements => empty send_list
506  // libmesh_assert (!send_list.empty());
507  libmesh_assert_less_equal (send_list.size(), solution->size());
508 
509  // Create current_local_solution from solution. This will
510  // put a local copy of solution into current_local_solution.
511  // Only the necessary values (specified by the send_list)
512  // are copied to minimize communication
513  solution->localize (*current_local_solution, send_list);
514 }
515 
516 
517 
519 {
520  parallel_object_only();
521 
522  // If this system is empty... don't do anything!
523  if (!this->n_vars())
524  return;
525 
526  const std::vector<dof_id_type> & send_list = this->get_dof_map().get_send_list ();
527 
528  // Check sizes
529  libmesh_assert_equal_to (current_local_solution->size(), solution->size());
530  // Not true with ghosted vectors
531  // libmesh_assert_equal_to (current_local_solution->local_size(), solution->size());
532  // libmesh_assert (!send_list.empty());
533  libmesh_assert_less_equal (send_list.size(), solution->size());
534 
535  // Create current_local_solution from solution. This will
536  // put a local copy of solution into current_local_solution.
537  solution->localize (*current_local_solution, send_list);
538 }
539 
540 
541 
543  const SubsetSolveMode /*subset_solve_mode*/)
544 {
545  if (subset != nullptr)
546  libmesh_not_implemented();
547 }
548 
549 
550 
552 {
553  // Log how long the user's assembly code takes
554  LOG_SCOPE("assemble()", "System");
555 
556  // Call the user-specified assembly function
557  this->user_assembly();
558 }
559 
560 
561 
562 void System::assemble_qoi (const QoISet & qoi_indices)
563 {
564  // Log how long the user's assembly code takes
565  LOG_SCOPE("assemble_qoi()", "System");
566 
567  // Call the user-specified quantity of interest function
568  this->user_QOI(qoi_indices);
569 }
570 
571 
572 
573 void System::assemble_qoi_derivative(const QoISet & qoi_indices,
574  bool include_liftfunc,
575  bool apply_constraints)
576 {
577  // Log how long the user's assembly code takes
578  LOG_SCOPE("assemble_qoi_derivative()", "System");
579 
580  // Call the user-specified quantity of interest function
581  this->user_QOI_derivative(qoi_indices, include_liftfunc,
582  apply_constraints);
583 }
584 
585 
586 
587 void System::qoi_parameter_sensitivity (const QoISet & qoi_indices,
588  const ParameterVector & parameters_vec,
589  SensitivityData & sensitivities)
590 {
591  // Forward sensitivities are more efficient for Nq > Np
592  if (qoi_indices.size(*this) > parameters_vec.size())
593  forward_qoi_parameter_sensitivity(qoi_indices, parameters_vec, sensitivities);
594  // Adjoint sensitivities are more efficient for Np > Nq,
595  // and an adjoint may be more reusable than a forward
596  // solution sensitivity in the Np == Nq case.
597  else
598  adjoint_qoi_parameter_sensitivity(qoi_indices, parameters_vec, sensitivities);
599 }
600 
601 
602 
603 bool System::compare (const System & other_system,
604  const Real threshold,
605  const bool verbose) const
606 {
607  // we do not care for matrices, but for vectors
609  libmesh_assert (other_system._is_initialized);
610 
611  if (verbose)
612  {
613  libMesh::out << " Systems \"" << _sys_name << "\"" << std::endl;
614  libMesh::out << " comparing matrices not supported." << std::endl;
615  libMesh::out << " comparing names...";
616  }
617 
618  // compare the name: 0 means identical
619  const int name_result = _sys_name.compare(other_system.name());
620  if (verbose)
621  {
622  if (name_result == 0)
623  libMesh::out << " identical." << std::endl;
624  else
625  libMesh::out << " names not identical." << std::endl;
626  libMesh::out << " comparing solution vector...";
627  }
628 
629 
630  // compare the solution: -1 means identical
631  const int solu_result = solution->compare (*other_system.solution.get(),
632  threshold);
633 
634  if (verbose)
635  {
636  if (solu_result == -1)
637  libMesh::out << " identical up to threshold." << std::endl;
638  else
639  libMesh::out << " first difference occurred at index = "
640  << solu_result << "." << std::endl;
641  }
642 
643 
644  // safety check, whether we handle at least the same number
645  // of vectors
646  std::vector<int> ov_result;
647 
648  if (this->n_vectors() != other_system.n_vectors())
649  {
650  if (verbose)
651  {
652  libMesh::out << " Fatal difference. This system handles "
653  << this->n_vectors() << " add'l vectors," << std::endl
654  << " while the other system handles "
655  << other_system.n_vectors()
656  << " add'l vectors." << std::endl
657  << " Aborting comparison." << std::endl;
658  }
659  return false;
660  }
661  else if (this->n_vectors() == 0)
662  {
663  // there are no additional vectors...
664  ov_result.clear ();
665  }
666  else
667  {
668  // compare other vectors
669  for (auto & [vec_name, vec] : _vectors)
670  {
671  if (verbose)
672  libMesh::out << " comparing vector \""
673  << vec_name << "\" ...";
674 
675  // assume they have the same name
676  const NumericVector<Number> & other_system_vector =
677  other_system.get_vector(vec_name);
678 
679  ov_result.push_back(vec->compare(other_system_vector, threshold));
680 
681  if (verbose)
682  {
683  if (ov_result[ov_result.size()-1] == -1)
684  libMesh::out << " identical up to threshold." << std::endl;
685  else
686  libMesh::out << " first difference occurred at" << std::endl
687  << " index = " << ov_result[ov_result.size()-1] << "." << std::endl;
688  }
689  }
690  } // finished comparing additional vectors
691 
692 
693  bool overall_result;
694 
695  // sum up the results
696  if ((name_result==0) && (solu_result==-1))
697  {
698  if (ov_result.size()==0)
699  overall_result = true;
700  else
701  {
702  bool ov_identical;
703  unsigned int n = 0;
704  do
705  {
706  ov_identical = (ov_result[n]==-1);
707  n++;
708  }
709  while (ov_identical && n<ov_result.size());
710  overall_result = ov_identical;
711  }
712  }
713  else
714  overall_result = false;
715 
716  if (verbose)
717  {
718  libMesh::out << " finished comparisons, ";
719  if (overall_result)
720  libMesh::out << "found no differences." << std::endl << std::endl;
721  else
722  libMesh::out << "found differences." << std::endl << std::endl;
723  }
724 
725  return overall_result;
726 }
727 
728 
729 
730 void System::update_global_solution (std::vector<Number> & global_soln) const
731 {
732  parallel_object_only();
733 
734  global_soln.resize (solution->size());
735 
736  solution->localize (global_soln);
737 }
738 
739 
740 
741 void System::update_global_solution (std::vector<Number> & global_soln,
742  const processor_id_type dest_proc) const
743 {
744  parallel_object_only();
745 
746  global_soln.resize (solution->size());
747 
748  solution->localize_to_one (global_soln, dest_proc);
749 }
750 
751 
752 
753 NumericVector<Number> & System::add_vector (std::string_view vec_name,
754  const bool projections,
755  const ParallelType type)
756 {
757  parallel_object_only();
758 
759  libmesh_assert(this->comm().verify(std::string(vec_name)));
760  libmesh_assert(this->comm().verify(int(type)));
761  libmesh_assert(this->comm().verify(projections));
762 
763  // Return the vector if it is already there.
764  if (auto it = this->_vectors.find(vec_name);
765  it != this->_vectors.end())
766  {
767  // If the projection setting has *upgraded*, change it.
768  if (projections) // only do expensive lookup if needed
769  libmesh_map_find(_vector_projections, vec_name) = projections;
770 
771  NumericVector<Number> & vec = *it->second;
772 
773  // If we're in serial, our vectors are effectively SERIAL, so
774  // we'll ignore any type setting. If we're in parallel, we
775  // might have a type change to deal with.
776 
777  if (this->n_processors() > 1)
778  {
779  // If the type setting has changed in a way we can't
780  // perceive as an upgrade or a downgrade, scream.
781  libmesh_assert_equal_to(type == SERIAL,
782  vec.type() == SERIAL);
783 
784  // If the type setting has *upgraded*, change it.
785  if (type == GHOSTED && vec.type() == PARALLEL)
786  {
787  // A *really* late upgrade is expensive, but better not
788  // to risk zeroing data.
789  if (vec.initialized())
790  {
791  if (!vec.closed())
792  vec.close();
793 
794  // Ideally we'd move parallel coefficients and then
795  // add ghosted coefficients, but copy and swap is
796  // simpler. If anyone actually ever uses this case
797  // for real we can look into optimizing it.
798  auto new_vec = NumericVector<Number>::build(this->comm());
799 #ifdef LIBMESH_ENABLE_GHOSTED
800  new_vec->init (this->n_dofs(), this->n_local_dofs(),
801  _dof_map->get_send_list(), /*fast=*/false,
802  GHOSTED);
803 #else
804  libmesh_error_msg("Cannot initialize ghosted vectors when they are not enabled.");
805 #endif
806 
807  *new_vec = vec;
808  vec.swap(*new_vec);
809  }
810  else
811  // The PARALLEL vec is not yet initialized, so we can
812  // just "upgrade" it to GHOSTED.
813  vec.set_type(type);
814  }
815  }
816 
817  // Any upgrades are done; we're happy here.
818  return vec;
819  }
820 
821  // Otherwise, build the vector. The following emplace() is
822  // guaranteed to succeed because, if we made it here, we don't
823  // already have a vector named "vec_name". We pass the user's
824  // requested ParallelType directly to NumericVector::build() so
825  // that, even if the vector is not initialized now, it will get the
826  // right type when it is initialized later.
827  auto pr =
828  _vectors.emplace(vec_name,
831  type));
832  auto buf = pr.first->second.get();
833  _vector_projections.emplace(vec_name, projections);
834 
835  // Vectors are primal by default
836  _vector_is_adjoint.emplace(vec_name, -1);
837 
838  // Initialize it if necessary
839  if (_is_initialized)
840  {
841  if (type == GHOSTED)
842  {
843 #ifdef LIBMESH_ENABLE_GHOSTED
844  buf->init (this->n_dofs(), this->n_local_dofs(),
845  _dof_map->get_send_list(), /*fast=*/false,
846  GHOSTED);
847 #else
848  libmesh_error_msg("Cannot initialize ghosted vectors when they are not enabled.");
849 #endif
850  }
851  else
852  buf->init (this->n_dofs(), this->n_local_dofs(), false, type);
853  }
854 
855  return *buf;
856 }
857 
858 void System::remove_vector (std::string_view vec_name)
859 {
860  parallel_object_only(); // Not strictly needed, but the only safe way to keep in sync
861 
862  if (const auto pos = _vectors.find(vec_name);
863  pos != _vectors.end())
864  {
865  _vectors.erase(pos);
866  auto proj_it = _vector_projections.find(vec_name);
867  libmesh_assert(proj_it != _vector_projections.end());
868  _vector_projections.erase(proj_it);
869 
870  auto adj_it = _vector_is_adjoint.find(vec_name);
871  libmesh_assert(adj_it != _vector_is_adjoint.end());
872  _vector_is_adjoint.erase(adj_it);
873  }
874 }
875 
876 const NumericVector<Number> * System::request_vector (std::string_view vec_name) const
877 {
878  if (const auto pos = _vectors.find(vec_name);
879  pos != _vectors.end())
880  return pos->second.get();
881 
882  // Otherwise, vec_name was not found
883  return nullptr;
884 }
885 
886 
887 
888 NumericVector<Number> * System::request_vector (std::string_view vec_name)
889 {
890  if (auto pos = _vectors.find(vec_name);
891  pos != _vectors.end())
892  return pos->second.get();
893 
894  // Otherwise, vec_name was not found
895  return nullptr;
896 }
897 
898 
899 
900 const NumericVector<Number> * System::request_vector (const unsigned int vec_num) const
901 {
902  // If we don't have that many vectors, return nullptr
903  if (vec_num >= _vectors.size())
904  return nullptr;
905 
906  // Otherwise return a pointer to the vec_num'th vector
907  auto it = vectors_begin();
908  std::advance(it, vec_num);
909  return it->second.get();
910 }
911 
912 
913 
914 NumericVector<Number> * System::request_vector (const unsigned int vec_num)
915 {
916  // If we don't have that many vectors, return nullptr
917  if (vec_num >= _vectors.size())
918  return nullptr;
919 
920  // Otherwise return a pointer to the vec_num'th vector
921  auto it = vectors_begin();
922  std::advance(it, vec_num);
923  return it->second.get();
924 }
925 
926 
927 
928 const NumericVector<Number> & System::get_vector (std::string_view vec_name) const
929 {
930  return *(libmesh_map_find(_vectors, vec_name));
931 }
932 
933 
934 
935 NumericVector<Number> & System::get_vector (std::string_view vec_name)
936 {
937  return *(libmesh_map_find(_vectors, vec_name));
938 }
939 
940 
941 
942 const NumericVector<Number> & System::get_vector (const unsigned int vec_num) const
943 {
944  // If we don't have that many vectors, throw an error
945  libmesh_assert_less(vec_num, _vectors.size());
946 
947  // Otherwise return a reference to the vec_num'th vector
948  auto it = vectors_begin();
949  std::advance(it, vec_num);
950  return *(it->second);
951 }
952 
953 
954 
955 NumericVector<Number> & System::get_vector (const unsigned int vec_num)
956 {
957  // If we don't have that many vectors, throw an error
958  libmesh_assert_less(vec_num, _vectors.size());
959 
960  // Otherwise return a reference to the vec_num'th vector
961  auto it = vectors_begin();
962  std::advance(it, vec_num);
963  return *(it->second);
964 }
965 
966 
967 
968 const std::string & System::vector_name (const unsigned int vec_num) const
969 {
970  // If we don't have that many vectors, throw an error
971  libmesh_assert_less(vec_num, _vectors.size());
972 
973  // Otherwise return a reference to the vec_num'th vector name
974  auto it = vectors_begin();
975  std::advance(it, vec_num);
976  return it->first;
977 }
978 
979 const std::string & System::vector_name (const NumericVector<Number> & vec_reference) const
980 {
981  // Linear search for a vector whose pointer matches vec_reference
982  auto it = std::find_if(vectors_begin(), vectors_end(),
983  [&vec_reference](const decltype(_vectors)::value_type & pr)
984  { return &vec_reference == pr.second.get(); });
985 
986  // Before returning, make sure we didn't loop till the end and not find any match
987  libmesh_assert (it != vectors_end());
988 
989  // Return the string associated with the current vector
990  return it->first;
991 }
992 
993 
994 
995 SparseMatrix<Number> & System::add_matrix (std::string_view mat_name,
996  const ParallelType type,
997  const MatrixBuildType mat_build_type)
998 {
999  parallel_object_only();
1000 
1001  libmesh_assert(this->comm().verify(std::string(mat_name)));
1002  libmesh_assert(this->comm().verify(int(type)));
1003  libmesh_assert(this->comm().verify(int(mat_build_type)));
1004 
1005  // Return the matrix if it is already there.
1006  if (auto it = this->_matrices.find(mat_name);
1007  it != this->_matrices.end())
1008  return *it->second;
1009 
1010  // Otherwise build the matrix to return.
1011  std::unique_ptr<SparseMatrix<Number>> matrix;
1012  if (this->has_static_condensation())
1013  {
1014  if (mat_build_type == MatrixBuildType::DIAGONAL)
1015  libmesh_error_msg(
1016  "We do not currently support static condensation of the diagonal matrix type");
1017  matrix = std::make_unique<StaticCondensation>(this->get_mesh(),
1018  *this,
1019  this->get_dof_map(),
1021  }
1022  else
1024  auto & mat = *matrix;
1025 
1026  _matrices.emplace(mat_name, std::move(matrix));
1027 
1028  _matrix_types.emplace(mat_name, type);
1029 
1030  // Initialize it first if we've already initialized the others.
1031  this->late_matrix_init(mat, type);
1032 
1033  return mat;
1034 }
1035 
1036 
1037 
1038 SparseMatrix<Number> & System::add_matrix (std::string_view mat_name,
1039  std::unique_ptr<SparseMatrix<Number>> matrix,
1040  const ParallelType type)
1041 {
1042  parallel_object_only();
1043 
1044  libmesh_assert(this->comm().verify(std::string(mat_name)));
1045  libmesh_assert(this->comm().verify(int(type)));
1046 
1047  auto [it, inserted] = _matrices.emplace(mat_name, std::move(matrix));
1048  libmesh_error_msg_if(!inserted,
1049  "Tried to add '" << mat_name << "' but the matrix already exists");
1050 
1051  _matrix_types.emplace(mat_name, type);
1052 
1053  SparseMatrix<Number> & mat = *(it->second);
1054 
1055  // Initialize it first if we've already initialized the others.
1056  this->late_matrix_init(mat, type);
1057 
1058  return mat;
1059 }
1060 
1062  ParallelType type)
1063 {
1065  {
1066  this->get_dof_map().attach_matrix(mat);
1067  mat.init(type);
1068  }
1069 }
1070 
1071 
1072 
1073 
1074 void System::remove_matrix (std::string_view mat_name)
1075 {
1076  parallel_object_only(); // Not strictly needed, but the only safe way to keep in sync
1077 
1078  if (const auto pos = _matrices.find(mat_name);
1079  pos != _matrices.end())
1080  _matrices.erase(pos); // erase()'d entries are destroyed
1081 }
1082 
1083 
1084 
1085 const SparseMatrix<Number> * System::request_matrix (std::string_view mat_name) const
1086 {
1087  if (const auto pos = _matrices.find(mat_name);
1088  pos != _matrices.end())
1089  return pos->second.get();
1090 
1091  // Otherwise, mat_name does not exist
1092  return nullptr;
1093 }
1094 
1095 
1096 
1097 SparseMatrix<Number> * System::request_matrix (std::string_view mat_name)
1098 {
1099  if (auto pos = _matrices.find(mat_name);
1100  pos != _matrices.end())
1101  return pos->second.get();
1102 
1103  // Otherwise, mat_name does not exist
1104  return nullptr;
1105 }
1106 
1107 
1108 
1109 const SparseMatrix<Number> & System::get_matrix (std::string_view mat_name) const
1110 {
1111  return *libmesh_map_find(_matrices, mat_name);
1112 }
1113 
1114 
1115 
1116 SparseMatrix<Number> & System::get_matrix (std::string_view mat_name)
1117 {
1118  return *libmesh_map_find(_matrices, mat_name);
1119 }
1120 
1121 
1122 
1123 void System::set_vector_preservation (const std::string & vec_name,
1124  bool preserve)
1125 {
1126  parallel_object_only(); // Not strictly needed, but the only safe way to keep in sync
1127 
1128  _vector_projections[vec_name] = preserve;
1129 }
1130 
1131 
1132 
1133 bool System::vector_preservation (std::string_view vec_name) const
1134 {
1135  if (auto it = _vector_projections.find(vec_name);
1136  it != _vector_projections.end())
1137  return it->second;
1138 
1139  // vec_name was not in the map, return false
1140  return false;
1141 }
1142 
1143 
1144 
1145 void System::set_vector_as_adjoint (const std::string & vec_name,
1146  int qoi_num)
1147 {
1148  parallel_object_only(); // Not strictly needed, but the only safe way to keep in sync
1149 
1150  // We reserve -1 for vectors which get primal constraints, -2 for
1151  // vectors which get no constraints
1152  libmesh_assert_greater_equal(qoi_num, -2);
1153  _vector_is_adjoint[vec_name] = qoi_num;
1154 }
1155 
1156 
1157 
1158 int System::vector_is_adjoint (std::string_view vec_name) const
1159 {
1160  const auto it = _vector_is_adjoint.find(vec_name);
1161  libmesh_assert(it != _vector_is_adjoint.end());
1162  return it->second;
1163 }
1164 
1165 
1166 
1168 {
1169  std::ostringstream sensitivity_name;
1170  sensitivity_name << "sensitivity_solution" << i;
1171 
1172  return this->add_vector(sensitivity_name.str());
1173 }
1174 
1175 
1176 
1178 {
1179  std::ostringstream sensitivity_name;
1180  sensitivity_name << "sensitivity_solution" << i;
1181 
1182  return this->get_vector(sensitivity_name.str());
1183 }
1184 
1185 
1186 
1188 {
1189  std::ostringstream sensitivity_name;
1190  sensitivity_name << "sensitivity_solution" << i;
1191 
1192  return this->get_vector(sensitivity_name.str());
1193 }
1194 
1195 
1196 
1198 {
1199  return this->add_vector("weighted_sensitivity_solution");
1200 }
1201 
1202 
1203 
1205 {
1206  return this->get_vector("weighted_sensitivity_solution");
1207 }
1208 
1209 
1210 
1212 {
1213  return this->get_vector("weighted_sensitivity_solution");
1214 }
1215 
1216 
1217 
1219 {
1220  std::ostringstream adjoint_name;
1221  adjoint_name << "adjoint_solution" << i;
1222 
1223  NumericVector<Number> & returnval = this->add_vector(adjoint_name.str());
1224  this->set_vector_as_adjoint(adjoint_name.str(), i);
1225  return returnval;
1226 }
1227 
1228 
1229 
1231 {
1232  std::ostringstream adjoint_name;
1233  adjoint_name << "adjoint_solution" << i;
1234 
1235  return this->get_vector(adjoint_name.str());
1236 }
1237 
1238 
1239 
1241 {
1242  std::ostringstream adjoint_name;
1243  adjoint_name << "adjoint_solution" << i;
1244 
1245  return this->get_vector(adjoint_name.str());
1246 }
1247 
1248 
1249 
1251 {
1252  std::ostringstream adjoint_name;
1253  adjoint_name << "weighted_sensitivity_adjoint_solution" << i;
1254 
1255  NumericVector<Number> & returnval = this->add_vector(adjoint_name.str());
1256  this->set_vector_as_adjoint(adjoint_name.str(), i);
1257  return returnval;
1258 }
1259 
1260 
1261 
1263 {
1264  std::ostringstream adjoint_name;
1265  adjoint_name << "weighted_sensitivity_adjoint_solution" << i;
1266 
1267  return this->get_vector(adjoint_name.str());
1268 }
1269 
1270 
1271 
1273 {
1274  std::ostringstream adjoint_name;
1275  adjoint_name << "weighted_sensitivity_adjoint_solution" << i;
1276 
1277  return this->get_vector(adjoint_name.str());
1278 }
1279 
1280 
1281 
1283 {
1284  std::ostringstream adjoint_rhs_name;
1285  adjoint_rhs_name << "adjoint_rhs" << i;
1286 
1287  return this->add_vector(adjoint_rhs_name.str(), false);
1288 }
1289 
1290 
1291 
1293 {
1294  std::ostringstream adjoint_rhs_name;
1295  adjoint_rhs_name << "adjoint_rhs" << i;
1296 
1297  return this->get_vector(adjoint_rhs_name.str());
1298 }
1299 
1300 
1301 
1302 const NumericVector<Number> & System::get_adjoint_rhs (unsigned int i) const
1303 {
1304  std::ostringstream adjoint_rhs_name;
1305  adjoint_rhs_name << "adjoint_rhs" << i;
1306 
1307  return this->get_vector(adjoint_rhs_name.str());
1308 }
1309 
1310 
1311 
1313 {
1314  std::ostringstream sensitivity_rhs_name;
1315  sensitivity_rhs_name << "sensitivity_rhs" << i;
1316 
1317  return this->add_vector(sensitivity_rhs_name.str(), false);
1318 }
1319 
1320 
1321 
1323 {
1324  std::ostringstream sensitivity_rhs_name;
1325  sensitivity_rhs_name << "sensitivity_rhs" << i;
1326 
1327  return this->get_vector(sensitivity_rhs_name.str());
1328 }
1329 
1330 
1331 
1333 {
1334  std::ostringstream sensitivity_rhs_name;
1335  sensitivity_rhs_name << "sensitivity_rhs" << i;
1336 
1337  return this->get_vector(sensitivity_rhs_name.str());
1338 }
1339 
1340 
1341 
1342 unsigned int System::add_variable (std::string_view var,
1343  const FEType & type,
1344  const std::set<subdomain_id_type> * const active_subdomains)
1345 {
1346  return this->get_dof_map().add_variable(*this, var, type, active_subdomains);
1347 }
1348 
1349 
1350 
1351 unsigned int System::add_variable (std::string_view var,
1352  const Order order,
1353  const FEFamily family,
1354  const std::set<subdomain_id_type> * const active_subdomains)
1355 {
1356  return this->add_variable(var,
1357  FEType(order, family),
1358  active_subdomains);
1359 }
1360 
1361 
1362 
1363 unsigned int System::add_variables (const std::vector<std::string> & vars,
1364  const FEType & type,
1365  const std::set<subdomain_id_type> * const active_subdomains)
1366 {
1367  return this->get_dof_map().add_variables(*this, vars, type, active_subdomains);
1368 }
1369 
1370 
1371 
1372 unsigned int System::add_variables (const std::vector<std::string> & vars,
1373  const Order order,
1374  const FEFamily family,
1375  const std::set<subdomain_id_type> * const active_subdomains)
1376 {
1377  return this->add_variables(vars,
1378  FEType(order, family),
1379  active_subdomains);
1380 }
1381 
1382 unsigned int System::add_variable_array (const std::vector<std::string> & vars,
1383  const FEType & type,
1384  const std::set<subdomain_id_type> * const active_subdomains)
1385 {
1386  return this->get_dof_map().add_variable_array(*this, vars, type, active_subdomains);
1387 }
1388 
1389 bool System::has_variable (std::string_view var) const
1390 {
1391  return this->get_dof_map().has_variable(var);
1392 }
1393 
1394 unsigned int System::variable_number (std::string_view var) const
1395 {
1396  return this->get_dof_map().variable_number(var);
1397 }
1398 
1399 void System::get_all_variable_numbers(std::vector<unsigned int> & all_variable_numbers) const
1400 {
1401  this->get_dof_map().get_all_variable_numbers(all_variable_numbers);
1402 }
1403 
1404 
1405 void System::local_dof_indices(const unsigned int var,
1406  std::set<dof_id_type> & var_indices) const
1407 {
1408  // Make sure the set is clear
1409  var_indices.clear();
1410 
1411  std::vector<dof_id_type> dof_indices;
1412 
1413  const dof_id_type
1414  first_local = this->get_dof_map().first_dof(),
1415  end_local = this->get_dof_map().end_dof();
1416 
1417  // Begin the loop over the elements
1418  for (const auto & elem : this->get_mesh().active_local_element_ptr_range())
1419  {
1420  this->get_dof_map().dof_indices (elem, dof_indices, var);
1421 
1422  for (dof_id_type dof : dof_indices)
1423  //If the dof is owned by the local processor
1424  if (first_local <= dof && dof < end_local)
1425  var_indices.insert(dof);
1426  }
1427 
1428  // we may have missed assigning DOFs to nodes that we own
1429  // but to which we have no connected elements matching our
1430  // variable restriction criterion. this will happen, for example,
1431  // if variable V is restricted to subdomain S. We may not own
1432  // any elements which live in S, but we may own nodes which are
1433  // *connected* to elements which do.
1434  for (const auto & node : this->get_mesh().local_node_ptr_range())
1435  {
1436  libmesh_assert(node);
1437  this->get_dof_map().dof_indices (node, dof_indices, var);
1438  for (auto dof : dof_indices)
1439  if (first_local <= dof && dof < end_local)
1440  var_indices.insert(dof);
1441  }
1442 }
1443 
1444 
1445 
1447  unsigned int var_num) const
1448 {
1449  /* Make sure the call makes sense. */
1450  libmesh_assert_less (var_num, this->n_vars());
1451 
1452  /* Get a reference to the mesh. */
1453  const MeshBase & mesh = this->get_mesh();
1454 
1455  /* Check which system we are. */
1456  const unsigned int sys_num = this->number();
1457 
1458  // Loop over nodes.
1459  for (const auto & node : mesh.local_node_ptr_range())
1460  {
1461  unsigned int n_comp = node->n_comp(sys_num,var_num);
1462  for (unsigned int i=0; i<n_comp; i++)
1463  {
1464  const dof_id_type index = node->dof_number(sys_num,var_num,i);
1465  v.set(index,0.0);
1466  }
1467  }
1468 
1469  // Loop over elements.
1470  for (const auto & elem : mesh.active_local_element_ptr_range())
1471  {
1472  unsigned int n_comp = elem->n_comp(sys_num,var_num);
1473  for (unsigned int i=0; i<n_comp; i++)
1474  {
1475  const dof_id_type index = elem->dof_number(sys_num,var_num,i);
1476  v.set(index,0.0);
1477  }
1478  }
1479 }
1480 
1481 
1482 
1484  unsigned int var,
1485  FEMNormType norm_type) const
1486 {
1487  std::set<dof_id_type> var_indices;
1488  local_dof_indices(var, var_indices);
1489 
1490  if (norm_type == DISCRETE_L1)
1491  return v.subset_l1_norm(var_indices);
1492  if (norm_type == DISCRETE_L2)
1493  return v.subset_l2_norm(var_indices);
1494  if (norm_type == DISCRETE_L_INF)
1495  return v.subset_linfty_norm(var_indices);
1496  else
1497  libmesh_error_msg("Invalid norm_type = " << Utility::enum_to_string(norm_type));
1498 }
1499 
1500 
1501 
1503  unsigned int var,
1504  FEMNormType norm_type,
1505  std::set<unsigned int> * skip_dimensions) const
1506 {
1507  //short circuit to save time
1508  if (norm_type == DISCRETE_L1 ||
1509  norm_type == DISCRETE_L2 ||
1510  norm_type == DISCRETE_L_INF)
1511  return discrete_var_norm(v,var,norm_type);
1512 
1513  // Not a discrete norm
1514  std::vector<FEMNormType> norms(this->n_vars(), L2);
1515  std::vector<Real> weights(this->n_vars(), 0.0);
1516  norms[var] = norm_type;
1517  weights[var] = 1.0;
1518  Real val = this->calculate_norm(v, SystemNorm(norms, weights), skip_dimensions);
1519  return val;
1520 }
1521 
1522 
1523 
1525  const SystemNorm & norm,
1526  std::set<unsigned int> * skip_dimensions) const
1527 {
1528  // This function must be run on all processors at once
1529  parallel_object_only();
1530 
1531  LOG_SCOPE ("calculate_norm()", "System");
1532 
1533  // Zero the norm before summation
1534  Real v_norm = 0.;
1535 
1536  if (norm.is_discrete())
1537  {
1538  //Check to see if all weights are 1.0 and all types are equal
1539  FEMNormType norm_type0 = norm.type(0);
1540  unsigned int check_var = 0, check_end = this->n_vars();
1541  for (; check_var != check_end; ++check_var)
1542  if ((norm.weight(check_var) != 1.0) || (norm.type(check_var) != norm_type0))
1543  break;
1544 
1545  //All weights were 1.0 so just do the full vector discrete norm
1546  if (check_var == this->n_vars())
1547  {
1548  if (norm_type0 == DISCRETE_L1)
1549  return v.l1_norm();
1550  if (norm_type0 == DISCRETE_L2)
1551  return v.l2_norm();
1552  if (norm_type0 == DISCRETE_L_INF)
1553  return v.linfty_norm();
1554  else
1555  libmesh_error_msg("Invalid norm_type0 = " << Utility::enum_to_string(norm_type0));
1556  }
1557 
1558  for (auto var : make_range(this->n_vars()))
1559  {
1560  // Skip any variables we don't need to integrate
1561  if (norm.weight(var) == 0.0)
1562  continue;
1563 
1564  v_norm += norm.weight(var) * discrete_var_norm(v, var, norm.type(var));
1565  }
1566 
1567  return v_norm;
1568  }
1569 
1570  // Localize the potentially parallel vector
1571  std::unique_ptr<NumericVector<Number>> local_v = NumericVector<Number>::build(this->comm());
1572  local_v->init(v.size(), v.local_size(), _dof_map->get_send_list(),
1573  true, GHOSTED);
1574  v.localize (*local_v, _dof_map->get_send_list());
1575 
1576  // I'm not sure how best to mix Hilbert norms on some variables (for
1577  // which we'll want to square then sum then square root) with norms
1578  // like L_inf (for which we'll just want to take an absolute value
1579  // and then sum).
1580  bool using_hilbert_norm = true,
1581  using_nonhilbert_norm = true;
1582 
1583  // Loop over all variables
1584  for (auto var : make_range(this->n_vars()))
1585  {
1586  // Skip any variables we don't need to integrate
1587  Real norm_weight_sq = norm.weight_sq(var);
1588  if (norm_weight_sq == 0.0)
1589  continue;
1590  Real norm_weight = norm.weight(var);
1591 
1592  // Check for unimplemented norms (rather than just returning 0).
1593  FEMNormType norm_type = norm.type(var);
1594  if ((norm_type==H1) ||
1595  (norm_type==H2) ||
1596  (norm_type==L2) ||
1597  (norm_type==H1_SEMINORM) ||
1598  (norm_type==H2_SEMINORM))
1599  {
1600  if (!using_hilbert_norm)
1601  libmesh_not_implemented();
1602  using_nonhilbert_norm = false;
1603  }
1604  else if ((norm_type==L1) ||
1605  (norm_type==L_INF) ||
1606  (norm_type==W1_INF_SEMINORM) ||
1607  (norm_type==W2_INF_SEMINORM))
1608  {
1609  if (!using_nonhilbert_norm)
1610  libmesh_not_implemented();
1611  using_hilbert_norm = false;
1612  }
1613  else
1614  libmesh_not_implemented();
1615 
1616  const FEType & fe_type = this->get_dof_map().variable_type(var);
1617 
1618  // Allow space for dims 0-3, and for both scalar and vector
1619  // elements, even if we don't use them all
1620  std::vector<std::unique_ptr<FEBase>> fe_ptrs(4);
1621  std::vector<std::unique_ptr<FEVectorBase>> vec_fe_ptrs(4);
1622  std::vector<std::unique_ptr<QBase>> q_rules(4);
1623 
1624  const std::set<unsigned char> & elem_dims = _mesh.elem_dimensions();
1625 
1626  // Prepare finite elements for each dimension present in the mesh
1627  for (const auto & dim : elem_dims)
1628  {
1629  if (skip_dimensions && skip_dimensions->find(dim) != skip_dimensions->end())
1630  continue;
1631 
1632  // Construct quadrature and finite element objects
1633  q_rules[dim] = fe_type.default_quadrature_rule (dim);
1634 
1635  const FEFieldType field_type = FEInterface::field_type(fe_type);
1636  if (field_type == TYPE_SCALAR)
1637  {
1638  fe_ptrs[dim] = FEBase::build(dim, fe_type);
1639  fe_ptrs[dim]->attach_quadrature_rule (q_rules[dim].get());
1640  }
1641  else
1642  {
1643  vec_fe_ptrs[dim] = FEVectorBase::build(dim, fe_type);
1644  vec_fe_ptrs[dim]->attach_quadrature_rule (q_rules[dim].get());
1645  libmesh_assert_equal_to(field_type, TYPE_VECTOR);
1646  }
1647 
1648  }
1649 
1650  std::vector<dof_id_type> dof_indices;
1651 
1652  // Begin the loop over the elements
1653  for (const auto & elem : this->get_mesh().active_local_element_ptr_range())
1654  {
1655  const unsigned int dim = elem->dim();
1656 
1657  // One way for implementing this would be to exchange the fe with the FEInterface- class.
1658  // However, it needs to be discussed whether integral-norms make sense for infinite elements.
1659  // or in which sense they could make sense.
1660  if (elem->infinite() )
1661  libmesh_not_implemented();
1662 
1663  if (skip_dimensions && skip_dimensions->find(dim) != skip_dimensions->end())
1664  continue;
1665 
1666  QBase * qrule = q_rules[dim].get();
1667  libmesh_assert(qrule);
1668 
1669  this->get_dof_map().dof_indices (elem, dof_indices, var);
1670 
1671  auto element_calculation = [&dof_indices, &elem,
1672  norm_type, norm_weight, norm_weight_sq, &qrule,
1673  &local_v, &v_norm](auto & fe) {
1674  typedef typename std::remove_reference<decltype(fe)>::type::OutputShape OutputShape;
1675  typedef typename TensorTools::MakeNumber<OutputShape>::type OutputNumberShape;
1676  typedef typename std::remove_reference<decltype(fe)>::type::OutputGradient OutputGradient;
1677  typedef typename TensorTools::MakeNumber<OutputGradient>::type OutputNumberGradient;
1678 
1679  const std::vector<Real> & JxW = fe.get_JxW();
1680  const std::vector<std::vector<OutputShape>> * phi = nullptr;
1681  if (norm_type == H1 ||
1682  norm_type == H2 ||
1683  norm_type == L2 ||
1684  norm_type == L1 ||
1685  norm_type == L_INF)
1686  phi = &(fe.get_phi());
1687 
1688  const std::vector<std::vector<OutputGradient>> * dphi = nullptr;
1689  if (norm_type == H1 ||
1690  norm_type == H2 ||
1691  norm_type == H1_SEMINORM ||
1692  norm_type == W1_INF_SEMINORM)
1693  dphi = &(fe.get_dphi());
1694 
1695 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1696  typedef typename std::remove_reference<decltype(fe)>::type::OutputTensor OutputTensor;
1697 
1698  const std::vector<std::vector<OutputTensor>> * d2phi = nullptr;
1699  if (norm_type == H2 ||
1700  norm_type == H2_SEMINORM ||
1701  norm_type == W2_INF_SEMINORM)
1702  d2phi = &(fe.get_d2phi());
1703 #endif
1704 
1705  fe.reinit (elem);
1706 
1707  const unsigned int n_qp = qrule->n_points();
1708 
1709  const unsigned int n_sf = cast_int<unsigned int>
1710  (dof_indices.size());
1711 
1712  // Begin the loop over the Quadrature points.
1713  for (unsigned int qp=0; qp<n_qp; qp++)
1714  {
1715  if (norm_type == L1)
1716  {
1717  OutputNumberShape u_h = 0.;
1718  for (unsigned int i=0; i != n_sf; ++i)
1719  u_h += (*phi)[i][qp] * (*local_v)(dof_indices[i]);
1720  v_norm += norm_weight *
1721  JxW[qp] * TensorTools::norm(u_h);
1722  }
1723 
1724  if (norm_type == L_INF)
1725  {
1726  OutputNumberShape u_h = 0.;
1727  for (unsigned int i=0; i != n_sf; ++i)
1728  u_h += (*phi)[i][qp] * (*local_v)(dof_indices[i]);
1729  v_norm = std::max(v_norm, norm_weight * TensorTools::norm(u_h));
1730  }
1731 
1732  if (norm_type == H1 ||
1733  norm_type == H2 ||
1734  norm_type == L2)
1735  {
1736  OutputNumberShape u_h = 0.;
1737  for (unsigned int i=0; i != n_sf; ++i)
1738  u_h += (*phi)[i][qp] * (*local_v)(dof_indices[i]);
1739  v_norm += norm_weight_sq *
1740  JxW[qp] * TensorTools::norm_sq(u_h);
1741  }
1742 
1743  if (norm_type == H1 ||
1744  norm_type == H2 ||
1745  norm_type == H1_SEMINORM)
1746  {
1747  OutputNumberGradient grad_u_h;
1748  for (unsigned int i=0; i != n_sf; ++i)
1749  grad_u_h.add_scaled((*dphi)[i][qp], (*local_v)(dof_indices[i]));
1750  v_norm += norm_weight_sq *
1751  JxW[qp] * grad_u_h.norm_sq();
1752  }
1753 
1754  if (norm_type == W1_INF_SEMINORM)
1755  {
1756  OutputNumberGradient grad_u_h;
1757  for (unsigned int i=0; i != n_sf; ++i)
1758  grad_u_h.add_scaled((*dphi)[i][qp], (*local_v)(dof_indices[i]));
1759  v_norm = std::max(v_norm, norm_weight * grad_u_h.norm());
1760  }
1761 
1762 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1763  typedef typename TensorTools::MakeNumber<OutputTensor>::type OutputNumberTensor;
1764 
1765  if (norm_type == H2 ||
1766  norm_type == H2_SEMINORM)
1767  {
1768  OutputNumberTensor hess_u_h;
1769  for (unsigned int i=0; i != n_sf; ++i)
1770  hess_u_h.add_scaled((*d2phi)[i][qp], (*local_v)(dof_indices[i]));
1771  v_norm += norm_weight_sq *
1772  JxW[qp] * hess_u_h.norm_sq();
1773  }
1774 
1775  if (norm_type == W2_INF_SEMINORM)
1776  {
1777  OutputNumberTensor hess_u_h;
1778  for (unsigned int i=0; i != n_sf; ++i)
1779  hess_u_h.add_scaled((*d2phi)[i][qp], (*local_v)(dof_indices[i]));
1780  v_norm = std::max(v_norm, norm_weight * hess_u_h.norm());
1781  }
1782 #endif
1783  }
1784  };
1785 
1786  FEBase * scalar_fe = fe_ptrs[dim].get();
1787  FEVectorBase * vec_fe = vec_fe_ptrs[dim].get();
1788 
1789  if (scalar_fe)
1790  {
1791  libmesh_assert(!vec_fe);
1792  element_calculation(*scalar_fe);
1793  }
1794 
1795  if (vec_fe)
1796  {
1797  libmesh_assert(!scalar_fe);
1798  element_calculation(*vec_fe);
1799  }
1800  }
1801  }
1802 
1803  if (using_hilbert_norm)
1804  {
1805  this->comm().sum(v_norm);
1806  v_norm = std::sqrt(v_norm);
1807  }
1808  else
1809  {
1810  this->comm().max(v_norm);
1811  }
1812 
1813  return v_norm;
1814 }
1815 
1816 
1817 
1818 std::string System::get_info() const
1819 {
1820  std::ostringstream oss;
1821 
1822 
1823  const std::string & sys_name = this->name();
1824 
1825  oss << " System #" << this->number() << ", \"" << sys_name << "\"\n"
1826  << " Type \"" << this->system_type() << "\"\n"
1827  << " Variables=";
1828 
1829  for (auto vg : make_range(this->n_variable_groups()))
1830  {
1831  const VariableGroup & vg_description (this->variable_group(vg));
1832 
1833  if (vg_description.n_variables() > 1) oss << "{ ";
1834  for (auto vn : make_range(vg_description.n_variables()))
1835  oss << "\"" << vg_description.name(vn) << "\" ";
1836  if (vg_description.n_variables() > 1) oss << "} ";
1837  }
1838 
1839  oss << '\n';
1840 
1841  oss << " Finite Element Types=";
1842 #ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS
1843  for (auto vg : make_range(this->n_variable_groups()))
1844  oss << "\""
1845  << Utility::enum_to_string<FEFamily>(this->get_dof_map().variable_group(vg).type().family)
1846  << "\" ";
1847 #else
1848  for (auto vg : make_range(this->n_variable_groups()))
1849  {
1850  oss << "\""
1851  << Utility::enum_to_string<FEFamily>(this->get_dof_map().variable_group(vg).type().family)
1852  << "\", \""
1853  << Utility::enum_to_string<FEFamily>(this->get_dof_map().variable_group(vg).type().radial_family)
1854  << "\" ";
1855  }
1856 
1857  oss << '\n' << " Infinite Element Mapping=";
1858  for (auto vg : make_range(this->n_variable_groups()))
1859  oss << "\""
1860  << Utility::enum_to_string<InfMapType>(this->get_dof_map().variable_group(vg).type().inf_map)
1861  << "\" ";
1862 #endif
1863 
1864  oss << '\n';
1865 
1866  oss << " Approximation Orders=";
1867  for (auto vg : make_range(this->n_variable_groups()))
1868  {
1869 #ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS
1870  oss << "\""
1871  << Utility::enum_to_string<Order>(this->get_dof_map().variable_group(vg).type().order)
1872  << "\" ";
1873 #else
1874  oss << "\""
1875  << Utility::enum_to_string<Order>(this->get_dof_map().variable_group(vg).type().order)
1876  << "\", \""
1877  << Utility::enum_to_string<Order>(this->get_dof_map().variable_group(vg).type().radial_order)
1878  << "\" ";
1879 #endif
1880  }
1881 
1882  oss << '\n';
1883 
1884  oss << " n_dofs()=" << this->n_dofs() << '\n';
1885  dof_id_type local_dofs = this->n_local_dofs();
1886  oss << " n_local_dofs()=" << local_dofs << '\n';
1887  this->comm().max(local_dofs);
1888  oss << " max(n_local_dofs())=" << local_dofs << '\n';
1889 #ifdef LIBMESH_ENABLE_CONSTRAINTS
1890  oss << " n_constrained_dofs()=" << this->n_constrained_dofs() << '\n';
1891  oss << " n_local_constrained_dofs()=" << this->n_local_constrained_dofs() << '\n';
1892  dof_id_type local_unconstrained_dofs = this->n_local_dofs() - this->n_local_constrained_dofs();
1893  this->comm().max(local_unconstrained_dofs);
1894  oss << " max(local unconstrained dofs)=" << local_unconstrained_dofs << '\n';
1895 #endif
1896 
1897  oss << " " << "n_vectors()=" << this->n_vectors() << '\n';
1898  oss << " " << "n_matrices()=" << this->n_matrices() << '\n';
1899  // oss << " " << "n_additional_matrices()=" << this->n_additional_matrices() << '\n';
1900 
1901  oss << this->get_dof_map().get_info();
1902 
1903  return oss.str();
1904 }
1905 
1906 
1907 
1909  const std::string & name))
1910 {
1912 
1913  if (_init_system_object != nullptr)
1914  {
1915  libmesh_warning("WARNING: Cannot specify both initialization function and object!");
1916 
1917  _init_system_object = nullptr;
1918  }
1919 
1921 }
1922 
1923 
1924 
1926 {
1927  if (_init_system_function != nullptr)
1928  {
1929  libmesh_warning("WARNING: Cannot specify both initialization object and function!");
1930 
1931  _init_system_function = nullptr;
1932  }
1933 
1934  _init_system_object = &init_in;
1935 }
1936 
1937 
1938 
1940  const std::string & name))
1941 {
1943 
1944  if (_assemble_system_object != nullptr)
1945  {
1946  libmesh_warning("WARNING: Cannot specify both assembly function and object!");
1947 
1948  _assemble_system_object = nullptr;
1949  }
1950 
1952 }
1953 
1954 
1955 
1957 {
1958  if (_assemble_system_function != nullptr)
1959  {
1960  libmesh_warning("WARNING: Cannot specify both assembly object and function!");
1961 
1962  _assemble_system_function = nullptr;
1963  }
1964 
1965  _assemble_system_object = &assemble_in;
1966 }
1967 
1968 
1969 
1971  const std::string & name))
1972 {
1974 
1975  if (_constrain_system_object != nullptr)
1976  {
1977  libmesh_warning("WARNING: Cannot specify both constraint function and object!");
1978 
1979  _constrain_system_object = nullptr;
1980  }
1981 
1983 }
1984 
1985 
1986 
1988 {
1989  if (_constrain_system_function != nullptr)
1990  {
1991  libmesh_warning("WARNING: Cannot specify both constraint object and function!");
1992 
1993  _constrain_system_function = nullptr;
1994  }
1995 
1996  _constrain_system_object = &constrain;
1997 }
1998 
2000 {
2001  return _constrain_system_object != nullptr;
2002 }
2003 
2005 {
2006  libmesh_assert_msg(_constrain_system_object,"No constraint object available.");
2007  return *_constrain_system_object;
2008 }
2009 
2010 
2011 
2013  const std::string &,
2014  const QoISet &))
2015 {
2017 
2018  if (_qoi_evaluate_object != nullptr)
2019  {
2020  libmesh_warning("WARNING: Cannot specify both QOI function and object!");
2021 
2022  _qoi_evaluate_object = nullptr;
2023  }
2024 
2026 }
2027 
2028 
2029 
2031 {
2032  if (_qoi_evaluate_function != nullptr)
2033  {
2034  libmesh_warning("WARNING: Cannot specify both QOI object and function!");
2035 
2036  _qoi_evaluate_function = nullptr;
2037  }
2038 
2039  _qoi_evaluate_object = &qoi_in;
2040 }
2041 
2042 
2043 
2044 void System::attach_QOI_derivative(void fptr(EquationSystems &, const std::string &,
2045  const QoISet &, bool, bool))
2046 {
2048 
2049  if (_qoi_evaluate_derivative_object != nullptr)
2050  {
2051  libmesh_warning("WARNING: Cannot specify both QOI derivative function and object!");
2052 
2054  }
2055 
2057 }
2058 
2059 
2060 
2062 {
2063  if (_qoi_evaluate_derivative_function != nullptr)
2064  {
2065  libmesh_warning("WARNING: Cannot specify both QOI derivative object and function!");
2066 
2068  }
2069 
2070  _qoi_evaluate_derivative_object = &qoi_derivative;
2071 }
2072 
2073 
2074 
2076 {
2077  // Call the user-provided initialization function,
2078  // if it was provided
2079  if (_init_system_function != nullptr)
2080  this->_init_system_function (_equation_systems, this->name());
2081 
2082  // ...or the user-provided initialization object.
2083  else if (_init_system_object != nullptr)
2085 }
2086 
2087 
2088 
2090 {
2091  // Call the user-provided assembly function,
2092  // if it was provided
2093  if (_assemble_system_function != nullptr)
2095 
2096  // ...or the user-provided assembly object.
2097  else if (_assemble_system_object != nullptr)
2099 }
2100 
2101 
2102 
2104 {
2105  // Call the user-provided constraint function,
2106  // if it was provided
2107  if (_constrain_system_function!= nullptr)
2109 
2110  // ...or the user-provided constraint object.
2111  else if (_constrain_system_object != nullptr)
2113 }
2114 
2115 
2116 
2117 void System::user_QOI (const QoISet & qoi_indices)
2118 {
2119  // Call the user-provided quantity of interest function,
2120  // if it was provided
2121  if (_qoi_evaluate_function != nullptr)
2122  this->_qoi_evaluate_function(_equation_systems, this->name(), qoi_indices);
2123 
2124  // ...or the user-provided QOI function object.
2125  else if (_qoi_evaluate_object != nullptr)
2126  this->_qoi_evaluate_object->qoi(qoi_indices);
2127 }
2128 
2129 
2130 
2131 void System::user_QOI_derivative(const QoISet & qoi_indices,
2132  bool include_liftfunc,
2133  bool apply_constraints)
2134 {
2135  // Call the user-provided quantity of interest derivative,
2136  // if it was provided
2137  if (_qoi_evaluate_derivative_function != nullptr)
2139  (_equation_systems, this->name(), qoi_indices, include_liftfunc,
2140  apply_constraints);
2141 
2142  // ...or the user-provided QOI derivative function object.
2143  else if (_qoi_evaluate_derivative_object != nullptr)
2145  (qoi_indices, include_liftfunc, apply_constraints);
2146 }
2147 
2148 
2149 void System::init_qois(unsigned int n_qois)
2150 {
2151  qoi.resize(n_qois);
2152  qoi_error_estimates.resize(n_qois);
2153 }
2154 
2155 
2156 void System::set_qoi(unsigned int qoi_index, Number qoi_value)
2157 {
2158  libmesh_assert(qoi_index < qoi.size());
2159 
2160  qoi[qoi_index] = qoi_value;
2161 }
2162 
2163 
2164 Number System::get_qoi_value(unsigned int qoi_index) const
2165 {
2166  libmesh_assert(qoi_index < qoi.size());
2167  return qoi[qoi_index];
2168 }
2169 
2170 
2171 std::vector<Number> System::get_qoi_values() const
2172 {
2173  return this->qoi;
2174 }
2175 
2176 
2177 void System::set_qoi(std::vector<Number> new_qoi)
2178 {
2179  libmesh_assert_equal_to(this->qoi.size(), new_qoi.size());
2180  this->qoi = std::move(new_qoi);
2181 }
2182 
2183 
2184 void System::set_qoi_error_estimate(unsigned int qoi_index, Number qoi_error_estimate)
2185 {
2186  libmesh_assert(qoi_index < qoi_error_estimates.size());
2187 
2188  qoi_error_estimates[qoi_index] = qoi_error_estimate;
2189 }
2190 
2191 Number System::get_qoi_error_estimate_value(unsigned int qoi_index) const
2192 {
2193  libmesh_assert(qoi_index < qoi_error_estimates.size());
2194  return qoi_error_estimates[qoi_index];
2195 }
2196 
2197 
2198 
2199 Number System::point_value(unsigned int var,
2200  const Point & p,
2201  const bool insist_on_success,
2202  const NumericVector<Number> *sol) const
2203 {
2204  // This function must be called on every processor; there's no
2205  // telling where in the partition p falls.
2206  parallel_object_only();
2207 
2208  // And every processor had better agree about which point we're
2209  // looking for
2210 #ifndef NDEBUG
2211  libmesh_assert(this->comm().verify(p(0)));
2212 #if LIBMESH_DIM > 1
2213  libmesh_assert(this->comm().verify(p(1)));
2214 #endif
2215 #if LIBMESH_DIM > 2
2216  libmesh_assert(this->comm().verify(p(2)));
2217 #endif
2218 #endif // NDEBUG
2219 
2220  // Get a reference to the mesh object associated with the system object that calls this function
2221  const MeshBase & mesh = this->get_mesh();
2222 
2223  // Use an existing PointLocator or create a new one
2224  std::unique_ptr<PointLocatorBase> locator_ptr = mesh.sub_point_locator();
2225  PointLocatorBase & locator = *locator_ptr;
2226 
2227  if (!insist_on_success || !mesh.is_serial())
2228  locator.enable_out_of_mesh_mode();
2229 
2230  // Get a pointer to an element that contains p and allows us to
2231  // evaluate var
2232  const std::set<subdomain_id_type> & raw_subdomains =
2233  this->variable(var).active_subdomains();
2234  const std::set<subdomain_id_type> * implicit_subdomains =
2235  raw_subdomains.empty() ? nullptr : &raw_subdomains;
2236  const Elem * e = locator(p, implicit_subdomains);
2237 
2238  Number u = 0;
2239 
2240  if (e && this->get_dof_map().is_evaluable(*e, var))
2241  u = point_value(var, p, *e, sol);
2242 
2243  // If I have an element containing p, then let's let everyone know
2244  processor_id_type lowest_owner =
2245  (e && (e->processor_id() == this->processor_id())) ?
2246  this->processor_id() : this->n_processors();
2247  this->comm().min(lowest_owner);
2248 
2249  // Everybody should get their value from a processor that was able
2250  // to compute it.
2251  // If nobody admits owning the point, we have a problem.
2252  if (lowest_owner != this->n_processors())
2253  this->comm().broadcast(u, lowest_owner);
2254  else
2255  libmesh_assert(!insist_on_success);
2256 
2257  return u;
2258 }
2259 
2260 Number System::point_value(unsigned int var,
2261  const Point & p,
2262  const Elem & e,
2263  const NumericVector<Number> *sol) const
2264 {
2265  // Ensuring that the given point is really in the element is an
2266  // expensive assert, but as long as debugging is turned on we might
2267  // as well try to catch a particularly nasty potential error
2269 
2270  if (!sol)
2271  sol = this->current_local_solution.get();
2272 
2273  // Get the dof map to get the proper indices for our computation
2274  const DofMap & dof_map = this->get_dof_map();
2275 
2276  // Make sure we can evaluate on this element.
2277  libmesh_assert (dof_map.is_evaluable(e, var));
2278 
2279  // Need dof_indices for phi[i][j]
2280  std::vector<dof_id_type> dof_indices;
2281 
2282  // Fill in the dof_indices for our element
2283  dof_map.dof_indices (&e, dof_indices, var);
2284 
2285  // Get the no of dofs associated with this point
2286  const unsigned int num_dofs = cast_int<unsigned int>
2287  (dof_indices.size());
2288 
2289  FEType fe_type = dof_map.variable_type(var);
2290 
2291  // Map the physical co-ordinates to the master co-ordinates
2292  Point coor = FEMap::inverse_map(e.dim(), &e, p);
2293 
2294  // get the shape function value via the FEInterface to also handle the case
2295  // of infinite elements correctly, the shape function is not fe->phi().
2296  FEComputeData fe_data(this->get_equation_systems(), coor);
2297  FEInterface::compute_data(e.dim(), fe_type, &e, fe_data);
2298 
2299  // Get ready to accumulate a value
2300  Number u = 0;
2301 
2302  for (unsigned int l=0; l<num_dofs; l++)
2303  {
2304  u += fe_data.shape[l] * (*sol)(dof_indices[l]);
2305  }
2306 
2307  return u;
2308 }
2309 
2310 
2311 
2312 Number System::point_value(unsigned int var, const Point & p, const Elem * e) const
2313 {
2314  libmesh_assert(e);
2315  return this->point_value(var, p, *e);
2316 }
2317 
2318 
2319 
2320 Number System::point_value(unsigned int var, const Point & p, const NumericVector<Number> * sol) const
2321 {
2322  return this->point_value(var, p, true, sol);
2323 }
2324 
2325 
2326 
2327 
2329  const Point & p,
2330  const bool insist_on_success,
2331  const NumericVector<Number> *sol) const
2332 {
2333  // This function must be called on every processor; there's no
2334  // telling where in the partition p falls.
2335  parallel_object_only();
2336 
2337  // And every processor had better agree about which point we're
2338  // looking for
2339 #ifndef NDEBUG
2340  libmesh_assert(this->comm().verify(p(0)));
2341 #if LIBMESH_DIM > 1
2342  libmesh_assert(this->comm().verify(p(1)));
2343 #endif
2344 #if LIBMESH_DIM > 2
2345  libmesh_assert(this->comm().verify(p(2)));
2346 #endif
2347 #endif // NDEBUG
2348 
2349  // Get a reference to the mesh object associated with the system object that calls this function
2350  const MeshBase & mesh = this->get_mesh();
2351 
2352  // Use an existing PointLocator or create a new one
2353  std::unique_ptr<PointLocatorBase> locator_ptr = mesh.sub_point_locator();
2354  PointLocatorBase & locator = *locator_ptr;
2355 
2356  if (!insist_on_success || !mesh.is_serial())
2357  locator.enable_out_of_mesh_mode();
2358 
2359  // Get a pointer to an element that contains p and allows us to
2360  // evaluate var
2361  const std::set<subdomain_id_type> & raw_subdomains =
2362  this->variable(var).active_subdomains();
2363  const std::set<subdomain_id_type> * implicit_subdomains =
2364  raw_subdomains.empty() ? nullptr : &raw_subdomains;
2365  const Elem * e = locator(p, implicit_subdomains);
2366 
2367  Gradient grad_u;
2368 
2369  if (e && this->get_dof_map().is_evaluable(*e, var))
2370  grad_u = point_gradient(var, p, *e, sol);
2371 
2372  // If I have an element containing p, then let's let everyone know
2373  processor_id_type lowest_owner =
2374  (e && (e->processor_id() == this->processor_id())) ?
2375  this->processor_id() : this->n_processors();
2376  this->comm().min(lowest_owner);
2377 
2378  // Everybody should get their value from a processor that was able
2379  // to compute it.
2380  // If nobody admits owning the point, we may have a problem.
2381  if (lowest_owner != this->n_processors())
2382  this->comm().broadcast(grad_u, lowest_owner);
2383  else
2384  libmesh_assert(!insist_on_success);
2385 
2386  return grad_u;
2387 }
2388 
2389 
2391  const Point & p,
2392  const Elem & e,
2393  const NumericVector<Number> *sol) const
2394 {
2395  // Ensuring that the given point is really in the element is an
2396  // expensive assert, but as long as debugging is turned on we might
2397  // as well try to catch a particularly nasty potential error
2399 
2400  if (!sol)
2401  sol = this->current_local_solution.get();
2402 
2403  // Get the dof map to get the proper indices for our computation
2404  const DofMap & dof_map = this->get_dof_map();
2405 
2406  // write the element dimension into a separate variable.
2407  const unsigned int dim = e.dim();
2408 
2409  // Make sure we can evaluate on this element.
2410  libmesh_assert (dof_map.is_evaluable(e, var));
2411 
2412  // Need dof_indices for phi[i][j]
2413  std::vector<dof_id_type> dof_indices;
2414 
2415  // Fill in the dof_indices for our element
2416  dof_map.dof_indices (&e, dof_indices, var);
2417 
2418  // Get the no of dofs associated with this point
2419  const unsigned int num_dofs = cast_int<unsigned int>
2420  (dof_indices.size());
2421 
2422  FEType fe_type = dof_map.variable_type(var);
2423 
2424  // Map the physical co-ordinates to the master co-ordinates
2425  Point coor = FEMap::inverse_map(dim, &e, p);
2426 
2427  // get the shape function value via the FEInterface to also handle the case
2428  // of infinite elements correctly, the shape function is not fe->phi().
2429  FEComputeData fe_data(this->get_equation_systems(), coor);
2430  fe_data.enable_derivative();
2431  FEInterface::compute_data(dim, fe_type, &e, fe_data);
2432 
2433  // Get ready to accumulate a gradient
2434  Gradient grad_u;
2435 
2436  for (unsigned int l=0; l<num_dofs; l++)
2437  {
2438  // Chartesian coordinates have always LIBMESH_DIM entries,
2439  // local coordinates have as many coordinates as the element has.
2440  for (std::size_t v=0; v<dim; v++)
2441  for (std::size_t xyz=0; xyz<LIBMESH_DIM; xyz++)
2442  {
2443  // FIXME: this needs better syntax: It is matrix-vector multiplication.
2444  grad_u(xyz) += fe_data.local_transform[v][xyz]
2445  * fe_data.dshape[l](v)
2446  * (*sol)(dof_indices[l]);
2447  }
2448  }
2449 
2450  return grad_u;
2451 }
2452 
2453 
2454 
2455 Gradient System::point_gradient(unsigned int var, const Point & p, const Elem * e) const
2456 {
2457  libmesh_assert(e);
2458  return this->point_gradient(var, p, *e);
2459 }
2460 
2461 
2462 
2463 Gradient System::point_gradient(unsigned int var, const Point & p, const NumericVector<Number> * sol) const
2464 {
2465  return this->point_gradient(var, p, true, sol);
2466 }
2467 
2468 
2469 
2470 // We can only accumulate a hessian with --enable-second
2471 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
2473  const Point & p,
2474  const bool insist_on_success,
2475  const NumericVector<Number> *sol) const
2476 {
2477  // This function must be called on every processor; there's no
2478  // telling where in the partition p falls.
2479  parallel_object_only();
2480 
2481  // And every processor had better agree about which point we're
2482  // looking for
2483 #ifndef NDEBUG
2484  libmesh_assert(this->comm().verify(p(0)));
2485 #if LIBMESH_DIM > 1
2486  libmesh_assert(this->comm().verify(p(1)));
2487 #endif
2488 #if LIBMESH_DIM > 2
2489  libmesh_assert(this->comm().verify(p(2)));
2490 #endif
2491 #endif // NDEBUG
2492 
2493  // Get a reference to the mesh object associated with the system object that calls this function
2494  const MeshBase & mesh = this->get_mesh();
2495 
2496  // Use an existing PointLocator or create a new one
2497  std::unique_ptr<PointLocatorBase> locator_ptr = mesh.sub_point_locator();
2498  PointLocatorBase & locator = *locator_ptr;
2499 
2500  if (!insist_on_success || !mesh.is_serial())
2501  locator.enable_out_of_mesh_mode();
2502 
2503  // Get a pointer to an element that contains p and allows us to
2504  // evaluate var
2505  const std::set<subdomain_id_type> & raw_subdomains =
2506  this->variable(var).active_subdomains();
2507  const std::set<subdomain_id_type> * implicit_subdomains =
2508  raw_subdomains.empty() ? nullptr : &raw_subdomains;
2509  const Elem * e = locator(p, implicit_subdomains);
2510 
2511  Tensor hess_u;
2512 
2513  if (e && this->get_dof_map().is_evaluable(*e, var))
2514  hess_u = point_hessian(var, p, *e, sol);
2515 
2516  // If I have an element containing p, then let's let everyone know
2517  processor_id_type lowest_owner =
2518  (e && (e->processor_id() == this->processor_id())) ?
2519  this->processor_id() : this->n_processors();
2520  this->comm().min(lowest_owner);
2521 
2522  // Everybody should get their value from a processor that was able
2523  // to compute it.
2524  // If nobody admits owning the point, we may have a problem.
2525  if (lowest_owner != this->n_processors())
2526  this->comm().broadcast(hess_u, lowest_owner);
2527  else
2528  libmesh_assert(!insist_on_success);
2529 
2530  return hess_u;
2531 }
2532 
2534  const Point & p,
2535  const Elem & e,
2536  const NumericVector<Number> *sol) const
2537 {
2538  // Ensuring that the given point is really in the element is an
2539  // expensive assert, but as long as debugging is turned on we might
2540  // as well try to catch a particularly nasty potential error
2542 
2543  if (!sol)
2544  sol = this->current_local_solution.get();
2545 
2546  if (e.infinite())
2547  libmesh_not_implemented();
2548 
2549  // Get the dof map to get the proper indices for our computation
2550  const DofMap & dof_map = this->get_dof_map();
2551 
2552  // Make sure we can evaluate on this element.
2553  libmesh_assert (dof_map.is_evaluable(e, var));
2554 
2555  // Need dof_indices for phi[i][j]
2556  std::vector<dof_id_type> dof_indices;
2557 
2558  // Fill in the dof_indices for our element
2559  dof_map.dof_indices (&e, dof_indices, var);
2560 
2561  // Get the no of dofs associated with this point
2562  const unsigned int num_dofs = cast_int<unsigned int>
2563  (dof_indices.size());
2564 
2565  FEType fe_type = dof_map.variable_type(var);
2566 
2567  // Build a FE again so we can calculate u(p)
2568  std::unique_ptr<FEBase> fe (FEBase::build(e.dim(), fe_type));
2569 
2570  // Map the physical co-ordinates to the master co-ordinates
2571  // Build a vector of point co-ordinates to send to reinit
2572  std::vector<Point> coor(1, FEMap::inverse_map(e.dim(), &e, p));
2573 
2574  // Get the values of the shape function derivatives
2575  const std::vector<std::vector<RealTensor>> & d2phi = fe->get_d2phi();
2576 
2577  // Reinitialize the element and compute the shape function values at coor
2578  fe->reinit (&e, &coor);
2579 
2580  // Get ready to accumulate a hessian
2581  Tensor hess_u;
2582 
2583  for (unsigned int l=0; l<num_dofs; l++)
2584  {
2585  hess_u.add_scaled (d2phi[l][0], (*sol)(dof_indices[l]));
2586  }
2587 
2588  return hess_u;
2589 }
2590 
2591 
2592 
2593 Tensor System::point_hessian(unsigned int var, const Point & p, const Elem * e) const
2594 {
2595  libmesh_assert(e);
2596  return this->point_hessian(var, p, *e);
2597 }
2598 
2599 
2600 
2601 Tensor System::point_hessian(unsigned int var, const Point & p, const NumericVector<Number> * sol) const
2602 {
2603  return this->point_hessian(var, p, true, sol);
2604 }
2605 
2606 #else
2607 
2608 Tensor System::point_hessian(unsigned int, const Point &, const bool,
2609  const NumericVector<Number> *) const
2610 {
2611  libmesh_error_msg("We can only accumulate a hessian with --enable-second");
2612 
2613  // Avoid compiler warnings
2614  return Tensor();
2615 }
2616 
2617 Tensor System::point_hessian(unsigned int, const Point &, const Elem &,
2618  const NumericVector<Number> *) const
2619 {
2620  libmesh_error_msg("We can only accumulate a hessian with --enable-second");
2621 
2622  // Avoid compiler warnings
2623  return Tensor();
2624 }
2625 
2626 Tensor System::point_hessian(unsigned int, const Point &, const Elem *) const
2627 {
2628  libmesh_error_msg("We can only accumulate a hessian with --enable-second");
2629 
2630  // Avoid compiler warnings
2631  return Tensor();
2632 }
2633 
2634 Tensor System::point_hessian(unsigned int, const Point &, const NumericVector<Number> *) const
2635 {
2636  libmesh_error_msg("We can only accumulate a hessian with --enable-second");
2637 
2638  // Avoid compiler warnings
2639  return Tensor();
2640 }
2641 
2642 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
2643 
2645 {
2646  this->get_dof_map().create_static_condensation(this->get_mesh(), *this);
2647 }
2648 
2650 {
2651  return this->get_dof_map().has_static_condensation();
2652 }
2653 
2654 unsigned int System::n_vars() const
2655 {
2656  return this->get_dof_map().n_vars();
2657 }
2658 
2659 const std::string & System::variable_name (const unsigned int i) const
2660 {
2661  return this->get_dof_map().variable_name(i);
2662 }
2663 
2665 {
2666  return this->get_dof_map().identify_variable_groups();
2667 }
2668 
2670 {
2672 }
2673 
2674 unsigned int System::n_components() const
2675 {
2676  return this->get_dof_map().n_components(this->get_mesh());
2677 }
2678 
2679 unsigned int System::n_variable_groups() const
2680 {
2681  return this->get_dof_map().n_variable_groups();
2682 }
2683 
2684 const Variable & System::variable (const unsigned int i) const
2685 {
2686  return this->get_dof_map().variable(i);
2687 }
2688 
2689 const VariableGroup & System::variable_group (const unsigned int vg) const
2690 {
2691  return this->get_dof_map().variable_group(vg);
2692 }
2693 
2694 unsigned int
2695 System::variable_scalar_number (unsigned int var_num,
2696  unsigned int component) const
2697 {
2698  return this->get_dof_map().variable_scalar_number(var_num, component);
2699 }
2700 
2701 const FEType & System::variable_type (const unsigned int i) const
2702 {
2703  return this->get_dof_map().variable_type(i);
2704 }
2705 
2706 const FEType & System::variable_type (std::string_view var) const
2707 {
2708  return this->get_dof_map().variable_type(var);
2709 }
2710 
2711 } // namespace libMesh
void set_vector_as_adjoint(const std::string &vec_name, int qoi_num)
Allows one to set the QoI index controlling whether the vector identified by vec_name represents a so...
Definition: system.C:1145
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
unsigned int add_variables(const std::vector< std::string > &vars, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variables vars to the list of variables for this system.
Definition: system.C:1363
vectors_iterator vectors_end()
End of vectors container.
Definition: system.h:2475
virtual bool initialized() const
bool is_initialized() const
Definition: system.h:2415
FEFamily family
The type of finite element.
Definition: fe_type.h:221
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Definition: system.h:1638
void local_dof_indices(const unsigned int var, std::set< dof_id_type > &var_indices) const
Fills the std::set with the degrees of freedom on the local processor corresponding the the variable ...
Definition: system.C:1405
bool closed()
Checks that the library has been closed.
Definition: libmesh.C:341
unsigned int n_variable_groups() const
Definition: dof_map.h:733
dof_id_type end_dof(const processor_id_type proc) const
Definition: dof_map_base.h:191
virtual Real subset_l2_norm(const std::set< numeric_index_type > &indices) const
This is the EquationSystems class.
unsigned int variable_scalar_number(std::string_view var, unsigned int component) const
Definition: system.h:2432
virtual void clear()
Clear all the data structures associated with the system.
Definition: system.C:174
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
virtual void create_static_condensation()
Request that static condensation be performed for this system.
Definition: system.C:2644
void update_global_solution(std::vector< Number > &global_soln) const
Fill the input vector global_soln so that it contains the global solution on all processors.
Definition: system.C:730
Abstract base class to be used for system initialization.
Definition: system.h:132
std::size_t size() const
const Variable & variable(unsigned int var) const
Return a constant reference to Variable var.
Definition: system.C:2684
Assembly * _assemble_system_object
Object that assembles the system.
Definition: system.h:2159
Gradient point_gradient(unsigned int var, const Point &p, const bool insist_on_success=true, const NumericVector< Number > *sol=nullptr) const
Definition: system.C:2328
bool _basic_system_only
Holds true if the components of more advanced system types (e.g.
Definition: system.h:2277
virtual void init_matrices()
Initializes the matrices associated with this system.
Definition: system.C:308
bool _is_initialized
true when additional vectors and variables do not require immediate initialization, false otherwise.
Definition: system.h:2283
void(* _constrain_system_function)(EquationSystems &es, const std::string &name)
Function to impose constraints.
Definition: system.h:2164
NumericVector< Number > & add_adjoint_solution(unsigned int i=0)
Definition: system.C:1218
virtual void forward_qoi_parameter_sensitivity(const QoISet &qoi_indices, const ParameterVector &parameters, SensitivityData &sensitivities)
Solves for parameter sensitivities using the forward method.
Definition: system.h:2571
virtual void get(const std::vector< numeric_index_type > &index, T *values) const
Access multiple components at once.
Data structure for specifying which Parameters should be independent variables in a parameter sensiti...
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Definition: dof_map.C:2197
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
Number get_qoi_value(unsigned int qoi_index) const
Definition: system.C:2164
virtual void reinit()
Reinitializes degrees of freedom and other required data on the current mesh.
Definition: system.C:439
unsigned int dim
bool has_variable(std::string_view var) const
Definition: dof_map.h:2975
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true, const bool extra_checks=true)
Definition: fe_map.C:1628
class FEComputeData hides arbitrary data to be passed to and from children of FEBase through the FEIn...
unsigned int n_qois() const
Number of currently active quantities of interest.
Definition: system.h:2520
unsigned int n_components() const
Definition: system.C:2674
virtual bool initialized() const
virtual void initialize()=0
Initialization function.
virtual numeric_index_type size() const =0
NumericVector< Number > & get_sensitivity_solution(unsigned int i=0)
Definition: system.C:1177
OrderWrapper radial_order
The approximation order in radial direction of the infinite element.
Definition: fe_type.h:254
unsigned int n_components(const MeshBase &mesh) const
Definition: dof_map.h:2952
virtual void assemble_qoi_derivative(const QoISet &qoi_indices=QoISet(), bool include_liftfunc=true, bool apply_constraints=true)
Calls user qoi derivative function.
Definition: system.C:573
NumericVector< Number > & get_sensitivity_rhs(unsigned int i=0)
Definition: system.C:1322
int vector_is_adjoint(std::string_view vec_name) const
Definition: system.C:1158
virtual void assemble()=0
Assembly function.
unsigned int add_variables(System &sys, const std::vector< std::string > &vars, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variables vars to the list of variables for this system.
Definition: dof_map.C:3249
Constraint * _constrain_system_object
Object that constrains the system.
Definition: system.h:2170
const EquationSystems & get_equation_systems() const
Definition: system.h:722
void sum(T &r) const
unsigned int n_variable_groups() const
Definition: system.C:2679
bool is_attached(SparseMatrix< Number > &matrix)
Matrices should not be attached more than once.
Definition: dof_map.C:295
virtual void assemble()
Prepares matrix and _dof_map for matrix assembly.
Definition: system.C:551
void attach_matrix(SparseMatrix< Number > &matrix)
Additional matrices may be attached to this DofMap.
Definition: dof_map.C:240
unsigned int add_variable_array(const std::vector< std::string > &vars, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds variables vars to the list of variables for this system.
Definition: system.C:1382
virtual void init_data()
Initializes the data for the system.
Definition: system.C:204
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
MeshBase & mesh
static FEFieldType field_type(const FEType &fe_type)
FEMNormType
defines an enum for norms defined on vectors of finite element coefficients
void get_all_variable_numbers(std::vector< unsigned int > &all_variable_numbers) const
Fills all_variable_numbers with all the variable numbers for the variables that have been added to th...
Definition: dof_map.C:3385
void init_qois(unsigned int n_qois)
Accessors for qoi and qoi_error_estimates vectors.
Definition: system.C:2149
Number point_value(unsigned int var, const Point &p, const bool insist_on_success=true, const NumericVector< Number > *sol=nullptr) const
Definition: system.C:2199
std::unique_ptr< DofMap > _dof_map
Data structure describing the relationship between nodes, variables, etc...
Definition: system.h:2202
virtual void user_initialization()
Calls user&#39;s attached initialization function, or is overridden by the user in derived classes...
Definition: system.C:2075
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
Definition: vector_fe_ex5.C:44
std::vector< Gradient > dshape
Storage for the computed shape derivative values.
const Parallel::Communicator & comm() const
std::map< std::string, std::unique_ptr< NumericVector< Number > >, std::less<> > _vectors
Some systems need an arbitrary number of vectors.
Definition: system.h:2237
NumericVector< Number > & add_vector(std::string_view vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Adds the additional vector vec_name to this system.
Definition: system.C:753
void get_all_variable_numbers(std::vector< unsigned int > &all_variable_numbers) const
Fills all_variable_numbers with all the variable numbers for the variables that have been added to th...
Definition: system.C:1399
OrderWrapper order
The approximation order of the element.
Definition: fe_type.h:215
This class defines a norm/seminorm to be applied to a NumericVector which contains coefficients in a ...
Definition: system_norm.h:49
const std::string & name(unsigned int v) const
Definition: variable.h:295
NumericVector< Number > & add_sensitivity_rhs(unsigned int i=0)
Definition: system.C:1312
void reinit_static_condensation()
Calls reinit on the static condensation map if it exists.
Definition: dof_map.C:3136
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
std::vector< std::vector< Real > > local_transform
Storage for local to global mapping at p.
bool has_static_condensation() const
Definition: system.C:2649
The libMesh namespace provides an interface to certain functionality in the library.
vectors_iterator vectors_begin()
Beginning of vectors container.
Definition: system.h:2463
virtual Real subset_linfty_norm(const std::set< numeric_index_type > &indices) const
void attach_QOI_derivative_object(QOIDerivative &qoi_derivative)
Register a user object for evaluating derivatives of a quantity of interest with respect to test func...
Definition: system.C:2061
virtual void adjoint_qoi_parameter_sensitivity(const QoISet &qoi_indices, const ParameterVector &parameters, SensitivityData &sensitivities)
Solves for parameter sensitivities using the adjoint method.
Definition: system.h:2562
dof_id_type n_local_dofs() const
Definition: system.C:156
Abstract base class to be used for system assembly.
Definition: system.h:156
NumericVector< Number > & add_weighted_sensitivity_solution()
Definition: system.C:1197
void(* _qoi_evaluate_function)(EquationSystems &es, const std::string &name, const QoISet &qoi_indices)
Function to evaluate quantity of interest.
Definition: system.h:2175
const SparseMatrix< Number > * request_matrix(std::string_view mat_name) const
Definition: system.C:1085
void init()
Initializes degrees of freedom on the current mesh.
Definition: system.C:193
const MeshBase & get_mesh() const
Definition: system.h:2359
NumericVector< Number > & get_weighted_sensitivity_solution()
Definition: system.C:1204
static std::unique_ptr< SparseMatrix< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package(), const MatrixBuildType matrix_build_type=MatrixBuildType::AUTOMATIC)
Builds a SparseMatrix<T> using the linear solver package specified by solver_package.
virtual Real subset_l1_norm(const std::set< numeric_index_type > &indices) const
std::string get_info() const
Gets summary info about the sparsity bandwidth and constraints.
Definition: dof_map.C:2981
dof_id_type n_dofs() const
Definition: system.C:119
uint8_t processor_id_type
Definition: id_types.h:104
This is the MeshBase class.
Definition: mesh_base.h:75
virtual void prolong_vectors()
Prolong vectors after the mesh has refined.
Definition: system.C:429
const Variable & variable(const unsigned int c) const override
Definition: dof_map.h:2356
virtual void user_QOI(const QoISet &qoi_indices)
Calls user&#39;s attached quantity of interest function, or is overridden by the user in derived classes...
Definition: system.C:2117
Number current_solution(const dof_id_type global_dof_number) const
Definition: system.C:163
Number fptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
Definition: projection.C:81
unsigned int variable_number(std::string_view var) const
Definition: system.C:1394
std::map< std::string, std::unique_ptr< SparseMatrix< Number > >, std::less<> > _matrices
Some systems need an arbitrary number of matrices.
Definition: system.h:2254
QOIDerivative * _qoi_evaluate_derivative_object
Object to compute derivatives of quantities of interest.
Definition: system.h:2196
NumericVector< Number > & add_weighted_sensitivity_adjoint_solution(unsigned int i=0)
Definition: system.C:1250
virtual void qoi_derivative(const QoISet &qoi_indices, bool include_liftfunc, bool apply_constraints)=0
Quantity of interest derivative function.
virtual void qoi_parameter_sensitivity(const QoISet &qoi_indices, const ParameterVector &parameters, SensitivityData &sensitivities)
Solves for the derivative of each of the system&#39;s quantities of interest q in qoi[qoi_indices] with r...
Definition: system.C:587
void attach_constraint_object(Constraint &constrain)
Register a user object for imposing constraints.
Definition: system.C:1987
SolverPackage default_solver_package()
Definition: libmesh.C:1188
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
std::unique_ptr< QBase > default_quadrature_rule(const unsigned int dim, const int extraorder=0) const
Definition: fe_type.C:34
void attach_init_object(Initialization &init)
Register a user class to use to initialize the system.
Definition: system.C:1925
unsigned int variable_scalar_number(unsigned int var_num, unsigned int component) const
Definition: dof_map.h:2963
bool has_variable(std::string_view var) const
Definition: system.C:1389
processor_id_type n_processors() const
void libmesh_ignore(const Args &...)
bool identify_variable_groups() const
Definition: dof_map.h:2940
unsigned int number() const
Definition: system.h:2351
void remove_vector(std::string_view vec_name)
Removes the additional vector vec_name from this system.
Definition: system.C:858
virtual bool contains_point(const Point &p, Real tol=TOLERANCE) const
Definition: elem.C:2755
This class defines the notion of a variable in the system.
Definition: variable.h:50
const std::set< subdomain_id_type > & active_subdomains() const
Definition: variable.h:181
bool has_static_condensation() const
Checks whether we have static condensation.
Definition: dof_map.h:1796
virtual Real l2_norm() const =0
std::map< std::string, int, std::less<> > _vector_is_adjoint
Holds non-negative if a vector by that name should be projected using adjoint constraints/BCs, -1 if primal.
Definition: system.h:2249
unsigned int add_variable(System &sys, std::string_view var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
Definition: dof_map.C:3142
Constraint & get_constraint_object()
Return the user object for imposing constraints.
Definition: system.C:2004
void min(const T &r, T &o, Request &req) const
Data structure for holding completed parameter sensitivity calculations.
bool _is_initialized
Flag that tells if init() has been called.
Definition: libmesh.C:315
StaticCondensationDofMap & get_static_condensation()
Definition: dof_map.h:2855
unsigned int n_variables() const
Definition: variable.h:266
std::vector< Number > qoi
Values of the quantities of interest.
Definition: system.h:1654
void remove_matrix(std::string_view mat_name)
Removes the additional matrix mat_name from this system.
Definition: system.C:1074
This is a base class for classes which represent subsets of the dofs of a System. ...
Definition: system_subset.h:42
unsigned int n_vectors() const
Definition: system.h:2457
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:97
void zero_variable(NumericVector< Number > &v, unsigned int var_num) const
Zeroes all dofs in v that correspond to variable number var_num.
Definition: system.C:1446
bool _prefer_hash_table_matrix_assembly
Whether to use hash table matrix assembly if the matrix sub-classes support it.
Definition: system.h:2325
virtual void add_matrices()
Insertion point for adding matrices in derived classes before init_matrices() is called.
Definition: system.h:1987
static void compute_data(const unsigned int dim, const FEType &fe_t, const Elem *elem, FEComputeData &data)
Lets the appropriate child of FEBase compute the requested data for the input specified in data...
bool has_constraint_object() const
Definition: system.C:1999
std::string get_info() const
Definition: system.C:1818
virtual void qoi(const QoISet &qoi_indices)=0
Quantity of interest function.
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1616
void attach_QOI_derivative(void fptr(EquationSystems &es, const std::string &name, const QoISet &qoi_indices, bool include_liftfunc, bool apply_constraints))
Register a user function for evaluating derivatives of a quantity of interest with respect to test fu...
Definition: system.C:2044
Real calculate_norm(const NumericVector< Number > &v, unsigned int var, FEMNormType norm_type, std::set< unsigned int > *skip_dimensions=nullptr) const
Definition: system.C:1502
void print_dof_constraints(std::ostream &os=libMesh::out, bool print_nonlocal=false) const
Prints (from processor 0) all DoF and Node constraints.
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
libmesh_assert(ctx)
QOI * _qoi_evaluate_object
Object to compute quantities of interest.
Definition: system.h:2182
Initialization * _init_system_object
Object that initializes the system.
Definition: system.h:2148
const FEType & variable_type(const unsigned int i) const
Definition: dof_map.h:2386
void attach_QOI_function(void fptr(EquationSystems &es, const std::string &name, const QoISet &qoi_indices))
Register a user function for evaluating the quantities of interest, whose values should be placed in ...
Definition: system.C:2012
const VariableGroup & variable_group(const unsigned int c) const
Definition: dof_map.h:2346
unsigned int add_variable(std::string_view var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
Definition: system.C:1342
const std::string & variable_name(const unsigned int i) const
Definition: system.C:2659
InfMapType inf_map
The coordinate mapping type of the infinite element.
Definition: fe_type.h:275
void(* _init_system_function)(EquationSystems &es, const std::string &name)
Function that initializes the system.
Definition: system.h:2142
const std::set< unsigned char > & elem_dimensions() const
Definition: mesh_base.h:282
std::map< std::string, ParallelType, std::less<> > _matrix_types
Holds the types of the matrices.
Definition: system.h:2259
This is the base class for point locators.
virtual void restrict_vectors()
Restrict vectors after the mesh has coarsened.
Definition: system.C:371
const std::string & variable_name(const unsigned int i) const
Definition: dof_map.h:2932
void create_static_condensation(MeshBase &mesh, System &system)
Add a static condensation class.
Definition: dof_map.C:3131
void attach_assemble_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in assembling the system matrix and RHS.
Definition: system.C:1939
auto norm(const T &a) -> decltype(std::abs(a))
Definition: tensor_tools.h:74
unsigned int n_matrices() const
Definition: system.h:2598
An object whose state is distributed along a set of processors.
This class defines a logically grouped set of variables in the system.
Definition: variable.h:203
std::size_t size(const System &sys) const
Definition: qoi_set.C:35
std::vector< Number > qoi_error_estimates
Vector to hold error estimates for qois, either from a steady state calculation, or from a single uns...
Definition: system.h:1662
void add_scaled(const TypeTensor< T2 > &, const T &)
Add a scaled tensor to this tensor without creating a temporary.
Definition: type_tensor.h:851
NumericVector< Number > & add_adjoint_rhs(unsigned int i=0)
Definition: system.C:1282
virtual Real l1_norm() const =0
virtual void restrict_solve_to(const SystemSubset *subset, const SubsetSolveMode subset_solve_mode=SUBSET_ZERO)
After calling this method, any solve will be restricted to the given subdomain.
Definition: system.C:542
virtual void user_QOI_derivative(const QoISet &qoi_indices=QoISet(), bool include_liftfunc=true, bool apply_constraints=true)
Calls user&#39;s attached quantity of interest derivative function, or is overridden by the user in deriv...
Definition: system.C:2131
void broadcast(T &data, const unsigned int root_id=0, const bool identical_sizes=false) const
virtual void reinit_constraints()
Reinitializes the constraints for this system.
Definition: system.C:480
bool identify_variable_groups() const
Definition: system.C:2664
virtual void close()=0
Calls the NumericVector&#39;s internal assembly routines, ensuring that the values are consistent across ...
virtual void reinit_mesh()
Reinitializes the system with a new mesh.
Definition: system.C:286
const NumericVector< Number > * request_vector(std::string_view vec_name) const
Definition: system.C:876
void attach_assemble_object(Assembly &assemble)
Register a user object to use in assembling the system matrix and RHS.
Definition: system.C:1956
FEFamily radial_family
The type of approximation in radial direction.
Definition: fe_type.h:267
bool _solution_projection
Holds true if the solution vector should be projected onto a changed grid, false if it should be zero...
Definition: system.h:2271
EquationSystems & _equation_systems
Constant reference to the EquationSystems object used for the simulation.
Definition: system.h:2208
std::string enum_to_string(const T e)
virtual bool closed() const
unsigned int n_vars() const
Definition: dof_map.h:2926
virtual std::string system_type() const
Definition: system.h:509
virtual void constrain()=0
Constraint function.
const std::string & vector_name(const unsigned int vec_num) const
Definition: system.C:968
void attach_QOI_object(QOI &qoi)
Register a user object for evaluating the quantities of interest, whose values should be placed in Sy...
Definition: system.C:2030
virtual void user_assembly()
Calls user&#39;s attached assembly function, or is overridden by the user in derived classes.
Definition: system.C:2089
ParallelType type() const
virtual void update()
Update the local values to reflect the solution on neighboring processors.
Definition: system.C:495
Real weight(unsigned int var) const
Definition: system_norm.C:134
const FEType & variable_type(const unsigned int i) const
Definition: system.C:2701
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void set_vector_preservation(const std::string &vec_name, bool preserve)
Allows one to set the boolean controlling whether the vector identified by vec_name should be "preser...
Definition: system.C:1123
void process_constraints(MeshBase &)
Postprocesses any constrained degrees of freedom to be constrained only in terms of unconstrained dof...
void set_qoi_error_estimate(unsigned int qoi_index, Number qoi_error_estimate)
Definition: system.C:2184
virtual void swap(NumericVector< T > &v)
Swaps the contents of this with v.
void max(const T &r, T &o, Request &req) const
NumberTensorValue Tensor
virtual unsigned short dim() const =0
bool _matrices_initialized
false when additional matrices being added require initialization, true otherwise.
Definition: system.h:2264
auto norm_sq(const T &a) -> decltype(std::norm(a))
Definition: tensor_tools.h:104
unsigned int add_variable_array(System &sys, const std::vector< std::string > &vars, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds variables vars to the list of variables for this system.
Definition: dof_map.C:3373
OStreamProxy out
bool _require_sparsity_pattern
Whether any of our matrices require an initial sparsity pattern computation in order to determine pre...
Definition: system.h:2330
const std::string _sys_name
A name associated with this system.
Definition: system.h:2219
virtual numeric_index_type local_size() const =0
Tensor point_hessian(unsigned int var, const Point &p, const bool insist_on_success=true, const NumericVector< Number > *sol=nullptr) const
Definition: system.C:2472
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, SolverPackage solver_package=libMesh::default_solver_package(), ParallelType parallel_type=AUTOMATIC)
Builds a NumericVector on the processors in communicator comm using the linear solver package specifi...
void clear_sparsity()
Clears the sparsity pattern.
Definition: dof_map.C:1977
NumericVector< Number > & get_adjoint_solution(unsigned int i=0)
Definition: system.C:1230
virtual void re_update()
Re-update the local values when the mesh has changed.
Definition: system.C:518
unsigned int variable_number(std::string_view var) const
Definition: dof_map.h:2980
Abstract base class to be used for system constraints.
Definition: system.h:180
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
Definition: system.h:1628
void project_vector(NumericVector< Number > &new_vector, FunctionBase< Number > *f, FunctionBase< Gradient > *g=nullptr, int is_adjoint=-1) const
Projects arbitrary functions onto a vector of degree of freedom values for the current system...
virtual bool compare(const System &other_system, const Real threshold, const bool verbose) const
Definition: system.C:603
void set_qoi(unsigned int qoi_index, Number qoi_value)
Definition: system.C:2156
dof_id_type n_local_constrained_dofs() const
Definition: system.C:141
void create_dof_constraints(const MeshBase &, Real time=0)
Rebuilds the raw degree of freedom and DofObject constraints, based on attached DirichletBoundary obj...
virtual void set(const numeric_index_type i, const T value)=0
Sets v(i) = value.
void prepare_send_list()
Takes the _send_list vector (which may have duplicate entries) and sorts it.
Definition: dof_map.C:1831
bool on_command_line(std::string arg)
Definition: libmesh.C:1058
virtual bool infinite() const =0
const std::string & name() const
Definition: system.h:2343
void attach_constraint_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function for imposing constraints.
Definition: system.C:1970
virtual void user_constrain()
Calls user&#39;s attached constraint function, or is overridden by the user in derived classes...
Definition: system.C:2103
Real discrete_var_norm(const NumericVector< Number > &v, unsigned int var, FEMNormType norm_type) const
Finds the discrete norm for the entries in the vector corresponding to Dofs associated with var...
Definition: system.C:1483
dof_id_type first_dof(const processor_id_type proc) const
Definition: dof_map_base.h:185
Number get_qoi_error_estimate_value(unsigned int qoi_index) const
Definition: system.C:2191
virtual void assemble_qoi(const QoISet &qoi_indices=QoISet())
Calls user qoi function.
Definition: system.C:562
NumericVector< Number > & add_sensitivity_solution(unsigned int i=0)
Definition: system.C:1167
FEFamily
defines an enum for finite element families.
unsigned int n_vars() const
Definition: system.C:2654
void(* _assemble_system_function)(EquationSystems &es, const std::string &name)
Function that assembles the system.
Definition: system.h:2153
MatrixBuildType
Defines an enum for matrix build types.
virtual ~System()
Definition: system.C:112
void compute_sparsity(const MeshBase &)
Computes the sparsity pattern for the matrices corresponding to proc_id and sends that data to Linear...
Definition: dof_map.C:1956
virtual void enable_out_of_mesh_mode()=0
Enables out-of-mesh mode.
processor_id_type processor_id() const
System(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
Definition: system.C:63
void enable_derivative()
Enable the computation of shape gradients (dshape).
bool vector_preservation(std::string_view vec_name) const
Definition: system.C:1133
SparseMatrix< Number > & add_matrix(std::string_view mat_name, ParallelType type=PARALLEL, MatrixBuildType mat_build_type=MatrixBuildType::AUTOMATIC)
Adds the additional matrix mat_name to this system.
Definition: system.C:995
SubsetSolveMode
defines an enum for the question what happens to the dofs outside the given subset when a system is s...
const DofMap & get_dof_map() const
Definition: system.h:2375
processor_id_type processor_id() const
Definition: dof_object.h:905
void late_matrix_init(SparseMatrix< Number > &mat, ParallelType type)
Helper function to keep DofMap forward declarable in system.h.
Definition: system.C:1061
const SparseMatrix< Number > & get_matrix(std::string_view mat_name) const
Definition: system.C:1109
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
MeshBase & _mesh
Constant reference to the mesh data structure used for the simulation.
Definition: system.h:2214
template class LIBMESH_EXPORT NumericVector< Number >
The QBase class provides the basic functionality from which various quadrature rules can be derived...
Definition: quadrature.h:61
virtual Real linfty_norm() const =0
const VariableGroup & variable_group(unsigned int vg) const
Return a constant reference to VariableGroup vg.
Definition: system.C:2689
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:533
void attach_init_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in initializing the system.
Definition: system.C:1908
void set_type(ParallelType t)
Allow the user to change the ParallelType of the NumericVector under some circumstances.
std::map< std::string, bool, std::less<> > _vector_projections
Holds true if a vector by that name should be projected onto a changed grid, false if it should be ze...
Definition: system.h:2243
This class forms the foundation from which generic finite elements may be derived.
Abstract base class to be used for quantities of interest.
Definition: system.h:204
bool is_evaluable(const DofObjectSubclass &obj, unsigned int var_num=libMesh::invalid_uint) const
Definition: dof_map.C:2668
const NumericVector< Number > & get_vector(std::string_view vec_name) const
Definition: system.C:928
NumericVector< Number > & get_weighted_sensitivity_adjoint_solution(unsigned int i=0)
Definition: system.C:1262
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
dof_id_type n_constrained_dofs() const
Definition: system.C:126
uint8_t dof_id_type
Definition: id_types.h:67
std::vector< Number > get_qoi_values() const
Returns a copy of qoi, not a reference.
Definition: system.C:2171
virtual void init(const numeric_index_type m, const numeric_index_type n, const numeric_index_type m_l, const numeric_index_type n_l, const numeric_index_type nnz=30, const numeric_index_type noz=10, const numeric_index_type blocksize=1)=0
Initialize SparseMatrix with the specified sizes.
const FEType & type() const
Definition: variable.h:144
void(* _qoi_evaluate_derivative_function)(EquationSystems &es, const std::string &name, const QoISet &qoi_indices, bool include_liftfunc, bool apply_constraints)
Function to evaluate quantity of interest derivative.
Definition: system.h:2187
ParallelType
Defines an enum for parallel data structure types.
Abstract base class to be used for derivatives of quantities of interest.
Definition: system.h:228
NumericVector< Number > & get_adjoint_rhs(unsigned int i=0)
Definition: system.C:1292
bool prefix_with_name() const
Definition: system.h:1955
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.
FEFieldType
defines an enum for finite element field types - i.e.