Ifpack Package Browser (Single Doxygen Collection) Development
Loading...
Searching...
No Matches
Ifpack_ex_Amesos_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_Map.h"
39#include "Epetra_BlockMap.h"
40#include "Epetra_LinearProblem.h"
41#include "Galeri_Maps.h"
42#include "Galeri_CrsMatrices.h"
43#include "Teuchos_ParameterList.hpp"
44#include "Teuchos_RefCountPtr.hpp"
45#include "AztecOO.h"
46#include "Ifpack.h"
47
48int main(int argc, char *argv[])
49{
50
51 // initialize MPI and Epetra communicator
52#ifdef HAVE_MPI
53 MPI_Init(&argc,&argv);
54 Epetra_MpiComm Comm( MPI_COMM_WORLD );
55#else
57#endif
58
59 Teuchos::ParameterList GaleriList;
60
61 // The problem is defined on a 2D grid, global size is nx * nx.
62 int nx = 30;
63 GaleriList.set("nx", nx);
64 GaleriList.set("ny", nx * Comm.NumProc());
65 GaleriList.set("mx", 1);
66 GaleriList.set("my", Comm.NumProc());
67 Teuchos::RefCountPtr<Epetra_Map> Map = Teuchos::rcp( Galeri::CreateMap64("Cartesian2D", Comm, GaleriList) );
68 Teuchos::RefCountPtr<Epetra_RowMatrix> A = Teuchos::rcp( Galeri::CreateCrsMatrix("Laplace2D", &*Map, GaleriList) );
69
70 // =============================================================== //
71 // 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 //
72 // =============================================================== //
73
74 Teuchos::ParameterList List;
75
76 // allocates an IFPACK factory. No data is associated
77 // to this object (only method Create()).
78 Ifpack Factory;
79
80 // create the preconditioner. For valid PrecType values,
81 // please check the documentation
82 std::string PrecType = "Amesos";
83 int OverlapLevel = 2; // must be >= 0. If Comm.NumProc() == 1,
84 // it is ignored.
85
86 Teuchos::RefCountPtr<Ifpack_Preconditioner> Prec = Teuchos::rcp( Factory.Create(PrecType, &*A, OverlapLevel) );
87 assert(Prec != Teuchos::null);
88
89 // specify the Amesos solver to be used.
90 // If the selected solver is not available,
91 // IFPACK will try to use Amesos' KLU (which is usually always
92 // compiled). Amesos' serial solvers are:
93 // "Amesos_Klu", "Amesos_Umfpack", "Amesos_Superlu"
94 List.set("amesos: solver type", "Amesos_Klu");
95
96 // sets the parameters
97 IFPACK_CHK_ERR(Prec->SetParameters(List));
98
99 // initialize the preconditioner. At this point the matrix must
100 // have been FillComplete()'d, but actual values are ignored.
101 // At this call, Amesos will perform the symbolic factorization.
102 IFPACK_CHK_ERR(Prec->Initialize());
103
104 // Builds the preconditioners, by looking for the values of
105 // the matrix. At this call, Amesos will perform the
106 // numeric factorization.
107 IFPACK_CHK_ERR(Prec->Compute());
108
109 // =================================================== //
110 // E N D O F I F P A C K C O N S T R U C T I O N //
111 // =================================================== //
112
113 // At this point, we need some additional objects
114 // to define and solve the linear system.
115
116 // defines LHS and RHS
117 Epetra_Vector LHS(A->OperatorDomainMap());
118 Epetra_Vector RHS(A->OperatorDomainMap());
119
120 // solution is constant
121 LHS.PutScalar(1.0);
122 // now build corresponding RHS
123 A->Apply(LHS,RHS);
124
125 // now randomize the solution
126 RHS.Random();
127
128 // need an Epetra_LinearProblem to define AztecOO solver
129 Epetra_LinearProblem Problem(&*A,&LHS,&RHS);
130
131 // now we can allocate the AztecOO solver
132 AztecOO Solver(Problem);
133
134 // specify solver
135 Solver.SetAztecOption(AZ_solver,AZ_gmres);
136 Solver.SetAztecOption(AZ_output,32);
137
138 // HERE WE SET THE IFPACK PRECONDITIONER
139 Solver.SetPrecOperator(&*Prec);
140
141 // .. and here we solve
142 // NOTE: with one process, the solver must converge in
143 // one iteration.
144 Solver.Iterate(1550,1e-8);
145
146#ifdef HAVE_MPI
147 MPI_Finalize() ;
148#endif
149
150 return(EXIT_SUCCESS);
151}
Solver
#define IFPACK_CHK_ERR(ifpack_err)
int main(int argc, char *argv[])
#define RHS(a)
Definition: MatGenFD.c:60
int NumProc() const
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