Ifpack Package Browser (Single Doxygen Collection) Development
Loading...
Searching...
No Matches
Ifpack_ex_ICT.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 "Epetra_Time.h"
53#include "Galeri_Maps.h"
54#include "Galeri_CrsMatrices.h"
55#include "Teuchos_ParameterList.hpp"
56#include "Teuchos_RefCountPtr.hpp"
57#include "AztecOO.h"
59#include "Ifpack_ICT.h"
60
61int main(int argc, char *argv[])
62{
63
64 // initialize MPI and Epetra communicator
65#ifdef HAVE_MPI
66 MPI_Init(&argc,&argv);
67 Epetra_MpiComm Comm( MPI_COMM_WORLD );
68#else
70#endif
71
72 Teuchos::ParameterList GaleriList;
73
74 // The problem is defined on a 2D grid, global size is nx * nx.
75 int nx = 30;
76 GaleriList.set("nx", nx);
77 GaleriList.set("ny", nx * Comm.NumProc());
78 GaleriList.set("mx", 1);
79 GaleriList.set("my", Comm.NumProc());
80 Teuchos::RefCountPtr<Epetra_Map> Map = Teuchos::rcp( Galeri::CreateMap("Cartesian2D", Comm, GaleriList) );
81 Teuchos::RefCountPtr<Epetra_RowMatrix> A = Teuchos::rcp( Galeri::CreateCrsMatrix("Laplace2D", &*Map, GaleriList) );
82
83 // =============================================================== //
84 // 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 //
85 // =============================================================== //
86
87 Teuchos::ParameterList List;
88
89 // builds an Ifpack_AdditiveSchwarz. This is templated with
90 // the local solvers, in this case Ifpack_ICT. Note that any
91 // other Ifpack_Preconditioner-derived class can be used
92 // instead of Ifpack_ICT.
93
94 // In this example the overlap is zero. Use
95 // Prec(A,OverlapLevel) for the general case.
96 Ifpack_AdditiveSchwarz<Ifpack_ICT> Prec(&*A);
97
98 // `1.0' means that the factorization should approximatively
99 // keep the same number of nonzeros per row of the original matrix.
100 List.set("fact: ict level-of-fill", 1.0);
101 // no modifications on the diagonal
102 List.set("fact: absolute threshold", 0.0);
103 List.set("fact: relative threshold", 1.0);
104 List.set("fact: relaxation value", 0.0);
105 // matrix `laplace_2d_bc' is not symmetric because of the way
106 // boundary conditions are imposed. We can filter the singletons,
107 // (that is, Dirichlet nodes) and end up with a symmetric
108 // matrix (as ICT requires).
109 List.set("schwarz: filter singletons", true);
110
111 // sets the parameters
112 IFPACK_CHK_ERR(Prec.SetParameters(List));
113
114 // initialize the preconditioner. At this point the matrix must
115 // have been FillComplete()'d, but actual values are ignored.
116 IFPACK_CHK_ERR(Prec.Initialize());
117
118 // Builds the preconditioners, by looking for the values of
119 // the matrix.
120 IFPACK_CHK_ERR(Prec.Compute());
121
122 // =================================================== //
123 // E N D O F I F P A C K C O N S T R U C T I O N //
124 // =================================================== //
125
126 // At this point, we need some additional objects
127 // to define and solve the linear system.
128
129 // defines LHS and RHS
130 Epetra_Vector LHS(A->OperatorDomainMap());
131 Epetra_Vector RHS(A->OperatorDomainMap());
132
133 LHS.PutScalar(0.0);
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_cg_condnum);
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 // NOTE: with one process, the solver must converge in
151 // one iteration.
152 Solver.Iterate(1550,1e-5);
153
154 // Prints out some information about the preconditioner
155 std::cout << Prec;
156
157#ifdef HAVE_MPI
158 MPI_Finalize();
159#endif
160
161 return (EXIT_SUCCESS);
162}
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)