ML
Version of the Day
|
#include "ml_include.h"
#include "Epetra_SerialDenseMatrix.h"
#include "Epetra_SerialDenseVector.h"
#include "Epetra_SerialDenseSolver.h"
#include "Epetra_Operator.h"
#include "Epetra_RowMatrix.h"
#include "Epetra_MultiVector.h"
#include "Epetra_MsrMatrix.h"
#include "Teuchos_ParameterList.hpp"
#include "EpetraExt_SolverMap_CrsMatrix.h"
#include "ml_epetra_utils.h"
Classes | |
class | ML_Epetra::MultiLevelPreconditioner |
MultiLevelPreconditioner: a class to define black-box multilevel preconditioners using aggregation methods. More... | |
class | MLVec< T > |
struct | wrappedCommStruct |
Namespaces | |
ML_Epetra | |
ML_Epetra: default namespace for all Epetra interfaces. | |
Enumerations | |
enum | ML_Epetra::AMGType { ML_SA_FAMILY, ML_MAXWELL, ML_COMPOSITE, ML_CLASSICAL_FAMILY } |
Enumerated type indicating the type of AMG solver to be used. | |
enum | free_s { useFree, useDelete, neither } |
enum | epetraOrMLOrMueLu { epetraType, mlType, mueluType } |
Functions | |
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. More... | |
int | ML_Epetra::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 | ML_Epetra::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 | ML_Epetra::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 | ML_Epetra::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 | ML_Epetra::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 | ML_Epetra::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 | ML_Epetra::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 | ML_Epetra::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 | ML_Epetra::ReadXML (const std::string &FileName, Teuchos::ParameterList &List, const Epetra_Comm &Comm) |
Reads in parameter list options from file. | |
template<class T > | |
int | nodalComm (MLVec< T > &vector, MLVec< int > &myLocalNodeIds, struct wrappedCommStruct &framework) |
template<class T > | |
int | dofComm (MLVec< T > &vector, struct wrappedCommStruct &framework) |
int | dofCommUsingMlNodalMatrix (double *data, void *widget) |
int | MLsortCols (const MLVec< int > &ARowPtr, MLVec< int > &ACols, MLVec< double > &AVals) |
int | MLnMyGhost (struct wrappedCommStruct &framework) |
int | MLfillNodalMaps (MLVec< int > &amalgRowMap, MLVec< int > &amalgColMap, MLVec< int > &myLocalNodeIds, int nLocalDofs, struct wrappedCommStruct &framework, int nLocalNodes, int nLocalPlusGhostNodes) |
int | MLcolGlobalIds (struct wrappedCommStruct &framework, MLVec< int > &myGids) |
int | MLassignGhostLocalNodeIds (MLVec< int > &myLocalNodeIds, int nLocalDofs, int nLocalPlusGhostDofs, struct wrappedCommStruct &framework, int &nLocalNodes, int &nLocalPlusGhostNodes) |
int | MLextractDiag (const MLVec< int > &rowPtr, const MLVec< int > &cols, const MLVec< double > &, MLVec< double > &diagonal, struct wrappedCommStruct &framework) |
int | MLfindDirichlets (const MLVec< int > &rowPtr, const MLVec< int > &cols, const MLVec< double > &vals, const MLVec< double > &diagonal, double tol, MLVec< bool > &dirOrNot, struct wrappedCommStruct &framework) |
int | MLrmDirichletCols (MLVec< int > &rowPtr, MLVec< int > &cols, MLVec< double > &vals, const MLVec< double > &diagonal, bool squeeze, MLVec< double > &solution, MLVec< double > &rhs, const MLVec< bool > &dirOrNot, struct wrappedCommStruct &framework) |
int | MLsqueezeOutNnzs (MLVec< int > &rowPtr, MLVec< int > &cols, MLVec< double > &vals, const MLVec< bool > &keep) |
int | MLbuildMap (const MLVec< bool > &dofPresent, MLVec< int > &map, int NDof) |
int | MLvariableDofAmalg (int nCols, const MLVec< int > &rowPtr, const MLVec< int > &cols, const MLVec< double > &vals, int nNodes, int maxDofPerNode, const MLVec< int > &map, const MLVec< double > &diag, double tol, MLVec< int > &amalgRowPtr, MLVec< int > &amalgCols, struct wrappedCommStruct &framework, MLVec< int > &myLocalNodeIds) |
int | MLrmDifferentDofsCrossings (const MLVec< bool > &dofPresent, int maxDofPerNode, MLVec< int > &rowPtr, MLVec< int > &cols, int nCols, struct wrappedCommStruct &framework, MLVec< int > &myLocalNodeIds) |
int | MLbuildLaplacian (const MLVec< int > &rowPtr, const MLVec< int > &cols, MLVec< double > &vals, const MLVec< double > &x, const MLVec< double > &y, const MLVec< double > &z) |
int | MLunamalgP (const MLVec< int > &amalgRowPtr, const MLVec< int > &amalgCols, const MLVec< double > &amalgVals, int maxDofPerNode, const MLVec< char > &status, bool fineIsPadded, MLVec< int > &rowPtr, MLVec< int > &cols, MLVec< double > &vals) |
Unamalgamate prolongator so that it is suitable for PDE systems. More... | |
int | MLfindEmptyRows (const MLVec< int > &rowPtr, const MLVec< int > &cols, const MLVec< double > &vals, MLVec< bool > &rowEmpty) |
int | MLreplaceEmptyByDirichlet (MLVec< int > &rowPtr, MLVec< int > &cols, MLVec< double > &vals, const MLVec< bool > &colEmpty) |
int | MLfineStatus (const MLVec< bool > &dofPresent, const MLVec< int > &map, const MLVec< bool > &dirOrNot, MLVec< char > &status) |
int | MLcoarseStatus (const MLVec< bool > &rowEmpty, const MLVec< bool > &dirOrNot, MLVec< char > &status) |
int | MLShove (ML_Operator *Mat, MLVec< int > &rowPtr, MLVec< int > &cols, MLVec< double > &vals, int invec_leng, int(*commfunc)(double *vec, void *data), struct wrappedCommStruct &framework, int nGhost) |
Take Crs std::vec's, make arrays and shove them into an ML Operator. More... | |
int | ZeroDist (MLVec< double > &xxx, MLVec< double > &yyy, MLVec< double > &zzz, MLVec< int > &rowPtr, MLVec< int > &cols, MLVec< double > &vals, MLVec< double > &diagonal, double tol, MLVec< int > &rowZeros, MLVec< int > &colZeros, double disttol) |
int | MergeShared (MLVec< int > &cols, MLVec< int > &rowZeros, MLVec< int > &colZeros, MLVec< int > &groupHead, MLVec< int > &groupNext) |
int | BuildNonSharedMap (MLVec< int > &newMap, MLVec< int > &groupHead, MLVec< int > &groupNext) |
int | buildCompressedA (MLVec< int > &inputRowPtr, MLVec< int > &inputCols, MLVec< double > &inputVals, MLVec< double > &diagonal, MLVec< int > &groupHead, MLVec< int > &groupNext, MLVec< int > &outputRowptr, MLVec< int > &outputCols, MLVec< int > &map, int newN, double tol) |
int MLShove | ( | ML_Operator * | Mat, |
MLVec< int > & | rowPtr, | ||
MLVec< int > & | cols, | ||
MLVec< double > & | vals, | ||
int | invec_leng, | ||
int(*)(double *vec, void *data) | commfunc, | ||
struct wrappedCommStruct & | framework, | ||
int | nGhost | ||
) |
Take Crs std::vec's, make arrays and shove them into an ML Operator.
This includes allocating double arrays and a widget to point to them and calling all the appropriate ML functions to clean out what used to be in the operator before inserting the new matrix.
[in] | rowPtr,cols,vals | Crs matrix std:vectors |
Mat | On input, a populated ML Operator. On output, old ML Operator is cleaned and a new Crs matrix is filled into the operator |
References ML_Operator_Struct::N_nonzeros.
Referenced by ML_Epetra::MultiLevelPreconditioner::MultiLevelPreconditioner().
int MLunamalgP | ( | const MLVec< int > & | amalgRowPtr, |
const MLVec< int > & | amalgCols, | ||
const MLVec< double > & | amalgVals, | ||
int | maxDofPerNode, | ||
const MLVec< char > & | status, | ||
bool | fineIsPadded, | ||
MLVec< int > & | rowPtr, | ||
MLVec< int > & | cols, | ||
MLVec< double > & | vals | ||
) |
Unamalgamate prolongator so that it is suitable for PDE systems.
Unamalgamate Pmat (represented as a CRS matrix via amalgRowPtr, amalgCols, amalgVals for two different situations: Case 1: Fine level matrix is padded In this case, we basically replicate the unamalgamated operator, except that padded dofs and Dirichlet points use injection. Thus, PUnAmalg is equivalent to something like ( Pmat 0 0 ) ( 0 Pmat 0 ) ( 0 0 Pmat ) where entries associated with fine level padded dofs or Dir. BCs are replaced by rows with a single nonzero (whose value equals one) and the entire matrix is instead ordered so that dofs associated with nodes are consecutive.
Case 2: Fine level matrix is not padded. Thus, PUnAmalg is equivalent to ( Pmat 0 0 ) ( 0 Pmat 0 ) ( 0 0 Pmat ) where rows associated with dofs not present on finest level (e.g., an air pressure dof within water region) are removed and rows associated with Dir. BCs are replaced by rows with a single nonzero (whose value equals one) and the entire matrix is instead ordered so dofs associated with nodes are consecutive.
In both of the above cases, the coarse level discretization matrix is assumed to be padded.
amalgNRows local number of rows for amalgamted prolongator amalgRowPtr CRS local row pointer for amalgamted prolongator amalgCols CRS local cols for amalgamted prolongator amalgVals CRS vals for amalgamted prolongator maxDofPerNode maximum number of degrees-of-freedom at any mesh node status status[i*maxDofPerNode+j] refers to the jth dof at the ith node. status()==s ==> standard element: present in fine operator and not a Dirichlet BC. status()==d ==> element corresponds to Dirichlet BC. status()==p ==> element not present or is padded dof. fineIsPadded Indicates whether fine grid matrix includes padded dofs for those not really present in the PDE system rowPtr CRS local row pointer for resulting unamalgamated prolongator cols CRS local cols for resulting unamalgamated prolongator vals CRS vals for resulting unamalgamated prolongator
Referenced by ML_Epetra::MultiLevelPreconditioner::MultiLevelPreconditioner().