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