Line data Source code
1 : // The libMesh Finite Element Library.
2 : // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 :
4 : // This library is free software; you can redistribute it and/or
5 : // modify it under the terms of the GNU Lesser General Public
6 : // License as published by the Free Software Foundation; either
7 : // version 2.1 of the License, or (at your option) any later version.
8 :
9 : // This library is distributed in the hope that it will be useful,
10 : // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 : // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 : // Lesser General Public License for more details.
13 :
14 : // You should have received a copy of the GNU Lesser General Public
15 : // License along with this library; if not, write to the Free Software
16 : // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 :
18 :
19 :
20 : #ifndef LIBMESH_DOF_MAP_H
21 : #define LIBMESH_DOF_MAP_H
22 :
23 : // Local Includes
24 : #include "libmesh/libmesh_common.h"
25 : #include "libmesh/reference_counted_object.h"
26 : #include "libmesh/libmesh.h" // libMesh::invalid_uint
27 : #include "libmesh/variable.h"
28 : #include "libmesh/threads.h"
29 : #include "libmesh/threads_allocators.h"
30 : #include "libmesh/elem_range.h"
31 : #include "libmesh/ghosting_functor.h"
32 : #include "libmesh/sparsity_pattern.h"
33 : #include "libmesh/parallel_object.h"
34 : #include "libmesh/point.h"
35 : #include "libmesh/utility.h"
36 : #include "libmesh/elem.h"
37 : #include "libmesh/fe_interface.h"
38 : #include "libmesh/libmesh_logging.h"
39 : #include "libmesh/enum_elem_type.h"
40 : #include "libmesh/mesh_subdivision_support.h"
41 : #include "libmesh/dof_map_base.h"
42 :
43 : // TIMPI includes
44 : #include "timpi/parallel_implementation.h"
45 : #include "timpi/parallel_sync.h"
46 :
47 : // C++ Includes
48 : #include <algorithm>
49 : #include <cstddef>
50 : #include <iterator>
51 : #include <map>
52 : #include <string>
53 : #include <vector>
54 : #include <memory>
55 :
56 : namespace libMesh
57 : {
58 :
59 : // Forward Declarations
60 : class CouplingMatrix;
61 : class DefaultCoupling;
62 : class DirichletBoundary;
63 : class DirichletBoundaries;
64 : class DofMap;
65 : class DofObject;
66 : class FEType;
67 : class MeshBase;
68 : class PeriodicBoundaryBase;
69 : class PeriodicBoundaries;
70 : class System;
71 : class NonlinearImplicitSystem;
72 : class StaticCondensationDofMap;
73 : template <typename T> class DenseVectorBase;
74 : template <typename T> class DenseVector;
75 : template <typename T> class DenseMatrix;
76 : template <typename T> class SparseMatrix;
77 : template <typename T> class NumericVector;
78 : enum Order : int;
79 :
80 :
81 :
82 : // ------------------------------------------------------------
83 : // Do we need constraints for anything?
84 :
85 : #if defined(LIBMESH_ENABLE_AMR) || \
86 : defined(LIBMESH_ENABLE_PERIODIC) || \
87 : defined(LIBMESH_ENABLE_DIRICHLET)
88 : # define LIBMESH_ENABLE_CONSTRAINTS 1
89 : #endif
90 :
91 : // ------------------------------------------------------------
92 : // AMR constraint matrix types
93 :
94 : #ifdef LIBMESH_ENABLE_CONSTRAINTS
95 : /**
96 : * A row of the Dof constraint matrix.
97 : */
98 : typedef std::map<dof_id_type, Real,
99 : std::less<dof_id_type>,
100 : Threads::scalable_allocator<std::pair<const dof_id_type, Real>>> DofConstraintRow;
101 :
102 : /**
103 : * The constraint matrix storage format.
104 : * We're using a class instead of a typedef to allow forward
105 : * declarations and future flexibility. Don't delete this from
106 : * a pointer-to-std-map; the destructor isn't virtual!
107 : */
108 232166 : class DofConstraints : public std::map<dof_id_type,
109 : DofConstraintRow,
110 : std::less<dof_id_type>,
111 : Threads::scalable_allocator<std::pair<const dof_id_type, DofConstraintRow>>>
112 : {
113 : };
114 :
115 : /**
116 : * Storage for DofConstraint right hand sides for a particular
117 : * problem. Each dof id with a non-zero constraint offset
118 : * stores it in such a structure.
119 : */
120 232030 : class DofConstraintValueMap :
121 : public std::map<dof_id_type, Number,
122 : std::less<dof_id_type>,
123 : Threads::scalable_allocator<std::pair<const dof_id_type, Number>>>
124 : {
125 : };
126 :
127 : /**
128 : * Storage for DofConstraint right hand sides for all adjoint
129 : * problems.
130 : */
131 232030 : class AdjointDofConstraintValues :
132 : public std::map<unsigned int, DofConstraintValueMap,
133 : std::less<unsigned int>,
134 : Threads::scalable_allocator
135 : <std::pair<const unsigned int, DofConstraintValueMap>>>
136 : {
137 : };
138 :
139 : #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
140 : /**
141 : * A row of the Node constraint mapping. Currently this just
142 : * stores the topology of the constrained Nodes, but for forward
143 : * compatibility we also include coefficients, so we could add
144 : * Lagrange-positioned-node constraints later.
145 : */
146 : typedef std::map<const Node *, Real,
147 : std::less<const Node *>,
148 : Threads::scalable_allocator<std::pair<const Node * const, Real>>> NodeConstraintRow;
149 :
150 : /**
151 : * The Node constraint storage format.
152 : * We're using a class instead of a typedef to allow forward
153 : * declarations and future flexibility. Don't delete this from
154 : * a pointer-to-std-map; the destructor isn't virtual!
155 : */
156 6816 : class NodeConstraints : public std::map<const Node *,
157 : std::pair<NodeConstraintRow,Point>,
158 : std::less<const Node *>,
159 : Threads::scalable_allocator<std::pair<const Node * const, std::pair<NodeConstraintRow,Point>>>>
160 : {
161 : };
162 : #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
163 :
164 : #endif // LIBMESH_ENABLE_CONSTRAINTS
165 :
166 :
167 :
168 : /**
169 : * This class handles the numbering of degrees of freedom on a mesh.
170 : * For systems of equations the class supports a fixed number of variables.
171 : * The degrees of freedom are numbered such that sequential, contiguous blocks
172 : * belong to distinct processors. This is so that the resulting data
173 : * structures will work well with parallel linear algebra packages.
174 : *
175 : * \author Benjamin S. Kirk
176 : * \date 2002-2007
177 : * \brief Manages the degrees of freedom (DOFs) in a simulation.
178 : */
179 : class DofMap : public DofMapBase,
180 : public ReferenceCountedObject<DofMap>
181 : {
182 : public:
183 :
184 : /**
185 : * Constructor. Requires the number of the system for which we
186 : * will be numbering degrees of freedom & the parent object
187 : * we are contained in, which defines our communication space.
188 : */
189 : explicit
190 : DofMap(const unsigned int sys_number,
191 : MeshBase & mesh);
192 :
193 : /**
194 : * Destructor.
195 : */
196 : ~DofMap();
197 :
198 : /**
199 : * Backwards compatibility for prior AugmentSparsityPattern users.
200 : */
201 : class AugmentSparsityPattern : public SparsityPattern::AugmentSparsityPattern
202 : {};
203 :
204 : /**
205 : * Abstract base class to be used to add user-defined parallel
206 : * degree of freedom couplings.
207 : */
208 : class AugmentSendList
209 : {
210 : public:
211 : virtual ~AugmentSendList () = default;
212 :
213 : /**
214 : * User-defined function to augment the send list.
215 : */
216 : virtual void augment_send_list (std::vector<dof_id_type> & send_list) = 0;
217 : };
218 :
219 : /**
220 : * Additional matrices may be attached to this \p DofMap.
221 : * They are initialized to the same sparsity structure as
222 : * the major matrix.
223 : */
224 : void attach_matrix (SparseMatrix<Number> & matrix);
225 :
226 : /**
227 : * Additional matrices may be be temporarily initialized by this \p
228 : * DofMap.
229 : * They are initialized to the same sparsity structure as
230 : * the major matrix.
231 : */
232 : void update_sparsity_pattern(SparseMatrix<Number> & matrix) const;
233 :
234 : /**
235 : * Matrices should not be attached more than once. We can test for
236 : * an already-attached matrix if necessary using \p is_attached
237 : */
238 : bool is_attached (SparseMatrix<Number> & matrix);
239 :
240 : /**
241 : * Distribute dofs on the current mesh. Also builds the send list for
242 : * processor \p proc_id, which defaults to 0 for ease of use in serial
243 : * applications.
244 : * \returns The total number of DOFs for the System, summed across all procs.
245 : */
246 : std::size_t distribute_dofs (MeshBase &);
247 :
248 : /**
249 : * Computes the sparsity pattern for the matrices corresponding to
250 : * \p proc_id and sends that data to Linear Algebra packages for
251 : * preallocation of sparse matrices.
252 : */
253 : void compute_sparsity (const MeshBase &);
254 :
255 : /**
256 : * Returns true iff a sparsity pattern has already been computed.
257 : */
258 : bool computed_sparsity_already () const;
259 :
260 : /**
261 : * Sets the current policy for constructing sparsity patterns: if
262 : * \p use_constraints is true (for robustness), we explicitly
263 : * account for sparsity entries created by constraint matrix pre-
264 : * and post- application. If \p use_constraints is false (for
265 : * speed), we calculate only the sparsity pattern of an
266 : * unconstrained matrix. This is false by default, because in
267 : * nearly all applications our constraints do not increase the
268 : * number of non-zeros required in a sparse matrix.
269 : */
270 : void set_constrained_sparsity_construction(bool use_constraints);
271 :
272 : /**
273 : * Sets need_full_sparsity_pattern to true regardless of the requirements by matrices
274 : */
275 : void full_sparsity_pattern_needed();
276 :
277 : /**
278 : * Returns true iff the current policy when constructing sparsity
279 : * patterns is to explicitly account for sparsity entries created by
280 : * constraint matrix pre- and post- application.
281 : */
282 : bool constrained_sparsity_construction();
283 :
284 : /**
285 : * Clears the sparsity pattern
286 : */
287 : void clear_sparsity();
288 :
289 : /**
290 : * Remove any default ghosting functor(s). User-added ghosting
291 : * functors will be unaffected.
292 : *
293 : * Unless user-added equivalent ghosting functors exist, removing
294 : * the default coupling functor is only safe for explicit solves,
295 : * and removing the default algebraic ghosting functor is only safe
296 : * for codes where no evaluations on neighbor cells (e.g. no jump
297 : * error estimators) are done.
298 : *
299 : * Defaults can be restored manually via add_default_ghosting(), or
300 : * automatically if clear() returns the DofMap to a default state.
301 : */
302 : void remove_default_ghosting();
303 :
304 : /**
305 : * Add the default functor(s) for coupling and algebraic ghosting.
306 : * User-added ghosting functors will be unaffected.
307 : */
308 : void add_default_ghosting();
309 :
310 : /**
311 : * Iterator type for coupling and algebraic ghosting functor ranges.
312 : * This has changed in the past and may change again; code should
313 : * use auto or the type here.
314 : */
315 : typedef std::vector<GhostingFunctor *>::const_iterator GhostingFunctorIterator;
316 :
317 : /**
318 : * Adds a functor which can specify coupling requirements for
319 : * creation of sparse matrices.
320 : * Degree of freedom pairs which match the elements and variables
321 : * returned by these functors will be added to the sparsity pattern,
322 : * and the degrees of freedom which live on other processors will be
323 : * added to the send_list for use on ghosted vectors, and the
324 : * elements which live on other processors will be ghosted on a
325 : * distributed mesh.
326 : *
327 : * GhostingFunctor memory must be managed by the code which calls
328 : * this function; the GhostingFunctor lifetime is expected to extend
329 : * until either the functor is removed or the DofMap is destructed.
330 : *
331 : * When \p to_mesh is true, the \p coupling_functor is also added to
332 : * our associated mesh, to ensure that coupled elements do not get
333 : * lost during mesh distribution. (if coupled elements were
334 : * *already* lost there's no getting them back after the fact,
335 : * sorry)
336 : *
337 : * If \p to_mesh is false, no change to mesh ghosting is made;
338 : * the Mesh must already have ghosting functor(s) specifying a
339 : * superset of \p coupling_functor or this is a horrible bug.
340 : */
341 : void add_coupling_functor(GhostingFunctor & coupling_functor,
342 : bool to_mesh = true);
343 :
344 : /**
345 : * Adds a functor which can specify coupling requirements for
346 : * creation of sparse matrices.
347 : *
348 : * GhostingFunctor memory when using this method is managed by the
349 : * shared_ptr mechanism.
350 : */
351 : void add_coupling_functor(std::shared_ptr<GhostingFunctor> coupling_functor,
352 : bool to_mesh = true)
353 : { _shared_functors[coupling_functor.get()] = coupling_functor;
354 : this->add_coupling_functor(*coupling_functor, to_mesh); }
355 :
356 : /**
357 : * Removes a functor which was previously added to the set of
358 : * coupling functors, from both this DofMap and from the underlying
359 : * mesh.
360 : */
361 : void remove_coupling_functor(GhostingFunctor & coupling_functor);
362 :
363 : /**
364 : * Beginning of range of coupling functors
365 : */
366 314003 : GhostingFunctorIterator coupling_functors_begin() const
367 637787 : { return _coupling_functors.begin(); }
368 :
369 : /**
370 : * End of range of coupling functors
371 : */
372 314003 : GhostingFunctorIterator coupling_functors_end() const
373 637787 : { return _coupling_functors.end(); }
374 :
375 : /**
376 : * Default coupling functor
377 : */
378 0 : DefaultCoupling & default_coupling() { return *_default_coupling; }
379 :
380 : /**
381 : * Adds a functor which can specify algebraic ghosting requirements
382 : * for use with distributed vectors. Degrees of freedom on other
383 : * processors which match the elements and variables returned by
384 : * these functors will be added to the send_list, and the elements
385 : * on other processors will be ghosted on a distributed mesh, so
386 : * that the elements can always be found and the solutions on them
387 : * will always be evaluable.
388 : *
389 : * GhostingFunctor memory must be managed by the code which calls
390 : * this function; the GhostingFunctor lifetime is expected to extend
391 : * until either the functor is removed or the DofMap is destructed.
392 : *
393 : * When \p to_mesh is true, the \p coupling_functor is also added to
394 : * our associated mesh, to ensure that evaluable elements do not get
395 : * lost during mesh distribution. (if evaluable elements were
396 : * *already* lost there's no getting them back after the fact,
397 : * sorry)
398 : *
399 : * If \p to_mesh is false, no change to mesh ghosting is made;
400 : * the Mesh must already have ghosting functor(s) specifying a
401 : * superset of \p evaluable_functor or this is a horrible bug.
402 : */
403 : void add_algebraic_ghosting_functor(GhostingFunctor & evaluable_functor,
404 : bool to_mesh = true);
405 :
406 : /**
407 : * Adds a functor which can specify algebraic ghosting requirements
408 : * for use with distributed vectors.
409 : *
410 : * GhostingFunctor memory when using this method is managed by the
411 : * shared_ptr mechanism.
412 : */
413 : void add_algebraic_ghosting_functor(std::shared_ptr<GhostingFunctor> evaluable_functor,
414 : bool to_mesh = true)
415 : { _shared_functors[evaluable_functor.get()] = evaluable_functor;
416 : this->add_algebraic_ghosting_functor(*evaluable_functor, to_mesh); }
417 :
418 : /**
419 : * Removes a functor which was previously added to the set of
420 : * algebraic ghosting functors, from both this DofMap and from the
421 : * underlying mesh.
422 : */
423 : void remove_algebraic_ghosting_functor(GhostingFunctor & evaluable_functor);
424 :
425 : /**
426 : * Beginning of range of algebraic ghosting functors
427 : */
428 1568 : GhostingFunctorIterator algebraic_ghosting_functors_begin() const
429 9222 : { return _algebraic_ghosting_functors.begin(); }
430 :
431 : /**
432 : * End of range of algebraic ghosting functors
433 : */
434 1568 : GhostingFunctorIterator algebraic_ghosting_functors_end() const
435 9222 : { return _algebraic_ghosting_functors.end(); }
436 :
437 : /**
438 : * Default algebraic ghosting functor
439 : */
440 0 : DefaultCoupling & default_algebraic_ghosting() { return *_default_evaluating; }
441 :
442 : /**
443 : * Attach an object to use to populate the
444 : * sparsity pattern with extra entries.
445 : *
446 : * Care must be taken that when adding entries they are sorted into the Rows
447 : *
448 : * Further, you _must_ modify n_nz and n_oz properly!
449 : *
450 : * This is an advanced function... use at your own peril!
451 : */
452 : void attach_extra_sparsity_object (SparsityPattern::AugmentSparsityPattern & asp)
453 : {
454 : _augment_sparsity_pattern = &asp;
455 : }
456 :
457 : /**
458 : * Attach a function pointer to use as a callback to populate the
459 : * sparsity pattern with extra entries.
460 : *
461 : * Care must be taken that when adding entries they are sorted into the Rows
462 : *
463 : * Further, you _must_ modify n_nz and n_oz properly!
464 : *
465 : * This is an advanced function... use at your own peril!
466 : */
467 : void attach_extra_sparsity_function(void (*func)(SparsityPattern::Graph & sparsity,
468 : std::vector<dof_id_type> & n_nz,
469 : std::vector<dof_id_type> & n_oz,
470 : void *),
471 : void * context = nullptr)
472 : { _extra_sparsity_function = func; _extra_sparsity_context = context; }
473 :
474 : /**
475 : * Attach an object to populate the send_list with extra entries.
476 : * This should only add to the send list, but no checking is done
477 : * to enforce this behavior.
478 : *
479 : * This is an advanced function... use at your own peril!
480 : */
481 : void attach_extra_send_list_object (DofMap::AugmentSendList & asl)
482 : {
483 : _augment_send_list = &asl;
484 : }
485 :
486 : /**
487 : * Attach a function pointer to use as a callback to populate the
488 : * send_list with extra entries.
489 : */
490 : void attach_extra_send_list_function(void (*func)(std::vector<dof_id_type> &, void *),
491 : void * context = nullptr)
492 : { _extra_send_list_function = func; _extra_send_list_context = context; }
493 :
494 : /**
495 : * Takes the \p _send_list vector (which may have duplicate entries)
496 : * and sorts it. The duplicate entries are then removed, resulting in
497 : * a sorted \p _send_list with unique entries. Also calls any user-provided
498 : * methods for adding to the send list.
499 : */
500 : void prepare_send_list ();
501 :
502 : /**
503 : * Clears the \p _send_list vector. This should be done in order to completely
504 : * rebuild the send_list from scratch rather than merely adding to the existing
505 : * send_list.
506 : */
507 14496 : void clear_send_list ()
508 : {
509 14496 : _send_list.clear();
510 14496 : }
511 :
512 : /**
513 : * Clears the \p _send_list vector and then rebuilds it. This may be needed
514 : * in special situations, for example when an algebraic coupling functor cannot
515 : * be added to the \p DofMap until after it is completely setup. Then this method
516 : * can be used to rebuild the send_list once the algebraic coupling functor is
517 : * added. Note that while this will recommunicate constraints with the updated
518 : * send_list, this does assume no new constraints have been added since the previous
519 : * reinit_constraints call.
520 : */
521 : void reinit_send_list (MeshBase & mesh);
522 :
523 :
524 : /**
525 : * \returns A constant reference to the \p _send_list for this processor.
526 : *
527 : * The \p _send_list contains the global indices of all the
528 : * variables in the global solution vector that influence the
529 : * current processor. This information can be used for gathers at
530 : * each solution step to retrieve solution values needed for
531 : * computation.
532 : */
533 2790612 : const std::vector<dof_id_type> & get_send_list() const { return _send_list; }
534 :
535 : /**
536 : * \returns A constant reference to the \p _n_nz list for this processor.
537 : *
538 : * The vector contains the bandwidth of the on-processor coupling for each
539 : * row of the global matrix that the current processor owns. This
540 : * information can be used to preallocate space for a parallel sparse matrix.
541 : */
542 : const std::vector<dof_id_type> & get_n_nz() const
543 : {
544 : libmesh_assert(_sp);
545 : return _sp->get_n_nz();
546 : }
547 :
548 : /**
549 : * \returns A constant reference to the \p _n_oz list for this processor.
550 : *
551 : * The vector contains the bandwidth of the off-processor coupling for each
552 : * row of the global matrix that the current processor owns. This
553 : * information can be used to preallocate space for a parallel sparse matrix.
554 : */
555 : const std::vector<dof_id_type> & get_n_oz() const
556 : {
557 : libmesh_assert(_sp);
558 : return _sp->get_n_oz();
559 : }
560 :
561 :
562 : /**
563 : * \returns A constant pointer to the sparsity pattern stored here,
564 : * once that has been computed. Returns null if no sparsity pattern
565 : * has yet been computed.
566 : *
567 : * If need_full_sparsity_pattern is false, the "sparsity pattern"
568 : * may only own n_nz and n_oz lists.
569 : */
570 664 : const SparsityPattern::Build * get_sparsity_pattern() const
571 : {
572 664 : return _sp.get();
573 : }
574 :
575 : /**
576 : * \returns The number of variables in the system
577 : */
578 : unsigned int n_vars() const;
579 :
580 : /**
581 : * \returns The name of variable \p i.
582 : */
583 : const std::string & variable_name(const unsigned int i) const;
584 :
585 : /**
586 : * \returns The total number of scalar components in the system's
587 : * variables. This will equal \p n_vars() in the case of all
588 : * scalar-valued variables. If vector variables are involved, we
589 : * will need to leverage the \p mesh
590 : */
591 : unsigned int n_components(const MeshBase & mesh) const;
592 :
593 : /**
594 : * \returns \p true when \p VariableGroup structures should be
595 : * automatically identified, \p false otherwise.
596 : */
597 : bool identify_variable_groups () const;
598 :
599 : /**
600 : * Toggle automatic \p VariableGroup identification.
601 : */
602 : void identify_variable_groups (const bool);
603 :
604 : /**
605 : * \returns An index, starting from 0 for the first component of the
606 : * first variable, and incrementing for each component of each
607 : * (potentially vector-valued) variable in the system in order.
608 : * For systems with only scalar-valued variables, this will be the
609 : * same as \p var_num
610 : *
611 : * Irony: currently our only non-scalar-valued variable type is
612 : * SCALAR.
613 : */
614 : unsigned int variable_scalar_number (unsigned int var_num,
615 : unsigned int component) const;
616 :
617 : /**
618 : * \returns The finite element type for variable number \p i.
619 : */
620 : const FEType & variable_type (const unsigned int i) const;
621 :
622 : /**
623 : * \returns The finite element type for variable \p var.
624 : */
625 : const FEType & variable_type (std::string_view var) const;
626 :
627 : /**
628 : * \returns The variable number associated with
629 : * the user-specified variable named \p var.
630 : */
631 : unsigned int variable_number (std::string_view var) const;
632 :
633 : /**
634 : * \returns \p true if a variable named \p var exists in this System
635 : */
636 : bool has_variable(std::string_view var) const;
637 :
638 : /**
639 : * Fills \p all_variable_numbers with all the variable numbers for the
640 : * variables that have been added to this system.
641 : */
642 : void get_all_variable_numbers(std::vector<unsigned int> & all_variable_numbers) const;
643 :
644 : /**
645 : * Adds the variable \p var to the list of variables
646 : * for this system. If \p active_subdomains is either \p nullptr
647 : * (the default) or points to an empty set, then it will be assumed that
648 : * \p var has no subdomain restrictions
649 : *
650 : * \returns The index number for the new variable.
651 : */
652 : unsigned int add_variable (System & sys,
653 : std::string_view var,
654 : const FEType & type,
655 : const std::set<subdomain_id_type> * const active_subdomains = nullptr);
656 :
657 : /**
658 : * Adds the variables \p vars to the list of variables
659 : * for this system. If \p active_subdomains is either \p nullptr
660 : * (the default) or points to an empty set, then it will be assumed that
661 : * the \p vars have no subdomain restrictions
662 : *
663 : * \returns The index number for the last of the new variables.
664 : */
665 : unsigned int add_variables (System & sys,
666 : const std::vector<std::string> & vars,
667 : const FEType & type,
668 : const std::set<subdomain_id_type> * const active_subdomains = nullptr);
669 :
670 : /**
671 : * Adds variables \p vars to the list of variables
672 : * for this system. If \p active_subdomains is either \p nullptr
673 : * (the default) or points to an empty set, then it will be assumed that
674 : * the \p vars have no subdomain restrictions. This API will end up
675 : * calling \p this->add_variables(). However, we will additionally store data
676 : * that can be leveraged by the \p DofMap to build degrees of freedom
677 : * containers corresponding to all the variables in this variable array
678 : *
679 : * An 'array variable' is simply a sequence
680 : * of contiguous variable numbers defined by pair where the first member of the pair
681 : * is the first number in the variable sequence and the second member of the pair is
682 : * the number of the last variable in the sequence plus one. Array variables may be
683 : * used in tandem with variable grouping by downstream code to build optimized physics
684 : * kernels since each variable in the array will have the same shape functions.
685 : *
686 : * \returns The index number for the last of the new variables.
687 : */
688 : unsigned int add_variable_array (System & sys,
689 : const std::vector<std::string> & vars,
690 : const FEType & type,
691 : const std::set<subdomain_id_type> * const active_subdomains = nullptr);
692 :
693 : /**
694 : * Specify whether or not we perform an extra (opt-mode enabled) check
695 : * for constraint loops. If a constraint loop is present then
696 : * the system constraints are not valid, so if \p error_on_constraint_loop
697 : * is true we will throw an error in this case.
698 : *
699 : * \note We previously referred to these types of constraints as
700 : * "cyclic" but that has now been deprecated, and these will now
701 : * instead be referred to as "constraint loops" in libMesh.
702 : */
703 : void set_error_on_cyclic_constraint(bool error_on_cyclic_constraint);
704 : void set_error_on_constraint_loop(bool error_on_constraint_loop);
705 :
706 : /**
707 : * \returns The \p VariableGroup description object for group \p g.
708 : */
709 : const VariableGroup & variable_group (const unsigned int c) const;
710 :
711 : const Variable & variable (const unsigned int c) const override;
712 :
713 : /**
714 : * \returns The approximation order for variable \p c.
715 : */
716 : Order variable_order (const unsigned int c) const;
717 :
718 : /**
719 : * \returns The approximation order for \p VariableGroup \p vg.
720 : */
721 : Order variable_group_order (const unsigned int vg) const;
722 :
723 : /**
724 : * \returns The finite element type for \p VariableGroup \p vg.
725 : */
726 : const FEType & variable_group_type (const unsigned int vg) const;
727 :
728 : /**
729 : * \returns The number of variables in the global solution vector. Defaults
730 : * to 1, should be 1 for a scalar equation, 3 for 2D incompressible Navier
731 : * Stokes (u,v,p), etc...
732 : */
733 11913620 : unsigned int n_variable_groups() const
734 23619439 : { return cast_int<unsigned int>(_variable_groups.size()); }
735 :
736 2732453 : unsigned int n_variables() const override
737 2783447 : { return cast_int<unsigned int>(_variables.size()); }
738 :
739 : /**
740 : * \returns The variable group number that the provided variable number belongs to
741 : */
742 : unsigned int var_group_from_var_number(unsigned int var_num) const;
743 :
744 : /**
745 : * \returns \p true if the variables are capable of being stored in a blocked
746 : * form. Presently, this means that there can only be one variable group,
747 : * and that the group has more than one variable.
748 : */
749 40776 : bool has_blocked_representation() const
750 : {
751 42116 : return ((this->n_variable_groups() == 1) && (this->n_variables() > 1));
752 : }
753 :
754 : /**
755 : * \returns The block size, if the variables are amenable to block storage.
756 : * Otherwise 1.
757 : * This routine was originally designed to enable a blocked storage, but
758 : * it turns out this information is still super useful for solvers even when
759 : * we do not use the blocked storage (e.g., MATMPIBAIJ in PETSc). For example (in PCHMG),
760 : * for a system of PDEs, to construct an efficient multilevel preconditioner, we coarsen
761 : * the matrix of one single PDE instead of the entire huge matrix. In order to
762 : * accomplish this, we need to know many PDEs we have. Another use case,
763 : * the fieldsplit preconditioner can be constructed in place with this info without
764 : * involving any user efforts.
765 : */
766 42116 : unsigned int block_size() const
767 : {
768 40892 : return (this->has_blocked_representation() ? this->n_variables() : 1);
769 : }
770 :
771 : using DofMapBase::n_dofs;
772 : /**
773 : * \returns The total number of degrees of freedom for a particular
774 : * variable \p vn.
775 : */
776 : dof_id_type n_dofs(const unsigned int vn) const
777 : {
778 : dof_id_type n = this->n_local_dofs(vn);
779 : this->comm().sum(n);
780 : return n;
781 : }
782 :
783 : /**
784 : * \returns The number of SCALAR dofs.
785 : */
786 256797 : dof_id_type n_SCALAR_dofs() const { return _n_SCALAR_dofs; }
787 :
788 : using DofMapBase::n_local_dofs;
789 : /**
790 : * \returns The number of degrees of freedom on this processor for a
791 : * particular variable \p vn. This is an O(N) operation on serialized or
792 : * O(N/Nproc) operation on distributed meshes.
793 : */
794 : dof_id_type n_local_dofs(const unsigned int vn) const
795 : {
796 : dof_id_type n;
797 : this->local_variable_indices(n, _mesh, vn);
798 : return n;
799 : }
800 :
801 : /**
802 : * \returns The number of degrees of freedom on each partition for a
803 : * particular variable \p vn.
804 : */
805 : std::vector<dof_id_type> n_dofs_per_processor(const unsigned int vn) const
806 : {
807 : std::vector<dof_id_type> n_local_dofs(this->n_processors(), 0);
808 : this->comm().allgather(this->n_local_dofs(vn), n_local_dofs);
809 : return n_local_dofs;
810 : }
811 :
812 : /**
813 : * \returns The processor id that owns the dof index \p dof
814 : */
815 94065 : processor_id_type dof_owner(const dof_id_type dof) const
816 : { std::vector<dof_id_type>::const_iterator ub =
817 94065 : std::upper_bound(_end_df.begin(), _end_df.end(), dof);
818 2602 : libmesh_assert (ub != _end_df.end());
819 94065 : return cast_int<processor_id_type>(ub - _end_df.begin());
820 : }
821 :
822 : void dof_indices (const Elem * const elem,
823 : std::vector<dof_id_type> & di) const;
824 :
825 : /**
826 : * Fills the vector \p di with the global degree of freedom indices
827 : * for the element. For one variable, and potentially for a
828 : * non-default element p refinement level
829 : */
830 : void dof_indices (const Elem * const elem,
831 : std::vector<dof_id_type> & di,
832 : const unsigned int vn,
833 : int p_level = -12345) const override;
834 :
835 : /**
836 : * Fills the vector \p di with the global degree of freedom indices
837 : * for the element. This will aggregate all the degrees of the freedom
838 : * from the variable array that \p vn is a member of, and potentially for a
839 : * non-default element p refinement level
840 : */
841 : void array_dof_indices(const Elem * const elem,
842 : std::vector<dof_id_type> & di,
843 : const unsigned int vn,
844 : int p_level = -12345) const;
845 :
846 : void array_dof_indices(const Node * const node,
847 : std::vector<dof_id_type> & di,
848 : const unsigned int vn) const;
849 :
850 : template <typename DofIndicesFunctor>
851 : void array_dof_indices(const DofIndicesFunctor & functor,
852 : std::vector<dof_id_type> & di,
853 : const unsigned int vn) const;
854 :
855 : /**
856 : * Retrieves degree of freedom indices for a given \p elem and then performs actions for these
857 : * indices defined by the user-provided functors \p scalar_dofs_functor and \p field_dofs_functor.
858 : * This API is useful when a user wants to do more than simply fill a degree of freedom container
859 : * @param elem The element to get degrees of freedom for
860 : * @param di A container for degrees of freedom. It is up to the provided functors how this gets
861 : * filled
862 : * @param vn The variable number to retrieve degrees of freedom for
863 : * @param scalar_dofs_functor The functor that acts on scalar degrees of freedom. This functor has
864 : * the interface:
865 : * void scalar_dofs_functor(const Elem & elem,
866 : * std::vector<dof_id_type> & di,
867 : * const std::vector<dof_id_type> & scalar_dof_indices)
868 : * where \p di is the degree of freedom container described above and
869 : * \p scalar_dof_indices are the scalar dof indices available to
870 : * \p elem
871 : * @param field_dofs_functor The functor that acts on "field" (e.g. non-scalar, non-global)
872 : * degrees of freedom. This functor has
873 : * the interface:
874 : * void field_dofs_functor(const Elem & elem,
875 : * const unsigned int node_num,
876 : * const unsigned int var_num,
877 : * std::vector<dof_id_type> & di,
878 : * const dof_id_type field_dof)
879 : * where \p field_dof represents a field degree of freedom to act on and
880 : * is associated with \p node_num and \p var_num. If the degree of
881 : * freedom is elemental than \p node_num will be \p invalid_uint. \p di
882 : * is again the degree of freedom container provided above
883 : */
884 : template <typename ScalarDofsFunctor, typename FieldDofsFunctor>
885 : void dof_indices(const Elem * const elem,
886 : std::vector<dof_id_type> & di,
887 : const unsigned int vn,
888 : ScalarDofsFunctor scalar_dofs_functor,
889 : FieldDofsFunctor field_dofs_functor,
890 : int p_level = -12345) const;
891 :
892 : /**
893 : * Fills the vector \p di with the global degree of freedom indices
894 : * for the \p node.
895 : */
896 : void dof_indices (const Node * const node,
897 : std::vector<dof_id_type> & di) const;
898 :
899 : /**
900 : * Fills the vector \p di with the global degree of freedom indices
901 : * for the \p node, for one variable \p vn.
902 : */
903 : void dof_indices (const Node * const node,
904 : std::vector<dof_id_type> & di,
905 : const unsigned int vn) const override;
906 :
907 : /**
908 : * Appends to the vector \p di the global degree of freedom indices
909 : * for \p elem.node_ref(n), for one variable \p vn. On hanging
910 : * nodes with both vertex and non-vertex DoFs, only those indices
911 : * which are directly supported on \p elem are included.
912 : */
913 : void dof_indices (const Elem & elem,
914 : unsigned int n,
915 : std::vector<dof_id_type> & di,
916 : const unsigned int vn) const;
917 :
918 : #ifdef LIBMESH_ENABLE_AMR
919 :
920 : /**
921 : * Appends to the vector \p di the old global degree of freedom
922 : * indices for \p elem.node_ref(n), for one variable \p vn. On
923 : * hanging nodes with both vertex and non-vertex DoFs, only those
924 : * indices which are directly supported on \p elem are included.
925 : */
926 : void old_dof_indices (const Elem & elem,
927 : unsigned int n,
928 : std::vector<dof_id_type> & di,
929 : const unsigned int vn) const;
930 :
931 : #endif // LIBMESH_ENABLE_AMR
932 :
933 : /**
934 : * Fills the vector \p di with the global degree of freedom indices
935 : * corresponding to the SCALAR variable vn. If old_dofs=true,
936 : * the old SCALAR dof indices are returned.
937 : *
938 : * \note We do not need to pass in an element since SCALARs are
939 : * global variables.
940 : */
941 : void SCALAR_dof_indices (std::vector<dof_id_type> & di,
942 : const unsigned int vn,
943 : const bool old_dofs=false) const;
944 :
945 : /**
946 : * \returns \p true if degree of freedom index \p dof_index
947 : * is either a local index or in the \p send_list.
948 : *
949 : * \note This is an O(logN) operation for a send_list of size N; we
950 : * don't cache enough information for O(1) right now.
951 : */
952 : bool semilocal_index (dof_id_type dof_index) const;
953 :
954 : /**
955 : * \returns \p true if all degree of freedom indices in \p
956 : * dof_indices are either local indices or in the \p send_list.
957 : *
958 : * \note This is an O(logN) operation for a send_list of size N; we
959 : * don't cache enough information for O(1) right now.
960 : */
961 : bool all_semilocal_indices (const std::vector<dof_id_type> & dof_indices) const;
962 :
963 : /**
964 : * \returns \p true if degree of freedom index \p dof_index
965 : * is a local index.
966 : */
967 73744774 : bool local_index (dof_id_type dof_index) const
968 79953030 : { return (dof_index >= this->first_dof()) && (dof_index < this->end_dof()); }
969 :
970 : /**
971 : * \returns \p true iff our solutions can be locally evaluated on
972 : * \p obj (which should be an Elem or a Node) for variable number \p
973 : * var_num (for all variables, if \p var_num is invalid_uint)
974 : */
975 : template <typename DofObjectSubclass>
976 : bool is_evaluable(const DofObjectSubclass & obj,
977 : unsigned int var_num = libMesh::invalid_uint) const;
978 :
979 : /**
980 : * Allow the implicit_neighbor_dofs flag to be set programmatically.
981 : * This overrides the --implicit_neighbor_dofs commandline option.
982 : * We can use this to set the implicit neighbor dofs option differently
983 : * for different systems, whereas the commandline option is the same
984 : * for all systems.
985 : */
986 : void set_implicit_neighbor_dofs(bool implicit_neighbor_dofs);
987 :
988 : /**
989 : * Set the _verify_dirichlet_bc_consistency flag.
990 : */
991 : void set_verify_dirichlet_bc_consistency(bool val);
992 :
993 : /**
994 : * Tells other library functions whether or not this problem
995 : * includes coupling between dofs in neighboring cells, as can
996 : * currently be specified on the command line or inferred from
997 : * the use of all discontinuous variables.
998 : */
999 : bool use_coupled_neighbor_dofs(const MeshBase & mesh) const;
1000 :
1001 : /**
1002 : * Builds the local element vector \p Ue from the global vector \p Ug,
1003 : * accounting for any constrained degrees of freedom. For an element
1004 : * without constrained degrees of freedom this is the trivial mapping
1005 : * \f$ Ue[i] = Ug[dof_indices[i]] \f$
1006 : *
1007 : * \note The user must ensure that the element vector \p Ue is
1008 : * properly sized when calling this method. This is because there
1009 : * is no \p resize() method in the \p DenseVectorBase<> class.
1010 : */
1011 : void extract_local_vector (const NumericVector<Number> & Ug,
1012 : const std::vector<dof_id_type> & dof_indices,
1013 : DenseVectorBase<Number> & Ue) const;
1014 :
1015 : /**
1016 : * If T == dof_id_type, counts, if T == std::vector<dof_id_type>, fills an
1017 : * array of, those dof indices which belong to the given variable number and
1018 : * live on the current processor.
1019 : */
1020 : template <typename T, std::enable_if_t<std::is_same_v<T, dof_id_type> ||
1021 : std::is_same_v<T, std::vector<dof_id_type>>, int> = 0>
1022 : void local_variable_indices(T & idx,
1023 : const MeshBase & mesh,
1024 : unsigned int var_num) const;
1025 :
1026 : /**
1027 : * If T == dof_id_type, counts, if T == std::vector<dof_id_type>, fills an
1028 : * array of, those dof indices which belong to the given variable number and
1029 : * live on the current processor.
1030 : */
1031 : template <typename T,
1032 : std::enable_if_t<std::is_same_v<T, dof_id_type> ||
1033 : std::is_same_v<T, std::vector<dof_id_type>>,
1034 : int> = 0>
1035 : void local_variable_indices(T & idx, unsigned int var_num) const
1036 : { this->local_variable_indices(idx, this->_mesh, var_num); }
1037 :
1038 : #ifdef LIBMESH_ENABLE_CONSTRAINTS
1039 :
1040 : //--------------------------------------------------------------------
1041 : // Constraint-specific methods
1042 : /**
1043 : * \returns The total number of constrained degrees of freedom
1044 : * in the problem.
1045 : */
1046 : dof_id_type n_constrained_dofs() const;
1047 :
1048 : /**
1049 : * \returns The number of constrained degrees of freedom
1050 : * on this processor.
1051 : */
1052 : dof_id_type n_local_constrained_dofs() const;
1053 :
1054 : #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
1055 : /**
1056 : * \returns The total number of constrained Nodes
1057 : * in the mesh.
1058 : */
1059 : dof_id_type n_constrained_nodes() const
1060 : { return cast_int<dof_id_type>(_node_constraints.size()); }
1061 : #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
1062 :
1063 : /**
1064 : * Rebuilds the raw degree of freedom and DofObject constraints,
1065 : * based on attached DirichletBoundary objects and on non-conforming
1066 : * interface in adapted meshes.
1067 : *
1068 : * A time is specified for use in building time-dependent Dirichlet
1069 : * constraints.
1070 : */
1071 : void create_dof_constraints (const MeshBase &, Real time=0);
1072 :
1073 : /**
1074 : * Gathers constraint equation dependencies from other processors
1075 : */
1076 : void allgather_recursive_constraints (MeshBase &);
1077 :
1078 : /**
1079 : * Sends constraint equations to constraining processors
1080 : */
1081 : void scatter_constraints (MeshBase &);
1082 :
1083 : /**
1084 : * Helper function for querying about constraint equations on other
1085 : * processors. If any id in \p requested_dof_ids is constrained on
1086 : * another processor, its constraint will be added on this processor
1087 : * as well. If \p look_for_constrainees is true, then constraints
1088 : * will also be returned if the id appears as a constraining value
1089 : * not just if it appears as a constrained value.
1090 : *
1091 : * This function operates recursively: if the constraint for a
1092 : * constrained dof is newly added locally, then any other dofs which
1093 : * constrain it are queried to see if they are in turn constrained,
1094 : * and so on.
1095 : */
1096 : void gather_constraints (MeshBase & mesh,
1097 : std::set<dof_id_type> & unexpanded_dofs,
1098 : bool look_for_constrainees);
1099 :
1100 : /**
1101 : * Postprocesses any constrained degrees of freedom
1102 : * to be constrained only in terms of unconstrained dofs, then adds
1103 : * unconstrained dofs to the send_list and prepares that for use.
1104 : * This should be run after both system (create_dof_constraints) and
1105 : * user constraints have all been added.
1106 : */
1107 : void process_constraints (MeshBase &);
1108 :
1109 : /**
1110 : * Throw an error if we detect any constraint loops, i.e.
1111 : * A -> B -> C -> A
1112 : * that is, "dof A is constrained in terms of dof B which is
1113 : * constrained in terms of dof C which is constrained in terms of
1114 : * dof A", since these are not supported by libMesh and give
1115 : * erroneous results if they are present.
1116 : *
1117 : * \note The original "cyclic constraint" terminology was
1118 : * unfortunate since the word cyclic is used by some software to
1119 : * indicate an actual type of rotational/angular constraint and not
1120 : * (as here) a cyclic graph. The former nomenclature will eventually
1121 : * be deprecated in favor of "constraint loop".
1122 : */
1123 : void check_for_cyclic_constraints();
1124 : void check_for_constraint_loops();
1125 :
1126 : /**
1127 : * Adds a copy of the user-defined row to the constraint matrix, using
1128 : * an inhomogeneous right-hand-side for the constraint equation.
1129 : */
1130 : void add_constraint_row (const dof_id_type dof_number,
1131 : const DofConstraintRow & constraint_row,
1132 : const Number constraint_rhs,
1133 : const bool forbid_constraint_overwrite);
1134 :
1135 : /**
1136 : * Adds a copy of the user-defined row to the constraint matrix,
1137 : * using an inhomogeneous right-hand-side for the adjoint constraint
1138 : * equation.
1139 : *
1140 : * \p forbid_constraint_overwrite here only tests for overwriting
1141 : * the rhs. This method should only be used when an equivalent
1142 : * constraint (with a potentially different rhs) already exists for
1143 : * the primal problem.
1144 : */
1145 : void add_adjoint_constraint_row (const unsigned int qoi_index,
1146 : const dof_id_type dof_number,
1147 : const DofConstraintRow & constraint_row,
1148 : const Number constraint_rhs,
1149 : const bool forbid_constraint_overwrite);
1150 :
1151 : /**
1152 : * Adds a copy of the user-defined row to the constraint matrix, using
1153 : * a homogeneous right-hand-side for the constraint equation.
1154 : * By default, produces an error if the DOF was already constrained.
1155 : */
1156 19370 : void add_constraint_row (const dof_id_type dof_number,
1157 : const DofConstraintRow & constraint_row,
1158 : const bool forbid_constraint_overwrite = true)
1159 232736 : { add_constraint_row(dof_number, constraint_row, 0., forbid_constraint_overwrite); }
1160 :
1161 : /**
1162 : * \returns An iterator pointing to the first DoF constraint row.
1163 : */
1164 : DofConstraints::const_iterator constraint_rows_begin() const
1165 : { return _dof_constraints.begin(); }
1166 :
1167 : /**
1168 : * \returns An iterator pointing just past the last DoF constraint row.
1169 : */
1170 : DofConstraints::const_iterator constraint_rows_end() const
1171 : { return _dof_constraints.end(); }
1172 :
1173 : /**
1174 : * Provide a const accessor to the DofConstraints map. This allows the user
1175 : * to quickly search the data structure rather than just iterating over it.
1176 : */
1177 2812 : const DofConstraints & get_dof_constraints() const { return _dof_constraints; }
1178 :
1179 : void stash_dof_constraints()
1180 : {
1181 : libmesh_assert(_stashed_dof_constraints.empty());
1182 : _dof_constraints.swap(_stashed_dof_constraints);
1183 : }
1184 :
1185 : void unstash_dof_constraints()
1186 : {
1187 : libmesh_assert(_dof_constraints.empty());
1188 : _dof_constraints.swap(_stashed_dof_constraints);
1189 : }
1190 :
1191 : /**
1192 : * Similar to the stash/unstash_dof_constraints() API, but swaps
1193 : * _dof_constraints and _stashed_dof_constraints without asserting
1194 : * that the source or destination is empty first.
1195 : *
1196 : * \note There is an implicit assumption that swapping between sets
1197 : * of Constraints does not change the sparsity pattern or expand the
1198 : * send_list, since the only thing changed is the DofConstraints
1199 : * themselves. This is intended to work for swapping between
1200 : * DofConstraints A and B, where A is used to define the send_list,
1201 : * and B is a subset of A.
1202 : */
1203 : void swap_dof_constraints()
1204 : {
1205 : _dof_constraints.swap(_stashed_dof_constraints);
1206 : }
1207 :
1208 : #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
1209 : /**
1210 : * \returns An iterator pointing to the first Node constraint row.
1211 : */
1212 : NodeConstraints::const_iterator node_constraint_rows_begin() const
1213 : { return _node_constraints.begin(); }
1214 :
1215 : /**
1216 : * \returns An iterator pointing just past the last Node constraint row.
1217 : */
1218 : NodeConstraints::const_iterator node_constraint_rows_end() const
1219 : { return _node_constraints.end(); }
1220 : #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
1221 :
1222 : /**
1223 : * \returns \p true if the degree of freedom dof is constrained,
1224 : * \p false otherwise.
1225 : */
1226 : bool is_constrained_dof (const dof_id_type dof) const;
1227 :
1228 : /**
1229 : * \returns \p true if the system has any heterogeneous constraints for
1230 : * adjoint solution \p qoi_num, \p false otherwise.
1231 : */
1232 : bool has_heterogeneous_adjoint_constraints (const unsigned int qoi_num) const;
1233 :
1234 : /**
1235 : * Backwards compatibility with misspelling.
1236 : */
1237 4796 : bool has_heterogenous_adjoint_constraints (const unsigned int qoi_num) const
1238 : {
1239 83110 : return this->has_heterogeneous_adjoint_constraints (qoi_num);
1240 : }
1241 :
1242 : /**
1243 : * \returns The heterogeneous constraint value if the degree of
1244 : * freedom \p dof has a heterogeneous constraint for adjoint solution
1245 : * \p qoi_num, zero otherwise.
1246 : */
1247 : Number has_heterogeneous_adjoint_constraint (const unsigned int qoi_num,
1248 : const dof_id_type dof) const;
1249 :
1250 : /**
1251 : * Backwards compatibility with misspelling.
1252 : */
1253 1094256 : Number has_heterogenous_adjoint_constraint (const unsigned int qoi_num,
1254 : const dof_id_type dof) const
1255 : {
1256 1415236 : return this->has_heterogeneous_adjoint_constraint (qoi_num, dof);
1257 : }
1258 :
1259 : /**
1260 : * \returns A reference to the set of right-hand-side values in
1261 : * primal constraint equations
1262 : */
1263 : DofConstraintValueMap & get_primal_constraint_values();
1264 :
1265 : /**
1266 : * \returns \p true if the Node is constrained,
1267 : * false otherwise.
1268 : */
1269 : bool is_constrained_node (const Node * node) const;
1270 :
1271 : /**
1272 : * Prints (from processor 0) all DoF and Node constraints. If \p
1273 : * print_nonlocal is true, then each constraint is printed once for
1274 : * each processor that knows about it, which may be useful for \p
1275 : * DistributedMesh debugging.
1276 : */
1277 : void print_dof_constraints(std::ostream & os=libMesh::out,
1278 : bool print_nonlocal=false) const;
1279 :
1280 : /**
1281 : * Gets a string reporting all DoF and Node constraints local to
1282 : * this processor. If \p print_nonlocal is true, then nonlocal
1283 : * constraints which are locally known are included.
1284 : */
1285 : std::string get_local_constraints(bool print_nonlocal=false) const;
1286 :
1287 :
1288 : /**
1289 : * Tests the constrained degrees of freedom on the numeric vector \p v, which
1290 : * represents a solution defined on the mesh, returning a pair whose first
1291 : * entry is the maximum absolute error on a constrained DoF and whose second
1292 : * entry is the maximum relative error. Useful for debugging purposes.
1293 : *
1294 : * If \p v == nullptr, the system solution vector is tested.
1295 : */
1296 : std::pair<Real, Real> max_constraint_error(const System & system,
1297 : NumericVector<Number> * v = nullptr) const;
1298 :
1299 : #endif // LIBMESH_ENABLE_CONSTRAINTS
1300 :
1301 : //--------------------------------------------------------------------
1302 : // Constraint-specific methods
1303 : // Some of these methods are enabled (but inlined away to nothing)
1304 : // when constraints are disabled at configure-time. This is to
1305 : // increase API compatibility of user code with different library
1306 : // builds.
1307 :
1308 : /**
1309 : * Constrains the element matrix. This method requires the
1310 : * element matrix to be square, in which case the elem_dofs
1311 : * correspond to the global DOF indices of both the rows and
1312 : * columns of the element matrix. For this case the rows
1313 : * and columns of the matrix necessarily correspond to variables
1314 : * of the same approximation order.
1315 : *
1316 : * If \p asymmetric_constraint_rows is set to true (as it is by
1317 : * default), constraint row equations will be reinforced in a way
1318 : * which breaks matrix symmetry but makes inexact linear solver
1319 : * solutions more likely to satisfy hanging node constraints.
1320 : */
1321 : void constrain_element_matrix (DenseMatrix<Number> & matrix,
1322 : std::vector<dof_id_type> & elem_dofs,
1323 : bool asymmetric_constraint_rows = true) const;
1324 :
1325 : /**
1326 : * Constrains the element matrix. This method allows the
1327 : * element matrix to be non-square, in which case the row_dofs
1328 : * and col_dofs may be of different size and correspond to
1329 : * variables approximated in different spaces.
1330 : */
1331 : void constrain_element_matrix (DenseMatrix<Number> & matrix,
1332 : std::vector<dof_id_type> & row_dofs,
1333 : std::vector<dof_id_type> & col_dofs,
1334 : bool asymmetric_constraint_rows = true) const;
1335 :
1336 : /**
1337 : * Constrains the element vector.
1338 : */
1339 : void constrain_element_vector (DenseVector<Number> & rhs,
1340 : std::vector<dof_id_type> & dofs,
1341 : bool asymmetric_constraint_rows = true) const;
1342 :
1343 : /**
1344 : * Constrains the element matrix and vector. This method requires
1345 : * the element matrix to be square, in which case the elem_dofs
1346 : * correspond to the global DOF indices of both the rows and
1347 : * columns of the element matrix. For this case the rows
1348 : * and columns of the matrix necessarily correspond to variables
1349 : * of the same approximation order.
1350 : */
1351 : void constrain_element_matrix_and_vector (DenseMatrix<Number> & matrix,
1352 : DenseVector<Number> & rhs,
1353 : std::vector<dof_id_type> & elem_dofs,
1354 : bool asymmetric_constraint_rows = true) const;
1355 :
1356 : /**
1357 : * Constrains the element matrix and vector. This method requires
1358 : * the element matrix to be square, in which case the elem_dofs
1359 : * correspond to the global DOF indices of both the rows and
1360 : * columns of the element matrix. For this case the rows
1361 : * and columns of the matrix necessarily correspond to variables
1362 : * of the same approximation order.
1363 : *
1364 : * The heterogeneous version of this method creates linear systems in
1365 : * which heterogeneously constrained degrees of freedom will solve to
1366 : * their correct offset values, as would be appropriate for finding
1367 : * a solution to a linear problem in a single algebraic solve. The
1368 : * non-heterogeneous version of this method creates linear systems in
1369 : * which even heterogeneously constrained degrees of freedom are
1370 : * solved without offset values taken into account, as would be
1371 : * appropriate for finding linearized updates to a solution in which
1372 : * heterogeneous constraints are already satisfied.
1373 : *
1374 : * By default, the constraints for the primal solution of this
1375 : * system are used. If a non-negative \p qoi_index is passed in,
1376 : * then the constraints for the corresponding adjoint solution are
1377 : * used instead.
1378 : */
1379 : void heterogeneously_constrain_element_matrix_and_vector (DenseMatrix<Number> & matrix,
1380 : DenseVector<Number> & rhs,
1381 : std::vector<dof_id_type> & elem_dofs,
1382 : bool asymmetric_constraint_rows = true,
1383 : int qoi_index = -1) const;
1384 :
1385 : /*
1386 : * Backwards compatibility with misspelling.
1387 : */
1388 0 : void heterogenously_constrain_element_matrix_and_vector (DenseMatrix<Number> & matrix,
1389 : DenseVector<Number> & rhs,
1390 : std::vector<dof_id_type> & elem_dofs,
1391 : bool asymmetric_constraint_rows = true,
1392 : int qoi_index = -1) const
1393 : {
1394 : return this->heterogeneously_constrain_element_matrix_and_vector
1395 15527 : (matrix, rhs, elem_dofs, asymmetric_constraint_rows, qoi_index);
1396 : }
1397 :
1398 : /**
1399 : * Constrains the element vector. This method requires
1400 : * the element matrix to be square and not-yet-constrained, in which
1401 : * case the elem_dofs correspond to the global DOF indices of both
1402 : * the rows and columns of the element matrix.
1403 : *
1404 : * The heterogeneous version of this method creates linear systems in
1405 : * which heterogeneously constrained degrees of freedom will solve to
1406 : * their correct offset values, as would be appropriate for finding
1407 : * a solution to a linear problem in a single algebraic solve. The
1408 : * non-heterogeneous version of this method creates linear systems in
1409 : * which even heterogeneously constrained degrees of freedom are
1410 : * solved without offset values taken into account, as would be
1411 : * appropriate for finding linearized updates to a solution in which
1412 : * heterogeneous constraints are already satisfied.
1413 : *
1414 : * Note the sign difference from the nonlinear heterogeneous constraint
1415 : * method: Solving u:=K\f has the opposite sign convention from
1416 : * u:=u_in-J\r, and we apply heterogeneous constraints accordingly.
1417 : *
1418 : * By default, the constraints for the primal solution of this
1419 : * system are used. If a non-negative \p qoi_index is passed in,
1420 : * then the constraints for the corresponding adjoint solution are
1421 : * used instead.
1422 : */
1423 : void heterogeneously_constrain_element_vector (const DenseMatrix<Number> & matrix,
1424 : DenseVector<Number> & rhs,
1425 : std::vector<dof_id_type> & elem_dofs,
1426 : bool asymmetric_constraint_rows = true,
1427 : int qoi_index = -1) const;
1428 :
1429 : /*
1430 : * Backwards compatibility with misspelling.
1431 : */
1432 140 : void heterogenously_constrain_element_vector (const DenseMatrix<Number> & matrix,
1433 : DenseVector<Number> & rhs,
1434 : std::vector<dof_id_type> & elem_dofs,
1435 : bool asymmetric_constraint_rows = true,
1436 : int qoi_index = -1) const
1437 : {
1438 : return this->heterogeneously_constrain_element_vector
1439 1680 : (matrix, rhs, elem_dofs, asymmetric_constraint_rows, qoi_index);
1440 : }
1441 :
1442 : /**
1443 : * Constrains the element Jacobian and residual. The element
1444 : * Jacobian is square, and the elem_dofs should correspond to the
1445 : * global DOF indices of both the rows and columns of the element
1446 : * matrix.
1447 : *
1448 : * The residual-constraining version of this method creates linear
1449 : * systems in which heterogeneously constrained degrees of freedom
1450 : * create non-zero residual terms when not at their correct offset
1451 : * values, as would be appropriate for finding a solution to a
1452 : * nonlinear problem in a quasi-Newton solve.
1453 : *
1454 : * Note the sign difference from the linear heterogeneous constraint
1455 : * method: Solving u:=u_in-J\r has the opposite sign convention from
1456 : * u:=K\f, and we apply heterogeneous constraints accordingly.
1457 : *
1458 : * The \p solution vector passed in should be a serialized or
1459 : * ghosted primal solution
1460 : */
1461 : void heterogeneously_constrain_element_jacobian_and_residual (DenseMatrix<Number> & matrix,
1462 : DenseVector<Number> & rhs,
1463 : std::vector<dof_id_type> & elem_dofs,
1464 : NumericVector<Number> & solution_local) const;
1465 :
1466 : /**
1467 : * Constrains the element residual. The element Jacobian is square,
1468 : * and the elem_dofs should correspond to the global DOF indices of
1469 : * both the rows and columns of the element matrix.
1470 : *
1471 : * The residual-constraining version of this method creates linear
1472 : * systems in which heterogeneously constrained degrees of freedom
1473 : * create non-zero residual terms when not at their correct offset
1474 : * values, as would be appropriate for finding a solution to a
1475 : * nonlinear problem in a quasi-Newton solve.
1476 : *
1477 : * The \p solution vector passed in should be a serialized or
1478 : * ghosted primal solution
1479 : */
1480 : void heterogeneously_constrain_element_residual (DenseVector<Number> & rhs,
1481 : std::vector<dof_id_type> & elem_dofs,
1482 : NumericVector<Number> & solution_local) const;
1483 :
1484 :
1485 : /**
1486 : * Constrains the element residual. The element Jacobian is square,
1487 : * and the elem_dofs should correspond to the global DOF indices of
1488 : * both the rows and columns of the element matrix, and the dof
1489 : * constraint should not include any heterogeneous terms.
1490 : *
1491 : * The residual-constraining version of this method creates linear
1492 : * systems in which heterogeneously constrained degrees of freedom
1493 : * create non-zero residual terms when not at their correct offset
1494 : * values, as would be appropriate for finding a solution to a
1495 : * nonlinear problem in a quasi-Newton solve.
1496 : *
1497 : * The \p solution vector passed in should be a serialized or
1498 : * ghosted primal solution
1499 : */
1500 : void constrain_element_residual (DenseVector<Number> & rhs,
1501 : std::vector<dof_id_type> & elem_dofs,
1502 : NumericVector<Number> & solution_local) const;
1503 :
1504 : /**
1505 : * Constrains a dyadic element matrix B = v w'. This method
1506 : * requires the element matrix to be square, in which case the
1507 : * elem_dofs correspond to the global DOF indices of both the rows
1508 : * and columns of the element matrix. For this case the rows and
1509 : * columns of the matrix necessarily correspond to variables of the
1510 : * same approximation order.
1511 : */
1512 : void constrain_element_dyad_matrix (DenseVector<Number> & v,
1513 : DenseVector<Number> & w,
1514 : std::vector<dof_id_type> & row_dofs,
1515 : bool asymmetric_constraint_rows = true) const;
1516 :
1517 : /**
1518 : * Does not actually constrain anything, but modifies \p dofs in the
1519 : * same way as any of the constrain functions would do, i.e. adds
1520 : * those dofs in terms of which any of the existing dofs is
1521 : * constrained.
1522 : */
1523 : void constrain_nothing (std::vector<dof_id_type> & dofs) const;
1524 :
1525 : /**
1526 : * Constrains the numeric vector \p v, which represents a solution defined on
1527 : * the mesh. This may need to be used after a linear solve, if your linear
1528 : * solver's solutions do not satisfy your DoF constraints to a tight enough
1529 : * tolerance.
1530 : *
1531 : * If \p v == nullptr, the system solution vector is constrained
1532 : *
1533 : * If \p homogeneous == true, heterogeneous constraints are enforced
1534 : * as if they were homogeneous. This might be appropriate for e.g. a
1535 : * vector representing a difference between two
1536 : * heterogeneously-constrained solutions.
1537 : */
1538 : void enforce_constraints_exactly (const System & system,
1539 : NumericVector<Number> * v = nullptr,
1540 : bool homogeneous = false) const;
1541 :
1542 : /**
1543 : * Heterogeneously constrains the numeric vector \p v, which
1544 : * represents an adjoint solution defined on the mesh for quantity
1545 : * fo interest \p q. For homogeneous constraints, use \p
1546 : * enforce_constraints_exactly instead
1547 : */
1548 : void enforce_adjoint_constraints_exactly (NumericVector<Number> & v,
1549 : unsigned int q) const;
1550 :
1551 : void enforce_constraints_on_residual (const NonlinearImplicitSystem & system,
1552 : NumericVector<Number> * rhs,
1553 : NumericVector<Number> const * solution,
1554 : bool homogeneous = true) const;
1555 :
1556 : void enforce_constraints_on_jacobian (const NonlinearImplicitSystem & system,
1557 : SparseMatrix<Number> * jac) const;
1558 :
1559 : #ifdef LIBMESH_ENABLE_PERIODIC
1560 :
1561 : //--------------------------------------------------------------------
1562 : // PeriodicBoundary-specific methods
1563 :
1564 : /**
1565 : * Adds a copy of the specified periodic boundary to the system.
1566 : */
1567 : void add_periodic_boundary (const PeriodicBoundaryBase & periodic_boundary);
1568 :
1569 : /**
1570 : * Add a periodic boundary pair
1571 : *
1572 : * \param boundary - primary boundary
1573 : * \param inverse_boundary - inverse boundary
1574 : */
1575 : void add_periodic_boundary (const PeriodicBoundaryBase & boundary, const PeriodicBoundaryBase & inverse_boundary);
1576 :
1577 : /**
1578 : * \returns \p true if the boundary given by \p boundaryid is periodic,
1579 : * false otherwise
1580 : */
1581 : bool is_periodic_boundary (const boundary_id_type boundaryid) const;
1582 :
1583 : PeriodicBoundaries * get_periodic_boundaries()
1584 : {
1585 : return _periodic_boundaries.get();
1586 : }
1587 :
1588 : const PeriodicBoundaries * get_periodic_boundaries() const
1589 : {
1590 : return _periodic_boundaries.get();
1591 : }
1592 :
1593 : #endif // LIBMESH_ENABLE_PERIODIC
1594 :
1595 :
1596 : #ifdef LIBMESH_ENABLE_DIRICHLET
1597 :
1598 : //--------------------------------------------------------------------
1599 : // DirichletBoundary-specific methods
1600 :
1601 : /**
1602 : * Adds a copy of the specified Dirichlet boundary to the system.
1603 : *
1604 : * The constraints implied by DirichletBoundary objects are imposed
1605 : * in the same order in which DirichletBoundary objects are added to
1606 : * the DofMap. When multiple DirichletBoundary objects would impose
1607 : * competing constraints on a given DOF, the *first*
1608 : * DirichletBoundary to constrain the DOF "wins". This distinction
1609 : * is important when e.g. two surfaces (sidesets) intersect. The
1610 : * nodes on the intersection will be constrained according to
1611 : * whichever sideset's DirichletBoundary object was added to the
1612 : * DofMap first.
1613 : */
1614 : void add_dirichlet_boundary (const DirichletBoundary & dirichlet_boundary);
1615 :
1616 : /**
1617 : * Adds a copy of the specified Dirichlet boundary to the system,
1618 : * corresponding to the adjoint problem defined by Quantity of
1619 : * Interest \p q.
1620 : */
1621 : void add_adjoint_dirichlet_boundary (const DirichletBoundary & dirichlet_boundary,
1622 : unsigned int q);
1623 :
1624 : /**
1625 : * Removes the specified Dirichlet boundary from the system.
1626 : */
1627 : void remove_dirichlet_boundary (const DirichletBoundary & dirichlet_boundary);
1628 :
1629 : /**
1630 : * Removes from the system the specified Dirichlet boundary for the
1631 : * adjoint equation defined by Quantity of interest index q
1632 : */
1633 : void remove_adjoint_dirichlet_boundary (const DirichletBoundary & dirichlet_boundary,
1634 : unsigned int q);
1635 :
1636 : const DirichletBoundaries * get_dirichlet_boundaries() const
1637 : {
1638 : return _dirichlet_boundaries.get();
1639 : }
1640 :
1641 24 : DirichletBoundaries * get_dirichlet_boundaries()
1642 : {
1643 24 : return _dirichlet_boundaries.get();
1644 : }
1645 :
1646 : bool has_adjoint_dirichlet_boundaries(unsigned int q) const;
1647 :
1648 : const DirichletBoundaries *
1649 : get_adjoint_dirichlet_boundaries(unsigned int q) const;
1650 :
1651 : DirichletBoundaries *
1652 : get_adjoint_dirichlet_boundaries(unsigned int q);
1653 :
1654 : /**
1655 : * Check that all the ids in dirichlet_bcids are actually present in the mesh.
1656 : * If not, this will throw an error.
1657 : */
1658 : void check_dirichlet_bcid_consistency (const MeshBase & mesh,
1659 : const DirichletBoundary & boundary) const;
1660 : #endif // LIBMESH_ENABLE_DIRICHLET
1661 :
1662 :
1663 : #ifdef LIBMESH_ENABLE_AMR
1664 :
1665 : //--------------------------------------------------------------------
1666 : // AMR-specific methods
1667 :
1668 : /**
1669 : * After a mesh is refined and repartitioned it is possible that the
1670 : * \p _send_list will need to be augmented. This is the case when an
1671 : * element is refined and its children end up on different processors
1672 : * than the parent. These children will need values from the parent
1673 : * when projecting the solution onto the refined mesh, hence the parent's
1674 : * DOF indices need to be included in the \p _send_list.
1675 : */
1676 : // void augment_send_list_for_projection(const MeshBase &);
1677 :
1678 : #ifdef LIBMESH_ENABLE_AMR
1679 :
1680 : /**
1681 : * Fills the vector di with the global degree of freedom indices
1682 : * for the element using the \p DofMap::old_dof_object.
1683 : * If no variable number is specified then all
1684 : * variables are returned.
1685 : */
1686 : void old_dof_indices (const Elem * const elem,
1687 : std::vector<dof_id_type> & di,
1688 : const unsigned int vn = libMesh::invalid_uint) const;
1689 :
1690 : #endif // LIBMESH_ENABLE_AMR
1691 :
1692 : /**
1693 : * Constrains degrees of freedom on side \p s of element \p elem which
1694 : * correspond to variable number \p var and to p refinement levels
1695 : * above \p p.
1696 : */
1697 : void constrain_p_dofs (unsigned int var,
1698 : const Elem * elem,
1699 : unsigned int s,
1700 : unsigned int p);
1701 :
1702 : #endif // LIBMESH_ENABLE_AMR
1703 :
1704 : /**
1705 : * Reinitialize the underlying data structures conformal to the current mesh.
1706 : */
1707 : void reinit
1708 : (MeshBase & mesh,
1709 : const std::map<const Node *, std::set<subdomain_id_type>> &
1710 : constraining_subdomains);
1711 :
1712 : /**
1713 : * Free all new memory associated with the object, but restore its
1714 : * original state, with the mesh pointer and any default ghosting.
1715 : */
1716 : virtual void clear () override;
1717 :
1718 : /**
1719 : * Prints summary info about the sparsity bandwidth and constraints.
1720 : */
1721 : void print_info(std::ostream & os=libMesh::out) const;
1722 :
1723 : /**
1724 : * Gets summary info about the sparsity bandwidth and constraints.
1725 : */
1726 : std::string get_info() const;
1727 :
1728 : /**
1729 : * Degree of freedom coupling. If left empty each DOF
1730 : * couples to all others. Can be used to reduce memory
1731 : * requirements for sparse matrices. DOF 0 might only
1732 : * couple to itself, in which case \p dof_coupling(0,0)
1733 : * should be 1 and \p dof_coupling(0,j) = 0 for j not equal
1734 : * to 0.
1735 : *
1736 : * This variable is named as though it were class private,
1737 : * but it is in the public interface. Also there are no
1738 : * public methods for accessing it... This typically means
1739 : * you should only use it if you know what you are doing.
1740 : */
1741 : CouplingMatrix * _dof_coupling;
1742 :
1743 : /**
1744 : * \returns The number of the system we are responsible for.
1745 : */
1746 : unsigned int sys_number() const;
1747 :
1748 : /**
1749 : * Builds a sparsity pattern for matrices using the current
1750 : * degree-of-freedom numbering and coupling.
1751 : *
1752 : * By default, ignores constraint equations, for build speed; this
1753 : * is valid for the combination of !need_full_sparsity_pattern and
1754 : * constraints which only come from periodic boundary conditions and
1755 : * adaptive mesh refinement, where matrix constraint adds some
1756 : * matrix entries but removes equally many (or more) other entries.
1757 : *
1758 : * Can be told to calculate sparsity for the constrained matrix,
1759 : * which may be necessary in the case of spline control node
1760 : * constraints or sufficiently many user constraints.
1761 : */
1762 : std::unique_ptr<SparsityPattern::Build> build_sparsity(const MeshBase & mesh,
1763 : bool calculate_constrained = false,
1764 : bool use_condensed_system = false) const;
1765 :
1766 : /**
1767 : * Describe whether the given variable group should be p-refined. If this API is not called with
1768 : * \p false, the default is to p-refine
1769 : */
1770 : void should_p_refine(unsigned int g, bool p_refine);
1771 :
1772 : /**
1773 : * Whether the given variable group should be p-refined
1774 : */
1775 : bool should_p_refine(unsigned int g) const;
1776 :
1777 : /**
1778 : * Whether the given variable should be p-refined
1779 : */
1780 : bool should_p_refine_var(unsigned int var) const;
1781 :
1782 : // Prevent bad user implicit conversions
1783 : void should_p_refine(FEFamily, bool) = delete;
1784 : void should_p_refine(Order, bool) = delete;
1785 : bool should_p_refine(FEFamily) const = delete;
1786 : bool should_p_refine(Order) const = delete;
1787 :
1788 : /**
1789 : * Add a static condensation class
1790 : */
1791 : void create_static_condensation(MeshBase & mesh, System & system);
1792 :
1793 : /**
1794 : * Checks whether we have static condensation
1795 : */
1796 52428 : bool has_static_condensation() const { return _sc.get(); }
1797 :
1798 : /**
1799 : * @returns the static condensation class. This should have been already added with a call to \p
1800 : * add_static_condensation()
1801 : */
1802 : StaticCondensationDofMap & get_static_condensation();
1803 :
1804 : /**
1805 : * Calls reinit on the static condensation map if it exists
1806 : */
1807 : void reinit_static_condensation();
1808 :
1809 : private:
1810 :
1811 : /**
1812 : * Retrieve the array variable bounds for a given variable \p vi. This variable may
1813 : * lie anywhere within an array variable range. An 'array variable' is simply a sequence
1814 : * of contiguous variable numbers defined by pair where the first member of the pair
1815 : * is the first number in the variable sequence and the second member of the pair is
1816 : * the number of the last variable in the sequence plus one. Array variables may be
1817 : * used in tandem with variable grouping by downstream code to build optimized physics
1818 : * kernels since each variable in the array will have the same shape functions.
1819 : *
1820 : * We note that we store array variables as a container of the above described pairs. Within
1821 : * this API we will do a binary search such that the complexity is O(log(N)) where N is the
1822 : * number of array variables present in \p this
1823 : */
1824 : const std::pair<unsigned int, unsigned int> &
1825 : get_variable_array(unsigned int vi) const;
1826 :
1827 : /**
1828 : * Helper function that gets the dof indices on the current element
1829 : * for a non-SCALAR type variable, where the variable is identified
1830 : * by its variable group number \p vg and its offset \p vig from the
1831 : * first variable in that group.
1832 : *
1833 : * In DEBUG mode, the tot_size parameter will add up the total
1834 : * number of dof indices that should have been added to di, and v
1835 : * will be the variable number corresponding to vg and vig.
1836 : */
1837 : void _dof_indices (const Elem & elem,
1838 : int p_level,
1839 : std::vector<dof_id_type> & di,
1840 : const unsigned int vg,
1841 : const unsigned int vig,
1842 : const Node * const * nodes,
1843 : unsigned int n_nodes,
1844 : const unsigned int v
1845 : #ifdef DEBUG
1846 : ,
1847 : std::size_t & tot_size
1848 : #endif
1849 : ) const;
1850 :
1851 : /**
1852 : * As above except a \p field_dofs_functor must be provided. This method is useful when the caller
1853 : * wants to do more than simply fill a degree of freedom container
1854 : * @param field_dofs_functor This functor has the interface:
1855 : * void field_dofs_functor(const Elem & elem,
1856 : * const unsigned int node_num,
1857 : * const unsigned int var_num,
1858 : * std::vector<dof_id_type> & di,
1859 : * const dof_id_type field_dof)
1860 : * where \p field_dof represents a field degree of freedom to act on and
1861 : * is associated with \p node_num and \p var_num. If the degree of
1862 : * freedom is elemental than \p node_num will be \p invalid_uint. \p di
1863 : * is the degree of freedom container provided to the \p _dof_indices
1864 : * method
1865 : */
1866 : template <typename FieldDofsFunctor>
1867 : void _dof_indices (const Elem & elem,
1868 : int p_level,
1869 : std::vector<dof_id_type> & di,
1870 : const unsigned int vg,
1871 : const unsigned int vig,
1872 : const Node * const * nodes,
1873 : unsigned int n_nodes,
1874 : const unsigned int v,
1875 : #ifdef DEBUG
1876 : std::size_t & tot_size,
1877 : #endif
1878 : FieldDofsFunctor field_dofs_functor) const;
1879 :
1880 : /**
1881 : * Helper function that implements the element-nodal versions of
1882 : * dof_indices and old_dof_indices
1883 : */
1884 : void _node_dof_indices (const Elem & elem,
1885 : unsigned int n,
1886 : const DofObject & obj,
1887 : std::vector<dof_id_type> & di,
1888 : const unsigned int vn) const;
1889 :
1890 : /**
1891 : * Invalidates all active DofObject dofs for this system
1892 : */
1893 : void invalidate_dofs(MeshBase & mesh) const;
1894 :
1895 : /**
1896 : * \returns The Node pointer with index \p i from the \p mesh.
1897 : */
1898 : DofObject * node_ptr(MeshBase & mesh, dof_id_type i) const;
1899 :
1900 : /**
1901 : * \returns The Elem pointer with index \p i from the \p mesh.
1902 : */
1903 : DofObject * elem_ptr(MeshBase & mesh, dof_id_type i) const;
1904 :
1905 : /**
1906 : * A member function type like \p node_ptr() or \p elem_ptr().
1907 : */
1908 : typedef DofObject * (DofMap::*dofobject_accessor)
1909 : (MeshBase & mesh, dof_id_type i) const;
1910 :
1911 : /**
1912 : * Helper function for distributing dofs in parallel
1913 : */
1914 : template<typename iterator_type>
1915 : void set_nonlocal_dof_objects(iterator_type objects_begin,
1916 : iterator_type objects_end,
1917 : MeshBase & mesh,
1918 : dofobject_accessor objects);
1919 :
1920 : /**
1921 : * We may have mesh constraint rows with dependent nodes in one
1922 : * subdomain but dependency nodes in another subdomain, and we may
1923 : * have variables whose subdomain restriction includes the dependent
1924 : * subdomain but not the dependency. In those cases we need to
1925 : * place degrees of freedom on dependency nodes anyway.
1926 : *
1927 : * The set value for node n will include all subdomain ids of
1928 : * elements with nodes in subdomains constrained by n.
1929 : *
1930 : * We use a map<set> rather than a multimap here because we expect
1931 : * to be inserting the same subdomain multiple times and we don't
1932 : * need duplicate values.
1933 : */
1934 : std::map<const Node *, std::set<subdomain_id_type>>
1935 : calculate_constraining_subdomains();
1936 :
1937 : /**
1938 : * Distributes the global degrees of freedom, for dofs on
1939 : * this processor. In this format the local
1940 : * degrees of freedom are in a contiguous block for each
1941 : * variable in the system.
1942 : * Starts at index next_free_dof, and increments it to
1943 : * the post-final index.
1944 : *
1945 : * Uses the provided constraining_subdomains map from
1946 : * calculate_constraining_subdomains() to ensure allocation of all
1947 : * DoFs on constraining nodes.
1948 : */
1949 : void distribute_local_dofs_var_major
1950 : (dof_id_type & next_free_dof,
1951 : MeshBase & mesh,
1952 : const std::map<const Node *, std::set<subdomain_id_type>> &
1953 : constraining_subdomains);
1954 :
1955 : /**
1956 : * Distributes the global degrees of freedom for dofs on this
1957 : * processor. In this format all the degrees of freedom at a
1958 : * node/element are in contiguous blocks. Starts at index \p
1959 : * next_free_dof, and increments it to the post-final index. If \p
1960 : * build_send_list is \p true, builds the send list. If \p false,
1961 : * clears and reserves the send list.
1962 : *
1963 : * Uses the provided constraining_subdomains map from
1964 : * calculate_constraining_subdomains() to ensure allocation of all
1965 : * DoFs on constraining nodes.
1966 : *
1967 : * \note The degrees of freedom for a given variable are not in
1968 : * contiguous blocks, as in the case of \p distribute_local_dofs_var_major.
1969 : */
1970 : void distribute_local_dofs_node_major
1971 : (dof_id_type & next_free_dof,
1972 : MeshBase & mesh,
1973 : const std::map<const Node *, std::set<subdomain_id_type>> &
1974 : constraining_subdomains);
1975 :
1976 : /*
1977 : * Helper method for the above two to count + distriubte SCALAR dofs
1978 : */
1979 : void distribute_scalar_dofs (dof_id_type & next_free_dof);
1980 :
1981 : #ifdef DEBUG
1982 : /*
1983 : * Internal assertions for distribute_local_dofs_*
1984 : */
1985 : void assert_no_nodes_missed(MeshBase & mesh);
1986 : #endif
1987 :
1988 : /*
1989 : * A utility method for obtaining a set of elements to ghost along
1990 : * with merged coupling matrices.
1991 : */
1992 : typedef std::set<std::unique_ptr<CouplingMatrix>, Utility::CompareUnderlying> CouplingMatricesSet;
1993 : static void
1994 : merge_ghost_functor_outputs (GhostingFunctor::map_type & elements_to_ghost,
1995 : CouplingMatricesSet & temporary_coupling_matrices,
1996 : const GhostingFunctorIterator & gf_begin,
1997 : const GhostingFunctorIterator & gf_end,
1998 : const MeshBase::const_element_iterator & elems_begin,
1999 : const MeshBase::const_element_iterator & elems_end,
2000 : processor_id_type p);
2001 :
2002 : /**
2003 : * Adds entries to the \p _send_list vector corresponding to DoFs
2004 : * on elements neighboring the current processor.
2005 : */
2006 : void add_neighbors_to_send_list(MeshBase & mesh);
2007 :
2008 : #ifdef LIBMESH_ENABLE_CONSTRAINTS
2009 :
2010 : /**
2011 : * Build the constraint matrix C associated with the element
2012 : * degree of freedom indices elem_dofs. The optional parameter
2013 : * \p called_recursively should be left at the default value
2014 : * \p false. This is used to handle the special case of
2015 : * an element's degrees of freedom being constrained in terms
2016 : * of other, local degrees of freedom. The usual case is
2017 : * for an elements DOFs to be constrained by some other,
2018 : * external DOFs.
2019 : */
2020 : void build_constraint_matrix (DenseMatrix<Number> & C,
2021 : std::vector<dof_id_type> & elem_dofs,
2022 : const bool called_recursively=false) const;
2023 :
2024 : /**
2025 : * Build the constraint matrix C and the forcing vector H
2026 : * associated with the element degree of freedom indices elem_dofs.
2027 : * The optional parameter \p called_recursively should be left at
2028 : * the default value \p false. This is used to handle the special
2029 : * case of an element's degrees of freedom being constrained in
2030 : * terms of other, local degrees of freedom. The usual case is for
2031 : * an elements DOFs to be constrained by some other, external DOFs
2032 : * and/or Dirichlet conditions.
2033 : *
2034 : * The forcing vector will depend on which solution's heterogeneous
2035 : * constraints are being applied. For the default \p qoi_index this
2036 : * will be the primal solution; for \p qoi_index >= 0 the
2037 : * corresponding adjoint solution's constraints will be used.
2038 : */
2039 : void build_constraint_matrix_and_vector (DenseMatrix<Number> & C,
2040 : DenseVector<Number> & H,
2041 : std::vector<dof_id_type> & elem_dofs,
2042 : int qoi_index = -1,
2043 : const bool called_recursively=false) const;
2044 :
2045 : /**
2046 : * Finds all the DOFS associated with the element DOFs elem_dofs.
2047 : * This will account for off-element couplings via hanging nodes.
2048 : */
2049 : void find_connected_dofs (std::vector<dof_id_type> & elem_dofs) const;
2050 :
2051 : /**
2052 : * Finds all the DofObjects associated with the set in \p objs.
2053 : * This will account for off-element couplings via hanging nodes.
2054 : */
2055 : void find_connected_dof_objects (std::vector<const DofObject *> & objs) const;
2056 :
2057 : /**
2058 : * Adds entries to the \p _send_list vector corresponding to DoFs
2059 : * which are dependencies for constraint equations on the current
2060 : * processor.
2061 : */
2062 : void add_constraints_to_send_list();
2063 :
2064 : /**
2065 : * Adds any spline constraints from the Mesh to our DoF constraints.
2066 : * If any Dirichlet constraints exist on spline-constrained nodes,
2067 : * l2-projects those constraints onto the spline basis.
2068 : */
2069 : void process_mesh_constraint_rows(const MeshBase & mesh);
2070 :
2071 : #endif // LIBMESH_ENABLE_CONSTRAINTS
2072 :
2073 : /**
2074 : * This flag indicates whether or not we do an opt-mode check for
2075 : * the presence of constraint loops, i.e. cases where the constraint
2076 : * graph is cyclic.
2077 : */
2078 : bool _error_on_constraint_loop;
2079 :
2080 : /**
2081 : * This flag indicates whether or not we explicitly take constraint
2082 : * equations into account when computing a sparsity pattern.
2083 : */
2084 : bool _constrained_sparsity_construction;
2085 :
2086 : /**
2087 : * The variables in this system/degree of freedom map
2088 : */
2089 : std::vector<Variable> _variables;
2090 :
2091 : /**
2092 : * The variable groups in this system/degree of freedom map
2093 : */
2094 : std::vector<VariableGroup> _variable_groups;
2095 :
2096 : /**
2097 : * The variable group number for each variable.
2098 : */
2099 : std::vector<unsigned int> _variable_group_numbers;
2100 :
2101 : /**
2102 : * A map from variable number to variable group number
2103 : */
2104 : std::unordered_map<unsigned int, unsigned int> _var_to_vg;
2105 :
2106 : /**
2107 : * The variable numbers corresponding to user-specified
2108 : * names, useful for name-based lookups.
2109 : */
2110 : std::map<std::string, unsigned int, std::less<>> _variable_numbers;
2111 :
2112 : /**
2113 : * Array variable information storage. For a given array "variable", the first member of the pair
2114 : * denotes the first variable number present in the array variable and the second member of the
2115 : * pair denotes the last variable number present in the array variables plus one
2116 : */
2117 : std::vector<std::pair<unsigned int, unsigned int>> _array_variables;
2118 :
2119 : /**
2120 : * \p true when \p VariableGroup structures should be automatically
2121 : * identified, \p false otherwise. Defaults to \p true.
2122 : */
2123 : bool _identify_variable_groups = true;
2124 :
2125 : /**
2126 : * The number of the system we manage DOFs for.
2127 : */
2128 : const unsigned int _sys_number;
2129 :
2130 : /**
2131 : * The mesh that system uses.
2132 : */
2133 : MeshBase & _mesh;
2134 :
2135 : /**
2136 : * Additional matrices handled by this object. These pointers do \e
2137 : * not handle the memory, instead, \p System, who
2138 : * told \p DofMap about them, owns them.
2139 : */
2140 : std::vector<SparseMatrix<Number> * > _matrices;
2141 :
2142 : /**
2143 : * First DOF index for SCALAR variable v, or garbage for non-SCALAR
2144 : * variable v
2145 : */
2146 : std::vector<dof_id_type> _first_scalar_df;
2147 :
2148 : /**
2149 : * A list containing all the global DOF indices that affect the
2150 : * solution on my processor.
2151 : */
2152 : std::vector<dof_id_type> _send_list;
2153 :
2154 : /**
2155 : * Function object to call to add extra entries to the sparsity pattern
2156 : */
2157 : SparsityPattern::AugmentSparsityPattern * _augment_sparsity_pattern;
2158 :
2159 : /**
2160 : * A function pointer to a function to call to add extra entries to the sparsity pattern
2161 : */
2162 : void (*_extra_sparsity_function)(SparsityPattern::Graph &,
2163 : std::vector<dof_id_type> & n_nz,
2164 : std::vector<dof_id_type> & n_oz,
2165 : void *);
2166 : /**
2167 : * A pointer associated with the extra sparsity that can optionally be passed in
2168 : */
2169 : void * _extra_sparsity_context;
2170 :
2171 : /**
2172 : * Function object to call to add extra entries to the send list
2173 : */
2174 : AugmentSendList * _augment_send_list;
2175 :
2176 : /**
2177 : * A function pointer to a function to call to add extra entries to the send list
2178 : */
2179 : void (*_extra_send_list_function)(std::vector<dof_id_type> &, void *);
2180 :
2181 : /**
2182 : * A pointer associated with the extra send list that can optionally be passed in
2183 : */
2184 : void * _extra_send_list_context;
2185 :
2186 : /**
2187 : * The default coupling GhostingFunctor, used to implement standard
2188 : * libMesh sparsity pattern construction.
2189 : *
2190 : * We use a std::unique_ptr here to reduce header dependencies.
2191 : */
2192 : std::unique_ptr<DefaultCoupling> _default_coupling;
2193 :
2194 : /**
2195 : * The default algebraic GhostingFunctor, used to implement standard
2196 : * libMesh send_list construction.
2197 : *
2198 : * We use a std::unique_ptr here to reduce header dependencies.
2199 : */
2200 : std::unique_ptr<DefaultCoupling> _default_evaluating;
2201 :
2202 : /**
2203 : * The list of all GhostingFunctor objects to be used when
2204 : * distributing ghosted vectors.
2205 : *
2206 : * The library should automatically refer these functors to the
2207 : * MeshBase, too, so any algebraically ghosted dofs will live on
2208 : * geometrically ghosted elements.
2209 : *
2210 : * Keep these in a vector so any parallel computation is done in the
2211 : * same order on all processors.
2212 : */
2213 : std::vector<GhostingFunctor *> _algebraic_ghosting_functors;
2214 :
2215 : /**
2216 : * The list of all GhostingFunctor objects to be used when
2217 : * coupling degrees of freedom in matrix sparsity patterns.
2218 : *
2219 : * These objects will *also* be used as algebraic ghosting functors,
2220 : * but not vice-versa.
2221 : *
2222 : * The library should automatically refer these functors to the
2223 : * MeshBase, too, so any dofs coupled to local dofs will live on
2224 : * geometrically ghosted elements.
2225 : */
2226 : std::vector<GhostingFunctor *> _coupling_functors;
2227 :
2228 : /**
2229 : * Hang on to references to any GhostingFunctor objects we were
2230 : * passed in shared_ptr form
2231 : */
2232 : std::map<GhostingFunctor *, std::shared_ptr<GhostingFunctor> > _shared_functors;
2233 :
2234 : /**
2235 : * Default false; set to true if any attached matrix requires a full
2236 : * sparsity pattern.
2237 : */
2238 : bool need_full_sparsity_pattern;
2239 :
2240 : /**
2241 : * The sparsity pattern of the global matrix. If
2242 : * need_full_sparsity_pattern is true, we save the entire sparse
2243 : * graph here. Otherwise we save just the n_nz and n_oz vectors.
2244 : */
2245 : std::unique_ptr<SparsityPattern::Build> _sp;
2246 :
2247 : /**
2248 : * The total number of SCALAR dofs associated to
2249 : * all SCALAR variables.
2250 : */
2251 : dof_id_type _n_SCALAR_dofs;
2252 :
2253 : #ifdef LIBMESH_ENABLE_AMR
2254 :
2255 : /**
2256 : * First old DOF index for SCALAR variable v, or garbage for
2257 : * non-SCALAR variable v
2258 : */
2259 : std::vector<dof_id_type> _first_old_scalar_df;
2260 :
2261 : /**
2262 : * A container of variable groups that we should not p-refine
2263 : */
2264 : std::unordered_set<unsigned int> _dont_p_refine;
2265 : #endif
2266 :
2267 : #ifdef LIBMESH_ENABLE_CONSTRAINTS
2268 : /**
2269 : * Data structure containing DOF constraints. The ith
2270 : * entry is the constraint matrix row for DOF i.
2271 : */
2272 : DofConstraints _dof_constraints, _stashed_dof_constraints;
2273 :
2274 : DofConstraintValueMap _primal_constraint_values;
2275 :
2276 : AdjointDofConstraintValues _adjoint_constraint_values;
2277 : #endif
2278 :
2279 : #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
2280 : /**
2281 : * Data structure containing DofObject constraints.
2282 : */
2283 : NodeConstraints _node_constraints;
2284 : #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
2285 :
2286 :
2287 : #ifdef LIBMESH_ENABLE_PERIODIC
2288 : /**
2289 : * Data structure containing periodic boundaries. The ith
2290 : * entry is the constraint matrix row for boundaryid i.
2291 : */
2292 : std::unique_ptr<PeriodicBoundaries> _periodic_boundaries;
2293 : #endif
2294 :
2295 : #ifdef LIBMESH_ENABLE_DIRICHLET
2296 : /**
2297 : * Data structure containing Dirichlet functions. The ith
2298 : * entry is the constraint matrix row for boundaryid i.
2299 : */
2300 : std::unique_ptr<DirichletBoundaries> _dirichlet_boundaries;
2301 :
2302 : /**
2303 : * Data structure containing Dirichlet functions. The ith
2304 : * entry is the constraint matrix row for boundaryid i.
2305 : */
2306 : std::vector<std::unique_ptr<DirichletBoundaries>> _adjoint_dirichlet_boundaries;
2307 : #endif
2308 :
2309 : friend class SparsityPattern::Build;
2310 :
2311 : /**
2312 : * Bools to indicate if we override the --implicit_neighbor_dofs
2313 : * commandline options.
2314 : */
2315 : bool _implicit_neighbor_dofs_initialized;
2316 : bool _implicit_neighbor_dofs;
2317 :
2318 : /**
2319 : * Flag which determines whether we should do some additional
2320 : * checking of the consistency of the DirichletBoundary objects
2321 : * added by the user. Defaults to true, but can be disabled in cases
2322 : * where you only want to add DirichletBoundary objects "locally"
2323 : * and can guarantee that no repartitioning will be done, since
2324 : * repartitioning could cause processors to own new boundary sides
2325 : * for which they no longer have the proper DirichletBoundary
2326 : * objects stored.
2327 : */
2328 : bool _verify_dirichlet_bc_consistency;
2329 :
2330 : /// Static condensation class
2331 : std::unique_ptr<StaticCondensationDofMap> _sc;
2332 : };
2333 :
2334 :
2335 : // ------------------------------------------------------------
2336 : // Dof Map inline member functions
2337 : inline
2338 43379400 : unsigned int DofMap::sys_number() const
2339 : {
2340 300508730 : return _sys_number;
2341 : }
2342 :
2343 :
2344 :
2345 : inline
2346 48696766 : const VariableGroup & DofMap::variable_group (const unsigned int g) const
2347 : {
2348 48696766 : libmesh_assert_less (g, _variable_groups.size());
2349 :
2350 308903113 : return _variable_groups[g];
2351 : }
2352 :
2353 :
2354 :
2355 : inline
2356 63782684 : const Variable & DofMap::variable (const unsigned int c) const
2357 : {
2358 5606761 : libmesh_assert_less (c, _variables.size());
2359 :
2360 68740831 : return _variables[c];
2361 : }
2362 :
2363 :
2364 :
2365 : inline
2366 : Order DofMap::variable_order (const unsigned int c) const
2367 : {
2368 : libmesh_assert_less (c, _variables.size());
2369 :
2370 : return _variables[c].type().order;
2371 : }
2372 :
2373 :
2374 :
2375 : inline
2376 : Order DofMap::variable_group_order (const unsigned int vg) const
2377 : {
2378 : libmesh_assert_less (vg, _variable_groups.size());
2379 :
2380 : return _variable_groups[vg].type().order;
2381 : }
2382 :
2383 :
2384 :
2385 : inline
2386 2041198 : const FEType & DofMap::variable_type (const unsigned int c) const
2387 : {
2388 2041198 : libmesh_assert_less (c, _variables.size());
2389 :
2390 14885259 : return _variables[c].type();
2391 : }
2392 :
2393 :
2394 :
2395 : inline
2396 : const FEType & DofMap::variable_group_type (const unsigned int vg) const
2397 : {
2398 : libmesh_assert_less (vg, _variable_groups.size());
2399 :
2400 : return _variable_groups[vg].type();
2401 : }
2402 :
2403 :
2404 : #ifdef LIBMESH_ENABLE_CONSTRAINTS
2405 :
2406 :
2407 : inline
2408 1176677 : bool DofMap::is_constrained_node (const Node *
2409 : #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
2410 : node
2411 : #endif
2412 : ) const
2413 : {
2414 : #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
2415 1176677 : if (_node_constraints.count(node))
2416 13072 : return true;
2417 : #endif
2418 :
2419 1163605 : return false;
2420 : }
2421 :
2422 :
2423 : inline
2424 50804441 : bool DofMap::is_constrained_dof (const dof_id_type dof) const
2425 : {
2426 50804441 : if (_dof_constraints.count(dof))
2427 5419003 : return true;
2428 :
2429 45385438 : return false;
2430 : }
2431 :
2432 :
2433 : inline
2434 83110 : bool DofMap::has_heterogeneous_adjoint_constraints (const unsigned int qoi_num) const
2435 : {
2436 : AdjointDofConstraintValues::const_iterator it =
2437 4796 : _adjoint_constraint_values.find(qoi_num);
2438 87906 : if (it == _adjoint_constraint_values.end())
2439 1892 : return false;
2440 29012 : if (it->second.empty())
2441 25920 : return false;
2442 :
2443 24 : return true;
2444 : }
2445 :
2446 :
2447 : inline
2448 1415236 : Number DofMap::has_heterogeneous_adjoint_constraint (const unsigned int qoi_num,
2449 : const dof_id_type dof) const
2450 : {
2451 : AdjointDofConstraintValues::const_iterator it =
2452 1094256 : _adjoint_constraint_values.find(qoi_num);
2453 1415236 : if (it != _adjoint_constraint_values.end())
2454 : {
2455 : DofConstraintValueMap::const_iterator rhsit =
2456 104428 : it->second.find(dof);
2457 425408 : if (rhsit == it->second.end())
2458 154302 : return 0;
2459 : else
2460 1820 : return rhsit->second;
2461 : }
2462 :
2463 989828 : return 0;
2464 : }
2465 :
2466 :
2467 :
2468 : inline
2469 : DofConstraintValueMap & DofMap::get_primal_constraint_values()
2470 : {
2471 : return _primal_constraint_values;
2472 : }
2473 :
2474 :
2475 :
2476 : #else
2477 :
2478 : //--------------------------------------------------------------------
2479 : // Constraint-specific methods get inlined into nothing if
2480 : // constraints are disabled, so there's no reason for users not to
2481 : // use them.
2482 :
2483 : inline void DofMap::constrain_element_matrix (DenseMatrix<Number> &,
2484 : std::vector<dof_id_type> &,
2485 : bool) const {}
2486 :
2487 : inline void DofMap::constrain_element_matrix (DenseMatrix<Number> &,
2488 : std::vector<dof_id_type> &,
2489 : std::vector<dof_id_type> &,
2490 : bool) const {}
2491 :
2492 : inline void DofMap::constrain_element_vector (DenseVector<Number> &,
2493 : std::vector<dof_id_type> &,
2494 : bool) const {}
2495 :
2496 : inline void DofMap::constrain_element_matrix_and_vector (DenseMatrix<Number> &,
2497 : DenseVector<Number> &,
2498 : std::vector<dof_id_type> &,
2499 : bool) const {}
2500 :
2501 : inline void DofMap::heterogeneously_constrain_element_matrix_and_vector
2502 : (DenseMatrix<Number> &, DenseVector<Number> &,
2503 : std::vector<dof_id_type> &, bool, int) const {}
2504 :
2505 : inline void DofMap::heterogeneously_constrain_element_vector
2506 : (const DenseMatrix<Number> &, DenseVector<Number> &,
2507 : std::vector<dof_id_type> &, bool, int) const {}
2508 :
2509 : inline void DofMap::constrain_element_dyad_matrix (DenseVector<Number> &,
2510 : DenseVector<Number> &,
2511 : std::vector<dof_id_type> &,
2512 : bool) const {}
2513 :
2514 : inline void DofMap::constrain_nothing (std::vector<dof_id_type> &) const {}
2515 :
2516 : inline void DofMap::enforce_constraints_exactly (const System &,
2517 : NumericVector<Number> *,
2518 : bool) const {}
2519 :
2520 : inline void DofMap::enforce_adjoint_constraints_exactly (NumericVector<Number> &,
2521 : unsigned int) const {}
2522 :
2523 :
2524 : inline void DofMap::enforce_constraints_on_residual
2525 : (const NonlinearImplicitSystem &,
2526 : NumericVector<Number> *,
2527 : NumericVector<Number> const *,
2528 : bool) const {}
2529 :
2530 : inline void DofMap::enforce_constraints_on_jacobian
2531 : (const NonlinearImplicitSystem &,
2532 : SparseMatrix<Number> *) const {}
2533 :
2534 : #endif // LIBMESH_ENABLE_CONSTRAINTS
2535 :
2536 :
2537 :
2538 : inline
2539 : void DofMap::set_constrained_sparsity_construction(bool use_constraints)
2540 : {
2541 : // This got only partly finished...
2542 : if (use_constraints)
2543 : libmesh_not_implemented();
2544 :
2545 : #ifdef LIBMESH_ENABLE_CONSTRAINTS
2546 : _constrained_sparsity_construction = use_constraints;
2547 : #endif
2548 : libmesh_ignore(use_constraints);
2549 : }
2550 :
2551 : inline
2552 : void DofMap::full_sparsity_pattern_needed()
2553 : {
2554 : need_full_sparsity_pattern = true;
2555 : }
2556 :
2557 : inline
2558 : bool DofMap::constrained_sparsity_construction()
2559 : {
2560 : #ifdef LIBMESH_ENABLE_CONSTRAINTS
2561 : return _constrained_sparsity_construction;
2562 : #else
2563 : return true;
2564 : #endif
2565 : }
2566 :
2567 : inline
2568 : void DofMap::should_p_refine(const unsigned int g, const bool p_refine)
2569 : {
2570 : #ifdef LIBMESH_ENABLE_AMR
2571 : if (p_refine)
2572 : {
2573 : auto it = _dont_p_refine.find(g);
2574 : if (it != _dont_p_refine.end())
2575 : _dont_p_refine.erase(it);
2576 : }
2577 : else
2578 : _dont_p_refine.insert(g);
2579 : #else
2580 : libmesh_ignore(g, p_refine);
2581 : #endif
2582 : }
2583 :
2584 : inline
2585 11654231 : bool DofMap::should_p_refine(const unsigned int g) const
2586 : {
2587 : #ifdef LIBMESH_ENABLE_AMR
2588 12503580 : return !_dont_p_refine.count(g);
2589 : #else
2590 : libmesh_ignore(g);
2591 : return false;
2592 : #endif
2593 : }
2594 :
2595 : inline
2596 1181490 : unsigned int DofMap::var_group_from_var_number(const unsigned int var_num) const
2597 : {
2598 1181490 : libmesh_assert(var_num < n_variables());
2599 20716490 : return libmesh_map_find(_var_to_vg, var_num);
2600 : }
2601 :
2602 : inline
2603 6210322 : bool DofMap::should_p_refine_var(const unsigned int var) const
2604 : {
2605 : #ifdef LIBMESH_ENABLE_AMR
2606 6210322 : const auto vg = this->var_group_from_var_number(var);
2607 6210322 : return !_dont_p_refine.count(vg);
2608 : #else
2609 : libmesh_ignore(var);
2610 : return false;
2611 : #endif
2612 : }
2613 :
2614 : template <typename FieldDofsFunctor>
2615 272260243 : void DofMap::_dof_indices (const Elem & elem,
2616 : int p_level,
2617 : std::vector<dof_id_type> & di,
2618 : const unsigned int vg,
2619 : const unsigned int vig,
2620 : const Node * const * nodes,
2621 : unsigned int n_nodes,
2622 : const unsigned int v,
2623 : #ifdef DEBUG
2624 : std::size_t & tot_size,
2625 : #endif
2626 : FieldDofsFunctor field_dofs_functor) const
2627 : {
2628 272260243 : const VariableGroup & var = this->variable_group(vg);
2629 :
2630 272260243 : if (var.active_on_subdomain(elem.subdomain_id()))
2631 : {
2632 272121978 : const ElemType type = elem.type();
2633 24123966 : const unsigned int sys_num = this->sys_number();
2634 : #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2635 32082856 : const bool is_inf = elem.infinite();
2636 : #endif
2637 :
2638 : const bool extra_hanging_dofs =
2639 272121978 : FEInterface::extra_hanging_dofs(var.type());
2640 :
2641 272121978 : FEType fe_type = var.type();
2642 :
2643 24123966 : const bool add_p_level =
2644 : #ifdef LIBMESH_ENABLE_AMR
2645 24122698 : !_dont_p_refine.count(vg);
2646 : #else
2647 : false;
2648 : #endif
2649 :
2650 : #ifdef DEBUG
2651 : // The number of dofs per element is non-static for subdivision FE
2652 24122698 : if (var.type().family == SUBDIVISION)
2653 3936 : tot_size += n_nodes;
2654 : else
2655 : // FIXME: Is the passed-in p_level just elem.p_level()? If so,
2656 : // this seems redundant.
2657 24118762 : tot_size += FEInterface::n_dofs(fe_type, add_p_level*p_level, &elem);
2658 : #endif
2659 :
2660 : // The total Order is not required when getting the function
2661 : // pointer, it is only needed when the function is called (see
2662 : // below).
2663 : const FEInterface::n_dofs_at_node_ptr ndan =
2664 272121978 : FEInterface::n_dofs_at_node_function(fe_type, &elem);
2665 :
2666 : // Get the node-based DOF numbers
2667 1635777873 : for (unsigned int n=0; n != n_nodes; n++)
2668 : {
2669 1363655895 : const Node & node = *nodes[n];
2670 :
2671 : // Cache the intermediate lookups that are common to every
2672 : // component
2673 : #ifdef DEBUG
2674 : const std::pair<unsigned int, unsigned int>
2675 124010137 : vg_and_offset = node.var_to_vg_and_offset(sys_num,v);
2676 124010137 : libmesh_assert_equal_to (vg, vg_and_offset.first);
2677 124010137 : libmesh_assert_equal_to (vig, vg_and_offset.second);
2678 : #endif
2679 1363655895 : const unsigned int n_comp = node.n_comp_group(sys_num,vg);
2680 :
2681 : // There is a potential problem with h refinement. Imagine a
2682 : // quad9 that has a linear FE on it. Then, on the hanging side,
2683 : // it can falsely identify a DOF at the mid-edge node. This is why
2684 : // we go through FEInterface instead of node.n_comp() directly.
2685 1363655895 : const unsigned int nc =
2686 : #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2687 63517710 : is_inf ?
2688 82442 : FEInterface::n_dofs_at_node(fe_type, add_p_level*p_level, &elem, n) :
2689 : #endif
2690 1363582177 : ndan (type, fe_type.order + add_p_level*p_level, n);
2691 :
2692 : // If this is a non-vertex on a hanging node with extra
2693 : // degrees of freedom, we use the non-vertex dofs (which
2694 : // come in reverse order starting from the end, to
2695 : // simplify p refinement)
2696 1363655895 : if (extra_hanging_dofs && !elem.is_vertex(n))
2697 : {
2698 24869929 : const int dof_offset = n_comp - nc;
2699 :
2700 : // We should never have fewer dofs than necessary on a
2701 : // node unless we're getting indices on a parent element,
2702 : // and we should never need the indices on such a node
2703 24869929 : if (dof_offset < 0)
2704 : {
2705 0 : libmesh_assert(!elem.active());
2706 0 : di.resize(di.size() + nc, DofObject::invalid_id);
2707 : }
2708 : else
2709 56899517 : for (int i=int(n_comp)-1; i>=dof_offset; i--)
2710 : {
2711 2677638 : const dof_id_type d =
2712 32029588 : node.dof_number(sys_num, vg, vig, i, n_comp);
2713 2677638 : libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2714 32029588 : field_dofs_functor(elem, n, v, di, d);
2715 : }
2716 : }
2717 : // If this is a vertex or an element without extra hanging
2718 : // dofs, our dofs come in forward order coming from the
2719 : // beginning
2720 : else
2721 : {
2722 : // We have a good component index only if it's being
2723 : // used on this FE type (nc) *and* it's available on
2724 : // this DofObject (n_comp).
2725 1338785966 : const unsigned int good_nc = std::min(n_comp, nc);
2726 2595525489 : for (unsigned int i=0; i!=good_nc; ++i)
2727 : {
2728 111892123 : const dof_id_type d =
2729 1256734639 : node.dof_number(sys_num, vg, vig, i, n_comp);
2730 111892123 : libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2731 111892123 : libmesh_assert_less (d, this->n_dofs());
2732 1256739523 : field_dofs_functor(elem, n, v, di, d);
2733 : }
2734 :
2735 : // With fewer good component indices than we need, e.g.
2736 : // due to subdomain expansion, the remaining expected
2737 : // indices are marked invalid.
2738 1338785966 : if (n_comp < nc)
2739 0 : for (unsigned int i=n_comp; i!=nc; ++i)
2740 0 : di.push_back(DofObject::invalid_id);
2741 : }
2742 : }
2743 :
2744 : // If there are any element-based DOF numbers, get them
2745 272121978 : const unsigned int nc = FEInterface::n_dofs_per_elem(fe_type, add_p_level*p_level, &elem);
2746 :
2747 : // We should never have fewer dofs than necessary on an
2748 : // element unless we're getting indices on a parent element
2749 : // (and we should never need those indices) or off-domain for a
2750 : // subdomain-restricted variable (where invalid_id is the
2751 : // correct thing to return)
2752 272121978 : if (nc != 0)
2753 : {
2754 14799003 : const unsigned int n_comp = elem.n_comp_group(sys_num,vg);
2755 14799003 : if (elem.n_systems() > sys_num && nc <= n_comp)
2756 : {
2757 90066464 : for (unsigned int i=0; i<nc; i++)
2758 : {
2759 6779827 : const dof_id_type d =
2760 75264773 : elem.dof_number(sys_num, vg, vig, i, n_comp);
2761 6779827 : libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2762 :
2763 75267461 : field_dofs_functor(elem, invalid_uint, v, di, d);
2764 : }
2765 : }
2766 : else
2767 : {
2768 0 : libmesh_assert(!elem.active() || fe_type.family == LAGRANGE || fe_type.family == SUBDIVISION);
2769 0 : di.resize(di.size() + nc, DofObject::invalid_id);
2770 : }
2771 : }
2772 : }
2773 272260243 : }
2774 :
2775 :
2776 :
2777 : template <typename ScalarDofsFunctor, typename FieldDofsFunctor>
2778 115082704 : void DofMap::dof_indices (const Elem * const elem,
2779 : std::vector<dof_id_type> & di,
2780 : const unsigned int vn,
2781 : ScalarDofsFunctor scalar_dofs_functor,
2782 : FieldDofsFunctor field_dofs_functor,
2783 : int p_level) const
2784 : {
2785 : // We now allow elem==nullptr to request just SCALAR dofs
2786 : // libmesh_assert(elem);
2787 :
2788 11091874 : LOG_SCOPE("dof_indices()", "DofMap");
2789 :
2790 : // Clear the DOF indices vector
2791 11091874 : di.clear();
2792 :
2793 : // Use the default p refinement level?
2794 115082704 : if (p_level == -12345)
2795 114457636 : p_level = elem ? elem->p_level() : 0;
2796 :
2797 115084004 : const unsigned int vg = this->_variable_group_numbers[vn];
2798 11091874 : const VariableGroup & var = this->variable_group(vg);
2799 115082704 : const unsigned int vig = vn - var.number();
2800 :
2801 : #ifdef DEBUG
2802 : // Check that sizes match in DEBUG mode
2803 11091874 : std::size_t tot_size = 0;
2804 : #endif
2805 :
2806 115082704 : if (elem && elem->type() == TRI3SUBDIVISION)
2807 : {
2808 : // Subdivision surface FE require the 1-ring around elem
2809 2712 : const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *>(elem);
2810 :
2811 : // Ghost subdivision elements have no real dofs
2812 29832 : if (!sd_elem->is_ghost())
2813 : {
2814 : // Determine the nodes contributing to element elem
2815 4608 : std::vector<const Node *> elem_nodes;
2816 25344 : MeshTools::Subdivision::find_one_ring(sd_elem, elem_nodes);
2817 :
2818 25344 : _dof_indices(*elem, p_level, di, vg, vig, elem_nodes.data(),
2819 : cast_int<unsigned int>(elem_nodes.size()), vn,
2820 : #ifdef DEBUG
2821 : tot_size,
2822 : #endif
2823 : field_dofs_functor);
2824 : }
2825 :
2826 29832 : return;
2827 : }
2828 :
2829 : // Get the dof numbers
2830 115150823 : if (var.type().family == SCALAR &&
2831 1027177 : (!elem ||
2832 1027209 : var.active_on_subdomain(elem->subdomain_id())))
2833 : {
2834 : #ifdef DEBUG
2835 97951 : tot_size += var.type().order;
2836 : #endif
2837 195902 : std::vector<dof_id_type> di_new;
2838 1027209 : this->SCALAR_dof_indices(di_new,vn);
2839 1027209 : scalar_dofs_functor(*elem, di, di_new);
2840 : }
2841 114025663 : else if (elem)
2842 114025663 : _dof_indices(*elem, p_level, di, vg, vig, elem->get_nodes(),
2843 114025663 : elem->n_nodes(), vn,
2844 : #ifdef DEBUG
2845 : tot_size,
2846 : #endif
2847 : field_dofs_functor);
2848 :
2849 : #ifdef DEBUG
2850 11089162 : libmesh_assert_equal_to (tot_size, di.size());
2851 : #endif
2852 : }
2853 :
2854 : inline
2855 10 : StaticCondensationDofMap & DofMap::get_static_condensation()
2856 : {
2857 10 : libmesh_assert(_sc);
2858 10 : return *_sc;
2859 : }
2860 :
2861 : inline const std::pair<unsigned int, unsigned int> &
2862 104 : DofMap::get_variable_array(const unsigned int vi) const
2863 : {
2864 88 : auto it = std::upper_bound(
2865 : _array_variables.begin(),
2866 : _array_variables.end(),
2867 : vi,
2868 136 : [](unsigned int value, const std::pair<unsigned int, unsigned int> & b) { return value < b.first; });
2869 :
2870 16 : libmesh_assert_msg(it != _array_variables.begin(),
2871 : "Passed in " << std::to_string(vi) << " is not in any of our array variables");
2872 16 : --it;
2873 16 : libmesh_assert_msg(vi < it->second,
2874 : "Passed in " << std::to_string(vi) << " is not in any of our array variables");
2875 120 : return *it;
2876 : }
2877 :
2878 : template <typename DofIndicesFunctor>
2879 120 : void DofMap::array_dof_indices(const DofIndicesFunctor & functor,
2880 : std::vector<dof_id_type> & di,
2881 : const unsigned int vn) const
2882 : {
2883 136 : const auto [begin, end] = this->get_variable_array(vn);
2884 104 : functor(di, begin);
2885 :
2886 120 : const unsigned int count = end - begin;
2887 : // We make count, which could be >> ntest, the inner index in hopes of vectorization
2888 120 : if (count > 1)
2889 : {
2890 40 : const dof_id_type component_size = di.size();
2891 120 : di.resize(count * component_size);
2892 :
2893 1512 : const auto pack_container = [&di,
2894 : component_size](const unsigned int j,
2895 : const std::vector<dof_id_type> & j_dof_indices,
2896 32 : const unsigned int stride) {
2897 32 : if (&j_dof_indices != &di)
2898 4 : libmesh_assert(j_dof_indices.size() == component_size);
2899 1568 : for (const auto i : make_range(component_size))
2900 1536 : di[j * component_size + i] = j_dof_indices[i] + stride * j;
2901 : };
2902 104 : pack_container(0, di, 0);
2903 :
2904 120 : const auto & fe_type = _variable_groups[libmesh_map_find(_var_to_vg, vn)].type();
2905 136 : if (const bool lagrange = fe_type.family == LAGRANGE;
2906 120 : lagrange || (FEInterface::get_continuity(fe_type) == DISCONTINUOUS))
2907 : {
2908 90 : const auto stride = lagrange ? 1 : component_size;
2909 180 : for (const auto j : make_range((unsigned int)1, count))
2910 78 : pack_container(j, di, stride);
2911 : }
2912 : else
2913 : {
2914 30 : static thread_local std::vector<dof_id_type> work_dof_indices;
2915 4 : unsigned int j = 1;
2916 60 : for (const auto i : make_range(begin + 1, end))
2917 : {
2918 26 : functor(work_dof_indices, i);
2919 30 : pack_container(j++, work_dof_indices, 0);
2920 : }
2921 : }
2922 : }
2923 120 : }
2924 :
2925 : inline
2926 10536551 : unsigned int DofMap::n_vars() const
2927 : {
2928 21056184 : return cast_int<unsigned int>(_variables.size());
2929 : }
2930 :
2931 : inline
2932 2503548 : const std::string & DofMap::variable_name (const unsigned int i) const
2933 : {
2934 2503548 : libmesh_assert_less (i, _variables.size());
2935 :
2936 30524975 : return _variables[i].name();
2937 : }
2938 :
2939 : inline
2940 1184 : bool DofMap::identify_variable_groups () const
2941 : {
2942 40981 : return _identify_variable_groups;
2943 : }
2944 :
2945 : inline
2946 0 : void DofMap::identify_variable_groups (const bool ivg)
2947 : {
2948 0 : _identify_variable_groups = ivg;
2949 0 : }
2950 :
2951 : inline
2952 241793 : unsigned int DofMap::n_components(const MeshBase & mesh) const
2953 : {
2954 248893 : if (_variables.empty())
2955 6786 : return 0;
2956 :
2957 314 : const Variable & last = _variables.back();
2958 11103 : return last.first_scalar_number() + last.n_components(mesh);
2959 : }
2960 :
2961 : inline
2962 : unsigned int
2963 497554 : DofMap::variable_scalar_number (unsigned int var_num,
2964 : unsigned int component) const
2965 : {
2966 6609760 : return _variables[var_num].first_scalar_number() + component;
2967 : }
2968 :
2969 : inline
2970 34491 : const FEType & DofMap::variable_type (std::string_view var) const
2971 : {
2972 35463 : return _variables[this->variable_number(var)].type();
2973 : }
2974 :
2975 960 : inline bool DofMap::has_variable(std::string_view var) const
2976 : {
2977 960 : return _variable_numbers.count(var);
2978 : }
2979 :
2980 385695 : inline unsigned int DofMap::variable_number(std::string_view var) const
2981 : {
2982 8865541 : auto var_num = libmesh_map_find(_variable_numbers, var);
2983 385695 : libmesh_assert_equal_to(_variables[var_num].name(), var);
2984 385695 : return var_num;
2985 : }
2986 :
2987 : } // namespace libMesh
2988 :
2989 : #endif // LIBMESH_DOF_MAP_H
|