Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30: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
fe_values.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2021 - 2023 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#ifndef dealii_non_matching_fe_values
16#define dealii_non_matching_fe_values
17
20
22
26
28
32
36
37#include <deque>
38#include <optional>
39
41
42namespace NonMatching
43{
79
80
142 template <int dim>
144 {
145 public:
147
169 template <typename VectorType>
171 const Quadrature<1> &quadrature,
174 const DoFHandler<dim> &dof_handler,
175 const VectorType &level_set,
176 const AdditionalData &additional_data = AdditionalData());
177
204 template <typename VectorType>
207 const hp::QCollection<dim> &q_collection,
208 const hp::QCollection<1> &q_collection_1d,
211 const DoFHandler<dim> &dof_handler,
212 const VectorType &level_set,
213 const AdditionalData &additional_data = AdditionalData());
214
223 template <bool level_dof_access>
224 void
225 reinit(
227
236 const std::optional<::FEValues<dim>> &
237 get_inside_fe_values() const;
238
247 const std::optional<::FEValues<dim>> &
248 get_outside_fe_values() const;
249
257 const std::optional<FEImmersedSurfaceValues<dim>> &
258 get_surface_fe_values() const;
259
260 private:
266 void
267 initialize(const hp::QCollection<dim> &q_collection);
268
273
278
284
289
293 unsigned int active_fe_index;
294
299
304
318 std::deque<std::optional<::FEValues<dim>>>
320
334 std::deque<std::optional<::FEValues<dim>>>
336
344 std::optional<::FEValues<dim>> fe_values_inside;
345
353 std::optional<::FEValues<dim>> fe_values_outside;
354
362 std::optional<NonMatching::FEImmersedSurfaceValues<dim>> fe_values_surface;
363
368 };
369
370
371
434 template <int dim>
436 {
437 public:
440
462 template <typename VectorType>
464 const Quadrature<1> &quadrature,
467 const DoFHandler<dim> &dof_handler,
468 const VectorType &level_set,
469 const AdditionalData &additional_data = AdditionalData());
470
497 template <typename VectorType>
500 const hp::QCollection<dim - 1> &q_collection,
501 const hp::QCollection<1> &q_collection_1d,
504 const DoFHandler<dim> &dof_handler,
505 const VectorType &level_set,
506 const AdditionalData &additional_data = AdditionalData());
507
521 template <typename CellIteratorType, typename CellNeighborIteratorType>
522 void
523 reinit(const CellIteratorType &cell,
524 const unsigned int face_no,
525 const unsigned int sub_face_no,
526 const CellNeighborIteratorType &cell_neighbor,
527 const unsigned int face_no_neighbor,
528 const unsigned int sub_face_no_neighbor);
529
530
538 template <typename CellIteratorType>
539 void
540 reinit(const CellIteratorType &cell, const unsigned int face_no);
541
551 const std::optional<::FEInterfaceValues<dim>> &
552 get_inside_fe_values() const;
553
563 const std::optional<::FEInterfaceValues<dim>> &
564 get_outside_fe_values() const;
565
566 private:
572 void
573 initialize(const hp::QCollection<dim - 1> &q_collection);
574
581 template <bool level_dof_access>
582 void
583 do_reinit(
585 const unsigned int face_no,
586 const std::function<void(::FEInterfaceValues<dim> &)> &call_reinit);
587
592
597
603
608
612 unsigned int active_fe_index;
613
618
623
639 std::deque<std::optional<::FEInterfaceValues<dim>>>
641
656 std::deque<std::optional<::FEInterfaceValues<dim>>>
658
666 std::optional<::FEInterfaceValues<dim>> fe_values_inside;
667
675 std::optional<::FEInterfaceValues<dim>> fe_values_outside;
676
681 };
682
683
684#ifndef DOXYGEN
685
686 /*---------------------- Inline functions ---------------------*/
687
688 template <int dim>
689 template <typename CellIteratorType>
690 inline void
691 FEInterfaceValues<dim>::reinit(const CellIteratorType &cell,
692 const unsigned int face_no)
693 {
694 // Lambda describing how we should call reinit on a single
695 // ::FEInterfaceValues object.
696 const auto reinit_operation =
697 [&cell, face_no](::FEInterfaceValues<dim> &fe_interface_values) {
698 fe_interface_values.reinit(cell, face_no);
699 };
700
701 do_reinit(cell, face_no, reinit_operation);
702 }
703
704
705
706 template <int dim>
707 template <typename CellIteratorType, typename CellNeighborIteratorType>
708 inline void
709 FEInterfaceValues<dim>::reinit(const CellIteratorType &cell,
710 const unsigned int face_no,
711 const unsigned int sub_face_no,
712 const CellNeighborIteratorType &cell_neighbor,
713 const unsigned int face_no_neighbor,
714 const unsigned int sub_face_no_neighbor)
715 {
717 Assert(sub_face_no_neighbor == numbers::invalid_unsigned_int,
719
720 // Lambda describing how we should call reinit on a single
721 // ::FEInterfaceValues object.
722 const auto reinit_operation =
723 [&cell,
724 face_no,
725 sub_face_no,
726 &cell_neighbor,
727 face_no_neighbor,
728 sub_face_no_neighbor](
729 ::FEInterfaceValues<dim> &fe_interface_values) {
730 fe_interface_values.reinit(cell,
731 face_no,
732 sub_face_no,
733 cell_neighbor,
734 face_no_neighbor,
735 sub_face_no_neighbor);
736 };
737
738 do_reinit(cell, face_no, reinit_operation);
739 }
740
741#endif
742
743} // namespace NonMatching
745
746#endif
const std::optional<::FEInterfaceValues< dim > > & get_outside_fe_values() const
Definition fe_values.cc:469
const SmartPointer< const hp::FECollection< dim > > fe_collection
Definition fe_values.h:596
void initialize(const hp::QCollection< dim - 1 > &q_collection)
Definition fe_values.cc:328
const SmartPointer< const hp::MappingCollection< dim > > mapping_collection
Definition fe_values.h:591
const hp::QCollection< 1 > q_collection_1D
Definition fe_values.h:602
typename FaceQuadratureGenerator< dim >::AdditionalData AdditionalData
Definition fe_values.h:439
const RegionUpdateFlags region_update_flags
Definition fe_values.h:617
LocationToLevelSet current_face_location
Definition fe_values.h:607
std::deque< std::optional<::FEInterfaceValues< dim > > > fe_values_inside_full_quadrature
Definition fe_values.h:640
const std::optional<::FEInterfaceValues< dim > > & get_inside_fe_values() const
Definition fe_values.cc:457
const SmartPointer< const MeshClassifier< dim > > mesh_classifier
Definition fe_values.h:622
void reinit(const CellIteratorType &cell, const unsigned int face_no, const unsigned int sub_face_no, const CellNeighborIteratorType &cell_neighbor, const unsigned int face_no_neighbor, const unsigned int sub_face_no_neighbor)
std::deque< std::optional<::FEInterfaceValues< dim > > > fe_values_outside_full_quadrature
Definition fe_values.h:657
std::optional<::FEInterfaceValues< dim > > fe_values_inside
Definition fe_values.h:666
void reinit(const CellIteratorType &cell, const unsigned int face_no)
DiscreteFaceQuadratureGenerator< dim > face_quadrature_generator
Definition fe_values.h:680
void do_reinit(const TriaIterator< DoFCellAccessor< dim, dim, level_dof_access > > &cell, const unsigned int face_no, const std::function< void(::FEInterfaceValues< dim > &)> &call_reinit)
Definition fe_values.cc:380
std::optional<::FEInterfaceValues< dim > > fe_values_outside
Definition fe_values.h:675
std::optional< NonMatching::FEImmersedSurfaceValues< dim > > fe_values_surface
Definition fe_values.h:362
std::deque< std::optional<::FEValues< dim > > > fe_values_inside_full_quadrature
Definition fe_values.h:319
std::optional<::FEValues< dim > > fe_values_outside
Definition fe_values.h:353
LocationToLevelSet current_cell_location
Definition fe_values.h:288
const SmartPointer< const hp::FECollection< dim > > fe_collection
Definition fe_values.h:277
const std::optional<::FEValues< dim > > & get_outside_fe_values() const
Definition fe_values.cc:249
void initialize(const hp::QCollection< dim > &q_collection)
Definition fe_values.cc:101
const std::optional< FEImmersedSurfaceValues< dim > > & get_surface_fe_values() const
Definition fe_values.cc:261
const SmartPointer< const MeshClassifier< dim > > mesh_classifier
Definition fe_values.h:303
const hp::QCollection< 1 > q_collection_1D
Definition fe_values.h:283
const RegionUpdateFlags region_update_flags
Definition fe_values.h:298
void reinit(const TriaIterator< DoFCellAccessor< dim, dim, level_dof_access > > &cell)
Definition fe_values.cc:150
std::optional<::FEValues< dim > > fe_values_inside
Definition fe_values.h:344
const SmartPointer< const hp::MappingCollection< dim > > mapping_collection
Definition fe_values.h:272
typename QuadratureGenerator< dim >::AdditionalData AdditionalData
Definition fe_values.h:146
DiscreteQuadratureGenerator< dim > quadrature_generator
Definition fe_values.h:367
std::deque< std::optional<::FEValues< dim > > > fe_values_outside_full_quadrature
Definition fe_values.h:335
unsigned int active_fe_index
Definition fe_values.h:293
const std::optional<::FEValues< dim > > & get_inside_fe_values() const
Definition fe_values.cc:237
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
UpdateFlags
static const unsigned int invalid_unsigned_int
Definition types.h:220