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 33770b7498SAlp 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); 146770b7498SAlp 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 } 159770b7498SAlp 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 274*2f75a4aaSAlp Dener /* Routine for estimating the active set */ 275*2f75a4aaSAlp Dener 276*2f75a4aaSAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao) 277*2f75a4aaSAlp Dener { 278*2f75a4aaSAlp Dener PetscErrorCode ierr; 279*2f75a4aaSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 280*2f75a4aaSAlp Dener 281*2f75a4aaSAlp Dener PetscFunctionBegin; 282*2f75a4aaSAlp Dener switch (bnk->as_type) { 283*2f75a4aaSAlp Dener case BNK_AS_NONE: 284*2f75a4aaSAlp Dener ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr); 285*2f75a4aaSAlp Dener ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr); 286*2f75a4aaSAlp Dener ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr); 287*2f75a4aaSAlp Dener ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr); 288*2f75a4aaSAlp Dener break; 289*2f75a4aaSAlp Dener 290*2f75a4aaSAlp Dener case BNK_AS_BERTSEKAS: 291*2f75a4aaSAlp Dener /* Compute the trial step vector with which we will estimate the active set at the next iteration */ 292*2f75a4aaSAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { 293*2f75a4aaSAlp Dener /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */ 294*2f75a4aaSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr); 295*2f75a4aaSAlp Dener } else { 296*2f75a4aaSAlp Dener /* BFGS preconditioner doesn't exist so let's invert the diagonal of the Hessian instead onto the gradient*/ 297*2f75a4aaSAlp Dener ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr); 298*2f75a4aaSAlp Dener ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr); 299*2f75a4aaSAlp Dener ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr); 300*2f75a4aaSAlp Dener } 301*2f75a4aaSAlp Dener ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr); 302*2f75a4aaSAlp Dener ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, &bnk->bound_tol, &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr); 303*2f75a4aaSAlp Dener 304*2f75a4aaSAlp Dener default: 305*2f75a4aaSAlp Dener break; 306*2f75a4aaSAlp Dener } 307*2f75a4aaSAlp Dener PetscFunctionReturn(0); 308*2f75a4aaSAlp Dener } 309*2f75a4aaSAlp Dener 310*2f75a4aaSAlp Dener /*------------------------------------------------------------*/ 311*2f75a4aaSAlp Dener 312*2f75a4aaSAlp Dener /* Routine for bounding the step direction */ 313*2f75a4aaSAlp Dener 314*2f75a4aaSAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, Vec step) 315*2f75a4aaSAlp Dener { 316*2f75a4aaSAlp Dener PetscErrorCode ierr; 317*2f75a4aaSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 318*2f75a4aaSAlp Dener 319*2f75a4aaSAlp Dener PetscFunctionBegin; 320*2f75a4aaSAlp Dener switch (bnk->as_type) { 321*2f75a4aaSAlp Dener case BNK_AS_NONE: 322*2f75a4aaSAlp Dener if (bnk->active_idx) { 323*2f75a4aaSAlp Dener ierr = VecGetSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr); 324*2f75a4aaSAlp Dener ierr = VecSet(bnk->active_work, 0.0);CHKERRQ(ierr); 325*2f75a4aaSAlp Dener ierr = VecRestoreSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr); 326*2f75a4aaSAlp Dener } 327*2f75a4aaSAlp Dener break; 328*2f75a4aaSAlp Dener 329*2f75a4aaSAlp Dener case BNK_AS_BERTSEKAS: 330*2f75a4aaSAlp Dener ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, step);CHKERRQ(ierr); 331*2f75a4aaSAlp Dener break; 332*2f75a4aaSAlp Dener 333*2f75a4aaSAlp Dener default: 334*2f75a4aaSAlp Dener break; 335*2f75a4aaSAlp Dener } 336*2f75a4aaSAlp Dener PetscFunctionReturn(0); 337*2f75a4aaSAlp Dener } 338*2f75a4aaSAlp Dener 339*2f75a4aaSAlp Dener /*------------------------------------------------------------*/ 340*2f75a4aaSAlp Dener 341df278d8fSAlp Dener /* Routine for computing the Newton step. 342df278d8fSAlp Dener 343df278d8fSAlp Dener If the safeguard is enabled, the Newton step is verified to be a 344df278d8fSAlp Dener descent direction, with fallbacks onto BFGS, scaled gradient, and unscaled 345df278d8fSAlp Dener gradient steps if/when necessary. 346df278d8fSAlp Dener 347df278d8fSAlp Dener The function reports back on which type of step has ultimately been stored 348df278d8fSAlp Dener under tao->stepdirection. 349df278d8fSAlp Dener */ 350df278d8fSAlp Dener 351fed79b8eSAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool safeguard, PetscInt *stepType) 352eb910715SAlp Dener { 353eb910715SAlp Dener PetscErrorCode ierr; 354eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 355eb910715SAlp Dener KSPConvergedReason ksp_reason; 356eb910715SAlp Dener 357080d2917SAlp Dener PetscReal gdx, delta, e_min; 358eb910715SAlp Dener PetscInt bfgsUpdates = 0; 359eb910715SAlp Dener PetscInt kspits; 360eb910715SAlp Dener 361eb910715SAlp Dener PetscFunctionBegin; 362eb910715SAlp Dener /* Shift the Hessian matrix */ 363eb910715SAlp Dener if (bnk->pert > 0) { 364eb910715SAlp Dener ierr = MatShift(tao->hessian, bnk->pert);CHKERRQ(ierr); 365eb910715SAlp Dener if (tao->hessian != tao->hessian_pre) { 366eb910715SAlp Dener ierr = MatShift(tao->hessian_pre, bnk->pert);CHKERRQ(ierr); 367eb910715SAlp Dener } 368eb910715SAlp Dener } 369eb910715SAlp Dener 370*2f75a4aaSAlp Dener /* Determine the active and inactive sets */ 371*2f75a4aaSAlp Dener ierr = TaoBNKEstimateActiveSet(tao);CHKERRQ(ierr); 37209164190SAlp Dener 373eb3ba6a7SAlp Dener /* Prepare masked matrices for the inactive set */ 374*2f75a4aaSAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr); } 375*2f75a4aaSAlp Dener if (bnk->inactive_idx) { 376eb3ba6a7SAlp Dener ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr); 377eb3ba6a7SAlp Dener if (tao->hessian == tao->hessian_pre) { 378eb3ba6a7SAlp Dener bnk->Hpre_inactive = bnk->H_inactive; 379eb3ba6a7SAlp Dener } else { 380eb3ba6a7SAlp Dener ierr = TaoMatGetSubMat(tao->hessian_pre, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->Hpre_inactive);CHKERRQ(ierr); 381eb3ba6a7SAlp Dener } 382*2f75a4aaSAlp Dener } else { 383*2f75a4aaSAlp Dener bnk->H_inactive = tao->hessian; 384*2f75a4aaSAlp Dener bnk->Hpre_inactive = tao->hessian_pre; 385*2f75a4aaSAlp Dener } 38609164190SAlp Dener 387eb910715SAlp Dener /* Solve the Newton system of equations */ 388*2f75a4aaSAlp Dener ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr); 38909164190SAlp Dener ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr); 390eb910715SAlp Dener if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 391fed79b8eSAlp Dener ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 392eb3ba6a7SAlp Dener ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 393eb910715SAlp Dener ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 394eb910715SAlp Dener tao->ksp_its+=kspits; 395eb910715SAlp Dener tao->ksp_tot_its+=kspits; 396080d2917SAlp Dener ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr); 397eb910715SAlp Dener 398eb910715SAlp Dener if (0.0 == tao->trust) { 399eb910715SAlp Dener /* Radius was uninitialized; use the norm of the direction */ 400080d2917SAlp Dener if (bnk->dnorm > 0.0) { 401080d2917SAlp Dener tao->trust = bnk->dnorm; 402eb910715SAlp Dener 403eb910715SAlp Dener /* Modify the radius if it is too large or small */ 404eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 405eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 406eb910715SAlp Dener } else { 407eb910715SAlp Dener /* The direction was bad; set radius to default value and re-solve 408eb910715SAlp Dener the trust-region subproblem to get a direction */ 409eb910715SAlp Dener tao->trust = tao->trust0; 410eb910715SAlp Dener 411eb910715SAlp Dener /* Modify the radius if it is too large or small */ 412eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 413eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 414eb910715SAlp Dener 415fed79b8eSAlp Dener ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 416eb3ba6a7SAlp Dener ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 417eb910715SAlp Dener ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 418eb910715SAlp Dener tao->ksp_its+=kspits; 419eb910715SAlp Dener tao->ksp_tot_its+=kspits; 420080d2917SAlp Dener ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr); 421eb910715SAlp Dener 422080d2917SAlp Dener if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 423eb910715SAlp Dener } 424eb910715SAlp Dener } 425eb910715SAlp Dener } else { 426eb3ba6a7SAlp Dener ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 427eb910715SAlp Dener ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr); 428eb910715SAlp Dener tao->ksp_its += kspits; 429eb910715SAlp Dener tao->ksp_tot_its+=kspits; 430eb910715SAlp Dener } 431770b7498SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 432fed79b8eSAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 433*2f75a4aaSAlp Dener ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr); 434770b7498SAlp Dener 435770b7498SAlp Dener /* Record convergence reasons */ 436fed79b8eSAlp Dener ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 437770b7498SAlp Dener if (KSP_CONVERGED_ATOL == ksp_reason) { 438770b7498SAlp Dener ++bnk->ksp_atol; 439770b7498SAlp Dener } else if (KSP_CONVERGED_RTOL == ksp_reason) { 440770b7498SAlp Dener ++bnk->ksp_rtol; 441770b7498SAlp Dener } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) { 442770b7498SAlp Dener ++bnk->ksp_ctol; 443770b7498SAlp Dener } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) { 444770b7498SAlp Dener ++bnk->ksp_negc; 445770b7498SAlp Dener } else if (KSP_DIVERGED_DTOL == ksp_reason) { 446770b7498SAlp Dener ++bnk->ksp_dtol; 447770b7498SAlp Dener } else if (KSP_DIVERGED_ITS == ksp_reason) { 448770b7498SAlp Dener ++bnk->ksp_iter; 449770b7498SAlp Dener } else { 450770b7498SAlp Dener ++bnk->ksp_othr; 451770b7498SAlp Dener } 452fed79b8eSAlp Dener 453fed79b8eSAlp Dener /* Make sure the BFGS preconditioner is healthy */ 454fed79b8eSAlp Dener if (bnk->pc_type == BNK_PC_BFGS) { 455fed79b8eSAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 456fed79b8eSAlp Dener if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (bfgsUpdates > 1)) { 457fed79b8eSAlp Dener /* Preconditioner is numerically indefinite; reset the approximation. */ 458eb910715SAlp Dener if (bnk->f != 0.0) { 459eb910715SAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 460eb910715SAlp Dener } else { 461eb910715SAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 462eb910715SAlp Dener } 463eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr); 464eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 46509164190SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 466eb910715SAlp Dener bfgsUpdates = 1; 467eb910715SAlp Dener } 468fed79b8eSAlp Dener } 469eb910715SAlp Dener 470eb910715SAlp Dener /* Check for success (descent direction) */ 471fed79b8eSAlp Dener if (safeguard) { 472080d2917SAlp Dener ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 473eb910715SAlp Dener if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 474eb910715SAlp Dener /* Newton step is not descent or direction produced Inf or NaN 475eb910715SAlp Dener Update the perturbation for next time */ 476eb910715SAlp Dener if (bnk->pert <= 0.0) { 477eb910715SAlp Dener /* Initialize the perturbation */ 478eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 479eb910715SAlp Dener if (bnk->is_gltr) { 480eb910715SAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 481eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 482eb910715SAlp Dener } 483eb910715SAlp Dener } else { 484eb910715SAlp Dener /* Increase the perturbation */ 485eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 486eb910715SAlp Dener } 487eb910715SAlp Dener 488eb910715SAlp Dener if (BNK_PC_BFGS != bnk->pc_type) { 489eb910715SAlp Dener /* We don't have the bfgs matrix around and updated 490eb910715SAlp Dener Must use gradient direction in this case */ 491080d2917SAlp Dener ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 492eb910715SAlp Dener ++bnk->grad; 493eb910715SAlp Dener *stepType = BNK_GRADIENT; 494eb910715SAlp Dener } else { 495eb910715SAlp Dener /* Attempt to use the BFGS direction */ 49609164190SAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 497eb910715SAlp Dener 4988d5ead36SAlp Dener /* Check for success (descent direction) 4998d5ead36SAlp Dener NOTE: Negative gdx here means not a descent direction because 5008d5ead36SAlp Dener the fall-back step is missing a negative sign. */ 501080d2917SAlp Dener ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 5028d5ead36SAlp Dener if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 503eb910715SAlp Dener /* BFGS direction is not descent or direction produced not a number 504eb910715SAlp Dener We can assert bfgsUpdates > 1 in this case because 505eb910715SAlp Dener the first solve produces the scaled gradient direction, 506eb910715SAlp Dener which is guaranteed to be descent */ 507eb910715SAlp Dener 508eb910715SAlp Dener /* Use steepest descent direction (scaled) */ 509eb910715SAlp Dener if (bnk->f != 0.0) { 510eb910715SAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 511eb910715SAlp Dener } else { 512eb910715SAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 513eb910715SAlp Dener } 514eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 515eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 51609164190SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 51709164190SAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 518eb910715SAlp Dener 519eb910715SAlp Dener bfgsUpdates = 1; 520eb910715SAlp Dener ++bnk->sgrad; 521eb910715SAlp Dener *stepType = BNK_SCALED_GRADIENT; 522eb910715SAlp Dener } else { 523770b7498SAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 524eb910715SAlp Dener if (1 == bfgsUpdates) { 525eb910715SAlp Dener /* The first BFGS direction is always the scaled gradient */ 526eb910715SAlp Dener ++bnk->sgrad; 527eb910715SAlp Dener *stepType = BNK_SCALED_GRADIENT; 528eb910715SAlp Dener } else { 529eb910715SAlp Dener ++bnk->bfgs; 530eb910715SAlp Dener *stepType = BNK_BFGS; 531eb910715SAlp Dener } 532eb910715SAlp Dener } 533eb910715SAlp Dener } 5348d5ead36SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 5358d5ead36SAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 536*2f75a4aaSAlp Dener ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr); 537eb910715SAlp Dener } else { 538eb910715SAlp Dener /* Computed Newton step is descent */ 539eb910715SAlp Dener switch (ksp_reason) { 540eb910715SAlp Dener case KSP_DIVERGED_NANORINF: 541eb910715SAlp Dener case KSP_DIVERGED_BREAKDOWN: 542eb910715SAlp Dener case KSP_DIVERGED_INDEFINITE_MAT: 543eb910715SAlp Dener case KSP_DIVERGED_INDEFINITE_PC: 544eb910715SAlp Dener case KSP_CONVERGED_CG_NEG_CURVE: 545eb910715SAlp Dener /* Matrix or preconditioner is indefinite; increase perturbation */ 546eb910715SAlp Dener if (bnk->pert <= 0.0) { 547eb910715SAlp Dener /* Initialize the perturbation */ 548eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 549eb910715SAlp Dener if (bnk->is_gltr) { 550eb910715SAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr); 551eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 552eb910715SAlp Dener } 553eb910715SAlp Dener } else { 554eb910715SAlp Dener /* Increase the perturbation */ 555eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 556eb910715SAlp Dener } 557eb910715SAlp Dener break; 558eb910715SAlp Dener 559eb910715SAlp Dener default: 560eb910715SAlp Dener /* Newton step computation is good; decrease perturbation */ 561eb910715SAlp Dener bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm); 562eb910715SAlp Dener if (bnk->pert < bnk->pmin) { 563eb910715SAlp Dener bnk->pert = 0.0; 564eb910715SAlp Dener } 565eb910715SAlp Dener break; 566eb910715SAlp Dener } 567eb910715SAlp Dener 568eb910715SAlp Dener ++bnk->newt; 569fed79b8eSAlp Dener *stepType = BNK_NEWTON; 570fed79b8eSAlp Dener } 571fed79b8eSAlp Dener } else { 572fed79b8eSAlp Dener ++bnk->newt; 573fed79b8eSAlp Dener bnk->pert = 0.0; 574fed79b8eSAlp Dener *stepType = BNK_NEWTON; 575eb910715SAlp Dener } 576eb910715SAlp Dener PetscFunctionReturn(0); 577eb910715SAlp Dener } 578eb910715SAlp Dener 579df278d8fSAlp Dener /*------------------------------------------------------------*/ 580df278d8fSAlp Dener 581df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search. 582df278d8fSAlp Dener 583df278d8fSAlp Dener Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the 584df278d8fSAlp Dener Newton step does not produce a valid step length. 585df278d8fSAlp Dener */ 586df278d8fSAlp Dener 587c14b763aSAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason) 588c14b763aSAlp Dener { 589c14b763aSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 590c14b763aSAlp Dener PetscErrorCode ierr; 591c14b763aSAlp Dener TaoLineSearchConvergedReason ls_reason; 592c14b763aSAlp Dener 593c14b763aSAlp Dener PetscReal e_min, gdx, delta; 594c14b763aSAlp Dener PetscInt bfgsUpdates; 595c14b763aSAlp Dener 596c14b763aSAlp Dener PetscFunctionBegin; 597c14b763aSAlp Dener /* Perform the linesearch */ 598770b7498SAlp Dener *steplen = 1.0; 599c14b763aSAlp Dener ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr); 600c14b763aSAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 601c14b763aSAlp Dener ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 602c14b763aSAlp Dener 603c14b763aSAlp Dener while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != BNK_GRADIENT) { 604c14b763aSAlp Dener /* Linesearch failed, revert solution */ 605c14b763aSAlp Dener bnk->f = bnk->fold; 606c14b763aSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 607c14b763aSAlp Dener ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 608c14b763aSAlp Dener ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 609c14b763aSAlp Dener 610c14b763aSAlp Dener switch(stepType) { 611c14b763aSAlp Dener case BNK_NEWTON: 6128d5ead36SAlp Dener /* Failed to obtain acceptable iterate with Newton step 613c14b763aSAlp Dener Update the perturbation for next time */ 614c14b763aSAlp Dener --bnk->newt; 615c14b763aSAlp Dener if (bnk->pert <= 0.0) { 616c14b763aSAlp Dener /* Initialize the perturbation */ 617c14b763aSAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 618c14b763aSAlp Dener if (bnk->is_gltr) { 619c14b763aSAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 620c14b763aSAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 621c14b763aSAlp Dener } 622c14b763aSAlp Dener } else { 623c14b763aSAlp Dener /* Increase the perturbation */ 624c14b763aSAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 625c14b763aSAlp Dener } 626c14b763aSAlp Dener 627c14b763aSAlp Dener if (BNK_PC_BFGS != bnk->pc_type) { 628c14b763aSAlp Dener /* We don't have the bfgs matrix around and being updated 629c14b763aSAlp Dener Must use gradient direction in this case */ 630c14b763aSAlp Dener ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 631c14b763aSAlp Dener ++bnk->grad; 632c14b763aSAlp Dener stepType = BNK_GRADIENT; 633c14b763aSAlp Dener } else { 634c14b763aSAlp Dener /* Attempt to use the BFGS direction */ 635c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 6368d5ead36SAlp Dener /* Check for success (descent direction) 6378d5ead36SAlp Dener NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */ 638c14b763aSAlp Dener ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 639c14b763aSAlp Dener if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 640c14b763aSAlp Dener /* BFGS direction is not descent or direction produced not a number 641c14b763aSAlp Dener We can assert bfgsUpdates > 1 in this case 642c14b763aSAlp Dener Use steepest descent direction (scaled) */ 643c14b763aSAlp Dener if (bnk->f != 0.0) { 644c14b763aSAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 645c14b763aSAlp Dener } else { 646c14b763aSAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 647c14b763aSAlp Dener } 648c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 649c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 650c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 651c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 652c14b763aSAlp Dener 653c14b763aSAlp Dener bfgsUpdates = 1; 654c14b763aSAlp Dener ++bnk->sgrad; 655c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 656c14b763aSAlp Dener } else { 657c14b763aSAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 658c14b763aSAlp Dener if (1 == bfgsUpdates) { 659c14b763aSAlp Dener /* The first BFGS direction is always the scaled gradient */ 660c14b763aSAlp Dener ++bnk->sgrad; 661c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 662c14b763aSAlp Dener } else { 663c14b763aSAlp Dener ++bnk->bfgs; 664c14b763aSAlp Dener stepType = BNK_BFGS; 665c14b763aSAlp Dener } 666c14b763aSAlp Dener } 667c14b763aSAlp Dener } 668c14b763aSAlp Dener break; 669c14b763aSAlp Dener 670c14b763aSAlp Dener case BNK_BFGS: 671c14b763aSAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS 672c14b763aSAlp Dener Failed to obtain acceptable iterate with BFGS step 673c14b763aSAlp Dener Attempt to use the scaled gradient direction */ 674c14b763aSAlp Dener if (bnk->f != 0.0) { 675c14b763aSAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 676c14b763aSAlp Dener } else { 677c14b763aSAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 678c14b763aSAlp Dener } 679c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 680c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 681c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 682c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 683c14b763aSAlp Dener 684c14b763aSAlp Dener bfgsUpdates = 1; 685c14b763aSAlp Dener ++bnk->sgrad; 686c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 687c14b763aSAlp Dener break; 688c14b763aSAlp Dener 689c14b763aSAlp Dener case BNK_SCALED_GRADIENT: 690c14b763aSAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS 691c14b763aSAlp Dener The scaled gradient step did not produce a new iterate; 692c14b763aSAlp Dener attemp to use the gradient direction. 693c14b763aSAlp Dener Need to make sure we are not using a different diagonal scaling */ 694c14b763aSAlp Dener ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr); 695c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr); 696c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 697c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 698c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 699c14b763aSAlp Dener 700c14b763aSAlp Dener bfgsUpdates = 1; 701c14b763aSAlp Dener ++bnk->grad; 702c14b763aSAlp Dener stepType = BNK_GRADIENT; 703c14b763aSAlp Dener break; 704c14b763aSAlp Dener } 7058d5ead36SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 7068d5ead36SAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 707*2f75a4aaSAlp Dener ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr); 708c14b763aSAlp Dener 7098d5ead36SAlp Dener /* Perform one last line search with the fall-back step */ 710770b7498SAlp Dener *steplen = 1.0; 711c14b763aSAlp Dener ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr); 712c14b763aSAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 713c14b763aSAlp Dener ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 714c14b763aSAlp Dener } 715c14b763aSAlp Dener *reason = ls_reason; 716c14b763aSAlp Dener PetscFunctionReturn(0); 717c14b763aSAlp Dener } 718c14b763aSAlp Dener 719df278d8fSAlp Dener /*------------------------------------------------------------*/ 720df278d8fSAlp Dener 721df278d8fSAlp Dener /* Routine for updating the trust radius. 722df278d8fSAlp Dener 723df278d8fSAlp Dener Function features three different update methods: 724df278d8fSAlp Dener 1) Line-search step length based 725df278d8fSAlp Dener 2) Predicted decrease on the CG quadratic model 726df278d8fSAlp Dener 3) Interpolation 727df278d8fSAlp Dener */ 728df278d8fSAlp Dener 729b1c2d0e3SAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt stepType, PetscBool *accept) 730080d2917SAlp Dener { 731080d2917SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 732080d2917SAlp Dener PetscErrorCode ierr; 733080d2917SAlp Dener 734b1c2d0e3SAlp Dener PetscReal step, kappa; 735080d2917SAlp Dener PetscReal gdx, tau_1, tau_2, tau_min, tau_max; 736080d2917SAlp Dener 737080d2917SAlp Dener PetscFunctionBegin; 738080d2917SAlp Dener /* Update trust region radius */ 739080d2917SAlp Dener *accept = PETSC_FALSE; 740080d2917SAlp Dener switch(bnk->update_type) { 741080d2917SAlp Dener case BNK_UPDATE_STEP: 742c14b763aSAlp Dener *accept = PETSC_TRUE; /* always accept here because line search succeeded */ 743080d2917SAlp Dener if (stepType == BNK_NEWTON) { 744080d2917SAlp Dener ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr); 745080d2917SAlp Dener if (step < bnk->nu1) { 746080d2917SAlp Dener /* Very bad step taken; reduce radius */ 747080d2917SAlp Dener tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 748080d2917SAlp Dener } else if (step < bnk->nu2) { 749080d2917SAlp Dener /* Reasonably bad step taken; reduce radius */ 750080d2917SAlp Dener tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust); 751080d2917SAlp Dener } else if (step < bnk->nu3) { 752080d2917SAlp Dener /* Reasonable step was taken; leave radius alone */ 753080d2917SAlp Dener if (bnk->omega3 < 1.0) { 754080d2917SAlp Dener tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust); 755080d2917SAlp Dener } else if (bnk->omega3 > 1.0) { 756080d2917SAlp Dener tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust); 757080d2917SAlp Dener } 758080d2917SAlp Dener } else if (step < bnk->nu4) { 759080d2917SAlp Dener /* Full step taken; increase the radius */ 760080d2917SAlp Dener tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust); 761080d2917SAlp Dener } else { 762080d2917SAlp Dener /* More than full step taken; increase the radius */ 763080d2917SAlp Dener tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust); 764080d2917SAlp Dener } 765080d2917SAlp Dener } else { 766080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 767080d2917SAlp Dener tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 768080d2917SAlp Dener } 769080d2917SAlp Dener break; 770080d2917SAlp Dener 771080d2917SAlp Dener case BNK_UPDATE_REDUCTION: 772080d2917SAlp Dener if (stepType == BNK_NEWTON) { 773b1c2d0e3SAlp Dener if (prered < 0.0) { 774fed79b8eSAlp Dener /* The predicted reduction has the wrong sign. This cannot 775fed79b8eSAlp Dener happen in infinite precision arithmetic. Step should 776fed79b8eSAlp Dener be rejected! */ 777080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 778fed79b8eSAlp Dener } 779fed79b8eSAlp Dener else { 780b1c2d0e3SAlp Dener if (PetscIsInfOrNanReal(actred)) { 781080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 782080d2917SAlp Dener } else { 783080d2917SAlp Dener if ((PetscAbsScalar(actred) <= bnk->epsilon) && 784080d2917SAlp Dener (PetscAbsScalar(prered) <= bnk->epsilon)) { 785080d2917SAlp Dener kappa = 1.0; 786fed79b8eSAlp Dener } 787fed79b8eSAlp Dener else { 788080d2917SAlp Dener kappa = actred / prered; 789080d2917SAlp Dener } 790fed79b8eSAlp Dener 791fed79b8eSAlp Dener /* Accept or reject the step and update radius */ 792080d2917SAlp Dener if (kappa < bnk->eta1) { 793fed79b8eSAlp Dener /* Reject the step */ 794080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 795fed79b8eSAlp Dener } 796fed79b8eSAlp Dener else { 797fed79b8eSAlp Dener /* Accept the step */ 798c133c014SAlp Dener *accept = PETSC_TRUE; 799c133c014SAlp Dener /* Update the trust region radius only if the computed step is at the trust radius boundary */ 8008d5ead36SAlp Dener if (bnk->dnorm == tao->trust) { 801080d2917SAlp Dener if (kappa < bnk->eta2) { 802080d2917SAlp Dener /* Marginal bad step */ 803c133c014SAlp Dener tao->trust = bnk->alpha2 * tao->trust; 804080d2917SAlp Dener } 805fed79b8eSAlp Dener else if (kappa < bnk->eta3) { 806fed79b8eSAlp Dener /* Reasonable step */ 807fed79b8eSAlp Dener tao->trust = bnk->alpha3 * tao->trust; 808fed79b8eSAlp Dener } 809fed79b8eSAlp Dener else if (kappa < bnk->eta4) { 810080d2917SAlp Dener /* Good step */ 811c133c014SAlp Dener tao->trust = bnk->alpha4 * tao->trust; 812fed79b8eSAlp Dener } 813fed79b8eSAlp Dener else { 814080d2917SAlp Dener /* Very good step */ 815c133c014SAlp Dener tao->trust = bnk->alpha5 * tao->trust; 816080d2917SAlp Dener } 817c133c014SAlp Dener } 818080d2917SAlp Dener } 819080d2917SAlp Dener } 820080d2917SAlp Dener } 821080d2917SAlp Dener } else { 822080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 823080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust); 824080d2917SAlp Dener } 825080d2917SAlp Dener break; 826080d2917SAlp Dener 827080d2917SAlp Dener default: 828080d2917SAlp Dener if (stepType == BNK_NEWTON) { 829b1c2d0e3SAlp Dener if (prered < 0.0) { 830080d2917SAlp Dener /* The predicted reduction has the wrong sign. This cannot */ 831080d2917SAlp Dener /* happen in infinite precision arithmetic. Step should */ 832080d2917SAlp Dener /* be rejected! */ 833080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 834080d2917SAlp Dener } else { 835b1c2d0e3SAlp Dener if (PetscIsInfOrNanReal(actred)) { 836080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 837080d2917SAlp Dener } else { 838080d2917SAlp Dener if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 839080d2917SAlp Dener kappa = 1.0; 840080d2917SAlp Dener } else { 841080d2917SAlp Dener kappa = actred / prered; 842080d2917SAlp Dener } 843080d2917SAlp Dener 844080d2917SAlp Dener ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 845080d2917SAlp Dener tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred); 846080d2917SAlp Dener tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred); 847080d2917SAlp Dener tau_min = PetscMin(tau_1, tau_2); 848080d2917SAlp Dener tau_max = PetscMax(tau_1, tau_2); 849080d2917SAlp Dener 850080d2917SAlp Dener if (kappa >= 1.0 - bnk->mu1) { 851080d2917SAlp Dener /* Great agreement */ 852080d2917SAlp Dener *accept = PETSC_TRUE; 853080d2917SAlp Dener if (tau_max < 1.0) { 854080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 855080d2917SAlp Dener } else if (tau_max > bnk->gamma4) { 856080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm); 857080d2917SAlp Dener } else { 858080d2917SAlp Dener tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 859080d2917SAlp Dener } 860080d2917SAlp Dener } else if (kappa >= 1.0 - bnk->mu2) { 861080d2917SAlp Dener /* Good agreement */ 862080d2917SAlp Dener *accept = PETSC_TRUE; 863080d2917SAlp Dener if (tau_max < bnk->gamma2) { 864080d2917SAlp Dener tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 865080d2917SAlp Dener } else if (tau_max > bnk->gamma3) { 866080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 867080d2917SAlp Dener } else if (tau_max < 1.0) { 868080d2917SAlp Dener tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 869080d2917SAlp Dener } else { 870080d2917SAlp Dener tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 871080d2917SAlp Dener } 872080d2917SAlp Dener } else { 873080d2917SAlp Dener /* Not good agreement */ 874080d2917SAlp Dener if (tau_min > 1.0) { 875080d2917SAlp Dener tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 876080d2917SAlp Dener } else if (tau_max < bnk->gamma1) { 877080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 878080d2917SAlp Dener } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) { 879080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 880080d2917SAlp Dener } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) { 881080d2917SAlp Dener tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm); 882080d2917SAlp Dener } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) { 883080d2917SAlp Dener tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm); 884080d2917SAlp Dener } else { 885080d2917SAlp Dener tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 886080d2917SAlp Dener } 887080d2917SAlp Dener } 888080d2917SAlp Dener } 889080d2917SAlp Dener } 890080d2917SAlp Dener } else { 891080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 892080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust); 893080d2917SAlp Dener } 894080d2917SAlp Dener } 895c133c014SAlp Dener /* Make sure the radius does not violate min and max settings */ 896c133c014SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 897fed79b8eSAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 898080d2917SAlp Dener PetscFunctionReturn(0); 899080d2917SAlp Dener } 900080d2917SAlp Dener 901eb910715SAlp Dener /* ---------------------------------------------------------- */ 902df278d8fSAlp Dener 903eb910715SAlp Dener static PetscErrorCode TaoSetUp_BNK(Tao tao) 904eb910715SAlp Dener { 905eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 906eb910715SAlp Dener PetscErrorCode ierr; 907eb910715SAlp Dener 908eb910715SAlp Dener PetscFunctionBegin; 909eb910715SAlp Dener if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);} 910eb910715SAlp Dener if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);} 911eb910715SAlp Dener if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);} 912eb910715SAlp Dener if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);} 913eb910715SAlp Dener if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);} 91409164190SAlp Dener if (!bnk->Xwork) {ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);} 91509164190SAlp Dener if (!bnk->Gwork) {ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);} 91609164190SAlp Dener if (!bnk->unprojected_gradient) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);} 91709164190SAlp Dener if (!bnk->unprojected_gradient_old) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);} 918eb910715SAlp Dener bnk->Diag = 0; 919eb910715SAlp Dener bnk->M = 0; 92009164190SAlp Dener bnk->H_inactive = 0; 92109164190SAlp Dener bnk->Hpre_inactive = 0; 922eb910715SAlp Dener PetscFunctionReturn(0); 923eb910715SAlp Dener } 924eb910715SAlp Dener 925eb910715SAlp Dener /*------------------------------------------------------------*/ 926df278d8fSAlp Dener 927eb910715SAlp Dener static PetscErrorCode TaoDestroy_BNK(Tao tao) 928eb910715SAlp Dener { 929eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 930eb910715SAlp Dener PetscErrorCode ierr; 931eb910715SAlp Dener 932eb910715SAlp Dener PetscFunctionBegin; 933eb910715SAlp Dener if (tao->setupcalled) { 934eb910715SAlp Dener ierr = VecDestroy(&bnk->W);CHKERRQ(ierr); 935eb910715SAlp Dener ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr); 936eb910715SAlp Dener ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr); 93709164190SAlp Dener ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr); 93809164190SAlp Dener ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr); 93909164190SAlp Dener ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr); 94009164190SAlp Dener ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr); 941eb910715SAlp Dener } 942eb910715SAlp Dener ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr); 943eb910715SAlp Dener ierr = MatDestroy(&bnk->M);CHKERRQ(ierr); 944eb910715SAlp Dener ierr = PetscFree(tao->data);CHKERRQ(ierr); 945eb910715SAlp Dener PetscFunctionReturn(0); 946eb910715SAlp Dener } 947eb910715SAlp Dener 948eb910715SAlp Dener /*------------------------------------------------------------*/ 949df278d8fSAlp Dener 950eb910715SAlp Dener static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao) 951eb910715SAlp Dener { 952eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 953eb910715SAlp Dener PetscErrorCode ierr; 954eb910715SAlp Dener 955eb910715SAlp Dener PetscFunctionBegin; 956eb910715SAlp Dener ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr); 9578d5ead36SAlp 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); 9588d5ead36SAlp 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); 9598d5ead36SAlp 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); 9608d5ead36SAlp 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); 961*2f75a4aaSAlp 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); 9628d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr); 9638d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr); 9648d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr); 9658d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr); 9668d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr); 9678d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr); 9688d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr); 9698d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr); 9708d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr); 9718d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr); 9728d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr); 9738d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr); 9748d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr); 9758d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr); 9768d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr); 9778d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr); 9788d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr); 9798d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr); 9808d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr); 9818d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr); 9828d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr); 9838d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr); 9848d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr); 9858d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr); 9868d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr); 9878d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr); 9888d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr); 9898d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr); 9908d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr); 9918d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr); 9928d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr); 9938d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr); 9948d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr); 9958d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr); 9968d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr); 9978d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr); 9988d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr); 9998d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr); 10008d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr); 10018d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr); 10028d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr); 10038d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr); 10048d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr); 10058d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr); 10068d5ead36SAlp Dener ierr = PetscOptionsReal("-tao_bnk_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr); 1007*2f75a4aaSAlp Dener ierr = PetscOptionsReal("-tao_bnk_bound_tol", "initial tolerance used when estimating actively bounded variables", "", bnk->bound_tol, &bnk->bound_tol,NULL);CHKERRQ(ierr); 1008eb910715SAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 1009eb910715SAlp Dener ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 1010eb910715SAlp Dener ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 1011eb910715SAlp Dener PetscFunctionReturn(0); 1012eb910715SAlp Dener } 1013eb910715SAlp Dener 1014eb910715SAlp Dener /*------------------------------------------------------------*/ 1015df278d8fSAlp Dener 1016eb910715SAlp Dener static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer) 1017eb910715SAlp Dener { 1018eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 1019eb910715SAlp Dener PetscInt nrejects; 1020eb910715SAlp Dener PetscBool isascii; 1021eb910715SAlp Dener PetscErrorCode ierr; 1022eb910715SAlp Dener 1023eb910715SAlp Dener PetscFunctionBegin; 1024eb910715SAlp Dener ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 1025eb910715SAlp Dener if (isascii) { 1026eb910715SAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1027eb910715SAlp Dener if (BNK_PC_BFGS == bnk->pc_type && bnk->M) { 1028eb910715SAlp Dener ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr); 1029eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr); 1030eb910715SAlp Dener } 1031eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr); 1032eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr); 1033eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr); 1034eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr); 1035eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr); 1036eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr); 1037eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr); 1038eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr); 1039eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr); 1040eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr); 1041eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr); 1042eb910715SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr); 1043eb910715SAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1044eb910715SAlp Dener } 1045eb910715SAlp Dener PetscFunctionReturn(0); 1046eb910715SAlp Dener } 1047eb910715SAlp Dener 1048eb910715SAlp Dener /* ---------------------------------------------------------- */ 1049df278d8fSAlp Dener 1050eb910715SAlp Dener /*MC 1051eb910715SAlp Dener TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms. 105266ed3702SAlp Dener At each iteration, the BNK methods solve the symmetric 1053eb910715SAlp Dener system of equations to obtain the step diretion dk: 1054eb910715SAlp Dener Hk dk = -gk 10552b97c8d8SAlp Dener for free variables only. The step can be globalized either through 10562b97c8d8SAlp Dener trust-region methods, or a line search, or a heuristic mixture of both. 1057eb910715SAlp Dener 1058eb910715SAlp Dener Options Database Keys: 10598d5ead36SAlp Dener + -tao_bnk_pc_type - "none","ahess","bfgs","petsc" 10608d5ead36SAlp Dener . -tao_bnk_bfgs_scale_type - "ahess","phess","bfgs" 10618d5ead36SAlp Dener . -tao_bnk_init_type - "constant","direction","interpolation" 10628d5ead36SAlp Dener . -tao_bnk_update_type - "step","direction","interpolation" 1063*2f75a4aaSAlp Dener . -tao_bnk_as_type - "none","bertsekas" 10648d5ead36SAlp Dener . -tao_bnk_sval - perturbation starting value 10658d5ead36SAlp Dener . -tao_bnk_imin - minimum initial perturbation 10668d5ead36SAlp Dener . -tao_bnk_imax - maximum initial perturbation 10678d5ead36SAlp Dener . -tao_bnk_pmin - minimum perturbation 10688d5ead36SAlp Dener . -tao_bnk_pmax - maximum perturbation 10698d5ead36SAlp Dener . -tao_bnk_pgfac - growth factor 10708d5ead36SAlp Dener . -tao_bnk_psfac - shrink factor 10718d5ead36SAlp Dener . -tao_bnk_imfac - initial merit factor 10728d5ead36SAlp Dener . -tao_bnk_pmgfac - merit growth factor 10738d5ead36SAlp Dener . -tao_bnk_pmsfac - merit shrink factor 10748d5ead36SAlp Dener . -tao_bnk_eta1 - poor steplength; reduce radius 10758d5ead36SAlp Dener . -tao_bnk_eta2 - reasonable steplength; leave radius 10768d5ead36SAlp Dener . -tao_bnk_eta3 - good steplength; increase readius 10778d5ead36SAlp Dener . -tao_bnk_eta4 - excellent steplength; greatly increase radius 10788d5ead36SAlp Dener . -tao_bnk_alpha1 - alpha1 reduction 10798d5ead36SAlp Dener . -tao_bnk_alpha2 - alpha2 reduction 10808d5ead36SAlp Dener . -tao_bnk_alpha3 - alpha3 reduction 10818d5ead36SAlp Dener . -tao_bnk_alpha4 - alpha4 reduction 10828d5ead36SAlp Dener . -tao_bnk_alpha - alpha5 reduction 10838d5ead36SAlp Dener . -tao_bnk_mu1 - mu1 interpolation update 10848d5ead36SAlp Dener . -tao_bnk_mu2 - mu2 interpolation update 10858d5ead36SAlp Dener . -tao_bnk_gamma1 - gamma1 interpolation update 10868d5ead36SAlp Dener . -tao_bnk_gamma2 - gamma2 interpolation update 10878d5ead36SAlp Dener . -tao_bnk_gamma3 - gamma3 interpolation update 10888d5ead36SAlp Dener . -tao_bnk_gamma4 - gamma4 interpolation update 10898d5ead36SAlp Dener . -tao_bnk_theta - theta interpolation update 10908d5ead36SAlp Dener . -tao_bnk_omega1 - omega1 step update 10918d5ead36SAlp Dener . -tao_bnk_omega2 - omega2 step update 10928d5ead36SAlp Dener . -tao_bnk_omega3 - omega3 step update 10938d5ead36SAlp Dener . -tao_bnk_omega4 - omega4 step update 10948d5ead36SAlp Dener . -tao_bnk_omega5 - omega5 step update 10958d5ead36SAlp Dener . -tao_bnk_mu1_i - mu1 interpolation init factor 10968d5ead36SAlp Dener . -tao_bnk_mu2_i - mu2 interpolation init factor 10978d5ead36SAlp Dener . -tao_bnk_gamma1_i - gamma1 interpolation init factor 10988d5ead36SAlp Dener . -tao_bnk_gamma2_i - gamma2 interpolation init factor 10998d5ead36SAlp Dener . -tao_bnk_gamma3_i - gamma3 interpolation init factor 11008d5ead36SAlp Dener . -tao_bnk_gamma4_i - gamma4 interpolation init factor 1101*2f75a4aaSAlp Dener . -tao_bnk_theta_i - theta interpolation init factor 1102*2f75a4aaSAlp Dener - -tao_bnk_bound_tol - initial tolerance used in estimating bounded active variables 1103eb910715SAlp Dener 1104eb910715SAlp Dener Level: beginner 1105eb910715SAlp Dener M*/ 1106eb910715SAlp Dener 1107eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao) 1108eb910715SAlp Dener { 1109eb910715SAlp Dener TAO_BNK *bnk; 1110eb910715SAlp Dener const char *morethuente_type = TAOLINESEARCHMT; 1111eb910715SAlp Dener PetscErrorCode ierr; 1112eb910715SAlp Dener 1113eb910715SAlp Dener PetscFunctionBegin; 1114eb910715SAlp Dener ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr); 1115eb910715SAlp Dener 1116eb910715SAlp Dener tao->ops->setup = TaoSetUp_BNK; 1117eb910715SAlp Dener tao->ops->view = TaoView_BNK; 1118eb910715SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BNK; 1119eb910715SAlp Dener tao->ops->destroy = TaoDestroy_BNK; 1120eb910715SAlp Dener 1121eb910715SAlp Dener /* Override default settings (unless already changed) */ 1122eb910715SAlp Dener if (!tao->max_it_changed) tao->max_it = 50; 1123eb910715SAlp Dener if (!tao->trust0_changed) tao->trust0 = 100.0; 1124eb910715SAlp Dener 1125eb910715SAlp Dener tao->data = (void*)bnk; 1126eb910715SAlp Dener 112766ed3702SAlp Dener /* Hessian shifting parameters */ 1128eb910715SAlp Dener bnk->sval = 0.0; 1129eb910715SAlp Dener bnk->imin = 1.0e-4; 1130eb910715SAlp Dener bnk->imax = 1.0e+2; 1131eb910715SAlp Dener bnk->imfac = 1.0e-1; 1132eb910715SAlp Dener 1133eb910715SAlp Dener bnk->pmin = 1.0e-12; 1134eb910715SAlp Dener bnk->pmax = 1.0e+2; 1135eb910715SAlp Dener bnk->pgfac = 1.0e+1; 1136eb910715SAlp Dener bnk->psfac = 4.0e-1; 1137eb910715SAlp Dener bnk->pmgfac = 1.0e-1; 1138eb910715SAlp Dener bnk->pmsfac = 1.0e-1; 1139eb910715SAlp Dener 1140eb910715SAlp Dener /* Default values for trust-region radius update based on steplength */ 1141eb910715SAlp Dener bnk->nu1 = 0.25; 1142eb910715SAlp Dener bnk->nu2 = 0.50; 1143eb910715SAlp Dener bnk->nu3 = 1.00; 1144eb910715SAlp Dener bnk->nu4 = 1.25; 1145eb910715SAlp Dener 1146eb910715SAlp Dener bnk->omega1 = 0.25; 1147eb910715SAlp Dener bnk->omega2 = 0.50; 1148eb910715SAlp Dener bnk->omega3 = 1.00; 1149eb910715SAlp Dener bnk->omega4 = 2.00; 1150eb910715SAlp Dener bnk->omega5 = 4.00; 1151eb910715SAlp Dener 1152eb910715SAlp Dener /* Default values for trust-region radius update based on reduction */ 1153eb910715SAlp Dener bnk->eta1 = 1.0e-4; 1154eb910715SAlp Dener bnk->eta2 = 0.25; 1155eb910715SAlp Dener bnk->eta3 = 0.50; 1156eb910715SAlp Dener bnk->eta4 = 0.90; 1157eb910715SAlp Dener 1158eb910715SAlp Dener bnk->alpha1 = 0.25; 1159eb910715SAlp Dener bnk->alpha2 = 0.50; 1160eb910715SAlp Dener bnk->alpha3 = 1.00; 1161eb910715SAlp Dener bnk->alpha4 = 2.00; 1162eb910715SAlp Dener bnk->alpha5 = 4.00; 1163eb910715SAlp Dener 1164eb910715SAlp Dener /* Default values for trust-region radius update based on interpolation */ 1165eb910715SAlp Dener bnk->mu1 = 0.10; 1166eb910715SAlp Dener bnk->mu2 = 0.50; 1167eb910715SAlp Dener 1168eb910715SAlp Dener bnk->gamma1 = 0.25; 1169eb910715SAlp Dener bnk->gamma2 = 0.50; 1170eb910715SAlp Dener bnk->gamma3 = 2.00; 1171eb910715SAlp Dener bnk->gamma4 = 4.00; 1172eb910715SAlp Dener 1173eb910715SAlp Dener bnk->theta = 0.05; 1174eb910715SAlp Dener 1175eb910715SAlp Dener /* Default values for trust region initialization based on interpolation */ 1176eb910715SAlp Dener bnk->mu1_i = 0.35; 1177eb910715SAlp Dener bnk->mu2_i = 0.50; 1178eb910715SAlp Dener 1179eb910715SAlp Dener bnk->gamma1_i = 0.0625; 1180eb910715SAlp Dener bnk->gamma2_i = 0.5; 1181eb910715SAlp Dener bnk->gamma3_i = 2.0; 1182eb910715SAlp Dener bnk->gamma4_i = 5.0; 1183eb910715SAlp Dener 1184eb910715SAlp Dener bnk->theta_i = 0.25; 1185eb910715SAlp Dener 1186eb910715SAlp Dener /* Remaining parameters */ 1187eb910715SAlp Dener bnk->min_radius = 1.0e-10; 1188eb910715SAlp Dener bnk->max_radius = 1.0e10; 1189770b7498SAlp Dener bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); 1190*2f75a4aaSAlp Dener bnk->bound_tol = 1.0e-3; 1191eb910715SAlp Dener 1192eb910715SAlp Dener bnk->pc_type = BNK_PC_BFGS; 1193eb910715SAlp Dener bnk->bfgs_scale_type = BFGS_SCALE_PHESS; 1194eb910715SAlp Dener bnk->init_type = BNK_INIT_INTERPOLATION; 1195fed79b8eSAlp Dener bnk->update_type = BNK_UPDATE_INTERPOLATION; 1196*2f75a4aaSAlp Dener bnk->as_type = BNK_AS_BERTSEKAS; 1197eb910715SAlp Dener 1198eb910715SAlp Dener ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 1199eb910715SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 1200eb910715SAlp Dener ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 1201eb910715SAlp Dener ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 1202eb910715SAlp Dener ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 1203eb910715SAlp Dener 1204eb910715SAlp Dener /* Set linear solver to default for symmetric matrices */ 1205eb910715SAlp Dener ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 1206eb910715SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 1207eb910715SAlp Dener ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr); 1208eb910715SAlp Dener ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr); 1209eb910715SAlp Dener PetscFunctionReturn(0); 1210eb910715SAlp Dener } 1211