Ifpack Package Browser (Single Doxygen Collection) Development
Loading...
Searching...
No Matches
Ifpack_ex_Factory.cpp
Go to the documentation of this file.
1//@HEADER
2// ***********************************************************************
3//
4// Ifpack: Object-Oriented Algebraic Preconditioner Package
5// Copyright (2002) 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// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38//
39// ***********************************************************************
40//@HEADER
41
42#include "Ifpack_ConfigDefs.h"
43
44#ifdef HAVE_MPI
45#include "Epetra_MpiComm.h"
46#else
47#include "Epetra_SerialComm.h"
48#endif
49#include "Epetra_CrsMatrix.h"
50#include "Epetra_MultiVector.h"
52#include "Galeri_Maps.h"
53#include "Galeri_CrsMatrices.h"
54#include "Teuchos_ParameterList.hpp"
55#include "Teuchos_RefCountPtr.hpp"
56#include "AztecOO.h"
57#include "Ifpack.h"
59
60int main(int argc, char *argv[])
61{
62
63#ifdef HAVE_MPI
64 MPI_Init(&argc,&argv);
65 Epetra_MpiComm Comm( MPI_COMM_WORLD );
66#else
68#endif
69
70 Teuchos::ParameterList GaleriList;
71
72 // The problem is defined on a 2D grid, global size is nx * nx.
73 int nx = 30;
74 GaleriList.set("n", nx * nx);
75 GaleriList.set("nx", nx);
76 GaleriList.set("ny", nx);
77 Teuchos::RefCountPtr<Epetra_Map> Map = Teuchos::rcp( Galeri::CreateMap("Linear", Comm, GaleriList) );
78 Teuchos::RefCountPtr<Epetra_RowMatrix> A = Teuchos::rcp( Galeri::CreateCrsMatrix("Laplace2D", &*Map, GaleriList) );
79
80 // =============================================================== //
81 // 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 //
82 // =============================================================== //
83
84 Teuchos::ParameterList List;
85
86 // allocates an IFPACK factory. No data is associated
87 // to this object (only method Create()).
88 Ifpack Factory;
89
90 // create the preconditioner. For valid PrecType values,
91 // please check the documentation
92 std::string PrecType = "ILU"; // incomplete LU
93 int OverlapLevel = 1; // must be >= 0. If Comm.NumProc() == 1,
94 // it is ignored.
95
96 Teuchos::RefCountPtr<Ifpack_Preconditioner> Prec = Teuchos::rcp( Factory.Create(PrecType, &*A, OverlapLevel) );
97 assert(Prec != Teuchos::null);
98
99 // specify parameters for ILU
100 List.set("fact: drop tolerance", 1e-9);
101 List.set("fact: level-of-fill", 1);
102 // the combine mode is on the following:
103 // "Add", "Zero", "Insert", "InsertAdd", "Average", "AbsMax"
104 // Their meaning is as defined in file Epetra_CombineMode.h
105 List.set("schwarz: combine mode", "Add");
106 // sets the parameters
107 IFPACK_CHK_ERR(Prec->SetParameters(List));
108
109 // initialize the preconditioner. At this point the matrix must
110 // have been FillComplete()'d, but actual values are ignored.
111 IFPACK_CHK_ERR(Prec->Initialize());
112
113 // Builds the preconditioners, by looking for the values of
114 // the matrix.
115 IFPACK_CHK_ERR(Prec->Compute());
116
117 // =================================================== //
118 // E N D O F I F P A C K C O N S T R U C T I O N //
119 // =================================================== //
120
121 // At this point, we need some additional objects
122 // to define and solve the linear system.
123
124 // defines LHS and RHS
125 Epetra_Vector LHS(A->OperatorDomainMap());
126 Epetra_Vector RHS(A->OperatorDomainMap());
127
128 // solution is constant
129 LHS.PutScalar(1.0);
130 // now build corresponding RHS
131 A->Apply(LHS,RHS);
132
133 // now randomize the solution
134 RHS.Random();
135
136 // need an Epetra_LinearProblem to define AztecOO solver
137 Epetra_LinearProblem Problem(&*A,&LHS,&RHS);
138
139 // now we can allocate the AztecOO solver
140 AztecOO Solver(Problem);
141
142 // specify solver
143 Solver.SetAztecOption(AZ_solver,AZ_gmres);
144 Solver.SetAztecOption(AZ_output,32);
145
146 // HERE WE SET THE IFPACK PRECONDITIONER
147 Solver.SetPrecOperator(&*Prec);
148
149 // .. and here we solve
150 Solver.Iterate(1550,1e-8);
151
152 std::cout << *Prec;
153
154#ifdef HAVE_MPI
155 MPI_Finalize() ;
156#endif
157
158 return(EXIT_SUCCESS);
159}
Solver
#define IFPACK_CHK_ERR(ifpack_err)
int main(int argc, char *argv[])
#define RHS(a)
Definition MatGenFD.c:60
int PutScalar(double ScalarConstant)