wigner_seitz Subroutine

private subroutine wigner_seitz(count_pts)

Uses

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

Calculates a grid of lattice vectors r that fall inside (and eventually on the surface of) the Wigner-Seitz supercell centered on the origin of the Bravais superlattice with primitive translations mp_grid(1)a_1, mp_grid(2)a_2, and mp_grid(3)*a_3

Arguments

Type IntentOptional AttributesName
logical, intent(in) :: count_pts

Called by

proc~~wigner_seitz~~CalledByGraph proc~wigner_seitz wigner_seitz proc~pw90common_wanint_setup pw90common_wanint_setup proc~pw90common_wanint_setup->proc~wigner_seitz program~postw90 postw90 program~postw90->proc~pw90common_wanint_setup

Contents

Source Code


Source Code

  subroutine wigner_seitz(count_pts)
    !================================!
    !! Calculates a grid of lattice vectors r that fall inside (and eventually
    !! on the surface of) the Wigner-Seitz supercell centered on the
    !! origin of the Bravais superlattice with primitive translations
    !! mp_grid(1)*a_1, mp_grid(2)*a_2, and mp_grid(3)*a_3
    !==========================================================================!

    use w90_constants, only: dp
    use w90_io, only: stdout, io_error, io_stopwatch
    use w90_parameters, only: mp_grid, real_metric, iprint, timing_level

    ! irvec(i,irpt)     The irpt-th Wigner-Seitz grid point has components
    !                   irvec(1:3,irpt) in the basis of the lattice vectors
    ! ndegen(irpt)      Weight of the irpt-th point is 1/ndegen(irpt)
    ! nrpts             number of Wigner-Seitz grid points

    logical, intent(in) :: count_pts

    integer       :: ndiff(3)
    real(kind=dp) :: dist(125), tot, dist_min
    integer       :: n1, n2, n3, i1, i2, i3, icnt, i, j, ir

    if (timing_level > 1 .and. on_root) &
      call io_stopwatch('postw90_common: wigner_seitz', 1)

    ! The Wannier functions live in a periodic supercell of the real space unit
    ! cell. This supercell is mp_grid(i) unit cells long along each primitive
    ! translation vector a_i of the unit cell
    !
    ! We loop over grid points r on a cell that is approx. 8 times
    ! larger than this "primitive supercell."
    !
    ! One of these points is in the W-S supercell if it is closer to R=0 than any
    ! of the other points R (where R are the translation vectors of the
    ! supercell). In practice it is sufficient to inspect only 125 R-points.

    ! In the end, nrpts contains the total number of grid points that have been
    ! found in the Wigner-Seitz cell

    nrpts = 0
    do n1 = -mp_grid(1), mp_grid(1)
      do n2 = -mp_grid(2), mp_grid(2)
        do n3 = -mp_grid(3), mp_grid(3)
          ! Loop over the 125 points R. R=0 corresponds to i1=i2=i3=0,
          ! or icnt=63
          icnt = 0
          do i1 = -2, 2
            do i2 = -2, 2
              do i3 = -2, 2
                icnt = icnt + 1
                ! Calculate distance squared |r-R|^2
                ndiff(1) = n1 - i1*mp_grid(1)
                ndiff(2) = n2 - i2*mp_grid(2)
                ndiff(3) = n3 - i3*mp_grid(3)
                dist(icnt) = 0.0_dp
                do i = 1, 3
                  do j = 1, 3
                    dist(icnt) = dist(icnt) + &
                                 real(ndiff(i), dp)*real_metric(i, j)*real(ndiff(j), dp)
                  enddo
                enddo
              enddo
            enddo
          enddo
          dist_min = minval(dist)
          if (abs(dist(63) - dist_min) .lt. 1.e-7_dp) then
            nrpts = nrpts + 1
            if (.not. count_pts) then
              ndegen(nrpts) = 0
              do i = 1, 125
                if (abs(dist(i) - dist_min) .lt. 1.e-7_dp) &
                  ndegen(nrpts) = ndegen(nrpts) + 1
              end do
              irvec(1, nrpts) = n1
              irvec(2, nrpts) = n2
              irvec(3, nrpts) = n3
              !
              ! Remember which grid point r is at the origin
              !
              if (n1 == 0 .and. n2 == 0 .and. n3 == 0) rpt_origin = nrpts
            endif
          end if

          !n3
        enddo
        !n2
      enddo
      !n1
    enddo
    !
    if (count_pts) then
      if (timing_level > 1 .and. on_root) &
        call io_stopwatch('postw90_common: wigner_seitz', 2)
      return
    end if

    if (iprint >= 3 .and. on_root) then
      write (stdout, '(1x,i4,a,/)') nrpts, &
        ' lattice points in Wigner-Seitz supercell:'
      do ir = 1, nrpts
        write (stdout, '(4x,a,3(i3,1x),a,i2)') '  vector ', irvec(1, ir), &
          irvec(2, ir), irvec(3, ir), '  degeneracy: ', ndegen(ir)
      enddo
    endif
    ! Check the "sum rule"
    tot = 0.0_dp
    do ir = 1, nrpts
      !
      ! Corrects weights in Fourier sums for R-vectors on the boundary of the
      ! W-S supercell
      !
      tot = tot + 1.0_dp/real(ndegen(ir), dp)
    enddo
    if (abs(tot - real(mp_grid(1)*mp_grid(2)*mp_grid(3), dp)) > 1.e-8_dp) &
      call io_error('ERROR in wigner_seitz: error in finding Wigner-Seitz points')

    if (timing_level > 1 .and. on_root) &
      call io_stopwatch('postw90_common: wigner_seitz', 2)

    return
  end subroutine wigner_seitz