#include "ml_config.h"
#include "ml_common.h"
#ifdef HAVE_ML_MLAPI
#include "MLAPI.h"
using namespace Teuchos;
using namespace MLAPI;
int main(int argc, char *argv[])
{
#ifdef HAVE_MPI
MPI_Init(&argc,&argv);
#endif
try {
for (int i = 0 ; i < MySpace.GetNumGlobalElements() ; ++i) {
if (i)
A(i, i - 1) = -1.0;
A(i, i + 1) = -1.0;
A(i, i) = 2.0;
}
}
std::cout << A;
x = 1.0;
for (int i = 0 ; i < y.GetMyLength() ; ++i)
y(i) = 1.0 * i + x(i);
z = x + y;
double norm = sqrt(z * z);
double Anorm = sqrt(z * (A * z));
x = x / (x * y);
std::cout << "2-norm of z = " << norm << std::endl;
std::cout << "A-norm of z = " << Anorm << std::endl;
}
C = A + B;
C = 1.0 * A - B * 2.0;
x = invA * (A * x) - x;
double NormX = sqrt(x * x);
std::cout << "Norm of inv(A) * A * x - x = " << NormX << std::endl;
std::cout << MySpace;
std::cout << x;
std::cout << A;
std::cout << "ER(" << MySpace(i) << ") = " << ER(i) << std::endl;
matlab << "% you can add any MATLAB command this way\n";
x.SetLabel("x");
matlab << x;
A.SetLabel("A");
matlab << A;
matlab << ER;
matlab << EI;
matlab << "plot(ER, EI, 'o')\n";
}
catch (const int e) {
std::cout << "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("Please check your configure line.");
#ifdef HAVE_MPI
MPI_Finalize();
#endif
return(0);
}
#endif // #if defined(HAVE_ML_MLAPI)