#include "Didasko_ConfigDefs.h"
#if defined(HAVE_DIDASKO_EPETRA) && defined(HAVE_DIDASKO_NOX)
#include <iostream>
#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_Map.h"
#include "Epetra_RowMatrix.h"
#include "Epetra_CrsMatrix.h"
#include "NOX.H"
#include "NOX_Epetra_Interface.H"
#include "NOX_Epetra_Group.H"
class SimpleProblemInterface : public NOX::Epetra::Interface {
public:
SimpleProblemInterface( Epetra_Vector & InitialGuess,
Epetra_Vector & ExactSolution )
{
InitialGuess_ = new Epetra_Vector(InitialGuess);
ExactSolution_ = new Epetra_Vector(ExactSolution);
};
~SimpleProblemInterface()
{
};
bool computeF(const Epetra_Vector & x, Epetra_Vector & f,
NOX::Epetra::Interface::FillType F )
{
f[0] = x[0]*x[0] + x[1]*x[1] - 1.0;
f[1] = x[1] - x[0]*x[0];
return true;
};
bool computeJacobian(const Epetra_Vector & x, Epetra_Operator & Jac)
{
Epetra_CrsMatrix * J;
J = dynamic_cast<Epetra_CrsMatrix*>(&Jac);
if (J == NULL) {
cout << "*ERR* Problem_Interface::computeJacobian() - The supplied" << endl;
cout << "*ERR* Epetra_Operator is NOT an Epetra_CrsMatrix!" << endl;
throw;
}
int NumMyElements = J->Map().NumMyElements();
int * MyGlobalElements = J->Map().MyGlobalElements();
int indices[2];
double values[2];
indices[0] = 0; indices[1] = 1;
for( int i=0 ; i<NumMyElements ; ++i ) {
switch( MyGlobalElements[i] ) {
case 0:
values[0] = 2.0 * x[0];
values[1] = 2.0 * x[1];
J->ReplaceGlobalValues(0, 2, values, indices );
break;
case 1:
values[0] = - 2.0 * x[0];
values[1] = 1.0;
J->ReplaceGlobalValues(1, 2, values, indices );
break;
default:
cerr << "*ERR*" << endl;
exit( EXIT_FAILURE );
}
}
return true;
}
bool computePrecMatrix(const Epetra_Vector & x, Epetra_RowMatrix & M)
{
cout << "*ERR* SimpleProblem::preconditionVector()\n";
cout << "*ERR* don't use explicit preconditioning" << endl;
exit( 0 );
throw 1;
}
bool computePreconditioner(const Epetra_Vector & x, Epetra_Operator & O)
{
cout << "*ERR* SimpleProblem::preconditionVector()\n";
cout << "*ERR* don't use explicit preconditioning" << endl;
exit( 0 );
throw 1;
}
private:
Epetra_Vector * InitialGuess_;
Epetra_Vector * ExactSolution_;
};
int main( int argc, char **argv )
{
#ifdef HAVE_MPI
MPI_Init(&argc, &argv);
Epetra_MpiComm Comm(MPI_COMM_WORLD);
#else
Epetra_SerialComm Comm;
#endif
if (Comm.NumProc() != 1) {
if (Comm.MyPID() == 0)
cerr << "This example can be run with one process only!" << endl;
#ifdef HAVE_MPI
MPI_Finalize();
#endif
exit(EXIT_SUCCESS);
}
Epetra_Map Map(2,0,Comm);
Epetra_Vector ExactSolution(Map);
ExactSolution[0] = sqrt(0.5*(sqrt(5.0)-1));
ExactSolution[1] = 0.5*(sqrt(5.0)-1);
Epetra_Vector InitialGuess(Map);
InitialGuess[0] = 0.5;
InitialGuess[1] = 0.5;
SimpleProblemInterface Interface(InitialGuess,ExactSolution);
NOX::Parameter::List nlParams;
nlParams.setParameter("Nonlinear Solver", "Line Search Based");
NOX::Parameter::List& printParams = nlParams.sublist("Printing");
printParams.setParameter("MyPID", Comm.MyPID());
printParams.setParameter("Output Precision", 3);
printParams.setParameter("Output Processor", 0);
printParams.setParameter("Output Information",
NOX::Utils::OuterIteration +
NOX::Utils::OuterIterationStatusTest +
NOX::Utils::InnerIteration +
NOX::Utils::Parameters +
NOX::Utils::Details +
NOX::Utils::Warning);
NOX::Parameter::List& searchParams = nlParams.sublist("Line Search");
searchParams.setParameter("Method", "Full Step");
NOX::Parameter::List& dirParams = nlParams.sublist("Direction");
dirParams.setParameter("Method", "Newton");
NOX::Parameter::List& newtonParams = dirParams.sublist("Newton");
newtonParams.setParameter("Forcing Term Method", "Constant");
NOX::Parameter::List& lsParams = newtonParams.sublist("Linear Solver");
lsParams.setParameter("Aztec Solver", "GMRES");
lsParams.setParameter("Max Iterations", 800);
lsParams.setParameter("Tolerance", 1e-4);
lsParams.setParameter("Output Frequency", 50);
lsParams.setParameter("Aztec Preconditioner", "ilu");
Epetra_CrsMatrix A(Copy,Map,2);
int NumMyElements = A.Map().NumMyElements();
int * MyGlobalElements = A.Map().MyGlobalElements();
int indices[2];
double values[2];
indices[0]=0; indices[1]=1;
for( int i=0 ; i<NumMyElements ; ++i ) {
switch( MyGlobalElements[i] ) {
case 0:
values[0] = 2.0 * InitialGuess[0];
values[1] = 2.0 * InitialGuess[1];
A.InsertGlobalValues(0, 2, values, indices );
break;
case 1:
values[0] = - 2.0 * InitialGuess[0];
values[1] = 1.0;
A.InsertGlobalValues(1, 2, values, indices );
break;
}
}
A.TransformToLocal();
NOX::Epetra::Group grp(printParams, lsParams, Interface, InitialGuess,
dynamic_cast<Epetra_RowMatrix&>(A));
NOX::StatusTest::NormF statusTestA(grp, 1.0e-4);
NOX::StatusTest::MaxIters statusTestB(20);
NOX::StatusTest::Combo statusTestsCombo(NOX::StatusTest::Combo::OR, statusTestA, statusTestB);
NOX::Parameter::List solverParameters;
solverParameters.sublist("Printing").setParameter("Output Information",
NOX::Utils::OuterIteration);
solverParameters.setParameter("Nonlinear Solver", "Line Search Based");
NOX::Parameter::List& lineSearchParameters = solverParameters.sublist("Line Search");
lineSearchParameters.setParameter("Method","More'-Thuente");
NOX::Solver::Manager solver(grp, statusTestsCombo, solverParameters);
NOX::StatusTest::StatusType status = solver.solve();
cout << "\n" << "-- Parameter List From Solver --" << "\n";
solver.getParameterList().print(cout);
const NOX::Epetra::Group & finalGroup =
dynamic_cast<const NOX::Epetra::Group&>(solver.getSolutionGroup());
const Epetra_Vector & finalSolution =
(dynamic_cast<const NOX::Epetra::Vector&>(finalGroup.getX())).getEpetraVector();
cout << "Computed solution : " << endl;
cout << finalSolution;
cout << "Exact solution : " << endl;
cout << ExactSolution;
#ifdef HAVE_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:");
puts("--enable-epetra");
puts("--enable-nox");
puts("--enable-nox-epetra");
#ifdef HAVE_MPI
MPI_Finalize();
#endif
return(EXIT_SUCCESS);
}
#endif