]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make VectorTools::project fit for other types of vectors as well.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Tue, 18 Apr 2006 21:52:12 +0000 (21:52 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Tue, 18 Apr 2006 21:52:12 +0000 (21:52 +0000)
git-svn-id: https://svn.dealii.org/trunk@12859 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
deal.II/deal.II/source/numerics/vectors.instance.h
deal.II/doc/news/changes.html
tests/deal.II/Makefile

index 7d17b68c68402cc020ca17387249dae0859c0e10..66a4e4190baf67a5b21ed987bd15209e0ddca005 100644 (file)
@@ -25,7 +25,7 @@ template <int dim> class Point;
 template <int dim> class Function;
 template <int dim> class FunctionMap;
 template <int dim> class Quadrature;
-template <int dim> class QGauss2;
+template <int dim> class QGauss;
 
 template <typename number> class Vector;
 template <typename number> class FullMatrix;
@@ -491,15 +491,15 @@ class VectorTools
                                      * See the general documentation of this
                                      * class for further information.
                                      */
-    template <int dim>
+    template <int dim, class VECTOR>
     static void project (const Mapping<dim>       &mapping,
                         const DoFHandler<dim>    &dof,
                         const ConstraintMatrix   &constraints,
                         const Quadrature<dim>    &quadrature,
                         const Function<dim>      &function,
-                        Vector<double>           &vec,
+                        VECTOR                   &vec,
                         const bool                enforce_zero_boundary = false,
-                        const Quadrature<dim-1>  &q_boundary = QGauss2<dim-1>(),
+                        const Quadrature<dim-1>  &q_boundary = QGauss<dim-1>(2),
                         const bool                project_to_boundary_first = false);
 
                                     /**
@@ -512,12 +512,13 @@ class VectorTools
                                      * quadrature formula is an invalid
                                      * object since it makes no sense in 1d.
                                      */
+    template <class VECTOR>
     static void project (const Mapping<1>         &mapping,
                         const DoFHandler<1>      &dof,
                         const ConstraintMatrix   &constraints,
                         const Quadrature<1>      &quadrature,
                         const Function<1>        &function,
-                        Vector<double>           &vec,
+                        VECTOR                   &vec,
                         const bool                enforce_zero_boundary = false,
                         const Quadrature<0>      &q_boundary = *invalid_face_quadrature,
                         const bool                project_to_boundary_first = false);
@@ -527,14 +528,14 @@ class VectorTools
                                      * function, see above, with
                                      * <tt>mapping=MappingQ1@<dim@>()</tt>.
                                      */
-    template <int dim>
+    template <int dim, class VECTOR>
     static void project (const DoFHandler<dim>    &dof,
                         const ConstraintMatrix   &constraints,
                         const Quadrature<dim>    &quadrature,
                         const Function<dim>      &function,
-                        Vector<double>           &vec,
+                        VECTOR                   &vec,
                         const bool                enforce_zero_boundary = false,
-                        const Quadrature<dim-1>  &q_boundary = QGauss2<dim-1>(),
+                        const Quadrature<dim-1>  &q_boundary = QGauss<dim-1>(2),
                         const bool                project_to_boundary_first = false);
 
                                     /**
index c1e57726a5b9bdf527e834de376d8901e678d417..116bfe9bc0740fc78d467d27bdad54f8d164f129 100644 (file)
@@ -273,22 +273,25 @@ VectorTools::interpolate (const DoFHandler<dim>           &dof_1,
 
 #if deal_II_dimension == 1
 
+template <class VECTOR>
 void VectorTools::project (const Mapping<1>       &,
                           const DoFHandler<1>    &,
                           const ConstraintMatrix &,
                           const Quadrature<1>    &,
                           const Function<1>      &,
-                          Vector<double>         &,
+                          VECTOR                 &,
                           const bool              ,
                           const Quadrature<0>    &,
                           const bool              )
 {
-                                  // this function should easily be implemented
-                                  // using the template below. However some
-                                  // changes have to be made since faces don't
-                                  // exist in 1D. Maybe integrate the creation of
-                                  // zero boundary values into the
-                                  // project_boundary_values function?
+                                  // this function should easily be
+                                  // implemented using the template
+                                  // below. However some changes have
+                                  // to be made since faces don't
+                                  // exist in 1D. Maybe integrate the
+                                  // creation of zero boundary values
+                                  // into the project_boundary_values
+                                  // function?
   Assert (false, ExcNotImplemented());
 }
 
@@ -296,19 +299,22 @@ void VectorTools::project (const Mapping<1>       &,
 #endif
 
 
-template <int dim>
+template <int dim, class VECTOR>
 void VectorTools::project (const Mapping<dim>       &mapping,
                           const DoFHandler<dim>    &dof,
                           const ConstraintMatrix   &constraints,
                           const Quadrature<dim>    &quadrature,
                           const Function<dim>      &function,
-                          Vector<double>           &vec,
+                          VECTOR                   &vec_result,
                           const bool                enforce_zero_boundary,
                           const Quadrature<dim-1>  &q_boundary,
                           const bool                project_to_boundary_first)
 {
   Assert (dof.get_fe().n_components() == function.n_components,
          ExcInvalidFE());
+
+  Assert (vec_result.size() == dof.n_dofs(),
+          ExcDimensionMismatch (vec_result.size(), dof.n_dofs()));
   
   const FiniteElement<dim> &fe = dof.get_fe();
 
@@ -377,7 +383,7 @@ void VectorTools::project (const Mapping<dim>       &mapping,
 
 
                                   // set up mass matrix and right hand side
-  vec.reinit (dof.n_dofs());
+  Vector<double> vec (dof.n_dofs());
   SparsityPattern sparsity(dof.n_dofs(),
                           dof.n_dofs(),
                           dof.max_couplings_between_dofs());
@@ -409,15 +415,24 @@ void VectorTools::project (const Mapping<dim>       &mapping,
   
                                   // distribute solution
   constraints.distribute (vec);
+
+                                   // copy vec into vec_result. we
+                                   // can't use ve_result itself
+                                   // above, since it may be of
+                                   // another type than Vector<double>
+                                   // and that wouldn't necessarily go
+                                   // together with the matrix and
+                                   // other functions
+  std::copy (vec.begin(), vec.end(), vec_result.begin());
 }
 
 
-template <int dim>
+template <int dim, class VECTOR>
 void VectorTools::project (const DoFHandler<dim>    &dof,
                           const ConstraintMatrix   &constraints,
                           const Quadrature<dim>    &quadrature,
                           const Function<dim>      &function,
-                          Vector<double>           &vec,
+                          VECTOR                   &vec,
                           const bool                enforce_zero_boundary,
                           const Quadrature<dim-1>  &q_boundary,
                           const bool                project_to_boundary_first)
index bf045fa4e7772c098ff3cb70419e79f2fc45ca09..a69fb010fbd3f6f294c2a42138def01b26432770 100644 (file)
@@ -2,7 +2,7 @@
 //    $Id$
 //    Version: $Name$
 //
-//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005 by the deal.II authors
+//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006 by the deal.II authors
 //
 //    This file is subject to QPL and may not be  distributed
 //    without copyright and license information. Please refer
@@ -88,16 +88,6 @@ namespace
 #include "vectors.instance.h"
 #undef VEC
 
-template
-void VectorTools::project<deal_II_dimension>
-(const DoFHandler<deal_II_dimension>   &,
- const ConstraintMatrix                &,
- const Quadrature<deal_II_dimension>   &,
- const Function<deal_II_dimension>     &,
- Vector<double>                        &,
- const bool,
- const Quadrature<deal_II_dimension-1> &,
- const bool);
 
 template
 void VectorTools::create_right_hand_side<deal_II_dimension>
@@ -193,17 +183,6 @@ void VectorTools::interpolate_boundary_values<deal_II_dimension>
  std::map<unsigned int,double>       &,
  const std::vector<bool>    &);
 
-template
-void VectorTools::project<deal_II_dimension>
-(const Mapping<deal_II_dimension>      &,
- const DoFHandler<deal_II_dimension>   &,
- const ConstraintMatrix                &,
- const Quadrature<deal_II_dimension>   &,
- const Function<deal_II_dimension>     &,
- Vector<double>                        &,
- const bool,
- const Quadrature<deal_II_dimension-1> &,
- const bool);
 #endif
 
 
index 4faf6c45a153ea321b31859d1c1cb64448f0676c..7178be51a14d86fa122254ca8eaf46139d35557e 100644 (file)
@@ -114,3 +114,26 @@ double VectorTools::compute_mean_value<deal_II_dimension>
  const Quadrature<deal_II_dimension>&,
  const VEC&,
  const unsigned int);
+
+template
+void VectorTools::project
+(const Mapping<deal_II_dimension>      &,
+ const DoFHandler<deal_II_dimension>   &,
+ const ConstraintMatrix                &,
+ const Quadrature<deal_II_dimension>   &,
+ const Function<deal_II_dimension>     &,
+ VEC                                   &,
+ const bool,
+ const Quadrature<deal_II_dimension-1> &,
+ const bool);
+
+template
+void VectorTools::project
+(const DoFHandler<deal_II_dimension>   &,
+ const ConstraintMatrix                &,
+ const Quadrature<deal_II_dimension>   &,
+ const Function<deal_II_dimension>     &,
+ VEC                                   &,
+ const bool,
+ const Quadrature<deal_II_dimension-1> &,
+ const bool);
index d6750f954dd3b4b7eac6805c4323dd3a13e7b6aa..021e7db38da7819806dfa4bf90b92e7ea6a1e8cb 100644 (file)
@@ -583,6 +583,14 @@ inconvenience this causes.
 <h3>deal.II</h3>
 
 <ol>
+  <li> <p>
+       New: The function <code>VectorTools::project</code> functions can now
+       also be used for vector arguments of type other than
+       <code>Vector&lt;double&gt;</code>.
+       <br> 
+       (WB 2006/04/17)
+       </p>
+
   <li> <p>
        New: The function <code>GridTools::get_finest_common_cells</code> can be
        used to find those cells that are shared by two triangulations that are
index 31c34477064c93141eea8025eae985c2605f671d..9d837e4ff2aa9fee52b447de8446330359daaf5b 100644 (file)
@@ -53,7 +53,8 @@ tests_x = block_matrices \
        subcelldata \
        interpolate_boundary_values_* \
        have_same_coarse_mesh_* \
-       get_finest_common_cells_* 
+       get_finest_common_cells_* \
+       project_*
 
 # from above list of regular expressions, generate the real set of
 # tests

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.