check_g2gf Subroutine

private subroutine check_g2gf()

FIXME : Add documentation

Arguments

None

Contents

Source Code


Source Code

  subroutine check_g2gf
    use file_utils, only: error_unit
    use mp, only: finish_mp, iproc, proc0
    use species, only : nspec
    use theta_grid, only: ntgrid
    use gs2_layouts, only: g_lo, gf_lo
    use gs2_layouts, only: ig_idx, ik_idx, it_idx, il_idx, ie_idx, is_idx
    use redistribute, only: gather, scatter, report_map_property

    implicit none

    integer :: iglo, ig, isgn, ik, it, il, ie, igf, is, ierr
    complex, dimension (:,:,:), allocatable :: gtmp
    complex, dimension (:,:,:,:,:,:), allocatable :: gftmp

    if (proc0) then
       ierr = error_unit()
    else
       ierr = 6
    end if

    allocate (gtmp(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc))
    allocate (gftmp(-ntgrid:ntgrid, 2, nspec, negrid, nlambda, gf_lo%llim_proc:gf_lo%ulim_alloc))

    ! gather check
    gtmp = 0.0
    do iglo = g_lo%llim_proc, g_lo%ulim_proc
       ik = ik_idx(g_lo,iglo)
       it = it_idx(g_lo,iglo)
       il = il_idx(g_lo,iglo)
       ie = ie_idx(g_lo,iglo)
       is = is_idx(g_lo,iglo)
       do isgn = 1, 2
          do ig = -ntgrid, ntgrid
             gtmp(ig,isgn,iglo) = rule(ig,isgn,ik,it,il,ie,is)
          end do
       end do
    end do
    call gather (g2gf, gtmp, gftmp, ntgrid)
    do igf = gf_lo%llim_proc, gf_lo%ulim_proc
       ik = ik_idx(gf_lo,igf)
       it = it_idx(gf_lo,igf)
       do is = 1, nspec
          do il = 1,nlambda
             do ie = 1,negrid
                do isgn = 1,2
                   do ig=-ntgrid,ntgrid
                      if (int(real(gftmp(ig,isgn,is,ie,il,igf))) /= rule(ig,isgn,ik,it,il,ie,is)) &
                           write (ierr,'(a,8i6)') 'ERROR: gather by g2gf broken!', iproc
                   end do
                end do
             end do
          end do
       end do
    end do
    if (proc0) write (ierr,'(a)') 'g2gf gather check done'

    ! scatter check
    gtmp = 0.0
    call scatter (g2gf, gftmp, gtmp, ntgrid)
    do iglo = g_lo%llim_proc, g_lo%ulim_proc
       ik = ik_idx(g_lo,iglo)
       it = it_idx(g_lo,iglo)
       il = il_idx(g_lo,iglo)
       ie = ie_idx(g_lo,iglo)
       is = is_idx(g_lo,iglo)
       do isgn = 1, 2
          do ig = -ntgrid, ntgrid
             if (gtmp(ig,isgn,iglo) /= rule(ig,isgn,ik,it,il,ie,is)) &
                  write (ierr,'(a,i6)') 'ERROR: scatter by g2gf broken!', iproc
          end do
       end do
    end do
    if (proc0) write (ierr,'(a)') 'g2gf scatter check done'

    deallocate (gtmp,gftmp)

!    call finish_mp
!    stop

  contains
    !> FIXME : Add documentation
    function rule (ig, isgn, ik, it, il, ie, is)
      integer, intent (in) :: ig, isgn, ik, it, il, ie, is
      integer :: rule
      rule = ig + isgn + ik + it + il + ie + is  ! make whatever you want
    end function rule
  end subroutine check_g2gf