00001 #include "NLPInterfacePack_NLPFirstOrderThyraModelEvaluator.hpp"
00002 #ifdef HAVE_MPI
00003 #include "EpetraExt_MultiPointModelEvaluator.h"
00004 #include "EpetraExt_MultiMpiComm.h"
00005 #endif
00006 #include "EpetraMultiPointModelEval4DOpt.hpp"
00007 #include "MoochoPack_MoochoThyraSolver.hpp"
00008 #include "Thyra_DefaultRealLinearSolverBuilder.hpp"
00009 #include "Thyra_EpetraModelEvaluator.hpp"
00010 #include "Teuchos_RCP.hpp"
00011 #include "Teuchos_GlobalMPISession.hpp"
00012 #include "Teuchos_VerboseObject.hpp"
00013 #include "Teuchos_CommandLineProcessor.hpp"
00014 #include "Teuchos_StandardCatchMacros.hpp"
00015
00016 int main( int argc, char* argv[] )
00017 #ifdef HAVE_MPI
00018 {
00019 using Teuchos::CommandLineProcessor;
00020 typedef AbstractLinAlgPack::value_type Scalar;
00021 using MoochoPack::MoochoSolver;
00022 using MoochoPack::MoochoThyraSolver;
00023
00024 bool dummySuccess = true;
00025
00026 Teuchos::GlobalMPISession mpiSession(&argc,&argv);
00027
00028 Teuchos::RCP<Teuchos::FancyOStream>
00029 out = Teuchos::VerboseObjectBase::getDefaultOStream();
00030
00031 try {
00032
00033 Thyra::DefaultRealLinearSolverBuilder lowsfCreator;
00034 MoochoThyraSolver solver;
00035
00036
00037
00038
00039
00040 Scalar xt0 = 1.0;
00041 Scalar xt1 = 1.0;
00042 Scalar pt0 = 2.0;
00043 Scalar pt1 = 0.0;
00044 Scalar d = 10.0;
00045 Scalar x00 = 1.0;
00046 Scalar x01 = 1.0;
00047 Scalar p00 = 2.0;
00048 Scalar p01 = 0.0;
00049 Scalar pL0 = -1e+50;
00050 Scalar pL1 = -1e+50;
00051 Scalar pU0 = +1e+50;
00052 Scalar pU1 = +1e+50;
00053 Scalar q0 = 0.0;
00054
00055 Scalar xL0 = -1e+50;
00056 Scalar xL1 = -1e+50;
00057 Scalar xU0 = +1e+50;
00058 Scalar xU1 = +1e+50;
00059
00060 CommandLineProcessor clp(false);
00061
00062 lowsfCreator.setupCLP(&clp);
00063 solver.setupCLP(&clp);
00064
00065 clp.setOption( "xt0", &xt0 );
00066 clp.setOption( "xt1", &xt1 );
00067 clp.setOption( "pt0", &pt0 );
00068 clp.setOption( "pt1", &pt1 );
00069 clp.setOption( "d", &d );
00070 clp.setOption( "x00", &x00 );
00071 clp.setOption( "x01", &x01 );
00072 clp.setOption( "p00", &p00 );
00073 clp.setOption( "q0", &q0 );
00074 clp.setOption( "p01", &p01 );
00075 clp.setOption( "pL0", &pL0 );
00076 clp.setOption( "pL1", &pL1 );
00077 clp.setOption( "pU0", &pU0 );
00078 clp.setOption( "pU1", &pU1 );
00079 clp.setOption( "xL0", &xL0 );
00080 clp.setOption( "xL1", &xL1 );
00081 clp.setOption( "xU0", &xU0 );
00082 clp.setOption( "xU1", &xU1 );
00083
00084 CommandLineProcessor::EParseCommandLineReturn
00085 parse_return = clp.parse(argc,argv,&std::cerr);
00086
00087 if( parse_return != CommandLineProcessor::PARSE_SUCCESSFUL )
00088 return parse_return;
00089
00090 lowsfCreator.readParameters(out.get());
00091 solver.readParameters(out.get());
00092
00093
00094
00095
00096
00097
00098 int numPoints = 20;
00099 Teuchos::RCP<EpetraExt::MultiMpiComm> globalComm =
00100 Teuchos::rcp(new EpetraExt::MultiMpiComm(MPI_COMM_WORLD, 1, numPoints));
00101 int myPoints = globalComm->NumTimeStepsOnDomain();
00102 int myFirstPoint = globalComm->FirstTimeStepOnDomain();
00103
00104 Teuchos::RCP<Epetra_MpiComm> epetra_comm = Teuchos::rcp(&globalComm->SubDomainComm(),false);
00105
00106
00107 Teuchos::RCP<EpetraMultiPointModelEval4DOpt>
00108 epetraModel = Teuchos::rcp(new EpetraMultiPointModelEval4DOpt(epetra_comm,xt0,xt1,pt0,pt1,d,x00,x01,p00,p01, q0));
00109 epetraModel->set_p_bounds(pL0,pL1,pU0,pU1);
00110 epetraModel->set_x_bounds(xL0,xL1,xU0,xU1);
00111
00112
00113 std::vector<Epetra_Vector*> multi_x_init(myPoints);
00114 Epetra_Vector init_vec(*(epetraModel->get_x_init().get()));
00115 for (int i=0; i<myPoints; i++) { multi_x_init[i] = &init_vec;}
00116
00117
00118 Teuchos::RCP<std::vector< Teuchos::RCP<Epetra_Vector> > > q_vec
00119 = Teuchos::rcp(new std::vector< Teuchos::RCP<Epetra_Vector> >(myPoints));
00120 for (int i=0; i<myPoints; i++) {
00121 q_vec->operator[](i) = Teuchos::rcp(new Epetra_Vector( *(epetraModel->get_p_map(1))));
00122 q_vec->operator[](i)->operator[](0) = 0.0 + 0.1*(i+myFirstPoint);
00123 }
00124
00125
00126 Teuchos::RCP<EpetraExt::MultiPointModelEvaluator>
00127 multiPointModel = rcp(new EpetraExt::MultiPointModelEvaluator(
00128 epetraModel, globalComm, multi_x_init, q_vec));
00129
00130
00131
00132 Teuchos::RCP<Thyra::LinearOpWithSolveFactoryBase<double> >
00133 lowsFactory = lowsfCreator.createLinearSolveStrategy("");
00134
00135 Teuchos::RCP<Thyra::EpetraModelEvaluator>
00136 epetraThyraModel = rcp(new Thyra::EpetraModelEvaluator());
00137
00138 epetraThyraModel->initialize(multiPointModel,lowsFactory);
00139
00140
00141
00142
00143
00144
00145 solver.setModel(epetraThyraModel);
00146
00147
00148 solver.readInitialGuess(out.get());
00149
00150
00151 const MoochoSolver::ESolutionStatus solution_status = solver.solve();
00152
00153
00154 lowsfCreator.writeParamsFile(*lowsFactory);
00155 solver.writeParamsFile();
00156
00157
00158
00159
00160
00161 return solution_status;
00162
00163 }
00164 TEUCHOS_STANDARD_CATCH_STATEMENTS(true,*out,dummySuccess)
00165
00166 return MoochoSolver::SOLVE_RETURN_EXCEPTION;
00167
00168 }
00169 #else
00170 {
00171 cout << "NLPThyraEpetraModelEval4DOpt does nothing when HAVE_MPI is not defined" << endl;
00172 return 0;
00173 }
00174 #endif