Reference documentation for deal.II version 8.5.1
data_postprocessor.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2007 - 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 #ifndef dealii__data_postprocessor_h
17 #define dealii__data_postprocessor_h
18 
19 
20 
21 #include <deal.II/base/subscriptor.h>
22 #include <deal.II/base/tensor.h>
23 #include <deal.II/base/point.h>
24 #include <deal.II/lac/vector.h>
25 #include <deal.II/fe/fe_update_flags.h>
26 #include <deal.II/numerics/data_component_interpretation.h>
27 
28 #include <boost/any.hpp>
29 
30 #include <vector>
31 #include <string>
32 
33 DEAL_II_NAMESPACE_OPEN
34 
35 
41 {
144  template <int spacedim>
146  {
163  std::vector<Tensor<1, spacedim> > normals;
164 
177  std::vector<Point<spacedim> > evaluation_points;
178 
187  template <typename DoFHandlerType>
188  void
189  set_cell (const typename DoFHandlerType::cell_iterator &cell);
190 
196  template <typename DoFHandlerType>
197  typename DoFHandlerType::cell_iterator
198  get_cell () const;
199 
200  private:
208  boost::any cell;
209  };
210 
224  template <int spacedim>
225  struct Scalar : public CommonInputs<spacedim>
226  {
232  std::vector<double> solution_values;
233 
247  std::vector<Tensor<1, spacedim> > solution_gradients;
248 
262  std::vector<Tensor<2, spacedim> > solution_hessians;
263  };
264 
265 
266 
280  template <int spacedim>
281  struct Vector : public CommonInputs<spacedim>
282  {
292  std::vector<::Vector<double> > solution_values;
293 
311  std::vector<std::vector<Tensor<1, spacedim> > > solution_gradients;
312 
330  std::vector<std::vector<Tensor<2, spacedim> > > solution_hessians;
331  };
332 
333 }
334 
335 
414 template <int dim>
416 {
417 public:
418 
424  virtual ~DataPostprocessor ();
425 
444  virtual
445  void
447  std::vector<Vector<double> > &computed_quantities) const;
448 
469  virtual
470  void
471  compute_derived_quantities_scalar (const std::vector<double> &solution_values,
472  const std::vector<Tensor<1,dim> > &solution_gradients,
473  const std::vector<Tensor<2,dim> > &solution_hessians,
474  const std::vector<Point<dim> > &normals,
475  const std::vector<Point<dim> > &evaluation_points,
476  std::vector<Vector<double> > &computed_quantities) const DEAL_II_DEPRECATED;
477 
483  virtual
484  void
486  std::vector<Vector<double> > &computed_quantities) const;
487 
508  virtual
509  void
510  compute_derived_quantities_vector (const std::vector<Vector<double> > &solution_values,
511  const std::vector<std::vector<Tensor<1,dim> > > &solution_gradients,
512  const std::vector<std::vector<Tensor<2,dim> > > &solution_hessians,
513  const std::vector<Point<dim> > &normals,
514  const std::vector<Point<dim> > &evaluation_points,
515  std::vector<Vector<double> > &computed_quantities) const DEAL_II_DEPRECATED;
516 
521  virtual std::vector<std::string> get_names () const = 0;
522 
545  virtual
546  std::vector<DataComponentInterpretation::DataComponentInterpretation>
548 
556  virtual UpdateFlags get_needed_update_flags () const = 0;
557 };
558 
559 
560 
580 template <int dim>
582 {
583 public:
596  DataPostprocessorScalar (const std::string &name,
597  const UpdateFlags update_flags);
598 
604  virtual std::vector<std::string> get_names () const;
605 
613  virtual
614  std::vector<DataComponentInterpretation::DataComponentInterpretation>
616 
622  virtual UpdateFlags get_needed_update_flags () const;
623 
624 private:
628  const std::string name;
629  const UpdateFlags update_flags;
630 };
631 
632 
633 
655 template <int dim>
657 {
658 public:
671  DataPostprocessorVector (const std::string &name,
672  const UpdateFlags update_flags);
673 
679  virtual std::vector<std::string> get_names () const;
680 
688  virtual
689  std::vector<DataComponentInterpretation::DataComponentInterpretation>
691 
697  virtual UpdateFlags get_needed_update_flags () const;
698 
699 private:
703  const std::string name;
704  const UpdateFlags update_flags;
705 };
706 
707 
708 
709 #ifndef DOXYGEN
710 // -------------------- template functions ----------------------
711 
712 namespace DataPostprocessorInputs
713 {
714 
715  template <int spacedim>
716  template <typename DoFHandlerType>
717  void
718  CommonInputs<spacedim>::set_cell (const typename DoFHandlerType::cell_iterator &new_cell)
719  {
720  // see if we had previously already stored a cell that has the same
721  // data type; if so, reuse the memory location and avoid calling 'new'
722  // inside boost::any
723  if (typename DoFHandlerType::cell_iterator * storage_location
724  = boost::any_cast<typename DoFHandlerType::cell_iterator>(&cell))
725  *storage_location = new_cell;
726  else
727  // if we had nothing stored before, or if we had stored a different
728  // data type, just let boost::any replace things
729  cell = new_cell;
730  }
731 
732 
733 
734  template <int spacedim>
735  template <typename DoFHandlerType>
736  typename DoFHandlerType::cell_iterator
738  {
739  Assert(cell.empty() == false,
740  ExcMessage ("You are trying to access the cell associated with a "
741  "DataPostprocessorInputs::Scalar object for which no cell has "
742  "been set."));
743  Assert(boost::any_cast<typename DoFHandlerType::cell_iterator>(&cell) != 0,
744  ExcMessage ("You are trying to access the cell associated with a "
745  "DataPostprocessorInputs::Scalar with a DoFHandler type that "
746  "is different from the type with which it has been set. For "
747  "example, if the cell for which output is currently being "
748  "generated belongs to a hp::DoFHandler<2,3> object, then you can "
749  "only call the current function with a template argument "
750  "equal to hp::DoFHandler<2,3>, but not with any other class "
751  "type or dimension template argument."));
752  return boost::any_cast<typename DoFHandlerType::cell_iterator>(cell);
753  }
754 }
755 
756 #endif
757 
758 DEAL_II_NAMESPACE_CLOSE
759 
760 #endif
virtual void evaluate_vector_field(const DataPostprocessorInputs::Vector< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
std::vector< std::vector< Tensor< 2, spacedim > > > solution_hessians
virtual void evaluate_scalar_field(const DataPostprocessorInputs::Scalar< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
std::vector< double > solution_values
DoFHandlerType::cell_iterator get_cell() const
virtual std::vector< DataComponentInterpretation::DataComponentInterpretation > get_data_component_interpretation() const
DataPostprocessorScalar(const std::string &name, const UpdateFlags update_flags)
virtual UpdateFlags get_needed_update_flags() const
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:313
UpdateFlags
std::vector< std::vector< Tensor< 1, spacedim > > > solution_gradients
void set_cell(const typename DoFHandlerType::cell_iterator &cell)
std::vector< Tensor< 2, spacedim > > solution_hessians
virtual std::vector< std::string > get_names() const
virtual void compute_derived_quantities_scalar(const std::vector< double > &solution_values, const std::vector< Tensor< 1, dim > > &solution_gradients, const std::vector< Tensor< 2, dim > > &solution_hessians, const std::vector< Point< dim > > &normals, const std::vector< Point< dim > > &evaluation_points, std::vector< Vector< double > > &computed_quantities) const 1
std::vector< Tensor< 1, spacedim > > solution_gradients
virtual std::vector< std::string > get_names() const
virtual UpdateFlags get_needed_update_flags() const
std::vector< Point< spacedim > > evaluation_points
virtual void compute_derived_quantities_vector(const std::vector< Vector< double > > &solution_values, const std::vector< std::vector< Tensor< 1, dim > > > &solution_gradients, const std::vector< std::vector< Tensor< 2, dim > > > &solution_hessians, const std::vector< Point< dim > > &normals, const std::vector< Point< dim > > &evaluation_points, std::vector< Vector< double > > &computed_quantities) const 1
DataPostprocessorVector(const std::string &name, const UpdateFlags update_flags)
virtual std::vector< DataComponentInterpretation::DataComponentInterpretation > get_data_component_interpretation() const
virtual std::vector< DataComponentInterpretation::DataComponentInterpretation > get_data_component_interpretation() const
std::vector< Tensor< 1, spacedim > > normals
virtual std::vector< std::string > get_names() const =0
std::vector<::Vector< double > > solution_values
virtual UpdateFlags get_needed_update_flags() const =0