libMesh
dtk_adapter.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 
19 
20 #include "libmesh/libmesh_config.h"
21 
22 #ifdef LIBMESH_TRILINOS_HAVE_DTK
23 
24 // libMesh includes
25 #include "libmesh/dtk_adapter.h"
26 #include "libmesh/dtk_evaluator.h"
27 #include "libmesh/mesh.h"
28 #include "libmesh/numeric_vector.h"
29 #include "libmesh/elem.h"
30 #include "libmesh/equation_systems.h"
31 
32 // Trilinos includes
33 #include "libmesh/ignore_warnings.h"
34 #include <DTK_MeshTypes.hpp>
35 #include <Teuchos_Comm.hpp>
36 #include "libmesh/restore_warnings.h"
37 
38 // C++ includes
39 #include <vector>
40 #include <numeric>
41 
42 namespace libMesh
43 {
44 
45 DTKAdapter::DTKAdapter(Teuchos::RCP<const Teuchos::Comm<int>> in_comm, EquationSystems & in_es):
46  comm(in_comm),
47  es(in_es),
48  mesh(in_es.get_mesh()),
49  dim(mesh.mesh_dimension())
50 {
51  std::set<unsigned int> semi_local_nodes;
52  get_semi_local_nodes(semi_local_nodes);
53 
54  num_local_nodes = semi_local_nodes.size();
55 
56  vertices.resize(num_local_nodes);
57  Teuchos::ArrayRCP<double> coordinates(num_local_nodes * dim);
58 
59  // Fill in the vertices and coordinates
60  {
61  unsigned int i = 0;
62 
63  for (const auto & id : semi_local_nodes)
64  {
65  const Node & node = mesh.node_ref(id);
66 
67  vertices[i] = node.id();
68 
69  for (unsigned int j=0; j<dim; j++)
70  coordinates[(j*num_local_nodes) + i] = node(j);
71 
72  i++;
73  }
74  }
75 
76  // Currently assuming all elements are the same!
77  DataTransferKit::DTK_ElementTopology element_topology = get_element_topology(mesh.elem_ptr(0));
78  unsigned int n_nodes_per_elem = mesh.elem_ptr(0)->n_nodes();
79 
80  unsigned int n_local_elem = mesh.n_local_elem();
81 
82  Teuchos::ArrayRCP<int> elements(n_local_elem);
83  Teuchos::ArrayRCP<int> connectivity(n_nodes_per_elem*n_local_elem);
84 
85  // Fill in the elements and connectivity
86  {
87  unsigned int i = 0;
88 
89  for (const auto & elem : as_range(mesh.local_elements_begin(),
91  {
92  elements[i] = elem->id();
93 
94  for (unsigned int j=0; j<n_nodes_per_elem; j++)
95  connectivity[(j*n_local_elem)+i] = elem->node_id(j);
96 
97  i++;
98  }
99  }
100 
101  Teuchos::ArrayRCP<int> permutation_list(n_nodes_per_elem);
102  std::iota(permutation_list.begin(), permutation_list.end(), 0);
103 
104  /*
105  if (this->processor_id() == 1)
106  sleep(1);
107 
108  libMesh::out<<"n_nodes_per_elem: "<<n_nodes_per_elem<<std::endl;
109 
110  libMesh::out<<"Dim: "<<dim<<std::endl;
111 
112  libMesh::err<<"Vertices size: "<<vertices.size()<<std::endl;
113  {
114  libMesh::err<<this->processor_id()<<" Vertices: ";
115 
116  for (std::size_t i=0; i<vertices.size(); i++)
117  libMesh::err<<vertices[i]<<" ";
118 
119  libMesh::err<<std::endl;
120  }
121 
122  libMesh::err<<"Coordinates size: "<<coordinates.size()<<std::endl;
123  {
124  libMesh::err<<this->processor_id()<<" Coordinates: ";
125 
126  for (std::size_t i=0; i<coordinates.size(); i++)
127  libMesh::err<<coordinates[i]<<" ";
128 
129  libMesh::err<<std::endl;
130  }
131 
132  libMesh::err<<"Connectivity size: "<<connectivity.size()<<std::endl;
133  {
134  libMesh::err<<this->processor_id()<<" Connectivity: ";
135 
136  for (std::size_t i=0; i<connectivity.size(); i++)
137  libMesh::err<<connectivity[i]<<" ";
138 
139  libMesh::err<<std::endl;
140  }
141 
142  libMesh::err<<"Permutation_List size: "<<permutation_list.size()<<std::endl;
143  {
144  libMesh::err<<this->processor_id()<<" Permutation_List: ";
145 
146  for (std::size_t i=0; i<permutation_list.size(); i++)
147  libMesh::err<<permutation_list[i]<<" ";
148 
149  libMesh::err<<std::endl;
150  }
151 
152  */
153  Teuchos::RCP<MeshContainerType>
154  mesh_container = Teuchos::rcp(new MeshContainerType(dim,
155  vertices,
156  coordinates,
157  element_topology,
158  n_nodes_per_elem,
159  elements,
160  connectivity,
161  permutation_list));
162 
163  // We only have 1 element topology in this grid so we make just one mesh block
164  Teuchos::ArrayRCP<Teuchos::RCP<MeshContainerType>> mesh_blocks(1);
165  mesh_blocks[0] = mesh_container;
166 
167  // Create the MeshManager
168  mesh_manager = Teuchos::rcp(new DataTransferKit::MeshManager<MeshContainerType>(mesh_blocks, comm, dim) );
169 
170  // Pack the coordinates into a field, this will be the positions we'll ask for other systems fields at
171  target_coords = Teuchos::rcp(new DataTransferKit::FieldManager<MeshContainerType>(mesh_container, comm));
172 }
173 
176 {
177  if (evaluators.find(var_name) == evaluators.end()) // We haven't created an evaluator for the variable yet
178  {
179  System * sys = find_sys(var_name);
180 
181  // Create the FieldEvaluator
182  evaluators[var_name] = Teuchos::rcp(new DTKEvaluator(*sys, var_name));
183  }
184 
185  return evaluators[var_name];
186 }
187 
188 Teuchos::RCP<DataTransferKit::FieldManager<DTKAdapter::FieldContainerType>>
189 DTKAdapter::get_values_to_fill(std::string var_name)
190 {
191  if (values_to_fill.find(var_name) == values_to_fill.end())
192  {
193  Teuchos::ArrayRCP<double> data_space(num_local_nodes);
194  Teuchos::RCP<FieldContainerType> field_container = Teuchos::rcp(new FieldContainerType(data_space, 1));
195  values_to_fill[var_name] = Teuchos::rcp(new DataTransferKit::FieldManager<FieldContainerType>(field_container, comm));
196  }
197 
198  return values_to_fill[var_name];
199 }
200 
201 void
203 {
204  System * sys = find_sys(var_name);
205  unsigned int var_num = sys->variable_number(var_name);
206 
207  Teuchos::RCP<FieldContainerType> values = values_to_fill[var_name]->field();
208 
209  unsigned int i=0;
210  // Loop over the values (one for each node) and assign the value of this variable at each node
211  for (const auto & value : *values)
212  {
213  unsigned int node_num = vertices[i];
214  const Node & node = mesh.node_ref(node_num);
215 
216  if (node.processor_id() == sys->processor_id())
217  {
218  // The 0 is for the component... this only works for LAGRANGE!
219  dof_id_type dof = node.dof_number(sys->number(), var_num, 0);
220  sys->solution->set(dof, value);
221  }
222 
223  i++;
224  }
225 
226  sys->solution->close();
227 }
228 
229 
235 System *
236 DTKAdapter::find_sys(std::string var_name)
237 {
238  System * sys = nullptr;
239 
240  // Find the system this variable is from
241  for (unsigned int i=0; i<es.n_systems(); i++)
242  {
243  if (es.get_system(i).has_variable(var_name))
244  {
245  sys = &es.get_system(i);
246  break;
247  }
248  }
249 
250  libmesh_assert(sys);
251 
252  return sys;
253 }
254 
255 DataTransferKit::DTK_ElementTopology
257 {
258  ElemType type = elem->type();
259 
260  if (type == EDGE2)
261  return DataTransferKit::DTK_LINE_SEGMENT;
262  else if (type == TRI3)
263  return DataTransferKit::DTK_TRIANGLE;
264  else if (type == QUAD4)
265  return DataTransferKit::DTK_QUADRILATERAL;
266  else if (type == TET4)
267  return DataTransferKit::DTK_TETRAHEDRON;
268  else if (type == HEX8)
269  return DataTransferKit::DTK_HEXAHEDRON;
270  else if (type == PYRAMID5)
271  return DataTransferKit::DTK_PYRAMID;
272 
273  libmesh_error_msg("Element type not supported by DTK!");
274 }
275 
276 void
277 DTKAdapter::get_semi_local_nodes(std::set<unsigned int> & semi_local_nodes)
278 {
279  for (const auto & elem : as_range(mesh.local_elements_begin(), mesh.local_elements_end()))
280  for (const Node & node : elem->node_ref_range())
281  semi_local_nodes.insert(node.id());
282 }
283 
284 } // namespace libMesh
285 
286 #endif // #ifdef LIBMESH_TRILINOS_HAVE_DTK
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::DTKAdapter::evaluators
std::map< std::string, RCP_Evaluator > evaluators
Map of variable names to RCP_Evaluator objects.
Definition: dtk_adapter.h:110
libMesh::Elem::n_nodes
virtual unsigned int n_nodes() const =0
libMesh::HEX8
Definition: enum_elem_type.h:47
libMesh::DTKAdapter::get_element_topology
DataTransferKit::DTK_ElementTopology get_element_topology(const Elem *elem)
Definition: dtk_adapter.C:256
libMesh::DTKAdapter::values_to_fill
std::map< std::string, Teuchos::RCP< DataTransferKit::FieldManager< FieldContainerType > > > values_to_fill
Map of variable names to arrays to be filled by a transfer.
Definition: dtk_adapter.h:107
libMesh::DTKAdapter::num_local_nodes
unsigned int num_local_nodes
Definition: dtk_adapter.h:99
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::DTKAdapter::es
EquationSystems & es
Definition: dtk_adapter.h:95
libMesh::DTKAdapter::FieldContainerType
DataTransferKit::FieldContainer< double > FieldContainerType
Definition: dtk_adapter.h:58
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
libMesh::DTKAdapter::update_variable_values
void update_variable_values(std::string var_name)
After computing values for a variable in this EquationSystems we need to take those values and put th...
Definition: dtk_adapter.C:202
libMesh::System::number
unsigned int number() const
Definition: system.h:2075
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::DTKAdapter::mesh_manager
Teuchos::RCP< DataTransferKit::MeshManager< MeshContainerType > > mesh_manager
Definition: dtk_adapter.h:102
libMesh::DofObject::processor_id
processor_id_type processor_id() const
Definition: dof_object.h:829
libMesh::MeshBase::elem_ptr
virtual const Elem * elem_ptr(const dof_id_type i) const =0
libMesh::DTKAdapter::comm
Teuchos::RCP< const Teuchos::Comm< int > > comm
Definition: dtk_adapter.h:94
libMesh::MeshBase::local_elements_begin
virtual element_iterator local_elements_begin()=0
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::TET4
Definition: enum_elem_type.h:45
libMesh::DTKAdapter::dim
unsigned int dim
Definition: dtk_adapter.h:97
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::libmesh_assert
libmesh_assert(ctx)
libMesh::DTKAdapter::vertices
Teuchos::ArrayRCP< int > vertices
Definition: dtk_adapter.h:100
libMesh::Utility::iota
void iota(ForwardIter first, ForwardIter last, T value)
Utility::iota is a duplication of the SGI STL extension std::iota.
Definition: utility.h:105
libMesh::ParallelObject::processor_id
processor_id_type processor_id() const
Definition: parallel_object.h:106
libMesh::QUAD4
Definition: enum_elem_type.h:41
libMesh::TRI3
Definition: enum_elem_type.h:39
libMesh::Node
A Node is like a Point, but with more information.
Definition: node.h:52
libMesh::DTKAdapter::mesh
const MeshBase & mesh
Definition: dtk_adapter.h:96
libMesh::as_range
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
Helper function that allows us to treat a homogenous pair as a range.
Definition: simple_range.h:57
libMesh::MeshBase::n_local_elem
dof_id_type n_local_elem() const
Definition: mesh_base.h:403
libMesh::EquationSystems::n_systems
unsigned int n_systems() const
Definition: equation_systems.h:652
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::DTKAdapter::target_coords
Teuchos::RCP< DataTransferKit::FieldManager< MeshContainerType > > target_coords
Definition: dtk_adapter.h:104
libMesh::MeshBase::node_ref
virtual const Node & node_ref(const dof_id_type i) const
Definition: mesh_base.h:451
libMesh::DTKAdapter::RCP_Evaluator
Teuchos::RCP< EvaluatorType > RCP_Evaluator
Definition: dtk_adapter.h:61
libMesh::System::solution
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1539
libMesh::MeshBase::local_elements_end
virtual element_iterator local_elements_end()=0
value
static const bool value
Definition: xdr_io.C:56
libMesh::PYRAMID5
Definition: enum_elem_type.h:53
libMesh::DTKAdapter::find_sys
System * find_sys(std::string var_name)
Small helper function for finding the system containing the variable.
Definition: dtk_adapter.C:236
libMesh::DofObject::id
dof_id_type id() const
Definition: dof_object.h:767
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::System::variable_number
unsigned short int variable_number(const std::string &var) const
Definition: system.C:1232
libMesh::DTKAdapter::get_values_to_fill
Teuchos::RCP< DataTransferKit::FieldManager< FieldContainerType > > get_values_to_fill(std::string var_name)
Definition: dtk_adapter.C:189
libMesh::DTKAdapter::get_variable_evaluator
RCP_Evaluator get_variable_evaluator(std::string var_name)
Definition: dtk_adapter.C:175
libMesh::DTKEvaluator
Implements the evaluate() function to compute FE solution values at points requested by DTK.
Definition: dtk_evaluator.h:60
libMesh::DTKAdapter::MeshContainerType
DataTransferKit::MeshContainer< int > MeshContainerType
Definition: dtk_adapter.h:57
libMesh::DTKAdapter::get_semi_local_nodes
void get_semi_local_nodes(std::set< unsigned int > &semi_local_nodes)
Helper function that fills the std::set with all of the node numbers of nodes connected to local elem...
Definition: dtk_adapter.C:277
libMesh::Elem::type
virtual ElemType type() const =0
libMesh::EDGE2
Definition: enum_elem_type.h:35
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
libMesh::DTKAdapter::DTKAdapter
DTKAdapter(Teuchos::RCP< const Teuchos::Comm< int >> in_comm, EquationSystems &in_es)
Definition: dtk_adapter.C:45