Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/Relaxation_LL/cxx_main.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 
43 #include "Ifpack_ConfigDefs.h"
44 
45 #ifdef HAVE_MPI
46 #include "Epetra_MpiComm.h"
47 #else
48 #include "Epetra_SerialComm.h"
49 #endif
50 #include "Epetra_CrsMatrix.h"
51 #include "Epetra_Vector.h"
52 #include "Epetra_LinearProblem.h"
53 #include "Epetra_Map.h"
54 #include "Galeri_Maps.h"
55 #include "Galeri_CrsMatrices.h"
56 #include "Galeri_Utils.h"
57 #include "Teuchos_ParameterList.hpp"
58 #include "Teuchos_RefCountPtr.hpp"
59 #include "Ifpack_PointRelaxation.h"
60 #include "Ifpack_BlockRelaxation.h"
61 #include "Ifpack_SparseContainer.h"
62 #include "Ifpack_Amesos.h"
63 #include "AztecOO.h"
64 
65 static bool verbose = false;
66 static bool SymmetricGallery = false;
67 static bool Solver = AZ_gmres;
68 const int NumVectors = 3;
69 
70 // ======================================================================
71 int CompareBlockOverlap(const Teuchos::RefCountPtr<Epetra_RowMatrix>& A, int Overlap)
72 {
73  Epetra_MultiVector LHS(A->RowMatrixRowMap(), NumVectors);
74  Epetra_MultiVector RHS(A->RowMatrixRowMap(), NumVectors);
75  LHS.PutScalar(0.0); RHS.Random();
76 
77  Epetra_LinearProblem Problem(&*A, &LHS, &RHS);
78 
80  List.set("relaxation: damping factor", 1.0);
81  List.set("relaxation: type", "Jacobi");
82  List.set("relaxation: sweeps",1);
83  List.set("partitioner: overlap", Overlap);
84  List.set("partitioner: type", "linear");
85  List.set("partitioner: local parts", 16);
86 
87  RHS.PutScalar(1.0);
88  LHS.PutScalar(0.0);
89 
91  Prec.SetParameters(List);
92  Prec.Compute();
93 
94  // set AztecOO solver object
95  AztecOO AztecOOSolver(Problem);
96  AztecOOSolver.SetAztecOption(AZ_solver,Solver);
97  if (verbose)
98  AztecOOSolver.SetAztecOption(AZ_output,32);
99  else
100  AztecOOSolver.SetAztecOption(AZ_output,AZ_none);
101  AztecOOSolver.SetPrecOperator(&Prec);
102 
103  AztecOOSolver.Iterate(1550,1e-5);
104 
105  return(AztecOOSolver.NumIters());
106 }
107 
108 // ======================================================================
109 int CompareBlockSizes(std::string PrecType, const Teuchos::RefCountPtr<Epetra_RowMatrix>& A, int NumParts)
110 {
111  Epetra_MultiVector RHS(A->RowMatrixRowMap(), NumVectors);
112  Epetra_MultiVector LHS(A->RowMatrixRowMap(), NumVectors);
113  LHS.PutScalar(0.0); RHS.Random();
114 
115  Epetra_LinearProblem Problem(&*A, &LHS, &RHS);
116 
118  List.set("relaxation: damping factor", 1.0);
119  List.set("relaxation: type", PrecType);
120  List.set("relaxation: sweeps",1);
121  List.set("partitioner: type", "linear");
122  List.set("partitioner: local parts", NumParts);
123 
124  RHS.PutScalar(1.0);
125  LHS.PutScalar(0.0);
126 
128  Prec.SetParameters(List);
129  Prec.Compute();
130 
131  // set AztecOO solver object
132  AztecOO AztecOOSolver(Problem);
133  AztecOOSolver.SetAztecOption(AZ_solver,Solver);
134  if (verbose)
135  AztecOOSolver.SetAztecOption(AZ_output,32);
136  else
137  AztecOOSolver.SetAztecOption(AZ_output,AZ_none);
138  AztecOOSolver.SetPrecOperator(&Prec);
139 
140  AztecOOSolver.Iterate(1550,1e-5);
141 
142  return(AztecOOSolver.NumIters());
143 }
144 
145 // ======================================================================
146 bool ComparePointAndBlock(std::string PrecType, const Teuchos::RefCountPtr<Epetra_RowMatrix>& A, int sweeps)
147 {
148  Epetra_MultiVector RHS(A->RowMatrixRowMap(), NumVectors);
149  Epetra_MultiVector LHS(A->RowMatrixRowMap(), NumVectors);
150  LHS.PutScalar(0.0); RHS.Random();
151 
152  Epetra_LinearProblem Problem(&*A, &LHS, &RHS);
153 
154  // Set up the list
156  List.set("relaxation: damping factor", 1.0);
157  List.set("relaxation: type", PrecType);
158  List.set("relaxation: sweeps",sweeps);
159  List.set("partitioner: type", "linear");
160  List.set("partitioner: local parts", A->NumMyRows());
161 
162  int ItersPoint, ItersBlock;
163 
164  // ================================================== //
165  // get the number of iterations with point relaxation //
166  // ================================================== //
167  {
168  RHS.PutScalar(1.0);
169  LHS.PutScalar(0.0);
170 
171  Ifpack_PointRelaxation Point(&*A);
172  Point.SetParameters(List);
173  Point.Compute();
174 
175  // set AztecOO solver object
176  AztecOO AztecOOSolver(Problem);
177  AztecOOSolver.SetAztecOption(AZ_solver,Solver);
178  if (verbose)
179  AztecOOSolver.SetAztecOption(AZ_output,32);
180  else
181  AztecOOSolver.SetAztecOption(AZ_output,AZ_none);
182  AztecOOSolver.SetPrecOperator(&Point);
183 
184  AztecOOSolver.Iterate(1550,1e-2);
185 
186  double TrueResidual = AztecOOSolver.TrueResidual();
187  ItersPoint = AztecOOSolver.NumIters();
188  // some output
189  if (verbose && Problem.GetMatrix()->Comm().MyPID() == 0) {
190  cout << "Iterations = " << ItersPoint << endl;
191  cout << "Norm of the true residual = " << TrueResidual << endl;
192  }
193  }
194 
195  // ================================================== //
196  // get the number of iterations with block relaxation //
197  // ================================================== //
198  {
199 
200  RHS.PutScalar(1.0);
201  LHS.PutScalar(0.0);
202 
204  Block.SetParameters(List);
205  Block.Compute();
206 
207  // set AztecOO solver object
208  AztecOO AztecOOSolver(Problem);
209  AztecOOSolver.SetAztecOption(AZ_solver,Solver);
210  if (verbose)
211  AztecOOSolver.SetAztecOption(AZ_output,32);
212  else
213  AztecOOSolver.SetAztecOption(AZ_output,AZ_none);
214  AztecOOSolver.SetPrecOperator(&Block);
215 
216  AztecOOSolver.Iterate(1550,1e-2);
217 
218  double TrueResidual = AztecOOSolver.TrueResidual();
219  ItersBlock = AztecOOSolver.NumIters();
220  // some output
221  if (verbose && Problem.GetMatrix()->Comm().MyPID() == 0) {
222  cout << "Iterations " << ItersBlock << endl;
223  cout << "Norm of the true residual = " << TrueResidual << endl;
224  }
225  }
226 
227  int diff = ItersPoint - ItersBlock;
228  if (diff < 0) diff = -diff;
229 
230  if (diff > 10)
231  {
232  if (verbose)
233  cout << "TEST FAILED!" << endl;
234  return(false);
235  }
236  else {
237  if (verbose)
238  cout << "TEST PASSED" << endl;
239  return(true);
240  }
241 }
242 
243 // ======================================================================
244 bool KrylovTest(std::string PrecType, const Teuchos::RefCountPtr<Epetra_RowMatrix>& A, bool backward)
245 {
246  Epetra_MultiVector LHS(A->RowMatrixRowMap(), NumVectors);
247  Epetra_MultiVector RHS(A->RowMatrixRowMap(), NumVectors);
248  LHS.PutScalar(0.0); RHS.Random();
249 
250  Epetra_LinearProblem Problem(&*A, &LHS, &RHS);
251 
252  // Set up the list
254  List.set("relaxation: damping factor", 1.0);
255  List.set("relaxation: type", PrecType);
256  if(backward) List.set("relaxation: backward mode",backward);
257 
258  int Iters1, Iters10;
259 
260  if (verbose) {
261  cout << "Krylov test: Using " << PrecType
262  << " with AztecOO" << endl;
263  }
264 
265  // ============================================== //
266  // get the number of iterations with 1 sweep only //
267  // ============================================== //
268  {
269 
270  List.set("relaxation: sweeps",1);
271  Ifpack_PointRelaxation Point(&*A);
272  Point.SetParameters(List);
273  Point.Compute();
274 
275  // set AztecOO solver object
276  AztecOO AztecOOSolver(Problem);
277  AztecOOSolver.SetAztecOption(AZ_solver,Solver);
278  AztecOOSolver.SetAztecOption(AZ_output,AZ_none);
279  AztecOOSolver.SetPrecOperator(&Point);
280 
281  AztecOOSolver.Iterate(1550,1e-5);
282 
283  double TrueResidual = AztecOOSolver.TrueResidual();
284  // some output
285  if (verbose && Problem.GetMatrix()->Comm().MyPID() == 0) {
286  cout << "Norm of the true residual = " << TrueResidual << endl;
287  }
288  Iters1 = AztecOOSolver.NumIters();
289  }
290 
291  // ======================================================== //
292  // now re-run with 10 sweeps, solver should converge faster
293  // ======================================================== //
294  {
295  List.set("relaxation: sweeps",10);
296  Ifpack_PointRelaxation Point(&*A);
297  Point.SetParameters(List);
298  Point.Compute();
299  LHS.PutScalar(0.0);
300 
301  // set AztecOO solver object
302  AztecOO AztecOOSolver(Problem);
303  AztecOOSolver.SetAztecOption(AZ_solver,Solver);
304  AztecOOSolver.SetAztecOption(AZ_output,AZ_none);
305  AztecOOSolver.SetPrecOperator(&Point);
306  AztecOOSolver.Iterate(1550,1e-5);
307 
308  double TrueResidual = AztecOOSolver.TrueResidual();
309  // some output
310  if (verbose && Problem.GetMatrix()->Comm().MyPID() == 0) {
311  cout << "Norm of the true residual = " << TrueResidual << endl;
312  }
313  Iters10 = AztecOOSolver.NumIters();
314  }
315 
316  if (verbose) {
317  cout << "Iters_1 = " << Iters1 << ", Iters_10 = " << Iters10 << endl;
318  cout << "(second number should be smaller than first one)" << endl;
319  }
320 
321  if (Iters10 > Iters1) {
322  if (verbose)
323  cout << "TEST FAILED!" << endl;
324  return(false);
325  }
326  else {
327  if (verbose)
328  cout << "TEST PASSED" << endl;
329  return(true);
330  }
331 }
332 
333 // ======================================================================
334 bool BasicTest(std::string PrecType, const Teuchos::RefCountPtr<Epetra_RowMatrix>& A,bool backward)
335 {
336  Epetra_MultiVector LHS(A->RowMatrixRowMap(), NumVectors);
337  Epetra_MultiVector RHS(A->RowMatrixRowMap(), NumVectors);
338  LHS.PutScalar(0.0); RHS.Random();
339 
340  double starting_residual = Galeri::ComputeNorm(&*A, &LHS, &RHS);
341  Epetra_LinearProblem Problem(&*A, &LHS, &RHS);
342 
343  // Set up the list
345  List.set("relaxation: damping factor", 1.0);
346  List.set("relaxation: sweeps",1550);
347  List.set("relaxation: type", PrecType);
348  if(backward) List.set("relaxation: backward mode",backward);
349 
350 
351  Ifpack_PointRelaxation Point(&*A);
352 
353  Point.SetParameters(List);
354  Point.Compute();
355  // use the preconditioner as solver, with 1550 iterations
356  Point.ApplyInverse(RHS,LHS);
357 
358  // compute the real residual
359 
360  double residual = Galeri::ComputeNorm(&*A, &LHS, &RHS);
361 
362  if (A->Comm().MyPID() == 0 && verbose)
363  cout << "||A * x - b||_2 (scaled) = " << residual / starting_residual << endl;
364 
365  // Jacobi is very slow to converge here
366  if (residual / starting_residual < 1e-2) {
367  if (verbose)
368  cout << "Test passed" << endl;
369  return(true);
370  }
371  else {
372  if (verbose)
373  cout << "Test failed!" << endl;
374  return(false);
375  }
376 }
377 
378 // ======================================================================
379 int main(int argc, char *argv[])
380 {
381 #ifdef HAVE_MPI
382  MPI_Init(&argc,&argv);
383  Epetra_MpiComm Comm( MPI_COMM_WORLD );
384 #else
385  Epetra_SerialComm Comm;
386 #endif
387 
388  verbose = (Comm.MyPID() == 0);
389 
390  for (int i = 1 ; i < argc ; ++i) {
391  if (strcmp(argv[i],"-s") == 0) {
392  SymmetricGallery = true;
393  Solver = AZ_cg;
394  }
395  }
396 
397  // size of the global matrix.
398  Teuchos::ParameterList GaleriList;
399  int nx = 30;
400  GaleriList.set("nx", nx);
401  GaleriList.set("ny", nx * Comm.NumProc());
402  GaleriList.set("mx", 1);
403  GaleriList.set("my", Comm.NumProc());
404  Teuchos::RefCountPtr<Epetra_Map> Map = Teuchos::rcp( Galeri::CreateMap64("Cartesian2D", Comm, GaleriList) );
405  Teuchos::RefCountPtr<Epetra_CrsMatrix> A;
406  if (SymmetricGallery)
407  A = Teuchos::rcp( Galeri::CreateCrsMatrix("Laplace2D", &*Map, GaleriList) );
408  else
409  A = Teuchos::rcp( Galeri::CreateCrsMatrix("Recirc2D", &*Map, GaleriList) );
410 
411  // test the preconditioner
412  int TestPassed = true;
413 
414  // ======================================== //
415  // first verify that we can get convergence //
416  // with all point relaxation methods //
417  // ======================================== //
418 
419  if(!BasicTest("Jacobi",A,false))
420  TestPassed = false;
421 
422  if(!BasicTest("symmetric Gauss-Seidel",A,false))
423  TestPassed = false;
424 
425  if (!SymmetricGallery) {
426  if(!BasicTest("Gauss-Seidel",A,false))
427  TestPassed = false;
428  if(!BasicTest("Gauss-Seidel",A,true))
429  TestPassed = false;
430 
431  }
432 
433  // ============================= //
434  // check uses as preconditioners //
435  // ============================= //
436 
437  if(!KrylovTest("symmetric Gauss-Seidel",A,false))
438  TestPassed = false;
439 
440  if (!SymmetricGallery) {
441  if(!KrylovTest("Gauss-Seidel",A,false))
442  TestPassed = false;
443  if(!KrylovTest("Gauss-Seidel",A,true))
444  TestPassed = false;
445 
446  }
447 
448  // ================================== //
449  // compare point and block relaxation //
450  // ================================== //
451 
452  //TestPassed = TestPassed &&
453  // ComparePointAndBlock("Jacobi",A,1);
454 
455  TestPassed = TestPassed &&
456  ComparePointAndBlock("Jacobi",A,10);
457 
458  //TestPassed = TestPassed &&
459  //ComparePointAndBlock("symmetric Gauss-Seidel",A,1);
460 
461  TestPassed = TestPassed &&
462  ComparePointAndBlock("symmetric Gauss-Seidel",A,10);
463 
464  if (!SymmetricGallery) {
465  //TestPassed = TestPassed &&
466  //ComparePointAndBlock("Gauss-Seidel",A,1);
467 
468  TestPassed = TestPassed &&
469  ComparePointAndBlock("Gauss-Seidel",A,10);
470  }
471 
472  // ============================ //
473  // verify effect of # of blocks //
474  // ============================ //
475 
476  {
477  int Iters4, Iters8, Iters16;
478  Iters4 = CompareBlockSizes("Jacobi",A,4);
479  Iters8 = CompareBlockSizes("Jacobi",A,8);
480  Iters16 = CompareBlockSizes("Jacobi",A,16);
481 
482  if ((Iters16 > Iters8) && (Iters8 > Iters4)) {
483  if (verbose)
484  cout << "Test passed" << endl;
485  }
486  else {
487  if (verbose)
488  cout << "TEST FAILED!" << endl;
489  TestPassed = TestPassed && false;
490  }
491  }
492 
493  // ================================== //
494  // verify effect of overlap in Jacobi //
495  // ================================== //
496 
497  {
498  int Iters0, Iters2, Iters4;
499  Iters0 = CompareBlockOverlap(A,0);
500  Iters2 = CompareBlockOverlap(A,2);
501  Iters4 = CompareBlockOverlap(A,4);
502  if ((Iters4 < Iters2) && (Iters2 < Iters0)) {
503  if (verbose)
504  cout << "Test passed" << endl;
505  }
506  else {
507  if (verbose)
508  cout << "TEST FAILED!" << endl;
509  TestPassed = TestPassed && false;
510  }
511  }
512 
513  // ============ //
514  // final output //
515  // ============ //
516 
517  if (!TestPassed) {
518  cout << "Test `TestRelaxation.exe' failed!" << endl;
519  exit(EXIT_FAILURE);
520  }
521 
522 #ifdef HAVE_MPI
523  MPI_Finalize();
524 #endif
525 
526  cout << endl;
527  cout << "Test `TestRelaxation.exe' passed!" << endl;
528  cout << endl;
529  return(EXIT_SUCCESS);
530 }
Ifpack_BlockRelaxation: a class to define block relaxation preconditioners of Epetra_RowMatrix&#39;s.
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
bool BasicTest(std::string PrecType, const Teuchos::RefCountPtr< Epetra_RowMatrix > &A, bool backward, bool reorder=false)
bool KrylovTest(std::string PrecType, const Teuchos::RefCountPtr< Epetra_RowMatrix > &A, bool backward, bool reorder=false)
const int NumVectors
Definition: performance.cpp:71
int MyPID() const
virtual int MyPID() const =0
virtual int SetParameters(Teuchos::ParameterList &List)
Sets all the parameters for the preconditioner.
static bool Solver
Definition: performance.cpp:70
virtual const Epetra_Comm & Comm() const =0
virtual int Compute()
Computes the preconditioner.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Ifpack_PointRelaxation: a class to define point relaxation preconditioners of for Epetra_RowMatrix&#39;s...
virtual int SetParameters(Teuchos::ParameterList &List)
Sets all the parameters for the preconditioner.
int main(int argc, char *argv[])
virtual int Compute()
Computes the preconditioners.
int NumProc() const
bool ComparePointAndBlock(std::string PrecType, const Teuchos::RefCountPtr< Epetra_RowMatrix > &A, int sweeps)
Definition: performance.cpp:74
Epetra_RowMatrix * GetMatrix() const
static bool SymmetricGallery
int CompareBlockOverlap(const Teuchos::RefCountPtr< Epetra_RowMatrix > &A, int Overlap)
virtual int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Applies the preconditioner to X, returns the result in Y.
int CompareBlockSizes(std::string PrecType, const Teuchos::RefCountPtr< Epetra_RowMatrix > &A, int NumParts)
#define RHS(a)
Definition: MatGenFD.c:60