Actual source code: ntr.c

  1: /*$Id$*/

  3: #include "ntr.h"

  5: #ifdef TAO_USE_PETSC
  6: #include "petscksp.h"
  7: #include "petscpc.h"
  8: #include "src/petsctao/linearsolver/taolinearsolver_petsc.h"
  9: #include "src/petsctao/vector/taovec_petsc.h"

 11: #include "private/kspimpl.h"
 12: #include "private/pcimpl.h"

 14: #define NTR_KSP_STCG    0
 15: #define NTR_KSP_GLTR    1
 16: #define NTR_KSP_TYPES   2

 18: #define NTR_PC_NONE        0
 19: #define NTR_PC_AHESS    1
 20: #define NTR_PC_BFGS     2
 21: #define NTR_PC_PETSC    3
 22: #define NTR_PC_TYPES    4

 24: #define BFGS_SCALE_AHESS   0
 25: #define BFGS_SCALE_BFGS    1
 26: #define BFGS_SCALE_TYPES   2

 28: #define NTR_INIT_CONSTANT          0
 29: #define NTR_INIT_DIRECTION          1
 30: #define NTR_INIT_INTERPOLATION          2
 31: #define NTR_INIT_TYPES                  3

 33: #define NTR_UPDATE_REDUCTION      0
 34: #define NTR_UPDATE_INTERPOLATION  1
 35: #define NTR_UPDATE_TYPES          2

 37: static const char *NTR_KSP[64] = {
 38:   "stcg", "gltr"
 39: };

 41: static const char *NTR_PC[64] = {
 42:   "none", "ahess", "bfgs", "petsc"
 43: };

 45: static const char *BFGS_SCALE[64] = {
 46:   "ahess", "bfgs"
 47: };

 49: static const char *NTR_INIT[64] = {
 50:   "constant", "direction", "interpolation"
 51: };

 53: static const char *NTR_UPDATE[64] = {
 54:   "reduction", "interpolation"
 55: };

 57: // Routine for BFGS preconditioner

 61: static PetscErrorCode bfgs_apply(void *ptr, Vec xin, Vec xout)
 62: {
 63:   TaoLMVMMat *M = (TaoLMVMMat *)ptr;
 64:   TaoVecPetsc Xin(xin);
 65:   TaoVecPetsc Xout(xout);
 66:   TaoTruth info2;
 67:   int info;

 70:   info = M->Solve(&Xin, &Xout, &info2); CHKERRQ(info);
 71:   return(0);
 72: }

 74: // TaoSolve_NTR - Implements Newton's Method with a trust region approach 
 75: // for solving unconstrained minimization problems.  
 76: //
 77: // The basic algorithm is taken from MINPACK-2 (dstrn).
 78: //
 79: // TaoSolve_NTR computes a local minimizer of a twice differentiable function
 80: // f by applying a trust region variant of Newton's method.  At each stage 
 81: // of the algorithm, we use the prconditioned conjugate gradient method to
 82: // determine an approximate minimizer of the quadratic equation
 83: //
 84: //      q(s) = <s, Hs + g>
 85: //
 86: // subject to the trust region constraint
 87: //
 88: //      || s ||_M <= radius,
 89: //
 90: // where radius is the trust region radius and M is a symmetric positive
 91: // definite matrix (the preconditioner).  Here g is the gradient and H 
 92: // is the Hessian matrix.
 93: //
 94: // Note:  TaoSolve_NTR MUST use the iterative solver KSPSTCG or KSPGLTR.
 95: //        Thus, we set KSPSTCG or KSPGLTR in this routine regardless of 
 96: //        what the user may have previously specified.

100: static int TaoSolve_NTR(TAO_SOLVER tao, void *solver)
101: {
102:   TAO_NTR *tr = (TAO_NTR *)solver;
103:   TaoVec *X, *G = tr->G, *D = tr->D, *W = tr->W;
104:   TaoVec *Diag = tr->Diag;
105:   TaoMat *H;
106:   TaoLMVMMat *M = tr->M;

108:   TaoLinearSolver *tls;
109:   TaoLinearSolverPetsc *pls;

111:   KSP pksp;
112:   PC ppc;

114:   TaoTerminateReason reason;
115:   TaoTruth success;

117:   double fmin, ftrial, prered, actred, kappa, sigma, beta;
118:   double tau, tau_1, tau_2, tau_max, tau_min, max_radius;
119:   double f, gnorm;

121:   double delta;
122:   double radius, norm_d;

124:   int iter = 0, info;
125:   int needH;

127:   int i_max = 5;
128:   int j_max = 1;
129:   int i, j;

131:   TaoFunctionBegin;

133:   // Get the initial trust-region radius
134:   info = TaoGetInitialTrustRegionRadius(tao, &radius); CHKERRQ(info);
135:   if (radius < 0.0) {
136:     SETERRQ(1, "Initial radius negative");
137:   }

139:   // Modify the radius if it is too large or small
140:   radius = TaoMax(radius, tr->min_radius);
141:   radius = TaoMin(radius, tr->max_radius);

143:   // Get vectors we will need
144:   info = TaoGetSolution(tao, &X); CHKERRQ(info);
145:   info = TaoGetHessian(tao, &H); CHKERRQ(info);

147:   if (NTR_PC_BFGS == tr->pc_type && !M) {
148:     tr->M = new TaoLMVMMat(X);
149:     M = tr->M;
150:   }

152:   // Check convergence criteria
153:   info = TaoComputeFunctionGradient(tao, X, &f, G); CHKERRQ(info);
154:   info = G->Norm2(&gnorm); CHKERRQ(info);
155:   if (TaoInfOrNaN(f) || TaoInfOrNaN(gnorm)) {
156:     SETERRQ(1, "User provided compute function generated Inf or NaN");
157:   }
158:   needH = 1;

160:   info = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason); CHKERRQ(info);
161:   if (reason != TAO_CONTINUE_ITERATING) {
162:     TaoFunctionReturn(0);
163:   }

165:   // Create vectors for the limited memory preconditioner
166:   if ((NTR_PC_BFGS == tr->pc_type) && 
167:       (BFGS_SCALE_BFGS != tr->bfgs_scale_type)) {
168:     if (!Diag) {
169:       info = X->Clone(&tr->Diag); CHKERRQ(info);
170:       Diag = tr->Diag;
171:     }
172:   }

174:   // Modify the linear solver to a conjugate gradient method
175:   info = TaoGetLinearSolver(tao, &tls); CHKERRQ(info);
176:   pls  = dynamic_cast <TaoLinearSolverPetsc *> (tls);

178:   pksp = pls->GetKSP();
179:   switch(tr->ksp_type) {
180:   case NTR_KSP_STCG:
181:     info = KSPSetType(pksp, KSPSTCG); CHKERRQ(info);
182:     if (pksp->ops->setfromoptions) {
183:       (*pksp->ops->setfromoptions)(pksp);
184:     }
185:     break;

187:   case NTR_KSP_GLTR:
188:     info = KSPSetType(pksp, KSPGLTR); CHKERRQ(info);
189:     if (pksp->ops->setfromoptions) {
190:       (*pksp->ops->setfromoptions)(pksp);
191:     }
192:     break;

194:   default:
195:     // Use the method set by the ksp_type
196:     break;
197:   }

199:   // Modify the preconditioner to use the bfgs approximation
200:   info = KSPGetPC(pksp, &ppc); CHKERRQ(info);
201:   switch(tr->pc_type) {
202:   case NTR_PC_NONE:
203:     info = PCSetType(ppc, PCNONE); CHKERRQ(info);
204:     if (ppc->ops->setfromoptions) {
205:       (*ppc->ops->setfromoptions)(ppc);
206:     }
207:     break;
208:  
209:   case NTR_PC_AHESS:
210:     info = PCSetType(ppc, PCJACOBI); CHKERRQ(info);
211:     if (ppc->ops->setfromoptions) {
212:       (*ppc->ops->setfromoptions)(ppc);
213:     }
214:     info = PCJacobiSetUseAbs(ppc); CHKERRQ(info);
215:     break;

217:   case NTR_PC_BFGS:
218:     info = PCSetType(ppc, PCSHELL); CHKERRQ(info);
219:     if (ppc->ops->setfromoptions) {
220:       (*ppc->ops->setfromoptions)(ppc);
221:     }
222:     info = PCShellSetName(ppc, "bfgs"); CHKERRQ(info);
223:     info = PCShellSetContext(ppc, M); CHKERRQ(info);
224:     info = PCShellSetApply(ppc, bfgs_apply); CHKERRQ(info);
225:     break;

227:   default:
228:     // Use the pc method set by pc_type
229:     break;
230:   }

232:   // Initialize trust-region radius
233:   switch(tr->init_type) {
234:   case NTR_INIT_CONSTANT:
235:     // Use the initial radius specified
236:     break;

238:   case NTR_INIT_INTERPOLATION:
239:     // Use the initial radius specified
240:     max_radius = 0.0;

242:     for (j = 0; j < j_max; ++j) {
243:       fmin = f;
244:       sigma = 0.0;

246:       if (needH) {
247:         info = TaoComputeHessian(tao, X, H); CHKERRQ(info);
248:         needH = 0;
249:       }

251:       for (i = 0; i < i_max; ++i) {
252:         info = W->Waxpby(1.0, X, -radius / gnorm, G); CHKERRQ(info);

254:         info = TaoComputeFunction(tao, W, &ftrial); CHKERRQ(info);
255:         if (TaoInfOrNaN(ftrial)) {
256:           tau = tr->gamma1_i;
257:         }
258:         else {
259:           if (ftrial < fmin) {
260:             fmin = ftrial;
261:             sigma = -radius / gnorm;
262:           }

264:           info = H->Multiply(G, D); CHKERRQ(info);
265:           info = D->Dot(G, &prered); CHKERRQ(info);

267:           prered = radius * (gnorm - 0.5 * radius * prered / (gnorm * gnorm));
268:           actred = f - ftrial;
269:           if ((fabs(actred) <= tr->epsilon) && 
270:               (fabs(prered) <= tr->epsilon)) {
271:             kappa = 1.0;
272:           }
273:           else {
274:             kappa = actred / prered;
275:           }

277:           tau_1 = tr->theta_i * gnorm * radius / (tr->theta_i * gnorm * radius + (1.0 - tr->theta_i) * prered - actred);
278:           tau_2 = tr->theta_i * gnorm * radius / (tr->theta_i * gnorm * radius - (1.0 + tr->theta_i) * prered + actred);
279:           tau_min = TaoMin(tau_1, tau_2);
280:           tau_max = TaoMax(tau_1, tau_2);

282:           if (fabs(kappa - 1.0) <= tr->mu1_i) {
283:             // Great agreement
284:             max_radius = TaoMax(max_radius, radius);

286:             if (tau_max < 1.0) {
287:               tau = tr->gamma3_i;
288:             }
289:             else if (tau_max > tr->gamma4_i) {
290:               tau = tr->gamma4_i;
291:             }
292:             else {
293:               tau = tau_max;
294:             }
295:           }
296:           else if (fabs(kappa - 1.0) <= tr->mu2_i) {
297:             // Good agreement
298:             max_radius = TaoMax(max_radius, radius);

300:             if (tau_max < tr->gamma2_i) {
301:               tau = tr->gamma2_i;
302:             }
303:             else if (tau_max > tr->gamma3_i) {
304:               tau = tr->gamma3_i;
305:             }
306:             else {
307:               tau = tau_max;
308:             }
309:           }
310:           else {
311:             // Not good agreement
312:             if (tau_min > 1.0) {
313:               tau = tr->gamma2_i;
314:             }
315:             else if (tau_max < tr->gamma1_i) {
316:               tau = tr->gamma1_i;
317:             }
318:             else if ((tau_min < tr->gamma1_i) && (tau_max >= 1.0)) {
319:               tau = tr->gamma1_i;
320:             }
321:             else if ((tau_1 >= tr->gamma1_i) && (tau_1 < 1.0) && 
322:                      ((tau_2 < tr->gamma1_i) || (tau_2 >= 1.0))) {
323:               tau = tau_1;
324:             }
325:             else if ((tau_2 >= tr->gamma1_i) && (tau_2 < 1.0) && 
326:                      ((tau_1 < tr->gamma1_i) || (tau_2 >= 1.0))) {
327:               tau = tau_2;
328:             }
329:             else {
330:               tau = tau_max;
331:             }
332:           }
333:         }
334:         radius = tau * radius;
335:       }
336:   
337:       if (fmin < f) {
338:         f = fmin;
339:         info = X->Waxpby(1.0, X, sigma, G); CHKERRQ(info);
340:         info = TaoComputeGradient(tao, X, G); CHKERRQ(info);

342:         info = G->Norm2(&gnorm); CHKERRQ(info);
343:         if (TaoInfOrNaN(f) || TaoInfOrNaN(gnorm)) {
344:           SETERRQ(1, "User provided compute function generated Inf or NaN");
345:         }
346:         needH = 1;

348:         info = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason); CHKERRQ(info);
349:         if (reason != TAO_CONTINUE_ITERATING) {
350:           TaoFunctionReturn(0);
351:         }
352:       }
353:     }
354:     radius = TaoMax(radius, max_radius);

356:     // Modify the radius if it is too large or small
357:     radius = TaoMax(radius, tr->min_radius);
358:     radius = TaoMin(radius, tr->max_radius);
359:     break;

361:   default:
362:     // Norm of the first direction will initialize radius
363:     radius = 0.0;
364:     break;
365:   }

367:   // Set initial scaling for the BFGS preconditioner 
368:   // This step is done after computing the initial trust-region radius
369:   // since the function value may have decreased
370:   if (NTR_PC_BFGS == tr->pc_type) {
371:     if (f != 0.0) {
372:       delta = 2.0 * TaoAbsDouble(f) / (gnorm*gnorm);
373:     }
374:     else {
375:       delta = 2.0 / (gnorm*gnorm);
376:     }
377:     info = M->SetDelta(delta); CHKERRQ(info);
378:   }

380:   // Have not converged; continue with Newton method
381:   while (reason == TAO_CONTINUE_ITERATING) {
382:     ++iter;

384:     // Compute the Hessian
385:     if (needH) {
386:       info = TaoComputeHessian(tao, X, H); CHKERRQ(info);
387:       needH = 0;
388:     }

390:     if (NTR_PC_BFGS == tr->pc_type) {
391:       if (BFGS_SCALE_AHESS == tr->bfgs_scale_type) {
392:         // Obtain diagonal for the bfgs preconditioner
393:         info = H->GetDiagonal(Diag); CHKERRQ(info);
394:         info = Diag->AbsoluteValue(); CHKERRQ(info);
395:         info = Diag->Reciprocal(); CHKERRQ(info);
396:         info = M->SetScale(Diag); CHKERRQ(info);
397:       }

399:       // Update the limited memory preconditioner
400:       info = M->Update(X, G); CHKERRQ(info);
401:     }

403:     while (reason == TAO_CONTINUE_ITERATING) {
404:       // Solve the trust region subproblem
405:       info = TaoPreLinearSolve(tao, H); CHKERRQ(info);
406:       info = TaoLinearSolveTrustRegion(tao, H, G, D, radius, &success); CHKERRQ(info);
407:       info = pls->GetNormDirection(&norm_d); CHKERRQ(info);
408:       if (0.0 == radius) {
409:         // Radius was uninitialized; use the norm of the direction
410:         if (norm_d > 0.0) {
411:           radius = norm_d;

413:           // Modify the radius if it is too large or small
414:           radius = TaoMax(radius, tr->min_radius);
415:           radius = TaoMin(radius, tr->max_radius);
416:         }
417:         else {
418:           // The direction was bad; set radius to default value and re-solve 
419:           // the trust-region subproblem to get a direction
420:           info = TaoGetInitialTrustRegionRadius(tao, &radius); CHKERRQ(info);

422:           // Modify the radius if it is too large or small
423:           radius = TaoMax(radius, tr->min_radius);
424:           radius = TaoMin(radius, tr->max_radius);

426:           info = TaoLinearSolveTrustRegion(tao, H, G, D, radius, &success); CHKERRQ(info);
427:           info = pls->GetNormDirection(&norm_d); CHKERRQ(info);
428:           if (norm_d == 0.0) {
429:             SETERRQ(1, "Initial direction zero");
430:           }
431:         }
432:       }
433:       info = D->Negate(); CHKERRQ(info);

435:       if (NTR_UPDATE_REDUCTION == tr->update_type) {
436:         // Get predicted reduction
437:         info = pls->GetObjFcn(&prered); CHKERRQ(info);
438:         
439:         if (prered >= 0.0) {
440:           // The predicted reduction has the wrong sign.  This cannot
441:           // happen in infinite precision arithmetic.  Step should
442:           // be rejected!
443:           radius = tr->alpha1 * TaoMin(radius, norm_d);
444:         }
445:         else {
446:           // Compute trial step and function value
447:           info = W->Waxpby(1.0, X, 1.0, D); CHKERRQ(info);
448:           info = TaoComputeFunction(tao, W, &ftrial); CHKERRQ(info);
449:           if (TaoInfOrNaN(ftrial)) {
450:             radius = tr->alpha1 * TaoMin(radius, norm_d);
451:           }
452:           else {
453:             // Compute and actual reduction
454:             actred = f - ftrial;
455:             prered = -prered;
456:             if ((fabs(actred) <= tr->epsilon) && 
457:                 (fabs(prered) <= tr->epsilon)) {
458:               kappa = 1.0;
459:             }
460:             else {
461:               kappa = actred / prered;
462:             }
463:             
464:             // Accept of reject the step and update radius
465:             if (kappa < tr->eta1) {
466:               // Reject the step
467:               radius = tr->alpha1 * TaoMin(radius, norm_d);
468:             } 
469:             else {
470:               // Accept the step
471:               if (kappa < tr->eta2) { 
472:                 // Marginal bad step
473:                 radius = tr->alpha2 * TaoMin(radius, norm_d);
474:               }
475:               else if (kappa < tr->eta3) {
476:                 // Reasonable step
477:                 radius = tr->alpha3 * radius;
478:               }
479:               else if (kappa < tr->eta4) { 
480:                 // Good step
481:                 radius = TaoMax(tr->alpha4 * norm_d, radius);
482:               }
483:               else {
484:                 // Very good step
485:                 radius = TaoMax(tr->alpha5 * norm_d, radius);
486:               }
487:               break;
488:             }
489:           }
490:         } 
491:       }
492:       else {
493:         // Get predicted reduction
494:         info = pls->GetObjFcn(&prered); CHKERRQ(info);

496:         if (prered >= 0.0) {
497:           // The predicted reduction has the wrong sign.  This cannot
498:           // happen in infinite precision arithmetic.  Step should
499:           // be rejected!
500:           radius = tr->gamma1 * TaoMin(radius, norm_d);
501:         }
502:         else {
503:           info = W->Waxpby(1.0, X, 1.0, D); CHKERRQ(info);
504:           info = TaoComputeFunction(tao, W, &ftrial); CHKERRQ(info);
505:           if (TaoInfOrNaN(ftrial)) {
506:             radius = tr->gamma1 * TaoMin(radius, norm_d);
507:           }
508:           else {
509:             info = D->Dot(G, &beta); CHKERRQ(info);

511:             actred = f - ftrial;
512:             prered = -prered;
513:             if ((fabs(actred) <= tr->epsilon) && 
514:                 (fabs(prered) <= tr->epsilon)) {
515:               kappa = 1.0;
516:             }
517:             else {
518:               kappa = actred / prered;
519:             }

521:             tau_1 = tr->theta * beta / (tr->theta * beta - (1.0 - tr->theta) * prered + actred);
522:             tau_2 = tr->theta * beta / (tr->theta * beta + (1.0 + tr->theta) * prered - actred);
523:             tau_min = TaoMin(tau_1, tau_2);
524:             tau_max = TaoMax(tau_1, tau_2);

526:             if (kappa >= 1.0 - tr->mu1) {
527:               // Great agreement; accept step and update radius
528:               if (tau_max < 1.0) {
529:                 radius = TaoMax(radius, tr->gamma3 * norm_d);
530:               }
531:               else if (tau_max > tr->gamma4) {
532:                 radius = TaoMax(radius, tr->gamma4 * norm_d);
533:               }
534:               else {
535:                 radius = TaoMax(radius, tau_max * norm_d);
536:               }
537:               break;
538:             }
539:             else if (kappa >= 1.0 - tr->mu2) {
540:               // Good agreement

542:               if (tau_max < tr->gamma2) {
543:                 radius = tr->gamma2 * TaoMin(radius, norm_d);
544:               }
545:               else if (tau_max > tr->gamma3) {
546:                 radius = TaoMax(radius, tr->gamma3 * norm_d);
547:               }
548:               else if (tau_max < 1.0) {
549:                 radius = tau_max * TaoMin(radius, norm_d);
550:               }
551:               else {
552:                 radius = TaoMax(radius, tau_max * norm_d);
553:               }
554:               break;
555:             }
556:             else {
557:               // Not good agreement
558:               if (tau_min > 1.0) {
559:                 radius = tr->gamma2 * TaoMin(radius, norm_d);
560:               }
561:               else if (tau_max < tr->gamma1) {
562:                 radius = tr->gamma1 * TaoMin(radius, norm_d);
563:               }
564:               else if ((tau_min < tr->gamma1) && (tau_max >= 1.0)) {
565:                 radius = tr->gamma1 * TaoMin(radius, norm_d);
566:               }
567:               else if ((tau_1 >= tr->gamma1) && (tau_1 < 1.0) && 
568:                        ((tau_2 < tr->gamma1) || (tau_2 >= 1.0))) {
569:                 radius = tau_1 * TaoMin(radius, norm_d);
570:               }
571:               else if ((tau_2 >= tr->gamma1) && (tau_2 < 1.0) && 
572:                        ((tau_1 < tr->gamma1) || (tau_2 >= 1.0))) {
573:                 radius = tau_2 * TaoMin(radius, norm_d);
574:               }
575:               else {
576:                 radius = tau_max * TaoMin(radius, norm_d);
577:               }
578:             }
579:           }
580:         }
581:       }

583:       // The step computed was not good and the radius was decreased.
584:       // Monitor the radius to terminate.
585:       info = TaoMonitor(tao, iter, f, gnorm, 0.0, radius, &reason); CHKERRQ(info);
586:     }

588:     // The radius may have been increased; modify if it is too large
589:     radius = TaoMin(radius, tr->max_radius);

591:     if (reason == TAO_CONTINUE_ITERATING) {
592:       info = X->CopyFrom(W); CHKERRQ(info);
593:       f = ftrial;
594:       info = TaoComputeGradient(tao, X, G); CHKERRQ(info);
595:       info = G->Norm2(&gnorm); CHKERRQ(info);
596:       if (TaoInfOrNaN(f) || TaoInfOrNaN(gnorm)) {
597:         SETERRQ(1, "User provided compute function generated Inf or NaN");
598:       }
599:       needH = 1;

601:       info = TaoMonitor(tao, iter, f, gnorm, 0.0, radius, &reason); CHKERRQ(info);
602:     }
603:   }
604:   TaoFunctionReturn(0);
605: }

607: /*------------------------------------------------------------*/
610: static int TaoSetUp_NTR(TAO_SOLVER tao, void *solver)
611: {
612:   TAO_NTR *tr = (TAO_NTR *)solver;
613:   TaoVec *X;
614:   TaoMat *H;
615:   int info;

617:   TaoFunctionBegin;

619:   info = TaoGetSolution(tao, &X); CHKERRQ(info);
620:   info = X->Clone(&tr->G); CHKERRQ(info);
621:   info = X->Clone(&tr->D); CHKERRQ(info);
622:   info = X->Clone(&tr->W); CHKERRQ(info);

624:   tr->Diag = 0;
625:   tr->M = 0;

627:   info = TaoSetLagrangianGradientVector(tao, tr->G); CHKERRQ(info);
628:   info = TaoSetStepDirectionVector(tao, tr->D); CHKERRQ(info);

630:   // Set linear solver to default for trust region
631:   info = TaoGetHessian(tao, &H); CHKERRQ(info);
632:   info = TaoCreateLinearSolver(tao, H, 220, 0); CHKERRQ(info); 

634:   // Check sizes for compatability
635:   info = TaoCheckFGH(tao); CHKERRQ(info);
636:   TaoFunctionReturn(0);
637: }

639: /*------------------------------------------------------------*/
642: static int TaoDestroy_NTR(TAO_SOLVER tao, void *solver)
643: {
644:   TAO_NTR *tr = (TAO_NTR *)solver;
645:   int info;

647:   TaoFunctionBegin;
648:   info = TaoVecDestroy(tr->G); CHKERRQ(info);
649:   info = TaoVecDestroy(tr->D); CHKERRQ(info);
650:   info = TaoVecDestroy(tr->W); CHKERRQ(info);

652:   info = TaoSetLagrangianGradientVector(tao, 0); CHKERRQ(info);
653:   info = TaoSetStepDirectionVector(tao, 0); CHKERRQ(info);

655:   if (tr->Diag) {
656:     info = TaoVecDestroy(tr->Diag); CHKERRQ(info);
657:     tr->Diag = 0;
658:   }

660:   if (tr->M) {
661:     info = TaoMatDestroy(tr->M); CHKERRQ(info);
662:     tr->M = 0;
663:   }
664:   TaoFunctionReturn(0);
665: }

667: /*------------------------------------------------------------*/
670: static int TaoSetOptions_NTR(TAO_SOLVER tao, void*solver)
671: {
672:   TAO_NTR *tr = (TAO_NTR *)solver;
673:   int info;

675:   TaoFunctionBegin;
676:   info = TaoOptionsHead("Newton trust region method for unconstrained optimization"); CHKERRQ(info);
677:   info = TaoOptionList("-tao_ntr_ksp_type", "ksp type", "", NTR_KSP, NTR_KSP_TYPES, NTR_KSP[tr->ksp_type], &tr->ksp_type, 0); CHKERRQ(info);
678:   info = TaoOptionList("-tao_ntr_pc_type", "pc type", "", NTR_PC, NTR_PC_TYPES, NTR_PC[tr->pc_type], &tr->pc_type, 0); CHKERRQ(info);
679:   info = TaoOptionList("-tao_ntr_bfgs_scale_type", "bfgs scale type", "", BFGS_SCALE, BFGS_SCALE_TYPES, BFGS_SCALE[tr->bfgs_scale_type], &tr->bfgs_scale_type, 0); CHKERRQ(info);
680:   info = TaoOptionList("-tao_ntr_init_type", "radius initialization type", "", NTR_INIT, NTR_INIT_TYPES, NTR_INIT[tr->init_type], &tr->init_type, 0); CHKERRQ(info);
681:   info = TaoOptionList("-tao_ntr_update_type", "radius update type", "", NTR_UPDATE, NTR_UPDATE_TYPES, NTR_UPDATE[tr->update_type], &tr->update_type, 0); CHKERRQ(info);
682:   info = TaoOptionDouble("-tao_ntr_eta1", "step is unsuccessful if actual reduction < eta1 * predicted reduction", "", tr->eta1, &tr->eta1, 0); CHKERRQ(info);
683:   info = TaoOptionDouble("-tao_ntr_eta2", "", "", tr->eta2, &tr->eta2, 0); CHKERRQ(info);
684:   info = TaoOptionDouble("-tao_ntr_eta3", "", "", tr->eta3, &tr->eta3, 0); CHKERRQ(info);
685:   info = TaoOptionDouble("-tao_ntr_eta4", "", "", tr->eta4, &tr->eta4, 0); CHKERRQ(info);
686:   info = TaoOptionDouble("-tao_ntr_alpha1", "", "", tr->alpha1, &tr->alpha1, 0); CHKERRQ(info);
687:   info = TaoOptionDouble("-tao_ntr_alpha2", "", "", tr->alpha2, &tr->alpha2, 0); CHKERRQ(info);
688:   info = TaoOptionDouble("-tao_ntr_alpha3", "", "", tr->alpha3, &tr->alpha3, 0); CHKERRQ(info);
689:   info = TaoOptionDouble("-tao_ntr_alpha4", "", "", tr->alpha4, &tr->alpha4, 0); CHKERRQ(info);
690:   info = TaoOptionDouble("-tao_ntr_alpha5", "", "", tr->alpha5, &tr->alpha5, 0); CHKERRQ(info);
691:   info = TaoOptionDouble("-tao_ntr_mu1", "", "", tr->mu1, &tr->mu1, 0); CHKERRQ(info);
692:   info = TaoOptionDouble("-tao_ntr_mu2", "", "", tr->mu2, &tr->mu2, 0); CHKERRQ(info);
693:   info = TaoOptionDouble("-tao_ntr_gamma1", "", "", tr->gamma1, &tr->gamma1, 0); CHKERRQ(info);
694:   info = TaoOptionDouble("-tao_ntr_gamma2", "", "", tr->gamma2, &tr->gamma2, 0); CHKERRQ(info);
695:   info = TaoOptionDouble("-tao_ntr_gamma3", "", "", tr->gamma3, &tr->gamma3, 0); CHKERRQ(info);
696:   info = TaoOptionDouble("-tao_ntr_gamma4", "", "", tr->gamma4, &tr->gamma4, 0); CHKERRQ(info);
697:   info = TaoOptionDouble("-tao_ntr_theta", "", "", tr->theta, &tr->theta, 0); CHKERRQ(info);
698:   info = TaoOptionDouble("-tao_ntr_mu1_i", "", "", tr->mu1_i, &tr->mu1_i, 0); CHKERRQ(info);
699:   info = TaoOptionDouble("-tao_ntr_mu2_i", "", "", tr->mu2_i, &tr->mu2_i, 0); CHKERRQ(info);
700:   info = TaoOptionDouble("-tao_ntr_gamma1_i", "", "", tr->gamma1_i, &tr->gamma1_i, 0); CHKERRQ(info);
701:   info = TaoOptionDouble("-tao_ntr_gamma2_i", "", "", tr->gamma2_i, &tr->gamma2_i, 0); CHKERRQ(info);
702:   info = TaoOptionDouble("-tao_ntr_gamma3_i", "", "", tr->gamma3_i, &tr->gamma3_i, 0); CHKERRQ(info);
703:   info = TaoOptionDouble("-tao_ntr_gamma4_i", "", "", tr->gamma4_i, &tr->gamma4_i, 0); CHKERRQ(info);
704:   info = TaoOptionDouble("-tao_ntr_theta_i", "", "", tr->theta_i, &tr->theta_i, 0); CHKERRQ(info);
705:   info = TaoOptionDouble("-tao_ntr_min_radius", "lower bound on initial trust-region radius", "", tr->min_radius, &tr->min_radius, 0); CHKERRQ(info);
706:   info = TaoOptionDouble("-tao_ntr_max_radius", "upper bound on trust-region radius", "", tr->max_radius, &tr->max_radius, 0); CHKERRQ(info);
707:   info = TaoOptionDouble("-tao_ntr_epsilon", "tolerance used when computing actual and predicted reduction", "", tr->epsilon, &tr->epsilon, 0); CHKERRQ(info);
708:   info = TaoOptionsTail(); CHKERRQ(info);
709:   TaoFunctionReturn(0);
710: }

712: /*------------------------------------------------------------*/
715: static int TaoView_NTR(TAO_SOLVER tao,void*solver)
716: {
717:   TAO_NTR *tr = (TAO_NTR *)solver;
718:   int info;

720:   TaoFunctionBegin;
721:   if (NTR_PC_BFGS == tr->pc_type && tr->M) {
722:     info = TaoPrintInt(tao, "  Rejected matrix updates: %d\n", tr->M->GetRejects()); CHKERRQ(info);
723:   }
724:   TaoFunctionReturn(0);
725: }

727: /*------------------------------------------------------------*/
731: int TaoCreate_NTR(TAO_SOLVER tao)
732: {
733:   TAO_NTR *tr;
734:   int info;

736:   TaoFunctionBegin;

738:   info = TaoNew(TAO_NTR, &tr); CHKERRQ(info);
739:   info = PetscLogObjectMemory(tao, sizeof(TAO_NTR)); CHKERRQ(info);

741:   info = TaoSetTaoSolveRoutine(tao, TaoSolve_NTR, (void *)tr); CHKERRQ(info);
742:   info = TaoSetTaoSetUpDownRoutines(tao, TaoSetUp_NTR, TaoDestroy_NTR); CHKERRQ(info);
743:   info = TaoSetTaoOptionsRoutine(tao, TaoSetOptions_NTR); CHKERRQ(info);
744:   info = TaoSetTaoViewRoutine(tao, TaoView_NTR); CHKERRQ(info);

746:   info = TaoSetMaximumIterates(tao, 50); CHKERRQ(info);
747:   info = TaoSetTolerances(tao, 1e-10, 1e-10, 0, 0); CHKERRQ(info);

749:   info = TaoSetTrustRegionRadius(tao, 100.0); CHKERRQ(info);
750:   info = TaoSetTrustRegionTolerance(tao, 1.0e-12); CHKERRQ(info);

752:   // Standard trust region update parameters
753:   tr->eta1 = 1.0e-4;
754:   tr->eta2 = 0.25;
755:   tr->eta3 = 0.50;
756:   tr->eta4 = 0.90;

758:   tr->alpha1 = 0.25;
759:   tr->alpha2 = 0.50;
760:   tr->alpha3 = 1.00;
761:   tr->alpha4 = 2.00;
762:   tr->alpha5 = 4.00;

764:   // Interpolation parameters
765:   tr->mu1_i = 0.35;
766:   tr->mu2_i = 0.50;

768:   tr->gamma1_i = 0.0625;
769:   tr->gamma2_i = 0.50;
770:   tr->gamma3_i = 2.00;
771:   tr->gamma4_i = 5.00;

773:   tr->theta_i = 0.25;

775:   // Interpolation trust region update parameters
776:   tr->mu1 = 0.10;
777:   tr->mu2 = 0.50;

779:   tr->gamma1 = 0.25;
780:   tr->gamma2 = 0.50;
781:   tr->gamma3 = 2.00;
782:   tr->gamma4 = 4.00;

784:   tr->theta = 0.05;

786:   tr->min_radius = 1.0e-10;
787:   tr->max_radius = 1.0e10;
788:   tr->epsilon = 1.0e-6;

790:   tr->ksp_type        = NTR_KSP_STCG;
791:   tr->pc_type         = NTR_PC_BFGS;
792:   tr->bfgs_scale_type = BFGS_SCALE_AHESS;
793:   tr->init_type              = NTR_INIT_INTERPOLATION;
794:   tr->update_type     = NTR_UPDATE_REDUCTION;
795:   TaoFunctionReturn(0);
796: }

799: #endif