#include "Didasko_ConfigDefs.h"
#if defined(HAVE_DIDASKO_EPETRA) && defined(HAVE_DIDASKO_TRIUTILS)
#include "Epetra_ConfigDefs.h"
#ifdef HAVE_MPI
#include "mpi.h"
#include "Epetra_MpiComm.h"
#else
#include "Epetra_SerialComm.h"
#endif
#include "AztecOO.h"
#include "Epetra_CrsMatrix.h"
#include "Trilinos_Util_CrsMatrixGallery.h"
using namespace Trilinos_Util;
int main(int argc, char *argv[])
{
#ifdef HAVE_MPI
MPI_Init(&argc,&argv);
#else
#endif
CrsMatrixGallery Gallery("laplace_2d", Comm, false);
Gallery.Set("problem_size",900);
AztecOO solver(*Problem);
solver.SetAztecOption(AZ_solver, AZ_cg_condnum);
solver.SetAztecOption(AZ_precond, AZ_dom_decomp);
solver.SetAztecOption(AZ_overlap,0);
solver.SetAztecOption(AZ_subdomain_solve, AZ_icc);
solver.Iterate(1550, 1e-12);
double status[AZ_STATUS_SIZE];
solver.GetAllAztecStatus(status);
double residual, diff;
Gallery.ComputeResidual(&residual);
Gallery.ComputeDiffBetweenStartingAndExactSolutions(&diff);
cout << "||b-Ax||_2 = " << residual << endl;
cout << "||x_exact - x||_2 = " << diff << endl;
}
#ifdef HAVE_MPI
MPI_Finalize();
#endif
return(EXIT_SUCCESS);
}
#else
#include <stdlib.h>
#include <stdio.h>
#ifdef HAVE_MPI
#include "mpi.h"
#endif
int main(int argc, char *argv[])
{
#ifdef HAVE_MPI
MPI_Init(&argc,&argv);
#endif
puts("Please configure Didasko with:\n"
"--enable-epetra\n"
"--enable-triutils\n"
"--enable-aztecoo\n");
#ifdef HAVE_MPI
MPI_Finalize();
#endif
return 0;
}
#endif