MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ExampleNLPDirectMain.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Moocho: Multi-functional Object-Oriented arCHitecture for Optimization
5 // Copyright (2003) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 //
42 
43 #include <assert.h>
44 
45 #include <fstream>
46 #include <fstream>
47 #include <iostream>
48 #include <iomanip>
49 
53 #include "Teuchos_Workspace.hpp"
55 
56 int main( int argc, char* argv[] ) {
57 
58  using std::endl;
59  using std::setw;
60  namespace mmp = MemMngPack;
61  using Teuchos::RCP;
64  namespace NLPIP = NLPInterfacePack;
65  namespace rsqp = MoochoPack;
66  using rsqp::MoochoSolver;
67 
69 
70  using Teuchos::Workspace;
72 
73  int err = 0;
74 
75 /*
76  // Print out the input arguments
77  printf("argc = %d\n",argc);
78  {for( int i = 0; i < argc; ++i) {
79  printf("argv[%d] = %s\n",i,argv[i]);
80  }}
81 */
82  // Get an idea of what processors we have.
83  MPI_Init(&argc,&argv);
84  int num_proc, proc_rank;
85  MPI_Comm_size( MPI_COMM_WORLD, &num_proc );
86  MPI_Comm_rank( MPI_COMM_WORLD, &proc_rank );
87 /*
88  // Print out the input arguments
89  printf("\nproc_rank = %d\n",proc_rank);
90  printf("argc = %d\n",argc);
91  {for( int i = 0; i < argc; ++i) {
92  printf("argv[%d] = %s\n",i,argv[i]);
93  }}
94 */
95 
96  // Define program return values
97  const int
98  PROG_SUCCESS = 0,
99  PROG_NLP_TEST_ERR = -1,
100  PROG_EXCEPTION = -2,
101  PROG_MAX_ITER_EXEEDED = -3,
102  PROG_MAX_TIME_EXEEDED = -4;
103 
104  int prog_return = PROG_SUCCESS;
105 
106  // Set the output stream
107  std::ostream &out = std::cout;
108  std::ostream &eout = std::cerr;
109  Teuchos::oblackholestream blackhole;
110 
111  try {
112 
113  //
114  // Initialize stuff
115  //
116 
117  size_type n;
118  value_type xo;
119  bool has_bounds;
120  bool dep_bounded;
121 
122  VectorSpace::space_ptr_t vec_space;
123  const int err = AbstractLinAlgPack::exampleNLPDiagSetup(argc,argv,MPI_COMM_WORLD,&vec_space,&n,&xo,&has_bounds,&dep_bounded);
124  if(err) return err;
125 
126  // Create and test the NLP using this vector space object
127  const MoochoSolver::ESolutionStatus
128  solve_return = NLPIP::ExampleNLPDirectRun(
129  *vec_space, xo, has_bounds, dep_bounded
130  ,proc_rank == 0 ? &out : &blackhole // console_out
131  ,proc_rank == 0 ? &eout : &blackhole // error_out
132  ,proc_rank == 0 ? false : true // throw_solve_exception
133  ,proc_rank == 0 ? NULL : &blackhole // algo_out
134  ,proc_rank == 0 ? NULL : &blackhole // summary_out
135  ,proc_rank == 0 ? NULL : &blackhole // journal_out
136  );
137 
138  switch(solve_return) {
139  case MoochoSolver::SOLVE_RETURN_SOLVED:
140  prog_return = PROG_SUCCESS;
141  break;
142  case MoochoSolver::SOLVE_RETURN_MAX_ITER:
143  prog_return = PROG_MAX_ITER_EXEEDED;
144  break;
145  case MoochoSolver::SOLVE_RETURN_MAX_RUN_TIME:
146  prog_return = PROG_MAX_TIME_EXEEDED;
147  break;
148  case MoochoSolver::SOLVE_RETURN_NLP_TEST_FAILED:
149  prog_return = PROG_NLP_TEST_ERR;
150  break;
151  case MoochoSolver::SOLVE_RETURN_EXCEPTION:
152  prog_return = PROG_EXCEPTION;
153  break;
154  default:
156  }
157 
158  } // end try
159  catch(const std::exception& excpt) {
160  eout << "\nCaught a std::exception on process " << proc_rank<< ": " << excpt.what() << endl;
161  prog_return = PROG_EXCEPTION;
162  }
163  catch(...) {
164  eout << "\nCaught an unknown exception on process " << proc_rank<< "\n";
165  prog_return = PROG_EXCEPTION;
166  }
167 
168  MPI_Finalize();
169 
170  return prog_return;
171 }
AbstractLinAlgPack::size_type size_type
MoochoPack::MoochoSolver::ESolutionStatus ExampleNLPDirectRun(const VectorSpace &vec_space, value_type xo, bool has_bounds, bool dep_bounded, std::ostream *console_out, std::ostream *error_out, bool throw_solve_exception=false, std::ostream *algo_out=NULL, std::ostream *summary_out=NULL, std::ostream *journal_out=NULL)
Function accepts a VectorSpace object and then uses it to define an example NLP and run MoochoPack::M...
int MPI_Comm_size(MPI_Comm comm, int *size)
Definition: RTOp_mpi.c:71
int MPI_Comm_rank(MPI_Comm comm, int *rank)
Definition: RTOp_mpi.c:77
int MPI_Init(int *argc, char ***argv)
Definition: RTOp_mpi.c:61
int main(int argc, char *argv[])
int exampleNLPDiagSetup(int argc, char *argv[], MPI_Comm comm, Teuchos::RCP< const VectorSpace > *vec_space, int *n, value_type *xo, bool *has_bounds, bool *dep_bounded)
Create a vector space given the input arguments argc, argv[] and an MPI communicator.
Abstract interface for objects that represent a space for mutable coordinate vectors.
std::ostream * out
AbstractLinAlgPack::value_type value_type
#define MPI_COMM_WORLD
Definition: RTOp_mpi.h:68
int n
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
TEUCHOSCORE_LIB_DLL_EXPORT Teuchos::RCP< WorkspaceStore > get_default_workspace_store()
int MPI_Finalize(void)
Definition: RTOp_mpi.c:66