22#ifndef INCLUDE_EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_QUAD_COUPLING_H_
23#define INCLUDE_EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_QUAD_COUPLING_H_
25#include <deal.II/matrix_free/fe_evaluation.h>
27#include <exadg/fluid_structure_interaction/precice/coupling_base.h>
40template<
int dim,
int data_dim,
typename VectorizedArrayType>
41class QuadCoupling :
public CouplingBase<dim, data_dim, VectorizedArrayType>
44 QuadCoupling(dealii::MatrixFree<dim, double, VectorizedArrayType>
const & data,
45#ifdef EXADG_WITH_PRECICE
46 std::shared_ptr<precice::SolverInterface> precice,
49 dealii::types::boundary_id
const surface_id,
50 int const mf_dof_index,
51 int const mf_quad_index);
56 using value_type =
typename CouplingBase<dim, data_dim, VectorizedArrayType>::value_type;
78 write_data(dealii::LinearAlgebra::distributed::Vector<double>
const & data_vector,
79 std::string
const & data_name)
override;
93 dealii::LinearAlgebra::distributed::Vector<double>
const & data_vector,
94 int const write_data_id,
95 dealii::EvaluationFlags::EvaluationFlags
const flags,
96 std::function<value_type(
FEFaceIntegrator &,
unsigned int)>
const & get_write_value);
99 std::vector<std::array<int, VectorizedArrayType::size()>> coupling_nodes_ids;
103 int const mf_dof_index;
104 int const mf_quad_index;
107 get_surface_type()
const override;
112template<
int dim,
int data_dim,
typename VectorizedArrayType>
113QuadCoupling<dim, data_dim, VectorizedArrayType>::QuadCoupling(
114 dealii::MatrixFree<dim, double, VectorizedArrayType>
const & data,
115#ifdef EXADG_WITH_PRECICE
116 std::shared_ptr<precice::SolverInterface> precice,
118 std::string
const mesh_name,
119 dealii::types::boundary_id
const surface_id,
120 int const mf_dof_index_,
121 int const mf_quad_index_)
122 :
CouplingBase<dim, data_dim, VectorizedArrayType>(data,
123#ifdef EXADG_WITH_PRECICE
128 mf_dof_index(mf_dof_index_),
129 mf_quad_index(mf_quad_index_)
135template<
int dim,
int data_dim,
typename VectorizedArrayType>
139 Assert(this->mesh_id != -1, dealii::ExcNotInitialized());
143 if(coupling_nodes_ids.size() > 0)
147 coupling_nodes_ids.reserve(this->
matrix_free.n_boundary_face_batches() / 2);
151 std::array<double, dim * VectorizedArrayType::size()> unrolled_vertices;
152 std::array<int, VectorizedArrayType::size()> node_ids;
153 unsigned int size = 0;
155 for(
unsigned int face = this->
matrix_free.n_inner_face_batches();
156 face < this->
matrix_free.n_boundary_face_batches() + this->matrix_free.n_inner_face_batches();
159 auto const boundary_id = this->
matrix_free.get_boundary_id(face);
162 if(boundary_id != this->dealii_boundary_surface_id)
166 int const active_faces = this->
matrix_free.n_active_entries_per_face_batch(face);
169 for(
unsigned int q = 0; q < phi.n_q_points; ++q)
171 auto const local_vertex = phi.quadrature_point(q);
177 for(
int d = 0; d < dim; ++d)
178 for(
unsigned int v = 0; v < VectorizedArrayType::size(); ++v)
179 unrolled_vertices[d + dim * v] = local_vertex[d][v];
181#ifdef EXADG_WITH_PRECICE
182 this->precice->setMeshVertices(this->mesh_id,
184 unrolled_vertices.data(),
189 coupling_nodes_ids.emplace_back(node_ids);
195 coupling_nodes_ids.resize(size);
197#ifdef EXADG_WITH_PRECICE
200 Assert(size * VectorizedArrayType::size() >=
201 static_cast<unsigned int>(this->precice->getMeshVertexSize(this->mesh_id)),
202 dealii::ExcInternalError());
204 if(this->read_data_map.size() > 0)
205 this->
print_info(
true, this->precice->getMeshVertexSize(this->mesh_id));
206 if(this->write_data_map.size() > 0)
207 this->
print_info(
false, this->precice->getMeshVertexSize(this->mesh_id));
212template<
int dim,
int data_dim,
typename VectorizedArrayType>
215 dealii::LinearAlgebra::distributed::Vector<double>
const & data_vector,
216 std::string
const & data_name)
218 int const write_data_id = this->write_data_map.at(data_name);
220 switch(this->write_data_type)
222 case WriteDataType::values_on_q_points:
223 write_data_factory(data_vector,
225 dealii::EvaluationFlags::values,
226 [](
auto & phi,
auto q_point) {
return phi.get_value(q_point); });
228 case WriteDataType::normal_gradients_on_q_points:
229 write_data_factory(data_vector,
231 dealii::EvaluationFlags::gradients,
232 [](
auto & phi,
auto q_point) {
233 return phi.get_normal_derivative(q_point);
237 AssertThrow(
false, dealii::ExcNotImplemented());
243template<
int dim,
int data_dim,
typename VectorizedArrayType>
245QuadCoupling<dim, data_dim, VectorizedArrayType>::write_data_factory(
246 dealii::LinearAlgebra::distributed::Vector<double>
const & data_vector,
247 int const write_data_id,
248 dealii::EvaluationFlags::EvaluationFlags
const flags,
249 std::function<value_type(FEFaceIntegrator &,
unsigned int)>
const & get_write_value)
251 Assert(write_data_id != -1, dealii::ExcNotInitialized());
252 Assert(coupling_nodes_ids.size() > 0, dealii::ExcNotInitialized());
254 FEFaceIntegrator phi(this->matrix_free,
true, mf_dof_index, mf_quad_index);
257 std::array<double, data_dim * VectorizedArrayType::size()> unrolled_local_data;
258 (void)unrolled_local_data;
260 auto index = coupling_nodes_ids.begin();
263 for(
unsigned int face = this->matrix_free.n_inner_face_batches();
264 face < this->matrix_free.n_boundary_face_batches() + this->matrix_free.n_inner_face_batches();
267 auto const boundary_id = this->matrix_free.get_boundary_id(face);
270 if(boundary_id != this->dealii_boundary_surface_id)
275 phi.read_dof_values_plain(data_vector);
277 int const active_faces = this->matrix_free.n_active_entries_per_face_batch(face);
279 for(
unsigned int q = 0; q < phi.n_q_points; ++q, ++index)
281 Assert(index != coupling_nodes_ids.end(), dealii::ExcInternalError());
282 auto const local_data = get_write_value(phi, q);
286 if constexpr(data_dim > 1)
290 for(
int d = 0; d < data_dim; ++d)
291 for(
unsigned int v = 0; v < VectorizedArrayType::size(); ++v)
292 unrolled_local_data[d + data_dim * v] = local_data[d][v];
294#ifdef EXADG_WITH_PRECICE
295 this->precice->writeBlockVectorData(write_data_id,
298 unrolled_local_data.data());
306#ifdef EXADG_WITH_PRECICE
307 this->precice->writeBlockScalarData(write_data_id,
319template<
int dim,
int data_dim,
typename VectorizedArrayType>
321QuadCoupling<dim, data_dim, VectorizedArrayType>::get_surface_type()
const
323 return "quadrature points using matrix-free quad index " +
324 dealii::Utilities::to_string(mf_quad_index);
Definition coupling_base.h:66
std::string const mesh_name
public precice solverinterface
Definition coupling_base.h:152
FaceIntegrator< dim, data_dim, double, VectorizedArrayType > FEFaceIntegrator
Alias for the face integrator.
Definition coupling_base.h:78
dealii::MatrixFree< dim, double, VectorizedArrayType > const & matrix_free
The dealii::MatrixFree object (preCICE can only handle double precision)
Definition coupling_base.h:144
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 quad_coupling.h:137
typename CouplingBase< dim, data_dim, VectorizedArrayType >::FEFaceIntegrator FEFaceIntegrator
Alias as defined in the base class.
Definition quad_coupling.h:54
virtual void write_data(dealii::LinearAlgebra::distributed::Vector< double > const &data_vector, std::string const &data_name) override
write_data Evaluates the given
Definition quad_coupling.h:214