Reference documentation for deal.II version GIT relicensing-487-ge9eb5ab491 2024-04-25 07:20: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_data.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2001 - 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_fe_data_h
16#define dealii_fe_data_h
17
18#include <deal.II/base/config.h>
19
21
23
25
26#include <vector>
27
29
30// Forward declarations:
31#ifndef DOXYGEN
32template <int dim>
34#endif
35
136
137namespace internal
138{
183 {
187 std::vector<std::vector<unsigned int>> dofs_per_object_exclusive;
188
192 std::vector<std::vector<unsigned int>> dofs_per_object_inclusive;
193
197 std::vector<std::vector<unsigned int>> object_index;
198
202 std::vector<std::vector<unsigned int>> first_object_index_on_face;
203
207 template <int dim>
210 };
211} // namespace internal
212
227template <int dim>
229{
230public:
262 {
266 unknown = 0x00,
267
271 L2 = 0x01,
272
277 Hcurl = 0x02,
278
283 Hdiv = 0x04,
284
289
294 H2 = 0x0e
295 };
296
301 static constexpr unsigned int dimension = dim;
302
303private:
308
315
320 const unsigned int number_unique_faces;
321
322public:
326 const unsigned int dofs_per_vertex;
327
332 const unsigned int dofs_per_line;
333
334private:
339 const std::vector<unsigned int> n_dofs_on_quad;
340
341public:
346 const unsigned int dofs_per_quad;
347
348private:
352 const unsigned int dofs_per_quad_max;
353
354public:
359 const unsigned int dofs_per_hex;
360
364 const unsigned int first_line_index;
365
366private:
372 const std::vector<unsigned int> first_index_of_quads;
373
374public:
378 const unsigned int first_quad_index;
379
383 const unsigned int first_hex_index;
384
385private:
389 const std::vector<unsigned int> first_line_index_of_faces;
390
391public:
395 const unsigned int first_face_line_index;
396
397private:
401 const std::vector<unsigned int> first_quad_index_of_faces;
402
403public:
407 const unsigned int first_face_quad_index;
408
409private:
414 const std::vector<unsigned int> n_dofs_on_face;
415
416public:
422 const unsigned int dofs_per_face;
423
424private:
428 const unsigned int dofs_per_face_max;
429
430public:
436 const unsigned int dofs_per_cell;
437
446 const unsigned int components;
447
452 const unsigned int degree;
453
458
465
504 FiniteElementData(const std::vector<unsigned int> &dofs_per_object,
505 const unsigned int n_components,
506 const unsigned int degree,
507 const Conformity conformity = unknown,
509
514 FiniteElementData(const std::vector<unsigned int> &dofs_per_object,
516 const unsigned int n_components,
517 const unsigned int degree,
518 const Conformity conformity = unknown,
520
529 const unsigned int n_components,
530 const unsigned int degree,
531 const Conformity conformity = unknown,
533
541
554 unsigned int
556
563 unsigned int
565
569 unsigned int
571
575 unsigned int
577
581 unsigned int
582 n_dofs_per_quad(unsigned int face_no = 0) const;
583
588 unsigned int
590
594 unsigned int
596
601 unsigned int
602 n_dofs_per_face(unsigned int face_no = 0, unsigned int child = 0) const;
603
608 unsigned int
610
615 unsigned int
617
626 template <int structdim>
627 unsigned int
628 n_dofs_per_object(const unsigned int i = 0) const;
629
635 unsigned int
637
643 unsigned int
644 n_blocks() const;
645
649 const BlockIndices &
651
658 unsigned int
660
667 bool
668 conforms(const Conformity) const;
669
673 bool
674 operator==(const FiniteElementData &) const;
675
679 unsigned int
681
685 unsigned int
686 get_first_quad_index(const unsigned int quad_no = 0) const;
687
691 unsigned int
693
697 unsigned int
698 get_first_face_line_index(const unsigned int face_no = 0) const;
699
703 unsigned int
704 get_first_face_quad_index(const unsigned int face_no = 0) const;
705};
706
707namespace internal
708{
714 expand(const unsigned int dim,
715 const std::vector<unsigned int> &dofs_per_object,
716 const ReferenceCell reference_cell);
717} // namespace internal
718
719
720
721// --------- inline and template functions ---------------
722
723
724#ifndef DOXYGEN
725
727{
728 inline Domination
729 operator&(const Domination d1, const Domination d2)
730 {
731 // go through the entire list of possibilities. note that if we were into
732 // speed, obfuscation and cared enough, we could implement this operator
733 // by doing a bitwise & (and) if we gave these values to the enum values:
734 // neither_element_dominates=0, this_element_dominates=1,
735 // other_element_dominates=2, either_element_can_dominate=3
736 // =this_element_dominates|other_element_dominates
737 switch (d1)
738 {
740 if ((d2 == this_element_dominates) ||
743 else
745
747 if ((d2 == other_element_dominates) ||
750 else
752
755
757 if (d2 == no_requirements)
759 else
760 return d2;
761
762 case no_requirements:
763 return d2;
764
765 default:
766 // shouldn't get here
768 }
769
771 }
772} // namespace FiniteElementDomination
773
774
775template <int dim>
776inline ReferenceCell
778{
779 return reference_cell_kind;
780}
781
782
783
784template <int dim>
785inline unsigned int
787{
788 return number_of_unique_2d_subobjects;
789}
790
791
792
793template <int dim>
794inline unsigned int
796{
797 return number_unique_faces;
798}
799
800
801
802template <int dim>
803inline unsigned int
805{
806 return dofs_per_vertex;
807}
808
809
810
811template <int dim>
812inline unsigned int
814{
815 return dofs_per_line;
816}
817
818
819
820template <int dim>
821inline unsigned int
822FiniteElementData<dim>::n_dofs_per_quad(unsigned int face_no) const
823{
824 return n_dofs_on_quad[n_dofs_on_quad.size() == 1 ? 0 : face_no];
825}
826
827
828
829template <int dim>
830inline unsigned int
832{
833 return dofs_per_quad_max;
834}
835
836
837
838template <int dim>
839inline unsigned int
841{
842 return dofs_per_hex;
843}
844
845
846
847template <int dim>
848inline unsigned int
849FiniteElementData<dim>::n_dofs_per_face(unsigned int face_no,
850 unsigned int child_no) const
851{
852 (void)child_no;
853
854 return n_dofs_on_face[n_dofs_on_face.size() == 1 ? 0 : face_no];
855}
856
857
858
859template <int dim>
860inline unsigned int
862{
863 return dofs_per_face_max;
864}
865
866
867
868template <int dim>
869inline unsigned int
871{
872 return dofs_per_cell;
873}
874
875
876
877template <int dim>
878template <int structdim>
879inline unsigned int
880FiniteElementData<dim>::n_dofs_per_object(const unsigned int i) const
881{
882 switch (structdim)
883 {
884 case 0:
885 return n_dofs_per_vertex();
886 case 1:
887 return n_dofs_per_line();
888 case 2:
889 return n_dofs_per_quad((structdim == 2 && dim == 3) ? i : 0);
890 case 3:
891 return n_dofs_per_hex();
892 default:
894 }
896}
897
898
899
900template <int dim>
901inline unsigned int
903{
904 return components;
905}
906
907
908
909template <int dim>
910inline const BlockIndices &
912{
913 return block_indices_data;
914}
915
916
917
918template <int dim>
919inline unsigned int
921{
922 return block_indices_data.size();
923}
924
925
926
927template <int dim>
928inline unsigned int
930{
931 return degree;
932}
933
934
935template <int dim>
936inline bool
937FiniteElementData<dim>::conforms(const Conformity space) const
938{
939 return ((space & conforming_space) == space);
940}
941
942
943
944template <int dim>
945unsigned int
947{
948 return first_line_index;
949}
950
951template <int dim>
952unsigned int
953FiniteElementData<dim>::get_first_quad_index(const unsigned int quad_no) const
954{
955 if (first_index_of_quads.size() == 1)
956 return first_index_of_quads[0] + quad_no * n_dofs_per_quad(0);
957 else
958 return first_index_of_quads[quad_no];
959}
960
961template <int dim>
962unsigned int
964{
965 return first_hex_index;
966}
967
968template <int dim>
969unsigned int
971 const unsigned int face_no) const
972{
973 return first_line_index_of_faces[first_line_index_of_faces.size() == 1 ?
974 0 :
975 face_no];
976}
977
978template <int dim>
979unsigned int
981 const unsigned int face_no) const
982{
983 return first_quad_index_of_faces[first_quad_index_of_faces.size() == 1 ?
984 0 :
985 face_no];
986}
987
988template <int dim>
991{
992 const auto reference_cell = fe.reference_cell();
993
995
996 result.dofs_per_object_exclusive.resize(4);
997 result.dofs_per_object_inclusive.resize(4);
998 result.object_index.resize(4);
999
1000 unsigned int counter = 0;
1001
1002 for (const unsigned int v : reference_cell.vertex_indices())
1003 {
1004 const auto c = fe.template n_dofs_per_object<0>(v);
1005
1006 result.dofs_per_object_exclusive[0].emplace_back(c);
1007 result.dofs_per_object_inclusive[0].emplace_back(c);
1008 result.object_index[0].emplace_back(counter);
1009
1010 counter += c;
1011 }
1012
1013 if (dim >= 2)
1014 for (const unsigned int l : reference_cell.line_indices())
1015 {
1016 const auto c = fe.template n_dofs_per_object<1>(l);
1017
1018 result.dofs_per_object_exclusive[1].emplace_back(c);
1019 result.dofs_per_object_inclusive[1].emplace_back(
1020 c + 2 * fe.template n_dofs_per_object<0>());
1021 result.object_index[1].emplace_back(counter);
1022
1023 counter += c;
1024 }
1025
1026 if (dim == 3)
1027 for (const unsigned int f : reference_cell.face_indices())
1028 {
1029 const auto c = fe.template n_dofs_per_object<2>(f);
1030
1031 result.dofs_per_object_exclusive[2].emplace_back(c);
1032 result.dofs_per_object_inclusive[2].emplace_back(fe.n_dofs_per_face(f));
1033 result.object_index[2].emplace_back(counter);
1034
1035 counter += c;
1036 }
1037
1038 {
1039 const auto c = fe.template n_dofs_per_object<dim>();
1040
1041 result.dofs_per_object_exclusive[dim].emplace_back(c);
1042 result.dofs_per_object_inclusive[dim].emplace_back(fe.n_dofs_per_cell());
1043 result.object_index[dim].emplace_back(counter);
1044
1045 counter += c;
1046 }
1047
1048 for (unsigned int d = dim + 1; d <= 3; ++d)
1049 {
1050 result.dofs_per_object_exclusive[d].emplace_back(0);
1051 result.dofs_per_object_inclusive[d].emplace_back(0);
1052 result.object_index[d].emplace_back(counter);
1053 }
1054
1055 result.first_object_index_on_face.resize(3);
1056 for (const unsigned int face_no : reference_cell.face_indices())
1057 {
1058 result.first_object_index_on_face[0].emplace_back(0);
1059
1060 result.first_object_index_on_face[1].emplace_back(
1061 fe.get_first_face_line_index(face_no));
1062
1063 result.first_object_index_on_face[2].emplace_back(
1064 fe.get_first_face_quad_index(face_no));
1065 }
1066
1067 return result;
1068}
1069
1070
1071#endif // DOXYGEN
1072
1073
1075
1076#endif
unsigned int size() const
const unsigned int dofs_per_quad
Definition fe_data.h:346
unsigned int get_first_line_index() const
unsigned int n_dofs_per_vertex() const
const unsigned int first_face_quad_index
Definition fe_data.h:407
const std::vector< unsigned int > first_index_of_quads
Definition fe_data.h:372
bool conforms(const Conformity) const
const std::vector< unsigned int > first_line_index_of_faces
Definition fe_data.h:389
const unsigned int components
Definition fe_data.h:446
const unsigned int degree
Definition fe_data.h:452
const std::vector< unsigned int > n_dofs_on_face
Definition fe_data.h:414
unsigned int n_dofs_per_cell() const
const unsigned int dofs_per_face
Definition fe_data.h:422
unsigned int n_dofs_per_line() const
unsigned int get_first_quad_index(const unsigned int quad_no=0) const
bool operator==(const FiniteElementData &) const
Definition fe_data.cc:195
const unsigned int dofs_per_line
Definition fe_data.h:332
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int n_blocks() const
unsigned int max_dofs_per_face() const
unsigned int tensor_degree() const
unsigned int n_components() const
ReferenceCell reference_cell() const
const unsigned int first_hex_index
Definition fe_data.h:383
unsigned int n_dofs_per_object(const unsigned int i=0) const
const unsigned int first_face_line_index
Definition fe_data.h:395
const unsigned int dofs_per_quad_max
Definition fe_data.h:352
unsigned int n_unique_2d_subobjects() const
unsigned int get_first_face_line_index(const unsigned int face_no=0) const
const unsigned int first_line_index
Definition fe_data.h:364
const unsigned int dofs_per_hex
Definition fe_data.h:359
const BlockIndices block_indices_data
Definition fe_data.h:464
unsigned int get_first_face_quad_index(const unsigned int face_no=0) const
const unsigned int number_unique_faces
Definition fe_data.h:320
unsigned int max_dofs_per_quad() const
const BlockIndices & block_indices() const
const std::vector< unsigned int > first_quad_index_of_faces
Definition fe_data.h:401
const ReferenceCell reference_cell_kind
Definition fe_data.h:307
const unsigned int number_of_unique_2d_subobjects
Definition fe_data.h:314
const std::vector< unsigned int > n_dofs_on_quad
Definition fe_data.h:339
const Conformity conforming_space
Definition fe_data.h:457
const unsigned int first_quad_index
Definition fe_data.h:378
unsigned int n_unique_faces() const
unsigned int n_dofs_per_quad(unsigned int face_no=0) const
const unsigned int dofs_per_vertex
Definition fe_data.h:326
const unsigned int dofs_per_cell
Definition fe_data.h:436
static constexpr unsigned int dimension
Definition fe_data.h:301
unsigned int n_dofs_per_hex() const
const unsigned int dofs_per_face_max
Definition fe_data.h:428
unsigned int get_first_hex_index() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
unsigned int vertex_indices[2]
#define DEAL_II_ASSERT_UNREACHABLE()
Domination operator&(const Domination d1, const Domination d2)
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
internal::GenericDoFsPerObject expand(const unsigned int dim, const std::vector< unsigned int > &dofs_per_object, const ReferenceCell reference_cell)
Definition fe_data.cc:24
static const unsigned int invalid_unsigned_int
Definition types.h:220
std::vector< std::vector< unsigned int > > object_index
Definition fe_data.h:197
std::vector< std::vector< unsigned int > > first_object_index_on_face
Definition fe_data.h:202
std::vector< std::vector< unsigned int > > dofs_per_object_inclusive
Definition fe_data.h:192
std::vector< std::vector< unsigned int > > dofs_per_object_exclusive
Definition fe_data.h:187
static GenericDoFsPerObject generate(const FiniteElementData< dim > &fe)