libMesh
|
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< int > | MeshContainerType |
typedef DataTransferKit::FieldContainer< double > | FieldContainerType |
typedef DataTransferKit::MeshTraits< MeshContainerType >::global_ordinal_type | GlobalOrdinal |
typedef DataTransferKit::FieldEvaluator< GlobalOrdinal, FieldContainerType > | EvaluatorType |
typedef Teuchos::RCP< EvaluatorType > | RCP_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 | |
System * | find_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 |
EquationSystems & | es |
const MeshBase & | mesh |
unsigned int | dim |
unsigned int | num_local_nodes |
Teuchos::ArrayRCP< int > | vertices |
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_Evaluator > | evaluators |
Map of variable names to RCP_Evaluator objects. More... | |
The DTKAdapter is used with the DTKSolutionTransfer object to adapt libmesh data to the DTK interface.
Definition at line 52 of file dtk_adapter.h.
typedef DataTransferKit::FieldEvaluator<GlobalOrdinal,FieldContainerType> libMesh::DTKAdapter::EvaluatorType |
Definition at line 60 of file dtk_adapter.h.
typedef DataTransferKit::FieldContainer<double> libMesh::DTKAdapter::FieldContainerType |
Definition at line 58 of file dtk_adapter.h.
typedef DataTransferKit::MeshTraits<MeshContainerType>::global_ordinal_type libMesh::DTKAdapter::GlobalOrdinal |
Definition at line 59 of file dtk_adapter.h.
typedef DataTransferKit::MeshContainer<int> libMesh::DTKAdapter::MeshContainerType |
Definition at line 57 of file dtk_adapter.h.
typedef Teuchos::RCP<EvaluatorType> libMesh::DTKAdapter::RCP_Evaluator |
Definition at line 61 of file dtk_adapter.h.
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.
|
protected |
Small helper function for finding the system containing the variable.
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().
|
protected |
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().
|
inline |
Definition at line 64 of file dtk_adapter.h.
References mesh_manager.
Referenced by libMesh::DTKSolutionTransfer::transfer().
|
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().
|
inline |
Definition at line 66 of file dtk_adapter.h.
References target_coords.
Referenced by libMesh::DTKSolutionTransfer::transfer().
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().
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().
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().
|
protected |
Definition at line 94 of file dtk_adapter.h.
Referenced by DTKAdapter(), and get_values_to_fill().
|
protected |
Definition at line 97 of file dtk_adapter.h.
Referenced by DTKAdapter().
|
protected |
Definition at line 95 of file dtk_adapter.h.
Referenced by find_sys().
|
protected |
Map of variable names to RCP_Evaluator objects.
Definition at line 110 of file dtk_adapter.h.
Referenced by get_variable_evaluator().
|
protected |
Definition at line 103 of file dtk_adapter.h.
|
protected |
Definition at line 96 of file dtk_adapter.h.
Referenced by DTKAdapter(), get_semi_local_nodes(), and update_variable_values().
|
protected |
Definition at line 102 of file dtk_adapter.h.
Referenced by DTKAdapter(), and get_mesh_manager().
|
protected |
Definition at line 99 of file dtk_adapter.h.
Referenced by DTKAdapter(), and get_values_to_fill().
|
protected |
Definition at line 104 of file dtk_adapter.h.
Referenced by DTKAdapter(), and get_target_coords().
|
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().
|
protected |
Definition at line 100 of file dtk_adapter.h.
Referenced by DTKAdapter(), and update_variable_values().