1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424
|
program zndrv1
c
c Example program to illustrate the idea of reverse communication
c for a standard complex nonsymmetric eigenvalue problem.
c
c We implement example one of ex-complex.doc in DOCUMENTS directory
c
c\Example-1
c ... Suppose we want to solve A*x = lambda*x in regular mode,
c where A is obtained from the standard central difference
c discretization of the convection-diffusion operator
c (Laplacian u) + rho*(du / dx)
c on the unit squre [0,1]x[0,1] with zero Dirichlet boundary
c condition.
c
c ... OP = A and B = I.
c
c ... Assume "call av (nx,x,y)" computes y = A*x
c
c ... Use mode 1 of ZNAUPD .
c
c\BeginLib
c
c\Routines called
c znaupd ARPACK reverse communication interface routine.
c zneupd ARPACK routine that returns Ritz values and (optionally)
c Ritz vectors.
c dlapy2 LAPACK routine to compute sqrt(x**2+y**2) carefully.
c dznrm2 Level 1 BLAS that computes the norm of a complex vector.
c zaxpy Level 1 BLAS that computes y <- alpha*x+y.
c av Matrix vector multiplication routine that computes A*x.
c tv Matrix vector multiplication routine that computes T*x,
c where T is a tridiagonal matrix. It is used in routine
c av.
c
c\Author
c Richard Lehoucq
c Danny Sorensen
c Chao Yang
c Dept. of Computational &
c Applied Mathematics
c Rice University
c Houston, Texas
c
c\SCCS Information: @(#)
c FILE: ndrv1.F SID: 2.4 DATE OF SID: 10/17/00 RELEASE: 2
c
c\Remarks
c 1. None
c
c\EndLib
c---------------------------------------------------------------------------
c
c %-----------------------------%
c | Define maximum dimensions |
c | for all arrays. |
c | MAXN: Maximum dimension |
c | of the A allowed. |
c | MAXNEV: Maximum NEV allowed |
c | MAXNCV: Maximum NCV allowed |
c %-----------------------------%
c
integer maxn, maxnev, maxncv, ldv
parameter (maxn=256, maxnev=12, maxncv=30, ldv=maxn)
c
c %--------------%
c | Local Arrays |
c %--------------%
c
integer iparam(11), ipntr(14)
logical select(maxncv)
Complex*16
& ax(maxn), d(maxncv),
& v(ldv,maxncv), workd(3*maxn),
& workev(3*maxncv), resid(maxn),
& workl(3*maxncv*maxncv+5*maxncv)
Double precision
& rwork(maxncv), rd(maxncv,3)
c
c %---------------%
c | Local Scalars |
c %---------------%
c
character bmat*1, which*2
integer ido, n, nx, nev, ncv, lworkl, info, j,
& ierr, nconv, maxitr, ishfts, mode
Complex*16
& sigma
Double precision
& tol
logical rvec
c
c %-----------------------------%
c | BLAS & LAPACK routines used |
c %-----------------------------%
c
Double precision
& dznrm2 , dlapy2
external dznrm2 , zaxpy , dlapy2
c
c %-----------------------%
c | Executable Statements |
c %-----------------------%
c
c %--------------------------------------------------%
c | The number NX is the number of interior points |
c | in the discretization of the 2-dimensional |
c | convection-diffusion operator on the unit |
c | square with zero Dirichlet boundary condition. |
c | The number N(=NX*NX) is the dimension of the |
c | matrix. A standard eigenvalue problem is |
c | solved (BMAT = 'I'). NEV is the number of |
c | eigenvalues to be approximated. The user can |
c | modify NX, NEV, NCV, WHICH to solve problems of |
c | different sizes, and to get different parts of |
c | the spectrum. However, The following |
c | conditions must be satisfied: |
c | N <= MAXN |
c | NEV <= MAXNEV |
c | NEV + 2 <= NCV <= MAXNCV |
c %--------------------------------------------------%
c
nx = 10
n = nx*nx
nev = 4
ncv = 20
if ( n .gt. maxn ) then
print *, ' ERROR with _NDRV1: N is greater than MAXN '
go to 9000
else if ( nev .gt. maxnev ) then
print *, ' ERROR with _NDRV1: NEV is greater than MAXNEV '
go to 9000
else if ( ncv .gt. maxncv ) then
print *, ' ERROR with _NDRV1: NCV is greater than MAXNCV '
go to 9000
end if
bmat = 'I'
which = 'LM'
c
c %---------------------------------------------------%
c | The work array WORKL is used in ZNAUPD as |
c | workspace. Its dimension LWORKL is set as |
c | illustrated below. The parameter TOL determines |
c | the stopping criterion. If TOL<=0, machine |
c | precision is used. The variable IDO is used for |
c | reverse communication, and is initially set to 0. |
c | Setting INFO=0 indicates that a random vector is |
c | generated to start the ARNOLDI iteration. |
c %---------------------------------------------------%
c
lworkl = 3*ncv**2+5*ncv
tol = 0.0
ido = 0
info = 0
c
c %---------------------------------------------------%
c | This program uses exact shift with respect to |
c | the current Hessenberg matrix (IPARAM(1) = 1). |
c | IPARAM(3) specifies the maximum number of Arnoldi |
c | iterations allowed. Mode 1 of ZNAUPD is used |
c | (IPARAM(7) = 1). All these options can be changed |
c | by the user. For details see the documentation in |
c | ZNAUPD . |
c %---------------------------------------------------%
c
ishfts = 1
maxitr = 300
mode = 1
c
iparam(1) = ishfts
iparam(3) = maxitr
iparam(7) = mode
c
c %-------------------------------------------%
c | M A I N L O O P (Reverse communication) |
c %-------------------------------------------%
c
10 continue
c
c %---------------------------------------------%
c | Repeatedly call the routine ZNAUPD and take |
c | actions indicated by parameter IDO until |
c | either convergence is indicated or maxitr |
c | has been exceeded. |
c %---------------------------------------------%
c
call znaupd ( ido, bmat, n, which, nev, tol, resid, ncv,
& v, ldv, iparam, ipntr, workd, workl, lworkl,
& rwork,info )
c
if (ido .eq. -1 .or. ido .eq. 1) then
c
c %-------------------------------------------%
c | Perform matrix vector multiplication |
c | y <--- OP*x |
c | The user should supply his/her own |
c | matrix vector multiplication routine here |
c | that takes workd(ipntr(1)) as the input |
c | vector, and return the matrix vector |
c | product to workd(ipntr(2)). |
c %-------------------------------------------%
c
call av (nx, workd(ipntr(1)), workd(ipntr(2)))
c
c %-----------------------------------------%
c | L O O P B A C K to call ZNAUPD again. |
c %-----------------------------------------%
c
go to 10
end if
c
c %----------------------------------------%
c | Either we have convergence or there is |
c | an error. |
c %----------------------------------------%
c
if ( info .lt. 0 ) then
c
c %--------------------------%
c | Error message, check the |
c | documentation in ZNAUPD |
c %--------------------------%
c
print *, ' '
print *, ' Error with _naupd, info = ', info
print *, ' Check the documentation of _naupd'
print *, ' '
c
else
c
c %-------------------------------------------%
c | No fatal errors occurred. |
c | Post-Process using ZNEUPD . |
c | |
c | Computed eigenvalues may be extracted. |
c | |
c | Eigenvectors may also be computed now if |
c | desired. (indicated by rvec = .true.) |
c %-------------------------------------------%
c
rvec = .true.
c
call zneupd (rvec, 'A', select, d, v, ldv, sigma,
& workev, bmat, n, which, nev, tol, resid, ncv,
& v, ldv, iparam, ipntr, workd, workl, lworkl,
& rwork, ierr)
c
c %----------------------------------------------%
c | Eigenvalues are returned in the one |
c | dimensional array D. The corresponding |
c | eigenvectors are returned in the first NCONV |
c | (=IPARAM(5)) columns of the two dimensional |
c | array V if requested. Otherwise, an |
c | orthogonal basis for the invariant subspace |
c | corresponding to the eigenvalues in D is |
c | returned in V. |
c %----------------------------------------------%
c
if ( ierr .ne. 0) then
c
c %------------------------------------%
c | Error condition: |
c | Check the documentation of ZNEUPD . |
c %------------------------------------%
c
print *, ' '
print *, ' Error with _neupd, info = ', ierr
print *, ' Check the documentation of _neupd. '
print *, ' '
c
else
c
nconv = iparam(5)
do 20 j=1, nconv
c
c %---------------------------%
c | Compute the residual norm |
c | |
c | || A*x - lambda*x || |
c | |
c | for the NCONV accurately |
c | computed eigenvalues and |
c | eigenvectors. (iparam(5) |
c | indicates how many are |
c | accurate to the requested |
c | tolerance) |
c %---------------------------%
c
call av(nx, v(1,j), ax)
call zaxpy (n, -d(j), v(1,j), 1, ax, 1)
rd(j,1) = dble (d(j))
rd(j,2) = dimag (d(j))
rd(j,3) = dznrm2 (n, ax, 1)
rd(j,3) = rd(j,3) / dlapy2 (rd(j,1),rd(j,2))
20 continue
c
c %-----------------------------%
c | Display computed residuals. |
c %-----------------------------%
c
call dmout (6, nconv, 3, rd, maxncv, -6,
& 'Ritz values (Real, Imag) and relative residuals')
end if
c
c %-------------------------------------------%
c | Print additional convergence information. |
c %-------------------------------------------%
c
if ( info .eq. 1) then
print *, ' '
print *, ' Maximum number of iterations reached.'
print *, ' '
else if ( info .eq. 3) then
print *, ' '
print *, ' No shifts could be applied during implicit',
& ' Arnoldi update, try increasing NCV.'
print *, ' '
end if
c
print *, ' '
print *, '_NDRV1'
print *, '====== '
print *, ' '
print *, ' Size of the matrix is ', n
print *, ' The number of Ritz values requested is ', nev
print *, ' The number of Arnoldi vectors generated',
& ' (NCV) is ', ncv
print *, ' What portion of the spectrum: ', which
print *, ' The number of converged Ritz values is ',
& nconv
print *, ' The number of Implicit Arnoldi update',
& ' iterations taken is ', iparam(3)
print *, ' The number of OP*x is ', iparam(9)
print *, ' The convergence criterion is ', tol
print *, ' '
c
end if
c
c %---------------------------%
c | Done with program zndrv1 . |
c %---------------------------%
c
9000 continue
c
end
c
c==========================================================================
c
c matrix vector subroutine
c
c The matrix used is the convection-diffusion operator
c discretized using centered difference.
c
subroutine av (nx, v, w)
integer nx, j, lo
Complex*16
& v(nx*nx), w(nx*nx), one, h2
parameter (one = (1.0D+0, 0.0D+0) )
external zaxpy , tv
c
c Computes w <--- OP*v, where OP is the nx*nx by nx*nx block
c tridiagonal matrix
c
c | T -I |
c |-I T -I |
c OP = | -I T |
c | ... -I|
c | -I T|
c
c derived from the standard central difference discretization
c of the convection-diffusion operator (Laplacian u) + rho*(du/dx)
c with zero boundary condition.
c
c The subroutine TV is called to computed y<---T*x.
c
c
h2 = one / dcmplx ((nx+1)*(nx+1))
c
call tv(nx,v(1),w(1))
call zaxpy (nx, -one/h2, v(nx+1), 1, w(1), 1)
c
do 10 j = 2, nx-1
lo = (j-1)*nx
call tv(nx, v(lo+1), w(lo+1))
call zaxpy (nx, -one/h2, v(lo-nx+1), 1, w(lo+1), 1)
call zaxpy (nx, -one/h2, v(lo+nx+1), 1, w(lo+1), 1)
10 continue
c
lo = (nx-1)*nx
call tv(nx, v(lo+1), w(lo+1))
call zaxpy (nx, -one/h2, v(lo-nx+1), 1, w(lo+1), 1)
c
return
end
c=========================================================================
subroutine tv (nx, x, y)
c
integer nx, j
Complex*16
& x(nx), y(nx), h, h2, dd, dl, du
c
Complex*16
& one, rho
parameter (one = (1.0D+0, 0.0D+0) ,
& rho = (1.0D+2, 0.0D+0) )
c
c Compute the matrix vector multiplication y<---T*x
c where T is a nx by nx tridiagonal matrix with DD on the
c diagonal, DL on the subdiagonal, and DU on the superdiagonal
c
h = one / dcmplx (nx+1)
h2 = h*h
dd = (4.0D+0, 0.0D+0) / h2
dl = -one/h2 - (5.0D-1, 0.0D+0) *rho/h
du = -one/h2 + (5.0D-1, 0.0D+0) *rho/h
c
y(1) = dd*x(1) + du*x(2)
do 10 j = 2,nx-1
y(j) = dl*x(j-1) + dd*x(j) + du*x(j+1)
10 continue
y(nx) = dl*x(nx-1) + dd*x(nx)
return
end
|