Go to the documentation of this file.
18 #include "libmesh/boundary_volume_solution_transfer.h"
19 #include "libmesh/elem.h"
20 #include "libmesh/numeric_vector.h"
21 #include "libmesh/dof_map.h"
37 if (from_dimension == to_dimension)
38 libmesh_error_msg(
"Error: Transfer must be from volume mesh to its boundary or vice-versa!");
40 if (from_dimension > to_dimension)
61 const unsigned int from_sys_number = from_sys->
number();
62 const unsigned int from_var_number = from_var.
number();
64 const unsigned int to_sys_number = to_sys->
number();
65 const unsigned int to_var_number = to_var.
number();
68 const unsigned int from_n_comp = from_var.
n_components();
72 libmesh_assert_equal_to(from_n_comp, to_n_comp);
75 typedef std::map<numeric_index_type, numeric_index_type> DofMapping;
76 DofMapping dof_mapping;
84 libmesh_error_msg(
"Error, transfer must be between a Mesh and its associated BoundaryMesh.");
87 for (
const Node & to_node : to_elem->node_ref_range())
92 const dof_id_type from_dof = from_node.dof_number(from_sys_number,
98 DofMapping::iterator it = dof_mapping.find(from_dof);
102 if (it == dof_mapping.end() &&
103 from_node.absolute_fuzzy_equals(to_node,
TOLERANCE))
106 const dof_id_type to_dof = to_node.dof_number(to_sys_number,
111 dof_mapping[from_dof] = to_dof;
119 std::vector<numeric_index_type> needed_indices;
120 needed_indices.reserve(dof_mapping.size());
123 it = dof_mapping.begin(),
124 end = dof_mapping.end();
126 for (; it!=
end; ++it)
127 needed_indices.push_back(it->first);
132 std::vector<Number> needed_values;
133 from_sys->
solution->localize(needed_values, needed_indices);
139 it = dof_mapping.begin(),
140 end = dof_mapping.end();
142 for (
unsigned idx=0; it!=
end; ++it, ++
idx)
143 to_sys->
solution->set(it->second, needed_values[
idx]);
162 const unsigned int to_sys_number = to_sys->
number();
163 const unsigned int to_var_number = to_var.
number();
164 const unsigned int from_var_number = from_var.
number();
168 std::vector<dof_id_type> from_dof_indices;
169 std::vector<Number>
value;
177 libmesh_error_msg(
"Error, transfer must be between a Mesh and its associated BoundaryMesh.");
180 dof_map.
dof_indices(from_elem, from_dof_indices, from_var_number);
186 for (
auto node : from_elem->node_index_range())
189 const Node * from_node = from_elem->node_ptr(node);
195 if (to_node.absolute_fuzzy_equals(*from_node,
TOLERANCE))
Manages consistently variables, degrees of freedom, and coefficient vectors.
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
void transfer_volume_boundary(const Variable &from_var, const Variable &to_var)
Transfer values from volume mesh to boundary mesh.
unsigned int number() const
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
The libMesh namespace provides an interface to certain functionality in the library.
unsigned int n_components() const
IterBase * end
Also have a polymorphic pointer to the end object, this prevents iterating past the end.
static const Real TOLERANCE
unsigned int number() const
virtual void transfer(const Variable &from_var, const Variable &to_var) override
Transfer values from a Variable in a System associated with a volume mesh to a Variable in a System a...
unsigned int mesh_dimension() const
void transfer_boundary_volume(const Variable &from_var, const Variable &to_var)
Transfer values from boundary mesh to volume mesh.
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
This is the MeshBase class.
SimpleRange< NodeRefIter > node_ref_range()
Returns a range with all nodes of an element, usable in range-based for loops.
const MeshBase & get_mesh() const
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
This class defines the notion of a variable in the system.
A Node is like a Point, but with more information.
const Elem * interior_parent() const
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
This class handles the numbering of degrees of freedom on a mesh.
This is the base class from which all geometric element types are derived.
const DofMap & get_dof_map() const