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