conv_write_chkpt_fmt Subroutine

public subroutine conv_write_chkpt_fmt()

Uses

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

Write formatted checkpoint file

Arguments

None

Calls

proc~~conv_write_chkpt_fmt~~CallsGraph proc~conv_write_chkpt_fmt conv_write_chkpt_fmt proc~io_file_unit io_file_unit proc~conv_write_chkpt_fmt->proc~io_file_unit

Called by

proc~~conv_write_chkpt_fmt~~CalledByGraph proc~conv_write_chkpt_fmt conv_write_chkpt_fmt program~w90chk2chk w90chk2chk program~w90chk2chk->proc~conv_write_chkpt_fmt

Contents

Source Code


Source Code

  subroutine conv_write_chkpt_fmt()
    !=======================================!
    !! Write formatted checkpoint file
    !=======================================!

    use w90_io, only: io_file_unit, io_date, seedname
    use w90_parameters

    implicit none

    integer :: chk_unit, nkp, i, j, k, l
    character(len=9) :: cdate, ctime

    write (stdout, '(/1x,3a)', advance='no') 'Writing formatted checkpoint file ', trim(seedname), '.chk.fmt...'

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

    write (chk_unit, '(A33)') header                                   ! Date and time from the read file
    write (chk_unit, '(I0)') num_bands                                ! Number of bands
    write (chk_unit, '(I0)') num_exclude_bands                        ! Number of excluded bands
    do i = 1, num_exclude_bands
      write (chk_unit, '(I0)') exclude_bands(i) ! Excluded bands
    end do
    write (chk_unit, '(9G25.17)') ((real_lattice(i, j), i=1, 3), j=1, 3)        ! Real lattice
    write (chk_unit, '(9G25.17)') ((recip_lattice(i, j), i=1, 3), j=1, 3)       ! Reciprocal lattice
    write (chk_unit, '(I0)') num_kpts                                 ! Number of k-points
    write (chk_unit, '(I0," ",I0," ",I0)') (mp_grid(i), i=1, 3)                       ! M-P grid
    do nkp = 1, num_kpts
      write (chk_unit, '(3G25.17)') (kpt_latt(i, nkp), i=1, 3) ! K-points
    end do
    write (chk_unit, '(I0)') nntot                                    ! Number of nearest k-point neighbours
    write (chk_unit, '(I0)') num_wann                                 ! Number of wannier functions
    write (chk_unit, '(A20)') checkpoint                               ! Position of checkpoint
    if (have_disentangled) then
      write (chk_unit, '(I1)') 1      ! Whether a disentanglement has been performed
    else
      write (chk_unit, '(I1)') 0      ! Whether a disentanglement has been performed
    end if
    if (have_disentangled) then
      write (chk_unit, '(G25.17)') omega_invariant     ! Omega invariant
      ! lwindow, ndimwin and U_matrix_opt
      do nkp = 1, num_kpts
        do i = 1, num_bands
          if (lwindow(i, nkp)) then
            write (chk_unit, '(I1)') 1
          else
            write (chk_unit, '(I1)') 0
          end if
        end do
      end do
      do nkp = 1, num_kpts
        write (chk_unit, '(I0)') ndimwin(nkp)
      end do
      do nkp = 1, num_kpts
        do j = 1, num_wann
          do i = 1, num_bands
            write (chk_unit, '(2G25.17)') u_matrix_opt(i, j, nkp)
          end do
        end do
      end do
    endif
    do k = 1, num_kpts
      do j = 1, num_wann
        do i = 1, num_wann
          write (chk_unit, '(2G25.17)') u_matrix(i, j, k)
        end do
      end do
    end do
    do l = 1, num_kpts
      do k = 1, nntot
        do j = 1, num_wann
          do i = 1, num_wann
            write (chk_unit, '(2G25.17)') m_matrix(i, j, k, l)
          end do
        end do
      end do
    end do
    do j = 1, num_wann
      write (chk_unit, '(3G25.17)') (wannier_centres(i, j), i=1, 3)
    end do
    do i = 1, num_wann
      write (chk_unit, '(G25.17)') wannier_spreads(i)
    end do
    close (chk_unit)

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

  end subroutine conv_write_chkpt_fmt