ML
Version of the Day
|
ML_Epetra: default namespace for all Epetra interfaces. More...
Classes | |
class | CrsGraphWrapper |
ML_Epetra::CrsGraphWrapper: a class to wrap an Epetra_CrsGraph as Epetra_RowMatrix. More... | |
class | Ifpack_ML |
Wraps an ML preconditioner as an Ifpack_Preconditioner. More... | |
class | MultiLevelOperator |
MultiLevelOperator: An implementation of the Epetra_Operator class. More... | |
class | MultiLevelPreconditioner |
MultiLevelPreconditioner: a class to define black-box multilevel preconditioners using aggregation methods. More... | |
class | RowMatrix |
Basic wrapper from ML_Operator to Epetra_RowMatrix. More... | |
class | MatrixFreePreconditioner |
MatrixFreePreconditioner: a class to define preconditioners for Epetra_Operator's. More... | |
class | EdgeMatrixFreePreconditioner |
class | ML_RMP |
class | ML_RefMaxwell_11_Operator |
class | CoordPack |
Enumerations | |
enum | AMGType { ML_SA_FAMILY, ML_MAXWELL, ML_COMPOSITE, ML_CLASSICAL_FAMILY } |
Enumerated type indicating the type of AMG solver to be used. | |
Functions | |
int | ML_Epetra_PtAP (const Epetra_CrsMatrix &A, const Epetra_CrsMatrix &P, Epetra_CrsMatrix *&Result, bool verbose=false) |
Does an P^TAP for Epetra_CrsMatrices using ML's kernels. | |
int | ML_Epetra_RAP (const Epetra_CrsMatrix &A, const Epetra_CrsMatrix &P, const Epetra_CrsMatrix &R, Epetra_CrsMatrix *&Result, bool verbose) |
Does an RAP for Epetra_CrsMatrices using ML's kernels. | |
int | Epetra_PtAP (const Epetra_CrsMatrix &A, const Epetra_CrsMatrix &P, Epetra_CrsMatrix *&Result, bool keep_zero_rows=false, bool verbose=false) |
Does an P^TAP for Epetra_CrsMatrices using EpetraExt's kernels. | |
int * | FindLocalDiricheltRowsFromOnesAndZeros (const Epetra_CrsMatrix &Matrix, int &numBCRows) |
Finds the Dirichlet rows in a square matrix that got the one-and-zeros. More... | |
void | Apply_BCsToMatrixColumns (const int *dirichletRows, int numBCRows, const Epetra_CrsMatrix &Matrix) |
Applies Dirichlet conditions to columns that rows already have. More... | |
void | Apply_BCsToMatrixRows (const int *dirichletRows, int numBCRows, const Epetra_CrsMatrix &Matrix) |
Applies Dirichlet conditions to matrix rows. | |
void | Apply_BCsToMatrixRowsNodal (const int *dirichletRows, int numBCRows, int dofsPerNode, const Epetra_CrsMatrix &Matrix) |
Applies Dirichlet conditions to matrix rows (marking all dofs on a node as Dirichlet if any one dof is) | |
void | Apply_BCsToMatrixColumns (const Epetra_RowMatrix &iBoundaryMatrix, const Epetra_RowMatrix &iMatrix) |
Applies Dirichlet conditions to columns that rows already have. More... | |
void | Apply_BCsToMatrixColumns (const Epetra_IntVector &dirichletColumns, const Epetra_CrsMatrix &Matrix) |
void | Apply_BCsToGradient (const Epetra_RowMatrix &EdgeMatrix, const Epetra_RowMatrix &T) |
Applies boundary conditions to gradient matrix. (Maxwell's equations) | |
void | Apply_OAZToMatrix (int *dirichletRows, int numBCRows, const Epetra_CrsMatrix &Matrix) |
Does Row/Column OAZ to a matrix. | |
Epetra_IntVector * | LocalRowstoColumns (int *Rows, int numRows, const Epetra_CrsMatrix &Matrix) |
Returns the local column numbers of the local rows passed in. | |
Epetra_IntVector * | FindLocalDirichletColumnsFromRows (const int *dirichletRows, int numBCRows, const Epetra_CrsMatrix &Matrix) |
Finds Dirichlet the local Dirichlet columns, given the local Dirichlet rows. | |
Epetra_IntVector * | FindLocalDirichletDomainsFromRows (const int *dirichletRows, int numBCRows, const Epetra_CrsMatrix &Matrix) |
Finds Dirichlet the local Dirichlet domains, given the local Dirichlet rows. | |
void | Remove_Zeroed_Rows (const Epetra_CrsMatrix &Matrix, double tol=0.0) |
Drops a 1 on the diagonal of zero'd our rows. | |
Epetra_RowMatrix * | ModifyEpetraMatrixColMap (const Epetra_RowMatrix &A, EpetraExt::CrsMatrix_SolverMap &transform, const char *matrixName=0, bool verbose=false) |
Transforms Epetra matrix column map (if necessary) to be compatible with. More... | |
int | UpdateList (Teuchos::ParameterList &source, Teuchos::ParameterList &dest, bool OverWrite) |
int | SetDefaults (std::string ProblemType, Teuchos::ParameterList &List, int *options=0, double *params=0, const bool OverWrite=true) |
Sets default parameters for aggregation-based preconditioners. More... | |
int | SetDefaultsDD (Teuchos::ParameterList &List, Teuchos::RCP< std::vector< int > > &options, Teuchos::RCP< std::vector< double > > ¶ms, bool Overwrite=true) |
Sets default parameters for aggregation-based 2-level domain decomposition preconditioners. More... | |
int | SetDefaultsDD_LU (Teuchos::ParameterList &List, Teuchos::RCP< std::vector< int > > &options, Teuchos::RCP< std::vector< double > > ¶ms, bool Overwrite=true) |
Sets default parameters for aggregation-based 2-level domain decomposition preconditioners, using LU on each subdomain. More... | |
int | SetDefaultsDD_3Levels (Teuchos::ParameterList &List, Teuchos::RCP< std::vector< int > > &options, Teuchos::RCP< std::vector< double > > ¶ms, bool Overwrite=true) |
Sets default parameters for aggregation-based 3-level domain decomposition preconditioners. More... | |
int | SetDefaultsDD_3Levels_LU (Teuchos::ParameterList &List, Teuchos::RCP< std::vector< int > > &options, Teuchos::RCP< std::vector< double > > ¶ms, bool Overwrite=true) |
Sets default parameters for aggregation-based 3-level domain decomposition preconditioners with LU. More... | |
int | SetDefaultsMaxwell (Teuchos::ParameterList &List, Teuchos::RCP< std::vector< int > > &options, Teuchos::RCP< std::vector< double > > ¶ms, bool Overwrite=true) |
Sets default parameters for the eddy current equations equations. More... | |
int | SetDefaultsSA (Teuchos::ParameterList &List, Teuchos::RCP< std::vector< int > > &options, Teuchos::RCP< std::vector< double > > ¶ms, bool Overwrite=true) |
Sets default parameters for classical smoothed aggregation. More... | |
int | SetDefaultsNSSA (Teuchos::ParameterList &List, Teuchos::RCP< std::vector< int > > &options, Teuchos::RCP< std::vector< double > > ¶ms, bool Overwrite=true) |
Sets defaults for energy minimization preconditioning for nonsymmetric problems. More... | |
int | SetDefaultsClassicalAMG (Teuchos::ParameterList &List, Teuchos::RCP< std::vector< int > > &options, Teuchos::RCP< std::vector< double > > ¶ms, bool Overwrite=true) |
Sets defaults for classical amg. | |
int | ReadXML (const std::string &FileName, Teuchos::ParameterList &List, const Epetra_Comm &Comm) |
Reads in parameter list options from file. | |
int | SetDefaultsRefMaxwell (Teuchos::ParameterList &inList, bool OverWrite=true) |
Sets default parameters for aggregation-based 2-level domain decomposition preconditioners. | |
void | IVOUT (const Epetra_IntVector &A, const char *of) |
int | CSR_getrow_ones (ML_Operator *data, int N_requested_rows, int requested_rows[], int allocated_space, int columns[], double values[], int row_lengths[]) |
int | RefMaxwell_SetupCoordinates (ML_Operator *A, const double *in_coordx, const double *in_coordy, const double *in_coordz, const double *in_material, double *&coordx, double *&coordy, double *&coordz, double *&material) |
int | RefMaxwell_Aggregate_Nodes (const Epetra_CrsMatrix &A, Teuchos::ParameterList &List, ML_Comm *ml_comm, std::string PrintMsg, ML_Aggregate_Struct *&MLAggr, ML_Operator *&P, int &NumAggregates, CoordPack &pack) |
Epetra_MultiVector * | Build_Edge_Nullspace (const Epetra_CrsMatrix &D0Clean_Matrix, const Teuchos::ArrayRCP< int > BCedges, Teuchos::ParameterList &List, bool verbose) |
Epetra_CrsMatrix * | Build_Pi_Matrix (const Epetra_CrsMatrix &D0Clean_Matrix, const Teuchos::ArrayRCP< int > BCedges, Teuchos::ParameterList &List, bool verbose) |
ML_Epetra: default namespace for all Epetra interfaces.
void ML_Epetra::Apply_BCsToMatrixColumns | ( | const int * | dirichletRows, |
int | numBCRows, | ||
const Epetra_CrsMatrix & | Matrix | ||
) |
Applies Dirichlet conditions to columns that rows already have.
Apply the Dirichlet BC's specified by the dirichletRows to remove all columns (across all processors) that have entries zero'd by the BC's. This will symmetrize a square matrix, though the routine can be run on non-square matrices.
void ML_Epetra::Apply_BCsToMatrixColumns | ( | const Epetra_RowMatrix & | iBoundaryMatrix, |
const Epetra_RowMatrix & | iMatrix | ||
) |
Applies Dirichlet conditions to columns that rows already have.
Apply the Dirichlet BC's specified by BoundaryMatrix to remove all columns (across all processors) that have entries zero'd by the BC's. This will symmetrize a square matrix, though the routine can be run on non-square matrices.
int* ML_Epetra::FindLocalDiricheltRowsFromOnesAndZeros | ( | const Epetra_CrsMatrix & | Matrix, |
int & | numBCRows | ||
) |
Finds the Dirichlet rows in a square matrix that got the one-and-zeros.
Returns the local Dirichlet rows for a square matrix that go the ones-and-zeros treatment for BCs.
Epetra_RowMatrix* ML_Epetra::ModifyEpetraMatrixColMap | ( | const Epetra_RowMatrix & | A, |
EpetraExt::CrsMatrix_SolverMap & | transform, | ||
const char * | matrixName = 0 , |
||
bool | verbose = false |
||
) |
Transforms Epetra matrix column map (if necessary) to be compatible with.
how ML handles column indices. Any matrix that cannot be dynamically cast to an Epetra_CrsMatrix will not be changed.
A | - (In) Matrix that is to be transformed. |
transform | - (In) EpetraExt widget that does the transformation. |
matrixName | - (In) Optional label for the incoming matrix. |
Referenced by ML_Epetra::MultiLevelPreconditioner::ComputePreconditioner().
int ML_Epetra::SetDefaults | ( | std::string | ProblemType, |
Teuchos::ParameterList & | List, | ||
int * | options = 0 , |
||
double * | params = 0 , |
||
const bool | OverWrite = true |
||
) |
Sets default parameters for aggregation-based preconditioners.
This function, defined in the namespace ML_Epetra, can be used to set default values in a user's defined Teuchos::ParameterList.
ProblemType | (In) : a std::string, whose possible values are:
|
List | (Out) : list which will populated by the default parameters |
options | (In/Out) : integer array, of size AZ_OPTIONS_SIZE , that will be populated with suitable values. A pointer to options will be stick into the parameters list. Note that this array is still required to apply the preconditioner! Do not delete options, nor let it go out of scope. The default value is 0, meaning that SetDefaults() will allocate the array. |
params | (In/Out) : double array, of size AZ_PARAMS_SIZE . See comments for options . |
OverWrite | (In) : boolean. If false, any pre-existing values in the parameter list will be preserved. Default value is true, i.e., any pre-existing values may be overwritten. |
References Teuchos::rcp(), SetDefaultsClassicalAMG(), SetDefaultsDD(), SetDefaultsDD_3Levels(), SetDefaultsDD_3Levels_LU(), SetDefaultsDD_LU(), SetDefaultsMaxwell(), SetDefaultsNSSA(), and SetDefaultsSA().
Referenced by ML_Epetra::MultiLevelPreconditioner::MultiLevelPreconditioner(), ReadXML(), ML_Epetra::Ifpack_ML::SetParameters(), and ML_Epetra::MultiLevelPreconditioner::TestSmoothers().
int ML_Epetra::SetDefaultsDD | ( | Teuchos::ParameterList & | List, |
Teuchos::RCP< std::vector< int > > & | options, | ||
Teuchos::RCP< std::vector< double > > & | params, | ||
bool | OverWrite = true |
||
) |
Sets default parameters for aggregation-based 2-level domain decomposition preconditioners.
Sets the following default values for "DD".
Note: This method should not be called directly. Instead, ML_Epetra::SetDefaults("DD",...) should be used.
"default values"
=
"DD"
"max levels"
=
2
"prec type"
=
"MGV"
"increasing or decreasing"
=
"increasing"
"aggregation: type"
=
"METIS"
"aggregation: local aggregates"
=
1
"aggregation: damping factor"
=
1.333
"eigen-analysis: type"
=
"power-method"
"eigen-analysis: iterations"
=
20
"smoother: sweeps"
=
1
"smoother: pre or post"
=
"both"
"smoother: type"
=
Aztec".
- \c "smoother:
Aztec options" \c = \c options
- \c "smoother: Aztec params" \c = \c params
- \c "smoother: Aztec as solver" =
false
"coarse: type"
=
"Amesos-KLU"
"coarse: max size"
=
128
References Teuchos::ParameterList::set(), and Teuchos::ParameterList::setName().
Referenced by SetDefaults().
int ML_Epetra::SetDefaultsDD_3Levels | ( | Teuchos::ParameterList & | List, |
Teuchos::RCP< std::vector< int > > & | options, | ||
Teuchos::RCP< std::vector< double > > & | params, | ||
bool | OverWrite = true |
||
) |
Sets default parameters for aggregation-based 3-level domain decomposition preconditioners.
Sets the following default values for "DD-ML".
This method should not be called directly. Instead, ML_Epetra::SetDefaults("DD-ML",...) should be used.
"default values"
=
"DD-ML"
"max levels"
=
3
"prec type"
=
"MGV"
"increasing or decreasing"
=
"increasing"
"aggregation: type"
=
"METIS"
"aggregation: nodes per aggregate"
=
512
"aggregation: next-level aggregates per process"
=
128
"aggregation: damping factor"
=
1.333
"eigen-analysis: type"
=
"power-method"
"eigen-analysis: iterations"
=
20
"smoother: sweeps"
=
1
"smoother: pre or post"
=
"both"
"smoother: type"
=
"Aztec"
"smoother: Aztec options"
=
options
"smoother: Aztec params"
=
params
"smoother: Aztec as solver"
=
false
"coarse: type"
=
"Amesos-KLU"
"coarse: max size"
=
128
References Teuchos::ParameterList::set(), and Teuchos::ParameterList::setName().
Referenced by SetDefaults().
int ML_Epetra::SetDefaultsDD_3Levels_LU | ( | Teuchos::ParameterList & | List, |
Teuchos::RCP< std::vector< int > > & | options, | ||
Teuchos::RCP< std::vector< double > > & | params, | ||
bool | OverWrite = true |
||
) |
Sets default parameters for aggregation-based 3-level domain decomposition preconditioners with LU.
Same as SetDefaultsDD_3Levels but with LU factorizations subdomains.
Note: This method should not be called directly. Instead, ML_Epetra::SetDefaults("DD-ML-LU",...) should be used.
"default values"
=
"DD-ML-LU"
"max levels"
=
3
"prec type"
=
"MGV"
"increasing or decreasing"
=
"increasing"
"aggregation: type"
=
"METIS"
"aggregation: nodes per aggregate"
=
512
"aggregation: next-level aggregates per process"
=
128
"aggregation: damping factor"
=
1.333
"eigen-analysis: type"
=
"power-method"
"eigen-analysis: iterations"
=
20
"smoother: sweeps"
=
1
"smoother: pre or post"
=
"both"
"smoother: type"
=
"Aztec"
"smoother: Aztec options"
=
options
"smoother: Aztec params"
=
params
"smoother: Aztec as solver"
=
false
"coarse: type"
=
"Amesos-KLU"
"coarse: max size"
=
128
References Teuchos::ParameterList::set(), and Teuchos::ParameterList::setName().
Referenced by SetDefaults().
int ML_Epetra::SetDefaultsDD_LU | ( | Teuchos::ParameterList & | List, |
Teuchos::RCP< std::vector< int > > & | options, | ||
Teuchos::RCP< std::vector< double > > & | params, | ||
bool | OverWrite = true |
||
) |
Sets default parameters for aggregation-based 2-level domain decomposition preconditioners, using LU on each subdomain.
Same as SetDefaultsDD(), but used exact LU decompositions on subdomains.
Note: This method should not be called directly. Instead, ML_Epetra::SetDefaults("DD-LU",...) should be used.
=
"DD-LU"
"max levels"
=
2
"prec type"
=
"MGV"
"increasing or decreasing"
=
"increasing"
"aggregation: type"
=
"METIS"
"aggregation: local aggregates"
=
1
"aggregation: damping factor"
=
1.333
"eigen-analysis: type"
=
"power-method"
"eigen-analysis: iterations"
=
20
"smoother: sweeps"
=
1
"smoother: pre or post"
=
"both"
"smoother: type"
=
"Aztec"
"smoother: Aztec options"
=
options
"smoother: Aztec params"
=
params
"smoother: Aztec as solver"
=
false
"coarse: type"
=
"Amesos-KLU"
"coarse: max size"
=
128
References Teuchos::ParameterList::set(), and Teuchos::ParameterList::setName().
Referenced by SetDefaults().
int ML_Epetra::SetDefaultsMaxwell | ( | Teuchos::ParameterList & | List, |
Teuchos::RCP< std::vector< int > > & | options, | ||
Teuchos::RCP< std::vector< double > > & | params, | ||
bool | OverWrite = true |
||
) |
Sets default parameters for the eddy current equations equations.
Set values for Maxwell:
Note: This method should not be called directly. Instead, ML_Epetra::SetDefaults("maxwell",...) should be used.
"default values"
=
"maxwell"
"max levels"
=
10
"prec type"
=
"MGV"
"increasing or decreasing"
=
"decreasing"
"aggregation: type"
=
"Uncoupled-MIS"
"aggregation: damping factor"
=
1.333
"eigen-analysis: type"
=
"cg"
"eigen-analysis: iterations"
=
10
"aggregation: edge prolongator drop threshold"
=
0.0
"smoother: sweeps"
=
1
"smoother: damping factor"
=
1.0
"smoother: pre or post"
=
"both"
"smoother: type"
=
"Hiptmair"
"smoother: Hiptmair efficient symmetric"
=
true
"subsmoother: type"
=
"Chebyshev"
"subsmoother: Chebyshev alpha"
=
20.0
"subsmoother: node sweeps"
=
4
"subsmoother: edge sweeps"
=
4
"coarse: type"
=
"Amesos-KLU"
"coarse: max size"
=
128
References Teuchos::ParameterList::set(), and Teuchos::ParameterList::setName().
Referenced by SetDefaults().
int ML_Epetra::SetDefaultsNSSA | ( | Teuchos::ParameterList & | List, |
Teuchos::RCP< std::vector< int > > & | options, | ||
Teuchos::RCP< std::vector< double > > & | params, | ||
bool | OverWrite = true |
||
) |
Sets defaults for energy minimization preconditioning for nonsymmetric problems.
Set default values for smoothed aggregation for nonsymmetric problems:
Note: This method should not be called directly. Instead, ML_Epetra::SetDefaults("NSSA",...) should be used.
"default values"
=
"NSSA"
"max levels"
=
10
"prec type"
=
"MGW"
"increasing or decreasing"
=
"increasing"
"aggregation: type"
=
"Uncoupled-MIS"
"energy minimization: enable"
=
true
"eigen-analysis: type"
=
"power-method"
"eigen-analysis: iterations"
=
20
"smoother: sweeps"
=
4
"smoother: damping factor"
=
0.67
"smoother: pre or post"
=
"post"
"smoother: type"
=
"symmetric Gauss-Seidel"
"coarse: type"
=
"Amesos-KLU"
"coarse: max size"
=
256
References Teuchos::ParameterList::set(), and Teuchos::ParameterList::setName().
Referenced by SetDefaults().
int ML_Epetra::SetDefaultsSA | ( | Teuchos::ParameterList & | List, |
Teuchos::RCP< std::vector< int > > & | options, | ||
Teuchos::RCP< std::vector< double > > & | params, | ||
bool | OverWrite = true |
||
) |
Sets default parameters for classical smoothed aggregation.
Set default values for classical smoothed aggregation.
This method should not be called directly. Instead, ML_Epetra::SetDefaults("SA",...) should be used.
"default values"
=
"SA"
"max levels"
=
10
"prec type"
=
"MGV"
"increasing or decreasing"
=
"increasing"
"aggregation: type"
=
"Uncoupled-MIS"
"aggregation: damping factor"
=
1.333
"eigen-analysis: type"
=
"cg"
"eigen-analysis: iterations"
=
10
"smoother: sweeps"
=
2
"smoother: damping factor"
=
1.0
"smoother: pre or post"
=
"both"
"smoother: type"
=
"symmetric Gauss-Seidel"
"coarse: type"
=
"Amesos-KLU"
"coarse: max size"
=
128
References Teuchos::ParameterList::set(), and Teuchos::ParameterList::setName().
Referenced by SetDefaults().