]> https://gitweb.dealii.org/ - dealii.git/commitdiff
base: Implement a Lazy<T> wrapper
authorMatthias Maier <tamiko@43-1.org>
Tue, 24 Oct 2023 16:33:04 +0000 (11:33 -0500)
committerMatthias Maier <tamiko@43-1.org>
Fri, 27 Oct 2023 07:01:36 +0000 (02:01 -0500)
Implement an optimized, thin wrapper that provides a convenient
mechanism for lazy initialization of the contained object on first use.
The class ensures that on-demand initialization of some expensive data
structure happens (a) in a thread-safe manner, and that (b) subsequent
checks in hot paths are cheap:

 - The class implements proper copy and move semantics inherited from
   std::optional. This is particularly desirable for our typical use
   case in deal.II where we want to initialize some payload exactly once
   (and simply copy or move it along with the container).

 - The "passive" call, i.e., when the object is already initialized is
   much cheaper than a call to std::call_once, or thread synchronization
   with a std::shared_mutex (because the initialize() or
   value_or_initialize() functions are inlined and consist of a mere
   check of a boolean followed by a conditional jump).

include/deal.II/base/lazy.h [new file with mode: 0644]

diff --git a/include/deal.II/base/lazy.h b/include/deal.II/base/lazy.h
new file mode 100644 (file)
index 0000000..e2c5f22
--- /dev/null
@@ -0,0 +1,378 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2023 - 2023 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.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii_lazy_h
+#define dealii_lazy_h
+
+
+#include <deal.II/base/config.h>
+
+#include <deal.II/base/exceptions.h>
+#include <deal.II/base/mutex.h>
+
+#include <atomic>
+#include <mutex>
+#include <optional>
+
+DEAL_II_NAMESPACE_OPEN
+
+/**
+ * @addtogroup threads
+ * @{
+ */
+
+/**
+ * This class is a wrapper that provides a convenient mechanism for lazy
+ * initialization of the contained object on first use. The class ensures
+ * that on-demand initialization of some expensive data structure happens
+ * (a) exactly once in a thread-safe manner, and that (b) subsequent checks
+ * in hot paths are cheap.
+ *
+ * Lazy<T> is closely modeled after the `std::optional` interface providing
+ * a `reset()` and `value()` method, but also and extending it with two
+ * methods: `ensure_initialized(creator)` which, as the name suggests,
+ * ensures that the contained object is properly initialized. If the
+ * `Lazy<T>` happens happens to contain no value yet, it initializes the
+ * wrapped object by calling the `creator()` function object and storing
+ * the return value. In addition a `value_or_initialize(creator)` function
+ * is provided that, similarly, ensures that the object is properly
+ * initialized and then returns a reference to the contained value.
+ *
+ * Example usage could look like the following, where the `FE` class stores a
+ * matrix that is expensive to compute and so we do not want to do it unless it
+ * is actually needed. As a consequence, rather than storing a matrix, we store
+ * a `Lazy<FullMatrix<double>>` that by default is empty; whenever the matrix is
+ * first requested, we create it and store it for later reuse:
+ * ```
+ * template<...>
+ * class FE
+ * {
+ * public:
+ *   const FullMatrix<double> & get_prolongation_matrix() const
+ *   {
+ *     prolongation_matrix.ensure_initialized([&](){
+ *       // Some expensive operation initializing the prolongation matrix
+ *       // that we only want to perform once and when necessary.
+ *       });
+ *     return prolongation_matrix.value();
+ *   }
+ *
+ * private:
+ *   Lazy<FullMatrix<double>> prolongation_matrix;
+ * };
+ * ```
+ */
+template <typename T>
+class Lazy
+{
+public:
+  /**
+   * Default Constructor.
+   */
+  Lazy();
+
+
+  /**
+   * Copy constructor.
+   */
+  Lazy(const Lazy &other);
+
+
+  /**
+   * Move constructor.
+   */
+  Lazy(Lazy &&other) noexcept;
+
+
+  /**
+   * Copy assignment.
+   */
+  Lazy &
+  operator=(const Lazy &other);
+
+
+  /**
+   * Move assignment.
+   */
+  Lazy &
+  operator=(Lazy &&other) noexcept;
+
+
+  /**
+   * Reset the Lazy<T> object to an uninitialized state.
+   */
+  void
+  reset() noexcept;
+
+
+  /**
+   * Initialize the wrapped object.
+   *
+   * If the contained object is already initialized this function simply
+   * returns and does nothing.
+   *
+   * If, instead, the object has not yet been initialized then the @p
+   * creator function object (oftentimes a lambda function) is called to
+   * initialize the contained object.
+   *
+   * This operation is thread safe: The ensure_initialized() method
+   * guarantees that the creator function object is only called once on one
+   * of the calling threads and that after completion the initialization
+   * result (which is stored in the std::optional) is visible on all
+   * threads.
+   */
+  template <typename Callable>
+  void
+  ensure_initialized(const Callable &creator) const;
+
+
+  /**
+   * Return a const reference to the contained object.
+   *
+   * @pre The object has been initialized with a call to
+   * ensure_initialized() or value_or_initialized().
+   */
+  const T &
+  value() const;
+
+
+  /**
+   * Return a reference to the contained object.
+   *
+   * @pre The object has been initialized with a call to
+   * ensure_initialized() or value_or_initialized().
+   */
+  T &
+  value();
+
+
+  /**
+   * If the underlying object is initialized the function simply returns a
+   * const reference to the contained value. Otherwise, the @p creator()
+   * function object is called to initialize the object first.
+   *
+   * This function mimics the syntax of the std::optional<T> interface and
+   * is functionally equivalent to calling ensure_initialized() followed by
+   * value().
+   *
+   * @note This method can be called from a context where the Lazy<T>
+   * wrapper itself is marked const. FIXME
+   *
+   * @post The underlying object is initialized, meaning, has_value()
+   * returns true.
+   */
+  template <typename Callable>
+  const T &
+  value_or_initialize(const Callable &creator) const;
+
+
+  /**
+   * Variant of above function that returns a non-const reference.
+   */
+  template <typename Callable>
+  DEAL_II_ALWAYS_INLINE inline T &
+  value_or_initialize(const Callable &creator);
+
+
+private:
+  /**
+   * The lazily initialized object stored as a std::optional<T>.
+   */
+  mutable std::optional<T> object;
+
+
+  /**
+   * An atomic bool used for checking whether the object is initialized in
+   * a thread-safe manner.
+   */
+  mutable std::atomic<bool> object_is_initialized;
+
+
+  /**
+   * A mutex used for protecting the initialization of the object.
+   */
+  mutable Threads::Mutex initialization_mutex;
+};
+
+/**
+ * @}
+ */
+
+
+// ------------------------------- inline functions --------------------------
+
+
+template <typename T>
+inline Lazy<T>::Lazy()
+  : object_is_initialized(false)
+{}
+
+
+template <typename T>
+inline Lazy<T>::Lazy(const Lazy &other)
+  : object(other.object)
+{
+  object_is_initialized.store(other.object_is_initialized.load());
+}
+
+
+template <typename T>
+inline Lazy<T>::Lazy(Lazy &&other) noexcept
+  : object(other.object)
+{
+  object_is_initialized.store(other.object_is_initialized.load());
+}
+
+
+template <typename T>
+inline Lazy<T> &
+Lazy<T>::operator=(const Lazy &other)
+{
+  object = other.object;
+  object_is_initialized.store(other.object_is_initialized.load());
+  return *this;
+}
+
+
+template <typename T>
+inline Lazy<T> &
+Lazy<T>::operator=(Lazy &&other) noexcept
+{
+  object = other.object;
+  object_is_initialized.store(other.object_is_initialized.load());
+  return *this;
+}
+
+
+template <typename T>
+inline void
+Lazy<T>::reset() noexcept
+{
+  object_is_initialized.store(false);
+  object.reset();
+}
+
+
+template <typename T>
+template <typename Callable>
+inline DEAL_II_ALWAYS_INLINE void
+Lazy<T>::ensure_initialized(const Callable &creator) const
+{
+  //
+  // Use Schmidt's double checking [1] for checking and initializing the
+  // object.
+  //
+  // [1] https://en.wikipedia.org/wiki/Double-checked_locking
+  //
+
+  //
+  // Check the object_is_initialized atomic with "acquire" semantics [1].
+  //
+  // This ensures that (a) all subsequent reads (of the object) are
+  // ordered after this check, and that (b) all writes to the object
+  // before the atomic bool was set to true with "release" semantics are
+  // visible on this thread.
+  //
+  // [1]
+  // https://en.cppreference.com/w/cpp/atomic/memory_order#Release-Acquire_ordering
+  //
+  if (!object_is_initialized.load(std::memory_order_acquire))
+#ifdef DEAL_II_HAVE_CXX20
+    [[unlikely]]
+#endif
+    {
+      std::lock_guard<std::mutex> lock(initialization_mutex);
+
+      //
+      // Check again. If this thread won the race to the lock then we
+      // initialize the object. Otherwise another thread has already
+      // initialized the object and flipped the object_is_initialized
+      // bit. (Here, the initialization_mutex ensures consistent ordering
+      // with a memory fence, so we will observe the updated bool without
+      // acquire semantics.)
+      //
+      if (!object_is_initialized.load(std::memory_order_relaxed))
+        {
+          object = std::move(creator());
+
+          //
+          // Flip the object_is_initialized boolean with "release"
+          // semantics [1].
+          //
+          // This ensures that the above move is visible on all threads
+          // before checking the atomic bool with acquire semantics.
+          //
+          object_is_initialized.store(true, std::memory_order_release);
+        }
+    }
+
+  Assert(
+    object.has_value(),
+    dealii::ExcMessage(
+      "The internal std::optional<T> object does not contain a valid object "
+      "even though we have just initialized it."));
+}
+
+
+template <typename T>
+inline DEAL_II_ALWAYS_INLINE const T &
+Lazy<T>::value() const
+{
+  Assert(
+    object.has_value(),
+    dealii::ExcMessage(
+      "value() has been called but the contained object has not been "
+      "initialized. Did you forget to call 'ensure_initialized()' first?"));
+
+  return object.value();
+}
+
+
+template <typename T>
+inline DEAL_II_ALWAYS_INLINE T &
+Lazy<T>::value()
+{
+  Assert(
+    object.has_value(),
+    dealii::ExcMessage(
+      "value() has been called but the contained object has not been "
+      "initialized. Did you forget to call 'ensure_initialized()' first?"));
+
+  return object.value();
+}
+
+
+template <typename T>
+template <typename Callable>
+inline DEAL_II_ALWAYS_INLINE const T &
+Lazy<T>::value_or_initialize(const Callable &creator) const
+{
+  ensure_initialized(creator);
+  return object.value();
+}
+
+
+template <typename T>
+template <typename Callable>
+inline DEAL_II_ALWAYS_INLINE T &
+Lazy<T>::value_or_initialize(const Callable &creator)
+{
+  ensure_initialized(creator);
+  return object.value();
+}
+
+
+DEAL_II_NAMESPACE_CLOSE
+#endif

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.