]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Implementation of a new FEFaceValuesBase::get_face_index() function which returns...
authorRalf Hartmann <Ralf.Hartmann@dlr.de>
Fri, 15 Feb 2008 15:02:33 +0000 (15:02 +0000)
committerRalf Hartmann <Ralf.Hartmann@dlr.de>
Fri, 15 Feb 2008 15:02:33 +0000 (15:02 +0000)
git-svn-id: https://svn.dealii.org/trunk@15743 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/fe/fe_values.h
deal.II/deal.II/source/fe/fe_values.cc

index b675d6e4d62df57761c5fbad153965e80421e810..f2e9e3421065b90fc5bdea735dd6ef2bed663b32 100644 (file)
@@ -2077,6 +2077,13 @@ class FEFaceValuesBase : public FEValuesBase<dim>
                                      */
     const std::vector<Tensor<1,dim> > & get_boundary_forms () const;
 
+                                    /**
+                                     * Return the index of the face
+                                     * selected the last time the
+                                     * reinit() function was called.
+                                     */
+    unsigned int get_face_index() const;
+
                                     /**
                                      * Return a reference to the copy
                                      * of the quadrature formula
@@ -2092,6 +2099,14 @@ class FEFaceValuesBase : public FEValuesBase<dim>
     unsigned int memory_consumption () const;
     
   protected:
+
+                                    /**
+                                     * Index of the face selected the
+                                     * last time the reinit()
+                                     * function was called.
+                                     */
+    unsigned int present_face_index;
+    
                                     /**
                                      * Store a copy of the quadrature
                                      * formula here.
@@ -2964,6 +2979,14 @@ FEFaceValuesBase<dim>::normal_vector (const unsigned int i) const
 }
 
 
+template <int dim>
+inline
+unsigned int
+FEFaceValuesBase<dim>::get_face_index () const
+{
+  return present_face_index;
+}
+
 
 /*------------------------ Inline functions: FE*FaceValues --------------------*/
 
index 8408e3574c881648912bf0af4f38df40e7d3e165..866c10c757ec10d93ee5b80d609d7a8615015e3b 100644 (file)
@@ -1540,6 +1540,9 @@ void FEFaceValues<dim>::reinit (const typename DoFHandler<dim>::cell_iterator &c
     (new typename FEValuesBase<dim>::template
      CellIterator<typename DoFHandler<dim>::cell_iterator> (cell));
 
+                                  // set the present face index
+  this->present_face_index=cell->face_index(face_no);
+  
                                    // this was the part of the work
                                    // that is dependent on the actual
                                    // data type of the iterator. now
@@ -1575,6 +1578,9 @@ void FEFaceValues<dim>::reinit (const typename hp::DoFHandler<dim>::cell_iterato
     (new typename FEValuesBase<dim>::template
      CellIterator<typename hp::DoFHandler<dim>::cell_iterator> (cell));
 
+                                  // set the present face index
+  this->present_face_index=cell->face_index(face_no);
+
                                    // this was the part of the work
                                    // that is dependent on the actual
                                    // data type of the iterator. now
@@ -1608,6 +1614,9 @@ void FEFaceValues<dim>::reinit (const typename MGDoFHandler<dim>::cell_iterator
     (new typename FEValuesBase<dim>::template
      CellIterator<typename MGDoFHandler<dim>::cell_iterator> (cell));
 
+                                  // set the present face index
+  this->present_face_index=cell->face_index(face_no);
+
                                    // this was the part of the work
                                    // that is dependent on the actual
                                    // data type of the iterator. now
@@ -1632,6 +1641,9 @@ void FEFaceValues<dim>::reinit (const typename Triangulation<dim>::cell_iterator
   this->present_cell.reset 
     (new typename FEValuesBase<dim>::TriaCellIterator (cell));
 
+                                  // set the present face index
+  this->present_face_index=cell->face_index(face_no);
+
                                    // this was the part of the work
                                    // that is dependent on the actual
                                    // data type of the iterator. now
@@ -1765,6 +1777,13 @@ void FESubfaceValues<dim>::reinit (const typename DoFHandler<dim>::cell_iterator
           ExcMessage ("You can't use subface data for cells that are "
                       "already refined. Iterate over their children "
                       "instead in these cases."));
+
+                                  // set the present face index
+  unsigned int real_subface_no=subface_no;
+  if (dim==3)
+    real_subface_no=GeometryInfo<dim>::standard_to_real_face_vertex(
+      subface_no, cell->face_orientation(face_no), cell->face_flip(face_no), cell->face_rotation(face_no));
+  this->present_face_index=cell->face(face_no)->child_index(real_subface_no);
   
                                    // set new cell. auto_ptr will take
                                    // care that old object gets
@@ -1816,6 +1835,13 @@ void FESubfaceValues<dim>::reinit (const typename hp::DoFHandler<dim>::cell_iter
     (new typename FEValuesBase<dim>::template
      CellIterator<typename hp::DoFHandler<dim>::cell_iterator> (cell));
 
+                                  // set the present face index
+  unsigned int real_subface_no=subface_no;
+  if (dim==3)
+    real_subface_no=GeometryInfo<dim>::standard_to_real_face_vertex(
+      subface_no, cell->face_orientation(face_no), cell->face_flip(face_no), cell->face_rotation(face_no));
+  this->present_face_index=cell->face(face_no)->child_index(real_subface_no);
+
                                    // this was the part of the work
                                    // that is dependent on the actual
                                    // data type of the iterator. now
@@ -1851,6 +1877,13 @@ void FESubfaceValues<dim>::reinit (const typename MGDoFHandler<dim>::cell_iterat
     (new typename FEValuesBase<dim>::template
      CellIterator<typename MGDoFHandler<dim>::cell_iterator> (cell));
 
+                                  // set the present face index
+  unsigned int real_subface_no=subface_no;
+  if (dim==3)
+    real_subface_no=GeometryInfo<dim>::standard_to_real_face_vertex(
+      subface_no, cell->face_orientation(face_no), cell->face_flip(face_no), cell->face_rotation(face_no));
+  this->present_face_index=cell->face(face_no)->child_index(real_subface_no);
+
                                    // this was the part of the work
                                    // that is dependent on the actual
                                    // data type of the iterator. now
@@ -1879,6 +1912,13 @@ void FESubfaceValues<dim>::reinit (const typename Triangulation<dim>::cell_itera
   this->present_cell.reset 
     (new typename FEValuesBase<dim>::TriaCellIterator (cell));
 
+                                  // set the present face index
+  unsigned int real_subface_no=subface_no;
+  if (dim==3)
+    real_subface_no=GeometryInfo<dim>::standard_to_real_face_vertex(
+      subface_no, cell->face_orientation(face_no), cell->face_flip(face_no), cell->face_rotation(face_no));
+  this->present_face_index=cell->face(face_no)->child_index(real_subface_no);
+
                                    // this was the part of the work
                                    // that is dependent on the actual
                                    // data type of the iterator. now

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.