Reference documentation for deal.II version 8.5.1
data_out_rotation.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2017 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/grid/tria.h>
21 #include <deal.II/dofs/dof_handler.h>
22 #include <deal.II/dofs/dof_accessor.h>
23 #include <deal.II/grid/tria_iterator.h>
24 #include <deal.II/fe/fe.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 #include <deal.II/numerics/data_out_rotation.h>
29 
30 #include <cmath>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 
35 //TODO: Update documentation
36 //TODO: Unify code for dimensions
37 
38 
39 //TODO: build_some_patches isn't going to work if first_cell/next_cell
40 //don't iterate over all cells and if cell data is requested. in that
41 //case, we need to calculate cell_number as in the DataOut class
42 
43 // Not implemented for 3D
44 
45 
46 namespace internal
47 {
48  namespace DataOutRotation
49  {
50  template <int dim, int spacedim>
51  ParallelData<dim,spacedim>::
52  ParallelData (const unsigned int n_datasets,
53  const unsigned int n_subdivisions,
54  const unsigned int n_patches_per_circle,
55  const std::vector<unsigned int> &n_postprocessor_outputs,
56  const Mapping<dim,spacedim> &mapping,
57  const std::vector<std_cxx11::shared_ptr<::hp::FECollection<dim,spacedim> > > &finite_elements,
58  const UpdateFlags update_flags)
59  :
60  internal::DataOut::
61  ParallelDataBase<dim,spacedim> (n_datasets,
62  n_subdivisions,
63  n_postprocessor_outputs,
64  mapping,
65  finite_elements,
66  update_flags,
67  false),
68  n_patches_per_circle (n_patches_per_circle)
69  {}
70 
71 
72 
77  template <int dim, int spacedim>
78  void
79  append_patch_to_list (const std::vector<DataOutBase::Patch<dim+1,spacedim+1> > &new_patches,
80  std::vector<DataOutBase::Patch<dim+1,spacedim+1> > &patches)
81  {
82  for (unsigned int i=0; i<new_patches.size(); ++i)
83  {
84  patches.push_back (new_patches[i]);
85  patches.back().patch_index = patches.size()-1;
86  }
87  }
88  }
89 }
90 
91 
92 
93 template <int dim, typename DoFHandlerType>
94 void
98  std::vector<DataOutBase::Patch<dimension+1,space_dimension+1> > &my_patches)
99 {
100  if (dim == 3)
101  {
102  // would this function make any sense after all? who would want to
103  // output/compute in four space dimensions?
104  Assert (false, ExcNotImplemented());
105  return;
106  }
107 
108  Assert ((*cell)->is_locally_owned(),
110 
111  const unsigned int n_patches_per_circle = data.n_patches_per_circle;
112 
113  // another abbreviation denoting the number of q_points in each direction
114  const unsigned int n_points = data.n_subdivisions+1;
115 
116  // set up an array that holds the directions in the plane of rotation in
117  // which we will put points in the whole domain (not the rotationally
118  // reduced one in which the computation took place. for simplicity add the
119  // initial direction at the end again
120  std::vector<Point<dimension+1> > angle_directions (n_patches_per_circle+1);
121  for (unsigned int i=0; i<=n_patches_per_circle; ++i)
122  {
123  angle_directions[i][dimension-1] = std::cos(2*numbers::PI *
124  i/n_patches_per_circle);
125  angle_directions[i][dimension] = std::sin(2*numbers::PI *
126  i/n_patches_per_circle);
127  }
128 
129  for (unsigned int angle=0; angle<n_patches_per_circle; ++angle)
130  {
131  // first compute the vertices of the patch. note that they will have to
132  // be computed from the vertices of the cell, which has one dimension
133  // less, however.
134  switch (dimension)
135  {
136  case 1:
137  {
138  const double r1 = (*cell)->vertex(0)(0),
139  r2 = (*cell)->vertex(1)(0);
140  Assert (r1 >= 0, ExcRadialVariableHasNegativeValues(r1));
141  Assert (r2 >= 0, ExcRadialVariableHasNegativeValues(r2));
142 
143  my_patches[angle].vertices[0] = r1*angle_directions[angle];
144  my_patches[angle].vertices[1] = r2*angle_directions[angle];
145  my_patches[angle].vertices[2] = r1*angle_directions[angle+1];
146  my_patches[angle].vertices[3] = r2*angle_directions[angle+1];
147 
148  break;
149  };
150 
151  case 2:
152  {
153  for (unsigned int vertex=0;
154  vertex<GeometryInfo<dimension>::vertices_per_cell;
155  ++vertex)
156  {
157  const Point<dimension> v = (*cell)->vertex(vertex);
158 
159  // make sure that the radial variable does attain negative
160  // values
161  Assert (v(0) >= 0, ExcRadialVariableHasNegativeValues(v(0)));
162 
163  // now set the vertices of the patch
164  my_patches[angle].vertices[vertex] = v(0) * angle_directions[angle];
165  my_patches[angle].vertices[vertex][0] = v(1);
166 
167  my_patches[angle].vertices[vertex+GeometryInfo<dimension>::vertices_per_cell]
168  = v(0) * angle_directions[angle+1];
169  my_patches[angle].vertices[vertex+GeometryInfo<dimension>::vertices_per_cell][0]
170  = v(1);
171  };
172 
173  break;
174  };
175 
176  default:
177  Assert (false, ExcNotImplemented());
178  };
179 
180  // then fill in data
181  if (data.n_datasets > 0)
182  {
183  unsigned int offset=0;
184 
185  data.reinit_all_fe_values(this->dof_data, *cell);
186  // first fill dof_data
187  for (unsigned int dataset=0; dataset<this->dof_data.size(); ++dataset)
188  {
189  const FEValuesBase<dimension> &fe_patch_values
190  = data.get_present_fe_values(dataset);
191  const unsigned int n_components
192  = fe_patch_values.get_fe().n_components();
193  const DataPostprocessor<dim> *postprocessor=this->dof_data[dataset]->postprocessor;
194  if (postprocessor != 0)
195  {
196  // we have to postprocess the
197  // data, so determine, which
198  // fields have to be updated
199  const UpdateFlags update_flags=postprocessor->get_needed_update_flags();
200 
201  if (n_components == 1)
202  {
203  // at each point there is
204  // only one component of
205  // value, gradient etc.
206  if (update_flags & update_values)
207  this->dof_data[dataset]->get_function_values (fe_patch_values,
208  data.patch_values_scalar.solution_values);
209  if (update_flags & update_gradients)
210  this->dof_data[dataset]->get_function_gradients (fe_patch_values,
211  data.patch_values_scalar.solution_gradients);
212  if (update_flags & update_hessians)
213  this->dof_data[dataset]->get_function_hessians (fe_patch_values,
214  data.patch_values_scalar.solution_hessians);
215 
216  if (update_flags & update_quadrature_points)
217  data.patch_values_scalar.evaluation_points = fe_patch_values.get_quadrature_points();
218 
219  const typename DoFHandlerType::active_cell_iterator dh_cell(&(*cell)->get_triangulation(),
220  (*cell)->level(),
221  (*cell)->index(),
222  this->dof_data[dataset]->dof_handler);
223  data.patch_values_scalar.template set_cell<DoFHandlerType> (dh_cell);
224 
225  postprocessor->
226  evaluate_scalar_field(data.patch_values_scalar,
227  data.postprocessed_values[dataset]);
228  }
229  else
230  {
231  data.resize_system_vectors(n_components);
232 
233  // at each point there is a vector valued function and
234  // its derivative...
235  if (update_flags & update_values)
236  this->dof_data[dataset]->get_function_values (fe_patch_values,
237  data.patch_values_system.solution_values);
238  if (update_flags & update_gradients)
239  this->dof_data[dataset]->get_function_gradients (fe_patch_values,
240  data.patch_values_system.solution_gradients);
241  if (update_flags & update_hessians)
242  this->dof_data[dataset]->get_function_hessians (fe_patch_values,
243  data.patch_values_system.solution_hessians);
244 
245  if (update_flags & update_quadrature_points)
246  data.patch_values_system.evaluation_points = fe_patch_values.get_quadrature_points();
247 
248  const typename DoFHandlerType::active_cell_iterator dh_cell(&(*cell)->get_triangulation(),
249  (*cell)->level(),
250  (*cell)->index(),
251  this->dof_data[dataset]->dof_handler);
252  data.patch_values_system.template set_cell<DoFHandlerType> (dh_cell);
253 
254  postprocessor->
255  evaluate_vector_field(data.patch_values_system,
256  data.postprocessed_values[dataset]);
257  }
258 
259  for (unsigned int component=0;
260  component<this->dof_data[dataset]->n_output_variables;
261  ++component)
262  {
263  switch (dimension)
264  {
265  case 1:
266  for (unsigned int x=0; x<n_points; ++x)
267  for (unsigned int y=0; y<n_points; ++y)
268  my_patches[angle].data(offset+component,
269  x*n_points + y)
270  = data.postprocessed_values[dataset][x](component);
271  break;
272 
273  case 2:
274  for (unsigned int x=0; x<n_points; ++x)
275  for (unsigned int y=0; y<n_points; ++y)
276  for (unsigned int z=0; z<n_points; ++z)
277  my_patches[angle].data(offset+component,
278  x*n_points*n_points +
279  y*n_points +
280  z)
281  = data.postprocessed_values[dataset][x*n_points+z](component);
282  break;
283 
284  default:
285  Assert (false, ExcNotImplemented());
286  }
287  }
288  }
289  else if (n_components == 1)
290  {
291  this->dof_data[dataset]->get_function_values (fe_patch_values,
292  data.patch_values_scalar.solution_values);
293 
294  switch (dimension)
295  {
296  case 1:
297  for (unsigned int x=0; x<n_points; ++x)
298  for (unsigned int y=0; y<n_points; ++y)
299  my_patches[angle].data(offset,
300  x*n_points + y)
301  = data.patch_values_scalar.solution_values[x];
302  break;
303 
304  case 2:
305  for (unsigned int x=0; x<n_points; ++x)
306  for (unsigned int y=0; y<n_points; ++y)
307  for (unsigned int z=0; z<n_points; ++z)
308  my_patches[angle].data(offset,
309  x*n_points*n_points +
310  y +
311  z*n_points)
312  = data.patch_values_scalar.solution_values[x*n_points+z];
313  break;
314 
315  default:
316  Assert (false, ExcNotImplemented());
317  }
318  }
319  else
320  // system of components
321  {
322  data.resize_system_vectors(n_components);
323  this->dof_data[dataset]->get_function_values (fe_patch_values,
324  data.patch_values_system.solution_values);
325 
326  for (unsigned int component=0; component<n_components;
327  ++component)
328  {
329  switch (dimension)
330  {
331  case 1:
332  for (unsigned int x=0; x<n_points; ++x)
333  for (unsigned int y=0; y<n_points; ++y)
334  my_patches[angle].data(offset+component,
335  x*n_points + y)
336  = data.patch_values_system.solution_values[x](component);
337  break;
338 
339  case 2:
340  for (unsigned int x=0; x<n_points; ++x)
341  for (unsigned int y=0; y<n_points; ++y)
342  for (unsigned int z=0; z<n_points; ++z)
343  my_patches[angle].data(offset+component,
344  x*n_points*n_points +
345  y*n_points +
346  z)
347  = data.patch_values_system.solution_values[x*n_points+z](component);
348  break;
349 
350  default:
351  Assert (false, ExcNotImplemented());
352  }
353  }
354  }
355  offset+=this->dof_data[dataset]->n_output_variables;
356  }
357 
358  // then do the cell data
359  for (unsigned int dataset=0; dataset<this->cell_data.size(); ++dataset)
360  {
361  // we need to get at the number of the cell to which this face
362  // belongs in order to access the cell data. this is not readily
363  // available, so choose the following rather inefficient way:
364  Assert ((*cell)->active(),
365  ExcMessage("Cell must be active for cell data"));
366  const unsigned int cell_number
367  = std::distance (this->triangulation->begin_active(),
369  const double value
370  = this->cell_data[dataset]->get_cell_data_value (cell_number);
371  switch (dimension)
372  {
373  case 1:
374  for (unsigned int x=0; x<n_points; ++x)
375  for (unsigned int y=0; y<n_points; ++y)
376  my_patches[angle].data(dataset+offset,
377  x*n_points +
378  y)
379  = value;
380  break;
381 
382  case 2:
383  for (unsigned int x=0; x<n_points; ++x)
384  for (unsigned int y=0; y<n_points; ++y)
385  for (unsigned int z=0; z<n_points; ++z)
386  my_patches[angle].data(dataset+offset,
387  x*n_points*n_points +
388  y*n_points +
389  z)
390  = value;
391  break;
392 
393  default:
394  Assert (false, ExcNotImplemented());
395  }
396  }
397  }
398  }
399 }
400 
401 
402 
403 template <int dim, typename DoFHandlerType>
404 void DataOutRotation<dim,DoFHandlerType>::build_patches (const unsigned int n_patches_per_circle,
405  const unsigned int nnnn_subdivisions)
406 {
407  // Check consistency of redundant
408  // template parameter
409  Assert (dim==dimension, ExcDimensionMismatch(dim, dimension));
410  Assert (this->triangulation != 0,
412 
413  const unsigned int n_subdivisions = (nnnn_subdivisions != 0)
414  ? nnnn_subdivisions
415  : this->default_subdivisions;
416  Assert (n_subdivisions >= 1,
418 
419  this->validate_dataset_names();
420 
421  unsigned int n_datasets=this->cell_data.size();
422  for (unsigned int i=0; i<this->dof_data.size(); ++i)
423  n_datasets+= this->dof_data[i]->n_output_variables;
424 
426  for (unsigned int i=0; i<this->dof_data.size(); ++i)
427  if (this->dof_data[i]->postprocessor)
428  update_flags |= this->dof_data[i]->postprocessor->get_needed_update_flags();
429  // perhaps update_normal_vectors is present,
430  // which would only be useful on faces, but
431  // we may not use it here.
432  Assert (!(update_flags & update_normal_vectors),
433  ExcMessage("The update of normal vectors may not be requested for "
434  "evaluation of data on cells via DataPostprocessor."));
435 
436  // first count the cells we want to
437  // create patches of and make sure
438  // there is enough memory for that
439  std::vector<cell_iterator> all_cells;
440  for (cell_iterator cell=first_cell(); cell != this->triangulation->end();
441  cell = next_cell(cell))
442  all_cells.push_back (cell);
443 
444  // then also take into account that
445  // we want more than one patch to
446  // come out of every cell, as they
447  // are repeated around the axis of
448  // rotation
449  this->patches.clear();
450  this->patches.reserve (all_cells.size() * n_patches_per_circle);
451 
452 
453  std::vector<unsigned int> n_postprocessor_outputs (this->dof_data.size());
454  for (unsigned int dataset=0; dataset<this->dof_data.size(); ++dataset)
455  if (this->dof_data[dataset]->postprocessor)
456  n_postprocessor_outputs[dataset] = this->dof_data[dataset]->n_output_variables;
457  else
458  n_postprocessor_outputs[dataset] = 0;
459 
461  thread_data (n_datasets,
462  n_subdivisions, n_patches_per_circle,
463  n_postprocessor_outputs,
465  this->get_finite_elements(),
466  update_flags);
467  std::vector<DataOutBase::Patch<dimension+1,space_dimension+1> >
468  new_patches (n_patches_per_circle);
469  for (unsigned int i=0; i<new_patches.size(); ++i)
470  {
471  new_patches[i].n_subdivisions = n_subdivisions;
472  new_patches[i].data.reinit (n_datasets,
473  Utilities::fixed_power<dimension+1>(n_subdivisions+1));
474  }
475 
476  // now build the patches in parallel
477  WorkStream::run (&all_cells[0],
478  &all_cells[0]+all_cells.size(),
480  this, std_cxx11::_1, std_cxx11::_2, std_cxx11::_3),
481  std_cxx11::bind(&internal::DataOutRotation
482  ::append_patch_to_list<dim,space_dimension>,
483  std_cxx11::_1, std_cxx11::ref(this->patches)),
484  thread_data,
485  new_patches);
486 }
487 
488 
489 
490 template <int dim, typename DoFHandlerType>
493 {
494  return this->triangulation->begin_active ();
495 }
496 
497 
498 template <int dim, typename DoFHandlerType>
501 {
502  // convert the iterator to an
503  // active_iterator and advance
504  // this to the next active cell
506  ++active_cell;
507  return active_cell;
508 }
509 
510 
511 
512 // explicit instantiations
513 #include "data_out_rotation.inst"
514 
515 
516 DEAL_II_NAMESPACE_CLOSE
TriaActiveIterator< CellAccessor< dim, spacedim > > active_cell_iterator
Definition: tria.h:1499
Shape function values.
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
virtual cell_iterator next_cell(const cell_iterator &cell)
const FiniteElement< dim, spacedim > & get_fe() const
Transformed quadrature points.
Definition: point.h:89
const std::vector< Point< spacedim > > & get_quadrature_points() const
static const double PI
Definition: numbers.h:94
void build_one_patch(const cell_iterator *cell, internal::DataOutRotation::ParallelData< dimension, space_dimension > &data, std::vector< DataOutBase::Patch< dimension+1, space_dimension+1 > > &my_patches)
static ::ExceptionBase & ExcMessage(std::string arg1)
#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
Second derivatives of shape functions.
virtual cell_iterator first_cell()
virtual void build_patches(const unsigned int n_patches_per_circle, const unsigned int n_subdivisions=0)
DataOut_DoFData< DoFHandlerType, dimension+1 >::cell_iterator cell_iterator
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()
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
static ::ExceptionBase & ExcNoTriangulationSelected()
virtual UpdateFlags get_needed_update_flags() const =0