FIXME : Add documentation
This routine depends on the values stored in xx as input
to the legendre_polynomials method. However, with Radau grids the
values in xx are not Legendre roots so this routine is probably
incorrect when [[le_grids_knobs::radau_gauss_grid]] is .true.
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
complex, | intent(in), | dimension (-ntgrid:,:,g_lo%llim_proc:) | :: | g | ||
complex, | intent(out), | dimension (0:,-ntgrid:,:,:,:) | :: | tote | ||
complex, | intent(out), | dimension (0:,-ntgrid:,:,:,:) | :: | totl | ||
complex, | intent(out), | optional, | dimension (0:,-ntgrid:,:,:,:) | :: | tott |
subroutine legendre_transform (g, tote, totl, tott)
use mp, only: nproc
use theta_grid, only: ntgrid, bmag, bmax
use species, only: nspec
use kt_grids, only: naky, ntheta0
use gs2_layouts, only: g_lo, idx, idx_local
use mp, only: sum_reduce
implicit none
complex, dimension (-ntgrid:,:,g_lo%llim_proc:), intent (in) :: g
complex, dimension (0:,-ntgrid:,:,:,:), intent (out) :: tote, totl
complex, dimension (0:,-ntgrid:,:,:,:), intent (out), optional :: tott
complex :: totfac
real :: ulim
integer :: is, il, ie, ik, it, iglo, ig, im, ntrap
integer, save :: lpesize
real, dimension (:), allocatable :: nodes
real, dimension (:,:), allocatable :: lpltmp, lpttmp
! real, dimension (:,:), allocatable, save :: lpe, lpl
! real, dimension (:,:,:,:), allocatable, save :: lpt
if (.not. allocated(lpl)) then
allocate(lpltmp(ng2,0:ng2-1))
allocate(lpl(nlambda,0:ng2-1))
lpesize = nesub
allocate(lpe(negrid,0:lpesize-1,nspec)) ; lpe = 0.0
! get value of first nesub legendre polynomials
! at each of the grid points on (0,vcut)
do is = 1,nspec
call legendre_polynomials (0.0,vcut,speed(:lpesize,is),lpe(:lpesize,:,is))
end do
! TEMP FOR TESTING -- MAB
! lpe = 2.*lpe/vcut
! get value of first ng2 legendre polynomials
! at each of the grid points on (0,1)
call legendre_polynomials (0.0,1.0,xx,lpltmp)
lpl = 0.0
lpl(1:ng2,:) = lpltmp
if (present(tott)) then
allocate (lpt(nlambda,0:2*(nlambda-ng2-1),-ntgrid:ntgrid,2))
lpt = 0.0
do ig = -ntgrid, ntgrid
ntrap = 1
if (jend(ig) > ng2+1) then
ntrap = jend(ig)-ng2
allocate (nodes(2*ntrap-1))
allocate (lpttmp(2*ntrap-1,0:2*(ntrap-1)))
do il = 1, ntrap
nodes(il) = -sqrt(max(0.0,1.0-al(ng2+il)*bmag(ig)))
end do
nodes(ntrap+1:) = -nodes(ntrap-1:1:-1)
!Can we remove this?
! TEMP FOR TESTING -- MAB
! nodes = nodes + sqrt(1.0-bmag(ig)/bmax)
! ulim = 2.*sqrt(1.0-bmag(ig)/bmax)
ulim = sqrt(1.0-bmag(ig)/bmax)
call legendre_polynomials (-ulim,ulim,nodes,lpttmp)
lpt(ng2+1:jend(ig),0:2*(ntrap-1),ig,2) = lpttmp(1:ntrap,:)
lpt(ng2+1:jend(ig)-1,0:2*(ntrap-1),ig,1) = lpttmp(2*ntrap-1:ntrap+1:-1,:)
!Can we remove this?
! lpt(ng2+1:jend(ig),0:2*(ntrap-1),ig,1) = lpttmp(2*ntrap-1:ntrap+1:-1,:)
! do ie = 0, 2*(ntrap-1)
! do il = 1, 2*ntrap-1
! if (proc0) write (*,*) 'lptrap', ig, ntrap, ulim, ie, il, nodes(il), lpttmp(il,ie)
! end do
! end do
! do ie = 0, 2*(ntrap-1)
! do il = ng2+1,jend(ig)
! write (*,*) 'lpt', ig, ie, il, lpt(il,ie,ig,1), lpt(il,ie,ig,2)
! end do
! end do
deallocate (nodes, lpttmp)
end if
end do
end if
deallocate (lpltmp)
end if
! carry out legendre transform to get coefficients of
! legendre polynomial expansion of g
totfac = 0. ; tote = 0. ; totl = 0.
if (present(tott)) tott = 0.
!Loop over all indices, note this loop is optimal only for layout 'xyles' (at least in terms of
!g memory access)
do is = 1, nspec
do ie = 1, negrid
do il = 1, nlambda
do ik = 1, naky !Swapped ik and it loop order.
do it = 1, ntheta0
iglo = idx (g_lo, ik, it, il, ie, is)
if (idx_local (g_lo, iglo)) then
do ig=-ntgrid,ntgrid
totfac = w(ie,is)*wl(ig,il)*(g(ig,1,iglo)+g(ig,2,iglo))
do im=0,lpesize-1
tote(im, ig, it, ik, is) = tote(im, ig, it, ik, is) + totfac*lpe(ie,im,is)*(2*im+1)
end do
do im=0,ng2-1
totl(im, ig, it, ik, is) = totl(im, ig, it, ik, is) + totfac*lpl(il,im)*(2*im+1)
end do
if (present(tott)) then
do im=0,2*(jend(ig)-ng2-1)
tott(im, ig, it, ik, is) = tott(im, ig, it, ik, is) + &
w(ie,is)*wl(ig,il)*(lpt(il,im,ig,1)*g(ig,1,iglo)+lpt(il,im,ig,2)*g(ig,2,iglo))*(2*im+1)
end do
end if
end do
end if
end do
end do
end do
end do
end do
!Do we really need this if?
if (nproc > 1) then
!Now complete velocity integral, bringing back results to proc0
call sum_reduce (tote, 0)
call sum_reduce (totl, 0)
if (present(tott)) call sum_reduce (tott, 0)
end if
end subroutine legendre_transform