]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Squashed changes for build working both on MSVC and other compilers. 289/head
authorLukas Korous <lukas.korous@gmail.com>
Fri, 12 Dec 2014 15:05:07 +0000 (16:05 +0100)
committerLukas Korous <lukas.korous@gmail.com>
Fri, 12 Dec 2014 15:05:07 +0000 (16:05 +0100)
include/deal.II/grid/grid_tools.h
include/deal.II/numerics/derivative_approximation.h
include/deal.II/numerics/fe_field_function.templates.h
source/grid/grid_tools.cc
source/grid/grid_tools.inst.in
source/numerics/derivative_approximation.cc

index 4e615cc1dd285b30cfc447bc919f7717273b7e0b..f9ce7420f9d88f67478fcbd590e9fbb455dce646 100644 (file)
@@ -40,6 +40,37 @@ namespace hp
 
 class SparsityPattern;
 
+namespace internal
+{
+  template<int dim, int spacedim, class Container>
+  class ActiveCellIterator
+  {
+  public:
+    typedef typename Container::active_cell_iterator type;
+  };
+
+  template<int dim, int spacedim>
+  class ActiveCellIterator<dim, spacedim, dealii::DoFHandler<dim, spacedim> >
+  {
+  public:
+#ifndef _MSC_VER
+    typedef typename dealii::DoFHandler<dim, spacedim>::active_cell_iterator type;
+#else
+    typedef TriaActiveIterator < dealii::DoFCellAccessor < dealii::DoFHandler<dim, spacedim>, false > > type;
+#endif
+  };
+
+  template<int dim, int spacedim>
+  class ActiveCellIterator<dim, spacedim, dealii::hp::DoFHandler<dim, spacedim> >
+  {
+  public:
+#ifndef _MSC_VER
+    typedef typename dealii::hp::DoFHandler<dim, spacedim>::active_cell_iterator type;
+#else
+    typedef TriaActiveIterator < dealii::DoFCellAccessor < dealii::hp::DoFHandler<dim, spacedim>, false > > type;
+#endif
+  };
+}
 
 /**
  * This namespace is a collection of algorithms working on triangulations,
@@ -426,7 +457,11 @@ namespace GridTools
    * to be checked for this case.
    */
   template<int dim, template <int, int> class Container, int spacedim>
-  std::vector<typename Container<dim,spacedim>::active_cell_iterator>
+#ifndef _MSC_VER
+  std::vector<typename Container<dim, spacedim>::active_cell_iterator>
+#else
+  std::vector<typename dealii::internal::ActiveCellIterator<dim, spacedim, Container<dim, spacedim> >::type>
+#endif
   find_cells_adjacent_to_vertex (const Container<dim,spacedim> &container,
                                  const unsigned int    vertex_index);
 
@@ -463,7 +498,11 @@ namespace GridTools
    * will have to decide what to do in that case.
    */
   template <int dim, template <int,int> class Container, int spacedim>
+#ifndef _MSC_VER
   typename Container<dim,spacedim>::active_cell_iterator
+#else
+  typename dealii::internal::ActiveCellIterator<dim, spacedim, Container<dim, spacedim> >::type
+#endif
   find_active_cell_around_point (const Container<dim,spacedim>  &container,
                                  const Point<spacedim> &p);
 
@@ -514,7 +553,11 @@ namespace GridTools
    * will have to decide what to do in that case.
    */
   template <int dim, template<int, int> class Container, int spacedim>
-  std::pair<typename Container<dim,spacedim>::active_cell_iterator, Point<dim> >
+#ifndef _MSC_VER
+  std::pair<typename Container<dim, spacedim>::active_cell_iterator, Point<dim> >
+#else
+  std::pair<typename dealii::internal::ActiveCellIterator<dim, spacedim, Container<dim, spacedim> >::type, Point<dim> >
+#endif
   find_active_cell_around_point (const Mapping<dim,spacedim>   &mapping,
                                  const Container<dim,spacedim> &container,
                                  const Point<spacedim>     &p);
@@ -540,7 +583,7 @@ namespace GridTools
    * will have to decide what to do in that case.
    */
   template <int dim, int spacedim>
-  std::pair<typename hp::DoFHandler<dim,spacedim>::active_cell_iterator, Point<dim> >
+  std::pair<typename hp::DoFHandler<dim, spacedim>::active_cell_iterator, Point<dim> >
   find_active_cell_around_point (const hp::MappingCollection<dim,spacedim>   &mapping,
                                  const hp::DoFHandler<dim,spacedim> &container,
                                  const Point<spacedim>     &p);
index 875cbcb10195e19a5220390e6f21c195128876e8..da68daee45e66394fb3f2f09b65d8024ea12afe2 100644 (file)
@@ -24,7 +24,9 @@
 #include <deal.II/fe/mapping.h>
 #include <deal.II/lac/vector.h>
 #include <deal.II/grid/filtered_iterator.h>
-
+#ifdef _MSC_VER
+#include <deal.II/dofs/dof_accessor.h>
+#endif
 #include <utility>
 
 DEAL_II_NAMESPACE_OPEN
@@ -275,9 +277,9 @@ namespace DerivativeApproximation
                                 const DH                      &dof,
                                 const InputVector                            &solution,
 #ifndef _MSC_VER
-                                const typename DH::active_cell_iterator      &cell,
+                                const typename DH::active_cell_iterator &cell,
 #else
-                                const TriaActiveIterator < dealii::DoFCellAccessor < DH, false > >      &cell,
+                                const TriaActiveIterator < dealii::DoFCellAccessor < DH, false > > &cell,
 #endif
                                 Tensor<order, dim>  &derivative,
                                 const unsigned int                            component = 0);
@@ -290,9 +292,9 @@ namespace DerivativeApproximation
   approximate_derivative_tensor(const DH                    &dof,
                                 const InputVector                            &solution,
 #ifndef _MSC_VER
-                                const typename DH::active_cell_iterator      &cell,
+                                const typename DH::active_cell_iterator &cell,
 #else
-                                const TriaActiveIterator < dealii::DoFCellAccessor < DH, false > >      &cell,
+                                const TriaActiveIterator < dealii::DoFCellAccessor < DH, false > > &cell,
 #endif
                                 Tensor<order, dim>                  &derivative,
                                 const unsigned int                            component = 0);
index c84f8e0a24d9828c82f41c37f927820c6b2ec57d..b13770f4100481396e0abc830f9910038a1de265 100644 (file)
@@ -70,7 +70,7 @@ namespace Functions
     qp = get_reference_coordinates (cell, p);
     if (!qp)
       {
-        const std::pair<typename DH::active_cell_iterator, Point<dim> > my_pair
+        const std::pair<typename dealii::internal::ActiveCellIterator<dim, dim, DH>::type, Point<dim> > my_pair
           = GridTools::find_active_cell_around_point (mapping, *dh, p);
         AssertThrow (my_pair.first->is_locally_owned(),
                      ExcPointNotAvailableHere());
@@ -121,7 +121,7 @@ namespace Functions
     qp = get_reference_coordinates (cell, p);
     if (!qp)
       {
-        const std::pair<typename DH::active_cell_iterator, Point<dim> > my_pair
+        const std::pair<typename dealii::internal::ActiveCellIterator<dim, dim, DH>::type, Point<dim> > my_pair
           = GridTools::find_active_cell_around_point (mapping, *dh, p);
         AssertThrow (my_pair.first->is_locally_owned(),
                      ExcPointNotAvailableHere());
@@ -174,7 +174,7 @@ namespace Functions
     qp = get_reference_coordinates (cell, p);
     if (!qp)
       {
-        const std::pair<typename DH::active_cell_iterator, Point<dim> > my_pair
+        const std::pair<typename dealii::internal::ActiveCellIterator<dim, dim, DH>::type, Point<dim> > my_pair
           = GridTools::find_active_cell_around_point (mapping, *dh, p);
         AssertThrow (my_pair.first->is_locally_owned(),
                      ExcPointNotAvailableHere());
@@ -441,7 +441,7 @@ namespace Functions
       qp = get_reference_coordinates (cell, points[0]);
       if (!qp)
         {
-          const std::pair<typename DH::active_cell_iterator, Point<dim> >
+          const std::pair<typename dealii::internal::ActiveCellIterator<dim, dim, DH>::type, Point<dim> >
           my_pair  = GridTools::find_active_cell_around_point
                      (mapping, *dh, points[0]);
           AssertThrow (my_pair.first->is_locally_owned(),
@@ -506,7 +506,7 @@ namespace Functions
         // the next cell
         if (left_over == true)
           {
-            const std::pair<typename DH::active_cell_iterator, Point<dim> > my_pair
+            const std::pair<typename dealii::internal::ActiveCellIterator<dim, dim, DH>::type, Point<dim> > my_pair
               = GridTools::find_active_cell_around_point (mapping, *dh, points[first_outside]);
             AssertThrow (my_pair.first->is_locally_owned(),
                          ExcPointNotAvailableHere());
index 72a11fedef54e79e857465b695bfccbbf909f7a2..95ee495eee291cfb825300ecf1c0a897f60433fe 100644 (file)
@@ -910,7 +910,11 @@ namespace GridTools
 
 
   template<int dim, template<int, int> class Container, int spacedim>
-  std::vector<typename Container<dim,spacedim>::active_cell_iterator>
+#ifndef _MSC_VER
+  std::vector<typename Container<dim, spacedim>::active_cell_iterator>
+#else
+  std::vector<typename dealii::internal::ActiveCellIterator<dim, spacedim, Container<dim, spacedim> >::type>
+#endif
   find_cells_adjacent_to_vertex(const Container<dim,spacedim> &container,
                                 const unsigned int    vertex)
   {
@@ -925,9 +929,9 @@ namespace GridTools
     // use a set instead of a vector
     // to ensure that cells are inserted only
     // once
-    std::set<typename Container<dim,spacedim>::active_cell_iterator> adjacent_cells;
+    std::set<typename dealii::internal::ActiveCellIterator<dim, spacedim, Container<dim, spacedim> >::type> adjacent_cells;
 
-    typename Container<dim,spacedim>::active_cell_iterator
+    typename dealii::internal::ActiveCellIterator<dim, spacedim, Container<dim, spacedim> >::type
     cell = container.begin_active(),
     endc = container.end();
 
@@ -1039,7 +1043,7 @@ next_cell:
 
     // return the result as a vector, rather than the set we built above
     return
-      std::vector<typename Container<dim,spacedim>::active_cell_iterator>
+      std::vector<typename dealii::internal::ActiveCellIterator<dim, spacedim, Container<dim, spacedim> >::type>
       (adjacent_cells.begin(), adjacent_cells.end());
   }
 
@@ -1049,10 +1053,19 @@ next_cell:
   {
     template <int dim, template<int, int> class Container, int spacedim>
     void find_active_cell_around_point_internal(const Container<dim,spacedim> &container,
-                                                std::set<typename Container<dim,spacedim>::active_cell_iterator> &searched_cells,
-                                                std::set<typename Container<dim,spacedim>::active_cell_iterator> &adjacent_cells)
+#ifndef _MSC_VER
+                                                std::set<typename Container<dim, spacedim>::active_cell_iterator> &searched_cells,
+                                                std::set<typename Container<dim, spacedim>::active_cell_iterator> &adjacent_cells)
+#else
+                                                std::set<typename dealii::internal::ActiveCellIterator<dim, spacedim, Container<dim, spacedim> >::type> &searched_cells,
+                                                std::set<typename dealii::internal::ActiveCellIterator<dim, spacedim, Container<dim, spacedim> >::type> &adjacent_cells)
+#endif
     {
-      typedef typename Container<dim,spacedim>::active_cell_iterator cell_iterator;
+#ifndef _MSC_VER
+      typedef typename Container<dim, spacedim>::active_cell_iterator cell_iterator;
+#else
+      typedef typename dealii::internal::ActiveCellIterator<dim, spacedim, Container<dim, spacedim> >::type cell_iterator;
+#endif
 
       // update the searched cells
       searched_cells.insert(adjacent_cells.begin(), adjacent_cells.end());
@@ -1094,7 +1107,11 @@ next_cell:
   }
 
   template <int dim, template<int, int> class Container, int spacedim>
-  typename Container<dim,spacedim>::active_cell_iterator
+#ifndef _MSC_VER
+  typename Container<dim, spacedim>::active_cell_iterator
+#else
+  typename dealii::internal::ActiveCellIterator<dim, spacedim, Container<dim, spacedim> >::type
+#endif
   find_active_cell_around_point (const Container<dim,spacedim>  &container,
                                  const Point<spacedim> &p)
   {
@@ -1106,12 +1123,16 @@ next_cell:
 
 
   template <int dim, template <int, int> class Container, int spacedim>
-  std::pair<typename Container<dim,spacedim>::active_cell_iterator, Point<dim> >
+#ifndef _MSC_VER
+  std::pair<typename Container<dim, spacedim>::active_cell_iterator, Point<dim> >
+#else
+  std::pair<typename dealii::internal::ActiveCellIterator<dim, spacedim, Container<dim, spacedim> >::type, Point<dim> >
+#endif
   find_active_cell_around_point (const Mapping<dim,spacedim>   &mapping,
                                  const Container<dim,spacedim> &container,
                                  const Point<spacedim>     &p)
   {
-    typedef typename Container<dim,spacedim>::active_cell_iterator active_cell_iterator;
+    typedef typename dealii::internal::ActiveCellIterator<dim, spacedim, Container<dim, spacedim> >::type active_cell_iterator;
 
     // The best distance is set to the
     // maximum allowable distance from
@@ -1367,7 +1388,7 @@ next_cell:
     std::map< std::pair<unsigned int,unsigned int>, unsigned int >
     indexmap;
     unsigned int index = 0;
-    for (typename Triangulation<dim,spacedim>::active_cell_iterator
+    for (typename dealii::internal::ActiveCellIterator<dim, spacedim, Triangulation<dim, spacedim> >::type
          cell = triangulation.begin_active();
          cell != triangulation.end(); ++cell, ++index)
       indexmap[std::pair<unsigned int,unsigned int>(cell->level(),cell->index())] = index;
@@ -1389,7 +1410,7 @@ next_cell:
     // add entries in both directions
     // for both cells
     index = 0;
-    for (typename Triangulation<dim,spacedim>::active_cell_iterator
+    for (typename dealii::internal::ActiveCellIterator<dim, spacedim, Triangulation<dim, spacedim> >::type
          cell = triangulation.begin_active();
          cell != triangulation.end(); ++cell, ++index)
       {
@@ -1428,7 +1449,7 @@ next_cell:
     // check for an easy return
     if (n_partitions == 1)
       {
-        for (typename Triangulation<dim,spacedim>::active_cell_iterator
+        for (typename dealii::internal::ActiveCellIterator<dim, spacedim, Triangulation<dim, spacedim> >::type
              cell = triangulation.begin_active();
              cell != triangulation.end(); ++cell)
           cell->set_subdomain_id (0);
@@ -1472,7 +1493,7 @@ next_cell:
     // check for an easy return
     if (n_partitions == 1)
       {
-        for (typename Triangulation<dim,spacedim>::active_cell_iterator
+        for (typename dealii::internal::ActiveCellIterator<dim, spacedim, Triangulation<dim, spacedim> >::type
              cell = triangulation.begin_active();
              cell != triangulation.end(); ++cell)
           cell->set_subdomain_id (0);
@@ -1489,7 +1510,7 @@ next_cell:
     // finally loop over all cells and set the
     // subdomain ids
     unsigned int index = 0;
-    for (typename Triangulation<dim,spacedim>::active_cell_iterator
+    for (typename dealii::internal::ActiveCellIterator<dim, spacedim, Triangulation<dim, spacedim> >::type
          cell = triangulation.begin_active();
          cell != triangulation.end(); ++cell, ++index)
       cell->set_subdomain_id (partition_indices[index]);
@@ -1547,7 +1568,7 @@ next_cell:
     if (const parallel::distributed::Triangulation<dim,spacedim> *tr
         = dynamic_cast<const parallel::distributed::Triangulation<dim,spacedim> *>
         (&triangulation))
-      for (typename Triangulation<dim,spacedim>::active_cell_iterator
+      for (typename dealii::internal::ActiveCellIterator<dim, spacedim, Triangulation<dim, spacedim> >::type
            cell = triangulation.begin_active();
            cell != triangulation.end(); ++cell)
         if (cell->is_artificial()
index 1c2ee8d9ec836b5890e7cef96d363db986aeb6f5..6749419546d6c172116842153d60b3f76eb2d569 100644 (file)
@@ -37,16 +37,15 @@ for (X : TRIANGULATION_AND_DOFHANDLERS; deal_II_dimension : DIMENSIONS ; deal_II
                         const Point<deal_II_space_dimension> &);
 
   template
-    std::vector<X::active_cell_iterator>
-    find_cells_adjacent_to_vertex(const X &,
-                                 const unsigned int);
+    std::vector<dealii::internal::ActiveCellIterator<deal_II_dimension, deal_II_space_dimension, X>::type>
+    find_cells_adjacent_to_vertex(const X &, const unsigned int);
+
   template
-    X::active_cell_iterator
-    find_active_cell_around_point (const X &,
-                                  const Point<deal_II_space_dimension> &p);
+    dealii::internal::ActiveCellIterator<deal_II_dimension, deal_II_space_dimension, X>::type
+    find_active_cell_around_point (const X &, const Point<deal_II_space_dimension> &p);
 
   template
-    std::pair<X::active_cell_iterator, Point<deal_II_dimension> >
+    std::pair<dealii::internal::ActiveCellIterator<deal_II_dimension, deal_II_space_dimension, X>::type, Point<deal_II_dimension> >
     find_active_cell_around_point (const Mapping<deal_II_dimension, deal_II_space_dimension> &,
                                   const X &,
                                   const Point<deal_II_space_dimension> &);
@@ -81,16 +80,16 @@ for (deal_II_dimension : DIMENSIONS ; deal_II_space_dimension : SPACE_DIMENSIONS
                         const Point<deal_II_space_dimension> &);
 
   template
-    std::vector<parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension>::active_cell_iterator>
+    std::vector<dealii::internal::ActiveCellIterator<deal_II_dimension, deal_II_space_dimension, parallel::distributed::Triangulation<deal_II_dimension, deal_II_space_dimension> >::type>
     find_cells_adjacent_to_vertex(const parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension> &,
                                  const unsigned int);
   template
-    parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension>::active_cell_iterator
+    dealii::internal::ActiveCellIterator<deal_II_dimension, deal_II_space_dimension, parallel::distributed::Triangulation<deal_II_dimension, deal_II_space_dimension> >::type
     find_active_cell_around_point (const parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension> &,
                                   const Point<deal_II_space_dimension> &p);
 
   template
-    std::pair<parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension>::active_cell_iterator, Point<deal_II_dimension> >
+    std::pair<dealii::internal::ActiveCellIterator<deal_II_dimension, deal_II_space_dimension, parallel::distributed::Triangulation<deal_II_dimension, deal_II_space_dimension> >::type, Point<deal_II_dimension> >
     find_active_cell_around_point (const Mapping<deal_II_dimension, deal_II_space_dimension> &,
                                   const parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension> &,
                                   const Point<deal_II_space_dimension> &);
@@ -116,12 +115,12 @@ for (deal_II_dimension : DIMENSIONS ; deal_II_space_dimension : SPACE_DIMENSIONS
 for (deal_II_space_dimension : SPACE_DIMENSIONS)
 {
 
-    parallel::distributed::Triangulation<deal_II_space_dimension>::active_cell_iterator
+    dealii::internal::ActiveCellIterator<deal_II_space_dimension, deal_II_space_dimension, parallel::distributed::Triangulation<deal_II_space_dimension, deal_II_space_dimension> >::type
     find_active_cell_around_point (const parallel::distributed::Triangulation<deal_II_space_dimension> &,
                                   const Point<deal_II_space_dimension> &p);
 
 
-    std::pair<parallel::distributed::Triangulation<deal_II_space_dimension>::active_cell_iterator, Point<deal_II_space_dimension> >
+    std::pair<dealii::internal::ActiveCellIterator<deal_II_space_dimension, deal_II_space_dimension, parallel::distributed::Triangulation<deal_II_space_dimension, deal_II_space_dimension> >::type, Point<deal_II_space_dimension> >
     find_active_cell_around_point (const Mapping<deal_II_space_dimension> &,
                                   const parallel::distributed::Triangulation<deal_II_space_dimension> &,
                                   const Point<deal_II_space_dimension> &);
index 23e568a545089dd5183764600a2845acd4cac3d5..4f3500a8002f45bfea8ea1f898cfa4b4f81a62c4 100644 (file)
@@ -753,17 +753,16 @@ namespace DerivativeApproximation
   namespace internal
   {
     /**
-     * Compute the derivative approximation on one cell. This computes the full
-     * derivative tensor.
-     */
-    template <class DerivativeDescription, int dim,
-              template <int, int> class DH, class InputVector, int spacedim>
+    * Compute the derivative approximation on one cell. This computes the full
+    * derivative tensor.
+    */
+    template <class DerivativeDescription, int dim, template <int, int> class DH, class InputVector, int spacedim>
     void
     approximate_cell (const Mapping<dim,spacedim>                   &mapping,
                       const DH<dim,spacedim>                        &dof_handler,
                       const InputVector                             &solution,
                       const unsigned int                             component,
-                      const typename DH<dim,spacedim>::active_cell_iterator  &cell,
+                      const TriaActiveIterator < dealii::DoFCellAccessor < DH < dim, spacedim >, false > >  &cell,
                       typename DerivativeDescription::Derivative    &derivative)
     {
       QMidpoint<dim> midpoint_rule;
@@ -792,7 +791,8 @@ namespace DerivativeApproximation
       // active neighbors of a cell
       // reserve the maximal number of
       // active neighbors
-      std::vector<typename DH<dim,spacedim>::active_cell_iterator> active_neighbors;
+      std::vector<TriaActiveIterator < dealii::DoFCellAccessor < DH < dim, spacedim >, false > > > active_neighbors;
+
       active_neighbors.reserve (GeometryInfo<dim>::faces_per_cell *
                                 GeometryInfo<dim>::max_children_per_face);
 
@@ -834,11 +834,11 @@ namespace DerivativeApproximation
       // now loop over all active
       // neighbors and collect the
       // data we need
-      typename std::vector<typename DH<dim,spacedim>::active_cell_iterator>::const_iterator
+      typename std::vector<TriaActiveIterator < dealii::DoFCellAccessor < DH < dim, spacedim >, false > > >::const_iterator
       neighbor_ptr = active_neighbors.begin();
       for (; neighbor_ptr!=active_neighbors.end(); ++neighbor_ptr)
         {
-          const typename DH<dim,spacedim>::active_cell_iterator
+          const TriaActiveIterator < dealii::DoFCellAccessor < DH < dim, spacedim >, false > >
           neighbor = *neighbor_ptr;
 
           // reinit fe values object...
@@ -929,11 +929,12 @@ namespace DerivativeApproximation
     template <class DerivativeDescription, int dim,
               template <int, int> class DH, class InputVector, int spacedim>
     void
-    approximate (SynchronousIterators<std_cxx11::tuple<typename DH<dim,spacedim>::active_cell_iterator,Vector<float>::iterator> > const &cell,
-                 const Mapping<dim,spacedim>                  &mapping,
-                 const DH<dim,spacedim>                       &dof_handler,
-                 const InputVector                            &solution,
-                 const unsigned int                            component)
+    approximate(
+      SynchronousIterators<std_cxx11::tuple<TriaActiveIterator < dealii::DoFCellAccessor < DH < dim, spacedim >, false > >, Vector<float>::iterator> > const &cell,
+      const Mapping<dim,spacedim>                  &mapping,
+      const DH<dim,spacedim>                       &dof_handler,
+      const InputVector                            &solution,
+      const unsigned int                            component)
     {
       // if the cell is not locally owned, then there is nothing to do
       if (std_cxx11::get<0>(cell.iterators)->is_locally_owned() == false)
@@ -943,8 +944,9 @@ namespace DerivativeApproximation
           typename DerivativeDescription::Derivative derivative;
           // call the function doing the actual
           // work on this cell
-          approximate_cell<DerivativeDescription,dim,DH,InputVector>
+          approximate_cell<DerivativeDescription,dim,DH,InputVector, spacedim>
           (mapping,dof_handler,solution,component,std_cxx11::get<0>(cell.iterators),derivative);
+
           // evaluate the norm and fill the vector
           //*derivative_norm_on_this_cell
           *std_cxx11::get<1>(cell.iterators) = DerivativeDescription::derivative_norm (derivative);
@@ -977,7 +979,7 @@ namespace DerivativeApproximation
       Assert (component < dof_handler.get_fe().n_components(),
               ExcIndexRange (component, 0, dof_handler.get_fe().n_components()));
 
-      typedef std_cxx11::tuple<typename DH<dim,spacedim>::active_cell_iterator,Vector<float>::iterator>
+      typedef std_cxx11::tuple<TriaActiveIterator < dealii::DoFCellAccessor < DH < dim, spacedim >, false > >, Vector<float>::iterator>
       Iterators;
       SynchronousIterators<Iterators> begin(Iterators(dof_handler.begin_active(),
                                                       derivative_norm.begin())),
@@ -1076,9 +1078,9 @@ namespace DerivativeApproximation
                                 const DH                    &dof,
                                 const InputVector                            &solution,
 #ifndef _MSC_VER
-                                const typename DH::active_cell_iterator      &cell,
+                                const typename DH::active_cell_iterator &cell,
 #else
-                                const TriaActiveIterator < dealii::DoFCellAccessor < DH, false > >      &cell,
+                                const TriaActiveIterator < dealii::DoFCellAccessor < DH, false > > &cell,
 #endif
                                 Tensor<order, dim>  &derivative,
                                 const unsigned int                            component)
@@ -1099,9 +1101,9 @@ namespace DerivativeApproximation
   approximate_derivative_tensor(const DH                     &dof,
                                 const InputVector                            &solution,
 #ifndef _MSC_VER
-                                const typename DH::active_cell_iterator      &cell,
+                                const typename DH::active_cell_iterator &cell,
 #else
-                                const TriaActiveIterator < dealii::DoFCellAccessor < DH, false > >      &cell,
+                                const TriaActiveIterator < dealii::DoFCellAccessor < DH, false > > &cell,
 #endif
                                 Tensor<order, dim>                  &derivative,
                                 const unsigned int                            component)

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.