Ifpack Package Browser (Single Doxygen Collection) Development
Loading...
Searching...
No Matches
Ifpack_ex_Factory_LL.cpp
Go to the documentation of this file.
1// @HEADER
2// ***********************************************************************
3//
4// IFPACK
5// Copyright (2004) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// This library is free software; you can redistribute it and/or modify
11// it under the terms of the GNU Lesser General Public License as
12// published by the Free Software Foundation; either version 2.1 of the
13// License, or (at your option) any later version.
14//
15// This library is distributed in the hope that it will be useful, but
16// WITHOUT ANY WARRANTY; without even the implied warranty of
17// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18// Lesser General Public License for more details.
19//
20// You should have received a copy of the GNU Lesser General Public
21// License along with this library; if not, write to the Free Software
22// Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301
23// USA
24// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
25//
26// ***********************************************************************
27// @HEADER
28
29#include "Ifpack_ConfigDefs.h"
30
31#ifdef HAVE_MPI
32#include "Epetra_MpiComm.h"
33#else
34#include "Epetra_SerialComm.h"
35#endif
36#include "Epetra_CrsMatrix.h"
37#include "Epetra_MultiVector.h"
38#include "Epetra_LinearProblem.h"
39#include "Galeri_Maps.h"
40#include "Galeri_CrsMatrices.h"
41#include "Teuchos_ParameterList.hpp"
42#include "Teuchos_RefCountPtr.hpp"
43#include "AztecOO.h"
44#include "Ifpack.h"
46
47int main(int argc, char *argv[])
48{
49
50#ifdef HAVE_MPI
51 MPI_Init(&argc,&argv);
52 Epetra_MpiComm Comm( MPI_COMM_WORLD );
53#else
55#endif
56
57 Teuchos::ParameterList GaleriList;
58
59 // The problem is defined on a 2D grid, global size is nx * nx.
60 int nx = 30;
61 GaleriList.set("n", nx * nx);
62 GaleriList.set("nx", nx);
63 GaleriList.set("ny", nx);
64 Teuchos::RefCountPtr<Epetra_Map> Map = Teuchos::rcp( Galeri::CreateMap64("Linear", Comm, GaleriList) );
65 Teuchos::RefCountPtr<Epetra_RowMatrix> A = Teuchos::rcp( Galeri::CreateCrsMatrix("Laplace2D", &*Map, GaleriList) );
66
67 // =============================================================== //
68 // B E G I N N I N G O F I F P A C K C O N S T R U C T I O N //
69 // =============================================================== //
70
71 Teuchos::ParameterList List;
72
73 // allocates an IFPACK factory. No data is associated
74 // to this object (only method Create()).
75 Ifpack Factory;
76
77 // create the preconditioner. For valid PrecType values,
78 // please check the documentation
79 std::string PrecType = "ILU"; // incomplete LU
80 int OverlapLevel = 1; // must be >= 0. If Comm.NumProc() == 1,
81 // it is ignored.
82
83 Teuchos::RefCountPtr<Ifpack_Preconditioner> Prec = Teuchos::rcp( Factory.Create(PrecType, &*A, OverlapLevel) );
84 assert(Prec != Teuchos::null);
85
86 // specify parameters for ILU
87 List.set("fact: drop tolerance", 1e-9);
88 List.set("fact: level-of-fill", 1);
89 // the combine mode is on the following:
90 // "Add", "Zero", "Insert", "InsertAdd", "Average", "AbsMax"
91 // Their meaning is as defined in file Epetra_CombineMode.h
92 List.set("schwarz: combine mode", "Add");
93 // sets the parameters
94 IFPACK_CHK_ERR(Prec->SetParameters(List));
95
96 // initialize the preconditioner. At this point the matrix must
97 // have been FillComplete()'d, but actual values are ignored.
98 IFPACK_CHK_ERR(Prec->Initialize());
99
100 // Builds the preconditioners, by looking for the values of
101 // the matrix.
102 IFPACK_CHK_ERR(Prec->Compute());
103
104 // =================================================== //
105 // E N D O F I F P A C K C O N S T R U C T I O N //
106 // =================================================== //
107
108 // At this point, we need some additional objects
109 // to define and solve the linear system.
110
111 // defines LHS and RHS
112 Epetra_Vector LHS(A->OperatorDomainMap());
113 Epetra_Vector RHS(A->OperatorDomainMap());
114
115 // solution is constant
116 LHS.PutScalar(1.0);
117 // now build corresponding RHS
118 A->Apply(LHS,RHS);
119
120 // now randomize the solution
121 RHS.Random();
122
123 // need an Epetra_LinearProblem to define AztecOO solver
124 Epetra_LinearProblem Problem(&*A,&LHS,&RHS);
125
126 // now we can allocate the AztecOO solver
127 AztecOO Solver(Problem);
128
129 // specify solver
130 Solver.SetAztecOption(AZ_solver,AZ_gmres);
131 Solver.SetAztecOption(AZ_output,32);
132
133 // HERE WE SET THE IFPACK PRECONDITIONER
134 Solver.SetPrecOperator(&*Prec);
135
136 // .. and here we solve
137 Solver.Iterate(1550,1e-8);
138
139 std::cout << *Prec;
140
141#ifdef HAVE_MPI
142 MPI_Finalize() ;
143#endif
144
145 return(EXIT_SUCCESS);
146}
Solver
#define IFPACK_CHK_ERR(ifpack_err)
int main(int argc, char *argv[])
#define RHS(a)
Definition: MatGenFD.c:60
int PutScalar(double ScalarConstant)
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...
Definition: Ifpack.cpp:253