Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/PointPreconditioner/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 #ifdef HAVE_IFPACK_AZTECOO
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 "Trilinos_Util_CrsMatrixGallery.h"
54 #include "Teuchos_ParameterList.hpp"
55 #include "Ifpack_PointRelaxation.h"
56 #include "Ifpack_BlockRelaxation.h"
57 #include "Ifpack_DenseContainer.h"
58 #include "Ifpack_AdditiveSchwarz.h"
59 #include "AztecOO.h"
60 
61 using namespace Trilinos_Util;
62 static bool verbose = false;
63 static bool SymmetricGallery = false;
64 static bool Solver = AZ_gmres;
65 
66 
67 // ======================================================================
68 int CompareBlockOverlap(CrsMatrixGallery& Gallery, int Overlap)
69 {
70 
71  Epetra_RowMatrix* A = Gallery.GetMatrix();
72  Epetra_LinearProblem* Problem = Gallery.GetLinearProblem();
73 
74  Epetra_MultiVector& RHS = *(Problem->GetRHS());
75  Epetra_MultiVector& LHS = *(Problem->GetLHS());
76 
78  List.set("relaxation: damping factor", 1.0);
79  List.set("relaxation: type", "Jacobi");
80  List.set("relaxation: sweeps",1);
81  List.set("partitioner: overlap", Overlap);
82  List.set("partitioner: type", "greedy");
83  List.set("partitioner: local parts", 16);
84 
85  RHS.PutScalar(1.0);
86  LHS.PutScalar(0.0);
87 
89  Prec.SetParameters(List);
90  Prec.Compute();
91 
92  // set AztecOO solver object
93  AztecOO AztecOOSolver(*Problem);
94  AztecOOSolver.SetAztecOption(AZ_solver,Solver);
95  if (verbose)
96  AztecOOSolver.SetAztecOption(AZ_output,32);
97  else
98  AztecOOSolver.SetAztecOption(AZ_output,AZ_none);
99  AztecOOSolver.SetPrecOperator(&Prec);
100 
101  AztecOOSolver.Iterate(1550,1e-8);
102 
103  return(AztecOOSolver.NumIters());
104 }
105 
106 // ======================================================================
107 int CompareBlockSizes(std::string PrecType,
108  CrsMatrixGallery& Gallery, int NumParts)
109 {
110 
111  Epetra_RowMatrix* A = Gallery.GetMatrix();
112  Epetra_LinearProblem* Problem = Gallery.GetLinearProblem();
113 
114  Epetra_MultiVector& RHS = *(Problem->GetRHS());
115  Epetra_MultiVector& LHS = *(Problem->GetLHS());
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", "greedy");
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-8);
141 
142  return(AztecOOSolver.NumIters());
143 }
144 
145 // ======================================================================
146 bool ComparePointAndBlock(std::string PrecType,
147  CrsMatrixGallery& Gallery, int sweeps)
148 {
149 
150  Epetra_RowMatrix* A = Gallery.GetMatrix();
151  Epetra_LinearProblem* Problem = Gallery.GetLinearProblem();
152 
153  Epetra_MultiVector& RHS = *(Problem->GetRHS());
154  Epetra_MultiVector& LHS = *(Problem->GetLHS());
155 
156  // Set up the list
158  List.set("relaxation: damping factor", 1.0);
159  List.set("relaxation: type", PrecType);
160  List.set("relaxation: sweeps",sweeps);
161  List.set("partitioner: type", "linear");
162  List.set("partitioner: local parts", A->NumMyRows());
163 
164  int ItersPoint, ItersBlock;
165 
166  // ================================================== //
167  // get the number of iterations with point relaxation //
168  // ================================================== //
169  {
170 
171  RHS.PutScalar(1.0);
172  LHS.PutScalar(0.0);
173 
175  Point.SetParameters(List);
176  Point.Compute();
177 
178  // set AztecOO solver object
179  AztecOO AztecOOSolver(*Problem);
180  AztecOOSolver.SetAztecOption(AZ_solver,Solver);
181  if (verbose)
182  AztecOOSolver.SetAztecOption(AZ_output,32);
183  else
184  AztecOOSolver.SetAztecOption(AZ_output,AZ_none);
185  AztecOOSolver.SetPrecOperator(&Point);
186 
187  AztecOOSolver.Iterate(1550,1e-8);
188 
189  double TrueResidual = AztecOOSolver.TrueResidual();
190  ItersPoint = AztecOOSolver.NumIters();
191  // some output
192  if (verbose && Problem->GetMatrix()->Comm().MyPID() == 0) {
193  cout << "Iterations = " << ItersPoint << endl;
194  cout << "Norm of the true residual = " << TrueResidual << endl;
195  }
196  }
197 
198  // ================================================== //
199  // get the number of iterations with block relaxation //
200  // ================================================== //
201  {
202 
203  RHS.PutScalar(1.0);
204  LHS.PutScalar(0.0);
205 
207  Block.SetParameters(List);
208  Block.Compute();
209 
210  // set AztecOO solver object
211  AztecOO AztecOOSolver(*Problem);
212  AztecOOSolver.SetAztecOption(AZ_solver,Solver);
213  if (verbose)
214  AztecOOSolver.SetAztecOption(AZ_output,32);
215  else
216  AztecOOSolver.SetAztecOption(AZ_output,AZ_none);
217  AztecOOSolver.SetPrecOperator(&Block);
218 
219  AztecOOSolver.Iterate(1550,1e-8);
220 
221  double TrueResidual = AztecOOSolver.TrueResidual();
222  ItersBlock = AztecOOSolver.NumIters();
223  // some output
224  if (verbose && Problem->GetMatrix()->Comm().MyPID() == 0) {
225  cout << "Iterations " << ItersBlock << endl;
226  cout << "Norm of the true residual = " << TrueResidual << endl;
227  }
228  }
229 
230  if (ItersPoint != ItersBlock) {
231  if (verbose)
232  cout << "TEST FAILED!" << endl;
233  return(false);
234  }
235  else {
236  if (verbose)
237  cout << "TEST PASSED" << endl;
238  return(true);
239  }
240 }
241 
242 // ======================================================================
243 bool KrylovTest(std::string PrecType,
244  CrsMatrixGallery& Gallery)
245 {
246 
247  // The following methods of CrsMatrixGallery are used to get pointers
248  // to internally stored Epetra_RowMatrix and Epetra_LinearProblem.
249  Epetra_RowMatrix* A = Gallery.GetMatrix();
250  Epetra_LinearProblem* Problem = Gallery.GetLinearProblem();
251 
252  Epetra_MultiVector& LHS = *(Problem->GetLHS());
253  LHS.PutScalar(0.0);
254 
255  // Set up the list
257  List.set("relaxation: damping factor", 1.0);
258  List.set("relaxation: type", PrecType);
259 
260  int Iters1, Iters10;
261 
262  // ============================================== //
263  // get the number of iterations with 1 sweep only //
264  // ============================================== //
265  {
266 
267  List.set("relaxation: sweeps",1);
269  Point.SetParameters(List);
270  Point.Compute();
271 
272  // set AztecOO solver object
273  AztecOO AztecOOSolver(*Problem);
274  AztecOOSolver.SetAztecOption(AZ_solver,Solver);
275  AztecOOSolver.SetAztecOption(AZ_output,AZ_none);
276  AztecOOSolver.SetPrecOperator(&Point);
277 
278  AztecOOSolver.Iterate(1550,1e-8);
279 
280  double TrueResidual = AztecOOSolver.TrueResidual();
281  // some output
282  if (verbose && Problem->GetMatrix()->Comm().MyPID() == 0) {
283  cout << "Norm of the true residual = " << TrueResidual << endl;
284  }
285  Iters1 = AztecOOSolver.NumIters();
286  }
287 
288  // ======================================================== //
289  // now re-run with 10 sweeps, solver should converge faster
290  // ======================================================== //
291  {
292  List.set("relaxation: sweeps",10);
294  Point.SetParameters(List);
295  Point.Compute();
296  LHS.PutScalar(0.0);
297 
298  // set AztecOO solver object
299  AztecOO AztecOOSolver(*Problem);
300  AztecOOSolver.SetAztecOption(AZ_solver,Solver);
301  AztecOOSolver.SetAztecOption(AZ_output,AZ_none);
302  AztecOOSolver.SetPrecOperator(&Point);
303  AztecOOSolver.Iterate(1550,1e-8);
304 
305  double TrueResidual = AztecOOSolver.TrueResidual();
306  // some output
307  if (verbose && Problem->GetMatrix()->Comm().MyPID() == 0) {
308  cout << "Norm of the true residual = " << TrueResidual << endl;
309  }
310  Iters10 = AztecOOSolver.NumIters();
311  }
312 
313  if (verbose) {
314  cout << "Iters_1 = " << Iters1 << ", Iters_10 = " << Iters10 << endl;
315  cout << "(second number should be smaller than first one)" << endl;
316  }
317 
318  if (Iters10 > Iters1) {
319  if (verbose)
320  cout << "TEST FAILED!" << endl;
321  return(false);
322  }
323  else {
324  if (verbose)
325  cout << "TEST PASSED" << endl;
326  return(true);
327  }
328 }
329 
330 // ======================================================================
331 bool BasicTest(std::string PrecType,
332  CrsMatrixGallery& Gallery)
333 {
334 
335  // The following methods of CrsMatrixGallery are used to get pointers
336  // to internally stored Epetra_RowMatrix and Epetra_LinearProblem.
337  Epetra_RowMatrix* A = Gallery.GetMatrix();
338  Epetra_LinearProblem* Problem = Gallery.GetLinearProblem();
339 
340  Epetra_MultiVector& RHS = *(Problem->GetRHS());
341  Epetra_MultiVector& LHS = *(Problem->GetLHS());
342 
343  LHS.PutScalar(0.0);
344 
345  // Set up the list
347  List.set("relaxation: damping factor", 1.0);
348  List.set("relaxation: sweeps",1550);
349  List.set("relaxation: type", PrecType);
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, diff;
361  Gallery.ComputeResidual(&residual);
362  Gallery.ComputeDiffBetweenStartingAndExactSolutions(&diff);
363 
364  if (verbose && A->Comm().MyPID()==0) {
365  cout << "||b-Ax||_2 = " << residual << endl;
366  cout << "||x_exact - x||_2 = " << diff << endl;
367  }
368 
369  // Jacobi is very slow to converge here
370  if (residual < 1e-2) {
371  if (verbose)
372  cout << "Test passed" << endl;
373  return(true);
374  }
375  else {
376  if (verbose)
377  cout << "Test failed!" << endl;
378  return(false);
379  }
380 }
381 
382 // ======================================================================
383 int main(int argc, char *argv[])
384 {
385 
386 #ifdef HAVE_MPI
387  MPI_Init(&argc,&argv);
388  Epetra_MpiComm Comm( MPI_COMM_WORLD );
389 #else
390  Epetra_SerialComm Comm;
391 #endif
392 
393  bool verbose = (Comm.MyPID() == 0);
394 
395  for (int i = 1 ; i < argc ; ++i) {
396  if (strcmp(argv[i],"-s") == 0) {
397  SymmetricGallery = true;
398  Solver = AZ_cg;
399  }
400  }
401 
402  // size of the global matrix.
403  const int NumPoints = 900;
404 
405  CrsMatrixGallery Gallery("", Comm);
406  Gallery.Set("problem_size", NumPoints);
407  Gallery.Set("map_type", "linear");
408  if (SymmetricGallery)
409  Gallery.Set("problem_type", "laplace_2d");
410  else
411  Gallery.Set("problem_type", "recirc_2d");
412 
413  // test the preconditioner
414  int TestPassed = true;
415 
416  // ======================================== //
417  // first verify that we can get convergence //
418  // with all point relaxation methods //
419  // ======================================== //
420  if(!BasicTest("Jacobi",Gallery))
421  TestPassed = false;
422 
423  if(!BasicTest("symmetric Gauss-Seidel",Gallery))
424  TestPassed = false;
425 
426  if (!SymmetricGallery) {
427  if(!BasicTest("Gauss-Seidel",Gallery))
428  TestPassed = false;
429  }
430 
431  // ============================= //
432  // check uses as preconditioners //
433  // ============================= //
434 
435  if(!KrylovTest("symmetric Gauss-Seidel",Gallery))
436  TestPassed = false;
437 
438  if (!SymmetricGallery) {
439  if(!KrylovTest("Gauss-Seidel",Gallery))
440  TestPassed = false;
441  }
442 
443  // ================================== //
444  // compare point and block relaxation //
445  // ================================== //
446 
447  TestPassed = TestPassed &&
448  ComparePointAndBlock("Jacobi",Gallery,1);
449 
450  TestPassed = TestPassed &&
451  ComparePointAndBlock("Jacobi",Gallery,10);
452 
453  TestPassed = TestPassed &&
454  ComparePointAndBlock("symmetric Gauss-Seidel",Gallery,1);
455 
456  TestPassed = TestPassed &&
457  ComparePointAndBlock("symmetric Gauss-Seidel",Gallery,10);
458 
459  if (!SymmetricGallery) {
460  TestPassed = TestPassed &&
461  ComparePointAndBlock("Gauss-Seidel",Gallery,1);
462 
463  TestPassed = TestPassed &&
464  ComparePointAndBlock("Gauss-Seidel",Gallery,10);
465  }
466 
467  // ============================ //
468  // verify effect of # of blocks //
469  // ============================ //
470 
471  {
472  int Iters4, Iters8, Iters16;
473  Iters4 = CompareBlockSizes("Jacobi",Gallery,4);
474  Iters8 = CompareBlockSizes("Jacobi",Gallery,8);
475  Iters16 = CompareBlockSizes("Jacobi",Gallery,16);
476  if ((Iters16 > Iters8) && (Iters8 > Iters4)) {
477  if (verbose)
478  cout << "Test passed" << endl;
479  }
480  else {
481  if (verbose)
482  cout << "TEST FAILED!" << endl;
483  TestPassed = TestPassed && false;
484  }
485  }
486 
487  // ================================== //
488  // verify effect of overlap in Jacobi //
489  // ================================== //
490 
491  {
492  int Iters0, Iters2, Iters4;
493  Iters0 = CompareBlockOverlap(Gallery,0);
494  Iters2 = CompareBlockOverlap(Gallery,2);
495  Iters4 = CompareBlockOverlap(Gallery,4);
496  if ((Iters4 < Iters2) && (Iters2 < Iters0)) {
497  if (verbose)
498  cout << "Test passed" << endl;
499  }
500  else {
501  if (verbose)
502  cout << "TEST FAILED!" << endl;
503  TestPassed = TestPassed && false;
504  }
505  }
506 
507  // ============ //
508  // final output //
509  // ============ //
510 
511  if (!TestPassed) {
512  cout << "Test `PointPreconditioner.exe' FAILED!" << endl;
513  exit(EXIT_FAILURE);
514  }
515 
516 #ifdef HAVE_MPI
517  MPI_Finalize();
518 #endif
519 
520  cout << "Test `PointPreconditioner.exe' passed!" << endl;
521  exit(EXIT_SUCCESS);
522 }
523 
524 #else
525 
526 #ifdef HAVE_MPI
527 #include "Epetra_MpiComm.h"
528 #else
529 #include "Epetra_SerialComm.h"
530 #endif
531 
532 int main(int argc, char *argv[])
533 {
534 
535 #ifdef HAVE_MPI
536  MPI_Init(&argc,&argv);
537  Epetra_MpiComm Comm( MPI_COMM_WORLD );
538 #else
539  Epetra_SerialComm Comm;
540 #endif
541 
542  puts("please configure IFPACK with --eanble-aztecoo --enable-teuchos");
543  puts("to run this test");
544 
545 #ifdef HAVE_MPI
546  MPI_Finalize() ;
547 #endif
548  return(EXIT_SUCCESS);
549 }
550 
551 #endif
Epetra_MultiVector * GetLHS() const
Epetra_MultiVector * GetRHS() const
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)
int MyPID() const
virtual int MyPID() const =0
static bool Solver
Definition: performance.cpp:70
virtual const Epetra_Comm & Comm() const =0
Ifpack_AdditiveSchwarz: a class to define Additive Schwarz preconditioners of Epetra_RowMatrix&#39;s.
Ifpack_PointRelaxation: a class to define point relaxation preconditioners of for Epetra_RowMatrix&#39;s...
virtual int NumMyRows() const =0
int main(int argc, char *argv[])
bool ComparePointAndBlock(std::string PrecType, const Teuchos::RefCountPtr< Epetra_RowMatrix > &A, int sweeps)
Definition: performance.cpp:74
Epetra_RowMatrix * GetMatrix() const
static bool SymmetricGallery
Definition: performance.cpp:69
int CompareBlockOverlap(const Teuchos::RefCountPtr< Epetra_RowMatrix > &A, int Overlap)
int CompareBlockSizes(std::string PrecType, const Teuchos::RefCountPtr< Epetra_RowMatrix > &A, int NumParts)
#define RHS(a)
Definition: MatGenFD.c:60