1aaa7dc30SBarry Smith #include <../src/tao/matrix/lmvmmat.h> 2aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/ntl/ntl.h> 3a7e14dcfSSatish Balay 4aaa7dc30SBarry Smith #include <petscksp.h> 5a7e14dcfSSatish Balay 6a7e14dcfSSatish Balay #define NTL_KSP_NASH 0 7a7e14dcfSSatish Balay #define NTL_KSP_STCG 1 8a7e14dcfSSatish Balay #define NTL_KSP_GLTR 2 9a7e14dcfSSatish Balay #define NTL_KSP_TYPES 3 10a7e14dcfSSatish Balay 11a7e14dcfSSatish Balay #define NTL_PC_NONE 0 12a7e14dcfSSatish Balay #define NTL_PC_AHESS 1 13a7e14dcfSSatish Balay #define NTL_PC_BFGS 2 14a7e14dcfSSatish Balay #define NTL_PC_PETSC 3 15a7e14dcfSSatish Balay #define NTL_PC_TYPES 4 16a7e14dcfSSatish Balay 17a7e14dcfSSatish Balay #define BFGS_SCALE_AHESS 0 18a7e14dcfSSatish Balay #define BFGS_SCALE_BFGS 1 19a7e14dcfSSatish Balay #define BFGS_SCALE_TYPES 2 20a7e14dcfSSatish Balay 21a7e14dcfSSatish Balay #define NTL_INIT_CONSTANT 0 22a7e14dcfSSatish Balay #define NTL_INIT_DIRECTION 1 23a7e14dcfSSatish Balay #define NTL_INIT_INTERPOLATION 2 24a7e14dcfSSatish Balay #define NTL_INIT_TYPES 3 25a7e14dcfSSatish Balay 26a7e14dcfSSatish Balay #define NTL_UPDATE_REDUCTION 0 27a7e14dcfSSatish Balay #define NTL_UPDATE_INTERPOLATION 1 28a7e14dcfSSatish Balay #define NTL_UPDATE_TYPES 2 29a7e14dcfSSatish Balay 3087f595a5SBarry Smith static const char *NTL_KSP[64] = {"nash", "stcg", "gltr"}; 31a7e14dcfSSatish Balay 3287f595a5SBarry Smith static const char *NTL_PC[64] = {"none", "ahess", "bfgs", "petsc"}; 33a7e14dcfSSatish Balay 3487f595a5SBarry Smith static const char *BFGS_SCALE[64] = {"ahess", "bfgs"}; 35a7e14dcfSSatish Balay 3687f595a5SBarry Smith static const char *NTL_INIT[64] = {"constant", "direction", "interpolation"}; 37a7e14dcfSSatish Balay 3887f595a5SBarry Smith static const char *NTL_UPDATE[64] = {"reduction", "interpolation"}; 39a7e14dcfSSatish Balay 40a7e14dcfSSatish Balay /* Routine for BFGS preconditioner */ 41a7e14dcfSSatish Balay 42a7e14dcfSSatish Balay #undef __FUNCT__ 43a7e14dcfSSatish Balay #define __FUNCT__ "MatLMVMSolveShell" 44a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x) 45a7e14dcfSSatish Balay { 46a7e14dcfSSatish Balay PetscErrorCode ierr; 47a7e14dcfSSatish Balay Mat M; 48a7e14dcfSSatish Balay 49a7e14dcfSSatish Balay PetscFunctionBegin; 50a7e14dcfSSatish Balay PetscValidHeaderSpecific(pc,PC_CLASSID,1); 51a7e14dcfSSatish Balay PetscValidHeaderSpecific(b,VEC_CLASSID,2); 52a7e14dcfSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,3); 53a7e14dcfSSatish Balay ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr); 54a7e14dcfSSatish Balay ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr); 55a7e14dcfSSatish Balay PetscFunctionReturn(0); 56a7e14dcfSSatish Balay } 57a7e14dcfSSatish Balay 58a7e14dcfSSatish Balay /* Implements Newton's Method with a trust-region, line-search approach for 59a7e14dcfSSatish Balay solving unconstrained minimization problems. A More'-Thuente line search 60a7e14dcfSSatish Balay is used to guarantee that the bfgs preconditioner remains positive 61a7e14dcfSSatish Balay definite. */ 62a7e14dcfSSatish Balay 63a7e14dcfSSatish Balay #define NTL_NEWTON 0 64a7e14dcfSSatish Balay #define NTL_BFGS 1 65a7e14dcfSSatish Balay #define NTL_SCALED_GRADIENT 2 66a7e14dcfSSatish Balay #define NTL_GRADIENT 3 67a7e14dcfSSatish Balay 68a7e14dcfSSatish Balay #undef __FUNCT__ 69a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_NTL" 70441846f8SBarry Smith static PetscErrorCode TaoSolve_NTL(Tao tao) 71a7e14dcfSSatish Balay { 72a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 73a7e14dcfSSatish Balay PC pc; 74a7e14dcfSSatish Balay KSPConvergedReason ksp_reason; 75e4cb33bbSBarry Smith TaoConvergedReason reason; 76e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 77a7e14dcfSSatish Balay 78a7e14dcfSSatish Balay PetscReal fmin, ftrial, prered, actred, kappa, sigma; 79a7e14dcfSSatish Balay PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 80a7e14dcfSSatish Balay PetscReal f, fold, gdx, gnorm; 81a7e14dcfSSatish Balay PetscReal step = 1.0; 82a7e14dcfSSatish Balay 83a7e14dcfSSatish Balay PetscReal delta; 84a7e14dcfSSatish Balay PetscReal norm_d = 0.0; 85a7e14dcfSSatish Balay PetscErrorCode ierr; 86a7e14dcfSSatish Balay PetscInt stepType; 878931d482SJason Sarich PetscInt its; 88a7e14dcfSSatish Balay 89a7e14dcfSSatish Balay PetscInt bfgsUpdates = 0; 90a7e14dcfSSatish Balay PetscInt needH; 91a7e14dcfSSatish Balay 92a7e14dcfSSatish Balay PetscInt i_max = 5; 93a7e14dcfSSatish Balay PetscInt j_max = 1; 94a7e14dcfSSatish Balay PetscInt i, j, n, N; 95a7e14dcfSSatish Balay 96a7e14dcfSSatish Balay PetscInt tr_reject; 97a7e14dcfSSatish Balay 98a7e14dcfSSatish Balay PetscFunctionBegin; 99a7e14dcfSSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 100a7e14dcfSSatish Balay ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by ntl algorithm\n");CHKERRQ(ierr); 101a7e14dcfSSatish Balay } 102a7e14dcfSSatish Balay 103a7e14dcfSSatish Balay /* Initialize trust-region radius */ 104a7e14dcfSSatish Balay tao->trust = tao->trust0; 105a7e14dcfSSatish Balay 106a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 107a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 108a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 109a7e14dcfSSatish Balay 110a7e14dcfSSatish Balay if (NTL_PC_BFGS == tl->pc_type && !tl->M) { 111a7e14dcfSSatish Balay ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 112a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 113a7e14dcfSSatish Balay ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&tl->M);CHKERRQ(ierr); 114a7e14dcfSSatish Balay ierr = MatLMVMAllocateVectors(tl->M,tao->solution);CHKERRQ(ierr); 115a7e14dcfSSatish Balay } 116a7e14dcfSSatish Balay 117a7e14dcfSSatish Balay /* Check convergence criteria */ 118a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 119a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 12053506e15SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 121a7e14dcfSSatish Balay needH = 1; 122a7e14dcfSSatish Balay 1238931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 12453506e15SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 125a7e14dcfSSatish Balay 126a7e14dcfSSatish Balay /* Create vectors for the limited memory preconditioner */ 12753506e15SBarry Smith if ((NTL_PC_BFGS == tl->pc_type) && (BFGS_SCALE_BFGS != tl->bfgs_scale_type)) { 128a7e14dcfSSatish Balay if (!tl->Diag) { 129a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tl->Diag);CHKERRQ(ierr); 130a7e14dcfSSatish Balay } 131a7e14dcfSSatish Balay } 132a7e14dcfSSatish Balay 133a7e14dcfSSatish Balay /* Modify the linear solver to a conjugate gradient method */ 134a7e14dcfSSatish Balay switch(tl->ksp_type) { 135a7e14dcfSSatish Balay case NTL_KSP_NASH: 136*ba7fe8fbSTodd Munson ierr = KSPSetType(tao->ksp, KSPCGNASH);CHKERRQ(ierr); 1371a1499c8SBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 138a7e14dcfSSatish Balay break; 139a7e14dcfSSatish Balay 140a7e14dcfSSatish Balay case NTL_KSP_STCG: 141*ba7fe8fbSTodd Munson ierr = KSPSetType(tao->ksp, KSPCGSTCG);CHKERRQ(ierr); 1421a1499c8SBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 143a7e14dcfSSatish Balay break; 144a7e14dcfSSatish Balay 145a7e14dcfSSatish Balay default: 146*ba7fe8fbSTodd Munson ierr = KSPSetType(tao->ksp, KSPCGGLTR);CHKERRQ(ierr); 1471a1499c8SBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 148a7e14dcfSSatish Balay break; 149a7e14dcfSSatish Balay } 150a7e14dcfSSatish Balay 151a7e14dcfSSatish Balay /* Modify the preconditioner to use the bfgs approximation */ 152a7e14dcfSSatish Balay ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 153a7e14dcfSSatish Balay switch(tl->pc_type) { 154a7e14dcfSSatish Balay case NTL_PC_NONE: 155a7e14dcfSSatish Balay ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr); 1561a1499c8SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 157a7e14dcfSSatish Balay break; 158a7e14dcfSSatish Balay 159a7e14dcfSSatish Balay case NTL_PC_AHESS: 160a7e14dcfSSatish Balay ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr); 1611a1499c8SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 162baa89ecbSBarry Smith ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 163a7e14dcfSSatish Balay break; 164a7e14dcfSSatish Balay 165a7e14dcfSSatish Balay case NTL_PC_BFGS: 166a7e14dcfSSatish Balay ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr); 1671a1499c8SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 168a7e14dcfSSatish Balay ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr); 169a7e14dcfSSatish Balay ierr = PCShellSetContext(pc, tl->M);CHKERRQ(ierr); 170a7e14dcfSSatish Balay ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr); 171a7e14dcfSSatish Balay break; 172a7e14dcfSSatish Balay 173a7e14dcfSSatish Balay default: 174a7e14dcfSSatish Balay /* Use the pc method set by pc_type */ 175a7e14dcfSSatish Balay break; 176a7e14dcfSSatish Balay } 177a7e14dcfSSatish Balay 178a7e14dcfSSatish Balay /* Initialize trust-region radius. The initialization is only performed 179a7e14dcfSSatish Balay when we are using Steihaug-Toint or the Generalized Lanczos method. */ 180a7e14dcfSSatish Balay switch(tl->init_type) { 181a7e14dcfSSatish Balay case NTL_INIT_CONSTANT: 182a7e14dcfSSatish Balay /* Use the initial radius specified */ 183a7e14dcfSSatish Balay break; 184a7e14dcfSSatish Balay 185a7e14dcfSSatish Balay case NTL_INIT_INTERPOLATION: 186a7e14dcfSSatish Balay /* Use the initial radius specified */ 187a7e14dcfSSatish Balay max_radius = 0.0; 188a7e14dcfSSatish Balay 189a7e14dcfSSatish Balay for (j = 0; j < j_max; ++j) { 190a7e14dcfSSatish Balay fmin = f; 191a7e14dcfSSatish Balay sigma = 0.0; 192a7e14dcfSSatish Balay 193a7e14dcfSSatish Balay if (needH) { 194ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 195a7e14dcfSSatish Balay needH = 0; 196a7e14dcfSSatish Balay } 197a7e14dcfSSatish Balay 198a7e14dcfSSatish Balay for (i = 0; i < i_max; ++i) { 199a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 200a7e14dcfSSatish Balay ierr = VecAXPY(tl->W, -tao->trust/gnorm, tao->gradient);CHKERRQ(ierr); 201a7e14dcfSSatish Balay 202a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 203a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 204a7e14dcfSSatish Balay tau = tl->gamma1_i; 20553506e15SBarry Smith } else { 206a7e14dcfSSatish Balay if (ftrial < fmin) { 207a7e14dcfSSatish Balay fmin = ftrial; 208a7e14dcfSSatish Balay sigma = -tao->trust / gnorm; 209a7e14dcfSSatish Balay } 210a7e14dcfSSatish Balay 211a7e14dcfSSatish Balay ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 212a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr); 213a7e14dcfSSatish Balay 214a7e14dcfSSatish Balay prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm)); 215a7e14dcfSSatish Balay actred = f - ftrial; 21653506e15SBarry Smith if ((PetscAbsScalar(actred) <= tl->epsilon) && (PetscAbsScalar(prered) <= tl->epsilon)) { 217a7e14dcfSSatish Balay kappa = 1.0; 21853506e15SBarry Smith } else { 219a7e14dcfSSatish Balay kappa = actred / prered; 220a7e14dcfSSatish Balay } 221a7e14dcfSSatish Balay 222a7e14dcfSSatish Balay tau_1 = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust + (1.0 - tl->theta_i) * prered - actred); 223a7e14dcfSSatish Balay tau_2 = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust - (1.0 + tl->theta_i) * prered + actred); 224a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 225a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 226a7e14dcfSSatish Balay 227a7e14dcfSSatish Balay if (PetscAbsScalar(kappa - 1.0) <= tl->mu1_i) { 228a7e14dcfSSatish Balay /* Great agreement */ 229a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 230a7e14dcfSSatish Balay 231a7e14dcfSSatish Balay if (tau_max < 1.0) { 232a7e14dcfSSatish Balay tau = tl->gamma3_i; 23353506e15SBarry Smith } else if (tau_max > tl->gamma4_i) { 234a7e14dcfSSatish Balay tau = tl->gamma4_i; 23553506e15SBarry Smith } else if (tau_1 >= 1.0 && tau_1 <= tl->gamma4_i && tau_2 < 1.0) { 236a7e14dcfSSatish Balay tau = tau_1; 23753506e15SBarry Smith } else if (tau_2 >= 1.0 && tau_2 <= tl->gamma4_i && tau_1 < 1.0) { 238a7e14dcfSSatish Balay tau = tau_2; 23953506e15SBarry Smith } else { 240a7e14dcfSSatish Balay tau = tau_max; 241a7e14dcfSSatish Balay } 24253506e15SBarry Smith } else if (PetscAbsScalar(kappa - 1.0) <= tl->mu2_i) { 243a7e14dcfSSatish Balay /* Good agreement */ 244a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 245a7e14dcfSSatish Balay 246a7e14dcfSSatish Balay if (tau_max < tl->gamma2_i) { 247a7e14dcfSSatish Balay tau = tl->gamma2_i; 24853506e15SBarry Smith } else if (tau_max > tl->gamma3_i) { 249a7e14dcfSSatish Balay tau = tl->gamma3_i; 25053506e15SBarry Smith } else { 251a7e14dcfSSatish Balay tau = tau_max; 252a7e14dcfSSatish Balay } 25353506e15SBarry Smith } else { 254a7e14dcfSSatish Balay /* Not good agreement */ 255a7e14dcfSSatish Balay if (tau_min > 1.0) { 256a7e14dcfSSatish Balay tau = tl->gamma2_i; 25753506e15SBarry Smith } else if (tau_max < tl->gamma1_i) { 258a7e14dcfSSatish Balay tau = tl->gamma1_i; 25953506e15SBarry Smith } else if ((tau_min < tl->gamma1_i) && (tau_max >= 1.0)) { 260a7e14dcfSSatish Balay tau = tl->gamma1_i; 26153506e15SBarry Smith } else if ((tau_1 >= tl->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < tl->gamma1_i) || (tau_2 >= 1.0))) { 262a7e14dcfSSatish Balay tau = tau_1; 26353506e15SBarry Smith } else if ((tau_2 >= tl->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < tl->gamma1_i) || (tau_2 >= 1.0))) { 264a7e14dcfSSatish Balay tau = tau_2; 26553506e15SBarry Smith } else { 266a7e14dcfSSatish Balay tau = tau_max; 267a7e14dcfSSatish Balay } 268a7e14dcfSSatish Balay } 269a7e14dcfSSatish Balay } 270a7e14dcfSSatish Balay tao->trust = tau * tao->trust; 271a7e14dcfSSatish Balay } 272a7e14dcfSSatish Balay 273a7e14dcfSSatish Balay if (fmin < f) { 274a7e14dcfSSatish Balay f = fmin; 275a7e14dcfSSatish Balay ierr = VecAXPY(tao->solution, sigma, tao->gradient);CHKERRQ(ierr); 276a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr); 277a7e14dcfSSatish Balay 278a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 27953506e15SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 280a7e14dcfSSatish Balay needH = 1; 281a7e14dcfSSatish Balay 2828931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 28353506e15SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 284a7e14dcfSSatish Balay } 285a7e14dcfSSatish Balay } 286a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, max_radius); 287a7e14dcfSSatish Balay 288a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 289a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 290a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 291a7e14dcfSSatish Balay break; 292a7e14dcfSSatish Balay 293a7e14dcfSSatish Balay default: 294a7e14dcfSSatish Balay /* Norm of the first direction will initialize radius */ 295a7e14dcfSSatish Balay tao->trust = 0.0; 296a7e14dcfSSatish Balay break; 297a7e14dcfSSatish Balay } 298a7e14dcfSSatish Balay 299a7e14dcfSSatish Balay /* Set initial scaling for the BFGS preconditioner 300a7e14dcfSSatish Balay This step is done after computing the initial trust-region radius 301a7e14dcfSSatish Balay since the function value may have decreased */ 302a7e14dcfSSatish Balay if (NTL_PC_BFGS == tl->pc_type) { 303a7e14dcfSSatish Balay if (f != 0.0) { 304a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 30553506e15SBarry Smith } else { 306a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 307a7e14dcfSSatish Balay } 308a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 309a7e14dcfSSatish Balay } 310a7e14dcfSSatish Balay 311a7e14dcfSSatish Balay /* Set counter for gradient/reset steps */ 312a7e14dcfSSatish Balay tl->ntrust = 0; 313a7e14dcfSSatish Balay tl->newt = 0; 314a7e14dcfSSatish Balay tl->bfgs = 0; 315a7e14dcfSSatish Balay tl->sgrad = 0; 316a7e14dcfSSatish Balay tl->grad = 0; 317a7e14dcfSSatish Balay 318a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 319a7e14dcfSSatish Balay while (reason == TAO_CONTINUE_ITERATING) { 3208931d482SJason Sarich ++tao->niter; 321ae93cb3cSJason Sarich tao->ksp_its=0; 322a7e14dcfSSatish Balay /* Compute the Hessian */ 323a7e14dcfSSatish Balay if (needH) { 324ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 325a7e14dcfSSatish Balay } 326a7e14dcfSSatish Balay 327a7e14dcfSSatish Balay if (NTL_PC_BFGS == tl->pc_type) { 328a7e14dcfSSatish Balay if (BFGS_SCALE_AHESS == tl->bfgs_scale_type) { 329a7e14dcfSSatish Balay /* Obtain diagonal for the bfgs preconditioner */ 330a7e14dcfSSatish Balay ierr = MatGetDiagonal(tao->hessian, tl->Diag);CHKERRQ(ierr); 331a7e14dcfSSatish Balay ierr = VecAbs(tl->Diag);CHKERRQ(ierr); 332a7e14dcfSSatish Balay ierr = VecReciprocal(tl->Diag);CHKERRQ(ierr); 333a7e14dcfSSatish Balay ierr = MatLMVMSetScale(tl->M, tl->Diag);CHKERRQ(ierr); 334a7e14dcfSSatish Balay } 335a7e14dcfSSatish Balay 336a7e14dcfSSatish Balay /* Update the limited memory preconditioner */ 337a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M,tao->solution, tao->gradient);CHKERRQ(ierr); 338a7e14dcfSSatish Balay ++bfgsUpdates; 339a7e14dcfSSatish Balay } 34023ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre);CHKERRQ(ierr); 341a7e14dcfSSatish Balay /* Solve the Newton system of equations */ 342*ba7fe8fbSTodd Munson ierr = KSPCGSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 343a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 344a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 345a7e14dcfSSatish Balay tao->ksp_its+=its; 346ae93cb3cSJason Sarich tao->ksp_tot_its+=its; 347*ba7fe8fbSTodd Munson ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 348a7e14dcfSSatish Balay 349a7e14dcfSSatish Balay if (0.0 == tao->trust) { 350a7e14dcfSSatish Balay /* Radius was uninitialized; use the norm of the direction */ 351a7e14dcfSSatish Balay if (norm_d > 0.0) { 352a7e14dcfSSatish Balay tao->trust = norm_d; 353a7e14dcfSSatish Balay 354a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 355a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 356a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 35753506e15SBarry Smith } else { 358a7e14dcfSSatish Balay /* The direction was bad; set radius to default value and re-solve 359a7e14dcfSSatish Balay the trust-region subproblem to get a direction */ 360a7e14dcfSSatish Balay tao->trust = tao->trust0; 361a7e14dcfSSatish Balay 362a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 363a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 364a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 365a7e14dcfSSatish Balay 366*ba7fe8fbSTodd Munson ierr = KSPCGSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 367a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 368a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 369a7e14dcfSSatish Balay tao->ksp_its+=its; 3702d9aa51bSJason Sarich tao->ksp_tot_its+=its; 371*ba7fe8fbSTodd Munson ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 372a7e14dcfSSatish Balay 37353506e15SBarry Smith if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 374a7e14dcfSSatish Balay } 375a7e14dcfSSatish Balay } 376a7e14dcfSSatish Balay 377a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 378a7e14dcfSSatish Balay ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 37953506e15SBarry Smith if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (NTL_PC_BFGS == tl->pc_type) && (bfgsUpdates > 1)) { 380a7e14dcfSSatish Balay /* Preconditioner is numerically indefinite; reset the 381a7e14dcfSSatish Balay approximate if using BFGS preconditioning. */ 382a7e14dcfSSatish Balay 383a7e14dcfSSatish Balay if (f != 0.0) { 384a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 38553506e15SBarry Smith } else { 386a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 387a7e14dcfSSatish Balay } 388a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 389a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 390a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 391a7e14dcfSSatish Balay bfgsUpdates = 1; 392a7e14dcfSSatish Balay } 393a7e14dcfSSatish Balay 394a7e14dcfSSatish Balay /* Check trust-region reduction conditions */ 395a7e14dcfSSatish Balay tr_reject = 0; 396a7e14dcfSSatish Balay if (NTL_UPDATE_REDUCTION == tl->update_type) { 397a7e14dcfSSatish Balay /* Get predicted reduction */ 398*ba7fe8fbSTodd Munson ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 399a7e14dcfSSatish Balay if (prered >= 0.0) { 400a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot 401a7e14dcfSSatish Balay happen in infinite precision arithmetic. Step should 402a7e14dcfSSatish Balay be rejected! */ 403a7e14dcfSSatish Balay tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 404a7e14dcfSSatish Balay tr_reject = 1; 40553506e15SBarry Smith } else { 406a7e14dcfSSatish Balay /* Compute trial step and function value */ 407a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 408a7e14dcfSSatish Balay ierr = VecAXPY(tl->W, 1.0, tao->stepdirection);CHKERRQ(ierr); 409a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 410a7e14dcfSSatish Balay 411a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 412a7e14dcfSSatish Balay tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 413a7e14dcfSSatish Balay tr_reject = 1; 41453506e15SBarry Smith } else { 415a7e14dcfSSatish Balay /* Compute and actual reduction */ 416a7e14dcfSSatish Balay actred = f - ftrial; 417a7e14dcfSSatish Balay prered = -prered; 418a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= tl->epsilon) && 419a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= tl->epsilon)) { 420a7e14dcfSSatish Balay kappa = 1.0; 42153506e15SBarry Smith } else { 422a7e14dcfSSatish Balay kappa = actred / prered; 423a7e14dcfSSatish Balay } 424a7e14dcfSSatish Balay 425a7e14dcfSSatish Balay /* Accept of reject the step and update radius */ 426a7e14dcfSSatish Balay if (kappa < tl->eta1) { 427a7e14dcfSSatish Balay /* Reject the step */ 428a7e14dcfSSatish Balay tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 429a7e14dcfSSatish Balay tr_reject = 1; 43053506e15SBarry Smith } else { 431a7e14dcfSSatish Balay /* Accept the step */ 432a7e14dcfSSatish Balay if (kappa < tl->eta2) { 433a7e14dcfSSatish Balay /* Marginal bad step */ 434a7e14dcfSSatish Balay tao->trust = tl->alpha2 * PetscMin(tao->trust, norm_d); 43553506e15SBarry Smith } else if (kappa < tl->eta3) { 436a7e14dcfSSatish Balay /* Reasonable step */ 437a7e14dcfSSatish Balay tao->trust = tl->alpha3 * tao->trust; 43853506e15SBarry Smith } else if (kappa < tl->eta4) { 439a7e14dcfSSatish Balay /* Good step */ 440a7e14dcfSSatish Balay tao->trust = PetscMax(tl->alpha4 * norm_d, tao->trust); 44153506e15SBarry Smith } else { 442a7e14dcfSSatish Balay /* Very good step */ 443a7e14dcfSSatish Balay tao->trust = PetscMax(tl->alpha5 * norm_d, tao->trust); 444a7e14dcfSSatish Balay } 445a7e14dcfSSatish Balay } 446a7e14dcfSSatish Balay } 447a7e14dcfSSatish Balay } 44853506e15SBarry Smith } else { 449a7e14dcfSSatish Balay /* Get predicted reduction */ 450*ba7fe8fbSTodd Munson ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 451a7e14dcfSSatish Balay if (prered >= 0.0) { 452a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot 453a7e14dcfSSatish Balay happen in infinite precision arithmetic. Step should 454a7e14dcfSSatish Balay be rejected! */ 455a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 456a7e14dcfSSatish Balay tr_reject = 1; 45753506e15SBarry Smith } else { 458a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 459a7e14dcfSSatish Balay ierr = VecAXPY(tl->W, 1.0, tao->stepdirection);CHKERRQ(ierr); 460a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 461a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 462a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 463a7e14dcfSSatish Balay tr_reject = 1; 46453506e15SBarry Smith } else { 465a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 466a7e14dcfSSatish Balay 467a7e14dcfSSatish Balay actred = f - ftrial; 468a7e14dcfSSatish Balay prered = -prered; 469a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= tl->epsilon) && 470a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= tl->epsilon)) { 471a7e14dcfSSatish Balay kappa = 1.0; 47253506e15SBarry Smith } else { 473a7e14dcfSSatish Balay kappa = actred / prered; 474a7e14dcfSSatish Balay } 475a7e14dcfSSatish Balay 476a7e14dcfSSatish Balay tau_1 = tl->theta * gdx / (tl->theta * gdx - (1.0 - tl->theta) * prered + actred); 477a7e14dcfSSatish Balay tau_2 = tl->theta * gdx / (tl->theta * gdx + (1.0 + tl->theta) * prered - actred); 478a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 479a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 480a7e14dcfSSatish Balay 481a7e14dcfSSatish Balay if (kappa >= 1.0 - tl->mu1) { 482a7e14dcfSSatish Balay /* Great agreement; accept step and update radius */ 483a7e14dcfSSatish Balay if (tau_max < 1.0) { 484a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d); 48553506e15SBarry Smith } else if (tau_max > tl->gamma4) { 486a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->gamma4 * norm_d); 48753506e15SBarry Smith } else { 488a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 489a7e14dcfSSatish Balay } 49053506e15SBarry Smith } else if (kappa >= 1.0 - tl->mu2) { 491a7e14dcfSSatish Balay /* Good agreement */ 492a7e14dcfSSatish Balay 493a7e14dcfSSatish Balay if (tau_max < tl->gamma2) { 494a7e14dcfSSatish Balay tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d); 49553506e15SBarry Smith } else if (tau_max > tl->gamma3) { 496a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d); 497a7e14dcfSSatish Balay } else if (tau_max < 1.0) { 498a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 49953506e15SBarry Smith } else { 500a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 501a7e14dcfSSatish Balay } 50253506e15SBarry Smith } else { 503a7e14dcfSSatish Balay /* Not good agreement */ 504a7e14dcfSSatish Balay if (tau_min > 1.0) { 505a7e14dcfSSatish Balay tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d); 50653506e15SBarry Smith } else if (tau_max < tl->gamma1) { 507a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 50853506e15SBarry Smith } else if ((tau_min < tl->gamma1) && (tau_max >= 1.0)) { 509a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 51053506e15SBarry Smith } else if ((tau_1 >= tl->gamma1) && (tau_1 < 1.0) && ((tau_2 < tl->gamma1) || (tau_2 >= 1.0))) { 511a7e14dcfSSatish Balay tao->trust = tau_1 * PetscMin(tao->trust, norm_d); 51253506e15SBarry Smith } else if ((tau_2 >= tl->gamma1) && (tau_2 < 1.0) && ((tau_1 < tl->gamma1) || (tau_2 >= 1.0))) { 513a7e14dcfSSatish Balay tao->trust = tau_2 * PetscMin(tao->trust, norm_d); 51453506e15SBarry Smith } else { 515a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 516a7e14dcfSSatish Balay } 517a7e14dcfSSatish Balay tr_reject = 1; 518a7e14dcfSSatish Balay } 519a7e14dcfSSatish Balay } 520a7e14dcfSSatish Balay } 521a7e14dcfSSatish Balay } 522a7e14dcfSSatish Balay 523a7e14dcfSSatish Balay if (tr_reject) { 524a7e14dcfSSatish Balay /* The trust-region constraints rejected the step. Apply a linesearch. 525a7e14dcfSSatish Balay Check for descent direction. */ 526a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 527a7e14dcfSSatish Balay if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 528a7e14dcfSSatish Balay /* Newton step is not descent or direction produced Inf or NaN */ 529a7e14dcfSSatish Balay 530a7e14dcfSSatish Balay if (NTL_PC_BFGS != tl->pc_type) { 531a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and updated 532a7e14dcfSSatish Balay Must use gradient direction in this case */ 533a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 534a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 535a7e14dcfSSatish Balay ++tl->grad; 536a7e14dcfSSatish Balay stepType = NTL_GRADIENT; 53753506e15SBarry Smith } else { 538a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 539a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 540a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 541a7e14dcfSSatish Balay 542a7e14dcfSSatish Balay /* Check for success (descent direction) */ 543a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 544a7e14dcfSSatish Balay if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) { 545a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced not a number 546a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case because 547a7e14dcfSSatish Balay the first solve produces the scaled gradient direction, 548a7e14dcfSSatish Balay which is guaranteed to be descent */ 549a7e14dcfSSatish Balay 550a7e14dcfSSatish Balay /* Use steepest descent direction (scaled) */ 551a7e14dcfSSatish Balay if (f != 0.0) { 552a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 55353506e15SBarry Smith } else { 554a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 555a7e14dcfSSatish Balay } 556a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 557a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 558a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 559a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 560a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 561a7e14dcfSSatish Balay 562a7e14dcfSSatish Balay bfgsUpdates = 1; 563a7e14dcfSSatish Balay ++tl->sgrad; 564a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 56553506e15SBarry Smith } else { 566a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 567a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 568a7e14dcfSSatish Balay ++tl->sgrad; 569a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 57053506e15SBarry Smith } else { 571a7e14dcfSSatish Balay ++tl->bfgs; 572a7e14dcfSSatish Balay stepType = NTL_BFGS; 573a7e14dcfSSatish Balay } 574a7e14dcfSSatish Balay } 575a7e14dcfSSatish Balay } 57653506e15SBarry Smith } else { 577a7e14dcfSSatish Balay /* Computed Newton step is descent */ 578a7e14dcfSSatish Balay ++tl->newt; 579a7e14dcfSSatish Balay stepType = NTL_NEWTON; 580a7e14dcfSSatish Balay } 581a7e14dcfSSatish Balay 582a7e14dcfSSatish Balay /* Perform the linesearch */ 583a7e14dcfSSatish Balay fold = f; 584a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->Xold);CHKERRQ(ierr); 585a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tl->Gold);CHKERRQ(ierr); 586a7e14dcfSSatish Balay 587a7e14dcfSSatish Balay step = 1.0; 588a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 589a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 590a7e14dcfSSatish Balay 59153506e15SBarry Smith while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NTL_GRADIENT) { /* Linesearch failed */ 592a7e14dcfSSatish Balay /* Linesearch failed */ 593a7e14dcfSSatish Balay f = fold; 594a7e14dcfSSatish Balay ierr = VecCopy(tl->Xold, tao->solution);CHKERRQ(ierr); 595a7e14dcfSSatish Balay ierr = VecCopy(tl->Gold, tao->gradient);CHKERRQ(ierr); 596a7e14dcfSSatish Balay 597a7e14dcfSSatish Balay switch(stepType) { 598a7e14dcfSSatish Balay case NTL_NEWTON: 599a7e14dcfSSatish Balay /* Failed to obtain acceptable iterate with Newton step */ 600a7e14dcfSSatish Balay 601a7e14dcfSSatish Balay if (NTL_PC_BFGS != tl->pc_type) { 602a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and being updated 603a7e14dcfSSatish Balay Must use gradient direction in this case */ 604a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 605a7e14dcfSSatish Balay ++tl->grad; 606a7e14dcfSSatish Balay stepType = NTL_GRADIENT; 60753506e15SBarry Smith } else { 608a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 609a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 610a7e14dcfSSatish Balay 611a7e14dcfSSatish Balay 612a7e14dcfSSatish Balay /* Check for success (descent direction) */ 613a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 614a7e14dcfSSatish Balay if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 615a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced 616a7e14dcfSSatish Balay not a number. We can assert bfgsUpdates > 1 in this case 617a7e14dcfSSatish Balay Use steepest descent direction (scaled) */ 618a7e14dcfSSatish Balay 619a7e14dcfSSatish Balay if (f != 0.0) { 620a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 62153506e15SBarry Smith } else { 622a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 623a7e14dcfSSatish Balay } 624a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 625a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 626a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 627a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 628a7e14dcfSSatish Balay 629a7e14dcfSSatish Balay bfgsUpdates = 1; 630a7e14dcfSSatish Balay ++tl->sgrad; 631a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 63253506e15SBarry Smith } else { 633a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 634a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 635a7e14dcfSSatish Balay ++tl->sgrad; 636a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 63753506e15SBarry Smith } else { 638a7e14dcfSSatish Balay ++tl->bfgs; 639a7e14dcfSSatish Balay stepType = NTL_BFGS; 640a7e14dcfSSatish Balay } 641a7e14dcfSSatish Balay } 642a7e14dcfSSatish Balay } 643a7e14dcfSSatish Balay break; 644a7e14dcfSSatish Balay 645a7e14dcfSSatish Balay case NTL_BFGS: 646a7e14dcfSSatish Balay /* Can only enter if pc_type == NTL_PC_BFGS 647a7e14dcfSSatish Balay Failed to obtain acceptable iterate with BFGS step 648a7e14dcfSSatish Balay Attempt to use the scaled gradient direction */ 649a7e14dcfSSatish Balay 650a7e14dcfSSatish Balay if (f != 0.0) { 651a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 65253506e15SBarry Smith } else { 653a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 654a7e14dcfSSatish Balay } 655a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 656a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 657a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 658a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 659a7e14dcfSSatish Balay 660a7e14dcfSSatish Balay bfgsUpdates = 1; 661a7e14dcfSSatish Balay ++tl->sgrad; 662a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 663a7e14dcfSSatish Balay break; 664a7e14dcfSSatish Balay 665a7e14dcfSSatish Balay case NTL_SCALED_GRADIENT: 666a7e14dcfSSatish Balay /* Can only enter if pc_type == NTL_PC_BFGS 667a7e14dcfSSatish Balay The scaled gradient step did not produce a new iterate; 668a7e14dcfSSatish Balay attemp to use the gradient direction. 669a7e14dcfSSatish Balay Need to make sure we are not using a different diagonal scaling */ 670a7e14dcfSSatish Balay ierr = MatLMVMSetScale(tl->M, tl->Diag);CHKERRQ(ierr); 671a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, 1.0);CHKERRQ(ierr); 672a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 673a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 674a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 675a7e14dcfSSatish Balay 676a7e14dcfSSatish Balay bfgsUpdates = 1; 677a7e14dcfSSatish Balay ++tl->grad; 678a7e14dcfSSatish Balay stepType = NTL_GRADIENT; 679a7e14dcfSSatish Balay break; 680a7e14dcfSSatish Balay } 681a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 682a7e14dcfSSatish Balay 683a7e14dcfSSatish Balay /* This may be incorrect; linesearch has values for stepmax and stepmin 684a7e14dcfSSatish Balay that should be reset. */ 685a7e14dcfSSatish Balay step = 1.0; 686a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 687a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 688a7e14dcfSSatish Balay } 689a7e14dcfSSatish Balay 69053506e15SBarry Smith if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 691a7e14dcfSSatish Balay /* Failed to find an improving point */ 692a7e14dcfSSatish Balay f = fold; 693a7e14dcfSSatish Balay ierr = VecCopy(tl->Xold, tao->solution);CHKERRQ(ierr); 694a7e14dcfSSatish Balay ierr = VecCopy(tl->Gold, tao->gradient);CHKERRQ(ierr); 695a7e14dcfSSatish Balay tao->trust = 0.0; 696a7e14dcfSSatish Balay step = 0.0; 697a7e14dcfSSatish Balay reason = TAO_DIVERGED_LS_FAILURE; 698a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 699a7e14dcfSSatish Balay break; 70053506e15SBarry Smith } else if (stepType == NTL_NEWTON) { 701a7e14dcfSSatish Balay if (step < tl->nu1) { 702a7e14dcfSSatish Balay /* Very bad step taken; reduce radius */ 703a7e14dcfSSatish Balay tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust); 70453506e15SBarry Smith } else if (step < tl->nu2) { 705a7e14dcfSSatish Balay /* Reasonably bad step taken; reduce radius */ 706a7e14dcfSSatish Balay tao->trust = tl->omega2 * PetscMin(norm_d, tao->trust); 70753506e15SBarry Smith } else if (step < tl->nu3) { 708a7e14dcfSSatish Balay /* Reasonable step was taken; leave radius alone */ 709a7e14dcfSSatish Balay if (tl->omega3 < 1.0) { 710a7e14dcfSSatish Balay tao->trust = tl->omega3 * PetscMin(norm_d, tao->trust); 71153506e15SBarry Smith } else if (tl->omega3 > 1.0) { 712a7e14dcfSSatish Balay tao->trust = PetscMax(tl->omega3 * norm_d, tao->trust); 713a7e14dcfSSatish Balay } 71453506e15SBarry Smith } else if (step < tl->nu4) { 715a7e14dcfSSatish Balay /* Full step taken; increase the radius */ 716a7e14dcfSSatish Balay tao->trust = PetscMax(tl->omega4 * norm_d, tao->trust); 71753506e15SBarry Smith } else { 718a7e14dcfSSatish Balay /* More than full step taken; increase the radius */ 719a7e14dcfSSatish Balay tao->trust = PetscMax(tl->omega5 * norm_d, tao->trust); 720a7e14dcfSSatish Balay } 72153506e15SBarry Smith } else { 722a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 723a7e14dcfSSatish Balay tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust); 724a7e14dcfSSatish Balay } 72553506e15SBarry Smith } else { 726a7e14dcfSSatish Balay /* Trust-region step is accepted */ 727a7e14dcfSSatish Balay ierr = VecCopy(tl->W, tao->solution);CHKERRQ(ierr); 728a7e14dcfSSatish Balay f = ftrial; 729a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr); 730a7e14dcfSSatish Balay ++tl->ntrust; 731a7e14dcfSSatish Balay } 732a7e14dcfSSatish Balay 733a7e14dcfSSatish Balay /* The radius may have been increased; modify if it is too large */ 734a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 735a7e14dcfSSatish Balay 736e4cb33bbSBarry Smith /* Check for converged */ 737a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 73853506e15SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 739a7e14dcfSSatish Balay needH = 1; 740a7e14dcfSSatish Balay 7418931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, tao->trust, &reason);CHKERRQ(ierr); 742a7e14dcfSSatish Balay } 743a7e14dcfSSatish Balay PetscFunctionReturn(0); 744a7e14dcfSSatish Balay } 745a7e14dcfSSatish Balay 746a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 747a7e14dcfSSatish Balay #undef __FUNCT__ 748a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_NTL" 749441846f8SBarry Smith static PetscErrorCode TaoSetUp_NTL(Tao tao) 750a7e14dcfSSatish Balay { 751a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 752a7e14dcfSSatish Balay PetscErrorCode ierr; 753a7e14dcfSSatish Balay 754a7e14dcfSSatish Balay PetscFunctionBegin; 755a7e14dcfSSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); } 756a7e14dcfSSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);} 757a7e14dcfSSatish Balay if (!tl->W) { ierr = VecDuplicate(tao->solution, &tl->W);CHKERRQ(ierr);} 758a7e14dcfSSatish Balay if (!tl->Xold) { ierr = VecDuplicate(tao->solution, &tl->Xold);CHKERRQ(ierr);} 759a7e14dcfSSatish Balay if (!tl->Gold) { ierr = VecDuplicate(tao->solution, &tl->Gold);CHKERRQ(ierr);} 760a7e14dcfSSatish Balay tl->Diag = 0; 761a7e14dcfSSatish Balay tl->M = 0; 762a7e14dcfSSatish Balay PetscFunctionReturn(0); 763a7e14dcfSSatish Balay } 764a7e14dcfSSatish Balay 765a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 766a7e14dcfSSatish Balay #undef __FUNCT__ 767a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_NTL" 768441846f8SBarry Smith static PetscErrorCode TaoDestroy_NTL(Tao tao) 769a7e14dcfSSatish Balay { 770a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 771a7e14dcfSSatish Balay PetscErrorCode ierr; 772a7e14dcfSSatish Balay 773a7e14dcfSSatish Balay PetscFunctionBegin; 774a7e14dcfSSatish Balay if (tao->setupcalled) { 775a7e14dcfSSatish Balay ierr = VecDestroy(&tl->W);CHKERRQ(ierr); 776a7e14dcfSSatish Balay ierr = VecDestroy(&tl->Xold);CHKERRQ(ierr); 777a7e14dcfSSatish Balay ierr = VecDestroy(&tl->Gold);CHKERRQ(ierr); 778a7e14dcfSSatish Balay } 779a7e14dcfSSatish Balay ierr = VecDestroy(&tl->Diag);CHKERRQ(ierr); 780a7e14dcfSSatish Balay ierr = MatDestroy(&tl->M);CHKERRQ(ierr); 781a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 782a7e14dcfSSatish Balay PetscFunctionReturn(0); 783a7e14dcfSSatish Balay } 784a7e14dcfSSatish Balay 785a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 786a7e14dcfSSatish Balay #undef __FUNCT__ 787a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_NTL" 7884416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_NTL(PetscOptionItems *PetscOptionsObject,Tao tao) 789a7e14dcfSSatish Balay { 790a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 791a7e14dcfSSatish Balay PetscErrorCode ierr; 792a7e14dcfSSatish Balay 793a7e14dcfSSatish Balay PetscFunctionBegin; 7941a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Newton trust region with line search method for unconstrained optimization");CHKERRQ(ierr); 79594ae4db5SBarry Smith ierr = PetscOptionsEList("-tao_ntl_ksp_type", "ksp type", "", NTL_KSP, NTL_KSP_TYPES, NTL_KSP[tl->ksp_type], &tl->ksp_type,NULL);CHKERRQ(ierr); 79694ae4db5SBarry Smith ierr = PetscOptionsEList("-tao_ntl_pc_type", "pc type", "", NTL_PC, NTL_PC_TYPES, NTL_PC[tl->pc_type], &tl->pc_type,NULL);CHKERRQ(ierr); 79794ae4db5SBarry Smith ierr = PetscOptionsEList("-tao_ntl_bfgs_scale_type", "bfgs scale type", "", BFGS_SCALE, BFGS_SCALE_TYPES, BFGS_SCALE[tl->bfgs_scale_type], &tl->bfgs_scale_type,NULL);CHKERRQ(ierr); 79894ae4db5SBarry 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); 79994ae4db5SBarry 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); 80094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_eta1", "poor steplength; reduce radius", "", tl->eta1, &tl->eta1,NULL);CHKERRQ(ierr); 80194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_eta2", "reasonable steplength; leave radius alone", "", tl->eta2, &tl->eta2,NULL);CHKERRQ(ierr); 80294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_eta3", "good steplength; increase radius", "", tl->eta3, &tl->eta3,NULL);CHKERRQ(ierr); 80394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_eta4", "excellent steplength; greatly increase radius", "", tl->eta4, &tl->eta4,NULL);CHKERRQ(ierr); 80494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha1", "", "", tl->alpha1, &tl->alpha1,NULL);CHKERRQ(ierr); 80594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha2", "", "", tl->alpha2, &tl->alpha2,NULL);CHKERRQ(ierr); 80694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha3", "", "", tl->alpha3, &tl->alpha3,NULL);CHKERRQ(ierr); 80794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha4", "", "", tl->alpha4, &tl->alpha4,NULL);CHKERRQ(ierr); 80894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha5", "", "", tl->alpha5, &tl->alpha5,NULL);CHKERRQ(ierr); 80994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_nu1", "poor steplength; reduce radius", "", tl->nu1, &tl->nu1,NULL);CHKERRQ(ierr); 81094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_nu2", "reasonable steplength; leave radius alone", "", tl->nu2, &tl->nu2,NULL);CHKERRQ(ierr); 81194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_nu3", "good steplength; increase radius", "", tl->nu3, &tl->nu3,NULL);CHKERRQ(ierr); 81294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_nu4", "excellent steplength; greatly increase radius", "", tl->nu4, &tl->nu4,NULL);CHKERRQ(ierr); 81394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega1", "", "", tl->omega1, &tl->omega1,NULL);CHKERRQ(ierr); 81494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega2", "", "", tl->omega2, &tl->omega2,NULL);CHKERRQ(ierr); 81594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega3", "", "", tl->omega3, &tl->omega3,NULL);CHKERRQ(ierr); 81694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega4", "", "", tl->omega4, &tl->omega4,NULL);CHKERRQ(ierr); 81794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega5", "", "", tl->omega5, &tl->omega5,NULL);CHKERRQ(ierr); 81894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_mu1_i", "", "", tl->mu1_i, &tl->mu1_i,NULL);CHKERRQ(ierr); 81994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_mu2_i", "", "", tl->mu2_i, &tl->mu2_i,NULL);CHKERRQ(ierr); 82094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma1_i", "", "", tl->gamma1_i, &tl->gamma1_i,NULL);CHKERRQ(ierr); 82194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma2_i", "", "", tl->gamma2_i, &tl->gamma2_i,NULL);CHKERRQ(ierr); 82294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma3_i", "", "", tl->gamma3_i, &tl->gamma3_i,NULL);CHKERRQ(ierr); 82394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma4_i", "", "", tl->gamma4_i, &tl->gamma4_i,NULL);CHKERRQ(ierr); 82494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_theta_i", "", "", tl->theta_i, &tl->theta_i,NULL);CHKERRQ(ierr); 82594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_mu1", "", "", tl->mu1, &tl->mu1,NULL);CHKERRQ(ierr); 82694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_mu2", "", "", tl->mu2, &tl->mu2,NULL);CHKERRQ(ierr); 82794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma1", "", "", tl->gamma1, &tl->gamma1,NULL);CHKERRQ(ierr); 82894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma2", "", "", tl->gamma2, &tl->gamma2,NULL);CHKERRQ(ierr); 82994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma3", "", "", tl->gamma3, &tl->gamma3,NULL);CHKERRQ(ierr); 83094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma4", "", "", tl->gamma4, &tl->gamma4,NULL);CHKERRQ(ierr); 83194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_theta", "", "", tl->theta, &tl->theta,NULL);CHKERRQ(ierr); 83294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_min_radius", "lower bound on initial radius", "", tl->min_radius, &tl->min_radius,NULL);CHKERRQ(ierr); 83394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_max_radius", "upper bound on radius", "", tl->max_radius, &tl->max_radius,NULL);CHKERRQ(ierr); 83494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_epsilon", "tolerance used when computing actual and predicted reduction", "", tl->epsilon, &tl->epsilon,NULL);CHKERRQ(ierr); 835a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 836a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 837a7e14dcfSSatish Balay ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 838a7e14dcfSSatish Balay PetscFunctionReturn(0); 839a7e14dcfSSatish Balay } 840a7e14dcfSSatish Balay 841a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 842a7e14dcfSSatish Balay #undef __FUNCT__ 843a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_NTL" 844441846f8SBarry Smith static PetscErrorCode TaoView_NTL(Tao tao, PetscViewer viewer) 845a7e14dcfSSatish Balay { 846a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 847a7e14dcfSSatish Balay PetscInt nrejects; 848a7e14dcfSSatish Balay PetscBool isascii; 849a7e14dcfSSatish Balay PetscErrorCode ierr; 850a7e14dcfSSatish Balay 851a7e14dcfSSatish Balay PetscFunctionBegin; 852a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 853a7e14dcfSSatish Balay if (isascii) { 854a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 855a7e14dcfSSatish Balay if (NTL_PC_BFGS == tl->pc_type && tl->M) { 856a7e14dcfSSatish Balay ierr = MatLMVMGetRejects(tl->M, &nrejects);CHKERRQ(ierr); 857a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n", nrejects);CHKERRQ(ierr); 858a7e14dcfSSatish Balay } 859a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Trust-region steps: %D\n", tl->ntrust);CHKERRQ(ierr); 860a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Newton search steps: %D\n", tl->newt);CHKERRQ(ierr); 861a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "BFGS search steps: %D\n", tl->bfgs);CHKERRQ(ierr); 862a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient search steps: %D\n", tl->sgrad);CHKERRQ(ierr); 863a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Gradient search steps: %D\n", tl->grad);CHKERRQ(ierr); 864a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 865a7e14dcfSSatish Balay } 866a7e14dcfSSatish Balay PetscFunctionReturn(0); 867a7e14dcfSSatish Balay } 868a7e14dcfSSatish Balay 869a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 8701522df2eSJason Sarich /*MC 8711522df2eSJason Sarich TAONTR - Newton's method with trust region and linesearch 8721522df2eSJason Sarich for unconstrained minimization. 8731522df2eSJason Sarich At each iteration, the Newton trust region method solves the system for d 8741522df2eSJason Sarich and performs a line search in the d direction: 8751522df2eSJason Sarich 8761522df2eSJason Sarich min_d .5 dT Hk d + gkT d, s.t. ||d|| < Delta_k 8771522df2eSJason Sarich 8781522df2eSJason Sarich Options Database Keys: 8791522df2eSJason Sarich + -tao_ntl_ksp_type - "nash","stcg","gltr" 8801522df2eSJason Sarich . -tao_ntl_pc_type - "none","ahess","bfgs","petsc" 8811522df2eSJason Sarich . -tao_ntl_bfgs_scale_type - type of scaling with bfgs pc, "ahess" or "bfgs" 8821522df2eSJason Sarich . -tao_ntl_init_type - "constant","direction","interpolation" 8831522df2eSJason Sarich . -tao_ntl_update_type - "reduction","interpolation" 8841522df2eSJason Sarich . -tao_ntl_min_radius - lower bound on trust region radius 8851522df2eSJason Sarich . -tao_ntl_max_radius - upper bound on trust region radius 8861522df2eSJason Sarich . -tao_ntl_epsilon - tolerance for accepting actual / predicted reduction 8871522df2eSJason Sarich . -tao_ntl_mu1_i - mu1 interpolation init factor 8881522df2eSJason Sarich . -tao_ntl_mu2_i - mu2 interpolation init factor 8891522df2eSJason Sarich . -tao_ntl_gamma1_i - gamma1 interpolation init factor 8901522df2eSJason Sarich . -tao_ntl_gamma2_i - gamma2 interpolation init factor 8911522df2eSJason Sarich . -tao_ntl_gamma3_i - gamma3 interpolation init factor 8921522df2eSJason Sarich . -tao_ntl_gamma4_i - gamma4 interpolation init factor 8931522df2eSJason Sarich . -tao_ntl_theta_i - thetha1 interpolation init factor 8941522df2eSJason Sarich . -tao_ntl_eta1 - eta1 reduction update factor 8951522df2eSJason Sarich . -tao_ntl_eta2 - eta2 reduction update factor 8961522df2eSJason Sarich . -tao_ntl_eta3 - eta3 reduction update factor 8971522df2eSJason Sarich . -tao_ntl_eta4 - eta4 reduction update factor 8981522df2eSJason Sarich . -tao_ntl_alpha1 - alpha1 reduction update factor 8991522df2eSJason Sarich . -tao_ntl_alpha2 - alpha2 reduction update factor 9001522df2eSJason Sarich . -tao_ntl_alpha3 - alpha3 reduction update factor 9011522df2eSJason Sarich . -tao_ntl_alpha4 - alpha4 reduction update factor 9021522df2eSJason Sarich . -tao_ntl_alpha4 - alpha4 reduction update factor 9031522df2eSJason Sarich . -tao_ntl_mu1 - mu1 interpolation update 9041522df2eSJason Sarich . -tao_ntl_mu2 - mu2 interpolation update 9051522df2eSJason Sarich . -tao_ntl_gamma1 - gamma1 interpolcation update 9061522df2eSJason Sarich . -tao_ntl_gamma2 - gamma2 interpolcation update 9071522df2eSJason Sarich . -tao_ntl_gamma3 - gamma3 interpolcation update 9081522df2eSJason Sarich . -tao_ntl_gamma4 - gamma4 interpolation update 9091522df2eSJason Sarich - -tao_ntl_theta - theta1 interpolation update 9101522df2eSJason Sarich 9111eb8069cSJason Sarich Level: beginner 9121522df2eSJason Sarich M*/ 9131522df2eSJason Sarich 914a7e14dcfSSatish Balay #undef __FUNCT__ 915a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_NTL" 916728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NTL(Tao tao) 917a7e14dcfSSatish Balay { 918a7e14dcfSSatish Balay TAO_NTL *tl; 919a7e14dcfSSatish Balay PetscErrorCode ierr; 9208caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 92153506e15SBarry Smith 922a7e14dcfSSatish Balay PetscFunctionBegin; 9233c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&tl);CHKERRQ(ierr); 924a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_NTL; 925a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_NTL; 926a7e14dcfSSatish Balay tao->ops->view = TaoView_NTL; 927a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_NTL; 928a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_NTL; 929a7e14dcfSSatish Balay 9306552cf8aSJason Sarich /* Override default settings (unless already changed) */ 9316552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 50; 9326552cf8aSJason Sarich if (!tao->trust0_changed) tao->trust0 = 100.0; 9336552cf8aSJason Sarich 934a7e14dcfSSatish Balay tao->data = (void*)tl; 935a7e14dcfSSatish Balay 936a7e14dcfSSatish Balay /* Default values for trust-region radius update based on steplength */ 937a7e14dcfSSatish Balay tl->nu1 = 0.25; 938a7e14dcfSSatish Balay tl->nu2 = 0.50; 939a7e14dcfSSatish Balay tl->nu3 = 1.00; 940a7e14dcfSSatish Balay tl->nu4 = 1.25; 941a7e14dcfSSatish Balay 942a7e14dcfSSatish Balay tl->omega1 = 0.25; 943a7e14dcfSSatish Balay tl->omega2 = 0.50; 944a7e14dcfSSatish Balay tl->omega3 = 1.00; 945a7e14dcfSSatish Balay tl->omega4 = 2.00; 946a7e14dcfSSatish Balay tl->omega5 = 4.00; 947a7e14dcfSSatish Balay 948a7e14dcfSSatish Balay /* Default values for trust-region radius update based on reduction */ 949a7e14dcfSSatish Balay tl->eta1 = 1.0e-4; 950a7e14dcfSSatish Balay tl->eta2 = 0.25; 951a7e14dcfSSatish Balay tl->eta3 = 0.50; 952a7e14dcfSSatish Balay tl->eta4 = 0.90; 953a7e14dcfSSatish Balay 954a7e14dcfSSatish Balay tl->alpha1 = 0.25; 955a7e14dcfSSatish Balay tl->alpha2 = 0.50; 956a7e14dcfSSatish Balay tl->alpha3 = 1.00; 957a7e14dcfSSatish Balay tl->alpha4 = 2.00; 958a7e14dcfSSatish Balay tl->alpha5 = 4.00; 959a7e14dcfSSatish Balay 960a7e14dcfSSatish Balay /* Default values for trust-region radius update based on interpolation */ 961a7e14dcfSSatish Balay tl->mu1 = 0.10; 962a7e14dcfSSatish Balay tl->mu2 = 0.50; 963a7e14dcfSSatish Balay 964a7e14dcfSSatish Balay tl->gamma1 = 0.25; 965a7e14dcfSSatish Balay tl->gamma2 = 0.50; 966a7e14dcfSSatish Balay tl->gamma3 = 2.00; 967a7e14dcfSSatish Balay tl->gamma4 = 4.00; 968a7e14dcfSSatish Balay 969a7e14dcfSSatish Balay tl->theta = 0.05; 970a7e14dcfSSatish Balay 971a7e14dcfSSatish Balay /* Default values for trust region initialization based on interpolation */ 972a7e14dcfSSatish Balay tl->mu1_i = 0.35; 973a7e14dcfSSatish Balay tl->mu2_i = 0.50; 974a7e14dcfSSatish Balay 975a7e14dcfSSatish Balay tl->gamma1_i = 0.0625; 976a7e14dcfSSatish Balay tl->gamma2_i = 0.5; 977a7e14dcfSSatish Balay tl->gamma3_i = 2.0; 978a7e14dcfSSatish Balay tl->gamma4_i = 5.0; 979a7e14dcfSSatish Balay 980a7e14dcfSSatish Balay tl->theta_i = 0.25; 981a7e14dcfSSatish Balay 982a7e14dcfSSatish Balay /* Remaining parameters */ 983a7e14dcfSSatish Balay tl->min_radius = 1.0e-10; 984a7e14dcfSSatish Balay tl->max_radius = 1.0e10; 985a7e14dcfSSatish Balay tl->epsilon = 1.0e-6; 986a7e14dcfSSatish Balay 987a7e14dcfSSatish Balay tl->ksp_type = NTL_KSP_STCG; 988a7e14dcfSSatish Balay tl->pc_type = NTL_PC_BFGS; 989a7e14dcfSSatish Balay tl->bfgs_scale_type = BFGS_SCALE_AHESS; 990a7e14dcfSSatish Balay tl->init_type = NTL_INIT_INTERPOLATION; 991a7e14dcfSSatish Balay tl->update_type = NTL_UPDATE_REDUCTION; 992a7e14dcfSSatish Balay 993a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 994a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 995441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch, tao);CHKERRQ(ierr); 9965d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 997a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr); 9985d527766SPatrick Farrell ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr); 999a7e14dcfSSatish Balay PetscFunctionReturn(0); 1000a7e14dcfSSatish Balay } 1001728e0ed0SBarry Smith 1002