Loading [MathJax]/extensions/tex2jax.js
libMesh
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
equation_systems.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 // Local Includes
20 #include "libmesh/default_coupling.h" // For downconversion
21 #include "libmesh/dof_map.h"
22 #include "libmesh/eigen_system.h"
23 #include "libmesh/elem.h"
24 #include "libmesh/explicit_system.h"
25 #include "libmesh/fe_interface.h"
26 #include "libmesh/frequency_system.h"
27 #include "libmesh/int_range.h"
28 #include "libmesh/libmesh_logging.h"
29 #include "libmesh/linear_implicit_system.h"
30 #include "libmesh/mesh_base.h"
31 #include "libmesh/mesh_refinement.h"
32 #include "libmesh/newmark_system.h"
33 #include "libmesh/nonlinear_implicit_system.h"
34 #include "libmesh/parallel.h"
35 #include "libmesh/rb_construction.h"
36 #include "libmesh/remote_elem.h"
37 #include "libmesh/transient_rb_construction.h"
38 #include "libmesh/transient_system.h"
39 
40 // System includes
41 #include <functional> // std::plus
42 #include <numeric> // std::iota
43 #include <sstream>
44 
45 // Include the systems before this one to avoid
46 // overlapping forward declarations.
47 #include "libmesh/equation_systems.h"
48 
49 namespace libMesh
50 {
51 
53  ParallelObject (m),
54  _mesh (m),
55  _refine_in_reinit(true),
56  _enable_default_ghosting(true)
57 {
58  // Set default parameters
59  this->parameters.set<Real> ("linear solver tolerance") = TOLERANCE * TOLERANCE;
60  this->parameters.set<unsigned int>("linear solver maximum iterations") = 5000;
61 }
62 
63 
64 
66 
67 
68 
70 {
71  // Clear any additional parameters
72  parameters.clear ();
73 
74  // Clear the systems.
75  _systems.clear();
76 }
77 
78 
79 
81 {
82 #ifndef NDEBUG
83  for (auto i : make_range(this->n_systems()))
84  libmesh_assert(!this->get_system(i).is_initialized());
85 #endif
86 
87  this->reinit_mesh();
88 }
89 
90 
91 
93 {
94  const bool mesh_changed = this->reinit_solutions();
95 
96  // If the mesh has changed, systems will need to reinitialize their
97  // own data on the new mesh.
98  if (mesh_changed)
99  this->reinit_systems();
100 }
101 
103 {
104  const unsigned int n_sys = this->n_systems();
105 
106  libmesh_assert_not_equal_to (n_sys, 0);
107 
108  // Tell all the \p DofObject entities how many systems
109  // there are.
110  for (auto & node : _mesh.node_ptr_range())
111  node->set_n_systems(n_sys);
112 
113  for (auto & elem : _mesh.element_ptr_range())
114  elem->set_n_systems(n_sys);
115 
116  //for (auto i : make_range(this->n_systems()))
117  //this->get_system(i).init();
118 
119 #ifdef LIBMESH_ENABLE_AMR
120  MeshRefinement mesh_refine(_mesh);
121  mesh_refine.clean_refinement_flags();
122 #endif
123 
124  // Now loop over all the systems belonging to this ES
125  // and call reinit_mesh for each system
126  for (auto i : make_range(this->n_systems()))
127  this->get_system(i).reinit_mesh();
128 
129 }
130 
132 {
133  parallel_object_only();
134 
135  const unsigned int n_sys = this->n_systems();
136  libmesh_assert_not_equal_to (n_sys, 0);
137 
138  // And any new systems will need initialization
139  for (unsigned int i=0; i != n_sys; ++i)
140  if (!this->get_system(i).is_initialized())
141  this->get_system(i).init();
142 
143  // We used to assert that all nodes and elements *already* had
144  // n_systems() properly set; however this is false in the case where
145  // user code has manually added nodes and/or elements to an
146  // already-initialized system.
147 
148  // Make sure all the \p DofObject entities know how many systems
149  // there are.
150  {
151  // All the nodes
152  for (auto & node : _mesh.node_ptr_range())
153  node->set_n_systems(n_sys);
154 
155  // All the elements
156  for (auto & elem : _mesh.element_ptr_range())
157  elem->set_n_systems(n_sys);
158  }
159 
160  // Localize each system's vectors
161  for (unsigned int i=0; i != n_sys; ++i)
162  this->get_system(i).re_update();
163 
164 #ifdef LIBMESH_ENABLE_AMR
165 
166  bool mesh_changed = false;
167 
168  // FIXME: For backwards compatibility, assume
169  // refine_and_coarsen_elements or refine_uniformly have already
170  // been called
171  {
172  for (unsigned int i=0; i != n_sys; ++i)
173  {
174  System & sys = this->get_system(i);
175 
176  // Even if the system doesn't have any variables in it we want
177  // consistent behavior; e.g. distribute_dofs should have the
178  // opportunity to count up zero dofs on each processor.
179  //
180  // Who's been adding zero-var systems anyway, outside of my
181  // unit tests? - RHS
182  // if (!sys.n_vars())
183  // continue;
184 
186 
187  // Recreate any user or internal constraints
188  sys.reinit_constraints();
189 
191 
192  sys.prolong_vectors();
193  }
194  mesh_changed = true;
195  }
196 
197  if (this->_refine_in_reinit)
198  {
199  // Don't override any user refinement settings
200  MeshRefinement mesh_refine(_mesh);
201  mesh_refine.face_level_mismatch_limit() = 0; // unlimited
202  mesh_refine.overrefined_boundary_limit() = -1; // unlimited
203  mesh_refine.underrefined_boundary_limit() = -1; // unlimited
204 
205  // Try to coarsen the mesh, then restrict each system's vectors
206  // if necessary
207  if (mesh_refine.coarsen_elements())
208  {
209  for (auto i : make_range(this->n_systems()))
210  {
211  System & sys = this->get_system(i);
213  sys.reinit_constraints();
215  sys.restrict_vectors();
216  }
217  mesh_changed = true;
218  }
219 
220  // Once vectors are all restricted, we can delete
221  // children of coarsened elements
222  if (mesh_changed)
223  this->get_mesh().contract();
224 
225  // Try to refine the mesh, then prolong each system's vectors
226  // if necessary
227  if (mesh_refine.refine_elements())
228  {
229  for (auto i : make_range(this->n_systems()))
230  {
231  System & sys = this->get_system(i);
233  sys.reinit_constraints();
235  sys.prolong_vectors();
236  }
237  mesh_changed = true;
238  }
239  }
240 
241  return mesh_changed;
242 
243 #endif // #ifdef LIBMESH_ENABLE_AMR
244 
245  return false;
246 }
247 
248 
249 
251 {
252  for (auto i : make_range(this->n_systems()))
253  this->get_system(i).reinit();
254 }
255 
256 
257 
259 {
260  // A serial mesh means nothing needs to be done
261  if (_mesh.is_serial())
262  return;
263 
264  const unsigned int n_sys = this->n_systems();
265 
266  libmesh_assert_not_equal_to (n_sys, 0);
267 
268  // Gather the mesh
269  _mesh.allgather();
270 
271  // Tell all the \p DofObject entities how many systems
272  // there are.
273  for (auto & node : _mesh.node_ptr_range())
274  node->set_n_systems(n_sys);
275 
276  for (auto & elem : _mesh.element_ptr_range())
277  elem->set_n_systems(n_sys);
278 
279  // And distribute each system's dofs
280  for (auto i : make_range(this->n_systems()))
281  {
282  System & sys = this->get_system(i);
283  DofMap & dof_map = sys.get_dof_map();
284  dof_map.distribute_dofs(_mesh);
285 
286  // The user probably won't need constraint equations or the
287  // send_list after an allgather, but let's keep it in consistent
288  // shape just in case.
289  sys.reinit_constraints();
290  dof_map.prepare_send_list();
291  }
292 }
293 
294 
295 
297 {
298  _enable_default_ghosting = enable;
299  MeshBase &mesh = this->get_mesh();
300 
301  if (enable)
302  mesh.add_ghosting_functor(mesh.default_ghosting());
303  else
304  mesh.remove_ghosting_functor(mesh.default_ghosting());
305 
306  for (auto i : make_range(this->n_systems()))
307  {
308  DofMap & dof_map = this->get_system(i).get_dof_map();
309  if (enable)
310  dof_map.add_default_ghosting();
311  else
312  dof_map.remove_default_ghosting();
313  }
314 }
315 
316 
317 
319 {
320  LOG_SCOPE("update()", "EquationSystems");
321 
322  // Localize each system's vectors
323  for (auto i : make_range(this->n_systems()))
324  this->get_system(i).update();
325 }
326 
327 
328 
329 System & EquationSystems::add_system (std::string_view sys_type,
330  std::string_view name)
331 {
332  // If the user already built a system with this name, we'll
333  // trust them and we'll use it. That way they can pre-add
334  // non-standard derived system classes, and if their restart file
335  // has some non-standard sys_type we won't throw an error.
336  if (_systems.count(name))
337  {
338  return this->get_system(name);
339  }
340  // Build a basic System
341  else if (sys_type == "Basic")
342  this->add_system<System> (name);
343 
344  // Build a Newmark system
345  else if (sys_type == "Newmark")
346  this->add_system<NewmarkSystem> (name);
347 
348  // Build an Explicit system
349  else if ((sys_type == "Explicit"))
350  this->add_system<ExplicitSystem> (name);
351 
352  // Build an Implicit system
353  else if ((sys_type == "Implicit") ||
354  (sys_type == "Steady" ))
355  this->add_system<ImplicitSystem> (name);
356 
357  // build a transient implicit linear system
358  else if ((sys_type == "Transient") ||
359  (sys_type == "TransientImplicit") ||
360  (sys_type == "TransientLinearImplicit"))
361  this->add_system<TransientLinearImplicitSystem> (name);
362 
363  // build a transient implicit nonlinear system
364  else if (sys_type == "TransientNonlinearImplicit")
365  this->add_system<TransientNonlinearImplicitSystem> (name);
366 
367  // build a transient explicit system
368  else if (sys_type == "TransientExplicit")
369  this->add_system<TransientExplicitSystem> (name);
370 
371  // build a linear implicit system
372  else if (sys_type == "LinearImplicit")
373  this->add_system<LinearImplicitSystem> (name);
374 
375  // build a nonlinear implicit system
376  else if (sys_type == "NonlinearImplicit")
377  this->add_system<NonlinearImplicitSystem> (name);
378 
379  // build a Reduced Basis Construction system
380  else if (sys_type == "RBConstruction")
381  this->add_system<RBConstruction> (name);
382 
383  // build a transient Reduced Basis Construction system
384  else if (sys_type == "TransientRBConstruction")
385  this->add_system<TransientRBConstruction> (name);
386 
387 #ifdef LIBMESH_HAVE_SLEPC
388  // build an eigen system
389  else if (sys_type == "Eigen")
390  this->add_system<EigenSystem> (name);
391  else if (sys_type == "TransientEigenSystem")
392  this->add_system<TransientEigenSystem> (name);
393 #endif
394 
395 #if defined(LIBMESH_USE_COMPLEX_NUMBERS)
396  // build a frequency system
397  else if (sys_type == "Frequency")
398  this->add_system<FrequencySystem> (name);
399 #endif
400 
401  else
402  libmesh_error_msg("ERROR: Unknown system type: " << sys_type);
403 
404  // Return a reference to the new system
405  //return (*this)(name);
406  return this->get_system(name);
407 }
408 
409 
410 
412 {
413  libmesh_assert (this->n_systems());
414 
415  for (auto i : make_range(this->n_systems()))
416  this->get_system(i).solve();
417 }
418 
419 
420 
422 {
423  libmesh_assert (this->n_systems());
424 
425  for (auto i : make_range(this->n_systems()))
426  this->get_system(i).sensitivity_solve(parameters_in);
427 }
428 
429 
430 
431 void EquationSystems::adjoint_solve (const QoISet & qoi_indices)
432 {
433  libmesh_assert (this->n_systems());
434 
435  for (unsigned int i=this->n_systems(); i != 0; --i)
436  this->get_system(i-1).adjoint_solve(qoi_indices);
437 }
438 
439 
440 
441 void EquationSystems::build_variable_names (std::vector<std::string> & var_names,
442  const FEType * type,
443  const std::set<std::string> * system_names) const
444 {
445  // start indexing at end of possibly non-empty vector of variable names to avoid overwriting them
446  unsigned int var_num = var_names.size();
447 
448  // We'll want to double-check that we don't have any naming
449  // conflicts; this API causes problems down the line if so.
450  std::unordered_multiset<std::string> seen_names;
451 
452  // Need to size var_names by scalar variables plus all the
453  // vector components for all the vector variables
454  //Could this be replaced by a/some convenience methods?[PB]
455  {
456  unsigned int n_scalar_vars = 0;
457  unsigned int n_vector_vars = 0;
458 
459  for (const auto & [sys_name, sys_ptr] : _systems)
460  {
461  // Check current system is listed in system_names, and skip pos if not
462  bool use_current_system = (system_names == nullptr);
463  if (!use_current_system)
464  use_current_system = system_names->count(sys_name);
465  if (!use_current_system || sys_ptr->hide_output())
466  {
467  for (auto vn : make_range(sys_ptr->n_vars()))
468  seen_names.insert(sys_ptr->variable_name(vn));
469  continue;
470  }
471 
472  for (auto vn : make_range(sys_ptr->n_vars()))
473  {
474  seen_names.insert(sys_ptr->variable_name(vn));
475  if (FEInterface::field_type(sys_ptr->variable_type(vn)) == TYPE_VECTOR)
476  n_vector_vars++;
477  else
478  n_scalar_vars++;
479  }
480  }
481 
482  // Here, we're assuming the number of vector components is the same
483  // as the mesh spatial dimension.
484  unsigned int dim = this->get_mesh().spatial_dimension();
485  unsigned int nv = n_scalar_vars + dim*n_vector_vars;
486 
487  // We'd better not have more than dim*this->n_vars() (all vector variables)
488  // Treat the NodeElem-only mesh case as dim=1
489  libmesh_assert_less_equal ( nv, (dim > 0 ? dim : 1)*this->n_vars() );
490 
491  // 'nv' represents the max possible number of output variables, so allocate enough memory for
492  // all variables in the system to be populated here. When this is called more than once on a
493  // single 'var_names' vector, different filters should be used such that duplicates don't occur.
494  var_names.resize( nv );
495  }
496 
497  for (const auto & [sys_name, sys_ptr] : _systems)
498  {
499  // Check current system is listed in system_names, and skip pos if not
500  bool use_current_system = (system_names == nullptr);
501  if (!use_current_system)
502  use_current_system = system_names->count(sys_name);
503  if (!use_current_system || sys_ptr->hide_output())
504  continue;
505 
506  for (auto vn : make_range(sys_ptr->n_vars()))
507  {
508  const std::string & var_name = sys_ptr->variable_name(vn);
509  const FEType & fe_type = sys_ptr->variable_type(vn);
510 
511  unsigned int n_vec_dim = FEInterface::n_vec_dim( sys_ptr->get_mesh(), fe_type);
512 
513  // Filter on the type if requested
514  if (type == nullptr || (type && *type == fe_type))
515  {
516  if (FEInterface::field_type(fe_type) == TYPE_VECTOR)
517  {
518  switch(n_vec_dim)
519  {
520  case 0:
521  case 1:
522  var_names[var_num++] = var_name;
523  libmesh_error_msg_if(seen_names.count(var_name) > 1,
524  "Duplicate variable name "+var_name);
525  break;
526  case 2:
527  var_names[var_num++] = var_name+"_x";
528  var_names[var_num++] = var_name+"_y";
529  libmesh_error_msg_if(seen_names.count(var_name+"_x"),
530  "Duplicate variable name "+var_name+"_x");
531  libmesh_error_msg_if(seen_names.count(var_name+"_y"),
532  "Duplicate variable name "+var_name+"_y");
533  break;
534  case 3:
535  var_names[var_num++] = var_name+"_x";
536  var_names[var_num++] = var_name+"_y";
537  var_names[var_num++] = var_name+"_z";
538  libmesh_error_msg_if(seen_names.count(var_name+"_x"),
539  "Duplicate variable name "+var_name+"_x");
540  libmesh_error_msg_if(seen_names.count(var_name+"_y"),
541  "Duplicate variable name "+var_name+"_y");
542  libmesh_error_msg_if(seen_names.count(var_name+"_z"),
543  "Duplicate variable name "+var_name+"_z");
544  break;
545  default:
546  libmesh_error_msg("Invalid dim in build_variable_names");
547  }
548  }
549  else
550  var_names[var_num++] = var_name;
551  }
552  }
553  }
554  // Now resize again in case we filtered any names
555  var_names.resize(var_num);
556 }
557 
558 
559 
560 void EquationSystems::build_solution_vector (std::vector<Number> &,
561  std::string_view,
562  std::string_view) const
563 {
564  // TODO:[BSK] re-implement this from the method below
565  libmesh_not_implemented();
566 }
567 
568 
569 
570 
571 std::unique_ptr<NumericVector<Number>>
572 EquationSystems::build_parallel_solution_vector(const std::set<std::string> * system_names,
573  bool add_sides) const
574 {
575  LOG_SCOPE("build_parallel_solution_vector()", "EquationSystems");
576 
577  // This function must be run on all processors at once
578  parallel_object_only();
579 
580  const unsigned int dim = _mesh.spatial_dimension();
581  const dof_id_type max_nn = _mesh.max_node_id();
582 
583  // allocate vector storage to hold
584  // (max_node_id)*(number_of_variables) entries.
585  //
586  // If node renumbering is disabled and adaptive coarsening has
587  // created gaps between node numbers, then this vector will be
588  // sparse.
589  //
590  // We have to differentiate between between scalar and vector
591  // variables. We intercept vector variables and treat each
592  // component as a scalar variable (consistently with build_solution_names).
593 
594  unsigned int nv = 0;
595 
596  //Could this be replaced by a/some convenience methods?[PB]
597  {
598  unsigned int n_scalar_vars = 0;
599  unsigned int n_vector_vars = 0;
600  for (const auto & [sys_name, sys_ptr] : _systems)
601  {
602  // Check current system is listed in system_names, and skip pos if not
603  bool use_current_system = (system_names == nullptr);
604  if (!use_current_system)
605  use_current_system = system_names->count(sys_name);
606  if (!use_current_system || sys_ptr->hide_output())
607  continue;
608 
609  for (auto vn : make_range(sys_ptr->n_vars()))
610  {
611  if (FEInterface::field_type(sys_ptr->variable_type(vn)) == TYPE_VECTOR)
612  n_vector_vars++;
613  else
614  n_scalar_vars++;
615  }
616  }
617  // Here, we're assuming the number of vector components is the same
618  // as the mesh spatial dimension.
619  nv = n_scalar_vars + dim*n_vector_vars;
620  }
621 
622  // Get the number of nodes to store locally.
623  dof_id_type n_local_nodes = cast_int<dof_id_type>
624  (std::distance(_mesh.local_nodes_begin(),
625  _mesh.local_nodes_end()));
626 
627  // If node renumbering has been disabled, nodes may not be numbered
628  // contiguously, and the number of nodes might not match the
629  // max_node_id. In this case we just do our best.
630  dof_id_type n_total_nodes = n_local_nodes;
631  _mesh.comm().sum(n_total_nodes);
632 
633  const processor_id_type n_proc = _mesh.comm().size();
634  const processor_id_type my_pid = _mesh.comm().rank();
635  const dof_id_type n_gaps = max_nn - n_total_nodes;
636  const dof_id_type gaps_per_processor = n_gaps / n_proc;
637  const dof_id_type remainder_gaps = n_gaps % n_proc;
638 
639  n_local_nodes = n_local_nodes + // Actual nodes
640  gaps_per_processor + // Our even share of gaps
641  (my_pid < remainder_gaps); // Leftovers
642 
643  // If we've been asked to build added sides' data, we need space to
644  // add it. Keep track of how much space.
645  dof_id_type local_added_side_nodes = 0,
646  added_side_nodes = 0;
647 
648  // others_added_side_nodes[p]: local_added_side_nodes on rank p
649  std::vector<dof_id_type> others_added_side_nodes;
650 
651  // A map of (element_id, side, side_node) pairs to the corresponding
652  // added side node index.
653  std::map<std::tuple<dof_id_type, unsigned short, unsigned short>,
654  dof_id_type> discontinuous_node_indices;
655 
656  // If we don't have any added side nodes, we'll have no offsets from
657  // them, and we won't care about which offsets apply to which node
658  // ids either.
659 
660  // Number of true nodes on processors [0,p]
661  std::vector<dof_id_type> true_node_offsets;
662  // Number of added (fake) nodes on processors [0,p)
663  std::vector<dof_id_type> added_node_offsets;
664 
665  auto node_id_to_vec_id =
666  [&true_node_offsets, &added_node_offsets]
667  (const dof_id_type node_id)
668  {
669  if (true_node_offsets.empty())
670  return node_id; // O(1) in the common !add_sides case
671 
672  // Find the processor id that has node_id in the parallel vec
673  const auto lb = std::upper_bound(true_node_offsets.begin(),
674  true_node_offsets.end(), node_id);
675  libmesh_assert(lb != true_node_offsets.end());
676  const processor_id_type p = lb - true_node_offsets.begin();
677 
678  return node_id + added_node_offsets[p];
679  };
680 
681  if (add_sides)
682  {
683  true_node_offsets.resize(n_proc);
684  added_node_offsets.resize(n_proc);
685 
686  // One loop to count everyone's new side nodes
687  for (const auto & elem : _mesh.active_element_ptr_range())
688  {
689  for (auto s : elem->side_index_range())
690  {
691  if (redundant_added_side(*elem,s))
692  continue;
693 
694  const std::vector<unsigned int> side_nodes =
695  elem->nodes_on_side(s);
696 
697  if (elem->processor_id() == this->processor_id())
698  local_added_side_nodes += side_nodes.size();
699  }
700  }
701 
702  others_added_side_nodes.resize(n_proc);
703  _mesh.comm().allgather(local_added_side_nodes,
704  others_added_side_nodes);
705 
706  added_side_nodes = std::accumulate(others_added_side_nodes.begin(),
707  others_added_side_nodes.end(), 0,
708  std::plus<>());
709 
710  _mesh.comm().allgather(n_local_nodes, true_node_offsets);
711  for (auto p : make_range(n_proc-1))
712  true_node_offsets[p+1] += true_node_offsets[p];
713  libmesh_assert_equal_to(true_node_offsets[n_proc-1], _mesh.max_node_id());
714 
715  // For nodes that exist in the mesh, we just need an offset to
716  // tell where to put their solutions.
717  added_node_offsets[0] = 0;
718  for (auto p : make_range(n_proc-1))
719  added_node_offsets[p+1] =
720  added_node_offsets[p] + others_added_side_nodes[p];
721 
722  // For added side nodes, we need to fill a map. Start after all
723  // the true node for our pid plus all the side nodes for
724  // previous pids
725  dof_id_type node_counter = true_node_offsets[my_pid];
726  for (auto p : make_range(my_pid))
727  node_counter += others_added_side_nodes[p];
728 
729  // One loop to figure out whose added side nodes get which index
730  for (const auto & elem : _mesh.active_local_element_ptr_range())
731  {
732  for (auto s : elem->side_index_range())
733  {
734  if (redundant_added_side(*elem,s))
735  continue;
736 
737  const std::vector<unsigned int> side_nodes =
738  elem->nodes_on_side(s);
739 
740  for (auto n : index_range(side_nodes))
741  discontinuous_node_indices
742  [std::make_tuple(elem->id(),s,n)] = node_counter++;
743  }
744  }
745  }
746 
747  const dof_id_type
748  n_global_vals = (max_nn + added_side_nodes) * nv,
749  n_local_vals = (n_local_nodes + local_added_side_nodes) * nv;
750 
751  // Create a NumericVector to hold the parallel solution
752  std::unique_ptr<NumericVector<Number>> parallel_soln_ptr = NumericVector<Number>::build(_communicator);
753  NumericVector<Number> & parallel_soln = *parallel_soln_ptr;
754  parallel_soln.init(n_global_vals, n_local_vals, false, PARALLEL);
755 
756  // Create a NumericVector to hold the "repeat_count" for each node - this is essentially
757  // the number of elements contributing to that node's value
758  std::unique_ptr<NumericVector<Number>> repeat_count_ptr = NumericVector<Number>::build(_communicator);
759  NumericVector<Number> & repeat_count = *repeat_count_ptr;
760  repeat_count.init(n_global_vals, n_local_vals, false, PARALLEL);
761 
762  repeat_count.close();
763 
764  unsigned int var_num=0;
765 
766  // For each system in this EquationSystems object,
767  // update the global solution and if we are on processor 0,
768  // loop over the elements and build the nodal solution
769  // from the element solution. Then insert this nodal solution
770  // into the vector passed to build_solution_vector.
771  for (const auto & [sys_name, sys_ptr] : _systems)
772  {
773  // Check current system is listed in system_names, and skip pos if not
774  bool use_current_system = (system_names == nullptr);
775  if (!use_current_system)
776  use_current_system = system_names->count(sys_name);
777  if (!use_current_system || sys_ptr->hide_output())
778  continue;
779 
780  const System & system = *sys_ptr;
781  const unsigned int nv_sys = system.n_vars();
782  const unsigned int sys_num = system.number();
783 
784  //Could this be replaced by a/some convenience methods?[PB]
785  unsigned int n_scalar_vars = 0;
786  unsigned int n_vector_vars = 0;
787  for (auto vn : make_range(sys_ptr->n_vars()))
788  {
789  if (FEInterface::field_type(sys_ptr->variable_type(vn)) == TYPE_VECTOR)
790  n_vector_vars++;
791  else
792  n_scalar_vars++;
793  }
794 
795  // Here, we're assuming the number of vector components is the same
796  // as the mesh spatial dimension.
797  unsigned int nv_sys_split = n_scalar_vars + dim*n_vector_vars;
798 
799  // Update the current_local_solution
800  {
801  System & non_const_sys = const_cast<System &>(system);
802  // We used to simply call non_const_sys.solution->close()
803  // here, but that is not allowed when the solution vector is
804  // locked read-only, for example when printing the solution
805  // during the middle of a solve... So try to be a bit
806  // more careful about calling close() unnecessarily.
807  libmesh_assert(this->comm().verify(non_const_sys.solution->closed()));
808  if (!non_const_sys.solution->closed())
809  non_const_sys.solution->close();
810  non_const_sys.update();
811  }
812 
813  NumericVector<Number> & sys_soln(*system.current_local_solution);
814 
815  const DofMap & dof_map = system.get_dof_map();
816 
817  std::vector<Number> elem_soln; // The finite element solution
818  std::vector<Number> nodal_soln; // The FE solution interpolated to the nodes
819  std::vector<dof_id_type> dof_indices; // The DOF indices for the finite element
820 
821  unsigned var_inc = 0;
822  for (unsigned int var=0; var<nv_sys; var++)
823  {
824  const FEType & fe_type = system.variable_type(var);
825  const Variable & var_description = system.variable(var);
826  unsigned int n_vec_dim = FEInterface::n_vec_dim( sys_ptr->get_mesh(), fe_type );
827  const auto vg = dof_map.var_group_from_var_number(var);
828  const bool add_p_level = dof_map.should_p_refine(vg);
829 
830  for (const auto & elem : _mesh.active_local_element_ptr_range())
831  {
832  if (var_description.active_on_subdomain(elem->subdomain_id()))
833  {
834  dof_map.dof_indices (elem, dof_indices, var);
835  sys_soln.get(dof_indices, elem_soln);
836 
837  FEInterface::nodal_soln (elem->dim(),
838  fe_type,
839  elem,
840  elem_soln,
841  nodal_soln,
842  add_p_level,
843  n_vec_dim);
844 
845  // infinite elements should be skipped...
846  if (!elem->infinite())
847  {
848  libmesh_assert_equal_to (nodal_soln.size(), n_vec_dim*elem->n_nodes());
849 
850  for (auto n : elem->node_index_range())
851  {
852  const Node & node = elem->node_ref(n);
853 
854  const dof_id_type node_idx =
855  nv * node_id_to_vec_id(node.id());
856 
857  for (unsigned int d=0; d < n_vec_dim; d++)
858  {
859  // For vector-valued elements, all components are in nodal_soln. For each
860  // node, the components are stored in order, i.e. node_0 -> s0_x, s0_y, s0_z
861  parallel_soln.add(node_idx + (var_inc+d + var_num), nodal_soln[n_vec_dim*n+d]);
862 
863  // Increment the repeat count for this position
864  repeat_count.add(node_idx + (var_inc+d + var_num), 1);
865  }
866  }
867 
868  if (add_sides)
869  {
870  for (auto s : elem->side_index_range())
871  {
872  if (redundant_added_side(*elem,s))
873  continue;
874 
875  // Compute the FE solution at all the
876  // side nodes
878  (fe_type, elem, s, elem_soln,
879  nodal_soln, add_p_level, n_vec_dim);
880 
881 #ifdef DEBUG
882  const std::vector<unsigned int> side_nodes =
883  elem->nodes_on_side(s);
884 
885  libmesh_assert_equal_to
886  (nodal_soln.size(),
887  side_nodes.size());
888 #endif
889 
890  for (auto n : index_range(nodal_soln))
891  {
892  // Retrieve index into global solution vector.
893  std::size_t node_index =
894  nv * libmesh_map_find(discontinuous_node_indices,
895  std::make_tuple(elem->id(), s, n));
896 
897  for (unsigned int d=0; d < n_vec_dim; d++)
898  {
899  parallel_soln.add(node_index + (var_inc+d + var_num), nodal_soln[n_vec_dim*n+d]);
900  repeat_count.add(node_index + (var_inc+d + var_num), 1);
901  }
902  }
903  }
904  }
905  }
906  }
907  else // If this variable doesn't exist on this subdomain we have to still increment repeat_count so that we won't divide by 0 later:
908  for (auto n : elem->node_index_range())
909  {
910  const Node & node = elem->node_ref(n);
911  // Only do this if this variable has NO DoFs at
912  // this node... it might have some from an
913  // adjoining element...
914  if (!node.n_dofs(sys_num, var))
915  {
916  const dof_id_type node_idx =
917  nv * node_id_to_vec_id(node.id());
918 
919  for (unsigned int d=0; d < n_vec_dim; d++)
920  repeat_count.add(node_idx + (var_inc+d + var_num), 1);
921  }
922  }
923 
924  } // end loop over elements
925  var_inc += n_vec_dim;
926  } // end loop on variables in this system
927 
928  var_num += nv_sys_split;
929  } // end loop over systems
930 
931  // Sum the nodal solution values and repeat counts.
932  parallel_soln.close();
933  repeat_count.close();
934 
935  // If there were gaps in the node numbering, there will be
936  // corresponding zeros in the parallel_soln and repeat_count
937  // vectors. We need to set those repeat_count entries to 1
938  // in order to avoid dividing by zero.
939  if (n_gaps)
940  {
941  for (numeric_index_type i=repeat_count.first_local_index();
942  i<repeat_count.last_local_index(); ++i)
943  {
944  // repeat_count entries are integral values but let's avoid a
945  // direct floating point comparison with 0 just in case some
946  // roundoff noise crept in during vector assembly?
947  if (std::abs(repeat_count(i)) < TOLERANCE)
948  repeat_count.set(i, 1.);
949  }
950 
951  // Make sure the repeat_count vector is up-to-date on all
952  // processors.
953  repeat_count.close();
954  }
955 
956  // Divide to get the average value at the nodes
957  parallel_soln /= repeat_count;
958 
959  return parallel_soln_ptr;
960 }
961 
962 
963 
964 void EquationSystems::build_solution_vector (std::vector<Number> & soln,
965  const std::set<std::string> * system_names,
966  bool add_sides) const
967 {
968  LOG_SCOPE("build_solution_vector()", "EquationSystems");
969 
970  // Call the parallel implementation
971  std::unique_ptr<NumericVector<Number>> parallel_soln =
972  this->build_parallel_solution_vector(system_names, add_sides);
973 
974  // Localize the NumericVector into the provided std::vector.
975  parallel_soln->localize_to_one(soln);
976 }
977 
978 
979 
980 void EquationSystems::get_vars_active_subdomains(const std::vector<std::string> & names,
981  std::vector<std::set<subdomain_id_type>> & vars_active_subdomains) const
982 {
983  vars_active_subdomains.clear();
984  vars_active_subdomains.resize(names.size());
985 
986  for (const auto & pr : _systems)
987  {
988  const auto & sys_ptr = pr.second;
989  for (auto vn : make_range(sys_ptr->n_vars()))
990  {
991  const std::string & var_name = sys_ptr->variable_name(vn);
992 
993  auto names_it = std::find(names.begin(), names.end(), var_name);
994  if(names_it != names.end())
995  {
996  const Variable & variable = sys_ptr->variable(vn);
997  const std::set<subdomain_id_type> & active_subdomains = variable.active_subdomains();
998  vars_active_subdomains[std::distance(names.begin(), names_it)] = active_subdomains;
999  }
1000  }
1001  }
1002 }
1003 
1004 
1005 
1006 #ifdef LIBMESH_ENABLE_DEPRECATED
1007 void EquationSystems::get_solution (std::vector<Number> & soln,
1008  std::vector<std::string> & names) const
1009 {
1010  libmesh_deprecated();
1011  this->build_elemental_solution_vector(soln, names);
1012 }
1013 #endif // LIBMESH_ENABLE_DEPRECATED
1014 
1015 
1016 
1017 void
1019  std::vector<std::string> & names) const
1020 {
1021  // Call the parallel version of this function
1022  std::unique_ptr<NumericVector<Number>> parallel_soln =
1024 
1025  // Localize into 'soln', provided that parallel_soln is not empty.
1026  // Note: parallel_soln will be empty in the event that none of the
1027  // input names were CONSTANT, MONOMIAL nor components of CONSTANT,
1028  // MONOMIAL_VEC variables, or there were simply none of these in
1029  // the EquationSystems object.
1030  soln.clear();
1031  if (parallel_soln)
1032  parallel_soln->localize_to_one(soln);
1033 }
1034 
1035 std::vector<std::pair<unsigned int, unsigned int>>
1037  (std::vector<std::string> & names, const FEType * type, const std::vector<FEType> * types) const
1038 {
1039  // This function must be run on all processors at once
1040  parallel_object_only();
1041 
1042  libmesh_assert (this->n_systems());
1043 
1044  // Resolve class of type input and assert that at least one of them is null
1045  libmesh_assert_msg(!type || !types,
1046  "Input 'type', 'types', or neither in find_variable_numbers, but not both.");
1047 
1048  std::vector<FEType> type_filter;
1049  if (type)
1050  type_filter.push_back(*type);
1051  else if (types)
1052  type_filter = *types;
1053 
1054  // Store a copy of the valid variable names, if any. The names vector will be repopulated with any
1055  // valid names (or all if 'is_names_empty') in the system that passes through the type filter. If
1056  // the variable is a vector, its name will be decomposed into its separate components in
1057  // accordance with build_variable_names().
1058  std::vector<std::string> name_filter = names;
1059  bool is_names_empty = name_filter.empty();
1060  names.clear();
1061 
1062  // initialize convenience variables
1063  FEType var_type;
1064  std::string name;
1065 
1066  const std::vector<std::string> component_suffix = {"_x", "_y", "_z"};
1067  unsigned int dim = _mesh.spatial_dimension();
1068  libmesh_error_msg_if(dim > 3, "Invalid dim in find_variable_numbers");
1069 
1070  // Now filter through the variables in each system and store the system index and their index
1071  // within that system. This way, we know where to find their data even after we sort them.
1072  std::vector<std::pair<unsigned int, unsigned int>> var_nums;
1073 
1074  for (const auto & pr : _systems)
1075  {
1076  const System & system = *(pr.second);
1077 
1078  for (auto var : make_range(system.n_vars()))
1079  {
1080  // apply the type filter
1081  var_type = system.variable_type(var);
1082  if (type_filter.size() &&
1083  std::find(type_filter.begin(), type_filter.end(), var_type) == type_filter.end())
1084  continue;
1085 
1086  // apply the name filter (note that all variables pass if it is empty)
1087  if (FEInterface::field_type(var_type) == TYPE_VECTOR)
1088  {
1089  std::vector<std::string> component_names;
1090  for (unsigned int comp = 0; comp < dim; ++comp)
1091  {
1092  name = system.variable_name(var) + component_suffix[comp];
1093  if (is_names_empty ||
1094  (std::find(name_filter.begin(), name_filter.end(), name) != name_filter.end()))
1095  component_names.push_back(name);
1096  }
1097 
1098  if (! component_names.empty())
1099  names.insert(names.end(), component_names.begin(), component_names.end());
1100  else
1101  continue;
1102  }
1103  else /*scalar-valued variable*/
1104  {
1105  name = system.variable_name(var);
1106  if (is_names_empty ||
1107  (std::find(name_filter.begin(), name_filter.end(), name) != name_filter.end()))
1108  names.push_back(name);
1109  else
1110  continue;
1111  }
1112 
1113  // if the variable made it through both filters get its system indices
1114  var_nums.emplace_back(system.number(), var);
1115  }
1116  }
1117 
1118  // Sort the var_nums vector pairs alphabetically based on the variable name
1119  std::vector<unsigned int> sort_index(var_nums.size());
1120  std::iota(sort_index.begin(), sort_index.end(), 0);
1121  std::sort(sort_index.begin(), sort_index.end(),
1122  [&](const unsigned int & lhs, const unsigned int & rhs)
1123  {return this->get_system(var_nums[lhs].first).variable_name(var_nums[lhs].second) <
1124  this->get_system(var_nums[rhs].first).variable_name(var_nums[rhs].second);});
1125 
1126  std::vector<std::pair<unsigned int, unsigned int>> var_nums_sorted(var_nums.size());
1127  for (auto i : index_range(var_nums_sorted))
1128  {
1129  var_nums_sorted[i].first = var_nums[sort_index[i]].first;
1130  var_nums_sorted[i].second = var_nums[sort_index[i]].second;
1131  }
1132 
1133  // Also sort the names vector
1134  std::sort(names.begin(), names.end());
1135 
1136  // Return the sorted vector pairs
1137  return var_nums_sorted;
1138 }
1139 
1140 
1141 std::unique_ptr<NumericVector<Number>>
1142 EquationSystems::build_parallel_elemental_solution_vector (std::vector<std::string> & names) const
1143 {
1144  // Filter any names that aren't elemental variables and get the system indices for those that are.
1145  // Note that it's probably fine if the names vector is empty since we'll still at least filter
1146  // out all non-monomials. If there are no monomials, then nothing is output here.
1147  std::vector<FEType> type = {FEType(CONSTANT, MONOMIAL), FEType(CONSTANT, MONOMIAL_VEC)};
1148  std::vector<std::pair<unsigned int, unsigned int>> var_nums =
1149  this->find_variable_numbers(names, /*type=*/nullptr, &type);
1150 
1151  const std::size_t nv = names.size(); /*total number of vars including vector components*/
1152  const dof_id_type ne = _mesh.n_elem();
1153  libmesh_assert_equal_to (ne, _mesh.max_elem_id());
1154 
1155  // If there are no variables to write out don't do anything...
1156  if (!nv)
1157  return std::unique_ptr<NumericVector<Number>>(nullptr);
1158 
1159  // We can handle the case where there are nullptrs in the Elem vector
1160  // by just having extra zeros in the solution vector.
1161  numeric_index_type parallel_soln_global_size = ne*nv;
1162 
1163  numeric_index_type div = parallel_soln_global_size / this->n_processors();
1164  numeric_index_type mod = parallel_soln_global_size % this->n_processors();
1165 
1166  // Initialize all processors to the average size.
1167  numeric_index_type parallel_soln_local_size = div;
1168 
1169  // The first "mod" processors get an extra entry.
1170  if (this->processor_id() < mod)
1171  parallel_soln_local_size = div+1;
1172 
1173  // Create a NumericVector to hold the parallel solution
1174  std::unique_ptr<NumericVector<Number>> parallel_soln_ptr = NumericVector<Number>::build(_communicator);
1175  NumericVector<Number> & parallel_soln = *parallel_soln_ptr;
1176  parallel_soln.init(parallel_soln_global_size,
1177  parallel_soln_local_size,
1178  /*fast=*/false,
1179  /*ParallelType=*/PARALLEL);
1180 
1181  unsigned int sys_ctr = 0;
1182  unsigned int var_ctr = 0;
1183  for (auto i : index_range(var_nums))
1184  {
1185  std::pair<unsigned int, unsigned int> var_num = var_nums[i];
1186  const System & system = this->get_system(var_num.first);
1187 
1188  // Update the current_local_solution if necessary
1189  if (sys_ctr != var_num.first)
1190  {
1191  System & non_const_sys = const_cast<System &>(system);
1192  // We used to simply call non_const_sys.solution->close()
1193  // here, but that is not allowed when the solution vector is
1194  // locked read-only, for example when printing the solution
1195  // during during the middle of a solve... So try to be a bit
1196  // more careful about calling close() unnecessarily.
1197  libmesh_assert(this->comm().verify(non_const_sys.solution->closed()));
1198  if (!non_const_sys.solution->closed())
1199  non_const_sys.solution->close();
1200  non_const_sys.update();
1201  sys_ctr = var_num.first;
1202  }
1203 
1204  NumericVector<Number> & sys_soln(*system.current_local_solution);
1205 
1206  // The DOF indices for the finite element
1207  std::vector<dof_id_type> dof_indices;
1208 
1209  const unsigned int var = var_num.second;
1210 
1211  const Variable & variable = system.variable(var);
1212  const DofMap & dof_map = system.get_dof_map();
1213 
1214  // We need to check if the constant monomial is a scalar or a vector and set the number of
1215  // components as the mesh spatial dimension for the latter as per es.find_variable_numbers().
1216  // Even for the case where a variable is not active on any subdomain belonging to the
1217  // processor, we still need to know this number to update 'var_ctr'.
1218  const unsigned int n_comps =
1219  (system.variable_type(var) == type[1]) ? _mesh.spatial_dimension() : 1;
1220 
1221  // Loop over all elements in the mesh and index all components of the variable if it's active
1222  for (const auto & elem : _mesh.active_local_element_ptr_range())
1223  if (variable.active_on_subdomain(elem->subdomain_id()))
1224  {
1225  dof_map.dof_indices(elem, dof_indices, var);
1226 
1227  // The number of DOF components needs to be equal to the expected number so that we know
1228  // where to store data to correctly correspond to variable names.
1229  libmesh_assert_equal_to(dof_indices.size(), n_comps);
1230 
1231  for (unsigned int comp = 0; comp < n_comps; comp++)
1232  parallel_soln.set(ne * (var_ctr + comp) + elem->id(), sys_soln(dof_indices[comp]));
1233  }
1234 
1235  var_ctr += n_comps;
1236  } // end loop over var_nums
1237 
1238  // NOTE: number of output names might not be equal to the number passed to this function. Any that
1239  // aren't CONSTANT MONOMIALS or components of CONSTANT MONOMIAL_VECS have been filtered out (see
1240  // EquationSystems::find_variable_numbers).
1241  //
1242  // But, if everything is accounted for properly, then names.size() == var_ctr
1243  libmesh_assert_equal_to(names.size(), var_ctr);
1244 
1245  parallel_soln.close();
1246  return parallel_soln_ptr;
1247 }
1248 
1249 
1250 
1251 void
1253 (std::vector<Number> & soln,
1254  const std::set<std::string> * system_names,
1255  const std::vector<std::string> * var_names,
1256  bool vertices_only,
1257  bool add_sides) const
1258 {
1259  LOG_SCOPE("build_discontinuous_solution_vector()", "EquationSystems");
1260 
1261  libmesh_assert (this->n_systems());
1262 
1263  // Get the number of variables (nv) by counting the number of variables
1264  // in each system listed in system_names
1265  unsigned int nv = 0;
1266 
1267  for (const auto & [sys_name, sys_ptr] : _systems)
1268  {
1269  // Check current system is listed in system_names, and skip pos if not
1270  bool use_current_system = (system_names == nullptr);
1271  if (!use_current_system)
1272  use_current_system = system_names->count(sys_name);
1273  if (!use_current_system || sys_ptr->hide_output())
1274  continue;
1275 
1276  // Loop over all variables in this System and check whether we
1277  // are supposed to use each one.
1278  for (auto var_id : make_range(sys_ptr->n_vars()))
1279  {
1280  bool use_current_var = (var_names == nullptr);
1281  if (!use_current_var)
1282  use_current_var = std::count(var_names->begin(),
1283  var_names->end(),
1284  sys_ptr->variable_name(var_id));
1285 
1286  // Only increment the total number of vars if we are
1287  // supposed to use this one.
1288  if (use_current_var)
1289  nv++;
1290  }
1291  }
1292 
1293  // get the total "weight" - the number of nodal values to write for
1294  // each variable.
1295  unsigned int tw=0;
1296  for (const auto & elem : _mesh.active_element_ptr_range())
1297  {
1298  tw += vertices_only ? elem->n_vertices() : elem->n_nodes();
1299 
1300  if (add_sides)
1301  {
1302  for (auto s : elem->side_index_range())
1303  {
1304  if (redundant_added_side(*elem,s))
1305  continue;
1306 
1307  const std::vector<unsigned int> side_nodes =
1308  elem->nodes_on_side(s);
1309 
1310  if (!vertices_only)
1311  tw += side_nodes.size();
1312  else
1313  for (auto n : index_range(side_nodes))
1314  if (elem->is_vertex(side_nodes[n]))
1315  ++tw;
1316  }
1317  }
1318  }
1319 
1320  // Only if we are on processor zero, allocate the storage
1321  // to hold (number_of_nodes)*(number_of_variables) entries.
1322  if (_mesh.processor_id() == 0)
1323  soln.resize(tw*nv);
1324 
1325  std::vector<Number> sys_soln;
1326 
1327  // Keep track of the variable "offset". This is used for indexing
1328  // into the global solution vector.
1329  unsigned int var_offset = 0;
1330 
1331  // For each system in this EquationSystems object,
1332  // update the global solution and if we are on processor 0,
1333  // loop over the elements and build the nodal solution
1334  // from the element solution. Then insert this nodal solution
1335  // into the vector passed to build_solution_vector.
1336  for (const auto & [sys_name, system] : _systems)
1337  {
1338  // Check current system is listed in system_names, and skip pos if not
1339  bool use_current_system = (system_names == nullptr);
1340  if (!use_current_system)
1341  use_current_system = system_names->count(sys_name);
1342  if (!use_current_system || system->hide_output())
1343  continue;
1344 
1345  const unsigned int nv_sys = system->n_vars();
1346  const auto & dof_map = system->get_dof_map();
1347 
1348  system->update_global_solution (sys_soln, 0);
1349 
1350  // Keep track of the number of vars actually written.
1351  unsigned int n_vars_written_current_system = 0;
1352 
1353  if (_mesh.processor_id() == 0)
1354  {
1355  std::vector<Number> soln_coeffs; // The finite element solution coeffs
1356  std::vector<Number> nodal_soln; // The FE solution interpolated to the nodes
1357  std::vector<dof_id_type> dof_indices; // The DOF indices for the finite element
1358 
1359  // For each variable, determine if we are supposed to
1360  // write it, then loop over the active elements, compute
1361  // the nodal_soln and store it to the "soln" vector. We
1362  // store zeros for subdomain-restricted variables on
1363  // elements where they are not active.
1364  for (unsigned int var=0; var<nv_sys; var++)
1365  {
1366  bool use_current_var = (var_names == nullptr);
1367  if (!use_current_var)
1368  use_current_var = std::count(var_names->begin(),
1369  var_names->end(),
1370  system->variable_name(var));
1371 
1372  // If we aren't supposed to write this var, go to the
1373  // next loop iteration.
1374  if (!use_current_var)
1375  continue;
1376 
1377  const FEType & fe_type = system->variable_type(var);
1378  const Variable & var_description = system->variable(var);
1379  const auto vg = dof_map.var_group_from_var_number(var);
1380  const bool add_p_level = dof_map.should_p_refine(vg);
1381 
1382  unsigned int nn=0;
1383 
1384  for (auto & elem : _mesh.active_element_ptr_range())
1385  {
1386  if (var_description.active_on_subdomain(elem->subdomain_id()))
1387  {
1388  system->get_dof_map().dof_indices (elem, dof_indices, var);
1389 
1390  soln_coeffs.resize(dof_indices.size());
1391 
1392  for (auto i : index_range(dof_indices))
1393  soln_coeffs[i] = sys_soln[dof_indices[i]];
1394 
1395  // Compute the FE solution at all the nodes, but
1396  // only use the first n_vertices() entries if
1397  // vertices_only == true.
1398  FEInterface::nodal_soln (elem->dim(),
1399  fe_type,
1400  elem,
1401  soln_coeffs,
1402  nodal_soln,
1403  add_p_level,
1404  FEInterface::n_vec_dim(_mesh, fe_type));
1405 
1406  // infinite elements should be skipped...
1407  if (!elem->infinite())
1408  {
1409  libmesh_assert_equal_to (nodal_soln.size(), elem->n_nodes());
1410 
1411  const unsigned int n_vals =
1412  vertices_only ? elem->n_vertices() : elem->n_nodes();
1413 
1414  for (unsigned int n=0; n<n_vals; n++)
1415  {
1416  // Compute index into global solution vector.
1417  std::size_t index =
1418  nv * (nn++) + (n_vars_written_current_system + var_offset);
1419 
1420  soln[index] += nodal_soln[n];
1421  }
1422  }
1423  }
1424  else
1425  nn += vertices_only ? elem->n_vertices() : elem->n_nodes();
1426  } // end loop over active elements writing interiors
1427 
1428  // Loop writing "fake" sides, if requested
1429  if (add_sides)
1430  {
1431  // We don't build discontinuous solution vectors in
1432  // parallel yet, but we'll do ordering of fake side
1433  // values as if we did, for consistency with the
1434  // parallel continuous ordering and for future
1435  // compatibility.
1436  std::vector<std::vector<const Elem *>>
1437  elems_by_pid(_mesh.n_processors());
1438 
1439  for (const auto & elem : _mesh.active_element_ptr_range())
1440  elems_by_pid[elem->processor_id()].push_back(elem);
1441 
1442  for (auto p : index_range(elems_by_pid))
1443  for (const Elem * elem : elems_by_pid[p])
1444  {
1445  if (var_description.active_on_subdomain(elem->subdomain_id()))
1446  {
1447  system->get_dof_map().dof_indices (elem, dof_indices, var);
1448 
1449  soln_coeffs.resize(dof_indices.size());
1450 
1451  for (auto i : index_range(dof_indices))
1452  soln_coeffs[i] = sys_soln[dof_indices[i]];
1453 
1454  for (auto s : elem->side_index_range())
1455  {
1456  if (redundant_added_side(*elem,s))
1457  continue;
1458 
1459  const std::vector<unsigned int> side_nodes =
1460  elem->nodes_on_side(s);
1461 
1462  // Compute the FE solution at all the
1463  // side nodes, but only use those for
1464  // which is_vertex() == true if
1465  // vertices_only == true.
1467  (fe_type, elem, s, soln_coeffs,
1468  nodal_soln, add_p_level,
1469  FEInterface::n_vec_dim(_mesh, fe_type));
1470 
1471  libmesh_assert_equal_to
1472  (nodal_soln.size(),
1473  side_nodes.size());
1474 
1475  // If we don't have a continuous FE
1476  // then we want to average between
1477  // sides, at least in the equal-level
1478  // case where it's easy. This is
1479  // analogous to our repeat_count
1480  // behavior elsewhere.
1481  const FEContinuity cont =
1482  FEInterface::get_continuity(fe_type);
1483  const Elem * const neigh = elem->neighbor_ptr(s);
1484 
1485  if ((cont == DISCONTINUOUS || cont == H_CURL || cont == H_DIV) &&
1486  neigh &&
1487  neigh->level() == elem->level() &&
1488  var_description.active_on_subdomain(neigh->subdomain_id()))
1489  {
1490  std::vector<dof_id_type> neigh_indices;
1491  system->get_dof_map().dof_indices (neigh, neigh_indices, var);
1492  std::vector<Number> neigh_coeffs(neigh_indices.size());
1493 
1494  for (auto i : index_range(neigh_indices))
1495  neigh_coeffs[i] = sys_soln[neigh_indices[i]];
1496 
1497  const unsigned int s_neigh =
1498  neigh->which_neighbor_am_i(elem);
1499  std::vector<Number> neigh_soln;
1501  (fe_type, neigh, s_neigh,
1502  neigh_coeffs, neigh_soln, add_p_level,
1503  FEInterface::n_vec_dim(_mesh, fe_type));
1504 
1505  const std::vector<unsigned int> neigh_nodes =
1506  neigh->nodes_on_side(s_neigh);
1507  for (auto n : index_range(side_nodes))
1508  for (auto neigh_n : index_range(neigh_nodes))
1509  if (neigh->node_ptr(neigh_nodes[neigh_n])
1510  == elem->node_ptr(side_nodes[n]))
1511  {
1512  nodal_soln[n] += neigh_soln[neigh_n];
1513  nodal_soln[n] /= 2;
1514  }
1515  }
1516 
1517  for (auto n : index_range(side_nodes))
1518  {
1519  if (vertices_only &&
1520  !elem->is_vertex(n))
1521  continue;
1522 
1523  // Compute index into global solution vector.
1524  std::size_t index =
1525  nv * (nn++) + (n_vars_written_current_system + var_offset);
1526 
1527  soln[index] += nodal_soln[n];
1528  }
1529  }
1530  }
1531  else
1532  {
1533  nn += vertices_only ? elem->n_vertices() : elem->n_nodes();
1534 
1535  for (auto s : elem->side_index_range())
1536  {
1537  if (redundant_added_side(*elem,s))
1538  continue;
1539 
1540  const std::vector<unsigned int> side_nodes =
1541  elem->nodes_on_side(s);
1542 
1543  for (auto n : index_range(side_nodes))
1544  {
1545  if (vertices_only &&
1546  !elem->is_vertex(n))
1547  continue;
1548  nn++;
1549  }
1550  }
1551  }
1552  } // end loop over active elements, writing "fake" sides
1553  }
1554  // If we made it here, we actually wrote a variable, so increment
1555  // the number of variables actually written for the current system.
1556  n_vars_written_current_system++;
1557 
1558  } // end loop over vars
1559  } // end if proc 0
1560 
1561  // Update offset for next loop iteration.
1562  var_offset += n_vars_written_current_system;
1563  } // end loop over systems
1564 }
1565 
1566 
1567 
1568 bool EquationSystems::redundant_added_side(const Elem & elem, unsigned int side)
1569 {
1570  libmesh_assert(elem.active());
1571 
1572  const Elem * neigh = elem.neighbor_ptr(side);
1573 
1574  // Write boundary sides.
1575  if (!neigh)
1576  return false;
1577 
1578  // Write ghost sides in Nemesis
1579  if (neigh == remote_elem)
1580  return false;
1581 
1582  // Don't write a coarser side if a finer side exists
1583  if (!neigh->active())
1584  return true;
1585 
1586  // Don't write a side redundantly from both of the
1587  // elements sharing it. We'll disambiguate with id().
1588  return (neigh->id() < elem.id());
1589 }
1590 
1591 
1592 
1594  const Real threshold,
1595  const bool verbose) const
1596 {
1597  // safety check, whether we handle at least the same number
1598  // of systems
1599  std::vector<bool> os_result;
1600 
1601  if (this->n_systems() != other_es.n_systems())
1602  {
1603  if (verbose)
1604  {
1605  libMesh::out << " Fatal difference. This system handles "
1606  << this->n_systems() << " systems," << std::endl
1607  << " while the other system handles "
1608  << other_es.n_systems()
1609  << " systems." << std::endl
1610  << " Aborting comparison." << std::endl;
1611  }
1612  return false;
1613  }
1614  else
1615  {
1616  // start comparing each system
1617  for (const auto & [sys_name, sys_ptr] : _systems)
1618  {
1619  // get the other system
1620  const System & other_system = other_es.get_system (sys_name);
1621 
1622  os_result.push_back (sys_ptr->compare (other_system, threshold, verbose));
1623 
1624  }
1625 
1626  }
1627 
1628 
1629  // sum up the results
1630  if (os_result.size()==0)
1631  return true;
1632  else
1633  {
1634  bool os_identical;
1635  unsigned int n = 0;
1636  do
1637  {
1638  os_identical = os_result[n];
1639  n++;
1640  }
1641  while (os_identical && n<os_result.size());
1642  return os_identical;
1643  }
1644 }
1645 
1646 
1647 
1648 std::string EquationSystems::get_info () const
1649 {
1650  std::ostringstream oss;
1651 
1652  oss << " EquationSystems\n"
1653  << " n_systems()=" << this->n_systems() << '\n';
1654 
1655  // Print the info for the individual systems
1656  for (const auto & pr : _systems)
1657  oss << pr.second->get_info();
1658 
1659 
1660  // // Possibly print the parameters
1661  // if (!this->parameters.empty())
1662  // {
1663  // oss << " n_parameters()=" << this->n_parameters() << '\n';
1664  // oss << " Parameters:\n";
1665 
1666  // for (const auto & [key, val] : _parameters)
1667  // oss << " "
1668  // << "\""
1669  // << key
1670  // << "\""
1671  // << "="
1672  // << val
1673  // << '\n';
1674  // }
1675 
1676  return oss.str();
1677 }
1678 
1679 
1680 
1681 void EquationSystems::print_info (std::ostream & os) const
1682 {
1683  os << this->get_info()
1684  << std::endl;
1685 }
1686 
1687 
1688 
1689 std::ostream & operator << (std::ostream & os,
1690  const EquationSystems & es)
1691 {
1692  es.print_info(os);
1693  return os;
1694 }
1695 
1696 
1697 
1698 unsigned int EquationSystems::n_vars () const
1699 {
1700  unsigned int tot=0;
1701 
1702  for (const auto & pr : _systems)
1703  tot += pr.second->n_vars();
1704 
1705  return tot;
1706 }
1707 
1708 
1709 
1710 std::size_t EquationSystems::n_dofs () const
1711 {
1712  std::size_t tot=0;
1713 
1714  for (const auto & pr : _systems)
1715  tot += pr.second->n_dofs();
1716 
1717  return tot;
1718 }
1719 
1720 
1721 
1722 
1724 {
1725  std::size_t tot=0;
1726 
1727  for (const auto & pr : _systems)
1728  tot += pr.second->n_active_dofs();
1729 
1730  return tot;
1731 }
1732 
1733 
1735 {
1736  // All the nodes
1737  for (auto & node : _mesh.node_ptr_range())
1738  node->add_system();
1739 
1740  // All the elements
1741  for (auto & elem : _mesh.element_ptr_range())
1742  elem->add_system();
1743 }
1744 
1746 {
1747  this->get_system(sys_num).get_dof_map().remove_default_ghosting();
1748 }
1749 
1750 } // namespace libMesh
EquationSystems(MeshBase &mesh)
Constructor.
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
void allgather(const T &send_data, std::vector< T, A > &recv_data) const
This is the EquationSystems class.
unsigned int n_vars() const
void build_variable_names(std::vector< std::string > &var_names, const FEType *type=nullptr, const std::set< std::string > *system_names=nullptr) const
Fill the input vector var_names with the names of the variables for each system.
A Node is like a Point, but with more information.
Definition: node.h:52
unsigned int n_systems() const
const Variable & variable(unsigned int var) const
Return a constant reference to Variable var.
Definition: system.h:2445
bool _enable_default_ghosting
Flag for whether to enable default ghosting on newly added Systems.
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
Fills the vector di with the global degree of freedom indices for the element.
Definition: dof_map.C:2184
std::unique_ptr< NumericVector< Number > > build_parallel_elemental_solution_vector(std::vector< std::string > &names) const
Builds a parallel vector of CONSTANT MONOMIAL and/or components of CONSTANT MONOMIAL_VEC solution val...
std::size_t distribute_dofs(MeshBase &)
Distribute dofs on the current mesh.
Definition: dof_map.C:975
static constexpr Real TOLERANCE
std::size_t n_dofs() const
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
unsigned int dim
void build_discontinuous_solution_vector(std::vector< Number > &soln, const std::set< std::string > *system_names=nullptr, const std::vector< std::string > *var_names=nullptr, bool vertices_only=false, bool add_sides=false) const
Fill the input vector soln with solution values.
void add_default_ghosting()
Add the default functor(s) for coupling and algebraic ghosting.
Definition: dof_map.C:2044
bool refine_elements()
Only refines the user-requested elements.
bool reinit_solutions()
Handle any mesh changes and project any solutions onto the updated mesh.
virtual void clear()
Restores the data structure to a pristine state.
virtual void allgather()
Gathers all elements and nodes of the mesh onto every processor.
Definition: mesh_base.h:240
virtual ~EquationSystems()
Destructor.
void sum(T &r) const
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
MeshBase & mesh
static FEFieldType field_type(const FEType &fe_type)
std::ostream & operator<<(std::ostream &os, const OrderWrapper &order)
Overload stream operators.
Definition: fe_type.h:182
processor_id_type rank() const
const Parallel::Communicator & comm() const
virtual void enable_default_ghosting(bool enable)
Enable or disable default ghosting functors on the Mesh and on all Systems.
void clean_refinement_flags()
Sets the refinement flag to Elem::DO_NOTHING for each element in the mesh.
void should_p_refine(unsigned int g, bool p_refine)
Describe whether the given variable group should be p-refined.
Definition: dof_map.h:2496
virtual void adjoint_solve(const QoISet &qoi_indices=QoISet())
Call adjoint_solve on all the individual equation systems.
const Parallel::Communicator & _communicator
The libMesh namespace provides an interface to certain functionality in the library.
const T_sys & get_system(std::string_view name) const
Real distance(const Point &p)
virtual void init(const numeric_index_type n, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC)=0
Change the dimension of the vector to n.
void iota(ForwardIter first, ForwardIter last, T value)
Utility::iota was created back when std::iota was just an SGI STL extension.
Definition: utility.h:229
static void side_nodal_soln(const FEType &fe_t, const Elem *elem, const unsigned int side, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln, const bool add_p_level=true, const unsigned int vdim=1)
Build the nodal soln on one side from the (full) element soln.
Definition: fe_interface.C:650
uint8_t processor_id_type
Definition: id_types.h:104
This is the MeshBase class.
Definition: mesh_base.h:75
virtual std::string get_info() const
virtual void prolong_vectors()
Prolong vectors after the mesh has refined.
Definition: system.C:435
MeshBase & _mesh
The mesh data structure.
bool coarsen_elements()
Only coarsens the user-requested elements.
unsigned int var_group_from_var_number(unsigned int var_num) const
Definition: dof_map.h:2524
signed char & underrefined_boundary_limit()
If underrefined_boundary_limit is set to a nonnegative value, then refinement and coarsening will pro...
unsigned int n_dofs(const unsigned int s, const unsigned int var=libMesh::invalid_uint) const
Definition: dof_object.h:806
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:178
void allgather()
Serializes a distributed mesh and its associated degree of freedom numbering for all systems...
void get_vars_active_subdomains(const std::vector< std::string > &names, std::vector< std::set< subdomain_id_type >> &vars_active_subdomains) const
Retrieve vars_active_subdomains, which indicates the active subdomains for each variable in names...
processor_id_type size() const
virtual void sensitivity_solve(const ParameterVector &parameters)
Call sensitivity_solve on all the individual equation systems.
virtual void reinit_systems()
Reinitialize all systems on the current mesh.
processor_id_type n_processors() const
virtual bool is_serial() const
Definition: mesh_base.h:211
unsigned int number() const
Definition: system.h:2337
This class defines the notion of a variable in the system.
Definition: variable.h:49
const std::set< subdomain_id_type > & active_subdomains() const
Definition: variable.h:171
Implements (adaptive) mesh refinement algorithms for a MeshBase.
dof_id_type id() const
Definition: dof_object.h:828
dof_id_type numeric_index_type
Definition: id_types.h:99
unsigned int which_neighbor_am_i(const Elem *e) const
This function tells you which neighbor e is.
Definition: elem.h:2916
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
virtual void reinit()
Handle any mesh changes and reinitialize all the systems on the updated mesh.
virtual dof_id_type max_elem_id() const =0
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1590
libmesh_assert(ctx)
const std::string & variable_name(const unsigned int i) const
Definition: system.h:2465
void build_solution_vector(std::vector< Number > &soln, std::string_view system_name, std::string_view variable_name="all_vars") const
Fill the input vector soln with the solution values for the system named name.
virtual void restrict_vectors()
Restrict vectors after the mesh has coarsened.
Definition: system.C:377
bool active_on_subdomain(subdomain_id_type sid) const
Definition: variable.h:157
bool _refine_in_reinit
Flag for whether to call coarsen/refine in reinit().
virtual void reinit_mesh()
Handle the association of a completely new mesh with the EquationSystem and all the Systems assigned ...
unsigned char & face_level_mismatch_limit()
If face_level_mismatch_limit is set to a nonzero value, then refinement and coarsening will produce m...
An object whose state is distributed along a set of processors.
static unsigned int n_vec_dim(const MeshBase &mesh, const FEType &fe_type)
virtual void reinit_constraints()
Reinitializes the constraints for this system.
Definition: system.C:482
virtual void close()=0
Calls the NumericVector&#39;s internal assembly routines, ensuring that the values are consistent across ...
static FEContinuity get_continuity(const FEType &fe_type)
Returns the input FEType&#39;s FEContinuity based on the underlying FEFamily and potentially the Order...
virtual void solve()
Call solve on all the individual equation systems.
virtual numeric_index_type first_local_index() const =0
const Elem * neighbor_ptr(unsigned int i) const
Definition: elem.h:2567
void get_solution(std::vector< Number > &soln, std::vector< std::string > &names) const
Retrieve the solution data for CONSTANT MONOMIALs and/or components of CONSTANT MONOMIAL_VECs.
virtual void update()
Update the local values to reflect the solution on neighboring processors.
Definition: system.C:497
unsigned int level() const
Definition: elem.h:3071
const FEType & variable_type(const unsigned int i) const
Definition: system.h:2495
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool contract()=0
Delete subactive (i.e.
subdomain_id_type subdomain_id() const
Definition: elem.h:2551
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:2476
unsigned int spatial_dimension() const
Definition: mesh_base.C:525
T & set(const std::string &)
Definition: parameters.h:469
std::unique_ptr< NumericVector< Number > > build_parallel_solution_vector(const std::set< std::string > *system_names=nullptr, bool add_sides=false) const
A version of build_solution_vector which is appropriate for "parallel" output formats like Nemesis...
OStreamProxy out
virtual void clear()
Clears internal data structures & frees any allocated memory.
Definition: parameters.h:324
signed char & overrefined_boundary_limit()
If overrefined_boundary_limit is set to a nonnegative value, then refinement and coarsening will prod...
void remove_default_ghosting()
Remove any default ghosting functor(s).
Definition: dof_map.C:2036
void _remove_default_ghosting(unsigned int sys_num)
This just calls DofMap::remove_default_ghosting() but using a shim lets us forward-declare DofMap...
const MeshBase & get_mesh() const
std::size_t n_active_dofs() const
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...
std::vector< std::pair< unsigned int, unsigned int > > find_variable_numbers(std::vector< std::string > &names, const FEType *type=nullptr, const std::vector< FEType > *types=nullptr) const
Finds system and variable numbers for any variables of &#39;type&#39; or of &#39;types&#39; corresponding to the entr...
Parameters parameters
Data structure holding arbitrary parameters.
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:1602
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:1884
virtual void init()
Initialize all the systems.
virtual bool compare(const EquationSystems &other_es, const Real threshold, const bool verbose) const
unsigned int n_vars() const
Definition: system.h:2417
virtual dof_id_type max_node_id() const =0
virtual dof_id_type n_elem() const =0
virtual void add(const numeric_index_type i, const T value)=0
Adds value to the vector entry specified by i.
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
processor_id_type processor_id() const
static void nodal_soln(const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln, const bool add_p_level=true, const unsigned int vdim=1)
Build the nodal soln from the element soln.
Definition: fe_interface.C:625
bool active() const
Definition: elem.h:2938
void _add_system_to_nodes_and_elems()
This function is used in the implementation of add_system, it loops over the nodes and elements of th...
const DofMap & get_dof_map() const
Definition: system.h:2361
virtual numeric_index_type last_local_index() const =0
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117
void update()
Updates local values for all the systems.
virtual std::vector< unsigned int > nodes_on_side(const unsigned int) const =0
void build_elemental_solution_vector(std::vector< Number > &soln, std::vector< std::string > &names) const
Retrieve the solution data for CONSTANT MONOMIALs and/or components of CONSTANT MONOMIAL_VECs.
uint8_t dof_id_type
Definition: id_types.h:67
std::map< std::string, std::unique_ptr< System >, std::less<> > _systems
Data structure holding the systems.
static bool redundant_added_side(const Elem &elem, unsigned int side)
const RemoteElem * remote_elem
Definition: remote_elem.C:57