49 #ifndef ROL_HELPERFUNCTIONS_HPP
50 #define ROL_HELPERFUNCTIONS_HPP
56 #include "Teuchos_SerialDenseMatrix.hpp"
57 #include "Teuchos_SerialDenseVector.hpp"
58 #include "Teuchos_LAPACK.hpp"
68 Teuchos::SerialDenseMatrix<int, Real> H(dim, dim);
70 Teuchos::RCP<Vector<Real> > e = x.
clone();
71 Teuchos::RCP<Vector<Real> > h = x.
clone();
73 for (
int i=0; i<dim; i++) {
76 for (
int j=0; j<dim; j++) {
91 Teuchos::SerialDenseMatrix<int, Real> M(dim, dim);
93 Teuchos::RCP<Vector<Real> > ei = x.
clone();
94 Teuchos::RCP<Vector<Real> > ej = x.
clone();
96 for (
int i=0; i<dim; i++) {
98 for (
int j=0; j<dim; j++) {
100 M(j,i) = ej->dot(*ei);
109 std::vector<std::vector<Real> >
computeEigenvalues(
const Teuchos::SerialDenseMatrix<int, Real> & mat) {
111 Teuchos::LAPACK<int, Real> lapack;
113 Teuchos::SerialDenseMatrix<int, Real> mymat(Teuchos::Copy, mat);
118 int n = mat.numRows();
120 std::vector<Real> real(n, 0);
121 std::vector<Real> imag(n, 0);
122 std::vector<std::vector<Real> > eigenvals;
132 std::vector<Real> work(lwork, 0);
136 lapack.GEEV(jobvl, jobvr, n, &mymat(0,0), n, &real[0], &imag[0], vl, ldvl, vr, ldvr, &work[0], lwork, &info);
138 eigenvals.push_back(real);
139 eigenvals.push_back(imag);
148 const Teuchos::SerialDenseMatrix<int, Real> & B) {
150 Teuchos::LAPACK<int, Real> lapack;
152 Teuchos::SerialDenseMatrix<int, Real> myA(Teuchos::Copy, A);
153 Teuchos::SerialDenseMatrix<int, Real> myB(Teuchos::Copy, B);
160 std::vector<Real> real(n, 0);
161 std::vector<Real> imag(n, 0);
162 std::vector<Real> beta(n, 0);
163 std::vector<std::vector<Real> > eigenvals;
173 std::vector<Real> work(lwork, 0);
177 lapack.GGEV(jobvl, jobvr, n, &myA(0,0), n, &myB(0,0), n, &real[0], &imag[0], &beta[0],
178 vl, ldvl, vr, ldvr, &work[0], lwork, &info);
180 for (
int i=0; i<n; i++) {
185 eigenvals.push_back(real);
186 eigenvals.push_back(imag);
194 Teuchos::SerialDenseMatrix<int, Real>
computeInverse(
const Teuchos::SerialDenseMatrix<int, Real> & mat) {
196 Teuchos::LAPACK<int, Real> lapack;
198 Teuchos::SerialDenseMatrix<int, Real> mymat(Teuchos::Copy, mat);
200 int n = mat.numRows();
202 std::vector<int> ipiv(n, 0);
206 std::vector<Real> work(lwork, 0);
210 lapack.GETRF(n, n, &mymat(0,0), n, &ipiv[0], &info);
211 lapack.GETRI(n, &mymat(0,0), n, &ipiv[0], &work[0], lwork, &info);
222 Teuchos::RCP<Objective<Real> >
obj_;
223 Teuchos::RCP<BoundConstraint<Real> >
con_;
232 bool useSecantPrecond =
false,
bool useSecantHessVec =
false, Real eps = 0.0 ) {
233 obj_ = Teuchos::rcp(&obj,
false);
234 con_ = Teuchos::rcp(&con,
false);
242 this->
obj_->update(x,flag,iter);
243 this->
con_->update(x,flag,iter);
247 return this->
obj_->value(x,tol);
251 this->
obj_->gradient(g,x,tol);
255 return this->
obj_->dirDeriv(x,d,tol);
260 this->
secant_->applyB( Hv, v, x );
263 this->
obj_->hessVec( Hv, v, x, tol );
272 this->
obj_->invHessVec(Hv,v,x,tol);
278 this->
secant_->applyH( Mv, v, x );
281 this->
obj_->precond( Mv, v, x, tol );
298 if ( this->
con_->isActivated() ) {
299 Teuchos::RCP<Vector<Real> > vnew = x.
clone();
303 this->
con_->pruneActive(*vnew,d,p,this->
eps_);
307 this->
con_->pruneActive(Hv,d,p,this->
eps_);
311 this->
con_->pruneInactive(*vnew,d,p,this->
eps_);
332 if ( this->
con_->isActivated() ) {
333 Teuchos::RCP<Vector<Real> > vnew = x.
clone();
337 this->
con_->pruneActive(*vnew,p,this->
eps_);
341 this->
con_->pruneActive(Hv,p,this->
eps_);
345 this->
con_->pruneInactive(*vnew,p,this->
eps_);
367 if ( this->
con_->isActivated() ) {
368 Teuchos::RCP<Vector<Real> > vnew = x.
clone();
372 this->
con_->pruneActive(*vnew,d,p,this->
eps_);
376 this->
con_->pruneActive(Hv,d,p,this->
eps_);
380 this->
con_->pruneInactive(*vnew,d,p,this->
eps_);
401 if ( this->
con_->isActivated() ) {
402 Teuchos::RCP<Vector<Real> > vnew = x.
clone();
406 this->
con_->pruneActive(*vnew,p,this->
eps_);
410 this->
con_->pruneActive(Hv,p,this->
eps_);
414 this->
con_->pruneInactive(*vnew,p,this->
eps_);
436 if ( this->
con_->isActivated() ) {
437 Teuchos::RCP<Vector<Real> > vnew = x.
clone();
441 this->
con_->pruneActive(*vnew,d,p,this->
eps_);
445 this->
con_->pruneActive(Mv,d,p,this->
eps_);
449 this->
con_->pruneInactive(*vnew,d,p,this->
eps_);
470 if ( this->
con_->isActivated() ) {
471 Teuchos::RCP<Vector<Real> > vnew = x.
clone();
475 this->
con_->pruneActive(*vnew,p,this->
eps_);
479 this->
con_->pruneActive(Mv,p,this->
eps_);
483 this->
con_->pruneInactive(*vnew,p,this->
eps_);
493 this->
con_->project(x);
497 this->
con_->pruneActive(v,g,x,this->
eps_);
501 this->
con_->pruneActive(v,x,this->
eps_);
505 this->
con_->pruneInactive(v,g,x,this->
eps_);
509 this->
con_->pruneInactive(v,x,this->
eps_);
513 return this->
con_->isFeasible(v);
517 return this->
con_->isActivated();
521 this->
con_->computeProjectedStep(v,x);
Provides the interface to evaluate objective functions.
Real value(const Vector< Real > &x, Real &tol)
Compute value.
bool isConActivated(void)
virtual int dimension() const
Return dimension of the vector space.
virtual void plus(const Vector &x)=0
Compute , where .
void reducedHessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &d, const Vector< Real > &x, Real &tol)
Apply the reduced Hessian to a vector, v. The reduced Hessian first removes elements of v correspondi...
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Real dirDeriv(const Vector< Real > &x, const Vector< Real > &d, Real &tol)
Compute directional derivative.
void pruneInactive(Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x)
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
Teuchos::RCP< BoundConstraint< Real > > con_
bool isFeasible(const Vector< Real > &v)
void pruneActive(Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x)
void reducedPrecond(Vector< Real > &Mv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &x, Real &tol)
Apply the reduced preconditioner to a vector, v. The reduced preconditioner first removes elements of...
Defines the linear algebra or vector space interface.
void hessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Teuchos::SerialDenseMatrix< int, Real > computeDenseHessian(Objective< Real > &obj, const Vector< Real > &x)
std::vector< std::vector< Real > > computeGenEigenvalues(const Teuchos::SerialDenseMatrix< int, Real > &A, const Teuchos::SerialDenseMatrix< int, Real > &B)
void reducedInvHessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &x, Real &tol)
Apply the reduced inverse Hessian to a vector, v. The reduced inverse Hessian first removes elements ...
ProjectedObjective(Objective< Real > &obj, BoundConstraint< Real > &con, Teuchos::RCP< Secant< Real > > &secant, bool useSecantPrecond=false, bool useSecantHessVec=false, Real eps=0.0)
Provides interface for and implements limited-memory secant operators.
Teuchos::RCP< Objective< Real > > obj_
std::vector< std::vector< Real > > computeEigenvalues(const Teuchos::SerialDenseMatrix< int, Real > &mat)
void reducedInvHessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &d, const Vector< Real > &x, Real &tol)
Apply the reduced inverse Hessian to a vector, v. The reduced inverse Hessian first removes elements ...
void pruneInactive(Vector< Real > &v, const Vector< Real > &x)
void reducedHessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &x, Real &tol)
Apply the reduced Hessian to a vector, v. The reduced Hessian first removes elements of v correspondi...
Provides the interface to apply upper and lower bound constraints.
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
Teuchos::SerialDenseMatrix< int, Real > computeDotMatrix(const Vector< Real > &x)
virtual Teuchos::RCP< Vector > basis(const int i) const
Return i-th basis vector.
void computeProjectedStep(Vector< Real > &v, const Vector< Real > &x)
void reducedPrecond(Vector< Real > &Mv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &d, const Vector< Real > &x, Real &tol)
Apply the reduced preconditioner to a vector, v. The reduced preconditioner first removes elements of...
void project(Vector< Real > &x)
void invHessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply inverse Hessian approximation to vector.
Teuchos::SerialDenseMatrix< int, Real > computeInverse(const Teuchos::SerialDenseMatrix< int, Real > &mat)
void precond(Vector< Real > &Mv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply preconditioner to vector.
void pruneActive(Vector< Real > &v, const Vector< Real > &x)
Teuchos::RCP< Secant< Real > > secant_
static const double ROL_EPSILON
Platform-dependent machine epsilon.