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