program main #ifdef INTEL use ifport, only : chdir, getcwd #endif use const_and_precisions, only : wp_ use logger, only : INFO, ERROR, WARNING, set_logger, log_message use utils, only : dirname use types, only : contour use gray_cli, only : cli_options, parse_cli_options, & parse_param_overrides use gray_core, only : gray_main use gray_errors, only : gray_error use gray_equil, only : abstract_equil, load_equil use gray_plasma, only : abstract_plasma, load_plasma use gray_params, only : gray_parameters, gray_results, & read_gray_params, read_gray_config, & make_gray_header implicit none integer :: status ! exit status no_save: block ! CLI options type(cli_options) :: opts ! gray_main subroutine arguments type(gray_parameters) :: params ! Inputs class(abstract_equil), allocatable :: equil ! class(abstract_plasma), allocatable :: plasma ! type(contour) :: limiter ! type(gray_results) :: results ! Outputs integer(kind=gray_error) :: sim_error ! ! Store the original working directory integer :: err character(len=256) :: cwd err = getcwd(cwd) ! Parse the command-line options call parse_cli_options(opts) ! Initialise logging if (opts%quiet) opts%verbose = ERROR call set_logger(level=opts%verbose) ! Load the parameters from file and move to its directory ! (all other filepaths are assumed relative to it) if (allocated(opts%config_file)) then ! Use the newer gray.ini configuration file call log_message(level=INFO, mod='main', msg='using GRAY INI config') call read_gray_config(opts%config_file, params, err) if (err /= 0) call exit(1) err = chdir(dirname(opts%config_file)) else ! Use the legacy gray_params.data file call log_message(level=INFO, mod='main', msg='using GRAY params file') call read_gray_params(opts%params_file, params, err) if (err /= 0) call exit(1) err = chdir(dirname(opts%params_file)) end if if (err /= 0) then call log_message(level=ERROR, mod='main', & msg='chdir to configuration file directory failed!') call exit(1) end if ! Initialise antenna parameters before parsing the ! cli arguments, so antenna%pos can be overriden call init_antenna(params%antenna, err) if (err /= 0) call exit(1) ! Apply CLI overrides to the parameters call parse_param_overrides(params) ! Do some checks on the inputs associate (p => params%raytracing) if (p%igrad == 1 .and. p%nrayr < 5) then p%igrad = 0 call log_message(level=WARNING, mod='main', & msg='igrad = 1 but nrayr < 5: disabling beamtracing') end if if (p%nrayr == 1) p%nrayth = 1 if (p%nrayr > 1) then p%rwmax = real(p%nrayr - 1) / max(nint((p%nrayr - 1)/p%rwmax), 1) end if end associate ! Read and initialise the equilibrium and limiter objects call load_equil(params%equilibrium, equil, limiter, err) if (err /= 0) call exit(1) ! Read and initialise the plasma state object call load_plasma(params, equil, plasma, err) if (err /= 0) call exit(1) ! Create a simple limiter, if necessary if (.not. allocated(limiter%R) .or. params%raytracing%ipass < 0) then call log_message('using fallback limiter', level=INFO, mod='main') params%raytracing%ipass = abs(params%raytracing%ipass) limiter = make_fallback_limiter(params, equil) end if ! Get back to the initial directory err = chdir(cwd) if (allocated(opts%sum_filelist)) then call log_message(level=INFO, mod='main', msg='summing profiles') call sum_profiles(params, equil, opts%sum_filelist, opts%output_dir, results) status = 0 else ! Activate the given output tables params%misc%active_tables = opts%tables ! Finally run the simulation call gray_main(params, equil, plasma, limiter, results, sim_error) status = 2 * sim_error end if print_res: block character(256) :: msg write(msg, '(a,g0.3," MW")') 'total absoption: P=', results%pabs call log_message(msg, level=INFO, mod='main') write(msg, '(a,g0.3," kA")') 'total current drive: I=', results%icd * 1.0e3_wp_ call log_message(msg, level=INFO, mod='main') end block print_res ! Change the current directory to output files there if (allocated(opts%output_dir)) then if (chdir(opts%output_dir) /= 0) then call log_message(level=ERROR, mod='main', & msg='chdir to output_dir ('//opts%output_dir//') failed!') call exit(1) end if end if ! Write results to file write_res: block integer :: i do i = 1, size(results%tables%iter) if (.not. associated(results%tables%iter(i)%ptr)) cycle associate (tbl => results%tables%iter(i)%ptr) call tbl%save(err, header=make_gray_header(params)) if (err /= 0) then call log_message('failed to save '//trim(tbl%title)//' table', & level=ERROR, mod='main') status = 1 end if end associate end do end block write_res end block no_save ! Exit status is >1 if any error occurred call exit(status) contains subroutine init_antenna(params, err) ! Reads the wave launcher file (containing the wave frequency, launcher ! position, direction and beam description) and initialises the respective ! GRAY parameters. use beams, only : read_beam0, read_beam1, read_beam2 use gray_params, only : antenna_parameters, BEAM_0D, BEAM_1D, BEAM_2D ! subroutine arguments type(antenna_parameters), intent(inout) :: params integer, intent(out) :: err ! Note: α, β are loaded from gray_params.data select case (params%ibeam) case (BEAM_2D) ! 2 degrees of freedom ! w(z, α, β), 1/R(z, α, β) ! FIXME: 1st beam is always selected, iox read from table call read_beam2(params, beamid=1, err=err, set_pol=.false.) case (BEAM_1D) ! 1 degree of freedom ! w(z, α), 1/R(z, α) call read_beam1(params, err) case (BEAM_0D) ! fixed w(z), 1/R(z) call read_beam0(params, err) end select end subroutine init_antenna function make_fallback_limiter(params, equil) result(limiter) ! Creates a simple rectangular limiter: ! ! (R_wall, z_launch+Δz)╔═════╗(R_max, z_launch+Δz) ! ║4 3║ ! ║ ║ ! ║1 2║ ! (R_wall, z_launch-Δz)╚═════╝(R_max, z_launch-Δz) ! ! Note: R_max = { +∞ for vacuum, ! { R₀+a for analytical model, ! { grid_r[-1] for numerical data. ! use logger, only : log_info use const_and_precisions, only : cm ! function arguments type(gray_parameters), intent(in) :: params class(abstract_equil), intent(in) :: equil type(contour) :: limiter ! local variables real(wp_) :: R_launch, z_launch, R_max, delta_z ! Launcher coordinates R_launch = norm2(params%antenna%pos(1:2)) * cm z_launch = params%antenna%pos(3) * cm ! Max vertical distance a ray can travel delta_z = abs(params%raytracing%ipass) * & params%raytracing%dst * params%raytracing%nstep * cm ! Avoid clipping out the launcher R_max = max(R_launch, equil%r_range(2)) ! Build a rectangle limiter = contour( & Rmin=params%misc%rwall, Rmax=R_max, & zmin=z_launch - delta_z, zmax=z_launch + delta_z) end function make_fallback_limiter subroutine sum_profiles(params, equil, filelist, output_dir, results) ! Combines the EC profiles obtained from multiple gray runs ! (of different beams) and recomputes the summary table use gray_params, only : gray_parameters use gray_equil, only : abstract_equil use pec, only : pec_init, postproc_profiles, dealloc_pec, & rhot_tab ! subroutine arguments type(gray_parameters), intent(inout) :: params class(abstract_equil), intent(in) :: equil character(len=*), intent(in) :: filelist, output_dir type(gray_results), intent(inout) :: results ! local variables real(wp_), dimension(:), allocatable :: jphi real(wp_), dimension(:), allocatable :: currins, pins, rtin, rpin real(wp_), dimension(:), allocatable :: extracol integer :: i, j, k, l, nbeams, ncases, nextracols, irt, err character(len=255) :: filename, line, extralabels, fmtstr real(wp_), dimension(5) :: sums real(wp_) :: jphip, dpdvp, rhotj, rhotjava, rhotp, rhotpav, & drhotjava, drhotpav, ratjamx, ratjbmx real(wp_) :: drhotp, drhotj, dpdvmx, jphimx ! unit numbers integer :: list, prof, summ integer, allocatable :: beams(:) ! Initialise the output profiles call pec_init(params%output, equil) associate(nrho =>params%output%nrho) allocate(jphi(nrho), currins(nrho), pins(nrho), rtin(nrho), rpin(nrho)) allocate(results%jcd(nrho), results%dpdv(nrho)) end associate ! Read the list of files to be summed open(newunit=list, file=filelist, action='read', status='old', iostat=err) if (err /= 0) then call log_message(level=ERROR, mod='main', & msg='opening file list ('//filelist//') failed!') call exit(1) end if ! Open combined output files open(newunit=prof, file=output_dir//'/'//'sum-ec-profiles.txt', action='write') open(newunit=summ, file=output_dir//'/'//'sum-summary.txt', action='write') ! Move to the filelist directory ! (the filepaths are assumed relative to it) if (chdir(dirname(filelist)) /= 0) then call log_message(level=ERROR, mod='main', & msg='chdir to file list directory failed!') call exit(1) end if ! Read the filelist header read(list, *) nbeams, ncases, nextracols allocate(beams(nbeams), extracol(nextracols)) ! Read the list of files do i = 1, nbeams read(list, *) filename ! Open the summand file open(newunit=beams(i), file=filename, & action='read', status='old', iostat=err) if (err /= 0) then call log_message(level=ERROR, mod='main', & msg='opening summand file ('//trim(filename)//') failed!') call exit(1) end if ! Copy the Gray header from the ec-profiles file do j = 1, 21 read(beams(i), '(a255)') line if (i == 1) write(prof, '(a)') trim(line) if (i == 1) write(summ, '(a)') trim(line) end do ! Copy the column labels ! Note: for the summary.7 we copy only the first columns read(beams(i), '(a255)') line if (i == 1) then extralabels = line(1:index(line, "ρ_p")-1) write(prof, '(a)') trim(line) write(summ, '(a,x,a,x,a)') trim(extralabels), & "I_cd P_abs J_φ_peak dPdV_peak ρ_max_J ρ_avg_J ρ_max_P", & "ρ_avg_P Δρ_avg_J Δρ_avg_P ratio_Ja_max ratio_Jb_max" end if end do close(list) do k = 1, ncases jphi = 0 currins = 0 pins = 0 results%jcd = 0 results%dpdv = 0 ! Sum the radial profiles do j = 1, params%output%nrho do i = 1, nbeams read(beams(i), *) & (extracol(l), l=1,nextracols), rpin(j), rtin(j), sums(1:5), irt jphi(j) = sums(1) + jphi(j) results%jcd(j) = sums(2) + results%jcd(j) results%dpdv(j) = sums(3) + results%dpdv(j) currins(j) = sums(4) + currins(j) pins(j) = sums(5) + pins(j) end do ! Store the sums write(fmtstr, '("(",i0,"(1x,e16.8e3),i5)")') nextracols + 7 write(prof, trim(fmtstr)) & (extracol(l), l=1,nextracols), rpin(j), rtin(j), jphi(j), & results%jcd(j), results%dpdv(j), currins(j), pins(j), irt end do ! Combine the total power, current results%pabs = pins(params%output%nrho) results%icd = currins(params%output%nrho) ! Recompute the summary values call postproc_profiles(equil, & results%Pabs, results%Icd, rhot_tab, results%dPdV, jphi, & rhotpav, drhotpav, rhotjava, drhotjava, dpdvp, jphip, & rhotp, drhotp, rhotj, drhotj, dpdvmx, jphimx, ratjamx, ratjbmx) write(fmtstr, '("(",i0,"(1x,e12.5),i5,7(1x,e12.5))")') nextracols+12 write(summ, trim(fmtstr)) & (extracol(l), l=1,nextracols), results%icd, results%pabs, & jphip, dpdvp, rhotj, rhotjava, rhotp, rhotpav, drhotjava, & drhotpav, ratjamx, ratjbmx ! Separate the case write(prof, *) end do ! Close all files close(prof); close(summ) do i = 1, nbeams close(beams(i)) end do ! Free memory call dealloc_pec end subroutine sum_profiles end program main