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