1 #include <petscsys.h> 2 #include <petscblaslapack.h> 3 4 static PetscErrorCode estsv(PetscInt n, PetscReal *r, PetscInt ldr, PetscReal *svmin, PetscReal *z) 5 { 6 PetscBLASInt blas1=1, blasn=n, blasnmi, blasj, blasldr = ldr; 7 PetscInt i,j; 8 PetscReal e,temp,w,wm,ynorm,znorm,s,sm; 9 10 PetscFunctionBegin; 11 for (i=0;i<n;i++) { 12 z[i]=0.0; 13 } 14 e = PetscAbs(r[0]); 15 if (e == 0.0) { 16 *svmin = 0.0; 17 z[0] = 1.0; 18 } else { 19 /* Solve R'*y = e */ 20 for (i=0;i<n;i++) { 21 /* Scale y. The scaling factor (0.01) reduces the number of scalings */ 22 if (z[i] >= 0.0) e =-PetscAbs(e); 23 else e = PetscAbs(e); 24 25 if (PetscAbs(e - z[i]) > PetscAbs(r[i + ldr*i])) { 26 temp = PetscMin(0.01,PetscAbs(r[i + ldr*i]))/PetscAbs(e-z[i]); 27 PetscStackCallBLAS("BLASscal",BLASscal_(&blasn, &temp, z, &blas1)); 28 e = temp*e; 29 } 30 31 /* Determine the two possible choices of y[i] */ 32 if (r[i + ldr*i] == 0.0) { 33 w = wm = 1.0; 34 } else { 35 w = (e - z[i]) / r[i + ldr*i]; 36 wm = - (e + z[i]) / r[i + ldr*i]; 37 } 38 39 /* Chose y[i] based on the predicted value of y[j] for j>i */ 40 s = PetscAbs(e - z[i]); 41 sm = PetscAbs(e + z[i]); 42 for (j=i+1;j<n;j++) { 43 sm += PetscAbs(z[j] + wm * r[i + ldr*j]); 44 } 45 if (i < n-1) { 46 blasnmi = n-i-1; 47 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasnmi, &w, &r[i + ldr*(i+1)], &blasldr, &z[i+1], &blas1)); 48 s += BLASasum_(&blasnmi, &z[i+1], &blas1); 49 } 50 if (s < sm) { 51 temp = wm - w; 52 w = wm; 53 if (i < n-1) { 54 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasnmi, &temp, &r[i + ldr*(i+1)], &blasldr, &z[i+1], &blas1)); 55 } 56 } 57 z[i] = w; 58 } 59 60 ynorm = BLASnrm2_(&blasn, z, &blas1); 61 62 /* Solve R*z = y */ 63 for (j=n-1; j>=0; j--) { 64 /* Scale z */ 65 if (PetscAbs(z[j]) > PetscAbs(r[j + ldr*j])) { 66 temp = PetscMin(0.01, PetscAbs(r[j + ldr*j] / z[j])); 67 PetscStackCallBLAS("BLASscal",BLASscal_(&blasn, &temp, z, &blas1)); 68 ynorm *=temp; 69 } 70 if (r[j + ldr*j] == 0) { 71 z[j] = 1.0; 72 } else { 73 z[j] = z[j] / r[j + ldr*j]; 74 } 75 temp = -z[j]; 76 blasj=j; 77 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasj,&temp,&r[0+ldr*j],&blas1,z,&blas1)); 78 } 79 80 /* Compute svmin and normalize z */ 81 znorm = 1.0 / BLASnrm2_(&blasn, z, &blas1); 82 *svmin = ynorm*znorm; 83 PetscStackCallBLAS("BLASscal",BLASscal_(&blasn, &znorm, z, &blas1)); 84 } 85 PetscFunctionReturn(0); 86 } 87 88 /* 89 c *********** 90 c 91 c Subroutine dgqt 92 c 93 c Given an n by n symmetric matrix A, an n-vector b, and a 94 c positive number delta, this subroutine determines a vector 95 c x which approximately minimizes the quadratic function 96 c 97 c f(x) = (1/2)*x'*A*x + b'*x 98 c 99 c subject to the Euclidean norm constraint 100 c 101 c norm(x) <= delta. 102 c 103 c This subroutine computes an approximation x and a Lagrange 104 c multiplier par such that either par is zero and 105 c 106 c norm(x) <= (1+rtol)*delta, 107 c 108 c or par is positive and 109 c 110 c abs(norm(x) - delta) <= rtol*delta. 111 c 112 c If xsol is the solution to the problem, the approximation x 113 c satisfies 114 c 115 c f(x) <= ((1 - rtol)**2)*f(xsol) 116 c 117 c The subroutine statement is 118 c 119 c subroutine dgqt(n,a,lda,b,delta,rtol,atol,itmax, 120 c par,f,x,info,z,wa1,wa2) 121 c 122 c where 123 c 124 c n is an integer variable. 125 c On entry n is the order of A. 126 c On exit n is unchanged. 127 c 128 c a is a double precision array of dimension (lda,n). 129 c On entry the full upper triangle of a must contain the 130 c full upper triangle of the symmetric matrix A. 131 c On exit the array contains the matrix A. 132 c 133 c lda is an integer variable. 134 c On entry lda is the leading dimension of the array a. 135 c On exit lda is unchanged. 136 c 137 c b is an double precision array of dimension n. 138 c On entry b specifies the linear term in the quadratic. 139 c On exit b is unchanged. 140 c 141 c delta is a double precision variable. 142 c On entry delta is a bound on the Euclidean norm of x. 143 c On exit delta is unchanged. 144 c 145 c rtol is a double precision variable. 146 c On entry rtol is the relative accuracy desired in the 147 c solution. Convergence occurs if 148 c 149 c f(x) <= ((1 - rtol)**2)*f(xsol) 150 c 151 c On exit rtol is unchanged. 152 c 153 c atol is a double precision variable. 154 c On entry atol is the absolute accuracy desired in the 155 c solution. Convergence occurs when 156 c 157 c norm(x) <= (1 + rtol)*delta 158 c 159 c max(-f(x),-f(xsol)) <= atol 160 c 161 c On exit atol is unchanged. 162 c 163 c itmax is an integer variable. 164 c On entry itmax specifies the maximum number of iterations. 165 c On exit itmax is unchanged. 166 c 167 c par is a double precision variable. 168 c On entry par is an initial estimate of the Lagrange 169 c multiplier for the constraint norm(x) <= delta. 170 c On exit par contains the final estimate of the multiplier. 171 c 172 c f is a double precision variable. 173 c On entry f need not be specified. 174 c On exit f is set to f(x) at the output x. 175 c 176 c x is a double precision array of dimension n. 177 c On entry x need not be specified. 178 c On exit x is set to the final estimate of the solution. 179 c 180 c info is an integer variable. 181 c On entry info need not be specified. 182 c On exit info is set as follows: 183 c 184 c info = 1 The function value f(x) has the relative 185 c accuracy specified by rtol. 186 c 187 c info = 2 The function value f(x) has the absolute 188 c accuracy specified by atol. 189 c 190 c info = 3 Rounding errors prevent further progress. 191 c On exit x is the best available approximation. 192 c 193 c info = 4 Failure to converge after itmax iterations. 194 c On exit x is the best available approximation. 195 c 196 c z is a double precision work array of dimension n. 197 c 198 c wa1 is a double precision work array of dimension n. 199 c 200 c wa2 is a double precision work array of dimension n. 201 c 202 c Subprograms called 203 c 204 c MINPACK-2 ...... destsv 205 c 206 c LAPACK ......... dpotrf 207 c 208 c Level 1 BLAS ... daxpy, dcopy, ddot, dnrm2, dscal 209 c 210 c Level 2 BLAS ... dtrmv, dtrsv 211 c 212 c MINPACK-2 Project. October 1993. 213 c Argonne National Laboratory and University of Minnesota. 214 c Brett M. Averick, Richard Carter, and Jorge J. More' 215 c 216 c *********** 217 */ 218 PetscErrorCode gqt(PetscInt n, PetscReal *a, PetscInt lda, PetscReal *b, 219 PetscReal delta, PetscReal rtol, PetscReal atol, 220 PetscInt itmax, PetscReal *retpar, PetscReal *retf, 221 PetscReal *x, PetscInt *retinfo, PetscInt *retits, 222 PetscReal *z, PetscReal *wa1, PetscReal *wa2) 223 { 224 PetscErrorCode ierr; 225 PetscReal f=0.0,p001=0.001,p5=0.5,minusone=-1,delta2=delta*delta; 226 PetscInt iter, j, rednc,info; 227 PetscBLASInt indef; 228 PetscBLASInt blas1=1, blasn=n, iblas, blaslda = lda,blasldap1=lda+1,blasinfo; 229 PetscReal alpha, anorm, bnorm, parc, parf, parl, pars, par=*retpar,paru, prod, rxnorm, rznorm=0.0, temp, xnorm; 230 231 PetscFunctionBegin; 232 parf = 0.0; 233 xnorm = 0.0; 234 rxnorm = 0.0; 235 rednc = 0; 236 for (j=0; j<n; j++) { 237 x[j] = 0.0; 238 z[j] = 0.0; 239 } 240 241 /* Copy the diagonal and save A in its lower triangle */ 242 PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasn,a,&blasldap1, wa1, &blas1)); 243 for (j=0;j<n-1;j++) { 244 iblas = n - j - 1; 245 PetscStackCallBLAS("BLAScopy",BLAScopy_(&iblas,&a[j + lda*(j+1)], &blaslda, &a[j+1 + lda*j], &blas1)); 246 } 247 248 /* Calculate the l1-norm of A, the Gershgorin row sums, and the 249 l2-norm of b */ 250 anorm = 0.0; 251 for (j=0;j<n;j++) { 252 wa2[j] = BLASasum_(&blasn, &a[0 + lda*j], &blas1); 253 CHKMEMQ; 254 anorm = PetscMax(anorm,wa2[j]); 255 } 256 for (j=0;j<n;j++) { 257 wa2[j] = wa2[j] - PetscAbs(wa1[j]); 258 } 259 bnorm = BLASnrm2_(&blasn,b,&blas1); 260 CHKMEMQ; 261 /* Calculate a lower bound, pars, for the domain of the problem. 262 Also calculate an upper bound, paru, and a lower bound, parl, 263 for the Lagrange multiplier. */ 264 pars = parl = paru = -anorm; 265 for (j=0;j<n;j++) { 266 pars = PetscMax(pars, -wa1[j]); 267 parl = PetscMax(parl, wa1[j] + wa2[j]); 268 paru = PetscMax(paru, -wa1[j] + wa2[j]); 269 } 270 parl = PetscMax(bnorm/delta - parl,pars); 271 parl = PetscMax(0.0,parl); 272 paru = PetscMax(0.0, bnorm/delta + paru); 273 274 /* If the input par lies outside of the interval (parl, paru), 275 set par to the closer endpoint. */ 276 277 par = PetscMax(par,parl); 278 par = PetscMin(par,paru); 279 280 /* Special case: parl == paru */ 281 paru = PetscMax(paru, (1.0 + rtol)*parl); 282 283 /* Beginning of an iteration */ 284 285 info = 0; 286 for (iter=1;iter<=itmax;iter++) { 287 /* Safeguard par */ 288 if (par <= pars && paru > 0) { 289 par = PetscMax(p001, PetscSqrtScalar(parl/paru)) * paru; 290 } 291 292 /* Copy the lower triangle of A into its upper triangle and 293 compute A + par*I */ 294 295 for (j=0;j<n-1;j++) { 296 iblas = n - j - 1; 297 PetscStackCallBLAS("BLAScopy",BLAScopy_(&iblas,&a[j+1 + j*lda], &blas1,&a[j + (j+1)*lda], &blaslda)); 298 } 299 for (j=0;j<n;j++) { 300 a[j + j*lda] = wa1[j] + par; 301 } 302 303 /* Attempt the Cholesky factorization of A without referencing 304 the lower triangular part. */ 305 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("U",&blasn,a,&blaslda,&indef)); 306 307 /* Case 1: A + par*I is pos. def. */ 308 if (indef == 0) { 309 310 /* Compute an approximate solution x and save the 311 last value of par with A + par*I pos. def. */ 312 313 parf = par; 314 PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasn, b, &blas1, wa2, &blas1)); 315 PetscStackCallBLAS("LAPACKtrtrs",LAPACKtrtrs_("U","T","N",&blasn,&blas1,a,&blaslda,wa2,&blasn,&blasinfo)); 316 rxnorm = BLASnrm2_(&blasn, wa2, &blas1); 317 PetscStackCallBLAS("LAPACKtrtrs",LAPACKtrtrs_("U","N","N",&blasn,&blas1,a,&blaslda,wa2,&blasn,&blasinfo)); 318 PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasn, wa2, &blas1, x, &blas1)); 319 PetscStackCallBLAS("BLASscal",BLASscal_(&blasn, &minusone, x, &blas1)); 320 xnorm = BLASnrm2_(&blasn, x, &blas1); 321 CHKMEMQ; 322 323 /* Test for convergence */ 324 if (PetscAbs(xnorm - delta) <= rtol*delta || 325 (par == 0 && xnorm <= (1.0+rtol)*delta)) { 326 info = 1; 327 } 328 329 /* Compute a direction of negative curvature and use this 330 information to improve pars. */ 331 332 iblas=blasn*blasn; 333 334 ierr = estsv(n,a,lda,&rznorm,z);CHKERRQ(ierr); 335 CHKMEMQ; 336 pars = PetscMax(pars, par-rznorm*rznorm); 337 338 /* Compute a negative curvature solution of the form 339 x + alpha*z, where norm(x+alpha*z)==delta */ 340 341 rednc = 0; 342 if (xnorm < delta) { 343 /* Compute alpha */ 344 prod = BLASdot_(&blasn, z, &blas1, x, &blas1) / delta; 345 temp = (delta - xnorm)*((delta + xnorm)/delta); 346 alpha = temp/(PetscAbs(prod) + PetscSqrtScalar(prod*prod + temp/delta)); 347 if (prod >= 0) alpha = PetscAbs(alpha); 348 else alpha =-PetscAbs(alpha); 349 350 /* Test to decide if the negative curvature step 351 produces a larger reduction than with z=0 */ 352 rznorm = PetscAbs(alpha) * rznorm; 353 if ((rznorm*rznorm + par*xnorm*xnorm)/(delta2) <= par) { 354 rednc = 1; 355 } 356 /* Test for convergence */ 357 if (p5 * rznorm*rznorm / delta2 <= rtol*(1.0-p5*rtol)*(par + rxnorm*rxnorm/delta2)) { 358 info = 1; 359 } else if (info == 0 && (p5*(par + rxnorm*rxnorm/delta2) <= atol/delta2)) { 360 info = 2; 361 } 362 } 363 364 /* Compute the Newton correction parc to par. */ 365 if (xnorm == 0) { 366 parc = -par; 367 } else { 368 PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasn, x, &blas1, wa2, &blas1)); 369 temp = 1.0/xnorm; 370 PetscStackCallBLAS("BLASscal",BLASscal_(&blasn, &temp, wa2, &blas1)); 371 PetscStackCallBLAS("LAPACKtrtrs",LAPACKtrtrs_("U","T","N",&blasn, &blas1, a, &blaslda, wa2, &blasn, &blasinfo)); 372 temp = BLASnrm2_(&blasn, wa2, &blas1); 373 parc = (xnorm - delta)/(delta*temp*temp); 374 } 375 376 /* update parl or paru */ 377 if (xnorm > delta) { 378 parl = PetscMax(parl, par); 379 } else if (xnorm < delta) { 380 paru = PetscMin(paru, par); 381 } 382 } else { 383 /* Case 2: A + par*I is not pos. def. */ 384 385 /* Use the rank information from the Cholesky 386 decomposition to update par. */ 387 388 if (indef > 1) { 389 /* Restore column indef to A + par*I. */ 390 iblas = indef - 1; 391 PetscStackCallBLAS("BLAScopy",BLAScopy_(&iblas,&a[indef-1 + 0*lda],&blaslda,&a[0 + (indef-1)*lda],&blas1)); 392 a[indef-1 + (indef-1)*lda] = wa1[indef-1] + par; 393 394 /* compute parc. */ 395 PetscStackCallBLAS("BLAScopy",BLAScopy_(&iblas,&a[0 + (indef-1)*lda], &blas1, wa2, &blas1)); 396 PetscStackCallBLAS("LAPACKtrtrs",LAPACKtrtrs_("U","T","N",&iblas,&blas1,a,&blaslda,wa2,&blasn,&blasinfo)); 397 PetscStackCallBLAS("BLAScopy",BLAScopy_(&iblas,wa2,&blas1,&a[0 + (indef-1)*lda],&blas1)); 398 temp = BLASnrm2_(&iblas,&a[0 + (indef-1)*lda],&blas1); 399 CHKMEMQ; 400 a[indef-1 + (indef-1)*lda] -= temp*temp; 401 PetscStackCallBLAS("LAPACKtrtrs",LAPACKtrtrs_("U","N","N",&iblas,&blas1,a,&blaslda,wa2,&blasn,&blasinfo)); 402 } 403 404 wa2[indef-1] = -1.0; 405 iblas = indef; 406 temp = BLASnrm2_(&iblas,wa2,&blas1); 407 parc = - a[indef-1 + (indef-1)*lda]/(temp*temp); 408 pars = PetscMax(pars,par+parc); 409 410 /* If necessary, increase paru slightly. 411 This is needed because in some exceptional situations 412 paru is the optimal value of par. */ 413 414 paru = PetscMax(paru, (1.0+rtol)*pars); 415 } 416 417 /* Use pars to update parl */ 418 parl = PetscMax(parl,pars); 419 420 /* Test for converged. */ 421 if (info == 0) { 422 if (iter == itmax) info=4; 423 if (paru <= (1.0+p5*rtol)*pars) info=3; 424 if (paru == 0.0) info = 2; 425 } 426 427 /* If exiting, store the best approximation and restore 428 the upper triangle of A. */ 429 430 if (info != 0) { 431 /* Compute the best current estimates for x and f. */ 432 par = parf; 433 f = -p5 * (rxnorm*rxnorm + par*xnorm*xnorm); 434 if (rednc) { 435 f = -p5 * (rxnorm*rxnorm + par*delta*delta - rznorm*rznorm); 436 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&blasn, &alpha, z, &blas1, x, &blas1)); 437 } 438 /* Restore the upper triangle of A */ 439 for (j = 0; j<n; j++) { 440 iblas = n - j - 1; 441 PetscStackCallBLAS("BLAScopy",BLAScopy_(&iblas,&a[j+1 + j*lda],&blas1, &a[j + (j+1)*lda],&blaslda)); 442 } 443 iblas = lda+1; 444 PetscStackCallBLAS("BLAScopy",BLAScopy_(&blasn,wa1,&blas1,a,&iblas)); 445 break; 446 } 447 par = PetscMax(parl,par+parc); 448 } 449 *retpar = par; 450 *retf = f; 451 *retinfo = info; 452 *retits = iter; 453 CHKMEMQ; 454 PetscFunctionReturn(0); 455 } 456