Reference documentation for deal.II version 9.3.3
\(\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\}}\)
fe_base.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2000 - 2021 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#ifndef dealii_fe_base_h
17#define dealii_fe_base_h
18
19#include <deal.II/base/config.h>
20
22
24
26
27#include <vector>
28
30
36{
86 {
107 };
108
109
128 inline Domination operator&(const Domination d1, const Domination d2);
129} // namespace FiniteElementDomination
130
131namespace internal
132{
174 {
178 std::vector<std::vector<unsigned int>> dofs_per_object_exclusive;
179
183 std::vector<std::vector<unsigned int>> dofs_per_object_inclusive;
184
188 std::vector<std::vector<unsigned int>> object_index;
189
193 std::vector<std::vector<unsigned int>> first_object_index_on_face;
194 };
195} // namespace internal
196
211template <int dim>
213{
214public:
246 {
250 unknown = 0x00,
251
255 L2 = 0x01,
256
261 Hcurl = 0x02,
262
267 Hdiv = 0x04,
268
273
278 H2 = 0x0e
279 };
280
285 static const unsigned int dimension = dim;
286
287private:
292
297 const unsigned int number_unique_quads;
298
303 const unsigned int number_unique_faces;
304
305public:
309 const unsigned int dofs_per_vertex;
310
315 const unsigned int dofs_per_line;
316
317private:
322 const std::vector<unsigned int> n_dofs_on_quad;
323
324public:
329 const unsigned int dofs_per_quad;
330
331private:
335 const unsigned int dofs_per_quad_max;
336
337public:
342 const unsigned int dofs_per_hex;
343
347 const unsigned int first_line_index;
348
349private:
355 const std::vector<unsigned int> first_index_of_quads;
356
357public:
361 const unsigned int first_quad_index;
362
366 const unsigned int first_hex_index;
367
368private:
372 const std::vector<unsigned int> first_line_index_of_faces;
373
374public:
378 const unsigned int first_face_line_index;
379
380private:
384 const std::vector<unsigned int> first_quad_index_of_faces;
385
386public:
390 const unsigned int first_face_quad_index;
391
392private:
397 const std::vector<unsigned int> n_dofs_on_face;
398
399public:
405 const unsigned int dofs_per_face;
406
407private:
411 const unsigned int dofs_per_face_max;
412
413public:
419 const unsigned int dofs_per_cell;
420
429 const unsigned int components;
430
435 const unsigned int degree;
436
441
448
487 FiniteElementData(const std::vector<unsigned int> &dofs_per_object,
488 const unsigned int n_components,
489 const unsigned int degree,
490 const Conformity conformity = unknown,
492
497 FiniteElementData(const std::vector<unsigned int> &dofs_per_object,
499 const unsigned int n_components,
500 const unsigned int degree,
501 const Conformity conformity = unknown,
503
512 const unsigned int n_components,
513 const unsigned int degree,
514 const Conformity conformity = unknown,
516
524
529 unsigned int
531
536 unsigned int
538
542 unsigned int
544
548 unsigned int
550
554 unsigned int
555 n_dofs_per_quad(unsigned int face_no = 0) const;
556
561 unsigned int
563
567 unsigned int
569
574 unsigned int
575 n_dofs_per_face(unsigned int face_no = 0, unsigned int child = 0) const;
576
581 unsigned int
583
588 unsigned int
590
599 template <int structdim>
600 unsigned int
601 n_dofs_per_object(const unsigned int i = 0) const;
602
608 unsigned int
610
616 unsigned int
617 n_blocks() const;
618
622 const BlockIndices &
624
631 unsigned int
633
640 bool
641 conforms(const Conformity) const;
642
646 bool
647 operator==(const FiniteElementData &) const;
648
652 unsigned int
654
658 unsigned int
659 get_first_quad_index(const unsigned int quad_no = 0) const;
660
664 unsigned int
666
670 unsigned int
671 get_first_face_line_index(const unsigned int face_no = 0) const;
672
676 unsigned int
677 get_first_face_quad_index(const unsigned int face_no = 0) const;
678};
679
680namespace internal
681{
687 expand(const unsigned int dim,
688 const std::vector<unsigned int> &dofs_per_object,
689 const ::ReferenceCell reference_cell);
690} // namespace internal
691
692
693
694// --------- inline and template functions ---------------
695
696
697#ifndef DOXYGEN
698
700{
701 inline Domination operator&(const Domination d1, const Domination d2)
702 {
703 // go through the entire list of possibilities. note that if we were into
704 // speed, obfuscation and cared enough, we could implement this operator
705 // by doing a bitwise & (and) if we gave these values to the enum values:
706 // neither_element_dominates=0, this_element_dominates=1,
707 // other_element_dominates=2, either_element_can_dominate=3
708 // =this_element_dominates|other_element_dominates
709 switch (d1)
710 {
712 if ((d2 == this_element_dominates) ||
715 else
717
719 if ((d2 == other_element_dominates) ||
722 else
724
727
729 if (d2 == no_requirements)
731 else
732 return d2;
733
734 case no_requirements:
735 return d2;
736
737 default:
738 // shouldn't get here
739 Assert(false, ExcInternalError());
740 }
741
743 }
744} // namespace FiniteElementDomination
745
746
747template <int dim>
748inline ReferenceCell
750{
751 return reference_cell_kind;
752}
753
754
755
756template <int dim>
757inline unsigned int
759{
760 return number_unique_quads;
761}
762
763
764
765template <int dim>
766inline unsigned int
768{
769 return number_unique_faces;
770}
771
772
773
774template <int dim>
775inline unsigned int
777{
778 return dofs_per_vertex;
779}
780
781
782
783template <int dim>
784inline unsigned int
786{
787 return dofs_per_line;
788}
789
790
791
792template <int dim>
793inline unsigned int
794FiniteElementData<dim>::n_dofs_per_quad(unsigned int face_no) const
795{
796 return n_dofs_on_quad[n_dofs_on_quad.size() == 1 ? 0 : face_no];
797}
798
799
800
801template <int dim>
802inline unsigned int
804{
805 return dofs_per_quad_max;
806}
807
808
809
810template <int dim>
811inline unsigned int
813{
814 return dofs_per_hex;
815}
816
817
818
819template <int dim>
820inline unsigned int
821FiniteElementData<dim>::n_dofs_per_face(unsigned int face_no,
822 unsigned int child_no) const
823{
824 (void)child_no;
825
826 return n_dofs_on_face[n_dofs_on_face.size() == 1 ? 0 : face_no];
827}
828
829
830
831template <int dim>
832inline unsigned int
834{
835 return dofs_per_face_max;
836}
837
838
839
840template <int dim>
841inline unsigned int
843{
844 return dofs_per_cell;
845}
846
847
848
849template <int dim>
850template <int structdim>
851inline unsigned int
852FiniteElementData<dim>::n_dofs_per_object(const unsigned int i) const
853{
854 switch (structdim)
855 {
856 case 0:
857 return n_dofs_per_vertex();
858 case 1:
859 return n_dofs_per_line();
860 case 2:
861 return n_dofs_per_quad((structdim == 2 && dim == 3) ? i : 0);
862 case 3:
863 return n_dofs_per_hex();
864 default:
865 Assert(false, ExcInternalError());
866 }
868}
869
870
871
872template <int dim>
873inline unsigned int
875{
876 return components;
877}
878
879
880
881template <int dim>
882inline const BlockIndices &
884{
885 return block_indices_data;
886}
887
888
889
890template <int dim>
891inline unsigned int
893{
894 return block_indices_data.size();
895}
896
897
898
899template <int dim>
900inline unsigned int
902{
903 return degree;
904}
905
906
907template <int dim>
908inline bool
909FiniteElementData<dim>::conforms(const Conformity space) const
910{
911 return ((space & conforming_space) == space);
912}
913
914
915
916template <int dim>
917unsigned int
919{
920 return first_line_index;
921}
922
923template <int dim>
924unsigned int
925FiniteElementData<dim>::get_first_quad_index(const unsigned int quad_no) const
926{
927 if (first_index_of_quads.size() == 1)
928 return first_index_of_quads[0] + quad_no * n_dofs_per_quad(0);
929 else
930 return first_index_of_quads[quad_no];
931}
932
933template <int dim>
934unsigned int
936{
937 return first_hex_index;
938}
939
940template <int dim>
941unsigned int
943 const unsigned int face_no) const
944{
945 return first_line_index_of_faces[first_line_index_of_faces.size() == 1 ?
946 0 :
947 face_no];
948}
949
950template <int dim>
951unsigned int
953 const unsigned int face_no) const
954{
955 return first_quad_index_of_faces[first_quad_index_of_faces.size() == 1 ?
956 0 :
957 face_no];
958}
959
960
961#endif // DOXYGEN
962
963
965
966#endif
unsigned int size() const
const unsigned int dofs_per_quad
Definition: fe_base.h:329
unsigned int get_first_line_index() const
unsigned int n_dofs_per_vertex() const
const unsigned int first_face_quad_index
Definition: fe_base.h:390
const std::vector< unsigned int > first_index_of_quads
Definition: fe_base.h:355
bool conforms(const Conformity) const
const std::vector< unsigned int > first_line_index_of_faces
Definition: fe_base.h:372
const unsigned int components
Definition: fe_base.h:429
const unsigned int degree
Definition: fe_base.h:435
const std::vector< unsigned int > n_dofs_on_face
Definition: fe_base.h:397
unsigned int n_dofs_per_cell() const
const unsigned int dofs_per_face
Definition: fe_base.h:405
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:180
const unsigned int dofs_per_line
Definition: fe_base.h:315
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
unsigned int n_unique_quads() const
ReferenceCell reference_cell() const
const unsigned int first_hex_index
Definition: fe_base.h:366
unsigned int n_dofs_per_object(const unsigned int i=0) const
const unsigned int first_face_line_index
Definition: fe_base.h:378
const unsigned int dofs_per_quad_max
Definition: fe_base.h:335
unsigned int get_first_face_line_index(const unsigned int face_no=0) const
const unsigned int first_line_index
Definition: fe_base.h:347
const unsigned int dofs_per_hex
Definition: fe_base.h:342
const BlockIndices block_indices_data
Definition: fe_base.h:447
unsigned int get_first_face_quad_index(const unsigned int face_no=0) const
const unsigned int number_unique_faces
Definition: fe_base.h:303
unsigned int max_dofs_per_quad() const
const BlockIndices & block_indices() const
const std::vector< unsigned int > first_quad_index_of_faces
Definition: fe_base.h:384
const ReferenceCell reference_cell_kind
Definition: fe_base.h:291
const unsigned int number_unique_quads
Definition: fe_base.h:297
const std::vector< unsigned int > n_dofs_on_quad
Definition: fe_base.h:322
const Conformity conforming_space
Definition: fe_base.h:440
const unsigned int first_quad_index
Definition: fe_base.h:361
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_base.h:309
const unsigned int dofs_per_cell
Definition: fe_base.h:419
unsigned int n_dofs_per_hex() const
FiniteElementData(const std::vector< unsigned int > &dofs_per_object, const unsigned int n_components, const unsigned int degree, const Conformity conformity=unknown, const BlockIndices &block_indices=BlockIndices())
Definition: fe_data.cc:102
static const unsigned int dimension
Definition: fe_base.h:285
const unsigned int dofs_per_face_max
Definition: fe_base.h:411
unsigned int get_first_hex_index() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcInternalError()
Domination operator&(const Domination d1, const Domination d2)
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
internal::GenericDoFsPerObject expand(const unsigned int dim, const std::vector< unsigned int > &dofs_per_object, const ::ReferenceCell reference_cell)
Definition: fe_data.cc:25
static const unsigned int invalid_unsigned_int
Definition: types.h:196
std::vector< std::vector< unsigned int > > object_index
Definition: fe_base.h:188
std::vector< std::vector< unsigned int > > first_object_index_on_face
Definition: fe_base.h:193
std::vector< std::vector< unsigned int > > dofs_per_object_inclusive
Definition: fe_base.h:183
std::vector< std::vector< unsigned int > > dofs_per_object_exclusive
Definition: fe_base.h:178