/usr/include/deal.II/lac/petsc_precondition.h is in libdeal.ii-dev 8.4.2-2+b1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 | // ---------------------------------------------------------------------
//
// Copyright (C) 2004 - 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.
//
// ---------------------------------------------------------------------
#ifndef dealii__petsc_precondition_h
#define dealii__petsc_precondition_h
#include <deal.II/base/config.h>
#ifdef DEAL_II_WITH_PETSC
# include <deal.II/lac/exceptions.h>
# include <petscpc.h>
DEAL_II_NAMESPACE_OPEN
namespace PETScWrappers
{
// forward declarations
class MatrixBase;
class VectorBase;
class SolverBase;
/**
* Base class for preconditioner classes using the PETSc functionality. The
* classes in this hierarchy don't do a whole lot, except for providing a
* function that sets the preconditioner and certain parameters on the
* preconditioning context of the solver. These classes are basically here
* only to allow a similar interface as already used for the deal.II solver
* and preconditioner classes.
*
* Note that derived classes only provide interfaces to the relevant
* functionality of PETSc. PETSc does not implement all preconditioners for
* all matrix types. In particular, some preconditioners are not going to
* work for parallel jobs, such as for example the ILU preconditioner.
*
* @ingroup PETScWrappers
* @author Wolfgang Bangerth, Timo Heister, 2004, 2011
*/
class PreconditionerBase
{
public:
/**
* Constructor.
*/
PreconditionerBase ();
/**
* Destructor.
*/
virtual ~PreconditionerBase ();
/**
* Apply the preconditioner once to the given src vector.
*/
void vmult (VectorBase &dst,
const VectorBase &src) const;
/**
* Gives access to the underlying PETSc object.
*/
const PC &get_pc () const;
protected:
/**
* the PETSc preconditioner object
*/
PC pc;
/**
* A pointer to the matrix that acts as a preconditioner.
*/
Mat matrix;
/**
* Internal function to create the PETSc preconditioner object. Fails if
* called twice.
*/
void create_pc ();
/**
* Conversion operator to get a representation of the matrix that
* represents this preconditioner. We use this inside the actual solver,
* where we need to pass this matrix to the PETSc solvers.
*/
operator Mat () const;
/**
* Make the solver class a friend, since it needs to call the conversion
* operator.
*/
friend class SolverBase;
};
/**
* A class that implements the interface to use the PETSc Jacobi
* preconditioner.
*
* See the comment in the base class
* @ref PreconditionerBase
* for when this preconditioner may or may not work.
*
* @ingroup PETScWrappers
* @author Wolfgang Bangerth, Timo Heister, 2004, 2011
*/
class PreconditionJacobi : 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.
*/
PreconditionJacobi ();
/**
* Constructor. Take the matrix which is used to form the preconditioner,
* and additional flags if there are any.
*/
PreconditionJacobi (const MatrixBase &matrix,
const AdditionalData &additional_data = AdditionalData());
/**
* Same as above but without setting a matrix to form the preconditioner.
* Intended to be used with SLEPc objects.
*/
PreconditionJacobi (const MPI_Comm communicator,
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());
protected:
/**
* Store a copy of the flags for this particular preconditioner.
*/
AdditionalData additional_data;
/**
* Initializes the preconditioner object without knowing a particular
* matrix. This function sets up appropriate parameters to the underlying
* PETSc object after it has been created.
*/
void initialize();
};
/**
* A class that implements the interface to use the PETSc Block Jacobi
* preconditioner. The blocking structure of the matrix is determined by the
* association of degrees of freedom to the individual processors in an MPI-
* parallel job. If you use this preconditioner on a sequential job (or an
* MPI job with only one process) then the entire matrix is the only block.
*
* By default, PETSc uses an ILU(0) decomposition of each diagonal block of
* the matrix for preconditioning. This can be changed, as is explained in
* the relevant section of the PETSc manual, but is not implemented here.
*
* See the comment in the base class
* @ref PreconditionerBase
* for when this preconditioner may or may not work.
*
* @ingroup PETScWrappers
* @author Wolfgang Bangerth, Timo Heister, 2004, 2011
*/
class PreconditionBlockJacobi : 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.
*/
PreconditionBlockJacobi ();
/**
* Constructor. Take the matrix which is used to form the preconditioner,
* and additional flags if there are any.
*/
PreconditionBlockJacobi (const MatrixBase &matrix,
const AdditionalData &additional_data = AdditionalData());
/**
* Same as above but without setting a matrix to form the preconditioner.
* Intended to be used with SLEPc objects.
*/
PreconditionBlockJacobi (const MPI_Comm communicator,
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());
protected:
/**
* Store a copy of the flags for this particular preconditioner.
*/
AdditionalData additional_data;
/**
* Initializes the preconditioner object without knowing a particular
* matrix. This function sets up appropriate parameters to the underlying
* PETSc object after it has been created.
*/
void initialize();
};
/**
* A class that implements the interface to use the PETSc SOR
* preconditioner.
*
* See the comment in the base class
* @ref PreconditionerBase
* for when this preconditioner may or may not work.
*
* @ingroup PETScWrappers
* @author Wolfgang Bangerth, Timo Heister, 2004, 2011
*/
class PreconditionSOR : public PreconditionerBase
{
public:
/**
* Standardized data struct to pipe additional flags to the
* preconditioner.
*/
struct AdditionalData
{
/**
* Constructor. By default, set the damping parameter to one.
*/
AdditionalData (const double omega = 1);
/**
* Relaxation parameter.
*/
double omega;
};
/**
* Empty Constructor. You need to call initialize() before using this
* object.
*/
PreconditionSOR ();
/**
* Constructor. Take the matrix which is used to form the preconditioner,
* and additional flags if there are any.
*/
PreconditionSOR (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());
protected:
/**
* Store a copy of the flags for this particular preconditioner.
*/
AdditionalData additional_data;
};
/**
* A class that implements the interface to use the PETSc SSOR
* preconditioner.
*
* See the comment in the base class
* @ref PreconditionerBase
* for when this preconditioner may or may not work.
*
* @ingroup PETScWrappers
* @author Wolfgang Bangerth, Timo Heister, 2004, 2011
*/
class PreconditionSSOR : public PreconditionerBase
{
public:
/**
* Standardized data struct to pipe additional flags to the
* preconditioner.
*/
struct AdditionalData
{
/**
* Constructor. By default, set the damping parameter to one.
*/
AdditionalData (const double omega = 1);
/**
* Relaxation parameter.
*/
double omega;
};
/**
* Empty Constructor. You need to call initialize() before using this
* object.
*/
PreconditionSSOR ();
/**
* Constructor. Take the matrix which is used to form the preconditioner,
* and additional flags if there are any.
*/
PreconditionSSOR (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());
protected:
/**
* Store a copy of the flags for this particular preconditioner.
*/
AdditionalData additional_data;
};
/**
* A class that implements the interface to use the PETSc Eisenstat
* preconditioner.
*
* See the comment in the base class
* @ref PreconditionerBase
* for when this preconditioner may or may not work.
*
* @ingroup PETScWrappers
* @author Wolfgang Bangerth, Timo Heister, 2004, 2011
*/
class PreconditionEisenstat : public PreconditionerBase
{
public:
/**
* Standardized data struct to pipe additional flags to the
* preconditioner.
*/
struct AdditionalData
{
/**
* Constructor. By default, set the damping parameter to one.
*/
AdditionalData (const double omega = 1);
/**
* Relaxation parameter.
*/
double omega;
};
/**
* Empty Constructor. You need to call initialize() before using this
* object.
*/
PreconditionEisenstat ();
/**
* Constructor. Take the matrix which is used to form the preconditioner,
* and additional flags if there are any.
*/
PreconditionEisenstat (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());
protected:
/**
* Store a copy of the flags for this particular preconditioner.
*/
AdditionalData additional_data;
};
/**
* A class that implements the interface to use the PETSc Incomplete
* Cholesky preconditioner.
*
* See the comment in the base class
* @ref PreconditionerBase
* for when this preconditioner may or may not work.
*
* @ingroup PETScWrappers
* @author Wolfgang Bangerth, Timo Heister, 2004, 2011
*/
class PreconditionICC : public PreconditionerBase
{
public:
/**
* Standardized data struct to pipe additional flags to the
* preconditioner.
*/
struct AdditionalData
{
/**
* Constructor. By default, set the fill-in parameter to zero.
*/
AdditionalData (const unsigned int levels = 0);
/**
* Fill-in parameter.
*/
unsigned int levels;
};
/**
* Empty Constructor. You need to call initialize() before using this
* object.
*/
PreconditionICC ();
/**
* Constructor. Take the matrix which is used to form the preconditioner,
* and additional flags if there are any.
*/
PreconditionICC (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());
protected:
/**
* Store a copy of the flags for this particular preconditioner.
*/
AdditionalData additional_data;
};
/**
* A class that implements the interface to use the PETSc ILU
* preconditioner.
*
* See the comment in the base class
* @ref PreconditionerBase
* for when this preconditioner may or may not work.
*
* @ingroup PETScWrappers
* @author Wolfgang Bangerth, Timo Heister, 2004, 2011
*/
class PreconditionILU : public PreconditionerBase
{
public:
/**
* Standardized data struct to pipe additional flags to the
* preconditioner.
*/
struct AdditionalData
{
/**
* Constructor. By default, set the fill-in parameter to zero.
*/
AdditionalData (const unsigned int levels = 0);
/**
* Fill-in parameter.
*/
unsigned int levels;
};
/**
* Empty Constructor. You need to call initialize() before using this
* object.
*/
PreconditionILU ();
/**
* Constructor. Take the matrix which is used to form the preconditioner,
* and additional flags if there are any.
*/
PreconditionILU (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());
protected:
/**
* Store a copy of the flags for this particular preconditioner.
*/
AdditionalData additional_data;
};
/**
* A class that implements the interface to use the PETSc LU preconditioner.
* The LU decomposition is only implemented for single processor machines.
* It should provide a convenient interface to another direct solver.
*
* See the comment in the base class
* @ref PreconditionerBase
* for when this preconditioner may or may not work.
*
* @ingroup PETScWrappers
* @author Oliver Kayser-Herold, 2004
*/
class PreconditionLU : public PreconditionerBase
{
public:
/**
* Standardized data struct to pipe additional flags to the
* preconditioner.
*/
struct AdditionalData
{
/**
* Constructor. (Default values taken from function PCCreate_LU of the
* PetSC lib.)
*/
AdditionalData (const double pivoting = 1.e-6,
const double zero_pivot = 1.e-12,
const double damping = 0.0);
/**
* Determines, when Pivoting is done during LU decomposition. 0.0
* indicates no pivoting, and 1.0 complete pivoting. Confer PetSC manual
* for more details.
*/
double pivoting;
/**
* Size at which smaller pivots are declared to be zero. Confer PetSC
* manual for more details.
*/
double zero_pivot;
/**
* This quantity is added to the diagonal of the matrix during
* factorisation.
*/
double damping;
};
/**
* Empty Constructor. You need to call initialize() before using this
* object.
*/
PreconditionLU ();
/**
* Constructor. Take the matrix which is used to form the preconditioner,
* and additional flags if there are any.
*/
PreconditionLU (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());
protected:
/**
* Store a copy of the flags for this particular preconditioner.
*/
AdditionalData additional_data;
};
/**
* A class that implements the interface to use the BoomerAMG algebraic
* multigrid preconditioner from the HYPRE suite. Note that PETSc has to be
* configured with HYPRE (e.g. with --download-hypre=1).
*
* The preconditioner does support parallel distributed computations. See
* step-40 for an example.
*
* @ingroup PETScWrappers
* @author Timo Heister, 2010
*/
class PreconditionBoomerAMG : public PreconditionerBase
{
public:
/**
* Standardized data struct to pipe additional flags to the
* preconditioner.
*/
struct AdditionalData
{
/**
* Constructor. Note that BoomerAMG offers a lot more options to set
* than what is exposed here.
*/
AdditionalData (
const bool symmetric_operator = false,
const double strong_threshold = 0.25,
const double max_row_sum = 0.9,
const unsigned int aggressive_coarsening_num_levels = 0,
const bool output_details = false
);
/**
* Set this flag to true if you have a symmetric system matrix and you
* want to use a solver which assumes a symmetric preconditioner like
* CG. The relaxation is done with SSOR/Jacobi when set to true and with
* SOR/Jacobi otherwise.
*/
bool symmetric_operator;
/**
* Threshold of when nodes are considered strongly connected. See
* HYPRE_BoomerAMGSetStrongThreshold(). Recommended values are 0.25 for
* 2d and 0.5 for 3d problems, but it is problem dependent.
*/
double strong_threshold;
/**
* If set to a value smaller than 1.0 then diagonally dominant parts of
* the matrix are treated as having no strongly connected nodes. If the
* row sum weighted by the diagonal entry is bigger than the given
* value, it is considered diagonally dominant. This feature is turned
* of by setting the value to 1.0. This is the default as some matrices
* can result in having only diagonally dominant entries and thus no
* multigrid levels are constructed. The default in BoomerAMG for this
* is 0.9. When you try this, check for a reasonable number of levels
* created.
*/
double max_row_sum;
/**
* Number of levels of aggressive coarsening. Increasing this value
* reduces the construction time and memory requirements but may
* decrease effectiveness.
*/
unsigned int aggressive_coarsening_num_levels;
/**
* Setting this flag to true produces debug output from HYPRE, when the
* preconditioner is constructed.
*/
bool output_details;
};
/**
* Empty Constructor. You need to call initialize() before using this
* object.
*/
PreconditionBoomerAMG ();
/**
* Constructor. Take the matrix which is used to form the preconditioner,
* and additional flags if there are any.
*/
PreconditionBoomerAMG (const MatrixBase &matrix,
const AdditionalData &additional_data = AdditionalData());
/**
* Same as above but without setting a matrix to form the preconditioner.
* Intended to be used with SLEPc objects.
*/
PreconditionBoomerAMG (const MPI_Comm communicator,
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());
protected:
/**
* Store a copy of the flags for this particular preconditioner.
*/
AdditionalData additional_data;
/**
* Initializes the preconditioner object without knowing a particular
* matrix. This function sets up appropriate parameters to the underlying
* PETSc object after it has been created.
*/
void initialize();
};
/**
* 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:
* <ul>
* <li> @p 0: nonsymmetric and/or indefinite problem, and nonsymmetric
* preconditioner
* <li> @p 1: SPD problem, and SPD (factored) preconditioner
* <li> @p 2: nonsymmetric, definite problem, and SPD (factored)
* preconditioner
* </ul>
* Default is <tt>symmetric = 1</tt>.
*/
unsigned int symmetric;
/**
* The sparsity pattern used for the approximate inverse is the pattern
* of a power <tt>B^m</tt> where <tt>B</tt> has been sparsified from the
* given matrix <tt>A</tt>, <tt>n_level</tt> is equal to <tt>m+1</tt>.
* Default value is <tt>n_levels = 1</tt>.
*/
unsigned int n_levels;
/**
* Sparsification is performed by dropping nonzeros which are smaller
* than <tt>thresh</tt> in magnitude. Lower values of <tt>thresh</tt>
* lead to more accurate, but also more expensive preconditioners.
* Default value is <tt>thresh = 0.1</tt>. Setting <tt>thresh < 0</tt> a
* threshold is selected automatically, such that <tt>-thresh</tt>
* represents the fraction of nonzero elements that are dropped. For
* example, if <tt>thresh = -0.9</tt>, then <tt>B</tt> will contain
* about ten percent of the nonzeros of the given matrix <tt>A</tt>.
*/
double threshold;
/**
* Filtering is a post-processing procedure, <tt>filter</tt> represents
* a fraction of nonzero elements that are dropped after creating the
* approximate inverse sparsity pattern. Default value is <tt>filter =
* 0.05</tt>. Setting <tt>filter < 0</tt> a value is selected
* automatically, such that <tt>-filter</tt> represents the fraction of
* nonzero elements that are dropped. For example, if <tt>thresh =
* -0.9</tt>, 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;
};
}
DEAL_II_NAMESPACE_CLOSE
#endif // DEAL_II_WITH_PETSC
/*---------------------------- petsc_precondition.h ---------------------------*/
#endif
/*---------------------------- petsc_precondition.h ---------------------------*/
|