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 1052d9aa51bSJason Sarich /* Number of times ksp stopped because of these reasons */ 106a7e14dcfSSatish Balay nlsP->ksp_atol = 0; 107a7e14dcfSSatish Balay nlsP->ksp_rtol = 0; 108a7e14dcfSSatish Balay nlsP->ksp_dtol = 0; 109a7e14dcfSSatish Balay nlsP->ksp_ctol = 0; 110a7e14dcfSSatish Balay nlsP->ksp_negc = 0; 111a7e14dcfSSatish Balay nlsP->ksp_iter = 0; 112a7e14dcfSSatish Balay nlsP->ksp_othr = 0; 113a7e14dcfSSatish Balay 114a7e14dcfSSatish Balay /* Modify the linear solver to a trust region method if desired */ 115a7e14dcfSSatish Balay switch(nlsP->ksp_type) { 116a7e14dcfSSatish Balay case NLS_KSP_CG: 117a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPCG);CHKERRQ(ierr); 11872b7fd4bSBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 119a7e14dcfSSatish Balay break; 120a7e14dcfSSatish Balay 121a7e14dcfSSatish Balay case NLS_KSP_NASH: 122a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPNASH);CHKERRQ(ierr); 12372b7fd4bSBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 124a7e14dcfSSatish Balay break; 125a7e14dcfSSatish Balay 126a7e14dcfSSatish Balay case NLS_KSP_STCG: 127a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPSTCG);CHKERRQ(ierr); 12872b7fd4bSBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 129a7e14dcfSSatish Balay break; 130a7e14dcfSSatish Balay 131a7e14dcfSSatish Balay case NLS_KSP_GLTR: 132a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp, KSPGLTR);CHKERRQ(ierr); 13372b7fd4bSBarry Smith ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 134a7e14dcfSSatish Balay break; 135a7e14dcfSSatish Balay 136a7e14dcfSSatish Balay default: 137a7e14dcfSSatish Balay /* Use the method set by the ksp_type */ 138a7e14dcfSSatish Balay break; 139a7e14dcfSSatish Balay } 140a7e14dcfSSatish Balay 141a7e14dcfSSatish Balay /* Initialize trust-region radius when using nash, stcg, or gltr 142a7e14dcfSSatish Balay Will be reset during the first iteration */ 143a7e14dcfSSatish Balay if (NLS_KSP_NASH == nlsP->ksp_type) { 144a7e14dcfSSatish Balay ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 145a7e14dcfSSatish Balay } else if (NLS_KSP_STCG == nlsP->ksp_type) { 146a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 147a7e14dcfSSatish Balay } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 148a7e14dcfSSatish Balay ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 149a7e14dcfSSatish Balay } 150a7e14dcfSSatish Balay 15187f595a5SBarry Smith if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) { 152a7e14dcfSSatish Balay tao->trust = tao->trust0; 153a7e14dcfSSatish Balay 15487f595a5SBarry Smith if (tao->trust < 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial radius negative"); 155a7e14dcfSSatish Balay 156a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 157a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 158a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 159a7e14dcfSSatish Balay } 160a7e14dcfSSatish Balay 161a7e14dcfSSatish Balay /* Get vectors we will need */ 162a7e14dcfSSatish Balay if (NLS_PC_BFGS == nlsP->pc_type && !nlsP->M) { 163a7e14dcfSSatish Balay ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 164a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 165a7e14dcfSSatish Balay ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&nlsP->M);CHKERRQ(ierr); 166a7e14dcfSSatish Balay ierr = MatLMVMAllocateVectors(nlsP->M,tao->solution);CHKERRQ(ierr); 167a7e14dcfSSatish Balay } 168a7e14dcfSSatish Balay 169a7e14dcfSSatish Balay /* Check convergence criteria */ 170a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 171a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 17287f595a5SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 173a7e14dcfSSatish Balay needH = 1; 174a7e14dcfSSatish Balay 175a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 17687f595a5SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 177a7e14dcfSSatish Balay 178a7e14dcfSSatish Balay /* create vectors for the limited memory preconditioner */ 17987f595a5SBarry Smith if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_BFGS != nlsP->bfgs_scale_type)) { 180a7e14dcfSSatish Balay if (!nlsP->Diag) { 181a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&nlsP->Diag);CHKERRQ(ierr); 182a7e14dcfSSatish Balay } 183a7e14dcfSSatish Balay } 184a7e14dcfSSatish Balay 185a7e14dcfSSatish Balay /* Modify the preconditioner to use the bfgs approximation */ 186a7e14dcfSSatish Balay ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 187a7e14dcfSSatish Balay switch(nlsP->pc_type) { 188a7e14dcfSSatish Balay case NLS_PC_NONE: 189a7e14dcfSSatish Balay ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr); 190*1a1499c8SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 191a7e14dcfSSatish Balay break; 192a7e14dcfSSatish Balay 193a7e14dcfSSatish Balay case NLS_PC_AHESS: 194a7e14dcfSSatish Balay ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr); 195*1a1499c8SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 196baa89ecbSBarry Smith ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 197a7e14dcfSSatish Balay break; 198a7e14dcfSSatish Balay 199a7e14dcfSSatish Balay case NLS_PC_BFGS: 200a7e14dcfSSatish Balay ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr); 201*1a1499c8SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 202a7e14dcfSSatish Balay ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr); 203a7e14dcfSSatish Balay ierr = PCShellSetContext(pc, nlsP->M);CHKERRQ(ierr); 204a7e14dcfSSatish Balay ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr); 205a7e14dcfSSatish Balay break; 206a7e14dcfSSatish Balay 207a7e14dcfSSatish Balay default: 208a7e14dcfSSatish Balay /* Use the pc method set by pc_type */ 209a7e14dcfSSatish Balay break; 210a7e14dcfSSatish Balay } 211a7e14dcfSSatish Balay 212a7e14dcfSSatish Balay /* Initialize trust-region radius. The initialization is only performed 213a7e14dcfSSatish Balay when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */ 21487f595a5SBarry Smith if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) { 215a7e14dcfSSatish Balay switch(nlsP->init_type) { 216a7e14dcfSSatish Balay case NLS_INIT_CONSTANT: 217a7e14dcfSSatish Balay /* Use the initial radius specified */ 218a7e14dcfSSatish Balay break; 219a7e14dcfSSatish Balay 220a7e14dcfSSatish Balay case NLS_INIT_INTERPOLATION: 221a7e14dcfSSatish Balay /* Use the initial radius specified */ 222a7e14dcfSSatish Balay max_radius = 0.0; 223a7e14dcfSSatish Balay 224a7e14dcfSSatish Balay for (j = 0; j < j_max; ++j) { 225a7e14dcfSSatish Balay fmin = f; 226a7e14dcfSSatish Balay sigma = 0.0; 227a7e14dcfSSatish Balay 228a7e14dcfSSatish Balay if (needH) { 229ffad9901SBarry Smith ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 230a7e14dcfSSatish Balay needH = 0; 231a7e14dcfSSatish Balay } 232a7e14dcfSSatish Balay 233a7e14dcfSSatish Balay for (i = 0; i < i_max; ++i) { 234a7e14dcfSSatish Balay ierr = VecCopy(tao->solution,nlsP->W);CHKERRQ(ierr); 235a7e14dcfSSatish Balay ierr = VecAXPY(nlsP->W,-tao->trust/gnorm,tao->gradient);CHKERRQ(ierr); 236a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao, nlsP->W, &ftrial);CHKERRQ(ierr); 237a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) { 238a7e14dcfSSatish Balay tau = nlsP->gamma1_i; 23987f595a5SBarry Smith } else { 240a7e14dcfSSatish Balay if (ftrial < fmin) { 241a7e14dcfSSatish Balay fmin = ftrial; 242a7e14dcfSSatish Balay sigma = -tao->trust / gnorm; 243a7e14dcfSSatish Balay } 244a7e14dcfSSatish Balay 245a7e14dcfSSatish Balay ierr = MatMult(tao->hessian, tao->gradient, nlsP->D);CHKERRQ(ierr); 246a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, nlsP->D, &prered);CHKERRQ(ierr); 247a7e14dcfSSatish Balay 248a7e14dcfSSatish Balay prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm)); 249a7e14dcfSSatish Balay actred = f - ftrial; 25087f595a5SBarry Smith if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) { 251a7e14dcfSSatish Balay kappa = 1.0; 25287f595a5SBarry Smith } else { 253a7e14dcfSSatish Balay kappa = actred / prered; 254a7e14dcfSSatish Balay } 255a7e14dcfSSatish Balay 256a7e14dcfSSatish Balay tau_1 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred); 257a7e14dcfSSatish Balay tau_2 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred); 258a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 259a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 260a7e14dcfSSatish Balay 261a7e14dcfSSatish Balay if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu1_i) { 262a7e14dcfSSatish Balay /* Great agreement */ 263a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 264a7e14dcfSSatish Balay 265a7e14dcfSSatish Balay if (tau_max < 1.0) { 266a7e14dcfSSatish Balay tau = nlsP->gamma3_i; 26787f595a5SBarry Smith } else if (tau_max > nlsP->gamma4_i) { 268a7e14dcfSSatish Balay tau = nlsP->gamma4_i; 26987f595a5SBarry Smith } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) { 270a7e14dcfSSatish Balay tau = tau_1; 27187f595a5SBarry Smith } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) { 272a7e14dcfSSatish Balay tau = tau_2; 27387f595a5SBarry Smith } else { 274a7e14dcfSSatish Balay tau = tau_max; 275a7e14dcfSSatish Balay } 27687f595a5SBarry Smith } else if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu2_i) { 277a7e14dcfSSatish Balay /* Good agreement */ 278a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust); 279a7e14dcfSSatish Balay 280a7e14dcfSSatish Balay if (tau_max < nlsP->gamma2_i) { 281a7e14dcfSSatish Balay tau = nlsP->gamma2_i; 28287f595a5SBarry Smith } else if (tau_max > nlsP->gamma3_i) { 283a7e14dcfSSatish Balay tau = nlsP->gamma3_i; 28487f595a5SBarry Smith } else { 285a7e14dcfSSatish Balay tau = tau_max; 286a7e14dcfSSatish Balay } 28787f595a5SBarry Smith } else { 288a7e14dcfSSatish Balay /* Not good agreement */ 289a7e14dcfSSatish Balay if (tau_min > 1.0) { 290a7e14dcfSSatish Balay tau = nlsP->gamma2_i; 29187f595a5SBarry Smith } else if (tau_max < nlsP->gamma1_i) { 292a7e14dcfSSatish Balay tau = nlsP->gamma1_i; 29387f595a5SBarry Smith } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) { 294a7e14dcfSSatish Balay tau = nlsP->gamma1_i; 29587f595a5SBarry Smith } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) { 296a7e14dcfSSatish Balay tau = tau_1; 29787f595a5SBarry Smith } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) { 298a7e14dcfSSatish Balay tau = tau_2; 29987f595a5SBarry Smith } else { 300a7e14dcfSSatish Balay tau = tau_max; 301a7e14dcfSSatish Balay } 302a7e14dcfSSatish Balay } 303a7e14dcfSSatish Balay } 304a7e14dcfSSatish Balay tao->trust = tau * tao->trust; 305a7e14dcfSSatish Balay } 306a7e14dcfSSatish Balay 307a7e14dcfSSatish Balay if (fmin < f) { 308a7e14dcfSSatish Balay f = fmin; 309a7e14dcfSSatish Balay ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr); 310a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao,tao->solution,tao->gradient);CHKERRQ(ierr); 311a7e14dcfSSatish Balay 312a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 31387f595a5SBarry Smith if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN"); 314a7e14dcfSSatish Balay needH = 1; 315a7e14dcfSSatish Balay 316a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 31787f595a5SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 318a7e14dcfSSatish Balay } 319a7e14dcfSSatish Balay } 320a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, max_radius); 321a7e14dcfSSatish Balay 322a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 323a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 324a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 325a7e14dcfSSatish Balay break; 326a7e14dcfSSatish Balay 327a7e14dcfSSatish Balay default: 328a7e14dcfSSatish Balay /* Norm of the first direction will initialize radius */ 329a7e14dcfSSatish Balay tao->trust = 0.0; 330a7e14dcfSSatish Balay break; 331a7e14dcfSSatish Balay } 332a7e14dcfSSatish Balay } 333a7e14dcfSSatish Balay 334a7e14dcfSSatish Balay /* Set initial scaling for the BFGS preconditioner 335a7e14dcfSSatish Balay This step is done after computing the initial trust-region radius 336a7e14dcfSSatish Balay since the function value may have decreased */ 337a7e14dcfSSatish Balay if (NLS_PC_BFGS == nlsP->pc_type) { 338a7e14dcfSSatish Balay if (f != 0.0) { 339a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 34087f595a5SBarry Smith } else { 341a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 342a7e14dcfSSatish Balay } 343a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr); 344a7e14dcfSSatish Balay } 345a7e14dcfSSatish Balay 346a7e14dcfSSatish Balay /* Set counter for gradient/reset steps*/ 347a7e14dcfSSatish Balay nlsP->newt = 0; 348a7e14dcfSSatish Balay nlsP->bfgs = 0; 349a7e14dcfSSatish Balay nlsP->sgrad = 0; 350a7e14dcfSSatish Balay nlsP->grad = 0; 351a7e14dcfSSatish Balay 352a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 353a7e14dcfSSatish Balay while (reason == TAO_CONTINUE_ITERATING) { 354a7e14dcfSSatish Balay ++iter; 355ae93cb3cSJason Sarich tao->ksp_its=0; 356a7e14dcfSSatish Balay 357a7e14dcfSSatish Balay /* Compute the Hessian */ 358a7e14dcfSSatish Balay if (needH) { 359ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 360a7e14dcfSSatish Balay needH = 0; 361a7e14dcfSSatish Balay } 362a7e14dcfSSatish Balay 36387f595a5SBarry Smith if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_AHESS == nlsP->bfgs_scale_type)) { 364a7e14dcfSSatish Balay /* Obtain diagonal for the bfgs preconditioner */ 365a7e14dcfSSatish Balay ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr); 366a7e14dcfSSatish Balay ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr); 367a7e14dcfSSatish Balay ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr); 368a7e14dcfSSatish Balay ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr); 369a7e14dcfSSatish Balay } 370a7e14dcfSSatish Balay 371a7e14dcfSSatish Balay /* Shift the Hessian matrix */ 372a7e14dcfSSatish Balay if (pert > 0) { 373a7e14dcfSSatish Balay ierr = MatShift(tao->hessian, pert); 374a7e14dcfSSatish Balay if (tao->hessian != tao->hessian_pre) { 375a7e14dcfSSatish Balay ierr = MatShift(tao->hessian_pre, pert);CHKERRQ(ierr); 376a7e14dcfSSatish Balay } 377a7e14dcfSSatish Balay } 378a7e14dcfSSatish Balay 379a7e14dcfSSatish Balay if (NLS_PC_BFGS == nlsP->pc_type) { 380a7e14dcfSSatish Balay if (BFGS_SCALE_PHESS == nlsP->bfgs_scale_type) { 381a7e14dcfSSatish Balay /* Obtain diagonal for the bfgs preconditioner */ 382a7e14dcfSSatish Balay ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr); 383a7e14dcfSSatish Balay ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr); 384a7e14dcfSSatish Balay ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr); 385a7e14dcfSSatish Balay ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr); 386a7e14dcfSSatish Balay } 387a7e14dcfSSatish Balay /* Update the limited memory preconditioner */ 388a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 389a7e14dcfSSatish Balay ++bfgsUpdates; 390a7e14dcfSSatish Balay } 391a7e14dcfSSatish Balay 392a7e14dcfSSatish Balay /* Solve the Newton system of equations */ 39323ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 39487f595a5SBarry Smith if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) { 395a7e14dcfSSatish Balay 396a7e14dcfSSatish Balay if (NLS_KSP_NASH == nlsP->ksp_type) { 397a7e14dcfSSatish Balay ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 398a7e14dcfSSatish Balay } else if (NLS_KSP_STCG == nlsP->ksp_type) { 399a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 400a7e14dcfSSatish Balay } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 401a7e14dcfSSatish Balay ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 402a7e14dcfSSatish Balay } 403a7e14dcfSSatish Balay 404a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr); 405a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 406a7e14dcfSSatish Balay tao->ksp_its+=kspits; 407ae93cb3cSJason Sarich tao->ksp_tot_its+=kspits; 408a7e14dcfSSatish Balay 409a7e14dcfSSatish Balay if (NLS_KSP_NASH == nlsP->ksp_type) { 410a7e14dcfSSatish Balay ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 411a7e14dcfSSatish Balay } else if (NLS_KSP_STCG == nlsP->ksp_type) { 412a7e14dcfSSatish Balay ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 413a7e14dcfSSatish Balay } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 414a7e14dcfSSatish Balay ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 415a7e14dcfSSatish Balay } 416a7e14dcfSSatish Balay 417a7e14dcfSSatish Balay if (0.0 == tao->trust) { 418a7e14dcfSSatish Balay /* Radius was uninitialized; use the norm of the direction */ 419a7e14dcfSSatish Balay if (norm_d > 0.0) { 420a7e14dcfSSatish Balay tao->trust = norm_d; 421a7e14dcfSSatish Balay 422a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 423a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 424a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 42587f595a5SBarry Smith } else { 426a7e14dcfSSatish Balay /* The direction was bad; set radius to default value and re-solve 427a7e14dcfSSatish Balay the trust-region subproblem to get a direction */ 428a7e14dcfSSatish Balay tao->trust = tao->trust0; 429a7e14dcfSSatish Balay 430a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */ 431a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->min_radius); 432a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 433a7e14dcfSSatish Balay 434a7e14dcfSSatish Balay if (NLS_KSP_NASH == nlsP->ksp_type) { 435a7e14dcfSSatish Balay ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 436a7e14dcfSSatish Balay } else if (NLS_KSP_STCG == nlsP->ksp_type) { 437a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 438a7e14dcfSSatish Balay } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 439a7e14dcfSSatish Balay ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 440a7e14dcfSSatish Balay } 441a7e14dcfSSatish Balay 442a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr); 443a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 444a7e14dcfSSatish Balay tao->ksp_its+=kspits; 445ae93cb3cSJason Sarich tao->ksp_tot_its+=kspits; 446a7e14dcfSSatish Balay if (NLS_KSP_NASH == nlsP->ksp_type) { 447a7e14dcfSSatish Balay ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 448a7e14dcfSSatish Balay } else if (NLS_KSP_STCG == nlsP->ksp_type) { 449a7e14dcfSSatish Balay ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 450a7e14dcfSSatish Balay } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 451a7e14dcfSSatish Balay ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 452a7e14dcfSSatish Balay } 453a7e14dcfSSatish Balay 45487f595a5SBarry Smith if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 455a7e14dcfSSatish Balay } 456a7e14dcfSSatish Balay } 45787f595a5SBarry Smith } else { 458a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr); 459a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr); 460a7e14dcfSSatish Balay tao->ksp_its += kspits; 461ae93cb3cSJason Sarich tao->ksp_tot_its+=kspits; 462a7e14dcfSSatish Balay } 463a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 464a7e14dcfSSatish Balay ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 46587f595a5SBarry Smith if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (NLS_PC_BFGS == nlsP->pc_type) && (bfgsUpdates > 1)) { 466a7e14dcfSSatish Balay /* Preconditioner is numerically indefinite; reset the 467a7e14dcfSSatish Balay approximate if using BFGS preconditioning. */ 468a7e14dcfSSatish Balay 469a7e14dcfSSatish Balay if (f != 0.0) { 470a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 47187f595a5SBarry Smith } else { 472a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 473a7e14dcfSSatish Balay } 474a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr); 475a7e14dcfSSatish Balay ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 476a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 477a7e14dcfSSatish Balay bfgsUpdates = 1; 478a7e14dcfSSatish Balay } 479a7e14dcfSSatish Balay 480a7e14dcfSSatish Balay if (KSP_CONVERGED_ATOL == ksp_reason) { 481a7e14dcfSSatish Balay ++nlsP->ksp_atol; 48287f595a5SBarry Smith } else if (KSP_CONVERGED_RTOL == ksp_reason) { 483a7e14dcfSSatish Balay ++nlsP->ksp_rtol; 48487f595a5SBarry Smith } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) { 485a7e14dcfSSatish Balay ++nlsP->ksp_ctol; 48687f595a5SBarry Smith } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) { 487a7e14dcfSSatish Balay ++nlsP->ksp_negc; 48887f595a5SBarry Smith } else if (KSP_DIVERGED_DTOL == ksp_reason) { 489a7e14dcfSSatish Balay ++nlsP->ksp_dtol; 49087f595a5SBarry Smith } else if (KSP_DIVERGED_ITS == ksp_reason) { 491a7e14dcfSSatish Balay ++nlsP->ksp_iter; 49287f595a5SBarry Smith } else { 493a7e14dcfSSatish Balay ++nlsP->ksp_othr; 494a7e14dcfSSatish Balay } 495a7e14dcfSSatish Balay 496a7e14dcfSSatish Balay /* Check for success (descent direction) */ 497a7e14dcfSSatish Balay ierr = VecDot(nlsP->D, tao->gradient, &gdx);CHKERRQ(ierr); 498a7e14dcfSSatish Balay if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 499a7e14dcfSSatish Balay /* Newton step is not descent or direction produced Inf or NaN 500a7e14dcfSSatish Balay Update the perturbation for next time */ 501a7e14dcfSSatish Balay if (pert <= 0.0) { 502a7e14dcfSSatish Balay /* Initialize the perturbation */ 503a7e14dcfSSatish Balay pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 504a7e14dcfSSatish Balay if (NLS_KSP_GLTR == nlsP->ksp_type) { 505a7e14dcfSSatish Balay ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 506a7e14dcfSSatish Balay pert = PetscMax(pert, -e_min); 507a7e14dcfSSatish Balay } 50887f595a5SBarry Smith } else { 509a7e14dcfSSatish Balay /* Increase the perturbation */ 510a7e14dcfSSatish Balay pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 511a7e14dcfSSatish Balay } 512a7e14dcfSSatish Balay 513a7e14dcfSSatish Balay if (NLS_PC_BFGS != nlsP->pc_type) { 514a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and updated 515a7e14dcfSSatish Balay Must use gradient direction in this case */ 516a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr); 517a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 518a7e14dcfSSatish Balay ++nlsP->grad; 519a7e14dcfSSatish Balay stepType = NLS_GRADIENT; 52087f595a5SBarry Smith } else { 521a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 522a7e14dcfSSatish Balay ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 523a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 524a7e14dcfSSatish Balay 525a7e14dcfSSatish Balay /* Check for success (descent direction) */ 526a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, nlsP->D, &gdx);CHKERRQ(ierr); 527a7e14dcfSSatish Balay if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) { 528a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced not a number 529a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case because 530a7e14dcfSSatish Balay the first solve produces the scaled gradient direction, 531a7e14dcfSSatish Balay which is guaranteed to be descent */ 532a7e14dcfSSatish Balay 533a7e14dcfSSatish Balay /* Use steepest descent direction (scaled) */ 534a7e14dcfSSatish Balay 535a7e14dcfSSatish Balay if (f != 0.0) { 536a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 53787f595a5SBarry Smith } else { 538a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 539a7e14dcfSSatish Balay } 540a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr); 541a7e14dcfSSatish Balay ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 542a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 543a7e14dcfSSatish Balay ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 544a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 545a7e14dcfSSatish Balay 546a7e14dcfSSatish Balay bfgsUpdates = 1; 547a7e14dcfSSatish Balay ++nlsP->sgrad; 548a7e14dcfSSatish Balay stepType = NLS_SCALED_GRADIENT; 54987f595a5SBarry Smith } else { 550a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 551a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 552a7e14dcfSSatish Balay ++nlsP->sgrad; 553a7e14dcfSSatish Balay stepType = NLS_SCALED_GRADIENT; 55487f595a5SBarry Smith } else { 555a7e14dcfSSatish Balay ++nlsP->bfgs; 556a7e14dcfSSatish Balay stepType = NLS_BFGS; 557a7e14dcfSSatish Balay } 558a7e14dcfSSatish Balay } 559a7e14dcfSSatish Balay } 56087f595a5SBarry Smith } else { 561a7e14dcfSSatish Balay /* Computed Newton step is descent */ 562a7e14dcfSSatish Balay switch (ksp_reason) { 563a7e14dcfSSatish Balay case KSP_DIVERGED_NANORINF: 564a7e14dcfSSatish Balay case KSP_DIVERGED_BREAKDOWN: 565a7e14dcfSSatish Balay case KSP_DIVERGED_INDEFINITE_MAT: 566a7e14dcfSSatish Balay case KSP_DIVERGED_INDEFINITE_PC: 567a7e14dcfSSatish Balay case KSP_CONVERGED_CG_NEG_CURVE: 568a7e14dcfSSatish Balay /* Matrix or preconditioner is indefinite; increase perturbation */ 569a7e14dcfSSatish Balay if (pert <= 0.0) { 570a7e14dcfSSatish Balay /* Initialize the perturbation */ 571a7e14dcfSSatish Balay pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 572a7e14dcfSSatish Balay if (NLS_KSP_GLTR == nlsP->ksp_type) { 573a7e14dcfSSatish Balay ierr = KSPGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr); 574a7e14dcfSSatish Balay pert = PetscMax(pert, -e_min); 575a7e14dcfSSatish Balay } 57687f595a5SBarry Smith } else { 577a7e14dcfSSatish Balay /* Increase the perturbation */ 578a7e14dcfSSatish Balay pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 579a7e14dcfSSatish Balay } 580a7e14dcfSSatish Balay break; 581a7e14dcfSSatish Balay 582a7e14dcfSSatish Balay default: 583a7e14dcfSSatish Balay /* Newton step computation is good; decrease perturbation */ 584a7e14dcfSSatish Balay pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm); 585a7e14dcfSSatish Balay if (pert < nlsP->pmin) { 586a7e14dcfSSatish Balay pert = 0.0; 587a7e14dcfSSatish Balay } 588a7e14dcfSSatish Balay break; 589a7e14dcfSSatish Balay } 590a7e14dcfSSatish Balay 591a7e14dcfSSatish Balay ++nlsP->newt; 592a7e14dcfSSatish Balay stepType = NLS_NEWTON; 593a7e14dcfSSatish Balay } 594a7e14dcfSSatish Balay 595a7e14dcfSSatish Balay /* Perform the linesearch */ 596a7e14dcfSSatish Balay fold = f; 597a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, nlsP->Xold);CHKERRQ(ierr); 598a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, nlsP->Gold);CHKERRQ(ierr); 599a7e14dcfSSatish Balay 600a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr); 601a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 602a7e14dcfSSatish Balay 60387f595a5SBarry Smith while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) { /* Linesearch failed */ 604a7e14dcfSSatish Balay f = fold; 605a7e14dcfSSatish Balay ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr); 606a7e14dcfSSatish Balay ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr); 607a7e14dcfSSatish Balay 608a7e14dcfSSatish Balay switch(stepType) { 609a7e14dcfSSatish Balay case NLS_NEWTON: 610a7e14dcfSSatish Balay /* Failed to obtain acceptable iterate with Newton 1step 611a7e14dcfSSatish Balay Update the perturbation for next time */ 612a7e14dcfSSatish Balay if (pert <= 0.0) { 613a7e14dcfSSatish Balay /* Initialize the perturbation */ 614a7e14dcfSSatish Balay pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 615a7e14dcfSSatish Balay if (NLS_KSP_GLTR == nlsP->ksp_type) { 616a7e14dcfSSatish Balay ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 617a7e14dcfSSatish Balay pert = PetscMax(pert, -e_min); 618a7e14dcfSSatish Balay } 61987f595a5SBarry Smith } else { 620a7e14dcfSSatish Balay /* Increase the perturbation */ 621a7e14dcfSSatish Balay pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 622a7e14dcfSSatish Balay } 623a7e14dcfSSatish Balay 624a7e14dcfSSatish Balay if (NLS_PC_BFGS != nlsP->pc_type) { 625a7e14dcfSSatish Balay /* We don't have the bfgs matrix around and being updated 626a7e14dcfSSatish Balay Must use gradient direction in this case */ 627a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr); 628a7e14dcfSSatish Balay ++nlsP->grad; 629a7e14dcfSSatish Balay stepType = NLS_GRADIENT; 63087f595a5SBarry Smith } else { 631a7e14dcfSSatish Balay /* Attempt to use the BFGS direction */ 632a7e14dcfSSatish Balay ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 633a7e14dcfSSatish Balay /* Check for success (descent direction) */ 634a7e14dcfSSatish Balay ierr = VecDot(tao->solution, nlsP->D, &gdx);CHKERRQ(ierr); 635a7e14dcfSSatish Balay if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 636a7e14dcfSSatish Balay /* BFGS direction is not descent or direction produced not a number 637a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case 638a7e14dcfSSatish Balay Use steepest descent direction (scaled) */ 639a7e14dcfSSatish Balay 640a7e14dcfSSatish Balay if (f != 0.0) { 641a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 64287f595a5SBarry Smith } else { 643a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 644a7e14dcfSSatish Balay } 645a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr); 646a7e14dcfSSatish Balay ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 647a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 648a7e14dcfSSatish Balay ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 649a7e14dcfSSatish Balay 650a7e14dcfSSatish Balay bfgsUpdates = 1; 651a7e14dcfSSatish Balay ++nlsP->sgrad; 652a7e14dcfSSatish Balay stepType = NLS_SCALED_GRADIENT; 65387f595a5SBarry Smith } else { 654a7e14dcfSSatish Balay if (1 == bfgsUpdates) { 655a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 656a7e14dcfSSatish Balay ++nlsP->sgrad; 657a7e14dcfSSatish Balay stepType = NLS_SCALED_GRADIENT; 65887f595a5SBarry Smith } else { 659a7e14dcfSSatish Balay ++nlsP->bfgs; 660a7e14dcfSSatish Balay stepType = NLS_BFGS; 661a7e14dcfSSatish Balay } 662a7e14dcfSSatish Balay } 663a7e14dcfSSatish Balay } 664a7e14dcfSSatish Balay break; 665a7e14dcfSSatish Balay 666a7e14dcfSSatish Balay case NLS_BFGS: 667a7e14dcfSSatish Balay /* Can only enter if pc_type == NLS_PC_BFGS 668a7e14dcfSSatish Balay Failed to obtain acceptable iterate with BFGS step 669a7e14dcfSSatish Balay Attempt to use the scaled gradient direction */ 670a7e14dcfSSatish Balay 671a7e14dcfSSatish Balay if (f != 0.0) { 672a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 67387f595a5SBarry Smith } else { 674a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 675a7e14dcfSSatish Balay } 676a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr); 677a7e14dcfSSatish Balay ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 678a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 679a7e14dcfSSatish Balay ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 680a7e14dcfSSatish Balay 681a7e14dcfSSatish Balay bfgsUpdates = 1; 682a7e14dcfSSatish Balay ++nlsP->sgrad; 683a7e14dcfSSatish Balay stepType = NLS_SCALED_GRADIENT; 684a7e14dcfSSatish Balay break; 685a7e14dcfSSatish Balay 686a7e14dcfSSatish Balay case NLS_SCALED_GRADIENT: 687a7e14dcfSSatish Balay /* Can only enter if pc_type == NLS_PC_BFGS 688a7e14dcfSSatish Balay The scaled gradient step did not produce a new iterate; 689a7e14dcfSSatish Balay attemp to use the gradient direction. 690a7e14dcfSSatish Balay Need to make sure we are not using a different diagonal scaling */ 691a7e14dcfSSatish Balay 692a7e14dcfSSatish Balay ierr = MatLMVMSetScale(nlsP->M,0);CHKERRQ(ierr); 693a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(nlsP->M,1.0);CHKERRQ(ierr); 694a7e14dcfSSatish Balay ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 695a7e14dcfSSatish Balay ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 696a7e14dcfSSatish Balay ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 697a7e14dcfSSatish Balay 698a7e14dcfSSatish Balay bfgsUpdates = 1; 699a7e14dcfSSatish Balay ++nlsP->grad; 700a7e14dcfSSatish Balay stepType = NLS_GRADIENT; 701a7e14dcfSSatish Balay break; 702a7e14dcfSSatish Balay } 703a7e14dcfSSatish Balay ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 704a7e14dcfSSatish Balay 705a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr); 706a7e14dcfSSatish Balay ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr); 707a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 708a7e14dcfSSatish Balay } 709a7e14dcfSSatish Balay 71087f595a5SBarry Smith if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 711a7e14dcfSSatish Balay /* Failed to find an improving point */ 712a7e14dcfSSatish Balay f = fold; 713a7e14dcfSSatish Balay ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr); 714a7e14dcfSSatish Balay ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr); 715a7e14dcfSSatish Balay step = 0.0; 716a7e14dcfSSatish Balay reason = TAO_DIVERGED_LS_FAILURE; 717a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 718a7e14dcfSSatish Balay break; 719a7e14dcfSSatish Balay } 720a7e14dcfSSatish Balay 721a7e14dcfSSatish Balay /* Update trust region radius */ 72287f595a5SBarry Smith if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) { 723a7e14dcfSSatish Balay switch(nlsP->update_type) { 724a7e14dcfSSatish Balay case NLS_UPDATE_STEP: 725a7e14dcfSSatish Balay if (stepType == NLS_NEWTON) { 726a7e14dcfSSatish Balay if (step < nlsP->nu1) { 727a7e14dcfSSatish Balay /* Very bad step taken; reduce radius */ 728a7e14dcfSSatish Balay tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust); 72987f595a5SBarry Smith } else if (step < nlsP->nu2) { 730a7e14dcfSSatish Balay /* Reasonably bad step taken; reduce radius */ 731a7e14dcfSSatish Balay tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust); 73287f595a5SBarry Smith } else if (step < nlsP->nu3) { 733a7e14dcfSSatish Balay /* Reasonable step was taken; leave radius alone */ 734a7e14dcfSSatish Balay if (nlsP->omega3 < 1.0) { 735a7e14dcfSSatish Balay tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust); 73687f595a5SBarry Smith } else if (nlsP->omega3 > 1.0) { 737a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust); 738a7e14dcfSSatish Balay } 73987f595a5SBarry Smith } else if (step < nlsP->nu4) { 740a7e14dcfSSatish Balay /* Full step taken; increase the radius */ 741a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust); 74287f595a5SBarry Smith } else { 743a7e14dcfSSatish Balay /* More than full step taken; increase the radius */ 744a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust); 745a7e14dcfSSatish Balay } 74687f595a5SBarry Smith } else { 747a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 748a7e14dcfSSatish Balay tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust); 749a7e14dcfSSatish Balay } 750a7e14dcfSSatish Balay break; 751a7e14dcfSSatish Balay 752a7e14dcfSSatish Balay case NLS_UPDATE_REDUCTION: 753a7e14dcfSSatish Balay if (stepType == NLS_NEWTON) { 754a7e14dcfSSatish Balay /* Get predicted reduction */ 755a7e14dcfSSatish Balay 756a7e14dcfSSatish Balay if (NLS_KSP_STCG == nlsP->ksp_type) { 757a7e14dcfSSatish Balay ierr = KSPSTCGGetObjFcn(tao->ksp,&prered); 758a7e14dcfSSatish Balay } else if (NLS_KSP_NASH == nlsP->ksp_type) { 759a7e14dcfSSatish Balay ierr = KSPNASHGetObjFcn(tao->ksp,&prered); 760a7e14dcfSSatish Balay } else { 761a7e14dcfSSatish Balay ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 762a7e14dcfSSatish Balay } 763a7e14dcfSSatish Balay 764a7e14dcfSSatish Balay if (prered >= 0.0) { 765a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot */ 766a7e14dcfSSatish Balay /* happen in infinite precision arithmetic. Step should */ 767a7e14dcfSSatish Balay /* be rejected! */ 768a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 76987f595a5SBarry Smith } else { 770a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(f_full)) { 771a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 77287f595a5SBarry Smith } else { 773a7e14dcfSSatish Balay /* Compute and actual reduction */ 774a7e14dcfSSatish Balay actred = fold - f_full; 775a7e14dcfSSatish Balay prered = -prered; 776a7e14dcfSSatish Balay if ((PetscAbsScalar(actred) <= nlsP->epsilon) && 777a7e14dcfSSatish Balay (PetscAbsScalar(prered) <= nlsP->epsilon)) { 778a7e14dcfSSatish Balay kappa = 1.0; 77987f595a5SBarry Smith } else { 780a7e14dcfSSatish Balay kappa = actred / prered; 781a7e14dcfSSatish Balay } 782a7e14dcfSSatish Balay 783a7e14dcfSSatish Balay /* Accept of reject the step and update radius */ 784a7e14dcfSSatish Balay if (kappa < nlsP->eta1) { 785a7e14dcfSSatish Balay /* Very bad step */ 786a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 78787f595a5SBarry Smith } else if (kappa < nlsP->eta2) { 788a7e14dcfSSatish Balay /* Marginal bad step */ 789a7e14dcfSSatish Balay tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d); 79087f595a5SBarry Smith } else if (kappa < nlsP->eta3) { 791a7e14dcfSSatish Balay /* Reasonable step */ 792a7e14dcfSSatish Balay if (nlsP->alpha3 < 1.0) { 793a7e14dcfSSatish Balay tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust); 79487f595a5SBarry Smith } else if (nlsP->alpha3 > 1.0) { 795a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust); 796a7e14dcfSSatish Balay } 79787f595a5SBarry Smith } else if (kappa < nlsP->eta4) { 798a7e14dcfSSatish Balay /* Good step */ 799a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust); 80087f595a5SBarry Smith } else { 801a7e14dcfSSatish Balay /* Very good step */ 802a7e14dcfSSatish Balay tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust); 803a7e14dcfSSatish Balay } 804a7e14dcfSSatish Balay } 805a7e14dcfSSatish Balay } 80687f595a5SBarry Smith } else { 807a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 808a7e14dcfSSatish Balay tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust); 809a7e14dcfSSatish Balay } 810a7e14dcfSSatish Balay break; 811a7e14dcfSSatish Balay 812a7e14dcfSSatish Balay default: 813a7e14dcfSSatish Balay if (stepType == NLS_NEWTON) { 814a7e14dcfSSatish Balay 815a7e14dcfSSatish Balay if (NLS_KSP_STCG == nlsP->ksp_type) { 816a7e14dcfSSatish Balay ierr = KSPSTCGGetObjFcn(tao->ksp,&prered); 817a7e14dcfSSatish Balay } else if (NLS_KSP_NASH == nlsP->ksp_type) { 818a7e14dcfSSatish Balay ierr = KSPNASHGetObjFcn(tao->ksp,&prered); 819a7e14dcfSSatish Balay } else { 820a7e14dcfSSatish Balay ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 821a7e14dcfSSatish Balay } 822a7e14dcfSSatish Balay if (prered >= 0.0) { 823a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot */ 824a7e14dcfSSatish Balay /* happen in infinite precision arithmetic. Step should */ 825a7e14dcfSSatish Balay /* be rejected! */ 826a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 82787f595a5SBarry Smith } else { 828a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(f_full)) { 829a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 83087f595a5SBarry Smith } else { 831a7e14dcfSSatish Balay actred = fold - f_full; 832a7e14dcfSSatish Balay prered = -prered; 83387f595a5SBarry Smith if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) { 834a7e14dcfSSatish Balay kappa = 1.0; 83587f595a5SBarry Smith } else { 836a7e14dcfSSatish Balay kappa = actred / prered; 837a7e14dcfSSatish Balay } 838a7e14dcfSSatish Balay 839a7e14dcfSSatish Balay tau_1 = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred); 840a7e14dcfSSatish Balay tau_2 = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred); 841a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2); 842a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2); 843a7e14dcfSSatish Balay 844a7e14dcfSSatish Balay if (kappa >= 1.0 - nlsP->mu1) { 845a7e14dcfSSatish Balay /* Great agreement */ 846a7e14dcfSSatish Balay if (tau_max < 1.0) { 847a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d); 84887f595a5SBarry Smith } else if (tau_max > nlsP->gamma4) { 849a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d); 85087f595a5SBarry Smith } else { 851a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 852a7e14dcfSSatish Balay } 85387f595a5SBarry Smith } else if (kappa >= 1.0 - nlsP->mu2) { 854a7e14dcfSSatish Balay /* Good agreement */ 855a7e14dcfSSatish Balay 856a7e14dcfSSatish Balay if (tau_max < nlsP->gamma2) { 857a7e14dcfSSatish Balay tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d); 85887f595a5SBarry Smith } else if (tau_max > nlsP->gamma3) { 859a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d); 86087f595a5SBarry Smith } else if (tau_max < 1.0) { 861a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 86287f595a5SBarry Smith } else { 863a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d); 864a7e14dcfSSatish Balay } 86587f595a5SBarry Smith } else { 866a7e14dcfSSatish Balay /* Not good agreement */ 867a7e14dcfSSatish Balay if (tau_min > 1.0) { 868a7e14dcfSSatish Balay tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d); 86987f595a5SBarry Smith } else if (tau_max < nlsP->gamma1) { 870a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 87187f595a5SBarry Smith } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) { 872a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 87387f595a5SBarry Smith } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) { 874a7e14dcfSSatish Balay tao->trust = tau_1 * PetscMin(tao->trust, norm_d); 87587f595a5SBarry Smith } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) { 876a7e14dcfSSatish Balay tao->trust = tau_2 * PetscMin(tao->trust, norm_d); 87787f595a5SBarry Smith } else { 878a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d); 879a7e14dcfSSatish Balay } 880a7e14dcfSSatish Balay } 881a7e14dcfSSatish Balay } 882a7e14dcfSSatish Balay } 88387f595a5SBarry Smith } else { 884a7e14dcfSSatish Balay /* Newton step was not good; reduce the radius */ 885a7e14dcfSSatish Balay tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust); 886a7e14dcfSSatish Balay } 887a7e14dcfSSatish Balay break; 888a7e14dcfSSatish Balay } 889a7e14dcfSSatish Balay 890a7e14dcfSSatish Balay /* The radius may have been increased; modify if it is too large */ 891a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, nlsP->max_radius); 892a7e14dcfSSatish Balay } 893a7e14dcfSSatish Balay 894a7e14dcfSSatish Balay /* Check for termination */ 895a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 89687f595a5SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 897a7e14dcfSSatish Balay needH = 1; 898a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr); 899a7e14dcfSSatish Balay } 900a7e14dcfSSatish Balay PetscFunctionReturn(0); 901a7e14dcfSSatish Balay } 902a7e14dcfSSatish Balay 903a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 904a7e14dcfSSatish Balay #undef __FUNCT__ 905a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_NLS" 906441846f8SBarry Smith static PetscErrorCode TaoSetUp_NLS(Tao tao) 907a7e14dcfSSatish Balay { 908a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 909a7e14dcfSSatish Balay PetscErrorCode ierr; 910a7e14dcfSSatish Balay 911a7e14dcfSSatish Balay PetscFunctionBegin; 912a7e14dcfSSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);} 913a7e14dcfSSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);} 914a7e14dcfSSatish Balay if (!nlsP->W) {ierr = VecDuplicate(tao->solution,&nlsP->W);CHKERRQ(ierr);} 915a7e14dcfSSatish Balay if (!nlsP->D) {ierr = VecDuplicate(tao->solution,&nlsP->D);CHKERRQ(ierr);} 916a7e14dcfSSatish Balay if (!nlsP->Xold) {ierr = VecDuplicate(tao->solution,&nlsP->Xold);CHKERRQ(ierr);} 917a7e14dcfSSatish Balay if (!nlsP->Gold) {ierr = VecDuplicate(tao->solution,&nlsP->Gold);CHKERRQ(ierr);} 918a7e14dcfSSatish Balay nlsP->Diag = 0; 919a7e14dcfSSatish Balay nlsP->M = 0; 920a7e14dcfSSatish Balay PetscFunctionReturn(0); 921a7e14dcfSSatish Balay } 922a7e14dcfSSatish Balay 923a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 924a7e14dcfSSatish Balay #undef __FUNCT__ 925a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_NLS" 926441846f8SBarry Smith static PetscErrorCode TaoDestroy_NLS(Tao tao) 927a7e14dcfSSatish Balay { 928a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 929a7e14dcfSSatish Balay PetscErrorCode ierr; 930a7e14dcfSSatish Balay 931a7e14dcfSSatish Balay PetscFunctionBegin; 932a7e14dcfSSatish Balay if (tao->setupcalled) { 933a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->D);CHKERRQ(ierr); 934a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->W);CHKERRQ(ierr); 935a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->Xold);CHKERRQ(ierr); 936a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->Gold);CHKERRQ(ierr); 937a7e14dcfSSatish Balay } 938a7e14dcfSSatish Balay ierr = VecDestroy(&nlsP->Diag);CHKERRQ(ierr); 939a7e14dcfSSatish Balay ierr = MatDestroy(&nlsP->M);CHKERRQ(ierr); 940a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 941a7e14dcfSSatish Balay PetscFunctionReturn(0); 942a7e14dcfSSatish Balay } 943a7e14dcfSSatish Balay 944a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 945a7e14dcfSSatish Balay #undef __FUNCT__ 946a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_NLS" 947*1a1499c8SBarry Smith static PetscErrorCode TaoSetFromOptions_NLS(PetscOptionsObjectType *PetscOptionsObject,Tao tao) 948a7e14dcfSSatish Balay { 949a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 950a7e14dcfSSatish Balay PetscErrorCode ierr; 951a7e14dcfSSatish Balay 952a7e14dcfSSatish Balay PetscFunctionBegin; 953*1a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr); 954a7e14dcfSSatish 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); 955a7e14dcfSSatish 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); 956a7e14dcfSSatish 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); 957a7e14dcfSSatish 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); 958a7e14dcfSSatish 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); 95994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval,NULL);CHKERRQ(ierr); 96094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin,NULL);CHKERRQ(ierr); 96194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax,NULL);CHKERRQ(ierr); 96294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac,NULL);CHKERRQ(ierr); 96394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin,NULL);CHKERRQ(ierr); 96494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax,NULL);CHKERRQ(ierr); 96594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac,NULL);CHKERRQ(ierr); 96694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac,NULL);CHKERRQ(ierr); 96794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac,NULL);CHKERRQ(ierr); 96894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac,NULL);CHKERRQ(ierr); 96994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1,NULL);CHKERRQ(ierr); 97094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2,NULL);CHKERRQ(ierr); 97194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3,NULL);CHKERRQ(ierr); 97294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4,NULL);CHKERRQ(ierr); 97394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1,NULL);CHKERRQ(ierr); 97494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2,NULL);CHKERRQ(ierr); 97594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3,NULL);CHKERRQ(ierr); 97694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4,NULL);CHKERRQ(ierr); 97794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5,NULL);CHKERRQ(ierr); 97894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1,NULL);CHKERRQ(ierr); 97994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2,NULL);CHKERRQ(ierr); 98094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3,NULL);CHKERRQ(ierr); 98194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4,NULL);CHKERRQ(ierr); 98294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1,NULL);CHKERRQ(ierr); 98394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2,NULL);CHKERRQ(ierr); 98494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3,NULL);CHKERRQ(ierr); 98594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4,NULL);CHKERRQ(ierr); 98694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5,NULL);CHKERRQ(ierr); 98794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i,NULL);CHKERRQ(ierr); 98894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i,NULL);CHKERRQ(ierr); 98994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i,NULL);CHKERRQ(ierr); 99094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i,NULL);CHKERRQ(ierr); 99194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i,NULL);CHKERRQ(ierr); 99294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i,NULL);CHKERRQ(ierr); 99394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i,NULL);CHKERRQ(ierr); 99494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1,NULL);CHKERRQ(ierr); 99594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2,NULL);CHKERRQ(ierr); 99694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1,NULL);CHKERRQ(ierr); 99794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2,NULL);CHKERRQ(ierr); 99894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3,NULL);CHKERRQ(ierr); 99994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4,NULL);CHKERRQ(ierr); 100094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta,NULL);CHKERRQ(ierr); 100194ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius,NULL);CHKERRQ(ierr); 100294ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius,NULL);CHKERRQ(ierr); 100394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon,NULL);CHKERRQ(ierr); 1004a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 1005a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 1006a7e14dcfSSatish Balay ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 1007a7e14dcfSSatish Balay PetscFunctionReturn(0); 1008a7e14dcfSSatish Balay } 1009a7e14dcfSSatish Balay 1010a7e14dcfSSatish Balay 1011a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 1012a7e14dcfSSatish Balay #undef __FUNCT__ 1013a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_NLS" 1014441846f8SBarry Smith static PetscErrorCode TaoView_NLS(Tao tao, PetscViewer viewer) 1015a7e14dcfSSatish Balay { 1016a7e14dcfSSatish Balay TAO_NLS *nlsP = (TAO_NLS *)tao->data; 1017a7e14dcfSSatish Balay PetscInt nrejects; 1018a7e14dcfSSatish Balay PetscBool isascii; 1019a7e14dcfSSatish Balay PetscErrorCode ierr; 1020a7e14dcfSSatish Balay 1021a7e14dcfSSatish Balay PetscFunctionBegin; 1022a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 1023a7e14dcfSSatish Balay if (isascii) { 1024a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1025a7e14dcfSSatish Balay if (NLS_PC_BFGS == nlsP->pc_type && nlsP->M) { 1026a7e14dcfSSatish Balay ierr = MatLMVMGetRejects(nlsP->M,&nrejects);CHKERRQ(ierr); 1027a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr); 1028a7e14dcfSSatish Balay } 1029a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", nlsP->newt);CHKERRQ(ierr); 1030a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", nlsP->bfgs);CHKERRQ(ierr); 1031a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", nlsP->sgrad);CHKERRQ(ierr); 1032a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", nlsP->grad);CHKERRQ(ierr); 1033a7e14dcfSSatish Balay 1034a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp atol: %D\n", nlsP->ksp_atol);CHKERRQ(ierr); 1035a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %D\n", nlsP->ksp_rtol);CHKERRQ(ierr); 1036a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %D\n", nlsP->ksp_ctol);CHKERRQ(ierr); 1037a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp negc: %D\n", nlsP->ksp_negc);CHKERRQ(ierr); 1038a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %D\n", nlsP->ksp_dtol);CHKERRQ(ierr); 1039a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp iter: %D\n", nlsP->ksp_iter);CHKERRQ(ierr); 1040a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "nls ksp othr: %D\n", nlsP->ksp_othr);CHKERRQ(ierr); 1041a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1042a7e14dcfSSatish Balay } 1043a7e14dcfSSatish Balay PetscFunctionReturn(0); 1044a7e14dcfSSatish Balay } 1045a7e14dcfSSatish Balay 1046a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 10474aa34175SJason Sarich /*MC 10484aa34175SJason Sarich TAONLS - Newton's method with linesearch for unconstrained minimization. 10494aa34175SJason Sarich At each iteration, the Newton line search method solves the symmetric 10504aa34175SJason Sarich system of equations to obtain the step diretion dk: 10514aa34175SJason Sarich Hk dk = -gk 10524aa34175SJason Sarich a More-Thuente line search is applied on the direction dk to approximately 10534aa34175SJason Sarich solve 10544aa34175SJason Sarich min_t f(xk + t d_k) 10554aa34175SJason Sarich 10564aa34175SJason Sarich Options Database Keys: 10574aa34175SJason Sarich + -tao_nls_ksp_type - "cg","nash","stcg","gltr","petsc" 10584aa34175SJason Sarich . -tao_nls_pc_type - "none","ahess","bfgs","petsc" 10594aa34175SJason Sarich . -tao_nls_bfgs_scale_type - "ahess","phess","bfgs" 10604aa34175SJason Sarich . -tao_nls_init_type - "constant","direction","interpolation" 10614aa34175SJason Sarich . -tao_nls_update_type - "step","direction","interpolation" 10624aa34175SJason Sarich . -tao_nls_sval - perturbation starting value 10634aa34175SJason Sarich . -tao_nls_imin - minimum initial perturbation 10644aa34175SJason Sarich . -tao_nls_imax - maximum initial perturbation 10654aa34175SJason Sarich . -tao_nls_pmin - minimum perturbation 10664aa34175SJason Sarich . -tao_nls_pmax - maximum perturbation 10674aa34175SJason Sarich . -tao_nls_pgfac - growth factor 10684aa34175SJason Sarich . -tao_nls_psfac - shrink factor 10694aa34175SJason Sarich . -tao_nls_imfac - initial merit factor 10704aa34175SJason Sarich . -tao_nls_pmgfac - merit growth factor 10714aa34175SJason Sarich . -tao_nls_pmsfac - merit shrink factor 10724aa34175SJason Sarich . -tao_nls_eta1 - poor steplength; reduce radius 10734aa34175SJason Sarich . -tao_nls_eta2 - reasonable steplength; leave radius 10744aa34175SJason Sarich . -tao_nls_eta3 - good steplength; increase readius 10754aa34175SJason Sarich . -tao_nls_eta4 - excellent steplength; greatly increase radius 10764aa34175SJason Sarich . -tao_nls_alpha1 - alpha1 reduction 10774aa34175SJason Sarich . -tao_nls_alpha2 - alpha2 reduction 10784aa34175SJason Sarich . -tao_nls_alpha3 - alpha3 reduction 10794aa34175SJason Sarich . -tao_nls_alpha4 - alpha4 reduction 10804aa34175SJason Sarich . -tao_nls_alpha - alpha5 reduction 10814aa34175SJason Sarich . -tao_nls_mu1 - mu1 interpolation update 10824aa34175SJason Sarich . -tao_nls_mu2 - mu2 interpolation update 10834aa34175SJason Sarich . -tao_nls_gamma1 - gamma1 interpolation update 10844aa34175SJason Sarich . -tao_nls_gamma2 - gamma2 interpolation update 10854aa34175SJason Sarich . -tao_nls_gamma3 - gamma3 interpolation update 10864aa34175SJason Sarich . -tao_nls_gamma4 - gamma4 interpolation update 10874aa34175SJason Sarich . -tao_nls_theta - theta interpolation update 10884aa34175SJason Sarich . -tao_nls_omega1 - omega1 step update 10894aa34175SJason Sarich . -tao_nls_omega2 - omega2 step update 10904aa34175SJason Sarich . -tao_nls_omega3 - omega3 step update 10914aa34175SJason Sarich . -tao_nls_omega4 - omega4 step update 10924aa34175SJason Sarich . -tao_nls_omega5 - omega5 step update 10931522df2eSJason Sarich . -tao_nls_mu1_i - mu1 interpolation init factor 10941522df2eSJason Sarich . -tao_nls_mu2_i - mu2 interpolation init factor 10951522df2eSJason Sarich . -tao_nls_gamma1_i - gamma1 interpolation init factor 10961522df2eSJason Sarich . -tao_nls_gamma2_i - gamma2 interpolation init factor 10971522df2eSJason Sarich . -tao_nls_gamma3_i - gamma3 interpolation init factor 10981522df2eSJason Sarich . -tao_nls_gamma4_i - gamma4 interpolation init factor 10991522df2eSJason Sarich - -tao_nls_theta_i - theta interpolation init factor 11001eb8069cSJason Sarich 11011eb8069cSJason Sarich Level: beginner 11021522df2eSJason Sarich M*/ 11034aa34175SJason Sarich 1104a7e14dcfSSatish Balay #undef __FUNCT__ 1105a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_NLS" 1106728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NLS(Tao tao) 1107a7e14dcfSSatish Balay { 1108a7e14dcfSSatish Balay TAO_NLS *nlsP; 11098caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 1110a7e14dcfSSatish Balay PetscErrorCode ierr; 1111a7e14dcfSSatish Balay 1112a7e14dcfSSatish Balay PetscFunctionBegin; 11133c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&nlsP);CHKERRQ(ierr); 1114a7e14dcfSSatish Balay 1115a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_NLS; 1116a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_NLS; 1117a7e14dcfSSatish Balay tao->ops->view = TaoView_NLS; 1118a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_NLS; 1119a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_NLS; 1120a7e14dcfSSatish Balay 1121a7e14dcfSSatish Balay tao->max_it = 50; 11226f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 11236f4723b1SBarry Smith tao->fatol = 1e-5; 11246f4723b1SBarry Smith tao->frtol = 1e-5; 11256f4723b1SBarry Smith #else 1126a7e14dcfSSatish Balay tao->fatol = 1e-10; 1127a7e14dcfSSatish Balay tao->frtol = 1e-10; 11286f4723b1SBarry Smith #endif 1129a7e14dcfSSatish Balay tao->data = (void*)nlsP; 1130a7e14dcfSSatish Balay tao->trust0 = 100.0; 1131a7e14dcfSSatish Balay 1132a7e14dcfSSatish Balay nlsP->sval = 0.0; 1133a7e14dcfSSatish Balay nlsP->imin = 1.0e-4; 1134a7e14dcfSSatish Balay nlsP->imax = 1.0e+2; 1135a7e14dcfSSatish Balay nlsP->imfac = 1.0e-1; 1136a7e14dcfSSatish Balay 1137a7e14dcfSSatish Balay nlsP->pmin = 1.0e-12; 1138a7e14dcfSSatish Balay nlsP->pmax = 1.0e+2; 1139a7e14dcfSSatish Balay nlsP->pgfac = 1.0e+1; 1140a7e14dcfSSatish Balay nlsP->psfac = 4.0e-1; 1141a7e14dcfSSatish Balay nlsP->pmgfac = 1.0e-1; 1142a7e14dcfSSatish Balay nlsP->pmsfac = 1.0e-1; 1143a7e14dcfSSatish Balay 1144a7e14dcfSSatish Balay /* Default values for trust-region radius update based on steplength */ 1145a7e14dcfSSatish Balay nlsP->nu1 = 0.25; 1146a7e14dcfSSatish Balay nlsP->nu2 = 0.50; 1147a7e14dcfSSatish Balay nlsP->nu3 = 1.00; 1148a7e14dcfSSatish Balay nlsP->nu4 = 1.25; 1149a7e14dcfSSatish Balay 1150a7e14dcfSSatish Balay nlsP->omega1 = 0.25; 1151a7e14dcfSSatish Balay nlsP->omega2 = 0.50; 1152a7e14dcfSSatish Balay nlsP->omega3 = 1.00; 1153a7e14dcfSSatish Balay nlsP->omega4 = 2.00; 1154a7e14dcfSSatish Balay nlsP->omega5 = 4.00; 1155a7e14dcfSSatish Balay 1156a7e14dcfSSatish Balay /* Default values for trust-region radius update based on reduction */ 1157a7e14dcfSSatish Balay nlsP->eta1 = 1.0e-4; 1158a7e14dcfSSatish Balay nlsP->eta2 = 0.25; 1159a7e14dcfSSatish Balay nlsP->eta3 = 0.50; 1160a7e14dcfSSatish Balay nlsP->eta4 = 0.90; 1161a7e14dcfSSatish Balay 1162a7e14dcfSSatish Balay nlsP->alpha1 = 0.25; 1163a7e14dcfSSatish Balay nlsP->alpha2 = 0.50; 1164a7e14dcfSSatish Balay nlsP->alpha3 = 1.00; 1165a7e14dcfSSatish Balay nlsP->alpha4 = 2.00; 1166a7e14dcfSSatish Balay nlsP->alpha5 = 4.00; 1167a7e14dcfSSatish Balay 1168a7e14dcfSSatish Balay /* Default values for trust-region radius update based on interpolation */ 1169a7e14dcfSSatish Balay nlsP->mu1 = 0.10; 1170a7e14dcfSSatish Balay nlsP->mu2 = 0.50; 1171a7e14dcfSSatish Balay 1172a7e14dcfSSatish Balay nlsP->gamma1 = 0.25; 1173a7e14dcfSSatish Balay nlsP->gamma2 = 0.50; 1174a7e14dcfSSatish Balay nlsP->gamma3 = 2.00; 1175a7e14dcfSSatish Balay nlsP->gamma4 = 4.00; 1176a7e14dcfSSatish Balay 1177a7e14dcfSSatish Balay nlsP->theta = 0.05; 1178a7e14dcfSSatish Balay 1179a7e14dcfSSatish Balay /* Default values for trust region initialization based on interpolation */ 1180a7e14dcfSSatish Balay nlsP->mu1_i = 0.35; 1181a7e14dcfSSatish Balay nlsP->mu2_i = 0.50; 1182a7e14dcfSSatish Balay 1183a7e14dcfSSatish Balay nlsP->gamma1_i = 0.0625; 1184a7e14dcfSSatish Balay nlsP->gamma2_i = 0.5; 1185a7e14dcfSSatish Balay nlsP->gamma3_i = 2.0; 1186a7e14dcfSSatish Balay nlsP->gamma4_i = 5.0; 1187a7e14dcfSSatish Balay 1188a7e14dcfSSatish Balay nlsP->theta_i = 0.25; 1189a7e14dcfSSatish Balay 1190a7e14dcfSSatish Balay /* Remaining parameters */ 1191a7e14dcfSSatish Balay nlsP->min_radius = 1.0e-10; 1192a7e14dcfSSatish Balay nlsP->max_radius = 1.0e10; 1193a7e14dcfSSatish Balay nlsP->epsilon = 1.0e-6; 1194a7e14dcfSSatish Balay 1195a7e14dcfSSatish Balay nlsP->ksp_type = NLS_KSP_STCG; 1196a7e14dcfSSatish Balay nlsP->pc_type = NLS_PC_BFGS; 1197a7e14dcfSSatish Balay nlsP->bfgs_scale_type = BFGS_SCALE_PHESS; 1198a7e14dcfSSatish Balay nlsP->init_type = NLS_INIT_INTERPOLATION; 1199a7e14dcfSSatish Balay nlsP->update_type = NLS_UPDATE_STEP; 1200a7e14dcfSSatish Balay 1201a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 1202a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 1203441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 1204a7e14dcfSSatish Balay 1205a7e14dcfSSatish Balay /* Set linear solver to default for symmetric matrices */ 1206a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 1207a7e14dcfSSatish Balay PetscFunctionReturn(0); 1208a7e14dcfSSatish Balay } 1209a7e14dcfSSatish Balay 1210a7e14dcfSSatish Balay #undef __FUNCT__ 1211a7e14dcfSSatish Balay #define __FUNCT__ "MatLMVMSolveShell" 1212a7e14dcfSSatish Balay static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x) 1213a7e14dcfSSatish Balay { 1214a7e14dcfSSatish Balay PetscErrorCode ierr; 1215a7e14dcfSSatish Balay Mat M; 121687f595a5SBarry Smith 1217a7e14dcfSSatish Balay PetscFunctionBegin; 1218a7e14dcfSSatish Balay PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1219a7e14dcfSSatish Balay PetscValidHeaderSpecific(b,VEC_CLASSID,2); 1220a7e14dcfSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,3); 1221a7e14dcfSSatish Balay ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr); 1222a7e14dcfSSatish Balay ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr); 1223a7e14dcfSSatish Balay PetscFunctionReturn(0); 1224a7e14dcfSSatish Balay } 1225