This page describes how to write the main file for using NOX's Epetra implementation of the Group and Vector.
This example can be found in Trilinos/packages/nox/test/epetra/1Dfem/1Dfem.C
Begin by including the NOX library header files and any other headers needed by your application:
#include "NOX.H"
#include "NOX_Epetra.H"
#ifdef HAVE_MPI
#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_Map.h"
#include "Epetra_LinearProblem.h"
#include "1DfemInterface.H"
#include "1DfemPrePostOperator.H"
#include "Teuchos_ParameterList.hpp"
For convenience define the following:
Begin by initializing MPI (if building a parallel version) and create the Epetra communicator for MPI. In this case autoconf defines the flag HAVE_MPI if we are building a parallel version.
int main(int argc, char *argv[])
{
#ifdef HAVE_MPI
MPI_Init(&argc,&argv);
#endif
#ifdef HAVE_MPI
Epetra_MpiComm Comm( MPI_COMM_WORLD );
#else
Epetra_SerialComm Comm;
#endif
Setup some initial values.
int MyPID = Comm.MyPID();
int NumProc = Comm.NumProc();
bool verbose = false;
if (argc > 1)
if (argv[1][0]=='-' && argv[1][1]=='v')
verbose = true;
int NumGlobalElements = 0;
if ((argc > 2) && (verbose))
NumGlobalElements = atoi(argv[2]) + 1;
else if ((argc > 1) && (!verbose))
NumGlobalElements = atoi(argv[1]) + 1;
else
NumGlobalElements = 101;
if (NumGlobalElements < NumProc) {
cout << "numGlobalBlocks = " << NumGlobalElements
<< " cannot be < number of processors = " << NumProc << std::endl;
cout << "Test failed!" << std::endl;
throw "NOX Error";
}
Create your object that derives from NOX::Epetra::Interface::Required so that nox has access to the set of nonlinear equations to solve.
Teuchos::RCP<Interface> interface =
Teuchos::rcp(new Interface(NumGlobalElements, Comm));
Grab the initial guess vector and initialize it and the interface object.
Teuchos::RCP<Epetra_Vector> soln = interface->getSolution();
Teuchos::RCP<NOX::Epetra::Vector> noxSoln =
interface->setPDEfactor(1000.0);
soln->PutScalar(1.0);
Create a parameter list and choose the options for the NOX solver. Note that the linear solver parameters in teh "Linear Solver" sublist are dependent upon what NOX::Epetra::LinearSystem object you are using. Currently NOX only comes with one concrete linear system implementation: NOX::EpetraLinearSystemAztecOO.
Teuchos::RCP<Teuchos::ParameterList> nlParamsPtr =
Teuchos::rcp(new Teuchos::ParameterList);
Teuchos::ParameterList& nlParams = *(nlParamsPtr.get());
nlParams.set("Nonlinear Solver", "Line Search Based");
Teuchos::ParameterList& printParams = nlParams.sublist("Printing");
printParams.set("MyPID", MyPID);
printParams.set("Output Precision", 3);
printParams.set("Output Processor", 0);
if (verbose)
printParams.set("Output Information",
else
Teuchos::ParameterList& searchParams = nlParams.sublist("Line Search");
searchParams.set("Method", "Full Step");
Teuchos::ParameterList& dirParams = nlParams.sublist("Direction");
dirParams.set("Method", "Newton");
Teuchos::ParameterList& newtonParams = dirParams.sublist("Newton");
newtonParams.set("Forcing Term Method", "Constant");
Teuchos::ParameterList& lsParams = newtonParams.sublist("Linear Solver");
lsParams.set("Aztec Solver", "GMRES");
lsParams.set("Max Iterations", 800);
lsParams.set("Tolerance", 1e-4);
lsParams.set("Preconditioner", "New Ifpack");
lsParams.set("Preconditioner Reuse Policy", "Reuse");
lsParams.set("Max Age Of Prec", 5);
Optionally the user can define methods that will be called before and after each nonlineaer iteration and before and after each nonlinear solve. This is done by creating an object derived from NOX::Abstract::PrePostOperator.
Teuchos::RCP<NOX::Abstract::PrePostOperator> ppo =
Teuchos::rcp(new UserPrePostOperator(printing));
nlParams.sublist("Solver Options").set("User Defined Pre/Post Operator",
ppo);
Set the status test check option.
nlParams.sublist("Solver Options").
Create a Jacobian-Free Newton-Krylov method by using a NOX::Epetra::MatrixFree operator for the Jacobian.
Teuchos::RCP<NOX::Epetra::MatrixFree> MF =
*noxSoln));
Create a Finite Difference operator to estimate the Jacobian matrix for preconditioning.
Teuchos::RCP<NOX::Epetra::FiniteDifference> FD =
*soln));
Create a linear system derived from NOX::Epetra::LinearSystem. NOX comes with a concrete implementation of the linear system class that uses AztecOO as the linear solver called NOX::Epetra::LinearSystemAztecOO. In this case we will use the constructor that specifies the Jacobian operator and a preconditioner matrix to be used with an internally constructed AztecOO preconditioner.
Teuchos::RCP<NOX::Epetra::Interface::Jacobian> iJac = MF;
Teuchos::RCP<NOX::Epetra::Interface::Preconditioner> iPrec = FD;
Teuchos::RCP<NOX::Epetra::LinearSystemAztecOO> linSys =
iJac, MF,
iPrec, FD,
*soln));
Create the group using the linear system.
Teuchos::RCP<NOX::Epetra::Group> grpPtr =
iReq,
initialGuess,
linSys));
Create the status tests to determine how to terminate the nonlinear solve.
Teuchos::RCP<NOX::StatusTest::NormF> absresid =
Teuchos::RCP<NOX::StatusTest::NormF> relresid =
Teuchos::RCP<NOX::StatusTest::NormUpdate> update =
Teuchos::RCP<NOX::StatusTest::NormWRMS> wrms =
Teuchos::RCP<NOX::StatusTest::Combo> converged =
converged->addStatusTest(absresid);
converged->addStatusTest(relresid);
converged->addStatusTest(wrms);
converged->addStatusTest(update);
Teuchos::RCP<NOX::StatusTest::MaxIters> maxiters =
Teuchos::RCP<NOX::StatusTest::FiniteValue> fv =
Teuchos::RCP<NOX::StatusTest::Combo> combo =
combo->addStatusTest(fv);
combo->addStatusTest(converged);
combo->addStatusTest(maxiters);
Create the solver and solve the problem
Teuchos::RCP<NOX::Solver::Generic> solver =
NOX::Solver::buildSolver(grpPtr, combo, nlParamsPtr);
Get the underlying solution vector in terms of an Epetra_Vector and do whatever you want with it.
const Epetra_Vector& finalSolution =
getEpetraVector();
Exit simulation.