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