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