Actual source code: pounders.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  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: }
 33: static PetscErrorCode pounders_feval(Tao tao, Vec x, Vec F, PetscReal *fsum)
 34: {
 36:   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;
 37:   PetscInt i,row,col;
 38:   PetscReal fr,fc;
 40:   TaoComputeSeparableObjective(tao,x,F);
 41:   if (tao->sep_weights_v) {
 42:     VecPointwiseMult(mfqP->workfvec,tao->sep_weights_v,F);
 43:     VecNorm(mfqP->workfvec,NORM_2,fsum);
 44:     *fsum*=(*fsum);
 45:   } else if (tao->sep_weights_w) {
 46:     *fsum=0;
 47:     for (i=0;i<tao->sep_weights_n;i++) {
 48:       row=tao->sep_weights_rows[i];
 49:       col=tao->sep_weights_cols[i];
 50:       VecGetValues(F,1,&row,&fr);
 51:       VecGetValues(F,1,&col,&fc);
 52:       *fsum += tao->sep_weights_w[i]*fc*fr;
 53:     }
 54:   } else {
 55:     VecNorm(F,NORM_2,fsum);
 56:     *fsum*=(*fsum);
 57:   }
 58:   if (PetscIsInfOrNanReal(*fsum)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
 59:   return(0);
 60: }

 64: PetscErrorCode gqtwrap(Tao tao,PetscReal *gnorm, PetscReal *qmin)
 65: {
 67: #if defined(PETSC_USE_REAL_SINGLE)
 68:   PetscReal      atol=1.0e-5;
 69: #else
 70:   PetscReal      atol=1.0e-10;
 71: #endif
 72:   PetscInt       info,its;
 73:   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;

 76:   if (! mfqP->usegqt) {
 77:     PetscReal maxval;
 78:     PetscInt  i,j;

 80:     VecSetValues(mfqP->subb,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);
 81:     VecAssemblyBegin(mfqP->subb);
 82:     VecAssemblyEnd(mfqP->subb);

 84:     VecSet(mfqP->subx,0.0);

 86:     VecSet(mfqP->subndel,-mfqP->delta);
 87:     VecSet(mfqP->subpdel,mfqP->delta);

 89:     for (i=0;i<mfqP->n;i++) {
 90:       for (j=i;j<mfqP->n;j++) {
 91:         mfqP->Hres[j+mfqP->n*i] = mfqP->Hres[mfqP->n*j+i];
 92:       }
 93:     }
 94:     MatSetValues(mfqP->subH,mfqP->n,mfqP->indices,mfqP->n,mfqP->indices,mfqP->Hres,INSERT_VALUES);
 95:     MatAssemblyBegin(mfqP->subH,MAT_FINAL_ASSEMBLY);
 96:     MatAssemblyEnd(mfqP->subH,MAT_FINAL_ASSEMBLY);

 98:     TaoResetStatistics(mfqP->subtao);
 99:     TaoSetTolerances(mfqP->subtao,*gnorm,*gnorm,PETSC_DEFAULT);
100:     /* enforce bound constraints -- experimental */
101:     if (tao->XU && tao->XL) {
102:       VecCopy(tao->XU,mfqP->subxu);
103:       VecAXPY(mfqP->subxu,-1.0,tao->solution);
104:       VecScale(mfqP->subxu,1.0/mfqP->delta);
105:       VecCopy(tao->XL,mfqP->subxl);
106:       VecAXPY(mfqP->subxl,-1.0,tao->solution);
107:       VecScale(mfqP->subxl,1.0/mfqP->delta);

109:       VecPointwiseMin(mfqP->subxu,mfqP->subxu,mfqP->subpdel);
110:       VecPointwiseMax(mfqP->subxl,mfqP->subxl,mfqP->subndel);
111:     } else {
112:       VecCopy(mfqP->subpdel,mfqP->subxu);
113:       VecCopy(mfqP->subndel,mfqP->subxl);
114:     }
115:     /* Make sure xu > xl */
116:     VecCopy(mfqP->subxl,mfqP->subpdel);
117:     VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);
118:     VecMax(mfqP->subpdel,NULL,&maxval);
119:     if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"upper bound < lower bound in subproblem");
120:     /* Make sure xu > tao->solution > xl */
121:     VecCopy(mfqP->subxl,mfqP->subpdel);
122:     VecAXPY(mfqP->subpdel,-1.0,mfqP->subx);
123:     VecMax(mfqP->subpdel,NULL,&maxval);
124:     if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"initial guess < lower bound in subproblem");

126:     VecCopy(mfqP->subx,mfqP->subpdel);
127:     VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);
128:     VecMax(mfqP->subpdel,NULL,&maxval);
129:     if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"initial guess > upper bound in subproblem");

131:     TaoSolve(mfqP->subtao);
132:     TaoGetSolutionStatus(mfqP->subtao,NULL,qmin,NULL,NULL,NULL,NULL);

134:     /* test bounds post-solution*/
135:     VecCopy(mfqP->subxl,mfqP->subpdel);
136:     VecAXPY(mfqP->subpdel,-1.0,mfqP->subx);
137:     VecMax(mfqP->subpdel,NULL,&maxval);
138:     if (maxval > 1e-5) {
139:       PetscInfo(tao,"subproblem solution < lower bound\n");
140:       tao->reason = TAO_DIVERGED_TR_REDUCTION;
141:     }

143:     VecCopy(mfqP->subx,mfqP->subpdel);
144:     VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);
145:     VecMax(mfqP->subpdel,NULL,&maxval);
146:     if (maxval > 1e-5) {
147:       PetscInfo(tao,"subproblem solution > upper bound\n");
148:       tao->reason = TAO_DIVERGED_TR_REDUCTION;
149:     }
150:   } else {
151:     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);
152:   }
153:   *qmin *= -1;
154:   return(0);
155: }

159: static PetscErrorCode pounders_update_res(Tao tao)
160: {
161:   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;
162:   PetscInt i,row,col;
163:   PetscBLASInt blasn=mfqP->n,blasn2=blasn*blasn,blasm=mfqP->m,ione=1;
164:   PetscReal zero=0.0,one=1.0,wii,factor;

168:   for (i=0;i<mfqP->n;i++) {
169:     mfqP->Gres[i]=0;
170:   }
171:   for (i=0;i<mfqP->n*mfqP->n;i++) {
172:     mfqP->Hres[i]=0;
173:   }

175:   /* Compute Gres= sum_ij[wij * (cjgi + cigj)] */
176:   if (tao->sep_weights_v) {
177:     /* Vector(diagonal) weights: gres = sum_i(wii*ci*gi) */
178:     for (i=0;i<mfqP->m;i++) {
179:       VecGetValues(tao->sep_weights_v,1,&i,&factor);
180:       factor=factor*mfqP->C[i];
181:       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasn,&factor,&mfqP->Fdiff[blasn*i],&ione,mfqP->Gres,&ione));
182:     }
183:   } else if (tao->sep_weights_n) {
184:     /* general case: .5 * Gres= sum_ij[wij * (cjgi + cigj)] */
185:     for (i=0;i<tao->sep_weights_n;i++) {
186:       row=tao->sep_weights_rows[i];
187:       col=tao->sep_weights_cols[i];

189:       factor = tao->sep_weights_w[i]*mfqP->C[col]/2.0;
190:       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasn,&factor,&mfqP->Fdiff[blasn*row],&ione,mfqP->Gres,&ione));
191:       factor = tao->sep_weights_w[i]*mfqP->C[row]/2.0;
192:       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasn,&factor,&mfqP->Fdiff[blasn*col],&ione,mfqP->Gres,&ione));
193:     }
194:   } else {
195:     /* default: Gres= sum_i[cigi] = G*c' */
196:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->C,&ione,&zero,mfqP->Gres,&ione));
197:   }

199:   /* compute Hres = sum_ij [wij * (*ci*Hj + cj*Hi + gi gj' + gj gi') ] */
200:   if (tao->sep_weights_v) {
201:     /* vector(diagonal weights) Hres = sum_i(wii*(ci*Hi + gi * gi' )*/
202:     for (i=0;i<mfqP->m;i++) {
203:       VecGetValues(tao->sep_weights_v,1,&i,&wii);
204:       if (tao->niter>1) {
205:         factor=wii*mfqP->C[i];
206:         /* add wii * ci * Hi */
207:         PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasn2,&factor,&mfqP->H[i],&blasn2,mfqP->Hres,&ione));
208:       }
209:       /* add wii * gi * gi' */
210:       PetscStackCallBLAS("BLASgemm_",BLASgemm_("N","T",&blasn,&blasn,&ione,&wii,&mfqP->Fdiff[blasn*i],&blasn,&mfqP->Fdiff[blasn*i],&blasn,&one,mfqP->Hres,&blasn));
211:     }
212:   } else if (tao->sep_weights_w) {
213:     /* .5 * sum_ij [wij * (*ci*Hj + cj*Hi + gi gj' + gj gi') ] */
214:     for (i=0;i<tao->sep_weights_n;i++) {
215:       row=tao->sep_weights_rows[i];
216:       col=tao->sep_weights_cols[i];
217:       factor=tao->sep_weights_w[i]/2.0;
218:       /* add wij * gi gj' + wij * gj gi' */
219:       PetscStackCallBLAS("BLASgemm_",BLASgemm_("N","T",&blasn,&blasn,&ione,&factor,&mfqP->Fdiff[blasn*row],&blasn,&mfqP->Fdiff[blasn*col],&blasn,&one,mfqP->Hres,&blasn));
220:       PetscStackCallBLAS("BLASgemm_",BLASgemm_("N","T",&blasn,&blasn,&ione,&factor,&mfqP->Fdiff[blasn*col],&blasn,&mfqP->Fdiff[blasn*row],&blasn,&one,mfqP->Hres,&blasn));
221:     }
222:     if (tao->niter > 1) {
223:       for (i=0;i<tao->sep_weights_n;i++) {
224:         row=tao->sep_weights_rows[i];
225:         col=tao->sep_weights_cols[i];

227:         /* add  wij*cj*Hi */
228:         factor = tao->sep_weights_w[i]*mfqP->C[col]/2.0;
229:         PetscStackCallBLAS("BLASaxpy_",BLASaxpy_(&blasn2,&factor,&mfqP->H[row],&blasn2,mfqP->Hres,&ione));

231:         /* add wij*ci*Hj */
232:         factor = tao->sep_weights_w[i]*mfqP->C[row]/2.0;
233:         PetscStackCallBLAS("BLASaxpy_",BLASaxpy_(&blasn2,&factor,&mfqP->H[col],&blasn2,mfqP->Hres,&ione));
234:       }
235:     }
236:   } else {
237:     /*  Hres = G*G' + 0.5 sum {F(xkin,i)*H(:,:,i)}  */
238:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&blasn,&blasn,&blasm,&one,mfqP->Fdiff, &blasn,mfqP->Fdiff, &blasn,&zero,mfqP->Hres,&blasn));

240:     /* sum(F(xkin,i)*H(:,:,i)) */
241:     if (tao->niter>1) {
242:       for (i=0;i<mfqP->m;i++) {
243:         factor = mfqP->C[i];
244:         PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasn2,&factor,&mfqP->H[i],&blasn2,mfqP->Hres,&ione));
245:       }
246:     }
247:   }

249:   return(0);
250: }
253: PetscErrorCode phi2eval(PetscReal *x, PetscInt n, PetscReal *phi)
254: {
255: /* 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] */
256:   PetscInt  i,j,k;
257:   PetscReal sqrt2 = PetscSqrtReal(2.0);

260:   j=0;
261:   for (i=0;i<n;i++) {
262:     phi[j] = 0.5 * x[i]*x[i];
263:     j++;
264:     for (k=i+1;k<n;k++) {
265:       phi[j]  = x[i]*x[k]/sqrt2;
266:       j++;
267:     }
268:   }
269:   return(0);
270: }

274: PetscErrorCode getquadpounders(TAO_POUNDERS *mfqP)
275: {
276: /* Computes the parameters of the quadratic Q(x) = c + g'*x + 0.5*x*G*x'
277:    that satisfies the interpolation conditions Q(X[:,j]) = f(j)
278:    for j=1,...,m and with a Hessian matrix of least Frobenius norm */

280:     /* NB --we are ignoring c */
281:   PetscInt     i,j,k,num,np = mfqP->nmodelpoints;
282:   PetscReal    one = 1.0,zero=0.0,negone=-1.0;
283:   PetscBLASInt blasnpmax = mfqP->npmax;
284:   PetscBLASInt blasnplus1 = mfqP->n+1;
285:   PetscBLASInt blasnp = np;
286:   PetscBLASInt blasint = mfqP->n*(mfqP->n+1) / 2;
287:   PetscBLASInt blasint2 = np - mfqP->n-1;
288:   PetscBLASInt info,ione=1;
289:   PetscReal    sqrt2 = PetscSqrtReal(2.0);

292:   for (i=0;i<mfqP->n*mfqP->m;i++) {
293:     mfqP->Gdel[i] = 0;
294:   }
295:   for (i=0;i<mfqP->n*mfqP->n*mfqP->m;i++) {
296:     mfqP->Hdel[i] = 0;
297:   }

299:     /* factor M */
300:   PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&blasnplus1,&blasnp,mfqP->M,&blasnplus1,mfqP->npmaxiwork,&info));
301:   if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine getrf returned with value %d\n",info);

303:   if (np == mfqP->n+1) {
304:     for (i=0;i<mfqP->npmax-mfqP->n-1;i++) {
305:       mfqP->omega[i]=0.0;
306:     }
307:     for (i=0;i<mfqP->n*(mfqP->n+1)/2;i++) {
308:       mfqP->beta[i]=0.0;
309:     }
310:   } else {
311:     /* Let Ltmp = (L'*L) */
312:     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));

314:     /* factor Ltmp */
315:     PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&blasint2,mfqP->L_tmp,&blasint,&info));
316:     if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine potrf returned with value %d\n",info);
317:   }

319:   for (k=0;k<mfqP->m;k++) {
320:     if (np != mfqP->n+1) {
321:       /* Solve L'*L*Omega = Z' * RESk*/
322:       PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&blasnp,&blasint2,&one,mfqP->Z,&blasnpmax,&mfqP->RES[mfqP->npmax*k],&ione,&zero,mfqP->omega,&ione));
323:       PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&blasint2,&ione,mfqP->L_tmp,&blasint,mfqP->omega,&blasint2,&info));
324:       if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine potrs returned with value %d\n",info);

326:       /* Beta = L*Omega */
327:       PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasint,&blasint2,&one,&mfqP->L[(mfqP->n+1)*blasint],&blasint,mfqP->omega,&ione,&zero,mfqP->beta,&ione));
328:     }

330:     /* solve M'*Alpha = RESk - N'*Beta */
331:     PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&blasint,&blasnp,&negone,mfqP->N,&blasint,mfqP->beta,&ione,&one,&mfqP->RES[mfqP->npmax*k],&ione));
332:     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&blasnplus1,&ione,mfqP->M,&blasnplus1,mfqP->npmaxiwork,&mfqP->RES[mfqP->npmax*k],&blasnplus1,&info));
333:     if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine getrs returned with value %d\n",info);

335:     /* Gdel(:,k) = Alpha(2:n+1) */
336:     for (i=0;i<mfqP->n;i++) {
337:       mfqP->Gdel[i + mfqP->n*k] = mfqP->RES[mfqP->npmax*k + i+1];
338:     }

340:     /* Set Hdels */
341:     num=0;
342:     for (i=0;i<mfqP->n;i++) {
343:       /* H[i,i,k] = Beta(num) */
344:       mfqP->Hdel[(i*mfqP->n + i)*mfqP->m + k] = mfqP->beta[num];
345:       num++;
346:       for (j=i+1;j<mfqP->n;j++) {
347:         /* H[i,j,k] = H[j,i,k] = Beta(num)/sqrt(2) */
348:         mfqP->Hdel[(j*mfqP->n + i)*mfqP->m + k] = mfqP->beta[num]/sqrt2;
349:         mfqP->Hdel[(i*mfqP->n + j)*mfqP->m + k] = mfqP->beta[num]/sqrt2;
350:         num++;
351:       }
352:     }
353:   }
354:   return(0);
355: }

359: PetscErrorCode morepoints(TAO_POUNDERS *mfqP)
360: {
361:   /* Assumes mfqP->model_indices[0]  is minimum index
362:    Finishes adding points to mfqP->model_indices (up to npmax)
363:    Computes L,Z,M,N
364:    np is actual number of points in model (should equal npmax?) */
365:   PetscInt        point,i,j,offset;
366:   PetscInt        reject;
367:   PetscBLASInt    blasn=mfqP->n,blasnpmax=mfqP->npmax,blasnplus1=mfqP->n+1,info,blasnmax=mfqP->nmax,blasint,blasint2,blasnp,blasmaxmn;
368:   const PetscReal *x;
369:   PetscReal       normd;
370:   PetscErrorCode  ierr;

373:   /* Initialize M,N */
374:   for (i=0;i<mfqP->n+1;i++) {
375:     VecGetArrayRead(mfqP->Xhist[mfqP->model_indices[i]],&x);
376:     mfqP->M[(mfqP->n+1)*i] = 1.0;
377:     for (j=0;j<mfqP->n;j++) {
378:       mfqP->M[j+1+((mfqP->n+1)*i)] = (x[j]  - mfqP->xmin[j]) / mfqP->delta;
379:     }
380:     VecRestoreArrayRead(mfqP->Xhist[mfqP->model_indices[i]],&x);
381:     phi2eval(&mfqP->M[1+((mfqP->n+1)*i)],mfqP->n,&mfqP->N[mfqP->n*(mfqP->n+1)/2 * i]);
382:   }

384:   /* Now we add points until we have npmax starting with the most recent ones */
385:   point = mfqP->nHist-1;
386:   mfqP->nmodelpoints = mfqP->n+1;
387:   while (mfqP->nmodelpoints < mfqP->npmax && point>=0) {
388:     /* Reject any points already in the model */
389:     reject = 0;
390:     for (j=0;j<mfqP->n+1;j++) {
391:       if (point == mfqP->model_indices[j]) {
392:         reject = 1;
393:         break;
394:       }
395:     }

397:     /* Reject if norm(d) >c2 */
398:     if (!reject) {
399:       VecCopy(mfqP->Xhist[point],mfqP->workxvec);
400:       VecAXPY(mfqP->workxvec,-1.0,mfqP->Xhist[mfqP->minindex]);
401:       VecNorm(mfqP->workxvec,NORM_2,&normd);
402:       normd /= mfqP->delta;
403:       if (normd > mfqP->c2) {
404:         reject =1;
405:       }
406:     }
407:     if (reject){
408:       point--;
409:       continue;
410:     }

412:     VecGetArrayRead(mfqP->Xhist[point],&x);
413:     mfqP->M[(mfqP->n+1)*mfqP->nmodelpoints] = 1.0;
414:     for (j=0;j<mfqP->n;j++) {
415:       mfqP->M[j+1+((mfqP->n+1)*mfqP->nmodelpoints)] = (x[j]  - mfqP->xmin[j]) / mfqP->delta;
416:     }
417:     VecRestoreArrayRead(mfqP->Xhist[point],&x);
418:     phi2eval(&mfqP->M[1+(mfqP->n+1)*mfqP->nmodelpoints],mfqP->n,&mfqP->N[mfqP->n*(mfqP->n+1)/2 * (mfqP->nmodelpoints)]);

420:     /* Update QR factorization */
421:     /* Copy M' to Q_tmp */
422:     for (i=0;i<mfqP->n+1;i++) {
423:       for (j=0;j<mfqP->npmax;j++) {
424:         mfqP->Q_tmp[j+mfqP->npmax*i] = mfqP->M[i+(mfqP->n+1)*j];
425:       }
426:     }
427:     blasnp = mfqP->nmodelpoints+1;
428:     /* Q_tmp,R = qr(M') */
429:     blasmaxmn=PetscMax(mfqP->m,mfqP->n+1);
430:     PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&blasnp,&blasnplus1,mfqP->Q_tmp,&blasnpmax,mfqP->tau_tmp,mfqP->mwork,&blasmaxmn,&info));
431:     if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine geqrf returned with value %d\n",info);

433:     /* Reject if min(svd(N*Q(:,n+2:np+1)) <= theta2 */
434:     /* L = N*Qtmp */
435:     blasint2 = mfqP->n * (mfqP->n+1) / 2;
436:     /* Copy N to L_tmp */
437:     for (i=0;i<mfqP->n*(mfqP->n+1)/2 * mfqP->npmax;i++) {
438:       mfqP->L_tmp[i]= mfqP->N[i];
439:     }
440:     /* Copy L_save to L_tmp */

442:     /* L_tmp = N*Qtmp' */
443:     PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("R","N",&blasint2,&blasnp,&blasnplus1,mfqP->Q_tmp,&blasnpmax,mfqP->tau_tmp,mfqP->L_tmp,&blasint2,mfqP->npmaxwork,&blasnmax,&info));
444:     if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine ormqr returned with value %d\n",info);

446:     /* Copy L_tmp to L_save */
447:     for (i=0;i<mfqP->npmax * mfqP->n*(mfqP->n+1)/2;i++) {
448:       mfqP->L_save[i] = mfqP->L_tmp[i];
449:     }

451:     /* Get svd for L_tmp(:,n+2:np+1) (L_tmp is modified in process) */
452:     blasint = mfqP->nmodelpoints - mfqP->n;
453:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
454:     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));
455:     PetscFPTrapPop();
456:     if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine gesvd returned with value %d\n",info);

458:     if (mfqP->beta[PetscMin(blasint,blasint2)-1] > mfqP->theta2) {
459:       /* accept point */
460:       mfqP->model_indices[mfqP->nmodelpoints] = point;
461:       /* Copy Q_tmp to Q */
462:       for (i=0;i<mfqP->npmax* mfqP->npmax;i++) {
463:         mfqP->Q[i] = mfqP->Q_tmp[i];
464:       }
465:       for (i=0;i<mfqP->npmax;i++){
466:         mfqP->tau[i] = mfqP->tau_tmp[i];
467:       }
468:       mfqP->nmodelpoints++;
469:       blasnp = mfqP->nmodelpoints;

471:       /* Copy L_save to L */
472:       for (i=0;i<mfqP->npmax * mfqP->n*(mfqP->n+1)/2;i++) {
473:         mfqP->L[i] = mfqP->L_save[i];
474:       }
475:     }
476:     point--;
477:   }

479:   blasnp = mfqP->nmodelpoints;
480:   /* Copy Q(:,n+2:np) to Z */
481:   /* First set Q_tmp to I */
482:   for (i=0;i<mfqP->npmax*mfqP->npmax;i++) {
483:     mfqP->Q_tmp[i] = 0.0;
484:   }
485:   for (i=0;i<mfqP->npmax;i++) {
486:     mfqP->Q_tmp[i + mfqP->npmax*i] = 1.0;
487:   }

489:   /* Q_tmp = I * Q */
490:   PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("R","N",&blasnp,&blasnp,&blasnplus1,mfqP->Q,&blasnpmax,mfqP->tau,mfqP->Q_tmp,&blasnpmax,mfqP->npmaxwork,&blasnmax,&info));
491:   if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine ormqr returned with value %d\n",info);

493:   /* Copy Q_tmp(:,n+2:np) to Z) */
494:   offset = mfqP->npmax * (mfqP->n+1);
495:   for (i=offset;i<mfqP->npmax*mfqP->npmax;i++) {
496:     mfqP->Z[i-offset] = mfqP->Q_tmp[i];
497:   }

499:   if (mfqP->nmodelpoints == mfqP->n + 1) {
500:     /* Set L to I_{n+1} */
501:     for (i=0;i<mfqP->npmax * mfqP->n*(mfqP->n+1)/2;i++) {
502:       mfqP->L[i] = 0.0;
503:     }
504:     for (i=0;i<mfqP->n;i++) {
505:       mfqP->L[(mfqP->n*(mfqP->n+1)/2)*i + i] = 1.0;
506:     }
507:   }
508:   return(0);
509: }

513: /* Only call from modelimprove, addpoint() needs ->Q_tmp and ->work to be set */
514: PetscErrorCode addpoint(Tao tao, TAO_POUNDERS *mfqP, PetscInt index)
515: {

519:   /* Create new vector in history: X[newidx] = X[mfqP->index] + delta*X[index]*/
520:   VecDuplicate(mfqP->Xhist[0],&mfqP->Xhist[mfqP->nHist]);
521:   VecSetValues(mfqP->Xhist[mfqP->nHist],mfqP->n,mfqP->indices,&mfqP->Q_tmp[index*mfqP->npmax],INSERT_VALUES);
522:   VecAssemblyBegin(mfqP->Xhist[mfqP->nHist]);
523:   VecAssemblyEnd(mfqP->Xhist[mfqP->nHist]);
524:   VecAYPX(mfqP->Xhist[mfqP->nHist],mfqP->delta,mfqP->Xhist[mfqP->minindex]);

526:   /* Project into feasible region */
527:   if (tao->XU && tao->XL) {
528:     VecMedian(mfqP->Xhist[mfqP->nHist], tao->XL, tao->XU, mfqP->Xhist[mfqP->nHist]);
529:   }

531:   /* Compute value of new vector */
532:   VecDuplicate(mfqP->Fhist[0],&mfqP->Fhist[mfqP->nHist]);
533:   CHKMEMQ;
534:   pounders_feval(tao,mfqP->Xhist[mfqP->nHist],mfqP->Fhist[mfqP->nHist],&mfqP->Fres[mfqP->nHist]);

536:   /* Add new vector to model */
537:   mfqP->model_indices[mfqP->nmodelpoints] = mfqP->nHist;
538:   mfqP->nmodelpoints++;
539:   mfqP->nHist++;
540:   return(0);
541: }

545: PetscErrorCode modelimprove(Tao tao, TAO_POUNDERS *mfqP, PetscInt addallpoints)
546: {
547:   /* modeld = Q(:,np+1:n)' */
549:   PetscInt       i,j,minindex=0;
550:   PetscReal      dp,half=0.5,one=1.0,minvalue=PETSC_INFINITY;
551:   PetscBLASInt   blasn=mfqP->n,  blasnpmax = mfqP->npmax, blask,info;
552:   PetscBLASInt   blas1=1,blasnmax = mfqP->nmax;

554:   blask = mfqP->nmodelpoints;
555:   /* Qtmp = I(n x n) */
556:   for (i=0;i<mfqP->n;i++) {
557:     for (j=0;j<mfqP->n;j++) {
558:       mfqP->Q_tmp[i + mfqP->npmax*j] = 0.0;
559:     }
560:   }
561:   for (j=0;j<mfqP->n;j++) {
562:     mfqP->Q_tmp[j + mfqP->npmax*j] = 1.0;
563:   }

565:   /* Qtmp = Q * I */
566:   PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("R","N",&blasn,&blasn,&blask,mfqP->Q,&blasnpmax,mfqP->tau, mfqP->Q_tmp, &blasnpmax, mfqP->npmaxwork,&blasnmax, &info));

568:   for (i=mfqP->nmodelpoints;i<mfqP->n;i++) {
569:     dp = BLASdot_(&blasn,&mfqP->Q_tmp[i*mfqP->npmax],&blas1,mfqP->Gres,&blas1);
570:     if (dp>0.0) { /* Model says use the other direction! */
571:       for (j=0;j<mfqP->n;j++) {
572:         mfqP->Q_tmp[i*mfqP->npmax+j] *= -1;
573:       }
574:     }
575:     /* mfqP->work[i] = Cres+Modeld(i,:)*(Gres+.5*Hres*Modeld(i,:)') */
576:     for (j=0;j<mfqP->n;j++) {
577:       mfqP->work2[j] = mfqP->Gres[j];
578:     }
579:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasn,&half,mfqP->Hres,&blasn,&mfqP->Q_tmp[i*mfqP->npmax], &blas1, &one, mfqP->work2,&blas1));
580:     mfqP->work[i] = BLASdot_(&blasn,&mfqP->Q_tmp[i*mfqP->npmax],&blas1,mfqP->work2,&blas1);
581:     if (i==mfqP->nmodelpoints || mfqP->work[i] < minvalue) {
582:       minindex=i;
583:       minvalue = mfqP->work[i];
584:     }
585:     if (addallpoints != 0) {
586:       addpoint(tao,mfqP,i);
587:     }
588:   }
589:   if (!addallpoints) {
590:     addpoint(tao,mfqP,minindex);
591:   }
592:   return(0);
593: }


598: PetscErrorCode affpoints(TAO_POUNDERS *mfqP, PetscReal *xmin,PetscReal c)
599: {
600:   PetscInt        i,j;
601:   PetscBLASInt    blasm=mfqP->m,blasj,blask,blasn=mfqP->n,ione=1,info;
602:   PetscBLASInt    blasnpmax = mfqP->npmax,blasmaxmn;
603:   PetscReal       proj,normd;
604:   const PetscReal *x;
605:   PetscErrorCode  ierr;

608:   for (i=mfqP->nHist-1;i>=0;i--) {
609:     VecGetArrayRead(mfqP->Xhist[i],&x);
610:     for (j=0;j<mfqP->n;j++) {
611:       mfqP->work[j] = (x[j] - xmin[j])/mfqP->delta;
612:     }
613:     VecRestoreArrayRead(mfqP->Xhist[i],&x);
614:     PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasn,mfqP->work,&ione,mfqP->work2,&ione));
615:     normd = BLASnrm2_(&blasn,mfqP->work,&ione);
616:     if (normd <= c*c) {
617:       blasj=PetscMax((mfqP->n - mfqP->nmodelpoints),0);
618:       if (!mfqP->q_is_I) {
619:         /* project D onto null */
620:         blask=(mfqP->nmodelpoints);
621:         PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("R","N",&ione,&blasn,&blask,mfqP->Q,&blasnpmax,mfqP->tau,mfqP->work2,&ione,mfqP->mwork,&blasm,&info));
622:         if (info < 0) SETERRQ1(PETSC_COMM_SELF,1,"ormqr returned value %d\n",info);
623:       }
624:       proj = BLASnrm2_(&blasj,&mfqP->work2[mfqP->nmodelpoints],&ione);

626:       if (proj >= mfqP->theta1) { /* add this index to model */
627:         mfqP->model_indices[mfqP->nmodelpoints]=i;
628:         mfqP->nmodelpoints++;
629:         PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasn,mfqP->work,&ione,&mfqP->Q_tmp[mfqP->npmax*(mfqP->nmodelpoints-1)],&ione));
630:         blask=mfqP->npmax*(mfqP->nmodelpoints);
631:         PetscStackCallBLAS("BLAScopy",BLAScopy_(&blask,mfqP->Q_tmp,&ione,mfqP->Q,&ione));
632:         blask = mfqP->nmodelpoints;
633:         blasmaxmn = PetscMax(mfqP->m,mfqP->n);
634:         PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&blasn,&blask,mfqP->Q,&blasnpmax,mfqP->tau,mfqP->mwork,&blasmaxmn,&info));
635:         if (info < 0) SETERRQ1(PETSC_COMM_SELF,1,"geqrf returned value %d\n",info);
636:         mfqP->q_is_I = 0;
637:       }
638:       if (mfqP->nmodelpoints == mfqP->n)  {
639:         break;
640:       }
641:     }
642:   }

644:   return(0);
645: }

649: static PetscErrorCode TaoSolve_POUNDERS(Tao tao)
650: {
651:   TAO_POUNDERS       *mfqP = (TAO_POUNDERS *)tao->data;
652:   PetscInt           i,ii,j,k,l;
653:   PetscReal          step=1.0;
654:   TaoConvergedReason reason = TAO_CONTINUE_ITERATING;
655:   PetscInt           low,high;
656:   PetscReal          minnorm;
657:   PetscReal          *x,*f;
658:   const PetscReal    *xmint,*fmin;
659:   PetscReal          cres,deltaold;
660:   PetscReal          gnorm;
661:   PetscBLASInt       info,ione=1,iblas;
662:   PetscBool          valid,same;
663:   PetscReal          mdec, rho, normxsp;
664:   PetscReal          one=1.0,zero=0.0,ratio;
665:   PetscBLASInt       blasm,blasn,blasncopy,blasnpmax;
666:   PetscErrorCode     ierr;
667:   static PetscBool   set = PETSC_FALSE;

669:   /* n = # of parameters
670:      m = dimension (components) of function  */
672:   PetscCitationsRegister("@article{UNEDF0,\n"
673:                                 "title = {Nuclear energy density optimization},\n"
674:                                 "author = {Kortelainen, M.  and Lesinski, T.  and Mor\'e, J.  and Nazarewicz, W.\n"
675:                                 "          and Sarich, J.  and Schunck, N.  and Stoitsov, M. V. and Wild, S. },\n"
676:                                 "journal = {Phys. Rev. C},\n"
677:                                 "volume = {82},\n"
678:                                 "number = {2},\n"
679:                                 "pages = {024313},\n"
680:                                 "numpages = {18},\n"
681:                                 "year = {2010},\n"
682:                                 "month = {Aug},\n"
683:                                 "doi = {10.1103/PhysRevC.82.024313}\n}\n",&set);
684:   tao->niter=0;
685:   if (tao->XL && tao->XU) {
686:     /* Check x0 <= XU */
687:     PetscReal val;
688:     VecCopy(tao->solution,mfqP->Xhist[0]);
689:     VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);
690:     VecMax(mfqP->Xhist[0],NULL,&val);
691:     if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 > upper bound");

693:     /* Check x0 >= xl */
694:     VecCopy(tao->XL,mfqP->Xhist[0]);
695:     VecAXPY(mfqP->Xhist[0],-1.0,tao->solution);
696:     VecMax(mfqP->Xhist[0],NULL,&val);
697:     if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 < lower bound");

699:     /* Check x0 + delta < XU  -- should be able to get around this eventually */

701:     VecSet(mfqP->Xhist[0],mfqP->delta);
702:     VecAXPY(mfqP->Xhist[0],1.0,tao->solution);
703:     VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);
704:     VecMax(mfqP->Xhist[0],NULL,&val);
705:     if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 + delta > upper bound");
706:   }

708:   CHKMEMQ;
709:   blasm = mfqP->m; blasn=mfqP->n; blasnpmax = mfqP->npmax;
710:   for (i=0;i<mfqP->n*mfqP->n*mfqP->m;i++) mfqP->H[i]=0;

712:   VecCopy(tao->solution,mfqP->Xhist[0]);
713:   CHKMEMQ;
714:   pounders_feval(tao,tao->solution,mfqP->Fhist[0],&mfqP->Fres[0]);
715:   mfqP->minindex = 0;
716:   minnorm = mfqP->Fres[0];
717:   TaoMonitor(tao, tao->niter, minnorm, PETSC_INFINITY, 0.0, step, &reason);
718:   tao->niter++;

720:   VecGetOwnershipRange(mfqP->Xhist[0],&low,&high);
721:   for (i=1;i<mfqP->n+1;i++) {
722:     VecCopy(tao->solution,mfqP->Xhist[i]);
723:     if (i-1 >= low && i-1 < high) {
724:       VecGetArray(mfqP->Xhist[i],&x);
725:       x[i-1-low] += mfqP->delta;
726:       VecRestoreArray(mfqP->Xhist[i],&x);
727:     }
728:     CHKMEMQ;
729:     pounders_feval(tao,mfqP->Xhist[i],mfqP->Fhist[i],&mfqP->Fres[i]);
730:     if (mfqP->Fres[i] < minnorm) {
731:       mfqP->minindex = i;
732:       minnorm = mfqP->Fres[i];
733:     }
734:   }
735:   VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);
736:   VecCopy(mfqP->Fhist[mfqP->minindex],tao->sep_objective);
737:   /* Gather mpi vecs to one big local vec */

739:   /* Begin serial code */

741:   /* Disp[i] = Xi-xmin, i=1,..,mfqP->minindex-1,mfqP->minindex+1,..,n */
742:   /* Fdiff[i] = (Fi-Fmin)', i=1,..,mfqP->minindex-1,mfqP->minindex+1,..,n */
743:   /* (Column oriented for blas calls) */
744:   ii=0;

746:   if (mfqP->size == 1) {
747:     VecGetArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);
748:     for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i];
749:     VecRestoreArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);
750:     VecGetArrayRead(mfqP->Fhist[mfqP->minindex],&fmin);
751:     for (i=0;i<mfqP->n+1;i++) {
752:       if (i == mfqP->minindex) continue;

754:       VecGetArray(mfqP->Xhist[i],&x);
755:       for (j=0;j<mfqP->n;j++) {
756:         mfqP->Disp[ii+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta;
757:       }
758:       VecRestoreArray(mfqP->Xhist[i],&x);

760:       VecGetArray(mfqP->Fhist[i],&f);
761:       for (j=0;j<mfqP->m;j++) {
762:         mfqP->Fdiff[ii+mfqP->n*j] = f[j] - fmin[j];
763:       }
764:       VecRestoreArray(mfqP->Fhist[i],&f);
765:       mfqP->model_indices[ii++] = i;

767:     }
768:     for (j=0;j<mfqP->m;j++) {
769:       mfqP->C[j] = fmin[j];
770:     }
771:     VecRestoreArrayRead(mfqP->Fhist[mfqP->minindex],&fmin);
772:   } else {
773:     VecScatterBegin(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);
774:     VecScatterEnd(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);
775:     VecGetArrayRead(mfqP->localxmin,&xmint);
776:     for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i];
777:     VecRestoreArrayRead(mfqP->localxmin,&xmint);

779:     VecScatterBegin(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);
780:     VecScatterEnd(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);
781:     VecGetArrayRead(mfqP->localfmin,&fmin);
782:     for (i=0;i<mfqP->n+1;i++) {
783:       if (i == mfqP->minindex) continue;

785:       mfqP->model_indices[ii++] = i;
786:       VecScatterBegin(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);
787:       VecScatterEnd(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);
788:       VecGetArray(mfqP->localx,&x);
789:       for (j=0;j<mfqP->n;j++) {
790:         mfqP->Disp[i+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta;
791:       }
792:       VecRestoreArray(mfqP->localx,&x);

794:       VecScatterBegin(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);
795:       VecScatterEnd(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);
796:       VecGetArray(mfqP->localf,&f);
797:       for (j=0;j<mfqP->m;j++) {
798:         mfqP->Fdiff[i*mfqP->n+j] = f[j] - fmin[j];
799:       }
800:       VecRestoreArray(mfqP->localf,&f);
801:     }
802:     for (j=0;j<mfqP->m;j++) {
803:       mfqP->C[j] = fmin[j];
804:     }
805:     VecRestoreArrayRead(mfqP->localfmin,&fmin);
806:   }

808:   /* Determine the initial quadratic models */
809:   /* G = D(ModelIn,:) \ (F(ModelIn,1:m)-repmat(F(xkin,1:m),n,1)); */
810:   /* D (nxn) Fdiff (nxm)  => G (nxm) */
811:   blasncopy = blasn;
812:   PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&blasn,&blasm,mfqP->Disp,&blasnpmax,mfqP->iwork,mfqP->Fdiff,&blasncopy,&info));
813:   PetscInfo1(tao,"gesv returned %d\n",info);

815:   cres = minnorm;
816:   pounders_update_res(tao);

818:   valid = PETSC_TRUE;

820:   VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);
821:   VecAssemblyBegin(tao->gradient);
822:   VecAssemblyEnd(tao->gradient);
823:   VecNorm(tao->gradient,NORM_2,&gnorm);
824:   gnorm *= mfqP->delta;
825:   VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);
826:   TaoMonitor(tao, tao->niter, minnorm, gnorm, 0.0, step, &reason);
827:   mfqP->nHist = mfqP->n+1;
828:   mfqP->nmodelpoints = mfqP->n+1;

830:   while (reason == TAO_CONTINUE_ITERATING) {
831:     PetscReal gnm = 1e-4;
832:     tao->niter++;
833:     /* Solve the subproblem min{Q(s): ||s|| <= delta} */
834:     gqtwrap(tao,&gnm,&mdec);
835:     /* Evaluate the function at the new point */

837:     for (i=0;i<mfqP->n;i++) {
838:         mfqP->work[i] = mfqP->Xsubproblem[i]*mfqP->delta + mfqP->xmin[i];
839:     }
840:     VecDuplicate(tao->solution,&mfqP->Xhist[mfqP->nHist]);
841:     VecDuplicate(tao->sep_objective,&mfqP->Fhist[mfqP->nHist]);
842:     VecSetValues(mfqP->Xhist[mfqP->nHist],mfqP->n,mfqP->indices,mfqP->work,INSERT_VALUES);
843:     VecAssemblyBegin(mfqP->Xhist[mfqP->nHist]);
844:     VecAssemblyEnd(mfqP->Xhist[mfqP->nHist]);

846:     pounders_feval(tao,mfqP->Xhist[mfqP->nHist],mfqP->Fhist[mfqP->nHist],&mfqP->Fres[mfqP->nHist]);

848:     rho = (mfqP->Fres[mfqP->minindex] - mfqP->Fres[mfqP->nHist]) / mdec;
849:     mfqP->nHist++;

851:     /* Update the center */
852:     if ((rho >= mfqP->eta1) || (rho > mfqP->eta0 && valid==PETSC_TRUE)) {
853:       /* Update model to reflect new base point */
854:       for (i=0;i<mfqP->n;i++) {
855:         mfqP->work[i] = (mfqP->work[i] - mfqP->xmin[i])/mfqP->delta;
856:       }
857:       for (j=0;j<mfqP->m;j++) {
858:         /* C(j) = C(j) + work*G(:,j) + .5*work*H(:,:,j)*work';
859:          G(:,j) = G(:,j) + H(:,:,j)*work' */
860:         for (k=0;k<mfqP->n;k++) {
861:           mfqP->work2[k]=0.0;
862:           for (l=0;l<mfqP->n;l++) {
863:             mfqP->work2[k]+=mfqP->H[j + mfqP->m*(k + l*mfqP->n)]*mfqP->work[l];
864:           }
865:         }
866:         for (i=0;i<mfqP->n;i++) {
867:           mfqP->C[j]+=mfqP->work[i]*(mfqP->Fdiff[i + mfqP->n* j] + 0.5*mfqP->work2[i]);
868:           mfqP->Fdiff[i+mfqP->n*j] +=mfqP-> work2[i];
869:         }
870:       }
871:       /* Cres += work*Gres + .5*work*Hres*work';
872:        Gres += Hres*work'; */

874:       PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasn,&one,mfqP->Hres,&blasn,mfqP->work,&ione,&zero,mfqP->work2,&ione));
875:       for (i=0;i<mfqP->n;i++) {
876:         cres += mfqP->work[i]*(mfqP->Gres[i]  + 0.5*mfqP->work2[i]);
877:         mfqP->Gres[i] += mfqP->work2[i];
878:       }
879:       mfqP->minindex = mfqP->nHist-1;
880:       minnorm = mfqP->Fres[mfqP->minindex];
881:       VecCopy(mfqP->Fhist[mfqP->minindex],tao->sep_objective);
882:       /* Change current center */
883:       VecGetArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);
884:       for (i=0;i<mfqP->n;i++) {
885:         mfqP->xmin[i] = xmint[i];
886:       }
887:       VecRestoreArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);
888:     }

890:     /* Evaluate at a model-improving point if necessary */
891:     if (valid == PETSC_FALSE) {
892:       mfqP->q_is_I = 1;
893:       mfqP->nmodelpoints = 0;
894:       affpoints(mfqP,mfqP->xmin,mfqP->c1);
895:       if (mfqP->nmodelpoints < mfqP->n) {
896:         PetscInfo(tao,"Model not valid -- model-improving\n");
897:         modelimprove(tao,mfqP,1);
898:       }
899:     }

901:     /* Update the trust region radius */
902:     deltaold = mfqP->delta;
903:     normxsp = 0;
904:     for (i=0;i<mfqP->n;i++) {
905:       normxsp += mfqP->Xsubproblem[i]*mfqP->Xsubproblem[i];
906:     }
907:     normxsp = PetscSqrtReal(normxsp);
908:     if (rho >= mfqP->eta1 && normxsp > 0.5*mfqP->delta) {
909:       mfqP->delta = PetscMin(mfqP->delta*mfqP->gamma1,mfqP->deltamax);
910:     } else if (valid == PETSC_TRUE) {
911:       mfqP->delta = PetscMax(mfqP->delta*mfqP->gamma0,mfqP->deltamin);
912:     }

914:     /* Compute the next interpolation set */
915:     mfqP->q_is_I = 1;
916:     mfqP->nmodelpoints=0;
917:     affpoints(mfqP,mfqP->xmin,mfqP->c1);
918:     if (mfqP->nmodelpoints == mfqP->n) {
919:       valid = PETSC_TRUE;
920:     } else {
921:       valid = PETSC_FALSE;
922:       affpoints(mfqP,mfqP->xmin,mfqP->c2);
923:       if (mfqP->n > mfqP->nmodelpoints) {
924:         PetscInfo(tao,"Model not valid -- adding geometry points\n");
925:         modelimprove(tao,mfqP,mfqP->n - mfqP->nmodelpoints);
926:       }
927:     }
928:     for (i=mfqP->nmodelpoints;i>0;i--) {
929:       mfqP->model_indices[i] = mfqP->model_indices[i-1];
930:     }
931:     mfqP->nmodelpoints++;
932:     mfqP->model_indices[0] = mfqP->minindex;
933:     morepoints(mfqP);
934:     for (i=0;i<mfqP->nmodelpoints;i++) {
935:       VecGetArray(mfqP->Xhist[mfqP->model_indices[i]],&x);
936:       for (j=0;j<mfqP->n;j++) {
937:         mfqP->Disp[i + mfqP->npmax*j] = (x[j]  - mfqP->xmin[j]) / deltaold;
938:       }
939:       VecRestoreArray(mfqP->Xhist[mfqP->model_indices[i]],&x);
940:       VecGetArray(mfqP->Fhist[mfqP->model_indices[i]],&f);
941:       for (j=0;j<mfqP->m;j++) {
942:         for (k=0;k<mfqP->n;k++)  {
943:           mfqP->work[k]=0.0;
944:           for (l=0;l<mfqP->n;l++) {
945:             mfqP->work[k] += mfqP->H[j + mfqP->m*(k + mfqP->n*l)] * mfqP->Disp[i + mfqP->npmax*l];
946:           }
947:         }
948:         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];
949:       }
950:       VecRestoreArray(mfqP->Fhist[mfqP->model_indices[i]],&f);
951:     }

953:     /* Update the quadratic model */
954:     getquadpounders(mfqP);
955:     VecGetArrayRead(mfqP->Fhist[mfqP->minindex],&fmin);
956:     PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasm,fmin,&ione,mfqP->C,&ione));
957:     /* G = G*(delta/deltaold) + Gdel */
958:     ratio = mfqP->delta/deltaold;
959:     iblas = blasm*blasn;
960:     PetscStackCallBLAS("BLASscal",BLASscal_(&iblas,&ratio,mfqP->Fdiff,&ione));
961:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&iblas,&one,mfqP->Gdel,&ione,mfqP->Fdiff,&ione));
962:     /* H = H*(delta/deltaold)^2 + Hdel */
963:     iblas = blasm*blasn*blasn;
964:     ratio *= ratio;
965:     PetscStackCallBLAS("BLASscal",BLASscal_(&iblas,&ratio,mfqP->H,&ione));
966:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&iblas,&one,mfqP->Hdel,&ione,mfqP->H,&ione));

968:     /* Get residuals */
969:     cres = mfqP->Fres[mfqP->minindex];
970:     pounders_update_res(tao);

972:     /* Export solution and gradient residual to TAO */
973:     VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);
974:     VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);
975:     VecAssemblyBegin(tao->gradient);
976:     VecAssemblyEnd(tao->gradient);
977:     VecNorm(tao->gradient,NORM_2,&gnorm);
978:     gnorm *= mfqP->delta;
979:     /*  final criticality test */
980:     TaoMonitor(tao, tao->niter, minnorm, gnorm, 0.0, step, &reason);
981:     /* test for repeated model */
982:     if (mfqP->nmodelpoints==mfqP->last_nmodelpoints) {
983:       same = PETSC_TRUE;
984:     } else {
985:       same = PETSC_FALSE;
986:     }
987:     for (i=0;i<mfqP->nmodelpoints;i++) {
988:       if (same) {
989:         if (mfqP->model_indices[i] == mfqP->last_model_indices[i]) {
990:           same = PETSC_TRUE;
991:         } else {
992:           same = PETSC_FALSE;
993:         }
994:       }
995:       mfqP->last_model_indices[i] = mfqP->model_indices[i];
996:     }
997:     mfqP->last_nmodelpoints = mfqP->nmodelpoints;
998:     if (same && mfqP->delta == deltaold) {
999:       PetscInfo(tao,"Identical model used in successive iterations\n");
1000:       reason = TAO_CONVERGED_STEPTOL;
1001:       tao->reason = TAO_CONVERGED_STEPTOL;
1002:     }
1003:   }
1004:   return(0);
1005: }

1009: static PetscErrorCode TaoSetUp_POUNDERS(Tao tao)
1010: {
1011:   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;
1012:   PetscInt       i,j;
1013:   IS             isfloc,isfglob,isxloc,isxglob;

1017:   if (!tao->gradient) {VecDuplicate(tao->solution,&tao->gradient);  }
1018:   if (!tao->stepdirection) {VecDuplicate(tao->solution,&tao->stepdirection);  }
1019:   VecGetSize(tao->solution,&mfqP->n);
1020:   VecGetSize(tao->sep_objective,&mfqP->m);
1021:   mfqP->c1 = PetscSqrtReal((PetscReal)mfqP->n);
1022:   if (mfqP->npmax == PETSC_DEFAULT) {
1023:     mfqP->npmax = 2*mfqP->n + 1;
1024:   }
1025:   mfqP->npmax = PetscMin((mfqP->n+1)*(mfqP->n+2)/2,mfqP->npmax);
1026:   mfqP->npmax = PetscMax(mfqP->npmax, mfqP->n+2);

1028:   PetscMalloc1(tao->max_funcs+10,&mfqP->Xhist);
1029:   PetscMalloc1(tao->max_funcs+10,&mfqP->Fhist);
1030:   for (i=0;i<mfqP->n +1;i++) {
1031:     VecDuplicate(tao->solution,&mfqP->Xhist[i]);
1032:     VecDuplicate(tao->sep_objective,&mfqP->Fhist[i]);
1033:   }
1034:   VecDuplicate(tao->solution,&mfqP->workxvec);
1035:   VecDuplicate(tao->sep_objective,&mfqP->workfvec);
1036:   mfqP->nHist = 0;

1038:   PetscMalloc1(tao->max_funcs+10,&mfqP->Fres);
1039:   PetscMalloc1(mfqP->npmax*mfqP->m,&mfqP->RES);
1040:   PetscMalloc1(mfqP->n,&mfqP->work);
1041:   PetscMalloc1(mfqP->n,&mfqP->work2);
1042:   PetscMalloc1(mfqP->n,&mfqP->work3);
1043:   PetscMalloc1(PetscMax(mfqP->m,mfqP->n+1),&mfqP->mwork);
1044:   PetscMalloc1(mfqP->npmax - mfqP->n - 1,&mfqP->omega);
1045:   PetscMalloc1(mfqP->n * (mfqP->n+1) / 2,&mfqP->beta);
1046:   PetscMalloc1(mfqP->n + 1 ,&mfqP->alpha);

1048:   PetscMalloc1(mfqP->n*mfqP->n*mfqP->m,&mfqP->H);
1049:   PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q);
1050:   PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q_tmp);
1051:   PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L);
1052:   PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_tmp);
1053:   PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_save);
1054:   PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->N);
1055:   PetscMalloc1(mfqP->npmax * (mfqP->n+1) ,&mfqP->M);
1056:   PetscMalloc1(mfqP->npmax * (mfqP->npmax - mfqP->n - 1) , &mfqP->Z);
1057:   PetscMalloc1(mfqP->npmax,&mfqP->tau);
1058:   PetscMalloc1(mfqP->npmax,&mfqP->tau_tmp);
1059:   mfqP->nmax = PetscMax(5*mfqP->npmax,mfqP->n*(mfqP->n+1)/2);
1060:   PetscMalloc1(mfqP->nmax,&mfqP->npmaxwork);
1061:   PetscMalloc1(mfqP->nmax,&mfqP->npmaxiwork);
1062:   PetscMalloc1(mfqP->n,&mfqP->xmin);
1063:   PetscMalloc1(mfqP->m,&mfqP->C);
1064:   PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Fdiff);
1065:   PetscMalloc1(mfqP->npmax*mfqP->n,&mfqP->Disp);
1066:   PetscMalloc1(mfqP->n,&mfqP->Gres);
1067:   PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Hres);
1068:   PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Gpoints);
1069:   PetscMalloc1(mfqP->npmax,&mfqP->model_indices);
1070:   PetscMalloc1(mfqP->npmax,&mfqP->last_model_indices);
1071:   PetscMalloc1(mfqP->n,&mfqP->Xsubproblem);
1072:   PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Gdel);
1073:   PetscMalloc1(mfqP->n*mfqP->n*mfqP->m, &mfqP->Hdel);
1074:   PetscMalloc1(PetscMax(mfqP->m,mfqP->n),&mfqP->indices);
1075:   PetscMalloc1(mfqP->n,&mfqP->iwork);
1076:   PetscMalloc1(mfqP->m*mfqP->m,&mfqP->w);
1077:   for (i=0;i<mfqP->m;i++) {
1078:     for (j=0;j<mfqP->m;j++) {
1079:       if (i==j) {
1080:         mfqP->w[i+mfqP->m*j]=1.0;
1081:       } else {
1082:         mfqP->w[i+mfqP->m*j]=0.0;
1083:       }
1084:     }
1085:   }
1086:   for (i=0;i<PetscMax(mfqP->m,mfqP->n);i++) {
1087:     mfqP->indices[i] = i;
1088:   }
1089:   MPI_Comm_size(((PetscObject)tao)->comm,&mfqP->size);
1090:   if (mfqP->size > 1) {
1091:     VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localx);
1092:     VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localxmin);
1093:     VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localf);
1094:     VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localfmin);
1095:     ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxloc);
1096:     ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxglob);
1097:     ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfloc);
1098:     ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfglob);


1101:     VecScatterCreate(tao->solution,isxglob,mfqP->localx,isxloc,&mfqP->scatterx);
1102:     VecScatterCreate(tao->sep_objective,isfglob,mfqP->localf,isfloc,&mfqP->scatterf);

1104:     ISDestroy(&isxloc);
1105:     ISDestroy(&isxglob);
1106:     ISDestroy(&isfloc);
1107:     ISDestroy(&isfglob);
1108:   }

1110:   if (!mfqP->usegqt) {
1111:     KSP       ksp;
1112:     PC        pc;
1113:     VecCreateSeqWithArray(PETSC_COMM_SELF,mfqP->n,mfqP->n,mfqP->Xsubproblem,&mfqP->subx);
1114:     VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->subxl);
1115:     VecDuplicate(mfqP->subxl,&mfqP->subb);
1116:     VecDuplicate(mfqP->subxl,&mfqP->subxu);
1117:     VecDuplicate(mfqP->subxl,&mfqP->subpdel);
1118:     VecDuplicate(mfqP->subxl,&mfqP->subndel);
1119:     TaoCreate(PETSC_COMM_SELF,&mfqP->subtao);
1120:     TaoSetType(mfqP->subtao,TAOTRON);
1121:     TaoSetOptionsPrefix(mfqP->subtao,"pounders_subsolver_");
1122:     TaoSetInitialVector(mfqP->subtao,mfqP->subx);
1123:     TaoSetObjectiveAndGradientRoutine(mfqP->subtao,pounders_fg,(void*)mfqP);
1124:     TaoSetMaximumIterations(mfqP->subtao,mfqP->gqt_maxits);
1125:     TaoSetFromOptions(mfqP->subtao);
1126:     TaoGetKSP(mfqP->subtao,&ksp);
1127:     if (ksp) {
1128:       KSPGetPC(ksp,&pc);
1129:       PCSetType(pc,PCNONE);
1130:     }
1131:     TaoSetVariableBounds(mfqP->subtao,mfqP->subxl,mfqP->subxu);
1132:     MatCreateSeqDense(PETSC_COMM_SELF,mfqP->n,mfqP->n,mfqP->Hres,&mfqP->subH);
1133:     TaoSetHessianRoutine(mfqP->subtao,mfqP->subH,mfqP->subH,pounders_h,(void*)mfqP);
1134:   }
1135:   return(0);
1136: }

1140: static PetscErrorCode TaoDestroy_POUNDERS(Tao tao)
1141: {
1142:   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;
1143:   PetscInt       i;

1147:   if (!mfqP->usegqt) {
1148:     TaoDestroy(&mfqP->subtao);
1149:     VecDestroy(&mfqP->subx);
1150:     VecDestroy(&mfqP->subxl);
1151:     VecDestroy(&mfqP->subxu);
1152:     VecDestroy(&mfqP->subb);
1153:     VecDestroy(&mfqP->subpdel);
1154:     VecDestroy(&mfqP->subndel);
1155:     MatDestroy(&mfqP->subH);
1156:   }
1157:   PetscFree(mfqP->Fres);
1158:   PetscFree(mfqP->RES);
1159:   PetscFree(mfqP->work);
1160:   PetscFree(mfqP->work2);
1161:   PetscFree(mfqP->work3);
1162:   PetscFree(mfqP->mwork);
1163:   PetscFree(mfqP->omega);
1164:   PetscFree(mfqP->beta);
1165:   PetscFree(mfqP->alpha);
1166:   PetscFree(mfqP->H);
1167:   PetscFree(mfqP->Q);
1168:   PetscFree(mfqP->Q_tmp);
1169:   PetscFree(mfqP->L);
1170:   PetscFree(mfqP->L_tmp);
1171:   PetscFree(mfqP->L_save);
1172:   PetscFree(mfqP->N);
1173:   PetscFree(mfqP->M);
1174:   PetscFree(mfqP->Z);
1175:   PetscFree(mfqP->tau);
1176:   PetscFree(mfqP->tau_tmp);
1177:   PetscFree(mfqP->npmaxwork);
1178:   PetscFree(mfqP->npmaxiwork);
1179:   PetscFree(mfqP->xmin);
1180:   PetscFree(mfqP->C);
1181:   PetscFree(mfqP->Fdiff);
1182:   PetscFree(mfqP->Disp);
1183:   PetscFree(mfqP->Gres);
1184:   PetscFree(mfqP->Hres);
1185:   PetscFree(mfqP->Gpoints);
1186:   PetscFree(mfqP->model_indices);
1187:   PetscFree(mfqP->last_model_indices);
1188:   PetscFree(mfqP->Xsubproblem);
1189:   PetscFree(mfqP->Gdel);
1190:   PetscFree(mfqP->Hdel);
1191:   PetscFree(mfqP->indices);
1192:   PetscFree(mfqP->iwork);
1193:   PetscFree(mfqP->w);
1194:   for (i=0;i<mfqP->nHist;i++) {
1195:     VecDestroy(&mfqP->Xhist[i]);
1196:     VecDestroy(&mfqP->Fhist[i]);
1197:   }
1198:   VecDestroy(&mfqP->workxvec);
1199:   VecDestroy(&mfqP->workfvec);
1200:   PetscFree(mfqP->Xhist);
1201:   PetscFree(mfqP->Fhist);

1203:   if (mfqP->size > 1) {
1204:     VecDestroy(&mfqP->localx);
1205:     VecDestroy(&mfqP->localxmin);
1206:     VecDestroy(&mfqP->localf);
1207:     VecDestroy(&mfqP->localfmin);
1208:   }
1209:   PetscFree(tao->data);
1210:   return(0);
1211: }

1215: static PetscErrorCode TaoSetFromOptions_POUNDERS(PetscOptionItems *PetscOptionsObject,Tao tao)
1216: {
1217:   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;

1221:   PetscOptionsHead(PetscOptionsObject,"POUNDERS method for least-squares optimization");
1222:   PetscOptionsReal("-tao_pounders_delta","initial delta","",mfqP->delta,&mfqP->delta0,NULL);
1223:   mfqP->delta = mfqP->delta0;
1224:   mfqP->npmax = PETSC_DEFAULT;
1225:   PetscOptionsInt("-tao_pounders_npmax","max number of points in model","",mfqP->npmax,&mfqP->npmax,NULL);
1226:   mfqP->usegqt = PETSC_FALSE;
1227:   PetscOptionsBool("-tao_pounders_gqt","use gqt algorithm for subproblem","",mfqP->usegqt,&mfqP->usegqt,NULL);
1228:   PetscOptionsTail();
1229:   return(0);
1230: }

1234: static PetscErrorCode TaoView_POUNDERS(Tao tao, PetscViewer viewer)
1235: {
1236:   TAO_POUNDERS   *mfqP = (TAO_POUNDERS *)tao->data;
1237:   PetscBool      isascii;
1238:   PetscInt       nits;

1242:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1243:   if (isascii) {
1244:     PetscViewerASCIIPushTab(viewer);
1245:     PetscViewerASCIIPrintf(viewer, "initial delta: %g\n",(double)mfqP->delta0);
1246:     PetscViewerASCIIPrintf(viewer, "final delta: %g\n",(double)mfqP->delta);
1247:     PetscViewerASCIIPrintf(viewer, "model points: %D\n",mfqP->nmodelpoints);
1248:     if (mfqP->usegqt) {
1249:       PetscViewerASCIIPrintf(viewer, "subproblem solver: gqt\n");
1250:     } else {
1251:       PetscViewerASCIIPrintf(viewer, "subproblem solver: %s\n",((PetscObject)mfqP->subtao)->type_name);
1252:       TaoGetTotalIterationNumber(mfqP->subtao,&nits);
1253:       PetscViewerASCIIPrintf(viewer, "total subproblem iterations: %D\n",nits);
1254:     }
1255:     PetscViewerASCIIPopTab(viewer);
1256:   }
1257:   return(0);
1258: }
1259: /*MC
1260:   TAOPOUNDERS - POUNDERS derivate-free model-based algorithm for nonlinear least squares

1262:   Options Database Keys:
1263: + -tao_pounders_delta - initial step length
1264: . -tao_pounders_npmax - maximum number of points in model
1265: - -tao_pounders_gqt - use gqt algorithm for subproblem instead of TRON

1267:   Level: beginner
1268:  
1269: M*/

1273: PETSC_EXTERN PetscErrorCode TaoCreate_POUNDERS(Tao tao)
1274: {
1275:   TAO_POUNDERS   *mfqP = (TAO_POUNDERS*)tao->data;

1279:   tao->ops->setup = TaoSetUp_POUNDERS;
1280:   tao->ops->solve = TaoSolve_POUNDERS;
1281:   tao->ops->view = TaoView_POUNDERS;
1282:   tao->ops->setfromoptions = TaoSetFromOptions_POUNDERS;
1283:   tao->ops->destroy = TaoDestroy_POUNDERS;

1285:   PetscNewLog(tao,&mfqP);
1286:   tao->data = (void*)mfqP;
1287:   /* Override default settings (unless already changed) */
1288:   if (!tao->max_it_changed) tao->max_it = 2000;
1289:   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
1290:   mfqP->delta0 = 0.1;
1291:   mfqP->delta = 0.1;
1292:   mfqP->deltamax=1e3;
1293:   mfqP->deltamin=1e-6;
1294:   mfqP->c2 = 100.0;
1295:   mfqP->theta1=1e-5;
1296:   mfqP->theta2=1e-4;
1297:   mfqP->gamma0=0.5;
1298:   mfqP->gamma1=2.0;
1299:   mfqP->eta0 = 0.0;
1300:   mfqP->eta1 = 0.1;
1301:   mfqP->gqt_rtol = 0.001;
1302:   mfqP->gqt_maxits = 50;
1303:   mfqP->workxvec = 0;
1304:   return(0);
1305: }