Reference documentation for deal.II version 8.5.1
geometry_info.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2016 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/geometry_info.h>
17 #include <deal.II/base/tensor.h>
18 
19 DEAL_II_NAMESPACE_OPEN
20 
21 
22 unsigned int
24 {
25  return 0;
26 }
27 
28 
29 
30 template <int dim> const unsigned int GeometryInfo<dim>::max_children_per_cell;
31 template <int dim> const unsigned int GeometryInfo<dim>::faces_per_cell;
32 template <int dim> const unsigned int GeometryInfo<dim>::max_children_per_face;
33 template <int dim> const unsigned int GeometryInfo<dim>::vertices_per_cell;
34 template <int dim> const unsigned int GeometryInfo<dim>::vertices_per_face;
35 template <int dim> const unsigned int GeometryInfo<dim>::lines_per_face;
36 template <int dim> const unsigned int GeometryInfo<dim>::quads_per_face;
37 template <int dim> const unsigned int GeometryInfo<dim>::lines_per_cell;
38 template <int dim> const unsigned int GeometryInfo<dim>::quads_per_cell;
39 template <int dim> const unsigned int GeometryInfo<dim>::hexes_per_cell;
40 
41 
42 using namespace numbers;
43 
44 template <>
45 const unsigned int
47  = { 0, 0 };
48 
49 template <>
50 const unsigned int
52  = { 0, 0, 1, 1 };
53 
54 template <>
55 const unsigned int
57  = { 0, 0, 1, 1, 2, 2 };
58 
59 template <>
60 const unsigned int
62  = { 0, 0, 1, 1, 2, 2, 3, 3 };
63 
64 
65 
66 template <>
67 const int
69  = { -1, 1 };
70 
71 template <>
72 const int
74  = { -1, 1, -1, 1 };
75 
76 template <>
77 const int
79  = { -1, 1, -1, 1, -1, 1 };
80 
81 template <>
82 const int
84  = { -1, 1, -1, 1, -1, 1, -1, 1 };
85 
86 
87 
88 template <>
89 const unsigned int
90 GeometryInfo<1>::opposite_face[faces_per_cell]
91  = { 1, 0 };
92 
93 template <>
94 const unsigned int
95 GeometryInfo<2>::opposite_face[faces_per_cell]
96  = { 1, 0, 3, 2 };
97 
98 template <>
99 const unsigned int
100 GeometryInfo<3>::opposite_face[faces_per_cell]
101  = { 1, 0, 3, 2, 5, 4 };
102 
103 template <>
104 const unsigned int
105 GeometryInfo<4>::opposite_face[faces_per_cell]
106  = { 1, 0, 3, 2, 5, 4, 7, 6 };
107 
108 
109 
110 template <>
112  = { 0, 1};
113 
114 template <>
116  = { 0, 1, 3, 2};
117 
118 template <>
120  = { 0, 1, 5, 4, 2, 3, 7, 6};
121 
122 template <>
140  };
141 
142 
143 template <>
145  = { 0, 1};
146 
147 template <>
149  = { 0, 2, 1, 3};
150 
151 template <>
153  = { 0, 4, 2, 6, 1, 5, 3, 7};
154 
155 template <>
173  };
174 
175 template <>
176 const unsigned int GeometryInfo<1>::vertex_to_face
178 = { { 0 },
179  { 1 }
180 };
181 
182 template <>
183 const unsigned int GeometryInfo<2>::vertex_to_face
185 = { { 0, 2 },
186  { 1, 2 },
187  { 0, 3 },
188  { 1, 3 }
189 };
190 
191 template <>
192 const unsigned int GeometryInfo<3>::vertex_to_face
194 = { { 0, 2, 4 },
195  { 1, 2, 4 },
196  { 0, 3, 4 },
197  { 1, 3, 4 },
198  { 0, 2, 5 },
199  { 1, 2, 5 },
200  { 0, 3, 5 },
201  { 1, 3, 5 }
202 };
203 
204 template <>
205 const unsigned int GeometryInfo<4>::vertex_to_face
223 };
224 
225 
226 template<int dim>
227 unsigned int
229 {
230  static const unsigned int n_children[RefinementCase<3>::cut_xyz+1]=
231  {0, 2, 2, 4, 2, 4, 4, 8};
232 
233  return n_children[ref_case];
234 }
235 
236 
237 template<>
238 unsigned int
240 {
241  Assert(false, ExcImpossibleInDim(1));
242  return 0;
243 }
244 
245 
246 
247 template<>
248 unsigned int
250 {
251  return (subface_case == internal::SubfaceCase<2>::case_x) ? 2 : 0;
252 }
253 
254 
255 
256 template<>
257 unsigned int
259 {
260  static const unsigned int nsubs[internal::SubfaceCase<3>::case_isotropic+1]=
261  {0, 2, 3, 3, 4, 2, 3, 3, 4, 4};
262  return nsubs[subface_case];
263 }
264 
265 
266 template<>
267 double
269  const unsigned int)
270 {
271  return 1;
272 }
273 
274 
275 template<>
276 double
278  const unsigned int)
279 {
280  const unsigned int dim=2;
281 
282  double ratio=1;
283  switch (subface_case)
284  {
286  // Here, an
287  // Assert(false,ExcInternalError())
288  // would be the right
289  // choice, but
290  // unfortunately the
291  // current function is
292  // also called for faces
293  // without children (see
294  // tests/fe/mapping.cc).
295 // Assert(false, ExcMessage("Face has no subfaces."));
296  // Furthermore, assign
297  // following value as
298  // otherwise the
299  // bits/volume_x tests
300  // break
302  break;
304  ratio=0.5;
305  break;
306  default:
307  // there should be no
308  // cases left
309  Assert(false, ExcInternalError());
310  break;
311  }
312 
313  return ratio;
314 }
315 
316 
317 template<>
318 double
320  const unsigned int subface_no)
321 {
322  const unsigned int dim=3;
323 
324  double ratio=1;
325  switch (subface_case)
326  {
328  // Here, an
329  // Assert(false,ExcInternalError())
330  // would be the right
331  // choice, but
332  // unfortunately the
333  // current function is
334  // also called for faces
335  // without children (see
336  // tests/bits/mesh_3d_16.cc). Add
337  // following switch to
338  // avoid diffs in
339  // tests/bits/mesh_3d_16
341  break;
344  ratio=0.5;
345  break;
349  ratio=0.25;
350  break;
353  if (subface_no<2)
354  ratio=0.25;
355  else
356  ratio=0.5;
357  break;
360  if (subface_no==0)
361  ratio=0.5;
362  else
363  ratio=0.25;
364  break;
365  default:
366  // there should be no
367  // cases left
368  Assert(false, ExcInternalError());
369  break;
370  }
371 
372  return ratio;
373 }
374 
375 
376 
377 template<>
380  const unsigned int,
381  const bool,
382  const bool,
383  const bool)
384 {
385  Assert(false, ExcImpossibleInDim(1));
386 
388 }
389 
390 
391 template<>
393 GeometryInfo<2>::face_refinement_case(const RefinementCase<2> &cell_refinement_case,
394  const unsigned int face_no,
395  const bool,
396  const bool,
397  const bool)
398 {
399  const unsigned int dim=2;
400  Assert(cell_refinement_case<RefinementCase<dim>::isotropic_refinement+1,
401  ExcIndexRange(cell_refinement_case, 0, RefinementCase<dim>::isotropic_refinement+1));
404 
405  static const RefinementCase<dim-1>
407  {
408  {
409  RefinementCase<dim-1>::no_refinement, // no_refinement
410  RefinementCase<dim-1>::no_refinement
411  },
412 
413  {
414  RefinementCase<dim-1>::no_refinement,
415  RefinementCase<dim-1>::cut_x
416  },
417 
418  {
419  RefinementCase<dim-1>::cut_x,
420  RefinementCase<dim-1>::no_refinement
421  },
422 
423  {
424  RefinementCase<dim-1>::cut_x, // cut_xy
425  RefinementCase<dim-1>::cut_x
426  }
427  };
428 
429  return ref_cases[cell_refinement_case][face_no/2];
430 }
431 
432 
433 template<>
435 GeometryInfo<3>::face_refinement_case(const RefinementCase<3> &cell_refinement_case,
436  const unsigned int face_no,
437  const bool face_orientation,
438  const bool /*face_flip*/,
439  const bool face_rotation)
440 {
441  const unsigned int dim=3;
442  Assert(cell_refinement_case<RefinementCase<dim>::isotropic_refinement+1,
443  ExcIndexRange(cell_refinement_case, 0, RefinementCase<dim>::isotropic_refinement+1));
446 
447  static const RefinementCase<dim-1>
449  {
450  {
451  RefinementCase<dim-1>::no_refinement, // no_refinement
452  RefinementCase<dim-1>::no_refinement,
453  RefinementCase<dim-1>::no_refinement
454  },
455 
456  {
457  RefinementCase<dim-1>::no_refinement, // cut_x
458  RefinementCase<dim-1>::cut_y,
459  RefinementCase<dim-1>::cut_x
460  },
461 
462  {
463  RefinementCase<dim-1>::cut_x, // cut_y
464  RefinementCase<dim-1>::no_refinement,
465  RefinementCase<dim-1>::cut_y
466  },
467 
468  {
469  RefinementCase<dim-1>::cut_x, // cut_xy
470  RefinementCase<dim-1>::cut_y,
471  RefinementCase<dim-1>::cut_xy
472  },
473 
474  {
475  RefinementCase<dim-1>::cut_y, // cut_z
476  RefinementCase<dim-1>::cut_x,
477  RefinementCase<dim-1>::no_refinement
478  },
479 
480  {
481  RefinementCase<dim-1>::cut_y, // cut_xz
482  RefinementCase<dim-1>::cut_xy,
483  RefinementCase<dim-1>::cut_x
484  },
485 
486  {
487  RefinementCase<dim-1>::cut_xy, // cut_yz
488  RefinementCase<dim-1>::cut_x,
489  RefinementCase<dim-1>::cut_y
490  },
491 
492  {
493  RefinementCase<dim-1>::cut_xy, // cut_xyz
494  RefinementCase<dim-1>::cut_xy,
495  RefinementCase<dim-1>::cut_xy
496  },
497  };
498 
499  const RefinementCase<dim-1> ref_case=ref_cases[cell_refinement_case][face_no/2];
500 
501  static const RefinementCase<dim-1> flip[4]=
502  {
503  RefinementCase<dim-1>::no_refinement,
504  RefinementCase<dim-1>::cut_y,
505  RefinementCase<dim-1>::cut_x,
506  RefinementCase<dim-1>::cut_xy
507  };
508 
509  // correct the ref_case for face_orientation
510  // and face_rotation. for face_orientation,
511  // 'true' is the default value whereas for
512  // face_rotation, 'false' is standard. If
513  // <tt>face_rotation==face_orientation</tt>,
514  // then one of them is non-standard and we
515  // have to swap cut_x and cut_y, otherwise no
516  // change is necessary. face_flip has no
517  // influence. however, in order to keep the
518  // interface consistent with other functions,
519  // we still include it as an argument to this
520  // function
521  return (face_orientation==face_rotation) ? flip[ref_case] : ref_case;
522 }
523 
524 
525 
526 template<>
528 GeometryInfo<1>::line_refinement_case(const RefinementCase<1> &cell_refinement_case,
529  const unsigned int line_no)
530 {
531  (void)line_no;
532  const unsigned int dim = 1;
533  (void)dim;
534  Assert(cell_refinement_case<RefinementCase<dim>::isotropic_refinement+1,
535  ExcIndexRange(cell_refinement_case, 0, RefinementCase<dim>::isotropic_refinement+1));
538 
539  return cell_refinement_case;
540 }
541 
542 
543 template<>
545 GeometryInfo<2>::line_refinement_case(const RefinementCase<2> &cell_refinement_case,
546  const unsigned int line_no)
547 {
548  // Assertions are in face_refinement_case()
549  return face_refinement_case(cell_refinement_case, line_no);
550 }
551 
552 
553 template<>
555 GeometryInfo<3>::line_refinement_case(const RefinementCase<3> &cell_refinement_case,
556  const unsigned int line_no)
557 {
558  const unsigned int dim=3;
559  Assert(cell_refinement_case<RefinementCase<dim>::isotropic_refinement+1,
560  ExcIndexRange(cell_refinement_case, 0, RefinementCase<dim>::isotropic_refinement+1));
563 
564  // array indicating, which simple refine
565  // case cuts a line in direction x, y or
566  // z. For example, cut_y and everything
567  // containing cut_y (cut_xy, cut_yz,
568  // cut_xyz) cuts lines, which are in y
569  // direction.
570  static const RefinementCase<dim>
571  cut_one[dim] =
572  {
576  };
577 
578  // order the direction of lines
579  // 0->x, 1->y, 2->z
580  static const unsigned int direction[lines_per_cell]=
581  {1,1,0,0,1,1,0,0,2,2,2,2};
582 
583  return ((cell_refinement_case & cut_one[direction[line_no]]) ?
585 }
586 
587 
588 
589 template<>
592  const unsigned int,
593  const bool,
594  const bool,
595  const bool)
596 {
597  const unsigned int dim = 1;
598  Assert(false, ExcImpossibleInDim(dim));
599 
601 }
602 
603 
604 template<>
607  const unsigned int face_no,
608  const bool,
609  const bool,
610  const bool)
611 {
612  const unsigned int dim = 2;
614  ExcIndexRange(face_refinement_case, 0, RefinementCase<dim-1>::isotropic_refinement+1));
617 
618  if (face_refinement_case==RefinementCase<dim>::cut_x)
620  else
622 }
623 
624 
625 template<>
628  const unsigned int face_no,
629  const bool face_orientation,
630  const bool /*face_flip*/,
631  const bool face_rotation)
632 {
633  const unsigned int dim=3;
635  ExcIndexRange(face_refinement_case, 0, RefinementCase<dim-1>::isotropic_refinement+1));
638 
639  static const RefinementCase<2> flip[4]=
640  {
645  };
646 
647  // correct the face_refinement_case for
648  // face_orientation and face_rotation. for
649  // face_orientation, 'true' is the default
650  // value whereas for face_rotation, 'false'
651  // is standard. If
652  // <tt>face_rotation==face_orientation</tt>,
653  // then one of them is non-standard and we
654  // have to swap cut_x and cut_y, otherwise no
655  // change is necessary. face_flip has no
656  // influence. however, in order to keep the
657  // interface consistent with other functions,
658  // we still include it as an argument to this
659  // function
660  const RefinementCase<dim-1> std_face_ref = (face_orientation==face_rotation) ? flip[face_refinement_case] : face_refinement_case;
661 
662  static const RefinementCase<dim> face_to_cell[3][4]=
663  {
664  {
665  RefinementCase<dim>::no_refinement, // faces 0 and 1
666  RefinementCase<dim>::cut_y, // cut_x in face 0 means cut_y for the cell
669  },
670 
671  {
672  RefinementCase<dim>::no_refinement, // faces 2 and 3 (note that x and y are "exchanged on faces 2 and 3")
676  },
677 
678  {
679  RefinementCase<dim>::no_refinement, // faces 4 and 5
683  }
684  };
685 
686  return face_to_cell[face_no/2][std_face_ref];
687 }
688 
689 
690 
691 template<>
694 {
695  (void)line_no;
696  Assert(line_no==0, ExcIndexRange(line_no,0,1));
697 
699 }
700 
701 
702 template<>
705 {
706  const unsigned int dim = 2;
707  (void)dim;
710 
711  return (line_no/2) ? RefinementCase<2>::cut_x : RefinementCase<2>::cut_y;
712 }
713 
714 
715 template<>
718 {
719  const unsigned int dim=3;
722 
723  static const RefinementCase<dim> ref_cases[6]=
724  {
725  RefinementCase<dim>::cut_y, // lines 0 and 1
726  RefinementCase<dim>::cut_x, // lines 2 and 3
727  RefinementCase<dim>::cut_y, // lines 4 and 5
728  RefinementCase<dim>::cut_x, // lines 6 and 7
729  RefinementCase<dim>::cut_z, // lines 8 and 9
731  }; // lines 10 and 11
732 
733  return ref_cases[line_no/2];
734 }
735 
736 
737 
738 template <>
739 unsigned int
740 GeometryInfo<3>::standard_to_real_face_vertex(const unsigned int vertex,
741  const bool face_orientation,
742  const bool face_flip,
743  const bool face_rotation)
744 {
747 
748  // set up a table to make sure that
749  // we handle non-standard faces correctly
750  //
751  // so set up a table that for each vertex (of
752  // a quad in standard position) describes
753  // which vertex to take
754  //
755  // first index: four vertices 0...3
756  //
757  // second index: face_orientation; 0:
758  // opposite normal, 1: standard
759  //
760  // third index: face_flip; 0: standard, 1:
761  // face rotated by 180 degrees
762  //
763  // forth index: face_rotation: 0: standard,
764  // 1: face rotated by 90 degrees
765 
766  static const unsigned int vertex_translation[4][2][2][2] =
767  {
768  { { { 0, 2 }, // vertex 0, face_orientation=false, face_flip=false, face_rotation=false and true
769  { 3, 1 }
770  }, // vertex 0, face_orientation=false, face_flip=true, face_rotation=false and true
771  { { 0, 2 }, // vertex 0, face_orientation=true, face_flip=false, face_rotation=false and true
772  { 3, 1 }
773  }
774  },// vertex 0, face_orientation=true, face_flip=true, face_rotation=false and true
775 
776  { { { 2, 3 }, // vertex 1 ...
777  { 1, 0 }
778  },
779  { { 1, 0 },
780  { 2, 3 }
781  }
782  },
783 
784  { { { 1, 0 }, // vertex 2 ...
785  { 2, 3 }
786  },
787  { { 2, 3 },
788  { 1, 0 }
789  }
790  },
791 
792  { { { 3, 1 }, // vertex 3 ...
793  { 0, 2 }
794  },
795  { { 3, 1 },
796  { 0, 2 }
797  }
798  }
799  };
800 
801  return vertex_translation[vertex][face_orientation][face_flip][face_rotation];
802 }
803 
804 
805 
806 template <int dim>
807 unsigned int
809  const bool,
810  const bool,
811  const bool)
812 {
813  Assert(dim>1, ExcImpossibleInDim(dim));
816  return vertex;
817 }
818 
819 
820 
821 template <>
822 unsigned int
823 GeometryInfo<3>::real_to_standard_face_vertex(const unsigned int vertex,
824  const bool face_orientation,
825  const bool face_flip,
826  const bool face_rotation)
827 {
830 
831  // set up a table to make sure that
832  // we handle non-standard faces correctly
833  //
834  // so set up a table that for each vertex (of
835  // a quad in standard position) describes
836  // which vertex to take
837  //
838  // first index: four vertices 0...3
839  //
840  // second index: face_orientation; 0:
841  // opposite normal, 1: standard
842  //
843  // third index: face_flip; 0: standard, 1:
844  // face rotated by 180 degrees
845  //
846  // forth index: face_rotation: 0: standard,
847  // 1: face rotated by 90 degrees
848 
849  static const unsigned int vertex_translation[4][2][2][2] =
850  {
851  { { { 0, 2 }, // vertex 0, face_orientation=false, face_flip=false, face_rotation=false and true
852  { 3, 1 }
853  }, // vertex 0, face_orientation=false, face_flip=true, face_rotation=false and true
854  { { 0, 1 }, // vertex 0, face_orientation=true, face_flip=false, face_rotation=false and true
855  { 3, 2 }
856  }
857  },// vertex 0, face_orientation=true, face_flip=true, face_rotation=false and true
858 
859  { { { 2, 3 }, // vertex 1 ...
860  { 1, 0 }
861  },
862  { { 1, 3 },
863  { 2, 0 }
864  }
865  },
866 
867  { { { 1, 0 }, // vertex 2 ...
868  { 2, 3 }
869  },
870  { { 2, 0 },
871  { 1, 3 }
872  }
873  },
874 
875  { { { 3, 1 }, // vertex 3 ...
876  { 0, 2 }
877  },
878  { { 3, 2 },
879  { 0, 1 }
880  }
881  }
882  };
883 
884  return vertex_translation[vertex][face_orientation][face_flip][face_rotation];
885 }
886 
887 
888 
889 template <int dim>
890 unsigned int
892  const bool,
893  const bool,
894  const bool)
895 {
896  Assert(dim>1, ExcImpossibleInDim(dim));
899  return vertex;
900 }
901 
902 
903 
904 template <>
905 unsigned int
906 GeometryInfo<3>::standard_to_real_face_line(const unsigned int line,
907  const bool face_orientation,
908  const bool face_flip,
909  const bool face_rotation)
910 {
913 
914 
915  // make sure we handle
916  // non-standard faces correctly
917  //
918  // so set up a table that for each line (of a
919  // quad) describes which line to take
920  //
921  // first index: four lines 0...3
922  //
923  // second index: face_orientation; 0:
924  // opposite normal, 1: standard
925  //
926  // third index: face_flip; 0: standard, 1:
927  // face rotated by 180 degrees
928  //
929  // forth index: face_rotation: 0: standard,
930  // 1: face rotated by 90 degrees
931 
932  static const unsigned int line_translation[4][2][2][2] =
933  {
934  { { { 2, 0 }, // line 0, face_orientation=false, face_flip=false, face_rotation=false and true
935  { 3, 1 }
936  }, // line 0, face_orientation=false, face_flip=true, face_rotation=false and true
937  { { 0, 3 }, // line 0, face_orientation=true, face_flip=false, face_rotation=false and true
938  { 1, 2 }
939  }
940  },// line 0, face_orientation=true, face_flip=true, face_rotation=false and true
941 
942  { { { 3, 1 }, // line 1 ...
943  { 2, 0 }
944  },
945  { { 1, 2 },
946  { 0, 3 }
947  }
948  },
949 
950  { { { 0, 3 }, // line 2 ...
951  { 1, 2 }
952  },
953  { { 2, 0 },
954  { 3, 1 }
955  }
956  },
957 
958  { { { 1, 2 }, // line 3 ...
959  { 0, 3 }
960  },
961  { { 3, 1 },
962  { 2, 0 }
963  }
964  }
965  };
966 
967  return line_translation[line][face_orientation][face_flip][face_rotation];
968 }
969 
970 
971 
972 template <int dim>
973 unsigned int
975  const bool,
976  const bool,
977  const bool)
978 {
979  Assert(false, ExcNotImplemented());
980  return line;
981 }
982 
983 
984 
985 template <>
986 unsigned int
987 GeometryInfo<3>::real_to_standard_face_line(const unsigned int line,
988  const bool face_orientation,
989  const bool face_flip,
990  const bool face_rotation)
991 {
994 
995 
996  // make sure we handle
997  // non-standard faces correctly
998  //
999  // so set up a table that for each line (of a
1000  // quad) describes which line to take
1001  //
1002  // first index: four lines 0...3
1003  //
1004  // second index: face_orientation; 0:
1005  // opposite normal, 1: standard
1006  //
1007  // third index: face_flip; 0: standard, 1:
1008  // face rotated by 180 degrees
1009  //
1010  // forth index: face_rotation: 0: standard,
1011  // 1: face rotated by 90 degrees
1012 
1013  static const unsigned int line_translation[4][2][2][2] =
1014  {
1015  { { { 2, 0 }, // line 0, face_orientation=false, face_flip=false, face_rotation=false and true
1016  { 3, 1 }
1017  }, // line 0, face_orientation=false, face_flip=true, face_rotation=false and true
1018  { { 0, 2 }, // line 0, face_orientation=true, face_flip=false, face_rotation=false and true
1019  { 1, 3 }
1020  }
1021  },// line 0, face_orientation=true, face_flip=true, face_rotation=false and true
1022 
1023  { { { 3, 1 }, // line 1 ...
1024  { 2, 0 }
1025  },
1026  { { 1, 3 },
1027  { 0, 2 }
1028  }
1029  },
1030 
1031  { { { 0, 3 }, // line 2 ...
1032  { 1, 2 }
1033  },
1034  { { 2, 1 },
1035  { 3, 0 }
1036  }
1037  },
1038 
1039  { { { 1, 2 }, // line 3 ...
1040  { 0, 3 }
1041  },
1042  { { 3, 0 },
1043  { 2, 1 }
1044  }
1045  }
1046  };
1047 
1048  return line_translation[line][face_orientation][face_flip][face_rotation];
1049 }
1050 
1051 
1052 
1053 template <int dim>
1054 unsigned int
1056  const bool,
1057  const bool,
1058  const bool)
1059 {
1060  Assert(false, ExcNotImplemented());
1061  return line;
1062 }
1063 
1064 
1065 
1066 template <>
1067 unsigned int
1069  const unsigned int face,
1070  const unsigned int subface,
1071  const bool, const bool, const bool,
1072  const RefinementCase<0> &)
1073 {
1074  (void)subface;
1075  Assert (face<faces_per_cell, ExcIndexRange(face, 0, faces_per_cell));
1076  Assert (subface<max_children_per_face,
1077  ExcIndexRange(subface, 0, max_children_per_face));
1078 
1079  return face;
1080 }
1081 
1082 
1083 
1084 template <>
1085 unsigned int
1087  const unsigned int face,
1088  const unsigned int subface,
1089  const bool /*face_orientation*/,
1090  const bool face_flip,
1091  const bool /*face_rotation*/,
1092  const RefinementCase<1> &)
1093 {
1094  Assert (face<faces_per_cell, ExcIndexRange(face, 0, faces_per_cell));
1095  Assert (subface<max_children_per_face,
1096  ExcIndexRange(subface, 0, max_children_per_face));
1097 
1098  // always return the child adjacent to the specified
1099  // subface. if the face of a cell is not refined, don't
1100  // throw an assertion but deliver the child adjacent to
1101  // the face nevertheless, i.e. deliver the child of
1102  // this cell adjacent to the subface of a possibly
1103  // refined neighbor. this simplifies setting neighbor
1104  // information in execute_refinement.
1105  static const unsigned int
1106  subcells[2][RefinementCase<2>::isotropic_refinement][faces_per_cell][max_children_per_face] =
1107  {
1108  {
1109  // Normal orientation (face_filp = false)
1110  {{0,0},{1,1},{0,1},{0,1}}, // cut_x
1111  {{0,1},{0,1},{0,0},{1,1}}, // cut_y
1112  {{0,2},{1,3},{0,1},{2,3}} // cut_z
1113  },
1114  {
1115  // Flipped orientation (face_flip = true)
1116  {{0,0},{1,1},{1,0},{1,0}}, // cut_x
1117  {{1,0},{1,0},{0,0},{1,1}}, // cut_y
1118  {{2,0},{3,1},{1,0},{3,2}} // cut_z
1119  }
1120  };
1121 
1122  return subcells[face_flip][ref_case-1][face][subface];
1123 }
1124 
1125 
1126 
1127 template <>
1128 unsigned int
1130  const unsigned int face,
1131  const unsigned int subface,
1132  const bool face_orientation,
1133  const bool face_flip,
1134  const bool face_rotation,
1135  const RefinementCase<2> &face_ref_case)
1136 {
1137  const unsigned int dim = 3;
1138 
1139  Assert (ref_case>RefinementCase<dim-1>::no_refinement, ExcMessage("Cell has no children."));
1140  Assert (face<faces_per_cell, ExcIndexRange(face, 0, faces_per_cell));
1141  Assert (subface<GeometryInfo<dim-1>::n_children(face_ref_case) ||
1142  (subface==0 && face_ref_case==RefinementCase<dim-1>::no_refinement),
1143  ExcIndexRange(subface, 0, GeometryInfo<2>::n_children(face_ref_case)));
1144 
1145  // invalid number used for invalid cases,
1146  // e.g. when the children are more refined at
1147  // a given face than the face itself
1148  static const unsigned int e=invalid_unsigned_int;
1149 
1150  // the whole process of finding a child cell
1151  // at a given subface considering the
1152  // possibly anisotropic refinement cases of
1153  // the cell and the face as well as
1154  // orientation, flip and rotation of the face
1155  // is quite complicated. thus, we break it
1156  // down into several steps.
1157 
1158  // first step: convert the given face refine
1159  // case to a face refine case concerning the
1160  // face in standard orientation (, flip and
1161  // rotation). This only affects cut_x and
1162  // cut_y
1163  static const RefinementCase<dim-1> flip[4]=
1164  {
1165  RefinementCase<dim-1>::no_refinement,
1166  RefinementCase<dim-1>::cut_y,
1167  RefinementCase<dim-1>::cut_x,
1168  RefinementCase<dim-1>::cut_xy
1169  };
1170  // for face_orientation, 'true' is the
1171  // default value whereas for face_rotation,
1172  // 'false' is standard. If
1173  // <tt>face_rotation==face_orientation</tt>,
1174  // then one of them is non-standard and we
1175  // have to swap cut_x and cut_y, otherwise no
1176  // change is necessary.
1177  const RefinementCase<dim-1> std_face_ref = (face_orientation==face_rotation) ? flip[face_ref_case] : face_ref_case;
1178 
1179  // second step: convert the given subface
1180  // index to the one for a standard face
1181  // respecting face_orientation, face_flip and
1182  // face_rotation
1183 
1184  // first index: face_ref_case
1185  // second index: face_orientation
1186  // third index: face_flip
1187  // forth index: face_rotation
1188  // fifth index: subface index
1189  static const unsigned int subface_exchange[4][2][2][2][4]=
1190  {
1191  // no_refinement (subface 0 stays 0,
1192  // all others are invalid)
1193  { { { {0,e,e,e},
1194  {0,e,e,e}
1195  },
1196  { {0,e,e,e},
1197  {0,e,e,e}
1198  }
1199  },
1200  { { {0,e,e,e},
1201  {0,e,e,e}
1202  },
1203  { {0,e,e,e},
1204  {0,e,e,e}
1205  }
1206  }
1207  },
1208  // cut_x (here, if the face is only
1209  // rotated OR only falsely oriented,
1210  // then subface 0 of the non-standard
1211  // face does NOT correspond to one of
1212  // the subfaces of a standard
1213  // face. Thus we indicate the subface
1214  // which is located at the lower left
1215  // corner (the origin of the face's
1216  // local coordinate system) with
1217  // '0'. The rest of this issue is
1218  // taken care of using the above
1219  // conversion to a 'standard face
1220  // refine case')
1221  { { { {0,1,e,e},
1222  {0,1,e,e}
1223  },
1224  { {1,0,e,e},
1225  {1,0,e,e}
1226  }
1227  },
1228  { { {0,1,e,e},
1229  {0,1,e,e}
1230  },
1231  { {1,0,e,e},
1232  {1,0,e,e}
1233  }
1234  }
1235  },
1236  // cut_y (the same applies as for
1237  // cut_x)
1238  { { { {0,1,e,e},
1239  {1,0,e,e}
1240  },
1241  { {1,0,e,e},
1242  {0,1,e,e}
1243  }
1244  },
1245  { { {0,1,e,e},
1246  {1,0,e,e}
1247  },
1248  { {1,0,e,e},
1249  {0,1,e,e}
1250  }
1251  }
1252  },
1253  // cut_xyz: this information is
1254  // identical to the information
1255  // returned by
1256  // GeometryInfo<3>::real_to_standard_face_vertex()
1257  { { { {0,2,1,3}, // face_orientation=false, face_flip=false, face_rotation=false, subfaces 0,1,2,3
1258  {2,3,0,1}
1259  }, // face_orientation=false, face_flip=false, face_rotation=true, subfaces 0,1,2,3
1260  { {3,1,2,0}, // face_orientation=false, face_flip=true, face_rotation=false, subfaces 0,1,2,3
1261  {1,0,3,2}
1262  }
1263  }, // face_orientation=false, face_flip=true, face_rotation=true, subfaces 0,1,2,3
1264  { { {0,1,2,3}, // face_orientation=true, face_flip=false, face_rotation=false, subfaces 0,1,2,3
1265  {1,3,0,2}
1266  }, // face_orientation=true, face_flip=false, face_rotation=true, subfaces 0,1,2,3
1267  { {3,2,1,0}, // face_orientation=true, face_flip=true, face_rotation=false, subfaces 0,1,2,3
1268  {2,0,3,1}
1269  }
1270  }
1271  }
1272  };// face_orientation=true, face_flip=true, face_rotation=true, subfaces 0,1,2,3
1273 
1274  const unsigned int std_subface=subface_exchange
1275  [face_ref_case]
1276  [face_orientation]
1277  [face_flip]
1278  [face_rotation]
1279  [subface];
1280  Assert (std_subface!=e, ExcInternalError());
1281 
1282  // third step: these are the children, which
1283  // can be found at the given subfaces of an
1284  // isotropically refined (standard) face
1285  //
1286  // first index: (refinement_case-1)
1287  // second index: face_index
1288  // third index: subface_index (isotropic refinement)
1289  static const unsigned int
1290  iso_children[RefinementCase<dim>::cut_xyz][faces_per_cell][max_children_per_face] =
1291  {
1292  // cut_x
1293  { {0, 0, 0, 0}, // face 0, subfaces 0,1,2,3
1294  {1, 1, 1, 1}, // face 1, subfaces 0,1,2,3
1295  {0, 0, 1, 1}, // face 2, subfaces 0,1,2,3
1296  {0, 0, 1, 1}, // face 3, subfaces 0,1,2,3
1297  {0, 1, 0, 1}, // face 4, subfaces 0,1,2,3
1298  {0, 1, 0, 1}
1299  }, // face 5, subfaces 0,1,2,3
1300  // cut_y
1301  { {0, 1, 0, 1},
1302  {0, 1, 0, 1},
1303  {0, 0, 0, 0},
1304  {1, 1, 1, 1},
1305  {0, 0, 1, 1},
1306  {0, 0, 1, 1}
1307  },
1308  // cut_xy
1309  { {0, 2, 0, 2},
1310  {1, 3, 1, 3},
1311  {0, 0, 1, 1},
1312  {2, 2, 3, 3},
1313  {0, 1, 2, 3},
1314  {0, 1, 2, 3}
1315  },
1316  // cut_z
1317  { {0, 0, 1, 1},
1318  {0, 0, 1, 1},
1319  {0, 1, 0, 1},
1320  {0, 1, 0, 1},
1321  {0, 0, 0, 0},
1322  {1, 1, 1, 1}
1323  },
1324  // cut_xz
1325  { {0, 0, 1, 1},
1326  {2, 2, 3, 3},
1327  {0, 1, 2, 3},
1328  {0, 1, 2, 3},
1329  {0, 2, 0, 2},
1330  {1, 3, 1, 3}
1331  },
1332  // cut_yz
1333  { {0, 1, 2, 3},
1334  {0, 1, 2, 3},
1335  {0, 2, 0, 2},
1336  {1, 3, 1, 3},
1337  {0, 0, 1, 1},
1338  {2, 2, 3, 3}
1339  },
1340  // cut_xyz
1341  { {0, 2, 4, 6},
1342  {1, 3, 5, 7},
1343  {0, 4, 1, 5},
1344  {2, 6, 3, 7},
1345  {0, 1, 2, 3},
1346  {4, 5, 6, 7}
1347  }
1348  };
1349 
1350  // forth step: check, whether the given face
1351  // refine case is valid for the given cell
1352  // refine case. this is the case, if the
1353  // given face refine case is at least as
1354  // refined as the face is for the given cell
1355  // refine case
1356 
1357  // note, that we are considering standard
1358  // face refinement cases here and thus must
1359  // not pass the given orientation, flip and
1360  // rotation flags
1361  if ((std_face_ref & face_refinement_case(ref_case, face))
1362  == face_refinement_case(ref_case, face))
1363  {
1364  // all is fine. for anisotropic face
1365  // refine cases, select one of the
1366  // isotropic subfaces which neighbors the
1367  // same child
1368 
1369  // first index: (standard) face refine case
1370  // second index: subface index
1371  static const unsigned int equivalent_iso_subface[4][4]=
1372  {
1373  {0,e,e,e}, // no_refinement
1374  {0,3,e,e}, // cut_x
1375  {0,3,e,e}, // cut_y
1376  {0,1,2,3}
1377  }; // cut_xy
1378 
1379  const unsigned int equ_std_subface
1380  =equivalent_iso_subface[std_face_ref][std_subface];
1381  Assert (equ_std_subface!=e, ExcInternalError());
1382 
1383  return iso_children[ref_case-1][face][equ_std_subface];
1384  }
1385  else
1386  {
1387  // the face_ref_case was too coarse,
1388  // throw an error
1389  Assert(false,
1390  ExcMessage("The face RefineCase is too coarse "
1391  "for the given cell RefineCase."));
1392  }
1393  // we only get here in case of an error
1394  return e;
1395 }
1396 
1397 
1398 
1399 template <>
1400 unsigned int
1402  const unsigned int,
1403  const unsigned int,
1404  const bool, const bool, const bool,
1405  const RefinementCase<3> &)
1406 {
1407  Assert(false, ExcNotImplemented());
1408  return invalid_unsigned_int;
1409 }
1410 
1411 
1412 
1413 template <>
1414 unsigned int
1415 GeometryInfo<1>::line_to_cell_vertices (const unsigned int line,
1416  const unsigned int vertex)
1417 {
1418  (void)line;
1419  Assert (line<lines_per_cell, ExcIndexRange(line, 0, lines_per_cell));
1420  Assert (vertex<2, ExcIndexRange(vertex, 0, 2));
1421 
1422  return vertex;
1423 }
1424 
1425 
1426 template <>
1427 unsigned int
1428 GeometryInfo<2>::line_to_cell_vertices (const unsigned int line,
1429  const unsigned int vertex)
1430 {
1431  return child_cell_on_face(RefinementCase<2>::isotropic_refinement, line, vertex);
1432 }
1433 
1434 
1435 
1436 template <>
1437 unsigned int
1438 GeometryInfo<3>::line_to_cell_vertices (const unsigned int line,
1439  const unsigned int vertex)
1440 {
1441  Assert (line<lines_per_cell, ExcIndexRange(line, 0, lines_per_cell));
1442  Assert (vertex<2, ExcIndexRange(vertex, 0, 2));
1443 
1444  static const unsigned
1445  vertices[lines_per_cell][2] = {{0, 2}, // bottom face
1446  {1, 3},
1447  {0, 1},
1448  {2, 3},
1449  {4, 6}, // top face
1450  {5, 7},
1451  {4, 5},
1452  {6, 7},
1453  {0, 4}, // connects of bottom
1454  {1, 5}, // top face
1455  {2, 6},
1456  {3, 7}
1457  };
1458  return vertices[line][vertex];
1459 }
1460 
1461 
1462 template <>
1463 unsigned int
1464 GeometryInfo<4>::line_to_cell_vertices (const unsigned int,
1465  const unsigned int)
1466 {
1467  Assert(false, ExcNotImplemented());
1468  return invalid_unsigned_int;
1469 }
1470 
1471 
1472 template <>
1473 unsigned int
1474 GeometryInfo<1>::face_to_cell_lines (const unsigned int face,
1475  const unsigned int line,
1476  const bool, const bool, const bool)
1477 {
1478  (void)face;
1479  (void)line;
1480  Assert (face+1<faces_per_cell+1, ExcIndexRange(face, 0, faces_per_cell));
1481  Assert (line+1<lines_per_face+1, ExcIndexRange(line, 0, lines_per_face));
1482 
1483  // There is only a single line, so
1484  // it must be this.
1485  return 0;
1486 }
1487 
1488 
1489 
1490 template <>
1491 unsigned int
1492 GeometryInfo<2>::face_to_cell_lines (const unsigned int face,
1493  const unsigned int line,
1494  const bool, const bool, const bool)
1495 {
1496  (void)line;
1497  Assert (face<faces_per_cell, ExcIndexRange(face, 0, faces_per_cell));
1498  Assert (line<lines_per_face, ExcIndexRange(line, 0, lines_per_face));
1499 
1500  // The face is a line itself.
1501  return face;
1502 }
1503 
1504 
1505 
1506 template <>
1507 unsigned int
1508 GeometryInfo<3>::face_to_cell_lines (const unsigned int face,
1509  const unsigned int line,
1510  const bool face_orientation,
1511  const bool face_flip,
1512  const bool face_rotation)
1513 {
1514  Assert (face<faces_per_cell, ExcIndexRange(face, 0, faces_per_cell));
1515  Assert (line<lines_per_face, ExcIndexRange(line, 0, lines_per_face));
1516 
1517  static const unsigned
1518  lines[faces_per_cell][lines_per_face] = {{8,10, 0, 4}, // left face
1519  {9,11, 1, 5}, // right face
1520  {2, 6, 8, 9}, // front face
1521  {3, 7,10,11}, // back face
1522  {0, 1, 2, 3}, // bottom face
1523  {4, 5, 6, 7}
1524  };// top face
1525  return lines[face][real_to_standard_face_line(line,
1526  face_orientation,
1527  face_flip,
1528  face_rotation)];
1529 }
1530 
1531 
1532 
1533 template<int dim>
1534 unsigned int
1536  const unsigned int,
1537  const bool, const bool, const bool)
1538 {
1539  Assert(false, ExcNotImplemented());
1540  return invalid_unsigned_int;
1541 }
1542 
1543 
1544 
1545 template <int dim>
1546 unsigned int
1548  const unsigned int vertex,
1549  const bool face_orientation,
1550  const bool face_flip,
1551  const bool face_rotation)
1552 {
1553  return child_cell_on_face(RefinementCase<dim>::isotropic_refinement, face, vertex,
1554  face_orientation, face_flip, face_rotation);
1555 }
1556 
1557 
1558 
1559 template <int dim>
1560 Point<dim>
1562 {
1563  Point<dim> p = q;
1564  for (unsigned int i=0; i<dim; i++)
1565  if (p[i] < 0.) p[i] = 0.;
1566  else if (p[i] > 1.) p[i] = 1.;
1567 
1568  return p;
1569 }
1570 
1571 
1572 
1573 template <int dim>
1574 double
1576 {
1577  double result = 0.0;
1578 
1579  for (unsigned int i=0; i<dim; i++)
1580  if ((-p[i]) > result)
1581  result = -p[i];
1582  else if ((p[i]-1.) > result)
1583  result = (p[i] - 1.);
1584 
1585  return result;
1586 }
1587 
1588 
1589 
1590 template <int dim>
1591 double
1594  const unsigned int i)
1595 {
1598 
1599  switch (dim)
1600  {
1601  case 1:
1602  {
1603  const double x = xi[0];
1604  switch (i)
1605  {
1606  case 0:
1607  return 1-x;
1608  case 1:
1609  return x;
1610  }
1611  }
1612 
1613  case 2:
1614  {
1615  const double x = xi[0];
1616  const double y = xi[1];
1617  switch (i)
1618  {
1619  case 0:
1620  return (1-x)*(1-y);
1621  case 1:
1622  return x*(1-y);
1623  case 2:
1624  return (1-x)*y;
1625  case 3:
1626  return x*y;
1627  }
1628  }
1629 
1630  case 3:
1631  {
1632  const double x = xi[0];
1633  const double y = xi[1];
1634  const double z = xi[2];
1635  switch (i)
1636  {
1637  case 0:
1638  return (1-x)*(1-y)*(1-z);
1639  case 1:
1640  return x*(1-y)*(1-z);
1641  case 2:
1642  return (1-x)*y*(1-z);
1643  case 3:
1644  return x*y*(1-z);
1645  case 4:
1646  return (1-x)*(1-y)*z;
1647  case 5:
1648  return x*(1-y)*z;
1649  case 6:
1650  return (1-x)*y*z;
1651  case 7:
1652  return x*y*z;
1653  }
1654  }
1655 
1656  default:
1657  Assert (false, ExcNotImplemented());
1658  }
1659  return -1e9;
1660 }
1661 
1662 
1663 
1664 template <>
1668  const unsigned int i)
1669 {
1672 
1673  switch (i)
1674  {
1675  case 0:
1676  return Point<1>(-1.);
1677  case 1:
1678  return Point<1>(1.);
1679  }
1680 
1681  return Point<1>(-1e9);
1682 }
1683 
1684 
1685 
1686 template <>
1690  const unsigned int i)
1691 {
1694 
1695  const double x = xi[0];
1696  const double y = xi[1];
1697  switch (i)
1698  {
1699  case 0:
1700  return Point<2>(-(1-y),-(1-x));
1701  case 1:
1702  return Point<2>(1-y,-x);
1703  case 2:
1704  return Point<2>(-y, 1-x);
1705  case 3:
1706  return Point<2>(y,x);
1707  }
1708  return Point<2> (-1e9, -1e9);
1709 }
1710 
1711 
1712 
1713 template <>
1717  const unsigned int i)
1718 {
1721 
1722  const double x = xi[0];
1723  const double y = xi[1];
1724  const double z = xi[2];
1725  switch (i)
1726  {
1727  case 0:
1728  return Point<3>(-(1-y)*(1-z),
1729  -(1-x)*(1-z),
1730  -(1-x)*(1-y));
1731  case 1:
1732  return Point<3>((1-y)*(1-z),
1733  -x*(1-z),
1734  -x*(1-y));
1735  case 2:
1736  return Point<3>(-y*(1-z),
1737  (1-x)*(1-z),
1738  -(1-x)*y);
1739  case 3:
1740  return Point<3>(y*(1-z),
1741  x*(1-z),
1742  -x*y);
1743  case 4:
1744  return Point<3>(-(1-y)*z,
1745  -(1-x)*z,
1746  (1-x)*(1-y));
1747  case 5:
1748  return Point<3>((1-y)*z,
1749  -x*z,
1750  x*(1-y));
1751  case 6:
1752  return Point<3>(-y*z,
1753  (1-x)*z,
1754  (1-x)*y);
1755  case 7:
1756  return Point<3>(y*z, x*z, x*y);
1757  }
1758 
1759  return Point<3> (-1e9, -1e9, -1e9);
1760 }
1761 
1762 
1763 
1764 template <int dim>
1768  const unsigned int)
1769 {
1770  Assert (false, ExcNotImplemented());
1771  return Tensor<1,dim>();
1772 }
1773 
1774 
1775 
1776 
1777 
1778 namespace internal
1779 {
1780  namespace GeometryInfo
1781  {
1782  // wedge product of a single
1783  // vector in 2d: we just have to
1784  // rotate it by 90 degrees to the
1785  // right
1786  inline
1787  Tensor<1,2>
1788  wedge_product (const Tensor<1,2> (&derivative)[1])
1789  {
1790  Tensor<1,2> result;
1791  result[0] = derivative[0][1];
1792  result[1] = -derivative[0][0];
1793 
1794  return result;
1795  }
1796 
1797 
1798  // wedge product of 2 vectors in
1799  // 3d is the cross product
1800  inline
1801  Tensor<1,3>
1802  wedge_product (const Tensor<1,3> (&derivative)[2])
1803  {
1804  return cross_product_3d (derivative[0], derivative[1]);
1805  }
1806 
1807 
1808  // wedge product of dim vectors
1809  // in dim-d: that's the
1810  // determinant of the matrix
1811  template <int dim>
1812  inline
1814  wedge_product (const Tensor<1,dim> (&derivative)[dim])
1815  {
1816  Tensor<2,dim> jacobian;
1817  for (unsigned int i=0; i<dim; ++i)
1818  jacobian[i] = derivative[i];
1819 
1820  return determinant (jacobian);
1821  }
1822  }
1823 }
1824 
1825 
1826 template <int dim>
1827 template <int spacedim>
1828 void
1830 alternating_form_at_vertices
1831 #ifndef DEAL_II_CONSTEXPR_BUG
1832 (const Point<spacedim> (&vertices)[vertices_per_cell],
1833  Tensor<spacedim-dim,spacedim> (&forms)[vertices_per_cell])
1834 #else
1835 (const Point<spacedim> *vertices,
1837 #endif
1838 {
1839  // for each of the vertices,
1840  // compute the alternating form
1841  // of the mapped unit
1842  // vectors. consider for
1843  // example the case of a quad
1844  // in spacedim==3: to do so, we
1845  // need to see how the
1846  // infinitesimal vectors
1847  // (d\xi_1,0) and (0,d\xi_2)
1848  // are transformed into
1849  // spacedim-dimensional space
1850  // and then form their cross
1851  // product (i.e. the wedge product
1852  // of two vectors). to this end, note
1853  // that
1854  // \vec x = sum_i \vec v_i phi_i(\vec xi)
1855  // so the transformed vectors are
1856  // [x(\xi+(d\xi_1,0))-x(\xi)]/d\xi_1
1857  // and
1858  // [x(\xi+(0,d\xi_2))-x(\xi)]/d\xi_2
1859  // which boils down to the columns
1860  // of the 3x2 matrix \grad_\xi x(\xi)
1861  //
1862  // a similar reasoning would
1863  // hold for all dim,spacedim
1864  // pairs -- we only have to
1865  // compute the wedge product of
1866  // the columns of the
1867  // derivatives
1868  for (unsigned int i=0; i<vertices_per_cell; ++i)
1869  {
1870  Tensor<1,spacedim> derivatives[dim];
1871 
1872  for (unsigned int j=0; j<vertices_per_cell; ++j)
1873  {
1874  const Tensor<1,dim> grad_phi_j
1875  = d_linear_shape_function_gradient (unit_cell_vertex(i),
1876  j);
1877  for (unsigned int l=0; l<dim; ++l)
1878  derivatives[l] += vertices[j] * grad_phi_j[l];
1879  }
1880 
1881  forms[i] = internal::GeometryInfo::wedge_product (derivatives);
1882  }
1883 }
1884 
1885 
1886 template struct GeometryInfo<1>;
1887 template struct GeometryInfo<2>;
1888 template struct GeometryInfo<3>;
1889 template struct GeometryInfo<4>;
1890 
1891 template
1892 void
1894 alternating_form_at_vertices
1895 #ifndef DEAL_II_CONSTEXPR_BUG
1896 (const Point<1> (&)[vertices_per_cell],
1897  Tensor<1-1,1> (&)[vertices_per_cell])
1898 #else
1899 (const Point<1> *, Tensor<1-1,1> *)
1900 #endif
1901 ;
1902 
1903 template
1904 void
1906 alternating_form_at_vertices
1907 #ifndef DEAL_II_CONSTEXPR_BUG
1908 (const Point<2> (&)[vertices_per_cell],
1909  Tensor<2-1,2> (&)[vertices_per_cell])
1910 #else
1911 (const Point<2> *, Tensor<2-1,2> *)
1912 #endif
1913 ;
1914 
1915 template
1916 void
1918 alternating_form_at_vertices
1919 #ifndef DEAL_II_CONSTEXPR_BUG
1920 (const Point<2> (&vertices)[vertices_per_cell],
1921  Tensor<2-2,2> (&forms)[vertices_per_cell])
1922 #else
1923 (const Point<2> *, Tensor<2-2,2> *)
1924 #endif
1925 ;
1926 
1927 template
1928 void
1930 alternating_form_at_vertices
1931 #ifndef DEAL_II_CONSTEXPR_BUG
1932 (const Point<3> (&vertices)[vertices_per_cell],
1933  Tensor<3-2,3> (&forms)[vertices_per_cell])
1934 #else
1935 (const Point<3> *, Tensor<3-2,3> *)
1936 #endif
1937 ;
1938 
1939 
1940 template
1941 void
1943 alternating_form_at_vertices
1944 #ifndef DEAL_II_CONSTEXPR_BUG
1945 (const Point<3> (&vertices)[vertices_per_cell],
1946  Tensor<3-3,3> (&forms)[vertices_per_cell])
1947 #else
1948 (const Point<3> *, Tensor<3-3,3> *)
1949 #endif
1950 ;
1951 
1952 
1953 DEAL_II_NAMESPACE_CLOSE
static unsigned int standard_to_real_face_vertex(const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static const unsigned int invalid_unsigned_int
Definition: types.h:170
static RefinementCase< 1 > line_refinement_case(const RefinementCase< dim > &cell_refinement_case, const unsigned int line_no)
static unsigned int real_to_standard_face_line(const unsigned int line, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
static unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static unsigned int line_to_cell_vertices(const unsigned int line, const unsigned int vertex)
static RefinementCase< dim-1 > face_refinement_case(const RefinementCase< dim > &cell_refinement_case, const unsigned int face_no, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
static unsigned int face_to_cell_lines(const unsigned int face, const unsigned int line, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static unsigned int child_cell_on_face(const RefinementCase< dim > &ref_case, const unsigned int face, const unsigned int subface, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const RefinementCase< dim-1 > &face_refinement_case=RefinementCase< dim-1 >::isotropic_refinement)
#define Assert(cond, exc)
Definition: exceptions.h:313
static double distance_to_unit_cell(const Point< dim > &p)
static unsigned int standard_to_real_face_line(const unsigned int line, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static RefinementCase< dim > min_cell_refinement_case_for_line_refinement(const unsigned int line_no)
Definition: mpi.h:41
static unsigned int n_subfaces(const internal::SubfaceCase< dim > &subface_case)
static unsigned int real_to_standard_face_vertex(const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static Point< dim > project_to_unit_cell(const Point< dim > &p)
static ::ExceptionBase & ExcNotImplemented()
static Tensor< 1, dim > d_linear_shape_function_gradient(const Point< dim > &xi, const unsigned int i)
static double subface_ratio(const internal::SubfaceCase< dim > &subface_case, const unsigned int subface_no)
static RefinementCase< dim > min_cell_refinement_case_for_face_refinement(const RefinementCase< dim-1 > &face_refinement_case, const unsigned int face_no, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
static ::ExceptionBase & ExcInternalError()