]> https://gitweb.dealii.org/ - dealii.git/commitdiff
fix compilation on intel 13.1
authorTimo Heister <timo.heister@gmail.com>
Sun, 3 Nov 2013 20:43:28 +0000 (20:43 +0000)
committerTimo Heister <timo.heister@gmail.com>
Sun, 3 Nov 2013 20:43:28 +0000 (20:43 +0000)
git-svn-id: https://svn.dealii.org/trunk@31537 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/source/dofs/dof_tools_constraints.cc
deal.II/source/numerics/derivative_approximation.cc
deal.II/source/numerics/matrix_tools.cc

index 9dabf622e954b96e785509cecba6a9cee17543a8..acc7ba92b13d191282b16940b0acd268e8d99560 100644 (file)
@@ -2344,12 +2344,22 @@ namespace DoFTools
         WorkStream::run(coarse_grid.begin_active(),coarse_grid.end(),
             static_cast<std_cxx1x::function<void (typename dealii::DoFHandler<dim,spacedim>
             ::active_cell_iterator const &,Assembler::Scratch const &,Assembler::CopyData<dim,spacedim> 
-            &)> > (std_cxx1x::bind(compute_intergrid_weights_3<dim,spacedim>,std_cxx1x::_1,std_cxx1x::_2,
-            std_cxx1x::_3,std_cxx1x::cref(coarse_grid.get_fe()),std_cxx1x::cref(coarse_to_fine_grid_map),
-            std_cxx1x::cref(parameter_dofs))), static_cast<std_cxx1x::function<void (Assembler
-            ::CopyData<dim,spacedim> const &)> > (std_cxx1x::bind(copy_intergrid_weights_3<dim,spacedim>,
-            std_cxx1x::_1,std_cxx1x::cref(coarse_grid.get_fe()),std_cxx1x::ref(weights))),scratch,
-            copy_data);
+            &)> > (
+                  std_cxx1x::bind(&compute_intergrid_weights_3<dim,spacedim>,
+                                  std_cxx1x::_1,
+                                  std_cxx1x::_2,
+                                  std_cxx1x::_3,
+                                  std_cxx1x::cref(coarse_grid.get_fe()),
+                                  std_cxx1x::cref(coarse_to_fine_grid_map),
+                                  std_cxx1x::cref(parameter_dofs))),
+                       static_cast<std_cxx1x::function<void (Assembler
+            ::CopyData<dim,spacedim> const &)> > (
+                                                 std_cxx1x::bind(&copy_intergrid_weights_3<dim,spacedim>,
+                                                                 std_cxx1x::_1,
+                                                                 std_cxx1x::cref(coarse_grid.get_fe()),
+                                                                 std_cxx1x::ref(weights))),
+                       scratch,
+                       copy_data);
       }
 
 
index a3311a49611d400036887ac68d4db9dc86d56680..e1358c3abbe939d7cd8a6cbc21ff11a43a3bd900 100644 (file)
@@ -663,7 +663,7 @@ approximate_derivative (const Mapping<dim,spacedim>    &mapping,
   WorkStream::run(begin,end,
       static_cast<std_cxx1x::function<void (SynchronousIterators<Iterators> const&,
         internal::Assembler::Scratch const&,internal::Assembler::CopyData &)> >
-      (std_cxx1x::bind(DerivativeApproximation::template approximate<DerivativeDescription,dim,DH,
+      (std_cxx1x::bind(&DerivativeApproximation::template approximate<DerivativeDescription,dim,DH,
                        InputVector,spacedim>,
                        std_cxx1x::_1,
                        std_cxx1x::cref(mapping),
index 054d62345234fc40a4f0797b5be7019d86d34eda..0423f1f79491d201df022c22cc3201e0238d5165 100644 (file)
@@ -1228,17 +1228,21 @@ namespace MatrixCreator
         static_cast<std_cxx1x::function<void (typename DoFHandler<dim,spacedim>::active_cell_iterator
           const &,MatrixCreator::internal::AssemblerBoundary::Scratch const &,
           MatrixCreator::internal::AssemblerBoundary::CopyData<DoFHandler<dim,spacedim> > &)> > 
-        (std_cxx1x::bind(create_boundary_mass_matrix_1<dim,spacedim>,std_cxx1x::_1,std_cxx1x::_2,
+        (std_cxx1x::bind(&create_boundary_mass_matrix_1<dim,spacedim>,std_cxx1x::_1,std_cxx1x::_2,
                          std_cxx1x::_3,
                          std_cxx1x::cref(mapping),std_cxx1x::cref(fe),std_cxx1x::cref(q),
                          std_cxx1x::cref(boundary_functions),coefficient,
                          std_cxx1x::cref(component_mapping))), 
         static_cast<std_cxx1x::function<void (MatrixCreator::internal::AssemblerBoundary
           ::CopyData<DoFHandler<dim,spacedim> > const &)> > (std_cxx1x::bind(
-              copy_boundary_mass_matrix_1<dim,spacedim>,std_cxx1x::_1,
-            std_cxx1x::cref(boundary_functions),std_cxx1x::cref(dof_to_boundary_mapping),
-            std_cxx1x::ref(matrix),std_cxx1x::ref(rhs_vector))),
-        scratch,copy_data);
+              &copy_boundary_mass_matrix_1<dim,spacedim>,
+             std_cxx1x::_1,
+             std_cxx1x::cref(boundary_functions),
+             std_cxx1x::cref(dof_to_boundary_mapping),
+             std_cxx1x::ref(matrix),
+             std_cxx1x::ref(rhs_vector))),
+                   scratch,
+                   copy_data);
   }
 
 
@@ -1639,17 +1643,21 @@ namespace MatrixCreator
         static_cast<std_cxx1x::function<void (typename hp::DoFHandler<dim,spacedim>::active_cell_iterator
           const &,MatrixCreator::internal::AssemblerBoundary::Scratch const &,
           MatrixCreator::internal::AssemblerBoundary::CopyData<hp::DoFHandler<dim,spacedim> > &)> > 
-        (std_cxx1x::bind(create_hp_boundary_mass_matrix_1<dim,spacedim>,std_cxx1x::_1,std_cxx1x::_2,
+        (std_cxx1x::bind( &create_hp_boundary_mass_matrix_1<dim,spacedim>,std_cxx1x::_1,std_cxx1x::_2,
                          std_cxx1x::_3,
                          std_cxx1x::cref(mapping),std_cxx1x::cref(fe_collection),std_cxx1x::cref(q),
                          std_cxx1x::cref(boundary_functions),coefficient,
                          std_cxx1x::cref(component_mapping))), 
-        static_cast<std_cxx1x::function<void (MatrixCreator::internal::AssemblerBoundary
-          ::CopyData<hp::DoFHandler<dim,spacedim> > const &)> > (std_cxx1x::bind(
-              copy_hp_boundary_mass_matrix_1<dim,spacedim>,std_cxx1x::_1,
-            std_cxx1x::cref(boundary_functions),std_cxx1x::cref(dof_to_boundary_mapping),
-            std_cxx1x::ref(matrix),std_cxx1x::ref(rhs_vector))),
-        scratch,copy_data);
+                   static_cast<std_cxx1x::function<void (MatrixCreator::internal::AssemblerBoundary
+                                                         ::CopyData<hp::DoFHandler<dim,spacedim> > const &)> > (
+                      std_cxx1x::bind( &copy_hp_boundary_mass_matrix_1<dim,spacedim>,
+                                       std_cxx1x::_1,
+                                       std_cxx1x::cref(boundary_functions),
+                                       std_cxx1x::cref(dof_to_boundary_mapping),
+                                       std_cxx1x::ref(matrix),
+                                       std_cxx1x::ref(rhs_vector))),
+                   scratch,
+                   copy_data);
   }
 
 
@@ -1700,17 +1708,14 @@ namespace MatrixCreator
     copy_data.dof_indices.resize (assembler_data.fe_collection.max_dofs_per_cell());
     copy_data.constraints = &constraints;
 
-    void (*copy_local_to_global) (const MatrixCreator::internal::AssemblerData::CopyData &,
-                                  SparseMatrix<double> *,
-                                  Vector<double> *)
-      = &MatrixCreator::internal::
-        copy_local_to_global<SparseMatrix<double>, Vector<double> >;
-
     WorkStream::run (dof.begin_active(),
                      static_cast<typename DoFHandler<dim,spacedim>::active_cell_iterator>(dof.end()),
                      &MatrixCreator::internal::laplace_assembler<dim, spacedim, typename DoFHandler<dim,spacedim>::active_cell_iterator>,
-                     std_cxx1x::bind (copy_local_to_global,
-                                      std_cxx1x::_1, &matrix, (Vector<double> *)0),
+                     std_cxx1x::bind (&MatrixCreator::internal::
+                                     copy_local_to_global<SparseMatrix<double>, Vector<double> >,
+                                      std_cxx1x::_1, 
+                                     &matrix, 
+                                     (Vector<double> *)NULL),
                      assembler_data,
                      copy_data);
   }
@@ -1760,17 +1765,14 @@ namespace MatrixCreator
     copy_data.dof_indices.resize (assembler_data.fe_collection.max_dofs_per_cell());
     copy_data.constraints = &constraints;
 
-    void (*copy_local_to_global) (const MatrixCreator::internal::AssemblerData::CopyData &,
-                                  SparseMatrix<double> *,
-                                  Vector<double> *)
-      = &MatrixCreator::internal::
-        copy_local_to_global<SparseMatrix<double>, Vector<double> >;
-
     WorkStream::run (dof.begin_active(),
                      static_cast<typename DoFHandler<dim,spacedim>::active_cell_iterator>(dof.end()),
                      &MatrixCreator::internal::laplace_assembler<dim, spacedim, typename DoFHandler<dim,spacedim>::active_cell_iterator>,
-                     std_cxx1x::bind (copy_local_to_global,
-                                      std_cxx1x::_1, &matrix, &rhs_vector),
+                     std_cxx1x::bind (&MatrixCreator::internal::
+                                     copy_local_to_global<SparseMatrix<double>, Vector<double> >,
+                                      std_cxx1x::_1, 
+                                     &matrix, 
+                                     &rhs_vector),
                      assembler_data,
                      copy_data);
   }
@@ -1818,17 +1820,14 @@ namespace MatrixCreator
     copy_data.dof_indices.resize (assembler_data.fe_collection.max_dofs_per_cell());
     copy_data.constraints = &constraints;
 
-    void (*copy_local_to_global) (const MatrixCreator::internal::AssemblerData::CopyData &,
-                                  SparseMatrix<double> *,
-                                  Vector<double> *)
-      = &MatrixCreator::internal::
-        copy_local_to_global<SparseMatrix<double>, Vector<double> >;
-
     WorkStream::run (dof.begin_active(),
                      static_cast<typename hp::DoFHandler<dim,spacedim>::active_cell_iterator>(dof.end()),
                      &MatrixCreator::internal::laplace_assembler<dim, spacedim, typename hp::DoFHandler<dim,spacedim>::active_cell_iterator>,
-                     std_cxx1x::bind (copy_local_to_global,
-                                      std_cxx1x::_1, &matrix, (Vector<double> *)0),
+                     std_cxx1x::bind (&MatrixCreator::internal::
+                                     copy_local_to_global<SparseMatrix<double>, Vector<double> >,
+                                      std_cxx1x::_1, 
+                                     &matrix, 
+                                     (Vector<double> *)0),
                      assembler_data,
                      copy_data);
   }
@@ -1875,17 +1874,14 @@ namespace MatrixCreator
     copy_data.dof_indices.resize (assembler_data.fe_collection.max_dofs_per_cell());
     copy_data.constraints = &constraints;
 
-    void (*copy_local_to_global) (const MatrixCreator::internal::AssemblerData::CopyData &,
-                                  SparseMatrix<double> *,
-                                  Vector<double> *)
-      = &MatrixCreator::internal::
-        copy_local_to_global<SparseMatrix<double>, Vector<double> >;
-
     WorkStream::run (dof.begin_active(),
                      static_cast<typename hp::DoFHandler<dim,spacedim>::active_cell_iterator>(dof.end()),
                      &MatrixCreator::internal::laplace_assembler<dim, spacedim, typename hp::DoFHandler<dim,spacedim>::active_cell_iterator>,
-                     std_cxx1x::bind (copy_local_to_global,
-                                      std_cxx1x::_1, &matrix, &rhs_vector),
+                     std_cxx1x::bind (&MatrixCreator::internal::
+                                     copy_local_to_global<SparseMatrix<double>, Vector<double> >,
+                                      std_cxx1x::_1, 
+                                     &matrix, 
+                                     &rhs_vector),
                      assembler_data,
                      copy_data);
   }

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.