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