#include "Galeri_ConfigDefs.h"
#include "Galeri_Utils.h"
#include "Galeri_FiniteElements.h"
#ifdef HAVE_MPI
#include "mpi.h"
#include "Epetra_MpiComm.h"
#else
#include "Epetra_SerialComm.h"
#endif
double Diffusion(const double& x, const double& y, const double& z)
{
  return (1.0);
}
double Source(const double& x, const double& y, const double& z)
{
  return (0.0);
}
double Force(const double& x, const double& y, const double& z)
{
  return (-6.0);
}
double BoundaryValue(const double& x, const double& y, 
                     const double& z, const int& PatchID)
{
  return(x * x + y * y + z * z);
}
int BoundaryType(const int& PatchID)
{
  return(Galeri::FiniteElements::GALERI_DIRICHLET);
}
int ExactSolution(double x, double y, double z, double* res)
{
  res[0] = x * x + y * y + z * z;
  res[1] = 2 * x;
  res[2] = 2 * y;
  res[3] = 2 * z;
  return(0);
}
using namespace Galeri;
using namespace Galeri::FiniteElements;
int main(int argc, char *argv[])
{
  
#ifdef HAVE_MPI
  MPI_Init(&argc,&argv);
  Epetra_MpiComm Comm(MPI_COMM_WORLD);
#else
  Epetra_SerialComm Comm;
#endif
  try {
    
    
    
    
    
    
    
    
    
    int nx = 4 * Comm.NumProc();
    int ny = 4; 
    int nz = 4;
    int mx = Comm.NumProc(), my = 1, mz = 1;
    
    
    
    
    
    
    
    Epetra_CrsMatrix A(Copy, Grid.RowMap(), 0);
    Epetra_Vector    LHS(Grid.RowMap());
    Epetra_Vector    RHS(Grid.RowMap());
    int NumQuadratureNodes = 1;
      Laplacian(NumQuadratureNodes, Diffusion, Source, Force, 
                BoundaryValue, BoundaryType);
    LinearProblem FiniteElementProblem(Grid, Laplacian, A, LHS, RHS); 
 
    FiniteElementProblem.Compute();
    
    
    
    
    
    
    
    
    
    
    
    Solve(&A, &LHS, &RHS);
    
    
    
    
    
    
    
    FiniteElementProblem.ComputeNorms(LHS, ExactSolution);
    
    MEDIT.Write(Grid, "Laplacian3D", LHS);
  }
  catch (Exception& rhs)
  {
    if (Comm.MyPID() == 0)
      rhs.Print();
  }
  catch (int e) {
    cerr << "Caught exception, value = " << e << endl;
  }
  catch (...) {
    cerr << "Caught generic exception" << endl;
  }
#ifdef HAVE_MPI
  MPI_Finalize();
#endif
  return(0);
}