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