xref: /petsc/src/snes/impls/al/al.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
197276fddSZach Atkins #include <../src/snes/impls/al/alimpl.h> /*I "petscsnes.h" I*/
297276fddSZach Atkins 
397276fddSZach Atkins /*
497276fddSZach Atkins      This file implements a truncated Newton method with arc length continuation,
597276fddSZach Atkins      for solving a system of nonlinear equations, using the KSP, Vec,
697276fddSZach Atkins      and Mat interfaces for linear solvers, vectors, and matrices,
797276fddSZach Atkins      respectively.
897276fddSZach Atkins */
997276fddSZach Atkins 
1097276fddSZach Atkins const char NewtonALExactCitation[]   = "@article{Ritto-CorreaCamotim2008,\n"
1197276fddSZach Atkins                                        "  title={On the arc-length and other quadratic control methods: Established, less known and new implementation procedures},\n"
1297276fddSZach Atkins                                        "  volume={86},\n"
1397276fddSZach Atkins                                        "  ISSN={0045-7949},\n"
1497276fddSZach Atkins                                        "  DOI={10.1016/j.compstruc.2007.08.003},\n"
1597276fddSZach Atkins                                        "  number={11},\n"
1697276fddSZach Atkins                                        "  journal={Computers & Structures},\n"
1797276fddSZach Atkins                                        "  author={Ritto-Corr{\\^{e}}a, Manuel and Camotim, Dinar},\n"
1897276fddSZach Atkins                                        "  year={2008},\n"
1997276fddSZach Atkins                                        "  month=jun,\n"
2097276fddSZach Atkins                                        "  pages={1353-1368},\n"
2197276fddSZach Atkins                                        "}\n";
2297276fddSZach Atkins PetscBool  NewtonALExactCitationSet  = PETSC_FALSE;
2397276fddSZach Atkins const char NewtonALNormalCitation[]  = "@article{LeonPaulinoPereiraMenezesLages_2011,\n"
2497276fddSZach Atkins                                        "  title={A Unified Library of Nonlinear Solution Schemes},\n"
2597276fddSZach Atkins                                        "  volume={64},\n"
2697276fddSZach Atkins                                        "  ISSN={0003-6900, 2379-0407},\n"
2797276fddSZach Atkins                                        "  DOI={10.1115/1.4006992},\n"
2897276fddSZach Atkins                                        "  number={4},\n"
2997276fddSZach Atkins                                        "  journal={Applied Mechanics Reviews},\n"
3097276fddSZach Atkins                                        "  author={Leon, Sofie E. and Paulino, Glaucio H. and Pereira, Anderson and Menezes, Ivan F. M. and Lages, Eduardo N.},\n"
3197276fddSZach Atkins                                        "  year={2011},\n"
3297276fddSZach Atkins                                        "  month=jul,\n"
3397276fddSZach Atkins                                        "  pages={040803},\n"
3497276fddSZach Atkins                                        "  language={en}\n"
3597276fddSZach Atkins                                        "}\n";
3697276fddSZach Atkins PetscBool  NewtonALNormalCitationSet = PETSC_FALSE;
3797276fddSZach Atkins 
3897276fddSZach Atkins const char *const SNESNewtonALCorrectionTypes[] = {"EXACT", "NORMAL", "SNESNewtonALCorrectionType", "SNES_NEWTONAL_CORRECTION_", NULL};
3997276fddSZach Atkins 
SNESNewtonALCheckArcLength(SNES snes,Vec XStep,PetscReal lambdaStep,PetscReal stepSize)4097276fddSZach Atkins static PetscErrorCode SNESNewtonALCheckArcLength(SNES snes, Vec XStep, PetscReal lambdaStep, PetscReal stepSize)
4197276fddSZach Atkins {
4297276fddSZach Atkins   PetscReal      arcLength, arcLengthError;
4397276fddSZach Atkins   SNES_NEWTONAL *al = (SNES_NEWTONAL *)snes->data;
4497276fddSZach Atkins 
4597276fddSZach Atkins   PetscFunctionBegin;
4697276fddSZach Atkins   PetscCall(VecNorm(XStep, NORM_2, &arcLength));
4797276fddSZach Atkins   arcLength      = PetscSqrtReal(PetscSqr(arcLength) + al->psisq * lambdaStep * lambdaStep);
4897276fddSZach Atkins   arcLengthError = PetscAbsReal(arcLength - stepSize);
4997276fddSZach Atkins 
5097276fddSZach 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));
5197276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
5297276fddSZach Atkins }
5397276fddSZach Atkins 
5497276fddSZach Atkins /* stable implementation of roots of a*x^2 + b*x + c = 0 */
PetscQuadraticRoots(PetscReal a,PetscReal b,PetscReal c,PetscReal * xm,PetscReal * xp)5597276fddSZach Atkins static inline void PetscQuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal *xm, PetscReal *xp)
5697276fddSZach Atkins {
5797276fddSZach Atkins   PetscReal temp = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(b * b - 4 * a * c));
5897276fddSZach Atkins   PetscReal x1   = temp / a;
5997276fddSZach Atkins   PetscReal x2   = c / temp;
6097276fddSZach Atkins   *xm            = PetscMin(x1, x2);
6197276fddSZach Atkins   *xp            = PetscMax(x1, x2);
6297276fddSZach Atkins }
6397276fddSZach Atkins 
SNESNewtonALSetCorrectionType_NEWTONAL(SNES snes,SNESNewtonALCorrectionType ctype)6497276fddSZach Atkins static PetscErrorCode SNESNewtonALSetCorrectionType_NEWTONAL(SNES snes, SNESNewtonALCorrectionType ctype)
6597276fddSZach Atkins {
6697276fddSZach Atkins   SNES_NEWTONAL *al = (SNES_NEWTONAL *)snes->data;
6797276fddSZach Atkins 
6897276fddSZach Atkins   PetscFunctionBegin;
6997276fddSZach Atkins   al->correction_type = ctype;
7097276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
7197276fddSZach Atkins }
7297276fddSZach Atkins 
7397276fddSZach Atkins /*@
7497276fddSZach Atkins   SNESNewtonALSetCorrectionType - Set the type of correction to use in the arc-length continuation method.
7597276fddSZach Atkins 
7697276fddSZach Atkins   Logically Collective
7797276fddSZach Atkins 
7897276fddSZach Atkins   Input Parameters:
7997276fddSZach Atkins + snes  - the nonlinear solver object
8097276fddSZach Atkins - ctype - the type of correction to use
8197276fddSZach Atkins 
8297276fddSZach Atkins   Options Database Key:
8397276fddSZach Atkins . -snes_newtonal_correction_type <type> - Set the type of correction to use; use -help for a list of available types
8497276fddSZach Atkins 
8597276fddSZach Atkins   Level: intermediate
8697276fddSZach Atkins 
8797276fddSZach Atkins .seealso: [](ch_snes), `SNES`, `SNESNEWTONAL`, `SNESNewtonALCorrectionType`
8897276fddSZach Atkins @*/
SNESNewtonALSetCorrectionType(SNES snes,SNESNewtonALCorrectionType ctype)8997276fddSZach Atkins PetscErrorCode SNESNewtonALSetCorrectionType(SNES snes, SNESNewtonALCorrectionType ctype)
9097276fddSZach Atkins {
9197276fddSZach Atkins   PetscFunctionBegin;
9297276fddSZach Atkins   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
9397276fddSZach Atkins   PetscValidLogicalCollectiveEnum(snes, ctype, 2);
9497276fddSZach Atkins   PetscTryMethod(snes, "SNESNewtonALSetCorrectionType_C", (SNES, SNESNewtonALCorrectionType), (snes, ctype));
9597276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
9697276fddSZach Atkins }
9797276fddSZach Atkins 
SNESNewtonALSetFunction_NEWTONAL(SNES snes,SNESFunctionFn * func,PetscCtx ctx)98*2a8381b2SBarry Smith static PetscErrorCode SNESNewtonALSetFunction_NEWTONAL(SNES snes, SNESFunctionFn *func, PetscCtx ctx)
9997276fddSZach Atkins {
10097276fddSZach Atkins   SNES_NEWTONAL *al = (SNES_NEWTONAL *)snes->data;
10197276fddSZach Atkins 
10297276fddSZach Atkins   PetscFunctionBegin;
10397276fddSZach Atkins   al->computealfunction = func;
10497276fddSZach Atkins   al->alctx             = ctx;
10597276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
10697276fddSZach Atkins }
10797276fddSZach Atkins 
10897276fddSZach Atkins /*@C
10997276fddSZach Atkins   SNESNewtonALSetFunction - Sets a user function that is called at each function evaluation to
11097276fddSZach Atkins   compute the tangent load vector for the arc-length continuation method.
11197276fddSZach Atkins 
11297276fddSZach Atkins   Logically Collective
11397276fddSZach Atkins 
11497276fddSZach Atkins   Input Parameters:
11597276fddSZach Atkins + snes - the nonlinear solver object
11697276fddSZach 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
11797276fddSZach Atkins - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
11897276fddSZach Atkins 
11997276fddSZach Atkins   Level: intermediate
12097276fddSZach Atkins 
12197276fddSZach Atkins   Notes:
12297276fddSZach Atkins   If the current value of the load parameter is needed in `func`, it can be obtained with `SNESNewtonALGetLoadParameter()`.
12397276fddSZach Atkins 
12497276fddSZach Atkins   The tangent load vector is the partial derivative of external load with respect to the load parameter.
12597276fddSZach Atkins   In the case of proportional loading, the tangent load vector is the full external load vector at the end of the load step.
12697276fddSZach Atkins 
12797276fddSZach Atkins .seealso: [](ch_snes), `SNES`, `SNESNEWTONAL`, `SNESNewtonALGetFunction()`, `SNESNewtonALGetLoadParameter()`
12897276fddSZach Atkins @*/
SNESNewtonALSetFunction(SNES snes,SNESFunctionFn * func,PetscCtx ctx)129*2a8381b2SBarry Smith PetscErrorCode SNESNewtonALSetFunction(SNES snes, SNESFunctionFn *func, PetscCtx ctx)
13097276fddSZach Atkins {
13197276fddSZach Atkins   PetscFunctionBegin;
13297276fddSZach Atkins   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
13397276fddSZach Atkins   PetscTryMethod(snes, "SNESNewtonALSetFunction_C", (SNES, SNESFunctionFn *, void *), (snes, func, ctx));
13497276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
13597276fddSZach Atkins }
13697276fddSZach Atkins 
SNESNewtonALGetFunction_NEWTONAL(SNES snes,SNESFunctionFn ** func,PetscCtxRt ctx)137*2a8381b2SBarry Smith static PetscErrorCode SNESNewtonALGetFunction_NEWTONAL(SNES snes, SNESFunctionFn **func, PetscCtxRt ctx)
13897276fddSZach Atkins {
13997276fddSZach Atkins   SNES_NEWTONAL *al = (SNES_NEWTONAL *)snes->data;
14097276fddSZach Atkins 
14197276fddSZach Atkins   PetscFunctionBegin;
14297276fddSZach Atkins   if (func) *func = al->computealfunction;
143*2a8381b2SBarry Smith   if (ctx) *(void **)ctx = al->alctx;
14497276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
14597276fddSZach Atkins }
14697276fddSZach Atkins 
14797276fddSZach Atkins /*@C
14897276fddSZach Atkins   SNESNewtonALGetFunction - Get the user function and context set with `SNESNewtonALSetFunction`
14997276fddSZach Atkins 
15097276fddSZach Atkins   Logically Collective
15197276fddSZach Atkins 
15297276fddSZach Atkins   Input Parameters:
15397276fddSZach Atkins + snes - the nonlinear solver object
15497276fddSZach Atkins . func - [optional] tangent load function evaluation routine, see `SNESNewtonALSetFunction()` for the call sequence
15597276fddSZach Atkins - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
15697276fddSZach Atkins 
15797276fddSZach Atkins   Level: intermediate
15897276fddSZach Atkins 
15997276fddSZach Atkins .seealso: [](ch_snes), `SNES`, `SNESNEWTONAL`, `SNESNewtonALSetFunction()`
16097276fddSZach Atkins @*/
SNESNewtonALGetFunction(SNES snes,SNESFunctionFn ** func,PetscCtxRt ctx)161*2a8381b2SBarry Smith PetscErrorCode SNESNewtonALGetFunction(SNES snes, SNESFunctionFn **func, PetscCtxRt ctx)
16297276fddSZach Atkins {
16397276fddSZach Atkins   PetscFunctionBegin;
16497276fddSZach Atkins   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
165*2a8381b2SBarry Smith   PetscUseMethod(snes, "SNESNewtonALGetFunction_C", (SNES, SNESFunctionFn **, PetscCtxRt), (snes, func, ctx));
16697276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
16797276fddSZach Atkins }
16897276fddSZach Atkins 
SNESNewtonALGetLoadParameter_NEWTONAL(SNES snes,PetscReal * lambda)16997276fddSZach Atkins static PetscErrorCode SNESNewtonALGetLoadParameter_NEWTONAL(SNES snes, PetscReal *lambda)
17097276fddSZach Atkins {
17197276fddSZach Atkins   SNES_NEWTONAL *al;
17297276fddSZach Atkins 
17397276fddSZach Atkins   PetscFunctionBeginHot;
17497276fddSZach Atkins   al      = (SNES_NEWTONAL *)snes->data;
17597276fddSZach Atkins   *lambda = al->lambda;
17697276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
17797276fddSZach Atkins }
17897276fddSZach Atkins 
17997276fddSZach Atkins /*@C
18097276fddSZach Atkins   SNESNewtonALGetLoadParameter - Get the value of the load parameter `lambda` for the arc-length continuation method.
18197276fddSZach Atkins 
18297276fddSZach Atkins   Logically Collective
18397276fddSZach Atkins 
18497276fddSZach Atkins   Input Parameter:
18597276fddSZach Atkins . snes - the nonlinear solver object
18697276fddSZach Atkins 
18797276fddSZach Atkins   Output Parameter:
18897276fddSZach Atkins . lambda - the arc-length parameter
18997276fddSZach Atkins 
19097276fddSZach Atkins   Level: intermediate
19197276fddSZach Atkins 
19297276fddSZach Atkins   Notes:
19397276fddSZach Atkins   This function should be used in the functions provided to `SNESSetFunction()` and `SNESNewtonALSetFunction()`
19497276fddSZach Atkins   to compute the residual and tangent load vectors for a given value of `lambda` (0 <= lambda <= 1).
19597276fddSZach Atkins 
19697276fddSZach Atkins   Usually, `lambda` is used to scale the external force vector in the residual function, i.e. proportional loading,
19797276fddSZach Atkins   in which case the tangent load vector is the full external force vector.
19897276fddSZach Atkins 
19997276fddSZach Atkins .seealso: [](ch_snes), `SNES`, `SNESNEWTONAL`, `SNESNewtonALSetFunction()`
20097276fddSZach Atkins @*/
SNESNewtonALGetLoadParameter(SNES snes,PetscReal * lambda)20197276fddSZach Atkins PetscErrorCode SNESNewtonALGetLoadParameter(SNES snes, PetscReal *lambda)
20297276fddSZach Atkins {
20397276fddSZach Atkins   PetscFunctionBeginHot;
20497276fddSZach Atkins   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
20597276fddSZach Atkins   PetscAssertPointer(lambda, 2);
20697276fddSZach Atkins   PetscUseMethod(snes, "SNESNewtonALGetLoadParameter_C", (SNES, PetscReal *), (snes, lambda));
20797276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
20897276fddSZach Atkins }
20997276fddSZach Atkins 
SNESNewtonALComputeFunction_NEWTONAL(SNES snes,Vec X,Vec Q)21097276fddSZach Atkins static PetscErrorCode SNESNewtonALComputeFunction_NEWTONAL(SNES snes, Vec X, Vec Q)
21197276fddSZach Atkins {
21297276fddSZach Atkins   void           *ctx               = NULL;
21397276fddSZach Atkins   SNESFunctionFn *computealfunction = NULL;
21497276fddSZach Atkins   SNES_NEWTONAL  *al;
21597276fddSZach Atkins 
21697276fddSZach Atkins   PetscFunctionBegin;
21797276fddSZach Atkins   al = (SNES_NEWTONAL *)snes->data;
21897276fddSZach Atkins   PetscCall(SNESNewtonALGetFunction(snes, &computealfunction, &ctx));
21997276fddSZach Atkins 
22097276fddSZach Atkins   PetscCall(VecZeroEntries(Q));
22197276fddSZach 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");
22297276fddSZach Atkins   if (computealfunction) {
22397276fddSZach Atkins     PetscCall(VecLockReadPush(X));
22497276fddSZach Atkins     PetscCallBack("SNES callback NewtonAL tangent load function", (*computealfunction)(snes, X, Q, ctx));
22597276fddSZach Atkins     PetscCall(VecLockReadPop(X));
22697276fddSZach Atkins   }
22797276fddSZach Atkins   if (al->scale_rhs && snes->vec_rhs) {
22897276fddSZach Atkins     /* Save original RHS vector values, then scale `snes->vec_rhs` by load parameter */
2293a7d0413SPierre Jolivet     if (!al->vec_rhs_orig) PetscCall(VecDuplicate(snes->vec_rhs, &al->vec_rhs_orig));
23097276fddSZach Atkins     if (!al->copied_rhs) {
23197276fddSZach Atkins       PetscCall(VecCopy(snes->vec_rhs, al->vec_rhs_orig));
23297276fddSZach Atkins       al->copied_rhs = PETSC_TRUE;
23397276fddSZach Atkins     }
23497276fddSZach Atkins     PetscCall(VecAXPBY(snes->vec_rhs, al->lambda, 0.0, al->vec_rhs_orig));
23597276fddSZach Atkins     PetscCall(VecAXPY(Q, 1, al->vec_rhs_orig));
23697276fddSZach Atkins   }
23797276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
23897276fddSZach Atkins }
23997276fddSZach Atkins 
24097276fddSZach Atkins /*@C
24197276fddSZach Atkins   SNESNewtonALComputeFunction - Calls the function that has been set with `SNESNewtonALSetFunction()`.
24297276fddSZach Atkins 
24397276fddSZach Atkins   Collective
24497276fddSZach Atkins 
24597276fddSZach Atkins   Input Parameters:
24697276fddSZach Atkins + snes - the `SNES` context
24797276fddSZach Atkins - X    - input vector
24897276fddSZach Atkins 
24997276fddSZach Atkins   Output Parameter:
25097276fddSZach Atkins . Q - tangent load vector, as set by `SNESNewtonALSetFunction()`
25197276fddSZach Atkins 
25297276fddSZach Atkins   Level: developer
25397276fddSZach Atkins 
25497276fddSZach Atkins   Note:
25597276fddSZach Atkins   `SNESNewtonALComputeFunction()` is typically used within nonlinear solvers
25697276fddSZach Atkins   implementations, so users would not generally call this routine themselves.
25797276fddSZach Atkins 
25897276fddSZach Atkins .seealso: [](ch_snes), `SNES`, `SNESNewtonALSetFunction()`, `SNESNewtonALGetFunction()`
25997276fddSZach Atkins @*/
SNESNewtonALComputeFunction(SNES snes,Vec X,Vec Q)26097276fddSZach Atkins PetscErrorCode SNESNewtonALComputeFunction(SNES snes, Vec X, Vec Q)
26197276fddSZach Atkins {
26297276fddSZach Atkins   PetscFunctionBegin;
26397276fddSZach Atkins   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
26497276fddSZach Atkins   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
26597276fddSZach Atkins   PetscValidHeaderSpecific(Q, VEC_CLASSID, 3);
26697276fddSZach Atkins   PetscCheckSameComm(snes, 1, X, 2);
26797276fddSZach Atkins   PetscCheckSameComm(snes, 1, Q, 3);
26897276fddSZach Atkins   PetscCall(VecValidValues_Internal(X, 2, PETSC_TRUE));
26997276fddSZach Atkins   PetscCall(PetscLogEventBegin(SNES_NewtonALEval, snes, X, Q, 0));
27097276fddSZach Atkins   PetscTryMethod(snes, "SNESNewtonALComputeFunction_C", (SNES, Vec, Vec), (snes, X, Q));
27197276fddSZach Atkins   PetscCall(PetscLogEventEnd(SNES_NewtonALEval, snes, X, Q, 0));
27297276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
27397276fddSZach Atkins }
27497276fddSZach Atkins 
27597276fddSZach Atkins /*
27697276fddSZach Atkins   SNESSolve_NEWTONAL - Solves a nonlinear system with Newton's method with arc length continuation.
27797276fddSZach Atkins 
27897276fddSZach Atkins   Input Parameter:
27997276fddSZach Atkins . snes - the `SNES` context
28097276fddSZach Atkins 
28197276fddSZach Atkins   Application Interface Routine: SNESSolve()
28297276fddSZach Atkins */
SNESSolve_NEWTONAL(SNES snes)28397276fddSZach Atkins static PetscErrorCode SNESSolve_NEWTONAL(SNES snes)
28497276fddSZach Atkins {
28597276fddSZach Atkins   SNES_NEWTONAL *data = (SNES_NEWTONAL *)snes->data;
28697276fddSZach Atkins   PetscInt       maxits, maxincs, lits;
28797276fddSZach Atkins   PetscReal      fnorm, xnorm, ynorm, stepSize;
28897276fddSZach Atkins   Vec            DeltaX, deltaX, X, R, Q, deltaX_Q, deltaX_R;
28997276fddSZach Atkins 
29097276fddSZach Atkins   PetscFunctionBegin;
29197276fddSZach 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);
29297276fddSZach Atkins 
29397276fddSZach Atkins   /* Register citations */
29497276fddSZach Atkins   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
29597276fddSZach Atkins   if (data->correction_type == SNES_NEWTONAL_CORRECTION_EXACT) {
29697276fddSZach Atkins     PetscCall(PetscCitationsRegister(NewtonALExactCitation, &NewtonALExactCitationSet));
29797276fddSZach Atkins   } else if (data->correction_type == SNES_NEWTONAL_CORRECTION_NORMAL) {
29897276fddSZach Atkins     PetscCall(PetscCitationsRegister(NewtonALNormalCitation, &NewtonALNormalCitationSet));
29997276fddSZach Atkins   }
30097276fddSZach Atkins 
30197276fddSZach Atkins   data->copied_rhs             = PETSC_FALSE;
30297276fddSZach Atkins   data->lambda_update          = 0.0;
30397276fddSZach Atkins   data->lambda                 = 0.0;
30497276fddSZach Atkins   snes->numFailures            = 0;
30597276fddSZach Atkins   snes->numLinearSolveFailures = 0;
30697276fddSZach Atkins   snes->reason                 = SNES_CONVERGED_ITERATING;
3079350dab8SStefano Zampini   snes->iter                   = 0;
30897276fddSZach Atkins 
30997276fddSZach Atkins   maxits   = snes->max_its;                /* maximum number of iterations */
31097276fddSZach Atkins   maxincs  = data->max_continuation_steps; /* maximum number of increments */
31197276fddSZach Atkins   X        = snes->vec_sol;                /* solution vector */
31297276fddSZach Atkins   R        = snes->vec_func;               /* residual vector */
31397276fddSZach Atkins   Q        = snes->work[0];                /* tangent load vector */
31497276fddSZach Atkins   deltaX_Q = snes->work[1];                /* variation of X with respect to lambda */
31597276fddSZach Atkins   deltaX_R = snes->work[2];                /* linearized error correction */
31697276fddSZach Atkins   DeltaX   = snes->work[3];                /* step from equilibrium */
31797276fddSZach Atkins   deltaX   = snes->vec_sol_update;         /* full newton step */
31897276fddSZach Atkins   stepSize = data->step_size;              /* initial step size */
31997276fddSZach Atkins 
32097276fddSZach Atkins   PetscCall(VecZeroEntries(DeltaX));
32197276fddSZach Atkins 
32297276fddSZach Atkins   /* set snes->max_its for convergence test */
32397276fddSZach Atkins   snes->max_its = maxits * maxincs;
32497276fddSZach Atkins 
32597276fddSZach Atkins   /* main incremental-iterative loop */
32697276fddSZach Atkins   for (PetscInt i = 0; i < maxincs || maxincs < 0; i++) {
32797276fddSZach Atkins     PetscReal deltaLambda;
32897276fddSZach Atkins 
32997276fddSZach Atkins     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
33097276fddSZach Atkins     snes->norm = 0.0;
33197276fddSZach Atkins     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
33297276fddSZach Atkins     PetscCall(SNESNewtonALComputeFunction(snes, X, Q));
33397276fddSZach Atkins     PetscCall(SNESComputeFunction(snes, X, R));
33497276fddSZach Atkins     PetscCall(VecAXPY(R, 1, Q));           /* R <- R + Q */
33597276fddSZach Atkins     PetscCall(VecNorm(R, NORM_2, &fnorm)); /* fnorm <- ||R|| */
33676c63389SBarry Smith     SNESCheckFunctionDomainError(snes, fnorm);
33797276fddSZach Atkins 
33897276fddSZach Atkins     /* Monitor convergence */
33997276fddSZach Atkins     PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm));
34097276fddSZach Atkins     PetscCall(SNESMonitor(snes, snes->iter, fnorm));
34197276fddSZach Atkins     if (i == 0 && snes->reason) break;
34297276fddSZach Atkins     for (PetscInt j = 0; j < maxits; j++) {
34397276fddSZach Atkins       PetscReal normsqX_Q, deltaS = 1;
34497276fddSZach Atkins 
34597276fddSZach Atkins       /* Call general purpose update function */
34697276fddSZach Atkins       PetscTryTypeMethod(snes, update, snes->iter);
34797276fddSZach Atkins 
34897276fddSZach Atkins       PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
34976c63389SBarry Smith       SNESCheckJacobianDomainError(snes);
35097276fddSZach Atkins       PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
35197276fddSZach Atkins       /* Solve J deltaX_Q = Q, where J is Jacobian matrix */
35297276fddSZach Atkins       PetscCall(KSPSolve(snes->ksp, Q, deltaX_Q));
35397276fddSZach Atkins       SNESCheckKSPSolve(snes);
35497276fddSZach Atkins       PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
35597276fddSZach Atkins       PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", tangent load linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
35697276fddSZach Atkins       /* Compute load parameter variation */
35797276fddSZach Atkins       PetscCall(VecNorm(deltaX_Q, NORM_2, &normsqX_Q));
35897276fddSZach Atkins       normsqX_Q *= normsqX_Q;
35997276fddSZach Atkins       /* On first iter, use predictor. This is the same regardless of corrector scheme. */
36097276fddSZach Atkins       if (j == 0) {
36197276fddSZach Atkins         PetscReal sign = 1.0;
36297276fddSZach Atkins         if (i > 0) {
36397276fddSZach Atkins           PetscCall(VecDotRealPart(DeltaX, deltaX_Q, &sign));
36497276fddSZach Atkins           sign += data->psisq * data->lambda_update;
36597276fddSZach Atkins           sign = sign >= 0 ? 1.0 : -1.0;
36697276fddSZach Atkins         }
36797276fddSZach Atkins         data->lambda_update = 0.0;
36897276fddSZach Atkins         PetscCall(VecZeroEntries(DeltaX));
36997276fddSZach Atkins         deltaLambda = sign * stepSize / PetscSqrtReal(normsqX_Q + data->psisq);
37097276fddSZach Atkins       } else {
37197276fddSZach Atkins         /* Solve J deltaX_R = -R */
37297276fddSZach Atkins         PetscCall(KSPSolve(snes->ksp, R, deltaX_R));
37397276fddSZach Atkins         SNESCheckKSPSolve(snes);
37497276fddSZach Atkins         PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
37597276fddSZach Atkins         PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", residual linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
37697276fddSZach Atkins         PetscCall(VecScale(deltaX_R, -1));
37797276fddSZach Atkins 
37897276fddSZach Atkins         if (data->correction_type == SNES_NEWTONAL_CORRECTION_NORMAL) {
37997276fddSZach Atkins           /*
38097276fddSZach Atkins             Take a step orthogonal to the current incremental update DeltaX.
38197276fddSZach Atkins             Note, this approach is cheaper than the exact correction, but may exhibit convergence
382d7c1f440SPierre Jolivet             issues due to the iterative trial points not being on the quadratic constraint surface.
38397276fddSZach Atkins             On the bright side, we always have a real and unique solution for deltaLambda.
38497276fddSZach Atkins           */
38597276fddSZach Atkins           PetscScalar coefs[2];
38697276fddSZach Atkins           Vec         rhs[] = {deltaX_R, deltaX_Q};
38797276fddSZach Atkins 
38897276fddSZach Atkins           PetscCall(VecMDot(DeltaX, 2, rhs, coefs));
38997276fddSZach Atkins           deltaLambda = -PetscRealPart(coefs[0]) / (PetscRealPart(coefs[1]) + data->psisq * data->lambda_update);
39097276fddSZach Atkins         } else {
39197276fddSZach Atkins           /*
39297276fddSZach Atkins             Solve
39397276fddSZach Atkins               a*deltaLambda^2 + b*deltaLambda + c = 0  (*)
39497276fddSZach Atkins             where
39597276fddSZach Atkins               a = a0
39697276fddSZach Atkins               b = b0 + b1*deltaS
39797276fddSZach Atkins               c = c0 + c1*deltaS + c2*deltaS^2
39897276fddSZach Atkins             and deltaS is either 1, or the largest value in (0, 1) that satisfies
39997276fddSZach Atkins               b^2 - 4*a*c = as*deltaS^2 + bs*deltaS + cs >= 0
40097276fddSZach Atkins             where
40197276fddSZach Atkins               as = b1^2 - 4*a0*c2
40297276fddSZach Atkins               bs = 2*b1*b0 - 4*a0*c1
40397276fddSZach Atkins               cs = b0^2 - 4*a0*c0
40497276fddSZach Atkins             These "partial corrections" prevent (*) from having complex roots.
40597276fddSZach Atkins           */
40697276fddSZach Atkins           PetscReal   psisqLambdaUpdate, discriminant;
40797276fddSZach Atkins           PetscReal   as, bs, cs;
40897276fddSZach Atkins           PetscReal   a0, b0, b1, c0, c1, c2;
40997276fddSZach Atkins           PetscScalar coefs1[3]; /* coefs[0] = deltaX_Q*DeltaX, coefs[1] = deltaX_R*DeltaX, coefs[2] = DeltaX*DeltaX */
41097276fddSZach Atkins           PetscScalar coefs2[2]; /* coefs[0] = deltaX_Q*deltaX_R, coefs[1] = deltaX_R*deltaX_R */
41197276fddSZach Atkins           const Vec   rhs1[3] = {deltaX_Q, deltaX_R, DeltaX};
41297276fddSZach Atkins           const Vec   rhs2[2] = {deltaX_Q, deltaX_R};
41397276fddSZach Atkins 
41497276fddSZach Atkins           psisqLambdaUpdate = data->psisq * data->lambda_update;
41597276fddSZach Atkins           PetscCall(VecMDotBegin(DeltaX, 3, rhs1, coefs1));
41697276fddSZach Atkins           PetscCall(VecMDotBegin(deltaX_R, 2, rhs2, coefs2));
41797276fddSZach Atkins           PetscCall(VecMDotEnd(DeltaX, 3, rhs1, coefs1));
41897276fddSZach Atkins           PetscCall(VecMDotEnd(deltaX_R, 2, rhs2, coefs2));
41997276fddSZach Atkins 
42097276fddSZach Atkins           a0 = normsqX_Q + data->psisq;
42197276fddSZach Atkins           b0 = 2 * (PetscRealPart(coefs1[0]) + psisqLambdaUpdate);
42297276fddSZach Atkins           b1 = 2 * PetscRealPart(coefs2[0]);
42397276fddSZach Atkins           c0 = PetscRealPart(coefs1[2]) + psisqLambdaUpdate * data->lambda_update - stepSize * stepSize;
42497276fddSZach Atkins           c1 = 2 * PetscRealPart(coefs1[1]);
42597276fddSZach Atkins           c2 = PetscRealPart(coefs2[1]);
42697276fddSZach Atkins 
42797276fddSZach Atkins           as = b1 * b1 - 4 * a0 * c2;
42897276fddSZach Atkins           bs = 2 * (b1 * b0 - 2 * a0 * c1);
42997276fddSZach Atkins           cs = b0 * b0 - 4 * a0 * c0;
43097276fddSZach Atkins 
43197276fddSZach Atkins           discriminant = cs + bs * deltaS + as * deltaS * deltaS;
43297276fddSZach Atkins 
43397276fddSZach Atkins           if (discriminant < 0) {
43497276fddSZach Atkins             /* Take deltaS < 1 with the unique root -b/(2*a) */
43597276fddSZach Atkins             PetscReal x1;
43697276fddSZach Atkins 
43797276fddSZach Atkins             /* Compute deltaS to be the largest root of (as * x^2 + bs * x + cs = 0) */
43897276fddSZach Atkins             PetscQuadraticRoots(as, bs, cs, &x1, &deltaS);
43997276fddSZach 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));
44097276fddSZach Atkins             deltaLambda = -0.5 * (b0 + b1 * deltaS) / a0;
44197276fddSZach Atkins           } else {
44297276fddSZach Atkins             /* Use deltaS = 1, pick root that is closest to the last point to prevent doubling back */
44397276fddSZach Atkins             PetscReal dlambda1, dlambda2;
44497276fddSZach Atkins 
44597276fddSZach Atkins             PetscQuadraticRoots(a0, b0 + b1, c0 + c1 + c2, &dlambda1, &dlambda2);
44697276fddSZach Atkins             deltaLambda = (b0 * dlambda1) > (b0 * dlambda2) ? dlambda1 : dlambda2;
44797276fddSZach Atkins           }
44897276fddSZach Atkins         }
44997276fddSZach Atkins       }
45097276fddSZach Atkins       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
45197276fddSZach Atkins       data->lambda = data->lambda + deltaLambda;
45297276fddSZach Atkins       if (data->lambda > data->lambda_max) {
45397276fddSZach 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. */
45497276fddSZach Atkins         deltaLambda  = deltaLambda - (data->lambda - data->lambda_max);
45597276fddSZach Atkins         data->lambda = data->lambda_max;
45697276fddSZach Atkins       }
45797276fddSZach Atkins       if (data->lambda < data->lambda_min) {
45897276fddSZach Atkins         // LCOV_EXCL_START
45997276fddSZach Atkins         /* Ensure that lambda >= lambda_min. This prevents some potential oscillatory behavior. */
46097276fddSZach Atkins         deltaLambda  = deltaLambda - (data->lambda - data->lambda_min);
46197276fddSZach Atkins         data->lambda = data->lambda_min;
46297276fddSZach Atkins         // LCOV_EXCL_STOP
46397276fddSZach Atkins       }
46497276fddSZach Atkins       data->lambda_update = data->lambda_update + deltaLambda;
46597276fddSZach Atkins       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
46697276fddSZach Atkins       PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", lambda=%18.16e, lambda_update=%18.16e\n", snes->iter, (double)data->lambda, (double)data->lambda_update));
46797276fddSZach Atkins       if (j == 0) {
46897276fddSZach Atkins         /* deltaX = deltaLambda*deltaX_Q */
46997276fddSZach Atkins         PetscCall(VecCopy(deltaX_Q, deltaX));
47097276fddSZach Atkins         PetscCall(VecScale(deltaX, deltaLambda));
47197276fddSZach Atkins       } else {
47297276fddSZach Atkins         /* deltaX = deltaS*deltaX_R + deltaLambda*deltaX_Q */
47397276fddSZach Atkins         PetscCall(VecAXPBYPCZ(deltaX, deltaS, deltaLambda, 0, deltaX_R, deltaX_Q));
47497276fddSZach Atkins       }
47597276fddSZach Atkins       PetscCall(VecAXPY(DeltaX, 1, deltaX));
47697276fddSZach Atkins       PetscCall(VecAXPY(X, 1, deltaX));
47797276fddSZach Atkins       /* Q = -dF/dlambda(X, lambda)*/
47897276fddSZach Atkins       PetscCall(SNESNewtonALComputeFunction(snes, X, Q));
47997276fddSZach Atkins       /* R = F(X, lambda) */
48097276fddSZach Atkins       PetscCall(SNESComputeFunction(snes, X, R));
48197276fddSZach Atkins       PetscCall(VecNormBegin(R, NORM_2, &fnorm));
48297276fddSZach Atkins       PetscCall(VecNormBegin(X, NORM_2, &xnorm));
48397276fddSZach Atkins       PetscCall(VecNormBegin(deltaX, NORM_2, &ynorm));
48497276fddSZach Atkins       PetscCall(VecNormEnd(R, NORM_2, &fnorm));
48597276fddSZach Atkins       PetscCall(VecNormEnd(X, NORM_2, &xnorm));
48697276fddSZach Atkins       PetscCall(VecNormEnd(deltaX, NORM_2, &ynorm));
48776c63389SBarry Smith       SNESCheckFunctionDomainError(snes, fnorm);
48897276fddSZach Atkins       if (PetscLogPrintInfo) PetscCall(SNESNewtonALCheckArcLength(snes, DeltaX, data->lambda_update, stepSize));
48997276fddSZach Atkins 
49097276fddSZach Atkins       /* Monitor convergence */
49197276fddSZach Atkins       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
49297276fddSZach Atkins       snes->iter++;
49397276fddSZach Atkins       snes->norm  = fnorm;
49497276fddSZach Atkins       snes->ynorm = ynorm;
49597276fddSZach Atkins       snes->xnorm = xnorm;
49697276fddSZach Atkins       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
49797276fddSZach Atkins       PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
49897276fddSZach Atkins       PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
49997276fddSZach Atkins       PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
5009350dab8SStefano Zampini       if (!snes->reason && j == maxits - 1) snes->reason = SNES_DIVERGED_MAX_IT;
50197276fddSZach Atkins       if (snes->reason) break;
50297276fddSZach Atkins     }
50397276fddSZach Atkins     if (snes->reason < 0) break;
50497276fddSZach Atkins     if (data->lambda >= data->lambda_max) {
50597276fddSZach Atkins       break;
50697276fddSZach Atkins     } else if (maxincs > 0 && i == maxincs - 1) {
50797276fddSZach Atkins       snes->reason = SNES_DIVERGED_MAX_IT;
50897276fddSZach Atkins       break;
50997276fddSZach Atkins     } else {
51097276fddSZach Atkins       snes->reason = SNES_CONVERGED_ITERATING;
51197276fddSZach Atkins     }
51297276fddSZach Atkins   }
51397276fddSZach Atkins   /* Reset RHS vector, if changed */
51497276fddSZach Atkins   if (data->copied_rhs) {
51597276fddSZach Atkins     PetscCall(VecCopy(data->vec_rhs_orig, snes->vec_rhs));
51697276fddSZach Atkins     data->copied_rhs = PETSC_FALSE;
51797276fddSZach Atkins   }
51897276fddSZach Atkins   snes->max_its = maxits; /* reset snes->max_its */
51997276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
52097276fddSZach Atkins }
52197276fddSZach Atkins 
52297276fddSZach Atkins /*
52397276fddSZach Atkins    SNESSetUp_NEWTONAL - Sets up the internal data structures for the later use
52497276fddSZach Atkins    of the SNESNEWTONAL nonlinear solver.
52597276fddSZach Atkins 
52697276fddSZach Atkins    Input Parameter:
52797276fddSZach Atkins .  snes - the SNES context
52897276fddSZach Atkins .  x - the solution vector
52997276fddSZach Atkins 
53097276fddSZach Atkins    Application Interface Routine: SNESSetUp()
53197276fddSZach Atkins  */
SNESSetUp_NEWTONAL(SNES snes)53297276fddSZach Atkins static PetscErrorCode SNESSetUp_NEWTONAL(SNES snes)
53397276fddSZach Atkins {
53497276fddSZach Atkins   PetscFunctionBegin;
53597276fddSZach Atkins   PetscCall(SNESSetWorkVecs(snes, 4));
53697276fddSZach Atkins   PetscCall(SNESSetUpMatrices(snes));
53797276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
53897276fddSZach Atkins }
53997276fddSZach Atkins 
54097276fddSZach Atkins /*
54197276fddSZach Atkins    SNESSetFromOptions_NEWTONAL - Sets various parameters for the SNESNEWTONAL method.
54297276fddSZach Atkins 
54397276fddSZach Atkins    Input Parameter:
54497276fddSZach Atkins .  snes - the SNES context
54597276fddSZach Atkins 
54697276fddSZach Atkins    Application Interface Routine: SNESSetFromOptions()
54797276fddSZach Atkins */
SNESSetFromOptions_NEWTONAL(SNES snes,PetscOptionItems PetscOptionsObject)548ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONAL(SNES snes, PetscOptionItems PetscOptionsObject)
54997276fddSZach Atkins {
55097276fddSZach Atkins   SNES_NEWTONAL             *data            = (SNES_NEWTONAL *)snes->data;
55197276fddSZach Atkins   SNESNewtonALCorrectionType correction_type = data->correction_type;
55297276fddSZach Atkins 
55397276fddSZach Atkins   PetscFunctionBegin;
55497276fddSZach Atkins   PetscOptionsHeadBegin(PetscOptionsObject, "SNES Newton Arc Length options");
55597276fddSZach Atkins   PetscCall(PetscOptionsReal("-snes_newtonal_step_size", "Initial arc length increment step size", "SNESNewtonAL", data->step_size, &data->step_size, NULL));
55697276fddSZach Atkins   PetscCall(PetscOptionsInt("-snes_newtonal_max_continuation_steps", "Maximum number of increment steps", "SNESNewtonAL", data->max_continuation_steps, &data->max_continuation_steps, NULL));
55797276fddSZach Atkins   PetscCall(PetscOptionsReal("-snes_newtonal_psisq", "Regularization parameter for arc length continuation, 0 for cylindrical", "SNESNewtonAL", data->psisq, &data->psisq, NULL));
55897276fddSZach Atkins   PetscCall(PetscOptionsReal("-snes_newtonal_lambda_min", "Minimum value of the load parameter lambda", "SNESNewtonAL", data->lambda_min, &data->lambda_min, NULL));
55997276fddSZach Atkins   PetscCall(PetscOptionsReal("-snes_newtonal_lambda_max", "Maximum value of the load parameter lambda", "SNESNewtonAL", data->lambda_max, &data->lambda_max, NULL));
56097276fddSZach 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));
56197276fddSZach 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));
56297276fddSZach Atkins   PetscCall(SNESNewtonALSetCorrectionType(snes, correction_type));
56397276fddSZach Atkins   PetscOptionsHeadEnd();
56497276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
56597276fddSZach Atkins }
56697276fddSZach Atkins 
SNESReset_NEWTONAL(SNES snes)56797276fddSZach Atkins static PetscErrorCode SNESReset_NEWTONAL(SNES snes)
56897276fddSZach Atkins {
56997276fddSZach Atkins   SNES_NEWTONAL *al = (SNES_NEWTONAL *)snes->data;
57097276fddSZach Atkins 
57197276fddSZach Atkins   PetscFunctionBegin;
57297276fddSZach Atkins   PetscCall(VecDestroy(&al->vec_rhs_orig));
57397276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
57497276fddSZach Atkins }
57597276fddSZach Atkins 
57697276fddSZach Atkins /*
57797276fddSZach Atkins    SNESDestroy_NEWTONAL - Destroys the private SNES_NEWTONAL context that was created
57897276fddSZach Atkins    with SNESCreate_NEWTONAL().
57997276fddSZach Atkins 
58097276fddSZach Atkins    Input Parameter:
58197276fddSZach Atkins .  snes - the SNES context
58297276fddSZach Atkins 
58397276fddSZach Atkins    Application Interface Routine: SNESDestroy()
58497276fddSZach Atkins  */
SNESDestroy_NEWTONAL(SNES snes)58597276fddSZach Atkins static PetscErrorCode SNESDestroy_NEWTONAL(SNES snes)
58697276fddSZach Atkins {
58797276fddSZach Atkins   PetscFunctionBegin;
58897276fddSZach Atkins   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALSetCorrectionType_C", NULL));
58997276fddSZach Atkins   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALGetLoadParameter_C", NULL));
59097276fddSZach Atkins   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALSetFunction_C", NULL));
59197276fddSZach Atkins   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALGetFunction_C", NULL));
59297276fddSZach Atkins   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALComputeFunction_C", NULL));
59397276fddSZach Atkins   PetscCall(PetscFree(snes->data));
59497276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
59597276fddSZach Atkins }
59697276fddSZach Atkins 
59797276fddSZach Atkins /*MC
59897276fddSZach Atkins   SNESNEWTONAL - Newton based nonlinear solver that uses a arc-length continuation method to solve the nonlinear system.
59997276fddSZach Atkins 
60097276fddSZach Atkins   Options Database Keys:
60197276fddSZach Atkins +   -snes_newtonal_step_size <1.0>              - Initial arc length increment step size
60297276fddSZach Atkins .   -snes_newtonal_max_continuation_steps <100> - Maximum number of continuation steps, or negative for no limit (not recommended)
60397276fddSZach Atkins .   -snes_newtonal_psisq <1.0>                  - Regularization parameter for arc length continuation, 0 for cylindrical. Larger values generally lead to more steps.
60497276fddSZach Atkins .   -snes_newtonal_lambda_min <0.0>             - Minimum value of the load parameter lambda
60597276fddSZach Atkins .   -snes_newtonal_lambda_max <1.0>             - Maximum value of the load parameter lambda
60697276fddSZach Atkins .   -snes_newtonal_scale_rhs <true>             - Scale the constant vector passed to `SNESSolve` by the load parameter lambda
60797276fddSZach Atkins -   -snes_newtonal_correction_type <exact>      - Type of correction to use in the arc-length continuation method, `exact` or `normal`
60897276fddSZach Atkins 
60997276fddSZach Atkins   Level: intermediate
61097276fddSZach Atkins 
61197276fddSZach Atkins   Note:
61297276fddSZach Atkins   The exact correction scheme with partial updates is detailed in {cite}`Ritto-CorreaCamotim2008` and the implementation of the
61397276fddSZach Atkins   normal correction scheme is based on {cite}`LeonPaulinoPereiraMenezesLages_2011`.
61497276fddSZach Atkins 
61597276fddSZach Atkins .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONAL`, `SNESNewtonALSetFunction()`, `SNESNewtonALGetFunction()`, `SNESNewtonALGetLoadParameter()`
61697276fddSZach Atkins M*/
SNESCreate_NEWTONAL(SNES snes)61797276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONAL(SNES snes)
61897276fddSZach Atkins {
61997276fddSZach Atkins   SNES_NEWTONAL *arclengthParameters;
62097276fddSZach Atkins 
62197276fddSZach Atkins   PetscFunctionBegin;
62297276fddSZach Atkins   snes->ops->setup          = SNESSetUp_NEWTONAL;
62397276fddSZach Atkins   snes->ops->solve          = SNESSolve_NEWTONAL;
62497276fddSZach Atkins   snes->ops->destroy        = SNESDestroy_NEWTONAL;
62597276fddSZach Atkins   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONAL;
62697276fddSZach Atkins   snes->ops->reset          = SNESReset_NEWTONAL;
62797276fddSZach Atkins 
62897276fddSZach Atkins   snes->usesksp = PETSC_TRUE;
62997276fddSZach Atkins   snes->usesnpc = PETSC_FALSE;
63097276fddSZach Atkins 
63197276fddSZach Atkins   PetscCall(SNESParametersInitialize(snes));
63297276fddSZach Atkins   PetscObjectParameterSetDefault(snes, max_funcs, 30000);
63397276fddSZach Atkins   PetscObjectParameterSetDefault(snes, max_its, 50);
63497276fddSZach Atkins 
63597276fddSZach Atkins   snes->alwayscomputesfinalresidual = PETSC_TRUE;
63697276fddSZach Atkins 
63797276fddSZach Atkins   PetscCall(PetscNew(&arclengthParameters));
63897276fddSZach Atkins   arclengthParameters->lambda                 = 0.0;
63997276fddSZach Atkins   arclengthParameters->lambda_update          = 0.0;
64097276fddSZach Atkins   arclengthParameters->step_size              = 1.0;
64197276fddSZach Atkins   arclengthParameters->max_continuation_steps = 100;
64297276fddSZach Atkins   arclengthParameters->psisq                  = 1.0;
64397276fddSZach Atkins   arclengthParameters->lambda_min             = 0.0;
64497276fddSZach Atkins   arclengthParameters->lambda_max             = 1.0;
64597276fddSZach Atkins   arclengthParameters->scale_rhs              = PETSC_TRUE;
64697276fddSZach Atkins   arclengthParameters->correction_type        = SNES_NEWTONAL_CORRECTION_EXACT;
64797276fddSZach Atkins   snes->data                                  = (void *)arclengthParameters;
64897276fddSZach Atkins 
64997276fddSZach Atkins   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALSetCorrectionType_C", SNESNewtonALSetCorrectionType_NEWTONAL));
65097276fddSZach Atkins   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALGetLoadParameter_C", SNESNewtonALGetLoadParameter_NEWTONAL));
65197276fddSZach Atkins   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALSetFunction_C", SNESNewtonALSetFunction_NEWTONAL));
65297276fddSZach Atkins   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALGetFunction_C", SNESNewtonALGetFunction_NEWTONAL));
65397276fddSZach Atkins   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonALComputeFunction_C", SNESNewtonALComputeFunction_NEWTONAL));
65497276fddSZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
65597276fddSZach Atkins }
656