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 TAO_NLS *nlsP = (TAO_NLS *)tao->data; 411daac38eSTodd Munson KSPType ksp_type; 420ad3a497SAlp Dener PetscBool is_nash,is_stcg,is_gltr,is_bfgs,is_jacobi,is_symmetric,sym_set; 43a7e14dcfSSatish Balay KSPConvergedReason ksp_reason; 441daac38eSTodd Munson PC pc; 451daac38eSTodd Munson TaoLineSearchConvergedReason ls_reason; 46a7e14dcfSSatish Balay 47a7e14dcfSSatish Balay PetscReal fmin, ftrial, f_full, prered, actred, kappa, sigma; 48a7e14dcfSSatish Balay PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 49a7e14dcfSSatish Balay PetscReal f, fold, gdx, gnorm, pert; 50a7e14dcfSSatish Balay PetscReal step = 1.0; 51a7e14dcfSSatish Balay PetscReal norm_d = 0.0, e_min; 52a7e14dcfSSatish Balay 53a7e14dcfSSatish Balay PetscInt stepType; 54a7e14dcfSSatish Balay PetscInt bfgsUpdates = 0; 55a7e14dcfSSatish Balay PetscInt n,N,kspits; 56b4690188SBarry Smith PetscInt needH = 1; 57a7e14dcfSSatish Balay 58a7e14dcfSSatish Balay PetscInt i_max = 5; 59a7e14dcfSSatish Balay PetscInt j_max = 1; 60a7e14dcfSSatish Balay PetscInt i, j; 61a7e14dcfSSatish Balay 62a7e14dcfSSatish Balay PetscFunctionBegin; 63a7e14dcfSSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 649566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao,"WARNING: Variable bounds have been set but will be ignored by nls algorithm\n")); 65a7e14dcfSSatish Balay } 66a7e14dcfSSatish Balay 67a7e14dcfSSatish Balay /* Initialized variables */ 68a7e14dcfSSatish Balay pert = nlsP->sval; 69a7e14dcfSSatish Balay 702d9aa51bSJason Sarich /* Number of times ksp stopped because of these reasons */ 71a7e14dcfSSatish Balay nlsP->ksp_atol = 0; 72a7e14dcfSSatish Balay nlsP->ksp_rtol = 0; 73a7e14dcfSSatish Balay nlsP->ksp_dtol = 0; 74a7e14dcfSSatish Balay nlsP->ksp_ctol = 0; 75a7e14dcfSSatish Balay nlsP->ksp_negc = 0; 76a7e14dcfSSatish Balay nlsP->ksp_iter = 0; 77a7e14dcfSSatish Balay nlsP->ksp_othr = 0; 78a7e14dcfSSatish Balay 79a7e14dcfSSatish Balay /* Initialize trust-region radius when using nash, stcg, or gltr 80ba7fe8fbSTodd Munson Command automatically ignored for other methods 81ba7fe8fbSTodd Munson Will be reset during the first iteration 82ba7fe8fbSTodd Munson */ 839566063dSJacob Faibussowitsch PetscCall(KSPGetType(tao->ksp,&ksp_type)); 849566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ksp_type,KSPNASH,&is_nash)); 859566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ksp_type,KSPSTCG,&is_stcg)); 869566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ksp_type,KSPGLTR,&is_gltr)); 871daac38eSTodd Munson 889566063dSJacob Faibussowitsch PetscCall(KSPCGSetRadius(tao->ksp,nlsP->max_radius)); 89a7e14dcfSSatish Balay 901daac38eSTodd Munson if (is_nash || is_stcg || is_gltr) { 913c859ba3SBarry Smith PetscCheck(tao->trust0 >= 0.0,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_OUTOFRANGE,"Initial radius negative"); 92a7e14dcfSSatish Balay tao->trust = tao->trust0; 93a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 94a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 95a7e14dcfSSatish Balay } 96a7e14dcfSSatish Balay 97a7e14dcfSSatish Balay /* Check convergence criteria */ 989566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient)); 999566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm)); 1003c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 101a7e14dcfSSatish Balay 1023ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 1039566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its)); 1049566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,step)); 1059566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 1063ecd9318SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 107a7e14dcfSSatish Balay 1080c51296cSAlp Dener /* Allocate the vectors needed for the BFGS approximation */ 1099566063dSJacob Faibussowitsch PetscCall(KSPGetPC(tao->ksp, &pc)); 1109566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs)); 1119566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi)); 1120c51296cSAlp Dener if (is_bfgs) { 1130c51296cSAlp Dener nlsP->bfgs_pre = pc; 1149566063dSJacob Faibussowitsch PetscCall(PCLMVMGetMatLMVM(nlsP->bfgs_pre, &nlsP->M)); 1159566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution, &n)); 1169566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &N)); 1179566063dSJacob Faibussowitsch PetscCall(MatSetSizes(nlsP->M, n, n, N, N)); 1189566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(nlsP->M, tao->solution, tao->gradient)); 1199566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(nlsP->M, &sym_set, &is_symmetric)); 1203c859ba3SBarry Smith PetscCheck(sym_set && is_symmetric,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric."); 1210c51296cSAlp Dener } else if (is_jacobi) { 1229566063dSJacob Faibussowitsch PetscCall(PCJacobiSetUseAbs(pc,PETSC_TRUE)); 123a7e14dcfSSatish Balay } 124a7e14dcfSSatish Balay 125a7e14dcfSSatish Balay /* Initialize trust-region radius. The initialization is only performed 126a7e14dcfSSatish Balay when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */ 1271daac38eSTodd Munson if (is_nash || is_stcg || is_gltr) { 128a7e14dcfSSatish Balay switch (nlsP->init_type) { 129a7e14dcfSSatish Balay case NLS_INIT_CONSTANT: 130a7e14dcfSSatish Balay /* Use the initial radius specified */ 131a7e14dcfSSatish Balay break; 132a7e14dcfSSatish Balay 133a7e14dcfSSatish Balay case NLS_INIT_INTERPOLATION: 134a7e14dcfSSatish Balay /* Use the initial radius specified */ 135a7e14dcfSSatish Balay max_radius = 0.0; 136a7e14dcfSSatish Balay 137a7e14dcfSSatish Balay for (j = 0; j < j_max; ++j) { 138a7e14dcfSSatish Balay fmin = f; 139a7e14dcfSSatish Balay sigma = 0.0; 140a7e14dcfSSatish Balay 141a7e14dcfSSatish Balay if (needH) { 1429566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre)); 143a7e14dcfSSatish Balay needH = 0; 144a7e14dcfSSatish Balay } 145a7e14dcfSSatish Balay 146a7e14dcfSSatish Balay for (i = 0; i < i_max; ++i) { 1479566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution,nlsP->W)); 1489566063dSJacob Faibussowitsch PetscCall(VecAXPY(nlsP->W,-tao->trust/gnorm,tao->gradient)); 1499566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao, nlsP->W, &ftrial)); 150a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 151a7e14dcfSSatish Balay tau = nlsP->gamma1_i; 15287f595a5SBarry Smith } else { 153a7e14dcfSSatish Balay if (ftrial < fmin) { 154a7e14dcfSSatish Balay fmin = ftrial; 155a7e14dcfSSatish Balay sigma = -tao->trust / gnorm; 156a7e14dcfSSatish Balay } 157a7e14dcfSSatish Balay 1589566063dSJacob Faibussowitsch PetscCall(MatMult(tao->hessian, tao->gradient, nlsP->D)); 1599566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, nlsP->D, &prered)); 160a7e14dcfSSatish Balay 161a7e14dcfSSatish Balay prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm)); 162a7e14dcfSSatish Balay actred = f - ftrial; 16387f595a5SBarry Smith if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) { 164a7e14dcfSSatish Balay kappa = 1.0; 16587f595a5SBarry Smith } else { 166a7e14dcfSSatish Balay kappa = actred / prered; 167a7e14dcfSSatish Balay } 168a7e14dcfSSatish Balay 169a7e14dcfSSatish Balay tau_1 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred); 170a7e14dcfSSatish Balay tau_2 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred); 171a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 172a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 173a7e14dcfSSatish Balay 17418cfbf8eSSatish Balay if (PetscAbsScalar(kappa - (PetscReal)1.0) <= nlsP->mu1_i) { 175a7e14dcfSSatish Balay /* Great agreement */ 176a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 177a7e14dcfSSatish Balay 178a7e14dcfSSatish Balay if (tau_max < 1.0) { 179a7e14dcfSSatish Balay tau = nlsP->gamma3_i; 18087f595a5SBarry Smith } else if (tau_max > nlsP->gamma4_i) { 181a7e14dcfSSatish Balay tau = nlsP->gamma4_i; 18287f595a5SBarry Smith } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) { 183a7e14dcfSSatish Balay tau = tau_1; 18487f595a5SBarry Smith } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) { 185a7e14dcfSSatish Balay tau = tau_2; 18687f595a5SBarry Smith } else { 187a7e14dcfSSatish Balay tau = tau_max; 188a7e14dcfSSatish Balay } 18918cfbf8eSSatish Balay } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= nlsP->mu2_i) { 190a7e14dcfSSatish Balay /* Good agreement */ 191a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 192a7e14dcfSSatish Balay 193a7e14dcfSSatish Balay if (tau_max < nlsP->gamma2_i) { 194a7e14dcfSSatish Balay tau = nlsP->gamma2_i; 19587f595a5SBarry Smith } else if (tau_max > nlsP->gamma3_i) { 196a7e14dcfSSatish Balay tau = nlsP->gamma3_i; 19787f595a5SBarry Smith } else { 198a7e14dcfSSatish Balay tau = tau_max; 199a7e14dcfSSatish Balay } 20087f595a5SBarry Smith } else { 201a7e14dcfSSatish Balay /* Not good agreement */ 202a7e14dcfSSatish Balay if (tau_min > 1.0) { 203a7e14dcfSSatish Balay tau = nlsP->gamma2_i; 20487f595a5SBarry Smith } else if (tau_max < nlsP->gamma1_i) { 205a7e14dcfSSatish Balay tau = nlsP->gamma1_i; 20687f595a5SBarry Smith } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) { 207a7e14dcfSSatish Balay tau = nlsP->gamma1_i; 20887f595a5SBarry Smith } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) { 209a7e14dcfSSatish Balay tau = tau_1; 21087f595a5SBarry Smith } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) { 211a7e14dcfSSatish Balay tau = tau_2; 21287f595a5SBarry Smith } else { 213a7e14dcfSSatish Balay tau = tau_max; 214a7e14dcfSSatish Balay } 215a7e14dcfSSatish Balay } 216a7e14dcfSSatish Balay } 217a7e14dcfSSatish Balay tao->trust = tau * tao->trust; 218a7e14dcfSSatish Balay } 219a7e14dcfSSatish Balay 220a7e14dcfSSatish Balay if (fmin < f) { 221a7e14dcfSSatish Balay f = fmin; 2229566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution,sigma,tao->gradient)); 2239566063dSJacob Faibussowitsch PetscCall(TaoComputeGradient(tao,tao->solution,tao->gradient)); 224a7e14dcfSSatish Balay 2259566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm)); 2263c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute gradient generated Inf or NaN"); 227a7e14dcfSSatish Balay needH = 1; 228a7e14dcfSSatish Balay 2299566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its)); 2309566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,step)); 2319566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 2323ecd9318SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 233a7e14dcfSSatish Balay } 234a7e14dcfSSatish Balay } 235a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, max_radius); 236a7e14dcfSSatish Balay 237a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 238a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 239a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 240a7e14dcfSSatish Balay break; 241a7e14dcfSSatish Balay 242a7e14dcfSSatish Balay default: 243a7e14dcfSSatish Balay /* Norm of the first direction will initialize radius */ 244a7e14dcfSSatish Balay tao->trust = 0.0; 245a7e14dcfSSatish Balay break; 246a7e14dcfSSatish Balay } 247a7e14dcfSSatish Balay } 248a7e14dcfSSatish Balay 249a7e14dcfSSatish Balay /* Set counter for gradient/reset steps*/ 250a7e14dcfSSatish Balay nlsP->newt = 0; 251a7e14dcfSSatish Balay nlsP->bfgs = 0; 252a7e14dcfSSatish Balay nlsP->grad = 0; 253a7e14dcfSSatish Balay 254a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 2553ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 256e1e80dc8SAlp Dener /* Call general purpose update function */ 257e1e80dc8SAlp Dener if (tao->ops->update) { 2589566063dSJacob Faibussowitsch PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update)); 259e1e80dc8SAlp Dener } 2608931d482SJason Sarich ++tao->niter; 261ae93cb3cSJason Sarich tao->ksp_its = 0; 262a7e14dcfSSatish Balay 263a7e14dcfSSatish Balay /* Compute the Hessian */ 264a7e14dcfSSatish Balay if (needH) { 2659566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre)); 266a7e14dcfSSatish Balay } 267a7e14dcfSSatish Balay 268a7e14dcfSSatish Balay /* Shift the Hessian matrix */ 269a7e14dcfSSatish Balay if (pert > 0) { 2709566063dSJacob Faibussowitsch PetscCall(MatShift(tao->hessian, pert)); 271a7e14dcfSSatish Balay if (tao->hessian != tao->hessian_pre) { 2729566063dSJacob Faibussowitsch PetscCall(MatShift(tao->hessian_pre, pert)); 273a7e14dcfSSatish Balay } 274a7e14dcfSSatish Balay } 275a7e14dcfSSatish Balay 2760c51296cSAlp Dener if (nlsP->bfgs_pre) { 2779566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient)); 278a7e14dcfSSatish Balay ++bfgsUpdates; 279a7e14dcfSSatish Balay } 280a7e14dcfSSatish Balay 281a7e14dcfSSatish Balay /* Solve the Newton system of equations */ 2829566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre)); 2831daac38eSTodd Munson if (is_nash || is_stcg || is_gltr) { 2849566063dSJacob Faibussowitsch PetscCall(KSPCGSetRadius(tao->ksp,nlsP->max_radius)); 2859566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, tao->gradient, nlsP->D)); 2869566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&kspits)); 287a7e14dcfSSatish Balay tao->ksp_its += kspits; 288ae93cb3cSJason Sarich tao->ksp_tot_its += kspits; 2899566063dSJacob Faibussowitsch PetscCall(KSPCGGetNormD(tao->ksp,&norm_d)); 290a7e14dcfSSatish Balay 291a7e14dcfSSatish Balay if (0.0 == tao->trust) { 292a7e14dcfSSatish Balay /* Radius was uninitialized; use the norm of the direction */ 293a7e14dcfSSatish Balay if (norm_d > 0.0) { 294a7e14dcfSSatish Balay tao->trust = norm_d; 295a7e14dcfSSatish Balay 296a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 297a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 298a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 29987f595a5SBarry Smith } else { 300a7e14dcfSSatish Balay /* The direction was bad; set radius to default value and re-solve 301a7e14dcfSSatish Balay the trust-region subproblem to get a direction */ 302a7e14dcfSSatish Balay tao->trust = tao->trust0; 303a7e14dcfSSatish Balay 304a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 305a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 306a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 307a7e14dcfSSatish Balay 3089566063dSJacob Faibussowitsch PetscCall(KSPCGSetRadius(tao->ksp,nlsP->max_radius)); 3099566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, tao->gradient, nlsP->D)); 3109566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&kspits)); 311a7e14dcfSSatish Balay tao->ksp_its += kspits; 312ae93cb3cSJason Sarich tao->ksp_tot_its += kspits; 3139566063dSJacob Faibussowitsch PetscCall(KSPCGGetNormD(tao->ksp,&norm_d)); 314a7e14dcfSSatish Balay 3153c859ba3SBarry Smith PetscCheck(norm_d != 0.0,PETSC_COMM_SELF,PETSC_ERR_PLIB, "Initial direction zero"); 316a7e14dcfSSatish Balay } 317a7e14dcfSSatish Balay } 31887f595a5SBarry Smith } else { 3199566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, tao->gradient, nlsP->D)); 3209566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &kspits)); 321a7e14dcfSSatish Balay tao->ksp_its += kspits; 322ae93cb3cSJason Sarich tao->ksp_tot_its += kspits; 323a7e14dcfSSatish Balay } 3249566063dSJacob Faibussowitsch PetscCall(VecScale(nlsP->D, -1.0)); 3259566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(tao->ksp, &ksp_reason)); 3260c51296cSAlp Dener if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (nlsP->bfgs_pre)) { 327a7e14dcfSSatish Balay /* Preconditioner is numerically indefinite; reset the 328a7e14dcfSSatish Balay approximate if using BFGS preconditioning. */ 3299566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE)); 3309566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient)); 331a7e14dcfSSatish Balay bfgsUpdates = 1; 332a7e14dcfSSatish Balay } 333a7e14dcfSSatish Balay 334a7e14dcfSSatish Balay if (KSP_CONVERGED_ATOL == ksp_reason) { 335a7e14dcfSSatish Balay ++nlsP->ksp_atol; 33687f595a5SBarry Smith } else if (KSP_CONVERGED_RTOL == ksp_reason) { 337a7e14dcfSSatish Balay ++nlsP->ksp_rtol; 33887f595a5SBarry Smith } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) { 339a7e14dcfSSatish Balay ++nlsP->ksp_ctol; 34087f595a5SBarry Smith } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) { 341a7e14dcfSSatish Balay ++nlsP->ksp_negc; 34287f595a5SBarry Smith } else if (KSP_DIVERGED_DTOL == ksp_reason) { 343a7e14dcfSSatish Balay ++nlsP->ksp_dtol; 34487f595a5SBarry Smith } else if (KSP_DIVERGED_ITS == ksp_reason) { 345a7e14dcfSSatish Balay ++nlsP->ksp_iter; 34687f595a5SBarry Smith } else { 347a7e14dcfSSatish Balay ++nlsP->ksp_othr; 348a7e14dcfSSatish Balay } 349a7e14dcfSSatish Balay 350a7e14dcfSSatish Balay /* Check for success (descent direction) */ 3519566063dSJacob Faibussowitsch PetscCall(VecDot(nlsP->D, tao->gradient, &gdx)); 352a7e14dcfSSatish Balay if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 353a7e14dcfSSatish Balay /* Newton step is not descent or direction produced Inf or NaN 354a7e14dcfSSatish Balay Update the perturbation for next time */ 355a7e14dcfSSatish Balay if (pert <= 0.0) { 356a7e14dcfSSatish Balay /* Initialize the perturbation */ 357a7e14dcfSSatish Balay pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 3581daac38eSTodd Munson if (is_gltr) { 3599566063dSJacob Faibussowitsch PetscCall(KSPGLTRGetMinEig(tao->ksp,&e_min)); 360a7e14dcfSSatish Balay pert = PetscMax(pert, -e_min); 361a7e14dcfSSatish Balay } 36287f595a5SBarry Smith } else { 363a7e14dcfSSatish Balay /* Increase the perturbation */ 364a7e14dcfSSatish Balay pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 365a7e14dcfSSatish Balay } 366a7e14dcfSSatish Balay 3670c51296cSAlp Dener if (!nlsP->bfgs_pre) { 368a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and updated 369a7e14dcfSSatish Balay Must use gradient direction in this case */ 3709566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, nlsP->D)); 3719566063dSJacob Faibussowitsch PetscCall(VecScale(nlsP->D, -1.0)); 372a7e14dcfSSatish Balay ++nlsP->grad; 373a7e14dcfSSatish Balay stepType = NLS_GRADIENT; 37487f595a5SBarry Smith } else { 375a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 3769566063dSJacob Faibussowitsch PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D)); 3779566063dSJacob Faibussowitsch PetscCall(VecScale(nlsP->D, -1.0)); 378a7e14dcfSSatish Balay 379a7e14dcfSSatish Balay /* Check for success (descent direction) */ 3809566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, nlsP->D, &gdx)); 381a7e14dcfSSatish Balay if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) { 382a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced not a number 383a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case because 384a7e14dcfSSatish Balay the first solve produces the scaled gradient direction, 385a7e14dcfSSatish Balay which is guaranteed to be descent */ 386a7e14dcfSSatish Balay 387a7e14dcfSSatish Balay /* Use steepest descent direction (scaled) */ 3889566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE)); 3899566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient)); 3909566063dSJacob Faibussowitsch PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D)); 3919566063dSJacob Faibussowitsch PetscCall(VecScale(nlsP->D, -1.0)); 392a7e14dcfSSatish Balay 393a7e14dcfSSatish Balay bfgsUpdates = 1; 3940c51296cSAlp Dener ++nlsP->grad; 3950c51296cSAlp Dener stepType = NLS_GRADIENT; 39687f595a5SBarry Smith } else { 3979566063dSJacob Faibussowitsch PetscCall(MatLMVMGetUpdateCount(nlsP->M, &bfgsUpdates)); 398a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 399a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 4000c51296cSAlp Dener ++nlsP->grad; 4010c51296cSAlp Dener stepType = NLS_GRADIENT; 40287f595a5SBarry Smith } else { 403a7e14dcfSSatish Balay ++nlsP->bfgs; 404a7e14dcfSSatish Balay stepType = NLS_BFGS; 405a7e14dcfSSatish Balay } 406a7e14dcfSSatish Balay } 407a7e14dcfSSatish Balay } 40887f595a5SBarry Smith } else { 409a7e14dcfSSatish Balay /* Computed Newton step is descent */ 410a7e14dcfSSatish Balay switch (ksp_reason) { 411a7e14dcfSSatish Balay case KSP_DIVERGED_NANORINF: 412a7e14dcfSSatish Balay case KSP_DIVERGED_BREAKDOWN: 413a7e14dcfSSatish Balay case KSP_DIVERGED_INDEFINITE_MAT: 414a7e14dcfSSatish Balay case KSP_DIVERGED_INDEFINITE_PC: 415a7e14dcfSSatish Balay case KSP_CONVERGED_CG_NEG_CURVE: 416a7e14dcfSSatish Balay /* Matrix or preconditioner is indefinite; increase perturbation */ 417a7e14dcfSSatish Balay if (pert <= 0.0) { 418a7e14dcfSSatish Balay /* Initialize the perturbation */ 419a7e14dcfSSatish Balay pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 4201daac38eSTodd Munson if (is_gltr) { 4219566063dSJacob Faibussowitsch PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min)); 422a7e14dcfSSatish Balay pert = PetscMax(pert, -e_min); 423a7e14dcfSSatish Balay } 42487f595a5SBarry Smith } else { 425a7e14dcfSSatish Balay /* Increase the perturbation */ 426a7e14dcfSSatish Balay pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 427a7e14dcfSSatish Balay } 428a7e14dcfSSatish Balay break; 429a7e14dcfSSatish Balay 430a7e14dcfSSatish Balay default: 431a7e14dcfSSatish Balay /* Newton step computation is good; decrease perturbation */ 432a7e14dcfSSatish Balay pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm); 433a7e14dcfSSatish Balay if (pert < nlsP->pmin) { 434a7e14dcfSSatish Balay pert = 0.0; 435a7e14dcfSSatish Balay } 436a7e14dcfSSatish Balay break; 437a7e14dcfSSatish Balay } 438a7e14dcfSSatish Balay 439a7e14dcfSSatish Balay ++nlsP->newt; 440a7e14dcfSSatish Balay stepType = NLS_NEWTON; 441a7e14dcfSSatish Balay } 442a7e14dcfSSatish Balay 443a7e14dcfSSatish Balay /* Perform the linesearch */ 444a7e14dcfSSatish Balay fold = f; 4459566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, nlsP->Xold)); 4469566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, nlsP->Gold)); 447a7e14dcfSSatish Balay 4489566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason)); 4499566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 450a7e14dcfSSatish Balay 45187f595a5SBarry Smith while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) { /* Linesearch failed */ 452a7e14dcfSSatish Balay f = fold; 4539566063dSJacob Faibussowitsch PetscCall(VecCopy(nlsP->Xold, tao->solution)); 4549566063dSJacob Faibussowitsch PetscCall(VecCopy(nlsP->Gold, tao->gradient)); 455a7e14dcfSSatish Balay 456a7e14dcfSSatish Balay switch(stepType) { 457a7e14dcfSSatish Balay case NLS_NEWTON: 458a7e14dcfSSatish Balay /* Failed to obtain acceptable iterate with Newton 1step 459a7e14dcfSSatish Balay Update the perturbation for next time */ 460a7e14dcfSSatish Balay if (pert <= 0.0) { 461a7e14dcfSSatish Balay /* Initialize the perturbation */ 462a7e14dcfSSatish Balay pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 4631daac38eSTodd Munson if (is_gltr) { 4649566063dSJacob Faibussowitsch PetscCall(KSPGLTRGetMinEig(tao->ksp,&e_min)); 465a7e14dcfSSatish Balay pert = PetscMax(pert, -e_min); 466a7e14dcfSSatish Balay } 46787f595a5SBarry Smith } else { 468a7e14dcfSSatish Balay /* Increase the perturbation */ 469a7e14dcfSSatish Balay pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 470a7e14dcfSSatish Balay } 471a7e14dcfSSatish Balay 4720c51296cSAlp Dener if (!nlsP->bfgs_pre) { 473a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and being updated 474a7e14dcfSSatish Balay Must use gradient direction in this case */ 4759566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, nlsP->D)); 476a7e14dcfSSatish Balay ++nlsP->grad; 477a7e14dcfSSatish Balay stepType = NLS_GRADIENT; 47887f595a5SBarry Smith } else { 479a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 4809566063dSJacob Faibussowitsch PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D)); 481a7e14dcfSSatish Balay /* Check for success (descent direction) */ 4829566063dSJacob Faibussowitsch PetscCall(VecDot(tao->solution, nlsP->D, &gdx)); 483a7e14dcfSSatish Balay if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 484a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced not a number 485a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case 486a7e14dcfSSatish Balay Use steepest descent direction (scaled) */ 4879566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE)); 4889566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient)); 4899566063dSJacob Faibussowitsch PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D)); 490a7e14dcfSSatish Balay 491a7e14dcfSSatish Balay bfgsUpdates = 1; 4920c51296cSAlp Dener ++nlsP->grad; 4930c51296cSAlp Dener stepType = NLS_GRADIENT; 49487f595a5SBarry Smith } else { 495a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 496a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 4970c51296cSAlp Dener ++nlsP->grad; 4980c51296cSAlp Dener stepType = NLS_GRADIENT; 49987f595a5SBarry Smith } else { 500a7e14dcfSSatish Balay ++nlsP->bfgs; 501a7e14dcfSSatish Balay stepType = NLS_BFGS; 502a7e14dcfSSatish Balay } 503a7e14dcfSSatish Balay } 504a7e14dcfSSatish Balay } 505a7e14dcfSSatish Balay break; 506a7e14dcfSSatish Balay 507a7e14dcfSSatish Balay case NLS_BFGS: 508a7e14dcfSSatish Balay /* Can only enter if pc_type == NLS_PC_BFGS 509a7e14dcfSSatish Balay Failed to obtain acceptable iterate with BFGS step 510a7e14dcfSSatish Balay Attempt to use the scaled gradient direction */ 5119566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE)); 5129566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient)); 5139566063dSJacob Faibussowitsch PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D)); 514a7e14dcfSSatish Balay 515a7e14dcfSSatish Balay bfgsUpdates = 1; 516a7e14dcfSSatish Balay ++nlsP->grad; 517a7e14dcfSSatish Balay stepType = NLS_GRADIENT; 518a7e14dcfSSatish Balay break; 519a7e14dcfSSatish Balay } 5209566063dSJacob Faibussowitsch PetscCall(VecScale(nlsP->D, -1.0)); 521a7e14dcfSSatish Balay 5229566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason)); 5239566063dSJacob Faibussowitsch PetscCall(TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full)); 5249566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 525a7e14dcfSSatish Balay } 526a7e14dcfSSatish Balay 52787f595a5SBarry Smith if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 528a7e14dcfSSatish Balay /* Failed to find an improving point */ 529a7e14dcfSSatish Balay f = fold; 5309566063dSJacob Faibussowitsch PetscCall(VecCopy(nlsP->Xold, tao->solution)); 5319566063dSJacob Faibussowitsch PetscCall(VecCopy(nlsP->Gold, tao->gradient)); 532a7e14dcfSSatish Balay step = 0.0; 533a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 534a7e14dcfSSatish Balay break; 535a7e14dcfSSatish Balay } 536a7e14dcfSSatish Balay 537a7e14dcfSSatish Balay /* Update trust region radius */ 5381daac38eSTodd Munson if (is_nash || is_stcg || is_gltr) { 539a7e14dcfSSatish Balay switch(nlsP->update_type) { 540a7e14dcfSSatish Balay case NLS_UPDATE_STEP: 541a7e14dcfSSatish Balay if (stepType == NLS_NEWTON) { 542a7e14dcfSSatish Balay if (step < nlsP->nu1) { 543a7e14dcfSSatish Balay /* Very bad step taken; reduce radius */ 544a7e14dcfSSatish Balay tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust); 54587f595a5SBarry Smith } else if (step < nlsP->nu2) { 546a7e14dcfSSatish Balay /* Reasonably bad step taken; reduce radius */ 547a7e14dcfSSatish Balay tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust); 54887f595a5SBarry Smith } else if (step < nlsP->nu3) { 549a7e14dcfSSatish Balay /* Reasonable step was taken; leave radius alone */ 550a7e14dcfSSatish Balay if (nlsP->omega3 < 1.0) { 551a7e14dcfSSatish Balay tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust); 55287f595a5SBarry Smith } else if (nlsP->omega3 > 1.0) { 553a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust); 554a7e14dcfSSatish Balay } 55587f595a5SBarry Smith } else if (step < nlsP->nu4) { 556a7e14dcfSSatish Balay /* Full step taken; increase the radius */ 557a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust); 55887f595a5SBarry Smith } else { 559a7e14dcfSSatish Balay /* More than full step taken; increase the radius */ 560a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust); 561a7e14dcfSSatish Balay } 56287f595a5SBarry Smith } else { 563a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 564a7e14dcfSSatish Balay tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust); 565a7e14dcfSSatish Balay } 566a7e14dcfSSatish Balay break; 567a7e14dcfSSatish Balay 568a7e14dcfSSatish Balay case NLS_UPDATE_REDUCTION: 569a7e14dcfSSatish Balay if (stepType == NLS_NEWTON) { 570a7e14dcfSSatish Balay /* Get predicted reduction */ 571a7e14dcfSSatish Balay 5729566063dSJacob Faibussowitsch PetscCall(KSPCGGetObjFcn(tao->ksp,&prered)); 573a7e14dcfSSatish Balay if (prered >= 0.0) { 574a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot */ 575a7e14dcfSSatish Balay /* happen in infinite precision arithmetic. Step should */ 576a7e14dcfSSatish Balay /* be rejected! */ 577a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 57887f595a5SBarry Smith } else { 579a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(f_full)) { 580a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 58187f595a5SBarry Smith } else { 582a7e14dcfSSatish Balay /* Compute and actual reduction */ 583a7e14dcfSSatish Balay actred = fold - f_full; 584a7e14dcfSSatish Balay prered = -prered; 585a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= nlsP->epsilon) && 586a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= nlsP->epsilon)) { 587a7e14dcfSSatish Balay kappa = 1.0; 58887f595a5SBarry Smith } else { 589a7e14dcfSSatish Balay kappa = actred / prered; 590a7e14dcfSSatish Balay } 591a7e14dcfSSatish Balay 592a7e14dcfSSatish Balay /* Accept of reject the step and update radius */ 593a7e14dcfSSatish Balay if (kappa < nlsP->eta1) { 594a7e14dcfSSatish Balay /* Very bad step */ 595a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 59687f595a5SBarry Smith } else if (kappa < nlsP->eta2) { 597a7e14dcfSSatish Balay /* Marginal bad step */ 598a7e14dcfSSatish Balay tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d); 59987f595a5SBarry Smith } else if (kappa < nlsP->eta3) { 600a7e14dcfSSatish Balay /* Reasonable step */ 601a7e14dcfSSatish Balay if (nlsP->alpha3 < 1.0) { 602a7e14dcfSSatish Balay tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust); 60387f595a5SBarry Smith } else if (nlsP->alpha3 > 1.0) { 604a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust); 605a7e14dcfSSatish Balay } 60687f595a5SBarry Smith } else if (kappa < nlsP->eta4) { 607a7e14dcfSSatish Balay /* Good step */ 608a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust); 60987f595a5SBarry Smith } else { 610a7e14dcfSSatish Balay /* Very good step */ 611a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust); 612a7e14dcfSSatish Balay } 613a7e14dcfSSatish Balay } 614a7e14dcfSSatish Balay } 61587f595a5SBarry Smith } else { 616a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 617a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust); 618a7e14dcfSSatish Balay } 619a7e14dcfSSatish Balay break; 620a7e14dcfSSatish Balay 621a7e14dcfSSatish Balay default: 622a7e14dcfSSatish Balay if (stepType == NLS_NEWTON) { 6239566063dSJacob Faibussowitsch PetscCall(KSPCGGetObjFcn(tao->ksp,&prered)); 624a7e14dcfSSatish Balay if (prered >= 0.0) { 625a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot */ 626a7e14dcfSSatish Balay /* happen in infinite precision arithmetic. Step should */ 627a7e14dcfSSatish Balay /* be rejected! */ 628a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 62987f595a5SBarry Smith } else { 630a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(f_full)) { 631a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 63287f595a5SBarry Smith } else { 633a7e14dcfSSatish Balay actred = fold - f_full; 634a7e14dcfSSatish Balay prered = -prered; 63587f595a5SBarry Smith if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) { 636a7e14dcfSSatish Balay kappa = 1.0; 63787f595a5SBarry Smith } else { 638a7e14dcfSSatish Balay kappa = actred / prered; 639a7e14dcfSSatish Balay } 640a7e14dcfSSatish Balay 641a7e14dcfSSatish Balay tau_1 = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred); 642a7e14dcfSSatish Balay tau_2 = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred); 643a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 644a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 645a7e14dcfSSatish Balay 646a7e14dcfSSatish Balay if (kappa >= 1.0 - nlsP->mu1) { 647a7e14dcfSSatish Balay /* Great agreement */ 648a7e14dcfSSatish Balay if (tau_max < 1.0) { 649a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d); 65087f595a5SBarry Smith } else if (tau_max > nlsP->gamma4) { 651a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d); 65287f595a5SBarry Smith } else { 653a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 654a7e14dcfSSatish Balay } 65587f595a5SBarry Smith } else if (kappa >= 1.0 - nlsP->mu2) { 656a7e14dcfSSatish Balay /* Good agreement */ 657a7e14dcfSSatish Balay 658a7e14dcfSSatish Balay if (tau_max < nlsP->gamma2) { 659a7e14dcfSSatish Balay tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d); 66087f595a5SBarry Smith } else if (tau_max > nlsP->gamma3) { 661a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d); 66287f595a5SBarry Smith } else if (tau_max < 1.0) { 663a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 66487f595a5SBarry Smith } else { 665a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 666a7e14dcfSSatish Balay } 66787f595a5SBarry Smith } else { 668a7e14dcfSSatish Balay /* Not good agreement */ 669a7e14dcfSSatish Balay if (tau_min > 1.0) { 670a7e14dcfSSatish Balay tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d); 67187f595a5SBarry Smith } else if (tau_max < nlsP->gamma1) { 672a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 67387f595a5SBarry Smith } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) { 674a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 67587f595a5SBarry Smith } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) { 676a7e14dcfSSatish Balay tao->trust = tau_1 * PetscMin(tao->trust, norm_d); 67787f595a5SBarry Smith } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) { 678a7e14dcfSSatish Balay tao->trust = tau_2 * PetscMin(tao->trust, norm_d); 67987f595a5SBarry Smith } else { 680a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 681a7e14dcfSSatish Balay } 682a7e14dcfSSatish Balay } 683a7e14dcfSSatish Balay } 684a7e14dcfSSatish Balay } 68587f595a5SBarry Smith } else { 686a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 687a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust); 688a7e14dcfSSatish Balay } 689a7e14dcfSSatish Balay break; 690a7e14dcfSSatish Balay } 691a7e14dcfSSatish Balay 692a7e14dcfSSatish Balay /* The radius may have been increased; modify if it is too large */ 693a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 694a7e14dcfSSatish Balay } 695a7e14dcfSSatish Balay 696a7e14dcfSSatish Balay /* Check for termination */ 6979566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm)); 6983c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER,"User provided compute function generated Not-a-Number"); 699a7e14dcfSSatish Balay needH = 1; 7009566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its)); 7019566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,step)); 7029566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 703a7e14dcfSSatish Balay } 704a7e14dcfSSatish Balay PetscFunctionReturn(0); 705a7e14dcfSSatish Balay } 706a7e14dcfSSatish Balay 707a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 708441846f8SBarry Smith static PetscErrorCode TaoSetUp_NLS(Tao tao) 709a7e14dcfSSatish Balay { 710a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 711a7e14dcfSSatish Balay 712a7e14dcfSSatish Balay PetscFunctionBegin; 7139566063dSJacob Faibussowitsch if (!tao->gradient) PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 7149566063dSJacob Faibussowitsch if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution,&tao->stepdirection)); 7159566063dSJacob Faibussowitsch if (!nlsP->W) PetscCall(VecDuplicate(tao->solution,&nlsP->W)); 7169566063dSJacob Faibussowitsch if (!nlsP->D) PetscCall(VecDuplicate(tao->solution,&nlsP->D)); 7179566063dSJacob Faibussowitsch if (!nlsP->Xold) PetscCall(VecDuplicate(tao->solution,&nlsP->Xold)); 7189566063dSJacob Faibussowitsch if (!nlsP->Gold) PetscCall(VecDuplicate(tao->solution,&nlsP->Gold)); 71983c8fe1dSLisandro Dalcin nlsP->M = NULL; 72083c8fe1dSLisandro Dalcin nlsP->bfgs_pre = NULL; 721a7e14dcfSSatish Balay PetscFunctionReturn(0); 722a7e14dcfSSatish Balay } 723a7e14dcfSSatish Balay 724a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 725441846f8SBarry Smith static PetscErrorCode TaoDestroy_NLS(Tao tao) 726a7e14dcfSSatish Balay { 727a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 728a7e14dcfSSatish Balay 729a7e14dcfSSatish Balay PetscFunctionBegin; 730a7e14dcfSSatish Balay if (tao->setupcalled) { 7319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&nlsP->D)); 7329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&nlsP->W)); 7339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&nlsP->Xold)); 7349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&nlsP->Gold)); 735a7e14dcfSSatish Balay } 736*a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 7379566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 738a7e14dcfSSatish Balay PetscFunctionReturn(0); 739a7e14dcfSSatish Balay } 740a7e14dcfSSatish Balay 741a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 7424416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_NLS(PetscOptionItems *PetscOptionsObject,Tao tao) 743a7e14dcfSSatish Balay { 744a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 745a7e14dcfSSatish Balay 746a7e14dcfSSatish Balay PetscFunctionBegin; 747d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Newton line search method for unconstrained optimization"); 7489566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-tao_nls_init_type", "radius initialization type", "", NLS_INIT, NLS_INIT_TYPES, NLS_INIT[nlsP->init_type], &nlsP->init_type, NULL)); 7499566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-tao_nls_update_type", "radius update type", "", NLS_UPDATE, NLS_UPDATE_TYPES, NLS_UPDATE[nlsP->update_type], &nlsP->update_type, NULL)); 7509566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval,NULL)); 7519566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin,NULL)); 7529566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax,NULL)); 7539566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac,NULL)); 7549566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin,NULL)); 7559566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax,NULL)); 7569566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac,NULL)); 7579566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac,NULL)); 7589566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac,NULL)); 7599566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac,NULL)); 7609566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1,NULL)); 7619566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2,NULL)); 7629566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3,NULL)); 7639566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4,NULL)); 7649566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1,NULL)); 7659566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2,NULL)); 7669566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3,NULL)); 7679566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4,NULL)); 7689566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5,NULL)); 7699566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1,NULL)); 7709566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2,NULL)); 7719566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3,NULL)); 7729566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4,NULL)); 7739566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1,NULL)); 7749566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2,NULL)); 7759566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3,NULL)); 7769566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4,NULL)); 7779566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5,NULL)); 7789566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i,NULL)); 7799566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i,NULL)); 7809566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i,NULL)); 7819566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i,NULL)); 7829566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i,NULL)); 7839566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i,NULL)); 7849566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i,NULL)); 7859566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1,NULL)); 7869566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2,NULL)); 7879566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1,NULL)); 7889566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2,NULL)); 7899566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3,NULL)); 7909566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4,NULL)); 7919566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta,NULL)); 7929566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius,NULL)); 7939566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius,NULL)); 7949566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon,NULL)); 795d0609cedSBarry Smith PetscOptionsHeadEnd(); 7969566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 7979566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp)); 798a7e14dcfSSatish Balay PetscFunctionReturn(0); 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 807a7e14dcfSSatish Balay PetscFunctionBegin; 8089566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); 809a7e14dcfSSatish Balay if (isascii) { 8109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 81163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Newton steps: %" PetscInt_FMT "\n", nlsP->newt)); 81263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS steps: %" PetscInt_FMT "\n", nlsP->bfgs)); 81363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", nlsP->grad)); 814a7e14dcfSSatish Balay 81563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp atol: %" PetscInt_FMT "\n", nlsP->ksp_atol)); 81663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %" PetscInt_FMT "\n", nlsP->ksp_rtol)); 81763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %" PetscInt_FMT "\n", nlsP->ksp_ctol)); 81863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp negc: %" PetscInt_FMT "\n", nlsP->ksp_negc)); 81963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %" PetscInt_FMT "\n", nlsP->ksp_dtol)); 82063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp iter: %" PetscInt_FMT "\n", nlsP->ksp_iter)); 82163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp othr: %" PetscInt_FMT "\n", nlsP->ksp_othr)); 8229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 823a7e14dcfSSatish Balay } 824a7e14dcfSSatish Balay PetscFunctionReturn(0); 825a7e14dcfSSatish Balay } 826a7e14dcfSSatish Balay 827a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 8284aa34175SJason Sarich /*MC 8294aa34175SJason Sarich TAONLS - Newton's method with linesearch for unconstrained minimization. 8304aa34175SJason Sarich At each iteration, the Newton line search method solves the symmetric 8314aa34175SJason Sarich system of equations to obtain the step diretion dk: 8324aa34175SJason Sarich Hk dk = -gk 8334aa34175SJason Sarich a More-Thuente line search is applied on the direction dk to approximately 8344aa34175SJason Sarich solve 8354aa34175SJason Sarich min_t f(xk + t d_k) 8364aa34175SJason Sarich 8374aa34175SJason Sarich Options Database Keys: 8389d0a60b2SAlp Dener + -tao_nls_init_type - "constant","direction","interpolation" 8394aa34175SJason Sarich . -tao_nls_update_type - "step","direction","interpolation" 8404aa34175SJason Sarich . -tao_nls_sval - perturbation starting value 8414aa34175SJason Sarich . -tao_nls_imin - minimum initial perturbation 8424aa34175SJason Sarich . -tao_nls_imax - maximum initial perturbation 8434aa34175SJason Sarich . -tao_nls_pmin - minimum perturbation 8444aa34175SJason Sarich . -tao_nls_pmax - maximum perturbation 8454aa34175SJason Sarich . -tao_nls_pgfac - growth factor 8464aa34175SJason Sarich . -tao_nls_psfac - shrink factor 8474aa34175SJason Sarich . -tao_nls_imfac - initial merit factor 8484aa34175SJason Sarich . -tao_nls_pmgfac - merit growth factor 8494aa34175SJason Sarich . -tao_nls_pmsfac - merit shrink factor 8504aa34175SJason Sarich . -tao_nls_eta1 - poor steplength; reduce radius 8514aa34175SJason Sarich . -tao_nls_eta2 - reasonable steplength; leave radius 8524aa34175SJason Sarich . -tao_nls_eta3 - good steplength; increase readius 8534aa34175SJason Sarich . -tao_nls_eta4 - excellent steplength; greatly increase radius 8544aa34175SJason Sarich . -tao_nls_alpha1 - alpha1 reduction 8554aa34175SJason Sarich . -tao_nls_alpha2 - alpha2 reduction 8564aa34175SJason Sarich . -tao_nls_alpha3 - alpha3 reduction 8574aa34175SJason Sarich . -tao_nls_alpha4 - alpha4 reduction 8584aa34175SJason Sarich . -tao_nls_alpha - alpha5 reduction 8594aa34175SJason Sarich . -tao_nls_mu1 - mu1 interpolation update 8604aa34175SJason Sarich . -tao_nls_mu2 - mu2 interpolation update 8614aa34175SJason Sarich . -tao_nls_gamma1 - gamma1 interpolation update 8624aa34175SJason Sarich . -tao_nls_gamma2 - gamma2 interpolation update 8634aa34175SJason Sarich . -tao_nls_gamma3 - gamma3 interpolation update 8644aa34175SJason Sarich . -tao_nls_gamma4 - gamma4 interpolation update 8654aa34175SJason Sarich . -tao_nls_theta - theta interpolation update 8664aa34175SJason Sarich . -tao_nls_omega1 - omega1 step update 8674aa34175SJason Sarich . -tao_nls_omega2 - omega2 step update 8684aa34175SJason Sarich . -tao_nls_omega3 - omega3 step update 8694aa34175SJason Sarich . -tao_nls_omega4 - omega4 step update 8704aa34175SJason Sarich . -tao_nls_omega5 - omega5 step update 8711522df2eSJason Sarich . -tao_nls_mu1_i - mu1 interpolation init factor 8721522df2eSJason Sarich . -tao_nls_mu2_i - mu2 interpolation init factor 8731522df2eSJason Sarich . -tao_nls_gamma1_i - gamma1 interpolation init factor 8741522df2eSJason Sarich . -tao_nls_gamma2_i - gamma2 interpolation init factor 8751522df2eSJason Sarich . -tao_nls_gamma3_i - gamma3 interpolation init factor 8761522df2eSJason Sarich . -tao_nls_gamma4_i - gamma4 interpolation init factor 8771522df2eSJason Sarich - -tao_nls_theta_i - theta interpolation init factor 8781eb8069cSJason Sarich 8791eb8069cSJason Sarich Level: beginner 8801522df2eSJason Sarich M*/ 8814aa34175SJason Sarich 882728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NLS(Tao tao) 883a7e14dcfSSatish Balay { 884a7e14dcfSSatish Balay TAO_NLS *nlsP; 8858caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 886a7e14dcfSSatish Balay 887a7e14dcfSSatish Balay PetscFunctionBegin; 8889566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&nlsP)); 889a7e14dcfSSatish Balay 890a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_NLS; 891a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_NLS; 892a7e14dcfSSatish Balay tao->ops->view = TaoView_NLS; 893a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_NLS; 894a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_NLS; 895a7e14dcfSSatish Balay 8966552cf8aSJason Sarich /* Override default settings (unless already changed) */ 8976552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 50; 8986552cf8aSJason Sarich if (!tao->trust0_changed) tao->trust0 = 100.0; 8996552cf8aSJason Sarich 900a7e14dcfSSatish Balay tao->data = (void*)nlsP; 901a7e14dcfSSatish Balay 902a7e14dcfSSatish Balay nlsP->sval = 0.0; 903a7e14dcfSSatish Balay nlsP->imin = 1.0e-4; 904a7e14dcfSSatish Balay nlsP->imax = 1.0e+2; 905a7e14dcfSSatish Balay nlsP->imfac = 1.0e-1; 906a7e14dcfSSatish Balay 907a7e14dcfSSatish Balay nlsP->pmin = 1.0e-12; 908a7e14dcfSSatish Balay nlsP->pmax = 1.0e+2; 909a7e14dcfSSatish Balay nlsP->pgfac = 1.0e+1; 910a7e14dcfSSatish Balay nlsP->psfac = 4.0e-1; 911a7e14dcfSSatish Balay nlsP->pmgfac = 1.0e-1; 912a7e14dcfSSatish Balay nlsP->pmsfac = 1.0e-1; 913a7e14dcfSSatish Balay 914a7e14dcfSSatish Balay /* Default values for trust-region radius update based on steplength */ 915a7e14dcfSSatish Balay nlsP->nu1 = 0.25; 916a7e14dcfSSatish Balay nlsP->nu2 = 0.50; 917a7e14dcfSSatish Balay nlsP->nu3 = 1.00; 918a7e14dcfSSatish Balay nlsP->nu4 = 1.25; 919a7e14dcfSSatish Balay 920a7e14dcfSSatish Balay nlsP->omega1 = 0.25; 921a7e14dcfSSatish Balay nlsP->omega2 = 0.50; 922a7e14dcfSSatish Balay nlsP->omega3 = 1.00; 923a7e14dcfSSatish Balay nlsP->omega4 = 2.00; 924a7e14dcfSSatish Balay nlsP->omega5 = 4.00; 925a7e14dcfSSatish Balay 926a7e14dcfSSatish Balay /* Default values for trust-region radius update based on reduction */ 927a7e14dcfSSatish Balay nlsP->eta1 = 1.0e-4; 928a7e14dcfSSatish Balay nlsP->eta2 = 0.25; 929a7e14dcfSSatish Balay nlsP->eta3 = 0.50; 930a7e14dcfSSatish Balay nlsP->eta4 = 0.90; 931a7e14dcfSSatish Balay 932a7e14dcfSSatish Balay nlsP->alpha1 = 0.25; 933a7e14dcfSSatish Balay nlsP->alpha2 = 0.50; 934a7e14dcfSSatish Balay nlsP->alpha3 = 1.00; 935a7e14dcfSSatish Balay nlsP->alpha4 = 2.00; 936a7e14dcfSSatish Balay nlsP->alpha5 = 4.00; 937a7e14dcfSSatish Balay 938a7e14dcfSSatish Balay /* Default values for trust-region radius update based on interpolation */ 939a7e14dcfSSatish Balay nlsP->mu1 = 0.10; 940a7e14dcfSSatish Balay nlsP->mu2 = 0.50; 941a7e14dcfSSatish Balay 942a7e14dcfSSatish Balay nlsP->gamma1 = 0.25; 943a7e14dcfSSatish Balay nlsP->gamma2 = 0.50; 944a7e14dcfSSatish Balay nlsP->gamma3 = 2.00; 945a7e14dcfSSatish Balay nlsP->gamma4 = 4.00; 946a7e14dcfSSatish Balay 947a7e14dcfSSatish Balay nlsP->theta = 0.05; 948a7e14dcfSSatish Balay 949a7e14dcfSSatish Balay /* Default values for trust region initialization based on interpolation */ 950a7e14dcfSSatish Balay nlsP->mu1_i = 0.35; 951a7e14dcfSSatish Balay nlsP->mu2_i = 0.50; 952a7e14dcfSSatish Balay 953a7e14dcfSSatish Balay nlsP->gamma1_i = 0.0625; 954a7e14dcfSSatish Balay nlsP->gamma2_i = 0.5; 955a7e14dcfSSatish Balay nlsP->gamma3_i = 2.0; 956a7e14dcfSSatish Balay nlsP->gamma4_i = 5.0; 957a7e14dcfSSatish Balay 958a7e14dcfSSatish Balay nlsP->theta_i = 0.25; 959a7e14dcfSSatish Balay 960a7e14dcfSSatish Balay /* Remaining parameters */ 961a7e14dcfSSatish Balay nlsP->min_radius = 1.0e-10; 962a7e14dcfSSatish Balay nlsP->max_radius = 1.0e10; 963a7e14dcfSSatish Balay nlsP->epsilon = 1.0e-6; 964a7e14dcfSSatish Balay 965a7e14dcfSSatish Balay nlsP->init_type = NLS_INIT_INTERPOLATION; 966a7e14dcfSSatish Balay nlsP->update_type = NLS_UPDATE_STEP; 967a7e14dcfSSatish Balay 9689566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch)); 9699566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch,(PetscObject)tao,1)); 9709566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch,morethuente_type)); 9719566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch,tao)); 9729566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix)); 973a7e14dcfSSatish Balay 974a7e14dcfSSatish Balay /* Set linear solver to default for symmetric matrices */ 9759566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm,&tao->ksp)); 9769566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp,(PetscObject)tao,1)); 9779566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix)); 9789566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(tao->ksp,"tao_nls_")); 9799566063dSJacob Faibussowitsch PetscCall(KSPSetType(tao->ksp,KSPSTCG)); 980a7e14dcfSSatish Balay PetscFunctionReturn(0); 981a7e14dcfSSatish Balay } 982