libMesh
boundary_volume_solution_transfer.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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/boundary_volume_solution_transfer.h"
19 #include "libmesh/elem.h"
20 #include "libmesh/numeric_vector.h"
21 #include "libmesh/dof_map.h"
22 
23 namespace libMesh {
24 
26  const Variable & to_var)
27 {
28  // Determine which direction the transfer is in
29  System * from_sys = from_var.system();
30  System * to_sys = to_var.system();
31 
32  unsigned int
33  from_dimension = from_sys->get_mesh().mesh_dimension(),
34  to_dimension = to_sys->get_mesh().mesh_dimension();
35 
36  // Sanity check
37  if (from_dimension == to_dimension)
38  libmesh_error_msg("Error: Transfer must be from volume mesh to its boundary or vice-versa!");
39 
40  if (from_dimension > to_dimension)
41  this->transfer_volume_boundary(from_var, to_var);
42  else
43  this->transfer_boundary_volume(from_var, to_var);
44 }
45 
46 
47 
49 transfer_volume_boundary(const Variable & from_var, const Variable & to_var)
50 {
51  // Get references to the Systems from the Variables
52  System * from_sys = from_var.system(); // volume system
53  System * to_sys = to_var.system(); // boundary system
54 
55  // Get reference to the BoundaryMesh. Note: we always loop over the
56  // BoundaryMesh since, by definition, it has fewer dofs than the
57  // volume mesh.
58  const MeshBase & to_mesh = to_sys->get_mesh();
59 
60  // Get system number and variable numbers
61  const unsigned int from_sys_number = from_sys->number();
62  const unsigned int from_var_number = from_var.number();
63 
64  const unsigned int to_sys_number = to_sys->number();
65  const unsigned int to_var_number = to_var.number();
66 
67  // Get a constant reference to variables, get their number of components
68  const unsigned int from_n_comp = from_var.n_components();
69  const unsigned int to_n_comp = to_var.n_components();
70 
71  // Sanity check that the variables have the same number of components
72  libmesh_assert_equal_to(from_n_comp, to_n_comp);
73 
74  // Construct map from "from" dofs to "to" dofs.
75  typedef std::map<numeric_index_type, numeric_index_type> DofMapping;
76  DofMapping dof_mapping;
77 
78  // Loop through all boundary elements.
79  for (const auto & to_elem : to_mesh.active_local_element_ptr_range())
80  {
81  const Elem * from_elem = to_elem->interior_parent();
82 
83  if (!from_elem)
84  libmesh_error_msg("Error, transfer must be between a Mesh and its associated BoundaryMesh.");
85 
86  // loop through all nodes in each boundary element.
87  for (const Node & to_node : to_elem->node_ref_range())
88  {
89  // Nodes in interior_parent element.
90  for (const Node & from_node : from_elem->node_ref_range())
91  {
92  const dof_id_type from_dof = from_node.dof_number(from_sys_number,
93  from_var_number,
94  from_n_comp - 1);
95 
96  // See if we've already encountered this DOF in the loop
97  // over boundary elements.
98  DofMapping::iterator it = dof_mapping.find(from_dof);
99 
100  // If we've already mapped this dof, we don't need to map
101  // it again or do floating point comparisons.
102  if (it == dof_mapping.end() &&
103  from_node.absolute_fuzzy_equals(to_node, TOLERANCE))
104  {
105  // Global dof_index for node in BoundaryMesh.
106  const dof_id_type to_dof = to_node.dof_number(to_sys_number,
107  to_var_number,
108  to_n_comp - 1);
109 
110  // Keep track of the volume system dof index which is needed.
111  dof_mapping[from_dof] = to_dof;
112  }
113  }
114  }
115  }
116 
117  // Construct a vector of the indices needed from the Volume system's
118  // global solution vector on this processor.
119  std::vector<numeric_index_type> needed_indices;
120  needed_indices.reserve(dof_mapping.size());
121  {
122  DofMapping::iterator
123  it = dof_mapping.begin(),
124  end = dof_mapping.end();
125 
126  for (; it!=end; ++it)
127  needed_indices.push_back(it->first);
128  }
129 
130  // Communicate just the required values without making a copy of the
131  // global solution vector.
132  std::vector<Number> needed_values;
133  from_sys->solution->localize(needed_values, needed_indices);
134 
135  // Loop over DofMapping again, assigning values in the
136  // Boundary System solution vector.
137  {
138  DofMapping::iterator
139  it = dof_mapping.begin(),
140  end = dof_mapping.end();
141 
142  for (unsigned idx=0; it!=end; ++it, ++idx)
143  to_sys->solution->set(it->second, needed_values[idx]);
144  }
145 }
146 
147 
149 transfer_boundary_volume(const Variable & from_var, const Variable & to_var)
150 {
151  // Get references to the Systems from the Variables
152  System * from_sys = from_var.system(); // boundary system
153  System * to_sys = to_var.system(); // volume system
154 
155  // Get reference to BoundaryMesh.
156  const MeshBase & from_mesh = from_sys->get_mesh();
157 
158  // DofMap for BoundaryMesh
159  const DofMap & dof_map = from_sys->get_dof_map();
160 
161  // Get system number and variable numbers
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();
165  const unsigned int to_n_comp = to_var.n_components();
166 
167  // In order to get solution vectors from BoundaryMesh
168  std::vector<dof_id_type> from_dof_indices;
169  std::vector<Number> value;
170 
171  // Loop through all boundary elements.
172  for (const auto & from_elem : from_mesh.active_local_element_ptr_range())
173  {
174  const Elem * to_elem = from_elem->interior_parent();
175 
176  if (!to_elem)
177  libmesh_error_msg("Error, transfer must be between a Mesh and its associated BoundaryMesh.");
178 
179  // Get dof indices for phi2 for all nodes on this boundary element
180  dof_map.dof_indices(from_elem, from_dof_indices, from_var_number);
181 
182  // Get values from BoundaryMesh for this element
183  from_sys->current_local_solution->get(from_dof_indices, value);
184 
185  // loop through all nodes in each boundary element.
186  for (auto node : from_elem->node_index_range())
187  {
188  // Node in boundary element.
189  const Node * from_node = from_elem->node_ptr(node);
190 
191  // Nodes in interior_parent element.
192  for (const Node & to_node : to_elem->node_ref_range())
193  {
194  // Match BoundaryNode & VolumeNode.
195  if (to_node.absolute_fuzzy_equals(*from_node, TOLERANCE))
196  {
197  // Global dof_index for node in VolumeMesh.
198  const dof_id_type to_dof = to_node.dof_number(to_sys_number,
199  to_var_number,
200  to_n_comp - 1);
201 
202  // Assign values to boundary in VolumeMesh.
203  to_sys->solution->set(to_dof, value[node]);
204  }
205  }
206  }
207  }
208 }
209 
210 } // namespace libMesh
libMesh::System
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:100
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::Variable::system
System * system() const
Definition: variable.h:92
libMesh::DofMap::dof_indices
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.
Definition: dof_map.C:1967
libMesh::BoundaryVolumeSolutionTransfer::transfer_volume_boundary
void transfer_volume_boundary(const Variable &from_var, const Variable &to_var)
Transfer values from volume mesh to boundary mesh.
Definition: boundary_volume_solution_transfer.C:49
libMesh::Variable::number
unsigned int number() const
Definition: variable.h:106
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::Variable::n_components
unsigned int n_components() const
Definition: variable.h:125
end
IterBase * end
Also have a polymorphic pointer to the end object, this prevents iterating past the end.
Definition: variant_filter_iterator.h:343
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::System::number
unsigned int number() const
Definition: system.h:2075
libMesh::BoundaryVolumeSolutionTransfer::transfer
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...
Definition: boundary_volume_solution_transfer.C:25
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::BoundaryVolumeSolutionTransfer::transfer_boundary_volume
void transfer_boundary_volume(const Variable &from_var, const Variable &to_var)
Transfer values from boundary mesh to volume mesh.
Definition: boundary_volume_solution_transfer.C:149
libMesh::DofObject::dof_number
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
Definition: dof_object.h:956
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::Elem::node_ref_range
SimpleRange< NodeRefIter > node_ref_range()
Returns a range with all nodes of an element, usable in range-based for loops.
Definition: elem.h:2152
libMesh::System::get_mesh
const MeshBase & get_mesh() const
Definition: system.h:2083
libMesh::System::current_local_solution
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
Definition: system.h:1551
libMesh::Variable
This class defines the notion of a variable in the system.
Definition: variable.h:49
libMesh::Node
A Node is like a Point, but with more information.
Definition: node.h:52
libMesh::MeshTools::Generation::Private::idx
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.
Definition: mesh_generation.C:72
libMesh::Elem::interior_parent
const Elem * interior_parent() const
Definition: elem.C:749
libMesh::System::solution
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1539
value
static const bool value
Definition: xdr_io.C:56
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099