#if defined(HAVE_ML_MLAPI) && defined(HAVE_ML_EPETRA) && defined(HAVE_ML_TEUCHOS) && defined(HAVE_ML_IFPACK) && defined(HAVE_ML_AMESOS) && defined(HAVE_ML_AZTECOO) && defined(HAVE_ML_TRIUTILS)
#include "ml_config.h"
#include "ml_common.h"
#include "Epetra_Map.h"
#include "Epetra_Vector.h"
#include "Epetra_CrsMatrix.h"
#include "Epetra_LinearProblem.h"
#include "Trilinos_Util_CrsMatrixGallery.h"
#include "AztecOO.h"
#include "MLAPI_Space.h"
#include "MLAPI_Operator.h"
#include "MLAPI_MultiLevelSA.h"
#include "MLAPI_EpetraBaseOperator.h"
#include "MLAPI_Workspace.h"
using namespace Teuchos;
using namespace MLAPI;
using namespace Trilinos_Util;
int main(int argc, char *argv[])
{
#ifdef HAVE_MPI
MPI_Init(&argc,&argv);
Epetra_MpiComm Comm(MPI_COMM_WORLD);
#else
Epetra_SerialComm Comm;
#endif
CrsMatrixGallery Gallery("laplace_2d", Comm);
Gallery.Set("problem_size", 10000);
Epetra_RowMatrix* A = Gallery.GetMatrix();
Epetra_LinearProblem* Problem = Gallery.GetLinearProblem();
AztecOO solver(*Problem);
Init();
try {
Space S(-1, A->NumMyRows(), A->RowMatrixRowMap().MyGlobalElements());
Operator A_MLAPI(S, S, A, false);
Teuchos::ParameterList MLList;
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");
MultiLevelSA Prec_MLAPI(A_MLAPI, MLList);
Epetra_Operator* Prec = new EpetraBaseOperator(A->RowMatrixRowMap(), Prec_MLAPI);
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;
}
Finalize();
double residual, diff;
Gallery.ComputeResidual(&residual);
Gallery.ComputeDiffBetweenStartingAndExactSolutions(&diff);
if( Comm.MyPID()==0 ) {
cout << "||b-Ax||_2 = " << residual << endl;
cout << "||x_exact - x||_2 = " << diff << endl;
}
if (residual > 1e-5)
exit(EXIT_FAILURE);
#ifdef HAVE_MPI
MPI_Finalize();
#endif
exit(EXIT_SUCCESS);
}
#else
#include "ml_include.h"
int main(int argc, char *argv[])
{
#ifdef HAVE_MPI
MPI_Init(&argc,&argv);
#endif
puts("The ML API requires the following configuration options:");
puts("\t--enable-epetra");
puts("\t--enable-teuchos");
puts("\t--enable-ifpack");
puts("\t--enable-amesos");
puts("Please check your configure line.");
#ifdef HAVE_MPI
MPI_Finalize();
#endif
return(0);
}
#endif // #if defined(HAVE_ML_MLAPI)