1 #include <petsc-private/linesearchimpl.h> 2 #include <petscsnes.h> 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "SNESLineSearchApply_CP" 6 static PetscErrorCode SNESLineSearchApply_CP(SNESLineSearch linesearch) 7 { 8 PetscBool changed_y, changed_w, domainerror; 9 PetscErrorCode ierr; 10 Vec X, Y, F, W; 11 SNES snes; 12 PetscReal xnorm, ynorm, gnorm, steptol, atol, rtol, ltol, maxstep; 13 14 PetscReal lambda, lambda_old, lambda_update, delLambda; 15 PetscScalar fty, fty_init, fty_old, fty_mid1, fty_mid2, s; 16 PetscInt i, max_its; 17 18 PetscViewer monitor; 19 20 PetscFunctionBegin; 21 ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, PETSC_NULL);CHKERRQ(ierr); 22 ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 23 ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 24 ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); 25 ierr = SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, &rtol, &atol, <ol, &max_its);CHKERRQ(ierr); 26 ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 27 ierr = SNESLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr); 28 29 /* precheck */ 30 ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr); 31 lambda_old = 0.0; 32 ierr = VecDot(F, Y, &fty_old);CHKERRQ(ierr); 33 fty_init = fty_old; 34 35 for (i = 0; i < max_its; i++) { 36 37 /* compute the norm at lambda */ 38 ierr = VecCopy(X, W);CHKERRQ(ierr); 39 ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); 40 if (linesearch->ops->viproject) { 41 ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 42 } 43 ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 44 45 ierr = VecDot(F, Y, &fty);CHKERRQ(ierr); 46 47 delLambda = lambda - lambda_old; 48 49 /* check for convergence */ 50 if (PetscAbsReal(delLambda) < steptol*lambda) break; 51 if (PetscAbsScalar(fty) / PetscAbsScalar(fty_init) < rtol) break; 52 if (PetscAbsScalar(fty) < atol && i > 0) break; 53 if (monitor) { 54 ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 55 ierr = PetscViewerASCIIPrintf(monitor," Line search: lambdas = [%g, %g], ftys = [%g, %g]\n",(double)lambda, (double)lambda_old, (double)PetscRealPart(fty), (double)PetscRealPart(fty_old));CHKERRQ(ierr); 56 ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 57 } 58 59 /* compute the search direction */ 60 if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) { 61 s = (fty - fty_old) / delLambda; 62 } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 63 ierr = VecCopy(X, W);CHKERRQ(ierr); 64 ierr = VecAXPY(W, -0.5*(lambda + lambda_old), Y);CHKERRQ(ierr); 65 if (linesearch->ops->viproject) { 66 ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 67 } 68 ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 69 ierr = VecDot(F, Y, &fty_mid1);CHKERRQ(ierr); 70 s = (3.*fty - 4.*fty_mid1 + fty_old) / delLambda; 71 } else { 72 ierr = VecCopy(X, W);CHKERRQ(ierr); 73 ierr = VecAXPY(W, -0.5*(lambda + lambda_old), Y);CHKERRQ(ierr); 74 if (linesearch->ops->viproject) { 75 ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 76 } 77 ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 78 ierr = VecDot(F, Y, &fty_mid1);CHKERRQ(ierr); 79 ierr = VecCopy(X, W);CHKERRQ(ierr); 80 ierr = VecAXPY(W, -(lambda + 0.5*(lambda - lambda_old)), Y);CHKERRQ(ierr); 81 if (linesearch->ops->viproject) { 82 ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 83 } 84 ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 85 ierr = VecDot(F, Y, &fty_mid2);CHKERRQ(ierr); 86 s = (2.*fty_mid2 + 3.*fty - 6.*fty_mid1 + fty_old) / (3.*delLambda); 87 } 88 /* if the solve is going in the wrong direction, fix it */ 89 if (PetscRealPart(s) > 0.) s = -s; 90 lambda_update = lambda - PetscRealPart(fty / s); 91 92 /* switch directions if we stepped out of bounds */ 93 if (lambda_update < steptol) { 94 lambda_update = lambda + PetscRealPart(fty / s); 95 } 96 97 if (PetscIsInfOrNanScalar(lambda_update)) break; 98 if (lambda_update > maxstep) { 99 break; 100 } 101 102 /* compute the new state of the line search */ 103 lambda_old = lambda; 104 lambda = lambda_update; 105 fty_old = fty; 106 } 107 /* construct the solution */ 108 ierr = VecCopy(X, W);CHKERRQ(ierr); 109 ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); 110 if (linesearch->ops->viproject) { 111 ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 112 } 113 /* postcheck */ 114 ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 115 if (changed_y) { 116 ierr = VecAXPY(X, -lambda, Y);CHKERRQ(ierr); 117 if (linesearch->ops->viproject) { 118 ierr = (*linesearch->ops->viproject)(snes, X);CHKERRQ(ierr); 119 } 120 } else { 121 ierr = VecCopy(W, X);CHKERRQ(ierr); 122 } 123 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 124 ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 125 if (domainerror) { 126 ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 127 PetscFunctionReturn(0); 128 } 129 130 ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr); 131 ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 132 133 ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 134 135 if (monitor) { 136 ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 137 ierr = PetscViewerASCIIPrintf(monitor," Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm);CHKERRQ(ierr); 138 ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 139 } 140 if (lambda <= steptol) { 141 ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 142 } 143 PetscFunctionReturn(0); 144 } 145 146 #undef __FUNCT__ 147 #define __FUNCT__ "SNESLineSearchCreate_CP" 148 /*MC 149 SNESLINESEARCHCP - Critical point line search. This line search assumes that there exists some 150 artificial G(x) for which the SNESFunction F(x) = grad G(x). Therefore, this line search seeks 151 to find roots of dot(F, Y) via a secant method. 152 153 Options Database Keys: 154 + -snes_linesearch_minlambda - the minimum acceptable lambda 155 . -snes_linesearch_damping - initial trial step length 156 - -snes_linesearch_max_it - the maximum number of secant steps performed. 157 158 Notes: 159 This method is the preferred line search for SNESQN and SNESNCG. 160 161 Level: advanced 162 163 .keywords: SNES, SNESLineSearch, damping 164 165 .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 166 M*/ 167 PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch) 168 { 169 PetscFunctionBegin; 170 linesearch->ops->apply = SNESLineSearchApply_CP; 171 linesearch->ops->destroy = PETSC_NULL; 172 linesearch->ops->setfromoptions = PETSC_NULL; 173 linesearch->ops->reset = PETSC_NULL; 174 linesearch->ops->view = PETSC_NULL; 175 linesearch->ops->setup = PETSC_NULL; 176 linesearch->order = SNES_LINESEARCH_ORDER_LINEAR; 177 178 linesearch->max_its = 1; 179 PetscFunctionReturn(0); 180 } 181