Reference documentation for deal.II version 8.5.1
integration_info.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 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 
17 #ifndef dealii__mesh_worker_integration_info_h
18 #define dealii__mesh_worker_integration_info_h
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/quadrature_lib.h>
22 #include <deal.II/base/std_cxx11/shared_ptr.h>
23 #include <deal.II/dofs/block_info.h>
24 #include <deal.II/fe/fe_values.h>
25 #include <deal.II/meshworker/local_results.h>
26 #include <deal.II/meshworker/dof_info.h>
27 #include <deal.II/meshworker/vector_selector.h>
28 
29 DEAL_II_NAMESPACE_OPEN
30 
31 namespace MeshWorker
32 {
72  template<int dim, int spacedim = dim>
74  {
75  private:
77  std::vector<std_cxx11::shared_ptr<FEValuesBase<dim, spacedim> > > fevalv;
78  public:
79  static const unsigned int dimension = dim;
80  static const unsigned int space_dimension = spacedim;
81 
86 
91 
110  template <class FEVALUES>
112  const Mapping<dim,spacedim> &mapping,
114  const UpdateFlags flags,
115  const BlockInfo *local_block_info = 0);
116 
120  void initialize_data(const std_cxx11::shared_ptr<VectorDataBase<dim,spacedim> > &data);
121 
125  void clear();
126 
132 
134  bool multigrid;
136 
141  const FEValuesBase<dim, spacedim> &fe_values () const;
142 
144 
148  const FEValuesBase<dim, spacedim> &fe_values (const unsigned int i) const;
149 
158  std::vector<std::vector<std::vector<double> > > values;
159 
168  std::vector<std::vector<std::vector<Tensor<1,dim> > > > gradients;
169 
178  std::vector<std::vector<std::vector<Tensor<2,dim> > > > hessians;
179 
183  template <typename number>
184  void reinit(const DoFInfo<dim, spacedim, number> &i);
185 
190  template<typename number>
191  void fill_local_data(const DoFInfo<dim, spacedim, number> &info, bool split_fevalues);
192 
197  std_cxx11::shared_ptr<VectorDataBase<dim, spacedim> > global_data;
198 
202  std::size_t memory_consumption () const;
203 
204  private:
209 
215  template <typename TYPE>
216  void fill_local_data(
217  std::vector<std::vector<std::vector<TYPE> > > &data,
218  VectorSelector &selector,
219  bool split_fevalues) const;
223  unsigned int n_components;
224  };
225 
280  template <int dim, int spacedim=dim>
282  {
283  public:
284 
289 
294 
303  const Mapping<dim, spacedim> &mapping,
304  const BlockInfo *block_info = 0);
305 
313  template <typename VectorType>
315  const Mapping<dim, spacedim> &mapping,
316  const AnyData &data,
317  const VectorType &dummy,
318  const BlockInfo *block_info = 0);
326  template <typename VectorType>
328  const Mapping<dim, spacedim> &mapping,
329  const AnyData &data,
330  const MGLevelObject<VectorType> &dummy,
331  const BlockInfo *block_info = 0);
335  /* @{ */
336 
347  void initialize_update_flags(bool neighbor_geometry = false);
348 
353  void add_update_flags_all (const UpdateFlags flags);
354 
358  void add_update_flags_cell(const UpdateFlags flags);
359 
363  void add_update_flags_boundary(const UpdateFlags flags);
364 
368  void add_update_flags_face(const UpdateFlags flags);
369 
376  void add_update_flags(const UpdateFlags flags,
377  const bool cell = true,
378  const bool boundary = true,
379  const bool face = true,
380  const bool neighbor = true);
381 
391  void initialize_gauss_quadrature(unsigned int n_cell_points,
392  unsigned int n_boundary_points,
393  unsigned int n_face_points,
394  const bool force = true);
395 
399  std::size_t memory_consumption () const;
400 
413 
420 
428 
433 
438 
443  /* @} */
444 
448  /* @{ */
449 
464 
470 
476 
477  std_cxx11::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim> > cell_data;
478  std_cxx11::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim> > boundary_data;
479  std_cxx11::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim> > face_data;
480  /* @} */
481 
485  /* @{ */
503  template <class DOFINFO>
504  void post_cell(const DoFInfoBox<dim, DOFINFO> &);
505 
522  template <class DOFINFO>
523  void post_faces(const DoFInfoBox<dim, DOFINFO> &);
524 
546 
547  /* @} */
548  };
549 
550 
551 //----------------------------------------------------------------------//
552 
553  template<int dim, int sdim>
554  inline
556  :
557  fevalv(0),
558  multigrid(false),
559  global_data(std_cxx11::shared_ptr<VectorDataBase<dim, sdim> >(new VectorDataBase<dim, sdim>)),
560  n_components(numbers::invalid_unsigned_int)
561  {}
562 
563 
564  template<int dim, int sdim>
565  inline
567  :
568  multigrid(other.multigrid),
569  values(other.values),
570  gradients(other.gradients),
571  hessians(other.hessians),
572  global_data(other.global_data),
573  fe_pointer(other.fe_pointer),
574  n_components(other.n_components)
575  {
576  fevalv.resize(other.fevalv.size());
577  for (unsigned int i=0; i<other.fevalv.size(); ++i)
578  {
579  const FEValuesBase<dim,sdim> &p = *other.fevalv[i];
580  const FEValues<dim,sdim> *pc = dynamic_cast<const FEValues<dim,sdim>*>(&p);
581  const FEFaceValues<dim,sdim> *pf = dynamic_cast<const FEFaceValues<dim,sdim>*>(&p);
582  const FESubfaceValues<dim,sdim> *ps = dynamic_cast<const FESubfaceValues<dim,sdim>*>(&p);
583 
584  if (pc != 0)
585  fevalv[i] = std_cxx11::shared_ptr<FEValuesBase<dim,sdim> > (
586  new FEValues<dim,sdim> (pc->get_mapping(), pc->get_fe(),
587  pc->get_quadrature(), pc->get_update_flags()));
588  else if (pf != 0)
589  fevalv[i] = std_cxx11::shared_ptr<FEValuesBase<dim,sdim> > (
590  new FEFaceValues<dim,sdim> (pf->get_mapping(), pf->get_fe(), pf->get_quadrature(), pf->get_update_flags()));
591  else if (ps != 0)
592  fevalv[i] = std_cxx11::shared_ptr<FEValuesBase<dim,sdim> > (
594  else
595  Assert(false, ExcInternalError());
596  }
597  }
598 
599 
600 
601  template<int dim, int sdim>
602  template <class FEVALUES>
603  inline void
605  const FiniteElement<dim,sdim> &el,
606  const Mapping<dim,sdim> &mapping,
608  const UpdateFlags flags,
609  const BlockInfo *block_info)
610  {
611  fe_pointer = &el;
612  if (block_info == 0 || block_info->local().size() == 0)
613  {
614  fevalv.resize(1);
615  fevalv[0] = std_cxx11::shared_ptr<FEValuesBase<dim,sdim> > (
616  new FEVALUES (mapping, el, quadrature, flags));
617  }
618  else
619  {
620  fevalv.resize(el.n_base_elements());
621  for (unsigned int i=0; i<fevalv.size(); ++i)
622  {
623  fevalv[i] = std_cxx11::shared_ptr<FEValuesBase<dim,sdim> > (
624  new FEVALUES (mapping, el.base_element(i), quadrature, flags));
625  }
626  }
627  n_components = el.n_components();
628  }
629 
630 
631  template <int dim, int spacedim>
632  inline const FiniteElement<dim, spacedim> &
634  {
635  Assert (fe_pointer !=0, ExcNotInitialized());
636  return *fe_pointer;
637  }
638 
639  template <int dim, int spacedim>
640  inline const FEValuesBase<dim, spacedim> &
642  {
643  AssertDimension(fevalv.size(), 1);
644  return *fevalv[0];
645  }
646 
647 
648  template <int dim, int spacedim>
649  inline const FEValuesBase<dim, spacedim> &
651  {
652  Assert (i<fevalv.size(), ExcIndexRange(i,0,fevalv.size()));
653  return *fevalv[i];
654  }
655 
656 
657  template <int dim, int spacedim>
658  template <typename number>
659  inline void
661  {
662  for (unsigned int i=0; i<fevalv.size(); ++i)
663  {
664  FEValuesBase<dim, spacedim> &febase = *fevalv[i];
666  {
667  // This is a subface
668  FESubfaceValues<dim> &fe = dynamic_cast<FESubfaceValues<dim>&> (febase);
669  fe.reinit(info.cell, info.face_number, info.sub_number);
670  }
672  {
673  // This is a face
674  FEFaceValues<dim> &fe = dynamic_cast<FEFaceValues<dim>&> (febase);
675  fe.reinit(info.cell, info.face_number);
676  }
677  else
678  {
679  // This is a cell
680  FEValues<dim,spacedim> &fe = dynamic_cast<FEValues<dim,spacedim>&> (febase);
681  fe.reinit(info.cell);
682  }
683  }
684 
685  const bool split_fevalues = info.block_info != 0;
686  if (!global_data->empty())
687  fill_local_data(info, split_fevalues);
688  }
689 
690 
691 
692 
693 //----------------------------------------------------------------------//
694 
695  template <int dim, int sdim>
696  inline
697  void
699  unsigned int cp,
700  unsigned int bp,
701  unsigned int fp,
702  bool force)
703  {
704  if (force || cell_quadrature.size() == 0)
705  cell_quadrature = QGauss<dim>(cp);
706  if (force || boundary_quadrature.size() == 0)
707  boundary_quadrature = QGauss<dim-1>(bp);
708  if (force || face_quadrature.size() == 0)
709  face_quadrature = QGauss<dim-1>(fp);
710  }
711 
712 
713  template <int dim, int sdim>
714  inline
715  void
717  {
718  add_update_flags(flags, true, true, true, true);
719  }
720 
721 
722  template <int dim, int sdim>
723  inline
724  void
726  {
727  add_update_flags(flags, true, false, false, false);
728  }
729 
730 
731  template <int dim, int sdim>
732  inline
733  void
735  {
736  add_update_flags(flags, false, true, false, false);
737  }
738 
739 
740  template <int dim, int sdim>
741  inline
742  void
744  {
745  add_update_flags(flags, false, false, true, true);
746  }
747 
748 
749  template <int dim, int sdim>
750  inline
751  void
753  const FiniteElement<dim,sdim> &el,
754  const Mapping<dim,sdim> &mapping,
755  const BlockInfo *block_info)
756  {
757  initialize_update_flags();
758  initialize_gauss_quadrature(
759  (cell_flags & update_values) ? (el.tensor_degree()+1) : el.tensor_degree(),
760  (boundary_flags & update_values) ? (el.tensor_degree()+1) : el.tensor_degree(),
761  (face_flags & update_values) ? (el.tensor_degree()+1) : el.tensor_degree(), false);
762 
763  cell.template initialize<FEValues<dim,sdim> >(el, mapping, cell_quadrature,
764  cell_flags, block_info);
765  boundary.template initialize<FEFaceValues<dim,sdim> >(el, mapping, boundary_quadrature,
766  boundary_flags, block_info);
767  face.template initialize<FEFaceValues<dim,sdim> >(el, mapping, face_quadrature,
768  face_flags, block_info);
769  subface.template initialize<FESubfaceValues<dim,sdim> >(el, mapping, face_quadrature,
770  face_flags, block_info);
771  neighbor.template initialize<FEFaceValues<dim,sdim> >(el, mapping, face_quadrature,
772  neighbor_flags, block_info);
773  }
774 
775 
776  template <int dim, int sdim>
777  template <typename VectorType>
778  void
781  const Mapping<dim,sdim> &mapping,
782  const AnyData &data,
783  const VectorType &,
784  const BlockInfo *block_info)
785  {
786  initialize(el, mapping, block_info);
787  std_cxx11::shared_ptr<VectorData<VectorType, dim, sdim> > p;
789 
790  p = std_cxx11::shared_ptr<VectorData<VectorType, dim, sdim> >(new VectorData<VectorType, dim, sdim> (cell_selector));
791  // Public member function of parent class was not found without
792  // explicit cast
793  pp = &*p;
794  pp->initialize(data);
795  cell_data = p;
796  cell.initialize_data(p);
797 
798  p = std_cxx11::shared_ptr<VectorData<VectorType, dim, sdim> >(new VectorData<VectorType, dim, sdim> (boundary_selector));
799  pp = &*p;
800  pp->initialize(data);
801  boundary_data = p;
802  boundary.initialize_data(p);
803 
804  p = std_cxx11::shared_ptr<VectorData<VectorType, dim, sdim> >(new VectorData<VectorType, dim, sdim> (face_selector));
805  pp = &*p;
806  pp->initialize(data);
807  face_data = p;
808  face.initialize_data(p);
809  subface.initialize_data(p);
810  neighbor.initialize_data(p);
811  }
812 
813  template <int dim, int sdim>
814  template <typename VectorType>
815  void
818  const Mapping<dim,sdim> &mapping,
819  const AnyData &data,
821  const BlockInfo *block_info)
822  {
823  initialize(el, mapping, block_info);
824  std_cxx11::shared_ptr<MGVectorData<VectorType, dim, sdim> > p;
826 
827  p = std_cxx11::shared_ptr<MGVectorData<VectorType, dim, sdim> >(new MGVectorData<VectorType, dim, sdim> (cell_selector));
828  // Public member function of parent class was not found without
829  // explicit cast
830  pp = &*p;
831  pp->initialize(data);
832  cell_data = p;
833  cell.initialize_data(p);
834 
835  p = std_cxx11::shared_ptr<MGVectorData<VectorType, dim, sdim> >(new MGVectorData<VectorType, dim, sdim> (boundary_selector));
836  pp = &*p;
837  pp->initialize(data);
838  boundary_data = p;
839  boundary.initialize_data(p);
840 
841  p = std_cxx11::shared_ptr<MGVectorData<VectorType, dim, sdim> >(new MGVectorData<VectorType, dim, sdim> (face_selector));
842  pp = &*p;
843  pp->initialize(data);
844  face_data = p;
845  face.initialize_data(p);
846  subface.initialize_data(p);
847  neighbor.initialize_data(p);
848  }
849 
850  template <int dim, int sdim>
851  template <class DOFINFO>
852  void
854  {}
855 
856 
857  template <int dim, int sdim>
858  template <class DOFINFO>
859  void
861  {}
862 
863 
864 }
865 
866 DEAL_II_NAMESPACE_CLOSE
867 
868 #endif
Shape function values.
std::vector< std::vector< std::vector< double > > > values
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
Definition: dof_info.h:72
static const unsigned int invalid_unsigned_int
Definition: types.h:170
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
void initialize(const FiniteElement< dim, spacedim > &el, const Mapping< dim, spacedim > &mapping, const BlockInfo *block_info=0)
void post_faces(const DoFInfoBox< dim, DOFINFO > &)
Quadrature< dim-1 > boundary_quadrature
void add_update_flags_cell(const UpdateFlags flags)
const FEValuesBase< dim, spacedim > & fe_values() const
Access to finite element.
std::vector< std::vector< std::vector< Tensor< 1, dim > > > > gradients
const BlockIndices & local() const
Definition: block_info.h:217
std::size_t memory_consumption() const
void add_update_flags_all(const UpdateFlags flags)
void fill_local_data(const DoFInfo< dim, spacedim, number > &info, bool split_fevalues)
const FiniteElement< dim, spacedim > & get_fe() const
UpdateFlags get_update_flags() const
unsigned int sub_number
Definition: dof_info.h:91
static ::ExceptionBase & ExcNotInitialized()
void initialize_update_flags(bool neighbor_geometry=false)
void initialize(const FiniteElement< dim, spacedim > &el, const Mapping< dim, spacedim > &mapping, const Quadrature< FEVALUES::integral_dimension > &quadrature, const UpdateFlags flags, const BlockInfo *local_block_info=0)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
void add_update_flags_face(const UpdateFlags flags)
const Mapping< dim, spacedim > & get_mapping() const
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access > > &cell, const unsigned int face_no, const unsigned int subface_no)
Definition: fe_values.cc:4093
const FiniteElement< dim, spacedim > & finite_element() const
void reinit(const DoFInfo< dim, spacedim, number > &i)
void initialize_gauss_quadrature(unsigned int n_cell_points, unsigned int n_boundary_points, unsigned int n_face_points, const bool force=true)
Quadrature< dim-1 > face_quadrature
std::size_t memory_consumption() const
MeshWorker::VectorSelector cell_selector
#define Assert(cond, exc)
Definition: exceptions.h:313
UpdateFlags
void post_cell(const DoFInfoBox< dim, DOFINFO > &)
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access > > &cell, const unsigned int face_no)
Definition: fe_values.cc:3927
Abstract base class for mapping classes.
Definition: dof_tools.h:46
const Quadrature< dim-1 > & get_quadrature() const
SmartPointer< const BlockInfo, DoFInfo< dim, spacedim > > block_info
The block structure of the system.
Definition: dof_info.h:160
void initialize_data(const std_cxx11::shared_ptr< VectorDataBase< dim, spacedim > > &data)
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
Definition: fe.cc:1251
std::vector< std_cxx11::shared_ptr< FEValuesBase< dim, spacedim > > > fevalv
vector of FEValues objects
void add_update_flags_boundary(const UpdateFlags flags)
const Quadrature< dim > & get_quadrature() const
MeshWorker::VectorSelector boundary_selector
bool multigrid
This is true if we are assembling for multigrid.
unsigned int n_components() const
void add_update_flags(const UpdateFlags flags, const bool cell=true, const bool boundary=true, const bool face=true, const bool neighbor=true)
Definition: fe.h:30
unsigned int face_number
Definition: dof_info.h:83
A small class collecting the different BlockIndices involved in global, multilevel and local computat...
Definition: block_info.h:93
IntegrationInfo< dim, spacedim > CellInfo
SmartPointer< const FiniteElement< dim, spacedim >, IntegrationInfo< dim, spacedim > > fe_pointer
std::vector< std::vector< std::vector< Tensor< 2, dim > > > > hessians
unsigned int n_base_elements() const
Definition: fe.h:2853
unsigned int size() const
void initialize(const AnyData &)
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access > > &cell)
Definition: fe_values.cc:3729
std_cxx11::shared_ptr< VectorDataBase< dim, spacedim > > global_data
MeshWorker::VectorSelector face_selector
unsigned int tensor_degree() const
static ::ExceptionBase & ExcInternalError()