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