xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision e4cb33bb7dbdbae9285fba102465ca0f1dcb3977)
1ba92ff59SBarry Smith #include <petsctaolinesearch.h>
2aaa7dc30SBarry Smith #include <../src/tao/matrix/lmvmmat.h>
3aaa7dc30SBarry Smith #include <../src/tao/bound/impls/blmvm/blmvm.h>
4a7e14dcfSSatish Balay 
5a7e14dcfSSatish Balay /*------------------------------------------------------------*/
6a7e14dcfSSatish Balay #undef __FUNCT__
7a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_BLMVM"
8441846f8SBarry Smith static PetscErrorCode TaoSolve_BLMVM(Tao tao)
9a7e14dcfSSatish Balay {
10a7e14dcfSSatish Balay   PetscErrorCode               ierr;
11a7e14dcfSSatish Balay   TAO_BLMVM                    *blmP = (TAO_BLMVM *)tao->data;
12*e4cb33bbSBarry Smith   TaoConvergedReason           reason = TAO_CONTINUE_ITERATING;
13*e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
14a7e14dcfSSatish Balay   PetscReal                    f, fold, gdx, gnorm;
15a7e14dcfSSatish Balay   PetscReal                    stepsize = 1.0,delta;
16a7e14dcfSSatish Balay   PetscInt                     iter = 0;
17a7e14dcfSSatish Balay 
18a7e14dcfSSatish Balay   PetscFunctionBegin;
19a7e14dcfSSatish Balay   /*  Project initial point onto bounds */
20a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
21a7e14dcfSSatish Balay   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
22a7e14dcfSSatish Balay   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
23a7e14dcfSSatish Balay 
24a7e14dcfSSatish Balay   /* Check convergence criteria */
25a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);CHKERRQ(ierr);
26a7e14dcfSSatish Balay   ierr = VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
27a7e14dcfSSatish Balay 
28a7e14dcfSSatish Balay   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
2953506e15SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf pr NaN");
30a7e14dcfSSatish Balay 
31a7e14dcfSSatish Balay   ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, stepsize, &reason);CHKERRQ(ierr);
3253506e15SBarry Smith   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
33a7e14dcfSSatish Balay 
34a7e14dcfSSatish Balay   /* Set initial scaling for the function */
35a7e14dcfSSatish Balay   if (f != 0.0) {
36a7e14dcfSSatish Balay     delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
3753506e15SBarry Smith   } else {
38a7e14dcfSSatish Balay     delta = 2.0 / (gnorm*gnorm);
39a7e14dcfSSatish Balay   }
40a7e14dcfSSatish Balay   ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
41a7e14dcfSSatish Balay 
42a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps */
43a7e14dcfSSatish Balay   blmP->grad = 0;
44a7e14dcfSSatish Balay   blmP->reset = 0;
45a7e14dcfSSatish Balay 
46a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
47a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
48a7e14dcfSSatish Balay     /* Compute direction */
49a7e14dcfSSatish Balay     ierr = MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
50a7e14dcfSSatish Balay     ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
51a7e14dcfSSatish Balay     ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
52a7e14dcfSSatish Balay 
53a7e14dcfSSatish Balay     /* Check for success (descent direction) */
54a7e14dcfSSatish Balay     ierr = VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);CHKERRQ(ierr);
55a7e14dcfSSatish Balay     if (gdx <= 0) {
56a7e14dcfSSatish Balay       /* Step is not descent or solve was not successful
57a7e14dcfSSatish Balay          Use steepest descent direction (scaled) */
58a7e14dcfSSatish Balay       ++blmP->grad;
59a7e14dcfSSatish Balay 
60a7e14dcfSSatish Balay       if (f != 0.0) {
61a7e14dcfSSatish Balay         delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
6253506e15SBarry Smith       } else {
63a7e14dcfSSatish Balay         delta = 2.0 / (gnorm*gnorm);
64a7e14dcfSSatish Balay       }
65a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
66a7e14dcfSSatish Balay       ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr);
67a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
68a7e14dcfSSatish Balay       ierr = MatLMVMSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
69a7e14dcfSSatish Balay     }
70a7e14dcfSSatish Balay     ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
71a7e14dcfSSatish Balay 
72a7e14dcfSSatish Balay     /* Perform the linesearch */
73a7e14dcfSSatish Balay     fold = f;
74a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, blmP->Xold);CHKERRQ(ierr);
75a7e14dcfSSatish Balay     ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold);CHKERRQ(ierr);
76a7e14dcfSSatish Balay     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
77a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr);
78a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
79a7e14dcfSSatish Balay 
80a7e14dcfSSatish Balay     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
81a7e14dcfSSatish Balay       /* Linesearch failed
82a7e14dcfSSatish Balay          Reset factors and use scaled (projected) gradient step */
83a7e14dcfSSatish Balay       ++blmP->reset;
84a7e14dcfSSatish Balay 
85a7e14dcfSSatish Balay       f = fold;
86a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
87a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
88a7e14dcfSSatish Balay 
89a7e14dcfSSatish Balay       if (f != 0.0) {
90a7e14dcfSSatish Balay         delta = 2.0* PetscAbsScalar(f) / (gnorm*gnorm);
9153506e15SBarry Smith       } else {
92a7e14dcfSSatish Balay         delta = 2.0/ (gnorm*gnorm);
93a7e14dcfSSatish Balay       }
94a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
95a7e14dcfSSatish Balay       ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr);
96a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
97a7e14dcfSSatish Balay       ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
98a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
99a7e14dcfSSatish Balay 
100a7e14dcfSSatish Balay       /* This may be incorrect; linesearch has values fo stepmax and stepmin
101a7e14dcfSSatish Balay          that should be reset. */
102a7e14dcfSSatish Balay       ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
103a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status);CHKERRQ(ierr);
104a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
105a7e14dcfSSatish Balay 
106a7e14dcfSSatish Balay       if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
107a7e14dcfSSatish Balay         tao->reason = TAO_DIVERGED_LS_FAILURE;
108a7e14dcfSSatish Balay         break;
109a7e14dcfSSatish Balay       }
110a7e14dcfSSatish Balay     }
111a7e14dcfSSatish Balay 
112*e4cb33bbSBarry Smith     /* Check for converged */
113a7e14dcfSSatish Balay     ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr);
114a7e14dcfSSatish Balay     ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
115a7e14dcfSSatish Balay 
116a7e14dcfSSatish Balay 
11753506e15SBarry Smith     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
118a7e14dcfSSatish Balay     iter++;
119a7e14dcfSSatish Balay     ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, stepsize, &reason);CHKERRQ(ierr);
120a7e14dcfSSatish Balay   }
121a7e14dcfSSatish Balay   PetscFunctionReturn(0);
122a7e14dcfSSatish Balay }
123a7e14dcfSSatish Balay 
124a7e14dcfSSatish Balay #undef __FUNCT__
125a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetup_BLMVM"
126441846f8SBarry Smith static PetscErrorCode TaoSetup_BLMVM(Tao tao)
127a7e14dcfSSatish Balay {
128a7e14dcfSSatish Balay   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
129a7e14dcfSSatish Balay   PetscInt       n,N;
130a7e14dcfSSatish Balay   PetscErrorCode ierr;
131a7e14dcfSSatish Balay 
132a7e14dcfSSatish Balay   PetscFunctionBegin;
133a7e14dcfSSatish Balay   /* Existence of tao->solution checked in TaoSetup() */
134a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr);
135a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr);
136a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);
137a7e14dcfSSatish Balay 
138a7e14dcfSSatish Balay   if (!tao->stepdirection) {
13953506e15SBarry Smith     ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
140a7e14dcfSSatish Balay   }
141a7e14dcfSSatish Balay   if (!tao->gradient) {
142a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
143a7e14dcfSSatish Balay   }
144a7e14dcfSSatish Balay   if (!tao->XL) {
145a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&tao->XL);CHKERRQ(ierr);
146e270355aSBarry Smith     ierr = VecSet(tao->XL,PETSC_NINFINITY);CHKERRQ(ierr);
147a7e14dcfSSatish Balay   }
148a7e14dcfSSatish Balay   if (!tao->XU) {
149a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&tao->XU);CHKERRQ(ierr);
150e270355aSBarry Smith     ierr = VecSet(tao->XU,PETSC_INFINITY);CHKERRQ(ierr);
151a7e14dcfSSatish Balay   }
152a7e14dcfSSatish Balay   /* Create matrix for the limited memory approximation */
153a7e14dcfSSatish Balay   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
154a7e14dcfSSatish Balay   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
155a7e14dcfSSatish Balay   ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&blmP->M);CHKERRQ(ierr);
156a7e14dcfSSatish Balay   ierr = MatLMVMAllocateVectors(blmP->M,tao->solution);CHKERRQ(ierr);
157a7e14dcfSSatish Balay   PetscFunctionReturn(0);
158a7e14dcfSSatish Balay }
159a7e14dcfSSatish Balay 
160a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
161a7e14dcfSSatish Balay #undef __FUNCT__
162a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_BLMVM"
163441846f8SBarry Smith static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
164a7e14dcfSSatish Balay {
165a7e14dcfSSatish Balay   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
166a7e14dcfSSatish Balay   PetscErrorCode ierr;
167a7e14dcfSSatish Balay 
168a7e14dcfSSatish Balay   PetscFunctionBegin;
169a7e14dcfSSatish Balay   if (tao->setupcalled) {
170a7e14dcfSSatish Balay     ierr = MatDestroy(&blmP->M);CHKERRQ(ierr);
171a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr);
172a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr);
173a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr);
174a7e14dcfSSatish Balay   }
175a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
176a7e14dcfSSatish Balay   PetscFunctionReturn(0);
177a7e14dcfSSatish Balay }
178a7e14dcfSSatish Balay 
179a7e14dcfSSatish Balay /*------------------------------------------------------------*/
180a7e14dcfSSatish Balay #undef __FUNCT__
181a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_BLMVM"
182441846f8SBarry Smith static PetscErrorCode TaoSetFromOptions_BLMVM(Tao tao)
183a7e14dcfSSatish Balay {
184a7e14dcfSSatish Balay   PetscErrorCode ierr;
185a7e14dcfSSatish Balay 
186a7e14dcfSSatish Balay   PetscFunctionBegin;
187a7e14dcfSSatish Balay   ierr = PetscOptionsHead("Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr);
188a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
189a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
190a7e14dcfSSatish Balay   PetscFunctionReturn(0);
191a7e14dcfSSatish Balay }
192a7e14dcfSSatish Balay 
193a7e14dcfSSatish Balay 
194a7e14dcfSSatish Balay /*------------------------------------------------------------*/
195a7e14dcfSSatish Balay #undef __FUNCT__
196a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_BLMVM"
197441846f8SBarry Smith static int TaoView_BLMVM(Tao tao, PetscViewer viewer)
198a7e14dcfSSatish Balay {
199a7e14dcfSSatish Balay   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
200a7e14dcfSSatish Balay   PetscBool      isascii;
201a7e14dcfSSatish Balay   PetscErrorCode ierr;
202a7e14dcfSSatish Balay 
203a7e14dcfSSatish Balay   PetscFunctionBegin;
204a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
205a7e14dcfSSatish Balay   if (isascii) {
206a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
207a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr);
208a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
209a7e14dcfSSatish Balay   }
210a7e14dcfSSatish Balay   PetscFunctionReturn(0);
211a7e14dcfSSatish Balay }
212a7e14dcfSSatish Balay 
213a7e14dcfSSatish Balay #undef __FUNCT__
214a7e14dcfSSatish Balay #define __FUNCT__ "TaoComputeDual_BLMVM"
215441846f8SBarry Smith static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
216a7e14dcfSSatish Balay {
217a7e14dcfSSatish Balay   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
218a7e14dcfSSatish Balay   PetscErrorCode ierr;
219a7e14dcfSSatish Balay 
220a7e14dcfSSatish Balay   PetscFunctionBegin;
221441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
222a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
223a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
22453506e15SBarry 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");
225a7e14dcfSSatish Balay 
226a7e14dcfSSatish Balay   ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr);
227a7e14dcfSSatish Balay   ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr);
228a7e14dcfSSatish Balay   ierr = VecSet(DXU,0.0);CHKERRQ(ierr);
229a7e14dcfSSatish Balay   ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr);
230a7e14dcfSSatish Balay 
231a7e14dcfSSatish Balay   ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr);
232a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr);
233a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr);
234a7e14dcfSSatish Balay   PetscFunctionReturn(0);
235a7e14dcfSSatish Balay }
236a7e14dcfSSatish Balay 
237a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
238a7e14dcfSSatish Balay EXTERN_C_BEGIN
239a7e14dcfSSatish Balay #undef __FUNCT__
240a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_BLMVM"
241441846f8SBarry Smith PetscErrorCode TaoCreate_BLMVM(Tao tao)
242a7e14dcfSSatish Balay {
243a7e14dcfSSatish Balay   TAO_BLMVM      *blmP;
244a7e14dcfSSatish Balay   const char     *morethuente_type = TAOLINESEARCH_MT;
245a7e14dcfSSatish Balay   PetscErrorCode ierr;
246a7e14dcfSSatish Balay 
247a7e14dcfSSatish Balay   PetscFunctionBegin;
248a7e14dcfSSatish Balay   tao->ops->setup = TaoSetup_BLMVM;
249a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_BLMVM;
250a7e14dcfSSatish Balay   tao->ops->view = TaoView_BLMVM;
251a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
252a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_BLMVM;
253a7e14dcfSSatish Balay   tao->ops->computedual = TaoComputeDual_BLMVM;
254a7e14dcfSSatish Balay 
2553c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr);
256a7e14dcfSSatish Balay   tao->data = (void*)blmP;
257a7e14dcfSSatish Balay   tao->max_it = 2000;
258a7e14dcfSSatish Balay   tao->max_funcs = 4000;
259a7e14dcfSSatish Balay   tao->fatol = 1e-4;
260a7e14dcfSSatish Balay   tao->frtol = 1e-4;
261a7e14dcfSSatish Balay 
262a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
263a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
264441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
265a7e14dcfSSatish Balay   PetscFunctionReturn(0);
266a7e14dcfSSatish Balay }
267a7e14dcfSSatish Balay EXTERN_C_END
268a7e14dcfSSatish Balay 
269