#include "Didasko_ConfigDefs.h"
#if defined(HAVE_DIDASKO_EPETRA) && defined(HAVE_DIDASKO_ML) && defined(HAVE_DIDASKO_TRIUTILS) && defined(HAVE_DIDASKO_AZTECOO) && defined(HAVE_DIDASKO_TEUCHOS)
#include "Epetra_ConfigDefs.h"
#ifdef HAVE_MPI
#include "mpi.h"
#include "Epetra_MpiComm.h"
#else
#include "Epetra_SerialComm.h"
#endif
#include "Epetra_Map.h"
#include "Epetra_Vector.h"
#include "Epetra_RowMatrix.h"
#include "Epetra_CrsMatrix.h"
#include "Epetra_LinearProblem.h"
#include "Epetra_Time.h"
#include "AztecOO.h"
#include "ml_epetra_preconditioner.h"
#include "Trilinos_Util_CrsMatrixGallery.h"
using namespace Teuchos;
using namespace Trilinos_Util;
#include <iostream>
int main(int argc, char *argv[])
{
#ifdef EPETRA_MPI
MPI_Init(&argc,&argv);
#else
#endif
CrsMatrixGallery Gallery("laplace_3d", Comm, false);
Gallery.Set("problem_size", 1000);
AztecOO solver(*Problem);
ML_Epetra::MultiLevelPreconditioner * MLPrec =
new ML_Epetra::MultiLevelPreconditioner(*A, true);
solver.SetPrecOperator(MLPrec);
solver.SetAztecOption(AZ_solver, AZ_gmres_condnum);
solver.SetAztecOption(AZ_output, 32);
solver.SetAztecOption(AZ_kspace, 160);
int Niters = 500;
solver.Iterate(Niters, 1e-12);
if( Comm.
MyPID() == 0 ) cout << MLPrec->GetOutputList();
delete MLPrec;
double residual, diff;
Gallery.ComputeResidual(&residual);
Gallery.ComputeDiffBetweenStartingAndExactSolutions(&diff);
cout << "||b-Ax||_2 = " << residual << endl;
cout << "||x_exact - x||_2 = " << diff << endl;
cout << "Total Time = " << Time.ElapsedTime() << endl;
}
if (residual > 1e-5)
exit(EXIT_FAILURE);
#ifdef EPETRA_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-teuchos\n"
"--enable-triutils\n"
"--enable-aztecoo\n"
"--enable-ml");
#ifdef HAVE_MPI
MPI_Finalize();
#endif
return(EXIT_SUCCESS);
}
#endif