Reference documentation for deal.II version 9.2.0
\(\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 - 2020 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 
34 namespace parallel
35 {
36  namespace distributed
37  {
127  template <int dim, int spacedim = dim, typename VectorType = Vector<double>>
129  {
130  private:
134  using value_type = typename VectorType::value_type;
135 
136  public:
149  {
155  DEAL_II_DEPRECATED static typename VectorType::value_type
158  const VectorType &input_vector)
159  {
160  const typename VectorType::value_type value =
161  input_vector[parent->child(0)->active_cell_index()];
162 
163  for (unsigned int child_index = 1; child_index < parent->n_children();
164  ++child_index)
165  Assert(
166  value ==
167  input_vector[parent->child(child_index)->active_cell_index()],
168  ExcMessage(
169  "Values on cells that will be coarsened are not equal!"));
170 
171  return value;
172  }
173 
179  DEAL_II_DEPRECATED static typename VectorType::value_type
181  cell_iterator & parent,
182  const VectorType &input_vector)
183  {
184  typename VectorType::value_type sum = 0;
185 
186  for (unsigned int child_index = 0; child_index < parent->n_children();
187  ++child_index)
188  sum +=
189  input_vector[parent->child(child_index)->active_cell_index()];
190 
191  return sum;
192  }
193 
199  DEAL_II_DEPRECATED static typename VectorType::value_type
200  mean(const typename parallel::distributed::
202  const VectorType & input_vector)
203  {
204  return sum(parent, input_vector) / parent->n_children();
205  }
206  };
207 
225  & triangulation,
226  const bool transfer_variable_size_data = false,
227  const std::function<std::vector<value_type>(
228  const typename ::Triangulation<dim, spacedim>::cell_iterator
229  & parent,
230  const value_type parent_value)> refinement_strategy =
232  preserve<dim, spacedim, value_type>,
233  const std::function<value_type(
234  const typename ::Triangulation<dim, spacedim>::cell_iterator
235  & parent,
236  const std::vector<value_type> &children_values)> coarsening_strategy =
238  check_equality<dim, spacedim, value_type>);
239 
257  & triangulation,
258  const bool transfer_variable_size_data,
259  const std::function<value_type(
260  const typename parallel::distributed::
262  const VectorType &input_vector)> coarsening_strategy);
263 
276  void
278 
282  void
284  const std::vector<const VectorType *> &all_in);
285 
296  void
298 
302  void
303  prepare_for_serialization(const std::vector<const VectorType *> &all_in);
304 
309  void
310  unpack(VectorType &out);
311 
315  void
316  unpack(std::vector<VectorType *> &all_out);
317 
322  void
323  deserialize(VectorType &out);
324 
328  void
329  deserialize(std::vector<VectorType *> &all_out);
330 
331  private:
338 
343 
348  const std::function<std::vector<value_type>(
349  const typename Triangulation<dim, spacedim>::cell_iterator &parent,
350  const value_type parent_value)>
352 
357  const std::function<value_type(
358  const typename Triangulation<dim, spacedim>::cell_iterator &parent,
359  const std::vector<value_type> &children_values)>
361 
366  std::vector<const VectorType *> input_vectors;
367 
372  unsigned int handle;
373 
378  void
380 
386  std::vector<char>
387  pack_callback(const typename parallel::distributed::
389  const typename parallel::distributed::
391 
397  void
400  cell_iterator &cell,
402  CellStatus status,
403  const boost::iterator_range<std::vector<char>::const_iterator>
404  & data_range,
405  std::vector<VectorType *> &all_out);
406  };
407  } // namespace distributed
408 } // namespace parallel
409 
410 
412 
413 #endif /* dealii_distributed_cell_data_transfer_h */
parallel::distributed::CellDataTransfer::CoarseningStrategies::check_equality
static VectorType::value_type check_equality(const typename parallel::distributed::Triangulation< dim, spacedim >::cell_iterator &parent, const VectorType &input_vector)
Definition: cell_data_transfer.h:156
parallel::distributed::CellDataTransfer::coarsening_strategy
const std::function< value_type(const typename Triangulation< dim, spacedim >::cell_iterator &parent, const std::vector< value_type > &children_values)> coarsening_strategy
Definition: cell_data_transfer.h:360
parallel::distributed::CellDataTransfer::triangulation
SmartPointer< const parallel::distributed::Triangulation< dim, spacedim >, CellDataTransfer< dim, spacedim, VectorType > > triangulation
Definition: cell_data_transfer.h:337
adaptation_strategies.h
parallel::distributed::CellDataTransfer::deserialize
void deserialize(VectorType &out)
parallel::distributed::CellDataTransfer::unpack
void unpack(VectorType &out)
parallel::distributed::CellDataTransfer::CellDataTransfer
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 >)
parallel::distributed::CellDataTransfer::transfer_variable_size_data
const bool transfer_variable_size_data
Definition: cell_data_transfer.h:342
VectorType
parallel::distributed::CellDataTransfer::input_vectors
std::vector< const VectorType * > input_vectors
Definition: cell_data_transfer.h:366
parallel::distributed::CellDataTransfer< dim, spacedim, std::vector< unsigned int > >::value_type
typename std::vector< unsigned int > ::value_type value_type
Definition: cell_data_transfer.h:134
parallel::distributed::CellDataTransfer::prepare_for_coarsening_and_refinement
void prepare_for_coarsening_and_refinement(const VectorType &in)
parallel::distributed::CellDataTransfer::unpack_callback
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)
Triangulation< dim, dim >::CellStatus
CellStatus
Definition: tria.h:1969
parallel::distributed::CellDataTransfer::register_data_attach
void register_data_attach()
AdaptationStrategies::Coarsening
Definition: adaptation_strategies.h:120
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
parallel::distributed::CellDataTransfer::pack_callback
std::vector< char > pack_callback(const typename parallel::distributed::Triangulation< dim, spacedim >::cell_iterator &cell, const typename parallel::distributed::Triangulation< dim, spacedim >::CellStatus status)
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
parallel::distributed::CellDataTransfer::CoarseningStrategies
Definition: cell_data_transfer.h:148
parallel::distributed::Triangulation
Definition: tria.h:242
AdaptationStrategies::Refinement
Definition: adaptation_strategies.h:55
DEAL_II_DEPRECATED
#define DEAL_II_DEPRECATED
Definition: config.h:98
Utilities::MPI::sum
T sum(const T &t, const MPI_Comm &mpi_communicator)
smartpointer.h
parallel::distributed::CellDataTransfer::prepare_for_serialization
void prepare_for_serialization(const VectorType &in)
parallel::distributed::CellDataTransfer::CoarseningStrategies::sum
static VectorType::value_type sum(const typename parallel::distributed::Triangulation< dim, spacedim >::cell_iterator &parent, const VectorType &input_vector)
Definition: cell_data_transfer.h:180
parallel::distributed
Definition: distributed.h:429
parallel::distributed::CellDataTransfer
Definition: cell_data_transfer.h:128
SmartPointer
Definition: smartpointer.h:68
value
static const bool value
Definition: dof_tools_constraints.cc:433
tria.h
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
parallel::distributed::CellDataTransfer::handle
unsigned int handle
Definition: cell_data_transfer.h:372
config.h
parallel::distributed::CellDataTransfer::refinement_strategy
const std::function< std::vector< value_type > const typename Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)> refinement_strategy
Definition: cell_data_transfer.h:351
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
value_type
parallel::distributed::CellDataTransfer::CoarseningStrategies::mean
static VectorType::value_type mean(const typename parallel::distributed::Triangulation< dim, spacedim >::cell_iterator &parent, const VectorType &input_vector)
Definition: cell_data_transfer.h:200
TriaIterator
Definition: tria_iterator.h:578
parallel
Definition: distributed.h:416