Reference documentation for deal.II version 8.4.2
data_out_faces.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 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/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 DataOutFaces
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_cxx11::shared_ptr<::hp::FECollection<dim,spacedim> > > &finite_elements,
44  const UpdateFlags update_flags)
45  :
46  internal::DataOut::
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,
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,
81  ExcNotImplemented());
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(),
94  ExcNotImplemented());
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 != 0)
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  // get normals, if needed. this is a geometrical information and
153  // thus does not depend on the number of components of the data
154  // vector
155  if (update_flags & update_normal_vectors)
156  {
157 //TODO: undo this copying when we can change the data type of
158 // data.patch_normals to Tensor<1,spacedim> as well
159  for (unsigned int q=0; q<this_fe_patch_values.n_quadrature_points; ++q)
160  data.patch_normals[q] = Point<dim>(this_fe_patch_values.get_all_normal_vectors()[q]);
161  }
162 
163  if (n_components == 1)
164  {
165  // at each point there is only one component of value,
166  // gradient etc.
167  if (update_flags & update_values)
168  this->dof_data[dataset]->get_function_values (this_fe_patch_values,
169  data.patch_values);
170  if (update_flags & update_gradients)
171  this->dof_data[dataset]->get_function_gradients (this_fe_patch_values,
172  data.patch_gradients);
173  if (update_flags & update_hessians)
174  this->dof_data[dataset]->get_function_hessians (this_fe_patch_values,
175  data.patch_hessians);
176 
177  if (update_flags & update_quadrature_points)
178  data.patch_evaluation_points = this_fe_patch_values.get_quadrature_points();
179 
180  postprocessor->
181  compute_derived_quantities_scalar(data.patch_values,
182  data.patch_gradients,
183  data.patch_hessians,
184  data.patch_normals,
185  data.patch_evaluation_points,
186  data.postprocessed_values[dataset]);
187  }
188  else
189  {
190  // at each point there is a vector valued function and its
191  // derivative...
192  data.resize_system_vectors(n_components);
193  if (update_flags & update_values)
194  this->dof_data[dataset]->get_function_values (this_fe_patch_values,
195  data.patch_values_system);
196  if (update_flags & update_gradients)
197  this->dof_data[dataset]->get_function_gradients (this_fe_patch_values,
198  data.patch_gradients_system);
199  if (update_flags & update_hessians)
200  this->dof_data[dataset]->get_function_hessians (this_fe_patch_values,
201  data.patch_hessians_system);
202 
203  if (update_flags & update_quadrature_points)
204  data.patch_evaluation_points = this_fe_patch_values.get_quadrature_points();
205 
206  postprocessor->
207  compute_derived_quantities_vector(data.patch_values_system,
208  data.patch_gradients_system,
209  data.patch_hessians_system,
210  data.patch_normals,
211  data.patch_evaluation_points,
212  data.postprocessed_values[dataset]);
213  }
214 
215  for (unsigned int q=0; q<n_q_points; ++q)
216  for (unsigned int component=0;
217  component<this->dof_data[dataset]->n_output_variables; ++component)
218  patch.data(offset+component,q)
219  = data.postprocessed_values[dataset][q](component);
220  }
221  else
222  // now we use the given data vector without modifications. again,
223  // we treat single component functions separately for efficiency
224  // reasons.
225  if (n_components == 1)
226  {
227  this->dof_data[dataset]->get_function_values (this_fe_patch_values,
228  data.patch_values);
229  for (unsigned int q=0; q<n_q_points; ++q)
230  patch.data(offset,q) = data.patch_values[q];
231  }
232  else
233  {
234  data.resize_system_vectors(n_components);
235  this->dof_data[dataset]->get_function_values (this_fe_patch_values,
236  data.patch_values_system);
237  for (unsigned int component=0; component<n_components;
238  ++component)
239  for (unsigned int q=0; q<n_q_points; ++q)
240  patch.data(offset+component,q) =
241  data.patch_values_system[q](component);
242  }
243  // increment the counter for the actual data record
244  offset+=this->dof_data[dataset]->n_output_variables;
245  }
246 
247  // then do the cell data
248  for (unsigned int dataset=0; dataset<this->cell_data.size(); ++dataset)
249  {
250  // we need to get at the number of the cell to which this face
251  // belongs in order to access the cell data. this is not readily
252  // available, so choose the following rather inefficient way:
253  Assert (cell_and_face->first->active(),
254  ExcMessage("The current function is trying to generate cell-data output "
255  "for a face that does not belong to an active cell. This is "
256  "not supported."));
257  const unsigned int cell_number
258  = std::distance (this->triangulation->begin_active(),
260 
261  const double value
262  = this->cell_data[dataset]->get_cell_data_value (cell_number);
263  for (unsigned int q=0; q<n_q_points; ++q)
264  patch.data(dataset+offset,q) = value;
265  }
266  }
267 }
268 
269 
270 
271 
272 template <int dim, typename DoFHandlerType>
273 void DataOutFaces<dim,DoFHandlerType>::build_patches (const unsigned int n_subdivisions_)
274 {
276 }
277 
278 
279 
280 template <int dim, typename DoFHandlerType>
282  const unsigned int n_subdivisions_)
283 {
284  // Check consistency of redundant template parameter
285  Assert (dim==dimension, ExcDimensionMismatch(dim, dimension));
286 
287  const unsigned int n_subdivisions = (n_subdivisions_ != 0)
288  ? n_subdivisions_
289  : this->default_subdivisions;
290 
291  Assert (n_subdivisions >= 1,
292  Exceptions::DataOut::ExcInvalidNumberOfSubdivisions(n_subdivisions));
293 
294  Assert (this->triangulation != 0,
295  Exceptions::DataOut::ExcNoTriangulationSelected());
296 
297  unsigned int n_datasets = this->cell_data.size();
298  for (unsigned int i=0; i<this->dof_data.size(); ++i)
299  n_datasets += this->dof_data[i]->n_output_variables;
300 
301  // first count the cells we want to create patches of and make sure there is
302  // enough memory for that
303  std::vector<FaceDescriptor> all_faces;
304  for (FaceDescriptor face=first_face();
305  face.first != this->triangulation->end();
306  face = next_face(face))
307  all_faces.push_back (face);
308 
309  // clear the patches array and allocate the right number of elements
310  this->patches.clear ();
311  this->patches.reserve (all_faces.size());
312  Assert (this->patches.size() == 0, ExcInternalError());
313 
314 
315  std::vector<unsigned int> n_postprocessor_outputs (this->dof_data.size());
316  for (unsigned int dataset=0; dataset<this->dof_data.size(); ++dataset)
317  if (this->dof_data[dataset]->postprocessor)
318  n_postprocessor_outputs[dataset] = this->dof_data[dataset]->n_output_variables;
319  else
320  n_postprocessor_outputs[dataset] = 0;
321 
322  UpdateFlags update_flags=update_values;
323  for (unsigned int i=0; i<this->dof_data.size(); ++i)
324  if (this->dof_data[i]->postprocessor)
325  update_flags |= this->dof_data[i]->postprocessor->get_needed_update_flags();
326  update_flags |= update_quadrature_points;
327 
329  thread_data (n_datasets,
330  n_subdivisions,
331  n_postprocessor_outputs,
332  mapping,
333  this->get_finite_elements(),
334  update_flags);
336  sample_patch.n_subdivisions = n_subdivisions;
337  sample_patch.data.reinit (n_datasets,
338  Utilities::fixed_power<dimension-1>(n_subdivisions+1));
339 
340  // now build the patches in parallel
341  WorkStream::run (&all_faces[0],
342  &all_faces[0]+all_faces.size(),
344  this, std_cxx11::_1, std_cxx11::_2, std_cxx11::_3),
345  std_cxx11::bind(&internal::DataOutFaces::
346  append_patch_to_list<dim,space_dimension>,
347  std_cxx11::_1, std_cxx11::ref(this->patches)),
348  thread_data,
349  sample_patch);
350 }
351 
352 
353 
354 template <int dim, typename DoFHandlerType>
357 {
358  // simply find first active cell with a face on the boundary
360  for (; cell != this->triangulation->end(); ++cell)
361  if (cell->is_locally_owned())
362  for (unsigned int f=0; f<GeometryInfo<dimension>::faces_per_cell; ++f)
363  if (!surface_only || cell->face(f)->at_boundary())
364  return FaceDescriptor(cell, f);
365 
366  // just return an invalid descriptor if we haven't found a locally
367  // owned face. this can happen in parallel where all boundary
368  // faces are owned by other processors
369  return FaceDescriptor();
370 }
371 
372 
373 
374 template <int dim, typename DoFHandlerType>
377 {
378  FaceDescriptor face = old_face;
379 
380  // first check whether the present cell has more faces on the boundary. since
381  // we started with this face, its cell must clearly be locally owned
382  Assert (face.first->is_locally_owned(), ExcInternalError());
383  for (unsigned int f=face.second+1; f<GeometryInfo<dimension>::faces_per_cell; ++f)
384  if (!surface_only || face.first->face(f)->at_boundary())
385  // yup, that is so, so return it
386  {
387  face.second = f;
388  return face;
389  }
390 
391  // otherwise find the next active cell that has a face on the boundary
392 
393  // convert the iterator to an active_iterator and advance this to the next
394  // active cell
395  typename Triangulation<dimension,space_dimension>::active_cell_iterator active_cell = face.first;
396 
397  // increase face pointer by one
398  ++active_cell;
399 
400  // while there are active cells
401  while (active_cell != this->triangulation->end())
402  {
403  // check all the faces of this active cell. but skip it altogether
404  // if it isn't locally owned
405  if (active_cell->is_locally_owned())
406  for (unsigned int f=0; f<GeometryInfo<dimension>::faces_per_cell; ++f)
407  if (!surface_only || active_cell->face(f)->at_boundary())
408  {
409  face.first = active_cell;
410  face.second = f;
411  return face;
412  }
413 
414  // the present cell had no faces on the boundary (or was not locally
415  // owned), so check next cell
416  ++active_cell;
417  }
418 
419  // we fell off the edge, so return with invalid pointer
420  face.first = this->triangulation->end();
421  face.second = 0;
422  return face;
423 }
424 
425 
426 
427 // explicit instantiations
428 #include "data_out_faces.inst"
429 
430 DEAL_II_NAMESPACE_CLOSE
std::vector< std_cxx11::shared_ptr< internal::DataOut::DataEntryBase< DoFHandlerType > > > cell_data
TriaActiveIterator< CellAccessor< dim, spacedim > > active_cell_iterator
Definition: tria.h:1410
Shape function values.
static const unsigned int space_dimension
::ExceptionBase & ExcMessage(std::string arg1)
const FiniteElement< dim, spacedim > & get_fe() const
void build_one_patch(const FaceDescriptor *cell_and_face, internal::DataOutFaces::ParallelData< dimension, dimension > &data, DataOutBase::Patch< dimension-1, space_dimension > &patch)
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: tria.cc:10397
Transformed quadrature points.
const bool surface_only
Point< spacedim > vertices[GeometryInfo< dim >::vertices_per_cell]
static const unsigned int dimension
virtual FaceDescriptor first_face()
cell_iterator end() const
Definition: tria.cc:10465
const std::vector< Point< spacedim > > & get_quadrature_points() const
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
Definition: exceptions.h:294
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:52
std::vector< Patch > patches
virtual FaceDescriptor next_face(const FaceDescriptor &face)
Second derivatives of shape functions.
const std::vector< Tensor< 1, spacedim > > & get_all_normal_vectors() const
Definition: fe_values.cc:3394
DataOutFaces(const bool surface_only=true)
const unsigned int n_quadrature_points
Definition: fe_values.h:1457
std::pair< cell_iterator, unsigned int > FaceDescriptor
std::vector< std_cxx11::shared_ptr<::hp::FECollection< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > > get_finite_elements() const
static const unsigned int space_dim
Shape function gradients.
Normal vectors.
SmartPointer< const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > triangulation
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:1119
std::vector< std_cxx11::shared_ptr< internal::DataOut::DataEntryBase< DoFHandlerType > > > dof_data
unsigned int size(const unsigned int i) const
virtual void build_patches(const unsigned int n_subdivisions=0)
Table< 2, float > data
virtual UpdateFlags get_needed_update_flags() const =0