22#ifndef INCLUDE_EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_PRECICE_ADAPTER_H_
23#define INCLUDE_EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_PRECICE_ADAPTER_H_
26#include <deal.II/base/exceptions.h>
27#include <deal.II/base/mpi.h>
28#include <deal.II/base/utilities.h>
30#include <deal.II/fe/fe.h>
31#include <deal.II/fe/mapping_q_generic.h>
33#include <deal.II/matrix_free/matrix_free.h>
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>
42#ifdef EXADG_WITH_PRECICE
43# include <precice/SolverInterface.hpp>
60 typename VectorizedArrayType = dealii::VectorizedArray<double>>
64 static unsigned int const rank = n_components_to_rank<data_dim, dim>();
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);
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);
116 add_read_surface(dealii::MatrixFree<dim, double, VectorizedArrayType>
const & data,
118 std::string
const & mesh_name,
119 std::vector<std::string>
const & read_data_name);
125 advance(
double const computed_timestep_length);
161 write_data(std::string
const & write_mesh_name,
162 std::string
const & write_data_name,
164 double const computed_timestep_length);
167 read_block_data(std::string
const & mesh_name,
const std::string & data_name)
const;
190#ifdef EXADG_WITH_PRECICE
191 std::shared_ptr<precice::SolverInterface> precice;
195 std::map<std::string, std::shared_ptr<CouplingBase<dim, data_dim, VectorizedArrayType>>> writer;
199 std::map<std::string, std::shared_ptr<ExaDGCoupling<dim, data_dim, VectorizedArrayType>>> reader;
202 std::vector<VectorType> old_state_data;
207template<
int dim,
int data_dim,
typename VectorType,
typename VectorizedArrayType>
208template<
typename ParameterClass>
213#ifdef EXADG_WITH_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));
220 AssertThrow(dim == precice->getDimensions(), dealii::ExcInternalError());
225 dealii::ExcMessage(
"EXADG_WITH_PRECICE has to be activated to use this code."));
228 AssertThrow(dim > 1, dealii::ExcNotImplemented());
233template<
int dim,
int data_dim,
typename VectorType,
typename VectorizedArrayType>
235Adapter<dim, data_dim, VectorType, VectorizedArrayType>::add_write_surface(
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)
245 auto const found_reader = reader.find(mesh_name);
247 if(found_reader != reader.end())
250 writer.insert({mesh_name, found_reader->second});
252 else if(write_data_type == WriteDataType::values_on_dofs)
256 std::make_shared<DoFCoupling<dim, data_dim, VectorizedArrayType>>(data,
257#ifdef EXADG_WITH_PRECICE
261 dealii_boundary_surface_id,
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>>(
272#ifdef EXADG_WITH_PRECICE
276 dealii_boundary_surface_id,
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);
286 writer.at(mesh_name)->set_write_data_type(write_data_type);
288 writer.at(mesh_name)->define_coupling_mesh();
293template<
int dim,
int data_dim,
typename VectorType,
typename VectorizedArrayType>
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)
303 std::make_shared<ExaDGCoupling<dim, data_dim, VectorizedArrayType>>(data,
304#ifdef EXADG_WITH_PRECICE
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();
317template<
int dim,
int data_dim,
typename VectorType,
typename VectorizedArrayType>
322#ifdef EXADG_WITH_PRECICE
327 precice->initialize();
334 if(precice->isActionRequired(precice::constants::actionWriteInitialData()))
336 writer[0]->write_data(dealii_to_precice,
"");
338 precice->markActionFulfilled(precice::constants::actionWriteInitialData());
340 precice->initializeData();
350 (void)dealii_to_precice;
354template<
int dim,
int data_dim,
typename VectorType,
typename VectorizedArrayType>
356Adapter<dim, data_dim, VectorType, VectorizedArrayType>::write_data(
357 std::string
const & write_mesh_name,
358 std::string
const & write_data_name,
360 double const computed_timestep_length)
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);
366 (void)write_mesh_name;
367 (void)write_data_name;
368 (void)dealii_to_precice;
369 (void)computed_timestep_length;
373template<
int dim,
int data_dim,
typename VectorType,
typename VectorizedArrayType>
376 double const computed_timestep_length)
378#ifdef EXADG_WITH_PRECICE
382 precice->advance(computed_timestep_length);
384 (void)computed_timestep_length;
390template<
int dim,
int data_dim,
typename VectorType,
typename VectorizedArrayType>
392Adapter<dim, data_dim, VectorType, VectorizedArrayType>::read_block_data(
393 std::string
const & mesh_name,
394 std::string
const & data_name)
const
396 reader.at(mesh_name)->read_block_data(data_name);
401template<
int dim,
int data_dim,
typename VectorType,
typename VectorizedArrayType>
404 std::function<
void()>
const & save_state)
406#ifdef EXADG_WITH_PRECICE
409 if(precice->isActionRequired(precice::constants::actionWriteIterationCheckpoint()))
412 precice->markActionFulfilled(precice::constants::actionWriteIterationCheckpoint());
421template<
int dim,
int data_dim,
typename VectorType,
typename VectorizedArrayType>
424 std::function<
void()>
const & reload_old_state)
426#ifdef EXADG_WITH_PRECICE
429 if(precice->isActionRequired(precice::constants::actionReadIterationCheckpoint()))
432 precice->markActionFulfilled(precice::constants::actionReadIterationCheckpoint());
435 (void)reload_old_state;
441template<
int dim,
int data_dim,
typename VectorType,
typename VectorizedArrayType>
445#ifdef EXADG_WITH_PRECICE
446 return precice->isCouplingOngoing();
454template<
int dim,
int data_dim,
typename VectorType,
typename VectorizedArrayType>
458#ifdef EXADG_WITH_PRECICE
459 return precice->isTimeWindowComplete();
Definition container_interface_data.h:45
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 ¶meters, 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