Reference documentation for deal.II version 9.4.1
\(\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\}}\)
Loading...
Searching...
No Matches
integration_info.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2006 - 2022 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
27
31
32#include <memory>
33
35
36namespace MeshWorker
37{
76 template <int dim, int spacedim = dim>
78 {
79 private:
81 std::vector<std::shared_ptr<FEValuesBase<dim, spacedim>>> fevalv;
82
83 public:
84 static constexpr unsigned int dimension = dim;
85 static constexpr unsigned int space_dimension = spacedim;
86
91
96
115 template <class FEVALUES>
116 void
118 const Mapping<dim, spacedim> & mapping,
120 const UpdateFlags flags,
121 const BlockInfo *local_block_info = nullptr);
122
126 void
127 initialize_data(const std::shared_ptr<VectorDataBase<dim, spacedim>> &data);
128
132 void
134
141
145
151 fe_values() const;
152
154
159 fe_values(const unsigned int i) const;
160
169 std::vector<std::vector<std::vector<double>>> values;
170
179 std::vector<std::vector<std::vector<Tensor<1, spacedim>>>> gradients;
180
189 std::vector<std::vector<std::vector<Tensor<2, spacedim>>>> hessians;
190
194 template <typename number>
195 void
197
202 template <typename number>
203 void
205 bool split_fevalues);
206
211 std::shared_ptr<VectorDataBase<dim, spacedim>> global_data;
212
216 std::size_t
218
219 private:
226
232 template <typename TYPE>
233 void
234 fill_local_data(std::vector<std::vector<std::vector<TYPE>>> &data,
235 VectorSelector & selector,
236 bool split_fevalues) const;
240 unsigned int n_components;
241 };
242
296 template <int dim, int spacedim = dim>
298 {
299 public:
304
309
317 void
319 const Mapping<dim, spacedim> & mapping,
320 const BlockInfo * block_info = nullptr);
321
329 template <typename VectorType>
330 void
332 const Mapping<dim, spacedim> & mapping,
333 const AnyData & data,
334 const VectorType & dummy,
335 const BlockInfo * block_info = nullptr);
343 template <typename VectorType>
344 void
346 const Mapping<dim, spacedim> & mapping,
347 const AnyData & data,
348 const MGLevelObject<VectorType> & dummy,
349 const BlockInfo * block_info = nullptr);
353 /* @{ */
354
365 void
366 initialize_update_flags(bool neighbor_geometry = false);
367
372 void
374
378 void
380
384 void
386
390 void
392
399 void
401 const bool cell = true,
402 const bool boundary = true,
403 const bool face = true,
404 const bool neighbor = true);
405
415 void
416 initialize_gauss_quadrature(unsigned int n_cell_points,
417 unsigned int n_boundary_points,
418 unsigned int n_face_points,
419 const bool force = true);
420
424 std::size_t
426
439
446
454
459
464
469 /* @} */
470
474 /* @{ */
475
490
496
502
503 std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim>> cell_data;
504 std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim>> boundary_data;
505 std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim>> face_data;
506 /* @} */
507
511 /* @{ */
529 template <class DOFINFO>
530 void
532
549 template <class DOFINFO>
550 void
552
574
575 /* @} */
576 };
577
578
579 //----------------------------------------------------------------------//
580
581 template <int dim, int sdim>
583 : fevalv(0)
584 , multigrid(false)
585 , global_data(std::make_shared<VectorDataBase<dim, sdim>>())
586 , n_components(numbers::invalid_unsigned_int)
587 {}
588
589
590 template <int dim, int sdim>
592 const IntegrationInfo<dim, sdim> &other)
593 : multigrid(other.multigrid)
594 , values(other.values)
595 , gradients(other.gradients)
596 , hessians(other.hessians)
597 , global_data(other.global_data)
598 , fe_pointer(other.fe_pointer)
599 , n_components(other.n_components)
600 {
601 fevalv.resize(other.fevalv.size());
602 for (unsigned int i = 0; i < other.fevalv.size(); ++i)
603 {
604 const FEValuesBase<dim, sdim> &p = *other.fevalv[i];
605 const FEValues<dim, sdim> * pc =
606 dynamic_cast<const FEValues<dim, sdim> *>(&p);
607 const FEFaceValues<dim, sdim> *pf =
608 dynamic_cast<const FEFaceValues<dim, sdim> *>(&p);
610 dynamic_cast<const FESubfaceValues<dim, sdim> *>(&p);
611
612 if (pc != nullptr)
613 fevalv[i] =
614 std::make_shared<FEValues<dim, sdim>>(pc->get_mapping(),
615 pc->get_fe(),
616 pc->get_quadrature(),
617 pc->get_update_flags());
618 else if (pf != nullptr)
619 fevalv[i] =
620 std::make_shared<FEFaceValues<dim, sdim>>(pf->get_mapping(),
621 pf->get_fe(),
622 pf->get_quadrature(),
623 pf->get_update_flags());
624 else if (ps != nullptr)
625 fevalv[i] = std::make_shared<FESubfaceValues<dim, sdim>>(
626 ps->get_mapping(),
627 ps->get_fe(),
628 ps->get_quadrature(),
629 ps->get_update_flags());
630 else
631 Assert(false, ExcInternalError());
632 }
633 }
634
635
636
637 template <int dim, int sdim>
638 template <class FEVALUES>
639 inline void
641 const FiniteElement<dim, sdim> & el,
642 const Mapping<dim, sdim> & mapping,
644 const UpdateFlags flags,
645 const BlockInfo * block_info)
646 {
647 fe_pointer = &el;
648 if (block_info == nullptr || block_info->local().size() == 0)
649 {
650 fevalv.resize(1);
651 fevalv[0] = std::make_shared<FEVALUES>(mapping, el, quadrature, flags);
652 }
653 else
654 {
655 fevalv.resize(el.n_base_elements());
656 for (unsigned int i = 0; i < fevalv.size(); ++i)
657 fevalv[i] = std::make_shared<FEVALUES>(mapping,
658 el.base_element(i),
659 quadrature,
660 flags);
661 }
662 n_components = el.n_components();
663 }
664
665
666 template <int dim, int spacedim>
667 inline const FiniteElement<dim, spacedim> &
669 {
670 Assert(fe_pointer != nullptr, ExcNotInitialized());
671 return *fe_pointer;
672 }
673
674 template <int dim, int spacedim>
675 inline const FEValuesBase<dim, spacedim> &
677 {
678 AssertDimension(fevalv.size(), 1);
679 return *fevalv[0];
680 }
681
682
683 template <int dim, int spacedim>
684 inline const FEValuesBase<dim, spacedim> &
686 {
687 AssertIndexRange(i, fevalv.size());
688 return *fevalv[i];
689 }
690
691
692 template <int dim, int spacedim>
693 template <typename number>
694 inline void
697 {
698 for (unsigned int i = 0; i < fevalv.size(); ++i)
699 {
700 FEValuesBase<dim, spacedim> &febase = *fevalv[i];
702 {
703 // This is a subface
705 dynamic_cast<FESubfaceValues<dim, spacedim> &>(febase);
706 fe.reinit(info.cell, info.face_number, info.sub_number);
707 }
709 {
710 // This is a face
712 dynamic_cast<FEFaceValues<dim, spacedim> &>(febase);
713 fe.reinit(info.cell, info.face_number);
714 }
715 else
716 {
717 // This is a cell
719 dynamic_cast<FEValues<dim, spacedim> &>(febase);
720 fe.reinit(info.cell);
721 }
722 }
723
724 const bool split_fevalues = info.block_info != nullptr;
725 if (!global_data->empty())
726 fill_local_data(info, split_fevalues);
727 }
728
729
730
731 //----------------------------------------------------------------------//
732
733 template <int dim, int sdim>
734 inline void
736 unsigned int bp,
737 unsigned int fp,
738 bool force)
739 {
740 if (force || cell_quadrature.size() == 0)
741 cell_quadrature = QGauss<dim>(cp);
742 if (force || boundary_quadrature.size() == 0)
743 boundary_quadrature = QGauss<dim - 1>(bp);
744 if (force || face_quadrature.size() == 0)
745 face_quadrature = QGauss<dim - 1>(fp);
746 }
747
748
749 template <int dim, int sdim>
750 inline void
752 {
753 add_update_flags(flags, true, true, true, true);
754 }
755
756
757 template <int dim, int sdim>
758 inline void
760 {
761 add_update_flags(flags, true, false, false, false);
762 }
763
764
765 template <int dim, int sdim>
766 inline void
768 const UpdateFlags flags)
769 {
770 add_update_flags(flags, false, true, false, false);
771 }
772
773
774 template <int dim, int sdim>
775 inline void
777 {
778 add_update_flags(flags, false, false, true, true);
779 }
780
781
782 template <int dim, int sdim>
783 inline void
785 const Mapping<dim, sdim> &mapping,
786 const BlockInfo *block_info)
787 {
788 initialize_update_flags();
789 initialize_gauss_quadrature((cell_flags & update_values) ?
790 (el.tensor_degree() + 1) :
791 el.tensor_degree(),
792 (boundary_flags & update_values) ?
793 (el.tensor_degree() + 1) :
794 el.tensor_degree(),
795 (face_flags & update_values) ?
796 (el.tensor_degree() + 1) :
797 el.tensor_degree(),
798 false);
799
800 cell.template initialize<FEValues<dim, sdim>>(
801 el, mapping, cell_quadrature, cell_flags, block_info);
802 boundary.template initialize<FEFaceValues<dim, sdim>>(
803 el, mapping, boundary_quadrature, boundary_flags, block_info);
804 face.template initialize<FEFaceValues<dim, sdim>>(
805 el, mapping, face_quadrature, face_flags, block_info);
806 subface.template initialize<FESubfaceValues<dim, sdim>>(
807 el, mapping, face_quadrature, face_flags, block_info);
808 neighbor.template initialize<FEFaceValues<dim, sdim>>(
809 el, mapping, face_quadrature, neighbor_flags, block_info);
810 }
811
812
813 template <int dim, int sdim>
814 template <typename VectorType>
815 void
817 const Mapping<dim, sdim> &mapping,
818 const AnyData & data,
819 const VectorType &,
820 const BlockInfo *block_info)
821 {
822 initialize(el, mapping, block_info);
823 std::shared_ptr<VectorData<VectorType, dim, sdim>> p;
825
826 p = std::make_shared<VectorData<VectorType, dim, sdim>>(cell_selector);
827 // Public member function of parent class was not found without
828 // explicit cast
829 pp = &*p;
830 pp->initialize(data);
831 cell_data = p;
832 cell.initialize_data(p);
833
834 p = std::make_shared<VectorData<VectorType, dim, sdim>>(boundary_selector);
835 pp = &*p;
836 pp->initialize(data);
837 boundary_data = p;
838 boundary.initialize_data(p);
839
840 p = std::make_shared<VectorData<VectorType, dim, sdim>>(face_selector);
841 pp = &*p;
842 pp->initialize(data);
843 face_data = p;
844 face.initialize_data(p);
845 subface.initialize_data(p);
846 neighbor.initialize_data(p);
847 }
848
849 template <int dim, int sdim>
850 template <typename VectorType>
851 void
853 const Mapping<dim, sdim> &mapping,
854 const AnyData & data,
856 const BlockInfo *block_info)
857 {
858 initialize(el, mapping, block_info);
859 std::shared_ptr<MGVectorData<VectorType, dim, sdim>> p;
861
862 p = std::make_shared<MGVectorData<VectorType, dim, sdim>>(cell_selector);
863 // Public member function of parent class was not found without
864 // explicit cast
865 pp = &*p;
866 pp->initialize(data);
867 cell_data = p;
868 cell.initialize_data(p);
869
870 p =
871 std::make_shared<MGVectorData<VectorType, dim, sdim>>(boundary_selector);
872 pp = &*p;
873 pp->initialize(data);
874 boundary_data = p;
875 boundary.initialize_data(p);
876
877 p = std::make_shared<MGVectorData<VectorType, dim, sdim>>(face_selector);
878 pp = &*p;
879 pp->initialize(data);
880 face_data = p;
881 face.initialize_data(p);
882 subface.initialize_data(p);
883 neighbor.initialize_data(p);
884 }
885
886 template <int dim, int sdim>
887 template <class DOFINFO>
888 void
890 {}
891
892
893 template <int dim, int sdim>
894 template <class DOFINFO>
895 void
897 {}
898
899
900} // namespace MeshWorker
901
903
904#endif
unsigned int size() const
A small class collecting the different BlockIndices involved in global, multilevel and local computat...
Definition: block_info.h:94
const BlockIndices & local() const
Definition: block_info.h:227
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.
Definition: mapping.h:311
unsigned int sub_number
Definition: dof_info.h:98
SmartPointer< const BlockInfo, DoFInfo< dim, spacedim > > block_info
The block structure of the system.
Definition: dof_info.h:172
unsigned int face_number
Definition: dof_info.h:90
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
Definition: dof_info.h:79
std::size_t memory_consumption() const
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
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)
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
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
SmartPointer< const FiniteElement< dim, spacedim >, IntegrationInfo< dim, spacedim > > fe_pointer
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
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
UpdateFlags
@ update_values
Shape function values.
#define Assert(cond, exc)
Definition: exceptions.h:1473
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1667
#define AssertIndexRange(index, range)
Definition: exceptions.h:1732
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcNotInitialized()
static const unsigned int invalid_unsigned_int
Definition: types.h:201
STL namespace.