1*97276fddSZach Atkins #include <../src/snes/impls/al/alimpl.h> /*I "petscsnes.h" I*/ 2*97276fddSZach Atkins 3*97276fddSZach Atkins /* 4*97276fddSZach Atkins This file implements a truncated Newton method with arc length continuation, 5*97276fddSZach Atkins for solving a system of nonlinear equations, using the KSP, Vec, 6*97276fddSZach Atkins and Mat interfaces for linear solvers, vectors, and matrices, 7*97276fddSZach Atkins respectively. 8*97276fddSZach Atkins */ 9*97276fddSZach Atkins 10*97276fddSZach Atkins const char NewtonALExactCitation[] = "@article{Ritto-CorreaCamotim2008,\n" 11*97276fddSZach Atkins " title={On the arc-length and other quadratic control methods: Established, less known and new implementation procedures},\n" 12*97276fddSZach Atkins " volume={86},\n" 13*97276fddSZach Atkins " ISSN={0045-7949},\n" 14*97276fddSZach Atkins " DOI={10.1016/j.compstruc.2007.08.003},\n" 15*97276fddSZach Atkins " number={11},\n" 16*97276fddSZach Atkins " journal={Computers & Structures},\n" 17*97276fddSZach Atkins " author={Ritto-Corr{\\^{e}}a, Manuel and Camotim, Dinar},\n" 18*97276fddSZach Atkins " year={2008},\n" 19*97276fddSZach Atkins " month=jun,\n" 20*97276fddSZach Atkins " pages={1353-1368},\n" 21*97276fddSZach Atkins "}\n"; 22*97276fddSZach Atkins PetscBool NewtonALExactCitationSet = PETSC_FALSE; 23*97276fddSZach Atkins const char NewtonALNormalCitation[] = "@article{LeonPaulinoPereiraMenezesLages_2011,\n" 24*97276fddSZach Atkins " title={A Unified Library of Nonlinear Solution Schemes},\n" 25*97276fddSZach Atkins " volume={64},\n" 26*97276fddSZach Atkins " ISSN={0003-6900, 2379-0407},\n" 27*97276fddSZach Atkins " DOI={10.1115/1.4006992},\n" 28*97276fddSZach Atkins " number={4},\n" 29*97276fddSZach Atkins " journal={Applied Mechanics Reviews},\n" 30*97276fddSZach Atkins " author={Leon, Sofie E. and Paulino, Glaucio H. and Pereira, Anderson and Menezes, Ivan F. M. and Lages, Eduardo N.},\n" 31*97276fddSZach Atkins " year={2011},\n" 32*97276fddSZach Atkins " month=jul,\n" 33*97276fddSZach Atkins " pages={040803},\n" 34*97276fddSZach Atkins " language={en}\n" 35*97276fddSZach Atkins "}\n"; 36*97276fddSZach Atkins PetscBool NewtonALNormalCitationSet = PETSC_FALSE; 37*97276fddSZach Atkins 38*97276fddSZach Atkins const char *const SNESNewtonALCorrectionTypes[] = {"EXACT", "NORMAL", "SNESNewtonALCorrectionType", "SNES_NEWTONAL_CORRECTION_", NULL}; 39*97276fddSZach Atkins 40*97276fddSZach Atkins static PetscErrorCode SNESNewtonALCheckArcLength(SNES snes, Vec XStep, PetscReal lambdaStep, PetscReal stepSize) 41*97276fddSZach Atkins { 42*97276fddSZach Atkins PetscReal arcLength, arcLengthError; 43*97276fddSZach Atkins SNES_NEWTONAL *al = (SNES_NEWTONAL *)snes->data; 44*97276fddSZach Atkins 45*97276fddSZach Atkins PetscFunctionBegin; 46*97276fddSZach Atkins PetscCall(VecNorm(XStep, NORM_2, &arcLength)); 47*97276fddSZach Atkins arcLength = PetscSqrtReal(PetscSqr(arcLength) + al->psisq * lambdaStep * lambdaStep); 48*97276fddSZach Atkins arcLengthError = PetscAbsReal(arcLength - stepSize); 49*97276fddSZach Atkins 50*97276fddSZach Atkins if (arcLengthError > PETSC_SQRT_MACHINE_EPSILON) PetscCall(PetscInfo(snes, "Arc length differs from specified step size: computed=%18.16e, expected=%18.16e, error=%18.16e \n", (double)arcLength, (double)stepSize, (double)arcLengthError)); 51*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 52*97276fddSZach Atkins } 53*97276fddSZach Atkins 54*97276fddSZach Atkins /* stable implementation of roots of a*x^2 + b*x + c = 0 */ 55*97276fddSZach Atkins static inline void PetscQuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal *xm, PetscReal *xp) 56*97276fddSZach Atkins { 57*97276fddSZach Atkins PetscReal temp = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(b * b - 4 * a * c)); 58*97276fddSZach Atkins PetscReal x1 = temp / a; 59*97276fddSZach Atkins PetscReal x2 = c / temp; 60*97276fddSZach Atkins *xm = PetscMin(x1, x2); 61*97276fddSZach Atkins *xp = PetscMax(x1, x2); 62*97276fddSZach Atkins } 63*97276fddSZach Atkins 64*97276fddSZach Atkins static PetscErrorCode SNESNewtonALSetCorrectionType_NEWTONAL(SNES snes, SNESNewtonALCorrectionType ctype) 65*97276fddSZach Atkins { 66*97276fddSZach Atkins SNES_NEWTONAL *al = (SNES_NEWTONAL *)snes->data; 67*97276fddSZach Atkins 68*97276fddSZach Atkins PetscFunctionBegin; 69*97276fddSZach Atkins al->correction_type = ctype; 70*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 71*97276fddSZach Atkins } 72*97276fddSZach Atkins 73*97276fddSZach Atkins /*@ 74*97276fddSZach Atkins SNESNewtonALSetCorrectionType - Set the type of correction to use in the arc-length continuation method. 75*97276fddSZach Atkins 76*97276fddSZach Atkins Logically Collective 77*97276fddSZach Atkins 78*97276fddSZach Atkins Input Parameters: 79*97276fddSZach Atkins + snes - the nonlinear solver object 80*97276fddSZach Atkins - ctype - the type of correction to use 81*97276fddSZach Atkins 82*97276fddSZach Atkins Options Database Key: 83*97276fddSZach Atkins . -snes_newtonal_correction_type <type> - Set the type of correction to use; use -help for a list of available types 84*97276fddSZach Atkins 85*97276fddSZach Atkins Level: intermediate 86*97276fddSZach Atkins 87*97276fddSZach Atkins .seealso: [](ch_snes), `SNES`, `SNESNEWTONAL`, `SNESNewtonALCorrectionType` 88*97276fddSZach Atkins @*/ 89*97276fddSZach Atkins PetscErrorCode SNESNewtonALSetCorrectionType(SNES snes, SNESNewtonALCorrectionType ctype) 90*97276fddSZach Atkins { 91*97276fddSZach Atkins PetscFunctionBegin; 92*97276fddSZach Atkins PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 93*97276fddSZach Atkins PetscValidLogicalCollectiveEnum(snes, ctype, 2); 94*97276fddSZach Atkins PetscTryMethod(snes, "SNESNewtonALSetCorrectionType_C", (SNES, SNESNewtonALCorrectionType), (snes, ctype)); 95*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 96*97276fddSZach Atkins } 97*97276fddSZach Atkins 98*97276fddSZach Atkins static PetscErrorCode SNESNewtonALSetFunction_NEWTONAL(SNES snes, SNESFunctionFn *func, void *ctx) 99*97276fddSZach Atkins { 100*97276fddSZach Atkins SNES_NEWTONAL *al = (SNES_NEWTONAL *)snes->data; 101*97276fddSZach Atkins 102*97276fddSZach Atkins PetscFunctionBegin; 103*97276fddSZach Atkins al->computealfunction = func; 104*97276fddSZach Atkins al->alctx = ctx; 105*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 106*97276fddSZach Atkins } 107*97276fddSZach Atkins 108*97276fddSZach Atkins /*@C 109*97276fddSZach Atkins SNESNewtonALSetFunction - Sets a user function that is called at each function evaluation to 110*97276fddSZach Atkins compute the tangent load vector for the arc-length continuation method. 111*97276fddSZach Atkins 112*97276fddSZach Atkins Logically Collective 113*97276fddSZach Atkins 114*97276fddSZach Atkins Input Parameters: 115*97276fddSZach Atkins + snes - the nonlinear solver object 116*97276fddSZach Atkins . func - [optional] tangent load function evaluation routine, see `SNESFunctionFn` for the calling sequence. `U` is the current solution vector, `Q` is the output tangent load vector 117*97276fddSZach Atkins - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 118*97276fddSZach Atkins 119*97276fddSZach Atkins Level: intermediate 120*97276fddSZach Atkins 121*97276fddSZach Atkins Notes: 122*97276fddSZach Atkins If the current value of the load parameter is needed in `func`, it can be obtained with `SNESNewtonALGetLoadParameter()`. 123*97276fddSZach Atkins 124*97276fddSZach Atkins The tangent load vector is the partial derivative of external load with respect to the load parameter. 125*97276fddSZach Atkins In the case of proportional loading, the tangent load vector is the full external load vector at the end of the load step. 126*97276fddSZach Atkins 127*97276fddSZach Atkins .seealso: [](ch_snes), `SNES`, `SNESNEWTONAL`, `SNESNewtonALGetFunction()`, `SNESNewtonALGetLoadParameter()` 128*97276fddSZach Atkins @*/ 129*97276fddSZach Atkins PetscErrorCode SNESNewtonALSetFunction(SNES snes, SNESFunctionFn *func, void *ctx) 130*97276fddSZach Atkins { 131*97276fddSZach Atkins PetscFunctionBegin; 132*97276fddSZach Atkins PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 133*97276fddSZach Atkins PetscTryMethod(snes, "SNESNewtonALSetFunction_C", (SNES, SNESFunctionFn *, void *), (snes, func, ctx)); 134*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 135*97276fddSZach Atkins } 136*97276fddSZach Atkins 137*97276fddSZach Atkins static PetscErrorCode SNESNewtonALGetFunction_NEWTONAL(SNES snes, SNESFunctionFn **func, void **ctx) 138*97276fddSZach Atkins { 139*97276fddSZach Atkins SNES_NEWTONAL *al = (SNES_NEWTONAL *)snes->data; 140*97276fddSZach Atkins 141*97276fddSZach Atkins PetscFunctionBegin; 142*97276fddSZach Atkins if (func) *func = al->computealfunction; 143*97276fddSZach Atkins if (ctx) *ctx = al->alctx; 144*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 145*97276fddSZach Atkins } 146*97276fddSZach Atkins 147*97276fddSZach Atkins /*@C 148*97276fddSZach Atkins SNESNewtonALGetFunction - Get the user function and context set with `SNESNewtonALSetFunction` 149*97276fddSZach Atkins 150*97276fddSZach Atkins Logically Collective 151*97276fddSZach Atkins 152*97276fddSZach Atkins Input Parameters: 153*97276fddSZach Atkins + snes - the nonlinear solver object 154*97276fddSZach Atkins . func - [optional] tangent load function evaluation routine, see `SNESNewtonALSetFunction()` for the call sequence 155*97276fddSZach Atkins - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 156*97276fddSZach Atkins 157*97276fddSZach Atkins Level: intermediate 158*97276fddSZach Atkins 159*97276fddSZach Atkins .seealso: [](ch_snes), `SNES`, `SNESNEWTONAL`, `SNESNewtonALSetFunction()` 160*97276fddSZach Atkins @*/ 161*97276fddSZach Atkins PetscErrorCode SNESNewtonALGetFunction(SNES snes, SNESFunctionFn **func, void **ctx) 162*97276fddSZach Atkins { 163*97276fddSZach Atkins PetscFunctionBegin; 164*97276fddSZach Atkins PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 165*97276fddSZach Atkins PetscUseMethod(snes, "SNESNewtonALGetFunction_C", (SNES, SNESFunctionFn **, void **), (snes, func, ctx)); 166*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 167*97276fddSZach Atkins } 168*97276fddSZach Atkins 169*97276fddSZach Atkins static PetscErrorCode SNESNewtonALGetLoadParameter_NEWTONAL(SNES snes, PetscReal *lambda) 170*97276fddSZach Atkins { 171*97276fddSZach Atkins SNES_NEWTONAL *al; 172*97276fddSZach Atkins 173*97276fddSZach Atkins PetscFunctionBeginHot; 174*97276fddSZach Atkins al = (SNES_NEWTONAL *)snes->data; 175*97276fddSZach Atkins *lambda = al->lambda; 176*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 177*97276fddSZach Atkins } 178*97276fddSZach Atkins 179*97276fddSZach Atkins /*@C 180*97276fddSZach Atkins SNESNewtonALGetLoadParameter - Get the value of the load parameter `lambda` for the arc-length continuation method. 181*97276fddSZach Atkins 182*97276fddSZach Atkins Logically Collective 183*97276fddSZach Atkins 184*97276fddSZach Atkins Input Parameter: 185*97276fddSZach Atkins . snes - the nonlinear solver object 186*97276fddSZach Atkins 187*97276fddSZach Atkins Output Parameter: 188*97276fddSZach Atkins . lambda - the arc-length parameter 189*97276fddSZach Atkins 190*97276fddSZach Atkins Level: intermediate 191*97276fddSZach Atkins 192*97276fddSZach Atkins Notes: 193*97276fddSZach Atkins This function should be used in the functions provided to `SNESSetFunction()` and `SNESNewtonALSetFunction()` 194*97276fddSZach Atkins to compute the residual and tangent load vectors for a given value of `lambda` (0 <= lambda <= 1). 195*97276fddSZach Atkins 196*97276fddSZach Atkins Usually, `lambda` is used to scale the external force vector in the residual function, i.e. proportional loading, 197*97276fddSZach Atkins in which case the tangent load vector is the full external force vector. 198*97276fddSZach Atkins 199*97276fddSZach Atkins .seealso: [](ch_snes), `SNES`, `SNESNEWTONAL`, `SNESNewtonALSetFunction()` 200*97276fddSZach Atkins @*/ 201*97276fddSZach Atkins PetscErrorCode SNESNewtonALGetLoadParameter(SNES snes, PetscReal *lambda) 202*97276fddSZach Atkins { 203*97276fddSZach Atkins PetscFunctionBeginHot; 204*97276fddSZach Atkins PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 205*97276fddSZach Atkins PetscAssertPointer(lambda, 2); 206*97276fddSZach Atkins PetscUseMethod(snes, "SNESNewtonALGetLoadParameter_C", (SNES, PetscReal *), (snes, lambda)); 207*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 208*97276fddSZach Atkins } 209*97276fddSZach Atkins 210*97276fddSZach Atkins static PetscErrorCode SNESNewtonALComputeFunction_NEWTONAL(SNES snes, Vec X, Vec Q) 211*97276fddSZach Atkins { 212*97276fddSZach Atkins void *ctx = NULL; 213*97276fddSZach Atkins SNESFunctionFn *computealfunction = NULL; 214*97276fddSZach Atkins SNES_NEWTONAL *al; 215*97276fddSZach Atkins 216*97276fddSZach Atkins PetscFunctionBegin; 217*97276fddSZach Atkins al = (SNES_NEWTONAL *)snes->data; 218*97276fddSZach Atkins PetscCall(SNESNewtonALGetFunction(snes, &computealfunction, &ctx)); 219*97276fddSZach Atkins 220*97276fddSZach Atkins PetscCall(VecZeroEntries(Q)); 221*97276fddSZach Atkins PetscCheck(computealfunction || (snes->vec_rhs && al->scale_rhs), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "No tangent load function or rhs vector has been set"); 222*97276fddSZach Atkins if (computealfunction) { 223*97276fddSZach Atkins PetscCall(VecLockReadPush(X)); 224*97276fddSZach Atkins PetscCallBack("SNES callback NewtonAL tangent load function", (*computealfunction)(snes, X, Q, ctx)); 225*97276fddSZach Atkins PetscCall(VecLockReadPop(X)); 226*97276fddSZach Atkins } 227*97276fddSZach Atkins if (al->scale_rhs && snes->vec_rhs) { 228*97276fddSZach Atkins /* Save original RHS vector values, then scale `snes->vec_rhs` by load parameter */ 229*97276fddSZach Atkins if (!al->vec_rhs_orig) { PetscCall(VecDuplicate(snes->vec_rhs, &al->vec_rhs_orig)); } 230*97276fddSZach Atkins if (!al->copied_rhs) { 231*97276fddSZach Atkins PetscCall(VecCopy(snes->vec_rhs, al->vec_rhs_orig)); 232*97276fddSZach Atkins al->copied_rhs = PETSC_TRUE; 233*97276fddSZach Atkins } 234*97276fddSZach Atkins PetscCall(VecAXPBY(snes->vec_rhs, al->lambda, 0.0, al->vec_rhs_orig)); 235*97276fddSZach Atkins PetscCall(VecAXPY(Q, 1, al->vec_rhs_orig)); 236*97276fddSZach Atkins } 237*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 238*97276fddSZach Atkins } 239*97276fddSZach Atkins 240*97276fddSZach Atkins /*@C 241*97276fddSZach Atkins SNESNewtonALComputeFunction - Calls the function that has been set with `SNESNewtonALSetFunction()`. 242*97276fddSZach Atkins 243*97276fddSZach Atkins Collective 244*97276fddSZach Atkins 245*97276fddSZach Atkins Input Parameters: 246*97276fddSZach Atkins + snes - the `SNES` context 247*97276fddSZach Atkins - X - input vector 248*97276fddSZach Atkins 249*97276fddSZach Atkins Output Parameter: 250*97276fddSZach Atkins . Q - tangent load vector, as set by `SNESNewtonALSetFunction()` 251*97276fddSZach Atkins 252*97276fddSZach Atkins Level: developer 253*97276fddSZach Atkins 254*97276fddSZach Atkins Note: 255*97276fddSZach Atkins `SNESNewtonALComputeFunction()` is typically used within nonlinear solvers 256*97276fddSZach Atkins implementations, so users would not generally call this routine themselves. 257*97276fddSZach Atkins 258*97276fddSZach Atkins .seealso: [](ch_snes), `SNES`, `SNESNewtonALSetFunction()`, `SNESNewtonALGetFunction()` 259*97276fddSZach Atkins @*/ 260*97276fddSZach Atkins PetscErrorCode SNESNewtonALComputeFunction(SNES snes, Vec X, Vec Q) 261*97276fddSZach Atkins { 262*97276fddSZach Atkins PetscFunctionBegin; 263*97276fddSZach Atkins PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 264*97276fddSZach Atkins PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 265*97276fddSZach Atkins PetscValidHeaderSpecific(Q, VEC_CLASSID, 3); 266*97276fddSZach Atkins PetscCheckSameComm(snes, 1, X, 2); 267*97276fddSZach Atkins PetscCheckSameComm(snes, 1, Q, 3); 268*97276fddSZach Atkins PetscCall(VecValidValues_Internal(X, 2, PETSC_TRUE)); 269*97276fddSZach Atkins PetscCall(PetscLogEventBegin(SNES_NewtonALEval, snes, X, Q, 0)); 270*97276fddSZach Atkins PetscTryMethod(snes, "SNESNewtonALComputeFunction_C", (SNES, Vec, Vec), (snes, X, Q)); 271*97276fddSZach Atkins PetscCall(PetscLogEventEnd(SNES_NewtonALEval, snes, X, Q, 0)); 272*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 273*97276fddSZach Atkins } 274*97276fddSZach Atkins 275*97276fddSZach Atkins /* 276*97276fddSZach Atkins SNESSolve_NEWTONAL - Solves a nonlinear system with Newton's method with arc length continuation. 277*97276fddSZach Atkins 278*97276fddSZach Atkins Input Parameter: 279*97276fddSZach Atkins . snes - the `SNES` context 280*97276fddSZach Atkins 281*97276fddSZach Atkins Application Interface Routine: SNESSolve() 282*97276fddSZach Atkins */ 283*97276fddSZach Atkins static PetscErrorCode SNESSolve_NEWTONAL(SNES snes) 284*97276fddSZach Atkins { 285*97276fddSZach Atkins SNES_NEWTONAL *data = (SNES_NEWTONAL *)snes->data; 286*97276fddSZach Atkins PetscInt maxits, maxincs, lits; 287*97276fddSZach Atkins PetscReal fnorm, xnorm, ynorm, stepSize; 288*97276fddSZach Atkins Vec DeltaX, deltaX, X, R, Q, deltaX_Q, deltaX_R; 289*97276fddSZach Atkins 290*97276fddSZach Atkins PetscFunctionBegin; 291*97276fddSZach Atkins PetscCheck(!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); 292*97276fddSZach Atkins 293*97276fddSZach Atkins /* Register citations */ 294*97276fddSZach Atkins PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 295*97276fddSZach Atkins if (data->correction_type == SNES_NEWTONAL_CORRECTION_EXACT) { 296*97276fddSZach Atkins PetscCall(PetscCitationsRegister(NewtonALExactCitation, &NewtonALExactCitationSet)); 297*97276fddSZach Atkins } else if (data->correction_type == SNES_NEWTONAL_CORRECTION_NORMAL) { 298*97276fddSZach Atkins PetscCall(PetscCitationsRegister(NewtonALNormalCitation, &NewtonALNormalCitationSet)); 299*97276fddSZach Atkins } 300*97276fddSZach Atkins 301*97276fddSZach Atkins data->copied_rhs = PETSC_FALSE; 302*97276fddSZach Atkins data->lambda_update = 0.0; 303*97276fddSZach Atkins data->lambda = 0.0; 304*97276fddSZach Atkins snes->numFailures = 0; 305*97276fddSZach Atkins snes->numLinearSolveFailures = 0; 306*97276fddSZach Atkins snes->reason = SNES_CONVERGED_ITERATING; 307*97276fddSZach Atkins 308*97276fddSZach Atkins maxits = snes->max_its; /* maximum number of iterations */ 309*97276fddSZach Atkins maxincs = data->max_continuation_steps; /* maximum number of increments */ 310*97276fddSZach Atkins X = snes->vec_sol; /* solution vector */ 311*97276fddSZach Atkins R = snes->vec_func; /* residual vector */ 312*97276fddSZach Atkins Q = snes->work[0]; /* tangent load vector */ 313*97276fddSZach Atkins deltaX_Q = snes->work[1]; /* variation of X with respect to lambda */ 314*97276fddSZach Atkins deltaX_R = snes->work[2]; /* linearized error correction */ 315*97276fddSZach Atkins DeltaX = snes->work[3]; /* step from equilibrium */ 316*97276fddSZach Atkins deltaX = snes->vec_sol_update; /* full newton step */ 317*97276fddSZach Atkins stepSize = data->step_size; /* initial step size */ 318*97276fddSZach Atkins 319*97276fddSZach Atkins PetscCall(VecZeroEntries(DeltaX)); 320*97276fddSZach Atkins 321*97276fddSZach Atkins /* set snes->max_its for convergence test */ 322*97276fddSZach Atkins snes->max_its = maxits * maxincs; 323*97276fddSZach Atkins 324*97276fddSZach Atkins /* main incremental-iterative loop */ 325*97276fddSZach Atkins for (PetscInt i = 0; i < maxincs || maxincs < 0; i++) { 326*97276fddSZach Atkins PetscReal deltaLambda; 327*97276fddSZach Atkins 328*97276fddSZach Atkins PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 329*97276fddSZach Atkins snes->norm = 0.0; 330*97276fddSZach Atkins PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 331*97276fddSZach Atkins PetscCall(SNESNewtonALComputeFunction(snes, X, Q)); 332*97276fddSZach Atkins PetscCall(SNESComputeFunction(snes, X, R)); 333*97276fddSZach Atkins PetscCall(VecAXPY(R, 1, Q)); /* R <- R + Q */ 334*97276fddSZach Atkins PetscCall(VecNorm(R, NORM_2, &fnorm)); /* fnorm <- ||R|| */ 335*97276fddSZach Atkins SNESCheckFunctionNorm(snes, fnorm); 336*97276fddSZach Atkins 337*97276fddSZach Atkins /* Monitor convergence */ 338*97276fddSZach Atkins PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm)); 339*97276fddSZach Atkins PetscCall(SNESMonitor(snes, snes->iter, fnorm)); 340*97276fddSZach Atkins if (i == 0 && snes->reason) break; 341*97276fddSZach Atkins for (PetscInt j = 0; j < maxits; j++) { 342*97276fddSZach Atkins PetscReal normsqX_Q, deltaS = 1; 343*97276fddSZach Atkins 344*97276fddSZach Atkins /* Call general purpose update function */ 345*97276fddSZach Atkins PetscTryTypeMethod(snes, update, snes->iter); 346*97276fddSZach Atkins 347*97276fddSZach Atkins PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 348*97276fddSZach Atkins SNESCheckJacobianDomainerror(snes); 349*97276fddSZach Atkins PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 350*97276fddSZach Atkins /* Solve J deltaX_Q = Q, where J is Jacobian matrix */ 351*97276fddSZach Atkins PetscCall(KSPSolve(snes->ksp, Q, deltaX_Q)); 352*97276fddSZach Atkins SNESCheckKSPSolve(snes); 353*97276fddSZach Atkins PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 354*97276fddSZach Atkins PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", tangent load linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 355*97276fddSZach Atkins /* Compute load parameter variation */ 356*97276fddSZach Atkins PetscCall(VecNorm(deltaX_Q, NORM_2, &normsqX_Q)); 357*97276fddSZach Atkins normsqX_Q *= normsqX_Q; 358*97276fddSZach Atkins /* On first iter, use predictor. This is the same regardless of corrector scheme. */ 359*97276fddSZach Atkins if (j == 0) { 360*97276fddSZach Atkins PetscReal sign = 1.0; 361*97276fddSZach Atkins if (i > 0) { 362*97276fddSZach Atkins PetscCall(VecDotRealPart(DeltaX, deltaX_Q, &sign)); 363*97276fddSZach Atkins sign += data->psisq * data->lambda_update; 364*97276fddSZach Atkins sign = sign >= 0 ? 1.0 : -1.0; 365*97276fddSZach Atkins } 366*97276fddSZach Atkins data->lambda_update = 0.0; 367*97276fddSZach Atkins PetscCall(VecZeroEntries(DeltaX)); 368*97276fddSZach Atkins deltaLambda = sign * stepSize / PetscSqrtReal(normsqX_Q + data->psisq); 369*97276fddSZach Atkins } else { 370*97276fddSZach Atkins /* Solve J deltaX_R = -R */ 371*97276fddSZach Atkins PetscCall(KSPSolve(snes->ksp, R, deltaX_R)); 372*97276fddSZach Atkins SNESCheckKSPSolve(snes); 373*97276fddSZach Atkins PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 374*97276fddSZach Atkins PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", residual linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 375*97276fddSZach Atkins PetscCall(VecScale(deltaX_R, -1)); 376*97276fddSZach Atkins 377*97276fddSZach Atkins if (data->correction_type == SNES_NEWTONAL_CORRECTION_NORMAL) { 378*97276fddSZach Atkins /* 379*97276fddSZach Atkins Take a step orthogonal to the current incremental update DeltaX. 380*97276fddSZach Atkins Note, this approach is cheaper than the exact correction, but may exhibit convergence 381*97276fddSZach Atkins issues due to the iterative trial points not being on the quadratic contraint surface. 382*97276fddSZach Atkins On the bright side, we always have a real and unique solution for deltaLambda. 383*97276fddSZach Atkins */ 384*97276fddSZach Atkins PetscScalar coefs[2]; 385*97276fddSZach Atkins Vec rhs[] = {deltaX_R, deltaX_Q}; 386*97276fddSZach Atkins 387*97276fddSZach Atkins PetscCall(VecMDot(DeltaX, 2, rhs, coefs)); 388*97276fddSZach Atkins deltaLambda = -PetscRealPart(coefs[0]) / (PetscRealPart(coefs[1]) + data->psisq * data->lambda_update); 389*97276fddSZach Atkins } else { 390*97276fddSZach Atkins /* 391*97276fddSZach Atkins Solve 392*97276fddSZach Atkins a*deltaLambda^2 + b*deltaLambda + c = 0 (*) 393*97276fddSZach Atkins where 394*97276fddSZach Atkins a = a0 395*97276fddSZach Atkins b = b0 + b1*deltaS 396*97276fddSZach Atkins c = c0 + c1*deltaS + c2*deltaS^2 397*97276fddSZach Atkins and deltaS is either 1, or the largest value in (0, 1) that satisfies 398*97276fddSZach Atkins b^2 - 4*a*c = as*deltaS^2 + bs*deltaS + cs >= 0 399*97276fddSZach Atkins where 400*97276fddSZach Atkins as = b1^2 - 4*a0*c2 401*97276fddSZach Atkins bs = 2*b1*b0 - 4*a0*c1 402*97276fddSZach Atkins cs = b0^2 - 4*a0*c0 403*97276fddSZach Atkins These "partial corrections" prevent (*) from having complex roots. 404*97276fddSZach Atkins */ 405*97276fddSZach Atkins PetscReal psisqLambdaUpdate, discriminant; 406*97276fddSZach Atkins PetscReal as, bs, cs; 407*97276fddSZach Atkins PetscReal a0, b0, b1, c0, c1, c2; 408*97276fddSZach Atkins PetscScalar coefs1[3]; /* coefs[0] = deltaX_Q*DeltaX, coefs[1] = deltaX_R*DeltaX, coefs[2] = DeltaX*DeltaX */ 409*97276fddSZach Atkins PetscScalar coefs2[2]; /* coefs[0] = deltaX_Q*deltaX_R, coefs[1] = deltaX_R*deltaX_R */ 410*97276fddSZach Atkins const Vec rhs1[3] = {deltaX_Q, deltaX_R, DeltaX}; 411*97276fddSZach Atkins const Vec rhs2[2] = {deltaX_Q, deltaX_R}; 412*97276fddSZach Atkins 413*97276fddSZach Atkins psisqLambdaUpdate = data->psisq * data->lambda_update; 414*97276fddSZach Atkins PetscCall(VecMDotBegin(DeltaX, 3, rhs1, coefs1)); 415*97276fddSZach Atkins PetscCall(VecMDotBegin(deltaX_R, 2, rhs2, coefs2)); 416*97276fddSZach Atkins PetscCall(VecMDotEnd(DeltaX, 3, rhs1, coefs1)); 417*97276fddSZach Atkins PetscCall(VecMDotEnd(deltaX_R, 2, rhs2, coefs2)); 418*97276fddSZach Atkins 419*97276fddSZach Atkins a0 = normsqX_Q + data->psisq; 420*97276fddSZach Atkins b0 = 2 * (PetscRealPart(coefs1[0]) + psisqLambdaUpdate); 421*97276fddSZach Atkins b1 = 2 * PetscRealPart(coefs2[0]); 422*97276fddSZach Atkins c0 = PetscRealPart(coefs1[2]) + psisqLambdaUpdate * data->lambda_update - stepSize * stepSize; 423*97276fddSZach Atkins c1 = 2 * PetscRealPart(coefs1[1]); 424*97276fddSZach Atkins c2 = PetscRealPart(coefs2[1]); 425*97276fddSZach Atkins 426*97276fddSZach Atkins as = b1 * b1 - 4 * a0 * c2; 427*97276fddSZach Atkins bs = 2 * (b1 * b0 - 2 * a0 * c1); 428*97276fddSZach Atkins cs = b0 * b0 - 4 * a0 * c0; 429*97276fddSZach Atkins 430*97276fddSZach Atkins discriminant = cs + bs * deltaS + as * deltaS * deltaS; 431*97276fddSZach Atkins 432*97276fddSZach Atkins if (discriminant < 0) { 433*97276fddSZach Atkins /* Take deltaS < 1 with the unique root -b/(2*a) */ 434*97276fddSZach Atkins PetscReal x1; 435*97276fddSZach Atkins 436*97276fddSZach Atkins /* Compute deltaS to be the largest root of (as * x^2 + bs * x + cs = 0) */ 437*97276fddSZach Atkins PetscQuadraticRoots(as, bs, cs, &x1, &deltaS); 438*97276fddSZach Atkins PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", discriminant=%18.16e < 0, shrinking residual update size to deltaS = %18.16e\n", snes->iter, (double)discriminant, (double)deltaS)); 439*97276fddSZach Atkins deltaLambda = -0.5 * (b0 + b1 * deltaS) / a0; 440*97276fddSZach Atkins } else { 441*97276fddSZach Atkins /* Use deltaS = 1, pick root that is closest to the last point to prevent doubling back */ 442*97276fddSZach Atkins PetscReal dlambda1, dlambda2; 443*97276fddSZach Atkins 444*97276fddSZach Atkins PetscQuadraticRoots(a0, b0 + b1, c0 + c1 + c2, &dlambda1, &dlambda2); 445*97276fddSZach Atkins deltaLambda = (b0 * dlambda1) > (b0 * dlambda2) ? dlambda1 : dlambda2; 446*97276fddSZach Atkins } 447*97276fddSZach Atkins } 448*97276fddSZach Atkins } 449*97276fddSZach Atkins PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 450*97276fddSZach Atkins data->lambda = data->lambda + deltaLambda; 451*97276fddSZach Atkins if (data->lambda > data->lambda_max) { 452*97276fddSZach Atkins /* Ensure that lambda = lambda_max exactly at the end of incremental process. This ensures the final solution matches the problem we want to solve. */ 453*97276fddSZach Atkins deltaLambda = deltaLambda - (data->lambda - data->lambda_max); 454*97276fddSZach Atkins data->lambda = data->lambda_max; 455*97276fddSZach Atkins } 456*97276fddSZach Atkins if (data->lambda < data->lambda_min) { 457*97276fddSZach Atkins // LCOV_EXCL_START 458*97276fddSZach Atkins /* Ensure that lambda >= lambda_min. This prevents some potential oscillatory behavior. */ 459*97276fddSZach Atkins deltaLambda = deltaLambda - (data->lambda - data->lambda_min); 460*97276fddSZach Atkins data->lambda = data->lambda_min; 461*97276fddSZach Atkins // LCOV_EXCL_STOP 462*97276fddSZach Atkins } 463*97276fddSZach Atkins data->lambda_update = data->lambda_update + deltaLambda; 464*97276fddSZach Atkins PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 465*97276fddSZach Atkins PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", lambda=%18.16e, lambda_update=%18.16e\n", snes->iter, (double)data->lambda, (double)data->lambda_update)); 466*97276fddSZach Atkins if (j == 0) { 467*97276fddSZach Atkins /* deltaX = deltaLambda*deltaX_Q */ 468*97276fddSZach Atkins PetscCall(VecCopy(deltaX_Q, deltaX)); 469*97276fddSZach Atkins PetscCall(VecScale(deltaX, deltaLambda)); 470*97276fddSZach Atkins } else { 471*97276fddSZach Atkins /* deltaX = deltaS*deltaX_R + deltaLambda*deltaX_Q */ 472*97276fddSZach Atkins PetscCall(VecAXPBYPCZ(deltaX, deltaS, deltaLambda, 0, deltaX_R, deltaX_Q)); 473*97276fddSZach Atkins } 474*97276fddSZach Atkins PetscCall(VecAXPY(DeltaX, 1, deltaX)); 475*97276fddSZach Atkins PetscCall(VecAXPY(X, 1, deltaX)); 476*97276fddSZach Atkins /* Q = -dF/dlambda(X, lambda)*/ 477*97276fddSZach Atkins PetscCall(SNESNewtonALComputeFunction(snes, X, Q)); 478*97276fddSZach Atkins /* R = F(X, lambda) */ 479*97276fddSZach Atkins PetscCall(SNESComputeFunction(snes, X, R)); 480*97276fddSZach Atkins PetscCall(VecNormBegin(R, NORM_2, &fnorm)); 481*97276fddSZach Atkins PetscCall(VecNormBegin(X, NORM_2, &xnorm)); 482*97276fddSZach Atkins PetscCall(VecNormBegin(deltaX, NORM_2, &ynorm)); 483*97276fddSZach Atkins PetscCall(VecNormEnd(R, NORM_2, &fnorm)); 484*97276fddSZach Atkins PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 485*97276fddSZach Atkins PetscCall(VecNormEnd(deltaX, NORM_2, &ynorm)); 486*97276fddSZach Atkins 487*97276fddSZach Atkins if (PetscLogPrintInfo) PetscCall(SNESNewtonALCheckArcLength(snes, DeltaX, data->lambda_update, stepSize)); 488*97276fddSZach Atkins 489*97276fddSZach Atkins /* Monitor convergence */ 490*97276fddSZach Atkins PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 491*97276fddSZach Atkins snes->iter++; 492*97276fddSZach Atkins snes->norm = fnorm; 493*97276fddSZach Atkins snes->ynorm = ynorm; 494*97276fddSZach Atkins snes->xnorm = xnorm; 495*97276fddSZach Atkins PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 496*97276fddSZach Atkins PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 497*97276fddSZach Atkins PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 498*97276fddSZach Atkins PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 499*97276fddSZach Atkins if (snes->reason) break; 500*97276fddSZach Atkins if (j == maxits - 1) { 501*97276fddSZach Atkins snes->reason = SNES_DIVERGED_MAX_IT; 502*97276fddSZach Atkins break; 503*97276fddSZach Atkins } 504*97276fddSZach Atkins } 505*97276fddSZach Atkins if (snes->reason < 0) break; 506*97276fddSZach Atkins if (data->lambda >= data->lambda_max) { 507*97276fddSZach Atkins break; 508*97276fddSZach Atkins } else if (maxincs > 0 && i == maxincs - 1) { 509*97276fddSZach Atkins snes->reason = SNES_DIVERGED_MAX_IT; 510*97276fddSZach Atkins break; 511*97276fddSZach Atkins } else { 512*97276fddSZach Atkins snes->reason = SNES_CONVERGED_ITERATING; 513*97276fddSZach Atkins } 514*97276fddSZach Atkins } 515*97276fddSZach Atkins /* Reset RHS vector, if changed */ 516*97276fddSZach Atkins if (data->copied_rhs) { 517*97276fddSZach Atkins PetscCall(VecCopy(data->vec_rhs_orig, snes->vec_rhs)); 518*97276fddSZach Atkins data->copied_rhs = PETSC_FALSE; 519*97276fddSZach Atkins } 520*97276fddSZach Atkins snes->max_its = maxits; /* reset snes->max_its */ 521*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 522*97276fddSZach Atkins } 523*97276fddSZach Atkins 524*97276fddSZach Atkins /* 525*97276fddSZach Atkins SNESSetUp_NEWTONAL - Sets up the internal data structures for the later use 526*97276fddSZach Atkins of the SNESNEWTONAL nonlinear solver. 527*97276fddSZach Atkins 528*97276fddSZach Atkins Input Parameter: 529*97276fddSZach Atkins . snes - the SNES context 530*97276fddSZach Atkins . x - the solution vector 531*97276fddSZach Atkins 532*97276fddSZach Atkins Application Interface Routine: SNESSetUp() 533*97276fddSZach Atkins */ 534*97276fddSZach Atkins static PetscErrorCode SNESSetUp_NEWTONAL(SNES snes) 535*97276fddSZach Atkins { 536*97276fddSZach Atkins PetscFunctionBegin; 537*97276fddSZach Atkins PetscCall(SNESSetWorkVecs(snes, 4)); 538*97276fddSZach Atkins PetscCall(SNESSetUpMatrices(snes)); 539*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 540*97276fddSZach Atkins } 541*97276fddSZach Atkins 542*97276fddSZach Atkins /* 543*97276fddSZach Atkins SNESSetFromOptions_NEWTONAL - Sets various parameters for the SNESNEWTONAL method. 544*97276fddSZach Atkins 545*97276fddSZach Atkins Input Parameter: 546*97276fddSZach Atkins . snes - the SNES context 547*97276fddSZach Atkins 548*97276fddSZach Atkins Application Interface Routine: SNESSetFromOptions() 549*97276fddSZach Atkins */ 550*97276fddSZach Atkins static PetscErrorCode SNESSetFromOptions_NEWTONAL(SNES snes, PetscOptionItems *PetscOptionsObject) 551*97276fddSZach Atkins { 552*97276fddSZach Atkins SNES_NEWTONAL *data = (SNES_NEWTONAL *)snes->data; 553*97276fddSZach Atkins SNESNewtonALCorrectionType correction_type = data->correction_type; 554*97276fddSZach Atkins 555*97276fddSZach Atkins PetscFunctionBegin; 556*97276fddSZach Atkins PetscOptionsHeadBegin(PetscOptionsObject, "SNES Newton Arc Length options"); 557*97276fddSZach Atkins PetscCall(PetscOptionsReal("-snes_newtonal_step_size", "Initial arc length increment step size", "SNESNewtonAL", data->step_size, &data->step_size, NULL)); 558*97276fddSZach Atkins PetscCall(PetscOptionsInt("-snes_newtonal_max_continuation_steps", "Maximum number of increment steps", "SNESNewtonAL", data->max_continuation_steps, &data->max_continuation_steps, NULL)); 559*97276fddSZach Atkins PetscCall(PetscOptionsReal("-snes_newtonal_psisq", "Regularization parameter for arc length continuation, 0 for cylindrical", "SNESNewtonAL", data->psisq, &data->psisq, NULL)); 560*97276fddSZach Atkins PetscCall(PetscOptionsReal("-snes_newtonal_lambda_min", "Minimum value of the load parameter lambda", "SNESNewtonAL", data->lambda_min, &data->lambda_min, NULL)); 561*97276fddSZach Atkins PetscCall(PetscOptionsReal("-snes_newtonal_lambda_max", "Maximum value of the load parameter lambda", "SNESNewtonAL", data->lambda_max, &data->lambda_max, NULL)); 562*97276fddSZach Atkins PetscCall(PetscOptionsBool("-snes_newtonal_scale_rhs", "Scale the constant vector passed to `SNESSolve` by the load parameter lambda", "SNESNewtonAL", data->scale_rhs, &data->scale_rhs, NULL)); 563*97276fddSZach Atkins PetscCall(PetscOptionsEnum("-snes_newtonal_correction_type", "Type of correction to use in the arc-length continuation method", "SNESNewtonALCorrectionType", SNESNewtonALCorrectionTypes, (PetscEnum)correction_type, (PetscEnum *)&correction_type, NULL)); 564*97276fddSZach Atkins PetscCall(SNESNewtonALSetCorrectionType(snes, correction_type)); 565*97276fddSZach Atkins PetscOptionsHeadEnd(); 566*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 567*97276fddSZach Atkins } 568*97276fddSZach Atkins 569*97276fddSZach Atkins static PetscErrorCode SNESReset_NEWTONAL(SNES snes) 570*97276fddSZach Atkins { 571*97276fddSZach Atkins SNES_NEWTONAL *al = (SNES_NEWTONAL *)snes->data; 572*97276fddSZach Atkins 573*97276fddSZach Atkins PetscFunctionBegin; 574*97276fddSZach Atkins PetscCall(VecDestroy(&al->vec_rhs_orig)); 575*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 576*97276fddSZach Atkins } 577*97276fddSZach Atkins 578*97276fddSZach Atkins /* 579*97276fddSZach Atkins SNESDestroy_NEWTONAL - Destroys the private SNES_NEWTONAL context that was created 580*97276fddSZach Atkins with SNESCreate_NEWTONAL(). 581*97276fddSZach Atkins 582*97276fddSZach Atkins Input Parameter: 583*97276fddSZach Atkins . snes - the SNES context 584*97276fddSZach Atkins 585*97276fddSZach Atkins Application Interface Routine: SNESDestroy() 586*97276fddSZach Atkins */ 587*97276fddSZach Atkins static PetscErrorCode SNESDestroy_NEWTONAL(SNES snes) 588*97276fddSZach Atkins { 589*97276fddSZach Atkins PetscFunctionBegin; 590*97276fddSZach Atkins PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALSetCorrectionType_C", NULL)); 591*97276fddSZach Atkins PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALGetLoadParameter_C", NULL)); 592*97276fddSZach Atkins PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALSetFunction_C", NULL)); 593*97276fddSZach Atkins PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALGetFunction_C", NULL)); 594*97276fddSZach Atkins PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALComputeFunction_C", NULL)); 595*97276fddSZach Atkins PetscCall(PetscFree(snes->data)); 596*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 597*97276fddSZach Atkins } 598*97276fddSZach Atkins 599*97276fddSZach Atkins /*MC 600*97276fddSZach Atkins SNESNEWTONAL - Newton based nonlinear solver that uses a arc-length continuation method to solve the nonlinear system. 601*97276fddSZach Atkins 602*97276fddSZach Atkins Options Database Keys: 603*97276fddSZach Atkins + -snes_newtonal_step_size <1.0> - Initial arc length increment step size 604*97276fddSZach Atkins . -snes_newtonal_max_continuation_steps <100> - Maximum number of continuation steps, or negative for no limit (not recommended) 605*97276fddSZach Atkins . -snes_newtonal_psisq <1.0> - Regularization parameter for arc length continuation, 0 for cylindrical. Larger values generally lead to more steps. 606*97276fddSZach Atkins . -snes_newtonal_lambda_min <0.0> - Minimum value of the load parameter lambda 607*97276fddSZach Atkins . -snes_newtonal_lambda_max <1.0> - Maximum value of the load parameter lambda 608*97276fddSZach Atkins . -snes_newtonal_scale_rhs <true> - Scale the constant vector passed to `SNESSolve` by the load parameter lambda 609*97276fddSZach Atkins - -snes_newtonal_correction_type <exact> - Type of correction to use in the arc-length continuation method, `exact` or `normal` 610*97276fddSZach Atkins 611*97276fddSZach Atkins Level: intermediate 612*97276fddSZach Atkins 613*97276fddSZach Atkins Note: 614*97276fddSZach Atkins The exact correction scheme with partial updates is detailed in {cite}`Ritto-CorreaCamotim2008` and the implementation of the 615*97276fddSZach Atkins normal correction scheme is based on {cite}`LeonPaulinoPereiraMenezesLages_2011`. 616*97276fddSZach Atkins 617*97276fddSZach Atkins .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONAL`, `SNESNewtonALSetFunction()`, `SNESNewtonALGetFunction()`, `SNESNewtonALGetLoadParameter()` 618*97276fddSZach Atkins M*/ 619*97276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONAL(SNES snes) 620*97276fddSZach Atkins { 621*97276fddSZach Atkins SNES_NEWTONAL *arclengthParameters; 622*97276fddSZach Atkins 623*97276fddSZach Atkins PetscFunctionBegin; 624*97276fddSZach Atkins snes->ops->setup = SNESSetUp_NEWTONAL; 625*97276fddSZach Atkins snes->ops->solve = SNESSolve_NEWTONAL; 626*97276fddSZach Atkins snes->ops->destroy = SNESDestroy_NEWTONAL; 627*97276fddSZach Atkins snes->ops->setfromoptions = SNESSetFromOptions_NEWTONAL; 628*97276fddSZach Atkins snes->ops->reset = SNESReset_NEWTONAL; 629*97276fddSZach Atkins 630*97276fddSZach Atkins snes->usesksp = PETSC_TRUE; 631*97276fddSZach Atkins snes->usesnpc = PETSC_FALSE; 632*97276fddSZach Atkins 633*97276fddSZach Atkins PetscCall(SNESParametersInitialize(snes)); 634*97276fddSZach Atkins PetscObjectParameterSetDefault(snes, max_funcs, 30000); 635*97276fddSZach Atkins PetscObjectParameterSetDefault(snes, max_its, 50); 636*97276fddSZach Atkins 637*97276fddSZach Atkins snes->alwayscomputesfinalresidual = PETSC_TRUE; 638*97276fddSZach Atkins 639*97276fddSZach Atkins PetscCall(PetscNew(&arclengthParameters)); 640*97276fddSZach Atkins arclengthParameters->lambda = 0.0; 641*97276fddSZach Atkins arclengthParameters->lambda_update = 0.0; 642*97276fddSZach Atkins arclengthParameters->step_size = 1.0; 643*97276fddSZach Atkins arclengthParameters->max_continuation_steps = 100; 644*97276fddSZach Atkins arclengthParameters->psisq = 1.0; 645*97276fddSZach Atkins arclengthParameters->lambda_min = 0.0; 646*97276fddSZach Atkins arclengthParameters->lambda_max = 1.0; 647*97276fddSZach Atkins arclengthParameters->scale_rhs = PETSC_TRUE; 648*97276fddSZach Atkins arclengthParameters->correction_type = SNES_NEWTONAL_CORRECTION_EXACT; 649*97276fddSZach Atkins snes->data = (void *)arclengthParameters; 650*97276fddSZach Atkins 651*97276fddSZach Atkins PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALSetCorrectionType_C", SNESNewtonALSetCorrectionType_NEWTONAL)); 652*97276fddSZach Atkins PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALGetLoadParameter_C", SNESNewtonALGetLoadParameter_NEWTONAL)); 653*97276fddSZach Atkins PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALSetFunction_C", SNESNewtonALSetFunction_NEWTONAL)); 654*97276fddSZach Atkins PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALGetFunction_C", SNESNewtonALGetFunction_NEWTONAL)); 655*97276fddSZach Atkins PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALComputeFunction_C", SNESNewtonALComputeFunction_NEWTONAL)); 656*97276fddSZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 657*97276fddSZach Atkins } 658