1fef7b6d8SPeter Brune #include <../src/snes/impls/ncg/snesncgimpl.h> 2fef7b6d8SPeter Brune 3fef7b6d8SPeter Brune #undef __FUNCT__ 4fef7b6d8SPeter Brune #define __FUNCT__ "SNESReset_NCG" 5fef7b6d8SPeter Brune PetscErrorCode SNESReset_NCG(SNES snes) 6fef7b6d8SPeter Brune { 7fef7b6d8SPeter Brune PetscErrorCode ierr; 8fef7b6d8SPeter Brune 9fef7b6d8SPeter Brune PetscFunctionBegin; 10fef7b6d8SPeter Brune if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 11fef7b6d8SPeter Brune PetscFunctionReturn(0); 12fef7b6d8SPeter Brune } 13fef7b6d8SPeter Brune 14fef7b6d8SPeter Brune /* 15fef7b6d8SPeter Brune SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG(). 16fef7b6d8SPeter Brune 17fef7b6d8SPeter Brune Input Parameter: 18fef7b6d8SPeter Brune . snes - the SNES context 19fef7b6d8SPeter Brune 20fef7b6d8SPeter Brune Application Interface Routine: SNESDestroy() 21fef7b6d8SPeter Brune */ 22fef7b6d8SPeter Brune #undef __FUNCT__ 23fef7b6d8SPeter Brune #define __FUNCT__ "SNESDestroy_NCG" 24fef7b6d8SPeter Brune PetscErrorCode SNESDestroy_NCG(SNES snes) 25fef7b6d8SPeter Brune { 26fef7b6d8SPeter Brune PetscErrorCode ierr; 27fef7b6d8SPeter Brune SNES_NCG *neP = (SNES_NCG*)snes->data; 28fef7b6d8SPeter Brune 29fef7b6d8SPeter Brune PetscFunctionBegin; 30fef7b6d8SPeter Brune ierr = SNESReset_NCG(snes);CHKERRQ(ierr); 31fef7b6d8SPeter Brune ierr = PetscViewerDestroy(&neP->monitor);CHKERRQ(ierr); 32fef7b6d8SPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 33fef7b6d8SPeter Brune PetscFunctionReturn(0); 34fef7b6d8SPeter Brune } 35fef7b6d8SPeter Brune 36fef7b6d8SPeter Brune /* 37fef7b6d8SPeter Brune SNESSetUp_NCG - Sets up the internal data structures for the later use 38fef7b6d8SPeter Brune of the SNESNCG nonlinear solver. 39fef7b6d8SPeter Brune 40fef7b6d8SPeter Brune Input Parameters: 41fef7b6d8SPeter Brune + snes - the SNES context 42fef7b6d8SPeter Brune - x - the solution vector 43fef7b6d8SPeter Brune 44fef7b6d8SPeter Brune Application Interface Routine: SNESSetUp() 45fef7b6d8SPeter Brune */ 46fef7b6d8SPeter Brune #undef __FUNCT__ 47fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetUp_NCG" 48fef7b6d8SPeter Brune PetscErrorCode SNESSetUp_NCG(SNES snes) 49fef7b6d8SPeter Brune { 50fef7b6d8SPeter Brune PetscErrorCode ierr; 51fef7b6d8SPeter Brune 52fef7b6d8SPeter Brune PetscFunctionBegin; 53fef7b6d8SPeter Brune ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr); 54fef7b6d8SPeter Brune PetscFunctionReturn(0); 55fef7b6d8SPeter Brune } 56fef7b6d8SPeter Brune 57fef7b6d8SPeter Brune EXTERN_C_BEGIN 58fef7b6d8SPeter Brune #undef __FUNCT__ 59fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor_NCG" 60fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchSetMonitor_NCG(SNES snes,PetscBool flg) 61fef7b6d8SPeter Brune { 62fef7b6d8SPeter Brune SNES_NCG *neP = (SNES_NCG*)snes->data; 63fef7b6d8SPeter Brune PetscErrorCode ierr; 64fef7b6d8SPeter Brune 65fef7b6d8SPeter Brune PetscFunctionBegin; 66fef7b6d8SPeter Brune if (flg && !neP->monitor) { 67fef7b6d8SPeter Brune ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,"stdout",&neP->monitor);CHKERRQ(ierr); 68fef7b6d8SPeter Brune } else if (!flg && neP->monitor) { 69fef7b6d8SPeter Brune ierr = PetscViewerDestroy(&neP->monitor);CHKERRQ(ierr); 70fef7b6d8SPeter Brune } 71fef7b6d8SPeter Brune PetscFunctionReturn(0); 72fef7b6d8SPeter Brune } 73fef7b6d8SPeter Brune EXTERN_C_END 74fef7b6d8SPeter Brune 75fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchNoNorms_NCG(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool*); 76fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchNo_NCG(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool*); 77fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchQuadratic_NCG(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool*); 78fef7b6d8SPeter Brune /* 79fef7b6d8SPeter Brune SNESSetFromOptions_NCG - Sets various parameters for the SNESLS method. 80fef7b6d8SPeter Brune 81fef7b6d8SPeter Brune Input Parameter: 82fef7b6d8SPeter Brune . snes - the SNES context 83fef7b6d8SPeter Brune 84fef7b6d8SPeter Brune Application Interface Routine: SNESSetFromOptions() 85fef7b6d8SPeter Brune */ 86fef7b6d8SPeter Brune #undef __FUNCT__ 87fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NCG" 88fef7b6d8SPeter Brune static PetscErrorCode SNESSetFromOptions_NCG(SNES snes) 89fef7b6d8SPeter Brune { 90fef7b6d8SPeter Brune SNES_NCG *ls = (SNES_NCG *)snes->data; 91fef7b6d8SPeter Brune SNESLineSearchType indx; 92fef7b6d8SPeter Brune PetscBool flg,set; 93fef7b6d8SPeter Brune PetscErrorCode ierr; 94fef7b6d8SPeter Brune 95fef7b6d8SPeter Brune PetscFunctionBegin; 96fef7b6d8SPeter Brune ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr); 97fef7b6d8SPeter Brune ierr = PetscOptionsEnum("-snes_ls","NCG line search type","SNESLineSearchSet",SNESLineSearchTypes,(PetscEnum)ls->type,(PetscEnum*)&indx,&flg);CHKERRQ(ierr); 98fef7b6d8SPeter Brune if (flg) { 99fef7b6d8SPeter Brune ls->type = indx; 100fef7b6d8SPeter Brune switch (indx) { 101fef7b6d8SPeter Brune case SNES_LS_BASIC: 102fef7b6d8SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo_NCG,PETSC_NULL);CHKERRQ(ierr); 103fef7b6d8SPeter Brune break; 104fef7b6d8SPeter Brune case SNES_LS_BASIC_NONORMS: 105fef7b6d8SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_NCG,PETSC_NULL);CHKERRQ(ierr); 106fef7b6d8SPeter Brune break; 107fef7b6d8SPeter Brune case SNES_LS_QUADRATIC: 108fef7b6d8SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_NCG,PETSC_NULL);CHKERRQ(ierr); 109fef7b6d8SPeter Brune break; 110fef7b6d8SPeter Brune case SNES_LS_CUBIC: 111fef7b6d8SPeter Brune SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_SUP,"No support for cubic search"); 112fef7b6d8SPeter Brune break; 113fef7b6d8SPeter Brune } 114fef7b6d8SPeter Brune } 115fef7b6d8SPeter Brune ierr = PetscOptionsReal("-snes_ls_damping","Damping parameter","SNES",ls->damping,&ls->damping,&flg);CHKERRQ(ierr); 116fef7b6d8SPeter Brune ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",ls->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 117fef7b6d8SPeter Brune if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 118fef7b6d8SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 119fef7b6d8SPeter Brune PetscFunctionReturn(0); 120fef7b6d8SPeter Brune } 121fef7b6d8SPeter Brune 122fef7b6d8SPeter Brune /* 123fef7b6d8SPeter Brune SNESView_NCG - Prints info from the SNESNCG data structure. 124fef7b6d8SPeter Brune 125fef7b6d8SPeter Brune Input Parameters: 126fef7b6d8SPeter Brune + SNES - the SNES context 127fef7b6d8SPeter Brune - viewer - visualization context 128fef7b6d8SPeter Brune 129fef7b6d8SPeter Brune Application Interface Routine: SNESView() 130fef7b6d8SPeter Brune */ 131fef7b6d8SPeter Brune #undef __FUNCT__ 132fef7b6d8SPeter Brune #define __FUNCT__ "SNESView_NCG" 133fef7b6d8SPeter Brune static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) 134fef7b6d8SPeter Brune { 135fef7b6d8SPeter Brune SNES_NCG *ls = (SNES_NCG *)snes->data; 136fef7b6d8SPeter Brune PetscBool iascii; 137fef7b6d8SPeter Brune PetscErrorCode ierr; 138fef7b6d8SPeter Brune 139fef7b6d8SPeter Brune PetscFunctionBegin; 140fef7b6d8SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 141fef7b6d8SPeter Brune if (iascii) { 142fef7b6d8SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," line search type variant: %s\n", SNESLineSearchTypeName(ls->type));CHKERRQ(ierr); 143fef7b6d8SPeter Brune } 144fef7b6d8SPeter Brune PetscFunctionReturn(0); 145fef7b6d8SPeter Brune } 146fef7b6d8SPeter Brune 147fef7b6d8SPeter Brune #undef __FUNCT__ 148fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchNo_NCG" 149fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchNo_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag) 150fef7b6d8SPeter Brune { 151fef7b6d8SPeter Brune PetscErrorCode ierr; 152fef7b6d8SPeter Brune SNES_NCG *neP = (SNES_NCG *) snes->data; 153fef7b6d8SPeter Brune 154fef7b6d8SPeter Brune PetscFunctionBegin; 155fef7b6d8SPeter Brune ierr = VecAXPY(X, neP->damping, Y);CHKERRQ(ierr); 156fef7b6d8SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 157fef7b6d8SPeter Brune ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr); 158fef7b6d8SPeter Brune if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated norm"); 159fef7b6d8SPeter Brune PetscFunctionReturn(0); 160fef7b6d8SPeter Brune } 161fef7b6d8SPeter Brune 162fef7b6d8SPeter Brune #undef __FUNCT__ 163fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchNoNorms_NCG" 164fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchNoNorms_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag) 165fef7b6d8SPeter Brune { 166fef7b6d8SPeter Brune PetscErrorCode ierr; 167fef7b6d8SPeter Brune SNES_NCG *neP = (SNES_NCG *) snes->data; 168fef7b6d8SPeter Brune 169fef7b6d8SPeter Brune PetscFunctionBegin; 170fef7b6d8SPeter Brune ierr = VecAXPY(X, neP->damping, Y);CHKERRQ(ierr); 171fef7b6d8SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 172fef7b6d8SPeter Brune PetscFunctionReturn(0); 173fef7b6d8SPeter Brune } 174fef7b6d8SPeter Brune 175fef7b6d8SPeter Brune #undef __FUNCT__ 176fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchQuadratic_NCG" 177fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchQuadratic_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec G,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag) 178fef7b6d8SPeter Brune { 179fef7b6d8SPeter Brune PetscInt i; 180fef7b6d8SPeter Brune PetscReal alphas[3] = {0.0, 0.5, 1.0}; 181fef7b6d8SPeter Brune PetscReal norms[3]; 182fef7b6d8SPeter Brune PetscReal alpha,a,b; 183fef7b6d8SPeter Brune PetscErrorCode ierr; 184fef7b6d8SPeter Brune SNES_NCG *neP = (SNES_NCG *) snes->data; 185fef7b6d8SPeter Brune 186fef7b6d8SPeter Brune PetscFunctionBegin; 187fef7b6d8SPeter Brune norms[0] = fnorm; 188fef7b6d8SPeter Brune for(i=1; i < 3; ++i) { 189fef7b6d8SPeter Brune ierr = VecWAXPY(W, alphas[i], Y, X);CHKERRQ(ierr); /* W = X^n - \alpha Y */ 190fef7b6d8SPeter Brune ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 191fef7b6d8SPeter Brune ierr = VecNorm(F, NORM_2, &norms[i]);CHKERRQ(ierr); 192fef7b6d8SPeter Brune } 193fef7b6d8SPeter Brune for(i = 0; i < 3; ++i) { 194fef7b6d8SPeter Brune norms[i] = PetscSqr(norms[i]); 195fef7b6d8SPeter Brune } 196fef7b6d8SPeter Brune /* Fit a quadratic: 197fef7b6d8SPeter Brune If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2} 198fef7b6d8SPeter Brune a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1) 199fef7b6d8SPeter Brune b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1) 200fef7b6d8SPeter Brune c = y_0 201fef7b6d8SPeter Brune x_min = -b/2a 202fef7b6d8SPeter Brune 203fef7b6d8SPeter Brune If we let x_{0,1,2} = 0, 0.5, 1.0 204fef7b6d8SPeter Brune a = 2 y_2 - 4 y_1 + 2 y_0 205fef7b6d8SPeter Brune b = -y_2 + 4 y_1 - 3 y_0 206fef7b6d8SPeter Brune c = y_0 207fef7b6d8SPeter Brune */ 208fef7b6d8SPeter Brune a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1])); 209fef7b6d8SPeter Brune b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]); 210fef7b6d8SPeter Brune /* Check for positive a (concave up) */ 211fef7b6d8SPeter Brune if (a >= 0.0) { 212fef7b6d8SPeter Brune alpha = -b/(2.0*a); 213fef7b6d8SPeter Brune alpha = PetscMin(alpha, alphas[2]); 214fef7b6d8SPeter Brune alpha = PetscMax(alpha, alphas[0]); 215fef7b6d8SPeter Brune } else { 216fef7b6d8SPeter Brune alpha = 1.0; 217fef7b6d8SPeter Brune } 218fef7b6d8SPeter Brune if (neP->monitor) { 219fef7b6d8SPeter Brune ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 220fef7b6d8SPeter Brune ierr = PetscViewerASCIIPrintf(neP->monitor," Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", sqrt(norms[0]),sqrt(norms[1]),sqrt(norms[2]),alpha);CHKERRQ(ierr); 221fef7b6d8SPeter Brune ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 222fef7b6d8SPeter Brune } 223fef7b6d8SPeter Brune ierr = VecAXPY(X, alpha, Y);CHKERRQ(ierr); 224fef7b6d8SPeter Brune if (alpha != 1.0) { 225fef7b6d8SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 226fef7b6d8SPeter Brune ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr); 227fef7b6d8SPeter Brune } else { 228fef7b6d8SPeter Brune *gnorm = PetscSqrtReal(norms[2]); 229fef7b6d8SPeter Brune } 230fef7b6d8SPeter Brune if (alpha == 0.0) *flag = PETSC_FALSE; 231fef7b6d8SPeter Brune else *flag = PETSC_TRUE; 232fef7b6d8SPeter Brune PetscFunctionReturn(0); 233fef7b6d8SPeter Brune } 234fef7b6d8SPeter Brune 235fef7b6d8SPeter Brune /* 236fef7b6d8SPeter Brune SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method. 237fef7b6d8SPeter Brune 238fef7b6d8SPeter Brune Input Parameters: 239fef7b6d8SPeter Brune . snes - the SNES context 240fef7b6d8SPeter Brune 241fef7b6d8SPeter Brune Output Parameter: 242fef7b6d8SPeter Brune . outits - number of iterations until termination 243fef7b6d8SPeter Brune 244fef7b6d8SPeter Brune Application Interface Routine: SNESSolve() 245fef7b6d8SPeter Brune */ 246fef7b6d8SPeter Brune #undef __FUNCT__ 247fef7b6d8SPeter Brune #define __FUNCT__ "SNESSolve_NCG" 248fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes) 249fef7b6d8SPeter Brune { 250fef7b6d8SPeter Brune SNES_NCG *neP = (SNES_NCG *) snes->data; 251302dbbaeSPeter Brune Vec X, dX, lX, F, W, B; 252fef7b6d8SPeter Brune PetscReal fnorm, dummyNorm; 253fef7b6d8SPeter Brune PetscScalar dXdot, dXdot_old; 254fef7b6d8SPeter Brune PetscInt maxits, i; 255fef7b6d8SPeter Brune PetscErrorCode ierr; 256fef7b6d8SPeter Brune SNESConvergedReason reason; 257fef7b6d8SPeter Brune PetscBool lsSuccess = PETSC_TRUE; 258fef7b6d8SPeter Brune PetscFunctionBegin; 259fef7b6d8SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 260fef7b6d8SPeter Brune 261fef7b6d8SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 262fef7b6d8SPeter Brune X = snes->vec_sol; /* X^n */ 263169e2be8SPeter Brune dX = snes->work[1]; /* the steepest direction */ 264169e2be8SPeter Brune lX = snes->vec_sol_update; /* the conjugate direction */ 265fef7b6d8SPeter Brune F = snes->vec_func; /* residual vector */ 266fef7b6d8SPeter Brune W = snes->work[0]; /* work vector for the line search */ 267302dbbaeSPeter Brune B = snes->vec_rhs; /* the right hand side */ 268fef7b6d8SPeter Brune 269fef7b6d8SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 270fef7b6d8SPeter Brune snes->iter = 0; 271fef7b6d8SPeter Brune snes->norm = 0.; 272fef7b6d8SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 273fef7b6d8SPeter Brune 274fef7b6d8SPeter Brune /* compute the initial function and preconditioned update delX */ 275fef7b6d8SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 276fef7b6d8SPeter Brune if (snes->domainerror) { 277fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 278fef7b6d8SPeter Brune PetscFunctionReturn(0); 279fef7b6d8SPeter Brune } 280fef7b6d8SPeter Brune /* convergence test */ 281fef7b6d8SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 282fef7b6d8SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 283fef7b6d8SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 284fef7b6d8SPeter Brune snes->norm = fnorm; 285fef7b6d8SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 286fef7b6d8SPeter Brune SNESLogConvHistory(snes,fnorm,0); 287fef7b6d8SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 288fef7b6d8SPeter Brune 289fef7b6d8SPeter Brune /* set parameter for default relative tolerance convergence test */ 290fef7b6d8SPeter Brune snes->ttol = fnorm*snes->rtol; 291fef7b6d8SPeter Brune /* test convergence */ 292fef7b6d8SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 293fef7b6d8SPeter Brune if (snes->reason) PetscFunctionReturn(0); 294fef7b6d8SPeter Brune 295fef7b6d8SPeter Brune /* Call general purpose update function */ 296fef7b6d8SPeter Brune if (snes->ops->update) { 297fef7b6d8SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 298fef7b6d8SPeter Brune } 299fef7b6d8SPeter Brune 300fef7b6d8SPeter Brune /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 301fef7b6d8SPeter Brune if (!snes->pc) { 302fef7b6d8SPeter Brune ierr = VecCopy(F, lX);CHKERRQ(ierr); 303fef7b6d8SPeter Brune ierr = VecScale(lX, -1.0);CHKERRQ(ierr); 304fef7b6d8SPeter Brune } else { 305fef7b6d8SPeter Brune ierr = VecCopy(X, lX);CHKERRQ(ierr); 306421d9b32SPeter Brune ierr = SNESSolve(snes->pc, 0, lX);CHKERRQ(ierr); 307fef7b6d8SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 308302dbbaeSPeter Brune ierr = VecAXPY(lX, -1.0, X);CHKERRQ(ierr); 309fef7b6d8SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 310fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_INNER; 311fef7b6d8SPeter Brune PetscFunctionReturn(0); 312fef7b6d8SPeter Brune } 313fef7b6d8SPeter Brune } 314fef7b6d8SPeter Brune ierr = VecDot(lX, lX, &dXdot);CHKERRQ(ierr); 315fef7b6d8SPeter Brune for(i = 1; i < maxits; i++) { 316fef7b6d8SPeter Brune lsSuccess = PETSC_TRUE; 317fef7b6d8SPeter Brune ierr = (*neP->LineSearch)(snes, neP->lsP, X, F, lX, fnorm, 0.0, 0, W, &dummyNorm, &fnorm, &lsSuccess);CHKERRQ(ierr); 318fef7b6d8SPeter Brune if (!lsSuccess) { 319fef7b6d8SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 320fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 321fef7b6d8SPeter Brune break; 322fef7b6d8SPeter Brune } 323fef7b6d8SPeter Brune } 324fef7b6d8SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 325fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 326fef7b6d8SPeter Brune break; 327fef7b6d8SPeter Brune } 328fef7b6d8SPeter Brune if (snes->domainerror) { 329fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 330fef7b6d8SPeter Brune PetscFunctionReturn(0); 331fef7b6d8SPeter Brune } 332302dbbaeSPeter Brune 333fef7b6d8SPeter Brune /* Monitor convergence */ 334fef7b6d8SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 335169e2be8SPeter Brune snes->iter = i; 336fef7b6d8SPeter Brune snes->norm = fnorm; 337fef7b6d8SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 338fef7b6d8SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 339fef7b6d8SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 340302dbbaeSPeter Brune 341fef7b6d8SPeter Brune /* Test for convergence */ 342fef7b6d8SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 343fef7b6d8SPeter Brune if (snes->reason) break; 344302dbbaeSPeter Brune 345302dbbaeSPeter Brune /* Call general purpose update function */ 346302dbbaeSPeter Brune if (snes->ops->update) { 347302dbbaeSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 348302dbbaeSPeter Brune } 349*ee78dd50SPeter Brune if (!snes->pc) { 350302dbbaeSPeter Brune ierr = VecCopy(F,dX);CHKERRQ(ierr); 351302dbbaeSPeter Brune ierr = VecScale(dX,-1.0);CHKERRQ(ierr); 352302dbbaeSPeter Brune } else { 353302dbbaeSPeter Brune ierr = VecCopy(X,dX);CHKERRQ(ierr); 354*ee78dd50SPeter Brune ierr = SNESSolve(snes->pc, snes->vec_rhs, dX);CHKERRQ(ierr); 355302dbbaeSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 356302dbbaeSPeter Brune if (reason < 0) { 357302dbbaeSPeter Brune snes->reason = SNES_DIVERGED_INNER; 358302dbbaeSPeter Brune PetscFunctionReturn(0); 359302dbbaeSPeter Brune } 360302dbbaeSPeter Brune ierr = VecAXPY(dX,-1.0,X);CHKERRQ(ierr); 361302dbbaeSPeter Brune } 362302dbbaeSPeter Brune 363169e2be8SPeter Brune /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 364302dbbaeSPeter Brune dXdot_old = dXdot; 365302dbbaeSPeter Brune ierr = VecDot(dX, dX, &dXdot);CHKERRQ(ierr); 366*ee78dd50SPeter Brune #if 0 367421d9b32SPeter Brune if (1) 368169e2be8SPeter Brune ierr = PetscPrintf(PETSC_COMM_WORLD, "beta %e = %e / %e \n", dXdot / dXdot_old, dXdot, dXdot_old);CHKERRQ(ierr); 369*ee78dd50SPeter Brune #endif 370302dbbaeSPeter Brune ierr = VecAYPX(lX, dXdot / dXdot_old, dX);CHKERRQ(ierr); 371302dbbaeSPeter Brune /* line search for the proper contribution of lX to the solution */ 372fef7b6d8SPeter Brune } 373fef7b6d8SPeter Brune if (i == maxits) { 374fef7b6d8SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 375fef7b6d8SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 376fef7b6d8SPeter Brune } 377fef7b6d8SPeter Brune PetscFunctionReturn(0); 378fef7b6d8SPeter Brune } 379fef7b6d8SPeter Brune 380fef7b6d8SPeter Brune typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,void*,PetscBool *); /* force argument to next function to not be extern C*/ 381fef7b6d8SPeter Brune EXTERN_C_BEGIN 382fef7b6d8SPeter Brune #undef __FUNCT__ 383fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck_NCG" 384fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchSetPreCheck_NCG(SNES snes, FCN1 func, void *checkctx) 385fef7b6d8SPeter Brune { 386fef7b6d8SPeter Brune PetscFunctionBegin; 387fef7b6d8SPeter Brune ((SNES_NCG *)(snes->data))->precheckstep = func; 388fef7b6d8SPeter Brune ((SNES_NCG *)(snes->data))->precheck = checkctx; 389fef7b6d8SPeter Brune PetscFunctionReturn(0); 390fef7b6d8SPeter Brune } 391fef7b6d8SPeter Brune EXTERN_C_END 392fef7b6d8SPeter Brune 393fef7b6d8SPeter Brune typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/ 394fef7b6d8SPeter Brune EXTERN_C_BEGIN 395fef7b6d8SPeter Brune #undef __FUNCT__ 396fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchSet_NCG" 397fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchSet_NCG(SNES snes, FCN2 func, void *lsctx) 398fef7b6d8SPeter Brune { 399fef7b6d8SPeter Brune PetscFunctionBegin; 400fef7b6d8SPeter Brune ((SNES_NCG *)(snes->data))->LineSearch = func; 401fef7b6d8SPeter Brune ((SNES_NCG *)(snes->data))->lsP = lsctx; 402fef7b6d8SPeter Brune PetscFunctionReturn(0); 403fef7b6d8SPeter Brune } 404fef7b6d8SPeter Brune EXTERN_C_END 405fef7b6d8SPeter Brune 406fef7b6d8SPeter Brune typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *); /* force argument to next function to not be extern C*/ 407fef7b6d8SPeter Brune EXTERN_C_BEGIN 408fef7b6d8SPeter Brune #undef __FUNCT__ 409fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck_NCG" 410fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchSetPostCheck_NCG(SNES snes, FCN3 func, void *checkctx) 411fef7b6d8SPeter Brune { 412fef7b6d8SPeter Brune PetscFunctionBegin; 413fef7b6d8SPeter Brune ((SNES_NCG *)(snes->data))->postcheckstep = func; 414fef7b6d8SPeter Brune ((SNES_NCG *)(snes->data))->postcheck = checkctx; 415fef7b6d8SPeter Brune PetscFunctionReturn(0); 416fef7b6d8SPeter Brune } 417fef7b6d8SPeter Brune EXTERN_C_END 418fef7b6d8SPeter Brune 419fef7b6d8SPeter Brune /*MC 420fef7b6d8SPeter Brune SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems. 421fef7b6d8SPeter Brune 422fef7b6d8SPeter Brune Level: beginner 423fef7b6d8SPeter Brune 424fef7b6d8SPeter Brune Options Database: 425fef7b6d8SPeter Brune + -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms) 426fef7b6d8SPeter Brune - -snes_ls <basic,basicnormnorms,quadratic> 427fef7b6d8SPeter Brune 428fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate 429fef7b6d8SPeter Brune gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 430fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x. 431fef7b6d8SPeter Brune 432fef7b6d8SPeter Brune 433fef7b6d8SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN 434fef7b6d8SPeter Brune M*/ 435fef7b6d8SPeter Brune EXTERN_C_BEGIN 436fef7b6d8SPeter Brune #undef __FUNCT__ 437fef7b6d8SPeter Brune #define __FUNCT__ "SNESCreate_NCG" 438fef7b6d8SPeter Brune PetscErrorCode SNESCreate_NCG(SNES snes) 439fef7b6d8SPeter Brune { 440fef7b6d8SPeter Brune SNES_NCG *neP; 441fef7b6d8SPeter Brune PetscErrorCode ierr; 442fef7b6d8SPeter Brune 443fef7b6d8SPeter Brune PetscFunctionBegin; 444fef7b6d8SPeter Brune snes->ops->destroy = SNESDestroy_NCG; 445fef7b6d8SPeter Brune snes->ops->setup = SNESSetUp_NCG; 446fef7b6d8SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NCG; 447fef7b6d8SPeter Brune snes->ops->view = SNESView_NCG; 448fef7b6d8SPeter Brune snes->ops->solve = SNESSolve_NCG; 449fef7b6d8SPeter Brune snes->ops->reset = SNESReset_NCG; 450fef7b6d8SPeter Brune 451fef7b6d8SPeter Brune snes->usesksp = PETSC_FALSE; 452fef7b6d8SPeter Brune snes->usespc = PETSC_TRUE; 453fef7b6d8SPeter Brune 454fef7b6d8SPeter Brune ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr); 455fef7b6d8SPeter Brune snes->data = (void*) neP; 456fef7b6d8SPeter Brune neP->maxstep = 1.e8; 457fef7b6d8SPeter Brune neP->steptol = 1.e-12; 458fef7b6d8SPeter Brune neP->type = SNES_LS_QUADRATIC; 459fef7b6d8SPeter Brune neP->damping = 1.0; 460fef7b6d8SPeter Brune neP->LineSearch = SNESLineSearchQuadratic_NCG; 461fef7b6d8SPeter Brune neP->lsP = PETSC_NULL; 462fef7b6d8SPeter Brune neP->postcheckstep = PETSC_NULL; 463fef7b6d8SPeter Brune neP->postcheck = PETSC_NULL; 464fef7b6d8SPeter Brune neP->precheckstep = PETSC_NULL; 465fef7b6d8SPeter Brune neP->precheck = PETSC_NULL; 466fef7b6d8SPeter Brune 467fef7b6d8SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_NCG",SNESLineSearchSetMonitor_NCG);CHKERRQ(ierr); 468fef7b6d8SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_NCG",SNESLineSearchSet_NCG);CHKERRQ(ierr); 469fef7b6d8SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_NCG",SNESLineSearchSetPostCheck_NCG);CHKERRQ(ierr); 470fef7b6d8SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_NCG",SNESLineSearchSetPreCheck_NCG);CHKERRQ(ierr); 471fef7b6d8SPeter Brune PetscFunctionReturn(0); 472fef7b6d8SPeter Brune } 473fef7b6d8SPeter Brune EXTERN_C_END 474