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> 6aaa7dc30SBarry Smith #include <petsc-private/kspimpl.h> 7aaa7dc30SBarry 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; 78*e4cb33bbSBarry Smith TaoConvergedReason reason; 79*e4cb33bbSBarry 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 MatStructure matflag; 89a7e14dcfSSatish Balay PetscErrorCode ierr; 90a7e14dcfSSatish Balay PetscInt stepType; 91a7e14dcfSSatish Balay PetscInt iter = 0,its; 92a7e14dcfSSatish Balay 93a7e14dcfSSatish Balay PetscInt bfgsUpdates = 0; 94a7e14dcfSSatish Balay PetscInt needH; 95a7e14dcfSSatish Balay 96a7e14dcfSSatish Balay PetscInt i_max = 5; 97a7e14dcfSSatish Balay PetscInt j_max = 1; 98a7e14dcfSSatish Balay PetscInt i, j, n, N; 99a7e14dcfSSatish Balay 100a7e14dcfSSatish Balay PetscInt tr_reject; 101a7e14dcfSSatish Balay 102a7e14dcfSSatish Balay PetscFunctionBegin; 103a7e14dcfSSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 104a7e14dcfSSatish Balay ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by ntl algorithm\n");CHKERRQ(ierr); 105a7e14dcfSSatish Balay } 106a7e14dcfSSatish Balay 107a7e14dcfSSatish Balay /* Initialize trust-region radius */ 108a7e14dcfSSatish Balay tao->trust = tao->trust0; 109a7e14dcfSSatish Balay 110a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 111a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 112a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 113a7e14dcfSSatish Balay 114a7e14dcfSSatish Balay if (NTL_PC_BFGS == tl->pc_type && !tl->M) { 115a7e14dcfSSatish Balay ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 116a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 117a7e14dcfSSatish Balay ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&tl->M);CHKERRQ(ierr); 118a7e14dcfSSatish Balay ierr = MatLMVMAllocateVectors(tl->M,tao->solution);CHKERRQ(ierr); 119a7e14dcfSSatish Balay } 120a7e14dcfSSatish Balay 121a7e14dcfSSatish Balay /* Check convergence criteria */ 122a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 123a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 12453506e15SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 125a7e14dcfSSatish Balay needH = 1; 126a7e14dcfSSatish Balay 127a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 12853506e15SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 129a7e14dcfSSatish Balay 130a7e14dcfSSatish Balay /* Create vectors for the limited memory preconditioner */ 13153506e15SBarry Smith if ((NTL_PC_BFGS == tl->pc_type) && (BFGS_SCALE_BFGS != tl->bfgs_scale_type)) { 132a7e14dcfSSatish Balay if (!tl->Diag) { 133a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tl->Diag);CHKERRQ(ierr); 134a7e14dcfSSatish Balay } 135a7e14dcfSSatish Balay } 136a7e14dcfSSatish Balay 137a7e14dcfSSatish Balay /* Modify the linear solver to a conjugate gradient method */ 138a7e14dcfSSatish Balay switch(tl->ksp_type) { 139a7e14dcfSSatish Balay case NTL_KSP_NASH: 140a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPNASH);CHKERRQ(ierr); 141a7e14dcfSSatish Balay if (tao->ksp->ops->setfromoptions) { 142a7e14dcfSSatish Balay (*tao->ksp->ops->setfromoptions)(tao->ksp); 143a7e14dcfSSatish Balay } 144a7e14dcfSSatish Balay break; 145a7e14dcfSSatish Balay 146a7e14dcfSSatish Balay case NTL_KSP_STCG: 147a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPSTCG);CHKERRQ(ierr); 148a7e14dcfSSatish Balay if (tao->ksp->ops->setfromoptions) { 149a7e14dcfSSatish Balay (*tao->ksp->ops->setfromoptions)(tao->ksp); 150a7e14dcfSSatish Balay } 151a7e14dcfSSatish Balay break; 152a7e14dcfSSatish Balay 153a7e14dcfSSatish Balay default: 154a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPGLTR);CHKERRQ(ierr); 155a7e14dcfSSatish Balay if (tao->ksp->ops->setfromoptions) { 156a7e14dcfSSatish Balay (*tao->ksp->ops->setfromoptions)(tao->ksp); 157a7e14dcfSSatish Balay } 158a7e14dcfSSatish Balay break; 159a7e14dcfSSatish Balay } 160a7e14dcfSSatish Balay 161a7e14dcfSSatish Balay /* Modify the preconditioner to use the bfgs approximation */ 162a7e14dcfSSatish Balay ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 163a7e14dcfSSatish Balay switch(tl->pc_type) { 164a7e14dcfSSatish Balay case NTL_PC_NONE: 165a7e14dcfSSatish Balay ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr); 166a7e14dcfSSatish Balay if (pc->ops->setfromoptions) { 167a7e14dcfSSatish Balay (*pc->ops->setfromoptions)(pc); 168a7e14dcfSSatish Balay } 169a7e14dcfSSatish Balay break; 170a7e14dcfSSatish Balay 171a7e14dcfSSatish Balay case NTL_PC_AHESS: 172a7e14dcfSSatish Balay ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr); 173a7e14dcfSSatish Balay if (pc->ops->setfromoptions) { 174a7e14dcfSSatish Balay (*pc->ops->setfromoptions)(pc); 175a7e14dcfSSatish Balay } 176a7e14dcfSSatish Balay ierr = PCJacobiSetUseAbs(pc);CHKERRQ(ierr); 177a7e14dcfSSatish Balay break; 178a7e14dcfSSatish Balay 179a7e14dcfSSatish Balay case NTL_PC_BFGS: 180a7e14dcfSSatish Balay ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr); 181a7e14dcfSSatish Balay if (pc->ops->setfromoptions) { 182a7e14dcfSSatish Balay (*pc->ops->setfromoptions)(pc); 183a7e14dcfSSatish Balay } 184a7e14dcfSSatish Balay ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr); 185a7e14dcfSSatish Balay ierr = PCShellSetContext(pc, tl->M);CHKERRQ(ierr); 186a7e14dcfSSatish Balay ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr); 187a7e14dcfSSatish Balay break; 188a7e14dcfSSatish Balay 189a7e14dcfSSatish Balay default: 190a7e14dcfSSatish Balay /* Use the pc method set by pc_type */ 191a7e14dcfSSatish Balay break; 192a7e14dcfSSatish Balay } 193a7e14dcfSSatish Balay 194a7e14dcfSSatish Balay /* Initialize trust-region radius. The initialization is only performed 195a7e14dcfSSatish Balay when we are using Steihaug-Toint or the Generalized Lanczos method. */ 196a7e14dcfSSatish Balay switch(tl->init_type) { 197a7e14dcfSSatish Balay case NTL_INIT_CONSTANT: 198a7e14dcfSSatish Balay /* Use the initial radius specified */ 199a7e14dcfSSatish Balay break; 200a7e14dcfSSatish Balay 201a7e14dcfSSatish Balay case NTL_INIT_INTERPOLATION: 202a7e14dcfSSatish Balay /* Use the initial radius specified */ 203a7e14dcfSSatish Balay max_radius = 0.0; 204a7e14dcfSSatish Balay 205a7e14dcfSSatish Balay for (j = 0; j < j_max; ++j) { 206a7e14dcfSSatish Balay fmin = f; 207a7e14dcfSSatish Balay sigma = 0.0; 208a7e14dcfSSatish Balay 209a7e14dcfSSatish Balay if (needH) { 210a7e14dcfSSatish Balay ierr = TaoComputeHessian(tao, tao->solution, &tao->hessian, &tao->hessian_pre, &matflag);CHKERRQ(ierr); 211a7e14dcfSSatish Balay needH = 0; 212a7e14dcfSSatish Balay } 213a7e14dcfSSatish Balay 214a7e14dcfSSatish Balay for (i = 0; i < i_max; ++i) { 215a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 216a7e14dcfSSatish Balay ierr = VecAXPY(tl->W, -tao->trust/gnorm, tao->gradient);CHKERRQ(ierr); 217a7e14dcfSSatish Balay 218a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 219a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 220a7e14dcfSSatish Balay tau = tl->gamma1_i; 22153506e15SBarry Smith } else { 222a7e14dcfSSatish Balay if (ftrial < fmin) { 223a7e14dcfSSatish Balay fmin = ftrial; 224a7e14dcfSSatish Balay sigma = -tao->trust / gnorm; 225a7e14dcfSSatish Balay } 226a7e14dcfSSatish Balay 227a7e14dcfSSatish Balay ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 228a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr); 229a7e14dcfSSatish Balay 230a7e14dcfSSatish Balay prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm)); 231a7e14dcfSSatish Balay actred = f - ftrial; 23253506e15SBarry Smith if ((PetscAbsScalar(actred) <= tl->epsilon) && (PetscAbsScalar(prered) <= tl->epsilon)) { 233a7e14dcfSSatish Balay kappa = 1.0; 23453506e15SBarry Smith } else { 235a7e14dcfSSatish Balay kappa = actred / prered; 236a7e14dcfSSatish Balay } 237a7e14dcfSSatish Balay 238a7e14dcfSSatish Balay tau_1 = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust + (1.0 - tl->theta_i) * prered - actred); 239a7e14dcfSSatish Balay tau_2 = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust - (1.0 + tl->theta_i) * prered + actred); 240a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 241a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 242a7e14dcfSSatish Balay 243a7e14dcfSSatish Balay if (PetscAbsScalar(kappa - 1.0) <= tl->mu1_i) { 244a7e14dcfSSatish Balay /* Great agreement */ 245a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 246a7e14dcfSSatish Balay 247a7e14dcfSSatish Balay if (tau_max < 1.0) { 248a7e14dcfSSatish Balay tau = tl->gamma3_i; 24953506e15SBarry Smith } else if (tau_max > tl->gamma4_i) { 250a7e14dcfSSatish Balay tau = tl->gamma4_i; 25153506e15SBarry Smith } else if (tau_1 >= 1.0 && tau_1 <= tl->gamma4_i && tau_2 < 1.0) { 252a7e14dcfSSatish Balay tau = tau_1; 25353506e15SBarry Smith } else if (tau_2 >= 1.0 && tau_2 <= tl->gamma4_i && tau_1 < 1.0) { 254a7e14dcfSSatish Balay tau = tau_2; 25553506e15SBarry Smith } else { 256a7e14dcfSSatish Balay tau = tau_max; 257a7e14dcfSSatish Balay } 25853506e15SBarry Smith } else if (PetscAbsScalar(kappa - 1.0) <= tl->mu2_i) { 259a7e14dcfSSatish Balay /* Good agreement */ 260a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 261a7e14dcfSSatish Balay 262a7e14dcfSSatish Balay if (tau_max < tl->gamma2_i) { 263a7e14dcfSSatish Balay tau = tl->gamma2_i; 26453506e15SBarry Smith } else if (tau_max > tl->gamma3_i) { 265a7e14dcfSSatish Balay tau = tl->gamma3_i; 26653506e15SBarry Smith } else { 267a7e14dcfSSatish Balay tau = tau_max; 268a7e14dcfSSatish Balay } 26953506e15SBarry Smith } else { 270a7e14dcfSSatish Balay /* Not good agreement */ 271a7e14dcfSSatish Balay if (tau_min > 1.0) { 272a7e14dcfSSatish Balay tau = tl->gamma2_i; 27353506e15SBarry Smith } else if (tau_max < tl->gamma1_i) { 274a7e14dcfSSatish Balay tau = tl->gamma1_i; 27553506e15SBarry Smith } else if ((tau_min < tl->gamma1_i) && (tau_max >= 1.0)) { 276a7e14dcfSSatish Balay tau = tl->gamma1_i; 27753506e15SBarry Smith } else if ((tau_1 >= tl->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < tl->gamma1_i) || (tau_2 >= 1.0))) { 278a7e14dcfSSatish Balay tau = tau_1; 27953506e15SBarry Smith } else if ((tau_2 >= tl->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < tl->gamma1_i) || (tau_2 >= 1.0))) { 280a7e14dcfSSatish Balay tau = tau_2; 28153506e15SBarry Smith } else { 282a7e14dcfSSatish Balay tau = tau_max; 283a7e14dcfSSatish Balay } 284a7e14dcfSSatish Balay } 285a7e14dcfSSatish Balay } 286a7e14dcfSSatish Balay tao->trust = tau * tao->trust; 287a7e14dcfSSatish Balay } 288a7e14dcfSSatish Balay 289a7e14dcfSSatish Balay if (fmin < f) { 290a7e14dcfSSatish Balay f = fmin; 291a7e14dcfSSatish Balay ierr = VecAXPY(tao->solution, sigma, tao->gradient);CHKERRQ(ierr); 292a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr); 293a7e14dcfSSatish Balay 294a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 29553506e15SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 296a7e14dcfSSatish Balay needH = 1; 297a7e14dcfSSatish Balay 298a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 29953506e15SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 300a7e14dcfSSatish Balay } 301a7e14dcfSSatish Balay } 302a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, max_radius); 303a7e14dcfSSatish Balay 304a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 305a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 306a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 307a7e14dcfSSatish Balay break; 308a7e14dcfSSatish Balay 309a7e14dcfSSatish Balay default: 310a7e14dcfSSatish Balay /* Norm of the first direction will initialize radius */ 311a7e14dcfSSatish Balay tao->trust = 0.0; 312a7e14dcfSSatish Balay break; 313a7e14dcfSSatish Balay } 314a7e14dcfSSatish Balay 315a7e14dcfSSatish Balay /* Set initial scaling for the BFGS preconditioner 316a7e14dcfSSatish Balay This step is done after computing the initial trust-region radius 317a7e14dcfSSatish Balay since the function value may have decreased */ 318a7e14dcfSSatish Balay if (NTL_PC_BFGS == tl->pc_type) { 319a7e14dcfSSatish Balay if (f != 0.0) { 320a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 32153506e15SBarry Smith } else { 322a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 323a7e14dcfSSatish Balay } 324a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 325a7e14dcfSSatish Balay } 326a7e14dcfSSatish Balay 327a7e14dcfSSatish Balay /* Set counter for gradient/reset steps */ 328a7e14dcfSSatish Balay tl->ntrust = 0; 329a7e14dcfSSatish Balay tl->newt = 0; 330a7e14dcfSSatish Balay tl->bfgs = 0; 331a7e14dcfSSatish Balay tl->sgrad = 0; 332a7e14dcfSSatish Balay tl->grad = 0; 333a7e14dcfSSatish Balay 334a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 335a7e14dcfSSatish Balay while (reason == TAO_CONTINUE_ITERATING) { 336a7e14dcfSSatish Balay ++iter; 337a7e14dcfSSatish Balay 338a7e14dcfSSatish Balay /* Compute the Hessian */ 339a7e14dcfSSatish Balay if (needH) { 340a7e14dcfSSatish Balay ierr = TaoComputeHessian(tao, tao->solution, &tao->hessian, &tao->hessian_pre, &matflag);CHKERRQ(ierr); 341a7e14dcfSSatish Balay needH = 0; 342a7e14dcfSSatish Balay } 343a7e14dcfSSatish Balay 344a7e14dcfSSatish Balay if (NTL_PC_BFGS == tl->pc_type) { 345a7e14dcfSSatish Balay if (BFGS_SCALE_AHESS == tl->bfgs_scale_type) { 346a7e14dcfSSatish Balay /* Obtain diagonal for the bfgs preconditioner */ 347a7e14dcfSSatish Balay ierr = MatGetDiagonal(tao->hessian, tl->Diag);CHKERRQ(ierr); 348a7e14dcfSSatish Balay ierr = VecAbs(tl->Diag);CHKERRQ(ierr); 349a7e14dcfSSatish Balay ierr = VecReciprocal(tl->Diag);CHKERRQ(ierr); 350a7e14dcfSSatish Balay ierr = MatLMVMSetScale(tl->M, tl->Diag);CHKERRQ(ierr); 351a7e14dcfSSatish Balay } 352a7e14dcfSSatish Balay 353a7e14dcfSSatish Balay /* Update the limited memory preconditioner */ 354a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M,tao->solution, tao->gradient);CHKERRQ(ierr); 355a7e14dcfSSatish Balay ++bfgsUpdates; 356a7e14dcfSSatish Balay } 357a7e14dcfSSatish Balay ierr = KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre, matflag);CHKERRQ(ierr); 358a7e14dcfSSatish Balay /* Solve the Newton system of equations */ 359a7e14dcfSSatish Balay if (NTL_KSP_NASH == tl->ksp_type) { 360a7e14dcfSSatish Balay ierr = KSPNASHSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 361a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 362a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 363a7e14dcfSSatish Balay tao->ksp_its+=its; 364a7e14dcfSSatish Balay ierr = KSPNASHGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 365a7e14dcfSSatish Balay } else if (NTL_KSP_STCG == tl->ksp_type) { 366a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(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; 370a7e14dcfSSatish Balay ierr = KSPSTCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 371a7e14dcfSSatish Balay } else { /* NTL_KSP_GLTR */ 372a7e14dcfSSatish Balay ierr = KSPGLTRSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 373a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 374a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 375a7e14dcfSSatish Balay tao->ksp_its+=its; 376a7e14dcfSSatish Balay ierr = KSPGLTRGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 377a7e14dcfSSatish Balay } 378a7e14dcfSSatish Balay 379a7e14dcfSSatish Balay if (0.0 == tao->trust) { 380a7e14dcfSSatish Balay /* Radius was uninitialized; use the norm of the direction */ 381a7e14dcfSSatish Balay if (norm_d > 0.0) { 382a7e14dcfSSatish Balay tao->trust = norm_d; 383a7e14dcfSSatish Balay 384a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 385a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 386a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 38753506e15SBarry Smith } else { 388a7e14dcfSSatish Balay /* The direction was bad; set radius to default value and re-solve 389a7e14dcfSSatish Balay the trust-region subproblem to get a direction */ 390a7e14dcfSSatish Balay tao->trust = tao->trust0; 391a7e14dcfSSatish Balay 392a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 393a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->min_radius); 394a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 395a7e14dcfSSatish Balay 396a7e14dcfSSatish Balay if (NTL_KSP_NASH == tl->ksp_type) { 397a7e14dcfSSatish Balay ierr = KSPNASHSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 398a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 399a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 400a7e14dcfSSatish Balay tao->ksp_its+=its; 401a7e14dcfSSatish Balay ierr = KSPNASHGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 402a7e14dcfSSatish Balay } else if (NTL_KSP_STCG == tl->ksp_type) { 403a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 404a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 405a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 406a7e14dcfSSatish Balay tao->ksp_its+=its; 407a7e14dcfSSatish Balay ierr = KSPSTCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 408a7e14dcfSSatish Balay } else { /* NTL_KSP_GLTR */ 409a7e14dcfSSatish Balay ierr = KSPGLTRSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 410a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 411a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 412a7e14dcfSSatish Balay tao->ksp_its+=its; 413a7e14dcfSSatish Balay ierr = KSPGLTRGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 414a7e14dcfSSatish Balay } 415a7e14dcfSSatish Balay 416a7e14dcfSSatish Balay 41753506e15SBarry Smith if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 418a7e14dcfSSatish Balay } 419a7e14dcfSSatish Balay } 420a7e14dcfSSatish Balay 421a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 422a7e14dcfSSatish Balay ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 42353506e15SBarry Smith if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (NTL_PC_BFGS == tl->pc_type) && (bfgsUpdates > 1)) { 424a7e14dcfSSatish Balay /* Preconditioner is numerically indefinite; reset the 425a7e14dcfSSatish Balay approximate if using BFGS preconditioning. */ 426a7e14dcfSSatish Balay 427a7e14dcfSSatish Balay if (f != 0.0) { 428a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 42953506e15SBarry Smith } else { 430a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 431a7e14dcfSSatish Balay } 432a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 433a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 434a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 435a7e14dcfSSatish Balay bfgsUpdates = 1; 436a7e14dcfSSatish Balay } 437a7e14dcfSSatish Balay 438a7e14dcfSSatish Balay /* Check trust-region reduction conditions */ 439a7e14dcfSSatish Balay tr_reject = 0; 440a7e14dcfSSatish Balay if (NTL_UPDATE_REDUCTION == tl->update_type) { 441a7e14dcfSSatish Balay /* Get predicted reduction */ 442a7e14dcfSSatish Balay if (NTL_KSP_NASH == tl->ksp_type) { 443a7e14dcfSSatish Balay ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 444a7e14dcfSSatish Balay } else if (NTL_KSP_STCG == tl->ksp_type) { 445a7e14dcfSSatish Balay ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 446a7e14dcfSSatish Balay } else { /* gltr */ 447a7e14dcfSSatish Balay ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 448a7e14dcfSSatish Balay } 449a7e14dcfSSatish Balay 450a7e14dcfSSatish Balay if (prered >= 0.0) { 451a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot 452a7e14dcfSSatish Balay happen in infinite precision arithmetic. Step should 453a7e14dcfSSatish Balay be rejected! */ 454a7e14dcfSSatish Balay tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 455a7e14dcfSSatish Balay tr_reject = 1; 45653506e15SBarry Smith } else { 457a7e14dcfSSatish Balay /* Compute trial step and function value */ 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 462a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 463a7e14dcfSSatish Balay tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 464a7e14dcfSSatish Balay tr_reject = 1; 46553506e15SBarry Smith } else { 466a7e14dcfSSatish Balay /* Compute and actual reduction */ 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 /* Accept of reject the step and update radius */ 477a7e14dcfSSatish Balay if (kappa < tl->eta1) { 478a7e14dcfSSatish Balay /* Reject the step */ 479a7e14dcfSSatish Balay tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 480a7e14dcfSSatish Balay tr_reject = 1; 48153506e15SBarry Smith } else { 482a7e14dcfSSatish Balay /* Accept the step */ 483a7e14dcfSSatish Balay if (kappa < tl->eta2) { 484a7e14dcfSSatish Balay /* Marginal bad step */ 485a7e14dcfSSatish Balay tao->trust = tl->alpha2 * PetscMin(tao->trust, norm_d); 48653506e15SBarry Smith } else if (kappa < tl->eta3) { 487a7e14dcfSSatish Balay /* Reasonable step */ 488a7e14dcfSSatish Balay tao->trust = tl->alpha3 * tao->trust; 48953506e15SBarry Smith } else if (kappa < tl->eta4) { 490a7e14dcfSSatish Balay /* Good step */ 491a7e14dcfSSatish Balay tao->trust = PetscMax(tl->alpha4 * norm_d, tao->trust); 49253506e15SBarry Smith } else { 493a7e14dcfSSatish Balay /* Very good step */ 494a7e14dcfSSatish Balay tao->trust = PetscMax(tl->alpha5 * norm_d, tao->trust); 495a7e14dcfSSatish Balay } 496a7e14dcfSSatish Balay } 497a7e14dcfSSatish Balay } 498a7e14dcfSSatish Balay } 49953506e15SBarry Smith } else { 500a7e14dcfSSatish Balay /* Get predicted reduction */ 501a7e14dcfSSatish Balay if (NTL_KSP_NASH == tl->ksp_type) { 502a7e14dcfSSatish Balay ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 503a7e14dcfSSatish Balay } else if (NTL_KSP_STCG == tl->ksp_type) { 504a7e14dcfSSatish Balay ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 505a7e14dcfSSatish Balay } else { /* gltr */ 506a7e14dcfSSatish Balay ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 507a7e14dcfSSatish Balay } 508a7e14dcfSSatish Balay 509a7e14dcfSSatish Balay if (prered >= 0.0) { 510a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot 511a7e14dcfSSatish Balay happen in infinite precision arithmetic. Step should 512a7e14dcfSSatish Balay be rejected! */ 513a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 514a7e14dcfSSatish Balay tr_reject = 1; 51553506e15SBarry Smith } else { 516a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 517a7e14dcfSSatish Balay ierr = VecAXPY(tl->W, 1.0, tao->stepdirection);CHKERRQ(ierr); 518a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 519a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 520a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 521a7e14dcfSSatish Balay tr_reject = 1; 52253506e15SBarry Smith } else { 523a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 524a7e14dcfSSatish Balay 525a7e14dcfSSatish Balay actred = f - ftrial; 526a7e14dcfSSatish Balay prered = -prered; 527a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= tl->epsilon) && 528a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= tl->epsilon)) { 529a7e14dcfSSatish Balay kappa = 1.0; 53053506e15SBarry Smith } else { 531a7e14dcfSSatish Balay kappa = actred / prered; 532a7e14dcfSSatish Balay } 533a7e14dcfSSatish Balay 534a7e14dcfSSatish Balay tau_1 = tl->theta * gdx / (tl->theta * gdx - (1.0 - tl->theta) * prered + actred); 535a7e14dcfSSatish Balay tau_2 = tl->theta * gdx / (tl->theta * gdx + (1.0 + tl->theta) * prered - actred); 536a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 537a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 538a7e14dcfSSatish Balay 539a7e14dcfSSatish Balay if (kappa >= 1.0 - tl->mu1) { 540a7e14dcfSSatish Balay /* Great agreement; accept step and update radius */ 541a7e14dcfSSatish Balay if (tau_max < 1.0) { 542a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d); 54353506e15SBarry Smith } else if (tau_max > tl->gamma4) { 544a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->gamma4 * norm_d); 54553506e15SBarry Smith } else { 546a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 547a7e14dcfSSatish Balay } 54853506e15SBarry Smith } else if (kappa >= 1.0 - tl->mu2) { 549a7e14dcfSSatish Balay /* Good agreement */ 550a7e14dcfSSatish Balay 551a7e14dcfSSatish Balay if (tau_max < tl->gamma2) { 552a7e14dcfSSatish Balay tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d); 55353506e15SBarry Smith } else if (tau_max > tl->gamma3) { 554a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d); 555a7e14dcfSSatish Balay } else if (tau_max < 1.0) { 556a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 55753506e15SBarry Smith } else { 558a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 559a7e14dcfSSatish Balay } 56053506e15SBarry Smith } else { 561a7e14dcfSSatish Balay /* Not good agreement */ 562a7e14dcfSSatish Balay if (tau_min > 1.0) { 563a7e14dcfSSatish Balay tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d); 56453506e15SBarry Smith } else if (tau_max < tl->gamma1) { 565a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 56653506e15SBarry Smith } else if ((tau_min < tl->gamma1) && (tau_max >= 1.0)) { 567a7e14dcfSSatish Balay tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 56853506e15SBarry Smith } else if ((tau_1 >= tl->gamma1) && (tau_1 < 1.0) && ((tau_2 < tl->gamma1) || (tau_2 >= 1.0))) { 569a7e14dcfSSatish Balay tao->trust = tau_1 * PetscMin(tao->trust, norm_d); 57053506e15SBarry Smith } else if ((tau_2 >= tl->gamma1) && (tau_2 < 1.0) && ((tau_1 < tl->gamma1) || (tau_2 >= 1.0))) { 571a7e14dcfSSatish Balay tao->trust = tau_2 * PetscMin(tao->trust, norm_d); 57253506e15SBarry Smith } else { 573a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 574a7e14dcfSSatish Balay } 575a7e14dcfSSatish Balay tr_reject = 1; 576a7e14dcfSSatish Balay } 577a7e14dcfSSatish Balay } 578a7e14dcfSSatish Balay } 579a7e14dcfSSatish Balay } 580a7e14dcfSSatish Balay 581a7e14dcfSSatish Balay if (tr_reject) { 582a7e14dcfSSatish Balay /* The trust-region constraints rejected the step. Apply a linesearch. 583a7e14dcfSSatish Balay Check for descent direction. */ 584a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 585a7e14dcfSSatish Balay if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 586a7e14dcfSSatish Balay /* Newton step is not descent or direction produced Inf or NaN */ 587a7e14dcfSSatish Balay 588a7e14dcfSSatish Balay if (NTL_PC_BFGS != tl->pc_type) { 589a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and updated 590a7e14dcfSSatish Balay Must use gradient direction in this case */ 591a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 592a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 593a7e14dcfSSatish Balay ++tl->grad; 594a7e14dcfSSatish Balay stepType = NTL_GRADIENT; 59553506e15SBarry Smith } else { 596a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 597a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 598a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 599a7e14dcfSSatish Balay 600a7e14dcfSSatish Balay /* Check for success (descent direction) */ 601a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 602a7e14dcfSSatish Balay if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) { 603a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced not a number 604a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case because 605a7e14dcfSSatish Balay the first solve produces the scaled gradient direction, 606a7e14dcfSSatish Balay which is guaranteed to be descent */ 607a7e14dcfSSatish Balay 608a7e14dcfSSatish Balay /* Use steepest descent direction (scaled) */ 609a7e14dcfSSatish Balay if (f != 0.0) { 610a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 61153506e15SBarry Smith } else { 612a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 613a7e14dcfSSatish Balay } 614a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 615a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 616a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 617a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 618a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 619a7e14dcfSSatish Balay 620a7e14dcfSSatish Balay bfgsUpdates = 1; 621a7e14dcfSSatish Balay ++tl->sgrad; 622a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 62353506e15SBarry Smith } else { 624a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 625a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 626a7e14dcfSSatish Balay ++tl->sgrad; 627a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 62853506e15SBarry Smith } else { 629a7e14dcfSSatish Balay ++tl->bfgs; 630a7e14dcfSSatish Balay stepType = NTL_BFGS; 631a7e14dcfSSatish Balay } 632a7e14dcfSSatish Balay } 633a7e14dcfSSatish Balay } 63453506e15SBarry Smith } else { 635a7e14dcfSSatish Balay /* Computed Newton step is descent */ 636a7e14dcfSSatish Balay ++tl->newt; 637a7e14dcfSSatish Balay stepType = NTL_NEWTON; 638a7e14dcfSSatish Balay } 639a7e14dcfSSatish Balay 640a7e14dcfSSatish Balay /* Perform the linesearch */ 641a7e14dcfSSatish Balay fold = f; 642a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tl->Xold);CHKERRQ(ierr); 643a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tl->Gold);CHKERRQ(ierr); 644a7e14dcfSSatish Balay 645a7e14dcfSSatish Balay step = 1.0; 646a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 647a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 648a7e14dcfSSatish Balay 64953506e15SBarry Smith while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NTL_GRADIENT) { /* Linesearch failed */ 650a7e14dcfSSatish Balay /* Linesearch failed */ 651a7e14dcfSSatish Balay f = fold; 652a7e14dcfSSatish Balay ierr = VecCopy(tl->Xold, tao->solution);CHKERRQ(ierr); 653a7e14dcfSSatish Balay ierr = VecCopy(tl->Gold, tao->gradient);CHKERRQ(ierr); 654a7e14dcfSSatish Balay 655a7e14dcfSSatish Balay switch(stepType) { 656a7e14dcfSSatish Balay case NTL_NEWTON: 657a7e14dcfSSatish Balay /* Failed to obtain acceptable iterate with Newton step */ 658a7e14dcfSSatish Balay 659a7e14dcfSSatish Balay if (NTL_PC_BFGS != tl->pc_type) { 660a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and being updated 661a7e14dcfSSatish Balay Must use gradient direction in this case */ 662a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 663a7e14dcfSSatish Balay ++tl->grad; 664a7e14dcfSSatish Balay stepType = NTL_GRADIENT; 66553506e15SBarry Smith } else { 666a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 667a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 668a7e14dcfSSatish Balay 669a7e14dcfSSatish Balay 670a7e14dcfSSatish Balay /* Check for success (descent direction) */ 671a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 672a7e14dcfSSatish Balay if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 673a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced 674a7e14dcfSSatish Balay not a number. We can assert bfgsUpdates > 1 in this case 675a7e14dcfSSatish Balay Use steepest descent direction (scaled) */ 676a7e14dcfSSatish Balay 677a7e14dcfSSatish Balay if (f != 0.0) { 678a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 67953506e15SBarry Smith } else { 680a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 681a7e14dcfSSatish Balay } 682a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 683a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 684a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 685a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 686a7e14dcfSSatish Balay 687a7e14dcfSSatish Balay bfgsUpdates = 1; 688a7e14dcfSSatish Balay ++tl->sgrad; 689a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 69053506e15SBarry Smith } else { 691a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 692a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 693a7e14dcfSSatish Balay ++tl->sgrad; 694a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 69553506e15SBarry Smith } else { 696a7e14dcfSSatish Balay ++tl->bfgs; 697a7e14dcfSSatish Balay stepType = NTL_BFGS; 698a7e14dcfSSatish Balay } 699a7e14dcfSSatish Balay } 700a7e14dcfSSatish Balay } 701a7e14dcfSSatish Balay break; 702a7e14dcfSSatish Balay 703a7e14dcfSSatish Balay case NTL_BFGS: 704a7e14dcfSSatish Balay /* Can only enter if pc_type == NTL_PC_BFGS 705a7e14dcfSSatish Balay Failed to obtain acceptable iterate with BFGS step 706a7e14dcfSSatish Balay Attempt to use the scaled gradient direction */ 707a7e14dcfSSatish Balay 708a7e14dcfSSatish Balay if (f != 0.0) { 709a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 71053506e15SBarry Smith } else { 711a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 712a7e14dcfSSatish Balay } 713a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, delta);CHKERRQ(ierr); 714a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 715a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 716a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 717a7e14dcfSSatish Balay 718a7e14dcfSSatish Balay bfgsUpdates = 1; 719a7e14dcfSSatish Balay ++tl->sgrad; 720a7e14dcfSSatish Balay stepType = NTL_SCALED_GRADIENT; 721a7e14dcfSSatish Balay break; 722a7e14dcfSSatish Balay 723a7e14dcfSSatish Balay case NTL_SCALED_GRADIENT: 724a7e14dcfSSatish Balay /* Can only enter if pc_type == NTL_PC_BFGS 725a7e14dcfSSatish Balay The scaled gradient step did not produce a new iterate; 726a7e14dcfSSatish Balay attemp to use the gradient direction. 727a7e14dcfSSatish Balay Need to make sure we are not using a different diagonal scaling */ 728a7e14dcfSSatish Balay ierr = MatLMVMSetScale(tl->M, tl->Diag);CHKERRQ(ierr); 729a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tl->M, 1.0);CHKERRQ(ierr); 730a7e14dcfSSatish Balay ierr = MatLMVMReset(tl->M);CHKERRQ(ierr); 731a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 732a7e14dcfSSatish Balay ierr = MatLMVMSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 733a7e14dcfSSatish Balay 734a7e14dcfSSatish Balay bfgsUpdates = 1; 735a7e14dcfSSatish Balay ++tl->grad; 736a7e14dcfSSatish Balay stepType = NTL_GRADIENT; 737a7e14dcfSSatish Balay break; 738a7e14dcfSSatish Balay } 739a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 740a7e14dcfSSatish Balay 741a7e14dcfSSatish Balay /* This may be incorrect; linesearch has values for stepmax and stepmin 742a7e14dcfSSatish Balay that should be reset. */ 743a7e14dcfSSatish Balay step = 1.0; 744a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 745a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 746a7e14dcfSSatish Balay } 747a7e14dcfSSatish Balay 74853506e15SBarry Smith if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 749a7e14dcfSSatish Balay /* Failed to find an improving point */ 750a7e14dcfSSatish Balay f = fold; 751a7e14dcfSSatish Balay ierr = VecCopy(tl->Xold, tao->solution);CHKERRQ(ierr); 752a7e14dcfSSatish Balay ierr = VecCopy(tl->Gold, tao->gradient);CHKERRQ(ierr); 753a7e14dcfSSatish Balay tao->trust = 0.0; 754a7e14dcfSSatish Balay step = 0.0; 755a7e14dcfSSatish Balay reason = TAO_DIVERGED_LS_FAILURE; 756a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 757a7e14dcfSSatish Balay break; 75853506e15SBarry Smith } else if (stepType == NTL_NEWTON) { 759a7e14dcfSSatish Balay if (step < tl->nu1) { 760a7e14dcfSSatish Balay /* Very bad step taken; reduce radius */ 761a7e14dcfSSatish Balay tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust); 76253506e15SBarry Smith } else if (step < tl->nu2) { 763a7e14dcfSSatish Balay /* Reasonably bad step taken; reduce radius */ 764a7e14dcfSSatish Balay tao->trust = tl->omega2 * PetscMin(norm_d, tao->trust); 76553506e15SBarry Smith } else if (step < tl->nu3) { 766a7e14dcfSSatish Balay /* Reasonable step was taken; leave radius alone */ 767a7e14dcfSSatish Balay if (tl->omega3 < 1.0) { 768a7e14dcfSSatish Balay tao->trust = tl->omega3 * PetscMin(norm_d, tao->trust); 76953506e15SBarry Smith } else if (tl->omega3 > 1.0) { 770a7e14dcfSSatish Balay tao->trust = PetscMax(tl->omega3 * norm_d, tao->trust); 771a7e14dcfSSatish Balay } 77253506e15SBarry Smith } else if (step < tl->nu4) { 773a7e14dcfSSatish Balay /* Full step taken; increase the radius */ 774a7e14dcfSSatish Balay tao->trust = PetscMax(tl->omega4 * norm_d, tao->trust); 77553506e15SBarry Smith } else { 776a7e14dcfSSatish Balay /* More than full step taken; increase the radius */ 777a7e14dcfSSatish Balay tao->trust = PetscMax(tl->omega5 * norm_d, tao->trust); 778a7e14dcfSSatish Balay } 77953506e15SBarry Smith } else { 780a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 781a7e14dcfSSatish Balay tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust); 782a7e14dcfSSatish Balay } 78353506e15SBarry Smith } else { 784a7e14dcfSSatish Balay /* Trust-region step is accepted */ 785a7e14dcfSSatish Balay ierr = VecCopy(tl->W, tao->solution);CHKERRQ(ierr); 786a7e14dcfSSatish Balay f = ftrial; 787a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr); 788a7e14dcfSSatish Balay ++tl->ntrust; 789a7e14dcfSSatish Balay } 790a7e14dcfSSatish Balay 791a7e14dcfSSatish Balay /* The radius may have been increased; modify if it is too large */ 792a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tl->max_radius); 793a7e14dcfSSatish Balay 794*e4cb33bbSBarry Smith /* Check for converged */ 795a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 79653506e15SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 797a7e14dcfSSatish Balay needH = 1; 798a7e14dcfSSatish Balay 799a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, tao->trust, &reason);CHKERRQ(ierr); 800a7e14dcfSSatish Balay } 801a7e14dcfSSatish Balay PetscFunctionReturn(0); 802a7e14dcfSSatish Balay } 803a7e14dcfSSatish Balay 804a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 805a7e14dcfSSatish Balay #undef __FUNCT__ 806a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_NTL" 807441846f8SBarry Smith static PetscErrorCode TaoSetUp_NTL(Tao tao) 808a7e14dcfSSatish Balay { 809a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 810a7e14dcfSSatish Balay PetscErrorCode ierr; 811a7e14dcfSSatish Balay 812a7e14dcfSSatish Balay PetscFunctionBegin; 813a7e14dcfSSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); } 814a7e14dcfSSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);} 815a7e14dcfSSatish Balay if (!tl->W) { ierr = VecDuplicate(tao->solution, &tl->W);CHKERRQ(ierr);} 816a7e14dcfSSatish Balay if (!tl->Xold) { ierr = VecDuplicate(tao->solution, &tl->Xold);CHKERRQ(ierr);} 817a7e14dcfSSatish Balay if (!tl->Gold) { ierr = VecDuplicate(tao->solution, &tl->Gold);CHKERRQ(ierr);} 818a7e14dcfSSatish Balay tl->Diag = 0; 819a7e14dcfSSatish Balay tl->M = 0; 820a7e14dcfSSatish Balay PetscFunctionReturn(0); 821a7e14dcfSSatish Balay } 822a7e14dcfSSatish Balay 823a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 824a7e14dcfSSatish Balay #undef __FUNCT__ 825a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_NTL" 826441846f8SBarry Smith static PetscErrorCode TaoDestroy_NTL(Tao tao) 827a7e14dcfSSatish Balay { 828a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 829a7e14dcfSSatish Balay PetscErrorCode ierr; 830a7e14dcfSSatish Balay 831a7e14dcfSSatish Balay PetscFunctionBegin; 832a7e14dcfSSatish Balay if (tao->setupcalled) { 833a7e14dcfSSatish Balay ierr = VecDestroy(&tl->W);CHKERRQ(ierr); 834a7e14dcfSSatish Balay ierr = VecDestroy(&tl->Xold);CHKERRQ(ierr); 835a7e14dcfSSatish Balay ierr = VecDestroy(&tl->Gold);CHKERRQ(ierr); 836a7e14dcfSSatish Balay } 837a7e14dcfSSatish Balay ierr = VecDestroy(&tl->Diag);CHKERRQ(ierr); 838a7e14dcfSSatish Balay ierr = MatDestroy(&tl->M);CHKERRQ(ierr); 839a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 840a7e14dcfSSatish Balay PetscFunctionReturn(0); 841a7e14dcfSSatish Balay } 842a7e14dcfSSatish Balay 843a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 844a7e14dcfSSatish Balay #undef __FUNCT__ 845a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_NTL" 846441846f8SBarry Smith static PetscErrorCode TaoSetFromOptions_NTL(Tao tao) 847a7e14dcfSSatish Balay { 848a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 849a7e14dcfSSatish Balay PetscErrorCode ierr; 850a7e14dcfSSatish Balay 851a7e14dcfSSatish Balay PetscFunctionBegin; 852a7e14dcfSSatish Balay ierr = PetscOptionsHead("Newton line search method for unconstrained optimization");CHKERRQ(ierr); 853a7e14dcfSSatish Balay ierr = PetscOptionsEList("-tao_ntl_ksp_type", "ksp type", "", NTL_KSP, NTL_KSP_TYPES, NTL_KSP[tl->ksp_type], &tl->ksp_type, 0);CHKERRQ(ierr); 854a7e14dcfSSatish Balay ierr = PetscOptionsEList("-tao_ntl_pc_type", "pc type", "", NTL_PC, NTL_PC_TYPES, NTL_PC[tl->pc_type], &tl->pc_type, 0);CHKERRQ(ierr); 855a7e14dcfSSatish Balay 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, 0);CHKERRQ(ierr); 856a7e14dcfSSatish Balay ierr = PetscOptionsEList("-tao_ntl_init_type", "radius initialization type", "", NTL_INIT, NTL_INIT_TYPES, NTL_INIT[tl->init_type], &tl->init_type, 0);CHKERRQ(ierr); 857a7e14dcfSSatish Balay ierr = PetscOptionsEList("-tao_ntl_update_type", "radius update type", "", NTL_UPDATE, NTL_UPDATE_TYPES, NTL_UPDATE[tl->update_type], &tl->update_type, 0);CHKERRQ(ierr); 858a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_eta1", "poor steplength; reduce radius", "", tl->eta1, &tl->eta1, 0);CHKERRQ(ierr); 859a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_eta2", "reasonable steplength; leave radius alone", "", tl->eta2, &tl->eta2, 0);CHKERRQ(ierr); 860a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_eta3", "good steplength; increase radius", "", tl->eta3, &tl->eta3, 0);CHKERRQ(ierr); 861a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_eta4", "excellent steplength; greatly increase radius", "", tl->eta4, &tl->eta4, 0);CHKERRQ(ierr); 862a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_alpha1", "", "", tl->alpha1, &tl->alpha1, 0);CHKERRQ(ierr); 863a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_alpha2", "", "", tl->alpha2, &tl->alpha2, 0);CHKERRQ(ierr); 864a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_alpha3", "", "", tl->alpha3, &tl->alpha3, 0);CHKERRQ(ierr); 865a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_alpha4", "", "", tl->alpha4, &tl->alpha4, 0);CHKERRQ(ierr); 866a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_alpha5", "", "", tl->alpha5, &tl->alpha5, 0);CHKERRQ(ierr); 867a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_nu1", "poor steplength; reduce radius", "", tl->nu1, &tl->nu1, 0);CHKERRQ(ierr); 868a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_nu2", "reasonable steplength; leave radius alone", "", tl->nu2, &tl->nu2, 0);CHKERRQ(ierr); 869a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_nu3", "good steplength; increase radius", "", tl->nu3, &tl->nu3, 0);CHKERRQ(ierr); 870a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_nu4", "excellent steplength; greatly increase radius", "", tl->nu4, &tl->nu4, 0);CHKERRQ(ierr); 871a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_omega1", "", "", tl->omega1, &tl->omega1, 0);CHKERRQ(ierr); 872a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_omega2", "", "", tl->omega2, &tl->omega2, 0);CHKERRQ(ierr); 873a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_omega3", "", "", tl->omega3, &tl->omega3, 0);CHKERRQ(ierr); 874a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_omega4", "", "", tl->omega4, &tl->omega4, 0);CHKERRQ(ierr); 875a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_omega5", "", "", tl->omega5, &tl->omega5, 0);CHKERRQ(ierr); 876a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_mu1_i", "", "", tl->mu1_i, &tl->mu1_i, 0);CHKERRQ(ierr); 877a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_mu2_i", "", "", tl->mu2_i, &tl->mu2_i, 0);CHKERRQ(ierr); 878a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_gamma1_i", "", "", tl->gamma1_i, &tl->gamma1_i, 0);CHKERRQ(ierr); 879a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_gamma2_i", "", "", tl->gamma2_i, &tl->gamma2_i, 0);CHKERRQ(ierr); 880a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_gamma3_i", "", "", tl->gamma3_i, &tl->gamma3_i, 0);CHKERRQ(ierr); 881a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_gamma4_i", "", "", tl->gamma4_i, &tl->gamma4_i, 0);CHKERRQ(ierr); 882a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_theta_i", "", "", tl->theta_i, &tl->theta_i, 0);CHKERRQ(ierr); 883a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_mu1", "", "", tl->mu1, &tl->mu1, 0);CHKERRQ(ierr); 884a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_mu2", "", "", tl->mu2, &tl->mu2, 0);CHKERRQ(ierr); 885a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_gamma1", "", "", tl->gamma1, &tl->gamma1, 0);CHKERRQ(ierr); 886a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_gamma2", "", "", tl->gamma2, &tl->gamma2, 0);CHKERRQ(ierr); 887a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_gamma3", "", "", tl->gamma3, &tl->gamma3, 0);CHKERRQ(ierr); 888a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_gamma4", "", "", tl->gamma4, &tl->gamma4, 0);CHKERRQ(ierr); 889a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_theta", "", "", tl->theta, &tl->theta, 0);CHKERRQ(ierr); 890a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_min_radius", "lower bound on initial radius", "", tl->min_radius, &tl->min_radius, 0);CHKERRQ(ierr); 891a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_max_radius", "upper bound on radius", "", tl->max_radius, &tl->max_radius, 0);CHKERRQ(ierr); 892a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_ntl_epsilon", "tolerance used when computing actual and predicted reduction", "", tl->epsilon, &tl->epsilon, 0);CHKERRQ(ierr); 893a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 894a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 895a7e14dcfSSatish Balay ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 896a7e14dcfSSatish Balay PetscFunctionReturn(0); 897a7e14dcfSSatish Balay } 898a7e14dcfSSatish Balay 899a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 900a7e14dcfSSatish Balay #undef __FUNCT__ 901a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_NTL" 902441846f8SBarry Smith static PetscErrorCode TaoView_NTL(Tao tao, PetscViewer viewer) 903a7e14dcfSSatish Balay { 904a7e14dcfSSatish Balay TAO_NTL *tl = (TAO_NTL *)tao->data; 905a7e14dcfSSatish Balay PetscInt nrejects; 906a7e14dcfSSatish Balay PetscBool isascii; 907a7e14dcfSSatish Balay PetscErrorCode ierr; 908a7e14dcfSSatish Balay 909a7e14dcfSSatish Balay PetscFunctionBegin; 910a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 911a7e14dcfSSatish Balay if (isascii) { 912a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 913a7e14dcfSSatish Balay if (NTL_PC_BFGS == tl->pc_type && tl->M) { 914a7e14dcfSSatish Balay ierr = MatLMVMGetRejects(tl->M, &nrejects);CHKERRQ(ierr); 915a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n", nrejects);CHKERRQ(ierr); 916a7e14dcfSSatish Balay } 917a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Trust-region steps: %D\n", tl->ntrust);CHKERRQ(ierr); 918a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Newton search steps: %D\n", tl->newt);CHKERRQ(ierr); 919a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "BFGS search steps: %D\n", tl->bfgs);CHKERRQ(ierr); 920a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient search steps: %D\n", tl->sgrad);CHKERRQ(ierr); 921a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Gradient search steps: %D\n", tl->grad);CHKERRQ(ierr); 922a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 923a7e14dcfSSatish Balay } 924a7e14dcfSSatish Balay PetscFunctionReturn(0); 925a7e14dcfSSatish Balay } 926a7e14dcfSSatish Balay 927a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 928a7e14dcfSSatish Balay EXTERN_C_BEGIN 929a7e14dcfSSatish Balay #undef __FUNCT__ 930a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_NTL" 931441846f8SBarry Smith PetscErrorCode TaoCreate_NTL(Tao tao) 932a7e14dcfSSatish Balay { 933a7e14dcfSSatish Balay TAO_NTL *tl; 934a7e14dcfSSatish Balay PetscErrorCode ierr; 935a7e14dcfSSatish Balay const char *morethuente_type = TAOLINESEARCH_MT; 93653506e15SBarry Smith 937a7e14dcfSSatish Balay PetscFunctionBegin; 9383c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&tl);CHKERRQ(ierr); 939a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_NTL; 940a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_NTL; 941a7e14dcfSSatish Balay tao->ops->view = TaoView_NTL; 942a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_NTL; 943a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_NTL; 944a7e14dcfSSatish Balay 945a7e14dcfSSatish Balay tao->max_it = 50; 9466f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 9476f4723b1SBarry Smith tao->fatol = 1e-5; 9486f4723b1SBarry Smith tao->frtol = 1e-5; 9496f4723b1SBarry Smith #else 950a7e14dcfSSatish Balay tao->fatol = 1e-10; 951a7e14dcfSSatish Balay tao->frtol = 1e-10; 9526f4723b1SBarry Smith #endif 953a7e14dcfSSatish Balay tao->data = (void*)tl; 954a7e14dcfSSatish Balay 955a7e14dcfSSatish Balay tao->trust0 = 100.0; 956a7e14dcfSSatish Balay 957a7e14dcfSSatish Balay 958a7e14dcfSSatish Balay /* Default values for trust-region radius update based on steplength */ 959a7e14dcfSSatish Balay tl->nu1 = 0.25; 960a7e14dcfSSatish Balay tl->nu2 = 0.50; 961a7e14dcfSSatish Balay tl->nu3 = 1.00; 962a7e14dcfSSatish Balay tl->nu4 = 1.25; 963a7e14dcfSSatish Balay 964a7e14dcfSSatish Balay tl->omega1 = 0.25; 965a7e14dcfSSatish Balay tl->omega2 = 0.50; 966a7e14dcfSSatish Balay tl->omega3 = 1.00; 967a7e14dcfSSatish Balay tl->omega4 = 2.00; 968a7e14dcfSSatish Balay tl->omega5 = 4.00; 969a7e14dcfSSatish Balay 970a7e14dcfSSatish Balay /* Default values for trust-region radius update based on reduction */ 971a7e14dcfSSatish Balay tl->eta1 = 1.0e-4; 972a7e14dcfSSatish Balay tl->eta2 = 0.25; 973a7e14dcfSSatish Balay tl->eta3 = 0.50; 974a7e14dcfSSatish Balay tl->eta4 = 0.90; 975a7e14dcfSSatish Balay 976a7e14dcfSSatish Balay tl->alpha1 = 0.25; 977a7e14dcfSSatish Balay tl->alpha2 = 0.50; 978a7e14dcfSSatish Balay tl->alpha3 = 1.00; 979a7e14dcfSSatish Balay tl->alpha4 = 2.00; 980a7e14dcfSSatish Balay tl->alpha5 = 4.00; 981a7e14dcfSSatish Balay 982a7e14dcfSSatish Balay /* Default values for trust-region radius update based on interpolation */ 983a7e14dcfSSatish Balay tl->mu1 = 0.10; 984a7e14dcfSSatish Balay tl->mu2 = 0.50; 985a7e14dcfSSatish Balay 986a7e14dcfSSatish Balay tl->gamma1 = 0.25; 987a7e14dcfSSatish Balay tl->gamma2 = 0.50; 988a7e14dcfSSatish Balay tl->gamma3 = 2.00; 989a7e14dcfSSatish Balay tl->gamma4 = 4.00; 990a7e14dcfSSatish Balay 991a7e14dcfSSatish Balay tl->theta = 0.05; 992a7e14dcfSSatish Balay 993a7e14dcfSSatish Balay /* Default values for trust region initialization based on interpolation */ 994a7e14dcfSSatish Balay tl->mu1_i = 0.35; 995a7e14dcfSSatish Balay tl->mu2_i = 0.50; 996a7e14dcfSSatish Balay 997a7e14dcfSSatish Balay tl->gamma1_i = 0.0625; 998a7e14dcfSSatish Balay tl->gamma2_i = 0.5; 999a7e14dcfSSatish Balay tl->gamma3_i = 2.0; 1000a7e14dcfSSatish Balay tl->gamma4_i = 5.0; 1001a7e14dcfSSatish Balay 1002a7e14dcfSSatish Balay tl->theta_i = 0.25; 1003a7e14dcfSSatish Balay 1004a7e14dcfSSatish Balay /* Remaining parameters */ 1005a7e14dcfSSatish Balay tl->min_radius = 1.0e-10; 1006a7e14dcfSSatish Balay tl->max_radius = 1.0e10; 1007a7e14dcfSSatish Balay tl->epsilon = 1.0e-6; 1008a7e14dcfSSatish Balay 1009a7e14dcfSSatish Balay tl->ksp_type = NTL_KSP_STCG; 1010a7e14dcfSSatish Balay tl->pc_type = NTL_PC_BFGS; 1011a7e14dcfSSatish Balay tl->bfgs_scale_type = BFGS_SCALE_AHESS; 1012a7e14dcfSSatish Balay tl->init_type = NTL_INIT_INTERPOLATION; 1013a7e14dcfSSatish Balay tl->update_type = NTL_UPDATE_REDUCTION; 1014a7e14dcfSSatish Balay 1015a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 1016a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 1017441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch, tao);CHKERRQ(ierr); 1018a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr); 1019a7e14dcfSSatish Balay PetscFunctionReturn(0); 1020a7e14dcfSSatish Balay } 1021a7e14dcfSSatish Balay EXTERN_C_END 1022a7e14dcfSSatish Balay 1023a7e14dcfSSatish Balay 1024a7e14dcfSSatish Balay 1025