Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
data_out.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2019 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 #include <deal.II/base/work_stream.h>
17 
18 #include <deal.II/dofs/dof_accessor.h>
19 #include <deal.II/dofs/dof_handler.h>
20 
21 #include <deal.II/fe/fe.h>
22 #include <deal.II/fe/fe_dgq.h>
23 #include <deal.II/fe/fe_values.h>
24 #include <deal.II/fe/mapping_q1.h>
25 
26 #include <deal.II/grid/tria.h>
27 #include <deal.II/grid/tria_iterator.h>
28 
29 #include <deal.II/hp/dof_handler.h>
30 #include <deal.II/hp/fe_values.h>
31 
32 #include <deal.II/numerics/data_out.h>
33 
34 #include <sstream>
35 
36 DEAL_II_NAMESPACE_OPEN
37 
38 
39 namespace internal
40 {
41  namespace DataOutImplementation
42  {
43  template <int dim, int spacedim>
44  ParallelData<dim, spacedim>::ParallelData(
45  const unsigned int n_datasets,
46  const unsigned int n_subdivisions,
47  const std::vector<unsigned int> &n_postprocessor_outputs,
48  const Mapping<dim, spacedim> & mapping,
49  const std::vector<
50  std::shared_ptr<::hp::FECollection<dim, spacedim>>>
51  & finite_elements,
52  const UpdateFlags update_flags,
53  const std::vector<std::vector<unsigned int>> &cell_to_patch_index_map)
54  : ParallelDataBase<dim, spacedim>(n_datasets,
55  n_subdivisions,
56  n_postprocessor_outputs,
57  mapping,
58  finite_elements,
59  update_flags,
60  false)
61  , cell_to_patch_index_map(&cell_to_patch_index_map)
62  {}
63  } // namespace DataOutImplementation
64 } // namespace internal
65 
66 
67 
68 template <int dim, typename DoFHandlerType>
69 void
71  const std::pair<cell_iterator, unsigned int> *cell_and_index,
72  internal::DataOutImplementation::ParallelData<DoFHandlerType::dimension,
73  DoFHandlerType::space_dimension>
74  & scratch_data,
75  const unsigned int n_subdivisions,
76  const CurvedCellRegion curved_cell_region)
77 {
78  // first create the output object that we will write into
79  ::DataOutBase::Patch<DoFHandlerType::dimension,
80  DoFHandlerType::space_dimension>
81  patch;
82  patch.n_subdivisions = n_subdivisions;
83 
84  // set the vertices of the patch. if the mapping does not preserve locations
85  // (e.g. MappingQEulerian), we need to compute the offset of the vertex for
86  // the graphical output. Otherwise, we can just use the vertex info.
87  for (unsigned int vertex = 0;
88  vertex < GeometryInfo<DoFHandlerType::dimension>::vertices_per_cell;
89  ++vertex)
90  if (scratch_data.mapping_collection[0].preserves_vertex_locations())
91  patch.vertices[vertex] = cell_and_index->first->vertex(vertex);
92  else
93  patch.vertices[vertex] =
94  scratch_data.mapping_collection[0].transform_unit_to_real_cell(
95  cell_and_index->first,
97 
98  // create DoFHandlerType::active_cell_iterator and initialize FEValues
99  scratch_data.reinit_all_fe_values(this->dof_data, cell_and_index->first);
100 
102  &fe_patch_values = scratch_data.get_present_fe_values(0);
103 
104  const unsigned int n_q_points = fe_patch_values.n_quadrature_points;
105 
106  // depending on the requested output of curved cells, if necessary
107  // append the quadrature points to the last rows of the patch.data
108  // member. This is the case if we want to produce curved cells at the
109  // boundary and this cell actually is at the boundary, or else if we
110  // want to produce curved cells everywhere
111  //
112  // note: a cell is *always* at the boundary if dim<spacedim
113  if (curved_cell_region == curved_inner_cells ||
114  (curved_cell_region == curved_boundary &&
115  (cell_and_index->first->at_boundary() ||
116  (DoFHandlerType::dimension != DoFHandlerType::space_dimension))))
117  {
118  Assert(patch.space_dim == DoFHandlerType::space_dimension,
119  ExcInternalError());
120 
121  // set the flag indicating that for this cell the points are
122  // explicitly given
123  patch.points_are_available = true;
124 
125  // then resize the patch.data member in order to have enough memory for
126  // the quadrature points as well, and copy the quadrature points there
127  const std::vector<Point<DoFHandlerType::space_dimension>> &q_points =
128  fe_patch_values.get_quadrature_points();
129 
130  patch.data.reinit(scratch_data.n_datasets +
131  DoFHandlerType::space_dimension,
132  n_q_points);
133  for (unsigned int i = 0; i < DoFHandlerType::space_dimension; ++i)
134  for (unsigned int q = 0; q < n_q_points; ++q)
135  patch.data(patch.data.size(0) - DoFHandlerType::space_dimension + i,
136  q) = q_points[q][i];
137  }
138  else
139  {
140  patch.data.reinit(scratch_data.n_datasets, n_q_points);
141  patch.points_are_available = false;
142  }
143 
144 
145  if (scratch_data.n_datasets > 0)
146  {
147  // counter for data records
148  unsigned int offset = 0;
149 
150  // first fill dof_data
151  for (unsigned int dataset = 0; dataset < this->dof_data.size(); ++dataset)
152  {
153  const FEValuesBase<DoFHandlerType::dimension,
154  DoFHandlerType::space_dimension>
155  &this_fe_patch_values = scratch_data.get_present_fe_values(dataset);
156  const unsigned int n_components =
157  this_fe_patch_values.get_fe().n_components();
158 
160  *postprocessor = this->dof_data[dataset]->postprocessor;
161 
162  if (postprocessor != nullptr)
163  {
164  // we have to postprocess the data, so determine, which fields
165  // have to be updated
166  const UpdateFlags update_flags =
167  postprocessor->get_needed_update_flags();
168  if (n_components == 1)
169  {
170  // at each point there is only one component of value,
171  // gradient etc.
172  if (update_flags & update_values)
173  this->dof_data[dataset]->get_function_values(
174  this_fe_patch_values,
175  internal::DataOutImplementation::ComponentExtractor::
176  real_part,
177  scratch_data.patch_values_scalar.solution_values);
178  if (update_flags & update_gradients)
179  this->dof_data[dataset]->get_function_gradients(
180  this_fe_patch_values,
181  internal::DataOutImplementation::ComponentExtractor::
182  real_part,
183  scratch_data.patch_values_scalar.solution_gradients);
184  if (update_flags & update_hessians)
185  this->dof_data[dataset]->get_function_hessians(
186  this_fe_patch_values,
187  internal::DataOutImplementation::ComponentExtractor::
188  real_part,
189  scratch_data.patch_values_scalar.solution_hessians);
190 
191  if (update_flags & update_quadrature_points)
192  scratch_data.patch_values_scalar.evaluation_points =
193  this_fe_patch_values.get_quadrature_points();
194 
195  const typename DoFHandlerType::active_cell_iterator dh_cell(
196  &cell_and_index->first->get_triangulation(),
197  cell_and_index->first->level(),
198  cell_and_index->first->index(),
199  this->dof_data[dataset]->dof_handler);
200  scratch_data.patch_values_scalar
201  .template set_cell<DoFHandlerType>(dh_cell);
202 
203  postprocessor->evaluate_scalar_field(
204  scratch_data.patch_values_scalar,
205  scratch_data.postprocessed_values[dataset]);
206  }
207  else
208  {
209  scratch_data.resize_system_vectors(n_components);
210 
211  // at each point there is a vector valued function and its
212  // derivative...
213  if (update_flags & update_values)
214  this->dof_data[dataset]->get_function_values(
215  this_fe_patch_values,
216  internal::DataOutImplementation::ComponentExtractor::
217  real_part,
218  scratch_data.patch_values_system.solution_values);
219  if (update_flags & update_gradients)
220  this->dof_data[dataset]->get_function_gradients(
221  this_fe_patch_values,
222  internal::DataOutImplementation::ComponentExtractor::
223  real_part,
224  scratch_data.patch_values_system.solution_gradients);
225  if (update_flags & update_hessians)
226  this->dof_data[dataset]->get_function_hessians(
227  this_fe_patch_values,
228  internal::DataOutImplementation::ComponentExtractor::
229  real_part,
230  scratch_data.patch_values_system.solution_hessians);
231 
232  if (update_flags & update_quadrature_points)
233  scratch_data.patch_values_system.evaluation_points =
234  this_fe_patch_values.get_quadrature_points();
235 
236  const typename DoFHandlerType::active_cell_iterator dh_cell(
237  &cell_and_index->first->get_triangulation(),
238  cell_and_index->first->level(),
239  cell_and_index->first->index(),
240  this->dof_data[dataset]->dof_handler);
241  scratch_data.patch_values_system
242  .template set_cell<DoFHandlerType>(dh_cell);
243 
244  postprocessor->evaluate_vector_field(
245  scratch_data.patch_values_system,
246  scratch_data.postprocessed_values[dataset]);
247  }
248 
249  for (unsigned int q = 0; q < n_q_points; ++q)
250  for (unsigned int component = 0;
251  component < this->dof_data[dataset]->n_output_variables;
252  ++component)
253  patch.data(offset + component, q) =
254  scratch_data.postprocessed_values[dataset][q](component);
255  }
256  else
257  // use the given data vector directly, without a postprocessor.
258  // again, we treat single component functions separately for
259  // efficiency reasons.
260  if (n_components == 1)
261  {
262  // first output the real part of the solution vector
263  this->dof_data[dataset]->get_function_values(
264  this_fe_patch_values,
265  internal::DataOutImplementation::ComponentExtractor::real_part,
266  scratch_data.patch_values_scalar.solution_values);
267  for (unsigned int q = 0; q < n_q_points; ++q)
268  patch.data(offset, q) =
269  scratch_data.patch_values_scalar.solution_values[q];
270 
271  // and if there is one, also output the imaginary part
272  if (this->dof_data[dataset]->is_complex_valued() == true)
273  {
274  this->dof_data[dataset]->get_function_values(
275  this_fe_patch_values,
276  internal::DataOutImplementation::ComponentExtractor::
277  imaginary_part,
278  scratch_data.patch_values_scalar.solution_values);
279  for (unsigned int q = 0; q < n_q_points; ++q)
280  patch.data(offset + 1, q) =
281  scratch_data.patch_values_scalar.solution_values[q];
282  }
283  }
284  else
285  {
286  scratch_data.resize_system_vectors(n_components);
287 
288  // same as above: first the real part
289  const unsigned int stride =
290  (this->dof_data[dataset]->is_complex_valued() ? 2 : 1);
291  this->dof_data[dataset]->get_function_values(
292  this_fe_patch_values,
293  internal::DataOutImplementation::ComponentExtractor::real_part,
294  scratch_data.patch_values_system.solution_values);
295  for (unsigned int component = 0; component < n_components;
296  ++component)
297  for (unsigned int q = 0; q < n_q_points; ++q)
298  patch.data(offset + component * stride, q) =
299  scratch_data.patch_values_system.solution_values[q](
300  component);
301 
302  // and if there is one, also output the imaginary part
303  if (this->dof_data[dataset]->is_complex_valued() == true)
304  {
305  this->dof_data[dataset]->get_function_values(
306  this_fe_patch_values,
307  internal::DataOutImplementation::ComponentExtractor::
308  imaginary_part,
309  scratch_data.patch_values_system.solution_values);
310  for (unsigned int component = 0; component < n_components;
311  ++component)
312  for (unsigned int q = 0; q < n_q_points; ++q)
313  patch.data(offset + component * stride + 1, q) =
314  scratch_data.patch_values_system.solution_values[q](
315  component);
316  }
317  }
318 
319  // increment the counter for the actual data record
320  offset += this->dof_data[dataset]->n_output_variables *
321  (this->dof_data[dataset]->is_complex_valued() ? 2 : 1);
322  }
323 
324  // then do the cell data. only compute the number of a cell if needed;
325  // also make sure that we only access cell data if the
326  // first_cell/next_cell functions only return active cells
327  if (this->cell_data.size() != 0)
328  {
329  Assert(!cell_and_index->first->has_children(), ExcNotImplemented());
330 
331  for (unsigned int dataset = 0; dataset < this->cell_data.size();
332  ++dataset)
333  {
334  // as above, first output the real part
335  {
336  const double value =
337  this->cell_data[dataset]->get_cell_data_value(
338  cell_and_index->second,
339  internal::DataOutImplementation::ComponentExtractor::
340  real_part);
341  for (unsigned int q = 0; q < n_q_points; ++q)
342  patch.data(offset, q) = value;
343  }
344 
345  // and if there is one, also output the imaginary part
346  if (this->cell_data[dataset]->is_complex_valued() == true)
347  {
348  const double value =
349  this->cell_data[dataset]->get_cell_data_value(
350  cell_and_index->second,
351  internal::DataOutImplementation::ComponentExtractor::
352  imaginary_part);
353  for (unsigned int q = 0; q < n_q_points; ++q)
354  patch.data(offset + 1, q) = value;
355  }
356 
357  offset += (this->cell_data[dataset]->is_complex_valued() ? 2 : 1);
358  }
359  }
360  }
361 
362 
363  for (unsigned int f = 0;
364  f < GeometryInfo<DoFHandlerType::dimension>::faces_per_cell;
365  ++f)
366  {
367  // let's look up whether the neighbor behind that face is noted in the
368  // table of cells which we treat. this can only happen if the neighbor
369  // exists, and is on the same level as this cell, but it may also happen
370  // that the neighbor is not a member of the range of cells over which we
371  // loop, in which case the respective entry in the
372  // cell_to_patch_index_map will have the value no_neighbor. (note that
373  // since we allocated only as much space in this array as the maximum
374  // index of the cells we loop over, not every neighbor may have its
375  // space in it, so we have to assume that it is extended by values
376  // no_neighbor)
377  if (cell_and_index->first->at_boundary(f) ||
378  (cell_and_index->first->neighbor(f)->level() !=
379  cell_and_index->first->level()))
380  {
381  patch.neighbors[f] = numbers::invalid_unsigned_int;
382  continue;
383  }
384 
385  const cell_iterator neighbor = cell_and_index->first->neighbor(f);
386  Assert(static_cast<unsigned int>(neighbor->level()) <
387  scratch_data.cell_to_patch_index_map->size(),
388  ExcInternalError());
389  if ((static_cast<unsigned int>(neighbor->index()) >=
390  (*scratch_data.cell_to_patch_index_map)[neighbor->level()].size()) ||
391  ((*scratch_data.cell_to_patch_index_map)[neighbor->level()]
392  [neighbor->index()] ==
394  {
395  patch.neighbors[f] = numbers::invalid_unsigned_int;
396  continue;
397  }
398 
399  // now, there is a neighbor, so get its patch number and set it for the
400  // neighbor index
401  patch.neighbors[f] =
402  (*scratch_data
403  .cell_to_patch_index_map)[neighbor->level()][neighbor->index()];
404  }
405 
406  const unsigned int patch_idx =
407  (*scratch_data.cell_to_patch_index_map)[cell_and_index->first->level()]
408  [cell_and_index->first->index()];
409  // did we mess up the indices?
410  Assert(patch_idx < this->patches.size(), ExcInternalError());
411  patch.patch_index = patch_idx;
412 
413  // Put the patch into the patches vector. instead of copying the data,
414  // simply swap the contents to avoid the penalty of writing into another
415  // processor's memory
416  this->patches[patch_idx].swap(patch);
417 }
418 
419 
420 
421 template <int dim, typename DoFHandlerType>
422 void
423 DataOut<dim, DoFHandlerType>::build_patches(const unsigned int n_subdivisions)
424 {
425  build_patches(StaticMappingQ1<DoFHandlerType::dimension,
426  DoFHandlerType::space_dimension>::mapping,
427  n_subdivisions,
428  no_curved_cells);
429 }
430 
431 
432 
433 template <int dim, typename DoFHandlerType>
434 void
437  & mapping,
438  const unsigned int n_subdivisions_,
439  const CurvedCellRegion curved_region)
440 {
441  // Check consistency of redundant template parameter
442  Assert(dim == DoFHandlerType::dimension,
443  ExcDimensionMismatch(dim, DoFHandlerType::dimension));
444 
445  Assert(this->triangulation != nullptr,
447 
448  const unsigned int n_subdivisions =
449  (n_subdivisions_ != 0) ? n_subdivisions_ : this->default_subdivisions;
450  Assert(n_subdivisions >= 1,
452  n_subdivisions));
453 
454  this->validate_dataset_names();
455 
456  // First count the cells we want to create patches of. Also fill the object
457  // that maps the cell indices to the patch numbers, as this will be needed
458  // for generation of neighborship information.
459  // Note, there is a confusing mess of different indices here at play:
460  // - patch_index: the index of a patch in all_cells
461  // - cell->index: only unique on each level, used in cell_to_patch_index_map
462  // - active_index: index for a cell when counting from begin_active() using
463  // ++cell (identical to cell->active_cell_index())
464  // - cell_index: unique index of a cell counted using
465  // next_locally_owned_cell() starting from first_locally_owned_cell()
466  //
467  // It turns out that we create one patch for each selected cell, so
468  // patch_index==cell_index.
469  //
470  // Now construct the map such that
471  // cell_to_patch_index_map[cell->level][cell->index] = patch_index
472  std::vector<std::vector<unsigned int>> cell_to_patch_index_map;
473  cell_to_patch_index_map.resize(this->triangulation->n_levels());
474  for (unsigned int l = 0; l < this->triangulation->n_levels(); ++l)
475  {
476  // max_index is the largest cell->index on level l
477  unsigned int max_index = 0;
478  for (cell_iterator cell = first_locally_owned_cell();
479  cell != this->triangulation->end();
480  cell = next_locally_owned_cell(cell))
481  if (static_cast<unsigned int>(cell->level()) == l)
482  max_index =
483  std::max(max_index, static_cast<unsigned int>(cell->index()));
484 
485  cell_to_patch_index_map[l].resize(
486  max_index + 1,
488  DoFHandlerType::dimension,
489  DoFHandlerType::space_dimension>::no_neighbor);
490  }
491 
492  // will be all_cells[patch_index] = pair(cell, active_index)
493  std::vector<std::pair<cell_iterator, unsigned int>> all_cells;
494  {
495  // important: we need to compute the active_index of the cell in the range
496  // 0..n_active_cells() because this is where we need to look up cell
497  // data from (cell data vectors do not have the length distance computed by
498  // first_locally_owned_cell/next_locally_owned_cell because this might skip
499  // some values (FilteredIterator).
500  active_cell_iterator active_cell = this->triangulation->begin_active();
501  unsigned int active_index = 0;
502  cell_iterator cell = first_locally_owned_cell();
503  for (; cell != this->triangulation->end();
504  cell = next_locally_owned_cell(cell))
505  {
506  // move forward until active_cell points at the cell (cell) we are
507  // looking at to compute the current active_index
508  while (active_cell != this->triangulation->end() && cell->active() &&
509  active_cell_iterator(cell) != active_cell)
510  {
511  ++active_cell;
512  ++active_index;
513  }
514 
515  Assert(static_cast<unsigned int>(cell->level()) <
516  cell_to_patch_index_map.size(),
517  ExcInternalError());
518  Assert(static_cast<unsigned int>(cell->index()) <
519  cell_to_patch_index_map[cell->level()].size(),
520  ExcInternalError());
521  Assert(active_index < this->triangulation->n_active_cells(),
522  ExcInternalError());
523  cell_to_patch_index_map[cell->level()][cell->index()] =
524  all_cells.size();
525 
526  all_cells.emplace_back(cell, active_index);
527  }
528  }
529 
530  this->patches.clear();
531  this->patches.resize(all_cells.size());
532 
533  // now create a default object for the WorkStream object to work with
534  unsigned int n_datasets = 0;
535  for (unsigned int i = 0; i < this->cell_data.size(); ++i)
536  n_datasets += (this->cell_data[i]->is_complex_valued() ? 2 : 1);
537  for (unsigned int i = 0; i < this->dof_data.size(); ++i)
538  n_datasets += (this->dof_data[i]->n_output_variables *
539  (this->dof_data[i]->is_complex_valued() ? 2 : 1));
540 
541  std::vector<unsigned int> n_postprocessor_outputs(this->dof_data.size());
542  for (unsigned int dataset = 0; dataset < this->dof_data.size(); ++dataset)
543  if (this->dof_data[dataset]->postprocessor)
544  n_postprocessor_outputs[dataset] =
545  this->dof_data[dataset]->n_output_variables;
546  else
547  n_postprocessor_outputs[dataset] = 0;
548 
549  const CurvedCellRegion curved_cell_region =
550  (n_subdivisions < 2 ? no_curved_cells : curved_region);
551 
552  UpdateFlags update_flags = update_values;
553  if (curved_cell_region != no_curved_cells)
554  update_flags |= update_quadrature_points;
555 
556  for (unsigned int i = 0; i < this->dof_data.size(); ++i)
557  if (this->dof_data[i]->postprocessor)
558  update_flags |=
559  this->dof_data[i]->postprocessor->get_needed_update_flags();
560  // perhaps update_normal_vectors is present, which would only be useful on
561  // faces, but we may not use it here.
562  Assert(
563  !(update_flags & update_normal_vectors),
564  ExcMessage(
565  "The update of normal vectors may not be requested for evaluation of "
566  "data on cells via DataPostprocessor."));
567 
568  internal::DataOutImplementation::ParallelData<DoFHandlerType::dimension,
569  DoFHandlerType::space_dimension>
570  thread_data(n_datasets,
571  n_subdivisions,
572  n_postprocessor_outputs,
573  mapping,
574  this->get_fes(),
575  update_flags,
576  cell_to_patch_index_map);
577 
578  // now build the patches in parallel
579  if (all_cells.size() > 0)
581  all_cells.data(),
582  all_cells.data() + all_cells.size(),
584  this,
585  std::placeholders::_1,
586  std::placeholders::_2,
587  /* no std::placeholders::_3, since this function doesn't
588  actually need a copy data object -- it just writes everything
589  right into the output array */
590  n_subdivisions,
591  curved_cell_region),
592  // no copy-local-to-global function needed here
593  std::function<void(const int)>(),
594  thread_data,
595  /* dummy CopyData object = */ 0,
596  // experimenting shows that we can make things run a bit
597  // faster if we increase the number of cells we work on
598  // per item (i.e., WorkStream's chunk_size argument,
599  // about 10% improvement) and the items in flight at any
600  // given time (another 5% on the testcase discussed in
601  // @ref workstream_paper, on 32 cores) and if
603  64);
604 }
605 
606 
607 
608 template <int dim, typename DoFHandlerType>
611 {
612  return this->triangulation->begin_active();
613 }
614 
615 
616 
617 template <int dim, typename DoFHandlerType>
620  const typename DataOut<dim, DoFHandlerType>::cell_iterator &cell)
621 {
622  // convert the iterator to an active_iterator and advance this to the next
623  // active cell
624  typename Triangulation<DoFHandlerType::dimension,
625  DoFHandlerType::space_dimension>::active_cell_iterator
626  active_cell = cell;
627  ++active_cell;
628  return active_cell;
629 }
630 
631 
632 
633 template <int dim, typename DoFHandlerType>
636 {
637  typename DataOut<dim, DoFHandlerType>::cell_iterator cell = first_cell();
638 
639  // skip cells if the current one has no children (is active) and is a ghost
640  // or artificial cell
641  while ((cell != this->triangulation->end()) &&
642  (cell->has_children() == false) && !cell->is_locally_owned())
643  cell = next_cell(cell);
644 
645  return cell;
646 }
647 
648 
649 
650 template <int dim, typename DoFHandlerType>
653  const typename DataOut<dim, DoFHandlerType>::cell_iterator &old_cell)
654 {
656  next_cell(old_cell);
657  while ((cell != this->triangulation->end()) &&
658  (cell->has_children() == false) && !cell->is_locally_owned())
659  cell = next_cell(cell);
660  return cell;
661 }
662 
663 
664 // explicit instantiations
665 #include "data_out.inst"
666 
667 DEAL_II_NAMESPACE_CLOSE
Shape function values.
static const unsigned int invalid_unsigned_int
Definition: types.h:173
static ::ExceptionBase & ExcNoTriangulationSelected()
void build_one_patch(const std::pair< cell_iterator, unsigned int > *cell_and_index, internal::DataOutImplementation::ParallelData< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &scratch_data, const unsigned int n_subdivisions, const CurvedCellRegion curved_cell_region)
Definition: data_out.cc:70
virtual cell_iterator next_cell(const cell_iterator &cell)
Definition: data_out.cc:619
typename DataOut_DoFData< DoFHandler< dim >, DoFHandler< dim > ::dimension, DoFHandler< dim > ::space_dimension >::cell_iterator cell_iterator
Definition: data_out.h:172
virtual void evaluate_vector_field(const DataPostprocessorInputs::Vector< dim > &input_data, std::vector< Vector< double >> &computed_quantities) const
Transformed quadrature points.
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
virtual void build_patches(const unsigned int n_subdivisions=0)
Definition: data_out.cc:423
const std::vector< Point< spacedim > > & get_quadrature_points() const
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual cell_iterator first_locally_owned_cell()
Definition: data_out.cc:635
#define Assert(cond, exc)
Definition: exceptions.h:1407
UpdateFlags
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
Definition: dof_tools.h:57
virtual cell_iterator first_cell()
Definition: data_out.cc:610
unsigned int n_subdivisions
Second derivatives of shape functions.
const unsigned int n_quadrature_points
Definition: fe_values.h:2099
Shape function gradients.
virtual void evaluate_scalar_field(const DataPostprocessorInputs::Scalar< dim > &input_data, std::vector< Vector< double >> &computed_quantities) const
Normal vectors.
static ::ExceptionBase & ExcNotImplemented()
void run(const std::vector< std::vector< Iterator >> &colored_iterators, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
Definition: work_stream.h:1167
static unsigned int n_threads()
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual cell_iterator next_locally_owned_cell(const cell_iterator &cell)
Definition: data_out.cc:652
static ::ExceptionBase & ExcInternalError()
virtual UpdateFlags get_needed_update_flags() const =0