xref: /petsc/src/tao/bound/impls/bqnk/bqnk.c (revision 4f4fdda464867e01e2f7b514a4f52c421e99e383)
1e0ed867bSAlp Dener #include <../src/tao/bound/impls/bqnk/bqnk.h>
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;
9e0ed867bSAlp Dener 
10e0ed867bSAlp Dener   PetscFunctionBegin;
11e0ed867bSAlp Dener   /* Alias the LMVM matrix into the TAO hessian */
12e0ed867bSAlp Dener   if (tao->hessian) {
13e0ed867bSAlp Dener     ierr = MatDestroy(&tao->hessian);CHKERRQ(ierr);
14e0ed867bSAlp Dener   }
15e0ed867bSAlp Dener   if (tao->hessian_pre) {
16e0ed867bSAlp Dener     ierr = MatDestroy(&tao->hessian_pre);CHKERRQ(ierr);
17e0ed867bSAlp Dener   }
18e0ed867bSAlp Dener   ierr = PetscObjectReference((PetscObject)bqnk->B);CHKERRQ(ierr);
19e0ed867bSAlp Dener   tao->hessian = bqnk->B;
20e0ed867bSAlp Dener   ierr = PetscObjectReference((PetscObject)bqnk->B);CHKERRQ(ierr);
21e0ed867bSAlp Dener   tao->hessian_pre = bqnk->B;
22e0ed867bSAlp Dener   /* Update the Hessian with the latest solution */
23e0ed867bSAlp Dener   ierr = MatLMVMUpdate(tao->hessian, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
24e0ed867bSAlp Dener   ierr = MatLMVMResetShift(tao->hessian);CHKERRQ(ierr);
25e0ed867bSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
26*4f4fdda4SAlp Dener   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
27e0ed867bSAlp Dener   if (bnk->active_idx) {
28e0ed867bSAlp Dener     ierr = MatCreateSubMatrixVirtual(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, &bnk->H_inactive);CHKERRQ(ierr);
29*4f4fdda4SAlp Dener     if (bqnk->pc) {
30*4f4fdda4SAlp Dener       ierr = PCLMVMSetIS(bqnk->pc, bnk->inactive_idx);CHKERRQ(ierr);
31*4f4fdda4SAlp Dener     }
32e0ed867bSAlp Dener   } else {
33e0ed867bSAlp Dener     ierr = PetscObjectReference((PetscObject)tao->hessian);CHKERRQ(ierr);
34e0ed867bSAlp Dener     bnk->H_inactive = tao->hessian;
35*4f4fdda4SAlp Dener     if (bqnk->pc) {
36*4f4fdda4SAlp Dener       ierr = PCLMVMClearIS(bqnk->pc);CHKERRQ(ierr);
37e0ed867bSAlp Dener     }
38*4f4fdda4SAlp Dener   }
39*4f4fdda4SAlp Dener   ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
40*4f4fdda4SAlp Dener   ierr = PetscObjectReference((PetscObject)bnk->H_inactive);CHKERRQ(ierr);
41e0ed867bSAlp Dener   bnk->Hpre_inactive = bnk->H_inactive;
42e0ed867bSAlp Dener   PetscFunctionReturn(0);
43e0ed867bSAlp Dener }
44e0ed867bSAlp Dener 
45e0ed867bSAlp Dener static PetscErrorCode TaoBQNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason)
46e0ed867bSAlp Dener {
47e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
48e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
49e0ed867bSAlp Dener   PetscErrorCode ierr;
50e0ed867bSAlp Dener 
51e0ed867bSAlp Dener   PetscFunctionBegin;
52e0ed867bSAlp Dener   ierr = TaoBNKComputeStep(tao, shift, ksp_reason);CHKERRQ(ierr);
53e0ed867bSAlp Dener   if (*ksp_reason < 0) {
54e0ed867bSAlp Dener     /* Krylov solver failed to converge so reset the LMVM matrix */
55e0ed867bSAlp Dener     ierr = MatLMVMReset(bqnk->B, PETSC_FALSE);CHKERRQ(ierr);
56e0ed867bSAlp Dener   }
57e0ed867bSAlp Dener   PetscFunctionReturn(0);
58e0ed867bSAlp Dener }
59e0ed867bSAlp Dener 
60*4f4fdda4SAlp Dener static PetscErrorCode TaoSetUp_BQNK(Tao tao)
61*4f4fdda4SAlp Dener {
62*4f4fdda4SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
63*4f4fdda4SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
64*4f4fdda4SAlp Dener   PetscErrorCode ierr;
65*4f4fdda4SAlp Dener   PetscInt       n, N;
66*4f4fdda4SAlp Dener   PetscBool      is_lmvm, is_sym, is_sr1;
67*4f4fdda4SAlp Dener 
68*4f4fdda4SAlp Dener   PetscFunctionBegin;
69*4f4fdda4SAlp Dener   ierr = TaoSetUp_BNK(tao);CHKERRQ(ierr);
70*4f4fdda4SAlp Dener   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
71*4f4fdda4SAlp Dener   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
72*4f4fdda4SAlp Dener   ierr = MatSetSizes(bqnk->B, n, n, N, N);CHKERRQ(ierr);
73*4f4fdda4SAlp Dener   ierr = MatLMVMAllocate(bqnk->B,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
74*4f4fdda4SAlp Dener   ierr = PetscObjectBaseTypeCompare((PetscObject)bqnk->B, MATLMVM, &is_lmvm);CHKERRQ(ierr);
75*4f4fdda4SAlp Dener   if (!is_lmvm) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Matrix must be an LMVM-type");
76*4f4fdda4SAlp Dener   ierr = MatGetOption(bqnk->B, MAT_SYMMETRIC, &is_sym);CHKERRQ(ierr);
77*4f4fdda4SAlp Dener   if (!is_sym) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric");
78*4f4fdda4SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)bqnk->B, MATLMVMSR1, &is_sr1);CHKERRQ(ierr);
79*4f4fdda4SAlp Dener   ierr = KSPGetPC(tao->ksp, &bqnk->pc);CHKERRQ(ierr);
80*4f4fdda4SAlp Dener   if (is_sr1) {
81*4f4fdda4SAlp Dener     ierr = PCSetType(bqnk->pc, PCNONE);CHKERRQ(ierr);
82*4f4fdda4SAlp Dener     bqnk->pc = NULL;
83*4f4fdda4SAlp Dener     bqnk->no_scale = PETSC_TRUE;
84*4f4fdda4SAlp Dener   } else {
85*4f4fdda4SAlp Dener     ierr = PCSetType(bqnk->pc, PCLMVM);CHKERRQ(ierr);
86*4f4fdda4SAlp Dener     ierr = PCLMVMSetMatLMVM(bqnk->pc, bqnk->B);CHKERRQ(ierr);
87*4f4fdda4SAlp Dener   }
88*4f4fdda4SAlp Dener   if (!bqnk->no_scale) {
89*4f4fdda4SAlp Dener     if (!bqnk->Bscale) {
90*4f4fdda4SAlp Dener       ierr = MatCreateLMVMDiagBrdn(PetscObjectComm((PetscObject)bqnk->B), n, N, &bqnk->Bscale);CHKERRQ(ierr);
91*4f4fdda4SAlp Dener       ierr = MatSetOptionsPrefix(bqnk->Bscale, "tao_bqnk_scale_");CHKERRQ(ierr);
92*4f4fdda4SAlp Dener       ierr = MatLMVMAllocate(bqnk->Bscale, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
93*4f4fdda4SAlp Dener     }
94*4f4fdda4SAlp Dener     ierr = MatLMVMSetJ0(bqnk->B, bqnk->Bscale);CHKERRQ(ierr);
95*4f4fdda4SAlp Dener   }
96*4f4fdda4SAlp Dener   PetscFunctionReturn(0);
97*4f4fdda4SAlp Dener }
98*4f4fdda4SAlp Dener 
99e0ed867bSAlp Dener static PetscErrorCode TaoSetFromOptions_BQNK(PetscOptionItems *PetscOptionsObject, Tao tao)
100e0ed867bSAlp Dener {
101e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
102e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
103e0ed867bSAlp Dener   PetscErrorCode ierr;
104e0ed867bSAlp Dener 
105e0ed867bSAlp Dener   PetscFunctionBegin;
106*4f4fdda4SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Quasi-Newton-Krylov method for bound constrained optimization");CHKERRQ(ierr);
107*4f4fdda4SAlp Dener   ierr = PetscOptionsBool("-tao_bqnk_no_scale","(developer) disable the diagonal Broyden scaling of the BFGS approximation","",bqnk->no_scale,&bqnk->no_scale,NULL);CHKERRQ(ierr);
108*4f4fdda4SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
109e0ed867bSAlp Dener   ierr = TaoSetFromOptions_BNK(PetscOptionsObject, tao);CHKERRQ(ierr);
110e0ed867bSAlp Dener   if (bnk->init_type == BNK_INIT_INTERPOLATION) bnk->init_type = BNK_INIT_DIRECTION;
111e0ed867bSAlp Dener   ierr = MatSetFromOptions(bqnk->B);CHKERRQ(ierr);
112e0ed867bSAlp Dener   PetscFunctionReturn(0);
113e0ed867bSAlp Dener }
114e0ed867bSAlp Dener 
115e0ed867bSAlp Dener static PetscErrorCode TaoView_BQNK(Tao tao, PetscViewer viewer)
116e0ed867bSAlp Dener {
117e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
118e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
119e0ed867bSAlp Dener   PetscErrorCode ierr;
120e0ed867bSAlp Dener   PetscBool      isascii;
121e0ed867bSAlp Dener 
122e0ed867bSAlp Dener   PetscFunctionBegin;
123e0ed867bSAlp Dener   ierr = TaoView_BNK(tao, viewer);CHKERRQ(ierr);
124e0ed867bSAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
125e0ed867bSAlp Dener   if (isascii) {
126e0ed867bSAlp Dener     ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
127e0ed867bSAlp Dener     ierr = MatView(bqnk->B, viewer);CHKERRQ(ierr);
128e0ed867bSAlp Dener     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
129e0ed867bSAlp Dener   }
130e0ed867bSAlp Dener   PetscFunctionReturn(0);
131e0ed867bSAlp Dener }
132e0ed867bSAlp Dener 
133e0ed867bSAlp Dener static PetscErrorCode TaoDestroy_BQNK(Tao tao)
134e0ed867bSAlp Dener {
135e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
136e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
137e0ed867bSAlp Dener   PetscErrorCode ierr;
138e0ed867bSAlp Dener 
139e0ed867bSAlp Dener   PetscFunctionBegin;
140*4f4fdda4SAlp Dener   ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
141cb384e1eSAlp Dener   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
142e0ed867bSAlp Dener   ierr = MatDestroy(&bqnk->B);CHKERRQ(ierr);
143*4f4fdda4SAlp Dener   if (!bqnk->no_scale) {
144*4f4fdda4SAlp Dener     ierr = MatDestroy(&bqnk->Bscale);CHKERRQ(ierr);
145*4f4fdda4SAlp Dener   }
146e0ed867bSAlp Dener   ierr = PetscFree(bnk->ctx);CHKERRQ(ierr);
147e0ed867bSAlp Dener   ierr = TaoDestroy_BNK(tao);CHKERRQ(ierr);
148e0ed867bSAlp Dener   PetscFunctionReturn(0);
149e0ed867bSAlp Dener }
150e0ed867bSAlp Dener 
151e0ed867bSAlp Dener PETSC_INTERN PetscErrorCode TaoCreate_BQNK(Tao tao)
152e0ed867bSAlp Dener {
153e0ed867bSAlp Dener   TAO_BNK        *bnk;
154e0ed867bSAlp Dener   TAO_BQNK       *bqnk;
155e0ed867bSAlp Dener   PetscErrorCode ierr;
156e0ed867bSAlp Dener 
157e0ed867bSAlp Dener   PetscFunctionBegin;
158e0ed867bSAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
159e0ed867bSAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BQNK;
160e0ed867bSAlp Dener   tao->ops->destroy = TaoDestroy_BQNK;
161e0ed867bSAlp Dener   tao->ops->view = TaoView_BQNK;
162*4f4fdda4SAlp Dener   tao->ops->setup = TaoSetUp_BQNK;
163e0ed867bSAlp Dener 
164e0ed867bSAlp Dener   bnk = (TAO_BNK *)tao->data;
165e0ed867bSAlp Dener   bnk->computehessian = TaoBQNKComputeHessian;
166e0ed867bSAlp Dener   bnk->computestep = TaoBQNKComputeStep;
167e0ed867bSAlp Dener   bnk->init_type = BNK_INIT_DIRECTION;
168e0ed867bSAlp Dener 
169e0ed867bSAlp Dener   ierr = PetscNewLog(tao,&bqnk);CHKERRQ(ierr);
170e0ed867bSAlp Dener   bnk->ctx = (void*)bqnk;
171*4f4fdda4SAlp Dener   bqnk->no_scale = PETSC_FALSE;
172e0ed867bSAlp Dener 
173e0ed867bSAlp Dener   ierr = MatCreate(PetscObjectComm((PetscObject)tao), &bqnk->B);CHKERRQ(ierr);
174e0ed867bSAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)bqnk->B, (PetscObject)tao, 1);CHKERRQ(ierr);
175e0ed867bSAlp Dener   ierr = MatSetOptionsPrefix(bqnk->B, "tao_bqnk_");CHKERRQ(ierr);
176e0ed867bSAlp Dener   ierr = MatSetType(bqnk->B, MATLMVMSR1);CHKERRQ(ierr);
177e0ed867bSAlp Dener   PetscFunctionReturn(0);
178e0ed867bSAlp Dener }