]> https://gitweb.dealii.org/ - dealii.git/commitdiff
MeshWorker::loop() uses WorkStream::run()
authorGuido Kanschat <dr.guido.kanschat@gmail.com>
Tue, 23 Feb 2010 23:05:01 +0000 (23:05 +0000)
committerGuido Kanschat <dr.guido.kanschat@gmail.com>
Tue, 23 Feb 2010 23:05:01 +0000 (23:05 +0000)
git-svn-id: https://svn.dealii.org/trunk@20679 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/numerics/mesh_worker_info.h
deal.II/deal.II/include/numerics/mesh_worker_info.templates.h
deal.II/deal.II/include/numerics/mesh_worker_loop.h

index 54c2e54564847de85b10ba966d2443cebc9da7a6..f48fe4b9e19b825048f6a6cd3297a2ada013d422 100644 (file)
@@ -359,7 +359,7 @@ namespace MeshWorker
  *
  * @author Guido Kanschat, 2010
  */
-  template < int dim, int spacedim=dim>
+  template <int dim, int spacedim=dim>
   class DoFInfoBox
   {
     public:
@@ -367,8 +367,15 @@ namespace MeshWorker
                                        * Constructor copying the seed
                                        * into all other objects.
                                        */
-      DoFInfoBox(const MeshWorker::DoFInfo<dim, spacedim>& seed);
+      DoFInfoBox(const DoFInfo<dim, spacedim>& seed);
 
+                                      /**
+                                       * Copy constructor, taking
+                                       * #cell and using it as a seed
+                                       * in the other constructor.
+                                       */
+      DoFInfoBox(const DoFInfoBox<dim, spacedim>&);
+      
                                       /**
                                        * Reset all the availability flags.
                                        */
@@ -390,15 +397,15 @@ namespace MeshWorker
                                       /**
                                        * The data for the cell.
                                        */
-      MeshWorker::DoFInfo<dim> cell;
+      DoFInfo<dim> cell;
                                       /**
                                        * The data for the faces from inside.
                                        */
-      MeshWorker::DoFInfo<dim> interior[GeometryInfo<dim>::faces_per_cell];
+      DoFInfo<dim> interior[GeometryInfo<dim>::faces_per_cell];
                                       /**
                                        * The data for the faces from outside.
                                        */
-      MeshWorker::DoFInfo<dim> exterior[GeometryInfo<dim>::faces_per_cell];
+      DoFInfo<dim> exterior[GeometryInfo<dim>::faces_per_cell];
 
                                       /**
                                        * A set of flags, indicating
@@ -474,11 +481,21 @@ namespace MeshWorker
                                       /// vector of FEValues objects
       std::vector<boost::shared_ptr<FEVALUESBASE> > fevalv;
     public:
+      static const unsigned int dimension = FEVALUESBASE::dimension;
+      static const unsigned int space_dimension = FEVALUESBASE::space_dimension;
+      
                                       /**
                                        * Constructor.
                                        */
       IntegrationInfo();
-
+      
+                                      /**
+                                       * Copy constructor, creating a
+                                       * clone to be used by
+                                       * WorksTream::run().
+                                       */
+      IntegrationInfo(const IntegrationInfo<dim, FEVALUESBASE, spacedim>& other);
+      
                                       /**
                                        * Build all internal
                                        * structures, in particular
@@ -977,6 +994,22 @@ namespace MeshWorker
   }
 
 
+  template < int dim, int spacedim>
+  inline
+  DoFInfoBox<dim, spacedim>::DoFInfoBox(const DoFInfoBox<dim, spacedim>& other)
+                 :
+                 cell(other.cell)
+  {
+    for (unsigned int i=0;i<GeometryInfo<dim>::faces_per_cell;++i)
+      {
+       exterior[i] = other.exterior[i];
+       interior[i] = other.interior[i];
+       interior_face_available[i] = false;
+       exterior_face_available[i] = false;
+      }
+  }
+
+
   template < int dim, int spacedim>
   inline void
   DoFInfoBox<dim, spacedim>::reset ()
index 3aa4d16b039ada79fe4f00b4821279e516b9f783..a5b1a956d2a53b659b12c95a26031094ef4e95ae 100644 (file)
@@ -75,6 +75,40 @@ namespace MeshWorker
                  global_data(boost::shared_ptr<VectorDataBase<dim, sdim> >(new VectorDataBase<dim, sdim>))
   {}
 
+  
+  template<int dim, class FVB, int sdim>
+  IntegrationInfo<dim,FVB,sdim>::IntegrationInfo(const IntegrationInfo<dim,FVB,sdim>& other)
+                 :
+                 multigrid(other.multigrid),
+                 values(other.values),
+                 gradients(other.gradients),
+                 hessians(other.hessians),
+                 global_data(other.global_data),
+                 n_components(other.n_components)
+  {
+    fevalv.resize(other.fevalv.size());
+    for (unsigned int i=0;i<other.fevalv.size();++i)
+      {
+       const FVB& p = *other.fevalv[i];
+       const FEValues<dim,sdim>* pc = dynamic_cast<const FEValues<dim,sdim>*>(&p);
+       const FEFaceValues<dim,sdim>* pf = dynamic_cast<const FEFaceValues<dim,sdim>*>(&p);
+       const FESubfaceValues<dim,sdim>* ps = dynamic_cast<const FESubfaceValues<dim,sdim>*>(&p);
+       
+       if (pc != 0)
+         fevalv[i] = typename boost::shared_ptr<FVB> (
+           reinterpret_cast<FVB*>(new FEValues<dim,sdim> (pc->get_mapping(), pc->get_fe(),
+                                                          pc->get_quadrature(), pc->get_update_flags())));
+       else if (pf != 0)
+         fevalv[i] = typename boost::shared_ptr<FVB> (
+           new FEFaceValues<dim,sdim> (pf->get_mapping(), pf->get_fe(), pf->get_quadrature(), pf->get_update_flags()));
+       else if (ps != 0)
+         fevalv[i] = typename boost::shared_ptr<FVB> (
+           new FESubfaceValues<dim,sdim> (ps->get_mapping(), ps->get_fe(), ps->get_quadrature(), ps->get_update_flags()));
+       else
+         Assert(false, ExcInternalError());
+      }
+  }
+  
 
   template<int dim, class FVB, int sdim>
   template <class FEVALUES>
index 05023e5e9ad776e56eebdebe5d92a91461f9df07..fe628ed2556ac194515a7a2ce530be93a15c8546 100644 (file)
@@ -15,6 +15,7 @@
 
 #include <base/config.h>
 #include <base/std_cxx1x/function.h>
+#include <base/work_stream.h>
 #include <base/template_constraints.h>
 #include <grid/tria.h>
 #include <numerics/mesh_worker_info.h>
@@ -41,6 +42,11 @@ namespace internal
     return true;
   }
 
+  template<int dim, int sdim, class A>
+  void assemble(const MeshWorker::DoFInfoBox<dim, sdim>& dinfo, A& assembler)
+  {
+    dinfo.assemble(assembler);
+  }
 }
 
 
@@ -242,11 +248,17 @@ namespace MeshWorker
       }
     
                                     // Loop over all cells
-    for (ITERATOR cell = begin; cell != end; ++cell)
-      {
-       cell_action(cell, dof_info, info, cell_worker, boundary_worker, face_worker, cells_first);
-       dof_info.assemble(assembler);
-      }
+    WorkStream::run(begin, end,
+                   std_cxx1x::bind(cell_action<INFOBOX, dim, spacedim, ITERATOR>, _1, _3, _2,
+                                   cell_worker, boundary_worker, face_worker, cells_first),
+                   std_cxx1x::bind(internal::assemble<dim,spacedim,ASSEMBLER>, _1, assembler),
+                   info, dof_info);
+    
+//     for (ITERATOR cell = begin; cell != end; ++cell)
+//       {
+//     cell_action(cell, dof_info, info, cell_worker, boundary_worker, face_worker, cells_first);
+//     dof_info.assemble(assembler);
+//       }
   }
 
 /**

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.