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