Actual source code: bnk.c
1: #include <petsctaolinesearch.h>
2: #include <../src/tao/bound/impls/bnk/bnk.h>
3: #include <petscksp.h>
5: static const char *BNK_INIT[64] = {"constant", "direction", "interpolation"};
6: static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"};
7: static const char *BNK_AS[64] = {"none", "bertsekas"};
9: /* Extracts from the full Hessian the part associated with the current bnk->inactive_idx and set the PCLMVM preconditioner */
11: static PetscErrorCode TaoBNKComputeSubHessian(Tao tao)
12: {
13: TAO_BNK *bnk = (TAO_BNK *)tao->data;
15: PetscFunctionBegin;
16: PetscCall(MatDestroy(&bnk->Hpre_inactive));
17: PetscCall(MatDestroy(&bnk->H_inactive));
18: if (bnk->active_idx) {
19: PetscCall(MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive));
20: if (tao->hessian == tao->hessian_pre) {
21: PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive));
22: bnk->Hpre_inactive = bnk->H_inactive;
23: } else {
24: PetscCall(MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive));
25: }
26: if (bnk->bfgs_pre) PetscCall(PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx));
27: } else {
28: PetscCall(PetscObjectReference((PetscObject)tao->hessian));
29: bnk->H_inactive = tao->hessian;
30: PetscCall(PetscObjectReference((PetscObject)tao->hessian_pre));
31: bnk->Hpre_inactive = tao->hessian_pre;
32: if (bnk->bfgs_pre) PetscCall(PCLMVMClearIS(bnk->bfgs_pre));
33: }
34: PetscFunctionReturn(PETSC_SUCCESS);
35: }
37: /* Initializes the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
39: PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH)
40: {
41: TAO_BNK *bnk = (TAO_BNK *)tao->data;
42: PC pc;
43: PetscReal f_min, ftrial, prered, actred, kappa, sigma, resnorm;
44: PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius;
45: PetscBool is_bfgs, is_jacobi, is_symmetric, sym_set;
46: PetscInt n, N, nDiff;
47: PetscInt i_max = 5;
48: PetscInt j_max = 1;
49: PetscBool kspTR;
51: PetscFunctionBegin;
52: /* Project the current point onto the feasible set */
53: PetscCall(TaoComputeVariableBounds(tao));
54: PetscCall(TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU));
55: if (tao->bounded) PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU));
57: /* Project the initial point onto the feasible region */
58: PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
60: /* Check convergence criteria */
61: PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient));
62: PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
63: PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
64: if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
65: PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
67: /* Test the initial point for convergence */
68: PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
69: PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
70: PetscCheck(!PetscIsInfOrNanReal(bnk->f) && !PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
71: PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
72: PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, 1.0));
73: PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
74: if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
76: /* Reset KSP stopping reason counters */
77: bnk->ksp_atol = 0;
78: bnk->ksp_rtol = 0;
79: bnk->ksp_dtol = 0;
80: bnk->ksp_ctol = 0;
81: bnk->ksp_negc = 0;
82: bnk->ksp_iter = 0;
83: bnk->ksp_othr = 0;
85: /* Reset accepted step type counters */
86: bnk->tot_cg_its = 0;
87: bnk->newt = 0;
88: bnk->bfgs = 0;
89: bnk->sgrad = 0;
90: bnk->grad = 0;
92: /* Initialize the Hessian perturbation */
93: bnk->pert = bnk->sval;
95: /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
96: PetscCall(VecSet(tao->stepdirection, 0.0));
98: /* Allocate the vectors needed for the BFGS approximation */
99: PetscCall(KSPGetPC(tao->ksp, &pc));
100: PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs));
101: PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi));
102: if (is_bfgs) {
103: bnk->bfgs_pre = pc;
104: PetscCall(PCLMVMGetMatLMVM(bnk->bfgs_pre, &bnk->M));
105: PetscCall(VecGetLocalSize(tao->solution, &n));
106: PetscCall(VecGetSize(tao->solution, &N));
107: PetscCall(MatSetSizes(bnk->M, n, n, N, N));
108: PetscCall(MatLMVMAllocate(bnk->M, tao->solution, bnk->unprojected_gradient));
109: PetscCall(MatIsSymmetricKnown(bnk->M, &sym_set, &is_symmetric));
110: PetscCheck(sym_set && is_symmetric, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
111: } else if (is_jacobi) PetscCall(PCJacobiSetUseAbs(pc, PETSC_TRUE));
113: /* Prepare the min/max vectors for safeguarding diagonal scales */
114: PetscCall(VecSet(bnk->Diag_min, bnk->dmin));
115: PetscCall(VecSet(bnk->Diag_max, bnk->dmax));
117: /* Initialize trust-region radius. The initialization is only performed
118: when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
119: *needH = PETSC_TRUE;
120: PetscCall(PetscObjectHasFunction((PetscObject)tao->ksp, "KSPCGSetRadius_C", &kspTR));
121: if (kspTR) {
122: switch (initType) {
123: case BNK_INIT_CONSTANT:
124: /* Use the initial radius specified */
125: tao->trust = tao->trust0;
126: break;
128: case BNK_INIT_INTERPOLATION:
129: /* Use interpolation based on the initial Hessian */
130: max_radius = 0.0;
131: tao->trust = tao->trust0;
132: for (PetscInt j = 0; j < j_max; ++j) {
133: f_min = bnk->f;
134: sigma = 0.0;
136: if (*needH) {
137: /* Compute the Hessian at the new step, and extract the inactive subsystem */
138: PetscCall((*bnk->computehessian)(tao));
139: PetscCall(TaoBNKEstimateActiveSet(tao, BNK_AS_NONE));
140: PetscCall(TaoBNKComputeSubHessian(tao));
141: *needH = PETSC_FALSE;
142: }
144: for (PetscInt i = 0; i < i_max; ++i) {
145: /* Take a steepest descent step and snap it to bounds */
146: PetscCall(VecCopy(tao->solution, bnk->Xold));
147: PetscCall(VecAXPY(tao->solution, -tao->trust / bnk->gnorm, tao->gradient));
148: PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
149: /* Compute the step we actually accepted */
150: PetscCall(VecCopy(tao->solution, bnk->W));
151: PetscCall(VecAXPY(bnk->W, -1.0, bnk->Xold));
152: /* Compute the objective at the trial */
153: PetscCall(TaoComputeObjective(tao, tao->solution, &ftrial));
154: PetscCheck(!PetscIsInfOrNanReal(bnk->f), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
155: PetscCall(VecCopy(bnk->Xold, tao->solution));
156: if (PetscIsInfOrNanReal(ftrial)) {
157: tau = bnk->gamma1_i;
158: } else {
159: if (ftrial < f_min) {
160: f_min = ftrial;
161: sigma = -tao->trust / bnk->gnorm;
162: }
164: /* Compute the predicted and actual reduction */
165: if (bnk->active_idx) {
166: PetscCall(VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive));
167: PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
168: } else {
169: bnk->X_inactive = bnk->W;
170: bnk->inactive_work = bnk->Xwork;
171: }
172: PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work));
173: PetscCall(VecDot(bnk->X_inactive, bnk->inactive_work, &prered));
174: if (bnk->active_idx) {
175: PetscCall(VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive));
176: PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
177: }
178: prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
179: actred = bnk->f - ftrial;
180: if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
181: kappa = 1.0;
182: } else {
183: kappa = actred / prered;
184: }
186: tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
187: tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
188: tau_min = PetscMin(tau_1, tau_2);
189: tau_max = PetscMax(tau_1, tau_2);
191: if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu1_i) {
192: /* Great agreement */
193: max_radius = PetscMax(max_radius, tao->trust);
195: if (tau_max < 1.0) {
196: tau = bnk->gamma3_i;
197: } else if (tau_max > bnk->gamma4_i) {
198: tau = bnk->gamma4_i;
199: } else {
200: tau = tau_max;
201: }
202: } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu2_i) {
203: /* Good agreement */
204: max_radius = PetscMax(max_radius, tao->trust);
206: if (tau_max < bnk->gamma2_i) {
207: tau = bnk->gamma2_i;
208: } else if (tau_max > bnk->gamma3_i) {
209: tau = bnk->gamma3_i;
210: } else {
211: tau = tau_max;
212: }
213: } else {
214: /* Not good agreement */
215: if (tau_min > 1.0) {
216: tau = bnk->gamma2_i;
217: } else if (tau_max < bnk->gamma1_i) {
218: tau = bnk->gamma1_i;
219: } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
220: tau = bnk->gamma1_i;
221: } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
222: tau = tau_1;
223: } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
224: tau = tau_2;
225: } else {
226: tau = tau_max;
227: }
228: }
229: }
230: tao->trust = tau * tao->trust;
231: }
233: if (f_min < bnk->f) {
234: /* We accidentally found a solution better than the initial, so accept it */
235: bnk->f = f_min;
236: PetscCall(VecCopy(tao->solution, bnk->Xold));
237: PetscCall(VecAXPY(tao->solution, sigma, tao->gradient));
238: PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
239: PetscCall(VecCopy(tao->solution, tao->stepdirection));
240: PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold));
241: PetscCall(TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient));
242: PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
243: PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
244: if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
245: /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
246: PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
247: *needH = PETSC_TRUE;
248: /* Test the new step for convergence */
249: PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
250: PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
251: PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
252: PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
253: PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, 1.0));
254: PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
255: if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
256: /* active BNCG recycling early because we have a stepdirection computed */
257: PetscCall(TaoSetRecycleHistory(bnk->bncg, PETSC_TRUE));
258: }
259: }
260: tao->trust = PetscMax(tao->trust, max_radius);
262: /* Ensure that the trust radius is within the limits */
263: tao->trust = PetscMax(tao->trust, bnk->min_radius);
264: tao->trust = PetscMin(tao->trust, bnk->max_radius);
265: break;
267: default:
268: /* Norm of the first direction will initialize radius */
269: tao->trust = 0.0;
270: break;
271: }
272: }
273: PetscFunctionReturn(PETSC_SUCCESS);
274: }
276: /* Computes the exact Hessian and extracts its subHessian */
278: PetscErrorCode TaoBNKComputeHessian(Tao tao)
279: {
280: TAO_BNK *bnk = (TAO_BNK *)tao->data;
282: PetscFunctionBegin;
283: /* Compute the Hessian */
284: PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
285: /* Add a correction to the BFGS preconditioner */
286: if (bnk->M) PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
287: /* Prepare the reduced sub-matrices for the inactive set */
288: PetscCall(TaoBNKComputeSubHessian(tao));
289: PetscFunctionReturn(PETSC_SUCCESS);
290: }
292: /* Routine for estimating the active set */
294: PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
295: {
296: TAO_BNK *bnk = (TAO_BNK *)tao->data;
297: PetscBool hessComputed, diagExists, hadactive;
299: PetscFunctionBegin;
300: hadactive = bnk->active_idx ? PETSC_TRUE : PETSC_FALSE;
301: switch (asType) {
302: case BNK_AS_NONE:
303: PetscCall(ISDestroy(&bnk->inactive_idx));
304: PetscCall(VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx));
305: PetscCall(ISDestroy(&bnk->active_idx));
306: PetscCall(ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx));
307: break;
309: case BNK_AS_BERTSEKAS:
310: /* Compute the trial step vector with which we will estimate the active set at the next iteration */
311: if (bnk->M) {
312: /* If the BFGS matrix is available, we will construct a trial step with it */
313: PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W));
314: } else {
315: hessComputed = diagExists = PETSC_FALSE;
316: if (tao->hessian) PetscCall(MatAssembled(tao->hessian, &hessComputed));
317: if (hessComputed) PetscCall(MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists));
318: if (diagExists) {
319: /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
320: PetscCall(MatGetDiagonal(tao->hessian, bnk->Xwork));
321: PetscCall(VecAbs(bnk->Xwork));
322: PetscCall(VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork));
323: PetscCall(VecReciprocal(bnk->Xwork));
324: PetscCall(VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient));
325: } else {
326: /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
327: PetscCall(VecCopy(bnk->unprojected_gradient, bnk->W));
328: }
329: }
330: PetscCall(VecScale(bnk->W, -1.0));
331: PetscCall(TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol, &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx));
332: break;
334: default:
335: break;
336: }
337: bnk->resetksp = (PetscBool)(bnk->active_idx || hadactive); /* inactive Hessian size may have changed, need to reset operators */
338: PetscFunctionReturn(PETSC_SUCCESS);
339: }
341: /* Routine for bounding the step direction */
343: PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
344: {
345: TAO_BNK *bnk = (TAO_BNK *)tao->data;
347: PetscFunctionBegin;
348: switch (asType) {
349: case BNK_AS_NONE:
350: if (bnk->active_idx) PetscCall(VecISSet(step, bnk->active_idx, 0.0));
351: break;
352: case BNK_AS_BERTSEKAS:
353: PetscCall(TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step));
354: break;
355: default:
356: break;
357: }
358: PetscFunctionReturn(PETSC_SUCCESS);
359: }
361: /* Routine for taking a finite number of BNCG iterations to
362: accelerate Newton convergence.
364: In practice, this approach simply trades off Hessian evaluations
365: for more gradient evaluations.
366: */
368: PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
369: {
370: TAO_BNK *bnk = (TAO_BNK *)tao->data;
372: PetscFunctionBegin;
373: *terminate = PETSC_FALSE;
374: if (bnk->max_cg_its > 0) {
375: /* Copy the current function value (important vectors are already shared) */
376: bnk->bncg_ctx->f = bnk->f;
377: /* Take some small finite number of BNCG iterations */
378: PetscCall(TaoSolve(bnk->bncg));
379: /* Add the number of gradient and function evaluations to the total *
380: * Note: nfuncs are not copied as tao and subsolvers share same TaoTerm */
381: bnk->tot_cg_its += bnk->bncg->niter;
382: /* Extract the BNCG function value out and save it into BNK */
383: bnk->f = bnk->bncg_ctx->f;
384: if (bnk->bncg->reason == TAO_CONVERGED_GATOL || bnk->bncg->reason == TAO_CONVERGED_GRTOL || bnk->bncg->reason == TAO_CONVERGED_GTTOL || bnk->bncg->reason == TAO_CONVERGED_MINF) {
385: *terminate = PETSC_TRUE;
386: } else {
387: PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
388: }
389: }
390: PetscFunctionReturn(PETSC_SUCCESS);
391: }
393: /* Routine for computing the Newton step. */
395: PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
396: {
397: TAO_BNK *bnk = (TAO_BNK *)tao->data;
398: PetscInt bfgsUpdates = 0;
399: PetscInt kspits;
400: PetscBool is_lmvm;
401: PetscBool kspTR;
403: PetscFunctionBegin;
404: /* If there are no inactive variables left, save some computation and return an adjusted zero step
405: that has (l-x) and (u-x) for lower and upper bounded variables. */
406: if (!bnk->inactive_idx) {
407: PetscCall(VecSet(tao->stepdirection, 0.0));
408: PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
409: PetscFunctionReturn(PETSC_SUCCESS);
410: }
412: /* Shift the reduced Hessian matrix */
413: if (shift && bnk->pert > 0) {
414: PetscCall(PetscObjectTypeCompare((PetscObject)tao->hessian, MATLMVM, &is_lmvm));
415: if (is_lmvm) {
416: PetscCall(MatShift(tao->hessian, bnk->pert));
417: } else {
418: PetscCall(MatShift(bnk->H_inactive, bnk->pert));
419: if (bnk->H_inactive != bnk->Hpre_inactive) PetscCall(MatShift(bnk->Hpre_inactive, bnk->pert));
420: }
421: }
423: /* Solve the Newton system of equations */
424: tao->ksp_its = 0;
425: PetscCall(VecSet(tao->stepdirection, 0.0));
426: if (bnk->resetksp) {
427: PetscCall(KSPReset(tao->ksp));
428: PetscCall(KSPResetFromOptions(tao->ksp));
429: bnk->resetksp = PETSC_FALSE;
430: }
431: PetscCall(KSPSetOperators(tao->ksp, bnk->H_inactive, bnk->Hpre_inactive));
432: PetscCall(VecCopy(bnk->unprojected_gradient, bnk->Gwork));
433: if (bnk->active_idx) {
434: PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
435: PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
436: } else {
437: bnk->G_inactive = bnk->unprojected_gradient;
438: bnk->X_inactive = tao->stepdirection;
439: }
440: PetscCall(KSPCGSetRadius(tao->ksp, tao->trust));
441: PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive));
442: PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
443: tao->ksp_its += kspits;
444: tao->ksp_tot_its += kspits;
445: PetscCall(PetscObjectHasFunction((PetscObject)tao->ksp, "KSPCGGetNormD_C", &kspTR));
446: if (kspTR) {
447: PetscCall(KSPCGGetNormD(tao->ksp, &bnk->dnorm));
449: if (0.0 == tao->trust) {
450: /* Radius was uninitialized; use the norm of the direction */
451: if (bnk->dnorm > 0.0) {
452: tao->trust = bnk->dnorm;
454: /* Modify the radius if it is too large or small */
455: tao->trust = PetscMax(tao->trust, bnk->min_radius);
456: tao->trust = PetscMin(tao->trust, bnk->max_radius);
457: } else {
458: /* The direction was bad; set radius to default value and re-solve
459: the trust-region subproblem to get a direction */
460: tao->trust = tao->trust0;
462: /* Modify the radius if it is too large or small */
463: tao->trust = PetscMax(tao->trust, bnk->min_radius);
464: tao->trust = PetscMin(tao->trust, bnk->max_radius);
466: PetscCall(KSPCGSetRadius(tao->ksp, tao->trust));
467: PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive));
468: PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
469: tao->ksp_its += kspits;
470: tao->ksp_tot_its += kspits;
471: PetscCall(KSPCGGetNormD(tao->ksp, &bnk->dnorm));
473: PetscCheck(bnk->dnorm != 0.0, PetscObjectComm((PetscObject)tao), PETSC_ERR_PLIB, "Initial direction zero");
474: }
475: }
476: }
477: /* Restore sub vectors back */
478: if (bnk->active_idx) {
479: PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
480: PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
481: }
482: /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
483: PetscCall(VecScale(tao->stepdirection, -1.0));
484: PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
486: /* Record convergence reasons */
487: PetscCall(KSPGetConvergedReason(tao->ksp, ksp_reason));
488: if (KSP_CONVERGED_ATOL == *ksp_reason) {
489: ++bnk->ksp_atol;
490: } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
491: ++bnk->ksp_rtol;
492: } else if (KSP_CONVERGED_STEP_LENGTH == *ksp_reason) {
493: ++bnk->ksp_ctol;
494: } else if (KSP_CONVERGED_NEG_CURVE == *ksp_reason) {
495: ++bnk->ksp_negc;
496: } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
497: ++bnk->ksp_dtol;
498: } else if (KSP_DIVERGED_ITS == *ksp_reason) {
499: ++bnk->ksp_iter;
500: } else {
501: ++bnk->ksp_othr;
502: }
504: /* Make sure the BFGS preconditioner is healthy */
505: if (bnk->M) {
506: PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
507: if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) {
508: /* Preconditioner is numerically indefinite; reset the approximation. */
509: PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
510: PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
511: }
512: }
513: *step_type = BNK_NEWTON;
514: PetscFunctionReturn(PETSC_SUCCESS);
515: }
517: /* Routine for recomputing the predicted reduction for a given step vector */
519: PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
520: {
521: TAO_BNK *bnk = (TAO_BNK *)tao->data;
523: PetscFunctionBegin;
524: /* Extract subvectors associated with the inactive set */
525: if (bnk->active_idx) {
526: PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
527: PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
528: PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
529: } else {
530: bnk->X_inactive = tao->stepdirection;
531: bnk->inactive_work = bnk->Xwork;
532: bnk->G_inactive = bnk->Gwork;
533: }
534: /* Recompute the predicted decrease based on the quadratic model */
535: PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work));
536: PetscCall(VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive));
537: PetscCall(VecDot(bnk->inactive_work, bnk->X_inactive, prered));
538: /* Restore the sub vectors */
539: if (bnk->active_idx) {
540: PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
541: PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
542: PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
543: }
544: PetscFunctionReturn(PETSC_SUCCESS);
545: }
547: /* Routine for ensuring that the Newton step is a descent direction.
549: The step direction falls back onto BFGS, scaled gradient and gradient steps
550: in the event that the Newton step fails the test.
551: */
553: PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
554: {
555: TAO_BNK *bnk = (TAO_BNK *)tao->data;
556: PetscReal gdx, e_min;
557: PetscInt bfgsUpdates;
559: PetscFunctionBegin;
560: switch (*stepType) {
561: case BNK_NEWTON:
562: PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
563: if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
564: /* Newton step is not descent or direction produced infinity or NaN
565: Update the perturbation for next time */
566: if (bnk->pert <= 0.0) {
567: PetscBool is_gltr;
569: /* Initialize the perturbation */
570: bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
571: PetscCall(PetscObjectTypeCompare((PetscObject)tao->ksp, KSPGLTR, &is_gltr));
572: if (is_gltr) {
573: PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
574: bnk->pert = PetscMax(bnk->pert, -e_min);
575: }
576: } else {
577: /* Increase the perturbation */
578: bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
579: }
581: if (!bnk->M) {
582: /* We don't have the bfgs matrix around and updated
583: Must use gradient direction in this case */
584: PetscCall(VecCopy(tao->gradient, tao->stepdirection));
585: *stepType = BNK_GRADIENT;
586: } else {
587: /* Attempt to use the BFGS direction */
588: PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
590: /* Check for success (descent direction)
591: NOTE: Negative gdx here means not a descent direction because
592: the fall-back step is missing a negative sign. */
593: PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
594: if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
595: /* BFGS direction is not descent or direction produced not a number
596: We can assert bfgsUpdates > 1 in this case because
597: the first solve produces the scaled gradient direction,
598: which is guaranteed to be descent */
600: /* Use steepest descent direction (scaled) */
601: PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
602: PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
603: PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
605: *stepType = BNK_SCALED_GRADIENT;
606: } else {
607: PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
608: if (1 == bfgsUpdates) {
609: /* The first BFGS direction is always the scaled gradient */
610: *stepType = BNK_SCALED_GRADIENT;
611: } else {
612: *stepType = BNK_BFGS;
613: }
614: }
615: }
616: /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
617: PetscCall(VecScale(tao->stepdirection, -1.0));
618: PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
619: } else {
620: /* Computed Newton step is descent */
621: switch (ksp_reason) {
622: case KSP_DIVERGED_NANORINF:
623: case KSP_DIVERGED_BREAKDOWN:
624: case KSP_DIVERGED_INDEFINITE_MAT:
625: case KSP_DIVERGED_INDEFINITE_PC:
626: case KSP_CONVERGED_NEG_CURVE:
627: /* Matrix or preconditioner is indefinite; increase perturbation */
628: if (bnk->pert <= 0.0) {
629: PetscBool is_gltr;
631: /* Initialize the perturbation */
632: bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
633: PetscCall(PetscObjectTypeCompare((PetscObject)tao->ksp, KSPGLTR, &is_gltr));
634: if (is_gltr) {
635: PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
636: bnk->pert = PetscMax(bnk->pert, -e_min);
637: }
638: } else {
639: /* Increase the perturbation */
640: bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
641: }
642: break;
644: default:
645: /* Newton step computation is good; decrease perturbation */
646: bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
647: if (bnk->pert < bnk->pmin) bnk->pert = 0.0;
648: break;
649: }
650: *stepType = BNK_NEWTON;
651: }
652: break;
654: case BNK_BFGS:
655: /* Check for success (descent direction) */
656: PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
657: if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) {
658: /* Step is not descent or solve was not successful
659: Use steepest descent direction (scaled) */
660: PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
661: PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
662: PetscCall(MatSolve(bnk->M, tao->gradient, tao->stepdirection));
663: PetscCall(VecScale(tao->stepdirection, -1.0));
664: PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
665: *stepType = BNK_SCALED_GRADIENT;
666: } else {
667: *stepType = BNK_BFGS;
668: }
669: break;
671: case BNK_SCALED_GRADIENT:
672: break;
674: default:
675: break;
676: }
677: PetscFunctionReturn(PETSC_SUCCESS);
678: }
680: /* Routine for performing a bound-projected More-Thuente line search.
682: Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
683: Newton step does not produce a valid step length.
684: */
686: PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
687: {
688: TAO_BNK *bnk = (TAO_BNK *)tao->data;
689: TaoLineSearchConvergedReason ls_reason;
690: PetscReal e_min, gdx;
691: PetscInt bfgsUpdates;
693: PetscFunctionBegin;
694: /* Perform the linesearch */
695: PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason));
696: PetscCall(TaoAddLineSearchCounts(tao));
698: while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) {
699: /* Linesearch failed, revert solution */
700: bnk->f = bnk->fold;
701: PetscCall(VecCopy(bnk->Xold, tao->solution));
702: PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient));
704: switch (*stepType) {
705: case BNK_NEWTON:
706: /* Failed to obtain acceptable iterate with Newton step
707: Update the perturbation for next time */
708: if (bnk->pert <= 0.0) {
709: PetscBool is_gltr;
711: /* Initialize the perturbation */
712: bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
713: PetscCall(PetscObjectTypeCompare((PetscObject)tao->ksp, KSPGLTR, &is_gltr));
714: if (is_gltr) {
715: PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
716: bnk->pert = PetscMax(bnk->pert, -e_min);
717: }
718: } else {
719: /* Increase the perturbation */
720: bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
721: }
723: if (!bnk->M) {
724: /* We don't have the bfgs matrix around and being updated
725: Must use gradient direction in this case */
726: PetscCall(VecCopy(bnk->unprojected_gradient, tao->stepdirection));
727: *stepType = BNK_GRADIENT;
728: } else {
729: /* Attempt to use the BFGS direction */
730: PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
731: /* Check for success (descent direction)
732: NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
733: PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
734: if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
735: /* BFGS direction is not descent or direction produced not a number
736: We can assert bfgsUpdates > 1 in this case
737: Use steepest descent direction (scaled) */
738: PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
739: PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
740: PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
742: bfgsUpdates = 1;
743: *stepType = BNK_SCALED_GRADIENT;
744: } else {
745: PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
746: if (1 == bfgsUpdates) {
747: /* The first BFGS direction is always the scaled gradient */
748: *stepType = BNK_SCALED_GRADIENT;
749: } else {
750: *stepType = BNK_BFGS;
751: }
752: }
753: }
754: break;
756: case BNK_BFGS:
757: /* Can only enter if pc_type == BNK_PC_BFGS
758: Failed to obtain acceptable iterate with BFGS step
759: Attempt to use the scaled gradient direction */
760: PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
761: PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
762: PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
764: bfgsUpdates = 1;
765: *stepType = BNK_SCALED_GRADIENT;
766: break;
767: }
768: /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
769: PetscCall(VecScale(tao->stepdirection, -1.0));
770: PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
772: /* Perform one last line search with the fall-back step */
773: PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason));
774: PetscCall(TaoAddLineSearchCounts(tao));
775: }
776: *reason = ls_reason;
777: PetscFunctionReturn(PETSC_SUCCESS);
778: }
780: /* Routine for updating the trust radius.
782: Function features three different update methods:
783: 1) Line-search step length based
784: 2) Predicted decrease on the CG quadratic model
785: 3) Interpolation
786: */
788: PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
789: {
790: TAO_BNK *bnk = (TAO_BNK *)tao->data;
792: PetscReal step, kappa;
793: PetscReal gdx, tau_1, tau_2, tau_min, tau_max;
795: PetscFunctionBegin;
796: /* Update trust region radius */
797: *accept = PETSC_FALSE;
798: switch (updateType) {
799: case BNK_UPDATE_STEP:
800: *accept = PETSC_TRUE; /* always accept here because line search succeeded */
801: if (stepType == BNK_NEWTON) {
802: PetscCall(TaoLineSearchGetStepLength(tao->linesearch, &step));
803: if (step < bnk->nu1) {
804: /* Very bad step taken; reduce radius */
805: tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
806: } else if (step < bnk->nu2) {
807: /* Reasonably bad step taken; reduce radius */
808: tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
809: } else if (step < bnk->nu3) {
810: /* Reasonable step was taken; leave radius alone */
811: if (bnk->omega3 < 1.0) {
812: tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
813: } else if (bnk->omega3 > 1.0) {
814: tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
815: }
816: } else if (step < bnk->nu4) {
817: /* Full step taken; increase the radius */
818: tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
819: } else {
820: /* More than full step taken; increase the radius */
821: tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
822: }
823: } else {
824: /* Newton step was not good; reduce the radius */
825: tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
826: }
827: break;
829: case BNK_UPDATE_REDUCTION:
830: if (stepType == BNK_NEWTON) {
831: if ((prered < 0.0) || PetscIsInfOrNanReal(prered)) {
832: /* The predicted reduction has the wrong sign. This cannot
833: happen in infinite precision arithmetic. Step should
834: be rejected! */
835: tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
836: } else {
837: if (PetscIsInfOrNanReal(actred)) {
838: tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
839: } else {
840: if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f)) * bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f)) * bnk->epsilon)) {
841: kappa = 1.0;
842: } else {
843: kappa = actred / prered;
844: }
845: /* Accept or reject the step and update radius */
846: if (kappa < bnk->eta1) {
847: /* Reject the step */
848: tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
849: } else {
850: /* Accept the step */
851: *accept = PETSC_TRUE;
852: /* Update the trust region radius only if the computed step is at the trust radius boundary */
853: if (bnk->dnorm == tao->trust) {
854: if (kappa < bnk->eta2) {
855: /* Marginal bad step */
856: tao->trust = bnk->alpha2 * tao->trust;
857: } else if (kappa < bnk->eta3) {
858: /* Reasonable step */
859: tao->trust = bnk->alpha3 * tao->trust;
860: } else if (kappa < bnk->eta4) {
861: /* Good step */
862: tao->trust = bnk->alpha4 * tao->trust;
863: } else {
864: /* Very good step */
865: tao->trust = bnk->alpha5 * tao->trust;
866: }
867: }
868: }
869: }
870: }
871: } else {
872: /* Newton step was not good; reduce the radius */
873: tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
874: }
875: break;
877: default:
878: if (stepType == BNK_NEWTON) {
879: if (prered < 0.0) {
880: /* The predicted reduction has the wrong sign. This cannot */
881: /* happen in infinite precision arithmetic. Step should */
882: /* be rejected! */
883: tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
884: } else {
885: if (PetscIsInfOrNanReal(actred)) {
886: tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
887: } else {
888: if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
889: kappa = 1.0;
890: } else {
891: kappa = actred / prered;
892: }
894: PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
895: tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
896: tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
897: tau_min = PetscMin(tau_1, tau_2);
898: tau_max = PetscMax(tau_1, tau_2);
900: if (kappa >= 1.0 - bnk->mu1) {
901: /* Great agreement */
902: *accept = PETSC_TRUE;
903: if (tau_max < 1.0) {
904: tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
905: } else if (tau_max > bnk->gamma4) {
906: tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
907: } else {
908: tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
909: }
910: } else if (kappa >= 1.0 - bnk->mu2) {
911: /* Good agreement */
912: *accept = PETSC_TRUE;
913: if (tau_max < bnk->gamma2) {
914: tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
915: } else if (tau_max > bnk->gamma3) {
916: tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
917: } else if (tau_max < 1.0) {
918: tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
919: } else {
920: tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
921: }
922: } else {
923: /* Not good agreement */
924: if (tau_min > 1.0) {
925: tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
926: } else if (tau_max < bnk->gamma1) {
927: tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
928: } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
929: tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
930: } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
931: tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
932: } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
933: tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
934: } else {
935: tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
936: }
937: }
938: }
939: }
940: } else {
941: /* Newton step was not good; reduce the radius */
942: tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
943: }
944: break;
945: }
946: /* Make sure the radius does not violate min and max settings */
947: tao->trust = PetscMin(tao->trust, bnk->max_radius);
948: tao->trust = PetscMax(tao->trust, bnk->min_radius);
949: PetscFunctionReturn(PETSC_SUCCESS);
950: }
952: PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
953: {
954: TAO_BNK *bnk = (TAO_BNK *)tao->data;
956: PetscFunctionBegin;
957: switch (stepType) {
958: case BNK_NEWTON:
959: ++bnk->newt;
960: break;
961: case BNK_BFGS:
962: ++bnk->bfgs;
963: break;
964: case BNK_SCALED_GRADIENT:
965: ++bnk->sgrad;
966: break;
967: case BNK_GRADIENT:
968: ++bnk->grad;
969: break;
970: default:
971: break;
972: }
973: PetscFunctionReturn(PETSC_SUCCESS);
974: }
976: PetscErrorCode TaoSetUp_BNK(Tao tao)
977: {
978: TAO_BNK *bnk = (TAO_BNK *)tao->data;
980: PetscFunctionBegin;
981: if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
982: if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
983: if (!bnk->W) PetscCall(VecDuplicate(tao->solution, &bnk->W));
984: if (!bnk->Xold) PetscCall(VecDuplicate(tao->solution, &bnk->Xold));
985: if (!bnk->Gold) PetscCall(VecDuplicate(tao->solution, &bnk->Gold));
986: if (!bnk->Xwork) PetscCall(VecDuplicate(tao->solution, &bnk->Xwork));
987: if (!bnk->Gwork) PetscCall(VecDuplicate(tao->solution, &bnk->Gwork));
988: if (!bnk->unprojected_gradient) PetscCall(VecDuplicate(tao->solution, &bnk->unprojected_gradient));
989: if (!bnk->unprojected_gradient_old) PetscCall(VecDuplicate(tao->solution, &bnk->unprojected_gradient_old));
990: if (!bnk->Diag_min) PetscCall(VecDuplicate(tao->solution, &bnk->Diag_min));
991: if (!bnk->Diag_max) PetscCall(VecDuplicate(tao->solution, &bnk->Diag_max));
992: PetscCall(TaoSetSolution(bnk->bncg, tao->solution));
993: if (bnk->max_cg_its > 0) {
994: /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
995: bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
996: PetscCall(PetscObjectReference((PetscObject)bnk->unprojected_gradient_old));
997: PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old));
998: bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
999: PetscCall(PetscObjectReference((PetscObject)bnk->unprojected_gradient));
1000: PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient));
1001: bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
1002: PetscCall(PetscObjectReference((PetscObject)bnk->Gold));
1003: PetscCall(VecDestroy(&bnk->bncg_ctx->G_old));
1004: bnk->bncg_ctx->G_old = bnk->Gold;
1005: PetscCall(PetscObjectReference((PetscObject)tao->gradient));
1006: PetscCall(VecDestroy(&bnk->bncg->gradient));
1007: bnk->bncg->gradient = tao->gradient;
1008: PetscCall(PetscObjectReference((PetscObject)tao->stepdirection));
1009: PetscCall(VecDestroy(&bnk->bncg->stepdirection));
1010: bnk->bncg->stepdirection = tao->stepdirection;
1011: /* Copy over some settings from BNK into BNCG */
1012: PetscCall(TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its));
1013: PetscCall(TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol));
1014: PetscCall(TaoSetFunctionLowerBound(bnk->bncg, tao->fmin));
1015: PetscCall(TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP));
1016: {
1017: TaoTerm term;
1018: Vec params;
1019: PetscReal scale;
1020: Mat map;
1022: // Note: tao->objective_term.term and bnk->bncg->objective_term.term will point to same address
1023: PetscCall(TaoGetTerm(tao, &scale, &term, ¶ms, &map));
1024: PetscCall(TaoAddTerm(bnk->bncg, NULL, scale, term, params, map));
1025: }
1026: PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)bnk->bncg));
1027: }
1028: bnk->X_inactive = NULL;
1029: bnk->G_inactive = NULL;
1030: bnk->inactive_work = NULL;
1031: bnk->active_work = NULL;
1032: bnk->inactive_idx = NULL;
1033: bnk->active_idx = NULL;
1034: bnk->active_lower = NULL;
1035: bnk->active_upper = NULL;
1036: bnk->active_fixed = NULL;
1037: bnk->M = NULL;
1038: bnk->H_inactive = NULL;
1039: bnk->Hpre_inactive = NULL;
1040: PetscFunctionReturn(PETSC_SUCCESS);
1041: }
1043: PetscErrorCode TaoDestroy_BNK(Tao tao)
1044: {
1045: TAO_BNK *bnk = (TAO_BNK *)tao->data;
1047: PetscFunctionBegin;
1048: PetscCall(VecDestroy(&bnk->W));
1049: PetscCall(VecDestroy(&bnk->Xold));
1050: PetscCall(VecDestroy(&bnk->Gold));
1051: PetscCall(VecDestroy(&bnk->Xwork));
1052: PetscCall(VecDestroy(&bnk->Gwork));
1053: PetscCall(VecDestroy(&bnk->unprojected_gradient));
1054: PetscCall(VecDestroy(&bnk->unprojected_gradient_old));
1055: PetscCall(VecDestroy(&bnk->Diag_min));
1056: PetscCall(VecDestroy(&bnk->Diag_max));
1057: PetscCall(ISDestroy(&bnk->active_lower));
1058: PetscCall(ISDestroy(&bnk->active_upper));
1059: PetscCall(ISDestroy(&bnk->active_fixed));
1060: PetscCall(ISDestroy(&bnk->active_idx));
1061: PetscCall(ISDestroy(&bnk->inactive_idx));
1062: PetscCall(MatDestroy(&bnk->Hpre_inactive));
1063: PetscCall(MatDestroy(&bnk->H_inactive));
1064: PetscCall(TaoDestroy(&bnk->bncg));
1065: PetscCall(KSPDestroy(&tao->ksp));
1066: PetscCall(PetscFree(tao->data));
1067: PetscFunctionReturn(PETSC_SUCCESS);
1068: }
1070: PetscErrorCode TaoSetFromOptions_BNK(Tao tao, PetscOptionItems PetscOptionsObject)
1071: {
1072: TAO_BNK *bnk = (TAO_BNK *)tao->data;
1074: PetscFunctionBegin;
1075: PetscOptionsHeadBegin(PetscOptionsObject, "Newton-Krylov method for bound constrained optimization");
1076: PetscCall(PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, NULL));
1077: PetscCall(PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, NULL));
1078: PetscCall(PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, NULL));
1079: PetscCall(PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval, NULL));
1080: PetscCall(PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin, NULL));
1081: PetscCall(PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax, NULL));
1082: PetscCall(PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac, NULL));
1083: PetscCall(PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin, NULL));
1084: PetscCall(PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax, NULL));
1085: PetscCall(PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac, NULL));
1086: PetscCall(PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac, NULL));
1087: PetscCall(PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac, NULL));
1088: PetscCall(PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac, NULL));
1089: PetscCall(PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1, NULL));
1090: PetscCall(PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2, NULL));
1091: PetscCall(PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3, NULL));
1092: PetscCall(PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4, NULL));
1093: PetscCall(PetscOptionsReal("-tao_bnk_alpha1", "(developer) radius reduction factor for rejected step (-tao_bnk_update_type reduction)", "", bnk->alpha1, &bnk->alpha1, NULL));
1094: PetscCall(PetscOptionsReal("-tao_bnk_alpha2", "(developer) radius reduction factor for marginally accepted bad step (-tao_bnk_update_type reduction)", "", bnk->alpha2, &bnk->alpha2, NULL));
1095: PetscCall(PetscOptionsReal("-tao_bnk_alpha3", "(developer) radius increase factor for reasonable accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha3, &bnk->alpha3, NULL));
1096: PetscCall(PetscOptionsReal("-tao_bnk_alpha4", "(developer) radius increase factor for good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha4, &bnk->alpha4, NULL));
1097: PetscCall(PetscOptionsReal("-tao_bnk_alpha5", "(developer) radius increase factor for very good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha5, &bnk->alpha5, NULL));
1098: PetscCall(PetscOptionsReal("-tao_bnk_nu1", "(developer) threshold for small line-search step length (-tao_bnk_update_type step)", "", bnk->nu1, &bnk->nu1, NULL));
1099: PetscCall(PetscOptionsReal("-tao_bnk_nu2", "(developer) threshold for reasonable line-search step length (-tao_bnk_update_type step)", "", bnk->nu2, &bnk->nu2, NULL));
1100: PetscCall(PetscOptionsReal("-tao_bnk_nu3", "(developer) threshold for large line-search step length (-tao_bnk_update_type step)", "", bnk->nu3, &bnk->nu3, NULL));
1101: PetscCall(PetscOptionsReal("-tao_bnk_nu4", "(developer) threshold for very large line-search step length (-tao_bnk_update_type step)", "", bnk->nu4, &bnk->nu4, NULL));
1102: PetscCall(PetscOptionsReal("-tao_bnk_omega1", "(developer) radius reduction factor for very small line-search step length (-tao_bnk_update_type step)", "", bnk->omega1, &bnk->omega1, NULL));
1103: PetscCall(PetscOptionsReal("-tao_bnk_omega2", "(developer) radius reduction factor for small line-search step length (-tao_bnk_update_type step)", "", bnk->omega2, &bnk->omega2, NULL));
1104: PetscCall(PetscOptionsReal("-tao_bnk_omega3", "(developer) radius factor for decent line-search step length (-tao_bnk_update_type step)", "", bnk->omega3, &bnk->omega3, NULL));
1105: PetscCall(PetscOptionsReal("-tao_bnk_omega4", "(developer) radius increase factor for large line-search step length (-tao_bnk_update_type step)", "", bnk->omega4, &bnk->omega4, NULL));
1106: PetscCall(PetscOptionsReal("-tao_bnk_omega5", "(developer) radius increase factor for very large line-search step length (-tao_bnk_update_type step)", "", bnk->omega5, &bnk->omega5, NULL));
1107: PetscCall(PetscOptionsReal("-tao_bnk_mu1_i", "(developer) threshold for accepting very good step (-tao_bnk_init_type interpolation)", "", bnk->mu1_i, &bnk->mu1_i, NULL));
1108: PetscCall(PetscOptionsReal("-tao_bnk_mu2_i", "(developer) threshold for accepting good step (-tao_bnk_init_type interpolation)", "", bnk->mu2_i, &bnk->mu2_i, NULL));
1109: PetscCall(PetscOptionsReal("-tao_bnk_gamma1_i", "(developer) radius reduction factor for rejected very bad step (-tao_bnk_init_type interpolation)", "", bnk->gamma1_i, &bnk->gamma1_i, NULL));
1110: PetscCall(PetscOptionsReal("-tao_bnk_gamma2_i", "(developer) radius reduction factor for rejected bad step (-tao_bnk_init_type interpolation)", "", bnk->gamma2_i, &bnk->gamma2_i, NULL));
1111: PetscCall(PetscOptionsReal("-tao_bnk_gamma3_i", "(developer) radius increase factor for accepted good step (-tao_bnk_init_type interpolation)", "", bnk->gamma3_i, &bnk->gamma3_i, NULL));
1112: PetscCall(PetscOptionsReal("-tao_bnk_gamma4_i", "(developer) radius increase factor for accepted very good step (-tao_bnk_init_type interpolation)", "", bnk->gamma4_i, &bnk->gamma4_i, NULL));
1113: PetscCall(PetscOptionsReal("-tao_bnk_theta_i", "(developer) trust region interpolation factor (-tao_bnk_init_type interpolation)", "", bnk->theta_i, &bnk->theta_i, NULL));
1114: PetscCall(PetscOptionsReal("-tao_bnk_mu1", "(developer) threshold for accepting very good step (-tao_bnk_update_type interpolation)", "", bnk->mu1, &bnk->mu1, NULL));
1115: PetscCall(PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2, NULL));
1116: PetscCall(PetscOptionsReal("-tao_bnk_gamma1", "(developer) radius reduction factor for rejected very bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma1, &bnk->gamma1, NULL));
1117: PetscCall(PetscOptionsReal("-tao_bnk_gamma2", "(developer) radius reduction factor for rejected bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma2, &bnk->gamma2, NULL));
1118: PetscCall(PetscOptionsReal("-tao_bnk_gamma3", "(developer) radius increase factor for accepted good step (-tao_bnk_update_type interpolation)", "", bnk->gamma3, &bnk->gamma3, NULL));
1119: PetscCall(PetscOptionsReal("-tao_bnk_gamma4", "(developer) radius increase factor for accepted very good step (-tao_bnk_update_type interpolation)", "", bnk->gamma4, &bnk->gamma4, NULL));
1120: PetscCall(PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta, NULL));
1121: PetscCall(PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius, NULL));
1122: PetscCall(PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius, NULL));
1123: PetscCall(PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon, NULL));
1124: PetscCall(PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol, NULL));
1125: PetscCall(PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step, NULL));
1126: PetscCall(PetscOptionsInt("-tao_bnk_max_cg_its", "number of BNCG iterations to take for each Newton step", "", bnk->max_cg_its, &bnk->max_cg_its, NULL));
1127: PetscOptionsHeadEnd();
1129: PetscCall(TaoSetOptionsPrefix(bnk->bncg, ((PetscObject)tao)->prefix));
1130: PetscCall(TaoAppendOptionsPrefix(bnk->bncg, "tao_bnk_cg_"));
1131: PetscCall(TaoSetFromOptions(bnk->bncg));
1133: PetscCall(KSPSetOptionsPrefix(tao->ksp, ((PetscObject)tao)->prefix));
1134: PetscCall(KSPAppendOptionsPrefix(tao->ksp, "tao_bnk_"));
1135: PetscCall(KSPSetFromOptions(tao->ksp));
1136: PetscFunctionReturn(PETSC_SUCCESS);
1137: }
1139: PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1140: {
1141: TAO_BNK *bnk = (TAO_BNK *)tao->data;
1142: PetscInt nrejects;
1143: PetscBool isascii;
1145: PetscFunctionBegin;
1146: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1147: if (isascii) {
1148: PetscCall(PetscViewerASCIIPushTab(viewer));
1149: PetscCall(TaoView(bnk->bncg, viewer));
1150: if (bnk->M) {
1151: PetscCall(MatLMVMGetRejectCount(bnk->M, &nrejects));
1152: PetscCall(PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %" PetscInt_FMT "\n", nrejects));
1153: }
1154: PetscCall(PetscViewerASCIIPrintf(viewer, "CG steps: %" PetscInt_FMT "\n", bnk->tot_cg_its));
1155: PetscCall(PetscViewerASCIIPrintf(viewer, "Newton steps: %" PetscInt_FMT "\n", bnk->newt));
1156: if (bnk->M) PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS steps: %" PetscInt_FMT "\n", bnk->bfgs));
1157: PetscCall(PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %" PetscInt_FMT "\n", bnk->sgrad));
1158: PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", bnk->grad));
1159: PetscCall(PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n"));
1160: PetscCall(PetscViewerASCIIPrintf(viewer, " atol: %" PetscInt_FMT "\n", bnk->ksp_atol));
1161: PetscCall(PetscViewerASCIIPrintf(viewer, " rtol: %" PetscInt_FMT "\n", bnk->ksp_rtol));
1162: PetscCall(PetscViewerASCIIPrintf(viewer, " ctol: %" PetscInt_FMT "\n", bnk->ksp_ctol));
1163: PetscCall(PetscViewerASCIIPrintf(viewer, " negc: %" PetscInt_FMT "\n", bnk->ksp_negc));
1164: PetscCall(PetscViewerASCIIPrintf(viewer, " dtol: %" PetscInt_FMT "\n", bnk->ksp_dtol));
1165: PetscCall(PetscViewerASCIIPrintf(viewer, " iter: %" PetscInt_FMT "\n", bnk->ksp_iter));
1166: PetscCall(PetscViewerASCIIPrintf(viewer, " othr: %" PetscInt_FMT "\n", bnk->ksp_othr));
1167: PetscCall(PetscViewerASCIIPopTab(viewer));
1168: }
1169: PetscFunctionReturn(PETSC_SUCCESS);
1170: }
1172: /*MC
1173: TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
1174: At each iteration, the BNK methods solve the symmetric
1175: system of equations to obtain the step direction $d_k$:
1176: $ H_k d_k = -g_k $
1177: for free variables only. The step can be globalized either through
1178: trust-region methods, or a line search, or a heuristic mixture of both.
1180: Options Database Keys:
1181: + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
1182: . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
1183: . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
1184: . -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
1185: . -tao_bnk_as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
1186: . -tao_bnk_as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
1187: . -tao_bnk_sval - (developer) Hessian perturbation starting value
1188: . -tao_bnk_imin - (developer) minimum initial Hessian perturbation
1189: . -tao_bnk_imax - (developer) maximum initial Hessian perturbation
1190: . -tao_bnk_pmin - (developer) minimum Hessian perturbation
1191: . -tao_bnk_pmax - (developer) aximum Hessian perturbation
1192: . -tao_bnk_pgfac - (developer) Hessian perturbation growth factor
1193: . -tao_bnk_psfac - (developer) Hessian perturbation shrink factor
1194: . -tao_bnk_imfac - (developer) initial merit factor for Hessian perturbation
1195: . -tao_bnk_pmgfac - (developer) merit growth factor for Hessian perturbation
1196: . -tao_bnk_pmsfac - (developer) merit shrink factor for Hessian perturbation
1197: . -tao_bnk_eta1 - (developer) threshold for rejecting step (-update_type reduction)
1198: . -tao_bnk_eta2 - (developer) threshold for accepting marginal step (-update_type reduction)
1199: . -tao_bnk_eta3 - (developer) threshold for accepting reasonable step (-update_type reduction)
1200: . -tao_bnk_eta4 - (developer) threshold for accepting good step (-update_type reduction)
1201: . -tao_bnk_alpha1 - (developer) radius reduction factor for rejected step (-update_type reduction)
1202: . -tao_bnk_alpha2 - (developer) radius reduction factor for marginally accepted bad step (-update_type reduction)
1203: . -tao_bnk_alpha3 - (developer) radius increase factor for reasonable accepted step (-update_type reduction)
1204: . -tao_bnk_alpha4 - (developer) radius increase factor for good accepted step (-update_type reduction)
1205: . -tao_bnk_alpha5 - (developer) radius increase factor for very good accepted step (-update_type reduction)
1206: . -tao_bnk_epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-update_type reduction)
1207: . -tao_bnk_mu1 - (developer) threshold for accepting very good step (-update_type interpolation)
1208: . -tao_bnk_mu2 - (developer) threshold for accepting good step (-update_type interpolation)
1209: . -tao_bnk_gamma1 - (developer) radius reduction factor for rejected very bad step (-update_type interpolation)
1210: . -tao_bnk_gamma2 - (developer) radius reduction factor for rejected bad step (-update_type interpolation)
1211: . -tao_bnk_gamma3 - (developer) radius increase factor for accepted good step (-update_type interpolation)
1212: . -tao_bnk_gamma4 - (developer) radius increase factor for accepted very good step (-update_type interpolation)
1213: . -tao_bnk_theta - (developer) trust region interpolation factor (-update_type interpolation)
1214: . -tao_bnk_nu1 - (developer) threshold for small line-search step length (-update_type step)
1215: . -tao_bnk_nu2 - (developer) threshold for reasonable line-search step length (-update_type step)
1216: . -tao_bnk_nu3 - (developer) threshold for large line-search step length (-update_type step)
1217: . -tao_bnk_nu4 - (developer) threshold for very large line-search step length (-update_type step)
1218: . -tao_bnk_omega1 - (developer) radius reduction factor for very small line-search step length (-update_type step)
1219: . -tao_bnk_omega2 - (developer) radius reduction factor for small line-search step length (-update_type step)
1220: . -tao_bnk_omega3 - (developer) radius factor for decent line-search step length (-update_type step)
1221: . -tao_bnk_omega4 - (developer) radius increase factor for large line-search step length (-update_type step)
1222: . -tao_bnk_omega5 - (developer) radius increase factor for very large line-search step length (-update_type step)
1223: . -tao_bnk_mu1_i - (developer) threshold for accepting very good step (-init_type interpolation)
1224: . -tao_bnk_mu2_i - (developer) threshold for accepting good step (-init_type interpolation)
1225: . -tao_bnk_gamma1_i - (developer) radius reduction factor for rejected very bad step (-init_type interpolation)
1226: . -tao_bnk_gamma2_i - (developer) radius reduction factor for rejected bad step (-init_type interpolation)
1227: . -tao_bnk_gamma3_i - (developer) radius increase factor for accepted good step (-init_type interpolation)
1228: . -tao_bnk_gamma4_i - (developer) radius increase factor for accepted very good step (-init_type interpolation)
1229: - -tao_bnk_theta_i - (developer) trust region interpolation factor (-init_type interpolation)
1231: Level: beginner
1233: The various algorithmic factors can only be supplied via the options database
1235: .seealso: `Tao`, `TAONLS`, `TAONTL`, `TAONM`, `TaoType`, `TaoCreate()`
1236: M*/
1238: PetscErrorCode TaoCreate_BNK(Tao tao)
1239: {
1240: TAO_BNK *bnk;
1241: PC pc;
1243: PetscFunctionBegin;
1244: PetscCall(PetscNew(&bnk));
1246: tao->ops->setup = TaoSetUp_BNK;
1247: tao->ops->view = TaoView_BNK;
1248: tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1249: tao->ops->destroy = TaoDestroy_BNK;
1250: tao->uses_gradient = PETSC_TRUE;
1251: tao->uses_hessian_matrices = PETSC_TRUE;
1253: /* Override default settings (unless already changed) */
1254: PetscCall(TaoParametersInitialize(tao));
1255: PetscObjectParameterSetDefault(tao, max_it, 50);
1256: PetscObjectParameterSetDefault(tao, trust0, 100.0);
1258: tao->data = (void *)bnk;
1260: /* Hessian shifting parameters */
1261: bnk->computehessian = TaoBNKComputeHessian;
1262: bnk->computestep = TaoBNKComputeStep;
1264: bnk->sval = 0.0;
1265: bnk->imin = 1.0e-4;
1266: bnk->imax = 1.0e+2;
1267: bnk->imfac = 1.0e-1;
1269: bnk->pmin = 1.0e-12;
1270: bnk->pmax = 1.0e+2;
1271: bnk->pgfac = 1.0e+1;
1272: bnk->psfac = 4.0e-1;
1273: bnk->pmgfac = 1.0e-1;
1274: bnk->pmsfac = 1.0e-1;
1276: /* Default values for trust-region radius update based on steplength */
1277: bnk->nu1 = 0.25;
1278: bnk->nu2 = 0.50;
1279: bnk->nu3 = 1.00;
1280: bnk->nu4 = 1.25;
1282: bnk->omega1 = 0.25;
1283: bnk->omega2 = 0.50;
1284: bnk->omega3 = 1.00;
1285: bnk->omega4 = 2.00;
1286: bnk->omega5 = 4.00;
1288: /* Default values for trust-region radius update based on reduction */
1289: bnk->eta1 = 1.0e-4;
1290: bnk->eta2 = 0.25;
1291: bnk->eta3 = 0.50;
1292: bnk->eta4 = 0.90;
1294: bnk->alpha1 = 0.25;
1295: bnk->alpha2 = 0.50;
1296: bnk->alpha3 = 1.00;
1297: bnk->alpha4 = 2.00;
1298: bnk->alpha5 = 4.00;
1300: /* Default values for trust-region radius update based on interpolation */
1301: bnk->mu1 = 0.10;
1302: bnk->mu2 = 0.50;
1304: bnk->gamma1 = 0.25;
1305: bnk->gamma2 = 0.50;
1306: bnk->gamma3 = 2.00;
1307: bnk->gamma4 = 4.00;
1309: bnk->theta = 0.05;
1311: /* Default values for trust region initialization based on interpolation */
1312: bnk->mu1_i = 0.35;
1313: bnk->mu2_i = 0.50;
1315: bnk->gamma1_i = 0.0625;
1316: bnk->gamma2_i = 0.5;
1317: bnk->gamma3_i = 2.0;
1318: bnk->gamma4_i = 5.0;
1320: bnk->theta_i = 0.25;
1322: /* Remaining parameters */
1323: bnk->max_cg_its = 0;
1324: bnk->min_radius = 1.0e-10;
1325: bnk->max_radius = 1.0e10;
1326: bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0 / 3.0);
1327: bnk->as_tol = 1.0e-3;
1328: bnk->as_step = 1.0e-3;
1329: bnk->dmin = 1.0e-6;
1330: bnk->dmax = 1.0e6;
1332: bnk->M = NULL;
1333: bnk->bfgs_pre = NULL;
1334: bnk->init_type = BNK_INIT_INTERPOLATION;
1335: bnk->update_type = BNK_UPDATE_REDUCTION;
1336: bnk->as_type = BNK_AS_BERTSEKAS;
1338: /* Create the embedded BNCG solver */
1339: PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg));
1340: PetscCall(PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1));
1341: PetscCall(TaoSetType(bnk->bncg, TAOBNCG));
1343: /* Create the line search */
1344: PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
1345: PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
1346: PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT));
1347: PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
1349: /* Set linear solver to default for symmetric matrices */
1350: PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
1351: PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
1352: PetscCall(KSPSetType(tao->ksp, KSPSTCG));
1353: PetscCall(KSPGetPC(tao->ksp, &pc));
1354: PetscCall(PCSetType(pc, PCLMVM));
1355: PetscFunctionReturn(PETSC_SUCCESS);
1356: }