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\}}\)
data_out.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 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 
17 
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>
28 
29 #include <deal.II/hp/dof_handler.h>
30 #include <deal.II/hp/fe_values.h>
31 
33 
34 #include <sstream>
35 
37 
38 
39 namespace internal
40 {
41  namespace DataOutImplementation
42  {
43  template <int dim, int spacedim>
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>
70 {
71  // For the moment, just call the existing virtual functions. This
72  // preserves backward compatibility. When these deprecated functions are
73  // removed, we'll just inline their functionality into the lambda
74  // functions created here.
75  set_cell_selection(
76  [this](const Triangulation<dim, spacedim> &) {
77  return this->first_locally_owned_cell();
78  },
79  [this](const Triangulation<dim, spacedim> &, const cell_iterator &cell) {
80  return this->next_locally_owned_cell(cell);
81  });
82 }
83 
84 
85 
86 template <int dim, typename DoFHandlerType>
87 void
89  const std::pair<cell_iterator, unsigned int> *cell_and_index,
90  internal::DataOutImplementation::ParallelData<DoFHandlerType::dimension,
91  DoFHandlerType::space_dimension>
92  & scratch_data,
93  const unsigned int n_subdivisions,
94  const CurvedCellRegion curved_cell_region)
95 {
96  // first create the output object that we will write into
97  ::DataOutBase::Patch<DoFHandlerType::dimension,
98  DoFHandlerType::space_dimension>
99  patch;
100  patch.n_subdivisions = n_subdivisions;
101 
102  // set the vertices of the patch. if the mapping does not preserve locations
103  // (e.g. MappingQEulerian), we need to compute the offset of the vertex for
104  // the graphical output. Otherwise, we can just use the vertex info.
105  for (const unsigned int vertex :
107  if (scratch_data.mapping_collection[0].preserves_vertex_locations())
108  patch.vertices[vertex] = cell_and_index->first->vertex(vertex);
109  else
110  patch.vertices[vertex] =
111  scratch_data.mapping_collection[0].transform_unit_to_real_cell(
112  cell_and_index->first,
114 
115  // initialize FEValues
116  scratch_data.reinit_all_fe_values(this->dof_data, cell_and_index->first);
117 
119  &fe_patch_values = scratch_data.get_present_fe_values(0);
120 
121  const unsigned int n_q_points = fe_patch_values.n_quadrature_points;
122 
123  // First fill the geometric information for the patch: Where are the
124  // nodes in question located.
125  //
126  // Depending on the requested output of curved cells, if necessary
127  // append the quadrature points to the last rows of the patch.data
128  // member. This is the case if we want to produce curved cells at the
129  // boundary and this cell actually is at the boundary, or else if we
130  // want to produce curved cells everywhere
131  //
132  // note: a cell is *always* at the boundary if dim<spacedim
133  if (curved_cell_region == curved_inner_cells ||
134  (curved_cell_region == curved_boundary &&
135  (cell_and_index->first->at_boundary() ||
136  (DoFHandlerType::dimension != DoFHandlerType::space_dimension))))
137  {
138  Assert(patch.space_dim == DoFHandlerType::space_dimension,
139  ExcInternalError());
140 
141  // set the flag indicating that for this cell the points are
142  // explicitly given
143  patch.points_are_available = true;
144 
145  // then resize the patch.data member in order to have enough memory for
146  // the quadrature points as well, and copy the quadrature points there
147  const std::vector<Point<DoFHandlerType::space_dimension>> &q_points =
148  fe_patch_values.get_quadrature_points();
149 
150  patch.data.reinit(scratch_data.n_datasets +
151  DoFHandlerType::space_dimension,
152  n_q_points);
153  for (unsigned int i = 0; i < DoFHandlerType::space_dimension; ++i)
154  for (unsigned int q = 0; q < n_q_points; ++q)
155  patch.data(patch.data.size(0) - DoFHandlerType::space_dimension + i,
156  q) = q_points[q][i];
157  }
158  else
159  {
160  patch.data.reinit(scratch_data.n_datasets, n_q_points);
161  patch.points_are_available = false;
162  }
163 
164 
165  // Next fill the information we get from DoF data
166  if (scratch_data.n_datasets > 0)
167  {
168  // counter for data records
169  unsigned int offset = 0;
170 
171  // first fill dof_data
172  unsigned int dataset_number = 0;
173  for (const auto &dataset : this->dof_data)
174  {
175  const FEValuesBase<DoFHandlerType::dimension,
176  DoFHandlerType::space_dimension>
177  &this_fe_patch_values =
178  scratch_data.get_present_fe_values(dataset_number);
179  const unsigned int n_components =
180  this_fe_patch_values.get_fe().n_components();
181 
183  *postprocessor = dataset->postprocessor;
184 
185  if (postprocessor != nullptr)
186  {
187  // We have to postprocess the data, so determine, which fields
188  // have to be updated
189  const UpdateFlags update_flags =
190  postprocessor->get_needed_update_flags();
191 
192  if ((n_components == 1) &&
193  (dataset->is_complex_valued() == false))
194  {
195  // At each point there is only one component of value,
196  // gradient etc. Based on the 'if' statement above, we
197  // know that the solution is scalar and real-valued, so
198  // we do not need to worry about getting any imaginary
199  // components to the postprocessor, and we can safely
200  // call the function that evaluates a scalar field
201  if (update_flags & update_values)
202  dataset->get_function_values(
203  this_fe_patch_values,
205  real_part,
206  scratch_data.patch_values_scalar.solution_values);
207 
208  if (update_flags & update_gradients)
209  dataset->get_function_gradients(
210  this_fe_patch_values,
212  real_part,
213  scratch_data.patch_values_scalar.solution_gradients);
214 
215  if (update_flags & update_hessians)
216  dataset->get_function_hessians(
217  this_fe_patch_values,
219  real_part,
220  scratch_data.patch_values_scalar.solution_hessians);
221 
222  // Also fill some of the other fields postprocessors may
223  // want to access.
224  if (update_flags & update_quadrature_points)
225  scratch_data.patch_values_scalar.evaluation_points =
226  this_fe_patch_values.get_quadrature_points();
227 
228  const typename DoFHandlerType::cell_iterator dh_cell(
229  &cell_and_index->first->get_triangulation(),
230  cell_and_index->first->level(),
231  cell_and_index->first->index(),
232  dataset->dof_handler);
233  scratch_data.patch_values_scalar
234  .template set_cell<DoFHandlerType>(dh_cell);
235 
236  // Finally call the postprocessor's function that
237  // deals with scalar inputs.
238  postprocessor->evaluate_scalar_field(
239  scratch_data.patch_values_scalar,
240  scratch_data.postprocessed_values[dataset_number]);
241  }
242  else
243  {
244  // At each point we now have to evaluate a vector valued
245  // function and its derivatives. It may be that the solution
246  // is scalar and complex-valued, but we treat this as a vector
247  // field with two components.
248 
249  // At this point, we need to ask how we fill the fields that
250  // we want to pass on to the postprocessor. If the field in
251  // question is real-valued, we'll just extract the (only)
252  // real component from the solution fields
253  if (dataset->is_complex_valued() == false)
254  {
255  scratch_data.resize_system_vectors(n_components);
256 
257  if (update_flags & update_values)
258  dataset->get_function_values(
259  this_fe_patch_values,
261  real_part,
262  scratch_data.patch_values_system.solution_values);
263 
264  if (update_flags & update_gradients)
265  dataset->get_function_gradients(
266  this_fe_patch_values,
268  real_part,
269  scratch_data.patch_values_system.solution_gradients);
270 
271  if (update_flags & update_hessians)
272  dataset->get_function_hessians(
273  this_fe_patch_values,
275  real_part,
276  scratch_data.patch_values_system.solution_hessians);
277  }
278  else
279  {
280  // The solution is complex-valued. Let's cover the scalar
281  // case first (i.e., one scalar but complex-valued field,
282  // which we will have to split into its real and imaginar
283  // parts).
284  if (n_components == 1)
285  {
286  scratch_data.resize_system_vectors(2);
287 
288  // First get the real component of the scalar solution
289  // and copy the data into the
290  // scratch_data.patch_values_system output fields
291  if (update_flags & update_values)
292  {
293  dataset->get_function_values(
294  this_fe_patch_values,
296  ComponentExtractor::real_part,
297  scratch_data.patch_values_scalar
298  .solution_values);
299 
300  for (unsigned int i = 0;
301  i < scratch_data.patch_values_scalar
302  .solution_values.size();
303  ++i)
304  {
306  scratch_data.patch_values_system
307  .solution_values[i]
308  .size(),
309  2);
310  scratch_data.patch_values_system
311  .solution_values[i][0] =
312  scratch_data.patch_values_scalar
313  .solution_values[i];
314  }
315  }
316 
317  if (update_flags & update_gradients)
318  {
319  dataset->get_function_gradients(
320  this_fe_patch_values,
322  ComponentExtractor::real_part,
323  scratch_data.patch_values_scalar
324  .solution_gradients);
325 
326  for (unsigned int i = 0;
327  i < scratch_data.patch_values_scalar
328  .solution_gradients.size();
329  ++i)
330  {
332  scratch_data.patch_values_system
333  .solution_values[i]
334  .size(),
335  2);
336  scratch_data.patch_values_system
337  .solution_gradients[i][0] =
338  scratch_data.patch_values_scalar
339  .solution_gradients[i];
340  }
341  }
342 
343  if (update_flags & update_hessians)
344  {
345  dataset->get_function_hessians(
346  this_fe_patch_values,
348  ComponentExtractor::real_part,
349  scratch_data.patch_values_scalar
350  .solution_hessians);
351 
352  for (unsigned int i = 0;
353  i < scratch_data.patch_values_scalar
354  .solution_hessians.size();
355  ++i)
356  {
358  scratch_data.patch_values_system
359  .solution_hessians[i]
360  .size(),
361  2);
362  scratch_data.patch_values_system
363  .solution_hessians[i][0] =
364  scratch_data.patch_values_scalar
365  .solution_hessians[i];
366  }
367  }
368 
369  // Now we also have to get the imaginary
370  // component of the scalar solution
371  // and copy the data into the
372  // scratch_data.patch_values_system output fields
373  // that follow the real one
374  if (update_flags & update_values)
375  {
376  dataset->get_function_values(
377  this_fe_patch_values,
379  ComponentExtractor::imaginary_part,
380  scratch_data.patch_values_scalar
381  .solution_values);
382 
383  for (unsigned int i = 0;
384  i < scratch_data.patch_values_scalar
385  .solution_values.size();
386  ++i)
387  {
389  scratch_data.patch_values_system
390  .solution_values[i]
391  .size(),
392  2);
393  scratch_data.patch_values_system
394  .solution_values[i][1] =
395  scratch_data.patch_values_scalar
396  .solution_values[i];
397  }
398  }
399 
400  if (update_flags & update_gradients)
401  {
402  dataset->get_function_gradients(
403  this_fe_patch_values,
405  ComponentExtractor::imaginary_part,
406  scratch_data.patch_values_scalar
407  .solution_gradients);
408 
409  for (unsigned int i = 0;
410  i < scratch_data.patch_values_scalar
411  .solution_gradients.size();
412  ++i)
413  {
415  scratch_data.patch_values_system
416  .solution_values[i]
417  .size(),
418  2);
419  scratch_data.patch_values_system
420  .solution_gradients[i][1] =
421  scratch_data.patch_values_scalar
422  .solution_gradients[i];
423  }
424  }
425 
426  if (update_flags & update_hessians)
427  {
428  dataset->get_function_hessians(
429  this_fe_patch_values,
431  ComponentExtractor::imaginary_part,
432  scratch_data.patch_values_scalar
433  .solution_hessians);
434 
435  for (unsigned int i = 0;
436  i < scratch_data.patch_values_scalar
437  .solution_hessians.size();
438  ++i)
439  {
441  scratch_data.patch_values_system
442  .solution_hessians[i]
443  .size(),
444  2);
445  scratch_data.patch_values_system
446  .solution_hessians[i][1] =
447  scratch_data.patch_values_scalar
448  .solution_hessians[i];
449  }
450  }
451  }
452  else
453  {
454  scratch_data.resize_system_vectors(2 * n_components);
455 
456  // This is the vector-valued, complex-valued case. In
457  // essence, we just need to do the same as above,
458  // i.e., call the functions in dataset
459  // to retrieve first the real and then the imaginary
460  // part of the solution, then copy them to the
461  // scratch_data.patch_values_system. The difference to
462  // the scalar case is that there, we could (ab)use the
463  // scratch_data.patch_values_scalar members for first
464  // the real part and then the imaginary part, copying
465  // them into the scratch_data.patch_values_system
466  // variable one after the other. We can't do this
467  // here because the solution is vector-valued, and
468  // so using the single
469  // scratch_data.patch_values_system object doesn't
470  // work.
471  //
472  // Rather, we need to come up with a temporary object
473  // for this (one for each the values, the gradients,
474  // and the hessians).
475  //
476  // Compared to the previous code path, we here
477  // first get the real and imaginary parts of the
478  // values, then the real and imaginary parts of the
479  // gradients, etc. This allows us to scope the
480  // temporary objects better
481  if (update_flags & update_values)
482  {
483  std::vector<Vector<double>> tmp(
484  scratch_data.patch_values_system.solution_values
485  .size(),
487 
488  // First get the real part into the tmp object
489  dataset->get_function_values(
490  this_fe_patch_values,
492  ComponentExtractor::real_part,
493  tmp);
494 
495  // Then copy these values into the first
496  // n_components slots of the output object.
497  for (unsigned int i = 0;
498  i < scratch_data.patch_values_system
499  .solution_values.size();
500  ++i)
501  {
503  scratch_data.patch_values_system
504  .solution_values[i]
505  .size(),
506  2 * n_components);
507  for (unsigned int j = 0; j < n_components;
508  ++j)
509  scratch_data.patch_values_system
510  .solution_values[i][j] = tmp[i][j];
511  }
512 
513  // Then do the same with the imaginary part,
514  // copying past the end of the previous set of
515  // values.
516  dataset->get_function_values(
517  this_fe_patch_values,
519  ComponentExtractor::imaginary_part,
520  tmp);
521 
522  for (unsigned int i = 0;
523  i < scratch_data.patch_values_system
524  .solution_values.size();
525  ++i)
526  {
527  for (unsigned int j = 0; j < n_components;
528  ++j)
529  scratch_data.patch_values_system
530  .solution_values[i][j + n_components] =
531  tmp[i][j];
532  }
533  }
534 
535  // Now do the exact same thing for the gradients
536  if (update_flags & update_gradients)
537  {
538  std::vector<std::vector<
540  tmp(
541  scratch_data.patch_values_system
542  .solution_gradients.size(),
543  std::vector<
545  n_components));
546 
547  // First the real part
548  dataset->get_function_gradients(
549  this_fe_patch_values,
551  ComponentExtractor::real_part,
552  tmp);
553 
554  for (unsigned int i = 0;
555  i < scratch_data.patch_values_system
556  .solution_gradients.size();
557  ++i)
558  {
560  scratch_data.patch_values_system
561  .solution_gradients[i]
562  .size(),
563  2 * n_components);
564  for (unsigned int j = 0; j < n_components;
565  ++j)
566  scratch_data.patch_values_system
567  .solution_gradients[i][j] = tmp[i][j];
568  }
569 
570  // Then the imaginary part
571  dataset->get_function_gradients(
572  this_fe_patch_values,
574  ComponentExtractor::imaginary_part,
575  tmp);
576 
577  for (unsigned int i = 0;
578  i < scratch_data.patch_values_system
579  .solution_gradients.size();
580  ++i)
581  {
582  for (unsigned int j = 0; j < n_components;
583  ++j)
584  scratch_data.patch_values_system
585  .solution_gradients[i][j + n_components] =
586  tmp[i][j];
587  }
588  }
589 
590  // And finally the Hessians. Same scheme as above.
591  if (update_flags & update_hessians)
592  {
593  std::vector<std::vector<
595  tmp(
596  scratch_data.patch_values_system
597  .solution_gradients.size(),
598  std::vector<
600  n_components));
601 
602  // First the real part
603  dataset->get_function_hessians(
604  this_fe_patch_values,
606  ComponentExtractor::real_part,
607  tmp);
608 
609  for (unsigned int i = 0;
610  i < scratch_data.patch_values_system
611  .solution_hessians.size();
612  ++i)
613  {
615  scratch_data.patch_values_system
616  .solution_hessians[i]
617  .size(),
618  2 * n_components);
619  for (unsigned int j = 0; j < n_components;
620  ++j)
621  scratch_data.patch_values_system
622  .solution_hessians[i][j] = tmp[i][j];
623  }
624 
625  // Then the imaginary part
626  dataset->get_function_hessians(
627  this_fe_patch_values,
629  ComponentExtractor::imaginary_part,
630  tmp);
631 
632  for (unsigned int i = 0;
633  i < scratch_data.patch_values_system
634  .solution_hessians.size();
635  ++i)
636  {
637  for (unsigned int j = 0; j < n_components;
638  ++j)
639  scratch_data.patch_values_system
640  .solution_hessians[i][j + n_components] =
641  tmp[i][j];
642  }
643  }
644  }
645  }
646 
647  // Now set other fields we may need
648  if (update_flags & update_quadrature_points)
649  scratch_data.patch_values_system.evaluation_points =
650  this_fe_patch_values.get_quadrature_points();
651 
652  const typename DoFHandlerType::cell_iterator dh_cell(
653  &cell_and_index->first->get_triangulation(),
654  cell_and_index->first->level(),
655  cell_and_index->first->index(),
656  dataset->dof_handler);
657  scratch_data.patch_values_system
658  .template set_cell<DoFHandlerType>(dh_cell);
659 
660  // Whether the solution was complex-scalar or
661  // complex-vector-valued doesn't matter -- we took it apart
662  // into several fields and so we have to call the
663  // evaluate_vector_field() function.
664  postprocessor->evaluate_vector_field(
665  scratch_data.patch_values_system,
666  scratch_data.postprocessed_values[dataset_number]);
667  }
668 
669  // Now we need to copy the result of the postprocessor to
670  // the Patch object where it can then be further processed
671  // by the functions in DataOutBase
672  for (unsigned int q = 0; q < n_q_points; ++q)
673  for (unsigned int component = 0;
674  component < dataset->n_output_variables;
675  ++component)
676  patch.data(offset + component, q) =
677  scratch_data.postprocessed_values[dataset_number][q](
678  component);
679 
680  // Move the counter for the output location forward as
681  // appropriate
682  offset += dataset->n_output_variables;
683  }
684  else
685  {
686  // use the given data vector directly, without a postprocessor.
687  // again, we treat single component functions separately for
688  // efficiency reasons.
689  if (n_components == 1)
690  {
691  Assert(dataset->n_output_variables == 1, ExcInternalError());
692 
693  // First output the real part of the solution vector
694  dataset->get_function_values(
695  this_fe_patch_values,
697  real_part,
698  scratch_data.patch_values_scalar.solution_values);
699  for (unsigned int q = 0; q < n_q_points; ++q)
700  patch.data(offset, q) =
701  scratch_data.patch_values_scalar.solution_values[q];
702  offset += 1;
703 
704  // And if there is one, also output the imaginary part. Note
705  // that the problem is scalar-valued, so we can freely add the
706  // imaginary part after the real part without having to worry
707  // that we are interleaving the real components of a vector
708  // with the imaginary components of the same vector.
709  if (dataset->is_complex_valued() == true)
710  {
711  dataset->get_function_values(
712  this_fe_patch_values,
714  imaginary_part,
715  scratch_data.patch_values_scalar.solution_values);
716  for (unsigned int q = 0; q < n_q_points; ++q)
717  patch.data(offset, q) =
718  scratch_data.patch_values_scalar.solution_values[q];
719  offset += 1;
720  }
721  }
722  else
723  {
724  scratch_data.resize_system_vectors(n_components);
725 
726  // So we have a multi-component DoFHandler here. That's more
727  // complicated. If the vector is real-valued, then we can just
728  // get everything at all quadrature points and copy them into
729  // the output array. In fact, we don't have to worry at all
730  // about the interpretation of the components.
731  if (dataset->is_complex_valued() == false)
732  {
733  dataset->get_function_values(
734  this_fe_patch_values,
736  real_part,
737  scratch_data.patch_values_system.solution_values);
738  for (unsigned int component = 0; component < n_components;
739  ++component)
740  for (unsigned int q = 0; q < n_q_points; ++q)
741  patch.data(offset + component, q) =
742  scratch_data.patch_values_system.solution_values[q](
743  component);
744 
745  // Increment the counter for the actual data record.
746  offset += dataset->n_output_variables;
747  }
748  else
749  // The situation is more complicated if the input vector is
750  // complex-valued. The easiest approach would have been to
751  // just have all real and then all imaginary components.
752  // This would have been conceptually easy, but it has the
753  // annoying downside that if you have a vector-valued
754  // problem (say, [u v]) then the output order would have
755  // been [u_re, v_re, u_im, v_im]. That's tolerable, but not
756  // quite so nice because one typically thinks of real and
757  // imaginary parts as belonging together. We would really
758  // like the output order to be [u_re, u_im, v_re, v_im].
759  // That, too, would have been easy to implement because one
760  // just has to interleave real and imaginary parts.
761  //
762  // But that's also not what we want. That's because if one
763  // were, for example, to solve a complex-valued Stokes
764  // problem (e.g., computing eigenfunctions of the Stokes
765  // operator), then one has solution components
766  // [[u v] p] and the proper output order is
767  // [[u_re v_re] [u_im v_im] p_re p_im].
768  // In other words, the order in which we want to output
769  // data depends on the *interpretation* of components.
770  //
771  // Doing this requires a bit more code, and also needs to
772  // be in sync with what we do in
773  // DataOut_DoFData::get_dataset_names() and
774  // DataOut_DoFData::get_nonscalar_data_ranges().
775  {
776  // Given this description, first get the real parts of
777  // all components:
778  dataset->get_function_values(
779  this_fe_patch_values,
781  real_part,
782  scratch_data.patch_values_system.solution_values);
783 
784  // Then we need to distribute them to the correct
785  // location. This requires knowledge of the interpretation
786  // of components as discussed above.
787  {
788  Assert(dataset->data_component_interpretation.size() ==
789  n_components,
790  ExcInternalError());
791 
792  unsigned int destination = offset;
793  for (unsigned int component = 0;
794  component < n_components;
795  /* component is updated below */)
796  {
797  switch (
798  dataset->data_component_interpretation[component])
799  {
802  {
803  // OK, a scalar component. Put all of the
804  // values into the current row
805  // ('destination'); then move 'component'
806  // forward by one (so we treat the next
807  // component) and 'destination' forward by
808  // two (because we're going to put the
809  // imaginary part of the current component
810  // into the next slot).
811  for (unsigned int q = 0; q < n_q_points;
812  ++q)
813  patch.data(destination, q) =
814  scratch_data.patch_values_system
815  .solution_values[q](component);
816 
817  ++component;
818  destination += 2;
819 
820  break;
821  }
822 
825  {
826  // A vector component. Put the
827  // DoFHandlerType::space_dimension
828  // components into the next set of
829  // contiguous rows
830  // ('destination+c'); then move 'component'
831  // forward by spacedim (so we get to the
832  // next component after the current vector)
833  // and 'destination' forward by two*spacedim
834  // (because we're going to put the imaginary
835  // part of the vector into the subsequent
836  // spacedim slots).
837  const unsigned int size =
838  DoFHandlerType::space_dimension;
839  for (unsigned int c = 0; c < size; ++c)
840  for (unsigned int q = 0; q < n_q_points;
841  ++q)
842  patch.data(destination + c, q) =
843  scratch_data.patch_values_system
844  .solution_values[q](component + c);
845 
846  component += size;
847  destination += 2 * size;
848 
849  break;
850  }
851 
854  {
855  // Same approach as for vectors above.
856  const unsigned int size =
857  DoFHandlerType::space_dimension *
858  DoFHandlerType::space_dimension;
859  for (unsigned int c = 0; c < size; ++c)
860  for (unsigned int q = 0; q < n_q_points;
861  ++q)
862  patch.data(destination + c, q) =
863  scratch_data.patch_values_system
864  .solution_values[q](component + c);
865 
866  component += size;
867  destination += 2 * size;
868 
869  break;
870  }
871 
872  default:
873  Assert(false, ExcNotImplemented());
874  }
875  }
876  }
877 
878  // And now we need to do the same thing again for the
879  // imaginary parts, starting at the top of the list of
880  // components/destinations again.
881  dataset->get_function_values(
882  this_fe_patch_values,
884  imaginary_part,
885  scratch_data.patch_values_system.solution_values);
886  {
887  unsigned int destination = offset;
888  for (unsigned int component = 0;
889  component < n_components;
890  /* component is updated below */)
891  {
892  switch (
893  dataset->data_component_interpretation[component])
894  {
897  {
898  // OK, a scalar component. Put all of the
899  // values into the row past the current one
900  // ('destination+1') since 'destination' is
901  // occupied by the real part.
902  for (unsigned int q = 0; q < n_q_points;
903  ++q)
904  patch.data(destination + 1, q) =
905  scratch_data.patch_values_system
906  .solution_values[q](component);
907 
908  ++component;
909  destination += 2;
910 
911  break;
912  }
913 
916  {
917  // A vector component. Put the
918  // DoFHandlerType::space_dimension
919  // components into the set of contiguous
920  // rows that follow the real parts
921  // ('destination+spacedim+c').
922  const unsigned int size =
923  DoFHandlerType::space_dimension;
924  for (unsigned int c = 0; c < size; ++c)
925  for (unsigned int q = 0; q < n_q_points;
926  ++q)
927  patch.data(destination + size + c, q) =
928  scratch_data.patch_values_system
929  .solution_values[q](component + c);
930 
931  component += size;
932  destination += 2 * size;
933 
934  break;
935  }
936 
939  {
940  // Same as for vectors.
941  const unsigned int size =
942  DoFHandlerType::space_dimension *
943  DoFHandlerType::space_dimension;
944  for (unsigned int c = 0; c < size; ++c)
945  for (unsigned int q = 0; q < n_q_points;
946  ++q)
947  patch.data(destination + size + c, q) =
948  scratch_data.patch_values_system
949  .solution_values[q](component + c);
950 
951  component += size;
952  destination += 2 * size;
953 
954  break;
955  }
956 
957  default:
958  Assert(false, ExcNotImplemented());
959  }
960  }
961  }
962 
963  // Increment the counter for the actual data record. We
964  // need to move it forward a number of positions equal to
965  // the number of components of this data set, times two
966  // because we dealt with a complex-valued input vector
967  offset += dataset->n_output_variables * 2;
968  }
969  }
970  }
971 
972  // Also update the dataset_number index that we carry along with the
973  // for-loop over all data sets.
974  ++dataset_number;
975  }
976 
977  // Then do the cell data. At least, we don't have to worry about
978  // complex-valued vectors/tensors since cell data is always scalar.
979  if (this->cell_data.size() != 0)
980  {
981  Assert(!cell_and_index->first->has_children(), ExcNotImplemented());
982 
983  for (const auto &dataset : this->cell_data)
984  {
985  // as above, first output the real part
986  {
987  const double value =
988  dataset->get_cell_data_value(cell_and_index->second,
990  ComponentExtractor::real_part);
991  for (unsigned int q = 0; q < n_q_points; ++q)
992  patch.data(offset, q) = value;
993  }
994 
995  // and if there is one, also output the imaginary part
996  if (dataset->is_complex_valued() == true)
997  {
998  const double value = dataset->get_cell_data_value(
999  cell_and_index->second,
1001  imaginary_part);
1002  for (unsigned int q = 0; q < n_q_points; ++q)
1003  patch.data(offset + 1, q) = value;
1004  }
1005 
1006  offset += (dataset->is_complex_valued() ? 2 : 1);
1007  }
1008  }
1009  }
1010 
1011 
1012  for (const unsigned int f :
1014  {
1015  // let's look up whether the neighbor behind that face is noted in the
1016  // table of cells which we treat. this can only happen if the neighbor
1017  // exists, and is on the same level as this cell, but it may also happen
1018  // that the neighbor is not a member of the range of cells over which we
1019  // loop, in which case the respective entry in the
1020  // cell_to_patch_index_map will have the value no_neighbor. (note that
1021  // since we allocated only as much space in this array as the maximum
1022  // index of the cells we loop over, not every neighbor may have its
1023  // space in it, so we have to assume that it is extended by values
1024  // no_neighbor)
1025  if (cell_and_index->first->at_boundary(f) ||
1026  (cell_and_index->first->neighbor(f)->level() !=
1027  cell_and_index->first->level()))
1028  {
1029  patch.neighbors[f] = numbers::invalid_unsigned_int;
1030  continue;
1031  }
1032 
1033  const cell_iterator neighbor = cell_and_index->first->neighbor(f);
1034  Assert(static_cast<unsigned int>(neighbor->level()) <
1035  scratch_data.cell_to_patch_index_map->size(),
1036  ExcInternalError());
1037  if ((static_cast<unsigned int>(neighbor->index()) >=
1038  (*scratch_data.cell_to_patch_index_map)[neighbor->level()].size()) ||
1039  ((*scratch_data.cell_to_patch_index_map)[neighbor->level()]
1040  [neighbor->index()] ==
1042  {
1043  patch.neighbors[f] = numbers::invalid_unsigned_int;
1044  continue;
1045  }
1046 
1047  // now, there is a neighbor, so get its patch number and set it for the
1048  // neighbor index
1049  patch.neighbors[f] =
1050  (*scratch_data
1051  .cell_to_patch_index_map)[neighbor->level()][neighbor->index()];
1052  }
1053 
1054  const unsigned int patch_idx =
1055  (*scratch_data.cell_to_patch_index_map)[cell_and_index->first->level()]
1056  [cell_and_index->first->index()];
1057  // did we mess up the indices?
1058  Assert(patch_idx < this->patches.size(), ExcInternalError());
1059  patch.patch_index = patch_idx;
1060 
1061  // Put the patch into the patches vector. instead of copying the data,
1062  // simply swap the contents to avoid the penalty of writing into another
1063  // processor's memory
1064  this->patches[patch_idx].swap(patch);
1065 }
1066 
1067 
1068 
1069 template <int dim, typename DoFHandlerType>
1070 void
1071 DataOut<dim, DoFHandlerType>::build_patches(const unsigned int n_subdivisions)
1072 {
1073  build_patches(StaticMappingQ1<DoFHandlerType::dimension,
1074  DoFHandlerType::space_dimension>::mapping,
1075  n_subdivisions,
1076  no_curved_cells);
1077 }
1078 
1079 
1080 
1081 template <int dim, typename DoFHandlerType>
1082 void
1085  & mapping,
1086  const unsigned int n_subdivisions_,
1087  const CurvedCellRegion curved_region)
1088 {
1089  // Check consistency of redundant template parameter
1090  Assert(dim == DoFHandlerType::dimension,
1091  ExcDimensionMismatch(dim, DoFHandlerType::dimension));
1092 
1093  Assert(this->triangulation != nullptr,
1095 
1096  const unsigned int n_subdivisions =
1097  (n_subdivisions_ != 0) ? n_subdivisions_ : this->default_subdivisions;
1098  Assert(n_subdivisions >= 1,
1100  n_subdivisions));
1101 
1102  this->validate_dataset_names();
1103 
1104  // First count the cells we want to create patches of. Also fill the object
1105  // that maps the cell indices to the patch numbers, as this will be needed
1106  // for generation of neighborship information.
1107  // Note, there is a confusing mess of different indices here at play:
1108  // - patch_index: the index of a patch in all_cells
1109  // - cell->index: only unique on each level, used in cell_to_patch_index_map
1110  // - active_index: index for a cell when counting from begin_active() using
1111  // ++cell (identical to cell->active_cell_index())
1112  // - cell_index: unique index of a cell counted using
1113  // next_cell_function() starting from first_cell_function()
1114  //
1115  // It turns out that we create one patch for each selected cell, so
1116  // patch_index==cell_index.
1117  //
1118  // Now construct the map such that
1119  // cell_to_patch_index_map[cell->level][cell->index] = patch_index
1120  std::vector<std::vector<unsigned int>> cell_to_patch_index_map;
1121  cell_to_patch_index_map.resize(this->triangulation->n_levels());
1122  for (unsigned int l = 0; l < this->triangulation->n_levels(); ++l)
1123  {
1124  // max_index is the largest cell->index on level l
1125  unsigned int max_index = 0;
1126  for (cell_iterator cell = first_cell_function(*this->triangulation);
1127  cell != this->triangulation->end();
1128  cell = next_cell_function(*this->triangulation, cell))
1129  if (static_cast<unsigned int>(cell->level()) == l)
1130  max_index =
1131  std::max(max_index, static_cast<unsigned int>(cell->index()));
1132 
1133  cell_to_patch_index_map[l].resize(
1134  max_index + 1,
1136  DoFHandlerType::dimension,
1137  DoFHandlerType::space_dimension>::no_neighbor);
1138  }
1139 
1140  // will be all_cells[patch_index] = pair(cell, active_index)
1141  std::vector<std::pair<cell_iterator, unsigned int>> all_cells;
1142  {
1143  // important: we need to compute the active_index of the cell in the range
1144  // 0..n_active_cells() because this is where we need to look up cell
1145  // data from (cell data vectors do not have the length distance computed by
1146  // first_cell_function/next_cell_function because this might skip
1147  // some values (FilteredIterator).
1148  auto active_cell = this->triangulation->begin_active();
1149  unsigned int active_index = 0;
1150  cell_iterator cell = first_cell_function(*this->triangulation);
1151  for (; cell != this->triangulation->end();
1152  cell = next_cell_function(*this->triangulation, cell))
1153  {
1154  // move forward until active_cell points at the cell (cell) we are
1155  // looking at to compute the current active_index
1156  while (active_cell != this->triangulation->end() && cell->is_active() &&
1157  decltype(active_cell)(cell) != active_cell)
1158  {
1159  ++active_cell;
1160  ++active_index;
1161  }
1162 
1163  Assert(static_cast<unsigned int>(cell->level()) <
1164  cell_to_patch_index_map.size(),
1165  ExcInternalError());
1166  Assert(static_cast<unsigned int>(cell->index()) <
1167  cell_to_patch_index_map[cell->level()].size(),
1168  ExcInternalError());
1169  Assert(active_index < this->triangulation->n_active_cells(),
1170  ExcInternalError());
1171  cell_to_patch_index_map[cell->level()][cell->index()] =
1172  all_cells.size();
1173 
1174  all_cells.emplace_back(cell, active_index);
1175  }
1176  }
1177 
1178  this->patches.clear();
1179  this->patches.resize(all_cells.size());
1180 
1181  // Now create a default object for the WorkStream object to work with. The
1182  // first step is to count how many output data sets there will be. This is,
1183  // in principle, just the number of components of each data set, but we
1184  // need to allocate two entries per component if there are
1185  // complex-valued input data (unless we use a postprocessor on this
1186  // output -- all postprocessor outputs are real-valued)
1187  unsigned int n_datasets = 0;
1188  for (unsigned int i = 0; i < this->cell_data.size(); ++i)
1189  n_datasets += (this->cell_data[i]->is_complex_valued() &&
1190  (this->cell_data[i]->postprocessor == nullptr) ?
1191  2 :
1192  1);
1193  for (unsigned int i = 0; i < this->dof_data.size(); ++i)
1194  n_datasets += (this->dof_data[i]->n_output_variables *
1195  (this->dof_data[i]->is_complex_valued() &&
1196  (this->dof_data[i]->postprocessor == nullptr) ?
1197  2 :
1198  1));
1199 
1200  std::vector<unsigned int> n_postprocessor_outputs(this->dof_data.size());
1201  for (unsigned int dataset = 0; dataset < this->dof_data.size(); ++dataset)
1202  if (this->dof_data[dataset]->postprocessor)
1203  n_postprocessor_outputs[dataset] =
1204  this->dof_data[dataset]->n_output_variables;
1205  else
1206  n_postprocessor_outputs[dataset] = 0;
1207 
1208  const CurvedCellRegion curved_cell_region =
1209  (n_subdivisions < 2 ? no_curved_cells : curved_region);
1210 
1211  UpdateFlags update_flags = update_values;
1212  if (curved_cell_region != no_curved_cells)
1213  update_flags |= update_quadrature_points;
1214 
1215  for (unsigned int i = 0; i < this->dof_data.size(); ++i)
1216  if (this->dof_data[i]->postprocessor)
1217  update_flags |=
1218  this->dof_data[i]->postprocessor->get_needed_update_flags();
1219  // perhaps update_normal_vectors is present, which would only be useful on
1220  // faces, but we may not use it here.
1221  Assert(
1222  !(update_flags & update_normal_vectors),
1223  ExcMessage(
1224  "The update of normal vectors may not be requested for evaluation of "
1225  "data on cells via DataPostprocessor."));
1226 
1227  internal::DataOutImplementation::ParallelData<DoFHandlerType::dimension,
1228  DoFHandlerType::space_dimension>
1229  thread_data(n_datasets,
1230  n_subdivisions,
1231  n_postprocessor_outputs,
1232  mapping,
1233  this->get_fes(),
1234  update_flags,
1235  cell_to_patch_index_map);
1236 
1237  auto worker = [this, n_subdivisions, curved_cell_region](
1238  const std::pair<cell_iterator, unsigned int> *cell_and_index,
1240  DoFHandlerType::dimension,
1241  DoFHandlerType::space_dimension> &scratch_data,
1242  // this function doesn't actually need a copy data object --
1243  // it just writes everything right into the output array
1244  int) {
1245  this->build_one_patch(cell_and_index,
1246  scratch_data,
1247  n_subdivisions,
1248  curved_cell_region);
1249  };
1250 
1251  // now build the patches in parallel
1252  if (all_cells.size() > 0)
1253  WorkStream::run(all_cells.data(),
1254  all_cells.data() + all_cells.size(),
1255  worker,
1256  // no copy-local-to-global function needed here
1257  std::function<void(const int)>(),
1258  thread_data,
1259  /* dummy CopyData object = */ 0,
1260  // experimenting shows that we can make things run a bit
1261  // faster if we increase the number of cells we work on
1262  // per item (i.e., WorkStream's chunk_size argument,
1263  // about 10% improvement) and the items in flight at any
1264  // given time (another 5% on the testcase discussed in
1265  // @ref workstream_paper, on 32 cores) and if
1267  64);
1268 }
1269 
1270 
1271 
1272 template <int dim, typename DoFHandlerType>
1273 void
1275  const std::function<cell_iterator(const Triangulation<dim, spacedim> &)>
1276  & first_cell,
1277  const std::function<cell_iterator(const Triangulation<dim, spacedim> &,
1278  const cell_iterator &)> &next_cell)
1279 {
1280  first_cell_function = first_cell;
1281  next_cell_function = next_cell;
1282 }
1283 
1284 
1285 
1286 template <int dim, typename DoFHandlerType>
1287 void
1289  const FilteredIterator<cell_iterator> &filtered_iterator)
1290 {
1291  const auto first_cell =
1292  [filtered_iterator](const Triangulation<dim, spacedim> &triangulation) {
1293  // Create a copy of the filtered iterator so that we can
1294  // call a non-const function -- though we are really only
1295  // interested in the return value of that function, not the
1296  // state of the object
1297  FilteredIterator<cell_iterator> x = filtered_iterator;
1298  return x.set_to_next_positive(triangulation.begin());
1299  };
1300 
1301 
1302  const auto next_cell =
1303  [filtered_iterator](const Triangulation<dim, spacedim> &,
1304  const cell_iterator &cell) {
1305  // Create a copy of the filtered iterator so that we can
1306  // call a non-const function -- though we are really only
1307  // interested in the return value of that function, not the
1308  // state of the object
1309  FilteredIterator<cell_iterator> x = filtered_iterator;
1310 
1311  // Set the iterator to 'cell'. Since 'cell' must satisfy the
1312  // predicate (that's how it was created), set_to_next_positive
1313  // simply sets the iterator to 'cell'.
1314  x.set_to_next_positive(cell);
1315 
1316  // Advance by one:
1317  ++x;
1318 
1319  return x;
1320  };
1321 
1322  set_cell_selection(first_cell, next_cell);
1323 }
1324 
1325 
1326 
1327 template <int dim, typename DoFHandlerType>
1328 const std::pair<typename DataOut<dim, DoFHandlerType>::FirstCellFunctionType,
1331 {
1332  return std::make_pair(first_cell_function, next_cell_function);
1333 }
1334 
1335 
1336 
1337 template <int dim, typename DoFHandlerType>
1340 {
1341  return this->triangulation->begin_active();
1342 }
1343 
1344 
1345 
1346 template <int dim, typename DoFHandlerType>
1349  const typename DataOut<dim, DoFHandlerType>::cell_iterator &cell)
1350 {
1351  // convert the iterator to an active_iterator and advance this to the next
1352  // active cell
1353  typename Triangulation<DoFHandlerType::dimension,
1354  DoFHandlerType::space_dimension>::active_cell_iterator
1355  active_cell = cell;
1356  ++active_cell;
1357  return active_cell;
1358 }
1359 
1360 
1361 
1362 template <int dim, typename DoFHandlerType>
1365 {
1366  typename DataOut<dim, DoFHandlerType>::cell_iterator cell = first_cell();
1367 
1368  // skip cells if the current one has no children (is active) and is a ghost
1369  // or artificial cell
1370  while ((cell != this->triangulation->end()) && cell->is_active() &&
1371  !cell->is_locally_owned())
1372  cell = next_cell(cell);
1373 
1374  return cell;
1375 }
1376 
1377 
1378 
1379 template <int dim, typename DoFHandlerType>
1382  const typename DataOut<dim, DoFHandlerType>::cell_iterator &old_cell)
1383 {
1385  next_cell(old_cell);
1386  while ((cell != this->triangulation->end()) && cell->is_active() &&
1387  !cell->is_locally_owned())
1388  cell = next_cell(cell);
1389  return cell;
1390 }
1391 
1392 
1393 // explicit instantiations
1394 #include "data_out.inst"
1395 
DataPostprocessor< DoFHandlerType::space_dimension >
fe_values.h
update_quadrature_points
@ update_quadrature_points
Transformed quadrature points.
Definition: fe_update_flags.h:122
StaticMappingQ1
Definition: mapping_q1.h:88
MultithreadInfo::n_threads
static unsigned int n_threads()
Definition: multithread_info.cc:169
DataComponentInterpretation::component_is_scalar
@ component_is_scalar
Definition: data_component_interpretation.h:55
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
internal::DataOutImplementation::ParallelDataBase
Definition: data_out_dof_data.h:405
Triangulation
Definition: tria.h:1109
tria.h
FilteredIterator::set_to_next_positive
FilteredIterator & set_to_next_positive(const BaseIterator &bi)
Definition: filtered_iterator.h:994
tria_iterator.h
mapping_q1.h
DataOut::next_locally_owned_cell
virtual cell_iterator next_locally_owned_cell(const cell_iterator &cell)
Definition: data_out.cc:1381
GeometryInfo
Definition: geometry_info.h:1224
DataOut::set_cell_selection
void set_cell_selection(const std::function< cell_iterator(const Triangulation< dim, spacedim > &)> &first_cell, const std::function< cell_iterator(const Triangulation< dim, spacedim > &, const cell_iterator &)> &next_cell)
Definition: data_out.cc:1274
DoFHandler::n_components
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
update_values
@ update_values
Shape function values.
Definition: fe_update_flags.h:78
DataOut::build_patches
virtual void build_patches(const unsigned int n_subdivisions=0)
Definition: data_out.cc:1071
WorkStream::run
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:1185
work_stream.h
DataOut::first_cell
virtual cell_iterator first_cell()
Definition: data_out.cc:1339
internal::DataOutImplementation::ComponentExtractor
ComponentExtractor
Definition: data_out_dof_data.h:197
DataPostprocessor::evaluate_vector_field
virtual void evaluate_vector_field(const DataPostprocessorInputs::Vector< dim > &input_data, std::vector< Vector< double >> &computed_quantities) const
Definition: data_postprocessor.cc:37
internal::DataOutImplementation
Definition: data_out.h:33
Mapping
Abstract base class for mapping classes.
Definition: mapping.h:302
Exceptions::DataOutImplementation::ExcInvalidNumberOfSubdivisions
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
DataOut::DataOut
DataOut()
Definition: data_out.cc:69
FilteredIterator
Definition: filtered_iterator.h:529
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
Tensor
Definition: tensor.h:450
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
DataOut::NextCellFunctionType
typename std::function< cell_iterator(const Triangulation< dim, spacedim > &, const cell_iterator &)> NextCellFunctionType
Definition: data_out.h:181
update_gradients
@ update_gradients
Shape function gradients.
Definition: fe_update_flags.h:84
DataOut< dim >::CurvedCellRegion
CurvedCellRegion
Definition: data_out.h:196
Physics::Elasticity::Kinematics::l
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
FEValuesBase::get_quadrature_points
const std::vector< Point< spacedim > > & get_quadrature_points() const
fe.h
DataOut< dim >::cell_iterator
typename DataOut_DoFData< DoFHandler< dim >, DoFHandler< dim > ::dimension, DoFHandler< dim > ::space_dimension >::cell_iterator cell_iterator
Definition: data_out.h:166
fe_dgq.h
internal::DataOutImplementation::ParallelData::ParallelData
ParallelData(const unsigned int n_datasets, const unsigned int n_subdivisions, const std::vector< unsigned int > &n_postprocessor_outputs, const Mapping< dim, spacedim > &mapping, const std::vector< std::shared_ptr<::hp::FECollection< dim, spacedim >>> &finite_elements, const UpdateFlags update_flags, const std::vector< std::vector< unsigned int >> &cell_to_patch_index_map)
Definition: data_out.cc:44
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
FEValuesBase
Definition: fe.h:36
UpdateFlags
UpdateFlags
Definition: fe_update_flags.h:66
value
static const bool value
Definition: dof_tools_constraints.cc:433
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
DataOut::get_cell_selection
const std::pair< FirstCellFunctionType, NextCellFunctionType > get_cell_selection() const
Definition: data_out.cc:1330
update_normal_vectors
@ update_normal_vectors
Normal vectors.
Definition: fe_update_flags.h:136
DataOut::next_cell
virtual cell_iterator next_cell(const cell_iterator &cell)
Definition: data_out.cc:1348
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
update_hessians
@ update_hessians
Second derivatives of shape functions.
Definition: fe_update_flags.h:90
DataOutBase::Patch::n_subdivisions
unsigned int n_subdivisions
Definition: data_out_base.h:287
dof_handler.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)
Exceptions::DataOutImplementation::ExcNoTriangulationSelected
static ::ExceptionBase & ExcNoTriangulationSelected()
FEValuesBase::n_quadrature_points
const unsigned int n_quadrature_points
Definition: fe_values.h:2101
DataPostprocessor::get_needed_update_flags
virtual UpdateFlags get_needed_update_flags() const =0
triangulation
const typename ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
Definition: p4est_wrappers.cc:69
internal
Definition: aligned_vector.h:369
DataOutBase::Patch
Definition: data_out_base.h:248
fe_values.h
DataComponentInterpretation::component_is_part_of_tensor
@ component_is_part_of_tensor
Definition: data_component_interpretation.h:67
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
data_out.h
Vector< double >
DataPostprocessor::evaluate_scalar_field
virtual void evaluate_scalar_field(const DataPostprocessorInputs::Scalar< dim > &input_data, std::vector< Vector< double >> &computed_quantities) const
Definition: data_postprocessor.cc:26
internal::DataOutImplementation::ParallelData
Definition: data_out.h:41
DataOut::first_locally_owned_cell
virtual cell_iterator first_locally_owned_cell()
Definition: data_out.cc:1364
Utilities::MPI::max
T max(const T &t, const MPI_Comm &mpi_communicator)
hp::FECollection
Definition: fe_collection.h:54
DataComponentInterpretation::component_is_part_of_vector
@ component_is_part_of_vector
Definition: data_component_interpretation.h:61
int
DataOut::build_one_patch
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:88
dof_handler.h