Write normalised to the value of at the outboard midplane
Writes to text file <runname>.eigenfunc
and/or netCDF
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | file_id |
NetCDF ID of the file to write to |
||
logical, | intent(in) | :: | write_text |
If true, write text file |
||
complex, | intent(inout), | dimension (:,:) | :: | phi0 |
The normalising factor for the fields that is actually used. See write_eigenfunc for more details |
subroutine do_write_eigenfunc(file_id, write_text, phi0)
use mp, only: proc0
use file_utils, only: open_output_file, close_output_file
use fields_arrays, only: phi, apar, bpar
use kt_grids, only: ntheta0, naky, theta0, aky
use theta_grid, only: theta, ntgrid
use gs2_io, only: nc_eigenfunc
use dist_fn, only: def_parity, even
use run_parameters, only: has_apar
implicit none
!> NetCDF ID of the file to write to
integer, intent(in) :: file_id
!> If true, write text file
logical, intent(in) :: write_text
!> The normalising factor for the fields that is actually used. See
!> [[gs2_diagnostics_knobs:write_eigenfunc]] for more details
complex, dimension (:,:), intent(inout) :: phi0
integer :: it, ik, ig, unit
if(.not.proc0) return
!Should this actually use igomega instead of 0?
!What if fphi==0? --> See where statements below
phi0 = phi(0,:,:)
!This looks like a hack for the case where we know we've forced phi(theta=0) to be 0
!this could probably be better addressed by the use of igomega above
if (def_parity .and. has_apar .and. (.not. even)) phi0 = apar(0, :, :)
!Address locations where phi0=0 by using next point
where (abs(phi0) < 10.0*epsilon(0.0))
phi0 = phi(1,:,:)/(theta(1)-theta(0))
end where
!Check again if any locations are 0, this could be true if fphi (or fapar)
!is zero.
where (abs(phi0) < 10.0*epsilon(0.0))
phi0 = 1.0
end where
if (write_text) then
call open_output_file (unit, ".eigenfunc")
do ik = 1, naky
do it = 1, ntheta0
do ig = -ntgrid, ntgrid
write (unit, "(9(1x,e12.5))") &
theta(ig), theta0(it,ik), aky(ik), &
phi(ig,it,ik)/phi0(it,ik), &
apar(ig,it,ik)/phi0(it,ik), &
bpar(ig,it,ik)/phi0(it,ik)
end do
write (unit, "()")
end do
end do
call close_output_file (unit)
end if
call nc_eigenfunc (file_id, phi0)
end subroutine do_write_eigenfunc