1eb910715SAlp Dener #include <petsctaolinesearch.h> 2eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h> 3eb910715SAlp Dener #include <petscksp.h> 4eb910715SAlp Dener 570a3f44bSAlp Dener static const char *BNK_INIT[64] = {"constant", "direction", "interpolation"}; 670a3f44bSAlp Dener static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"}; 770a3f44bSAlp Dener static const char *BNK_AS[64] = {"none", "bertsekas"}; 870a3f44bSAlp Dener 9e031d6f5SAlp Dener /*------------------------------------------------------------*/ 10e031d6f5SAlp Dener 11df278d8fSAlp Dener /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */ 12df278d8fSAlp Dener 13c0f10754SAlp Dener PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH) 14eb910715SAlp Dener { 15eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 16eb910715SAlp Dener PC pc; 1789da521bSAlp Dener PetscReal f_min, ftrial, prered, actred, kappa, sigma, resnorm; 18eb910715SAlp Dener PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 190ad3a497SAlp Dener PetscBool is_bfgs, is_jacobi, is_symmetric, sym_set; 20c4b75bccSAlp Dener PetscInt n, N, nDiff; 21eb910715SAlp Dener PetscInt i_max = 5; 22eb910715SAlp Dener PetscInt j_max = 1; 23eb910715SAlp Dener PetscInt i, j; 242e6e4ca1SStefano Zampini PetscVoidFunction kspTR; 25eb910715SAlp Dener 26eb910715SAlp Dener PetscFunctionBegin; 2728017e9fSAlp Dener /* Project the current point onto the feasible set */ 289566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 299566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU)); 301baa6e33SBarry Smith if (tao->bounded) PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU)); 3128017e9fSAlp Dener 3228017e9fSAlp Dener /* Project the initial point onto the feasible region */ 339566063dSJacob Faibussowitsch PetscCall(TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution)); 3428017e9fSAlp Dener 3528017e9fSAlp Dener /* Check convergence criteria */ 369566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient)); 379566063dSJacob Faibussowitsch PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type)); 389566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient)); 399566063dSJacob Faibussowitsch PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0)); 409566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm)); 4128017e9fSAlp Dener 42c0f10754SAlp Dener /* Test the initial point for convergence */ 439566063dSJacob Faibussowitsch PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W)); 449566063dSJacob Faibussowitsch PetscCall(VecNorm(bnk->W, NORM_2, &resnorm)); 453c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(bnk->f) && !PetscIsInfOrNanReal(resnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 469566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its)); 479566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0)); 48*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest ,tao->cnvP); 49c0f10754SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 50c0f10754SAlp Dener 51e031d6f5SAlp Dener /* Reset KSP stopping reason counters */ 52eb910715SAlp Dener bnk->ksp_atol = 0; 53eb910715SAlp Dener bnk->ksp_rtol = 0; 54eb910715SAlp Dener bnk->ksp_dtol = 0; 55eb910715SAlp Dener bnk->ksp_ctol = 0; 56eb910715SAlp Dener bnk->ksp_negc = 0; 57eb910715SAlp Dener bnk->ksp_iter = 0; 58eb910715SAlp Dener bnk->ksp_othr = 0; 59eb910715SAlp Dener 60e031d6f5SAlp Dener /* Reset accepted step type counters */ 61e031d6f5SAlp Dener bnk->tot_cg_its = 0; 62e031d6f5SAlp Dener bnk->newt = 0; 63e031d6f5SAlp Dener bnk->bfgs = 0; 64e031d6f5SAlp Dener bnk->sgrad = 0; 65e031d6f5SAlp Dener bnk->grad = 0; 66e031d6f5SAlp Dener 67fed79b8eSAlp Dener /* Initialize the Hessian perturbation */ 68fed79b8eSAlp Dener bnk->pert = bnk->sval; 69fed79b8eSAlp Dener 70937a31a1SAlp Dener /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */ 719566063dSJacob Faibussowitsch PetscCall(VecSet(tao->stepdirection, 0.0)); 72937a31a1SAlp Dener 73e031d6f5SAlp Dener /* Allocate the vectors needed for the BFGS approximation */ 749566063dSJacob Faibussowitsch PetscCall(KSPGetPC(tao->ksp, &pc)); 759566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs)); 769566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi)); 77b9ac7092SAlp Dener if (is_bfgs) { 78b9ac7092SAlp Dener bnk->bfgs_pre = pc; 799566063dSJacob Faibussowitsch PetscCall(PCLMVMGetMatLMVM(bnk->bfgs_pre, &bnk->M)); 809566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution, &n)); 819566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &N)); 829566063dSJacob Faibussowitsch PetscCall(MatSetSizes(bnk->M, n, n, N, N)); 839566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(bnk->M, tao->solution, bnk->unprojected_gradient)); 849566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(bnk->M, &sym_set, &is_symmetric)); 853c859ba3SBarry Smith PetscCheck(sym_set && is_symmetric,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric."); 861baa6e33SBarry Smith } else if (is_jacobi) PetscCall(PCJacobiSetUseAbs(pc,PETSC_TRUE)); 87e031d6f5SAlp Dener 88e031d6f5SAlp Dener /* Prepare the min/max vectors for safeguarding diagonal scales */ 899566063dSJacob Faibussowitsch PetscCall(VecSet(bnk->Diag_min, bnk->dmin)); 909566063dSJacob Faibussowitsch PetscCall(VecSet(bnk->Diag_max, bnk->dmax)); 91eb910715SAlp Dener 92eb910715SAlp Dener /* Initialize trust-region radius. The initialization is only performed 93eb910715SAlp Dener when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */ 94c0f10754SAlp Dener *needH = PETSC_TRUE; 959566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)tao->ksp,"KSPCGSetRadius_C",&kspTR)); 962e6e4ca1SStefano Zampini if (kspTR) { 9762675beeSAlp Dener switch (initType) { 98eb910715SAlp Dener case BNK_INIT_CONSTANT: 99eb910715SAlp Dener /* Use the initial radius specified */ 100c0f10754SAlp Dener tao->trust = tao->trust0; 101eb910715SAlp Dener break; 102eb910715SAlp Dener 103eb910715SAlp Dener case BNK_INIT_INTERPOLATION: 104c0f10754SAlp Dener /* Use interpolation based on the initial Hessian */ 105eb910715SAlp Dener max_radius = 0.0; 10608752603SAlp Dener tao->trust = tao->trust0; 107eb910715SAlp Dener for (j = 0; j < j_max; ++j) { 1080a4511e9SAlp Dener f_min = bnk->f; 109eb910715SAlp Dener sigma = 0.0; 110eb910715SAlp Dener 111c0f10754SAlp Dener if (*needH) { 11262602cfbSAlp Dener /* Compute the Hessian at the new step, and extract the inactive subsystem */ 1139566063dSJacob Faibussowitsch PetscCall((*bnk->computehessian)(tao)); 1149566063dSJacob Faibussowitsch PetscCall(TaoBNKEstimateActiveSet(tao, BNK_AS_NONE)); 1159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bnk->H_inactive)); 11689da521bSAlp Dener if (bnk->active_idx) { 1179566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive)); 11828017e9fSAlp Dener } else { 1199566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)tao->hessian)); 120c5e9d94cSAlp Dener bnk->H_inactive = tao->hessian; 12128017e9fSAlp Dener } 122c0f10754SAlp Dener *needH = PETSC_FALSE; 123eb910715SAlp Dener } 124eb910715SAlp Dener 125eb910715SAlp Dener for (i = 0; i < i_max; ++i) { 12662602cfbSAlp Dener /* Take a steepest descent step and snap it to bounds */ 1279566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, bnk->Xold)); 1289566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient)); 1299566063dSJacob Faibussowitsch PetscCall(TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution)); 13089da521bSAlp Dener /* Compute the step we actually accepted */ 1319566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, bnk->W)); 1329566063dSJacob Faibussowitsch PetscCall(VecAXPY(bnk->W, -1.0, bnk->Xold)); 13362602cfbSAlp Dener /* Compute the objective at the trial */ 1349566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao, tao->solution, &ftrial)); 1353c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(bnk->f),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 1369566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->Xold, tao->solution)); 137eb910715SAlp Dener if (PetscIsInfOrNanReal(ftrial)) { 138eb910715SAlp Dener tau = bnk->gamma1_i; 139eb910715SAlp Dener } else { 1400a4511e9SAlp Dener if (ftrial < f_min) { 1410a4511e9SAlp Dener f_min = ftrial; 142eb910715SAlp Dener sigma = -tao->trust / bnk->gnorm; 143eb910715SAlp Dener } 14408752603SAlp Dener 145770b7498SAlp Dener /* Compute the predicted and actual reduction */ 14689da521bSAlp Dener if (bnk->active_idx) { 1479566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive)); 1489566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work)); 1492ab2a32cSAlp Dener } else { 15008752603SAlp Dener bnk->X_inactive = bnk->W; 15108752603SAlp Dener bnk->inactive_work = bnk->Xwork; 1522ab2a32cSAlp Dener } 1539566063dSJacob Faibussowitsch PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work)); 1549566063dSJacob Faibussowitsch PetscCall(VecDot(bnk->X_inactive, bnk->inactive_work, &prered)); 15589da521bSAlp Dener if (bnk->active_idx) { 1569566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive)); 1579566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work)); 1582ab2a32cSAlp Dener } 159eb910715SAlp Dener prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm)); 160eb910715SAlp Dener actred = bnk->f - ftrial; 1613105154fSTodd Munson if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 162eb910715SAlp Dener kappa = 1.0; 1633105154fSTodd Munson } else { 164eb910715SAlp Dener kappa = actred / prered; 165eb910715SAlp Dener } 166eb910715SAlp Dener 167eb910715SAlp Dener tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred); 168eb910715SAlp Dener tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred); 169eb910715SAlp Dener tau_min = PetscMin(tau_1, tau_2); 170eb910715SAlp Dener tau_max = PetscMax(tau_1, tau_2); 171eb910715SAlp Dener 17218cfbf8eSSatish Balay if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu1_i) { 173eb910715SAlp Dener /* Great agreement */ 174eb910715SAlp Dener max_radius = PetscMax(max_radius, tao->trust); 175eb910715SAlp Dener 176eb910715SAlp Dener if (tau_max < 1.0) { 177eb910715SAlp Dener tau = bnk->gamma3_i; 1783105154fSTodd Munson } else if (tau_max > bnk->gamma4_i) { 179eb910715SAlp Dener tau = bnk->gamma4_i; 1803105154fSTodd Munson } else { 181eb910715SAlp Dener tau = tau_max; 182eb910715SAlp Dener } 18318cfbf8eSSatish Balay } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu2_i) { 184eb910715SAlp Dener /* Good agreement */ 185eb910715SAlp Dener max_radius = PetscMax(max_radius, tao->trust); 186eb910715SAlp Dener 187eb910715SAlp Dener if (tau_max < bnk->gamma2_i) { 188eb910715SAlp Dener tau = bnk->gamma2_i; 189eb910715SAlp Dener } else if (tau_max > bnk->gamma3_i) { 190eb910715SAlp Dener tau = bnk->gamma3_i; 191eb910715SAlp Dener } else { 192eb910715SAlp Dener tau = tau_max; 193eb910715SAlp Dener } 1948f8a4e06SAlp Dener } else { 195eb910715SAlp Dener /* Not good agreement */ 196eb910715SAlp Dener if (tau_min > 1.0) { 197eb910715SAlp Dener tau = bnk->gamma2_i; 198eb910715SAlp Dener } else if (tau_max < bnk->gamma1_i) { 199eb910715SAlp Dener tau = bnk->gamma1_i; 200eb910715SAlp Dener } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) { 201eb910715SAlp Dener tau = bnk->gamma1_i; 2023105154fSTodd Munson } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 203eb910715SAlp Dener tau = tau_1; 2043105154fSTodd Munson } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 205eb910715SAlp Dener tau = tau_2; 206eb910715SAlp Dener } else { 207eb910715SAlp Dener tau = tau_max; 208eb910715SAlp Dener } 209eb910715SAlp Dener } 210eb910715SAlp Dener } 211eb910715SAlp Dener tao->trust = tau * tao->trust; 212eb910715SAlp Dener } 213eb910715SAlp Dener 2140a4511e9SAlp Dener if (f_min < bnk->f) { 215937a31a1SAlp Dener /* We accidentally found a solution better than the initial, so accept it */ 2160a4511e9SAlp Dener bnk->f = f_min; 2179566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, bnk->Xold)); 2189566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution,sigma,tao->gradient)); 2199566063dSJacob Faibussowitsch PetscCall(TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution)); 2209566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, tao->stepdirection)); 2219566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold)); 2229566063dSJacob Faibussowitsch PetscCall(TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient)); 2239566063dSJacob Faibussowitsch PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type)); 2249566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient)); 2259566063dSJacob Faibussowitsch PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0)); 226937a31a1SAlp Dener /* Compute gradient at the new iterate and flip switch to compute the Hessian later */ 2279566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm)); 228c0f10754SAlp Dener *needH = PETSC_TRUE; 229937a31a1SAlp Dener /* Test the new step for convergence */ 2309566063dSJacob Faibussowitsch PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W)); 2319566063dSJacob Faibussowitsch PetscCall(VecNorm(bnk->W, NORM_2, &resnorm)); 2323c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(resnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 2339566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its)); 2349566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0)); 235*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest ,tao->cnvP); 236eb910715SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 237937a31a1SAlp Dener /* active BNCG recycling early because we have a stepdirection computed */ 2389566063dSJacob Faibussowitsch PetscCall(TaoSetRecycleHistory(bnk->bncg, PETSC_TRUE)); 239eb910715SAlp Dener } 240eb910715SAlp Dener } 241eb910715SAlp Dener tao->trust = PetscMax(tao->trust, max_radius); 242e031d6f5SAlp Dener 243e031d6f5SAlp Dener /* Ensure that the trust radius is within the limits */ 244e031d6f5SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 245e031d6f5SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 246eb910715SAlp Dener break; 247eb910715SAlp Dener 248eb910715SAlp Dener default: 249eb910715SAlp Dener /* Norm of the first direction will initialize radius */ 250eb910715SAlp Dener tao->trust = 0.0; 251eb910715SAlp Dener break; 252eb910715SAlp Dener } 253eb910715SAlp Dener } 254eb910715SAlp Dener PetscFunctionReturn(0); 255eb910715SAlp Dener } 256eb910715SAlp Dener 257df278d8fSAlp Dener /*------------------------------------------------------------*/ 258df278d8fSAlp Dener 259e0ed867bSAlp Dener /* Routine for computing the exact Hessian and preparing the preconditioner at the new iterate */ 26062675beeSAlp Dener 26162675beeSAlp Dener PetscErrorCode TaoBNKComputeHessian(Tao tao) 26262675beeSAlp Dener { 26362675beeSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 26462675beeSAlp Dener 26562675beeSAlp Dener PetscFunctionBegin; 26662675beeSAlp Dener /* Compute the Hessian */ 2679566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre)); 26862675beeSAlp Dener /* Add a correction to the BFGS preconditioner */ 2691baa6e33SBarry Smith if (bnk->M) PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient)); 270e0ed867bSAlp Dener /* Prepare the reduced sub-matrices for the inactive set */ 2719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bnk->Hpre_inactive)); 2729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bnk->H_inactive)); 273f5766c09SAlp Dener if (bnk->active_idx) { 2749566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive)); 275e0ed867bSAlp Dener if (tao->hessian == tao->hessian_pre) { 2769566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive)); 277e0ed867bSAlp Dener bnk->Hpre_inactive = bnk->H_inactive; 278e0ed867bSAlp Dener } else { 2799566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive)); 280e0ed867bSAlp Dener } 2811baa6e33SBarry Smith if (bnk->bfgs_pre) PetscCall(PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx)); 282e0ed867bSAlp Dener } else { 2839566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)tao->hessian)); 284c5e9d94cSAlp Dener bnk->H_inactive = tao->hessian; 285e0ed867bSAlp Dener if (tao->hessian == tao->hessian_pre) { 2869566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive)); 287e0ed867bSAlp Dener bnk->Hpre_inactive = bnk->H_inactive; 288e0ed867bSAlp Dener } else { 2899566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)tao->hessian_pre)); 290c5e9d94cSAlp Dener bnk->Hpre_inactive = tao->hessian_pre; 291e0ed867bSAlp Dener } 2921baa6e33SBarry Smith if (bnk->bfgs_pre) PetscCall(PCLMVMClearIS(bnk->bfgs_pre)); 293e0ed867bSAlp Dener } 29462675beeSAlp Dener PetscFunctionReturn(0); 29562675beeSAlp Dener } 29662675beeSAlp Dener 29762675beeSAlp Dener /*------------------------------------------------------------*/ 29862675beeSAlp Dener 2992f75a4aaSAlp Dener /* Routine for estimating the active set */ 3002f75a4aaSAlp Dener 30108752603SAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType) 3022f75a4aaSAlp Dener { 3032f75a4aaSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 304f4db9bf7SStefano Zampini PetscBool hessComputed, diagExists, hadactive; 3052f75a4aaSAlp Dener 3062f75a4aaSAlp Dener PetscFunctionBegin; 307f4db9bf7SStefano Zampini hadactive = bnk->active_idx ? PETSC_TRUE : PETSC_FALSE; 30808752603SAlp Dener switch (asType) { 3092f75a4aaSAlp Dener case BNK_AS_NONE: 3109566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->inactive_idx)); 3119566063dSJacob Faibussowitsch PetscCall(VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx)); 3129566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->active_idx)); 3139566063dSJacob Faibussowitsch PetscCall(ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx)); 3142f75a4aaSAlp Dener break; 3152f75a4aaSAlp Dener 3162f75a4aaSAlp Dener case BNK_AS_BERTSEKAS: 3172f75a4aaSAlp Dener /* Compute the trial step vector with which we will estimate the active set at the next iteration */ 318b9ac7092SAlp Dener if (bnk->M) { 3192f75a4aaSAlp Dener /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */ 3209566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W)); 3212f75a4aaSAlp Dener } else { 322fc5ca067SStefano Zampini hessComputed = diagExists = PETSC_FALSE; 323f5766c09SAlp Dener if (tao->hessian) { 3249566063dSJacob Faibussowitsch PetscCall(MatAssembled(tao->hessian, &hessComputed)); 325f5766c09SAlp Dener } 326fc5ca067SStefano Zampini if (hessComputed) { 3279566063dSJacob Faibussowitsch PetscCall(MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists)); 328fc5ca067SStefano Zampini } 329fc5ca067SStefano Zampini if (diagExists) { 3309b6ef848SAlp Dener /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */ 3319566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(tao->hessian, bnk->Xwork)); 3329566063dSJacob Faibussowitsch PetscCall(VecAbs(bnk->Xwork)); 3339566063dSJacob Faibussowitsch PetscCall(VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork)); 3349566063dSJacob Faibussowitsch PetscCall(VecReciprocal(bnk->Xwork)); 3359566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient)); 33661be54a6SAlp Dener } else { 337c4b75bccSAlp Dener /* If the Hessian or its diagonal does not exist, we will simply use gradient step */ 3389566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, bnk->W)); 33961be54a6SAlp Dener } 3402f75a4aaSAlp Dener } 3419566063dSJacob Faibussowitsch PetscCall(VecScale(bnk->W, -1.0)); 342d0609cedSBarry Smith PetscCall(TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol, 343d0609cedSBarry Smith &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx)); 344c4b75bccSAlp Dener break; 3452f75a4aaSAlp Dener 3462f75a4aaSAlp Dener default: 3472f75a4aaSAlp Dener break; 3482f75a4aaSAlp Dener } 349f4db9bf7SStefano Zampini bnk->resetksp = (PetscBool)(bnk->active_idx || hadactive); /* inactive Hessian size may have changed, need to reset operators */ 3502f75a4aaSAlp Dener PetscFunctionReturn(0); 3512f75a4aaSAlp Dener } 3522f75a4aaSAlp Dener 3532f75a4aaSAlp Dener /*------------------------------------------------------------*/ 3542f75a4aaSAlp Dener 3552f75a4aaSAlp Dener /* Routine for bounding the step direction */ 3562f75a4aaSAlp Dener 357a1318120SAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step) 3582f75a4aaSAlp Dener { 3592f75a4aaSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 3602f75a4aaSAlp Dener 3612f75a4aaSAlp Dener PetscFunctionBegin; 362a1318120SAlp Dener switch (asType) { 3632f75a4aaSAlp Dener case BNK_AS_NONE: 3649566063dSJacob Faibussowitsch PetscCall(VecISSet(step, bnk->active_idx, 0.0)); 3652f75a4aaSAlp Dener break; 3662f75a4aaSAlp Dener 3672f75a4aaSAlp Dener case BNK_AS_BERTSEKAS: 3689566063dSJacob Faibussowitsch PetscCall(TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step)); 3692f75a4aaSAlp Dener break; 3702f75a4aaSAlp Dener 3712f75a4aaSAlp Dener default: 3722f75a4aaSAlp Dener break; 3732f75a4aaSAlp Dener } 3742f75a4aaSAlp Dener PetscFunctionReturn(0); 3752f75a4aaSAlp Dener } 3762f75a4aaSAlp Dener 377e031d6f5SAlp Dener /*------------------------------------------------------------*/ 378e031d6f5SAlp Dener 379e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to 380e031d6f5SAlp Dener accelerate Newton convergence. 381e031d6f5SAlp Dener 382e031d6f5SAlp Dener In practice, this approach simply trades off Hessian evaluations 383e031d6f5SAlp Dener for more gradient evaluations. 384e031d6f5SAlp Dener */ 385e031d6f5SAlp Dener 386c0f10754SAlp Dener PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate) 387c0f10754SAlp Dener { 388c0f10754SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 389c0f10754SAlp Dener 390c0f10754SAlp Dener PetscFunctionBegin; 391c0f10754SAlp Dener *terminate = PETSC_FALSE; 392c0f10754SAlp Dener if (bnk->max_cg_its > 0) { 393c4b75bccSAlp Dener /* Copy the current function value (important vectors are already shared) */ 394c0f10754SAlp Dener bnk->bncg_ctx->f = bnk->f; 395c0f10754SAlp Dener /* Take some small finite number of BNCG iterations */ 3969566063dSJacob Faibussowitsch PetscCall(TaoSolve(bnk->bncg)); 397c0f10754SAlp Dener /* Add the number of gradient and function evaluations to the total */ 398c0f10754SAlp Dener tao->nfuncs += bnk->bncg->nfuncs; 399c0f10754SAlp Dener tao->nfuncgrads += bnk->bncg->nfuncgrads; 400c0f10754SAlp Dener tao->ngrads += bnk->bncg->ngrads; 401c0f10754SAlp Dener tao->nhess += bnk->bncg->nhess; 402e031d6f5SAlp Dener bnk->tot_cg_its += bnk->bncg->niter; 403c4b75bccSAlp Dener /* Extract the BNCG function value out and save it into BNK */ 404c0f10754SAlp Dener bnk->f = bnk->bncg_ctx->f; 405c0f10754SAlp Dener if (bnk->bncg->reason == TAO_CONVERGED_GATOL || bnk->bncg->reason == TAO_CONVERGED_GRTOL || bnk->bncg->reason == TAO_CONVERGED_GTTOL || bnk->bncg->reason == TAO_CONVERGED_MINF) { 406c0f10754SAlp Dener *terminate = PETSC_TRUE; 40761be54a6SAlp Dener } else { 4089566063dSJacob Faibussowitsch PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type)); 409c0f10754SAlp Dener } 410c0f10754SAlp Dener } 411c0f10754SAlp Dener PetscFunctionReturn(0); 412c0f10754SAlp Dener } 413c0f10754SAlp Dener 4142f75a4aaSAlp Dener /*------------------------------------------------------------*/ 4152f75a4aaSAlp Dener 416c0f10754SAlp Dener /* Routine for computing the Newton step. */ 417df278d8fSAlp Dener 4186b591159SAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type) 419eb910715SAlp Dener { 420eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 421eb910715SAlp Dener PetscInt bfgsUpdates = 0; 422eb910715SAlp Dener PetscInt kspits; 423bddd1ffdSAlp Dener PetscBool is_lmvm; 4242e6e4ca1SStefano Zampini PetscVoidFunction kspTR; 425eb910715SAlp Dener 426eb910715SAlp Dener PetscFunctionBegin; 42789da521bSAlp Dener /* If there are no inactive variables left, save some computation and return an adjusted zero step 42889da521bSAlp Dener that has (l-x) and (u-x) for lower and upper bounded variables. */ 42989da521bSAlp Dener if (!bnk->inactive_idx) { 4309566063dSJacob Faibussowitsch PetscCall(VecSet(tao->stepdirection, 0.0)); 4319566063dSJacob Faibussowitsch PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection)); 43289da521bSAlp Dener PetscFunctionReturn(0); 43389da521bSAlp Dener } 43489da521bSAlp Dener 43562675beeSAlp Dener /* Shift the reduced Hessian matrix */ 436e831869dSStefano Zampini if (shift && bnk->pert > 0) { 4379566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao->hessian, MATLMVM, &is_lmvm)); 438f7bf01afSAlp Dener if (is_lmvm) { 4399566063dSJacob Faibussowitsch PetscCall(MatShift(tao->hessian, bnk->pert)); 440f7bf01afSAlp Dener } else { 4419566063dSJacob Faibussowitsch PetscCall(MatShift(bnk->H_inactive, bnk->pert)); 44262675beeSAlp Dener if (bnk->H_inactive != bnk->Hpre_inactive) { 4439566063dSJacob Faibussowitsch PetscCall(MatShift(bnk->Hpre_inactive, bnk->pert)); 44462675beeSAlp Dener } 44562675beeSAlp Dener } 446f7bf01afSAlp Dener } 44762675beeSAlp Dener 448eb910715SAlp Dener /* Solve the Newton system of equations */ 449937a31a1SAlp Dener tao->ksp_its = 0; 4509566063dSJacob Faibussowitsch PetscCall(VecSet(tao->stepdirection, 0.0)); 451f4db9bf7SStefano Zampini if (bnk->resetksp) { 4529566063dSJacob Faibussowitsch PetscCall(KSPReset(tao->ksp)); 4539566063dSJacob Faibussowitsch PetscCall(KSPResetFromOptions(tao->ksp)); 454f4db9bf7SStefano Zampini bnk->resetksp = PETSC_FALSE; 455f4db9bf7SStefano Zampini } 4569566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive)); 4579566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, bnk->Gwork)); 45889da521bSAlp Dener if (bnk->active_idx) { 4599566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive)); 4609566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive)); 4615e9b73cbSAlp Dener } else { 4625e9b73cbSAlp Dener bnk->G_inactive = bnk->unprojected_gradient; 4635e9b73cbSAlp Dener bnk->X_inactive = tao->stepdirection; 46428017e9fSAlp Dener } 4659566063dSJacob Faibussowitsch PetscCall(KSPCGSetRadius(tao->ksp,tao->trust)); 4669566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive)); 4679566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&kspits)); 468eb910715SAlp Dener tao->ksp_its += kspits; 469eb910715SAlp Dener tao->ksp_tot_its += kspits; 470f4db9bf7SStefano Zampini PetscCall(PetscObjectQueryFunction((PetscObject)tao->ksp,"KSPCGGetNormD_C",&kspTR)); 471f4db9bf7SStefano Zampini if (kspTR) { 4729566063dSJacob Faibussowitsch PetscCall(KSPCGGetNormD(tao->ksp,&bnk->dnorm)); 473eb910715SAlp Dener 474eb910715SAlp Dener if (0.0 == tao->trust) { 475eb910715SAlp Dener /* Radius was uninitialized; use the norm of the direction */ 476080d2917SAlp Dener if (bnk->dnorm > 0.0) { 477080d2917SAlp Dener tao->trust = bnk->dnorm; 478eb910715SAlp Dener 479eb910715SAlp Dener /* Modify the radius if it is too large or small */ 480eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 481eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 482eb910715SAlp Dener } else { 483eb910715SAlp Dener /* The direction was bad; set radius to default value and re-solve 484eb910715SAlp Dener the trust-region subproblem to get a direction */ 485eb910715SAlp Dener tao->trust = tao->trust0; 486eb910715SAlp Dener 487eb910715SAlp Dener /* Modify the radius if it is too large or small */ 488eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 489eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 490eb910715SAlp Dener 4919566063dSJacob Faibussowitsch PetscCall(KSPCGSetRadius(tao->ksp,tao->trust)); 4929566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive)); 4939566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&kspits)); 494eb910715SAlp Dener tao->ksp_its += kspits; 495eb910715SAlp Dener tao->ksp_tot_its += kspits; 4969566063dSJacob Faibussowitsch PetscCall(KSPCGGetNormD(tao->ksp,&bnk->dnorm)); 497eb910715SAlp Dener 4983c859ba3SBarry Smith PetscCheck(bnk->dnorm != 0.0,PetscObjectComm((PetscObject)tao),PETSC_ERR_PLIB, "Initial direction zero"); 499eb910715SAlp Dener } 500eb910715SAlp Dener } 501eb910715SAlp Dener } 5025e9b73cbSAlp Dener /* Restore sub vectors back */ 50389da521bSAlp Dener if (bnk->active_idx) { 5049566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive)); 5059566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive)); 5065e9b73cbSAlp Dener } 507770b7498SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 5089566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0)); 5099566063dSJacob Faibussowitsch PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection)); 510770b7498SAlp Dener 511770b7498SAlp Dener /* Record convergence reasons */ 5129566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(tao->ksp, ksp_reason)); 513e465cd6fSAlp Dener if (KSP_CONVERGED_ATOL == *ksp_reason) { 514770b7498SAlp Dener ++bnk->ksp_atol; 515e465cd6fSAlp Dener } else if (KSP_CONVERGED_RTOL == *ksp_reason) { 516770b7498SAlp Dener ++bnk->ksp_rtol; 517e465cd6fSAlp Dener } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) { 518770b7498SAlp Dener ++bnk->ksp_ctol; 519e465cd6fSAlp Dener } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) { 520770b7498SAlp Dener ++bnk->ksp_negc; 521e465cd6fSAlp Dener } else if (KSP_DIVERGED_DTOL == *ksp_reason) { 522770b7498SAlp Dener ++bnk->ksp_dtol; 523e465cd6fSAlp Dener } else if (KSP_DIVERGED_ITS == *ksp_reason) { 524770b7498SAlp Dener ++bnk->ksp_iter; 525770b7498SAlp Dener } else { 526770b7498SAlp Dener ++bnk->ksp_othr; 527770b7498SAlp Dener } 528fed79b8eSAlp Dener 529fed79b8eSAlp Dener /* Make sure the BFGS preconditioner is healthy */ 530b9ac7092SAlp Dener if (bnk->M) { 5319566063dSJacob Faibussowitsch PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates)); 532b2d8c577SAlp Dener if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) { 533fed79b8eSAlp Dener /* Preconditioner is numerically indefinite; reset the approximation. */ 5349566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE)); 5359566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient)); 536eb910715SAlp Dener } 537fed79b8eSAlp Dener } 5386b591159SAlp Dener *step_type = BNK_NEWTON; 539e465cd6fSAlp Dener PetscFunctionReturn(0); 540e465cd6fSAlp Dener } 541eb910715SAlp Dener 54262675beeSAlp Dener /*------------------------------------------------------------*/ 54362675beeSAlp Dener 5445e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */ 5455e9b73cbSAlp Dener 5465e9b73cbSAlp Dener PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered) 5475e9b73cbSAlp Dener { 5485e9b73cbSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 5495e9b73cbSAlp Dener 5505e9b73cbSAlp Dener PetscFunctionBegin; 5515e9b73cbSAlp Dener /* Extract subvectors associated with the inactive set */ 55289da521bSAlp Dener if (bnk->active_idx) { 5539566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive)); 5549566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work)); 5559566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive)); 5565e9b73cbSAlp Dener } else { 5575e9b73cbSAlp Dener bnk->X_inactive = tao->stepdirection; 5585e9b73cbSAlp Dener bnk->inactive_work = bnk->Xwork; 5595e9b73cbSAlp Dener bnk->G_inactive = bnk->Gwork; 5605e9b73cbSAlp Dener } 5615e9b73cbSAlp Dener /* Recompute the predicted decrease based on the quadratic model */ 5629566063dSJacob Faibussowitsch PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work)); 5639566063dSJacob Faibussowitsch PetscCall(VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive)); 5649566063dSJacob Faibussowitsch PetscCall(VecDot(bnk->inactive_work, bnk->X_inactive, prered)); 5655e9b73cbSAlp Dener /* Restore the sub vectors */ 56689da521bSAlp Dener if (bnk->active_idx) { 5679566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive)); 5689566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work)); 5699566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive)); 5705e9b73cbSAlp Dener } 5715e9b73cbSAlp Dener PetscFunctionReturn(0); 5725e9b73cbSAlp Dener } 5735e9b73cbSAlp Dener 5745e9b73cbSAlp Dener /*------------------------------------------------------------*/ 5755e9b73cbSAlp Dener 57662675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction. 57762675beeSAlp Dener 57862675beeSAlp Dener The step direction falls back onto BFGS, scaled gradient and gradient steps 57962675beeSAlp Dener in the event that the Newton step fails the test. 58062675beeSAlp Dener */ 58162675beeSAlp Dener 582e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType) 583e465cd6fSAlp Dener { 584e465cd6fSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 585b2d8c577SAlp Dener PetscReal gdx, e_min; 586e465cd6fSAlp Dener PetscInt bfgsUpdates; 587e465cd6fSAlp Dener 588e465cd6fSAlp Dener PetscFunctionBegin; 5896b591159SAlp Dener switch (*stepType) { 5906b591159SAlp Dener case BNK_NEWTON: 5919566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx)); 592eb910715SAlp Dener if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 593eb910715SAlp Dener /* Newton step is not descent or direction produced Inf or NaN 594eb910715SAlp Dener Update the perturbation for next time */ 595eb910715SAlp Dener if (bnk->pert <= 0.0) { 5962e6e4ca1SStefano Zampini PetscBool is_gltr; 5972e6e4ca1SStefano Zampini 598eb910715SAlp Dener /* Initialize the perturbation */ 599eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 6009566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp),KSPGLTR,&is_gltr)); 6012e6e4ca1SStefano Zampini if (is_gltr) { 6029566063dSJacob Faibussowitsch PetscCall(KSPGLTRGetMinEig(tao->ksp,&e_min)); 603eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 604eb910715SAlp Dener } 605eb910715SAlp Dener } else { 606eb910715SAlp Dener /* Increase the perturbation */ 607eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 608eb910715SAlp Dener } 609eb910715SAlp Dener 6100ad3a497SAlp Dener if (!bnk->M) { 611eb910715SAlp Dener /* We don't have the bfgs matrix around and updated 612eb910715SAlp Dener Must use gradient direction in this case */ 6139566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, tao->stepdirection)); 614eb910715SAlp Dener *stepType = BNK_GRADIENT; 615eb910715SAlp Dener } else { 616eb910715SAlp Dener /* Attempt to use the BFGS direction */ 6179566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection)); 618eb910715SAlp Dener 6198d5ead36SAlp Dener /* Check for success (descent direction) 6208d5ead36SAlp Dener NOTE: Negative gdx here means not a descent direction because 6218d5ead36SAlp Dener the fall-back step is missing a negative sign. */ 6229566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx)); 6233105154fSTodd Munson if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) { 624eb910715SAlp Dener /* BFGS direction is not descent or direction produced not a number 625eb910715SAlp Dener We can assert bfgsUpdates > 1 in this case because 626eb910715SAlp Dener the first solve produces the scaled gradient direction, 627eb910715SAlp Dener which is guaranteed to be descent */ 628eb910715SAlp Dener 629eb910715SAlp Dener /* Use steepest descent direction (scaled) */ 6309566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE)); 6319566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient)); 6329566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection)); 633eb910715SAlp Dener 634eb910715SAlp Dener *stepType = BNK_SCALED_GRADIENT; 635eb910715SAlp Dener } else { 6369566063dSJacob Faibussowitsch PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates)); 637eb910715SAlp Dener if (1 == bfgsUpdates) { 638eb910715SAlp Dener /* The first BFGS direction is always the scaled gradient */ 639eb910715SAlp Dener *stepType = BNK_SCALED_GRADIENT; 640eb910715SAlp Dener } else { 641eb910715SAlp Dener *stepType = BNK_BFGS; 642eb910715SAlp Dener } 643eb910715SAlp Dener } 644eb910715SAlp Dener } 6458d5ead36SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 6469566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0)); 6479566063dSJacob Faibussowitsch PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection)); 648eb910715SAlp Dener } else { 649eb910715SAlp Dener /* Computed Newton step is descent */ 650eb910715SAlp Dener switch (ksp_reason) { 651eb910715SAlp Dener case KSP_DIVERGED_NANORINF: 652eb910715SAlp Dener case KSP_DIVERGED_BREAKDOWN: 653eb910715SAlp Dener case KSP_DIVERGED_INDEFINITE_MAT: 654eb910715SAlp Dener case KSP_DIVERGED_INDEFINITE_PC: 655eb910715SAlp Dener case KSP_CONVERGED_CG_NEG_CURVE: 656eb910715SAlp Dener /* Matrix or preconditioner is indefinite; increase perturbation */ 657eb910715SAlp Dener if (bnk->pert <= 0.0) { 6582e6e4ca1SStefano Zampini PetscBool is_gltr; 6592e6e4ca1SStefano Zampini 660eb910715SAlp Dener /* Initialize the perturbation */ 661eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 6629566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp),KSPGLTR,&is_gltr)); 6632e6e4ca1SStefano Zampini if (is_gltr) { 6649566063dSJacob Faibussowitsch PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min)); 665eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 666eb910715SAlp Dener } 667eb910715SAlp Dener } else { 668eb910715SAlp Dener /* Increase the perturbation */ 669eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 670eb910715SAlp Dener } 671eb910715SAlp Dener break; 672eb910715SAlp Dener 673eb910715SAlp Dener default: 674eb910715SAlp Dener /* Newton step computation is good; decrease perturbation */ 675eb910715SAlp Dener bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm); 676eb910715SAlp Dener if (bnk->pert < bnk->pmin) { 677eb910715SAlp Dener bnk->pert = 0.0; 678eb910715SAlp Dener } 679eb910715SAlp Dener break; 680eb910715SAlp Dener } 681fed79b8eSAlp Dener *stepType = BNK_NEWTON; 682eb910715SAlp Dener } 6836b591159SAlp Dener break; 6846b591159SAlp Dener 6856b591159SAlp Dener case BNK_BFGS: 6866b591159SAlp Dener /* Check for success (descent direction) */ 6879566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx)); 6886b591159SAlp Dener if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) { 6896b591159SAlp Dener /* Step is not descent or solve was not successful 6906b591159SAlp Dener Use steepest descent direction (scaled) */ 6919566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE)); 6929566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient)); 6939566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, tao->gradient, tao->stepdirection)); 6949566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection,-1.0)); 6959566063dSJacob Faibussowitsch PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection)); 6966b591159SAlp Dener *stepType = BNK_SCALED_GRADIENT; 6976b591159SAlp Dener } else { 6986b591159SAlp Dener *stepType = BNK_BFGS; 6996b591159SAlp Dener } 7006b591159SAlp Dener break; 7016b591159SAlp Dener 7026b591159SAlp Dener case BNK_SCALED_GRADIENT: 7036b591159SAlp Dener break; 7046b591159SAlp Dener 7056b591159SAlp Dener default: 7066b591159SAlp Dener break; 7076b591159SAlp Dener } 7086b591159SAlp Dener 709eb910715SAlp Dener PetscFunctionReturn(0); 710eb910715SAlp Dener } 711eb910715SAlp Dener 712df278d8fSAlp Dener /*------------------------------------------------------------*/ 713df278d8fSAlp Dener 714df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search. 715df278d8fSAlp Dener 716df278d8fSAlp Dener Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the 717df278d8fSAlp Dener Newton step does not produce a valid step length. 718df278d8fSAlp Dener */ 719df278d8fSAlp Dener 720937a31a1SAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason) 721c14b763aSAlp Dener { 722c14b763aSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 723c14b763aSAlp Dener TaoLineSearchConvergedReason ls_reason; 724b2d8c577SAlp Dener PetscReal e_min, gdx; 725c14b763aSAlp Dener PetscInt bfgsUpdates; 726c14b763aSAlp Dener 727c14b763aSAlp Dener PetscFunctionBegin; 728c14b763aSAlp Dener /* Perform the linesearch */ 7299566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason)); 7309566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 731c14b763aSAlp Dener 732b2d8c577SAlp Dener while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) { 733c14b763aSAlp Dener /* Linesearch failed, revert solution */ 734c14b763aSAlp Dener bnk->f = bnk->fold; 7359566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->Xold, tao->solution)); 7369566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient)); 737c14b763aSAlp Dener 738937a31a1SAlp Dener switch(*stepType) { 739c14b763aSAlp Dener case BNK_NEWTON: 7408d5ead36SAlp Dener /* Failed to obtain acceptable iterate with Newton step 741c14b763aSAlp Dener Update the perturbation for next time */ 742c14b763aSAlp Dener if (bnk->pert <= 0.0) { 7432e6e4ca1SStefano Zampini PetscBool is_gltr; 7442e6e4ca1SStefano Zampini 745c14b763aSAlp Dener /* Initialize the perturbation */ 746c14b763aSAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 7479566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp),KSPGLTR,&is_gltr)); 7482e6e4ca1SStefano Zampini if (is_gltr) { 7499566063dSJacob Faibussowitsch PetscCall(KSPGLTRGetMinEig(tao->ksp,&e_min)); 750c14b763aSAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 751c14b763aSAlp Dener } 752c14b763aSAlp Dener } else { 753c14b763aSAlp Dener /* Increase the perturbation */ 754c14b763aSAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 755c14b763aSAlp Dener } 756c14b763aSAlp Dener 7570ad3a497SAlp Dener if (!bnk->M) { 758c14b763aSAlp Dener /* We don't have the bfgs matrix around and being updated 759c14b763aSAlp Dener Must use gradient direction in this case */ 7609566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, tao->stepdirection)); 761937a31a1SAlp Dener *stepType = BNK_GRADIENT; 762c14b763aSAlp Dener } else { 763c14b763aSAlp Dener /* Attempt to use the BFGS direction */ 7649566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection)); 7658d5ead36SAlp Dener /* Check for success (descent direction) 7668d5ead36SAlp Dener NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */ 7679566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx)); 7683105154fSTodd Munson if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) { 769c14b763aSAlp Dener /* BFGS direction is not descent or direction produced not a number 770c14b763aSAlp Dener We can assert bfgsUpdates > 1 in this case 771c14b763aSAlp Dener Use steepest descent direction (scaled) */ 7729566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE)); 7739566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient)); 7749566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection)); 775c14b763aSAlp Dener 776c14b763aSAlp Dener bfgsUpdates = 1; 777937a31a1SAlp Dener *stepType = BNK_SCALED_GRADIENT; 778c14b763aSAlp Dener } else { 7799566063dSJacob Faibussowitsch PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates)); 780c14b763aSAlp Dener if (1 == bfgsUpdates) { 781c14b763aSAlp Dener /* The first BFGS direction is always the scaled gradient */ 782937a31a1SAlp Dener *stepType = BNK_SCALED_GRADIENT; 783c14b763aSAlp Dener } else { 784937a31a1SAlp Dener *stepType = BNK_BFGS; 785c14b763aSAlp Dener } 786c14b763aSAlp Dener } 787c14b763aSAlp Dener } 788c14b763aSAlp Dener break; 789c14b763aSAlp Dener 790c14b763aSAlp Dener case BNK_BFGS: 791c14b763aSAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS 792c14b763aSAlp Dener Failed to obtain acceptable iterate with BFGS step 793c14b763aSAlp Dener Attempt to use the scaled gradient direction */ 7949566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE)); 7959566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient)); 7969566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection)); 797c14b763aSAlp Dener 798c14b763aSAlp Dener bfgsUpdates = 1; 799937a31a1SAlp Dener *stepType = BNK_SCALED_GRADIENT; 800c14b763aSAlp Dener break; 801c14b763aSAlp Dener } 8028d5ead36SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 8039566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0)); 8049566063dSJacob Faibussowitsch PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection)); 805c14b763aSAlp Dener 8068d5ead36SAlp Dener /* Perform one last line search with the fall-back step */ 8079566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason)); 8089566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 809c14b763aSAlp Dener } 810c14b763aSAlp Dener *reason = ls_reason; 811c14b763aSAlp Dener PetscFunctionReturn(0); 812c14b763aSAlp Dener } 813c14b763aSAlp Dener 814df278d8fSAlp Dener /*------------------------------------------------------------*/ 815df278d8fSAlp Dener 816df278d8fSAlp Dener /* Routine for updating the trust radius. 817df278d8fSAlp Dener 818df278d8fSAlp Dener Function features three different update methods: 819df278d8fSAlp Dener 1) Line-search step length based 820df278d8fSAlp Dener 2) Predicted decrease on the CG quadratic model 821df278d8fSAlp Dener 3) Interpolation 822df278d8fSAlp Dener */ 823df278d8fSAlp Dener 82428017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept) 825080d2917SAlp Dener { 826080d2917SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 827080d2917SAlp Dener 828b1c2d0e3SAlp Dener PetscReal step, kappa; 829080d2917SAlp Dener PetscReal gdx, tau_1, tau_2, tau_min, tau_max; 830080d2917SAlp Dener 831080d2917SAlp Dener PetscFunctionBegin; 832080d2917SAlp Dener /* Update trust region radius */ 833080d2917SAlp Dener *accept = PETSC_FALSE; 83428017e9fSAlp Dener switch(updateType) { 835080d2917SAlp Dener case BNK_UPDATE_STEP: 836c14b763aSAlp Dener *accept = PETSC_TRUE; /* always accept here because line search succeeded */ 837080d2917SAlp Dener if (stepType == BNK_NEWTON) { 8389566063dSJacob Faibussowitsch PetscCall(TaoLineSearchGetStepLength(tao->linesearch, &step)); 839080d2917SAlp Dener if (step < bnk->nu1) { 840080d2917SAlp Dener /* Very bad step taken; reduce radius */ 841080d2917SAlp Dener tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 842080d2917SAlp Dener } else if (step < bnk->nu2) { 843080d2917SAlp Dener /* Reasonably bad step taken; reduce radius */ 844080d2917SAlp Dener tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust); 845080d2917SAlp Dener } else if (step < bnk->nu3) { 846080d2917SAlp Dener /* Reasonable step was taken; leave radius alone */ 847080d2917SAlp Dener if (bnk->omega3 < 1.0) { 848080d2917SAlp Dener tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust); 849080d2917SAlp Dener } else if (bnk->omega3 > 1.0) { 850080d2917SAlp Dener tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust); 851080d2917SAlp Dener } 852080d2917SAlp Dener } else if (step < bnk->nu4) { 853080d2917SAlp Dener /* Full step taken; increase the radius */ 854080d2917SAlp Dener tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust); 855080d2917SAlp Dener } else { 856080d2917SAlp Dener /* More than full step taken; increase the radius */ 857080d2917SAlp Dener tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust); 858080d2917SAlp Dener } 859080d2917SAlp Dener } else { 860080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 861080d2917SAlp Dener tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 862080d2917SAlp Dener } 863080d2917SAlp Dener break; 864080d2917SAlp Dener 865080d2917SAlp Dener case BNK_UPDATE_REDUCTION: 866080d2917SAlp Dener if (stepType == BNK_NEWTON) { 867e0ed867bSAlp Dener if ((prered < 0.0) || PetscIsInfOrNanReal(prered)) { 868fed79b8eSAlp Dener /* The predicted reduction has the wrong sign. This cannot 869fed79b8eSAlp Dener happen in infinite precision arithmetic. Step should 870fed79b8eSAlp Dener be rejected! */ 871080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 8723105154fSTodd Munson } else { 873b1c2d0e3SAlp Dener if (PetscIsInfOrNanReal(actred)) { 874080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 875080d2917SAlp Dener } else { 8763105154fSTodd Munson if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) { 877080d2917SAlp Dener kappa = 1.0; 8783105154fSTodd Munson } else { 879080d2917SAlp Dener kappa = actred / prered; 880080d2917SAlp Dener } 881fed79b8eSAlp Dener /* Accept or reject the step and update radius */ 882080d2917SAlp Dener if (kappa < bnk->eta1) { 883fed79b8eSAlp Dener /* Reject the step */ 884080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 8853105154fSTodd Munson } else { 886fed79b8eSAlp Dener /* Accept the step */ 887c133c014SAlp Dener *accept = PETSC_TRUE; 888c133c014SAlp Dener /* Update the trust region radius only if the computed step is at the trust radius boundary */ 8898d5ead36SAlp Dener if (bnk->dnorm == tao->trust) { 890080d2917SAlp Dener if (kappa < bnk->eta2) { 891080d2917SAlp Dener /* Marginal bad step */ 892c133c014SAlp Dener tao->trust = bnk->alpha2 * tao->trust; 8933105154fSTodd Munson } else if (kappa < bnk->eta3) { 894fed79b8eSAlp Dener /* Reasonable step */ 895fed79b8eSAlp Dener tao->trust = bnk->alpha3 * tao->trust; 8963105154fSTodd Munson } else if (kappa < bnk->eta4) { 897080d2917SAlp Dener /* Good step */ 898c133c014SAlp Dener tao->trust = bnk->alpha4 * tao->trust; 8993105154fSTodd Munson } else { 900080d2917SAlp Dener /* Very good step */ 901c133c014SAlp Dener tao->trust = bnk->alpha5 * tao->trust; 902080d2917SAlp Dener } 903c133c014SAlp Dener } 904080d2917SAlp Dener } 905080d2917SAlp Dener } 906080d2917SAlp Dener } 907080d2917SAlp Dener } else { 908080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 909080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust); 910080d2917SAlp Dener } 911080d2917SAlp Dener break; 912080d2917SAlp Dener 913080d2917SAlp Dener default: 914080d2917SAlp Dener if (stepType == BNK_NEWTON) { 915b1c2d0e3SAlp Dener if (prered < 0.0) { 916080d2917SAlp Dener /* The predicted reduction has the wrong sign. This cannot */ 917080d2917SAlp Dener /* happen in infinite precision arithmetic. Step should */ 918080d2917SAlp Dener /* be rejected! */ 919080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 920080d2917SAlp Dener } else { 921b1c2d0e3SAlp Dener if (PetscIsInfOrNanReal(actred)) { 922080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 923080d2917SAlp Dener } else { 924080d2917SAlp Dener if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 925080d2917SAlp Dener kappa = 1.0; 926080d2917SAlp Dener } else { 927080d2917SAlp Dener kappa = actred / prered; 928080d2917SAlp Dener } 929080d2917SAlp Dener 9309566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx)); 931080d2917SAlp Dener tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred); 932080d2917SAlp Dener tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred); 933080d2917SAlp Dener tau_min = PetscMin(tau_1, tau_2); 934080d2917SAlp Dener tau_max = PetscMax(tau_1, tau_2); 935080d2917SAlp Dener 936080d2917SAlp Dener if (kappa >= 1.0 - bnk->mu1) { 937080d2917SAlp Dener /* Great agreement */ 938080d2917SAlp Dener *accept = PETSC_TRUE; 939080d2917SAlp Dener if (tau_max < 1.0) { 940080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 941080d2917SAlp Dener } else if (tau_max > bnk->gamma4) { 942080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm); 943080d2917SAlp Dener } else { 944080d2917SAlp Dener tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 945080d2917SAlp Dener } 946080d2917SAlp Dener } else if (kappa >= 1.0 - bnk->mu2) { 947080d2917SAlp Dener /* Good agreement */ 948080d2917SAlp Dener *accept = PETSC_TRUE; 949080d2917SAlp Dener if (tau_max < bnk->gamma2) { 950080d2917SAlp Dener tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 951080d2917SAlp Dener } else if (tau_max > bnk->gamma3) { 952080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 953080d2917SAlp Dener } else if (tau_max < 1.0) { 954080d2917SAlp Dener tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 955080d2917SAlp Dener } else { 956080d2917SAlp Dener tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 957080d2917SAlp Dener } 958080d2917SAlp Dener } else { 959080d2917SAlp Dener /* Not good agreement */ 960080d2917SAlp Dener if (tau_min > 1.0) { 961080d2917SAlp Dener tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 962080d2917SAlp Dener } else if (tau_max < bnk->gamma1) { 963080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 964080d2917SAlp Dener } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) { 965080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 966080d2917SAlp Dener } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) { 967080d2917SAlp Dener tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm); 968080d2917SAlp Dener } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) { 969080d2917SAlp Dener tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm); 970080d2917SAlp Dener } else { 971080d2917SAlp Dener tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 972080d2917SAlp Dener } 973080d2917SAlp Dener } 974080d2917SAlp Dener } 975080d2917SAlp Dener } 976080d2917SAlp Dener } else { 977080d2917SAlp Dener /* Newton step was not good; reduce the radius */ 978080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust); 979080d2917SAlp Dener } 98028017e9fSAlp Dener break; 981080d2917SAlp Dener } 982c133c014SAlp Dener /* Make sure the radius does not violate min and max settings */ 983c133c014SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius); 984fed79b8eSAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius); 985080d2917SAlp Dener PetscFunctionReturn(0); 986080d2917SAlp Dener } 987080d2917SAlp Dener 988eb910715SAlp Dener /* ---------------------------------------------------------- */ 989df278d8fSAlp Dener 99062675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType) 99162675beeSAlp Dener { 99262675beeSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 99362675beeSAlp Dener 99462675beeSAlp Dener PetscFunctionBegin; 99562675beeSAlp Dener switch (stepType) { 99662675beeSAlp Dener case BNK_NEWTON: 99762675beeSAlp Dener ++bnk->newt; 99862675beeSAlp Dener break; 99962675beeSAlp Dener case BNK_BFGS: 100062675beeSAlp Dener ++bnk->bfgs; 100162675beeSAlp Dener break; 100262675beeSAlp Dener case BNK_SCALED_GRADIENT: 100362675beeSAlp Dener ++bnk->sgrad; 100462675beeSAlp Dener break; 100562675beeSAlp Dener case BNK_GRADIENT: 100662675beeSAlp Dener ++bnk->grad; 100762675beeSAlp Dener break; 100862675beeSAlp Dener default: 100962675beeSAlp Dener break; 101062675beeSAlp Dener } 101162675beeSAlp Dener PetscFunctionReturn(0); 101262675beeSAlp Dener } 101362675beeSAlp Dener 101462675beeSAlp Dener /* ---------------------------------------------------------- */ 101562675beeSAlp Dener 10169b6ef848SAlp Dener PetscErrorCode TaoSetUp_BNK(Tao tao) 1017eb910715SAlp Dener { 1018eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 1019e031d6f5SAlp Dener PetscInt i; 1020eb910715SAlp Dener 1021eb910715SAlp Dener PetscFunctionBegin; 1022c4b75bccSAlp Dener if (!tao->gradient) { 10239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 1024c4b75bccSAlp Dener } 1025c4b75bccSAlp Dener if (!tao->stepdirection) { 10269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->stepdirection)); 1027c4b75bccSAlp Dener } 1028c4b75bccSAlp Dener if (!bnk->W) { 10299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&bnk->W)); 1030c4b75bccSAlp Dener } 1031c4b75bccSAlp Dener if (!bnk->Xold) { 10329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&bnk->Xold)); 1033c4b75bccSAlp Dener } 1034c4b75bccSAlp Dener if (!bnk->Gold) { 10359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&bnk->Gold)); 1036c4b75bccSAlp Dener } 1037c4b75bccSAlp Dener if (!bnk->Xwork) { 10389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&bnk->Xwork)); 1039c4b75bccSAlp Dener } 1040c4b75bccSAlp Dener if (!bnk->Gwork) { 10419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&bnk->Gwork)); 1042c4b75bccSAlp Dener } 1043c4b75bccSAlp Dener if (!bnk->unprojected_gradient) { 10449566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&bnk->unprojected_gradient)); 1045c4b75bccSAlp Dener } 1046c4b75bccSAlp Dener if (!bnk->unprojected_gradient_old) { 10479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&bnk->unprojected_gradient_old)); 1048c4b75bccSAlp Dener } 1049c4b75bccSAlp Dener if (!bnk->Diag_min) { 10509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&bnk->Diag_min)); 1051c4b75bccSAlp Dener } 1052c4b75bccSAlp Dener if (!bnk->Diag_max) { 10539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&bnk->Diag_max)); 1054c4b75bccSAlp Dener } 1055e031d6f5SAlp Dener if (bnk->max_cg_its > 0) { 1056c4b75bccSAlp Dener /* Ensure that the important common vectors are shared between BNK and embedded BNCG */ 1057c4b75bccSAlp Dener bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data; 10589566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old))); 10599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old)); 106089da521bSAlp Dener bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old; 10619566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(bnk->unprojected_gradient))); 10629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient)); 1063c4b75bccSAlp Dener bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient; 10649566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(bnk->Gold))); 10659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->bncg_ctx->G_old)); 1066c4b75bccSAlp Dener bnk->bncg_ctx->G_old = bnk->Gold; 10679566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(tao->gradient))); 10689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->bncg->gradient)); 1069c4b75bccSAlp Dener bnk->bncg->gradient = tao->gradient; 10709566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(tao->stepdirection))); 10719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->bncg->stepdirection)); 1072c4b75bccSAlp Dener bnk->bncg->stepdirection = tao->stepdirection; 10739566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(bnk->bncg, tao->solution)); 1074c4b75bccSAlp Dener /* Copy over some settings from BNK into BNCG */ 10759566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its)); 10769566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol)); 10779566063dSJacob Faibussowitsch PetscCall(TaoSetFunctionLowerBound(bnk->bncg, tao->fmin)); 10789566063dSJacob Faibussowitsch PetscCall(TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP)); 10799566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(bnk->bncg, tao->ops->computeobjective, tao->user_objP)); 10809566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(bnk->bncg, NULL, tao->ops->computegradient, tao->user_gradP)); 10819566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(bnk->bncg, NULL, tao->ops->computeobjectiveandgradient, tao->user_objgradP)); 10829566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg))); 1083c4b75bccSAlp Dener for (i=0; i<tao->numbermonitors; ++i) { 10849566063dSJacob Faibussowitsch PetscCall(TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i])); 10859566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(tao->monitorcontext[i]))); 1086e031d6f5SAlp Dener } 1087e031d6f5SAlp Dener } 108883c8fe1dSLisandro Dalcin bnk->X_inactive = NULL; 108983c8fe1dSLisandro Dalcin bnk->G_inactive = NULL; 109083c8fe1dSLisandro Dalcin bnk->inactive_work = NULL; 109183c8fe1dSLisandro Dalcin bnk->active_work = NULL; 109283c8fe1dSLisandro Dalcin bnk->inactive_idx = NULL; 109383c8fe1dSLisandro Dalcin bnk->active_idx = NULL; 109483c8fe1dSLisandro Dalcin bnk->active_lower = NULL; 109583c8fe1dSLisandro Dalcin bnk->active_upper = NULL; 109683c8fe1dSLisandro Dalcin bnk->active_fixed = NULL; 109783c8fe1dSLisandro Dalcin bnk->M = NULL; 109883c8fe1dSLisandro Dalcin bnk->H_inactive = NULL; 109983c8fe1dSLisandro Dalcin bnk->Hpre_inactive = NULL; 1100eb910715SAlp Dener PetscFunctionReturn(0); 1101eb910715SAlp Dener } 1102eb910715SAlp Dener 1103eb910715SAlp Dener /*------------------------------------------------------------*/ 1104df278d8fSAlp Dener 1105e0ed867bSAlp Dener PetscErrorCode TaoDestroy_BNK(Tao tao) 1106eb910715SAlp Dener { 1107eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 1108eb910715SAlp Dener 1109eb910715SAlp Dener PetscFunctionBegin; 11109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->W)); 11119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->Xold)); 11129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->Gold)); 11139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->Xwork)); 11149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->Gwork)); 11159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->unprojected_gradient)); 11169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->unprojected_gradient_old)); 11179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->Diag_min)); 11189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->Diag_max)); 11199566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->active_lower)); 11209566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->active_upper)); 11219566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->active_fixed)); 11229566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->active_idx)); 11239566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->inactive_idx)); 11249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bnk->Hpre_inactive)); 11259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bnk->H_inactive)); 11269566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&bnk->bncg)); 1127a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 11289566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 1129eb910715SAlp Dener PetscFunctionReturn(0); 1130eb910715SAlp Dener } 1131eb910715SAlp Dener 1132eb910715SAlp Dener /*------------------------------------------------------------*/ 1133df278d8fSAlp Dener 1134*dbbe0bcdSBarry Smith PetscErrorCode TaoSetFromOptions_BNK(Tao tao,PetscOptionItems *PetscOptionsObject) 1135eb910715SAlp Dener { 1136eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 1137eb910715SAlp Dener 1138eb910715SAlp Dener PetscFunctionBegin; 1139d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Newton-Krylov method for bound constrained optimization"); 11409566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, NULL)); 11419566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, NULL)); 11429566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, NULL)); 11439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval,NULL)); 11449566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin,NULL)); 11459566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax,NULL)); 11469566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac,NULL)); 11479566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin,NULL)); 11489566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax,NULL)); 11499566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac,NULL)); 11509566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac,NULL)); 11519566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac,NULL)); 11529566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac,NULL)); 11539566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1,NULL)); 11549566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2,NULL)); 11559566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3,NULL)); 11569566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4,NULL)); 11579566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_alpha1", "(developer) radius reduction factor for rejected step (-tao_bnk_update_type reduction)", "", bnk->alpha1, &bnk->alpha1,NULL)); 11589566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_alpha2", "(developer) radius reduction factor for marginally accepted bad step (-tao_bnk_update_type reduction)", "", bnk->alpha2, &bnk->alpha2,NULL)); 11599566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_alpha3", "(developer) radius increase factor for reasonable accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha3, &bnk->alpha3,NULL)); 11609566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_alpha4", "(developer) radius increase factor for good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha4, &bnk->alpha4,NULL)); 11619566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_alpha5", "(developer) radius increase factor for very good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha5, &bnk->alpha5,NULL)); 11629566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_nu1", "(developer) threshold for small line-search step length (-tao_bnk_update_type step)", "", bnk->nu1, &bnk->nu1,NULL)); 11639566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_nu2", "(developer) threshold for reasonable line-search step length (-tao_bnk_update_type step)", "", bnk->nu2, &bnk->nu2,NULL)); 11649566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_nu3", "(developer) threshold for large line-search step length (-tao_bnk_update_type step)", "", bnk->nu3, &bnk->nu3,NULL)); 11659566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_nu4", "(developer) threshold for very large line-search step length (-tao_bnk_update_type step)", "", bnk->nu4, &bnk->nu4,NULL)); 11669566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_omega1", "(developer) radius reduction factor for very small line-search step length (-tao_bnk_update_type step)", "", bnk->omega1, &bnk->omega1,NULL)); 11679566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_omega2", "(developer) radius reduction factor for small line-search step length (-tao_bnk_update_type step)", "", bnk->omega2, &bnk->omega2,NULL)); 11689566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_omega3", "(developer) radius factor for decent line-search step length (-tao_bnk_update_type step)", "", bnk->omega3, &bnk->omega3,NULL)); 11699566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_omega4", "(developer) radius increase factor for large line-search step length (-tao_bnk_update_type step)", "", bnk->omega4, &bnk->omega4,NULL)); 11709566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_omega5", "(developer) radius increase factor for very large line-search step length (-tao_bnk_update_type step)", "", bnk->omega5, &bnk->omega5,NULL)); 11719566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_mu1_i", "(developer) threshold for accepting very good step (-tao_bnk_init_type interpolation)", "", bnk->mu1_i, &bnk->mu1_i,NULL)); 11729566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_mu2_i", "(developer) threshold for accepting good step (-tao_bnk_init_type interpolation)", "", bnk->mu2_i, &bnk->mu2_i,NULL)); 11739566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_gamma1_i", "(developer) radius reduction factor for rejected very bad step (-tao_bnk_init_type interpolation)", "", bnk->gamma1_i, &bnk->gamma1_i,NULL)); 11749566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_gamma2_i", "(developer) radius reduction factor for rejected bad step (-tao_bnk_init_type interpolation)", "", bnk->gamma2_i, &bnk->gamma2_i,NULL)); 11759566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_gamma3_i", "(developer) radius increase factor for accepted good step (-tao_bnk_init_type interpolation)", "", bnk->gamma3_i, &bnk->gamma3_i,NULL)); 11769566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_gamma4_i", "(developer) radius increase factor for accepted very good step (-tao_bnk_init_type interpolation)", "", bnk->gamma4_i, &bnk->gamma4_i,NULL)); 11779566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_theta_i", "(developer) trust region interpolation factor (-tao_bnk_init_type interpolation)", "", bnk->theta_i, &bnk->theta_i,NULL)); 11789566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_mu1", "(developer) threshold for accepting very good step (-tao_bnk_update_type interpolation)", "", bnk->mu1, &bnk->mu1,NULL)); 11799566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2,NULL)); 11809566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_gamma1", "(developer) radius reduction factor for rejected very bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma1, &bnk->gamma1,NULL)); 11819566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_gamma2", "(developer) radius reduction factor for rejected bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma2, &bnk->gamma2,NULL)); 11829566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_gamma3", "(developer) radius increase factor for accepted good step (-tao_bnk_update_type interpolation)", "", bnk->gamma3, &bnk->gamma3,NULL)); 11839566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_gamma4", "(developer) radius increase factor for accepted very good step (-tao_bnk_update_type interpolation)", "", bnk->gamma4, &bnk->gamma4,NULL)); 11849566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta,NULL)); 11859566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL)); 11869566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL)); 11879566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL)); 11889566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL)); 11899566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL)); 11909566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_bnk_max_cg_its", "number of BNCG iterations to take for each Newton step", "", bnk->max_cg_its, &bnk->max_cg_its,NULL)); 1191d0609cedSBarry Smith PetscOptionsHeadEnd(); 11928ebe3e4eSStefano Zampini 11939566063dSJacob Faibussowitsch PetscCall(TaoSetOptionsPrefix(bnk->bncg,((PetscObject)(tao))->prefix)); 11949566063dSJacob Faibussowitsch PetscCall(TaoAppendOptionsPrefix(bnk->bncg,"tao_bnk_cg_")); 11959566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(bnk->bncg)); 11968ebe3e4eSStefano Zampini 11979566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp,((PetscObject)(tao))->prefix)); 11989566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(tao->ksp,"tao_bnk_")); 11999566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp)); 1200eb910715SAlp Dener PetscFunctionReturn(0); 1201eb910715SAlp Dener } 1202eb910715SAlp Dener 1203eb910715SAlp Dener /*------------------------------------------------------------*/ 1204df278d8fSAlp Dener 1205e0ed867bSAlp Dener PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer) 1206eb910715SAlp Dener { 1207eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 1208eb910715SAlp Dener PetscInt nrejects; 1209eb910715SAlp Dener PetscBool isascii; 1210eb910715SAlp Dener 1211eb910715SAlp Dener PetscFunctionBegin; 12129566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); 1213eb910715SAlp Dener if (isascii) { 12149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1215b9ac7092SAlp Dener if (bnk->M) { 12169566063dSJacob Faibussowitsch PetscCall(MatLMVMGetRejectCount(bnk->M,&nrejects)); 121763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %" PetscInt_FMT "\n",nrejects)); 1218eb910715SAlp Dener } 121963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "CG steps: %" PetscInt_FMT "\n", bnk->tot_cg_its)); 122063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Newton steps: %" PetscInt_FMT "\n", bnk->newt)); 1221b9ac7092SAlp Dener if (bnk->M) { 122263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS steps: %" PetscInt_FMT "\n", bnk->bfgs)); 1223b9ac7092SAlp Dener } 122463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %" PetscInt_FMT "\n", bnk->sgrad)); 122563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", bnk->grad)); 12269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n")); 122763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " atol: %" PetscInt_FMT "\n", bnk->ksp_atol)); 122863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " rtol: %" PetscInt_FMT "\n", bnk->ksp_rtol)); 122963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ctol: %" PetscInt_FMT "\n", bnk->ksp_ctol)); 123063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " negc: %" PetscInt_FMT "\n", bnk->ksp_negc)); 123163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dtol: %" PetscInt_FMT "\n", bnk->ksp_dtol)); 123263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " iter: %" PetscInt_FMT "\n", bnk->ksp_iter)); 123363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " othr: %" PetscInt_FMT "\n", bnk->ksp_othr)); 12349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1235eb910715SAlp Dener } 1236eb910715SAlp Dener PetscFunctionReturn(0); 1237eb910715SAlp Dener } 1238eb910715SAlp Dener 1239eb910715SAlp Dener /* ---------------------------------------------------------- */ 1240df278d8fSAlp Dener 1241eb910715SAlp Dener /*MC 1242eb910715SAlp Dener TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms. 124366ed3702SAlp Dener At each iteration, the BNK methods solve the symmetric 1244eb910715SAlp Dener system of equations to obtain the step diretion dk: 1245eb910715SAlp Dener Hk dk = -gk 12462b97c8d8SAlp Dener for free variables only. The step can be globalized either through 12472b97c8d8SAlp Dener trust-region methods, or a line search, or a heuristic mixture of both. 1248eb910715SAlp Dener 1249eb910715SAlp Dener Options Database Keys: 12509fa2b5dcSStefano Zampini + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop 12519fa2b5dcSStefano Zampini . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation") 12529fa2b5dcSStefano Zampini . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation") 12539fa2b5dcSStefano Zampini . -tao_bnk_as_type - active-set estimation method ("none", "bertsekas") 12549fa2b5dcSStefano Zampini . -tao_bnk_as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas) 12559fa2b5dcSStefano Zampini . -tao_bnk_as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas) 12569fa2b5dcSStefano Zampini . -tao_bnk_sval - (developer) Hessian perturbation starting value 12579fa2b5dcSStefano Zampini . -tao_bnk_imin - (developer) minimum initial Hessian perturbation 12589fa2b5dcSStefano Zampini . -tao_bnk_imax - (developer) maximum initial Hessian perturbation 12599fa2b5dcSStefano Zampini . -tao_bnk_pmin - (developer) minimum Hessian perturbation 12609fa2b5dcSStefano Zampini . -tao_bnk_pmax - (developer) aximum Hessian perturbation 12619fa2b5dcSStefano Zampini . -tao_bnk_pgfac - (developer) Hessian perturbation growth factor 12629fa2b5dcSStefano Zampini . -tao_bnk_psfac - (developer) Hessian perturbation shrink factor 12639fa2b5dcSStefano Zampini . -tao_bnk_imfac - (developer) initial merit factor for Hessian perturbation 12649fa2b5dcSStefano Zampini . -tao_bnk_pmgfac - (developer) merit growth factor for Hessian perturbation 12659fa2b5dcSStefano Zampini . -tao_bnk_pmsfac - (developer) merit shrink factor for Hessian perturbation 12669fa2b5dcSStefano Zampini . -tao_bnk_eta1 - (developer) threshold for rejecting step (-update_type reduction) 12679fa2b5dcSStefano Zampini . -tao_bnk_eta2 - (developer) threshold for accepting marginal step (-update_type reduction) 12689fa2b5dcSStefano Zampini . -tao_bnk_eta3 - (developer) threshold for accepting reasonable step (-update_type reduction) 12699fa2b5dcSStefano Zampini . -tao_bnk_eta4 - (developer) threshold for accepting good step (-update_type reduction) 12709fa2b5dcSStefano Zampini . -tao_bnk_alpha1 - (developer) radius reduction factor for rejected step (-update_type reduction) 12719fa2b5dcSStefano Zampini . -tao_bnk_alpha2 - (developer) radius reduction factor for marginally accepted bad step (-update_type reduction) 12729fa2b5dcSStefano Zampini . -tao_bnk_alpha3 - (developer) radius increase factor for reasonable accepted step (-update_type reduction) 12739fa2b5dcSStefano Zampini . -tao_bnk_alpha4 - (developer) radius increase factor for good accepted step (-update_type reduction) 12749fa2b5dcSStefano Zampini . -tao_bnk_alpha5 - (developer) radius increase factor for very good accepted step (-update_type reduction) 12759fa2b5dcSStefano Zampini . -tao_bnk_epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-update_type reduction) 12769fa2b5dcSStefano Zampini . -tao_bnk_mu1 - (developer) threshold for accepting very good step (-update_type interpolation) 12779fa2b5dcSStefano Zampini . -tao_bnk_mu2 - (developer) threshold for accepting good step (-update_type interpolation) 12789fa2b5dcSStefano Zampini . -tao_bnk_gamma1 - (developer) radius reduction factor for rejected very bad step (-update_type interpolation) 12799fa2b5dcSStefano Zampini . -tao_bnk_gamma2 - (developer) radius reduction factor for rejected bad step (-update_type interpolation) 12809fa2b5dcSStefano Zampini . -tao_bnk_gamma3 - (developer) radius increase factor for accepted good step (-update_type interpolation) 12819fa2b5dcSStefano Zampini . -tao_bnk_gamma4 - (developer) radius increase factor for accepted very good step (-update_type interpolation) 12829fa2b5dcSStefano Zampini . -tao_bnk_theta - (developer) trust region interpolation factor (-update_type interpolation) 12839fa2b5dcSStefano Zampini . -tao_bnk_nu1 - (developer) threshold for small line-search step length (-update_type step) 12849fa2b5dcSStefano Zampini . -tao_bnk_nu2 - (developer) threshold for reasonable line-search step length (-update_type step) 12859fa2b5dcSStefano Zampini . -tao_bnk_nu3 - (developer) threshold for large line-search step length (-update_type step) 12869fa2b5dcSStefano Zampini . -tao_bnk_nu4 - (developer) threshold for very large line-search step length (-update_type step) 12879fa2b5dcSStefano Zampini . -tao_bnk_omega1 - (developer) radius reduction factor for very small line-search step length (-update_type step) 12889fa2b5dcSStefano Zampini . -tao_bnk_omega2 - (developer) radius reduction factor for small line-search step length (-update_type step) 12899fa2b5dcSStefano Zampini . -tao_bnk_omega3 - (developer) radius factor for decent line-search step length (-update_type step) 12909fa2b5dcSStefano Zampini . -tao_bnk_omega4 - (developer) radius increase factor for large line-search step length (-update_type step) 12919fa2b5dcSStefano Zampini . -tao_bnk_omega5 - (developer) radius increase factor for very large line-search step length (-update_type step) 12929fa2b5dcSStefano Zampini . -tao_bnk_mu1_i - (developer) threshold for accepting very good step (-init_type interpolation) 12939fa2b5dcSStefano Zampini . -tao_bnk_mu2_i - (developer) threshold for accepting good step (-init_type interpolation) 12949fa2b5dcSStefano Zampini . -tao_bnk_gamma1_i - (developer) radius reduction factor for rejected very bad step (-init_type interpolation) 12959fa2b5dcSStefano Zampini . -tao_bnk_gamma2_i - (developer) radius reduction factor for rejected bad step (-init_type interpolation) 12969fa2b5dcSStefano Zampini . -tao_bnk_gamma3_i - (developer) radius increase factor for accepted good step (-init_type interpolation) 12979fa2b5dcSStefano Zampini . -tao_bnk_gamma4_i - (developer) radius increase factor for accepted very good step (-init_type interpolation) 12989fa2b5dcSStefano Zampini - -tao_bnk_theta_i - (developer) trust region interpolation factor (-init_type interpolation) 1299eb910715SAlp Dener 1300eb910715SAlp Dener Level: beginner 1301eb910715SAlp Dener M*/ 1302eb910715SAlp Dener 1303eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao) 1304eb910715SAlp Dener { 1305eb910715SAlp Dener TAO_BNK *bnk; 1306b9ac7092SAlp Dener PC pc; 1307eb910715SAlp Dener 1308eb910715SAlp Dener PetscFunctionBegin; 13099566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&bnk)); 1310eb910715SAlp Dener 1311eb910715SAlp Dener tao->ops->setup = TaoSetUp_BNK; 1312eb910715SAlp Dener tao->ops->view = TaoView_BNK; 1313eb910715SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BNK; 1314eb910715SAlp Dener tao->ops->destroy = TaoDestroy_BNK; 1315eb910715SAlp Dener 1316eb910715SAlp Dener /* Override default settings (unless already changed) */ 1317eb910715SAlp Dener if (!tao->max_it_changed) tao->max_it = 50; 1318eb910715SAlp Dener if (!tao->trust0_changed) tao->trust0 = 100.0; 1319eb910715SAlp Dener 1320eb910715SAlp Dener tao->data = (void*)bnk; 1321eb910715SAlp Dener 132266ed3702SAlp Dener /* Hessian shifting parameters */ 1323e0ed867bSAlp Dener bnk->computehessian = TaoBNKComputeHessian; 1324e0ed867bSAlp Dener bnk->computestep = TaoBNKComputeStep; 1325e0ed867bSAlp Dener 1326eb910715SAlp Dener bnk->sval = 0.0; 1327eb910715SAlp Dener bnk->imin = 1.0e-4; 1328eb910715SAlp Dener bnk->imax = 1.0e+2; 1329eb910715SAlp Dener bnk->imfac = 1.0e-1; 1330eb910715SAlp Dener 1331eb910715SAlp Dener bnk->pmin = 1.0e-12; 1332eb910715SAlp Dener bnk->pmax = 1.0e+2; 1333eb910715SAlp Dener bnk->pgfac = 1.0e+1; 1334eb910715SAlp Dener bnk->psfac = 4.0e-1; 1335eb910715SAlp Dener bnk->pmgfac = 1.0e-1; 1336eb910715SAlp Dener bnk->pmsfac = 1.0e-1; 1337eb910715SAlp Dener 1338eb910715SAlp Dener /* Default values for trust-region radius update based on steplength */ 1339eb910715SAlp Dener bnk->nu1 = 0.25; 1340eb910715SAlp Dener bnk->nu2 = 0.50; 1341eb910715SAlp Dener bnk->nu3 = 1.00; 1342eb910715SAlp Dener bnk->nu4 = 1.25; 1343eb910715SAlp Dener 1344eb910715SAlp Dener bnk->omega1 = 0.25; 1345eb910715SAlp Dener bnk->omega2 = 0.50; 1346eb910715SAlp Dener bnk->omega3 = 1.00; 1347eb910715SAlp Dener bnk->omega4 = 2.00; 1348eb910715SAlp Dener bnk->omega5 = 4.00; 1349eb910715SAlp Dener 1350eb910715SAlp Dener /* Default values for trust-region radius update based on reduction */ 1351eb910715SAlp Dener bnk->eta1 = 1.0e-4; 1352eb910715SAlp Dener bnk->eta2 = 0.25; 1353eb910715SAlp Dener bnk->eta3 = 0.50; 1354eb910715SAlp Dener bnk->eta4 = 0.90; 1355eb910715SAlp Dener 1356eb910715SAlp Dener bnk->alpha1 = 0.25; 1357eb910715SAlp Dener bnk->alpha2 = 0.50; 1358eb910715SAlp Dener bnk->alpha3 = 1.00; 1359eb910715SAlp Dener bnk->alpha4 = 2.00; 1360eb910715SAlp Dener bnk->alpha5 = 4.00; 1361eb910715SAlp Dener 1362eb910715SAlp Dener /* Default values for trust-region radius update based on interpolation */ 1363eb910715SAlp Dener bnk->mu1 = 0.10; 1364eb910715SAlp Dener bnk->mu2 = 0.50; 1365eb910715SAlp Dener 1366eb910715SAlp Dener bnk->gamma1 = 0.25; 1367eb910715SAlp Dener bnk->gamma2 = 0.50; 1368eb910715SAlp Dener bnk->gamma3 = 2.00; 1369eb910715SAlp Dener bnk->gamma4 = 4.00; 1370eb910715SAlp Dener 1371eb910715SAlp Dener bnk->theta = 0.05; 1372eb910715SAlp Dener 1373eb910715SAlp Dener /* Default values for trust region initialization based on interpolation */ 1374eb910715SAlp Dener bnk->mu1_i = 0.35; 1375eb910715SAlp Dener bnk->mu2_i = 0.50; 1376eb910715SAlp Dener 1377eb910715SAlp Dener bnk->gamma1_i = 0.0625; 1378eb910715SAlp Dener bnk->gamma2_i = 0.5; 1379eb910715SAlp Dener bnk->gamma3_i = 2.0; 1380eb910715SAlp Dener bnk->gamma4_i = 5.0; 1381eb910715SAlp Dener 1382eb910715SAlp Dener bnk->theta_i = 0.25; 1383eb910715SAlp Dener 1384eb910715SAlp Dener /* Remaining parameters */ 1385c0f10754SAlp Dener bnk->max_cg_its = 0; 1386eb910715SAlp Dener bnk->min_radius = 1.0e-10; 1387eb910715SAlp Dener bnk->max_radius = 1.0e10; 1388770b7498SAlp Dener bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); 13890a4511e9SAlp Dener bnk->as_tol = 1.0e-3; 13900a4511e9SAlp Dener bnk->as_step = 1.0e-3; 139162675beeSAlp Dener bnk->dmin = 1.0e-6; 139262675beeSAlp Dener bnk->dmax = 1.0e6; 1393eb910715SAlp Dener 139483c8fe1dSLisandro Dalcin bnk->M = NULL; 139583c8fe1dSLisandro Dalcin bnk->bfgs_pre = NULL; 1396eb910715SAlp Dener bnk->init_type = BNK_INIT_INTERPOLATION; 13977b1c7716SAlp Dener bnk->update_type = BNK_UPDATE_REDUCTION; 13982f75a4aaSAlp Dener bnk->as_type = BNK_AS_BERTSEKAS; 1399eb910715SAlp Dener 1400e031d6f5SAlp Dener /* Create the embedded BNCG solver */ 14019566063dSJacob Faibussowitsch PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&bnk->bncg)); 14029566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)bnk->bncg,(PetscObject)tao,1)); 14039566063dSJacob Faibussowitsch PetscCall(TaoSetType(bnk->bncg,TAOBNCG)); 1404e031d6f5SAlp Dener 1405c0f10754SAlp Dener /* Create the line search */ 14069566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch)); 14079566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch,(PetscObject)tao,1)); 1408f4db9bf7SStefano Zampini PetscCall(TaoLineSearchSetType(tao->linesearch,TAOLINESEARCHMT)); 14099566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch,tao)); 1410eb910715SAlp Dener 1411eb910715SAlp Dener /* Set linear solver to default for symmetric matrices */ 14129566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm,&tao->ksp)); 14139566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp,(PetscObject)tao,1)); 14149566063dSJacob Faibussowitsch PetscCall(KSPSetType(tao->ksp,KSPSTCG)); 14159566063dSJacob Faibussowitsch PetscCall(KSPGetPC(tao->ksp,&pc)); 14169566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCLMVM)); 1417eb910715SAlp Dener PetscFunctionReturn(0); 1418eb910715SAlp Dener } 1419