Go to the documentation of this file.
   20 #include "libmesh/libmesh_config.h" 
   22 #ifdef LIBMESH_TRILINOS_HAVE_DTK 
   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" 
   33 #include "libmesh/ignore_warnings.h" 
   34 #include <DTK_MeshTypes.hpp> 
   35 #include <Teuchos_Comm.hpp> 
   36 #include "libmesh/restore_warnings.h" 
   48   mesh(in_es.get_mesh()),
 
   51   std::set<unsigned int> semi_local_nodes;
 
   63     for (
const auto & 
id : semi_local_nodes)
 
   69         for (
unsigned int j=0; j<
dim; j++)
 
   82   Teuchos::ArrayRCP<int> elements(n_local_elem);
 
   83   Teuchos::ArrayRCP<int> connectivity(n_nodes_per_elem*n_local_elem);
 
   92         elements[i] = elem->id();
 
   94         for (
unsigned int j=0; j<n_nodes_per_elem; j++)
 
   95           connectivity[(j*n_local_elem)+i] = elem->node_id(j);
 
  101   Teuchos::ArrayRCP<int> permutation_list(n_nodes_per_elem);
 
  102   std::iota(permutation_list.begin(), permutation_list.end(), 0);
 
  153   Teuchos::RCP<MeshContainerType>
 
  164   Teuchos::ArrayRCP<Teuchos::RCP<MeshContainerType>> mesh_blocks(1);
 
  165   mesh_blocks[0] = mesh_container;
 
  168   mesh_manager = Teuchos::rcp(
new DataTransferKit::MeshManager<MeshContainerType>(mesh_blocks, 
comm, 
dim) );
 
  171   target_coords = Teuchos::rcp(
new DataTransferKit::FieldManager<MeshContainerType>(mesh_container, 
comm));
 
  188 Teuchos::RCP<DataTransferKit::FieldManager<DTKAdapter::FieldContainerType>>
 
  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));
 
  207   Teuchos::RCP<FieldContainerType> values = 
values_to_fill[var_name]->field();
 
  211   for (
const auto & 
value : *values)
 
  213       unsigned int node_num = 
vertices[i];
 
  255 DataTransferKit::DTK_ElementTopology
 
  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;
 
  271     return DataTransferKit::DTK_PYRAMID;
 
  273   libmesh_error_msg(
"Element type not supported by DTK!");
 
  280     for (
const Node & node : elem->node_ref_range())
 
  281       semi_local_nodes.insert(node.id());
 
  286 #endif // #ifdef LIBMESH_TRILINOS_HAVE_DTK 
  
Manages consistently variables, degrees of freedom, and coefficient vectors.
 
std::map< std::string, RCP_Evaluator > evaluators
Map of variable names to RCP_Evaluator objects.
 
virtual unsigned int n_nodes() const =0
 
DataTransferKit::DTK_ElementTopology get_element_topology(const Elem *elem)
 
std::map< std::string, Teuchos::RCP< DataTransferKit::FieldManager< FieldContainerType > > > values_to_fill
Map of variable names to arrays to be filled by a transfer.
 
unsigned int num_local_nodes
 
The libMesh namespace provides an interface to certain functionality in the library.
 
DataTransferKit::FieldContainer< double > FieldContainerType
 
const T_sys & get_system(const std::string &name) const
 
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...
 
unsigned int number() const
 
Teuchos::RCP< DataTransferKit::MeshManager< MeshContainerType > > mesh_manager
 
processor_id_type processor_id() const
 
virtual const Elem * elem_ptr(const dof_id_type i) const =0
 
Teuchos::RCP< const Teuchos::Comm< int > > comm
 
virtual element_iterator local_elements_begin()=0
 
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
 
Teuchos::ArrayRCP< int > vertices
 
void iota(ForwardIter first, ForwardIter last, T value)
Utility::iota is a duplication of the SGI STL extension std::iota.
 
processor_id_type processor_id() const
 
A Node is like a Point, but with more information.
 
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
Helper function that allows us to treat a homogenous pair as a range.
 
dof_id_type n_local_elem() const
 
unsigned int n_systems() const
 
This is the EquationSystems class.
 
Teuchos::RCP< DataTransferKit::FieldManager< MeshContainerType > > target_coords
 
virtual const Node & node_ref(const dof_id_type i) const
 
Teuchos::RCP< EvaluatorType > RCP_Evaluator
 
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
 
virtual element_iterator local_elements_end()=0
 
System * find_sys(std::string var_name)
Small helper function for finding the system containing the variable.
 
This is the base class from which all geometric element types are derived.
 
unsigned short int variable_number(const std::string &var) const
 
Teuchos::RCP< DataTransferKit::FieldManager< FieldContainerType > > get_values_to_fill(std::string var_name)
 
RCP_Evaluator get_variable_evaluator(std::string var_name)
 
Implements the evaluate() function to compute FE solution values at points requested by DTK.
 
DataTransferKit::MeshContainer< int > MeshContainerType
 
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...
 
virtual ElemType type() const =0
 
ElemType
Defines an enum for geometric element types.
 
DTKAdapter(Teuchos::RCP< const Teuchos::Comm< int >> in_comm, EquationSystems &in_es)