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