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