Bug in dlasd4 in 3.4.2

Post here if you want to report a bug to the LAPACK team

Bug in dlasd4 in 3.4.2

Postby andreasnoackjensen » Wed Feb 20, 2013 3:51 am

This small program gives either sigma=NaN or the wrong answer 0.5 in LAPACK 3.4.2 but not in previous versions. Also, the answer changes if I comment out the write statements.

program test

implicit none

double precision :: v(5), z(5), delta(5), work(5), sigma, rho, dnrm2
integer :: info

v = (/0.1d0, 0.2d0, 0.3d0, 0.4d0, 0.5d0/)
z = (/0.1d0, 0.2d0, 0.3d0, 0.4d0, 0.5d0/)
rho = dnrm2(5_8, z, 1_8)
write(*,*) "z", z
z = z / rho
rho = rho**2
write(*,*) "v", v
write(*,*) "z", z
write(*,*) "rho", rho
call dlasd4(5_8, 5_8, v, z, delta, rho, sigma, work, info)
write(*,*) "sigma", sigma
write(*,*) "info", info
end program
andreasnoackjensen
 
Posts: 1
Joined: Wed Feb 20, 2013 3:40 am

Return to Bug report

Who is online

Users browsing this forum: No registered users and 1 guest