46 #include "Epetra_ConfigDefs.h"
49 #include "Epetra_MpiComm.h"
51 #include "Epetra_SerialComm.h"
53 #include "Epetra_Comm.h"
54 #include "Epetra_Map.h"
55 #include "Epetra_Time.h"
56 #include "Epetra_BlockMap.h"
57 #include "Epetra_MultiVector.h"
58 #include "Epetra_Vector.h"
59 #include "Epetra_Export.h"
61 #include "Galeri_Maps.h"
62 #include "Galeri_CrsMatrices.h"
65 #include "Teuchos_RefCountPtr.hpp"
72 return std::string(s);
78 return std::string(s);
83 int main(
int argc,
char *argv[]) {
86 MPI_Init(&argc,&argv);
95 GaleriList.
set(
"nx", nx);
97 GaleriList.
set(
"mx", 1);
99 Teuchos::RefCountPtr<Epetra_Map> Map =
Teuchos::rcp( Galeri::CreateMap64(
"Cartesian2D", Comm, GaleriList) );
100 Teuchos::RefCountPtr<Epetra_CrsMatrix>
A =
Teuchos::rcp( Galeri::CreateCrsMatrix(
"Laplace2D", &*Map, GaleriList) );
103 LHS->PutScalar(0.0); RHS->Random();
112 double DropTol = 0.3333;
113 double Athresh = 0.0123;
114 double Rthresh = 0.9876;
118 Teuchos::RefCountPtr<Ifpack_IlukGraph>
Graph;
119 Teuchos::RefCountPtr<Ifpack_CrsRiluk> RILU;
123 ierr = Graph->ConstructFilledGraph();
127 RILU->SetAbsoluteThreshold(Athresh);
128 RILU->SetRelativeThreshold(Rthresh);
129 RILU->SetRelaxValue(Relax);
130 int initerr = RILU->InitValues(*A);
131 if (initerr!=0) cout << Comm <<
"*ERR* InitValues = " << initerr;
136 std::string label =
"Ifpack_CrsRiluk Preconditioner: LevelFill = " +
toString(LevelFill) +
146 solver.SetUserMatrix(&*A);
147 solver.SetLHS(&*LHS);
148 solver.SetRHS(&*RHS);
149 solver.SetAztecOption(AZ_solver,AZ_gmres);
150 solver.SetPrecOperator(&*RILU);
151 solver.SetAztecOption(AZ_output, 16);
152 solver.Iterate(Niters, 5.0e-5);
154 int OldIters = solver.NumIters();
160 Teuchos::RefCountPtr<Ifpack_Preconditioner> Prec =
Teuchos::rcp( Factory.
Create(
"ILU", &*A, Overlap) );
163 List.
get(
"fact: level-of-fill", LevelFill);
164 List.get(
"fact: drop tolerance", DropTol);
165 List.get(
"fact: absolute threshold", Athresh);
166 List.get(
"fact: relative threshold", Rthresh);
167 List.get(
"fact: relax value", Relax);
175 solver.SetUserMatrix(&*A);
176 solver.SetLHS(&*LHS);
177 solver.SetRHS(&*RHS);
178 solver.SetAztecOption(AZ_solver,AZ_gmres);
179 solver.SetPrecOperator(&*Prec);
180 solver.SetAztecOption(AZ_output, 16);
181 solver.Iterate(Niters, 5.0e-5);
183 int NewIters = solver.NumIters();
185 if (OldIters != NewIters)
189 #ifdef HAVE_IFPACK_SUPERLU
191 Teuchos::RefCountPtr<Ifpack_Preconditioner> Prec2 =
Teuchos::rcp( Factory.
Create(
"SILU", &*A,0) );
193 SList.
set(
"fact: drop tolerance",1e-4);
194 SList.
set(
"fact: zero pivot threshold",.1);
195 SList.
set(
"fact: maximum fill factor",10.0);
198 SList.
set(
"fact: silu drop rule",9);
204 solver.SetPrecOperator(&*Prec2);
205 solver.Iterate(Niters, 5.0e-5);
215 return(EXIT_SUCCESS);
Ifpack_CrsRiluk: A class for constructing and using an incomplete lower/upper (ILU) factorization of ...
T & get(ParameterList &l, const std::string &name)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
std::string toString(const int &x)
Ifpack_IlukGraph: A class for constructing level filled graphs for use with ILU(k) class precondition...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
int main(int argc, char *argv[])
adjacency_list< vecS, vecS, undirectedS, no_property, property< edge_weight_t, double > > Graph
Ifpack: a function class to define Ifpack preconditioners.
static Ifpack_Preconditioner * Create(EPrecType PrecType, Epetra_RowMatrix *Matrix, const int overlap=0, bool overrideSerialDefault=false)
Creates an instance of Ifpack_Preconditioner given the enum value of the preconditioner type (can not...
#define IFPACK_CHK_ERR(ifpack_err)