Actual source code: trimpl.h
1: /*
2: Context for a Newton trust region method for solving a system
3: of nonlinear equations
4: */
6: #pragma once
7: #include <petsc/private/snesimpl.h>
9: typedef struct {
10: PetscReal delta; /* trust region radius */
12: PetscObjectParameterDeclare(PetscReal, delta0); /* initial radius for trust region */
13: PetscObjectParameterDeclare(PetscReal, deltaM); /* maximum radius for trust region */
14: PetscObjectParameterDeclare(PetscReal, deltam); /* minimum radius for trust region */
16: PetscReal kmdc; /* sufficient decrease parameter */
18: /*
19: Given rho = (fk - fkp1) / (m(0) - m(pk))
21: The radius is modified as:
22: rho < eta2 -> delta *= t1
23: rho > eta3 -> delta *= t2
24: delta = min(delta,deltaM)
26: The step is accepted if rho > eta1
27: The iterative process is halted if delta < deltam
28: */
29: PetscObjectParameterDeclare(PetscReal, eta1);
30: PetscObjectParameterDeclare(PetscReal, eta2);
31: PetscObjectParameterDeclare(PetscReal, eta3);
32: PetscObjectParameterDeclare(PetscReal, t1);
33: PetscObjectParameterDeclare(PetscReal, t2);
35: /* Use quasi-Newton models for J and (possibly different) Jp */
36: SNESNewtonTRQNType qn;
37: Mat qnB;
38: Mat qnB_pre;
40: /* The type of norm for the trust region */
41: NormType norm;
43: SNESNewtonTRFallbackType fallback; /* enum to distinguish fallback in case Newton step is outside of the trust region */
45: PetscErrorCode (*precheck)(SNES, Vec, Vec, PetscBool *, void *);
46: void *precheckctx;
47: PetscErrorCode (*postcheck)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *);
48: void *postcheckctx;
49: } SNES_NEWTONTR;