param_write_chkpt Subroutine

public subroutine param_write_chkpt(chkpt)

Uses

  • proc~~param_write_chkpt~~UsesGraph proc~param_write_chkpt param_write_chkpt module~w90_io w90_io proc~param_write_chkpt->module~w90_io module~w90_constants w90_constants module~w90_io->module~w90_constants

Write checkpoint file IMPORTANT! If you change the chkpt format, adapt accordingly also the w90chk2chk.x utility! Also, note that this routine writes the u_matrix and the m_matrix - in parallel mode these are however stored in distributed form in, e.g., u_matrix_loc only, so if you are changing the u_matrix, remember to gather it from u_matrix_loc first!

Arguments

Type IntentOptional AttributesName
character(len=*), intent(in) :: chkpt

Calls

proc~~param_write_chkpt~~CallsGraph proc~param_write_chkpt param_write_chkpt proc~io_date io_date proc~param_write_chkpt->proc~io_date proc~io_file_unit io_file_unit proc~param_write_chkpt->proc~io_file_unit

Called by

proc~~param_write_chkpt~~CalledByGraph proc~param_write_chkpt param_write_chkpt proc~wann_main_gamma wann_main_gamma proc~wann_main_gamma->proc~param_write_chkpt proc~wannier_run wannier_run proc~wannier_run->proc~param_write_chkpt proc~wannier_run->proc~wann_main_gamma program~wannier wannier program~wannier->proc~wann_main_gamma

Contents

Source Code


Source Code

  subroutine param_write_chkpt(chkpt)
    !=================================================!
    !! Write checkpoint file
    !! IMPORTANT! If you change the chkpt format, adapt
    !! accordingly also the w90chk2chk.x utility!
    !! Also, note that this routine writes the u_matrix and the m_matrix - in parallel
    !! mode these are however stored in distributed form in, e.g., u_matrix_loc only, so
    !! if you are changing the u_matrix, remember to gather it from u_matrix_loc first!
    !=================================================!

    use w90_io, only: io_file_unit, io_date, seedname

    implicit none

    character(len=*), intent(in) :: chkpt

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

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

    call io_date(cdate, ctime)
    header = 'written on '//cdate//' at '//ctime

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

    write (chk_unit) header                                   ! Date and time
    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
    chkpt1 = adjustl(trim(chkpt))
    write (chk_unit) chkpt1                 ! 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'

    return

  end subroutine param_write_chkpt