Integrate species using gf_lo data format
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
complex, | intent(in), | dimension (-ntgrid:,:,g_lo%llim_proc:) | :: | g | ||
real, | intent(inout), | dimension (:) | :: | weights | ||
complex, | intent(out), | dimension (-ntgrid:ntgrid,gf_lo%ntheta0,gf_lo%naky) | :: | total |
subroutine integrate_species_gf (g, weights, total)
use species, only : nspec
use theta_grid, only: ntgrid
use gs2_layouts, only: gf_lo, g_lo
use gs2_layouts, only: is_idx, ik_idx, it_idx, ie_idx, il_idx
use mp, only: sum_allreduce_sub, sum_allreduce
use kt_grids, only: kwork_filter
use redistribute, only: gather
use species, only: spec, tracer_species
use array_utils, only: zero_array
implicit none
complex, dimension (-ntgrid:,:,g_lo%llim_proc:), intent (in) :: g
complex, dimension (:, :, :, :, :, :), allocatable :: gf
real, dimension (:), intent (in out) :: weights
complex, dimension (-ntgrid:ntgrid,gf_lo%ntheta0,gf_lo%naky), intent (out) :: total
integer :: is, il, ie, igf, it, ik
call zero_array(total)
where (spec%type == tracer_species) weights = 0
allocate(gf(-ntgrid:ntgrid,2,nspec,negrid,nlambda,gf_lo%llim_proc:gf_lo%ulim_alloc))
call gather(g2gf, g, gf, ntgrid)
!Performed integral (weighted sum) over local velocity space and species
if(any(kwork_filter)) then
do igf = gf_lo%llim_proc,gf_lo%ulim_proc
it = it_idx(gf_lo,igf)
ik = ik_idx(gf_lo,igf)
if(kwork_filter(it,ik)) cycle
do il = 1,gf_lo%nlambda
do ie = 1,gf_lo%negrid
do is = 1,gf_lo%nspec
total(:,it,ik) = total(:,it,ik) + &
(weights(is)*w(ie,is))*wl(:,il)*(gf(:,1,is,ie,il,igf)+gf(:,2,is,ie,il,igf))
end do
end do
end do
end do
else
do igf = gf_lo%llim_proc,gf_lo%ulim_proc
it = it_idx(gf_lo,igf)
ik = ik_idx(gf_lo,igf)
do il = 1,gf_lo%nlambda
do ie = 1,gf_lo%negrid
do is = 1,gf_lo%nspec
total(:,it,ik) = total(:,it,ik) + &
(weights(is)*w(ie,is))*wl(:,il)*(gf(:,1,is,ie,il,igf)+gf(:,2,is,ie,il,igf))
end do
end do
end do
end do
end if
deallocate(gf)
end subroutine integrate_species_gf