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 46 of file dtk_adapter.C.

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

46  :
47  comm(in_comm),
48  es(in_es),
49  mesh(in_es.get_mesh()),
51 {
52  std::set<unsigned int> semi_local_nodes;
53  get_semi_local_nodes(semi_local_nodes);
54 
55  num_local_nodes = semi_local_nodes.size();
56 
57  vertices.resize(num_local_nodes);
58  Teuchos::ArrayRCP<double> coordinates(num_local_nodes * dim);
59 
60  // Fill in the vertices and coordinates
61  {
62  unsigned int i = 0;
63 
64  for (const auto & id : semi_local_nodes)
65  {
66  const Node & node = mesh.node_ref(id);
67 
68  vertices[i] = node.id();
69 
70  for (unsigned int j=0; j<dim; j++)
71  coordinates[(j*num_local_nodes) + i] = node(j);
72 
73  i++;
74  }
75  }
76 
77  // Currently assuming all elements are the same!
78  DataTransferKit::DTK_ElementTopology element_topology = get_element_topology(mesh.elem_ptr(0));
79  unsigned int n_nodes_per_elem = mesh.elem_ptr(0)->n_nodes();
80 
81  unsigned int n_local_elem = mesh.n_local_elem();
82 
83  Teuchos::ArrayRCP<int> elements(n_local_elem);
84  Teuchos::ArrayRCP<int> connectivity(n_nodes_per_elem*n_local_elem);
85 
86  // Fill in the elements and connectivity
87  {
88  unsigned int i = 0;
89 
90  for (const auto & elem : as_range(mesh.local_elements_begin(),
91  mesh.local_elements_end()))
92  {
93  elements[i] = elem->id();
94 
95  for (unsigned int j=0; j<n_nodes_per_elem; j++)
96  connectivity[(j*n_local_elem)+i] = elem->node_id(j);
97 
98  i++;
99  }
100  }
101 
102  Teuchos::ArrayRCP<int> permutation_list(n_nodes_per_elem);
103  std::iota(permutation_list.begin(), permutation_list.end(), 0);
104 
105  /*
106  if (this->processor_id() == 1)
107  sleep(1);
108 
109  libMesh::out<<"n_nodes_per_elem: "<<n_nodes_per_elem<<std::endl;
110 
111  libMesh::out<<"Dim: "<<dim<<std::endl;
112 
113  libMesh::err<<"Vertices size: "<<vertices.size()<<std::endl;
114  {
115  libMesh::err<<this->processor_id()<<" Vertices: ";
116 
117  for (std::size_t i=0; i<vertices.size(); i++)
118  libMesh::err<<vertices[i]<<" ";
119 
120  libMesh::err<<std::endl;
121  }
122 
123  libMesh::err<<"Coordinates size: "<<coordinates.size()<<std::endl;
124  {
125  libMesh::err<<this->processor_id()<<" Coordinates: ";
126 
127  for (std::size_t i=0; i<coordinates.size(); i++)
128  libMesh::err<<coordinates[i]<<" ";
129 
130  libMesh::err<<std::endl;
131  }
132 
133  libMesh::err<<"Connectivity size: "<<connectivity.size()<<std::endl;
134  {
135  libMesh::err<<this->processor_id()<<" Connectivity: ";
136 
137  for (std::size_t i=0; i<connectivity.size(); i++)
138  libMesh::err<<connectivity[i]<<" ";
139 
140  libMesh::err<<std::endl;
141  }
142 
143  libMesh::err<<"Permutation_List size: "<<permutation_list.size()<<std::endl;
144  {
145  libMesh::err<<this->processor_id()<<" Permutation_List: ";
146 
147  for (std::size_t i=0; i<permutation_list.size(); i++)
148  libMesh::err<<permutation_list[i]<<" ";
149 
150  libMesh::err<<std::endl;
151  }
152 
153  */
154  Teuchos::RCP<MeshContainerType>
155  mesh_container = Teuchos::rcp(new MeshContainerType(dim,
156  vertices,
157  coordinates,
158  element_topology,
159  n_nodes_per_elem,
160  elements,
161  connectivity,
162  permutation_list));
163 
164  // We only have 1 element topology in this grid so we make just one mesh block
165  Teuchos::ArrayRCP<Teuchos::RCP<MeshContainerType>> mesh_blocks(1);
166  mesh_blocks[0] = mesh_container;
167 
168  // Create the MeshManager
169  mesh_manager = Teuchos::rcp(new DataTransferKit::MeshManager<MeshContainerType>(mesh_blocks, comm, dim) );
170 
171  // Pack the coordinates into a field, this will be the positions we'll ask for other systems fields at
172  target_coords = Teuchos::rcp(new DataTransferKit::FieldManager<MeshContainerType>(mesh_container, comm));
173 }
Teuchos::RCP< const Teuchos::Comm< int > > comm
Definition: dtk_adapter.h:94
unsigned int dim
Definition: dtk_adapter.h:97
unsigned int num_local_nodes
Definition: dtk_adapter.h:99
DataTransferKit::MeshContainer< int > MeshContainerType
Definition: dtk_adapter.h:57
const MeshBase & mesh
Definition: dtk_adapter.h:96
dof_id_type n_local_elem() const
Definition: mesh_base.h:527
void iota(ForwardIter first, ForwardIter last, T value)
Utility::iota was created back when std::iota was just an SGI STL extension.
Definition: utility.h:229
DataTransferKit::DTK_ElementTopology get_element_topology(const Elem *elem)
Definition: dtk_adapter.C:257
EquationSystems & es
Definition: dtk_adapter.h:95
virtual unsigned int n_nodes() const =0
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
virtual const Elem * elem_ptr(const dof_id_type i) const =0
unsigned int mesh_dimension() const
Definition: mesh_base.C:324
Teuchos::RCP< DataTransferKit::FieldManager< MeshContainerType > > target_coords
Definition: dtk_adapter.h:104
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:278
virtual const Node & node_ref(const dof_id_type i) const
Definition: mesh_base.h:575
Teuchos::RCP< DataTransferKit::MeshManager< MeshContainerType > > mesh_manager
Definition: dtk_adapter.h:102
Teuchos::ArrayRCP< int > vertices
Definition: dtk_adapter.h:100

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 237 of file dtk_adapter.C.

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

Referenced by get_variable_evaluator(), and update_variable_values().

238 {
239  System * sys = nullptr;
240 
241  // Find the system this variable is from
242  for (auto i : make_range(es.n_systems()))
243  {
244  if (es.get_system(i).has_variable(var_name))
245  {
246  sys = &es.get_system(i);
247  break;
248  }
249  }
250 
251  libmesh_assert(sys);
252 
253  return sys;
254 }
unsigned int n_systems() const
const T_sys & get_system(std::string_view name) const
EquationSystems & es
Definition: dtk_adapter.h:95
libmesh_assert(ctx)
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...
Definition: int_range.h:134

◆ 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 257 of file dtk_adapter.C.

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

Referenced by DTKAdapter().

258 {
259  ElemType type = elem->type();
260 
261  if (type == EDGE2)
262  return DataTransferKit::DTK_LINE_SEGMENT;
263  else if (type == TRI3)
264  return DataTransferKit::DTK_TRIANGLE;
265  else if (type == QUAD4)
266  return DataTransferKit::DTK_QUADRILATERAL;
267  else if (type == TET4)
268  return DataTransferKit::DTK_TETRAHEDRON;
269  else if (type == HEX8)
270  return DataTransferKit::DTK_HEXAHEDRON;
271  else if (type == PYRAMID5)
272  return DataTransferKit::DTK_PYRAMID;
273 
274  libmesh_error_msg("Element type not supported by DTK!");
275 }
ElemType
Defines an enum for geometric element types.

◆ get_mesh_manager()

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

Definition at line 64 of file dtk_adapter.h.

References mesh_manager.

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

64 { return mesh_manager; }
Teuchos::RCP< DataTransferKit::MeshManager< MeshContainerType > > mesh_manager
Definition: dtk_adapter.h:102

◆ 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 278 of file dtk_adapter.C.

References libMesh::as_range(), and mesh.

Referenced by DTKAdapter().

279 {
280  for (const auto & elem : as_range(mesh.local_elements_begin(), mesh.local_elements_end()))
281  for (const Node & node : elem->node_ref_range())
282  semi_local_nodes.insert(node.id());
283 }
const MeshBase & mesh
Definition: dtk_adapter.h:96
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

◆ get_target_coords()

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

Definition at line 66 of file dtk_adapter.h.

References target_coords.

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

66 { return target_coords; }
Teuchos::RCP< DataTransferKit::FieldManager< MeshContainerType > > target_coords
Definition: dtk_adapter.h:104

◆ get_values_to_fill()

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

Definition at line 190 of file dtk_adapter.C.

References comm, num_local_nodes, and values_to_fill.

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

191 {
192  if (values_to_fill.find(var_name) == values_to_fill.end())
193  {
194  Teuchos::ArrayRCP<double> data_space(num_local_nodes);
195  Teuchos::RCP<FieldContainerType> field_container = Teuchos::rcp(new FieldContainerType(data_space, 1));
196  values_to_fill[var_name] = Teuchos::rcp(new DataTransferKit::FieldManager<FieldContainerType>(field_container, comm));
197  }
198 
199  return values_to_fill[var_name];
200 }
Teuchos::RCP< const Teuchos::Comm< int > > comm
Definition: dtk_adapter.h:94
unsigned int num_local_nodes
Definition: dtk_adapter.h:99
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
DataTransferKit::FieldContainer< double > FieldContainerType
Definition: dtk_adapter.h:58

◆ get_variable_evaluator()

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

Definition at line 176 of file dtk_adapter.C.

References evaluators, and find_sys().

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

177 {
178  if (evaluators.find(var_name) == evaluators.end()) // We haven't created an evaluator for the variable yet
179  {
180  System * sys = find_sys(var_name);
181 
182  // Create the FieldEvaluator
183  evaluators[var_name] = Teuchos::rcp(new DTKEvaluator(*sys, var_name));
184  }
185 
186  return evaluators[var_name];
187 }
std::map< std::string, RCP_Evaluator > evaluators
Map of variable names to RCP_Evaluator objects.
Definition: dtk_adapter.h:110
System * find_sys(std::string var_name)
Small helper function for finding the system containing the variable.
Definition: dtk_adapter.C:237

◆ 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 203 of file dtk_adapter.C.

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().

204 {
205  System * sys = find_sys(var_name);
206  unsigned int var_num = sys->variable_number(var_name);
207 
208  Teuchos::RCP<FieldContainerType> values = values_to_fill[var_name]->field();
209 
210  unsigned int i=0;
211  // Loop over the values (one for each node) and assign the value of this variable at each node
212  for (const auto & value : *values)
213  {
214  unsigned int node_num = vertices[i];
215  const Node & node = mesh.node_ref(node_num);
216 
217  if (node.processor_id() == sys->processor_id())
218  {
219  // The 0 is for the component... this only works for LAGRANGE!
220  dof_id_type dof = node.dof_number(sys->number(), var_num, 0);
221  sys->solution->set(dof, value);
222  }
223 
224  i++;
225  }
226 
227  sys->solution->close();
228 }
const MeshBase & mesh
Definition: dtk_adapter.h:96
static const bool value
Definition: xdr_io.C:54
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
virtual const Node & node_ref(const dof_id_type i) const
Definition: mesh_base.h:575
Teuchos::ArrayRCP< int > vertices
Definition: dtk_adapter.h:100
System * find_sys(std::string var_name)
Small helper function for finding the system containing the variable.
Definition: dtk_adapter.C:237
uint8_t dof_id_type
Definition: id_types.h:67

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: