|
ROL
|
#include <ROL_StochasticProblem.hpp>
Inheritance diagram for ROL::StochasticProblem< Real >:Public Member Functions | |
| StochasticProblem (const Ptr< Objective< Real >> &obj, const Ptr< Vector< Real >> &x, const Ptr< Vector< Real >> &g=nullPtr) | |
| Default constructor for StochasticProblem. More... | |
| StochasticProblem (const Problem< Real > &problem) | |
| void | makeObjectiveStochastic (ParameterList &list, const Ptr< SampleGenerator< Real >> &fsampler, const Ptr< SampleGenerator< Real >> &gsampler=nullPtr, const Ptr< SampleGenerator< Real >> &hsampler=nullPtr) |
| void | makeObjectiveStochastic (const Ptr< RandVarFunctional< Real >> &rvf, ParameterList &list, const Ptr< SampleGenerator< Real >> &fsampler, const Ptr< SampleGenerator< Real >> &gsampler=nullPtr, const Ptr< SampleGenerator< Real >> &hsampler=nullPtr) |
| void | makeConstraintStochastic (std::string name, ParameterList &list, const Ptr< SampleGenerator< Real >> &sampler, const Ptr< BatchManager< Real >> &bman=nullPtr) |
| void | makeLinearConstraintStochastic (std::string name, ParameterList &list, const Ptr< SampleGenerator< Real >> &sampler, const Ptr< BatchManager< Real >> &bman=nullPtr) |
| void | resetStochasticObjective (void) |
| void | resetStochasticConstraint (std::string name) |
| void | resetStochasticLinearConstraint (std::string name) |
| void | resetStochastic (void) |
| std::vector< Real > | getObjectiveStatistic (void) const |
| std::vector< Real > | getConstraintStatistic (std::string name) const |
| Real | getSolutionStatistic (int comp=0, std::string name="") const |
| void | finalize (bool lumpConstraints=false, bool printToStream=false, std::ostream &outStream=std::cout) override |
| Tranform user-supplied constraints to consist of only bounds and equalities. Optimization problem cannot be modified after finalize has been called without calling the edit function. More... | |
| void | edit (void) override |
| Resume editting optimization problem after finalize has been called. More... | |
Public Member Functions inherited from ROL::Problem< Real > | |
| virtual | ~Problem () |
| Problem (const Ptr< Objective< Real >> &obj, const Ptr< Vector< Real >> &x, const Ptr< Vector< Real >> &g=nullPtr) | |
| Default constructor for OptimizationProblem. More... | |
| Problem (const Problem &problem) | |
| Copy constructor for OptimizationProblem. More... | |
| void | addBoundConstraint (const Ptr< BoundConstraint< Real >> &bnd) |
| Add a bound constraint. More... | |
| void | removeBoundConstraint () |
| Remove an existing bound constraint. More... | |
| void | addProxObjective (const Ptr< ProxObjective< Real >> &prox) |
| Add a prox objective. More... | |
| void | removeProxObjective () |
| Remove an existing prox objective. More... | |
| void | addConstraint (std::string name, const Ptr< Constraint< Real >> &econ, const Ptr< Vector< Real >> &emul, const Ptr< Vector< Real >> &eres=nullPtr, bool reset=false) |
| Add an equality constraint. More... | |
| void | addConstraint (std::string name, const Ptr< Constraint< Real >> &icon, const Ptr< Vector< Real >> &imul, const Ptr< BoundConstraint< Real >> &ibnd, const Ptr< Vector< Real >> &ires=nullPtr, bool reset=false) |
| Add an inequality constraint. More... | |
| void | removeConstraint (std::string name) |
| Remove an existing constraint. More... | |
| void | addLinearConstraint (std::string name, const Ptr< Constraint< Real >> &linear_econ, const Ptr< Vector< Real >> &linear_emul, const Ptr< Vector< Real >> &linear_eres=nullPtr, bool reset=false) |
| Add a linear equality constraint. More... | |
| void | addLinearConstraint (std::string name, const Ptr< Constraint< Real >> &linear_icon, const Ptr< Vector< Real >> &linear_imul, const Ptr< BoundConstraint< Real >> &linear_ibnd, const Ptr< Vector< Real >> &linear_ires=nullPtr, bool reset=false) |
| Add a linear inequality constraint. More... | |
| void | removeLinearConstraint (std::string name) |
| Remove an existing linear constraint. More... | |
| void | setProjectionAlgorithm (ParameterList &parlist) |
| Set polyhedral projection algorithm. More... | |
| const Ptr< Objective< Real > > & | getObjective () |
| Get the objective function. More... | |
| const Ptr< Vector< Real > > & | getPrimalOptimizationVector () |
| Get the primal optimization space vector. More... | |
| const Ptr< Vector< Real > > & | getDualOptimizationVector () |
| Get the dual optimization space vector. More... | |
| const Ptr< BoundConstraint < Real > > & | getBoundConstraint () |
| Get the bound constraint. More... | |
| const Ptr< Constraint< Real > > & | getConstraint () |
| Get the equality constraint. More... | |
| const Ptr< Vector< Real > > & | getMultiplierVector () |
| Get the dual constraint space vector. More... | |
| const Ptr< Vector< Real > > & | getResidualVector () |
| Get the primal constraint space vector. More... | |
| const Ptr < PolyhedralProjection< Real > > & | getPolyhedralProjection () |
| Get the polyhedral projection object. This is a null pointer if no linear constraints and/or bounds are present. More... | |
| EProblem | getProblemType () |
| Get the optimization problem type (U, B, E, or G). More... | |
| Real | checkLinearity (bool printToStream=false, std::ostream &outStream=std::cout) const |
| Check if user-supplied linear constraints are affine. More... | |
| void | checkVectors (bool printToStream=false, std::ostream &outStream=std::cout) const |
| Run vector checks for user-supplied vectors. More... | |
| void | checkDerivatives (bool printToStream=false, std::ostream &outStream=std::cout) const |
| Run derivative checks for user-supplied objective function and constraints. More... | |
| void | check (bool printToStream=false, std::ostream &outStream=std::cout) const |
| Run vector, linearity and derivative checks for user-supplied vectors, objective function and constraints. More... | |
| bool | isFinalized () const |
| Indicate whether or no finalize has been called. More... | |
| void | finalizeIteration () |
| Transform the optimization variables to the native parameterization after an optimization algorithm has finished. More... | |
| virtual | ~Problem () |
| Problem (const Ptr< Objective< Real >> &obj, const Ptr< Vector< Real >> &x, const Ptr< Vector< Real >> &g=nullPtr) | |
| Default constructor for OptimizationProblem. More... | |
| Problem (const Problem &problem) | |
| Copy constructor for OptimizationProblem. More... | |
| void | addBoundConstraint (const Ptr< BoundConstraint< Real >> &bnd) |
| Add a bound constraint. More... | |
| void | removeBoundConstraint () |
| Remove an existing bound constraint. More... | |
| void | addConstraint (std::string name, const Ptr< Constraint< Real >> &econ, const Ptr< Vector< Real >> &emul, const Ptr< Vector< Real >> &eres=nullPtr, bool reset=false) |
| Add an equality constraint. More... | |
| void | addConstraint (std::string name, const Ptr< Constraint< Real >> &icon, const Ptr< Vector< Real >> &imul, const Ptr< BoundConstraint< Real >> &ibnd, const Ptr< Vector< Real >> &ires=nullPtr, bool reset=false) |
| Add an inequality constraint. More... | |
| void | removeConstraint (std::string name) |
| Remove an existing constraint. More... | |
| void | addLinearConstraint (std::string name, const Ptr< Constraint< Real >> &linear_econ, const Ptr< Vector< Real >> &linear_emul, const Ptr< Vector< Real >> &linear_eres=nullPtr, bool reset=false) |
| Add a linear equality constraint. More... | |
| void | addLinearConstraint (std::string name, const Ptr< Constraint< Real >> &linear_icon, const Ptr< Vector< Real >> &linear_imul, const Ptr< BoundConstraint< Real >> &linear_ibnd, const Ptr< Vector< Real >> &linear_ires=nullPtr, bool reset=false) |
| Add a linear inequality constraint. More... | |
| void | removeLinearConstraint (std::string name) |
| Remove an existing linear constraint. More... | |
| void | setProjectionAlgorithm (ParameterList &parlist) |
| Set polyhedral projection algorithm. More... | |
| void | addProximableObjective (const Ptr< Objective< Real >> &nobj) |
| void | removeProximableObjective () |
| const Ptr< Objective< Real > > & | getObjective () |
| Get the objective function. More... | |
| const Ptr< Objective< Real > > & | getProximableObjective () |
| Get proximable objective. More... | |
| const Ptr< Vector< Real > > & | getPrimalOptimizationVector () |
| Get the primal optimization space vector. More... | |
| const Ptr< Vector< Real > > & | getDualOptimizationVector () |
| Get the dual optimization space vector. More... | |
| const Ptr< BoundConstraint < Real > > & | getBoundConstraint () |
| Get the bound constraint. More... | |
| const Ptr< Constraint< Real > > & | getConstraint () |
| Get the equality constraint. More... | |
| const Ptr< Vector< Real > > & | getMultiplierVector () |
| Get the dual constraint space vector. More... | |
| const Ptr< Vector< Real > > & | getResidualVector () |
| Get the primal constraint space vector. More... | |
| const Ptr < PolyhedralProjection< Real > > & | getPolyhedralProjection () |
| Get the polyhedral projection object. This is a null pointer if no linear constraints and/or bounds are present. More... | |
| EProblem | getProblemType () |
| Get the optimization problem type (U, B, E, G, or P). More... | |
| Real | checkLinearity (bool printToStream=false, std::ostream &outStream=std::cout) const |
| Check if user-supplied linear constraints are affine. More... | |
| void | checkVectors (bool printToStream=false, std::ostream &outStream=std::cout) const |
| Run vector checks for user-supplied vectors. More... | |
| void | checkDerivatives (bool printToStream=false, std::ostream &outStream=std::cout, const Ptr< Vector< Real >> &x0=nullPtr, Real scale=Real(1)) const |
| Run derivative checks for user-supplied objective function and constraints. More... | |
| void | check (bool printToStream=false, std::ostream &outStream=std::cout, const Ptr< Vector< Real >> &x0=nullPtr, Real scale=Real(1)) const |
| Run vector, linearity and derivative checks for user-supplied vectors, objective function and constraints. More... | |
| bool | isFinalized () const |
| Indicate whether or no finalize has been called. More... | |
| void | finalizeIteration () |
| Transform the optimization variables to the native parameterization after an optimization algorithm has finished. More... | |
Private Attributes | |
| Ptr< Objective< Real > > | ORIGINAL_obj_ |
| Ptr< Objective< Real > > | ORIGINAL_nobj_ |
| Ptr< Vector< Real > > | ORIGINAL_xprim_ |
| Ptr< Vector< Real > > | ORIGINAL_xdual_ |
| Ptr< BoundConstraint< Real > > | ORIGINAL_bnd_ |
| std::unordered_map < std::string, ConstraintData < Real > > | ORIGINAL_con_ |
| std::unordered_map < std::string, ConstraintData < Real > > | ORIGINAL_linear_con_ |
| bool | needRiskLessObj_ |
| bool | risk_ |
| std::vector< bool > | needRiskLessCon_ |
| Ptr< ParameterList > | objList_ |
| std::unordered_map < std::string, std::pair< Ptr < ParameterList >, bool > > | conList_ |
| std::unordered_map < std::string, size_t > | statMap_ |
Additional Inherited Members | |
Protected Attributes inherited from ROL::Problem< Real > | |
| Ptr< Objective< Real > > | INPUT_obj_ |
| Ptr< Vector< Real > > | INPUT_xprim_ |
| Ptr< Vector< Real > > | INPUT_xdual_ |
| Ptr< BoundConstraint< Real > > | INPUT_bnd_ |
| std::unordered_map < std::string, ConstraintData < Real > > | INPUT_con_ |
| std::unordered_map < std::string, ConstraintData < Real > > | INPUT_linear_con_ |
| Ptr< ProxObjective< Real > > | INPUT_prox_ |
| Ptr< Objective< Real > > | INPUT_nobj_ |
Definition at line 33 of file ROL_StochasticProblem.hpp.
| ROL::StochasticProblem< Real >::StochasticProblem | ( | const Ptr< Objective< Real >> & | obj, |
| const Ptr< Vector< Real >> & | x, | ||
| const Ptr< Vector< Real >> & | g = nullPtr |
||
| ) |
Default constructor for StochasticProblem.
| [in] | obj | objective function object |
| [in] | x | primal optimization space vector |
| [in] | g | dual optimization space vector |
Definition at line 16 of file ROL_StochasticProblem_Def.hpp.
|
inline |
Definition at line 69 of file ROL_StochasticProblem.hpp.
| void ROL::StochasticProblem< Real >::makeObjectiveStochastic | ( | ParameterList & | list, |
| const Ptr< SampleGenerator< Real >> & | fsampler, | ||
| const Ptr< SampleGenerator< Real >> & | gsampler = nullPtr, |
||
| const Ptr< SampleGenerator< Real >> & | hsampler = nullPtr |
||
| ) |
Definition at line 22 of file ROL_StochasticProblem_Def.hpp.
| void ROL::StochasticProblem< Real >::makeObjectiveStochastic | ( | const Ptr< RandVarFunctional< Real >> & | rvf, |
| ParameterList & | list, | ||
| const Ptr< SampleGenerator< Real >> & | fsampler, | ||
| const Ptr< SampleGenerator< Real >> & | gsampler = nullPtr, |
||
| const Ptr< SampleGenerator< Real >> & | hsampler = nullPtr |
||
| ) |
Definition at line 65 of file ROL_StochasticProblem_Def.hpp.
| void ROL::StochasticProblem< Real >::makeConstraintStochastic | ( | std::string | name, |
| ParameterList & | list, | ||
| const Ptr< SampleGenerator< Real >> & | sampler, | ||
| const Ptr< BatchManager< Real >> & | bman = nullPtr |
||
| ) |
Definition at line 94 of file ROL_StochasticProblem_Def.hpp.
References ROL::Problem< Real >::addConstraint(), and ROL::Problem< Real >::removeConstraint().
| void ROL::StochasticProblem< Real >::makeLinearConstraintStochastic | ( | std::string | name, |
| ParameterList & | list, | ||
| const Ptr< SampleGenerator< Real >> & | sampler, | ||
| const Ptr< BatchManager< Real >> & | bman = nullPtr |
||
| ) |
Definition at line 160 of file ROL_StochasticProblem_Def.hpp.
References ROL::Problem< Real >::addLinearConstraint(), and ROL::Problem< Real >::removeLinearConstraint().
| void ROL::StochasticProblem< Real >::resetStochasticObjective | ( | void | ) |
Definition at line 214 of file ROL_StochasticProblem_Def.hpp.
| void ROL::StochasticProblem< Real >::resetStochasticConstraint | ( | std::string | name | ) |
Definition at line 226 of file ROL_StochasticProblem_Def.hpp.
References ROL::Problem< Real >::addConstraint(), and ROL::Problem< Real >::removeConstraint().
| void ROL::StochasticProblem< Real >::resetStochasticLinearConstraint | ( | std::string | name | ) |
Definition at line 244 of file ROL_StochasticProblem_Def.hpp.
References ROL::Problem< Real >::addLinearConstraint(), and ROL::Problem< Real >::removeLinearConstraint().
| void ROL::StochasticProblem< Real >::resetStochastic | ( | void | ) |
Definition at line 261 of file ROL_StochasticProblem_Def.hpp.
| std::vector< Real > ROL::StochasticProblem< Real >::getObjectiveStatistic | ( | void | ) | const |
Definition at line 300 of file ROL_StochasticProblem_Def.hpp.
| std::vector< Real > ROL::StochasticProblem< Real >::getConstraintStatistic | ( | std::string | name | ) | const |
Definition at line 315 of file ROL_StochasticProblem_Def.hpp.
| Real ROL::StochasticProblem< Real >::getSolutionStatistic | ( | int | comp = 0, |
| std::string | name = "" |
||
| ) | const |
Definition at line 333 of file ROL_StochasticProblem_Def.hpp.
|
overridevirtual |
Tranform user-supplied constraints to consist of only bounds and equalities. Optimization problem cannot be modified after finalize has been called without calling the edit function.
| [in] | lumpConstraints | combine both linear and nonlinear constraints |
| [in] | printToStream | determines whether to print to the supplied std::ostream |
| [in,out] | outStream | user supplied std::ostream |
Reimplemented from ROL::Problem< Real >.
Definition at line 365 of file ROL_StochasticProblem_Def.hpp.
References ROL::Problem< Real >::addConstraint(), ROL::Problem< Real >::addLinearConstraint(), ROL::Problem< Real >::finalize(), ROL::Problem< Real >::removeConstraint(), and ROL::Problem< Real >::removeLinearConstraint().
|
overridevirtual |
Resume editting optimization problem after finalize has been called.
Reimplemented from ROL::Problem< Real >.
Definition at line 468 of file ROL_StochasticProblem_Def.hpp.
References ROL::Problem< Real >::addConstraint(), ROL::Problem< Real >::addLinearConstraint(), ROL::Problem< Real >::edit(), ROL::Problem< Real >::removeConstraint(), and ROL::Problem< Real >::removeLinearConstraint().
|
private |
Definition at line 35 of file ROL_StochasticProblem.hpp.
|
private |
Definition at line 36 of file ROL_StochasticProblem.hpp.
|
private |
Definition at line 37 of file ROL_StochasticProblem.hpp.
|
private |
Definition at line 38 of file ROL_StochasticProblem.hpp.
|
private |
Definition at line 39 of file ROL_StochasticProblem.hpp.
|
private |
Definition at line 40 of file ROL_StochasticProblem.hpp.
|
private |
Definition at line 41 of file ROL_StochasticProblem.hpp.
|
private |
Definition at line 43 of file ROL_StochasticProblem.hpp.
|
private |
Definition at line 43 of file ROL_StochasticProblem.hpp.
|
private |
Definition at line 44 of file ROL_StochasticProblem.hpp.
|
private |
Definition at line 45 of file ROL_StochasticProblem.hpp.
|
private |
Definition at line 46 of file ROL_StochasticProblem.hpp.
|
private |
Definition at line 47 of file ROL_StochasticProblem.hpp.
1.8.5