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> 5aaa7dc30SBarry Smith #include <petscpc.h> 6af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 7af0996ceSBarry Smith #include <petsc/private/pcimpl.h> 8a7e14dcfSSatish Balay 9a7e14dcfSSatish Balay #define NTL_KSP_NASH 0 10a7e14dcfSSatish Balay #define NTL_KSP_STCG 1 11a7e14dcfSSatish Balay #define NTL_KSP_GLTR 2 12a7e14dcfSSatish Balay #define NTL_KSP_TYPES 3 13a7e14dcfSSatish Balay 14a7e14dcfSSatish Balay #define NTL_PC_NONE 0 15a7e14dcfSSatish Balay #define NTL_PC_AHESS 1 16a7e14dcfSSatish Balay #define NTL_PC_BFGS 2 17a7e14dcfSSatish Balay #define NTL_PC_PETSC 3 18a7e14dcfSSatish Balay #define NTL_PC_TYPES 4 19a7e14dcfSSatish Balay 20a7e14dcfSSatish Balay #define BFGS_SCALE_AHESS 0 21a7e14dcfSSatish Balay #define BFGS_SCALE_BFGS 1 22a7e14dcfSSatish Balay #define BFGS_SCALE_TYPES 2 23a7e14dcfSSatish Balay 24a7e14dcfSSatish Balay #define NTL_INIT_CONSTANT 0 25a7e14dcfSSatish Balay #define NTL_INIT_DIRECTION 1 26a7e14dcfSSatish Balay #define NTL_INIT_INTERPOLATION 2 27a7e14dcfSSatish Balay #define NTL_INIT_TYPES 3 28a7e14dcfSSatish Balay 29a7e14dcfSSatish Balay #define NTL_UPDATE_REDUCTION 0 30a7e14dcfSSatish Balay #define NTL_UPDATE_INTERPOLATION 1 31a7e14dcfSSatish Balay #define NTL_UPDATE_TYPES 2 32a7e14dcfSSatish Balay 3387f595a5SBarry Smith static const char *NTL_KSP[64] = {"nash", "stcg", "gltr"}; 34a7e14dcfSSatish Balay 3587f595a5SBarry Smith static const char *NTL_PC[64] = {"none", "ahess", "bfgs", "petsc"}; 36a7e14dcfSSatish Balay 3787f595a5SBarry Smith static const char *BFGS_SCALE[64] = {"ahess", "bfgs"}; 38a7e14dcfSSatish Balay 3987f595a5SBarry Smith static const char *NTL_INIT[64] = {"constant", "direction", "interpolation"}; 40a7e14dcfSSatish Balay 4187f595a5SBarry Smith static const char *NTL_UPDATE[64] = {"reduction", "interpolation"}; 42a7e14dcfSSatish Balay 43a7e14dcfSSatish Balay /* Routine for BFGS preconditioner */ 44a7e14dcfSSatish Balay 45a7e14dcfSSatish Balay #undef __FUNCT__ 46a7e14dcfSSatish Balay #define __FUNCT__ "MatLMVMSolveShell" 47a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x) 48a7e14dcfSSatish Balay { 49a7e14dcfSSatish Balay PetscErrorCode ierr; 50a7e14dcfSSatish Balay Mat M; 51a7e14dcfSSatish Balay 52a7e14dcfSSatish Balay PetscFunctionBegin; 53a7e14dcfSSatish Balay PetscValidHeaderSpecific(pc,PC_CLASSID,1); 54a7e14dcfSSatish Balay PetscValidHeaderSpecific(b,VEC_CLASSID,2); 55a7e14dcfSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,3); 56a7e14dcfSSatish Balay ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr); 57a7e14dcfSSatish Balay ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr); 58a7e14dcfSSatish Balay PetscFunctionReturn(0); 59a7e14dcfSSatish Balay } 60a7e14dcfSSatish Balay 61a7e14dcfSSatish Balay /* Implements Newton's Method with a trust-region, line-search approach for 62a7e14dcfSSatish Balay solving unconstrained minimization problems. A More'-Thuente line search 63a7e14dcfSSatish Balay is used to guarantee that the bfgs preconditioner remains positive 64a7e14dcfSSatish Balay definite. */ 65a7e14dcfSSatish Balay 66a7e14dcfSSatish Balay #define NTL_NEWTON 0 67a7e14dcfSSatish Balay #define NTL_BFGS 1 68a7e14dcfSSatish Balay #define NTL_SCALED_GRADIENT 2 69a7e14dcfSSatish Balay #define NTL_GRADIENT 3 70a7e14dcfSSatish Balay 71a7e14dcfSSatish Balay #undef __FUNCT__ 72a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_NTL" 73441846f8SBarry Smith static PetscErrorCode TaoSolve_NTL(Tao tao) 74a7e14dcfSSatish Balay { 75a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 76a7e14dcfSSatish Balay PC pc; 77a7e14dcfSSatish Balay KSPConvergedReason ksp_reason; 78e4cb33bbSBarry Smith TaoConvergedReason reason; 79e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 80a7e14dcfSSatish Balay 81a7e14dcfSSatish Balay PetscReal fmin, ftrial, prered, actred, kappa, sigma; 82a7e14dcfSSatish Balay PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 83a7e14dcfSSatish Balay PetscReal f, fold, gdx, gnorm; 84a7e14dcfSSatish Balay PetscReal step = 1.0; 85a7e14dcfSSatish Balay 86a7e14dcfSSatish Balay PetscReal delta; 87a7e14dcfSSatish Balay PetscReal norm_d = 0.0; 88a7e14dcfSSatish Balay PetscErrorCode ierr; 89a7e14dcfSSatish Balay PetscInt stepType; 908931d482SJason Sarich PetscInt its; 91a7e14dcfSSatish Balay 92a7e14dcfSSatish Balay PetscInt bfgsUpdates = 0; 93a7e14dcfSSatish Balay PetscInt needH; 94a7e14dcfSSatish Balay 95a7e14dcfSSatish Balay PetscInt i_max = 5; 96a7e14dcfSSatish Balay PetscInt j_max = 1; 97a7e14dcfSSatish Balay PetscInt i, j, n, N; 98a7e14dcfSSatish Balay 99a7e14dcfSSatish Balay PetscInt tr_reject; 100a7e14dcfSSatish Balay 101a7e14dcfSSatish Balay PetscFunctionBegin; 102a7e14dcfSSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 103a7e14dcfSSatish Balay ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by ntl algorithm\n");CHKERRQ(ierr); 104a7e14dcfSSatish Balay } 105a7e14dcfSSatish Balay 106a7e14dcfSSatish Balay /* Initialize trust-region radius */ 107a7e14dcfSSatish Balay tao->trust = tao->trust0; 108a7e14dcfSSatish Balay 109a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 110a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 111a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 112a7e14dcfSSatish Balay 113a7e14dcfSSatish Balay if (NTL_PC_BFGS == tl->pc_type && !tl->M) { 114a7e14dcfSSatish Balay ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 115a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 116a7e14dcfSSatish Balay ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&tl->M);CHKERRQ(ierr); 117a7e14dcfSSatish Balay ierr = MatLMVMAllocateVectors(tl->M,tao->solution);CHKERRQ(ierr); 118a7e14dcfSSatish Balay } 119a7e14dcfSSatish Balay 120a7e14dcfSSatish Balay /* Check convergence criteria */ 121a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 122a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 12353506e15SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 124a7e14dcfSSatish Balay needH = 1; 125a7e14dcfSSatish Balay 1268931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 12753506e15SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 128a7e14dcfSSatish Balay 129a7e14dcfSSatish Balay /* Create vectors for the limited memory preconditioner */ 13053506e15SBarry Smith if ((NTL_PC_BFGS == tl->pc_type) && (BFGS_SCALE_BFGS != tl->bfgs_scale_type)) { 131a7e14dcfSSatish Balay if (!tl->Diag) { 132a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tl->Diag);CHKERRQ(ierr); 133a7e14dcfSSatish Balay } 134a7e14dcfSSatish Balay } 135a7e14dcfSSatish Balay 136a7e14dcfSSatish Balay /* Modify the linear solver to a conjugate gradient method */ 137a7e14dcfSSatish Balay switch(tl->ksp_type) { 138a7e14dcfSSatish Balay case NTL_KSP_NASH: 139a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPNASH);CHKERRQ(ierr); 1401a1499c8SBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 141a7e14dcfSSatish Balay break; 142a7e14dcfSSatish Balay 143a7e14dcfSSatish Balay case NTL_KSP_STCG: 144a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPSTCG);CHKERRQ(ierr); 1451a1499c8SBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 146a7e14dcfSSatish Balay break; 147a7e14dcfSSatish Balay 148a7e14dcfSSatish Balay default: 149a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPGLTR);CHKERRQ(ierr); 1501a1499c8SBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 151a7e14dcfSSatish Balay break; 152a7e14dcfSSatish Balay } 153a7e14dcfSSatish Balay 154a7e14dcfSSatish Balay /* Modify the preconditioner to use the bfgs approximation */ 155a7e14dcfSSatish Balay ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 156a7e14dcfSSatish Balay switch(tl->pc_type) { 157a7e14dcfSSatish Balay case NTL_PC_NONE: 158a7e14dcfSSatish Balay ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr); 1591a1499c8SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 160a7e14dcfSSatish Balay break; 161a7e14dcfSSatish Balay 162a7e14dcfSSatish Balay case NTL_PC_AHESS: 163a7e14dcfSSatish Balay ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr); 1641a1499c8SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 165baa89ecbSBarry Smith ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 166a7e14dcfSSatish Balay break; 167a7e14dcfSSatish Balay 168a7e14dcfSSatish Balay case NTL_PC_BFGS: 169a7e14dcfSSatish Balay ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr); 1701a1499c8SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 171a7e14dcfSSatish Balay ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr); 172a7e14dcfSSatish Balay ierr = PCShellSetContext(pc, tl->M);CHKERRQ(ierr); 173a7e14dcfSSatish Balay ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr); 174a7e14dcfSSatish Balay break; 175a7e14dcfSSatish Balay 176a7e14dcfSSatish Balay default: 177a7e14dcfSSatish Balay /* Use the pc method set by pc_type */ 178a7e14dcfSSatish Balay break; 179a7e14dcfSSatish Balay } 180a7e14dcfSSatish Balay 181a7e14dcfSSatish Balay /* Initialize trust-region radius. The initialization is only performed 182a7e14dcfSSatish Balay when we are using Steihaug-Toint or the Generalized Lanczos method. */ 183a7e14dcfSSatish Balay switch(tl->init_type) { 184a7e14dcfSSatish Balay case NTL_INIT_CONSTANT: 185a7e14dcfSSatish Balay /* Use the initial radius specified */ 186a7e14dcfSSatish Balay break; 187a7e14dcfSSatish Balay 188a7e14dcfSSatish Balay case NTL_INIT_INTERPOLATION: 189a7e14dcfSSatish Balay /* Use the initial radius specified */ 190a7e14dcfSSatish Balay max_radius = 0.0; 191a7e14dcfSSatish Balay 192a7e14dcfSSatish Balay for (j = 0; j < j_max; ++j) { 193a7e14dcfSSatish Balay fmin = f; 194a7e14dcfSSatish Balay sigma = 0.0; 195a7e14dcfSSatish Balay 196a7e14dcfSSatish Balay if (needH) { 197ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 198a7e14dcfSSatish Balay needH = 0; 199a7e14dcfSSatish Balay } 200a7e14dcfSSatish Balay 201a7e14dcfSSatish Balay for (i = 0; i < i_max; ++i) { 202a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 203a7e14dcfSSatish Balay ierr = VecAXPY(tl->W, -tao->trust/gnorm, tao->gradient);CHKERRQ(ierr); 204a7e14dcfSSatish Balay 205a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 206a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 207a7e14dcfSSatish Balay tau = tl->gamma1_i; 20853506e15SBarry Smith } else { 209a7e14dcfSSatish Balay if (ftrial < fmin) { 210a7e14dcfSSatish Balay fmin = ftrial; 211a7e14dcfSSatish Balay sigma = -tao->trust / gnorm; 212a7e14dcfSSatish Balay } 213a7e14dcfSSatish Balay 214a7e14dcfSSatish Balay ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 215a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr); 216a7e14dcfSSatish Balay 217a7e14dcfSSatish Balay prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm)); 218a7e14dcfSSatish Balay actred = f - ftrial; 21953506e15SBarry Smith if ((PetscAbsScalar(actred) <= tl->epsilon) && (PetscAbsScalar(prered) <= tl->epsilon)) { 220a7e14dcfSSatish Balay kappa = 1.0; 22153506e15SBarry Smith } else { 222a7e14dcfSSatish Balay kappa = actred / prered; 223a7e14dcfSSatish Balay } 224a7e14dcfSSatish Balay 225a7e14dcfSSatish Balay tau_1 = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust + (1.0 - tl->theta_i) * prered - actred); 226a7e14dcfSSatish Balay tau_2 = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust - (1.0 + tl->theta_i) * prered + actred); 227a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 228a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 229a7e14dcfSSatish Balay 230a7e14dcfSSatish Balay if (PetscAbsScalar(kappa - 1.0) <= tl->mu1_i) { 231a7e14dcfSSatish Balay /* Great agreement */ 232a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 233a7e14dcfSSatish Balay 234a7e14dcfSSatish Balay if (tau_max < 1.0) { 235a7e14dcfSSatish Balay tau = tl->gamma3_i; 23653506e15SBarry Smith } else if (tau_max > tl->gamma4_i) { 237a7e14dcfSSatish Balay tau = tl->gamma4_i; 23853506e15SBarry Smith } else if (tau_1 >= 1.0 && tau_1 <= tl->gamma4_i && tau_2 < 1.0) { 239a7e14dcfSSatish Balay tau = tau_1; 24053506e15SBarry Smith } else if (tau_2 >= 1.0 && tau_2 <= tl->gamma4_i && tau_1 < 1.0) { 241a7e14dcfSSatish Balay tau = tau_2; 24253506e15SBarry Smith } else { 243a7e14dcfSSatish Balay tau = tau_max; 244a7e14dcfSSatish Balay } 24553506e15SBarry Smith } else if (PetscAbsScalar(kappa - 1.0) <= tl->mu2_i) { 246a7e14dcfSSatish Balay /* Good agreement */ 247a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 248a7e14dcfSSatish Balay 249a7e14dcfSSatish Balay if (tau_max < tl->gamma2_i) { 250a7e14dcfSSatish Balay tau = tl->gamma2_i; 25153506e15SBarry Smith } else if (tau_max > tl->gamma3_i) { 252a7e14dcfSSatish Balay tau = tl->gamma3_i; 25353506e15SBarry Smith } else { 254a7e14dcfSSatish Balay tau = tau_max; 255a7e14dcfSSatish Balay } 25653506e15SBarry Smith } else { 257a7e14dcfSSatish Balay /* Not good agreement */ 258a7e14dcfSSatish Balay if (tau_min > 1.0) { 259a7e14dcfSSatish Balay tau = tl->gamma2_i; 26053506e15SBarry Smith } else if (tau_max < tl->gamma1_i) { 261a7e14dcfSSatish Balay tau = tl->gamma1_i; 26253506e15SBarry Smith } else if ((tau_min < tl->gamma1_i) && (tau_max >= 1.0)) { 263a7e14dcfSSatish Balay tau = tl->gamma1_i; 26453506e15SBarry Smith } else if ((tau_1 >= tl->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < tl->gamma1_i) || (tau_2 >= 1.0))) { 265a7e14dcfSSatish Balay tau = tau_1; 26653506e15SBarry Smith } else if ((tau_2 >= tl->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < tl->gamma1_i) || (tau_2 >= 1.0))) { 267a7e14dcfSSatish Balay tau = tau_2; 26853506e15SBarry Smith } else { 269a7e14dcfSSatish Balay tau = tau_max; 270a7e14dcfSSatish Balay } 271a7e14dcfSSatish Balay } 272a7e14dcfSSatish Balay } 273a7e14dcfSSatish Balay tao->trust = tau * tao->trust; 274a7e14dcfSSatish Balay } 275a7e14dcfSSatish Balay 276a7e14dcfSSatish Balay if (fmin < f) { 277a7e14dcfSSatish Balay f = fmin; 278a7e14dcfSSatish Balay ierr = VecAXPY(tao->solution, sigma, tao->gradient);CHKERRQ(ierr); 279a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr); 280a7e14dcfSSatish Balay 281a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 28253506e15SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 283a7e14dcfSSatish Balay needH = 1; 284a7e14dcfSSatish Balay 2858931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 28653506e15SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 287a7e14dcfSSatish Balay } 288a7e14dcfSSatish Balay } 289a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, max_radius); 290a7e14dcfSSatish Balay 291a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 292a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 293a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 294a7e14dcfSSatish Balay break; 295a7e14dcfSSatish Balay 296a7e14dcfSSatish Balay default: 297a7e14dcfSSatish Balay /* Norm of the first direction will initialize radius */ 298a7e14dcfSSatish Balay tao->trust = 0.0; 299a7e14dcfSSatish Balay break; 300a7e14dcfSSatish Balay } 301a7e14dcfSSatish Balay 302a7e14dcfSSatish Balay /* Set initial scaling for the BFGS preconditioner 303a7e14dcfSSatish Balay This step is done after computing the initial trust-region radius 304a7e14dcfSSatish Balay since the function value may have decreased */ 305a7e14dcfSSatish Balay if (NTL_PC_BFGS == tl->pc_type) { 306a7e14dcfSSatish Balay if (f != 0.0) { 307a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 30853506e15SBarry Smith } else { 309a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 310a7e14dcfSSatish Balay } 311a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 312a7e14dcfSSatish Balay } 313a7e14dcfSSatish Balay 314a7e14dcfSSatish Balay /* Set counter for gradient/reset steps */ 315a7e14dcfSSatish Balay tl->ntrust = 0; 316a7e14dcfSSatish Balay tl->newt = 0; 317a7e14dcfSSatish Balay tl->bfgs = 0; 318a7e14dcfSSatish Balay tl->sgrad = 0; 319a7e14dcfSSatish Balay tl->grad = 0; 320a7e14dcfSSatish Balay 321a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 322a7e14dcfSSatish Balay while (reason == TAO_CONTINUE_ITERATING) { 3238931d482SJason Sarich ++tao->niter; 324ae93cb3cSJason Sarich tao->ksp_its=0; 325a7e14dcfSSatish Balay /* Compute the Hessian */ 326a7e14dcfSSatish Balay if (needH) { 327ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 328a7e14dcfSSatish Balay needH = 0; 329a7e14dcfSSatish Balay } 330a7e14dcfSSatish Balay 331a7e14dcfSSatish Balay if (NTL_PC_BFGS == tl->pc_type) { 332a7e14dcfSSatish Balay if (BFGS_SCALE_AHESS == tl->bfgs_scale_type) { 333a7e14dcfSSatish Balay /* Obtain diagonal for the bfgs preconditioner */ 334a7e14dcfSSatish Balay ierr = MatGetDiagonal(tao->hessian, tl->Diag);CHKERRQ(ierr); 335a7e14dcfSSatish Balay ierr = VecAbs(tl->Diag);CHKERRQ(ierr); 336a7e14dcfSSatish Balay ierr = VecReciprocal(tl->Diag);CHKERRQ(ierr); 337a7e14dcfSSatish Balay ierr = MatLMVMSetScale(tl->M, tl->Diag);CHKERRQ(ierr); 338a7e14dcfSSatish Balay } 339a7e14dcfSSatish Balay 340a7e14dcfSSatish Balay /* Update the limited memory preconditioner */ 341a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M,tao->solution, tao->gradient);CHKERRQ(ierr); 342a7e14dcfSSatish Balay ++bfgsUpdates; 343a7e14dcfSSatish Balay } 34423ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre);CHKERRQ(ierr); 345a7e14dcfSSatish Balay /* Solve the Newton system of equations */ 346a7e14dcfSSatish Balay if (NTL_KSP_NASH == tl->ksp_type) { 347a7e14dcfSSatish Balay ierr = KSPNASHSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 348a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 349a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 350a7e14dcfSSatish Balay tao->ksp_its+=its; 351ae93cb3cSJason Sarich tao->ksp_tot_its+=its; 352a7e14dcfSSatish Balay ierr = KSPNASHGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 353a7e14dcfSSatish Balay } else if (NTL_KSP_STCG == tl->ksp_type) { 354a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 355a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 356a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 357a7e14dcfSSatish Balay tao->ksp_its+=its; 358ae93cb3cSJason Sarich tao->ksp_tot_its+=its; 359a7e14dcfSSatish Balay ierr = KSPSTCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 360a7e14dcfSSatish Balay } else { /* NTL_KSP_GLTR */ 361a7e14dcfSSatish Balay ierr = KSPGLTRSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 362a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 363a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 364a7e14dcfSSatish Balay tao->ksp_its+=its; 3652d9aa51bSJason Sarich tao->ksp_tot_its+=its; 366a7e14dcfSSatish Balay ierr = KSPGLTRGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 367a7e14dcfSSatish Balay } 368a7e14dcfSSatish Balay 369a7e14dcfSSatish Balay if (0.0 == tao->trust) { 370a7e14dcfSSatish Balay /* Radius was uninitialized; use the norm of the direction */ 371a7e14dcfSSatish Balay if (norm_d > 0.0) { 372a7e14dcfSSatish Balay tao->trust = norm_d; 373a7e14dcfSSatish Balay 374a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 375a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 376a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 37753506e15SBarry Smith } else { 378a7e14dcfSSatish Balay /* The direction was bad; set radius to default value and re-solve 379a7e14dcfSSatish Balay the trust-region subproblem to get a direction */ 380a7e14dcfSSatish Balay tao->trust = tao->trust0; 381a7e14dcfSSatish Balay 382a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 383a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 384a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 385a7e14dcfSSatish Balay 386a7e14dcfSSatish Balay if (NTL_KSP_NASH == tl->ksp_type) { 387a7e14dcfSSatish Balay ierr = KSPNASHSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 388a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 389a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 390a7e14dcfSSatish Balay tao->ksp_its+=its; 3912d9aa51bSJason Sarich tao->ksp_tot_its+=its; 392a7e14dcfSSatish Balay ierr = KSPNASHGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 393a7e14dcfSSatish Balay } else if (NTL_KSP_STCG == tl->ksp_type) { 394a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 395a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 396a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 397a7e14dcfSSatish Balay tao->ksp_its+=its; 3982d9aa51bSJason Sarich tao->ksp_tot_its+=its; 399a7e14dcfSSatish Balay ierr = KSPSTCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 400a7e14dcfSSatish Balay } else { /* NTL_KSP_GLTR */ 401a7e14dcfSSatish Balay ierr = KSPGLTRSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 402a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 403a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 404a7e14dcfSSatish Balay tao->ksp_its+=its; 4052d9aa51bSJason Sarich tao->ksp_tot_its+=its; 406a7e14dcfSSatish Balay ierr = KSPGLTRGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 407a7e14dcfSSatish Balay } 408a7e14dcfSSatish Balay 409a7e14dcfSSatish Balay 41053506e15SBarry Smith if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 411a7e14dcfSSatish Balay } 412a7e14dcfSSatish Balay } 413a7e14dcfSSatish Balay 414a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 415a7e14dcfSSatish Balay ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 41653506e15SBarry Smith if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (NTL_PC_BFGS == tl->pc_type) && (bfgsUpdates > 1)) { 417a7e14dcfSSatish Balay /* Preconditioner is numerically indefinite; reset the 418a7e14dcfSSatish Balay approximate if using BFGS preconditioning. */ 419a7e14dcfSSatish Balay 420a7e14dcfSSatish Balay if (f != 0.0) { 421a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 42253506e15SBarry Smith } else { 423a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 424a7e14dcfSSatish Balay } 425a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 426a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 427a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 428a7e14dcfSSatish Balay bfgsUpdates = 1; 429a7e14dcfSSatish Balay } 430a7e14dcfSSatish Balay 431a7e14dcfSSatish Balay /* Check trust-region reduction conditions */ 432a7e14dcfSSatish Balay tr_reject = 0; 433a7e14dcfSSatish Balay if (NTL_UPDATE_REDUCTION == tl->update_type) { 434a7e14dcfSSatish Balay /* Get predicted reduction */ 435a7e14dcfSSatish Balay if (NTL_KSP_NASH == tl->ksp_type) { 436a7e14dcfSSatish Balay ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 437a7e14dcfSSatish Balay } else if (NTL_KSP_STCG == tl->ksp_type) { 438a7e14dcfSSatish Balay ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 439a7e14dcfSSatish Balay } else { /* gltr */ 440a7e14dcfSSatish Balay ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 441a7e14dcfSSatish Balay } 442a7e14dcfSSatish Balay 443a7e14dcfSSatish Balay if (prered >= 0.0) { 444a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot 445a7e14dcfSSatish Balay happen in infinite precision arithmetic. Step should 446a7e14dcfSSatish Balay be rejected! */ 447a7e14dcfSSatish Balay tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 448a7e14dcfSSatish Balay tr_reject = 1; 44953506e15SBarry Smith } else { 450a7e14dcfSSatish Balay /* Compute trial step and function value */ 451a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 452a7e14dcfSSatish Balay ierr = VecAXPY(tl->W, 1.0, tao->stepdirection);CHKERRQ(ierr); 453a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 454a7e14dcfSSatish Balay 455a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 456a7e14dcfSSatish Balay tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 457a7e14dcfSSatish Balay tr_reject = 1; 45853506e15SBarry Smith } else { 459a7e14dcfSSatish Balay /* Compute and actual reduction */ 460a7e14dcfSSatish Balay actred = f - ftrial; 461a7e14dcfSSatish Balay prered = -prered; 462a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= tl->epsilon) && 463a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= tl->epsilon)) { 464a7e14dcfSSatish Balay kappa = 1.0; 46553506e15SBarry Smith } else { 466a7e14dcfSSatish Balay kappa = actred / prered; 467a7e14dcfSSatish Balay } 468a7e14dcfSSatish Balay 469a7e14dcfSSatish Balay /* Accept of reject the step and update radius */ 470a7e14dcfSSatish Balay if (kappa < tl->eta1) { 471a7e14dcfSSatish Balay /* Reject the step */ 472a7e14dcfSSatish Balay tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 473a7e14dcfSSatish Balay tr_reject = 1; 47453506e15SBarry Smith } else { 475a7e14dcfSSatish Balay /* Accept the step */ 476a7e14dcfSSatish Balay if (kappa < tl->eta2) { 477a7e14dcfSSatish Balay /* Marginal bad step */ 478a7e14dcfSSatish Balay tao->trust = tl->alpha2 * PetscMin(tao->trust, norm_d); 47953506e15SBarry Smith } else if (kappa < tl->eta3) { 480a7e14dcfSSatish Balay /* Reasonable step */ 481a7e14dcfSSatish Balay tao->trust = tl->alpha3 * tao->trust; 48253506e15SBarry Smith } else if (kappa < tl->eta4) { 483a7e14dcfSSatish Balay /* Good step */ 484a7e14dcfSSatish Balay tao->trust = PetscMax(tl->alpha4 * norm_d, tao->trust); 48553506e15SBarry Smith } else { 486a7e14dcfSSatish Balay /* Very good step */ 487a7e14dcfSSatish Balay tao->trust = PetscMax(tl->alpha5 * norm_d, tao->trust); 488a7e14dcfSSatish Balay } 489a7e14dcfSSatish Balay } 490a7e14dcfSSatish Balay } 491a7e14dcfSSatish Balay } 49253506e15SBarry Smith } else { 493a7e14dcfSSatish Balay /* Get predicted reduction */ 494a7e14dcfSSatish Balay if (NTL_KSP_NASH == tl->ksp_type) { 495a7e14dcfSSatish Balay ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 496a7e14dcfSSatish Balay } else if (NTL_KSP_STCG == tl->ksp_type) { 497a7e14dcfSSatish Balay ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 498a7e14dcfSSatish Balay } else { /* gltr */ 499a7e14dcfSSatish Balay ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 500a7e14dcfSSatish Balay } 501a7e14dcfSSatish Balay 502a7e14dcfSSatish Balay if (prered >= 0.0) { 503a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot 504a7e14dcfSSatish Balay happen in infinite precision arithmetic. Step should 505a7e14dcfSSatish Balay be rejected! */ 506a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 507a7e14dcfSSatish Balay tr_reject = 1; 50853506e15SBarry Smith } else { 509a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 510a7e14dcfSSatish Balay ierr = VecAXPY(tl->W, 1.0, tao->stepdirection);CHKERRQ(ierr); 511a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 512a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 513a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 514a7e14dcfSSatish Balay tr_reject = 1; 51553506e15SBarry Smith } else { 516a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 517a7e14dcfSSatish Balay 518a7e14dcfSSatish Balay actred = f - ftrial; 519a7e14dcfSSatish Balay prered = -prered; 520a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= tl->epsilon) && 521a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= tl->epsilon)) { 522a7e14dcfSSatish Balay kappa = 1.0; 52353506e15SBarry Smith } else { 524a7e14dcfSSatish Balay kappa = actred / prered; 525a7e14dcfSSatish Balay } 526a7e14dcfSSatish Balay 527a7e14dcfSSatish Balay tau_1 = tl->theta * gdx / (tl->theta * gdx - (1.0 - tl->theta) * prered + actred); 528a7e14dcfSSatish Balay tau_2 = tl->theta * gdx / (tl->theta * gdx + (1.0 + tl->theta) * prered - actred); 529a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 530a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 531a7e14dcfSSatish Balay 532a7e14dcfSSatish Balay if (kappa >= 1.0 - tl->mu1) { 533a7e14dcfSSatish Balay /* Great agreement; accept step and update radius */ 534a7e14dcfSSatish Balay if (tau_max < 1.0) { 535a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d); 53653506e15SBarry Smith } else if (tau_max > tl->gamma4) { 537a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->gamma4 * norm_d); 53853506e15SBarry Smith } else { 539a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 540a7e14dcfSSatish Balay } 54153506e15SBarry Smith } else if (kappa >= 1.0 - tl->mu2) { 542a7e14dcfSSatish Balay /* Good agreement */ 543a7e14dcfSSatish Balay 544a7e14dcfSSatish Balay if (tau_max < tl->gamma2) { 545a7e14dcfSSatish Balay tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d); 54653506e15SBarry Smith } else if (tau_max > tl->gamma3) { 547a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d); 548a7e14dcfSSatish Balay } else if (tau_max < 1.0) { 549a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 55053506e15SBarry Smith } else { 551a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 552a7e14dcfSSatish Balay } 55353506e15SBarry Smith } else { 554a7e14dcfSSatish Balay /* Not good agreement */ 555a7e14dcfSSatish Balay if (tau_min > 1.0) { 556a7e14dcfSSatish Balay tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d); 55753506e15SBarry Smith } else if (tau_max < tl->gamma1) { 558a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 55953506e15SBarry Smith } else if ((tau_min < tl->gamma1) && (tau_max >= 1.0)) { 560a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 56153506e15SBarry Smith } else if ((tau_1 >= tl->gamma1) && (tau_1 < 1.0) && ((tau_2 < tl->gamma1) || (tau_2 >= 1.0))) { 562a7e14dcfSSatish Balay tao->trust = tau_1 * PetscMin(tao->trust, norm_d); 56353506e15SBarry Smith } else if ((tau_2 >= tl->gamma1) && (tau_2 < 1.0) && ((tau_1 < tl->gamma1) || (tau_2 >= 1.0))) { 564a7e14dcfSSatish Balay tao->trust = tau_2 * PetscMin(tao->trust, norm_d); 56553506e15SBarry Smith } else { 566a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 567a7e14dcfSSatish Balay } 568a7e14dcfSSatish Balay tr_reject = 1; 569a7e14dcfSSatish Balay } 570a7e14dcfSSatish Balay } 571a7e14dcfSSatish Balay } 572a7e14dcfSSatish Balay } 573a7e14dcfSSatish Balay 574a7e14dcfSSatish Balay if (tr_reject) { 575a7e14dcfSSatish Balay /* The trust-region constraints rejected the step. Apply a linesearch. 576a7e14dcfSSatish Balay Check for descent direction. */ 577a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 578a7e14dcfSSatish Balay if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 579a7e14dcfSSatish Balay /* Newton step is not descent or direction produced Inf or NaN */ 580a7e14dcfSSatish Balay 581a7e14dcfSSatish Balay if (NTL_PC_BFGS != tl->pc_type) { 582a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and updated 583a7e14dcfSSatish Balay Must use gradient direction in this case */ 584a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 585a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 586a7e14dcfSSatish Balay ++tl->grad; 587a7e14dcfSSatish Balay stepType = NTL_GRADIENT; 58853506e15SBarry Smith } else { 589a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 590a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 591a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 592a7e14dcfSSatish Balay 593a7e14dcfSSatish Balay /* Check for success (descent direction) */ 594a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 595a7e14dcfSSatish Balay if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) { 596a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced not a number 597a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case because 598a7e14dcfSSatish Balay the first solve produces the scaled gradient direction, 599a7e14dcfSSatish Balay which is guaranteed to be descent */ 600a7e14dcfSSatish Balay 601a7e14dcfSSatish Balay /* Use steepest descent direction (scaled) */ 602a7e14dcfSSatish Balay if (f != 0.0) { 603a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 60453506e15SBarry Smith } else { 605a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 606a7e14dcfSSatish Balay } 607a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 608a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 609a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 610a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 611a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 612a7e14dcfSSatish Balay 613a7e14dcfSSatish Balay bfgsUpdates = 1; 614a7e14dcfSSatish Balay ++tl->sgrad; 615a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 61653506e15SBarry Smith } else { 617a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 618a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 619a7e14dcfSSatish Balay ++tl->sgrad; 620a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 62153506e15SBarry Smith } else { 622a7e14dcfSSatish Balay ++tl->bfgs; 623a7e14dcfSSatish Balay stepType = NTL_BFGS; 624a7e14dcfSSatish Balay } 625a7e14dcfSSatish Balay } 626a7e14dcfSSatish Balay } 62753506e15SBarry Smith } else { 628a7e14dcfSSatish Balay /* Computed Newton step is descent */ 629a7e14dcfSSatish Balay ++tl->newt; 630a7e14dcfSSatish Balay stepType = NTL_NEWTON; 631a7e14dcfSSatish Balay } 632a7e14dcfSSatish Balay 633a7e14dcfSSatish Balay /* Perform the linesearch */ 634a7e14dcfSSatish Balay fold = f; 635a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->Xold);CHKERRQ(ierr); 636a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tl->Gold);CHKERRQ(ierr); 637a7e14dcfSSatish Balay 638a7e14dcfSSatish Balay step = 1.0; 639a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 640a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 641a7e14dcfSSatish Balay 64253506e15SBarry Smith while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NTL_GRADIENT) { /* Linesearch failed */ 643a7e14dcfSSatish Balay /* Linesearch failed */ 644a7e14dcfSSatish Balay f = fold; 645a7e14dcfSSatish Balay ierr = VecCopy(tl->Xold, tao->solution);CHKERRQ(ierr); 646a7e14dcfSSatish Balay ierr = VecCopy(tl->Gold, tao->gradient);CHKERRQ(ierr); 647a7e14dcfSSatish Balay 648a7e14dcfSSatish Balay switch(stepType) { 649a7e14dcfSSatish Balay case NTL_NEWTON: 650a7e14dcfSSatish Balay /* Failed to obtain acceptable iterate with Newton step */ 651a7e14dcfSSatish Balay 652a7e14dcfSSatish Balay if (NTL_PC_BFGS != tl->pc_type) { 653a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and being updated 654a7e14dcfSSatish Balay Must use gradient direction in this case */ 655a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 656a7e14dcfSSatish Balay ++tl->grad; 657a7e14dcfSSatish Balay stepType = NTL_GRADIENT; 65853506e15SBarry Smith } else { 659a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 660a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 661a7e14dcfSSatish Balay 662a7e14dcfSSatish Balay 663a7e14dcfSSatish Balay /* Check for success (descent direction) */ 664a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 665a7e14dcfSSatish Balay if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 666a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced 667a7e14dcfSSatish Balay not a number. We can assert bfgsUpdates > 1 in this case 668a7e14dcfSSatish Balay Use steepest descent direction (scaled) */ 669a7e14dcfSSatish Balay 670a7e14dcfSSatish Balay if (f != 0.0) { 671a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 67253506e15SBarry Smith } else { 673a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 674a7e14dcfSSatish Balay } 675a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 676a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 677a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 678a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 679a7e14dcfSSatish Balay 680a7e14dcfSSatish Balay bfgsUpdates = 1; 681a7e14dcfSSatish Balay ++tl->sgrad; 682a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 68353506e15SBarry Smith } else { 684a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 685a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 686a7e14dcfSSatish Balay ++tl->sgrad; 687a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 68853506e15SBarry Smith } else { 689a7e14dcfSSatish Balay ++tl->bfgs; 690a7e14dcfSSatish Balay stepType = NTL_BFGS; 691a7e14dcfSSatish Balay } 692a7e14dcfSSatish Balay } 693a7e14dcfSSatish Balay } 694a7e14dcfSSatish Balay break; 695a7e14dcfSSatish Balay 696a7e14dcfSSatish Balay case NTL_BFGS: 697a7e14dcfSSatish Balay /* Can only enter if pc_type == NTL_PC_BFGS 698a7e14dcfSSatish Balay Failed to obtain acceptable iterate with BFGS step 699a7e14dcfSSatish Balay Attempt to use the scaled gradient direction */ 700a7e14dcfSSatish Balay 701a7e14dcfSSatish Balay if (f != 0.0) { 702a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 70353506e15SBarry Smith } else { 704a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 705a7e14dcfSSatish Balay } 706a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 707a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 708a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 709a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 710a7e14dcfSSatish Balay 711a7e14dcfSSatish Balay bfgsUpdates = 1; 712a7e14dcfSSatish Balay ++tl->sgrad; 713a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 714a7e14dcfSSatish Balay break; 715a7e14dcfSSatish Balay 716a7e14dcfSSatish Balay case NTL_SCALED_GRADIENT: 717a7e14dcfSSatish Balay /* Can only enter if pc_type == NTL_PC_BFGS 718a7e14dcfSSatish Balay The scaled gradient step did not produce a new iterate; 719a7e14dcfSSatish Balay attemp to use the gradient direction. 720a7e14dcfSSatish Balay Need to make sure we are not using a different diagonal scaling */ 721a7e14dcfSSatish Balay ierr = MatLMVMSetScale(tl->M, tl->Diag);CHKERRQ(ierr); 722a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, 1.0);CHKERRQ(ierr); 723a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 724a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 725a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 726a7e14dcfSSatish Balay 727a7e14dcfSSatish Balay bfgsUpdates = 1; 728a7e14dcfSSatish Balay ++tl->grad; 729a7e14dcfSSatish Balay stepType = NTL_GRADIENT; 730a7e14dcfSSatish Balay break; 731a7e14dcfSSatish Balay } 732a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 733a7e14dcfSSatish Balay 734a7e14dcfSSatish Balay /* This may be incorrect; linesearch has values for stepmax and stepmin 735a7e14dcfSSatish Balay that should be reset. */ 736a7e14dcfSSatish Balay step = 1.0; 737a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 738a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 739a7e14dcfSSatish Balay } 740a7e14dcfSSatish Balay 74153506e15SBarry Smith if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 742a7e14dcfSSatish Balay /* Failed to find an improving point */ 743a7e14dcfSSatish Balay f = fold; 744a7e14dcfSSatish Balay ierr = VecCopy(tl->Xold, tao->solution);CHKERRQ(ierr); 745a7e14dcfSSatish Balay ierr = VecCopy(tl->Gold, tao->gradient);CHKERRQ(ierr); 746a7e14dcfSSatish Balay tao->trust = 0.0; 747a7e14dcfSSatish Balay step = 0.0; 748a7e14dcfSSatish Balay reason = TAO_DIVERGED_LS_FAILURE; 749a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 750a7e14dcfSSatish Balay break; 75153506e15SBarry Smith } else if (stepType == NTL_NEWTON) { 752a7e14dcfSSatish Balay if (step < tl->nu1) { 753a7e14dcfSSatish Balay /* Very bad step taken; reduce radius */ 754a7e14dcfSSatish Balay tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust); 75553506e15SBarry Smith } else if (step < tl->nu2) { 756a7e14dcfSSatish Balay /* Reasonably bad step taken; reduce radius */ 757a7e14dcfSSatish Balay tao->trust = tl->omega2 * PetscMin(norm_d, tao->trust); 75853506e15SBarry Smith } else if (step < tl->nu3) { 759a7e14dcfSSatish Balay /* Reasonable step was taken; leave radius alone */ 760a7e14dcfSSatish Balay if (tl->omega3 < 1.0) { 761a7e14dcfSSatish Balay tao->trust = tl->omega3 * PetscMin(norm_d, tao->trust); 76253506e15SBarry Smith } else if (tl->omega3 > 1.0) { 763a7e14dcfSSatish Balay tao->trust = PetscMax(tl->omega3 * norm_d, tao->trust); 764a7e14dcfSSatish Balay } 76553506e15SBarry Smith } else if (step < tl->nu4) { 766a7e14dcfSSatish Balay /* Full step taken; increase the radius */ 767a7e14dcfSSatish Balay tao->trust = PetscMax(tl->omega4 * norm_d, tao->trust); 76853506e15SBarry Smith } else { 769a7e14dcfSSatish Balay /* More than full step taken; increase the radius */ 770a7e14dcfSSatish Balay tao->trust = PetscMax(tl->omega5 * norm_d, tao->trust); 771a7e14dcfSSatish Balay } 77253506e15SBarry Smith } else { 773a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 774a7e14dcfSSatish Balay tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust); 775a7e14dcfSSatish Balay } 77653506e15SBarry Smith } else { 777a7e14dcfSSatish Balay /* Trust-region step is accepted */ 778a7e14dcfSSatish Balay ierr = VecCopy(tl->W, tao->solution);CHKERRQ(ierr); 779a7e14dcfSSatish Balay f = ftrial; 780a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr); 781a7e14dcfSSatish Balay ++tl->ntrust; 782a7e14dcfSSatish Balay } 783a7e14dcfSSatish Balay 784a7e14dcfSSatish Balay /* The radius may have been increased; modify if it is too large */ 785a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 786a7e14dcfSSatish Balay 787e4cb33bbSBarry Smith /* Check for converged */ 788a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 78953506e15SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 790a7e14dcfSSatish Balay needH = 1; 791a7e14dcfSSatish Balay 7928931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, tao->trust, &reason);CHKERRQ(ierr); 793a7e14dcfSSatish Balay } 794a7e14dcfSSatish Balay PetscFunctionReturn(0); 795a7e14dcfSSatish Balay } 796a7e14dcfSSatish Balay 797a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 798a7e14dcfSSatish Balay #undef __FUNCT__ 799a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_NTL" 800441846f8SBarry Smith static PetscErrorCode TaoSetUp_NTL(Tao tao) 801a7e14dcfSSatish Balay { 802a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 803a7e14dcfSSatish Balay PetscErrorCode ierr; 804a7e14dcfSSatish Balay 805a7e14dcfSSatish Balay PetscFunctionBegin; 806a7e14dcfSSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); } 807a7e14dcfSSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);} 808a7e14dcfSSatish Balay if (!tl->W) { ierr = VecDuplicate(tao->solution, &tl->W);CHKERRQ(ierr);} 809a7e14dcfSSatish Balay if (!tl->Xold) { ierr = VecDuplicate(tao->solution, &tl->Xold);CHKERRQ(ierr);} 810a7e14dcfSSatish Balay if (!tl->Gold) { ierr = VecDuplicate(tao->solution, &tl->Gold);CHKERRQ(ierr);} 811a7e14dcfSSatish Balay tl->Diag = 0; 812a7e14dcfSSatish Balay tl->M = 0; 813a7e14dcfSSatish Balay PetscFunctionReturn(0); 814a7e14dcfSSatish Balay } 815a7e14dcfSSatish Balay 816a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 817a7e14dcfSSatish Balay #undef __FUNCT__ 818a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_NTL" 819441846f8SBarry Smith static PetscErrorCode TaoDestroy_NTL(Tao tao) 820a7e14dcfSSatish Balay { 821a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 822a7e14dcfSSatish Balay PetscErrorCode ierr; 823a7e14dcfSSatish Balay 824a7e14dcfSSatish Balay PetscFunctionBegin; 825a7e14dcfSSatish Balay if (tao->setupcalled) { 826a7e14dcfSSatish Balay ierr = VecDestroy(&tl->W);CHKERRQ(ierr); 827a7e14dcfSSatish Balay ierr = VecDestroy(&tl->Xold);CHKERRQ(ierr); 828a7e14dcfSSatish Balay ierr = VecDestroy(&tl->Gold);CHKERRQ(ierr); 829a7e14dcfSSatish Balay } 830a7e14dcfSSatish Balay ierr = VecDestroy(&tl->Diag);CHKERRQ(ierr); 831a7e14dcfSSatish Balay ierr = MatDestroy(&tl->M);CHKERRQ(ierr); 832a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 833a7e14dcfSSatish Balay PetscFunctionReturn(0); 834a7e14dcfSSatish Balay } 835a7e14dcfSSatish Balay 836a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 837a7e14dcfSSatish Balay #undef __FUNCT__ 838a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_NTL" 8398c34d3f5SBarry Smith static PetscErrorCode TaoSetFromOptions_NTL(PetscOptions *PetscOptionsObject,Tao tao) 840a7e14dcfSSatish Balay { 841a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 842a7e14dcfSSatish Balay PetscErrorCode ierr; 843a7e14dcfSSatish Balay 844a7e14dcfSSatish Balay PetscFunctionBegin; 8451a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Newton trust region with line search method for unconstrained optimization");CHKERRQ(ierr); 84694ae4db5SBarry 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); 84794ae4db5SBarry 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); 84894ae4db5SBarry 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); 84994ae4db5SBarry 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); 85094ae4db5SBarry 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); 85194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_eta1", "poor steplength; reduce radius", "", tl->eta1, &tl->eta1,NULL);CHKERRQ(ierr); 85294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_eta2", "reasonable steplength; leave radius alone", "", tl->eta2, &tl->eta2,NULL);CHKERRQ(ierr); 85394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_eta3", "good steplength; increase radius", "", tl->eta3, &tl->eta3,NULL);CHKERRQ(ierr); 85494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_eta4", "excellent steplength; greatly increase radius", "", tl->eta4, &tl->eta4,NULL);CHKERRQ(ierr); 85594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha1", "", "", tl->alpha1, &tl->alpha1,NULL);CHKERRQ(ierr); 85694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha2", "", "", tl->alpha2, &tl->alpha2,NULL);CHKERRQ(ierr); 85794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha3", "", "", tl->alpha3, &tl->alpha3,NULL);CHKERRQ(ierr); 85894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha4", "", "", tl->alpha4, &tl->alpha4,NULL);CHKERRQ(ierr); 85994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_alpha5", "", "", tl->alpha5, &tl->alpha5,NULL);CHKERRQ(ierr); 86094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_nu1", "poor steplength; reduce radius", "", tl->nu1, &tl->nu1,NULL);CHKERRQ(ierr); 86194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_nu2", "reasonable steplength; leave radius alone", "", tl->nu2, &tl->nu2,NULL);CHKERRQ(ierr); 86294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_nu3", "good steplength; increase radius", "", tl->nu3, &tl->nu3,NULL);CHKERRQ(ierr); 86394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_nu4", "excellent steplength; greatly increase radius", "", tl->nu4, &tl->nu4,NULL);CHKERRQ(ierr); 86494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega1", "", "", tl->omega1, &tl->omega1,NULL);CHKERRQ(ierr); 86594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega2", "", "", tl->omega2, &tl->omega2,NULL);CHKERRQ(ierr); 86694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega3", "", "", tl->omega3, &tl->omega3,NULL);CHKERRQ(ierr); 86794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega4", "", "", tl->omega4, &tl->omega4,NULL);CHKERRQ(ierr); 86894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_omega5", "", "", tl->omega5, &tl->omega5,NULL);CHKERRQ(ierr); 86994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_mu1_i", "", "", tl->mu1_i, &tl->mu1_i,NULL);CHKERRQ(ierr); 87094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_mu2_i", "", "", tl->mu2_i, &tl->mu2_i,NULL);CHKERRQ(ierr); 87194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma1_i", "", "", tl->gamma1_i, &tl->gamma1_i,NULL);CHKERRQ(ierr); 87294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma2_i", "", "", tl->gamma2_i, &tl->gamma2_i,NULL);CHKERRQ(ierr); 87394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma3_i", "", "", tl->gamma3_i, &tl->gamma3_i,NULL);CHKERRQ(ierr); 87494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma4_i", "", "", tl->gamma4_i, &tl->gamma4_i,NULL);CHKERRQ(ierr); 87594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_theta_i", "", "", tl->theta_i, &tl->theta_i,NULL);CHKERRQ(ierr); 87694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_mu1", "", "", tl->mu1, &tl->mu1,NULL);CHKERRQ(ierr); 87794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_mu2", "", "", tl->mu2, &tl->mu2,NULL);CHKERRQ(ierr); 87894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma1", "", "", tl->gamma1, &tl->gamma1,NULL);CHKERRQ(ierr); 87994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma2", "", "", tl->gamma2, &tl->gamma2,NULL);CHKERRQ(ierr); 88094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma3", "", "", tl->gamma3, &tl->gamma3,NULL);CHKERRQ(ierr); 88194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_gamma4", "", "", tl->gamma4, &tl->gamma4,NULL);CHKERRQ(ierr); 88294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_theta", "", "", tl->theta, &tl->theta,NULL);CHKERRQ(ierr); 88394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_min_radius", "lower bound on initial radius", "", tl->min_radius, &tl->min_radius,NULL);CHKERRQ(ierr); 88494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_max_radius", "upper bound on radius", "", tl->max_radius, &tl->max_radius,NULL);CHKERRQ(ierr); 88594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntl_epsilon", "tolerance used when computing actual and predicted reduction", "", tl->epsilon, &tl->epsilon,NULL);CHKERRQ(ierr); 886a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 887a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 888a7e14dcfSSatish Balay ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 889a7e14dcfSSatish Balay PetscFunctionReturn(0); 890a7e14dcfSSatish Balay } 891a7e14dcfSSatish Balay 892a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 893a7e14dcfSSatish Balay #undef __FUNCT__ 894a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_NTL" 895441846f8SBarry Smith static PetscErrorCode TaoView_NTL(Tao tao, PetscViewer viewer) 896a7e14dcfSSatish Balay { 897a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 898a7e14dcfSSatish Balay PetscInt nrejects; 899a7e14dcfSSatish Balay PetscBool isascii; 900a7e14dcfSSatish Balay PetscErrorCode ierr; 901a7e14dcfSSatish Balay 902a7e14dcfSSatish Balay PetscFunctionBegin; 903a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 904a7e14dcfSSatish Balay if (isascii) { 905a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 906a7e14dcfSSatish Balay if (NTL_PC_BFGS == tl->pc_type && tl->M) { 907a7e14dcfSSatish Balay ierr = MatLMVMGetRejects(tl->M, &nrejects);CHKERRQ(ierr); 908a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n", nrejects);CHKERRQ(ierr); 909a7e14dcfSSatish Balay } 910a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Trust-region steps: %D\n", tl->ntrust);CHKERRQ(ierr); 911a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Newton search steps: %D\n", tl->newt);CHKERRQ(ierr); 912a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "BFGS search steps: %D\n", tl->bfgs);CHKERRQ(ierr); 913a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient search steps: %D\n", tl->sgrad);CHKERRQ(ierr); 914a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Gradient search steps: %D\n", tl->grad);CHKERRQ(ierr); 915a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 916a7e14dcfSSatish Balay } 917a7e14dcfSSatish Balay PetscFunctionReturn(0); 918a7e14dcfSSatish Balay } 919a7e14dcfSSatish Balay 920a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 9211522df2eSJason Sarich /*MC 9221522df2eSJason Sarich TAONTR - Newton's method with trust region and linesearch 9231522df2eSJason Sarich for unconstrained minimization. 9241522df2eSJason Sarich At each iteration, the Newton trust region method solves the system for d 9251522df2eSJason Sarich and performs a line search in the d direction: 9261522df2eSJason Sarich 9271522df2eSJason Sarich min_d .5 dT Hk d + gkT d, s.t. ||d|| < Delta_k 9281522df2eSJason Sarich 9291522df2eSJason Sarich Options Database Keys: 9301522df2eSJason Sarich + -tao_ntl_ksp_type - "nash","stcg","gltr" 9311522df2eSJason Sarich . -tao_ntl_pc_type - "none","ahess","bfgs","petsc" 9321522df2eSJason Sarich . -tao_ntl_bfgs_scale_type - type of scaling with bfgs pc, "ahess" or "bfgs" 9331522df2eSJason Sarich . -tao_ntl_init_type - "constant","direction","interpolation" 9341522df2eSJason Sarich . -tao_ntl_update_type - "reduction","interpolation" 9351522df2eSJason Sarich . -tao_ntl_min_radius - lower bound on trust region radius 9361522df2eSJason Sarich . -tao_ntl_max_radius - upper bound on trust region radius 9371522df2eSJason Sarich . -tao_ntl_epsilon - tolerance for accepting actual / predicted reduction 9381522df2eSJason Sarich . -tao_ntl_mu1_i - mu1 interpolation init factor 9391522df2eSJason Sarich . -tao_ntl_mu2_i - mu2 interpolation init factor 9401522df2eSJason Sarich . -tao_ntl_gamma1_i - gamma1 interpolation init factor 9411522df2eSJason Sarich . -tao_ntl_gamma2_i - gamma2 interpolation init factor 9421522df2eSJason Sarich . -tao_ntl_gamma3_i - gamma3 interpolation init factor 9431522df2eSJason Sarich . -tao_ntl_gamma4_i - gamma4 interpolation init factor 9441522df2eSJason Sarich . -tao_ntl_theta_i - thetha1 interpolation init factor 9451522df2eSJason Sarich . -tao_ntl_eta1 - eta1 reduction update factor 9461522df2eSJason Sarich . -tao_ntl_eta2 - eta2 reduction update factor 9471522df2eSJason Sarich . -tao_ntl_eta3 - eta3 reduction update factor 9481522df2eSJason Sarich . -tao_ntl_eta4 - eta4 reduction update factor 9491522df2eSJason Sarich . -tao_ntl_alpha1 - alpha1 reduction update factor 9501522df2eSJason Sarich . -tao_ntl_alpha2 - alpha2 reduction update factor 9511522df2eSJason Sarich . -tao_ntl_alpha3 - alpha3 reduction update factor 9521522df2eSJason Sarich . -tao_ntl_alpha4 - alpha4 reduction update factor 9531522df2eSJason Sarich . -tao_ntl_alpha4 - alpha4 reduction update factor 9541522df2eSJason Sarich . -tao_ntl_mu1 - mu1 interpolation update 9551522df2eSJason Sarich . -tao_ntl_mu2 - mu2 interpolation update 9561522df2eSJason Sarich . -tao_ntl_gamma1 - gamma1 interpolcation update 9571522df2eSJason Sarich . -tao_ntl_gamma2 - gamma2 interpolcation update 9581522df2eSJason Sarich . -tao_ntl_gamma3 - gamma3 interpolcation update 9591522df2eSJason Sarich . -tao_ntl_gamma4 - gamma4 interpolation update 9601522df2eSJason Sarich - -tao_ntl_theta - theta1 interpolation update 9611522df2eSJason Sarich 9621eb8069cSJason Sarich Level: beginner 9631522df2eSJason Sarich M*/ 9641522df2eSJason Sarich 965a7e14dcfSSatish Balay #undef __FUNCT__ 966a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_NTL" 967728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NTL(Tao tao) 968a7e14dcfSSatish Balay { 969a7e14dcfSSatish Balay TAO_NTL *tl; 970a7e14dcfSSatish Balay PetscErrorCode ierr; 9718caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 97253506e15SBarry Smith 973a7e14dcfSSatish Balay PetscFunctionBegin; 9743c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&tl);CHKERRQ(ierr); 975a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_NTL; 976a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_NTL; 977a7e14dcfSSatish Balay tao->ops->view = TaoView_NTL; 978a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_NTL; 979a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_NTL; 980a7e14dcfSSatish Balay 981a7e14dcfSSatish Balay tao->max_it = 50; 9826f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 9836f4723b1SBarry Smith tao->fatol = 1e-5; 9846f4723b1SBarry Smith tao->frtol = 1e-5; 9856f4723b1SBarry Smith #else 986a7e14dcfSSatish Balay tao->fatol = 1e-10; 987a7e14dcfSSatish Balay tao->frtol = 1e-10; 9886f4723b1SBarry Smith #endif 989a7e14dcfSSatish Balay tao->data = (void*)tl; 990a7e14dcfSSatish Balay 991a7e14dcfSSatish Balay tao->trust0 = 100.0; 992a7e14dcfSSatish Balay 993a7e14dcfSSatish Balay 994a7e14dcfSSatish Balay /* Default values for trust-region radius update based on steplength */ 995a7e14dcfSSatish Balay tl->nu1 = 0.25; 996a7e14dcfSSatish Balay tl->nu2 = 0.50; 997a7e14dcfSSatish Balay tl->nu3 = 1.00; 998a7e14dcfSSatish Balay tl->nu4 = 1.25; 999a7e14dcfSSatish Balay 1000a7e14dcfSSatish Balay tl->omega1 = 0.25; 1001a7e14dcfSSatish Balay tl->omega2 = 0.50; 1002a7e14dcfSSatish Balay tl->omega3 = 1.00; 1003a7e14dcfSSatish Balay tl->omega4 = 2.00; 1004a7e14dcfSSatish Balay tl->omega5 = 4.00; 1005a7e14dcfSSatish Balay 1006a7e14dcfSSatish Balay /* Default values for trust-region radius update based on reduction */ 1007a7e14dcfSSatish Balay tl->eta1 = 1.0e-4; 1008a7e14dcfSSatish Balay tl->eta2 = 0.25; 1009a7e14dcfSSatish Balay tl->eta3 = 0.50; 1010a7e14dcfSSatish Balay tl->eta4 = 0.90; 1011a7e14dcfSSatish Balay 1012a7e14dcfSSatish Balay tl->alpha1 = 0.25; 1013a7e14dcfSSatish Balay tl->alpha2 = 0.50; 1014a7e14dcfSSatish Balay tl->alpha3 = 1.00; 1015a7e14dcfSSatish Balay tl->alpha4 = 2.00; 1016a7e14dcfSSatish Balay tl->alpha5 = 4.00; 1017a7e14dcfSSatish Balay 1018a7e14dcfSSatish Balay /* Default values for trust-region radius update based on interpolation */ 1019a7e14dcfSSatish Balay tl->mu1 = 0.10; 1020a7e14dcfSSatish Balay tl->mu2 = 0.50; 1021a7e14dcfSSatish Balay 1022a7e14dcfSSatish Balay tl->gamma1 = 0.25; 1023a7e14dcfSSatish Balay tl->gamma2 = 0.50; 1024a7e14dcfSSatish Balay tl->gamma3 = 2.00; 1025a7e14dcfSSatish Balay tl->gamma4 = 4.00; 1026a7e14dcfSSatish Balay 1027a7e14dcfSSatish Balay tl->theta = 0.05; 1028a7e14dcfSSatish Balay 1029a7e14dcfSSatish Balay /* Default values for trust region initialization based on interpolation */ 1030a7e14dcfSSatish Balay tl->mu1_i = 0.35; 1031a7e14dcfSSatish Balay tl->mu2_i = 0.50; 1032a7e14dcfSSatish Balay 1033a7e14dcfSSatish Balay tl->gamma1_i = 0.0625; 1034a7e14dcfSSatish Balay tl->gamma2_i = 0.5; 1035a7e14dcfSSatish Balay tl->gamma3_i = 2.0; 1036a7e14dcfSSatish Balay tl->gamma4_i = 5.0; 1037a7e14dcfSSatish Balay 1038a7e14dcfSSatish Balay tl->theta_i = 0.25; 1039a7e14dcfSSatish Balay 1040a7e14dcfSSatish Balay /* Remaining parameters */ 1041a7e14dcfSSatish Balay tl->min_radius = 1.0e-10; 1042a7e14dcfSSatish Balay tl->max_radius = 1.0e10; 1043a7e14dcfSSatish Balay tl->epsilon = 1.0e-6; 1044a7e14dcfSSatish Balay 1045a7e14dcfSSatish Balay tl->ksp_type = NTL_KSP_STCG; 1046a7e14dcfSSatish Balay tl->pc_type = NTL_PC_BFGS; 1047a7e14dcfSSatish Balay tl->bfgs_scale_type = BFGS_SCALE_AHESS; 1048a7e14dcfSSatish Balay tl->init_type = NTL_INIT_INTERPOLATION; 1049a7e14dcfSSatish Balay tl->update_type = NTL_UPDATE_REDUCTION; 1050a7e14dcfSSatish Balay 1051a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 1052a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 1053441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch, tao);CHKERRQ(ierr); 1054*5d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 1055a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr); 1056*5d527766SPatrick Farrell ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr); 1057a7e14dcfSSatish Balay PetscFunctionReturn(0); 1058a7e14dcfSSatish Balay } 1059728e0ed0SBarry Smith 1060a7e14dcfSSatish Balay 1061a7e14dcfSSatish Balay 1062a7e14dcfSSatish Balay 1063