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