Millepede-II  V04-01-00
 All Classes Files Functions Variables Enumerator Pages
minres.f90
Go to the documentation of this file.
1 
2 ! Code converted using TO_F90 by Alan Miller
3 ! Date: 2012-03-16 Time: 11:06:43
4 
7 
288 
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 )
292  USE mpdef
293 
294  IMPLICIT NONE
295  EXTERNAL aprod, msolve
296  INTEGER(mpi) :: n, nout, itnlim, istop, itn
297  LOGICAL :: checka, precon
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)
300 
301  EXTERNAL ddot , dnrm2
302  REAL(mpd) :: ddot , dnrm2
303 
304  ! Local variables
305 
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
310  INTEGER(mpi) :: i
311  LOGICAL :: debug, prnt
312 
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 )
316 
317  CHARACTER (LEN=16) :: enter, exit
318  CHARACTER (LEN=52) :: msg(-1:8)
319 
320  DATA enter /' Enter MINRES. '/, exit /' Exit MINRES. '/
321 
322  DATA msg &
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' /
332  ! ------------------------------------------------------------------
333 
334  debug = .false.
335 
336  ! ------------------------------------------------------------------
337  ! Compute eps, the machine precision. The call to daxpy is
338  ! intended to fool compilers that use extra-length registers.
339  ! 31 May 1999: Hardwire eps so the debugger can step thru easily.
340  ! ------------------------------------------------------------------
341  eps = epsilon(eps) ! Set eps = zero here if you want it computed.
342 
343  eps = 0.0d0 !!!!!!!!!!!!!!!
344  gmin = 0.0d0
345  gmax = 0.0d0
346 
347  IF (eps <= zero) THEN
348  eps = two**(-12)
349  DO
350  eps = eps / two
351  x(1) = eps
352  y(1) = one
353  CALL daxpy( 1, one, x, 1, y, 1 )
354  IF (y(1) <= one) exit
355  END DO
356  eps = eps * two
357  END IF
358  ! WRITE(*,*) 'computed epsilon is ',eps !!!!!!!!!!!!!!!
359 
360  ! ------------------------------------------------------------------
361  ! Print heading and initialize.
362  ! ------------------------------------------------------------------
363  IF (nout > 0) THEN
364  WRITE(nout, 1000) enter, n, checka, precon, itnlim, rtol, shift
365  END IF
366  istop = 0
367  itn = 0
368  anorm = zero
369  acond = zero
370  rnorm = zero
371  ynorm = zero
372  CALL dload2( n, zero, x )
373 
374  ! ------------------------------------------------------------------
375  ! Set up y and v for the first Lanczos vector v1.
376  ! y = beta1 P' v1, where P = C**(-1).
377  ! v is really P' v1.
378  ! ------------------------------------------------------------------
379  CALL dcopy( n, b, 1, y , 1 ) ! y = b
380  CALL dcopy( n, b, 1, r1, 1 ) ! r1 = b
381  IF ( precon ) CALL msolve( n, b, y )
382  beta1 = ddot( n, b, 1, y, 1 )
383 
384  IF (beta1 < zero) THEN ! m must be indefinite.
385  istop = 8
386  go to 900
387  END IF
388 
389  IF (beta1 == zero) THEN ! b = 0 exactly. Stop with x = 0.
390  istop = 0
391  go to 900
392  END IF
393 
394  beta1 = sqrt( beta1 ) ! Normalize y to get v1 later.
395 
396  ! ------------------------------------------------------------------
397  ! See if Msolve is symmetric.
398  ! ------------------------------------------------------------------
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 )
403  z = abs( s - t )
404  epsa = (s + eps) * eps**0.33333_mpd
405  IF (z > epsa) THEN
406  istop = 7
407  go to 900
408  END IF
409  END IF
410 
411  ! ------------------------------------------------------------------
412  ! See if Aprod is symmetric.
413  ! ------------------------------------------------------------------
414  IF (checka) THEN
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 )
419  z = abs( s - t )
420  epsa = (s + eps) * eps**0.33333_mpd
421  IF (z > epsa) THEN
422  istop = 6
423  go to 900
424  END IF
425  END IF
426 
427  ! ------------------------------------------------------------------
428  ! Initialize other quantities.
429  ! ------------------------------------------------------------------
430  oldb = zero
431  beta = beta1
432  dbar = zero
433  epsln = zero
434  qrnorm = beta1
435  phibar = beta1
436  rhs1 = beta1
437  rhs2 = zero
438  tnorm2 = zero
439  ynorm2 = zero
440  cs = - one
441  sn = zero
442  CALL dload2( n, zero, w ) ! w = 0
443  CALL dload2( n, zero, w2 ) ! w2 = 0
444  CALL dcopy( n, r1, 1, r2, 1 ) ! r2 = r1
445 
446  IF (debug) THEN
447  WRITE(*,*) ' '
448  WRITE(*,*) 'b ', b
449  WRITE(*,*) 'beta ', beta
450  WRITE(*,*) ' '
451  END IF
452 
453  ! ------------------------------------------------------------------
454  ! Main iteration loop.
455  ! ------------------------------------------------------------------
456  DO
457  itn = itn + 1 ! k = itn = 1 first time through
458  IF (istop /= 0) exit
459 
460  !-----------------------------------------------------------------
461  ! Obtain quantities for the next Lanczos vector vk+1, k = 1, 2,...
462  ! The general iteration is similar to the case k = 1 with v0 = 0:
463  !
464  ! p1 = Operator * v1 - beta1 * v0,
465  ! alpha1 = v1'p1,
466  ! q2 = p2 - alpha1 * v1,
467  ! beta2^2 = q2'q2,
468  ! v2 = (1/beta2) q2.
469  !
470  ! Again, y = betak P vk, where P = C**(-1).
471  ! .... more description needed.
472  !-----------------------------------------------------------------
473  s = one / beta ! Normalize previous vector (in y).
474  CALL dscal2( n, s, y, v ) ! v = vk if P = I
475 
476  CALL aprod( n, v, y )
477  CALL daxpy( n, (- shift), v, 1, y, 1 )
478  IF (itn >= 2) THEN
479  CALL daxpy( n, (- beta/oldb), r1, 1, y, 1 )
480  END IF
481 
482  alfa = ddot( n, v, 1, y, 1 ) ! alphak
483 
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 )
488 
489  oldb = beta ! oldb = betak
490  beta = ddot( n, r2, 1, y, 1 ) ! beta = betak+1^2
491  IF (beta < zero) THEN
492  istop = 6
493  WRITE(*,*) 'unfortunately beta < 0 ',beta,zero
494  exit
495  END IF
496 
497  beta = sqrt( beta ) ! beta = betak+1
498  tnorm2 = tnorm2 + alfa**2 + oldb**2 + beta**2
499 
500  IF (itn == 1) THEN ! Initialize a few things.
501  IF (beta/beta1 <= ten*eps) THEN ! beta2 = 0 or ~ 0.
502  istop = -1 ! Terminate later.
503  END IF
504  !tnorm2 = alfa**2
505  gmax = abs( alfa ) ! alpha1
506  gmin = gmax ! alpha1
507  END IF
508 
509  ! Apply previous rotation Qk-1 to get
510  ! [deltak epslnk+1] = [cs sn][dbark 0 ]
511  ! [gbar k dbar k+1] [sn -cs][alfak betak+1].
512 
513  oldeps = epsln
514  delta = cs * dbar + sn * alfa ! delta1 = 0 deltak
515  gbar = sn * dbar - cs * alfa ! gbar 1 = alfa1 gbar k
516  epsln = sn * beta ! epsln2 = 0 epslnk+1
517  dbar = - cs * beta ! dbar 2 = beta2 dbar k+1
518 
519  ! Compute the next plane rotation Qk
520 
521  agamma = sqrt( gbar**2 + beta**2 ) ! gammak
522  cs = gbar / agamma ! ck
523  sn = beta / agamma ! sk
524  phi = cs * phibar ! phik
525  phibar = sn * phibar ! phibark+1
526 
527  IF (debug) THEN
528  WRITE(*,*) ' '
529  WRITE(*,*) 'v ', v
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
539  WRITE(*,*) ' '
540  END IF
541 
542  ! Update x.
543 
544  denom = one/agamma
545 
546  DO i = 1, n
547  w1(i) = w2(i)
548  w2(i) = w(i)
549  w(i) = ( v(i) - oldeps*w1(i) - delta*w2(i) ) * denom
550  x(i) = x(i) + phi * w(i)
551  END DO
552 
553  ! Go round again.
554 
555  gmax = max( gmax, agamma)
556  gmin = min( gmin, agamma)
557  z = rhs1 / agamma
558  ynorm2 = z**2 + ynorm2
559  rhs1 = rhs2 - delta * z
560  rhs2 = - epsln * z
561 
562  ! Estimate various norms and test for convergence.
563 
564  anorm = sqrt( tnorm2 )
565  ynorm = sqrt( ynorm2 )
566  epsa = anorm * eps
567  epsx = anorm * ynorm * eps
568  epsr = anorm * ynorm * rtol
569  diag = gbar
570  IF (diag == zero) diag = epsa
571 
572  qrnorm = phibar
573  rnorm = qrnorm
574  relarnorml = sqrt(gbar**2+dbar**2)/anorm
575  print *, ' iter ', itn, anorm, ynorm, qrnorm, relarnorml, rtol
576 
577  ! Estimate cond(A).
578  ! In this version we look at the diagonals of R in the
579  ! factorization of the lower Hessenberg matrix, Q * H = R,
580  ! where H is the tridiagonal matrix from Lanczos with one
581  ! extra row, beta(k+1) e_k^T.
582 
583  acond = gmax / gmin
584 
585  ! See if any of the stopping criteria are satisfied.
586  ! In rare cases, istop is already -1 from above (Abar = const*I).
587 
588  IF (istop == 0) THEN
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
594  END IF
595 
596 
597  ! See if it is time to print something.
598 
599  IF (nout > 0) THEN
600  prnt = .false.
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.
609 
610  IF ( prnt ) THEN
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)
614  END IF
615  END IF
616 
617  END DO
618  ! ------------------------------------------------------------------
619  ! End of main iteration loop.
620  ! ------------------------------------------------------------------
621 
622  ! Display final status.
623 
624 900 IF (nout > 0) THEN
625  WRITE(nout, 2000) exit, istop, itn, exit, anorm, acond, &
626  exit, rnorm, ynorm
627  WRITE(nout, 3000) exit, msg(istop)
628  END IF
629 
630  return
631 
632 
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)
639 1500 format(1x)
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 )
644 
645 END SUBROUTINE minres