Actual source code: pounders.c
petsc-dev 2014-02-02
1: #include <../src/tao/leastsquares/impls/pounders/pounders.h>
5: static PetscErrorCode pounders_h(Tao subtao, Vec v, Mat *H, Mat *Hpre, MatStructure *flag, void *ctx)
6: {
8: *flag = SAME_NONZERO_PATTERN;
9: return(0);
10: }
13: static PetscErrorCode pounders_fg(Tao subtao, Vec x, PetscReal *f, Vec g, void *ctx)
14: {
15: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)ctx;
16: PetscReal d1,d2;
20: /* g = A*x (add b later)*/
21: MatMult(mfqP->subH,x,g);
23: /* f = 1/2 * x'*(Ax) + b'*x */
24: VecDot(x,g,&d1);
25: VecDot(mfqP->subb,x,&d2);
26: *f = 0.5 *d1 + d2;
28: /* now g = g + b */
29: VecAXPY(g, 1.0, mfqP->subb);
30: return(0);
31: }
35: PetscErrorCode gqtwrap(Tao tao,PetscReal *gnorm, PetscReal *qmin)
36: {
38: #if defined(PETSC_USE_REAL_SINGLE)
39: PetscReal atol=1.0e-5;
40: #else
41: PetscReal atol=1.0e-10;
42: #endif
43: PetscInt info,its;
44: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
47: if (! mfqP->usegqt) {
48: PetscReal maxval;
49: PetscInt i,j;
51: VecSetValues(mfqP->subb,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);
52: VecAssemblyBegin(mfqP->subb);
53: VecAssemblyEnd(mfqP->subb);
55: VecSet(mfqP->subx,0.0);
57: VecSet(mfqP->subndel,-mfqP->delta);
58: VecSet(mfqP->subpdel,mfqP->delta);
60: for (i=0;i<mfqP->n;i++) {
61: for (j=i;j<mfqP->n;j++) {
62: mfqP->Hres[j+mfqP->n*i] = mfqP->Hres[mfqP->n*j+i];
63: }
64: }
65: MatSetValues(mfqP->subH,mfqP->n,mfqP->indices,mfqP->n,mfqP->indices,mfqP->Hres,INSERT_VALUES);
66: MatAssemblyBegin(mfqP->subH,MAT_FINAL_ASSEMBLY);
67: MatAssemblyEnd(mfqP->subH,MAT_FINAL_ASSEMBLY);
69: TaoResetStatistics(mfqP->subtao);
70: TaoSetTolerances(mfqP->subtao,PETSC_DEFAULT,PETSC_DEFAULT,*gnorm,*gnorm,PETSC_DEFAULT);
71: /* enforce bound constraints -- experimental */
72: if (tao->XU && tao->XL) {
73: VecCopy(tao->XU,mfqP->subxu);
74: VecAXPY(mfqP->subxu,-1.0,tao->solution);
75: VecScale(mfqP->subxu,1.0/mfqP->delta);
76: VecCopy(tao->XL,mfqP->subxl);
77: VecAXPY(mfqP->subxl,-1.0,tao->solution);
78: VecScale(mfqP->subxl,1.0/mfqP->delta);
80: VecPointwiseMin(mfqP->subxu,mfqP->subxu,mfqP->subpdel);
81: VecPointwiseMax(mfqP->subxl,mfqP->subxl,mfqP->subndel);
82: } else {
83: VecCopy(mfqP->subpdel,mfqP->subxu);
84: VecCopy(mfqP->subndel,mfqP->subxl);
85: }
86: /* Make sure xu > xl */
87: VecCopy(mfqP->subxl,mfqP->subpdel);
88: VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);
89: VecMax(mfqP->subpdel,NULL,&maxval);
90: if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"upper bound < lower bound in subproblem");
91: /* Make sure xu > tao->solution > xl */
92: VecCopy(mfqP->subxl,mfqP->subpdel);
93: VecAXPY(mfqP->subpdel,-1.0,mfqP->subx);
94: VecMax(mfqP->subpdel,NULL,&maxval);
95: if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"initial guess < lower bound in subproblem");
97: VecCopy(mfqP->subx,mfqP->subpdel);
98: VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);
99: VecMax(mfqP->subpdel,NULL,&maxval);
100: if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"initial guess > upper bound in subproblem");
102: TaoSolve(mfqP->subtao);
103: TaoGetSolutionStatus(mfqP->subtao,NULL,qmin,NULL,NULL,NULL,NULL);
105: /* test bounds post-solution*/
106: VecCopy(mfqP->subxl,mfqP->subpdel);
107: VecAXPY(mfqP->subpdel,-1.0,mfqP->subx);
108: VecMax(mfqP->subpdel,NULL,&maxval);
109: if (maxval > 1e-5) {
110: PetscInfo(tao,"subproblem solution < lower bound");
111: tao->reason = TAO_DIVERGED_TR_REDUCTION;
112: }
114: VecCopy(mfqP->subx,mfqP->subpdel);
115: VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);
116: VecMax(mfqP->subpdel,NULL,&maxval);
117: if (maxval > 1e-5) {
118: PetscInfo(tao,"subproblem solution > upper bound");
119: tao->reason = TAO_DIVERGED_TR_REDUCTION;
120: }
121: } else {
122: 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);
123: }
124: *qmin *= -1;
125: return(0);
126: }
130: PetscErrorCode phi2eval(PetscReal *x, PetscInt n, PetscReal *phi)
131: {
132: /* 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] */
133: PetscInt i,j,k;
134: PetscReal sqrt2 = PetscSqrtReal(2.0);
137: j=0;
138: for (i=0;i<n;i++) {
139: phi[j] = 0.5 * x[i]*x[i];
140: j++;
141: for (k=i+1;k<n;k++) {
142: phi[j] = x[i]*x[k]/sqrt2;
143: j++;
144: }
145: }
146: return(0);
147: }
151: PetscErrorCode getquadpounders(TAO_POUNDERS *mfqP)
152: {
153: /* Computes the parameters of the quadratic Q(x) = c + g'*x + 0.5*x*G*x'
154: that satisfies the interpolation conditions Q(X[:,j]) = f(j)
155: for j=1,...,m and with a Hessian matrix of least Frobenius norm */
157: /* NB --we are ignoring c */
158: PetscInt i,j,k,num,np = mfqP->nmodelpoints;
159: PetscReal one = 1.0,zero=0.0,negone=-1.0;
160: PetscBLASInt blasnpmax = mfqP->npmax;
161: PetscBLASInt blasnplus1 = mfqP->n+1;
162: PetscBLASInt blasnp = np;
163: PetscBLASInt blasint = mfqP->n*(mfqP->n+1) / 2;
164: PetscBLASInt blasint2 = np - mfqP->n-1;
165: PetscBLASInt info,ione=1;
166: PetscReal sqrt2 = PetscSqrtReal(2.0);
169: for (i=0;i<mfqP->n*mfqP->m;i++) {
170: mfqP->Gdel[i] = 0;
171: }
172: for (i=0;i<mfqP->n*mfqP->n*mfqP->m;i++) {
173: mfqP->Hdel[i] = 0;
174: }
176: /* factor M */
177: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&blasnplus1,&blasnpmax,mfqP->M,&blasnplus1,mfqP->npmaxiwork,&info));
178: if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine getrf returned with value %d\n",info);
180: if (np == mfqP->n+1) {
181: for (i=0;i<mfqP->npmax-mfqP->n-1;i++) {
182: mfqP->omega[i]=0.0;
183: }
184: for (i=0;i<mfqP->n*(mfqP->n+1)/2;i++) {
185: mfqP->beta[i]=0.0;
186: }
187: } else {
188: /* Let Ltmp = (L'*L) */
189: 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));
191: /* factor Ltmp */
192: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&blasint2,mfqP->L_tmp,&blasint,&info));
193: if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine potrf returned with value %d\n",info);
194: }
196: for (k=0;k<mfqP->m;k++) {
197: if (np != mfqP->n+1) {
198: /* Solve L'*L*Omega = Z' * RESk*/
199: PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&blasnp,&blasint2,&one,mfqP->Z,&blasnpmax,&mfqP->RES[mfqP->npmax*k],&ione,&zero,mfqP->omega,&ione));
200: PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&blasint2,&ione,mfqP->L_tmp,&blasint,mfqP->omega,&blasint2,&info));
201: if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine potrs returned with value %d\n",info);
203: /* Beta = L*Omega */
204: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasint,&blasint2,&one,&mfqP->L[(mfqP->n+1)*blasint],&blasint,mfqP->omega,&ione,&zero,mfqP->beta,&ione));
205: }
207: /* solve M'*Alpha = RESk - N'*Beta */
208: PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&blasint,&blasnp,&negone,mfqP->N,&blasint,mfqP->beta,&ione,&one,&mfqP->RES[mfqP->npmax*k],&ione));
209: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&blasnplus1,&ione,mfqP->M,&blasnplus1,mfqP->npmaxiwork,&mfqP->RES[mfqP->npmax*k],&blasnplus1,&info));
210: if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine getrs returned with value %d\n",info);
212: /* Gdel(:,k) = Alpha(2:n+1) */
213: for (i=0;i<mfqP->n;i++) {
214: mfqP->Gdel[i + mfqP->n*k] = mfqP->RES[mfqP->npmax*k + i+1];
215: }
217: /* Set Hdels */
218: num=0;
219: for (i=0;i<mfqP->n;i++) {
220: /* H[i,i,k] = Beta(num) */
221: mfqP->Hdel[(i*mfqP->n + i)*mfqP->m + k] = mfqP->beta[num];
222: num++;
223: for (j=i+1;j<mfqP->n;j++) {
224: /* H[i,j,k] = H[j,i,k] = Beta(num)/sqrt(2) */
225: mfqP->Hdel[(j*mfqP->n + i)*mfqP->m + k] = mfqP->beta[num]/sqrt2;
226: mfqP->Hdel[(i*mfqP->n + j)*mfqP->m + k] = mfqP->beta[num]/sqrt2;
227: num++;
228: }
229: }
230: }
231: return(0);
232: }
236: PetscErrorCode morepoints(TAO_POUNDERS *mfqP)
237: {
238: /* Assumes mfqP->model_indices[0] is minimum index
239: Finishes adding points to mfqP->model_indices (up to npmax)
240: Computes L,Z,M,N
241: np is actual number of points in model (should equal npmax?) */
242: PetscInt point,i,j,offset;
243: PetscInt reject;
244: PetscBLASInt blasn=mfqP->n,blasnpmax=mfqP->npmax,blasnplus1=mfqP->n+1,info,blasnmax=mfqP->nmax,blasint,blasint2,blasnp,blasmaxmn;
245: PetscReal *x,normd;
249: /* Initialize M,N */
250: for (i=0;i<mfqP->n+1;i++) {
251: VecGetArray(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: VecRestoreArray(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: VecGetArray(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: VecRestoreArray(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+1;
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: PetscReal *x;
485: for (i=mfqP->nHist-1;i>=0;i--) {
486: VecGetArray(mfqP->Xhist[i],&x);
487: for (j=0;j<mfqP->n;j++) {
488: mfqP->work[j] = (x[j] - xmin[j])/mfqP->delta;
489: }
490: VecRestoreArray(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: }
520: return(0);
521: }
525: static PetscErrorCode TaoSolve_POUNDERS(Tao tao)
526: {
527: TAO_POUNDERS *mfqP = (TAO_POUNDERS *)tao->data;
528: PetscInt i,ii,j,k,l,iter=0;
529: PetscReal step=1.0;
530: TaoTerminationReason reason = TAO_CONTINUE_ITERATING;
531: PetscInt low,high;
532: PetscReal minnorm;
533: PetscReal *x,*f,*fmin,*xmint;
534: PetscReal cres,deltaold;
535: PetscReal gnorm;
536: PetscBLASInt info,ione=1,iblas;
537: PetscBool valid;
538: PetscReal mdec, rho, normxsp;
539: PetscReal one=1.0,zero=0.0,ratio;
540: PetscBLASInt blasm,blasn,blasnpmax,blasn2;
541: PetscErrorCode ierr;
543: /* n = # of parameters
544: m = dimension (components) of function */
546: if (tao->XL && tao->XU) {
547: /* Check x0 <= XU */
548: PetscReal val;
549: VecCopy(tao->solution,mfqP->Xhist[0]);
550: VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);
551: VecMax(mfqP->Xhist[0],NULL,&val);
552: if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 > upper bound");
554: /* Check x0 >= xl */
555: VecCopy(tao->XL,mfqP->Xhist[0]);
556: VecAXPY(mfqP->Xhist[0],-1.0,tao->solution);
557: VecMax(mfqP->Xhist[0],NULL,&val);
558: if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 < lower bound");
560: /* Check x0 + delta < XU -- should be able to get around this eventually */
562: VecSet(mfqP->Xhist[0],mfqP->delta);
563: VecAXPY(mfqP->Xhist[0],1.0,tao->solution);
564: VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);
565: VecMax(mfqP->Xhist[0],NULL,&val);
566: if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 + delta > upper bound");
567: }
569: CHKMEMQ;
570: blasm = mfqP->m; blasn=mfqP->n; blasnpmax = mfqP->npmax;
571: for (i=0;i<mfqP->n*mfqP->n*mfqP->m;i++) mfqP->H[i]=0;
573: VecCopy(tao->solution,mfqP->Xhist[0]);
574: CHKMEMQ;
575: TaoComputeSeparableObjective(tao,tao->solution,mfqP->Fhist[0]);
577: VecNorm(mfqP->Fhist[0],NORM_2,&mfqP->Fres[0]);
578: if (PetscIsInfOrNanReal(mfqP->Fres[0])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
579: mfqP->Fres[0]*=mfqP->Fres[0];
580: mfqP->minindex = 0;
581: minnorm = mfqP->Fres[0];
582: VecGetOwnershipRange(mfqP->Xhist[0],&low,&high);
583: for (i=1;i<mfqP->n+1;i++) {
584: VecCopy(tao->solution,mfqP->Xhist[i]);
585: if (i-1 >= low && i-1 < high) {
586: VecGetArray(mfqP->Xhist[i],&x);
587: x[i-1-low] += mfqP->delta;
588: VecRestoreArray(mfqP->Xhist[i],&x);
589: }
590: CHKMEMQ;
591: TaoComputeSeparableObjective(tao,mfqP->Xhist[i],mfqP->Fhist[i]);
592: VecNorm(mfqP->Fhist[i],NORM_2,&mfqP->Fres[i]);
593: if (PetscIsInfOrNanReal(mfqP->Fres[i])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
594: mfqP->Fres[i]*=mfqP->Fres[i];
595: if (mfqP->Fres[i] < minnorm) {
596: mfqP->minindex = i;
597: minnorm = mfqP->Fres[i];
598: }
599: }
601: VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);
602: VecCopy(mfqP->Fhist[mfqP->minindex],tao->sep_objective);
603: /* Gather mpi vecs to one big local vec */
605: /* Begin serial code */
607: /* Disp[i] = Xi-xmin, i=1,..,mfqP->minindex-1,mfqP->minindex+1,..,n */
608: /* Fdiff[i] = (Fi-Fmin)', i=1,..,mfqP->minindex-1,mfqP->minindex+1,..,n */
609: /* (Column oriented for blas calls) */
610: ii=0;
612: if (mfqP->size == 1) {
613: VecGetArray(mfqP->Xhist[mfqP->minindex],&xmint);
614: for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i];
615: VecRestoreArray(mfqP->Xhist[mfqP->minindex],&xmint);
616: VecGetArray(mfqP->Fhist[mfqP->minindex],&fmin);
617: for (i=0;i<mfqP->n+1;i++) {
618: if (i == mfqP->minindex) continue;
620: VecGetArray(mfqP->Xhist[i],&x);
621: for (j=0;j<mfqP->n;j++) {
622: mfqP->Disp[ii+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta;
623: }
624: VecRestoreArray(mfqP->Xhist[i],&x);
626: VecGetArray(mfqP->Fhist[i],&f);
627: for (j=0;j<mfqP->m;j++) {
628: mfqP->Fdiff[ii+mfqP->n*j] = f[j] - fmin[j];
629: }
630: VecRestoreArray(mfqP->Fhist[i],&f);
631: mfqP->model_indices[ii++] = i;
633: }
634: for (j=0;j<mfqP->m;j++) {
635: mfqP->C[j] = fmin[j];
636: }
637: VecRestoreArray(mfqP->Fhist[mfqP->minindex],&fmin);
638: } else {
639: VecScatterBegin(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);
640: VecScatterEnd(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);
641: VecGetArray(mfqP->localxmin,&xmint);
642: for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i];
643: VecRestoreArray(mfqP->localxmin,&mfqP->xmin);
645: VecScatterBegin(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);
646: VecScatterEnd(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);
647: VecGetArray(mfqP->localfmin,&fmin);
648: for (i=0;i<mfqP->n+1;i++) {
649: if (i == mfqP->minindex) continue;
651: mfqP->model_indices[ii++] = i;
652: VecScatterBegin(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);
653: VecScatterEnd(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);
654: VecGetArray(mfqP->localx,&x);
655: for (j=0;j<mfqP->n;j++) {
656: mfqP->Disp[i+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta;
657: }
658: VecRestoreArray(mfqP->localx,&x);
660: VecScatterBegin(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);
661: VecScatterEnd(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);
662: VecGetArray(mfqP->localf,&f);
663: for (j=0;j<mfqP->m;j++) {
664: mfqP->Fdiff[i*mfqP->n+j] = f[j] - fmin[j];
665: }
666: VecRestoreArray(mfqP->localf,&f);
667: }
668: for (j=0;j<mfqP->m;j++) {
669: mfqP->C[j] = fmin[j];
670: }
671: VecRestoreArray(mfqP->localfmin,&fmin);
672: }
674: /* Determine the initial quadratic models */
675: /* G = D(ModelIn,:) \ (F(ModelIn,1:m)-repmat(F(xkin,1:m),n,1)); */
676: /* D (nxn) Fdiff (nxm) => G (nxm) */
677: blasn2 = blasn;
678: PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&blasn,&blasm,mfqP->Disp,&blasnpmax,mfqP->iwork,mfqP->Fdiff,&blasn2,&info));
679: PetscInfo1(tao,"gesv returned %d\n",info);
681: cres = minnorm;
682: /* Gres = G*F(xkin,1:m)' G (nxm) Fk (m) */
683: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->C,&ione,&zero,mfqP->Gres,&ione));
685: /* Hres = G*G' */
686: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&blasn,&blasn,&blasm,&one,mfqP->Fdiff, &blasn,mfqP->Fdiff,&blasn,&zero,mfqP->Hres,&blasn));
688: valid = PETSC_TRUE;
690: VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);
691: VecAssemblyBegin(tao->gradient);
692: VecAssemblyEnd(tao->gradient);
693: VecNorm(tao->gradient,NORM_2,&gnorm);
694: gnorm *= mfqP->delta;
695: VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);
696: TaoMonitor(tao, iter, minnorm, gnorm, 0.0, step, &reason);
697: mfqP->nHist = mfqP->n+1;
698: mfqP->nmodelpoints = mfqP->n+1;
700: while (reason == TAO_CONTINUE_ITERATING) {
701: PetscReal gnm = 1e-4;
702: iter++;
703: /* Solve the subproblem min{Q(s): ||s|| <= delta} */
704: gqtwrap(tao,&gnm,&mdec);
705: /* Evaluate the function at the new point */
707: for (i=0;i<mfqP->n;i++) {
708: mfqP->work[i] = mfqP->Xsubproblem[i]*mfqP->delta + mfqP->xmin[i];
709: }
710: VecDuplicate(tao->solution,&mfqP->Xhist[mfqP->nHist]);
711: VecDuplicate(tao->sep_objective,&mfqP->Fhist[mfqP->nHist]);
712: VecSetValues(mfqP->Xhist[mfqP->nHist],mfqP->n,mfqP->indices,mfqP->work,INSERT_VALUES);
713: VecAssemblyBegin(mfqP->Xhist[mfqP->nHist]);
714: VecAssemblyEnd(mfqP->Xhist[mfqP->nHist]);
716: TaoComputeSeparableObjective(tao,mfqP->Xhist[mfqP->nHist],mfqP->Fhist[mfqP->nHist]);
717: VecNorm(mfqP->Fhist[mfqP->nHist],NORM_2,&mfqP->Fres[mfqP->nHist]);
718: if (PetscIsInfOrNanReal(mfqP->Fres[mfqP->nHist])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
719: mfqP->Fres[mfqP->nHist]*=mfqP->Fres[mfqP->nHist];
720: rho = (mfqP->Fres[mfqP->minindex] - mfqP->Fres[mfqP->nHist]) / mdec;
721: mfqP->nHist++;
723: /* Update the center */
724: if ((rho >= mfqP->eta1) || (rho > mfqP->eta0 && valid==PETSC_TRUE)) {
725: /* Update model to reflect new base point */
726: for (i=0;i<mfqP->n;i++) {
727: mfqP->work[i] = (mfqP->work[i] - mfqP->xmin[i])/mfqP->delta;
728: }
729: for (j=0;j<mfqP->m;j++) {
730: /* C(j) = C(j) + work*G(:,j) + .5*work*H(:,:,j)*work';
731: G(:,j) = G(:,j) + H(:,:,j)*work' */
732: for (k=0;k<mfqP->n;k++) {
733: mfqP->work2[k]=0.0;
734: for (l=0;l<mfqP->n;l++) {
735: mfqP->work2[k]+=mfqP->H[j + mfqP->m*(k + l*mfqP->n)]*mfqP->work[l];
736: }
737: }
738: for (i=0;i<mfqP->n;i++) {
739: mfqP->C[j]+=mfqP->work[i]*(mfqP->Fdiff[i + mfqP->n* j] + 0.5*mfqP->work2[i]);
740: mfqP->Fdiff[i+mfqP->n*j] +=mfqP-> work2[i];
741: }
742: }
743: /* Cres += work*Gres + .5*work*Hres*work';
744: Gres += Hres*work'; */
746: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasn,&one,mfqP->Hres,&blasn,mfqP->work,&ione,&zero,mfqP->work2,&ione));
747: for (i=0;j<mfqP->n;j++) {
748: cres += mfqP->work[i]*(mfqP->Gres[i] + 0.5*mfqP->work2[i]);
749: mfqP->Gres[i] += mfqP->work2[i];
750: }
751: mfqP->minindex = mfqP->nHist-1;
752: minnorm = mfqP->Fres[mfqP->minindex];
753: /* Change current center */
754: VecGetArray(mfqP->Xhist[mfqP->minindex],&xmint);
755: for (i=0;i<mfqP->n;i++) {
756: mfqP->xmin[i] = xmint[i];
757: }
758: VecRestoreArray(mfqP->Xhist[mfqP->minindex],&xmint);
759: }
761: /* Evaluate at a model-improving point if necessary */
762: if (valid == PETSC_FALSE) {
763: mfqP->q_is_I = 1;
764: mfqP->nmodelpoints = 0;
765: affpoints(mfqP,mfqP->xmin,mfqP->c1);
766: if (mfqP->nmodelpoints < mfqP->n) {
767: PetscInfo(tao,"Model not valid -- model-improving");
768: modelimprove(tao,mfqP,1);
769: }
770: }
772: /* Update the trust region radius */
773: deltaold = mfqP->delta;
774: normxsp = 0;
775: for (i=0;i<mfqP->n;i++) {
776: normxsp += mfqP->Xsubproblem[i]*mfqP->Xsubproblem[i];
777: }
778: normxsp = PetscSqrtReal(normxsp);
779: if (rho >= mfqP->eta1 && normxsp > 0.5*mfqP->delta) {
780: mfqP->delta = PetscMin(mfqP->delta*mfqP->gamma1,mfqP->deltamax);
781: } else if (valid == PETSC_TRUE) {
782: mfqP->delta = PetscMax(mfqP->delta*mfqP->gamma0,mfqP->deltamin);
783: }
785: /* Compute the next interpolation set */
786: mfqP->q_is_I = 1;
787: mfqP->nmodelpoints=0;
788: affpoints(mfqP,mfqP->xmin,mfqP->c1);
789: if (mfqP->nmodelpoints == mfqP->n) {
790: valid = PETSC_TRUE;
791: } else {
792: valid = PETSC_FALSE;
793: affpoints(mfqP,mfqP->xmin,mfqP->c2);
794: if (mfqP->n > mfqP->nmodelpoints) {
795: PetscInfo(tao,"Model not valid -- adding geometry points");
796: modelimprove(tao,mfqP,mfqP->n - mfqP->nmodelpoints);
797: }
798: }
799: for (i=mfqP->nmodelpoints;i>0;i--) {
800: mfqP->model_indices[i] = mfqP->model_indices[i-1];
801: }
802: mfqP->nmodelpoints++;
803: mfqP->model_indices[0] = mfqP->minindex;
804: morepoints(mfqP);
805: for (i=0;i<mfqP->nmodelpoints;i++) {
806: VecGetArray(mfqP->Xhist[mfqP->model_indices[i]],&x);
807: for (j=0;j<mfqP->n;j++) {
808: mfqP->Disp[i + mfqP->npmax*j] = (x[j] - mfqP->xmin[j]) / deltaold;
809: }
810: VecRestoreArray(mfqP->Xhist[mfqP->model_indices[i]],&x);
811: VecGetArray(mfqP->Fhist[mfqP->model_indices[i]],&f);
812: for (j=0;j<mfqP->m;j++) {
813: for (k=0;k<mfqP->n;k++) {
814: mfqP->work[k]=0.0;
815: for (l=0;l<mfqP->n;l++) {
816: mfqP->work[k] += mfqP->H[j + mfqP->m*(k + mfqP->n*l)] * mfqP->Disp[i + mfqP->npmax*l];
817: }
818: }
819: 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];
820: }
821: VecRestoreArray(mfqP->Fhist[mfqP->model_indices[i]],&f);
822: }
824: /* Update the quadratic model */
825: getquadpounders(mfqP);
826: VecGetArray(mfqP->Fhist[mfqP->minindex],&fmin);
827: PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasm,fmin,&ione,mfqP->C,&ione));
828: /* G = G*(delta/deltaold) + Gdel */
829: ratio = mfqP->delta/deltaold;
830: iblas = blasm*blasn;
831: PetscStackCallBLAS("BLASscal",BLASscal_(&iblas,&ratio,mfqP->Fdiff,&ione));
832: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&iblas,&one,mfqP->Gdel,&ione,mfqP->Fdiff,&ione));
833: /* H = H*(delta/deltaold) + Hdel */
834: iblas = blasm*blasn*blasn;
835: ratio *= ratio;
836: PetscStackCallBLAS("BLASscal",BLASscal_(&iblas,&ratio,mfqP->H,&ione));
837: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&iblas,&one,mfqP->Hdel,&ione,mfqP->H,&ione));
839: /* Get residuals */
840: cres = mfqP->Fres[mfqP->minindex];
841: /* Gres = G*F(xkin,1:m)' */
842: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->C,&ione,&zero,mfqP->Gres,&ione));
843: /* Hres = sum i=1..m {F(xkin,i)*H(:,:,i)} + G*G' */
844: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&blasn,&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->Fdiff,&blasn,&zero,mfqP->Hres,&blasn));
846: iblas = mfqP->n*mfqP->n;
848: for (j=0;j<mfqP->m;j++) {
849: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&iblas,&fmin[j],&mfqP->H[j],&blasm,mfqP->Hres,&ione));
850: }
852: /* Export solution and gradient residual to TAO */
853: VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);
854: VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);
855: VecAssemblyBegin(tao->gradient);
856: VecAssemblyEnd(tao->gradient);
857: VecNorm(tao->gradient,NORM_2,&gnorm);
858: gnorm *= mfqP->delta;
859: /* final criticality test */
860: TaoMonitor(tao, iter, minnorm, gnorm, 0.0, step, &reason);
861: }
862: return(0);
863: }
867: static PetscErrorCode TaoSetUp_POUNDERS(Tao tao)
868: {
869: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
870: PetscInt i;
871: IS isfloc,isfglob,isxloc,isxglob;
875: if (!tao->gradient) {VecDuplicate(tao->solution,&tao->gradient); }
876: if (!tao->stepdirection) {VecDuplicate(tao->solution,&tao->stepdirection); }
877: VecGetSize(tao->solution,&mfqP->n);
878: VecGetSize(tao->sep_objective,&mfqP->m);
879: mfqP->c1 = PetscSqrtReal((PetscReal)mfqP->n);
880: if (mfqP->npmax == PETSC_DEFAULT) {
881: mfqP->npmax = 2*mfqP->n + 1;
882: }
883: mfqP->npmax = PetscMin((mfqP->n+1)*(mfqP->n+2)/2,mfqP->npmax);
884: mfqP->npmax = PetscMax(mfqP->npmax, mfqP->n+2);
886: PetscMalloc1((tao->max_funcs+10),&mfqP->Xhist);
887: PetscMalloc1((tao->max_funcs+10),&mfqP->Fhist);
888: for (i=0;i<mfqP->n +1;i++) {
889: VecDuplicate(tao->solution,&mfqP->Xhist[i]);
890: VecDuplicate(tao->sep_objective,&mfqP->Fhist[i]);
891: }
892: VecDuplicate(tao->solution,&mfqP->workxvec);
893: mfqP->nHist = 0;
895: PetscMalloc1((tao->max_funcs+10),&mfqP->Fres);
896: PetscMalloc1(mfqP->npmax*mfqP->m,&mfqP->RES);
897: PetscMalloc1(mfqP->n,&mfqP->work);
898: PetscMalloc1(mfqP->n,&mfqP->work2);
899: PetscMalloc1(mfqP->n,&mfqP->work3);
900: PetscMalloc1(PetscMax(mfqP->m,mfqP->n+1),&mfqP->mwork);
901: PetscMalloc1((mfqP->npmax - mfqP->n - 1),&mfqP->omega);
902: PetscMalloc1((mfqP->n * (mfqP->n+1) / 2),&mfqP->beta);
903: PetscMalloc1((mfqP->n + 1) ,&mfqP->alpha);
905: PetscMalloc1(mfqP->n*mfqP->n*mfqP->m,&mfqP->H);
906: PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q);
907: PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q_tmp);
908: PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L);
909: PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_tmp);
910: PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_save);
911: PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->N);
912: PetscMalloc1(mfqP->npmax * (mfqP->n+1) ,&mfqP->M);
913: PetscMalloc1(mfqP->npmax * (mfqP->npmax - mfqP->n - 1) , &mfqP->Z);
914: PetscMalloc1(mfqP->npmax,&mfqP->tau);
915: PetscMalloc1(mfqP->npmax,&mfqP->tau_tmp);
916: mfqP->nmax = PetscMax(5*mfqP->npmax,mfqP->n*(mfqP->n+1)/2);
917: PetscMalloc1(mfqP->nmax,&mfqP->npmaxwork);
918: PetscMalloc1(mfqP->nmax,&mfqP->npmaxiwork);
919: PetscMalloc1(mfqP->n,&mfqP->xmin);
920: PetscMalloc1(mfqP->m,&mfqP->C);
921: PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Fdiff);
922: PetscMalloc1(mfqP->npmax*mfqP->n,&mfqP->Disp);
923: PetscMalloc1(mfqP->n,&mfqP->Gres);
924: PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Hres);
925: PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Gpoints);
926: PetscMalloc1(mfqP->npmax,&mfqP->model_indices);
927: PetscMalloc1(mfqP->n,&mfqP->Xsubproblem);
928: PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Gdel);
929: PetscMalloc1(mfqP->n*mfqP->n*mfqP->m, &mfqP->Hdel);
930: PetscMalloc1(PetscMax(mfqP->m,mfqP->n),&mfqP->indices);
931: PetscMalloc1(mfqP->n,&mfqP->iwork);
932: for (i=0;i<PetscMax(mfqP->m,mfqP->n);i++) {
933: mfqP->indices[i] = i;
934: }
935: MPI_Comm_size(((PetscObject)tao)->comm,&mfqP->size);
936: if (mfqP->size > 1) {
937: VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localx);
938: VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localxmin);
939: VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localf);
940: VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localfmin);
941: ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxloc);
942: ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxglob);
943: ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfloc);
944: ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfglob);
947: VecScatterCreate(tao->solution,isxglob,mfqP->localx,isxloc,&mfqP->scatterx);
948: VecScatterCreate(tao->sep_objective,isfglob,mfqP->localf,isfloc,&mfqP->scatterf);
950: ISDestroy(&isxloc);
951: ISDestroy(&isxglob);
952: ISDestroy(&isfloc);
953: ISDestroy(&isfglob);
954: }
956: if (!mfqP->usegqt) {
957: KSP ksp;
958: PC pc;
959: VecCreateSeqWithArray(PETSC_COMM_SELF,mfqP->n,mfqP->n,mfqP->Xsubproblem,&mfqP->subx);
960: VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->subxl);
961: VecDuplicate(mfqP->subxl,&mfqP->subb);
962: VecDuplicate(mfqP->subxl,&mfqP->subxu);
963: VecDuplicate(mfqP->subxl,&mfqP->subpdel);
964: VecDuplicate(mfqP->subxl,&mfqP->subndel);
965: TaoCreate(PETSC_COMM_SELF,&mfqP->subtao);
966: /* TaoSetType(mfqP->subtao,"tao_bqpip"); */
967: TaoSetType(mfqP->subtao,"tao_tron");
968: TaoSetOptionsPrefix(mfqP->subtao,"pounders_subsolver_");
969: TaoSetInitialVector(mfqP->subtao,mfqP->subx);
970: TaoSetObjectiveAndGradientRoutine(mfqP->subtao,pounders_fg,(void*)mfqP);
971: TaoSetMaximumIterations(mfqP->subtao,mfqP->gqt_maxits);
972: TaoGetKSP(mfqP->subtao,&ksp);
973: if (ksp) {
974: KSPGetPC(ksp,&pc);
975: PCSetType(pc,PCNONE);
976: }
977: TaoSetVariableBounds(mfqP->subtao,mfqP->subxl,mfqP->subxu);
978: TaoSetFromOptions(mfqP->subtao);
979: MatCreateSeqDense(PETSC_COMM_SELF,mfqP->n,mfqP->n,mfqP->Hres,&mfqP->subH);
980: TaoSetHessianRoutine(mfqP->subtao,mfqP->subH,mfqP->subH,pounders_h,(void*)mfqP);
981: }
982: return(0);
983: }
987: static PetscErrorCode TaoDestroy_POUNDERS(Tao tao)
988: {
989: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
990: PetscInt i;
994: if (!mfqP->usegqt) {
995: TaoDestroy(&mfqP->subtao);
996: VecDestroy(&mfqP->subx);
997: VecDestroy(&mfqP->subxl);
998: VecDestroy(&mfqP->subxu);
999: VecDestroy(&mfqP->subb);
1000: VecDestroy(&mfqP->subpdel);
1001: VecDestroy(&mfqP->subndel);
1002: MatDestroy(&mfqP->subH);
1003: }
1004: PetscFree(mfqP->Fres);
1005: PetscFree(mfqP->RES);
1006: PetscFree(mfqP->work);
1007: PetscFree(mfqP->work2);
1008: PetscFree(mfqP->work3);
1009: PetscFree(mfqP->mwork);
1010: PetscFree(mfqP->omega);
1011: PetscFree(mfqP->beta);
1012: PetscFree(mfqP->alpha);
1013: PetscFree(mfqP->H);
1014: PetscFree(mfqP->Q);
1015: PetscFree(mfqP->Q_tmp);
1016: PetscFree(mfqP->L);
1017: PetscFree(mfqP->L_tmp);
1018: PetscFree(mfqP->L_save);
1019: PetscFree(mfqP->N);
1020: PetscFree(mfqP->M);
1021: PetscFree(mfqP->Z);
1022: PetscFree(mfqP->tau);
1023: PetscFree(mfqP->tau_tmp);
1024: PetscFree(mfqP->npmaxwork);
1025: PetscFree(mfqP->npmaxiwork);
1026: PetscFree(mfqP->xmin);
1027: PetscFree(mfqP->C);
1028: PetscFree(mfqP->Fdiff);
1029: PetscFree(mfqP->Disp);
1030: PetscFree(mfqP->Gres);
1031: PetscFree(mfqP->Hres);
1032: PetscFree(mfqP->Gpoints);
1033: PetscFree(mfqP->model_indices);
1034: PetscFree(mfqP->Xsubproblem);
1035: PetscFree(mfqP->Gdel);
1036: PetscFree(mfqP->Hdel);
1037: PetscFree(mfqP->indices);
1038: PetscFree(mfqP->iwork);
1040: for (i=0;i<mfqP->nHist;i++) {
1041: VecDestroy(&mfqP->Xhist[i]);
1042: VecDestroy(&mfqP->Fhist[i]);
1043: }
1044: if (mfqP->workxvec) {
1045: VecDestroy(&mfqP->workxvec);
1046: }
1047: PetscFree(mfqP->Xhist);
1048: PetscFree(mfqP->Fhist);
1050: if (mfqP->size > 1) {
1051: VecDestroy(&mfqP->localx);
1052: VecDestroy(&mfqP->localxmin);
1053: VecDestroy(&mfqP->localf);
1054: VecDestroy(&mfqP->localfmin);
1055: }
1056: PetscFree(tao->data);
1057: return(0);
1058: }
1062: static PetscErrorCode TaoSetFromOptions_POUNDERS(Tao tao)
1063: {
1064: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
1068: PetscOptionsHead("POUNDERS method for least-squares optimization");
1069: PetscOptionsReal("-tao_pounders_delta","initial delta","",mfqP->delta,&mfqP->delta,0);
1070: mfqP->npmax = PETSC_DEFAULT;
1071: PetscOptionsInt("-tao_pounders_npmax","max number of points in model","",mfqP->npmax,&mfqP->npmax,0);
1072: mfqP->usegqt = PETSC_FALSE;
1073: PetscOptionsBool("-tao_pounders_gqt","use gqt algorithm for subproblem","",mfqP->usegqt,&mfqP->usegqt,0);
1074: PetscOptionsTail();
1075: return(0);
1076: }
1080: static PetscErrorCode TaoView_POUNDERS(Tao tao, PetscViewer viewer)
1081: {
1082: TAO_POUNDERS *mfqP = (TAO_POUNDERS *)tao->data;
1083: PetscBool isascii;
1087: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1088: if (isascii) {
1089: PetscViewerASCIIPushTab(viewer);
1090: if (mfqP->usegqt) {
1091: PetscViewerASCIIPrintf(viewer, "Subproblem solver: gqt\n");
1092: } else {
1093: PetscViewerASCIIPrintf(viewer, "Subproblem solver: tron\n");
1094: }
1095: PetscViewerASCIIPopTab(viewer);
1096: }
1097: return(0);
1098: }
1100: EXTERN_C_BEGIN
1103: PetscErrorCode TaoCreate_POUNDERS(Tao tao)
1104: {
1105: TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data;
1109: tao->ops->setup = TaoSetUp_POUNDERS;
1110: tao->ops->solve = TaoSolve_POUNDERS;
1111: tao->ops->view = TaoView_POUNDERS;
1112: tao->ops->setfromoptions = TaoSetFromOptions_POUNDERS;
1113: tao->ops->destroy = TaoDestroy_POUNDERS;
1115: PetscNewLog(tao,&mfqP);
1116: tao->data = (void*)mfqP;
1117: tao->max_it = 2000;
1118: tao->max_funcs = 4000;
1119: #if defined(PETSC_USE_REAL_SINGLE)
1120: tao->fatol = 1e-4;
1121: tao->frtol = 1e-4;
1122: mfqP->deltamin=1e-3;
1123: #else
1124: tao->fatol = 1e-8;
1125: tao->frtol = 1e-8;
1126: mfqP->deltamin=1e-6;
1127: #endif
1128: mfqP->delta = 0.1;
1129: mfqP->deltamax=1e3;
1130: mfqP->c2 = 100.0;
1131: mfqP->theta1=1e-5;
1132: mfqP->theta2=1e-4;
1133: mfqP->gamma0=0.5;
1134: mfqP->gamma1=2.0;
1135: mfqP->eta0 = 0.0;
1136: mfqP->eta1 = 0.1;
1137: mfqP->gqt_rtol = 0.001;
1138: mfqP->gqt_maxits = 50;
1139: mfqP->workxvec = 0;
1140: return(0);
1141: }
1142: EXTERN_C_END