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