155119615STodd Munson #include <../src/tao/unconstrained/impls/ntl/ntlimpl.h> 2a7e14dcfSSatish Balay 3aaa7dc30SBarry Smith #include <petscksp.h> 4a7e14dcfSSatish Balay 5a7e14dcfSSatish Balay #define NTL_INIT_CONSTANT 0 6a7e14dcfSSatish Balay #define NTL_INIT_DIRECTION 1 7a7e14dcfSSatish Balay #define NTL_INIT_INTERPOLATION 2 8a7e14dcfSSatish Balay #define NTL_INIT_TYPES 3 9a7e14dcfSSatish Balay 10a7e14dcfSSatish Balay #define NTL_UPDATE_REDUCTION 0 11a7e14dcfSSatish Balay #define NTL_UPDATE_INTERPOLATION 1 12a7e14dcfSSatish Balay #define NTL_UPDATE_TYPES 2 13a7e14dcfSSatish Balay 1487f595a5SBarry Smith static const char *NTL_INIT[64] = {"constant","direction","interpolation"}; 15a7e14dcfSSatish Balay 1687f595a5SBarry Smith static const char *NTL_UPDATE[64] = {"reduction","interpolation"}; 17a7e14dcfSSatish Balay 18a7e14dcfSSatish Balay /* Implements Newton's Method with a trust-region, line-search approach for 19a7e14dcfSSatish Balay solving unconstrained minimization problems. A More'-Thuente line search 20a7e14dcfSSatish Balay is used to guarantee that the bfgs preconditioner remains positive 21a7e14dcfSSatish Balay definite. */ 22a7e14dcfSSatish Balay 23a7e14dcfSSatish Balay #define NTL_NEWTON 0 24a7e14dcfSSatish Balay #define NTL_BFGS 1 25a7e14dcfSSatish Balay #define NTL_SCALED_GRADIENT 2 26a7e14dcfSSatish Balay #define NTL_GRADIENT 3 27a7e14dcfSSatish Balay 28441846f8SBarry Smith static PetscErrorCode TaoSolve_NTL(Tao tao) 29a7e14dcfSSatish Balay { 30a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 3155119615STodd Munson KSPType ksp_type; 320ad3a497SAlp Dener PetscBool is_nash,is_stcg,is_gltr,is_bfgs,is_jacobi,is_symmetric,sym_set; 33a7e14dcfSSatish Balay KSPConvergedReason ksp_reason; 3455119615STodd Munson PC pc; 35e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 36a7e14dcfSSatish Balay 37a7e14dcfSSatish Balay PetscReal fmin, ftrial, prered, actred, kappa, sigma; 38a7e14dcfSSatish Balay PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 39a7e14dcfSSatish Balay PetscReal f, fold, gdx, gnorm; 40a7e14dcfSSatish Balay PetscReal step = 1.0; 41a7e14dcfSSatish Balay 42a7e14dcfSSatish Balay PetscReal norm_d = 0.0; 43a7e14dcfSSatish Balay PetscErrorCode ierr; 44a7e14dcfSSatish Balay PetscInt stepType; 458931d482SJason Sarich PetscInt its; 46a7e14dcfSSatish Balay 47a7e14dcfSSatish Balay PetscInt bfgsUpdates = 0; 48a7e14dcfSSatish Balay PetscInt needH; 49a7e14dcfSSatish Balay 50a7e14dcfSSatish Balay PetscInt i_max = 5; 51a7e14dcfSSatish Balay PetscInt j_max = 1; 52a7e14dcfSSatish Balay PetscInt i, j, n, N; 53a7e14dcfSSatish Balay 54a7e14dcfSSatish Balay PetscInt tr_reject; 55a7e14dcfSSatish Balay 56a7e14dcfSSatish Balay PetscFunctionBegin; 57a7e14dcfSSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 589dcef436SStefano Zampini ierr = PetscInfo(tao,"WARNING: Variable bounds have been set but will be ignored by ntl algorithm\n");CHKERRQ(ierr); 59a7e14dcfSSatish Balay } 60a7e14dcfSSatish Balay 6155119615STodd Munson ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr); 6205de396fSBarry Smith ierr = PetscStrcmp(ksp_type,KSPNASH,&is_nash);CHKERRQ(ierr); 6305de396fSBarry Smith ierr = PetscStrcmp(ksp_type,KSPSTCG,&is_stcg);CHKERRQ(ierr); 6405de396fSBarry Smith ierr = PetscStrcmp(ksp_type,KSPGLTR,&is_gltr);CHKERRQ(ierr); 65*3c859ba3SBarry Smith PetscCheck(is_nash || is_stcg || is_gltr,PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"TAO_NTR requires nash, stcg, or gltr for the KSP"); 66a7e14dcfSSatish Balay 6755119615STodd Munson /* Initialize the radius and modify if it is too large or small */ 6855119615STodd Munson tao->trust = tao->trust0; 69a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 70a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 71a7e14dcfSSatish Balay 720c51296cSAlp Dener /* Allocate the vectors needed for the BFGS approximation */ 730c51296cSAlp Dener ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 740c51296cSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs);CHKERRQ(ierr); 750c51296cSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi);CHKERRQ(ierr); 760c51296cSAlp Dener if (is_bfgs) { 770c51296cSAlp Dener tl->bfgs_pre = pc; 780c51296cSAlp Dener ierr = PCLMVMGetMatLMVM(tl->bfgs_pre, &tl->M);CHKERRQ(ierr); 79a7e14dcfSSatish Balay ierr = VecGetLocalSize(tao->solution, &n);CHKERRQ(ierr); 80a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution, &N);CHKERRQ(ierr); 810c51296cSAlp Dener ierr = MatSetSizes(tl->M, n, n, N, N);CHKERRQ(ierr); 82cd929ea3SAlp Dener ierr = MatLMVMAllocate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 830ad3a497SAlp Dener ierr = MatIsSymmetricKnown(tl->M, &sym_set, &is_symmetric);CHKERRQ(ierr); 84*3c859ba3SBarry Smith PetscCheck(sym_set && is_symmetric,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric."); 850c51296cSAlp Dener } else if (is_jacobi) { 860c51296cSAlp Dener ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 87a7e14dcfSSatish Balay } 88a7e14dcfSSatish Balay 89a7e14dcfSSatish Balay /* Check convergence criteria */ 90a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 91a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 92*3c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 93a7e14dcfSSatish Balay needH = 1; 94a7e14dcfSSatish Balay 953ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 963ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 973ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 983ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 993ecd9318SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 100a7e14dcfSSatish Balay 10155119615STodd Munson /* Initialize trust-region radius */ 102a7e14dcfSSatish Balay switch(tl->init_type) { 103a7e14dcfSSatish Balay case NTL_INIT_CONSTANT: 104a7e14dcfSSatish Balay /* Use the initial radius specified */ 105a7e14dcfSSatish Balay break; 106a7e14dcfSSatish Balay 107a7e14dcfSSatish Balay case NTL_INIT_INTERPOLATION: 108a7e14dcfSSatish Balay /* Use the initial radius specified */ 109a7e14dcfSSatish Balay max_radius = 0.0; 110a7e14dcfSSatish Balay 111a7e14dcfSSatish Balay for (j = 0; j < j_max; ++j) { 112a7e14dcfSSatish Balay fmin = f; 113a7e14dcfSSatish Balay sigma = 0.0; 114a7e14dcfSSatish Balay 115a7e14dcfSSatish Balay if (needH) { 116ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 117a7e14dcfSSatish Balay needH = 0; 118a7e14dcfSSatish Balay } 119a7e14dcfSSatish Balay 120a7e14dcfSSatish Balay for (i = 0; i < i_max; ++i) { 121a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 122a7e14dcfSSatish Balay ierr = VecAXPY(tl->W, -tao->trust/gnorm, tao->gradient);CHKERRQ(ierr); 123a7e14dcfSSatish Balay 124a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 125a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 126a7e14dcfSSatish Balay tau = tl->gamma1_i; 12753506e15SBarry Smith } else { 128a7e14dcfSSatish Balay if (ftrial < fmin) { 129a7e14dcfSSatish Balay fmin = ftrial; 130a7e14dcfSSatish Balay sigma = -tao->trust / gnorm; 131a7e14dcfSSatish Balay } 132a7e14dcfSSatish Balay 133a7e14dcfSSatish Balay ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 134a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr); 135a7e14dcfSSatish Balay 136a7e14dcfSSatish Balay prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm)); 137a7e14dcfSSatish Balay actred = f - ftrial; 13853506e15SBarry Smith if ((PetscAbsScalar(actred) <= tl->epsilon) && (PetscAbsScalar(prered) <= tl->epsilon)) { 139a7e14dcfSSatish Balay kappa = 1.0; 14053506e15SBarry Smith } else { 141a7e14dcfSSatish Balay kappa = actred / prered; 142a7e14dcfSSatish Balay } 143a7e14dcfSSatish Balay 144a7e14dcfSSatish Balay tau_1 = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust + (1.0 - tl->theta_i) * prered - actred); 145a7e14dcfSSatish Balay tau_2 = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust - (1.0 + tl->theta_i) * prered + actred); 146a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 147a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 148a7e14dcfSSatish Balay 14918cfbf8eSSatish Balay if (PetscAbsScalar(kappa - (PetscReal)1.0) <= tl->mu1_i) { 150a7e14dcfSSatish Balay /* Great agreement */ 151a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 152a7e14dcfSSatish Balay 153a7e14dcfSSatish Balay if (tau_max < 1.0) { 154a7e14dcfSSatish Balay tau = tl->gamma3_i; 15553506e15SBarry Smith } else if (tau_max > tl->gamma4_i) { 156a7e14dcfSSatish Balay tau = tl->gamma4_i; 15753506e15SBarry Smith } else if (tau_1 >= 1.0 && tau_1 <= tl->gamma4_i && tau_2 < 1.0) { 158a7e14dcfSSatish Balay tau = tau_1; 15953506e15SBarry Smith } else if (tau_2 >= 1.0 && tau_2 <= tl->gamma4_i && tau_1 < 1.0) { 160a7e14dcfSSatish Balay tau = tau_2; 16153506e15SBarry Smith } else { 162a7e14dcfSSatish Balay tau = tau_max; 163a7e14dcfSSatish Balay } 16418cfbf8eSSatish Balay } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= tl->mu2_i) { 165a7e14dcfSSatish Balay /* Good agreement */ 166a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 167a7e14dcfSSatish Balay 168a7e14dcfSSatish Balay if (tau_max < tl->gamma2_i) { 169a7e14dcfSSatish Balay tau = tl->gamma2_i; 17053506e15SBarry Smith } else if (tau_max > tl->gamma3_i) { 171a7e14dcfSSatish Balay tau = tl->gamma3_i; 17253506e15SBarry Smith } else { 173a7e14dcfSSatish Balay tau = tau_max; 174a7e14dcfSSatish Balay } 17553506e15SBarry Smith } else { 176a7e14dcfSSatish Balay /* Not good agreement */ 177a7e14dcfSSatish Balay if (tau_min > 1.0) { 178a7e14dcfSSatish Balay tau = tl->gamma2_i; 17953506e15SBarry Smith } else if (tau_max < tl->gamma1_i) { 180a7e14dcfSSatish Balay tau = tl->gamma1_i; 18153506e15SBarry Smith } else if ((tau_min < tl->gamma1_i) && (tau_max >= 1.0)) { 182a7e14dcfSSatish Balay tau = tl->gamma1_i; 18353506e15SBarry Smith } else if ((tau_1 >= tl->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < tl->gamma1_i) || (tau_2 >= 1.0))) { 184a7e14dcfSSatish Balay tau = tau_1; 18553506e15SBarry Smith } else if ((tau_2 >= tl->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < tl->gamma1_i) || (tau_2 >= 1.0))) { 186a7e14dcfSSatish Balay tau = tau_2; 18753506e15SBarry Smith } else { 188a7e14dcfSSatish Balay tau = tau_max; 189a7e14dcfSSatish Balay } 190a7e14dcfSSatish Balay } 191a7e14dcfSSatish Balay } 192a7e14dcfSSatish Balay tao->trust = tau * tao->trust; 193a7e14dcfSSatish Balay } 194a7e14dcfSSatish Balay 195a7e14dcfSSatish Balay if (fmin < f) { 196a7e14dcfSSatish Balay f = fmin; 197a7e14dcfSSatish Balay ierr = VecAXPY(tao->solution, sigma, tao->gradient);CHKERRQ(ierr); 198a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr); 199a7e14dcfSSatish Balay 200a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 201*3c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 202a7e14dcfSSatish Balay needH = 1; 203a7e14dcfSSatish Balay 2043ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 2053ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 2063ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 2073ecd9318SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 208a7e14dcfSSatish Balay } 209a7e14dcfSSatish Balay } 210a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, max_radius); 211a7e14dcfSSatish Balay 212a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 213a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 214a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 215a7e14dcfSSatish Balay break; 216a7e14dcfSSatish Balay 217a7e14dcfSSatish Balay default: 218a7e14dcfSSatish Balay /* Norm of the first direction will initialize radius */ 219a7e14dcfSSatish Balay tao->trust = 0.0; 220a7e14dcfSSatish Balay break; 221a7e14dcfSSatish Balay } 222a7e14dcfSSatish Balay 223a7e14dcfSSatish Balay /* Set counter for gradient/reset steps */ 224a7e14dcfSSatish Balay tl->ntrust = 0; 225a7e14dcfSSatish Balay tl->newt = 0; 226a7e14dcfSSatish Balay tl->bfgs = 0; 227a7e14dcfSSatish Balay tl->grad = 0; 228a7e14dcfSSatish Balay 229a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 2303ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 231e1e80dc8SAlp Dener /* Call general purpose update function */ 232e1e80dc8SAlp Dener if (tao->ops->update) { 2338fcddce6SStefano Zampini ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr); 234e1e80dc8SAlp Dener } 2358931d482SJason Sarich ++tao->niter; 236ae93cb3cSJason Sarich tao->ksp_its=0; 237a7e14dcfSSatish Balay /* Compute the Hessian */ 238a7e14dcfSSatish Balay if (needH) { 239ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 240a7e14dcfSSatish Balay } 241a7e14dcfSSatish Balay 2420c51296cSAlp Dener if (tl->bfgs_pre) { 243a7e14dcfSSatish Balay /* Update the limited memory preconditioner */ 244a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M,tao->solution, tao->gradient);CHKERRQ(ierr); 245a7e14dcfSSatish Balay ++bfgsUpdates; 246a7e14dcfSSatish Balay } 24723ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre);CHKERRQ(ierr); 248a7e14dcfSSatish Balay /* Solve the Newton system of equations */ 249ba7fe8fbSTodd Munson ierr = KSPCGSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 250a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 251a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 252a7e14dcfSSatish Balay tao->ksp_its+=its; 253ae93cb3cSJason Sarich tao->ksp_tot_its+=its; 254ba7fe8fbSTodd Munson ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 255a7e14dcfSSatish Balay 256a7e14dcfSSatish Balay if (0.0 == tao->trust) { 257a7e14dcfSSatish Balay /* Radius was uninitialized; use the norm of the direction */ 258a7e14dcfSSatish Balay if (norm_d > 0.0) { 259a7e14dcfSSatish Balay tao->trust = norm_d; 260a7e14dcfSSatish Balay 261a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 262a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 263a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 26453506e15SBarry Smith } else { 265a7e14dcfSSatish Balay /* The direction was bad; set radius to default value and re-solve 266a7e14dcfSSatish Balay the trust-region subproblem to get a direction */ 267a7e14dcfSSatish Balay tao->trust = tao->trust0; 268a7e14dcfSSatish Balay 269a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 270a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 271a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 272a7e14dcfSSatish Balay 273ba7fe8fbSTodd Munson ierr = KSPCGSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 274a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 275a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 276a7e14dcfSSatish Balay tao->ksp_its+=its; 2772d9aa51bSJason Sarich tao->ksp_tot_its+=its; 278ba7fe8fbSTodd Munson ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 279a7e14dcfSSatish Balay 280*3c859ba3SBarry Smith PetscCheck(norm_d != 0.0,PetscObjectComm((PetscObject)tao),PETSC_ERR_PLIB, "Initial direction zero"); 281a7e14dcfSSatish Balay } 282a7e14dcfSSatish Balay } 283a7e14dcfSSatish Balay 284a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 285a7e14dcfSSatish Balay ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 2860c51296cSAlp Dener if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (tl->bfgs_pre)) { 287a7e14dcfSSatish Balay /* Preconditioner is numerically indefinite; reset the 288a7e14dcfSSatish Balay approximate if using BFGS preconditioning. */ 289cd929ea3SAlp Dener ierr = MatLMVMReset(tl->M, PETSC_FALSE);CHKERRQ(ierr); 290a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 291a7e14dcfSSatish Balay bfgsUpdates = 1; 292a7e14dcfSSatish Balay } 293a7e14dcfSSatish Balay 294a7e14dcfSSatish Balay /* Check trust-region reduction conditions */ 295a7e14dcfSSatish Balay tr_reject = 0; 296a7e14dcfSSatish Balay if (NTL_UPDATE_REDUCTION == tl->update_type) { 297a7e14dcfSSatish Balay /* Get predicted reduction */ 298ba7fe8fbSTodd Munson ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 299a7e14dcfSSatish Balay if (prered >= 0.0) { 300a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot 301a7e14dcfSSatish Balay happen in infinite precision arithmetic. Step should 302a7e14dcfSSatish Balay be rejected! */ 303a7e14dcfSSatish Balay tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 304a7e14dcfSSatish Balay tr_reject = 1; 30553506e15SBarry Smith } else { 306a7e14dcfSSatish Balay /* Compute trial step and function value */ 307a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 308a7e14dcfSSatish Balay ierr = VecAXPY(tl->W, 1.0, tao->stepdirection);CHKERRQ(ierr); 309a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 310a7e14dcfSSatish Balay 311a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 312a7e14dcfSSatish Balay tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 313a7e14dcfSSatish Balay tr_reject = 1; 31453506e15SBarry Smith } else { 315a7e14dcfSSatish Balay /* Compute and actual reduction */ 316a7e14dcfSSatish Balay actred = f - ftrial; 317a7e14dcfSSatish Balay prered = -prered; 318a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= tl->epsilon) && 319a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= tl->epsilon)) { 320a7e14dcfSSatish Balay kappa = 1.0; 32153506e15SBarry Smith } else { 322a7e14dcfSSatish Balay kappa = actred / prered; 323a7e14dcfSSatish Balay } 324a7e14dcfSSatish Balay 325a7e14dcfSSatish Balay /* Accept of reject the step and update radius */ 326a7e14dcfSSatish Balay if (kappa < tl->eta1) { 327a7e14dcfSSatish Balay /* Reject the step */ 328a7e14dcfSSatish Balay tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 329a7e14dcfSSatish Balay tr_reject = 1; 33053506e15SBarry Smith } else { 331a7e14dcfSSatish Balay /* Accept the step */ 332a7e14dcfSSatish Balay if (kappa < tl->eta2) { 333a7e14dcfSSatish Balay /* Marginal bad step */ 334a7e14dcfSSatish Balay tao->trust = tl->alpha2 * PetscMin(tao->trust, norm_d); 33553506e15SBarry Smith } else if (kappa < tl->eta3) { 336a7e14dcfSSatish Balay /* Reasonable step */ 337a7e14dcfSSatish Balay tao->trust = tl->alpha3 * tao->trust; 33853506e15SBarry Smith } else if (kappa < tl->eta4) { 339a7e14dcfSSatish Balay /* Good step */ 340a7e14dcfSSatish Balay tao->trust = PetscMax(tl->alpha4 * norm_d, tao->trust); 34153506e15SBarry Smith } else { 342a7e14dcfSSatish Balay /* Very good step */ 343a7e14dcfSSatish Balay tao->trust = PetscMax(tl->alpha5 * norm_d, tao->trust); 344a7e14dcfSSatish Balay } 345a7e14dcfSSatish Balay } 346a7e14dcfSSatish Balay } 347a7e14dcfSSatish Balay } 34853506e15SBarry Smith } else { 349a7e14dcfSSatish Balay /* Get predicted reduction */ 350ba7fe8fbSTodd Munson ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 351a7e14dcfSSatish Balay if (prered >= 0.0) { 352a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot 353a7e14dcfSSatish Balay happen in infinite precision arithmetic. Step should 354a7e14dcfSSatish Balay be rejected! */ 355a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 356a7e14dcfSSatish Balay tr_reject = 1; 35753506e15SBarry Smith } else { 358a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 359a7e14dcfSSatish Balay ierr = VecAXPY(tl->W, 1.0, tao->stepdirection);CHKERRQ(ierr); 360a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 361a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 362a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 363a7e14dcfSSatish Balay tr_reject = 1; 36453506e15SBarry Smith } else { 365a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 366a7e14dcfSSatish Balay 367a7e14dcfSSatish Balay actred = f - ftrial; 368a7e14dcfSSatish Balay prered = -prered; 369a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= tl->epsilon) && 370a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= tl->epsilon)) { 371a7e14dcfSSatish Balay kappa = 1.0; 37253506e15SBarry Smith } else { 373a7e14dcfSSatish Balay kappa = actred / prered; 374a7e14dcfSSatish Balay } 375a7e14dcfSSatish Balay 376a7e14dcfSSatish Balay tau_1 = tl->theta * gdx / (tl->theta * gdx - (1.0 - tl->theta) * prered + actred); 377a7e14dcfSSatish Balay tau_2 = tl->theta * gdx / (tl->theta * gdx + (1.0 + tl->theta) * prered - actred); 378a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 379a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 380a7e14dcfSSatish Balay 381a7e14dcfSSatish Balay if (kappa >= 1.0 - tl->mu1) { 382a7e14dcfSSatish Balay /* Great agreement; accept step and update radius */ 383a7e14dcfSSatish Balay if (tau_max < 1.0) { 384a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d); 38553506e15SBarry Smith } else if (tau_max > tl->gamma4) { 386a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->gamma4 * norm_d); 38753506e15SBarry Smith } else { 388a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 389a7e14dcfSSatish Balay } 39053506e15SBarry Smith } else if (kappa >= 1.0 - tl->mu2) { 391a7e14dcfSSatish Balay /* Good agreement */ 392a7e14dcfSSatish Balay 393a7e14dcfSSatish Balay if (tau_max < tl->gamma2) { 394a7e14dcfSSatish Balay tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d); 39553506e15SBarry Smith } else if (tau_max > tl->gamma3) { 396a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d); 397a7e14dcfSSatish Balay } else if (tau_max < 1.0) { 398a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 39953506e15SBarry Smith } else { 400a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 401a7e14dcfSSatish Balay } 40253506e15SBarry Smith } else { 403a7e14dcfSSatish Balay /* Not good agreement */ 404a7e14dcfSSatish Balay if (tau_min > 1.0) { 405a7e14dcfSSatish Balay tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d); 40653506e15SBarry Smith } else if (tau_max < tl->gamma1) { 407a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 40853506e15SBarry Smith } else if ((tau_min < tl->gamma1) && (tau_max >= 1.0)) { 409a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 41053506e15SBarry Smith } else if ((tau_1 >= tl->gamma1) && (tau_1 < 1.0) && ((tau_2 < tl->gamma1) || (tau_2 >= 1.0))) { 411a7e14dcfSSatish Balay tao->trust = tau_1 * PetscMin(tao->trust, norm_d); 41253506e15SBarry Smith } else if ((tau_2 >= tl->gamma1) && (tau_2 < 1.0) && ((tau_1 < tl->gamma1) || (tau_2 >= 1.0))) { 413a7e14dcfSSatish Balay tao->trust = tau_2 * PetscMin(tao->trust, norm_d); 41453506e15SBarry Smith } else { 415a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 416a7e14dcfSSatish Balay } 417a7e14dcfSSatish Balay tr_reject = 1; 418a7e14dcfSSatish Balay } 419a7e14dcfSSatish Balay } 420a7e14dcfSSatish Balay } 421a7e14dcfSSatish Balay } 422a7e14dcfSSatish Balay 423a7e14dcfSSatish Balay if (tr_reject) { 424a7e14dcfSSatish Balay /* The trust-region constraints rejected the step. Apply a linesearch. 425a7e14dcfSSatish Balay Check for descent direction. */ 426a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 427a7e14dcfSSatish Balay if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 428a7e14dcfSSatish Balay /* Newton step is not descent or direction produced Inf or NaN */ 429a7e14dcfSSatish Balay 4300c51296cSAlp Dener if (!tl->bfgs_pre) { 431a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and updated 432a7e14dcfSSatish Balay Must use gradient direction in this case */ 433a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 434a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 435a7e14dcfSSatish Balay ++tl->grad; 436a7e14dcfSSatish Balay stepType = NTL_GRADIENT; 43753506e15SBarry Smith } else { 438a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 439cd929ea3SAlp Dener ierr = MatSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 440a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 441a7e14dcfSSatish Balay 442a7e14dcfSSatish Balay /* Check for success (descent direction) */ 443a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 444a7e14dcfSSatish Balay if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) { 445a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced not a number 446a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case because 447a7e14dcfSSatish Balay the first solve produces the scaled gradient direction, 448a7e14dcfSSatish Balay which is guaranteed to be descent */ 449a7e14dcfSSatish Balay 450a7e14dcfSSatish Balay /* Use steepest descent direction (scaled) */ 451cd929ea3SAlp Dener ierr = MatLMVMReset(tl->M, PETSC_FALSE);CHKERRQ(ierr); 452a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 453cd929ea3SAlp Dener ierr = MatSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 454a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 455a7e14dcfSSatish Balay 456a7e14dcfSSatish Balay bfgsUpdates = 1; 4570c51296cSAlp Dener ++tl->grad; 4580c51296cSAlp Dener stepType = NTL_GRADIENT; 45953506e15SBarry Smith } else { 4600c51296cSAlp Dener ierr = MatLMVMGetUpdateCount(tl->M, &bfgsUpdates);CHKERRQ(ierr); 461a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 462a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 4630c51296cSAlp Dener ++tl->grad; 4640c51296cSAlp Dener stepType = NTL_GRADIENT; 46553506e15SBarry Smith } else { 466a7e14dcfSSatish Balay ++tl->bfgs; 467a7e14dcfSSatish Balay stepType = NTL_BFGS; 468a7e14dcfSSatish Balay } 469a7e14dcfSSatish Balay } 470a7e14dcfSSatish Balay } 47153506e15SBarry Smith } else { 472a7e14dcfSSatish Balay /* Computed Newton step is descent */ 473a7e14dcfSSatish Balay ++tl->newt; 474a7e14dcfSSatish Balay stepType = NTL_NEWTON; 475a7e14dcfSSatish Balay } 476a7e14dcfSSatish Balay 477a7e14dcfSSatish Balay /* Perform the linesearch */ 478a7e14dcfSSatish Balay fold = f; 479a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->Xold);CHKERRQ(ierr); 480a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tl->Gold);CHKERRQ(ierr); 481a7e14dcfSSatish Balay 482a7e14dcfSSatish Balay step = 1.0; 483a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 484a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 485a7e14dcfSSatish Balay 48653506e15SBarry Smith while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NTL_GRADIENT) { /* Linesearch failed */ 487a7e14dcfSSatish Balay /* Linesearch failed */ 488a7e14dcfSSatish Balay f = fold; 489a7e14dcfSSatish Balay ierr = VecCopy(tl->Xold, tao->solution);CHKERRQ(ierr); 490a7e14dcfSSatish Balay ierr = VecCopy(tl->Gold, tao->gradient);CHKERRQ(ierr); 491a7e14dcfSSatish Balay 492a7e14dcfSSatish Balay switch(stepType) { 493a7e14dcfSSatish Balay case NTL_NEWTON: 494a7e14dcfSSatish Balay /* Failed to obtain acceptable iterate with Newton step */ 495a7e14dcfSSatish Balay 4960c51296cSAlp Dener if (tl->bfgs_pre) { 497a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and being updated 498a7e14dcfSSatish Balay Must use gradient direction in this case */ 499a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 500a7e14dcfSSatish Balay ++tl->grad; 501a7e14dcfSSatish Balay stepType = NTL_GRADIENT; 50253506e15SBarry Smith } else { 503a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 504cd929ea3SAlp Dener ierr = MatSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 505a7e14dcfSSatish Balay 506a7e14dcfSSatish Balay /* Check for success (descent direction) */ 507a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 508a7e14dcfSSatish Balay if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 509a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced 510a7e14dcfSSatish Balay not a number. We can assert bfgsUpdates > 1 in this case 511a7e14dcfSSatish Balay Use steepest descent direction (scaled) */ 512cd929ea3SAlp Dener ierr = MatLMVMReset(tl->M, PETSC_FALSE);CHKERRQ(ierr); 513a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 514cd929ea3SAlp Dener ierr = MatSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 515a7e14dcfSSatish Balay 516a7e14dcfSSatish Balay bfgsUpdates = 1; 5170c51296cSAlp Dener ++tl->grad; 5180c51296cSAlp Dener stepType = NTL_GRADIENT; 51953506e15SBarry Smith } else { 5200c51296cSAlp Dener ierr = MatLMVMGetUpdateCount(tl->M, &bfgsUpdates);CHKERRQ(ierr); 521a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 522a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 5230c51296cSAlp Dener ++tl->grad; 5240c51296cSAlp Dener stepType = NTL_GRADIENT; 52553506e15SBarry Smith } else { 526a7e14dcfSSatish Balay ++tl->bfgs; 527a7e14dcfSSatish Balay stepType = NTL_BFGS; 528a7e14dcfSSatish Balay } 529a7e14dcfSSatish Balay } 530a7e14dcfSSatish Balay } 531a7e14dcfSSatish Balay break; 532a7e14dcfSSatish Balay 533a7e14dcfSSatish Balay case NTL_BFGS: 534a7e14dcfSSatish Balay /* Can only enter if pc_type == NTL_PC_BFGS 535a7e14dcfSSatish Balay Failed to obtain acceptable iterate with BFGS step 536a7e14dcfSSatish Balay Attempt to use the scaled gradient direction */ 537cd929ea3SAlp Dener ierr = MatLMVMReset(tl->M, PETSC_FALSE);CHKERRQ(ierr); 538a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 539cd929ea3SAlp Dener ierr = MatSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 540a7e14dcfSSatish Balay 541a7e14dcfSSatish Balay bfgsUpdates = 1; 542a7e14dcfSSatish Balay ++tl->grad; 543a7e14dcfSSatish Balay stepType = NTL_GRADIENT; 544a7e14dcfSSatish Balay break; 545a7e14dcfSSatish Balay } 546a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 547a7e14dcfSSatish Balay 548a7e14dcfSSatish Balay /* This may be incorrect; linesearch has values for stepmax and stepmin 549a7e14dcfSSatish Balay that should be reset. */ 550a7e14dcfSSatish Balay step = 1.0; 551a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 552a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 553a7e14dcfSSatish Balay } 554a7e14dcfSSatish Balay 55553506e15SBarry Smith if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 556a7e14dcfSSatish Balay /* Failed to find an improving point */ 557a7e14dcfSSatish Balay f = fold; 558a7e14dcfSSatish Balay ierr = VecCopy(tl->Xold, tao->solution);CHKERRQ(ierr); 559a7e14dcfSSatish Balay ierr = VecCopy(tl->Gold, tao->gradient);CHKERRQ(ierr); 560a7e14dcfSSatish Balay tao->trust = 0.0; 561a7e14dcfSSatish Balay step = 0.0; 562a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 563a7e14dcfSSatish Balay break; 56453506e15SBarry Smith } else if (stepType == NTL_NEWTON) { 565a7e14dcfSSatish Balay if (step < tl->nu1) { 566a7e14dcfSSatish Balay /* Very bad step taken; reduce radius */ 567a7e14dcfSSatish Balay tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust); 56853506e15SBarry Smith } else if (step < tl->nu2) { 569a7e14dcfSSatish Balay /* Reasonably bad step taken; reduce radius */ 570a7e14dcfSSatish Balay tao->trust = tl->omega2 * PetscMin(norm_d, tao->trust); 57153506e15SBarry Smith } else if (step < tl->nu3) { 572a7e14dcfSSatish Balay /* Reasonable step was taken; leave radius alone */ 573a7e14dcfSSatish Balay if (tl->omega3 < 1.0) { 574a7e14dcfSSatish Balay tao->trust = tl->omega3 * PetscMin(norm_d, tao->trust); 57553506e15SBarry Smith } else if (tl->omega3 > 1.0) { 576a7e14dcfSSatish Balay tao->trust = PetscMax(tl->omega3 * norm_d, tao->trust); 577a7e14dcfSSatish Balay } 57853506e15SBarry Smith } else if (step < tl->nu4) { 579a7e14dcfSSatish Balay /* Full step taken; increase the radius */ 580a7e14dcfSSatish Balay tao->trust = PetscMax(tl->omega4 * norm_d, tao->trust); 58153506e15SBarry Smith } else { 582a7e14dcfSSatish Balay /* More than full step taken; increase the radius */ 583a7e14dcfSSatish Balay tao->trust = PetscMax(tl->omega5 * norm_d, tao->trust); 584a7e14dcfSSatish Balay } 58553506e15SBarry Smith } else { 586a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 587a7e14dcfSSatish Balay tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust); 588a7e14dcfSSatish Balay } 58953506e15SBarry Smith } else { 590a7e14dcfSSatish Balay /* Trust-region step is accepted */ 591a7e14dcfSSatish Balay ierr = VecCopy(tl->W, tao->solution);CHKERRQ(ierr); 592a7e14dcfSSatish Balay f = ftrial; 593a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr); 594a7e14dcfSSatish Balay ++tl->ntrust; 595a7e14dcfSSatish Balay } 596a7e14dcfSSatish Balay 597a7e14dcfSSatish Balay /* The radius may have been increased; modify if it is too large */ 598a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 599a7e14dcfSSatish Balay 600e4cb33bbSBarry Smith /* Check for converged */ 601a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 602*3c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER,"User provided compute function generated Not-a-Number"); 603a7e14dcfSSatish Balay needH = 1; 604a7e14dcfSSatish Balay 6053ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 6063ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 6073ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 608a7e14dcfSSatish Balay } 609a7e14dcfSSatish Balay PetscFunctionReturn(0); 610a7e14dcfSSatish Balay } 611a7e14dcfSSatish Balay 612a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 613441846f8SBarry Smith static PetscErrorCode TaoSetUp_NTL(Tao tao) 614a7e14dcfSSatish Balay { 615a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 616a7e14dcfSSatish Balay PetscErrorCode ierr; 617a7e14dcfSSatish Balay 618a7e14dcfSSatish Balay PetscFunctionBegin; 619a7e14dcfSSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); } 620a7e14dcfSSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);} 621a7e14dcfSSatish Balay if (!tl->W) { ierr = VecDuplicate(tao->solution, &tl->W);CHKERRQ(ierr);} 622a7e14dcfSSatish Balay if (!tl->Xold) { ierr = VecDuplicate(tao->solution, &tl->Xold);CHKERRQ(ierr);} 623a7e14dcfSSatish Balay if (!tl->Gold) { ierr = VecDuplicate(tao->solution, &tl->Gold);CHKERRQ(ierr);} 62483c8fe1dSLisandro Dalcin tl->bfgs_pre = NULL; 62583c8fe1dSLisandro Dalcin tl->M = NULL; 626a7e14dcfSSatish Balay PetscFunctionReturn(0); 627a7e14dcfSSatish Balay } 628a7e14dcfSSatish Balay 629a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 630441846f8SBarry Smith static PetscErrorCode TaoDestroy_NTL(Tao tao) 631a7e14dcfSSatish Balay { 632a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 633a7e14dcfSSatish Balay PetscErrorCode ierr; 634a7e14dcfSSatish Balay 635a7e14dcfSSatish Balay PetscFunctionBegin; 636a7e14dcfSSatish Balay if (tao->setupcalled) { 637a7e14dcfSSatish Balay ierr = VecDestroy(&tl->W);CHKERRQ(ierr); 638a7e14dcfSSatish Balay ierr = VecDestroy(&tl->Xold);CHKERRQ(ierr); 639a7e14dcfSSatish Balay ierr = VecDestroy(&tl->Gold);CHKERRQ(ierr); 640a7e14dcfSSatish Balay } 641a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 642a7e14dcfSSatish Balay PetscFunctionReturn(0); 643a7e14dcfSSatish Balay } 644a7e14dcfSSatish Balay 645a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 6464416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_NTL(PetscOptionItems *PetscOptionsObject,Tao tao) 647a7e14dcfSSatish Balay { 648a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 649a7e14dcfSSatish Balay PetscErrorCode ierr; 650a7e14dcfSSatish Balay 651a7e14dcfSSatish Balay PetscFunctionBegin; 6521a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Newton trust region with line search method for unconstrained optimization");CHKERRQ(ierr); 65394ae4db5SBarry Smith ierr = PetscOptionsEList("-tao_ntl_init_type", "radius initialization type", "", NTL_INIT, NTL_INIT_TYPES, NTL_INIT[tl->init_type], &tl->init_type,NULL);CHKERRQ(ierr); 65494ae4db5SBarry Smith ierr = PetscOptionsEList("-tao_ntl_update_type", "radius update type", "", NTL_UPDATE, NTL_UPDATE_TYPES, NTL_UPDATE[tl->update_type], &tl->update_type,NULL);CHKERRQ(ierr); 65594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_eta1", "poor steplength; reduce radius", "", tl->eta1, &tl->eta1,NULL);CHKERRQ(ierr); 65694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_eta2", "reasonable steplength; leave radius alone", "", tl->eta2, &tl->eta2,NULL);CHKERRQ(ierr); 65794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_eta3", "good steplength; increase radius", "", tl->eta3, &tl->eta3,NULL);CHKERRQ(ierr); 65894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_eta4", "excellent steplength; greatly increase radius", "", tl->eta4, &tl->eta4,NULL);CHKERRQ(ierr); 65994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha1", "", "", tl->alpha1, &tl->alpha1,NULL);CHKERRQ(ierr); 66094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha2", "", "", tl->alpha2, &tl->alpha2,NULL);CHKERRQ(ierr); 66194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha3", "", "", tl->alpha3, &tl->alpha3,NULL);CHKERRQ(ierr); 66294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha4", "", "", tl->alpha4, &tl->alpha4,NULL);CHKERRQ(ierr); 66394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha5", "", "", tl->alpha5, &tl->alpha5,NULL);CHKERRQ(ierr); 66494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_nu1", "poor steplength; reduce radius", "", tl->nu1, &tl->nu1,NULL);CHKERRQ(ierr); 66594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_nu2", "reasonable steplength; leave radius alone", "", tl->nu2, &tl->nu2,NULL);CHKERRQ(ierr); 66694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_nu3", "good steplength; increase radius", "", tl->nu3, &tl->nu3,NULL);CHKERRQ(ierr); 66794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_nu4", "excellent steplength; greatly increase radius", "", tl->nu4, &tl->nu4,NULL);CHKERRQ(ierr); 66894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega1", "", "", tl->omega1, &tl->omega1,NULL);CHKERRQ(ierr); 66994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega2", "", "", tl->omega2, &tl->omega2,NULL);CHKERRQ(ierr); 67094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega3", "", "", tl->omega3, &tl->omega3,NULL);CHKERRQ(ierr); 67194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega4", "", "", tl->omega4, &tl->omega4,NULL);CHKERRQ(ierr); 67294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega5", "", "", tl->omega5, &tl->omega5,NULL);CHKERRQ(ierr); 67394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_mu1_i", "", "", tl->mu1_i, &tl->mu1_i,NULL);CHKERRQ(ierr); 67494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_mu2_i", "", "", tl->mu2_i, &tl->mu2_i,NULL);CHKERRQ(ierr); 67594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma1_i", "", "", tl->gamma1_i, &tl->gamma1_i,NULL);CHKERRQ(ierr); 67694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma2_i", "", "", tl->gamma2_i, &tl->gamma2_i,NULL);CHKERRQ(ierr); 67794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma3_i", "", "", tl->gamma3_i, &tl->gamma3_i,NULL);CHKERRQ(ierr); 67894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma4_i", "", "", tl->gamma4_i, &tl->gamma4_i,NULL);CHKERRQ(ierr); 67994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_theta_i", "", "", tl->theta_i, &tl->theta_i,NULL);CHKERRQ(ierr); 68094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_mu1", "", "", tl->mu1, &tl->mu1,NULL);CHKERRQ(ierr); 68194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_mu2", "", "", tl->mu2, &tl->mu2,NULL);CHKERRQ(ierr); 68294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma1", "", "", tl->gamma1, &tl->gamma1,NULL);CHKERRQ(ierr); 68394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma2", "", "", tl->gamma2, &tl->gamma2,NULL);CHKERRQ(ierr); 68494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma3", "", "", tl->gamma3, &tl->gamma3,NULL);CHKERRQ(ierr); 68594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma4", "", "", tl->gamma4, &tl->gamma4,NULL);CHKERRQ(ierr); 68694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_theta", "", "", tl->theta, &tl->theta,NULL);CHKERRQ(ierr); 68794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_min_radius", "lower bound on initial radius", "", tl->min_radius, &tl->min_radius,NULL);CHKERRQ(ierr); 68894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_max_radius", "upper bound on radius", "", tl->max_radius, &tl->max_radius,NULL);CHKERRQ(ierr); 68994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_epsilon", "tolerance used when computing actual and predicted reduction", "", tl->epsilon, &tl->epsilon,NULL);CHKERRQ(ierr); 690a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 691a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 692a7e14dcfSSatish Balay ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 693a7e14dcfSSatish Balay PetscFunctionReturn(0); 694a7e14dcfSSatish Balay } 695a7e14dcfSSatish Balay 696a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 697441846f8SBarry Smith static PetscErrorCode TaoView_NTL(Tao tao, PetscViewer viewer) 698a7e14dcfSSatish Balay { 699a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 700a7e14dcfSSatish Balay PetscBool isascii; 701a7e14dcfSSatish Balay PetscErrorCode ierr; 702a7e14dcfSSatish Balay 703a7e14dcfSSatish Balay PetscFunctionBegin; 704a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 705a7e14dcfSSatish Balay if (isascii) { 706a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 707a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Trust-region steps: %D\n", tl->ntrust);CHKERRQ(ierr); 708a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Newton search steps: %D\n", tl->newt);CHKERRQ(ierr); 709a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "BFGS search steps: %D\n", tl->bfgs);CHKERRQ(ierr); 710a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Gradient search steps: %D\n", tl->grad);CHKERRQ(ierr); 711a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 712a7e14dcfSSatish Balay } 713a7e14dcfSSatish Balay PetscFunctionReturn(0); 714a7e14dcfSSatish Balay } 715a7e14dcfSSatish Balay 716a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 7171522df2eSJason Sarich /*MC 7183850be85SAlp Dener TAONTL - Newton's method with trust region globalization and line search fallback. 7191522df2eSJason Sarich At each iteration, the Newton trust region method solves the system for d 7201522df2eSJason Sarich and performs a line search in the d direction: 7211522df2eSJason Sarich 7221522df2eSJason Sarich min_d .5 dT Hk d + gkT d, s.t. ||d|| < Delta_k 7231522df2eSJason Sarich 7241522df2eSJason Sarich Options Database Keys: 7259d0a60b2SAlp Dener + -tao_ntl_init_type - "constant","direction","interpolation" 7261522df2eSJason Sarich . -tao_ntl_update_type - "reduction","interpolation" 7271522df2eSJason Sarich . -tao_ntl_min_radius - lower bound on trust region radius 7281522df2eSJason Sarich . -tao_ntl_max_radius - upper bound on trust region radius 7291522df2eSJason Sarich . -tao_ntl_epsilon - tolerance for accepting actual / predicted reduction 7301522df2eSJason Sarich . -tao_ntl_mu1_i - mu1 interpolation init factor 7311522df2eSJason Sarich . -tao_ntl_mu2_i - mu2 interpolation init factor 7321522df2eSJason Sarich . -tao_ntl_gamma1_i - gamma1 interpolation init factor 7331522df2eSJason Sarich . -tao_ntl_gamma2_i - gamma2 interpolation init factor 7341522df2eSJason Sarich . -tao_ntl_gamma3_i - gamma3 interpolation init factor 7351522df2eSJason Sarich . -tao_ntl_gamma4_i - gamma4 interpolation init factor 7368966356dSPierre Jolivet . -tao_ntl_theta_i - theta1 interpolation init factor 7371522df2eSJason Sarich . -tao_ntl_eta1 - eta1 reduction update factor 7381522df2eSJason Sarich . -tao_ntl_eta2 - eta2 reduction update factor 7391522df2eSJason Sarich . -tao_ntl_eta3 - eta3 reduction update factor 7401522df2eSJason Sarich . -tao_ntl_eta4 - eta4 reduction update factor 7411522df2eSJason Sarich . -tao_ntl_alpha1 - alpha1 reduction update factor 7421522df2eSJason Sarich . -tao_ntl_alpha2 - alpha2 reduction update factor 7431522df2eSJason Sarich . -tao_ntl_alpha3 - alpha3 reduction update factor 7441522df2eSJason Sarich . -tao_ntl_alpha4 - alpha4 reduction update factor 7451522df2eSJason Sarich . -tao_ntl_alpha4 - alpha4 reduction update factor 7461522df2eSJason Sarich . -tao_ntl_mu1 - mu1 interpolation update 7471522df2eSJason Sarich . -tao_ntl_mu2 - mu2 interpolation update 7481522df2eSJason Sarich . -tao_ntl_gamma1 - gamma1 interpolcation update 7491522df2eSJason Sarich . -tao_ntl_gamma2 - gamma2 interpolcation update 7501522df2eSJason Sarich . -tao_ntl_gamma3 - gamma3 interpolcation update 7511522df2eSJason Sarich . -tao_ntl_gamma4 - gamma4 interpolation update 7521522df2eSJason Sarich - -tao_ntl_theta - theta1 interpolation update 7531522df2eSJason Sarich 7541eb8069cSJason Sarich Level: beginner 7551522df2eSJason Sarich M*/ 756728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NTL(Tao tao) 757a7e14dcfSSatish Balay { 758a7e14dcfSSatish Balay TAO_NTL *tl; 759a7e14dcfSSatish Balay PetscErrorCode ierr; 7608caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 76153506e15SBarry Smith 762a7e14dcfSSatish Balay PetscFunctionBegin; 7633c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&tl);CHKERRQ(ierr); 764a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_NTL; 765a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_NTL; 766a7e14dcfSSatish Balay tao->ops->view = TaoView_NTL; 767a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_NTL; 768a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_NTL; 769a7e14dcfSSatish Balay 7706552cf8aSJason Sarich /* Override default settings (unless already changed) */ 7716552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 50; 7726552cf8aSJason Sarich if (!tao->trust0_changed) tao->trust0 = 100.0; 7736552cf8aSJason Sarich 774a7e14dcfSSatish Balay tao->data = (void*)tl; 775a7e14dcfSSatish Balay 776a7e14dcfSSatish Balay /* Default values for trust-region radius update based on steplength */ 777a7e14dcfSSatish Balay tl->nu1 = 0.25; 778a7e14dcfSSatish Balay tl->nu2 = 0.50; 779a7e14dcfSSatish Balay tl->nu3 = 1.00; 780a7e14dcfSSatish Balay tl->nu4 = 1.25; 781a7e14dcfSSatish Balay 782a7e14dcfSSatish Balay tl->omega1 = 0.25; 783a7e14dcfSSatish Balay tl->omega2 = 0.50; 784a7e14dcfSSatish Balay tl->omega3 = 1.00; 785a7e14dcfSSatish Balay tl->omega4 = 2.00; 786a7e14dcfSSatish Balay tl->omega5 = 4.00; 787a7e14dcfSSatish Balay 788a7e14dcfSSatish Balay /* Default values for trust-region radius update based on reduction */ 789a7e14dcfSSatish Balay tl->eta1 = 1.0e-4; 790a7e14dcfSSatish Balay tl->eta2 = 0.25; 791a7e14dcfSSatish Balay tl->eta3 = 0.50; 792a7e14dcfSSatish Balay tl->eta4 = 0.90; 793a7e14dcfSSatish Balay 794a7e14dcfSSatish Balay tl->alpha1 = 0.25; 795a7e14dcfSSatish Balay tl->alpha2 = 0.50; 796a7e14dcfSSatish Balay tl->alpha3 = 1.00; 797a7e14dcfSSatish Balay tl->alpha4 = 2.00; 798a7e14dcfSSatish Balay tl->alpha5 = 4.00; 799a7e14dcfSSatish Balay 800a7e14dcfSSatish Balay /* Default values for trust-region radius update based on interpolation */ 801a7e14dcfSSatish Balay tl->mu1 = 0.10; 802a7e14dcfSSatish Balay tl->mu2 = 0.50; 803a7e14dcfSSatish Balay 804a7e14dcfSSatish Balay tl->gamma1 = 0.25; 805a7e14dcfSSatish Balay tl->gamma2 = 0.50; 806a7e14dcfSSatish Balay tl->gamma3 = 2.00; 807a7e14dcfSSatish Balay tl->gamma4 = 4.00; 808a7e14dcfSSatish Balay 809a7e14dcfSSatish Balay tl->theta = 0.05; 810a7e14dcfSSatish Balay 811a7e14dcfSSatish Balay /* Default values for trust region initialization based on interpolation */ 812a7e14dcfSSatish Balay tl->mu1_i = 0.35; 813a7e14dcfSSatish Balay tl->mu2_i = 0.50; 814a7e14dcfSSatish Balay 815a7e14dcfSSatish Balay tl->gamma1_i = 0.0625; 816a7e14dcfSSatish Balay tl->gamma2_i = 0.5; 817a7e14dcfSSatish Balay tl->gamma3_i = 2.0; 818a7e14dcfSSatish Balay tl->gamma4_i = 5.0; 819a7e14dcfSSatish Balay 820a7e14dcfSSatish Balay tl->theta_i = 0.25; 821a7e14dcfSSatish Balay 822a7e14dcfSSatish Balay /* Remaining parameters */ 823a7e14dcfSSatish Balay tl->min_radius = 1.0e-10; 824a7e14dcfSSatish Balay tl->max_radius = 1.0e10; 825a7e14dcfSSatish Balay tl->epsilon = 1.0e-6; 826a7e14dcfSSatish Balay 827a7e14dcfSSatish Balay tl->init_type = NTL_INIT_INTERPOLATION; 828a7e14dcfSSatish Balay tl->update_type = NTL_UPDATE_REDUCTION; 829a7e14dcfSSatish Balay 830a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 83163b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch,(PetscObject)tao,1);CHKERRQ(ierr); 832a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 833441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch, tao);CHKERRQ(ierr); 8345d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 835a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 83663b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp,(PetscObject)tao,1);CHKERRQ(ierr); 8375d527766SPatrick Farrell ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr); 838cbf034f8SAlp Dener ierr = KSPAppendOptionsPrefix(tao->ksp,"tao_ntl_");CHKERRQ(ierr); 83905de396fSBarry Smith ierr = KSPSetType(tao->ksp,KSPSTCG);CHKERRQ(ierr); 840a7e14dcfSSatish Balay PetscFunctionReturn(0); 841a7e14dcfSSatish Balay } 842