]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
avoid expensive computation in MGTools::reinit_vector
authorguido <guido@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 27 Aug 2004 12:38:27 +0000 (12:38 +0000)
committerguido <guido@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 27 Aug 2004 12:38:27 +0000 (12:38 +0000)
git-svn-id: https://svn.dealii.org/trunk@9592 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/multigrid/mg_dof_tools.h
deal.II/deal.II/include/multigrid/mg_transfer.h
deal.II/deal.II/include/multigrid/mg_transfer.templates.h
deal.II/deal.II/source/multigrid/mg_dof_tools.cc
deal.II/deal.II/source/multigrid/mg_transfer_block.cc

index fd12cf46e41eba2655cea0a844069577c8623b79..5b1271df7cc6a4f271d88ccb11c5dee884a21012 100644 (file)
@@ -180,10 +180,10 @@ class MGTools
                                      */
     template <int dim, typename number>
       static void
-      reinit_vector (const MGDoFHandler<dim> &mg_dof,
-                    MGLevelObject<BlockVector<number> > &v,
-                    const std::vector<bool> &selected = std::vector<bool>(),
-                    const std::vector<unsigned int> &target_component
+      reinit_vector (const MGDoFHandler<dim>mg_dof,
+                    MGLevelObject<BlockVector<number> >v,
+                    const std::vector<bool>selected = std::vector<bool>(),
+                    const std::vector<unsigned int>target_component
                     = std::vector<unsigned int>());
                                     /**
                                      * Adjust vectors on all levels
@@ -216,7 +216,9 @@ class MGTools
       reinit_vector (const MGDoFHandler<dim> &mg_dof,
                     MGLevelObject<Vector<number> > &v,
                     const std::vector<bool> &selected,
-                    const std::vector<unsigned int> &target_component);
+                    const std::vector<unsigned int> &target_component,
+                    std::vector<std::vector<unsigned int> >& cached_sizes
+                    = std::vector<std::vector<unsigned int> >());
 };
 
 
index ff4f49e383c4ab9586fa71c015034146ac1d29d7..1e34181d9422288cf0c413895cbf0ee8c8a0439b 100644 (file)
@@ -299,7 +299,7 @@ class MGTransferBlockBase
                                   /**
                                    * Sizes of the multi-level vectors.
                                    */
-    std::vector<std::vector<unsigned int> > sizes;
+    mutable std::vector<std::vector<unsigned int> > sizes;
   
                                   /**
                                    * Start index of each component.
@@ -347,7 +347,12 @@ class MGTransferBlockBase
     std::vector<boost::shared_ptr<BlockSparseMatrix<double> > > prolongation_matrices;
 #endif
 
-    std::vector<std::vector<std::pair<unsigned int, unsigned int> > >
+                                    /**
+                                     * Unused now, but intended to
+                                     * hold the mapping for the
+                                     * <tt>copy_to/from_mg</tt>-functions.
+                                     */
+    std::vector<std::map<unsigned int, unsigned int> >
     copy_to_and_from_indices;
 };
 
index 309f822c6c776ef70d5615b0f6e6e82c08ae10cb..db3e7a5b3ce43744946064a8f4f78738cd51cdba 100644 (file)
@@ -356,7 +356,7 @@ MGTransferSelect<number>::do_copy_to_mg (
         ExcMatricesNotBuilt());
 
   MGTools::reinit_vector(mg_dof_handler, dst,
-                        mg_selected, mg_target_component);
+                        mg_selected, mg_target_component, sizes);
   
   std::vector<unsigned int> global_dof_indices (dofs_per_cell);
   std::vector<unsigned int> level_dof_indices  (dofs_per_cell);
@@ -386,7 +386,8 @@ MGTransferSelect<number>::do_copy_to_mg (
                                   // already have built.
   for (int level=maxlevel; level>=static_cast<signed int>(minlevel); --level)
     {
-      typename MGDoFHandler<dim>::active_cell_iterator
+
+  typename MGDoFHandler<dim>::active_cell_iterator
        level_cell = mg_dof_handler.begin_active(level);
       const typename MGDoFHandler<dim>::active_cell_iterator
        level_end  = mg_dof_handler.end_active(level);
@@ -423,6 +424,7 @@ MGTransferSelect<number>::do_copy_to_mg (
                = src(global_dof_indices[*i] - offset);
            }
        }
+
                                       // for that part of the level
                                       // which is further refined:
                                       // get the defect by
index efffe1984c8525cc6f99f5122b31126aa59eef6a..268613979e4a737b8b5ca1c1ffb997162bf4bcab 100644 (file)
@@ -13,6 +13,7 @@
 
 
 #include <base/multithread_info.h>
+#include <base/logstream.h>
 #include <base/thread_management.h>
 #include <lac/sparsity_pattern.h>
 #include <lac/block_sparsity_pattern.h>
@@ -353,7 +354,7 @@ MGTools::count_dofs_per_component (
                                                   dofs_in_component[i]);
            };
          threads.join_all();
-         
+
                                           // next count what we got
          for (unsigned int i=0; i<n_components; ++i)
              result[l][target_component[i]]
@@ -462,7 +463,8 @@ void
 MGTools::reinit_vector (const MGDoFHandler<dim> &mg_dof,
                         MGLevelObject<Vector<number> > &v,
                         const std::vector<bool> &selected,
-                       const std::vector<unsigned int> &target_component)
+                       const std::vector<unsigned int> &target_component,
+                       std::vector<std::vector<unsigned int> >& ndofs)
 {
   Assert (selected.size() == target_component.size(),
          ExcDimensionMismatch(selected.size(), target_component.size()));
@@ -479,12 +481,15 @@ MGTools::reinit_vector (const MGDoFHandler<dim> &mg_dof,
   unsigned int selected_block = 0;
   while (!selected[selected_block])
     ++selected_block;
-  
-  std::vector<std::vector<unsigned int> >
-    ndofs(mg_dof.get_tria().n_levels(),
-         std::vector<unsigned int>(target_component.size()));
 
-  count_dofs_per_component (mg_dof, ndofs, target_component);
+  if (ndofs.size() == 0)
+    {
+      std::vector<std::vector<unsigned int> >
+       new_dofs(mg_dof.get_tria().n_levels(),
+                std::vector<unsigned int>(target_component.size()));
+      std::swap(ndofs, new_dofs);
+      count_dofs_per_component (mg_dof, ndofs, target_component);
+    }
   
   for (unsigned int level=v.get_minlevel();
        level<=v.get_maxlevel();++level)
@@ -610,12 +615,14 @@ template void MGTools::reinit_vector<deal_II_dimension> (
   const MGDoFHandler<deal_II_dimension>&,
   MGLevelObject<Vector<double> >&,
   const std::vector<bool>&,
-  const std::vector<unsigned int>&);
+  const std::vector<unsigned int>&,
+  std::vector<std::vector<unsigned int> >&);
 template void MGTools::reinit_vector<deal_II_dimension> (
   const MGDoFHandler<deal_II_dimension>&,
   MGLevelObject<Vector<float> >&,
   const std::vector<bool>&,
-  const std::vector<unsigned int>&);
+  const std::vector<unsigned int>&,
+  std::vector<std::vector<unsigned int> >&);
 
 
 template void MGTools::count_dofs_per_component<deal_II_dimension> (
index defa69125bc9d25ffc0692696609fde7e17474bf..04cea80a775d6b3aa3c4e3921c1464ff7e8cb0ce 100644 (file)
@@ -31,7 +31,7 @@
 
 template <int dim>
 void MGTransferBlockBase::build_matrices (
-  const DoFHandler<dim>& dof,
+  const DoFHandler<dim>&,
   const MGDoFHandler<dim>& mg_dof)
 {
                                   // Fill target component with
@@ -127,6 +127,8 @@ void MGTransferBlockBase::build_matrices (
                                   // be prebuilt, since the
                                   // get_dof_indices functions are
                                   // too slow
+
+  copy_to_and_from_indices.resize(n_levels);
   
 // Building the prolongation matrices starts here!
   
@@ -356,6 +358,7 @@ void MGTransferSelect<number>::build_matrices (
   const std::vector<unsigned int>& t_component,
   const std::vector<unsigned int>& mg_t_component)
 {
+  const FiniteElement<dim>& fe = mg_dof.get_fe();
   unsigned int ncomp = mg_dof.get_fe().n_components();
   
   target_component = t_component;
@@ -389,6 +392,45 @@ void MGTransferSelect<number>::build_matrices (
        }
     }    
   MGTransferBlockBase::build_matrices (dof, mg_dof);
+//TODO:[GK] Rewrite this to suit MGTransferBlock and move to base class
+  
+  std::vector<unsigned int> global_dof_indices (fe.dofs_per_cell);
+  std::vector<unsigned int> level_dof_indices  (fe.dofs_per_cell);
+  for (int level=dof.get_tria().n_levels()-1; level>=0; --level)
+    {
+      typename MGDoFHandler<dim>::active_cell_iterator
+       level_cell = mg_dof.begin_active(level);
+      const typename MGDoFHandler<dim>::active_cell_iterator
+       level_end  = mg_dof.end_active(level);
+      
+                                      // Compute coarse level right hand side
+                                      // by restricting from fine level.
+      for (; level_cell!=level_end; ++level_cell)
+       {
+         DoFObjectAccessor<dim, dim>& global_cell = *level_cell;
+                                          // get the dof numbers of
+                                          // this cell for the global
+                                          // and the level-wise
+                                          // numbering
+         global_cell.get_dof_indices(global_dof_indices);
+         level_cell->get_mg_dof_indices (level_dof_indices);
+         
+         for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
+           {
+             const unsigned int component
+               = mg_target_component[fe.system_to_component_index(i).first];
+             if (mg_selected[component])
+               {
+                 const unsigned int level_start
+                   = mg_component_start[level][component];
+                 copy_to_and_from_indices[level].insert(
+                   std::make_pair(global_dof_indices[i]
+                                  - component_start[selected_component],
+                                  level_dof_indices[i]-level_start));
+               }
+           }
+       }
+    }
 }
 
 

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.