]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add a functor interface to ReadWriteVector. 2657/head
authorBruno Turcksin <bruno.turcksin@gmail.com>
Tue, 31 May 2016 17:33:15 +0000 (13:33 -0400)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Tue, 31 May 2016 21:31:30 +0000 (17:31 -0400)
include/deal.II/grid/filtered_iterator.h
include/deal.II/lac/read_write_vector.h
include/deal.II/lac/read_write_vector.templates.h
tests/lac/readwritevector_functor.cc [new file with mode: 0644]
tests/lac/readwritevector_functor.output [new file with mode: 0644]

index 59fa12fb837173170daf8c8d5fcc5d1141b91514..ab6ba6122954f09183ff468b4a3c71bb09347613 100644 (file)
@@ -697,7 +697,7 @@ private:
 
     /**
      * Abstract function which in derived classes denotes the evaluation of
-     * the predicate on the give iterator.
+     * the predicate on the given iterator.
      */
     virtual bool operator () (const BaseIterator &bi) const = 0;
 
index 1841f4a2d17605532e4ea11f0aaabcfa0c00c170..00c1d7d31d14c18fbf0b8a7bbd7861ccdfd0826c 100644 (file)
@@ -195,6 +195,24 @@ namespace LinearAlgebra
     void reinit (const IndexSet &locally_stored_indices,
                  const bool      omit_zeroing_entries = false);
 
+#ifdef DEAL_II_WITH_CXX11
+    /**
+     * Apply the functor @p func to each element of the vector. The functor
+     * should look like
+     * <code>
+     * struct Functor
+     * {
+     *    void operator() (Number &value);
+     * };
+     * </code>
+     *
+     * @note This function requires C++11 and read_write_vector.templates.h
+     * needs to be included.
+     */
+    template <typename Functor>
+    void apply(const Functor &func);
+#endif
+
     /**
      * Swap the contents of this vector and the other vector @p v. One could
      * do this operation with a temporary variable and copying over the data
@@ -517,6 +535,42 @@ namespace LinearAlgebra
      * Make all other ReadWriteVector types friends.
      */
     template <typename Number2> friend class ReadWriteVector;
+
+#ifdef DEAL_II_WITH_CXX11
+  private:
+    /**
+     * This class provides a wrapper around a Functor which acts on
+     * single elements of the vector. This is necessary to use
+     * tbb::parallel_for which requires a TBBForFunctor.
+     */
+    template <typename Functor>
+    class FunctorTemplate
+    {
+    public:
+      /**
+       * Constructor. Take a functor and store a copy of it.
+       */
+      FunctorTemplate(ReadWriteVector<Number> &parent,
+                      const Functor &functor);
+
+      /**
+       * Evaluate the element with the stored copy of the functor.
+       */
+      virtual void operator() (const size_type begin,
+                               const size_type end);
+
+    private:
+      /**
+       * Alias to the ReadWriteVector object that owns the FunctorTemplate.
+       */
+      ReadWriteVector &parent;
+
+      /**
+       * Copy of the functor.
+       */
+      const Functor &functor;
+    };
+#endif
   };
 
   /*@}*/
@@ -790,6 +844,33 @@ namespace LinearAlgebra
       }
   }
 
+
+
+#ifdef DEAL_II_WITH_CXX11
+  template <typename Number>
+  template <typename Functor>
+  inline
+  ReadWriteVector<Number>::FunctorTemplate<Functor>::FunctorTemplate(
+    ReadWriteVector<Number> &parent,
+    const Functor &functor)
+    :
+    parent(parent),
+    functor(functor)
+  {}
+
+
+
+  template <typename Number>
+  template <typename Functor>
+  void
+  ReadWriteVector<Number>::FunctorTemplate<Functor>::operator() (const size_type begin,
+      const size_type end)
+  {
+    for (size_type i=begin; i<end; ++i)
+      functor(parent.val[i]);
+  }
+#endif
+
 #endif  // ifndef DOXYGEN
 
 } // end of namespace LinearAlgebra
index aa679f412d61c0a7c03f76fe7dc67048a4408097..11d5dfbe5540ddc3ec467926a587ef54975d9ee6 100644 (file)
@@ -124,6 +124,19 @@ namespace LinearAlgebra
 
 
 
+#ifdef DEAL_II_WITH_CXX11
+  template <typename Number>
+  template <typename Functor>
+  void
+  ReadWriteVector<Number>::apply(const Functor &func)
+  {
+    FunctorTemplate<Functor> functor(*this, func);
+    internal::parallel_for(functor, n_elements(), thread_loop_partitioner);
+  }
+#endif
+
+
+
   template <typename Number>
   ReadWriteVector<Number> &
   ReadWriteVector<Number>::operator= (const ReadWriteVector<Number> &in_vector)
diff --git a/tests/lac/readwritevector_functor.cc b/tests/lac/readwritevector_functor.cc
new file mode 100644 (file)
index 0000000..5c33474
--- /dev/null
@@ -0,0 +1,55 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 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.
+//
+// ---------------------------------------------------------------------
+
+// Check the functor interface
+
+
+#include "../tests.h"
+#include <deal.II/base/index_set.h>
+#include <deal.II/base/logstream.h>
+#include <deal.II/lac/read_write_vector.h>
+#include <deal.II/lac/read_write_vector.templates.h>
+
+#include <fstream>
+
+struct Functor
+{
+  void operator() (double &value) const
+  {
+    value *= 2.;
+  }
+};
+
+void test()
+{
+  const unsigned int size = 25;
+  LinearAlgebra::ReadWriteVector<double> vector(size);
+  for (unsigned int i=0; i<size; ++i)
+    vector[i] = i;
+
+  Functor functor;
+  vector.apply(functor);
+  for (unsigned int i=0; i<size; ++i)
+    deallog << vector[i] << std::endl;
+}
+
+int main()
+{
+  std::ofstream logfile("output");
+  deallog.attach(logfile);
+  test();
+
+  return 0;
+}
diff --git a/tests/lac/readwritevector_functor.output b/tests/lac/readwritevector_functor.output
new file mode 100644 (file)
index 0000000..55ccd41
--- /dev/null
@@ -0,0 +1,26 @@
+
+DEAL::0.00000
+DEAL::2.00000
+DEAL::4.00000
+DEAL::6.00000
+DEAL::8.00000
+DEAL::10.0000
+DEAL::12.0000
+DEAL::14.0000
+DEAL::16.0000
+DEAL::18.0000
+DEAL::20.0000
+DEAL::22.0000
+DEAL::24.0000
+DEAL::26.0000
+DEAL::28.0000
+DEAL::30.0000
+DEAL::32.0000
+DEAL::34.0000
+DEAL::36.0000
+DEAL::38.0000
+DEAL::40.0000
+DEAL::42.0000
+DEAL::44.0000
+DEAL::46.0000
+DEAL::48.0000

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.