]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use WorkStream instead of Thread in DerivativeApproximation::approximate_derivative.
authorBruno Turcksin <bruno.turcksin@gmail.com>
Thu, 10 Oct 2013 16:56:19 +0000 (16:56 +0000)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Thu, 10 Oct 2013 16:56:19 +0000 (16:56 +0000)
git-svn-id: https://svn.dealii.org/trunk@31195 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/numerics/derivative_approximation.h
deal.II/source/numerics/derivative_approximation.cc

index 34ebddb3815e6f793cb781651bc4705811448521..3fd49856c484eddff2fe64a7d55f4f8ca7260821 100644 (file)
@@ -665,21 +665,19 @@ private:
 
   /**
    * Compute the derivative
-   * approximation on the cells in
-   * the range given by the third
-   * parameter.
+   * approximation on a given cell.
    * Fill the @p derivative_norm vector with
    * the norm of the computed derivative
-   * tensors on each cell.
+   * tensors on the cell.
    */
   template <class DerivativeDescription, int dim,
            template <int, int> class DH, class InputVector, int spacedim>
   static void
-  approximate (const Mapping<dim,spacedim>    &mapping,
+  approximate (const typename DH<dim,spacedim>::active_cell_iterator &cell,
+               const Mapping<dim,spacedim>    &mapping,
                const DH<dim,spacedim>         &dof,
                const InputVector     &solution,
                const unsigned int     component,
-               const IndexInterval   &index_interval,
                Vector<float>         &derivative_norm);
 
   /**
index ce8f37fde0d245d4a79977f4aba604793a04f09a..397e688c60d98f74201d95f84856167b24e4ee8c 100644 (file)
@@ -15,8 +15,8 @@
 // ---------------------------------------------------------------------
 
 #include <deal.II/base/quadrature_lib.h>
-#include <deal.II/base/thread_management.h>
 #include <deal.II/base/multithread_info.h>
+#include <deal.II/base/work_stream.h>
 #include <deal.II/lac/vector.h>
 #include <deal.II/lac/block_vector.h>
 #include <deal.II/lac/parallel_vector.h>
@@ -27,6 +27,7 @@
 #include <deal.II/lac/trilinos_block_vector.h>
 #include <deal.II/grid/tria_iterator.h>
 #include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/filtered_iterator.h>
 #include <deal.II/dofs/dof_accessor.h>
 #include <deal.II/dofs/dof_handler.h>
 #include <deal.II/fe/fe.h>
@@ -64,6 +65,27 @@ const UpdateFlags DerivativeApproximation::ThirdDerivative<dim>::update_flags =
 
 
 
+// Dummy structures and dummy function used for WorkStream
+namespace internal
+{
+  namespace Assembler
+  {
+    struct Scratch
+    {
+      Scratch() {}
+    };
+
+    struct CopyData
+    {
+      CopyData() {}
+    };
+
+    void copier(CopyData const &) {}
+  }
+}
+
+
+
 template <int dim>
 template <class InputVector, int spacedim>
 inline
@@ -626,30 +648,25 @@ approximate_derivative (const Mapping<dim,spacedim>    &mapping,
   Assert (component < dof_handler.get_fe().n_components(),
           ExcIndexRange (component, 0, dof_handler.get_fe().n_components()));
 
-  const unsigned int n_threads = multithread_info.n_threads();
-  std::vector<IndexInterval> index_intervals
-    = Threads::split_interval (0, dof_handler.get_tria().n_active_cells(),
-                               n_threads);
-
-  typedef void (*FunPtr) (const Mapping<dim,spacedim> &,
-                          const DH<dim,spacedim> &,
-                          const InputVector &,
-                          const unsigned int     ,
-                          const IndexInterval &,
-                          Vector<float> &);
-  FunPtr fun_ptr = &DerivativeApproximation::
-                   template approximate<DerivativeDescription,dim,
-                   DH,InputVector>;
-
-//TODO: Use WorkStream here
-  Threads::TaskGroup<> tasks;
-  for (unsigned int i=0; i<n_threads; ++i)
-    tasks += Threads::new_task (fun_ptr,
-                                mapping, dof_handler,
-                                solution, component,
-                                index_intervals[i],
-                                derivative_norm);
-  tasks.join_all ();
+  // Only act on the locally owned cells
+  typedef FilteredIterator<typename DH<dim,spacedim>::active_cell_iterator> CellFilter;
+
+  // There is no need for a copier because there is no conflict between threads
+  // to write in derivative_norm. Scratch and CopyData are also useless.
+  WorkStream::run(CellFilter(IteratorFilters::LocallyOwnedCell(),dof_handler.begin_active()),
+      CellFilter(IteratorFilters::LocallyOwnedCell(),dof_handler.end()),
+      static_cast<std_cxx1x::function<void (typename DH<dim,spacedim>::active_cell_iterator const&,
+        internal::Assembler::Scratch const&,internal::Assembler::CopyData &)> >
+      (std_cxx1x::bind(DerivativeApproximation::template approximate<DerivativeDescription,dim,DH,
+                       InputVector,spacedim>,
+                       std_cxx1x::_1,
+                       std_cxx1x::cref(mapping),
+                       std_cxx1x::cref(dof_handler),
+                       std_cxx1x::cref(solution),component,
+                       std_cxx1x::ref(derivative_norm))),
+    static_cast<std_cxx1x::function<void (internal::Assembler::CopyData const &)> >
+      (internal::Assembler::copier),internal::Assembler::Scratch (),
+      internal::Assembler::CopyData ());
 }
 
 
@@ -657,48 +674,31 @@ approximate_derivative (const Mapping<dim,spacedim>    &mapping,
 template <class DerivativeDescription, int dim,
          template <int, int> class DH, class InputVector, int spacedim>
 void
-DerivativeApproximation::approximate (const Mapping<dim,spacedim>    &mapping,
-                                      const DH<dim,spacedim>         &dof_handler,
-                                      const InputVector     &solution,
-                                      const unsigned int     component,
-                                      const IndexInterval   &index_interval,
-                                      Vector<float>         &derivative_norm)
+DerivativeApproximation::approximate (const typename DH<dim,spacedim>::active_cell_iterator &cell,
+                                      const Mapping<dim,spacedim>                  &mapping,
+                                      const DH<dim,spacedim>                       &dof_handler,
+                                      const InputVector                            &solution,
+                                      const unsigned int                            component,
+                                      Vector<float>                                &derivative_norm)
 {
   // iterators over all cells and the
   // respective entries in the output
   // vector:
-  Vector<float>::iterator
-  derivative_norm_on_this_cell
-    = derivative_norm.begin() + index_interval.first;
-
-  typename DH<dim,spacedim>::active_cell_iterator cell, endc;
-  cell = endc = dof_handler.begin_active();
-  // (static_cast to avoid warnings
-  // about unsigned int always >=0)
-  std::advance (cell, static_cast<int>(index_interval.first));
-  std::advance (endc, static_cast<int>(index_interval.second));
-
-  for (; cell!=endc; ++cell, ++derivative_norm_on_this_cell)
-    if (cell->is_locally_owned())
-      {
-        typename DerivativeDescription::Derivative derivative;
-        // call the function doing the actual
-        // work on this cell
-        DerivativeApproximation::
-        template approximate_cell<DerivativeDescription,dim,DH,InputVector>
-        (mapping,
-         dof_handler,
-         solution,
-         component,
-         cell,
-         derivative);
-        // evaluate the norm and fill the vector
-        *derivative_norm_on_this_cell
-          = DerivativeDescription::derivative_norm (derivative);
-      }
+  Vector<float>::iterator derivative_norm_on_this_cell = derivative_norm.begin() + 
+    std::distance(dof_handler.begin_active(),cell);
+
+  typename DerivativeDescription::Derivative derivative;
+  // call the function doing the actual
+  // work on this cell
+  DerivativeApproximation::template approximate_cell<DerivativeDescription,dim,DH,InputVector>
+    (mapping,dof_handler,solution,component,cell,derivative);
+  // evaluate the norm and fill the vector
+  *derivative_norm_on_this_cell
+    = DerivativeDescription::derivative_norm (derivative);
 }
 
 
+
 template <class DerivativeDescription, int dim,
          template <int, int> class DH, class InputVector, int spacedim>
 void

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.