45class DoFCoupling :
public CouplingBase<dim, data_dim, VectorizedArrayType>
48 DoFCoupling(dealii::MatrixFree<dim, double, VectorizedArrayType>
const & data,
49#ifdef EXADG_WITH_PRECICE
50 std::shared_ptr<precice::SolverInterface> precice,
53 dealii::types::boundary_id
const surface_id,
54 int const mf_dof_index);
71 write_data(dealii::LinearAlgebra::distributed::Vector<double>
const & data_vector,
72 std::string
const & data_name)
override;
77 std::vector<int> coupling_nodes_ids;
79 std::vector<std::array<dealii::types::global_dof_index, data_dim>> global_indices;
83 int const mf_dof_index;
86 get_surface_type()
const override;
116 Assert(this->mesh_id != -1, dealii::ExcNotInitialized());
120 if(coupling_nodes_ids.size() > 0)
124 auto const get_component_dofs = [&](
int const component) {
126 dealii::ComponentMask component_mask(data_dim,
false);
127 component_mask.set(component,
true);
133 dealii::IndexSet
const indices =
134 (dealii::DoFTools::extract_boundary_dofs(this->
matrix_free.get_dof_handler(mf_dof_index),
136 std::set<dealii::types::boundary_id>{
137 this->dealii_boundary_surface_id}) &
138 this->
matrix_free.get_dof_handler(mf_dof_index).locally_owned_dofs());
140 Assert(indices.n_elements() * data_dim ==
141 this->matrix_free.get_dof_handler(mf_dof_index)
143 std::set<dealii::types::boundary_id>{this->dealii_boundary_surface_id}),
144 dealii::ExcInternalError());
148 global_indices.resize(indices.n_elements());
150 dealii::types::global_dof_index iterator = 0;
151 for(
auto const dof : indices)
153 global_indices[iterator][component] = dof;
159 for(
int d = 0; d < data_dim; ++d)
160 get_component_dofs(d);
163 std::map<dealii::types::global_dof_index, dealii::Point<dim>> support_points;
164 dealii::ComponentMask component_mask(data_dim,
false);
165 component_mask.set(0,
true);
167 dealii::DoFTools::map_boundary_dofs_to_support_points(
169 this->matrix_free.get_dof_handler(mf_dof_index),
172 this->dealii_boundary_surface_id);
176 coupling_nodes_ids.reserve(global_indices.size());
177 std::array<double, dim> nodes_position;
178 for(std::size_t i = 0; i < global_indices.size(); ++i)
181 auto const element = global_indices[i][0];
182 for(
int d = 0; d < dim; ++d)
183 nodes_position[d] = support_points[element][d];
186#ifdef EXADG_WITH_PRECICE
187 int const precice_id = this->precice->setMeshVertex(this->mesh_id, nodes_position.data());
189 int const precice_id = 0;
191 coupling_nodes_ids.emplace_back(precice_id);
194#ifdef EXADG_WITH_PRECICE
195 if(this->read_data_map.size() > 0)
196 this->
print_info(
true, this->precice->getMeshVertexSize(this->mesh_id));
197 if(this->write_data_map.size() > 0)
198 this->
print_info(
false, this->precice->getMeshVertexSize(this->mesh_id));
207 dealii::LinearAlgebra::distributed::Vector<double>
const & data_vector,
208 std::string
const & data_name)
210 int const write_data_id = this->write_data_map.at(data_name);
211 Assert(write_data_id != -1, dealii::ExcNotInitialized());
212 Assert(coupling_nodes_ids.size() > 0, dealii::ExcNotInitialized());
215 for(std::size_t i = 0; i < global_indices.size(); ++i)
218 for(
unsigned int d = 0; d < data_dim; ++d)
220 auto const element = global_indices[i][d];
224#ifdef EXADG_WITH_PRECICE
226 if constexpr(data_dim > 1)
228 this->precice->writeVectorData(write_data_id, coupling_nodes_ids[i],
write_data.data());
232 this->precice->writeScalarData(write_data_id, coupling_nodes_ids[i],
write_data[0]);
dealii::MatrixFree< dim, double, VectorizedArrayType > const & matrix_free
The dealii::MatrixFree object (preCICE can only handle double precision)
Definition coupling_base.h:143
void print_info(bool const reader, unsigned int const local_size) const
Print information of the current setup.
Definition coupling_base.h:255
virtual void write_data(dealii::LinearAlgebra::distributed::Vector< double > const &data_vector, std::string const &data_name) override
write_data Evaluates the given
Definition dof_coupling.h:206