xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision cd929ea3f739fd9f7b6394f772cb40b9d4e6d97c)
1ba92ff59SBarry Smith #include <petsctaolinesearch.h>
2a9603a14SPatrick Farrell #include <../src/tao/unconstrained/impls/lmvm/lmvm.h>
3aaa7dc30SBarry Smith #include <../src/tao/bound/impls/blmvm/blmvm.h>
4a7e14dcfSSatish Balay 
5*cd929ea3SAlp Dener #define BLMVM_STEP_BFGS     0
6*cd929ea3SAlp Dener #define BLMVM_STEP_GRAD     1
7*cd929ea3SAlp Dener 
8*cd929ea3SAlp Dener #define BLMVM_AS_NONE       0
9*cd929ea3SAlp Dener #define BLMVM_AS_BERTSEKAS  1
10*cd929ea3SAlp Dener #define BLMVM_AS_SIZE       2
11*cd929ea3SAlp Dener 
12*cd929ea3SAlp Dener static const char *BLMVM_AS_TYPE[64] = {"none", "bertsekas"};
13*cd929ea3SAlp Dener 
14*cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoBLMVMEstimateActiveSet(Tao tao, PetscInt asType)
15*cd929ea3SAlp Dener {
16*cd929ea3SAlp Dener   PetscErrorCode               ierr;
17*cd929ea3SAlp Dener   TAO_BLMVM                     *blmP = (TAO_BLMVM *)tao->data;
18*cd929ea3SAlp Dener 
19*cd929ea3SAlp Dener   PetscFunctionBegin;
20*cd929ea3SAlp Dener   if (!tao->bounded) PetscFunctionReturn(0);
21*cd929ea3SAlp Dener   switch (asType) {
22*cd929ea3SAlp Dener   case BLMVM_AS_NONE:
23*cd929ea3SAlp Dener     ierr = ISDestroy(&blmP->inactive_idx);CHKERRQ(ierr);
24*cd929ea3SAlp Dener     ierr = VecWhichInactive(tao->XL, tao->solution, blmP->unprojected_gradient, tao->XU, PETSC_TRUE, &blmP->inactive_idx);CHKERRQ(ierr);
25*cd929ea3SAlp Dener     ierr = ISDestroy(&blmP->active_idx);CHKERRQ(ierr);
26*cd929ea3SAlp Dener     ierr = ISComplementVec(blmP->inactive_idx, tao->solution, &blmP->active_idx);CHKERRQ(ierr);
27*cd929ea3SAlp Dener     break;
28*cd929ea3SAlp Dener 
29*cd929ea3SAlp Dener   case BLMVM_AS_BERTSEKAS:
30*cd929ea3SAlp Dener     /* Use gradient descent to estimate the active set */
31*cd929ea3SAlp Dener     ierr = VecCopy(blmP->unprojected_gradient, blmP->W);CHKERRQ(ierr);
32*cd929ea3SAlp Dener     ierr = VecScale(blmP->W, -1.0);CHKERRQ(ierr);
33*cd929ea3SAlp Dener     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, blmP->unprojected_gradient, blmP->W, blmP->work, blmP->as_step, &blmP->as_tol,
34*cd929ea3SAlp Dener                                    &blmP->active_lower, &blmP->active_upper, &blmP->active_fixed, &blmP->active_idx, &blmP->inactive_idx);CHKERRQ(ierr);
35*cd929ea3SAlp Dener     break;
36*cd929ea3SAlp Dener 
37*cd929ea3SAlp Dener   default:
38*cd929ea3SAlp Dener     break;
39*cd929ea3SAlp Dener   }
40*cd929ea3SAlp Dener   PetscFunctionReturn(0);
41*cd929ea3SAlp Dener }
42*cd929ea3SAlp Dener 
43*cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoBLMVMBoundStep(Tao tao, PetscInt asType, Vec step)
44*cd929ea3SAlp Dener {
45*cd929ea3SAlp Dener   PetscErrorCode               ierr;
46*cd929ea3SAlp Dener   TAO_BLMVM                     *blmP = (TAO_BLMVM *)tao->data;
47*cd929ea3SAlp Dener 
48*cd929ea3SAlp Dener   PetscFunctionBegin;
49*cd929ea3SAlp Dener   if (!tao->bounded) PetscFunctionReturn(0);
50*cd929ea3SAlp Dener   switch (asType) {
51*cd929ea3SAlp Dener   case BLMVM_AS_NONE:
52*cd929ea3SAlp Dener     ierr = VecISSet(step, blmP->active_idx, 0.0);CHKERRQ(ierr);
53*cd929ea3SAlp Dener     break;
54*cd929ea3SAlp Dener 
55*cd929ea3SAlp Dener   case BLMVM_AS_BERTSEKAS:
56*cd929ea3SAlp Dener     ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, blmP->active_lower, blmP->active_upper, blmP->active_fixed, 1.0, step);CHKERRQ(ierr);
57*cd929ea3SAlp Dener     break;
58*cd929ea3SAlp Dener 
59*cd929ea3SAlp Dener   default:
60*cd929ea3SAlp Dener     break;
61*cd929ea3SAlp Dener   }
62*cd929ea3SAlp Dener   PetscFunctionReturn(0);
63*cd929ea3SAlp Dener }
64*cd929ea3SAlp Dener 
65a7e14dcfSSatish Balay /*------------------------------------------------------------*/
66*cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoSolve_BLMVM(Tao tao)
67a7e14dcfSSatish Balay {
68a7e14dcfSSatish Balay   PetscErrorCode               ierr;
69a7e14dcfSSatish Balay   TAO_BLMVM                    *blmP = (TAO_BLMVM *)tao->data;
70e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
71*cd929ea3SAlp Dener   PetscReal                    f, fold, gdx, gnorm, resnorm;
72a7e14dcfSSatish Balay   PetscReal                    stepsize = 1.0,delta;
73*cd929ea3SAlp Dener   PetscInt                     nDiff, nupdates, stepType = BLMVM_STEP_GRAD;
74a7e14dcfSSatish Balay 
75a7e14dcfSSatish Balay   PetscFunctionBegin;
76a7e14dcfSSatish Balay   /*  Project initial point onto bounds */
77a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
78*cd929ea3SAlp Dener   ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
79*cd929ea3SAlp Dener   if (tao->bounded) {
80a7e14dcfSSatish Balay     ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
81*cd929ea3SAlp Dener   }
82a9603a14SPatrick Farrell 
83a7e14dcfSSatish Balay   /* Check convergence criteria */
84a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);CHKERRQ(ierr);
85a7e14dcfSSatish Balay   ierr = VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
86a7e14dcfSSatish Balay 
87*cd929ea3SAlp Dener   ierr = TaoGradientNorm(tao, blmP->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr);
883fcc4fd3STodd Munson   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
89a7e14dcfSSatish Balay 
903ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
91*cd929ea3SAlp Dener   ierr = VecFischer(tao->solution, blmP->unprojected_gradient, tao->XL, tao->XU, blmP->W);CHKERRQ(ierr);
92*cd929ea3SAlp Dener   ierr = VecNorm(blmP->W, NORM_2, &resnorm);CHKERRQ(ierr);
93*cd929ea3SAlp Dener   ierr = TaoLogConvergenceHistory(tao,f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
94*cd929ea3SAlp Dener   ierr = TaoMonitor(tao,tao->niter,f,resnorm,0.0,stepsize);CHKERRQ(ierr);
953ecd9318SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
963ecd9318SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
97a7e14dcfSSatish Balay 
98a7e14dcfSSatish Balay   /* Set initial scaling for the function */
99*cd929ea3SAlp Dener   if (!blmP->H0) {
100*cd929ea3SAlp Dener     delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm);
101*cd929ea3SAlp Dener     ierr = MatLMVMSetJ0Scale(blmP->M, delta);CHKERRQ(ierr);
102a7e14dcfSSatish Balay   }
103*cd929ea3SAlp Dener 
104a7e14dcfSSatish Balay 
105a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps */
106*cd929ea3SAlp Dener   if (!blmP->recycle) {
107*cd929ea3SAlp Dener     blmP->bfgs = 0;
108a7e14dcfSSatish Balay     blmP->grad = 0;
109*cd929ea3SAlp Dener     ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
110*cd929ea3SAlp Dener   }
111a7e14dcfSSatish Balay 
112a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
1133ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
114*cd929ea3SAlp Dener     /* Estimate active set at the current iterate */
115*cd929ea3SAlp Dener     ierr = TaoBLMVMEstimateActiveSet(tao, blmP->as_type);CHKERRQ(ierr);
116*cd929ea3SAlp Dener 
117a7e14dcfSSatish Balay     /* Compute direction */
118*cd929ea3SAlp Dener     if (blmP->H0) {
119*cd929ea3SAlp Dener       ierr = MatLMVMSetJ0(blmP->M, blmP->H0);CHKERRQ(ierr);
120*cd929ea3SAlp Dener       stepType = BLMVM_STEP_BFGS;
121*cd929ea3SAlp Dener     }
122a7e14dcfSSatish Balay     ierr = MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
123*cd929ea3SAlp Dener     ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
124*cd929ea3SAlp Dener     ierr = MatLMVMGetUpdateCount(blmP->M, &nupdates);CHKERRQ(ierr);
125*cd929ea3SAlp Dener     if (nupdates > 0) stepType = BLMVM_STEP_BFGS;
126a7e14dcfSSatish Balay 
127a7e14dcfSSatish Balay     /* Check for success (descent direction) */
128*cd929ea3SAlp Dener     ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
129*cd929ea3SAlp Dener     if (gdx <= 0 || PetscIsInfOrNanReal(gdx)) {
130a7e14dcfSSatish Balay       /* Step is not descent or solve was not successful
131a7e14dcfSSatish Balay          Use steepest descent direction (scaled) */
132a7e14dcfSSatish Balay 
133*cd929ea3SAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm);
134*cd929ea3SAlp Dener       ierr = MatLMVMSetJ0Scale(blmP->M, delta);CHKERRQ(ierr);
135*cd929ea3SAlp Dener       ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
136a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
137*cd929ea3SAlp Dener       ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
138*cd929ea3SAlp Dener       stepType = BLMVM_STEP_GRAD;
139a7e14dcfSSatish Balay     }
140a7e14dcfSSatish Balay     ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
141*cd929ea3SAlp Dener     ierr = TaoBLMVMBoundStep(tao, blmP->as_type, tao->stepdirection);CHKERRQ(ierr);
142a7e14dcfSSatish Balay 
143a7e14dcfSSatish Balay     /* Perform the linesearch */
144a7e14dcfSSatish Balay     fold = f;
145a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, blmP->Xold);CHKERRQ(ierr);
146a7e14dcfSSatish Balay     ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold);CHKERRQ(ierr);
147a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr);
148a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
149a7e14dcfSSatish Balay 
150*cd929ea3SAlp Dener     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && stepType != BLMVM_STEP_GRAD) {
151a7e14dcfSSatish Balay       /* Linesearch failed
152a7e14dcfSSatish Balay          Reset factors and use scaled (projected) gradient step */
153a7e14dcfSSatish Balay       f = fold;
154a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
155a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
156a7e14dcfSSatish Balay 
157*cd929ea3SAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm);
158*cd929ea3SAlp Dener       ierr = MatLMVMSetJ0Scale(blmP->M, delta);CHKERRQ(ierr);
159*cd929ea3SAlp Dener       ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
160a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
161*cd929ea3SAlp Dener       ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
162a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
163*cd929ea3SAlp Dener       ierr = TaoBLMVMBoundStep(tao, blmP->as_type, tao->stepdirection);CHKERRQ(ierr);
164*cd929ea3SAlp Dener       stepType = BLMVM_STEP_GRAD;
165a7e14dcfSSatish Balay 
166a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status);CHKERRQ(ierr);
167a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
168*cd929ea3SAlp Dener     }
169a7e14dcfSSatish Balay 
170a7e14dcfSSatish Balay     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
171*cd929ea3SAlp Dener       /* Line search failed on a gradient step, so just mark reason for divergence */
172*cd929ea3SAlp Dener       f = fold;
173*cd929ea3SAlp Dener       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
174*cd929ea3SAlp Dener       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
175a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
176*cd929ea3SAlp Dener     } else {
177*cd929ea3SAlp Dener       /* LS found valid step, so tally the step type and compute projected gradient */
178*cd929ea3SAlp Dener       switch (stepType) {
179*cd929ea3SAlp Dener       case BLMVM_STEP_BFGS:
180*cd929ea3SAlp Dener         ++blmP->bfgs;
181*cd929ea3SAlp Dener         break;
182*cd929ea3SAlp Dener       case BLMVM_STEP_GRAD:
183*cd929ea3SAlp Dener         ++blmP->grad;
184*cd929ea3SAlp Dener         break;
185*cd929ea3SAlp Dener       default:
186a7e14dcfSSatish Balay         break;
187a7e14dcfSSatish Balay       }
188*cd929ea3SAlp Dener       ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr);
189*cd929ea3SAlp Dener       ierr = TaoGradientNorm(tao, blmP->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr);
190*cd929ea3SAlp Dener       if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
191a7e14dcfSSatish Balay     }
192a7e14dcfSSatish Balay 
193e4cb33bbSBarry Smith     /* Check for converged */
1948931d482SJason Sarich     tao->niter++;
195*cd929ea3SAlp Dener     ierr = VecFischer(tao->solution, blmP->unprojected_gradient, tao->XL, tao->XU, blmP->W);CHKERRQ(ierr);
196*cd929ea3SAlp Dener     ierr = VecNorm(blmP->W, NORM_2, &resnorm);CHKERRQ(ierr);
197*cd929ea3SAlp Dener     ierr = TaoLogConvergenceHistory(tao,f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
198*cd929ea3SAlp Dener     ierr = TaoMonitor(tao,tao->niter,f,resnorm,0.0,stepsize);CHKERRQ(ierr);
1993ecd9318SAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
200a7e14dcfSSatish Balay   }
201a7e14dcfSSatish Balay   PetscFunctionReturn(0);
202a7e14dcfSSatish Balay }
203a7e14dcfSSatish Balay 
204*cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoSetup_BLMVM(Tao tao)
205a7e14dcfSSatish Balay {
206a7e14dcfSSatish Balay   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
207a7e14dcfSSatish Balay   PetscInt       n,N;
208a7e14dcfSSatish Balay   PetscErrorCode ierr;
209a7e14dcfSSatish Balay 
210a7e14dcfSSatish Balay   PetscFunctionBegin;
211a7e14dcfSSatish Balay   /* Existence of tao->solution checked in TaoSetup() */
212a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr);
213a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr);
214302440fdSBarry Smith   ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr);
215*cd929ea3SAlp Dener   ierr = VecDuplicate(tao->solution, &blmP->W);CHKERRQ(ierr);
216*cd929ea3SAlp Dener   ierr = VecDuplicate(tao->solution, &blmP->work);CHKERRQ(ierr);
217a7e14dcfSSatish Balay 
218a7e14dcfSSatish Balay   if (!tao->stepdirection) {
21953506e15SBarry Smith     ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
220a7e14dcfSSatish Balay   }
221a7e14dcfSSatish Balay   if (!tao->gradient) {
222a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
223a7e14dcfSSatish Balay   }
224a7e14dcfSSatish Balay   /* Create matrix for the limited memory approximation */
225a7e14dcfSSatish Balay   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
226a7e14dcfSSatish Balay   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
227*cd929ea3SAlp Dener   ierr = MatSetSizes(blmP->M, n, n, N, N);CHKERRQ(ierr);
228*cd929ea3SAlp Dener   ierr = MatLMVMAllocate(blmP->M,tao->solution,blmP->unprojected_gradient);CHKERRQ(ierr);
229a9603a14SPatrick Farrell 
230a9603a14SPatrick Farrell   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
231a9603a14SPatrick Farrell   if (blmP->H0) {
232*cd929ea3SAlp Dener     ierr = MatLMVMSetJ0(blmP->M, blmP->H0);CHKERRQ(ierr);
233a9603a14SPatrick Farrell   }
234a7e14dcfSSatish Balay   PetscFunctionReturn(0);
235a7e14dcfSSatish Balay }
236a7e14dcfSSatish Balay 
237a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
238*cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoDestroy_BLMVM(Tao tao)
239a7e14dcfSSatish Balay {
240a7e14dcfSSatish Balay   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
241a7e14dcfSSatish Balay   PetscErrorCode ierr;
242a7e14dcfSSatish Balay 
243a7e14dcfSSatish Balay   PetscFunctionBegin;
244a7e14dcfSSatish Balay   if (tao->setupcalled) {
245a7e14dcfSSatish Balay     ierr = MatDestroy(&blmP->M);CHKERRQ(ierr);
246a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr);
247a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr);
248a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr);
249*cd929ea3SAlp Dener     ierr = VecDestroy(&blmP->W);CHKERRQ(ierr);
250*cd929ea3SAlp Dener     ierr = VecDestroy(&blmP->work);CHKERRQ(ierr);
251a7e14dcfSSatish Balay   }
252*cd929ea3SAlp Dener   ierr = ISDestroy(&blmP->active_lower);CHKERRQ(ierr);
253*cd929ea3SAlp Dener   ierr = ISDestroy(&blmP->active_upper);CHKERRQ(ierr);
254*cd929ea3SAlp Dener   ierr = ISDestroy(&blmP->active_fixed);CHKERRQ(ierr);
255*cd929ea3SAlp Dener   ierr = ISDestroy(&blmP->active_idx);CHKERRQ(ierr);
256*cd929ea3SAlp Dener   ierr = ISDestroy(&blmP->inactive_idx);CHKERRQ(ierr);
257a9603a14SPatrick Farrell   if (blmP->H0) {
258a9603a14SPatrick Farrell     PetscObjectDereference((PetscObject)blmP->H0);
259a9603a14SPatrick Farrell   }
260a9603a14SPatrick Farrell 
261a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
262a7e14dcfSSatish Balay   PetscFunctionReturn(0);
263a7e14dcfSSatish Balay }
264a7e14dcfSSatish Balay 
265a7e14dcfSSatish Balay /*------------------------------------------------------------*/
266*cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
267a7e14dcfSSatish Balay {
268*cd929ea3SAlp Dener   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
269a7e14dcfSSatish Balay   PetscErrorCode ierr;
270a7e14dcfSSatish Balay 
271a7e14dcfSSatish Balay   PetscFunctionBegin;
2721a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr);
273*cd929ea3SAlp Dener   ierr = PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL);CHKERRQ(ierr);
274*cd929ea3SAlp Dener   ierr = PetscOptionsEList("-tao_blmvm_as_type","active set estimation method", "", BLMVM_AS_TYPE, BLMVM_AS_SIZE, BLMVM_AS_TYPE[blmP->as_type], &blmP->as_type,NULL);CHKERRQ(ierr);
275*cd929ea3SAlp Dener   ierr = PetscOptionsReal("-tao_blmvm_as_tol", "initial tolerance used when estimating actively bounded variables","",blmP->as_tol,&blmP->as_tol,NULL);CHKERRQ(ierr);
276*cd929ea3SAlp Dener   ierr = PetscOptionsReal("-tao_blmvm_as_step", "step length used when estimating actively bounded variables","",blmP->as_step,&blmP->as_step,NULL);CHKERRQ(ierr);
277a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
278a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
279a7e14dcfSSatish Balay   PetscFunctionReturn(0);
280a7e14dcfSSatish Balay }
281a7e14dcfSSatish Balay 
282a7e14dcfSSatish Balay 
283a7e14dcfSSatish Balay /*------------------------------------------------------------*/
284*cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer)
285a7e14dcfSSatish Balay {
286a7e14dcfSSatish Balay   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
287a7e14dcfSSatish Balay   PetscBool      isascii;
288a7e14dcfSSatish Balay   PetscErrorCode ierr;
289*cd929ea3SAlp Dener   PetscInt       recycled_its;
290a7e14dcfSSatish Balay 
291a7e14dcfSSatish Balay   PetscFunctionBegin;
292a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
293a7e14dcfSSatish Balay   if (isascii) {
294a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "  Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr);
295*cd929ea3SAlp Dener     if (lmP->recycle) {
296*cd929ea3SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "  Recycle: on\n");CHKERRQ(ierr);
297*cd929ea3SAlp Dener       recycled_its = lmP->bfgs + lmP->grad;
298*cd929ea3SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "  Total recycled iterations: %D\n", recycled_its);CHKERRQ(ierr);
299*cd929ea3SAlp Dener     }
300a7e14dcfSSatish Balay   }
301a7e14dcfSSatish Balay   PetscFunctionReturn(0);
302a7e14dcfSSatish Balay }
303a7e14dcfSSatish Balay 
304*cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
305a7e14dcfSSatish Balay {
306a7e14dcfSSatish Balay   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
307a7e14dcfSSatish Balay   PetscErrorCode ierr;
308a7e14dcfSSatish Balay 
309a7e14dcfSSatish Balay   PetscFunctionBegin;
310441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
311a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
312a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
31353506e15SBarry Smith   if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
314a7e14dcfSSatish Balay 
315a7e14dcfSSatish Balay   ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr);
316a7e14dcfSSatish Balay   ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr);
317a7e14dcfSSatish Balay   ierr = VecSet(DXU,0.0);CHKERRQ(ierr);
318a7e14dcfSSatish Balay   ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr);
319a7e14dcfSSatish Balay 
320a7e14dcfSSatish Balay   ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr);
321a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr);
322a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr);
323a7e14dcfSSatish Balay   PetscFunctionReturn(0);
324a7e14dcfSSatish Balay }
325a7e14dcfSSatish Balay 
326a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
3271522df2eSJason Sarich /*MC
3281522df2eSJason Sarich   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
3291522df2eSJason Sarich          for nonlinear minimization with bound constraints. It is an extension
3301522df2eSJason Sarich          of TAOLMVM
3311522df2eSJason Sarich 
3321eb8069cSJason Sarich   Level: beginner
3331522df2eSJason Sarich M*/
334728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
335a7e14dcfSSatish Balay {
336a7e14dcfSSatish Balay   TAO_BLMVM      *blmP;
3378caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
338a7e14dcfSSatish Balay   PetscErrorCode ierr;
339a7e14dcfSSatish Balay 
340a7e14dcfSSatish Balay   PetscFunctionBegin;
341a7e14dcfSSatish Balay   tao->ops->setup = TaoSetup_BLMVM;
342a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_BLMVM;
343a7e14dcfSSatish Balay   tao->ops->view = TaoView_BLMVM;
344a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
345a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_BLMVM;
346a7e14dcfSSatish Balay   tao->ops->computedual = TaoComputeDual_BLMVM;
347a7e14dcfSSatish Balay 
3483c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr);
349a9603a14SPatrick Farrell   blmP->H0 = NULL;
350*cd929ea3SAlp Dener   blmP->as_step = 0.001;
351*cd929ea3SAlp Dener   blmP->as_tol = 0.001;
352*cd929ea3SAlp Dener   blmP->as_type = BLMVM_AS_BERTSEKAS;
353*cd929ea3SAlp Dener   blmP->recycle = PETSC_FALSE;
354*cd929ea3SAlp Dener 
355a7e14dcfSSatish Balay   tao->data = (void*)blmP;
3566552cf8aSJason Sarich 
3576552cf8aSJason Sarich   /* Override default settings (unless already changed) */
3586552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 2000;
3596552cf8aSJason Sarich   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
360a7e14dcfSSatish Balay 
361a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
36263b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
363a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
364441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
3655d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
366*cd929ea3SAlp Dener 
367*cd929ea3SAlp Dener   ierr = KSPInitializePackage();CHKERRQ(ierr);
368*cd929ea3SAlp Dener   ierr = MatCreate(((PetscObject)tao)->comm, &blmP->M);CHKERRQ(ierr);
369*cd929ea3SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1);CHKERRQ(ierr);
370*cd929ea3SAlp Dener   ierr = MatSetType(blmP->M, MATLBFGS);CHKERRQ(ierr);
371*cd929ea3SAlp Dener   ierr = MatSetOptionsPrefix(blmP->M, "tao_blmvm_");CHKERRQ(ierr);
372a7e14dcfSSatish Balay   PetscFunctionReturn(0);
373a7e14dcfSSatish Balay }
374a7e14dcfSSatish Balay 
375*cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode TaoBLMVMSetH0(Tao tao, Mat H0)
376a9603a14SPatrick Farrell {
377a9603a14SPatrick Farrell   TAO_LMVM       *lmP;
378a9603a14SPatrick Farrell   TAO_BLMVM      *blmP;
379b625d6c7SJed Brown   TaoType        type;
380a9603a14SPatrick Farrell   PetscBool      is_lmvm, is_blmvm;
381a9603a14SPatrick Farrell   PetscErrorCode ierr;
382a9603a14SPatrick Farrell 
383a9603a14SPatrick Farrell   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
384a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
385a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
386a9603a14SPatrick Farrell 
387a9603a14SPatrick Farrell   if (is_lmvm) {
388a9603a14SPatrick Farrell     lmP = (TAO_LMVM *)tao->data;
389a9603a14SPatrick Farrell     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
390a9603a14SPatrick Farrell     lmP->H0 = H0;
391a9603a14SPatrick Farrell   } else if (is_blmvm) {
392a9603a14SPatrick Farrell     blmP = (TAO_BLMVM *)tao->data;
393a9603a14SPatrick Farrell     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
394a9603a14SPatrick Farrell     blmP->H0 = H0;
39574c66251SPatrick Farrell   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
396a9603a14SPatrick Farrell   PetscFunctionReturn(0);
397a9603a14SPatrick Farrell }
398a9603a14SPatrick Farrell 
399*cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode TaoBLMVMGetH0(Tao tao, Mat *H0)
400a9603a14SPatrick Farrell {
401a9603a14SPatrick Farrell   TAO_LMVM       *lmP;
402a9603a14SPatrick Farrell   TAO_BLMVM      *blmP;
403b625d6c7SJed Brown   TaoType        type;
404a9603a14SPatrick Farrell   PetscBool      is_lmvm, is_blmvm;
405a9603a14SPatrick Farrell   Mat            M;
406a9603a14SPatrick Farrell 
407a9603a14SPatrick Farrell   PetscErrorCode ierr;
408a9603a14SPatrick Farrell 
409a9603a14SPatrick Farrell   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
410a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
411a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
412a9603a14SPatrick Farrell 
413a9603a14SPatrick Farrell   if (is_lmvm) {
414a9603a14SPatrick Farrell     lmP = (TAO_LMVM *)tao->data;
415a9603a14SPatrick Farrell     M = lmP->M;
416a9603a14SPatrick Farrell   } else if (is_blmvm) {
417a9603a14SPatrick Farrell     blmP = (TAO_BLMVM *)tao->data;
418a9603a14SPatrick Farrell     M = blmP->M;
41974c66251SPatrick Farrell   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
420*cd929ea3SAlp Dener   ierr = MatLMVMGetJ0(M, H0);CHKERRQ(ierr);
421a9603a14SPatrick Farrell   PetscFunctionReturn(0);
422a9603a14SPatrick Farrell }
423a9603a14SPatrick Farrell 
424*cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode TaoBLMVMGetH0KSP(Tao tao, KSP *ksp)
425a9603a14SPatrick Farrell {
426a9603a14SPatrick Farrell   TAO_LMVM       *lmP;
427a9603a14SPatrick Farrell   TAO_BLMVM      *blmP;
428b625d6c7SJed Brown   TaoType        type;
429a9603a14SPatrick Farrell   PetscBool      is_lmvm, is_blmvm;
430a9603a14SPatrick Farrell   Mat            M;
431a9603a14SPatrick Farrell   PetscErrorCode ierr;
432a9603a14SPatrick Farrell 
433a9603a14SPatrick Farrell   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
434a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
435a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
436a9603a14SPatrick Farrell 
437a9603a14SPatrick Farrell   if (is_lmvm) {
438a9603a14SPatrick Farrell     lmP = (TAO_LMVM *)tao->data;
439a9603a14SPatrick Farrell     M = lmP->M;
440a9603a14SPatrick Farrell   } else if (is_blmvm) {
441a9603a14SPatrick Farrell     blmP = (TAO_BLMVM *)tao->data;
442a9603a14SPatrick Farrell     M = blmP->M;
44374c66251SPatrick Farrell   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
444*cd929ea3SAlp Dener   ierr = MatLMVMGetJ0KSP(M, ksp);CHKERRQ(ierr);
445a9603a14SPatrick Farrell   PetscFunctionReturn(0);
446a9603a14SPatrick Farrell }
447