]> https://gitweb.dealii.org/ - dealii.git/commitdiff
update namespace and documentation
authorGuido Kanschat <dr.guido.kanschat@gmail.com>
Fri, 16 Nov 2012 14:37:57 +0000 (14:37 +0000)
committerGuido Kanschat <dr.guido.kanschat@gmail.com>
Fri, 16 Nov 2012 14:37:57 +0000 (14:37 +0000)
git-svn-id: https://svn.dealii.org/trunk@27546 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/dofs/dof_renumbering.h
deal.II/include/deal.II/meshworker/integration_info.h
deal.II/include/deal.II/meshworker/vector_info.h

index 341d1da45e2a19ddc2c1496d3a1b14e92b240e46..1cbc12da7fc771db9c4a7d78bde6bf3af8f0a98d 100644 (file)
@@ -938,8 +938,7 @@ namespace DoFRenumbering
   hierarchical (DoFHandler<dim> &dof_handler);
 
                                    /**
-                                    * Cell-wise renumbering for DG
-                                    * elements.  This function takes
+                                    * Cell-wise renumbering. This function takes
                                     * the ordered set of cells in
                                     * <tt>cell_order</tt>, and makes
                                     * sure that all degrees of
@@ -949,12 +948,6 @@ namespace DoFRenumbering
                                     * lower index. The order inside
                                     * a cell block will be the same
                                     * as before this renumbering.
-                                    *
-                                    * This function only works with
-                                    * Discontinuous Galerkin Finite
-                                    * Elements, i.e. all degrees of
-                                    * freedom have to be associated
-                                    * with the interior of the cell.
                                     */
   template <class DH>
   void
@@ -1264,7 +1257,7 @@ namespace DoFRenumbering
                 const unsigned int level,
                 const Point<dim> &center,
                 const bool counter = false);
-
+  
                                    /**
                                     * Computes the renumbering
                                     * vector needed by the
index 6e715e36ec292a7842ded6d8480123d7692c38bb..c321c2b648461842d8e9c36b0b04296ef0174017 100644 (file)
@@ -30,18 +30,11 @@ namespace MeshWorker
  * Class for objects handed to local integration functions.
  *
  * Objects of this class contain one or more objects of type FEValues,
- * FEFaceValues or FESubfacevalues to be used in local
+ * FEFaceValues or FESubfaceValues to be used in local
  * integration. They are stored in an array of pointers to the base
- * classes FEValuesBase for cells and FEFaceValuesBase for faces and
- * subfaces, respectively. The template parameter VECTOR allows the
+ * classes FEValuesBase. The template parameter VECTOR allows the
  * use of different data types for the global system.
  *
- * The @p FEVALUESBASE template parameter should be either
- * FEValuesBase or FEFaceValuesBase, depending on whether the object
- * is used to integrate over cells or faces. The actual type of @p
- * FEVALUES object is fixed in the constructor and only used to
- * initialize the pointers in #fevalv.
- *
  * Additionally, this function contains space to store the values of
  * finite element functions stored in #global_data in the
  * quadrature points. These vectors are initialized automatically on
index 5f5920e5edfda58edcdde1a31d2debff2b815e55..4b3257ad4ed2d02f7eb663fd1fb377a7d41baa34 100644 (file)
@@ -1,7 +1,7 @@
 //---------------------------------------------------------------------------
 //    $Id: mesh_worker_info.h 23936 2011-07-09 17:02:27Z kanschat $
 //
-//    Copyright (C) 2012 by the deal.II authors
+//    Copyright (C) 2006, 2007, 2008, 2009, 2010, 2011, 2012 by the deal.II authors
 //
 //    This file is subject to QPL and may not be  distributed
 //    without copyright and license information. Please refer
 #include <deal.II/lac/block_indices.h>
 #include <deal.II/numerics/mesh_worker_info.h>
 
-//TODO[GK]: Remove the using directive
-using namespace dealii;
+DEAL_II_NAMESPACE_OPEN
 
-//TODO[GK]: Provide any kind of documentation for this class and the one below...
-class VectorInfo
+namespace MeshWorker
 {
-  public:
-                                    /**
-                                     * Initialize the data
-                                     * vector and cache the
-                                     * selector.
-                                     */
-    void initialize_data(const NamedData<BlockVector<double>*>&data);
-
-    template<int dim, int spacedim>
-    void reinit(const MeshWorker::DoFInfo<dim, spacedim>& i);
-
-    BlockVector<double>& operator() (unsigned int i);
+/**
+ * This class is a simple object that can be used in
+ * MeshWorker::loop(). For all vectors provided in the argument to
+ * initialize_data(), this class generates the local data vector on a
+ * cell specified by reinit().
+ *
+ * This objects of class as part of VectorInfoBox conform to the
+ * interface for INFOBOX laid out in the documentation of MeshWorker.
+ *
+ * See IntegrationInfo for an alternative.
+ *
+ * @author Guido Kanschat
+ * @date 2012
+ */
+  class VectorInfo
+  {
+    public:
+                                      /**
+                                       * Initialize the data
+                                       * vector and cache the
+                                       * selector.
+                                       */
+      void initialize_data(const NamedData<BlockVector<double>*>&data);
 
-  private:
-    std::vector<BlockVector<double> > local_data;
+                                      /**
+                                       * Reinitialize the local data
+                                       * vectors to represent the
+                                       * information on the given cell.
+                                       */
+      template<int dim, int spacedim>
+      void reinit(const MeshWorker::DoFInfo<dim, spacedim>& cell);
+
+                                      /**
+                                       * Return local data vector
+                                       * <i>i</i>, which selects the
+                                       * local data from vector
+                                       * <i>i</i> of the
+                                       * <code>data</code> object
+                                       * used in initialize_data().
+                                       */
+      BlockVector<double>& operator() (unsigned int i);
+      
+    private:
+      std::vector<BlockVector<double> > local_data;
                                       /**
                                        * The global data vector
                                        * used to compute function
                                        * values in quadrature
                                        * points.
                                        */
-    SmartPointer<const NamedData<BlockVector<double>*> > global_data;
-};
-
-
-class VectorInfoBox
- {
-   public:
-     typedef VectorInfo CellInfo;
-
-     void initialize_data(const NamedData<BlockVector<double>*>&data);
-
-     template <int dim, class DOFINFO>
-     void post_cell(const MeshWorker::DoFInfoBox<dim, DOFINFO>&)
-       {}
-
-     template <int dim, class DOFINFO>
-     void post_faces(const MeshWorker::DoFInfoBox<dim, DOFINFO>&)
-       {}
-
-     VectorInfo cell;
-     VectorInfo boundary;
-     VectorInfo face;
-     VectorInfo subface;
-     VectorInfo neighbor;
- };
-
-
-inline
-void
-VectorInfo::initialize_data(const NamedData<BlockVector<double>*>&data)
-{
-  global_data = &data;
-  local_data.resize(global_data->size());
-}
-
-
-template<int dim, int spacedim>
-inline
-void
-VectorInfo::reinit(const MeshWorker::DoFInfo<dim, spacedim>& i)
-{
-  const NamedData<BlockVector<double>*>& gd = *global_data;
-
-  for (unsigned int k=0;k<local_data.size();++k)
-    {
-      const BlockVector<double>& v = *gd(k);
-
-      local_data[k].reinit(i.block_info->local());
-      for (unsigned int j=0;j<local_data[k].size();++j)
-       local_data[k](j) = v(i.indices[j]);
-    }
-}
+      SmartPointer<const NamedData<BlockVector<double>*> > global_data;    
+  };
+
+
+/**
+ * A class conforming to the INFOBOX interface in the documentation of
+ * MeshWorker. This class provides local operators in
+ * MeshWorker::loop() with the local values of the degrees of freedom
+ * of global vectors, for instance for calculations based on cell
+ * matrices.
+ *
+ * For an alternative providing the values in quadrature points see
+ * IntegrationInfoBox.
+ *
+ * @author Guido Kanschat
+ * @date 2012
+ */
+  class VectorInfoBox
+  {
+    public:
+      typedef VectorInfo CellInfo;
+      
+      void initialize_data(const NamedData<BlockVector<double>*>&data);
+      
+      template <int dim, class DOFINFO>
+      void post_cell(const MeshWorker::DoFInfoBox<dim, DOFINFO>&)
+       {}
+      
+      template <int dim, class DOFINFO>
+      void post_faces(const MeshWorker::DoFInfoBox<dim, DOFINFO>&)
+       {}
+      
+      VectorInfo cell;
+      VectorInfo boundary;
+      VectorInfo face;
+      VectorInfo subface;
+      VectorInfo neighbor;
+  };
+  
+  
+  inline
+  void
+  VectorInfo::initialize_data(const NamedData<BlockVector<double>*>&data)
+  {
+    global_data = &data;
+    local_data.resize(global_data->size());
+  }
+
+
+  template<int dim, int spacedim>
+  inline
+  void
+  VectorInfo::reinit(const MeshWorker::DoFInfo<dim, spacedim>& i)
+  {
+    const NamedData<BlockVector<double>*>& gd = *global_data;
+
+    for (unsigned int k=0;k<local_data.size();++k)
+      {
+       const BlockVector<double>& v = *gd(k);
+
+       local_data[k].reinit(i.block_info->local());
+       for (unsigned int j=0;j<local_data[k].size();++j)
+         local_data[k](j) = v(i.indices[j]);
+      }
+  }
+
+
+  inline
+  BlockVector<double>&
+  VectorInfo::operator() (unsigned int i)
+  {
+    AssertIndexRange(i, local_data.size());
+    return local_data[i];
+  }
+
+
+  inline
+  void
+  VectorInfoBox::initialize_data(const NamedData<BlockVector<double>*>&data)
+  {
+    cell.initialize_data(data);
+    boundary.initialize_data(data);
+    face.initialize_data(data);
+    subface.initialize_data(data);
+    neighbor.initialize_data(data);
+  }
 
-
-inline
-BlockVector<double>&
-VectorInfo::operator() (unsigned int i)
-{
-  AssertIndexRange(i, local_data.size());
-  return local_data[i];
 }
 
-
-inline
-void
-VectorInfoBox::initialize_data(const NamedData<BlockVector<double>*>&data)
-{
-  cell.initialize_data(data);
-  boundary.initialize_data(data);
-  face.initialize_data(data);
-  subface.initialize_data(data);
-  neighbor.initialize_data(data);
-}
-
-
+DEAL_II_NAMESPACE_CLOSE
 
 #endif

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.