16#ifndef dealii_mesh_worker_integration_info_h
17#define dealii_mesh_worker_integration_info_h
75 template <
int dim,
int spacedim = dim>
80 std::vector<std::shared_ptr<FEValuesBase<dim, spacedim>>>
fevalv;
114 template <
class FEVALUES>
120 const BlockInfo *local_block_info =
nullptr);
168 std::vector<std::vector<std::vector<double>>>
values;
178 std::vector<std::vector<std::vector<Tensor<1, spacedim>>>>
gradients;
188 std::vector<std::vector<std::vector<Tensor<2, spacedim>>>>
hessians;
193 template <
typename number>
201 template <
typename number>
204 bool split_fevalues);
231 template <
typename TYPE>
235 bool split_fevalues)
const;
295 template <
int dim,
int spacedim = dim>
328 template <
typename VectorType>
333 const VectorType &dummy,
342 template <
typename VectorType>
400 const bool cell =
true,
402 const bool face =
true,
416 unsigned int n_boundary_points,
417 unsigned int n_face_points,
418 const bool force =
true);
502 std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim>>
cell_data;
504 std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim>>
face_data;
528 template <
class DOFINFO>
548 template <
class DOFINFO>
580 template <
int dim,
int sdim>
585 , n_components(
numbers::invalid_unsigned_int)
589 template <
int dim,
int sdim>
592 : multigrid(other.multigrid)
593 , values(other.values)
594 , gradients(other.gradients)
595 , hessians(other.hessians)
596 , global_data(other.global_data)
597 , fe_pointer(other.fe_pointer)
598 , n_components(other.n_components)
601 for (
unsigned int i = 0; i < other.
fevalv.size(); ++i)
613 std::make_shared<FEValues<dim, sdim>>(pc->
get_mapping(),
617 else if (pf !=
nullptr)
619 std::make_shared<FEFaceValues<dim, sdim>>(pf->
get_mapping(),
623 else if (ps !=
nullptr)
624 fevalv[i] = std::make_shared<FESubfaceValues<dim, sdim>>(
637 template <
int dim,
int sdim>
638 template <
class FEVALUES>
648 if (block_info ==
nullptr || block_info->
local().
size() == 0)
651 fevalv[0] = std::make_shared<FEVALUES>(mapping, el, quadrature, flags);
656 for (
unsigned int i = 0; i < fevalv.size(); ++i)
657 fevalv[i] = std::make_shared<FEVALUES>(mapping,
666 template <
int dim,
int spacedim>
674 template <
int dim,
int spacedim>
683 template <
int dim,
int spacedim>
692 template <
int dim,
int spacedim>
693 template <
typename number>
698 for (
unsigned int i = 0; i < fevalv.size(); ++i)
724 const bool split_fevalues = info.
block_info !=
nullptr;
725 if (!global_data->empty())
726 fill_local_data(info, split_fevalues);
733 template <
int dim,
int sdim>
740 if (force || cell_quadrature.empty())
742 if (force || boundary_quadrature.empty())
743 boundary_quadrature =
QGauss<dim - 1>(bp);
744 if (force || face_quadrature.empty())
749 template <
int dim,
int sdim>
753 add_update_flags(flags,
true,
true,
true,
true);
757 template <
int dim,
int sdim>
761 add_update_flags(flags,
true,
false,
false,
false);
765 template <
int dim,
int sdim>
770 add_update_flags(flags,
false,
true,
false,
false);
774 template <
int dim,
int sdim>
778 add_update_flags(flags,
false,
false,
true,
true);
783 template <
int dim,
int sdim>
789 initialize_update_flags();
794 (el.tensor_degree() + 1) :
797 (el.tensor_degree() + 1) :
801 cell.template initialize<FEValues<dim, sdim>>(
802 el, mapping, cell_quadrature, cell_flags, block_info);
803 boundary.template initialize<FEFaceValues<dim, sdim>>(
804 el, mapping, boundary_quadrature, boundary_flags, block_info);
805 face.template initialize<FEFaceValues<dim, sdim>>(
806 el, mapping, face_quadrature, face_flags, block_info);
807 subface.template initialize<FESubfaceValues<dim, sdim>>(
808 el, mapping, face_quadrature, face_flags, block_info);
809 neighbor.template initialize<FEFaceValues<dim, sdim>>(
810 el, mapping, face_quadrature, neighbor_flags, block_info);
815 template <
int dim,
int sdim>
816 template <
typename VectorType>
824 initialize(el, mapping, block_info);
825 std::shared_ptr<VectorData<VectorType, dim, sdim>> p;
828 p = std::make_shared<VectorData<VectorType, dim, sdim>>(cell_selector);
834 cell.initialize_data(p);
836 p = std::make_shared<VectorData<VectorType, dim, sdim>>(boundary_selector);
840 boundary.initialize_data(p);
842 p = std::make_shared<VectorData<VectorType, dim, sdim>>(face_selector);
846 face.initialize_data(p);
847 subface.initialize_data(p);
848 neighbor.initialize_data(p);
851 template <
int dim,
int sdim>
852 template <
typename VectorType>
860 initialize(el, mapping, block_info);
861 std::shared_ptr<MGVectorData<VectorType, dim, sdim>> p;
864 p = std::make_shared<MGVectorData<VectorType, dim, sdim>>(cell_selector);
870 cell.initialize_data(p);
873 std::make_shared<MGVectorData<VectorType, dim, sdim>>(boundary_selector);
877 boundary.initialize_data(p);
879 p = std::make_shared<MGVectorData<VectorType, dim, sdim>>(face_selector);
883 face.initialize_data(p);
884 subface.initialize_data(p);
885 neighbor.initialize_data(p);
888 template <
int dim,
int sdim>
889 template <
class DOFINFO>
895 template <
int dim,
int sdim>
896 template <
class DOFINFO>
unsigned int size() const
A small class collecting the different BlockIndices involved in global, multilevel and local computat...
const BlockIndices & local() const
const Quadrature< dim - 1 > & get_quadrature() const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const unsigned int face_no)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const unsigned int face_no, const unsigned int subface_no)
UpdateFlags get_update_flags() const
const Mapping< dim, spacedim > & get_mapping() const
const FiniteElement< dim, spacedim > & get_fe() const
const Quadrature< dim > & get_quadrature() const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
unsigned int tensor_degree() const
unsigned int n_components() const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
unsigned int n_base_elements() const
Abstract base class for mapping classes.
ObserverPointer< const BlockInfo, DoFInfo< dim, spacedim > > block_info
The block structure of the system.
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
std::size_t memory_consumption() const
Quadrature< dim > cell_quadrature
void initialize_update_flags(bool neighbor_geometry=false)
void initialize_gauss_quadrature(unsigned int n_cell_points, unsigned int n_boundary_points, unsigned int n_face_points, const bool force=true)
void post_cell(const DoFInfoBox< dim, DOFINFO > &)
void initialize(const FiniteElement< dim, spacedim > &el, const Mapping< dim, spacedim > &mapping, const BlockInfo *block_info=nullptr)
void add_update_flags_all(const UpdateFlags flags)
MeshWorker::VectorSelector boundary_selector
void add_update_flags(const UpdateFlags flags, const bool cell=true, const bool boundary=true, const bool face=true, const bool neighbor=true)
std::shared_ptr< MeshWorker::VectorDataBase< dim, spacedim > > face_data
std::shared_ptr< MeshWorker::VectorDataBase< dim, spacedim > > boundary_data
UpdateFlags neighbor_flags
std::shared_ptr< MeshWorker::VectorDataBase< dim, spacedim > > cell_data
void initialize(const FiniteElement< dim, spacedim > &el, const Mapping< dim, spacedim > &mapping, const AnyData &data, const MGLevelObject< VectorType > &dummy, const BlockInfo *block_info=nullptr)
Quadrature< dim - 1 > boundary_quadrature
void initialize(const FiniteElement< dim, spacedim > &el, const Mapping< dim, spacedim > &mapping, const AnyData &data, const VectorType &dummy, const BlockInfo *block_info=nullptr)
void add_update_flags_cell(const UpdateFlags flags)
Quadrature< dim - 1 > face_quadrature
MeshWorker::VectorSelector face_selector
void post_faces(const DoFInfoBox< dim, DOFINFO > &)
MeshWorker::VectorSelector cell_selector
void add_update_flags_face(const UpdateFlags flags)
UpdateFlags boundary_flags
void add_update_flags_boundary(const UpdateFlags flags)
static constexpr unsigned int space_dimension
static constexpr unsigned int dimension
std::vector< std::vector< std::vector< Tensor< 2, spacedim > > > > hessians
std::size_t memory_consumption() const
const FEValuesBase< dim, spacedim > & fe_values(const unsigned int i) const
Access to finite elements.
const FiniteElement< dim, spacedim > & finite_element() const
IntegrationInfo(const IntegrationInfo< dim, spacedim > &other)
const FEValuesBase< dim, spacedim > & fe_values() const
Access to finite element.
void fill_local_data(std::vector< std::vector< std::vector< TYPE > > > &data, VectorSelector &selector, bool split_fevalues) const
bool multigrid
This is true if we are assembling for multigrid.
std::shared_ptr< VectorDataBase< dim, spacedim > > global_data
ObserverPointer< const FiniteElement< dim, spacedim >, IntegrationInfo< dim, spacedim > > fe_pointer
void reinit(const DoFInfo< dim, spacedim, number > &i)
void fill_local_data(const DoFInfo< dim, spacedim, number > &info, bool split_fevalues)
std::vector< std::vector< std::vector< Tensor< 1, spacedim > > > > gradients
std::vector< std::shared_ptr< FEValuesBase< dim, spacedim > > > fevalv
vector of FEValues objects
unsigned int n_components
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)
void initialize_data(const std::shared_ptr< VectorDataBase< dim, spacedim > > &data)
std::vector< std::vector< std::vector< double > > > values
void initialize(const AnyData &)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcNotInitialized()
@ update_values
Shape function values.
#define DEAL_II_ASSERT_UNREACHABLE()
std::vector< index_type > data
static const unsigned int invalid_unsigned_int