xref: /petsc/src/tao/bound/impls/bnk/bnls.c (revision 62675bee48cd89cd8757a6de6d7b87f2edd3afcc)
1eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
2eb910715SAlp Dener #include <petscksp.h>
3eb910715SAlp Dener 
4eb910715SAlp Dener /*
5eb910715SAlp Dener  Implements Newton's Method with a line search approach for solving
62b97c8d8SAlp Dener  bound constrained minimization problems. A projected More'-Thuente line
7770b7498SAlp Dener  search is used to guarantee that the BFGS preconditioner remains positive
8eb910715SAlp Dener  definite.
9eb910715SAlp Dener 
10eb910715SAlp Dener  The method can shift the Hessian matrix. The shifting procedure is
11eb910715SAlp Dener  adapted from the PATH algorithm for solving complementarity
12eb910715SAlp Dener  problems.
13eb910715SAlp Dener 
14eb910715SAlp Dener  The linear system solve should be done with a conjugate gradient
15eb910715SAlp Dener  method, although any method can be used.
16eb910715SAlp Dener */
17eb910715SAlp Dener 
18eb910715SAlp Dener static PetscErrorCode TaoSolve_BNLS(Tao tao)
19eb910715SAlp Dener {
20eb910715SAlp Dener   PetscErrorCode               ierr;
21eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
22e465cd6fSAlp Dener   KSPConvergedReason           ksp_reason;
23eb910715SAlp Dener   TaoLineSearchConvergedReason ls_reason;
24eb910715SAlp Dener 
25c14b763aSAlp Dener   PetscReal                    steplen = 1.0;
26*62675beeSAlp Dener   PetscBool                    shift = PETSC_TRUE;
27eb910715SAlp Dener   PetscInt                     stepType;
28eb910715SAlp Dener 
29eb910715SAlp Dener   PetscFunctionBegin;
3028017e9fSAlp Dener   /* Initialize the preconditioner, KSP solver and trust radius/line search */
31eb910715SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
32*62675beeSAlp Dener   ierr = TaoBNKInitialize(tao, BNK_INIT_CONSTANT);CHKERRQ(ierr);
3328017e9fSAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
34eb910715SAlp Dener 
35eb910715SAlp Dener   /* Have not converged; continue with Newton method */
36eb910715SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
37eb910715SAlp Dener     ++tao->niter;
38eb910715SAlp Dener     tao->ksp_its=0;
39eb910715SAlp Dener 
40*62675beeSAlp Dener     /* Compute the hessian and update the BFGS preconditioner at the new iterate*/
41*62675beeSAlp Dener     ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
42fed79b8eSAlp Dener 
438d5ead36SAlp Dener     /* Use the common BNK kernel to compute the safeguarded Newton step (for inactive variables only) */
4428017e9fSAlp Dener     tao->trust = bnk->max_radius;
45*62675beeSAlp Dener     ierr = TaoBNKComputeStep(tao, shift, &ksp_reason);CHKERRQ(ierr);
46e465cd6fSAlp Dener     ierr = TaoBNKSafeguardStep(tao, ksp_reason, &stepType);CHKERRQ(ierr);
47eb910715SAlp Dener 
48080d2917SAlp Dener     /* Store current solution before it changes */
49080d2917SAlp Dener     bnk->fold = bnk->f;
50eb910715SAlp Dener     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
51eb910715SAlp Dener     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
5209164190SAlp Dener     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
53eb910715SAlp Dener 
54c14b763aSAlp Dener     /* Trigger the line search */
55c14b763aSAlp Dener     ierr = TaoBNKPerformLineSearch(tao, stepType, &steplen, &ls_reason);CHKERRQ(ierr);
56eb910715SAlp Dener 
57eb910715SAlp Dener     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
58eb910715SAlp Dener       /* Failed to find an improving point */
59080d2917SAlp Dener       bnk->f = bnk->fold;
60eb910715SAlp Dener       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
61eb910715SAlp Dener       ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
6209164190SAlp Dener       ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
63c14b763aSAlp Dener       steplen = 0.0;
64eb910715SAlp Dener       tao->reason = TAO_DIVERGED_LS_FAILURE;
65eb910715SAlp Dener       break;
66e465cd6fSAlp Dener     } else {
67*62675beeSAlp Dener       /* count the accepted step type */
68*62675beeSAlp Dener       ierr = TaoBNKAddStepCounts(tao, stepType);CHKERRQ(ierr);
69eb910715SAlp Dener     }
70eb910715SAlp Dener 
71eb910715SAlp Dener     /*  Check for termination */
72eb910715SAlp Dener     ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
73eb910715SAlp Dener     if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
74eb910715SAlp Dener     ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
75c14b763aSAlp Dener     ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,steplen);CHKERRQ(ierr);
76eb910715SAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
77eb910715SAlp Dener   }
78eb910715SAlp Dener   PetscFunctionReturn(0);
79eb910715SAlp Dener }
80eb910715SAlp Dener 
81df278d8fSAlp Dener /*------------------------------------------------------------*/
82df278d8fSAlp Dener 
83eb910715SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNLS(Tao tao)
84eb910715SAlp Dener {
85fed79b8eSAlp Dener   TAO_BNK        *bnk;
86eb910715SAlp Dener   PetscErrorCode ierr;
87eb910715SAlp Dener 
88eb910715SAlp Dener   PetscFunctionBegin;
89eb910715SAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
90eb910715SAlp Dener   tao->ops->solve = TaoSolve_BNLS;
91fed79b8eSAlp Dener 
92fed79b8eSAlp Dener   bnk = (TAO_BNK *)tao->data;
9328017e9fSAlp Dener   bnk->init_type = BNK_INIT_CONSTANT;
9466ed3702SAlp Dener   bnk->update_type = BNK_UPDATE_STEP; /* trust region updates based on line search step length */
95eb910715SAlp Dener   PetscFunctionReturn(0);
96eb910715SAlp Dener }