]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add functionality to FETools::compute_block_renumbering
authorGuido Kanschat <dr.guido.kanschat@gmail.com>
Wed, 26 Aug 2009 18:08:58 +0000 (18:08 +0000)
committerGuido Kanschat <dr.guido.kanschat@gmail.com>
Wed, 26 Aug 2009 18:08:58 +0000 (18:08 +0000)
git-svn-id: https://svn.dealii.org/trunk@19346 0785d39b-7218-0410-832d-ea1e28bc413d

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

index 6fae20c4f43bb7f253dcf802727af84e4db42273..dccdc7c1dc9cd7c706beb8fa1df664a89c94d97a 100644 (file)
@@ -163,10 +163,17 @@ class FETools
                                     /**
                                      * Compute the vector required to
                                      * renumber the dofs of a cell by
-                                     * block. Furthermore,
-                                     * compute the vector storing the
-                                     * start indices of each local block
-                                     * vector.
+                                     * block. Furthermore, compute
+                                     * the vector storing either the
+                                     * start indices or the size of
+                                     * each local block vector.
+                                     *
+                                     * If the @p bool parameter is
+                                     * true, @p block_data is filled
+                                     * with the start indices of each
+                                     * local block. If it is false,
+                                     * then the block sizes are
+                                     * returned.
                                      *
                                      * @todo Which way does this
                                      * vector map the numbers?
@@ -175,7 +182,8 @@ class FETools
     static void compute_block_renumbering (
       const FiniteElement<dim,spacedim>&  fe,
       std::vector<unsigned int>& renumbering,
-      std::vector<unsigned int>& start_indices);
+      std::vector<unsigned int>& block_data,
+      bool return_start_indices = true);
     
                                     /**
                                      * @name Generation of local matrices
index dd34c061ca04a6db49106b8f62c4a0e9df7aa242..d9ecc24fa88670b3f7d9b8eac63416da6f2bb69d 100644 (file)
@@ -289,13 +289,14 @@ template<int dim, int spacedim>
 void FETools::compute_block_renumbering (
   const FiniteElement<dim,spacedim>& element,
   std::vector<unsigned int>& renumbering,
-  std::vector<unsigned int>& start_indices)
+  std::vector<unsigned int>& block_data,
+  bool return_start_indices)
 {
   Assert(renumbering.size() == element.dofs_per_cell,
         ExcDimensionMismatch(renumbering.size(),
                              element.dofs_per_cell));
-  Assert(start_indices.size() == element.n_blocks(),
-        ExcDimensionMismatch(start_indices.size(),
+  Assert(block_data.size() == element.n_blocks(),
+        ExcDimensionMismatch(block_data.size(),
                              element.n_blocks()));
   
   unsigned int k=0;
@@ -303,10 +304,24 @@ void FETools::compute_block_renumbering (
   for (unsigned int b=0;b<element.n_base_elements();++b)
     for (unsigned int m=0;m<element.element_multiplicity(b);++m)
       {
-       start_indices[i++] = k;
+       block_data[i++] = (return_start_indices)
+                            ? k
+                            : (element.base_element(b).n_dofs_per_cell());
        k += element.base_element(b).n_dofs_per_cell();
       }
   Assert (i == element.n_blocks(), ExcInternalError());
+  
+  std::vector<unsigned int> start_indices(block_data.size());
+  k = 0;
+  for (unsigned int i=0;i<block_data.size();++i)
+    if (return_start_indices)
+      start_indices[i] = block_data[i];
+    else
+      {
+       start_indices[i] = k;
+       k += block_data[i];
+      }
+  
 //TODO:[GK] This does not work for a single RT  
   for (unsigned int i=0;i<element.dofs_per_cell;++i)
     {
@@ -2037,7 +2052,7 @@ void FETools::compute_component_wise(
 template
 void FETools::compute_block_renumbering (
   const FiniteElement<deal_II_dimension>& element,
-  std::vector<unsigned int>&, std::vector<unsigned int>&_indices);
+  std::vector<unsigned int>&, std::vector<unsigned int>&_indices, bool);
 template
 void FETools::get_interpolation_matrix<deal_II_dimension>
 (const FiniteElement<deal_II_dimension> &,

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.