conv_write_chkpt Subroutine

public subroutine conv_write_chkpt()

Uses

  • proc~~conv_write_chkpt~~UsesGraph proc~conv_write_chkpt conv_write_chkpt module~w90_parameters w90_parameters proc~conv_write_chkpt->module~w90_parameters module~w90_io w90_io proc~conv_write_chkpt->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~~CallsGraph proc~conv_write_chkpt conv_write_chkpt proc~io_file_unit io_file_unit proc~conv_write_chkpt->proc~io_file_unit

Called by

proc~~conv_write_chkpt~~CalledByGraph proc~conv_write_chkpt conv_write_chkpt program~w90chk2chk w90chk2chk program~w90chk2chk->proc~conv_write_chkpt

Contents

Source Code


Source Code

  subroutine conv_write_chkpt()
    !=======================================!
    !! 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 checkpoint file ', trim(seedname), '.chk...'

    chk_unit = io_file_unit()
    open (unit=chk_unit, file=trim(seedname)//'.chk', form='unformatted')

    write (chk_unit) header                                   ! Date and time from the read file
    write (chk_unit) num_bands                                ! Number of bands
    write (chk_unit) num_exclude_bands                        ! Number of excluded bands
    write (chk_unit) (exclude_bands(i), i=1, num_exclude_bands) ! Excluded bands
    write (chk_unit) ((real_lattice(i, j), i=1, 3), j=1, 3)        ! Real lattice
    write (chk_unit) ((recip_lattice(i, j), i=1, 3), j=1, 3)       ! Reciprocal lattice
    write (chk_unit) num_kpts                                 ! Number of k-points
    write (chk_unit) (mp_grid(i), i=1, 3)                       ! M-P grid
    write (chk_unit) ((kpt_latt(i, nkp), i=1, 3), nkp=1, num_kpts) ! K-points
    write (chk_unit) nntot                                    ! Number of nearest k-point neighbours
    write (chk_unit) num_wann                                 ! Number of wannier functions
    ! Next is correct: it always print out 20 characters
    write (chk_unit) checkpoint                               ! Position of checkpoint
    write (chk_unit) have_disentangled      ! Whether a disentanglement has been performed
    if (have_disentangled) then
      write (chk_unit) omega_invariant     ! Omega invariant
      ! lwindow, ndimwin and U_matrix_opt
      write (chk_unit) ((lwindow(i, nkp), i=1, num_bands), nkp=1, num_kpts)
      write (chk_unit) (ndimwin(nkp), nkp=1, num_kpts)
      write (chk_unit) (((u_matrix_opt(i, j, nkp), i=1, num_bands), j=1, num_wann), nkp=1, num_kpts)
    endif
    write (chk_unit) (((u_matrix(i, j, k), i=1, num_wann), j=1, num_wann), k=1, num_kpts)               ! U_matrix
    write (chk_unit) ((((m_matrix(i, j, k, l), i=1, num_wann), j=1, num_wann), k=1, nntot), l=1, num_kpts) ! M_matrix
    write (chk_unit) ((wannier_centres(i, j), i=1, 3), j=1, num_wann)
    write (chk_unit) (wannier_spreads(i), i=1, num_wann)
    close (chk_unit)

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

  end subroutine conv_write_chkpt