ROL
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ROL::TypeG::Algorithm< Real > Class Template Referenceabstract

Provides an interface to run general constrained optimization algorithms. More...

#include <ROL_TypeG_Algorithm.hpp>

+ Inheritance diagram for ROL::TypeG::Algorithm< Real >:

Public Member Functions

virtual ~Algorithm ()
 
 Algorithm ()
 Constructor, given a step and a status test. More...
 
void setStatusTest (const Ptr< StatusTest< Real >> &status, bool combineStatus=false)
 
virtual void run (Problem< Real > &problem, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &econ, Vector< Real > &emul, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, Constraint< Real > &econ, Vector< Real > &emul, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &econ, Vector< Real > &emul, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &econ, Vector< Real > &emul, const Vector< Real > &eres, std::ostream &outStream=std::cout)=0
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, const Vector< Real > &ires, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, const Vector< Real > &ires, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, Constraint< Real > &econ, Vector< Real > &emul, const Vector< Real > &eres, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, const Vector< Real > &ires, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &econ, Vector< Real > &emul, const Vector< Real > &eres, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, const Vector< Real > &ires, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &econ, Vector< Real > &emul, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, Constraint< Real > &econ, Vector< Real > &emul, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &econ, Vector< Real > &emul, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &econ, Vector< Real > &emul, const Vector< Real > &eres, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, const Vector< Real > &linear_eres, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, const Vector< Real > &ires, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, const Vector< Real > &linear_eres, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, const Vector< Real > &ires, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, const Vector< Real > &linear_eres, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, Constraint< Real > &econ, Vector< Real > &emul, const Vector< Real > &eres, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, const Vector< Real > &ires, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, const Vector< Real > &linear_eres, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &econ, Vector< Real > &emul, const Vector< Real > &eres, Constraint< Real > &icon, Vector< Real > &imul, BoundConstraint< Real > &ibnd, const Vector< Real > &ires, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, const Vector< Real > &linear_eres, std::ostream &outStream=std::cout)
 Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface. More...
 
virtual void writeHeader (std::ostream &os) const
 Print iterate header. More...
 
virtual void writeName (std::ostream &os) const
 Print step name. More...
 
virtual void writeOutput (std::ostream &os, const bool write_header=false) const
 Print iterate status. More...
 
virtual void writeExitStatus (std::ostream &os) const
 
Ptr< const AlgorithmState< Real > > getState () const
 
void reset ()
 

Protected Member Functions

void initialize (const Vector< Real > &x, const Vector< Real > &g, const Vector< Real > &mul, const Vector< Real > &c)
 

Protected Attributes

const Ptr< CombinedStatusTest
< Real > > 
status_
 
const Ptr< AlgorithmState< Real > > state_
 
Ptr< PolyhedralProjection< Real > > proj_
 

Detailed Description

template<typename Real>
class ROL::TypeG::Algorithm< Real >

Provides an interface to run general constrained optimization algorithms.

Definition at line 56 of file ROL_TypeG_Algorithm.hpp.

Constructor & Destructor Documentation

template<typename Real >
virtual ROL::TypeG::Algorithm< Real >::~Algorithm ( )
inlinevirtual

Definition at line 69 of file ROL_TypeG_Algorithm.hpp.

template<typename Real >
ROL::TypeG::Algorithm< Real >::Algorithm ( )

Constructor, given a step and a status test.

Definition at line 23 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::TypeG::Algorithm< Real >::status_.

Member Function Documentation

template<typename Real >
void ROL::TypeG::Algorithm< Real >::initialize ( const Vector< Real > &  x,
const Vector< Real > &  g,
const Vector< Real > &  mul,
const Vector< Real > &  c 
)
protected
template<typename Real >
void ROL::TypeG::Algorithm< Real >::setStatusTest ( const Ptr< StatusTest< Real >> &  status,
bool  combineStatus = false 
)

Definition at line 65 of file ROL_TypeG_Algorithm_Def.hpp.

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Problem< Real > &  problem,
std::ostream &  outStream = std::cout 
)
virtual
template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
Constraint< Real > &  econ,
Vector< Real > &  emul,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 94 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addBoundConstraint(), ROL::Problem< Real >::addConstraint(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
Objective< Real > &  obj,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 111 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addConstraint(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 127 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addBoundConstraint(), ROL::Problem< Real >::addConstraint(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
Objective< Real > &  obj,
Constraint< Real > &  econ,
Vector< Real > &  emul,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 145 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addConstraint(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
Constraint< Real > &  econ,
Vector< Real > &  emul,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 165 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addBoundConstraint(), ROL::Problem< Real >::addConstraint(), and ROL::Problem< Real >::finalize().

template<typename Real >
virtual void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
const Vector< Real > &  g,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
Constraint< Real > &  econ,
Vector< Real > &  emul,
const Vector< Real > &  eres,
std::ostream &  outStream = std::cout 
)
pure virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Implemented in ROL::TypeG::StabilizedLCLAlgorithm< Real >, ROL::TypeG::AugmentedLagrangianAlgorithm< Real >, ROL::TypeG::InteriorPointAlgorithm< Real >, and ROL::TypeG::MoreauYosidaAlgorithm< Real >.

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
const Vector< Real > &  g,
Objective< Real > &  obj,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
const Vector< Real > &  ires,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 189 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addConstraint(), ROL::Vector< Real >::clone(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
const Vector< Real > &  g,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
const Vector< Real > &  ires,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 216 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addBoundConstraint(), ROL::Problem< Real >::addConstraint(), ROL::Vector< Real >::clone(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
const Vector< Real > &  g,
Objective< Real > &  obj,
Constraint< Real > &  econ,
Vector< Real > &  emul,
const Vector< Real > &  eres,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
const Vector< Real > &  ires,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 246 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addConstraint(), ROL::Vector< Real >::clone(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
const Vector< Real > &  g,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
Constraint< Real > &  econ,
Vector< Real > &  emul,
const Vector< Real > &  eres,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
const Vector< Real > &  ires,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 283 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addBoundConstraint(), ROL::Problem< Real >::addConstraint(), ROL::Vector< Real >::clone(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
Constraint< Real > &  econ,
Vector< Real > &  emul,
Constraint< Real > &  linear_econ,
Vector< Real > &  linear_emul,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 324 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addBoundConstraint(), ROL::Problem< Real >::addConstraint(), ROL::Problem< Real >::addLinearConstraint(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
Objective< Real > &  obj,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
Constraint< Real > &  linear_econ,
Vector< Real > &  linear_emul,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 347 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addConstraint(), ROL::Problem< Real >::addLinearConstraint(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
Constraint< Real > &  linear_econ,
Vector< Real > &  linear_emul,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 369 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addBoundConstraint(), ROL::Problem< Real >::addConstraint(), ROL::Problem< Real >::addLinearConstraint(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
Objective< Real > &  obj,
Constraint< Real > &  econ,
Vector< Real > &  emul,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
Constraint< Real > &  linear_econ,
Vector< Real > &  linear_emul,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 393 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addConstraint(), ROL::Problem< Real >::addLinearConstraint(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
Constraint< Real > &  econ,
Vector< Real > &  emul,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
Constraint< Real > &  linear_econ,
Vector< Real > &  linear_emul,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 419 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addBoundConstraint(), ROL::Problem< Real >::addConstraint(), ROL::Problem< Real >::addLinearConstraint(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
const Vector< Real > &  g,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
Constraint< Real > &  econ,
Vector< Real > &  emul,
const Vector< Real > &  eres,
Constraint< Real > &  linear_econ,
Vector< Real > &  linear_emul,
const Vector< Real > &  linear_eres,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 449 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addBoundConstraint(), ROL::Problem< Real >::addConstraint(), ROL::Problem< Real >::addLinearConstraint(), ROL::Vector< Real >::clone(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
const Vector< Real > &  g,
Objective< Real > &  obj,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
const Vector< Real > &  ires,
Constraint< Real > &  linear_econ,
Vector< Real > &  linear_emul,
const Vector< Real > &  linear_eres,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 478 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addConstraint(), ROL::Problem< Real >::addLinearConstraint(), ROL::Vector< Real >::clone(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
const Vector< Real > &  g,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
const Vector< Real > &  ires,
Constraint< Real > &  linear_econ,
Vector< Real > &  linear_emul,
const Vector< Real > &  linear_eres,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 512 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addBoundConstraint(), ROL::Problem< Real >::addConstraint(), ROL::Problem< Real >::addLinearConstraint(), ROL::Vector< Real >::clone(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
const Vector< Real > &  g,
Objective< Real > &  obj,
Constraint< Real > &  econ,
Vector< Real > &  emul,
const Vector< Real > &  eres,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
const Vector< Real > &  ires,
Constraint< Real > &  linear_econ,
Vector< Real > &  linear_emul,
const Vector< Real > &  linear_eres,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 551 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addConstraint(), ROL::Problem< Real >::addLinearConstraint(), ROL::Vector< Real >::clone(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::run ( Vector< Real > &  x,
const Vector< Real > &  g,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
Constraint< Real > &  econ,
Vector< Real > &  emul,
const Vector< Real > &  eres,
Constraint< Real > &  icon,
Vector< Real > &  imul,
BoundConstraint< Real > &  ibnd,
const Vector< Real > &  ires,
Constraint< Real > &  linear_econ,
Vector< Real > &  linear_emul,
const Vector< Real > &  linear_eres,
std::ostream &  outStream = std::cout 
)
virtual

Run algorithm on general constrained problems (Type-G). This is the primary Type-G interface.

Definition at line 591 of file ROL_TypeG_Algorithm_Def.hpp.

References ROL::Problem< Real >::addBoundConstraint(), ROL::Problem< Real >::addConstraint(), ROL::Problem< Real >::addLinearConstraint(), ROL::Vector< Real >::clone(), and ROL::Problem< Real >::finalize().

template<typename Real >
void ROL::TypeG::Algorithm< Real >::writeHeader ( std::ostream &  os) const
virtual
template<typename Real >
void ROL::TypeG::Algorithm< Real >::writeName ( std::ostream &  os) const
virtual
template<typename Real >
void ROL::TypeG::Algorithm< Real >::writeOutput ( std::ostream &  os,
const bool  write_header = false 
) const
virtual
template<typename Real >
void ROL::TypeG::Algorithm< Real >::writeExitStatus ( std::ostream &  os) const
virtual
template<typename Real >
Ptr< const AlgorithmState< Real > > ROL::TypeG::Algorithm< Real >::getState ( void  ) const

Definition at line 691 of file ROL_TypeG_Algorithm_Def.hpp.

template<typename Real >
void ROL::TypeG::Algorithm< Real >::reset ( void  )

Definition at line 697 of file ROL_TypeG_Algorithm_Def.hpp.

Member Data Documentation

template<typename Real >
const Ptr<CombinedStatusTest<Real> > ROL::TypeG::Algorithm< Real >::status_
protected
template<typename Real >
const Ptr<AlgorithmState<Real> > ROL::TypeG::Algorithm< Real >::state_
protected
template<typename Real >
Ptr<PolyhedralProjection<Real> > ROL::TypeG::Algorithm< Real >::proj_
protected

Definition at line 60 of file ROL_TypeG_Algorithm.hpp.


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