1954e39ddSJose E. Roman #include <../src/tao/bound/impls/bqnk/bqnk.h> /*I "petsctao.h" I*/
2e0ed867bSAlp Dener #include <petscksp.h>
3e0ed867bSAlp Dener
TaoBQNKComputeHessian(Tao tao)4d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoBQNKComputeHessian(Tao tao)
5d71ae5a4SJacob Faibussowitsch {
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 */
1248a46eb9SPierre Jolivet if (tao->hessian) PetscCall(MatDestroy(&tao->hessian));
1348a46eb9SPierre Jolivet if (tao->hessian_pre) PetscCall(MatDestroy(&tao->hessian_pre));
149566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)bqnk->B));
15e0ed867bSAlp Dener tao->hessian = bqnk->B;
169566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)bqnk->B));
17e0ed867bSAlp Dener tao->hessian_pre = bqnk->B;
18e0ed867bSAlp Dener /* Update the Hessian with the latest solution */
19f5766c09SAlp Dener if (bqnk->is_spd) {
2065f5217aSAlp Dener gnorm2 = bnk->gnorm * bnk->gnorm;
218cabe928SAlp Dener if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON;
228cabe928SAlp Dener if (bnk->f == 0.0) {
238cabe928SAlp Dener delta = 2.0 / gnorm2;
248cabe928SAlp Dener } else {
258cabe928SAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / gnorm2;
268cabe928SAlp Dener }
279566063dSJacob Faibussowitsch PetscCall(MatLMVMSymBroydenSetDelta(bqnk->B, delta));
28f5766c09SAlp Dener }
299566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(tao->hessian, tao->solution, bnk->unprojected_gradient));
309566063dSJacob Faibussowitsch PetscCall(MatLMVMResetShift(tao->hessian));
31e0ed867bSAlp Dener /* Prepare the reduced sub-matrices for the inactive set */
329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bnk->H_inactive));
33e0ed867bSAlp Dener if (bnk->active_idx) {
349566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrixVirtual(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, &bnk->H_inactive));
359566063dSJacob Faibussowitsch PetscCall(PCLMVMSetIS(bqnk->pc, bnk->inactive_idx));
36e0ed867bSAlp Dener } else {
379566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)tao->hessian));
38e0ed867bSAlp Dener bnk->H_inactive = tao->hessian;
399566063dSJacob Faibussowitsch PetscCall(PCLMVMClearIS(bqnk->pc));
40e0ed867bSAlp Dener }
419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bnk->Hpre_inactive));
429566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive));
43e0ed867bSAlp Dener bnk->Hpre_inactive = bnk->H_inactive;
443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
45e0ed867bSAlp Dener }
46e0ed867bSAlp Dener
TaoBQNKComputeStep(Tao tao,PetscBool shift,KSPConvergedReason * ksp_reason,PetscInt * step_type)47d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoBQNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
48d71ae5a4SJacob Faibussowitsch {
49e0ed867bSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
50e0ed867bSAlp Dener TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
51e0ed867bSAlp Dener
52e0ed867bSAlp Dener PetscFunctionBegin;
539566063dSJacob Faibussowitsch PetscCall(TaoBNKComputeStep(tao, shift, ksp_reason, step_type));
54e0ed867bSAlp Dener if (*ksp_reason < 0) {
55e0ed867bSAlp Dener /* Krylov solver failed to converge so reset the LMVM matrix */
569566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(bqnk->B, PETSC_FALSE));
579566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(bqnk->B, tao->solution, bnk->unprojected_gradient));
58e0ed867bSAlp Dener }
593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
60e0ed867bSAlp Dener }
61e0ed867bSAlp Dener
TaoSolve_BQNK(Tao tao)62d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSolve_BQNK(Tao tao)
63d71ae5a4SJacob Faibussowitsch {
64414d97d3SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
65414d97d3SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
66414d97d3SAlp Dener Mat_LMVM *lmvm = (Mat_LMVM *)bqnk->B->data;
67414d97d3SAlp Dener Mat_LMVM *J0;
68414d97d3SAlp Dener PetscBool flg = PETSC_FALSE;
69414d97d3SAlp Dener
70414d97d3SAlp Dener PetscFunctionBegin;
71414d97d3SAlp Dener if (!tao->recycle) {
729566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(bqnk->B, PETSC_FALSE));
73414d97d3SAlp Dener lmvm->nresets = 0;
74414d97d3SAlp Dener if (lmvm->J0) {
759566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)lmvm->J0, MATLMVM, &flg));
76414d97d3SAlp Dener if (flg) {
77414d97d3SAlp Dener J0 = (Mat_LMVM *)lmvm->J0->data;
78414d97d3SAlp Dener J0->nresets = 0;
79414d97d3SAlp Dener }
80414d97d3SAlp Dener }
81414d97d3SAlp Dener }
829566063dSJacob Faibussowitsch PetscCall((*bqnk->solve)(tao));
833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
84414d97d3SAlp Dener }
85414d97d3SAlp Dener
TaoSetUp_BQNK(Tao tao)86d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetUp_BQNK(Tao tao)
87d71ae5a4SJacob Faibussowitsch {
884f4fdda4SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
894f4fdda4SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
904f4fdda4SAlp Dener PetscInt n, N;
91b94d7dedSBarry Smith PetscBool is_lmvm, is_set, is_sym;
924f4fdda4SAlp Dener
934f4fdda4SAlp Dener PetscFunctionBegin;
949566063dSJacob Faibussowitsch PetscCall(TaoSetUp_BNK(tao));
959566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution, &n));
969566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &N));
979566063dSJacob Faibussowitsch PetscCall(MatSetSizes(bqnk->B, n, n, N, N));
989566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(bqnk->B, tao->solution, bnk->unprojected_gradient));
999566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)bqnk->B, MATLMVM, &is_lmvm));
1003c859ba3SBarry Smith PetscCheck(is_lmvm, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Matrix must be an LMVM-type");
101b94d7dedSBarry Smith PetscCall(MatIsSymmetricKnown(bqnk->B, &is_set, &is_sym));
102b94d7dedSBarry Smith PetscCheck(is_set && is_sym, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric");
1039566063dSJacob Faibussowitsch PetscCall(KSPGetPC(tao->ksp, &bqnk->pc));
1049566063dSJacob Faibussowitsch PetscCall(PCSetType(bqnk->pc, PCLMVM));
1059566063dSJacob Faibussowitsch PetscCall(PCLMVMSetMatLMVM(bqnk->pc, bqnk->B));
1063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1074f4fdda4SAlp Dener }
1084f4fdda4SAlp Dener
TaoSetFromOptions_BQNK(Tao tao,PetscOptionItems PetscOptionsObject)109*ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_BQNK(Tao tao, PetscOptionItems PetscOptionsObject)
110d71ae5a4SJacob Faibussowitsch {
111e0ed867bSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
112e0ed867bSAlp Dener TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
113b94d7dedSBarry Smith PetscBool is_set;
114e0ed867bSAlp Dener
115e0ed867bSAlp Dener PetscFunctionBegin;
116dbbe0bcdSBarry Smith PetscCall(TaoSetFromOptions_BNK(tao, PetscOptionsObject));
117e0ed867bSAlp Dener if (bnk->init_type == BNK_INIT_INTERPOLATION) bnk->init_type = BNK_INIT_DIRECTION;
1189566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(bqnk->B, ((PetscObject)tao)->prefix));
1199566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(bqnk->B, "tao_bqnk_"));
1209566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(bqnk->B));
121b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(bqnk->B, &is_set, &bqnk->is_spd));
122b94d7dedSBarry Smith if (!is_set) bqnk->is_spd = PETSC_FALSE;
1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
124e0ed867bSAlp Dener }
125e0ed867bSAlp Dener
TaoView_BQNK(Tao tao,PetscViewer viewer)126d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_BQNK(Tao tao, PetscViewer viewer)
127d71ae5a4SJacob Faibussowitsch {
128e0ed867bSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
129e0ed867bSAlp Dener TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
130e0ed867bSAlp Dener PetscBool isascii;
131e0ed867bSAlp Dener
132e0ed867bSAlp Dener PetscFunctionBegin;
1339566063dSJacob Faibussowitsch PetscCall(TaoView_BNK(tao, viewer));
1349566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
135e0ed867bSAlp Dener if (isascii) {
1369566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
1379566063dSJacob Faibussowitsch PetscCall(MatView(bqnk->B, viewer));
1389566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer));
139e0ed867bSAlp Dener }
1403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
141e0ed867bSAlp Dener }
142e0ed867bSAlp Dener
TaoDestroy_BQNK(Tao tao)143d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_BQNK(Tao tao)
144d71ae5a4SJacob Faibussowitsch {
145e0ed867bSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
146e0ed867bSAlp Dener TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
147e0ed867bSAlp Dener
148e0ed867bSAlp Dener PetscFunctionBegin;
1499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bnk->Hpre_inactive));
1509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bnk->H_inactive));
1519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bqnk->B));
1529566063dSJacob Faibussowitsch PetscCall(PetscFree(bnk->ctx));
1539566063dSJacob Faibussowitsch PetscCall(TaoDestroy_BNK(tao));
1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
155e0ed867bSAlp Dener }
156e0ed867bSAlp Dener
TaoCreate_BQNK(Tao tao)157d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode TaoCreate_BQNK(Tao tao)
158d71ae5a4SJacob Faibussowitsch {
159e0ed867bSAlp Dener TAO_BNK *bnk;
160e0ed867bSAlp Dener TAO_BQNK *bqnk;
161e0ed867bSAlp Dener
162e0ed867bSAlp Dener PetscFunctionBegin;
1639566063dSJacob Faibussowitsch PetscCall(TaoCreate_BNK(tao));
164414d97d3SAlp Dener tao->ops->solve = TaoSolve_BQNK;
165e0ed867bSAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BQNK;
166e0ed867bSAlp Dener tao->ops->destroy = TaoDestroy_BQNK;
167e0ed867bSAlp Dener tao->ops->view = TaoView_BQNK;
1684f4fdda4SAlp Dener tao->ops->setup = TaoSetUp_BQNK;
169e0ed867bSAlp Dener
170e0ed867bSAlp Dener bnk = (TAO_BNK *)tao->data;
171e0ed867bSAlp Dener bnk->computehessian = TaoBQNKComputeHessian;
172e0ed867bSAlp Dener bnk->computestep = TaoBQNKComputeStep;
173e0ed867bSAlp Dener bnk->init_type = BNK_INIT_DIRECTION;
174e0ed867bSAlp Dener
1754dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&bqnk));
176e0ed867bSAlp Dener bnk->ctx = (void *)bqnk;
177f5766c09SAlp Dener bqnk->is_spd = PETSC_TRUE;
178e0ed867bSAlp Dener
1799566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &bqnk->B));
1809566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)bqnk->B, (PetscObject)tao, 1));
1819566063dSJacob Faibussowitsch PetscCall(MatSetType(bqnk->B, MATLMVMSR1));
1823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
183e0ed867bSAlp Dener }
184f5766c09SAlp Dener
185414d97d3SAlp Dener /*@
186414d97d3SAlp Dener TaoGetLMVMMatrix - Returns a pointer to the internal LMVM matrix. Valid
187414d97d3SAlp Dener only for quasi-Newton family of methods.
188414d97d3SAlp Dener
1892fe279fdSBarry Smith Input Parameter:
1902fe279fdSBarry Smith . tao - `Tao` solver context
191414d97d3SAlp Dener
1922fe279fdSBarry Smith Output Parameter:
193414d97d3SAlp Dener . B - LMVM matrix
194414d97d3SAlp Dener
195414d97d3SAlp Dener Level: advanced
196414d97d3SAlp Dener
197db781477SPatrick Sanan .seealso: `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTL`, `TAOBQNKTR`, `MATLMVM`, `TaoSetLMVMMatrix()`
198414d97d3SAlp Dener @*/
TaoGetLMVMMatrix(Tao tao,Mat * B)199d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetLMVMMatrix(Tao tao, Mat *B)
200d71ae5a4SJacob Faibussowitsch {
201f5766c09SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
202f5766c09SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
203414d97d3SAlp Dener PetscBool flg = PETSC_FALSE;
204f5766c09SAlp Dener
205f5766c09SAlp Dener PetscFunctionBegin;
2069566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)tao, &flg, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, ""));
2073c859ba3SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM Matrix only exists for quasi-Newton algorithms");
208f5766c09SAlp Dener *B = bqnk->B;
2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
210f5766c09SAlp Dener }
211414d97d3SAlp Dener
212414d97d3SAlp Dener /*@
213414d97d3SAlp Dener TaoSetLMVMMatrix - Sets an external LMVM matrix into the Tao solver. Valid
214414d97d3SAlp Dener only for quasi-Newton family of methods.
215414d97d3SAlp Dener
216414d97d3SAlp Dener QN family of methods create their own LMVM matrices and users who wish to
217414d97d3SAlp Dener manipulate this matrix should use TaoGetLMVMMatrix() instead.
218414d97d3SAlp Dener
219414d97d3SAlp Dener Input Parameters:
220414d97d3SAlp Dener + tao - Tao solver context
221414d97d3SAlp Dener - B - LMVM matrix
222414d97d3SAlp Dener
223414d97d3SAlp Dener Level: advanced
224414d97d3SAlp Dener
225db781477SPatrick Sanan .seealso: `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTL`, `TAOBQNKTR`, `MATLMVM`, `TaoGetLMVMMatrix()`
226414d97d3SAlp Dener @*/
TaoSetLMVMMatrix(Tao tao,Mat B)227d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetLMVMMatrix(Tao tao, Mat B)
228d71ae5a4SJacob Faibussowitsch {
229414d97d3SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
230414d97d3SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
231414d97d3SAlp Dener PetscBool flg = PETSC_FALSE;
232414d97d3SAlp Dener
233414d97d3SAlp Dener PetscFunctionBegin;
2349566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)tao, &flg, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, ""));
2353c859ba3SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM Matrix only exists for quasi-Newton algorithms");
2369566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)B, MATLMVM, &flg));
2373c859ba3SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given matrix is not an LMVM matrix");
23848a46eb9SPierre Jolivet if (bqnk->B) PetscCall(MatDestroy(&bqnk->B));
2399566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)B));
240414d97d3SAlp Dener bqnk->B = B;
2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
242414d97d3SAlp Dener }
243