#include "Didasko_ConfigDefs.h"
#if defined(HAVE_DIDASKO_EPETRA) && defined(HAVE_DIDASKO_TEUCHOS) && defined(HAVE_DIDASKO_TRIUTILS)
#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_IntVector.h"
#include "Epetra_SerialDenseVector.h"
#include "Epetra_Vector.h"
#include "Epetra_CrsMatrix.h"
#include "Epetra_Time.h"
#include "AztecOO.h"
#include "ml_include.h"
#include "Epetra_LinearProblem.h"
#include "ml_MultiLevelOperator.h"
#include "ml_epetra_utils.h"
#include "Trilinos_Util_CommandLineParser.h"
#include "Trilinos_Util_CrsMatrixGallery.h"
using namespace Trilinos_Util;
int main(int argc, char *argv[])
{
#ifdef EPETRA_MPI
MPI_Init(&argc,&argv);
#else
#endif
CommandLineParser CLP(argc,argv);
CrsMatrixGallery Gallery("", Comm, false);
if( CLP.Has("-problem_type") == false ) CLP.Add("-problem_type", "laplace_2d" );
if( CLP.Has("-problem_size") == false ) CLP.Add("-problem_size", "100" );
Gallery.Set(CLP);
AztecOO solver(*Problem);
solver.SetAztecOption(AZ_solver, AZ_cg);
ML *ml_handle;
int N_levels = 10;
ML_Set_PrintLevel(3);
ML_Create(&ml_handle,N_levels);
EpetraMatrix2MLMatrix(ml_handle, 0, Matrix);
ML_Aggregate *agg_object;
ML_Aggregate_Create(&agg_object);
ML_Aggregate_Set_MaxCoarseSize(agg_object,1);
N_levels = ML_Gen_MGHierarchy_UsingAggregation(ml_handle, 0,
ML_INCREASING, agg_object);
ML_Gen_Smoother_SymGaussSeidel(ml_handle, ML_ALL_LEVELS,
ML_BOTH, 1, ML_DEFAULT);
ML_Gen_Solver (ml_handle, ML_MGV, 0, N_levels-1);
ML_Epetra::MultiLevelOperator MLop(ml_handle,Comm,*Map,*Map);
solver.SetPrecOperator(&MLop);
solver.Iterate(1550, 1e-12);
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 EPETRA_MPI
MPI_Finalize() ;
#endif
exit(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-ml\n"
"--enable-triutils");
#ifdef HAVE_MPI
MPI_Finalize();
#endif
return 0;
}
#endif