--- /dev/null
+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
#ifndef dealii_sparse_vanka_h
#define dealii_sparse_vanka_h
-
-
#include <deal.II/base/config.h>
#include <deal.II/base/multithread_info.h>
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$.
/**
* 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
* 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;
/**
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.
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
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>
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(),
}
// 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 =
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();
}
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>;
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