Perform some basic checking of the diagnostic input parameters, and write the results to file
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | report_unit |
Unit of an open file to write to |
subroutine check_gs2_diagnostics(report_unit)
use file_utils, only: run_name
use nonlinear_terms, only: nonlin
use dist_fn, only : def_parity, even
use kt_grids, only : gridopt_switch, gridopt_box
use init_g, only : restart_file
use gs2_save, only: restart_writable
implicit none
!> Unit of an open file to write to
integer, intent(in) :: report_unit
logical :: writable
write (report_unit, *)
write (report_unit, fmt="('------------------------------------------------------------')")
write (report_unit, *)
write (report_unit, fmt="('Diagnostic control section.')")
if (print_line) then
write (report_unit, fmt="('print_line = T: Estimated frequencies &
& output to the screen every ',i4,' steps.')") nwrite
else
! nothing
end if
if (write_line) then
if (write_ascii) then
write (report_unit, fmt="('write_line = T: Estimated frequencies output to ',a,' every ',i4,' steps.')") &
& trim(run_name)//'.out', nwrite
end if
write (report_unit, fmt="('write_line = T: Estimated frequencies output to ',a,' every ',i4,' steps.')") &
& trim(run_name)//'.out.nc', nwrite
else
! nothing
end if
if (print_flux_line) then
write (report_unit, fmt="('print_flux_line = T: Instantaneous fluxes output to screen every ', &
& i4,' steps.')") nwrite
else
! nothing
end if
if (write_flux_line) then
if (write_ascii) then
write (report_unit, fmt="('write_flux_line = T: Instantaneous fluxes output to ',a,' every ',i4,' steps.')") &
& trim(run_name)//'.out', nwrite
end if
write (report_unit, fmt="('write_flux_line = T: Instantaneous fluxes output to ',a,' every ',i4,' steps.')") &
& trim(run_name)//'.out.nc', nwrite
else
! nothing
end if
if (write_omega) then
if (write_ascii) then
write (report_unit, fmt="('write_omega = T: Instantaneous frequencies written to ',a)") trim(run_name)//'.out'
else
write (report_unit, fmt="('write_omega = T: No effect.')")
end if
write (report_unit, fmt="(' Frequencies calculated at igomega = ',i4)") igomega
if (def_parity .and. .not. even) then
write (report_unit, fmt="('################# WARNING #######################')")
write (report_unit, fmt="(' You probably want igomega /= 0 for odd parity modes.')")
write (report_unit, fmt="('################# WARNING #######################')")
write (report_unit, *)
end if
end if
if (write_omavg) then
if (write_ascii) then
write (report_unit, fmt="('write_omavg = T: Time-averaged frequencies written to ',a)") trim(run_name)//'.out'
write (report_unit, fmt="(' Averages taken over ',i4,' timesteps.')") navg
else
write (report_unit, fmt="('write_omavg = T: No effect.')")
end if
end if
if (write_ascii) then
write (report_unit, fmt="('write_ascii = T: Write some data to ',a)") trim(run_name)//'.out'
end if
if (write_eigenfunc) then
if (write_ascii) then
write (report_unit, fmt="('write_eigenfunc = T: Normalized Phi(theta) written to ',a)") trim(run_name)//'.eigenfunc'
end if
write (report_unit, fmt="('write_eigenfunc = T: Normalized Phi(theta) written to ',a)") trim(run_name)//'.out.nc'
end if
if (write_final_fields) then
if (write_ascii) then
write (report_unit, fmt="('write_final_fields = T: Phi(theta), etc. written to ',a)") trim(run_name)//'.fields'
end if
write (report_unit, fmt="('write_final_fields = T: Phi(theta), etc. written to ',a)") trim(run_name)//'.out.nc'
end if
if (write_final_antot) then
if (write_ascii) then
write (report_unit, fmt="('write_final_antot = T: Sources for Maxwell eqns. written to ',a)") &
& trim(run_name)//'.antot'
end if
write (report_unit, fmt="('write_final_antot = T: Sources for Maxwell eqns. written to ',a)") &
& trim(run_name)//'.out.nc'
end if
if (write_final_moments) then
if (write_ascii) then
write (report_unit, fmt="('write_final_moments = T: Low-order moments of g written to ',a)") &
& trim(run_name)//'.moments'
write (report_unit, fmt="('write_final_moments = T: int dl/B average of low-order moments of g written to ',a)") &
& trim(run_name)//'.amoments'
end if
write (report_unit, fmt="('write_final_moments = T: Low-order moments of g written to ',a)") &
& trim(run_name)//'.out.nc'
write (report_unit, fmt="('write_final_moments = T: int dl/B average of low-order moments of g written to ',a)") &
& trim(run_name)//'.out.nc'
end if
if (write_avg_moments) then
if (gridopt_switch /= gridopt_box) then
write (report_unit, *)
write (report_unit, fmt="('################# WARNING #######################')")
write (report_unit, fmt="('write_avg_moments = T: Ignored unless grid_option=box')")
write (report_unit, fmt="('################# WARNING #######################')")
write (report_unit, *)
else
if (write_ascii) then
write (report_unit, fmt="('write_avg_moments = T: Flux surface averaged low-order moments of g written to ',a)") &
& trim(run_name)//'.moments'
end if
write (report_unit, fmt="('write_avg_moments = T: Flux surface averaged low-order moments of g written to ',a)") &
& trim(run_name)//'.out.nc'
end if
end if
if (write_final_epar) then
if (write_ascii) then
write (report_unit, fmt="('write_final_epar = T: E_parallel(theta) written to ',a)") trim(run_name)//'.epar'
end if
write (report_unit, fmt="('write_final_epar = T: E_parallel(theta) written to ',a)") trim(run_name)//'.out.nc'
end if
if (write_fluxes) then
if (write_ascii) then
write (report_unit, fmt="('write_fluxes = T: Phi**2(kx, ky) written to ',a)") trim(run_name)//'.out'
end if
else
write (report_unit, fmt="('write_fluxes = F: Phi**2(kx, ky) NOT written to ',a)") trim(run_name)//'.out'
end if
if (dump_check1) then
write (report_unit, fmt="('dump_check1 = T: Field-line avg of Phi written to ',a)") 'dump.check1'
write (report_unit, fmt="('This option is usually used for Rosenbluth-Hinton calculations.')")
write (report_unit, fmt="('THIS IS PROBABLY AN ERROR.')")
end if
if (dump_check2) then
write (report_unit, fmt="('dump_check2 = T: Apar(kx, ky, igomega) written to ',a)") trim(run_name)//'.dc2'
write (report_unit, fmt="('THIS IS PROBABLY AN ERROR.')")
end if
if (dump_fields_periodically) then
write (report_unit, fmt="('dump_fields_periodically = T: Phi, Apar, Bpar written to ',a)") 'dump.fields.t=(time)'
write (report_unit, fmt="('THIS IS PROBABLY AN ERROR. IT IS EXPENSIVE.')")
end if
if (save_for_restart) then
write (report_unit, fmt="('save_for_restart = T: Restart files written to ',a)") trim(restart_file)//'.(PE)'
else
if (nonlin) then
write (report_unit, *)
write (report_unit, fmt="('################# WARNING #######################')")
write (report_unit, fmt="('save_for_restart = F: This run cannot be continued.')")
write (report_unit, fmt="('THIS IS PROBABLY AN ERROR.')")
write (report_unit, fmt="('################# WARNING #######################')")
write (report_unit, *)
end if
end if
!Verify restart file can be written
if((save_for_restart.or.save_distfn).and.(file_safety_check))then
!Can we write file?
writable=restart_writable()
!If we can't write the restart file then we should probably quit
if((.not.writable))then
if(save_for_restart)then
write (report_unit, *)
write (report_unit, fmt="('################# WARNING #######################')")
write (report_unit, fmt="('save_for_restart = T: But we cannot write to a test file like ',A,'.')") trim(restart_file)
write (report_unit, fmt="('THIS IS PROBABLY AN ERROR --> Check restart_dir.')")
write (report_unit, fmt="('################# WARNING #######################')")
write (report_unit, *)
endif
if(save_distfn)then
write (report_unit, *)
write (report_unit, fmt="('################# WARNING #######################')")
write (report_unit, fmt="('save_distfn = T: But we cannot write to a test file like ',A,'.')") trim(restart_file)
write (report_unit, fmt="('THIS IS PROBABLY AN ERROR --> Check restart_dir.')")
write (report_unit, fmt="('################# WARNING #######################')")
write (report_unit, *)
endif
endif
endif
end subroutine check_gs2_diagnostics