conv_read_chkpt_fmt Subroutine

public subroutine conv_read_chkpt_fmt()

Uses

  • proc~~conv_read_chkpt_fmt~~UsesGraph proc~conv_read_chkpt_fmt conv_read_chkpt_fmt module~w90_constants w90_constants proc~conv_read_chkpt_fmt->module~w90_constants module~w90_parameters w90_parameters proc~conv_read_chkpt_fmt->module~w90_parameters module~w90_io w90_io proc~conv_read_chkpt_fmt->module~w90_io module~w90_parameters->module~w90_constants module~w90_parameters->module~w90_io module~w90_io->module~w90_constants

Read formatted checkpoint file

Arguments

None

Calls

proc~~conv_read_chkpt_fmt~~CallsGraph proc~conv_read_chkpt_fmt conv_read_chkpt_fmt proc~io_error io_error proc~conv_read_chkpt_fmt->proc~io_error proc~io_file_unit io_file_unit proc~conv_read_chkpt_fmt->proc~io_file_unit

Called by

proc~~conv_read_chkpt_fmt~~CalledByGraph proc~conv_read_chkpt_fmt conv_read_chkpt_fmt program~w90chk2chk w90chk2chk program~w90chk2chk->proc~conv_read_chkpt_fmt

Contents

Source Code


Source Code

  subroutine conv_read_chkpt_fmt()
    !=======================================!
    !! Read formatted checkpoint file
    !=======================================!

    use w90_constants, only: eps6
    use w90_io, only: io_error, io_file_unit, stdout, seedname
    use w90_parameters

    implicit none

    integer :: chk_unit, i, j, k, l, nkp, ierr, idum
    character(len=30) :: cdum
    real(kind=dp) :: rreal, rimag

    write (stdout, '(1x,3a)') 'Reading information from formatted file ', trim(seedname), '.chk.fmt :'

    chk_unit = io_file_unit()
    open (unit=chk_unit, file=trim(seedname)//'.chk.fmt', status='old', position='rewind', form='formatted', err=121)

    ! Read comment line
    read (chk_unit, '(A)') header
    write (stdout, '(1x,a)') trim(header)

    ! Consistency checks
    read (chk_unit, *) num_bands                           ! Number of bands
    write (stdout, '(a,i0)') "Number of bands: ", num_bands
    read (chk_unit, *) num_exclude_bands                   ! Number of excluded bands
    if (num_exclude_bands < 0) then
      call io_error('Invalid value for num_exclude_bands')
    endif
    allocate (exclude_bands(num_exclude_bands), stat=ierr)
    if (ierr /= 0) call io_error('Error allocating exclude_bands in conv_read_chkpt_fmt')
    do i = 1, num_exclude_bands
      read (chk_unit, *) exclude_bands(i) ! Excluded bands
    end do
    write (stdout, '(a)', advance='no') "Excluded bands: "
    if (num_exclude_bands == 0) then
      write (stdout, '(a)') "none."
    else
      do i = 1, num_exclude_bands - 1
        write (stdout, '(I0,a)', advance='no') exclude_bands(i), ','
      end do
      write (stdout, '(I0,a)') exclude_bands(num_exclude_bands), '.'
    end if
    read (chk_unit, *) ((real_lattice(i, j), i=1, 3), j=1, 3)  ! Real lattice
    write (stdout, '(a)') "Real lattice: read."
    read (chk_unit, *) ((recip_lattice(i, j), i=1, 3), j=1, 3)  ! Reciprocal lattice
    write (stdout, '(a)') "Reciprocal lattice: read."
    read (chk_unit, *) num_kpts                ! K-points
    write (stdout, '(a,I0)') "Num kpts:", num_kpts
    read (chk_unit, *) (mp_grid(i), i=1, 3)         ! M-P grid
    write (stdout, '(a)') "mp_grid: read."
    if (.not. allocated(kpt_latt)) then
      allocate (kpt_latt(3, num_kpts), stat=ierr)
      if (ierr /= 0) call io_error('Error allocating kpt_latt in conv_read_chkpt_fmt')
    endif
    do nkp = 1, num_kpts
      read (chk_unit, *, err=115) (kpt_latt(i, nkp), i=1, 3)
    end do
    write (stdout, '(a)') "kpt_latt: read."
    read (chk_unit, *) nntot                ! nntot
    write (stdout, '(a,I0)') "nntot:", nntot
    read (chk_unit, *) num_wann                ! num_wann
    write (stdout, '(a,I0)') "num_wann:", num_wann

    read (chk_unit, *) checkpoint             ! checkpoint
    checkpoint = adjustl(trim(checkpoint))
    write (stdout, '(a,I0)') "checkpoint: "//trim(checkpoint)

    read (chk_unit, *) idum
    if (idum == 1) then
      have_disentangled = .true.
    elseif (idum == 0) then
      have_disentangled = .false.
    else
      write (cdum, '(I0)') idum
      call io_error('Error reading formatted chk: have_distenangled should be 0 or 1, it is instead '//cdum)
    end if

    if (have_disentangled) then
      write (stdout, '(a)') "have_disentangled: TRUE"

      read (chk_unit, *) omega_invariant     ! omega invariant
      write (stdout, '(a)') "omega_invariant: read."

      ! lwindow
      if (.not. allocated(lwindow)) then
        allocate (lwindow(num_bands, num_kpts), stat=ierr)
        if (ierr /= 0) call io_error('Error allocating lwindow in conv_read_chkpt_fmt')
      endif
      do nkp = 1, num_kpts
        do i = 1, num_bands
          read (chk_unit, *) idum
          if (idum == 1) then
            lwindow(i, nkp) = .true.
          elseif (idum == 0) then
            lwindow(i, nkp) = .false.
          else
            write (cdum, '(I0)') idum
            call io_error('Error reading formatted chk: lwindow(i,nkp) should be 0 or 1, it is instead '//cdum)
          end if
        end do
      end do
      write (stdout, '(a)') "lwindow: read."

      ! ndimwin
      if (.not. allocated(ndimwin)) then
        allocate (ndimwin(num_kpts), stat=ierr)
        if (ierr /= 0) call io_error('Error allocating ndimwin in conv_read_chkpt_fmt')
      endif
      do nkp = 1, num_kpts
        read (chk_unit, *, err=123) ndimwin(nkp)
      end do
      write (stdout, '(a)') "ndimwin: read."

      ! U_matrix_opt
      if (.not. allocated(u_matrix_opt)) then
        allocate (u_matrix_opt(num_bands, num_wann, num_kpts), stat=ierr)
        if (ierr /= 0) call io_error('Error allocating u_matrix_opt in conv_read_chkpt_fmt')
      endif
      do nkp = 1, num_kpts
        do j = 1, num_wann
          do i = 1, num_bands
            read (chk_unit, *, err=124) rreal, rimag
            u_matrix_opt(i, j, nkp) = cmplx(rreal, rimag, kind=dp)
          end do
        end do
      end do
      write (stdout, '(a)') "U_matrix_opt: read."

    else
      write (stdout, '(a)') "have_disentangled: FALSE"
    endif

    ! U_matrix
    if (.not. allocated(u_matrix)) then
      allocate (u_matrix(num_wann, num_wann, num_kpts), stat=ierr)
      if (ierr /= 0) call io_error('Error allocating u_matrix in conv_read_chkpt_fmt')
    endif
    do k = 1, num_kpts
      do j = 1, num_wann
        do i = 1, num_wann
          read (chk_unit, *, err=124) rreal, rimag
          u_matrix(i, j, k) = cmplx(rreal, rimag, kind=dp)
        end do
      end do
    end do
    write (stdout, '(a)') "U_matrix: read."

    ! M_matrix
    if (.not. allocated(m_matrix)) then
      allocate (m_matrix(num_wann, num_wann, nntot, num_kpts), stat=ierr)
      if (ierr /= 0) call io_error('Error allocating m_matrix in conv_read_chkpt_fmt')
    endif
    do l = 1, num_kpts
      do k = 1, nntot
        do j = 1, num_wann
          do i = 1, num_wann
            read (chk_unit, *, err=124) rreal, rimag
            m_matrix(i, j, k, l) = cmplx(rreal, rimag, kind=dp)
          end do
        end do
      end do
    end do
    write (stdout, '(a)') "M_matrix: read."

    ! wannier_centres
    if (.not. allocated(wannier_centres)) then
      allocate (wannier_centres(3, num_wann), stat=ierr)
      if (ierr /= 0) call io_error('Error allocating wannier_centres in conv_read_chkpt_fmt')
    end if
    do j = 1, num_wann
      read (chk_unit, *, err=127) (wannier_centres(i, j), i=1, 3)
    end do
    write (stdout, '(a)') "wannier_centres: read."

    ! wannier spreads
    if (.not. allocated(wannier_spreads)) then
      allocate (wannier_spreads(num_wann), stat=ierr)
      if (ierr /= 0) call io_error('Error allocating wannier_centres in conv_read_chkpt_fmt')
    end if
    do i = 1, num_wann
      read (chk_unit, *, err=128) wannier_spreads(i)
    end do
    write (stdout, '(a)') "wannier_spreads: read."

    close (chk_unit)

    write (stdout, '(a/)') ' ... done'

    return

115 call io_error('Error reading variable from '//trim(seedname)//'.chk.fmt in conv_read_chkpt_fmt')
121 call io_error('Error opening '//trim(seedname)//'.chk.fmt in conv_read_chkpt_fmt')
122 call io_error('Error reading lwindow from '//trim(seedname)//'.chk.fmt in conv_read_chkpt_fmt')
123 call io_error('Error reading ndimwin from '//trim(seedname)//'.chk.fmt in conv_read_chkpt_fmt')
124 call io_error('Error reading u_matrix_opt from '//trim(seedname)//'.chk.fmt in conv_read_chkpt_fmt')
125 call io_error('Error reading u_matrix from '//trim(seedname)//'.chk.fmt in conv_read_chkpt_fmt')
126 call io_error('Error reading m_matrix from '//trim(seedname)//'.chk.fmt in conv_read_chkpt_fmt')
127 call io_error('Error reading wannier_centres from '//trim(seedname)//'.chk.fmt in conv_read_chkpt_fmt')
128 call io_error('Error reading wannier_spreads from '//trim(seedname)//'.chk.fmt in conv_read_chkpt_fmt')

  end subroutine conv_read_chkpt_fmt