ws_write_vec Subroutine

public subroutine ws_write_vec(nrpts, irvec)

Uses

  • proc~~ws_write_vec~~UsesGraph proc~ws_write_vec ws_write_vec module~w90_parameters w90_parameters proc~ws_write_vec->module~w90_parameters module~w90_io w90_io proc~ws_write_vec->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 to file the lattice vectors of the superlattice to be added to R vector in seedname_hr.dat, seedname_rmn.dat, etc. in order to have the second Wannier function inside the WS cell of the first one.

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: nrpts
integer, intent(in) :: irvec(3,nrpts)

Calls

proc~~ws_write_vec~~CallsGraph proc~ws_write_vec ws_write_vec proc~io_date io_date proc~ws_write_vec->proc~io_date proc~io_file_unit io_file_unit proc~ws_write_vec->proc~io_file_unit

Called by

proc~~ws_write_vec~~CalledByGraph proc~ws_write_vec ws_write_vec proc~plot_main plot_main proc~plot_main->proc~ws_write_vec program~wannier wannier program~wannier->proc~plot_main proc~wannier_run wannier_run proc~wannier_run->proc~plot_main

Contents

Source Code


Source Code

  subroutine ws_write_vec(nrpts, irvec)
    !! Write to file the lattice vectors of the superlattice
    !! to be added to R vector in seedname_hr.dat, seedname_rmn.dat, etc.
    !! in order to have the second Wannier function inside the WS cell
    !! of the first one.

    use w90_io, only: io_error, io_stopwatch, io_file_unit, &
      seedname, io_date
    use w90_parameters, only: num_wann

    implicit none

    integer, intent(in) :: nrpts
    integer, intent(in) :: irvec(3, nrpts)
    integer:: irpt, iw, jw, ideg, file_unit
    character(len=100) :: header
    character(len=9)  :: cdate, ctime

    file_unit = io_file_unit()
    call io_date(cdate, ctime)

    open (file_unit, file=trim(seedname)//'_wsvec.dat', form='formatted', &
          status='unknown', err=101)

    if (use_ws_distance) then
      header = '## written on '//cdate//' at '//ctime//' with use_ws_distance=.true.'
      write (file_unit, '(A)') trim(header)

      do irpt = 1, nrpts
        do iw = 1, num_wann
          do jw = 1, num_wann
            write (file_unit, '(5I5)') irvec(:, irpt), iw, jw
            write (file_unit, '(I5)') wdist_ndeg(iw, jw, irpt)
            do ideg = 1, wdist_ndeg(iw, jw, irpt)
              write (file_unit, '(5I5,2F12.6,I5)') irdist_ws(:, ideg, iw, jw, irpt) - &
                irvec(:, irpt)
            end do
          end do
        end do
      end do
    else
      header = '## written on '//cdate//' at '//ctime//' with use_ws_distance=.false.'
      write (file_unit, '(A)') trim(header)

      do irpt = 1, nrpts
        do iw = 1, num_wann
          do jw = 1, num_wann
            write (file_unit, '(5I5)') irvec(:, irpt), &
              iw, jw
            write (file_unit, '(I5)') 1
            write (file_unit, '(3I5)') 0, 0, 0
          end do
        end do
      end do
    end if

    close (file_unit)
    return

101 call io_error('Error: ws_write_vec: problem opening file '//trim(seedname)//'_ws_vec.dat')
    !====================================================!
  end subroutine ws_write_vec