1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2020 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 #ifndef DTKADAPTER_H
21 #define DTKADAPTER_H
22 
23 #include "libmesh/libmesh_config.h"
24 
25 #ifdef LIBMESH_TRILINOS_HAVE_DTK
26 
27 #include "libmesh/dtk_evaluator.h"
28 
29 #include "libmesh/ignore_warnings.h"
30 #include <DTK_MeshManager.hpp>
31 #include <DTK_MeshContainer.hpp>
32 #include <DTK_MeshTraits.hpp>
33 #include <DTK_MeshTraitsFieldAdapter.hpp>
34 #include <DTK_FieldManager.hpp>
35 #include <DTK_FieldContainer.hpp>
36 #include <DTK_FieldEvaluator.hpp>
37 
38 #include <Teuchos_RCP.hpp>
39 #include <Teuchos_ArrayRCP.hpp>
40 #include "libmesh/restore_warnings.h"
41 
42 namespace libMesh
43 {
44 
45 /**
46  * The DTKAdapter is used with the DTKSolutionTransfer object to adapt
47  * libmesh data to the DTK interface.
48  *
49  * \author Derek Gaston
50  * \date 2013
51  */
52 class DTKAdapter
53 {
54 public:
55   DTKAdapter(Teuchos::RCP<const Teuchos::Comm<int>> in_comm, EquationSystems & in_es);
56 
57   typedef DataTransferKit::MeshContainer<int>                                  MeshContainerType;
58   typedef DataTransferKit::FieldContainer<double>                              FieldContainerType;
59   typedef DataTransferKit::MeshTraits<MeshContainerType>::global_ordinal_type  GlobalOrdinal;
60   typedef DataTransferKit::FieldEvaluator<GlobalOrdinal,FieldContainerType>    EvaluatorType;
61   typedef Teuchos::RCP<EvaluatorType>                                          RCP_Evaluator;
62 
63 
get_mesh_manager()64   Teuchos::RCP<DataTransferKit::MeshManager<MeshContainerType>> get_mesh_manager() { return mesh_manager; }
65   RCP_Evaluator get_variable_evaluator(std::string var_name);
get_target_coords()66   Teuchos::RCP<DataTransferKit::FieldManager<MeshContainerType>> get_target_coords() { return target_coords; }
67   Teuchos::RCP<DataTransferKit::FieldManager<FieldContainerType>> get_values_to_fill(std::string var_name);
68 
69   /**
70    * After computing values for a variable in this EquationSystems
71    * we need to take those values and put them in the actual solution vector.
72    */
73   void update_variable_values(std::string var_name);
74 
75 protected:
76   /**
77    * Small helper function for finding the system containing the variable.
78    *
79    * \note This implies that variable names are unique across all systems!
80    */
81   System * find_sys(std::string var_name);
82 
83   /**
84    * \returns The DTK ElementTopology for a given Elem.
85    */
86   DataTransferKit::DTK_ElementTopology get_element_topology(const Elem * elem);
87 
88   /**
89    * Helper function that fills the std::set with all of the node numbers of
90    * nodes connected to local elements.
91    */
92   void get_semi_local_nodes(std::set<unsigned int> & semi_local_nodes);
93 
94   Teuchos::RCP<const Teuchos::Comm<int>> comm;
95   EquationSystems & es;
96   const MeshBase & mesh;
97   unsigned int dim;
98 
99   unsigned int num_local_nodes;
100   Teuchos::ArrayRCP<int> vertices;
101 
102   Teuchos::RCP<DataTransferKit::MeshManager<MeshContainerType>> mesh_manager;
103   RCP_Evaluator field_evaluator;
104   Teuchos::RCP<DataTransferKit::FieldManager<MeshContainerType>> target_coords;
105 
106   /// Map of variable names to arrays to be filled by a transfer
107   std::map<std::string, Teuchos::RCP<DataTransferKit::FieldManager<FieldContainerType>>> values_to_fill;
108 
109   /// Map of variable names to RCP_Evaluator objects
110   std::map<std::string, RCP_Evaluator> evaluators;
111 };
112 
113 } // namespace libMesh
114 
115 #endif // #ifdef LIBMESH_TRILINOS_HAVE_DTK
116 
117 #endif // #define DTKADAPTER_H
118