xref: /petsc/src/tao/bound/impls/bqnk/bqnk.c (revision 3c859ba3a04a72e8efcb87bc7ffc046a6cbab413)
1414d97d3SAlp Dener #include <../src/tao/bound/impls/bqnk/bqnk.h> /*I "petsctao.h" I*/ /*I "petscmat.h" I*/
2e0ed867bSAlp Dener #include <petscksp.h>
3e0ed867bSAlp Dener 
4e0ed867bSAlp Dener static PetscErrorCode TaoBQNKComputeHessian(Tao tao)
5e0ed867bSAlp Dener {
6e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
7e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
8e0ed867bSAlp Dener   PetscErrorCode ierr;
965f5217aSAlp Dener   PetscReal      gnorm2, delta;
10e0ed867bSAlp Dener 
11e0ed867bSAlp Dener   PetscFunctionBegin;
12e0ed867bSAlp Dener   /* Alias the LMVM matrix into the TAO hessian */
13e0ed867bSAlp Dener   if (tao->hessian) {
14e0ed867bSAlp Dener     ierr = MatDestroy(&tao->hessian);CHKERRQ(ierr);
15e0ed867bSAlp Dener   }
16e0ed867bSAlp Dener   if (tao->hessian_pre) {
17e0ed867bSAlp Dener     ierr = MatDestroy(&tao->hessian_pre);CHKERRQ(ierr);
18e0ed867bSAlp Dener   }
19e0ed867bSAlp Dener   ierr = PetscObjectReference((PetscObject)bqnk->B);CHKERRQ(ierr);
20e0ed867bSAlp Dener   tao->hessian = bqnk->B;
21e0ed867bSAlp Dener   ierr = PetscObjectReference((PetscObject)bqnk->B);CHKERRQ(ierr);
22e0ed867bSAlp Dener   tao->hessian_pre = bqnk->B;
23e0ed867bSAlp Dener   /* Update the Hessian with the latest solution */
24f5766c09SAlp Dener   if (bqnk->is_spd) {
2565f5217aSAlp Dener     gnorm2 = bnk->gnorm*bnk->gnorm;
268cabe928SAlp Dener     if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON;
278cabe928SAlp Dener     if (bnk->f == 0.0) {
288cabe928SAlp Dener       delta = 2.0 / gnorm2;
298cabe928SAlp Dener     } else {
308cabe928SAlp Dener       delta = 2.0 * PetscAbsScalar(bnk->f) / gnorm2;
318cabe928SAlp Dener     }
32864588a7SAlp Dener     ierr = MatLMVMSymBroydenSetDelta(bqnk->B, delta);CHKERRQ(ierr);
33f5766c09SAlp Dener   }
34e0ed867bSAlp Dener   ierr = MatLMVMUpdate(tao->hessian, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
35e0ed867bSAlp Dener   ierr = MatLMVMResetShift(tao->hessian);CHKERRQ(ierr);
36e0ed867bSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
374f4fdda4SAlp Dener   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
38e0ed867bSAlp Dener   if (bnk->active_idx) {
39e0ed867bSAlp Dener     ierr = MatCreateSubMatrixVirtual(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, &bnk->H_inactive);CHKERRQ(ierr);
404f4fdda4SAlp Dener     ierr = PCLMVMSetIS(bqnk->pc, bnk->inactive_idx);CHKERRQ(ierr);
41e0ed867bSAlp Dener   } else {
42e0ed867bSAlp Dener     ierr = PetscObjectReference((PetscObject)tao->hessian);CHKERRQ(ierr);
43e0ed867bSAlp Dener     bnk->H_inactive = tao->hessian;
444f4fdda4SAlp Dener     ierr = PCLMVMClearIS(bqnk->pc);CHKERRQ(ierr);
45e0ed867bSAlp Dener   }
464f4fdda4SAlp Dener   ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
474f4fdda4SAlp Dener   ierr = PetscObjectReference((PetscObject)bnk->H_inactive);CHKERRQ(ierr);
48e0ed867bSAlp Dener   bnk->Hpre_inactive = bnk->H_inactive;
49e0ed867bSAlp Dener   PetscFunctionReturn(0);
50e0ed867bSAlp Dener }
51e0ed867bSAlp Dener 
526b591159SAlp Dener static PetscErrorCode TaoBQNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
53e0ed867bSAlp Dener {
54e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
55e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
56e0ed867bSAlp Dener   PetscErrorCode ierr;
57e0ed867bSAlp Dener 
58e0ed867bSAlp Dener   PetscFunctionBegin;
596b591159SAlp Dener   ierr = TaoBNKComputeStep(tao, shift, ksp_reason, step_type);CHKERRQ(ierr);
60e0ed867bSAlp Dener   if (*ksp_reason < 0) {
61e0ed867bSAlp Dener     /* Krylov solver failed to converge so reset the LMVM matrix */
62e0ed867bSAlp Dener     ierr = MatLMVMReset(bqnk->B, PETSC_FALSE);CHKERRQ(ierr);
636b591159SAlp Dener     ierr = MatLMVMUpdate(bqnk->B, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
64e0ed867bSAlp Dener   }
65e0ed867bSAlp Dener   PetscFunctionReturn(0);
66e0ed867bSAlp Dener }
67e0ed867bSAlp Dener 
68414d97d3SAlp Dener PetscErrorCode TaoSolve_BQNK(Tao tao)
69414d97d3SAlp Dener {
70414d97d3SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
71414d97d3SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
72414d97d3SAlp Dener   Mat_LMVM       *lmvm = (Mat_LMVM*)bqnk->B->data;
73414d97d3SAlp Dener   Mat_LMVM       *J0;
74414d97d3SAlp Dener   Mat_SymBrdn    *diag_ctx;
75414d97d3SAlp Dener   PetscBool      flg = PETSC_FALSE;
76414d97d3SAlp Dener   PetscErrorCode ierr;
77414d97d3SAlp Dener 
78414d97d3SAlp Dener   PetscFunctionBegin;
79414d97d3SAlp Dener   if (!tao->recycle) {
80414d97d3SAlp Dener     ierr = MatLMVMReset(bqnk->B, PETSC_FALSE);CHKERRQ(ierr);
81414d97d3SAlp Dener     lmvm->nresets = 0;
82414d97d3SAlp Dener     if (lmvm->J0) {
83414d97d3SAlp Dener       ierr = PetscObjectBaseTypeCompare((PetscObject)lmvm->J0, MATLMVM, &flg);CHKERRQ(ierr);
84414d97d3SAlp Dener       if (flg) {
85414d97d3SAlp Dener         J0 = (Mat_LMVM*)lmvm->J0->data;
86414d97d3SAlp Dener         J0->nresets = 0;
87414d97d3SAlp Dener       }
88414d97d3SAlp Dener     }
89414d97d3SAlp Dener     flg = PETSC_FALSE;
901e1ea65dSPierre Jolivet     ierr = PetscObjectTypeCompareAny((PetscObject)bqnk->B, &flg, MATLMVMSYMBROYDEN, MATLMVMSYMBADBROYDEN, MATLMVMBFGS, MATLMVMDFP, "");CHKERRQ(ierr);
91414d97d3SAlp Dener     if (flg) {
92414d97d3SAlp Dener       diag_ctx = (Mat_SymBrdn*)lmvm->ctx;
93414d97d3SAlp Dener       J0 = (Mat_LMVM*)diag_ctx->D->data;
94414d97d3SAlp Dener       J0->nresets = 0;
95414d97d3SAlp Dener     }
96414d97d3SAlp Dener   }
97414d97d3SAlp Dener   ierr = (*bqnk->solve)(tao);CHKERRQ(ierr);
98414d97d3SAlp Dener   PetscFunctionReturn(0);
99414d97d3SAlp Dener }
100414d97d3SAlp Dener 
1016b591159SAlp Dener PetscErrorCode TaoSetUp_BQNK(Tao tao)
1024f4fdda4SAlp Dener {
1034f4fdda4SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1044f4fdda4SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
1054f4fdda4SAlp Dener   PetscErrorCode ierr;
1064f4fdda4SAlp Dener   PetscInt       n, N;
1076b591159SAlp Dener   PetscBool      is_lmvm, is_sym, is_spd;
1084f4fdda4SAlp Dener 
1094f4fdda4SAlp Dener   PetscFunctionBegin;
1104f4fdda4SAlp Dener   ierr = TaoSetUp_BNK(tao);CHKERRQ(ierr);
1114f4fdda4SAlp Dener   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
1124f4fdda4SAlp Dener   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
1134f4fdda4SAlp Dener   ierr = MatSetSizes(bqnk->B, n, n, N, N);CHKERRQ(ierr);
1144f4fdda4SAlp Dener   ierr = MatLMVMAllocate(bqnk->B,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
1154f4fdda4SAlp Dener   ierr = PetscObjectBaseTypeCompare((PetscObject)bqnk->B, MATLMVM, &is_lmvm);CHKERRQ(ierr);
116*3c859ba3SBarry Smith   PetscCheck(is_lmvm,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Matrix must be an LMVM-type");
1174f4fdda4SAlp Dener   ierr = MatGetOption(bqnk->B, MAT_SYMMETRIC, &is_sym);CHKERRQ(ierr);
118*3c859ba3SBarry Smith   PetscCheck(is_sym,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric");
1196b591159SAlp Dener   ierr = MatGetOption(bqnk->B, MAT_SPD, &is_spd);CHKERRQ(ierr);
1204f4fdda4SAlp Dener   ierr = KSPGetPC(tao->ksp, &bqnk->pc);CHKERRQ(ierr);
1214f4fdda4SAlp Dener   ierr = PCSetType(bqnk->pc, PCLMVM);CHKERRQ(ierr);
1224f4fdda4SAlp Dener   ierr = PCLMVMSetMatLMVM(bqnk->pc, bqnk->B);CHKERRQ(ierr);
1234f4fdda4SAlp Dener   PetscFunctionReturn(0);
1244f4fdda4SAlp Dener }
1254f4fdda4SAlp Dener 
126e0ed867bSAlp Dener static PetscErrorCode TaoSetFromOptions_BQNK(PetscOptionItems *PetscOptionsObject,Tao tao)
127e0ed867bSAlp Dener {
128e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
129e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
130e0ed867bSAlp Dener   PetscErrorCode ierr;
131e0ed867bSAlp Dener 
132e0ed867bSAlp Dener   PetscFunctionBegin;
1339fa2b5dcSStefano Zampini   ierr = TaoSetFromOptions_BNK(PetscOptionsObject,tao);CHKERRQ(ierr);
134e0ed867bSAlp Dener   if (bnk->init_type == BNK_INIT_INTERPOLATION) bnk->init_type = BNK_INIT_DIRECTION;
1358ebe3e4eSStefano Zampini   ierr = MatSetOptionsPrefix(bqnk->B, ((PetscObject)tao)->prefix);CHKERRQ(ierr);
1368ebe3e4eSStefano Zampini   ierr = MatAppendOptionsPrefix(bqnk->B, "tao_bqnk_");CHKERRQ(ierr);
137e0ed867bSAlp Dener   ierr = MatSetFromOptions(bqnk->B);CHKERRQ(ierr);
138f5766c09SAlp Dener   ierr = MatGetOption(bqnk->B, MAT_SPD, &bqnk->is_spd);CHKERRQ(ierr);
139e0ed867bSAlp Dener   PetscFunctionReturn(0);
140e0ed867bSAlp Dener }
141e0ed867bSAlp Dener 
142e0ed867bSAlp Dener static PetscErrorCode TaoView_BQNK(Tao tao, PetscViewer viewer)
143e0ed867bSAlp Dener {
144e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
145e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
146e0ed867bSAlp Dener   PetscErrorCode ierr;
147e0ed867bSAlp Dener   PetscBool      isascii;
148e0ed867bSAlp Dener 
149e0ed867bSAlp Dener   PetscFunctionBegin;
150e0ed867bSAlp Dener   ierr = TaoView_BNK(tao, viewer);CHKERRQ(ierr);
151e0ed867bSAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
152e0ed867bSAlp Dener   if (isascii) {
153e0ed867bSAlp Dener     ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
154e0ed867bSAlp Dener     ierr = MatView(bqnk->B, viewer);CHKERRQ(ierr);
155e0ed867bSAlp Dener     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
156e0ed867bSAlp Dener   }
157e0ed867bSAlp Dener   PetscFunctionReturn(0);
158e0ed867bSAlp Dener }
159e0ed867bSAlp Dener 
160e0ed867bSAlp Dener static PetscErrorCode TaoDestroy_BQNK(Tao tao)
161e0ed867bSAlp Dener {
162e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
163e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
164e0ed867bSAlp Dener   PetscErrorCode ierr;
165e0ed867bSAlp Dener 
166e0ed867bSAlp Dener   PetscFunctionBegin;
1674f4fdda4SAlp Dener   ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
168cb384e1eSAlp Dener   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
169e0ed867bSAlp Dener   ierr = MatDestroy(&bqnk->B);CHKERRQ(ierr);
170e0ed867bSAlp Dener   ierr = PetscFree(bnk->ctx);CHKERRQ(ierr);
171e0ed867bSAlp Dener   ierr = TaoDestroy_BNK(tao);CHKERRQ(ierr);
172e0ed867bSAlp Dener   PetscFunctionReturn(0);
173e0ed867bSAlp Dener }
174e0ed867bSAlp Dener 
175e0ed867bSAlp Dener PETSC_INTERN PetscErrorCode TaoCreate_BQNK(Tao tao)
176e0ed867bSAlp Dener {
177e0ed867bSAlp Dener   TAO_BNK        *bnk;
178e0ed867bSAlp Dener   TAO_BQNK       *bqnk;
179e0ed867bSAlp Dener   PetscErrorCode ierr;
180e0ed867bSAlp Dener 
181e0ed867bSAlp Dener   PetscFunctionBegin;
182e0ed867bSAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
183414d97d3SAlp Dener   tao->ops->solve = TaoSolve_BQNK;
184e0ed867bSAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BQNK;
185e0ed867bSAlp Dener   tao->ops->destroy = TaoDestroy_BQNK;
186e0ed867bSAlp Dener   tao->ops->view = TaoView_BQNK;
1874f4fdda4SAlp Dener   tao->ops->setup = TaoSetUp_BQNK;
188e0ed867bSAlp Dener 
189e0ed867bSAlp Dener   bnk = (TAO_BNK *)tao->data;
190e0ed867bSAlp Dener   bnk->computehessian = TaoBQNKComputeHessian;
191e0ed867bSAlp Dener   bnk->computestep = TaoBQNKComputeStep;
192e0ed867bSAlp Dener   bnk->init_type = BNK_INIT_DIRECTION;
193e0ed867bSAlp Dener 
194e0ed867bSAlp Dener   ierr = PetscNewLog(tao,&bqnk);CHKERRQ(ierr);
195e0ed867bSAlp Dener   bnk->ctx = (void*)bqnk;
196f5766c09SAlp Dener   bqnk->is_spd = PETSC_TRUE;
197e0ed867bSAlp Dener 
198e0ed867bSAlp Dener   ierr = MatCreate(PetscObjectComm((PetscObject)tao), &bqnk->B);CHKERRQ(ierr);
199e0ed867bSAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)bqnk->B, (PetscObject)tao, 1);CHKERRQ(ierr);
200e0ed867bSAlp Dener   ierr = MatSetType(bqnk->B, MATLMVMSR1);CHKERRQ(ierr);
201e0ed867bSAlp Dener   PetscFunctionReturn(0);
202e0ed867bSAlp Dener }
203f5766c09SAlp Dener 
204414d97d3SAlp Dener /*@
205414d97d3SAlp Dener    TaoGetLMVMMatrix - Returns a pointer to the internal LMVM matrix. Valid
206414d97d3SAlp Dener    only for quasi-Newton family of methods.
207414d97d3SAlp Dener 
208414d97d3SAlp Dener    Input Parameters:
209414d97d3SAlp Dener .  tao - Tao solver context
210414d97d3SAlp Dener 
211414d97d3SAlp Dener    Output Parameters:
212414d97d3SAlp Dener .  B - LMVM matrix
213414d97d3SAlp Dener 
214414d97d3SAlp Dener    Level: advanced
215414d97d3SAlp Dener 
216414d97d3SAlp Dener .seealso: TAOBQNLS, TAOBQNKLS, TAOBQNKTL, TAOBQNKTR, MATLMVM, TaoSetLMVMMatrix()
217414d97d3SAlp Dener @*/
218f5766c09SAlp Dener PetscErrorCode TaoGetLMVMMatrix(Tao tao, Mat *B)
219f5766c09SAlp Dener {
220f5766c09SAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
221f5766c09SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
222f5766c09SAlp Dener   PetscErrorCode ierr;
223414d97d3SAlp Dener   PetscBool      flg = PETSC_FALSE;
224f5766c09SAlp Dener 
225f5766c09SAlp Dener   PetscFunctionBegin;
226414d97d3SAlp Dener   ierr = PetscObjectTypeCompareAny((PetscObject)tao, &flg, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, "");CHKERRQ(ierr);
227*3c859ba3SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM Matrix only exists for quasi-Newton algorithms");
228f5766c09SAlp Dener   *B = bqnk->B;
229f5766c09SAlp Dener   PetscFunctionReturn(0);
230f5766c09SAlp Dener }
231414d97d3SAlp Dener 
232414d97d3SAlp Dener /*@
233414d97d3SAlp Dener    TaoSetLMVMMatrix - Sets an external LMVM matrix into the Tao solver. Valid
234414d97d3SAlp Dener    only for quasi-Newton family of methods.
235414d97d3SAlp Dener 
236414d97d3SAlp Dener    QN family of methods create their own LMVM matrices and users who wish to
237414d97d3SAlp Dener    manipulate this matrix should use TaoGetLMVMMatrix() instead.
238414d97d3SAlp Dener 
239414d97d3SAlp Dener    Input Parameters:
240414d97d3SAlp Dener +  tao - Tao solver context
241414d97d3SAlp Dener -  B - LMVM matrix
242414d97d3SAlp Dener 
243414d97d3SAlp Dener    Level: advanced
244414d97d3SAlp Dener 
245414d97d3SAlp Dener .seealso: TAOBQNLS, TAOBQNKLS, TAOBQNKTL, TAOBQNKTR, MATLMVM, TaoGetLMVMMatrix()
246414d97d3SAlp Dener @*/
247414d97d3SAlp Dener PetscErrorCode TaoSetLMVMMatrix(Tao tao, Mat B)
248414d97d3SAlp Dener {
249414d97d3SAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
250414d97d3SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
251414d97d3SAlp Dener   PetscErrorCode ierr;
252414d97d3SAlp Dener   PetscBool      flg = PETSC_FALSE;
253414d97d3SAlp Dener 
254414d97d3SAlp Dener   PetscFunctionBegin;
255414d97d3SAlp Dener   ierr = PetscObjectTypeCompareAny((PetscObject)tao, &flg, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, "");CHKERRQ(ierr);
256*3c859ba3SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM Matrix only exists for quasi-Newton algorithms");
257414d97d3SAlp Dener   ierr = PetscObjectBaseTypeCompare((PetscObject)B, MATLMVM, &flg);CHKERRQ(ierr);
258*3c859ba3SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given matrix is not an LMVM matrix");
259414d97d3SAlp Dener   if (bqnk->B) {
260414d97d3SAlp Dener     ierr = MatDestroy(&bqnk->B);CHKERRQ(ierr);
261414d97d3SAlp Dener   }
262414d97d3SAlp Dener   ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
263414d97d3SAlp Dener   bqnk->B = B;
264414d97d3SAlp Dener   PetscFunctionReturn(0);
265414d97d3SAlp Dener }
266