xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
1 #include <petsctaolinesearch.h>
2 #include <../src/tao/matrix/lmvmmat.h>
3 #include <../src/tao/bound/impls/blmvm/blmvm.h>
4 
5 /*------------------------------------------------------------*/
6 #undef __FUNCT__
7 #define __FUNCT__ "TaoSolve_BLMVM"
8 static PetscErrorCode TaoSolve_BLMVM(Tao tao)
9 {
10   PetscErrorCode               ierr;
11   TAO_BLMVM                    *blmP = (TAO_BLMVM *)tao->data;
12   TaoConvergedReason           reason = TAO_CONTINUE_ITERATING;
13   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
14   PetscReal                    f, fold, gdx, gnorm;
15   PetscReal                    stepsize = 1.0,delta;
16 
17   PetscFunctionBegin;
18   /*  Project initial point onto bounds */
19   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
20   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
21   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
22 
23   /* Check convergence criteria */
24   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);CHKERRQ(ierr);
25   ierr = VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
26 
27   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
28   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf pr NaN");
29 
30   ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize, &reason);CHKERRQ(ierr);
31   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
32 
33   /* Set initial scaling for the function */
34   if (f != 0.0) {
35     delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
36   } else {
37     delta = 2.0 / (gnorm*gnorm);
38   }
39   ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
40 
41   /* Set counter for gradient/reset steps */
42   blmP->grad = 0;
43   blmP->reset = 0;
44 
45   /* Have not converged; continue with Newton method */
46   while (reason == TAO_CONTINUE_ITERATING) {
47     /* Compute direction */
48     ierr = MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
49     ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
50     ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
51 
52     /* Check for success (descent direction) */
53     ierr = VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);CHKERRQ(ierr);
54     if (gdx <= 0) {
55       /* Step is not descent or solve was not successful
56          Use steepest descent direction (scaled) */
57       ++blmP->grad;
58 
59       if (f != 0.0) {
60         delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
61       } else {
62         delta = 2.0 / (gnorm*gnorm);
63       }
64       ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
65       ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr);
66       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
67       ierr = MatLMVMSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
68     }
69     ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
70 
71     /* Perform the linesearch */
72     fold = f;
73     ierr = VecCopy(tao->solution, blmP->Xold);CHKERRQ(ierr);
74     ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold);CHKERRQ(ierr);
75     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
76     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr);
77     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
78 
79     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
80       /* Linesearch failed
81          Reset factors and use scaled (projected) gradient step */
82       ++blmP->reset;
83 
84       f = fold;
85       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
86       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
87 
88       if (f != 0.0) {
89         delta = 2.0* PetscAbsScalar(f) / (gnorm*gnorm);
90       } else {
91         delta = 2.0/ (gnorm*gnorm);
92       }
93       ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
94       ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr);
95       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
96       ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
97       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
98 
99       /* This may be incorrect; linesearch has values fo stepmax and stepmin
100          that should be reset. */
101       ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
102       ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status);CHKERRQ(ierr);
103       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
104 
105       if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
106         tao->reason = TAO_DIVERGED_LS_FAILURE;
107         break;
108       }
109     }
110 
111     /* Check for converged */
112     ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr);
113     ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
114 
115 
116     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
117     tao->niter++;
118     ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize, &reason);CHKERRQ(ierr);
119   }
120   PetscFunctionReturn(0);
121 }
122 
123 #undef __FUNCT__
124 #define __FUNCT__ "TaoSetup_BLMVM"
125 static PetscErrorCode TaoSetup_BLMVM(Tao tao)
126 {
127   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
128   PetscInt       n,N;
129   PetscErrorCode ierr;
130 
131   PetscFunctionBegin;
132   /* Existence of tao->solution checked in TaoSetup() */
133   ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr);
134   ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr);
135   ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr);
136 
137   if (!tao->stepdirection) {
138     ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
139   }
140   if (!tao->gradient) {
141     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
142   }
143   if (!tao->XL) {
144     ierr = VecDuplicate(tao->solution,&tao->XL);CHKERRQ(ierr);
145     ierr = VecSet(tao->XL,PETSC_NINFINITY);CHKERRQ(ierr);
146   }
147   if (!tao->XU) {
148     ierr = VecDuplicate(tao->solution,&tao->XU);CHKERRQ(ierr);
149     ierr = VecSet(tao->XU,PETSC_INFINITY);CHKERRQ(ierr);
150   }
151   /* Create matrix for the limited memory approximation */
152   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
153   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
154   ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&blmP->M);CHKERRQ(ierr);
155   ierr = MatLMVMAllocateVectors(blmP->M,tao->solution);CHKERRQ(ierr);
156   PetscFunctionReturn(0);
157 }
158 
159 /* ---------------------------------------------------------- */
160 #undef __FUNCT__
161 #define __FUNCT__ "TaoDestroy_BLMVM"
162 static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
163 {
164   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
165   PetscErrorCode ierr;
166 
167   PetscFunctionBegin;
168   if (tao->setupcalled) {
169     ierr = MatDestroy(&blmP->M);CHKERRQ(ierr);
170     ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr);
171     ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr);
172     ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr);
173   }
174   ierr = PetscFree(tao->data);CHKERRQ(ierr);
175   PetscFunctionReturn(0);
176 }
177 
178 /*------------------------------------------------------------*/
179 #undef __FUNCT__
180 #define __FUNCT__ "TaoSetFromOptions_BLMVM"
181 static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptions* PetscOptionsObject,Tao tao)
182 {
183   PetscErrorCode ierr;
184 
185   PetscFunctionBegin;
186   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr);
187   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
188   ierr = PetscOptionsTail();CHKERRQ(ierr);
189   PetscFunctionReturn(0);
190 }
191 
192 
193 /*------------------------------------------------------------*/
194 #undef __FUNCT__
195 #define __FUNCT__ "TaoView_BLMVM"
196 static int TaoView_BLMVM(Tao tao, PetscViewer viewer)
197 {
198   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
199   PetscBool      isascii;
200   PetscErrorCode ierr;
201 
202   PetscFunctionBegin;
203   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
204   if (isascii) {
205     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
206     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr);
207     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
208   }
209   PetscFunctionReturn(0);
210 }
211 
212 #undef __FUNCT__
213 #define __FUNCT__ "TaoComputeDual_BLMVM"
214 static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
215 {
216   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
217   PetscErrorCode ierr;
218 
219   PetscFunctionBegin;
220   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
221   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
222   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
223   if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
224 
225   ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr);
226   ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr);
227   ierr = VecSet(DXU,0.0);CHKERRQ(ierr);
228   ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr);
229 
230   ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr);
231   ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr);
232   ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr);
233   PetscFunctionReturn(0);
234 }
235 
236 /* ---------------------------------------------------------- */
237 /*MC
238   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
239          for nonlinear minimization with bound constraints. It is an extension
240          of TAOLMVM
241 
242   Options Database Keys:
243 +     -tao_lmm_vectors - number of vectors to use for approximation
244 .     -tao_lmm_scale_type - "none","scalar","broyden"
245 .     -tao_lmm_limit_type - "none","average","relative","absolute"
246 .     -tao_lmm_rescale_type - "none","scalar","gl"
247 .     -tao_lmm_limit_mu - mu limiting factor
248 .     -tao_lmm_limit_nu - nu limiting factor
249 .     -tao_lmm_delta_min - minimum delta value
250 .     -tao_lmm_delta_max - maximum delta value
251 .     -tao_lmm_broyden_phi - phi factor for Broyden scaling
252 .     -tao_lmm_scalar_alpha - alpha factor for scalar scaling
253 .     -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal
254 .     -tao_lmm_rescale_beta - beta factor for rescaling diagonal
255 .     -tao_lmm_scalar_history - amount of history for scalar scaling
256 .     -tao_lmm_rescale_history - amount of history for rescaling diagonal
257 -     -tao_lmm_eps - rejection tolerance
258 
259   Level: beginner
260 M*/
261 #undef __FUNCT__
262 #define __FUNCT__ "TaoCreate_BLMVM"
263 PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
264 {
265   TAO_BLMVM      *blmP;
266   const char     *morethuente_type = TAOLINESEARCHMT;
267   PetscErrorCode ierr;
268 
269   PetscFunctionBegin;
270   tao->ops->setup = TaoSetup_BLMVM;
271   tao->ops->solve = TaoSolve_BLMVM;
272   tao->ops->view = TaoView_BLMVM;
273   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
274   tao->ops->destroy = TaoDestroy_BLMVM;
275   tao->ops->computedual = TaoComputeDual_BLMVM;
276 
277   ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr);
278   tao->data = (void*)blmP;
279 
280   /* Override default settings (unless already changed) */
281   if (!tao->max_it_changed) tao->max_it = 2000;
282   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
283   if (!tao->fatol_changed) tao->fatol = 1.0e-4;
284   if (!tao->frtol_changed) tao->frtol = 1.0e-4;
285 
286   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
287   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
288   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
289   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
290   PetscFunctionReturn(0);
291 }
292 
293