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