1ba92ff59SBarry Smith #include <petsctaolinesearch.h> 2aaa7dc30SBarry Smith #include <../src/tao/matrix/lmvmmat.h> 3aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/nls/nls.h> 4a7e14dcfSSatish Balay 5aaa7dc30SBarry Smith #include <petscksp.h> 6aaa7dc30SBarry Smith #include <petscpc.h> 7aaa7dc30SBarry Smith #include <petsc-private/kspimpl.h> 8aaa7dc30SBarry Smith #include <petsc-private/pcimpl.h> 9a7e14dcfSSatish Balay 10a7e14dcfSSatish Balay #define NLS_KSP_CG 0 11a7e14dcfSSatish Balay #define NLS_KSP_NASH 1 12a7e14dcfSSatish Balay #define NLS_KSP_STCG 2 13a7e14dcfSSatish Balay #define NLS_KSP_GLTR 3 14a7e14dcfSSatish Balay #define NLS_KSP_PETSC 4 15a7e14dcfSSatish Balay #define NLS_KSP_TYPES 5 16a7e14dcfSSatish Balay 17a7e14dcfSSatish Balay #define NLS_PC_NONE 0 18a7e14dcfSSatish Balay #define NLS_PC_AHESS 1 19a7e14dcfSSatish Balay #define NLS_PC_BFGS 2 20a7e14dcfSSatish Balay #define NLS_PC_PETSC 3 21a7e14dcfSSatish Balay #define NLS_PC_TYPES 4 22a7e14dcfSSatish Balay 23a7e14dcfSSatish Balay #define BFGS_SCALE_AHESS 0 24a7e14dcfSSatish Balay #define BFGS_SCALE_PHESS 1 25a7e14dcfSSatish Balay #define BFGS_SCALE_BFGS 2 26a7e14dcfSSatish Balay #define BFGS_SCALE_TYPES 3 27a7e14dcfSSatish Balay 28a7e14dcfSSatish Balay #define NLS_INIT_CONSTANT 0 29a7e14dcfSSatish Balay #define NLS_INIT_DIRECTION 1 30a7e14dcfSSatish Balay #define NLS_INIT_INTERPOLATION 2 31a7e14dcfSSatish Balay #define NLS_INIT_TYPES 3 32a7e14dcfSSatish Balay 33a7e14dcfSSatish Balay #define NLS_UPDATE_STEP 0 34a7e14dcfSSatish Balay #define NLS_UPDATE_REDUCTION 1 35a7e14dcfSSatish Balay #define NLS_UPDATE_INTERPOLATION 2 36a7e14dcfSSatish Balay #define NLS_UPDATE_TYPES 3 37a7e14dcfSSatish Balay 3887f595a5SBarry Smith static const char *NLS_KSP[64] = {"cg", "nash", "stcg", "gltr", "petsc"}; 39a7e14dcfSSatish Balay 4087f595a5SBarry Smith static const char *NLS_PC[64] = {"none", "ahess", "bfgs", "petsc"}; 41a7e14dcfSSatish Balay 4287f595a5SBarry Smith static const char *BFGS_SCALE[64] = {"ahess", "phess", "bfgs"}; 43a7e14dcfSSatish Balay 4487f595a5SBarry Smith static const char *NLS_INIT[64] = {"constant", "direction", "interpolation"}; 45a7e14dcfSSatish Balay 4687f595a5SBarry Smith static const char *NLS_UPDATE[64] = {"step", "reduction", "interpolation"}; 47a7e14dcfSSatish Balay 48a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x); 49a7e14dcfSSatish Balay /* Routine for BFGS preconditioner 50a7e14dcfSSatish Balay 51a7e14dcfSSatish Balay 52a7e14dcfSSatish Balay Implements Newton's Method with a line search approach for solving 53a7e14dcfSSatish Balay unconstrained minimization problems. A More'-Thuente line search 54a7e14dcfSSatish Balay is used to guarantee that the bfgs preconditioner remains positive 55a7e14dcfSSatish Balay definite. 56a7e14dcfSSatish Balay 57a7e14dcfSSatish Balay The method can shift the Hessian matrix. The shifting procedure is 58a7e14dcfSSatish Balay adapted from the PATH algorithm for solving complementarity 59a7e14dcfSSatish Balay problems. 60a7e14dcfSSatish Balay 61a7e14dcfSSatish Balay The linear system solve should be done with a conjugate gradient 62a7e14dcfSSatish Balay method, although any method can be used. */ 63a7e14dcfSSatish Balay 64a7e14dcfSSatish Balay #define NLS_NEWTON 0 65a7e14dcfSSatish Balay #define NLS_BFGS 1 66a7e14dcfSSatish Balay #define NLS_SCALED_GRADIENT 2 67a7e14dcfSSatish Balay #define NLS_GRADIENT 3 68a7e14dcfSSatish Balay 69a7e14dcfSSatish Balay #undef __FUNCT__ 70a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_NLS" 71441846f8SBarry Smith static PetscErrorCode TaoSolve_NLS(Tao tao) 72a7e14dcfSSatish Balay { 73a7e14dcfSSatish Balay PetscErrorCode ierr; 74a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 75a7e14dcfSSatish Balay PC pc; 76a7e14dcfSSatish Balay KSPConvergedReason ksp_reason; 77e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 78e4cb33bbSBarry Smith TaoConvergedReason reason; 79a7e14dcfSSatish Balay 80a7e14dcfSSatish Balay PetscReal fmin, ftrial, f_full, prered, actred, kappa, sigma; 81a7e14dcfSSatish Balay PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 82a7e14dcfSSatish Balay PetscReal f, fold, gdx, gnorm, pert; 83a7e14dcfSSatish Balay PetscReal step = 1.0; 84a7e14dcfSSatish Balay PetscReal delta; 85a7e14dcfSSatish Balay PetscReal norm_d = 0.0, e_min; 86a7e14dcfSSatish Balay 87a7e14dcfSSatish Balay PetscInt stepType; 88a7e14dcfSSatish Balay PetscInt iter = 0; 89a7e14dcfSSatish Balay PetscInt bfgsUpdates = 0; 90a7e14dcfSSatish Balay PetscInt n,N,kspits; 91a7e14dcfSSatish Balay PetscInt needH; 92a7e14dcfSSatish Balay 93a7e14dcfSSatish Balay PetscInt i_max = 5; 94a7e14dcfSSatish Balay PetscInt j_max = 1; 95a7e14dcfSSatish Balay PetscInt i, j; 96a7e14dcfSSatish Balay 97a7e14dcfSSatish Balay PetscFunctionBegin; 98a7e14dcfSSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 99a7e14dcfSSatish Balay ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by nls algorithm\n");CHKERRQ(ierr); 100a7e14dcfSSatish Balay } 101a7e14dcfSSatish Balay 102a7e14dcfSSatish Balay /* Initialized variables */ 103a7e14dcfSSatish Balay pert = nlsP->sval; 104a7e14dcfSSatish Balay 105a7e14dcfSSatish Balay nlsP->ksp_atol = 0; 106a7e14dcfSSatish Balay nlsP->ksp_rtol = 0; 107a7e14dcfSSatish Balay nlsP->ksp_dtol = 0; 108a7e14dcfSSatish Balay nlsP->ksp_ctol = 0; 109a7e14dcfSSatish Balay nlsP->ksp_negc = 0; 110a7e14dcfSSatish Balay nlsP->ksp_iter = 0; 111a7e14dcfSSatish Balay nlsP->ksp_othr = 0; 112a7e14dcfSSatish Balay 113a7e14dcfSSatish Balay /* Modify the linear solver to a trust region method if desired */ 114a7e14dcfSSatish Balay switch(nlsP->ksp_type) { 115a7e14dcfSSatish Balay case NLS_KSP_CG: 116a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPCG);CHKERRQ(ierr); 11772b7fd4bSBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 118a7e14dcfSSatish Balay break; 119a7e14dcfSSatish Balay 120a7e14dcfSSatish Balay case NLS_KSP_NASH: 121a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPNASH);CHKERRQ(ierr); 12272b7fd4bSBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 123a7e14dcfSSatish Balay break; 124a7e14dcfSSatish Balay 125a7e14dcfSSatish Balay case NLS_KSP_STCG: 126a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPSTCG);CHKERRQ(ierr); 12772b7fd4bSBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 128a7e14dcfSSatish Balay break; 129a7e14dcfSSatish Balay 130a7e14dcfSSatish Balay case NLS_KSP_GLTR: 131a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPGLTR);CHKERRQ(ierr); 13272b7fd4bSBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 133a7e14dcfSSatish Balay break; 134a7e14dcfSSatish Balay 135a7e14dcfSSatish Balay default: 136a7e14dcfSSatish Balay /* Use the method set by the ksp_type */ 137a7e14dcfSSatish Balay break; 138a7e14dcfSSatish Balay } 139a7e14dcfSSatish Balay 140a7e14dcfSSatish Balay /* Initialize trust-region radius when using nash, stcg, or gltr 141a7e14dcfSSatish Balay Will be reset during the first iteration */ 142a7e14dcfSSatish Balay if (NLS_KSP_NASH == nlsP->ksp_type) { 143a7e14dcfSSatish Balay ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 144a7e14dcfSSatish Balay } else if (NLS_KSP_STCG == nlsP->ksp_type) { 145a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 146a7e14dcfSSatish Balay } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 147a7e14dcfSSatish Balay ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 148a7e14dcfSSatish Balay } 149a7e14dcfSSatish Balay 15087f595a5SBarry Smith if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) { 151a7e14dcfSSatish Balay tao->trust = tao->trust0; 152a7e14dcfSSatish Balay 15387f595a5SBarry Smith if (tao->trust < 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial radius negative"); 154a7e14dcfSSatish Balay 155a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 156a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 157a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 158a7e14dcfSSatish Balay } 159a7e14dcfSSatish Balay 160a7e14dcfSSatish Balay /* Get vectors we will need */ 161a7e14dcfSSatish Balay if (NLS_PC_BFGS == nlsP->pc_type && !nlsP->M) { 162a7e14dcfSSatish Balay ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 163a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 164a7e14dcfSSatish Balay ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&nlsP->M);CHKERRQ(ierr); 165a7e14dcfSSatish Balay ierr = MatLMVMAllocateVectors(nlsP->M,tao->solution);CHKERRQ(ierr); 166a7e14dcfSSatish Balay } 167a7e14dcfSSatish Balay 168a7e14dcfSSatish Balay /* Check convergence criteria */ 169a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 170a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 17187f595a5SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 172a7e14dcfSSatish Balay needH = 1; 173a7e14dcfSSatish Balay 174a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 17587f595a5SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 176a7e14dcfSSatish Balay 177a7e14dcfSSatish Balay /* create vectors for the limited memory preconditioner */ 17887f595a5SBarry Smith if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_BFGS != nlsP->bfgs_scale_type)) { 179a7e14dcfSSatish Balay if (!nlsP->Diag) { 180a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&nlsP->Diag);CHKERRQ(ierr); 181a7e14dcfSSatish Balay } 182a7e14dcfSSatish Balay } 183a7e14dcfSSatish Balay 184a7e14dcfSSatish Balay /* Modify the preconditioner to use the bfgs approximation */ 185a7e14dcfSSatish Balay ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 186a7e14dcfSSatish Balay switch(nlsP->pc_type) { 187a7e14dcfSSatish Balay case NLS_PC_NONE: 188a7e14dcfSSatish Balay ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr); 189a7e14dcfSSatish Balay if (pc->ops->setfromoptions) { 190a7e14dcfSSatish Balay (*pc->ops->setfromoptions)(pc); 191a7e14dcfSSatish Balay } 192a7e14dcfSSatish Balay break; 193a7e14dcfSSatish Balay 194a7e14dcfSSatish Balay case NLS_PC_AHESS: 195a7e14dcfSSatish Balay ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr); 196a7e14dcfSSatish Balay if (pc->ops->setfromoptions) { 197a7e14dcfSSatish Balay (*pc->ops->setfromoptions)(pc); 198a7e14dcfSSatish Balay } 199a7e14dcfSSatish Balay ierr = PCJacobiSetUseAbs(pc);CHKERRQ(ierr); 200a7e14dcfSSatish Balay break; 201a7e14dcfSSatish Balay 202a7e14dcfSSatish Balay case NLS_PC_BFGS: 203a7e14dcfSSatish Balay ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr); 204a7e14dcfSSatish Balay if (pc->ops->setfromoptions) { 205a7e14dcfSSatish Balay (*pc->ops->setfromoptions)(pc); 206a7e14dcfSSatish Balay } 207a7e14dcfSSatish Balay ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr); 208a7e14dcfSSatish Balay ierr = PCShellSetContext(pc, nlsP->M);CHKERRQ(ierr); 209a7e14dcfSSatish Balay ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr); 210a7e14dcfSSatish Balay break; 211a7e14dcfSSatish Balay 212a7e14dcfSSatish Balay default: 213a7e14dcfSSatish Balay /* Use the pc method set by pc_type */ 214a7e14dcfSSatish Balay break; 215a7e14dcfSSatish Balay } 216a7e14dcfSSatish Balay 217a7e14dcfSSatish Balay /* Initialize trust-region radius. The initialization is only performed 218a7e14dcfSSatish Balay when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */ 21987f595a5SBarry Smith if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) { 220a7e14dcfSSatish Balay switch(nlsP->init_type) { 221a7e14dcfSSatish Balay case NLS_INIT_CONSTANT: 222a7e14dcfSSatish Balay /* Use the initial radius specified */ 223a7e14dcfSSatish Balay break; 224a7e14dcfSSatish Balay 225a7e14dcfSSatish Balay case NLS_INIT_INTERPOLATION: 226a7e14dcfSSatish Balay /* Use the initial radius specified */ 227a7e14dcfSSatish Balay max_radius = 0.0; 228a7e14dcfSSatish Balay 229a7e14dcfSSatish Balay for (j = 0; j < j_max; ++j) { 230a7e14dcfSSatish Balay fmin = f; 231a7e14dcfSSatish Balay sigma = 0.0; 232a7e14dcfSSatish Balay 233a7e14dcfSSatish Balay if (needH) { 234*ffad9901SBarry Smith ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 235a7e14dcfSSatish Balay needH = 0; 236a7e14dcfSSatish Balay } 237a7e14dcfSSatish Balay 238a7e14dcfSSatish Balay for (i = 0; i < i_max; ++i) { 239a7e14dcfSSatish Balay ierr = VecCopy(tao->solution,nlsP->W);CHKERRQ(ierr); 240a7e14dcfSSatish Balay ierr = VecAXPY(nlsP->W,-tao->trust/gnorm,tao->gradient);CHKERRQ(ierr); 241a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, nlsP->W, &ftrial);CHKERRQ(ierr); 242a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 243a7e14dcfSSatish Balay tau = nlsP->gamma1_i; 24487f595a5SBarry Smith } else { 245a7e14dcfSSatish Balay if (ftrial < fmin) { 246a7e14dcfSSatish Balay fmin = ftrial; 247a7e14dcfSSatish Balay sigma = -tao->trust / gnorm; 248a7e14dcfSSatish Balay } 249a7e14dcfSSatish Balay 250a7e14dcfSSatish Balay ierr = MatMult(tao->hessian, tao->gradient, nlsP->D);CHKERRQ(ierr); 251a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, nlsP->D, &prered);CHKERRQ(ierr); 252a7e14dcfSSatish Balay 253a7e14dcfSSatish Balay prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm)); 254a7e14dcfSSatish Balay actred = f - ftrial; 25587f595a5SBarry Smith if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) { 256a7e14dcfSSatish Balay kappa = 1.0; 25787f595a5SBarry Smith } else { 258a7e14dcfSSatish Balay kappa = actred / prered; 259a7e14dcfSSatish Balay } 260a7e14dcfSSatish Balay 261a7e14dcfSSatish Balay tau_1 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred); 262a7e14dcfSSatish Balay tau_2 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred); 263a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 264a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 265a7e14dcfSSatish Balay 266a7e14dcfSSatish Balay if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu1_i) { 267a7e14dcfSSatish Balay /* Great agreement */ 268a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 269a7e14dcfSSatish Balay 270a7e14dcfSSatish Balay if (tau_max < 1.0) { 271a7e14dcfSSatish Balay tau = nlsP->gamma3_i; 27287f595a5SBarry Smith } else if (tau_max > nlsP->gamma4_i) { 273a7e14dcfSSatish Balay tau = nlsP->gamma4_i; 27487f595a5SBarry Smith } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) { 275a7e14dcfSSatish Balay tau = tau_1; 27687f595a5SBarry Smith } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) { 277a7e14dcfSSatish Balay tau = tau_2; 27887f595a5SBarry Smith } else { 279a7e14dcfSSatish Balay tau = tau_max; 280a7e14dcfSSatish Balay } 28187f595a5SBarry Smith } else if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu2_i) { 282a7e14dcfSSatish Balay /* Good agreement */ 283a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 284a7e14dcfSSatish Balay 285a7e14dcfSSatish Balay if (tau_max < nlsP->gamma2_i) { 286a7e14dcfSSatish Balay tau = nlsP->gamma2_i; 28787f595a5SBarry Smith } else if (tau_max > nlsP->gamma3_i) { 288a7e14dcfSSatish Balay tau = nlsP->gamma3_i; 28987f595a5SBarry Smith } else { 290a7e14dcfSSatish Balay tau = tau_max; 291a7e14dcfSSatish Balay } 29287f595a5SBarry Smith } else { 293a7e14dcfSSatish Balay /* Not good agreement */ 294a7e14dcfSSatish Balay if (tau_min > 1.0) { 295a7e14dcfSSatish Balay tau = nlsP->gamma2_i; 29687f595a5SBarry Smith } else if (tau_max < nlsP->gamma1_i) { 297a7e14dcfSSatish Balay tau = nlsP->gamma1_i; 29887f595a5SBarry Smith } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) { 299a7e14dcfSSatish Balay tau = nlsP->gamma1_i; 30087f595a5SBarry Smith } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) { 301a7e14dcfSSatish Balay tau = tau_1; 30287f595a5SBarry Smith } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) { 303a7e14dcfSSatish Balay tau = tau_2; 30487f595a5SBarry Smith } else { 305a7e14dcfSSatish Balay tau = tau_max; 306a7e14dcfSSatish Balay } 307a7e14dcfSSatish Balay } 308a7e14dcfSSatish Balay } 309a7e14dcfSSatish Balay tao->trust = tau * tao->trust; 310a7e14dcfSSatish Balay } 311a7e14dcfSSatish Balay 312a7e14dcfSSatish Balay if (fmin < f) { 313a7e14dcfSSatish Balay f = fmin; 314a7e14dcfSSatish Balay ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr); 315a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao,tao->solution,tao->gradient);CHKERRQ(ierr); 316a7e14dcfSSatish Balay 317a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 31887f595a5SBarry Smith if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN"); 319a7e14dcfSSatish Balay needH = 1; 320a7e14dcfSSatish Balay 321a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 32287f595a5SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 323a7e14dcfSSatish Balay } 324a7e14dcfSSatish Balay } 325a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, max_radius); 326a7e14dcfSSatish Balay 327a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 328a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 329a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 330a7e14dcfSSatish Balay break; 331a7e14dcfSSatish Balay 332a7e14dcfSSatish Balay default: 333a7e14dcfSSatish Balay /* Norm of the first direction will initialize radius */ 334a7e14dcfSSatish Balay tao->trust = 0.0; 335a7e14dcfSSatish Balay break; 336a7e14dcfSSatish Balay } 337a7e14dcfSSatish Balay } 338a7e14dcfSSatish Balay 339a7e14dcfSSatish Balay /* Set initial scaling for the BFGS preconditioner 340a7e14dcfSSatish Balay This step is done after computing the initial trust-region radius 341a7e14dcfSSatish Balay since the function value may have decreased */ 342a7e14dcfSSatish Balay if (NLS_PC_BFGS == nlsP->pc_type) { 343a7e14dcfSSatish Balay if (f != 0.0) { 344a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 34587f595a5SBarry Smith } else { 346a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 347a7e14dcfSSatish Balay } 348a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr); 349a7e14dcfSSatish Balay } 350a7e14dcfSSatish Balay 351a7e14dcfSSatish Balay /* Set counter for gradient/reset steps*/ 352a7e14dcfSSatish Balay nlsP->newt = 0; 353a7e14dcfSSatish Balay nlsP->bfgs = 0; 354a7e14dcfSSatish Balay nlsP->sgrad = 0; 355a7e14dcfSSatish Balay nlsP->grad = 0; 356a7e14dcfSSatish Balay 357a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 358a7e14dcfSSatish Balay while (reason == TAO_CONTINUE_ITERATING) { 359a7e14dcfSSatish Balay ++iter; 360a7e14dcfSSatish Balay 361a7e14dcfSSatish Balay /* Compute the Hessian */ 362a7e14dcfSSatish Balay if (needH) { 363*ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 364a7e14dcfSSatish Balay needH = 0; 365a7e14dcfSSatish Balay } 366a7e14dcfSSatish Balay 36787f595a5SBarry Smith if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_AHESS == nlsP->bfgs_scale_type)) { 368a7e14dcfSSatish Balay /* Obtain diagonal for the bfgs preconditioner */ 369a7e14dcfSSatish Balay ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr); 370a7e14dcfSSatish Balay ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr); 371a7e14dcfSSatish Balay ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr); 372a7e14dcfSSatish Balay ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr); 373a7e14dcfSSatish Balay } 374a7e14dcfSSatish Balay 375a7e14dcfSSatish Balay /* Shift the Hessian matrix */ 376a7e14dcfSSatish Balay if (pert > 0) { 377a7e14dcfSSatish Balay ierr = MatShift(tao->hessian, pert); 378a7e14dcfSSatish Balay if (tao->hessian != tao->hessian_pre) { 379a7e14dcfSSatish Balay ierr = MatShift(tao->hessian_pre, pert);CHKERRQ(ierr); 380a7e14dcfSSatish Balay } 381a7e14dcfSSatish Balay } 382a7e14dcfSSatish Balay 383a7e14dcfSSatish Balay if (NLS_PC_BFGS == nlsP->pc_type) { 384a7e14dcfSSatish Balay if (BFGS_SCALE_PHESS == nlsP->bfgs_scale_type) { 385a7e14dcfSSatish Balay /* Obtain diagonal for the bfgs preconditioner */ 386a7e14dcfSSatish Balay ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr); 387a7e14dcfSSatish Balay ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr); 388a7e14dcfSSatish Balay ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr); 389a7e14dcfSSatish Balay ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr); 390a7e14dcfSSatish Balay } 391a7e14dcfSSatish Balay /* Update the limited memory preconditioner */ 392a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 393a7e14dcfSSatish Balay ++bfgsUpdates; 394a7e14dcfSSatish Balay } 395a7e14dcfSSatish Balay 396a7e14dcfSSatish Balay /* Solve the Newton system of equations */ 39723ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 39887f595a5SBarry Smith if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) { 399a7e14dcfSSatish Balay 400a7e14dcfSSatish Balay if (NLS_KSP_NASH == nlsP->ksp_type) { 401a7e14dcfSSatish Balay ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 402a7e14dcfSSatish Balay } else if (NLS_KSP_STCG == nlsP->ksp_type) { 403a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 404a7e14dcfSSatish Balay } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 405a7e14dcfSSatish Balay ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 406a7e14dcfSSatish Balay } 407a7e14dcfSSatish Balay 408a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr); 409a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 410a7e14dcfSSatish Balay tao->ksp_its+=kspits; 411a7e14dcfSSatish Balay 412a7e14dcfSSatish Balay if (NLS_KSP_NASH == nlsP->ksp_type) { 413a7e14dcfSSatish Balay ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 414a7e14dcfSSatish Balay } else if (NLS_KSP_STCG == nlsP->ksp_type) { 415a7e14dcfSSatish Balay ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 416a7e14dcfSSatish Balay } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 417a7e14dcfSSatish Balay ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 418a7e14dcfSSatish Balay } 419a7e14dcfSSatish Balay 420a7e14dcfSSatish Balay if (0.0 == tao->trust) { 421a7e14dcfSSatish Balay /* Radius was uninitialized; use the norm of the direction */ 422a7e14dcfSSatish Balay if (norm_d > 0.0) { 423a7e14dcfSSatish Balay tao->trust = norm_d; 424a7e14dcfSSatish Balay 425a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 426a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 427a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 42887f595a5SBarry Smith } else { 429a7e14dcfSSatish Balay /* The direction was bad; set radius to default value and re-solve 430a7e14dcfSSatish Balay the trust-region subproblem to get a direction */ 431a7e14dcfSSatish Balay tao->trust = tao->trust0; 432a7e14dcfSSatish Balay 433a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 434a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 435a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 436a7e14dcfSSatish Balay 437a7e14dcfSSatish Balay if (NLS_KSP_NASH == nlsP->ksp_type) { 438a7e14dcfSSatish Balay ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 439a7e14dcfSSatish Balay } else if (NLS_KSP_STCG == nlsP->ksp_type) { 440a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 441a7e14dcfSSatish Balay } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 442a7e14dcfSSatish Balay ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 443a7e14dcfSSatish Balay } 444a7e14dcfSSatish Balay 445a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr); 446a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 447a7e14dcfSSatish Balay tao->ksp_its+=kspits; 448a7e14dcfSSatish Balay if (NLS_KSP_NASH == nlsP->ksp_type) { 449a7e14dcfSSatish Balay ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 450a7e14dcfSSatish Balay } else if (NLS_KSP_STCG == nlsP->ksp_type) { 451a7e14dcfSSatish Balay ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 452a7e14dcfSSatish Balay } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 453a7e14dcfSSatish Balay ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 454a7e14dcfSSatish Balay } 455a7e14dcfSSatish Balay 45687f595a5SBarry Smith if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 457a7e14dcfSSatish Balay } 458a7e14dcfSSatish Balay } 45987f595a5SBarry Smith } else { 460a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr); 461a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr); 462a7e14dcfSSatish Balay tao->ksp_its += kspits; 463a7e14dcfSSatish Balay } 464a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 465a7e14dcfSSatish Balay ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 46687f595a5SBarry Smith if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (NLS_PC_BFGS == nlsP->pc_type) && (bfgsUpdates > 1)) { 467a7e14dcfSSatish Balay /* Preconditioner is numerically indefinite; reset the 468a7e14dcfSSatish Balay approximate if using BFGS preconditioning. */ 469a7e14dcfSSatish Balay 470a7e14dcfSSatish Balay if (f != 0.0) { 471a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 47287f595a5SBarry Smith } else { 473a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 474a7e14dcfSSatish Balay } 475a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr); 476a7e14dcfSSatish Balay ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 477a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 478a7e14dcfSSatish Balay bfgsUpdates = 1; 479a7e14dcfSSatish Balay } 480a7e14dcfSSatish Balay 481a7e14dcfSSatish Balay if (KSP_CONVERGED_ATOL == ksp_reason) { 482a7e14dcfSSatish Balay ++nlsP->ksp_atol; 48387f595a5SBarry Smith } else if (KSP_CONVERGED_RTOL == ksp_reason) { 484a7e14dcfSSatish Balay ++nlsP->ksp_rtol; 48587f595a5SBarry Smith } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) { 486a7e14dcfSSatish Balay ++nlsP->ksp_ctol; 48787f595a5SBarry Smith } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) { 488a7e14dcfSSatish Balay ++nlsP->ksp_negc; 48987f595a5SBarry Smith } else if (KSP_DIVERGED_DTOL == ksp_reason) { 490a7e14dcfSSatish Balay ++nlsP->ksp_dtol; 49187f595a5SBarry Smith } else if (KSP_DIVERGED_ITS == ksp_reason) { 492a7e14dcfSSatish Balay ++nlsP->ksp_iter; 49387f595a5SBarry Smith } else { 494a7e14dcfSSatish Balay ++nlsP->ksp_othr; 495a7e14dcfSSatish Balay } 496a7e14dcfSSatish Balay 497a7e14dcfSSatish Balay /* Check for success (descent direction) */ 498a7e14dcfSSatish Balay ierr = VecDot(nlsP->D, tao->gradient, &gdx);CHKERRQ(ierr); 499a7e14dcfSSatish Balay if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 500a7e14dcfSSatish Balay /* Newton step is not descent or direction produced Inf or NaN 501a7e14dcfSSatish Balay Update the perturbation for next time */ 502a7e14dcfSSatish Balay if (pert <= 0.0) { 503a7e14dcfSSatish Balay /* Initialize the perturbation */ 504a7e14dcfSSatish Balay pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 505a7e14dcfSSatish Balay if (NLS_KSP_GLTR == nlsP->ksp_type) { 506a7e14dcfSSatish Balay ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 507a7e14dcfSSatish Balay pert = PetscMax(pert, -e_min); 508a7e14dcfSSatish Balay } 50987f595a5SBarry Smith } else { 510a7e14dcfSSatish Balay /* Increase the perturbation */ 511a7e14dcfSSatish Balay pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 512a7e14dcfSSatish Balay } 513a7e14dcfSSatish Balay 514a7e14dcfSSatish Balay if (NLS_PC_BFGS != nlsP->pc_type) { 515a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and updated 516a7e14dcfSSatish Balay Must use gradient direction in this case */ 517a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr); 518a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 519a7e14dcfSSatish Balay ++nlsP->grad; 520a7e14dcfSSatish Balay stepType = NLS_GRADIENT; 52187f595a5SBarry Smith } else { 522a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 523a7e14dcfSSatish Balay ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 524a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 525a7e14dcfSSatish Balay 526a7e14dcfSSatish Balay /* Check for success (descent direction) */ 527a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, nlsP->D, &gdx);CHKERRQ(ierr); 528a7e14dcfSSatish Balay if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) { 529a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced not a number 530a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case because 531a7e14dcfSSatish Balay the first solve produces the scaled gradient direction, 532a7e14dcfSSatish Balay which is guaranteed to be descent */ 533a7e14dcfSSatish Balay 534a7e14dcfSSatish Balay /* Use steepest descent direction (scaled) */ 535a7e14dcfSSatish Balay 536a7e14dcfSSatish Balay if (f != 0.0) { 537a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 53887f595a5SBarry Smith } else { 539a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 540a7e14dcfSSatish Balay } 541a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr); 542a7e14dcfSSatish Balay ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 543a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 544a7e14dcfSSatish Balay ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 545a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 546a7e14dcfSSatish Balay 547a7e14dcfSSatish Balay bfgsUpdates = 1; 548a7e14dcfSSatish Balay ++nlsP->sgrad; 549a7e14dcfSSatish Balay stepType = NLS_SCALED_GRADIENT; 55087f595a5SBarry Smith } else { 551a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 552a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 553a7e14dcfSSatish Balay ++nlsP->sgrad; 554a7e14dcfSSatish Balay stepType = NLS_SCALED_GRADIENT; 55587f595a5SBarry Smith } else { 556a7e14dcfSSatish Balay ++nlsP->bfgs; 557a7e14dcfSSatish Balay stepType = NLS_BFGS; 558a7e14dcfSSatish Balay } 559a7e14dcfSSatish Balay } 560a7e14dcfSSatish Balay } 56187f595a5SBarry Smith } else { 562a7e14dcfSSatish Balay /* Computed Newton step is descent */ 563a7e14dcfSSatish Balay switch (ksp_reason) { 564a7e14dcfSSatish Balay case KSP_DIVERGED_NANORINF: 565a7e14dcfSSatish Balay case KSP_DIVERGED_BREAKDOWN: 566a7e14dcfSSatish Balay case KSP_DIVERGED_INDEFINITE_MAT: 567a7e14dcfSSatish Balay case KSP_DIVERGED_INDEFINITE_PC: 568a7e14dcfSSatish Balay case KSP_CONVERGED_CG_NEG_CURVE: 569a7e14dcfSSatish Balay /* Matrix or preconditioner is indefinite; increase perturbation */ 570a7e14dcfSSatish Balay if (pert <= 0.0) { 571a7e14dcfSSatish Balay /* Initialize the perturbation */ 572a7e14dcfSSatish Balay pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 573a7e14dcfSSatish Balay if (NLS_KSP_GLTR == nlsP->ksp_type) { 574a7e14dcfSSatish Balay ierr = KSPGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr); 575a7e14dcfSSatish Balay pert = PetscMax(pert, -e_min); 576a7e14dcfSSatish Balay } 57787f595a5SBarry Smith } else { 578a7e14dcfSSatish Balay /* Increase the perturbation */ 579a7e14dcfSSatish Balay pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 580a7e14dcfSSatish Balay } 581a7e14dcfSSatish Balay break; 582a7e14dcfSSatish Balay 583a7e14dcfSSatish Balay default: 584a7e14dcfSSatish Balay /* Newton step computation is good; decrease perturbation */ 585a7e14dcfSSatish Balay pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm); 586a7e14dcfSSatish Balay if (pert < nlsP->pmin) { 587a7e14dcfSSatish Balay pert = 0.0; 588a7e14dcfSSatish Balay } 589a7e14dcfSSatish Balay break; 590a7e14dcfSSatish Balay } 591a7e14dcfSSatish Balay 592a7e14dcfSSatish Balay ++nlsP->newt; 593a7e14dcfSSatish Balay stepType = NLS_NEWTON; 594a7e14dcfSSatish Balay } 595a7e14dcfSSatish Balay 596a7e14dcfSSatish Balay /* Perform the linesearch */ 597a7e14dcfSSatish Balay fold = f; 598a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, nlsP->Xold);CHKERRQ(ierr); 599a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, nlsP->Gold);CHKERRQ(ierr); 600a7e14dcfSSatish Balay 601a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr); 602a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 603a7e14dcfSSatish Balay 60487f595a5SBarry Smith while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) { /* Linesearch failed */ 605a7e14dcfSSatish Balay f = fold; 606a7e14dcfSSatish Balay ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr); 607a7e14dcfSSatish Balay ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr); 608a7e14dcfSSatish Balay 609a7e14dcfSSatish Balay switch(stepType) { 610a7e14dcfSSatish Balay case NLS_NEWTON: 611a7e14dcfSSatish Balay /* Failed to obtain acceptable iterate with Newton 1step 612a7e14dcfSSatish Balay Update the perturbation for next time */ 613a7e14dcfSSatish Balay if (pert <= 0.0) { 614a7e14dcfSSatish Balay /* Initialize the perturbation */ 615a7e14dcfSSatish Balay pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 616a7e14dcfSSatish Balay if (NLS_KSP_GLTR == nlsP->ksp_type) { 617a7e14dcfSSatish Balay ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 618a7e14dcfSSatish Balay pert = PetscMax(pert, -e_min); 619a7e14dcfSSatish Balay } 62087f595a5SBarry Smith } else { 621a7e14dcfSSatish Balay /* Increase the perturbation */ 622a7e14dcfSSatish Balay pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 623a7e14dcfSSatish Balay } 624a7e14dcfSSatish Balay 625a7e14dcfSSatish Balay if (NLS_PC_BFGS != nlsP->pc_type) { 626a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and being updated 627a7e14dcfSSatish Balay Must use gradient direction in this case */ 628a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr); 629a7e14dcfSSatish Balay ++nlsP->grad; 630a7e14dcfSSatish Balay stepType = NLS_GRADIENT; 63187f595a5SBarry Smith } else { 632a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 633a7e14dcfSSatish Balay ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 634a7e14dcfSSatish Balay /* Check for success (descent direction) */ 635a7e14dcfSSatish Balay ierr = VecDot(tao->solution, nlsP->D, &gdx);CHKERRQ(ierr); 636a7e14dcfSSatish Balay if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 637a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced not a number 638a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case 639a7e14dcfSSatish Balay Use steepest descent direction (scaled) */ 640a7e14dcfSSatish Balay 641a7e14dcfSSatish Balay if (f != 0.0) { 642a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 64387f595a5SBarry Smith } else { 644a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 645a7e14dcfSSatish Balay } 646a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr); 647a7e14dcfSSatish Balay ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 648a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 649a7e14dcfSSatish Balay ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 650a7e14dcfSSatish Balay 651a7e14dcfSSatish Balay bfgsUpdates = 1; 652a7e14dcfSSatish Balay ++nlsP->sgrad; 653a7e14dcfSSatish Balay stepType = NLS_SCALED_GRADIENT; 65487f595a5SBarry Smith } else { 655a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 656a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 657a7e14dcfSSatish Balay ++nlsP->sgrad; 658a7e14dcfSSatish Balay stepType = NLS_SCALED_GRADIENT; 65987f595a5SBarry Smith } else { 660a7e14dcfSSatish Balay ++nlsP->bfgs; 661a7e14dcfSSatish Balay stepType = NLS_BFGS; 662a7e14dcfSSatish Balay } 663a7e14dcfSSatish Balay } 664a7e14dcfSSatish Balay } 665a7e14dcfSSatish Balay break; 666a7e14dcfSSatish Balay 667a7e14dcfSSatish Balay case NLS_BFGS: 668a7e14dcfSSatish Balay /* Can only enter if pc_type == NLS_PC_BFGS 669a7e14dcfSSatish Balay Failed to obtain acceptable iterate with BFGS step 670a7e14dcfSSatish Balay Attempt to use the scaled gradient direction */ 671a7e14dcfSSatish Balay 672a7e14dcfSSatish Balay if (f != 0.0) { 673a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 67487f595a5SBarry Smith } else { 675a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 676a7e14dcfSSatish Balay } 677a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr); 678a7e14dcfSSatish Balay ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 679a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 680a7e14dcfSSatish Balay ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 681a7e14dcfSSatish Balay 682a7e14dcfSSatish Balay bfgsUpdates = 1; 683a7e14dcfSSatish Balay ++nlsP->sgrad; 684a7e14dcfSSatish Balay stepType = NLS_SCALED_GRADIENT; 685a7e14dcfSSatish Balay break; 686a7e14dcfSSatish Balay 687a7e14dcfSSatish Balay case NLS_SCALED_GRADIENT: 688a7e14dcfSSatish Balay /* Can only enter if pc_type == NLS_PC_BFGS 689a7e14dcfSSatish Balay The scaled gradient step did not produce a new iterate; 690a7e14dcfSSatish Balay attemp to use the gradient direction. 691a7e14dcfSSatish Balay Need to make sure we are not using a different diagonal scaling */ 692a7e14dcfSSatish Balay 693a7e14dcfSSatish Balay ierr = MatLMVMSetScale(nlsP->M,0);CHKERRQ(ierr); 694a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(nlsP->M,1.0);CHKERRQ(ierr); 695a7e14dcfSSatish Balay ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 696a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 697a7e14dcfSSatish Balay ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 698a7e14dcfSSatish Balay 699a7e14dcfSSatish Balay bfgsUpdates = 1; 700a7e14dcfSSatish Balay ++nlsP->grad; 701a7e14dcfSSatish Balay stepType = NLS_GRADIENT; 702a7e14dcfSSatish Balay break; 703a7e14dcfSSatish Balay } 704a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 705a7e14dcfSSatish Balay 706a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr); 707a7e14dcfSSatish Balay ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr); 708a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 709a7e14dcfSSatish Balay } 710a7e14dcfSSatish Balay 71187f595a5SBarry Smith if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 712a7e14dcfSSatish Balay /* Failed to find an improving point */ 713a7e14dcfSSatish Balay f = fold; 714a7e14dcfSSatish Balay ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr); 715a7e14dcfSSatish Balay ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr); 716a7e14dcfSSatish Balay step = 0.0; 717a7e14dcfSSatish Balay reason = TAO_DIVERGED_LS_FAILURE; 718a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 719a7e14dcfSSatish Balay break; 720a7e14dcfSSatish Balay } 721a7e14dcfSSatish Balay 722a7e14dcfSSatish Balay /* Update trust region radius */ 72387f595a5SBarry Smith if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) { 724a7e14dcfSSatish Balay switch(nlsP->update_type) { 725a7e14dcfSSatish Balay case NLS_UPDATE_STEP: 726a7e14dcfSSatish Balay if (stepType == NLS_NEWTON) { 727a7e14dcfSSatish Balay if (step < nlsP->nu1) { 728a7e14dcfSSatish Balay /* Very bad step taken; reduce radius */ 729a7e14dcfSSatish Balay tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust); 73087f595a5SBarry Smith } else if (step < nlsP->nu2) { 731a7e14dcfSSatish Balay /* Reasonably bad step taken; reduce radius */ 732a7e14dcfSSatish Balay tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust); 73387f595a5SBarry Smith } else if (step < nlsP->nu3) { 734a7e14dcfSSatish Balay /* Reasonable step was taken; leave radius alone */ 735a7e14dcfSSatish Balay if (nlsP->omega3 < 1.0) { 736a7e14dcfSSatish Balay tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust); 73787f595a5SBarry Smith } else if (nlsP->omega3 > 1.0) { 738a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust); 739a7e14dcfSSatish Balay } 74087f595a5SBarry Smith } else if (step < nlsP->nu4) { 741a7e14dcfSSatish Balay /* Full step taken; increase the radius */ 742a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust); 74387f595a5SBarry Smith } else { 744a7e14dcfSSatish Balay /* More than full step taken; increase the radius */ 745a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust); 746a7e14dcfSSatish Balay } 74787f595a5SBarry Smith } else { 748a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 749a7e14dcfSSatish Balay tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust); 750a7e14dcfSSatish Balay } 751a7e14dcfSSatish Balay break; 752a7e14dcfSSatish Balay 753a7e14dcfSSatish Balay case NLS_UPDATE_REDUCTION: 754a7e14dcfSSatish Balay if (stepType == NLS_NEWTON) { 755a7e14dcfSSatish Balay /* Get predicted reduction */ 756a7e14dcfSSatish Balay 757a7e14dcfSSatish Balay if (NLS_KSP_STCG == nlsP->ksp_type) { 758a7e14dcfSSatish Balay ierr = KSPSTCGGetObjFcn(tao->ksp,&prered); 759a7e14dcfSSatish Balay } else if (NLS_KSP_NASH == nlsP->ksp_type) { 760a7e14dcfSSatish Balay ierr = KSPNASHGetObjFcn(tao->ksp,&prered); 761a7e14dcfSSatish Balay } else { 762a7e14dcfSSatish Balay ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 763a7e14dcfSSatish Balay } 764a7e14dcfSSatish Balay 765a7e14dcfSSatish Balay if (prered >= 0.0) { 766a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot */ 767a7e14dcfSSatish Balay /* happen in infinite precision arithmetic. Step should */ 768a7e14dcfSSatish Balay /* be rejected! */ 769a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 77087f595a5SBarry Smith } else { 771a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(f_full)) { 772a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 77387f595a5SBarry Smith } else { 774a7e14dcfSSatish Balay /* Compute and actual reduction */ 775a7e14dcfSSatish Balay actred = fold - f_full; 776a7e14dcfSSatish Balay prered = -prered; 777a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= nlsP->epsilon) && 778a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= nlsP->epsilon)) { 779a7e14dcfSSatish Balay kappa = 1.0; 78087f595a5SBarry Smith } else { 781a7e14dcfSSatish Balay kappa = actred / prered; 782a7e14dcfSSatish Balay } 783a7e14dcfSSatish Balay 784a7e14dcfSSatish Balay /* Accept of reject the step and update radius */ 785a7e14dcfSSatish Balay if (kappa < nlsP->eta1) { 786a7e14dcfSSatish Balay /* Very bad step */ 787a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 78887f595a5SBarry Smith } else if (kappa < nlsP->eta2) { 789a7e14dcfSSatish Balay /* Marginal bad step */ 790a7e14dcfSSatish Balay tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d); 79187f595a5SBarry Smith } else if (kappa < nlsP->eta3) { 792a7e14dcfSSatish Balay /* Reasonable step */ 793a7e14dcfSSatish Balay if (nlsP->alpha3 < 1.0) { 794a7e14dcfSSatish Balay tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust); 79587f595a5SBarry Smith } else if (nlsP->alpha3 > 1.0) { 796a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust); 797a7e14dcfSSatish Balay } 79887f595a5SBarry Smith } else if (kappa < nlsP->eta4) { 799a7e14dcfSSatish Balay /* Good step */ 800a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust); 80187f595a5SBarry Smith } else { 802a7e14dcfSSatish Balay /* Very good step */ 803a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust); 804a7e14dcfSSatish Balay } 805a7e14dcfSSatish Balay } 806a7e14dcfSSatish Balay } 80787f595a5SBarry Smith } else { 808a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 809a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust); 810a7e14dcfSSatish Balay } 811a7e14dcfSSatish Balay break; 812a7e14dcfSSatish Balay 813a7e14dcfSSatish Balay default: 814a7e14dcfSSatish Balay if (stepType == NLS_NEWTON) { 815a7e14dcfSSatish Balay 816a7e14dcfSSatish Balay if (NLS_KSP_STCG == nlsP->ksp_type) { 817a7e14dcfSSatish Balay ierr = KSPSTCGGetObjFcn(tao->ksp,&prered); 818a7e14dcfSSatish Balay } else if (NLS_KSP_NASH == nlsP->ksp_type) { 819a7e14dcfSSatish Balay ierr = KSPNASHGetObjFcn(tao->ksp,&prered); 820a7e14dcfSSatish Balay } else { 821a7e14dcfSSatish Balay ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 822a7e14dcfSSatish Balay } 823a7e14dcfSSatish Balay if (prered >= 0.0) { 824a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot */ 825a7e14dcfSSatish Balay /* happen in infinite precision arithmetic. Step should */ 826a7e14dcfSSatish Balay /* be rejected! */ 827a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 82887f595a5SBarry Smith } else { 829a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(f_full)) { 830a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 83187f595a5SBarry Smith } else { 832a7e14dcfSSatish Balay actred = fold - f_full; 833a7e14dcfSSatish Balay prered = -prered; 83487f595a5SBarry Smith if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) { 835a7e14dcfSSatish Balay kappa = 1.0; 83687f595a5SBarry Smith } else { 837a7e14dcfSSatish Balay kappa = actred / prered; 838a7e14dcfSSatish Balay } 839a7e14dcfSSatish Balay 840a7e14dcfSSatish Balay tau_1 = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred); 841a7e14dcfSSatish Balay tau_2 = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred); 842a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 843a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 844a7e14dcfSSatish Balay 845a7e14dcfSSatish Balay if (kappa >= 1.0 - nlsP->mu1) { 846a7e14dcfSSatish Balay /* Great agreement */ 847a7e14dcfSSatish Balay if (tau_max < 1.0) { 848a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d); 84987f595a5SBarry Smith } else if (tau_max > nlsP->gamma4) { 850a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d); 85187f595a5SBarry Smith } else { 852a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 853a7e14dcfSSatish Balay } 85487f595a5SBarry Smith } else if (kappa >= 1.0 - nlsP->mu2) { 855a7e14dcfSSatish Balay /* Good agreement */ 856a7e14dcfSSatish Balay 857a7e14dcfSSatish Balay if (tau_max < nlsP->gamma2) { 858a7e14dcfSSatish Balay tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d); 85987f595a5SBarry Smith } else if (tau_max > nlsP->gamma3) { 860a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d); 86187f595a5SBarry Smith } else if (tau_max < 1.0) { 862a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 86387f595a5SBarry Smith } else { 864a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 865a7e14dcfSSatish Balay } 86687f595a5SBarry Smith } else { 867a7e14dcfSSatish Balay /* Not good agreement */ 868a7e14dcfSSatish Balay if (tau_min > 1.0) { 869a7e14dcfSSatish Balay tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d); 87087f595a5SBarry Smith } else if (tau_max < nlsP->gamma1) { 871a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 87287f595a5SBarry Smith } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) { 873a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 87487f595a5SBarry Smith } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) { 875a7e14dcfSSatish Balay tao->trust = tau_1 * PetscMin(tao->trust, norm_d); 87687f595a5SBarry Smith } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) { 877a7e14dcfSSatish Balay tao->trust = tau_2 * PetscMin(tao->trust, norm_d); 87887f595a5SBarry Smith } else { 879a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 880a7e14dcfSSatish Balay } 881a7e14dcfSSatish Balay } 882a7e14dcfSSatish Balay } 883a7e14dcfSSatish Balay } 88487f595a5SBarry Smith } else { 885a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 886a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust); 887a7e14dcfSSatish Balay } 888a7e14dcfSSatish Balay break; 889a7e14dcfSSatish Balay } 890a7e14dcfSSatish Balay 891a7e14dcfSSatish Balay /* The radius may have been increased; modify if it is too large */ 892a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 893a7e14dcfSSatish Balay } 894a7e14dcfSSatish Balay 895a7e14dcfSSatish Balay /* Check for termination */ 896a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 89787f595a5SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 898a7e14dcfSSatish Balay needH = 1; 899a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr); 900a7e14dcfSSatish Balay } 901a7e14dcfSSatish Balay PetscFunctionReturn(0); 902a7e14dcfSSatish Balay } 903a7e14dcfSSatish Balay 904a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 905a7e14dcfSSatish Balay #undef __FUNCT__ 906a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_NLS" 907441846f8SBarry Smith static PetscErrorCode TaoSetUp_NLS(Tao tao) 908a7e14dcfSSatish Balay { 909a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 910a7e14dcfSSatish Balay PetscErrorCode ierr; 911a7e14dcfSSatish Balay 912a7e14dcfSSatish Balay PetscFunctionBegin; 913a7e14dcfSSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);} 914a7e14dcfSSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);} 915a7e14dcfSSatish Balay if (!nlsP->W) {ierr = VecDuplicate(tao->solution,&nlsP->W);CHKERRQ(ierr);} 916a7e14dcfSSatish Balay if (!nlsP->D) {ierr = VecDuplicate(tao->solution,&nlsP->D);CHKERRQ(ierr);} 917a7e14dcfSSatish Balay if (!nlsP->Xold) {ierr = VecDuplicate(tao->solution,&nlsP->Xold);CHKERRQ(ierr);} 918a7e14dcfSSatish Balay if (!nlsP->Gold) {ierr = VecDuplicate(tao->solution,&nlsP->Gold);CHKERRQ(ierr);} 919a7e14dcfSSatish Balay nlsP->Diag = 0; 920a7e14dcfSSatish Balay nlsP->M = 0; 921a7e14dcfSSatish Balay PetscFunctionReturn(0); 922a7e14dcfSSatish Balay } 923a7e14dcfSSatish Balay 924a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 925a7e14dcfSSatish Balay #undef __FUNCT__ 926a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_NLS" 927441846f8SBarry Smith static PetscErrorCode TaoDestroy_NLS(Tao tao) 928a7e14dcfSSatish Balay { 929a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 930a7e14dcfSSatish Balay PetscErrorCode ierr; 931a7e14dcfSSatish Balay 932a7e14dcfSSatish Balay PetscFunctionBegin; 933a7e14dcfSSatish Balay if (tao->setupcalled) { 934a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->D);CHKERRQ(ierr); 935a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->W);CHKERRQ(ierr); 936a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->Xold);CHKERRQ(ierr); 937a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->Gold);CHKERRQ(ierr); 938a7e14dcfSSatish Balay } 939a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->Diag);CHKERRQ(ierr); 940a7e14dcfSSatish Balay ierr = MatDestroy(&nlsP->M);CHKERRQ(ierr); 941a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 942a7e14dcfSSatish Balay PetscFunctionReturn(0); 943a7e14dcfSSatish Balay } 944a7e14dcfSSatish Balay 945a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 946a7e14dcfSSatish Balay #undef __FUNCT__ 947a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_NLS" 948441846f8SBarry Smith static PetscErrorCode TaoSetFromOptions_NLS(Tao tao) 949a7e14dcfSSatish Balay { 950a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 951a7e14dcfSSatish Balay PetscErrorCode ierr; 952a7e14dcfSSatish Balay 953a7e14dcfSSatish Balay PetscFunctionBegin; 954a7e14dcfSSatish Balay ierr = PetscOptionsHead("Newton line search method for unconstrained optimization");CHKERRQ(ierr); 955a7e14dcfSSatish Balay ierr = PetscOptionsEList("-tao_nls_ksp_type", "ksp type", "", NLS_KSP, NLS_KSP_TYPES, NLS_KSP[nlsP->ksp_type], &nlsP->ksp_type, 0);CHKERRQ(ierr); 956a7e14dcfSSatish Balay ierr = PetscOptionsEList("-tao_nls_pc_type", "pc type", "", NLS_PC, NLS_PC_TYPES, NLS_PC[nlsP->pc_type], &nlsP->pc_type, 0);CHKERRQ(ierr); 957a7e14dcfSSatish Balay ierr = PetscOptionsEList("-tao_nls_bfgs_scale_type", "bfgs scale type", "", BFGS_SCALE, BFGS_SCALE_TYPES, BFGS_SCALE[nlsP->bfgs_scale_type], &nlsP->bfgs_scale_type, 0);CHKERRQ(ierr); 958a7e14dcfSSatish 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); 959a7e14dcfSSatish 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); 960a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval, 0);CHKERRQ(ierr); 961a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin, 0);CHKERRQ(ierr); 962a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax, 0);CHKERRQ(ierr); 963a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac, 0);CHKERRQ(ierr); 964a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin, 0);CHKERRQ(ierr); 965a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax, 0);CHKERRQ(ierr); 966a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac, 0);CHKERRQ(ierr); 967a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac, 0);CHKERRQ(ierr); 968a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac, 0);CHKERRQ(ierr); 969a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac, 0);CHKERRQ(ierr); 970a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1, 0);CHKERRQ(ierr); 971a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2, 0);CHKERRQ(ierr); 972a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3, 0);CHKERRQ(ierr); 973a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4, 0);CHKERRQ(ierr); 974a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1, 0);CHKERRQ(ierr); 975a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2, 0);CHKERRQ(ierr); 976a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3, 0);CHKERRQ(ierr); 977a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4, 0);CHKERRQ(ierr); 978a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5, 0);CHKERRQ(ierr); 979a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1, 0);CHKERRQ(ierr); 980a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2, 0);CHKERRQ(ierr); 981a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3, 0);CHKERRQ(ierr); 982a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4, 0);CHKERRQ(ierr); 983a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1, 0);CHKERRQ(ierr); 984a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2, 0);CHKERRQ(ierr); 985a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3, 0);CHKERRQ(ierr); 986a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4, 0);CHKERRQ(ierr); 987a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5, 0);CHKERRQ(ierr); 988a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i, 0);CHKERRQ(ierr); 989a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i, 0);CHKERRQ(ierr); 990a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i, 0);CHKERRQ(ierr); 991a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i, 0);CHKERRQ(ierr); 992a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i, 0);CHKERRQ(ierr); 993a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i, 0);CHKERRQ(ierr); 994a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i, 0);CHKERRQ(ierr); 995a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1, 0);CHKERRQ(ierr); 996a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2, 0);CHKERRQ(ierr); 997a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1, 0);CHKERRQ(ierr); 998a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2, 0);CHKERRQ(ierr); 999a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3, 0);CHKERRQ(ierr); 1000a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4, 0);CHKERRQ(ierr); 1001a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta, 0);CHKERRQ(ierr); 1002a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius, 0);CHKERRQ(ierr); 1003a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius, 0);CHKERRQ(ierr); 1004a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon, 0);CHKERRQ(ierr); 1005a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 1006a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 1007a7e14dcfSSatish Balay ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 1008a7e14dcfSSatish Balay PetscFunctionReturn(0); 1009a7e14dcfSSatish Balay } 1010a7e14dcfSSatish Balay 1011a7e14dcfSSatish Balay 1012a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 1013a7e14dcfSSatish Balay #undef __FUNCT__ 1014a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_NLS" 1015441846f8SBarry Smith static PetscErrorCode TaoView_NLS(Tao tao, PetscViewer viewer) 1016a7e14dcfSSatish Balay { 1017a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 1018a7e14dcfSSatish Balay PetscInt nrejects; 1019a7e14dcfSSatish Balay PetscBool isascii; 1020a7e14dcfSSatish Balay PetscErrorCode ierr; 1021a7e14dcfSSatish Balay 1022a7e14dcfSSatish Balay PetscFunctionBegin; 1023a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 1024a7e14dcfSSatish Balay if (isascii) { 1025a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1026a7e14dcfSSatish Balay if (NLS_PC_BFGS == nlsP->pc_type && nlsP->M) { 1027a7e14dcfSSatish Balay ierr = MatLMVMGetRejects(nlsP->M,&nrejects);CHKERRQ(ierr); 1028a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr); 1029a7e14dcfSSatish Balay } 1030a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", nlsP->newt);CHKERRQ(ierr); 1031a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", nlsP->bfgs);CHKERRQ(ierr); 1032a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", nlsP->sgrad);CHKERRQ(ierr); 1033a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", nlsP->grad);CHKERRQ(ierr); 1034a7e14dcfSSatish Balay 1035a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp atol: %D\n", nlsP->ksp_atol);CHKERRQ(ierr); 1036a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %D\n", nlsP->ksp_rtol);CHKERRQ(ierr); 1037a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %D\n", nlsP->ksp_ctol);CHKERRQ(ierr); 1038a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp negc: %D\n", nlsP->ksp_negc);CHKERRQ(ierr); 1039a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %D\n", nlsP->ksp_dtol);CHKERRQ(ierr); 1040a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp iter: %D\n", nlsP->ksp_iter);CHKERRQ(ierr); 1041a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp othr: %D\n", nlsP->ksp_othr);CHKERRQ(ierr); 1042a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1043a7e14dcfSSatish Balay } 1044a7e14dcfSSatish Balay PetscFunctionReturn(0); 1045a7e14dcfSSatish Balay } 1046a7e14dcfSSatish Balay 1047a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 1048a7e14dcfSSatish Balay EXTERN_C_BEGIN 1049a7e14dcfSSatish Balay #undef __FUNCT__ 1050a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_NLS" 1051441846f8SBarry Smith PetscErrorCode TaoCreate_NLS(Tao tao) 1052a7e14dcfSSatish Balay { 1053a7e14dcfSSatish Balay TAO_NLS *nlsP; 10548caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 1055a7e14dcfSSatish Balay PetscErrorCode ierr; 1056a7e14dcfSSatish Balay 1057a7e14dcfSSatish Balay PetscFunctionBegin; 10583c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&nlsP);CHKERRQ(ierr); 1059a7e14dcfSSatish Balay 1060a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_NLS; 1061a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_NLS; 1062a7e14dcfSSatish Balay tao->ops->view = TaoView_NLS; 1063a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_NLS; 1064a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_NLS; 1065a7e14dcfSSatish Balay 1066a7e14dcfSSatish Balay tao->max_it = 50; 10676f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 10686f4723b1SBarry Smith tao->fatol = 1e-5; 10696f4723b1SBarry Smith tao->frtol = 1e-5; 10706f4723b1SBarry Smith #else 1071a7e14dcfSSatish Balay tao->fatol = 1e-10; 1072a7e14dcfSSatish Balay tao->frtol = 1e-10; 10736f4723b1SBarry Smith #endif 1074a7e14dcfSSatish Balay tao->data = (void*)nlsP; 1075a7e14dcfSSatish Balay tao->trust0 = 100.0; 1076a7e14dcfSSatish Balay 1077a7e14dcfSSatish Balay nlsP->sval = 0.0; 1078a7e14dcfSSatish Balay nlsP->imin = 1.0e-4; 1079a7e14dcfSSatish Balay nlsP->imax = 1.0e+2; 1080a7e14dcfSSatish Balay nlsP->imfac = 1.0e-1; 1081a7e14dcfSSatish Balay 1082a7e14dcfSSatish Balay nlsP->pmin = 1.0e-12; 1083a7e14dcfSSatish Balay nlsP->pmax = 1.0e+2; 1084a7e14dcfSSatish Balay nlsP->pgfac = 1.0e+1; 1085a7e14dcfSSatish Balay nlsP->psfac = 4.0e-1; 1086a7e14dcfSSatish Balay nlsP->pmgfac = 1.0e-1; 1087a7e14dcfSSatish Balay nlsP->pmsfac = 1.0e-1; 1088a7e14dcfSSatish Balay 1089a7e14dcfSSatish Balay /* Default values for trust-region radius update based on steplength */ 1090a7e14dcfSSatish Balay nlsP->nu1 = 0.25; 1091a7e14dcfSSatish Balay nlsP->nu2 = 0.50; 1092a7e14dcfSSatish Balay nlsP->nu3 = 1.00; 1093a7e14dcfSSatish Balay nlsP->nu4 = 1.25; 1094a7e14dcfSSatish Balay 1095a7e14dcfSSatish Balay nlsP->omega1 = 0.25; 1096a7e14dcfSSatish Balay nlsP->omega2 = 0.50; 1097a7e14dcfSSatish Balay nlsP->omega3 = 1.00; 1098a7e14dcfSSatish Balay nlsP->omega4 = 2.00; 1099a7e14dcfSSatish Balay nlsP->omega5 = 4.00; 1100a7e14dcfSSatish Balay 1101a7e14dcfSSatish Balay /* Default values for trust-region radius update based on reduction */ 1102a7e14dcfSSatish Balay nlsP->eta1 = 1.0e-4; 1103a7e14dcfSSatish Balay nlsP->eta2 = 0.25; 1104a7e14dcfSSatish Balay nlsP->eta3 = 0.50; 1105a7e14dcfSSatish Balay nlsP->eta4 = 0.90; 1106a7e14dcfSSatish Balay 1107a7e14dcfSSatish Balay nlsP->alpha1 = 0.25; 1108a7e14dcfSSatish Balay nlsP->alpha2 = 0.50; 1109a7e14dcfSSatish Balay nlsP->alpha3 = 1.00; 1110a7e14dcfSSatish Balay nlsP->alpha4 = 2.00; 1111a7e14dcfSSatish Balay nlsP->alpha5 = 4.00; 1112a7e14dcfSSatish Balay 1113a7e14dcfSSatish Balay /* Default values for trust-region radius update based on interpolation */ 1114a7e14dcfSSatish Balay nlsP->mu1 = 0.10; 1115a7e14dcfSSatish Balay nlsP->mu2 = 0.50; 1116a7e14dcfSSatish Balay 1117a7e14dcfSSatish Balay nlsP->gamma1 = 0.25; 1118a7e14dcfSSatish Balay nlsP->gamma2 = 0.50; 1119a7e14dcfSSatish Balay nlsP->gamma3 = 2.00; 1120a7e14dcfSSatish Balay nlsP->gamma4 = 4.00; 1121a7e14dcfSSatish Balay 1122a7e14dcfSSatish Balay nlsP->theta = 0.05; 1123a7e14dcfSSatish Balay 1124a7e14dcfSSatish Balay /* Default values for trust region initialization based on interpolation */ 1125a7e14dcfSSatish Balay nlsP->mu1_i = 0.35; 1126a7e14dcfSSatish Balay nlsP->mu2_i = 0.50; 1127a7e14dcfSSatish Balay 1128a7e14dcfSSatish Balay nlsP->gamma1_i = 0.0625; 1129a7e14dcfSSatish Balay nlsP->gamma2_i = 0.5; 1130a7e14dcfSSatish Balay nlsP->gamma3_i = 2.0; 1131a7e14dcfSSatish Balay nlsP->gamma4_i = 5.0; 1132a7e14dcfSSatish Balay 1133a7e14dcfSSatish Balay nlsP->theta_i = 0.25; 1134a7e14dcfSSatish Balay 1135a7e14dcfSSatish Balay /* Remaining parameters */ 1136a7e14dcfSSatish Balay nlsP->min_radius = 1.0e-10; 1137a7e14dcfSSatish Balay nlsP->max_radius = 1.0e10; 1138a7e14dcfSSatish Balay nlsP->epsilon = 1.0e-6; 1139a7e14dcfSSatish Balay 1140a7e14dcfSSatish Balay nlsP->ksp_type = NLS_KSP_STCG; 1141a7e14dcfSSatish Balay nlsP->pc_type = NLS_PC_BFGS; 1142a7e14dcfSSatish Balay nlsP->bfgs_scale_type = BFGS_SCALE_PHESS; 1143a7e14dcfSSatish Balay nlsP->init_type = NLS_INIT_INTERPOLATION; 1144a7e14dcfSSatish Balay nlsP->update_type = NLS_UPDATE_STEP; 1145a7e14dcfSSatish Balay 1146a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 1147a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 1148441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 1149a7e14dcfSSatish Balay 1150a7e14dcfSSatish Balay /* Set linear solver to default for symmetric matrices */ 1151a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 1152a7e14dcfSSatish Balay PetscFunctionReturn(0); 1153a7e14dcfSSatish Balay } 1154a7e14dcfSSatish Balay EXTERN_C_END 1155a7e14dcfSSatish Balay 1156a7e14dcfSSatish Balay #undef __FUNCT__ 1157a7e14dcfSSatish Balay #define __FUNCT__ "MatLMVMSolveShell" 1158a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x) 1159a7e14dcfSSatish Balay { 1160a7e14dcfSSatish Balay PetscErrorCode ierr; 1161a7e14dcfSSatish Balay Mat M; 116287f595a5SBarry Smith 1163a7e14dcfSSatish Balay PetscFunctionBegin; 1164a7e14dcfSSatish Balay PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1165a7e14dcfSSatish Balay PetscValidHeaderSpecific(b,VEC_CLASSID,2); 1166a7e14dcfSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,3); 1167a7e14dcfSSatish Balay ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr); 1168a7e14dcfSSatish Balay ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr); 1169a7e14dcfSSatish Balay PetscFunctionReturn(0); 1170a7e14dcfSSatish Balay } 1171