]> https://gitweb.dealii.org/ - dealii.git/commitdiff
DerivativeApproximation can also work on BlockVector arguments.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 11 Apr 2003 20:03:57 +0000 (20:03 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 11 Apr 2003 20:03:57 +0000 (20:03 +0000)
git-svn-id: https://svn.dealii.org/trunk@7391 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/numerics/derivative_approximation.h
deal.II/deal.II/source/numerics/derivative_approximation.cc
deal.II/doc/news/2002/c-3-4.html

index 3ece87276ec95b895045ee4feffc5fd7a66cac07..208d95155a6a3c458435c657eeaea64b60bdcc70 100644 (file)
@@ -180,11 +180,11 @@ class DerivativeApproximation
                                      * computed. It defaults to the
                                      * first component.
                                      */
-    template <int dim>
+    template <int dim, class InputVector>
     static void
     approximate_gradient (const Mapping<dim>    &mapping,
                          const DoFHandler<dim> &dof,
-                         const Vector<double>  &solution,
+                         const InputVector     &solution,
                          Vector<float>         &derivative_norm,
                          const unsigned int     component = 0);
 
@@ -193,10 +193,10 @@ class DerivativeApproximation
                                      * function, see above, with
                                      * @p{mapping=MappingQ1<dim>()}.
                                      */
-    template <int dim>
+    template <int dim, class InputVector>
     static void
     approximate_gradient (const DoFHandler<dim> &dof,
-                         const Vector<double>  &solution,
+                         const InputVector     &solution,
                          Vector<float>         &derivative_norm,
                          const unsigned int     component = 0);
     
@@ -222,11 +222,11 @@ class DerivativeApproximation
                                      * computed. It defaults to the
                                      * first component.
                                      */
-    template <int dim>
+    template <int dim, class InputVector>
     static void
     approximate_second_derivative (const Mapping<dim>    &mapping,
                                   const DoFHandler<dim> &dof,
-                                  const Vector<double>  &solution,
+                                  const InputVector     &solution,
                                   Vector<float>         &derivative_norm,
                                   const unsigned int     component = 0);
     
@@ -235,10 +235,10 @@ class DerivativeApproximation
                                      * function, see above, with
                                      * @p{mapping=MappingQ1<dim>()}.
                                      */
-    template <int dim>
+    template <int dim, class InputVector>
     static void
     approximate_second_derivative (const DoFHandler<dim> &dof,
-                                  const Vector<double>  &solution,
+                                  const InputVector     &solution,
                                   Vector<float>         &derivative_norm,
                                   const unsigned int     component = 0);
 
@@ -306,9 +306,10 @@ class DerivativeApproximation
                                          * field at the center of each
                                          * cell).
                                          */
+       template <class InputVector>
        static ProjectedDerivative
        get_projected_derivative (const FEValues<dim>  &fe_values,
-                                 const Vector<double> &solution,
+                                 const InputVector    &solution,
                                  const unsigned int    component);
     
                                         /**
@@ -386,9 +387,10 @@ class DerivativeApproximation
                                          * field at the center of each
                                          * cell).
                                          */
+       template <class InputVector>
        static ProjectedDerivative
        get_projected_derivative (const FEValues<dim>  &fe_values,
-                                 const Vector<double> &solution,
+                                 const InputVector    &solution,
                                  const unsigned int    component);
        
                                         /**
@@ -452,11 +454,11 @@ class DerivativeApproximation
                                      * solution vector we are to work
                                      * on.
                                      */
-    template <class DerivativeDescription, int dim>
+    template <class DerivativeDescription, int dim, class InputVector>
     static void
     approximate_derivative (const Mapping<dim>    &mapping,
                            const DoFHandler<dim> &dof,
-                           const Vector<double>  &solution,
+                           const InputVector     &solution,
                            const unsigned int     component,
                            Vector<float>         &derivative_norm);
 
@@ -466,11 +468,11 @@ class DerivativeApproximation
                                      * the range given by the third
                                      * parameter.
                                      */
-    template <class DerivativeDescription, int dim>
+    template <class DerivativeDescription, int dim, class InputVector>
     static void
     approximate (const Mapping<dim>    &mapping,
                 const DoFHandler<dim> &dof,
-                const Vector<double>  &solution,
+                const InputVector     &solution,
                 const unsigned int     component,
                 const IndexInterval   &index_interval,
                 Vector<float>         &derivative_norm);    
index 7549c9b07975f308a798a402cdfc376001d1a96b..09dcd4deb6143456efa470faae693f5b5240da58 100644 (file)
@@ -16,6 +16,7 @@
 #include <base/thread_management.h>
 #include <base/multithread_info.h>
 #include <lac/vector.h>
+#include <lac/block_vector.h>
 #include <grid/tria_iterator.h>
 #include <dofs/dof_accessor.h>
 #include <dofs/dof_handler.h>
@@ -46,11 +47,12 @@ const UpdateFlags DerivativeApproximation::SecondDerivative<dim>::update_flags =
 
 
 template <int dim>
+template <class InputVector>
 inline
 typename DerivativeApproximation::Gradient<dim>::ProjectedDerivative
 DerivativeApproximation::Gradient<dim>::
 get_projected_derivative (const FEValues<dim>  &fe_values,
-                         const Vector<double> &solution,
+                         const InputVector    &solution,
                          const unsigned int    component)
 {
   if (fe_values.get_fe().n_components() == 1)
@@ -94,11 +96,12 @@ DerivativeApproximation::Gradient<dim>::symmetrize (Derivative &)
 
 
 template <int dim>
+template <class InputVector>
 inline
 typename DerivativeApproximation::SecondDerivative<dim>::ProjectedDerivative
 DerivativeApproximation::SecondDerivative<dim>::
 get_projected_derivative (const FEValues<dim>  &fe_values,
-                         const Vector<double> &solution,
+                         const InputVector    &solution,
                          const unsigned int    component)
 {
   if (fe_values.get_fe().n_components() == 1)
@@ -381,12 +384,12 @@ DerivativeApproximation::SecondDerivative<dim>::symmetrize (Derivative &d)
 
 
 
-template <int dim>
+template <int dim, class InputVector>
 void 
 DerivativeApproximation::
 approximate_gradient (const Mapping<dim>    &mapping,
                      const DoFHandler<dim> &dof_handler,
-                     const Vector<double>  &solution,
+                     const InputVector     &solution,
                      Vector<float>         &derivative_norm,
                      const unsigned int     component)
 {
@@ -398,11 +401,11 @@ approximate_gradient (const Mapping<dim>    &mapping,
 }
 
 
-template <int dim>
+template <int dim, class InputVector>
 void 
 DerivativeApproximation::
 approximate_gradient (const DoFHandler<dim> &dof_handler,
-                     const Vector<double>  &solution,
+                     const InputVector     &solution,
                      Vector<float>         &derivative_norm,
                      const unsigned int     component)
 {
@@ -416,12 +419,12 @@ approximate_gradient (const DoFHandler<dim> &dof_handler,
 }
 
 
-template <int dim>
+template <int dim, class InputVector>
 void 
 DerivativeApproximation::
 approximate_second_derivative (const Mapping<dim>    &mapping,
                               const DoFHandler<dim> &dof_handler,
-                              const Vector<double>  &solution,
+                              const InputVector     &solution,
                               Vector<float>         &derivative_norm,
                               const unsigned int     component)
 {
@@ -433,11 +436,11 @@ approximate_second_derivative (const Mapping<dim>    &mapping,
 }
 
 
-template <int dim>
+template <int dim, class InputVector>
 void 
 DerivativeApproximation::
 approximate_second_derivative (const DoFHandler<dim> &dof_handler,
-                              const Vector<double>  &solution,
+                              const InputVector     &solution,
                               Vector<float>         &derivative_norm,
                               const unsigned int     component)
 {
@@ -451,12 +454,12 @@ approximate_second_derivative (const DoFHandler<dim> &dof_handler,
 }
 
 
-template <class DerivativeDescription, int dim>
+template <class DerivativeDescription, int dim, class InputVector>
 void 
 DerivativeApproximation::
 approximate_derivative (const Mapping<dim>    &mapping,
                        const DoFHandler<dim> &dof_handler,
-                       const Vector<double>  &solution,
+                       const InputVector     &solution,
                        const unsigned int     component,
                        Vector<float>         &derivative_norm)
 {
@@ -471,11 +474,11 @@ approximate_derivative (const Mapping<dim>    &mapping,
   Threads::ThreadGroup<> threads;
   void (*fun_ptr) (const Mapping<dim>    &,
                   const DoFHandler<dim> &,
-                  const Vector<double>  &,
+                  const InputVector     &,
                   const unsigned int     ,
                   const IndexInterval   &,
                   Vector<float>         &)
-    = &DerivativeApproximation::template approximate<DerivativeDescription,dim>;
+    = &DerivativeApproximation::template approximate<DerivativeDescription,dim,InputVector>;
   for (unsigned int i=0; i<n_threads; ++i)
     threads += Threads::spawn (fun_ptr)(mapping, dof_handler, solution, component,
                                         index_intervals[i],
@@ -485,11 +488,11 @@ approximate_derivative (const Mapping<dim>    &mapping,
 
 
 
-template <class DerivativeDescription, int dim>
+template <class DerivativeDescription, int dim, class InputVector>
 void 
 DerivativeApproximation::approximate (const Mapping<dim>    &mapping,
                                      const DoFHandler<dim> &dof_handler,
-                                     const Vector<double>  &solution,
+                                     const InputVector     &solution,
                                      const unsigned int     component,
                                      const IndexInterval   &index_interval,
                                      Vector<float>         &derivative_norm)
@@ -719,11 +722,30 @@ template
 void 
 DerivativeApproximation::
 approximate_gradient<deal_II_dimension>
+(const Mapping<deal_II_dimension> &mapping,
+ const DoFHandler<deal_II_dimension> &dof_handler,
+ const BlockVector<double>  &solution,
+ Vector<float>         &derivative_norm,
+ const unsigned int     component);
+
+template
+void 
+DerivativeApproximation::
+approximate_gradient<deal_II_dimension>
 (const DoFHandler<deal_II_dimension> &dof_handler,
  const Vector<double>  &solution,
  Vector<float>         &derivative_norm,
  const unsigned int     component);
 
+template
+void 
+DerivativeApproximation::
+approximate_gradient<deal_II_dimension>
+(const DoFHandler<deal_II_dimension> &dof_handler,
+ const BlockVector<double>  &solution,
+ Vector<float>         &derivative_norm,
+ const unsigned int     component);
+
 template
 void 
 DerivativeApproximation::
@@ -733,6 +755,17 @@ approximate_second_derivative<deal_II_dimension>
  const Vector<double>  &solution,
  Vector<float>         &derivative_norm,
  const unsigned int     component);
+
+template
+void 
+DerivativeApproximation::
+approximate_second_derivative<deal_II_dimension>
+(const Mapping<deal_II_dimension>    &mapping,
+ const DoFHandler<deal_II_dimension> &dof_handler,
+ const BlockVector<double>  &solution,
+ Vector<float>         &derivative_norm,
+ const unsigned int     component);
+
 template
 void 
 DerivativeApproximation::
@@ -742,6 +775,16 @@ approximate_second_derivative<deal_II_dimension>
  Vector<float>         &derivative_norm,
  const unsigned int     component);
 
+template
+void 
+DerivativeApproximation::
+approximate_second_derivative<deal_II_dimension>
+(const DoFHandler<deal_II_dimension> &dof_handler,
+ const BlockVector<double>  &solution,
+ Vector<float>         &derivative_norm,
+ const unsigned int     component);
+
+
 // static variables
 // 
 // on AIX, the linker is unhappy about some missing symbols. they
index 1cda9498faac6849b5898c9c73b0ff25931d4ad1..3c53bdc8b0d9d85014f6e5cb88e16433d6aba76d 100644 (file)
@@ -747,6 +747,14 @@ contributor's names are abbreviated by WB (Wolfgang Bangerth), GK
 <h3>deal.II</h3>
 
 <ol>
+  <li> <p>
+       New: The functions in the <code
+       class="class">DerivativeApproximation</code> class can now also
+       work on <code class="class">BlockVector</code>s/
+       <br>
+       (WB 2003/04/11)
+       </p>
+
   <li> <p>
        Changed: The cell argument to <code
        class="member">Mapping::transform_unit_to_real_cell</code> (and its

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.