xref: /petsc/src/ksp/ksp/impls/cg/gltr/gltrimpl.h (revision ed0a315dc690bca58568ab127a9f322718ff62b7)
1 /*
2    Context for using preconditioned conjugate gradient method to minimized a
3    quadratic function subject to a trust region constraint.  If the matrix
4    is indefinite, a direction of negative curvature may be encountered.  If
5    a direction of negative curvature is found, we continue to build the
6    tridiagonal Lanczos matrix for a fixed number of iterations.  After this
7    matrix is computed, we compute a global solution to solve the trust-
8    region problem with the tridiagonal approximation by using a variant of
9    the More'-Sorenson algorithm.  The direction is then constructed from
10    this solution.
11 
12    This method is described in:
13      N. Gould, S. Lucidi, M. Roma, and Ph. Toint, "Solving the Trust-Region
14        Subproblem using the Lanczos Method", SIAM Journal on Optimization,
15        9, pages 504-525, 1999.
16 */
17 
18 #pragma once
19 
20 #include <petsc/private/kspimpl.h>
21 
22 typedef struct {
23   PetscReal *diag;   /* Diagonal part of Lanczos matrix           */
24   PetscReal *offd;   /* Off-diagonal part of Lanczos matrix       */
25   PetscReal *alpha;  /* Record of alpha values from CG            */
26   PetscReal *beta;   /* Record of beta values from CG             */
27   PetscReal *norm_r; /* Record of residual values from CG         */
28 
29   PetscReal    *rwork; /* Real workspace for solver computations    */
30   PetscBLASInt *iwork; /* Integer workspace for solver computations */
31 
32   PetscReal radius;
33   PetscReal norm_d;
34   PetscReal e_min;
35   PetscReal o_fcn;
36   PetscReal lambda;
37 
38   PetscReal init_pert;  /* Initial perturbation for solve            */
39   PetscReal eigen_tol;  /* Tolerance used when computing eigenvalue  */
40   PetscReal newton_tol; /* Tolerance used for newton method          */
41 
42   PetscInt alloced;    /* Size of workspace vectors allocated       */
43   PetscInt init_alloc; /* Initial size for workspace vectors        */
44 
45   PetscInt max_lanczos_its; /* Maximum lanczos iterations                */
46   PetscInt max_newton_its;  /* Maximum newton iterations                 */
47   PetscInt dtype;           /* Method used to measure the norm of step   */
48 } KSPCG_GLTR;
49