r1mpyq Subroutine

public subroutine r1mpyq(m, n, a, Lda, v, w)

given an m by n matrix a, this subroutine computes aq where q is the product of 2(n - 1) transformations

        gv(n-1)*...*gv(1)*gw(1)*...*gw(n-1)

and gv(i), gw(i) are givens rotations in the (i,n) plane which eliminate elements in the i-th and n-th planes, respectively. q itself is not given, rather the information to recover the gv, gw rotations is supplied.

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: m

a positive integer input variable set to the number of rows of a.

integer, intent(in) :: n

a positive integer input variable set to the number of columns of a.

real(kind=wp), intent(inout) :: a(Lda,n)

an m by n array. on input a must contain the matrix to be postmultiplied by the orthogonal matrix q described above. on output a*q has replaced a.

integer, intent(in) :: Lda

a positive integer input variable not less than m which specifies the leading dimension of the array a.

real(kind=wp), intent(in) :: v(n)

an input array of length n. v(i) must contain the information necessary to recover the givens rotation gv(i) described above.

real(kind=wp), intent(in) :: w(n)

an input array of length n. w(i) must contain the information necessary to recover the givens rotation gw(i) described above.


Called by

proc~~r1mpyq~~CalledByGraph proc~r1mpyq r1mpyq proc~hybrd hybrd proc~hybrd->proc~r1mpyq proc~hybrj hybrj proc~hybrj->proc~r1mpyq proc~minpack_hybrd minpack_hybrd proc~minpack_hybrd->proc~hybrd proc~hybrj1 hybrj1 proc~hybrj1->proc~hybrj proc~hybrd1 hybrd1 proc~hybrd1->proc~hybrd proc~minpack_hybrj minpack_hybrj proc~minpack_hybrj->proc~hybrj proc~minpack_hybrj1 minpack_hybrj1 proc~minpack_hybrj1->proc~hybrj1 proc~minpack_hybrd1 minpack_hybrd1 proc~minpack_hybrd1->proc~hybrd1

Contents

Source Code


Source Code

    subroutine r1mpyq(m, n, a, Lda, v, w)
        implicit none

        integer, intent(in) :: m !! a positive integer input variable set to the number
                                !! of rows of a.
        integer, intent(in) :: n !! a positive integer input variable set to the number
                                !! of columns of a.
        integer, intent(in) :: Lda !! a positive integer input variable not less than m
                                  !! which specifies the leading dimension of the array a.
        real(wp), intent(inout) :: a(Lda, n) !! an m by n array. on input a must contain the matrix
                                            !! to be postmultiplied by the orthogonal matrix q
                                            !! described above. on output a*q has replaced a.
        real(wp), intent(in) :: v(n) !! an input array of length n. v(i) must contain the
                                    !! information necessary to recover the givens rotation gv(i)
                                    !! described above.
        real(wp), intent(in) :: w(n) !! an input array of length n. w(i) must contain the
                                    !! information necessary to recover the givens rotation gw(i)
                                    !! described above.

        integer :: i, j, nmj, nm1
        real(wp) :: cos, sin, temp

        ! apply the first set of givens rotations to a.

        nm1 = n - 1
        if (nm1 >= 1) then
            do nmj = 1, nm1
                j = n - nmj
                if (abs(v(j)) > one) then
                    cos = one/v(j)
                    sin = sqrt(one - cos**2)
                else
                    sin = v(j)
                    cos = sqrt(one - sin**2)
                end if
                do i = 1, m
                    temp = cos*a(i, j) - sin*a(i, n)
                    a(i, n) = sin*a(i, j) + cos*a(i, n)
                    a(i, j) = temp
                end do
            end do

            ! apply the second set of givens rotations to a.

            do j = 1, nm1
                if (abs(w(j)) > one) cos = one/w(j)
                if (abs(w(j)) > one) sin = sqrt(one - cos**2)
                if (abs(w(j)) <= one) sin = w(j)
                if (abs(w(j)) <= one) cos = sqrt(one - sin**2)
                do i = 1, m
                    temp = cos*a(i, j) + sin*a(i, n)
                    a(i, n) = -sin*a(i, j) + cos*a(i, n)
                    a(i, j) = temp
                end do
            end do
        end if

    end subroutine r1mpyq