ML
Version of the Day
|
MultiLevelPreconditioner: a class to define black-box multilevel preconditioners using aggregation methods. More...
#include <ml_MultiLevelPreconditioner.h>
Public Member Functions | |
int | DestroyPreconditioner () |
Destroys all structures allocated in ComputePreconditioner() if the preconditioner has been computed. | |
const Epetra_RowMatrix & | RowMatrix () const |
Returns a reference to the internally stored RowMatrix. | |
const Epetra_BlockMap & | Map () const |
Returns a reference to RowMatrix->Map(). | |
int | NumGlobalRows () const |
Returns the global number of rows in the matrix. | |
int | NumGlobalCols () const |
Returns the global number of columns in the matrix. | |
int | NumMyRows () const |
Returns the local number of rows in the matrix. | |
int | NumMyCols () const |
Returns the local number of columns in the matrix. | |
int | PrintStencil2D (const int nx, const int ny, int NodeID=-1, const int EquationID=0) |
Prints the computational stencil for the specified row and equation (for 2D Cartesian grids only) More... | |
int | AnalyzeHierarchy (const bool AnalyzeMatrices, const int PreCycles, const int PostCycles, const int MLCycles) |
Cheap analysis of each level matrix. | |
int | AnalyzeSmoothers (const int NumPreCycles=1, const int NumPostCycles=1) |
Analyze the effect of each level's smoother on a random std::vector. | |
int | AnalyzeCoarse () |
Analyze the effect of the coarse solver on a random std::vector. | |
int | AnalyzeCycle (const int NumCycles=1) |
Analyze the effect of the ML cycle on a random std::vector. | |
int | TestSmoothers (Teuchos::ParameterList &InputList, const bool IsSymmetric=false) |
Test several smoothers on fine-level matrix. | |
int | TestSmoothers (const bool IsSymmetric=false) |
Test several smoothers on fine-level matrix using the current parameters. | |
const ML * | GetML (const int WhichML=-1) const |
Returns a pointer to the internally stored ml pointer. | |
const ML_Aggregate * | GetML_Aggregate () const |
Returns a pointer to the internally stored agg pointer. | |
int | Visualize (bool VizAggre, bool VizPreSmoother, bool VizPostSmoother, bool VizCycle, int NumApplPreSmoother, int NumApplPostSmoother, int NumCycleSmoother) |
Generic interface to visualization methods. | |
int | VisualizeAggregates () |
Visualizes the shape of the aggregates. | |
int | VisualizeSmoothers (int NumPrecCycles=1, int NumPostCycles=1) |
Visualizes the effect of smoothers on a random std::vector. | |
int | VisualizeCycle (int NumCycles=1) |
Visualizes the effect of the ML cycle on a random std::vector. | |
int | CreateLabel () |
void | ReportTime () |
void | Complexities (double &complexity, double &fineNnz) |
Return operator complexity and #nonzeros in fine grid matrix. | |
MultiLevelPreconditioner (const Epetra_RowMatrix &RowMatrix, const bool ComputePrec=true) | |
Constructs a MultiLevelPreconditioner with default values. | |
MultiLevelPreconditioner (const Epetra_RowMatrix &RowMatrix, const Teuchos::ParameterList &List, const bool ComputePrec=true) | |
Constructs a MultiLevelPreconditioner. Retrieves parameters from List . | |
MultiLevelPreconditioner (ML_Operator *Operator, const Teuchos::ParameterList &List, const bool ComputePrec=true) | |
Constructs a MultiLevelPreconditioner from an ML_Operator. Retrieves parameters from List . | |
MultiLevelPreconditioner (ML_Operator *Operator, const Teuchos::ParameterList &List, Epetra_RowMatrix **DiagOperators, Teuchos::ParameterList *DiagLists, int NBlocks=1, const bool ComputePrec=true) | |
Constructs a MultiLevelPreconditioner which is actually a composite AMG hierarchy using an array of ML_Operator's and an array of parameter lists. | |
MultiLevelPreconditioner (const Epetra_RowMatrix &EdgeMatrix, const Epetra_RowMatrix &GradMatrix, const Epetra_RowMatrix &NodeMatrix, const Teuchos::ParameterList &List, const bool ComputePrec=true, const bool UseNodeMatrixForSmoother=false) | |
MultiLevelPreconditioner constructor for Maxwell's equations. More... | |
MultiLevelPreconditioner (const Epetra_RowMatrix &CurlCurlMatrix, const Epetra_RowMatrix &MassMatrix, const Epetra_RowMatrix &TMatrix, const Epetra_RowMatrix &NodeMatrix, const Teuchos::ParameterList &List, const bool ComputePrec=true) | |
MultiLevelPreconditioner constructor for Maxwell's equations. More... | |
MultiLevelPreconditioner (Epetra_RowMatrix &RowMatrix, const Teuchos::ParameterList &List, const int &nNodes, const int &maxDofPerNode, bool *dofPresent, Epetra_Vector &Lhs, Epetra_Vector &Rhs, const bool rhsAndsolProvided, const bool ComputePrec=true) | |
Constructs a MultiLevelPreconditioner for multiphysics with variable dofs per node. More... | |
MultiLevelPreconditioner (Epetra_RowMatrix &RowMatrix, const Teuchos::ParameterList &List, const double distTol, const double tol, Epetra_Vector &Lhs, Epetra_Vector &Rhs, const bool rhsAndsolProvided, const bool ComputePrec=true) | |
MultiLevelPreconditioner (const Epetra_MsrMatrix &EdgeMatrix, ML_Operator *GradMatrix, AZ_MATRIX *NodeMatrix, int *proc_config, const Teuchos::ParameterList &List, const bool ComputePrec=true) | |
MultiLevelPreconditioner constructor for Maxwell's equations. More... | |
virtual | ~MultiLevelPreconditioner () |
Destroys the preconditioner. | |
const char * | Label () const |
Prints label associated to this object. | |
void | PrintUnused () const |
Prints unused parameters in the input ParameterList on standard output. | |
void | PrintUnused (std::ostream &os) const |
Prints unused parameters in the input ParameterList on the specified stream. | |
void | PrintUnused (const int MyPID) const |
Prints unused parameters in the input ParameterList to std::cout on proc MyPID . More... | |
Teuchos::ParameterList & | GetList () |
Gets a reference to the internally stored parameters' list. | |
Teuchos::ParameterList | GetOutputList () |
void | PrintList () |
Prints on std::cout the values of the internally stored parameter list. | |
int | SetParameterList (const Teuchos::ParameterList &List) |
Copies List into the internally stored parameter list object. | |
int | Apply (const Epetra_MultiVector &, Epetra_MultiVector &) const |
Apply the inverse of the preconditioner to an Epetra_MultiVector (NOT AVAILABLE) | |
int | ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
Apply the preconditioner to an Epetra_MultiVector X, puts the result in Y. | |
int | ComputePreconditioner (const bool CheckFiltering=false) |
Computes the multilevel hierarchy. More... | |
int | ReComputePreconditioner (bool keepFineLevelSmoother=false) |
Recompute the preconditioner (not implemented for Maxwell). More... | |
void | Print (int level=-2) |
Print the individual operators in the multigrid hierarchy. More... | |
int | ComputeAdaptivePreconditioner (int TentativeNullSpaceSize, double *TentativeNullSpace) |
int | IsPreconditionerComputed () const |
Queries whether multilevel hierarchy has been computed or not. | |
int | SetOwnership (bool ownership) |
Sets ownership. | |
int | SetUseTranspose (bool) |
Sets use transpose (not implemented). | |
double | NormInf () const |
Returns the infinity norm (not implemented). | |
bool | UseTranspose () const |
Returns the current UseTranspose setting. | |
bool | HasNormInf () const |
Returns true if the this object can provide an approximate Inf-norm, false otherwise. | |
const Epetra_Comm & | Comm () const |
Returns a pointer to the Epetra_Comm communicator associated with this operator. | |
const Epetra_Map & | OperatorDomainMap () const |
Returns the Epetra_Map object associated with the domain of this operator. | |
const Epetra_Map & | OperatorRangeMap () const |
Returns the Epetra_Map object associated with the range of this operator. | |
MultiLevelPreconditioner: a class to define black-box multilevel preconditioners using aggregation methods.
Class ML_Epetra::MultiLevelPreconditioner defined black-box algebraic multilevel preconditioners of matrices defined as Epetra_RowMatrix derived objects. The resulting preconditioner can be used in AztecOO, and in any other solver that accepts Epetra_Operator derived objects, and apply the action of the given Epetra_Operator using ApplyInverse().
Please refer to the user's guide for a detailed introduction to this class, examples, and description of input parameters.
This file requires ML to be configured with the following options:
–enable-epetra
–enable-teuchos
The following option is suggested:
–enable-amesos
–enable-ifpack
Some part of this class needs the following options:
–enable-aztecoo
–enable-anasazi
It is important to note that ML is more restrictive than Epetra for the definition of maps. It is required that RowMatrixRowMap() is equal to OperatorRangeMap(). This is because ML needs to perform matrix-std::vector product, as well as getrow() functions, on the same data distribution.
Also, for square matrices, OperatorDomainMap() must be as OperatorRangeMap().
Several examples are provided in the examples
subdirectories:
Defaults parameters can be specified using function SetDefaults().
MultiLevelPreconditioner::MultiLevelPreconditioner | ( | const Epetra_RowMatrix & | EdgeMatrix, |
const Epetra_RowMatrix & | TMatrix, | ||
const Epetra_RowMatrix & | NodeMatrix, | ||
const Teuchos::ParameterList & | List, | ||
const bool | ComputePrec = true , |
||
const bool | UseNodeMatrixForSmoother = false |
||
) |
MultiLevelPreconditioner constructor for Maxwell's equations.
Takes the stiffness and mass terms of the matrix combined.
EdgeMatrix | - (In) Linear matrix to be solved. |
GradMatrix | - (In) Node-to-edge connectivity matrix, a.k.a, topological gradient |
NodeMatrix | - (In) Auxiliary nodal finite element matrix |
List | - (In) Teuchos parameter list containing solver options. |
ComputePrec | - (In) Optional argument that specifies whether to create preconditioner immediately. Default is true. |
UseNodeMatrixForSmoother | - (In) Use the nodal matrix for the nodal portion of the Hipmair smoother (if used). |
The constructor for the Maxwell equations requires three Epetra_RowMatrices. Two conditions are required on their maps:
References ComputePreconditioner(), Epetra_Operator::OperatorDomainMap(), Epetra_Operator::OperatorRangeMap(), and Epetra_BlockMap::SameAs().
MultiLevelPreconditioner::MultiLevelPreconditioner | ( | const Epetra_RowMatrix & | CurlCurlMatrix, |
const Epetra_RowMatrix & | MassMatrix, | ||
const Epetra_RowMatrix & | TMatrix, | ||
const Epetra_RowMatrix & | NodeMatrix, | ||
const Teuchos::ParameterList & | List, | ||
const bool | ComputePrec = true |
||
) |
MultiLevelPreconditioner constructor for Maxwell's equations.
Takes the stiffness and mass terms of the matrix separately.
CurlCurlMatrix | - (In) The curl-curl (stiffness) term of the matrix to be solved. |
MassMatrix | - (In) The mass term of the matrix to be solved. |
GradMatrix | - (In) Node-to-edge connectivity matrix, a.k.a, topological gradient |
NodeMatrix | - (In) Auxiliary nodal finite element matrix |
List | - (In) Teuchos parameter list containing solver options. |
ComputePrec | - (In) Optional argument that specifies whether to create preconditioner immediately. Default is true. |
Constructor for the Maxwell equations. This version takes the stiffness and edge mass matrices separately. Two conditions are required on their maps:
References ComputePreconditioner(), Epetra_Operator::OperatorDomainMap(), Epetra_Operator::OperatorRangeMap(), and Epetra_BlockMap::SameAs().
MultiLevelPreconditioner::MultiLevelPreconditioner | ( | Epetra_RowMatrix & | RowMatrix, |
const Teuchos::ParameterList & | List, | ||
const int & | nNodes, | ||
const int & | maxDofPerNode, | ||
bool * | NotMLVecDofPresent, | ||
Epetra_Vector & | Lhs, | ||
Epetra_Vector & | Rhs, | ||
const bool | rhsAndsolProvided, | ||
const bool | ComputePrec = true |
||
) |
Constructs a MultiLevelPreconditioner for multiphysics with variable dofs per node.
Constructor for multiphysics problems with variable dofs per node. This version uses a discrete laplacian and padding on coarse grids to make a hierarchy. It also allows for the removal of column nonzeros associated with Dirichlet points. To use this option the rhs and initial guess must be provided.
References Epetra_Operator::Comm(), ML_Struct::comm, ComputePreconditioner(), Copy, CreateLabel(), Epetra_CrsMatrix::ExtractCrsDataPointers(), Epetra_Vector::ExtractView(), Teuchos::ParameterList::get(), Epetra_Comm::MaxAll(), ML_Struct::ML_num_actual_levels, MLShove(), MLunamalgP(), Epetra_Comm::MyPID(), Epetra_CrsMatrix::NumMyCols(), Epetra_CrsMatrix::NumMyRows(), Epetra_Comm::NumProc(), Epetra_Operator::OperatorDomainMap(), Epetra_Operator::OperatorRangeMap(), Epetra_CrsMatrix::RowMap(), RowMatrix(), Epetra_Comm::SumAll(), and TEUCHOS_TEST_FOR_EXCEPT_MSG.
MultiLevelPreconditioner::MultiLevelPreconditioner | ( | Epetra_RowMatrix & | RowMatrix, |
const Teuchos::ParameterList & | List, | ||
const double | distTol, | ||
const double | tol, | ||
Epetra_Vector & | Lhs, | ||
Epetra_Vector & | Rhs, | ||
const bool | rhsAndsolProvided, | ||
const bool | ComputePrec = true |
||
) |
Constructor for scalar PDE problems based on applying AMG to the distance Laplacian operator when constructing grid transfers. The main unique feature is that there may be some dofs that correspond to the same node location. These shared dofs fall into two categories. If these dofs are strongly connected to each other (as determined by tol), they are explicitly elminated from the Laplacian (merged into a supernode). Once a P is obtained, this P is then expanded to account for shared nodes by simply duplicating the supernodes row of P for each of the individual vertices that contribute to the supernode. If share dofs are weakly connected (or not connected at all), nothing special is done (other than the ususal ignoring of weak connections). One last trick is employed, connections between supernodes and non-supernodes (i.e., regular nodes) are always assumed to be weak. Shared nodes are often used to capture interfaces or other features. By breaking these connections, the AMG can better maintain these features throughout the hierarchy. Finally, the constructor also allows for * the removal of column nonzeros associated with Dirichlet points. To use this option the rhs and initial guess must be provided. Modification of the matrix, rhs, and initial guess must be allowable to use this option.
References Epetra_CrsMatrix::ColMap(), Epetra_Operator::Comm(), ML_Struct::comm, ComputePreconditioner(), Copy, CreateLabel(), Epetra_CrsMatrix::ExtractCrsDataPointers(), Epetra_Vector::ExtractView(), Teuchos::ParameterList::get(), ML_Struct::ML_num_actual_levels, Epetra_Comm::MyPID(), ML_Operator_Struct::N_nonzeros, Epetra_CrsMatrix::NumMyCols(), Epetra_BlockMap::NumMyElements(), Epetra_CrsMatrix::NumMyRows(), Epetra_Comm::NumProc(), Epetra_Operator::OperatorDomainMap(), Epetra_Operator::OperatorRangeMap(), Epetra_CrsMatrix::RowMap(), RowMatrix(), and TEUCHOS_TEST_FOR_EXCEPT_MSG.
MultiLevelPreconditioner::MultiLevelPreconditioner | ( | const Epetra_MsrMatrix & | EdgeMatrix, |
ML_Operator * | ML_TMatrix, | ||
AZ_MATRIX * | AZ_NodeMatrix, | ||
int * | proc_config, | ||
const Teuchos::ParameterList & | List, | ||
const bool | ComputePrec = true |
||
) |
MultiLevelPreconditioner constructor for Maxwell's equations.
Takes the stiffness and mass terms of the matrix combined. The edge matrix is of type Epetra_Msr, a light-weight wrapper for old-style Aztec MSR matrices. This is intended as transition code for Aztec users.
EdgeMatrix | - (In) Linear matrix to be solved. |
GradMatrix | - (In) Node-to-edge connectivity matrix, a.k.a, topological gradient |
NodeMatrix | - (In) Auxiliary nodal finite element matrix |
proc_config | - (In) Aztec array specifying processor layout. |
List | - (In) Teuchos parameter list containing solver options. |
ComputePrec | - (In) Optional argument that specifies whether to create preconditioner immediately. Default is true. |
Another constructor for Maxwell equations that takes an Epetra_MsrMatrix, an ML_Operator, and an AZ_MATRIX type. The Epetra_MsrMatrix type is defined in aztecoo. Two conditions are required on their maps:
References Epetra_MsrMatrix::Comm(), Epetra_MpiComm::Comm(), ComputePreconditioner(), ML_Operator2EpetraCrsMatrix(), Epetra_MsrMatrix::OperatorDomainMap(), Epetra_CrsMatrix::OperatorDomainMap(), Epetra_CrsMatrix::OperatorRangeMap(), and Epetra_BlockMap::SameAs().
int MultiLevelPreconditioner::ComputePreconditioner | ( | const bool | CheckFiltering = false | ) |
Computes the multilevel hierarchy.
Computes the multilevel hierarchy. This function retrives the user's defines parameters (as specified in the input ParameterList), or takes default values otherwise, and creates the ML objects for aggregation and hierarchy. Allocated data can be freed used DestroyPreconditioner(), or by the destructor,
In a Newton-type procedure, several linear systems have to be solved, Often, these systems are not too different. In this case, it might be convenient to keep the already computed preconditioner (with hierarchy, coarse solver, smoothers), and use it to precondition the next linear system. ML offers a way to determine whether the already available preconditioner is "good enough" for the next linear system. The user should proceed as follows:
"reuse: enable"
== true
ComputePreconditioner(true)
It is supposed that the pointer to the Epetra_RowMatrix remains constant. Currently, it is not possible to modify this pointer (other than creating a new preconditioner) References ML_Epetra::Apply_BCsToGradient(), Epetra_MpiComm::Comm(), Epetra_Time::ElapsedTime(), Teuchos::ParameterList::get(), ML_BreakForDebugger(), ML_Epetra_comm_wrapper(), ML_Epetra_getrow(), ML_Epetra_matvec(), ML_Struct::ML_num_actual_levels, ML_Operator_WrapEpetraMatrix(), ML_Epetra::ModifyEpetraMatrixColMap(), ML_Epetra::ReadXML(), Epetra_Time::ResetStartTime(), Teuchos::ParameterList::sublist(), and TEUCHOS_TEST_FOR_EXCEPT_MSG.
Referenced by MultiLevelPreconditioner().
int MultiLevelPreconditioner::CreateLabel | ( | ) |
Creates label for this object (printed out by AztecOO). This does not allocate/reallocate any memory.
Referenced by MultiLevelPreconditioner().
void MultiLevelPreconditioner::Print | ( | int | level = -2 | ) |
Print the individual operators in the multigrid hierarchy.
Print the individual operators in the multigrid hierarchy.
level | (In) By default, this method prints the entire multigrid hierarchy. If you desire only the operators associated with a particular level, pass in the level number. The fine level is 0, coarser levels are positive integers. |
int MultiLevelPreconditioner::PrintStencil2D | ( | const int | nx, |
const int | ny, | ||
int | NodeID = -1 , |
||
const int | EquationID = 0 |
||
) |
Prints the computational stencil for the specified row and equation (for 2D Cartesian grids only)
For problems defined on 2D Cartesian grids (with node numbering increasing along the x-axis), this function prints out the stencil in an intelligible form.
nx | (In) : number of nodes along the X-axis |
ny | (In) : number of nodes along the Y-axis |
NodeID | (In) : (local) ID of node that will be used to print the stencil. If set to -1, the code will automatically chose an internal node. Default: -1. |
EquationID | (In) : ID of the equation that will be used to print the stencil (default = 0) |
void MultiLevelPreconditioner::PrintUnused | ( | const int | MyPID | ) | const |
Prints unused parameters in the input ParameterList to std::cout on proc MyPID
.
Mispelled parameters are simply ignored. Therefore, it is often the best choice to print out the parameters that have not been used in the construction phase.
MyPID | (In) : ID of process that should print the unused parameters. |
int MultiLevelPreconditioner::ReComputePreconditioner | ( | bool | keepFineLevelSmoother = false | ) |
Recompute the preconditioner (not implemented for Maxwell).
keepFineLevelSmoother | (In) : If true, the fine level smoother is not recomputed. This is useful if the smoother is expensive to create, e.g., an incomplete factorization, and the fine level matrix has not changed. |
References Epetra_Time::ElapsedTime(), and Epetra_Time::ResetStartTime().