289 SUBROUTINE minres( n, b, r1, r2, v, w, w1, w2, x, y, &
290 aprod, msolve, checka,
precon, shift, nout , itnlim, rtol, &
291 istop, itn, anorm, acond, rnorm, ynorm )
295 EXTERNAL aprod, msolve
296 INTEGER(mpi) :: n, nout, itnlim, istop, itn
298 REAL(mpd) :: shift, rtol, anorm, acond, rnorm, ynorm, &
299 b(n), r1(n), r2(n), v(n), w(n), w1(n), w2(n), x(n), y(n)
306 REAL(mpd) :: alfa , beta , beta1 , cs , &
307 dbar , delta , denom , diag , eps , epsa , epsln , epsr , epsx , &
308 agamma, gbar , gmax , gmin , oldb , oldeps, qrnorm, phi , phibar, &
309 rhs1 , rhs2 , s , sn , t , tnorm2, ynorm2, z, relarnorml
311 LOGICAL :: debug, prnt
313 REAL(mpd) :: zero, one, two, ten
314 parameter( zero = 0.0_mpd, one = 1.0_mpd, &
315 two = 2.0_mpd, ten = 10.0_mpd )
317 CHARACTER (LEN=16) :: enter, exit
318 CHARACTER (LEN=52) :: msg(-1:8)
320 DATA enter /
' Enter MINRES. '/, exit /
' Exit MINRES. '/
323 /
'beta2 = 0. If M = I, b and x are eigenvectors of A', &
324 'beta1 = 0. The exact solution is x = 0', &
325 'Requested accuracy achieved, as determined by rtol', &
326 'Reasonable accuracy achieved, given eps', &
327 'x has converged to an eigenvector',
'Acond has exceeded 0.1/eps', &
328 'The iteration limit was reached', &
329 'Aprod does not define a symmetric matrix', &
330 'Msolve does not define a symmetric matrix', &
331 'Msolve does not define a pos-def preconditioner' /
347 IF (eps <= zero)
THEN
353 CALL
daxpy( 1, one, x, 1, y, 1 )
354 IF (y(1) <= one) exit
364 WRITE(nout, 1000) enter, n, checka,
precon, itnlim, rtol, shift
379 CALL
dcopy( n, b, 1, y , 1 )
380 CALL
dcopy( n, b, 1, r1, 1 )
381 IF (
precon ) CALL msolve( n, b, y )
382 beta1 =
ddot( n, b, 1, y, 1 )
384 IF (beta1 < zero)
THEN
389 IF (beta1 == zero)
THEN
394 beta1 = sqrt( beta1 )
399 IF (checka .AND.
precon)
THEN
400 CALL msolve( n, y, r2 )
401 s =
ddot( n, y, 1, y, 1 )
402 t =
ddot( n,r1, 1,r2, 1 )
404 epsa = (s + eps) * eps**0.33333_mpd
415 CALL aprod( n, y, w )
416 CALL aprod( n, w, r2 )
417 s =
ddot( n, w, 1, w, 1 )
418 t =
ddot( n, y, 1,r2, 1 )
420 epsa = (s + eps) * eps**0.33333_mpd
443 CALL
dload2( n, zero, w2 )
444 CALL
dcopy( n, r1, 1, r2, 1 )
449 WRITE(*,*)
'beta ', beta
476 CALL aprod( n, v, y )
477 CALL
daxpy( n, (- shift), v, 1, y, 1 )
479 CALL
daxpy( n, (- beta/oldb), r1, 1, y, 1 )
482 alfa =
ddot( n, v, 1, y, 1 )
484 CALL
daxpy( n, (- alfa/beta), r2, 1, y, 1 )
485 CALL
dcopy( n, r2, 1, r1, 1 )
486 CALL
dcopy( n, y, 1, r2, 1 )
487 IF (
precon ) CALL msolve( n, r2, y )
490 beta =
ddot( n, r2, 1, y, 1 )
491 IF (beta < zero)
THEN
493 WRITE(*,*)
'unfortunately beta < 0 ',beta,zero
498 tnorm2 = tnorm2 + alfa**2 + oldb**2 + beta**2
501 IF (beta/beta1 <= ten*eps)
THEN
514 delta = cs * dbar + sn * alfa
515 gbar = sn * dbar - cs * alfa
521 agamma = sqrt( gbar**2 + beta**2 )
530 WRITE(*,*)
'alfa ', alfa
531 WRITE(*,*)
'beta ', beta
532 WRITE(*,*)
'gamma', agamma
533 WRITE(*,*)
'delta', delta
534 WRITE(*,*)
'gbar ', gbar
535 WRITE(*,*)
'epsln', epsln
536 WRITE(*,*)
'dbar ', dbar
537 WRITE(*,*)
'phi ', phi
538 WRITE(*,*)
'phiba', phibar
549 w(i) = ( v(i) - oldeps*w1(i) - delta*w2(i) ) * denom
550 x(i) = x(i) + phi * w(i)
555 gmax = max( gmax, agamma)
556 gmin = min( gmin, agamma)
558 ynorm2 = z**2 + ynorm2
559 rhs1 = rhs2 - delta * z
564 anorm = sqrt( tnorm2 )
565 ynorm = sqrt( ynorm2 )
567 epsx = anorm * ynorm * eps
568 epsr = anorm * ynorm * rtol
570 IF (diag == zero) diag = epsa
574 relarnorml = sqrt(gbar**2+dbar**2)/anorm
575 print *,
' iter ', itn, anorm, ynorm, qrnorm, relarnorml, rtol
589 IF (itn >= itnlim ) istop = 5
590 IF (acond >= 0.1_mpd/eps) istop = 4
591 IF (epsx >= beta1 ) istop = 3
592 IF (qrnorm <= epsx ) istop = 2
593 IF (qrnorm <= epsr ) istop = 1
601 IF (n <= 40 ) prnt = .true.
602 IF (itn <= 10 ) prnt = .true.
603 IF (itn >= itnlim - 10) prnt = .true.
604 IF (mod(itn,10) == 0) prnt = .true.
605 IF (qrnorm <= ten * epsx) prnt = .true.
606 IF (qrnorm <= ten * epsr) prnt = .true.
607 IF (acond >= 1.0e-2_mpd/eps ) prnt = .true.
608 IF (istop /= 0 ) prnt = .true.
611 IF ( itn == 1)
WRITE(nout, 1200)
612 WRITE(nout, 1300) itn, x(1), qrnorm, anorm, acond
613 IF (mod(itn,10) == 0)
WRITE(nout, 1500)
624 900
IF (nout > 0)
THEN
625 WRITE(nout, 2000) exit, istop, itn, exit, anorm, acond, &
627 WRITE(nout, 3000) exit, msg(istop)
633 1000 format(// 1p, a, 5x,
'Solution of symmetric Ax = b' &
634 /
' n =', i7, 5x,
'checkA =', l4, 12x,
'precon =', l4 &
635 /
' itnlim =', i7, 5x,
'rtol =', e11.2, 5x,
'shift =', e23.14)
636 1200 format(// 5x,
'itn', 8x,
'x(1)', 10x, &
637 'norm(r)', 3x,
'norm(A)', 3x,
'cond(A)')
638 1300 format(1p, i8, e19.10, 3e10.2)
640 2000 format(/ 1p, a, 5x,
'istop =', i3, 14x,
'itn =', i8 &
641 / a, 5x,
'Anorm =', e12.4, 5x,
'Acond =', e12.4 &
642 / a, 5x,
'rnorm =', e12.4, 5x,
'ynorm =', e12.4)
643 3000 format( a, 5x, a )