MoochoPack : Framework for Large-Scale Optimization Algorithms  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
List of all members
MoochoPack::MoochoSolver Class Reference

Universal interface to a MOOCHO solver. More...

#include <MoochoPack_MoochoSolver.hpp>

enum  EOutputToBlackHole
 
enum  EConfigOptions
 
enum  ESolutionStatus
 
typedef RCP
< NLPInterfacePack::NLP
nlp_ptr_t
 
typedef RCP
< IterationPack::AlgorithmTracker
track_ptr_t
 
typedef RCP< NLPAlgoConfigconfig_ptr_t
 
typedef RCP
< OptionsFromStreamPack::OptionsFromStream
options_ptr_t
 
typedef RCP< std::ostream > ostream_ptr_t
 

Initialization and algorithm configuration

 MoochoSolver (const std::string &options_file_name="Moocho.opt", const std::string &extra_options_str="")
 Constructs to uninitialized. More...
 
OptionsFromStreamPack::CommandLineOptionsFromStreamProcessorcommandLineOptionsFromStreamProcessor ()
 
const
OptionsFromStreamPack::CommandLineOptionsFromStreamProcessor
commandLineOptionsFromStreamProcessor () const
 
void setup_commandline_processor (Teuchos::CommandLineProcessor *clp)
 Setup the commandline processor to process commandline options. More...
 
void set_nlp (const nlp_ptr_t &nlp)
 Set the NLP to be solved. More...
 
const nlp_ptr_tget_nlp () const
 Get the non-const smart pointer to the set NLP object. More...
 
void set_track (const track_ptr_t &track)
 
const track_ptr_tget_track () const
 Get the non-const smart pointer to the set AlgorithmTracker object. More...
 
void set_config (const config_ptr_t &config)
 Set the algorithm configuration object. More...
 
const config_ptr_tget_config () const
 Return the configuration object being used. More...
 
void set_options (const options_ptr_t &options)
 Set the various options to use. More...
 
const options_ptr_tget_options () const
 Get the OptionsFromStream object being used to extract the options from. More...
 

Exception handling

void set_error_handling (bool throw_exceptions, const ostream_ptr_t &error_out)
 Set the error output and whether exceptions will be thrown from these functions or not. More...
 
bool throw_exceptions () const
 Return if exceptions will be thrown out of this->solve_nlp(). More...
 
const ostream_ptr_terror_out () const
 Return the std::ostream object used for error reporting on exceptions. More...
 

Collective outputting control

void set_output_context (const std::string &file_context_postfix, EOutputToBlackHole output_to_black_hole=OUTPUT_TO_BLACK_HOLE_DEFAULT, const int procRank=-1, const int numProcs=-1)
 Setup the context for outputting. More...
 

Individual outputting control

void set_output_file_tag (const std::string &)
 Set a tag for output file names for all file names that are created internally. More...
 
void do_console_outputting (bool)
 Turn on and off console outputting. More...
 
bool do_console_outputting () const
 Return if console outputting is performed or not. More...
 
void set_console_out (const ostream_ptr_t &console_out)
 Set the std::ostream object to use for console output by a MoochoTrackerConsoleStd object. More...
 
const ostream_ptr_tget_console_out () const
 Get the non-const smart pointer to the set output stream for console outputting. More...
 
void do_summary_outputting (bool)
 Turn on and off summary outputting. More...
 
bool do_summary_outputting () const
 Return if summary outputting is performed or not. More...
 
void set_summary_out (const ostream_ptr_t &summary_out)
 Set the std::ostream object to use for summary output. More...
 
const ostream_ptr_tget_summary_out () const
 Get the non-const smart pointer to the set output stream for summary outputting. More...
 
void do_journal_outputting (bool)
 Turn on and off journal outputting. More...
 
bool do_journal_outputting () const
 Return if journal outputting is performed or not. More...
 
void set_journal_out (const ostream_ptr_t &journal_out)
 Set the std::ostream object to use for journal output by the MOOCHO step objects. More...
 
const ostream_ptr_tget_journal_out () const
 Get the non-const smart pointer to the set output stream for journal outputting. More...
 
void do_algo_outputting (bool)
 Turn on and off algo outputting. More...
 
bool do_algo_outputting () const
 Return if algo outputting is performed or not. More...
 
void generate_stats_file (bool)
 Turn on and off the generation of a statistics file. More...
 
bool generate_stats_file () const
 Return if a statistics file will be generated or not. More...
 
void set_algo_out (const ostream_ptr_t &algo_out)
 Set the std::ostream object to use for algorithm output. More...
 
const ostream_ptr_tget_algo_out () const
 Get the non-const smart pointer to the set output stream for algo outputting. More...
 
RCP< std::ostream > generate_output_file (const std::string &fileNameBase) const
 Generate an output file given a base file name. More...
 

Solve the NLP

void update_solver () const
 Setup the state of the solver and get ready for a solve. More...
 
ESolutionStatus solve_nlp () const
 Solve the NLP. More...
 

Get the underlying solver object

NLPSolverClientInterfaceget_solver ()
 Get the underlying NLPSolverClientInterface object. More...
 
const NLPSolverClientInterfaceget_solver () const
 

Detailed Description

Universal interface to a MOOCHO solver.

This class is designed to act as a simple encapsulation to several other smaller components needed to solve an NLP. This class is an instance of the popular "Facade" design pattern (Design Patterns, 1995). This class has a defualt implementation based on NLPAlgoConfigMamaJama but the client can set different algorithm configuration objects (see the requirments/specifications section below).

There are two distinct activities associated with using a MoochoSolver object:

  1. Algorithm configuration (i.e. configure an encapsulated NLPAlgoClientInterface object with an NLP and other objects).
  2. NLP solution (i.e. call NLPSolverClientInterface::find_min() on the encapuslaited solver object).

In the algorithm configuration phase, the client must, at a minimum, set the NLP object for the NLP to be solved using this->set_nlp(). The NLP object is needed so that the algorithm configuration object can adapt the MOOCHO algorithm to the NLP in the best way possible. The configuration phase can also include setting a user-defined track object(s) and a user-defined NLPAlgoConfig object. An NLP is solved by calling the method this->solve_nlp() which returns an enum stating what happended and reporting the final point to the NLP object.

This class encapsulates an NLPAlgoClientInterface object and takes over some of the drudgery of working with this interface. In most cases all the options that can be set to configuration object and other algorithmic objects can be set using an OptionsFromStreamPack::OptionsFromStream object by calling this->set_options().

Options specific to this class and the configuration object (down to the lower algorithmic objects that it creates) can be set through an OptionsFromStreamPack::OptionsFromStream object by passing it to this->set_options(). The files Moocho.opt.MoochoSolver, Moocho.opt.DecompositionSystemStateStepBuilderStd and Moocho.opt.NLPAlgoConfigMamaJama conatain the listing of these options as well as some documentation. An options file Moocho.opt can be generated automatically using the shell script generate_opt_file.pl.

Requirements / Specifications

The requirements and specifications for this class are stated below. More detailed scenarios are shown elsewhere (??? where ???).

  1. Base default implementation on NLPAlgoConfigMamaJama and require minimal effort to quickly solve an NLP. This includes setting up standard IterationPack::AlgorithmTracker objects and taking care of exceptions etc.
    Enabler: The client can simply not call this->set_config() or call this->set_config(NULL) which will result in the default configuration object being used.
  2. Allow clients the ability to insert a customized AlgorithmTracker object, in addition to the other standard track objects.
    Enabler: An extra user defined track object can be set using this->set_extra_track() and can be unset using this->set_extra_track(NULL). Multiple track objects can be handled using the IterationPack::AlgorithmTrackerComposite subclass.
  3. Allow clients to set the targets for standard output streams (i.e. summary_out, journal_out and console_out) at all times (i.e. between successive solves) or just use default output files.
    Enabler: Default output targets can be used by doing nothing. User specified output streams can be set using this->set_console_out(), this->set_summary_out() and this->set_journal_out(). The same output files can be appended to for successive NLP solves by doing nothing. The output files can be changed between NLP runs using this->set_console_out(), this->set_summary_out() and this->set_journal_out(). Default output files can be overwritten between successive NLP solves by calling this->set_console_out(NULL), this->set_summary_out(NULL) and this->set_journal_out(NULL).
  4. Allow clients to set any special options in NLPAlgoConfigMamaJama beyond what can be set through the OptionsFromStream object (i.e. set a specialized BasisSystem object).
    Enabler: This can be accomplished by having the client create a NLPAlgoConfigMamaJama object itself and then configure it using the published interface before explicitly setting it using this->set_config().
  5. Allow clients to precisly control how an algorithm is configured and initialized beyond NLPAlgoConfigMamaJama and between NLP solves.
    Enabler: This can be accomplised by allowing the client to set a customized NLPAlgoConfig object. Clients can simply modify algorithms created by NLPAlgoConfigMamaJama using delegation or subclassing (delegation is to be prefered).
  6. Allow clients to solve the same NLP (i.e. same dimensions, same structure etc) multiple times with the same configured MOOCHO algorithm.
    Enabler: This can be done By simply calling this->get_nlp() (if needed to access the NLP that was set using this->set_nlp()), modifying the NLP object in some way (i.e. a new initial point) and then calling this->solve_nlp().
  7. Allow clients to configure a new MOOCHO algorithm with a potentially new NLP object (i.e. different dimensions, different structure etc).
    Enabler: The client can just call this->set_uninitialized() which is equivalent to setting the state of the object after the default constructor. In this case the client will have to go through the entire reinitialization phase again. Or, in order to use the same NLP, track and configuration objects but start off with a fresh algorithm configuration the client can just call this->set_nlp().

ToDo: Finish documentation!

Definition at line 157 of file MoochoPack_MoochoSolver.hpp.

Member Typedef Documentation

Public types

Definition at line 164 of file MoochoPack_MoochoSolver.hpp.

Definition at line 166 of file MoochoPack_MoochoSolver.hpp.

Definition at line 168 of file MoochoPack_MoochoSolver.hpp.

Definition at line 170 of file MoochoPack_MoochoSolver.hpp.

Definition at line 172 of file MoochoPack_MoochoSolver.hpp.

Member Enumeration Documentation

Definition at line 178 of file MoochoPack_MoochoSolver.hpp.

Definition at line 184 of file MoochoPack_MoochoSolver.hpp.

Definition at line 189 of file MoochoPack_MoochoSolver.hpp.

Constructor & Destructor Documentation

MoochoPack::MoochoSolver::MoochoSolver ( const std::string &  options_file_name = "Moocho.opt",
const std::string &  extra_options_str = "" 
)

Constructs to uninitialized.

Postconditions:

Definition at line 75 of file MoochoPack_MoochoSolver.cpp.

Member Function Documentation

OptionsFromStreamPack::CommandLineOptionsFromStreamProcessor & MoochoPack::MoochoSolver::commandLineOptionsFromStreamProcessor ( )

Definition at line 122 of file MoochoPack_MoochoSolver.cpp.

const OptionsFromStreamPack::CommandLineOptionsFromStreamProcessor & MoochoPack::MoochoSolver::commandLineOptionsFromStreamProcessor ( ) const

Definition at line 128 of file MoochoPack_MoochoSolver.cpp.

void MoochoPack::MoochoSolver::setup_commandline_processor ( Teuchos::CommandLineProcessor clp)

Setup the commandline processor to process commandline options.

Definition at line 133 of file MoochoPack_MoochoSolver.cpp.

void MoochoPack::MoochoSolver::set_nlp ( const nlp_ptr_t nlp)

Set the NLP to be solved.

Preconditions:

  • nlp.get() != NULL (throw std::invalid_argument).
  • This NLP object must be ready to have nlp->initialize() called but nlp->is_initialized() can be false.

Postconditions:

  • this->get_nlp().get() == nlp.get()
  • This will cause the MOOCHO algorithm to be reconfigured before the NLP is solved again.

Definition at line 191 of file MoochoPack_MoochoSolver.cpp.

const MoochoSolver::nlp_ptr_t & MoochoPack::MoochoSolver::get_nlp ( ) const

Get the non-const smart pointer to the set NLP object.

Returns
Returns the nlp object pointer set by this->set_nlp().

Definition at line 198 of file MoochoPack_MoochoSolver.cpp.

void MoochoPack::MoochoSolver::set_track ( const track_ptr_t track)

Set an extra user defined AlgorithmTracker object to monitor the algorithm.

Postconditions:

Definition at line 203 of file MoochoPack_MoochoSolver.cpp.

const MoochoSolver::track_ptr_t & MoochoPack::MoochoSolver::get_track ( ) const

Get the non-const smart pointer to the set AlgorithmTracker object.

Returns
Returns the track object pointer set by this->set_track().

Definition at line 210 of file MoochoPack_MoochoSolver.cpp.

void MoochoPack::MoochoSolver::set_config ( const config_ptr_t config)

Set the algorithm configuration object.

Parameters
config[in] The algorithm configuration object to use the next time this->do_config_algo() or this->solve_nlp() are called.

Postconditions:

  • [config.get() == NULL] A NLPAlgoConfigMamaJama object will be used to configure the MOOCHO algorithm the next time that this->do_config_algo() or this->solve_nlp() are called.
  • [config.get() != NULL] The object *config will be used to configure the MOOCHO algorithm the next time that this->do_config_algo() or this->solve_nlp() are called.
  • A reconfiguration of the MOOCHO algorithm will be forced the next time that this->do_config_algo() or this->solve_nlp() are called.
  • this->get_config().get() == config.get()

Definition at line 215 of file MoochoPack_MoochoSolver.cpp.

const MoochoSolver::config_ptr_t & MoochoPack::MoochoSolver::get_config ( ) const

Return the configuration object being used.

Returns
Returns return.get() == config.get() passed to the last call to this->set_config(config).

Definition at line 223 of file MoochoPack_MoochoSolver.cpp.

void MoochoPack::MoochoSolver::set_options ( const options_ptr_t options)

Set the various options to use.

Parameters
options[in] Smart pointer to the OptionsFromStream object to extract the options for this as well as for the configuration object from. If options.get() != NULL then this object must not be destoryed until this->set_options(NULL) is called or this is destroyed.

Postconditions:

  • [options.get() == NULL] The file "Moocho.opt" will be looked for and opened in the current directory. If this file does not exist, then a default set of options will be used. If this file does exist then the options will be read from this file.
  • [options.get() != NULL] The options will be read from *options.
  • A reconfiguration of the MOOCHO algorithm will be forced the next time that this->do_config_algo() or this->solve_nlp() are called.
  • this->get_options().get() == options.get()

Definition at line 228 of file MoochoPack_MoochoSolver.cpp.

const MoochoSolver::options_ptr_t & MoochoPack::MoochoSolver::get_options ( ) const

Get the OptionsFromStream object being used to extract the options from.

ToDo: Finish documentation.

Definition at line 241 of file MoochoPack_MoochoSolver.cpp.

void MoochoPack::MoochoSolver::set_error_handling ( bool  throw_exceptions,
const ostream_ptr_t error_out 
)

Set the error output and whether exceptions will be thrown from these functions or not.

Parameters
throw_exceptions[in] If true, then after printing the error message (see error_out) the exception will be rethrown out of this->solve_nlp().
error_out[in] If error_out.get() != NULL, then the error messages from any thrown std::exception will be printed to *error_out. Otherwise, they will be printed to std::cerr.

Postconditions:

Definition at line 246 of file MoochoPack_MoochoSolver.cpp.

bool MoochoPack::MoochoSolver::throw_exceptions ( ) const

Return if exceptions will be thrown out of this->solve_nlp().

Definition at line 265 of file MoochoPack_MoochoSolver.cpp.

const MoochoSolver::ostream_ptr_t & MoochoPack::MoochoSolver::error_out ( ) const

Return the std::ostream object used for error reporting on exceptions.

If return.get() == NULL then this means that std::cerr wil be written to.

Definition at line 271 of file MoochoPack_MoochoSolver.cpp.

void MoochoPack::MoochoSolver::set_output_context ( const std::string &  file_context_postfix,
EOutputToBlackHole  output_to_black_hole = OUTPUT_TO_BLACK_HOLE_DEFAULT,
const int  procRank = -1,
const int  numProcs = -1 
)

Setup the context for outputting.

Parameters
file_context_postfix[in] Prefix to attach to every file name for default output files.
output_to_black_hole[in] Determines if output will be sent to black-hole streams or not.
procRank[in] If procRank >= 0 then this this will be considered to be the rank of this process. Otherwise, procRank will be extracted from Teuchos::GlobalMPISession::getRank(). This integer is used to build file names for output streams that get created.
numProcs[in] If numProcs > 0 then this this will be considered to be the number of processes. Otherwise, numProcs will be extracted from Teuchos::GlobalMPISession::getNProc(). This integer is used to build file names for output streams that get created.

Note that this function will wipe out all storage for all of the output streams that are set. Therefore, you should call this function before you set any streams manually.

This information affects the file names created by the generate_output_file() function which is used both by external clients and also internally to create various output files.

Definition at line 140 of file MoochoPack_MoochoSolver.cpp.

void MoochoPack::MoochoSolver::set_output_file_tag ( const std::string &  output_file_tag)
inline

Set a tag for output file names for all file names that are created internally.

Definition at line 835 of file MoochoPack_MoochoSolver.hpp.

void MoochoPack::MoochoSolver::do_console_outputting ( bool  do_console_outputting)
inline

Turn on and off console outputting.

Definition at line 841 of file MoochoPack_MoochoSolver.hpp.

bool MoochoPack::MoochoSolver::do_console_outputting ( ) const
inline

Return if console outputting is performed or not.

Definition at line 847 of file MoochoPack_MoochoSolver.hpp.

void MoochoPack::MoochoSolver::set_console_out ( const ostream_ptr_t console_out)

Set the std::ostream object to use for console output by a MoochoTrackerConsoleStd object.

Parameters
console_out[in] Smart pointer to an std::ostream object that output appropriate for the console is sent to. An obvious choice for this output stream is std::cout can can be set using this->set_console_out(rcp(&std::cout,false)).

Postconditions:

ToDo: Discuss exactly what is printed to this output stream.

Definition at line 276 of file MoochoPack_MoochoSolver.cpp.

const MoochoSolver::ostream_ptr_t & MoochoPack::MoochoSolver::get_console_out ( ) const

Get the non-const smart pointer to the set output stream for console outputting.

Returns
Returns the console_out smart pointer object pointer set by the last call to this->set_console_out(). Not that if console_out.get() == NULL on in the last call to this->set_console_out(console_out) this this method returns return.get() == NULL which is a flag that the stream std::cout is being written to and this function does not provide access to that std::ofstream object (the client can access that stream themselves).

Definition at line 284 of file MoochoPack_MoochoSolver.cpp.

void MoochoPack::MoochoSolver::do_summary_outputting ( bool  do_summary_outputting)
inline

Turn on and off summary outputting.

Definition at line 853 of file MoochoPack_MoochoSolver.hpp.

bool MoochoPack::MoochoSolver::do_summary_outputting ( ) const
inline

Return if summary outputting is performed or not.

Definition at line 859 of file MoochoPack_MoochoSolver.hpp.

void MoochoPack::MoochoSolver::set_summary_out ( const ostream_ptr_t summary_out)

Set the std::ostream object to use for summary output.

Parameters
summay_out[in] Smart pointer to std::ostream object that summary output is sent to.

Postconditions:

ToDo: Discuss exactly what is printed to this output stream.

Definition at line 289 of file MoochoPack_MoochoSolver.cpp.

const MoochoSolver::ostream_ptr_t & MoochoPack::MoochoSolver::get_summary_out ( ) const

Get the non-const smart pointer to the set output stream for summary outputting.

Postconditions:

  • return.get() == summary_out.get() where summary_out was the value last sent to this->set_summary_out(summary_out).
Returns
Returns the summary_out smart pointer object pointer set by the last call to this->set_summary_out(). Not that if summary_out.get() == NULL on in the last call to this->set_summary_out(summary_out) this this method returns return.get() == NULL which is a flag that the file "MoochoSummary.out" is being written to and this function does not provide access to that std::ofstream object.

Definition at line 297 of file MoochoPack_MoochoSolver.cpp.

void MoochoPack::MoochoSolver::do_journal_outputting ( bool  do_journal_outputting)
inline

Turn on and off journal outputting.

Definition at line 865 of file MoochoPack_MoochoSolver.hpp.

bool MoochoPack::MoochoSolver::do_journal_outputting ( ) const
inline

Return if journal outputting is performed or not.

Definition at line 871 of file MoochoPack_MoochoSolver.hpp.

void MoochoPack::MoochoSolver::set_journal_out ( const ostream_ptr_t journal_out)

Set the std::ostream object to use for journal output by the MOOCHO step objects.

Parameters
journal_out[in] Smart pointer to an std::ostream object that journal output will be set to.

Note that if the option NLPSolverClientInterface::journal_output_level == PRINT_NOTHING in the last call to this->set_options(options) then no output is sent to this stream at all.

Postconditions:

  • [journal_out.get() == NULL] The file "MoochoJournal.out" will be overwritten with journal output. object the next time this->solve_nlp() is called.
  • [journal_out.get() != NULL] The stream *journal_out will be written to with journal output the next time this->solve_nlp() is called.
  • this->get_journal_out().get() == journal_out.get()

Definition at line 302 of file MoochoPack_MoochoSolver.cpp.

const MoochoSolver::ostream_ptr_t & MoochoPack::MoochoSolver::get_journal_out ( ) const

Get the non-const smart pointer to the set output stream for journal outputting.

Postconditions:

  • return.get() == journal_out.get() where journal_out was the value last sent to this->set_journal_out(journal_out).
Returns
Returns the journal_out smart pointer object pointer set by the last call to this->set_journal_out(). Not that if journal_out.get() == NULL on in the last call to this->set_journal_out(journal_out) this this method returns return.get() == NULL which is a flag that the file "MoochoJournal.out" is being written to and this function does not provide access to that std::ofstream object.

Definition at line 310 of file MoochoPack_MoochoSolver.cpp.

void MoochoPack::MoochoSolver::do_algo_outputting ( bool  do_algo_outputting)
inline

Turn on and off algo outputting.

Definition at line 877 of file MoochoPack_MoochoSolver.hpp.

bool MoochoPack::MoochoSolver::do_algo_outputting ( ) const
inline

Return if algo outputting is performed or not.

Definition at line 883 of file MoochoPack_MoochoSolver.hpp.

void MoochoPack::MoochoSolver::generate_stats_file ( bool  generate_stats_file)
inline

Turn on and off the generation of a statistics file.

Definition at line 889 of file MoochoPack_MoochoSolver.hpp.

bool MoochoPack::MoochoSolver::generate_stats_file ( ) const
inline

Return if a statistics file will be generated or not.

Definition at line 895 of file MoochoPack_MoochoSolver.hpp.

void MoochoPack::MoochoSolver::set_algo_out ( const ostream_ptr_t algo_out)

Set the std::ostream object to use for algorithm output.

Parameters
algo_out[in] Smart pointer to an std::ostream object that static algorithm output will be set to.

Postconditions:

  • [algo_out.get() == NULL] The file "MoochoAlgo.out" will be overwritten with algorithm info. object the next time this->solve_nlp() is called.
  • [algo_out.get() != NULL] The stream *algo_out will be written to with algorithm info the next time this->solve_nlp() is called.
  • this->get_algo_out().get() == algo_out.get()

Note that if the option NLPSolverClientInterface::print_algo == false in the last call to this->set_options(options) then no output is sent to this stream at all.

ToDo: Discuss exactly what is printed to this output stream.

Definition at line 315 of file MoochoPack_MoochoSolver.cpp.

const MoochoSolver::ostream_ptr_t & MoochoPack::MoochoSolver::get_algo_out ( ) const

Get the non-const smart pointer to the set output stream for algo outputting.

Postconditions:

  • return.get() == algo_out.get() where algo_out was the value last sent to this->set_algo_out(algo_out).
Returns
Returns the algo_out smart pointer object pointer set by the last call to this->set_algo_out(). Not that if algo_out.get() == NULL on in the last call to this->set_algo_out(algo_out) this this method returns return.get() == NULL which is a flag that the file "MoochoAlgo.out" is being written to and this function does not provide access to that std::ofstream object.

Definition at line 322 of file MoochoPack_MoochoSolver.cpp.

Teuchos::RCP< std::ostream > MoochoPack::MoochoSolver::generate_output_file ( const std::string &  fileNameBase) const

Generate an output file given a base file name.

Note that this will typically only create a ofsteam object on the root process and a oblackholestream object on all other processes.

Definition at line 328 of file MoochoPack_MoochoSolver.cpp.

void MoochoPack::MoochoSolver::update_solver ( ) const

Setup the state of the solver and get ready for a solve.

This function gets called already by solve_nlp() so it does not have to be called manually. However, after calling this function all stream objects will be set up in order to be used prior to calling solve_nlp().

Definition at line 345 of file MoochoPack_MoochoSolver.cpp.

MoochoSolver::ESolutionStatus MoochoPack::MoochoSolver::solve_nlp ( ) const

Solve the NLP.

Preconditions:

  • this->get_nlp() != NULL (throw std::logic_error)

Algorithm configuration:
The NLPAlgoConfig object used to configure an optimization algorithm is specified by *this->get_config(). If this->get_config().get() == NULL then the default configuratiion class NLPAlgoConfigMamaJama is used. This default configuration class is full featured and should have a good shot at solving any NLP thrown at it.

Specifying options:
The options used by the configuration object to configure the optimization algorithm as well as solver tolerances, maximum number of iterations etc are taken from the OptionsFromStreamPack::OptionsFromStream object returned from *this->get_options(). If this->get_options().get() == NULL then an attempt is made to open the file 'Moocho.opt' in the current directory. If this file does not exist, then a default set of options is used which will be acceptable for most NLPs. The files Moocho.opt.MoochoSolver, Moocho.opt.DecompositionSystemStateStepBuilderStd and Moocho.opt.NLPAlgoConfigMamaJama show which options can be used with this solver interface and a NLPAlgoConfigMamaJama configuration object respectively. Other configuration classes will use a different set of options. See the documentation for those configuration classes for details.

Outputting to streams:

If this->throw_exceptions() == false then any exceptions that may be thown internally will be caught, std::exception::what() will be printed to *this->error_out() (or std::cerr if this->error_out().get() == NULL) and this method will return SOLVE_RETURN_EXCEPTION. If this->throw_exceptions() == true, then after the error has been reported, the exception will be rethrown out for the caller to deal with!

Even if no exception is thrown, then a short one-line summary message will be printed to *this->error_out() (or std::cerr if this->error_out().get() == NULL) stating whether the NLP was solved or not.

ToDo: Finish documentation!

Returns
The solution status:
  • ToDo: Fill these in and discuss them!

Definition at line 668 of file MoochoPack_MoochoSolver.cpp.

NLPSolverClientInterface & MoochoPack::MoochoSolver::get_solver ( )

Get the underlying NLPSolverClientInterface object.

If the algorithm has not already been configured it will be here using whatever NLPAlgoConfig object that it has (set using this->set_config*() or using the default). Whatever options where set (i.e. using this->set_options()) will be used when this algorithm is configured and the NLP is solved.

Preconditions:

  • this->get_nlp() != NULL (throw std::logic_error)

Warning! Do not try to modify the underlying solver object using this returned reference to the solver object. Instead, use the above public interface functions.

ToDo: Finish documentation!

Definition at line 918 of file MoochoPack_MoochoSolver.cpp.

const NLPSolverClientInterface & MoochoPack::MoochoSolver::get_solver ( ) const

Definition at line 924 of file MoochoPack_MoochoSolver.cpp.


The documentation for this class was generated from the following files: