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