xref: /petsc/src/tao/leastsquares/impls/pounders/pounders.c (revision 5b0d99aef51158554e42172a250c96a448fd1e2b)
1 #include "pounders.h"
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "pounders_h"
5 static PetscErrorCode pounders_h(TaoSolver subtao, Vec v, Mat *H, Mat *Hpre, MatStructure *flag, void *ctx)
6 {
7   PetscFunctionBegin;
8   *flag = SAME_NONZERO_PATTERN;
9   PetscFunctionReturn(0);
10 }
11 #undef __FUNCT__
12 #define __FUNCT__ "pounders_fg"
13 static PetscErrorCode  pounders_fg(TaoSolver subtao, Vec x, PetscReal *f, Vec g, void *ctx)
14 {
15   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)ctx;
16   PetscReal      d1,d2;
17   PetscErrorCode ierr;
18 
19   PetscFunctionBegin;
20   /* g = A*x  (add b later)*/
21   ierr = MatMult(mfqP->subH,x,g);CHKERRQ(ierr);
22 
23   /* f = 1/2 * x'*(Ax) + b'*x  */
24   ierr = VecDot(x,g,&d1);CHKERRQ(ierr);
25   ierr = VecDot(mfqP->subb,x,&d2);CHKERRQ(ierr);
26   *f = 0.5 *d1 + d2;
27 
28   /* now  g = g + b */
29   ierr = VecAXPY(g, 1.0, mfqP->subb);CHKERRQ(ierr);
30   PetscFunctionReturn(0);
31 }
32 
33 #undef __FUNCT__
34 #define __FUNCT__ "gqtwrap"
35 PetscErrorCode gqtwrap(TaoSolver tao,PetscReal *gnorm, PetscReal *qmin)
36 {
37   PetscErrorCode ierr;
38   PetscReal      atol=1.0e-10;
39   PetscInt       info,its;
40   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;
41 
42   PetscFunctionBegin;
43   if (! mfqP->usegqt) {
44     PetscReal maxval;
45     PetscInt  i,j;
46 
47     ierr = VecSetValues(mfqP->subb,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);CHKERRQ(ierr);
48     ierr = VecAssemblyBegin(mfqP->subb);CHKERRQ(ierr);
49     ierr = VecAssemblyEnd(mfqP->subb);CHKERRQ(ierr);
50 
51     ierr = VecSet(mfqP->subx,0.0);CHKERRQ(ierr);
52 
53     ierr = VecSet(mfqP->subndel,-mfqP->delta);CHKERRQ(ierr);
54     ierr = VecSet(mfqP->subpdel,mfqP->delta);CHKERRQ(ierr);
55 
56     for (i=0;i<mfqP->n;i++) {
57       for (j=i;j<mfqP->n;j++) {
58         mfqP->Hres[j+mfqP->n*i] = mfqP->Hres[mfqP->n*j+i];
59       }
60     }
61     ierr = MatSetValues(mfqP->subH,mfqP->n,mfqP->indices,mfqP->n,mfqP->indices,mfqP->Hres,INSERT_VALUES);
62     ierr = MatAssemblyBegin(mfqP->subH,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
63     ierr = MatAssemblyEnd(mfqP->subH,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
64 
65     ierr = TaoResetStatistics(mfqP->subtao);CHKERRQ(ierr);
66     ierr = TaoSetTolerances(mfqP->subtao,PETSC_DEFAULT,PETSC_DEFAULT,*gnorm,*gnorm,PETSC_DEFAULT);CHKERRQ(ierr);
67     /* enforce bound constraints -- experimental */
68     if (tao->XU && tao->XL) {
69       ierr = VecCopy(tao->XU,mfqP->subxu);CHKERRQ(ierr);
70       ierr = VecAXPY(mfqP->subxu,-1.0,tao->solution);CHKERRQ(ierr);
71       ierr = VecScale(mfqP->subxu,1.0/mfqP->delta);CHKERRQ(ierr);
72       ierr = VecCopy(tao->XL,mfqP->subxl);CHKERRQ(ierr);
73       ierr = VecAXPY(mfqP->subxl,-1.0,tao->solution);CHKERRQ(ierr);
74       ierr = VecScale(mfqP->subxl,1.0/mfqP->delta);CHKERRQ(ierr);
75 
76       ierr = VecPointwiseMin(mfqP->subxu,mfqP->subxu,mfqP->subpdel);CHKERRQ(ierr);
77       ierr = VecPointwiseMax(mfqP->subxl,mfqP->subxl,mfqP->subndel);CHKERRQ(ierr);
78     } else {
79       ierr = VecCopy(mfqP->subpdel,mfqP->subxu);CHKERRQ(ierr);
80       ierr = VecCopy(mfqP->subndel,mfqP->subxl);CHKERRQ(ierr);
81     }
82     /* Make sure xu > xl */
83     ierr = VecCopy(mfqP->subxl,mfqP->subpdel);CHKERRQ(ierr);
84     ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);CHKERRQ(ierr);
85     ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr);
86     if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"upper bound < lower bound in subproblem");
87     /* Make sure xu > tao->solution > xl */
88     ierr = VecCopy(mfqP->subxl,mfqP->subpdel);CHKERRQ(ierr);
89     ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subx);CHKERRQ(ierr);
90     ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr);
91     if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"initial guess < lower bound in subproblem");
92 
93     ierr = VecCopy(mfqP->subx,mfqP->subpdel);CHKERRQ(ierr);
94     ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);CHKERRQ(ierr);
95     ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr);
96     if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"initial guess > upper bound in subproblem");
97 
98     ierr = TaoSolve(mfqP->subtao);CHKERRQ(ierr);
99     ierr = TaoGetSolutionStatus(mfqP->subtao,NULL,qmin,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
100 
101     /* test bounds post-solution*/
102     ierr = VecCopy(mfqP->subxl,mfqP->subpdel);CHKERRQ(ierr);
103     ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subx);CHKERRQ(ierr);
104     ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr);
105     if (maxval > 1e-5) {
106       ierr = PetscInfo(tao,"subproblem solution < lower bound");CHKERRQ(ierr);
107       tao->reason = TAO_DIVERGED_TR_REDUCTION;
108     }
109 
110     ierr = VecCopy(mfqP->subx,mfqP->subpdel);CHKERRQ(ierr);
111     ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);CHKERRQ(ierr);
112     ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr);
113     if (maxval > 1e-5) {
114       ierr = PetscInfo(tao,"subproblem solution > upper bound");
115       tao->reason = TAO_DIVERGED_TR_REDUCTION;
116     }
117   } else {
118     gqt(mfqP->n,mfqP->Hres,mfqP->n,mfqP->Gres,1.0,mfqP->gqt_rtol,atol,mfqP->gqt_maxits,gnorm,qmin,mfqP->Xsubproblem,&info,&its,mfqP->work,mfqP->work2, mfqP->work3);
119   }
120   *qmin *= -1;
121   PetscFunctionReturn(0);
122 }
123 
124 #undef __FUNCT__
125 #define __FUNCT__ "phi2eval"
126 PetscErrorCode phi2eval(PetscReal *x, PetscInt n, PetscReal *phi)
127 {
128 /* Phi = .5*[x(1)^2  sqrt(2)*x(1)*x(2) ... sqrt(2)*x(1)*x(n) ... x(2)^2 sqrt(2)*x(2)*x(3) .. x(n)^2] */
129   PetscInt  i,j,k;
130   PetscReal sqrt2 = PetscSqrtReal(2.0);
131 
132   PetscFunctionBegin;
133   j=0;
134   for (i=0;i<n;i++) {
135     phi[j] = 0.5 * x[i]*x[i];
136     j++;
137     for (k=i+1;k<n;k++) {
138       phi[j]  = x[i]*x[k]/sqrt2;
139       j++;
140     }
141   }
142   PetscFunctionReturn(0);
143 }
144 
145 #undef __FUNCT__
146 #define __FUNCT__ "getquadpounders"
147 PetscErrorCode getquadpounders(TAO_POUNDERS *mfqP)
148 {
149 /* Computes the parameters of the quadratic Q(x) = c + g'*x + 0.5*x*G*x'
150    that satisfies the interpolation conditions Q(X[:,j]) = f(j)
151    for j=1,...,m and with a Hessian matrix of least Frobenius norm */
152 
153     /* NB --we are ignoring c */
154   PetscInt     i,j,k,num,np = mfqP->nmodelpoints;
155   PetscReal    one = 1.0,zero=0.0,negone=-1.0;
156   PetscBLASInt blasnpmax = mfqP->npmax;
157   PetscBLASInt blasnplus1 = mfqP->n+1;
158   PetscBLASInt blasnp = np;
159   PetscBLASInt blasint = mfqP->n*(mfqP->n+1) / 2;
160   PetscBLASInt blasint2 = np - mfqP->n-1;
161   PetscBLASInt info,ione=1;
162   PetscReal    sqrt2 = PetscSqrtReal(2.0);
163 
164   PetscFunctionBegin;
165   for (i=0;i<mfqP->n*mfqP->m;i++) {
166     mfqP->Gdel[i] = 0;
167   }
168   for (i=0;i<mfqP->n*mfqP->n*mfqP->m;i++) {
169     mfqP->Hdel[i] = 0;
170   }
171 
172     /* factor M */
173   LAPACKgetrf_(&blasnplus1,&blasnpmax,mfqP->M,&blasnplus1,mfqP->npmaxiwork,&info);
174   if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine getrf returned with value %d\n",info);
175 
176   if (np == mfqP->n+1) {
177     for (i=0;i<mfqP->npmax-mfqP->n-1;i++) {
178       mfqP->omega[i]=0.0;
179     }
180     for (i=0;i<mfqP->n*(mfqP->n+1)/2;i++) {
181       mfqP->beta[i]=0.0;
182     }
183   } else {
184     /* Let Ltmp = (L'*L) */
185     BLASgemm_("T","N",&blasint2,&blasint2,&blasint,&one,&mfqP->L[(mfqP->n+1)*blasint],&blasint,&mfqP->L[(mfqP->n+1)*blasint],&blasint,&zero,mfqP->L_tmp,&blasint);
186 
187     /* factor Ltmp */
188     LAPACKpotrf_("L",&blasint2,mfqP->L_tmp,&blasint,&info);
189     if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine potrf returned with value %d\n",info);
190   }
191 
192   for (k=0;k<mfqP->m;k++) {
193     if (np != mfqP->n+1) {
194       /* Solve L'*L*Omega = Z' * RESk*/
195       BLASgemv_("T",&blasnp,&blasint2,&one,mfqP->Z,&blasnpmax,&mfqP->RES[mfqP->npmax*k],&ione,&zero,mfqP->omega,&ione);
196       LAPACKpotrs_("L",&blasint2,&ione,mfqP->L_tmp,&blasint,mfqP->omega,&blasint2,&info);
197       if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine potrs returned with value %d\n",info);
198 
199       /* Beta = L*Omega */
200       BLASgemv_("N",&blasint,&blasint2,&one,&mfqP->L[(mfqP->n+1)*blasint],&blasint,mfqP->omega,&ione,&zero,mfqP->beta,&ione);
201     }
202 
203     /* solve M'*Alpha = RESk - N'*Beta */
204     BLASgemv_("T",&blasint,&blasnp,&negone,mfqP->N,&blasint,mfqP->beta,&ione,&one,&mfqP->RES[mfqP->npmax*k],&ione);
205     LAPACKgetrs_("T",&blasnplus1,&ione,mfqP->M,&blasnplus1,mfqP->npmaxiwork,&mfqP->RES[mfqP->npmax*k],&blasnplus1,&info);
206     if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine getrs returned with value %d\n",info);
207 
208     /* Gdel(:,k) = Alpha(2:n+1) */
209     for (i=0;i<mfqP->n;i++) {
210       mfqP->Gdel[i + mfqP->n*k] = mfqP->RES[mfqP->npmax*k + i+1];
211     }
212 
213     /* Set Hdels */
214     num=0;
215     for (i=0;i<mfqP->n;i++) {
216       /* H[i,i,k] = Beta(num) */
217       mfqP->Hdel[(i*mfqP->n + i)*mfqP->m + k] = mfqP->beta[num];
218       num++;
219       for (j=i+1;j<mfqP->n;j++) {
220         /* H[i,j,k] = H[j,i,k] = Beta(num)/sqrt(2) */
221         mfqP->Hdel[(j*mfqP->n + i)*mfqP->m + k] = mfqP->beta[num]/sqrt2;
222         mfqP->Hdel[(i*mfqP->n + j)*mfqP->m + k] = mfqP->beta[num]/sqrt2;
223         num++;
224       }
225     }
226   }
227   PetscFunctionReturn(0);
228 }
229 
230 #undef __FUNCT__
231 #define __FUNCT__ "morepoints"
232 PetscErrorCode morepoints(TAO_POUNDERS *mfqP)
233 {
234   /* Assumes mfqP->model_indices[0]  is minimum index
235    Finishes adding points to mfqP->model_indices (up to npmax)
236    Computes L,Z,M,N
237    np is actual number of points in model (should equal npmax?) */
238   PetscInt       point,i,j,offset;
239   PetscInt       reject;
240   PetscBLASInt   blasn=mfqP->n,blasnpmax=mfqP->npmax,blasnplus1=mfqP->n+1,info,blasnmax=mfqP->nmax,blasint,blasint2,blasnp,blasmaxmn;
241   PetscReal      *x,normd;
242   PetscErrorCode ierr;
243 
244   PetscFunctionBegin;
245   /* Initialize M,N */
246   for (i=0;i<mfqP->n+1;i++) {
247     ierr = VecGetArray(mfqP->Xhist[mfqP->model_indices[i]],&x);CHKERRQ(ierr);
248     mfqP->M[(mfqP->n+1)*i] = 1.0;
249     for (j=0;j<mfqP->n;j++) {
250       mfqP->M[j+1+((mfqP->n+1)*i)] = (x[j]  - mfqP->xmin[j]) / mfqP->delta;
251     }
252     ierr = VecRestoreArray(mfqP->Xhist[mfqP->model_indices[i]],&x);CHKERRQ(ierr);
253     ierr = phi2eval(&mfqP->M[1+((mfqP->n+1)*i)],mfqP->n,&mfqP->N[mfqP->n*(mfqP->n+1)/2 * i]);CHKERRQ(ierr);
254   }
255 
256   /* Now we add points until we have npmax starting with the most recent ones */
257   point = mfqP->nHist-1;
258   mfqP->nmodelpoints = mfqP->n+1;
259   while (mfqP->nmodelpoints < mfqP->npmax && point>=0) {
260     /* Reject any points already in the model */
261     reject = 0;
262     for (j=0;j<mfqP->n+1;j++) {
263       if (point == mfqP->model_indices[j]) {
264         reject = 1;
265         break;
266       }
267     }
268 
269     /* Reject if norm(d) >c2 */
270     if (!reject) {
271       ierr = VecCopy(mfqP->Xhist[point],mfqP->workxvec);CHKERRQ(ierr);
272       ierr = VecAXPY(mfqP->workxvec,-1.0,mfqP->Xhist[mfqP->minindex]);CHKERRQ(ierr);
273       ierr = VecNorm(mfqP->workxvec,NORM_2,&normd);CHKERRQ(ierr);
274       normd /= mfqP->delta;
275       if (normd > mfqP->c2) {
276         reject =1;
277       }
278     }
279     if (reject){
280       point--;
281       continue;
282     }
283 
284     ierr = VecGetArray(mfqP->Xhist[point],&x);CHKERRQ(ierr);
285     mfqP->M[(mfqP->n+1)*mfqP->nmodelpoints] = 1.0;
286     for (j=0;j<mfqP->n;j++) {
287       mfqP->M[j+1+((mfqP->n+1)*mfqP->nmodelpoints)] = (x[j]  - mfqP->xmin[j]) / mfqP->delta;
288     }
289     ierr = VecRestoreArray(mfqP->Xhist[point],&x);CHKERRQ(ierr);
290     ierr = phi2eval(&mfqP->M[1+(mfqP->n+1)*mfqP->nmodelpoints],mfqP->n,&mfqP->N[mfqP->n*(mfqP->n+1)/2 * (mfqP->nmodelpoints)]);CHKERRQ(ierr);
291 
292     /* Update QR factorization */
293     /* Copy M' to Q_tmp */
294     for (i=0;i<mfqP->n+1;i++) {
295       for (j=0;j<mfqP->npmax;j++) {
296         mfqP->Q_tmp[j+mfqP->npmax*i] = mfqP->M[i+(mfqP->n+1)*j];
297       }
298     }
299     blasnp = mfqP->nmodelpoints+1;
300     /* Q_tmp,R = qr(M') */
301     blasmaxmn=PetscMax(mfqP->m,mfqP->n+1);
302     LAPACKgeqrf_(&blasnp,&blasnplus1,mfqP->Q_tmp,&blasnpmax,mfqP->tau_tmp,mfqP->mwork,&blasmaxmn,&info);
303     if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine geqrf returned with value %d\n",info);
304 
305     /* Reject if min(svd(N*Q(:,n+2:np+1)) <= theta2 */
306     /* L = N*Qtmp */
307     blasint2 = mfqP->n * (mfqP->n+1) / 2;
308     /* Copy N to L_tmp */
309     for (i=0;i<mfqP->n*(mfqP->n+1)/2 * mfqP->npmax;i++) {
310       mfqP->L_tmp[i]= mfqP->N[i];
311     }
312     /* Copy L_save to L_tmp */
313 
314     /* L_tmp = N*Qtmp' */
315     LAPACKormqr_("R","N",&blasint2,&blasnp,&blasnplus1,mfqP->Q_tmp,&blasnpmax,mfqP->tau_tmp,mfqP->L_tmp,&blasint2,mfqP->npmaxwork,&blasnmax,&info);
316     if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine ormqr returned with value %d\n",info);
317 
318     /* Copy L_tmp to L_save */
319     for (i=0;i<mfqP->npmax * mfqP->n*(mfqP->n+1)/2;i++) {
320       mfqP->L_save[i] = mfqP->L_tmp[i];
321     }
322 
323     /* Get svd for L_tmp(:,n+2:np+1) (L_tmp is modified in process) */
324     blasint = mfqP->nmodelpoints - mfqP->n;
325     LAPACKgesvd_("N","N",&blasint2,&blasint,&mfqP->L_tmp[(mfqP->n+1)*blasint2],&blasint2,mfqP->beta,mfqP->work,&blasn,mfqP->work,&blasn,mfqP->npmaxwork,&blasnmax,&info);
326     if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine gesvd returned with value %d\n",info);
327 
328     if (mfqP->beta[PetscMin(blasint,blasint2)-1] > mfqP->theta2) {
329       /* accept point */
330       mfqP->model_indices[mfqP->nmodelpoints] = point;
331       /* Copy Q_tmp to Q */
332       for (i=0;i<mfqP->npmax* mfqP->npmax;i++) {
333         mfqP->Q[i] = mfqP->Q_tmp[i];
334       }
335       for (i=0;i<mfqP->npmax;i++){
336         mfqP->tau[i] = mfqP->tau_tmp[i];
337       }
338       mfqP->nmodelpoints++;
339       blasnp = mfqP->nmodelpoints+1;
340 
341       /* Copy L_save to L */
342       for (i=0;i<mfqP->npmax * mfqP->n*(mfqP->n+1)/2;i++) {
343         mfqP->L[i] = mfqP->L_save[i];
344       }
345     }
346     point--;
347   }
348 
349   blasnp = mfqP->nmodelpoints;
350   /* Copy Q(:,n+2:np) to Z */
351   /* First set Q_tmp to I */
352   for (i=0;i<mfqP->npmax*mfqP->npmax;i++) {
353     mfqP->Q_tmp[i] = 0.0;
354   }
355   for (i=0;i<mfqP->npmax;i++) {
356     mfqP->Q_tmp[i + mfqP->npmax*i] = 1.0;
357   }
358 
359   /* Q_tmp = I * Q */
360   LAPACKormqr_("R","N",&blasnp,&blasnp,&blasnplus1,mfqP->Q,&blasnpmax,mfqP->tau,mfqP->Q_tmp,&blasnpmax,mfqP->npmaxwork,&blasnmax,&info);
361   if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine ormqr returned with value %d\n",info);
362 
363   /* Copy Q_tmp(:,n+2:np) to Z) */
364   offset = mfqP->npmax * (mfqP->n+1);
365   for (i=offset;i<mfqP->npmax*mfqP->npmax;i++) {
366     mfqP->Z[i-offset] = mfqP->Q_tmp[i];
367   }
368 
369   if (mfqP->nmodelpoints == mfqP->n + 1) {
370     /* Set L to I_{n+1} */
371     for (i=0;i<mfqP->npmax * mfqP->n*(mfqP->n+1)/2;i++) {
372       mfqP->L[i] = 0.0;
373     }
374     for (i=0;i<mfqP->n;i++) {
375       mfqP->L[(mfqP->n*(mfqP->n+1)/2)*i + i] = 1.0;
376     }
377   }
378   PetscFunctionReturn(0);
379 }
380 
381 #undef __FUNCT__
382 #define __FUNCT__ "addpoint"
383 /* Only call from modelimprove, addpoint() needs ->Q_tmp and ->work to be set */
384 PetscErrorCode addpoint(TaoSolver tao, TAO_POUNDERS *mfqP, PetscInt index)
385 {
386   PetscErrorCode ierr;
387 
388   PetscFunctionBegin;
389   /* Create new vector in history: X[newidx] = X[mfqP->index] + delta*X[index]*/
390   ierr = VecDuplicate(mfqP->Xhist[0],&mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr);
391   ierr = VecSetValues(mfqP->Xhist[mfqP->nHist],mfqP->n,mfqP->indices,&mfqP->Q_tmp[index*mfqP->npmax],INSERT_VALUES);CHKERRQ(ierr);
392   ierr = VecAssemblyBegin(mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr);
393   ierr = VecAssemblyEnd(mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr);
394   ierr = VecAYPX(mfqP->Xhist[mfqP->nHist],mfqP->delta,mfqP->Xhist[mfqP->minindex]);CHKERRQ(ierr);
395 
396   /* Project into feasible region */
397   if (tao->XU && tao->XL) {
398     ierr = VecMedian(mfqP->Xhist[mfqP->nHist], tao->XL, tao->XU, mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr);
399   }
400 
401   /* Compute value of new vector */
402   ierr = VecDuplicate(mfqP->Fhist[0],&mfqP->Fhist[mfqP->nHist]);CHKERRQ(ierr);
403   CHKMEMQ;
404   ierr = TaoComputeSeparableObjective(tao,mfqP->Xhist[mfqP->nHist],mfqP->Fhist[mfqP->nHist]);CHKERRQ(ierr);
405   ierr = VecNorm(mfqP->Fhist[mfqP->nHist],NORM_2,&mfqP->Fres[mfqP->nHist]);CHKERRQ(ierr);
406   if (PetscIsInfOrNanReal(mfqP->Fres[mfqP->nHist])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
407   mfqP->Fres[mfqP->nHist]*=mfqP->Fres[mfqP->nHist];
408 
409   /* Add new vector to model */
410   mfqP->model_indices[mfqP->nmodelpoints] = mfqP->nHist;
411   mfqP->nmodelpoints++;
412   mfqP->nHist++;
413   PetscFunctionReturn(0);
414 }
415 
416 #undef __FUNCT__
417 #define __FUNCT__ "modelimprove"
418 PetscErrorCode modelimprove(TaoSolver tao, TAO_POUNDERS *mfqP, PetscInt addallpoints)
419 {
420   /* modeld = Q(:,np+1:n)' */
421   PetscErrorCode ierr;
422   PetscInt       i,j,minindex=0;
423   PetscReal      dp,half=0.5,one=1.0,minvalue=TAO_INFINITY;
424   PetscBLASInt   blasn=mfqP->n,  blasnpmax = mfqP->npmax, blask,info;
425   PetscBLASInt   blas1=1,blasnmax = mfqP->nmax;
426 
427   blask = mfqP->nmodelpoints;
428   /* Qtmp = I(n x n) */
429   for (i=0;i<mfqP->n;i++) {
430     for (j=0;j<mfqP->n;j++) {
431       mfqP->Q_tmp[i + mfqP->npmax*j] = 0.0;
432     }
433   }
434   for (j=0;j<mfqP->n;j++) {
435     mfqP->Q_tmp[j + mfqP->npmax*j] = 1.0;
436   }
437 
438   /* Qtmp = Q * I */
439   LAPACKormqr_("R","N",&blasn,&blasn,&blask,mfqP->Q,&blasnpmax,mfqP->tau, mfqP->Q_tmp, &blasnpmax, mfqP->npmaxwork,&blasnmax, &info);
440 
441   for (i=mfqP->nmodelpoints;i<mfqP->n;i++) {
442     dp = BLASdot_(&blasn,&mfqP->Q_tmp[i*mfqP->npmax],&blas1,mfqP->Gres,&blas1);
443     if (dp>0.0) { /* Model says use the other direction! */
444       for (j=0;j<mfqP->n;j++) {
445         mfqP->Q_tmp[i*mfqP->npmax+j] *= -1;
446       }
447     }
448     /* mfqP->work[i] = Cres+Modeld(i,:)*(Gres+.5*Hres*Modeld(i,:)') */
449     for (j=0;j<mfqP->n;j++) {
450       mfqP->work2[j] = mfqP->Gres[j];
451     }
452     BLASgemv_("N",&blasn,&blasn,&half,mfqP->Hres,&blasn,&mfqP->Q_tmp[i*mfqP->npmax], &blas1, &one, mfqP->work2,&blas1);
453     mfqP->work[i] = BLASdot_(&blasn,&mfqP->Q_tmp[i*mfqP->npmax],&blas1,mfqP->work2,&blas1);
454     if (i==mfqP->nmodelpoints || mfqP->work[i] < minvalue) {
455       minindex=i;
456       minvalue = mfqP->work[i];
457     }
458     if (addallpoints != 0) {
459       ierr = addpoint(tao,mfqP,i);CHKERRQ(ierr);
460     }
461   }
462   if (!addallpoints) {
463     ierr = addpoint(tao,mfqP,minindex);CHKERRQ(ierr);
464   }
465   PetscFunctionReturn(0);
466 }
467 
468 
469 #undef __FUNCT__
470 #define __FUNCT__ "affpoints"
471 PetscErrorCode affpoints(TAO_POUNDERS *mfqP, PetscReal *xmin,PetscReal c)
472 {
473   PetscInt       i,j;
474   PetscBLASInt   blasm=mfqP->m,blasj,blask,blasn=mfqP->n,ione=1,info;
475   PetscBLASInt   blasnpmax = mfqP->npmax,blasmaxmn;
476   PetscReal      proj,normd;
477   PetscReal      *x;
478   PetscErrorCode ierr;
479 
480   PetscFunctionBegin;
481   for (i=mfqP->nHist-1;i>=0;i--) {
482     ierr = VecGetArray(mfqP->Xhist[i],&x);CHKERRQ(ierr);
483     for (j=0;j<mfqP->n;j++) {
484       mfqP->work[j] = (x[j] - xmin[j])/mfqP->delta;
485     }
486     ierr = VecRestoreArray(mfqP->Xhist[i],&x);CHKERRQ(ierr);
487     BLAScopy_(&blasn,mfqP->work,&ione,mfqP->work2,&ione);
488     normd = BLASnrm2_(&blasn,mfqP->work,&ione);
489     if (normd <= c*c) {
490       blasj=PetscMax((mfqP->n - mfqP->nmodelpoints),0);
491       if (!mfqP->q_is_I) {
492         /* project D onto null */
493         blask=(mfqP->nmodelpoints);
494         LAPACKormqr_("R","N",&ione,&blasn,&blask,mfqP->Q,&blasnpmax,mfqP->tau,mfqP->work2,&ione,mfqP->mwork,&blasm,&info);
495         if (info < 0) SETERRQ1(PETSC_COMM_SELF,1,"ormqr returned value %d\n",info);
496       }
497       proj = BLASnrm2_(&blasj,&mfqP->work2[mfqP->nmodelpoints],&ione);
498 
499       if (proj >= mfqP->theta1) { /* add this index to model */
500         mfqP->model_indices[mfqP->nmodelpoints]=i;
501         mfqP->nmodelpoints++;
502         BLAScopy_(&blasn,mfqP->work,&ione,&mfqP->Q_tmp[mfqP->npmax*(mfqP->nmodelpoints-1)],&ione);
503         blask=mfqP->npmax*(mfqP->nmodelpoints);
504         BLAScopy_(&blask,mfqP->Q_tmp,&ione,mfqP->Q,&ione);
505         blask = mfqP->nmodelpoints;
506         blasmaxmn = PetscMax(mfqP->m,mfqP->n);
507         LAPACKgeqrf_(&blasn,&blask,mfqP->Q,&blasnpmax,mfqP->tau,mfqP->mwork,&blasmaxmn,&info);
508         if (info < 0) SETERRQ1(PETSC_COMM_SELF,1,"geqrf returned value %d\n",info);
509         mfqP->q_is_I = 0;
510       }
511       if (mfqP->nmodelpoints == mfqP->n)  {
512         break;
513       }
514     }
515   }
516   PetscFunctionReturn(0);
517 }
518 
519 #undef __FUNCT__
520 #define __FUNCT__ "TaoSolve_POUNDERS"
521 static PetscErrorCode TaoSolve_POUNDERS(TaoSolver tao)
522 {
523   TAO_POUNDERS               *mfqP = (TAO_POUNDERS *)tao->data;
524   PetscInt                   i,ii,j,k,l,iter=0;
525   PetscReal                  step=1.0;
526   TaoSolverTerminationReason reason = TAO_CONTINUE_ITERATING;
527   PetscInt                   low,high;
528   PetscReal                  minnorm;
529   PetscReal                  *x,*f,*fmin,*xmint;
530   PetscReal                  cres,deltaold;
531   PetscReal                  gnorm;
532   PetscBLASInt               info,ione=1,iblas;
533   PetscBool                  valid;
534   PetscReal                  mdec, rho, normxsp;
535   PetscReal                  one=1.0,zero=0.0,ratio;
536   PetscBLASInt               blasm,blasn,blasnpmax,blasn2;
537   PetscErrorCode             ierr;
538 
539   /* n = # of parameters
540      m = dimension (components) of function  */
541   PetscFunctionBegin;
542   if (tao->XL && tao->XU) {
543     /* Check x0 <= XU */
544     PetscReal val;
545     ierr = VecCopy(tao->solution,mfqP->Xhist[0]);CHKERRQ(ierr);
546     ierr = VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);CHKERRQ(ierr);
547     ierr = VecMax(mfqP->Xhist[0],NULL,&val);CHKERRQ(ierr);
548     if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 > upper bound");
549 
550     /* Check x0 >= xl */
551     ierr = VecCopy(tao->XL,mfqP->Xhist[0]);CHKERRQ(ierr);
552     ierr = VecAXPY(mfqP->Xhist[0],-1.0,tao->solution);CHKERRQ(ierr);
553     ierr = VecMax(mfqP->Xhist[0],NULL,&val);CHKERRQ(ierr);
554     if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 < lower bound");
555 
556     /* Check x0 + delta < XU  -- should be able to get around this eventually */
557 
558     ierr = VecSet(mfqP->Xhist[0],mfqP->delta);CHKERRQ(ierr);
559     ierr = VecAXPY(mfqP->Xhist[0],1.0,tao->solution);CHKERRQ(ierr);
560     ierr = VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);CHKERRQ(ierr);
561     ierr = VecMax(mfqP->Xhist[0],NULL,&val);CHKERRQ(ierr);
562     if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 + delta > upper bound");
563   }
564 
565   CHKMEMQ;
566   blasm = mfqP->m; blasn=mfqP->n; blasnpmax = mfqP->npmax;
567   for (i=0;i<mfqP->n*mfqP->n*mfqP->m;i++) mfqP->H[i]=0;
568 
569   ierr = VecCopy(tao->solution,mfqP->Xhist[0]);CHKERRQ(ierr);
570   CHKMEMQ;
571   ierr = TaoComputeSeparableObjective(tao,tao->solution,mfqP->Fhist[0]);CHKERRQ(ierr);
572 
573   ierr = VecNorm(mfqP->Fhist[0],NORM_2,&mfqP->Fres[0]);CHKERRQ(ierr);
574   if (PetscIsInfOrNanReal(mfqP->Fres[0])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
575   mfqP->Fres[0]*=mfqP->Fres[0];
576   mfqP->minindex = 0;
577   minnorm = mfqP->Fres[0];
578   ierr = VecGetOwnershipRange(mfqP->Xhist[0],&low,&high);CHKERRQ(ierr);
579   for (i=1;i<mfqP->n+1;i++) {
580     ierr = VecCopy(tao->solution,mfqP->Xhist[i]);CHKERRQ(ierr);
581     if (i-1 >= low && i-1 < high) {
582       ierr = VecGetArray(mfqP->Xhist[i],&x);CHKERRQ(ierr);
583       x[i-1-low] += mfqP->delta;
584       ierr = VecRestoreArray(mfqP->Xhist[i],&x);CHKERRQ(ierr);
585     }
586     CHKMEMQ;
587     ierr = TaoComputeSeparableObjective(tao,mfqP->Xhist[i],mfqP->Fhist[i]);CHKERRQ(ierr);
588     ierr = VecNorm(mfqP->Fhist[i],NORM_2,&mfqP->Fres[i]);CHKERRQ(ierr);
589     if (PetscIsInfOrNanReal(mfqP->Fres[i])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
590     mfqP->Fres[i]*=mfqP->Fres[i];
591     if (mfqP->Fres[i] < minnorm) {
592       mfqP->minindex = i;
593       minnorm = mfqP->Fres[i];
594     }
595   }
596 
597   ierr = VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);CHKERRQ(ierr);
598   ierr = VecCopy(mfqP->Fhist[mfqP->minindex],tao->sep_objective);CHKERRQ(ierr);
599   /* Gather mpi vecs to one big local vec */
600 
601   /* Begin serial code */
602 
603   /* Disp[i] = Xi-xmin, i=1,..,mfqP->minindex-1,mfqP->minindex+1,..,n */
604   /* Fdiff[i] = (Fi-Fmin)', i=1,..,mfqP->minindex-1,mfqP->minindex+1,..,n */
605   /* (Column oriented for blas calls) */
606   ii=0;
607 
608   if (mfqP->size == 1) {
609     ierr = VecGetArray(mfqP->Xhist[mfqP->minindex],&xmint);CHKERRQ(ierr);
610     for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i];
611     ierr = VecRestoreArray(mfqP->Xhist[mfqP->minindex],&xmint);CHKERRQ(ierr);
612     ierr = VecGetArray(mfqP->Fhist[mfqP->minindex],&fmin);CHKERRQ(ierr);
613     for (i=0;i<mfqP->n+1;i++) {
614       if (i == mfqP->minindex) continue;
615 
616       ierr = VecGetArray(mfqP->Xhist[i],&x);CHKERRQ(ierr);
617       for (j=0;j<mfqP->n;j++) {
618         mfqP->Disp[ii+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta;
619       }
620       ierr = VecRestoreArray(mfqP->Xhist[i],&x);CHKERRQ(ierr);
621 
622       ierr = VecGetArray(mfqP->Fhist[i],&f);CHKERRQ(ierr);
623       for (j=0;j<mfqP->m;j++) {
624         mfqP->Fdiff[ii+mfqP->n*j] = f[j] - fmin[j];
625       }
626       ierr = VecRestoreArray(mfqP->Fhist[i],&f);CHKERRQ(ierr);
627       mfqP->model_indices[ii++] = i;
628 
629     }
630     for (j=0;j<mfqP->m;j++) {
631       mfqP->C[j] = fmin[j];
632     }
633     ierr = VecRestoreArray(mfqP->Fhist[mfqP->minindex],&fmin);CHKERRQ(ierr);
634   } else {
635     ierr = VecScatterBegin(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
636     ierr = VecScatterEnd(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
637     ierr = VecGetArray(mfqP->localxmin,&xmint);CHKERRQ(ierr);
638     for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i];
639     ierr = VecRestoreArray(mfqP->localxmin,&mfqP->xmin);CHKERRQ(ierr);
640 
641     ierr = VecScatterBegin(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
642     ierr = VecScatterEnd(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
643     ierr = VecGetArray(mfqP->localfmin,&fmin);CHKERRQ(ierr);
644     for (i=0;i<mfqP->n+1;i++) {
645       if (i == mfqP->minindex) continue;
646 
647       mfqP->model_indices[ii++] = i;
648       ierr = VecScatterBegin(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
649       ierr = VecScatterEnd(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
650       ierr = VecGetArray(mfqP->localx,&x);CHKERRQ(ierr);
651       for (j=0;j<mfqP->n;j++) {
652         mfqP->Disp[i+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta;
653       }
654       ierr = VecRestoreArray(mfqP->localx,&x);CHKERRQ(ierr);
655 
656       ierr = VecScatterBegin(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
657       ierr = VecScatterEnd(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
658       ierr = VecGetArray(mfqP->localf,&f);CHKERRQ(ierr);
659       for (j=0;j<mfqP->m;j++) {
660         mfqP->Fdiff[i*mfqP->n+j] = f[j] - fmin[j];
661       }
662       ierr = VecRestoreArray(mfqP->localf,&f);CHKERRQ(ierr);
663     }
664     for (j=0;j<mfqP->m;j++) {
665       mfqP->C[j] = fmin[j];
666     }
667     ierr = VecRestoreArray(mfqP->localfmin,&fmin);CHKERRQ(ierr);
668   }
669 
670   /* Determine the initial quadratic models */
671   /* G = D(ModelIn,:) \ (F(ModelIn,1:m)-repmat(F(xkin,1:m),n,1)); */
672   /* D (nxn) Fdiff (nxm)  => G (nxm) */
673   blasn2 = blasn;
674   LAPACKgesv_(&blasn,&blasm,mfqP->Disp,&blasnpmax,mfqP->iwork,mfqP->Fdiff,&blasn2,&info);
675   ierr = PetscInfo1(tao,"gesv returned %d\n",info);CHKERRQ(ierr);
676 
677   cres = minnorm;
678   /* Gres = G*F(xkin,1:m)'  G (nxm)   Fk (m)   */
679   BLASgemv_("N",&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->C,&ione,&zero,mfqP->Gres,&ione);
680 
681   /*  Hres = G*G'  */
682   BLASgemm_("N","T",&blasn,&blasn,&blasm,&one,mfqP->Fdiff, &blasn,mfqP->Fdiff,&blasn,&zero,mfqP->Hres,&blasn);
683 
684   valid = PETSC_TRUE;
685 
686   ierr = VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);CHKERRQ(ierr);
687   ierr = VecAssemblyBegin(tao->gradient);
688   ierr = VecAssemblyEnd(tao->gradient);
689   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
690   gnorm *= mfqP->delta;
691   ierr = VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);CHKERRQ(ierr);
692   ierr = TaoMonitor(tao, iter, minnorm, gnorm, 0.0, step, &reason);CHKERRQ(ierr);
693   mfqP->nHist = mfqP->n+1;
694   mfqP->nmodelpoints = mfqP->n+1;
695 
696   while (reason == TAO_CONTINUE_ITERATING) {
697     PetscReal gnm = 1e-4;
698     iter++;
699     /* Solve the subproblem min{Q(s): ||s|| <= delta} */
700     ierr = gqtwrap(tao,&gnm,&mdec);CHKERRQ(ierr);
701     /* Evaluate the function at the new point */
702 
703     for (i=0;i<mfqP->n;i++) {
704         mfqP->work[i] = mfqP->Xsubproblem[i]*mfqP->delta + mfqP->xmin[i];
705     }
706     ierr = VecDuplicate(tao->solution,&mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr);
707     ierr = VecDuplicate(tao->sep_objective,&mfqP->Fhist[mfqP->nHist]);CHKERRQ(ierr);
708     ierr = VecSetValues(mfqP->Xhist[mfqP->nHist],mfqP->n,mfqP->indices,mfqP->work,INSERT_VALUES);CHKERRQ(ierr);
709     ierr = VecAssemblyBegin(mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr);
710     ierr = VecAssemblyEnd(mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr);
711 
712     ierr = TaoComputeSeparableObjective(tao,mfqP->Xhist[mfqP->nHist],mfqP->Fhist[mfqP->nHist]);CHKERRQ(ierr);
713     ierr = VecNorm(mfqP->Fhist[mfqP->nHist],NORM_2,&mfqP->Fres[mfqP->nHist]);CHKERRQ(ierr);
714     if (PetscIsInfOrNanReal(mfqP->Fres[mfqP->nHist])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
715     mfqP->Fres[mfqP->nHist]*=mfqP->Fres[mfqP->nHist];
716     rho = (mfqP->Fres[mfqP->minindex] - mfqP->Fres[mfqP->nHist]) / mdec;
717     mfqP->nHist++;
718 
719     /* Update the center */
720     if ((rho >= mfqP->eta1) || (rho > mfqP->eta0 && valid==PETSC_TRUE)) {
721       /* Update model to reflect new base point */
722       for (i=0;i<mfqP->n;i++) {
723         mfqP->work[i] = (mfqP->work[i] - mfqP->xmin[i])/mfqP->delta;
724       }
725       for (j=0;j<mfqP->m;j++) {
726         /* C(j) = C(j) + work*G(:,j) + .5*work*H(:,:,j)*work';
727          G(:,j) = G(:,j) + H(:,:,j)*work' */
728         for (k=0;k<mfqP->n;k++) {
729           mfqP->work2[k]=0.0;
730           for (l=0;l<mfqP->n;l++) {
731             mfqP->work2[k]+=mfqP->H[j + mfqP->m*(k + l*mfqP->n)]*mfqP->work[l];
732           }
733         }
734         for (i=0;i<mfqP->n;i++) {
735           mfqP->C[j]+=mfqP->work[i]*(mfqP->Fdiff[i + mfqP->n* j] + 0.5*mfqP->work2[i]);
736           mfqP->Fdiff[i+mfqP->n*j] +=mfqP-> work2[i];
737         }
738       }
739       /* Cres += work*Gres + .5*work*Hres*work';
740        Gres += Hres*work'; */
741 
742       BLASgemv_("N",&blasn,&blasn,&one,mfqP->Hres,&blasn,mfqP->work,&ione,&zero,mfqP->work2,&ione);
743       for (i=0;j<mfqP->n;j++) {
744         cres += mfqP->work[i]*(mfqP->Gres[i]  + 0.5*mfqP->work2[i]);
745         mfqP->Gres[i] += mfqP->work2[i];
746       }
747       mfqP->minindex = mfqP->nHist-1;
748       minnorm = mfqP->Fres[mfqP->minindex];
749       /* Change current center */
750       ierr = VecGetArray(mfqP->Xhist[mfqP->minindex],&xmint);CHKERRQ(ierr);
751       for (i=0;i<mfqP->n;i++) {
752         mfqP->xmin[i] = xmint[i];
753       }
754       ierr = VecRestoreArray(mfqP->Xhist[mfqP->minindex],&xmint);CHKERRQ(ierr);
755     }
756 
757     /* Evaluate at a model-improving point if necessary */
758     if (valid == PETSC_FALSE) {
759       mfqP->q_is_I = 1;
760       ierr = affpoints(mfqP,mfqP->xmin,mfqP->c1);CHKERRQ(ierr);
761       if (mfqP->nmodelpoints < mfqP->n) {
762         ierr = PetscInfo(tao,"Model not valid -- model-improving");
763         ierr = modelimprove(tao,mfqP,1);CHKERRQ(ierr);
764       }
765     }
766 
767     /* Update the trust region radius */
768     deltaold = mfqP->delta;
769     normxsp = 0;
770     for (i=0;i<mfqP->n;i++) {
771       normxsp += mfqP->Xsubproblem[i]*mfqP->Xsubproblem[i];
772     }
773     normxsp = PetscSqrtReal(normxsp);
774     if (rho >= mfqP->eta1 && normxsp > 0.5*mfqP->delta) {
775       mfqP->delta = PetscMin(mfqP->delta*mfqP->gamma1,mfqP->deltamax);
776     } else if (valid == PETSC_TRUE) {
777       mfqP->delta = PetscMax(mfqP->delta*mfqP->gamma0,mfqP->deltamin);
778     }
779 
780     /* Compute the next interpolation set */
781     mfqP->q_is_I = 1;
782     mfqP->nmodelpoints=0;
783     ierr = affpoints(mfqP,mfqP->xmin,mfqP->c1);CHKERRQ(ierr);
784     if (mfqP->nmodelpoints == mfqP->n) {
785       valid = PETSC_TRUE;
786     } else {
787       valid = PETSC_FALSE;
788       ierr = affpoints(mfqP,mfqP->xmin,mfqP->c2);CHKERRQ(ierr);
789       if (mfqP->n > mfqP->nmodelpoints) {
790         ierr = PetscInfo(tao,"Model not valid -- adding geometry points");
791         ierr = modelimprove(tao,mfqP,mfqP->n - mfqP->nmodelpoints);CHKERRQ(ierr);
792       }
793     }
794     for (i=mfqP->nmodelpoints;i>0;i--) {
795       mfqP->model_indices[i] = mfqP->model_indices[i-1];
796     }
797     mfqP->nmodelpoints++;
798     mfqP->model_indices[0] = mfqP->minindex;
799     ierr = morepoints(mfqP);CHKERRQ(ierr);
800     for (i=0;i<mfqP->nmodelpoints;i++) {
801       ierr = VecGetArray(mfqP->Xhist[mfqP->model_indices[i]],&x);CHKERRQ(ierr);
802       for (j=0;j<mfqP->n;j++) {
803         mfqP->Disp[i + mfqP->npmax*j] = (x[j]  - mfqP->xmin[j]) / deltaold;
804       }
805       ierr = VecRestoreArray(mfqP->Xhist[mfqP->model_indices[i]],&x);CHKERRQ(ierr);
806       ierr = VecGetArray(mfqP->Fhist[mfqP->model_indices[i]],&f);CHKERRQ(ierr);
807       for (j=0;j<mfqP->m;j++) {
808         for (k=0;k<mfqP->n;k++)  {
809           mfqP->work[k]=0.0;
810           for (l=0;l<mfqP->n;l++) {
811             mfqP->work[k] += mfqP->H[j + mfqP->m*(k + mfqP->n*l)] * mfqP->Disp[i + mfqP->npmax*l];
812           }
813         }
814         mfqP->RES[j*mfqP->npmax + i] = -mfqP->C[j] - BLASdot_(&blasn,&mfqP->Fdiff[j*mfqP->n],&ione,&mfqP->Disp[i],&blasnpmax) - 0.5*BLASdot_(&blasn,mfqP->work,&ione,&mfqP->Disp[i],&blasnpmax) + f[j];
815       }
816       ierr = VecRestoreArray(mfqP->Fhist[mfqP->model_indices[i]],&f);CHKERRQ(ierr);
817     }
818 
819     /* Update the quadratic model */
820     ierr = getquadpounders(mfqP);CHKERRQ(ierr);
821     ierr = VecGetArray(mfqP->Fhist[mfqP->minindex],&fmin);CHKERRQ(ierr);
822     BLAScopy_(&blasm,fmin,&ione,mfqP->C,&ione);
823     /* G = G*(delta/deltaold) + Gdel */
824     ratio = mfqP->delta/deltaold;
825     iblas = blasm*blasn;
826     BLASscal_(&iblas,&ratio,mfqP->Fdiff,&ione);
827     BLASaxpy_(&iblas,&one,mfqP->Gdel,&ione,mfqP->Fdiff,&ione);
828     /* H = H*(delta/deltaold) + Hdel */
829     iblas = blasm*blasn*blasn;
830     ratio *= ratio;
831     BLASscal_(&iblas,&ratio,mfqP->H,&ione);
832     BLASaxpy_(&iblas,&one,mfqP->Hdel,&ione,mfqP->H,&ione);
833 
834     /* Get residuals */
835     cres = mfqP->Fres[mfqP->minindex];
836     /* Gres = G*F(xkin,1:m)' */
837     BLASgemv_("N",&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->C,&ione,&zero,mfqP->Gres,&ione);
838     /* Hres = sum i=1..m {F(xkin,i)*H(:,:,i)}   + G*G' */
839     BLASgemm_("N","T",&blasn,&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->Fdiff,&blasn,&zero,mfqP->Hres,&blasn);
840 
841     iblas = mfqP->n*mfqP->n;
842 
843     for (j=0;j<mfqP->m;j++) {
844         BLASaxpy_(&iblas,&fmin[j],&mfqP->H[j],&blasm,mfqP->Hres,&ione);
845     }
846 
847     /* Export solution and gradient residual to TAO */
848     ierr = VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);CHKERRQ(ierr);
849     ierr = VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);CHKERRQ(ierr);
850     ierr = VecAssemblyBegin(tao->gradient);
851     ierr = VecAssemblyEnd(tao->gradient);
852     ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
853     gnorm *= mfqP->delta;
854     /*  final criticality test */
855     ierr = TaoMonitor(tao, iter, minnorm, gnorm, 0.0, step, &reason);CHKERRQ(ierr);
856   }
857   PetscFunctionReturn(0);
858 }
859 
860 #undef __FUNCT__
861 #define __FUNCT__ "TaoSetUp_POUNDERS"
862 static PetscErrorCode TaoSetUp_POUNDERS(TaoSolver tao)
863 {
864   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;
865   PetscInt       i;
866   IS             isfloc,isfglob,isxloc,isxglob;
867   PetscErrorCode ierr;
868 
869   PetscFunctionBegin;
870   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);  }
871   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);  }
872   ierr = VecGetSize(tao->solution,&mfqP->n);CHKERRQ(ierr);
873   ierr = VecGetSize(tao->sep_objective,&mfqP->m);CHKERRQ(ierr);
874   mfqP->c1 = PetscSqrtReal((PetscReal)mfqP->n);
875   if (mfqP->npmax == PETSC_DEFAULT) {
876     mfqP->npmax = 2*mfqP->n + 1;
877   }
878   mfqP->npmax = PetscMin((mfqP->n+1)*(mfqP->n+2)/2,mfqP->npmax);
879   mfqP->npmax = PetscMax(mfqP->npmax, mfqP->n+2);
880 
881   ierr = PetscMalloc1((tao->max_funcs+10),&mfqP->Xhist);CHKERRQ(ierr);
882   ierr = PetscMalloc1((tao->max_funcs+10),&mfqP->Fhist);CHKERRQ(ierr);
883   for (i=0;i<mfqP->n +1;i++) {
884     ierr = VecDuplicate(tao->solution,&mfqP->Xhist[i]);CHKERRQ(ierr);
885     ierr = VecDuplicate(tao->sep_objective,&mfqP->Fhist[i]);CHKERRQ(ierr);
886   }
887   ierr = VecDuplicate(tao->solution,&mfqP->workxvec);CHKERRQ(ierr);
888   mfqP->nHist = 0;
889 
890   ierr = PetscMalloc1((tao->max_funcs+10),&mfqP->Fres);CHKERRQ(ierr);
891   ierr = PetscMalloc1(mfqP->npmax*mfqP->m,&mfqP->RES);CHKERRQ(ierr);
892   ierr = PetscMalloc1(mfqP->n,&mfqP->work);CHKERRQ(ierr);
893   ierr = PetscMalloc1(mfqP->n,&mfqP->work2);CHKERRQ(ierr);
894   ierr = PetscMalloc1(mfqP->n,&mfqP->work3);CHKERRQ(ierr);
895   ierr = PetscMalloc1(PetscMax(mfqP->m,mfqP->n+1),&mfqP->mwork);CHKERRQ(ierr);
896   ierr = PetscMalloc1((mfqP->npmax - mfqP->n - 1),&mfqP->omega);CHKERRQ(ierr);
897   ierr = PetscMalloc1((mfqP->n * (mfqP->n+1) / 2),&mfqP->beta);CHKERRQ(ierr);
898   ierr = PetscMalloc1((mfqP->n + 1) ,&mfqP->alpha);CHKERRQ(ierr);
899 
900   ierr = PetscMalloc1(mfqP->n*mfqP->n*mfqP->m,&mfqP->H);CHKERRQ(ierr);
901   ierr = PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q);CHKERRQ(ierr);
902   ierr = PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q_tmp);CHKERRQ(ierr);
903   ierr = PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L);CHKERRQ(ierr);
904   ierr = PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_tmp);CHKERRQ(ierr);
905   ierr = PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_save);CHKERRQ(ierr);
906   ierr = PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->N);CHKERRQ(ierr);
907   ierr = PetscMalloc1(mfqP->npmax * (mfqP->n+1) ,&mfqP->M);CHKERRQ(ierr);
908   ierr = PetscMalloc1(mfqP->npmax * (mfqP->npmax - mfqP->n - 1) , &mfqP->Z);CHKERRQ(ierr);
909   ierr = PetscMalloc1(mfqP->npmax,&mfqP->tau);CHKERRQ(ierr);
910   ierr = PetscMalloc1(mfqP->npmax,&mfqP->tau_tmp);CHKERRQ(ierr);
911   mfqP->nmax = PetscMax(5*mfqP->npmax,mfqP->n*(mfqP->n+1)/2);
912   ierr = PetscMalloc1(mfqP->nmax,&mfqP->npmaxwork);CHKERRQ(ierr);
913   ierr = PetscMalloc1(mfqP->nmax,&mfqP->npmaxiwork);CHKERRQ(ierr);
914   ierr = PetscMalloc1(mfqP->n,&mfqP->xmin);CHKERRQ(ierr);
915   ierr = PetscMalloc1(mfqP->m,&mfqP->C);CHKERRQ(ierr);
916   ierr = PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Fdiff);CHKERRQ(ierr);
917   ierr = PetscMalloc1(mfqP->npmax*mfqP->n,&mfqP->Disp);CHKERRQ(ierr);
918   ierr = PetscMalloc1(mfqP->n,&mfqP->Gres);CHKERRQ(ierr);
919   ierr = PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Hres);CHKERRQ(ierr);
920   ierr = PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Gpoints);CHKERRQ(ierr);
921   ierr = PetscMalloc1(mfqP->npmax,&mfqP->model_indices);CHKERRQ(ierr);
922   ierr = PetscMalloc1(mfqP->n,&mfqP->Xsubproblem);CHKERRQ(ierr);
923   ierr = PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Gdel);CHKERRQ(ierr);
924   ierr = PetscMalloc1(mfqP->n*mfqP->n*mfqP->m, &mfqP->Hdel);CHKERRQ(ierr);
925   ierr = PetscMalloc1(PetscMax(mfqP->m,mfqP->n),&mfqP->indices);CHKERRQ(ierr);
926   ierr = PetscMalloc1(mfqP->n,&mfqP->iwork);CHKERRQ(ierr);
927   for (i=0;i<PetscMax(mfqP->m,mfqP->n);i++) {
928     mfqP->indices[i] = i;
929   }
930   ierr = MPI_Comm_size(((PetscObject)tao)->comm,&mfqP->size);CHKERRQ(ierr);
931   if (mfqP->size > 1) {
932     VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localx);CHKERRQ(ierr);
933     VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localxmin);CHKERRQ(ierr);
934     VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localf);CHKERRQ(ierr);
935     VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localfmin);CHKERRQ(ierr);
936     ierr = ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxloc);CHKERRQ(ierr);
937     ierr = ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxglob);CHKERRQ(ierr);
938     ierr = ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfloc);CHKERRQ(ierr);
939     ierr = ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfglob);CHKERRQ(ierr);
940 
941 
942     ierr = VecScatterCreate(tao->solution,isxglob,mfqP->localx,isxloc,&mfqP->scatterx);CHKERRQ(ierr);
943     ierr = VecScatterCreate(tao->sep_objective,isfglob,mfqP->localf,isfloc,&mfqP->scatterf);CHKERRQ(ierr);
944 
945     ierr = ISDestroy(&isxloc);CHKERRQ(ierr);
946     ierr = ISDestroy(&isxglob);CHKERRQ(ierr);
947     ierr = ISDestroy(&isfloc);CHKERRQ(ierr);
948     ierr = ISDestroy(&isfglob);CHKERRQ(ierr);
949   }
950 
951   if (!mfqP->usegqt) {
952     KSP       ksp;
953     PC        pc;
954     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,mfqP->n,mfqP->n,mfqP->Xsubproblem,&mfqP->subx);CHKERRQ(ierr);
955     ierr = VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->subxl);CHKERRQ(ierr);
956     ierr = VecDuplicate(mfqP->subxl,&mfqP->subb);CHKERRQ(ierr);
957     ierr = VecDuplicate(mfqP->subxl,&mfqP->subxu);CHKERRQ(ierr);
958     ierr = VecDuplicate(mfqP->subxl,&mfqP->subpdel);CHKERRQ(ierr);
959     ierr = VecDuplicate(mfqP->subxl,&mfqP->subndel);CHKERRQ(ierr);
960     ierr = TaoCreate(PETSC_COMM_SELF,&mfqP->subtao);CHKERRQ(ierr);
961     /* ierr = TaoSetType(mfqP->subtao,"tao_bqpip");CHKERRQ(ierr); */
962     ierr = TaoSetType(mfqP->subtao,"tao_tron");CHKERRQ(ierr);
963     ierr = TaoSetOptionsPrefix(mfqP->subtao,"pounders_subsolver_");CHKERRQ(ierr);
964     ierr = TaoSetInitialVector(mfqP->subtao,mfqP->subx);CHKERRQ(ierr);
965     ierr = TaoSetObjectiveAndGradientRoutine(mfqP->subtao,pounders_fg,(void*)mfqP);CHKERRQ(ierr);
966     ierr = TaoSetMaximumIterations(mfqP->subtao,mfqP->gqt_maxits);CHKERRQ(ierr);
967     ierr = TaoGetKSP(mfqP->subtao,&ksp);CHKERRQ(ierr);
968     if (ksp) {
969       ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
970       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
971     }
972     ierr = TaoSetVariableBounds(mfqP->subtao,mfqP->subxl,mfqP->subxu);CHKERRQ(ierr);
973     ierr = TaoSetFromOptions(mfqP->subtao);CHKERRQ(ierr);
974     ierr = MatCreateSeqDense(PETSC_COMM_SELF,mfqP->n,mfqP->n,mfqP->Hres,&mfqP->subH);CHKERRQ(ierr);
975     ierr = TaoSetHessianRoutine(mfqP->subtao,mfqP->subH,mfqP->subH,pounders_h,(void*)mfqP);CHKERRQ(ierr);
976   }
977   PetscFunctionReturn(0);
978 }
979 
980 #undef __FUNCT__
981 #define __FUNCT__ "TaoDestroy_POUNDERS"
982 static PetscErrorCode TaoDestroy_POUNDERS(TaoSolver tao)
983 {
984   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;
985   PetscInt       i;
986   PetscErrorCode ierr;
987 
988   PetscFunctionBegin;
989   if (!mfqP->usegqt) {
990     ierr = TaoDestroy(&mfqP->subtao);CHKERRQ(ierr);
991     ierr = VecDestroy(&mfqP->subx);CHKERRQ(ierr);
992     ierr = VecDestroy(&mfqP->subxl);CHKERRQ(ierr);
993     ierr = VecDestroy(&mfqP->subxu);CHKERRQ(ierr);
994     ierr = VecDestroy(&mfqP->subb);CHKERRQ(ierr);
995     ierr = VecDestroy(&mfqP->subpdel);CHKERRQ(ierr);
996     ierr = VecDestroy(&mfqP->subndel);CHKERRQ(ierr);
997     ierr = MatDestroy(&mfqP->subH);CHKERRQ(ierr);
998   }
999   ierr = PetscFree(mfqP->Fres);CHKERRQ(ierr);
1000   ierr = PetscFree(mfqP->RES);CHKERRQ(ierr);
1001   ierr = PetscFree(mfqP->work);CHKERRQ(ierr);
1002   ierr = PetscFree(mfqP->work2);CHKERRQ(ierr);
1003   ierr = PetscFree(mfqP->work3);CHKERRQ(ierr);
1004   ierr = PetscFree(mfqP->mwork);CHKERRQ(ierr);
1005   ierr = PetscFree(mfqP->omega);CHKERRQ(ierr);
1006   ierr = PetscFree(mfqP->beta);CHKERRQ(ierr);
1007   ierr = PetscFree(mfqP->alpha);CHKERRQ(ierr);
1008   ierr = PetscFree(mfqP->H);CHKERRQ(ierr);
1009   ierr = PetscFree(mfqP->Q);CHKERRQ(ierr);
1010   ierr = PetscFree(mfqP->Q_tmp);CHKERRQ(ierr);
1011   ierr = PetscFree(mfqP->L);CHKERRQ(ierr);
1012   ierr = PetscFree(mfqP->L_tmp);CHKERRQ(ierr);
1013   ierr = PetscFree(mfqP->L_save);CHKERRQ(ierr);
1014   ierr = PetscFree(mfqP->N);CHKERRQ(ierr);
1015   ierr = PetscFree(mfqP->M);CHKERRQ(ierr);
1016   ierr = PetscFree(mfqP->Z);CHKERRQ(ierr);
1017   ierr = PetscFree(mfqP->tau);CHKERRQ(ierr);
1018   ierr = PetscFree(mfqP->tau_tmp);CHKERRQ(ierr);
1019   ierr = PetscFree(mfqP->npmaxwork);CHKERRQ(ierr);
1020   ierr = PetscFree(mfqP->npmaxiwork);CHKERRQ(ierr);
1021   ierr = PetscFree(mfqP->xmin);CHKERRQ(ierr);
1022   ierr = PetscFree(mfqP->C);CHKERRQ(ierr);
1023   ierr = PetscFree(mfqP->Fdiff);CHKERRQ(ierr);
1024   ierr = PetscFree(mfqP->Disp);CHKERRQ(ierr);
1025   ierr = PetscFree(mfqP->Gres);CHKERRQ(ierr);
1026   ierr = PetscFree(mfqP->Hres);CHKERRQ(ierr);
1027   ierr = PetscFree(mfqP->Gpoints);CHKERRQ(ierr);
1028   ierr = PetscFree(mfqP->model_indices);CHKERRQ(ierr);
1029   ierr = PetscFree(mfqP->Xsubproblem);CHKERRQ(ierr);
1030   ierr = PetscFree(mfqP->Gdel);CHKERRQ(ierr);
1031   ierr = PetscFree(mfqP->Hdel);CHKERRQ(ierr);
1032   ierr = PetscFree(mfqP->indices);CHKERRQ(ierr);
1033   ierr = PetscFree(mfqP->iwork);CHKERRQ(ierr);
1034 
1035   for (i=0;i<mfqP->nHist;i++) {
1036     ierr = VecDestroy(&mfqP->Xhist[i]);CHKERRQ(ierr);
1037     ierr = VecDestroy(&mfqP->Fhist[i]);CHKERRQ(ierr);
1038   }
1039   if (mfqP->workxvec) {
1040     ierr = VecDestroy(&mfqP->workxvec);CHKERRQ(ierr);
1041   }
1042   ierr = PetscFree(mfqP->Xhist);CHKERRQ(ierr);
1043   ierr = PetscFree(mfqP->Fhist);CHKERRQ(ierr);
1044 
1045   if (mfqP->size > 1) {
1046     ierr = VecDestroy(&mfqP->localx);CHKERRQ(ierr);
1047     ierr = VecDestroy(&mfqP->localxmin);CHKERRQ(ierr);
1048     ierr = VecDestroy(&mfqP->localf);CHKERRQ(ierr);
1049     ierr = VecDestroy(&mfqP->localfmin);CHKERRQ(ierr);
1050   }
1051   ierr = PetscFree(tao->data);CHKERRQ(ierr);
1052   PetscFunctionReturn(0);
1053 }
1054 
1055 #undef __FUNCT__
1056 #define __FUNCT__ "TaoSetFromOptions_POUNDERS"
1057 static PetscErrorCode TaoSetFromOptions_POUNDERS(TaoSolver tao)
1058 {
1059   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;
1060   PetscErrorCode ierr;
1061 
1062   PetscFunctionBegin;
1063   ierr = PetscOptionsHead("POUNDERS method for least-squares optimization");CHKERRQ(ierr);
1064   ierr = PetscOptionsReal("-tao_pounders_delta","initial delta","",mfqP->delta,&mfqP->delta,0);CHKERRQ(ierr);
1065   mfqP->npmax = PETSC_DEFAULT;
1066   ierr = PetscOptionsInt("-tao_pounders_npmax","max number of points in model","",mfqP->npmax,&mfqP->npmax,0);CHKERRQ(ierr);
1067   mfqP->usegqt = PETSC_FALSE;
1068   ierr = PetscOptionsBool("-tao_pounders_gqt","use gqt algorithm for subproblem","",mfqP->usegqt,&mfqP->usegqt,0);CHKERRQ(ierr);
1069   ierr = PetscOptionsTail();CHKERRQ(ierr);
1070   PetscFunctionReturn(0);
1071 }
1072 
1073 #undef __FUNCT__
1074 #define __FUNCT__ "TaoView_POUNDERS"
1075 static PetscErrorCode TaoView_POUNDERS(TaoSolver tao, PetscViewer viewer)
1076 {
1077   TAO_POUNDERS   *mfqP = (TAO_POUNDERS *)tao->data;
1078   PetscBool      isascii;
1079   PetscErrorCode ierr;
1080 
1081   PetscFunctionBegin;
1082   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1083   if (isascii) {
1084     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1085     if (mfqP->usegqt) {
1086       ierr = PetscViewerASCIIPrintf(viewer, "Subproblem solver: gqt\n");CHKERRQ(ierr);
1087     } else {
1088       ierr = PetscViewerASCIIPrintf(viewer, "Subproblem solver: tron\n");CHKERRQ(ierr);
1089     }
1090     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1091   }
1092   PetscFunctionReturn(0);
1093 }
1094 
1095 EXTERN_C_BEGIN
1096 #undef __FUNCT__
1097 #define __FUNCT__ "TaoCreate_POUNDERS"
1098 PetscErrorCode TaoCreate_POUNDERS(TaoSolver tao)
1099 {
1100   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;
1101   PetscErrorCode ierr;
1102 
1103   PetscFunctionBegin;
1104   tao->ops->setup = TaoSetUp_POUNDERS;
1105   tao->ops->solve = TaoSolve_POUNDERS;
1106   tao->ops->view = TaoView_POUNDERS;
1107   tao->ops->setfromoptions = TaoSetFromOptions_POUNDERS;
1108   tao->ops->destroy = TaoDestroy_POUNDERS;
1109 
1110   ierr = PetscNewLog(tao,&mfqP);CHKERRQ(ierr);
1111   tao->data = (void*)mfqP;
1112   tao->max_it = 2000;
1113   tao->max_funcs = 4000;
1114   tao->fatol = 1e-8;
1115   tao->frtol = 1e-8;
1116   mfqP->delta = 0.1;
1117   mfqP->deltamax=1e3;
1118   mfqP->deltamin=1e-6;
1119   mfqP->c2 = 100.0;
1120   mfqP->theta1=1e-5;
1121   mfqP->theta2=1e-4;
1122   mfqP->gamma0=0.5;
1123   mfqP->gamma1=2.0;
1124   mfqP->eta0 = 0.0;
1125   mfqP->eta1 = 0.1;
1126   mfqP->gqt_rtol = 0.001;
1127   mfqP->gqt_maxits = 50;
1128   mfqP->workxvec = 0;
1129   PetscFunctionReturn(0);
1130 }
1131 EXTERN_C_END
1132 
1133