ExaDG
Loading...
Searching...
No Matches
dof_coupling.h
1
2/* ______________________________________________________________________
3 *
4 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
5 *
6 * Copyright (C) 2022 by the ExaDG authors
7 *
8 * This program is free software: you can redistribute it and/or modify
9 * it under the terms of the GNU General Public License as published by
10 * the Free Software Foundation, either version 3 of the License, or
11 * (at your option) any later version.
12 *
13 * This program is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the GNU General Public License
19 * along with this program. If not, see <https://www.gnu.org/licenses/>.
20 * ______________________________________________________________________
21 */
22
23#ifndef INCLUDE_EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_DOF_COUPLING_H_
24#define INCLUDE_EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_DOF_COUPLING_H_
25
26#include <deal.II/dofs/dof_tools.h>
27
28#include <exadg/fluid_structure_interaction/precice/coupling_base.h>
29#include <exadg/fluid_structure_interaction/precice/dof_tools_extension.h>
30
31namespace ExaDG
32{
33namespace preCICE
34{
42template<int dim, int data_dim, typename VectorizedArrayType>
43class DoFCoupling : public CouplingBase<dim, data_dim, VectorizedArrayType>
44{
45public:
46 DoFCoupling(dealii::MatrixFree<dim, double, VectorizedArrayType> const & data,
47#ifdef EXADG_WITH_PRECICE
48 std::shared_ptr<precice::SolverInterface> precice,
49#endif
50 std::string const mesh_name,
51 dealii::types::boundary_id const surface_id,
52 int const mf_dof_index);
57 virtual void
58 define_coupling_mesh() override;
59
68 virtual void
69 write_data(dealii::LinearAlgebra::distributed::Vector<double> const & data_vector,
70 std::string const & data_name) override;
71
72
73private:
75 std::vector<int> coupling_nodes_ids;
77 std::vector<std::array<dealii::types::global_dof_index, data_dim>> global_indices;
78
81 int const mf_dof_index;
82
83 virtual std::string
84 get_surface_type() const override;
85};
86
87
88
89template<int dim, int data_dim, typename VectorizedArrayType>
91 dealii::MatrixFree<dim, double, VectorizedArrayType> const & data,
92#ifdef EXADG_WITH_PRECICE
93 std::shared_ptr<precice::SolverInterface> precice,
94#endif
95 std::string const mesh_name,
96 dealii::types::boundary_id const surface_id,
97 int const mf_dof_index)
98 : CouplingBase<dim, data_dim, VectorizedArrayType>(data,
99#ifdef EXADG_WITH_PRECICE
100 precice,
101#endif
102 mesh_name,
103 surface_id),
104 mf_dof_index(mf_dof_index)
105{
106}
107
108
109
110template<int dim, int data_dim, typename VectorizedArrayType>
111void
113{
114 Assert(this->mesh_id != -1, dealii::ExcNotInitialized());
115
116 // In order to avoid that we define the surface multiple times when reader
117 // and writer refer to the same object
118 if(coupling_nodes_ids.size() > 0)
119 return;
120
121 // Get and sort the global dof indices
122 auto const get_component_dofs = [&](int const component) {
123 // Get a component mask of the vector component
124 dealii::ComponentMask component_mask(data_dim, false);
125 component_mask.set(component, true);
126
127 // Get the global DoF indices of the component
128 // Compute the intersection with locally owned dofs
129 // TODO: This is super inefficient, have a look at the
130 // dof_handler.n_boundary_dofs implementation for a proper version
131 dealii::IndexSet const indices =
132 (dealii::DoFTools::extract_boundary_dofs(this->matrix_free.get_dof_handler(mf_dof_index),
133 component_mask,
134 std::set<dealii::types::boundary_id>{
135 this->dealii_boundary_surface_id}) &
136 this->matrix_free.get_dof_handler(mf_dof_index).locally_owned_dofs());
137
138 Assert(indices.n_elements() * data_dim ==
139 this->matrix_free.get_dof_handler(mf_dof_index)
140 .n_boundary_dofs(
141 std::set<dealii::types::boundary_id>{this->dealii_boundary_surface_id}),
142 dealii::ExcInternalError());
143 // Resize the global dof index conatiner in case we call this lambda for
144 // the first time
145 if(component == 0)
146 global_indices.resize(indices.n_elements());
147 // fill the first array entry with the respective component
148 dealii::types::global_dof_index iterator = 0;
149 for(auto const dof : indices)
150 {
151 global_indices[iterator][component] = dof;
152 ++iterator;
153 }
154 };
155
156 // Fill the indices object this class works on
157 for(int d = 0; d < data_dim; ++d)
158 get_component_dofs(d);
159 // Compute the coordinates of the indices (only one component required here)
160 // We select the zeroth component
161 std::map<dealii::types::global_dof_index, dealii::Point<dim>> support_points;
162 dealii::ComponentMask component_mask(data_dim, false);
163 component_mask.set(0, true);
164
165 dealii::DoFTools::map_boundary_dofs_to_support_points(
166 *(this->matrix_free.get_mapping_info().mapping),
167 this->matrix_free.get_dof_handler(mf_dof_index),
168 support_points,
169 component_mask,
170 this->dealii_boundary_surface_id);
171
172
173 // Set size of the preCICE ID vector
174 coupling_nodes_ids.reserve(global_indices.size());
175 std::array<double, dim> nodes_position;
176 for(std::size_t i = 0; i < global_indices.size(); ++i)
177 {
178 // Get index of the zeroth component
179 auto const element = global_indices[i][0];
180 for(int d = 0; d < dim; ++d)
181 nodes_position[d] = support_points[element][d];
182
183 // pass node coordinates to precice
184#ifdef EXADG_WITH_PRECICE
185 int const precice_id = this->precice->setMeshVertex(this->mesh_id, nodes_position.data());
186#else
187 int const precice_id = 0;
188#endif
189 coupling_nodes_ids.emplace_back(precice_id);
190 }
191
192#ifdef EXADG_WITH_PRECICE
193 if(this->read_data_map.size() > 0)
194 this->print_info(true, this->precice->getMeshVertexSize(this->mesh_id));
195 if(this->write_data_map.size() > 0)
196 this->print_info(false, this->precice->getMeshVertexSize(this->mesh_id));
197#endif
198}
199
200
201
202template<int dim, int data_dim, typename VectorizedArrayType>
203void
205 dealii::LinearAlgebra::distributed::Vector<double> const & data_vector,
206 std::string const & data_name)
207{
208 int const write_data_id = this->write_data_map.at(data_name);
209 Assert(write_data_id != -1, dealii::ExcNotInitialized());
210 Assert(coupling_nodes_ids.size() > 0, dealii::ExcNotInitialized());
211
212 std::array<double, data_dim> write_data;
213 for(std::size_t i = 0; i < global_indices.size(); ++i)
214 {
215 // Extract relevant elements from global vector
216 for(unsigned int d = 0; d < data_dim; ++d)
217 {
218 auto const element = global_indices[i][d];
219 write_data[d] = data_vector[element];
220 }
221
222#ifdef EXADG_WITH_PRECICE
223 // and pass them to preCICE
224 if constexpr(data_dim > 1)
225 {
226 this->precice->writeVectorData(write_data_id, coupling_nodes_ids[i], write_data.data());
227 }
228 else
229 {
230 this->precice->writeScalarData(write_data_id, coupling_nodes_ids[i], write_data[0]);
231 }
232#else
233 (void)write_data_id;
234#endif
235 }
236}
237
238
239
240template<int dim, int data_dim, typename VectorizedArrayType>
241std::string
243{
244 return "DoF support points using matrix-free dof index " +
245 dealii::Utilities::to_string(mf_dof_index);
246}
247
248} // namespace preCICE
249} // namespace ExaDG
250
251#endif
Definition coupling_base.h:66
std::string const mesh_name
public precice solverinterface
Definition coupling_base.h:152
Definition dof_coupling.h:44
virtual void define_coupling_mesh() override
define_mesh_vertices Define a vertex coupling mesh for preCICE coupling the classical preCICE way
Definition dof_coupling.h:112
virtual void write_data(dealii::LinearAlgebra::distributed::Vector< double > const &data_vector, std::string const &data_name) override
write_data Evaluates the given
Definition dof_coupling.h:204
Definition driver.cpp:33