Loading [MathJax]/extensions/tex2jax.js
libMesh
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
dtk_solution_transfer.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #include "libmesh/libmesh_config.h"
21 
22 #ifdef LIBMESH_TRILINOS_HAVE_DTK
23 
24 #include "libmesh/dtk_solution_transfer.h"
25 #include "libmesh/system.h"
26 
27 #include "libmesh/ignore_warnings.h"
28 
29 // Trilinos Includes
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>
35 
36 // DTK Includes
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>
44 
45 #include "libmesh/restore_warnings.h"
46 
47 namespace libMesh
48 {
49 
51  SolutionTransfer(comm)
52 {
53  //comm_default = Teuchos::DefaultComm<int>::getComm();
54  comm_default = Teuchos::rcp(new Teuchos::MpiComm<int>(Teuchos::rcp(new Teuchos::OpaqueWrapper<MPI_Comm>(comm.get()))));
55 }
56 
58 
59 void
61  const Variable & to_var)
62 {
63  libmesh_experimental();
64 
65  EquationSystems * from_es = &from_var.system()->get_equation_systems();
66  EquationSystems * to_es = &to_var.system()->get_equation_systems();
67 
68  // Possibly make an Adapter for from_es
69  if (!adapters.count(from_es))
70  adapters[from_es] = std::make_unique<DTKAdapter>(comm_default, *from_es);
71 
72  // Possibly make an Adapter for to_es
73  if (!adapters.count(to_es))
74  adapters[to_es] = std::make_unique<DTKAdapter>(comm_default, *to_es);
75 
76  DTKAdapter * from_adapter = adapters[from_es].get();
77  DTKAdapter * to_adapter = adapters[to_es].get();
78 
79  // We first try emplacing a nullptr into the "dtk_maps" map with the desired key
80  auto [it, emplaced] = dtk_maps.emplace(std::make_pair(from_es, to_es), nullptr);
81 
82  // If the emplace succeeded, that means it was a new entry in the
83  // map, so we need to actually construct the object.
84  if (emplaced)
85  {
86  libmesh_assert(from_es->get_mesh().mesh_dimension() == to_es->get_mesh().mesh_dimension());
87 
88  it->second = std::make_unique<shared_domain_map_type>(comm_default, from_es->get_mesh().mesh_dimension(), true);
89 
90  // The tolerance here is for the "contains_point()" implementation in DTK. Set a larger value for a looser tolerance...
91  it->second->setup(from_adapter->get_mesh_manager(), to_adapter->get_target_coords(), 30*Teuchos::ScalarTraits<double>::eps());
92  }
93 
94  DTKAdapter::RCP_Evaluator from_evaluator = from_adapter->get_variable_evaluator(from_var.name());
95  Teuchos::RCP<DataTransferKit::FieldManager<DTKAdapter::FieldContainerType>> to_values = to_adapter->get_values_to_fill(to_var.name());
96 
97  it->second->apply(from_evaluator, to_values);
98 
99  if (it->second->getMissedTargetPoints().size())
100  libMesh::out << "Warning: Some points were missed in the transfer of " << from_var.name() << " to " << to_var.name() << "!" << std::endl;
101 
102  to_adapter->update_variable_values(to_var.name());
103 }
104 
105 } // namespace libMesh
106 
107 #endif // #ifdef LIBMESH_TRILINOS_HAVE_DTK
This is the EquationSystems class.
Teuchos::RCP< DataTransferKit::FieldManager< MeshContainerType > > get_target_coords()
Definition: dtk_adapter.h:66
std::map< std::pair< EquationSystems *, EquationSystems * >, std::unique_ptr< shared_domain_map_type > > dtk_maps
The dtk shared domain maps for pairs of EquationSystems (from, to)
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...
Definition: dtk_adapter.C:211
const EquationSystems & get_equation_systems() const
Definition: system.h:730
DTKSolutionTransfer(const libMesh::Parallel::Communicator &comm)
const Parallel::Communicator & comm() const
The DTKAdapter is used with the DTKSolutionTransfer object to adapt libmesh data to the DTK interface...
Definition: dtk_adapter.h:52
std::map< EquationSystems *, std::unique_ptr< DTKAdapter > > adapters
The DTKAdapter associated with each EquationSystems.
The libMesh namespace provides an interface to certain functionality in the library.
This class defines the notion of a variable in the system.
Definition: variable.h:49
RCP_Evaluator get_variable_evaluator(std::string var_name)
Definition: dtk_adapter.C:176
libmesh_assert(ctx)
virtual void transfer(const Variable &from_var, const Variable &to_var) override
Transfer the values of a variable to another.
Teuchos::RCP< const Teuchos::Comm< int > > comm_default
COMM_WORLD for now.
communicator & get()
Teuchos::RCP< EvaluatorType > RCP_Evaluator
Definition: dtk_adapter.h:61
System * system() const
Definition: variable.h:113
Teuchos::RCP< DataTransferKit::MeshManager< MeshContainerType > > get_mesh_manager()
Definition: dtk_adapter.h:64
OStreamProxy out
const MeshBase & get_mesh() const
unsigned int mesh_dimension() const
Definition: mesh_base.C:354
Base class for objects that allow transferring variable values between different systems with differe...
const std::string & name() const
Definition: variable.h:121
Teuchos::RCP< DataTransferKit::FieldManager< FieldContainerType > > get_values_to_fill(std::string var_name)
Definition: dtk_adapter.C:193