141ba4c6cSHeeho Park 241ba4c6cSHeeho Park #include <../src/snes/impls/ntrdc/ntrdcimpl.h> /*I "petscsnes.h" I*/ 341ba4c6cSHeeho Park 441ba4c6cSHeeho Park typedef struct { 541ba4c6cSHeeho Park SNES snes; 641ba4c6cSHeeho Park /* Information on the regular SNES convergence test; which may have been user provided 741ba4c6cSHeeho Park Copied from tr.c (maybe able to disposed, but this is a private function) - Heeho 841ba4c6cSHeeho Park Same with SNESTR_KSPConverged_Private, SNESTR_KSPConverged_Destroy, and SNESTR_Converged_Private 941ba4c6cSHeeho Park */ 1041ba4c6cSHeeho Park 1141ba4c6cSHeeho Park PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 1241ba4c6cSHeeho Park PetscErrorCode (*convdestroy)(void*); 1341ba4c6cSHeeho Park void *convctx; 1441ba4c6cSHeeho Park } SNES_TRDC_KSPConverged_Ctx; 1541ba4c6cSHeeho Park 1641ba4c6cSHeeho Park static PetscErrorCode SNESTRDC_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx) 1741ba4c6cSHeeho Park { 1841ba4c6cSHeeho Park SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx*)cctx; 1941ba4c6cSHeeho Park SNES snes = ctx->snes; 2041ba4c6cSHeeho Park SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC*)snes->data; 2141ba4c6cSHeeho Park Vec x; 2241ba4c6cSHeeho Park PetscReal nrm; 2341ba4c6cSHeeho Park PetscErrorCode ierr; 2441ba4c6cSHeeho Park 2541ba4c6cSHeeho Park PetscFunctionBegin; 2641ba4c6cSHeeho Park ierr = (*ctx->convtest)(ksp,n,rnorm,reason,ctx->convctx);CHKERRQ(ierr); 2741ba4c6cSHeeho Park if (*reason) { 287d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n",n,(double)rnorm);CHKERRQ(ierr); 2941ba4c6cSHeeho Park } 3041ba4c6cSHeeho Park /* Determine norm of solution */ 3141ba4c6cSHeeho Park ierr = KSPBuildSolution(ksp,NULL,&x);CHKERRQ(ierr); 3241ba4c6cSHeeho Park ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr); 3341ba4c6cSHeeho Park if (nrm >= neP->delta) { 347d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm);CHKERRQ(ierr); 3541ba4c6cSHeeho Park *reason = KSP_CONVERGED_STEP_LENGTH; 3641ba4c6cSHeeho Park } 3741ba4c6cSHeeho Park PetscFunctionReturn(0); 3841ba4c6cSHeeho Park } 3941ba4c6cSHeeho Park 4041ba4c6cSHeeho Park static PetscErrorCode SNESTRDC_KSPConverged_Destroy(void *cctx) 4141ba4c6cSHeeho Park { 4241ba4c6cSHeeho Park SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx*)cctx; 4341ba4c6cSHeeho Park PetscErrorCode ierr; 4441ba4c6cSHeeho Park 4541ba4c6cSHeeho Park PetscFunctionBegin; 4641ba4c6cSHeeho Park ierr = (*ctx->convdestroy)(ctx->convctx);CHKERRQ(ierr); 4741ba4c6cSHeeho Park ierr = PetscFree(ctx);CHKERRQ(ierr); 4841ba4c6cSHeeho Park 4941ba4c6cSHeeho Park PetscFunctionReturn(0); 5041ba4c6cSHeeho Park } 5141ba4c6cSHeeho Park 5241ba4c6cSHeeho Park /* ---------------------------------------------------------------- */ 5341ba4c6cSHeeho Park /* 5441ba4c6cSHeeho Park SNESTRDC_Converged_Private -test convergence JUST for 5541ba4c6cSHeeho Park the trust region tolerance. 5641ba4c6cSHeeho Park 5741ba4c6cSHeeho Park */ 5841ba4c6cSHeeho Park static PetscErrorCode SNESTRDC_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 5941ba4c6cSHeeho Park { 6041ba4c6cSHeeho Park SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC*)snes->data; 6141ba4c6cSHeeho Park PetscErrorCode ierr; 6241ba4c6cSHeeho Park 6341ba4c6cSHeeho Park PetscFunctionBegin; 6441ba4c6cSHeeho Park *reason = SNES_CONVERGED_ITERATING; 6541ba4c6cSHeeho Park if (neP->delta < xnorm * snes->deltatol) { 667d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"Diverged due to too small a trust region %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol);CHKERRQ(ierr); 6741ba4c6cSHeeho Park *reason = SNES_DIVERGED_TR_DELTA; 6841ba4c6cSHeeho Park } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 697d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n",snes->max_funcs);CHKERRQ(ierr); 7041ba4c6cSHeeho Park *reason = SNES_DIVERGED_FUNCTION_COUNT; 7141ba4c6cSHeeho Park } 7241ba4c6cSHeeho Park PetscFunctionReturn(0); 7341ba4c6cSHeeho Park } 7441ba4c6cSHeeho Park 7541ba4c6cSHeeho Park /*@ 7641ba4c6cSHeeho Park SNESNewtonTRDCGetRhoFlag - Get whether the solution update is within the trust-region. 7741ba4c6cSHeeho Park 7841ba4c6cSHeeho Park Input Parameters: 7941ba4c6cSHeeho Park . snes - the nonlinear solver object 8041ba4c6cSHeeho Park 8141ba4c6cSHeeho Park Output Parameters: 8241ba4c6cSHeeho Park . rho_flag: PETSC_TRUE if the solution update is in the trust-region; otherwise, PETSC_FALSE 8341ba4c6cSHeeho Park 8441ba4c6cSHeeho Park Level: developer 8541ba4c6cSHeeho Park 8641ba4c6cSHeeho Park @*/ 8741ba4c6cSHeeho Park PetscErrorCode SNESNewtonTRDCGetRhoFlag(SNES snes,PetscBool *rho_flag) 8841ba4c6cSHeeho Park { 8941ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 9041ba4c6cSHeeho Park 9141ba4c6cSHeeho Park PetscFunctionBegin; 9241ba4c6cSHeeho Park PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 9341ba4c6cSHeeho Park PetscValidBoolPointer(rho_flag,2); 9441ba4c6cSHeeho Park *rho_flag = tr->rho_satisfied; 9541ba4c6cSHeeho Park PetscFunctionReturn(0); 9641ba4c6cSHeeho Park } 9741ba4c6cSHeeho Park 9841ba4c6cSHeeho Park /*@C 9941ba4c6cSHeeho Park SNESNewtonTRDCSetPreCheck - Sets a user function that is called before the search step has been determined. 10041ba4c6cSHeeho Park Allows the user a chance to change or override the trust region decision. 10141ba4c6cSHeeho Park 10241ba4c6cSHeeho Park Logically Collective on snes 10341ba4c6cSHeeho Park 10441ba4c6cSHeeho Park Input Parameters: 10541ba4c6cSHeeho Park + snes - the nonlinear solver object 10641ba4c6cSHeeho Park . func - [optional] function evaluation routine, see SNESNewtonTRDCPreCheck() for the calling sequence 10741ba4c6cSHeeho Park - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 10841ba4c6cSHeeho Park 10941ba4c6cSHeeho Park Level: intermediate 11041ba4c6cSHeeho Park 11141ba4c6cSHeeho Park Note: This function is called BEFORE the function evaluation within the SNESNEWTONTRDC solver. 11241ba4c6cSHeeho Park 11341ba4c6cSHeeho Park .seealso: SNESNewtonTRDCPreCheck(), SNESNewtonTRDCGetPreCheck(), SNESNewtonTRDCSetPostCheck(), SNESNewtonTRDCGetPostCheck() 11441ba4c6cSHeeho Park @*/ 11541ba4c6cSHeeho Park PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES,Vec,Vec,PetscBool*,void*),void *ctx) 11641ba4c6cSHeeho Park { 11741ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 11841ba4c6cSHeeho Park 11941ba4c6cSHeeho Park PetscFunctionBegin; 12041ba4c6cSHeeho Park PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 12141ba4c6cSHeeho Park if (func) tr->precheck = func; 12241ba4c6cSHeeho Park if (ctx) tr->precheckctx = ctx; 12341ba4c6cSHeeho Park PetscFunctionReturn(0); 12441ba4c6cSHeeho Park } 12541ba4c6cSHeeho Park 12641ba4c6cSHeeho Park /*@C 12741ba4c6cSHeeho Park SNESNewtonTRDCGetPreCheck - Gets the pre-check function 12841ba4c6cSHeeho Park 12941ba4c6cSHeeho Park Not collective 13041ba4c6cSHeeho Park 13141ba4c6cSHeeho Park Input Parameter: 13241ba4c6cSHeeho Park . snes - the nonlinear solver context 13341ba4c6cSHeeho Park 13441ba4c6cSHeeho Park Output Parameters: 13541ba4c6cSHeeho Park + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRDCPreCheck() 13641ba4c6cSHeeho Park - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 13741ba4c6cSHeeho Park 13841ba4c6cSHeeho Park Level: intermediate 13941ba4c6cSHeeho Park 14041ba4c6cSHeeho Park .seealso: SNESNewtonTRDCSetPreCheck(), SNESNewtonTRDCPreCheck() 14141ba4c6cSHeeho Park @*/ 14241ba4c6cSHeeho Park PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES,Vec,Vec,PetscBool*,void*),void **ctx) 14341ba4c6cSHeeho Park { 14441ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 14541ba4c6cSHeeho Park 14641ba4c6cSHeeho Park PetscFunctionBegin; 14741ba4c6cSHeeho Park PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 14841ba4c6cSHeeho Park if (func) *func = tr->precheck; 14941ba4c6cSHeeho Park if (ctx) *ctx = tr->precheckctx; 15041ba4c6cSHeeho Park PetscFunctionReturn(0); 15141ba4c6cSHeeho Park } 15241ba4c6cSHeeho Park 15341ba4c6cSHeeho Park /*@C 15441ba4c6cSHeeho Park SNESNewtonTRDCSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 15541ba4c6cSHeeho Park function evaluation. Allows the user a chance to change or override the decision of the line search routine 15641ba4c6cSHeeho Park 15741ba4c6cSHeeho Park Logically Collective on snes 15841ba4c6cSHeeho Park 15941ba4c6cSHeeho Park Input Parameters: 16041ba4c6cSHeeho Park + snes - the nonlinear solver object 16141ba4c6cSHeeho Park . func - [optional] function evaluation routine, see SNESNewtonTRDCPostCheck() for the calling sequence 16241ba4c6cSHeeho Park - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 16341ba4c6cSHeeho Park 16441ba4c6cSHeeho Park Level: intermediate 16541ba4c6cSHeeho Park 16641ba4c6cSHeeho Park Note: This function is called BEFORE the function evaluation within the SNESNEWTONTRDC solver while the function set in 16741ba4c6cSHeeho Park SNESLineSearchSetPostCheck() is called AFTER the function evaluation. 16841ba4c6cSHeeho Park 16941ba4c6cSHeeho Park .seealso: SNESNewtonTRDCPostCheck(), SNESNewtonTRDCGetPostCheck() 17041ba4c6cSHeeho Park @*/ 17141ba4c6cSHeeho Park PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx) 17241ba4c6cSHeeho Park { 17341ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 17441ba4c6cSHeeho Park 17541ba4c6cSHeeho Park PetscFunctionBegin; 17641ba4c6cSHeeho Park PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 17741ba4c6cSHeeho Park if (func) tr->postcheck = func; 17841ba4c6cSHeeho Park if (ctx) tr->postcheckctx = ctx; 17941ba4c6cSHeeho Park PetscFunctionReturn(0); 18041ba4c6cSHeeho Park } 18141ba4c6cSHeeho Park 18241ba4c6cSHeeho Park /*@C 18341ba4c6cSHeeho Park SNESNewtonTRDCGetPostCheck - Gets the post-check function 18441ba4c6cSHeeho Park 18541ba4c6cSHeeho Park Not collective 18641ba4c6cSHeeho Park 18741ba4c6cSHeeho Park Input Parameter: 18841ba4c6cSHeeho Park . snes - the nonlinear solver context 18941ba4c6cSHeeho Park 19041ba4c6cSHeeho Park Output Parameters: 19141ba4c6cSHeeho Park + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRDCPostCheck() 19241ba4c6cSHeeho Park - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 19341ba4c6cSHeeho Park 19441ba4c6cSHeeho Park Level: intermediate 19541ba4c6cSHeeho Park 19641ba4c6cSHeeho Park .seealso: SNESNewtonTRDCSetPostCheck(), SNESNewtonTRDCPostCheck() 19741ba4c6cSHeeho Park @*/ 19841ba4c6cSHeeho Park PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES snes,PetscErrorCode (**func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx) 19941ba4c6cSHeeho Park { 20041ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 20141ba4c6cSHeeho Park 20241ba4c6cSHeeho Park PetscFunctionBegin; 20341ba4c6cSHeeho Park PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 20441ba4c6cSHeeho Park if (func) *func = tr->postcheck; 20541ba4c6cSHeeho Park if (ctx) *ctx = tr->postcheckctx; 20641ba4c6cSHeeho Park PetscFunctionReturn(0); 20741ba4c6cSHeeho Park } 20841ba4c6cSHeeho Park 20941ba4c6cSHeeho Park /*@C 21041ba4c6cSHeeho Park SNESNewtonTRDCPreCheck - Called before the step has been determined in SNESNEWTONTRDC 21141ba4c6cSHeeho Park 21241ba4c6cSHeeho Park Logically Collective on snes 21341ba4c6cSHeeho Park 21441ba4c6cSHeeho Park Input Parameters: 21541ba4c6cSHeeho Park + snes - the solver 21641ba4c6cSHeeho Park . X - The last solution 21741ba4c6cSHeeho Park - Y - The step direction 21841ba4c6cSHeeho Park 21941ba4c6cSHeeho Park Output Parameters: 22041ba4c6cSHeeho Park . changed_Y - Indicator that the step direction Y has been changed. 22141ba4c6cSHeeho Park 22241ba4c6cSHeeho Park Level: developer 22341ba4c6cSHeeho Park 22441ba4c6cSHeeho Park .seealso: SNESNewtonTRDCSetPreCheck(), SNESNewtonTRDCGetPreCheck() 22541ba4c6cSHeeho Park @*/ 22641ba4c6cSHeeho Park static PetscErrorCode SNESNewtonTRDCPreCheck(SNES snes,Vec X,Vec Y,PetscBool *changed_Y) 22741ba4c6cSHeeho Park { 22841ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 22941ba4c6cSHeeho Park PetscErrorCode ierr; 23041ba4c6cSHeeho Park 23141ba4c6cSHeeho Park PetscFunctionBegin; 23241ba4c6cSHeeho Park *changed_Y = PETSC_FALSE; 23341ba4c6cSHeeho Park if (tr->precheck) { 23441ba4c6cSHeeho Park ierr = (*tr->precheck)(snes,X,Y,changed_Y,tr->precheckctx);CHKERRQ(ierr); 23541ba4c6cSHeeho Park PetscValidLogicalCollectiveBool(snes,*changed_Y,4); 23641ba4c6cSHeeho Park } 23741ba4c6cSHeeho Park PetscFunctionReturn(0); 23841ba4c6cSHeeho Park } 23941ba4c6cSHeeho Park 24041ba4c6cSHeeho Park /*@C 24141ba4c6cSHeeho Park SNESNewtonTRDCPostCheck - Called after the step has been determined in SNESNEWTONTRDC but before the function evaluation 24241ba4c6cSHeeho Park 24341ba4c6cSHeeho Park Logically Collective on snes 24441ba4c6cSHeeho Park 24541ba4c6cSHeeho Park Input Parameters: 246*f1a722f8SMatthew G. Knepley + snes - the solver 247*f1a722f8SMatthew G. Knepley . X - The last solution 24841ba4c6cSHeeho Park . Y - The full step direction 24941ba4c6cSHeeho Park - W - The updated solution, W = X - Y 25041ba4c6cSHeeho Park 25141ba4c6cSHeeho Park Output Parameters: 25241ba4c6cSHeeho Park + changed_Y - indicator if step has been changed 25341ba4c6cSHeeho Park - changed_W - Indicator if the new candidate solution W has been changed. 25441ba4c6cSHeeho Park 25541ba4c6cSHeeho Park Notes: 25641ba4c6cSHeeho Park If Y is changed then W is recomputed as X - Y 25741ba4c6cSHeeho Park 25841ba4c6cSHeeho Park Level: developer 25941ba4c6cSHeeho Park 26041ba4c6cSHeeho Park .seealso: SNESNewtonTRDCSetPostCheck(), SNESNewtonTRDCGetPostCheck() 26141ba4c6cSHeeho Park @*/ 26241ba4c6cSHeeho Park static PetscErrorCode SNESNewtonTRDCPostCheck(SNES snes,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W) 26341ba4c6cSHeeho Park { 26441ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 26541ba4c6cSHeeho Park PetscErrorCode ierr; 26641ba4c6cSHeeho Park 26741ba4c6cSHeeho Park PetscFunctionBegin; 26841ba4c6cSHeeho Park *changed_Y = PETSC_FALSE; 26941ba4c6cSHeeho Park *changed_W = PETSC_FALSE; 27041ba4c6cSHeeho Park if (tr->postcheck) { 27141ba4c6cSHeeho Park ierr = (*tr->postcheck)(snes,X,Y,W,changed_Y,changed_W,tr->postcheckctx);CHKERRQ(ierr); 27241ba4c6cSHeeho Park PetscValidLogicalCollectiveBool(snes,*changed_Y,5); 27341ba4c6cSHeeho Park PetscValidLogicalCollectiveBool(snes,*changed_W,6); 27441ba4c6cSHeeho Park } 27541ba4c6cSHeeho Park PetscFunctionReturn(0); 27641ba4c6cSHeeho Park } 27741ba4c6cSHeeho Park 27841ba4c6cSHeeho Park /* 27941ba4c6cSHeeho Park SNESSolve_NEWTONTRDC - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy 28041ba4c6cSHeeho Park (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of 28141ba4c6cSHeeho Park nonlinear equations 28241ba4c6cSHeeho Park 28341ba4c6cSHeeho Park */ 28441ba4c6cSHeeho Park static PetscErrorCode SNESSolve_NEWTONTRDC(SNES snes) 28541ba4c6cSHeeho Park { 28641ba4c6cSHeeho Park SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC*)snes->data; 28741ba4c6cSHeeho Park Vec X,F,Y,G,W,GradF,YNtmp; 28841ba4c6cSHeeho Park Vec YCtmp; 28941ba4c6cSHeeho Park Mat jac; 29041ba4c6cSHeeho Park PetscErrorCode ierr; 29141ba4c6cSHeeho Park PetscInt maxits,i,j,lits,inner_count,bs; 29241ba4c6cSHeeho Park PetscReal rho,fnorm,gnorm,xnorm=0,delta,ynorm,temp_xnorm,temp_ynorm; /* TRDC inner iteration */ 29341ba4c6cSHeeho Park PetscReal inorms[99]; /* need to make it dynamic eventually, fixed max block size of 99 for now */ 29441ba4c6cSHeeho Park PetscReal deltaM,ynnorm,f0,mp,gTy,g,yTHy; /* rho calculation */ 29541ba4c6cSHeeho Park PetscReal auk,gfnorm,ycnorm,c0,c1,c2,tau,tau_pos,tau_neg,gTBg; /* Cauchy Point */ 29641ba4c6cSHeeho Park KSP ksp; 29741ba4c6cSHeeho Park SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 29841ba4c6cSHeeho Park PetscBool breakout = PETSC_FALSE; 29941ba4c6cSHeeho Park SNES_TRDC_KSPConverged_Ctx *ctx; 30041ba4c6cSHeeho Park PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),(*convdestroy)(void*); 30141ba4c6cSHeeho Park void *convctx; 30241ba4c6cSHeeho Park 30341ba4c6cSHeeho Park PetscFunctionBegin; 30441ba4c6cSHeeho Park maxits = snes->max_its; /* maximum number of iterations */ 30541ba4c6cSHeeho Park X = snes->vec_sol; /* solution vector */ 30641ba4c6cSHeeho Park F = snes->vec_func; /* residual vector */ 30741ba4c6cSHeeho Park Y = snes->work[0]; /* update vector */ 30841ba4c6cSHeeho Park G = snes->work[1]; /* updated residual */ 30941ba4c6cSHeeho Park W = snes->work[2]; /* temporary vector */ 31041ba4c6cSHeeho Park GradF = snes->work[3]; /* grad f = J^T F */ 31141ba4c6cSHeeho Park YNtmp = snes->work[4]; /* Newton solution */ 31241ba4c6cSHeeho Park YCtmp = snes->work[5]; /* Cauchy solution */ 31341ba4c6cSHeeho Park 3142c71b3e2SJacob Faibussowitsch PetscCheckFalse(snes->xl || snes->xu || snes->ops->computevariablebounds,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 31541ba4c6cSHeeho Park 31641ba4c6cSHeeho Park ierr = VecGetBlockSize(YNtmp,&bs);CHKERRQ(ierr); 31741ba4c6cSHeeho Park 31841ba4c6cSHeeho Park ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 31941ba4c6cSHeeho Park snes->iter = 0; 32041ba4c6cSHeeho Park ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 32141ba4c6cSHeeho Park 32241ba4c6cSHeeho Park /* Set the linear stopping criteria to use the More' trick. From tr.c */ 32341ba4c6cSHeeho Park ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 32441ba4c6cSHeeho Park ierr = KSPGetConvergenceTest(ksp,&convtest,&convctx,&convdestroy);CHKERRQ(ierr); 32541ba4c6cSHeeho Park if (convtest != SNESTRDC_KSPConverged_Private) { 32641ba4c6cSHeeho Park ierr = PetscNew(&ctx);CHKERRQ(ierr); 32741ba4c6cSHeeho Park ctx->snes = snes; 32841ba4c6cSHeeho Park ierr = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr); 32941ba4c6cSHeeho Park ierr = KSPSetConvergenceTest(ksp,SNESTRDC_KSPConverged_Private,ctx,SNESTRDC_KSPConverged_Destroy);CHKERRQ(ierr); 33041ba4c6cSHeeho Park ierr = PetscInfo(snes,"Using Krylov convergence test SNESTRDC_KSPConverged_Private\n");CHKERRQ(ierr); 33141ba4c6cSHeeho Park } 33241ba4c6cSHeeho Park 33341ba4c6cSHeeho Park if (!snes->vec_func_init_set) { 33441ba4c6cSHeeho Park ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 33541ba4c6cSHeeho Park } else snes->vec_func_init_set = PETSC_FALSE; 33641ba4c6cSHeeho Park 33741ba4c6cSHeeho Park ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 33841ba4c6cSHeeho Park SNESCheckFunctionNorm(snes,fnorm); 33941ba4c6cSHeeho Park ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- || X || */ 34041ba4c6cSHeeho Park 34141ba4c6cSHeeho Park ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 34241ba4c6cSHeeho Park snes->norm = fnorm; 34341ba4c6cSHeeho Park ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 34441ba4c6cSHeeho Park delta = xnorm ? neP->delta0*xnorm : neP->delta0; /* initial trust region size scaled by xnorm */ 34541ba4c6cSHeeho Park deltaM = xnorm ? neP->deltaM*xnorm : neP->deltaM; /* maximum trust region size scaled by xnorm */ 34641ba4c6cSHeeho Park neP->delta = delta; 34741ba4c6cSHeeho Park ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 34841ba4c6cSHeeho Park ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 34941ba4c6cSHeeho Park 35041ba4c6cSHeeho Park neP->rho_satisfied = PETSC_FALSE; 35141ba4c6cSHeeho Park 35241ba4c6cSHeeho Park /* test convergence */ 35341ba4c6cSHeeho Park ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 35441ba4c6cSHeeho Park if (snes->reason) PetscFunctionReturn(0); 35541ba4c6cSHeeho Park 35641ba4c6cSHeeho Park for (i=0; i<maxits; i++) { 35741ba4c6cSHeeho Park PetscBool changed_y; 35841ba4c6cSHeeho Park PetscBool changed_w; 35941ba4c6cSHeeho Park 36041ba4c6cSHeeho Park /* dogleg method */ 36141ba4c6cSHeeho Park ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 36241ba4c6cSHeeho Park SNESCheckJacobianDomainerror(snes); 36341ba4c6cSHeeho Park ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian);CHKERRQ(ierr); 36441ba4c6cSHeeho Park ierr = KSPSolve(snes->ksp,F,YNtmp);CHKERRQ(ierr); /* Quasi Newton Solution */ 36541ba4c6cSHeeho Park SNESCheckKSPSolve(snes); /* this is necessary but old tr.c did not have it*/ 36641ba4c6cSHeeho Park ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 36741ba4c6cSHeeho Park ierr = SNESGetJacobian(snes,&jac,NULL,NULL,NULL);CHKERRQ(ierr); 36841ba4c6cSHeeho Park 36941ba4c6cSHeeho Park /* rescale Jacobian, Newton solution update, and re-calculate delta for multiphase (multivariable) 37041ba4c6cSHeeho Park for inner iteration and Cauchy direction calculation 37141ba4c6cSHeeho Park */ 37241ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 37341ba4c6cSHeeho Park ierr = VecStrideNormAll(YNtmp,NORM_INFINITY,inorms);CHKERRQ(ierr); 37441ba4c6cSHeeho Park for (j=0; j<bs; j++) { 37541ba4c6cSHeeho Park if (neP->auto_scale_max > 1.0) { 37641ba4c6cSHeeho Park if (inorms[j] < 1.0/neP->auto_scale_max) { 37741ba4c6cSHeeho Park inorms[j] = 1.0/neP->auto_scale_max; 37841ba4c6cSHeeho Park } 37941ba4c6cSHeeho Park } 38041ba4c6cSHeeho Park ierr = VecStrideSet(W,j,inorms[j]);CHKERRQ(ierr); 38141ba4c6cSHeeho Park ierr = VecStrideScale(YNtmp,j,1.0/inorms[j]); 38241ba4c6cSHeeho Park ierr = VecStrideScale(X,j,1.0/inorms[j]); 38341ba4c6cSHeeho Park } 38441ba4c6cSHeeho Park ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); 38541ba4c6cSHeeho Park if (i == 0) { 38641ba4c6cSHeeho Park delta = neP->delta0*xnorm; 38741ba4c6cSHeeho Park } else { 38841ba4c6cSHeeho Park delta = neP->delta*xnorm; 38941ba4c6cSHeeho Park } 39041ba4c6cSHeeho Park deltaM = neP->deltaM*xnorm; 39141ba4c6cSHeeho Park ierr = MatDiagonalScale(jac,PETSC_NULL,W);CHKERRQ(ierr); 39241ba4c6cSHeeho Park } 39341ba4c6cSHeeho Park 39441ba4c6cSHeeho Park /* calculating GradF of minimization function */ 39541ba4c6cSHeeho Park ierr = MatMultTranspose(jac,F,GradF);CHKERRQ(ierr); /* grad f = J^T F */ 39641ba4c6cSHeeho Park ierr = VecNorm(YNtmp,NORM_2,&ynnorm);CHKERRQ(ierr); /* ynnorm <- || Y_newton || */ 39741ba4c6cSHeeho Park 39841ba4c6cSHeeho Park inner_count = 0; 39941ba4c6cSHeeho Park neP->rho_satisfied = PETSC_FALSE; 40041ba4c6cSHeeho Park while (1) { 40141ba4c6cSHeeho Park if (ynnorm <= delta) { /* see if the Newton solution is within the trust region */ 40241ba4c6cSHeeho Park ierr = VecCopy(YNtmp,Y);CHKERRQ(ierr); 40341ba4c6cSHeeho Park } else if (neP->use_cauchy) { /* use Cauchy direction if enabled */ 40441ba4c6cSHeeho Park ierr = MatMult(jac,GradF,W);CHKERRQ(ierr); 40541ba4c6cSHeeho Park ierr = VecDotRealPart(W,W,&gTBg);CHKERRQ(ierr); /* completes GradF^T J^T J GradF */ 40641ba4c6cSHeeho Park ierr = VecNorm(GradF,NORM_2,&gfnorm);CHKERRQ(ierr); /* grad f norm <- || grad f || */ 40741ba4c6cSHeeho Park if (gTBg <= 0.0) { 40841ba4c6cSHeeho Park auk = PETSC_MAX_REAL; 40941ba4c6cSHeeho Park } else { 41041ba4c6cSHeeho Park auk = PetscSqr(gfnorm)/gTBg; 41141ba4c6cSHeeho Park } 41241ba4c6cSHeeho Park auk = PetscMin(delta/gfnorm,auk); 41341ba4c6cSHeeho Park ierr = VecCopy(GradF,YCtmp);CHKERRQ(ierr); /* this could be improved */ 41441ba4c6cSHeeho Park ierr = VecScale(YCtmp,auk);CHKERRQ(ierr); /* YCtmp, Cauchy solution*/ 41541ba4c6cSHeeho Park ierr = VecNorm(YCtmp,NORM_2,&ycnorm);CHKERRQ(ierr); /* ycnorm <- || Y_cauchy || */ 41641ba4c6cSHeeho Park if (ycnorm >= delta) { /* see if the Cauchy solution meets the criteria */ 41741ba4c6cSHeeho Park ierr = VecCopy(YCtmp,Y);CHKERRQ(ierr); 4187d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"DL evaluated. delta: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n",(double)delta,(double)ynnorm,(double)ycnorm);CHKERRQ(ierr); 41941ba4c6cSHeeho Park } else { /* take ratio, tau, of Cauchy and Newton direction and step */ 42041ba4c6cSHeeho Park ierr = VecAXPY(YNtmp,-1.0,YCtmp);CHKERRQ(ierr); /* YCtmp = A, YNtmp = B */ 42141ba4c6cSHeeho Park ierr = VecNorm(YNtmp,NORM_2,&c0);CHKERRQ(ierr); /* this could be improved */ 42241ba4c6cSHeeho Park c0 = PetscSqr(c0); 42341ba4c6cSHeeho Park ierr = VecDotRealPart(YCtmp,YNtmp,&c1);CHKERRQ(ierr); 42441ba4c6cSHeeho Park c1 = 2.0*c1; 42541ba4c6cSHeeho Park ierr = VecNorm(YCtmp,NORM_2,&c2);CHKERRQ(ierr); /* this could be improved */ 42641ba4c6cSHeeho Park c2 = PetscSqr(c2) - PetscSqr(delta); 42741ba4c6cSHeeho Park tau_pos = (c1 + PetscSqrtReal(PetscSqr(c1) - 4.*c0*c2))/(2.*c0); /* quadratic formula */ 42841ba4c6cSHeeho Park tau_neg = (c1 - PetscSqrtReal(PetscSqr(c1) - 4.*c0*c2))/(2.*c0); 42941ba4c6cSHeeho Park tau = PetscMax(tau_pos, tau_neg); /* can tau_neg > tau_pos? I don't think so, but just in case. */ 4307d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"DL evaluated. tau: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n",(double)tau,(double)ynnorm,(double)ycnorm);CHKERRQ(ierr); 43141ba4c6cSHeeho Park ierr = VecWAXPY(W,tau,YNtmp,YCtmp);CHKERRQ(ierr); 43241ba4c6cSHeeho Park ierr = VecAXPY(W,-tau,YCtmp);CHKERRQ(ierr); 43341ba4c6cSHeeho Park ierr = VecCopy(W, Y);CHKERRQ(ierr); /* this could be improved */ 43441ba4c6cSHeeho Park } 43541ba4c6cSHeeho Park } else { 43641ba4c6cSHeeho Park /* if Cauchy is disabled, only use Newton direction */ 43741ba4c6cSHeeho Park auk = delta/ynnorm; 43841ba4c6cSHeeho Park ierr = VecScale(YNtmp,auk);CHKERRQ(ierr); 43941ba4c6cSHeeho Park ierr = VecCopy(YNtmp,Y);CHKERRQ(ierr); /* this could be improved (many VecCopy, VecNorm)*/ 44041ba4c6cSHeeho Park } 44141ba4c6cSHeeho Park 44241ba4c6cSHeeho Park ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); /* compute the final ynorm */ 44341ba4c6cSHeeho Park f0 = 0.5*PetscSqr(fnorm); /* minimizing function f(X) */ 44441ba4c6cSHeeho Park ierr = MatMult(jac,Y,W);CHKERRQ(ierr); 44541ba4c6cSHeeho Park ierr = VecDotRealPart(W,W,&yTHy);CHKERRQ(ierr); /* completes GradY^T J^T J GradY */ 44641ba4c6cSHeeho Park ierr = VecDotRealPart(GradF,Y,&gTy);CHKERRQ(ierr); 44741ba4c6cSHeeho Park mp = f0 - gTy + 0.5*yTHy; /* quadratic model to satisfy, -gTy because our update is X-Y*/ 44841ba4c6cSHeeho Park 44941ba4c6cSHeeho Park /* scale back solution update */ 45041ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 45141ba4c6cSHeeho Park for (j=0; j<bs; j++) { 45241ba4c6cSHeeho Park ierr = VecStrideScale(Y,j,inorms[j]); 45341ba4c6cSHeeho Park if (inner_count == 0) { 45441ba4c6cSHeeho Park /* TRDC inner algorithm does not need scaled X after calculating delta in the outer iteration */ 45541ba4c6cSHeeho Park /* need to scale back X to match Y and provide proper update to the external code */ 45641ba4c6cSHeeho Park ierr = VecStrideScale(X,j,inorms[j]); 45741ba4c6cSHeeho Park } 45841ba4c6cSHeeho Park } 45941ba4c6cSHeeho Park if (inner_count == 0) {ierr = VecNorm(X,NORM_2,&temp_xnorm);CHKERRQ(ierr);} /* only in the first iteration */ 46041ba4c6cSHeeho Park ierr = VecNorm(Y,NORM_2,&temp_ynorm);CHKERRQ(ierr); 46141ba4c6cSHeeho Park } else { 46241ba4c6cSHeeho Park temp_xnorm = xnorm; 46341ba4c6cSHeeho Park temp_ynorm = ynorm; 46441ba4c6cSHeeho Park } 46541ba4c6cSHeeho Park inner_count++; 46641ba4c6cSHeeho Park 46741ba4c6cSHeeho Park /* Evaluate the solution to meet the improvement ratio criteria */ 46841ba4c6cSHeeho Park ierr = SNESNewtonTRDCPreCheck(snes,X,Y,&changed_y);CHKERRQ(ierr); 46941ba4c6cSHeeho Park ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr); 47041ba4c6cSHeeho Park ierr = SNESNewtonTRDCPostCheck(snes,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 47141ba4c6cSHeeho Park if (changed_y) {ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr);} 47241ba4c6cSHeeho Park ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 47341ba4c6cSHeeho Park ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); /* F(X-Y) = G */ 47441ba4c6cSHeeho Park ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 47541ba4c6cSHeeho Park SNESCheckFunctionNorm(snes,gnorm); 47641ba4c6cSHeeho Park g = 0.5*PetscSqr(gnorm); /* minimizing function g(W) */ 47741ba4c6cSHeeho Park if (f0 == mp) rho = 0.0; 47841ba4c6cSHeeho Park else rho = (f0 - g)/(f0 - mp); /* actual improvement over predicted improvement */ 47941ba4c6cSHeeho Park 48041ba4c6cSHeeho Park if (rho < neP->eta2) { 48141ba4c6cSHeeho Park delta *= neP->t1; /* shrink the region */ 48241ba4c6cSHeeho Park } else if (rho > neP->eta3) { 48341ba4c6cSHeeho Park delta = PetscMin(neP->t2*delta,deltaM); /* expand the region, but not greater than deltaM */ 48441ba4c6cSHeeho Park } 48541ba4c6cSHeeho Park 48641ba4c6cSHeeho Park neP->delta = delta; 48741ba4c6cSHeeho Park if (rho >= neP->eta1) { 48841ba4c6cSHeeho Park /* unscale delta and xnorm before going to the next outer iteration */ 48941ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 49041ba4c6cSHeeho Park neP->delta = delta/xnorm; 49141ba4c6cSHeeho Park xnorm = temp_xnorm; 49241ba4c6cSHeeho Park ynorm = temp_ynorm; 49341ba4c6cSHeeho Park } 49441ba4c6cSHeeho Park neP->rho_satisfied = PETSC_TRUE; 49541ba4c6cSHeeho Park break; /* the improvement ratio is satisfactory */ 49641ba4c6cSHeeho Park } 49741ba4c6cSHeeho Park ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr); 49841ba4c6cSHeeho Park 49941ba4c6cSHeeho Park /* check to see if progress is hopeless */ 50041ba4c6cSHeeho Park neP->itflag = PETSC_FALSE; 50141ba4c6cSHeeho Park /* both delta, ynorm, and xnorm are either scaled or unscaled */ 50241ba4c6cSHeeho Park ierr = SNESTRDC_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 50341ba4c6cSHeeho Park if (!reason) { 50441ba4c6cSHeeho Park /* temp_xnorm, temp_ynorm is always unscaled */ 50541ba4c6cSHeeho Park /* also the inner iteration already calculated the Jacobian and solved the matrix */ 50641ba4c6cSHeeho Park /* therefore, it should be passing iteration number of iter+1 instead of iter+0 in the first iteration and after */ 50741ba4c6cSHeeho Park ierr = (*snes->ops->converged)(snes,snes->iter+1,temp_xnorm,temp_ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 50841ba4c6cSHeeho Park } 50941ba4c6cSHeeho Park /* if multiphase state changes, break out inner iteration */ 51041ba4c6cSHeeho Park if (reason == SNES_BREAKOUT_INNER_ITER) { 51141ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 51241ba4c6cSHeeho Park /* unscale delta and xnorm before going to the next outer iteration */ 51341ba4c6cSHeeho Park neP->delta = delta/xnorm; 51441ba4c6cSHeeho Park xnorm = temp_xnorm; 51541ba4c6cSHeeho Park ynorm = temp_ynorm; 51641ba4c6cSHeeho Park } 51741ba4c6cSHeeho Park reason = SNES_CONVERGED_ITERATING; 51841ba4c6cSHeeho Park break; 51941ba4c6cSHeeho Park } 52041ba4c6cSHeeho Park if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER; 52141ba4c6cSHeeho Park if (reason) { 52241ba4c6cSHeeho Park if (reason < 0) { 52341ba4c6cSHeeho Park /* We're not progressing, so return with the current iterate */ 52441ba4c6cSHeeho Park ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr); 52541ba4c6cSHeeho Park breakout = PETSC_TRUE; 52641ba4c6cSHeeho Park break; 52741ba4c6cSHeeho Park } else if (reason > 0) { 52841ba4c6cSHeeho Park /* We're converged, so return with the current iterate and update solution */ 52941ba4c6cSHeeho Park ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr); 53041ba4c6cSHeeho Park breakout = PETSC_FALSE; 53141ba4c6cSHeeho Park break; 53241ba4c6cSHeeho Park } 53341ba4c6cSHeeho Park } 53441ba4c6cSHeeho Park snes->numFailures++; 53541ba4c6cSHeeho Park } 53641ba4c6cSHeeho Park if (!breakout) { 53741ba4c6cSHeeho Park /* Update function and solution vectors */ 53841ba4c6cSHeeho Park fnorm = gnorm; 53941ba4c6cSHeeho Park ierr = VecCopy(G,F);CHKERRQ(ierr); 54041ba4c6cSHeeho Park ierr = VecCopy(W,X);CHKERRQ(ierr); 54141ba4c6cSHeeho Park /* Monitor convergence */ 54241ba4c6cSHeeho Park ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 54341ba4c6cSHeeho Park snes->iter = i+1; 54441ba4c6cSHeeho Park snes->norm = fnorm; 54541ba4c6cSHeeho Park snes->xnorm = xnorm; 54641ba4c6cSHeeho Park snes->ynorm = ynorm; 54741ba4c6cSHeeho Park ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 54841ba4c6cSHeeho Park ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 54941ba4c6cSHeeho Park ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 55041ba4c6cSHeeho Park /* Test for convergence, xnorm = || X || */ 55141ba4c6cSHeeho Park neP->itflag = PETSC_TRUE; 55241ba4c6cSHeeho Park if (snes->ops->converged != SNESConvergedSkip) {ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);} 55341ba4c6cSHeeho Park ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 55441ba4c6cSHeeho Park if (reason) break; 55541ba4c6cSHeeho Park } else break; 55641ba4c6cSHeeho Park } 55741ba4c6cSHeeho Park 55841ba4c6cSHeeho Park /* ierr = PetscFree(inorms);CHKERRQ(ierr); */ 55941ba4c6cSHeeho Park if (i == maxits) { 5607d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"Maximum number of iterations has been reached: %" PetscInt_FMT "\n",maxits);CHKERRQ(ierr); 56141ba4c6cSHeeho Park if (!reason) reason = SNES_DIVERGED_MAX_IT; 56241ba4c6cSHeeho Park } 56341ba4c6cSHeeho Park ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 56441ba4c6cSHeeho Park snes->reason = reason; 56541ba4c6cSHeeho Park ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 56641ba4c6cSHeeho Park if (convtest != SNESTRDC_KSPConverged_Private) { 56741ba4c6cSHeeho Park ierr = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr); 56841ba4c6cSHeeho Park ierr = PetscFree(ctx);CHKERRQ(ierr); 56941ba4c6cSHeeho Park ierr = KSPSetConvergenceTest(ksp,convtest,convctx,convdestroy);CHKERRQ(ierr); 57041ba4c6cSHeeho Park } 57141ba4c6cSHeeho Park PetscFunctionReturn(0); 57241ba4c6cSHeeho Park } 57341ba4c6cSHeeho Park 57441ba4c6cSHeeho Park /*------------------------------------------------------------*/ 57541ba4c6cSHeeho Park static PetscErrorCode SNESSetUp_NEWTONTRDC(SNES snes) 57641ba4c6cSHeeho Park { 57741ba4c6cSHeeho Park PetscErrorCode ierr; 57841ba4c6cSHeeho Park 57941ba4c6cSHeeho Park PetscFunctionBegin; 58041ba4c6cSHeeho Park ierr = SNESSetWorkVecs(snes,6);CHKERRQ(ierr); 58141ba4c6cSHeeho Park ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 58241ba4c6cSHeeho Park PetscFunctionReturn(0); 58341ba4c6cSHeeho Park } 58441ba4c6cSHeeho Park 58541ba4c6cSHeeho Park PetscErrorCode SNESReset_NEWTONTRDC(SNES snes) 58641ba4c6cSHeeho Park { 58741ba4c6cSHeeho Park 58841ba4c6cSHeeho Park PetscFunctionBegin; 58941ba4c6cSHeeho Park PetscFunctionReturn(0); 59041ba4c6cSHeeho Park } 59141ba4c6cSHeeho Park 59241ba4c6cSHeeho Park static PetscErrorCode SNESDestroy_NEWTONTRDC(SNES snes) 59341ba4c6cSHeeho Park { 59441ba4c6cSHeeho Park PetscErrorCode ierr; 59541ba4c6cSHeeho Park 59641ba4c6cSHeeho Park PetscFunctionBegin; 59741ba4c6cSHeeho Park ierr = SNESReset_NEWTONTRDC(snes);CHKERRQ(ierr); 59841ba4c6cSHeeho Park ierr = PetscFree(snes->data);CHKERRQ(ierr); 59941ba4c6cSHeeho Park PetscFunctionReturn(0); 60041ba4c6cSHeeho Park } 60141ba4c6cSHeeho Park /*------------------------------------------------------------*/ 60241ba4c6cSHeeho Park 60341ba4c6cSHeeho Park static PetscErrorCode SNESSetFromOptions_NEWTONTRDC(PetscOptionItems *PetscOptionsObject,SNES snes) 60441ba4c6cSHeeho Park { 60541ba4c6cSHeeho Park SNES_NEWTONTRDC *ctx = (SNES_NEWTONTRDC*)snes->data; 60641ba4c6cSHeeho Park PetscErrorCode ierr; 60741ba4c6cSHeeho Park 60841ba4c6cSHeeho Park PetscFunctionBegin; 60941ba4c6cSHeeho Park ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr); 61041ba4c6cSHeeho Park ierr = PetscOptionsReal("-snes_trdc_tol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr); 61141ba4c6cSHeeho Park ierr = PetscOptionsReal("-snes_trdc_eta1","eta1","None",ctx->eta1,&ctx->eta1,NULL);CHKERRQ(ierr); 61241ba4c6cSHeeho Park ierr = PetscOptionsReal("-snes_trdc_eta2","eta2","None",ctx->eta2,&ctx->eta2,NULL);CHKERRQ(ierr); 61341ba4c6cSHeeho Park ierr = PetscOptionsReal("-snes_trdc_eta3","eta3","None",ctx->eta3,&ctx->eta3,NULL);CHKERRQ(ierr); 61441ba4c6cSHeeho Park ierr = PetscOptionsReal("-snes_trdc_t1","t1","None",ctx->t1,&ctx->t1,NULL);CHKERRQ(ierr); 61541ba4c6cSHeeho Park ierr = PetscOptionsReal("-snes_trdc_t2","t2","None",ctx->t2,&ctx->t2,NULL);CHKERRQ(ierr); 61641ba4c6cSHeeho Park ierr = PetscOptionsReal("-snes_trdc_deltaM","deltaM","None",ctx->deltaM,&ctx->deltaM,NULL);CHKERRQ(ierr); 61741ba4c6cSHeeho Park ierr = PetscOptionsReal("-snes_trdc_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr); 61841ba4c6cSHeeho Park ierr = PetscOptionsReal("-snes_trdc_auto_scale_max","auto_scale_max","None",ctx->auto_scale_max,&ctx->auto_scale_max,NULL);CHKERRQ(ierr); 61941ba4c6cSHeeho Park ierr = PetscOptionsBool("-snes_trdc_use_cauchy","use_cauchy","use Cauchy step and direction",ctx->use_cauchy,&ctx->use_cauchy,NULL);CHKERRQ(ierr); 62041ba4c6cSHeeho Park ierr = PetscOptionsBool("-snes_trdc_auto_scale_multiphase","auto_scale_multiphase","Auto scaling for proper cauchy direction",ctx->auto_scale_multiphase,&ctx->auto_scale_multiphase,NULL);CHKERRQ(ierr); 62141ba4c6cSHeeho Park ierr = PetscOptionsTail();CHKERRQ(ierr); 62241ba4c6cSHeeho Park PetscFunctionReturn(0); 62341ba4c6cSHeeho Park } 62441ba4c6cSHeeho Park 62541ba4c6cSHeeho Park static PetscErrorCode SNESView_NEWTONTRDC(SNES snes,PetscViewer viewer) 62641ba4c6cSHeeho Park { 62741ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 62841ba4c6cSHeeho Park PetscErrorCode ierr; 62941ba4c6cSHeeho Park PetscBool iascii; 63041ba4c6cSHeeho Park 63141ba4c6cSHeeho Park PetscFunctionBegin; 63241ba4c6cSHeeho Park ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 63341ba4c6cSHeeho Park if (iascii) { 63441ba4c6cSHeeho Park ierr = PetscViewerASCIIPrintf(viewer," Trust region tolerance (-snes_trtol)\n",(double)snes->deltatol);CHKERRQ(ierr); 63541ba4c6cSHeeho Park ierr = PetscViewerASCIIPrintf(viewer," eta1=%g, eta2=%g, eta3=%g\n",(double)tr->eta1,(double)tr->eta2,(double)tr->eta3);CHKERRQ(ierr); 63641ba4c6cSHeeho Park ierr = PetscViewerASCIIPrintf(viewer," delta0=%g, t1=%g, t2=%g, deltaM=%g\n",(double)tr->delta0,(double)tr->t1,(double)tr->t2,(double)tr->deltaM);CHKERRQ(ierr); 63741ba4c6cSHeeho Park } 63841ba4c6cSHeeho Park PetscFunctionReturn(0); 63941ba4c6cSHeeho Park } 64041ba4c6cSHeeho Park /* ------------------------------------------------------------ */ 64141ba4c6cSHeeho Park /*MC 64241ba4c6cSHeeho Park SNESNEWTONTRDC - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction 64341ba4c6cSHeeho Park 64441ba4c6cSHeeho Park Options Database: 64541ba4c6cSHeeho Park + -snes_trdc_tol <tol> - trust region tolerance 64641ba4c6cSHeeho Park . -snes_trdc_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001) 64741ba4c6cSHeeho Park . -snes_trdc_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25) 64841ba4c6cSHeeho Park . -snes_trdc_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75) 64941ba4c6cSHeeho Park . -snes_trdc_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25) 65041ba4c6cSHeeho Park . -snes_trdc_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0) 65141ba4c6cSHeeho Park . -snes_trdc_deltaM <deltaM> - trust region parameter, max size of trust region, deltaM*norm2(x) (default: 0.5) 65241ba4c6cSHeeho Park . -snes_trdc_delta0 <delta0> - trust region parameter, initial size of trust region, delta0*norm2(x) (default: 0.1) 65341ba4c6cSHeeho Park . -snes_trdc_auto_scale_max <auto_scale_max> - used with auto_scale_multiphase, caps the maximum auto-scaling factor 65441ba4c6cSHeeho Park . -snes_trdc_use_cauchy <use_cauchy> - True uses dogleg Cauchy (Steepest Descent direction) step & direction in the trust region algorithm 65541ba4c6cSHeeho Park - -snes_trdc_auto_scale_multiphase <auto_scale_multiphase> - True turns on auto-scaling for multivariable block matrix for Cauchy and trust region 65641ba4c6cSHeeho Park 65741ba4c6cSHeeho Park Notes: 65841ba4c6cSHeeho Park The algorithm is taken from "Linear and Nonlinear Solvers for Simulating Multiphase Flow 65941ba4c6cSHeeho Park within Large-Scale Engineered Subsurface Systems" by Heeho D. Park, Glenn E. Hammond, 66041ba4c6cSHeeho Park Albert J. Valocchi, Tara LaForce. 66141ba4c6cSHeeho Park 66241ba4c6cSHeeho Park Level: intermediate 66341ba4c6cSHeeho Park 66441ba4c6cSHeeho Park .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance(), SNESNEWTONTRDC 66541ba4c6cSHeeho Park 66641ba4c6cSHeeho Park M*/ 66741ba4c6cSHeeho Park PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTRDC(SNES snes) 66841ba4c6cSHeeho Park { 66941ba4c6cSHeeho Park SNES_NEWTONTRDC *neP; 67041ba4c6cSHeeho Park PetscErrorCode ierr; 67141ba4c6cSHeeho Park 67241ba4c6cSHeeho Park PetscFunctionBegin; 67341ba4c6cSHeeho Park snes->ops->setup = SNESSetUp_NEWTONTRDC; 67441ba4c6cSHeeho Park snes->ops->solve = SNESSolve_NEWTONTRDC; 67541ba4c6cSHeeho Park snes->ops->destroy = SNESDestroy_NEWTONTRDC; 67641ba4c6cSHeeho Park snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTRDC; 67741ba4c6cSHeeho Park snes->ops->view = SNESView_NEWTONTRDC; 67841ba4c6cSHeeho Park snes->ops->reset = SNESReset_NEWTONTRDC; 67941ba4c6cSHeeho Park 68041ba4c6cSHeeho Park snes->usesksp = PETSC_TRUE; 68141ba4c6cSHeeho Park snes->usesnpc = PETSC_FALSE; 68241ba4c6cSHeeho Park 68341ba4c6cSHeeho Park snes->alwayscomputesfinalresidual = PETSC_TRUE; 68441ba4c6cSHeeho Park 68541ba4c6cSHeeho Park ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 68641ba4c6cSHeeho Park snes->data = (void*)neP; 68741ba4c6cSHeeho Park neP->delta = 0.0; 68841ba4c6cSHeeho Park neP->delta0 = 0.1; 68941ba4c6cSHeeho Park neP->eta1 = 0.001; 69041ba4c6cSHeeho Park neP->eta2 = 0.25; 69141ba4c6cSHeeho Park neP->eta3 = 0.75; 69241ba4c6cSHeeho Park neP->t1 = 0.25; 69341ba4c6cSHeeho Park neP->t2 = 2.0; 69441ba4c6cSHeeho Park neP->deltaM = 0.5; 69541ba4c6cSHeeho Park neP->sigma = 0.0001; 69641ba4c6cSHeeho Park neP->itflag = PETSC_FALSE; 69741ba4c6cSHeeho Park neP->rnorm0 = 0.0; 69841ba4c6cSHeeho Park neP->ttol = 0.0; 69941ba4c6cSHeeho Park neP->use_cauchy = PETSC_TRUE; 70041ba4c6cSHeeho Park neP->auto_scale_multiphase = PETSC_FALSE; 70141ba4c6cSHeeho Park neP->auto_scale_max = -1.0; 70241ba4c6cSHeeho Park neP->rho_satisfied = PETSC_FALSE; 70341ba4c6cSHeeho Park snes->deltatol = 1.e-12; 70441ba4c6cSHeeho Park 70541ba4c6cSHeeho Park /* for multiphase (multivariable) scaling */ 70641ba4c6cSHeeho Park /* may be used for dynamic allocation of inorms, but it fails snes_tutorials-ex3_13 70741ba4c6cSHeeho Park on test forced DIVERGED_JACOBIAN_DOMAIN test. I will use static array for now. 70841ba4c6cSHeeho Park ierr = VecGetBlockSize(snes->work[0],&neP->bs);CHKERRQ(ierr); 70941ba4c6cSHeeho Park ierr = PetscCalloc1(neP->bs,&neP->inorms);CHKERRQ(ierr); 71041ba4c6cSHeeho Park */ 71141ba4c6cSHeeho Park 71241ba4c6cSHeeho Park PetscFunctionReturn(0); 71341ba4c6cSHeeho Park } 714