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 26eb910715SAlp Dener PetscErrorCode TaoBNKInitialize(Tao tao) 27eb910715SAlp Dener { 28eb910715SAlp Dener PetscErrorCode ierr; 29eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 30eb910715SAlp Dener KSPType ksp_type; 31eb910715SAlp Dener PC pc; 32eb910715SAlp Dener 33*770b7498SAlp Dener PetscReal fmin, ftrial, prered, actred, kappa, sigma; 34eb910715SAlp Dener PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 35eb910715SAlp Dener PetscReal delta, step = 1.0; 36eb910715SAlp Dener 37eb910715SAlp Dener PetscInt n,N,needH = 1; 38eb910715SAlp Dener 39eb910715SAlp Dener PetscInt i_max = 5; 40eb910715SAlp Dener PetscInt j_max = 1; 41eb910715SAlp Dener PetscInt i, j; 42eb910715SAlp Dener 43eb910715SAlp Dener PetscFunctionBegin; 44eb910715SAlp Dener /* Number of times ksp stopped because of these reasons */ 45eb910715SAlp Dener bnk->ksp_atol = 0; 46eb910715SAlp Dener bnk->ksp_rtol = 0; 47eb910715SAlp Dener bnk->ksp_dtol = 0; 48eb910715SAlp Dener bnk->ksp_ctol = 0; 49eb910715SAlp Dener bnk->ksp_negc = 0; 50eb910715SAlp Dener bnk->ksp_iter = 0; 51eb910715SAlp Dener bnk->ksp_othr = 0; 52eb910715SAlp Dener 53fed79b8eSAlp Dener /* Initialize the Hessian perturbation */ 54fed79b8eSAlp Dener bnk->pert = bnk->sval; 55fed79b8eSAlp Dener 56eb910715SAlp Dener /* Initialize trust-region radius when using nash, stcg, or gltr 57eb910715SAlp Dener Command automatically ignored for other methods 58eb910715SAlp Dener Will be reset during the first iteration 59eb910715SAlp Dener */ 60eb910715SAlp Dener ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr); 61eb910715SAlp Dener ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr); 62eb910715SAlp Dener ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr); 63eb910715SAlp Dener ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr); 64eb910715SAlp Dener 65eb910715SAlp Dener ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr); 66eb910715SAlp Dener 67eb910715SAlp Dener if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 68eb910715SAlp Dener if (tao->trust0 < 0.0) SETERRQ(PETSC_COMM_SELF,1,"Initial radius negative"); 69eb910715SAlp Dener tao->trust = tao->trust0; 70eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 71eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 72eb910715SAlp Dener } 73eb910715SAlp Dener 74eb910715SAlp Dener /* Get vectors we will need */ 75eb910715SAlp Dener if (BNK_PC_BFGS == bnk->pc_type && !bnk->M) { 76eb910715SAlp Dener ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 77eb910715SAlp Dener ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 78eb910715SAlp Dener ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr); 79eb910715SAlp Dener ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr); 80eb910715SAlp Dener } 81eb910715SAlp Dener 82eb910715SAlp Dener /* create vectors for the limited memory preconditioner */ 83eb910715SAlp Dener if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_BFGS != bnk->bfgs_scale_type)) { 84eb910715SAlp Dener if (!bnk->Diag) { 85eb910715SAlp Dener ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr); 86eb910715SAlp Dener } 87eb910715SAlp Dener } 88eb910715SAlp Dener 89eb910715SAlp Dener /* Modify the preconditioner to use the bfgs approximation */ 90eb910715SAlp Dener ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 91eb910715SAlp Dener switch(bnk->pc_type) { 92eb910715SAlp Dener case BNK_PC_NONE: 93eb910715SAlp Dener ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr); 94eb910715SAlp Dener ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 95eb910715SAlp Dener break; 96eb910715SAlp Dener 97eb910715SAlp Dener case BNK_PC_AHESS: 98eb910715SAlp Dener ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr); 99eb910715SAlp Dener ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 100eb910715SAlp Dener ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 101eb910715SAlp Dener break; 102eb910715SAlp Dener 103eb910715SAlp Dener case BNK_PC_BFGS: 104eb910715SAlp Dener ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr); 105eb910715SAlp Dener ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 106eb910715SAlp Dener ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr); 107eb910715SAlp Dener ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr); 108eb910715SAlp Dener ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr); 109eb910715SAlp Dener break; 110eb910715SAlp Dener 111eb910715SAlp Dener default: 112eb910715SAlp Dener /* Use the pc method set by pc_type */ 113eb910715SAlp Dener break; 114eb910715SAlp Dener } 115eb910715SAlp Dener 116eb910715SAlp Dener /* Initialize trust-region radius. The initialization is only performed 117eb910715SAlp Dener when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */ 118eb910715SAlp Dener if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 119eb910715SAlp Dener switch(bnk->init_type) { 120eb910715SAlp Dener case BNK_INIT_CONSTANT: 121eb910715SAlp Dener /* Use the initial radius specified */ 122eb910715SAlp Dener break; 123eb910715SAlp Dener 124eb910715SAlp Dener case BNK_INIT_INTERPOLATION: 125eb910715SAlp Dener /* Use the initial radius specified */ 126eb910715SAlp Dener max_radius = 0.0; 127eb910715SAlp Dener 128eb910715SAlp Dener for (j = 0; j < j_max; ++j) { 129eb910715SAlp Dener fmin = bnk->f; 130eb910715SAlp Dener sigma = 0.0; 131eb910715SAlp Dener 132eb910715SAlp Dener if (needH) { 13362602cfbSAlp Dener /* Compute the Hessian at the new step, and extract the inactive subsystem */ 134eb910715SAlp Dener ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 13562602cfbSAlp Dener ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr); 13662602cfbSAlp Dener ierr = VecWhichInactive(tao->XL,tao->solution,bnk->unprojected_gradient,tao->XU,PETSC_TRUE,&bnk->inactive_idx);CHKERRQ(ierr); 13762602cfbSAlp Dener ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr); 138eb910715SAlp Dener needH = 0; 139eb910715SAlp Dener } 140eb910715SAlp Dener 141eb910715SAlp Dener for (i = 0; i < i_max; ++i) { 14262602cfbSAlp Dener /* Take a steepest descent step and snap it to bounds */ 14362602cfbSAlp Dener ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 14462602cfbSAlp Dener ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr); 14562602cfbSAlp Dener ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 146*770b7498SAlp Dener /* Recompute the step after bound snapping so that it can be used in predicted decrease calculation later */ 147eb910715SAlp Dener ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr); 14862602cfbSAlp Dener ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr); 14962602cfbSAlp Dener /* Compute the objective at the trial */ 15062602cfbSAlp Dener ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr); 15162602cfbSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 152eb910715SAlp Dener if (PetscIsInfOrNanReal(ftrial)) { 153eb910715SAlp Dener tau = bnk->gamma1_i; 154eb910715SAlp Dener } else { 155eb910715SAlp Dener if (ftrial < fmin) { 156eb910715SAlp Dener fmin = ftrial; 157eb910715SAlp Dener sigma = -tao->trust / bnk->gnorm; 158eb910715SAlp Dener } 159*770b7498SAlp Dener /* Compute the predicted and actual reduction */ 16062602cfbSAlp Dener ierr = MatMult(bnk->H_inactive, tao->gradient, bnk->W);CHKERRQ(ierr); 16162602cfbSAlp Dener ierr = VecDot(tao->gradient, bnk->W, &prered);CHKERRQ(ierr); 162eb910715SAlp Dener prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm)); 163eb910715SAlp Dener actred = bnk->f - ftrial; 164eb910715SAlp Dener if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 165eb910715SAlp Dener kappa = 1.0; 166eb910715SAlp Dener } else { 167eb910715SAlp Dener kappa = actred / prered; 168eb910715SAlp Dener } 169eb910715SAlp Dener 170eb910715SAlp Dener tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred); 171eb910715SAlp Dener tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred); 172eb910715SAlp Dener tau_min = PetscMin(tau_1, tau_2); 173eb910715SAlp Dener tau_max = PetscMax(tau_1, tau_2); 174eb910715SAlp Dener 175eb910715SAlp Dener if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) { 176eb910715SAlp Dener /* Great agreement */ 177eb910715SAlp Dener max_radius = PetscMax(max_radius, tao->trust); 178eb910715SAlp Dener 179eb910715SAlp Dener if (tau_max < 1.0) { 180eb910715SAlp Dener tau = bnk->gamma3_i; 181eb910715SAlp Dener } else if (tau_max > bnk->gamma4_i) { 182eb910715SAlp Dener tau = bnk->gamma4_i; 183eb910715SAlp Dener } else if (tau_1 >= 1.0 && tau_1 <= bnk->gamma4_i && tau_2 < 1.0) { 184eb910715SAlp Dener tau = tau_1; 185eb910715SAlp Dener } else if (tau_2 >= 1.0 && tau_2 <= bnk->gamma4_i && tau_1 < 1.0) { 186eb910715SAlp Dener tau = tau_2; 187eb910715SAlp Dener } else { 188eb910715SAlp Dener tau = tau_max; 189eb910715SAlp Dener } 190eb910715SAlp Dener } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) { 191eb910715SAlp Dener /* Good agreement */ 192eb910715SAlp Dener max_radius = PetscMax(max_radius, tao->trust); 193eb910715SAlp Dener 194eb910715SAlp Dener if (tau_max < bnk->gamma2_i) { 195eb910715SAlp Dener tau = bnk->gamma2_i; 196eb910715SAlp Dener } else if (tau_max > bnk->gamma3_i) { 197eb910715SAlp Dener tau = bnk->gamma3_i; 198eb910715SAlp Dener } else { 199eb910715SAlp Dener tau = tau_max; 200eb910715SAlp Dener } 201eb910715SAlp Dener } else { 202eb910715SAlp Dener /* Not good agreement */ 203eb910715SAlp Dener if (tau_min > 1.0) { 204eb910715SAlp Dener tau = bnk->gamma2_i; 205eb910715SAlp Dener } else if (tau_max < bnk->gamma1_i) { 206eb910715SAlp Dener tau = bnk->gamma1_i; 207eb910715SAlp Dener } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) { 208eb910715SAlp Dener tau = bnk->gamma1_i; 209eb910715SAlp Dener } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 210eb910715SAlp Dener tau = tau_1; 211eb910715SAlp Dener } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 212eb910715SAlp Dener tau = tau_2; 213eb910715SAlp Dener } else { 214eb910715SAlp Dener tau = tau_max; 215eb910715SAlp Dener } 216eb910715SAlp Dener } 217eb910715SAlp Dener } 218eb910715SAlp Dener tao->trust = tau * tao->trust; 219eb910715SAlp Dener } 220eb910715SAlp Dener 221eb910715SAlp Dener if (fmin < bnk->f) { 222eb910715SAlp Dener bnk->f = fmin; 223eb910715SAlp Dener ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr); 22487602d1eSAlp Dener ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 22509164190SAlp Dener ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr); 22609164190SAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 227eb910715SAlp Dener 228eb910715SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 229eb910715SAlp Dener if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN"); 230eb910715SAlp Dener needH = 1; 231eb910715SAlp Dener 232eb910715SAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 233eb910715SAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr); 234eb910715SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 235eb910715SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 236eb910715SAlp Dener } 237eb910715SAlp Dener } 238eb910715SAlp Dener tao->trust = PetscMax(tao->trust, max_radius); 239eb910715SAlp Dener 240eb910715SAlp Dener /* Modify the radius if it is too large or small */ 241eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 242eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 243eb910715SAlp Dener break; 244eb910715SAlp Dener 245eb910715SAlp Dener default: 246eb910715SAlp Dener /* Norm of the first direction will initialize radius */ 247eb910715SAlp Dener tao->trust = 0.0; 248eb910715SAlp Dener break; 249eb910715SAlp Dener } 250eb910715SAlp Dener } 251eb910715SAlp Dener 252eb910715SAlp Dener /* Set initial scaling for the BFGS preconditioner 253eb910715SAlp Dener This step is done after computing the initial trust-region radius 254eb910715SAlp Dener since the function value may have decreased */ 255eb910715SAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { 256eb910715SAlp Dener if (bnk->f != 0.0) { 257eb910715SAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 258eb910715SAlp Dener } else { 259eb910715SAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 260eb910715SAlp Dener } 261eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr); 262eb910715SAlp Dener } 263eb910715SAlp Dener 264eb910715SAlp Dener /* Set counter for gradient/reset steps*/ 265eb910715SAlp Dener bnk->newt = 0; 266eb910715SAlp Dener bnk->bfgs = 0; 267eb910715SAlp Dener bnk->sgrad = 0; 268eb910715SAlp Dener bnk->grad = 0; 269eb910715SAlp Dener PetscFunctionReturn(0); 270eb910715SAlp Dener } 271eb910715SAlp Dener 272df278d8fSAlp Dener /*------------------------------------------------------------*/ 273df278d8fSAlp Dener 274df278d8fSAlp Dener /* Routine for computing the Newton step. 275df278d8fSAlp Dener 276df278d8fSAlp Dener If the safeguard is enabled, the Newton step is verified to be a 277df278d8fSAlp Dener descent direction, with fallbacks onto BFGS, scaled gradient, and unscaled 278df278d8fSAlp Dener gradient steps if/when necessary. 279df278d8fSAlp Dener 280df278d8fSAlp Dener The function reports back on which type of step has ultimately been stored 281df278d8fSAlp Dener under tao->stepdirection. 282df278d8fSAlp Dener */ 283df278d8fSAlp Dener 284fed79b8eSAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool safeguard, PetscInt *stepType) 285eb910715SAlp Dener { 286eb910715SAlp Dener PetscErrorCode ierr; 287eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 288eb910715SAlp Dener KSPConvergedReason ksp_reason; 289eb910715SAlp Dener 290080d2917SAlp Dener PetscReal gdx, delta, e_min; 291eb910715SAlp Dener PetscInt bfgsUpdates = 0; 292eb910715SAlp Dener PetscInt kspits; 293eb910715SAlp Dener 294eb910715SAlp Dener PetscFunctionBegin; 295eb910715SAlp Dener /* Shift the Hessian matrix */ 296eb910715SAlp Dener if (bnk->pert > 0) { 297eb910715SAlp Dener ierr = MatShift(tao->hessian, bnk->pert);CHKERRQ(ierr); 298eb910715SAlp Dener if (tao->hessian != tao->hessian_pre) { 299eb910715SAlp Dener ierr = MatShift(tao->hessian_pre, bnk->pert);CHKERRQ(ierr); 300eb910715SAlp Dener } 301eb910715SAlp Dener } 302eb910715SAlp Dener 3038d5ead36SAlp Dener /* Determine the inactive and active sets */ 304*770b7498SAlp Dener ierr = ISDestroy(&bnk->inactive_old);CHKERRQ(ierr); 305*770b7498SAlp Dener ierr = ISDuplicate(bnk->inactive_idx, &bnk->inactive_old);CHKERRQ(ierr); 306*770b7498SAlp Dener ierr = ISCopy(bnk->inactive_idx, bnk->inactive_old);CHKERRQ(ierr); 30709164190SAlp Dener ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr); 30809164190SAlp Dener ierr = VecWhichInactive(tao->XL,tao->solution,bnk->unprojected_gradient,tao->XU,PETSC_TRUE,&bnk->inactive_idx);CHKERRQ(ierr); 30909164190SAlp Dener 310eb3ba6a7SAlp Dener /* Prepare masked matrices for the inactive set */ 31109164190SAlp Dener ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr); 31209164190SAlp Dener ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr); 313eb3ba6a7SAlp Dener ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr); 314eb3ba6a7SAlp Dener if (tao->hessian == tao->hessian_pre) { 315eb3ba6a7SAlp Dener bnk->Hpre_inactive = bnk->H_inactive; 316eb3ba6a7SAlp Dener } else { 317eb3ba6a7SAlp Dener ierr = TaoMatGetSubMat(tao->hessian_pre, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->Hpre_inactive);CHKERRQ(ierr); 318eb3ba6a7SAlp Dener } 31909164190SAlp Dener 320eb910715SAlp Dener /* Solve the Newton system of equations */ 32109164190SAlp Dener ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr); 322eb910715SAlp Dener if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 323fed79b8eSAlp Dener ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 324eb3ba6a7SAlp Dener ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 325eb910715SAlp Dener ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 326eb910715SAlp Dener tao->ksp_its+=kspits; 327eb910715SAlp Dener tao->ksp_tot_its+=kspits; 328080d2917SAlp Dener ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr); 329eb910715SAlp Dener 330eb910715SAlp Dener if (0.0 == tao->trust) { 331eb910715SAlp Dener /* Radius was uninitialized; use the norm of the direction */ 332080d2917SAlp Dener if (bnk->dnorm > 0.0) { 333080d2917SAlp Dener tao->trust = bnk->dnorm; 334eb910715SAlp Dener 335eb910715SAlp Dener /* Modify the radius if it is too large or small */ 336eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 337eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 338eb910715SAlp Dener } else { 339eb910715SAlp Dener /* The direction was bad; set radius to default value and re-solve 340eb910715SAlp Dener the trust-region subproblem to get a direction */ 341eb910715SAlp Dener tao->trust = tao->trust0; 342eb910715SAlp Dener 343eb910715SAlp Dener /* Modify the radius if it is too large or small */ 344eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 345eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 346eb910715SAlp Dener 347fed79b8eSAlp Dener ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 348eb3ba6a7SAlp Dener ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 349eb910715SAlp Dener ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 350eb910715SAlp Dener tao->ksp_its+=kspits; 351eb910715SAlp Dener tao->ksp_tot_its+=kspits; 352080d2917SAlp Dener ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr); 353eb910715SAlp Dener 354080d2917SAlp Dener if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 355eb910715SAlp Dener } 356eb910715SAlp Dener } 357eb910715SAlp Dener } else { 358eb3ba6a7SAlp Dener ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 359eb910715SAlp Dener ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr); 360eb910715SAlp Dener tao->ksp_its += kspits; 361eb910715SAlp Dener tao->ksp_tot_its+=kspits; 362eb910715SAlp Dener } 363*770b7498SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 364*770b7498SAlp Dener ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 365fed79b8eSAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 366*770b7498SAlp Dener 367*770b7498SAlp Dener /* Record convergence reasons */ 368fed79b8eSAlp Dener ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 369*770b7498SAlp Dener if (KSP_CONVERGED_ATOL == ksp_reason) { 370*770b7498SAlp Dener ++bnk->ksp_atol; 371*770b7498SAlp Dener } else if (KSP_CONVERGED_RTOL == ksp_reason) { 372*770b7498SAlp Dener ++bnk->ksp_rtol; 373*770b7498SAlp Dener } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) { 374*770b7498SAlp Dener ++bnk->ksp_ctol; 375*770b7498SAlp Dener } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) { 376*770b7498SAlp Dener ++bnk->ksp_negc; 377*770b7498SAlp Dener } else if (KSP_DIVERGED_DTOL == ksp_reason) { 378*770b7498SAlp Dener ++bnk->ksp_dtol; 379*770b7498SAlp Dener } else if (KSP_DIVERGED_ITS == ksp_reason) { 380*770b7498SAlp Dener ++bnk->ksp_iter; 381*770b7498SAlp Dener } else { 382*770b7498SAlp Dener ++bnk->ksp_othr; 383*770b7498SAlp Dener } 384fed79b8eSAlp Dener 385fed79b8eSAlp Dener /* Make sure the BFGS preconditioner is healthy */ 386fed79b8eSAlp Dener if (bnk->pc_type == BNK_PC_BFGS) { 387fed79b8eSAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 388fed79b8eSAlp Dener if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (bfgsUpdates > 1)) { 389fed79b8eSAlp Dener /* Preconditioner is numerically indefinite; reset the approximation. */ 390eb910715SAlp Dener if (bnk->f != 0.0) { 391eb910715SAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 392eb910715SAlp Dener } else { 393eb910715SAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 394eb910715SAlp Dener } 395eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr); 396eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 39709164190SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 398eb910715SAlp Dener bfgsUpdates = 1; 399eb910715SAlp Dener } 400fed79b8eSAlp Dener } 401eb910715SAlp Dener 402eb910715SAlp Dener /* Check for success (descent direction) */ 403fed79b8eSAlp Dener if (safeguard) { 404080d2917SAlp Dener ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 405eb910715SAlp Dener if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 406eb910715SAlp Dener /* Newton step is not descent or direction produced Inf or NaN 407eb910715SAlp Dener Update the perturbation for next time */ 408eb910715SAlp Dener if (bnk->pert <= 0.0) { 409eb910715SAlp Dener /* Initialize the perturbation */ 410eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 411eb910715SAlp Dener if (bnk->is_gltr) { 412eb910715SAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 413eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 414eb910715SAlp Dener } 415eb910715SAlp Dener } else { 416eb910715SAlp Dener /* Increase the perturbation */ 417eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 418eb910715SAlp Dener } 419eb910715SAlp Dener 420eb910715SAlp Dener if (BNK_PC_BFGS != bnk->pc_type) { 421eb910715SAlp Dener /* We don't have the bfgs matrix around and updated 422eb910715SAlp Dener Must use gradient direction in this case */ 423080d2917SAlp Dener ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 424eb910715SAlp Dener ++bnk->grad; 425eb910715SAlp Dener *stepType = BNK_GRADIENT; 426eb910715SAlp Dener } else { 427eb910715SAlp Dener /* Attempt to use the BFGS direction */ 42809164190SAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 429eb910715SAlp Dener 4308d5ead36SAlp Dener /* Check for success (descent direction) 4318d5ead36SAlp Dener NOTE: Negative gdx here means not a descent direction because 4328d5ead36SAlp Dener the fall-back step is missing a negative sign. */ 433080d2917SAlp Dener ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 4348d5ead36SAlp Dener if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 435eb910715SAlp Dener /* BFGS direction is not descent or direction produced not a number 436eb910715SAlp Dener We can assert bfgsUpdates > 1 in this case because 437eb910715SAlp Dener the first solve produces the scaled gradient direction, 438eb910715SAlp Dener which is guaranteed to be descent */ 439eb910715SAlp Dener 440eb910715SAlp Dener /* Use steepest descent direction (scaled) */ 441eb910715SAlp Dener if (bnk->f != 0.0) { 442eb910715SAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 443eb910715SAlp Dener } else { 444eb910715SAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 445eb910715SAlp Dener } 446eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 447eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 44809164190SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 44909164190SAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 450eb910715SAlp Dener 451eb910715SAlp Dener bfgsUpdates = 1; 452eb910715SAlp Dener ++bnk->sgrad; 453eb910715SAlp Dener *stepType = BNK_SCALED_GRADIENT; 454eb910715SAlp Dener } else { 455*770b7498SAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 456eb910715SAlp Dener if (1 == bfgsUpdates) { 457eb910715SAlp Dener /* The first BFGS direction is always the scaled gradient */ 458eb910715SAlp Dener ++bnk->sgrad; 459eb910715SAlp Dener *stepType = BNK_SCALED_GRADIENT; 460eb910715SAlp Dener } else { 461eb910715SAlp Dener ++bnk->bfgs; 462eb910715SAlp Dener *stepType = BNK_BFGS; 463eb910715SAlp Dener } 464eb910715SAlp Dener } 465eb910715SAlp Dener } 4668d5ead36SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 467*770b7498SAlp Dener ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 4688d5ead36SAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 469eb910715SAlp Dener } else { 470eb910715SAlp Dener /* Computed Newton step is descent */ 471eb910715SAlp Dener switch (ksp_reason) { 472eb910715SAlp Dener case KSP_DIVERGED_NANORINF: 473eb910715SAlp Dener case KSP_DIVERGED_BREAKDOWN: 474eb910715SAlp Dener case KSP_DIVERGED_INDEFINITE_MAT: 475eb910715SAlp Dener case KSP_DIVERGED_INDEFINITE_PC: 476eb910715SAlp Dener case KSP_CONVERGED_CG_NEG_CURVE: 477eb910715SAlp Dener /* Matrix or preconditioner is indefinite; increase perturbation */ 478eb910715SAlp Dener if (bnk->pert <= 0.0) { 479eb910715SAlp Dener /* Initialize the perturbation */ 480eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 481eb910715SAlp Dener if (bnk->is_gltr) { 482eb910715SAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr); 483eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 484eb910715SAlp Dener } 485eb910715SAlp Dener } else { 486eb910715SAlp Dener /* Increase the perturbation */ 487eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 488eb910715SAlp Dener } 489eb910715SAlp Dener break; 490eb910715SAlp Dener 491eb910715SAlp Dener default: 492eb910715SAlp Dener /* Newton step computation is good; decrease perturbation */ 493eb910715SAlp Dener bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm); 494eb910715SAlp Dener if (bnk->pert < bnk->pmin) { 495eb910715SAlp Dener bnk->pert = 0.0; 496eb910715SAlp Dener } 497eb910715SAlp Dener break; 498eb910715SAlp Dener } 499eb910715SAlp Dener 500eb910715SAlp Dener ++bnk->newt; 501fed79b8eSAlp Dener *stepType = BNK_NEWTON; 502fed79b8eSAlp Dener } 503fed79b8eSAlp Dener } else { 504fed79b8eSAlp Dener ++bnk->newt; 505fed79b8eSAlp Dener bnk->pert = 0.0; 506fed79b8eSAlp Dener *stepType = BNK_NEWTON; 507eb910715SAlp Dener } 508eb910715SAlp Dener PetscFunctionReturn(0); 509eb910715SAlp Dener } 510eb910715SAlp Dener 511df278d8fSAlp Dener /*------------------------------------------------------------*/ 512df278d8fSAlp Dener 513df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search. 514df278d8fSAlp Dener 515df278d8fSAlp Dener Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the 516df278d8fSAlp Dener Newton step does not produce a valid step length. 517df278d8fSAlp Dener */ 518df278d8fSAlp Dener 519c14b763aSAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason) 520c14b763aSAlp Dener { 521c14b763aSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 522c14b763aSAlp Dener PetscErrorCode ierr; 523c14b763aSAlp Dener TaoLineSearchConvergedReason ls_reason; 524c14b763aSAlp Dener 525c14b763aSAlp Dener PetscReal e_min, gdx, delta; 526c14b763aSAlp Dener PetscInt bfgsUpdates; 527c14b763aSAlp Dener 528c14b763aSAlp Dener PetscFunctionBegin; 529c14b763aSAlp Dener /* Perform the linesearch */ 530*770b7498SAlp Dener *steplen = 1.0; 531c14b763aSAlp Dener ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr); 532c14b763aSAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 533c14b763aSAlp Dener ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 534c14b763aSAlp Dener 535c14b763aSAlp Dener while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != BNK_GRADIENT) { 536c14b763aSAlp Dener /* Linesearch failed, revert solution */ 537c14b763aSAlp Dener bnk->f = bnk->fold; 538c14b763aSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 539c14b763aSAlp Dener ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 540c14b763aSAlp Dener ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 541c14b763aSAlp Dener 542c14b763aSAlp Dener switch(stepType) { 543c14b763aSAlp Dener case BNK_NEWTON: 5448d5ead36SAlp Dener /* Failed to obtain acceptable iterate with Newton step 545c14b763aSAlp Dener Update the perturbation for next time */ 546c14b763aSAlp Dener --bnk->newt; 547c14b763aSAlp Dener if (bnk->pert <= 0.0) { 548c14b763aSAlp Dener /* Initialize the perturbation */ 549c14b763aSAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 550c14b763aSAlp Dener if (bnk->is_gltr) { 551c14b763aSAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 552c14b763aSAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 553c14b763aSAlp Dener } 554c14b763aSAlp Dener } else { 555c14b763aSAlp Dener /* Increase the perturbation */ 556c14b763aSAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 557c14b763aSAlp Dener } 558c14b763aSAlp Dener 559c14b763aSAlp Dener if (BNK_PC_BFGS != bnk->pc_type) { 560c14b763aSAlp Dener /* We don't have the bfgs matrix around and being updated 561c14b763aSAlp Dener Must use gradient direction in this case */ 562c14b763aSAlp Dener ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 563c14b763aSAlp Dener ++bnk->grad; 564c14b763aSAlp Dener stepType = BNK_GRADIENT; 565c14b763aSAlp Dener } else { 566c14b763aSAlp Dener /* Attempt to use the BFGS direction */ 567c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 5688d5ead36SAlp Dener /* Check for success (descent direction) 5698d5ead36SAlp Dener NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */ 570c14b763aSAlp Dener ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 571c14b763aSAlp Dener if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 572c14b763aSAlp Dener /* BFGS direction is not descent or direction produced not a number 573c14b763aSAlp Dener We can assert bfgsUpdates > 1 in this case 574c14b763aSAlp Dener Use steepest descent direction (scaled) */ 575c14b763aSAlp Dener if (bnk->f != 0.0) { 576c14b763aSAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 577c14b763aSAlp Dener } else { 578c14b763aSAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 579c14b763aSAlp Dener } 580c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 581c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 582c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 583c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 584c14b763aSAlp Dener 585c14b763aSAlp Dener bfgsUpdates = 1; 586c14b763aSAlp Dener ++bnk->sgrad; 587c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 588c14b763aSAlp Dener } else { 589c14b763aSAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 590c14b763aSAlp Dener if (1 == bfgsUpdates) { 591c14b763aSAlp Dener /* The first BFGS direction is always the scaled gradient */ 592c14b763aSAlp Dener ++bnk->sgrad; 593c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 594c14b763aSAlp Dener } else { 595c14b763aSAlp Dener ++bnk->bfgs; 596c14b763aSAlp Dener stepType = BNK_BFGS; 597c14b763aSAlp Dener } 598c14b763aSAlp Dener } 599c14b763aSAlp Dener } 600c14b763aSAlp Dener break; 601c14b763aSAlp Dener 602c14b763aSAlp Dener case BNK_BFGS: 603c14b763aSAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS 604c14b763aSAlp Dener Failed to obtain acceptable iterate with BFGS step 605c14b763aSAlp Dener Attempt to use the scaled gradient direction */ 606c14b763aSAlp Dener if (bnk->f != 0.0) { 607c14b763aSAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 608c14b763aSAlp Dener } else { 609c14b763aSAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 610c14b763aSAlp Dener } 611c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 612c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 613c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 614c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 615c14b763aSAlp Dener 616c14b763aSAlp Dener bfgsUpdates = 1; 617c14b763aSAlp Dener ++bnk->sgrad; 618c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 619c14b763aSAlp Dener break; 620c14b763aSAlp Dener 621c14b763aSAlp Dener case BNK_SCALED_GRADIENT: 622c14b763aSAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS 623c14b763aSAlp Dener The scaled gradient step did not produce a new iterate; 624c14b763aSAlp Dener attemp to use the gradient direction. 625c14b763aSAlp Dener Need to make sure we are not using a different diagonal scaling */ 626c14b763aSAlp Dener ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr); 627c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr); 628c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 629c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 630c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 631c14b763aSAlp Dener 632c14b763aSAlp Dener bfgsUpdates = 1; 633c14b763aSAlp Dener ++bnk->grad; 634c14b763aSAlp Dener stepType = BNK_GRADIENT; 635c14b763aSAlp Dener break; 636c14b763aSAlp Dener } 6378d5ead36SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 638*770b7498SAlp Dener ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 6398d5ead36SAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 640c14b763aSAlp Dener 6418d5ead36SAlp Dener /* Perform one last line search with the fall-back step */ 642*770b7498SAlp Dener *steplen = 1.0; 643c14b763aSAlp Dener ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr); 644c14b763aSAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 645c14b763aSAlp Dener ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 646c14b763aSAlp Dener } 647c14b763aSAlp Dener *reason = ls_reason; 648c14b763aSAlp Dener PetscFunctionReturn(0); 649c14b763aSAlp Dener } 650c14b763aSAlp Dener 651df278d8fSAlp Dener /*------------------------------------------------------------*/ 652df278d8fSAlp Dener 653df278d8fSAlp Dener /* Routine for updating the trust radius. 654df278d8fSAlp Dener 655df278d8fSAlp Dener Function features three different update methods: 656df278d8fSAlp Dener 1) Line-search step length based 657df278d8fSAlp Dener 2) Predicted decrease on the CG quadratic model 658df278d8fSAlp Dener 3) Interpolation 659df278d8fSAlp Dener */ 660df278d8fSAlp Dener 661b1c2d0e3SAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt stepType, PetscBool *accept) 662080d2917SAlp Dener { 663080d2917SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 664080d2917SAlp Dener PetscErrorCode ierr; 665080d2917SAlp Dener 666b1c2d0e3SAlp Dener PetscReal step, kappa; 667080d2917SAlp Dener PetscReal gdx, tau_1, tau_2, tau_min, tau_max; 668080d2917SAlp Dener 669080d2917SAlp Dener PetscFunctionBegin; 670080d2917SAlp Dener /* Update trust region radius */ 671080d2917SAlp Dener *accept = PETSC_FALSE; 672080d2917SAlp Dener switch(bnk->update_type) { 673080d2917SAlp Dener case BNK_UPDATE_STEP: 674c14b763aSAlp Dener *accept = PETSC_TRUE; /* always accept here because line search succeeded */ 675080d2917SAlp Dener if (stepType == BNK_NEWTON) { 676080d2917SAlp Dener ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr); 677080d2917SAlp Dener if (step < bnk->nu1) { 678080d2917SAlp Dener /* Very bad step taken; reduce radius */ 679080d2917SAlp Dener tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 680080d2917SAlp Dener } else if (step < bnk->nu2) { 681080d2917SAlp Dener /* Reasonably bad step taken; reduce radius */ 682080d2917SAlp Dener tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust); 683080d2917SAlp Dener } else if (step < bnk->nu3) { 684080d2917SAlp Dener /* Reasonable step was taken; leave radius alone */ 685080d2917SAlp Dener if (bnk->omega3 < 1.0) { 686080d2917SAlp Dener tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust); 687080d2917SAlp Dener } else if (bnk->omega3 > 1.0) { 688080d2917SAlp Dener tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust); 689080d2917SAlp Dener } 690080d2917SAlp Dener } else if (step < bnk->nu4) { 691080d2917SAlp Dener /* Full step taken; increase the radius */ 692080d2917SAlp Dener tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust); 693080d2917SAlp Dener } else { 694080d2917SAlp Dener /* More than full step taken; increase the radius */ 695080d2917SAlp Dener tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust); 696080d2917SAlp Dener } 697080d2917SAlp Dener } else { 698080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 699080d2917SAlp Dener tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 700080d2917SAlp Dener } 701080d2917SAlp Dener break; 702080d2917SAlp Dener 703080d2917SAlp Dener case BNK_UPDATE_REDUCTION: 704080d2917SAlp Dener if (stepType == BNK_NEWTON) { 705b1c2d0e3SAlp Dener if (prered < 0.0) { 706fed79b8eSAlp Dener /* The predicted reduction has the wrong sign. This cannot 707fed79b8eSAlp Dener happen in infinite precision arithmetic. Step should 708fed79b8eSAlp Dener be rejected! */ 709080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 710fed79b8eSAlp Dener } 711fed79b8eSAlp Dener else { 712b1c2d0e3SAlp Dener if (PetscIsInfOrNanReal(actred)) { 713080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 714080d2917SAlp Dener } else { 715080d2917SAlp Dener if ((PetscAbsScalar(actred) <= bnk->epsilon) && 716080d2917SAlp Dener (PetscAbsScalar(prered) <= bnk->epsilon)) { 717080d2917SAlp Dener kappa = 1.0; 718fed79b8eSAlp Dener } 719fed79b8eSAlp Dener else { 720080d2917SAlp Dener kappa = actred / prered; 721080d2917SAlp Dener } 722fed79b8eSAlp Dener 723fed79b8eSAlp Dener /* Accept or reject the step and update radius */ 724080d2917SAlp Dener if (kappa < bnk->eta1) { 725fed79b8eSAlp Dener /* Reject the step */ 726080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 727fed79b8eSAlp Dener } 728fed79b8eSAlp Dener else { 729fed79b8eSAlp Dener /* Accept the step */ 730c133c014SAlp Dener *accept = PETSC_TRUE; 731c133c014SAlp Dener /* Update the trust region radius only if the computed step is at the trust radius boundary */ 7328d5ead36SAlp Dener if (bnk->dnorm == tao->trust) { 733080d2917SAlp Dener if (kappa < bnk->eta2) { 734080d2917SAlp Dener /* Marginal bad step */ 735c133c014SAlp Dener tao->trust = bnk->alpha2 * tao->trust; 736080d2917SAlp Dener } 737fed79b8eSAlp Dener else if (kappa < bnk->eta3) { 738fed79b8eSAlp Dener /* Reasonable step */ 739fed79b8eSAlp Dener tao->trust = bnk->alpha3 * tao->trust; 740fed79b8eSAlp Dener } 741fed79b8eSAlp Dener else if (kappa < bnk->eta4) { 742080d2917SAlp Dener /* Good step */ 743c133c014SAlp Dener tao->trust = bnk->alpha4 * tao->trust; 744fed79b8eSAlp Dener } 745fed79b8eSAlp Dener else { 746080d2917SAlp Dener /* Very good step */ 747c133c014SAlp Dener tao->trust = bnk->alpha5 * tao->trust; 748080d2917SAlp Dener } 749c133c014SAlp Dener } 750080d2917SAlp Dener } 751080d2917SAlp Dener } 752080d2917SAlp Dener } 753080d2917SAlp Dener } else { 754080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 755080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust); 756080d2917SAlp Dener } 757080d2917SAlp Dener break; 758080d2917SAlp Dener 759080d2917SAlp Dener default: 760080d2917SAlp Dener if (stepType == BNK_NEWTON) { 761b1c2d0e3SAlp Dener if (prered < 0.0) { 762080d2917SAlp Dener /* The predicted reduction has the wrong sign. This cannot */ 763080d2917SAlp Dener /* happen in infinite precision arithmetic. Step should */ 764080d2917SAlp Dener /* be rejected! */ 765080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 766080d2917SAlp Dener } else { 767b1c2d0e3SAlp Dener if (PetscIsInfOrNanReal(actred)) { 768080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 769080d2917SAlp Dener } else { 770080d2917SAlp Dener if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 771080d2917SAlp Dener kappa = 1.0; 772080d2917SAlp Dener } else { 773080d2917SAlp Dener kappa = actred / prered; 774080d2917SAlp Dener } 775080d2917SAlp Dener 776080d2917SAlp Dener ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 777080d2917SAlp Dener tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred); 778080d2917SAlp Dener tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred); 779080d2917SAlp Dener tau_min = PetscMin(tau_1, tau_2); 780080d2917SAlp Dener tau_max = PetscMax(tau_1, tau_2); 781080d2917SAlp Dener 782080d2917SAlp Dener if (kappa >= 1.0 - bnk->mu1) { 783080d2917SAlp Dener /* Great agreement */ 784080d2917SAlp Dener *accept = PETSC_TRUE; 785080d2917SAlp Dener if (tau_max < 1.0) { 786080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 787080d2917SAlp Dener } else if (tau_max > bnk->gamma4) { 788080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm); 789080d2917SAlp Dener } else { 790080d2917SAlp Dener tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 791080d2917SAlp Dener } 792080d2917SAlp Dener } else if (kappa >= 1.0 - bnk->mu2) { 793080d2917SAlp Dener /* Good agreement */ 794080d2917SAlp Dener *accept = PETSC_TRUE; 795080d2917SAlp Dener if (tau_max < bnk->gamma2) { 796080d2917SAlp Dener tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 797080d2917SAlp Dener } else if (tau_max > bnk->gamma3) { 798080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 799080d2917SAlp Dener } else if (tau_max < 1.0) { 800080d2917SAlp Dener tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 801080d2917SAlp Dener } else { 802080d2917SAlp Dener tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 803080d2917SAlp Dener } 804080d2917SAlp Dener } else { 805080d2917SAlp Dener /* Not good agreement */ 806080d2917SAlp Dener if (tau_min > 1.0) { 807080d2917SAlp Dener tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 808080d2917SAlp Dener } else if (tau_max < bnk->gamma1) { 809080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 810080d2917SAlp Dener } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) { 811080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 812080d2917SAlp Dener } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) { 813080d2917SAlp Dener tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm); 814080d2917SAlp Dener } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) { 815080d2917SAlp Dener tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm); 816080d2917SAlp Dener } else { 817080d2917SAlp Dener tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 818080d2917SAlp Dener } 819080d2917SAlp Dener } 820080d2917SAlp Dener } 821080d2917SAlp Dener } 822080d2917SAlp Dener } else { 823080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 824080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust); 825080d2917SAlp Dener } 826080d2917SAlp Dener } 827c133c014SAlp Dener /* Make sure the radius does not violate min and max settings */ 828c133c014SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 829fed79b8eSAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 830080d2917SAlp Dener PetscFunctionReturn(0); 831080d2917SAlp Dener } 832080d2917SAlp Dener 833eb910715SAlp Dener /* ---------------------------------------------------------- */ 834df278d8fSAlp Dener 835eb910715SAlp Dener static PetscErrorCode TaoSetUp_BNK(Tao tao) 836eb910715SAlp Dener { 837eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 838eb910715SAlp Dener PetscErrorCode ierr; 839eb910715SAlp Dener 840eb910715SAlp Dener PetscFunctionBegin; 841eb910715SAlp Dener if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);} 842eb910715SAlp Dener if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);} 843eb910715SAlp Dener if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);} 844eb910715SAlp Dener if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);} 845eb910715SAlp Dener if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);} 84609164190SAlp Dener if (!bnk->Xwork) {ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);} 84709164190SAlp Dener if (!bnk->Gwork) {ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);} 84809164190SAlp Dener if (!bnk->unprojected_gradient) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);} 84909164190SAlp Dener if (!bnk->unprojected_gradient_old) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);} 850eb910715SAlp Dener bnk->Diag = 0; 851eb910715SAlp Dener bnk->M = 0; 85209164190SAlp Dener bnk->H_inactive = 0; 85309164190SAlp Dener bnk->Hpre_inactive = 0; 854eb910715SAlp Dener PetscFunctionReturn(0); 855eb910715SAlp Dener } 856eb910715SAlp Dener 857eb910715SAlp Dener /*------------------------------------------------------------*/ 858df278d8fSAlp Dener 859eb910715SAlp Dener static PetscErrorCode TaoDestroy_BNK(Tao tao) 860eb910715SAlp Dener { 861eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 862eb910715SAlp Dener PetscErrorCode ierr; 863eb910715SAlp Dener 864eb910715SAlp Dener PetscFunctionBegin; 865eb910715SAlp Dener if (tao->setupcalled) { 866eb910715SAlp Dener ierr = VecDestroy(&bnk->W);CHKERRQ(ierr); 867eb910715SAlp Dener ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr); 868eb910715SAlp Dener ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr); 86909164190SAlp Dener ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr); 87009164190SAlp Dener ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr); 87109164190SAlp Dener ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr); 87209164190SAlp Dener ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr); 873eb910715SAlp Dener } 874eb910715SAlp Dener ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr); 875eb910715SAlp Dener ierr = MatDestroy(&bnk->M);CHKERRQ(ierr); 876eb910715SAlp Dener ierr = PetscFree(tao->data);CHKERRQ(ierr); 877eb910715SAlp Dener PetscFunctionReturn(0); 878eb910715SAlp Dener } 879eb910715SAlp Dener 880eb910715SAlp Dener /*------------------------------------------------------------*/ 881df278d8fSAlp Dener 882eb910715SAlp Dener static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao) 883eb910715SAlp Dener { 884eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 885eb910715SAlp Dener PetscErrorCode ierr; 886eb910715SAlp Dener 887eb910715SAlp Dener PetscFunctionBegin; 888eb910715SAlp Dener ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr); 8898d5ead36SAlp 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); 8908d5ead36SAlp 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); 8918d5ead36SAlp 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); 8928d5ead36SAlp 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); 8938d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr); 8948d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr); 8958d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr); 8968d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr); 8978d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr); 8988d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr); 8998d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr); 9008d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr); 9018d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr); 9028d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr); 9038d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr); 9048d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr); 9058d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr); 9068d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr); 9078d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr); 9088d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr); 9098d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr); 9108d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr); 9118d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr); 9128d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr); 9138d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr); 9148d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr); 9158d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr); 9168d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr); 9178d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr); 9188d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr); 9198d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr); 9208d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr); 9218d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr); 9228d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr); 9238d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr); 9248d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr); 9258d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr); 9268d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr); 9278d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr); 9288d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr); 9298d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr); 9308d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr); 9318d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr); 9328d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr); 9338d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr); 9348d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr); 9358d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr); 9368d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr); 9378d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr); 938eb910715SAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 939eb910715SAlp Dener ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 940eb910715SAlp Dener ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 941eb910715SAlp Dener PetscFunctionReturn(0); 942eb910715SAlp Dener } 943eb910715SAlp Dener 944eb910715SAlp Dener /*------------------------------------------------------------*/ 945df278d8fSAlp Dener 946eb910715SAlp Dener static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer) 947eb910715SAlp Dener { 948eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 949eb910715SAlp Dener PetscInt nrejects; 950eb910715SAlp Dener PetscBool isascii; 951eb910715SAlp Dener PetscErrorCode ierr; 952eb910715SAlp Dener 953eb910715SAlp Dener PetscFunctionBegin; 954eb910715SAlp Dener ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 955eb910715SAlp Dener if (isascii) { 956eb910715SAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 957eb910715SAlp Dener if (BNK_PC_BFGS == bnk->pc_type && bnk->M) { 958eb910715SAlp Dener ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr); 959eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr); 960eb910715SAlp Dener } 961eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr); 962eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr); 963eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr); 964eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr); 965eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr); 966eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr); 967eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr); 968eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr); 969eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr); 970eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr); 971eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr); 972eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr); 973eb910715SAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 974eb910715SAlp Dener } 975eb910715SAlp Dener PetscFunctionReturn(0); 976eb910715SAlp Dener } 977eb910715SAlp Dener 978eb910715SAlp Dener /* ---------------------------------------------------------- */ 979df278d8fSAlp Dener 980eb910715SAlp Dener /*MC 981eb910715SAlp Dener TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms. 98266ed3702SAlp Dener At each iteration, the BNK methods solve the symmetric 983eb910715SAlp Dener system of equations to obtain the step diretion dk: 984eb910715SAlp Dener Hk dk = -gk 9852b97c8d8SAlp Dener for free variables only. The step can be globalized either through 9862b97c8d8SAlp Dener trust-region methods, or a line search, or a heuristic mixture of both. 987eb910715SAlp Dener 988eb910715SAlp Dener Options Database Keys: 9898d5ead36SAlp Dener + -tao_bnk_pc_type - "none","ahess","bfgs","petsc" 9908d5ead36SAlp Dener . -tao_bnk_bfgs_scale_type - "ahess","phess","bfgs" 9918d5ead36SAlp Dener . -tao_bnk_init_type - "constant","direction","interpolation" 9928d5ead36SAlp Dener . -tao_bnk_update_type - "step","direction","interpolation" 9938d5ead36SAlp Dener . -tao_bnk_sval - perturbation starting value 9948d5ead36SAlp Dener . -tao_bnk_imin - minimum initial perturbation 9958d5ead36SAlp Dener . -tao_bnk_imax - maximum initial perturbation 9968d5ead36SAlp Dener . -tao_bnk_pmin - minimum perturbation 9978d5ead36SAlp Dener . -tao_bnk_pmax - maximum perturbation 9988d5ead36SAlp Dener . -tao_bnk_pgfac - growth factor 9998d5ead36SAlp Dener . -tao_bnk_psfac - shrink factor 10008d5ead36SAlp Dener . -tao_bnk_imfac - initial merit factor 10018d5ead36SAlp Dener . -tao_bnk_pmgfac - merit growth factor 10028d5ead36SAlp Dener . -tao_bnk_pmsfac - merit shrink factor 10038d5ead36SAlp Dener . -tao_bnk_eta1 - poor steplength; reduce radius 10048d5ead36SAlp Dener . -tao_bnk_eta2 - reasonable steplength; leave radius 10058d5ead36SAlp Dener . -tao_bnk_eta3 - good steplength; increase readius 10068d5ead36SAlp Dener . -tao_bnk_eta4 - excellent steplength; greatly increase radius 10078d5ead36SAlp Dener . -tao_bnk_alpha1 - alpha1 reduction 10088d5ead36SAlp Dener . -tao_bnk_alpha2 - alpha2 reduction 10098d5ead36SAlp Dener . -tao_bnk_alpha3 - alpha3 reduction 10108d5ead36SAlp Dener . -tao_bnk_alpha4 - alpha4 reduction 10118d5ead36SAlp Dener . -tao_bnk_alpha - alpha5 reduction 10128d5ead36SAlp Dener . -tao_bnk_mu1 - mu1 interpolation update 10138d5ead36SAlp Dener . -tao_bnk_mu2 - mu2 interpolation update 10148d5ead36SAlp Dener . -tao_bnk_gamma1 - gamma1 interpolation update 10158d5ead36SAlp Dener . -tao_bnk_gamma2 - gamma2 interpolation update 10168d5ead36SAlp Dener . -tao_bnk_gamma3 - gamma3 interpolation update 10178d5ead36SAlp Dener . -tao_bnk_gamma4 - gamma4 interpolation update 10188d5ead36SAlp Dener . -tao_bnk_theta - theta interpolation update 10198d5ead36SAlp Dener . -tao_bnk_omega1 - omega1 step update 10208d5ead36SAlp Dener . -tao_bnk_omega2 - omega2 step update 10218d5ead36SAlp Dener . -tao_bnk_omega3 - omega3 step update 10228d5ead36SAlp Dener . -tao_bnk_omega4 - omega4 step update 10238d5ead36SAlp Dener . -tao_bnk_omega5 - omega5 step update 10248d5ead36SAlp Dener . -tao_bnk_mu1_i - mu1 interpolation init factor 10258d5ead36SAlp Dener . -tao_bnk_mu2_i - mu2 interpolation init factor 10268d5ead36SAlp Dener . -tao_bnk_gamma1_i - gamma1 interpolation init factor 10278d5ead36SAlp Dener . -tao_bnk_gamma2_i - gamma2 interpolation init factor 10288d5ead36SAlp Dener . -tao_bnk_gamma3_i - gamma3 interpolation init factor 10298d5ead36SAlp Dener . -tao_bnk_gamma4_i - gamma4 interpolation init factor 10308d5ead36SAlp Dener - -tao_bnk_theta_i - theta interpolation init factor 1031eb910715SAlp Dener 1032eb910715SAlp Dener Level: beginner 1033eb910715SAlp Dener M*/ 1034eb910715SAlp Dener 1035eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao) 1036eb910715SAlp Dener { 1037eb910715SAlp Dener TAO_BNK *bnk; 1038eb910715SAlp Dener const char *morethuente_type = TAOLINESEARCHMT; 1039eb910715SAlp Dener PetscErrorCode ierr; 1040eb910715SAlp Dener 1041eb910715SAlp Dener PetscFunctionBegin; 1042eb910715SAlp Dener ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr); 1043eb910715SAlp Dener 1044eb910715SAlp Dener tao->ops->setup = TaoSetUp_BNK; 1045eb910715SAlp Dener tao->ops->view = TaoView_BNK; 1046eb910715SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BNK; 1047eb910715SAlp Dener tao->ops->destroy = TaoDestroy_BNK; 1048eb910715SAlp Dener 1049eb910715SAlp Dener /* Override default settings (unless already changed) */ 1050eb910715SAlp Dener if (!tao->max_it_changed) tao->max_it = 50; 1051eb910715SAlp Dener if (!tao->trust0_changed) tao->trust0 = 100.0; 1052eb910715SAlp Dener 1053eb910715SAlp Dener tao->data = (void*)bnk; 1054eb910715SAlp Dener 105566ed3702SAlp Dener /* Hessian shifting parameters */ 1056eb910715SAlp Dener bnk->sval = 0.0; 1057eb910715SAlp Dener bnk->imin = 1.0e-4; 1058eb910715SAlp Dener bnk->imax = 1.0e+2; 1059eb910715SAlp Dener bnk->imfac = 1.0e-1; 1060eb910715SAlp Dener 1061eb910715SAlp Dener bnk->pmin = 1.0e-12; 1062eb910715SAlp Dener bnk->pmax = 1.0e+2; 1063eb910715SAlp Dener bnk->pgfac = 1.0e+1; 1064eb910715SAlp Dener bnk->psfac = 4.0e-1; 1065eb910715SAlp Dener bnk->pmgfac = 1.0e-1; 1066eb910715SAlp Dener bnk->pmsfac = 1.0e-1; 1067eb910715SAlp Dener 1068eb910715SAlp Dener /* Default values for trust-region radius update based on steplength */ 1069eb910715SAlp Dener bnk->nu1 = 0.25; 1070eb910715SAlp Dener bnk->nu2 = 0.50; 1071eb910715SAlp Dener bnk->nu3 = 1.00; 1072eb910715SAlp Dener bnk->nu4 = 1.25; 1073eb910715SAlp Dener 1074eb910715SAlp Dener bnk->omega1 = 0.25; 1075eb910715SAlp Dener bnk->omega2 = 0.50; 1076eb910715SAlp Dener bnk->omega3 = 1.00; 1077eb910715SAlp Dener bnk->omega4 = 2.00; 1078eb910715SAlp Dener bnk->omega5 = 4.00; 1079eb910715SAlp Dener 1080eb910715SAlp Dener /* Default values for trust-region radius update based on reduction */ 1081eb910715SAlp Dener bnk->eta1 = 1.0e-4; 1082eb910715SAlp Dener bnk->eta2 = 0.25; 1083eb910715SAlp Dener bnk->eta3 = 0.50; 1084eb910715SAlp Dener bnk->eta4 = 0.90; 1085eb910715SAlp Dener 1086eb910715SAlp Dener bnk->alpha1 = 0.25; 1087eb910715SAlp Dener bnk->alpha2 = 0.50; 1088eb910715SAlp Dener bnk->alpha3 = 1.00; 1089eb910715SAlp Dener bnk->alpha4 = 2.00; 1090eb910715SAlp Dener bnk->alpha5 = 4.00; 1091eb910715SAlp Dener 1092eb910715SAlp Dener /* Default values for trust-region radius update based on interpolation */ 1093eb910715SAlp Dener bnk->mu1 = 0.10; 1094eb910715SAlp Dener bnk->mu2 = 0.50; 1095eb910715SAlp Dener 1096eb910715SAlp Dener bnk->gamma1 = 0.25; 1097eb910715SAlp Dener bnk->gamma2 = 0.50; 1098eb910715SAlp Dener bnk->gamma3 = 2.00; 1099eb910715SAlp Dener bnk->gamma4 = 4.00; 1100eb910715SAlp Dener 1101eb910715SAlp Dener bnk->theta = 0.05; 1102eb910715SAlp Dener 1103eb910715SAlp Dener /* Default values for trust region initialization based on interpolation */ 1104eb910715SAlp Dener bnk->mu1_i = 0.35; 1105eb910715SAlp Dener bnk->mu2_i = 0.50; 1106eb910715SAlp Dener 1107eb910715SAlp Dener bnk->gamma1_i = 0.0625; 1108eb910715SAlp Dener bnk->gamma2_i = 0.5; 1109eb910715SAlp Dener bnk->gamma3_i = 2.0; 1110eb910715SAlp Dener bnk->gamma4_i = 5.0; 1111eb910715SAlp Dener 1112eb910715SAlp Dener bnk->theta_i = 0.25; 1113eb910715SAlp Dener 1114eb910715SAlp Dener /* Remaining parameters */ 1115eb910715SAlp Dener bnk->min_radius = 1.0e-10; 1116eb910715SAlp Dener bnk->max_radius = 1.0e10; 1117*770b7498SAlp Dener bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); 1118eb910715SAlp Dener 1119eb910715SAlp Dener bnk->pc_type = BNK_PC_BFGS; 1120eb910715SAlp Dener bnk->bfgs_scale_type = BFGS_SCALE_PHESS; 1121eb910715SAlp Dener bnk->init_type = BNK_INIT_INTERPOLATION; 1122fed79b8eSAlp Dener bnk->update_type = BNK_UPDATE_INTERPOLATION; 1123eb910715SAlp Dener 1124eb910715SAlp Dener ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 1125eb910715SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 1126eb910715SAlp Dener ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 1127eb910715SAlp Dener ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 1128eb910715SAlp Dener ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 1129eb910715SAlp Dener 1130eb910715SAlp Dener /* Set linear solver to default for symmetric matrices */ 1131eb910715SAlp Dener ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 1132eb910715SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 1133eb910715SAlp Dener ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr); 1134eb910715SAlp Dener ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr); 1135eb910715SAlp Dener PetscFunctionReturn(0); 1136eb910715SAlp Dener } 1137