Reference documentation for deal.II version GIT relicensing-224-gc660c0d696 2024-03-28 18:40:02+00:00
\(\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// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2009 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
16#ifndef dealii_mesh_worker_integration_info_h
17#define dealii_mesh_worker_integration_info_h
18
19#include <deal.II/base/config.h>
20
22
24
26
30
31#include <memory>
32
34
35namespace MeshWorker
36{
75 template <int dim, int spacedim = dim>
77 {
78 private:
80 std::vector<std::shared_ptr<FEValuesBase<dim, spacedim>>> fevalv;
81
82 public:
83 static constexpr unsigned int dimension = dim;
84 static constexpr unsigned int space_dimension = spacedim;
85
90
95
114 template <class FEVALUES>
115 void
117 const Mapping<dim, spacedim> &mapping,
119 const UpdateFlags flags,
120 const BlockInfo *local_block_info = nullptr);
121
125 void
126 initialize_data(const std::shared_ptr<VectorDataBase<dim, spacedim>> &data);
127
131 void
133
140
144
150 fe_values() const;
151
153
158 fe_values(const unsigned int i) const;
159
168 std::vector<std::vector<std::vector<double>>> values;
169
178 std::vector<std::vector<std::vector<Tensor<1, spacedim>>>> gradients;
179
188 std::vector<std::vector<std::vector<Tensor<2, spacedim>>>> hessians;
189
193 template <typename number>
194 void
196
201 template <typename number>
202 void
204 bool split_fevalues);
205
210 std::shared_ptr<VectorDataBase<dim, spacedim>> global_data;
211
215 std::size_t
217
218 private:
225
231 template <typename TYPE>
232 void
233 fill_local_data(std::vector<std::vector<std::vector<TYPE>>> &data,
234 VectorSelector &selector,
235 bool split_fevalues) const;
239 unsigned int n_components;
240 };
241
295 template <int dim, int spacedim = dim>
297 {
298 public:
303
308
316 void
318 const Mapping<dim, spacedim> &mapping,
319 const BlockInfo *block_info = nullptr);
320
328 template <typename VectorType>
329 void
331 const Mapping<dim, spacedim> &mapping,
332 const AnyData &data,
333 const VectorType &dummy,
334 const BlockInfo *block_info = nullptr);
342 template <typename VectorType>
343 void
345 const Mapping<dim, spacedim> &mapping,
346 const AnyData &data,
347 const MGLevelObject<VectorType> &dummy,
348 const BlockInfo *block_info = nullptr);
364 void
365 initialize_update_flags(bool neighbor_geometry = false);
366
371 void
373
377 void
379
383 void
385
389 void
391
398 void
400 const bool cell = true,
401 const bool boundary = true,
402 const bool face = true,
403 const bool neighbor = true);
404
414 void
415 initialize_gauss_quadrature(unsigned int n_cell_points,
416 unsigned int n_boundary_points,
417 unsigned int n_face_points,
418 const bool force = true);
419
423 std::size_t
425
438
445
453
458
463
489
495
501
502 std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim>> cell_data;
503 std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim>> boundary_data;
504 std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim>> face_data;
528 template <class DOFINFO>
529 void
531
548 template <class DOFINFO>
549 void
551
573
575 };
576
577
578 //----------------------------------------------------------------------//
579
580 template <int dim, int sdim>
582 : fevalv(0)
583 , multigrid(false)
584 , global_data(std::make_shared<VectorDataBase<dim, sdim>>())
585 , n_components(numbers::invalid_unsigned_int)
586 {}
587
588#ifndef DOXYGEN
589 template <int dim, int sdim>
591 const IntegrationInfo<dim, sdim> &other)
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)
599 {
600 fevalv.resize(other.fevalv.size());
601 for (unsigned int i = 0; i < other.fevalv.size(); ++i)
602 {
603 const FEValuesBase<dim, sdim> &p = *other.fevalv[i];
604 const FEValues<dim, sdim> *pc =
605 dynamic_cast<const FEValues<dim, sdim> *>(&p);
606 const FEFaceValues<dim, sdim> *pf =
607 dynamic_cast<const FEFaceValues<dim, sdim> *>(&p);
609 dynamic_cast<const FESubfaceValues<dim, sdim> *>(&p);
610
611 if (pc != nullptr)
612 fevalv[i] =
613 std::make_shared<FEValues<dim, sdim>>(pc->get_mapping(),
614 pc->get_fe(),
615 pc->get_quadrature(),
616 pc->get_update_flags());
617 else if (pf != nullptr)
618 fevalv[i] =
619 std::make_shared<FEFaceValues<dim, sdim>>(pf->get_mapping(),
620 pf->get_fe(),
621 pf->get_quadrature(),
622 pf->get_update_flags());
623 else if (ps != nullptr)
624 fevalv[i] = std::make_shared<FESubfaceValues<dim, sdim>>(
625 ps->get_mapping(),
626 ps->get_fe(),
627 ps->get_quadrature(),
628 ps->get_update_flags());
629 else
631 }
632 }
633#endif
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.empty())
741 cell_quadrature = QGauss<dim>(cp);
742 if (force || boundary_quadrature.empty())
743 boundary_quadrature = QGauss<dim - 1>(bp);
744 if (force || face_quadrature.empty())
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#ifndef DOXYGEN
783 template <int dim, int sdim>
784 inline void
786 const Mapping<dim, sdim> &mapping,
787 const BlockInfo *block_info)
788 {
789 initialize_update_flags();
790 initialize_gauss_quadrature((cell_flags & update_values) ?
791 (el.tensor_degree() + 1) :
792 el.tensor_degree(),
793 (boundary_flags & update_values) ?
794 (el.tensor_degree() + 1) :
795 el.tensor_degree(),
796 (face_flags & update_values) ?
797 (el.tensor_degree() + 1) :
798 el.tensor_degree(),
799 false);
800
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);
811 }
812#endif
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
893
894
895 template <int dim, int sdim>
896 template <class DOFINFO>
897 void
900
901
902} // namespace MeshWorker
903
905
906#endif
unsigned int size() const
A small class collecting the different BlockIndices involved in global, multilevel and local computat...
Definition block_info.h:95
const BlockIndices & local() const
Definition block_info.h:228
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:316
unsigned int sub_number
Definition dof_info.h:97
SmartPointer< const BlockInfo, DoFInfo< dim, spacedim > > block_info
The block structure of the system.
Definition dof_info.h:171
unsigned int face_number
Definition dof_info.h:89
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
Definition dof_info.h:78
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:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcNotInitialized()
UpdateFlags
@ update_values
Shape function values.
#define DEAL_II_ASSERT_UNREACHABLE()
static const unsigned int invalid_unsigned_int
Definition types.h:220
STL namespace.