deal.II version GIT relicensing-2287-g6548a49e0a 2024-12-20 18:30:00+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) 2004 - 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#ifndef dealii_hp_fe_values_h
16#define dealii_hp_fe_values_h
17
18#include <deal.II/base/config.h>
19
21
22#include <deal.II/fe/fe.h>
24
25#include <deal.II/grid/tria.h>
27
31
32#include <memory>
33
35
36// Forward declaration
37#ifndef DOXYGEN
38template <int dim, int spacedim>
39class FiniteElement;
40#endif
41
42
43namespace hp
44{
66 template <int dim, int q_dim, typename FEValuesType>
68 {
69 public:
80
90 const std::vector<QCollection<q_dim>> &q_collection,
92
102
110 const std::vector<QCollection<q_dim>> &q_collection,
112
117
123 operator=(const FEValuesBase &) = delete;
124
134 void
135 precalculate_fe_values(const std::vector<unsigned int> &fe_indices,
136 const std::vector<unsigned int> &mapping_indices,
137 const std::vector<unsigned int> &q_indices);
138
151 void
153
160
166
170 const QCollection<q_dim> &
172
178
185 const FEValuesType &
187
188 protected:
196 FEValuesType &
197 select_fe_values(const unsigned int fe_index,
198 const unsigned int mapping_index,
199 const unsigned int q_index);
200
201 protected:
205 const ObserverPointer<
209
213 const ObserverPointer<
217
222
231 const std::vector<QCollection<q_dim>> q_collections;
232
233 private:
246
252
257 };
258
259} // namespace hp
260
261
262namespace hp
263{
313 template <int dim, int spacedim = dim>
315 : public hp::FEValuesBase<dim, dim, ::FEValues<dim, spacedim>>
316 {
317 public:
318 static constexpr unsigned int dimension = dim;
319
320 static constexpr unsigned int space_dimension = spacedim;
321
329
330
339
340
386 template <bool lda>
387 void
389 const unsigned int q_index = numbers::invalid_unsigned_int,
390 const unsigned int mapping_index = numbers::invalid_unsigned_int,
391 const unsigned int fe_index = numbers::invalid_unsigned_int);
392
405 void
407 const unsigned int q_index = numbers::invalid_unsigned_int,
408 const unsigned int mapping_index = numbers::invalid_unsigned_int,
409 const unsigned int fe_index = numbers::invalid_unsigned_int);
410 };
411
412
413
437 template <int dim, int spacedim = dim>
439 : public hp::FEValuesBase<dim, dim - 1, ::FEFaceValues<dim, spacedim>>
440 {
441 public:
449
461 const std::vector<hp::QCollection<dim - 1>> &q_collections,
463
464
473
484 const std::vector<hp::QCollection<dim - 1>> &q_collections,
486
535 template <bool lda>
536 void
538 const unsigned int face_no,
539 const unsigned int q_index = numbers::invalid_unsigned_int,
540 const unsigned int mapping_index = numbers::invalid_unsigned_int,
541 const unsigned int fe_index = numbers::invalid_unsigned_int);
542
548 template <bool lda>
549 void
552 const unsigned int q_index = numbers::invalid_unsigned_int,
553 const unsigned int mapping_index = numbers::invalid_unsigned_int,
554 const unsigned int fe_index = numbers::invalid_unsigned_int);
555
568 void
570 const unsigned int face_no,
571 const unsigned int q_index = numbers::invalid_unsigned_int,
572 const unsigned int mapping_index = numbers::invalid_unsigned_int,
573 const unsigned int fe_index = numbers::invalid_unsigned_int);
574
580 void
583 const unsigned int q_index = numbers::invalid_unsigned_int,
584 const unsigned int mapping_index = numbers::invalid_unsigned_int,
585 const unsigned int fe_index = numbers::invalid_unsigned_int);
586 };
587
588
589
596 template <int dim, int spacedim = dim>
598 : public hp::
599 FEValuesBase<dim, dim - 1, ::FESubfaceValues<dim, spacedim>>
600 {
601 public:
606 const hp::MappingCollection<dim, spacedim> &mapping_collection,
607 const hp::FECollection<dim, spacedim> &fe_collection,
608 const hp::QCollection<dim - 1> &q_collection,
609 const UpdateFlags update_flags);
610
611
618 const hp::QCollection<dim - 1> &q_collection,
619 const UpdateFlags update_flags);
620
656 template <bool lda>
657 void
659 const unsigned int face_no,
660 const unsigned int subface_no,
661 const unsigned int q_index = numbers::invalid_unsigned_int,
662 const unsigned int mapping_index = numbers::invalid_unsigned_int,
663 const unsigned int fe_index = numbers::invalid_unsigned_int);
664
677 void
679 const unsigned int face_no,
680 const unsigned int subface_no,
681 const unsigned int q_index = numbers::invalid_unsigned_int,
682 const unsigned int mapping_index = numbers::invalid_unsigned_int,
683 const unsigned int fe_index = numbers::invalid_unsigned_int);
684 };
685
686} // namespace hp
687
688
689// -------------- inline and template functions --------------
690
691namespace hp
692{
693 template <int dim, int q_dim, typename FEValuesType>
694 inline const FEValuesType &
696 {
697 Assert(
698 present_fe_values_index != TableIndices<3>(numbers::invalid_unsigned_int,
702 "The indices of the present FEValues object are all invalid. The reason"
703 " is likely that you have forgotten to call the reinit() function."));
704
705 return *fe_values_table(present_fe_values_index);
706 }
707
708
709
710 template <int dim, int q_dim, typename FEValuesType>
713 {
714 return *fe_collection;
715 }
716
717
718
719 template <int dim, int q_dim, typename FEValuesType>
722 {
723 return *mapping_collection;
724 }
725
726
727
728 template <int dim, int q_dim, typename FEValuesType>
729 inline const QCollection<q_dim> &
731 {
732 return q_collection;
733 }
734
735
736
737 template <int dim, int q_dim, typename FEValuesType>
738 inline UpdateFlags
740 {
741 return update_flags;
742 }
743} // namespace hp
744
746
747#endif
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, lda > > &cell, const unsigned int face_no, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
Definition fe_values.cc:423
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, lda > > &cell, const unsigned int face_no, const unsigned int subface_no, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
Definition fe_values.cc:563
const QCollection< q_dim > & get_quadrature_collection() const
Definition fe_values.h:730
const ObserverPointer< const FECollection< dim, FEValuesType::space_dimension >, FEValuesBase< dim, q_dim, FEValuesType > > fe_collection
Definition fe_values.h:208
FEValuesBase & operator=(const FEValuesBase &)=delete
const FEValuesType & get_present_fe_values() const
Definition fe_values.h:695
const std::vector< QCollection< q_dim > > q_collections
Definition fe_values.h:231
FEValuesType & select_fe_values(const unsigned int fe_index, const unsigned int mapping_index, const unsigned int q_index)
Definition fe_values.cc:175
const QCollection< q_dim > q_collection
Definition fe_values.h:221
UpdateFlags get_update_flags() const
Definition fe_values.h:739
const MappingCollection< dim, FEValuesType::space_dimension > & get_mapping_collection() const
Definition fe_values.h:721
const ObserverPointer< const MappingCollection< dim, FEValuesType::space_dimension >, FEValuesBase< dim, q_dim, FEValuesType > > mapping_collection
Definition fe_values.h:216
Table< 3, std::unique_ptr< FEValuesType > > fe_values_table
Definition fe_values.h:245
const FECollection< dim, FEValuesType::space_dimension > & get_fe_collection() const
Definition fe_values.h:712
TableIndices< 3 > present_fe_values_index
Definition fe_values.h:251
void precalculate_fe_values()
Definition fe_values.cc:242
const UpdateFlags update_flags
Definition fe_values.h:256
static constexpr unsigned int dimension
Definition fe_values.h:318
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, lda > > &cell, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
Definition fe_values.cc:293
static constexpr unsigned int space_dimension
Definition fe_values.h:320
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
#define Assert(cond, exc)
static ::ExceptionBase & ExcMessage(std::string arg1)
UpdateFlags
Definition hp.h:117
static const unsigned int invalid_unsigned_int
Definition types.h:220