Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
cell_data_transfer.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2018 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_distributed_cell_data_transfer_h
17#define dealii_distributed_cell_data_transfer_h
18
19#include <deal.II/base/config.h>
20
22
24
26
27#include <boost/range/iterator_range.hpp>
28
29#include <functional>
30
31
33
34namespace parallel
35{
36 namespace distributed
37 {
126 template <int dim, int spacedim = dim, typename VectorType = Vector<double>>
128 {
129 private:
133 using value_type = typename VectorType::value_type;
134
135 public:
148 {
154 DEAL_II_DEPRECATED static typename VectorType::value_type
157 const VectorType &input_vector)
158 {
159 const typename VectorType::value_type value =
160 input_vector[parent->child(0)->active_cell_index()];
161
162 for (unsigned int child_index = 1; child_index < parent->n_children();
163 ++child_index)
164 Assert(
165 value ==
166 input_vector[parent->child(child_index)->active_cell_index()],
168 "Values on cells that will be coarsened are not equal!"));
169
170 return value;
171 }
172
178 DEAL_II_DEPRECATED static typename VectorType::value_type
180 cell_iterator & parent,
181 const VectorType &input_vector)
182 {
183 typename VectorType::value_type sum = 0;
184
185 for (unsigned int child_index = 0; child_index < parent->n_children();
186 ++child_index)
187 sum +=
188 input_vector[parent->child(child_index)->active_cell_index()];
189
190 return sum;
191 }
192
198 DEAL_II_DEPRECATED static typename VectorType::value_type
201 const VectorType & input_vector)
202 {
203 return sum(parent, input_vector) / parent->n_children();
204 }
205 };
206
225 const bool transfer_variable_size_data = false,
226 const std::function<std::vector<value_type>(
227 const typename ::Triangulation<dim, spacedim>::cell_iterator
228 & parent,
229 const value_type parent_value)> refinement_strategy =
230 &::AdaptationStrategies::Refinement::
231 preserve<dim, spacedim, value_type>,
232 const std::function<value_type(
233 const typename ::Triangulation<dim, spacedim>::cell_iterator
234 & parent,
235 const std::vector<value_type> &children_values)> coarsening_strategy =
236 &::AdaptationStrategies::Coarsening::
237 check_equality<dim, spacedim, value_type>);
238
258 const std::function<value_type(
259 const typename parallel::distributed::
261 const VectorType &input_vector)> coarsening_strategy);
262
275 void
277
281 void
283 const std::vector<const VectorType *> &all_in);
284
295 void
296 prepare_for_serialization(const VectorType &in);
297
301 void
302 prepare_for_serialization(const std::vector<const VectorType *> &all_in);
303
308 void
309 unpack(VectorType &out);
310
314 void
315 unpack(std::vector<VectorType *> &all_out);
316
321 void
322 deserialize(VectorType &out);
323
327 void
328 deserialize(std::vector<VectorType *> &all_out);
329
330 private:
337
342
347 const std::function<std::vector<value_type>(
348 const typename Triangulation<dim, spacedim>::cell_iterator &parent,
349 const value_type parent_value)>
351
356 const std::function<value_type(
357 const typename Triangulation<dim, spacedim>::cell_iterator &parent,
358 const std::vector<value_type> &children_values)>
360
365 std::vector<const VectorType *> input_vectors;
366
371 unsigned int handle;
372
377 void
379
385 std::vector<char>
388 const typename parallel::distributed::
390
396 void
399 cell_iterator &cell,
401 CellStatus status,
402 const boost::iterator_range<std::vector<char>::const_iterator>
403 & data_range,
404 std::vector<VectorType *> &all_out);
405 };
406 } // namespace distributed
407} // namespace parallel
408
409
411
412#endif /* dealii_distributed_cell_data_transfer_h */
void prepare_for_serialization(const VectorType &in)
void prepare_for_serialization(const std::vector< const VectorType * > &all_in)
SmartPointer< const parallel::distributed::Triangulation< dim, spacedim >, CellDataTransfer< dim, spacedim, VectorType > > triangulation
CellDataTransfer(const parallel::distributed::Triangulation< dim, spacedim > &triangulation, const bool transfer_variable_size_data, const std::function< value_type(const typename parallel::distributed::Triangulation< dim, spacedim >::cell_iterator &parent, const VectorType &input_vector)> coarsening_strategy)
void unpack(std::vector< VectorType * > &all_out)
void prepare_for_coarsening_and_refinement(const VectorType &in)
CellDataTransfer(const parallel::distributed::Triangulation< dim, spacedim > &triangulation, const bool transfer_variable_size_data=false, const std::function< std::vector< value_type >(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)> refinement_strategy=&::AdaptationStrategies::Refinement::preserve< dim, spacedim, value_type >, const std::function< value_type(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const std::vector< value_type > &children_values)> coarsening_strategy=&::AdaptationStrategies::Coarsening::check_equality< dim, spacedim, value_type >)
void deserialize(std::vector< VectorType * > &all_out)
void prepare_for_coarsening_and_refinement(const std::vector< const VectorType * > &all_in)
const std::function< value_type(const typename Triangulation< dim, spacedim >::cell_iterator &parent, const std::vector< value_type > &children_values)> coarsening_strategy
const std::function< std::vector< value_type >(const typename Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)> refinement_strategy
std::vector< char > pack_callback(const typename parallel::distributed::Triangulation< dim, spacedim >::cell_iterator &cell, const typename parallel::distributed::Triangulation< dim, spacedim >::CellStatus status)
std::vector< const VectorType * > input_vectors
void unpack_callback(const typename parallel::distributed::Triangulation< dim, spacedim >::cell_iterator &cell, const typename parallel::distributed::Triangulation< dim, spacedim >::CellStatus status, const boost::iterator_range< std::vector< char >::const_iterator > &data_range, std::vector< VectorType * > &all_out)
typename VectorType::value_type value_type
#define DEAL_II_DEPRECATED
Definition: config.h:162
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcMessage(std::string arg1)
T sum(const T &t, const MPI_Comm &mpi_communicator)
static VectorType::value_type sum(const typename parallel::distributed::Triangulation< dim, spacedim >::cell_iterator &parent, const VectorType &input_vector)
static VectorType::value_type mean(const typename parallel::distributed::Triangulation< dim, spacedim >::cell_iterator &parent, const VectorType &input_vector)
static VectorType::value_type check_equality(const typename parallel::distributed::Triangulation< dim, spacedim >::cell_iterator &parent, const VectorType &input_vector)