]> https://gitweb.dealii.org/ - dealii.git/commitdiff
added clear() and Tvmult for SparseVanka 17383/head
authorChayapol Chaoveeraprasit <chayapol.com@gmail.com>
Wed, 24 Jul 2024 14:08:00 +0000 (17:08 +0300)
committerChayapol Chaoveeraprasit <chayapol.com@gmail.com>
Tue, 20 Aug 2024 12:39:27 +0000 (20:39 +0800)
undo copyright notices & fix wording issues

Add a changelog entry

fix formatting of sparse_vanka.h

fix typos

add spacing

doc/news/changes/minor/20240814chaycha [new file with mode: 0644]
include/deal.II/lac/sparse_vanka.h
include/deal.II/lac/sparse_vanka.templates.h
source/lac/sparse_vanka.cc

diff --git a/doc/news/changes/minor/20240814chaycha b/doc/news/changes/minor/20240814chaycha
new file mode 100644 (file)
index 0000000..e945e34
--- /dev/null
@@ -0,0 +1,4 @@
+New: Added SparseVanka::Tvmult() and SparseVanka::clear()
+SparseVanka can now be passed to MGSmootherPrecondition to be used as a multigrid smoother.
+<br>
+(Chayapol Chaoveeraprasit, 2024/08/14)
\ No newline at end of file
index 2d6a3e2b0d76eca46c0ff47077b103405221d76b..63ede4c36f6a0244963eb45a73eb99eab56a85bc 100644 (file)
@@ -15,8 +15,6 @@
 #ifndef dealii_sparse_vanka_h
 #define dealii_sparse_vanka_h
 
-
-
 #include <deal.II/base/config.h>
 
 #include <deal.II/base/multithread_info.h>
@@ -216,6 +214,12 @@ public:
   void
   Tvmult(Vector<number2> &dst, const Vector<number2> &src) const;
 
+  /**
+   * Clear all memory.
+   */
+  void
+  clear();
+
   /**
    * Return the dimension of the codomain (or range) space. Note that the
    * matrix is of dimension $m \times n$.
@@ -240,9 +244,11 @@ protected:
   /**
    * Apply the inverses corresponding to those degrees of freedom that have a
    * @p true value in @p dof_mask, to the @p src vector and move the result
-   * into @p dst. Actually, only values for allowed indices are written to @p
-   * dst, so the application of this function only does what is announced in
-   * the general documentation if the given mask sets all values to zero
+   * into @p dst. The transpose of the inverse is applied instead if
+   * @p transpose equals true. Actually, only values for allowed indices are
+   * written to @p dst, so the application of this function only does what is
+   * announced in the general documentation if the given mask sets all values
+   * to zero.
    *
    * The reason for providing the mask anyway is that in derived classes we
    * may want to apply the preconditioner to parts of the matrix only, in
@@ -258,10 +264,12 @@ protected:
    * The @p vmult of this class of course calls this function with a null
    * pointer
    */
+
   template <typename number2>
   void
   apply_preconditioner(Vector<number2>               &dst,
                        const Vector<number2>         &src,
+                       const bool                     transpose = false,
                        const std::vector<bool> *const dof_mask = nullptr) const;
 
   /**
@@ -511,6 +519,14 @@ public:
   void
   vmult(Vector<number2> &dst, const Vector<number2> &src) const;
 
+
+  /**
+   * Apply the transpose preconditioner.
+   */
+  template <typename number2>
+  void
+  Tvmult(Vector<number2> &dst, const Vector<number2> &src) const;
+
   /**
    * Determine an estimate for the memory consumption (in bytes) of this
    * object.
@@ -564,15 +580,6 @@ SparseVanka<number>::n() const
   return _n;
 }
 
-template <typename number>
-template <typename number2>
-inline void
-SparseVanka<number>::Tvmult(Vector<number2> & /*dst*/,
-                            const Vector<number2> & /*src*/) const
-{
-  AssertThrow(false, ExcNotImplemented());
-}
-
 #endif // DOXYGEN
 
 DEAL_II_NAMESPACE_CLOSE
index 168419546ef80ab3fd87e292241a05922b8e690a..34fda24b17ad89d835a8d48ae21d6a6d7dd9018a 100644 (file)
@@ -215,9 +215,42 @@ SparseVanka<number>::vmult(Vector<number2>       &dst,
   dst = 0;
   // then pass on to the function
   // that actually does the work
-  apply_preconditioner(dst, src);
+  apply_preconditioner(dst, src, false);
 }
 
+template <typename number>
+template <typename number2>
+void
+SparseVanka<number>::Tvmult(Vector<number2>       &dst,
+                            const Vector<number2> &src) const
+{
+  Assert(matrix != nullptr, ExcNotInitialized());
+  Assert(selected != nullptr, ExcNotInitialized());
+
+  // first set output vector to zero
+  dst = 0;
+  // then pass on to the function
+  // that actually does the work
+  apply_preconditioner(dst, src, true);
+}
+
+
+template <typename number>
+void
+SparseVanka<number>::clear()
+{ // Clear the inverses vector and deallocate memory
+  inverses.clear();
+
+  // Reset the matrix pointer
+  matrix = nullptr;
+
+  // Reset the selected pointer
+  selected = nullptr;
+
+  // Reset the sizes
+  _m = 0;
+  _n = 0;
+}
 
 template <typename number>
 template <typename number2>
@@ -225,6 +258,7 @@ void
 SparseVanka<number>::apply_preconditioner(
   Vector<number2>               &dst,
   const Vector<number2>         &src,
+  const bool                     transpose,
   const std::vector<bool> *const dof_mask) const
 {
   Assert(dst.size() == src.size(),
@@ -336,7 +370,10 @@ SparseVanka<number>::apply_preconditioner(
           }
 
         // apply preconditioner
-        inverses[row]->vmult(x, b);
+        if (transpose)
+          inverses[row]->Tvmult(x, b);
+        else
+          inverses[row]->vmult(x, b);
 
         // Distribute new values
         for (std::map<size_type, size_type>::const_iterator is =
@@ -567,8 +604,24 @@ SparseBlockVanka<number>::vmult(Vector<number2>       &dst,
 
   Threads::TaskGroup<> tasks;
   for (unsigned int block = 0; block < n_blocks; ++block)
-    tasks += Threads::new_task([&, block]() {
-      this->apply_preconditioner(dst, src, &dof_masks[block]);
+    tasks += Threads::new_task([&, block] {
+      this->apply_preconditioner(dst, src, false, &dof_masks[block]);
+    });
+  tasks.join_all();
+}
+
+template <typename number>
+template <typename number2>
+void
+SparseBlockVanka<number>::Tvmult(Vector<number2>       &dst,
+                                 const Vector<number2> &src) const
+{
+  dst = 0;
+
+  Threads::TaskGroup<> tasks;
+  for (unsigned int block = 0; block < n_blocks; ++block)
+    tasks += Threads::new_task([&, block] {
+      this->apply_preconditioner(dst, src, true, &dof_masks[block]);
     });
   tasks.join_all();
 }
index 9cb74bb9cbece2640edbe8645904ac821efa6d60..d49ba83de85497f815b9d5628ba7ddf43e103e10 100644 (file)
@@ -27,6 +27,12 @@ SparseVanka<double>::vmult<float>(Vector<float>       &dst,
 template void
 SparseVanka<double>::vmult<double>(Vector<double>       &dst,
                                    const Vector<double> &src) const;
+template void
+SparseVanka<double>::Tvmult<float>(Vector<float>       &dst,
+                                   const Vector<float> &src) const;
+template void
+SparseVanka<double>::Tvmult<double>(Vector<double>       &dst,
+                                    const Vector<double> &src) const;
 
 
 template class SparseBlockVanka<float>;
@@ -38,5 +44,11 @@ SparseBlockVanka<double>::vmult<float>(Vector<float>       &dst,
 template void
 SparseBlockVanka<double>::vmult<double>(Vector<double>       &dst,
                                         const Vector<double> &src) const;
+template void
+SparseBlockVanka<double>::Tvmult<float>(Vector<float>       &dst,
+                                        const Vector<float> &src) const;
+template void
+SparseBlockVanka<double>::Tvmult<double>(Vector<double>       &dst,
+                                         const Vector<double> &src) const;
 
 DEAL_II_NAMESPACE_CLOSE

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.