]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use VectorType, not Vector, as a template type.
authorDavid Wells <wellsd2@rpi.edu>
Mon, 19 Oct 2015 01:29:18 +0000 (21:29 -0400)
committerDavid Wells <wellsd2@rpi.edu>
Mon, 2 Nov 2015 16:21:24 +0000 (11:21 -0500)
Writing something like

template<typename Vector>

is ambiguous because Vector is also a class.

include/deal.II/base/index_set.h
include/deal.II/distributed/grid_refinement.h
include/deal.II/grid/grid_refinement.h
include/deal.II/numerics/vector_tools.templates.h
source/distributed/grid_refinement.cc
source/grid/grid_refinement.cc
tests/gla/extract_subvector_to.cc
tests/gla/extract_subvector_to_parallel.cc
tests/lac/block_matrices_04.cc
tests/lapack/solver_cg.cc

index 74de38b1a21204b47268d8b6457ad3796757d9d5..b0dc5a4016fa05a4ad2d81f99cf00eb9d3195eee 100644 (file)
@@ -295,8 +295,8 @@ public:
    * Vector, BlockVector, but also std::vector@<bool@>, std::vector@<int@>,
    * and std::vector@<double@>.
    */
-  template <typename Vector>
-  void fill_binary_vector (Vector &vector) const;
+  template <typename VectorType>
+  void fill_binary_vector (VectorType &vector) const;
 
   /**
    * Outputs a text representation of this IndexSet to the given stream. Used
index e26a8f0ad89f2a6d28644e830c2b4b24b1828719..ff774b1bc965847af26741106a672eb38dc4ed83 100644 (file)
@@ -65,14 +65,14 @@ namespace parallel
        *
        * The same is true for the fraction of cells that is coarsened.
        */
-      template <int dim, class Vector, int spacedim>
+      template <int dim, class VectorType, int spacedim>
       void
-      refine_and_coarsen_fixed_number (
-        parallel::distributed::Triangulation<dim,spacedim> &tria,
-        const Vector                &criteria,
-        const double                 top_fraction_of_cells,
-        const double                 bottom_fraction_of_cells,
-        const unsigned int           max_n_cells = std::numeric_limits<unsigned int>::max());
+      refine_and_coarsen_fixed_number
+      (parallel::distributed::Triangulation<dim,spacedim> &tria,
+       const VectorType                                   &criteria,
+       const double                                       top_fraction_of_cells,
+       const double                                       bottom_fraction_of_cells,
+       const unsigned int                                 max_n_cells = std::numeric_limits<unsigned int>::max());
 
       /**
        * Like dealii::GridRefinement::refine_and_coarsen_fixed_fraction, but
@@ -94,13 +94,13 @@ namespace parallel
        *
        * The same is true for the fraction of cells that is coarsened.
        */
-      template <int dim, class Vector, int spacedim>
+      template <int dim, class VectorType, int spacedim>
       void
-      refine_and_coarsen_fixed_fraction (
-        parallel::distributed::Triangulation<dim,spacedim> &tria,
-        const Vector                &criteria,
-        const double                top_fraction_of_error,
-        const double                bottom_fraction_of_error);
+      refine_and_coarsen_fixed_fraction
+      (parallel::distributed::Triangulation<dim,spacedim> &tria,
+       const VectorType                                   &criteria,
+       const double                                       top_fraction_of_error,
+       const double                                       bottom_fraction_of_error);
     }
   }
 }
index 318cbb785ea04f330e1136a754ef09253a9ff1e1..3d47f705bbc75776bd6b838b100ad7289f76f747 100644 (file)
@@ -146,14 +146,14 @@ namespace GridRefinement
    * indicator. The default value of this argument is to impose no limit on
    * the number of cells.
    */
-  template <int dim, class Vector, int spacedim>
+  template <int dim, class VectorType, int spacedim>
   void
-  refine_and_coarsen_fixed_number (
-    Triangulation<dim,spacedim> &tria,
-    const Vector                &criteria,
-    const double                top_fraction_of_cells,
-    const double                bottom_fraction_of_cells,
-    const unsigned int          max_n_cells = std::numeric_limits<unsigned int>::max());
+  refine_and_coarsen_fixed_number
+  (Triangulation<dim,spacedim> &tria,
+   const VectorType            &criteria,
+   const double                top_fraction_of_cells,
+   const double                bottom_fraction_of_cells,
+   const unsigned int          max_n_cells = std::numeric_limits<unsigned int>::max());
 
   /**
    * This function provides a refinement strategy controlling the reduction of
@@ -207,14 +207,14 @@ namespace GridRefinement
    * indicator. The default value of this argument is to impose no limit on
    * the number of cells.
    */
-  template <int dim, class Vector, int spacedim>
+  template <int dim, class VectorType, int spacedim>
   void
-  refine_and_coarsen_fixed_fraction (
-    Triangulation<dim,spacedim> &tria,
-    const Vector                &criteria,
-    const double                top_fraction,
-    const double                bottom_fraction,
-    const unsigned int          max_n_cells = std::numeric_limits<unsigned int>::max());
+  refine_and_coarsen_fixed_fraction
+  (Triangulation<dim,spacedim> &tria,
+   const VectorType            &criteria,
+   const double                top_fraction,
+   const double                bottom_fraction,
+   const unsigned int          max_n_cells = std::numeric_limits<unsigned int>::max());
 
 
 
@@ -290,11 +290,11 @@ namespace GridRefinement
    * thesis, University of Heidelberg, 2005. See in particular Section 4.3,
    * pp. 42-43.
    */
-  template <int dim, class Vector, int spacedim>
+  template <int dim, class VectorType, int spacedim>
   void
   refine_and_coarsen_optimize (Triangulation<dim,spacedim> &tria,
-                               const Vector                &criteria,
-                               const unsigned int           order=2);
+                               const VectorType            &criteria,
+                               const unsigned int          order=2);
 
   /**
    * Flag all mesh cells for which the value in @p criteria exceeds @p
@@ -310,11 +310,11 @@ namespace GridRefinement
    * This function does not implement a refinement strategy, it is more a
    * helper function for the actual strategies.
    */
-  template <int dim, class Vector, int spacedim>
+  template <int dim, class VectorType, int spacedim>
   void refine (Triangulation<dim,spacedim> &tria,
-               const Vector                &criteria,
+               const VectorType            &criteria,
                const double                threshold,
-               const unsigned int max_to_mark = numbers::invalid_unsigned_int);
+               const unsigned int          max_to_mark = numbers::invalid_unsigned_int);
 
   /**
    * Flag all mesh cells for which the value in @p criteria is less than @p
@@ -330,9 +330,9 @@ namespace GridRefinement
    * This function does not implement a refinement strategy, it is more a
    * helper function for the actual strategies.
    */
-  template <int dim, class Vector, int spacedim>
+  template <int dim, class VectorType, int spacedim>
   void coarsen (Triangulation<dim,spacedim> &tria,
-                const Vector                &criteria,
+                const VectorType            &criteria,
                 const double                threshold);
 
   /**
index eb8995710122d37b84ea9e23f94d140333dd6b8f..928a9b3acc8bbeb81b104bc24b8f004ee66bb2bf 100644 (file)
@@ -543,12 +543,12 @@ namespace VectorTools
 
   template <int dim, int spacedim,
             template <int,int> class DH,
-            class Vector>
+            class VectorType>
   void
   interpolate_to_different_mesh (const DH<dim, spacedim> &dof1,
-                                 const Vector            &u1,
+                                 const VectorType        &u1,
                                  const DH<dim, spacedim> &dof2,
-                                 Vector                  &u2)
+                                 VectorType              &u2)
   {
     Assert(GridTools::have_same_coarse_mesh(dof1, dof2),
            ExcMessage ("The two containers must represent triangulations that "
@@ -567,13 +567,13 @@ namespace VectorTools
 
   template <int dim, int spacedim,
             template <int,int> class DH,
-            class Vector>
+            class VectorType>
   void
   interpolate_to_different_mesh (const DH<dim, spacedim> &dof1,
-                                 const Vector            &u1,
+                                 const VectorType        &u1,
                                  const DH<dim, spacedim> &dof2,
                                  const ConstraintMatrix  &constraints,
-                                 Vector                  &u2)
+                                 VectorType              &u2)
   {
     Assert(GridTools::have_same_coarse_mesh(dof1, dof2),
            ExcMessage ("The two containers must represent triangulations that "
@@ -607,12 +607,12 @@ namespace VectorTools
 
   template <int dim, int spacedim,
             template <int,int> class DH,
-            class Vector>
+            class VectorType>
   void
   interpolate_to_different_mesh (const InterGridMap<DH<dim, spacedim> > &intergridmap,
-                                 const Vector                           &u1,
+                                 const VectorType                       &u1,
                                  const ConstraintMatrix                 &constraints,
-                                 Vector                                 &u2)
+                                 VectorType                             &u2)
   {
     const DH<dim, spacedim> &dof1 = intergridmap.get_source_grid();
     const DH<dim, spacedim> &dof2 = intergridmap.get_destination_grid();
@@ -750,7 +750,7 @@ namespace VectorTools
      * Generic implementation of the project() function
      */
     template <int dim, int spacedim,
-              class Vector,
+              class VectorType,
               template <int,int> class DH,
               template <int,int> class M_or_MC,
               template <int> class Q_or_QC>
@@ -759,7 +759,7 @@ namespace VectorTools
                      const ConstraintMatrix   &constraints,
                      const Q_or_QC<dim>       &quadrature,
                      const Function<spacedim> &function,
-                     Vector                   &vec_result,
+                     VectorType               &vec_result,
                      const bool                enforce_zero_boundary,
                      const Q_or_QC<dim-1>  &q_boundary,
                      const bool                project_to_boundary_first)
index 7325bc9e0ed8a797f80aa30f842c03e7659335f1..81c6f0f6d20f85036b4116fd28029fdd816c5ec7 100644 (file)
@@ -149,11 +149,11 @@ namespace
    * or not), extract those that pertain to
    * locally owned cells.
    */
-  template <int dim, int spacedim, class Vector>
+  template <int dim, int spacedim, typename VectorType>
   void
   get_locally_owned_indicators (const parallel::distributed::Triangulation<dim,spacedim> &tria,
-                                const Vector &criteria,
-                                dealii::Vector<typename Vector::value_type> &locally_owned_indicators)
+                                const VectorType &criteria,
+                                dealii::Vector<typename VectorType::value_type> &locally_owned_indicators)
   {
     Assert (locally_owned_indicators.size() == tria.n_locally_owned_active_cells(),
             ExcInternalError());
@@ -222,12 +222,12 @@ namespace
    * locally own as appropriate for
    * coarsening or refinement.
    */
-  template <int dim, int spacedim, class Vector>
+  template <int dim, int spacedim, typename VectorType>
   void
   mark_cells (parallel::distributed::Triangulation<dim,spacedim> &tria,
-              const Vector &criteria,
-              const double top_threshold,
-              const double bottom_threshold)
+              const VectorType                                   &criteria,
+              const double                                        top_threshold,
+              const double                                        bottom_threshold)
   {
     dealii::GridRefinement::refine (tria, criteria, top_threshold);
     dealii::GridRefinement::coarsen (tria, criteria, bottom_threshold);
@@ -451,14 +451,14 @@ namespace parallel
   {
     namespace GridRefinement
     {
-      template <int dim, class Vector, int spacedim>
+      template <int dim, typename VectorType, int spacedim>
       void
-      refine_and_coarsen_fixed_number (
-        parallel::distributed::Triangulation<dim,spacedim> &tria,
-        const Vector                &criteria,
-        const double                 top_fraction_of_cells,
-        const double                 bottom_fraction_of_cells,
-        const unsigned int           max_n_cells)
+      refine_and_coarsen_fixed_number
+      (parallel::distributed::Triangulation<dim,spacedim> &tria,
+       const VectorType                                   &criteria,
+       const double                                        top_fraction_of_cells,
+       const double                                        bottom_fraction_of_cells,
+       const unsigned int                                  max_n_cells)
       {
         Assert (criteria.size() == tria.n_active_cells(),
                 ExcDimensionMismatch (criteria.size(), tria.n_active_cells()));
@@ -482,7 +482,7 @@ namespace parallel
         // vector of indicators the
         // ones that correspond to
         // cells that we locally own
-        dealii::Vector<typename Vector::value_type>
+        dealii::Vector<typename VectorType::value_type>
         locally_owned_indicators (tria.n_locally_owned_active_cells());
         get_locally_owned_indicators (tria,
                                       criteria,
@@ -496,7 +496,7 @@ namespace parallel
         // need it here, but it's a
         // collective communication
         // call
-        const std::pair<typename Vector::value_type,typename Vector::value_type> global_min_and_max
+        const std::pair<typename VectorType::value_type,typename VectorType::value_type> global_min_and_max
           = compute_global_min_and_max_at_root (locally_owned_indicators,
                                                 mpi_communicator);
 
@@ -538,13 +538,13 @@ namespace parallel
       }
 
 
-      template <int dim, class Vector, int spacedim>
+      template <int dim, typename VectorType, int spacedim>
       void
-      refine_and_coarsen_fixed_fraction (
-        parallel::distributed::Triangulation<dim,spacedim> &tria,
-        const Vector                &criteria,
-        const double                top_fraction_of_error,
-        const double                bottom_fraction_of_error)
+      refine_and_coarsen_fixed_fraction
+      (parallel::distributed::Triangulation<dim,spacedim> &tria,
+       const VectorType                                   &criteria,
+       const double                                        top_fraction_of_error,
+       const double                                        bottom_fraction_of_error)
       {
         Assert (criteria.size() == tria.n_active_cells(),
                 ExcDimensionMismatch (criteria.size(), tria.n_active_cells()));
@@ -561,7 +561,7 @@ namespace parallel
         // vector of indicators the
         // ones that correspond to
         // cells that we locally own
-        dealii::Vector<typename Vector::value_type>
+        dealii::Vector<typename VectorType::value_type>
         locally_owned_indicators (tria.n_locally_owned_active_cells());
         get_locally_owned_indicators (tria,
                                       criteria,
index 8505d66fe9a6b1a751b9fc7b96dad4aa5a8b7361..f93af05734836c69eaed0030dc0421027c221ad1 100644 (file)
@@ -118,30 +118,30 @@ namespace
   } /* namespace internal */
 
 
-  template <typename Vector>
-  typename constraint_and_return_value<!IsBlockVector<Vector>::value,
-           typename Vector::value_type>::type
-           min_element (const Vector &criteria)
+  template <typename VectorType>
+  typename constraint_and_return_value<!IsBlockVector<VectorType>::value,
+           typename VectorType::value_type>::type
+           min_element (const VectorType &criteria)
   {
     return internal::min_element (criteria);
   }
 
 
-  template <typename Vector>
-  typename constraint_and_return_value<!IsBlockVector<Vector>::value,
-           typename Vector::value_type>::type
-           max_element (const Vector &criteria)
+  template <typename VectorType>
+  typename constraint_and_return_value<!IsBlockVector<VectorType>::value,
+           typename VectorType::value_type>::type
+           max_element (const VectorType &criteria)
   {
     return internal::max_element (criteria);
   }
 
 
-  template <typename Vector>
-  typename constraint_and_return_value<IsBlockVector<Vector>::value,
-           typename Vector::value_type>::type
-           min_element (const Vector &criteria)
+  template <typename VectorType>
+  typename constraint_and_return_value<IsBlockVector<VectorType>::value,
+           typename VectorType::value_type>::type
+           min_element (const VectorType &criteria)
   {
-    typename Vector::value_type t = internal::min_element(criteria.block(0));
+    typename VectorType::value_type t = internal::min_element(criteria.block(0));
     for (unsigned int b=1; b<criteria.n_blocks(); ++b)
       t = std::min (t, internal::min_element(criteria.block(b)));
 
@@ -149,12 +149,12 @@ namespace
   }
 
 
-  template <typename Vector>
-  typename constraint_and_return_value<IsBlockVector<Vector>::value,
-           typename Vector::value_type>::type
-           max_element (const Vector &criteria)
+  template <typename VectorType>
+  typename constraint_and_return_value<IsBlockVector<VectorType>::value,
+           typename VectorType::value_type>::type
+           max_element (const VectorType &criteria)
   {
-    typename Vector::value_type t = internal::max_element(criteria.block(0));
+    typename VectorType::value_type t = internal::max_element(criteria.block(0));
     for (unsigned int b=1; b<criteria.n_blocks(); ++b)
       t = std::max (t, internal::max_element(criteria.block(b)));
 
@@ -172,14 +172,14 @@ namespace
    * library version and is needed in @p refine_and_coarsen_optimize.
    */
 
-  template <class Vector>
-  void qsort_index (const Vector              &a,
+  template <class VectorType>
+  void qsort_index (const VectorType          &a,
                     std::vector<unsigned int> &ind,
                     int                        l,
                     int                        r)
   {
     int i,j;
-    typename Vector::value_type v;
+    typename VectorType::value_type v;
 
     if (r<=l)
       return;
@@ -214,9 +214,9 @@ namespace
 
 
 
-template <int dim, class Vector, int spacedim>
+template <int dim, class VectorType, int spacedim>
 void GridRefinement::refine (Triangulation<dim,spacedim> &tria,
-                             const Vector       &criteria,
+                             const VectorType   &criteria,
                              const double        threshold,
                              const unsigned int max_to_mark)
 {
@@ -261,10 +261,10 @@ void GridRefinement::refine (Triangulation<dim,spacedim> &tria,
 
 
 
-template <int dim, class Vector, int spacedim>
+template <int dim, class VectorType, int spacedim>
 void GridRefinement::coarsen (Triangulation<dim,spacedim> &tria,
-                              const Vector       &criteria,
-                              const double        threshold)
+                              const VectorType            &criteria,
+                              const double                 threshold)
 {
   Assert (criteria.size() == tria.n_active_cells(),
           ExcDimensionMismatch(criteria.size(), tria.n_active_cells()));
@@ -359,13 +359,13 @@ GridRefinement::adjust_refine_and_coarsen_number_fraction (const unsigned int  c
   return (adjusted_fractions);
 }
 
-template <int dim, class Vector, int spacedim>
+template <int dim, class VectorType, int spacedim>
 void
 GridRefinement::refine_and_coarsen_fixed_number (Triangulation<dim,spacedim> &tria,
-                                                 const Vector       &criteria,
-                                                 const double        top_fraction,
-                                                 const double        bottom_fraction,
-                                                 const unsigned int  max_n_cells)
+                                                 const VectorType            &criteria,
+                                                 const double                 top_fraction,
+                                                 const double                 bottom_fraction,
+                                                 const unsigned int           max_n_cells)
 {
   // correct number of cells is
   // checked in @p{refine}
@@ -385,7 +385,7 @@ GridRefinement::refine_and_coarsen_fixed_number (Triangulation<dim,spacedim> &tr
 
   if (refine_cells || coarsen_cells)
     {
-      dealii::Vector<typename Vector::value_type> tmp (criteria);
+      dealii::Vector<typename VectorType::value_type> tmp (criteria);
       if (refine_cells)
         {
           std::nth_element (tmp.begin(), tmp.begin()+refine_cells,
@@ -407,10 +407,10 @@ GridRefinement::refine_and_coarsen_fixed_number (Triangulation<dim,spacedim> &tr
 
 
 
-template <int dim, class Vector, int spacedim>
+template <int dim, typename VectorType, int spacedim>
 void
 GridRefinement::refine_and_coarsen_fixed_fraction (Triangulation<dim,spacedim> &tria,
-                                                   const Vector       &criteria,
+                                                   const VectorType   &criteria,
                                                    const double        top_fraction,
                                                    const double        bottom_fraction,
                                                    const unsigned int  max_n_cells)
@@ -426,7 +426,7 @@ GridRefinement::refine_and_coarsen_fixed_fraction (Triangulation<dim,spacedim> &
   // error, which is what we have to sum
   // up and compare with
   // @p{fraction_of_error*total_error}.
-  dealii::Vector<typename Vector::value_type> tmp;
+  dealii::Vector<typename VectorType::value_type> tmp;
   tmp = criteria;
   const double total_error = tmp.l1_norm();
 
@@ -435,7 +435,7 @@ GridRefinement::refine_and_coarsen_fixed_fraction (Triangulation<dim,spacedim> &
   std::sort (tmp.begin(), tmp.end(), std::greater<double>());
 
   // compute thresholds
-  typename dealii::Vector<typename Vector::value_type>::const_iterator
+  typename dealii::Vector<typename VectorType::value_type>::const_iterator
   pp=tmp.begin();
   for (double sum=0;
        (sum<top_fraction*total_error) && (pp!=(tmp.end()-1));
@@ -444,7 +444,7 @@ GridRefinement::refine_and_coarsen_fixed_fraction (Triangulation<dim,spacedim> &
   double top_threshold = ( pp != tmp.begin () ?
                            (*pp+*(pp-1))/2 :
                            *pp );
-  typename dealii::Vector<typename Vector::value_type>::const_iterator
+  typename dealii::Vector<typename VectorType::value_type>::const_iterator
   qq=(tmp.end()-1);
   for (double sum=0;
        (sum<bottom_fraction*total_error) && (qq!=tmp.begin());
@@ -537,11 +537,11 @@ GridRefinement::refine_and_coarsen_fixed_fraction (Triangulation<dim,spacedim> &
 
 
 
-template <int dim, class Vector, int spacedim>
+template <int dim, typename VectorType, int spacedim>
 void
 GridRefinement::refine_and_coarsen_optimize (Triangulation<dim,spacedim> &tria,
-                                             const Vector       &criteria,
-                                             const unsigned int  order)
+                                             const VectorType            &criteria,
+                                             const unsigned int           order)
 {
   Assert (criteria.size() == tria.n_active_cells(),
           ExcDimensionMismatch(criteria.size(), tria.n_active_cells()));
@@ -586,4 +586,3 @@ GridRefinement::refine_and_coarsen_optimize (Triangulation<dim,spacedim> &tria,
 #include "grid_refinement.inst"
 
 DEAL_II_NAMESPACE_CLOSE
-
index 5ebf4ec93ae2299f297f06cd3de0dfa875bfbeda..743c15619892983985f3ed3ffac22db898c8d91b 100644 (file)
@@ -29,8 +29,8 @@
 #include <vector>
 
 
-template <class Vector>
-void test (Vector &vector)
+template <typename VectorType>
+void test (VectorType &vector)
 {
   const unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
 
@@ -38,19 +38,19 @@ void test (Vector &vector)
     vector(i) = i;
 
   // select every other element
-  std::vector<typename Vector::size_type> indices;
+  std::vector<typename VectorType::size_type> indices;
   for (unsigned int j=0; j<vector.size()/2; ++j)
     indices.push_back (2*j);
 
   // do the extraction with the function that takes indices, then
   // assert correctness
-  std::vector<typename Vector::value_type> values1 (indices.size());
+  std::vector<typename VectorType::value_type> values1 (indices.size());
   vector.extract_subvector_to (indices, values1);
   for (unsigned int j=0; j<vector.size()/2; ++j)
     AssertThrow (values1[j] == 2*j, ExcInternalError());
 
   // do the same with the version of the function that takes iterators
-  std::vector<typename Vector::value_type> values2 (indices.size());
+  std::vector<typename VectorType::value_type> values2 (indices.size());
   vector.extract_subvector_to (indices.begin(),
                                indices.end(),
                                values2.begin());
index eddef20bda8888ebf8fa9652469ee1a21e260725..a231ba956fe03c41cee55accd03a5c7eb1d5d81b 100644 (file)
@@ -28,8 +28,8 @@
 #include <vector>
 
 
-template <class Vector>
-void set (Vector &vector)
+template <typename VectorType>
+void set (VectorType &vector)
 {
   for (unsigned int i=0; i<vector.size(); ++i)
     if (vector.locally_owned_elements().is_element(i))
@@ -38,25 +38,25 @@ void set (Vector &vector)
 }
 
 
-template <class Vector>
-void test (Vector &vector)
+template <typename VectorType>
+void test (VectorType &vector)
 {
   const unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
 
   // select every other element
-  std::vector<typename Vector::size_type> indices;
+  std::vector<typename VectorType::size_type> indices;
   for (unsigned int j=0; j<vector.size()/2; ++j)
     indices.push_back (2*j);
 
   // do the extraction with the function that takes indices, then
   // assert correctness
-  std::vector<typename Vector::value_type> values1 (indices.size());
+  std::vector<typename VectorType::value_type> values1 (indices.size());
   vector.extract_subvector_to (indices, values1);
   for (unsigned int j=0; j<vector.size()/2; ++j)
     AssertThrow (values1[j] == 2*j, ExcInternalError());
 
   // do the same with the version of the function that takes iterators
-  std::vector<typename Vector::value_type> values2 (indices.size());
+  std::vector<typename VectorType::value_type> values2 (indices.size());
   vector.extract_subvector_to (indices.begin(),
                                indices.end(),
                                values2.begin());
index c036f373195047267fa63f786c410bf74553d4fd..5aad98e6ab3025adb18d5c6c06a5dc13a0890da0 100644 (file)
@@ -54,7 +54,7 @@
 #include <fstream>
 
 
-template <class Vector, class Matrix, class Sparsity>
+template <typename VectorType, class Matrix, class Sparsity>
 class LaplaceProblem
 {
 public:
@@ -64,7 +64,7 @@ public:
   void reinit_sparsity ();
   void reinit_vectors ();
 
-  Vector       solution;
+  VectorType solution;
 
 private:
   void make_grid_and_dofs ();
@@ -73,19 +73,19 @@ private:
 
   const unsigned int n_blocks;
 
-  Triangulation<2>     triangulation;
-  FE_Q<2>              fe;
-  DoFHandler<2>        dof_handler;
+  Triangulation<2>   triangulation;
+  FE_Q<2>            fe;
+  DoFHandler<2>      dof_handler;
 
-  Sparsity      sparsity_pattern;
-  Matrix        system_matrix;
+  Sparsity           sparsity_pattern;
+  Matrix             system_matrix;
 
-  Vector       system_rhs;
+  VectorType         system_rhs;
 };
 
 
-template <class Vector, class Matrix, class Sparsity>
-LaplaceProblem<Vector,Matrix,Sparsity>::LaplaceProblem (const unsigned int n_blocks) :
+template <typename VectorType, class Matrix, class Sparsity>
+LaplaceProblem<VectorType,Matrix,Sparsity>::LaplaceProblem (const unsigned int n_blocks) :
   n_blocks (n_blocks),
   fe(1),
   dof_handler (triangulation)
@@ -113,8 +113,8 @@ LaplaceProblem<Vector<float>,SparseMatrix<float>,SparsityPattern>::LaplaceProble
 
 
 
-template <class Vector, class Matrix, class Sparsity>
-void LaplaceProblem<Vector,Matrix,Sparsity>::make_grid_and_dofs ()
+template <typename VectorType, class Matrix, class Sparsity>
+void LaplaceProblem<VectorType,Matrix,Sparsity>::make_grid_and_dofs ()
 {
   GridGenerator::hyper_cube (triangulation, -1, 1);
   triangulation.refine_global (3);
@@ -256,8 +256,8 @@ void LaplaceProblem<BlockVector<double>,BlockSparseMatrix<double>,BlockSparsityP
 
 
 
-template <class Vector, class Matrix, class Sparsity>
-void LaplaceProblem<Vector,Matrix,Sparsity>::assemble_system ()
+template <typename VectorType, class Matrix, class Sparsity>
+void LaplaceProblem<VectorType,Matrix,Sparsity>::assemble_system ()
 {
   QGauss<2>  quadrature_formula(2);
   FEValues<2> fe_values (fe, quadrature_formula,
@@ -320,12 +320,12 @@ void LaplaceProblem<Vector,Matrix,Sparsity>::assemble_system ()
 }
 
 
-template <class Vector, class Matrix, class Sparsity>
-void LaplaceProblem<Vector,Matrix,Sparsity>::solve ()
+template <typename VectorType, class Matrix, class Sparsity>
+void LaplaceProblem<VectorType,Matrix,Sparsity>::solve ()
 {
-  SolverControl           solver_control (1000, 1e-12, false, false);
-  PrimitiveVectorMemory<Vector> vector_memory;
-  SolverCG<Vector>        cg (solver_control, vector_memory);
+  SolverControl                     solver_control (1000, 1e-12, false, false);
+  PrimitiveVectorMemory<VectorType> vector_memory;
+  SolverCG<VectorType>              cg (solver_control, vector_memory);
 
   PreconditionJacobi<Matrix> preconditioner;
   preconditioner.initialize (system_matrix, 0.8);
@@ -335,8 +335,8 @@ void LaplaceProblem<Vector,Matrix,Sparsity>::solve ()
 }
 
 
-template <class Vector, class Matrix, class Sparsity>
-void LaplaceProblem<Vector,Matrix,Sparsity>::run ()
+template <typename VectorType, class Matrix, class Sparsity>
+void LaplaceProblem<VectorType,Matrix,Sparsity>::run ()
 {
   make_grid_and_dofs ();
   assemble_system ();
@@ -344,7 +344,7 @@ void LaplaceProblem<Vector,Matrix,Sparsity>::run ()
 
   for (unsigned int i=0; i<solution.size(); ++i)
     deallog
-    //<< typeid(Vector).name ()
+    //<< typeid(VectorType).name ()
     //<< ' '
     //<< typeid(Matrix).name ()
     //<< '-'
index 3d2643da1649b239bc236ce13ae376fe02685db6..e0ff8b84ae733e9785e3168057c1c8735a693e0b 100644 (file)
 #include <deal.II/lac/solver_cg.h>
 #include <deal.II/lac/precondition.h>
 
-template<class SOLVER, class MATRIX, class VECTOR, class PRECONDITION>
+template<class SOLVER, class MATRIX, typename VectorType, class PRECONDITION>
 void
-check_solve( SOLVER &solver, const MATRIX &A,
-             VECTOR &u, VECTOR &f, const PRECONDITION &P)
+check_solve (SOLVER             &solver,
+             const MATRIX       &A,
+             VectorType         &u,
+             VectorType         &f,
+             const PRECONDITION &P)
 {
   u = 0.;
   f = 1.;
@@ -48,10 +51,13 @@ check_solve( SOLVER &solver, const MATRIX &A,
     }
 }
 
-template<class SOLVER, class MATRIX, class VECTOR, class PRECONDITION>
+template<class SOLVER, class MATRIX, typename VectorType, class PRECONDITION>
 void
-check_Tsolve(SOLVER &solver, const MATRIX &A,
-             VECTOR &u, VECTOR &f, const PRECONDITION &P)
+check_Tsolve (SOLVER             &solver,
+              const MATRIX       &A,
+              VectorType         &u,
+              VectorType         &f,
+              const PRECONDITION &P)
 {
   u = 0.;
   f = 1.;
@@ -136,4 +142,3 @@ int main()
         }
     }
 }
-

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.