ExaDG
Loading...
Searching...
No Matches
coupling_base.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 EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_COUPLING_BASE_H_
23#define EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_COUPLING_BASE_H_
24
25// deal.II
26#include <deal.II/base/conditional_ostream.h>
27#include <deal.II/base/vectorization.h>
28#include <deal.II/matrix_free/matrix_free.h>
29
30// ExaDG
31#include <exadg/matrix_free/integrators.h>
32
33// preCICE
34#ifdef EXADG_WITH_PRECICE
35# include <precice/SolverInterface.hpp>
36#endif
37
38namespace ExaDG
39{
40namespace preCICE
41{
45enum class WriteDataType
46{
47 undefined,
48 values_on_dofs,
49 values_on_other_mesh,
50 gradients_on_other_mesh,
51 values_on_q_points,
52 normal_gradients_on_q_points
53};
54
63template<int dim, int data_dim, typename VectorizedArrayType>
64class CouplingBase
65{
66public:
67 CouplingBase(dealii::MatrixFree<dim, double, VectorizedArrayType> const & data,
68#ifdef EXADG_WITH_PRECICE
69 std::shared_ptr<precice::SolverInterface> precice,
70#endif
71 std::string const mesh_name,
72 dealii::types::boundary_id const surface_id);
73
74 virtual ~CouplingBase() = default;
75
77 using FEFaceIntegrator = FaceIntegrator<dim, data_dim, double, VectorizedArrayType>;
78 using value_type = typename FEFaceIntegrator::value_type;
83 virtual void
85
92 virtual void
94
102 virtual void
103 write_data(dealii::LinearAlgebra::distributed::Vector<double> const & data_vector,
104 std::string const & data_name) = 0;
105
106
107 virtual void
108 read_block_data(std::string const & data_name) const;
109
114 void
115 add_read_data(std::string const & read_data_name);
116
121 void
122 add_write_data(std::string const & write_data_name);
123
128 void
129 set_write_data_type(WriteDataType write_data_specification);
130
131protected:
139 void
140 print_info(bool const reader, unsigned int const local_size) const;
141
143 dealii::MatrixFree<dim, double, VectorizedArrayType> const & matrix_free;
144
146#ifdef EXADG_WITH_PRECICE
147 std::shared_ptr<precice::SolverInterface> precice;
148#endif
149
151 std::string const mesh_name;
152 int mesh_id;
153 // Map between data ID (preCICE) and the data name
154 std::map<std::string, int> read_data_map;
155 std::map<std::string, int> write_data_map;
156
157 dealii::types::boundary_id const dealii_boundary_surface_id;
158
159 WriteDataType write_data_type;
160
161 virtual std::string
162 get_surface_type() const = 0;
163};
164
165
166
167template<int dim, int data_dim, typename VectorizedArrayType>
168CouplingBase<dim, data_dim, VectorizedArrayType>::CouplingBase(
169 dealii::MatrixFree<dim, double, VectorizedArrayType> const & matrix_free_,
170#ifdef EXADG_WITH_PRECICE
171 std::shared_ptr<precice::SolverInterface> precice,
172#endif
173 std::string const mesh_name,
174 dealii::types::boundary_id const surface_id)
175 : matrix_free(matrix_free_),
176#ifdef EXADG_WITH_PRECICE
177 precice(precice),
178#endif
179 mesh_name(mesh_name),
180 dealii_boundary_surface_id(surface_id),
181 write_data_type(WriteDataType::undefined)
182{
183#ifdef EXADG_WITH_PRECICE
184 Assert(precice.get() != nullptr, dealii::ExcNotInitialized());
185
186 // Ask preCICE already in the constructor for the IDs
187 mesh_id = precice->getMeshID(mesh_name);
188#else
189 AssertThrow(false,
190 dealii::ExcMessage("EXADG_WITH_PRECICE has to be activated to use this code."));
191 mesh_id = 0;
192#endif
193}
194
195
196template<int dim, int data_dim, typename VectorizedArrayType>
197void
199{
200 Assert(mesh_id != -1, dealii::ExcNotInitialized());
201#ifdef EXADG_WITH_PRECICE
202 int const read_data_id = precice->getDataID(read_data_name, mesh_id);
203#else
204 int const read_data_id = 0;
205#endif
206 read_data_map.insert({read_data_name, read_data_id});
207}
208
209
210
211template<int dim, int data_dim, typename VectorizedArrayType>
212void
214 std::string const & write_data_name)
215{
216 Assert(mesh_id != -1, dealii::ExcNotInitialized());
217#ifdef EXADG_WITH_PRECICE
218 int const write_data_id = precice->getDataID(write_data_name, mesh_id);
219#else
220 int const write_data_id = 0;
221#endif
222 write_data_map.insert({write_data_name, write_data_id});
223}
224
225
226template<int dim, int data_dim, typename VectorizedArrayType>
227void
229 WriteDataType write_data_specification)
230{
231 write_data_type = write_data_specification;
232}
233
234
235
236template<int dim, int data_dim, typename VectorizedArrayType>
237void
242
243
244
245template<int dim, int data_dim, typename VectorizedArrayType>
246void
247CouplingBase<dim, data_dim, VectorizedArrayType>::read_block_data(std::string const &) const
248{
249 AssertThrow(false, dealii::ExcNotImplemented());
250}
251
252
253template<int dim, int data_dim, typename VectorizedArrayType>
254void
256 unsigned int const local_size) const
257{
258 dealii::ConditionalOStream pcout(std::cout,
259 dealii::Utilities::MPI::this_mpi_process(
260 matrix_free.get_dof_handler().get_mpi_communicator()) == 0);
261 auto const map = (reader ? read_data_map : write_data_map);
262
263 auto names = map.begin();
264 std::string data_names = names->first;
265 ++names;
266 for(; names != map.end(); ++names)
267 data_names += std::string(", ") + names->first;
268
269 pcout << "-- Data " << (reader ? "reading" : "writing") << ":\n"
270 << "-- . data name(s): " << data_names << "\n"
271 << "-- . associated mesh: " << mesh_name << "\n"
272 << "-- . Number of coupling nodes: "
273 << dealii::Utilities::MPI::sum(local_size,
274 matrix_free.get_dof_handler().get_mpi_communicator())
275 << "\n"
276 << "-- . Node location: " << get_surface_type() << "\n"
277 << std::endl;
278}
279
280} // namespace preCICE
281} // namespace ExaDG
282
283#endif /* EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_COUPLING_BASE_H_ */
void add_write_data(std::string const &write_data_name)
Queries data IDs from preCICE for the given write data name.
Definition coupling_base.h:213
std::string const mesh_name
public precice solverinterface
Definition coupling_base.h:151
virtual void process_coupling_mesh()
process_coupling_mesh (optional) Handle post-preCICE-initialization steps, e.g. do computations on re...
Definition coupling_base.h:238
FaceIntegrator< dim, data_dim, double, VectorizedArrayType > FEFaceIntegrator
Alias for the face integrator.
Definition coupling_base.h:77
virtual void write_data(dealii::LinearAlgebra::distributed::Vector< double > const &data_vector, std::string const &data_name)=0
write_data Write the data associated to the defined vertices to preCICE
dealii::MatrixFree< dim, double, VectorizedArrayType > const & matrix_free
The dealii::MatrixFree object (preCICE can only handle double precision)
Definition coupling_base.h:143
void set_write_data_type(WriteDataType write_data_specification)
Set the WriteDataType in this class which determines the location of the write data (e....
Definition coupling_base.h:228
virtual void define_coupling_mesh()=0
define_coupling_mesh Define the coupling mesh associated to the data points
void add_read_data(std::string const &read_data_name)
Queries data IDs from preCICE for the given read data name.
Definition coupling_base.h:198
void print_info(bool const reader, unsigned int const local_size) const
Print information of the current setup.
Definition coupling_base.h:255
Definition driver.cpp:33