Actual source code: pounders.c
petsc-3.8.4 2018-03-24
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: TaoComputeSeparableObjective(tao,x,F);
37: if (tao->sep_weights_v) {
38: VecPointwiseMult(mfqP->workfvec,tao->sep_weights_v,F);
39: VecNorm(mfqP->workfvec,NORM_2,fsum);
40: *fsum*=(*fsum);
41: } else if (tao->sep_weights_w) {
42: *fsum=0;
43: for (i=0;i<tao->sep_weights_n;i++) {
44: row=tao->sep_weights_rows[i];
45: col=tao->sep_weights_cols[i];
46: VecGetValues(F,1,&row,&fr);
47: VecGetValues(F,1,&col,&fc);
48: *fsum += tao->sep_weights_w[i]*fc*fr;
49: }
50: } else {
51: VecNorm(F,NORM_2,fsum);
52: *fsum*=(*fsum);
53: }
54: if (PetscIsInfOrNanReal(*fsum)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
55: return(0);
56: }
58: PetscErrorCode gqtwrap(Tao tao,PetscReal *gnorm, PetscReal *qmin)
59: {
61: #if defined(PETSC_USE_REAL_SINGLE)
62: PetscReal atol=1.0e-5;
63: #else
64: PetscReal atol=1.0e-10;
65: #endif
66: PetscInt info,its;
67: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
70: if (!mfqP->usegqt) {
71: PetscReal maxval;
72: PetscInt i,j;
74: VecSetValues(mfqP->subb,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);
75: VecAssemblyBegin(mfqP->subb);
76: VecAssemblyEnd(mfqP->subb);
78: VecSet(mfqP->subx,0.0);
80: VecSet(mfqP->subndel,-mfqP->delta);
81: VecSet(mfqP->subpdel,mfqP->delta);
83: for (i=0;i<mfqP->n;i++) {
84: for (j=i;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->sep_weights_v) {
169: /* Vector(diagonal) weights: gres = sum_i(wii*ci*gi) */
170: for (i=0;i<mfqP->m;i++) {
171: VecGetValues(tao->sep_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: }
175: } else if (tao->sep_weights_n) {
176: /* general case: .5 * Gres= sum_ij[wij * (cjgi + cigj)] */
177: for (i=0;i<tao->sep_weights_n;i++) {
178: row=tao->sep_weights_rows[i];
179: col=tao->sep_weights_cols[i];
181: factor = tao->sep_weights_w[i]*mfqP->C[col]/2.0;
182: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasn,&factor,&mfqP->Fdiff[blasn*row],&ione,mfqP->Gres,&ione));
183: factor = tao->sep_weights_w[i]*mfqP->C[row]/2.0;
184: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasn,&factor,&mfqP->Fdiff[blasn*col],&ione,mfqP->Gres,&ione));
185: }
186: } else {
187: /* default: Gres= sum_i[cigi] = G*c' */
188: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->C,&ione,&zero,mfqP->Gres,&ione));
189: }
191: /* compute Hres = sum_ij [wij * (*ci*Hj + cj*Hi + gi gj' + gj gi') ] */
192: if (tao->sep_weights_v) {
193: /* vector(diagonal weights) Hres = sum_i(wii*(ci*Hi + gi * gi' )*/
194: for (i=0;i<mfqP->m;i++) {
195: VecGetValues(tao->sep_weights_v,1,&i,&wii);
196: if (tao->niter>1) {
197: factor=wii*mfqP->C[i];
198: /* add wii * ci * Hi */
199: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasn2,&factor,&mfqP->H[i],&blasn2,mfqP->Hres,&ione));
200: }
201: /* add wii * gi * gi' */
202: PetscStackCallBLAS("BLASgemm_",BLASgemm_("N","T",&blasn,&blasn,&ione,&wii,&mfqP->Fdiff[blasn*i],&blasn,&mfqP->Fdiff[blasn*i],&blasn,&one,mfqP->Hres,&blasn));
203: }
204: } else if (tao->sep_weights_w) {
205: /* .5 * sum_ij [wij * (*ci*Hj + cj*Hi + gi gj' + gj gi') ] */
206: for (i=0;i<tao->sep_weights_n;i++) {
207: row=tao->sep_weights_rows[i];
208: col=tao->sep_weights_cols[i];
209: factor=tao->sep_weights_w[i]/2.0;
210: /* add wij * gi gj' + wij * gj gi' */
211: PetscStackCallBLAS("BLASgemm_",BLASgemm_("N","T",&blasn,&blasn,&ione,&factor,&mfqP->Fdiff[blasn*row],&blasn,&mfqP->Fdiff[blasn*col],&blasn,&one,mfqP->Hres,&blasn));
212: PetscStackCallBLAS("BLASgemm_",BLASgemm_("N","T",&blasn,&blasn,&ione,&factor,&mfqP->Fdiff[blasn*col],&blasn,&mfqP->Fdiff[blasn*row],&blasn,&one,mfqP->Hres,&blasn));
213: }
214: if (tao->niter > 1) {
215: for (i=0;i<tao->sep_weights_n;i++) {
216: row=tao->sep_weights_rows[i];
217: col=tao->sep_weights_cols[i];
219: /* add wij*cj*Hi */
220: factor = tao->sep_weights_w[i]*mfqP->C[col]/2.0;
221: PetscStackCallBLAS("BLASaxpy_",BLASaxpy_(&blasn2,&factor,&mfqP->H[row],&blasn2,mfqP->Hres,&ione));
223: /* add wij*ci*Hj */
224: factor = tao->sep_weights_w[i]*mfqP->C[row]/2.0;
225: PetscStackCallBLAS("BLASaxpy_",BLASaxpy_(&blasn2,&factor,&mfqP->H[col],&blasn2,mfqP->Hres,&ione));
226: }
227: }
228: } else {
229: /* Hres = G*G' + 0.5 sum {F(xkin,i)*H(:,:,i)} */
230: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&blasn,&blasn,&blasm,&one,mfqP->Fdiff, &blasn,mfqP->Fdiff, &blasn,&zero,mfqP->Hres,&blasn));
232: /* sum(F(xkin,i)*H(:,:,i)) */
233: if (tao->niter>1) {
234: for (i=0;i<mfqP->m;i++) {
235: factor = mfqP->C[i];
236: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasn2,&factor,&mfqP->H[i],&blasn2,mfqP->Hres,&ione));
237: }
238: }
239: }
240: return(0);
241: }
243: PetscErrorCode phi2eval(PetscReal *x, PetscInt n, PetscReal *phi)
244: {
245: /* 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] */
246: PetscInt i,j,k;
247: PetscReal sqrt2 = PetscSqrtReal(2.0);
250: j=0;
251: for (i=0;i<n;i++) {
252: phi[j] = 0.5 * x[i]*x[i];
253: j++;
254: for (k=i+1;k<n;k++) {
255: phi[j] = x[i]*x[k]/sqrt2;
256: j++;
257: }
258: }
259: return(0);
260: }
262: PetscErrorCode getquadpounders(TAO_POUNDERS *mfqP)
263: {
264: /* Computes the parameters of the quadratic Q(x) = c + g'*x + 0.5*x*G*x'
265: that satisfies the interpolation conditions Q(X[:,j]) = f(j)
266: for j=1,...,m and with a Hessian matrix of least Frobenius norm */
268: /* NB --we are ignoring c */
269: PetscInt i,j,k,num,np = mfqP->nmodelpoints;
270: PetscReal one = 1.0,zero=0.0,negone=-1.0;
271: PetscBLASInt blasnpmax = mfqP->npmax;
272: PetscBLASInt blasnplus1 = mfqP->n+1;
273: PetscBLASInt blasnp = np;
274: PetscBLASInt blasint = mfqP->n*(mfqP->n+1) / 2;
275: PetscBLASInt blasint2 = np - mfqP->n-1;
276: PetscBLASInt info,ione=1;
277: PetscReal sqrt2 = PetscSqrtReal(2.0);
280: for (i=0;i<mfqP->n*mfqP->m;i++) {
281: mfqP->Gdel[i] = 0;
282: }
283: for (i=0;i<mfqP->n*mfqP->n*mfqP->m;i++) {
284: mfqP->Hdel[i] = 0;
285: }
287: /* factor M */
288: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&blasnplus1,&blasnp,mfqP->M,&blasnplus1,mfqP->npmaxiwork,&info));
289: if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine getrf returned with value %d\n",info);
291: if (np == mfqP->n+1) {
292: for (i=0;i<mfqP->npmax-mfqP->n-1;i++) {
293: mfqP->omega[i]=0.0;
294: }
295: for (i=0;i<mfqP->n*(mfqP->n+1)/2;i++) {
296: mfqP->beta[i]=0.0;
297: }
298: } else {
299: /* Let Ltmp = (L'*L) */
300: PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&blasint2,&blasint2,&blasint,&one,&mfqP->L[(mfqP->n+1)*blasint],&blasint,&mfqP->L[(mfqP->n+1)*blasint],&blasint,&zero,mfqP->L_tmp,&blasint));
302: /* factor Ltmp */
303: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&blasint2,mfqP->L_tmp,&blasint,&info));
304: if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine potrf returned with value %d\n",info);
305: }
307: for (k=0;k<mfqP->m;k++) {
308: if (np != mfqP->n+1) {
309: /* Solve L'*L*Omega = Z' * RESk*/
310: PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&blasnp,&blasint2,&one,mfqP->Z,&blasnpmax,&mfqP->RES[mfqP->npmax*k],&ione,&zero,mfqP->omega,&ione));
311: PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&blasint2,&ione,mfqP->L_tmp,&blasint,mfqP->omega,&blasint2,&info));
312: if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine potrs returned with value %d\n",info);
314: /* Beta = L*Omega */
315: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasint,&blasint2,&one,&mfqP->L[(mfqP->n+1)*blasint],&blasint,mfqP->omega,&ione,&zero,mfqP->beta,&ione));
316: }
318: /* solve M'*Alpha = RESk - N'*Beta */
319: PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&blasint,&blasnp,&negone,mfqP->N,&blasint,mfqP->beta,&ione,&one,&mfqP->RES[mfqP->npmax*k],&ione));
320: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&blasnplus1,&ione,mfqP->M,&blasnplus1,mfqP->npmaxiwork,&mfqP->RES[mfqP->npmax*k],&blasnplus1,&info));
321: if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine getrs returned with value %d\n",info);
323: /* Gdel(:,k) = Alpha(2:n+1) */
324: for (i=0;i<mfqP->n;i++) {
325: mfqP->Gdel[i + mfqP->n*k] = mfqP->RES[mfqP->npmax*k + i+1];
326: }
328: /* Set Hdels */
329: num=0;
330: for (i=0;i<mfqP->n;i++) {
331: /* H[i,i,k] = Beta(num) */
332: mfqP->Hdel[(i*mfqP->n + i)*mfqP->m + k] = mfqP->beta[num];
333: num++;
334: for (j=i+1;j<mfqP->n;j++) {
335: /* H[i,j,k] = H[j,i,k] = Beta(num)/sqrt(2) */
336: mfqP->Hdel[(j*mfqP->n + i)*mfqP->m + k] = mfqP->beta[num]/sqrt2;
337: mfqP->Hdel[(i*mfqP->n + j)*mfqP->m + k] = mfqP->beta[num]/sqrt2;
338: num++;
339: }
340: }
341: }
342: return(0);
343: }
345: PetscErrorCode morepoints(TAO_POUNDERS *mfqP)
346: {
347: /* Assumes mfqP->model_indices[0] is minimum index
348: Finishes adding points to mfqP->model_indices (up to npmax)
349: Computes L,Z,M,N
350: np is actual number of points in model (should equal npmax?) */
351: PetscInt point,i,j,offset;
352: PetscInt reject;
353: PetscBLASInt blasn=mfqP->n,blasnpmax=mfqP->npmax,blasnplus1=mfqP->n+1,info,blasnmax=mfqP->nmax,blasint,blasint2,blasnp,blasmaxmn;
354: const PetscReal *x;
355: PetscReal normd;
356: PetscErrorCode ierr;
359: /* Initialize M,N */
360: for (i=0;i<mfqP->n+1;i++) {
361: VecGetArrayRead(mfqP->Xhist[mfqP->model_indices[i]],&x);
362: mfqP->M[(mfqP->n+1)*i] = 1.0;
363: for (j=0;j<mfqP->n;j++) {
364: mfqP->M[j+1+((mfqP->n+1)*i)] = (x[j] - mfqP->xmin[j]) / mfqP->delta;
365: }
366: VecRestoreArrayRead(mfqP->Xhist[mfqP->model_indices[i]],&x);
367: phi2eval(&mfqP->M[1+((mfqP->n+1)*i)],mfqP->n,&mfqP->N[mfqP->n*(mfqP->n+1)/2 * i]);
368: }
370: /* Now we add points until we have npmax starting with the most recent ones */
371: point = mfqP->nHist-1;
372: mfqP->nmodelpoints = mfqP->n+1;
373: while (mfqP->nmodelpoints < mfqP->npmax && point>=0) {
374: /* Reject any points already in the model */
375: reject = 0;
376: for (j=0;j<mfqP->n+1;j++) {
377: if (point == mfqP->model_indices[j]) {
378: reject = 1;
379: break;
380: }
381: }
383: /* Reject if norm(d) >c2 */
384: if (!reject) {
385: VecCopy(mfqP->Xhist[point],mfqP->workxvec);
386: VecAXPY(mfqP->workxvec,-1.0,mfqP->Xhist[mfqP->minindex]);
387: VecNorm(mfqP->workxvec,NORM_2,&normd);
388: normd /= mfqP->delta;
389: if (normd > mfqP->c2) {
390: reject =1;
391: }
392: }
393: if (reject){
394: point--;
395: continue;
396: }
398: VecGetArrayRead(mfqP->Xhist[point],&x);
399: mfqP->M[(mfqP->n+1)*mfqP->nmodelpoints] = 1.0;
400: for (j=0;j<mfqP->n;j++) {
401: mfqP->M[j+1+((mfqP->n+1)*mfqP->nmodelpoints)] = (x[j] - mfqP->xmin[j]) / mfqP->delta;
402: }
403: VecRestoreArrayRead(mfqP->Xhist[point],&x);
404: phi2eval(&mfqP->M[1+(mfqP->n+1)*mfqP->nmodelpoints],mfqP->n,&mfqP->N[mfqP->n*(mfqP->n+1)/2 * (mfqP->nmodelpoints)]);
406: /* Update QR factorization */
407: /* Copy M' to Q_tmp */
408: for (i=0;i<mfqP->n+1;i++) {
409: for (j=0;j<mfqP->npmax;j++) {
410: mfqP->Q_tmp[j+mfqP->npmax*i] = mfqP->M[i+(mfqP->n+1)*j];
411: }
412: }
413: blasnp = mfqP->nmodelpoints+1;
414: /* Q_tmp,R = qr(M') */
415: blasmaxmn=PetscMax(mfqP->m,mfqP->n+1);
416: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&blasnp,&blasnplus1,mfqP->Q_tmp,&blasnpmax,mfqP->tau_tmp,mfqP->mwork,&blasmaxmn,&info));
417: if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine geqrf returned with value %d\n",info);
419: /* Reject if min(svd(N*Q(:,n+2:np+1)) <= theta2 */
420: /* L = N*Qtmp */
421: blasint2 = mfqP->n * (mfqP->n+1) / 2;
422: /* Copy N to L_tmp */
423: for (i=0;i<mfqP->n*(mfqP->n+1)/2 * mfqP->npmax;i++) {
424: mfqP->L_tmp[i]= mfqP->N[i];
425: }
426: /* Copy L_save to L_tmp */
428: /* L_tmp = N*Qtmp' */
429: PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("R","N",&blasint2,&blasnp,&blasnplus1,mfqP->Q_tmp,&blasnpmax,mfqP->tau_tmp,mfqP->L_tmp,&blasint2,mfqP->npmaxwork,&blasnmax,&info));
430: if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine ormqr returned with value %d\n",info);
432: /* Copy L_tmp to L_save */
433: for (i=0;i<mfqP->npmax * mfqP->n*(mfqP->n+1)/2;i++) {
434: mfqP->L_save[i] = mfqP->L_tmp[i];
435: }
437: /* Get svd for L_tmp(:,n+2:np+1) (L_tmp is modified in process) */
438: blasint = mfqP->nmodelpoints - mfqP->n;
439: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
440: 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));
441: PetscFPTrapPop();
442: if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine gesvd returned with value %d\n",info);
444: if (mfqP->beta[PetscMin(blasint,blasint2)-1] > mfqP->theta2) {
445: /* accept point */
446: mfqP->model_indices[mfqP->nmodelpoints] = point;
447: /* Copy Q_tmp to Q */
448: for (i=0;i<mfqP->npmax* mfqP->npmax;i++) {
449: mfqP->Q[i] = mfqP->Q_tmp[i];
450: }
451: for (i=0;i<mfqP->npmax;i++){
452: mfqP->tau[i] = mfqP->tau_tmp[i];
453: }
454: mfqP->nmodelpoints++;
455: blasnp = mfqP->nmodelpoints;
457: /* Copy L_save to L */
458: for (i=0;i<mfqP->npmax * mfqP->n*(mfqP->n+1)/2;i++) {
459: mfqP->L[i] = mfqP->L_save[i];
460: }
461: }
462: point--;
463: }
465: blasnp = mfqP->nmodelpoints;
466: /* Copy Q(:,n+2:np) to Z */
467: /* First set Q_tmp to I */
468: for (i=0;i<mfqP->npmax*mfqP->npmax;i++) {
469: mfqP->Q_tmp[i] = 0.0;
470: }
471: for (i=0;i<mfqP->npmax;i++) {
472: mfqP->Q_tmp[i + mfqP->npmax*i] = 1.0;
473: }
475: /* Q_tmp = I * Q */
476: PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("R","N",&blasnp,&blasnp,&blasnplus1,mfqP->Q,&blasnpmax,mfqP->tau,mfqP->Q_tmp,&blasnpmax,mfqP->npmaxwork,&blasnmax,&info));
477: if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine ormqr returned with value %d\n",info);
479: /* Copy Q_tmp(:,n+2:np) to Z) */
480: offset = mfqP->npmax * (mfqP->n+1);
481: for (i=offset;i<mfqP->npmax*mfqP->npmax;i++) {
482: mfqP->Z[i-offset] = mfqP->Q_tmp[i];
483: }
485: if (mfqP->nmodelpoints == mfqP->n + 1) {
486: /* Set L to I_{n+1} */
487: for (i=0;i<mfqP->npmax * mfqP->n*(mfqP->n+1)/2;i++) {
488: mfqP->L[i] = 0.0;
489: }
490: for (i=0;i<mfqP->n;i++) {
491: mfqP->L[(mfqP->n*(mfqP->n+1)/2)*i + i] = 1.0;
492: }
493: }
494: return(0);
495: }
497: /* Only call from modelimprove, addpoint() needs ->Q_tmp and ->work to be set */
498: PetscErrorCode addpoint(Tao tao, TAO_POUNDERS *mfqP, PetscInt index)
499: {
503: /* Create new vector in history: X[newidx] = X[mfqP->index] + delta*X[index]*/
504: VecDuplicate(mfqP->Xhist[0],&mfqP->Xhist[mfqP->nHist]);
505: VecSetValues(mfqP->Xhist[mfqP->nHist],mfqP->n,mfqP->indices,&mfqP->Q_tmp[index*mfqP->npmax],INSERT_VALUES);
506: VecAssemblyBegin(mfqP->Xhist[mfqP->nHist]);
507: VecAssemblyEnd(mfqP->Xhist[mfqP->nHist]);
508: VecAYPX(mfqP->Xhist[mfqP->nHist],mfqP->delta,mfqP->Xhist[mfqP->minindex]);
510: /* Project into feasible region */
511: if (tao->XU && tao->XL) {
512: VecMedian(mfqP->Xhist[mfqP->nHist], tao->XL, tao->XU, mfqP->Xhist[mfqP->nHist]);
513: }
515: /* Compute value of new vector */
516: VecDuplicate(mfqP->Fhist[0],&mfqP->Fhist[mfqP->nHist]);
517: CHKMEMQ;
518: pounders_feval(tao,mfqP->Xhist[mfqP->nHist],mfqP->Fhist[mfqP->nHist],&mfqP->Fres[mfqP->nHist]);
520: /* Add new vector to model */
521: mfqP->model_indices[mfqP->nmodelpoints] = mfqP->nHist;
522: mfqP->nmodelpoints++;
523: mfqP->nHist++;
524: return(0);
525: }
527: PetscErrorCode modelimprove(Tao tao, TAO_POUNDERS *mfqP, PetscInt addallpoints)
528: {
529: /* modeld = Q(:,np+1:n)' */
531: PetscInt i,j,minindex=0;
532: PetscReal dp,half=0.5,one=1.0,minvalue=PETSC_INFINITY;
533: PetscBLASInt blasn=mfqP->n, blasnpmax = mfqP->npmax, blask,info;
534: PetscBLASInt blas1=1,blasnmax = mfqP->nmax;
536: blask = mfqP->nmodelpoints;
537: /* Qtmp = I(n x n) */
538: for (i=0;i<mfqP->n;i++) {
539: for (j=0;j<mfqP->n;j++) {
540: mfqP->Q_tmp[i + mfqP->npmax*j] = 0.0;
541: }
542: }
543: for (j=0;j<mfqP->n;j++) {
544: mfqP->Q_tmp[j + mfqP->npmax*j] = 1.0;
545: }
547: /* Qtmp = Q * I */
548: PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("R","N",&blasn,&blasn,&blask,mfqP->Q,&blasnpmax,mfqP->tau, mfqP->Q_tmp, &blasnpmax, mfqP->npmaxwork,&blasnmax, &info));
550: for (i=mfqP->nmodelpoints;i<mfqP->n;i++) {
551: dp = BLASdot_(&blasn,&mfqP->Q_tmp[i*mfqP->npmax],&blas1,mfqP->Gres,&blas1);
552: if (dp>0.0) { /* Model says use the other direction! */
553: for (j=0;j<mfqP->n;j++) {
554: mfqP->Q_tmp[i*mfqP->npmax+j] *= -1;
555: }
556: }
557: /* mfqP->work[i] = Cres+Modeld(i,:)*(Gres+.5*Hres*Modeld(i,:)') */
558: for (j=0;j<mfqP->n;j++) {
559: mfqP->work2[j] = mfqP->Gres[j];
560: }
561: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasn,&half,mfqP->Hres,&blasn,&mfqP->Q_tmp[i*mfqP->npmax], &blas1, &one, mfqP->work2,&blas1));
562: mfqP->work[i] = BLASdot_(&blasn,&mfqP->Q_tmp[i*mfqP->npmax],&blas1,mfqP->work2,&blas1);
563: if (i==mfqP->nmodelpoints || mfqP->work[i] < minvalue) {
564: minindex=i;
565: minvalue = mfqP->work[i];
566: }
567: if (addallpoints != 0) {
568: addpoint(tao,mfqP,i);
569: }
570: }
571: if (!addallpoints) {
572: addpoint(tao,mfqP,minindex);
573: }
574: return(0);
575: }
578: PetscErrorCode affpoints(TAO_POUNDERS *mfqP, PetscReal *xmin,PetscReal c)
579: {
580: PetscInt i,j;
581: PetscBLASInt blasm=mfqP->m,blasj,blask,blasn=mfqP->n,ione=1,info;
582: PetscBLASInt blasnpmax = mfqP->npmax,blasmaxmn;
583: PetscReal proj,normd;
584: const PetscReal *x;
585: PetscErrorCode ierr;
588: for (i=mfqP->nHist-1;i>=0;i--) {
589: VecGetArrayRead(mfqP->Xhist[i],&x);
590: for (j=0;j<mfqP->n;j++) {
591: mfqP->work[j] = (x[j] - xmin[j])/mfqP->delta;
592: }
593: VecRestoreArrayRead(mfqP->Xhist[i],&x);
594: PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasn,mfqP->work,&ione,mfqP->work2,&ione));
595: normd = BLASnrm2_(&blasn,mfqP->work,&ione);
596: if (normd <= c*c) {
597: blasj=PetscMax((mfqP->n - mfqP->nmodelpoints),0);
598: if (!mfqP->q_is_I) {
599: /* project D onto null */
600: blask=(mfqP->nmodelpoints);
601: PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("R","N",&ione,&blasn,&blask,mfqP->Q,&blasnpmax,mfqP->tau,mfqP->work2,&ione,mfqP->mwork,&blasm,&info));
602: if (info < 0) SETERRQ1(PETSC_COMM_SELF,1,"ormqr returned value %d\n",info);
603: }
604: proj = BLASnrm2_(&blasj,&mfqP->work2[mfqP->nmodelpoints],&ione);
606: if (proj >= mfqP->theta1) { /* add this index to model */
607: mfqP->model_indices[mfqP->nmodelpoints]=i;
608: mfqP->nmodelpoints++;
609: PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasn,mfqP->work,&ione,&mfqP->Q_tmp[mfqP->npmax*(mfqP->nmodelpoints-1)],&ione));
610: blask=mfqP->npmax*(mfqP->nmodelpoints);
611: PetscStackCallBLAS("BLAScopy",BLAScopy_(&blask,mfqP->Q_tmp,&ione,mfqP->Q,&ione));
612: blask = mfqP->nmodelpoints;
613: blasmaxmn = PetscMax(mfqP->m,mfqP->n);
614: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&blasn,&blask,mfqP->Q,&blasnpmax,mfqP->tau,mfqP->mwork,&blasmaxmn,&info));
615: if (info < 0) SETERRQ1(PETSC_COMM_SELF,1,"geqrf returned value %d\n",info);
616: mfqP->q_is_I = 0;
617: }
618: if (mfqP->nmodelpoints == mfqP->n) {
619: break;
620: }
621: }
622: }
624: return(0);
625: }
627: static PetscErrorCode TaoSolve_POUNDERS(Tao tao)
628: {
629: TAO_POUNDERS *mfqP = (TAO_POUNDERS *)tao->data;
630: PetscInt i,ii,j,k,l;
631: PetscReal step=1.0;
632: TaoConvergedReason reason = TAO_CONTINUE_ITERATING;
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;
666: VecCopy(tao->solution,mfqP->Xhist[0]);
667: VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);
668: VecMax(mfqP->Xhist[0],NULL,&val);
669: if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 > upper bound");
671: /* Check x0 >= xl */
672: VecCopy(tao->XL,mfqP->Xhist[0]);
673: VecAXPY(mfqP->Xhist[0],-1.0,tao->solution);
674: VecMax(mfqP->Xhist[0],NULL,&val);
675: if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 < lower bound");
677: /* Check x0 + delta < XU -- should be able to get around this eventually */
679: VecSet(mfqP->Xhist[0],mfqP->delta);
680: VecAXPY(mfqP->Xhist[0],1.0,tao->solution);
681: VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);
682: VecMax(mfqP->Xhist[0],NULL,&val);
683: if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 + delta > upper bound");
684: }
686: CHKMEMQ;
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]);
691: CHKMEMQ;
692: pounders_feval(tao,tao->solution,mfqP->Fhist[0],&mfqP->Fres[0]);
693: mfqP->minindex = 0;
694: minnorm = mfqP->Fres[0];
695: TaoMonitor(tao, tao->niter, minnorm, PETSC_INFINITY, 0.0, step, &reason);
696: tao->niter++;
698: VecGetOwnershipRange(mfqP->Xhist[0],&low,&high);
699: for (i=1;i<mfqP->n+1;i++) {
700: VecCopy(tao->solution,mfqP->Xhist[i]);
701: if (i-1 >= low && i-1 < high) {
702: VecGetArray(mfqP->Xhist[i],&x);
703: x[i-1-low] += mfqP->delta;
704: VecRestoreArray(mfqP->Xhist[i],&x);
705: }
706: CHKMEMQ;
707: pounders_feval(tao,mfqP->Xhist[i],mfqP->Fhist[i],&mfqP->Fres[i]);
708: if (mfqP->Fres[i] < minnorm) {
709: mfqP->minindex = i;
710: minnorm = mfqP->Fres[i];
711: }
712: }
713: VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);
714: VecCopy(mfqP->Fhist[mfqP->minindex],tao->sep_objective);
715: /* Gather mpi vecs to one big local vec */
717: /* Begin serial code */
719: /* Disp[i] = Xi-xmin, i=1,..,mfqP->minindex-1,mfqP->minindex+1,..,n */
720: /* Fdiff[i] = (Fi-Fmin)', i=1,..,mfqP->minindex-1,mfqP->minindex+1,..,n */
721: /* (Column oriented for blas calls) */
722: ii=0;
724: if (mfqP->size == 1) {
725: VecGetArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);
726: for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i];
727: VecRestoreArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);
728: VecGetArrayRead(mfqP->Fhist[mfqP->minindex],&fmin);
729: for (i=0;i<mfqP->n+1;i++) {
730: if (i == mfqP->minindex) continue;
732: VecGetArray(mfqP->Xhist[i],&x);
733: for (j=0;j<mfqP->n;j++) {
734: mfqP->Disp[ii+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta;
735: }
736: VecRestoreArray(mfqP->Xhist[i],&x);
738: VecGetArray(mfqP->Fhist[i],&f);
739: for (j=0;j<mfqP->m;j++) {
740: mfqP->Fdiff[ii+mfqP->n*j] = f[j] - fmin[j];
741: }
742: VecRestoreArray(mfqP->Fhist[i],&f);
743: mfqP->model_indices[ii++] = i;
745: }
746: for (j=0;j<mfqP->m;j++) {
747: mfqP->C[j] = fmin[j];
748: }
749: VecRestoreArrayRead(mfqP->Fhist[mfqP->minindex],&fmin);
750: } else {
751: VecScatterBegin(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);
752: VecScatterEnd(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);
753: VecGetArrayRead(mfqP->localxmin,&xmint);
754: for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i];
755: VecRestoreArrayRead(mfqP->localxmin,&xmint);
757: VecScatterBegin(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);
758: VecScatterEnd(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);
759: VecGetArrayRead(mfqP->localfmin,&fmin);
760: for (i=0;i<mfqP->n+1;i++) {
761: if (i == mfqP->minindex) continue;
763: mfqP->model_indices[ii++] = i;
764: VecScatterBegin(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);
765: VecScatterEnd(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);
766: VecGetArray(mfqP->localx,&x);
767: for (j=0;j<mfqP->n;j++) {
768: mfqP->Disp[i+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta;
769: }
770: VecRestoreArray(mfqP->localx,&x);
772: VecScatterBegin(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);
773: VecScatterEnd(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);
774: VecGetArray(mfqP->localf,&f);
775: for (j=0;j<mfqP->m;j++) {
776: mfqP->Fdiff[i*mfqP->n+j] = f[j] - fmin[j];
777: }
778: VecRestoreArray(mfqP->localf,&f);
779: }
780: for (j=0;j<mfqP->m;j++) {
781: mfqP->C[j] = fmin[j];
782: }
783: VecRestoreArrayRead(mfqP->localfmin,&fmin);
784: }
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: PetscInfo1(tao,"gesv returned %d\n",info);
793: cres = minnorm;
794: pounders_update_res(tao);
796: valid = PETSC_TRUE;
798: VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);
799: VecAssemblyBegin(tao->gradient);
800: VecAssemblyEnd(tao->gradient);
801: VecNorm(tao->gradient,NORM_2,&gnorm);
802: gnorm *= mfqP->delta;
803: VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);
804: TaoMonitor(tao, tao->niter, minnorm, gnorm, 0.0, step, &reason);
805: mfqP->nHist = mfqP->n+1;
806: mfqP->nmodelpoints = mfqP->n+1;
808: while (reason == TAO_CONTINUE_ITERATING) {
809: PetscReal gnm = 1e-4;
810: tao->niter++;
811: /* Solve the subproblem min{Q(s): ||s|| <= delta} */
812: gqtwrap(tao,&gnm,&mdec);
813: /* Evaluate the function at the new point */
815: for (i=0;i<mfqP->n;i++) {
816: mfqP->work[i] = mfqP->Xsubproblem[i]*mfqP->delta + mfqP->xmin[i];
817: }
818: VecDuplicate(tao->solution,&mfqP->Xhist[mfqP->nHist]);
819: VecDuplicate(tao->sep_objective,&mfqP->Fhist[mfqP->nHist]);
820: VecSetValues(mfqP->Xhist[mfqP->nHist],mfqP->n,mfqP->indices,mfqP->work,INSERT_VALUES);
821: VecAssemblyBegin(mfqP->Xhist[mfqP->nHist]);
822: VecAssemblyEnd(mfqP->Xhist[mfqP->nHist]);
824: pounders_feval(tao,mfqP->Xhist[mfqP->nHist],mfqP->Fhist[mfqP->nHist],&mfqP->Fres[mfqP->nHist]);
826: rho = (mfqP->Fres[mfqP->minindex] - mfqP->Fres[mfqP->nHist]) / mdec;
827: mfqP->nHist++;
829: /* Update the center */
830: if ((rho >= mfqP->eta1) || (rho > mfqP->eta0 && valid==PETSC_TRUE)) {
831: /* Update model to reflect new base point */
832: for (i=0;i<mfqP->n;i++) {
833: mfqP->work[i] = (mfqP->work[i] - mfqP->xmin[i])/mfqP->delta;
834: }
835: for (j=0;j<mfqP->m;j++) {
836: /* C(j) = C(j) + work*G(:,j) + .5*work*H(:,:,j)*work';
837: G(:,j) = G(:,j) + H(:,:,j)*work' */
838: for (k=0;k<mfqP->n;k++) {
839: mfqP->work2[k]=0.0;
840: for (l=0;l<mfqP->n;l++) {
841: mfqP->work2[k]+=mfqP->H[j + mfqP->m*(k + l*mfqP->n)]*mfqP->work[l];
842: }
843: }
844: for (i=0;i<mfqP->n;i++) {
845: mfqP->C[j]+=mfqP->work[i]*(mfqP->Fdiff[i + mfqP->n* j] + 0.5*mfqP->work2[i]);
846: mfqP->Fdiff[i+mfqP->n*j] +=mfqP-> work2[i];
847: }
848: }
849: /* Cres += work*Gres + .5*work*Hres*work';
850: Gres += Hres*work'; */
852: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasn,&one,mfqP->Hres,&blasn,mfqP->work,&ione,&zero,mfqP->work2,&ione));
853: for (i=0;i<mfqP->n;i++) {
854: cres += mfqP->work[i]*(mfqP->Gres[i] + 0.5*mfqP->work2[i]);
855: mfqP->Gres[i] += mfqP->work2[i];
856: }
857: mfqP->minindex = mfqP->nHist-1;
858: minnorm = mfqP->Fres[mfqP->minindex];
859: VecCopy(mfqP->Fhist[mfqP->minindex],tao->sep_objective);
860: /* Change current center */
861: VecGetArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);
862: for (i=0;i<mfqP->n;i++) {
863: mfqP->xmin[i] = xmint[i];
864: }
865: VecRestoreArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);
866: }
868: /* Evaluate at a model-improving point if necessary */
869: if (valid == PETSC_FALSE) {
870: mfqP->q_is_I = 1;
871: mfqP->nmodelpoints = 0;
872: affpoints(mfqP,mfqP->xmin,mfqP->c1);
873: if (mfqP->nmodelpoints < mfqP->n) {
874: PetscInfo(tao,"Model not valid -- model-improving\n");
875: modelimprove(tao,mfqP,1);
876: }
877: }
879: /* Update the trust region radius */
880: deltaold = mfqP->delta;
881: normxsp = 0;
882: for (i=0;i<mfqP->n;i++) {
883: normxsp += mfqP->Xsubproblem[i]*mfqP->Xsubproblem[i];
884: }
885: normxsp = PetscSqrtReal(normxsp);
886: if (rho >= mfqP->eta1 && normxsp > 0.5*mfqP->delta) {
887: mfqP->delta = PetscMin(mfqP->delta*mfqP->gamma1,mfqP->deltamax);
888: } else if (valid == PETSC_TRUE) {
889: mfqP->delta = PetscMax(mfqP->delta*mfqP->gamma0,mfqP->deltamin);
890: }
892: /* Compute the next interpolation set */
893: mfqP->q_is_I = 1;
894: mfqP->nmodelpoints=0;
895: affpoints(mfqP,mfqP->xmin,mfqP->c1);
896: if (mfqP->nmodelpoints == mfqP->n) {
897: valid = PETSC_TRUE;
898: } else {
899: valid = PETSC_FALSE;
900: affpoints(mfqP,mfqP->xmin,mfqP->c2);
901: if (mfqP->n > mfqP->nmodelpoints) {
902: PetscInfo(tao,"Model not valid -- adding geometry points\n");
903: modelimprove(tao,mfqP,mfqP->n - mfqP->nmodelpoints);
904: }
905: }
906: for (i=mfqP->nmodelpoints;i>0;i--) {
907: mfqP->model_indices[i] = mfqP->model_indices[i-1];
908: }
909: mfqP->nmodelpoints++;
910: mfqP->model_indices[0] = mfqP->minindex;
911: morepoints(mfqP);
912: for (i=0;i<mfqP->nmodelpoints;i++) {
913: VecGetArray(mfqP->Xhist[mfqP->model_indices[i]],&x);
914: for (j=0;j<mfqP->n;j++) {
915: mfqP->Disp[i + mfqP->npmax*j] = (x[j] - mfqP->xmin[j]) / deltaold;
916: }
917: VecRestoreArray(mfqP->Xhist[mfqP->model_indices[i]],&x);
918: VecGetArray(mfqP->Fhist[mfqP->model_indices[i]],&f);
919: for (j=0;j<mfqP->m;j++) {
920: for (k=0;k<mfqP->n;k++) {
921: mfqP->work[k]=0.0;
922: for (l=0;l<mfqP->n;l++) {
923: mfqP->work[k] += mfqP->H[j + mfqP->m*(k + mfqP->n*l)] * mfqP->Disp[i + mfqP->npmax*l];
924: }
925: }
926: 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];
927: }
928: VecRestoreArray(mfqP->Fhist[mfqP->model_indices[i]],&f);
929: }
931: /* Update the quadratic model */
932: getquadpounders(mfqP);
933: VecGetArrayRead(mfqP->Fhist[mfqP->minindex],&fmin);
934: PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasm,fmin,&ione,mfqP->C,&ione));
935: /* G = G*(delta/deltaold) + Gdel */
936: ratio = mfqP->delta/deltaold;
937: iblas = blasm*blasn;
938: PetscStackCallBLAS("BLASscal",BLASscal_(&iblas,&ratio,mfqP->Fdiff,&ione));
939: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&iblas,&one,mfqP->Gdel,&ione,mfqP->Fdiff,&ione));
940: /* H = H*(delta/deltaold)^2 + Hdel */
941: iblas = blasm*blasn*blasn;
942: ratio *= ratio;
943: PetscStackCallBLAS("BLASscal",BLASscal_(&iblas,&ratio,mfqP->H,&ione));
944: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&iblas,&one,mfqP->Hdel,&ione,mfqP->H,&ione));
946: /* Get residuals */
947: cres = mfqP->Fres[mfqP->minindex];
948: pounders_update_res(tao);
950: /* Export solution and gradient residual to TAO */
951: VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);
952: VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);
953: VecAssemblyBegin(tao->gradient);
954: VecAssemblyEnd(tao->gradient);
955: VecNorm(tao->gradient,NORM_2,&gnorm);
956: gnorm *= mfqP->delta;
957: /* final criticality test */
958: TaoMonitor(tao, tao->niter, minnorm, gnorm, 0.0, step, &reason);
959: /* test for repeated model */
960: if (mfqP->nmodelpoints==mfqP->last_nmodelpoints) {
961: same = PETSC_TRUE;
962: } else {
963: same = PETSC_FALSE;
964: }
965: for (i=0;i<mfqP->nmodelpoints;i++) {
966: if (same) {
967: if (mfqP->model_indices[i] == mfqP->last_model_indices[i]) {
968: same = PETSC_TRUE;
969: } else {
970: same = PETSC_FALSE;
971: }
972: }
973: mfqP->last_model_indices[i] = mfqP->model_indices[i];
974: }
975: mfqP->last_nmodelpoints = mfqP->nmodelpoints;
976: if (same && mfqP->delta == deltaold) {
977: PetscInfo(tao,"Identical model used in successive iterations\n");
978: reason = TAO_CONVERGED_STEPTOL;
979: tao->reason = TAO_CONVERGED_STEPTOL;
980: }
981: }
982: return(0);
983: }
985: static PetscErrorCode TaoSetUp_POUNDERS(Tao tao)
986: {
987: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
988: PetscInt i,j;
989: IS isfloc,isfglob,isxloc,isxglob;
993: if (!tao->gradient) {VecDuplicate(tao->solution,&tao->gradient); }
994: if (!tao->stepdirection) {VecDuplicate(tao->solution,&tao->stepdirection); }
995: VecGetSize(tao->solution,&mfqP->n);
996: VecGetSize(tao->sep_objective,&mfqP->m);
997: mfqP->c1 = PetscSqrtReal((PetscReal)mfqP->n);
998: if (mfqP->npmax == PETSC_DEFAULT) {
999: mfqP->npmax = 2*mfqP->n + 1;
1000: }
1001: mfqP->npmax = PetscMin((mfqP->n+1)*(mfqP->n+2)/2,mfqP->npmax);
1002: mfqP->npmax = PetscMax(mfqP->npmax, mfqP->n+2);
1004: PetscMalloc1(tao->max_funcs+10,&mfqP->Xhist);
1005: PetscMalloc1(tao->max_funcs+10,&mfqP->Fhist);
1006: for (i=0;i<mfqP->n +1;i++) {
1007: VecDuplicate(tao->solution,&mfqP->Xhist[i]);
1008: VecDuplicate(tao->sep_objective,&mfqP->Fhist[i]);
1009: }
1010: VecDuplicate(tao->solution,&mfqP->workxvec);
1011: VecDuplicate(tao->sep_objective,&mfqP->workfvec);
1012: mfqP->nHist = 0;
1014: PetscMalloc1(tao->max_funcs+10,&mfqP->Fres);
1015: PetscMalloc1(mfqP->npmax*mfqP->m,&mfqP->RES);
1016: PetscMalloc1(mfqP->n,&mfqP->work);
1017: PetscMalloc1(mfqP->n,&mfqP->work2);
1018: PetscMalloc1(mfqP->n,&mfqP->work3);
1019: PetscMalloc1(PetscMax(mfqP->m,mfqP->n+1),&mfqP->mwork);
1020: PetscMalloc1(mfqP->npmax - mfqP->n - 1,&mfqP->omega);
1021: PetscMalloc1(mfqP->n * (mfqP->n+1) / 2,&mfqP->beta);
1022: PetscMalloc1(mfqP->n + 1 ,&mfqP->alpha);
1024: PetscMalloc1(mfqP->n*mfqP->n*mfqP->m,&mfqP->H);
1025: PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q);
1026: PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q_tmp);
1027: PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L);
1028: PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_tmp);
1029: PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_save);
1030: PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->N);
1031: PetscMalloc1(mfqP->npmax * (mfqP->n+1) ,&mfqP->M);
1032: PetscMalloc1(mfqP->npmax * (mfqP->npmax - mfqP->n - 1) , &mfqP->Z);
1033: PetscMalloc1(mfqP->npmax,&mfqP->tau);
1034: PetscMalloc1(mfqP->npmax,&mfqP->tau_tmp);
1035: mfqP->nmax = PetscMax(5*mfqP->npmax,mfqP->n*(mfqP->n+1)/2);
1036: PetscMalloc1(mfqP->nmax,&mfqP->npmaxwork);
1037: PetscMalloc1(mfqP->nmax,&mfqP->npmaxiwork);
1038: PetscMalloc1(mfqP->n,&mfqP->xmin);
1039: PetscMalloc1(mfqP->m,&mfqP->C);
1040: PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Fdiff);
1041: PetscMalloc1(mfqP->npmax*mfqP->n,&mfqP->Disp);
1042: PetscMalloc1(mfqP->n,&mfqP->Gres);
1043: PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Hres);
1044: PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Gpoints);
1045: PetscMalloc1(mfqP->npmax,&mfqP->model_indices);
1046: PetscMalloc1(mfqP->npmax,&mfqP->last_model_indices);
1047: PetscMalloc1(mfqP->n,&mfqP->Xsubproblem);
1048: PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Gdel);
1049: PetscMalloc1(mfqP->n*mfqP->n*mfqP->m, &mfqP->Hdel);
1050: PetscMalloc1(PetscMax(mfqP->m,mfqP->n),&mfqP->indices);
1051: PetscMalloc1(mfqP->n,&mfqP->iwork);
1052: PetscMalloc1(mfqP->m*mfqP->m,&mfqP->w);
1053: for (i=0;i<mfqP->m;i++) {
1054: for (j=0;j<mfqP->m;j++) {
1055: if (i==j) {
1056: mfqP->w[i+mfqP->m*j]=1.0;
1057: } else {
1058: mfqP->w[i+mfqP->m*j]=0.0;
1059: }
1060: }
1061: }
1062: for (i=0;i<PetscMax(mfqP->m,mfqP->n);i++) {
1063: mfqP->indices[i] = i;
1064: }
1065: MPI_Comm_size(((PetscObject)tao)->comm,&mfqP->size);
1066: if (mfqP->size > 1) {
1067: VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localx);
1068: VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localxmin);
1069: VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localf);
1070: VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localfmin);
1071: ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxloc);
1072: ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxglob);
1073: ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfloc);
1074: ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfglob);
1077: VecScatterCreate(tao->solution,isxglob,mfqP->localx,isxloc,&mfqP->scatterx);
1078: VecScatterCreate(tao->sep_objective,isfglob,mfqP->localf,isfloc,&mfqP->scatterf);
1080: ISDestroy(&isxloc);
1081: ISDestroy(&isxglob);
1082: ISDestroy(&isfloc);
1083: ISDestroy(&isfglob);
1084: }
1086: if (!mfqP->usegqt) {
1087: KSP ksp;
1088: PC pc;
1089: VecCreateSeqWithArray(PETSC_COMM_SELF,mfqP->n,mfqP->n,mfqP->Xsubproblem,&mfqP->subx);
1090: VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->subxl);
1091: VecDuplicate(mfqP->subxl,&mfqP->subb);
1092: VecDuplicate(mfqP->subxl,&mfqP->subxu);
1093: VecDuplicate(mfqP->subxl,&mfqP->subpdel);
1094: VecDuplicate(mfqP->subxl,&mfqP->subndel);
1095: TaoCreate(PETSC_COMM_SELF,&mfqP->subtao);
1096: TaoSetType(mfqP->subtao,TAOTRON);
1097: TaoSetOptionsPrefix(mfqP->subtao,"pounders_subsolver_");
1098: TaoSetInitialVector(mfqP->subtao,mfqP->subx);
1099: TaoSetObjectiveAndGradientRoutine(mfqP->subtao,pounders_fg,(void*)mfqP);
1100: TaoSetMaximumIterations(mfqP->subtao,mfqP->gqt_maxits);
1101: TaoSetFromOptions(mfqP->subtao);
1102: TaoGetKSP(mfqP->subtao,&ksp);
1103: if (ksp) {
1104: KSPGetPC(ksp,&pc);
1105: PCSetType(pc,PCNONE);
1106: }
1107: TaoSetVariableBounds(mfqP->subtao,mfqP->subxl,mfqP->subxu);
1108: MatCreateSeqDense(PETSC_COMM_SELF,mfqP->n,mfqP->n,mfqP->Hres,&mfqP->subH);
1109: TaoSetHessianRoutine(mfqP->subtao,mfqP->subH,mfqP->subH,pounders_h,(void*)mfqP);
1110: }
1111: return(0);
1112: }
1114: static PetscErrorCode TaoDestroy_POUNDERS(Tao tao)
1115: {
1116: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
1117: PetscInt i;
1121: if (!mfqP->usegqt) {
1122: TaoDestroy(&mfqP->subtao);
1123: VecDestroy(&mfqP->subx);
1124: VecDestroy(&mfqP->subxl);
1125: VecDestroy(&mfqP->subxu);
1126: VecDestroy(&mfqP->subb);
1127: VecDestroy(&mfqP->subpdel);
1128: VecDestroy(&mfqP->subndel);
1129: MatDestroy(&mfqP->subH);
1130: }
1131: PetscFree(mfqP->Fres);
1132: PetscFree(mfqP->RES);
1133: PetscFree(mfqP->work);
1134: PetscFree(mfqP->work2);
1135: PetscFree(mfqP->work3);
1136: PetscFree(mfqP->mwork);
1137: PetscFree(mfqP->omega);
1138: PetscFree(mfqP->beta);
1139: PetscFree(mfqP->alpha);
1140: PetscFree(mfqP->H);
1141: PetscFree(mfqP->Q);
1142: PetscFree(mfqP->Q_tmp);
1143: PetscFree(mfqP->L);
1144: PetscFree(mfqP->L_tmp);
1145: PetscFree(mfqP->L_save);
1146: PetscFree(mfqP->N);
1147: PetscFree(mfqP->M);
1148: PetscFree(mfqP->Z);
1149: PetscFree(mfqP->tau);
1150: PetscFree(mfqP->tau_tmp);
1151: PetscFree(mfqP->npmaxwork);
1152: PetscFree(mfqP->npmaxiwork);
1153: PetscFree(mfqP->xmin);
1154: PetscFree(mfqP->C);
1155: PetscFree(mfqP->Fdiff);
1156: PetscFree(mfqP->Disp);
1157: PetscFree(mfqP->Gres);
1158: PetscFree(mfqP->Hres);
1159: PetscFree(mfqP->Gpoints);
1160: PetscFree(mfqP->model_indices);
1161: PetscFree(mfqP->last_model_indices);
1162: PetscFree(mfqP->Xsubproblem);
1163: PetscFree(mfqP->Gdel);
1164: PetscFree(mfqP->Hdel);
1165: PetscFree(mfqP->indices);
1166: PetscFree(mfqP->iwork);
1167: PetscFree(mfqP->w);
1168: for (i=0;i<mfqP->nHist;i++) {
1169: VecDestroy(&mfqP->Xhist[i]);
1170: VecDestroy(&mfqP->Fhist[i]);
1171: }
1172: VecDestroy(&mfqP->workxvec);
1173: VecDestroy(&mfqP->workfvec);
1174: PetscFree(mfqP->Xhist);
1175: PetscFree(mfqP->Fhist);
1177: if (mfqP->size > 1) {
1178: VecDestroy(&mfqP->localx);
1179: VecDestroy(&mfqP->localxmin);
1180: VecDestroy(&mfqP->localf);
1181: VecDestroy(&mfqP->localfmin);
1182: }
1183: PetscFree(tao->data);
1184: return(0);
1185: }
1187: static PetscErrorCode TaoSetFromOptions_POUNDERS(PetscOptionItems *PetscOptionsObject,Tao tao)
1188: {
1189: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
1193: PetscOptionsHead(PetscOptionsObject,"POUNDERS method for least-squares optimization");
1194: PetscOptionsReal("-tao_pounders_delta","initial delta","",mfqP->delta,&mfqP->delta0,NULL);
1195: mfqP->delta = mfqP->delta0;
1196: mfqP->npmax = PETSC_DEFAULT;
1197: PetscOptionsInt("-tao_pounders_npmax","max number of points in model","",mfqP->npmax,&mfqP->npmax,NULL);
1198: mfqP->usegqt = PETSC_FALSE;
1199: PetscOptionsBool("-tao_pounders_gqt","use gqt algorithm for subproblem","",mfqP->usegqt,&mfqP->usegqt,NULL);
1200: PetscOptionsTail();
1201: return(0);
1202: }
1204: static PetscErrorCode TaoView_POUNDERS(Tao tao, PetscViewer viewer)
1205: {
1206: TAO_POUNDERS *mfqP = (TAO_POUNDERS *)tao->data;
1207: PetscBool isascii;
1208: PetscInt nits;
1212: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1213: if (isascii) {
1214: PetscViewerASCIIPushTab(viewer);
1215: PetscViewerASCIIPrintf(viewer, "initial delta: %g\n",(double)mfqP->delta0);
1216: PetscViewerASCIIPrintf(viewer, "final delta: %g\n",(double)mfqP->delta);
1217: PetscViewerASCIIPrintf(viewer, "model points: %D\n",mfqP->nmodelpoints);
1218: if (mfqP->usegqt) {
1219: PetscViewerASCIIPrintf(viewer, "subproblem solver: gqt\n");
1220: } else {
1221: PetscViewerASCIIPrintf(viewer, "subproblem solver: %s\n",((PetscObject)mfqP->subtao)->type_name);
1222: TaoGetTotalIterationNumber(mfqP->subtao,&nits);
1223: PetscViewerASCIIPrintf(viewer, "total subproblem iterations: %D\n",nits);
1224: }
1225: PetscViewerASCIIPopTab(viewer);
1226: }
1227: return(0);
1228: }
1229: /*MC
1230: TAOPOUNDERS - POUNDERS derivate-free model-based algorithm for nonlinear least squares
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
1237: Level: beginner
1238:
1239: M*/
1241: PETSC_EXTERN PetscErrorCode TaoCreate_POUNDERS(Tao tao)
1242: {
1243: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
1247: tao->ops->setup = TaoSetUp_POUNDERS;
1248: tao->ops->solve = TaoSolve_POUNDERS;
1249: tao->ops->view = TaoView_POUNDERS;
1250: tao->ops->setfromoptions = TaoSetFromOptions_POUNDERS;
1251: tao->ops->destroy = TaoDestroy_POUNDERS;
1253: PetscNewLog(tao,&mfqP);
1254: tao->data = (void*)mfqP;
1255: /* Override default settings (unless already changed) */
1256: if (!tao->max_it_changed) tao->max_it = 2000;
1257: if (!tao->max_funcs_changed) tao->max_funcs = 4000;
1258: mfqP->delta0 = 0.1;
1259: mfqP->delta = 0.1;
1260: mfqP->deltamax=1e3;
1261: mfqP->deltamin=1e-6;
1262: mfqP->c2 = 100.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->gqt_rtol = 0.001;
1270: mfqP->gqt_maxits = 50;
1271: mfqP->workxvec = 0;
1272: return(0);
1273: }