xref: /petsc/src/tao/bound/impls/bqnk/bqnk.c (revision b94d7ded0a05f1bbd5e48daa6f92b28259c75b44)
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;
865f5217aSAlp Dener   PetscReal      gnorm2, delta;
9e0ed867bSAlp Dener 
10e0ed867bSAlp Dener   PetscFunctionBegin;
11e0ed867bSAlp Dener   /* Alias the LMVM matrix into the TAO hessian */
12e0ed867bSAlp Dener   if (tao->hessian) {
139566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&tao->hessian));
14e0ed867bSAlp Dener   }
15e0ed867bSAlp Dener   if (tao->hessian_pre) {
169566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&tao->hessian_pre));
17e0ed867bSAlp Dener   }
189566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)bqnk->B));
19e0ed867bSAlp Dener   tao->hessian = bqnk->B;
209566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)bqnk->B));
21e0ed867bSAlp Dener   tao->hessian_pre = bqnk->B;
22e0ed867bSAlp Dener   /* Update the Hessian with the latest solution */
23f5766c09SAlp Dener   if (bqnk->is_spd) {
2465f5217aSAlp Dener     gnorm2 = bnk->gnorm*bnk->gnorm;
258cabe928SAlp Dener     if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON;
268cabe928SAlp Dener     if (bnk->f == 0.0) {
278cabe928SAlp Dener       delta = 2.0 / gnorm2;
288cabe928SAlp Dener     } else {
298cabe928SAlp Dener       delta = 2.0 * PetscAbsScalar(bnk->f) / gnorm2;
308cabe928SAlp Dener     }
319566063dSJacob Faibussowitsch     PetscCall(MatLMVMSymBroydenSetDelta(bqnk->B, delta));
32f5766c09SAlp Dener   }
339566063dSJacob Faibussowitsch   PetscCall(MatLMVMUpdate(tao->hessian, tao->solution, bnk->unprojected_gradient));
349566063dSJacob Faibussowitsch   PetscCall(MatLMVMResetShift(tao->hessian));
35e0ed867bSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
369566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->H_inactive));
37e0ed867bSAlp Dener   if (bnk->active_idx) {
389566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrixVirtual(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, &bnk->H_inactive));
399566063dSJacob Faibussowitsch     PetscCall(PCLMVMSetIS(bqnk->pc, bnk->inactive_idx));
40e0ed867bSAlp Dener   } else {
419566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)tao->hessian));
42e0ed867bSAlp Dener     bnk->H_inactive = tao->hessian;
439566063dSJacob Faibussowitsch     PetscCall(PCLMVMClearIS(bqnk->pc));
44e0ed867bSAlp Dener   }
459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->Hpre_inactive));
469566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive));
47e0ed867bSAlp Dener   bnk->Hpre_inactive = bnk->H_inactive;
48e0ed867bSAlp Dener   PetscFunctionReturn(0);
49e0ed867bSAlp Dener }
50e0ed867bSAlp Dener 
516b591159SAlp Dener static PetscErrorCode TaoBQNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
52e0ed867bSAlp Dener {
53e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
54e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
55e0ed867bSAlp Dener 
56e0ed867bSAlp Dener   PetscFunctionBegin;
579566063dSJacob Faibussowitsch   PetscCall(TaoBNKComputeStep(tao, shift, ksp_reason, step_type));
58e0ed867bSAlp Dener   if (*ksp_reason < 0) {
59e0ed867bSAlp Dener     /* Krylov solver failed to converge so reset the LMVM matrix */
609566063dSJacob Faibussowitsch     PetscCall(MatLMVMReset(bqnk->B, PETSC_FALSE));
619566063dSJacob Faibussowitsch     PetscCall(MatLMVMUpdate(bqnk->B, tao->solution, bnk->unprojected_gradient));
62e0ed867bSAlp Dener   }
63e0ed867bSAlp Dener   PetscFunctionReturn(0);
64e0ed867bSAlp Dener }
65e0ed867bSAlp Dener 
66414d97d3SAlp Dener PetscErrorCode TaoSolve_BQNK(Tao tao)
67414d97d3SAlp Dener {
68414d97d3SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
69414d97d3SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
70414d97d3SAlp Dener   Mat_LMVM       *lmvm = (Mat_LMVM*)bqnk->B->data;
71414d97d3SAlp Dener   Mat_LMVM       *J0;
72414d97d3SAlp Dener   Mat_SymBrdn    *diag_ctx;
73414d97d3SAlp Dener   PetscBool      flg = PETSC_FALSE;
74414d97d3SAlp Dener 
75414d97d3SAlp Dener   PetscFunctionBegin;
76414d97d3SAlp Dener   if (!tao->recycle) {
779566063dSJacob Faibussowitsch     PetscCall(MatLMVMReset(bqnk->B, PETSC_FALSE));
78414d97d3SAlp Dener     lmvm->nresets = 0;
79414d97d3SAlp Dener     if (lmvm->J0) {
809566063dSJacob Faibussowitsch       PetscCall(PetscObjectBaseTypeCompare((PetscObject)lmvm->J0, MATLMVM, &flg));
81414d97d3SAlp Dener       if (flg) {
82414d97d3SAlp Dener         J0 = (Mat_LMVM*)lmvm->J0->data;
83414d97d3SAlp Dener         J0->nresets = 0;
84414d97d3SAlp Dener       }
85414d97d3SAlp Dener     }
86414d97d3SAlp Dener     flg = PETSC_FALSE;
879566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompareAny((PetscObject)bqnk->B, &flg, MATLMVMSYMBROYDEN, MATLMVMSYMBADBROYDEN, MATLMVMBFGS, MATLMVMDFP, ""));
88414d97d3SAlp Dener     if (flg) {
89414d97d3SAlp Dener       diag_ctx = (Mat_SymBrdn*)lmvm->ctx;
90414d97d3SAlp Dener       J0 = (Mat_LMVM*)diag_ctx->D->data;
91414d97d3SAlp Dener       J0->nresets = 0;
92414d97d3SAlp Dener     }
93414d97d3SAlp Dener   }
949566063dSJacob Faibussowitsch   PetscCall((*bqnk->solve)(tao));
95414d97d3SAlp Dener   PetscFunctionReturn(0);
96414d97d3SAlp Dener }
97414d97d3SAlp Dener 
986b591159SAlp Dener PetscErrorCode TaoSetUp_BQNK(Tao tao)
994f4fdda4SAlp Dener {
1004f4fdda4SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1014f4fdda4SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
1024f4fdda4SAlp Dener   PetscInt       n, N;
103*b94d7dedSBarry Smith   PetscBool      is_lmvm, is_set,is_sym;
1044f4fdda4SAlp Dener 
1054f4fdda4SAlp Dener   PetscFunctionBegin;
1069566063dSJacob Faibussowitsch   PetscCall(TaoSetUp_BNK(tao));
1079566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(tao->solution,&n));
1089566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tao->solution,&N));
1099566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(bqnk->B, n, n, N, N));
1109566063dSJacob Faibussowitsch   PetscCall(MatLMVMAllocate(bqnk->B,tao->solution,bnk->unprojected_gradient));
1119566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)bqnk->B, MATLMVM, &is_lmvm));
1123c859ba3SBarry Smith   PetscCheck(is_lmvm,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Matrix must be an LMVM-type");
113*b94d7dedSBarry Smith   PetscCall(MatIsSymmetricKnown(bqnk->B, &is_set,&is_sym));
114*b94d7dedSBarry Smith   PetscCheck(is_set && is_sym,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric");
1159566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(tao->ksp, &bqnk->pc));
1169566063dSJacob Faibussowitsch   PetscCall(PCSetType(bqnk->pc, PCLMVM));
1179566063dSJacob Faibussowitsch   PetscCall(PCLMVMSetMatLMVM(bqnk->pc, bqnk->B));
1184f4fdda4SAlp Dener   PetscFunctionReturn(0);
1194f4fdda4SAlp Dener }
1204f4fdda4SAlp Dener 
121e0ed867bSAlp Dener static PetscErrorCode TaoSetFromOptions_BQNK(PetscOptionItems *PetscOptionsObject,Tao tao)
122e0ed867bSAlp Dener {
123e0ed867bSAlp Dener   TAO_BNK   *bnk = (TAO_BNK *)tao->data;
124e0ed867bSAlp Dener   TAO_BQNK  *bqnk = (TAO_BQNK*)bnk->ctx;
125*b94d7dedSBarry Smith   PetscBool is_set;
126e0ed867bSAlp Dener 
127e0ed867bSAlp Dener   PetscFunctionBegin;
1289566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions_BNK(PetscOptionsObject,tao));
129e0ed867bSAlp Dener   if (bnk->init_type == BNK_INIT_INTERPOLATION) bnk->init_type = BNK_INIT_DIRECTION;
1309566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(bqnk->B, ((PetscObject)tao)->prefix));
1319566063dSJacob Faibussowitsch   PetscCall(MatAppendOptionsPrefix(bqnk->B, "tao_bqnk_"));
1329566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(bqnk->B));
133*b94d7dedSBarry Smith   PetscCall(MatIsSPDKnown(bqnk->B, &is_set, &bqnk->is_spd));
134*b94d7dedSBarry Smith   if (!is_set) bqnk->is_spd = PETSC_FALSE;
135e0ed867bSAlp Dener   PetscFunctionReturn(0);
136e0ed867bSAlp Dener }
137e0ed867bSAlp Dener 
138e0ed867bSAlp Dener static PetscErrorCode TaoView_BQNK(Tao tao, PetscViewer viewer)
139e0ed867bSAlp Dener {
140e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
141e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
142e0ed867bSAlp Dener   PetscBool      isascii;
143e0ed867bSAlp Dener 
144e0ed867bSAlp Dener   PetscFunctionBegin;
1459566063dSJacob Faibussowitsch   PetscCall(TaoView_BNK(tao, viewer));
1469566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
147e0ed867bSAlp Dener   if (isascii) {
1489566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
1499566063dSJacob Faibussowitsch     PetscCall(MatView(bqnk->B, viewer));
1509566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
151e0ed867bSAlp Dener   }
152e0ed867bSAlp Dener   PetscFunctionReturn(0);
153e0ed867bSAlp Dener }
154e0ed867bSAlp Dener 
155e0ed867bSAlp Dener static PetscErrorCode TaoDestroy_BQNK(Tao tao)
156e0ed867bSAlp Dener {
157e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
158e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
159e0ed867bSAlp Dener 
160e0ed867bSAlp Dener   PetscFunctionBegin;
1619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->Hpre_inactive));
1629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->H_inactive));
1639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bqnk->B));
1649566063dSJacob Faibussowitsch   PetscCall(PetscFree(bnk->ctx));
1659566063dSJacob Faibussowitsch   PetscCall(TaoDestroy_BNK(tao));
166e0ed867bSAlp Dener   PetscFunctionReturn(0);
167e0ed867bSAlp Dener }
168e0ed867bSAlp Dener 
169e0ed867bSAlp Dener PETSC_INTERN PetscErrorCode TaoCreate_BQNK(Tao tao)
170e0ed867bSAlp Dener {
171e0ed867bSAlp Dener   TAO_BNK        *bnk;
172e0ed867bSAlp Dener   TAO_BQNK       *bqnk;
173e0ed867bSAlp Dener 
174e0ed867bSAlp Dener   PetscFunctionBegin;
1759566063dSJacob Faibussowitsch   PetscCall(TaoCreate_BNK(tao));
176414d97d3SAlp Dener   tao->ops->solve = TaoSolve_BQNK;
177e0ed867bSAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BQNK;
178e0ed867bSAlp Dener   tao->ops->destroy = TaoDestroy_BQNK;
179e0ed867bSAlp Dener   tao->ops->view = TaoView_BQNK;
1804f4fdda4SAlp Dener   tao->ops->setup = TaoSetUp_BQNK;
181e0ed867bSAlp Dener 
182e0ed867bSAlp Dener   bnk = (TAO_BNK *)tao->data;
183e0ed867bSAlp Dener   bnk->computehessian = TaoBQNKComputeHessian;
184e0ed867bSAlp Dener   bnk->computestep = TaoBQNKComputeStep;
185e0ed867bSAlp Dener   bnk->init_type = BNK_INIT_DIRECTION;
186e0ed867bSAlp Dener 
1879566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tao,&bqnk));
188e0ed867bSAlp Dener   bnk->ctx = (void*)bqnk;
189f5766c09SAlp Dener   bqnk->is_spd = PETSC_TRUE;
190e0ed867bSAlp Dener 
1919566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &bqnk->B));
1929566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)bqnk->B, (PetscObject)tao, 1));
1939566063dSJacob Faibussowitsch   PetscCall(MatSetType(bqnk->B, MATLMVMSR1));
194e0ed867bSAlp Dener   PetscFunctionReturn(0);
195e0ed867bSAlp Dener }
196f5766c09SAlp Dener 
197414d97d3SAlp Dener /*@
198414d97d3SAlp Dener    TaoGetLMVMMatrix - Returns a pointer to the internal LMVM matrix. Valid
199414d97d3SAlp Dener    only for quasi-Newton family of methods.
200414d97d3SAlp Dener 
201414d97d3SAlp Dener    Input Parameters:
202414d97d3SAlp Dener .  tao - Tao solver context
203414d97d3SAlp Dener 
204414d97d3SAlp Dener    Output Parameters:
205414d97d3SAlp Dener .  B - LMVM matrix
206414d97d3SAlp Dener 
207414d97d3SAlp Dener    Level: advanced
208414d97d3SAlp Dener 
209db781477SPatrick Sanan .seealso: `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTL`, `TAOBQNKTR`, `MATLMVM`, `TaoSetLMVMMatrix()`
210414d97d3SAlp Dener @*/
211f5766c09SAlp Dener PetscErrorCode TaoGetLMVMMatrix(Tao tao, Mat *B)
212f5766c09SAlp Dener {
213f5766c09SAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
214f5766c09SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
215414d97d3SAlp Dener   PetscBool      flg = PETSC_FALSE;
216f5766c09SAlp Dener 
217f5766c09SAlp Dener   PetscFunctionBegin;
2189566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)tao, &flg, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, ""));
2193c859ba3SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM Matrix only exists for quasi-Newton algorithms");
220f5766c09SAlp Dener   *B = bqnk->B;
221f5766c09SAlp Dener   PetscFunctionReturn(0);
222f5766c09SAlp Dener }
223414d97d3SAlp Dener 
224414d97d3SAlp Dener /*@
225414d97d3SAlp Dener    TaoSetLMVMMatrix - Sets an external LMVM matrix into the Tao solver. Valid
226414d97d3SAlp Dener    only for quasi-Newton family of methods.
227414d97d3SAlp Dener 
228414d97d3SAlp Dener    QN family of methods create their own LMVM matrices and users who wish to
229414d97d3SAlp Dener    manipulate this matrix should use TaoGetLMVMMatrix() instead.
230414d97d3SAlp Dener 
231414d97d3SAlp Dener    Input Parameters:
232414d97d3SAlp Dener +  tao - Tao solver context
233414d97d3SAlp Dener -  B - LMVM matrix
234414d97d3SAlp Dener 
235414d97d3SAlp Dener    Level: advanced
236414d97d3SAlp Dener 
237db781477SPatrick Sanan .seealso: `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTL`, `TAOBQNKTR`, `MATLMVM`, `TaoGetLMVMMatrix()`
238414d97d3SAlp Dener @*/
239414d97d3SAlp Dener PetscErrorCode TaoSetLMVMMatrix(Tao tao, Mat B)
240414d97d3SAlp Dener {
241414d97d3SAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
242414d97d3SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
243414d97d3SAlp Dener   PetscBool      flg = PETSC_FALSE;
244414d97d3SAlp Dener 
245414d97d3SAlp Dener   PetscFunctionBegin;
2469566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)tao, &flg, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, ""));
2473c859ba3SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM Matrix only exists for quasi-Newton algorithms");
2489566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)B, MATLMVM, &flg));
2493c859ba3SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given matrix is not an LMVM matrix");
250414d97d3SAlp Dener   if (bqnk->B) {
2519566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&bqnk->B));
252414d97d3SAlp Dener   }
2539566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)B));
254414d97d3SAlp Dener   bqnk->B = B;
255414d97d3SAlp Dener   PetscFunctionReturn(0);
256414d97d3SAlp Dener }
257