Actual source code: alimpl.h

  1: /*
  2:    Private context for a Newton arc length continuation
  3:    method for solving systems of nonlinear equations
  4: */

  6: #pragma once
  7: #include <petsc/private/snesimpl.h>
  8: #include <petscsnes.h>

 10: typedef struct {
 11:   PetscInt                   max_continuation_steps; /* maximum number of continuation steps */
 12:   PetscReal                  step_size;              /* radius of quadratic constraint surface */
 13:   PetscReal                  psisq;                  /* load factor regularization parameter */
 14:   PetscReal                  lambda_update;          /* accumulated update to lambda over the current increment */
 15:   PetscReal                  lambda;                 /* load parameter */
 16:   PetscReal                  lambda_max;             /* maximum value of the load parameter */
 17:   PetscReal                  lambda_min;             /* minimum value of the load parameter */
 18:   PetscBool                  scale_rhs;              /* should the RHS vector be scaled by the load parameter? */
 19:   SNESNewtonALCorrectionType correction_type;        /* type of correction scheme to use */
 20:   PetscBool                  copied_rhs;             /* has the right-hand side vector been copied? */

 22:   Vec vec_rhs_orig; /* original right-hand side vector, used if `scale_rhs == PETSC_TRUE` */

 24:   SNESFunctionFn *computealfunction; /* user-provided function to compute the tangent load vector  */
 25:   void           *alctx;             /* user-provided context for the tangent load vector computation */
 26: } SNES_NEWTONAL;

 28: PETSC_INTERN const char NewtonALExactCitation[];
 29: PETSC_INTERN PetscBool  NewtonALExactCitationSet;
 30: PETSC_INTERN const char NewtonALNormalCitation[];
 31: PETSC_INTERN PetscBool  NewtonALNormalCitationSet;