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 EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_DOF_COUPLING_H_
24#define EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_DOF_COUPLING_H_
25
26// deal.II
27#include <deal.II/dofs/dof_tools.h>
28
29// ExaDG
30#include <exadg/fluid_structure_interaction/precice/coupling_base.h>
31#include <exadg/fluid_structure_interaction/precice/dof_tools_extension.h>
32
33namespace ExaDG
34{
35namespace preCICE
36{
44template<int dim, int data_dim, typename VectorizedArrayType>
45class DoFCoupling : public CouplingBase<dim, data_dim, VectorizedArrayType>
46{
47public:
48 DoFCoupling(dealii::MatrixFree<dim, double, VectorizedArrayType> const & data,
49#ifdef EXADG_WITH_PRECICE
50 std::shared_ptr<precice::SolverInterface> precice,
51#endif
52 std::string const mesh_name,
53 dealii::types::boundary_id const surface_id,
54 int const mf_dof_index);
59 virtual void
60 define_coupling_mesh() override;
61
70 virtual void
71 write_data(dealii::LinearAlgebra::distributed::Vector<double> const & data_vector,
72 std::string const & data_name) override;
73
74
75private:
77 std::vector<int> coupling_nodes_ids;
79 std::vector<std::array<dealii::types::global_dof_index, data_dim>> global_indices;
80
83 int const mf_dof_index;
84
85 virtual std::string
86 get_surface_type() const override;
87};
88
89
90
91template<int dim, int data_dim, typename VectorizedArrayType>
92DoFCoupling<dim, data_dim, VectorizedArrayType>::DoFCoupling(
93 dealii::MatrixFree<dim, double, VectorizedArrayType> const & data,
94#ifdef EXADG_WITH_PRECICE
95 std::shared_ptr<precice::SolverInterface> precice,
96#endif
97 std::string const mesh_name,
98 dealii::types::boundary_id const surface_id,
99 int const mf_dof_index)
100 : CouplingBase<dim, data_dim, VectorizedArrayType>(data,
101#ifdef EXADG_WITH_PRECICE
102 precice,
103#endif
104 mesh_name,
105 surface_id),
106 mf_dof_index(mf_dof_index)
107{
108}
109
110
111
112template<int dim, int data_dim, typename VectorizedArrayType>
113void
115{
116 Assert(this->mesh_id != -1, dealii::ExcNotInitialized());
117
118 // In order to avoid that we define the surface multiple times when reader
119 // and writer refer to the same object
120 if(coupling_nodes_ids.size() > 0)
121 return;
122
123 // Get and sort the global dof indices
124 auto const get_component_dofs = [&](int const component) {
125 // Get a component mask of the vector component
126 dealii::ComponentMask component_mask(data_dim, false);
127 component_mask.set(component, true);
128
129 // Get the global DoF indices of the component
130 // Compute the intersection with locally owned dofs
131 // TODO: This is super inefficient, have a look at the
132 // dof_handler.n_boundary_dofs implementation for a proper version
133 dealii::IndexSet const indices =
134 (dealii::DoFTools::extract_boundary_dofs(this->matrix_free.get_dof_handler(mf_dof_index),
135 component_mask,
136 std::set<dealii::types::boundary_id>{
137 this->dealii_boundary_surface_id}) &
138 this->matrix_free.get_dof_handler(mf_dof_index).locally_owned_dofs());
139
140 Assert(indices.n_elements() * data_dim ==
141 this->matrix_free.get_dof_handler(mf_dof_index)
142 .n_boundary_dofs(
143 std::set<dealii::types::boundary_id>{this->dealii_boundary_surface_id}),
144 dealii::ExcInternalError());
145 // Resize the global dof index conatiner in case we call this lambda for
146 // the first time
147 if(component == 0)
148 global_indices.resize(indices.n_elements());
149 // fill the first array entry with the respective component
150 dealii::types::global_dof_index iterator = 0;
151 for(auto const dof : indices)
152 {
153 global_indices[iterator][component] = dof;
154 ++iterator;
155 }
156 };
157
158 // Fill the indices object this class works on
159 for(int d = 0; d < data_dim; ++d)
160 get_component_dofs(d);
161 // Compute the coordinates of the indices (only one component required here)
162 // We select the zeroth component
163 std::map<dealii::types::global_dof_index, dealii::Point<dim>> support_points;
164 dealii::ComponentMask component_mask(data_dim, false);
165 component_mask.set(0, true);
166
167 dealii::DoFTools::map_boundary_dofs_to_support_points(
168 *(this->matrix_free.get_mapping_info().mapping),
169 this->matrix_free.get_dof_handler(mf_dof_index),
170 support_points,
171 component_mask,
172 this->dealii_boundary_surface_id);
173
174
175 // Set size of the preCICE ID vector
176 coupling_nodes_ids.reserve(global_indices.size());
177 std::array<double, dim> nodes_position;
178 for(std::size_t i = 0; i < global_indices.size(); ++i)
179 {
180 // Get index of the zeroth component
181 auto const element = global_indices[i][0];
182 for(int d = 0; d < dim; ++d)
183 nodes_position[d] = support_points[element][d];
184
185 // pass node coordinates to precice
186#ifdef EXADG_WITH_PRECICE
187 int const precice_id = this->precice->setMeshVertex(this->mesh_id, nodes_position.data());
188#else
189 int const precice_id = 0;
190#endif
191 coupling_nodes_ids.emplace_back(precice_id);
192 }
193
194#ifdef EXADG_WITH_PRECICE
195 if(this->read_data_map.size() > 0)
196 this->print_info(true, this->precice->getMeshVertexSize(this->mesh_id));
197 if(this->write_data_map.size() > 0)
198 this->print_info(false, this->precice->getMeshVertexSize(this->mesh_id));
199#endif
200}
201
202
203
204template<int dim, int data_dim, typename VectorizedArrayType>
205void
207 dealii::LinearAlgebra::distributed::Vector<double> const & data_vector,
208 std::string const & data_name)
209{
210 int const write_data_id = this->write_data_map.at(data_name);
211 Assert(write_data_id != -1, dealii::ExcNotInitialized());
212 Assert(coupling_nodes_ids.size() > 0, dealii::ExcNotInitialized());
213
214 std::array<double, data_dim> write_data;
215 for(std::size_t i = 0; i < global_indices.size(); ++i)
216 {
217 // Extract relevant elements from global vector
218 for(unsigned int d = 0; d < data_dim; ++d)
219 {
220 auto const element = global_indices[i][d];
221 write_data[d] = data_vector[element];
222 }
223
224#ifdef EXADG_WITH_PRECICE
225 // and pass them to preCICE
226 if constexpr(data_dim > 1)
227 {
228 this->precice->writeVectorData(write_data_id, coupling_nodes_ids[i], write_data.data());
229 }
230 else
231 {
232 this->precice->writeScalarData(write_data_id, coupling_nodes_ids[i], write_data[0]);
233 }
234#else
235 (void)write_data_id;
236#endif
237 }
238}
239
240
241
242template<int dim, int data_dim, typename VectorizedArrayType>
243std::string
244DoFCoupling<dim, data_dim, VectorizedArrayType>::get_surface_type() const
245{
246 return "DoF support points using matrix-free dof index " +
247 dealii::Utilities::to_string(mf_dof_index);
248}
249
250} // namespace preCICE
251} // namespace ExaDG
252
253#endif /* EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_DOF_COUPLING_H_ */
Definition coupling_base.h:65
std::string const mesh_name
public precice solverinterface
Definition coupling_base.h:151
dealii::MatrixFree< dim, double, VectorizedArrayType > const & matrix_free
The dealii::MatrixFree object (preCICE can only handle double precision)
Definition coupling_base.h:143
void print_info(bool const reader, unsigned int const local_size) const
Print information of the current setup.
Definition coupling_base.h:255
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:114
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:206
Definition driver.cpp:33