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