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