Go to the documentation of this file.
   20 #include "libmesh/libmesh_config.h" 
   22 #ifdef LIBMESH_TRILINOS_HAVE_DTK 
   24 #include "libmesh/dtk_solution_transfer.h" 
   25 #include "libmesh/system.h" 
   27 #include "libmesh/ignore_warnings.h" 
   30 #include <Teuchos_RCP.hpp> 
   31 #include <Teuchos_GlobalMPISession.hpp> 
   32 #include <Teuchos_Ptr.hpp> 
   33 #include <Teuchos_DefaultMpiComm.hpp> 
   34 #include <Teuchos_OpaqueWrapper.hpp> 
   37 #include <DTK_MeshManager.hpp> 
   38 #include <DTK_MeshContainer.hpp> 
   39 #include <DTK_FieldManager.hpp> 
   40 #include <DTK_FieldContainer.hpp> 
   41 #include <DTK_FieldTools.hpp> 
   42 #include <DTK_CommTools.hpp> 
   43 #include <DTK_CommIndexer.hpp> 
   45 #include "libmesh/restore_warnings.h" 
   54   comm_default = Teuchos::rcp(
new Teuchos::MpiComm<int>(Teuchos::rcp(
new Teuchos::OpaqueWrapper<MPI_Comm>(
comm.get()))));
 
   70   libmesh_experimental();
 
   86   std::pair<EquationSystems *, EquationSystems *> from_to(from_es, to_es);
 
  101   Teuchos::RCP<DataTransferKit::FieldManager<DTKAdapter::FieldContainerType>> to_values = to_adapter->
get_values_to_fill(to_var.
name());
 
  103   dtk_maps[from_to]->apply(from_evaluator, to_values);
 
  105   if (
dtk_maps[from_to]->getMissedTargetPoints().size())
 
  106     libMesh::out<<
"Warning: Some points were missed in the transfer of "<<from_var.
name()<<
" to "<<to_var.
name()<<
"!"<<std::endl;
 
  113 #endif // #ifdef LIBMESH_TRILINOS_HAVE_DTK 
  
const EquationSystems & get_equation_systems() const
 
const MeshBase & get_mesh() const
 
Base class for objects that allow transferring variable values between different systems with differe...
 
The libMesh namespace provides an interface to certain functionality in the library.
 
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...
 
const Parallel::Communicator & comm() const
 
unsigned int mesh_dimension() const
 
Teuchos::RCP< DataTransferKit::MeshManager< MeshContainerType > > get_mesh_manager()
 
std::map< EquationSystems *, DTKAdapter * > adapters
The DTKAdapter associated with each EquationSystems.
 
std::map< std::pair< EquationSystems *, EquationSystems * >, shared_domain_map_type * > dtk_maps
The dtk shared domain maps for pairs of EquationSystems (from, to)
 
const std::string & name() const
 
DTKSolutionTransfer(const libMesh::Parallel::Communicator &comm)
 
virtual ~DTKSolutionTransfer()
 
This class defines the notion of a variable in the system.
 
This is the EquationSystems class.
 
virtual void transfer(const Variable &from_var, const Variable &to_var) override
Transfer the values of a variable to another.
 
Teuchos::RCP< EvaluatorType > RCP_Evaluator
 
Teuchos::RCP< const Teuchos::Comm< int > > comm_default
COMM_WORLD for now.
 
Teuchos::RCP< DataTransferKit::FieldManager< MeshContainerType > > get_target_coords()
 
DataTransferKit::SharedDomainMap< DTKAdapter::MeshContainerType, DTKAdapter::MeshContainerType > shared_domain_map_type
 
Teuchos::RCP< DataTransferKit::FieldManager< FieldContainerType > > get_values_to_fill(std::string var_name)
 
RCP_Evaluator get_variable_evaluator(std::string var_name)
 
The DTKAdapter is used with the DTKSolutionTransfer object to adapt libmesh data to the DTK interface...