45 #include "Epetra_MpiComm.h" 47 #include "Epetra_SerialComm.h" 49 #include "Epetra_CrsMatrix.h" 50 #include "Epetra_MultiVector.h" 51 #include "Epetra_LinearProblem.h" 52 #include "Galeri_Maps.h" 53 #include "Galeri_CrsMatrices.h" 54 #include "Teuchos_ParameterList.hpp" 55 #include "Teuchos_RefCountPtr.hpp" 63 int main(
int argc,
char *argv[])
67 MPI_Init(&argc,&argv);
68 Epetra_MpiComm Comm( MPI_COMM_WORLD );
70 Epetra_SerialComm Comm;
73 Teuchos::ParameterList GaleriList;
77 GaleriList.set(
"nx", nx);
78 GaleriList.set(
"ny", nx * Comm.NumProc());
79 GaleriList.set(
"mx", 1);
80 GaleriList.set(
"my", Comm.NumProc());
81 Teuchos::RefCountPtr<Epetra_Map> Map = Teuchos::rcp( Galeri::CreateMap(
"Cartesian2D", Comm, GaleriList) );
82 Teuchos::RefCountPtr<Epetra_RowMatrix> A = Teuchos::rcp( Galeri::CreateCrsMatrix(
"Laplace2D", &*Map, GaleriList) );
88 Teuchos::ParameterList List;
108 int OverlapProcs = 2;
109 int OverlapBlocks = 0;
118 List.set(
"relaxation: type",
"symmetric Gauss-Seidel");
119 List.set(
"partitioner: overlap", OverlapBlocks);
120 #ifdef HAVE_IFPACK_METIS 123 List.set(
"partitioner: type",
"metis");
126 List.set(
"partitioner: type",
"greedy");
131 List.set(
"partitioner: local parts", 4);
150 Epetra_Vector LHS(A->OperatorDomainMap());
151 Epetra_Vector
RHS(A->OperatorDomainMap());
157 Epetra_LinearProblem Problem(&*A,&LHS,&
RHS);
163 Solver.SetAztecOption(AZ_solver,AZ_cg);
164 Solver.SetAztecOption(AZ_output,32);
167 Solver.SetPrecOperator(&Prec);
172 Solver.Iterate(1550,1e-5);
178 return(EXIT_SUCCESS);
int main(int argc, char *argv[])
virtual int Initialize()
Initialized the preconditioner.
Ifpack_AdditiveSchwarz: a class to define Additive Schwarz preconditioners of Epetra_RowMatrix's.
virtual int SetParameters(Teuchos::ParameterList &List)
Sets the parameters.
virtual int Compute()
Computes the preconditioner.
#define IFPACK_CHK_ERR(ifpack_err)