18 #include "libmesh/static_condensation_dof_map.h" 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" 27 #include <unordered_set> 38 _sc_is_initialized(false)
40 libmesh_experimental();
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,
51 nonlocal_uncondensed_dofs,
52 std::vector<dof_id_type> & elem_uncondensed_dofs,
54 std::unordered_set<dof_id_type> & constraint_dofs)
56 if (uncondensed_global_to_local_map.count(full_dof_number))
61 local_uncondensed_dofs_set.insert(full_dof_number);
63 nonlocal_uncondensed_dofs[
_dof_map.
dof_owner(full_dof_number)].insert(full_dof_number);
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);
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,
77 nonlocal_uncondensed_dofs,
78 std::vector<dof_id_type> & elem_uncondensed_dofs,
80 std::unordered_set<dof_id_type> & constraint_dofs)
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;
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,
96 for (
const auto [full_constraining_dof,
weight] : it->second)
102 uncondensed_global_to_local_map,
103 local_uncondensed_dofs_set,
104 nonlocal_uncondensed_dofs,
105 elem_uncondensed_dofs,
106 uncondensed_local_dof_number,
116 std::vector<dof_id_type> elem_dofs, elem_uncondensed_dofs;
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;
129 std::vector<dof_id_type> scalar_dof_indices;
131 for (
const auto vg_vn :
make_range(vg_description.n_variables()))
133 const auto vn = vg_description.number(vg_vn);
135 if (this->
comm().rank() == last_pid)
136 local_uncondensed_dofs_set.insert(scalar_dof_indices.begin(),
137 scalar_dof_indices.end());
139 nonlocal_uncondensed_dofs[last_pid].insert(scalar_dof_indices.begin(),
140 scalar_dof_indices.end());
144 auto scalar_dofs_functor =
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 & ,
153 const std::vector<dof_id_type> & scalar_dof_indices) {
155 for (
const auto global_dof : scalar_dof_indices)
158 *uncondensed_global_to_local_map,
159 local_uncondensed_dofs_set,
160 nonlocal_uncondensed_dofs,
161 elem_uncondensed_dofs,
162 uncondensed_local_dof_number,
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,
181 bool uncondensed_dof =
false;
186 "Users should not be providing continuous FEM variables to the uncondensed vars API");
187 uncondensed_dof =
true;
190 if (node_num !=
invalid_uint && !elem.is_internal(node_num))
191 uncondensed_dof =
true;
196 *uncondensed_global_to_local_map,
197 local_uncondensed_dofs_set,
198 nonlocal_uncondensed_dofs,
199 elem_uncondensed_dofs,
200 uncondensed_local_dof_number,
203 (*condensed_global_to_local_map)[field_dof] = condensed_local_dof_number++;
206 for (
auto elem :
_mesh.active_local_element_ptr_range())
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;
214 const auto sub_id = elem->subdomain_id();
218 if (!var_group.active_on_subdomain(sub_id))
221 for (
const auto v :
make_range(var_group.n_variables()))
223 const auto var_num = var_group.
number(v);
224 dof_data.reduced_space_indices.resize(var_num + 1);
225 elem_uncondensed_dofs.clear();
232 if (!elem_uncondensed_dofs.empty())
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);
245 local_uncondensed_dofs_set.end());
246 local_uncondensed_dofs_set.clear();
261 std::unordered_map<dof_id_type, dof_id_type> full_dof_to_reduced_dof;
268 this->
_mesh,
false,
true);
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)
294 auto & vec = nonlocal_uncondensed_dofs_mapvec[pid];
295 vec.assign(
set.begin(),
set.end());
298 nonlocal_uncondensed_dofs.clear();
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());
305 reduced_dof_ids[i] = libmesh_map_find(full_dof_to_reduced_dof, full_dof_ids[i]);
308 auto action_functor =
310 const std::vector<dof_id_type> & full_dof_ids,
311 const std::vector<dof_id_type> & reduced_dof_ids) {
315 full_dof_to_reduced_dof[full_dof_ids[i]] = reduced_dof_ids[i];
320 nonlocal_uncondensed_dofs_mapvec,
324 nonlocal_uncondensed_dofs_mapvec.clear();
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))
332 const auto full_var_num = *
std::next(full_vars_present_in_reduced_sys.begin(), i);
336 cast_int<unsigned int>(i),
343 std::vector<dof_id_type> var_full_dof_indices;
347 auto & reduced_space_indices = dof_data.reduced_space_indices;
349 if (reduced_space_indices.size())
350 for (
typename std::vector<dof_id_type>::difference_type i =
351 reduced_space_indices.size() - 1;
354 if (!full_vars_present_in_reduced_sys.count(i))
355 reduced_space_indices.erase(reduced_space_indices.begin() + i);
359 libmesh_assert(reduced_space_indices.size() <= full_vars_present_in_reduced_sys.size());
361 for (
auto & var_dof_indices : reduced_space_indices)
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));
371 for (
const auto full_dof : constraint_dofs)
373 libmesh_map_find(full_dof_to_reduced_dof, full_dof);
374 constraint_dofs.clear();
378 for (
auto & elem : nvpg)
387 for (
auto *
const nd :
_mesh.active_node_ptr_range())
406 std::vector<dof_id_type> & di,
407 const unsigned int vn,
415 std::vector<dof_id_type> &,
416 const unsigned int)
const 418 libmesh_error_msg(
"StaticCondensationDofMap dof indices are only meant to be queried with " 419 "elements, not nodes");
bool is_initialized() const
FEFamily family
The type of finite element.
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
A Node is like a Point, but with more information.
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...
processor_id_type dof_owner(const dof_id_type dof) const
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
std::vector< dof_id_type > _first_df
First DOF index on processor p.
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...
const EquationSystems & get_equation_systems() 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.
const Parallel::Communicator & comm() const
virtual ~StaticCondensationDofMap()
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.
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...
uint8_t processor_id_type
This is the MeshBase class.
System * _reduced_system
A dummyish system to help with DofObjects.
const Variable & variable(const unsigned int c) const override
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.
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
This class defines the notion of a variable in the system.
StaticCondensationDofMap(MeshBase &mesh, System &system, const DofMap &dof_map)
bool _sc_is_initialized
Whether our object has been initialized.
Manages consistently variables, degrees of freedom, and coefficient vectors.
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...
static const dof_id_type invalid_id
An invalid id to distinguish an uninitialized DofObject.
const VariableGroup & variable_group(const unsigned int c) const
void reinit()
Build the element global to local index maps.
const DofConstraints & get_dof_constraints() const
Provide a const accessor to the DofConstraints map.
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...
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
std::unique_ptr< SparsityPattern::Build > _reduced_sp
Owned storage of the reduced system sparsity pattern.
const std::string & name() const
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
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...
virtual void clear() override
bool local_index(dof_id_type dof_index) const
const FEType & type() const
void set_union(T &data, const unsigned int root_id) const