xref: /petsc/src/tao/leastsquares/impls/pounders/pounders.h (revision 9dd11ecf0918283bb567d8b33a92f53ac4ea7840)
1 #pragma once
2 #include <petsc/private/taoimpl.h>
3 #include <petscblaslapack.h>
4 
5 typedef struct {
6   PetscInt      npmax; /* Max number of interpolation points (>n+1) (def: 2n+1) */
7   PetscInt      nmax;  /* Max(n*(n+1)/2, 5*npmax) */
8   PetscInt      m, n;
9   Vec          *Xhist;
10   Vec          *Fhist;
11   PetscReal    *Fres;        /* (nfmax) */
12   PetscReal    *RES;         /* npxm */
13   PetscReal    *work;        /* (n) */
14   PetscReal    *work2;       /* (n) */
15   PetscReal    *work3;       /* (n) */
16   PetscReal    *xmin;        /* (n) */
17   PetscReal    *mwork;       /* (m) */
18   PetscReal    *Disp;        /* nxn */
19   PetscReal    *Fdiff;       /* nxm */
20   PetscReal    *H;           /* model hessians (mxnxn) */
21   PetscReal    *Hres;        /* nxn */
22   PetscReal    *Gres;        /* n */
23   PetscReal    *Gdel;        /* mxn */
24   PetscReal    *Hdel;        /* mxnxn */
25   PetscReal    *Gpoints;     /* nxn */
26   PetscReal    *C;           /* m */
27   PetscReal    *Xsubproblem; /* n */
28   PetscInt     *indices;     /* 1,2,3...m */
29   PetscInt      minindex;
30   PetscInt      nmodelpoints;
31   PetscInt     *model_indices; /* n */
32   PetscInt      last_nmodelpoints;
33   PetscInt     *last_model_indices; /* n */
34   PetscInt     *interp_indices;     /* n */
35   PetscBLASInt *iwork;              /* n */
36   PetscReal    *w;                  /* nxn */
37   PetscInt      nHist;
38   VecScatter    scatterf, scatterx;
39   Vec           localf, localx, localfmin, localxmin;
40   Vec           workxvec, workfvec;
41   PetscMPIInt   size;
42 
43   PetscReal delta; /* Trust region radius (>0) */
44   PetscReal delta0;
45   PetscBool usegqt;
46   Mat       Hs;
47   Vec       b;
48 
49   PetscReal deltamax;
50   PetscReal deltamin;
51   PetscReal c1;         /* Factor for checking validity */
52   PetscReal c2;         /* Factor for linear poisedness */
53   PetscReal theta1;     /* Pivot threshold for validity */
54   PetscReal theta2;     /* Pivot threshold for additional points */
55   PetscReal gamma0;     /* parameter for shrinking trust region (<1) */
56   PetscReal gamma1;     /* parameter for enlarging trust region (>2) */
57   PetscReal eta0;       /* parameter 1 for accepting point (0 <= eta0 < eta1)*/
58   PetscReal eta1;       /* parameter 2 for accepting point (eta0 < eta1 < 1)*/
59   PetscReal gqt_rtol;   /* parameter used by gqt */
60   PetscInt  gqt_maxits; /* parameter used by gqt */
61 
62   /* QR factorization data */
63   PetscInt      q_is_I;
64   PetscReal    *Q;          /* npmax x npmax */
65   PetscReal    *Q_tmp;      /* npmax x npmax */
66   PetscReal    *tau;        /* scalar factors of H(i) */
67   PetscReal    *tau_tmp;    /* scalar factors of H(i) */
68   PetscReal    *npmaxwork;  /* work vector of length npmax */
69   PetscBLASInt *npmaxiwork; /* integer work vector of length npmax */
70   /* morepoints and getquadnlsmfq */
71   PetscReal *L;      /* n*(n+1)/2 x npmax */
72   PetscReal *L_tmp;  /* n*(n+1)/2 x npmax */
73   PetscReal *L_save; /* n*(n+1)/2 x npmax */
74   PetscReal *Z;      /* npmax x npmax-(n+1) */
75   PetscReal *M;      /* npmax x n+1 */
76   PetscReal *N;      /* npmax x n*(n+1)/2  */
77   PetscReal *alpha;  /* n+1 */
78   PetscReal *beta;   /*  r(n+1)/2 */
79   PetscReal *omega;  /* npmax - np - 1 */
80 
81   Tao subtao;
82   Vec subxl, subxu, subx, subpdel, subndel, subb;
83   Mat subH;
84 
85 } TAO_POUNDERS;
86 
87 PetscErrorCode gqt(PetscInt, PetscReal *, PetscInt, PetscReal *, PetscReal, PetscReal, PetscReal, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscInt *, PetscInt *, PetscReal *, PetscReal *, PetscReal *);
88