1eb910715SAlp Dener #include <petsctaolinesearch.h> 2eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h> 3eb910715SAlp Dener 4eb910715SAlp Dener #include <petscksp.h> 5eb910715SAlp Dener 6eb910715SAlp Dener /* Routine for BFGS preconditioner */ 7eb910715SAlp Dener 8eb910715SAlp Dener PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x) 9eb910715SAlp Dener { 10eb910715SAlp Dener PetscErrorCode ierr; 11eb910715SAlp Dener Mat M; 12eb910715SAlp Dener 13eb910715SAlp Dener PetscFunctionBegin; 14eb910715SAlp Dener PetscValidHeaderSpecific(pc,PC_CLASSID,1); 15eb910715SAlp Dener PetscValidHeaderSpecific(b,VEC_CLASSID,2); 16eb910715SAlp Dener PetscValidHeaderSpecific(x,VEC_CLASSID,3); 17eb910715SAlp Dener ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr); 1809164190SAlp Dener ierr = MatLMVMSolveInactive(M, b, x);CHKERRQ(ierr); 19eb910715SAlp Dener PetscFunctionReturn(0); 20eb910715SAlp Dener } 21eb910715SAlp Dener 22df278d8fSAlp Dener /*------------------------------------------------------------*/ 23df278d8fSAlp Dener 24df278d8fSAlp Dener /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */ 25df278d8fSAlp Dener 2662675beeSAlp Dener PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType) 27eb910715SAlp Dener { 28eb910715SAlp Dener PetscErrorCode ierr; 29eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 30eb910715SAlp Dener PC pc; 31eb910715SAlp Dener 320a4511e9SAlp Dener PetscReal f_min, ftrial, prered, actred, kappa, sigma; 33eb910715SAlp Dener PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 3428017e9fSAlp Dener PetscReal delta; 35eb910715SAlp Dener 36eb910715SAlp Dener PetscInt n,N,needH = 1; 37eb910715SAlp Dener 38eb910715SAlp Dener PetscInt i_max = 5; 39eb910715SAlp Dener PetscInt j_max = 1; 40eb910715SAlp Dener PetscInt i, j; 41eb910715SAlp Dener 42eb910715SAlp Dener PetscFunctionBegin; 4328017e9fSAlp Dener /* Project the current point onto the feasible set */ 4428017e9fSAlp Dener ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 4528017e9fSAlp Dener ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 4628017e9fSAlp Dener 4728017e9fSAlp Dener /* Project the initial point onto the feasible region */ 4828017e9fSAlp Dener ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 4928017e9fSAlp Dener 5028017e9fSAlp Dener /* Check convergence criteria */ 5128017e9fSAlp Dener ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr); 5228017e9fSAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 5328017e9fSAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 5428017e9fSAlp Dener if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 5528017e9fSAlp Dener 56eb910715SAlp Dener /* Number of times ksp stopped because of these reasons */ 57eb910715SAlp Dener bnk->ksp_atol = 0; 58eb910715SAlp Dener bnk->ksp_rtol = 0; 59eb910715SAlp Dener bnk->ksp_dtol = 0; 60eb910715SAlp Dener bnk->ksp_ctol = 0; 61eb910715SAlp Dener bnk->ksp_negc = 0; 62eb910715SAlp Dener bnk->ksp_iter = 0; 63eb910715SAlp Dener bnk->ksp_othr = 0; 64eb910715SAlp Dener 65fed79b8eSAlp Dener /* Initialize the Hessian perturbation */ 66fed79b8eSAlp Dener bnk->pert = bnk->sval; 67fed79b8eSAlp Dener 68eb910715SAlp Dener /* Initialize trust-region radius when using nash, stcg, or gltr 69eb910715SAlp Dener Command automatically ignored for other methods 70eb910715SAlp Dener Will be reset during the first iteration 71eb910715SAlp Dener */ 72eb910715SAlp Dener ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr); 73eb910715SAlp Dener if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 74eb910715SAlp Dener if (tao->trust0 < 0.0) SETERRQ(PETSC_COMM_SELF,1,"Initial radius negative"); 75eb910715SAlp Dener tao->trust = tao->trust0; 76eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 77eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 78eb910715SAlp Dener } 79eb910715SAlp Dener 8028017e9fSAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 8128017e9fSAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,1.0);CHKERRQ(ierr); 8228017e9fSAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 8328017e9fSAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 8428017e9fSAlp Dener 85eb910715SAlp Dener /* Get vectors we will need */ 86eb910715SAlp Dener if (BNK_PC_BFGS == bnk->pc_type && !bnk->M) { 87eb910715SAlp Dener ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 88eb910715SAlp Dener ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 89eb910715SAlp Dener ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr); 90eb910715SAlp Dener ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr); 91eb910715SAlp Dener } 92eb910715SAlp Dener 93eb910715SAlp Dener /* create vectors for the limited memory preconditioner */ 94eb910715SAlp Dener if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_BFGS != bnk->bfgs_scale_type)) { 9562675beeSAlp Dener if (!bnk->Diag) {ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr);} 9662675beeSAlp Dener if (!bnk->Diag_min) {ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);} 9762675beeSAlp Dener if (!bnk->Diag_max) {ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);} 9862675beeSAlp Dener ierr = VecSet(bnk->Diag_min, bnk->dmin);CHKERRQ(ierr); 9962675beeSAlp Dener ierr = VecSet(bnk->Diag_max, bnk->dmax);CHKERRQ(ierr); 100eb910715SAlp Dener } 101eb910715SAlp Dener 102eb910715SAlp Dener /* Modify the preconditioner to use the bfgs approximation */ 103eb910715SAlp Dener ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 104eb910715SAlp Dener switch(bnk->pc_type) { 105eb910715SAlp Dener case BNK_PC_NONE: 106eb910715SAlp Dener ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr); 107eb910715SAlp Dener ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 108eb910715SAlp Dener break; 109eb910715SAlp Dener 110eb910715SAlp Dener case BNK_PC_AHESS: 111eb910715SAlp Dener ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr); 112eb910715SAlp Dener ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 113eb910715SAlp Dener ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 114eb910715SAlp Dener break; 115eb910715SAlp Dener 116eb910715SAlp Dener case BNK_PC_BFGS: 117eb910715SAlp Dener ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr); 118eb910715SAlp Dener ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 119eb910715SAlp Dener ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr); 120eb910715SAlp Dener ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr); 121eb910715SAlp Dener ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr); 122eb910715SAlp Dener break; 123eb910715SAlp Dener 124eb910715SAlp Dener default: 125eb910715SAlp Dener /* Use the pc method set by pc_type */ 126eb910715SAlp Dener break; 127eb910715SAlp Dener } 128eb910715SAlp Dener 129eb910715SAlp Dener /* Initialize trust-region radius. The initialization is only performed 130eb910715SAlp Dener when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */ 131eb910715SAlp Dener if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 13262675beeSAlp Dener switch(initType) { 133eb910715SAlp Dener case BNK_INIT_CONSTANT: 134eb910715SAlp Dener /* Use the initial radius specified */ 135eb910715SAlp Dener break; 136eb910715SAlp Dener 137eb910715SAlp Dener case BNK_INIT_INTERPOLATION: 138eb910715SAlp Dener /* Use the initial radius specified */ 139eb910715SAlp Dener max_radius = 0.0; 140eb910715SAlp Dener 141eb910715SAlp Dener for (j = 0; j < j_max; ++j) { 1420a4511e9SAlp Dener f_min = bnk->f; 143eb910715SAlp Dener sigma = 0.0; 144eb910715SAlp Dener 145eb910715SAlp Dener if (needH) { 14662602cfbSAlp Dener /* Compute the Hessian at the new step, and extract the inactive subsystem */ 147eb910715SAlp Dener ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 14828017e9fSAlp Dener if (bnk->inactive_idx) { 14962602cfbSAlp Dener ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr); 15062602cfbSAlp Dener ierr = VecWhichInactive(tao->XL,tao->solution,bnk->unprojected_gradient,tao->XU,PETSC_TRUE,&bnk->inactive_idx);CHKERRQ(ierr); 15162602cfbSAlp Dener ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr); 15228017e9fSAlp Dener } else { 15362675beeSAlp Dener ierr = MatDestroy(&bnk->H_inactive); 15462675beeSAlp Dener ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive); 15528017e9fSAlp Dener } 156eb910715SAlp Dener needH = 0; 157eb910715SAlp Dener } 158eb910715SAlp Dener 159eb910715SAlp Dener for (i = 0; i < i_max; ++i) { 16062602cfbSAlp Dener /* Take a steepest descent step and snap it to bounds */ 16162602cfbSAlp Dener ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 16262602cfbSAlp Dener ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr); 16362602cfbSAlp Dener ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 164770b7498SAlp Dener /* Recompute the step after bound snapping so that it can be used in predicted decrease calculation later */ 165eb910715SAlp Dener ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr); 16662602cfbSAlp Dener ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr); 16762602cfbSAlp Dener /* Compute the objective at the trial */ 16862602cfbSAlp Dener ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr); 16962602cfbSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 170eb910715SAlp Dener if (PetscIsInfOrNanReal(ftrial)) { 171eb910715SAlp Dener tau = bnk->gamma1_i; 172eb910715SAlp Dener } else { 1730a4511e9SAlp Dener if (ftrial < f_min) { 1740a4511e9SAlp Dener f_min = ftrial; 175eb910715SAlp Dener sigma = -tao->trust / bnk->gnorm; 176eb910715SAlp Dener } 177770b7498SAlp Dener /* Compute the predicted and actual reduction */ 17862602cfbSAlp Dener ierr = MatMult(bnk->H_inactive, tao->gradient, bnk->W);CHKERRQ(ierr); 17962602cfbSAlp Dener ierr = VecDot(tao->gradient, bnk->W, &prered);CHKERRQ(ierr); 180eb910715SAlp Dener prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm)); 181eb910715SAlp Dener actred = bnk->f - ftrial; 182eb910715SAlp Dener if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 183eb910715SAlp Dener kappa = 1.0; 184eb910715SAlp Dener } else { 185eb910715SAlp Dener kappa = actred / prered; 186eb910715SAlp Dener } 187eb910715SAlp Dener 188eb910715SAlp Dener tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred); 189eb910715SAlp Dener tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred); 190eb910715SAlp Dener tau_min = PetscMin(tau_1, tau_2); 191eb910715SAlp Dener tau_max = PetscMax(tau_1, tau_2); 192eb910715SAlp Dener 193eb910715SAlp Dener if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) { 194eb910715SAlp Dener /* Great agreement */ 195eb910715SAlp Dener max_radius = PetscMax(max_radius, tao->trust); 196eb910715SAlp Dener 197eb910715SAlp Dener if (tau_max < 1.0) { 198eb910715SAlp Dener tau = bnk->gamma3_i; 199eb910715SAlp Dener } else if (tau_max > bnk->gamma4_i) { 200eb910715SAlp Dener tau = bnk->gamma4_i; 201eb910715SAlp Dener } else if (tau_1 >= 1.0 && tau_1 <= bnk->gamma4_i && tau_2 < 1.0) { 202eb910715SAlp Dener tau = tau_1; 203eb910715SAlp Dener } else if (tau_2 >= 1.0 && tau_2 <= bnk->gamma4_i && tau_1 < 1.0) { 204eb910715SAlp Dener tau = tau_2; 205eb910715SAlp Dener } else { 206eb910715SAlp Dener tau = tau_max; 207eb910715SAlp Dener } 208eb910715SAlp Dener } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) { 209eb910715SAlp Dener /* Good agreement */ 210eb910715SAlp Dener max_radius = PetscMax(max_radius, tao->trust); 211eb910715SAlp Dener 212eb910715SAlp Dener if (tau_max < bnk->gamma2_i) { 213eb910715SAlp Dener tau = bnk->gamma2_i; 214eb910715SAlp Dener } else if (tau_max > bnk->gamma3_i) { 215eb910715SAlp Dener tau = bnk->gamma3_i; 216eb910715SAlp Dener } else { 217eb910715SAlp Dener tau = tau_max; 218eb910715SAlp Dener } 219eb910715SAlp Dener } else { 220eb910715SAlp Dener /* Not good agreement */ 221eb910715SAlp Dener if (tau_min > 1.0) { 222eb910715SAlp Dener tau = bnk->gamma2_i; 223eb910715SAlp Dener } else if (tau_max < bnk->gamma1_i) { 224eb910715SAlp Dener tau = bnk->gamma1_i; 225eb910715SAlp Dener } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) { 226eb910715SAlp Dener tau = bnk->gamma1_i; 227eb910715SAlp Dener } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 228eb910715SAlp Dener tau = tau_1; 229eb910715SAlp Dener } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 230eb910715SAlp Dener tau = tau_2; 231eb910715SAlp Dener } else { 232eb910715SAlp Dener tau = tau_max; 233eb910715SAlp Dener } 234eb910715SAlp Dener } 235eb910715SAlp Dener } 236eb910715SAlp Dener tao->trust = tau * tao->trust; 237eb910715SAlp Dener } 238eb910715SAlp Dener 2390a4511e9SAlp Dener if (f_min < bnk->f) { 2400a4511e9SAlp Dener bnk->f = f_min; 241eb910715SAlp Dener ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr); 24287602d1eSAlp Dener ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 24309164190SAlp Dener ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr); 24409164190SAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 245eb910715SAlp Dener 246eb910715SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 247eb910715SAlp Dener if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN"); 248eb910715SAlp Dener needH = 1; 249eb910715SAlp Dener 250eb910715SAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 25128017e9fSAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,1.0);CHKERRQ(ierr); 252eb910715SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 253eb910715SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 254eb910715SAlp Dener } 255eb910715SAlp Dener } 256eb910715SAlp Dener tao->trust = PetscMax(tao->trust, max_radius); 257eb910715SAlp Dener 258eb910715SAlp Dener /* Modify the radius if it is too large or small */ 259eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 260eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 261eb910715SAlp Dener break; 262eb910715SAlp Dener 263eb910715SAlp Dener default: 264eb910715SAlp Dener /* Norm of the first direction will initialize radius */ 265eb910715SAlp Dener tao->trust = 0.0; 266eb910715SAlp Dener break; 267eb910715SAlp Dener } 268eb910715SAlp Dener } 269eb910715SAlp Dener 270eb910715SAlp Dener /* Set initial scaling for the BFGS preconditioner 271eb910715SAlp Dener This step is done after computing the initial trust-region radius 272eb910715SAlp Dener since the function value may have decreased */ 273eb910715SAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { 274*9b6ef848SAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 275eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr); 276eb910715SAlp Dener } 277eb910715SAlp Dener 278eb910715SAlp Dener /* Set counter for gradient/reset steps*/ 279eb910715SAlp Dener bnk->newt = 0; 280eb910715SAlp Dener bnk->bfgs = 0; 281eb910715SAlp Dener bnk->sgrad = 0; 282eb910715SAlp Dener bnk->grad = 0; 283eb910715SAlp Dener PetscFunctionReturn(0); 284eb910715SAlp Dener } 285eb910715SAlp Dener 286df278d8fSAlp Dener /*------------------------------------------------------------*/ 287df278d8fSAlp Dener 28862675beeSAlp Dener /* Routine for computing the Hessian and preparing the preconditioner at the new iterate */ 28962675beeSAlp Dener 29062675beeSAlp Dener PetscErrorCode TaoBNKComputeHessian(Tao tao) 29162675beeSAlp Dener { 29262675beeSAlp Dener PetscErrorCode ierr; 29362675beeSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 29462675beeSAlp Dener 29562675beeSAlp Dener PetscFunctionBegin; 29662675beeSAlp Dener /* Compute the Hessian */ 29762675beeSAlp Dener ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 29862675beeSAlp Dener /* Add a correction to the BFGS preconditioner */ 29962675beeSAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { 30062675beeSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 30162675beeSAlp Dener /* Update the BFGS diagonal scaling */ 30262675beeSAlp Dener if (BFGS_SCALE_AHESS == bnk->bfgs_scale_type) { 30362675beeSAlp Dener ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr); 30462675beeSAlp Dener ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 30562675beeSAlp Dener ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr); 30662675beeSAlp Dener ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 30762675beeSAlp Dener ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 30862675beeSAlp Dener } 30962675beeSAlp Dener } 31062675beeSAlp Dener PetscFunctionReturn(0); 31162675beeSAlp Dener } 31262675beeSAlp Dener 31362675beeSAlp Dener /*------------------------------------------------------------*/ 31462675beeSAlp Dener 3152f75a4aaSAlp Dener /* Routine for estimating the active set */ 3162f75a4aaSAlp Dener 3172f75a4aaSAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao) 3182f75a4aaSAlp Dener { 3192f75a4aaSAlp Dener PetscErrorCode ierr; 3202f75a4aaSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 3212f75a4aaSAlp Dener 3222f75a4aaSAlp Dener PetscFunctionBegin; 3232f75a4aaSAlp Dener switch (bnk->as_type) { 3242f75a4aaSAlp Dener case BNK_AS_NONE: 3252f75a4aaSAlp Dener ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr); 3262f75a4aaSAlp Dener ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr); 3272f75a4aaSAlp Dener ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr); 3282f75a4aaSAlp Dener ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr); 3292f75a4aaSAlp Dener break; 3302f75a4aaSAlp Dener 3312f75a4aaSAlp Dener case BNK_AS_BERTSEKAS: 3322f75a4aaSAlp Dener /* Compute the trial step vector with which we will estimate the active set at the next iteration */ 3332f75a4aaSAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { 3342f75a4aaSAlp Dener /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */ 3352f75a4aaSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr); 3362f75a4aaSAlp Dener } else { 337*9b6ef848SAlp Dener /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */ 3382f75a4aaSAlp Dener ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr); 339*9b6ef848SAlp Dener ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr); 340*9b6ef848SAlp Dener ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr); 3412f75a4aaSAlp Dener ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr); 3422f75a4aaSAlp Dener ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr); 3432f75a4aaSAlp Dener } 3442f75a4aaSAlp Dener ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr); 3450a4511e9SAlp Dener ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->as_step, &bnk->as_tol, &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr); 3462f75a4aaSAlp Dener 3472f75a4aaSAlp Dener default: 3482f75a4aaSAlp Dener break; 3492f75a4aaSAlp Dener } 3502f75a4aaSAlp Dener PetscFunctionReturn(0); 3512f75a4aaSAlp Dener } 3522f75a4aaSAlp Dener 3532f75a4aaSAlp Dener /*------------------------------------------------------------*/ 3542f75a4aaSAlp Dener 3552f75a4aaSAlp Dener /* Routine for bounding the step direction */ 3562f75a4aaSAlp Dener 3572f75a4aaSAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, Vec step) 3582f75a4aaSAlp Dener { 3592f75a4aaSAlp Dener PetscErrorCode ierr; 3602f75a4aaSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 3612f75a4aaSAlp Dener 3622f75a4aaSAlp Dener PetscFunctionBegin; 36328017e9fSAlp Dener if (bnk->active_idx) { 3642f75a4aaSAlp Dener switch (bnk->as_type) { 3652f75a4aaSAlp Dener case BNK_AS_NONE: 3662f75a4aaSAlp Dener if (bnk->active_idx) { 3672f75a4aaSAlp Dener ierr = VecGetSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr); 3682f75a4aaSAlp Dener ierr = VecSet(bnk->active_work, 0.0);CHKERRQ(ierr); 3692f75a4aaSAlp Dener ierr = VecRestoreSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr); 3702f75a4aaSAlp Dener } 3712f75a4aaSAlp Dener break; 3722f75a4aaSAlp Dener 3732f75a4aaSAlp Dener case BNK_AS_BERTSEKAS: 3742f75a4aaSAlp Dener ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, step);CHKERRQ(ierr); 3752f75a4aaSAlp Dener break; 3762f75a4aaSAlp Dener 3772f75a4aaSAlp Dener default: 3782f75a4aaSAlp Dener break; 3792f75a4aaSAlp Dener } 380e465cd6fSAlp Dener } 3812f75a4aaSAlp Dener PetscFunctionReturn(0); 3822f75a4aaSAlp Dener } 3832f75a4aaSAlp Dener 3842f75a4aaSAlp Dener /*------------------------------------------------------------*/ 3852f75a4aaSAlp Dener 386df278d8fSAlp Dener /* Routine for computing the Newton step. 387df278d8fSAlp Dener 388df278d8fSAlp Dener If the safeguard is enabled, the Newton step is verified to be a 389df278d8fSAlp Dener descent direction, with fallbacks onto BFGS, scaled gradient, and unscaled 390df278d8fSAlp Dener gradient steps if/when necessary. 391df278d8fSAlp Dener 392df278d8fSAlp Dener The function reports back on which type of step has ultimately been stored 393df278d8fSAlp Dener under tao->stepdirection. 394df278d8fSAlp Dener */ 395df278d8fSAlp Dener 39662675beeSAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason) 397eb910715SAlp Dener { 398eb910715SAlp Dener PetscErrorCode ierr; 399eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 400eb910715SAlp Dener 401e465cd6fSAlp Dener PetscReal delta; 402eb910715SAlp Dener PetscInt bfgsUpdates = 0; 403eb910715SAlp Dener PetscInt kspits; 404eb910715SAlp Dener 405eb910715SAlp Dener PetscFunctionBegin; 4062f75a4aaSAlp Dener /* Determine the active and inactive sets */ 4072f75a4aaSAlp Dener ierr = TaoBNKEstimateActiveSet(tao);CHKERRQ(ierr); 40809164190SAlp Dener 409eb3ba6a7SAlp Dener /* Prepare masked matrices for the inactive set */ 4102f75a4aaSAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr); } 4112f75a4aaSAlp Dener if (bnk->inactive_idx) { 412eb3ba6a7SAlp Dener ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr); 413eb3ba6a7SAlp Dener if (tao->hessian == tao->hessian_pre) { 414eb3ba6a7SAlp Dener bnk->Hpre_inactive = bnk->H_inactive; 415eb3ba6a7SAlp Dener } else { 416eb3ba6a7SAlp Dener ierr = TaoMatGetSubMat(tao->hessian_pre, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->Hpre_inactive);CHKERRQ(ierr); 417eb3ba6a7SAlp Dener } 4182f75a4aaSAlp Dener } else { 41962675beeSAlp Dener ierr = MatDestroy(&bnk->H_inactive); 42062675beeSAlp Dener ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive); 42162675beeSAlp Dener if (tao->hessian == tao->hessian_pre) { 42262675beeSAlp Dener bnk->Hpre_inactive = bnk->H_inactive; 42362675beeSAlp Dener } else { 42462675beeSAlp Dener ierr = MatDestroy(&bnk->Hpre_inactive); 42562675beeSAlp Dener ierr = MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive); 42662675beeSAlp Dener } 42762675beeSAlp Dener } 42862675beeSAlp Dener 42962675beeSAlp Dener /* Shift the reduced Hessian matrix */ 43062675beeSAlp Dener if ((shift) && (bnk->pert > 0)) { 43162675beeSAlp Dener ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr); 43262675beeSAlp Dener if (bnk->H_inactive != bnk->Hpre_inactive) { 43362675beeSAlp Dener ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr); 43462675beeSAlp Dener } 43562675beeSAlp Dener } 43662675beeSAlp Dener 43762675beeSAlp Dener /* Update the diagonal scaling for the BFGS preconditioner, this time with the Hessian perturbation */ 43862675beeSAlp Dener if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_PHESS == bnk->bfgs_scale_type)) { 43962675beeSAlp Dener /* Obtain diagonal for the bfgs preconditioner */ 44062675beeSAlp Dener ierr = MatGetDiagonal(bnk->H_inactive, bnk->Diag);CHKERRQ(ierr); 44162675beeSAlp Dener ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 44262675beeSAlp Dener ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr); 44362675beeSAlp Dener ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 44462675beeSAlp Dener ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 4452f75a4aaSAlp Dener } 44609164190SAlp Dener 447eb910715SAlp Dener /* Solve the Newton system of equations */ 4482f75a4aaSAlp Dener ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr); 44909164190SAlp Dener ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr); 450198282dbSAlp Dener ierr = VecCopy(bnk->unprojected_gradient, bnk->G_inactive);CHKERRQ(ierr); 45128017e9fSAlp Dener if (bnk->active_idx) { 452198282dbSAlp Dener ierr = VecGetSubVector(bnk->G_inactive, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr); 45328017e9fSAlp Dener ierr = VecSet(bnk->active_work, 0.0);CHKERRQ(ierr); 454198282dbSAlp Dener ierr = VecRestoreSubVector(bnk->G_inactive, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr); 45528017e9fSAlp Dener } 456eb910715SAlp Dener if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 457fed79b8eSAlp Dener ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 458198282dbSAlp Dener ierr = KSPSolve(tao->ksp, bnk->G_inactive, tao->stepdirection);CHKERRQ(ierr); 459eb910715SAlp Dener ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 460eb910715SAlp Dener tao->ksp_its+=kspits; 461eb910715SAlp Dener tao->ksp_tot_its+=kspits; 462080d2917SAlp Dener ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr); 463eb910715SAlp Dener 464eb910715SAlp Dener if (0.0 == tao->trust) { 465eb910715SAlp Dener /* Radius was uninitialized; use the norm of the direction */ 466080d2917SAlp Dener if (bnk->dnorm > 0.0) { 467080d2917SAlp Dener tao->trust = bnk->dnorm; 468eb910715SAlp Dener 469eb910715SAlp Dener /* Modify the radius if it is too large or small */ 470eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 471eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 472eb910715SAlp Dener } else { 473eb910715SAlp Dener /* The direction was bad; set radius to default value and re-solve 474eb910715SAlp Dener the trust-region subproblem to get a direction */ 475eb910715SAlp Dener tao->trust = tao->trust0; 476eb910715SAlp Dener 477eb910715SAlp Dener /* Modify the radius if it is too large or small */ 478eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 479eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 480eb910715SAlp Dener 481fed79b8eSAlp Dener ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 482198282dbSAlp Dener ierr = KSPSolve(tao->ksp, bnk->G_inactive, tao->stepdirection);CHKERRQ(ierr); 483eb910715SAlp Dener ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 484eb910715SAlp Dener tao->ksp_its+=kspits; 485eb910715SAlp Dener tao->ksp_tot_its+=kspits; 486080d2917SAlp Dener ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr); 487eb910715SAlp Dener 488080d2917SAlp Dener if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 489eb910715SAlp Dener } 490eb910715SAlp Dener } 491eb910715SAlp Dener } else { 492198282dbSAlp Dener ierr = KSPSolve(tao->ksp, bnk->G_inactive, tao->stepdirection);CHKERRQ(ierr); 493eb910715SAlp Dener ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr); 494eb910715SAlp Dener tao->ksp_its += kspits; 495eb910715SAlp Dener tao->ksp_tot_its+=kspits; 496eb910715SAlp Dener } 497770b7498SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 498fed79b8eSAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 4992f75a4aaSAlp Dener ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr); 500770b7498SAlp Dener 501770b7498SAlp Dener /* Record convergence reasons */ 502e465cd6fSAlp Dener ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr); 503e465cd6fSAlp Dener if (KSP_CONVERGED_ATOL == *ksp_reason) { 504770b7498SAlp Dener ++bnk->ksp_atol; 505e465cd6fSAlp Dener } else if (KSP_CONVERGED_RTOL == *ksp_reason) { 506770b7498SAlp Dener ++bnk->ksp_rtol; 507e465cd6fSAlp Dener } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) { 508770b7498SAlp Dener ++bnk->ksp_ctol; 509e465cd6fSAlp Dener } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) { 510770b7498SAlp Dener ++bnk->ksp_negc; 511e465cd6fSAlp Dener } else if (KSP_DIVERGED_DTOL == *ksp_reason) { 512770b7498SAlp Dener ++bnk->ksp_dtol; 513e465cd6fSAlp Dener } else if (KSP_DIVERGED_ITS == *ksp_reason) { 514770b7498SAlp Dener ++bnk->ksp_iter; 515770b7498SAlp Dener } else { 516770b7498SAlp Dener ++bnk->ksp_othr; 517770b7498SAlp Dener } 518fed79b8eSAlp Dener 519fed79b8eSAlp Dener /* Make sure the BFGS preconditioner is healthy */ 520fed79b8eSAlp Dener if (bnk->pc_type == BNK_PC_BFGS) { 521fed79b8eSAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 522e465cd6fSAlp Dener if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 1)) { 523fed79b8eSAlp Dener /* Preconditioner is numerically indefinite; reset the approximation. */ 524*9b6ef848SAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 525eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr); 526eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 52709164190SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 528eb910715SAlp Dener } 529fed79b8eSAlp Dener } 530e465cd6fSAlp Dener PetscFunctionReturn(0); 531e465cd6fSAlp Dener } 532eb910715SAlp Dener 53362675beeSAlp Dener /*------------------------------------------------------------*/ 53462675beeSAlp Dener 53562675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction. 53662675beeSAlp Dener 53762675beeSAlp Dener The step direction falls back onto BFGS, scaled gradient and gradient steps 53862675beeSAlp Dener in the event that the Newton step fails the test. 53962675beeSAlp Dener */ 54062675beeSAlp Dener 541e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType) 542e465cd6fSAlp Dener { 543e465cd6fSAlp Dener PetscErrorCode ierr; 544e465cd6fSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 545e465cd6fSAlp Dener 546e465cd6fSAlp Dener PetscReal gdx, delta, e_min; 547e465cd6fSAlp Dener PetscInt bfgsUpdates; 548e465cd6fSAlp Dener 549e465cd6fSAlp Dener PetscFunctionBegin; 550080d2917SAlp Dener ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 551eb910715SAlp Dener if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 552eb910715SAlp Dener /* Newton step is not descent or direction produced Inf or NaN 553eb910715SAlp Dener Update the perturbation for next time */ 554eb910715SAlp Dener if (bnk->pert <= 0.0) { 555eb910715SAlp Dener /* Initialize the perturbation */ 556eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 557eb910715SAlp Dener if (bnk->is_gltr) { 558eb910715SAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 559eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 560eb910715SAlp Dener } 561eb910715SAlp Dener } else { 562eb910715SAlp Dener /* Increase the perturbation */ 563eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 564eb910715SAlp Dener } 565eb910715SAlp Dener 566eb910715SAlp Dener if (BNK_PC_BFGS != bnk->pc_type) { 567eb910715SAlp Dener /* We don't have the bfgs matrix around and updated 568eb910715SAlp Dener Must use gradient direction in this case */ 569080d2917SAlp Dener ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 570eb910715SAlp Dener *stepType = BNK_GRADIENT; 571eb910715SAlp Dener } else { 572eb910715SAlp Dener /* Attempt to use the BFGS direction */ 57309164190SAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 574eb910715SAlp Dener 5758d5ead36SAlp Dener /* Check for success (descent direction) 5768d5ead36SAlp Dener NOTE: Negative gdx here means not a descent direction because 5778d5ead36SAlp Dener the fall-back step is missing a negative sign. */ 578080d2917SAlp Dener ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 5798d5ead36SAlp Dener if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 580eb910715SAlp Dener /* BFGS direction is not descent or direction produced not a number 581eb910715SAlp Dener We can assert bfgsUpdates > 1 in this case because 582eb910715SAlp Dener the first solve produces the scaled gradient direction, 583eb910715SAlp Dener which is guaranteed to be descent */ 584eb910715SAlp Dener 585eb910715SAlp Dener /* Use steepest descent direction (scaled) */ 586*9b6ef848SAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 587eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 588eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 58909164190SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 59009164190SAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 591eb910715SAlp Dener 592eb910715SAlp Dener *stepType = BNK_SCALED_GRADIENT; 593eb910715SAlp Dener } else { 594770b7498SAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 595eb910715SAlp Dener if (1 == bfgsUpdates) { 596eb910715SAlp Dener /* The first BFGS direction is always the scaled gradient */ 597eb910715SAlp Dener *stepType = BNK_SCALED_GRADIENT; 598eb910715SAlp Dener } else { 599eb910715SAlp Dener *stepType = BNK_BFGS; 600eb910715SAlp Dener } 601eb910715SAlp Dener } 602eb910715SAlp Dener } 6038d5ead36SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 6048d5ead36SAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 6052f75a4aaSAlp Dener ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr); 606eb910715SAlp Dener } else { 607eb910715SAlp Dener /* Computed Newton step is descent */ 608eb910715SAlp Dener switch (ksp_reason) { 609eb910715SAlp Dener case KSP_DIVERGED_NANORINF: 610eb910715SAlp Dener case KSP_DIVERGED_BREAKDOWN: 611eb910715SAlp Dener case KSP_DIVERGED_INDEFINITE_MAT: 612eb910715SAlp Dener case KSP_DIVERGED_INDEFINITE_PC: 613eb910715SAlp Dener case KSP_CONVERGED_CG_NEG_CURVE: 614eb910715SAlp Dener /* Matrix or preconditioner is indefinite; increase perturbation */ 615eb910715SAlp Dener if (bnk->pert <= 0.0) { 616eb910715SAlp Dener /* Initialize the perturbation */ 617eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 618eb910715SAlp Dener if (bnk->is_gltr) { 619eb910715SAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr); 620eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 621eb910715SAlp Dener } 622eb910715SAlp Dener } else { 623eb910715SAlp Dener /* Increase the perturbation */ 624eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 625eb910715SAlp Dener } 626eb910715SAlp Dener break; 627eb910715SAlp Dener 628eb910715SAlp Dener default: 629eb910715SAlp Dener /* Newton step computation is good; decrease perturbation */ 630eb910715SAlp Dener bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm); 631eb910715SAlp Dener if (bnk->pert < bnk->pmin) { 632eb910715SAlp Dener bnk->pert = 0.0; 633eb910715SAlp Dener } 634eb910715SAlp Dener break; 635eb910715SAlp Dener } 636fed79b8eSAlp Dener *stepType = BNK_NEWTON; 637eb910715SAlp Dener } 638eb910715SAlp Dener PetscFunctionReturn(0); 639eb910715SAlp Dener } 640eb910715SAlp Dener 641df278d8fSAlp Dener /*------------------------------------------------------------*/ 642df278d8fSAlp Dener 643df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search. 644df278d8fSAlp Dener 645df278d8fSAlp Dener Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the 646df278d8fSAlp Dener Newton step does not produce a valid step length. 647df278d8fSAlp Dener */ 648df278d8fSAlp Dener 649c14b763aSAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason) 650c14b763aSAlp Dener { 651c14b763aSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 652c14b763aSAlp Dener PetscErrorCode ierr; 653c14b763aSAlp Dener TaoLineSearchConvergedReason ls_reason; 654c14b763aSAlp Dener 655c14b763aSAlp Dener PetscReal e_min, gdx, delta; 656c14b763aSAlp Dener PetscInt bfgsUpdates; 657c14b763aSAlp Dener 658c14b763aSAlp Dener PetscFunctionBegin; 659c14b763aSAlp Dener /* Perform the linesearch */ 660c14b763aSAlp Dener ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr); 661c14b763aSAlp Dener ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 662c14b763aSAlp Dener 663*9b6ef848SAlp Dener while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && (stepType != BNK_GRADIENT || stepType !=BNK_SCALED_GRADIENT)) { 664c14b763aSAlp Dener /* Linesearch failed, revert solution */ 665c14b763aSAlp Dener bnk->f = bnk->fold; 666c14b763aSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 667c14b763aSAlp Dener ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 668c14b763aSAlp Dener 669c14b763aSAlp Dener switch(stepType) { 670c14b763aSAlp Dener case BNK_NEWTON: 6718d5ead36SAlp Dener /* Failed to obtain acceptable iterate with Newton step 672c14b763aSAlp Dener Update the perturbation for next time */ 673c14b763aSAlp Dener if (bnk->pert <= 0.0) { 674c14b763aSAlp Dener /* Initialize the perturbation */ 675c14b763aSAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 676c14b763aSAlp Dener if (bnk->is_gltr) { 677c14b763aSAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 678c14b763aSAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 679c14b763aSAlp Dener } 680c14b763aSAlp Dener } else { 681c14b763aSAlp Dener /* Increase the perturbation */ 682c14b763aSAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 683c14b763aSAlp Dener } 684c14b763aSAlp Dener 685c14b763aSAlp Dener if (BNK_PC_BFGS != bnk->pc_type) { 686c14b763aSAlp Dener /* We don't have the bfgs matrix around and being updated 687c14b763aSAlp Dener Must use gradient direction in this case */ 688c14b763aSAlp Dener ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 689c14b763aSAlp Dener stepType = BNK_GRADIENT; 690c14b763aSAlp Dener } else { 691c14b763aSAlp Dener /* Attempt to use the BFGS direction */ 692c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 6938d5ead36SAlp Dener /* Check for success (descent direction) 6948d5ead36SAlp Dener NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */ 695c14b763aSAlp Dener ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 696c14b763aSAlp Dener if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 697c14b763aSAlp Dener /* BFGS direction is not descent or direction produced not a number 698c14b763aSAlp Dener We can assert bfgsUpdates > 1 in this case 699c14b763aSAlp Dener Use steepest descent direction (scaled) */ 700*9b6ef848SAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 701c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 702c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 703c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 704c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 705c14b763aSAlp Dener 706c14b763aSAlp Dener bfgsUpdates = 1; 707c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 708c14b763aSAlp Dener } else { 709c14b763aSAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 710c14b763aSAlp Dener if (1 == bfgsUpdates) { 711c14b763aSAlp Dener /* The first BFGS direction is always the scaled gradient */ 712c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 713c14b763aSAlp Dener } else { 714c14b763aSAlp Dener stepType = BNK_BFGS; 715c14b763aSAlp Dener } 716c14b763aSAlp Dener } 717c14b763aSAlp Dener } 718c14b763aSAlp Dener break; 719c14b763aSAlp Dener 720c14b763aSAlp Dener case BNK_BFGS: 721c14b763aSAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS 722c14b763aSAlp Dener Failed to obtain acceptable iterate with BFGS step 723c14b763aSAlp Dener Attempt to use the scaled gradient direction */ 724*9b6ef848SAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 725c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 726c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 727c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 728c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 729c14b763aSAlp Dener 730c14b763aSAlp Dener bfgsUpdates = 1; 731c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 732c14b763aSAlp Dener break; 733c14b763aSAlp Dener 734c14b763aSAlp Dener case BNK_SCALED_GRADIENT: 735c14b763aSAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS 736c14b763aSAlp Dener The scaled gradient step did not produce a new iterate; 737c14b763aSAlp Dener attemp to use the gradient direction. 738c14b763aSAlp Dener Need to make sure we are not using a different diagonal scaling */ 739c14b763aSAlp Dener ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr); 740c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr); 741c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 742c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 743c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 744c14b763aSAlp Dener 745c14b763aSAlp Dener bfgsUpdates = 1; 746c14b763aSAlp Dener stepType = BNK_GRADIENT; 747c14b763aSAlp Dener break; 748c14b763aSAlp Dener } 7498d5ead36SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 7508d5ead36SAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 7512f75a4aaSAlp Dener ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr); 752c14b763aSAlp Dener 7538d5ead36SAlp Dener /* Perform one last line search with the fall-back step */ 754c14b763aSAlp Dener ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr); 755c14b763aSAlp Dener ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 756c14b763aSAlp Dener } 757c14b763aSAlp Dener *reason = ls_reason; 758c14b763aSAlp Dener PetscFunctionReturn(0); 759c14b763aSAlp Dener } 760c14b763aSAlp Dener 761df278d8fSAlp Dener /*------------------------------------------------------------*/ 762df278d8fSAlp Dener 763df278d8fSAlp Dener /* Routine for updating the trust radius. 764df278d8fSAlp Dener 765df278d8fSAlp Dener Function features three different update methods: 766df278d8fSAlp Dener 1) Line-search step length based 767df278d8fSAlp Dener 2) Predicted decrease on the CG quadratic model 768df278d8fSAlp Dener 3) Interpolation 769df278d8fSAlp Dener */ 770df278d8fSAlp Dener 77128017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept) 772080d2917SAlp Dener { 773080d2917SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 774080d2917SAlp Dener PetscErrorCode ierr; 775080d2917SAlp Dener 776b1c2d0e3SAlp Dener PetscReal step, kappa; 777080d2917SAlp Dener PetscReal gdx, tau_1, tau_2, tau_min, tau_max; 778080d2917SAlp Dener 779080d2917SAlp Dener PetscFunctionBegin; 780080d2917SAlp Dener /* Update trust region radius */ 781080d2917SAlp Dener *accept = PETSC_FALSE; 78228017e9fSAlp Dener switch(updateType) { 783080d2917SAlp Dener case BNK_UPDATE_STEP: 784c14b763aSAlp Dener *accept = PETSC_TRUE; /* always accept here because line search succeeded */ 785080d2917SAlp Dener if (stepType == BNK_NEWTON) { 786080d2917SAlp Dener ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr); 787080d2917SAlp Dener if (step < bnk->nu1) { 788080d2917SAlp Dener /* Very bad step taken; reduce radius */ 789080d2917SAlp Dener tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 790080d2917SAlp Dener } else if (step < bnk->nu2) { 791080d2917SAlp Dener /* Reasonably bad step taken; reduce radius */ 792080d2917SAlp Dener tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust); 793080d2917SAlp Dener } else if (step < bnk->nu3) { 794080d2917SAlp Dener /* Reasonable step was taken; leave radius alone */ 795080d2917SAlp Dener if (bnk->omega3 < 1.0) { 796080d2917SAlp Dener tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust); 797080d2917SAlp Dener } else if (bnk->omega3 > 1.0) { 798080d2917SAlp Dener tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust); 799080d2917SAlp Dener } 800080d2917SAlp Dener } else if (step < bnk->nu4) { 801080d2917SAlp Dener /* Full step taken; increase the radius */ 802080d2917SAlp Dener tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust); 803080d2917SAlp Dener } else { 804080d2917SAlp Dener /* More than full step taken; increase the radius */ 805080d2917SAlp Dener tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust); 806080d2917SAlp Dener } 807080d2917SAlp Dener } else { 808080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 809080d2917SAlp Dener tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 810080d2917SAlp Dener } 811080d2917SAlp Dener break; 812080d2917SAlp Dener 813080d2917SAlp Dener case BNK_UPDATE_REDUCTION: 814080d2917SAlp Dener if (stepType == BNK_NEWTON) { 815b1c2d0e3SAlp Dener if (prered < 0.0) { 816fed79b8eSAlp Dener /* The predicted reduction has the wrong sign. This cannot 817fed79b8eSAlp Dener happen in infinite precision arithmetic. Step should 818fed79b8eSAlp Dener be rejected! */ 819080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 820fed79b8eSAlp Dener } 821fed79b8eSAlp Dener else { 822b1c2d0e3SAlp Dener if (PetscIsInfOrNanReal(actred)) { 823080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 824080d2917SAlp Dener } else { 8250a4511e9SAlp Dener if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && 8260a4511e9SAlp Dener (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) { 827080d2917SAlp Dener kappa = 1.0; 828fed79b8eSAlp Dener } 829fed79b8eSAlp Dener else { 830080d2917SAlp Dener kappa = actred / prered; 831080d2917SAlp Dener } 832fed79b8eSAlp Dener 833fed79b8eSAlp Dener /* Accept or reject the step and update radius */ 834080d2917SAlp Dener if (kappa < bnk->eta1) { 835fed79b8eSAlp Dener /* Reject the step */ 836080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 837fed79b8eSAlp Dener } 838fed79b8eSAlp Dener else { 839fed79b8eSAlp Dener /* Accept the step */ 840c133c014SAlp Dener *accept = PETSC_TRUE; 841c133c014SAlp Dener /* Update the trust region radius only if the computed step is at the trust radius boundary */ 8428d5ead36SAlp Dener if (bnk->dnorm == tao->trust) { 843080d2917SAlp Dener if (kappa < bnk->eta2) { 844080d2917SAlp Dener /* Marginal bad step */ 845c133c014SAlp Dener tao->trust = bnk->alpha2 * tao->trust; 846080d2917SAlp Dener } 847fed79b8eSAlp Dener else if (kappa < bnk->eta3) { 848fed79b8eSAlp Dener /* Reasonable step */ 849fed79b8eSAlp Dener tao->trust = bnk->alpha3 * tao->trust; 850fed79b8eSAlp Dener } 851fed79b8eSAlp Dener else if (kappa < bnk->eta4) { 852080d2917SAlp Dener /* Good step */ 853c133c014SAlp Dener tao->trust = bnk->alpha4 * tao->trust; 854fed79b8eSAlp Dener } 855fed79b8eSAlp Dener else { 856080d2917SAlp Dener /* Very good step */ 857c133c014SAlp Dener tao->trust = bnk->alpha5 * tao->trust; 858080d2917SAlp Dener } 859c133c014SAlp Dener } 860080d2917SAlp Dener } 861080d2917SAlp Dener } 862080d2917SAlp Dener } 863080d2917SAlp Dener } else { 864080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 865080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust); 866080d2917SAlp Dener } 867080d2917SAlp Dener break; 868080d2917SAlp Dener 869080d2917SAlp Dener default: 870080d2917SAlp Dener if (stepType == BNK_NEWTON) { 871b1c2d0e3SAlp Dener if (prered < 0.0) { 872080d2917SAlp Dener /* The predicted reduction has the wrong sign. This cannot */ 873080d2917SAlp Dener /* happen in infinite precision arithmetic. Step should */ 874080d2917SAlp Dener /* be rejected! */ 875080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 876080d2917SAlp Dener } else { 877b1c2d0e3SAlp Dener if (PetscIsInfOrNanReal(actred)) { 878080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 879080d2917SAlp Dener } else { 880080d2917SAlp Dener if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 881080d2917SAlp Dener kappa = 1.0; 882080d2917SAlp Dener } else { 883080d2917SAlp Dener kappa = actred / prered; 884080d2917SAlp Dener } 885080d2917SAlp Dener 886080d2917SAlp Dener ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 887080d2917SAlp Dener tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred); 888080d2917SAlp Dener tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred); 889080d2917SAlp Dener tau_min = PetscMin(tau_1, tau_2); 890080d2917SAlp Dener tau_max = PetscMax(tau_1, tau_2); 891080d2917SAlp Dener 892080d2917SAlp Dener if (kappa >= 1.0 - bnk->mu1) { 893080d2917SAlp Dener /* Great agreement */ 894080d2917SAlp Dener *accept = PETSC_TRUE; 895080d2917SAlp Dener if (tau_max < 1.0) { 896080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 897080d2917SAlp Dener } else if (tau_max > bnk->gamma4) { 898080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm); 899080d2917SAlp Dener } else { 900080d2917SAlp Dener tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 901080d2917SAlp Dener } 902080d2917SAlp Dener } else if (kappa >= 1.0 - bnk->mu2) { 903080d2917SAlp Dener /* Good agreement */ 904080d2917SAlp Dener *accept = PETSC_TRUE; 905080d2917SAlp Dener if (tau_max < bnk->gamma2) { 906080d2917SAlp Dener tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 907080d2917SAlp Dener } else if (tau_max > bnk->gamma3) { 908080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 909080d2917SAlp Dener } else if (tau_max < 1.0) { 910080d2917SAlp Dener tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 911080d2917SAlp Dener } else { 912080d2917SAlp Dener tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 913080d2917SAlp Dener } 914080d2917SAlp Dener } else { 915080d2917SAlp Dener /* Not good agreement */ 916080d2917SAlp Dener if (tau_min > 1.0) { 917080d2917SAlp Dener tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 918080d2917SAlp Dener } else if (tau_max < bnk->gamma1) { 919080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 920080d2917SAlp Dener } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) { 921080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 922080d2917SAlp Dener } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) { 923080d2917SAlp Dener tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm); 924080d2917SAlp Dener } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) { 925080d2917SAlp Dener tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm); 926080d2917SAlp Dener } else { 927080d2917SAlp Dener tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 928080d2917SAlp Dener } 929080d2917SAlp Dener } 930080d2917SAlp Dener } 931080d2917SAlp Dener } 932080d2917SAlp Dener } else { 933080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 934080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust); 935080d2917SAlp Dener } 93628017e9fSAlp Dener break; 937080d2917SAlp Dener } 938c133c014SAlp Dener /* Make sure the radius does not violate min and max settings */ 939c133c014SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 940fed79b8eSAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 941080d2917SAlp Dener PetscFunctionReturn(0); 942080d2917SAlp Dener } 943080d2917SAlp Dener 944eb910715SAlp Dener /* ---------------------------------------------------------- */ 945df278d8fSAlp Dener 94662675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType) 94762675beeSAlp Dener { 94862675beeSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 94962675beeSAlp Dener 95062675beeSAlp Dener PetscFunctionBegin; 95162675beeSAlp Dener switch (stepType) { 95262675beeSAlp Dener case BNK_NEWTON: 95362675beeSAlp Dener ++bnk->newt; 95462675beeSAlp Dener break; 95562675beeSAlp Dener case BNK_BFGS: 95662675beeSAlp Dener ++bnk->bfgs; 95762675beeSAlp Dener break; 95862675beeSAlp Dener case BNK_SCALED_GRADIENT: 95962675beeSAlp Dener ++bnk->sgrad; 96062675beeSAlp Dener break; 96162675beeSAlp Dener case BNK_GRADIENT: 96262675beeSAlp Dener ++bnk->grad; 96362675beeSAlp Dener break; 96462675beeSAlp Dener default: 96562675beeSAlp Dener break; 96662675beeSAlp Dener } 96762675beeSAlp Dener PetscFunctionReturn(0); 96862675beeSAlp Dener } 96962675beeSAlp Dener 97062675beeSAlp Dener /* ---------------------------------------------------------- */ 97162675beeSAlp Dener 972*9b6ef848SAlp Dener PetscErrorCode TaoSetUp_BNK(Tao tao) 973eb910715SAlp Dener { 974eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 975eb910715SAlp Dener PetscErrorCode ierr; 976*9b6ef848SAlp Dener KSPType ksp_type; 977eb910715SAlp Dener 978eb910715SAlp Dener PetscFunctionBegin; 979eb910715SAlp Dener if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);} 980eb910715SAlp Dener if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);} 981eb910715SAlp Dener if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);} 982eb910715SAlp Dener if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);} 983eb910715SAlp Dener if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);} 98409164190SAlp Dener if (!bnk->Xwork) {ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);} 98509164190SAlp Dener if (!bnk->Gwork) {ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);} 98609164190SAlp Dener if (!bnk->unprojected_gradient) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);} 98709164190SAlp Dener if (!bnk->unprojected_gradient_old) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);} 988198282dbSAlp Dener if (!bnk->G_inactive) {ierr = VecDuplicate(tao->solution,&bnk->G_inactive);CHKERRQ(ierr);} 989eb910715SAlp Dener bnk->Diag = 0; 99062675beeSAlp Dener bnk->Diag_min = 0; 99162675beeSAlp Dener bnk->Diag_max = 0; 99262675beeSAlp Dener bnk->inactive_work = 0; 99362675beeSAlp Dener bnk->active_work = 0; 99462675beeSAlp Dener bnk->inactive_idx = 0; 99562675beeSAlp Dener bnk->active_idx = 0; 99662675beeSAlp Dener bnk->active_lower = 0; 99762675beeSAlp Dener bnk->active_upper = 0; 99862675beeSAlp Dener bnk->active_fixed = 0; 999eb910715SAlp Dener bnk->M = 0; 100009164190SAlp Dener bnk->H_inactive = 0; 100109164190SAlp Dener bnk->Hpre_inactive = 0; 1002*9b6ef848SAlp Dener ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr); 1003*9b6ef848SAlp Dener ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr); 1004*9b6ef848SAlp Dener ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr); 1005*9b6ef848SAlp Dener ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr); 1006eb910715SAlp Dener PetscFunctionReturn(0); 1007eb910715SAlp Dener } 1008eb910715SAlp Dener 1009eb910715SAlp Dener /*------------------------------------------------------------*/ 1010df278d8fSAlp Dener 1011eb910715SAlp Dener static PetscErrorCode TaoDestroy_BNK(Tao tao) 1012eb910715SAlp Dener { 1013eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 1014eb910715SAlp Dener PetscErrorCode ierr; 1015eb910715SAlp Dener 1016eb910715SAlp Dener PetscFunctionBegin; 1017eb910715SAlp Dener if (tao->setupcalled) { 1018eb910715SAlp Dener ierr = VecDestroy(&bnk->W);CHKERRQ(ierr); 1019eb910715SAlp Dener ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr); 1020eb910715SAlp Dener ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr); 102109164190SAlp Dener ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr); 102209164190SAlp Dener ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr); 102309164190SAlp Dener ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr); 102409164190SAlp Dener ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr); 1025198282dbSAlp Dener ierr = VecDestroy(&bnk->G_inactive);CHKERRQ(ierr); 1026eb910715SAlp Dener } 1027eb910715SAlp Dener ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr); 102862675beeSAlp Dener ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr); 102962675beeSAlp Dener ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr); 1030eb910715SAlp Dener ierr = MatDestroy(&bnk->M);CHKERRQ(ierr); 103162675beeSAlp Dener if (bnk->Hpre_inactive != bnk->H_inactive) {ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);} 103262675beeSAlp Dener ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr); 1033eb910715SAlp Dener ierr = PetscFree(tao->data);CHKERRQ(ierr); 1034eb910715SAlp Dener PetscFunctionReturn(0); 1035eb910715SAlp Dener } 1036eb910715SAlp Dener 1037eb910715SAlp Dener /*------------------------------------------------------------*/ 1038df278d8fSAlp Dener 1039eb910715SAlp Dener static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao) 1040eb910715SAlp Dener { 1041eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 1042eb910715SAlp Dener PetscErrorCode ierr; 1043eb910715SAlp Dener 1044eb910715SAlp Dener PetscFunctionBegin; 1045eb910715SAlp Dener ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr); 10468d5ead36SAlp Dener ierr = PetscOptionsEList("-tao_bnk_pc_type", "pc type", "", BNK_PC, BNK_PC_TYPES, BNK_PC[bnk->pc_type], &bnk->pc_type, 0);CHKERRQ(ierr); 10478d5ead36SAlp Dener ierr = PetscOptionsEList("-tao_bnk_bfgs_scale_type", "bfgs scale type", "", BFGS_SCALE, BFGS_SCALE_TYPES, BFGS_SCALE[bnk->bfgs_scale_type], &bnk->bfgs_scale_type, 0);CHKERRQ(ierr); 10488d5ead36SAlp Dener ierr = PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, 0);CHKERRQ(ierr); 10498d5ead36SAlp Dener ierr = PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, 0);CHKERRQ(ierr); 10502f75a4aaSAlp Dener ierr = PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, 0);CHKERRQ(ierr); 10518d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr); 10528d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr); 10538d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr); 10548d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr); 10558d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr); 10568d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr); 10578d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr); 10588d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr); 10598d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr); 10608d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr); 10618d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr); 10628d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr); 10638d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr); 10648d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr); 10658d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr); 10668d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr); 10678d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr); 10688d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr); 10698d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr); 10708d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr); 10718d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr); 10728d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr); 10738d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr); 10748d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr); 10758d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr); 10768d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr); 10778d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr); 10788d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr); 10798d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr); 10808d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr); 10818d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr); 10828d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr); 10838d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr); 10848d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr); 10858d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr); 10868d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr); 10878d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr); 10888d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr); 10898d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr); 10908d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr); 10918d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr); 10928d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr); 10938d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr); 10948d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr); 10958d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr); 10960a4511e9SAlp Dener ierr = PetscOptionsReal("-tao_bnk_as_tol", "initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr); 10970a4511e9SAlp Dener ierr = PetscOptionsReal("-tao_bnk_as_step", "step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr); 1098eb910715SAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 1099eb910715SAlp Dener ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 1100eb910715SAlp Dener ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 1101eb910715SAlp Dener PetscFunctionReturn(0); 1102eb910715SAlp Dener } 1103eb910715SAlp Dener 1104eb910715SAlp Dener /*------------------------------------------------------------*/ 1105df278d8fSAlp Dener 1106eb910715SAlp Dener static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer) 1107eb910715SAlp Dener { 1108eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 1109eb910715SAlp Dener PetscInt nrejects; 1110eb910715SAlp Dener PetscBool isascii; 1111eb910715SAlp Dener PetscErrorCode ierr; 1112eb910715SAlp Dener 1113eb910715SAlp Dener PetscFunctionBegin; 1114eb910715SAlp Dener ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 1115eb910715SAlp Dener if (isascii) { 1116eb910715SAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1117eb910715SAlp Dener if (BNK_PC_BFGS == bnk->pc_type && bnk->M) { 1118eb910715SAlp Dener ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr); 1119eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr); 1120eb910715SAlp Dener } 1121eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr); 1122eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr); 1123eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr); 1124eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr); 1125eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr); 1126eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr); 1127eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr); 1128eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr); 1129eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr); 1130eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr); 1131eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr); 1132eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr); 1133eb910715SAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1134eb910715SAlp Dener } 1135eb910715SAlp Dener PetscFunctionReturn(0); 1136eb910715SAlp Dener } 1137eb910715SAlp Dener 1138eb910715SAlp Dener /* ---------------------------------------------------------- */ 1139df278d8fSAlp Dener 1140eb910715SAlp Dener /*MC 1141eb910715SAlp Dener TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms. 114266ed3702SAlp Dener At each iteration, the BNK methods solve the symmetric 1143eb910715SAlp Dener system of equations to obtain the step diretion dk: 1144eb910715SAlp Dener Hk dk = -gk 11452b97c8d8SAlp Dener for free variables only. The step can be globalized either through 11462b97c8d8SAlp Dener trust-region methods, or a line search, or a heuristic mixture of both. 1147eb910715SAlp Dener 1148eb910715SAlp Dener Options Database Keys: 11498d5ead36SAlp Dener + -tao_bnk_pc_type - "none","ahess","bfgs","petsc" 11508d5ead36SAlp Dener . -tao_bnk_bfgs_scale_type - "ahess","phess","bfgs" 11518d5ead36SAlp Dener . -tao_bnk_init_type - "constant","direction","interpolation" 11528d5ead36SAlp Dener . -tao_bnk_update_type - "step","direction","interpolation" 11532f75a4aaSAlp Dener . -tao_bnk_as_type - "none","bertsekas" 11548d5ead36SAlp Dener . -tao_bnk_sval - perturbation starting value 11558d5ead36SAlp Dener . -tao_bnk_imin - minimum initial perturbation 11568d5ead36SAlp Dener . -tao_bnk_imax - maximum initial perturbation 11578d5ead36SAlp Dener . -tao_bnk_pmin - minimum perturbation 11588d5ead36SAlp Dener . -tao_bnk_pmax - maximum perturbation 11598d5ead36SAlp Dener . -tao_bnk_pgfac - growth factor 11608d5ead36SAlp Dener . -tao_bnk_psfac - shrink factor 11618d5ead36SAlp Dener . -tao_bnk_imfac - initial merit factor 11628d5ead36SAlp Dener . -tao_bnk_pmgfac - merit growth factor 11638d5ead36SAlp Dener . -tao_bnk_pmsfac - merit shrink factor 11648d5ead36SAlp Dener . -tao_bnk_eta1 - poor steplength; reduce radius 11658d5ead36SAlp Dener . -tao_bnk_eta2 - reasonable steplength; leave radius 11668d5ead36SAlp Dener . -tao_bnk_eta3 - good steplength; increase readius 11678d5ead36SAlp Dener . -tao_bnk_eta4 - excellent steplength; greatly increase radius 11688d5ead36SAlp Dener . -tao_bnk_alpha1 - alpha1 reduction 11698d5ead36SAlp Dener . -tao_bnk_alpha2 - alpha2 reduction 11708d5ead36SAlp Dener . -tao_bnk_alpha3 - alpha3 reduction 11718d5ead36SAlp Dener . -tao_bnk_alpha4 - alpha4 reduction 11728d5ead36SAlp Dener . -tao_bnk_alpha - alpha5 reduction 11738d5ead36SAlp Dener . -tao_bnk_mu1 - mu1 interpolation update 11748d5ead36SAlp Dener . -tao_bnk_mu2 - mu2 interpolation update 11758d5ead36SAlp Dener . -tao_bnk_gamma1 - gamma1 interpolation update 11768d5ead36SAlp Dener . -tao_bnk_gamma2 - gamma2 interpolation update 11778d5ead36SAlp Dener . -tao_bnk_gamma3 - gamma3 interpolation update 11788d5ead36SAlp Dener . -tao_bnk_gamma4 - gamma4 interpolation update 11798d5ead36SAlp Dener . -tao_bnk_theta - theta interpolation update 11808d5ead36SAlp Dener . -tao_bnk_omega1 - omega1 step update 11818d5ead36SAlp Dener . -tao_bnk_omega2 - omega2 step update 11828d5ead36SAlp Dener . -tao_bnk_omega3 - omega3 step update 11838d5ead36SAlp Dener . -tao_bnk_omega4 - omega4 step update 11848d5ead36SAlp Dener . -tao_bnk_omega5 - omega5 step update 11858d5ead36SAlp Dener . -tao_bnk_mu1_i - mu1 interpolation init factor 11868d5ead36SAlp Dener . -tao_bnk_mu2_i - mu2 interpolation init factor 11878d5ead36SAlp Dener . -tao_bnk_gamma1_i - gamma1 interpolation init factor 11888d5ead36SAlp Dener . -tao_bnk_gamma2_i - gamma2 interpolation init factor 11898d5ead36SAlp Dener . -tao_bnk_gamma3_i - gamma3 interpolation init factor 11908d5ead36SAlp Dener . -tao_bnk_gamma4_i - gamma4 interpolation init factor 11912f75a4aaSAlp Dener . -tao_bnk_theta_i - theta interpolation init factor 11922f75a4aaSAlp Dener - -tao_bnk_bound_tol - initial tolerance used in estimating bounded active variables 1193eb910715SAlp Dener 1194eb910715SAlp Dener Level: beginner 1195eb910715SAlp Dener M*/ 1196eb910715SAlp Dener 1197eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao) 1198eb910715SAlp Dener { 1199eb910715SAlp Dener TAO_BNK *bnk; 1200eb910715SAlp Dener const char *morethuente_type = TAOLINESEARCHMT; 1201eb910715SAlp Dener PetscErrorCode ierr; 1202eb910715SAlp Dener 1203eb910715SAlp Dener PetscFunctionBegin; 1204eb910715SAlp Dener ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr); 1205eb910715SAlp Dener 1206eb910715SAlp Dener tao->ops->setup = TaoSetUp_BNK; 1207eb910715SAlp Dener tao->ops->view = TaoView_BNK; 1208eb910715SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BNK; 1209eb910715SAlp Dener tao->ops->destroy = TaoDestroy_BNK; 1210eb910715SAlp Dener 1211eb910715SAlp Dener /* Override default settings (unless already changed) */ 1212eb910715SAlp Dener if (!tao->max_it_changed) tao->max_it = 50; 1213eb910715SAlp Dener if (!tao->trust0_changed) tao->trust0 = 100.0; 1214eb910715SAlp Dener 1215eb910715SAlp Dener tao->data = (void*)bnk; 1216eb910715SAlp Dener 121766ed3702SAlp Dener /* Hessian shifting parameters */ 1218eb910715SAlp Dener bnk->sval = 0.0; 1219eb910715SAlp Dener bnk->imin = 1.0e-4; 1220eb910715SAlp Dener bnk->imax = 1.0e+2; 1221eb910715SAlp Dener bnk->imfac = 1.0e-1; 1222eb910715SAlp Dener 1223eb910715SAlp Dener bnk->pmin = 1.0e-12; 1224eb910715SAlp Dener bnk->pmax = 1.0e+2; 1225eb910715SAlp Dener bnk->pgfac = 1.0e+1; 1226eb910715SAlp Dener bnk->psfac = 4.0e-1; 1227eb910715SAlp Dener bnk->pmgfac = 1.0e-1; 1228eb910715SAlp Dener bnk->pmsfac = 1.0e-1; 1229eb910715SAlp Dener 1230eb910715SAlp Dener /* Default values for trust-region radius update based on steplength */ 1231eb910715SAlp Dener bnk->nu1 = 0.25; 1232eb910715SAlp Dener bnk->nu2 = 0.50; 1233eb910715SAlp Dener bnk->nu3 = 1.00; 1234eb910715SAlp Dener bnk->nu4 = 1.25; 1235eb910715SAlp Dener 1236eb910715SAlp Dener bnk->omega1 = 0.25; 1237eb910715SAlp Dener bnk->omega2 = 0.50; 1238eb910715SAlp Dener bnk->omega3 = 1.00; 1239eb910715SAlp Dener bnk->omega4 = 2.00; 1240eb910715SAlp Dener bnk->omega5 = 4.00; 1241eb910715SAlp Dener 1242eb910715SAlp Dener /* Default values for trust-region radius update based on reduction */ 1243eb910715SAlp Dener bnk->eta1 = 1.0e-4; 1244eb910715SAlp Dener bnk->eta2 = 0.25; 1245eb910715SAlp Dener bnk->eta3 = 0.50; 1246eb910715SAlp Dener bnk->eta4 = 0.90; 1247eb910715SAlp Dener 1248eb910715SAlp Dener bnk->alpha1 = 0.25; 1249eb910715SAlp Dener bnk->alpha2 = 0.50; 1250eb910715SAlp Dener bnk->alpha3 = 1.00; 1251eb910715SAlp Dener bnk->alpha4 = 2.00; 1252eb910715SAlp Dener bnk->alpha5 = 4.00; 1253eb910715SAlp Dener 1254eb910715SAlp Dener /* Default values for trust-region radius update based on interpolation */ 1255eb910715SAlp Dener bnk->mu1 = 0.10; 1256eb910715SAlp Dener bnk->mu2 = 0.50; 1257eb910715SAlp Dener 1258eb910715SAlp Dener bnk->gamma1 = 0.25; 1259eb910715SAlp Dener bnk->gamma2 = 0.50; 1260eb910715SAlp Dener bnk->gamma3 = 2.00; 1261eb910715SAlp Dener bnk->gamma4 = 4.00; 1262eb910715SAlp Dener 1263eb910715SAlp Dener bnk->theta = 0.05; 1264eb910715SAlp Dener 1265eb910715SAlp Dener /* Default values for trust region initialization based on interpolation */ 1266eb910715SAlp Dener bnk->mu1_i = 0.35; 1267eb910715SAlp Dener bnk->mu2_i = 0.50; 1268eb910715SAlp Dener 1269eb910715SAlp Dener bnk->gamma1_i = 0.0625; 1270eb910715SAlp Dener bnk->gamma2_i = 0.5; 1271eb910715SAlp Dener bnk->gamma3_i = 2.0; 1272eb910715SAlp Dener bnk->gamma4_i = 5.0; 1273eb910715SAlp Dener 1274eb910715SAlp Dener bnk->theta_i = 0.25; 1275eb910715SAlp Dener 1276eb910715SAlp Dener /* Remaining parameters */ 1277eb910715SAlp Dener bnk->min_radius = 1.0e-10; 1278eb910715SAlp Dener bnk->max_radius = 1.0e10; 1279770b7498SAlp Dener bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); 12800a4511e9SAlp Dener bnk->as_tol = 1.0e-3; 12810a4511e9SAlp Dener bnk->as_step = 1.0e-3; 128262675beeSAlp Dener bnk->dmin = 1.0e-6; 128362675beeSAlp Dener bnk->dmax = 1.0e6; 1284eb910715SAlp Dener 1285eb910715SAlp Dener bnk->pc_type = BNK_PC_BFGS; 1286eb910715SAlp Dener bnk->bfgs_scale_type = BFGS_SCALE_PHESS; 1287eb910715SAlp Dener bnk->init_type = BNK_INIT_INTERPOLATION; 1288fed79b8eSAlp Dener bnk->update_type = BNK_UPDATE_INTERPOLATION; 12892f75a4aaSAlp Dener bnk->as_type = BNK_AS_BERTSEKAS; 1290eb910715SAlp Dener 1291eb910715SAlp Dener ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 1292eb910715SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 1293eb910715SAlp Dener ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 1294eb910715SAlp Dener ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 1295eb910715SAlp Dener ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 1296eb910715SAlp Dener 1297eb910715SAlp Dener /* Set linear solver to default for symmetric matrices */ 1298eb910715SAlp Dener ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 1299eb910715SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 1300eb910715SAlp Dener ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr); 1301eb910715SAlp Dener ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr); 1302eb910715SAlp Dener PetscFunctionReturn(0); 1303eb910715SAlp Dener } 1304