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\}}\)
error_predictor.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2019 - 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 
17 #include <deal.II/base/config.h>
18 
19 #ifdef DEAL_II_WITH_P4EST
20 
23 
25 # include <deal.II/dofs/dof_tools.h>
26 
29 
30 # include <deal.II/hp/dof_handler.h>
31 
39 # include <deal.II/lac/vector.h>
40 
41 # include <functional>
42 # include <numeric>
43 
44 
46 
47 
48 namespace parallel
49 {
50  namespace distributed
51  {
52  template <int dim, int spacedim>
55  : dof_handler(&dof, typeid(*this).name())
56  , handle(numbers::invalid_unsigned_int)
57  , gamma_p(0.)
58  , gamma_h(0.)
59  , gamma_n(0.)
60  {
61  Assert(
63  *>(&dof_handler->get_triangulation()) != nullptr),
64  ExcMessage(
65  "parallel::distributed::ErrorPredictor requires a parallel::distributed::Triangulation object."));
66  }
67 
68 
69 
70  template <int dim, int spacedim>
71  void
73  const std::vector<const Vector<float> *> &all_in,
74  const double gamma_p,
75  const double gamma_h,
76  const double gamma_n)
77  {
78  error_indicators = all_in;
79  this->gamma_p = gamma_p;
80  this->gamma_h = gamma_h;
81  this->gamma_n = gamma_n;
82  register_data_attach();
83  }
84 
85 
86 
87  template <int dim, int spacedim>
88  void
90  {
91  // TODO: casting away constness is bad
94  const_cast<::Triangulation<dim, spacedim> *>(
95  &dof_handler->get_triangulation())));
96  Assert(tria != nullptr, ExcInternalError());
97 
98  handle = tria->register_data_attach(
99  [this](const typename Triangulation<dim, spacedim>::cell_iterator &cell,
100  const typename Triangulation<dim, spacedim>::CellStatus status) {
101  return this->pack_callback(cell, status);
102  },
103  /*returns_variable_size_data=*/false);
104  }
105 
106 
107 
108  template <int dim, int spacedim>
109  void
111  const Vector<float> &in,
112  const double gamma_p,
113  const double gamma_h,
114  const double gamma_n)
115  {
116  std::vector<const Vector<float> *> all_in(1, &in);
117  prepare_for_coarsening_and_refinement(all_in, gamma_p, gamma_h, gamma_n);
118  }
119 
120 
121 
122  template <int dim, int spacedim>
123  void
125  {
126  Assert(error_indicators.size() == all_out.size(),
127  ExcDimensionMismatch(error_indicators.size(), all_out.size()));
128 
129  // TODO: casting away constness is bad
132  const_cast<::Triangulation<dim, spacedim> *>(
133  &dof_handler->get_triangulation())));
134  Assert(tria != nullptr, ExcInternalError());
135 
137  handle,
138  [this, &all_out](
139  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
140  const typename Triangulation<dim, spacedim>::CellStatus status,
141  const boost::iterator_range<std::vector<char>::const_iterator>
142  &data_range) {
143  this->unpack_callback(cell, status, data_range, all_out);
144  });
145 
146  for (const auto &out : all_out)
147  out->compress(::VectorOperation::insert);
148 
149  error_indicators.clear();
150  }
151 
152 
153 
154  template <int dim, int spacedim>
155  void
157  {
158  std::vector<Vector<float> *> all_out(1, &out);
159  unpack(all_out);
160  }
161 
162 
163 
164  template <int dim, int spacedim>
165  std::vector<char>
167  const typename Triangulation<dim, spacedim>::cell_iterator &cell_,
168  const typename Triangulation<dim, spacedim>::CellStatus status)
169  {
171  dof_handler);
172 
173  // create buffer for each individual input vector
174  std::vector<float> predicted_errors(error_indicators.size());
175 
176  auto predicted_error_it = predicted_errors.begin();
177  auto estimated_error_it = error_indicators.cbegin();
178  for (; estimated_error_it != error_indicators.cend();
179  ++predicted_error_it, ++estimated_error_it)
180  switch (status)
181  {
183  spacedim>::CELL_PERSIST:
184  {
185  if (cell->future_fe_index_set())
186  {
187  const int degree_difference =
188  dof_handler->get_fe_collection()[cell->future_fe_index()]
189  .degree -
190  cell->get_fe().degree;
191 
192  *predicted_error_it =
193  (**estimated_error_it)[cell->active_cell_index()] *
194  std::pow(gamma_p, degree_difference);
195  }
196  else
197  {
198  *predicted_error_it =
199  (**estimated_error_it)[cell->active_cell_index()] *
200  gamma_n;
201  }
202  break;
203  }
204 
206  spacedim>::CELL_REFINE:
207  {
208  // Determine the exponent by the finite element degree on the
209  // adapted mesh.
210  const unsigned int future_fe_degree =
211  dof_handler->get_fe_collection()[cell->future_fe_index()]
212  .degree;
213 
214  *predicted_error_it =
215  (**estimated_error_it)[cell->active_cell_index()] *
216  (gamma_h * std::pow(.5, future_fe_degree));
217 
218  // If the future fe index differs from the active one, also take
219  // into account p-adaptation.
220  if (cell->future_fe_index_set())
221  *predicted_error_it *=
222  std::pow(gamma_p,
223  static_cast<int>(future_fe_degree -
224  cell->get_fe().degree));
225 
226  break;
227  }
228 
230  spacedim>::CELL_COARSEN:
231  {
232  // First figure out which finite element will be assigned to the
233  // parent cell after h-adaptation analogously to
234  // ::internal::hp::DoFHandlerImplementation::Implementation::
235  // collect_fe_indices_on_cells_to_be_refined()
236  std::set<unsigned int> fe_indices_children;
237  for (unsigned int child_index = 0;
238  child_index < cell->n_children();
239  ++child_index)
240  {
241  const auto child = cell->child(child_index);
242  Assert(child->is_active() && child->coarsen_flag_set(),
244  dim>::ExcInconsistentCoarseningFlags());
245 
246  fe_indices_children.insert(child->future_fe_index());
247  }
248 
249  const unsigned int future_fe_index =
250  dof_handler->get_fe_collection().find_dominated_fe_extended(
251  fe_indices_children, /*codim=*/0);
252 
253  const unsigned int future_fe_degree =
254  dof_handler->get_fe_collection()[future_fe_index].degree;
255 
256  // Then determine the actually contirbution to the predicted
257  // error of every single cells that is about to be coarsened.
258  float sqrsum_of_predicted_errors = 0.;
259  float predicted_error = 0.;
260  int degree_difference = 0;
261  for (unsigned int child_index = 0;
262  child_index < cell->n_children();
263  ++child_index)
264  {
265  const auto child = cell->child(child_index);
266 
267  predicted_error =
268  (**estimated_error_it)[child->active_cell_index()] /
269  (gamma_h * std::pow(.5, future_fe_degree));
270 
271  degree_difference =
272  future_fe_degree - child->get_fe().degree;
273 
274  if (degree_difference != 0)
275  predicted_error *= std::pow(gamma_p, degree_difference);
276 
277  sqrsum_of_predicted_errors +=
278  predicted_error * predicted_error;
279  }
280  *predicted_error_it = std::sqrt(sqrsum_of_predicted_errors);
281 
282  break;
283  }
284 
285  default:
286  Assert(false, ExcInternalError());
287  break;
288  }
289 
290  // We don't have to pack the whole container if there is just one entry.
291  if (error_indicators.size() == 1)
292  return Utilities::pack(predicted_errors[0],
293  /*allow_compression=*/false);
294  else
295  return Utilities::pack(predicted_errors, /*allow_compression=*/false);
296  }
297 
298 
299 
300  template <int dim, int spacedim>
301  void
303  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
304  const typename Triangulation<dim, spacedim>::CellStatus status,
305  const boost::iterator_range<std::vector<char>::const_iterator>
306  & data_range,
307  std::vector<Vector<float> *> &all_out)
308  {
309  std::vector<float> predicted_errors;
310 
311  if (all_out.size() == 1)
312  predicted_errors.push_back(
313  Utilities::unpack<float>(data_range.begin(),
314  data_range.end(),
315  /*allow_compression=*/false));
316  else
317  predicted_errors =
318  Utilities::unpack<std::vector<float>>(data_range.begin(),
319  data_range.end(),
320  /*allow_compression=*/false);
321 
322  Assert(predicted_errors.size() == all_out.size(), ExcInternalError());
323 
324  auto it_input = predicted_errors.cbegin();
325  auto it_output = all_out.begin();
326  for (; it_input != predicted_errors.cend(); ++it_input, ++it_output)
327  switch (status)
328  {
330  spacedim>::CELL_PERSIST:
332  spacedim>::CELL_COARSEN:
333  (**it_output)[cell->active_cell_index()] = *it_input;
334  break;
335 
336 
338  spacedim>::CELL_REFINE:
339  for (unsigned int child_index = 0;
340  child_index < cell->n_children();
341  ++child_index)
342  (**it_output)[cell->child(child_index)->active_cell_index()] =
343  (*it_input) / std::sqrt(cell->n_children());
344  break;
345 
346  default:
347  Assert(false, ExcInternalError());
348  break;
349  }
350  }
351  } // namespace distributed
352 } // namespace parallel
353 
354 
355 // explicit instantiations
356 # include "error_predictor.inst"
357 
359 
360 #endif /* DEAL_II_WITH_P4EST */
parallel::distributed::ErrorPredictor::register_data_attach
void register_data_attach()
Definition: error_predictor.cc:89
error_predictor.h
tria_accessor.h
VectorOperation::insert
@ insert
Definition: vector_operation.h:49
LinearAlgebra::distributed::Vector< Number >
parallel::distributed::ErrorPredictor::dof_handler
SmartPointer< const hp::DoFHandler< dim, spacedim >, ErrorPredictor< dim, spacedim > > dof_handler
Definition: error_predictor.h:197
tria_iterator.h
parallel::distributed::ErrorPredictor::pack_callback
std::vector< char > pack_callback(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status)
Definition: error_predictor.cc:166
parallel::distributed::Triangulation::notify_ready_to_unpack
void notify_ready_to_unpack(const unsigned int handle, const std::function< void(const cell_iterator &, const CellStatus, const boost::iterator_range< std::vector< char >::const_iterator > &)> &unpack_callback)
Definition: tria.cc:4378
Utilities::pack
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
Definition: utilities.h:1209
hp::DoFHandler::get_fe
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index) const
parallel::distributed::ErrorPredictor::unpack
void unpack(std::vector< Vector< float > * > &all_out)
Definition: error_predictor.cc:124
hp::DoFHandler
Definition: dof_handler.h:203
la_parallel_vector.h
Triangulation< dim, dim >::CellStatus
CellStatus
Definition: tria.h:1969
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
block_vector.h
parallel::distributed::ErrorPredictor::ErrorPredictor
ErrorPredictor(const hp::DoFHandler< dim, spacedim > &dof)
Definition: error_predictor.cc:53
la_parallel_block_vector.h
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
Utilities::unpack
T unpack(const std::vector< char > &buffer, const bool allow_compression=true)
Definition: utilities.h:1353
parallel::distributed::Triangulation
Definition: tria.h:242
petsc_block_vector.h
numbers
Definition: numbers.h:207
hp::DoFHandler::cell_iterator
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:325
parallel::Triangulation
TriangulationBase< dim, spacedim > Triangulation
Definition: tria_base.h:302
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
tria.h
dof_tools.h
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
vector.h
dof_accessor.h
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
StandardExceptions::ExcDimensionMismatch
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
petsc_vector.h
config.h
parallel::distributed::ErrorPredictor::unpack_callback
void unpack_callback(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status, const boost::iterator_range< std::vector< char >::const_iterator > &data_range, std::vector< Vector< float > * > &all_out)
Definition: error_predictor.cc:302
parallel::distributed::ErrorPredictor::prepare_for_coarsening_and_refinement
void prepare_for_coarsening_and_refinement(const std::vector< const Vector< float > * > &all_in, const double gamma_p=std::sqrt(0.4), const double gamma_h=2., const double gamma_n=1.)
Definition: error_predictor.cc:72
parallel::distributed::Triangulation::register_data_attach
unsigned int register_data_attach(const std::function< std::vector< char >(const cell_iterator &, const CellStatus)> &pack_callback, const bool returns_variable_size_data)
Definition: tria.cc:4348
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
trilinos_vector.h
TriaIterator
Definition: tria_iterator.h:578
parallel
Definition: distributed.h:416
trilinos_parallel_block_vector.h
dof_handler.h