Reference documentation for deal.II version GIT 0980a66d4b 2023-03-23 20:20:03+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\}}\)
fe_data.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2022 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_data_h
17 #define dealii_fe_data_h
18 
19 #include <deal.II/base/config.h>
20 
22 
24 
26 
27 #include <vector>
28 
30 
31 // Forward declarations:
32 #ifndef DOXYGEN
33 template <int dim>
34 class FiniteElementData;
35 #endif
36 
42 {
92  {
113  };
114 
115 
134  inline Domination
135  operator&(const Domination d1, const Domination d2);
136 } // namespace FiniteElementDomination
137 
138 namespace internal
139 {
181  {
185  std::vector<std::vector<unsigned int>> dofs_per_object_exclusive;
186 
190  std::vector<std::vector<unsigned int>> dofs_per_object_inclusive;
191 
195  std::vector<std::vector<unsigned int>> object_index;
196 
200  std::vector<std::vector<unsigned int>> first_object_index_on_face;
201 
205  template <int dim>
206  static GenericDoFsPerObject
208  };
209 } // namespace internal
210 
225 template <int dim>
227 {
228 public:
260  {
264  unknown = 0x00,
265 
269  L2 = 0x01,
270 
275  Hcurl = 0x02,
276 
281  Hdiv = 0x04,
282 
286  H1 = Hcurl | Hdiv,
287 
292  H2 = 0x0e
293  };
294 
299  static constexpr unsigned int dimension = dim;
300 
301 private:
306 
311  const unsigned int number_unique_quads;
312 
317  const unsigned int number_unique_faces;
318 
319 public:
323  const unsigned int dofs_per_vertex;
324 
329  const unsigned int dofs_per_line;
330 
331 private:
336  const std::vector<unsigned int> n_dofs_on_quad;
337 
338 public:
343  const unsigned int dofs_per_quad;
344 
345 private:
349  const unsigned int dofs_per_quad_max;
350 
351 public:
356  const unsigned int dofs_per_hex;
357 
361  const unsigned int first_line_index;
362 
363 private:
369  const std::vector<unsigned int> first_index_of_quads;
370 
371 public:
375  const unsigned int first_quad_index;
376 
380  const unsigned int first_hex_index;
381 
382 private:
386  const std::vector<unsigned int> first_line_index_of_faces;
387 
388 public:
392  const unsigned int first_face_line_index;
393 
394 private:
398  const std::vector<unsigned int> first_quad_index_of_faces;
399 
400 public:
404  const unsigned int first_face_quad_index;
405 
406 private:
411  const std::vector<unsigned int> n_dofs_on_face;
412 
413 public:
419  const unsigned int dofs_per_face;
420 
421 private:
425  const unsigned int dofs_per_face_max;
426 
427 public:
433  const unsigned int dofs_per_cell;
434 
443  const unsigned int components;
444 
449  const unsigned int degree;
450 
455 
462 
501  FiniteElementData(const std::vector<unsigned int> &dofs_per_object,
502  const unsigned int n_components,
503  const unsigned int degree,
504  const Conformity conformity = unknown,
506 
511  FiniteElementData(const std::vector<unsigned int> &dofs_per_object,
513  const unsigned int n_components,
514  const unsigned int degree,
515  const Conformity conformity = unknown,
517 
526  const unsigned int n_components,
527  const unsigned int degree,
528  const Conformity conformity = unknown,
530 
537  reference_cell() const;
538 
543  unsigned int
544  n_unique_quads() const;
545 
550  unsigned int
551  n_unique_faces() const;
552 
556  unsigned int
558 
562  unsigned int
564 
568  unsigned int
569  n_dofs_per_quad(unsigned int face_no = 0) const;
570 
575  unsigned int
577 
581  unsigned int
582  n_dofs_per_hex() const;
583 
588  unsigned int
589  n_dofs_per_face(unsigned int face_no = 0, unsigned int child = 0) const;
590 
595  unsigned int
597 
602  unsigned int
604 
613  template <int structdim>
614  unsigned int
615  n_dofs_per_object(const unsigned int i = 0) const;
616 
622  unsigned int
623  n_components() const;
624 
630  unsigned int
631  n_blocks() const;
632 
636  const BlockIndices &
637  block_indices() const;
638 
645  unsigned int
646  tensor_degree() const;
647 
654  bool
655  conforms(const Conformity) const;
656 
660  bool
661  operator==(const FiniteElementData &) const;
662 
666  unsigned int
668 
672  unsigned int
673  get_first_quad_index(const unsigned int quad_no = 0) const;
674 
678  unsigned int
680 
684  unsigned int
685  get_first_face_line_index(const unsigned int face_no = 0) const;
686 
690  unsigned int
691  get_first_face_quad_index(const unsigned int face_no = 0) const;
692 };
693 
694 namespace internal
695 {
701  expand(const unsigned int dim,
702  const std::vector<unsigned int> &dofs_per_object,
704 } // namespace internal
705 
706 
707 
708 // --------- inline and template functions ---------------
709 
710 
711 #ifndef DOXYGEN
712 
713 namespace FiniteElementDomination
714 {
715  inline Domination
716  operator&(const Domination d1, const Domination d2)
717  {
718  // go through the entire list of possibilities. note that if we were into
719  // speed, obfuscation and cared enough, we could implement this operator
720  // by doing a bitwise & (and) if we gave these values to the enum values:
721  // neither_element_dominates=0, this_element_dominates=1,
722  // other_element_dominates=2, either_element_can_dominate=3
723  // =this_element_dominates|other_element_dominates
724  switch (d1)
725  {
727  if ((d2 == this_element_dominates) ||
728  (d2 == either_element_can_dominate) || (d2 == no_requirements))
729  return this_element_dominates;
730  else
732 
734  if ((d2 == other_element_dominates) ||
735  (d2 == either_element_can_dominate) || (d2 == no_requirements))
737  else
739 
742 
744  if (d2 == no_requirements)
746  else
747  return d2;
748 
749  case no_requirements:
750  return d2;
751 
752  default:
753  // shouldn't get here
754  Assert(false, ExcInternalError());
755  }
756 
758  }
759 } // namespace FiniteElementDomination
760 
761 
762 template <int dim>
763 inline ReferenceCell
765 {
766  return reference_cell_kind;
767 }
768 
769 
770 
771 template <int dim>
772 inline unsigned int
774 {
775  return number_unique_quads;
776 }
777 
778 
779 
780 template <int dim>
781 inline unsigned int
783 {
784  return number_unique_faces;
785 }
786 
787 
788 
789 template <int dim>
790 inline unsigned int
792 {
793  return dofs_per_vertex;
794 }
795 
796 
797 
798 template <int dim>
799 inline unsigned int
801 {
802  return dofs_per_line;
803 }
804 
805 
806 
807 template <int dim>
808 inline unsigned int
809 FiniteElementData<dim>::n_dofs_per_quad(unsigned int face_no) const
810 {
811  return n_dofs_on_quad[n_dofs_on_quad.size() == 1 ? 0 : face_no];
812 }
813 
814 
815 
816 template <int dim>
817 inline unsigned int
819 {
820  return dofs_per_quad_max;
821 }
822 
823 
824 
825 template <int dim>
826 inline unsigned int
828 {
829  return dofs_per_hex;
830 }
831 
832 
833 
834 template <int dim>
835 inline unsigned int
836 FiniteElementData<dim>::n_dofs_per_face(unsigned int face_no,
837  unsigned int child_no) const
838 {
839  (void)child_no;
840 
841  return n_dofs_on_face[n_dofs_on_face.size() == 1 ? 0 : face_no];
842 }
843 
844 
845 
846 template <int dim>
847 inline unsigned int
849 {
850  return dofs_per_face_max;
851 }
852 
853 
854 
855 template <int dim>
856 inline unsigned int
858 {
859  return dofs_per_cell;
860 }
861 
862 
863 
864 template <int dim>
865 template <int structdim>
866 inline unsigned int
867 FiniteElementData<dim>::n_dofs_per_object(const unsigned int i) const
868 {
869  switch (structdim)
870  {
871  case 0:
872  return n_dofs_per_vertex();
873  case 1:
874  return n_dofs_per_line();
875  case 2:
876  return n_dofs_per_quad((structdim == 2 && dim == 3) ? i : 0);
877  case 3:
878  return n_dofs_per_hex();
879  default:
880  Assert(false, ExcInternalError());
881  }
883 }
884 
885 
886 
887 template <int dim>
888 inline unsigned int
890 {
891  return components;
892 }
893 
894 
895 
896 template <int dim>
897 inline const BlockIndices &
899 {
900  return block_indices_data;
901 }
902 
903 
904 
905 template <int dim>
906 inline unsigned int
908 {
909  return block_indices_data.size();
910 }
911 
912 
913 
914 template <int dim>
915 inline unsigned int
917 {
918  return degree;
919 }
920 
921 
922 template <int dim>
923 inline bool
924 FiniteElementData<dim>::conforms(const Conformity space) const
925 {
926  return ((space & conforming_space) == space);
927 }
928 
929 
930 
931 template <int dim>
932 unsigned int
934 {
935  return first_line_index;
936 }
937 
938 template <int dim>
939 unsigned int
940 FiniteElementData<dim>::get_first_quad_index(const unsigned int quad_no) const
941 {
942  if (first_index_of_quads.size() == 1)
943  return first_index_of_quads[0] + quad_no * n_dofs_per_quad(0);
944  else
945  return first_index_of_quads[quad_no];
946 }
947 
948 template <int dim>
949 unsigned int
951 {
952  return first_hex_index;
953 }
954 
955 template <int dim>
956 unsigned int
958  const unsigned int face_no) const
959 {
960  return first_line_index_of_faces[first_line_index_of_faces.size() == 1 ?
961  0 :
962  face_no];
963 }
964 
965 template <int dim>
966 unsigned int
968  const unsigned int face_no) const
969 {
970  return first_quad_index_of_faces[first_quad_index_of_faces.size() == 1 ?
971  0 :
972  face_no];
973 }
974 
975 template <int dim>
978 {
979  const auto reference_cell = fe.reference_cell();
980 
982 
983  result.dofs_per_object_exclusive.resize(4);
984  result.dofs_per_object_inclusive.resize(4);
985  result.object_index.resize(4);
986 
987  unsigned int counter = 0;
988 
989  for (unsigned int v : reference_cell.vertex_indices())
990  {
991  const auto c = fe.template n_dofs_per_object<0>(v);
992 
993  result.dofs_per_object_exclusive[0].emplace_back(c);
994  result.dofs_per_object_inclusive[0].emplace_back(c);
995  result.object_index[0].emplace_back(counter);
996 
997  counter += c;
998  }
999 
1000  if (dim >= 2)
1001  for (unsigned int l : reference_cell.line_indices())
1002  {
1003  const auto c = fe.template n_dofs_per_object<1>(l);
1004 
1005  result.dofs_per_object_exclusive[1].emplace_back(c);
1006  result.dofs_per_object_inclusive[1].emplace_back(
1007  c + 2 * fe.template n_dofs_per_object<0>());
1008  result.object_index[1].emplace_back(counter);
1009 
1010  counter += c;
1011  }
1012 
1013  if (dim == 3)
1014  for (unsigned int f : reference_cell.face_indices())
1015  {
1016  const auto c = fe.template n_dofs_per_object<2>(f);
1017 
1018  result.dofs_per_object_exclusive[2].emplace_back(c);
1019  result.dofs_per_object_inclusive[2].emplace_back(fe.n_dofs_per_face(f));
1020  result.object_index[2].emplace_back(counter);
1021 
1022  counter += c;
1023  }
1024 
1025  {
1026  const auto c = fe.template n_dofs_per_object<dim>();
1027 
1028  result.dofs_per_object_exclusive[dim].emplace_back(c);
1029  result.dofs_per_object_inclusive[dim].emplace_back(fe.n_dofs_per_cell());
1030  result.object_index[dim].emplace_back(counter);
1031 
1032  counter += c;
1033  }
1034 
1035  for (unsigned int d = dim + 1; d <= 3; ++d)
1036  {
1037  result.dofs_per_object_exclusive[d].emplace_back(0);
1038  result.dofs_per_object_inclusive[d].emplace_back(0);
1039  result.object_index[d].emplace_back(counter);
1040  }
1041 
1042  result.first_object_index_on_face.resize(3);
1043  for (unsigned int face_no : reference_cell.face_indices())
1044  {
1045  result.first_object_index_on_face[0].emplace_back(0);
1046 
1047  result.first_object_index_on_face[1].emplace_back(
1048  fe.get_first_face_line_index(face_no));
1049 
1050  result.first_object_index_on_face[2].emplace_back(
1051  fe.get_first_face_quad_index(face_no));
1052  }
1053 
1054  return result;
1055 }
1056 
1057 
1058 #endif // DOXYGEN
1059 
1060 
1062 
1063 #endif
unsigned int size() const
const unsigned int dofs_per_quad
Definition: fe_data.h:343
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:404
const std::vector< unsigned int > first_index_of_quads
Definition: fe_data.h:369
bool conforms(const Conformity) const
const std::vector< unsigned int > first_line_index_of_faces
Definition: fe_data.h:386
const unsigned int components
Definition: fe_data.h:443
const unsigned int degree
Definition: fe_data.h:449
const std::vector< unsigned int > n_dofs_on_face
Definition: fe_data.h:411
unsigned int n_dofs_per_cell() const
const unsigned int dofs_per_face
Definition: fe_data.h:419
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:196
const unsigned int dofs_per_line
Definition: fe_data.h:329
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_data.h:380
unsigned int n_dofs_per_object(const unsigned int i=0) const
const unsigned int first_face_line_index
Definition: fe_data.h:392
const unsigned int dofs_per_quad_max
Definition: fe_data.h:349
unsigned int get_first_face_line_index(const unsigned int face_no=0) const
const unsigned int first_line_index
Definition: fe_data.h:361
const unsigned int dofs_per_hex
Definition: fe_data.h:356
const BlockIndices block_indices_data
Definition: fe_data.h:461
unsigned int get_first_face_quad_index(const unsigned int face_no=0) const
const unsigned int number_unique_faces
Definition: fe_data.h:317
unsigned int max_dofs_per_quad() const
const std::vector< unsigned int > first_quad_index_of_faces
Definition: fe_data.h:398
const ReferenceCell reference_cell_kind
Definition: fe_data.h:305
const BlockIndices & block_indices() const
const unsigned int number_unique_quads
Definition: fe_data.h:311
const std::vector< unsigned int > n_dofs_on_quad
Definition: fe_data.h:336
const Conformity conforming_space
Definition: fe_data.h:454
const unsigned int first_quad_index
Definition: fe_data.h:375
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:323
const unsigned int dofs_per_cell
Definition: fe_data.h:433
static constexpr unsigned int dimension
Definition: fe_data.h:299
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:115
const unsigned int dofs_per_face_max
Definition: fe_data.h:425
unsigned int get_first_hex_index() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:474
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:475
static ::ExceptionBase & ExcInternalError()
#define Assert(cond, exc)
Definition: exceptions.h:1586
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)
Tensor< 2, dim, Number > l(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:25
static const unsigned int invalid_unsigned_int
Definition: types.h:213
std::vector< std::vector< unsigned int > > object_index
Definition: fe_data.h:195
std::vector< std::vector< unsigned int > > first_object_index_on_face
Definition: fe_data.h:200
std::vector< std::vector< unsigned int > > dofs_per_object_inclusive
Definition: fe_data.h:190
std::vector< std::vector< unsigned int > > dofs_per_object_exclusive
Definition: fe_data.h:185
static GenericDoFsPerObject generate(const FiniteElementData< dim > &fe)