libMesh
dof_map.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #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
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 
98 typedef std::map<dof_id_type, Real,
99  std::less<dof_id_type>,
101 
108 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 
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 
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 
146 typedef std::map<const Node *, Real,
147  std::less<const Node *>,
149 
156 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 
179 class DofMap : public DofMapBase,
180  public ReferenceCountedObject<DofMap>
181 {
182 public:
183 
189  explicit
190  DofMap(const unsigned int sys_number,
191  MeshBase & mesh);
192 
196  ~DofMap();
197 
202  {};
203 
209  {
210  public:
211  virtual ~AugmentSendList () = default;
212 
216  virtual void augment_send_list (std::vector<dof_id_type> & send_list) = 0;
217  };
218 
224  void attach_matrix (SparseMatrix<Number> & matrix);
225 
232  void update_sparsity_pattern(SparseMatrix<Number> & matrix) const;
233 
238  bool is_attached (SparseMatrix<Number> & matrix);
239 
246  std::size_t distribute_dofs (MeshBase &);
247 
253  void compute_sparsity (const MeshBase &);
254 
258  bool computed_sparsity_already () const;
259 
270  void set_constrained_sparsity_construction(bool use_constraints);
271 
276 
283 
287  void clear_sparsity();
288 
303 
308  void add_default_ghosting();
309 
315  typedef std::vector<GhostingFunctor *>::const_iterator GhostingFunctorIterator;
316 
341  void add_coupling_functor(GhostingFunctor & coupling_functor,
342  bool to_mesh = true);
343 
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 
361  void remove_coupling_functor(GhostingFunctor & coupling_functor);
362 
367  { return _coupling_functors.begin(); }
368 
373  { return _coupling_functors.end(); }
374 
379 
403  void add_algebraic_ghosting_functor(GhostingFunctor & evaluable_functor,
404  bool to_mesh = true);
405 
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 
423  void remove_algebraic_ghosting_functor(GhostingFunctor & evaluable_functor);
424 
429  { return _algebraic_ghosting_functors.begin(); }
430 
435  { return _algebraic_ghosting_functors.end(); }
436 
441 
453  {
455  }
456 
468  std::vector<dof_id_type> & n_nz,
469  std::vector<dof_id_type> & n_oz,
470  void *),
471  void * context = nullptr)
473 
482  {
483  _augment_send_list = &asl;
484  }
485 
490  void attach_extra_send_list_function(void (*func)(std::vector<dof_id_type> &, void *),
491  void * context = nullptr)
493 
500  void prepare_send_list ();
501 
508  {
509  _send_list.clear();
510  }
511 
521  void reinit_send_list (MeshBase & mesh);
522 
523 
533  const std::vector<dof_id_type> & get_send_list() const { return _send_list; }
534 
542  const std::vector<dof_id_type> & get_n_nz() const
543  {
545  return _sp->get_n_nz();
546  }
547 
555  const std::vector<dof_id_type> & get_n_oz() const
556  {
558  return _sp->get_n_oz();
559  }
560 
561 
571  {
572  return _sp.get();
573  }
574 
575  // /**
576  // * Add an unknown of order \p order and finite element type
577  // * \p type to the system of equations.
578  // */
579  // void add_variable (const Variable & var);
580 
585  void add_variable_group (VariableGroup var_group);
586 
597  void set_error_on_cyclic_constraint(bool error_on_cyclic_constraint);
598  void set_error_on_constraint_loop(bool error_on_constraint_loop);
599 
603  const VariableGroup & variable_group (const unsigned int c) const;
604 
605  const Variable & variable (const unsigned int c) const override;
606 
610  Order variable_order (const unsigned int c) const;
611 
615  Order variable_group_order (const unsigned int vg) const;
616 
620  const FEType & variable_type (const unsigned int c) const;
621 
625  const FEType & variable_group_type (const unsigned int vg) const;
626 
632  unsigned int n_variable_groups() const
633  { return cast_int<unsigned int>(_variable_groups.size()); }
634 
635  unsigned int n_variables() const override
636  { return cast_int<unsigned int>(_variables.size()); }
637 
641  unsigned int var_group_from_var_number(unsigned int var_num) const;
642 
649  {
650  return ((this->n_variable_groups() == 1) && (this->n_variables() > 1));
651  }
652 
665  unsigned int block_size() const
666  {
667  return (this->has_blocked_representation() ? this->n_variables() : 1);
668  }
669 
670  using DofMapBase::n_dofs;
675  dof_id_type n_dofs(const unsigned int vn) const
676  {
677  dof_id_type n = this->n_local_dofs(vn);
678  this->comm().sum(n);
679  return n;
680  }
681 
686 
693  dof_id_type n_local_dofs(const unsigned int vn) const
694  {
695  dof_id_type n;
696  this->local_variable_indices(n, _mesh, vn);
697  return n;
698  }
699 
704  std::vector<dof_id_type> n_dofs_per_processor(const unsigned int vn) const
705  {
706  std::vector<dof_id_type> n_local_dofs(this->n_processors(), 0);
707  this->comm().allgather(this->n_local_dofs(vn), n_local_dofs);
708  return n_local_dofs;
709  }
710 
715  { std::vector<dof_id_type>::const_iterator ub =
716  std::upper_bound(_end_df.begin(), _end_df.end(), dof);
717  libmesh_assert (ub != _end_df.end());
718  return cast_int<processor_id_type>(ub - _end_df.begin());
719  }
720 
721  void dof_indices (const Elem * const elem,
722  std::vector<dof_id_type> & di) const;
723 
729  void dof_indices (const Elem * const elem,
730  std::vector<dof_id_type> & di,
731  const unsigned int vn,
732  int p_level = -12345) const override;
733 
763  template <typename ScalarDofsFunctor, typename FieldDofsFunctor>
764  void dof_indices(const Elem * const elem,
765  std::vector<dof_id_type> & di,
766  const unsigned int vn,
767  ScalarDofsFunctor scalar_dofs_functor,
768  FieldDofsFunctor field_dofs_functor,
769  int p_level = -12345) const;
770 
775  void dof_indices (const Node * const node,
776  std::vector<dof_id_type> & di) const;
777 
782  void dof_indices (const Node * const node,
783  std::vector<dof_id_type> & di,
784  const unsigned int vn) const override;
785 
792  void dof_indices (const Elem & elem,
793  unsigned int n,
794  std::vector<dof_id_type> & di,
795  const unsigned int vn) const;
796 
797 #ifdef LIBMESH_ENABLE_AMR
798 
805  void old_dof_indices (const Elem & elem,
806  unsigned int n,
807  std::vector<dof_id_type> & di,
808  const unsigned int vn) const;
809 
810 #endif // LIBMESH_ENABLE_AMR
811 
820  void SCALAR_dof_indices (std::vector<dof_id_type> & di,
821  const unsigned int vn,
822  const bool old_dofs=false) const;
823 
831  bool semilocal_index (dof_id_type dof_index) const;
832 
840  bool all_semilocal_indices (const std::vector<dof_id_type> & dof_indices) const;
841 
846  bool local_index (dof_id_type dof_index) const
847  { return (dof_index >= this->first_dof()) && (dof_index < this->end_dof()); }
848 
854  template <typename DofObjectSubclass>
855  bool is_evaluable(const DofObjectSubclass & obj,
856  unsigned int var_num = libMesh::invalid_uint) const;
857 
865  void set_implicit_neighbor_dofs(bool implicit_neighbor_dofs);
866 
871 
878  bool use_coupled_neighbor_dofs(const MeshBase & mesh) const;
879 
891  const std::vector<dof_id_type> & dof_indices,
892  DenseVectorBase<Number> & Ue) const;
893 
899  template <typename T, std::enable_if_t<std::is_same_v<T, dof_id_type> ||
900  std::is_same_v<T, std::vector<dof_id_type>>, int> = 0>
901  void local_variable_indices(T & idx,
902  const MeshBase & mesh,
903  unsigned int var_num) const;
904 
910  template <typename T,
911  std::enable_if_t<std::is_same_v<T, dof_id_type> ||
912  std::is_same_v<T, std::vector<dof_id_type>>,
913  int> = 0>
914  void local_variable_indices(T & idx, unsigned int var_num) const
915  { this->local_variable_indices(idx, this->_mesh, var_num); }
916 
917 #ifdef LIBMESH_ENABLE_CONSTRAINTS
918 
919  //--------------------------------------------------------------------
920  // Constraint-specific methods
926 
932 
933 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
934 
939  { return cast_int<dof_id_type>(_node_constraints.size()); }
940 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
941 
950  void create_dof_constraints (const MeshBase &, Real time=0);
951 
956 
960  void scatter_constraints (MeshBase &);
961 
976  std::set<dof_id_type> & unexpanded_dofs,
977  bool look_for_constrainees);
978 
986  void process_constraints (MeshBase &);
987 
1004 
1009  void add_constraint_row (const dof_id_type dof_number,
1010  const DofConstraintRow & constraint_row,
1011  const Number constraint_rhs,
1012  const bool forbid_constraint_overwrite);
1013 
1024  void add_adjoint_constraint_row (const unsigned int qoi_index,
1025  const dof_id_type dof_number,
1026  const DofConstraintRow & constraint_row,
1027  const Number constraint_rhs,
1028  const bool forbid_constraint_overwrite);
1029 
1035  void add_constraint_row (const dof_id_type dof_number,
1036  const DofConstraintRow & constraint_row,
1037  const bool forbid_constraint_overwrite = true)
1038  { add_constraint_row(dof_number, constraint_row, 0., forbid_constraint_overwrite); }
1039 
1043  DofConstraints::const_iterator constraint_rows_begin() const
1044  { return _dof_constraints.begin(); }
1045 
1049  DofConstraints::const_iterator constraint_rows_end() const
1050  { return _dof_constraints.end(); }
1051 
1057 
1059  {
1062  }
1063 
1065  {
1068  }
1069 
1083  {
1085  }
1086 
1087 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
1088 
1091  NodeConstraints::const_iterator node_constraint_rows_begin() const
1092  { return _node_constraints.begin(); }
1093 
1097  NodeConstraints::const_iterator node_constraint_rows_end() const
1098  { return _node_constraints.end(); }
1099 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
1100 
1105  bool is_constrained_dof (const dof_id_type dof) const;
1106 
1111  bool has_heterogeneous_adjoint_constraints (const unsigned int qoi_num) const;
1112 
1116  bool has_heterogenous_adjoint_constraints (const unsigned int qoi_num) const
1117  {
1118  return this->has_heterogeneous_adjoint_constraints (qoi_num);
1119  }
1120 
1126  Number has_heterogeneous_adjoint_constraint (const unsigned int qoi_num,
1127  const dof_id_type dof) const;
1128 
1132  Number has_heterogenous_adjoint_constraint (const unsigned int qoi_num,
1133  const dof_id_type dof) const
1134  {
1135  return this->has_heterogeneous_adjoint_constraint (qoi_num, dof);
1136  }
1137 
1143 
1148  bool is_constrained_node (const Node * node) const;
1149 
1156  void print_dof_constraints(std::ostream & os=libMesh::out,
1157  bool print_nonlocal=false) const;
1158 
1164  std::string get_local_constraints(bool print_nonlocal=false) const;
1165 
1166 
1175  std::pair<Real, Real> max_constraint_error(const System & system,
1176  NumericVector<Number> * v = nullptr) const;
1177 
1178 #endif // LIBMESH_ENABLE_CONSTRAINTS
1179 
1180  //--------------------------------------------------------------------
1181  // Constraint-specific methods
1182  // Some of these methods are enabled (but inlined away to nothing)
1183  // when constraints are disabled at configure-time. This is to
1184  // increase API compatibility of user code with different library
1185  // builds.
1186 
1201  std::vector<dof_id_type> & elem_dofs,
1202  bool asymmetric_constraint_rows = true) const;
1203 
1211  std::vector<dof_id_type> & row_dofs,
1212  std::vector<dof_id_type> & col_dofs,
1213  bool asymmetric_constraint_rows = true) const;
1214 
1219  std::vector<dof_id_type> & dofs,
1220  bool asymmetric_constraint_rows = true) const;
1221 
1231  DenseVector<Number> & rhs,
1232  std::vector<dof_id_type> & elem_dofs,
1233  bool asymmetric_constraint_rows = true) const;
1234 
1259  DenseVector<Number> & rhs,
1260  std::vector<dof_id_type> & elem_dofs,
1261  bool asymmetric_constraint_rows = true,
1262  int qoi_index = -1) const;
1263 
1264  /*
1265  * Backwards compatibility with misspelling.
1266  */
1268  DenseVector<Number> & rhs,
1269  std::vector<dof_id_type> & elem_dofs,
1270  bool asymmetric_constraint_rows = true,
1271  int qoi_index = -1) const
1272  {
1274  (matrix, rhs, elem_dofs, asymmetric_constraint_rows, qoi_index);
1275  }
1276 
1303  DenseVector<Number> & rhs,
1304  std::vector<dof_id_type> & elem_dofs,
1305  bool asymmetric_constraint_rows = true,
1306  int qoi_index = -1) const;
1307 
1308  /*
1309  * Backwards compatibility with misspelling.
1310  */
1312  DenseVector<Number> & rhs,
1313  std::vector<dof_id_type> & elem_dofs,
1314  bool asymmetric_constraint_rows = true,
1315  int qoi_index = -1) const
1316  {
1318  (matrix, rhs, elem_dofs, asymmetric_constraint_rows, qoi_index);
1319  }
1320 
1341  DenseVector<Number> & rhs,
1342  std::vector<dof_id_type> & elem_dofs,
1343  NumericVector<Number> & solution_local) const;
1344 
1360  std::vector<dof_id_type> & elem_dofs,
1361  NumericVector<Number> & solution_local) const;
1362 
1363 
1380  std::vector<dof_id_type> & elem_dofs,
1381  NumericVector<Number> & solution_local) const;
1382 
1392  DenseVector<Number> & w,
1393  std::vector<dof_id_type> & row_dofs,
1394  bool asymmetric_constraint_rows = true) const;
1395 
1402  void constrain_nothing (std::vector<dof_id_type> & dofs) const;
1403 
1417  void enforce_constraints_exactly (const System & system,
1418  NumericVector<Number> * v = nullptr,
1419  bool homogeneous = false) const;
1420 
1428  unsigned int q) const;
1429 
1431  NumericVector<Number> * rhs,
1432  NumericVector<Number> const * solution,
1433  bool homogeneous = true) const;
1434 
1436  SparseMatrix<Number> * jac) const;
1437 
1438 #ifdef LIBMESH_ENABLE_PERIODIC
1439 
1440  //--------------------------------------------------------------------
1441  // PeriodicBoundary-specific methods
1442 
1446  void add_periodic_boundary (const PeriodicBoundaryBase & periodic_boundary);
1447 
1454  void add_periodic_boundary (const PeriodicBoundaryBase & boundary, const PeriodicBoundaryBase & inverse_boundary);
1455 
1460  bool is_periodic_boundary (const boundary_id_type boundaryid) const;
1461 
1463  {
1464  return _periodic_boundaries.get();
1465  }
1466 
1468  {
1469  return _periodic_boundaries.get();
1470  }
1471 
1472 #endif // LIBMESH_ENABLE_PERIODIC
1473 
1474 
1475 #ifdef LIBMESH_ENABLE_DIRICHLET
1476 
1477  //--------------------------------------------------------------------
1478  // DirichletBoundary-specific methods
1479 
1493  void add_dirichlet_boundary (const DirichletBoundary & dirichlet_boundary);
1494 
1500  void add_adjoint_dirichlet_boundary (const DirichletBoundary & dirichlet_boundary,
1501  unsigned int q);
1502 
1506  void remove_dirichlet_boundary (const DirichletBoundary & dirichlet_boundary);
1507 
1512  void remove_adjoint_dirichlet_boundary (const DirichletBoundary & dirichlet_boundary,
1513  unsigned int q);
1514 
1516  {
1517  return _dirichlet_boundaries.get();
1518  }
1519 
1521  {
1522  return _dirichlet_boundaries.get();
1523  }
1524 
1525  bool has_adjoint_dirichlet_boundaries(unsigned int q) const;
1526 
1527  const DirichletBoundaries *
1528  get_adjoint_dirichlet_boundaries(unsigned int q) const;
1529 
1531  get_adjoint_dirichlet_boundaries(unsigned int q);
1532 
1538  const DirichletBoundary & boundary) const;
1539 #endif // LIBMESH_ENABLE_DIRICHLET
1540 
1541 
1542 #ifdef LIBMESH_ENABLE_AMR
1543 
1544  //--------------------------------------------------------------------
1545  // AMR-specific methods
1546 
1555  // void augment_send_list_for_projection(const MeshBase &);
1556 
1557 #ifdef LIBMESH_ENABLE_AMR
1558 
1565  void old_dof_indices (const Elem * const elem,
1566  std::vector<dof_id_type> & di,
1567  const unsigned int vn = libMesh::invalid_uint) const;
1568 
1569 #endif // LIBMESH_ENABLE_AMR
1570 
1576  void constrain_p_dofs (unsigned int var,
1577  const Elem * elem,
1578  unsigned int s,
1579  unsigned int p);
1580 
1581 #endif // LIBMESH_ENABLE_AMR
1582 
1586  void reinit
1587  (MeshBase & mesh,
1588  const std::map<const Node *, std::set<subdomain_id_type>> &
1589  constraining_subdomains);
1590 
1595  virtual void clear () override;
1596 
1600  void print_info(std::ostream & os=libMesh::out) const;
1601 
1605  std::string get_info() const;
1606 
1621 
1625  unsigned int sys_number() const;
1626 
1641  std::unique_ptr<SparsityPattern::Build> build_sparsity(const MeshBase & mesh,
1642  bool calculate_constrained = false,
1643  bool use_condensed_system = false) const;
1644 
1649  void should_p_refine(unsigned int g, bool p_refine);
1650 
1654  bool should_p_refine(unsigned int g) const;
1655 
1659  bool should_p_refine_var(unsigned int var) const;
1660 
1661  // Prevent bad user implicit conversions
1662  void should_p_refine(FEFamily, bool) = delete;
1663  void should_p_refine(Order, bool) = delete;
1664  bool should_p_refine(FEFamily) const = delete;
1665  bool should_p_refine(Order) const = delete;
1666 
1670  void create_static_condensation(MeshBase & mesh, System & system);
1671 
1675  bool has_static_condensation() const { return _sc.get(); }
1676 
1682 
1687 
1688 private:
1689 
1700  void _dof_indices (const Elem & elem,
1701  int p_level,
1702  std::vector<dof_id_type> & di,
1703  const unsigned int vg,
1704  const unsigned int vig,
1705  const Node * const * nodes,
1706  unsigned int n_nodes,
1707  const unsigned int v
1708 #ifdef DEBUG
1709  ,
1710  std::size_t & tot_size
1711 #endif
1712  ) const;
1713 
1729  template <typename FieldDofsFunctor>
1730  void _dof_indices (const Elem & elem,
1731  int p_level,
1732  std::vector<dof_id_type> & di,
1733  const unsigned int vg,
1734  const unsigned int vig,
1735  const Node * const * nodes,
1736  unsigned int n_nodes,
1737  const unsigned int v,
1738 #ifdef DEBUG
1739  std::size_t & tot_size,
1740 #endif
1741  FieldDofsFunctor field_dofs_functor) const;
1742 
1747  void _node_dof_indices (const Elem & elem,
1748  unsigned int n,
1749  const DofObject & obj,
1750  std::vector<dof_id_type> & di,
1751  const unsigned int vn) const;
1752 
1756  void invalidate_dofs(MeshBase & mesh) const;
1757 
1761  DofObject * node_ptr(MeshBase & mesh, dof_id_type i) const;
1762 
1766  DofObject * elem_ptr(MeshBase & mesh, dof_id_type i) const;
1767 
1771  typedef DofObject * (DofMap::*dofobject_accessor)
1772  (MeshBase & mesh, dof_id_type i) const;
1773 
1777  template<typename iterator_type>
1778  void set_nonlocal_dof_objects(iterator_type objects_begin,
1779  iterator_type objects_end,
1780  MeshBase & mesh,
1781  dofobject_accessor objects);
1782 
1797  std::map<const Node *, std::set<subdomain_id_type>>
1799 
1813  (dof_id_type & next_free_dof,
1814  MeshBase & mesh,
1815  const std::map<const Node *, std::set<subdomain_id_type>> &
1816  constraining_subdomains);
1817 
1834  (dof_id_type & next_free_dof,
1835  MeshBase & mesh,
1836  const std::map<const Node *, std::set<subdomain_id_type>> &
1837  constraining_subdomains);
1838 
1839  /*
1840  * Helper method for the above two to count + distriubte SCALAR dofs
1841  */
1842  void distribute_scalar_dofs (dof_id_type & next_free_dof);
1843 
1844 #ifdef DEBUG
1845  /*
1846  * Internal assertions for distribute_local_dofs_*
1847  */
1849 #endif
1850 
1851  /*
1852  * A utility method for obtaining a set of elements to ghost along
1853  * with merged coupling matrices.
1854  */
1855  typedef std::set<std::unique_ptr<CouplingMatrix>, Utility::CompareUnderlying> CouplingMatricesSet;
1856  static void
1858  CouplingMatricesSet & temporary_coupling_matrices,
1859  const GhostingFunctorIterator & gf_begin,
1860  const GhostingFunctorIterator & gf_end,
1861  const MeshBase::const_element_iterator & elems_begin,
1862  const MeshBase::const_element_iterator & elems_end,
1863  processor_id_type p);
1864 
1870 
1871 #ifdef LIBMESH_ENABLE_CONSTRAINTS
1872 
1884  std::vector<dof_id_type> & elem_dofs,
1885  const bool called_recursively=false) const;
1886 
1903  DenseVector<Number> & H,
1904  std::vector<dof_id_type> & elem_dofs,
1905  int qoi_index = -1,
1906  const bool called_recursively=false) const;
1907 
1912  void find_connected_dofs (std::vector<dof_id_type> & elem_dofs) const;
1913 
1918  void find_connected_dof_objects (std::vector<const DofObject *> & objs) const;
1919 
1926 
1933 
1934 #endif // LIBMESH_ENABLE_CONSTRAINTS
1935 
1942 
1948 
1952  std::vector<Variable> _variables;
1953 
1957  std::vector<VariableGroup> _variable_groups;
1958 
1962  std::vector<unsigned int> _variable_group_numbers;
1963 
1967  std::unordered_map<unsigned int, unsigned int> _var_to_vg;
1968 
1972  const unsigned int _sys_number;
1973 
1978 
1984  std::vector<SparseMatrix<Number> * > _matrices;
1985 
1990  std::vector<dof_id_type> _first_scalar_df;
1991 
1996  std::vector<dof_id_type> _send_list;
1997 
2002 
2007  std::vector<dof_id_type> & n_nz,
2008  std::vector<dof_id_type> & n_oz,
2009  void *);
2014 
2019 
2023  void (*_extra_send_list_function)(std::vector<dof_id_type> &, void *);
2024 
2029 
2036  std::unique_ptr<DefaultCoupling> _default_coupling;
2037 
2044  std::unique_ptr<DefaultCoupling> _default_evaluating;
2045 
2057  std::vector<GhostingFunctor *> _algebraic_ghosting_functors;
2058 
2070  std::vector<GhostingFunctor *> _coupling_functors;
2071 
2076  std::map<GhostingFunctor *, std::shared_ptr<GhostingFunctor> > _shared_functors;
2077 
2083 
2089  std::unique_ptr<SparsityPattern::Build> _sp;
2090 
2096 
2097 #ifdef LIBMESH_ENABLE_AMR
2098 
2103  std::vector<dof_id_type> _first_old_scalar_df;
2104 
2108  std::unordered_set<unsigned int> _dont_p_refine;
2109 #endif
2110 
2111 #ifdef LIBMESH_ENABLE_CONSTRAINTS
2112 
2117 
2119 
2121 #endif
2122 
2123 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
2124 
2128 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
2129 
2130 
2131 #ifdef LIBMESH_ENABLE_PERIODIC
2132 
2136  std::unique_ptr<PeriodicBoundaries> _periodic_boundaries;
2137 #endif
2138 
2139 #ifdef LIBMESH_ENABLE_DIRICHLET
2140 
2144  std::unique_ptr<DirichletBoundaries> _dirichlet_boundaries;
2145 
2150  std::vector<std::unique_ptr<DirichletBoundaries>> _adjoint_dirichlet_boundaries;
2151 #endif
2152 
2154 
2161 
2173 
2175  std::unique_ptr<StaticCondensationDofMap> _sc;
2176 };
2177 
2178 
2179 // ------------------------------------------------------------
2180 // Dof Map inline member functions
2181 inline
2182 unsigned int DofMap::sys_number() const
2183 {
2184  return _sys_number;
2185 }
2186 
2187 
2188 
2189 inline
2190 const VariableGroup & DofMap::variable_group (const unsigned int g) const
2191 {
2192  libmesh_assert_less (g, _variable_groups.size());
2193 
2194  return _variable_groups[g];
2195 }
2196 
2197 
2198 
2199 inline
2200 const Variable & DofMap::variable (const unsigned int c) const
2201 {
2202  libmesh_assert_less (c, _variables.size());
2203 
2204  return _variables[c];
2205 }
2206 
2207 
2208 
2209 inline
2210 Order DofMap::variable_order (const unsigned int c) const
2211 {
2212  libmesh_assert_less (c, _variables.size());
2213 
2214  return _variables[c].type().order;
2215 }
2216 
2217 
2218 
2219 inline
2220 Order DofMap::variable_group_order (const unsigned int vg) const
2221 {
2222  libmesh_assert_less (vg, _variable_groups.size());
2223 
2224  return _variable_groups[vg].type().order;
2225 }
2226 
2227 
2228 
2229 inline
2230 const FEType & DofMap::variable_type (const unsigned int c) const
2231 {
2232  libmesh_assert_less (c, _variables.size());
2233 
2234  return _variables[c].type();
2235 }
2236 
2237 
2238 
2239 inline
2240 const FEType & DofMap::variable_group_type (const unsigned int vg) const
2241 {
2242  libmesh_assert_less (vg, _variable_groups.size());
2243 
2244  return _variable_groups[vg].type();
2245 }
2246 
2247 
2248 #ifdef LIBMESH_ENABLE_CONSTRAINTS
2249 
2250 
2251 inline
2253 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
2254  node
2255 #endif
2256  ) const
2257 {
2258 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
2259  if (_node_constraints.count(node))
2260  return true;
2261 #endif
2262 
2263  return false;
2264 }
2265 
2266 
2267 inline
2269 {
2270  if (_dof_constraints.count(dof))
2271  return true;
2272 
2273  return false;
2274 }
2275 
2276 
2277 inline
2278 bool DofMap::has_heterogeneous_adjoint_constraints (const unsigned int qoi_num) const
2279 {
2280  AdjointDofConstraintValues::const_iterator it =
2281  _adjoint_constraint_values.find(qoi_num);
2282  if (it == _adjoint_constraint_values.end())
2283  return false;
2284  if (it->second.empty())
2285  return false;
2286 
2287  return true;
2288 }
2289 
2290 
2291 inline
2293  const dof_id_type dof) const
2294 {
2295  AdjointDofConstraintValues::const_iterator it =
2296  _adjoint_constraint_values.find(qoi_num);
2297  if (it != _adjoint_constraint_values.end())
2298  {
2299  DofConstraintValueMap::const_iterator rhsit =
2300  it->second.find(dof);
2301  if (rhsit == it->second.end())
2302  return 0;
2303  else
2304  return rhsit->second;
2305  }
2306 
2307  return 0;
2308 }
2309 
2310 
2311 
2312 inline
2314 {
2316 }
2317 
2318 
2319 
2320 #else
2321 
2322 //--------------------------------------------------------------------
2323 // Constraint-specific methods get inlined into nothing if
2324 // constraints are disabled, so there's no reason for users not to
2325 // use them.
2326 
2328  std::vector<dof_id_type> &,
2329  bool) const {}
2330 
2332  std::vector<dof_id_type> &,
2333  std::vector<dof_id_type> &,
2334  bool) const {}
2335 
2337  std::vector<dof_id_type> &,
2338  bool) const {}
2339 
2342  std::vector<dof_id_type> &,
2343  bool) const {}
2344 
2347  std::vector<dof_id_type> &, bool, int) const {}
2348 
2351  std::vector<dof_id_type> &, bool, int) const {}
2352 
2355  std::vector<dof_id_type> &,
2356  bool) const {}
2357 
2358 inline void DofMap::constrain_nothing (std::vector<dof_id_type> &) const {}
2359 
2362  bool) const {}
2363 
2365  unsigned int) const {}
2366 
2367 
2371  NumericVector<Number> const *,
2372  bool) const {}
2373 
2376  SparseMatrix<Number> *) const {}
2377 
2378 #endif // LIBMESH_ENABLE_CONSTRAINTS
2379 
2380 
2381 
2382 inline
2384 {
2385  // This got only partly finished...
2386  if (use_constraints)
2387  libmesh_not_implemented();
2388 
2389 #ifdef LIBMESH_ENABLE_CONSTRAINTS
2390  _constrained_sparsity_construction = use_constraints;
2391 #endif
2392  libmesh_ignore(use_constraints);
2393 }
2394 
2395 inline
2397 {
2399 }
2400 
2401 inline
2403 {
2404 #ifdef LIBMESH_ENABLE_CONSTRAINTS
2406 #else
2407  return true;
2408 #endif
2409 }
2410 
2411 inline
2412 void DofMap::should_p_refine(const unsigned int g, const bool p_refine)
2413 {
2414 #ifdef LIBMESH_ENABLE_AMR
2415  if (p_refine)
2416  {
2417  auto it = _dont_p_refine.find(g);
2418  if (it != _dont_p_refine.end())
2419  _dont_p_refine.erase(it);
2420  }
2421  else
2422  _dont_p_refine.insert(g);
2423 #else
2424  libmesh_ignore(g, p_refine);
2425 #endif
2426 }
2427 
2428 inline
2429 bool DofMap::should_p_refine(const unsigned int g) const
2430 {
2431 #ifdef LIBMESH_ENABLE_AMR
2432  return !_dont_p_refine.count(g);
2433 #else
2434  libmesh_ignore(g);
2435  return false;
2436 #endif
2437 }
2438 
2439 inline
2440 unsigned int DofMap::var_group_from_var_number(const unsigned int var_num) const
2441 {
2442  libmesh_assert(var_num < n_variables());
2443  return libmesh_map_find(_var_to_vg, var_num);
2444 }
2445 
2446 inline
2447 bool DofMap::should_p_refine_var(const unsigned int var) const
2448 {
2449 #ifdef LIBMESH_ENABLE_AMR
2450  const auto vg = this->var_group_from_var_number(var);
2451  return !_dont_p_refine.count(vg);
2452 #else
2453  libmesh_ignore(var);
2454  return false;
2455 #endif
2456 }
2457 
2458 template <typename FieldDofsFunctor>
2459 void DofMap::_dof_indices (const Elem & elem,
2460  int p_level,
2461  std::vector<dof_id_type> & di,
2462  const unsigned int vg,
2463  const unsigned int vig,
2464  const Node * const * nodes,
2465  unsigned int n_nodes,
2466  const unsigned int v,
2467 #ifdef DEBUG
2468  std::size_t & tot_size,
2469 #endif
2470  FieldDofsFunctor field_dofs_functor) const
2471 {
2472  const VariableGroup & var = this->variable_group(vg);
2473 
2474  if (var.active_on_subdomain(elem.subdomain_id()))
2475  {
2476  const ElemType type = elem.type();
2477  const unsigned int sys_num = this->sys_number();
2478 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2479  const bool is_inf = elem.infinite();
2480 #endif
2481 
2482  const bool extra_hanging_dofs =
2484 
2485  FEType fe_type = var.type();
2486 
2487  const bool add_p_level =
2488 #ifdef LIBMESH_ENABLE_AMR
2489  !_dont_p_refine.count(vg);
2490 #else
2491  false;
2492 #endif
2493 
2494 #ifdef DEBUG
2495  // The number of dofs per element is non-static for subdivision FE
2496  if (var.type().family == SUBDIVISION)
2497  tot_size += n_nodes;
2498  else
2499  // FIXME: Is the passed-in p_level just elem.p_level()? If so,
2500  // this seems redundant.
2501  tot_size += FEInterface::n_dofs(fe_type, add_p_level*p_level, &elem);
2502 #endif
2503 
2504  // The total Order is not required when getting the function
2505  // pointer, it is only needed when the function is called (see
2506  // below).
2507  const FEInterface::n_dofs_at_node_ptr ndan =
2508  FEInterface::n_dofs_at_node_function(fe_type, &elem);
2509 
2510  // Get the node-based DOF numbers
2511  for (unsigned int n=0; n != n_nodes; n++)
2512  {
2513  const Node & node = *nodes[n];
2514 
2515  // Cache the intermediate lookups that are common to every
2516  // component
2517 #ifdef DEBUG
2518  const std::pair<unsigned int, unsigned int>
2519  vg_and_offset = node.var_to_vg_and_offset(sys_num,v);
2520  libmesh_assert_equal_to (vg, vg_and_offset.first);
2521  libmesh_assert_equal_to (vig, vg_and_offset.second);
2522 #endif
2523  const unsigned int n_comp = node.n_comp_group(sys_num,vg);
2524 
2525  // There is a potential problem with h refinement. Imagine a
2526  // quad9 that has a linear FE on it. Then, on the hanging side,
2527  // it can falsely identify a DOF at the mid-edge node. This is why
2528  // we go through FEInterface instead of node.n_comp() directly.
2529  const unsigned int nc =
2530 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2531  is_inf ?
2532  FEInterface::n_dofs_at_node(fe_type, add_p_level*p_level, &elem, n) :
2533 #endif
2534  ndan (type, fe_type.order + add_p_level*p_level, n);
2535 
2536  // If this is a non-vertex on a hanging node with extra
2537  // degrees of freedom, we use the non-vertex dofs (which
2538  // come in reverse order starting from the end, to
2539  // simplify p refinement)
2540  if (extra_hanging_dofs && !elem.is_vertex(n))
2541  {
2542  const int dof_offset = n_comp - nc;
2543 
2544  // We should never have fewer dofs than necessary on a
2545  // node unless we're getting indices on a parent element,
2546  // and we should never need the indices on such a node
2547  if (dof_offset < 0)
2548  {
2549  libmesh_assert(!elem.active());
2550  di.resize(di.size() + nc, DofObject::invalid_id);
2551  }
2552  else
2553  for (int i=int(n_comp)-1; i>=dof_offset; i--)
2554  {
2555  const dof_id_type d =
2556  node.dof_number(sys_num, vg, vig, i, n_comp);
2557  libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2558  field_dofs_functor(elem, n, v, di, d);
2559  }
2560  }
2561  // If this is a vertex or an element without extra hanging
2562  // dofs, our dofs come in forward order coming from the
2563  // beginning
2564  else
2565  {
2566  // We have a good component index only if it's being
2567  // used on this FE type (nc) *and* it's available on
2568  // this DofObject (n_comp).
2569  const unsigned int good_nc = std::min(n_comp, nc);
2570  for (unsigned int i=0; i!=good_nc; ++i)
2571  {
2572  const dof_id_type d =
2573  node.dof_number(sys_num, vg, vig, i, n_comp);
2574  libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2575  libmesh_assert_less (d, this->n_dofs());
2576  field_dofs_functor(elem, n, v, di, d);
2577  }
2578 
2579  // With fewer good component indices than we need, e.g.
2580  // due to subdomain expansion, the remaining expected
2581  // indices are marked invalid.
2582  if (n_comp < nc)
2583  for (unsigned int i=n_comp; i!=nc; ++i)
2584  di.push_back(DofObject::invalid_id);
2585  }
2586  }
2587 
2588  // If there are any element-based DOF numbers, get them
2589  const unsigned int nc = FEInterface::n_dofs_per_elem(fe_type, add_p_level*p_level, &elem);
2590 
2591  // We should never have fewer dofs than necessary on an
2592  // element unless we're getting indices on a parent element
2593  // (and we should never need those indices) or off-domain for a
2594  // subdomain-restricted variable (where invalid_id is the
2595  // correct thing to return)
2596  if (nc != 0)
2597  {
2598  const unsigned int n_comp = elem.n_comp_group(sys_num,vg);
2599  if (elem.n_systems() > sys_num && nc <= n_comp)
2600  {
2601  for (unsigned int i=0; i<nc; i++)
2602  {
2603  const dof_id_type d =
2604  elem.dof_number(sys_num, vg, vig, i, n_comp);
2605  libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2606 
2607  field_dofs_functor(elem, invalid_uint, v, di, d);
2608  }
2609  }
2610  else
2611  {
2612  libmesh_assert(!elem.active() || fe_type.family == LAGRANGE || fe_type.family == SUBDIVISION);
2613  di.resize(di.size() + nc, DofObject::invalid_id);
2614  }
2615  }
2616  }
2617 }
2618 
2619 
2620 
2621 template <typename ScalarDofsFunctor, typename FieldDofsFunctor>
2622 void DofMap::dof_indices (const Elem * const elem,
2623  std::vector<dof_id_type> & di,
2624  const unsigned int vn,
2625  ScalarDofsFunctor scalar_dofs_functor,
2626  FieldDofsFunctor field_dofs_functor,
2627  int p_level) const
2628 {
2629  // We now allow elem==nullptr to request just SCALAR dofs
2630  // libmesh_assert(elem);
2631 
2632  LOG_SCOPE("dof_indices()", "DofMap");
2633 
2634  // Clear the DOF indices vector
2635  di.clear();
2636 
2637  // Use the default p refinement level?
2638  if (p_level == -12345)
2639  p_level = elem ? elem->p_level() : 0;
2640 
2641  const unsigned int vg = this->_variable_group_numbers[vn];
2642  const VariableGroup & var = this->variable_group(vg);
2643  const unsigned int vig = vn - var.number();
2644 
2645 #ifdef DEBUG
2646  // Check that sizes match in DEBUG mode
2647  std::size_t tot_size = 0;
2648 #endif
2649 
2650  if (elem && elem->type() == TRI3SUBDIVISION)
2651  {
2652  // Subdivision surface FE require the 1-ring around elem
2653  const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *>(elem);
2654 
2655  // Ghost subdivision elements have no real dofs
2656  if (!sd_elem->is_ghost())
2657  {
2658  // Determine the nodes contributing to element elem
2659  std::vector<const Node *> elem_nodes;
2660  MeshTools::Subdivision::find_one_ring(sd_elem, elem_nodes);
2661 
2662  _dof_indices(*elem, p_level, di, vg, vig, elem_nodes.data(),
2663  cast_int<unsigned int>(elem_nodes.size()), vn,
2664 #ifdef DEBUG
2665  tot_size,
2666 #endif
2667  field_dofs_functor);
2668  }
2669 
2670  return;
2671  }
2672 
2673  // Get the dof numbers
2674  if (var.type().family == SCALAR &&
2675  (!elem ||
2676  var.active_on_subdomain(elem->subdomain_id())))
2677  {
2678 #ifdef DEBUG
2679  tot_size += var.type().order;
2680 #endif
2681  std::vector<dof_id_type> di_new;
2682  this->SCALAR_dof_indices(di_new,vn);
2683  scalar_dofs_functor(*elem, di, di_new);
2684  }
2685  else if (elem)
2686  _dof_indices(*elem, p_level, di, vg, vig, elem->get_nodes(),
2687  elem->n_nodes(), vn,
2688 #ifdef DEBUG
2689  tot_size,
2690 #endif
2691  field_dofs_functor);
2692 
2693 #ifdef DEBUG
2694  libmesh_assert_equal_to (tot_size, di.size());
2695 #endif
2696 }
2697 
2698 inline
2700 {
2702  return *_sc;
2703 }
2704 
2705 } // namespace libMesh
2706 
2707 #endif // LIBMESH_DOF_MAP_H
std::vector< VariableGroup > _variable_groups
The finite element type for each variable group.
Definition: dof_map.h:1957
void remove_adjoint_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary, unsigned int q)
Removes from the system the specified Dirichlet boundary for the adjoint equation defined by Quantity...
void find_connected_dofs(std::vector< dof_id_type > &elem_dofs) const
Finds all the DOFS associated with the element DOFs elem_dofs.
Definition: dof_map.C:2920
std::unique_ptr< SparsityPattern::Build > _sp
The sparsity pattern of the global matrix.
Definition: dof_map.h:2089
static unsigned int n_dofs_per_elem(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:530
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
std::string get_local_constraints(bool print_nonlocal=false) const
Gets a string reporting all DoF and Node constraints local to this processor.
FEFamily family
The type of finite element.
Definition: fe_type.h:221
void constrain_element_dyad_matrix(DenseVector< Number > &v, DenseVector< Number > &w, std::vector< dof_id_type > &row_dofs, bool asymmetric_constraint_rows=true) const
Constrains a dyadic element matrix B = v w&#39;.
Definition: dof_map.h:2353
A class holding degree of freedom information pertinent to static condensation.
bool _implicit_neighbor_dofs_initialized
Bools to indicate if we override the –implicit_neighbor_dofs commandline options.
Definition: dof_map.h:2159
ElemType
Defines an enum for geometric element types.
bool is_constrained_node(const Node *node) const
Definition: dof_map.h:2252
std::vector< GhostingFunctor * >::const_iterator GhostingFunctorIterator
Iterator type for coupling and algebraic ghosting functor ranges.
Definition: dof_map.h:315
void allgather(const T &send_data, std::vector< T, A > &recv_data) const
unsigned int n_variable_groups() const
Definition: dof_map.h:632
void add_adjoint_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary, unsigned int q)
Adds a copy of the specified Dirichlet boundary to the system, corresponding to the adjoint problem d...
const unsigned int _sys_number
The number of the system we manage DOFs for.
Definition: dof_map.h:1972
bool _implicit_neighbor_dofs
Definition: dof_map.h:2160
DofConstraintValueMap & get_primal_constraint_values()
Definition: dof_map.h:2313
void check_dirichlet_bcid_consistency(const MeshBase &mesh, const DirichletBoundary &boundary) const
Check that all the ids in dirichlet_bcids are actually present in the mesh.
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
Definition: dof_object.h:1032
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
DefaultCoupling & default_coupling()
Default coupling functor.
Definition: dof_map.h:378
const std::vector< dof_id_type > & get_n_oz() const
Definition: dof_map.h:555
A Node is like a Point, but with more information.
Definition: node.h:52
static unsigned int n_dofs(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:355
This abstract base class defines the interface by which library code and user code can report associa...
dof_id_type n_SCALAR_dofs() const
Definition: dof_map.h:685
void build_constraint_matrix_and_vector(DenseMatrix< Number > &C, DenseVector< Number > &H, std::vector< dof_id_type > &elem_dofs, int qoi_index=-1, const bool called_recursively=false) const
Build the constraint matrix C and the forcing vector H associated with the element degree of freedom ...
~DofMap()
Destructor.
Definition: dof_map.C:202
This helper class can be called on multiple threads to compute the sparsity pattern (or graph) of the...
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
Definition: libmesh.h:310
void local_variable_indices(T &idx, unsigned int var_num) const
If T == dof_id_type, counts, if T == std::vector<dof_id_type>, fills an array of, those dof indices w...
Definition: dof_map.h:914
processor_id_type dof_owner(const dof_id_type dof) const
Definition: dof_map.h:714
void scatter_constraints(MeshBase &)
Sends constraint equations to constraining processors.
void reinit(MeshBase &mesh, const std::map< const Node *, std::set< subdomain_id_type >> &constraining_subdomains)
Reinitialize the underlying data structures conformal to the current mesh.
Definition: dof_map.C:503
void add_periodic_boundary(const PeriodicBoundaryBase &periodic_boundary)
Adds a copy of the specified periodic boundary to the system.
bool _error_on_constraint_loop
This flag indicates whether or not we do an opt-mode check for the presence of constraint loops...
Definition: dof_map.h:1941
void add_adjoint_constraint_row(const unsigned int qoi_index, const dof_id_type dof_number, const DofConstraintRow &constraint_row, const Number constraint_rhs, const bool forbid_constraint_overwrite)
Adds a copy of the user-defined row to the constraint matrix, using an inhomogeneous right-hand-side ...
void process_mesh_constraint_rows(const MeshBase &mesh)
Adds any spline constraints from the Mesh to our DoF constraints.
void * _extra_sparsity_context
A pointer associated with the extra sparsity that can optionally be passed in.
Definition: dof_map.h:2013
void extract_local_vector(const NumericVector< Number > &Ug, const std::vector< dof_id_type > &dof_indices, DenseVectorBase< Number > &Ue) const
Builds the local element vector Ue from the global vector Ug, accounting for any constrained degrees ...
Definition: dof_map.C:2147
void constrain_element_matrix_and_vector(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true) const
Constrains the element matrix and vector.
Definition: dof_map.h:2340
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Definition: dof_map.C:2229
dof_id_type n_dofs() const
Definition: dof_map_base.h:105
The definition of the const_element_iterator struct.
Definition: mesh_base.h:2222
std::size_t distribute_dofs(MeshBase &)
Distribute dofs on the current mesh.
Definition: dof_map.C:978
void set_implicit_neighbor_dofs(bool implicit_neighbor_dofs)
Allow the implicit_neighbor_dofs flag to be set programmatically.
Definition: dof_map.C:1920
void add_default_ghosting()
Add the default functor(s) for coupling and algebraic ghosting.
Definition: dof_map.C:2024
We&#39;re using a class instead of a typedef to allow forward declarations and future flexibility...
void local_variable_indices(T &idx, const MeshBase &mesh, unsigned int var_num) const
If T == dof_id_type, counts, if T == std::vector<dof_id_type>, fills an array of, those dof indices w...
Definition: dof_map.C:1151
bool is_periodic_boundary(const boundary_id_type boundaryid) const
Definition: dof_map.C:215
Number has_heterogenous_adjoint_constraint(const unsigned int qoi_num, const dof_id_type dof) const
Backwards compatibility with misspelling.
Definition: dof_map.h:1132
GhostingFunctorIterator algebraic_ghosting_functors_begin() const
Beginning of range of algebraic ghosting functors.
Definition: dof_map.h:428
dof_id_type n_local_constrained_dofs() const
void remove_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Removes the specified Dirichlet boundary from the system.
const FEType & variable_group_type(const unsigned int vg) const
Definition: dof_map.h:2240
std::unique_ptr< StaticCondensationDofMap > _sc
Static condensation class.
Definition: dof_map.h:2175
NodeConstraints::const_iterator node_constraint_rows_begin() const
Definition: dof_map.h:1091
std::vector< dof_id_type > _send_list
A list containing all the global DOF indices that affect the solution on my processor.
Definition: dof_map.h:1996
const FEType & variable_type(const unsigned int c) const
Definition: dof_map.h:2230
std::unique_ptr< SparsityPattern::Build > build_sparsity(const MeshBase &mesh, bool calculate_constrained=false, bool use_condensed_system=false) const
Builds a sparsity pattern for matrices using the current degree-of-freedom numbering and coupling...
Definition: dof_map.C:61
void set_verify_dirichlet_bc_consistency(bool val)
Set the _verify_dirichlet_bc_consistency flag.
Definition: dof_map.C:1926
unsigned int block_size() const
Definition: dof_map.h:665
void sum(T &r) const
bool is_attached(SparseMatrix< Number > &matrix)
Matrices should not be attached more than once.
Definition: dof_map.C:329
void attach_matrix(SparseMatrix< Number > &matrix)
Additional matrices may be attached to this DofMap.
Definition: dof_map.C:274
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
void clear_send_list()
Clears the _send_list vector.
Definition: dof_map.h:507
MeshBase & mesh
void gather_constraints(MeshBase &mesh, std::set< dof_id_type > &unexpanded_dofs, bool look_for_constrainees)
Helper function for querying about constraint equations on other processors.
DefaultCoupling & default_algebraic_ghosting()
Default algebraic ghosting functor.
Definition: dof_map.h:440
std::unique_ptr< DirichletBoundaries > _dirichlet_boundaries
Data structure containing Dirichlet functions.
Definition: dof_map.h:2144
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
Definition: vector_fe_ex5.C:44
const Parallel::Communicator & comm() const
The Node constraint storage format.
Definition: dof_map.h:156
unsigned int p_level() const
Definition: elem.h:3108
OrderWrapper order
The approximation order of the element.
Definition: fe_type.h:215
bool use_coupled_neighbor_dofs(const MeshBase &mesh) const
Tells other library functions whether or not this problem includes coupling between dofs in neighbori...
Definition: dof_map.C:1932
std::map< const Elem *, const CouplingMatrix *, CompareDofObjectsByPIDAndThenID > map_type
What elements do we care about and what variables do we care about on each element?
void reinit_static_condensation()
Calls reinit on the static condensation map if it exists.
Definition: dof_map.C:3145
void should_p_refine(unsigned int g, bool p_refine)
Describe whether the given variable group should be p-refined.
Definition: dof_map.h:2412
dof_id_type n_dofs(const unsigned int vn) const
Definition: dof_map.h:675
GhostingFunctorIterator algebraic_ghosting_functors_end() const
End of range of algebraic ghosting functors.
Definition: dof_map.h:434
virtual void augment_send_list(std::vector< dof_id_type > &send_list)=0
User-defined function to augment the send list.
The libMesh namespace provides an interface to certain functionality in the library.
void add_constraint_row(const dof_id_type dof_number, const DofConstraintRow &constraint_row, const bool forbid_constraint_overwrite=true)
Adds a copy of the user-defined row to the constraint matrix, using a homogeneous right-hand-side for...
Definition: dof_map.h:1035
dof_id_type n_local_dofs(const unsigned int vn) const
Definition: dof_map.h:693
void set_error_on_constraint_loop(bool error_on_constraint_loop)
Definition: dof_map.C:241
void add_coupling_functor(GhostingFunctor &coupling_functor, bool to_mesh=true)
Adds a functor which can specify coupling requirements for creation of sparse matrices.
Definition: dof_map.C:2033
const SparsityPattern::Build * get_sparsity_pattern() const
Definition: dof_map.h:570
std::string get_info() const
Gets summary info about the sparsity bandwidth and constraints.
Definition: dof_map.C:2990
unsigned int sys_number() const
Definition: dof_map.h:2182
void SCALAR_dof_indices(std::vector< dof_id_type > &di, const unsigned int vn, const bool old_dofs=false) const
Fills the vector di with the global degree of freedom indices corresponding to the SCALAR variable vn...
Definition: dof_map.C:2612
uint8_t processor_id_type
Definition: id_types.h:104
This is the MeshBase class.
Definition: mesh_base.h:75
std::vector< dof_id_type > _first_scalar_df
First DOF index for SCALAR variable v, or garbage for non-SCALAR variable v.
Definition: dof_map.h:1990
const Variable & variable(const unsigned int c) const override
Definition: dof_map.h:2200
This class implements the default algebraic coupling in libMesh: elements couple to themselves...
void enforce_adjoint_constraints_exactly(NumericVector< Number > &v, unsigned int q) const
Heterogeneously constrains the numeric vector v, which represents an adjoint solution defined on the ...
Definition: dof_map.h:2364
AdjointDofConstraintValues _adjoint_constraint_values
Definition: dof_map.h:2120
bool _constrained_sparsity_construction
This flag indicates whether or not we explicitly take constraint equations into account when computin...
Definition: dof_map.h:1947
bool should_p_refine_var(unsigned int var) const
Whether the given variable should be p-refined.
Definition: dof_map.h:2447
unsigned int var_group_from_var_number(unsigned int var_num) const
Definition: dof_map.h:2440
dof_id_type end_dof() const
Definition: dof_map_base.h:83
AugmentSendList * _augment_send_list
Function object to call to add extra entries to the send list.
Definition: dof_map.h:2018
void heterogenously_constrain_element_vector(const DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true, int qoi_index=-1) const
Definition: dof_map.h:1311
bool need_full_sparsity_pattern
Default false; set to true if any attached matrix requires a full sparsity pattern.
Definition: dof_map.h:2082
bool constrained_sparsity_construction()
Returns true iff the current policy when constructing sparsity patterns is to explicitly account for ...
Definition: dof_map.h:2402
void distribute_local_dofs_var_major(dof_id_type &next_free_dof, MeshBase &mesh, const std::map< const Node *, std::set< subdomain_id_type >> &constraining_subdomains)
Distributes the global degrees of freedom, for dofs on this processor.
Definition: dof_map.C:1447
Generic sparse matrix.
Definition: vector_fe_ex5.C:46
void enforce_constraints_on_residual(const NonlinearImplicitSystem &system, NumericVector< Number > *rhs, NumericVector< Number > const *solution, bool homogeneous=true) const
Definition: dof_map.h:2369
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
virtual void clear() override
Free all new memory associated with the object, but restore its original state, with the mesh pointer...
Definition: dof_map.C:901
void _node_dof_indices(const Elem &elem, unsigned int n, const DofObject &obj, std::vector< dof_id_type > &di, const unsigned int vn) const
Helper function that implements the element-nodal versions of dof_indices and old_dof_indices.
Definition: dof_map.C:2496
void add_algebraic_ghosting_functor(std::shared_ptr< GhostingFunctor > evaluable_functor, bool to_mesh=true)
Adds a functor which can specify algebraic ghosting requirements for use with distributed vectors...
Definition: dof_map.h:413
void find_one_ring(const Tri3Subdivision *elem, std::vector< const Node *> &nodes)
Determines the 1-ring of element elem, and writes it to the nodes vector.
void reinit_send_list(MeshBase &mesh)
Clears the _send_list vector and then rebuilds it.
Definition: dof_map.C:1906
processor_id_type n_processors() const
void libmesh_ignore(const Args &...)
const dof_id_type n_nodes
Definition: tecplot_io.C:67
dof_id_type first_dof() const
Definition: dof_map_base.h:73
std::unordered_map< unsigned int, unsigned int > _var_to_vg
A map from variable number to variable group number.
Definition: dof_map.h:1967
This class defines the notion of a variable in the system.
Definition: variable.h:50
bool has_static_condensation() const
Checks whether we have static condensation.
Definition: dof_map.h:1675
void add_neighbors_to_send_list(MeshBase &mesh)
Adds entries to the _send_list vector corresponding to DoFs on elements neighboring the current proce...
Definition: dof_map.C:1716
std::vector< dof_id_type > _first_old_scalar_df
First old DOF index for SCALAR variable v, or garbage for non-SCALAR variable v.
Definition: dof_map.h:2103
int8_t boundary_id_type
Definition: id_types.h:51
bool has_adjoint_dirichlet_boundaries(unsigned int q) const
void(* _extra_sparsity_function)(SparsityPattern::Graph &, std::vector< dof_id_type > &n_nz, std::vector< dof_id_type > &n_oz, void *)
A function pointer to a function to call to add extra entries to the sparsity pattern.
Definition: dof_map.h:2006
dof_id_type _n_SCALAR_dofs
The total number of SCALAR dofs associated to all SCALAR variables.
Definition: dof_map.h:2095
void assert_no_nodes_missed(MeshBase &mesh)
Definition: dof_map.C:1592
PeriodicBoundaries * get_periodic_boundaries()
Definition: dof_map.h:1462
StaticCondensationDofMap & get_static_condensation()
Definition: dof_map.h:2699
void set_nonlocal_dof_objects(iterator_type objects_begin, iterator_type objects_end, MeshBase &mesh, dofobject_accessor objects)
Helper function for distributing dofs in parallel.
Definition: dof_map.C:352
void add_constraint_row(const dof_id_type dof_number, const DofConstraintRow &constraint_row, const Number constraint_rhs, const bool forbid_constraint_overwrite)
Adds a copy of the user-defined row to the constraint matrix, using an inhomogeneous right-hand-side ...
static bool extra_hanging_dofs(const FEType &fe_t)
virtual unsigned int n_nodes() const =0
std::set< std::unique_ptr< CouplingMatrix >, Utility::CompareUnderlying > CouplingMatricesSet
Definition: dof_map.h:1855
DofConstraints::const_iterator constraint_rows_begin() const
Definition: dof_map.h:1043
void heterogeneously_constrain_element_vector(const DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true, int qoi_index=-1) const
Constrains the element vector.
Definition: dof_map.h:2350
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
unsigned int n_variables() const override
Definition: dof_map.h:635
unsigned int n_systems() const
Definition: dof_object.h:937
void add_variable_group(VariableGroup var_group)
Add an unknown of order order and finite element type type to the system of equations.
Definition: dof_map.C:248
void swap_dof_constraints()
Similar to the stash/unstash_dof_constraints() API, but swaps _dof_constraints and _stashed_dof_const...
Definition: dof_map.h:1082
NodeConstraints::const_iterator node_constraint_rows_end() const
Definition: dof_map.h:1097
This base class provides a minimal set of interfaces for satisfying user requests for...
Definition: dof_map_base.h:51
dof_id_type n_constrained_nodes() const
Definition: dof_map.h:938
const Node *const * get_nodes() const
Definition: elem.h:2499
We&#39;re using a class instead of a typedef to allow forward declarations and future flexibility...
void print_dof_constraints(std::ostream &os=libMesh::out, bool print_nonlocal=false) const
Prints (from processor 0) all DoF and Node constraints.
libmesh_assert(ctx)
void attach_extra_sparsity_object(SparsityPattern::AugmentSparsityPattern &asp)
Attach an object to use to populate the sparsity pattern with extra entries.
Definition: dof_map.h:452
DirichletBoundaries * get_dirichlet_boundaries()
Definition: dof_map.h:1520
bool is_constrained_dof(const dof_id_type dof) const
Definition: dof_map.h:2268
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and non-linear solv...
The Tri3Subdivision element is a three-noded subdivision surface shell element used in mechanics calc...
static const dof_id_type invalid_id
An invalid id to distinguish an uninitialized DofObject.
Definition: dof_object.h:482
void allgather_recursive_constraints(MeshBase &)
Gathers constraint equation dependencies from other processors.
const VariableGroup & variable_group(const unsigned int c) const
Definition: dof_map.h:2190
std::vector< GhostingFunctor * > _coupling_functors
The list of all GhostingFunctor objects to be used when coupling degrees of freedom in matrix sparsit...
Definition: dof_map.h:2070
DofConstraints _dof_constraints
Data structure containing DOF constraints.
Definition: dof_map.h:2116
void add_coupling_functor(std::shared_ptr< GhostingFunctor > coupling_functor, bool to_mesh=true)
Adds a functor which can specify coupling requirements for creation of sparse matrices.
Definition: dof_map.h:351
DofMap(const unsigned int sys_number, MeshBase &mesh)
Constructor.
Definition: dof_map.C:136
CouplingMatrix * _dof_coupling
Degree of freedom coupling.
Definition: dof_map.h:1620
static void merge_ghost_functor_outputs(GhostingFunctor::map_type &elements_to_ghost, CouplingMatricesSet &temporary_coupling_matrices, const GhostingFunctorIterator &gf_begin, const GhostingFunctorIterator &gf_end, const MeshBase::const_element_iterator &elems_begin, const MeshBase::const_element_iterator &elems_end, processor_id_type p)
Definition: dof_map.C:1614
bool active_on_subdomain(subdomain_id_type sid) const
Definition: variable.h:167
void create_static_condensation(MeshBase &mesh, System &system)
Add a static condensation class.
Definition: dof_map.C:3140
void print_info(std::ostream &os=libMesh::out) const
Prints summary info about the sparsity bandwidth and constraints.
Definition: dof_map.C:2983
void * _extra_send_list_context
A pointer associated with the extra send list that can optionally be passed in.
Definition: dof_map.h:2028
Number has_heterogeneous_adjoint_constraint(const unsigned int qoi_num, const dof_id_type dof) const
Definition: dof_map.h:2292
This class implements reference counting.
const DirichletBoundaries * get_dirichlet_boundaries() const
Definition: dof_map.h:1515
bool computed_sparsity_already() const
Returns true iff a sparsity pattern has already been computed.
Definition: dof_map.C:292
This class defines a logically grouped set of variables in the system.
Definition: variable.h:203
void update_sparsity_pattern(SparseMatrix< Number > &matrix) const
Additional matrices may be be temporarily initialized by this DofMap.
Definition: dof_map.C:303
std::pair< Real, Real > max_constraint_error(const System &system, NumericVector< Number > *v=nullptr) const
Tests the constrained degrees of freedom on the numeric vector v, which represents a solution defined...
static n_dofs_at_node_ptr n_dofs_at_node_function(const unsigned int dim, const FEType &fe_t)
Definition: fe_interface.C:458
Storage for DofConstraint right hand sides for a particular problem.
Definition: dof_map.h:120
std::vector< SparseMatrix< Number > *> _matrices
Additional matrices handled by this object.
Definition: dof_map.h:1984
GhostingFunctorIterator coupling_functors_end() const
End of range of coupling functors.
Definition: dof_map.h:372
DofObject * elem_ptr(MeshBase &mesh, dof_id_type i) const
Definition: dof_map.C:344
Storage for DofConstraint right hand sides for all adjoint problems.
Definition: dof_map.h:131
std::vector< GhostingFunctor * > _algebraic_ghosting_functors
The list of all GhostingFunctor objects to be used when distributing ghosted vectors.
Definition: dof_map.h:2057
std::vector< std::unique_ptr< DirichletBoundaries > > _adjoint_dirichlet_boundaries
Data structure containing Dirichlet functions.
Definition: dof_map.h:2150
static unsigned int n_dofs_at_node(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)
Definition: fe_interface.C:436
void distribute_scalar_dofs(dof_id_type &next_free_dof)
Definition: dof_map.C:1568
DofConstraints::const_iterator constraint_rows_end() const
Definition: dof_map.h:1049
std::unique_ptr< PeriodicBoundaries > _periodic_boundaries
Data structure containing periodic boundaries.
Definition: dof_map.h:2136
void heterogenously_constrain_element_matrix_and_vector(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true, int qoi_index=-1) const
Definition: dof_map.h:1267
bool has_heterogenous_adjoint_constraints(const unsigned int qoi_num) const
Backwards compatibility with misspelling.
Definition: dof_map.h:1116
void unstash_dof_constraints()
Definition: dof_map.h:1064
void remove_coupling_functor(GhostingFunctor &coupling_functor)
Removes a functor which was previously added to the set of coupling functors, from both this DofMap a...
Definition: dof_map.C:2060
unsigned int(* n_dofs_at_node_ptr)(const ElemType, const Order, const unsigned int)
Definition: fe_interface.h:151
DofObject *(DofMap::* dofobject_accessor)(MeshBase &mesh, dof_id_type i) const
A member function type like node_ptr() or elem_ptr().
Definition: dof_map.h:1772
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Struct which defines a custom comparison object that can be used with std::sets of std::unique_ptrs...
Definition: utility.h:517
void process_constraints(MeshBase &)
Postprocesses any constrained degrees of freedom to be constrained only in terms of unconstrained dof...
subdomain_id_type subdomain_id() const
Definition: elem.h:2582
void stash_dof_constraints()
Definition: dof_map.h:1058
dof_id_type n_constrained_dofs() const
void constrain_element_vector(DenseVector< Number > &rhs, std::vector< dof_id_type > &dofs, bool asymmetric_constraint_rows=true) const
Constrains the element vector.
Definition: dof_map.h:2336
const DofConstraints & get_dof_constraints() const
Provide a const accessor to the DofConstraints map.
Definition: dof_map.h:1056
const std::vector< dof_id_type > & get_n_nz() const
Definition: dof_map.h:542
void _dof_indices(const Elem &elem, int p_level, std::vector< dof_id_type > &di, const unsigned int vg, const unsigned int vig, const Node *const *nodes, unsigned int n_nodes, const unsigned int v #ifdef DEBUG, std::size_t &tot_size #endif) const
Helper function that gets the dof indices on the current element for a non-SCALAR type variable...
Definition: dof_map.C:2580
void invalidate_dofs(MeshBase &mesh) const
Invalidates all active DofObject dofs for this system.
Definition: dof_map.C:886
const DirichletBoundaries * get_adjoint_dirichlet_boundaries(unsigned int q) const
void attach_extra_send_list_function(void(*func)(std::vector< dof_id_type > &, void *), void *context=nullptr)
Attach a function pointer to use as a callback to populate the send_list with extra entries...
Definition: dof_map.h:490
virtual bool is_vertex(const unsigned int i) const =0
bool _verify_dirichlet_bc_consistency
Flag which determines whether we should do some additional checking of the consistency of the Dirichl...
Definition: dof_map.h:2172
std::unique_ptr< DefaultCoupling > _default_coupling
The default coupling GhostingFunctor, used to implement standard libMesh sparsity pattern constructio...
Definition: dof_map.h:2036
void add_constraints_to_send_list()
Adds entries to the _send_list vector corresponding to DoFs which are dependencies for constraint equ...
OStreamProxy out
DofConstraintValueMap _primal_constraint_values
Definition: dof_map.h:2118
void heterogeneously_constrain_element_matrix_and_vector(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true, int qoi_index=-1) const
Constrains the element matrix and vector.
Definition: dof_map.h:2346
std::vector< Variable > _variables
The finite element type for each variable.
Definition: dof_map.h:1952
DofConstraints _stashed_dof_constraints
Definition: dof_map.h:2116
void remove_default_ghosting()
Remove any default ghosting functor(s).
Definition: dof_map.C:2016
void attach_extra_sparsity_function(void(*func)(SparsityPattern::Graph &sparsity, std::vector< dof_id_type > &n_nz, std::vector< dof_id_type > &n_oz, void *), void *context=nullptr)
Attach a function pointer to use as a callback to populate the sparsity pattern with extra entries...
Definition: dof_map.h:467
SparsityPattern::AugmentSparsityPattern * _augment_sparsity_pattern
Function object to call to add extra entries to the sparsity pattern.
Definition: dof_map.h:2001
void find_connected_dof_objects(std::vector< const DofObject *> &objs) const
Finds all the DofObjects associated with the set in objs.
std::map< dof_id_type, Real, std::less< dof_id_type >, Threads::scalable_allocator< std::pair< const dof_id_type, Real > > > DofConstraintRow
A row of the Dof constraint matrix.
Definition: dof_map.h:100
The DofObject defines an abstract base class for objects that have degrees of freedom associated with...
Definition: dof_object.h:54
void clear_sparsity()
Clears the sparsity pattern.
Definition: dof_map.C:2009
void set_constrained_sparsity_construction(bool use_constraints)
Sets the current policy for constructing sparsity patterns: if use_constraints is true (for robustnes...
Definition: dof_map.h:2383
Order variable_order(const unsigned int c) const
Definition: dof_map.h:2210
dof_id_type n_local_dofs() const
Definition: dof_map_base.h:115
std::unordered_set< unsigned int > _dont_p_refine
A container of variable groups that we should not p-refine.
Definition: dof_map.h:2108
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
unsigned int number(unsigned int v) const
Definition: variable.h:304
void build_constraint_matrix(DenseMatrix< Number > &C, std::vector< dof_id_type > &elem_dofs, const bool called_recursively=false) const
Build the constraint matrix C associated with the element degree of freedom indices elem_dofs...
std::pair< unsigned int, unsigned int > var_to_vg_and_offset(const unsigned int s, const unsigned int var) const
Definition: dof_object.h:1200
void(* _extra_send_list_function)(std::vector< dof_id_type > &, void *)
A function pointer to a function to call to add extra entries to the send list.
Definition: dof_map.h:2023
Scalable allocator to be used in multithreaded code chunks which allocate a lot of dynamic memory...
DofObject * node_ptr(MeshBase &mesh, dof_id_type i) const
Definition: dof_map.C:337
void heterogeneously_constrain_element_residual(DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, NumericVector< Number > &solution_local) const
Constrains the element residual.
The base class for defining periodic boundaries.
void create_dof_constraints(const MeshBase &, Real time=0)
Rebuilds the raw degree of freedom and DofObject constraints, based on attached DirichletBoundary obj...
Defines an abstract dense vector base class for use in Finite Element-type computations.
Definition: dof_map.h:73
std::map< const Node *, Real, std::less< const Node * >, Threads::scalable_allocator< std::pair< const Node *const, Real > > > NodeConstraintRow
A row of the Node constraint mapping.
Definition: dof_map.h:148
MeshBase & _mesh
The mesh that system uses.
Definition: dof_map.h:1977
void prepare_send_list()
Takes the _send_list vector (which may have duplicate entries) and sorts it.
Definition: dof_map.C:1864
bool all_semilocal_indices(const std::vector< dof_id_type > &dof_indices) const
Definition: dof_map.C:2664
Defines a dense vector for use in Finite Element-type computations.
Definition: dof_map.h:74
virtual bool infinite() const =0
std::unique_ptr< DefaultCoupling > _default_evaluating
The default algebraic GhostingFunctor, used to implement standard libMesh send_list construction...
Definition: dof_map.h:2044
void full_sparsity_pattern_needed()
Sets need_full_sparsity_pattern to true regardless of the requirements by matrices.
Definition: dof_map.h:2396
std::vector< unsigned int > _variable_group_numbers
The variable group number for each variable.
Definition: dof_map.h:1962
std::map< GhostingFunctor *, std::shared_ptr< GhostingFunctor > > _shared_functors
Hang on to references to any GhostingFunctor objects we were passed in shared_ptr form...
Definition: dof_map.h:2076
void heterogeneously_constrain_element_jacobian_and_residual(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, NumericVector< Number > &solution_local) const
Constrains the element Jacobian and residual.
FEFamily
defines an enum for finite element families.
void remove_algebraic_ghosting_functor(GhostingFunctor &evaluable_functor)
Removes a functor which was previously added to the set of algebraic ghosting functors, from both this DofMap and from the underlying mesh.
Definition: dof_map.C:2117
virtual ~AugmentSendList()=default
Backwards compatibility for prior AugmentSparsityPattern users.
Definition: dof_map.h:201
void attach_extra_send_list_object(DofMap::AugmentSendList &asl)
Attach an object to populate the send_list with extra entries.
Definition: dof_map.h:481
Order variable_group_order(const unsigned int vg) const
Definition: dof_map.h:2220
void compute_sparsity(const MeshBase &)
Computes the sparsity pattern for the matrices corresponding to proc_id and sends that data to Linear...
Definition: dof_map.C:1988
void add_algebraic_ghosting_functor(GhostingFunctor &evaluable_functor, bool to_mesh=true)
Adds a functor which can specify algebraic ghosting requirements for use with distributed vectors...
Definition: dof_map.C:2090
unsigned int n_comp_group(const unsigned int s, const unsigned int vg) const
Definition: dof_object.h:1015
bool has_heterogeneous_adjoint_constraints(const unsigned int qoi_num) const
Definition: dof_map.h:2278
void constrain_element_residual(DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, NumericVector< Number > &solution_local) const
Constrains the element residual.
const PeriodicBoundaries * get_periodic_boundaries() const
Definition: dof_map.h:1467
std::map< const Node *, std::set< subdomain_id_type > > calculate_constraining_subdomains()
We may have mesh constraint rows with dependent nodes in one subdomain but dependency nodes in anothe...
Definition: dof_map.C:1285
Defines a dense matrix for use in Finite Element-type computations.
Definition: dof_map.h:75
bool active() const
Definition: elem.h:2941
void set_error_on_cyclic_constraint(bool error_on_cyclic_constraint)
Specify whether or not we perform an extra (opt-mode enabled) check for constraint loops...
Definition: dof_map.C:234
virtual ElemType type() const =0
Abstract base class to be used to add user-defined implicit degree of freedom couplings.
bool has_blocked_representation() const
Definition: dof_map.h:648
The constraint matrix storage format.
Definition: dof_map.h:108
void ErrorVector unsigned int
Definition: adjoints_ex3.C:360
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:533
std::vector< dof_id_type > n_dofs_per_processor(const unsigned int vn) const
Definition: dof_map.h:704
Abstract base class to be used to add user-defined parallel degree of freedom couplings.
Definition: dof_map.h:208
void constrain_nothing(std::vector< dof_id_type > &dofs) const
Does not actually constrain anything, but modifies dofs in the same way as any of the constrain funct...
Definition: dof_map.h:2358
bool is_evaluable(const DofObjectSubclass &obj, unsigned int var_num=libMesh::invalid_uint) const
Definition: dof_map.C:2677
bool semilocal_index(dof_id_type dof_index) const
Definition: dof_map.C:2648
void distribute_local_dofs_node_major(dof_id_type &next_free_dof, MeshBase &mesh, const std::map< const Node *, std::set< subdomain_id_type >> &constraining_subdomains)
Distributes the global degrees of freedom for dofs on this processor.
Definition: dof_map.C:1319
void constrain_p_dofs(unsigned int var, const Elem *elem, unsigned int s, unsigned int p)
Constrains degrees of freedom on side s of element elem which correspond to variable number var and t...
void check_for_cyclic_constraints()
Throw an error if we detect any constraint loops, i.e.
unsigned int idx(const ElemType type, const unsigned int nx, const unsigned int i, const unsigned int j)
A useful inline function which replaces the macros used previously.
void old_dof_indices(const Elem &elem, unsigned int n, std::vector< dof_id_type > &di, const unsigned int vn) const
Appends to the vector di the old global degree of freedom indices for elem.node_ref(n), for one variable vn.
Definition: dof_map.C:2483
GhostingFunctorIterator coupling_functors_begin() const
Beginning of range of coupling functors.
Definition: dof_map.h:366
uint8_t dof_id_type
Definition: id_types.h:67
std::vector< dof_id_type > _end_df
Last DOF index (plus 1) on processor p.
Definition: dof_map_base.h:159
bool local_index(dof_id_type dof_index) const
Definition: dof_map.h:846
const FEType & type() const
Definition: variable.h:144
NodeConstraints _node_constraints
Data structure containing DofObject constraints.
Definition: dof_map.h:2127
void enforce_constraints_exactly(const System &system, NumericVector< Number > *v=nullptr, bool homogeneous=false) const
Constrains the numeric vector v, which represents a solution defined on the mesh. ...
Definition: dof_map.h:2360
This class defines a coupling matrix.
void constrain_element_matrix(DenseMatrix< Number > &matrix, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true) const
Constrains the element matrix.
Definition: dof_map.h:2327
void enforce_constraints_on_jacobian(const NonlinearImplicitSystem &system, SparseMatrix< Number > *jac) const
Definition: dof_map.h:2375