22#ifndef INCLUDE_EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_EXADG_COUPLING_H_
23#define INCLUDE_EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_EXADG_COUPLING_H_
26#include <exadg/fluid_structure_interaction/precice/coupling_base.h>
27#include <exadg/functions_and_boundary_conditions/interface_coupling.h>
28#include <exadg/utilities/tensor_utilities.h>
39template<
int dim,
int data_dim,
typename VectorizedArrayType>
40class ExaDGCoupling :
public CouplingBase<dim, data_dim, VectorizedArrayType>
43 static unsigned int const rank = n_components_to_rank<data_dim, dim>();
46 dealii::MatrixFree<dim, double, VectorizedArrayType>
const & data,
47#ifdef EXADG_WITH_PRECICE
48 std::shared_ptr<precice::SolverInterface> precice,
52 dealii::types::boundary_id
const surface_id = dealii::numbers::invalid_unsigned_int);
70 write_data(dealii::LinearAlgebra::distributed::Vector<double>
const & data_vector,
71 std::string
const & data_name)
override;
74 read_block_data(std::string
const & data_name)
const override;
78 std::shared_ptr<ContainerInterfaceData<rank, dim, double>> interface_data;
81 std::vector<int> coupling_nodes_ids;
84 get_surface_type()
const override;
40class ExaDGCoupling :
public CouplingBase<dim, data_dim, VectorizedArrayType> {
…};
89template<
int dim,
int data_dim,
typename VectorizedArrayType>
90ExaDGCoupling<dim, data_dim, VectorizedArrayType>::ExaDGCoupling(
91 dealii::MatrixFree<dim, double, VectorizedArrayType>
const & data,
92#ifdef EXADG_WITH_PRECICE
93 std::shared_ptr<precice::SolverInterface> precice,
95 std::string
const mesh_name,
97 dealii::types::boundary_id
const surface_id)
99#ifdef EXADG_WITH_PRECICE
104 interface_data(interface_data_)
110template<
int dim,
int data_dim,
typename VectorizedArrayType>
114 Assert(this->mesh_id != -1, dealii::ExcNotInitialized());
115 Assert(interface_data.get() !=
nullptr, dealii::ExcNotInitialized());
119 if(coupling_nodes_ids.size() > 0)
122 for(
auto quad_index : interface_data->get_quad_indices())
125 auto const & points = interface_data->get_array_q_points(quad_index);
128 auto const start_index = coupling_nodes_ids.size();
131 coupling_nodes_ids.resize(start_index + points.size());
134#ifdef EXADG_WITH_PRECICE
135 this->precice->setMeshVertices(this->mesh_id,
138 &coupling_nodes_ids[start_index]);
142#ifdef EXADG_WITH_PRECICE
143 if(this->read_data_map.size() > 0)
144 this->
print_info(
true, this->precice->getMeshVertexSize(this->mesh_id));
145 if(this->write_data_map.size() > 0)
146 this->
print_info(
false, this->precice->getMeshVertexSize(this->mesh_id));
151template<
int dim,
int data_dim,
typename VectorizedArrayType>
153ExaDGCoupling<dim, data_dim, VectorizedArrayType>::read_block_data(
154 std::string
const & data_name)
const
156 Assert(interface_data.get() !=
nullptr, dealii::ExcNotInitialized());
158 int const read_data_id = this->read_data_map.at(data_name);
161 unsigned int start_index = 0;
163 for(
auto quad_index : interface_data->get_quad_indices())
165 if constexpr(data_dim > 1)
167 auto & array_solution = interface_data->get_array_solution(quad_index);
168 auto const array_size = array_solution.size();
170 AssertIndexRange(start_index, coupling_nodes_ids.size());
171#ifdef EXADG_WITH_PRECICE
172 this->precice->readBlockVectorData(read_data_id,
174 &coupling_nodes_ids[start_index],
175 &array_solution[0][0]);
179 start_index += array_size;
183 AssertThrow(
false, dealii::ExcNotImplemented());
190template<
int dim,
int data_dim,
typename VectorizedArrayType>
193 dealii::LinearAlgebra::distributed::Vector<double>
const &,
196 AssertThrow(
false, dealii::ExcNotImplemented());
200template<
int dim,
int data_dim,
typename VectorizedArrayType>
202ExaDGCoupling<dim, data_dim, VectorizedArrayType>::get_surface_type()
const
204 return "exadg shallow wrapper";
Definition container_interface_data.h:45
Definition coupling_base.h:66
std::string const mesh_name
public precice solverinterface
Definition coupling_base.h:152
void print_info(bool const reader, unsigned int const local_size) const
Print information of the current setup.
Definition coupling_base.h:256
virtual void define_coupling_mesh() override
define_mesh_vertices Define a vertex coupling mesh for preCICE coupling the classical preCICE way
Definition exadg_coupling.h:112
virtual void write_data(dealii::LinearAlgebra::distributed::Vector< double > const &data_vector, std::string const &data_name) override
write_data
Definition exadg_coupling.h:192