45 #ifndef ROL_EQUALITYCONSTRAINT_DEF_H
46 #define ROL_EQUALITYCONSTRAINT_DEF_H
61 Real h = std::max(1.0,x.
norm()/v.
norm())*tol;
65 Teuchos::RCP<Vector<Real> > c = jv.
clone();
66 this->value(*c,x,ctol);
69 Teuchos::RCP<Vector<Real> > xnew = x.
clone();
76 this->value(jv,*xnew,ctol);
89 applyAdjointJacobian(ajv,v,x,v.
dual(),tol);
110 Teuchos::RCP<Vector<Real> > xnew = x.
clone();
111 Teuchos::RCP<Vector<Real> > ex = x.
clone();
112 Teuchos::RCP<Vector<Real> > eajv = ajv.
clone();
113 Teuchos::RCP<Vector<Real> > cnew = dualv.
clone();
114 Teuchos::RCP<Vector<Real> > c0 = dualv.
clone();
115 this->value(*c0,x,ctol);
117 for (
int i = 0; i < ajv.
dimension(); i++ ) {
120 h = std::max(1.0,x.
norm()/ex->norm())*tol;
124 this->value(*cnew,*xnew,ctol);
125 cnew->axpy(-1.0,*c0);
127 ajv.
axpy(cnew->dot(v.
dual()),*eajv);
164 template <
class Real>
173 template <
class Real>
191 tol = std::sqrt(b1.
dot(b1)+b2.
dot(b2))*1e-8;
197 Teuchos::RCP<Vector<Real> > r1 = b1.
clone();
198 Teuchos::RCP<Vector<Real> > r2 = b2.
clone();
199 Teuchos::RCP<Vector<Real> > z1 = v1.
clone();
200 Teuchos::RCP<Vector<Real> > z2 = v2.
clone();
201 Teuchos::RCP<Vector<Real> > w1 = b1.
clone();
202 Teuchos::RCP<Vector<Real> > w2 = b2.
clone();
203 std::vector<Teuchos::RCP<Vector<Real> > > V1;
204 std::vector<Teuchos::RCP<Vector<Real> > > V2;
205 Teuchos::RCP<Vector<Real> > V2temp = b2.
clone();
206 std::vector<Teuchos::RCP<Vector<Real> > > Z1;
207 std::vector<Teuchos::RCP<Vector<Real> > > Z2;
208 Teuchos::RCP<Vector<Real> > w1temp = b1.
clone();
209 Teuchos::RCP<Vector<Real> > Z2temp = v2.
clone();
211 std::vector<Real> res(m+1, zero);
212 Teuchos::SerialDenseMatrix<int, Real> H(m+1,m);
213 Teuchos::SerialDenseVector<int, Real> cs(m);
214 Teuchos::SerialDenseVector<int, Real> sn(m);
215 Teuchos::SerialDenseVector<int, Real> s(m+1);
216 Teuchos::SerialDenseVector<int, Real> y(m+1);
217 Teuchos::SerialDenseVector<int, Real> cnorm(m);
218 Teuchos::LAPACK<int, Real> lapack;
221 applyAdjointJacobian(*r1, v2, x, zerotol);
222 r1->scale(-one); r1->axpy(-one, v1.
dual()); r1->plus(b1);
223 applyJacobian(*r2, v1, x, zerotol);
224 r2->scale(-one); r2->plus(b2);
225 res[0] = std::sqrt(r1->dot(*r1) + r2->dot(*r2));
228 if (res[0] == zero) {
233 V1.push_back(b1.
clone()); (V1[0])->set(*r1); (V1[0])->scale(one/res[0]);
234 V2.push_back(b2.
clone()); (V2[0])->set(*r2); (V2[0])->scale(one/res[0]);
238 for (i=0; i<m; i++) {
241 V2temp->set(*(V2[i]));
242 applyPreconditioner(*Z2temp, *V2temp, x, b1, zerotol);
243 Z2.push_back(v2.
clone()); (Z2[i])->set(*Z2temp);
244 Z1.push_back(v1.
clone()); (Z1[i])->set((V1[i])->dual());
247 applyJacobian(*w2, *(Z1[i]), x, zerotol);
248 applyAdjointJacobian(*w1temp, *Z2temp, x, zerotol);
249 w1->set(*(V1[i])); w1->plus(*w1temp);
252 for (k=0; k<=i; k++) {
253 H(k,i) = w1->dot(*(V1[k])) + w2->dot(*(V2[k]));
254 w1->axpy(-H(k,i), *(V1[k]));
255 w2->axpy(-H(k,i), *(V2[k]));
257 H(i+1,i) = std::sqrt(w1->dot(*w1) + w2->dot(*w2));
259 V1.push_back(b1.
clone()); (V1[i+1])->set(*w1); (V1[i+1])->scale(one/H(i+1,i));
260 V2.push_back(b2.
clone()); (V2[i+1])->set(*w2); (V2[i+1])->scale(one/H(i+1,i));
263 for (k=0; k<=i-1; k++) {
264 temp = cs(k)*H(k,i) + sn(k)*H(k+1,i);
265 H(k+1,i) = -sn(k)*H(k,i) + cs(k)*H(k+1,i);
270 if ( H(i+1,i) == zero ) {
274 else if ( std::abs(H(i+1,i)) > std::abs(H(i,i)) ) {
275 temp = H(i,i) / H(i+1,i);
276 sn(i) = one / std::sqrt( one + temp*temp );
277 cs(i) = temp * sn(i);
280 temp = H(i+1,i) / H(i,i);
281 cs(i) = one / std::sqrt( one + temp*temp );
282 sn(i) = temp * cs(i);
287 s(i+1) = -sn(i)*s(i);
289 H(i,i) = cs(i)*H(i,i) + sn(i)*H(i+1,i);
291 resnrm = std::abs(s(i+1));
295 const char uplo =
'U';
296 const char trans =
'N';
297 const char diag =
'N';
298 const char normin =
'N';
302 lapack.LATRS(uplo, trans, diag, normin, i+1, H.values(), m+1, y.values(), &scaling, cnorm.values(), &info);
305 for (k=0; k<=i; k++) {
306 z1->axpy(y(k), *(Z1[k]));
307 z2->axpy(y(k), *(Z2[k]));
313 if (res[i+1] <= tol) {
341 template <
class Real>
345 const bool printToStream,
346 std::ostream & outStream,
349 std::vector<Real> steps(numSteps);
350 for(
int i=0;i<numSteps;++i) {
351 steps[i] = pow(10,-i);
354 return checkApplyJacobian(x,v,jv,steps,printToStream,outStream,order);
360 template <
class Real>
364 const std::vector<Real> &steps,
365 const bool printToStream,
366 std::ostream & outStream,
369 TEUCHOS_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
370 "Error: finite difference order must be 1,2,3, or 4" );
377 int numSteps = steps.size();
379 std::vector<Real> tmp(numVals);
380 std::vector<std::vector<Real> > jvCheck(numSteps, tmp);
383 Teuchos::oblackholestream oldFormatState;
384 oldFormatState.copyfmt(outStream);
387 Teuchos::RCP<Vector<Real> > c = jv.
clone();
388 this->value(*c, x, tol);
391 Teuchos::RCP<Vector<Real> > Jv = jv.
clone();
392 this->applyJacobian(*Jv, v, x, tol);
393 Real normJv = Jv->norm();
396 Teuchos::RCP<Vector<Real> > cdif = jv.
clone();
397 Teuchos::RCP<Vector<Real> > cnew = jv.
clone();
398 Teuchos::RCP<Vector<Real> > xnew = x.
clone();
400 for (
int i=0; i<numSteps; i++) {
407 cdif->scale(
weights[order-1][0]);
409 for(
int j=0; j<order; ++j) {
411 xnew->axpy(eta*
shifts[order-1][j], v);
413 if(
weights[order-1][j+1] != 0 ) {
415 this->value(*cnew,*xnew,tol);
416 cdif->axpy(
weights[order-1][j+1],*cnew);
421 cdif->scale(1.0/eta);
425 jvCheck[i][1] = normJv;
426 jvCheck[i][2] = cdif->norm();
427 cdif->axpy(-1.0, *Jv);
428 jvCheck[i][3] = cdif->norm();
431 std::stringstream hist;
434 << std::setw(20) <<
"Step size"
435 << std::setw(20) <<
"norm(Jac*vec)"
436 << std::setw(20) <<
"norm(FD approx)"
437 << std::setw(20) <<
"norm(abs error)"
439 << std::setw(20) <<
"---------"
440 << std::setw(20) <<
"-------------"
441 << std::setw(20) <<
"---------------"
442 << std::setw(20) <<
"---------------"
445 hist << std::scientific << std::setprecision(11) << std::right
446 << std::setw(20) << jvCheck[i][0]
447 << std::setw(20) << jvCheck[i][1]
448 << std::setw(20) << jvCheck[i][2]
449 << std::setw(20) << jvCheck[i][3]
451 outStream << hist.str();
457 outStream.copyfmt(oldFormatState);
463 template <
class Real>
468 const bool printToStream,
469 std::ostream & outStream,
470 const int numSteps) {
474 std::vector<Real> tmp(numVals);
475 std::vector<std::vector<Real> > ajvCheck(numSteps, tmp);
476 Real eta_factor = 1e-1;
480 Teuchos::RCP<Vector<Real> > c0 = c.
clone();
481 Teuchos::RCP<Vector<Real> > cnew = c.
clone();
482 Teuchos::RCP<Vector<Real> > xnew = x.
clone();
483 Teuchos::RCP<Vector<Real> > ajv0 = ajv.
clone();
484 Teuchos::RCP<Vector<Real> > ajv1 = ajv.
clone();
485 Teuchos::RCP<Vector<Real> > ex = x.
clone();
486 Teuchos::RCP<Vector<Real> > eajv = ajv.
clone();
489 Teuchos::oblackholestream oldFormatState;
490 oldFormatState.copyfmt(outStream);
493 this->value(*c0, x, tol);
496 this->applyAdjointJacobian(*ajv0, v, x, tol);
497 Real normAJv = ajv0->norm();
499 for (
int i=0; i<numSteps; i++) {
503 for (
int j = 0; j < ajv.
dimension(); j++ ) {
509 this->value(*cnew,*xnew,tol);
510 cnew->axpy(-1.0,*c0);
511 cnew->scale(1.0/eta);
512 ajv1->axpy(cnew->dot(v.
dual()),*eajv);
516 ajvCheck[i][0] = eta;
517 ajvCheck[i][1] = normAJv;
518 ajvCheck[i][2] = ajv1->norm();
519 ajv1->axpy(-1.0, *ajv0);
520 ajvCheck[i][3] = ajv1->norm();
523 std::stringstream hist;
526 << std::setw(20) <<
"Step size"
527 << std::setw(20) <<
"norm(adj(Jac)*vec)"
528 << std::setw(20) <<
"norm(FD approx)"
529 << std::setw(20) <<
"norm(abs error)"
531 << std::setw(20) <<
"---------"
532 << std::setw(20) <<
"------------------"
533 << std::setw(20) <<
"---------------"
534 << std::setw(20) <<
"---------------"
537 hist << std::scientific << std::setprecision(11) << std::right
538 << std::setw(20) << ajvCheck[i][0]
539 << std::setw(20) << ajvCheck[i][1]
540 << std::setw(20) << ajvCheck[i][2]
541 << std::setw(20) << ajvCheck[i][3]
543 outStream << hist.str();
547 eta = eta*eta_factor;
551 outStream.copyfmt(oldFormatState);
556 template <
class Real>
562 const bool printToStream,
563 std::ostream & outStream) {
566 Teuchos::RCP<Vector<Real> > Jv = dualw.
clone();
567 Teuchos::RCP<Vector<Real> > Jw = dualv.
clone();
569 applyJacobian(*Jv,v,x,tol);
570 applyAdjointJacobian(*Jw,w,x,tol);
572 Real vJw = v.
dot(Jw->dual());
573 Real wJv = w.
dot(Jv->dual());
575 Real diff = std::abs(wJv-vJw);
577 if ( printToStream ) {
578 std::stringstream hist;
579 hist << std::scientific << std::setprecision(8);
580 hist <<
"\nTest Consistency of Jacobian and its adjoint: \n |<w,Jv> - <adj(J)w,v>| = "
582 hist <<
" |<w,Jv>| = " << std::abs(wJv) <<
"\n";
583 hist <<
" Relative Error = " << diff / (std::abs(wJv)+
ROL_UNDERFLOW) <<
"\n";
584 outStream << hist.str();
589 template <
class Real>
594 const bool printToStream,
595 std::ostream & outStream,
598 std::vector<Real> steps(numSteps);
599 for(
int i=0;i<numSteps;++i) {
600 steps[i] = pow(10,-i);
603 return checkApplyAdjointHessian(x,u,v,hv,steps,printToStream,outStream,order);
608 template <
class Real>
613 const std::vector<Real> &steps,
614 const bool printToStream,
615 std::ostream & outStream,
622 int numSteps = steps.size();
624 std::vector<Real> tmp(numVals);
625 std::vector<std::vector<Real> > ahuvCheck(numSteps, tmp);
628 Teuchos::RCP<Vector<Real> > AJdif = hv.
clone();
629 Teuchos::RCP<Vector<Real> > AJu = hv.
clone();
630 Teuchos::RCP<Vector<Real> > AHuv = hv.
clone();
631 Teuchos::RCP<Vector<Real> > AJnew = hv.
clone();
632 Teuchos::RCP<Vector<Real> > xnew = x.
clone();
635 Teuchos::oblackholestream oldFormatState;
636 oldFormatState.copyfmt(outStream);
640 this->applyAdjointJacobian(*AJu, u, x, tol);
643 this->applyAdjointHessian(*AHuv, u, v, x, tol);
644 Real normAHuv = AHuv->norm();
646 for (
int i=0; i<numSteps; i++) {
654 AJdif->scale(
weights[order-1][0]);
656 for(
int j=0; j<order; ++j) {
658 xnew->axpy(eta*
shifts[order-1][j],v);
660 if(
weights[order-1][j+1] != 0 ) {
662 this->applyAdjointJacobian(*AJnew, u, *xnew, tol);
663 AJdif->axpy(
weights[order-1][j+1],*AJnew);
667 AJdif->scale(1.0/eta);
670 ahuvCheck[i][0] = eta;
671 ahuvCheck[i][1] = normAHuv;
672 ahuvCheck[i][2] = AJdif->norm();
673 AJdif->axpy(-1.0, *AHuv);
674 ahuvCheck[i][3] = AJdif->norm();
677 std::stringstream hist;
680 << std::setw(20) <<
"Step size"
681 << std::setw(20) <<
"norm(adj(H)(u,v))"
682 << std::setw(20) <<
"norm(FD approx)"
683 << std::setw(20) <<
"norm(abs error)"
685 << std::setw(20) <<
"---------"
686 << std::setw(20) <<
"-----------------"
687 << std::setw(20) <<
"---------------"
688 << std::setw(20) <<
"---------------"
691 hist << std::scientific << std::setprecision(11) << std::right
692 << std::setw(20) << ahuvCheck[i][0]
693 << std::setw(20) << ahuvCheck[i][1]
694 << std::setw(20) << ahuvCheck[i][2]
695 << std::setw(20) << ahuvCheck[i][3]
697 outStream << hist.str();
703 outStream.copyfmt(oldFormatState);
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
virtual void scale(const Real alpha)=0
Compute where .
virtual int dimension() const
Return dimension of the vector space.
virtual void plus(const Vector &x)=0
Compute , where .
const double weights[4][5]
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
virtual Real checkAdjointConsistencyJacobian(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &x, const bool printToStream=true, std::ostream &outStream=std::cout)
virtual void applyAdjointHessian(Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the derivative of the adjoint of the constraint Jacobian at to vector in direction ...
virtual std::vector< std::vector< Real > > checkApplyAdjointJacobian(const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &c, const Vector< Real > &ajv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS)
Finite-difference check for the application of the adjoint of constraint Jacobian.
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
virtual Real dot(const Vector &x) const =0
Compute where .
virtual std::vector< std::vector< Real > > checkApplyJacobian(const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &jv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
Finite-difference check for the constraint Jacobian application.
virtual std::vector< Real > solveAugmentedSystem(Vector< Real > &v1, Vector< Real > &v2, const Vector< Real > &b1, const Vector< Real > &b2, const Vector< Real > &x, Real &tol)
Approximately solves the augmented system where , , , , is an identity operator, and is a zero operator.
virtual void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the adjoint of the the constraint Jacobian at , , to vector .
virtual void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector .
virtual Teuchos::RCP< Vector > basis(const int i) const
Return i-th basis vector.
virtual Real norm() const =0
Returns where .
virtual std::vector< std::vector< Real > > checkApplyAdjointHessian(const Vector< Real > &x, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &step, const bool printToScreen=true, std::ostream &outStream=std::cout, const int order=1)
Finite-difference check for the application of the adjoint of constraint Hessian. ...
static const double ROL_UNDERFLOW
Platform-dependent minimum double.
static const double ROL_EPSILON
Platform-dependent machine epsilon.