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