22#ifndef EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_COUPLING_BASE_H_
23#define EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_COUPLING_BASE_H_
26#include <deal.II/base/conditional_ostream.h>
27#include <deal.II/base/vectorization.h>
28#include <deal.II/matrix_free/matrix_free.h>
31#include <exadg/matrix_free/integrators.h>
34#ifdef EXADG_WITH_PRECICE
35# include <precice/SolverInterface.hpp>
45enum class WriteDataType
50 gradients_on_other_mesh,
52 normal_gradients_on_q_points
63template<
int dim,
int data_dim,
typename VectorizedArrayType>
67 CouplingBase(dealii::MatrixFree<dim, double, VectorizedArrayType>
const & data,
68#ifdef EXADG_WITH_PRECICE
69 std::shared_ptr<precice::SolverInterface> precice,
72 dealii::types::boundary_id
const surface_id);
74 virtual ~CouplingBase() =
default;
78 using value_type =
typename FEFaceIntegrator::value_type;
103 write_data(dealii::LinearAlgebra::distributed::Vector<double>
const & data_vector,
104 std::string
const & data_name) = 0;
108 read_block_data(std::string
const & data_name)
const;
140 print_info(
bool const reader,
unsigned int const local_size)
const;
143 dealii::MatrixFree<dim, double, VectorizedArrayType>
const &
matrix_free;
146#ifdef EXADG_WITH_PRECICE
147 std::shared_ptr<precice::SolverInterface> precice;
154 std::map<std::string, int> read_data_map;
155 std::map<std::string, int> write_data_map;
157 dealii::types::boundary_id
const dealii_boundary_surface_id;
159 WriteDataType write_data_type;
162 get_surface_type()
const = 0;
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,
173 std::string
const mesh_name,
174 dealii::types::boundary_id
const surface_id)
175 : matrix_free(matrix_free_),
176#ifdef EXADG_WITH_PRECICE
179 mesh_name(mesh_name),
180 dealii_boundary_surface_id(surface_id),
181 write_data_type(WriteDataType::undefined)
183#ifdef EXADG_WITH_PRECICE
184 Assert(precice.get() !=
nullptr, dealii::ExcNotInitialized());
190 dealii::ExcMessage(
"EXADG_WITH_PRECICE has to be activated to use this code."));
196template<
int dim,
int data_dim,
typename VectorizedArrayType>
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);
204 int const read_data_id = 0;
206 read_data_map.insert({read_data_name, read_data_id});
211template<
int dim,
int data_dim,
typename VectorizedArrayType>
214 std::string
const & write_data_name)
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);
220 int const write_data_id = 0;
222 write_data_map.insert({write_data_name, write_data_id});
226template<
int dim,
int data_dim,
typename VectorizedArrayType>
229 WriteDataType write_data_specification)
231 write_data_type = write_data_specification;
236template<
int dim,
int data_dim,
typename VectorizedArrayType>
245template<
int dim,
int data_dim,
typename VectorizedArrayType>
247CouplingBase<dim, data_dim, VectorizedArrayType>::read_block_data(std::string
const &)
const
249 AssertThrow(
false, dealii::ExcNotImplemented());
253template<
int dim,
int data_dim,
typename VectorizedArrayType>
256 unsigned int const local_size)
const
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);
263 auto names = map.begin();
264 std::string data_names = names->first;
266 for(; names != map.end(); ++names)
267 data_names += std::string(
", ") + names->first;
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())
276 <<
"-- . Node location: " << get_surface_type() <<
"\n"
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