1 #include <../src/tao/leastsquares/impls/pounders/pounders.h> 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "pounders_h" 5 static PetscErrorCode pounders_h(Tao subtao, Vec v, Mat H, Mat Hpre, void *ctx) 6 { 7 PetscFunctionBegin; 8 PetscFunctionReturn(0); 9 } 10 #undef __FUNCT__ 11 #define __FUNCT__ "pounders_fg" 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; 16 PetscErrorCode ierr; 17 18 PetscFunctionBegin; 19 /* g = A*x (add b later)*/ 20 ierr = MatMult(mfqP->subH,x,g);CHKERRQ(ierr); 21 22 /* f = 1/2 * x'*(Ax) + b'*x */ 23 ierr = VecDot(x,g,&d1);CHKERRQ(ierr); 24 ierr = VecDot(mfqP->subb,x,&d2);CHKERRQ(ierr); 25 *f = 0.5 *d1 + d2; 26 27 /* now g = g + b */ 28 ierr = VecAXPY(g, 1.0, mfqP->subb);CHKERRQ(ierr); 29 PetscFunctionReturn(0); 30 } 31 32 #undef __FUNCT__ 33 #define __FUNCT__ "gqtwrap" 34 PetscErrorCode gqtwrap(Tao tao,PetscReal *gnorm, PetscReal *qmin) 35 { 36 PetscErrorCode ierr; 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; 44 45 PetscFunctionBegin; 46 if (! mfqP->usegqt) { 47 PetscReal maxval; 48 PetscInt i,j; 49 50 ierr = VecSetValues(mfqP->subb,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);CHKERRQ(ierr); 51 ierr = VecAssemblyBegin(mfqP->subb);CHKERRQ(ierr); 52 ierr = VecAssemblyEnd(mfqP->subb);CHKERRQ(ierr); 53 54 ierr = VecSet(mfqP->subx,0.0);CHKERRQ(ierr); 55 56 ierr = VecSet(mfqP->subndel,-mfqP->delta);CHKERRQ(ierr); 57 ierr = VecSet(mfqP->subpdel,mfqP->delta);CHKERRQ(ierr); 58 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 ierr = MatSetValues(mfqP->subH,mfqP->n,mfqP->indices,mfqP->n,mfqP->indices,mfqP->Hres,INSERT_VALUES);CHKERRQ(ierr); 65 ierr = MatAssemblyBegin(mfqP->subH,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 66 ierr = MatAssemblyEnd(mfqP->subH,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 67 68 ierr = TaoResetStatistics(mfqP->subtao);CHKERRQ(ierr); 69 ierr = TaoSetTolerances(mfqP->subtao,*gnorm,*gnorm,PETSC_DEFAULT);CHKERRQ(ierr); 70 /* enforce bound constraints -- experimental */ 71 if (tao->XU && tao->XL) { 72 ierr = VecCopy(tao->XU,mfqP->subxu);CHKERRQ(ierr); 73 ierr = VecAXPY(mfqP->subxu,-1.0,tao->solution);CHKERRQ(ierr); 74 ierr = VecScale(mfqP->subxu,1.0/mfqP->delta);CHKERRQ(ierr); 75 ierr = VecCopy(tao->XL,mfqP->subxl);CHKERRQ(ierr); 76 ierr = VecAXPY(mfqP->subxl,-1.0,tao->solution);CHKERRQ(ierr); 77 ierr = VecScale(mfqP->subxl,1.0/mfqP->delta);CHKERRQ(ierr); 78 79 ierr = VecPointwiseMin(mfqP->subxu,mfqP->subxu,mfqP->subpdel);CHKERRQ(ierr); 80 ierr = VecPointwiseMax(mfqP->subxl,mfqP->subxl,mfqP->subndel);CHKERRQ(ierr); 81 } else { 82 ierr = VecCopy(mfqP->subpdel,mfqP->subxu);CHKERRQ(ierr); 83 ierr = VecCopy(mfqP->subndel,mfqP->subxl);CHKERRQ(ierr); 84 } 85 /* Make sure xu > xl */ 86 ierr = VecCopy(mfqP->subxl,mfqP->subpdel);CHKERRQ(ierr); 87 ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);CHKERRQ(ierr); 88 ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr); 89 if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"upper bound < lower bound in subproblem"); 90 /* Make sure xu > tao->solution > xl */ 91 ierr = VecCopy(mfqP->subxl,mfqP->subpdel);CHKERRQ(ierr); 92 ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subx);CHKERRQ(ierr); 93 ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr); 94 if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"initial guess < lower bound in subproblem"); 95 96 ierr = VecCopy(mfqP->subx,mfqP->subpdel);CHKERRQ(ierr); 97 ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);CHKERRQ(ierr); 98 ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr); 99 if (maxval > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"initial guess > upper bound in subproblem"); 100 101 ierr = TaoSolve(mfqP->subtao);CHKERRQ(ierr); 102 ierr = TaoGetSolutionStatus(mfqP->subtao,NULL,qmin,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 103 104 /* test bounds post-solution*/ 105 ierr = VecCopy(mfqP->subxl,mfqP->subpdel);CHKERRQ(ierr); 106 ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subx);CHKERRQ(ierr); 107 ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr); 108 if (maxval > 1e-5) { 109 ierr = PetscInfo(tao,"subproblem solution < lower bound\n");CHKERRQ(ierr); 110 tao->reason = TAO_DIVERGED_TR_REDUCTION; 111 } 112 113 ierr = VecCopy(mfqP->subx,mfqP->subpdel);CHKERRQ(ierr); 114 ierr = VecAXPY(mfqP->subpdel,-1.0,mfqP->subxu);CHKERRQ(ierr); 115 ierr = VecMax(mfqP->subpdel,NULL,&maxval);CHKERRQ(ierr); 116 if (maxval > 1e-5) { 117 ierr = PetscInfo(tao,"subproblem solution > upper bound\n");CHKERRQ(ierr); 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 PetscFunctionReturn(0); 125 } 126 127 #undef __FUNCT__ 128 #define __FUNCT__ "phi2eval" 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); 134 135 PetscFunctionBegin; 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 PetscFunctionReturn(0); 146 } 147 148 #undef __FUNCT__ 149 #define __FUNCT__ "getquadpounders" 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 */ 155 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); 166 167 PetscFunctionBegin; 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 } 174 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); 178 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)); 189 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 } 194 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); 201 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 } 205 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); 210 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 } 215 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 PetscFunctionReturn(0); 231 } 232 233 #undef __FUNCT__ 234 #define __FUNCT__ "morepoints" 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; 247 248 PetscFunctionBegin; 249 /* Initialize M,N */ 250 for (i=0;i<mfqP->n+1;i++) { 251 ierr = VecGetArrayRead(mfqP->Xhist[mfqP->model_indices[i]],&x);CHKERRQ(ierr); 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 ierr = VecRestoreArrayRead(mfqP->Xhist[mfqP->model_indices[i]],&x);CHKERRQ(ierr); 257 ierr = phi2eval(&mfqP->M[1+((mfqP->n+1)*i)],mfqP->n,&mfqP->N[mfqP->n*(mfqP->n+1)/2 * i]);CHKERRQ(ierr); 258 } 259 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 } 272 273 /* Reject if norm(d) >c2 */ 274 if (!reject) { 275 ierr = VecCopy(mfqP->Xhist[point],mfqP->workxvec);CHKERRQ(ierr); 276 ierr = VecAXPY(mfqP->workxvec,-1.0,mfqP->Xhist[mfqP->minindex]);CHKERRQ(ierr); 277 ierr = VecNorm(mfqP->workxvec,NORM_2,&normd);CHKERRQ(ierr); 278 normd /= mfqP->delta; 279 if (normd > mfqP->c2) { 280 reject =1; 281 } 282 } 283 if (reject){ 284 point--; 285 continue; 286 } 287 288 ierr = VecGetArrayRead(mfqP->Xhist[point],&x);CHKERRQ(ierr); 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 ierr = VecRestoreArrayRead(mfqP->Xhist[point],&x);CHKERRQ(ierr); 294 ierr = phi2eval(&mfqP->M[1+(mfqP->n+1)*mfqP->nmodelpoints],mfqP->n,&mfqP->N[mfqP->n*(mfqP->n+1)/2 * (mfqP->nmodelpoints)]);CHKERRQ(ierr); 295 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); 308 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 */ 317 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); 321 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 } 326 327 /* Get svd for L_tmp(:,n+2:np+1) (L_tmp is modified in process) */ 328 blasint = mfqP->nmodelpoints - mfqP->n; 329 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 330 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)); 331 ierr = PetscFPTrapPop();CHKERRQ(ierr); 332 if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine gesvd returned with value %d\n",info); 333 334 if (mfqP->beta[PetscMin(blasint,blasint2)-1] > mfqP->theta2) { 335 /* accept point */ 336 mfqP->model_indices[mfqP->nmodelpoints] = point; 337 /* Copy Q_tmp to Q */ 338 for (i=0;i<mfqP->npmax* mfqP->npmax;i++) { 339 mfqP->Q[i] = mfqP->Q_tmp[i]; 340 } 341 for (i=0;i<mfqP->npmax;i++){ 342 mfqP->tau[i] = mfqP->tau_tmp[i]; 343 } 344 mfqP->nmodelpoints++; 345 blasnp = mfqP->nmodelpoints; 346 347 /* Copy L_save to L */ 348 for (i=0;i<mfqP->npmax * mfqP->n*(mfqP->n+1)/2;i++) { 349 mfqP->L[i] = mfqP->L_save[i]; 350 } 351 } 352 point--; 353 } 354 355 blasnp = mfqP->nmodelpoints; 356 /* Copy Q(:,n+2:np) to Z */ 357 /* First set Q_tmp to I */ 358 for (i=0;i<mfqP->npmax*mfqP->npmax;i++) { 359 mfqP->Q_tmp[i] = 0.0; 360 } 361 for (i=0;i<mfqP->npmax;i++) { 362 mfqP->Q_tmp[i + mfqP->npmax*i] = 1.0; 363 } 364 365 /* Q_tmp = I * Q */ 366 PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("R","N",&blasnp,&blasnp,&blasnplus1,mfqP->Q,&blasnpmax,mfqP->tau,mfqP->Q_tmp,&blasnpmax,mfqP->npmaxwork,&blasnmax,&info)); 367 if (info != 0) SETERRQ1(PETSC_COMM_SELF,1,"LAPACK routine ormqr returned with value %d\n",info); 368 369 /* Copy Q_tmp(:,n+2:np) to Z) */ 370 offset = mfqP->npmax * (mfqP->n+1); 371 for (i=offset;i<mfqP->npmax*mfqP->npmax;i++) { 372 mfqP->Z[i-offset] = mfqP->Q_tmp[i]; 373 } 374 375 if (mfqP->nmodelpoints == mfqP->n + 1) { 376 /* Set L to I_{n+1} */ 377 for (i=0;i<mfqP->npmax * mfqP->n*(mfqP->n+1)/2;i++) { 378 mfqP->L[i] = 0.0; 379 } 380 for (i=0;i<mfqP->n;i++) { 381 mfqP->L[(mfqP->n*(mfqP->n+1)/2)*i + i] = 1.0; 382 } 383 } 384 PetscFunctionReturn(0); 385 } 386 387 #undef __FUNCT__ 388 #define __FUNCT__ "addpoint" 389 /* Only call from modelimprove, addpoint() needs ->Q_tmp and ->work to be set */ 390 PetscErrorCode addpoint(Tao tao, TAO_POUNDERS *mfqP, PetscInt index) 391 { 392 PetscErrorCode ierr; 393 394 PetscFunctionBegin; 395 /* Create new vector in history: X[newidx] = X[mfqP->index] + delta*X[index]*/ 396 ierr = VecDuplicate(mfqP->Xhist[0],&mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr); 397 ierr = VecSetValues(mfqP->Xhist[mfqP->nHist],mfqP->n,mfqP->indices,&mfqP->Q_tmp[index*mfqP->npmax],INSERT_VALUES);CHKERRQ(ierr); 398 ierr = VecAssemblyBegin(mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr); 399 ierr = VecAssemblyEnd(mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr); 400 ierr = VecAYPX(mfqP->Xhist[mfqP->nHist],mfqP->delta,mfqP->Xhist[mfqP->minindex]);CHKERRQ(ierr); 401 402 /* Project into feasible region */ 403 if (tao->XU && tao->XL) { 404 ierr = VecMedian(mfqP->Xhist[mfqP->nHist], tao->XL, tao->XU, mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr); 405 } 406 407 /* Compute value of new vector */ 408 ierr = VecDuplicate(mfqP->Fhist[0],&mfqP->Fhist[mfqP->nHist]);CHKERRQ(ierr); 409 CHKMEMQ; 410 ierr = TaoComputeSeparableObjective(tao,mfqP->Xhist[mfqP->nHist],mfqP->Fhist[mfqP->nHist]);CHKERRQ(ierr); 411 ierr = VecNorm(mfqP->Fhist[mfqP->nHist],NORM_2,&mfqP->Fres[mfqP->nHist]);CHKERRQ(ierr); 412 if (PetscIsInfOrNanReal(mfqP->Fres[mfqP->nHist])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 413 mfqP->Fres[mfqP->nHist]*=mfqP->Fres[mfqP->nHist]; 414 415 /* Add new vector to model */ 416 mfqP->model_indices[mfqP->nmodelpoints] = mfqP->nHist; 417 mfqP->nmodelpoints++; 418 mfqP->nHist++; 419 PetscFunctionReturn(0); 420 } 421 422 #undef __FUNCT__ 423 #define __FUNCT__ "modelimprove" 424 PetscErrorCode modelimprove(Tao tao, TAO_POUNDERS *mfqP, PetscInt addallpoints) 425 { 426 /* modeld = Q(:,np+1:n)' */ 427 PetscErrorCode ierr; 428 PetscInt i,j,minindex=0; 429 PetscReal dp,half=0.5,one=1.0,minvalue=PETSC_INFINITY; 430 PetscBLASInt blasn=mfqP->n, blasnpmax = mfqP->npmax, blask,info; 431 PetscBLASInt blas1=1,blasnmax = mfqP->nmax; 432 433 blask = mfqP->nmodelpoints; 434 /* Qtmp = I(n x n) */ 435 for (i=0;i<mfqP->n;i++) { 436 for (j=0;j<mfqP->n;j++) { 437 mfqP->Q_tmp[i + mfqP->npmax*j] = 0.0; 438 } 439 } 440 for (j=0;j<mfqP->n;j++) { 441 mfqP->Q_tmp[j + mfqP->npmax*j] = 1.0; 442 } 443 444 /* Qtmp = Q * I */ 445 PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("R","N",&blasn,&blasn,&blask,mfqP->Q,&blasnpmax,mfqP->tau, mfqP->Q_tmp, &blasnpmax, mfqP->npmaxwork,&blasnmax, &info)); 446 447 for (i=mfqP->nmodelpoints;i<mfqP->n;i++) { 448 dp = BLASdot_(&blasn,&mfqP->Q_tmp[i*mfqP->npmax],&blas1,mfqP->Gres,&blas1); 449 if (dp>0.0) { /* Model says use the other direction! */ 450 for (j=0;j<mfqP->n;j++) { 451 mfqP->Q_tmp[i*mfqP->npmax+j] *= -1; 452 } 453 } 454 /* mfqP->work[i] = Cres+Modeld(i,:)*(Gres+.5*Hres*Modeld(i,:)') */ 455 for (j=0;j<mfqP->n;j++) { 456 mfqP->work2[j] = mfqP->Gres[j]; 457 } 458 PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasn,&half,mfqP->Hres,&blasn,&mfqP->Q_tmp[i*mfqP->npmax], &blas1, &one, mfqP->work2,&blas1)); 459 mfqP->work[i] = BLASdot_(&blasn,&mfqP->Q_tmp[i*mfqP->npmax],&blas1,mfqP->work2,&blas1); 460 if (i==mfqP->nmodelpoints || mfqP->work[i] < minvalue) { 461 minindex=i; 462 minvalue = mfqP->work[i]; 463 } 464 if (addallpoints != 0) { 465 ierr = addpoint(tao,mfqP,i);CHKERRQ(ierr); 466 } 467 } 468 if (!addallpoints) { 469 ierr = addpoint(tao,mfqP,minindex);CHKERRQ(ierr); 470 } 471 PetscFunctionReturn(0); 472 } 473 474 475 #undef __FUNCT__ 476 #define __FUNCT__ "affpoints" 477 PetscErrorCode affpoints(TAO_POUNDERS *mfqP, PetscReal *xmin,PetscReal c) 478 { 479 PetscInt i,j; 480 PetscBLASInt blasm=mfqP->m,blasj,blask,blasn=mfqP->n,ione=1,info; 481 PetscBLASInt blasnpmax = mfqP->npmax,blasmaxmn; 482 PetscReal proj,normd; 483 const PetscReal *x; 484 PetscErrorCode ierr; 485 486 PetscFunctionBegin; 487 for (i=mfqP->nHist-1;i>=0;i--) { 488 ierr = VecGetArrayRead(mfqP->Xhist[i],&x);CHKERRQ(ierr); 489 for (j=0;j<mfqP->n;j++) { 490 mfqP->work[j] = (x[j] - xmin[j])/mfqP->delta; 491 } 492 ierr = VecRestoreArrayRead(mfqP->Xhist[i],&x);CHKERRQ(ierr); 493 PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasn,mfqP->work,&ione,mfqP->work2,&ione)); 494 normd = BLASnrm2_(&blasn,mfqP->work,&ione); 495 if (normd <= c*c) { 496 blasj=PetscMax((mfqP->n - mfqP->nmodelpoints),0); 497 if (!mfqP->q_is_I) { 498 /* project D onto null */ 499 blask=(mfqP->nmodelpoints); 500 PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("R","N",&ione,&blasn,&blask,mfqP->Q,&blasnpmax,mfqP->tau,mfqP->work2,&ione,mfqP->mwork,&blasm,&info)); 501 if (info < 0) SETERRQ1(PETSC_COMM_SELF,1,"ormqr returned value %d\n",info); 502 } 503 proj = BLASnrm2_(&blasj,&mfqP->work2[mfqP->nmodelpoints],&ione); 504 505 if (proj >= mfqP->theta1) { /* add this index to model */ 506 mfqP->model_indices[mfqP->nmodelpoints]=i; 507 mfqP->nmodelpoints++; 508 PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasn,mfqP->work,&ione,&mfqP->Q_tmp[mfqP->npmax*(mfqP->nmodelpoints-1)],&ione)); 509 blask=mfqP->npmax*(mfqP->nmodelpoints); 510 PetscStackCallBLAS("BLAScopy",BLAScopy_(&blask,mfqP->Q_tmp,&ione,mfqP->Q,&ione)); 511 blask = mfqP->nmodelpoints; 512 blasmaxmn = PetscMax(mfqP->m,mfqP->n); 513 PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&blasn,&blask,mfqP->Q,&blasnpmax,mfqP->tau,mfqP->mwork,&blasmaxmn,&info)); 514 if (info < 0) SETERRQ1(PETSC_COMM_SELF,1,"geqrf returned value %d\n",info); 515 mfqP->q_is_I = 0; 516 } 517 if (mfqP->nmodelpoints == mfqP->n) { 518 break; 519 } 520 } 521 } 522 523 PetscFunctionReturn(0); 524 } 525 526 #undef __FUNCT__ 527 #define __FUNCT__ "TaoSolve_POUNDERS" 528 static PetscErrorCode TaoSolve_POUNDERS(Tao tao) 529 { 530 TAO_POUNDERS *mfqP = (TAO_POUNDERS *)tao->data; 531 PetscInt i,ii,j,k,l; 532 PetscReal step=1.0; 533 TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 534 PetscInt low,high; 535 PetscReal minnorm; 536 PetscReal *x,*f; 537 const PetscReal *xmint,*fmin; 538 PetscReal cres,deltaold; 539 PetscReal gnorm; 540 PetscBLASInt info,ione=1,iblas; 541 PetscBool valid,same; 542 PetscReal mdec, rho, normxsp; 543 PetscReal one=1.0,zero=0.0,ratio; 544 PetscBLASInt blasm,blasn,blasnpmax,blasn2; 545 PetscErrorCode ierr; 546 static PetscBool set = PETSC_FALSE; 547 548 /* n = # of parameters 549 m = dimension (components) of function */ 550 PetscFunctionBegin; 551 ierr = PetscCitationsRegister("@article{UNEDF0,\n" 552 "title = {Nuclear energy density optimization},\n" 553 "author = {Kortelainen, M. and Lesinski, T. and Mor\'e, J. and Nazarewicz, W.\n" 554 " and Sarich, J. and Schunck, N. and Stoitsov, M. V. and Wild, S. },\n" 555 "journal = {Phys. Rev. C},\n" 556 "volume = {82},\n" 557 "number = {2},\n" 558 "pages = {024313},\n" 559 "numpages = {18},\n" 560 "year = {2010},\n" 561 "month = {Aug},\n" 562 "doi = {10.1103/PhysRevC.82.024313}\n}\n",&set);CHKERRQ(ierr); 563 tao->niter=0; 564 if (tao->XL && tao->XU) { 565 /* Check x0 <= XU */ 566 PetscReal val; 567 ierr = VecCopy(tao->solution,mfqP->Xhist[0]);CHKERRQ(ierr); 568 ierr = VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);CHKERRQ(ierr); 569 ierr = VecMax(mfqP->Xhist[0],NULL,&val);CHKERRQ(ierr); 570 if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 > upper bound"); 571 572 /* Check x0 >= xl */ 573 ierr = VecCopy(tao->XL,mfqP->Xhist[0]);CHKERRQ(ierr); 574 ierr = VecAXPY(mfqP->Xhist[0],-1.0,tao->solution);CHKERRQ(ierr); 575 ierr = VecMax(mfqP->Xhist[0],NULL,&val);CHKERRQ(ierr); 576 if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 < lower bound"); 577 578 /* Check x0 + delta < XU -- should be able to get around this eventually */ 579 580 ierr = VecSet(mfqP->Xhist[0],mfqP->delta);CHKERRQ(ierr); 581 ierr = VecAXPY(mfqP->Xhist[0],1.0,tao->solution);CHKERRQ(ierr); 582 ierr = VecAXPY(mfqP->Xhist[0],-1.0,tao->XU);CHKERRQ(ierr); 583 ierr = VecMax(mfqP->Xhist[0],NULL,&val);CHKERRQ(ierr); 584 if (val > 1e-10) SETERRQ(PETSC_COMM_WORLD,1,"X0 + delta > upper bound"); 585 } 586 587 CHKMEMQ; 588 blasm = mfqP->m; blasn=mfqP->n; blasnpmax = mfqP->npmax; 589 for (i=0;i<mfqP->n*mfqP->n*mfqP->m;i++) mfqP->H[i]=0; 590 591 ierr = VecCopy(tao->solution,mfqP->Xhist[0]);CHKERRQ(ierr); 592 CHKMEMQ; 593 ierr = TaoComputeSeparableObjective(tao,tao->solution,mfqP->Fhist[0]);CHKERRQ(ierr); 594 595 ierr = VecNorm(mfqP->Fhist[0],NORM_2,&mfqP->Fres[0]);CHKERRQ(ierr); 596 if (PetscIsInfOrNanReal(mfqP->Fres[0])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 597 mfqP->Fres[0]*=mfqP->Fres[0]; 598 mfqP->minindex = 0; 599 minnorm = mfqP->Fres[0]; 600 ierr = TaoMonitor(tao, tao->niter, minnorm, PETSC_INFINITY, 0.0, step, &reason);CHKERRQ(ierr); 601 tao->niter++; 602 603 ierr = VecGetOwnershipRange(mfqP->Xhist[0],&low,&high);CHKERRQ(ierr); 604 for (i=1;i<mfqP->n+1;i++) { 605 ierr = VecCopy(tao->solution,mfqP->Xhist[i]);CHKERRQ(ierr); 606 if (i-1 >= low && i-1 < high) { 607 ierr = VecGetArray(mfqP->Xhist[i],&x);CHKERRQ(ierr); 608 x[i-1-low] += mfqP->delta; 609 ierr = VecRestoreArray(mfqP->Xhist[i],&x);CHKERRQ(ierr); 610 } 611 CHKMEMQ; 612 ierr = TaoComputeSeparableObjective(tao,mfqP->Xhist[i],mfqP->Fhist[i]);CHKERRQ(ierr); 613 ierr = VecNorm(mfqP->Fhist[i],NORM_2,&mfqP->Fres[i]);CHKERRQ(ierr); 614 if (PetscIsInfOrNanReal(mfqP->Fres[i])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 615 mfqP->Fres[i]*=mfqP->Fres[i]; 616 if (mfqP->Fres[i] < minnorm) { 617 mfqP->minindex = i; 618 minnorm = mfqP->Fres[i]; 619 } 620 } 621 ierr = VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);CHKERRQ(ierr); 622 ierr = VecCopy(mfqP->Fhist[mfqP->minindex],tao->sep_objective);CHKERRQ(ierr); 623 /* Gather mpi vecs to one big local vec */ 624 625 /* Begin serial code */ 626 627 /* Disp[i] = Xi-xmin, i=1,..,mfqP->minindex-1,mfqP->minindex+1,..,n */ 628 /* Fdiff[i] = (Fi-Fmin)', i=1,..,mfqP->minindex-1,mfqP->minindex+1,..,n */ 629 /* (Column oriented for blas calls) */ 630 ii=0; 631 632 if (mfqP->size == 1) { 633 ierr = VecGetArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);CHKERRQ(ierr); 634 for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i]; 635 ierr = VecRestoreArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);CHKERRQ(ierr); 636 ierr = VecGetArrayRead(mfqP->Fhist[mfqP->minindex],&fmin);CHKERRQ(ierr); 637 for (i=0;i<mfqP->n+1;i++) { 638 if (i == mfqP->minindex) continue; 639 640 ierr = VecGetArray(mfqP->Xhist[i],&x);CHKERRQ(ierr); 641 for (j=0;j<mfqP->n;j++) { 642 mfqP->Disp[ii+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta; 643 } 644 ierr = VecRestoreArray(mfqP->Xhist[i],&x);CHKERRQ(ierr); 645 646 ierr = VecGetArray(mfqP->Fhist[i],&f);CHKERRQ(ierr); 647 for (j=0;j<mfqP->m;j++) { 648 mfqP->Fdiff[ii+mfqP->n*j] = f[j] - fmin[j]; 649 } 650 ierr = VecRestoreArray(mfqP->Fhist[i],&f);CHKERRQ(ierr); 651 mfqP->model_indices[ii++] = i; 652 653 } 654 for (j=0;j<mfqP->m;j++) { 655 mfqP->C[j] = fmin[j]; 656 } 657 ierr = VecRestoreArrayRead(mfqP->Fhist[mfqP->minindex],&fmin);CHKERRQ(ierr); 658 } else { 659 ierr = VecScatterBegin(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 660 ierr = VecScatterEnd(mfqP->scatterx,mfqP->Xhist[mfqP->minindex],mfqP->localxmin,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 661 ierr = VecGetArrayRead(mfqP->localxmin,&xmint);CHKERRQ(ierr); 662 for (i=0;i<mfqP->n;i++) mfqP->xmin[i] = xmint[i]; 663 ierr = VecRestoreArrayRead(mfqP->localxmin,&xmint);CHKERRQ(ierr); 664 665 ierr = VecScatterBegin(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 666 ierr = VecScatterEnd(mfqP->scatterf,mfqP->Fhist[mfqP->minindex],mfqP->localfmin,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 667 ierr = VecGetArrayRead(mfqP->localfmin,&fmin);CHKERRQ(ierr); 668 for (i=0;i<mfqP->n+1;i++) { 669 if (i == mfqP->minindex) continue; 670 671 mfqP->model_indices[ii++] = i; 672 ierr = VecScatterBegin(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 673 ierr = VecScatterEnd(mfqP->scatterx,mfqP->Xhist[ii],mfqP->localx,INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 674 ierr = VecGetArray(mfqP->localx,&x);CHKERRQ(ierr); 675 for (j=0;j<mfqP->n;j++) { 676 mfqP->Disp[i+mfqP->npmax*j] = (x[j] - mfqP->xmin[j])/mfqP->delta; 677 } 678 ierr = VecRestoreArray(mfqP->localx,&x);CHKERRQ(ierr); 679 680 ierr = VecScatterBegin(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 681 ierr = VecScatterEnd(mfqP->scatterf,mfqP->Fhist[ii],mfqP->localf,INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 682 ierr = VecGetArray(mfqP->localf,&f);CHKERRQ(ierr); 683 for (j=0;j<mfqP->m;j++) { 684 mfqP->Fdiff[i*mfqP->n+j] = f[j] - fmin[j]; 685 } 686 ierr = VecRestoreArray(mfqP->localf,&f);CHKERRQ(ierr); 687 } 688 for (j=0;j<mfqP->m;j++) { 689 mfqP->C[j] = fmin[j]; 690 } 691 ierr = VecRestoreArrayRead(mfqP->localfmin,&fmin);CHKERRQ(ierr); 692 } 693 694 /* Determine the initial quadratic models */ 695 /* G = D(ModelIn,:) \ (F(ModelIn,1:m)-repmat(F(xkin,1:m),n,1)); */ 696 /* D (nxn) Fdiff (nxm) => G (nxm) */ 697 blasn2 = blasn; 698 PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&blasn,&blasm,mfqP->Disp,&blasnpmax,mfqP->iwork,mfqP->Fdiff,&blasn2,&info)); 699 ierr = PetscInfo1(tao,"gesv returned %d\n",info);CHKERRQ(ierr); 700 701 cres = minnorm; 702 /* Gres = G*F(xkin,1:m)' G (nxm) Fk (m) */ 703 PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->C,&ione,&zero,mfqP->Gres,&ione)); 704 705 /* Hres = G*G' */ 706 PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&blasn,&blasn,&blasm,&one,mfqP->Fdiff, &blasn,mfqP->Fdiff,&blasn,&zero,mfqP->Hres,&blasn)); 707 708 valid = PETSC_TRUE; 709 710 ierr = VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);CHKERRQ(ierr); 711 ierr = VecAssemblyBegin(tao->gradient);CHKERRQ(ierr); 712 ierr = VecAssemblyEnd(tao->gradient);CHKERRQ(ierr); 713 ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 714 gnorm *= mfqP->delta; 715 ierr = VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);CHKERRQ(ierr); 716 ierr = TaoMonitor(tao, tao->niter, minnorm, gnorm, 0.0, step, &reason);CHKERRQ(ierr); 717 mfqP->nHist = mfqP->n+1; 718 mfqP->nmodelpoints = mfqP->n+1; 719 720 while (reason == TAO_CONTINUE_ITERATING) { 721 PetscReal gnm = 1e-4; 722 tao->niter++; 723 /* Solve the subproblem min{Q(s): ||s|| <= delta} */ 724 ierr = gqtwrap(tao,&gnm,&mdec);CHKERRQ(ierr); 725 /* Evaluate the function at the new point */ 726 727 for (i=0;i<mfqP->n;i++) { 728 mfqP->work[i] = mfqP->Xsubproblem[i]*mfqP->delta + mfqP->xmin[i]; 729 } 730 ierr = VecDuplicate(tao->solution,&mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr); 731 ierr = VecDuplicate(tao->sep_objective,&mfqP->Fhist[mfqP->nHist]);CHKERRQ(ierr); 732 ierr = VecSetValues(mfqP->Xhist[mfqP->nHist],mfqP->n,mfqP->indices,mfqP->work,INSERT_VALUES);CHKERRQ(ierr); 733 ierr = VecAssemblyBegin(mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr); 734 ierr = VecAssemblyEnd(mfqP->Xhist[mfqP->nHist]);CHKERRQ(ierr); 735 736 ierr = TaoComputeSeparableObjective(tao,mfqP->Xhist[mfqP->nHist],mfqP->Fhist[mfqP->nHist]);CHKERRQ(ierr); 737 ierr = VecNorm(mfqP->Fhist[mfqP->nHist],NORM_2,&mfqP->Fres[mfqP->nHist]);CHKERRQ(ierr); 738 if (PetscIsInfOrNanReal(mfqP->Fres[mfqP->nHist])) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 739 mfqP->Fres[mfqP->nHist]*=mfqP->Fres[mfqP->nHist]; 740 rho = (mfqP->Fres[mfqP->minindex] - mfqP->Fres[mfqP->nHist]) / mdec; 741 mfqP->nHist++; 742 743 /* Update the center */ 744 if ((rho >= mfqP->eta1) || (rho > mfqP->eta0 && valid==PETSC_TRUE)) { 745 /* Update model to reflect new base point */ 746 for (i=0;i<mfqP->n;i++) { 747 mfqP->work[i] = (mfqP->work[i] - mfqP->xmin[i])/mfqP->delta; 748 } 749 for (j=0;j<mfqP->m;j++) { 750 /* C(j) = C(j) + work*G(:,j) + .5*work*H(:,:,j)*work'; 751 G(:,j) = G(:,j) + H(:,:,j)*work' */ 752 for (k=0;k<mfqP->n;k++) { 753 mfqP->work2[k]=0.0; 754 for (l=0;l<mfqP->n;l++) { 755 mfqP->work2[k]+=mfqP->H[j + mfqP->m*(k + l*mfqP->n)]*mfqP->work[l]; 756 } 757 } 758 for (i=0;i<mfqP->n;i++) { 759 mfqP->C[j]+=mfqP->work[i]*(mfqP->Fdiff[i + mfqP->n* j] + 0.5*mfqP->work2[i]); 760 mfqP->Fdiff[i+mfqP->n*j] +=mfqP-> work2[i]; 761 } 762 } 763 /* Cres += work*Gres + .5*work*Hres*work'; 764 Gres += Hres*work'; */ 765 766 PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasn,&one,mfqP->Hres,&blasn,mfqP->work,&ione,&zero,mfqP->work2,&ione)); 767 for (i=0;i<mfqP->n;i++) { 768 cres += mfqP->work[i]*(mfqP->Gres[i] + 0.5*mfqP->work2[i]); 769 mfqP->Gres[i] += mfqP->work2[i]; 770 } 771 mfqP->minindex = mfqP->nHist-1; 772 minnorm = mfqP->Fres[mfqP->minindex]; 773 ierr = VecCopy(mfqP->Fhist[mfqP->minindex],tao->sep_objective);CHKERRQ(ierr); 774 /* Change current center */ 775 ierr = VecGetArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);CHKERRQ(ierr); 776 for (i=0;i<mfqP->n;i++) { 777 mfqP->xmin[i] = xmint[i]; 778 } 779 ierr = VecRestoreArrayRead(mfqP->Xhist[mfqP->minindex],&xmint);CHKERRQ(ierr); 780 } 781 782 /* Evaluate at a model-improving point if necessary */ 783 if (valid == PETSC_FALSE) { 784 mfqP->q_is_I = 1; 785 mfqP->nmodelpoints = 0; 786 ierr = affpoints(mfqP,mfqP->xmin,mfqP->c1);CHKERRQ(ierr); 787 if (mfqP->nmodelpoints < mfqP->n) { 788 ierr = PetscInfo(tao,"Model not valid -- model-improving\n");CHKERRQ(ierr); 789 ierr = modelimprove(tao,mfqP,1);CHKERRQ(ierr); 790 } 791 } 792 793 /* Update the trust region radius */ 794 deltaold = mfqP->delta; 795 normxsp = 0; 796 for (i=0;i<mfqP->n;i++) { 797 normxsp += mfqP->Xsubproblem[i]*mfqP->Xsubproblem[i]; 798 } 799 normxsp = PetscSqrtReal(normxsp); 800 if (rho >= mfqP->eta1 && normxsp > 0.5*mfqP->delta) { 801 mfqP->delta = PetscMin(mfqP->delta*mfqP->gamma1,mfqP->deltamax); 802 } else if (valid == PETSC_TRUE) { 803 mfqP->delta = PetscMax(mfqP->delta*mfqP->gamma0,mfqP->deltamin); 804 } 805 806 /* Compute the next interpolation set */ 807 mfqP->q_is_I = 1; 808 mfqP->nmodelpoints=0; 809 ierr = affpoints(mfqP,mfqP->xmin,mfqP->c1);CHKERRQ(ierr); 810 if (mfqP->nmodelpoints == mfqP->n) { 811 valid = PETSC_TRUE; 812 } else { 813 valid = PETSC_FALSE; 814 ierr = affpoints(mfqP,mfqP->xmin,mfqP->c2);CHKERRQ(ierr); 815 if (mfqP->n > mfqP->nmodelpoints) { 816 ierr = PetscInfo(tao,"Model not valid -- adding geometry points\n");CHKERRQ(ierr); 817 ierr = modelimprove(tao,mfqP,mfqP->n - mfqP->nmodelpoints);CHKERRQ(ierr); 818 } 819 } 820 for (i=mfqP->nmodelpoints;i>0;i--) { 821 mfqP->model_indices[i] = mfqP->model_indices[i-1]; 822 } 823 mfqP->nmodelpoints++; 824 mfqP->model_indices[0] = mfqP->minindex; 825 ierr = morepoints(mfqP);CHKERRQ(ierr); 826 for (i=0;i<mfqP->nmodelpoints;i++) { 827 ierr = VecGetArray(mfqP->Xhist[mfqP->model_indices[i]],&x);CHKERRQ(ierr); 828 for (j=0;j<mfqP->n;j++) { 829 mfqP->Disp[i + mfqP->npmax*j] = (x[j] - mfqP->xmin[j]) / deltaold; 830 } 831 ierr = VecRestoreArray(mfqP->Xhist[mfqP->model_indices[i]],&x);CHKERRQ(ierr); 832 ierr = VecGetArray(mfqP->Fhist[mfqP->model_indices[i]],&f);CHKERRQ(ierr); 833 for (j=0;j<mfqP->m;j++) { 834 for (k=0;k<mfqP->n;k++) { 835 mfqP->work[k]=0.0; 836 for (l=0;l<mfqP->n;l++) { 837 mfqP->work[k] += mfqP->H[j + mfqP->m*(k + mfqP->n*l)] * mfqP->Disp[i + mfqP->npmax*l]; 838 } 839 } 840 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]; 841 } 842 ierr = VecRestoreArray(mfqP->Fhist[mfqP->model_indices[i]],&f);CHKERRQ(ierr); 843 } 844 845 /* Update the quadratic model */ 846 ierr = getquadpounders(mfqP);CHKERRQ(ierr); 847 ierr = VecGetArrayRead(mfqP->Fhist[mfqP->minindex],&fmin);CHKERRQ(ierr); 848 PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasm,fmin,&ione,mfqP->C,&ione)); 849 /* G = G*(delta/deltaold) + Gdel */ 850 ratio = mfqP->delta/deltaold; 851 iblas = blasm*blasn; 852 PetscStackCallBLAS("BLASscal",BLASscal_(&iblas,&ratio,mfqP->Fdiff,&ione)); 853 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&iblas,&one,mfqP->Gdel,&ione,mfqP->Fdiff,&ione)); 854 /* H = H*(delta/deltaold) + Hdel */ 855 iblas = blasm*blasn*blasn; 856 ratio *= ratio; 857 PetscStackCallBLAS("BLASscal",BLASscal_(&iblas,&ratio,mfqP->H,&ione)); 858 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&iblas,&one,mfqP->Hdel,&ione,mfqP->H,&ione)); 859 860 /* Get residuals */ 861 cres = mfqP->Fres[mfqP->minindex]; 862 /* Gres = G*F(xkin,1:m)' */ 863 PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->C,&ione,&zero,mfqP->Gres,&ione)); 864 /* Hres = sum i=1..m {F(xkin,i)*H(:,:,i)} + G*G' */ 865 PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&blasn,&blasn,&blasm,&one,mfqP->Fdiff,&blasn,mfqP->Fdiff,&blasn,&zero,mfqP->Hres,&blasn)); 866 867 iblas = mfqP->n*mfqP->n; 868 869 for (j=0;j<mfqP->m;j++) { 870 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&iblas,&fmin[j],&mfqP->H[j],&blasm,mfqP->Hres,&ione)); 871 } 872 873 /* Export solution and gradient residual to TAO */ 874 ierr = VecCopy(mfqP->Xhist[mfqP->minindex],tao->solution);CHKERRQ(ierr); 875 ierr = VecSetValues(tao->gradient,mfqP->n,mfqP->indices,mfqP->Gres,INSERT_VALUES);CHKERRQ(ierr); 876 ierr = VecAssemblyBegin(tao->gradient);CHKERRQ(ierr); 877 ierr = VecAssemblyEnd(tao->gradient);CHKERRQ(ierr); 878 ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 879 gnorm *= mfqP->delta; 880 /* final criticality test */ 881 ierr = TaoMonitor(tao, tao->niter, minnorm, gnorm, 0.0, step, &reason);CHKERRQ(ierr); 882 /* test for repeated model */ 883 if (mfqP->nmodelpoints==mfqP->last_nmodelpoints) { 884 same = PETSC_TRUE; 885 } else { 886 same = PETSC_FALSE; 887 } 888 for (i=0;i<mfqP->nmodelpoints;i++) { 889 if (same) { 890 if (mfqP->model_indices[i] == mfqP->last_model_indices[i]) { 891 same = PETSC_TRUE; 892 } else { 893 same = PETSC_FALSE; 894 } 895 } 896 mfqP->last_model_indices[i] = mfqP->model_indices[i]; 897 } 898 mfqP->last_nmodelpoints = mfqP->nmodelpoints; 899 if (same && mfqP->delta == deltaold) { 900 ierr = PetscInfo(tao,"Identical model used in successive iterations\n");CHKERRQ(ierr); 901 reason = TAO_CONVERGED_STEPTOL; 902 tao->reason = TAO_CONVERGED_STEPTOL; 903 } 904 } 905 PetscFunctionReturn(0); 906 } 907 908 #undef __FUNCT__ 909 #define __FUNCT__ "TaoSetUp_POUNDERS" 910 static PetscErrorCode TaoSetUp_POUNDERS(Tao tao) 911 { 912 TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data; 913 PetscInt i; 914 IS isfloc,isfglob,isxloc,isxglob; 915 PetscErrorCode ierr; 916 917 PetscFunctionBegin; 918 if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); } 919 if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); } 920 ierr = VecGetSize(tao->solution,&mfqP->n);CHKERRQ(ierr); 921 ierr = VecGetSize(tao->sep_objective,&mfqP->m);CHKERRQ(ierr); 922 mfqP->c1 = PetscSqrtReal((PetscReal)mfqP->n); 923 if (mfqP->npmax == PETSC_DEFAULT) { 924 mfqP->npmax = 2*mfqP->n + 1; 925 } 926 mfqP->npmax = PetscMin((mfqP->n+1)*(mfqP->n+2)/2,mfqP->npmax); 927 mfqP->npmax = PetscMax(mfqP->npmax, mfqP->n+2); 928 929 ierr = PetscMalloc1(tao->max_funcs+10,&mfqP->Xhist);CHKERRQ(ierr); 930 ierr = PetscMalloc1(tao->max_funcs+10,&mfqP->Fhist);CHKERRQ(ierr); 931 for (i=0;i<mfqP->n +1;i++) { 932 ierr = VecDuplicate(tao->solution,&mfqP->Xhist[i]);CHKERRQ(ierr); 933 ierr = VecDuplicate(tao->sep_objective,&mfqP->Fhist[i]);CHKERRQ(ierr); 934 } 935 ierr = VecDuplicate(tao->solution,&mfqP->workxvec);CHKERRQ(ierr); 936 mfqP->nHist = 0; 937 938 ierr = PetscMalloc1(tao->max_funcs+10,&mfqP->Fres);CHKERRQ(ierr); 939 ierr = PetscMalloc1(mfqP->npmax*mfqP->m,&mfqP->RES);CHKERRQ(ierr); 940 ierr = PetscMalloc1(mfqP->n,&mfqP->work);CHKERRQ(ierr); 941 ierr = PetscMalloc1(mfqP->n,&mfqP->work2);CHKERRQ(ierr); 942 ierr = PetscMalloc1(mfqP->n,&mfqP->work3);CHKERRQ(ierr); 943 ierr = PetscMalloc1(PetscMax(mfqP->m,mfqP->n+1),&mfqP->mwork);CHKERRQ(ierr); 944 ierr = PetscMalloc1(mfqP->npmax - mfqP->n - 1,&mfqP->omega);CHKERRQ(ierr); 945 ierr = PetscMalloc1(mfqP->n * (mfqP->n+1) / 2,&mfqP->beta);CHKERRQ(ierr); 946 ierr = PetscMalloc1(mfqP->n + 1 ,&mfqP->alpha);CHKERRQ(ierr); 947 948 ierr = PetscMalloc1(mfqP->n*mfqP->n*mfqP->m,&mfqP->H);CHKERRQ(ierr); 949 ierr = PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q);CHKERRQ(ierr); 950 ierr = PetscMalloc1(mfqP->npmax*mfqP->npmax,&mfqP->Q_tmp);CHKERRQ(ierr); 951 ierr = PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L);CHKERRQ(ierr); 952 ierr = PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_tmp);CHKERRQ(ierr); 953 ierr = PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->L_save);CHKERRQ(ierr); 954 ierr = PetscMalloc1(mfqP->n*(mfqP->n+1)/2*(mfqP->npmax),&mfqP->N);CHKERRQ(ierr); 955 ierr = PetscMalloc1(mfqP->npmax * (mfqP->n+1) ,&mfqP->M);CHKERRQ(ierr); 956 ierr = PetscMalloc1(mfqP->npmax * (mfqP->npmax - mfqP->n - 1) , &mfqP->Z);CHKERRQ(ierr); 957 ierr = PetscMalloc1(mfqP->npmax,&mfqP->tau);CHKERRQ(ierr); 958 ierr = PetscMalloc1(mfqP->npmax,&mfqP->tau_tmp);CHKERRQ(ierr); 959 mfqP->nmax = PetscMax(5*mfqP->npmax,mfqP->n*(mfqP->n+1)/2); 960 ierr = PetscMalloc1(mfqP->nmax,&mfqP->npmaxwork);CHKERRQ(ierr); 961 ierr = PetscMalloc1(mfqP->nmax,&mfqP->npmaxiwork);CHKERRQ(ierr); 962 ierr = PetscMalloc1(mfqP->n,&mfqP->xmin);CHKERRQ(ierr); 963 ierr = PetscMalloc1(mfqP->m,&mfqP->C);CHKERRQ(ierr); 964 ierr = PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Fdiff);CHKERRQ(ierr); 965 ierr = PetscMalloc1(mfqP->npmax*mfqP->n,&mfqP->Disp);CHKERRQ(ierr); 966 ierr = PetscMalloc1(mfqP->n,&mfqP->Gres);CHKERRQ(ierr); 967 ierr = PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Hres);CHKERRQ(ierr); 968 ierr = PetscMalloc1(mfqP->n*mfqP->n,&mfqP->Gpoints);CHKERRQ(ierr); 969 ierr = PetscMalloc1(mfqP->npmax,&mfqP->model_indices);CHKERRQ(ierr); 970 ierr = PetscMalloc1(mfqP->npmax,&mfqP->last_model_indices);CHKERRQ(ierr); 971 ierr = PetscMalloc1(mfqP->n,&mfqP->Xsubproblem);CHKERRQ(ierr); 972 ierr = PetscMalloc1(mfqP->m*mfqP->n,&mfqP->Gdel);CHKERRQ(ierr); 973 ierr = PetscMalloc1(mfqP->n*mfqP->n*mfqP->m, &mfqP->Hdel);CHKERRQ(ierr); 974 ierr = PetscMalloc1(PetscMax(mfqP->m,mfqP->n),&mfqP->indices);CHKERRQ(ierr); 975 ierr = PetscMalloc1(mfqP->n,&mfqP->iwork);CHKERRQ(ierr); 976 for (i=0;i<PetscMax(mfqP->m,mfqP->n);i++) { 977 mfqP->indices[i] = i; 978 } 979 ierr = MPI_Comm_size(((PetscObject)tao)->comm,&mfqP->size);CHKERRQ(ierr); 980 if (mfqP->size > 1) { 981 VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localx);CHKERRQ(ierr); 982 VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->localxmin);CHKERRQ(ierr); 983 VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localf);CHKERRQ(ierr); 984 VecCreateSeq(PETSC_COMM_SELF,mfqP->m,&mfqP->localfmin);CHKERRQ(ierr); 985 ierr = ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxloc);CHKERRQ(ierr); 986 ierr = ISCreateStride(MPI_COMM_SELF,mfqP->n,0,1,&isxglob);CHKERRQ(ierr); 987 ierr = ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfloc);CHKERRQ(ierr); 988 ierr = ISCreateStride(MPI_COMM_SELF,mfqP->m,0,1,&isfglob);CHKERRQ(ierr); 989 990 991 ierr = VecScatterCreate(tao->solution,isxglob,mfqP->localx,isxloc,&mfqP->scatterx);CHKERRQ(ierr); 992 ierr = VecScatterCreate(tao->sep_objective,isfglob,mfqP->localf,isfloc,&mfqP->scatterf);CHKERRQ(ierr); 993 994 ierr = ISDestroy(&isxloc);CHKERRQ(ierr); 995 ierr = ISDestroy(&isxglob);CHKERRQ(ierr); 996 ierr = ISDestroy(&isfloc);CHKERRQ(ierr); 997 ierr = ISDestroy(&isfglob);CHKERRQ(ierr); 998 } 999 1000 if (!mfqP->usegqt) { 1001 KSP ksp; 1002 PC pc; 1003 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,mfqP->n,mfqP->n,mfqP->Xsubproblem,&mfqP->subx);CHKERRQ(ierr); 1004 ierr = VecCreateSeq(PETSC_COMM_SELF,mfqP->n,&mfqP->subxl);CHKERRQ(ierr); 1005 ierr = VecDuplicate(mfqP->subxl,&mfqP->subb);CHKERRQ(ierr); 1006 ierr = VecDuplicate(mfqP->subxl,&mfqP->subxu);CHKERRQ(ierr); 1007 ierr = VecDuplicate(mfqP->subxl,&mfqP->subpdel);CHKERRQ(ierr); 1008 ierr = VecDuplicate(mfqP->subxl,&mfqP->subndel);CHKERRQ(ierr); 1009 ierr = TaoCreate(PETSC_COMM_SELF,&mfqP->subtao);CHKERRQ(ierr); 1010 ierr = TaoSetType(mfqP->subtao,TAOTRON);CHKERRQ(ierr); 1011 ierr = TaoSetOptionsPrefix(mfqP->subtao,"pounders_subsolver_");CHKERRQ(ierr); 1012 ierr = TaoSetInitialVector(mfqP->subtao,mfqP->subx);CHKERRQ(ierr); 1013 ierr = TaoSetObjectiveAndGradientRoutine(mfqP->subtao,pounders_fg,(void*)mfqP);CHKERRQ(ierr); 1014 ierr = TaoSetMaximumIterations(mfqP->subtao,mfqP->gqt_maxits);CHKERRQ(ierr); 1015 ierr = TaoSetFromOptions(mfqP->subtao);CHKERRQ(ierr); 1016 ierr = TaoGetKSP(mfqP->subtao,&ksp);CHKERRQ(ierr); 1017 if (ksp) { 1018 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1019 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1020 } 1021 ierr = TaoSetVariableBounds(mfqP->subtao,mfqP->subxl,mfqP->subxu);CHKERRQ(ierr); 1022 ierr = MatCreateSeqDense(PETSC_COMM_SELF,mfqP->n,mfqP->n,mfqP->Hres,&mfqP->subH);CHKERRQ(ierr); 1023 ierr = TaoSetHessianRoutine(mfqP->subtao,mfqP->subH,mfqP->subH,pounders_h,(void*)mfqP);CHKERRQ(ierr); 1024 } 1025 PetscFunctionReturn(0); 1026 } 1027 1028 #undef __FUNCT__ 1029 #define __FUNCT__ "TaoDestroy_POUNDERS" 1030 static PetscErrorCode TaoDestroy_POUNDERS(Tao tao) 1031 { 1032 TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data; 1033 PetscInt i; 1034 PetscErrorCode ierr; 1035 1036 PetscFunctionBegin; 1037 if (!mfqP->usegqt) { 1038 ierr = TaoDestroy(&mfqP->subtao);CHKERRQ(ierr); 1039 ierr = VecDestroy(&mfqP->subx);CHKERRQ(ierr); 1040 ierr = VecDestroy(&mfqP->subxl);CHKERRQ(ierr); 1041 ierr = VecDestroy(&mfqP->subxu);CHKERRQ(ierr); 1042 ierr = VecDestroy(&mfqP->subb);CHKERRQ(ierr); 1043 ierr = VecDestroy(&mfqP->subpdel);CHKERRQ(ierr); 1044 ierr = VecDestroy(&mfqP->subndel);CHKERRQ(ierr); 1045 ierr = MatDestroy(&mfqP->subH);CHKERRQ(ierr); 1046 } 1047 ierr = PetscFree(mfqP->Fres);CHKERRQ(ierr); 1048 ierr = PetscFree(mfqP->RES);CHKERRQ(ierr); 1049 ierr = PetscFree(mfqP->work);CHKERRQ(ierr); 1050 ierr = PetscFree(mfqP->work2);CHKERRQ(ierr); 1051 ierr = PetscFree(mfqP->work3);CHKERRQ(ierr); 1052 ierr = PetscFree(mfqP->mwork);CHKERRQ(ierr); 1053 ierr = PetscFree(mfqP->omega);CHKERRQ(ierr); 1054 ierr = PetscFree(mfqP->beta);CHKERRQ(ierr); 1055 ierr = PetscFree(mfqP->alpha);CHKERRQ(ierr); 1056 ierr = PetscFree(mfqP->H);CHKERRQ(ierr); 1057 ierr = PetscFree(mfqP->Q);CHKERRQ(ierr); 1058 ierr = PetscFree(mfqP->Q_tmp);CHKERRQ(ierr); 1059 ierr = PetscFree(mfqP->L);CHKERRQ(ierr); 1060 ierr = PetscFree(mfqP->L_tmp);CHKERRQ(ierr); 1061 ierr = PetscFree(mfqP->L_save);CHKERRQ(ierr); 1062 ierr = PetscFree(mfqP->N);CHKERRQ(ierr); 1063 ierr = PetscFree(mfqP->M);CHKERRQ(ierr); 1064 ierr = PetscFree(mfqP->Z);CHKERRQ(ierr); 1065 ierr = PetscFree(mfqP->tau);CHKERRQ(ierr); 1066 ierr = PetscFree(mfqP->tau_tmp);CHKERRQ(ierr); 1067 ierr = PetscFree(mfqP->npmaxwork);CHKERRQ(ierr); 1068 ierr = PetscFree(mfqP->npmaxiwork);CHKERRQ(ierr); 1069 ierr = PetscFree(mfqP->xmin);CHKERRQ(ierr); 1070 ierr = PetscFree(mfqP->C);CHKERRQ(ierr); 1071 ierr = PetscFree(mfqP->Fdiff);CHKERRQ(ierr); 1072 ierr = PetscFree(mfqP->Disp);CHKERRQ(ierr); 1073 ierr = PetscFree(mfqP->Gres);CHKERRQ(ierr); 1074 ierr = PetscFree(mfqP->Hres);CHKERRQ(ierr); 1075 ierr = PetscFree(mfqP->Gpoints);CHKERRQ(ierr); 1076 ierr = PetscFree(mfqP->model_indices);CHKERRQ(ierr); 1077 ierr = PetscFree(mfqP->last_model_indices);CHKERRQ(ierr); 1078 ierr = PetscFree(mfqP->Xsubproblem);CHKERRQ(ierr); 1079 ierr = PetscFree(mfqP->Gdel);CHKERRQ(ierr); 1080 ierr = PetscFree(mfqP->Hdel);CHKERRQ(ierr); 1081 ierr = PetscFree(mfqP->indices);CHKERRQ(ierr); 1082 ierr = PetscFree(mfqP->iwork);CHKERRQ(ierr); 1083 1084 for (i=0;i<mfqP->nHist;i++) { 1085 ierr = VecDestroy(&mfqP->Xhist[i]);CHKERRQ(ierr); 1086 ierr = VecDestroy(&mfqP->Fhist[i]);CHKERRQ(ierr); 1087 } 1088 if (mfqP->workxvec) { 1089 ierr = VecDestroy(&mfqP->workxvec);CHKERRQ(ierr); 1090 } 1091 ierr = PetscFree(mfqP->Xhist);CHKERRQ(ierr); 1092 ierr = PetscFree(mfqP->Fhist);CHKERRQ(ierr); 1093 1094 if (mfqP->size > 1) { 1095 ierr = VecDestroy(&mfqP->localx);CHKERRQ(ierr); 1096 ierr = VecDestroy(&mfqP->localxmin);CHKERRQ(ierr); 1097 ierr = VecDestroy(&mfqP->localf);CHKERRQ(ierr); 1098 ierr = VecDestroy(&mfqP->localfmin);CHKERRQ(ierr); 1099 } 1100 ierr = PetscFree(tao->data);CHKERRQ(ierr); 1101 PetscFunctionReturn(0); 1102 } 1103 1104 #undef __FUNCT__ 1105 #define __FUNCT__ "TaoSetFromOptions_POUNDERS" 1106 static PetscErrorCode TaoSetFromOptions_POUNDERS(PetscOptionItems *PetscOptionsObject,Tao tao) 1107 { 1108 TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data; 1109 PetscErrorCode ierr; 1110 1111 PetscFunctionBegin; 1112 ierr = PetscOptionsHead(PetscOptionsObject,"POUNDERS method for least-squares optimization");CHKERRQ(ierr); 1113 ierr = PetscOptionsReal("-tao_pounders_delta","initial delta","",mfqP->delta,&mfqP->delta0,NULL);CHKERRQ(ierr); 1114 mfqP->delta = mfqP->delta0; 1115 mfqP->npmax = PETSC_DEFAULT; 1116 ierr = PetscOptionsInt("-tao_pounders_npmax","max number of points in model","",mfqP->npmax,&mfqP->npmax,NULL);CHKERRQ(ierr); 1117 mfqP->usegqt = PETSC_FALSE; 1118 ierr = PetscOptionsBool("-tao_pounders_gqt","use gqt algorithm for subproblem","",mfqP->usegqt,&mfqP->usegqt,NULL);CHKERRQ(ierr); 1119 ierr = PetscOptionsTail();CHKERRQ(ierr); 1120 PetscFunctionReturn(0); 1121 } 1122 1123 #undef __FUNCT__ 1124 #define __FUNCT__ "TaoView_POUNDERS" 1125 static PetscErrorCode TaoView_POUNDERS(Tao tao, PetscViewer viewer) 1126 { 1127 TAO_POUNDERS *mfqP = (TAO_POUNDERS *)tao->data; 1128 PetscBool isascii; 1129 PetscInt nits; 1130 PetscErrorCode ierr; 1131 1132 PetscFunctionBegin; 1133 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 1134 if (isascii) { 1135 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1136 ierr = PetscViewerASCIIPrintf(viewer, "initial delta: %g\n",(double)mfqP->delta0);CHKERRQ(ierr); 1137 ierr = PetscViewerASCIIPrintf(viewer, "final delta: %g\n",(double)mfqP->delta);CHKERRQ(ierr); 1138 ierr = PetscViewerASCIIPrintf(viewer, "model points: %D\n",mfqP->nmodelpoints);CHKERRQ(ierr); 1139 if (mfqP->usegqt) { 1140 ierr = PetscViewerASCIIPrintf(viewer, "subproblem solver: gqt\n");CHKERRQ(ierr); 1141 } else { 1142 ierr = PetscViewerASCIIPrintf(viewer, "subproblem solver: %s\n",((PetscObject)mfqP->subtao)->type_name);CHKERRQ(ierr); 1143 ierr = TaoGetTotalIterationNumber(mfqP->subtao,&nits);CHKERRQ(ierr); 1144 ierr = PetscViewerASCIIPrintf(viewer, "total subproblem iterations: %D\n",nits);CHKERRQ(ierr); 1145 } 1146 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1147 } 1148 PetscFunctionReturn(0); 1149 } 1150 /*MC 1151 TAOPOUNDERS - POUNDERS derivate-free model-based algorithm for nonlinear least squares 1152 1153 Options Database Keys: 1154 + -tao_pounders_delta - initial step length 1155 . -tao_pounders_npmax - maximum number of points in model 1156 - -tao_pounders_gqt - use gqt algorithm for subproblem instead of TRON 1157 1158 Level: beginner 1159 1160 M*/ 1161 1162 #undef __FUNCT__ 1163 #define __FUNCT__ "TaoCreate_POUNDERS" 1164 PETSC_EXTERN PetscErrorCode TaoCreate_POUNDERS(Tao tao) 1165 { 1166 TAO_POUNDERS *mfqP = (TAO_POUNDERS*)tao->data; 1167 PetscErrorCode ierr; 1168 1169 PetscFunctionBegin; 1170 tao->ops->setup = TaoSetUp_POUNDERS; 1171 tao->ops->solve = TaoSolve_POUNDERS; 1172 tao->ops->view = TaoView_POUNDERS; 1173 tao->ops->setfromoptions = TaoSetFromOptions_POUNDERS; 1174 tao->ops->destroy = TaoDestroy_POUNDERS; 1175 1176 ierr = PetscNewLog(tao,&mfqP);CHKERRQ(ierr); 1177 tao->data = (void*)mfqP; 1178 /* Override default settings (unless already changed) */ 1179 if (!tao->max_it_changed) tao->max_it = 2000; 1180 if (!tao->max_funcs_changed) tao->max_funcs = 4000; 1181 mfqP->delta0 = 0.1; 1182 mfqP->delta = 0.1; 1183 mfqP->deltamax=1e3; 1184 mfqP->c2 = 100.0; 1185 mfqP->theta1=1e-5; 1186 mfqP->theta2=1e-4; 1187 mfqP->gamma0=0.5; 1188 mfqP->gamma1=2.0; 1189 mfqP->eta0 = 0.0; 1190 mfqP->eta1 = 0.1; 1191 mfqP->gqt_rtol = 0.001; 1192 mfqP->gqt_maxits = 50; 1193 mfqP->workxvec = 0; 1194 PetscFunctionReturn(0); 1195 } 1196 1197 1198