#include "ml_common.h"
#if defined(HAVE_ML_MLAPI) && defined(HAVE_ML_GALERI) && defined(HAVE_ML_AZTECOO)
#include "Epetra_Map.h"
#include "Epetra_Vector.h"
#include "Epetra_CrsMatrix.h"
#include "Epetra_LinearProblem.h"
#include "Galeri_Maps.h"
#include "Galeri_CrsMatrices.h"
#include "AztecOO.h"
using namespace Teuchos;
using namespace MLAPI;
using namespace Galeri;
int main(int argc, char *argv[])
{
#ifdef HAVE_MPI
  MPI_Init(&argc,&argv);
#else
#endif
  
  int nx = 8;
  GaleriList.
set(
"nx", nx);
  GaleriList.
set(
"ny", ny);
  Epetra_Map* Map = CreateMap(
"Cartesian2D", Comm, GaleriList);
 
  AztecOO solver(Problem);
  try {
    MLList.
set(
"max levels",3);
    MLList.
set(
"increasing or decreasing",
"increasing");
    MLList.
set(
"aggregation: type", 
"Uncoupled");
    MLList.
set(
"aggregation: damping factor", 0.0);
    MLList.
set(
"smoother: type",
"symmetric Gauss-Seidel");
    MLList.
set(
"smoother: sweeps",1);
    MLList.
set(
"smoother: damping factor",1.0);
    MLList.
set(
"coarse: max size",3);
    MLList.
set(
"smoother: pre or post", 
"both");
    MLList.
set(
"coarse: type",
"Amesos-KLU");
    solver.SetPrecOperator(Prec);
    solver.SetAztecOption(AZ_solver, AZ_cg_condnum);
    solver.SetAztecOption(AZ_output, 32);
    solver.Iterate(500, 1e-12);
    
    delete Prec;
  }
  catch (const int e) {
    cerr << "Caught exception, code = " << e << endl;
  }
  catch (...) {
    cerr << "Caught exception..." << endl;
  }
  
  double residual;
  LHS.Norm2(&residual);
    cout << "||b-Ax||_2 = " << residual << endl;
  }
  delete A;
  delete Map;
  
  if (residual > 1e-5)
    exit(EXIT_FAILURE);
#ifdef HAVE_MPI
  MPI_Finalize();
#endif
  return(EXIT_SUCCESS);
}
#else
#include "ml_include.h"
int main(int argc, char *argv[])
{
#ifdef HAVE_MPI
  MPI_Init(&argc,&argv);
#endif
  puts("This MLAPI example requires the following configuration options:");
  puts("\t--enable-epetra");
  puts("\t--enable-teuchos");
  puts("\t--enable-ifpack");
  puts("\t--enable-amesos");
  puts("\t--enable-galeri");
  puts("Please check your configure line.");
#ifdef HAVE_MPI
  MPI_Finalize();
#endif
  return(EXIT_SUCCESS);
}
#endif