Actual source code: pounders.h

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1: #ifndef __TAO_MFQNLS_H
  3:  #include <petsc/private/taoimpl.h>
  4:  #include <petscblaslapack.h>

  6: typedef struct {
  7:   PetscInt npmax;  /* Max number of interpolation points (>n+1) (def: 2n+1) */
  8:   PetscInt nmax; /* Max(n*(n+1)/2, 5*npmax) */
  9:   PetscInt m,n;
 10:   Vec *Xhist;
 11:   Vec *Fhist;
 12:   PetscReal *Fres; /* (nfmax) */
 13:   PetscReal *RES; /* npxm */
 14:   PetscReal *work; /* (n) */
 15:   PetscReal *work2; /* (n) */
 16:   PetscReal *work3; /* (n) */
 17:   PetscReal *xmin; /* (n) */
 18:   PetscReal *mwork; /* (m) */
 19:   PetscReal *Disp; /* nxn */
 20:   PetscReal *Fdiff;/* nxm */
 21:   PetscReal *H; /* model hessians (mxnxn) */
 22:   PetscReal *Hres;  /* nxn */
 23:   PetscReal *Gres;  /* n */
 24:   PetscReal *Gdel; /* mxn */
 25:   PetscReal *Hdel; /* mxnxn */
 26:   PetscReal *Gpoints; /* nxn */
 27:   PetscReal *C; /* m */
 28:   PetscReal *Xsubproblem; /* n */
 29:   PetscInt *indices; /* 1,2,3...m */
 30:   PetscInt minindex;
 31:   PetscInt nmodelpoints;
 32:   PetscInt *model_indices; /* n */
 33:   PetscInt last_nmodelpoints;
 34:   PetscInt *last_model_indices; /* n */
 35:   PetscInt *interp_indices; /* n */
 36:   PetscBLASInt *iwork; /* n */
 37:   PetscReal *w; /* nxn */
 38:   PetscInt nHist;
 39:   VecScatter scatterf,scatterx;
 40:   Vec localf, localx, localfmin, localxmin;
 41:   Vec workxvec,workfvec;
 42:   PetscMPIInt size;


 45:   PetscReal delta; /* Trust region radius (>0) */
 46:   PetscReal delta0;
 47:   PetscBool usegqt;
 48:   Mat Hs;
 49:   Vec b;

 51:   PetscReal deltamax;
 52:   PetscReal deltamin;
 53:   PetscReal c1; /* Factor for checking validity */
 54:   PetscReal c2; /* Factor for linear poisedness */
 55:   PetscReal theta1; /* Pivot threshold for validity */
 56:   PetscReal theta2; /* Pivot threshold for additional points */
 57:   PetscReal gamma0; /* parameter for shrinking trust region (<1) */
 58:   PetscReal gamma1; /* parameter for enlarging trust region (>2) */
 59:   PetscReal eta0;   /* parameter 1 for accepting point (0 <= eta0 < eta1)*/
 60:   PetscReal eta1;   /* parameter 2 for accepting point (eta0 < eta1 < 1)*/
 61:   PetscReal gqt_rtol;   /* parameter used by gqt */
 62:   PetscInt gqt_maxits; /* parameter used by gqt */
 63:   /* QR factorization data */
 64:   PetscInt q_is_I;
 65:   PetscReal *Q; /* npmax x npmax */
 66:   PetscReal *Q_tmp; /* npmax x npmax */
 67:   PetscReal *tau; /* scalar factors of H(i) */
 68:   PetscReal *tau_tmp; /* scalar factors of H(i) */
 69:   PetscReal *npmaxwork; /* work vector of length npmax */
 70:   PetscBLASInt *npmaxiwork; /* integer work vector of length npmax */
 71:   /* morepoints and getquadnlsmfq */
 72:   PetscReal *L;   /* n*(n+1)/2 x npmax */
 73:   PetscReal *L_tmp;   /* n*(n+1)/2 x npmax */
 74:   PetscReal *L_save;   /* n*(n+1)/2 x npmax */
 75:   PetscReal *Z;   /* npmax x npmax-(n+1) */
 76:   PetscReal *M;   /* npmax x n+1 */
 77:   PetscReal *N;   /* npmax x n*(n+1)/2  */
 78:   PetscReal *alpha; /* n+1 */
 79:   PetscReal *beta; /*  r(n+1)/2 */
 80:   PetscReal *omega; /* npmax - np - 1 */

 82:   Tao subtao;
 83:   Vec       subxl,subxu,subx,subpdel,subndel,subb;
 84:   Mat       subH;

 86: } TAO_POUNDERS;


 89: 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: PetscErrorCode gqtwrap(Tao tao,PetscReal *gnorm, PetscReal *qmin);
 92: PetscErrorCode phi2eval(PetscReal *x, PetscInt n, PetscReal *phi);
 93: PetscErrorCode getquadpounders(TAO_POUNDERS *mfqP);
 94: PetscErrorCode morepoints(TAO_POUNDERS *mfqP);
 95: PetscErrorCode addpoint(Tao tao, TAO_POUNDERS *mfqP, PetscInt index);
 96: PetscErrorCode modelimprove(Tao tao, TAO_POUNDERS *mfqP, PetscInt addallpoints);
 97: PetscErrorCode affpoints(TAO_POUNDERS *mfqP, PetscReal *xmin, PetscReal c);

 99: EXTERN_C_BEGIN
100: 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);
101: EXTERN_C_END
102: #endif /* ifndef __TAO_MFQNLS */