Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Ifpack_ex_ICT_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 "Epetra_Time.h"
40 #include "Galeri_Maps.h"
41 #include "Galeri_CrsMatrices.h"
42 #include "Teuchos_ParameterList.hpp"
43 #include "Teuchos_RefCountPtr.hpp"
44 #include "AztecOO.h"
45 #include "Ifpack_AdditiveSchwarz.h"
46 #include "Ifpack_ICT.h"
47 
48 int 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
56  Epetra_SerialComm Comm;
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 
75 
76  // builds an Ifpack_AdditiveSchwarz. This is templated with
77  // the local solvers, in this case Ifpack_ICT. Note that any
78  // other Ifpack_Preconditioner-derived class can be used
79  // instead of Ifpack_ICT.
80 
81  // In this example the overlap is zero. Use
82  // Prec(A,OverlapLevel) for the general case.
84 
85  // `1.0' means that the factorization should approximatively
86  // keep the same number of nonzeros per row of the original matrix.
87  List.set("fact: ict level-of-fill", 1.0);
88  // no modifications on the diagonal
89  List.set("fact: absolute threshold", 0.0);
90  List.set("fact: relative threshold", 1.0);
91  List.set("fact: relaxation value", 0.0);
92  // matrix `laplace_2d_bc' is not symmetric because of the way
93  // boundary conditions are imposed. We can filter the singletons,
94  // (that is, Dirichlet nodes) and end up with a symmetric
95  // matrix (as ICT requires).
96  List.set("schwarz: filter singletons", true);
97 
98  // sets the parameters
99  IFPACK_CHK_ERR(Prec.SetParameters(List));
100 
101  // initialize the preconditioner. At this point the matrix must
102  // have been FillComplete()'d, but actual values are ignored.
103  IFPACK_CHK_ERR(Prec.Initialize());
104 
105  // Builds the preconditioners, by looking for the values of
106  // the matrix.
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  LHS.PutScalar(0.0);
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_cg_condnum);
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  // NOTE: with one process, the solver must converge in
138  // one iteration.
139  Solver.Iterate(1550,1e-5);
140 
141  // Prints out some information about the preconditioner
142  std::cout << Prec;
143 
144 #ifdef HAVE_MPI
145  MPI_Finalize();
146 #endif
147 
148  return (EXIT_SUCCESS);
149 }
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
static bool Solver
Definition: performance.cpp:70
virtual int Initialize()
Initialized the preconditioner.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Ifpack_AdditiveSchwarz: a class to define Additive Schwarz preconditioners of Epetra_RowMatrix&#39;s.
int main(int argc, char *argv[])
int NumProc() const
virtual int SetParameters(Teuchos::ParameterList &List)
Sets the parameters.
virtual int Compute()
Computes the preconditioner.
#define IFPACK_CHK_ERR(ifpack_err)
#define RHS(a)
Definition: MatGenFD.c:60