1aaa7dc30SBarry Smith #include <../src/tao/matrix/lmvmmat.h> 2aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/ntr/ntr.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 NTR_KSP_NASH 0 10a7e14dcfSSatish Balay #define NTR_KSP_STCG 1 11a7e14dcfSSatish Balay #define NTR_KSP_GLTR 2 12a7e14dcfSSatish Balay #define NTR_KSP_TYPES 3 13a7e14dcfSSatish Balay 14a7e14dcfSSatish Balay #define NTR_PC_NONE 0 15a7e14dcfSSatish Balay #define NTR_PC_AHESS 1 16a7e14dcfSSatish Balay #define NTR_PC_BFGS 2 17a7e14dcfSSatish Balay #define NTR_PC_PETSC 3 18a7e14dcfSSatish Balay #define NTR_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 NTR_INIT_CONSTANT 0 25a7e14dcfSSatish Balay #define NTR_INIT_DIRECTION 1 26a7e14dcfSSatish Balay #define NTR_INIT_INTERPOLATION 2 27a7e14dcfSSatish Balay #define NTR_INIT_TYPES 3 28a7e14dcfSSatish Balay 29a7e14dcfSSatish Balay #define NTR_UPDATE_REDUCTION 0 30a7e14dcfSSatish Balay #define NTR_UPDATE_INTERPOLATION 1 31a7e14dcfSSatish Balay #define NTR_UPDATE_TYPES 2 32a7e14dcfSSatish Balay 3353506e15SBarry Smith static const char *NTR_KSP[64] = { "nash", "stcg", "gltr"}; 34a7e14dcfSSatish Balay 3553506e15SBarry Smith static const char *NTR_PC[64] = { "none", "ahess", "bfgs", "petsc"}; 36a7e14dcfSSatish Balay 3753506e15SBarry Smith static const char *BFGS_SCALE[64] = { "ahess", "bfgs"}; 38a7e14dcfSSatish Balay 3953506e15SBarry Smith static const char *NTR_INIT[64] = { "constant", "direction", "interpolation"}; 40a7e14dcfSSatish Balay 4153506e15SBarry Smith static const char *NTR_UPDATE[64] = { "reduction", "interpolation"}; 42a7e14dcfSSatish Balay 43a7e14dcfSSatish Balay /* Routine for BFGS preconditioner */ 44a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec xin, Vec xout); 45a7e14dcfSSatish Balay 46a7e14dcfSSatish Balay /* 47a7e14dcfSSatish Balay TaoSolve_NTR - Implements Newton's Method with a trust region approach 48a7e14dcfSSatish Balay for solving unconstrained minimization problems. 49a7e14dcfSSatish Balay 50a7e14dcfSSatish Balay The basic algorithm is taken from MINPACK-2 (dstrn). 51a7e14dcfSSatish Balay 52a7e14dcfSSatish Balay TaoSolve_NTR computes a local minimizer of a twice differentiable function 53a7e14dcfSSatish Balay f by applying a trust region variant of Newton's method. At each stage 54a7e14dcfSSatish Balay of the algorithm, we use the prconditioned conjugate gradient method to 55a7e14dcfSSatish Balay determine an approximate minimizer of the quadratic equation 56a7e14dcfSSatish Balay 57a7e14dcfSSatish Balay q(s) = <s, Hs + g> 58a7e14dcfSSatish Balay 59a7e14dcfSSatish Balay subject to the trust region constraint 60a7e14dcfSSatish Balay 61a7e14dcfSSatish Balay || s ||_M <= radius, 62a7e14dcfSSatish Balay 63a7e14dcfSSatish Balay where radius is the trust region radius and M is a symmetric positive 64a7e14dcfSSatish Balay definite matrix (the preconditioner). Here g is the gradient and H 65a7e14dcfSSatish Balay is the Hessian matrix. 66a7e14dcfSSatish Balay 67a7e14dcfSSatish Balay Note: TaoSolve_NTR MUST use the iterative solver KSPNASH, KSPSTCG, 68a7e14dcfSSatish Balay or KSPGLTR. Thus, we set KSPNASH, KSPSTCG, or KSPGLTR in this 69a7e14dcfSSatish Balay routine regardless of what the user may have previously specified. 70a7e14dcfSSatish Balay */ 71a7e14dcfSSatish Balay #undef __FUNCT__ 72a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_NTR" 73441846f8SBarry Smith static PetscErrorCode TaoSolve_NTR(Tao tao) 74a7e14dcfSSatish Balay { 75a7e14dcfSSatish Balay TAO_NTR *tr = (TAO_NTR *)tao->data; 76a7e14dcfSSatish Balay PC pc; 77a7e14dcfSSatish Balay KSPConvergedReason ksp_reason; 78e4cb33bbSBarry Smith TaoConvergedReason reason; 79a7e14dcfSSatish Balay PetscReal fmin, ftrial, prered, actred, kappa, sigma, beta; 80a7e14dcfSSatish Balay PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 81a7e14dcfSSatish Balay PetscReal f, gnorm; 82a7e14dcfSSatish Balay 83a7e14dcfSSatish Balay PetscReal delta; 84a7e14dcfSSatish Balay PetscReal norm_d; 85a7e14dcfSSatish Balay PetscErrorCode ierr; 86a7e14dcfSSatish Balay PetscInt iter = 0; 87a7e14dcfSSatish Balay PetscInt bfgsUpdates = 0; 88a7e14dcfSSatish Balay PetscInt needH; 89a7e14dcfSSatish Balay 90a7e14dcfSSatish Balay PetscInt i_max = 5; 91a7e14dcfSSatish Balay PetscInt j_max = 1; 92a7e14dcfSSatish Balay PetscInt i, j, N, n, its; 93a7e14dcfSSatish Balay 94a7e14dcfSSatish Balay PetscFunctionBegin; 95a7e14dcfSSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 96a7e14dcfSSatish Balay ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by ntr algorithm\n");CHKERRQ(ierr); 97a7e14dcfSSatish Balay } 98a7e14dcfSSatish Balay 99a7e14dcfSSatish Balay tao->trust = tao->trust0; 100a7e14dcfSSatish Balay 101a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 102a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->min_radius); 103a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tr->max_radius); 104a7e14dcfSSatish Balay 105a7e14dcfSSatish Balay 106a7e14dcfSSatish Balay if (NTR_PC_BFGS == tr->pc_type && !tr->M) { 107a7e14dcfSSatish Balay ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 108a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 109a7e14dcfSSatish Balay ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&tr->M);CHKERRQ(ierr); 110a7e14dcfSSatish Balay ierr = MatLMVMAllocateVectors(tr->M,tao->solution);CHKERRQ(ierr); 111a7e14dcfSSatish Balay } 112a7e14dcfSSatish Balay 113a7e14dcfSSatish Balay /* Check convergence criteria */ 114a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 115a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 11653506e15SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 117a7e14dcfSSatish Balay needH = 1; 118a7e14dcfSSatish Balay 119a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 12053506e15SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 121a7e14dcfSSatish Balay 122a7e14dcfSSatish Balay /* Create vectors for the limited memory preconditioner */ 123a7e14dcfSSatish Balay if ((NTR_PC_BFGS == tr->pc_type) && 124a7e14dcfSSatish Balay (BFGS_SCALE_BFGS != tr->bfgs_scale_type)) { 125a7e14dcfSSatish Balay if (!tr->Diag) { 126a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tr->Diag);CHKERRQ(ierr); 127a7e14dcfSSatish Balay } 128a7e14dcfSSatish Balay } 129a7e14dcfSSatish Balay 130a7e14dcfSSatish Balay switch(tr->ksp_type) { 131a7e14dcfSSatish Balay case NTR_KSP_NASH: 132a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPNASH);CHKERRQ(ierr); 133*1a1499c8SBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 134a7e14dcfSSatish Balay break; 135a7e14dcfSSatish Balay 136a7e14dcfSSatish Balay case NTR_KSP_STCG: 137a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPSTCG);CHKERRQ(ierr); 138*1a1499c8SBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 139a7e14dcfSSatish Balay break; 140a7e14dcfSSatish Balay 141a7e14dcfSSatish Balay default: 142a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPGLTR);CHKERRQ(ierr); 143*1a1499c8SBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 144a7e14dcfSSatish Balay break; 145a7e14dcfSSatish Balay } 146a7e14dcfSSatish Balay 147a7e14dcfSSatish Balay /* Modify the preconditioner to use the bfgs approximation */ 148a7e14dcfSSatish Balay ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 149a7e14dcfSSatish Balay switch(tr->pc_type) { 150a7e14dcfSSatish Balay case NTR_PC_NONE: 151a7e14dcfSSatish Balay ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr); 152*1a1499c8SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 153a7e14dcfSSatish Balay break; 154a7e14dcfSSatish Balay 155a7e14dcfSSatish Balay case NTR_PC_AHESS: 156a7e14dcfSSatish Balay ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr); 157*1a1499c8SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 158baa89ecbSBarry Smith ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 159a7e14dcfSSatish Balay break; 160a7e14dcfSSatish Balay 161a7e14dcfSSatish Balay case NTR_PC_BFGS: 162a7e14dcfSSatish Balay ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr); 163*1a1499c8SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 164a7e14dcfSSatish Balay ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr); 165a7e14dcfSSatish Balay ierr = PCShellSetContext(pc, tr->M);CHKERRQ(ierr); 166a7e14dcfSSatish Balay ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr); 167a7e14dcfSSatish Balay break; 168a7e14dcfSSatish Balay 169a7e14dcfSSatish Balay default: 170a7e14dcfSSatish Balay /* Use the pc method set by pc_type */ 171a7e14dcfSSatish Balay break; 172a7e14dcfSSatish Balay } 173a7e14dcfSSatish Balay 174a7e14dcfSSatish Balay /* Initialize trust-region radius */ 175a7e14dcfSSatish Balay switch(tr->init_type) { 176a7e14dcfSSatish Balay case NTR_INIT_CONSTANT: 177a7e14dcfSSatish Balay /* Use the initial radius specified */ 178a7e14dcfSSatish Balay break; 179a7e14dcfSSatish Balay 180a7e14dcfSSatish Balay case NTR_INIT_INTERPOLATION: 181a7e14dcfSSatish Balay /* Use the initial radius specified */ 182a7e14dcfSSatish Balay max_radius = 0.0; 183a7e14dcfSSatish Balay 184a7e14dcfSSatish Balay for (j = 0; j < j_max; ++j) { 185a7e14dcfSSatish Balay fmin = f; 186a7e14dcfSSatish Balay sigma = 0.0; 187a7e14dcfSSatish Balay 188a7e14dcfSSatish Balay if (needH) { 189ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 190a7e14dcfSSatish Balay needH = 0; 191a7e14dcfSSatish Balay } 192a7e14dcfSSatish Balay 193a7e14dcfSSatish Balay for (i = 0; i < i_max; ++i) { 194a7e14dcfSSatish Balay 195a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr); 196a7e14dcfSSatish Balay ierr = VecAXPY(tr->W, -tao->trust/gnorm, tao->gradient);CHKERRQ(ierr); 197a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr); 198a7e14dcfSSatish Balay 199a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 200a7e14dcfSSatish Balay tau = tr->gamma1_i; 201a7e14dcfSSatish Balay } 202a7e14dcfSSatish Balay else { 203a7e14dcfSSatish Balay if (ftrial < fmin) { 204a7e14dcfSSatish Balay fmin = ftrial; 205a7e14dcfSSatish Balay sigma = -tao->trust / gnorm; 206a7e14dcfSSatish Balay } 207a7e14dcfSSatish Balay 208a7e14dcfSSatish Balay ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 209a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr); 210a7e14dcfSSatish Balay 211a7e14dcfSSatish Balay prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm)); 212a7e14dcfSSatish Balay actred = f - ftrial; 213a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= tr->epsilon) && 214a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= tr->epsilon)) { 215a7e14dcfSSatish Balay kappa = 1.0; 216a7e14dcfSSatish Balay } 217a7e14dcfSSatish Balay else { 218a7e14dcfSSatish Balay kappa = actred / prered; 219a7e14dcfSSatish Balay } 220a7e14dcfSSatish Balay 221a7e14dcfSSatish Balay tau_1 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust + (1.0 - tr->theta_i) * prered - actred); 222a7e14dcfSSatish Balay tau_2 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust - (1.0 + tr->theta_i) * prered + actred); 223a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 224a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 225a7e14dcfSSatish Balay 226a7e14dcfSSatish Balay if (PetscAbsScalar(kappa - 1.0) <= tr->mu1_i) { 227a7e14dcfSSatish Balay /* Great agreement */ 228a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 229a7e14dcfSSatish Balay 230a7e14dcfSSatish Balay if (tau_max < 1.0) { 231a7e14dcfSSatish Balay tau = tr->gamma3_i; 232a7e14dcfSSatish Balay } 233a7e14dcfSSatish Balay else if (tau_max > tr->gamma4_i) { 234a7e14dcfSSatish Balay tau = tr->gamma4_i; 235a7e14dcfSSatish Balay } 236a7e14dcfSSatish Balay else { 237a7e14dcfSSatish Balay tau = tau_max; 238a7e14dcfSSatish Balay } 239a7e14dcfSSatish Balay } 240a7e14dcfSSatish Balay else if (PetscAbsScalar(kappa - 1.0) <= tr->mu2_i) { 241a7e14dcfSSatish Balay /* Good agreement */ 242a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 243a7e14dcfSSatish Balay 244a7e14dcfSSatish Balay if (tau_max < tr->gamma2_i) { 245a7e14dcfSSatish Balay tau = tr->gamma2_i; 246a7e14dcfSSatish Balay } 247a7e14dcfSSatish Balay else if (tau_max > tr->gamma3_i) { 248a7e14dcfSSatish Balay tau = tr->gamma3_i; 249a7e14dcfSSatish Balay } 250a7e14dcfSSatish Balay else { 251a7e14dcfSSatish Balay tau = tau_max; 252a7e14dcfSSatish Balay } 253a7e14dcfSSatish Balay } 254a7e14dcfSSatish Balay else { 255a7e14dcfSSatish Balay /* Not good agreement */ 256a7e14dcfSSatish Balay if (tau_min > 1.0) { 257a7e14dcfSSatish Balay tau = tr->gamma2_i; 258a7e14dcfSSatish Balay } 259a7e14dcfSSatish Balay else if (tau_max < tr->gamma1_i) { 260a7e14dcfSSatish Balay tau = tr->gamma1_i; 261a7e14dcfSSatish Balay } 262a7e14dcfSSatish Balay else if ((tau_min < tr->gamma1_i) && (tau_max >= 1.0)) { 263a7e14dcfSSatish Balay tau = tr->gamma1_i; 264a7e14dcfSSatish Balay } 265a7e14dcfSSatish Balay else if ((tau_1 >= tr->gamma1_i) && (tau_1 < 1.0) && 266a7e14dcfSSatish Balay ((tau_2 < tr->gamma1_i) || (tau_2 >= 1.0))) { 267a7e14dcfSSatish Balay tau = tau_1; 268a7e14dcfSSatish Balay } 269a7e14dcfSSatish Balay else if ((tau_2 >= tr->gamma1_i) && (tau_2 < 1.0) && 270a7e14dcfSSatish Balay ((tau_1 < tr->gamma1_i) || (tau_2 >= 1.0))) { 271a7e14dcfSSatish Balay tau = tau_2; 272a7e14dcfSSatish Balay } 273a7e14dcfSSatish Balay else { 274a7e14dcfSSatish Balay tau = tau_max; 275a7e14dcfSSatish Balay } 276a7e14dcfSSatish Balay } 277a7e14dcfSSatish Balay } 278a7e14dcfSSatish Balay tao->trust = tau * tao->trust; 279a7e14dcfSSatish Balay } 280a7e14dcfSSatish Balay 281a7e14dcfSSatish Balay if (fmin < f) { 282a7e14dcfSSatish Balay f = fmin; 283a7e14dcfSSatish Balay ierr = VecAXPY(tao->solution, sigma, tao->gradient);CHKERRQ(ierr); 284a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao,tao->solution, tao->gradient);CHKERRQ(ierr); 285a7e14dcfSSatish Balay 286a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 287a7e14dcfSSatish Balay 28853506e15SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 289a7e14dcfSSatish Balay needH = 1; 290a7e14dcfSSatish Balay 291a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 292a7e14dcfSSatish Balay if (reason != TAO_CONTINUE_ITERATING) { 293a7e14dcfSSatish Balay PetscFunctionReturn(0); 294a7e14dcfSSatish Balay } 295a7e14dcfSSatish Balay } 296a7e14dcfSSatish Balay } 297a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, max_radius); 298a7e14dcfSSatish Balay 299a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 300a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->min_radius); 301a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tr->max_radius); 302a7e14dcfSSatish Balay break; 303a7e14dcfSSatish Balay 304a7e14dcfSSatish Balay default: 305a7e14dcfSSatish Balay /* Norm of the first direction will initialize radius */ 306a7e14dcfSSatish Balay tao->trust = 0.0; 307a7e14dcfSSatish Balay break; 308a7e14dcfSSatish Balay } 309a7e14dcfSSatish Balay 310a7e14dcfSSatish Balay /* Set initial scaling for the BFGS preconditioner 311a7e14dcfSSatish Balay This step is done after computing the initial trust-region radius 312a7e14dcfSSatish Balay since the function value may have decreased */ 313a7e14dcfSSatish Balay if (NTR_PC_BFGS == tr->pc_type) { 314a7e14dcfSSatish Balay if (f != 0.0) { 315a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 316a7e14dcfSSatish Balay } 317a7e14dcfSSatish Balay else { 318a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 319a7e14dcfSSatish Balay } 320a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tr->M,delta);CHKERRQ(ierr); 321a7e14dcfSSatish Balay } 322a7e14dcfSSatish Balay 323a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 324a7e14dcfSSatish Balay while (reason == TAO_CONTINUE_ITERATING) { 325a7e14dcfSSatish Balay ++iter; 326ae93cb3cSJason Sarich tao->ksp_its=0; 327a7e14dcfSSatish Balay /* Compute the Hessian */ 328a7e14dcfSSatish Balay if (needH) { 329ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 330a7e14dcfSSatish Balay needH = 0; 331a7e14dcfSSatish Balay } 332a7e14dcfSSatish Balay 333a7e14dcfSSatish Balay if (NTR_PC_BFGS == tr->pc_type) { 334a7e14dcfSSatish Balay if (BFGS_SCALE_AHESS == tr->bfgs_scale_type) { 335a7e14dcfSSatish Balay /* Obtain diagonal for the bfgs preconditioner */ 336a7e14dcfSSatish Balay ierr = MatGetDiagonal(tao->hessian, tr->Diag);CHKERRQ(ierr); 337a7e14dcfSSatish Balay ierr = VecAbs(tr->Diag);CHKERRQ(ierr); 338a7e14dcfSSatish Balay ierr = VecReciprocal(tr->Diag);CHKERRQ(ierr); 339a7e14dcfSSatish Balay ierr = MatLMVMSetScale(tr->M,tr->Diag);CHKERRQ(ierr); 340a7e14dcfSSatish Balay } 341a7e14dcfSSatish Balay 342a7e14dcfSSatish Balay /* Update the limited memory preconditioner */ 343a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr); 344a7e14dcfSSatish Balay ++bfgsUpdates; 345a7e14dcfSSatish Balay } 346a7e14dcfSSatish Balay 347a7e14dcfSSatish Balay while (reason == TAO_CONTINUE_ITERATING) { 34823ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre);CHKERRQ(ierr); 349a7e14dcfSSatish Balay 350a7e14dcfSSatish Balay /* Solve the trust region subproblem */ 351a7e14dcfSSatish Balay if (NTR_KSP_NASH == tr->ksp_type) { 352a7e14dcfSSatish Balay ierr = KSPNASHSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 353a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 354a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 355a7e14dcfSSatish Balay tao->ksp_its+=its; 356ae93cb3cSJason Sarich tao->ksp_tot_its+=its; 357a7e14dcfSSatish Balay ierr = KSPNASHGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 358a7e14dcfSSatish Balay } else if (NTR_KSP_STCG == tr->ksp_type) { 359a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 360a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 361a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 362a7e14dcfSSatish Balay tao->ksp_its+=its; 363ae93cb3cSJason Sarich tao->ksp_tot_its+=its; 364a7e14dcfSSatish Balay ierr = KSPSTCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 365a7e14dcfSSatish Balay } else { /* NTR_KSP_GLTR */ 366a7e14dcfSSatish Balay ierr = KSPGLTRSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 367a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 368a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 369a7e14dcfSSatish Balay tao->ksp_its+=its; 3702d9aa51bSJason Sarich tao->ksp_tot_its+=its; 371a7e14dcfSSatish Balay ierr = KSPGLTRGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 372a7e14dcfSSatish Balay } 373a7e14dcfSSatish Balay 374a7e14dcfSSatish Balay if (0.0 == tao->trust) { 375a7e14dcfSSatish Balay /* Radius was uninitialized; use the norm of the direction */ 376a7e14dcfSSatish Balay if (norm_d > 0.0) { 377a7e14dcfSSatish Balay tao->trust = norm_d; 378a7e14dcfSSatish Balay 379a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 380a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->min_radius); 381a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tr->max_radius); 382a7e14dcfSSatish Balay } 383a7e14dcfSSatish Balay else { 384a7e14dcfSSatish Balay /* The direction was bad; set radius to default value and re-solve 385a7e14dcfSSatish Balay the trust-region subproblem to get a direction */ 386a7e14dcfSSatish Balay tao->trust = tao->trust0; 387a7e14dcfSSatish Balay 388a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 389a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->min_radius); 390a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tr->max_radius); 391a7e14dcfSSatish Balay 392a7e14dcfSSatish Balay if (NTR_KSP_NASH == tr->ksp_type) { 393a7e14dcfSSatish Balay ierr = KSPNASHSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 394a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 395a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 396a7e14dcfSSatish Balay tao->ksp_its+=its; 3972d9aa51bSJason Sarich tao->ksp_tot_its+=its; 398a7e14dcfSSatish Balay ierr = KSPNASHGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 399a7e14dcfSSatish Balay } else if (NTR_KSP_STCG == tr->ksp_type) { 400a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 401a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 402a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 403a7e14dcfSSatish Balay tao->ksp_its+=its; 4042d9aa51bSJason Sarich tao->ksp_tot_its+=its; 405a7e14dcfSSatish Balay ierr = KSPSTCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 406a7e14dcfSSatish Balay } else { /* NTR_KSP_GLTR */ 407a7e14dcfSSatish Balay ierr = KSPGLTRSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 408a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 409a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 410a7e14dcfSSatish Balay tao->ksp_its+=its; 4112d9aa51bSJason Sarich tao->ksp_tot_its+=its; 412a7e14dcfSSatish Balay ierr = KSPGLTRGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 413a7e14dcfSSatish Balay } 414a7e14dcfSSatish Balay 41553506e15SBarry Smith if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 416a7e14dcfSSatish Balay } 417a7e14dcfSSatish Balay } 418a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 419a7e14dcfSSatish Balay ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 420a7e14dcfSSatish Balay if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && 421a7e14dcfSSatish Balay (NTR_PC_BFGS == tr->pc_type) && (bfgsUpdates > 1)) { 422a7e14dcfSSatish Balay /* Preconditioner is numerically indefinite; reset the 423a7e14dcfSSatish Balay approximate if using BFGS preconditioning. */ 424a7e14dcfSSatish Balay 425a7e14dcfSSatish Balay if (f != 0.0) { 426a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 427a7e14dcfSSatish Balay } 428a7e14dcfSSatish Balay else { 429a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 430a7e14dcfSSatish Balay } 431a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(tr->M, delta);CHKERRQ(ierr); 432a7e14dcfSSatish Balay ierr = MatLMVMReset(tr->M);CHKERRQ(ierr); 433a7e14dcfSSatish Balay ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr); 434a7e14dcfSSatish Balay bfgsUpdates = 1; 435a7e14dcfSSatish Balay } 436a7e14dcfSSatish Balay 437a7e14dcfSSatish Balay if (NTR_UPDATE_REDUCTION == tr->update_type) { 438a7e14dcfSSatish Balay /* Get predicted reduction */ 439a7e14dcfSSatish Balay if (NTR_KSP_NASH == tr->ksp_type) { 440a7e14dcfSSatish Balay ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 441a7e14dcfSSatish Balay } else if (NTR_KSP_STCG == tr->ksp_type) { 442a7e14dcfSSatish Balay ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 443a7e14dcfSSatish Balay } else { /* gltr */ 444a7e14dcfSSatish Balay ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 445a7e14dcfSSatish Balay } 446a7e14dcfSSatish Balay 447a7e14dcfSSatish Balay if (prered >= 0.0) { 448a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot 449a7e14dcfSSatish Balay happen in infinite precision arithmetic. Step should 450a7e14dcfSSatish Balay be rejected! */ 451a7e14dcfSSatish Balay tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d); 452a7e14dcfSSatish Balay } 453a7e14dcfSSatish Balay else { 454a7e14dcfSSatish Balay /* Compute trial step and function value */ 455a7e14dcfSSatish Balay ierr = VecCopy(tao->solution,tr->W);CHKERRQ(ierr); 456a7e14dcfSSatish Balay ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr); 457a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr); 458a7e14dcfSSatish Balay 459a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 460a7e14dcfSSatish Balay tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d); 461a7e14dcfSSatish Balay } else { 462a7e14dcfSSatish Balay /* Compute and actual reduction */ 463a7e14dcfSSatish Balay actred = f - ftrial; 464a7e14dcfSSatish Balay prered = -prered; 465a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= tr->epsilon) && 466a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= tr->epsilon)) { 467a7e14dcfSSatish Balay kappa = 1.0; 468a7e14dcfSSatish Balay } 469a7e14dcfSSatish Balay else { 470a7e14dcfSSatish Balay kappa = actred / prered; 471a7e14dcfSSatish Balay } 472a7e14dcfSSatish Balay 473a7e14dcfSSatish Balay /* Accept or reject the step and update radius */ 474a7e14dcfSSatish Balay if (kappa < tr->eta1) { 475a7e14dcfSSatish Balay /* Reject the step */ 476a7e14dcfSSatish Balay tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d); 477a7e14dcfSSatish Balay } 478a7e14dcfSSatish Balay else { 479a7e14dcfSSatish Balay /* Accept the step */ 480a7e14dcfSSatish Balay if (kappa < tr->eta2) { 481a7e14dcfSSatish Balay /* Marginal bad step */ 482a7e14dcfSSatish Balay tao->trust = tr->alpha2 * PetscMin(tao->trust, norm_d); 483a7e14dcfSSatish Balay } 484a7e14dcfSSatish Balay else if (kappa < tr->eta3) { 485a7e14dcfSSatish Balay /* Reasonable step */ 486a7e14dcfSSatish Balay tao->trust = tr->alpha3 * tao->trust; 487a7e14dcfSSatish Balay } 488a7e14dcfSSatish Balay else if (kappa < tr->eta4) { 489a7e14dcfSSatish Balay /* Good step */ 490a7e14dcfSSatish Balay tao->trust = PetscMax(tr->alpha4 * norm_d, tao->trust); 491a7e14dcfSSatish Balay } 492a7e14dcfSSatish Balay else { 493a7e14dcfSSatish Balay /* Very good step */ 494a7e14dcfSSatish Balay tao->trust = PetscMax(tr->alpha5 * norm_d, tao->trust); 495a7e14dcfSSatish Balay } 496a7e14dcfSSatish Balay break; 497a7e14dcfSSatish Balay } 498a7e14dcfSSatish Balay } 499a7e14dcfSSatish Balay } 500a7e14dcfSSatish Balay } 501a7e14dcfSSatish Balay else { 502a7e14dcfSSatish Balay /* Get predicted reduction */ 503a7e14dcfSSatish Balay if (NTR_KSP_NASH == tr->ksp_type) { 504a7e14dcfSSatish Balay ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 505a7e14dcfSSatish Balay } else if (NTR_KSP_STCG == tr->ksp_type) { 506a7e14dcfSSatish Balay ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 507a7e14dcfSSatish Balay } else { /* gltr */ 508a7e14dcfSSatish Balay ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 509a7e14dcfSSatish Balay } 510a7e14dcfSSatish Balay 511a7e14dcfSSatish Balay if (prered >= 0.0) { 512a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot 513a7e14dcfSSatish Balay happen in infinite precision arithmetic. Step should 514a7e14dcfSSatish Balay be rejected! */ 515a7e14dcfSSatish Balay tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d); 516a7e14dcfSSatish Balay } 517a7e14dcfSSatish Balay else { 518a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr); 519a7e14dcfSSatish Balay ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr); 520a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr); 521a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 522a7e14dcfSSatish Balay tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d); 523a7e14dcfSSatish Balay } 524a7e14dcfSSatish Balay else { 525a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &beta);CHKERRQ(ierr); 526a7e14dcfSSatish Balay actred = f - ftrial; 527a7e14dcfSSatish Balay prered = -prered; 528a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= tr->epsilon) && 529a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= tr->epsilon)) { 530a7e14dcfSSatish Balay kappa = 1.0; 531a7e14dcfSSatish Balay } 532a7e14dcfSSatish Balay else { 533a7e14dcfSSatish Balay kappa = actred / prered; 534a7e14dcfSSatish Balay } 535a7e14dcfSSatish Balay 536a7e14dcfSSatish Balay tau_1 = tr->theta * beta / (tr->theta * beta - (1.0 - tr->theta) * prered + actred); 537a7e14dcfSSatish Balay tau_2 = tr->theta * beta / (tr->theta * beta + (1.0 + tr->theta) * prered - actred); 538a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 539a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 540a7e14dcfSSatish Balay 541a7e14dcfSSatish Balay if (kappa >= 1.0 - tr->mu1) { 542a7e14dcfSSatish Balay /* Great agreement; accept step and update radius */ 543a7e14dcfSSatish Balay if (tau_max < 1.0) { 544a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d); 545a7e14dcfSSatish Balay } 546a7e14dcfSSatish Balay else if (tau_max > tr->gamma4) { 547a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->gamma4 * norm_d); 548a7e14dcfSSatish Balay } 549a7e14dcfSSatish Balay else { 550a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 551a7e14dcfSSatish Balay } 552a7e14dcfSSatish Balay break; 553a7e14dcfSSatish Balay } 554a7e14dcfSSatish Balay else if (kappa >= 1.0 - tr->mu2) { 555a7e14dcfSSatish Balay /* Good agreement */ 556a7e14dcfSSatish Balay 557a7e14dcfSSatish Balay if (tau_max < tr->gamma2) { 558a7e14dcfSSatish Balay tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d); 559a7e14dcfSSatish Balay } 560a7e14dcfSSatish Balay else if (tau_max > tr->gamma3) { 561a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d); 562a7e14dcfSSatish Balay } 563a7e14dcfSSatish Balay else if (tau_max < 1.0) { 564a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 565a7e14dcfSSatish Balay } 566a7e14dcfSSatish Balay else { 567a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 568a7e14dcfSSatish Balay } 569a7e14dcfSSatish Balay break; 570a7e14dcfSSatish Balay } 571a7e14dcfSSatish Balay else { 572a7e14dcfSSatish Balay /* Not good agreement */ 573a7e14dcfSSatish Balay if (tau_min > 1.0) { 574a7e14dcfSSatish Balay tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d); 575a7e14dcfSSatish Balay } 576a7e14dcfSSatish Balay else if (tau_max < tr->gamma1) { 577a7e14dcfSSatish Balay tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d); 578a7e14dcfSSatish Balay } 579a7e14dcfSSatish Balay else if ((tau_min < tr->gamma1) && (tau_max >= 1.0)) { 580a7e14dcfSSatish Balay tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d); 581a7e14dcfSSatish Balay } 582a7e14dcfSSatish Balay else if ((tau_1 >= tr->gamma1) && (tau_1 < 1.0) && 583a7e14dcfSSatish Balay ((tau_2 < tr->gamma1) || (tau_2 >= 1.0))) { 584a7e14dcfSSatish Balay tao->trust = tau_1 * PetscMin(tao->trust, norm_d); 585a7e14dcfSSatish Balay } 586a7e14dcfSSatish Balay else if ((tau_2 >= tr->gamma1) && (tau_2 < 1.0) && 587a7e14dcfSSatish Balay ((tau_1 < tr->gamma1) || (tau_2 >= 1.0))) { 588a7e14dcfSSatish Balay tao->trust = tau_2 * PetscMin(tao->trust, norm_d); 589a7e14dcfSSatish Balay } 590a7e14dcfSSatish Balay else { 591a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 592a7e14dcfSSatish Balay } 593a7e14dcfSSatish Balay } 594a7e14dcfSSatish Balay } 595a7e14dcfSSatish Balay } 596a7e14dcfSSatish Balay } 597a7e14dcfSSatish Balay 598a7e14dcfSSatish Balay /* The step computed was not good and the radius was decreased. 599a7e14dcfSSatish Balay Monitor the radius to terminate. */ 600a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, tao->trust, &reason);CHKERRQ(ierr); 601a7e14dcfSSatish Balay } 602a7e14dcfSSatish Balay 603a7e14dcfSSatish Balay /* The radius may have been increased; modify if it is too large */ 604a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tr->max_radius); 605a7e14dcfSSatish Balay 606a7e14dcfSSatish Balay if (reason == TAO_CONTINUE_ITERATING) { 607a7e14dcfSSatish Balay ierr = VecCopy(tr->W, tao->solution);CHKERRQ(ierr); 608a7e14dcfSSatish Balay f = ftrial; 609a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao, tao->solution, tao->gradient); 610a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 61153506e15SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 612a7e14dcfSSatish Balay needH = 1; 613a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, tao->trust, &reason);CHKERRQ(ierr); 614a7e14dcfSSatish Balay } 615a7e14dcfSSatish Balay } 616a7e14dcfSSatish Balay PetscFunctionReturn(0); 617a7e14dcfSSatish Balay } 618a7e14dcfSSatish Balay 619a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 620a7e14dcfSSatish Balay #undef __FUNCT__ 621a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_NTR" 622441846f8SBarry Smith static PetscErrorCode TaoSetUp_NTR(Tao tao) 623a7e14dcfSSatish Balay { 624a7e14dcfSSatish Balay TAO_NTR *tr = (TAO_NTR *)tao->data; 625a7e14dcfSSatish Balay PetscErrorCode ierr; 626a7e14dcfSSatish Balay 627a7e14dcfSSatish Balay PetscFunctionBegin; 628a7e14dcfSSatish Balay 629a7e14dcfSSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);} 630a7e14dcfSSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);} 631a7e14dcfSSatish Balay if (!tr->W) {ierr = VecDuplicate(tao->solution, &tr->W);CHKERRQ(ierr);} 632a7e14dcfSSatish Balay 633a7e14dcfSSatish Balay tr->Diag = 0; 634a7e14dcfSSatish Balay tr->M = 0; 635a7e14dcfSSatish Balay 636a7e14dcfSSatish Balay 637a7e14dcfSSatish Balay PetscFunctionReturn(0); 638a7e14dcfSSatish Balay } 639a7e14dcfSSatish Balay 640a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 641a7e14dcfSSatish Balay #undef __FUNCT__ 642a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_NTR" 643441846f8SBarry Smith static PetscErrorCode TaoDestroy_NTR(Tao tao) 644a7e14dcfSSatish Balay { 645a7e14dcfSSatish Balay TAO_NTR *tr = (TAO_NTR *)tao->data; 646a7e14dcfSSatish Balay PetscErrorCode ierr; 647a7e14dcfSSatish Balay 648a7e14dcfSSatish Balay PetscFunctionBegin; 649a7e14dcfSSatish Balay if (tao->setupcalled) { 650a7e14dcfSSatish Balay ierr = VecDestroy(&tr->W);CHKERRQ(ierr); 651a7e14dcfSSatish Balay } 652a7e14dcfSSatish Balay ierr = MatDestroy(&tr->M);CHKERRQ(ierr); 653a7e14dcfSSatish Balay ierr = VecDestroy(&tr->Diag);CHKERRQ(ierr); 654a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 655a7e14dcfSSatish Balay PetscFunctionReturn(0); 656a7e14dcfSSatish Balay } 657a7e14dcfSSatish Balay 658a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 659a7e14dcfSSatish Balay #undef __FUNCT__ 660a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_NTR" 661*1a1499c8SBarry Smith static PetscErrorCode TaoSetFromOptions_NTR(PetscOptionsObjectType *PetscOptionsObject,Tao tao) 662a7e14dcfSSatish Balay { 663a7e14dcfSSatish Balay TAO_NTR *tr = (TAO_NTR *)tao->data; 664a7e14dcfSSatish Balay PetscErrorCode ierr; 665a7e14dcfSSatish Balay 666a7e14dcfSSatish Balay PetscFunctionBegin; 667*1a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Newton trust region method for unconstrained optimization");CHKERRQ(ierr); 66894ae4db5SBarry Smith ierr = PetscOptionsEList("-tao_ntr_ksp_type", "ksp type", "", NTR_KSP, NTR_KSP_TYPES, NTR_KSP[tr->ksp_type], &tr->ksp_type,NULL);CHKERRQ(ierr); 66994ae4db5SBarry Smith ierr = PetscOptionsEList("-tao_ntr_pc_type", "pc type", "", NTR_PC, NTR_PC_TYPES, NTR_PC[tr->pc_type], &tr->pc_type,NULL);CHKERRQ(ierr); 67094ae4db5SBarry Smith ierr = PetscOptionsEList("-tao_ntr_bfgs_scale_type", "bfgs scale type", "", BFGS_SCALE, BFGS_SCALE_TYPES, BFGS_SCALE[tr->bfgs_scale_type], &tr->bfgs_scale_type,NULL);CHKERRQ(ierr); 67194ae4db5SBarry Smith ierr = PetscOptionsEList("-tao_ntr_init_type", "tao->trust initialization type", "", NTR_INIT, NTR_INIT_TYPES, NTR_INIT[tr->init_type], &tr->init_type,NULL);CHKERRQ(ierr); 67294ae4db5SBarry Smith ierr = PetscOptionsEList("-tao_ntr_update_type", "radius update type", "", NTR_UPDATE, NTR_UPDATE_TYPES, NTR_UPDATE[tr->update_type], &tr->update_type,NULL);CHKERRQ(ierr); 67394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_eta1", "step is unsuccessful if actual reduction < eta1 * predicted reduction", "", tr->eta1, &tr->eta1,NULL);CHKERRQ(ierr); 67494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_eta2", "", "", tr->eta2, &tr->eta2,NULL);CHKERRQ(ierr); 67594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_eta3", "", "", tr->eta3, &tr->eta3,NULL);CHKERRQ(ierr); 67694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_eta4", "", "", tr->eta4, &tr->eta4,NULL);CHKERRQ(ierr); 67794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_alpha1", "", "", tr->alpha1, &tr->alpha1,NULL);CHKERRQ(ierr); 67894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_alpha2", "", "", tr->alpha2, &tr->alpha2,NULL);CHKERRQ(ierr); 67994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_alpha3", "", "", tr->alpha3, &tr->alpha3,NULL);CHKERRQ(ierr); 68094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_alpha4", "", "", tr->alpha4, &tr->alpha4,NULL);CHKERRQ(ierr); 68194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_alpha5", "", "", tr->alpha5, &tr->alpha5,NULL);CHKERRQ(ierr); 68294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_mu1", "", "", tr->mu1, &tr->mu1,NULL);CHKERRQ(ierr); 68394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_mu2", "", "", tr->mu2, &tr->mu2,NULL);CHKERRQ(ierr); 68494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_gamma1", "", "", tr->gamma1, &tr->gamma1,NULL);CHKERRQ(ierr); 68594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_gamma2", "", "", tr->gamma2, &tr->gamma2,NULL);CHKERRQ(ierr); 68694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_gamma3", "", "", tr->gamma3, &tr->gamma3,NULL);CHKERRQ(ierr); 68794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_gamma4", "", "", tr->gamma4, &tr->gamma4,NULL);CHKERRQ(ierr); 68894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_theta", "", "", tr->theta, &tr->theta,NULL);CHKERRQ(ierr); 68994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_mu1_i", "", "", tr->mu1_i, &tr->mu1_i,NULL);CHKERRQ(ierr); 69094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_mu2_i", "", "", tr->mu2_i, &tr->mu2_i,NULL);CHKERRQ(ierr); 69194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_gamma1_i", "", "", tr->gamma1_i, &tr->gamma1_i,NULL);CHKERRQ(ierr); 69294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_gamma2_i", "", "", tr->gamma2_i, &tr->gamma2_i,NULL);CHKERRQ(ierr); 69394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_gamma3_i", "", "", tr->gamma3_i, &tr->gamma3_i,NULL);CHKERRQ(ierr); 69494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_gamma4_i", "", "", tr->gamma4_i, &tr->gamma4_i,NULL);CHKERRQ(ierr); 69594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_theta_i", "", "", tr->theta_i, &tr->theta_i,NULL);CHKERRQ(ierr); 69694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_min_radius", "lower bound on initial trust-region radius", "", tr->min_radius, &tr->min_radius,NULL);CHKERRQ(ierr); 69794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_max_radius", "upper bound on trust-region radius", "", tr->max_radius, &tr->max_radius,NULL);CHKERRQ(ierr); 69894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ntr_epsilon", "tolerance used when computing actual and predicted reduction", "", tr->epsilon, &tr->epsilon,NULL);CHKERRQ(ierr); 699a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 700a7e14dcfSSatish Balay ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 701a7e14dcfSSatish Balay PetscFunctionReturn(0); 702a7e14dcfSSatish Balay } 703a7e14dcfSSatish Balay 704a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 705a7e14dcfSSatish Balay #undef __FUNCT__ 706a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_NTR" 707441846f8SBarry Smith static PetscErrorCode TaoView_NTR(Tao tao, PetscViewer viewer) 708a7e14dcfSSatish Balay { 709a7e14dcfSSatish Balay TAO_NTR *tr = (TAO_NTR *)tao->data; 710a7e14dcfSSatish Balay PetscErrorCode ierr; 711a7e14dcfSSatish Balay PetscInt nrejects; 712a7e14dcfSSatish Balay PetscBool isascii; 71353506e15SBarry Smith 714a7e14dcfSSatish Balay PetscFunctionBegin; 715a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 716a7e14dcfSSatish Balay if (isascii) { 717a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 718a7e14dcfSSatish Balay if (NTR_PC_BFGS == tr->pc_type && tr->M) { 719a7e14dcfSSatish Balay ierr = MatLMVMGetRejects(tr->M, &nrejects);CHKERRQ(ierr); 720a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n", nrejects);CHKERRQ(ierr); 721a7e14dcfSSatish Balay } 722a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 723a7e14dcfSSatish Balay } 724a7e14dcfSSatish Balay PetscFunctionReturn(0); 725a7e14dcfSSatish Balay } 726a7e14dcfSSatish Balay 727a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 7281522df2eSJason Sarich /*MC 7291522df2eSJason Sarich TAONTR - Newton's method with trust region for unconstrained minimization. 7301522df2eSJason Sarich At each iteration, the Newton trust region method solves the system. 7311522df2eSJason Sarich NTR expects a KSP solver with a trust region radius. 7321522df2eSJason Sarich min_d .5 dT Hk d + gkT d, s.t. ||d|| < Delta_k 7331522df2eSJason Sarich 7341522df2eSJason Sarich Options Database Keys: 7351522df2eSJason Sarich + -tao_ntr_ksp_type - "nash","stcg","gltr" 7361522df2eSJason Sarich . -tao_ntr_pc_type - "none","ahess","bfgs","petsc" 7371522df2eSJason Sarich . -tao_ntr_bfgs_scale_type - type of scaling with bfgs pc, "ahess" or "bfgs" 7381522df2eSJason Sarich . -tao_ntr_init_type - "constant","direction","interpolation" 7391522df2eSJason Sarich . -tao_ntr_update_type - "reduction","interpolation" 7401522df2eSJason Sarich . -tao_ntr_min_radius - lower bound on trust region radius 7411522df2eSJason Sarich . -tao_ntr_max_radius - upper bound on trust region radius 7421522df2eSJason Sarich . -tao_ntr_epsilon - tolerance for accepting actual / predicted reduction 7431522df2eSJason Sarich . -tao_ntr_mu1_i - mu1 interpolation init factor 7441522df2eSJason Sarich . -tao_ntr_mu2_i - mu2 interpolation init factor 7451522df2eSJason Sarich . -tao_ntr_gamma1_i - gamma1 interpolation init factor 7461522df2eSJason Sarich . -tao_ntr_gamma2_i - gamma2 interpolation init factor 7471522df2eSJason Sarich . -tao_ntr_gamma3_i - gamma3 interpolation init factor 7481522df2eSJason Sarich . -tao_ntr_gamma4_i - gamma4 interpolation init factor 7491522df2eSJason Sarich . -tao_ntr_theta_i - thetha1 interpolation init factor 7501522df2eSJason Sarich . -tao_ntr_eta1 - eta1 reduction update factor 7511522df2eSJason Sarich . -tao_ntr_eta2 - eta2 reduction update factor 7521522df2eSJason Sarich . -tao_ntr_eta3 - eta3 reduction update factor 7531522df2eSJason Sarich . -tao_ntr_eta4 - eta4 reduction update factor 7541522df2eSJason Sarich . -tao_ntr_alpha1 - alpha1 reduction update factor 7551522df2eSJason Sarich . -tao_ntr_alpha2 - alpha2 reduction update factor 7561522df2eSJason Sarich . -tao_ntr_alpha3 - alpha3 reduction update factor 7571522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor 7581522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor 7591522df2eSJason Sarich . -tao_ntr_mu1 - mu1 interpolation update 7601522df2eSJason Sarich . -tao_ntr_mu2 - mu2 interpolation update 7611522df2eSJason Sarich . -tao_ntr_gamma1 - gamma1 interpolcation update 7621522df2eSJason Sarich . -tao_ntr_gamma2 - gamma2 interpolcation update 7631522df2eSJason Sarich . -tao_ntr_gamma3 - gamma3 interpolcation update 7641522df2eSJason Sarich . -tao_ntr_gamma4 - gamma4 interpolation update 7651522df2eSJason Sarich - -tao_ntr_theta - theta interpolation update 7661522df2eSJason Sarich 7671eb8069cSJason Sarich Level: beginner 7681522df2eSJason Sarich M*/ 7691522df2eSJason Sarich 770a7e14dcfSSatish Balay #undef __FUNCT__ 771a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_NTR" 772728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NTR(Tao tao) 773a7e14dcfSSatish Balay { 774a7e14dcfSSatish Balay TAO_NTR *tr; 775a7e14dcfSSatish Balay PetscErrorCode ierr; 776a7e14dcfSSatish Balay 777a7e14dcfSSatish Balay PetscFunctionBegin; 778a7e14dcfSSatish Balay 7793c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&tr);CHKERRQ(ierr); 780a7e14dcfSSatish Balay 781a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_NTR; 782a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_NTR; 783a7e14dcfSSatish Balay tao->ops->view = TaoView_NTR; 784a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_NTR; 785a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_NTR; 786a7e14dcfSSatish Balay 787a7e14dcfSSatish Balay tao->max_it = 50; 7886f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 7896f4723b1SBarry Smith tao->fatol = 1e-5; 7906f4723b1SBarry Smith tao->frtol = 1e-5; 7916f4723b1SBarry Smith #else 792a7e14dcfSSatish Balay tao->fatol = 1e-10; 793a7e14dcfSSatish Balay tao->frtol = 1e-10; 7946f4723b1SBarry Smith #endif 795a7e14dcfSSatish Balay tao->data = (void*)tr; 796a7e14dcfSSatish Balay 797a7e14dcfSSatish Balay tao->trust0 = 100.0; 798a7e14dcfSSatish Balay 799a7e14dcfSSatish Balay /* Standard trust region update parameters */ 800a7e14dcfSSatish Balay tr->eta1 = 1.0e-4; 801a7e14dcfSSatish Balay tr->eta2 = 0.25; 802a7e14dcfSSatish Balay tr->eta3 = 0.50; 803a7e14dcfSSatish Balay tr->eta4 = 0.90; 804a7e14dcfSSatish Balay 805a7e14dcfSSatish Balay tr->alpha1 = 0.25; 806a7e14dcfSSatish Balay tr->alpha2 = 0.50; 807a7e14dcfSSatish Balay tr->alpha3 = 1.00; 808a7e14dcfSSatish Balay tr->alpha4 = 2.00; 809a7e14dcfSSatish Balay tr->alpha5 = 4.00; 810a7e14dcfSSatish Balay 811a7e14dcfSSatish Balay /* Interpolation parameters */ 812a7e14dcfSSatish Balay tr->mu1_i = 0.35; 813a7e14dcfSSatish Balay tr->mu2_i = 0.50; 814a7e14dcfSSatish Balay 815a7e14dcfSSatish Balay tr->gamma1_i = 0.0625; 816a7e14dcfSSatish Balay tr->gamma2_i = 0.50; 817a7e14dcfSSatish Balay tr->gamma3_i = 2.00; 818a7e14dcfSSatish Balay tr->gamma4_i = 5.00; 819a7e14dcfSSatish Balay 820a7e14dcfSSatish Balay tr->theta_i = 0.25; 821a7e14dcfSSatish Balay 822a7e14dcfSSatish Balay /* Interpolation trust region update parameters */ 823a7e14dcfSSatish Balay tr->mu1 = 0.10; 824a7e14dcfSSatish Balay tr->mu2 = 0.50; 825a7e14dcfSSatish Balay 826a7e14dcfSSatish Balay tr->gamma1 = 0.25; 827a7e14dcfSSatish Balay tr->gamma2 = 0.50; 828a7e14dcfSSatish Balay tr->gamma3 = 2.00; 829a7e14dcfSSatish Balay tr->gamma4 = 4.00; 830a7e14dcfSSatish Balay 831a7e14dcfSSatish Balay tr->theta = 0.05; 832a7e14dcfSSatish Balay 833a7e14dcfSSatish Balay tr->min_radius = 1.0e-10; 834a7e14dcfSSatish Balay tr->max_radius = 1.0e10; 835a7e14dcfSSatish Balay tr->epsilon = 1.0e-6; 836a7e14dcfSSatish Balay 837a7e14dcfSSatish Balay tr->ksp_type = NTR_KSP_STCG; 838a7e14dcfSSatish Balay tr->pc_type = NTR_PC_BFGS; 839a7e14dcfSSatish Balay tr->bfgs_scale_type = BFGS_SCALE_AHESS; 840a7e14dcfSSatish Balay tr->init_type = NTR_INIT_INTERPOLATION; 841a7e14dcfSSatish Balay tr->update_type = NTR_UPDATE_REDUCTION; 842a7e14dcfSSatish Balay 843a7e14dcfSSatish Balay 844a7e14dcfSSatish Balay /* Set linear solver to default for trust region */ 845a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr); 846a7e14dcfSSatish Balay PetscFunctionReturn(0); 847a7e14dcfSSatish Balay } 848a7e14dcfSSatish Balay 849a7e14dcfSSatish Balay 850a7e14dcfSSatish Balay #undef __FUNCT__ 851a7e14dcfSSatish Balay #define __FUNCT__ "MatLMVMSolveShell" 852a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x) 853a7e14dcfSSatish Balay { 854a7e14dcfSSatish Balay PetscErrorCode ierr; 855a7e14dcfSSatish Balay Mat M; 856a7e14dcfSSatish Balay PetscFunctionBegin; 857a7e14dcfSSatish Balay PetscValidHeaderSpecific(pc,PC_CLASSID,1); 858a7e14dcfSSatish Balay PetscValidHeaderSpecific(b,VEC_CLASSID,2); 859a7e14dcfSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,3); 860a7e14dcfSSatish Balay ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr); 861a7e14dcfSSatish Balay ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr); 862a7e14dcfSSatish Balay PetscFunctionReturn(0); 863a7e14dcfSSatish Balay } 864