libMesh
equation_systems.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 // System includes
20 #include <sstream>
21 
22 // Local Includes
23 #include "libmesh/default_coupling.h" // For downconversion
24 #include "libmesh/explicit_system.h"
25 #include "libmesh/fe_interface.h"
26 #include "libmesh/frequency_system.h"
27 #include "libmesh/linear_implicit_system.h"
28 #include "libmesh/mesh_refinement.h"
29 #include "libmesh/newmark_system.h"
30 #include "libmesh/nonlinear_implicit_system.h"
31 #include "libmesh/rb_construction.h"
32 #include "libmesh/transient_rb_construction.h"
33 #include "libmesh/eigen_system.h"
34 #include "libmesh/parallel.h"
35 #include "libmesh/transient_system.h"
36 #include "libmesh/dof_map.h"
37 #include "libmesh/mesh_base.h"
38 #include "libmesh/elem.h"
39 #include "libmesh/libmesh_logging.h"
40 
41 // Include the systems before this one to avoid
42 // overlapping forward declarations.
43 #include "libmesh/equation_systems.h"
44 
45 namespace libMesh
46 {
47 
48 // Forward Declarations
49 
50 
51 
52 
53 // ------------------------------------------------------------
54 // EquationSystems class implementation
56  ParallelObject (m),
57  _mesh (m),
58  _refine_in_reinit(true),
59  _enable_default_ghosting(true)
60 {
61  // Set default parameters
62  this->parameters.set<Real> ("linear solver tolerance") = TOLERANCE * TOLERANCE;
63  this->parameters.set<unsigned int>("linear solver maximum iterations") = 5000;
64 }
65 
66 
67 
69 {
70  this->clear ();
71 }
72 
73 
74 
76 {
77  // Clear any additional parameters
78  parameters.clear ();
79 
80  // clear the systems. We must delete them
81  // since we newed them!
82  while (!_systems.empty())
83  {
84  system_iterator pos = _systems.begin();
85 
86  System * sys = pos->second;
87  delete sys;
88  sys = nullptr;
89 
90  _systems.erase (pos);
91  }
92 }
93 
94 
95 
97 {
98  const unsigned int n_sys = this->n_systems();
99 
100  libmesh_assert_not_equal_to (n_sys, 0);
101 
102  // Tell all the \p DofObject entities how many systems
103  // there are.
104  for (auto & node : _mesh.node_ptr_range())
105  node->set_n_systems(n_sys);
106 
107  for (auto & elem : _mesh.element_ptr_range())
108  elem->set_n_systems(n_sys);
109 
110  for (unsigned int i=0; i != this->n_systems(); ++i)
111  this->get_system(i).init();
112 
113 #ifdef LIBMESH_ENABLE_AMR
114  MeshRefinement mesh_refine(_mesh);
115  mesh_refine.clean_refinement_flags();
116 #endif
117 }
118 
119 
120 
122 {
123  const bool mesh_changed = this->reinit_solutions();
124 
125  // If the mesh has changed, systems will need to reinitialize their
126  // own data on the new mesh.
127  if (mesh_changed)
128  this->reinit_systems();
129 }
130 
131 
132 
134 {
135  parallel_object_only();
136 
137  const unsigned int n_sys = this->n_systems();
138  libmesh_assert_not_equal_to (n_sys, 0);
139 
140  // We may have added new systems since our last
141  // EquationSystems::(re)init call
142  bool _added_new_systems = false;
143  for (unsigned int i=0; i != n_sys; ++i)
144  if (!this->get_system(i).is_initialized())
145  _added_new_systems = true;
146 
147  if (_added_new_systems)
148  {
149  // Our DofObjects will need space for the additional systems
150  for (auto & node : _mesh.node_ptr_range())
151  node->set_n_systems(n_sys);
152 
153  for (auto & elem : _mesh.element_ptr_range())
154  elem->set_n_systems(n_sys);
155 
156  // And any new systems will need initialization
157  for (unsigned int i=0; i != n_sys; ++i)
158  if (!this->get_system(i).is_initialized())
159  this->get_system(i).init();
160  }
161 
162 
163  // We used to assert that all nodes and elements *already* had
164  // n_systems() properly set; however this is false in the case where
165  // user code has manually added nodes and/or elements to an
166  // already-initialized system.
167 
168  // Make sure all the \p DofObject entities know how many systems
169  // there are.
170  {
171  // All the nodes
172  for (auto & node : _mesh.node_ptr_range())
173  node->set_n_systems(this->n_systems());
174 
175  // All the elements
176  for (auto & elem : _mesh.element_ptr_range())
177  elem->set_n_systems(this->n_systems());
178  }
179 
180  // Localize each system's vectors
181  for (unsigned int i=0; i != this->n_systems(); ++i)
182  this->get_system(i).re_update();
183 
184 #ifdef LIBMESH_ENABLE_AMR
185 
186  bool mesh_changed = false;
187 
188  // FIXME: For backwards compatibility, assume
189  // refine_and_coarsen_elements or refine_uniformly have already
190  // been called
191  {
192  for (unsigned int i=0; i != this->n_systems(); ++i)
193  {
194  System & sys = this->get_system(i);
195 
196  // Even if the system doesn't have any variables in it we want
197  // consistent behavior; e.g. distribute_dofs should have the
198  // opportunity to count up zero dofs on each processor.
199  //
200  // Who's been adding zero-var systems anyway, outside of my
201  // unit tests? - RHS
202  // if (!sys.n_vars())
203  // continue;
204 
206 
207  // Recreate any user or internal constraints
208  sys.reinit_constraints();
209 
210  sys.prolong_vectors();
211  }
212  mesh_changed = true;
213  }
214 
215  if (this->_refine_in_reinit)
216  {
217  // Don't override any user refinement settings
218  MeshRefinement mesh_refine(_mesh);
219  mesh_refine.face_level_mismatch_limit() = 0; // unlimited
220  mesh_refine.overrefined_boundary_limit() = -1; // unlimited
221  mesh_refine.underrefined_boundary_limit() = -1; // unlimited
222 
223  // Try to coarsen the mesh, then restrict each system's vectors
224  // if necessary
225  if (mesh_refine.coarsen_elements())
226  {
227  for (unsigned int i=0; i != this->n_systems(); ++i)
228  {
229  System & sys = this->get_system(i);
231  sys.reinit_constraints();
232  sys.restrict_vectors();
233  }
234  mesh_changed = true;
235  }
236 
237  // Once vectors are all restricted, we can delete
238  // children of coarsened elements
239  if (mesh_changed)
240  this->get_mesh().contract();
241 
242  // Try to refine the mesh, then prolong each system's vectors
243  // if necessary
244  if (mesh_refine.refine_elements())
245  {
246  for (unsigned int i=0; i != this->n_systems(); ++i)
247  {
248  System & sys = this->get_system(i);
250  sys.reinit_constraints();
251  sys.prolong_vectors();
252  }
253  mesh_changed = true;
254  }
255  }
256 
257  return mesh_changed;
258 
259 #endif // #ifdef LIBMESH_ENABLE_AMR
260 
261  return false;
262 }
263 
264 
265 
267 {
268  for (unsigned int i=0; i != this->n_systems(); ++i)
269  this->get_system(i).reinit();
270 }
271 
272 
273 
275 {
276  // A serial mesh means nothing needs to be done
277  if (_mesh.is_serial())
278  return;
279 
280  const unsigned int n_sys = this->n_systems();
281 
282  libmesh_assert_not_equal_to (n_sys, 0);
283 
284  // Gather the mesh
285  _mesh.allgather();
286 
287  // Tell all the \p DofObject entities how many systems
288  // there are.
289  for (auto & node : _mesh.node_ptr_range())
290  node->set_n_systems(n_sys);
291 
292  for (auto & elem : _mesh.element_ptr_range())
293  elem->set_n_systems(n_sys);
294 
295  // And distribute each system's dofs
296  for (unsigned int i=0; i != this->n_systems(); ++i)
297  {
298  System & sys = this->get_system(i);
299  DofMap & dof_map = sys.get_dof_map();
300  dof_map.distribute_dofs(_mesh);
301 
302  // The user probably won't need constraint equations or the
303  // send_list after an allgather, but let's keep it in consistent
304  // shape just in case.
305  sys.reinit_constraints();
306  dof_map.prepare_send_list();
307  }
308 }
309 
310 
311 
313 {
314  _enable_default_ghosting = enable;
315  MeshBase &mesh = this->get_mesh();
316 
317  if (enable)
318  mesh.add_ghosting_functor(mesh.default_ghosting());
319  else
320  mesh.remove_ghosting_functor(mesh.default_ghosting());
321 
322  for (unsigned int i=0; i != this->n_systems(); ++i)
323  {
324  DofMap & dof_map = this->get_system(i).get_dof_map();
325  if (enable)
326  dof_map.add_default_ghosting();
327  else
328  dof_map.remove_default_ghosting();
329  }
330 }
331 
332 
333 
335 {
336  LOG_SCOPE("update()", "EquationSystems");
337 
338  // Localize each system's vectors
339  for (unsigned int i=0; i != this->n_systems(); ++i)
340  this->get_system(i).update();
341 }
342 
343 
344 
345 System & EquationSystems::add_system (const std::string & sys_type,
346  const std::string & name)
347 {
348  // If the user already built a system with this name, we'll
349  // trust them and we'll use it. That way they can pre-add
350  // non-standard derived system classes, and if their restart file
351  // has some non-standard sys_type we won't throw an error.
352  if (_systems.count(name))
353  {
354  return this->get_system(name);
355  }
356  // Build a basic System
357  else if (sys_type == "Basic")
358  this->add_system<System> (name);
359 
360  // Build a Newmark system
361  else if (sys_type == "Newmark")
362  this->add_system<NewmarkSystem> (name);
363 
364  // Build an Explicit system
365  else if ((sys_type == "Explicit"))
366  this->add_system<ExplicitSystem> (name);
367 
368  // Build an Implicit system
369  else if ((sys_type == "Implicit") ||
370  (sys_type == "Steady" ))
371  this->add_system<ImplicitSystem> (name);
372 
373  // build a transient implicit linear system
374  else if ((sys_type == "Transient") ||
375  (sys_type == "TransientImplicit") ||
376  (sys_type == "TransientLinearImplicit"))
377  this->add_system<TransientLinearImplicitSystem> (name);
378 
379  // build a transient implicit nonlinear system
380  else if (sys_type == "TransientNonlinearImplicit")
381  this->add_system<TransientNonlinearImplicitSystem> (name);
382 
383  // build a transient explicit system
384  else if (sys_type == "TransientExplicit")
385  this->add_system<TransientExplicitSystem> (name);
386 
387  // build a linear implicit system
388  else if (sys_type == "LinearImplicit")
389  this->add_system<LinearImplicitSystem> (name);
390 
391  // build a nonlinear implicit system
392  else if (sys_type == "NonlinearImplicit")
393  this->add_system<NonlinearImplicitSystem> (name);
394 
395  // build a Reduced Basis Construction system
396  else if (sys_type == "RBConstruction")
397  this->add_system<RBConstruction> (name);
398 
399  // build a transient Reduced Basis Construction system
400  else if (sys_type == "TransientRBConstruction")
401  this->add_system<TransientRBConstruction> (name);
402 
403 #ifdef LIBMESH_HAVE_SLEPC
404  // build an eigen system
405  else if (sys_type == "Eigen")
406  this->add_system<EigenSystem> (name);
407  else if (sys_type == "TransientEigenSystem")
408  this->add_system<TransientEigenSystem> (name);
409 #endif
410 
411 #if defined(LIBMESH_USE_COMPLEX_NUMBERS)
412  // build a frequency system
413  else if (sys_type == "Frequency")
414  this->add_system<FrequencySystem> (name);
415 #endif
416 
417  else
418  libmesh_error_msg("ERROR: Unknown system type: " << sys_type);
419 
420  // Return a reference to the new system
421  //return (*this)(name);
422  return this->get_system(name);
423 }
424 
425 
426 
427 
428 
429 
430 #ifdef LIBMESH_ENABLE_DEPRECATED
431 void EquationSystems::delete_system (const std::string & name)
432 {
433  libmesh_deprecated();
434 
435  if (!_systems.count(name))
436  libmesh_error_msg("ERROR: no system named " << name);
437 
438  delete _systems[name];
439 
440  _systems.erase (name);
441 }
442 #endif
443 
444 
445 
447 {
448  libmesh_assert (this->n_systems());
449 
450  for (unsigned int i=0; i != this->n_systems(); ++i)
451  this->get_system(i).solve();
452 }
453 
454 
455 
457 {
458  libmesh_assert (this->n_systems());
459 
460  for (unsigned int i=0; i != this->n_systems(); ++i)
461  this->get_system(i).sensitivity_solve(parameters_in);
462 }
463 
464 
465 
466 void EquationSystems::adjoint_solve (const QoISet & qoi_indices)
467 {
468  libmesh_assert (this->n_systems());
469 
470  for (unsigned int i=this->n_systems(); i != 0; --i)
471  this->get_system(i-1).adjoint_solve(qoi_indices);
472 }
473 
474 
475 
476 void EquationSystems::build_variable_names (std::vector<std::string> & var_names,
477  const FEType * type,
478  const std::set<std::string> * system_names) const
479 {
480  unsigned int var_num=0;
481 
482  const_system_iterator pos = _systems.begin();
483  const const_system_iterator end = _systems.end();
484 
485  // Need to size var_names by scalar variables plus all the
486  // vector components for all the vector variables
487  //Could this be replaced by a/some convenience methods?[PB]
488  {
489  unsigned int n_scalar_vars = 0;
490  unsigned int n_vector_vars = 0;
491 
492  for (; pos != end; ++pos)
493  {
494  // Check current system is listed in system_names, and skip pos if not
495  bool use_current_system = (system_names == nullptr);
496  if (!use_current_system)
497  use_current_system = system_names->count(pos->first);
498  if (!use_current_system || pos->second->hide_output())
499  continue;
500 
501  for (auto vn : IntRange<unsigned int>(0, pos->second->n_vars()))
502  {
503  if (FEInterface::field_type(pos->second->variable_type(vn)) == TYPE_VECTOR)
504  n_vector_vars++;
505  else
506  n_scalar_vars++;
507  }
508  }
509 
510  // Here, we're assuming the number of vector components is the same
511  // as the mesh dimension. Will break for mixed dimension meshes.
512  unsigned int dim = this->get_mesh().mesh_dimension();
513  unsigned int nv = n_scalar_vars + dim*n_vector_vars;
514 
515  // We'd better not have more than dim*his->n_vars() (all vector variables)
516  libmesh_assert_less_equal ( nv, dim*this->n_vars() );
517 
518  // Here, we're assuming the number of vector components is the same
519  // as the mesh dimension. Will break for mixed dimension meshes.
520 
521  var_names.resize( nv );
522  }
523 
524  // reset
525  pos = _systems.begin();
526 
527  for (; pos != end; ++pos)
528  {
529  // Check current system is listed in system_names, and skip pos if not
530  bool use_current_system = (system_names == nullptr);
531  if (!use_current_system)
532  use_current_system = system_names->count(pos->first);
533  if (!use_current_system || pos->second->hide_output())
534  continue;
535 
536  for (auto vn : IntRange<unsigned int>(0, pos->second->n_vars()))
537  {
538  const std::string & var_name = pos->second->variable_name(vn);
539  const FEType & fe_type = pos->second->variable_type(vn);
540 
541  unsigned int n_vec_dim = FEInterface::n_vec_dim( pos->second->get_mesh(), fe_type);
542 
543  // Filter on the type if requested
544  if (type == nullptr || (type && *type == fe_type))
545  {
546  if (FEInterface::field_type(fe_type) == TYPE_VECTOR)
547  {
548  switch(n_vec_dim)
549  {
550  case 0:
551  case 1:
552  var_names[var_num++] = var_name;
553  break;
554  case 2:
555  var_names[var_num++] = var_name+"_x";
556  var_names[var_num++] = var_name+"_y";
557  break;
558  case 3:
559  var_names[var_num++] = var_name+"_x";
560  var_names[var_num++] = var_name+"_y";
561  var_names[var_num++] = var_name+"_z";
562  break;
563  default:
564  libmesh_error_msg("Invalid dim in build_variable_names");
565  }
566  }
567  else
568  var_names[var_num++] = var_name;
569  }
570  }
571  }
572  // Now resize again in case we filtered any names
573  var_names.resize(var_num);
574 }
575 
576 
577 
578 void EquationSystems::build_solution_vector (std::vector<Number> &,
579  const std::string &,
580  const std::string &) const
581 {
582  // TODO:[BSK] re-implement this from the method below
583  libmesh_not_implemented();
584 }
585 
586 
587 
588 
589 std::unique_ptr<NumericVector<Number>>
590 EquationSystems::build_parallel_solution_vector(const std::set<std::string> * system_names) const
591 {
592  LOG_SCOPE("build_parallel_solution_vector()", "EquationSystems");
593 
594  // This function must be run on all processors at once
595  parallel_object_only();
596 
597  const unsigned int dim = _mesh.mesh_dimension();
598  const dof_id_type max_nn = _mesh.max_node_id();
599 
600  // allocate vector storage to hold
601  // (max_node_id)*(number_of_variables) entries.
602  //
603  // If node renumbering is disabled and adaptive coarsening has
604  // created gaps between node numbers, then this vector will be
605  // sparse.
606  //
607  // We have to differentiate between between scalar and vector
608  // variables. We intercept vector variables and treat each
609  // component as a scalar variable (consistently with build_solution_names).
610 
611  unsigned int nv = 0;
612 
613  //Could this be replaced by a/some convenience methods?[PB]
614  {
615  unsigned int n_scalar_vars = 0;
616  unsigned int n_vector_vars = 0;
617  const_system_iterator pos = _systems.begin();
618  const const_system_iterator end = _systems.end();
619 
620  for (; pos != end; ++pos)
621  {
622  // Check current system is listed in system_names, and skip pos if not
623  bool use_current_system = (system_names == nullptr);
624  if (!use_current_system)
625  use_current_system = system_names->count(pos->first);
626  if (!use_current_system || pos->second->hide_output())
627  continue;
628 
629  for (auto vn : IntRange<unsigned int>(0, pos->second->n_vars()))
630  {
631  if (FEInterface::field_type(pos->second->variable_type(vn)) == TYPE_VECTOR)
632  n_vector_vars++;
633  else
634  n_scalar_vars++;
635  }
636  }
637  // Here, we're assuming the number of vector components is the same
638  // as the mesh dimension. Will break for mixed dimension meshes.
639  nv = n_scalar_vars + dim*n_vector_vars;
640  }
641 
642  // Get the number of nodes to store locally.
643  dof_id_type n_local_nodes = cast_int<dof_id_type>
646 
647  // If node renumbering has been disabled, nodes may not be numbered
648  // contiguously, and the number of nodes might not match the
649  // max_node_id. In this case we just do our best.
650  dof_id_type n_total_nodes = n_local_nodes;
651  _mesh.comm().sum(n_total_nodes);
652 
653  const dof_id_type n_gaps = max_nn - n_total_nodes;
654  const dof_id_type gaps_per_processor = n_gaps / _mesh.comm().size();
655  const dof_id_type remainder_gaps = n_gaps % _mesh.comm().size();
656 
657  n_local_nodes = n_local_nodes + // Actual nodes
658  gaps_per_processor + // Our even share of gaps
659  (_mesh.comm().rank() < remainder_gaps); // Leftovers
660 
661  // Create a NumericVector to hold the parallel solution
662  std::unique_ptr<NumericVector<Number>> parallel_soln_ptr = NumericVector<Number>::build(_communicator);
663  NumericVector<Number> & parallel_soln = *parallel_soln_ptr;
664  parallel_soln.init(max_nn*nv, n_local_nodes*nv, false, PARALLEL);
665 
666  // Create a NumericVector to hold the "repeat_count" for each node - this is essentially
667  // the number of elements contributing to that node's value
668  std::unique_ptr<NumericVector<Number>> repeat_count_ptr = NumericVector<Number>::build(_communicator);
669  NumericVector<Number> & repeat_count = *repeat_count_ptr;
670  repeat_count.init(max_nn*nv, n_local_nodes*nv, false, PARALLEL);
671 
672  repeat_count.close();
673 
674  unsigned int var_num=0;
675 
676  // For each system in this EquationSystems object,
677  // update the global solution and if we are on processor 0,
678  // loop over the elements and build the nodal solution
679  // from the element solution. Then insert this nodal solution
680  // into the vector passed to build_solution_vector.
681  const_system_iterator pos = _systems.begin();
682  const const_system_iterator end = _systems.end();
683 
684  for (; pos != end; ++pos)
685  {
686  // Check current system is listed in system_names, and skip pos if not
687  bool use_current_system = (system_names == nullptr);
688  if (!use_current_system)
689  use_current_system = system_names->count(pos->first);
690  if (!use_current_system || pos->second->hide_output())
691  continue;
692 
693  const System & system = *(pos->second);
694  const unsigned int nv_sys = system.n_vars();
695  const unsigned int sys_num = system.number();
696 
697  //Could this be replaced by a/some convenience methods?[PB]
698  unsigned int n_scalar_vars = 0;
699  unsigned int n_vector_vars = 0;
700  for (auto vn : IntRange<unsigned int>(0, pos->second->n_vars()))
701  {
702  if (FEInterface::field_type(pos->second->variable_type(vn)) == TYPE_VECTOR)
703  n_vector_vars++;
704  else
705  n_scalar_vars++;
706  }
707 
708  // Here, we're assuming the number of vector components is the same
709  // as the mesh dimension. Will break for mixed dimension meshes.
710  unsigned int nv_sys_split = n_scalar_vars + dim*n_vector_vars;
711 
712  // Update the current_local_solution
713  {
714  System & non_const_sys = const_cast<System &>(system);
715  // We used to simply call non_const_sys.solution->close()
716  // here, but that is not allowed when the solution vector is
717  // locked read-only, for example when printing the solution
718  // during the middle of a solve... So try to be a bit
719  // more careful about calling close() unnecessarily.
720  libmesh_assert(this->comm().verify(non_const_sys.solution->closed()));
721  if (!non_const_sys.solution->closed())
722  non_const_sys.solution->close();
723  non_const_sys.update();
724  }
725 
726  NumericVector<Number> & sys_soln(*system.current_local_solution);
727 
728  std::vector<Number> elem_soln; // The finite element solution
729  std::vector<Number> nodal_soln; // The FE solution interpolated to the nodes
730  std::vector<dof_id_type> dof_indices; // The DOF indices for the finite element
731 
732  unsigned var_inc = 0;
733  for (unsigned int var=0; var<nv_sys; var++)
734  {
735  const FEType & fe_type = system.variable_type(var);
736  const Variable & var_description = system.variable(var);
737  const DofMap & dof_map = system.get_dof_map();
738 
739  unsigned int n_vec_dim = FEInterface::n_vec_dim( pos->second->get_mesh(), fe_type );
740 
741  for (const auto & elem : _mesh.active_local_element_ptr_range())
742  {
743  if (var_description.active_on_subdomain(elem->subdomain_id()))
744  {
745  dof_map.dof_indices (elem, dof_indices, var);
746 
747  elem_soln.resize(dof_indices.size());
748 
749  for (auto i : index_range(dof_indices))
750  elem_soln[i] = sys_soln(dof_indices[i]);
751 
753  fe_type,
754  elem,
755  elem_soln,
756  nodal_soln);
757 
758 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
759  // infinite elements should be skipped...
760  if (!elem->infinite())
761 #endif
762  {
763  libmesh_assert_equal_to (nodal_soln.size(), n_vec_dim*elem->n_nodes());
764 
765  for (auto n : elem->node_index_range())
766  {
767  for (unsigned int d=0; d < n_vec_dim; d++)
768  {
769  // For vector-valued elements, all components are in nodal_soln. For each
770  // node, the components are stored in order, i.e. node_0 -> s0_x, s0_y, s0_z
771  parallel_soln.add(nv*(elem->node_id(n)) + (var_inc+d + var_num), nodal_soln[n_vec_dim*n+d]);
772 
773  // Increment the repeat count for this position
774  repeat_count.add(nv*(elem->node_id(n)) + (var_inc+d + var_num), 1);
775  }
776  }
777  }
778  }
779  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:
780  for (const Node & node : elem->node_ref_range())
781  // Only do this if this variable has NO DoFs at this node... it might have some from an adjoining element...
782  if (!node.n_dofs(sys_num, var))
783  for (unsigned int d=0; d < n_vec_dim; d++)
784  repeat_count.add(nv*node.id() + (var_inc+d + var_num), 1);
785 
786  } // end loop over elements
787  var_inc += n_vec_dim;
788  } // end loop on variables in this system
789 
790  var_num += nv_sys_split;
791  } // end loop over systems
792 
793  // Sum the nodal solution values and repeat counts.
794  parallel_soln.close();
795  repeat_count.close();
796 
797  // If there were gaps in the node numbering, there will be
798  // corresponding zeros in the parallel_soln and repeat_count
799  // vectors. We need to set those repeat_count entries to 1
800  // in order to avoid dividing by zero.
801  if (n_gaps)
802  {
803  for (numeric_index_type i=repeat_count.first_local_index();
804  i<repeat_count.last_local_index(); ++i)
805  {
806  // repeat_count entries are integral values but let's avoid a
807  // direct floating point comparison with 0 just in case some
808  // roundoff noise crept in during vector assembly?
809  if (std::abs(repeat_count(i)) < TOLERANCE)
810  repeat_count.set(i, 1.);
811  }
812 
813  // Make sure the repeat_count vector is up-to-date on all
814  // processors.
815  repeat_count.close();
816  }
817 
818  // Divide to get the average value at the nodes
819  parallel_soln /= repeat_count;
820 
821  return std::unique_ptr<NumericVector<Number>>(parallel_soln_ptr.release());
822 }
823 
824 
825 
826 void EquationSystems::build_solution_vector (std::vector<Number> & soln,
827  const std::set<std::string> * system_names) const
828 {
829  LOG_SCOPE("build_solution_vector()", "EquationSystems");
830 
831  // Call the parallel implementation
832  std::unique_ptr<NumericVector<Number>> parallel_soln =
833  this->build_parallel_solution_vector(system_names);
834 
835  // Localize the NumericVector into the provided std::vector.
836  parallel_soln->localize_to_one(soln);
837 }
838 
839 
840 
841 void EquationSystems::get_vars_active_subdomains(const std::vector<std::string> & names,
842  std::vector<std::set<subdomain_id_type>> & vars_active_subdomains) const
843 {
844  unsigned int var_num=0;
845 
846  vars_active_subdomains.clear();
847  vars_active_subdomains.resize(names.size());
848 
849  const_system_iterator pos = _systems.begin();
850  const const_system_iterator end = _systems.end();
851 
852  for (; pos != end; ++pos)
853  {
854  for (auto vn : IntRange<unsigned int>(0, pos->second->n_vars()))
855  {
856  const std::string & var_name = pos->second->variable_name(vn);
857 
858  auto names_it = std::find(names.begin(), names.end(), var_name);
859  if(names_it != names.end())
860  {
861  const Variable & variable = pos->second->variable(vn);
862  const std::set<subdomain_id_type> & active_subdomains = variable.active_subdomains();
863  vars_active_subdomains[var_num++] = active_subdomains;
864  }
865  }
866  }
867 
868  libmesh_assert_equal_to(var_num, names.size());
869 }
870 
871 
872 
873 void EquationSystems::get_solution (std::vector<Number> & soln,
874  std::vector<std::string> & names) const
875 {
876  libmesh_deprecated();
877  this->build_elemental_solution_vector(soln, names);
878 }
879 
880 
881 
882 void
884  std::vector<std::string> & names) const
885 {
886  // Call the parallel version of this function
887  std::unique_ptr<NumericVector<Number>> parallel_soln =
889 
890  // Localize into 'soln', provided that parallel_soln is not empty.
891  // Note: parallel_soln will be empty in the event that none of the
892  // input names were CONSTANT, MONOMIAL variables or there were
893  // simply no CONSTANT, MONOMIAL variables in the EquationSystems
894  // object.
895  soln.clear();
896  if (parallel_soln)
897  parallel_soln->localize_to_one(soln);
898 }
899 
900 
901 
902 std::vector<std::pair<unsigned int, unsigned int>>
904  (std::vector<std::string> & names, const FEType * type) const
905 {
906  // This function must be run on all processors at once
907  parallel_object_only();
908 
909  libmesh_assert (this->n_systems());
910 
911  // If the names vector has entries, we will only populate the soln vector
912  // with names included in that list. Note: The names vector may be
913  // reordered upon exiting this function
914  std::vector<std::pair<unsigned int, unsigned int>> var_nums;
915  std::vector<std::string> filter_names = names;
916  bool is_filter_names = !filter_names.empty();
917 
918  names.clear();
919 
920  const_system_iterator pos = _systems.begin();
921  const const_system_iterator end = _systems.end();
922  unsigned sys_ctr = 0;
923 
924  for (; pos != end; ++pos, ++sys_ctr)
925  {
926  const System & system = *(pos->second);
927  const unsigned int nv_sys = system.n_vars();
928 
929  for (unsigned int var=0; var < nv_sys; ++var)
930  {
931  const std::string & name = system.variable_name(var);
932  if ((type && system.variable_type(var) != *type) ||
933  (is_filter_names && std::find(filter_names.begin(), filter_names.end(), name) == filter_names.end()))
934  continue;
935 
936  // Otherwise, this variable should be output
937  var_nums.push_back
938  (std::make_pair(system.number(), var));
939  }
940  }
941 
942  std::sort(var_nums.begin(), var_nums.end());
943 
944  for (const auto & var_num : var_nums)
945  {
946  const std::string & name =
947  this->get_system(var_num.first).variable_name(var_num.second);
948  if (names.empty() || names.back() != name)
949  names.push_back(name);
950  }
951 
952  return var_nums;
953 }
954 
955 
956 std::unique_ptr<NumericVector<Number>>
957 EquationSystems::build_parallel_elemental_solution_vector (std::vector<std::string> & names) const
958 {
959  FEType type(CONSTANT, MONOMIAL);
960  std::vector<std::pair<unsigned int, unsigned int>> var_nums =
961  this->find_variable_numbers(names, &type);
962 
963  const std::size_t nv = var_nums.size();
964  const dof_id_type ne = _mesh.n_elem();
965  libmesh_assert_equal_to (ne, _mesh.max_elem_id());
966 
967  // If there are no variables to write out don't do anything...
968  if (!nv)
969  return std::unique_ptr<NumericVector<Number>>(nullptr);
970 
971  // We can handle the case where there are nullptrs in the Elem vector
972  // by just having extra zeros in the solution vector.
973  numeric_index_type parallel_soln_global_size = ne*nv;
974 
975  numeric_index_type div = parallel_soln_global_size / this->n_processors();
976  numeric_index_type mod = parallel_soln_global_size % this->n_processors();
977 
978  // Initialize all processors to the average size.
979  numeric_index_type parallel_soln_local_size = div;
980 
981  // The first "mod" processors get an extra entry.
982  if (this->processor_id() < mod)
983  parallel_soln_local_size = div+1;
984 
985  // Create a NumericVector to hold the parallel solution
986  std::unique_ptr<NumericVector<Number>> parallel_soln_ptr = NumericVector<Number>::build(_communicator);
987  NumericVector<Number> & parallel_soln = *parallel_soln_ptr;
988  parallel_soln.init(parallel_soln_global_size,
989  parallel_soln_local_size,
990  /*fast=*/false,
991  /*ParallelType=*/PARALLEL);
992 
993  unsigned int sys_ctr = 0;
994 
995  // For each system in this EquationSystems object,
996  // update the global solution and collect the
997  // CONSTANT MONOMIALs. The entries are in variable-major
998  // format.
999  for (auto i : index_range(var_nums))
1000  {
1001  std::pair<unsigned int, unsigned int> var_num = var_nums[i];
1002  const System & system = this->get_system(var_num.first);
1003 
1004  // Update the current_local_solution if necessary
1005  if (sys_ctr != var_num.first)
1006  {
1007  System & non_const_sys = const_cast<System &>(system);
1008  // We used to simply call non_const_sys.solution->close()
1009  // here, but that is not allowed when the solution vector is
1010  // locked read-only, for example when printing the solution
1011  // during during the middle of a solve... So try to be a bit
1012  // more careful about calling close() unnecessarily.
1013  libmesh_assert(this->comm().verify(non_const_sys.solution->closed()));
1014  if (!non_const_sys.solution->closed())
1015  non_const_sys.solution->close();
1016  non_const_sys.update();
1017  sys_ctr = var_num.first;
1018  }
1019 
1020  NumericVector<Number> & sys_soln(*system.current_local_solution);
1021 
1022  // The DOF indices for the finite element
1023  std::vector<dof_id_type> dof_indices;
1024 
1025  const unsigned int var = var_num.second;
1026 
1027  const Variable & variable = system.variable(var);
1028  const DofMap & dof_map = system.get_dof_map();
1029 
1030  for (const auto & elem : _mesh.active_local_element_ptr_range())
1031  if (variable.active_on_subdomain(elem->subdomain_id()))
1032  {
1033  dof_map.dof_indices (elem, dof_indices, var);
1034 
1035  libmesh_assert_equal_to (1, dof_indices.size());
1036 
1037  parallel_soln.set((ne*i)+elem->id(), sys_soln(dof_indices[0]));
1038  }
1039  } // end loop over var_nums
1040 
1041  parallel_soln.close();
1042  return std::unique_ptr<NumericVector<Number>>(parallel_soln_ptr.release());
1043 }
1044 
1045 
1046 
1047 void
1049 (std::vector<Number> & soln,
1050  const std::set<std::string> * system_names,
1051  const std::vector<std::string> * var_names,
1052  bool vertices_only) const
1053 {
1054  LOG_SCOPE("build_discontinuous_solution_vector()", "EquationSystems");
1055 
1056  libmesh_assert (this->n_systems());
1057 
1058  const unsigned int dim = _mesh.mesh_dimension();
1059 
1060  // Get the number of variables (nv) by counting the number of variables
1061  // in each system listed in system_names
1062  unsigned int nv = 0;
1063 
1064  for (const auto & pr : _systems)
1065  {
1066  // Check current system is listed in system_names, and skip pos if not
1067  bool use_current_system = (system_names == nullptr);
1068  if (!use_current_system)
1069  use_current_system = system_names->count(pr.first);
1070  if (!use_current_system || pr.second->hide_output())
1071  continue;
1072 
1073  const System * system = pr.second;
1074 
1075  // Loop over all variables in this System and check whether we
1076  // are supposed to use each one.
1077  for (auto var_id : IntRange<unsigned int>(0, system->n_vars()))
1078  {
1079  bool use_current_var = (var_names == nullptr);
1080  if (!use_current_var)
1081  use_current_var = std::count(var_names->begin(),
1082  var_names->end(),
1083  system->variable_name(var_id));
1084 
1085  // Only increment the total number of vars if we are
1086  // supposed to use this one.
1087  if (use_current_var)
1088  nv++;
1089  }
1090  }
1091 
1092  // get the total weight
1093  unsigned int tw=0;
1094  for (const auto & elem : _mesh.active_element_ptr_range())
1095  tw += vertices_only ? elem->n_vertices() : elem->n_nodes();
1096 
1097  // Only if we are on processor zero, allocate the storage
1098  // to hold (number_of_nodes)*(number_of_variables) entries.
1099  if (_mesh.processor_id() == 0)
1100  soln.resize(tw*nv);
1101 
1102  std::vector<Number> sys_soln;
1103 
1104  // Keep track of the variable "offset". This is used for indexing
1105  // into the global solution vector.
1106  unsigned int var_offset = 0;
1107 
1108  // For each system in this EquationSystems object,
1109  // update the global solution and if we are on processor 0,
1110  // loop over the elements and build the nodal solution
1111  // from the element solution. Then insert this nodal solution
1112  // into the vector passed to build_solution_vector.
1113  for (const auto & pr : _systems)
1114  {
1115  // Check current system is listed in system_names, and skip pos if not
1116  bool use_current_system = (system_names == nullptr);
1117  if (!use_current_system)
1118  use_current_system = system_names->count(pr.first);
1119  if (!use_current_system || pr.second->hide_output())
1120  continue;
1121 
1122  const System * system = pr.second;
1123  const unsigned int nv_sys = system->n_vars();
1124 
1125  system->update_global_solution (sys_soln, 0);
1126 
1127  // Keep track of the number of vars actually written.
1128  unsigned int n_vars_written_current_system = 0;
1129 
1130  if (_mesh.processor_id() == 0)
1131  {
1132  std::vector<Number> soln_coeffs; // The finite element solution coeffs
1133  std::vector<Number> nodal_soln; // The FE solution interpolated to the nodes
1134  std::vector<dof_id_type> dof_indices; // The DOF indices for the finite element
1135 
1136  // For each variable, determine if we are supposed to
1137  // write it, then loop over the active elements, compute
1138  // the nodal_soln and store it to the "soln" vector. We
1139  // store zeros for subdomain-restricted variables on
1140  // elements where they are not active.
1141  for (unsigned int var=0; var<nv_sys; var++)
1142  {
1143  bool use_current_var = (var_names == nullptr);
1144  if (!use_current_var)
1145  use_current_var = std::count(var_names->begin(),
1146  var_names->end(),
1147  system->variable_name(var));
1148 
1149  // If we aren't supposed to write this var, go to the
1150  // next loop iteration.
1151  if (!use_current_var)
1152  continue;
1153 
1154  const FEType & fe_type = system->variable_type(var);
1155  const Variable & var_description = system->variable(var);
1156 
1157  unsigned int nn=0;
1158 
1159  for (auto & elem : _mesh.active_element_ptr_range())
1160  {
1161  if (var_description.active_on_subdomain(elem->subdomain_id()))
1162  {
1163  system->get_dof_map().dof_indices (elem, dof_indices, var);
1164 
1165  soln_coeffs.resize(dof_indices.size());
1166 
1167  for (auto i : index_range(dof_indices))
1168  soln_coeffs[i] = sys_soln[dof_indices[i]];
1169 
1170  // Compute the FE solution at all the nodes, but
1171  // only use the first n_vertices() entries if
1172  // vertices_only == true.
1174  fe_type,
1175  elem,
1176  soln_coeffs,
1177  nodal_soln);
1178 
1179 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1180  // infinite elements should be skipped...
1181  if (!elem->infinite())
1182 #endif
1183  {
1184  libmesh_assert_equal_to (nodal_soln.size(), elem->n_nodes());
1185 
1186  const unsigned int n_vals =
1187  vertices_only ? elem->n_vertices() : elem->n_nodes();
1188 
1189  for (unsigned int n=0; n<n_vals; n++)
1190  {
1191  // Compute index into global solution vector.
1192  std::size_t index =
1193  nv * (nn++) + (n_vars_written_current_system + var_offset);
1194 
1195  soln[index] += nodal_soln[n];
1196  }
1197  }
1198  }
1199  else
1200  nn += vertices_only ? elem->n_vertices() : elem->n_nodes();
1201  } // end loop over active elements
1202 
1203  // If we made it here, we actually wrote a variable, so increment
1204  // the number of variables actually written for the current system.
1205  n_vars_written_current_system++;
1206 
1207  } // end loop over vars
1208  } // end if proc 0
1209 
1210  // Update offset for next loop iteration.
1211  var_offset += n_vars_written_current_system;
1212  } // end loop over systems
1213 }
1214 
1215 
1216 
1218  const Real threshold,
1219  const bool verbose) const
1220 {
1221  // safety check, whether we handle at least the same number
1222  // of systems
1223  std::vector<bool> os_result;
1224 
1225  if (this->n_systems() != other_es.n_systems())
1226  {
1227  if (verbose)
1228  {
1229  libMesh::out << " Fatal difference. This system handles "
1230  << this->n_systems() << " systems," << std::endl
1231  << " while the other system handles "
1232  << other_es.n_systems()
1233  << " systems." << std::endl
1234  << " Aborting comparison." << std::endl;
1235  }
1236  return false;
1237  }
1238  else
1239  {
1240  // start comparing each system
1241  const_system_iterator pos = _systems.begin();
1242  const const_system_iterator end = _systems.end();
1243 
1244  for (; pos != end; ++pos)
1245  {
1246  const std::string & sys_name = pos->first;
1247  const System & system = *(pos->second);
1248 
1249  // get the other system
1250  const System & other_system = other_es.get_system (sys_name);
1251 
1252  os_result.push_back (system.compare (other_system, threshold, verbose));
1253 
1254  }
1255 
1256  }
1257 
1258 
1259  // sum up the results
1260  if (os_result.size()==0)
1261  return true;
1262  else
1263  {
1264  bool os_identical;
1265  unsigned int n = 0;
1266  do
1267  {
1268  os_identical = os_result[n];
1269  n++;
1270  }
1271  while (os_identical && n<os_result.size());
1272  return os_identical;
1273  }
1274 }
1275 
1276 
1277 
1278 std::string EquationSystems::get_info () const
1279 {
1280  std::ostringstream oss;
1281 
1282  oss << " EquationSystems\n"
1283  << " n_systems()=" << this->n_systems() << '\n';
1284 
1285  // Print the info for the individual systems
1286  const_system_iterator pos = _systems.begin();
1287  const const_system_iterator end = _systems.end();
1288 
1289  for (; pos != end; ++pos)
1290  oss << pos->second->get_info();
1291 
1292 
1293  // // Possibly print the parameters
1294  // if (!this->parameters.empty())
1295  // {
1296  // oss << " n_parameters()=" << this->n_parameters() << '\n';
1297  // oss << " Parameters:\n";
1298 
1299  // for (const auto & pr : _parameters)
1300  // oss << " "
1301  // << "\""
1302  // << pr.first
1303  // << "\""
1304  // << "="
1305  // << pr.second
1306  // << '\n';
1307  // }
1308 
1309  return oss.str();
1310 }
1311 
1312 
1313 
1314 void EquationSystems::print_info (std::ostream & os) const
1315 {
1316  os << this->get_info()
1317  << std::endl;
1318 }
1319 
1320 
1321 
1322 std::ostream & operator << (std::ostream & os,
1323  const EquationSystems & es)
1324 {
1325  es.print_info(os);
1326  return os;
1327 }
1328 
1329 
1330 
1331 unsigned int EquationSystems::n_vars () const
1332 {
1333  unsigned int tot=0;
1334 
1335  const_system_iterator pos = _systems.begin();
1336  const const_system_iterator end = _systems.end();
1337 
1338  for (; pos != end; ++pos)
1339  tot += pos->second->n_vars();
1340 
1341  return tot;
1342 }
1343 
1344 
1345 
1346 std::size_t EquationSystems::n_dofs () const
1347 {
1348  std::size_t tot=0;
1349 
1350  const_system_iterator pos = _systems.begin();
1351  const const_system_iterator end = _systems.end();
1352 
1353  for (; pos != end; ++pos)
1354  tot += pos->second->n_dofs();
1355 
1356  return tot;
1357 }
1358 
1359 
1360 
1361 
1363 {
1364  std::size_t tot=0;
1365 
1366  const_system_iterator pos = _systems.begin();
1367  const const_system_iterator end = _systems.end();
1368 
1369  for (; pos != end; ++pos)
1370  tot += pos->second->n_active_dofs();
1371 
1372  return tot;
1373 }
1374 
1375 
1377 {
1378  // All the nodes
1379  for (auto & node : _mesh.node_ptr_range())
1380  node->add_system();
1381 
1382  // All the elements
1383  for (auto & elem : _mesh.element_ptr_range())
1384  elem->add_system();
1385 }
1386 
1388 {
1389  this->get_system(sys_num).get_dof_map().remove_default_ghosting();
1390 }
1391 
1392 } // namespace libMesh
libMesh::System
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:100
libMesh::EquationSystems::solve
virtual void solve()
Call solve on all the individual equation systems.
Definition: equation_systems.C:446
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::EquationSystems::build_elemental_solution_vector
void build_elemental_solution_vector(std::vector< Number > &soln, std::vector< std::string > &names) const
Retrieve the solution data for CONSTANT MONOMIALs.
Definition: equation_systems.C:883
libMesh::System::n_vars
unsigned int n_vars() const
Definition: system.h:2155
libMesh::NumericVector::add
virtual void add(const numeric_index_type i, const T value)=0
Adds value to each entry of the vector.
libMesh::EquationSystems::add_system
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
Definition: equation_systems.C:345
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
libMesh::Parameters::clear
virtual void clear()
Clears internal data structures & frees any allocated memory.
Definition: parameters.h:313
libMesh::DofMap::prepare_send_list
void prepare_send_list()
Takes the _send_list vector (which may have duplicate entries) and sorts it.
Definition: dof_map.C:1651
libMesh::MeshRefinement::refine_elements
bool refine_elements()
Only refines the user-requested elements.
Definition: mesh_refinement.C:683
libMesh::FEInterface::n_vec_dim
static unsigned int n_vec_dim(const MeshBase &mesh, const FEType &fe_type)
Definition: fe_interface.C:1701
libMesh::DofMap::dof_indices
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:1967
libMesh::MeshBase::is_serial
virtual bool is_serial() const
Definition: mesh_base.h:159
libMesh::System::compare
virtual bool compare(const System &other_system, const Real threshold, const bool verbose) const
Definition: system.C:514
libMesh::PARALLEL
Definition: enum_parallel_type.h:36
libMesh::EquationSystems::allgather
void allgather()
Serializes a distributed mesh and its associated degree of freedom numbering for all systems.
Definition: equation_systems.C:274
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::MeshBase::local_nodes_begin
virtual node_iterator local_nodes_begin()=0
Iterate over local nodes (nodes whose processor_id() matches the current processor).
libMesh::MeshBase::n_elem
virtual dof_id_type n_elem() const =0
libMesh::System::variable_name
const std::string & variable_name(const unsigned int i) const
Definition: system.h:2203
libMesh::NumericVector::last_local_index
virtual numeric_index_type last_local_index() const =0
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh::MeshBase::max_elem_id
virtual dof_id_type max_elem_id() const =0
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::MeshRefinement::overrefined_boundary_limit
signed char & overrefined_boundary_limit()
If overrefined_boundary_limit is set to a nonnegative value, then refinement and coarsening will prod...
Definition: mesh_refinement.h:930
libMesh::operator<<
std::ostream & operator<<(std::ostream &os, const OrderWrapper &order)
Overload stream operators.
Definition: fe_type.h:164
libMesh::System::restrict_vectors
virtual void restrict_vectors()
Restrict vectors after the mesh has coarsened.
Definition: system.C:324
libMesh::NumericVector::close
virtual void close()=0
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
libMesh::NumericVector::init
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.
libMesh::EquationSystems::build_solution_vector
void build_solution_vector(std::vector< Number > &soln, const std::string &system_name, const std::string &variable_name="all_vars") const
Fill the input vector soln with the solution values for the system named name.
Definition: equation_systems.C:578
end
IterBase * end
Also have a polymorphic pointer to the end object, this prevents iterating past the end.
Definition: variant_filter_iterator.h:343
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::ParallelObject::comm
const Parallel::Communicator & comm() const
Definition: parallel_object.h:94
libMesh::EquationSystems::sensitivity_solve
virtual void sensitivity_solve(const ParameterVector &parameters)
Call sensitivity_solve on all the individual equation systems.
Definition: equation_systems.C:456
libMesh::System::number
unsigned int number() const
Definition: system.h:2075
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::MeshBase::max_node_id
virtual dof_id_type max_node_id() const =0
libMesh::System::variable
const Variable & variable(unsigned int var) const
Return a constant reference to Variable var.
Definition: system.h:2183
libMesh::Variable::active_subdomains
const std::set< subdomain_id_type > & active_subdomains() const
Definition: variable.h:150
libMesh::EquationSystems::EquationSystems
EquationSystems(MeshBase &mesh)
Constructor.
Definition: equation_systems.C:55
libMesh::MeshRefinement::face_level_mismatch_limit
unsigned char & face_level_mismatch_limit()
If face_level_mismatch_limit is set to a nonzero value, then refinement and coarsening will produce m...
Definition: mesh_refinement.h:915
libMesh::EquationSystems::_remove_default_ghosting
void _remove_default_ghosting(unsigned int sys_num)
This just calls DofMap::remove_default_ghosting() but using a shim lets us forward-declare DofMap.
Definition: equation_systems.C:1387
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::NumericVector::build
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
Builds a NumericVector on the processors in communicator comm using the linear solver package specifi...
Definition: numeric_vector.C:49
libMesh::FEInterface::nodal_soln
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)
Build the nodal soln from the element soln.
Definition: fe_interface.C:580
libMesh::NumericVector< Number >
libMesh::DofMap::add_default_ghosting
void add_default_ghosting()
Add the default functor(s) for coupling and algebraic ghosting.
Definition: dof_map.C:1829
libMesh::MeshBase::element_ptr_range
virtual SimpleRange< element_iterator > element_ptr_range()=0
libMesh::MeshRefinement
Implements (adaptive) mesh refinement algorithms for a MeshBase.
Definition: mesh_refinement.h:61
libMesh::EquationSystems::find_variable_numbers
std::vector< std::pair< unsigned int, unsigned int > > find_variable_numbers(std::vector< std::string > &names, const FEType *type=nullptr) const
Finds system and variable numbers for any variables of type corresponding to the entries in the input...
Definition: equation_systems.C:904
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::IntRange
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
libMesh::DofMap::distribute_dofs
void distribute_dofs(MeshBase &)
Distribute dofs on the current mesh.
Definition: dof_map.C:910
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::MeshRefinement::underrefined_boundary_limit
signed char & underrefined_boundary_limit()
If underrefined_boundary_limit is set to a nonnegative value, then refinement and coarsening will pro...
Definition: mesh_refinement.h:935
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::System::reinit_constraints
virtual void reinit_constraints()
Reinitializes the constraints for this system.
Definition: system.C:397
libMesh::EquationSystems::n_active_dofs
std::size_t n_active_dofs() const
Definition: equation_systems.C:1362
libMesh::ParallelObject::n_processors
processor_id_type n_processors() const
Definition: parallel_object.h:100
libMesh::MeshBase::node_ptr_range
virtual SimpleRange< node_iterator > node_ptr_range()=0
libMesh::QoISet
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
libMesh::EquationSystems::reinit_systems
virtual void reinit_systems()
Reinitialize all systems on the current mesh.
Definition: equation_systems.C:266
libMesh::EquationSystems::build_variable_names
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.
Definition: equation_systems.C:476
libMesh::System::prolong_vectors
virtual void prolong_vectors()
Prolong vectors after the mesh has refined.
Definition: system.C:380
libMesh::EquationSystems::reinit_solutions
bool reinit_solutions()
Handle any mesh changes and project any solutions onto the updated mesh.
Definition: equation_systems.C:133
libMesh::ParallelObject::processor_id
processor_id_type processor_id() const
Definition: parallel_object.h:106
libMesh::EquationSystems::const_system_iterator
std::map< std::string, System * >::const_iterator const_system_iterator
Typedef for constant system iterators.
Definition: equation_systems.h:581
libMesh::System::update_global_solution
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:642
libMesh::EquationSystems::update
void update()
Updates local values for all the systems.
Definition: equation_systems.C:334
libMesh::System::current_local_solution
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:1551
libMesh::ParallelObject::_communicator
const Parallel::Communicator & _communicator
Definition: parallel_object.h:112
libMesh::EquationSystems::init
virtual void init()
Initialize all the systems.
Definition: equation_systems.C:96
libMesh::Variable
This class defines the notion of a variable in the system.
Definition: variable.h:49
libMesh::EquationSystems::n_vars
unsigned int n_vars() const
Definition: equation_systems.C:1331
libMesh::ParameterVector
Data structure for specifying which Parameters should be independent variables in a parameter sensiti...
Definition: parameter_vector.h:44
libMesh::EquationSystems::_enable_default_ghosting
bool _enable_default_ghosting
Flag for whether to enable default ghosting on newly added Systems.
Definition: equation_systems.h:593
libMesh::Node
A Node is like a Point, but with more information.
Definition: node.h:52
libMesh::MeshRefinement::clean_refinement_flags
void clean_refinement_flags()
Sets the refinement flag to Elem::DO_NOTHING for each element in the mesh.
Definition: mesh_refinement_flagging.C:682
libMesh::EquationSystems::system_iterator
std::map< std::string, System * >::iterator system_iterator
Typedef for system iterators.
Definition: equation_systems.h:576
libMesh::CONSTANT
Definition: enum_order.h:41
libMesh::EquationSystems::enable_default_ghosting
virtual void enable_default_ghosting(bool enable)
Enable or disable default ghosting functors on the Mesh and on all Systems.
Definition: equation_systems.C:312
libMesh::EquationSystems::n_systems
unsigned int n_systems() const
Definition: equation_systems.h:652
libMesh::EquationSystems::_mesh
MeshBase & _mesh
The mesh data structure.
Definition: equation_systems.h:566
libMesh::EquationSystems::n_dofs
std::size_t n_dofs() const
Definition: equation_systems.C:1346
libMesh::EquationSystems::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
Definition: equation_systems.C:1314
libMesh::numeric_index_type
dof_id_type numeric_index_type
Definition: id_types.h:99
libMesh::MONOMIAL
Definition: enum_fe_family.h:39
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::Parameters::set
T & set(const std::string &)
Definition: parameters.h:460
libMesh::MeshBase::allgather
virtual void allgather()
Gathers all elements and nodes of the mesh onto every processor.
Definition: mesh_base.h:188
libMesh::System::variable_type
const FEType & variable_type(const unsigned int i) const
Definition: system.h:2233
distance
Real distance(const Point &p)
Definition: subdomains_ex3.C:50
libMesh::EquationSystems::reinit
virtual void reinit()
Handle any mesh changes and reinitialize all the systems on the updated mesh.
Definition: equation_systems.C:121
libMesh::EquationSystems::compare
virtual bool compare(const EquationSystems &other_es, const Real threshold, const bool verbose) const
Definition: equation_systems.C:1217
libMesh::DofMap::remove_default_ghosting
void remove_default_ghosting()
Remove any default ghosting functor(s).
Definition: dof_map.C:1821
libMesh::EquationSystems::_systems
std::map< std::string, System * > _systems
Data structure holding the systems.
Definition: equation_systems.h:571
libMesh::System::solution
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1539
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::EquationSystems::get_info
virtual std::string get_info() const
Definition: equation_systems.C:1278
libMesh::NumericVector::set
virtual void set(const numeric_index_type i, const T value)=0
Sets v(i) = value.
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::FEInterface::field_type
static FEFieldType field_type(const FEType &fe_type)
Definition: fe_interface.C:1683
libMesh::EquationSystems::build_discontinuous_solution_vector
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) const
Fill the input vector soln with solution values.
Definition: equation_systems.C:1049
libMesh::TYPE_VECTOR
Definition: enum_fe_family.h:94
libMesh::EquationSystems::clear
virtual void clear()
Restores the data structure to a pristine state.
Definition: equation_systems.C:75
libMesh::EquationSystems::_refine_in_reinit
bool _refine_in_reinit
Flag for whether to call coarsen/refine in reinit().
Definition: equation_systems.h:587
libMesh::EquationSystems::adjoint_solve
virtual void adjoint_solve(const QoISet &qoi_indices=QoISet())
Call adjoint_solve on all the individual equation systems.
Definition: equation_systems.C:466
libMesh::EquationSystems::~EquationSystems
virtual ~EquationSystems()
Destructor.
Definition: equation_systems.C:68
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::MeshRefinement::coarsen_elements
bool coarsen_elements()
Only coarsens the user-requested elements.
Definition: mesh_refinement.C:608
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::Variable::active_on_subdomain
bool active_on_subdomain(subdomain_id_type sid) const
Definition: variable.h:136
libMesh::ParallelObject
An object whose state is distributed along a set of processors.
Definition: parallel_object.h:55
libMesh::EquationSystems::build_parallel_elemental_solution_vector
std::unique_ptr< NumericVector< Number > > build_parallel_elemental_solution_vector(std::vector< std::string > &names) const
Builds a parallel vector of CONSTANT MONOMIAL solution values corresponding to the entries in the inp...
Definition: equation_systems.C:957
libMesh::EquationSystems::build_parallel_solution_vector
std::unique_ptr< NumericVector< Number > > build_parallel_solution_vector(const std::set< std::string > *system_names=nullptr) const
A version of build_solution_vector which is appropriate for "parallel" output formats like Nemesis.
Definition: equation_systems.C:590
libMesh::EquationSystems::delete_system
void delete_system(const std::string &name)
Remove the system named name from the systems array.
Definition: equation_systems.C:431
libMesh::out
OStreamProxy out
libMesh::MeshBase::local_nodes_end
virtual node_iterator local_nodes_end()=0
libMesh::EquationSystems::get_solution
void get_solution(std::vector< Number > &soln, std::vector< std::string > &names) const
Retrieve the solution data for CONSTANT MONOMIALs.
Definition: equation_systems.C:873
libMesh::System::update
virtual void update()
Update the local values to reflect the solution on neighboring processors.
Definition: system.C:408
libMesh::EquationSystems::_add_system_to_nodes_and_elems
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...
Definition: equation_systems.C:1376
libMesh::MeshBase::contract
virtual bool contract()=0
Delete subactive (i.e.
libMesh::Quality::name
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
libMesh::EquationSystems::get_vars_active_subdomains
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.
Definition: equation_systems.C:841
libMesh::NumericVector::first_local_index
virtual numeric_index_type first_local_index() const =0
libMesh::EquationSystems::parameters
Parameters parameters
Data structure holding arbitrary parameters.
Definition: equation_systems.h:557