Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
integration_info.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 2020 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 
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 
23 
25 
26 #include <deal.II/fe/fe_values.h>
27 
31 
32 #include <memory>
33 
35 
36 namespace MeshWorker
37 {
77  template <int dim, int spacedim = dim>
79  {
80  private:
82  std::vector<std::shared_ptr<FEValuesBase<dim, spacedim>>> fevalv;
83 
84  public:
85  static const unsigned int dimension = dim;
86  static const unsigned int space_dimension = spacedim;
87 
92 
97 
116  template <class FEVALUES>
117  void
119  const Mapping<dim, spacedim> & mapping,
121  const UpdateFlags flags,
122  const BlockInfo *local_block_info = nullptr);
123 
127  void
128  initialize_data(const std::shared_ptr<VectorDataBase<dim, spacedim>> &data);
129 
133  void
134  clear();
135 
141  finite_element() const;
142 
144  bool multigrid;
146 
152  fe_values() const;
153 
155 
160  fe_values(const unsigned int i) const;
161 
170  std::vector<std::vector<std::vector<double>>> values;
171 
180  std::vector<std::vector<std::vector<Tensor<1, spacedim>>>> gradients;
181 
190  std::vector<std::vector<std::vector<Tensor<2, spacedim>>>> hessians;
191 
195  template <typename number>
196  void
198 
203  template <typename number>
204  void
206  bool split_fevalues);
207 
212  std::shared_ptr<VectorDataBase<dim, spacedim>> global_data;
213 
217  std::size_t
218  memory_consumption() const;
219 
220  private:
227 
233  template <typename TYPE>
234  void
235  fill_local_data(std::vector<std::vector<std::vector<TYPE>>> &data,
236  VectorSelector & selector,
237  bool split_fevalues) const;
241  unsigned int n_components;
242  };
243 
298  template <int dim, int spacedim = dim>
300  {
301  public:
306 
311 
319  void
321  const Mapping<dim, spacedim> & mapping,
322  const BlockInfo * block_info = nullptr);
323 
331  template <typename VectorType>
332  void
334  const Mapping<dim, spacedim> & mapping,
335  const AnyData & data,
336  const VectorType & dummy,
337  const BlockInfo * block_info = nullptr);
345  template <typename VectorType>
346  void
348  const Mapping<dim, spacedim> & mapping,
349  const AnyData & data,
351  const BlockInfo * block_info = nullptr);
355  /* @{ */
356 
367  void
368  initialize_update_flags(bool neighbor_geometry = false);
369 
374  void
375  add_update_flags_all(const UpdateFlags flags);
376 
380  void
381  add_update_flags_cell(const UpdateFlags flags);
382 
386  void
388 
392  void
393  add_update_flags_face(const UpdateFlags flags);
394 
401  void
402  add_update_flags(const UpdateFlags flags,
403  const bool cell = true,
404  const bool boundary = true,
405  const bool face = true,
406  const bool neighbor = true);
407 
417  void
418  initialize_gauss_quadrature(unsigned int n_cell_points,
419  unsigned int n_boundary_points,
420  unsigned int n_face_points,
421  const bool force = true);
422 
426  std::size_t
427  memory_consumption() const;
428 
441 
448 
456 
461 
466 
471  /* @} */
472 
476  /* @{ */
477 
492 
498 
504 
505  std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim>> cell_data;
506  std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim>> boundary_data;
507  std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim>> face_data;
508  /* @} */
509 
513  /* @{ */
531  template <class DOFINFO>
532  void
534 
551  template <class DOFINFO>
552  void
554 
576 
577  /* @} */
578  };
579 
580 
581  //----------------------------------------------------------------------//
582 
583  template <int dim, int sdim>
585  : fevalv(0)
586  , multigrid(false)
587  , global_data(std::make_shared<VectorDataBase<dim, sdim>>())
589  {}
590 
591 
592  template <int dim, int sdim>
594  const IntegrationInfo<dim, sdim> &other)
595  : multigrid(other.multigrid)
596  , values(other.values)
597  , gradients(other.gradients)
598  , hessians(other.hessians)
599  , global_data(other.global_data)
600  , fe_pointer(other.fe_pointer)
601  , n_components(other.n_components)
602  {
603  fevalv.resize(other.fevalv.size());
604  for (unsigned int i = 0; i < other.fevalv.size(); ++i)
605  {
606  const FEValuesBase<dim, sdim> &p = *other.fevalv[i];
607  const FEValues<dim, sdim> * pc =
608  dynamic_cast<const FEValues<dim, sdim> *>(&p);
609  const FEFaceValues<dim, sdim> *pf =
610  dynamic_cast<const FEFaceValues<dim, sdim> *>(&p);
611  const FESubfaceValues<dim, sdim> *ps =
612  dynamic_cast<const FESubfaceValues<dim, sdim> *>(&p);
613 
614  if (pc != nullptr)
615  fevalv[i] =
616  std::make_shared<FEValues<dim, sdim>>(pc->get_mapping(),
617  pc->get_fe(),
618  pc->get_quadrature(),
619  pc->get_update_flags());
620  else if (pf != nullptr)
621  fevalv[i] =
622  std::make_shared<FEFaceValues<dim, sdim>>(pf->get_mapping(),
623  pf->get_fe(),
624  pf->get_quadrature(),
625  pf->get_update_flags());
626  else if (ps != nullptr)
627  fevalv[i] = std::make_shared<FESubfaceValues<dim, sdim>>(
628  ps->get_mapping(),
629  ps->get_fe(),
630  ps->get_quadrature(),
631  ps->get_update_flags());
632  else
633  Assert(false, ExcInternalError());
634  }
635  }
636 
637 
638 
639  template <int dim, int sdim>
640  template <class FEVALUES>
641  inline void
643  const FiniteElement<dim, sdim> & el,
644  const Mapping<dim, sdim> & mapping,
646  const UpdateFlags flags,
647  const BlockInfo * block_info)
648  {
649  fe_pointer = &el;
650  if (block_info == nullptr || block_info->local().size() == 0)
651  {
652  fevalv.resize(1);
653  fevalv[0] = std::make_shared<FEVALUES>(mapping, el, quadrature, flags);
654  }
655  else
656  {
657  fevalv.resize(el.n_base_elements());
658  for (unsigned int i = 0; i < fevalv.size(); ++i)
659  fevalv[i] = std::make_shared<FEVALUES>(mapping,
660  el.base_element(i),
661  quadrature,
662  flags);
663  }
664  n_components = el.n_components();
665  }
666 
667 
668  template <int dim, int spacedim>
669  inline const FiniteElement<dim, spacedim> &
671  {
672  Assert(fe_pointer != nullptr, ExcNotInitialized());
673  return *fe_pointer;
674  }
675 
676  template <int dim, int spacedim>
677  inline const FEValuesBase<dim, spacedim> &
679  {
680  AssertDimension(fevalv.size(), 1);
681  return *fevalv[0];
682  }
683 
684 
685  template <int dim, int spacedim>
686  inline const FEValuesBase<dim, spacedim> &
688  {
689  AssertIndexRange(i, fevalv.size());
690  return *fevalv[i];
691  }
692 
693 
694  template <int dim, int spacedim>
695  template <typename number>
696  inline void
698  const DoFInfo<dim, spacedim, number> &info)
699  {
700  for (unsigned int i = 0; i < fevalv.size(); ++i)
701  {
702  FEValuesBase<dim, spacedim> &febase = *fevalv[i];
704  {
705  // This is a subface
707  dynamic_cast<FESubfaceValues<dim, spacedim> &>(febase);
708  fe.reinit(info.cell, info.face_number, info.sub_number);
709  }
711  {
712  // This is a face
714  dynamic_cast<FEFaceValues<dim, spacedim> &>(febase);
715  fe.reinit(info.cell, info.face_number);
716  }
717  else
718  {
719  // This is a cell
721  dynamic_cast<FEValues<dim, spacedim> &>(febase);
722  fe.reinit(info.cell);
723  }
724  }
725 
726  const bool split_fevalues = info.block_info != nullptr;
727  if (!global_data->empty())
728  fill_local_data(info, split_fevalues);
729  }
730 
731 
732 
733  //----------------------------------------------------------------------//
734 
735  template <int dim, int sdim>
736  inline void
738  unsigned int bp,
739  unsigned int fp,
740  bool force)
741  {
742  if (force || cell_quadrature.size() == 0)
743  cell_quadrature = QGauss<dim>(cp);
744  if (force || boundary_quadrature.size() == 0)
745  boundary_quadrature = QGauss<dim - 1>(bp);
746  if (force || face_quadrature.size() == 0)
747  face_quadrature = QGauss<dim - 1>(fp);
748  }
749 
750 
751  template <int dim, int sdim>
752  inline void
754  {
755  add_update_flags(flags, true, true, true, true);
756  }
757 
758 
759  template <int dim, int sdim>
760  inline void
762  {
763  add_update_flags(flags, true, false, false, false);
764  }
765 
766 
767  template <int dim, int sdim>
768  inline void
770  const UpdateFlags flags)
771  {
772  add_update_flags(flags, false, true, false, false);
773  }
774 
775 
776  template <int dim, int sdim>
777  inline void
779  {
780  add_update_flags(flags, false, false, true, true);
781  }
782 
783 
784  template <int dim, int sdim>
785  inline void
787  const Mapping<dim, sdim> &mapping,
788  const BlockInfo *block_info)
789  {
790  initialize_update_flags();
791  initialize_gauss_quadrature((cell_flags & update_values) ?
792  (el.tensor_degree() + 1) :
793  el.tensor_degree(),
794  (boundary_flags & update_values) ?
795  (el.tensor_degree() + 1) :
796  el.tensor_degree(),
797  (face_flags & update_values) ?
798  (el.tensor_degree() + 1) :
799  el.tensor_degree(),
800  false);
801 
802  cell.template initialize<FEValues<dim, sdim>>(
803  el, mapping, cell_quadrature, cell_flags, block_info);
804  boundary.template initialize<FEFaceValues<dim, sdim>>(
805  el, mapping, boundary_quadrature, boundary_flags, block_info);
806  face.template initialize<FEFaceValues<dim, sdim>>(
807  el, mapping, face_quadrature, face_flags, block_info);
808  subface.template initialize<FESubfaceValues<dim, sdim>>(
809  el, mapping, face_quadrature, face_flags, block_info);
810  neighbor.template initialize<FEFaceValues<dim, sdim>>(
811  el, mapping, face_quadrature, neighbor_flags, block_info);
812  }
813 
814 
815  template <int dim, int sdim>
816  template <typename VectorType>
817  void
819  const Mapping<dim, sdim> &mapping,
820  const AnyData & data,
821  const VectorType &,
822  const BlockInfo *block_info)
823  {
824  initialize(el, mapping, block_info);
825  std::shared_ptr<VectorData<VectorType, dim, sdim>> p;
827 
828  p = std::make_shared<VectorData<VectorType, dim, sdim>>(cell_selector);
829  // Public member function of parent class was not found without
830  // explicit cast
831  pp = &*p;
832  pp->initialize(data);
833  cell_data = p;
834  cell.initialize_data(p);
835 
836  p = std::make_shared<VectorData<VectorType, dim, sdim>>(boundary_selector);
837  pp = &*p;
838  pp->initialize(data);
839  boundary_data = p;
840  boundary.initialize_data(p);
841 
842  p = std::make_shared<VectorData<VectorType, dim, sdim>>(face_selector);
843  pp = &*p;
844  pp->initialize(data);
845  face_data = p;
846  face.initialize_data(p);
847  subface.initialize_data(p);
848  neighbor.initialize_data(p);
849  }
850 
851  template <int dim, int sdim>
852  template <typename VectorType>
853  void
855  const Mapping<dim, sdim> &mapping,
856  const AnyData & data,
858  const BlockInfo *block_info)
859  {
860  initialize(el, mapping, block_info);
861  std::shared_ptr<MGVectorData<VectorType, dim, sdim>> p;
863 
864  p = std::make_shared<MGVectorData<VectorType, dim, sdim>>(cell_selector);
865  // Public member function of parent class was not found without
866  // explicit cast
867  pp = &*p;
868  pp->initialize(data);
869  cell_data = p;
870  cell.initialize_data(p);
871 
872  p =
873  std::make_shared<MGVectorData<VectorType, dim, sdim>>(boundary_selector);
874  pp = &*p;
875  pp->initialize(data);
876  boundary_data = p;
877  boundary.initialize_data(p);
878 
879  p = std::make_shared<MGVectorData<VectorType, dim, sdim>>(face_selector);
880  pp = &*p;
881  pp->initialize(data);
882  face_data = p;
883  face.initialize_data(p);
884  subface.initialize_data(p);
885  neighbor.initialize_data(p);
886  }
887 
888  template <int dim, int sdim>
889  template <class DOFINFO>
890  void
892  {}
893 
894 
895  template <int dim, int sdim>
896  template <class DOFINFO>
897  void
899  {}
900 
901 
902 } // namespace MeshWorker
903 
905 
906 #endif
MeshWorker::IntegrationInfoBox::initialize
void initialize(const FiniteElement< dim, spacedim > &el, const Mapping< dim, spacedim > &mapping, const BlockInfo *block_info=nullptr)
AnyData
Definition: any_data.h:37
fe_values.h
MeshWorker::IntegrationInfoBox::face
CellInfo face
Definition: integration_info.h:566
MeshWorker::IntegrationInfoBox::cell
CellInfo cell
Definition: integration_info.h:558
MeshWorker::IntegrationInfo::fill_local_data
void fill_local_data(const DoFInfo< dim, spacedim, number > &info, bool split_fevalues)
MeshWorker::IntegrationInfo::IntegrationInfo
IntegrationInfo()
Definition: integration_info.h:584
MeshWorker::IntegrationInfoBox::subface
CellInfo subface
Definition: integration_info.h:571
BlockIndices::size
unsigned int size() const
Definition: block_indices.h:356
vector_selector.h
MeshWorker::IntegrationInfoBox::initialize_gauss_quadrature
void initialize_gauss_quadrature(unsigned int n_cell_points, unsigned int n_boundary_points, unsigned int n_face_points, const bool force=true)
Definition: integration_info.h:737
MeshWorker::IntegrationInfoBox::post_cell
void post_cell(const DoFInfoBox< dim, DOFINFO > &)
Definition: integration_info.h:891
FiniteElementData::tensor_degree
unsigned int tensor_degree() const
MeshWorker::IntegrationInfoBox::add_update_flags_boundary
void add_update_flags_boundary(const UpdateFlags flags)
Definition: integration_info.h:769
MeshWorker::IntegrationInfoBox::post_faces
void post_faces(const DoFInfoBox< dim, DOFINFO > &)
Definition: integration_info.h:898
MeshWorker
Definition: assemble_flags.h:30
MeshWorker::IntegrationInfo::values
std::vector< std::vector< std::vector< double > > > values
Definition: integration_info.h:170
MeshWorker::IntegrationInfo::clear
void clear()
MeshWorker::VectorSelector
Definition: vector_selector.h:52
MeshWorker::IntegrationInfo::fe_values
const FEValuesBase< dim, spacedim > & fe_values() const
Access to finite element.
Definition: integration_info.h:678
MeshWorker::IntegrationInfoBox::add_update_flags
void add_update_flags(const UpdateFlags flags, const bool cell=true, const bool boundary=true, const bool face=true, const bool neighbor=true)
MeshWorker::IntegrationInfo::global_data
std::shared_ptr< VectorDataBase< dim, spacedim > > global_data
Definition: integration_info.h:212
local_results.h
VectorType
MeshWorker::DoFInfo::block_info
SmartPointer< const BlockInfo, DoFInfo< dim, spacedim > > block_info
The block structure of the system.
Definition: dof_info.h:173
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
FiniteElement::n_base_elements
unsigned int n_base_elements() const
Definition: fe.h:3096
MeshWorker::IntegrationInfoBox::neighbor
CellInfo neighbor
Definition: integration_info.h:575
DoFHandler::n_components
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
MeshWorker::IntegrationInfo::fevalv
std::vector< std::shared_ptr< FEValuesBase< dim, spacedim > > > fevalv
vector of FEValues objects
Definition: integration_info.h:82
FESubfaceValues::reinit
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access >> &cell, const unsigned int face_no, const unsigned int subface_no)
MeshWorker::DoFInfo::sub_number
unsigned int sub_number
Definition: dof_info.h:99
update_values
@ update_values
Shape function values.
Definition: fe_update_flags.h:78
DoFTools::n_components
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
MeshWorker::IntegrationInfo
Definition: integration_info.h:78
MeshWorker::IntegrationInfo::finite_element
const FiniteElement< dim, spacedim > & finite_element() const
Definition: integration_info.h:670
quadrature_lib.h
MeshWorker::IntegrationInfoBox::initialize_update_flags
void initialize_update_flags(bool neighbor_geometry=false)
FEValues
Definition: fe.h:38
BlockInfo::local
const BlockIndices & local() const
Definition: block_info.h:252
MeshWorker::IntegrationInfo::gradients
std::vector< std::vector< std::vector< Tensor< 1, spacedim > > > > gradients
Definition: integration_info.h:180
MeshWorker::IntegrationInfoBox::face_selector
MeshWorker::VectorSelector face_selector
Definition: integration_info.h:503
MeshWorker::IntegrationInfo::reinit
void reinit(const DoFInfo< dim, spacedim, number > &i)
Definition: integration_info.h:697
Mapping
Abstract base class for mapping classes.
Definition: mapping.h:302
MeshWorker::IntegrationInfo::initialize_data
void initialize_data(const std::shared_ptr< VectorDataBase< dim, spacedim >> &data)
MeshWorker::IntegrationInfo::memory_consumption
std::size_t memory_consumption() const
MeshWorker::IntegrationInfoBox::cell_quadrature
Quadrature< dim > cell_quadrature
Definition: integration_info.h:460
internal::reinit
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:621
FiniteElement::base_element
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
Definition: fe.cc:1252
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
MeshWorker::IntegrationInfoBox::add_update_flags_face
void add_update_flags_face(const UpdateFlags flags)
Definition: integration_info.h:778
FEValues::reinit
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access >> &cell)
MeshWorker::IntegrationInfoBox::face_quadrature
Quadrature< dim - 1 > face_quadrature
Definition: integration_info.h:470
MeshWorker::VectorDataBase< dim, spacedim >
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
StandardExceptions::ExcNotInitialized
static ::ExceptionBase & ExcNotInitialized()
MeshWorker::IntegrationInfoBox::boundary_quadrature
Quadrature< dim - 1 > boundary_quadrature
Definition: integration_info.h:465
numbers
Definition: numbers.h:207
dof_info.h
MeshWorker::IntegrationInfoBox::face_flags
UpdateFlags face_flags
Definition: integration_info.h:447
FiniteElement
Definition: fe.h:648
FEValuesBase
Definition: fe.h:36
UpdateFlags
UpdateFlags
Definition: fe_update_flags.h:66
FEFaceValues::reinit
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access >> &cell, const unsigned int face_no)
QGauss
Definition: quadrature_lib.h:40
SmartPointer
Definition: smartpointer.h:68
MeshWorker::IntegrationInfoBox::boundary_data
std::shared_ptr< MeshWorker::VectorDataBase< dim, spacedim > > boundary_data
Definition: integration_info.h:506
internal::dummy
const types::global_dof_index * dummy()
Definition: dof_handler.cc:59
MeshWorker::IntegrationInfo::dimension
static const unsigned int dimension
Definition: integration_info.h:85
FESubfaceValues
Definition: fe.h:42
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
MeshWorker::IntegrationInfo::initialize
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=nullptr)
MeshWorker::IntegrationInfo::multigrid
bool multigrid
This is true if we are assembling for multigrid.
Definition: integration_info.h:144
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
MeshWorker::IntegrationInfoBox::IntegrationInfoBox
IntegrationInfoBox()
MeshWorker::IntegrationInfo::fe_pointer
SmartPointer< const FiniteElement< dim, spacedim >, IntegrationInfo< dim, spacedim > > fe_pointer
Definition: integration_info.h:226
MGLevelObject< VectorType >
MeshWorker::IntegrationInfoBox
Definition: integration_info.h:299
MeshWorker::IntegrationInfoBox::add_update_flags_all
void add_update_flags_all(const UpdateFlags flags)
Definition: integration_info.h:753
block_info.h
MeshWorker::DoFInfo::face_number
unsigned int face_number
Definition: dof_info.h:91
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
MeshWorker::IntegrationInfo::n_components
unsigned int n_components
Definition: integration_info.h:241
MeshWorker::IntegrationInfoBox::cell_data
std::shared_ptr< MeshWorker::VectorDataBase< dim, spacedim > > cell_data
Definition: integration_info.h:505
MeshWorker::IntegrationInfoBox::boundary
CellInfo boundary
Definition: integration_info.h:562
FiniteElementData::n_components
unsigned int n_components() const
MeshWorker::IntegrationInfoBox::cell_selector
MeshWorker::VectorSelector cell_selector
Definition: integration_info.h:491
config.h
MeshWorker::IntegrationInfoBox::memory_consumption
std::size_t memory_consumption() const
BlockInfo
A small class collecting the different BlockIndices involved in global, multilevel and local computat...
Definition: block_info.h:99
FEFaceValues
Definition: fe.h:40
MeshWorker::IntegrationInfoBox::boundary_selector
MeshWorker::VectorSelector boundary_selector
Definition: integration_info.h:497
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
MeshWorker::IntegrationInfoBox::cell_flags
UpdateFlags cell_flags
Definition: integration_info.h:434
Quadrature
Definition: quadrature.h:85
MeshWorker::IntegrationInfo::hessians
std::vector< std::vector< std::vector< Tensor< 2, spacedim > > > > hessians
Definition: integration_info.h:190
MeshWorker::IntegrationInfoBox::boundary_flags
UpdateFlags boundary_flags
Definition: integration_info.h:440
MeshWorker::DoFInfo
Definition: dof_info.h:76
FEValues::get_quadrature
const Quadrature< dim > & get_quadrature() const
MeshWorker::IntegrationInfo::space_dimension
static const unsigned int space_dimension
Definition: integration_info.h:86
MeshWorker::DoFInfo::cell
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
Definition: dof_info.h:80
MeshWorker::DoFInfoBox
Definition: dof_info.h:223
MeshWorker::IntegrationInfoBox::add_update_flags_cell
void add_update_flags_cell(const UpdateFlags flags)
Definition: integration_info.h:761
MeshWorker::IntegrationInfoBox::neighbor_flags
UpdateFlags neighbor_flags
Definition: integration_info.h:455
MeshWorker::IntegrationInfoBox::face_data
std::shared_ptr< MeshWorker::VectorDataBase< dim, spacedim > > face_data
Definition: integration_info.h:507