]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Move functions from paralution_vector.h to paralution_vector.cc (not inlined anymore...
authorturcksin <turcksin@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 3 Jun 2014 19:58:58 +0000 (19:58 +0000)
committerturcksin <turcksin@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 3 Jun 2014 19:58:58 +0000 (19:58 +0000)
git-svn-id: https://svn.dealii.org/branches/branch_paralution@33011 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/lac/paralution_vector.h
deal.II/source/lac/CMakeLists.txt
deal.II/source/lac/paralution_vector.cc [new file with mode: 0644]
tests/lac/paralution_vector_01.cc
tests/lac/paralution_vector_01.with_paralution=yes.output

index c8136a3866b1aa6953105f618b2de8029d01a139..16f553bf96c9d3411352dfc45413f6cdfc4e713f 100644 (file)
@@ -1,7 +1,7 @@
 // ---------------------------------------------------------------------
 // $Id: paralution_vector.h 30040 2013-07-18 17:06:48Z maier $
 //
-// Copyright (C) 2013 by the deal.II authors
+// Copyright (C) 2013, 2014 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -52,7 +52,7 @@ namespace ParalutionWrappers
    *
    * @ingroup ParalutionWrappers
    * @ingroup Vectors
-   * @author Bruno Turcksin, 2013
+   * @author Bruno Turcksin, 2013, 2014
    */
   //TODO: lots of functions are missing
   template <typename Number>
@@ -92,14 +92,12 @@ namespace ParalutionWrappers
     Vector();
 
     /**
-     * Copy constructor. Sets the dimension to that of the given vector, and
-     * copies all the elements. The copied vector stays on the host/device
-     * where it is create.
+     * Copy constructor. Set the dimension to that of the given vector, and
+     * copies all the elements. If copy_backend is set to false, the copied
+     * vector stays on the host/device where it is created. Otherwise the
+     * copied vector is moved to the host/device of the given vector.
      */
-    //TODO: Look to use choose between CopyFrom and Clone. Difference between
-    //copy and clone: copy the vector stays on his host/device, with clone
-    //the vector goes on the same host/device.
-    Vector(const Vector<Number> &v);
+    Vector(const Vector<Number> &v, bool copy_backend = false);
 
     /**
      * Constructor. Set dimensionto @p n and initialize all the elements
@@ -113,16 +111,6 @@ namespace ParalutionWrappers
      */
     explicit Vector(const size_type n);
 
-    /**
-     * Initialize the vector with a given range of values pointed to by the
-     * iterators. This function is there in anlogy to the @p std::vector
-     * class.
-     */
-    //TODO
-    // template <typename InputIterator>
-    // Vector (const InputIterator first,
-    //         const InputIterator last);
-
     /**
      * Destructor.
      */
@@ -133,7 +121,7 @@ namespace ParalutionWrappers
      * @p PETScWrappers::Vector class.
      *
      * For the PETSC vector wrapper class, this function compresses the
-     * underlying representation of the PETSc object, i.e. flusehes the
+     * underlying representation of the PETSc object, i.e. flushes the
      * buffers of the vector obkect if it has any. This function is necessary
      * after writing into a vector element-by-element and before anything else
      * can be done on it.
@@ -144,7 +132,6 @@ namespace ParalutionWrappers
     void compress (::dealii::VectorOperation::values operation
                    =::dealii::VectorOperation::unknown) const;
 
-
     /**
      * Change the dimension of the vector to @p N. The vector is filled with
      * zeros.
@@ -158,6 +145,11 @@ namespace ParalutionWrappers
      */
     void reinit(const Vector<Number> &v);
 
+    /**
+     * Free the vector.
+     */
+    void clear();
+
     /**
      * Return dimension of the vector.
      */
@@ -243,16 +235,6 @@ namespace ParalutionWrappers
      * index set.
      */
     IndexSet locally_owned_elements() const;
-
-    /**
-     * Move the Vector to the accelerator.
-     */
-    void move_to_accelerator();
-
-    /**
-     * Move the Vector to the host.
-     */
-    void move_to_host();
     //@}
 
     /**
@@ -280,7 +262,7 @@ namespace ParalutionWrappers
     /**
      * Access the @p ith component as a writeable reference.
      *
-     * Exactly thte asame as operator().
+     * Exactly thsame as operator().
      */
     Number &operator[] (const size_type i);
 
@@ -309,7 +291,7 @@ namespace ParalutionWrappers
     const paralution::LocalVector<Number>& paralution_vector() const;
 
     /**
-     * Return a (modifyable) reference to the underlying Paralution
+     * Return a (modifiable) reference to the underlying Paralution
      * LocalVector class.
      */
     paralution::LocalVector<Number>& paralution_vector();
@@ -351,9 +333,9 @@ namespace ParalutionWrappers
      * other two <tt>add()</tt> functions above.
      */
     template <typename Number2>
-    void add (const size_type n_elements,
+    void add (const size_type  n_elements,
               const size_type *indices,
-              const Number2 *values);
+              const Number2   *values);
 
     /**
      * Addition of @p s to all components. Note that @p s is a scalar and
@@ -429,17 +411,52 @@ namespace ParalutionWrappers
     void equ (const Number a, const Vector<Number> &u,
               const Number b, const Vector<Number> &v,
               const Number c, const Vector<Number> &w);
+    /**
+     * Return the scalar product <tt>*this^T *x</tt>.
+     */
+    Number scalar_product(const Vector<Number> &x) const;
+
     //@}
 
     /**
      * @name 4: Mixed stuff
      */
     //@{
+
+    /**
+     * Move the Vector to the accelerator.
+     */
+    void move_to_accelerator();
+
+    /**
+     * Move the Vector to the host.
+     */
+    void move_to_host();
+
+    /**
+     * Move the Vector to the accelerator. The function returns immediately and
+     * performs the asynchronous transfer in the background.
+     */
+    void move_to_accelerator_async();
+
+    /**
+     * Move the Vector to the host. The function returns immediately and
+     * performs the asynchronous transfer in the background.
+     */
+    void move_to_host_async();
+
+    /**
+     * Synchronize the code when move_to_host_async or move_to_accelerator_async
+     * is used.
+     */
+    void sync();
+
     /**
      * Determine an estimate for the memory consumption (in bytes) of this
      * object.
      */
     std::size_t memory_consumption () const;
+
     //@}
 
   private :
@@ -459,14 +476,6 @@ namespace ParalutionWrappers
 
 
 
-  template <typename Number>
-  inline Vector<Number>::Vector(const Vector <Number> &v)
-  {
-    local_vector.CopyFrom(v.paralution_vector());
-  }
-
-
-
   template <typename Number>
   inline Vector<Number>::Vector(const size_type n)
   {
@@ -481,10 +490,14 @@ namespace ParalutionWrappers
     local_vector.Clear();
   }
 
+
+
   template <typename Number>
-  inline void Vector<Number>:: compress(::dealii::VectorOperation::values operation) const
+  inline void Vector<Number>::compress(::dealii::VectorOperation::values operation) const
   {}
 
+
+
   template <typename Number>
   void Vector<Number>::reinit(const size_type n)
   {
@@ -495,7 +508,7 @@ namespace ParalutionWrappers
 
 
   template <typename Number>
-  void Vector<Number>::reinit(const Vector<Number> &v)
+  inline void Vector<Number>::reinit(const Vector<Number> &v)
   {
     local_vector.Clear();
     local_vector.Allocate("deal_ii_vector",v.size());
@@ -503,6 +516,14 @@ namespace ParalutionWrappers
 
 
 
+  template <typename Number>
+  inline void Vector<Number>::clear()
+  {
+    local_vector.Clear();
+  }
+
+
+
   template <typename Number>
   inline std::size_t Vector<Number>::size() const
   {
@@ -565,20 +586,6 @@ namespace ParalutionWrappers
 
 
 
-  template <typename Number>
-  inline Number Vector<Number>::mean_value() const
-  {
-    Number mean(0.);
-    unsigned int i_max(size());
-    for (unsigned int i=0; i<i_max; ++i)
-      mean += local_vector[i];
-    mean /= i_max;
-
-    return mean;
-  }
-
-
-
   template <typename Number>
   inline Number Vector<Number>::l2_norm() const
   {
@@ -593,6 +600,8 @@ namespace ParalutionWrappers
     return false;
   }
 
+
+
   template <typename Number>
   inline bool Vector<Number>::in_local_range(const size_type global_index) const
   {
@@ -609,22 +618,6 @@ namespace ParalutionWrappers
 
 
 
-  template <typename Number>
-  inline void Vector<Number>::move_to_accelerator()
-  {
-    local_vector.MoveToAccelerator();
-  }
-
-
-
-  template <typename Number>
-  inline void Vector<Number>::move_to_host()
-  {
-    local_vector.MoveToHost();
-  }
-
-
-
   template <typename Number>
   inline Number Vector<Number>::operator() (const size_type i) const
   {
@@ -665,21 +658,13 @@ namespace ParalutionWrappers
 
 
 
-  template <typename Number>
-  inline void Vector<Number>::extract_subvector_to (const std::vector<size_type> &indices,
-                                                    std::vector<Number>          &values) const
-  {
-    for (size_type i=0; i<indices.size(); ++i)
-      values[i] = operator()(indices[i]);
-  }
-
-
-
+  // This function is not inlined but it needs to be in the header file because
+  // of ForwardIterator and OutputIterator.
   template <typename Number>
   template <typename ForwardIterator, typename OutputIterator>
-  inline void Vector<Number>::extract_subvector_to (ForwardIterator       indices_begin,
-                                                    const ForwardIterator indices_end,
-                                                    OutputIterator        values_begin) const
+  void Vector<Number>::extract_subvector_to (ForwardIterator       indices_begin,
+                                             const ForwardIterator indices_end,
+                                             OutputIterator        values_begin) const
   {
     while (indices_begin != indices_end)
       {
@@ -732,214 +717,76 @@ namespace ParalutionWrappers
 
 
   template <typename Number>
-  template <typename Number2>
-  inline void Vector<Number>::add (const std::vector<size_type> &indices,
-                                   const std::vector<Number2>   &values)
-  {
-    Assert(indices.size()==values.size(),
-           ExcDimensionMismatch(indices.size(),values.size()));
-    add(indices.size(),&indices[0],&values[0]);
-  }
-
-
-
-  template <typename Number>
-  template <typename Number2>
-  inline void Vector<Number>::add (const std::vector<size_type>    &indices,
-                                   const ::dealii::Vector<Number2> &values)
-  {
-    Assert(indices.size()==values.size(),
-           ExcDimensionMismatch(indices.size(),values.size()));
-    add(indices.size(),&indices[0],values.begin());
-  }
-
-
-
-  template <typename Number>
-  template <typename Number2>
-  inline void Vector<Number>::add (const size_type n_indices,const size_type *indices,
-                                   const Number2 *values)
-  {
-    for (size_type i=0; i<n_indices; ++i)
-      {
-        Assert(indices[i]<size(),ExcIndexRange(indices[i],0,size()));
-        Assert(numbers::is_finite(values[i]),ExcNumberNotFinite());
-
-        local_vector[indices[i]] += values[i];
-      }
-  }
-
-
-
-  template <typename Number>
-  inline void Vector<Number>::add(const Number s)
-  {
-    Assert(numbers::is_finite(s),ExcNumberNotFinite());
-
-    size_type size = local_vector.get_size();
-    for (size_type i=0; i<size; ++i)
-      local_vector[i] += s;
-  }
-
-
-  template <typename Number>
-  inline void Vector<Number>::add(const Vector<Number> &V)
-  {
-    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
-    local_vector.AddScale(V.paralution_vector(),1.);
-  }
-
-
-
-  template <typename Number>
-  inline void Vector<Number>::add(const Number a,const Vector<Number> &V)
-  {
-    Assert(numbers::is_finite(a),ExcNumberNotFinite());
-    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
-
-    local_vector.ScaleAddScale(1.,V.paralution_vector(),a);
-  }
-
-
-
-  template <typename Number>
-  inline void Vector<Number>::add(const Number a,const Vector<Number> &V,
-                                  const Number b,const Vector<Number> &W)
+  inline Vector<Number>& Vector<Number>::operator*= (const Number factor)
   {
-    Assert(numbers::is_finite(a),ExcNumberNotFinite());
-    Assert(numbers::is_finite(b),ExcNumberNotFinite());
-    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
-    Assert(size()==W.size(),ExcDimensionMismatch(size(),W.size()));
-
-    local_vector.ScaleAdd2(1.,V.paralution_vector(),a,W.paralution_vector(),b);
-  }
-
-
+    Assert(numbers::is_finite(factor),ExcNumberNotFinite());
 
-  template <typename Number>
-  inline void Vector<Number>::sadd(const Number s, const Vector<Number> &V)
-  {
-    Assert(numbers::is_finite(s),ExcNumberNotFinite());
-    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
+    local_vector.Scale(factor);
 
-    local_vector.ScaleAdd(s,V.paralution_vector());
+    return *this;
   }
 
 
 
   template <typename Number>
-  inline void Vector<Number>::sadd(const Number s, const Number a, const Vector<Number> &V)
+  inline Vector<Number>& Vector<Number>::operator/= (const Number factor)
   {
-    Assert(numbers::is_finite(s),ExcNumberNotFinite());
-    Assert(numbers::is_finite(a),ExcNumberNotFinite());
-    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
-
-    local_vector.ScaleAddScale(s,V.paralution_vector(),a);
-  }
+    Assert(numbers::is_finite(factor),ExcNumberNotFinite());
 
+    const Number inv_factor(1./factor);
 
+    Assert(numbers::is_finite(inv_factor),ExcNumberNotFinite());
 
-  template <typename Number>
-  inline void Vector<Number>::sadd(const Number s, const Number a, const Vector<Number> &V,
-                                   const Number b, const Vector<Number> &W)
-  {
-    Assert(numbers::is_finite(s),ExcNumberNotFinite());
-    Assert(numbers::is_finite(a),ExcNumberNotFinite());
-    Assert(numbers::is_finite(b),ExcNumberNotFinite());
-    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
-    Assert(size()==W.size(),ExcDimensionMismatch(size(),W.size()));
+    local_vector.Scale(inv_factor);
 
-    local_vector.ScaleAdd2(s,V.paralution_vector(),a,W.paralution_vector(),b);
+    return *this;
   }
 
 
 
   template <typename Number>
-  inline void Vector<Number>::sadd(const Number s, const Number a, const Vector<Number> &V,
-                                   const Number b, const Vector<Number> &W, const Number c,
-                                   const Vector<Number> &X)
+  inline Number Vector<Number>::scalar_product(const Vector<Number> &x) const
   {
-    Assert(numbers::is_finite(s),ExcNumberNotFinite());
-    Assert(numbers::is_finite(a),ExcNumberNotFinite());
-    Assert(numbers::is_finite(b),ExcNumberNotFinite());
-    Assert(numbers::is_finite(c),ExcNumberNotFinite());
-    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
-    Assert(size()==W.size(),ExcDimensionMismatch(size(),W.size()));
-    Assert(size()==X.size(),ExcDimensionMismatch(size(),X.size()));
-
-    local_vector.ScaleAdd2(s,V.paralution_vector(),a,W.paralution_vector(),b);
-    local_vector.AddScale(X.paralution_vector(),c);
+    return local_vector.Dot(x.paralution_vector());
   }
 
 
 
   template <typename Number>
-  inline Vector<Number>& Vector<Number>::operator*= (const Number factor)
+  inline void Vector<Number>::move_to_accelerator()
   {
-    Assert(numbers::is_finite(factor),ExcNumberNotFinite());
-
-    local_vector.Scale(factor);
-
-    return *this;
+    local_vector.MoveToAccelerator();
   }
 
 
 
   template <typename Number>
-  inline Vector<Number>& Vector<Number>::operator/= (const Number factor)
+  inline void Vector<Number>::move_to_host()
   {
-    Assert(numbers::is_finite(factor),ExcNumberNotFinite());
-
-    const Number inv_factor(1./factor);
-
-    Assert(numbers::is_finite(inv_factor),ExcNumberNotFinite());
-
-    local_vector.Scale(inv_factor);
-
-    return *this;
+    local_vector.MoveToHost();
   }
 
 
 
   template <typename Number>
-  inline void Vector<Number>::equ (const Number a, const Vector<Number> &u)
+  inline void Vector<Number>::move_to_accelerator_async()
   {
-    Assert(numbers::is_finite(a), ExcNumberNotFinite());
-    Assert(size()==u.size(),ExcDimensionMismatch(size(),u.size()));
-
-    local_vector.ScaleAddScale(0.,u.paralution_vector(),a);
+    local_vector.MoveToAcceleratorAsync();
   }
 
 
-
   template <typename Number>
-  inline void Vector<Number>::equ (const Number a, const Vector<Number> &u,
-                                   const Number b, const Vector<Number> &v)
+  inline void Vector<Number>::move_to_host_async()
   {
-    Assert(numbers::is_finite(a), ExcNumberNotFinite());
-    Assert(numbers::is_finite(b), ExcNumberNotFinite());
-    Assert(size()==u.size(),ExcDimensionMismatch(size(),u.size()));
-    Assert(size()==v.size(),ExcDimensionMismatch(size(),v.size()));
-
-    local_vector.ScaleAdd2(0.,u.paralution_vector(),a,v.paralution_vector(),b);
+    local_vector.MoveToHostAsync();
   }
 
 
 
   template <typename Number>
-  inline void Vector<Number>::equ (const Number a, const Vector<Number> &u,
-                                   const Number b, const Vector<Number> &v,
-                                   const Number c, const Vector<Number> &w)
+  inline void Vector<Number>::sync()
   {
-    Assert(numbers::is_finite(a), ExcNumberNotFinite());
-    Assert(numbers::is_finite(b), ExcNumberNotFinite());
-    Assert(numbers::is_finite(c), ExcNumberNotFinite());
-    Assert(size()==u.size(),ExcDimensionMismatch(size(),u.size()));
-    Assert(size()==v.size(),ExcDimensionMismatch(size(),v.size()));
-    Assert(size()==w.size(),ExcDimensionMismatch(size(),w.size()));
-
-    local_vector.ScaleAdd2(0.,u.paralution_vector(),a,v.paralution_vector(),b);
-    local_vector.AddScale(w.paralution_vector(),c);
+    local_vector.Sync();
   }
 
 
index 259004d6c2c8d6caaaf10adf8183c5ff8239fbaa..1c811b69505378068c75eab52137b0ac786e708b 100644 (file)
@@ -36,6 +36,7 @@ SET(_src
   paralution_precondition.cc
   paralution_sparse_matrix.cc
   paralution_solver.cc
+  paralution_vector.cc
   petsc_block_sparse_matrix.cc
   petsc_full_matrix.cc
   petsc_matrix_base.cc
diff --git a/deal.II/source/lac/paralution_vector.cc b/deal.II/source/lac/paralution_vector.cc
new file mode 100644 (file)
index 0000000..87b7217
--- /dev/null
@@ -0,0 +1,282 @@
+// ---------------------------------------------------------------------
+// $Id: paralution_vector.cc 31349 2013-10-20 19:07:06Z maier $
+//
+// Copyright (C) 2014 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+#include <deal.II/lac/paralution_vector.h>
+
+#ifdef DEAL_II_WITH_PARALUTION
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace ParalutionWrappers
+{
+
+  template <typename Number>
+  Vector<Number>::Vector(const Vector <Number> &v, bool copy_backend)
+  {
+    if (copy_backend)
+      local_vector.CloneFrom(v.paralution_vector());
+    else
+      local_vector.CopyFrom(v.paralution_vector());
+  }
+
+
+  template <typename Number>
+  Number Vector<Number>::mean_value() const
+  {
+    Number mean(0.);
+    unsigned int i_max(size());
+    for (unsigned int i=0; i<i_max; ++i)
+      mean += local_vector[i];
+    mean /= i_max;
+
+    return mean;
+  }
+
+
+
+  template <typename Number>
+  void Vector<Number>::extract_subvector_to (const std::vector<size_type> &indices,
+                                             std::vector<Number>          &values) const
+  {
+    for (size_type i=0; i<indices.size(); ++i)
+      values[i] = operator()(indices[i]);
+  }
+
+
+
+  template <typename Number>
+  template <typename Number2>
+  void Vector<Number>::add (const std::vector<size_type> &indices,
+                            const std::vector<Number2>   &values)
+  {
+    Assert(indices.size()==values.size(),
+           ExcDimensionMismatch(indices.size(),values.size()));
+    add(indices.size(),&indices[0],&values[0]);
+  }
+
+
+
+  template <typename Number>
+  template <typename Number2>
+  void Vector<Number>::add (const std::vector<size_type>    &indices,
+                            const ::dealii::Vector<Number2> &values)
+  {
+    Assert(indices.size()==values.size(),
+           ExcDimensionMismatch(indices.size(),values.size()));
+    add(indices.size(),&indices[0],values.begin());
+  }
+
+
+
+  template <typename Number>
+  template <typename Number2>
+  void Vector<Number>::add (const size_type n_indices,const size_type *indices,
+                            const Number2 *values)
+  {
+    for (size_type i=0; i<n_indices; ++i)
+      {
+        Assert(indices[i]<size(),ExcIndexRange(indices[i],0,size()));
+        Assert(numbers::is_finite(values[i]),ExcNumberNotFinite());
+
+        local_vector[indices[i]] += values[i];
+      }
+  }
+
+
+
+  template <typename Number>
+  void Vector<Number>::add(const Number s)
+  {
+    Assert(numbers::is_finite(s),ExcNumberNotFinite());
+
+    size_type size = local_vector.get_size();
+    for (size_type i=0; i<size; ++i)
+      local_vector[i] += s;
+  }
+
+
+
+  template <typename Number>
+  void Vector<Number>::add(const Vector<Number> &V)
+  {
+    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
+    local_vector.AddScale(V.paralution_vector(),1.);
+  }
+
+
+
+  template <typename Number>
+  void Vector<Number>::add(const Number a,const Vector<Number> &V)
+  {
+    Assert(numbers::is_finite(a),ExcNumberNotFinite());
+    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
+
+    local_vector.ScaleAddScale(1.,V.paralution_vector(),a);
+  }
+
+
+
+  template <typename Number>
+  void Vector<Number>::add(const Number a,const Vector<Number> &V,
+                           const Number b,const Vector<Number> &W)
+  {
+    Assert(numbers::is_finite(a),ExcNumberNotFinite());
+    Assert(numbers::is_finite(b),ExcNumberNotFinite());
+    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
+    Assert(size()==W.size(),ExcDimensionMismatch(size(),W.size()));
+
+    local_vector.ScaleAdd2(1.,V.paralution_vector(),a,W.paralution_vector(),b);
+  }
+
+
+
+  template <typename Number>
+  void Vector<Number>::sadd(const Number s, const Vector<Number> &V)
+  {
+    Assert(numbers::is_finite(s),ExcNumberNotFinite());
+    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
+
+    local_vector.ScaleAdd(s,V.paralution_vector());
+  }
+
+
+
+  template <typename Number>
+  void Vector<Number>::sadd(const Number s, const Number a, const Vector<Number> &V)
+  {
+    Assert(numbers::is_finite(s),ExcNumberNotFinite());
+    Assert(numbers::is_finite(a),ExcNumberNotFinite());
+    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
+
+    local_vector.ScaleAddScale(s,V.paralution_vector(),a);
+  }
+
+
+
+  template <typename Number>
+  void Vector<Number>::sadd(const Number s, const Number a, const Vector<Number> &V,
+                            const Number b, const Vector<Number> &W)
+  {
+    Assert(numbers::is_finite(s),ExcNumberNotFinite());
+    Assert(numbers::is_finite(a),ExcNumberNotFinite());
+    Assert(numbers::is_finite(b),ExcNumberNotFinite());
+    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
+    Assert(size()==W.size(),ExcDimensionMismatch(size(),W.size()));
+
+    local_vector.ScaleAdd2(s,V.paralution_vector(),a,W.paralution_vector(),b);
+  }
+
+
+
+  template <typename Number>
+  void Vector<Number>::sadd(const Number s, const Number a, const Vector<Number> &V,
+                            const Number b, const Vector<Number> &W, const Number c,
+                            const Vector<Number> &X)
+  {
+    Assert(numbers::is_finite(s),ExcNumberNotFinite());
+    Assert(numbers::is_finite(a),ExcNumberNotFinite());
+    Assert(numbers::is_finite(b),ExcNumberNotFinite());
+    Assert(numbers::is_finite(c),ExcNumberNotFinite());
+    Assert(size()==V.size(),ExcDimensionMismatch(size(),V.size()));
+    Assert(size()==W.size(),ExcDimensionMismatch(size(),W.size()));
+    Assert(size()==X.size(),ExcDimensionMismatch(size(),X.size()));
+
+    local_vector.ScaleAdd2(s,V.paralution_vector(),a,W.paralution_vector(),b);
+    local_vector.AddScale(X.paralution_vector(),c);
+  }
+
+
+
+  template <typename Number>
+  void Vector<Number>::equ (const Number a, const Vector<Number> &u)
+  {
+    Assert(numbers::is_finite(a), ExcNumberNotFinite());
+    Assert(size()==u.size(),ExcDimensionMismatch(size(),u.size()));
+
+    local_vector.ScaleAddScale(0.,u.paralution_vector(),a);
+  }
+
+
+
+  template <typename Number>
+  void Vector<Number>::equ (const Number a, const Vector<Number> &u,
+                            const Number b, const Vector<Number> &v)
+  {
+    Assert(numbers::is_finite(a), ExcNumberNotFinite());
+    Assert(numbers::is_finite(b), ExcNumberNotFinite());
+    Assert(size()==u.size(),ExcDimensionMismatch(size(),u.size()));
+    Assert(size()==v.size(),ExcDimensionMismatch(size(),v.size()));
+
+    local_vector.ScaleAdd2(0.,u.paralution_vector(),a,v.paralution_vector(),b);
+  }
+
+
+
+  template <typename Number>
+  void Vector<Number>::equ (const Number a, const Vector<Number> &u,
+                            const Number b, const Vector<Number> &v,
+                            const Number c, const Vector<Number> &w)
+  {
+    Assert(numbers::is_finite(a), ExcNumberNotFinite());
+    Assert(numbers::is_finite(b), ExcNumberNotFinite());
+    Assert(numbers::is_finite(c), ExcNumberNotFinite());
+    Assert(size()==u.size(),ExcDimensionMismatch(size(),u.size()));
+    Assert(size()==v.size(),ExcDimensionMismatch(size(),v.size()));
+    Assert(size()==w.size(),ExcDimensionMismatch(size(),w.size()));
+
+    local_vector.ScaleAdd2(0.,u.paralution_vector(),a,v.paralution_vector(),b);
+    local_vector.AddScale(w.paralution_vector(),c);
+  }
+}
+
+// Explicit instantiations
+namespace ParalutionWrappers
+{
+  template class Vector<float>;
+  template class Vector<double>;
+  template void Vector<float>::add<float> (const std::vector<size_type> &indices,
+                                           const std::vector<float>     &values);
+  template void Vector<float>::add<float> (const std::vector<size_type>  &indices,
+                                           const ::dealii::Vector<float> &values);
+  template void Vector<float>::add<float> (const size_type  n_elements,
+                                           const size_type *indices,
+                                           const float     *values);
+  template void Vector<float>::add<double> (const std::vector<size_type> &indices,
+                                            const std::vector<double>    &values);
+  template void Vector<float>::add<double> (const std::vector<size_type>   &indices,
+                                            const ::dealii::Vector<double> &values);
+  template void Vector<float>::add<double> (const size_type  n_elements,
+                                            const size_type *indices,
+                                            const double    *values);
+  template void Vector<double>::add<float> (const std::vector<size_type> &indices,
+                                            const std::vector<float>     &values);
+  template void Vector<double>::add<float> (const std::vector<size_type>  &indices,
+                                            const ::dealii::Vector<float> &values);
+  template void Vector<double>::add<float> (const size_type  n_elements,
+                                            const size_type *indices,
+                                            const float     *values);
+  template void Vector<double>::add<double> (const std::vector<size_type> &indices,
+                                             const std::vector<double>    &values);
+  template void Vector<double>::add<double> (const std::vector<size_type>  &indices,
+                                             const ::dealii::Vector<double> &values);
+  template void Vector<double>::add<double> (const size_type  n_elements,
+                                             const size_type *indices,
+                                             const double    *values);
+}
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif // DEAL_II_WITH_PARALUTION
index fdf42a3285e83eb1ecfe328dfe2cc3d0b6e70b93..73311620d68b27f4caa7ed08257fc80815de008b 100644 (file)
@@ -27,8 +27,10 @@ void check()
   ParalutionWrappers::Vector<Number> vector_1;
   ParalutionWrappers::Vector<Number> vector_2(100);
 
-  deallog << vector_1.size() <<std::endl;
-  deallog << vector_2.size() <<std::endl;
+  deallog << vector_1.size() << std::endl;
+  deallog << vector_2.size() << std::endl;
+  vector_2.clear();
+  deallog << vector_2.size() << std::endl;
 }
 
 int main()
index 18cb6fab4e1cda00d0caa3ee55beedf117483aae..23bbdba40a8f72cac0c46ec52524267bd076305d 100644 (file)
@@ -2,6 +2,9 @@
 DEAL::0
 DEAL::100
 DEAL::0
+DEAL::0
 DEAL::100
 DEAL::0
+DEAL::0
 DEAL::100
+DEAL::0

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.