xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision 5aff1b4e483068a430fabe11c633beda7572be55)
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 
5cd929ea3SAlp Dener #define BLMVM_STEP_BFGS     0
6cd929ea3SAlp Dener #define BLMVM_STEP_GRAD     1
7cd929ea3SAlp Dener 
8cd929ea3SAlp Dener #define BLMVM_AS_NONE       0
9cd929ea3SAlp Dener #define BLMVM_AS_BERTSEKAS  1
10cd929ea3SAlp Dener #define BLMVM_AS_SIZE       2
11cd929ea3SAlp Dener 
12cd929ea3SAlp Dener static const char *BLMVM_AS_TYPE[64] = {"none", "bertsekas"};
13cd929ea3SAlp Dener 
14cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoBLMVMEstimateActiveSet(Tao tao, PetscInt asType)
15cd929ea3SAlp Dener {
16cd929ea3SAlp Dener   PetscErrorCode               ierr;
17cd929ea3SAlp Dener   TAO_BLMVM                     *blmP = (TAO_BLMVM *)tao->data;
18cd929ea3SAlp Dener 
19cd929ea3SAlp Dener   PetscFunctionBegin;
20cd929ea3SAlp Dener   if (!tao->bounded) PetscFunctionReturn(0);
21cd929ea3SAlp Dener   switch (asType) {
22cd929ea3SAlp Dener   case BLMVM_AS_NONE:
23cd929ea3SAlp Dener     ierr = ISDestroy(&blmP->inactive_idx);CHKERRQ(ierr);
24cd929ea3SAlp Dener     ierr = VecWhichInactive(tao->XL, tao->solution, blmP->unprojected_gradient, tao->XU, PETSC_TRUE, &blmP->inactive_idx);CHKERRQ(ierr);
25cd929ea3SAlp Dener     ierr = ISDestroy(&blmP->active_idx);CHKERRQ(ierr);
26cd929ea3SAlp Dener     ierr = ISComplementVec(blmP->inactive_idx, tao->solution, &blmP->active_idx);CHKERRQ(ierr);
27cd929ea3SAlp Dener     break;
28cd929ea3SAlp Dener 
29cd929ea3SAlp Dener   case BLMVM_AS_BERTSEKAS:
30cd929ea3SAlp Dener     /* Use gradient descent to estimate the active set */
31cd929ea3SAlp Dener     ierr = VecCopy(blmP->unprojected_gradient, blmP->W);CHKERRQ(ierr);
32cd929ea3SAlp Dener     ierr = VecScale(blmP->W, -1.0);CHKERRQ(ierr);
33cd929ea3SAlp Dener     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, blmP->unprojected_gradient, blmP->W, blmP->work, blmP->as_step, &blmP->as_tol,
34cd929ea3SAlp Dener                                    &blmP->active_lower, &blmP->active_upper, &blmP->active_fixed, &blmP->active_idx, &blmP->inactive_idx);CHKERRQ(ierr);
35cd929ea3SAlp Dener     break;
36cd929ea3SAlp Dener 
37cd929ea3SAlp Dener   default:
38cd929ea3SAlp Dener     break;
39cd929ea3SAlp Dener   }
40cd929ea3SAlp Dener   PetscFunctionReturn(0);
41cd929ea3SAlp Dener }
42cd929ea3SAlp Dener 
43cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoBLMVMBoundStep(Tao tao, PetscInt asType, Vec step)
44cd929ea3SAlp Dener {
45cd929ea3SAlp Dener   PetscErrorCode               ierr;
46cd929ea3SAlp Dener   TAO_BLMVM                     *blmP = (TAO_BLMVM *)tao->data;
47cd929ea3SAlp Dener 
48cd929ea3SAlp Dener   PetscFunctionBegin;
49cd929ea3SAlp Dener   if (!tao->bounded) PetscFunctionReturn(0);
50cd929ea3SAlp Dener   switch (asType) {
51cd929ea3SAlp Dener   case BLMVM_AS_NONE:
52cd929ea3SAlp Dener     ierr = VecISSet(step, blmP->active_idx, 0.0);CHKERRQ(ierr);
53cd929ea3SAlp Dener     break;
54cd929ea3SAlp Dener 
55cd929ea3SAlp Dener   case BLMVM_AS_BERTSEKAS:
56cd929ea3SAlp Dener     ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, blmP->active_lower, blmP->active_upper, blmP->active_fixed, 1.0, step);CHKERRQ(ierr);
57cd929ea3SAlp Dener     break;
58cd929ea3SAlp Dener 
59cd929ea3SAlp Dener   default:
60cd929ea3SAlp Dener     break;
61cd929ea3SAlp Dener   }
62cd929ea3SAlp Dener   PetscFunctionReturn(0);
63cd929ea3SAlp Dener }
64cd929ea3SAlp Dener 
65a7e14dcfSSatish Balay /*------------------------------------------------------------*/
66cd929ea3SAlp 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;
71cd929ea3SAlp Dener   PetscReal                    f, fold, gdx, gnorm, resnorm;
72b2d8c577SAlp Dener   PetscReal                    stepsize = 1.0;
73cd929ea3SAlp 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);
78cd929ea3SAlp Dener   ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
79cd929ea3SAlp Dener   if (tao->bounded) {
80a7e14dcfSSatish Balay     ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
81cd929ea3SAlp 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 
87cd929ea3SAlp 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;
91cd929ea3SAlp Dener   ierr = VecFischer(tao->solution, blmP->unprojected_gradient, tao->XL, tao->XU, blmP->W);CHKERRQ(ierr);
92cd929ea3SAlp Dener   ierr = VecNorm(blmP->W, NORM_2, &resnorm);CHKERRQ(ierr);
93cd929ea3SAlp Dener   ierr = TaoLogConvergenceHistory(tao,f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
94cd929ea3SAlp 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 counter for gradient/reset steps */
99cd929ea3SAlp Dener   if (!blmP->recycle) {
100cd929ea3SAlp Dener     blmP->bfgs = 0;
101a7e14dcfSSatish Balay     blmP->grad = 0;
102cd929ea3SAlp Dener     ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
103cd929ea3SAlp Dener   }
104a7e14dcfSSatish Balay 
105a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
1063ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
107cd929ea3SAlp Dener     /* Estimate active set at the current iterate */
108cd929ea3SAlp Dener     ierr = TaoBLMVMEstimateActiveSet(tao, blmP->as_type);CHKERRQ(ierr);
109cd929ea3SAlp Dener 
110a7e14dcfSSatish Balay     /* Compute direction */
111cd929ea3SAlp Dener     if (blmP->H0) {
112cd929ea3SAlp Dener       ierr = MatLMVMSetJ0(blmP->M, blmP->H0);CHKERRQ(ierr);
113cd929ea3SAlp Dener       stepType = BLMVM_STEP_BFGS;
1140ad3a497SAlp Dener     } else if (!blmP->no_scale) {
1150ad3a497SAlp Dener       ierr = MatLMVMSetJ0(blmP->M, blmP->Mscale);CHKERRQ(ierr);
116cd929ea3SAlp Dener     }
117b2d8c577SAlp Dener     ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
118b2d8c577SAlp Dener     ierr = VecCopy(blmP->unprojected_gradient, blmP->red_grad);CHKERRQ(ierr);
119b2d8c577SAlp Dener     ierr = VecISSet(blmP->red_grad, blmP->active_idx, 0.0);CHKERRQ(ierr);
120b2d8c577SAlp Dener     ierr = MatSolve(blmP->M, blmP->red_grad, tao->stepdirection);CHKERRQ(ierr);
121b2d8c577SAlp Dener     ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
122b2d8c577SAlp Dener     ierr = TaoBLMVMBoundStep(tao, blmP->as_type, tao->stepdirection);CHKERRQ(ierr);
123cd929ea3SAlp Dener     ierr = MatLMVMGetUpdateCount(blmP->M, &nupdates);CHKERRQ(ierr);
124cd929ea3SAlp Dener     if (nupdates > 0) stepType = BLMVM_STEP_BFGS;
125a7e14dcfSSatish Balay 
126a7e14dcfSSatish Balay     /* Check for success (descent direction) */
127b2d8c577SAlp Dener     ierr = VecDot(tao->stepdirection, blmP->red_grad, &gdx);CHKERRQ(ierr);
128b2d8c577SAlp Dener     if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) {
129a7e14dcfSSatish Balay       /* Step is not descent or solve was not successful
130a7e14dcfSSatish Balay          Use steepest descent direction (scaled) */
131cd929ea3SAlp Dener       ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
1320ad3a497SAlp Dener       ierr = MatLMVMClearJ0(blmP->M);CHKERRQ(ierr);
133a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
134cd929ea3SAlp Dener       ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
135a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
136cd929ea3SAlp Dener       ierr = TaoBLMVMBoundStep(tao, blmP->as_type, tao->stepdirection);CHKERRQ(ierr);
137b2d8c577SAlp Dener       stepType = BLMVM_STEP_GRAD;
138b2d8c577SAlp Dener     }
139a7e14dcfSSatish Balay 
140a7e14dcfSSatish Balay     /* Perform the linesearch */
141a7e14dcfSSatish Balay     fold = f;
142a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, blmP->Xold);CHKERRQ(ierr);
143a7e14dcfSSatish Balay     ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold);CHKERRQ(ierr);
144a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr);
145a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
146a7e14dcfSSatish Balay 
147cd929ea3SAlp Dener     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && stepType != BLMVM_STEP_GRAD) {
148a7e14dcfSSatish Balay       /* Linesearch failed
149a7e14dcfSSatish Balay          Reset factors and use scaled (projected) gradient step */
150a7e14dcfSSatish Balay       f = fold;
151a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
152a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
153a7e14dcfSSatish Balay 
154cd929ea3SAlp Dener       ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
1550ad3a497SAlp Dener       ierr = MatLMVMClearJ0(blmP->M);CHKERRQ(ierr);
156a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
157cd929ea3SAlp Dener       ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
158a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
159cd929ea3SAlp Dener       ierr = TaoBLMVMBoundStep(tao, blmP->as_type, tao->stepdirection);CHKERRQ(ierr);
160cd929ea3SAlp Dener       stepType = BLMVM_STEP_GRAD;
161a7e14dcfSSatish Balay 
162a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status);CHKERRQ(ierr);
163a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
164cd929ea3SAlp Dener     }
165a7e14dcfSSatish Balay 
166a7e14dcfSSatish Balay     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
167cd929ea3SAlp Dener       /* Line search failed on a gradient step, so just mark reason for divergence */
168cd929ea3SAlp Dener       f = fold;
169cd929ea3SAlp Dener       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
170cd929ea3SAlp Dener       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
171a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
172cd929ea3SAlp Dener     } else {
173cd929ea3SAlp Dener       /* LS found valid step, so tally the step type and compute projected gradient */
174cd929ea3SAlp Dener       switch (stepType) {
175cd929ea3SAlp Dener       case BLMVM_STEP_BFGS:
176cd929ea3SAlp Dener         ++blmP->bfgs;
177cd929ea3SAlp Dener         break;
178cd929ea3SAlp Dener       case BLMVM_STEP_GRAD:
179cd929ea3SAlp Dener         ++blmP->grad;
180cd929ea3SAlp Dener         break;
181cd929ea3SAlp Dener       default:
182a7e14dcfSSatish Balay         break;
183a7e14dcfSSatish Balay       }
184cd929ea3SAlp Dener       ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr);
185cd929ea3SAlp Dener       ierr = TaoGradientNorm(tao, blmP->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr);
186cd929ea3SAlp Dener       if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
187a7e14dcfSSatish Balay     }
188a7e14dcfSSatish Balay 
189e4cb33bbSBarry Smith     /* Check for converged */
1908931d482SJason Sarich     tao->niter++;
191cd929ea3SAlp Dener     ierr = VecFischer(tao->solution, blmP->unprojected_gradient, tao->XL, tao->XU, blmP->W);CHKERRQ(ierr);
192cd929ea3SAlp Dener     ierr = VecNorm(blmP->W, NORM_2, &resnorm);CHKERRQ(ierr);
193cd929ea3SAlp Dener     ierr = TaoLogConvergenceHistory(tao,f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
194cd929ea3SAlp Dener     ierr = TaoMonitor(tao,tao->niter,f,resnorm,0.0,stepsize);CHKERRQ(ierr);
1953ecd9318SAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
196a7e14dcfSSatish Balay   }
197a7e14dcfSSatish Balay   PetscFunctionReturn(0);
198a7e14dcfSSatish Balay }
199a7e14dcfSSatish Balay 
200cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoSetup_BLMVM(Tao tao)
201a7e14dcfSSatish Balay {
202a7e14dcfSSatish Balay   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
203a7e14dcfSSatish Balay   PetscInt       n,N;
204a7e14dcfSSatish Balay   PetscErrorCode ierr;
205a7e14dcfSSatish Balay 
206a7e14dcfSSatish Balay   PetscFunctionBegin;
207a7e14dcfSSatish Balay   /* Existence of tao->solution checked in TaoSetup() */
208a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr);
209a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr);
210302440fdSBarry Smith   ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr);
211cd929ea3SAlp Dener   ierr = VecDuplicate(tao->solution, &blmP->W);CHKERRQ(ierr);
212cd929ea3SAlp Dener   ierr = VecDuplicate(tao->solution, &blmP->work);CHKERRQ(ierr);
213b2d8c577SAlp Dener   ierr = VecDuplicate(tao->solution, &blmP->red_grad);CHKERRQ(ierr);
214a7e14dcfSSatish Balay 
215a7e14dcfSSatish Balay   if (!tao->stepdirection) {
21653506e15SBarry Smith     ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
217a7e14dcfSSatish Balay   }
218a7e14dcfSSatish Balay   if (!tao->gradient) {
219a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
220a7e14dcfSSatish Balay   }
221a7e14dcfSSatish Balay   /* Create matrix for the limited memory approximation */
222a7e14dcfSSatish Balay   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
223a7e14dcfSSatish Balay   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
224cd929ea3SAlp Dener   ierr = MatSetSizes(blmP->M, n, n, N, N);CHKERRQ(ierr);
225cd929ea3SAlp Dener   ierr = MatLMVMAllocate(blmP->M,tao->solution,blmP->unprojected_gradient);CHKERRQ(ierr);
226a9603a14SPatrick Farrell 
227a9603a14SPatrick Farrell   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
228a9603a14SPatrick Farrell   if (blmP->H0) {
229cd929ea3SAlp Dener     ierr = MatLMVMSetJ0(blmP->M, blmP->H0);CHKERRQ(ierr);
2300ad3a497SAlp Dener   } else if (!blmP->no_scale) {
2310ad3a497SAlp Dener     if (!blmP->Mscale) {
232*5aff1b4eSAlp Dener       ierr = MatCreate(PetscObjectComm((PetscObject)blmP->M), &blmP->Mscale);CHKERRQ(ierr);
233*5aff1b4eSAlp Dener       ierr = MatSetType(blmP->Mscale, MATLMVMDIAGBRDN);CHKERRQ(ierr);
2340ad3a497SAlp Dener       ierr = MatSetOptionsPrefix(blmP->Mscale, "tao_blmvm_scale_");CHKERRQ(ierr);
235*5aff1b4eSAlp Dener       ierr = MatSetFromOptions(blmP->Mscale);CHKERRQ(ierr);
2360ad3a497SAlp Dener       ierr = MatLMVMAllocate(blmP->Mscale, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
2370ad3a497SAlp Dener     }
2380ad3a497SAlp Dener     ierr = MatLMVMSetJ0(blmP->M, blmP->Mscale);CHKERRQ(ierr);
239a9603a14SPatrick Farrell   }
240a7e14dcfSSatish Balay   PetscFunctionReturn(0);
241a7e14dcfSSatish Balay }
242a7e14dcfSSatish Balay 
243a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
244cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoDestroy_BLMVM(Tao tao)
245a7e14dcfSSatish Balay {
246a7e14dcfSSatish Balay   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
247a7e14dcfSSatish Balay   PetscErrorCode ierr;
248a7e14dcfSSatish Balay 
249a7e14dcfSSatish Balay   PetscFunctionBegin;
250a7e14dcfSSatish Balay   if (tao->setupcalled) {
251a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr);
252a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr);
253a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr);
254cd929ea3SAlp Dener     ierr = VecDestroy(&blmP->W);CHKERRQ(ierr);
255cd929ea3SAlp Dener     ierr = VecDestroy(&blmP->work);CHKERRQ(ierr);
256b2d8c577SAlp Dener     ierr = VecDestroy(&blmP->red_grad);CHKERRQ(ierr);
257a7e14dcfSSatish Balay   }
258cd929ea3SAlp Dener   ierr = ISDestroy(&blmP->active_lower);CHKERRQ(ierr);
259cd929ea3SAlp Dener   ierr = ISDestroy(&blmP->active_upper);CHKERRQ(ierr);
260cd929ea3SAlp Dener   ierr = ISDestroy(&blmP->active_fixed);CHKERRQ(ierr);
261cd929ea3SAlp Dener   ierr = ISDestroy(&blmP->active_idx);CHKERRQ(ierr);
262cd929ea3SAlp Dener   ierr = ISDestroy(&blmP->inactive_idx);CHKERRQ(ierr);
263b2d8c577SAlp Dener   ierr = MatDestroy(&blmP->M);CHKERRQ(ierr);
264a9603a14SPatrick Farrell   if (blmP->H0) {
265a9603a14SPatrick Farrell     PetscObjectDereference((PetscObject)blmP->H0);
266a9603a14SPatrick Farrell   }
2670ad3a497SAlp Dener   if (blmP->Mscale) {
2680ad3a497SAlp Dener     ierr = MatDestroy(&blmP->Mscale);CHKERRQ(ierr);
2690ad3a497SAlp Dener   }
270a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
271a7e14dcfSSatish Balay   PetscFunctionReturn(0);
272a7e14dcfSSatish Balay }
273a7e14dcfSSatish Balay 
274a7e14dcfSSatish Balay /*------------------------------------------------------------*/
275cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
276a7e14dcfSSatish Balay {
277cd929ea3SAlp Dener   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
278a7e14dcfSSatish Balay   PetscErrorCode ierr;
2797b1c7716SAlp Dener   PetscBool      is_lmvm, is_spd;
280a7e14dcfSSatish Balay 
281a7e14dcfSSatish Balay   PetscFunctionBegin;
2821a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr);
283cd929ea3SAlp Dener   ierr = PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL);CHKERRQ(ierr);
2840ad3a497SAlp Dener   ierr = PetscOptionsBool("-tao_blmvm_no_scale","(developer) disable the diagonal Broyden scaling of the BFGS approximation","",blmP->no_scale,&blmP->no_scale,NULL);CHKERRQ(ierr);
285cd929ea3SAlp 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);
286cd929ea3SAlp Dener   ierr = PetscOptionsReal("-tao_blmvm_as_tol", "initial tolerance used when estimating actively bounded variables","",blmP->as_tol,&blmP->as_tol,NULL);CHKERRQ(ierr);
287cd929ea3SAlp Dener   ierr = PetscOptionsReal("-tao_blmvm_as_step", "step length used when estimating actively bounded variables","",blmP->as_step,&blmP->as_step,NULL);CHKERRQ(ierr);
2887b1c7716SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
289a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
290b2d8c577SAlp Dener   ierr = MatSetFromOptions(blmP->M);CHKERRQ(ierr);
291*5aff1b4eSAlp Dener   if (!blmP->no_scale) {
292*5aff1b4eSAlp Dener   }
2937b1c7716SAlp Dener   ierr = PetscObjectBaseTypeCompare((PetscObject)blmP->M, MATLMVM, &is_lmvm);CHKERRQ(ierr);
2947b1c7716SAlp Dener   if (!is_lmvm) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Matrix must be an LMVM-type");
2957b1c7716SAlp Dener   ierr = MatGetOption(blmP->M, MAT_SPD, &is_spd);CHKERRQ(ierr);
2967b1c7716SAlp Dener   if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite");
297a7e14dcfSSatish Balay   PetscFunctionReturn(0);
298a7e14dcfSSatish Balay }
299a7e14dcfSSatish Balay 
300a7e14dcfSSatish Balay /*------------------------------------------------------------*/
301cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer)
302a7e14dcfSSatish Balay {
303a7e14dcfSSatish Balay   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
304a7e14dcfSSatish Balay   PetscBool      isascii;
305a7e14dcfSSatish Balay   PetscErrorCode ierr;
306cd929ea3SAlp Dener   PetscInt       recycled_its;
307a7e14dcfSSatish Balay 
308a7e14dcfSSatish Balay   PetscFunctionBegin;
309a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
310a7e14dcfSSatish Balay   if (isascii) {
311a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "  Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr);
312cd929ea3SAlp Dener     if (lmP->recycle) {
313cd929ea3SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "  Recycle: on\n");CHKERRQ(ierr);
314cd929ea3SAlp Dener       recycled_its = lmP->bfgs + lmP->grad;
315cd929ea3SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "  Total recycled iterations: %D\n", recycled_its);CHKERRQ(ierr);
316cd929ea3SAlp Dener     }
317b2d8c577SAlp Dener     ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
318b2d8c577SAlp Dener     ierr = MatView(lmP->M, viewer);CHKERRQ(ierr);
319b2d8c577SAlp Dener     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
320a7e14dcfSSatish Balay   }
321a7e14dcfSSatish Balay   PetscFunctionReturn(0);
322a7e14dcfSSatish Balay }
323a7e14dcfSSatish Balay 
324cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
325a7e14dcfSSatish Balay {
326a7e14dcfSSatish Balay   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
327a7e14dcfSSatish Balay   PetscErrorCode ierr;
328a7e14dcfSSatish Balay 
329a7e14dcfSSatish Balay   PetscFunctionBegin;
330441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
331a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
332a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
33353506e15SBarry 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");
334a7e14dcfSSatish Balay 
335a7e14dcfSSatish Balay   ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr);
336a7e14dcfSSatish Balay   ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr);
337a7e14dcfSSatish Balay   ierr = VecSet(DXU,0.0);CHKERRQ(ierr);
338a7e14dcfSSatish Balay   ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr);
339a7e14dcfSSatish Balay 
340a7e14dcfSSatish Balay   ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr);
341a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr);
342a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr);
343a7e14dcfSSatish Balay   PetscFunctionReturn(0);
344a7e14dcfSSatish Balay }
345a7e14dcfSSatish Balay 
346a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
3471522df2eSJason Sarich /*MC
3481522df2eSJason Sarich   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
3491522df2eSJason Sarich          for nonlinear minimization with bound constraints. It is an extension
3501522df2eSJason Sarich          of TAOLMVM
3511522df2eSJason Sarich 
3520ad3a497SAlp Dener   Options Database Keys:
3530ad3a497SAlp Dener .   -tao_blmvm_recycle - enable recycling LMVM updates between TaoSolve() calls
3540ad3a497SAlp Dener .   -tao_blmvm_no_scale - (developer) disables diagonal Broyden scaling on the LMVM approximation
3550ad3a497SAlp Dener .   -tao_blmvm_as_type - (developer) active set estimation method <none,bertsekas>
3560ad3a497SAlp Dener .   -tao_blmvm_as_tol - (developer) initial distance tolerance used when estimating the active set
3570ad3a497SAlp Dener .   -tao_blmvm_as_step - (developer) trial step length used when estimating the active set
3580ad3a497SAlp Dener 
3591eb8069cSJason Sarich   Level: beginner
3601522df2eSJason Sarich M*/
361728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
362a7e14dcfSSatish Balay {
363a7e14dcfSSatish Balay   TAO_BLMVM      *blmP;
3648caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
365a7e14dcfSSatish Balay   PetscErrorCode ierr;
366a7e14dcfSSatish Balay 
367a7e14dcfSSatish Balay   PetscFunctionBegin;
368a7e14dcfSSatish Balay   tao->ops->setup = TaoSetup_BLMVM;
369a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_BLMVM;
370a7e14dcfSSatish Balay   tao->ops->view = TaoView_BLMVM;
371a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
372a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_BLMVM;
373a7e14dcfSSatish Balay   tao->ops->computedual = TaoComputeDual_BLMVM;
374a7e14dcfSSatish Balay 
3753c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr);
376a9603a14SPatrick Farrell   blmP->H0 = NULL;
377cd929ea3SAlp Dener   blmP->as_step = 0.001;
378cd929ea3SAlp Dener   blmP->as_tol = 0.001;
379cd929ea3SAlp Dener   blmP->as_type = BLMVM_AS_BERTSEKAS;
380cd929ea3SAlp Dener   blmP->recycle = PETSC_FALSE;
3810ad3a497SAlp Dener   blmP->no_scale = PETSC_FALSE;
382cd929ea3SAlp Dener 
383a7e14dcfSSatish Balay   tao->data = (void*)blmP;
3846552cf8aSJason Sarich 
3856552cf8aSJason Sarich   /* Override default settings (unless already changed) */
3866552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 2000;
3876552cf8aSJason Sarich   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
388a7e14dcfSSatish Balay 
389a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
39063b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
391a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
392441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
3935d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
394cd929ea3SAlp Dener 
395cd929ea3SAlp Dener   ierr = KSPInitializePackage();CHKERRQ(ierr);
396cd929ea3SAlp Dener   ierr = MatCreate(((PetscObject)tao)->comm, &blmP->M);CHKERRQ(ierr);
397cd929ea3SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1);CHKERRQ(ierr);
39878e4361aSAlp Dener   ierr = MatSetType(blmP->M, MATLMVMBFGS);CHKERRQ(ierr);
399cd929ea3SAlp Dener   ierr = MatSetOptionsPrefix(blmP->M, "tao_blmvm_");CHKERRQ(ierr);
400a7e14dcfSSatish Balay   PetscFunctionReturn(0);
401a7e14dcfSSatish Balay }
402a7e14dcfSSatish Balay 
403cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode TaoBLMVMSetH0(Tao tao, Mat H0)
404a9603a14SPatrick Farrell {
405a9603a14SPatrick Farrell   TAO_LMVM       *lmP;
406a9603a14SPatrick Farrell   TAO_BLMVM      *blmP;
407b625d6c7SJed Brown   TaoType        type;
408a9603a14SPatrick Farrell   PetscBool      is_lmvm, is_blmvm;
409a9603a14SPatrick Farrell   PetscErrorCode ierr;
410a9603a14SPatrick Farrell 
411a9603a14SPatrick Farrell   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
412a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
413a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
414a9603a14SPatrick Farrell 
415a9603a14SPatrick Farrell   if (is_lmvm) {
416a9603a14SPatrick Farrell     lmP = (TAO_LMVM *)tao->data;
417a9603a14SPatrick Farrell     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
418a9603a14SPatrick Farrell     lmP->H0 = H0;
419a9603a14SPatrick Farrell   } else if (is_blmvm) {
420a9603a14SPatrick Farrell     blmP = (TAO_BLMVM *)tao->data;
421a9603a14SPatrick Farrell     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
422a9603a14SPatrick Farrell     blmP->H0 = H0;
42374c66251SPatrick Farrell   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
424a9603a14SPatrick Farrell   PetscFunctionReturn(0);
425a9603a14SPatrick Farrell }
426a9603a14SPatrick Farrell 
427cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode TaoBLMVMGetH0(Tao tao, Mat *H0)
428a9603a14SPatrick Farrell {
429a9603a14SPatrick Farrell   TAO_LMVM       *lmP;
430a9603a14SPatrick Farrell   TAO_BLMVM      *blmP;
431b625d6c7SJed Brown   TaoType        type;
432a9603a14SPatrick Farrell   PetscBool      is_lmvm, is_blmvm;
433a9603a14SPatrick Farrell   Mat            M;
434a9603a14SPatrick Farrell 
435a9603a14SPatrick Farrell   PetscErrorCode ierr;
436a9603a14SPatrick Farrell 
437a9603a14SPatrick Farrell   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
438a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
439a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
440a9603a14SPatrick Farrell 
441a9603a14SPatrick Farrell   if (is_lmvm) {
442a9603a14SPatrick Farrell     lmP = (TAO_LMVM *)tao->data;
443a9603a14SPatrick Farrell     M = lmP->M;
444a9603a14SPatrick Farrell   } else if (is_blmvm) {
445a9603a14SPatrick Farrell     blmP = (TAO_BLMVM *)tao->data;
446a9603a14SPatrick Farrell     M = blmP->M;
44774c66251SPatrick Farrell   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
448cd929ea3SAlp Dener   ierr = MatLMVMGetJ0(M, H0);CHKERRQ(ierr);
449a9603a14SPatrick Farrell   PetscFunctionReturn(0);
450a9603a14SPatrick Farrell }
451a9603a14SPatrick Farrell 
452cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode TaoBLMVMGetH0KSP(Tao tao, KSP *ksp)
453a9603a14SPatrick Farrell {
454a9603a14SPatrick Farrell   TAO_LMVM       *lmP;
455a9603a14SPatrick Farrell   TAO_BLMVM      *blmP;
456b625d6c7SJed Brown   TaoType        type;
457a9603a14SPatrick Farrell   PetscBool      is_lmvm, is_blmvm;
458a9603a14SPatrick Farrell   Mat            M;
459a9603a14SPatrick Farrell   PetscErrorCode ierr;
460a9603a14SPatrick Farrell 
461a9603a14SPatrick Farrell   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
462a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
463a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
464a9603a14SPatrick Farrell 
465a9603a14SPatrick Farrell   if (is_lmvm) {
466a9603a14SPatrick Farrell     lmP = (TAO_LMVM *)tao->data;
467a9603a14SPatrick Farrell     M = lmP->M;
468a9603a14SPatrick Farrell   } else if (is_blmvm) {
469a9603a14SPatrick Farrell     blmP = (TAO_BLMVM *)tao->data;
470a9603a14SPatrick Farrell     M = blmP->M;
47174c66251SPatrick Farrell   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
472cd929ea3SAlp Dener   ierr = MatLMVMGetJ0KSP(M, ksp);CHKERRQ(ierr);
473a9603a14SPatrick Farrell   PetscFunctionReturn(0);
474a9603a14SPatrick Farrell }
475