]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
remove unneeded local indices from other simple assemblers
authorkanschat <kanschat@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 7 Mar 2013 23:20:54 +0000 (23:20 +0000)
committerkanschat <kanschat@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 7 Mar 2013 23:20:54 +0000 (23:20 +0000)
git-svn-id: https://svn.dealii.org/trunk@28806 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/meshworker/simple.h

index 8142d182dadcf210091b601f714b214cb738866a..3858b436ad75124c0100950c42f4e5fba31cfa44 100644 (file)
@@ -63,6 +63,9 @@ namespace MeshWorker
       void initialize(const ConstraintMatrix &constraints);
 
       /**
+       * @deprecated This function is of no effect. Only the block info
+       * structure in DoFInfo is being used.
+       *
        * Store information on the local block structure. If the
        * assembler is inititialized with this function,
        * initialize_info() will generate one local matrix for each
@@ -72,14 +75,9 @@ namespace MeshWorker
        * In spite of using local block structure, all blocks will be
        * enteres into the same global matrix, disregarding any global
        * block structure.
-       *
-       * @note The argument of this function will be copied into the
-       * member object #local_indices. Thus, every subsequent change
-       * in the block structure must be initialzied or will not be
-       * used by the assembler.
        */
 
-      void initialize_local_blocks(const BlockIndices &local_indices);
+      void initialize_local_blocks(const BlockIndices &);
 
       /**
        * Initialize the local data in the DoFInfo object used later
@@ -108,12 +106,6 @@ namespace MeshWorker
       void assemble(const DOFINFO &info1,
                     const DOFINFO &info2);
     private:
-      /**
-       * The object containing the local block structure. Set by
-       * initialize_local_blocks() and used by assembling functions.
-       */
-      BlockIndices local_indices;
-
       /**
        * The global residal vectors filled by assemble().
        */
@@ -144,11 +136,9 @@ namespace MeshWorker
      * by default, the obvious choice of taking a single local matrix with
      * dimensions equal to the number of degrees of freedom of the
      * cell. Alternatively, a local block structure can be initialized
-     * with initialize_local_blocks(). After this, the local data will be
+     * in DoFInfo. After this, the local data will be
      * arranged as an array of n by n FullMatrix blocks, which are
-     * ordered lexicographically in DoFInfo. Note that
-     * initialize_local_blocks() has to be called before initialize_info()
-     * to take the desired effect.
+     * ordered lexicographically in DoFInfo.
      *
      * @ingroup MeshWorker
      * @author Guido Kanschat, 2009
@@ -189,14 +179,9 @@ namespace MeshWorker
        * In spite of using local block structure, all blocks will be
        * enteres into the same global matrix, disregarding any global
        * block structure.
-       *
-       * @note The argument of this function will be copied into the
-       * member object #local_indices. Thus, every subsequent change
-       * in the block structure must be initialzied or will not be
-       * used by the assembler.
        */
 
-      void initialize_local_blocks(const BlockIndices &local_indices);
+      void initialize_local_blocks(const BlockIndices &);
 
       /**
        * Initialize the local data in the DoFInfo object used later
@@ -287,6 +272,9 @@ namespace MeshWorker
       void initialize(const MGConstrainedDoFs &mg_constrained_dofs);
 
       /**
+       * @deprecated This function is of no effect. Only the block info
+       * structure in DoFInfo is being used.
+       *
        * Store information on the local block structure. If the
        * assembler is inititialized with this function,
        * initialize_info() will generate one local matrix for each
@@ -296,13 +284,8 @@ namespace MeshWorker
        * In spite of using local block structure, all blocks will be
        * enteres into the same global matrix, disregarding any global
        * block structure.
-       *
-       * @note The argument of this function will be copied into the
-       * member object #local_indices. Thus, every subsequent change
-       * in the block structure must be initialzied or will not be
-       * used by the assembler.
        */
-      void initialize_local_blocks(const BlockIndices &local_indices);
+      void initialize_local_blocks(const BlockIndices &);
 
       /**
        * Initialize the matrices #flux_up and #flux_down used for
@@ -431,12 +414,6 @@ namespace MeshWorker
        */
       SmartPointer<const MGConstrainedDoFs,MGMatrixSimple<MATRIX> > mg_constrained_dofs;
 
-      /**
-       * The object containing the local block structure. Set by
-       * initialize_local_blocks() and used by assembling functions.
-       */
-      BlockIndices local_indices;
-
       /**
        * The smallest positive number that will be entered into the
        * global matrix. All smaller absolute values will be treated as
@@ -531,10 +508,8 @@ namespace MeshWorker
 
     template <class MATRIX>
     inline void
-    ResidualSimple<MATRIX>::initialize_local_blocks(const BlockIndices &b)
-    {
-      local_indices = b;
-    }
+    ResidualSimple<MATRIX>::initialize_local_blocks(const BlockIndices &)
+    {}
     
     
     template <class VECTOR>
@@ -560,7 +535,7 @@ namespace MeshWorker
             }
           else
             {
-             if (local_indices.size() == 0)
+             if (info.indices_by_block.size() == 0)
                constraints->distribute_local_to_global(info.vector(k).block(0), info.indices, (*residuals(k)));
              else
                for (unsigned int i=0;i != info.vector(k).n_blocks();++i)
@@ -586,14 +561,14 @@ namespace MeshWorker
             }
           else
             {
-             if (local_indices.size() == 0)
+             if (info1.indices_by_block.size() == 0 && info2.indices_by_block.size() == 0)
                {
                  constraints->distribute_local_to_global
                    (info1.vector(k).block(0), info1.indices, (*residuals(k)));
                  constraints->distribute_local_to_global
                    (info2.vector(k).block(0), info2.indices, (*residuals(k)));
                }
-             else
+             else if (info1.indices_by_block.size() != 0 && info2.indices_by_block.size() != 0)
                {
                  for (unsigned int i=0;i<info1.vector(k).n_blocks();++i)
                    {
@@ -603,6 +578,10 @@ namespace MeshWorker
                        (info2.vector(k).block(i), info2.indices_by_block[i], (*residuals(k)));
                    }
                }
+             else
+               {
+                 Assert(false, ExcNotImplemented());
+               }
             }
         }
     }
@@ -774,10 +753,8 @@ namespace MeshWorker
 
     template <class MATRIX>
     inline void
-    MGMatrixSimple<MATRIX>::initialize_local_blocks(const BlockIndices &b)
-    {
-      local_indices = b;
-    }
+    MGMatrixSimple<MATRIX>::initialize_local_blocks(const BlockIndices &)
+    {}
 
 
     template <class MATRIX>
@@ -805,7 +782,7 @@ namespace MeshWorker
     inline void
     MGMatrixSimple<MATRIX>::initialize_info(DOFINFO &info, bool face) const
     {
-      const unsigned int n = local_indices.size();
+      const unsigned int n = info.indices_by_block.size();
 
       if (n == 0)
         info.initialize_matrices(1, face);
@@ -1067,7 +1044,7 @@ namespace MeshWorker
       Assert(info.level_cell, ExcMessage("Cell must access level dofs"));
       const unsigned int level = info.cell->level();
 
-      if (local_indices.size() == 0)
+      if (info.indices_by_block.size() == 0)
         {
           assemble((*matrix)[level], info.matrix(0,false).matrix,
                    info.indices, info.indices, level);
@@ -1110,7 +1087,7 @@ namespace MeshWorker
       const unsigned int level1 = info1.cell->level();
       const unsigned int level2 = info2.cell->level();
 
-      if (local_indices.size() == 0)
+      if (info1.indices_by_block.size() == 0)
         {
           if (level1 == level2)
             {

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.