From f69389b993dca6797f4afa24fb62433d3e422e20 Mon Sep 17 00:00:00 2001 From: steigemann Date: Fri, 26 Oct 2012 16:11:16 +0000 Subject: [PATCH] Interface to PETSc::PreconditionNone and ParaSails Preconditioner from the Hypre Suite git-svn-id: https://svn.dealii.org/trunk@27208 0785d39b-7218-0410-832d-ea1e28bc413d --- .../include/deal.II/lac/petsc_precondition.h | 218 +++++++++++++++++- 1 file changed, 216 insertions(+), 2 deletions(-) diff --git a/deal.II/include/deal.II/lac/petsc_precondition.h b/deal.II/include/deal.II/lac/petsc_precondition.h index af3e05acab..6e4db17c10 100644 --- a/deal.II/include/deal.II/lac/petsc_precondition.h +++ b/deal.II/include/deal.II/lac/petsc_precondition.h @@ -804,8 +804,6 @@ namespace PETScWrappers bool output_details; }; - - /** * Empty Constructor. You need to call * initialize() before using this @@ -843,6 +841,222 @@ namespace PETScWrappers */ AdditionalData additional_data; }; + + + +/** + * A class that implements the interface to use the ParaSails sparse + * approximate inverse preconditioner from the HYPRE suite. Note that + * PETSc has to be configured with HYPRE (e.g. with --download-hypre=1). + * + * ParaSails uses least-squares minimization to compute a sparse + * approximate inverse. The sparsity pattern used is the pattern + * of a power of a sparsified matrix. ParaSails also uses a post-filtering + * technique to reduce the cost of applying the preconditioner. + * + * ParaSails solves symmetric positive definite (SPD) problems + * using a factorized SPD preconditioner and can also solve + * general (nonsymmetric and/or indefinite) problems with a + * nonfactorized preconditioner. The problem type has to be + * set in @p AdditionalData. + * + * The preconditioner does support parallel distributed computations. + * + * @ingroup PETScWrappers + * @author Martin Steigemann, 2012 + */ + class PreconditionParaSails : public PreconditionerBase + { + public: + /** + * Standardized data struct to + * pipe additional flags to the + * preconditioner. + */ + struct AdditionalData + { + /** + * Constructor. + */ + AdditionalData ( + const unsigned int symmetric = 1, + const unsigned int n_levels = 1, + const double threshold = 0.1, + const double filter = 0.05, + const bool output_details = false + ); + + /** + * This parameter specifies the + * type of problem to solve: + * + * Default is symmetric = 1. + */ + unsigned int symmetric; + + /** + * Patterns of powers of sparsified + * matrices can be approximate the + * pattern of large entries in the + * inverse of a matrix. Powers of a + * + * The sparsity pattern used for the + * apprximate inverse is the pattern + * of a power B^m where B + * has been sparsified from the given + * matrix A, n_level + * is equal to m+1. Default + * value is n_levels = 1. + */ + unsigned int n_levels; + + /** + * Sparsification is performed by + * dropping nonzeros which are smaller + * than thresh in magnitude. + * Lower values of thresh + * lead to more accurate, but also more + * expensive preconditioners. Default + * value is thresh = 0.1. Setting + * thresh < 0 a threshold + * is selected automatically, such that + * -thresh represents the + * fraction of nonzero elements that are + * dropped. For example, if thresh = -0.9, + * then B will contain about + * ten percent of the nonzeros of the + * given matrix A. + */ + double threshold; + + /** + * Filtering is a post-processing procedure, + * filter represents a fraction + * of nonzero elements that are dropped + * after creating the approximate inverse + * sparsity pattern. Default value is + * filter = 0.05. Setting filter < 0 + * a value is selected automatically, such + * that -filter represents the + * fraction of nonzero elements that are + * dropped. For example, if thresh = -0.9, + * then about 90 percent of the entries + * in the computed approximate inverse are + * dropped. + */ + double filter; + + /** + * Setting this flag to true + * produces output from HYPRE, + * when the preconditioner + * is constructed. + */ + bool output_details; + }; + + + + /** + * Empty Constructor. You need to call + * initialize() before using this + * object. + */ + PreconditionParaSails (); + + /** + * Constructor. Take the matrix which + * is used to form the preconditioner, + * and additional flags if there are + * any. + */ + PreconditionParaSails (const MatrixBase &matrix, + const AdditionalData &additional_data = AdditionalData()); + + /** + * Initializes the preconditioner + * object and calculate all data that + * is necessary for applying it in a + * solver. This function is + * automatically called when calling + * the constructor with the same + * arguments and is only used if you + * create the preconditioner without + * arguments. + */ + void initialize (const MatrixBase &matrix, + const AdditionalData &additional_data = AdditionalData()); + + private: + /** + * Store a copy of the flags for this + * particular preconditioner. + */ + AdditionalData additional_data; + }; + + + +/** + * A class that implements a non-preconditioned method. + * + * @ingroup PETScWrappers + * @author Martin Steigemann, 2012 + */ + class PreconditionNone : public PreconditionerBase + { + public: + /** + * Standardized data struct to + * pipe additional flags to the + * preconditioner. + */ + struct AdditionalData + {}; + + /** + * Empty Constructor. You need to call + * initialize() before using this + * object. + */ + PreconditionNone (); + + /** + * Constructor. Take the matrix which + * is used to form the preconditioner, + * and additional flags if there are + * any. The matrix is completely + * ignored in computations. + */ + PreconditionNone (const MatrixBase &matrix, + const AdditionalData &additional_data = AdditionalData()); + + /** + * Initializes the preconditioner + * object and calculate all data that + * is necessary for applying it in a + * solver. This function is + * automatically called when calling + * the constructor with the same + * arguments and is only used if you + * create the preconditioner without + * arguments. The matrix is completely + * ignored in computations. + */ + void initialize (const MatrixBase &matrix, + const AdditionalData &additional_data = AdditionalData()); + + private: + /** + * Store a copy of the flags for this + * particular preconditioner. + */ + AdditionalData additional_data; + }; } -- 2.39.5