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); 18*5e9b73cbSAlp Dener ierr = MatLMVMSolve(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);} 96*5e9b73cbSAlp Dener } 9762675beeSAlp Dener ierr = VecSet(bnk->Diag_min, bnk->dmin);CHKERRQ(ierr); 9862675beeSAlp Dener ierr = VecSet(bnk->Diag_max, bnk->dmax);CHKERRQ(ierr); 99eb910715SAlp Dener 100eb910715SAlp Dener /* Modify the preconditioner to use the bfgs approximation */ 101eb910715SAlp Dener ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 102eb910715SAlp Dener switch(bnk->pc_type) { 103eb910715SAlp Dener case BNK_PC_NONE: 104eb910715SAlp Dener ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr); 105eb910715SAlp Dener ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 106eb910715SAlp Dener break; 107eb910715SAlp Dener 108eb910715SAlp Dener case BNK_PC_AHESS: 109eb910715SAlp Dener ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr); 110eb910715SAlp Dener ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 111eb910715SAlp Dener ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 112eb910715SAlp Dener break; 113eb910715SAlp Dener 114eb910715SAlp Dener case BNK_PC_BFGS: 115eb910715SAlp Dener ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr); 116eb910715SAlp Dener ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 117eb910715SAlp Dener ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr); 118eb910715SAlp Dener ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr); 119eb910715SAlp Dener ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr); 120eb910715SAlp Dener break; 121eb910715SAlp Dener 122eb910715SAlp Dener default: 123eb910715SAlp Dener /* Use the pc method set by pc_type */ 124eb910715SAlp Dener break; 125eb910715SAlp Dener } 126eb910715SAlp Dener 127eb910715SAlp Dener /* Initialize trust-region radius. The initialization is only performed 128eb910715SAlp Dener when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */ 129eb910715SAlp Dener if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 13062675beeSAlp Dener switch(initType) { 131eb910715SAlp Dener case BNK_INIT_CONSTANT: 132eb910715SAlp Dener /* Use the initial radius specified */ 133eb910715SAlp Dener break; 134eb910715SAlp Dener 135eb910715SAlp Dener case BNK_INIT_INTERPOLATION: 136eb910715SAlp Dener /* Use the initial radius specified */ 137eb910715SAlp Dener max_radius = 0.0; 138eb910715SAlp Dener 139eb910715SAlp Dener for (j = 0; j < j_max; ++j) { 1400a4511e9SAlp Dener f_min = bnk->f; 141eb910715SAlp Dener sigma = 0.0; 142eb910715SAlp Dener 143eb910715SAlp Dener if (needH) { 14462602cfbSAlp Dener /* Compute the Hessian at the new step, and extract the inactive subsystem */ 145eb910715SAlp Dener ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 14628017e9fSAlp Dener if (bnk->inactive_idx) { 14762602cfbSAlp Dener ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr); 14862602cfbSAlp Dener ierr = VecWhichInactive(tao->XL,tao->solution,bnk->unprojected_gradient,tao->XU,PETSC_TRUE,&bnk->inactive_idx);CHKERRQ(ierr); 14962602cfbSAlp Dener ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr); 15028017e9fSAlp Dener } else { 15162675beeSAlp Dener ierr = MatDestroy(&bnk->H_inactive); 15262675beeSAlp Dener ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive); 15328017e9fSAlp Dener } 154eb910715SAlp Dener needH = 0; 155eb910715SAlp Dener } 156eb910715SAlp Dener 157eb910715SAlp Dener for (i = 0; i < i_max; ++i) { 15862602cfbSAlp Dener /* Take a steepest descent step and snap it to bounds */ 15962602cfbSAlp Dener ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 16062602cfbSAlp Dener ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr); 16162602cfbSAlp Dener ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 162770b7498SAlp Dener /* Recompute the step after bound snapping so that it can be used in predicted decrease calculation later */ 163eb910715SAlp Dener ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr); 16462602cfbSAlp Dener ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr); 16562602cfbSAlp Dener /* Compute the objective at the trial */ 16662602cfbSAlp Dener ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr); 16762602cfbSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 168eb910715SAlp Dener if (PetscIsInfOrNanReal(ftrial)) { 169eb910715SAlp Dener tau = bnk->gamma1_i; 170eb910715SAlp Dener } else { 1710a4511e9SAlp Dener if (ftrial < f_min) { 1720a4511e9SAlp Dener f_min = ftrial; 173eb910715SAlp Dener sigma = -tao->trust / bnk->gnorm; 174eb910715SAlp Dener } 175770b7498SAlp Dener /* Compute the predicted and actual reduction */ 17662602cfbSAlp Dener ierr = MatMult(bnk->H_inactive, tao->gradient, bnk->W);CHKERRQ(ierr); 17762602cfbSAlp Dener ierr = VecDot(tao->gradient, bnk->W, &prered);CHKERRQ(ierr); 178eb910715SAlp Dener prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm)); 179eb910715SAlp Dener actred = bnk->f - ftrial; 180eb910715SAlp Dener if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 181eb910715SAlp Dener kappa = 1.0; 182eb910715SAlp Dener } else { 183eb910715SAlp Dener kappa = actred / prered; 184eb910715SAlp Dener } 185eb910715SAlp Dener 186eb910715SAlp Dener tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred); 187eb910715SAlp Dener tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred); 188eb910715SAlp Dener tau_min = PetscMin(tau_1, tau_2); 189eb910715SAlp Dener tau_max = PetscMax(tau_1, tau_2); 190eb910715SAlp Dener 191eb910715SAlp Dener if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) { 192eb910715SAlp Dener /* Great agreement */ 193eb910715SAlp Dener max_radius = PetscMax(max_radius, tao->trust); 194eb910715SAlp Dener 195eb910715SAlp Dener if (tau_max < 1.0) { 196eb910715SAlp Dener tau = bnk->gamma3_i; 197eb910715SAlp Dener } else if (tau_max > bnk->gamma4_i) { 198eb910715SAlp Dener tau = bnk->gamma4_i; 199eb910715SAlp Dener } else if (tau_1 >= 1.0 && tau_1 <= bnk->gamma4_i && tau_2 < 1.0) { 200eb910715SAlp Dener tau = tau_1; 201eb910715SAlp Dener } else if (tau_2 >= 1.0 && tau_2 <= bnk->gamma4_i && tau_1 < 1.0) { 202eb910715SAlp Dener tau = tau_2; 203eb910715SAlp Dener } else { 204eb910715SAlp Dener tau = tau_max; 205eb910715SAlp Dener } 206eb910715SAlp Dener } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) { 207eb910715SAlp Dener /* Good agreement */ 208eb910715SAlp Dener max_radius = PetscMax(max_radius, tao->trust); 209eb910715SAlp Dener 210eb910715SAlp Dener if (tau_max < bnk->gamma2_i) { 211eb910715SAlp Dener tau = bnk->gamma2_i; 212eb910715SAlp Dener } else if (tau_max > bnk->gamma3_i) { 213eb910715SAlp Dener tau = bnk->gamma3_i; 214eb910715SAlp Dener } else { 215eb910715SAlp Dener tau = tau_max; 216eb910715SAlp Dener } 217eb910715SAlp Dener } else { 218eb910715SAlp Dener /* Not good agreement */ 219eb910715SAlp Dener if (tau_min > 1.0) { 220eb910715SAlp Dener tau = bnk->gamma2_i; 221eb910715SAlp Dener } else if (tau_max < bnk->gamma1_i) { 222eb910715SAlp Dener tau = bnk->gamma1_i; 223eb910715SAlp Dener } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) { 224eb910715SAlp Dener tau = bnk->gamma1_i; 225eb910715SAlp Dener } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 226eb910715SAlp Dener tau = tau_1; 227eb910715SAlp Dener } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 228eb910715SAlp Dener tau = tau_2; 229eb910715SAlp Dener } else { 230eb910715SAlp Dener tau = tau_max; 231eb910715SAlp Dener } 232eb910715SAlp Dener } 233eb910715SAlp Dener } 234eb910715SAlp Dener tao->trust = tau * tao->trust; 235eb910715SAlp Dener } 236eb910715SAlp Dener 2370a4511e9SAlp Dener if (f_min < bnk->f) { 2380a4511e9SAlp Dener bnk->f = f_min; 239eb910715SAlp Dener ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr); 24087602d1eSAlp Dener ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 24109164190SAlp Dener ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr); 24209164190SAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 243eb910715SAlp Dener 244eb910715SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 245eb910715SAlp Dener if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN"); 246eb910715SAlp Dener needH = 1; 247eb910715SAlp Dener 248eb910715SAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 24928017e9fSAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,1.0);CHKERRQ(ierr); 250eb910715SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 251eb910715SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 252eb910715SAlp Dener } 253eb910715SAlp Dener } 254eb910715SAlp Dener tao->trust = PetscMax(tao->trust, max_radius); 255eb910715SAlp Dener 256eb910715SAlp Dener /* Modify the radius if it is too large or small */ 257eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 258eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 259eb910715SAlp Dener break; 260eb910715SAlp Dener 261eb910715SAlp Dener default: 262eb910715SAlp Dener /* Norm of the first direction will initialize radius */ 263eb910715SAlp Dener tao->trust = 0.0; 264eb910715SAlp Dener break; 265eb910715SAlp Dener } 266eb910715SAlp Dener } 267eb910715SAlp Dener 268eb910715SAlp Dener /* Set initial scaling for the BFGS preconditioner 269eb910715SAlp Dener This step is done after computing the initial trust-region radius 270eb910715SAlp Dener since the function value may have decreased */ 271eb910715SAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { 2729b6ef848SAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 273eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr); 274eb910715SAlp Dener } 275eb910715SAlp Dener 276eb910715SAlp Dener /* Set counter for gradient/reset steps*/ 277eb910715SAlp Dener bnk->newt = 0; 278eb910715SAlp Dener bnk->bfgs = 0; 279eb910715SAlp Dener bnk->sgrad = 0; 280eb910715SAlp Dener bnk->grad = 0; 281eb910715SAlp Dener PetscFunctionReturn(0); 282eb910715SAlp Dener } 283eb910715SAlp Dener 284df278d8fSAlp Dener /*------------------------------------------------------------*/ 285df278d8fSAlp Dener 28662675beeSAlp Dener /* Routine for computing the Hessian and preparing the preconditioner at the new iterate */ 28762675beeSAlp Dener 28862675beeSAlp Dener PetscErrorCode TaoBNKComputeHessian(Tao tao) 28962675beeSAlp Dener { 29062675beeSAlp Dener PetscErrorCode ierr; 29162675beeSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 29262675beeSAlp Dener 29362675beeSAlp Dener PetscFunctionBegin; 29462675beeSAlp Dener /* Compute the Hessian */ 29562675beeSAlp Dener ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 29662675beeSAlp Dener /* Add a correction to the BFGS preconditioner */ 29762675beeSAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { 29862675beeSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 29962675beeSAlp Dener /* Update the BFGS diagonal scaling */ 30062675beeSAlp Dener if (BFGS_SCALE_AHESS == bnk->bfgs_scale_type) { 30162675beeSAlp Dener ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr); 30262675beeSAlp Dener ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 30362675beeSAlp Dener ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr); 30462675beeSAlp Dener ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 30562675beeSAlp Dener ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 30662675beeSAlp Dener } 30762675beeSAlp Dener } 30862675beeSAlp Dener PetscFunctionReturn(0); 30962675beeSAlp Dener } 31062675beeSAlp Dener 31162675beeSAlp Dener /*------------------------------------------------------------*/ 31262675beeSAlp Dener 3132f75a4aaSAlp Dener /* Routine for estimating the active set */ 3142f75a4aaSAlp Dener 3152f75a4aaSAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao) 3162f75a4aaSAlp Dener { 3172f75a4aaSAlp Dener PetscErrorCode ierr; 3182f75a4aaSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 3192f75a4aaSAlp Dener 3202f75a4aaSAlp Dener PetscFunctionBegin; 3212f75a4aaSAlp Dener switch (bnk->as_type) { 3222f75a4aaSAlp Dener case BNK_AS_NONE: 3232f75a4aaSAlp Dener ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr); 3242f75a4aaSAlp Dener ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr); 3252f75a4aaSAlp Dener ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr); 3262f75a4aaSAlp Dener ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr); 3272f75a4aaSAlp Dener break; 3282f75a4aaSAlp Dener 3292f75a4aaSAlp Dener case BNK_AS_BERTSEKAS: 3302f75a4aaSAlp Dener /* Compute the trial step vector with which we will estimate the active set at the next iteration */ 3312f75a4aaSAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { 3322f75a4aaSAlp Dener /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */ 333*5e9b73cbSAlp Dener ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr); 3342f75a4aaSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr); 3352f75a4aaSAlp Dener } else { 3369b6ef848SAlp Dener /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */ 3372f75a4aaSAlp Dener ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr); 3389b6ef848SAlp Dener ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr); 3399b6ef848SAlp Dener ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr); 3402f75a4aaSAlp Dener ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr); 3412f75a4aaSAlp Dener ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr); 3422f75a4aaSAlp Dener } 3432f75a4aaSAlp Dener ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr); 3440a4511e9SAlp 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); 3452f75a4aaSAlp Dener 3462f75a4aaSAlp Dener default: 3472f75a4aaSAlp Dener break; 3482f75a4aaSAlp Dener } 3492f75a4aaSAlp Dener PetscFunctionReturn(0); 3502f75a4aaSAlp Dener } 3512f75a4aaSAlp Dener 3522f75a4aaSAlp Dener /*------------------------------------------------------------*/ 3532f75a4aaSAlp Dener 3542f75a4aaSAlp Dener /* Routine for bounding the step direction */ 3552f75a4aaSAlp Dener 3562f75a4aaSAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, Vec step) 3572f75a4aaSAlp Dener { 3582f75a4aaSAlp Dener PetscErrorCode ierr; 3592f75a4aaSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 3602f75a4aaSAlp Dener 3612f75a4aaSAlp Dener PetscFunctionBegin; 36228017e9fSAlp Dener if (bnk->active_idx) { 3632f75a4aaSAlp Dener switch (bnk->as_type) { 3642f75a4aaSAlp Dener case BNK_AS_NONE: 3652f75a4aaSAlp Dener if (bnk->active_idx) { 3662f75a4aaSAlp Dener ierr = VecGetSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr); 3672f75a4aaSAlp Dener ierr = VecSet(bnk->active_work, 0.0);CHKERRQ(ierr); 3682f75a4aaSAlp Dener ierr = VecRestoreSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr); 3692f75a4aaSAlp Dener } 3702f75a4aaSAlp Dener break; 3712f75a4aaSAlp Dener 3722f75a4aaSAlp Dener case BNK_AS_BERTSEKAS: 3732f75a4aaSAlp Dener ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, step);CHKERRQ(ierr); 3742f75a4aaSAlp Dener break; 3752f75a4aaSAlp Dener 3762f75a4aaSAlp Dener default: 3772f75a4aaSAlp Dener break; 3782f75a4aaSAlp Dener } 379e465cd6fSAlp Dener } 3802f75a4aaSAlp Dener PetscFunctionReturn(0); 3812f75a4aaSAlp Dener } 3822f75a4aaSAlp Dener 3832f75a4aaSAlp Dener /*------------------------------------------------------------*/ 3842f75a4aaSAlp Dener 385df278d8fSAlp Dener /* Routine for computing the Newton step. 386df278d8fSAlp Dener 387df278d8fSAlp Dener If the safeguard is enabled, the Newton step is verified to be a 388df278d8fSAlp Dener descent direction, with fallbacks onto BFGS, scaled gradient, and unscaled 389df278d8fSAlp Dener gradient steps if/when necessary. 390df278d8fSAlp Dener 391df278d8fSAlp Dener The function reports back on which type of step has ultimately been stored 392df278d8fSAlp Dener under tao->stepdirection. 393df278d8fSAlp Dener */ 394df278d8fSAlp Dener 39562675beeSAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason) 396eb910715SAlp Dener { 397eb910715SAlp Dener PetscErrorCode ierr; 398eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 399eb910715SAlp Dener 400e465cd6fSAlp Dener PetscReal delta; 401eb910715SAlp Dener PetscInt bfgsUpdates = 0; 402eb910715SAlp Dener PetscInt kspits; 403eb910715SAlp Dener 404eb910715SAlp Dener PetscFunctionBegin; 4052f75a4aaSAlp Dener /* Determine the active and inactive sets */ 4062f75a4aaSAlp Dener ierr = TaoBNKEstimateActiveSet(tao);CHKERRQ(ierr); 40709164190SAlp Dener 408*5e9b73cbSAlp Dener /* Prepare the reduced sub-matrices for the inactive set */ 4092f75a4aaSAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr); } 4102f75a4aaSAlp Dener if (bnk->inactive_idx) { 411*5e9b73cbSAlp Dener ierr = MatDestroy(&bnk->H_inactive); 412*5e9b73cbSAlp Dener ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr); 413eb3ba6a7SAlp Dener if (tao->hessian == tao->hessian_pre) { 414eb3ba6a7SAlp Dener bnk->Hpre_inactive = bnk->H_inactive; 415eb3ba6a7SAlp Dener } else { 416*5e9b73cbSAlp Dener ierr = MatDestroy(&bnk->Hpre_inactive); 417*5e9b73cbSAlp Dener ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr); 418eb3ba6a7SAlp Dener } 4192f75a4aaSAlp Dener } else { 42062675beeSAlp Dener ierr = MatDestroy(&bnk->H_inactive); 42162675beeSAlp Dener ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive); 42262675beeSAlp Dener if (tao->hessian == tao->hessian_pre) { 42362675beeSAlp Dener bnk->Hpre_inactive = bnk->H_inactive; 42462675beeSAlp Dener } else { 42562675beeSAlp Dener ierr = MatDestroy(&bnk->Hpre_inactive); 42662675beeSAlp Dener ierr = MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive); 42762675beeSAlp Dener } 42862675beeSAlp Dener } 42962675beeSAlp Dener 43062675beeSAlp Dener /* Shift the reduced Hessian matrix */ 43162675beeSAlp Dener if ((shift) && (bnk->pert > 0)) { 43262675beeSAlp Dener ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr); 43362675beeSAlp Dener if (bnk->H_inactive != bnk->Hpre_inactive) { 43462675beeSAlp Dener ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr); 43562675beeSAlp Dener } 43662675beeSAlp Dener } 43762675beeSAlp Dener 43862675beeSAlp Dener /* Update the diagonal scaling for the BFGS preconditioner, this time with the Hessian perturbation */ 43962675beeSAlp Dener if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_PHESS == bnk->bfgs_scale_type)) { 44062675beeSAlp Dener /* Obtain diagonal for the bfgs preconditioner */ 441*5e9b73cbSAlp Dener ierr = VecSet(bnk->Diag, 1.0);CHKERRQ(ierr); 442*5e9b73cbSAlp Dener if (bnk->inactive_idx) { 443*5e9b73cbSAlp Dener ierr = VecGetSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr); 444*5e9b73cbSAlp Dener } else { 445*5e9b73cbSAlp Dener bnk->Diag_red = bnk->Diag; 446*5e9b73cbSAlp Dener } 447*5e9b73cbSAlp Dener ierr = MatGetDiagonal(bnk->H_inactive, bnk->Diag_red);CHKERRQ(ierr); 448*5e9b73cbSAlp Dener if (bnk->inactive_idx) { 449*5e9b73cbSAlp Dener ierr = VecRestoreSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr); 450*5e9b73cbSAlp Dener } 45162675beeSAlp Dener ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 45262675beeSAlp Dener ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr); 45362675beeSAlp Dener ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 45462675beeSAlp Dener ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 4552f75a4aaSAlp Dener } 45609164190SAlp Dener 457eb910715SAlp Dener /* Solve the Newton system of equations */ 4582f75a4aaSAlp Dener ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr); 459*5e9b73cbSAlp Dener ierr = KSPReset(tao->ksp);CHKERRQ(ierr); 46009164190SAlp Dener ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr); 461*5e9b73cbSAlp Dener ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr); 462*5e9b73cbSAlp Dener if (bnk->inactive_idx) { 463*5e9b73cbSAlp Dener ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr); 464*5e9b73cbSAlp Dener ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr); 465*5e9b73cbSAlp Dener } else { 466*5e9b73cbSAlp Dener bnk->G_inactive = bnk->unprojected_gradient; 467*5e9b73cbSAlp Dener bnk->X_inactive = tao->stepdirection; 46828017e9fSAlp Dener } 469eb910715SAlp Dener if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 470fed79b8eSAlp Dener ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 471*5e9b73cbSAlp Dener ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr); 472eb910715SAlp Dener ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 473eb910715SAlp Dener tao->ksp_its+=kspits; 474eb910715SAlp Dener tao->ksp_tot_its+=kspits; 475080d2917SAlp Dener ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr); 476eb910715SAlp Dener 477eb910715SAlp Dener if (0.0 == tao->trust) { 478eb910715SAlp Dener /* Radius was uninitialized; use the norm of the direction */ 479080d2917SAlp Dener if (bnk->dnorm > 0.0) { 480080d2917SAlp Dener tao->trust = bnk->dnorm; 481eb910715SAlp Dener 482eb910715SAlp Dener /* Modify the radius if it is too large or small */ 483eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 484eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 485eb910715SAlp Dener } else { 486eb910715SAlp Dener /* The direction was bad; set radius to default value and re-solve 487eb910715SAlp Dener the trust-region subproblem to get a direction */ 488eb910715SAlp Dener tao->trust = tao->trust0; 489eb910715SAlp Dener 490eb910715SAlp Dener /* Modify the radius if it is too large or small */ 491eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 492eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 493eb910715SAlp Dener 494fed79b8eSAlp Dener ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 495*5e9b73cbSAlp Dener ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr); 496eb910715SAlp Dener ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 497eb910715SAlp Dener tao->ksp_its+=kspits; 498eb910715SAlp Dener tao->ksp_tot_its+=kspits; 499080d2917SAlp Dener ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr); 500eb910715SAlp Dener 501080d2917SAlp Dener if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 502eb910715SAlp Dener } 503eb910715SAlp Dener } 504eb910715SAlp Dener } else { 505*5e9b73cbSAlp Dener ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr); 506eb910715SAlp Dener ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr); 507eb910715SAlp Dener tao->ksp_its += kspits; 508eb910715SAlp Dener tao->ksp_tot_its+=kspits; 509eb910715SAlp Dener } 510*5e9b73cbSAlp Dener /* Restore sub vectors back */ 511*5e9b73cbSAlp Dener if (bnk->inactive_idx) { 512*5e9b73cbSAlp Dener ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr); 513*5e9b73cbSAlp Dener ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr); 514*5e9b73cbSAlp Dener } 515770b7498SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 516fed79b8eSAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 5172f75a4aaSAlp Dener ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr); 518770b7498SAlp Dener 519770b7498SAlp Dener /* Record convergence reasons */ 520e465cd6fSAlp Dener ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr); 521e465cd6fSAlp Dener if (KSP_CONVERGED_ATOL == *ksp_reason) { 522770b7498SAlp Dener ++bnk->ksp_atol; 523e465cd6fSAlp Dener } else if (KSP_CONVERGED_RTOL == *ksp_reason) { 524770b7498SAlp Dener ++bnk->ksp_rtol; 525e465cd6fSAlp Dener } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) { 526770b7498SAlp Dener ++bnk->ksp_ctol; 527e465cd6fSAlp Dener } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) { 528770b7498SAlp Dener ++bnk->ksp_negc; 529e465cd6fSAlp Dener } else if (KSP_DIVERGED_DTOL == *ksp_reason) { 530770b7498SAlp Dener ++bnk->ksp_dtol; 531e465cd6fSAlp Dener } else if (KSP_DIVERGED_ITS == *ksp_reason) { 532770b7498SAlp Dener ++bnk->ksp_iter; 533770b7498SAlp Dener } else { 534770b7498SAlp Dener ++bnk->ksp_othr; 535770b7498SAlp Dener } 536fed79b8eSAlp Dener 537fed79b8eSAlp Dener /* Make sure the BFGS preconditioner is healthy */ 538fed79b8eSAlp Dener if (bnk->pc_type == BNK_PC_BFGS) { 539fed79b8eSAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 540e465cd6fSAlp Dener if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 1)) { 541fed79b8eSAlp Dener /* Preconditioner is numerically indefinite; reset the approximation. */ 5429b6ef848SAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 543eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr); 544eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 54509164190SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 546eb910715SAlp Dener } 547fed79b8eSAlp Dener } 548e465cd6fSAlp Dener PetscFunctionReturn(0); 549e465cd6fSAlp Dener } 550eb910715SAlp Dener 55162675beeSAlp Dener /*------------------------------------------------------------*/ 55262675beeSAlp Dener 553*5e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */ 554*5e9b73cbSAlp Dener 555*5e9b73cbSAlp Dener PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered) 556*5e9b73cbSAlp Dener { 557*5e9b73cbSAlp Dener PetscErrorCode ierr; 558*5e9b73cbSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 559*5e9b73cbSAlp Dener 560*5e9b73cbSAlp Dener PetscFunctionBegin; 561*5e9b73cbSAlp Dener /* Extract subvectors associated with the inactive set */ 562*5e9b73cbSAlp Dener if (bnk->inactive_idx){ 563*5e9b73cbSAlp Dener ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr); 564*5e9b73cbSAlp Dener ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr); 565*5e9b73cbSAlp Dener ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr); 566*5e9b73cbSAlp Dener } else { 567*5e9b73cbSAlp Dener bnk->X_inactive = tao->stepdirection; 568*5e9b73cbSAlp Dener bnk->inactive_work = bnk->Xwork; 569*5e9b73cbSAlp Dener bnk->G_inactive = bnk->Gwork; 570*5e9b73cbSAlp Dener } 571*5e9b73cbSAlp Dener /* Recompute the predicted decrease based on the quadratic model */ 572*5e9b73cbSAlp Dener ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr); 573*5e9b73cbSAlp Dener ierr = VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);CHKERRQ(ierr); 574*5e9b73cbSAlp Dener ierr = VecDot(bnk->inactive_work, bnk->X_inactive, prered); 575*5e9b73cbSAlp Dener /* Restore the sub vectors */ 576*5e9b73cbSAlp Dener if (bnk->inactive_idx){ 577*5e9b73cbSAlp Dener ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr); 578*5e9b73cbSAlp Dener ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr); 579*5e9b73cbSAlp Dener ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr); 580*5e9b73cbSAlp Dener } 581*5e9b73cbSAlp Dener PetscFunctionReturn(0); 582*5e9b73cbSAlp Dener } 583*5e9b73cbSAlp Dener 584*5e9b73cbSAlp Dener /*------------------------------------------------------------*/ 585*5e9b73cbSAlp Dener 58662675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction. 58762675beeSAlp Dener 58862675beeSAlp Dener The step direction falls back onto BFGS, scaled gradient and gradient steps 58962675beeSAlp Dener in the event that the Newton step fails the test. 59062675beeSAlp Dener */ 59162675beeSAlp Dener 592e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType) 593e465cd6fSAlp Dener { 594e465cd6fSAlp Dener PetscErrorCode ierr; 595e465cd6fSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 596e465cd6fSAlp Dener 597e465cd6fSAlp Dener PetscReal gdx, delta, e_min; 598e465cd6fSAlp Dener PetscInt bfgsUpdates; 599e465cd6fSAlp Dener 600e465cd6fSAlp Dener PetscFunctionBegin; 601080d2917SAlp Dener ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 602eb910715SAlp Dener if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 603eb910715SAlp Dener /* Newton step is not descent or direction produced Inf or NaN 604eb910715SAlp Dener Update the perturbation for next time */ 605eb910715SAlp Dener if (bnk->pert <= 0.0) { 606eb910715SAlp Dener /* Initialize the perturbation */ 607eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 608eb910715SAlp Dener if (bnk->is_gltr) { 609eb910715SAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 610eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 611eb910715SAlp Dener } 612eb910715SAlp Dener } else { 613eb910715SAlp Dener /* Increase the perturbation */ 614eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 615eb910715SAlp Dener } 616eb910715SAlp Dener 617eb910715SAlp Dener if (BNK_PC_BFGS != bnk->pc_type) { 618eb910715SAlp Dener /* We don't have the bfgs matrix around and updated 619eb910715SAlp Dener Must use gradient direction in this case */ 620080d2917SAlp Dener ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 621eb910715SAlp Dener *stepType = BNK_GRADIENT; 622eb910715SAlp Dener } else { 623eb910715SAlp Dener /* Attempt to use the BFGS direction */ 62409164190SAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 625eb910715SAlp Dener 6268d5ead36SAlp Dener /* Check for success (descent direction) 6278d5ead36SAlp Dener NOTE: Negative gdx here means not a descent direction because 6288d5ead36SAlp Dener the fall-back step is missing a negative sign. */ 629080d2917SAlp Dener ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 6308d5ead36SAlp Dener if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 631eb910715SAlp Dener /* BFGS direction is not descent or direction produced not a number 632eb910715SAlp Dener We can assert bfgsUpdates > 1 in this case because 633eb910715SAlp Dener the first solve produces the scaled gradient direction, 634eb910715SAlp Dener which is guaranteed to be descent */ 635eb910715SAlp Dener 636eb910715SAlp Dener /* Use steepest descent direction (scaled) */ 6379b6ef848SAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 638eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 639eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 64009164190SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 64109164190SAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 642eb910715SAlp Dener 643eb910715SAlp Dener *stepType = BNK_SCALED_GRADIENT; 644eb910715SAlp Dener } else { 645770b7498SAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 646eb910715SAlp Dener if (1 == bfgsUpdates) { 647eb910715SAlp Dener /* The first BFGS direction is always the scaled gradient */ 648eb910715SAlp Dener *stepType = BNK_SCALED_GRADIENT; 649eb910715SAlp Dener } else { 650eb910715SAlp Dener *stepType = BNK_BFGS; 651eb910715SAlp Dener } 652eb910715SAlp Dener } 653eb910715SAlp Dener } 6548d5ead36SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 6558d5ead36SAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 6562f75a4aaSAlp Dener ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr); 657eb910715SAlp Dener } else { 658eb910715SAlp Dener /* Computed Newton step is descent */ 659eb910715SAlp Dener switch (ksp_reason) { 660eb910715SAlp Dener case KSP_DIVERGED_NANORINF: 661eb910715SAlp Dener case KSP_DIVERGED_BREAKDOWN: 662eb910715SAlp Dener case KSP_DIVERGED_INDEFINITE_MAT: 663eb910715SAlp Dener case KSP_DIVERGED_INDEFINITE_PC: 664eb910715SAlp Dener case KSP_CONVERGED_CG_NEG_CURVE: 665eb910715SAlp Dener /* Matrix or preconditioner is indefinite; increase perturbation */ 666eb910715SAlp Dener if (bnk->pert <= 0.0) { 667eb910715SAlp Dener /* Initialize the perturbation */ 668eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 669eb910715SAlp Dener if (bnk->is_gltr) { 670eb910715SAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr); 671eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 672eb910715SAlp Dener } 673eb910715SAlp Dener } else { 674eb910715SAlp Dener /* Increase the perturbation */ 675eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 676eb910715SAlp Dener } 677eb910715SAlp Dener break; 678eb910715SAlp Dener 679eb910715SAlp Dener default: 680eb910715SAlp Dener /* Newton step computation is good; decrease perturbation */ 681eb910715SAlp Dener bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm); 682eb910715SAlp Dener if (bnk->pert < bnk->pmin) { 683eb910715SAlp Dener bnk->pert = 0.0; 684eb910715SAlp Dener } 685eb910715SAlp Dener break; 686eb910715SAlp Dener } 687fed79b8eSAlp Dener *stepType = BNK_NEWTON; 688eb910715SAlp Dener } 689eb910715SAlp Dener PetscFunctionReturn(0); 690eb910715SAlp Dener } 691eb910715SAlp Dener 692df278d8fSAlp Dener /*------------------------------------------------------------*/ 693df278d8fSAlp Dener 694df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search. 695df278d8fSAlp Dener 696df278d8fSAlp Dener Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the 697df278d8fSAlp Dener Newton step does not produce a valid step length. 698df278d8fSAlp Dener */ 699df278d8fSAlp Dener 700c14b763aSAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason) 701c14b763aSAlp Dener { 702c14b763aSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 703c14b763aSAlp Dener PetscErrorCode ierr; 704c14b763aSAlp Dener TaoLineSearchConvergedReason ls_reason; 705c14b763aSAlp Dener 706c14b763aSAlp Dener PetscReal e_min, gdx, delta; 707c14b763aSAlp Dener PetscInt bfgsUpdates; 708c14b763aSAlp Dener 709c14b763aSAlp Dener PetscFunctionBegin; 710c14b763aSAlp Dener /* Perform the linesearch */ 711c14b763aSAlp Dener ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr); 712c14b763aSAlp Dener ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 713c14b763aSAlp Dener 7149b6ef848SAlp Dener while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && (stepType != BNK_GRADIENT || stepType !=BNK_SCALED_GRADIENT)) { 715c14b763aSAlp Dener /* Linesearch failed, revert solution */ 716c14b763aSAlp Dener bnk->f = bnk->fold; 717c14b763aSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 718c14b763aSAlp Dener ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 719c14b763aSAlp Dener 720c14b763aSAlp Dener switch(stepType) { 721c14b763aSAlp Dener case BNK_NEWTON: 7228d5ead36SAlp Dener /* Failed to obtain acceptable iterate with Newton step 723c14b763aSAlp Dener Update the perturbation for next time */ 724c14b763aSAlp Dener if (bnk->pert <= 0.0) { 725c14b763aSAlp Dener /* Initialize the perturbation */ 726c14b763aSAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 727c14b763aSAlp Dener if (bnk->is_gltr) { 728c14b763aSAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 729c14b763aSAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 730c14b763aSAlp Dener } 731c14b763aSAlp Dener } else { 732c14b763aSAlp Dener /* Increase the perturbation */ 733c14b763aSAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 734c14b763aSAlp Dener } 735c14b763aSAlp Dener 736c14b763aSAlp Dener if (BNK_PC_BFGS != bnk->pc_type) { 737c14b763aSAlp Dener /* We don't have the bfgs matrix around and being updated 738c14b763aSAlp Dener Must use gradient direction in this case */ 739c14b763aSAlp Dener ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 740c14b763aSAlp Dener stepType = BNK_GRADIENT; 741c14b763aSAlp Dener } else { 742c14b763aSAlp Dener /* Attempt to use the BFGS direction */ 743c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 7448d5ead36SAlp Dener /* Check for success (descent direction) 7458d5ead36SAlp Dener NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */ 746c14b763aSAlp Dener ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 747c14b763aSAlp Dener if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 748c14b763aSAlp Dener /* BFGS direction is not descent or direction produced not a number 749c14b763aSAlp Dener We can assert bfgsUpdates > 1 in this case 750c14b763aSAlp Dener Use steepest descent direction (scaled) */ 7519b6ef848SAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 752c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 753c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 754c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 755c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 756c14b763aSAlp Dener 757c14b763aSAlp Dener bfgsUpdates = 1; 758c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 759c14b763aSAlp Dener } else { 760c14b763aSAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 761c14b763aSAlp Dener if (1 == bfgsUpdates) { 762c14b763aSAlp Dener /* The first BFGS direction is always the scaled gradient */ 763c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 764c14b763aSAlp Dener } else { 765c14b763aSAlp Dener stepType = BNK_BFGS; 766c14b763aSAlp Dener } 767c14b763aSAlp Dener } 768c14b763aSAlp Dener } 769c14b763aSAlp Dener break; 770c14b763aSAlp Dener 771c14b763aSAlp Dener case BNK_BFGS: 772c14b763aSAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS 773c14b763aSAlp Dener Failed to obtain acceptable iterate with BFGS step 774c14b763aSAlp Dener Attempt to use the scaled gradient direction */ 7759b6ef848SAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 776c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 777c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 778c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 779c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 780c14b763aSAlp Dener 781c14b763aSAlp Dener bfgsUpdates = 1; 782c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 783c14b763aSAlp Dener break; 784c14b763aSAlp Dener 785c14b763aSAlp Dener case BNK_SCALED_GRADIENT: 786c14b763aSAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS 787c14b763aSAlp Dener The scaled gradient step did not produce a new iterate; 788c14b763aSAlp Dener attemp to use the gradient direction. 789c14b763aSAlp Dener Need to make sure we are not using a different diagonal scaling */ 790c14b763aSAlp Dener ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr); 791c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr); 792c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 793c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 794c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 795c14b763aSAlp Dener 796c14b763aSAlp Dener bfgsUpdates = 1; 797c14b763aSAlp Dener stepType = BNK_GRADIENT; 798c14b763aSAlp Dener break; 799c14b763aSAlp Dener } 8008d5ead36SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 8018d5ead36SAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 8022f75a4aaSAlp Dener ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr); 803c14b763aSAlp Dener 8048d5ead36SAlp Dener /* Perform one last line search with the fall-back step */ 805c14b763aSAlp Dener ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr); 806c14b763aSAlp Dener ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 807c14b763aSAlp Dener } 808c14b763aSAlp Dener *reason = ls_reason; 809c14b763aSAlp Dener PetscFunctionReturn(0); 810c14b763aSAlp Dener } 811c14b763aSAlp Dener 812df278d8fSAlp Dener /*------------------------------------------------------------*/ 813df278d8fSAlp Dener 814df278d8fSAlp Dener /* Routine for updating the trust radius. 815df278d8fSAlp Dener 816df278d8fSAlp Dener Function features three different update methods: 817df278d8fSAlp Dener 1) Line-search step length based 818df278d8fSAlp Dener 2) Predicted decrease on the CG quadratic model 819df278d8fSAlp Dener 3) Interpolation 820df278d8fSAlp Dener */ 821df278d8fSAlp Dener 82228017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept) 823080d2917SAlp Dener { 824080d2917SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 825080d2917SAlp Dener PetscErrorCode ierr; 826080d2917SAlp Dener 827b1c2d0e3SAlp Dener PetscReal step, kappa; 828080d2917SAlp Dener PetscReal gdx, tau_1, tau_2, tau_min, tau_max; 829080d2917SAlp Dener 830080d2917SAlp Dener PetscFunctionBegin; 831080d2917SAlp Dener /* Update trust region radius */ 832080d2917SAlp Dener *accept = PETSC_FALSE; 83328017e9fSAlp Dener switch(updateType) { 834080d2917SAlp Dener case BNK_UPDATE_STEP: 835c14b763aSAlp Dener *accept = PETSC_TRUE; /* always accept here because line search succeeded */ 836080d2917SAlp Dener if (stepType == BNK_NEWTON) { 837080d2917SAlp Dener ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr); 838080d2917SAlp Dener if (step < bnk->nu1) { 839080d2917SAlp Dener /* Very bad step taken; reduce radius */ 840080d2917SAlp Dener tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 841080d2917SAlp Dener } else if (step < bnk->nu2) { 842080d2917SAlp Dener /* Reasonably bad step taken; reduce radius */ 843080d2917SAlp Dener tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust); 844080d2917SAlp Dener } else if (step < bnk->nu3) { 845080d2917SAlp Dener /* Reasonable step was taken; leave radius alone */ 846080d2917SAlp Dener if (bnk->omega3 < 1.0) { 847080d2917SAlp Dener tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust); 848080d2917SAlp Dener } else if (bnk->omega3 > 1.0) { 849080d2917SAlp Dener tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust); 850080d2917SAlp Dener } 851080d2917SAlp Dener } else if (step < bnk->nu4) { 852080d2917SAlp Dener /* Full step taken; increase the radius */ 853080d2917SAlp Dener tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust); 854080d2917SAlp Dener } else { 855080d2917SAlp Dener /* More than full step taken; increase the radius */ 856080d2917SAlp Dener tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust); 857080d2917SAlp Dener } 858080d2917SAlp Dener } else { 859080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 860080d2917SAlp Dener tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 861080d2917SAlp Dener } 862080d2917SAlp Dener break; 863080d2917SAlp Dener 864080d2917SAlp Dener case BNK_UPDATE_REDUCTION: 865080d2917SAlp Dener if (stepType == BNK_NEWTON) { 866b1c2d0e3SAlp Dener if (prered < 0.0) { 867fed79b8eSAlp Dener /* The predicted reduction has the wrong sign. This cannot 868fed79b8eSAlp Dener happen in infinite precision arithmetic. Step should 869fed79b8eSAlp Dener be rejected! */ 870080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 871fed79b8eSAlp Dener } 872fed79b8eSAlp Dener else { 873b1c2d0e3SAlp Dener if (PetscIsInfOrNanReal(actred)) { 874080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 875080d2917SAlp Dener } else { 8760a4511e9SAlp Dener if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && 8770a4511e9SAlp Dener (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) { 878080d2917SAlp Dener kappa = 1.0; 879fed79b8eSAlp Dener } 880fed79b8eSAlp Dener else { 881080d2917SAlp Dener kappa = actred / prered; 882080d2917SAlp Dener } 883fed79b8eSAlp Dener 884fed79b8eSAlp Dener /* Accept or reject the step and update radius */ 885080d2917SAlp Dener if (kappa < bnk->eta1) { 886fed79b8eSAlp Dener /* Reject the step */ 887080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 888fed79b8eSAlp Dener } 889fed79b8eSAlp Dener else { 890fed79b8eSAlp Dener /* Accept the step */ 891c133c014SAlp Dener *accept = PETSC_TRUE; 892c133c014SAlp Dener /* Update the trust region radius only if the computed step is at the trust radius boundary */ 8938d5ead36SAlp Dener if (bnk->dnorm == tao->trust) { 894080d2917SAlp Dener if (kappa < bnk->eta2) { 895080d2917SAlp Dener /* Marginal bad step */ 896c133c014SAlp Dener tao->trust = bnk->alpha2 * tao->trust; 897080d2917SAlp Dener } 898fed79b8eSAlp Dener else if (kappa < bnk->eta3) { 899fed79b8eSAlp Dener /* Reasonable step */ 900fed79b8eSAlp Dener tao->trust = bnk->alpha3 * tao->trust; 901fed79b8eSAlp Dener } 902fed79b8eSAlp Dener else if (kappa < bnk->eta4) { 903080d2917SAlp Dener /* Good step */ 904c133c014SAlp Dener tao->trust = bnk->alpha4 * tao->trust; 905fed79b8eSAlp Dener } 906fed79b8eSAlp Dener else { 907080d2917SAlp Dener /* Very good step */ 908c133c014SAlp Dener tao->trust = bnk->alpha5 * tao->trust; 909080d2917SAlp Dener } 910c133c014SAlp Dener } 911080d2917SAlp Dener } 912080d2917SAlp Dener } 913080d2917SAlp Dener } 914080d2917SAlp Dener } else { 915080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 916080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust); 917080d2917SAlp Dener } 918080d2917SAlp Dener break; 919080d2917SAlp Dener 920080d2917SAlp Dener default: 921080d2917SAlp Dener if (stepType == BNK_NEWTON) { 922b1c2d0e3SAlp Dener if (prered < 0.0) { 923080d2917SAlp Dener /* The predicted reduction has the wrong sign. This cannot */ 924080d2917SAlp Dener /* happen in infinite precision arithmetic. Step should */ 925080d2917SAlp Dener /* be rejected! */ 926080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 927080d2917SAlp Dener } else { 928b1c2d0e3SAlp Dener if (PetscIsInfOrNanReal(actred)) { 929080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 930080d2917SAlp Dener } else { 931080d2917SAlp Dener if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 932080d2917SAlp Dener kappa = 1.0; 933080d2917SAlp Dener } else { 934080d2917SAlp Dener kappa = actred / prered; 935080d2917SAlp Dener } 936080d2917SAlp Dener 937080d2917SAlp Dener ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 938080d2917SAlp Dener tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred); 939080d2917SAlp Dener tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred); 940080d2917SAlp Dener tau_min = PetscMin(tau_1, tau_2); 941080d2917SAlp Dener tau_max = PetscMax(tau_1, tau_2); 942080d2917SAlp Dener 943080d2917SAlp Dener if (kappa >= 1.0 - bnk->mu1) { 944080d2917SAlp Dener /* Great agreement */ 945080d2917SAlp Dener *accept = PETSC_TRUE; 946080d2917SAlp Dener if (tau_max < 1.0) { 947080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 948080d2917SAlp Dener } else if (tau_max > bnk->gamma4) { 949080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm); 950080d2917SAlp Dener } else { 951080d2917SAlp Dener tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 952080d2917SAlp Dener } 953080d2917SAlp Dener } else if (kappa >= 1.0 - bnk->mu2) { 954080d2917SAlp Dener /* Good agreement */ 955080d2917SAlp Dener *accept = PETSC_TRUE; 956080d2917SAlp Dener if (tau_max < bnk->gamma2) { 957080d2917SAlp Dener tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 958080d2917SAlp Dener } else if (tau_max > bnk->gamma3) { 959080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 960080d2917SAlp Dener } else if (tau_max < 1.0) { 961080d2917SAlp Dener tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 962080d2917SAlp Dener } else { 963080d2917SAlp Dener tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 964080d2917SAlp Dener } 965080d2917SAlp Dener } else { 966080d2917SAlp Dener /* Not good agreement */ 967080d2917SAlp Dener if (tau_min > 1.0) { 968080d2917SAlp Dener tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 969080d2917SAlp Dener } else if (tau_max < bnk->gamma1) { 970080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 971080d2917SAlp Dener } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) { 972080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 973080d2917SAlp Dener } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) { 974080d2917SAlp Dener tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm); 975080d2917SAlp Dener } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) { 976080d2917SAlp Dener tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm); 977080d2917SAlp Dener } else { 978080d2917SAlp Dener tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 979080d2917SAlp Dener } 980080d2917SAlp Dener } 981080d2917SAlp Dener } 982080d2917SAlp Dener } 983080d2917SAlp Dener } else { 984080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 985080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust); 986080d2917SAlp Dener } 98728017e9fSAlp Dener break; 988080d2917SAlp Dener } 989c133c014SAlp Dener /* Make sure the radius does not violate min and max settings */ 990c133c014SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 991fed79b8eSAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 992080d2917SAlp Dener PetscFunctionReturn(0); 993080d2917SAlp Dener } 994080d2917SAlp Dener 995eb910715SAlp Dener /* ---------------------------------------------------------- */ 996df278d8fSAlp Dener 99762675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType) 99862675beeSAlp Dener { 99962675beeSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 100062675beeSAlp Dener 100162675beeSAlp Dener PetscFunctionBegin; 100262675beeSAlp Dener switch (stepType) { 100362675beeSAlp Dener case BNK_NEWTON: 100462675beeSAlp Dener ++bnk->newt; 100562675beeSAlp Dener break; 100662675beeSAlp Dener case BNK_BFGS: 100762675beeSAlp Dener ++bnk->bfgs; 100862675beeSAlp Dener break; 100962675beeSAlp Dener case BNK_SCALED_GRADIENT: 101062675beeSAlp Dener ++bnk->sgrad; 101162675beeSAlp Dener break; 101262675beeSAlp Dener case BNK_GRADIENT: 101362675beeSAlp Dener ++bnk->grad; 101462675beeSAlp Dener break; 101562675beeSAlp Dener default: 101662675beeSAlp Dener break; 101762675beeSAlp Dener } 101862675beeSAlp Dener PetscFunctionReturn(0); 101962675beeSAlp Dener } 102062675beeSAlp Dener 102162675beeSAlp Dener /* ---------------------------------------------------------- */ 102262675beeSAlp Dener 10239b6ef848SAlp Dener PetscErrorCode TaoSetUp_BNK(Tao tao) 1024eb910715SAlp Dener { 1025eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 1026eb910715SAlp Dener PetscErrorCode ierr; 10279b6ef848SAlp Dener KSPType ksp_type; 1028eb910715SAlp Dener 1029eb910715SAlp Dener PetscFunctionBegin; 1030eb910715SAlp Dener if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);} 1031eb910715SAlp Dener if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);} 1032eb910715SAlp Dener if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);} 1033eb910715SAlp Dener if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);} 1034eb910715SAlp Dener if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);} 103509164190SAlp Dener if (!bnk->Xwork) {ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);} 103609164190SAlp Dener if (!bnk->Gwork) {ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);} 103709164190SAlp Dener if (!bnk->unprojected_gradient) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);} 103809164190SAlp Dener if (!bnk->unprojected_gradient_old) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);} 1039198282dbSAlp Dener if (!bnk->G_inactive) {ierr = VecDuplicate(tao->solution,&bnk->G_inactive);CHKERRQ(ierr);} 1040*5e9b73cbSAlp Dener if (!bnk->Diag_min) {ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);} 1041*5e9b73cbSAlp Dener if (!bnk->Diag_max) {ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);} 1042eb910715SAlp Dener bnk->Diag = 0; 104362675beeSAlp Dener bnk->inactive_work = 0; 104462675beeSAlp Dener bnk->active_work = 0; 104562675beeSAlp Dener bnk->inactive_idx = 0; 104662675beeSAlp Dener bnk->active_idx = 0; 104762675beeSAlp Dener bnk->active_lower = 0; 104862675beeSAlp Dener bnk->active_upper = 0; 104962675beeSAlp Dener bnk->active_fixed = 0; 1050eb910715SAlp Dener bnk->M = 0; 105109164190SAlp Dener bnk->H_inactive = 0; 105209164190SAlp Dener bnk->Hpre_inactive = 0; 10539b6ef848SAlp Dener ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr); 10549b6ef848SAlp Dener ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr); 10559b6ef848SAlp Dener ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr); 10569b6ef848SAlp Dener ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr); 1057eb910715SAlp Dener PetscFunctionReturn(0); 1058eb910715SAlp Dener } 1059eb910715SAlp Dener 1060eb910715SAlp Dener /*------------------------------------------------------------*/ 1061df278d8fSAlp Dener 1062eb910715SAlp Dener static PetscErrorCode TaoDestroy_BNK(Tao tao) 1063eb910715SAlp Dener { 1064eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 1065eb910715SAlp Dener PetscErrorCode ierr; 1066eb910715SAlp Dener 1067eb910715SAlp Dener PetscFunctionBegin; 1068eb910715SAlp Dener if (tao->setupcalled) { 1069eb910715SAlp Dener ierr = VecDestroy(&bnk->W);CHKERRQ(ierr); 1070eb910715SAlp Dener ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr); 1071eb910715SAlp Dener ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr); 107209164190SAlp Dener ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr); 107309164190SAlp Dener ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr); 107409164190SAlp Dener ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr); 107509164190SAlp Dener ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr); 1076198282dbSAlp Dener ierr = VecDestroy(&bnk->G_inactive);CHKERRQ(ierr); 107762675beeSAlp Dener ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr); 107862675beeSAlp Dener ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr); 1079*5e9b73cbSAlp Dener } 1080*5e9b73cbSAlp Dener ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr); 1081eb910715SAlp Dener ierr = MatDestroy(&bnk->M);CHKERRQ(ierr); 108262675beeSAlp Dener if (bnk->Hpre_inactive != bnk->H_inactive) {ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);} 108362675beeSAlp Dener ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr); 1084eb910715SAlp Dener ierr = PetscFree(tao->data);CHKERRQ(ierr); 1085eb910715SAlp Dener PetscFunctionReturn(0); 1086eb910715SAlp Dener } 1087eb910715SAlp Dener 1088eb910715SAlp Dener /*------------------------------------------------------------*/ 1089df278d8fSAlp Dener 1090eb910715SAlp Dener static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao) 1091eb910715SAlp Dener { 1092eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 1093eb910715SAlp Dener PetscErrorCode ierr; 1094eb910715SAlp Dener 1095eb910715SAlp Dener PetscFunctionBegin; 1096eb910715SAlp Dener ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr); 10978d5ead36SAlp 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); 10988d5ead36SAlp 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); 10998d5ead36SAlp 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); 11008d5ead36SAlp 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); 11012f75a4aaSAlp 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); 11028d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr); 11038d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr); 11048d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr); 11058d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr); 11068d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr); 11078d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr); 11088d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr); 11098d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr); 11108d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr); 11118d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr); 11128d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr); 11138d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr); 11148d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr); 11158d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr); 11168d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr); 11178d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr); 11188d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr); 11198d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr); 11208d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr); 11218d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr); 11228d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr); 11238d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr); 11248d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr); 11258d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr); 11268d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr); 11278d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr); 11288d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr); 11298d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr); 11308d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr); 11318d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr); 11328d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr); 11338d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr); 11348d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr); 11358d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr); 11368d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr); 11378d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr); 11388d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr); 11398d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr); 11408d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr); 11418d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr); 11428d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr); 11438d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr); 11448d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr); 11458d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr); 11468d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr); 11470a4511e9SAlp Dener ierr = PetscOptionsReal("-tao_bnk_as_tol", "initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr); 11480a4511e9SAlp Dener ierr = PetscOptionsReal("-tao_bnk_as_step", "step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr); 1149eb910715SAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 1150eb910715SAlp Dener ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 1151eb910715SAlp Dener ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 1152eb910715SAlp Dener PetscFunctionReturn(0); 1153eb910715SAlp Dener } 1154eb910715SAlp Dener 1155eb910715SAlp Dener /*------------------------------------------------------------*/ 1156df278d8fSAlp Dener 1157eb910715SAlp Dener static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer) 1158eb910715SAlp Dener { 1159eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 1160eb910715SAlp Dener PetscInt nrejects; 1161eb910715SAlp Dener PetscBool isascii; 1162eb910715SAlp Dener PetscErrorCode ierr; 1163eb910715SAlp Dener 1164eb910715SAlp Dener PetscFunctionBegin; 1165eb910715SAlp Dener ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 1166eb910715SAlp Dener if (isascii) { 1167eb910715SAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1168eb910715SAlp Dener if (BNK_PC_BFGS == bnk->pc_type && bnk->M) { 1169eb910715SAlp Dener ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr); 1170eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr); 1171eb910715SAlp Dener } 1172eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr); 1173eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr); 1174eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr); 1175eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr); 1176eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr); 1177eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr); 1178eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr); 1179eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr); 1180eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr); 1181eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr); 1182eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr); 1183eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr); 1184eb910715SAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1185eb910715SAlp Dener } 1186eb910715SAlp Dener PetscFunctionReturn(0); 1187eb910715SAlp Dener } 1188eb910715SAlp Dener 1189eb910715SAlp Dener /* ---------------------------------------------------------- */ 1190df278d8fSAlp Dener 1191eb910715SAlp Dener /*MC 1192eb910715SAlp Dener TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms. 119366ed3702SAlp Dener At each iteration, the BNK methods solve the symmetric 1194eb910715SAlp Dener system of equations to obtain the step diretion dk: 1195eb910715SAlp Dener Hk dk = -gk 11962b97c8d8SAlp Dener for free variables only. The step can be globalized either through 11972b97c8d8SAlp Dener trust-region methods, or a line search, or a heuristic mixture of both. 1198eb910715SAlp Dener 1199eb910715SAlp Dener Options Database Keys: 12008d5ead36SAlp Dener + -tao_bnk_pc_type - "none","ahess","bfgs","petsc" 12018d5ead36SAlp Dener . -tao_bnk_bfgs_scale_type - "ahess","phess","bfgs" 12028d5ead36SAlp Dener . -tao_bnk_init_type - "constant","direction","interpolation" 12038d5ead36SAlp Dener . -tao_bnk_update_type - "step","direction","interpolation" 12042f75a4aaSAlp Dener . -tao_bnk_as_type - "none","bertsekas" 12058d5ead36SAlp Dener . -tao_bnk_sval - perturbation starting value 12068d5ead36SAlp Dener . -tao_bnk_imin - minimum initial perturbation 12078d5ead36SAlp Dener . -tao_bnk_imax - maximum initial perturbation 12088d5ead36SAlp Dener . -tao_bnk_pmin - minimum perturbation 12098d5ead36SAlp Dener . -tao_bnk_pmax - maximum perturbation 12108d5ead36SAlp Dener . -tao_bnk_pgfac - growth factor 12118d5ead36SAlp Dener . -tao_bnk_psfac - shrink factor 12128d5ead36SAlp Dener . -tao_bnk_imfac - initial merit factor 12138d5ead36SAlp Dener . -tao_bnk_pmgfac - merit growth factor 12148d5ead36SAlp Dener . -tao_bnk_pmsfac - merit shrink factor 12158d5ead36SAlp Dener . -tao_bnk_eta1 - poor steplength; reduce radius 12168d5ead36SAlp Dener . -tao_bnk_eta2 - reasonable steplength; leave radius 12178d5ead36SAlp Dener . -tao_bnk_eta3 - good steplength; increase readius 12188d5ead36SAlp Dener . -tao_bnk_eta4 - excellent steplength; greatly increase radius 12198d5ead36SAlp Dener . -tao_bnk_alpha1 - alpha1 reduction 12208d5ead36SAlp Dener . -tao_bnk_alpha2 - alpha2 reduction 12218d5ead36SAlp Dener . -tao_bnk_alpha3 - alpha3 reduction 12228d5ead36SAlp Dener . -tao_bnk_alpha4 - alpha4 reduction 12238d5ead36SAlp Dener . -tao_bnk_alpha - alpha5 reduction 12248d5ead36SAlp Dener . -tao_bnk_mu1 - mu1 interpolation update 12258d5ead36SAlp Dener . -tao_bnk_mu2 - mu2 interpolation update 12268d5ead36SAlp Dener . -tao_bnk_gamma1 - gamma1 interpolation update 12278d5ead36SAlp Dener . -tao_bnk_gamma2 - gamma2 interpolation update 12288d5ead36SAlp Dener . -tao_bnk_gamma3 - gamma3 interpolation update 12298d5ead36SAlp Dener . -tao_bnk_gamma4 - gamma4 interpolation update 12308d5ead36SAlp Dener . -tao_bnk_theta - theta interpolation update 12318d5ead36SAlp Dener . -tao_bnk_omega1 - omega1 step update 12328d5ead36SAlp Dener . -tao_bnk_omega2 - omega2 step update 12338d5ead36SAlp Dener . -tao_bnk_omega3 - omega3 step update 12348d5ead36SAlp Dener . -tao_bnk_omega4 - omega4 step update 12358d5ead36SAlp Dener . -tao_bnk_omega5 - omega5 step update 12368d5ead36SAlp Dener . -tao_bnk_mu1_i - mu1 interpolation init factor 12378d5ead36SAlp Dener . -tao_bnk_mu2_i - mu2 interpolation init factor 12388d5ead36SAlp Dener . -tao_bnk_gamma1_i - gamma1 interpolation init factor 12398d5ead36SAlp Dener . -tao_bnk_gamma2_i - gamma2 interpolation init factor 12408d5ead36SAlp Dener . -tao_bnk_gamma3_i - gamma3 interpolation init factor 12418d5ead36SAlp Dener . -tao_bnk_gamma4_i - gamma4 interpolation init factor 12422f75a4aaSAlp Dener . -tao_bnk_theta_i - theta interpolation init factor 12432f75a4aaSAlp Dener - -tao_bnk_bound_tol - initial tolerance used in estimating bounded active variables 1244eb910715SAlp Dener 1245eb910715SAlp Dener Level: beginner 1246eb910715SAlp Dener M*/ 1247eb910715SAlp Dener 1248eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao) 1249eb910715SAlp Dener { 1250eb910715SAlp Dener TAO_BNK *bnk; 1251eb910715SAlp Dener const char *morethuente_type = TAOLINESEARCHMT; 1252eb910715SAlp Dener PetscErrorCode ierr; 1253eb910715SAlp Dener 1254eb910715SAlp Dener PetscFunctionBegin; 1255eb910715SAlp Dener ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr); 1256eb910715SAlp Dener 1257eb910715SAlp Dener tao->ops->setup = TaoSetUp_BNK; 1258eb910715SAlp Dener tao->ops->view = TaoView_BNK; 1259eb910715SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BNK; 1260eb910715SAlp Dener tao->ops->destroy = TaoDestroy_BNK; 1261eb910715SAlp Dener 1262eb910715SAlp Dener /* Override default settings (unless already changed) */ 1263eb910715SAlp Dener if (!tao->max_it_changed) tao->max_it = 50; 1264eb910715SAlp Dener if (!tao->trust0_changed) tao->trust0 = 100.0; 1265eb910715SAlp Dener 1266eb910715SAlp Dener tao->data = (void*)bnk; 1267eb910715SAlp Dener 126866ed3702SAlp Dener /* Hessian shifting parameters */ 1269eb910715SAlp Dener bnk->sval = 0.0; 1270eb910715SAlp Dener bnk->imin = 1.0e-4; 1271eb910715SAlp Dener bnk->imax = 1.0e+2; 1272eb910715SAlp Dener bnk->imfac = 1.0e-1; 1273eb910715SAlp Dener 1274eb910715SAlp Dener bnk->pmin = 1.0e-12; 1275eb910715SAlp Dener bnk->pmax = 1.0e+2; 1276eb910715SAlp Dener bnk->pgfac = 1.0e+1; 1277eb910715SAlp Dener bnk->psfac = 4.0e-1; 1278eb910715SAlp Dener bnk->pmgfac = 1.0e-1; 1279eb910715SAlp Dener bnk->pmsfac = 1.0e-1; 1280eb910715SAlp Dener 1281eb910715SAlp Dener /* Default values for trust-region radius update based on steplength */ 1282eb910715SAlp Dener bnk->nu1 = 0.25; 1283eb910715SAlp Dener bnk->nu2 = 0.50; 1284eb910715SAlp Dener bnk->nu3 = 1.00; 1285eb910715SAlp Dener bnk->nu4 = 1.25; 1286eb910715SAlp Dener 1287eb910715SAlp Dener bnk->omega1 = 0.25; 1288eb910715SAlp Dener bnk->omega2 = 0.50; 1289eb910715SAlp Dener bnk->omega3 = 1.00; 1290eb910715SAlp Dener bnk->omega4 = 2.00; 1291eb910715SAlp Dener bnk->omega5 = 4.00; 1292eb910715SAlp Dener 1293eb910715SAlp Dener /* Default values for trust-region radius update based on reduction */ 1294eb910715SAlp Dener bnk->eta1 = 1.0e-4; 1295eb910715SAlp Dener bnk->eta2 = 0.25; 1296eb910715SAlp Dener bnk->eta3 = 0.50; 1297eb910715SAlp Dener bnk->eta4 = 0.90; 1298eb910715SAlp Dener 1299eb910715SAlp Dener bnk->alpha1 = 0.25; 1300eb910715SAlp Dener bnk->alpha2 = 0.50; 1301eb910715SAlp Dener bnk->alpha3 = 1.00; 1302eb910715SAlp Dener bnk->alpha4 = 2.00; 1303eb910715SAlp Dener bnk->alpha5 = 4.00; 1304eb910715SAlp Dener 1305eb910715SAlp Dener /* Default values for trust-region radius update based on interpolation */ 1306eb910715SAlp Dener bnk->mu1 = 0.10; 1307eb910715SAlp Dener bnk->mu2 = 0.50; 1308eb910715SAlp Dener 1309eb910715SAlp Dener bnk->gamma1 = 0.25; 1310eb910715SAlp Dener bnk->gamma2 = 0.50; 1311eb910715SAlp Dener bnk->gamma3 = 2.00; 1312eb910715SAlp Dener bnk->gamma4 = 4.00; 1313eb910715SAlp Dener 1314eb910715SAlp Dener bnk->theta = 0.05; 1315eb910715SAlp Dener 1316eb910715SAlp Dener /* Default values for trust region initialization based on interpolation */ 1317eb910715SAlp Dener bnk->mu1_i = 0.35; 1318eb910715SAlp Dener bnk->mu2_i = 0.50; 1319eb910715SAlp Dener 1320eb910715SAlp Dener bnk->gamma1_i = 0.0625; 1321eb910715SAlp Dener bnk->gamma2_i = 0.5; 1322eb910715SAlp Dener bnk->gamma3_i = 2.0; 1323eb910715SAlp Dener bnk->gamma4_i = 5.0; 1324eb910715SAlp Dener 1325eb910715SAlp Dener bnk->theta_i = 0.25; 1326eb910715SAlp Dener 1327eb910715SAlp Dener /* Remaining parameters */ 1328eb910715SAlp Dener bnk->min_radius = 1.0e-10; 1329eb910715SAlp Dener bnk->max_radius = 1.0e10; 1330770b7498SAlp Dener bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); 13310a4511e9SAlp Dener bnk->as_tol = 1.0e-3; 13320a4511e9SAlp Dener bnk->as_step = 1.0e-3; 133362675beeSAlp Dener bnk->dmin = 1.0e-6; 133462675beeSAlp Dener bnk->dmax = 1.0e6; 1335eb910715SAlp Dener 1336eb910715SAlp Dener bnk->pc_type = BNK_PC_BFGS; 1337eb910715SAlp Dener bnk->bfgs_scale_type = BFGS_SCALE_PHESS; 1338eb910715SAlp Dener bnk->init_type = BNK_INIT_INTERPOLATION; 1339fed79b8eSAlp Dener bnk->update_type = BNK_UPDATE_INTERPOLATION; 13402f75a4aaSAlp Dener bnk->as_type = BNK_AS_BERTSEKAS; 1341eb910715SAlp Dener 1342eb910715SAlp Dener ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 1343eb910715SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 1344eb910715SAlp Dener ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 1345eb910715SAlp Dener ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 1346eb910715SAlp Dener ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 1347eb910715SAlp Dener 1348eb910715SAlp Dener /* Set linear solver to default for symmetric matrices */ 1349eb910715SAlp Dener ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 1350eb910715SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 1351eb910715SAlp Dener ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr); 1352eb910715SAlp Dener ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr); 1353eb910715SAlp Dener PetscFunctionReturn(0); 1354eb910715SAlp Dener } 1355