Read formatted spn file
subroutine conv_read_spn_fmt()
!=======================================!
!! Read formatted spn file
!=======================================!
use w90_constants, only: eps6, dp
use w90_io, only: io_error, io_file_unit, stdout, seedname
use w90_parameters, only: num_bands, num_kpts
implicit none
integer :: spn_unit, m, n, ik, ierr
real(kind=dp) :: s_real, s_img
write (stdout, '(3a)') 'Reading information from formatted file ', trim(seedname), '.spn.fmt :'
spn_unit = io_file_unit()
open (unit=spn_unit, file=trim(seedname)//'.spn.fmt', status='old', position='rewind', form='formatted', err=109)
! Read comment line
read (spn_unit, '(a)') header
header = ADJUSTL(header)
write (stdout, '(1x,a)') trim(header)
! Consistency checks
read (spn_unit, *, err=110, end=110) num_bands, num_kpts
write (stdout, '(1x,a,i0)') "Number of bands: ", num_bands
write (stdout, '(1x,a,i0)') "Number of k-points: ", num_kpts
allocate (spn_o(num_bands, num_bands, num_kpts, 3), stat=ierr)
if (ierr /= 0) call io_error('Error in allocating spn_o in conv_read_spn_fmt')
do ik = 1, num_kpts
do m = 1, num_bands
do n = 1, m
read (spn_unit, *, err=110, end=110) s_real, s_img
spn_o(n, m, ik, 1) = cmplx(s_real, s_img, dp)
read (spn_unit, *, err=110, end=110) s_real, s_img
spn_o(n, m, ik, 2) = cmplx(s_real, s_img, dp)
read (spn_unit, *, err=110, end=110) s_real, s_img
spn_o(n, m, ik, 3) = cmplx(s_real, s_img, dp)
! Although each diagonal element of spn_o should be a real number,
! actually it has a very small imaginary part.
! We skip the conjugation on the diagonal elements so that
! the file after formatted <==> unformatted conversions is exactly
! the same as the original file, otherwise the diagonal elements
! are the conjugations of those of the original file.
if (m == n) cycle
! Read upper-triangular part, now build the rest
spn_o(m, n, ik, 1) = conjg(spn_o(n, m, ik, 1))
spn_o(m, n, ik, 2) = conjg(spn_o(n, m, ik, 2))
spn_o(m, n, ik, 3) = conjg(spn_o(n, m, ik, 3))
end do
end do
enddo
close (spn_unit)
write (stdout, '(1x,a)') 'read done.'
return
109 call io_error('Error opening '//trim(seedname)//'.spn.fmt in conv_read_spn_fmt')
110 call io_error('Error reading '//trim(seedname)//'.spn.fmt in conv_read_spn_fmt')
end subroutine conv_read_spn_fmt