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