libMesh
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
libMesh::DTKAdapter Class Reference

The DTKAdapter is used with the DTKSolutionTransfer object to adapt libmesh data to the DTK interface. More...

#include <dtk_adapter.h>

Public Types

typedef DataTransferKit::MeshContainer< intMeshContainerType
 
typedef DataTransferKit::FieldContainer< double > FieldContainerType
 
typedef DataTransferKit::MeshTraits< MeshContainerType >::global_ordinal_type GlobalOrdinal
 
typedef DataTransferKit::FieldEvaluator< GlobalOrdinal, FieldContainerTypeEvaluatorType
 
typedef Teuchos::RCP< EvaluatorTypeRCP_Evaluator
 

Public Member Functions

 DTKAdapter (Teuchos::RCP< const Teuchos::Comm< int >> in_comm, EquationSystems &in_es)
 
Teuchos::RCP< DataTransferKit::MeshManager< MeshContainerType > > get_mesh_manager ()
 
RCP_Evaluator get_variable_evaluator (std::string var_name)
 
Teuchos::RCP< DataTransferKit::FieldManager< MeshContainerType > > get_target_coords ()
 
Teuchos::RCP< DataTransferKit::FieldManager< FieldContainerType > > get_values_to_fill (std::string var_name)
 
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 them in the actual solution vector. More...
 

Protected Member Functions

Systemfind_sys (std::string var_name)
 Small helper function for finding the system containing the variable. More...
 
DataTransferKit::DTK_ElementTopology get_element_topology (const Elem *elem)
 
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 elements. More...
 

Protected Attributes

Teuchos::RCP< const Teuchos::Comm< int > > comm
 
EquationSystemses
 
const MeshBasemesh
 
unsigned int dim
 
unsigned int num_local_nodes
 
Teuchos::ArrayRCP< intvertices
 
Teuchos::RCP< DataTransferKit::MeshManager< MeshContainerType > > mesh_manager
 
RCP_Evaluator field_evaluator
 
Teuchos::RCP< DataTransferKit::FieldManager< MeshContainerType > > target_coords
 
std::map< std::string, Teuchos::RCP< DataTransferKit::FieldManager< FieldContainerType > > > values_to_fill
 Map of variable names to arrays to be filled by a transfer. More...
 
std::map< std::string, RCP_Evaluatorevaluators
 Map of variable names to RCP_Evaluator objects. More...
 

Detailed Description

The DTKAdapter is used with the DTKSolutionTransfer object to adapt libmesh data to the DTK interface.

Author
Derek Gaston
Date
2013

Definition at line 52 of file dtk_adapter.h.

Member Typedef Documentation

◆ EvaluatorType

typedef DataTransferKit::FieldEvaluator<GlobalOrdinal,FieldContainerType> libMesh::DTKAdapter::EvaluatorType

Definition at line 60 of file dtk_adapter.h.

◆ FieldContainerType

typedef DataTransferKit::FieldContainer<double> libMesh::DTKAdapter::FieldContainerType

Definition at line 58 of file dtk_adapter.h.

◆ GlobalOrdinal

typedef DataTransferKit::MeshTraits<MeshContainerType>::global_ordinal_type libMesh::DTKAdapter::GlobalOrdinal

Definition at line 59 of file dtk_adapter.h.

◆ MeshContainerType

typedef DataTransferKit::MeshContainer<int> libMesh::DTKAdapter::MeshContainerType

Definition at line 57 of file dtk_adapter.h.

◆ RCP_Evaluator

Definition at line 61 of file dtk_adapter.h.

Constructor & Destructor Documentation

◆ DTKAdapter()

libMesh::DTKAdapter::DTKAdapter ( Teuchos::RCP< const Teuchos::Comm< int >>  in_comm,
EquationSystems in_es 
)

Definition at line 45 of file dtk_adapter.C.

45  :
46  comm(in_comm),
47  es(in_es),
48  mesh(in_es.get_mesh()),
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 }

References libMesh::as_range(), comm, dim, libMesh::MeshBase::elem_ptr(), get_element_topology(), get_semi_local_nodes(), libMesh::DofObject::id(), libMesh::Utility::iota(), libMesh::MeshBase::local_elements_begin(), libMesh::MeshBase::local_elements_end(), mesh, mesh_manager, libMesh::MeshBase::n_local_elem(), libMesh::Elem::n_nodes(), libMesh::MeshBase::node_ref(), num_local_nodes, target_coords, and vertices.

Member Function Documentation

◆ find_sys()

System * libMesh::DTKAdapter::find_sys ( std::string  var_name)
protected

Small helper function for finding the system containing the variable.

Note
This implies that variable names are unique across all systems!

Note that this implies that variable names are unique across all systems!

Definition at line 236 of file dtk_adapter.C.

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 }

References es, libMesh::EquationSystems::get_system(), libMesh::libmesh_assert(), and libMesh::EquationSystems::n_systems().

Referenced by get_variable_evaluator(), and update_variable_values().

◆ get_element_topology()

DataTransferKit::DTK_ElementTopology libMesh::DTKAdapter::get_element_topology ( const Elem elem)
protected
Returns
The DTK ElementTopology for a given Elem.

Definition at line 256 of file dtk_adapter.C.

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 }

References libMesh::EDGE2, libMesh::HEX8, libMesh::PYRAMID5, libMesh::QUAD4, libMesh::TET4, libMesh::TRI3, and libMesh::Elem::type().

Referenced by DTKAdapter().

◆ get_mesh_manager()

Teuchos::RCP<DataTransferKit::MeshManager<MeshContainerType> > libMesh::DTKAdapter::get_mesh_manager ( )
inline

Definition at line 64 of file dtk_adapter.h.

64 { return mesh_manager; }

References mesh_manager.

Referenced by libMesh::DTKSolutionTransfer::transfer().

◆ get_semi_local_nodes()

void libMesh::DTKAdapter::get_semi_local_nodes ( std::set< unsigned int > &  semi_local_nodes)
protected

Helper function that fills the std::set with all of the node numbers of nodes connected to local elements.

Definition at line 277 of file dtk_adapter.C.

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 }

References libMesh::as_range(), libMesh::MeshBase::local_elements_begin(), libMesh::MeshBase::local_elements_end(), and mesh.

Referenced by DTKAdapter().

◆ get_target_coords()

Teuchos::RCP<DataTransferKit::FieldManager<MeshContainerType> > libMesh::DTKAdapter::get_target_coords ( )
inline

Definition at line 66 of file dtk_adapter.h.

66 { return target_coords; }

References target_coords.

Referenced by libMesh::DTKSolutionTransfer::transfer().

◆ get_values_to_fill()

Teuchos::RCP< DataTransferKit::FieldManager< DTKAdapter::FieldContainerType > > libMesh::DTKAdapter::get_values_to_fill ( std::string  var_name)

Definition at line 189 of file dtk_adapter.C.

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 }

References comm, num_local_nodes, and values_to_fill.

Referenced by libMesh::DTKSolutionTransfer::transfer().

◆ get_variable_evaluator()

DTKAdapter::RCP_Evaluator libMesh::DTKAdapter::get_variable_evaluator ( std::string  var_name)

Definition at line 175 of file dtk_adapter.C.

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 }

References evaluators, and find_sys().

Referenced by libMesh::DTKSolutionTransfer::transfer().

◆ update_variable_values()

void libMesh::DTKAdapter::update_variable_values ( std::string  var_name)

After computing values for a variable in this EquationSystems we need to take those values and put them in the actual solution vector.

Definition at line 202 of file dtk_adapter.C.

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 }

References libMesh::DofObject::dof_number(), find_sys(), mesh, libMesh::MeshBase::node_ref(), libMesh::System::number(), libMesh::ParallelObject::processor_id(), libMesh::DofObject::processor_id(), libMesh::System::solution, value, values_to_fill, libMesh::System::variable_number(), and vertices.

Referenced by libMesh::DTKSolutionTransfer::transfer().

Member Data Documentation

◆ comm

Teuchos::RCP<const Teuchos::Comm<int> > libMesh::DTKAdapter::comm
protected

Definition at line 94 of file dtk_adapter.h.

Referenced by DTKAdapter(), and get_values_to_fill().

◆ dim

unsigned int libMesh::DTKAdapter::dim
protected

Definition at line 97 of file dtk_adapter.h.

Referenced by DTKAdapter().

◆ es

EquationSystems& libMesh::DTKAdapter::es
protected

Definition at line 95 of file dtk_adapter.h.

Referenced by find_sys().

◆ evaluators

std::map<std::string, RCP_Evaluator> libMesh::DTKAdapter::evaluators
protected

Map of variable names to RCP_Evaluator objects.

Definition at line 110 of file dtk_adapter.h.

Referenced by get_variable_evaluator().

◆ field_evaluator

RCP_Evaluator libMesh::DTKAdapter::field_evaluator
protected

Definition at line 103 of file dtk_adapter.h.

◆ mesh

const MeshBase& libMesh::DTKAdapter::mesh
protected

Definition at line 96 of file dtk_adapter.h.

Referenced by DTKAdapter(), get_semi_local_nodes(), and update_variable_values().

◆ mesh_manager

Teuchos::RCP<DataTransferKit::MeshManager<MeshContainerType> > libMesh::DTKAdapter::mesh_manager
protected

Definition at line 102 of file dtk_adapter.h.

Referenced by DTKAdapter(), and get_mesh_manager().

◆ num_local_nodes

unsigned int libMesh::DTKAdapter::num_local_nodes
protected

Definition at line 99 of file dtk_adapter.h.

Referenced by DTKAdapter(), and get_values_to_fill().

◆ target_coords

Teuchos::RCP<DataTransferKit::FieldManager<MeshContainerType> > libMesh::DTKAdapter::target_coords
protected

Definition at line 104 of file dtk_adapter.h.

Referenced by DTKAdapter(), and get_target_coords().

◆ values_to_fill

std::map<std::string, Teuchos::RCP<DataTransferKit::FieldManager<FieldContainerType> > > libMesh::DTKAdapter::values_to_fill
protected

Map of variable names to arrays to be filled by a transfer.

Definition at line 107 of file dtk_adapter.h.

Referenced by get_values_to_fill(), and update_variable_values().

◆ vertices

Teuchos::ArrayRCP<int> libMesh::DTKAdapter::vertices
protected

Definition at line 100 of file dtk_adapter.h.

Referenced by DTKAdapter(), and update_variable_values().


The documentation for this class was generated from the following files:
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::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::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::DTKAdapter::mesh_manager
Teuchos::RCP< DataTransferKit::MeshManager< MeshContainerType > > mesh_manager
Definition: dtk_adapter.h:102
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
libMesh::TET4
Definition: enum_elem_type.h:45
libMesh::DTKAdapter::dim
unsigned int dim
Definition: dtk_adapter.h:97
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::QUAD4
Definition: enum_elem_type.h:41
libMesh::TRI3
Definition: enum_elem_type.h:39
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::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::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::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::EDGE2
Definition: enum_elem_type.h:35
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33