xref: /petsc/src/tao/leastsquares/impls/pounders/pounders.c (revision babe810ccf84b2c80abb81cf6f7b7ea2cb9ec828)
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, 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(Tao 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(Tao tao,PetscReal *gnorm, PetscReal *qmin)
36 {
37   PetscErrorCode ierr;
38 #if defined(PETSC_USE_REAL_SINGLE)
39   PetscReal      atol=1.0e-5;
40 #else
41   PetscReal      atol=1.0e-10;
42 #endif
43   PetscInt       info,its;
44   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;
45 
46   PetscFunctionBegin;
47   if (! mfqP->usegqt) {
48     PetscReal maxval;
49     PetscInt  i,j;
50 
51     ierr = VecSetValues(mfqP->subb,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);CHKERRQ(ierr);
52     ierr = VecAssemblyBegin(mfqP->subb);CHKERRQ(ierr);
53     ierr = VecAssemblyEnd(mfqP->subb);CHKERRQ(ierr);
54 
55     ierr = VecSet(mfqP->subx,0.0);CHKERRQ(ierr);
56 
57     ierr = VecSet(mfqP->subndel,-mfqP->delta);CHKERRQ(ierr);
58     ierr = VecSet(mfqP->subpdel,mfqP->delta);CHKERRQ(ierr);
59 
60     for (i=0;i<mfqP->n;i++) {
61       for (j=i;j<mfqP->n;j++) {
62         mfqP->Hres[j+mfqP->n*i] = mfqP->Hres[mfqP->n*j+i];
63       }
64     }
65     ierr = MatSetValues(mfqP->subH,mfqP->n,mfqP->indices,mfqP->n,mfqP->indices,mfqP->Hres,INSERT_VALUES);
66     ierr = MatAssemblyBegin(mfqP->subH,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
67     ierr = MatAssemblyEnd(mfqP->subH,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68 
69     ierr = TaoResetStatistics(mfqP->subtao);CHKERRQ(ierr);
70     ierr = TaoSetTolerances(mfqP->subtao,PETSC_DEFAULT,PETSC_DEFAULT,*gnorm,*gnorm,PETSC_DEFAULT);CHKERRQ(ierr);
71     /* enforce bound constraints -- experimental */
72     if (tao->XU && tao->XL) {
73       ierr = VecCopy(tao->XU,mfqP->subxu);CHKERRQ(ierr);
74       ierr = VecAXPY(mfqP->subxu,-1.0,tao->solution);CHKERRQ(ierr);
75       ierr = VecScale(mfqP->subxu,1.0/mfqP->delta);CHKERRQ(ierr);
76       ierr = VecCopy(tao->XL,mfqP->subxl);CHKERRQ(ierr);
77       ierr = VecAXPY(mfqP->subxl,-1.0,tao->solution);CHKERRQ(ierr);
78       ierr = VecScale(mfqP->subxl,1.0/mfqP->delta);CHKERRQ(ierr);
79 
80       ierr = VecPointwiseMin(mfqP->subxu,mfqP->subxu,mfqP->subpdel);CHKERRQ(ierr);
81       ierr = VecPointwiseMax(mfqP->subxl,mfqP->subxl,mfqP->subndel);CHKERRQ(ierr);
82     } else {
83       ierr = VecCopy(mfqP->subpdel,mfqP->subxu);CHKERRQ(ierr);
84       ierr = VecCopy(mfqP->subndel,mfqP->subxl);CHKERRQ(ierr);
85     }
86     /* Make sure xu > xl */
87     ierr = VecCopy(mfqP->subxl,mfqP->subpdel);CHKERRQ(ierr);
88     ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);CHKERRQ(ierr);
89     ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr);
90     if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"upper bound < lower bound in subproblem");
91     /* Make sure xu > tao->solution > xl */
92     ierr = VecCopy(mfqP->subxl,mfqP->subpdel);CHKERRQ(ierr);
93     ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subx);CHKERRQ(ierr);
94     ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr);
95     if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"initial guess < lower bound in subproblem");
96 
97     ierr = VecCopy(mfqP->subx,mfqP->subpdel);CHKERRQ(ierr);
98     ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);CHKERRQ(ierr);
99     ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr);
100     if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"initial guess > upper bound in subproblem");
101 
102     ierr = TaoSolve(mfqP->subtao);CHKERRQ(ierr);
103     ierr = TaoGetSolutionStatus(mfqP->subtao,NULL,qmin,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
104 
105     /* test bounds post-solution*/
106     ierr = VecCopy(mfqP->subxl,mfqP->subpdel);CHKERRQ(ierr);
107     ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subx);CHKERRQ(ierr);
108     ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr);
109     if (maxval > 1e-5) {
110       ierr = PetscInfo(tao,"subproblem solution < lower bound");CHKERRQ(ierr);
111       tao->reason = TAO_DIVERGED_TR_REDUCTION;
112     }
113 
114     ierr = VecCopy(mfqP->subx,mfqP->subpdel);CHKERRQ(ierr);
115     ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);CHKERRQ(ierr);
116     ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr);
117     if (maxval > 1e-5) {
118       ierr = PetscInfo(tao,"subproblem solution > upper bound");
119       tao->reason = TAO_DIVERGED_TR_REDUCTION;
120     }
121   } else {
122     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);
123   }
124   *qmin *= -1;
125   PetscFunctionReturn(0);
126 }
127 
128 #undef __FUNCT__
129 #define __FUNCT__ "phi2eval"
130 PetscErrorCode phi2eval(PetscReal *x, PetscInt n, PetscReal *phi)
131 {
132 /* 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] */
133   PetscInt  i,j,k;
134   PetscReal sqrt2 = PetscSqrtReal(2.0);
135 
136   PetscFunctionBegin;
137   j=0;
138   for (i=0;i<n;i++) {
139     phi[j] = 0.5 * x[i]*x[i];
140     j++;
141     for (k=i+1;k<n;k++) {
142       phi[j]  = x[i]*x[k]/sqrt2;
143       j++;
144     }
145   }
146   PetscFunctionReturn(0);
147 }
148 
149 #undef __FUNCT__
150 #define __FUNCT__ "getquadpounders"
151 PetscErrorCode getquadpounders(TAO_POUNDERS *mfqP)
152 {
153 /* Computes the parameters of the quadratic Q(x) = c + g'*x + 0.5*x*G*x'
154    that satisfies the interpolation conditions Q(X[:,j]) = f(j)
155    for j=1,...,m and with a Hessian matrix of least Frobenius norm */
156 
157     /* NB --we are ignoring c */
158   PetscInt     i,j,k,num,np = mfqP->nmodelpoints;
159   PetscReal    one = 1.0,zero=0.0,negone=-1.0;
160   PetscBLASInt blasnpmax = mfqP->npmax;
161   PetscBLASInt blasnplus1 = mfqP->n+1;
162   PetscBLASInt blasnp = np;
163   PetscBLASInt blasint = mfqP->n*(mfqP->n+1) / 2;
164   PetscBLASInt blasint2 = np - mfqP->n-1;
165   PetscBLASInt info,ione=1;
166   PetscReal    sqrt2 = PetscSqrtReal(2.0);
167 
168   PetscFunctionBegin;
169   for (i=0;i<mfqP->n*mfqP->m;i++) {
170     mfqP->Gdel[i] = 0;
171   }
172   for (i=0;i<mfqP->n*mfqP->n*mfqP->m;i++) {
173     mfqP->Hdel[i] = 0;
174   }
175 
176     /* factor M */
177   PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&blasnplus1,&blasnpmax,mfqP->M,&blasnplus1,mfqP->npmaxiwork,&info));
178   if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine getrf returned with value %d\n",info);
179 
180   if (np == mfqP->n+1) {
181     for (i=0;i<mfqP->npmax-mfqP->n-1;i++) {
182       mfqP->omega[i]=0.0;
183     }
184     for (i=0;i<mfqP->n*(mfqP->n+1)/2;i++) {
185       mfqP->beta[i]=0.0;
186     }
187   } else {
188     /* Let Ltmp = (L'*L) */
189     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));
190 
191     /* factor Ltmp */
192     PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&blasint2,mfqP->L_tmp,&blasint,&info));
193     if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine potrf returned with value %d\n",info);
194   }
195 
196   for (k=0;k<mfqP->m;k++) {
197     if (np != mfqP->n+1) {
198       /* Solve L'*L*Omega = Z' * RESk*/
199       PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&blasnp,&blasint2,&one,mfqP->Z,&blasnpmax,&mfqP->RES[mfqP->npmax*k],&ione,&zero,mfqP->omega,&ione));
200       PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&blasint2,&ione,mfqP->L_tmp,&blasint,mfqP->omega,&blasint2,&info));
201       if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine potrs returned with value %d\n",info);
202 
203       /* Beta = L*Omega */
204       PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasint,&blasint2,&one,&mfqP->L[(mfqP->n+1)*blasint],&blasint,mfqP->omega,&ione,&zero,mfqP->beta,&ione));
205     }
206 
207     /* solve M'*Alpha = RESk - N'*Beta */
208     PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&blasint,&blasnp,&negone,mfqP->N,&blasint,mfqP->beta,&ione,&one,&mfqP->RES[mfqP->npmax*k],&ione));
209     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&blasnplus1,&ione,mfqP->M,&blasnplus1,mfqP->npmaxiwork,&mfqP->RES[mfqP->npmax*k],&blasnplus1,&info));
210     if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine getrs returned with value %d\n",info);
211 
212     /* Gdel(:,k) = Alpha(2:n+1) */
213     for (i=0;i<mfqP->n;i++) {
214       mfqP->Gdel[i + mfqP->n*k] = mfqP->RES[mfqP->npmax*k + i+1];
215     }
216 
217     /* Set Hdels */
218     num=0;
219     for (i=0;i<mfqP->n;i++) {
220       /* H[i,i,k] = Beta(num) */
221       mfqP->Hdel[(i*mfqP->n + i)*mfqP->m + k] = mfqP->beta[num];
222       num++;
223       for (j=i+1;j<mfqP->n;j++) {
224         /* H[i,j,k] = H[j,i,k] = Beta(num)/sqrt(2) */
225         mfqP->Hdel[(j*mfqP->n + i)*mfqP->m + k] = mfqP->beta[num]/sqrt2;
226         mfqP->Hdel[(i*mfqP->n + j)*mfqP->m + k] = mfqP->beta[num]/sqrt2;
227         num++;
228       }
229     }
230   }
231   PetscFunctionReturn(0);
232 }
233 
234 #undef __FUNCT__
235 #define __FUNCT__ "morepoints"
236 PetscErrorCode morepoints(TAO_POUNDERS *mfqP)
237 {
238   /* Assumes mfqP->model_indices[0]  is minimum index
239    Finishes adding points to mfqP->model_indices (up to npmax)
240    Computes L,Z,M,N
241    np is actual number of points in model (should equal npmax?) */
242   PetscInt       point,i,j,offset;
243   PetscInt       reject;
244   PetscBLASInt   blasn=mfqP->n,blasnpmax=mfqP->npmax,blasnplus1=mfqP->n+1,info,blasnmax=mfqP->nmax,blasint,blasint2,blasnp,blasmaxmn;
245   PetscReal      *x,normd;
246   PetscErrorCode ierr;
247 
248   PetscFunctionBegin;
249   /* Initialize M,N */
250   for (i=0;i<mfqP->n+1;i++) {
251     ierr = VecGetArray(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 = VecRestoreArray(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 = VecGetArray(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 = VecRestoreArray(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   PetscReal      *x;
482   PetscErrorCode ierr;
483 
484   PetscFunctionBegin;
485   for (i=mfqP->nHist-1;i>=0;i--) {
486     ierr = VecGetArray(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 = VecRestoreArray(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   TaoTerminationReason reason = TAO_CONTINUE_ITERATING;
532   PetscInt                   low,high;
533   PetscReal                  minnorm;
534   PetscReal                  *x,*f,*fmin,*xmint;
535   PetscReal                  cres,deltaold;
536   PetscReal                  gnorm;
537   PetscBLASInt               info,ione=1,iblas;
538   PetscBool                  valid;
539   PetscReal                  mdec, rho, normxsp;
540   PetscReal                  one=1.0,zero=0.0,ratio;
541   PetscBLASInt               blasm,blasn,blasnpmax,blasn2;
542   PetscErrorCode             ierr;
543 
544   /* n = # of parameters
545      m = dimension (components) of function  */
546   PetscFunctionBegin;
547   if (tao->XL && tao->XU) {
548     /* Check x0 <= XU */
549     PetscReal val;
550     ierr = VecCopy(tao->solution,mfqP->Xhist[0]);CHKERRQ(ierr);
551     ierr = VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);CHKERRQ(ierr);
552     ierr = VecMax(mfqP->Xhist[0],NULL,&val);CHKERRQ(ierr);
553     if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 > upper bound");
554 
555     /* Check x0 >= xl */
556     ierr = VecCopy(tao->XL,mfqP->Xhist[0]);CHKERRQ(ierr);
557     ierr = VecAXPY(mfqP->Xhist[0],-1.0,tao->solution);CHKERRQ(ierr);
558     ierr = VecMax(mfqP->Xhist[0],NULL,&val);CHKERRQ(ierr);
559     if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 < lower bound");
560 
561     /* Check x0 + delta < XU  -- should be able to get around this eventually */
562 
563     ierr = VecSet(mfqP->Xhist[0],mfqP->delta);CHKERRQ(ierr);
564     ierr = VecAXPY(mfqP->Xhist[0],1.0,tao->solution);CHKERRQ(ierr);
565     ierr = VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);CHKERRQ(ierr);
566     ierr = VecMax(mfqP->Xhist[0],NULL,&val);CHKERRQ(ierr);
567     if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 + delta > upper bound");
568   }
569 
570   CHKMEMQ;
571   blasm = mfqP->m; blasn=mfqP->n; blasnpmax = mfqP->npmax;
572   for (i=0;i<mfqP->n*mfqP->n*mfqP->m;i++) mfqP->H[i]=0;
573 
574   ierr = VecCopy(tao->solution,mfqP->Xhist[0]);CHKERRQ(ierr);
575   CHKMEMQ;
576   ierr = TaoComputeSeparableObjective(tao,tao->solution,mfqP->Fhist[0]);CHKERRQ(ierr);
577 
578   ierr = VecNorm(mfqP->Fhist[0],NORM_2,&mfqP->Fres[0]);CHKERRQ(ierr);
579   if (PetscIsInfOrNanReal(mfqP->Fres[0])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
580   mfqP->Fres[0]*=mfqP->Fres[0];
581   mfqP->minindex = 0;
582   minnorm = mfqP->Fres[0];
583   ierr = VecGetOwnershipRange(mfqP->Xhist[0],&low,&high);CHKERRQ(ierr);
584   for (i=1;i<mfqP->n+1;i++) {
585     ierr = VecCopy(tao->solution,mfqP->Xhist[i]);CHKERRQ(ierr);
586     if (i-1 >= low && i-1 < high) {
587       ierr = VecGetArray(mfqP->Xhist[i],&x);CHKERRQ(ierr);
588       x[i-1-low] += mfqP->delta;
589       ierr = VecRestoreArray(mfqP->Xhist[i],&x);CHKERRQ(ierr);
590     }
591     CHKMEMQ;
592     ierr = TaoComputeSeparableObjective(tao,mfqP->Xhist[i],mfqP->Fhist[i]);CHKERRQ(ierr);
593     ierr = VecNorm(mfqP->Fhist[i],NORM_2,&mfqP->Fres[i]);CHKERRQ(ierr);
594     if (PetscIsInfOrNanReal(mfqP->Fres[i])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
595     mfqP->Fres[i]*=mfqP->Fres[i];
596     if (mfqP->Fres[i] < minnorm) {
597       mfqP->minindex = i;
598       minnorm = mfqP->Fres[i];
599     }
600   }
601 
602   ierr = VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);CHKERRQ(ierr);
603   ierr = VecCopy(mfqP->Fhist[mfqP->minindex],tao->sep_objective);CHKERRQ(ierr);
604   /* Gather mpi vecs to one big local vec */
605 
606   /* Begin serial code */
607 
608   /* Disp[i] = Xi-xmin, i=1,..,mfqP->minindex-1,mfqP->minindex+1,..,n */
609   /* Fdiff[i] = (Fi-Fmin)', i=1,..,mfqP->minindex-1,mfqP->minindex+1,..,n */
610   /* (Column oriented for blas calls) */
611   ii=0;
612 
613   if (mfqP->size == 1) {
614     ierr = VecGetArray(mfqP->Xhist[mfqP->minindex],&xmint);CHKERRQ(ierr);
615     for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i];
616     ierr = VecRestoreArray(mfqP->Xhist[mfqP->minindex],&xmint);CHKERRQ(ierr);
617     ierr = VecGetArray(mfqP->Fhist[mfqP->minindex],&fmin);CHKERRQ(ierr);
618     for (i=0;i<mfqP->n+1;i++) {
619       if (i == mfqP->minindex) continue;
620 
621       ierr = VecGetArray(mfqP->Xhist[i],&x);CHKERRQ(ierr);
622       for (j=0;j<mfqP->n;j++) {
623         mfqP->Disp[ii+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta;
624       }
625       ierr = VecRestoreArray(mfqP->Xhist[i],&x);CHKERRQ(ierr);
626 
627       ierr = VecGetArray(mfqP->Fhist[i],&f);CHKERRQ(ierr);
628       for (j=0;j<mfqP->m;j++) {
629         mfqP->Fdiff[ii+mfqP->n*j] = f[j] - fmin[j];
630       }
631       ierr = VecRestoreArray(mfqP->Fhist[i],&f);CHKERRQ(ierr);
632       mfqP->model_indices[ii++] = i;
633 
634     }
635     for (j=0;j<mfqP->m;j++) {
636       mfqP->C[j] = fmin[j];
637     }
638     ierr = VecRestoreArray(mfqP->Fhist[mfqP->minindex],&fmin);CHKERRQ(ierr);
639   } else {
640     ierr = VecScatterBegin(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
641     ierr = VecScatterEnd(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
642     ierr = VecGetArray(mfqP->localxmin,&xmint);CHKERRQ(ierr);
643     for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i];
644     ierr = VecRestoreArray(mfqP->localxmin,&mfqP->xmin);CHKERRQ(ierr);
645 
646     ierr = VecScatterBegin(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
647     ierr = VecScatterEnd(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
648     ierr = VecGetArray(mfqP->localfmin,&fmin);CHKERRQ(ierr);
649     for (i=0;i<mfqP->n+1;i++) {
650       if (i == mfqP->minindex) continue;
651 
652       mfqP->model_indices[ii++] = i;
653       ierr = VecScatterBegin(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
654       ierr = VecScatterEnd(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
655       ierr = VecGetArray(mfqP->localx,&x);CHKERRQ(ierr);
656       for (j=0;j<mfqP->n;j++) {
657         mfqP->Disp[i+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta;
658       }
659       ierr = VecRestoreArray(mfqP->localx,&x);CHKERRQ(ierr);
660 
661       ierr = VecScatterBegin(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
662       ierr = VecScatterEnd(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
663       ierr = VecGetArray(mfqP->localf,&f);CHKERRQ(ierr);
664       for (j=0;j<mfqP->m;j++) {
665         mfqP->Fdiff[i*mfqP->n+j] = f[j] - fmin[j];
666       }
667       ierr = VecRestoreArray(mfqP->localf,&f);CHKERRQ(ierr);
668     }
669     for (j=0;j<mfqP->m;j++) {
670       mfqP->C[j] = fmin[j];
671     }
672     ierr = VecRestoreArray(mfqP->localfmin,&fmin);CHKERRQ(ierr);
673   }
674 
675   /* Determine the initial quadratic models */
676   /* G = D(ModelIn,:) \ (F(ModelIn,1:m)-repmat(F(xkin,1:m),n,1)); */
677   /* D (nxn) Fdiff (nxm)  => G (nxm) */
678   blasn2 = blasn;
679   PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&blasn,&blasm,mfqP->Disp,&blasnpmax,mfqP->iwork,mfqP->Fdiff,&blasn2,&info));
680   ierr = PetscInfo1(tao,"gesv returned %d\n",info);CHKERRQ(ierr);
681 
682   cres = minnorm;
683   /* Gres = G*F(xkin,1:m)'  G (nxm)   Fk (m)   */
684   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->C,&ione,&zero,mfqP->Gres,&ione));
685 
686   /*  Hres = G*G'  */
687   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&blasn,&blasn,&blasm,&one,mfqP->Fdiff, &blasn,mfqP->Fdiff,&blasn,&zero,mfqP->Hres,&blasn));
688 
689   valid = PETSC_TRUE;
690 
691   ierr = VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);CHKERRQ(ierr);
692   ierr = VecAssemblyBegin(tao->gradient);
693   ierr = VecAssemblyEnd(tao->gradient);
694   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
695   gnorm *= mfqP->delta;
696   ierr = VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);CHKERRQ(ierr);
697   ierr = TaoMonitor(tao, iter, minnorm, gnorm, 0.0, step, &reason);CHKERRQ(ierr);
698   mfqP->nHist = mfqP->n+1;
699   mfqP->nmodelpoints = mfqP->n+1;
700 
701   while (reason == TAO_CONTINUE_ITERATING) {
702     PetscReal gnm = 1e-4;
703     iter++;
704     /* Solve the subproblem min{Q(s): ||s|| <= delta} */
705     ierr = gqtwrap(tao,&gnm,&mdec);CHKERRQ(ierr);
706     /* Evaluate the function at the new point */
707 
708     for (i=0;i<mfqP->n;i++) {
709         mfqP->work[i] = mfqP->Xsubproblem[i]*mfqP->delta + mfqP->xmin[i];
710     }
711     ierr = VecDuplicate(tao->solution,&mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr);
712     ierr = VecDuplicate(tao->sep_objective,&mfqP->Fhist[mfqP->nHist]);CHKERRQ(ierr);
713     ierr = VecSetValues(mfqP->Xhist[mfqP->nHist],mfqP->n,mfqP->indices,mfqP->work,INSERT_VALUES);CHKERRQ(ierr);
714     ierr = VecAssemblyBegin(mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr);
715     ierr = VecAssemblyEnd(mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr);
716 
717     ierr = TaoComputeSeparableObjective(tao,mfqP->Xhist[mfqP->nHist],mfqP->Fhist[mfqP->nHist]);CHKERRQ(ierr);
718     ierr = VecNorm(mfqP->Fhist[mfqP->nHist],NORM_2,&mfqP->Fres[mfqP->nHist]);CHKERRQ(ierr);
719     if (PetscIsInfOrNanReal(mfqP->Fres[mfqP->nHist])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
720     mfqP->Fres[mfqP->nHist]*=mfqP->Fres[mfqP->nHist];
721     rho = (mfqP->Fres[mfqP->minindex] - mfqP->Fres[mfqP->nHist]) / mdec;
722     mfqP->nHist++;
723 
724     /* Update the center */
725     if ((rho >= mfqP->eta1) || (rho > mfqP->eta0 && valid==PETSC_TRUE)) {
726       /* Update model to reflect new base point */
727       for (i=0;i<mfqP->n;i++) {
728         mfqP->work[i] = (mfqP->work[i] - mfqP->xmin[i])/mfqP->delta;
729       }
730       for (j=0;j<mfqP->m;j++) {
731         /* C(j) = C(j) + work*G(:,j) + .5*work*H(:,:,j)*work';
732          G(:,j) = G(:,j) + H(:,:,j)*work' */
733         for (k=0;k<mfqP->n;k++) {
734           mfqP->work2[k]=0.0;
735           for (l=0;l<mfqP->n;l++) {
736             mfqP->work2[k]+=mfqP->H[j + mfqP->m*(k + l*mfqP->n)]*mfqP->work[l];
737           }
738         }
739         for (i=0;i<mfqP->n;i++) {
740           mfqP->C[j]+=mfqP->work[i]*(mfqP->Fdiff[i + mfqP->n* j] + 0.5*mfqP->work2[i]);
741           mfqP->Fdiff[i+mfqP->n*j] +=mfqP-> work2[i];
742         }
743       }
744       /* Cres += work*Gres + .5*work*Hres*work';
745        Gres += Hres*work'; */
746 
747       PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasn,&one,mfqP->Hres,&blasn,mfqP->work,&ione,&zero,mfqP->work2,&ione));
748       for (i=0;i<mfqP->n;i++) {
749         cres += mfqP->work[i]*(mfqP->Gres[i]  + 0.5*mfqP->work2[i]);
750         mfqP->Gres[i] += mfqP->work2[i];
751       }
752       mfqP->minindex = mfqP->nHist-1;
753       minnorm = mfqP->Fres[mfqP->minindex];
754       /* Change current center */
755       ierr = VecGetArray(mfqP->Xhist[mfqP->minindex],&xmint);CHKERRQ(ierr);
756       for (i=0;i<mfqP->n;i++) {
757         mfqP->xmin[i] = xmint[i];
758       }
759       ierr = VecRestoreArray(mfqP->Xhist[mfqP->minindex],&xmint);CHKERRQ(ierr);
760     }
761 
762     /* Evaluate at a model-improving point if necessary */
763     if (valid == PETSC_FALSE) {
764       mfqP->q_is_I = 1;
765       mfqP->nmodelpoints = 0;
766       ierr = affpoints(mfqP,mfqP->xmin,mfqP->c1);CHKERRQ(ierr);
767       if (mfqP->nmodelpoints < mfqP->n) {
768         ierr = PetscInfo(tao,"Model not valid -- model-improving");
769         ierr = modelimprove(tao,mfqP,1);CHKERRQ(ierr);
770       }
771     }
772 
773     /* Update the trust region radius */
774     deltaold = mfqP->delta;
775     normxsp = 0;
776     for (i=0;i<mfqP->n;i++) {
777       normxsp += mfqP->Xsubproblem[i]*mfqP->Xsubproblem[i];
778     }
779     normxsp = PetscSqrtReal(normxsp);
780     if (rho >= mfqP->eta1 && normxsp > 0.5*mfqP->delta) {
781       mfqP->delta = PetscMin(mfqP->delta*mfqP->gamma1,mfqP->deltamax);
782     } else if (valid == PETSC_TRUE) {
783       mfqP->delta = PetscMax(mfqP->delta*mfqP->gamma0,mfqP->deltamin);
784     }
785 
786     /* Compute the next interpolation set */
787     mfqP->q_is_I = 1;
788     mfqP->nmodelpoints=0;
789     ierr = affpoints(mfqP,mfqP->xmin,mfqP->c1);CHKERRQ(ierr);
790     if (mfqP->nmodelpoints == mfqP->n) {
791       valid = PETSC_TRUE;
792     } else {
793       valid = PETSC_FALSE;
794       ierr = affpoints(mfqP,mfqP->xmin,mfqP->c2);CHKERRQ(ierr);
795       if (mfqP->n > mfqP->nmodelpoints) {
796         ierr = PetscInfo(tao,"Model not valid -- adding geometry points");
797         ierr = modelimprove(tao,mfqP,mfqP->n - mfqP->nmodelpoints);CHKERRQ(ierr);
798       }
799     }
800     for (i=mfqP->nmodelpoints;i>0;i--) {
801       mfqP->model_indices[i] = mfqP->model_indices[i-1];
802     }
803     mfqP->nmodelpoints++;
804     mfqP->model_indices[0] = mfqP->minindex;
805     ierr = morepoints(mfqP);CHKERRQ(ierr);
806     for (i=0;i<mfqP->nmodelpoints;i++) {
807       ierr = VecGetArray(mfqP->Xhist[mfqP->model_indices[i]],&x);CHKERRQ(ierr);
808       for (j=0;j<mfqP->n;j++) {
809         mfqP->Disp[i + mfqP->npmax*j] = (x[j]  - mfqP->xmin[j]) / deltaold;
810       }
811       ierr = VecRestoreArray(mfqP->Xhist[mfqP->model_indices[i]],&x);CHKERRQ(ierr);
812       ierr = VecGetArray(mfqP->Fhist[mfqP->model_indices[i]],&f);CHKERRQ(ierr);
813       for (j=0;j<mfqP->m;j++) {
814         for (k=0;k<mfqP->n;k++)  {
815           mfqP->work[k]=0.0;
816           for (l=0;l<mfqP->n;l++) {
817             mfqP->work[k] += mfqP->H[j + mfqP->m*(k + mfqP->n*l)] * mfqP->Disp[i + mfqP->npmax*l];
818           }
819         }
820         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];
821       }
822       ierr = VecRestoreArray(mfqP->Fhist[mfqP->model_indices[i]],&f);CHKERRQ(ierr);
823     }
824 
825     /* Update the quadratic model */
826     ierr = getquadpounders(mfqP);CHKERRQ(ierr);
827     ierr = VecGetArray(mfqP->Fhist[mfqP->minindex],&fmin);CHKERRQ(ierr);
828     PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasm,fmin,&ione,mfqP->C,&ione));
829     /* G = G*(delta/deltaold) + Gdel */
830     ratio = mfqP->delta/deltaold;
831     iblas = blasm*blasn;
832     PetscStackCallBLAS("BLASscal",BLASscal_(&iblas,&ratio,mfqP->Fdiff,&ione));
833     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&iblas,&one,mfqP->Gdel,&ione,mfqP->Fdiff,&ione));
834     /* H = H*(delta/deltaold) + Hdel */
835     iblas = blasm*blasn*blasn;
836     ratio *= ratio;
837     PetscStackCallBLAS("BLASscal",BLASscal_(&iblas,&ratio,mfqP->H,&ione));
838     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&iblas,&one,mfqP->Hdel,&ione,mfqP->H,&ione));
839 
840     /* Get residuals */
841     cres = mfqP->Fres[mfqP->minindex];
842     /* Gres = G*F(xkin,1:m)' */
843     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->C,&ione,&zero,mfqP->Gres,&ione));
844     /* Hres = sum i=1..m {F(xkin,i)*H(:,:,i)}   + G*G' */
845     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&blasn,&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->Fdiff,&blasn,&zero,mfqP->Hres,&blasn));
846 
847     iblas = mfqP->n*mfqP->n;
848 
849     for (j=0;j<mfqP->m;j++) {
850       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&iblas,&fmin[j],&mfqP->H[j],&blasm,mfqP->Hres,&ione));
851     }
852 
853     /* Export solution and gradient residual to TAO */
854     ierr = VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);CHKERRQ(ierr);
855     ierr = VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);CHKERRQ(ierr);
856     ierr = VecAssemblyBegin(tao->gradient);
857     ierr = VecAssemblyEnd(tao->gradient);
858     ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
859     gnorm *= mfqP->delta;
860     /*  final criticality test */
861     ierr = TaoMonitor(tao, iter, minnorm, gnorm, 0.0, step, &reason);CHKERRQ(ierr);
862   }
863   PetscFunctionReturn(0);
864 }
865 
866 #undef __FUNCT__
867 #define __FUNCT__ "TaoSetUp_POUNDERS"
868 static PetscErrorCode TaoSetUp_POUNDERS(Tao tao)
869 {
870   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;
871   PetscInt       i;
872   IS             isfloc,isfglob,isxloc,isxglob;
873   PetscErrorCode ierr;
874 
875   PetscFunctionBegin;
876   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);  }
877   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);  }
878   ierr = VecGetSize(tao->solution,&mfqP->n);CHKERRQ(ierr);
879   ierr = VecGetSize(tao->sep_objective,&mfqP->m);CHKERRQ(ierr);
880   mfqP->c1 = PetscSqrtReal((PetscReal)mfqP->n);
881   if (mfqP->npmax == PETSC_DEFAULT) {
882     mfqP->npmax = 2*mfqP->n + 1;
883   }
884   mfqP->npmax = PetscMin((mfqP->n+1)*(mfqP->n+2)/2,mfqP->npmax);
885   mfqP->npmax = PetscMax(mfqP->npmax, mfqP->n+2);
886 
887   ierr = PetscMalloc1((tao->max_funcs+10),&mfqP->Xhist);CHKERRQ(ierr);
888   ierr = PetscMalloc1((tao->max_funcs+10),&mfqP->Fhist);CHKERRQ(ierr);
889   for (i=0;i<mfqP->n +1;i++) {
890     ierr = VecDuplicate(tao->solution,&mfqP->Xhist[i]);CHKERRQ(ierr);
891     ierr = VecDuplicate(tao->sep_objective,&mfqP->Fhist[i]);CHKERRQ(ierr);
892   }
893   ierr = VecDuplicate(tao->solution,&mfqP->workxvec);CHKERRQ(ierr);
894   mfqP->nHist = 0;
895 
896   ierr = PetscMalloc1((tao->max_funcs+10),&mfqP->Fres);CHKERRQ(ierr);
897   ierr = PetscMalloc1(mfqP->npmax*mfqP->m,&mfqP->RES);CHKERRQ(ierr);
898   ierr = PetscMalloc1(mfqP->n,&mfqP->work);CHKERRQ(ierr);
899   ierr = PetscMalloc1(mfqP->n,&mfqP->work2);CHKERRQ(ierr);
900   ierr = PetscMalloc1(mfqP->n,&mfqP->work3);CHKERRQ(ierr);
901   ierr = PetscMalloc1(PetscMax(mfqP->m,mfqP->n+1),&mfqP->mwork);CHKERRQ(ierr);
902   ierr = PetscMalloc1((mfqP->npmax - mfqP->n - 1),&mfqP->omega);CHKERRQ(ierr);
903   ierr = PetscMalloc1((mfqP->n * (mfqP->n+1) / 2),&mfqP->beta);CHKERRQ(ierr);
904   ierr = PetscMalloc1((mfqP->n + 1) ,&mfqP->alpha);CHKERRQ(ierr);
905 
906   ierr = PetscMalloc1(mfqP->n*mfqP->n*mfqP->m,&mfqP->H);CHKERRQ(ierr);
907   ierr = PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q);CHKERRQ(ierr);
908   ierr = PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q_tmp);CHKERRQ(ierr);
909   ierr = PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L);CHKERRQ(ierr);
910   ierr = PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_tmp);CHKERRQ(ierr);
911   ierr = PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_save);CHKERRQ(ierr);
912   ierr = PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->N);CHKERRQ(ierr);
913   ierr = PetscMalloc1(mfqP->npmax * (mfqP->n+1) ,&mfqP->M);CHKERRQ(ierr);
914   ierr = PetscMalloc1(mfqP->npmax * (mfqP->npmax - mfqP->n - 1) , &mfqP->Z);CHKERRQ(ierr);
915   ierr = PetscMalloc1(mfqP->npmax,&mfqP->tau);CHKERRQ(ierr);
916   ierr = PetscMalloc1(mfqP->npmax,&mfqP->tau_tmp);CHKERRQ(ierr);
917   mfqP->nmax = PetscMax(5*mfqP->npmax,mfqP->n*(mfqP->n+1)/2);
918   ierr = PetscMalloc1(mfqP->nmax,&mfqP->npmaxwork);CHKERRQ(ierr);
919   ierr = PetscMalloc1(mfqP->nmax,&mfqP->npmaxiwork);CHKERRQ(ierr);
920   ierr = PetscMalloc1(mfqP->n,&mfqP->xmin);CHKERRQ(ierr);
921   ierr = PetscMalloc1(mfqP->m,&mfqP->C);CHKERRQ(ierr);
922   ierr = PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Fdiff);CHKERRQ(ierr);
923   ierr = PetscMalloc1(mfqP->npmax*mfqP->n,&mfqP->Disp);CHKERRQ(ierr);
924   ierr = PetscMalloc1(mfqP->n,&mfqP->Gres);CHKERRQ(ierr);
925   ierr = PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Hres);CHKERRQ(ierr);
926   ierr = PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Gpoints);CHKERRQ(ierr);
927   ierr = PetscMalloc1(mfqP->npmax,&mfqP->model_indices);CHKERRQ(ierr);
928   ierr = PetscMalloc1(mfqP->n,&mfqP->Xsubproblem);CHKERRQ(ierr);
929   ierr = PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Gdel);CHKERRQ(ierr);
930   ierr = PetscMalloc1(mfqP->n*mfqP->n*mfqP->m, &mfqP->Hdel);CHKERRQ(ierr);
931   ierr = PetscMalloc1(PetscMax(mfqP->m,mfqP->n),&mfqP->indices);CHKERRQ(ierr);
932   ierr = PetscMalloc1(mfqP->n,&mfqP->iwork);CHKERRQ(ierr);
933   for (i=0;i<PetscMax(mfqP->m,mfqP->n);i++) {
934     mfqP->indices[i] = i;
935   }
936   ierr = MPI_Comm_size(((PetscObject)tao)->comm,&mfqP->size);CHKERRQ(ierr);
937   if (mfqP->size > 1) {
938     VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localx);CHKERRQ(ierr);
939     VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localxmin);CHKERRQ(ierr);
940     VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localf);CHKERRQ(ierr);
941     VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localfmin);CHKERRQ(ierr);
942     ierr = ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxloc);CHKERRQ(ierr);
943     ierr = ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxglob);CHKERRQ(ierr);
944     ierr = ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfloc);CHKERRQ(ierr);
945     ierr = ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfglob);CHKERRQ(ierr);
946 
947 
948     ierr = VecScatterCreate(tao->solution,isxglob,mfqP->localx,isxloc,&mfqP->scatterx);CHKERRQ(ierr);
949     ierr = VecScatterCreate(tao->sep_objective,isfglob,mfqP->localf,isfloc,&mfqP->scatterf);CHKERRQ(ierr);
950 
951     ierr = ISDestroy(&isxloc);CHKERRQ(ierr);
952     ierr = ISDestroy(&isxglob);CHKERRQ(ierr);
953     ierr = ISDestroy(&isfloc);CHKERRQ(ierr);
954     ierr = ISDestroy(&isfglob);CHKERRQ(ierr);
955   }
956 
957   if (!mfqP->usegqt) {
958     KSP       ksp;
959     PC        pc;
960     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,mfqP->n,mfqP->n,mfqP->Xsubproblem,&mfqP->subx);CHKERRQ(ierr);
961     ierr = VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->subxl);CHKERRQ(ierr);
962     ierr = VecDuplicate(mfqP->subxl,&mfqP->subb);CHKERRQ(ierr);
963     ierr = VecDuplicate(mfqP->subxl,&mfqP->subxu);CHKERRQ(ierr);
964     ierr = VecDuplicate(mfqP->subxl,&mfqP->subpdel);CHKERRQ(ierr);
965     ierr = VecDuplicate(mfqP->subxl,&mfqP->subndel);CHKERRQ(ierr);
966     ierr = TaoCreate(PETSC_COMM_SELF,&mfqP->subtao);CHKERRQ(ierr);
967     /* ierr = TaoSetType(mfqP->subtao,"tao_bqpip");CHKERRQ(ierr); */
968     ierr = TaoSetType(mfqP->subtao,"tao_tron");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