#include "ml_config.h"
#include "ml_common.h"
#if defined(HAVE_ML_MLAPI) && defined(HAVE_ML_GALERI)
using namespace Teuchos;
using namespace MLAPI;
int main(int argc, char *argv[])
{
#ifdef HAVE_MPI
  MPI_Init(&argc,&argv);
#endif
  try {
    
    
    int NumGlobalElements = 10000;
    
    Space S(NumGlobalElements);
 
    
    
    
    
    
    MLList.
set(
"max levels",3);
    MLList.
set(
"aggregation: type", 
"Uncoupled");
    MLList.
set(
"aggregation: damping factor", 1.333);
    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(
"coarse: type",
"Amesos-KLU");
    
    
    
    
    x_ex.Random();
    b = A * x_ex;
    x = 0.0;
    double OldNorm   = 1.0;
    double Tolerance = 1e-13;
    int    MaxIters  = 30;
    
    
    
    for (int i = 0 ; i < MaxIters ; ++i) {
      r = b - A * x; 
      z = P * r;     
      x = x + z;     
      
      double NewNorm = sqrt((x - x_ex) * (A * (x - x_ex)));
        std::cout << "||x - x_ex||_A = ";
        std::cout.width(15);
        std::cout << NewNorm << ", ";
        std::cout << "reduction = ";
        std::cout.width(15);
        std::cout << NewNorm / OldNorm << std::endl;
      }
      if (NewNorm < Tolerance)
        break;
      OldNorm = NewNorm;
    }
    
  }
  catch (const int e) {
    std::cout << "Caught integer exception, code = " << e << std::endl;
  }
  catch (...) {
    std::cout << "problems here..." << std::endl;
  }
#ifdef HAVE_MPI
  MPI_Finalize();
#endif
  return(0);
}
#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(0);
}
#endif // if defined(HAVE_ML_MLAPI)