Ifpack Package Browser (Single Doxygen Collection) Development
Loading...
Searching...
No Matches
Ifpack_ex_BlockRelaxation_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"
39#include "Galeri_Maps.h"
40#include "Galeri_CrsMatrices.h"
41#include "Teuchos_ParameterList.hpp"
42#include "Teuchos_RefCountPtr.hpp"
43#include "AztecOO.h"
48#include "Ifpack_Amesos.h"
49
50int main(int argc, char *argv[])
51{
52 // initialize MPI and Epetra communicator
53#ifdef HAVE_MPI
54 MPI_Init(&argc,&argv);
55 Epetra_MpiComm Comm( MPI_COMM_WORLD );
56#else
58#endif
59
60 Teuchos::ParameterList GaleriList;
61
62 // The problem is defined on a 2D grid, global size is nx * nx.
63 int nx = 30;
64 GaleriList.set("nx", nx);
65 GaleriList.set("ny", nx * Comm.NumProc());
66 GaleriList.set("mx", 1);
67 GaleriList.set("my", Comm.NumProc());
68 Teuchos::RefCountPtr<Epetra_Map> Map = Teuchos::rcp( Galeri::CreateMap64("Cartesian2D", Comm, GaleriList) );
69 Teuchos::RefCountPtr<Epetra_RowMatrix> A = Teuchos::rcp( Galeri::CreateCrsMatrix("Laplace2D", &*Map, GaleriList) );
70
71 // =============================================================== //
72 // 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 //
73 // =============================================================== //
74
75 Teuchos::ParameterList List;
76
77 // builds an Ifpack_AdditiveSchwarz. This is templated with
78 // the local solvers, in this case Ifpack_BlockRelaxation.
79 // Ifpack_BlockRelaxation requires as a templated a container
80 // class. A container defines
81 // how to store the diagonal blocks. Two choices are available:
82 // Ifpack_DenseContainer (to store them as dense block,
83 // than use LAPACK' factorization to apply the inverse of
84 // each block), of Ifpack_SparseContainer (to store
85 // the diagonal block as Epetra_CrsMatrix's).
86 //
87 // Here, we use Ifpack_SparseContainer, which in turn is
88 // templated with the class to use to apply the inverse
89 // of each block. For example, we can use Ifpack_Amesos.
90
91 // We still have to decide the overlap among the processes,
92 // and the overlap among the blocks. The two values
93 // can be different. The overlap among the blocks is
94 // considered only if block Jacobi is used.
95 int OverlapProcs = 2;
96 int OverlapBlocks = 0;
97
98 // define the block below to use dense containers
99#if 0
100 Ifpack_AdditiveSchwarz<Ifpack_BlockRelaxation<Ifpack_DenseContainer> > Prec(A, OverlapProcs);
101#else
102 Ifpack_AdditiveSchwarz<Ifpack_BlockRelaxation<Ifpack_SparseContainer<Ifpack_Amesos> > > Prec(&*A, OverlapProcs);
103#endif
104
105 List.set("relaxation: type", "symmetric Gauss-Seidel");
106 List.set("partitioner: overlap", OverlapBlocks);
107#ifdef HAVE_IFPACK_METIS
108 // use METIS to create the blocks. This requires --enable-ifpack-metis.
109 // If METIS is not installed, the user may select "linear".
110 List.set("partitioner: type", "metis");
111#else
112 // or a simple greedy algorithm is METIS is not enabled
113 List.set("partitioner: type", "greedy");
114#endif
115 // defines here the number of local blocks. If 1,
116 // and only one process is used in the computation, then
117 // the preconditioner must converge in one iteration.
118 List.set("partitioner: local parts", 4);
119
120 // sets the parameters
121 IFPACK_CHK_ERR(Prec.SetParameters(List));
122
123 // initialize the preconditioner.
124 IFPACK_CHK_ERR(Prec.Initialize());
125
126 // Builds the preconditioners.
127 IFPACK_CHK_ERR(Prec.Compute());
128
129 // =================================================== //
130 // E N D O F I F P A C K C O N S T R U C T I O N //
131 // =================================================== //
132
133 // At this point, we need some additional objects
134 // to define and solve the linear system.
135
136 // defines LHS and RHS
137 Epetra_Vector LHS(A->OperatorDomainMap());
138 Epetra_Vector RHS(A->OperatorDomainMap());
139
140 LHS.PutScalar(0.0);
141 RHS.Random();
142
143 // need an Epetra_LinearProblem to define AztecOO solver
144 Epetra_LinearProblem Problem(&*A,&LHS,&RHS);
145
146 // now we can allocate the AztecOO solver
147 AztecOO Solver(Problem);
148
149 // specify solver
150 Solver.SetAztecOption(AZ_solver,AZ_cg);
151 Solver.SetAztecOption(AZ_output,32);
152
153 // HERE WE SET THE IFPACK PRECONDITIONER
154 Solver.SetPrecOperator(&Prec);
155
156 // .. and here we solve
157 // NOTE: with one process, the solver must converge in
158 // one iteration.
159 Solver.Iterate(1550,1e-5);
160
161#ifdef HAVE_MPI
162 MPI_Finalize() ;
163#endif
164
165 return(EXIT_SUCCESS);
166}
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)