ExaDG
Loading...
Searching...
No Matches
precice_adapter.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_PRECICE_ADAPTER_H_
23#define INCLUDE_EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_PRECICE_ADAPTER_H_
24
25// deal.II
26#include <deal.II/base/exceptions.h>
27#include <deal.II/base/mpi.h>
28#include <deal.II/base/utilities.h>
29
30#include <deal.II/fe/fe.h>
31#include <deal.II/fe/mapping_q_generic.h>
32
33#include <deal.II/matrix_free/matrix_free.h>
34
35// ExaDG
36#include <exadg/fluid_structure_interaction/precice/dof_coupling.h>
37#include <exadg/fluid_structure_interaction/precice/exadg_coupling.h>
38#include <exadg/fluid_structure_interaction/precice/quad_coupling.h>
39#include <exadg/utilities/tensor_utilities.h>
40
41// preCICE
42#ifdef EXADG_WITH_PRECICE
43# include <precice/SolverInterface.hpp>
44#endif
45
46#include <ostream>
47
48namespace ExaDG
49{
50namespace preCICE
51{
57template<int dim,
58 int data_dim,
59 typename VectorType,
60 typename VectorizedArrayType = dealii::VectorizedArray<double>>
62{
63public:
64 static unsigned int const rank = n_components_to_rank<data_dim, dim>();
65
66 using value_type = typename CouplingBase<dim, data_dim, VectorizedArrayType>::value_type;
88 template<typename ParameterClass>
89 Adapter(ParameterClass const & parameters, MPI_Comm mpi_comm);
90
91
102 void
103 initialize_precice(VectorType const & dealii_to_precice);
104
105 void
106 add_write_surface(dealii::types::boundary_id const surface_id,
107 std::string const & mesh_name,
108 std::vector<std::string> const & write_data_names,
109 WriteDataType const write_data_type,
110 dealii::MatrixFree<dim, double, VectorizedArrayType> const & data,
111 unsigned int const dof_index,
112 unsigned int const write_quad_index);
113
114
115 void
116 add_read_surface(dealii::MatrixFree<dim, double, VectorizedArrayType> const & data,
117 std::shared_ptr<ContainerInterfaceData<rank, dim, double>> interface_data,
118 std::string const & mesh_name,
119 std::vector<std::string> const & read_data_name);
120
124 void
125 advance(double const computed_timestep_length);
126
139 void
140 save_current_state_if_required(std::function<void()> const & save_state);
141
156 void
157 reload_old_state_if_required(std::function<void()> const & reload_old_state);
158
159
160 void
161 write_data(std::string const & write_mesh_name,
162 std::string const & write_data_name,
163 VectorType const & write_data,
164 double const computed_timestep_length);
165
166 void
167 read_block_data(std::string const & mesh_name, const std::string & data_name) const;
168
174 bool
175 is_coupling_ongoing() const;
176
183 bool
185
186
187private:
188 // public precice solverinterface, needed in order to steer the time loop
189 // inside the solver.
190#ifdef EXADG_WITH_PRECICE
191 std::shared_ptr<precice::SolverInterface> precice;
192#endif
193
195 std::map<std::string, std::shared_ptr<CouplingBase<dim, data_dim, VectorizedArrayType>>> writer;
196
197 // We restrict the reader to be of type ExaDGCoupling for the moment, as all other choices don't
198 // make sense
199 std::map<std::string, std::shared_ptr<ExaDGCoupling<dim, data_dim, VectorizedArrayType>>> reader;
200
201 // Container to store time dependent data in case of an implicit coupling
202 std::vector<VectorType> old_state_data;
203};
204
205
206
207template<int dim, int data_dim, typename VectorType, typename VectorizedArrayType>
208template<typename ParameterClass>
210 MPI_Comm mpi_comm)
211
212{
213#ifdef EXADG_WITH_PRECICE
214 precice =
215 std::make_shared<precice::SolverInterface>(parameters.participant_name,
216 parameters.config_file,
217 dealii::Utilities::MPI::this_mpi_process(mpi_comm),
218 dealii::Utilities::MPI::n_mpi_processes(mpi_comm));
219
220 AssertThrow(dim == precice->getDimensions(), dealii::ExcInternalError());
221#else
222 (void)parameters;
223 (void)mpi_comm;
224 AssertThrow(false,
225 dealii::ExcMessage("EXADG_WITH_PRECICE has to be activated to use this code."));
226#endif
227
228 AssertThrow(dim > 1, dealii::ExcNotImplemented());
229}
230
231
232
233template<int dim, int data_dim, typename VectorType, typename VectorizedArrayType>
234void
236 dealii::types::boundary_id const dealii_boundary_surface_id,
237 std::string const & mesh_name,
238 std::vector<std::string> const & write_data_names,
239 WriteDataType const write_data_type,
240 dealii::MatrixFree<dim, double, VectorizedArrayType> const & data,
241 unsigned int const dof_index,
242 unsigned int const quad_index)
243{
244 // Check, if we already have such an interface
245 auto const found_reader = reader.find(mesh_name);
246
247 if(found_reader != reader.end())
248 {
249 // Duplicate the shared pointer
250 writer.insert({mesh_name, found_reader->second});
251 }
252 else if(write_data_type == WriteDataType::values_on_dofs)
253 {
254 writer.insert(
255 {mesh_name,
256 std::make_shared<DoFCoupling<dim, data_dim, VectorizedArrayType>>(data,
257#ifdef EXADG_WITH_PRECICE
258 precice,
259#endif
260 mesh_name,
261 dealii_boundary_surface_id,
262 dof_index)});
263 }
264 else
265 {
266 Assert(write_data_type == WriteDataType::values_on_q_points or
267 write_data_type == WriteDataType::normal_gradients_on_q_points,
268 dealii::ExcNotImplemented());
269 writer.insert({mesh_name,
270 std::make_shared<QuadCoupling<dim, data_dim, VectorizedArrayType>>(
271 data,
272#ifdef EXADG_WITH_PRECICE
273 precice,
274#endif
275 mesh_name,
276 dealii_boundary_surface_id,
277 dof_index,
278 quad_index)});
279 }
280
281 // Register the write data
282 std::vector<dealii::Point<dim>> const points;
283 for(auto const & data_name : write_data_names)
284 writer.at(mesh_name)->add_write_data(data_name);
285
286 writer.at(mesh_name)->set_write_data_type(write_data_type);
287 // Finally initialize the surface
288 writer.at(mesh_name)->define_coupling_mesh();
289}
290
291
292
293template<int dim, int data_dim, typename VectorType, typename VectorizedArrayType>
294void
295Adapter<dim, data_dim, VectorType, VectorizedArrayType>::add_read_surface(
296 dealii::MatrixFree<dim, double, VectorizedArrayType> const & data,
297 std::shared_ptr<ContainerInterfaceData<rank, dim, double>> interface_data,
298 std::string const & mesh_name,
299 std::vector<std::string> const & read_data_names)
300{
301 reader.insert(
302 {mesh_name,
303 std::make_shared<ExaDGCoupling<dim, data_dim, VectorizedArrayType>>(data,
304#ifdef EXADG_WITH_PRECICE
305 precice,
306#endif
307 mesh_name,
308 interface_data)});
309
310 for(auto const & data_name : read_data_names)
311 reader.at(mesh_name)->add_read_data(data_name);
312 reader.at(mesh_name)->define_coupling_mesh();
313}
314
315
316
317template<int dim, int data_dim, typename VectorType, typename VectorizedArrayType>
318void
320 VectorType const & dealii_to_precice)
321{
322#ifdef EXADG_WITH_PRECICE
323 // if(not dealii_to_precice.has_ghost_elements())
324 // dealii_to_precice.update_ghost_values();
325
326 // Initialize preCICE internally
327 precice->initialize();
328
329 // Only the writer needs potentially to process the coupling mesh, if the
330 // mapping is carried out in the solver
331 // writer->process_coupling_mesh();
332
333 // write initial writeData to preCICE if required
334 if(precice->isActionRequired(precice::constants::actionWriteInitialData()))
335 {
336 writer[0]->write_data(dealii_to_precice, "");
337
338 precice->markActionFulfilled(precice::constants::actionWriteInitialData());
339 }
340 precice->initializeData();
341
342 // Maybe, read block-wise and work with an AlignedVector since the read data
343 // (forces) is multiple times required during the Newton iteration
344 // if (shared_memory_parallel and precice->isReadDataAvailable())
345 // precice->readBlockVectorData(read_data_id,
346 // read_nodes_ids.size(),
347 // read_nodes_ids.data(),
348 // read_data.data());
349#else
350 (void)dealii_to_precice;
351#endif
352}
353
354template<int dim, int data_dim, typename VectorType, typename VectorizedArrayType>
355void
357 std::string const & write_mesh_name,
358 std::string const & write_data_name,
359 VectorType const & dealii_to_precice,
360 double const computed_timestep_length)
361{
362#ifdef EXADG_WITH_PRECICE
363 if(precice->isWriteDataRequired(computed_timestep_length))
364 writer.at(write_mesh_name)->write_data(dealii_to_precice, write_data_name);
365#else
366 (void)write_mesh_name;
367 (void)write_data_name;
368 (void)dealii_to_precice;
369 (void)computed_timestep_length;
370#endif
371}
372
373template<int dim, int data_dim, typename VectorType, typename VectorizedArrayType>
374void
376 double const computed_timestep_length)
377{
378#ifdef EXADG_WITH_PRECICE
379 // Here, we need to specify the computed time step length and pass it to
380 // preCICE
381 // TODO: The function returns the available time-step size which is required for time-step sync
382 precice->advance(computed_timestep_length);
383#else
384 (void)computed_timestep_length;
385#endif
386}
387
388
389
390template<int dim, int data_dim, typename VectorType, typename VectorizedArrayType>
391void
393 std::string const & mesh_name,
394 std::string const & data_name) const
395{
396 reader.at(mesh_name)->read_block_data(data_name);
397}
398
399
400
401template<int dim, int data_dim, typename VectorType, typename VectorizedArrayType>
402inline void
404 std::function<void()> const & save_state)
405{
406#ifdef EXADG_WITH_PRECICE
407 // First, we let preCICE check, whether we need to store the variables.
408 // Then, the data is stored in the class
409 if(precice->isActionRequired(precice::constants::actionWriteIterationCheckpoint()))
410 {
411 save_state();
412 precice->markActionFulfilled(precice::constants::actionWriteIterationCheckpoint());
413 }
414#else
415 (void)save_state;
416#endif
417}
418
419
420
421template<int dim, int data_dim, typename VectorType, typename VectorizedArrayType>
422inline void
424 std::function<void()> const & reload_old_state)
425{
426#ifdef EXADG_WITH_PRECICE
427 // In case we need to reload a state, we just take the internally stored
428 // data vectors and write then in to the input data
429 if(precice->isActionRequired(precice::constants::actionReadIterationCheckpoint()))
430 {
431 reload_old_state();
432 precice->markActionFulfilled(precice::constants::actionReadIterationCheckpoint());
433 }
434#else
435 (void)reload_old_state;
436#endif
437}
438
439
440
441template<int dim, int data_dim, typename VectorType, typename VectorizedArrayType>
442inline bool
444{
445#ifdef EXADG_WITH_PRECICE
446 return precice->isCouplingOngoing();
447#else
448 return true;
449#endif
450}
451
452
453
454template<int dim, int data_dim, typename VectorType, typename VectorizedArrayType>
455inline bool
457{
458#ifdef EXADG_WITH_PRECICE
459 return precice->isTimeWindowComplete();
460#else
461 return true;
462#endif
463}
464
465
466} // namespace preCICE
467} // namespace ExaDG
468
469#endif
Definition container_interface_data.h:45
Definition precice_adapter.h:62
void reload_old_state_if_required(std::function< void()> const &reload_old_state)
Reloads the previously stored variables in case of an implicit coupling. The current implementation s...
Definition precice_adapter.h:423
void initialize_precice(VectorType const &dealii_to_precice)
Initializes preCICE and passes all relevant data to preCICE.
Definition precice_adapter.h:319
Adapter(ParameterClass const &parameters, MPI_Comm mpi_comm)
Constructor, which sets up the precice Solverinterface.
Definition precice_adapter.h:209
void save_current_state_if_required(std::function< void()> const &save_state)
Saves current state of time dependent variables in case of an implicit coupling.
Definition precice_adapter.h:403
bool is_coupling_ongoing() const
is_coupling_ongoing Calls the preCICE API function isCouplingOnGoing
Definition precice_adapter.h:443
void advance(double const computed_timestep_length)
Advances preCICE after every timestep.
Definition precice_adapter.h:375
bool is_time_window_complete() const
is_time_window_complete Calls the preCICE API function isTimeWindowComplete
Definition precice_adapter.h:456
Definition driver.cpp:33