]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
simplify declarations by removing specializations
authorkanschat <kanschat@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 25 Jul 2007 07:38:02 +0000 (07:38 +0000)
committerkanschat <kanschat@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 25 Jul 2007 07:38:02 +0000 (07:38 +0000)
git-svn-id: https://svn.dealii.org/trunk@14869 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/numerics/vectors.h
deal.II/deal.II/include/numerics/vectors.templates.h
deal.II/deal.II/source/numerics/vectors.cc

index 9f07bf08cb45740307f4f6837749016dc7ee20fe..7f927e2497cc79677c44b9316747c94f7a2b6f38 100644 (file)
@@ -620,32 +620,9 @@ class VectorTools
                                  const std::vector<bool>       &component_mask = std::vector<bool>());
 
                                     /**
-                                     * Declaration of specialization
-                                     * of the previous function for
-                                     * 1d and a DoFHandler.
-                                     */
-    static
-    void
-    interpolate_boundary_values (const Mapping<1>              &mapping,
-                                 const DoFHandler<1>           &dof,
-                                 const FunctionMap<1>::type    &function_map,
-                                 std::map<unsigned int,double> &boundary_values,
-                                 const std::vector<bool>       &component_mask = std::vector<bool>());
-
-                                    /**
-                                     * Declaration of specialization
-                                     * of the previous function for
-                                     * 1d and a hp::DoFHandler().
-                                     */
-    static
-    void
-    interpolate_boundary_values (const Mapping<1>              &mapping,
-                                 const hp::DoFHandler<1>         &dof,
-                                 const FunctionMap<1>::type    &function_map,
-                                 std::map<unsigned int,double> &boundary_values,
-                                 const std::vector<bool>       &component_mask = std::vector<bool>());
-
-                                    /**
+                                     * @deprecated This function is there mainly
+                                     * for backward compatibility.
+                                     *
                                      * Same function as above, but
                                      * taking only one pair of
                                      * boundary indicator and
@@ -654,8 +631,6 @@ class VectorTools
                                      * function with remapped
                                      * arguments.
                                      *
-                                     * This function is there mainly
-                                     * for backward compatibility.
                                      */
     template <int dim, template <int> class DH>
     static
@@ -667,22 +642,6 @@ class VectorTools
                                  std::map<unsigned int,double> &boundary_values,
                                  const std::vector<bool>       &component_mask = std::vector<bool>());
 
-                                    /**
-                                     * Declaration of specialization
-                                     * of the previous function for
-                                     * 1d.
-                                     */
-    template <template <int> class DH>
-    static
-    void
-    interpolate_boundary_values (const Mapping<1>              &mapping,
-                                 const DH<1>                   &dof,
-                                 const unsigned char            boundary_component,
-                                 const Function<1>             &boundary_function,
-                                 std::map<unsigned int,double> &boundary_values,
-                                 const std::vector<bool>       &component_mask = std::vector<bool>());
-
-    
                                     /**
                                      * Calls the other
                                      * @p interpolate_boundary_values
@@ -730,6 +689,11 @@ class VectorTools
                                      * matches that of the finite
                                      * element used by @p dof.
                                      *
+                                     * In 1d, projection equals
+                                     * interpolation. Therefore,
+                                     * interpolate_boundary_values is
+                                     * called.
+                                     *
                                      * See the general documentation of this
                                      * class for further information.
                                      */
@@ -740,20 +704,6 @@ class VectorTools
                                         const Quadrature<dim-1>  &q,
                                         std::map<unsigned int,double> &boundary_values);
 
-                                    /**
-                                     * Declaration of specialization
-                                     * of the previous function for
-                                     * 1d. Since in 1d projection
-                                     * equals interpolation, the
-                                     * interpolation function is
-                                     * called.
-                                     */
-    static void project_boundary_values (const Mapping<1>           &mapping,
-                                        const DoFHandler<1>        &dof,
-                                        const FunctionMap<1>::type &boundary_functions,
-                                        const Quadrature<0>        &q,
-                                        std::map<unsigned int,double> &boundary_values);
-
                                     /**
                                      * Calls the @p project_boundary_values
                                      * function, see above, with
@@ -888,20 +838,6 @@ class VectorTools
                                                 Vector<double>          &rhs_vector,
                                                 const std::set<unsigned char> &boundary_indicators = std::set<unsigned char>());
 
-                                    /**
-                                     * Specialization of above
-                                     * function for 1d. Since the
-                                     * computation is not useful in
-                                     * 1d, this function simply
-                                     * throws an exception.
-                                     */
-    static void create_boundary_right_hand_side (const Mapping<1>    &mapping,
-                                                const DoFHandler<1> &dof,
-                                                const Quadrature<0> &q,
-                                                const Function<1>   &rhs,
-                                                Vector<double>      &rhs_vector,
-                                                const std::set<unsigned char> &boundary_indicators = std::set<unsigned char>());
-
                                     /**
                                      * Calls the
                                      * @p create_boundary_right_hand_side
@@ -927,20 +863,6 @@ class VectorTools
                                                 Vector<double>          &rhs_vector,
                                                 const std::set<unsigned char> &boundary_indicators = std::set<unsigned char>());
 
-                                    /**
-                                     * Specialization of above
-                                     * function for 1d. Since the
-                                     * computation is not useful in
-                                     * 1d, this function simply
-                                     * throws an exception.
-                                     */
-    static void create_boundary_right_hand_side (const hp::MappingCollection<1>    &mapping,
-                                                const hp::DoFHandler<1> &dof,
-                                                const hp::QCollection<0> &q,
-                                                const Function<1>   &rhs,
-                                                Vector<double>      &rhs_vector,
-                                                const std::set<unsigned char> &boundary_indicators = std::set<unsigned char>());
-
                                     /**
                                      * Calls the @p
                                      * create_boundary_right_hand_side
index 13e5071a7edc52524688910f75a4dcdb8ca5f045..bc5f5d8be8ccd29e13873a744197c8fc6da8744f 100644 (file)
@@ -980,15 +980,17 @@ VectorTools::create_boundary_right_hand_side (const Mapping<dim>      &mapping,
 
 #else
 
+// Implementation for 1D
+template <int dim>
 void
-VectorTools::create_boundary_right_hand_side (const Mapping<1>    &,
-                                             const DoFHandler<1> &,
-                                             const Quadrature<0> &,
-                                             const Function<1>   &,
+VectorTools::create_boundary_right_hand_side (const Mapping<dim>    &,
+                                             const DoFHandler<dim> &,
+                                             const Quadrature<dim-1> &,
+                                             const Function<dim>   &,
                                              Vector<double>      &,
                                              const std::set<unsigned char> &)
 {
-  Assert (false, ExcImpossibleInDim(1));
+  Assert (false, ExcImpossibleInDim(dim));
 }
 
 #endif
@@ -1145,11 +1147,13 @@ VectorTools::create_boundary_right_hand_side (const hp::MappingCollection<dim>
 
 #else
 
+// Implementation for 1D
+template <int dim>
 void
-VectorTools::create_boundary_right_hand_side (const hp::MappingCollection<1>    &,
-                                             const hp::DoFHandler<1> &,
-                                             const hp::QCollection<0> &,
-                                             const Function<1>   &,
+VectorTools::create_boundary_right_hand_side (const hp::MappingCollection<dim>    &,
+                                             const hp::DoFHandler<dim> &,
+                                             const hp::QCollection<dim-1> &,
+                                             const Function<dim>   &,
                                              Vector<double>      &,
                                              const std::set<unsigned char> &)
 {
@@ -1180,12 +1184,12 @@ VectorTools::create_boundary_right_hand_side (const hp::DoFHandler<dim>   &dof_h
 
 //TODO[?] Actually the Mapping object should be a MappingCollection object for the
 // hp::DoFHandler.
-template <template <int> class DH>
+template <int dim, template <int> class DH>
 void
-VectorTools::interpolate_boundary_values (const Mapping<1>         &,
-                                         const DH<1>              &dof,
+VectorTools::interpolate_boundary_values (const Mapping<dim>         &,
+                                         const DH<dim>              &dof,
                                          const unsigned char       boundary_component,
-                                         const Function<1>        &boundary_function,
+                                         const Function<dim>        &boundary_function,
                                          std::map<unsigned int,double> &boundary_values,
                                          const std::vector<bool>       &component_mask_)
 {
@@ -1207,7 +1211,7 @@ VectorTools::interpolate_boundary_values (const Mapping<1>         &,
                                   // cell by first traversing the coarse
                                   // grid to its end and then going
                                   // to the children
-  typename DH<1>::cell_iterator outermost_cell = dof.begin(0);
+  typename DH<dim>::cell_iterator outermost_cell = dof.begin(0);
   while (outermost_cell->neighbor(direction).state() == IteratorState::valid)
     outermost_cell = outermost_cell->neighbor(direction);
   
@@ -1216,7 +1220,7 @@ VectorTools::interpolate_boundary_values (const Mapping<1>         &,
 
                                    // get the FE corresponding to this
                                    // cell
-  const FiniteElement<1> &fe = outermost_cell->get_fe();
+  const FiniteElement<dim> &fe = outermost_cell->get_fe();
   Assert (fe.n_components() == boundary_function.n_components,
          ExcDimensionMismatch(fe.n_components(), boundary_function.n_components));
 
@@ -1257,14 +1261,17 @@ VectorTools::interpolate_boundary_values (const Mapping<1>         &,
 
 //TODO[?] Actually the Mapping object should be a MappingCollection object for the
 // hp::DoFHandler.
+
+// Implementation for 1D
+template <int dim, template <int> class DH>
 void
-VectorTools::interpolate_boundary_values (const Mapping<1>              &mapping,
-                                         const DoFHandler<1>           &dof,
-                                         const FunctionMap<1>::type    &function_map,
+VectorTools::interpolate_boundary_values (const Mapping<dim>              &mapping,
+                                         const DH<dim>           &dof,
+                                         const typename FunctionMap<dim>::type    &function_map,
                                          std::map<unsigned int,double> &boundary_values,
                                          const std::vector<bool>       &component_mask)
 {
-  for (FunctionMap<1>::type::const_iterator i=function_map.begin();
+  for (typename FunctionMap<dim>::type::const_iterator i=function_map.begin();
        i!=function_map.end(); ++i)
     interpolate_boundary_values (mapping, dof, i->first, *i->second,
                                 boundary_values, component_mask);
@@ -1273,21 +1280,9 @@ VectorTools::interpolate_boundary_values (const Mapping<1>              &mapping
 
 //TODO[?] Actually the Mapping object should be a MappingCollection object for the
 // hp::DoFHandler.
-void
-VectorTools::interpolate_boundary_values (const Mapping<1>              &mapping,
-                                         const hp::DoFHandler<1>           &dof,
-                                         const FunctionMap<1>::type    &function_map,
-                                         std::map<unsigned int,double> &boundary_values,
-                                         const std::vector<bool>       &component_mask)
-{
-  for (FunctionMap<1>::type::const_iterator i=function_map.begin();
-       i!=function_map.end(); ++i)
-    interpolate_boundary_values (mapping, dof, i->first, *i->second,
-                                boundary_values, component_mask);
-}
 
 
-#endif
+#else
 
 
 template <int dim, template <int> class DH>
@@ -1615,7 +1610,7 @@ VectorTools::interpolate_boundary_values (const Mapping<dim>            &mapping
                               component_mask);
 }
 
-
+#endif
   
 //TODO[?] Change for real hp::DoFHandler
 // This function might not work anymore if the real hp::DoFHandler is available.
@@ -1653,11 +1648,13 @@ VectorTools::interpolate_boundary_values (const DH<dim>                 &dof,
 
 #if deal_II_dimension == 1
 
+// Implementation for 1D
+template <int dim>
 void
-VectorTools::project_boundary_values (const Mapping<1>       &mapping,
-                                     const DoFHandler<1>    &dof,
-                                     const FunctionMap<1>::type &boundary_functions,
-                                     const Quadrature<0>  &,
+VectorTools::project_boundary_values (const Mapping<dim>       &mapping,
+                                     const DoFHandler<dim>    &dof,
+                                     const typename FunctionMap<dim>::type &boundary_functions,
+                                     const Quadrature<dim-1>  &,
                                      std::map<unsigned int,double> &boundary_values)
 {
                                   // projection in 1d is equivalent
@@ -1666,7 +1663,7 @@ VectorTools::project_boundary_values (const Mapping<1>       &mapping,
                               boundary_values, std::vector<bool>());
 }
 
-#endif
+#else
 
 
 template <int dim>
@@ -1768,6 +1765,7 @@ VectorTools::project_boundary_values (const Mapping<dim>       &mapping,
       boundary_values[i] = boundary_projection(dof_to_boundary_mapping[i]);
 }
 
+#endif
 
 template <int dim>
 void
index 860a756bd1abfeab6dd91af5c452916af7970e38..6d37f31ba587bc27c047fa925d62f498ba7dcad1 100644 (file)
@@ -89,7 +89,6 @@ void VectorTools::create_point_source_vector<deal_II_dimension>
  const Point<deal_II_dimension>      &,
  Vector<double>                      &);
 
-#if deal_II_dimension != 1
 template
 void
 VectorTools::create_boundary_right_hand_side<deal_II_dimension>
@@ -99,7 +98,6 @@ VectorTools::create_boundary_right_hand_side<deal_II_dimension>
  const Function<deal_II_dimension>   &,
  Vector<double>                      &,
  const std::set<unsigned char> &);
-#endif
 
 template
 void
@@ -110,7 +108,6 @@ VectorTools::create_boundary_right_hand_side<deal_II_dimension>
  Vector<double>                      &,
  const std::set<unsigned char> &);
 
-#if deal_II_dimension != 1
 template
 void
 VectorTools::create_boundary_right_hand_side<deal_II_dimension>
@@ -120,7 +117,6 @@ VectorTools::create_boundary_right_hand_side<deal_II_dimension>
  const Function<deal_II_dimension>   &,
  Vector<double>                      &,
  const std::set<unsigned char> &);
-#endif
 
 template
 void
@@ -155,7 +151,6 @@ void VectorTools::interpolate_boundary_values<deal_II_dimension> (
   std::map<unsigned int,double>       &,
   const std::vector<bool>    &);
 
-#if deal_II_dimension != 1
 template
 void VectorTools::project_boundary_values<deal_II_dimension>
 (const Mapping<deal_II_dimension>     &,
@@ -163,7 +158,6 @@ void VectorTools::project_boundary_values<deal_II_dimension>
  const FunctionMap<deal_II_dimension>::type &,
  const Quadrature<deal_II_dimension-1>&,
  std::map<unsigned int,double>        &);
-#endif
 
 template
 void VectorTools::project_boundary_values<deal_II_dimension>
@@ -174,22 +168,22 @@ void VectorTools::project_boundary_values<deal_II_dimension>
 
 
 
-// Due to introducing the DoFHandler as a template parameter,
-// the following instantiations are required in 1d
-#if deal_II_dimension == 1
-template
-void VectorTools::interpolate_boundary_values<deal_II_dimension> 
-(const Mapping<1>         &,
- const DoFHandler<1>      &,
- const unsigned char,
- const Function<1>        &,
- std::map<unsigned int,double> &,
- const std::vector<bool>       &);
-#endif
+// // Due to introducing the DoFHandler as a template parameter,
+// // the following instantiations are required in 1d
+// #if deal_II_dimension == 1
+// template
+// void VectorTools::interpolate_boundary_values<deal_II_dimension> 
+// (const Mapping<1>         &,
+//  const DoFHandler<1>      &,
+//  const unsigned char,
+//  const Function<1>        &,
+//  std::map<unsigned int,double> &,
+//  const std::vector<bool>       &);
+// #endif
 
 // the following two functions are not derived from a template in 1d
 // and thus need no explicit instantiation
-#if deal_II_dimension > 1
+//#if deal_II_dimension > 1
 template
 void VectorTools::interpolate_boundary_values<deal_II_dimension>
 (const Mapping<deal_II_dimension>    &,
@@ -207,7 +201,7 @@ void VectorTools::interpolate_boundary_values<deal_II_dimension>
  std::map<unsigned int,double>       &,
  const std::vector<bool>    &);
 
-#endif
+//#endif
 
 
 template

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.