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