#include <ml_MultiLevelPreconditioner.h>
Inheritance diagram for ML_Epetra::MultiLevelPreconditioner:
Public Member Functions | |
| int | DestroyPreconditioner () |
Destroys all structures allocated in ComputePreconditioner() if the preconditioner has been computed. | |
| const Epetra_BlockMap & | Map () const |
| Returns a reference to RowMatrix->Map(). | |
| int | NumGlobalRows () const |
| int | NumGlobalCols () const |
| int | NumMyRows () const |
| int | NumMyCols () const |
Constructors. | |
| MultiLevelPreconditioner (const Epetra_RowMatrix &RowMatrix, const bool ComputePrec) | |
| Constructs an MultiLevelPreconditioner with default values. | |
| MultiLevelPreconditioner (const Epetra_RowMatrix &RowMatrix, const Teuchos::ParameterList &List, const bool ComputePrec=true) | |
Constructs an MultiLevelPreconditioner. Retrives parameters from List. | |
| MultiLevelPreconditioner (ML_Operator *Operator, const Teuchos::ParameterList &List, const bool ComputePrec=true) | |
Constructs an MultiLevelPreconditioner from an ML_Operator. Retrives parameters from List. | |
| MultiLevelPreconditioner (const Epetra_RowMatrix &EdgeMatrix, const Epetra_RowMatrix &TMatrix, const Epetra_RowMatrix &NodeMatrix, const Teuchos::ParameterList &List, const bool ComputePrec=true) | |
Constructs an MultiLevelPreconditioner for Maxwell equations. Retrives parameters from List. | |
| MultiLevelPreconditioner (const Epetra_MsrMatrix &EdgeMatrix, ML_Operator *ML_TMatrix, AZ_MATRIX *AZ_NodeMatrix, int *proc_config, const Teuchos::ParameterList &List, const bool ComputePrec) | |
Constructs an MultiLevelPreconditioner for Maxwell equations. Retrives parameters from List. | |
Destructor. | |
| ~MultiLevelPreconditioner () | |
| Destroys the preconditioner. | |
Query functions | |
| 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 (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 cout on proc MyPID. | |
| Teuchos::ParameterList & | GetList () |
| Gets a reference to the internally stored parameters' list. | |
| Teuchos::ParameterList | GetOutputList () |
| void | PrintList (int MyPID) |
Prints on cout the values of the internally stored parameter list for processor MyPID. | |
| int | SetParameterList (const Teuchos::ParameterList &List) |
Copies List into the internally stored parameter list object. | |
Mathematical functions. | |
| int | Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) 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. | |
Atribute access functions | |
| int | ComputePreconditioner (const bool CheckFiltering=false) |
| Computes the multilevel hierarchy. | |
| 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 UseTranspose) |
| 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. | |
debugging and other utilities | |
| int | BreakForDebugger () |
| Stops the code, waiting for a debugger to attach. | |
| 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). | |
| int | AnalyzeSmoothersSparse (const int NumPreCycles=1, const int NumPostCycles=1) |
| Analyze the effect of each level's smoother on a random vector. | |
| int | AnalyzeSmoothersDense (const int NumPreCycles=1, const int NumPostCycles=1, const int MaxSize=1024) |
| Analyze the effect of each level's smoother on a random vector. | |
| int | AnalyzeMatrixCheap () |
| Cheap analysis of each level matrix. | |
| int | AnalyzeMatrixEigenvaluesSparse (char *MatVec, bool IsSymmetric=false) |
| Compute the lowest and largest magniture eigenvalues of fine-level matrix. | |
| int | AnalyzeMatrixEigenvaluesDense (char *MatVec, bool IsSymmetric=false) |
| Compute the lowest and largest magniture eigenvalues of fine-level matrix. | |
| int | AnalyzeCycle (const int NumCycles=1) |
| Analyze the effect of the ML cycle on a random 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 |
| 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) |
| int | VisualizeAggregates () |
| Visualizes the shape of the aggregates. | |
| int | VisualizeSmoothers (int NumPrecCycles=1, int NumPostCycles=1) |
| Visualizes the effect of smoothers on a random vector. | |
| int | VisualizeCycle (int NumCycles=1) |
| Visualizes the effect of the ML cycle on a random vector. | |
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-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:
|
||||||||||||||||||||||||
|
Constructs an MultiLevelPreconditioner for Maxwell equations. Retrives parameters from Constructs an MultiLevelPreconditioner for Maxwell equations. The constructor requires the edge matrix, the connectivity matrix T, the nodal matrix. |
|
||||||||||||||||||||||||||||
|
Constructs an MultiLevelPreconditioner for Maxwell equations. Retrives parameters from Constructs an MultiLevelPreconditioner for Maxwell equations. The constructor requires the edge matrix, the connectivity matrix T, the nodal matrix. |
|
|
Stops the code, waiting for a debugger to attach. BreakForDebugger() is useful when the user wants to attach to the running process(es). This is a very easy task for serial runs -- just run gdb. Parallel runs may result more problematic. In this case, one can proceed as follows:
|
|
|
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:
|
|
||||||||||||||||||||
|
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.
|
|
|
Prints unused parameters in the input ParameterList to cout on proc 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.
|
1.3.9.1