Reference documentation for deal.II version 9.1.0-pre
data_out_stack.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/memory_consumption.h>
17 #include <deal.II/base/quadrature_lib.h>
18 
19 #include <deal.II/dofs/dof_accessor.h>
20 #include <deal.II/dofs/dof_handler.h>
21 
22 #include <deal.II/fe/fe.h>
23 #include <deal.II/fe/fe_values.h>
24 #include <deal.II/fe/mapping_q1.h>
25 
26 #include <deal.II/grid/tria_iterator.h>
27 
28 #include <deal.II/hp/fe_values.h>
29 
30 #include <deal.II/lac/block_vector.h>
31 #include <deal.II/lac/vector.h>
32 
33 #include <deal.II/numerics/data_out_stack.h>
34 
35 #include <sstream>
36 
37 DEAL_II_NAMESPACE_OPEN
38 
39 
40 template <int dim, int spacedim, typename DoFHandlerType>
41 std::size_t
43  const
44 {
47 }
48 
49 
50 
51 template <int dim, int spacedim, typename DoFHandlerType>
52 void
54  const double p,
55  const double dp)
56 {
57  parameter = p;
58  parameter_step = dp;
59 
60  // check whether the user called finish_parameter_value() at the end of the
61  // previous parameter step
62  //
63  // this is to prevent serious waste of memory
64  for (typename std::vector<DataVector>::const_iterator i = dof_data.begin();
65  i != dof_data.end();
66  ++i)
67  Assert(i->data.size() == 0, ExcDataNotCleared());
68  for (typename std::vector<DataVector>::const_iterator i = cell_data.begin();
69  i != cell_data.end();
70  ++i)
71  Assert(i->data.size() == 0, ExcDataNotCleared());
72 }
73 
74 
75 template <int dim, int spacedim, typename DoFHandlerType>
76 void
78  const DoFHandlerType &dof)
79 {
80  // Check consistency of redundant
81  // template parameter
82  Assert(dim == DoFHandlerType::dimension,
83  ExcDimensionMismatch(dim, DoFHandlerType::dimension));
84 
85  dof_handler = &dof;
86 }
87 
88 
89 template <int dim, int spacedim, typename DoFHandlerType>
90 void
92  const std::string &name,
93  const VectorType vector_type)
94 {
95  std::vector<std::string> names;
96  names.push_back(name);
97  declare_data_vector(names, vector_type);
98 }
99 
100 
101 template <int dim, int spacedim, typename DoFHandlerType>
102 void
104  const std::vector<std::string> &names,
105  const VectorType vector_type)
106 {
107  // make sure this function is
108  // not called after some parameter
109  // values have already been
110  // processed
111  Assert(patches.size() == 0, ExcDataAlreadyAdded());
112 
113  // also make sure that no name is
114  // used twice
115  for (std::vector<std::string>::const_iterator name = names.begin();
116  name != names.end();
117  ++name)
118  {
119  for (typename std::vector<DataVector>::const_iterator data_set =
120  dof_data.begin();
121  data_set != dof_data.end();
122  ++data_set)
123  for (unsigned int i = 0; i < data_set->names.size(); ++i)
124  Assert(*name != data_set->names[i], ExcNameAlreadyUsed(*name));
125 
126  for (typename std::vector<DataVector>::const_iterator data_set =
127  cell_data.begin();
128  data_set != cell_data.end();
129  ++data_set)
130  for (unsigned int i = 0; i < data_set->names.size(); ++i)
131  Assert(*name != data_set->names[i], ExcNameAlreadyUsed(*name));
132  };
133 
134  switch (vector_type)
135  {
136  case dof_vector:
137  dof_data.emplace_back();
138  dof_data.back().names = names;
139  break;
140 
141  case cell_vector:
142  cell_data.emplace_back();
143  cell_data.back().names = names;
144  break;
145  };
146 }
147 
148 
149 template <int dim, int spacedim, typename DoFHandlerType>
150 template <typename number>
151 void
153  const Vector<number> &vec,
154  const std::string & name)
155 {
156  const unsigned int n_components = dof_handler->get_fe(0).n_components();
157 
158  std::vector<std::string> names;
159  // if only one component or vector
160  // is cell vector: we only need one
161  // name
162  if ((n_components == 1) ||
163  (vec.size() == dof_handler->get_triangulation().n_active_cells()))
164  {
165  names.resize(1, name);
166  }
167  else
168  // otherwise append _i to the
169  // given name
170  {
171  names.resize(n_components);
172  for (unsigned int i = 0; i < n_components; ++i)
173  {
174  std::ostringstream namebuf;
175  namebuf << '_' << i;
176  names[i] = name + namebuf.str();
177  }
178  }
179 
180  add_data_vector(vec, names);
181 }
182 
183 
184 template <int dim, int spacedim, typename DoFHandlerType>
185 template <typename number>
186 void
188  const Vector<number> & vec,
189  const std::vector<std::string> &names)
190 {
191  Assert(dof_handler != nullptr,
193  // either cell data and one name,
194  // or dof data and n_components names
195  Assert(((vec.size() == dof_handler->get_triangulation().n_active_cells()) &&
196  (names.size() == 1)) ||
197  ((vec.size() == dof_handler->n_dofs()) &&
198  (names.size() == dof_handler->get_fe(0).n_components())),
200  names.size(), dof_handler->get_fe(0).n_components()));
201  for (unsigned int i = 0; i < names.size(); ++i)
202  Assert(names[i].find_first_not_of("abcdefghijklmnopqrstuvwxyz"
203  "ABCDEFGHIJKLMNOPQRSTUVWXYZ"
204  "0123456789_<>()") == std::string::npos,
206  names[i],
207  names[i].find_first_not_of("abcdefghijklmnopqrstuvwxyz"
208  "ABCDEFGHIJKLMNOPQRSTUVWXYZ"
209  "0123456789_<>()")));
210 
211  if (vec.size() == dof_handler->n_dofs())
212  {
213  typename std::vector<DataVector>::iterator data_vector = dof_data.begin();
214  for (; data_vector != dof_data.end(); ++data_vector)
215  if (data_vector->names == names)
216  {
217  data_vector->data.reinit(vec.size());
218  std::copy(vec.begin(), vec.end(), data_vector->data.begin());
219  return;
220  };
221 
222  // ok. not found. there is a
223  // slight chance that
224  // n_dofs==n_cells, so only
225  // bomb out if the next if
226  // statement will not be run
227  if (dof_handler->n_dofs() !=
228  dof_handler->get_triangulation().n_active_cells())
229  Assert(false, ExcVectorNotDeclared(names[0]));
230  }
231 
232  // search cell data
233  if ((vec.size() != dof_handler->n_dofs()) ||
234  (dof_handler->n_dofs() ==
235  dof_handler->get_triangulation().n_active_cells()))
236  {
237  typename std::vector<DataVector>::iterator data_vector =
238  cell_data.begin();
239  for (; data_vector != cell_data.end(); ++data_vector)
240  if (data_vector->names == names)
241  {
242  data_vector->data.reinit(vec.size());
243  std::copy(vec.begin(), vec.end(), data_vector->data.begin());
244  return;
245  };
246  Assert(false, ExcVectorNotDeclared(names[0]));
247  };
248 
249  // we have either return or Assert
250  // statements above, so shouldn't
251  // get here!
252  Assert(false, ExcInternalError());
253 }
254 
255 
256 template <int dim, int spacedim, typename DoFHandlerType>
257 void
259  const unsigned int nnnn_subdivisions)
260 {
261  // this is mostly copied from the
262  // DataOut class
263  unsigned int n_subdivisions =
264  (nnnn_subdivisions != 0) ? nnnn_subdivisions : this->default_subdivisions;
265 
266  Assert(n_subdivisions >= 1,
268  n_subdivisions));
269  Assert(dof_handler != nullptr,
271 
272  this->validate_dataset_names();
273 
274  const unsigned int n_components = dof_handler->get_fe(0).n_components();
275  const unsigned int n_datasets =
276  dof_data.size() * n_components + cell_data.size();
277 
278  // first count the cells we want to
279  // create patches of and make sure
280  // there is enough memory for that
281  unsigned int n_patches = 0;
282  for (typename DoFHandlerType::active_cell_iterator cell =
283  dof_handler->begin_active();
284  cell != dof_handler->end();
285  ++cell)
286  ++n_patches;
287 
288 
289  // before we start the loop:
290  // create a quadrature rule that
291  // actually has the points on this
292  // patch, and an object that
293  // extracts the data on each
294  // cell to these points
295  QTrapez<1> q_trapez;
296  QIterated<dim> patch_points(q_trapez, n_subdivisions);
297 
298  // create collection objects from
299  // single quadratures,
300  // and finite elements. if we have
301  // an hp DoFHandler,
302  // dof_handler.get_fe() returns a
303  // collection of which we do a
304  // shallow copy instead
305  const hp::QCollection<dim> q_collection(patch_points);
306  const hp::FECollection<dim> &fe_collection = dof_handler->get_fe_collection();
307 
308  hp::FEValues<dim> x_fe_patch_values(fe_collection,
309  q_collection,
310  update_values);
311 
312  const unsigned int n_q_points = patch_points.size();
313  std::vector<double> patch_values(n_q_points);
314  std::vector<Vector<double>> patch_values_system(n_q_points,
315  Vector<double>(n_components));
316 
317  // add the required number of
318  // patches. first initialize a template
319  // patch with n_q_points (in the plane
320  // of the cells) times n_subdivisions+1 (in
321  // the time direction) points
323  default_patch.n_subdivisions = n_subdivisions;
324  default_patch.data.reinit(n_datasets, n_q_points * (n_subdivisions + 1));
325  patches.insert(patches.end(), n_patches, default_patch);
326 
327  // now loop over all cells and
328  // actually create the patches
329  typename std::vector<::DataOutBase::Patch<dim + 1, dim + 1>>::iterator
330  patch = patches.begin() + (patches.size() - n_patches);
331  unsigned int cell_number = 0;
332  for (typename DoFHandlerType::active_cell_iterator cell =
333  dof_handler->begin_active();
334  cell != dof_handler->end();
335  ++cell, ++patch, ++cell_number)
336  {
337  Assert(cell->is_locally_owned(), ExcNotImplemented());
338 
339  Assert(patch != patches.end(), ExcInternalError());
340 
341  // first fill in the vertices of the patch
342 
343  // Patches are organized such
344  // that the parameter direction
345  // is the last
346  // coordinate. Thus, vertices
347  // are two copies of the space
348  // patch, one at parameter-step
349  // and one at parameter.
350  switch (dim)
351  {
352  case 1:
353  patch->vertices[0] =
354  Point<dim + 1>(cell->vertex(0)(0), parameter - parameter_step);
355  patch->vertices[1] =
356  Point<dim + 1>(cell->vertex(1)(0), parameter - parameter_step);
357  patch->vertices[2] = Point<dim + 1>(cell->vertex(0)(0), parameter);
358  patch->vertices[3] = Point<dim + 1>(cell->vertex(1)(0), parameter);
359  break;
360 
361  case 2:
362  patch->vertices[0] = Point<dim + 1>(cell->vertex(0)(0),
363  cell->vertex(0)(1),
365  patch->vertices[1] = Point<dim + 1>(cell->vertex(1)(0),
366  cell->vertex(1)(1),
368  patch->vertices[2] = Point<dim + 1>(cell->vertex(2)(0),
369  cell->vertex(2)(1),
371  patch->vertices[3] = Point<dim + 1>(cell->vertex(3)(0),
372  cell->vertex(3)(1),
374  patch->vertices[4] =
375  Point<dim + 1>(cell->vertex(0)(0), cell->vertex(0)(1), parameter);
376  patch->vertices[5] =
377  Point<dim + 1>(cell->vertex(1)(0), cell->vertex(1)(1), parameter);
378  patch->vertices[6] =
379  Point<dim + 1>(cell->vertex(2)(0), cell->vertex(2)(1), parameter);
380  patch->vertices[7] =
381  Point<dim + 1>(cell->vertex(3)(0), cell->vertex(3)(1), parameter);
382  break;
383 
384  default:
385  Assert(false, ExcNotImplemented());
386  };
387 
388 
389  // now fill in the data values.
390  // note that the required order is
391  // with highest coordinate running
392  // fastest, we need to enter each
393  // value (n_subdivisions+1) times
394  // in succession
395  if (n_datasets > 0)
396  {
397  x_fe_patch_values.reinit(cell);
398  const FEValues<dim> &fe_patch_values =
399  x_fe_patch_values.get_present_fe_values();
400 
401  // first fill dof_data
402  for (unsigned int dataset = 0; dataset < dof_data.size(); ++dataset)
403  {
404  if (n_components == 1)
405  {
406  fe_patch_values.get_function_values(dof_data[dataset].data,
407  patch_values);
408  for (unsigned int i = 0; i < n_subdivisions + 1; ++i)
409  for (unsigned int q = 0; q < n_q_points; ++q)
410  patch->data(dataset, q + n_q_points * i) =
411  patch_values[q];
412  }
413  else
414  // system of components
415  {
416  fe_patch_values.get_function_values(dof_data[dataset].data,
417  patch_values_system);
418  for (unsigned int component = 0; component < n_components;
419  ++component)
420  for (unsigned int i = 0; i < n_subdivisions + 1; ++i)
421  for (unsigned int q = 0; q < n_q_points; ++q)
422  patch->data(dataset * n_components + component,
423  q + n_q_points * i) =
424  patch_values_system[q](component);
425  }
426  }
427 
428  // then do the cell data
429  for (unsigned int dataset = 0; dataset < cell_data.size(); ++dataset)
430  {
431  const double value = cell_data[dataset].data(cell_number);
432  for (unsigned int q = 0; q < n_q_points; ++q)
433  for (unsigned int i = 0; i < n_subdivisions + 1; ++i)
434  patch->data(dataset + dof_data.size() * n_components,
435  q * (n_subdivisions + 1) + i) = value;
436  }
437  }
438  }
439 }
440 
441 
442 template <int dim, int spacedim, typename DoFHandlerType>
443 void
445 {
446  // release lock on dof handler
447  dof_handler = nullptr;
448  for (typename std::vector<DataVector>::iterator i = dof_data.begin();
449  i != dof_data.end();
450  ++i)
451  i->data.reinit(0);
452 
453  for (typename std::vector<DataVector>::iterator i = cell_data.begin();
454  i != cell_data.end();
455  ++i)
456  i->data.reinit(0);
457 }
458 
459 
460 
461 template <int dim, int spacedim, typename DoFHandlerType>
462 std::size_t
464 {
472 }
473 
474 
475 
476 template <int dim, int spacedim, typename DoFHandlerType>
477 const std::vector<::DataOutBase::Patch<dim + 1, dim + 1>> &
479 {
480  return patches;
481 }
482 
483 
484 
485 template <int dim, int spacedim, typename DoFHandlerType>
486 std::vector<std::string>
488 {
489  std::vector<std::string> names;
490  for (typename std::vector<DataVector>::const_iterator dataset =
491  dof_data.begin();
492  dataset != dof_data.end();
493  ++dataset)
494  names.insert(names.end(), dataset->names.begin(), dataset->names.end());
495  for (typename std::vector<DataVector>::const_iterator dataset =
496  cell_data.begin();
497  dataset != cell_data.end();
498  ++dataset)
499  names.insert(names.end(), dataset->names.begin(), dataset->names.end());
500 
501  return names;
502 }
503 
504 
505 
506 // explicit instantiations
507 #include "data_out_stack.inst"
508 
509 
510 DEAL_II_NAMESPACE_CLOSE
static::ExceptionBase & ExcNameAlreadyUsed(std::string arg1)
void get_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
Definition: fe_values.cc:3657
Shape function values.
std::vector< DataVector > cell_data
size_type size() const
static::ExceptionBase & ExcDataAlreadyAdded()
void attach_dof_handler(const DoFHandlerType &dof_handler)
std::vector< DataVector > dof_data
unsigned int default_subdivisions
static::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
iterator end()
Definition: point.h:106
virtual std::vector< std::string > get_dataset_names() const override
static::ExceptionBase & ExcDataNotCleared()
SmartPointer< const DoFHandlerType, DataOutStack< dim, spacedim, DoFHandlerType > > dof_handler
double parameter_step
std::vector< std::string > names
unsigned int size() const
std::size_t memory_consumption() const
static::ExceptionBase & ExcNoDoFHandlerSelected()
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
void declare_data_vector(const std::string &name, const VectorType vector_type)
#define Assert(cond, exc)
Definition: exceptions.h:1227
static::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static::ExceptionBase & ExcInvalidNumberOfNames(int arg1, int arg2)
virtual const std::vector<::DataOutBase::Patch< dim+1, dim+1 > > & get_patches() const override
static::ExceptionBase & ExcInvalidCharacter(std::string arg1, size_t arg2)
unsigned int n_subdivisions
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType, lda >> cell, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
Definition: fe_values.cc:142
iterator begin()
static::ExceptionBase & ExcVectorNotDeclared(std::string arg1)
void build_patches(const unsigned int n_subdivisions=0)
std::vector<::DataOutBase::Patch< dim+1, dim+1 > > patches
static::ExceptionBase & ExcNotImplemented()
const ::FEValues< dim, spacedim > & get_present_fe_values() const
void finish_parameter_value()
Table< 2, float > data
void validate_dataset_names() const
void new_parameter_value(const double parameter_value, const double parameter_step)
void add_data_vector(const Vector< number > &vec, const std::string &name)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static::ExceptionBase & ExcInternalError()
std::size_t memory_consumption() const