xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision 2ec5c1acdc8744b9aa587f7e7afc34821d593ce8)
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 
5a7e14dcfSSatish Balay /*------------------------------------------------------------*/
6*2ec5c1acSAlp Dener static PetscErrorCode TaoSolve_BLMVM(Tao tao)
7a7e14dcfSSatish Balay {
8a7e14dcfSSatish Balay   PetscErrorCode               ierr;
9a7e14dcfSSatish Balay   TAO_BLMVM                    *blmP = (TAO_BLMVM *)tao->data;
10e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
11*2ec5c1acSAlp Dener   PetscReal                    f, fold, gdx, gnorm;
12*2ec5c1acSAlp Dener   PetscReal                    stepsize = 1.0,delta;
13*2ec5c1acSAlp Dener   PetscInt                     stepType = BLMVM_STEP_GRAD;
14a7e14dcfSSatish Balay 
15a7e14dcfSSatish Balay   PetscFunctionBegin;
16a7e14dcfSSatish Balay   /*  Project initial point onto bounds */
17a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
18*2ec5c1acSAlp Dener   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
19cd929ea3SAlp Dener   if (tao->bounded) {
20a7e14dcfSSatish Balay     ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
21cd929ea3SAlp Dener   }
22a9603a14SPatrick Farrell 
23a7e14dcfSSatish Balay   /* Check convergence criteria */
24a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);CHKERRQ(ierr);
25a7e14dcfSSatish Balay   ierr = VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
26a7e14dcfSSatish Balay 
27*2ec5c1acSAlp Dener   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
283fcc4fd3STodd Munson   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
29a7e14dcfSSatish Balay 
303ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
31*2ec5c1acSAlp Dener   ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
32*2ec5c1acSAlp Dener   ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize);CHKERRQ(ierr);
333ecd9318SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
343ecd9318SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
35a7e14dcfSSatish Balay 
36*2ec5c1acSAlp Dener   /* Set initial scaling for the function */
37*2ec5c1acSAlp Dener   if (blmP->Mscale) {
38*2ec5c1acSAlp Dener     if (f != 0.0) {
39*2ec5c1acSAlp Dener       delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
40*2ec5c1acSAlp Dener     } else {
41*2ec5c1acSAlp Dener       delta = 2.0 / (gnorm*gnorm);
42cd929ea3SAlp Dener     }
43*2ec5c1acSAlp Dener     ierr = MatLMVMSetJ0Scale(blmP->Mscale, 1.0/delta);CHKERRQ(ierr);
44*2ec5c1acSAlp Dener   }
45*2ec5c1acSAlp Dener   ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
46*2ec5c1acSAlp Dener 
47*2ec5c1acSAlp Dener   /* Set counter for gradient/reset steps */
48*2ec5c1acSAlp Dener   blmP->grad = 0;
49*2ec5c1acSAlp Dener   blmP->bfgs = 0;
50a7e14dcfSSatish Balay 
51a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
523ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
53a7e14dcfSSatish Balay     /* Compute direction */
54*2ec5c1acSAlp Dener     ierr = MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
55*2ec5c1acSAlp Dener     ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
56*2ec5c1acSAlp Dener     ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
57a7e14dcfSSatish Balay 
58a7e14dcfSSatish Balay     /* Check for success (descent direction) */
59*2ec5c1acSAlp Dener     ierr = VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);CHKERRQ(ierr);
60*2ec5c1acSAlp Dener     if (gdx <= 0) {
61a7e14dcfSSatish Balay       /* Step is not descent or solve was not successful
62a7e14dcfSSatish Balay          Use steepest descent direction (scaled) */
63*2ec5c1acSAlp Dener       stepType = BLMVM_STEP_GRAD;
64*2ec5c1acSAlp Dener 
65*2ec5c1acSAlp Dener       if (blmP->Mscale) {
66*2ec5c1acSAlp Dener         if (f != 0.0) {
67*2ec5c1acSAlp Dener           delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
68*2ec5c1acSAlp Dener         } else {
69*2ec5c1acSAlp Dener           delta = 2.0 / (gnorm*gnorm);
70*2ec5c1acSAlp Dener         }
71*2ec5c1acSAlp Dener         ierr = MatLMVMSetJ0Scale(blmP->Mscale, 1.0/delta);CHKERRQ(ierr);
72*2ec5c1acSAlp Dener       }
73cd929ea3SAlp Dener       ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
74a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
75cd929ea3SAlp Dener       ierr = MatSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
76b2d8c577SAlp Dener     }
77*2ec5c1acSAlp Dener     ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
78a7e14dcfSSatish Balay 
79a7e14dcfSSatish Balay     /* Perform the linesearch */
80a7e14dcfSSatish Balay     fold = f;
81a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, blmP->Xold);CHKERRQ(ierr);
82a7e14dcfSSatish Balay     ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold);CHKERRQ(ierr);
83*2ec5c1acSAlp Dener     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
84a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr);
85a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
86a7e14dcfSSatish Balay 
87cd929ea3SAlp Dener     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && stepType != BLMVM_STEP_GRAD) {
88a7e14dcfSSatish Balay       /* Linesearch failed
89a7e14dcfSSatish Balay          Reset factors and use scaled (projected) gradient step */
90a7e14dcfSSatish Balay       f = fold;
91a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
92a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
93a7e14dcfSSatish Balay 
94*2ec5c1acSAlp Dener       if (blmP->Mscale) {
95*2ec5c1acSAlp Dener         if (f != 0.0) {
96*2ec5c1acSAlp Dener           delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
97*2ec5c1acSAlp Dener         } else {
98*2ec5c1acSAlp Dener           delta = 2.0 / (gnorm*gnorm);
99*2ec5c1acSAlp Dener         }
100*2ec5c1acSAlp Dener         ierr = MatLMVMSetJ0Scale(blmP->Mscale, 1.0/delta);CHKERRQ(ierr);
101*2ec5c1acSAlp Dener       }
102cd929ea3SAlp Dener       ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
103a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
104cd929ea3SAlp Dener       ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
105a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
106a7e14dcfSSatish Balay 
107*2ec5c1acSAlp Dener       /* This may be incorrect; linesearch has values for stepmax and stepmin
108*2ec5c1acSAlp Dener          that should be reset. */
109*2ec5c1acSAlp Dener       ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
110a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status);CHKERRQ(ierr);
111a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
112cd929ea3SAlp Dener     }
113a7e14dcfSSatish Balay 
114a7e14dcfSSatish Balay     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
115a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
116cd929ea3SAlp Dener       break;
117a7e14dcfSSatish Balay     }
118a7e14dcfSSatish Balay 
119e4cb33bbSBarry Smith     /* Check for converged */
120*2ec5c1acSAlp Dener     ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr);
121*2ec5c1acSAlp Dener     ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
122*2ec5c1acSAlp Dener 
123*2ec5c1acSAlp Dener 
124*2ec5c1acSAlp Dener     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
1258931d482SJason Sarich     tao->niter++;
126*2ec5c1acSAlp Dener     ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
127*2ec5c1acSAlp Dener     ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize);CHKERRQ(ierr);
1283ecd9318SAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
129a7e14dcfSSatish Balay   }
130a7e14dcfSSatish Balay   PetscFunctionReturn(0);
131a7e14dcfSSatish Balay }
132a7e14dcfSSatish Balay 
133*2ec5c1acSAlp Dener static PetscErrorCode TaoSetup_BLMVM(Tao tao)
134a7e14dcfSSatish Balay {
135a7e14dcfSSatish Balay   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
136*2ec5c1acSAlp Dener   const char     *prefix;
137a7e14dcfSSatish Balay   PetscErrorCode ierr;
138a7e14dcfSSatish Balay 
139a7e14dcfSSatish Balay   PetscFunctionBegin;
140*2ec5c1acSAlp Dener   if (!blmP->Xold) {
141a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr);
142*2ec5c1acSAlp Dener   }
143*2ec5c1acSAlp Dener   if (!blmP->Gold) {
144a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr);
145*2ec5c1acSAlp Dener   }
146*2ec5c1acSAlp Dener   if (!blmP->unprojected_gradient) {
147302440fdSBarry Smith     ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr);
148*2ec5c1acSAlp Dener   }
149a7e14dcfSSatish Balay   if (!tao->stepdirection) {
15053506e15SBarry Smith     ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
151a7e14dcfSSatish Balay   }
152a7e14dcfSSatish Balay   if (!tao->gradient) {
153a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
154a7e14dcfSSatish Balay   }
155a7e14dcfSSatish Balay   /* Create matrix for the limited memory approximation */
156cd929ea3SAlp Dener   ierr = MatLMVMAllocate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
157a9603a14SPatrick Farrell 
158a9603a14SPatrick Farrell   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
159a9603a14SPatrick Farrell   if (blmP->H0) {
160cd929ea3SAlp Dener     ierr = MatLMVMSetJ0(blmP->M, blmP->H0);CHKERRQ(ierr);
1610ad3a497SAlp Dener   } else if (!blmP->no_scale) {
162*2ec5c1acSAlp Dener     ierr = MatCreate(((PetscObject)tao)->comm, &blmP->Mscale);CHKERRQ(ierr);
1635aff1b4eSAlp Dener     ierr = MatSetType(blmP->Mscale, MATLMVMDIAGBRDN);CHKERRQ(ierr);
164*2ec5c1acSAlp Dener     ierr = TaoGetOptionsPrefix(tao, &prefix);CHKERRQ(ierr);
165*2ec5c1acSAlp Dener     ierr = MatSetOptionsPrefix(blmP->Mscale, prefix);CHKERRQ(ierr);
166*2ec5c1acSAlp Dener     ierr = MatAppendOptionsPrefix(blmP->Mscale, "tao_blmvm_scale_");CHKERRQ(ierr);
1675aff1b4eSAlp Dener     ierr = MatSetFromOptions(blmP->Mscale);CHKERRQ(ierr);
1680ad3a497SAlp Dener     ierr = MatLMVMAllocate(blmP->Mscale, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
1690ad3a497SAlp Dener     ierr = MatLMVMSetJ0(blmP->M, blmP->Mscale);CHKERRQ(ierr);
170a9603a14SPatrick Farrell   }
171a7e14dcfSSatish Balay   PetscFunctionReturn(0);
172a7e14dcfSSatish Balay }
173a7e14dcfSSatish Balay 
174a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
175*2ec5c1acSAlp Dener static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
176a7e14dcfSSatish Balay {
177a7e14dcfSSatish Balay   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
178a7e14dcfSSatish Balay   PetscErrorCode ierr;
179a7e14dcfSSatish Balay 
180a7e14dcfSSatish Balay   PetscFunctionBegin;
181a7e14dcfSSatish Balay   if (tao->setupcalled) {
182b2d8c577SAlp Dener     ierr = MatDestroy(&blmP->M);CHKERRQ(ierr);
1830ad3a497SAlp Dener     if (blmP->Mscale) {
1840ad3a497SAlp Dener       ierr = MatDestroy(&blmP->Mscale);CHKERRQ(ierr);
1850ad3a497SAlp Dener     }
186*2ec5c1acSAlp Dener     ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr);
187*2ec5c1acSAlp Dener     ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr);
188*2ec5c1acSAlp Dener     ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr);
189*2ec5c1acSAlp Dener   }
190*2ec5c1acSAlp Dener 
191*2ec5c1acSAlp Dener   if (blmP->H0) {
192*2ec5c1acSAlp Dener     PetscObjectDereference((PetscObject)blmP->H0);
193*2ec5c1acSAlp Dener   }
194*2ec5c1acSAlp Dener 
195a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
196a7e14dcfSSatish Balay   PetscFunctionReturn(0);
197a7e14dcfSSatish Balay }
198a7e14dcfSSatish Balay 
199a7e14dcfSSatish Balay /*------------------------------------------------------------*/
200*2ec5c1acSAlp Dener static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
201a7e14dcfSSatish Balay {
202cd929ea3SAlp Dener   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
203a7e14dcfSSatish Balay   PetscErrorCode ierr;
204a7e14dcfSSatish Balay 
205a7e14dcfSSatish Balay   PetscFunctionBegin;
2061a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr);
207*2ec5c1acSAlp Dener   ierr = PetscOptionsBool("-tao_blmvm_no_scale","(developer) disable the diagonal Broyden scaling of the QN approximation","",blmP->no_scale,&blmP->no_scale,NULL);CHKERRQ(ierr);
2087b1c7716SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
209a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
210b2d8c577SAlp Dener   ierr = MatSetFromOptions(blmP->M);CHKERRQ(ierr);
211a7e14dcfSSatish Balay   PetscFunctionReturn(0);
212a7e14dcfSSatish Balay }
213a7e14dcfSSatish Balay 
214*2ec5c1acSAlp Dener 
215a7e14dcfSSatish Balay /*------------------------------------------------------------*/
216*2ec5c1acSAlp Dener static int TaoView_BLMVM(Tao tao, PetscViewer viewer)
217a7e14dcfSSatish Balay {
218a7e14dcfSSatish Balay   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
219a7e14dcfSSatish Balay   PetscBool      isascii;
220a7e14dcfSSatish Balay   PetscErrorCode ierr;
221a7e14dcfSSatish Balay 
222a7e14dcfSSatish Balay   PetscFunctionBegin;
223a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
224a7e14dcfSSatish Balay   if (isascii) {
225a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr);
226a7e14dcfSSatish Balay   }
227a7e14dcfSSatish Balay   PetscFunctionReturn(0);
228a7e14dcfSSatish Balay }
229a7e14dcfSSatish Balay 
230*2ec5c1acSAlp Dener static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
231a7e14dcfSSatish Balay {
232a7e14dcfSSatish Balay   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
233a7e14dcfSSatish Balay   PetscErrorCode ierr;
234a7e14dcfSSatish Balay 
235a7e14dcfSSatish Balay   PetscFunctionBegin;
236441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
237a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
238a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
23953506e15SBarry 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");
240a7e14dcfSSatish Balay 
241a7e14dcfSSatish Balay   ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr);
242a7e14dcfSSatish Balay   ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr);
243a7e14dcfSSatish Balay   ierr = VecSet(DXU,0.0);CHKERRQ(ierr);
244a7e14dcfSSatish Balay   ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr);
245a7e14dcfSSatish Balay 
246a7e14dcfSSatish Balay   ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr);
247a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr);
248a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr);
249a7e14dcfSSatish Balay   PetscFunctionReturn(0);
250a7e14dcfSSatish Balay }
251a7e14dcfSSatish Balay 
252a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
2531522df2eSJason Sarich /*MC
2541522df2eSJason Sarich   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
2551522df2eSJason Sarich          for nonlinear minimization with bound constraints. It is an extension
2561522df2eSJason Sarich          of TAOLMVM
2571522df2eSJason Sarich 
2580ad3a497SAlp Dener   Options Database Keys:
259*2ec5c1acSAlp Dener +     -tao_lmm_vectors - number of vectors to use for approximation
260*2ec5c1acSAlp Dener .     -tao_lmm_scale_type - "none","scalar","broyden"
261*2ec5c1acSAlp Dener .     -tao_lmm_limit_type - "none","average","relative","absolute"
262*2ec5c1acSAlp Dener .     -tao_lmm_rescale_type - "none","scalar","gl"
263*2ec5c1acSAlp Dener .     -tao_lmm_limit_mu - mu limiting factor
264*2ec5c1acSAlp Dener .     -tao_lmm_limit_nu - nu limiting factor
265*2ec5c1acSAlp Dener .     -tao_lmm_delta_min - minimum delta value
266*2ec5c1acSAlp Dener .     -tao_lmm_delta_max - maximum delta value
267*2ec5c1acSAlp Dener .     -tao_lmm_broyden_phi - phi factor for Broyden scaling
268*2ec5c1acSAlp Dener .     -tao_lmm_scalar_alpha - alpha factor for scalar scaling
269*2ec5c1acSAlp Dener .     -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal
270*2ec5c1acSAlp Dener .     -tao_lmm_rescale_beta - beta factor for rescaling diagonal
271*2ec5c1acSAlp Dener .     -tao_lmm_scalar_history - amount of history for scalar scaling
272*2ec5c1acSAlp Dener .     -tao_lmm_rescale_history - amount of history for rescaling diagonal
273*2ec5c1acSAlp Dener -     -tao_lmm_eps - rejection tolerance
2740ad3a497SAlp Dener 
2751eb8069cSJason Sarich   Level: beginner
2761522df2eSJason Sarich M*/
277728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
278a7e14dcfSSatish Balay {
279a7e14dcfSSatish Balay   TAO_BLMVM      *blmP;
280*2ec5c1acSAlp Dener   const char     *prefix;
2818caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
282a7e14dcfSSatish Balay   PetscErrorCode ierr;
283a7e14dcfSSatish Balay 
284a7e14dcfSSatish Balay   PetscFunctionBegin;
285a7e14dcfSSatish Balay   tao->ops->setup = TaoSetup_BLMVM;
286a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_BLMVM;
287a7e14dcfSSatish Balay   tao->ops->view = TaoView_BLMVM;
288a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
289a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_BLMVM;
290a7e14dcfSSatish Balay   tao->ops->computedual = TaoComputeDual_BLMVM;
291a7e14dcfSSatish Balay 
2923c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr);
293a9603a14SPatrick Farrell   blmP->H0 = NULL;
2940ad3a497SAlp Dener   blmP->no_scale = PETSC_FALSE;
295a7e14dcfSSatish Balay   tao->data = (void*)blmP;
2966552cf8aSJason Sarich 
2976552cf8aSJason Sarich   /* Override default settings (unless already changed) */
2986552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 2000;
2996552cf8aSJason Sarich   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
300a7e14dcfSSatish Balay 
301a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
30263b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
303a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
304441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
3055d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
306cd929ea3SAlp Dener 
307cd929ea3SAlp Dener   ierr = MatCreate(((PetscObject)tao)->comm, &blmP->M);CHKERRQ(ierr);
30878e4361aSAlp Dener   ierr = MatSetType(blmP->M, MATLMVMBFGS);CHKERRQ(ierr);
309*2ec5c1acSAlp Dener   ierr = TaoGetOptionsPrefix(tao, &prefix);CHKERRQ(ierr);
310*2ec5c1acSAlp Dener   ierr = MatSetOptionsPrefix(blmP->M, prefix);CHKERRQ(ierr);
311*2ec5c1acSAlp Dener   ierr = MatAppendOptionsPrefix(blmP->M, "tao_blmvm_");CHKERRQ(ierr);
312a7e14dcfSSatish Balay   PetscFunctionReturn(0);
313a7e14dcfSSatish Balay }
314a7e14dcfSSatish Balay 
315cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode TaoBLMVMSetH0(Tao tao, Mat H0)
316a9603a14SPatrick Farrell {
317a9603a14SPatrick Farrell   TAO_LMVM       *lmP;
318a9603a14SPatrick Farrell   TAO_BLMVM      *blmP;
319b625d6c7SJed Brown   TaoType        type;
320a9603a14SPatrick Farrell   PetscBool      is_lmvm, is_blmvm;
321a9603a14SPatrick Farrell   PetscErrorCode ierr;
322a9603a14SPatrick Farrell 
323a9603a14SPatrick Farrell   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
324a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
325a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
326a9603a14SPatrick Farrell 
327a9603a14SPatrick Farrell   if (is_lmvm) {
328a9603a14SPatrick Farrell     lmP = (TAO_LMVM *)tao->data;
329a9603a14SPatrick Farrell     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
330a9603a14SPatrick Farrell     lmP->H0 = H0;
331a9603a14SPatrick Farrell   } else if (is_blmvm) {
332a9603a14SPatrick Farrell     blmP = (TAO_BLMVM *)tao->data;
333a9603a14SPatrick Farrell     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
334a9603a14SPatrick Farrell     blmP->H0 = H0;
33574c66251SPatrick Farrell   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
336a9603a14SPatrick Farrell   PetscFunctionReturn(0);
337a9603a14SPatrick Farrell }
338a9603a14SPatrick Farrell 
339cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode TaoBLMVMGetH0(Tao tao, Mat *H0)
340a9603a14SPatrick Farrell {
341a9603a14SPatrick Farrell   TAO_LMVM       *lmP;
342a9603a14SPatrick Farrell   TAO_BLMVM      *blmP;
343b625d6c7SJed Brown   TaoType        type;
344a9603a14SPatrick Farrell   PetscBool      is_lmvm, is_blmvm;
345a9603a14SPatrick Farrell   Mat            M;
346a9603a14SPatrick Farrell 
347a9603a14SPatrick Farrell   PetscErrorCode ierr;
348a9603a14SPatrick Farrell 
349a9603a14SPatrick Farrell   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
350a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
351a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
352a9603a14SPatrick Farrell 
353a9603a14SPatrick Farrell   if (is_lmvm) {
354a9603a14SPatrick Farrell     lmP = (TAO_LMVM *)tao->data;
355a9603a14SPatrick Farrell     M = lmP->M;
356a9603a14SPatrick Farrell   } else if (is_blmvm) {
357a9603a14SPatrick Farrell     blmP = (TAO_BLMVM *)tao->data;
358a9603a14SPatrick Farrell     M = blmP->M;
35974c66251SPatrick Farrell   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
360cd929ea3SAlp Dener   ierr = MatLMVMGetJ0(M, H0);CHKERRQ(ierr);
361a9603a14SPatrick Farrell   PetscFunctionReturn(0);
362a9603a14SPatrick Farrell }
363a9603a14SPatrick Farrell 
364cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode TaoBLMVMGetH0KSP(Tao tao, KSP *ksp)
365a9603a14SPatrick Farrell {
366a9603a14SPatrick Farrell   TAO_LMVM       *lmP;
367a9603a14SPatrick Farrell   TAO_BLMVM      *blmP;
368b625d6c7SJed Brown   TaoType        type;
369a9603a14SPatrick Farrell   PetscBool      is_lmvm, is_blmvm;
370a9603a14SPatrick Farrell   Mat            M;
371a9603a14SPatrick Farrell   PetscErrorCode ierr;
372a9603a14SPatrick Farrell 
373a9603a14SPatrick Farrell   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
374a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
375a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
376a9603a14SPatrick Farrell 
377a9603a14SPatrick Farrell   if (is_lmvm) {
378a9603a14SPatrick Farrell     lmP = (TAO_LMVM *)tao->data;
379a9603a14SPatrick Farrell     M = lmP->M;
380a9603a14SPatrick Farrell   } else if (is_blmvm) {
381a9603a14SPatrick Farrell     blmP = (TAO_BLMVM *)tao->data;
382a9603a14SPatrick Farrell     M = blmP->M;
38374c66251SPatrick Farrell   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
384cd929ea3SAlp Dener   ierr = MatLMVMGetJ0KSP(M, ksp);CHKERRQ(ierr);
385a9603a14SPatrick Farrell   PetscFunctionReturn(0);
386a9603a14SPatrick Farrell }