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