Actual source code: pounders.c
petsc-3.11.4 2019-09-28
1: #include <../src/tao/leastsquares/impls/pounders/pounders.h>
3: static PetscErrorCode pounders_h(Tao subtao, Vec v, Mat H, Mat Hpre, void *ctx)
4: {
6: return(0);
7: }
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;
16: /* g = A*x (add b later)*/
17: MatMult(mfqP->subH,x,g);
19: /* f = 1/2 * x'*(Ax) + b'*x */
20: VecDot(x,g,&d1);
21: VecDot(mfqP->subb,x,&d2);
22: *f = 0.5 *d1 + d2;
24: /* now g = g + b */
25: VecAXPY(g, 1.0, mfqP->subb);
26: return(0);
27: }
29: static PetscErrorCode pounders_feval(Tao tao, Vec x, Vec F, PetscReal *fsum)
30: {
32: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
33: PetscInt i,row,col;
34: PetscReal fr,fc;
36: TaoComputeResidual(tao,x,F);
37: if (tao->res_weights_v) {
38: VecPointwiseMult(mfqP->workfvec,tao->res_weights_v,F);
39: VecDot(mfqP->workfvec,mfqP->workfvec,fsum);
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: VecGetValues(F,1,&row,&fr);
46: VecGetValues(F,1,&col,&fc);
47: *fsum += tao->res_weights_w[i]*fc*fr;
48: }
49: } else {
50: VecDot(F,F,fsum);
51: }
52: PetscInfo1(tao,"Least-squares residual norm: %20.19e\n",(double)*fsum);
53: if (PetscIsInfOrNanReal(*fsum)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
54: return(0);
55: }
57: PetscErrorCode gqtwrap(Tao tao,PetscReal *gnorm, PetscReal *qmin)
58: {
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;
69: if (!mfqP->usegqt) {
70: PetscReal maxval;
71: PetscInt i,j;
73: VecSetValues(mfqP->subb,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);
74: VecAssemblyBegin(mfqP->subb);
75: VecAssemblyEnd(mfqP->subb);
77: VecSet(mfqP->subx,0.0);
79: VecSet(mfqP->subndel,-1.0);
80: VecSet(mfqP->subpdel,+1.0);
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: MatSetValues(mfqP->subH,mfqP->n,mfqP->indices,mfqP->n,mfqP->indices,mfqP->Hres,INSERT_VALUES);
89: MatAssemblyBegin(mfqP->subH,MAT_FINAL_ASSEMBLY);
90: MatAssemblyEnd(mfqP->subH,MAT_FINAL_ASSEMBLY);
92: TaoResetStatistics(mfqP->subtao);
93: /* TaoSetTolerances(mfqP->subtao,*gnorm,*gnorm,PETSC_DEFAULT); */
94: /* enforce bound constraints -- experimental */
95: if (tao->XU && tao->XL) {
96: VecCopy(tao->XU,mfqP->subxu);
97: VecAXPY(mfqP->subxu,-1.0,tao->solution);
98: VecScale(mfqP->subxu,1.0/mfqP->delta);
99: VecCopy(tao->XL,mfqP->subxl);
100: VecAXPY(mfqP->subxl,-1.0,tao->solution);
101: VecScale(mfqP->subxl,1.0/mfqP->delta);
103: VecPointwiseMin(mfqP->subxu,mfqP->subxu,mfqP->subpdel);
104: VecPointwiseMax(mfqP->subxl,mfqP->subxl,mfqP->subndel);
105: } else {
106: VecCopy(mfqP->subpdel,mfqP->subxu);
107: VecCopy(mfqP->subndel,mfqP->subxl);
108: }
109: /* Make sure xu > xl */
110: VecCopy(mfqP->subxl,mfqP->subpdel);
111: VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);
112: VecMax(mfqP->subpdel,NULL,&maxval);
113: if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"upper bound < lower bound in subproblem");
114: /* Make sure xu > tao->solution > xl */
115: VecCopy(mfqP->subxl,mfqP->subpdel);
116: VecAXPY(mfqP->subpdel,-1.0,mfqP->subx);
117: VecMax(mfqP->subpdel,NULL,&maxval);
118: if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"initial guess < lower bound in subproblem");
120: VecCopy(mfqP->subx,mfqP->subpdel);
121: VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);
122: VecMax(mfqP->subpdel,NULL,&maxval);
123: if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"initial guess > upper bound in subproblem");
125: TaoSolve(mfqP->subtao);
126: TaoGetSolutionStatus(mfqP->subtao,NULL,qmin,NULL,NULL,NULL,NULL);
128: /* test bounds post-solution*/
129: VecCopy(mfqP->subxl,mfqP->subpdel);
130: VecAXPY(mfqP->subpdel,-1.0,mfqP->subx);
131: VecMax(mfqP->subpdel,NULL,&maxval);
132: if (maxval > 1e-5) {
133: PetscInfo(tao,"subproblem solution < lower bound\n");
134: tao->reason = TAO_DIVERGED_TR_REDUCTION;
135: }
137: VecCopy(mfqP->subx,mfqP->subpdel);
138: VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);
139: VecMax(mfqP->subpdel,NULL,&maxval);
140: if (maxval > 1e-5) {
141: PetscInfo(tao,"subproblem solution > upper bound\n");
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: return(0);
149: }
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;
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: }
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: VecGetValues(tao->res_weights_v,1,&i,&factor);
172: factor=factor*mfqP->C[i];
173: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasn,&factor,&mfqP->Fdiff[blasn*i],&ione,mfqP->Gres,&ione));
174: }
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: VecGetValues(tao->res_weights_v,1,&i,&wii);
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];
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: }
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];
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));
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: PetscInfo(tao,"Identity weights\n");
227: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->C,&ione,&zero,mfqP->Gres,&ione));
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));
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: return(0);
242: }
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);
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: return(0);
261: }
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 */
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);
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: }
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);
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));
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: }
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);
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: }
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);
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: }
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: return(0);
344: }
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;
360: /* Initialize M,N */
361: for (i=0;i<mfqP->n+1;i++) {
362: VecGetArrayRead(mfqP->Xhist[mfqP->model_indices[i]],&x);
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: VecRestoreArrayRead(mfqP->Xhist[mfqP->model_indices[i]],&x);
368: phi2eval(&mfqP->M[1+((mfqP->n+1)*i)],mfqP->n,&mfqP->N[mfqP->n*(mfqP->n+1)/2 * i]);
369: }
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: }
384: /* Reject if norm(d) >c2 */
385: if (!reject) {
386: VecCopy(mfqP->Xhist[point],mfqP->workxvec);
387: VecAXPY(mfqP->workxvec,-1.0,mfqP->Xhist[mfqP->minindex]);
388: VecNorm(mfqP->workxvec,NORM_2,&normd);
389: normd /= mfqP->delta;
390: if (normd > mfqP->c2) {
391: reject =1;
392: }
393: }
394: if (reject){
395: point--;
396: continue;
397: }
399: VecGetArrayRead(mfqP->Xhist[point],&x);
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: VecRestoreArrayRead(mfqP->Xhist[point],&x);
405: phi2eval(&mfqP->M[1+(mfqP->n+1)*mfqP->nmodelpoints],mfqP->n,&mfqP->N[mfqP->n*(mfqP->n+1)/2 * (mfqP->nmodelpoints)]);
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);
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 */
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);
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: }
438: /* Get svd for L_tmp(:,n+2:np+1) (L_tmp is modified in process) */
439: blasint = mfqP->nmodelpoints - mfqP->n;
440: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
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: PetscFPTrapPop();
443: if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine gesvd returned with value %d\n",info);
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;
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: }
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: }
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);
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: }
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: return(0);
496: }
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: {
504: /* Create new vector in history: X[newidx] = X[mfqP->index] + delta*X[index]*/
505: VecDuplicate(mfqP->Xhist[0],&mfqP->Xhist[mfqP->nHist]);
506: VecSetValues(mfqP->Xhist[mfqP->nHist],mfqP->n,mfqP->indices,&mfqP->Q_tmp[index*mfqP->npmax],INSERT_VALUES);
507: VecAssemblyBegin(mfqP->Xhist[mfqP->nHist]);
508: VecAssemblyEnd(mfqP->Xhist[mfqP->nHist]);
509: VecAYPX(mfqP->Xhist[mfqP->nHist],mfqP->delta,mfqP->Xhist[mfqP->minindex]);
511: /* Project into feasible region */
512: if (tao->XU && tao->XL) {
513: VecMedian(mfqP->Xhist[mfqP->nHist], tao->XL, tao->XU, mfqP->Xhist[mfqP->nHist]);
514: }
516: /* Compute value of new vector */
517: VecDuplicate(mfqP->Fhist[0],&mfqP->Fhist[mfqP->nHist]);
518: CHKMEMQ;
519: pounders_feval(tao,mfqP->Xhist[mfqP->nHist],mfqP->Fhist[mfqP->nHist],&mfqP->Fres[mfqP->nHist]);
521: /* Add new vector to model */
522: mfqP->model_indices[mfqP->nmodelpoints] = mfqP->nHist;
523: mfqP->nmodelpoints++;
524: mfqP->nHist++;
525: return(0);
526: }
528: PetscErrorCode modelimprove(Tao tao, TAO_POUNDERS *mfqP, PetscInt addallpoints)
529: {
530: /* modeld = Q(:,np+1:n)' */
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;
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: }
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));
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: addpoint(tao,mfqP,i);
570: }
571: }
572: if (!addallpoints) {
573: addpoint(tao,mfqP,minindex);
574: }
575: return(0);
576: }
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;
589: for (i=mfqP->nHist-1;i>=0;i--) {
590: VecGetArrayRead(mfqP->Xhist[i],&x);
591: for (j=0;j<mfqP->n;j++) {
592: mfqP->work[j] = (x[j] - xmin[j])/mfqP->delta;
593: }
594: VecRestoreArrayRead(mfqP->Xhist[i],&x);
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);
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: }
625: return(0);
626: }
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;
647: /* n = # of parameters
648: m = dimension (components) of function */
650: 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);
662: tao->niter=0;
663: if (tao->XL && tao->XU) {
664: /* Check x0 <= XU */
665: PetscReal val;
667: VecCopy(tao->solution,mfqP->Xhist[0]);
668: VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);
669: VecMax(mfqP->Xhist[0],NULL,&val);
670: if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 > upper bound");
672: /* Check x0 >= xl */
673: VecCopy(tao->XL,mfqP->Xhist[0]);
674: VecAXPY(mfqP->Xhist[0],-1.0,tao->solution);
675: VecMax(mfqP->Xhist[0],NULL,&val);
676: if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 < lower bound");
678: /* Check x0 + delta < XU -- should be able to get around this eventually */
680: VecSet(mfqP->Xhist[0],mfqP->delta);
681: VecAXPY(mfqP->Xhist[0],1.0,tao->solution);
682: VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);
683: VecMax(mfqP->Xhist[0],NULL,&val);
684: if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 + delta > upper bound");
685: }
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;
690: VecCopy(tao->solution,mfqP->Xhist[0]);
692: /* This provides enough information to approximate the gradient of the objective */
693: /* using a forward difference scheme. */
695: PetscInfo1(tao,"Initialize simplex; delta = %10.9e\n",(double)mfqP->delta);
696: pounders_feval(tao,mfqP->Xhist[0],mfqP->Fhist[0],&mfqP->Fres[0]);
697: mfqP->minindex = 0;
698: minnorm = mfqP->Fres[0];
700: VecGetOwnershipRange(mfqP->Xhist[0],&low,&high);
701: for (i=1;i<mfqP->n+1;++i) {
702: VecCopy(mfqP->Xhist[0],mfqP->Xhist[i]);
704: if (i-1 >= low && i-1 < high) {
705: VecGetArray(mfqP->Xhist[i],&x);
706: x[i-1-low] += mfqP->delta;
707: VecRestoreArray(mfqP->Xhist[i],&x);
708: }
709: CHKMEMQ;
710: pounders_feval(tao,mfqP->Xhist[i],mfqP->Fhist[i],&mfqP->Fres[i]);
711: if (mfqP->Fres[i] < minnorm) {
712: mfqP->minindex = i;
713: minnorm = mfqP->Fres[i];
714: }
715: }
716: VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);
717: VecCopy(mfqP->Fhist[mfqP->minindex],tao->ls_res);
718: PetscInfo1(tao,"Finalize simplex; minnorm = %10.9e\n",(double)minnorm);
720: /* Gather mpi vecs to one big local vec */
722: /* Begin serial code */
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;
729: PetscInfo1(tao,"Build matrix: %D\n",(PetscInt)mfqP->size);
730: if (1 == mfqP->size) {
731: VecGetArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);
732: for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i];
733: VecRestoreArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);
734: VecGetArrayRead(mfqP->Fhist[mfqP->minindex],&fmin);
735: for (i=0;i<mfqP->n+1;i++) {
736: if (i == mfqP->minindex) continue;
738: VecGetArray(mfqP->Xhist[i],&x);
739: for (j=0;j<mfqP->n;j++) {
740: mfqP->Disp[ii+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta;
741: }
742: VecRestoreArray(mfqP->Xhist[i],&x);
744: VecGetArray(mfqP->Fhist[i],&f);
745: for (j=0;j<mfqP->m;j++) {
746: mfqP->Fdiff[ii+mfqP->n*j] = f[j] - fmin[j];
747: }
748: VecRestoreArray(mfqP->Fhist[i],&f);
750: mfqP->model_indices[ii++] = i;
751: }
752: for (j=0;j<mfqP->m;j++) {
753: mfqP->C[j] = fmin[j];
754: }
755: VecRestoreArrayRead(mfqP->Fhist[mfqP->minindex],&fmin);
756: } else {
757: VecSet(mfqP->localxmin,0);
758: VecScatterBegin(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);
759: VecScatterEnd(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);
761: VecGetArrayRead(mfqP->localxmin,&xmint);
762: for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i];
763: VecRestoreArrayRead(mfqP->localxmin,&xmint);
765: VecScatterBegin(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);
766: VecScatterEnd(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);
767: VecGetArrayRead(mfqP->localfmin,&fmin);
768: for (i=0;i<mfqP->n+1;i++) {
769: if (i == mfqP->minindex) continue;
771: VecScatterBegin(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);
772: VecScatterEnd(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);
773: VecGetArray(mfqP->localx,&x);
774: for (j=0;j<mfqP->n;j++) {
775: mfqP->Disp[ii+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta;
776: }
777: VecRestoreArray(mfqP->localx,&x);
779: VecScatterBegin(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);
780: VecScatterEnd(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);
781: VecGetArray(mfqP->localf,&f);
782: for (j=0;j<mfqP->m;j++) {
783: mfqP->Fdiff[ii+mfqP->n*j] = f[j] - fmin[j];
784: }
785: VecRestoreArray(mfqP->localf,&f);
787: mfqP->model_indices[ii++] = i;
788: }
789: for (j=0;j<mfqP->m;j++) {
790: mfqP->C[j] = fmin[j];
791: }
792: VecRestoreArrayRead(mfqP->localfmin,&fmin);
793: }
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: PetscInfo1(tao,"Linear solve return: %D\n",(PetscInt)info);
802: cres = minnorm;
803: pounders_update_res(tao);
805: valid = PETSC_TRUE;
807: VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);
808: VecAssemblyBegin(tao->gradient);
809: VecAssemblyEnd(tao->gradient);
810: VecNorm(tao->gradient,NORM_2,&gnorm);
811: gnorm *= mfqP->delta;
812: VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);
813:
814: tao->reason = TAO_CONTINUE_ITERATING;
815: TaoLogConvergenceHistory(tao,minnorm,gnorm,0.0,tao->ksp_its);
816: TaoMonitor(tao,tao->niter,minnorm,gnorm,0.0,step);
817: (*tao->ops->convergencetest)(tao,tao->cnvP);
818:
819: mfqP->nHist = mfqP->n+1;
820: mfqP->nmodelpoints = mfqP->n+1;
821: PetscInfo1(tao,"Initial gradient: %20.19e\n",(double)gnorm);
823: while (tao->reason == TAO_CONTINUE_ITERATING) {
824: PetscReal gnm = 1e-4;
825: /* Call general purpose update function */
826: if (tao->ops->update) {
827: (*tao->ops->update)(tao, tao->niter, tao->user_update);
828: }
829: tao->niter++;
830: /* Solve the subproblem min{Q(s): ||s|| <= 1.0} */
831: gqtwrap(tao,&gnm,&mdec);
832: /* Evaluate the function at the new point */
834: for (i=0;i<mfqP->n;i++) {
835: mfqP->work[i] = mfqP->Xsubproblem[i]*mfqP->delta + mfqP->xmin[i];
836: }
837: VecDuplicate(tao->solution,&mfqP->Xhist[mfqP->nHist]);
838: VecDuplicate(tao->ls_res,&mfqP->Fhist[mfqP->nHist]);
839: VecSetValues(mfqP->Xhist[mfqP->nHist],mfqP->n,mfqP->indices,mfqP->work,INSERT_VALUES);
840: VecAssemblyBegin(mfqP->Xhist[mfqP->nHist]);
841: VecAssemblyEnd(mfqP->Xhist[mfqP->nHist]);
843: pounders_feval(tao,mfqP->Xhist[mfqP->nHist],mfqP->Fhist[mfqP->nHist],&mfqP->Fres[mfqP->nHist]);
845: rho = (mfqP->Fres[mfqP->minindex] - mfqP->Fres[mfqP->nHist]) / mdec;
846: mfqP->nHist++;
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'; */
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: VecCopy(mfqP->Fhist[mfqP->minindex],tao->ls_res);
879: /* Change current center */
880: VecGetArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);
881: for (i=0;i<mfqP->n;i++) {
882: mfqP->xmin[i] = xmint[i];
883: }
884: VecRestoreArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);
885: }
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: affpoints(mfqP,mfqP->xmin,mfqP->c1);
892: if (mfqP->nmodelpoints < mfqP->n) {
893: PetscInfo(tao,"Model not valid -- model-improving\n");
894: modelimprove(tao,mfqP,1);
895: }
896: }
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: }
911: /* Compute the next interpolation set */
912: mfqP->q_is_I = 1;
913: mfqP->nmodelpoints=0;
914: PetscInfo2(tao,"Affine Points: xmin = %20.19e, c1 = %20.19e\n",(double)*mfqP->xmin,(double)mfqP->c1);
915: affpoints(mfqP,mfqP->xmin,mfqP->c1);
916: if (mfqP->nmodelpoints == mfqP->n) {
917: valid = PETSC_TRUE;
918: } else {
919: valid = PETSC_FALSE;
920: PetscInfo2(tao,"Affine Points: xmin = %20.19e, c2 = %20.19e\n",(double)*mfqP->xmin,(double)mfqP->c2);
921: affpoints(mfqP,mfqP->xmin,mfqP->c2);
922: if (mfqP->n > mfqP->nmodelpoints) {
923: PetscInfo(tao,"Model not valid -- adding geometry points\n");
924: modelimprove(tao,mfqP,mfqP->n - mfqP->nmodelpoints);
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: morepoints(mfqP);
933: for (i=0;i<mfqP->nmodelpoints;i++) {
934: VecGetArray(mfqP->Xhist[mfqP->model_indices[i]],&x);
935: for (j=0;j<mfqP->n;j++) {
936: mfqP->Disp[i + mfqP->npmax*j] = (x[j] - mfqP->xmin[j]) / deltaold;
937: }
938: VecRestoreArray(mfqP->Xhist[mfqP->model_indices[i]],&x);
939: VecGetArray(mfqP->Fhist[mfqP->model_indices[i]],&f);
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: VecRestoreArray(mfqP->Fhist[mfqP->model_indices[i]],&f);
950: }
952: /* Update the quadratic model */
953: PetscInfo2(tao,"Get Quad, size: %D, points: %D\n",mfqP->n,mfqP->nmodelpoints);
954: getquadpounders(mfqP);
955: VecGetArrayRead(mfqP->Fhist[mfqP->minindex],&fmin);
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));
968: /* Get residuals */
969: cres = mfqP->Fres[mfqP->minindex];
970: pounders_update_res(tao);
972: /* Export solution and gradient residual to TAO */
973: VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);
974: VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);
975: VecAssemblyBegin(tao->gradient);
976: VecAssemblyEnd(tao->gradient);
977: VecNorm(tao->gradient,NORM_2,&gnorm);
978: gnorm *= mfqP->delta;
979: /* final criticality test */
980: TaoLogConvergenceHistory(tao,minnorm,gnorm,0.0,tao->ksp_its);
981: TaoMonitor(tao,tao->niter,minnorm,gnorm,0.0,step);
982: (*tao->ops->convergencetest)(tao,tao->cnvP);
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: PetscInfo(tao,"Identical model used in successive iterations\n");
1002: tao->reason = TAO_CONVERGED_STEPTOL;
1003: }
1004: }
1005: return(0);
1006: }
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;
1016: if (!tao->gradient) {VecDuplicate(tao->solution,&tao->gradient); }
1017: if (!tao->stepdirection) {VecDuplicate(tao->solution,&tao->stepdirection); }
1018: VecGetSize(tao->solution,&mfqP->n);
1019: VecGetSize(tao->ls_res,&mfqP->m);
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);
1027: PetscMalloc1(tao->max_funcs+100,&mfqP->Xhist);
1028: PetscMalloc1(tao->max_funcs+100,&mfqP->Fhist);
1029: for (i=0;i<mfqP->n+1;i++) {
1030: VecDuplicate(tao->solution,&mfqP->Xhist[i]);
1031: VecDuplicate(tao->ls_res,&mfqP->Fhist[i]);
1032: }
1033: VecDuplicate(tao->solution,&mfqP->workxvec);
1034: VecDuplicate(tao->ls_res,&mfqP->workfvec);
1035: mfqP->nHist = 0;
1037: PetscMalloc1(tao->max_funcs+100,&mfqP->Fres);
1038: PetscMalloc1(mfqP->npmax*mfqP->m,&mfqP->RES);
1039: PetscMalloc1(mfqP->n,&mfqP->work);
1040: PetscMalloc1(mfqP->n,&mfqP->work2);
1041: PetscMalloc1(mfqP->n,&mfqP->work3);
1042: PetscMalloc1(PetscMax(mfqP->m,mfqP->n+1),&mfqP->mwork);
1043: PetscMalloc1(mfqP->npmax - mfqP->n - 1,&mfqP->omega);
1044: PetscMalloc1(mfqP->n * (mfqP->n+1) / 2,&mfqP->beta);
1045: PetscMalloc1(mfqP->n + 1 ,&mfqP->alpha);
1047: PetscMalloc1(mfqP->n*mfqP->n*mfqP->m,&mfqP->H);
1048: PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q);
1049: PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q_tmp);
1050: PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L);
1051: PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_tmp);
1052: PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_save);
1053: PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->N);
1054: PetscMalloc1(mfqP->npmax * (mfqP->n+1) ,&mfqP->M);
1055: PetscMalloc1(mfqP->npmax * (mfqP->npmax - mfqP->n - 1) , &mfqP->Z);
1056: PetscMalloc1(mfqP->npmax,&mfqP->tau);
1057: PetscMalloc1(mfqP->npmax,&mfqP->tau_tmp);
1058: mfqP->nmax = PetscMax(5*mfqP->npmax,mfqP->n*(mfqP->n+1)/2);
1059: PetscMalloc1(mfqP->nmax,&mfqP->npmaxwork);
1060: PetscMalloc1(mfqP->nmax,&mfqP->npmaxiwork);
1061: PetscMalloc1(mfqP->n,&mfqP->xmin);
1062: PetscMalloc1(mfqP->m,&mfqP->C);
1063: PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Fdiff);
1064: PetscMalloc1(mfqP->npmax*mfqP->n,&mfqP->Disp);
1065: PetscMalloc1(mfqP->n,&mfqP->Gres);
1066: PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Hres);
1067: PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Gpoints);
1068: PetscMalloc1(mfqP->npmax,&mfqP->model_indices);
1069: PetscMalloc1(mfqP->npmax,&mfqP->last_model_indices);
1070: PetscMalloc1(mfqP->n,&mfqP->Xsubproblem);
1071: PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Gdel);
1072: PetscMalloc1(mfqP->n*mfqP->n*mfqP->m, &mfqP->Hdel);
1073: PetscMalloc1(PetscMax(mfqP->m,mfqP->n),&mfqP->indices);
1074: PetscMalloc1(mfqP->n,&mfqP->iwork);
1075: PetscMalloc1(mfqP->m*mfqP->m,&mfqP->w);
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: MPI_Comm_size(((PetscObject)tao)->comm,&mfqP->size);
1089: if (mfqP->size > 1) {
1090: VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localx);
1091: VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localxmin);
1092: VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localf);
1093: VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localfmin);
1094: ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxloc);
1095: ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxglob);
1096: ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfloc);
1097: ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfglob);
1100: VecScatterCreate(tao->solution,isxglob,mfqP->localx,isxloc,&mfqP->scatterx);
1101: VecScatterCreate(tao->ls_res,isfglob,mfqP->localf,isfloc,&mfqP->scatterf);
1103: ISDestroy(&isxloc);
1104: ISDestroy(&isxglob);
1105: ISDestroy(&isfloc);
1106: ISDestroy(&isfglob);
1107: }
1109: if (!mfqP->usegqt) {
1110: KSP ksp;
1111: PC pc;
1112: VecCreateSeqWithArray(PETSC_COMM_SELF,mfqP->n,mfqP->n,mfqP->Xsubproblem,&mfqP->subx);
1113: VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->subxl);
1114: VecDuplicate(mfqP->subxl,&mfqP->subb);
1115: VecDuplicate(mfqP->subxl,&mfqP->subxu);
1116: VecDuplicate(mfqP->subxl,&mfqP->subpdel);
1117: VecDuplicate(mfqP->subxl,&mfqP->subndel);
1118: TaoCreate(PETSC_COMM_SELF,&mfqP->subtao);
1119: PetscObjectIncrementTabLevel((PetscObject)mfqP->subtao, (PetscObject)tao, 1);
1120: TaoSetType(mfqP->subtao,TAOBNTR);
1121: TaoSetOptionsPrefix(mfqP->subtao,"pounders_subsolver_");
1122: TaoSetInitialVector(mfqP->subtao,mfqP->subx);
1123: TaoSetObjectiveAndGradientRoutine(mfqP->subtao,pounders_fg,(void*)mfqP);
1124: TaoSetMaximumIterations(mfqP->subtao,mfqP->gqt_maxits);
1125: TaoSetFromOptions(mfqP->subtao);
1126: TaoGetKSP(mfqP->subtao,&ksp);
1127: if (ksp) {
1128: KSPGetPC(ksp,&pc);
1129: PCSetType(pc,PCNONE);
1130: }
1131: TaoSetVariableBounds(mfqP->subtao,mfqP->subxl,mfqP->subxu);
1132: MatCreateSeqDense(PETSC_COMM_SELF,mfqP->n,mfqP->n,mfqP->Hres,&mfqP->subH);
1133: TaoSetHessianRoutine(mfqP->subtao,mfqP->subH,mfqP->subH,pounders_h,(void*)mfqP);
1134: }
1135: return(0);
1136: }
1138: static PetscErrorCode TaoDestroy_POUNDERS(Tao tao)
1139: {
1140: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
1141: PetscInt i;
1145: if (!mfqP->usegqt) {
1146: TaoDestroy(&mfqP->subtao);
1147: VecDestroy(&mfqP->subx);
1148: VecDestroy(&mfqP->subxl);
1149: VecDestroy(&mfqP->subxu);
1150: VecDestroy(&mfqP->subb);
1151: VecDestroy(&mfqP->subpdel);
1152: VecDestroy(&mfqP->subndel);
1153: MatDestroy(&mfqP->subH);
1154: }
1155: PetscFree(mfqP->Fres);
1156: PetscFree(mfqP->RES);
1157: PetscFree(mfqP->work);
1158: PetscFree(mfqP->work2);
1159: PetscFree(mfqP->work3);
1160: PetscFree(mfqP->mwork);
1161: PetscFree(mfqP->omega);
1162: PetscFree(mfqP->beta);
1163: PetscFree(mfqP->alpha);
1164: PetscFree(mfqP->H);
1165: PetscFree(mfqP->Q);
1166: PetscFree(mfqP->Q_tmp);
1167: PetscFree(mfqP->L);
1168: PetscFree(mfqP->L_tmp);
1169: PetscFree(mfqP->L_save);
1170: PetscFree(mfqP->N);
1171: PetscFree(mfqP->M);
1172: PetscFree(mfqP->Z);
1173: PetscFree(mfqP->tau);
1174: PetscFree(mfqP->tau_tmp);
1175: PetscFree(mfqP->npmaxwork);
1176: PetscFree(mfqP->npmaxiwork);
1177: PetscFree(mfqP->xmin);
1178: PetscFree(mfqP->C);
1179: PetscFree(mfqP->Fdiff);
1180: PetscFree(mfqP->Disp);
1181: PetscFree(mfqP->Gres);
1182: PetscFree(mfqP->Hres);
1183: PetscFree(mfqP->Gpoints);
1184: PetscFree(mfqP->model_indices);
1185: PetscFree(mfqP->last_model_indices);
1186: PetscFree(mfqP->Xsubproblem);
1187: PetscFree(mfqP->Gdel);
1188: PetscFree(mfqP->Hdel);
1189: PetscFree(mfqP->indices);
1190: PetscFree(mfqP->iwork);
1191: PetscFree(mfqP->w);
1192: for (i=0;i<mfqP->nHist;i++) {
1193: VecDestroy(&mfqP->Xhist[i]);
1194: VecDestroy(&mfqP->Fhist[i]);
1195: }
1196: VecDestroy(&mfqP->workxvec);
1197: VecDestroy(&mfqP->workfvec);
1198: PetscFree(mfqP->Xhist);
1199: PetscFree(mfqP->Fhist);
1201: if (mfqP->size > 1) {
1202: VecDestroy(&mfqP->localx);
1203: VecDestroy(&mfqP->localxmin);
1204: VecDestroy(&mfqP->localf);
1205: VecDestroy(&mfqP->localfmin);
1206: }
1207: PetscFree(tao->data);
1208: return(0);
1209: }
1211: static PetscErrorCode TaoSetFromOptions_POUNDERS(PetscOptionItems *PetscOptionsObject,Tao tao)
1212: {
1213: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
1217: PetscOptionsHead(PetscOptionsObject,"POUNDERS method for least-squares optimization");
1218: PetscOptionsReal("-tao_pounders_delta","initial delta","",mfqP->delta,&mfqP->delta0,NULL);
1219: mfqP->delta = mfqP->delta0;
1220: PetscOptionsInt("-tao_pounders_npmax","max number of points in model","",mfqP->npmax,&mfqP->npmax,NULL);
1221: PetscOptionsBool("-tao_pounders_gqt","use gqt algorithm for subproblem","",mfqP->usegqt,&mfqP->usegqt,NULL);
1222: PetscOptionsTail();
1223: return(0);
1224: }
1226: static PetscErrorCode TaoView_POUNDERS(Tao tao, PetscViewer viewer)
1227: {
1228: TAO_POUNDERS *mfqP = (TAO_POUNDERS *)tao->data;
1229: PetscBool isascii;
1233: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1234: if (isascii) {
1235: PetscViewerASCIIPrintf(viewer, "initial delta: %g\n",(double)mfqP->delta0);
1236: PetscViewerASCIIPrintf(viewer, "final delta: %g\n",(double)mfqP->delta);
1237: PetscViewerASCIIPrintf(viewer, "model points: %D\n",mfqP->nmodelpoints);
1238: if (mfqP->usegqt) {
1239: PetscViewerASCIIPrintf(viewer, "subproblem solver: gqt\n");
1240: } else {
1241: TaoView(mfqP->subtao, viewer);
1242: }
1243: }
1244: return(0);
1245: }
1246: /*MC
1247: TAOPOUNDERS - POUNDERS derivate-free model-based algorithm for nonlinear least squares
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
1254: Level: beginner
1256: M*/
1258: PETSC_EXTERN PetscErrorCode TaoCreate_POUNDERS(Tao tao)
1259: {
1260: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
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;
1270: PetscNewLog(tao,&mfqP);
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: return(0);
1292: }