libMesh
static_condensation_dof_map.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 #include "libmesh/static_condensation_dof_map.h"
19 
20 #include "libmesh/mesh_base.h"
21 #include "libmesh/dof_map.h"
22 #include "libmesh/elem.h"
23 #include "libmesh/int_range.h"
24 #include "libmesh/system.h"
25 #include "libmesh/equation_systems.h"
26 #include "timpi/parallel_sync.h"
27 #include <unordered_set>
28 
29 namespace libMesh
30 {
32  System & system,
33  const DofMap & dof_map)
34  : DofMapBase(dof_map.comm()),
35  _mesh(mesh),
36  _system(system),
37  _dof_map(dof_map),
38  _sc_is_initialized(false)
39 {
40  libmesh_experimental();
41 }
42 
44 
46  const dof_id_type full_dof_number,
47  const bool involved_in_constraints,
48  std::unordered_map<dof_id_type, dof_id_type> & uncondensed_global_to_local_map,
49  std::unordered_set<dof_id_type> & local_uncondensed_dofs_set,
50  std::unordered_map<processor_id_type, std::unordered_set<dof_id_type>> &
51  nonlocal_uncondensed_dofs,
52  std::vector<dof_id_type> & elem_uncondensed_dofs,
53  dof_id_type & uncondensed_local_dof_number,
54  std::unordered_set<dof_id_type> & constraint_dofs)
55 {
56  if (uncondensed_global_to_local_map.count(full_dof_number))
57  // We've already seen this dof on this element
58  return;
59 
60  if (_dof_map.local_index(full_dof_number))
61  local_uncondensed_dofs_set.insert(full_dof_number);
62  else
63  nonlocal_uncondensed_dofs[_dof_map.dof_owner(full_dof_number)].insert(full_dof_number);
64 
65  elem_uncondensed_dofs.push_back(full_dof_number);
66  (uncondensed_global_to_local_map)[full_dof_number] = uncondensed_local_dof_number++;
67  if (involved_in_constraints)
68  constraint_dofs.insert(full_dof_number);
69 }
70 
72  const dof_id_type full_dof_number,
73  bool involved_in_constraints,
74  std::unordered_map<dof_id_type, dof_id_type> & uncondensed_global_to_local_map,
75  std::unordered_set<dof_id_type> & local_uncondensed_dofs_set,
76  std::unordered_map<processor_id_type, std::unordered_set<dof_id_type>> &
77  nonlocal_uncondensed_dofs,
78  std::vector<dof_id_type> & elem_uncondensed_dofs,
79  dof_id_type & uncondensed_local_dof_number,
80  std::unordered_set<dof_id_type> & constraint_dofs)
81 {
82  const auto & full_dof_constraints = _dof_map.get_dof_constraints();
83  auto it = full_dof_constraints.find(full_dof_number);
84  const bool is_constrained = it != full_dof_constraints.end();
85  involved_in_constraints = involved_in_constraints || is_constrained;
86 
87  this->add_uncondensed_dof(full_dof_number,
88  involved_in_constraints,
89  uncondensed_global_to_local_map,
90  local_uncondensed_dofs_set,
91  nonlocal_uncondensed_dofs,
92  elem_uncondensed_dofs,
93  uncondensed_local_dof_number,
94  constraint_dofs);
95  if (is_constrained)
96  for (const auto [full_constraining_dof, weight] : it->second)
97  {
99  // Our constraining dofs may themselves be constrained
100  this->add_uncondensed_dof_plus_constraint_dofs(full_constraining_dof,
101  /*involved_in_constraints=*/true,
102  uncondensed_global_to_local_map,
103  local_uncondensed_dofs_set,
104  nonlocal_uncondensed_dofs,
105  elem_uncondensed_dofs,
106  uncondensed_local_dof_number,
107  constraint_dofs);
108  }
109 }
110 
112 {
113  if (this->initialized())
114  this->clear();
115 
116  std::vector<dof_id_type> elem_dofs, elem_uncondensed_dofs; // only used to satisfy API
117  dof_id_type condensed_local_dof_number = 0, uncondensed_local_dof_number = 0;
118  std::unordered_map<dof_id_type, dof_id_type> *condensed_global_to_local_map = nullptr,
119  *uncondensed_global_to_local_map = nullptr;
120  std::set<unsigned int> full_vars_present_in_reduced_sys;
121  std::unordered_set<dof_id_type> local_uncondensed_dofs_set, constraint_dofs;
122  std::unordered_map<processor_id_type, std::unordered_set<dof_id_type>> nonlocal_uncondensed_dofs;
123 
124  // Handle SCALAR dofs
125  for (const auto vg : make_range(_dof_map.n_variable_groups()))
126  if (const auto & vg_description = _dof_map.variable_group(vg);
127  vg_description.type().family == SCALAR)
128  {
129  std::vector<dof_id_type> scalar_dof_indices;
130  const processor_id_type last_pid = this->comm().size() - 1;
131  for (const auto vg_vn : make_range(vg_description.n_variables()))
132  {
133  const auto vn = vg_description.number(vg_vn);
134  _dof_map.SCALAR_dof_indices(scalar_dof_indices, vn);
135  if (this->comm().rank() == last_pid)
136  local_uncondensed_dofs_set.insert(scalar_dof_indices.begin(),
137  scalar_dof_indices.end());
138  else
139  nonlocal_uncondensed_dofs[last_pid].insert(scalar_dof_indices.begin(),
140  scalar_dof_indices.end());
141  }
142  }
143 
144  auto scalar_dofs_functor =
145  [this,
146  &uncondensed_global_to_local_map,
147  &local_uncondensed_dofs_set,
148  &nonlocal_uncondensed_dofs,
149  &elem_uncondensed_dofs,
150  &uncondensed_local_dof_number,
151  &constraint_dofs](const Elem & /*elem*/,
152  std::vector<dof_id_type> & dof_indices,
153  const std::vector<dof_id_type> & scalar_dof_indices) {
154  dof_indices.insert(dof_indices.end(), scalar_dof_indices.begin(), scalar_dof_indices.end());
155  for (const auto global_dof : scalar_dof_indices)
157  false,
158  *uncondensed_global_to_local_map,
159  local_uncondensed_dofs_set,
160  nonlocal_uncondensed_dofs,
161  elem_uncondensed_dofs,
162  uncondensed_local_dof_number,
163  constraint_dofs);
164  };
165 
166  auto field_dofs_functor = [this,
167  &condensed_local_dof_number,
168  &condensed_global_to_local_map,
169  &uncondensed_global_to_local_map,
170  &local_uncondensed_dofs_set,
171  &nonlocal_uncondensed_dofs,
172  &elem_uncondensed_dofs,
173  &uncondensed_local_dof_number,
174  &constraint_dofs](const Elem & elem,
175  const unsigned int node_num,
176  const unsigned int var_num,
177  std::vector<dof_id_type> & dof_indices,
178  const dof_id_type field_dof) {
179  dof_indices.push_back(field_dof);
180 
181  bool uncondensed_dof = false;
182  if (_uncondensed_vars.count(var_num))
183  {
184  libmesh_assert_msg(
185  node_num == invalid_uint,
186  "Users should not be providing continuous FEM variables to the uncondensed vars API");
187  uncondensed_dof = true;
188  }
189 
190  if (node_num != invalid_uint && !elem.is_internal(node_num))
191  uncondensed_dof = true;
192 
193  if (uncondensed_dof)
195  false,
196  *uncondensed_global_to_local_map,
197  local_uncondensed_dofs_set,
198  nonlocal_uncondensed_dofs,
199  elem_uncondensed_dofs,
200  uncondensed_local_dof_number,
201  constraint_dofs);
202  else
203  (*condensed_global_to_local_map)[field_dof] = condensed_local_dof_number++;
204  };
205 
206  for (auto elem : _mesh.active_local_element_ptr_range())
207  {
208  auto & dof_data = _elem_to_dof_data[elem->id()];
209  condensed_local_dof_number = 0;
210  uncondensed_local_dof_number = 0;
211  condensed_global_to_local_map = &dof_data.condensed_global_to_local_map;
212  uncondensed_global_to_local_map = &dof_data.uncondensed_global_to_local_map;
213 
214  const auto sub_id = elem->subdomain_id();
215  for (const auto vg : make_range(_dof_map.n_variable_groups()))
216  {
217  const auto & var_group = _dof_map.variable_group(vg);
218  if (!var_group.active_on_subdomain(sub_id))
219  continue;
220 
221  for (const auto v : make_range(var_group.n_variables()))
222  {
223  const auto var_num = var_group.number(v);
224  dof_data.reduced_space_indices.resize(var_num + 1);
225  elem_uncondensed_dofs.clear();
226  _dof_map.dof_indices(elem,
227  elem_dofs,
228  var_num,
229  scalar_dofs_functor,
230  field_dofs_functor,
231  elem->p_level());
232  if (!elem_uncondensed_dofs.empty())
233  {
234  auto & var_reduced_space_indices = dof_data.reduced_space_indices[var_num];
235  var_reduced_space_indices.insert(var_reduced_space_indices.end(),
236  elem_uncondensed_dofs.begin(),
237  elem_uncondensed_dofs.end());
238  full_vars_present_in_reduced_sys.insert(var_num);
239  }
240  }
241  }
242  }
243 
244  _local_uncondensed_dofs.assign(local_uncondensed_dofs_set.begin(),
245  local_uncondensed_dofs_set.end());
246  local_uncondensed_dofs_set.clear();
247 
248  //
249  // Build the reduced system data
250  //
251 
252  const dof_id_type n_local = _local_uncondensed_dofs.size();
253  dof_id_type n = n_local;
254  this->comm().sum(n);
255 
256  // Get DOF counts on all processors
257  this->compute_dof_info(n_local);
258 
259  // Build a map from the full size problem uncondensed dof indices to the reduced problem
260  // (uncondensed) dof indices
261  std::unordered_map<dof_id_type, dof_id_type> full_dof_to_reduced_dof;
262  const auto local_start = _first_df[this->processor_id()];
263  for (const auto i : index_range(_local_uncondensed_dofs))
264  full_dof_to_reduced_dof[_local_uncondensed_dofs[i]] = i + local_start;
265 
266  // Build the condensed system sparsity pattern
268  this->_mesh, /*calculate_constrained=*/false, /*use_condensed_system=*/true);
269  const auto & nnz = _reduced_sp->get_n_nz();
270  const auto & noz = _reduced_sp->get_n_oz();
271  libmesh_assert(nnz.size() == noz.size());
272 
273  // Optimization for PETSc. This is critical for problems in which there are SCALAR dofs that
274  // introduce dense rows to avoid allocating a dense matrix
275  _reduced_nnz.resize(_local_uncondensed_dofs.size());
276  _reduced_noz.resize(_local_uncondensed_dofs.size());
277  for (const dof_id_type local_reduced_i : index_range(_local_uncondensed_dofs))
278  {
279  const dof_id_type full_i = _local_uncondensed_dofs[local_reduced_i];
280  const dof_id_type local_full_i = full_i - _dof_map.first_dof();
281  libmesh_assert(local_full_i < nnz.size());
282  _reduced_nnz[local_reduced_i] = nnz[local_full_i];
283  _reduced_noz[local_reduced_i] = noz[local_full_i];
284  }
285 
286  //
287  // Now we need to pull our nonlocal data
288  //
289 
290  // build our queries
291  std::unordered_map<processor_id_type, std::vector<dof_id_type>> nonlocal_uncondensed_dofs_mapvec;
292  for (const auto & [pid, set] : nonlocal_uncondensed_dofs)
293  {
294  auto & vec = nonlocal_uncondensed_dofs_mapvec[pid];
295  vec.assign(set.begin(), set.end());
296  }
297  // clear no longer needed memory
298  nonlocal_uncondensed_dofs.clear();
299 
300  auto gather_functor = [&full_dof_to_reduced_dof](processor_id_type,
301  const std::vector<dof_id_type> & full_dof_ids,
302  std::vector<dof_id_type> & reduced_dof_ids) {
303  reduced_dof_ids.resize(full_dof_ids.size());
304  for (const auto i : index_range(full_dof_ids))
305  reduced_dof_ids[i] = libmesh_map_find(full_dof_to_reduced_dof, full_dof_ids[i]);
306  };
307 
308  auto action_functor =
309  [&full_dof_to_reduced_dof](processor_id_type,
310  const std::vector<dof_id_type> & full_dof_ids,
311  const std::vector<dof_id_type> & reduced_dof_ids) {
312  for (const auto i : index_range(full_dof_ids))
313  {
314  libmesh_assert(!full_dof_to_reduced_dof.count(full_dof_ids[i]));
315  full_dof_to_reduced_dof[full_dof_ids[i]] = reduced_dof_ids[i];
316  }
317  };
318 
320  nonlocal_uncondensed_dofs_mapvec,
321  gather_functor,
322  action_functor,
324  nonlocal_uncondensed_dofs_mapvec.clear();
325 
326  // Determine the variables with any degrees of freedom present in the reduced system
327  _communicator.set_union(full_vars_present_in_reduced_sys);
328  _reduced_vars.reserve(full_vars_present_in_reduced_sys.size());
329  unsigned int first_local_number = 0;
330  for (const auto i : index_range(full_vars_present_in_reduced_sys))
331  {
332  const auto full_var_num = *std::next(full_vars_present_in_reduced_sys.begin(), i);
333  const auto & full_var = _dof_map.variable(full_var_num);
334  _reduced_vars.push_back(Variable{nullptr,
335  full_var.name(),
336  cast_int<unsigned int>(i),
337  first_local_number,
338  full_var.type()});
339  first_local_number += _reduced_vars.back().n_components(_mesh);
340  }
341 
342  // Now we can finally set our element reduced dof indices
343  std::vector<dof_id_type> var_full_dof_indices;
344  for (auto & [elem, dof_data] : _elem_to_dof_data)
345  {
346  libmesh_ignore(elem);
347  auto & reduced_space_indices = dof_data.reduced_space_indices;
348  // Start erasing from the back to reduce the number of copy assignment operations
349  if (reduced_space_indices.size())
350  for (typename std::vector<dof_id_type>::difference_type i =
351  reduced_space_indices.size() - 1;
352  i >= 0;
353  --i)
354  if (!full_vars_present_in_reduced_sys.count(i))
355  reduced_space_indices.erase(reduced_space_indices.begin() + i);
356  // It is theoretically possible that we have an element that doesn't have dofs for one of the
357  // variables present in our reduced system, which is why the assertion below is not an
358  // equality assertion
359  libmesh_assert(reduced_space_indices.size() <= full_vars_present_in_reduced_sys.size());
360 
361  for (auto & var_dof_indices : reduced_space_indices)
362  {
363  var_full_dof_indices = var_dof_indices;
364  var_dof_indices.clear();
365  for (const auto full_dof : var_full_dof_indices)
366  var_dof_indices.push_back(libmesh_map_find(full_dof_to_reduced_dof, full_dof));
367  }
368  }
369 
370  // Build our dof constraints map
371  for (const auto full_dof : constraint_dofs)
373  libmesh_map_find(full_dof_to_reduced_dof, full_dof);
374  constraint_dofs.clear();
375 
376  // Prevent querying Nodes for dof indices
377  std::vector<unsigned int> nvpg(_reduced_vars.size());
378  for (auto & elem : nvpg)
379  elem = 1;
380 
381  // add_system returns a system if it already exists instead of erroring so there's no harm if
382  // we do this multiple times
383  _reduced_system = &_system.get_equation_systems().add_system("Basic", "reduced");
385  {
387  for (auto * const nd : _mesh.active_node_ptr_range())
388  {
389  nd->set_n_vars_per_group(_reduced_system->number(), nvpg);
390  for (const auto g : index_range(nvpg))
391  nd->set_n_comp_group(_reduced_system->number(), g, 0);
392  }
393  }
394 
395  _sc_is_initialized = true;
396 }
397 
398 unsigned int StaticCondensationDofMap::n_variables() const { return _reduced_vars.size(); }
399 
400 const Variable & StaticCondensationDofMap::variable(const unsigned int c) const
401 {
402  return _reduced_vars[c];
403 }
404 
406  std::vector<dof_id_type> & di,
407  const unsigned int vn,
408  int /*p_level*/) const
409 {
410  di.clear();
411  di = libmesh_map_find(_elem_to_dof_data, elem->id()).reduced_space_indices[vn];
412 }
413 
415  std::vector<dof_id_type> &,
416  const unsigned int) const
417 {
418  libmesh_error_msg("StaticCondensationDofMap dof indices are only meant to be queried with "
419  "elements, not nodes");
420 }
421 
423 {
425  _elem_to_dof_data.clear();
426  _uncondensed_vars.clear();
427  _reduced_vars.clear();
428  _reduced_sp = nullptr;
429  _reduced_nnz.clear();
430  _reduced_noz.clear();
431  _sc_is_initialized = false;
432 }
433 }
virtual void clear()
Definition: dof_map_base.C:71
bool is_initialized() const
Definition: system.h:2414
FEFamily family
The type of finite element.
Definition: fe_type.h:221
void pull_parallel_vector_data(const Communicator &comm, const MapToVectors &queries, GatherFunctor &gather_data, const ActionFunctor &act_on_data, const datum *example)
unsigned int n_variable_groups() const
Definition: dof_map.h:625
A Node is like a Point, but with more information.
Definition: node.h:52
std::unordered_set< unsigned int > _uncondensed_vars
Variables for which we will keep all dofs.
virtual const Variable & variable(const unsigned int c) const override
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
processor_id_type dof_owner(const dof_id_type dof) const
Definition: dof_map.h:707
std::unordered_map< dof_id_type, dof_id_type > _full_to_reduced_constraint_dofs
A small map from full system degrees of freedom to reduced/condensed system degrees of freedom involv...
std::vector< dof_id_type > _reduced_nnz
Number of on-diagonal nonzeros per row in the reduced system.
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Definition: dof_map.C:2164
std::vector< dof_id_type > _first_df
First DOF index on processor p.
Definition: dof_map_base.h:154
std::vector< dof_id_type > _reduced_noz
Number of off-diagonal nonzeros per row in the reduced system.
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
const EquationSystems & get_equation_systems() const
Definition: system.h:721
void sum(T &r) const
std::vector< Variable > _reduced_vars
The variables in the reduced system.
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
MeshBase & mesh
const Parallel::Communicator & comm() const
const Parallel::Communicator & _communicator
The libMesh namespace provides an interface to certain functionality in the library.
void init()
Initializes degrees of freedom on the current mesh.
Definition: system.C:197
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:2547
uint8_t processor_id_type
Definition: id_types.h:104
This is the MeshBase class.
Definition: mesh_base.h:75
System * _reduced_system
A dummyish system to help with DofObjects.
const Variable & variable(const unsigned int c) const override
Definition: dof_map.h:2190
std::unordered_map< dof_id_type, DofData > _elem_to_dof_data
A map from element ID to Schur complement data.
virtual unsigned int n_variables() const override
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
processor_id_type size() const
void add_uncondensed_dof(dof_id_type full_dof_number, bool involved_in_constraints, std::unordered_map< dof_id_type, dof_id_type > &uncondensed_global_to_local_map, std::unordered_set< dof_id_type > &local_uncondensed_dofs_set, std::unordered_map< processor_id_type, std::unordered_set< dof_id_type >> &nonlocal_uncondensed_dofs, std::vector< dof_id_type > &elem_uncondensed_dofs, dof_id_type &uncondensed_local_dof_number, std::unordered_set< dof_id_type > &constraint_dofs)
Add an uncondensed dof.
void libmesh_ignore(const Args &...)
unsigned int number() const
Definition: system.h:2350
This class defines the notion of a variable in the system.
Definition: variable.h:50
StaticCondensationDofMap(MeshBase &mesh, System &system, const DofMap &dof_map)
dof_id_type id() const
Definition: dof_object.h:828
bool _sc_is_initialized
Whether our object has been initialized.
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:437
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
std::vector< dof_id_type > _local_uncondensed_dofs
All the uncondensed degrees of freedom (numbered in the "full" uncondensed + condensed space)...
This base class provides a minimal set of interfaces for satisfying user requests for...
Definition: dof_map_base.h:51
libmesh_assert(ctx)
static const dof_id_type invalid_id
An invalid id to distinguish an uninitialized DofObject.
Definition: dof_object.h:482
const VariableGroup & variable_group(const unsigned int c) const
Definition: dof_map.h:2180
void reinit()
Build the element global to local index maps.
static const unsigned int next[3]
A lookup table for the increment modulo 3 operation, for iterating through the three nodes per elemen...
const DofConstraints & get_dof_constraints() const
Provide a const accessor to the DofConstraints map.
Definition: dof_map.h:1049
virtual void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di, const unsigned int vn, int p_level=-12345) const override
Fills the vector di with the global degree of freedom indices for the element.
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
std::size_t compute_dof_info(dof_id_type n_local_dofs)
compute the key degree of freedom information given the local number of degrees of freedom on this pr...
unsigned int number(unsigned int v) const
Definition: variable.h:304
std::unique_ptr< SparsityPattern::Build > _reduced_sp
Owned storage of the reduced system sparsity pattern.
const std::string & name() const
Definition: variable.h:122
bool initialized() const
Whether we are initialized.
void add_uncondensed_dof_plus_constraint_dofs(dof_id_type full_dof_number, bool involved_in_constraints, std::unordered_map< dof_id_type, dof_id_type > &uncondensed_global_to_local_map, std::unordered_set< dof_id_type > &local_uncondensed_dofs_set, std::unordered_map< processor_id_type, std::unordered_set< dof_id_type >> &nonlocal_uncondensed_dofs, std::vector< dof_id_type > &elem_uncondensed_dofs, dof_id_type &uncondensed_local_dof_number, std::unordered_set< dof_id_type > &constraint_dofs)
Add an uncondensed dof potentially along with constraining dofs which themselves must/will also be un...
dof_id_type first_dof(const processor_id_type proc) const
Definition: dof_map_base.h:185
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
processor_id_type processor_id() const
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117
uint8_t dof_id_type
Definition: id_types.h:67
bool local_index(dof_id_type dof_index) const
Definition: dof_map.h:839
const FEType & type() const
Definition: variable.h:144
void set_union(T &data, const unsigned int root_id) const