Output the poloidal distributions of the fluxes of particles, perpendicular momentum, parallel momentum, and energy ! JRB
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | file_id |
NetCDF ID of the file to write to |
||
integer, | intent(in) | :: | nout |
Current timestep |
||
real, | intent(in), | dimension (-ntgrid:,:) | :: | part_fluxes_dist | ||
real, | intent(in), | dimension (-ntgrid:,:) | :: | mom_fluxes_par_dist | ||
real, | intent(in), | dimension (-ntgrid:,:) | :: | mom_fluxes_perp_dist | ||
real, | intent(in), | dimension (-ntgrid:,:) | :: | heat_fluxes_dist | ||
real, | intent(in), | dimension (2,-ntgrid:ntgrid) | :: | phi_dist |
subroutine nc_loop_dist (file_id, nout, part_fluxes_dist, mom_fluxes_par_dist, mom_fluxes_perp_dist, heat_fluxes_dist, phi_dist)
use run_parameters, only: has_phi
use theta_grid, only: ntgrid
# ifdef NETCDF
use neasyf, only: neasyf_write
# endif
!> NetCDF ID of the file to write to
integer, intent(in) :: file_id
!> Current timestep
integer, intent (in) :: nout
real, dimension (-ntgrid:,:), intent (in) :: part_fluxes_dist, mom_fluxes_par_dist, mom_fluxes_perp_dist, heat_fluxes_dist
real, dimension (2,-ntgrid:ntgrid), intent (in) :: phi_dist
# ifdef NETCDF
character(*), dimension(*), parameter :: flux_dist_dim = [theta_dim, species_dim, time_dim]
character(*), dimension(*), parameter :: phi_dist_dim = [complex_dim, theta_dim, time_dim]
if (has_phi) then
call neasyf_write(file_id, "es_part_flux_dist", part_fluxes_dist, &
dim_names=flux_dist_dim, start=starts(3, nout), &
long_name= "Particle flux density, see Ball et al PPCF 58 045023 2016")
call neasyf_write(file_id, "es_mom_flux_par_dist", mom_fluxes_par_dist, &
dim_names=flux_dist_dim, start=starts(3, nout), &
long_name= "Parallel momentum flux density, see Ball et al PPCF 58 045023 2016")
call neasyf_write(file_id, "es_mom_flux_perp_dist", mom_fluxes_perp_dist, &
dim_names=flux_dist_dim, start=starts(3, nout), &
long_name= "Perpendicular momentum flux density, see Ball et al PPCF 58 045023 2016")
call neasyf_write(file_id, "es_heat_flux_dist", heat_fluxes_dist, &
dim_names=flux_dist_dim, start=starts(3, nout), &
long_name= "Heat flux density, see Ball et al PPCF 58 045023 2016")
call neasyf_write(file_id, "es_phi_dist", phi_dist, &
dim_names=phi_dist_dim, start=starts(3, nout))
end if
# endif
end subroutine nc_loop_dist