kmesh_write Subroutine

public subroutine kmesh_write()

Uses

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

Writes nnkp file (list of overlaps needed)

Arguments

None

Calls

proc~~kmesh_write~~CallsGraph proc~kmesh_write kmesh_write proc~io_date io_date proc~kmesh_write->proc~io_date proc~io_file_unit io_file_unit proc~kmesh_write->proc~io_file_unit

Called by

proc~~kmesh_write~~CalledByGraph proc~kmesh_write kmesh_write proc~wannier_setup wannier_setup proc~wannier_setup->proc~kmesh_write

Contents

Source Code


Source Code

  subroutine kmesh_write()
    !==================================================================!
    !                                                                  !
    !! Writes nnkp file (list of overlaps needed)
    !                                                                  !
    ! Note that the format is different to (and more compact than)     !
    ! that used by the old f77 code.                                   !
    !                                                                  !
    ! The file consists of num_kpts blocks of data, one block for each !
    ! k-point of the mesh. Each block consists of nntot+1 lines,       !
    ! where nntot is the (integer) number of nearest neighbours        !
    ! belonging to k-point nkp.                                        !
    !                                                                  !
    ! The first line in each block is just nntot.                      !
    !                                                                  !
    ! The second line consists of 5 integers. The first is the k-point !
    ! nkp. The second to the fifth specify it's nearest neighbours     !
    ! k+b: the second integer points to the k-point that is the        !
    ! periodic image of k+b that we want; the last three integers give !
    ! the G-vector, in reciprocal lattice units, that brings the       !
    ! k-point specified by the second integer (which is in the first   !
    ! BZ) to the actual k+b that we need.                              !
    !                                                                  !
    ! So wannier.nnkp specifies the nearest neighbours of each         !
    ! k-point, and therefore provides the information required to      !
    ! calculate the M_mn(k,b) matrix elements -- Marzari & Vanderbilt  !
    ! PRB 56, 12847 (1997) Eq. (25) -- for each pair of band indices   !
    ! m and n.                                                         !
    !===================================================================
    use w90_io, only: io_file_unit, seedname, io_date, io_stopwatch

    implicit none

    integer           :: i, nkp, nn, nnkpout
    character(len=9) :: cdate, ctime

    if (timing_level > 0) call io_stopwatch('kmesh: write', 1)

    nnkpout = io_file_unit()
    open (unit=nnkpout, file=trim(seedname)//'.nnkp', form='formatted')

    ! Date and time
    call io_date(cdate, ctime)
    write (nnkpout, '(4(a),/)') 'File written on ', cdate, ' at ', ctime

    ! Calc_only_A
    write (nnkpout, '(a,l2,/)') 'calc_only_A  : ', calc_only_A

    ! Real lattice
    write (nnkpout, '(a)') 'begin real_lattice'
    write (nnkpout, '(3(f12.7))') (real_lattice(1, i), i=1, 3)
    write (nnkpout, '(3(f12.7))') (real_lattice(2, i), i=1, 3)
    write (nnkpout, '(3(f12.7))') (real_lattice(3, i), i=1, 3)
    write (nnkpout, '(a/)') 'end real_lattice'

    ! Reciprocal lattice
    write (nnkpout, '(a)') 'begin recip_lattice'
    write (nnkpout, '(3f12.7)') (recip_lattice(1, i), i=1, 3)
    write (nnkpout, '(3f12.7)') (recip_lattice(2, i), i=1, 3)
    write (nnkpout, '(3f12.7)') (recip_lattice(3, i), i=1, 3)
    write (nnkpout, '(a/)') 'end recip_lattice'

    ! K-points
    write (nnkpout, '(a)') 'begin kpoints'
    write (nnkpout, '(i6)') num_kpts
    do nkp = 1, num_kpts
      write (nnkpout, '(3f14.8)') (kpt_latt(i, nkp), i=1, 3)
    enddo
    write (nnkpout, '(a/)') 'end kpoints'

    if (spinors) then
      ! Projections
      write (nnkpout, '(a)') 'begin spinor_projections'
      if (allocated(input_proj_site)) then
        write (nnkpout, '(i6)') num_proj
        do i = 1, num_proj
          write (nnkpout, '(3(f10.5,1x),2x,3i3)') &
            input_proj_site(1, i), input_proj_site(2, i), input_proj_site(3, i), &
            input_proj_l(i), input_proj_m(i), input_proj_radial(i)
!~           write(nnkpout,'(3x,3f7.3,1x,3f7.3,1x,f7.2)') &
          write (nnkpout, '(2x,3f11.7,1x,3f11.7,1x,f7.2)') &
            input_proj_z(1, i), input_proj_z(2, i), input_proj_z(3, i), &
            input_proj_x(1, i), input_proj_x(2, i), input_proj_x(3, i), &
            input_proj_zona(i)
          write (nnkpout, '(2x,1i3,1x,3f11.7)') &
            input_proj_s(i), &
            input_proj_s_qaxis(1, i), input_proj_s_qaxis(2, i), input_proj_s_qaxis(3, i)
        enddo
      else
        ! No projections
        write (nnkpout, '(i6)') 0
      end if
      write (nnkpout, '(a/)') 'end spinor_projections'
    else
      ! Projections
      write (nnkpout, '(a)') 'begin projections'
      if (allocated(input_proj_site)) then
        write (nnkpout, '(i6)') num_proj
        do i = 1, num_proj
          write (nnkpout, '(3(f10.5,1x),2x,3i3)') &
            input_proj_site(1, i), input_proj_site(2, i), input_proj_site(3, i), &
            input_proj_l(i), input_proj_m(i), input_proj_radial(i)
!~           write(nnkpout,'(3x,3f7.3,1x,3f7.3,1x,f7.2)') &
          write (nnkpout, '(2x,3f11.7,1x,3f11.7,1x,f7.2)') &
            input_proj_z(1, i), input_proj_z(2, i), input_proj_z(3, i), &
            input_proj_x(1, i), input_proj_x(2, i), input_proj_x(3, i), &
            input_proj_zona(i)
        enddo
      else
        ! No projections
        write (nnkpout, '(i6)') 0
      end if
      write (nnkpout, '(a/)') 'end projections'
    endif

    ! Info for automatic generation of projections
    if (auto_projections) then
      write (nnkpout, '(a)') 'begin auto_projections'
      write (nnkpout, '(i6)') num_proj
      write (nnkpout, '(i6)') 0
      write (nnkpout, '(a/)') 'end auto_projections'
    end if

    ! Nearest neighbour k-points
    write (nnkpout, '(a)') 'begin nnkpts'
    write (nnkpout, '(i4)') nntot
    do nkp = 1, num_kpts
      do nn = 1, nntot
        write (nnkpout, '(2i6,3x,3i4)') &
          nkp, nnlist(nkp, nn), (nncell(i, nkp, nn), i=1, 3)
      end do
    end do
    write (nnkpout, '(a/)') 'end nnkpts'

    !states to exclude
    write (nnkpout, '(a)') 'begin exclude_bands'
    write (nnkpout, '(i4)') num_exclude_bands
    if (num_exclude_bands > 0) then
      do i = 1, num_exclude_bands
        write (nnkpout, '(i4)') exclude_bands(i)
      end do
    endif
    write (nnkpout, '(a)') 'end exclude_bands'

    close (nnkpout)

    if (timing_level > 0) call io_stopwatch('kmesh: write', 2)

    return

  end subroutine kmesh_write