Reference documentation for deal.II version 9.0.0
data_out_faces.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2018 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/quadrature_lib.h>
17 #include <deal.II/base/work_stream.h>
18 #include <deal.II/lac/vector.h>
19 #include <deal.II/lac/block_vector.h>
20 #include <deal.II/numerics/data_out_faces.h>
21 #include <deal.II/grid/tria.h>
22 #include <deal.II/dofs/dof_handler.h>
23 #include <deal.II/dofs/dof_accessor.h>
24 #include <deal.II/grid/tria_iterator.h>
25 #include <deal.II/fe/fe.h>
26 #include <deal.II/fe/fe_values.h>
27 #include <deal.II/hp/fe_values.h>
28 #include <deal.II/fe/mapping_q1.h>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 
33 namespace internal
34 {
35  namespace DataOutFacesImplementation
36  {
37  template <int dim, int spacedim>
38  ParallelData<dim,spacedim>::
39  ParallelData (const unsigned int n_datasets,
40  const unsigned int n_subdivisions,
41  const std::vector<unsigned int> &n_postprocessor_outputs,
42  const Mapping<dim,spacedim> &mapping,
43  const std::vector<std::shared_ptr<::hp::FECollection<dim,spacedim> > > &finite_elements,
44  const UpdateFlags update_flags)
45  :
46  internal::DataOutImplementation::
47  ParallelDataBase<dim,spacedim> (n_datasets,
48  n_subdivisions,
49  n_postprocessor_outputs,
50  mapping,
51  finite_elements,
52  update_flags,
53  true)
54  {}
55 
56 
57 
62  template <int dim, int spacedim>
63  void
64  append_patch_to_list (const DataOutBase::Patch<dim-1,spacedim> &patch,
65  std::vector<DataOutBase::Patch<dim-1,spacedim> > &patches)
66  {
67  patches.push_back (patch);
68  patches.back().patch_index = patches.size()-1;
69  }
70  }
71 }
72 
73 
74 
75 template <int dim, typename DoFHandlerType>
77  :
78  surface_only(so)
79 {
80  Assert (dim == DoFHandlerType::dimension,
82 }
83 
84 
85 
86 template <int dim, typename DoFHandlerType>
87 void
89 build_one_patch (const FaceDescriptor *cell_and_face,
92 {
93  Assert (cell_and_face->first->is_locally_owned(),
95 
96  // we use the mapping to transform the vertices. However, the mapping works
97  // on cells, not faces, so transform the face vertex to a cell vertex, that
98  // to a unit cell vertex and then, finally, that to the mapped vertex. In
99  // most cases this complicated procedure will be the identity.
100  for (unsigned int vertex=0; vertex<GeometryInfo<dimension-1>::vertices_per_cell; ++vertex)
101  patch.vertices[vertex] = data.mapping_collection[0].transform_unit_to_real_cell
102  (cell_and_face->first,
105  (cell_and_face->second,
106  vertex,
107  cell_and_face->first->face_orientation(cell_and_face->second),
108  cell_and_face->first->face_flip(cell_and_face->second),
109  cell_and_face->first->face_rotation(cell_and_face->second))));
110 
111  if (data.n_datasets > 0)
112  {
113  data.reinit_all_fe_values(this->dof_data, cell_and_face->first,
114  cell_and_face->second);
115  const FEValuesBase<dimension> &fe_patch_values
116  = data.get_present_fe_values (0);
117 
118  const unsigned int n_q_points = fe_patch_values.n_quadrature_points;
119 
120  // store the intermediate points
121  Assert(patch.space_dim==dimension, ExcInternalError());
122  const std::vector<Point<dimension> > &q_points=fe_patch_values.get_quadrature_points();
123  // resize the patch.data member in order to have enough memory for the
124  // quadrature points as well
125  patch.data.reinit(data.n_datasets+dimension,
126  patch.data.size(1));
127  // set the flag indicating that for this cell the points are explicitly
128  // given
129  patch.points_are_available=true;
130  // copy points to patch.data
131  for (unsigned int i=0; i<dimension; ++i)
132  for (unsigned int q=0; q<n_q_points; ++q)
133  patch.data(patch.data.size(0)-dimension+i,q)=q_points[q][i];
134 
135  // counter for data records
136  unsigned int offset=0;
137 
138  // first fill dof_data
139  for (unsigned int dataset=0; dataset<this->dof_data.size(); ++dataset)
140  {
141  const FEValuesBase<dimension> &this_fe_patch_values
142  = data.get_present_fe_values (dataset);
143  const unsigned int n_components
144  = this_fe_patch_values.get_fe().n_components();
145  const DataPostprocessor<dim> *postprocessor=this->dof_data[dataset]->postprocessor;
146  if (postprocessor != nullptr)
147  {
148  // we have to postprocess the data, so determine, which fields
149  // have to be updated
150  const UpdateFlags update_flags=postprocessor->get_needed_update_flags();
151 
152  if (n_components == 1)
153  {
154  // at each point there is only one component of value,
155  // gradient etc.
156  if (update_flags & update_values)
157  this->dof_data[dataset]->get_function_values (this_fe_patch_values,
158  internal::DataOutImplementation::ComponentExtractor::real_part,
159  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  internal::DataOutImplementation::ComponentExtractor::real_part,
163  data.patch_values_scalar.solution_gradients);
164  if (update_flags & update_hessians)
165  this->dof_data[dataset]->get_function_hessians (this_fe_patch_values,
166  internal::DataOutImplementation::ComponentExtractor::real_part,
167  data.patch_values_scalar.solution_hessians);
168 
169  if (update_flags & update_quadrature_points)
170  data.patch_values_scalar.evaluation_points = this_fe_patch_values.get_quadrature_points();
171 
172  if (update_flags & update_normal_vectors)
173  data.patch_values_scalar.normals = this_fe_patch_values.get_all_normal_vectors();
174 
175  const typename DoFHandlerType::active_cell_iterator dh_cell(&cell_and_face->first->get_triangulation(),
176  cell_and_face->first->level(),
177  cell_and_face->first->index(),
178  this->dof_data[dataset]->dof_handler);
179  data.patch_values_scalar.template set_cell<DoFHandlerType> (dh_cell);
180 
181  postprocessor->
182  evaluate_scalar_field(data.patch_values_scalar,
183  data.postprocessed_values[dataset]);
184  }
185  else
186  {
187  // at each point there is a vector valued function and its
188  // derivative...
189  data.resize_system_vectors(n_components);
190  if (update_flags & update_values)
191  this->dof_data[dataset]->get_function_values (this_fe_patch_values,
192  internal::DataOutImplementation::ComponentExtractor::real_part,
193  data.patch_values_system.solution_values);
194  if (update_flags & update_gradients)
195  this->dof_data[dataset]->get_function_gradients (this_fe_patch_values,
196  internal::DataOutImplementation::ComponentExtractor::real_part,
197  data.patch_values_system.solution_gradients);
198  if (update_flags & update_hessians)
199  this->dof_data[dataset]->get_function_hessians (this_fe_patch_values,
200  internal::DataOutImplementation::ComponentExtractor::real_part,
201  data.patch_values_system.solution_hessians);
202 
203  if (update_flags & update_quadrature_points)
204  data.patch_values_system.evaluation_points = this_fe_patch_values.get_quadrature_points();
205 
206  if (update_flags & update_normal_vectors)
207  data.patch_values_system.normals = this_fe_patch_values.get_all_normal_vectors();
208 
209  const typename DoFHandlerType::active_cell_iterator dh_cell(&cell_and_face->first->get_triangulation(),
210  cell_and_face->first->level(),
211  cell_and_face->first->index(),
212  this->dof_data[dataset]->dof_handler);
213  data.patch_values_system.template set_cell<DoFHandlerType> (dh_cell);
214 
215  postprocessor->
216  evaluate_vector_field(data.patch_values_system,
217  data.postprocessed_values[dataset]);
218  }
219 
220  for (unsigned int q=0; q<n_q_points; ++q)
221  for (unsigned int component=0;
222  component<this->dof_data[dataset]->n_output_variables; ++component)
223  patch.data(offset+component,q)
224  = data.postprocessed_values[dataset][q](component);
225  }
226  else
227  // now we use the given data vector without modifications. again,
228  // we treat single component functions separately for efficiency
229  // reasons.
230  if (n_components == 1)
231  {
232  this->dof_data[dataset]->get_function_values (this_fe_patch_values,
233  internal::DataOutImplementation::ComponentExtractor::real_part,
234  data.patch_values_scalar.solution_values);
235  for (unsigned int q=0; q<n_q_points; ++q)
236  patch.data(offset,q) = data.patch_values_scalar.solution_values[q];
237  }
238  else
239  {
240  data.resize_system_vectors(n_components);
241  this->dof_data[dataset]->get_function_values (this_fe_patch_values,
242  internal::DataOutImplementation::ComponentExtractor::real_part,
243  data.patch_values_system.solution_values);
244  for (unsigned int component=0; component<n_components;
245  ++component)
246  for (unsigned int q=0; q<n_q_points; ++q)
247  patch.data(offset+component,q) =
248  data.patch_values_system.solution_values[q](component);
249  }
250  // increment the counter for the actual data record
251  offset += this->dof_data[dataset]->n_output_variables;
252  }
253 
254  // then do the cell data
255  for (unsigned int dataset=0; dataset<this->cell_data.size(); ++dataset)
256  {
257  // we need to get at the number of the cell to which this face
258  // belongs in order to access the cell data. this is not readily
259  // available, so choose the following rather inefficient way:
260  Assert (cell_and_face->first->active(),
261  ExcMessage("The current function is trying to generate cell-data output "
262  "for a face that does not belong to an active cell. This is "
263  "not supported."));
264  const unsigned int cell_number
265  = std::distance (this->triangulation->begin_active(),
267 
268  const double value
269  = this->cell_data[dataset]->get_cell_data_value (cell_number,
270  internal::DataOutImplementation::ComponentExtractor::real_part);
271  for (unsigned int q=0; q<n_q_points; ++q)
272  patch.data(dataset+offset,q) = value;
273  }
274  }
275 }
276 
277 
278 
279 
280 template <int dim, typename DoFHandlerType>
281 void DataOutFaces<dim,DoFHandlerType>::build_patches (const unsigned int n_subdivisions_)
282 {
283  build_patches (StaticMappingQ1<dimension>::mapping, n_subdivisions_);
284 }
285 
286 
287 
288 template <int dim, typename DoFHandlerType>
290  const unsigned int n_subdivisions_)
291 {
292  // Check consistency of redundant template parameter
293  Assert (dim==dimension, ExcDimensionMismatch(dim, dimension));
294 
295  const unsigned int n_subdivisions = (n_subdivisions_ != 0)
296  ? n_subdivisions_
297  : this->default_subdivisions;
298 
299  Assert (n_subdivisions >= 1,
301 
302  Assert (this->triangulation != nullptr,
304 
305  this->validate_dataset_names();
306 
307  unsigned int n_datasets = this->cell_data.size();
308  for (unsigned int i=0; i<this->dof_data.size(); ++i)
309  n_datasets += this->dof_data[i]->n_output_variables;
310 
311  // first collect the cells we want to create patches of; we will
312  // then iterate over them. the end-condition of the loop needs to
313  // test that next_face() returns an end iterator, as well as for the
314  // case that first_face() returns an invalid FaceDescriptor object
315  std::vector<FaceDescriptor> all_faces;
316  for (FaceDescriptor face=first_face();
317  ((face.first != this->triangulation->end())
318  &&
319  (face != FaceDescriptor()));
320  face = next_face(face))
321  all_faces.push_back (face);
322 
323  // clear the patches array and allocate the right number of elements
324  this->patches.clear ();
325  this->patches.reserve (all_faces.size());
326  Assert (this->patches.size() == 0, ExcInternalError());
327 
328 
329  std::vector<unsigned int> n_postprocessor_outputs (this->dof_data.size());
330  for (unsigned int dataset=0; dataset<this->dof_data.size(); ++dataset)
331  if (this->dof_data[dataset]->postprocessor)
332  n_postprocessor_outputs[dataset] = this->dof_data[dataset]->n_output_variables;
333  else
334  n_postprocessor_outputs[dataset] = 0;
335 
336  UpdateFlags update_flags=update_values;
337  for (unsigned int i=0; i<this->dof_data.size(); ++i)
338  if (this->dof_data[i]->postprocessor)
339  update_flags |= this->dof_data[i]->postprocessor->get_needed_update_flags();
340  update_flags |= update_quadrature_points;
341 
343  thread_data (n_datasets,
344  n_subdivisions,
345  n_postprocessor_outputs,
346  mapping,
347  this->get_fes(),
348  update_flags);
349  DataOutBase::Patch<dimension-1,space_dimension> sample_patch;
350  sample_patch.n_subdivisions = n_subdivisions;
351  sample_patch.data.reinit (n_datasets,
352  Utilities::fixed_power<dimension-1>(n_subdivisions+1));
353 
354  // now build the patches in parallel
355  WorkStream::run (&all_faces[0],
356  &all_faces[0]+all_faces.size(),
358  this, std::placeholders::_1, std::placeholders::_2, std::placeholders::_3),
360  append_patch_to_list<dim,space_dimension>,
361  std::placeholders::_1, std::ref(this->patches)),
362  thread_data,
363  sample_patch);
364 }
365 
366 
367 
368 template <int dim, typename DoFHandlerType>
371 {
372  // simply find first active cell with a face on the boundary
373  typename Triangulation<dimension,space_dimension>::active_cell_iterator cell = this->triangulation->begin_active();
374  for (; cell != this->triangulation->end(); ++cell)
375  if (cell->is_locally_owned())
376  for (unsigned int f=0; f<GeometryInfo<dimension>::faces_per_cell; ++f)
377  if (!surface_only || cell->face(f)->at_boundary())
378  return FaceDescriptor(cell, f);
379 
380  // just return an invalid descriptor if we haven't found a locally
381  // owned face. this can happen in parallel where all boundary
382  // faces are owned by other processors
383  return FaceDescriptor();
384 }
385 
386 
387 
388 template <int dim, typename DoFHandlerType>
391 {
392  FaceDescriptor face = old_face;
393 
394  // first check whether the present cell has more faces on the boundary. since
395  // we started with this face, its cell must clearly be locally owned
396  Assert (face.first->is_locally_owned(), ExcInternalError());
397  for (unsigned int f=face.second+1; f<GeometryInfo<dimension>::faces_per_cell; ++f)
398  if (!surface_only || face.first->face(f)->at_boundary())
399  // yup, that is so, so return it
400  {
401  face.second = f;
402  return face;
403  }
404 
405  // otherwise find the next active cell that has a face on the boundary
406 
407  // convert the iterator to an active_iterator and advance this to the next
408  // active cell
409  typename Triangulation<dimension,space_dimension>::active_cell_iterator active_cell = face.first;
410 
411  // increase face pointer by one
412  ++active_cell;
413 
414  // while there are active cells
415  while (active_cell != this->triangulation->end())
416  {
417  // check all the faces of this active cell. but skip it altogether
418  // if it isn't locally owned
419  if (active_cell->is_locally_owned())
420  for (unsigned int f=0; f<GeometryInfo<dimension>::faces_per_cell; ++f)
421  if (!surface_only || active_cell->face(f)->at_boundary())
422  {
423  face.first = active_cell;
424  face.second = f;
425  return face;
426  }
427 
428  // the present cell had no faces on the boundary (or was not locally
429  // owned), so check next cell
430  ++active_cell;
431  }
432 
433  // we fell off the edge, so return with invalid pointer
434  face.first = this->triangulation->end();
435  face.second = 0;
436  return face;
437 }
438 
439 
440 
441 // explicit instantiations
442 #include "data_out_faces.inst"
443 
444 DEAL_II_NAMESPACE_CLOSE
TriaActiveIterator< CellAccessor< dim, spacedim > > active_cell_iterator
Definition: tria.h:1527
Shape function values.
static ::ExceptionBase & ExcNoTriangulationSelected()
void build_one_patch(const FaceDescriptor *cell_and_face, internal::DataOutFacesImplementation::ParallelData< dimension, dimension > &data, DataOutBase::Patch< dimension-1, space_dimension > &patch)
const FiniteElement< dim, spacedim > & get_fe() const
Transformed quadrature points.
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
Point< spacedim > vertices[GeometryInfo< dim >::vertices_per_cell]
virtual FaceDescriptor first_face()
const std::vector< Point< spacedim > > & get_quadrature_points() const
static ::ExceptionBase & ExcMessage(std::string arg1)
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
Definition: exceptions.h:1142
UpdateFlags
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
Definition: dof_tools.h:46
virtual FaceDescriptor next_face(const FaceDescriptor &face)
unsigned int n_subdivisions
Second derivatives of shape functions.
const std::vector< Tensor< 1, spacedim > > & get_all_normal_vectors() const
Definition: fe_values.cc:3762
DataOutFaces(const bool surface_only=true)
const unsigned int n_quadrature_points
Definition: fe_values.h:1840
std::pair< cell_iterator, unsigned int > FaceDescriptor
size_type size(const unsigned int i) const
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:1106
static ::ExceptionBase & ExcNotImplemented()
virtual void build_patches(const unsigned int n_subdivisions=0)
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
Table< 2, float > data
static ::ExceptionBase & ExcInternalError()
virtual UpdateFlags get_needed_update_flags() const =0