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