ExaDG
Loading...
Searching...
No Matches
quad_coupling.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2022 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef INCLUDE_EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_QUAD_COUPLING_H_
23#define INCLUDE_EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_QUAD_COUPLING_H_
24
25#include <deal.II/matrix_free/fe_evaluation.h>
26
27#include <exadg/fluid_structure_interaction/precice/coupling_base.h>
28
29namespace ExaDG
30{
31namespace preCICE
32{
40template<int dim, int data_dim, typename VectorizedArrayType>
41class QuadCoupling : public CouplingBase<dim, data_dim, VectorizedArrayType>
42{
43public:
44 QuadCoupling(dealii::MatrixFree<dim, double, VectorizedArrayType> const & data,
45#ifdef EXADG_WITH_PRECICE
46 std::shared_ptr<precice::SolverInterface> precice,
47#endif
48 std::string const mesh_name,
49 dealii::types::boundary_id const surface_id,
50 int const mf_dof_index,
51 int const mf_quad_index);
52
56 using value_type = typename CouplingBase<dim, data_dim, VectorizedArrayType>::value_type;
57
62 virtual void
63 define_coupling_mesh() override;
64
77 virtual void
78 write_data(dealii::LinearAlgebra::distributed::Vector<double> const & data_vector,
79 std::string const & data_name) override;
80
81private:
91 void
92 write_data_factory(
93 dealii::LinearAlgebra::distributed::Vector<double> const & data_vector,
94 int const write_data_id,
95 dealii::EvaluationFlags::EvaluationFlags const flags,
96 std::function<value_type(FEFaceIntegrator &, unsigned int)> const & get_write_value);
97
99 std::vector<std::array<int, VectorizedArrayType::size()>> coupling_nodes_ids;
100
103 int const mf_dof_index;
104 int const mf_quad_index;
105
106 virtual std::string
107 get_surface_type() const override;
108};
109
110
111
112template<int dim, int data_dim, typename VectorizedArrayType>
114 dealii::MatrixFree<dim, double, VectorizedArrayType> const & data,
115#ifdef EXADG_WITH_PRECICE
116 std::shared_ptr<precice::SolverInterface> precice,
117#endif
118 std::string const mesh_name,
119 dealii::types::boundary_id const surface_id,
120 int const mf_dof_index_,
121 int const mf_quad_index_)
122 : CouplingBase<dim, data_dim, VectorizedArrayType>(data,
123#ifdef EXADG_WITH_PRECICE
124 precice,
125#endif
126 mesh_name,
127 surface_id),
128 mf_dof_index(mf_dof_index_),
129 mf_quad_index(mf_quad_index_)
130{
131}
132
133
134
135template<int dim, int data_dim, typename VectorizedArrayType>
136void
138{
139 Assert(this->mesh_id != -1, dealii::ExcNotInitialized());
140
141 // In order to avoid that we define the surface multiple times when reader
142 // and writer refer to the same object
143 if(coupling_nodes_ids.size() > 0)
144 return;
145
146 // Initial guess: half of the boundary is part of the coupling surface
147 coupling_nodes_ids.reserve(this->matrix_free.n_boundary_face_batches() / 2);
148
149 // Set up data structures
150 FEFaceIntegrator phi(this->matrix_free, true, mf_dof_index, mf_quad_index);
151 std::array<double, dim * VectorizedArrayType::size()> unrolled_vertices;
152 std::array<int, VectorizedArrayType::size()> node_ids;
153 unsigned int size = 0;
154 // Loop over all boundary faces
155 for(unsigned int face = this->matrix_free.n_inner_face_batches();
156 face < this->matrix_free.n_boundary_face_batches() + this->matrix_free.n_inner_face_batches();
157 ++face)
158 {
159 auto const boundary_id = this->matrix_free.get_boundary_id(face);
160
161 // Only for interface nodes
162 if(boundary_id != this->dealii_boundary_surface_id)
163 continue;
164
165 phi.reinit(face);
166 int const active_faces = this->matrix_free.n_active_entries_per_face_batch(face);
167
168 // Loop over all quadrature points and pass the vertices to preCICE
169 for(unsigned int q = 0; q < phi.n_q_points; ++q)
170 {
171 auto const local_vertex = phi.quadrature_point(q);
172
173 // Transform dealii::Point<Vectorized> into preCICE conform format
174 // We store here also the potential 'dummy'/empty lanes (not only
175 // active_faces), but it allows us to use a static loop as well as a
176 // static array for the indices
177 for(int d = 0; d < dim; ++d)
178 for(unsigned int v = 0; v < VectorizedArrayType::size(); ++v)
179 unrolled_vertices[d + dim * v] = local_vertex[d][v];
180
181#ifdef EXADG_WITH_PRECICE
182 this->precice->setMeshVertices(this->mesh_id,
183 active_faces,
184 unrolled_vertices.data(),
185 node_ids.data());
186#else
187 (void)active_faces;
188#endif
189 coupling_nodes_ids.emplace_back(node_ids);
190 ++size;
191 }
192 }
193
194 // resize the IDs in case the initial guess was too large
195 coupling_nodes_ids.resize(size);
196
197#ifdef EXADG_WITH_PRECICE
198 // Consistency check: the number of IDs we stored is equal or greater than
199 // the IDs preCICE knows
200 Assert(size * VectorizedArrayType::size() >=
201 static_cast<unsigned int>(this->precice->getMeshVertexSize(this->mesh_id)),
202 dealii::ExcInternalError());
203
204 if(this->read_data_map.size() > 0)
205 this->print_info(true, this->precice->getMeshVertexSize(this->mesh_id));
206 if(this->write_data_map.size() > 0)
207 this->print_info(false, this->precice->getMeshVertexSize(this->mesh_id));
208#endif
209}
210
211
212template<int dim, int data_dim, typename VectorizedArrayType>
213void
215 dealii::LinearAlgebra::distributed::Vector<double> const & data_vector,
216 std::string const & data_name)
217{
218 int const write_data_id = this->write_data_map.at(data_name);
219
220 switch(this->write_data_type)
221 {
222 case WriteDataType::values_on_q_points:
223 write_data_factory(data_vector,
224 write_data_id,
225 dealii::EvaluationFlags::values,
226 [](auto & phi, auto q_point) { return phi.get_value(q_point); });
227 break;
228 case WriteDataType::normal_gradients_on_q_points:
229 write_data_factory(data_vector,
230 write_data_id,
231 dealii::EvaluationFlags::gradients,
232 [](auto & phi, auto q_point) {
233 return phi.get_normal_derivative(q_point);
234 });
235 break;
236 default:
237 AssertThrow(false, dealii::ExcNotImplemented());
238 }
239}
240
241
242
243template<int dim, int data_dim, typename VectorizedArrayType>
244void
246 dealii::LinearAlgebra::distributed::Vector<double> const & data_vector,
247 int const write_data_id,
248 dealii::EvaluationFlags::EvaluationFlags const flags,
249 std::function<value_type(FEFaceIntegrator &, unsigned int)> const & get_write_value)
250{
251 Assert(write_data_id != -1, dealii::ExcNotInitialized());
252 Assert(coupling_nodes_ids.size() > 0, dealii::ExcNotInitialized());
253 // Similar as in define_coupling_mesh
254 FEFaceIntegrator phi(this->matrix_free, true, mf_dof_index, mf_quad_index);
255
256 // In order to unroll the vectorization
257 std::array<double, data_dim * VectorizedArrayType::size()> unrolled_local_data;
258 (void)unrolled_local_data;
259
260 auto index = coupling_nodes_ids.begin();
261
262 // Loop over all faces
263 for(unsigned int face = this->matrix_free.n_inner_face_batches();
264 face < this->matrix_free.n_boundary_face_batches() + this->matrix_free.n_inner_face_batches();
265 ++face)
266 {
267 auto const boundary_id = this->matrix_free.get_boundary_id(face);
268
269 // Only for interface nodes
270 if(boundary_id != this->dealii_boundary_surface_id)
271 continue;
272
273 // Read and interpolate
274 phi.reinit(face);
275 phi.read_dof_values_plain(data_vector);
276 phi.evaluate(flags);
277 int const active_faces = this->matrix_free.n_active_entries_per_face_batch(face);
278
279 for(unsigned int q = 0; q < phi.n_q_points; ++q, ++index)
280 {
281 Assert(index != coupling_nodes_ids.end(), dealii::ExcInternalError());
282 auto const local_data = get_write_value(phi, q);
283
284 // Constexpr evaluation required in order to comply with the
285 // compiler here
286 if constexpr(data_dim > 1)
287 {
288 // Transform Tensor<1,dim,VectorizedArrayType> into preCICE
289 // conform format
290 for(int d = 0; d < data_dim; ++d)
291 for(unsigned int v = 0; v < VectorizedArrayType::size(); ++v)
292 unrolled_local_data[d + data_dim * v] = local_data[d][v];
293
294#ifdef EXADG_WITH_PRECICE
295 this->precice->writeBlockVectorData(write_data_id,
296 active_faces,
297 index->data(),
298 unrolled_local_data.data());
299#else
300 (void)active_faces;
301 (void)write_data_id;
302#endif
303 }
304 else
305 {
306#ifdef EXADG_WITH_PRECICE
307 this->precice->writeBlockScalarData(write_data_id,
308 active_faces,
309 index->data(),
310 &local_data[0]);
311#endif
312 }
313 }
314 }
315}
316
317
318
319template<int dim, int data_dim, typename VectorizedArrayType>
320std::string
321QuadCoupling<dim, data_dim, VectorizedArrayType>::get_surface_type() const
322{
323 return "quadrature points using matrix-free quad index " +
324 dealii::Utilities::to_string(mf_quad_index);
325}
326
327} // namespace preCICE
328} // namespace ExaDG
329
330#endif
Definition coupling_base.h:66
std::string const mesh_name
public precice solverinterface
Definition coupling_base.h:152
FaceIntegrator< dim, data_dim, double, VectorizedArrayType > FEFaceIntegrator
Alias for the face integrator.
Definition coupling_base.h:78
Definition quad_coupling.h:42
virtual void define_coupling_mesh() override
define_mesh_vertices Define a vertex coupling mesh for preCICE coupling the classical preCICE way
Definition quad_coupling.h:137
typename CouplingBase< dim, data_dim, VectorizedArrayType >::FEFaceIntegrator FEFaceIntegrator
Alias as defined in the base class.
Definition quad_coupling.h:54
virtual void write_data(dealii::LinearAlgebra::distributed::Vector< double > const &data_vector, std::string const &data_name) override
write_data Evaluates the given
Definition quad_coupling.h:214
Definition driver.cpp:33