Accumulates the stabilized elementary similarity transformations used in the reduction of a real general matrix to upper Hessenberg form by ELMHES.
This subroutine is a translation of the ALGOL procedure ELMTRANS, NUM. MATH. 16, 181-204(1970) by Peters and Wilkinson. HANDBOOK FOR AUTO. COMP., VOL.II-LINEAR ALGEBRA, 372-395(1971).
This subroutine accumulates the stabilized elementary similarity transformations used in the reduction of a REAL GENERAL matrix to upper Hessenberg form by ELMHES.
NM must be set to the row dimension of the two-dimensional
array parameters, A and Z, as declared in the calling
program dimension statement. NM is an INTEGER variable.
N is the order of the matrix A. N is an INTEGER variable.
N must be less than or equal to NM.
LOW and IGH are two INTEGER variables determined by the
balancing subroutine BALANC. If BALANC has not been
used, set LOW=1 and IGH equal to the order of the matrix, N.
A contains the multipliers which were used in the reduction
by ELMHES in its lower triangle below the subdiagonal.
A is a two-dimensional REAL array, dimensioned A(NM,IGH).
INT contains information on the rows and columns interchanged
in the reduction by ELMHES. Only elements LOW through IGH
are used. INT is a one-dimensional INTEGER array,
dimensioned INT(IGH).
Z contains the transformation matrix produced in the reduction
by ELMHES. Z is a two-dimensional REAL array, dimensioned
Z(NM,N).
Questions and comments should be directed to B. S. Garbow,
APPLIED MATHEMATICS DIVISION, ARGONNE NATIONAL LABORATORY
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer | :: | Nm | ||||
integer | :: | n | ||||
integer | :: | Low | ||||
integer | :: | Igh | ||||
real(kind=wp) | :: | a(Nm,*) | ||||
integer | :: | Int(*) | ||||
real(kind=wp) | :: | z(Nm,*) |
subroutine eltran(Nm, n, Low, Igh, a, Int, z) implicit none integer i, j, n, kl, mm, mp, Nm, Igh, Low, mp1 real(wp) a(Nm, *), z(Nm, *) integer Int(*) do i = 1, n do j = 1, n z(i, j) = 0.0_wp enddo z(i, i) = 1.0_wp enddo kl = Igh - Low - 1 if (kl >= 1) then ! for mp=igh-1 step -1 until low+1 do -- do mm = 1, kl mp = Igh - mm mp1 = mp + 1 do i = mp1, Igh z(i, mp) = a(i, mp - 1) enddo i = Int(mp) if (i /= mp) then do j = mp, Igh z(mp, j) = z(i, j) z(i, j) = 0.0_wp enddo z(i, mp) = 1.0_wp endif enddo endif end subroutine eltran