1ba92ff59SBarry Smith #include <petsctaolinesearch.h> 21daac38eSTodd Munson #include <../src/tao/unconstrained/impls/nls/nlsimpl.h> 3a7e14dcfSSatish Balay 4aaa7dc30SBarry Smith #include <petscksp.h> 5a7e14dcfSSatish Balay 6a7e14dcfSSatish Balay #define NLS_INIT_CONSTANT 0 7a7e14dcfSSatish Balay #define NLS_INIT_DIRECTION 1 8a7e14dcfSSatish Balay #define NLS_INIT_INTERPOLATION 2 9a7e14dcfSSatish Balay #define NLS_INIT_TYPES 3 10a7e14dcfSSatish Balay 11a7e14dcfSSatish Balay #define NLS_UPDATE_STEP 0 12a7e14dcfSSatish Balay #define NLS_UPDATE_REDUCTION 1 13a7e14dcfSSatish Balay #define NLS_UPDATE_INTERPOLATION 2 14a7e14dcfSSatish Balay #define NLS_UPDATE_TYPES 3 15a7e14dcfSSatish Balay 1687f595a5SBarry Smith static const char *NLS_INIT[64] = {"constant", "direction", "interpolation"}; 17a7e14dcfSSatish Balay 1887f595a5SBarry Smith static const char *NLS_UPDATE[64] = {"step", "reduction", "interpolation"}; 19a7e14dcfSSatish Balay 201daac38eSTodd Munson /* 21a7e14dcfSSatish Balay Implements Newton's Method with a line search approach for solving 22a7e14dcfSSatish Balay unconstrained minimization problems. A More'-Thuente line search 23a7e14dcfSSatish Balay is used to guarantee that the bfgs preconditioner remains positive 24a7e14dcfSSatish Balay definite. 25a7e14dcfSSatish Balay 26a7e14dcfSSatish Balay The method can shift the Hessian matrix. The shifting procedure is 27a7e14dcfSSatish Balay adapted from the PATH algorithm for solving complementarity 28a7e14dcfSSatish Balay problems. 29a7e14dcfSSatish Balay 30a7e14dcfSSatish Balay The linear system solve should be done with a conjugate gradient 311daac38eSTodd Munson method, although any method can be used. 321daac38eSTodd Munson */ 33a7e14dcfSSatish Balay 34a7e14dcfSSatish Balay #define NLS_NEWTON 0 35a7e14dcfSSatish Balay #define NLS_BFGS 1 360c51296cSAlp Dener #define NLS_GRADIENT 2 37a7e14dcfSSatish Balay 38441846f8SBarry Smith static PetscErrorCode TaoSolve_NLS(Tao tao) 39a7e14dcfSSatish Balay { 40a7e14dcfSSatish Balay PetscErrorCode ierr; 41a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 421daac38eSTodd Munson KSPType ksp_type; 430ad3a497SAlp Dener PetscBool is_nash,is_stcg,is_gltr,is_bfgs,is_jacobi,is_symmetric,sym_set; 44a7e14dcfSSatish Balay KSPConvergedReason ksp_reason; 451daac38eSTodd Munson PC pc; 461daac38eSTodd Munson TaoLineSearchConvergedReason ls_reason; 47a7e14dcfSSatish Balay 48a7e14dcfSSatish Balay PetscReal fmin, ftrial, f_full, prered, actred, kappa, sigma; 49a7e14dcfSSatish Balay PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 50a7e14dcfSSatish Balay PetscReal f, fold, gdx, gnorm, pert; 51a7e14dcfSSatish Balay PetscReal step = 1.0; 52a7e14dcfSSatish Balay PetscReal norm_d = 0.0, e_min; 53a7e14dcfSSatish Balay 54a7e14dcfSSatish Balay PetscInt stepType; 55a7e14dcfSSatish Balay PetscInt bfgsUpdates = 0; 56a7e14dcfSSatish Balay PetscInt n,N,kspits; 57b4690188SBarry Smith PetscInt needH = 1; 58a7e14dcfSSatish Balay 59a7e14dcfSSatish Balay PetscInt i_max = 5; 60a7e14dcfSSatish Balay PetscInt j_max = 1; 61a7e14dcfSSatish Balay PetscInt i, j; 62a7e14dcfSSatish Balay 63a7e14dcfSSatish Balay PetscFunctionBegin; 64a7e14dcfSSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 65a7e14dcfSSatish Balay ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by nls algorithm\n");CHKERRQ(ierr); 66a7e14dcfSSatish Balay } 67a7e14dcfSSatish Balay 68a7e14dcfSSatish Balay /* Initialized variables */ 69a7e14dcfSSatish Balay pert = nlsP->sval; 70a7e14dcfSSatish Balay 712d9aa51bSJason Sarich /* Number of times ksp stopped because of these reasons */ 72a7e14dcfSSatish Balay nlsP->ksp_atol = 0; 73a7e14dcfSSatish Balay nlsP->ksp_rtol = 0; 74a7e14dcfSSatish Balay nlsP->ksp_dtol = 0; 75a7e14dcfSSatish Balay nlsP->ksp_ctol = 0; 76a7e14dcfSSatish Balay nlsP->ksp_negc = 0; 77a7e14dcfSSatish Balay nlsP->ksp_iter = 0; 78a7e14dcfSSatish Balay nlsP->ksp_othr = 0; 79a7e14dcfSSatish Balay 80a7e14dcfSSatish Balay /* Initialize trust-region radius when using nash, stcg, or gltr 81ba7fe8fbSTodd Munson Command automatically ignored for other methods 82ba7fe8fbSTodd Munson Will be reset during the first iteration 83ba7fe8fbSTodd Munson */ 841daac38eSTodd Munson ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr); 851daac38eSTodd Munson ierr = PetscStrcmp(ksp_type,KSPCGNASH,&is_nash);CHKERRQ(ierr); 861daac38eSTodd Munson ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&is_stcg);CHKERRQ(ierr); 871daac38eSTodd Munson ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&is_gltr);CHKERRQ(ierr); 881daac38eSTodd Munson 89ba7fe8fbSTodd Munson ierr = KSPCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 90a7e14dcfSSatish Balay 911daac38eSTodd Munson if (is_nash || is_stcg || is_gltr) { 921daac38eSTodd Munson if (tao->trust0 < 0.0) SETERRQ(PETSC_COMM_SELF,1,"Initial radius negative"); 93a7e14dcfSSatish Balay tao->trust = tao->trust0; 94a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 95a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 96a7e14dcfSSatish Balay } 97a7e14dcfSSatish Balay 98a7e14dcfSSatish Balay /* Check convergence criteria */ 99a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 100a9603a14SPatrick Farrell ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 10187f595a5SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 102a7e14dcfSSatish Balay 1033ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 1043ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 1053ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 1063ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 1073ecd9318SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 108a7e14dcfSSatish Balay 1090c51296cSAlp Dener /* Allocate the vectors needed for the BFGS approximation */ 110a7e14dcfSSatish Balay ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 1110c51296cSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs);CHKERRQ(ierr); 1120c51296cSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi);CHKERRQ(ierr); 1130c51296cSAlp Dener if (is_bfgs) { 1140c51296cSAlp Dener nlsP->bfgs_pre = pc; 1150c51296cSAlp Dener ierr = PCLMVMGetMatLMVM(nlsP->bfgs_pre, &nlsP->M);CHKERRQ(ierr); 1160c51296cSAlp Dener ierr = VecGetLocalSize(tao->solution, &n);CHKERRQ(ierr); 1170c51296cSAlp Dener ierr = VecGetSize(tao->solution, &N);CHKERRQ(ierr); 1180c51296cSAlp Dener ierr = MatSetSizes(nlsP->M, n, n, N, N);CHKERRQ(ierr); 1190c51296cSAlp Dener ierr = MatLMVMAllocate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 1200ad3a497SAlp Dener ierr = MatIsSymmetricKnown(nlsP->M, &sym_set, &is_symmetric);CHKERRQ(ierr); 1210ad3a497SAlp Dener if (!sym_set || !is_symmetric) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric."); 1220c51296cSAlp Dener } else if (is_jacobi) { 123baa89ecbSBarry Smith ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 124a7e14dcfSSatish Balay } 125a7e14dcfSSatish Balay 126a7e14dcfSSatish Balay /* Initialize trust-region radius. The initialization is only performed 127a7e14dcfSSatish Balay when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */ 1281daac38eSTodd Munson if (is_nash || is_stcg || is_gltr) { 129a7e14dcfSSatish Balay switch(nlsP->init_type) { 130a7e14dcfSSatish Balay case NLS_INIT_CONSTANT: 131a7e14dcfSSatish Balay /* Use the initial radius specified */ 132a7e14dcfSSatish Balay break; 133a7e14dcfSSatish Balay 134a7e14dcfSSatish Balay case NLS_INIT_INTERPOLATION: 135a7e14dcfSSatish Balay /* Use the initial radius specified */ 136a7e14dcfSSatish Balay max_radius = 0.0; 137a7e14dcfSSatish Balay 138a7e14dcfSSatish Balay for (j = 0; j < j_max; ++j) { 139a7e14dcfSSatish Balay fmin = f; 140a7e14dcfSSatish Balay sigma = 0.0; 141a7e14dcfSSatish Balay 142a7e14dcfSSatish Balay if (needH) { 143ffad9901SBarry Smith ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 144a7e14dcfSSatish Balay needH = 0; 145a7e14dcfSSatish Balay } 146a7e14dcfSSatish Balay 147a7e14dcfSSatish Balay for (i = 0; i < i_max; ++i) { 148a7e14dcfSSatish Balay ierr = VecCopy(tao->solution,nlsP->W);CHKERRQ(ierr); 149a7e14dcfSSatish Balay ierr = VecAXPY(nlsP->W,-tao->trust/gnorm,tao->gradient);CHKERRQ(ierr); 150a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, nlsP->W, &ftrial);CHKERRQ(ierr); 151a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 152a7e14dcfSSatish Balay tau = nlsP->gamma1_i; 15387f595a5SBarry Smith } else { 154a7e14dcfSSatish Balay if (ftrial < fmin) { 155a7e14dcfSSatish Balay fmin = ftrial; 156a7e14dcfSSatish Balay sigma = -tao->trust / gnorm; 157a7e14dcfSSatish Balay } 158a7e14dcfSSatish Balay 159a7e14dcfSSatish Balay ierr = MatMult(tao->hessian, tao->gradient, nlsP->D);CHKERRQ(ierr); 160a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, nlsP->D, &prered);CHKERRQ(ierr); 161a7e14dcfSSatish Balay 162a7e14dcfSSatish Balay prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm)); 163a7e14dcfSSatish Balay actred = f - ftrial; 16487f595a5SBarry Smith if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) { 165a7e14dcfSSatish Balay kappa = 1.0; 16687f595a5SBarry Smith } else { 167a7e14dcfSSatish Balay kappa = actred / prered; 168a7e14dcfSSatish Balay } 169a7e14dcfSSatish Balay 170a7e14dcfSSatish Balay tau_1 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred); 171a7e14dcfSSatish Balay tau_2 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred); 172a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 173a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 174a7e14dcfSSatish Balay 175a7e14dcfSSatish Balay if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu1_i) { 176a7e14dcfSSatish Balay /* Great agreement */ 177a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 178a7e14dcfSSatish Balay 179a7e14dcfSSatish Balay if (tau_max < 1.0) { 180a7e14dcfSSatish Balay tau = nlsP->gamma3_i; 18187f595a5SBarry Smith } else if (tau_max > nlsP->gamma4_i) { 182a7e14dcfSSatish Balay tau = nlsP->gamma4_i; 18387f595a5SBarry Smith } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) { 184a7e14dcfSSatish Balay tau = tau_1; 18587f595a5SBarry Smith } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) { 186a7e14dcfSSatish Balay tau = tau_2; 18787f595a5SBarry Smith } else { 188a7e14dcfSSatish Balay tau = tau_max; 189a7e14dcfSSatish Balay } 19087f595a5SBarry Smith } else if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu2_i) { 191a7e14dcfSSatish Balay /* Good agreement */ 192a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 193a7e14dcfSSatish Balay 194a7e14dcfSSatish Balay if (tau_max < nlsP->gamma2_i) { 195a7e14dcfSSatish Balay tau = nlsP->gamma2_i; 19687f595a5SBarry Smith } else if (tau_max > nlsP->gamma3_i) { 197a7e14dcfSSatish Balay tau = nlsP->gamma3_i; 19887f595a5SBarry Smith } else { 199a7e14dcfSSatish Balay tau = tau_max; 200a7e14dcfSSatish Balay } 20187f595a5SBarry Smith } else { 202a7e14dcfSSatish Balay /* Not good agreement */ 203a7e14dcfSSatish Balay if (tau_min > 1.0) { 204a7e14dcfSSatish Balay tau = nlsP->gamma2_i; 20587f595a5SBarry Smith } else if (tau_max < nlsP->gamma1_i) { 206a7e14dcfSSatish Balay tau = nlsP->gamma1_i; 20787f595a5SBarry Smith } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) { 208a7e14dcfSSatish Balay tau = nlsP->gamma1_i; 20987f595a5SBarry Smith } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) { 210a7e14dcfSSatish Balay tau = tau_1; 21187f595a5SBarry Smith } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) { 212a7e14dcfSSatish Balay tau = tau_2; 21387f595a5SBarry Smith } else { 214a7e14dcfSSatish Balay tau = tau_max; 215a7e14dcfSSatish Balay } 216a7e14dcfSSatish Balay } 217a7e14dcfSSatish Balay } 218a7e14dcfSSatish Balay tao->trust = tau * tao->trust; 219a7e14dcfSSatish Balay } 220a7e14dcfSSatish Balay 221a7e14dcfSSatish Balay if (fmin < f) { 222a7e14dcfSSatish Balay f = fmin; 223a7e14dcfSSatish Balay ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr); 224a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao,tao->solution,tao->gradient);CHKERRQ(ierr); 225a7e14dcfSSatish Balay 226a9603a14SPatrick Farrell ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 22787f595a5SBarry Smith if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN"); 228a7e14dcfSSatish Balay needH = 1; 229a7e14dcfSSatish Balay 2303ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 2313ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 2323ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 2333ecd9318SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 234a7e14dcfSSatish Balay } 235a7e14dcfSSatish Balay } 236a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, max_radius); 237a7e14dcfSSatish Balay 238a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 239a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 240a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 241a7e14dcfSSatish Balay break; 242a7e14dcfSSatish Balay 243a7e14dcfSSatish Balay default: 244a7e14dcfSSatish Balay /* Norm of the first direction will initialize radius */ 245a7e14dcfSSatish Balay tao->trust = 0.0; 246a7e14dcfSSatish Balay break; 247a7e14dcfSSatish Balay } 248a7e14dcfSSatish Balay } 249a7e14dcfSSatish Balay 250a7e14dcfSSatish Balay /* Set counter for gradient/reset steps*/ 251a7e14dcfSSatish Balay nlsP->newt = 0; 252a7e14dcfSSatish Balay nlsP->bfgs = 0; 253a7e14dcfSSatish Balay nlsP->grad = 0; 254a7e14dcfSSatish Balay 255a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 2563ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 2578931d482SJason Sarich ++tao->niter; 258ae93cb3cSJason Sarich tao->ksp_its=0; 259a7e14dcfSSatish Balay 260a7e14dcfSSatish Balay /* Compute the Hessian */ 261a7e14dcfSSatish Balay if (needH) { 262ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 263a7e14dcfSSatish Balay } 264a7e14dcfSSatish Balay 265a7e14dcfSSatish Balay /* Shift the Hessian matrix */ 266a7e14dcfSSatish Balay if (pert > 0) { 267302440fdSBarry Smith ierr = MatShift(tao->hessian, pert);CHKERRQ(ierr); 268a7e14dcfSSatish Balay if (tao->hessian != tao->hessian_pre) { 269a7e14dcfSSatish Balay ierr = MatShift(tao->hessian_pre, pert);CHKERRQ(ierr); 270a7e14dcfSSatish Balay } 271a7e14dcfSSatish Balay } 272a7e14dcfSSatish Balay 2730c51296cSAlp Dener if (nlsP->bfgs_pre) { 274a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 275a7e14dcfSSatish Balay ++bfgsUpdates; 276a7e14dcfSSatish Balay } 277a7e14dcfSSatish Balay 278a7e14dcfSSatish Balay /* Solve the Newton system of equations */ 27923ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 2801daac38eSTodd Munson if (is_nash || is_stcg || is_gltr) { 281ba7fe8fbSTodd Munson ierr = KSPCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 282a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr); 283a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 284a7e14dcfSSatish Balay tao->ksp_its+=kspits; 285ae93cb3cSJason Sarich tao->ksp_tot_its+=kspits; 286ba7fe8fbSTodd Munson ierr = KSPCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 287a7e14dcfSSatish Balay 288a7e14dcfSSatish Balay if (0.0 == tao->trust) { 289a7e14dcfSSatish Balay /* Radius was uninitialized; use the norm of the direction */ 290a7e14dcfSSatish Balay if (norm_d > 0.0) { 291a7e14dcfSSatish Balay tao->trust = norm_d; 292a7e14dcfSSatish Balay 293a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 294a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 295a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 29687f595a5SBarry Smith } else { 297a7e14dcfSSatish Balay /* The direction was bad; set radius to default value and re-solve 298a7e14dcfSSatish Balay the trust-region subproblem to get a direction */ 299a7e14dcfSSatish Balay tao->trust = tao->trust0; 300a7e14dcfSSatish Balay 301a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 302a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 303a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 304a7e14dcfSSatish Balay 305ba7fe8fbSTodd Munson ierr = KSPCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 306a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr); 307a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 308a7e14dcfSSatish Balay tao->ksp_its+=kspits; 309ae93cb3cSJason Sarich tao->ksp_tot_its+=kspits; 310ba7fe8fbSTodd Munson ierr = KSPCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 311a7e14dcfSSatish Balay 31287f595a5SBarry Smith if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 313a7e14dcfSSatish Balay } 314a7e14dcfSSatish Balay } 31587f595a5SBarry Smith } else { 316a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr); 317a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr); 318a7e14dcfSSatish Balay tao->ksp_its += kspits; 319ae93cb3cSJason Sarich tao->ksp_tot_its+=kspits; 320a7e14dcfSSatish Balay } 321a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 322a7e14dcfSSatish Balay ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 3230c51296cSAlp Dener if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (nlsP->bfgs_pre)) { 324a7e14dcfSSatish Balay /* Preconditioner is numerically indefinite; reset the 325a7e14dcfSSatish Balay approximate if using BFGS preconditioning. */ 326cd929ea3SAlp Dener ierr = MatLMVMReset(nlsP->M, PETSC_FALSE);CHKERRQ(ierr); 327a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 328a7e14dcfSSatish Balay bfgsUpdates = 1; 329a7e14dcfSSatish Balay } 330a7e14dcfSSatish Balay 331a7e14dcfSSatish Balay if (KSP_CONVERGED_ATOL == ksp_reason) { 332a7e14dcfSSatish Balay ++nlsP->ksp_atol; 33387f595a5SBarry Smith } else if (KSP_CONVERGED_RTOL == ksp_reason) { 334a7e14dcfSSatish Balay ++nlsP->ksp_rtol; 33587f595a5SBarry Smith } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) { 336a7e14dcfSSatish Balay ++nlsP->ksp_ctol; 33787f595a5SBarry Smith } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) { 338a7e14dcfSSatish Balay ++nlsP->ksp_negc; 33987f595a5SBarry Smith } else if (KSP_DIVERGED_DTOL == ksp_reason) { 340a7e14dcfSSatish Balay ++nlsP->ksp_dtol; 34187f595a5SBarry Smith } else if (KSP_DIVERGED_ITS == ksp_reason) { 342a7e14dcfSSatish Balay ++nlsP->ksp_iter; 34387f595a5SBarry Smith } else { 344a7e14dcfSSatish Balay ++nlsP->ksp_othr; 345a7e14dcfSSatish Balay } 346a7e14dcfSSatish Balay 347a7e14dcfSSatish Balay /* Check for success (descent direction) */ 348a7e14dcfSSatish Balay ierr = VecDot(nlsP->D, tao->gradient, &gdx);CHKERRQ(ierr); 349a7e14dcfSSatish Balay if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 350a7e14dcfSSatish Balay /* Newton step is not descent or direction produced Inf or NaN 351a7e14dcfSSatish Balay Update the perturbation for next time */ 352a7e14dcfSSatish Balay if (pert <= 0.0) { 353a7e14dcfSSatish Balay /* Initialize the perturbation */ 354a7e14dcfSSatish Balay pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 3551daac38eSTodd Munson if (is_gltr) { 356ba7fe8fbSTodd Munson ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 357a7e14dcfSSatish Balay pert = PetscMax(pert, -e_min); 358a7e14dcfSSatish Balay } 35987f595a5SBarry Smith } else { 360a7e14dcfSSatish Balay /* Increase the perturbation */ 361a7e14dcfSSatish Balay pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 362a7e14dcfSSatish Balay } 363a7e14dcfSSatish Balay 3640c51296cSAlp Dener if (!nlsP->bfgs_pre) { 365a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and updated 366a7e14dcfSSatish Balay Must use gradient direction in this case */ 367a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr); 368a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 369a7e14dcfSSatish Balay ++nlsP->grad; 370a7e14dcfSSatish Balay stepType = NLS_GRADIENT; 37187f595a5SBarry Smith } else { 372a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 373cd929ea3SAlp Dener ierr = MatSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 374a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 375a7e14dcfSSatish Balay 376a7e14dcfSSatish Balay /* Check for success (descent direction) */ 377a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, nlsP->D, &gdx);CHKERRQ(ierr); 378a7e14dcfSSatish Balay if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) { 379a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced not a number 380a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case because 381a7e14dcfSSatish Balay the first solve produces the scaled gradient direction, 382a7e14dcfSSatish Balay which is guaranteed to be descent */ 383a7e14dcfSSatish Balay 384a7e14dcfSSatish Balay /* Use steepest descent direction (scaled) */ 385cd929ea3SAlp Dener ierr = MatLMVMReset(nlsP->M, PETSC_FALSE);CHKERRQ(ierr); 386a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 387cd929ea3SAlp Dener ierr = MatSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 388a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 389a7e14dcfSSatish Balay 390a7e14dcfSSatish Balay bfgsUpdates = 1; 3910c51296cSAlp Dener ++nlsP->grad; 3920c51296cSAlp Dener stepType = NLS_GRADIENT; 39387f595a5SBarry Smith } else { 3940c51296cSAlp Dener ierr = MatLMVMGetUpdateCount(nlsP->M, &bfgsUpdates);CHKERRQ(ierr); 395a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 396a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 3970c51296cSAlp Dener ++nlsP->grad; 3980c51296cSAlp Dener stepType = NLS_GRADIENT; 39987f595a5SBarry Smith } else { 400a7e14dcfSSatish Balay ++nlsP->bfgs; 401a7e14dcfSSatish Balay stepType = NLS_BFGS; 402a7e14dcfSSatish Balay } 403a7e14dcfSSatish Balay } 404a7e14dcfSSatish Balay } 40587f595a5SBarry Smith } else { 406a7e14dcfSSatish Balay /* Computed Newton step is descent */ 407a7e14dcfSSatish Balay switch (ksp_reason) { 408a7e14dcfSSatish Balay case KSP_DIVERGED_NANORINF: 409a7e14dcfSSatish Balay case KSP_DIVERGED_BREAKDOWN: 410a7e14dcfSSatish Balay case KSP_DIVERGED_INDEFINITE_MAT: 411a7e14dcfSSatish Balay case KSP_DIVERGED_INDEFINITE_PC: 412a7e14dcfSSatish Balay case KSP_CONVERGED_CG_NEG_CURVE: 413a7e14dcfSSatish Balay /* Matrix or preconditioner is indefinite; increase perturbation */ 414a7e14dcfSSatish Balay if (pert <= 0.0) { 415a7e14dcfSSatish Balay /* Initialize the perturbation */ 416a7e14dcfSSatish Balay pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 4171daac38eSTodd Munson if (is_gltr) { 418ba7fe8fbSTodd Munson ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr); 419a7e14dcfSSatish Balay pert = PetscMax(pert, -e_min); 420a7e14dcfSSatish Balay } 42187f595a5SBarry Smith } else { 422a7e14dcfSSatish Balay /* Increase the perturbation */ 423a7e14dcfSSatish Balay pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 424a7e14dcfSSatish Balay } 425a7e14dcfSSatish Balay break; 426a7e14dcfSSatish Balay 427a7e14dcfSSatish Balay default: 428a7e14dcfSSatish Balay /* Newton step computation is good; decrease perturbation */ 429a7e14dcfSSatish Balay pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm); 430a7e14dcfSSatish Balay if (pert < nlsP->pmin) { 431a7e14dcfSSatish Balay pert = 0.0; 432a7e14dcfSSatish Balay } 433a7e14dcfSSatish Balay break; 434a7e14dcfSSatish Balay } 435a7e14dcfSSatish Balay 436a7e14dcfSSatish Balay ++nlsP->newt; 437a7e14dcfSSatish Balay stepType = NLS_NEWTON; 438a7e14dcfSSatish Balay } 439a7e14dcfSSatish Balay 440a7e14dcfSSatish Balay /* Perform the linesearch */ 441a7e14dcfSSatish Balay fold = f; 442a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, nlsP->Xold);CHKERRQ(ierr); 443a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, nlsP->Gold);CHKERRQ(ierr); 444a7e14dcfSSatish Balay 445a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr); 446a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 447a7e14dcfSSatish Balay 44887f595a5SBarry Smith while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) { /* Linesearch failed */ 449a7e14dcfSSatish Balay f = fold; 450a7e14dcfSSatish Balay ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr); 451a7e14dcfSSatish Balay ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr); 452a7e14dcfSSatish Balay 453a7e14dcfSSatish Balay switch(stepType) { 454a7e14dcfSSatish Balay case NLS_NEWTON: 455a7e14dcfSSatish Balay /* Failed to obtain acceptable iterate with Newton 1step 456a7e14dcfSSatish Balay Update the perturbation for next time */ 457a7e14dcfSSatish Balay if (pert <= 0.0) { 458a7e14dcfSSatish Balay /* Initialize the perturbation */ 459a7e14dcfSSatish Balay pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 4601daac38eSTodd Munson if (is_gltr) { 461ba7fe8fbSTodd Munson ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 462a7e14dcfSSatish Balay pert = PetscMax(pert, -e_min); 463a7e14dcfSSatish Balay } 46487f595a5SBarry Smith } else { 465a7e14dcfSSatish Balay /* Increase the perturbation */ 466a7e14dcfSSatish Balay pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 467a7e14dcfSSatish Balay } 468a7e14dcfSSatish Balay 4690c51296cSAlp Dener if (!nlsP->bfgs_pre) { 470a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and being updated 471a7e14dcfSSatish Balay Must use gradient direction in this case */ 472a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr); 473a7e14dcfSSatish Balay ++nlsP->grad; 474a7e14dcfSSatish Balay stepType = NLS_GRADIENT; 47587f595a5SBarry Smith } else { 476a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 477cd929ea3SAlp Dener ierr = MatSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 478a7e14dcfSSatish Balay /* Check for success (descent direction) */ 479a7e14dcfSSatish Balay ierr = VecDot(tao->solution, nlsP->D, &gdx);CHKERRQ(ierr); 480a7e14dcfSSatish Balay if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 481a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced not a number 482a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case 483a7e14dcfSSatish Balay Use steepest descent direction (scaled) */ 484cd929ea3SAlp Dener ierr = MatLMVMReset(nlsP->M, PETSC_FALSE);CHKERRQ(ierr); 485a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 486cd929ea3SAlp Dener ierr = MatSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 487a7e14dcfSSatish Balay 488a7e14dcfSSatish Balay bfgsUpdates = 1; 4890c51296cSAlp Dener ++nlsP->grad; 4900c51296cSAlp Dener stepType = NLS_GRADIENT; 49187f595a5SBarry Smith } else { 492a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 493a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 4940c51296cSAlp Dener ++nlsP->grad; 4950c51296cSAlp Dener stepType = NLS_GRADIENT; 49687f595a5SBarry Smith } else { 497a7e14dcfSSatish Balay ++nlsP->bfgs; 498a7e14dcfSSatish Balay stepType = NLS_BFGS; 499a7e14dcfSSatish Balay } 500a7e14dcfSSatish Balay } 501a7e14dcfSSatish Balay } 502a7e14dcfSSatish Balay break; 503a7e14dcfSSatish Balay 504a7e14dcfSSatish Balay case NLS_BFGS: 505a7e14dcfSSatish Balay /* Can only enter if pc_type == NLS_PC_BFGS 506a7e14dcfSSatish Balay Failed to obtain acceptable iterate with BFGS step 507a7e14dcfSSatish Balay Attempt to use the scaled gradient direction */ 508cd929ea3SAlp Dener ierr = MatLMVMReset(nlsP->M, PETSC_FALSE);CHKERRQ(ierr); 509a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 510cd929ea3SAlp Dener ierr = MatSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 511a7e14dcfSSatish Balay 512a7e14dcfSSatish Balay bfgsUpdates = 1; 513a7e14dcfSSatish Balay ++nlsP->grad; 514a7e14dcfSSatish Balay stepType = NLS_GRADIENT; 515a7e14dcfSSatish Balay break; 516a7e14dcfSSatish Balay } 517a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 518a7e14dcfSSatish Balay 519a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr); 520a7e14dcfSSatish Balay ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr); 521a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 522a7e14dcfSSatish Balay } 523a7e14dcfSSatish Balay 52487f595a5SBarry Smith if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 525a7e14dcfSSatish Balay /* Failed to find an improving point */ 526a7e14dcfSSatish Balay f = fold; 527a7e14dcfSSatish Balay ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr); 528a7e14dcfSSatish Balay ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr); 529a7e14dcfSSatish Balay step = 0.0; 530a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 531a7e14dcfSSatish Balay break; 532a7e14dcfSSatish Balay } 533a7e14dcfSSatish Balay 534a7e14dcfSSatish Balay /* Update trust region radius */ 5351daac38eSTodd Munson if (is_nash || is_stcg || is_gltr) { 536a7e14dcfSSatish Balay switch(nlsP->update_type) { 537a7e14dcfSSatish Balay case NLS_UPDATE_STEP: 538a7e14dcfSSatish Balay if (stepType == NLS_NEWTON) { 539a7e14dcfSSatish Balay if (step < nlsP->nu1) { 540a7e14dcfSSatish Balay /* Very bad step taken; reduce radius */ 541a7e14dcfSSatish Balay tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust); 54287f595a5SBarry Smith } else if (step < nlsP->nu2) { 543a7e14dcfSSatish Balay /* Reasonably bad step taken; reduce radius */ 544a7e14dcfSSatish Balay tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust); 54587f595a5SBarry Smith } else if (step < nlsP->nu3) { 546a7e14dcfSSatish Balay /* Reasonable step was taken; leave radius alone */ 547a7e14dcfSSatish Balay if (nlsP->omega3 < 1.0) { 548a7e14dcfSSatish Balay tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust); 54987f595a5SBarry Smith } else if (nlsP->omega3 > 1.0) { 550a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust); 551a7e14dcfSSatish Balay } 55287f595a5SBarry Smith } else if (step < nlsP->nu4) { 553a7e14dcfSSatish Balay /* Full step taken; increase the radius */ 554a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust); 55587f595a5SBarry Smith } else { 556a7e14dcfSSatish Balay /* More than full step taken; increase the radius */ 557a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust); 558a7e14dcfSSatish Balay } 55987f595a5SBarry Smith } else { 560a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 561a7e14dcfSSatish Balay tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust); 562a7e14dcfSSatish Balay } 563a7e14dcfSSatish Balay break; 564a7e14dcfSSatish Balay 565a7e14dcfSSatish Balay case NLS_UPDATE_REDUCTION: 566a7e14dcfSSatish Balay if (stepType == NLS_NEWTON) { 567a7e14dcfSSatish Balay /* Get predicted reduction */ 568a7e14dcfSSatish Balay 569ba7fe8fbSTodd Munson ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 570a7e14dcfSSatish Balay if (prered >= 0.0) { 571a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot */ 572a7e14dcfSSatish Balay /* happen in infinite precision arithmetic. Step should */ 573a7e14dcfSSatish Balay /* be rejected! */ 574a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 57587f595a5SBarry Smith } else { 576a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(f_full)) { 577a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 57887f595a5SBarry Smith } else { 579a7e14dcfSSatish Balay /* Compute and actual reduction */ 580a7e14dcfSSatish Balay actred = fold - f_full; 581a7e14dcfSSatish Balay prered = -prered; 582a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= nlsP->epsilon) && 583a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= nlsP->epsilon)) { 584a7e14dcfSSatish Balay kappa = 1.0; 58587f595a5SBarry Smith } else { 586a7e14dcfSSatish Balay kappa = actred / prered; 587a7e14dcfSSatish Balay } 588a7e14dcfSSatish Balay 589a7e14dcfSSatish Balay /* Accept of reject the step and update radius */ 590a7e14dcfSSatish Balay if (kappa < nlsP->eta1) { 591a7e14dcfSSatish Balay /* Very bad step */ 592a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 59387f595a5SBarry Smith } else if (kappa < nlsP->eta2) { 594a7e14dcfSSatish Balay /* Marginal bad step */ 595a7e14dcfSSatish Balay tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d); 59687f595a5SBarry Smith } else if (kappa < nlsP->eta3) { 597a7e14dcfSSatish Balay /* Reasonable step */ 598a7e14dcfSSatish Balay if (nlsP->alpha3 < 1.0) { 599a7e14dcfSSatish Balay tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust); 60087f595a5SBarry Smith } else if (nlsP->alpha3 > 1.0) { 601a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust); 602a7e14dcfSSatish Balay } 60387f595a5SBarry Smith } else if (kappa < nlsP->eta4) { 604a7e14dcfSSatish Balay /* Good step */ 605a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust); 60687f595a5SBarry Smith } else { 607a7e14dcfSSatish Balay /* Very good step */ 608a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust); 609a7e14dcfSSatish Balay } 610a7e14dcfSSatish Balay } 611a7e14dcfSSatish Balay } 61287f595a5SBarry Smith } else { 613a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 614a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust); 615a7e14dcfSSatish Balay } 616a7e14dcfSSatish Balay break; 617a7e14dcfSSatish Balay 618a7e14dcfSSatish Balay default: 619a7e14dcfSSatish Balay if (stepType == NLS_NEWTON) { 620ba7fe8fbSTodd Munson ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 621a7e14dcfSSatish Balay if (prered >= 0.0) { 622a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot */ 623a7e14dcfSSatish Balay /* happen in infinite precision arithmetic. Step should */ 624a7e14dcfSSatish Balay /* be rejected! */ 625a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 62687f595a5SBarry Smith } else { 627a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(f_full)) { 628a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 62987f595a5SBarry Smith } else { 630a7e14dcfSSatish Balay actred = fold - f_full; 631a7e14dcfSSatish Balay prered = -prered; 63287f595a5SBarry Smith if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) { 633a7e14dcfSSatish Balay kappa = 1.0; 63487f595a5SBarry Smith } else { 635a7e14dcfSSatish Balay kappa = actred / prered; 636a7e14dcfSSatish Balay } 637a7e14dcfSSatish Balay 638a7e14dcfSSatish Balay tau_1 = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred); 639a7e14dcfSSatish Balay tau_2 = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred); 640a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 641a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 642a7e14dcfSSatish Balay 643a7e14dcfSSatish Balay if (kappa >= 1.0 - nlsP->mu1) { 644a7e14dcfSSatish Balay /* Great agreement */ 645a7e14dcfSSatish Balay if (tau_max < 1.0) { 646a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d); 64787f595a5SBarry Smith } else if (tau_max > nlsP->gamma4) { 648a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d); 64987f595a5SBarry Smith } else { 650a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 651a7e14dcfSSatish Balay } 65287f595a5SBarry Smith } else if (kappa >= 1.0 - nlsP->mu2) { 653a7e14dcfSSatish Balay /* Good agreement */ 654a7e14dcfSSatish Balay 655a7e14dcfSSatish Balay if (tau_max < nlsP->gamma2) { 656a7e14dcfSSatish Balay tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d); 65787f595a5SBarry Smith } else if (tau_max > nlsP->gamma3) { 658a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d); 65987f595a5SBarry Smith } else if (tau_max < 1.0) { 660a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 66187f595a5SBarry Smith } else { 662a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 663a7e14dcfSSatish Balay } 66487f595a5SBarry Smith } else { 665a7e14dcfSSatish Balay /* Not good agreement */ 666a7e14dcfSSatish Balay if (tau_min > 1.0) { 667a7e14dcfSSatish Balay tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d); 66887f595a5SBarry Smith } else if (tau_max < nlsP->gamma1) { 669a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 67087f595a5SBarry Smith } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) { 671a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 67287f595a5SBarry Smith } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) { 673a7e14dcfSSatish Balay tao->trust = tau_1 * PetscMin(tao->trust, norm_d); 67487f595a5SBarry Smith } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) { 675a7e14dcfSSatish Balay tao->trust = tau_2 * PetscMin(tao->trust, norm_d); 67687f595a5SBarry Smith } else { 677a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 678a7e14dcfSSatish Balay } 679a7e14dcfSSatish Balay } 680a7e14dcfSSatish Balay } 681a7e14dcfSSatish Balay } 68287f595a5SBarry Smith } else { 683a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 684a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust); 685a7e14dcfSSatish Balay } 686a7e14dcfSSatish Balay break; 687a7e14dcfSSatish Balay } 688a7e14dcfSSatish Balay 689a7e14dcfSSatish Balay /* The radius may have been increased; modify if it is too large */ 690a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 691a7e14dcfSSatish Balay } 692a7e14dcfSSatish Balay 693a7e14dcfSSatish Balay /* Check for termination */ 694a9603a14SPatrick Farrell ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 69587f595a5SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 696a7e14dcfSSatish Balay needH = 1; 6973ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 6983ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 6993ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 700a7e14dcfSSatish Balay } 701a7e14dcfSSatish Balay PetscFunctionReturn(0); 702a7e14dcfSSatish Balay } 703a7e14dcfSSatish Balay 704a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 705441846f8SBarry Smith static PetscErrorCode TaoSetUp_NLS(Tao tao) 706a7e14dcfSSatish Balay { 707a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 708a7e14dcfSSatish Balay PetscErrorCode ierr; 709a7e14dcfSSatish Balay 710a7e14dcfSSatish Balay PetscFunctionBegin; 711a7e14dcfSSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);} 712a7e14dcfSSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);} 713a7e14dcfSSatish Balay if (!nlsP->W) {ierr = VecDuplicate(tao->solution,&nlsP->W);CHKERRQ(ierr);} 714a7e14dcfSSatish Balay if (!nlsP->D) {ierr = VecDuplicate(tao->solution,&nlsP->D);CHKERRQ(ierr);} 715a7e14dcfSSatish Balay if (!nlsP->Xold) {ierr = VecDuplicate(tao->solution,&nlsP->Xold);CHKERRQ(ierr);} 716a7e14dcfSSatish Balay if (!nlsP->Gold) {ierr = VecDuplicate(tao->solution,&nlsP->Gold);CHKERRQ(ierr);} 717a7e14dcfSSatish Balay nlsP->M = 0; 7180c51296cSAlp Dener nlsP->bfgs_pre = 0; 719a7e14dcfSSatish Balay PetscFunctionReturn(0); 720a7e14dcfSSatish Balay } 721a7e14dcfSSatish Balay 722a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 723441846f8SBarry Smith static PetscErrorCode TaoDestroy_NLS(Tao tao) 724a7e14dcfSSatish Balay { 725a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 726a7e14dcfSSatish Balay PetscErrorCode ierr; 727a7e14dcfSSatish Balay 728a7e14dcfSSatish Balay PetscFunctionBegin; 729a7e14dcfSSatish Balay if (tao->setupcalled) { 730a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->D);CHKERRQ(ierr); 731a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->W);CHKERRQ(ierr); 732a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->Xold);CHKERRQ(ierr); 733a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->Gold);CHKERRQ(ierr); 734a7e14dcfSSatish Balay } 735a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 736a7e14dcfSSatish Balay PetscFunctionReturn(0); 737a7e14dcfSSatish Balay } 738a7e14dcfSSatish Balay 739a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 7404416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_NLS(PetscOptionItems *PetscOptionsObject,Tao tao) 741a7e14dcfSSatish Balay { 742a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 743a7e14dcfSSatish Balay PetscErrorCode ierr; 744a7e14dcfSSatish Balay 745a7e14dcfSSatish Balay PetscFunctionBegin; 7461a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr); 747a7e14dcfSSatish Balay ierr = PetscOptionsEList("-tao_nls_init_type", "radius initialization type", "", NLS_INIT, NLS_INIT_TYPES, NLS_INIT[nlsP->init_type], &nlsP->init_type, 0);CHKERRQ(ierr); 748a7e14dcfSSatish Balay ierr = PetscOptionsEList("-tao_nls_update_type", "radius update type", "", NLS_UPDATE, NLS_UPDATE_TYPES, NLS_UPDATE[nlsP->update_type], &nlsP->update_type, 0);CHKERRQ(ierr); 74994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval,NULL);CHKERRQ(ierr); 75094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin,NULL);CHKERRQ(ierr); 75194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax,NULL);CHKERRQ(ierr); 75294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac,NULL);CHKERRQ(ierr); 75394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin,NULL);CHKERRQ(ierr); 75494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax,NULL);CHKERRQ(ierr); 75594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac,NULL);CHKERRQ(ierr); 75694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac,NULL);CHKERRQ(ierr); 75794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac,NULL);CHKERRQ(ierr); 75894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac,NULL);CHKERRQ(ierr); 75994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1,NULL);CHKERRQ(ierr); 76094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2,NULL);CHKERRQ(ierr); 76194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3,NULL);CHKERRQ(ierr); 76294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4,NULL);CHKERRQ(ierr); 76394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1,NULL);CHKERRQ(ierr); 76494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2,NULL);CHKERRQ(ierr); 76594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3,NULL);CHKERRQ(ierr); 76694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4,NULL);CHKERRQ(ierr); 76794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5,NULL);CHKERRQ(ierr); 76894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1,NULL);CHKERRQ(ierr); 76994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2,NULL);CHKERRQ(ierr); 77094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3,NULL);CHKERRQ(ierr); 77194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4,NULL);CHKERRQ(ierr); 77294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1,NULL);CHKERRQ(ierr); 77394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2,NULL);CHKERRQ(ierr); 77494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3,NULL);CHKERRQ(ierr); 77594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4,NULL);CHKERRQ(ierr); 77694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5,NULL);CHKERRQ(ierr); 77794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i,NULL);CHKERRQ(ierr); 77894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i,NULL);CHKERRQ(ierr); 77994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i,NULL);CHKERRQ(ierr); 78094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i,NULL);CHKERRQ(ierr); 78194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i,NULL);CHKERRQ(ierr); 78294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i,NULL);CHKERRQ(ierr); 78394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i,NULL);CHKERRQ(ierr); 78494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1,NULL);CHKERRQ(ierr); 78594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2,NULL);CHKERRQ(ierr); 78694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1,NULL);CHKERRQ(ierr); 78794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2,NULL);CHKERRQ(ierr); 78894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3,NULL);CHKERRQ(ierr); 78994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4,NULL);CHKERRQ(ierr); 79094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta,NULL);CHKERRQ(ierr); 79194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius,NULL);CHKERRQ(ierr); 79294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius,NULL);CHKERRQ(ierr); 79394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon,NULL);CHKERRQ(ierr); 794a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 795a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 796a7e14dcfSSatish Balay ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 797a7e14dcfSSatish Balay PetscFunctionReturn(0); 798a7e14dcfSSatish Balay } 799a7e14dcfSSatish Balay 800a7e14dcfSSatish Balay 801a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 802441846f8SBarry Smith static PetscErrorCode TaoView_NLS(Tao tao, PetscViewer viewer) 803a7e14dcfSSatish Balay { 804a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 805a7e14dcfSSatish Balay PetscBool isascii; 806a7e14dcfSSatish Balay PetscErrorCode ierr; 807a7e14dcfSSatish Balay 808a7e14dcfSSatish Balay PetscFunctionBegin; 809a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 810a7e14dcfSSatish Balay if (isascii) { 811a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 812a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", nlsP->newt);CHKERRQ(ierr); 813a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", nlsP->bfgs);CHKERRQ(ierr); 814a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", nlsP->grad);CHKERRQ(ierr); 815a7e14dcfSSatish Balay 816a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp atol: %D\n", nlsP->ksp_atol);CHKERRQ(ierr); 817a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %D\n", nlsP->ksp_rtol);CHKERRQ(ierr); 818a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %D\n", nlsP->ksp_ctol);CHKERRQ(ierr); 819a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp negc: %D\n", nlsP->ksp_negc);CHKERRQ(ierr); 820a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %D\n", nlsP->ksp_dtol);CHKERRQ(ierr); 821a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp iter: %D\n", nlsP->ksp_iter);CHKERRQ(ierr); 822a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp othr: %D\n", nlsP->ksp_othr);CHKERRQ(ierr); 823a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 824a7e14dcfSSatish Balay } 825a7e14dcfSSatish Balay PetscFunctionReturn(0); 826a7e14dcfSSatish Balay } 827a7e14dcfSSatish Balay 828a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 8294aa34175SJason Sarich /*MC 8304aa34175SJason Sarich TAONLS - Newton's method with linesearch for unconstrained minimization. 8314aa34175SJason Sarich At each iteration, the Newton line search method solves the symmetric 8324aa34175SJason Sarich system of equations to obtain the step diretion dk: 8334aa34175SJason Sarich Hk dk = -gk 8344aa34175SJason Sarich a More-Thuente line search is applied on the direction dk to approximately 8354aa34175SJason Sarich solve 8364aa34175SJason Sarich min_t f(xk + t d_k) 8374aa34175SJason Sarich 8384aa34175SJason Sarich Options Database Keys: 8391daac38eSTodd Munson + -tao_nls_pc_type - "none","ahess","bfgs","petsc" 8404aa34175SJason Sarich . -tao_nls_bfgs_scale_type - "ahess","phess","bfgs" 8414aa34175SJason Sarich . -tao_nls_init_type - "constant","direction","interpolation" 8424aa34175SJason Sarich . -tao_nls_update_type - "step","direction","interpolation" 8434aa34175SJason Sarich . -tao_nls_sval - perturbation starting value 8444aa34175SJason Sarich . -tao_nls_imin - minimum initial perturbation 8454aa34175SJason Sarich . -tao_nls_imax - maximum initial perturbation 8464aa34175SJason Sarich . -tao_nls_pmin - minimum perturbation 8474aa34175SJason Sarich . -tao_nls_pmax - maximum perturbation 8484aa34175SJason Sarich . -tao_nls_pgfac - growth factor 8494aa34175SJason Sarich . -tao_nls_psfac - shrink factor 8504aa34175SJason Sarich . -tao_nls_imfac - initial merit factor 8514aa34175SJason Sarich . -tao_nls_pmgfac - merit growth factor 8524aa34175SJason Sarich . -tao_nls_pmsfac - merit shrink factor 8534aa34175SJason Sarich . -tao_nls_eta1 - poor steplength; reduce radius 8544aa34175SJason Sarich . -tao_nls_eta2 - reasonable steplength; leave radius 8554aa34175SJason Sarich . -tao_nls_eta3 - good steplength; increase readius 8564aa34175SJason Sarich . -tao_nls_eta4 - excellent steplength; greatly increase radius 8574aa34175SJason Sarich . -tao_nls_alpha1 - alpha1 reduction 8584aa34175SJason Sarich . -tao_nls_alpha2 - alpha2 reduction 8594aa34175SJason Sarich . -tao_nls_alpha3 - alpha3 reduction 8604aa34175SJason Sarich . -tao_nls_alpha4 - alpha4 reduction 8614aa34175SJason Sarich . -tao_nls_alpha - alpha5 reduction 8624aa34175SJason Sarich . -tao_nls_mu1 - mu1 interpolation update 8634aa34175SJason Sarich . -tao_nls_mu2 - mu2 interpolation update 8644aa34175SJason Sarich . -tao_nls_gamma1 - gamma1 interpolation update 8654aa34175SJason Sarich . -tao_nls_gamma2 - gamma2 interpolation update 8664aa34175SJason Sarich . -tao_nls_gamma3 - gamma3 interpolation update 8674aa34175SJason Sarich . -tao_nls_gamma4 - gamma4 interpolation update 8684aa34175SJason Sarich . -tao_nls_theta - theta interpolation update 8694aa34175SJason Sarich . -tao_nls_omega1 - omega1 step update 8704aa34175SJason Sarich . -tao_nls_omega2 - omega2 step update 8714aa34175SJason Sarich . -tao_nls_omega3 - omega3 step update 8724aa34175SJason Sarich . -tao_nls_omega4 - omega4 step update 8734aa34175SJason Sarich . -tao_nls_omega5 - omega5 step update 8741522df2eSJason Sarich . -tao_nls_mu1_i - mu1 interpolation init factor 8751522df2eSJason Sarich . -tao_nls_mu2_i - mu2 interpolation init factor 8761522df2eSJason Sarich . -tao_nls_gamma1_i - gamma1 interpolation init factor 8771522df2eSJason Sarich . -tao_nls_gamma2_i - gamma2 interpolation init factor 8781522df2eSJason Sarich . -tao_nls_gamma3_i - gamma3 interpolation init factor 8791522df2eSJason Sarich . -tao_nls_gamma4_i - gamma4 interpolation init factor 8801522df2eSJason Sarich - -tao_nls_theta_i - theta interpolation init factor 8811eb8069cSJason Sarich 8821eb8069cSJason Sarich Level: beginner 8831522df2eSJason Sarich M*/ 8844aa34175SJason Sarich 885728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NLS(Tao tao) 886a7e14dcfSSatish Balay { 887a7e14dcfSSatish Balay TAO_NLS *nlsP; 8888caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 889a7e14dcfSSatish Balay PetscErrorCode ierr; 890a7e14dcfSSatish Balay 891a7e14dcfSSatish Balay PetscFunctionBegin; 8923c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&nlsP);CHKERRQ(ierr); 893a7e14dcfSSatish Balay 894a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_NLS; 895a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_NLS; 896a7e14dcfSSatish Balay tao->ops->view = TaoView_NLS; 897a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_NLS; 898a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_NLS; 899a7e14dcfSSatish Balay 9006552cf8aSJason Sarich /* Override default settings (unless already changed) */ 9016552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 50; 9026552cf8aSJason Sarich if (!tao->trust0_changed) tao->trust0 = 100.0; 9036552cf8aSJason Sarich 904a7e14dcfSSatish Balay tao->data = (void*)nlsP; 905a7e14dcfSSatish Balay 906a7e14dcfSSatish Balay nlsP->sval = 0.0; 907a7e14dcfSSatish Balay nlsP->imin = 1.0e-4; 908a7e14dcfSSatish Balay nlsP->imax = 1.0e+2; 909a7e14dcfSSatish Balay nlsP->imfac = 1.0e-1; 910a7e14dcfSSatish Balay 911a7e14dcfSSatish Balay nlsP->pmin = 1.0e-12; 912a7e14dcfSSatish Balay nlsP->pmax = 1.0e+2; 913a7e14dcfSSatish Balay nlsP->pgfac = 1.0e+1; 914a7e14dcfSSatish Balay nlsP->psfac = 4.0e-1; 915a7e14dcfSSatish Balay nlsP->pmgfac = 1.0e-1; 916a7e14dcfSSatish Balay nlsP->pmsfac = 1.0e-1; 917a7e14dcfSSatish Balay 918a7e14dcfSSatish Balay /* Default values for trust-region radius update based on steplength */ 919a7e14dcfSSatish Balay nlsP->nu1 = 0.25; 920a7e14dcfSSatish Balay nlsP->nu2 = 0.50; 921a7e14dcfSSatish Balay nlsP->nu3 = 1.00; 922a7e14dcfSSatish Balay nlsP->nu4 = 1.25; 923a7e14dcfSSatish Balay 924a7e14dcfSSatish Balay nlsP->omega1 = 0.25; 925a7e14dcfSSatish Balay nlsP->omega2 = 0.50; 926a7e14dcfSSatish Balay nlsP->omega3 = 1.00; 927a7e14dcfSSatish Balay nlsP->omega4 = 2.00; 928a7e14dcfSSatish Balay nlsP->omega5 = 4.00; 929a7e14dcfSSatish Balay 930a7e14dcfSSatish Balay /* Default values for trust-region radius update based on reduction */ 931a7e14dcfSSatish Balay nlsP->eta1 = 1.0e-4; 932a7e14dcfSSatish Balay nlsP->eta2 = 0.25; 933a7e14dcfSSatish Balay nlsP->eta3 = 0.50; 934a7e14dcfSSatish Balay nlsP->eta4 = 0.90; 935a7e14dcfSSatish Balay 936a7e14dcfSSatish Balay nlsP->alpha1 = 0.25; 937a7e14dcfSSatish Balay nlsP->alpha2 = 0.50; 938a7e14dcfSSatish Balay nlsP->alpha3 = 1.00; 939a7e14dcfSSatish Balay nlsP->alpha4 = 2.00; 940a7e14dcfSSatish Balay nlsP->alpha5 = 4.00; 941a7e14dcfSSatish Balay 942a7e14dcfSSatish Balay /* Default values for trust-region radius update based on interpolation */ 943a7e14dcfSSatish Balay nlsP->mu1 = 0.10; 944a7e14dcfSSatish Balay nlsP->mu2 = 0.50; 945a7e14dcfSSatish Balay 946a7e14dcfSSatish Balay nlsP->gamma1 = 0.25; 947a7e14dcfSSatish Balay nlsP->gamma2 = 0.50; 948a7e14dcfSSatish Balay nlsP->gamma3 = 2.00; 949a7e14dcfSSatish Balay nlsP->gamma4 = 4.00; 950a7e14dcfSSatish Balay 951a7e14dcfSSatish Balay nlsP->theta = 0.05; 952a7e14dcfSSatish Balay 953a7e14dcfSSatish Balay /* Default values for trust region initialization based on interpolation */ 954a7e14dcfSSatish Balay nlsP->mu1_i = 0.35; 955a7e14dcfSSatish Balay nlsP->mu2_i = 0.50; 956a7e14dcfSSatish Balay 957a7e14dcfSSatish Balay nlsP->gamma1_i = 0.0625; 958a7e14dcfSSatish Balay nlsP->gamma2_i = 0.5; 959a7e14dcfSSatish Balay nlsP->gamma3_i = 2.0; 960a7e14dcfSSatish Balay nlsP->gamma4_i = 5.0; 961a7e14dcfSSatish Balay 962a7e14dcfSSatish Balay nlsP->theta_i = 0.25; 963a7e14dcfSSatish Balay 964a7e14dcfSSatish Balay /* Remaining parameters */ 965a7e14dcfSSatish Balay nlsP->min_radius = 1.0e-10; 966a7e14dcfSSatish Balay nlsP->max_radius = 1.0e10; 967a7e14dcfSSatish Balay nlsP->epsilon = 1.0e-6; 968a7e14dcfSSatish Balay 969a7e14dcfSSatish Balay nlsP->init_type = NLS_INIT_INTERPOLATION; 970a7e14dcfSSatish Balay nlsP->update_type = NLS_UPDATE_STEP; 971a7e14dcfSSatish Balay 972a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 97363b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 974a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 975441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 9765d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 977a7e14dcfSSatish Balay 978a7e14dcfSSatish Balay /* Set linear solver to default for symmetric matrices */ 979a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 98063b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 981*cbf034f8SAlp Dener ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr); 982*cbf034f8SAlp Dener ierr = KSPAppendOptionsPrefix(tao->ksp, "tao_nls_");CHKERRQ(ierr); 9831daac38eSTodd Munson ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr); 984a7e14dcfSSatish Balay PetscFunctionReturn(0); 985a7e14dcfSSatish Balay } 986