Calculate the poloidal distributions of the fluxes of particles, parallel momentum, perpendicular momentum, and energy (see section 3.1 and appendix A of Ball et al. PPCF 58 (2016) 045023 as well as section 5 of "GS2 analytic geometry specification")
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | file_id |
NetCDF ID of the file to write to |
||
integer, | intent(in) | :: | nout |
Current timestep |
subroutine do_write_nl_flux_dist(file_id, nout)
use diagnostics_fluxes, only: flux_dist
use theta_grid, only: ntgrid
use species, only: nspec, spec
use mp, only: proc0
use gs2_io, only: nc_loop_dist
use convert, only: c2r
use kt_grids, only: naky, ntheta0
use fields_arrays, only: phinew, bparnew
use run_parameters, only: has_phi
implicit none
!> NetCDF ID of the file to write to
integer, intent(in) :: file_id
!> Current timestep
integer, intent(in) :: nout
real, dimension (:,:), allocatable :: part_fluxes_dist, mom_fluxes_par_dist, mom_fluxes_perp_dist, heat_fluxes_dist
real, dimension (:,:), allocatable :: phi_dist
real, dimension (:,:,:,:), allocatable :: pflux_dist, vflux_par_dist, vflux_perp_dist, qflux_dist
real, dimension (:,:,:,:), allocatable :: phi_byMode_dist
integer :: ig, is
allocate (part_fluxes_dist(-ntgrid:ntgrid,nspec))
allocate (mom_fluxes_par_dist(-ntgrid:ntgrid,nspec))
allocate (mom_fluxes_perp_dist(-ntgrid:ntgrid,nspec))
allocate (heat_fluxes_dist(-ntgrid:ntgrid,nspec))
allocate (phi_dist(2,-ntgrid:ntgrid))
allocate (pflux_dist(-ntgrid:ntgrid,ntheta0,naky,nspec))
allocate (vflux_par_dist(-ntgrid:ntgrid,ntheta0,naky,nspec))
allocate (vflux_perp_dist(-ntgrid:ntgrid,ntheta0,naky,nspec))
allocate (qflux_dist(-ntgrid:ntgrid,ntheta0,naky,nspec))
allocate (phi_byMode_dist(2,-ntgrid:ntgrid,ntheta0,naky))
call c2r(phinew,phi_byMode_dist)
call flux_dist (phinew, bparnew, pflux_dist, vflux_par_dist, vflux_perp_dist, qflux_dist)
if (proc0) then
if (has_phi) then
do ig = -ntgrid, ntgrid
do is = 1, nspec
pflux_dist(ig,:,:,is) = pflux_dist(ig,:,:,is) * spec(is)%dens
call get_volume_average (pflux_dist(ig,:,:,is), part_fluxes_dist(ig,is))
vflux_par_dist(ig,:,:,is) = vflux_par_dist(ig,:,:,is) * spec(is)%dens*sqrt(spec(is)%mass*spec(is)%temp)
call get_volume_average (vflux_par_dist(ig,:,:,is), mom_fluxes_par_dist(ig,is))
vflux_perp_dist(ig,:,:,is) = vflux_perp_dist(ig,:,:,is) * spec(is)%dens*sqrt(spec(is)%mass*spec(is)%temp)
call get_volume_average (vflux_perp_dist(ig,:,:,is), mom_fluxes_perp_dist(ig,is))
qflux_dist(ig,:,:,is) = qflux_dist(ig,:,:,is) * spec(is)%dens*spec(is)%temp
call get_volume_average (qflux_dist(ig,:,:,is), heat_fluxes_dist(ig,is))
end do
call get_volume_average (phi_byMode_dist(1,ig,:,:), phi_dist(1,ig))
call get_volume_average (phi_byMode_dist(2,ig,:,:), phi_dist(2,ig))
end do
end if
end if
if (proc0) call nc_loop_dist (file_id, nout, part_fluxes_dist, mom_fluxes_par_dist, &
mom_fluxes_perp_dist, heat_fluxes_dist, phi_dist)
deallocate (part_fluxes_dist)
deallocate (mom_fluxes_par_dist)
deallocate (mom_fluxes_perp_dist)
deallocate (heat_fluxes_dist)
deallocate (phi_dist)
deallocate (pflux_dist)
deallocate (vflux_par_dist)
deallocate (vflux_perp_dist)
deallocate (qflux_dist)
deallocate (phi_byMode_dist)
end subroutine do_write_nl_flux_dist