// ====================================================================== 
bool BasicTest(string PrecType, const Teuchos::RefCountPtr<Epetra_RowMatrix>& A,bool backward, bool reorder=false)
{
  Epetra_MultiVector LHS(A->RowMatrixRowMap(), NumVectors);
  Epetra_MultiVector RHS(A->RowMatrixRowMap(), NumVectors);
  LHS.PutScalar(0.0); RHS.Random();

  double starting_residual = Galeri::ComputeNorm(&*A, &LHS, &RHS);
  Epetra_LinearProblem Problem(&*A, &LHS, &RHS);

  // Set up the list
  Teuchos::ParameterList List;
  List.set("relaxation: damping factor", 1.0);
  List.set("relaxation: sweeps",2550);
  List.set("relaxation: type", PrecType);
  if(backward) List.set("relaxation: backward mode",backward);

  // Reordering if needed
  int NumRows=A->NumMyRows();
  std::vector<int> RowList(NumRows);
  if(reorder) {
    for(int i=0; i<NumRows; i++)
      RowList[i]=i;
    List.set("relaxation: number of local smoothing indices",NumRows);
    List.set("relaxation: local smoothing indices",RowList.size()>0? &RowList[0] : (int*)0);
  }

  Ifpack_PointRelaxation Point(&*A);

  Point.SetParameters(List);
  Point.Compute();
  // use the preconditioner as solver, with 1550 iterations
  Point.ApplyInverse(RHS,LHS);

  // compute the real residual

  double residual = Galeri::ComputeNorm(&*A, &LHS, &RHS);
  
  if (A->Comm().MyPID() == 0 && verbose)
    cout << "||A * x - b||_2 (scaled) = " << residual / starting_residual << endl;
  
  // Jacobi is very slow to converge here
  if (residual / starting_residual < 1e-2) {
    if (verbose)
      cout << "BasicTest Test passed" << endl;
    return(true);
  }
  else {
    if (verbose)
      cout << "BasicTest Test failed!" << endl;
    return(false);
  }
}
// ====================================================================== 
int main(int argc, char *argv[])
{
#ifdef HAVE_MPI
  MPI_Init(&argc,&argv);
  Epetra_MpiComm Comm( MPI_COMM_WORLD );
#else
  Epetra_SerialComm Comm;
#endif

  verbose = (Comm.MyPID() == 0);

  for (int i = 1 ; i < argc ; ++i) {
    if (strcmp(argv[i],"-s") == 0) {
      SymmetricGallery = true;
      Solver = AZ_cg;
    }
  }

  // size of the global matrix. 
  Teuchos::ParameterList GaleriList;
  int nx = 30; 
  GaleriList.set("nx", nx);
  GaleriList.set("ny", nx * Comm.NumProc());
  GaleriList.set("mx", 1);
  GaleriList.set("my", Comm.NumProc());
  Teuchos::RefCountPtr<Epetra_Map> Map = Teuchos::rcp( Galeri::CreateMap("Cartesian2D", Comm, GaleriList) );
  Teuchos::RefCountPtr<Epetra_CrsMatrix> A;
  if (SymmetricGallery)
    A = Teuchos::rcp( Galeri::CreateCrsMatrix("Laplace2D", &*Map, GaleriList) );
  else
    A = Teuchos::rcp( Galeri::CreateCrsMatrix("Recirc2D", &*Map, GaleriList) );

  // coordinates
  Teuchos::RCP<Epetra_MultiVector> coord = Teuchos::rcp( Galeri::CreateCartesianCoordinates("2D",&*Map,GaleriList));

  // test the preconditioner
  int TestPassed = true;

  // ======================================== //
  // first verify that we can get convergence //
  // with all point relaxation methods        //
  // ======================================== //

  if(!BasicTest("Jacobi",A,false))
    TestPassed = false;

  if(!BasicTest("symmetric Gauss-Seidel",A,false))
    TestPassed = false;

  if(!BasicTest("symmetric Gauss-Seidel",A,false,true))
    TestPassed = false;

  if (!SymmetricGallery) {
    if(!BasicTest("Gauss-Seidel",A,false))
      TestPassed = false;
    if(!BasicTest("Gauss-Seidel",A,true))
      TestPassed = false;  

    if(!BasicTest("Gauss-Seidel",A,false,true))
      TestPassed = false;
    if(!BasicTest("Gauss-Seidel",A,true,true))
      TestPassed = false;  

  }

  // ============================= //
  // check uses as preconditioners //
  // ============================= //
  
  if(!KrylovTest("symmetric Gauss-Seidel",A,false))
    TestPassed = false;

  if(!KrylovTest("symmetric Gauss-Seidel",A,false,true))
    TestPassed = false;


  if (!SymmetricGallery) {
    if(!KrylovTest("Gauss-Seidel",A,false))
      TestPassed = false;
    if(!KrylovTest("Gauss-Seidel",A,true))
      TestPassed = false;

    if(!KrylovTest("Gauss-Seidel",A,false,true))
      TestPassed = false;
    if(!KrylovTest("Gauss-Seidel",A,true,true))
      TestPassed = false;

  }

  // ================================== //
  // compare point and block relaxation //
  // ================================== //

  //TestPassed = TestPassed && 
   // ComparePointAndBlock("Jacobi",A,1);

  TestPassed = TestPassed && 
    ComparePointAndBlock("Jacobi",A,10);

  //TestPassed = TestPassed && 
    //ComparePointAndBlock("symmetric Gauss-Seidel",A,1);

  TestPassed = TestPassed && 
    ComparePointAndBlock("symmetric Gauss-Seidel",A,10);

  if (!SymmetricGallery) {
    //TestPassed = TestPassed && 
      //ComparePointAndBlock("Gauss-Seidel",A,1);

    TestPassed = TestPassed && 
      ComparePointAndBlock("Gauss-Seidel",A,10);
  }

  // ============================ //
  // verify effect of # of blocks //
  // ============================ //
  
  {
    int Iters4, Iters8, Iters16;
    Iters4 = CompareBlockSizes("Jacobi",A,4);
    Iters8 = CompareBlockSizes("Jacobi",A,8);
    Iters16 = CompareBlockSizes("Jacobi",A,16);

    if ((Iters16 > Iters8) && (Iters8 > Iters4)) {
      if (verbose)
        cout << "CompareBlockSizes Test passed" << endl;
    }
    else {
      if (verbose) 
        cout << "CompareBlockSizes TEST FAILED!" << endl;
      TestPassed = TestPassed && false;
    }
  }

  // ================================== //
  // verify effect of overlap in Jacobi //
  // ================================== //

  {
    int Iters0, Iters2, Iters4;
    Iters0 = CompareBlockOverlap(A,0);
    Iters2 = CompareBlockOverlap(A,2);
    Iters4 = CompareBlockOverlap(A,4);
    if ((Iters4 < Iters2) && (Iters2 < Iters0)) {
      if (verbose)
        cout << "CompareBlockOverlap Test passed" << endl;
    }
    else {
      if (verbose) 
        cout << "CompareBlockOverlap TEST FAILED!" << endl;
      TestPassed = TestPassed && false;
    }
  }

  // ================================== //
  // check if line smoothing works      //
  // ================================== //
  {
    int Iters1=
    CompareLineSmoother(A,coord);    
    printf(" comparelinesmoother iters %d \n",Iters1);
  }				
 // ================================== //
  // check if All singleton version of CompareLineSmoother    //
  // ================================== //
  {

    AllSingle(A,coord);    

  }				

  // ================================== //
  // test variable blocking             //
  // ================================== //
  {
    TestPassed = TestPassed && TestVariableBlocking(A->Comm());
  }

  // ================================== //
  // test variable blocking             //
  // ================================== //
  {
    TestPassed = TestPassed && TestTriDiVariableBlocking(A->Comm());
  }


  // ============ //
  // final output //
  // ============ //

  if (!TestPassed) {
    cout << "Test `TestRelaxation.exe' failed!" << endl;
    exit(EXIT_FAILURE);
  }
  
#ifdef HAVE_MPI
  MPI_Finalize(); 
#endif

  cout << endl;
  cout << "Test `TestRelaxation.exe' passed!" << endl;
  cout << endl;
  return(EXIT_SUCCESS);
}
Esempio n. 3
0
// ======================================================================
bool TestContainer(std::string Type, const Teuchos::RefCountPtr<Epetra_RowMatrix>& A)
{
  using std::cout;
  using std::endl;

  int NumVectors = 3;
  int NumMyRows = A->NumMyRows();

  Epetra_MultiVector LHS_exact(A->RowMatrixRowMap(), NumVectors);
  Epetra_MultiVector LHS(A->RowMatrixRowMap(), NumVectors);
  Epetra_MultiVector RHS(A->RowMatrixRowMap(), NumVectors);
  LHS_exact.Random(); LHS.PutScalar(0.0);
  A->Multiply(false, LHS_exact, RHS);

  Epetra_LinearProblem Problem(&*A, &LHS, &RHS);

  if (verbose) {
    cout << "Container type = " << Type << endl;
    cout << "NumMyRows = " << NumMyRows << ", NumVectors = " << NumVectors << endl;
  }
  LHS.PutScalar(0.0);

  Teuchos::RefCountPtr<Ifpack_Container> Container;

  if (Type == "dense")
    Container = Teuchos::rcp( new Ifpack_DenseContainer(A->NumMyRows(), NumVectors) );
  else
    Container = Teuchos::rcp( new Ifpack_SparseContainer<Ifpack_Amesos>(A->NumMyRows(), NumVectors) );

  assert (Container != Teuchos::null);

  IFPACK_CHK_ERR(Container->Initialize());
  // set as ID all the local rows of A
  for (int i = 0 ; i < A->NumMyRows() ; ++i)
    Container->ID(i) = i;

  // extract submatrix (in this case, the entire matrix)
  // and complete setup
  IFPACK_CHK_ERR(Container->Compute(*A));

  // set the RHS and LHS
  for (int i = 0 ; i < A->NumMyRows() ; ++i)
    for (int j = 0 ; j < NumVectors ; ++j) {
      Container->RHS(i,j) = RHS[j][i];
      Container->LHS(i,j) = LHS[j][i];
    }

  // set parameters (empty for dense containers)
  Teuchos::ParameterList List;
  List.set("amesos: solver type", Type);
  IFPACK_CHK_ERR(Container->SetParameters(List));

  // solve the linear system
  IFPACK_CHK_ERR(Container->ApplyInverse());

  // get the computed solution, store it in LHS
  for (int i = 0 ; i < A->NumMyRows() ; ++i)
    for (int j = 0 ; j < NumVectors ; ++j) {
       LHS[j][i] = Container->LHS(i,j);
    }

  double residual = Galeri::ComputeNorm(&LHS, &LHS_exact);

  if (A->Comm().MyPID() == 0 && verbose) {
    cout << "||x_exact - x||_2 = " << residual << endl;
    cout << *Container;
  }

  bool passed = false;
  if (residual < 1e-5)
    passed = true;

  return(passed);
}
Esempio n. 4
0
bool CompareWithAztecOO(Epetra_LinearProblem& Problem, const std::string what,
                       int Overlap, int ival)
{
  using std::cout;
  using std::endl;

  AztecOO AztecOOSolver(Problem);
  AztecOOSolver.SetAztecOption(AZ_solver,AZ_gmres);
  AztecOOSolver.SetAztecOption(AZ_output,AZ_none);
  AztecOOSolver.SetAztecOption(AZ_overlap,Overlap);
  AztecOOSolver.SetAztecOption(AZ_graph_fill,ival);
  AztecOOSolver.SetAztecOption(AZ_poly_ord, ival);
  AztecOOSolver.SetAztecParam(AZ_drop, 0.0);
  AztecOOSolver.SetAztecParam(AZ_athresh, 0.0);
  AztecOOSolver.SetAztecParam(AZ_rthresh, 0.0);

  Epetra_MultiVector& RHS = *(Problem.GetRHS());
  Epetra_MultiVector& LHS = *(Problem.GetLHS());
  Teuchos::RefCountPtr<Epetra_RowMatrix> A = Teuchos::rcp(Problem.GetMatrix(), false);

  LHS.Random();
  A->Multiply(false,LHS,RHS);

  Teuchos::ParameterList List;
  List.set("fact: level-of-fill", ival);
  List.set("relaxation: sweeps", ival);
  List.set("relaxation: damping factor", 1.0);
  List.set("relaxation: zero starting solution", true);

  //default combine mode is as for AztecOO
  List.set("schwarz: combine mode", Zero);

  Epetra_Time Time(A->Comm());

  Teuchos::RefCountPtr<Ifpack_Preconditioner> Prec;

  if (what == "Jacobi") {
    Prec = Teuchos::rcp( new Ifpack_PointRelaxation(&*A) );
    List.set("relaxation: type", "Jacobi");
    AztecOOSolver.SetAztecOption(AZ_precond,AZ_Jacobi);
    AztecOOSolver.SetAztecOption(AZ_reorder,0);
  }
  else if (what == "IC no reord") {
    Prec = Teuchos::rcp( new Ifpack_AdditiveSchwarz<Ifpack_IC>(&*A,Overlap) );
    AztecOOSolver.SetAztecOption(AZ_precond,AZ_dom_decomp);
    AztecOOSolver.SetAztecOption(AZ_subdomain_solve,AZ_icc);
    AztecOOSolver.SetAztecOption(AZ_reorder,0);
  }
  else if (what == "IC reord") {
    Prec = Teuchos::rcp( new Ifpack_AdditiveSchwarz<Ifpack_IC>(&*A,Overlap) );
    List.set("schwarz: use reordering", true);
    AztecOOSolver.SetAztecOption(AZ_precond,AZ_dom_decomp);
    AztecOOSolver.SetAztecOption(AZ_subdomain_solve,AZ_icc);
    AztecOOSolver.SetAztecOption(AZ_reorder,1);
  }
  else if (what == "ILU no reord") {
    Prec = Teuchos::rcp( new Ifpack_AdditiveSchwarz<Ifpack_ILU>(&*A,Overlap) );
    AztecOOSolver.SetAztecOption(AZ_precond,AZ_dom_decomp);
    AztecOOSolver.SetAztecOption(AZ_subdomain_solve,AZ_ilu);
    AztecOOSolver.SetAztecOption(AZ_reorder,0);
  }
  else if (what == "ILU reord") {
    Prec = Teuchos::rcp( new Ifpack_AdditiveSchwarz<Ifpack_ILU>(&*A,Overlap) );
    List.set("schwarz: use reordering", true);
    AztecOOSolver.SetAztecOption(AZ_precond,AZ_dom_decomp);
    AztecOOSolver.SetAztecOption(AZ_subdomain_solve,AZ_ilu);
    AztecOOSolver.SetAztecOption(AZ_reorder,1);
  }
#ifdef HAVE_IFPACK_AMESOS
  else if (what == "LU") {
    Prec = Teuchos::rcp( new Ifpack_AdditiveSchwarz<Ifpack_Amesos>(&*A,Overlap) );
    List.set("amesos: solver type", "Klu");
    AztecOOSolver.SetAztecOption(AZ_precond,AZ_dom_decomp);
    AztecOOSolver.SetAztecOption(AZ_subdomain_solve,AZ_lu);
  }
#endif
  else {
    cerr << "Option not recognized" << endl;
    exit(EXIT_FAILURE);
  }

  // ==================================== //
  // Solve with AztecOO's preconditioners //
  // ==================================== //

  LHS.PutScalar(0.0);

  Time.ResetStartTime();
  AztecOOSolver.Iterate(150,1e-5);

  if (verbose) {
    cout << endl;
    cout << "==================================================" << endl;
    cout << "Testing `" << what << "', Overlap = "
         << Overlap << ", ival = " << ival << endl;
    cout << endl;
    cout << "[AztecOO] Total time = " << Time.ElapsedTime() << " (s)" << endl;
    cout << "[AztecOO] Residual   = " << AztecOOSolver.TrueResidual() << " (s)" << endl;
    cout << "[AztecOO] Iterations = " << AztecOOSolver.NumIters() << endl;
    cout << endl;
  }

  int AztecOOPrecIters = AztecOOSolver.NumIters();

  // =========================================== //
  // Create the IFPACK preconditioner and solver //
  // =========================================== //

  Epetra_Time Time2(A->Comm());
  assert(Prec != Teuchos::null);
  IFPACK_CHK_ERR(Prec->SetParameters(List));

  Time.ResetStartTime();
  IFPACK_CHK_ERR(Prec->Initialize());
  if (verbose)
    cout << "[IFPACK] Time for Initialize() = "
         << Time.ElapsedTime() << " (s)" << endl;

  Time.ResetStartTime();
  IFPACK_CHK_ERR(Prec->Compute());
  if (verbose)
    cout << "[IFPACK] Time for Compute() = "
         << Time.ElapsedTime() << " (s)" << endl;


  AztecOOSolver.SetPrecOperator(&*Prec);

  LHS.PutScalar(0.0);

  Time.ResetStartTime();
  AztecOOSolver.Iterate(150,1e-5);

  if (verbose) {
    cout << "[IFPACK] Total time = " << Time2.ElapsedTime() << " (s)" << endl;
    cout << "[IFPACK] Residual   = " << AztecOOSolver.TrueResidual() << " (s)" << endl;
    cout << "[IFPACK] Iterations = " << AztecOOSolver.NumIters() << endl;
    cout << endl;
  }

  int IFPACKPrecIters = AztecOOSolver.NumIters();

  if (IFPACK_ABS(AztecOOPrecIters - IFPACKPrecIters) > 3) {
    cerr << "TEST FAILED (" << AztecOOPrecIters << " != "
         << IFPACKPrecIters << ")" << endl;
    return(false);
  }
  else
    return(true);

}