Routine to select first scan case author: J Morris, UKAEA, Culham Science Centre
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | nwp | |||
real(kind=dp), | intent(in), | dimension(:) | :: | swp | ||
integer, | intent(in) | :: | iscn | |||
character(len=25), | intent(out) | :: | vlab | |||
character(len=25), | intent(out) | :: | xlab |
subroutine scan_select(nwp, swp, iscn, vlab, xlab)
!! Routine to select first scan case
!! author: J Morris, UKAEA, Culham Science Centre
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
use build_variables, only: blnkoth, shldith, scrapli, scraplo, ohcth
use constraint_variables, only: fiooic, walalw, bmxlim, fqval, taulimit, &
gammax, t_burn_min, tbrmin, fjprot, pnetelin, powfmax
use cost_variables, only: cfactr, iavail, fkind, startupratio
use current_drive_variables, only: bootstrap_current_fraction_max, etaech
use divertor_variables, only: hldivlim
use error_handling, only: idiags, report_error
use fwbs_variables, only: inlet_temp_liq, outlet_temp_liq, blpressure_liq, &
n_liq_recirc, bz_channel_conduct_liq, pnuc_fw_ratio_dcll, f_nuc_pow_bz_struct, pitch
use impurity_radiation_module, only: fimp, coreradius, impurity_arr_frac
use physics_variables, only: kappa, dnbeta, te, aspect, ftar, bt, &
rad_fraction_sol, triang, rmajor, beamfus0, hfact
use numerics, only: epsvmc, boundu, boundl
use tfcoil_variables, only: tmargmin_tf, sig_tf_case_max, n_pancake, oacdcp, &
n_layer, b_crit_upper_nbti, sig_tf_wp_max, fcoolcp, n_tf_turn
use heat_transport_variables, only: crypmw_max, etath
use rebco_variables, only: copperaoh_m2_max
use pfcoil_variables, only: coheof, ohhghf, oh_steel_frac
use CS_fatigue_variables, only: n_cycle_min, t_crack_vertical
implicit none
! Arguments
integer, intent(in) :: nwp, iscn
real(dp), intent(in), dimension(:) :: swp
character(len=25), intent(out) :: vlab, xlab
select case (nwp)
! Use underscores instead of spaces in xlabel
! MDK Remove the "=" from vlabel, to make it easier to compare with
! list of iteration variables
case (1)
aspect = swp(iscn)
vlab = 'aspect' ; xlab = 'Aspect_ratio'
case (2)
hldivlim = swp(iscn)
vlab = 'hldivlim' ; xlab = 'Div_heat_limit_(MW/m2)'
case (3)
pnetelin = swp(iscn)
vlab = 'pnetelin' ; xlab = 'Net_electric_power_(MW)'
case (4)
hfact = swp(iscn)
vlab = 'hfact' ; xlab = 'Confinement_H_factor'
case (5)
oacdcp = swp(iscn)
vlab = 'oacdcp' ; xlab = 'TF_inboard_leg_J_(MA/m2)'
case (6)
walalw = swp(iscn)
vlab = 'walalw' ; xlab = 'Allow._wall_load_(MW/m2)'
case (7)
beamfus0 = swp(iscn)
vlab = 'beamfus0' ; xlab = 'Beam_bkgrd_multiplier'
case (8)
fqval = swp(iscn)
vlab = 'fqval' ; xlab = 'Big_Q_f-value'
case (9)
te = swp(iscn)
vlab = 'te' ; xlab = 'Electron_temperature_keV'
case (10)
boundu(15) = swp(iscn)
vlab = 'boundu(15)' ; xlab = 'Volt-second_upper_bound'
case (11)
dnbeta = swp(iscn)
vlab = 'dnbeta' ; xlab = 'Beta_coefficient'
case (12)
bootstrap_current_fraction_max = swp(iscn)
vlab = 'bootstrap_current_fraction_max' ; xlab = 'Bootstrap_fraction'
case (13)
boundu(10) = swp(iscn)
vlab = 'boundu(10)' ; xlab = 'H_factor_upper_bound'
case (14)
fiooic = swp(iscn)
vlab = 'fiooic' ; xlab = 'TFC_Iop_/_Icrit_f-value'
case (15)
fjprot = swp(iscn)
vlab = 'fjprot' ; xlab = 'TFC_Jprot_limit_f-value'
case (16)
rmajor = swp(iscn)
vlab = 'rmajor' ; xlab = 'Plasma_major_radius_(m)'
case (17)
bmxlim = swp(iscn)
vlab = 'bmxlim' ; xlab = 'Max_toroidal_field_(T)'
case (18)
gammax = swp(iscn)
vlab = 'gammax' ; xlab = 'Maximum_CD_gamma'
case (19)
boundl(16) = swp(iscn)
vlab = 'boundl(16)' ; xlab = 'CS_thickness_lower_bound'
case (20)
t_burn_min = swp(iscn)
vlab = 't_burn_min' ; xlab = 'Minimum_burn_time_(s)'
case (21)
! sigpfalw = swp(iscn)
vlab = 'obsolete' ; xlab = 'obsolete'
case (22)
if (iavail == 1) call report_error(95)
cfactr = swp(iscn)
vlab = 'cfactr' ; xlab = 'Plant_availability_factor'
case (23)
boundu(72) = swp(iscn)
vlab = 'boundu(72)' ; xlab = 'Ip/Irod_upper_bound'
case (24)
powfmax = swp(iscn)
vlab = 'powfmax' ; xlab = 'Fusion_power_limit_(MW)'
case (25)
kappa = swp(iscn)
vlab = 'kappa' ; xlab = 'Plasma_elongation'
case (26)
triang = swp(iscn)
vlab = 'triang' ; xlab = 'Plasma_triangularity'
case (27)
tbrmin = swp(iscn)
vlab = 'tbrmin' ; xlab = 'Min_tritium_breed._ratio'
case (28)
bt = swp(iscn)
vlab = 'bt' ; xlab = 'Tor._field_on_axis_(T)'
case (29)
coreradius = swp(iscn)
vlab = 'coreradius' ; xlab = 'Core_radius'
case (30)
!fimpvar = swp(iscn)
vlab = 'OBSOLETE' ; xlab = 'OBSOLETE'
case (31)
taulimit = swp(iscn)
vlab = 'taulimit' ; xlab = 'Taup/taueff_lower_limit'
case (32)
epsvmc = swp(iscn)
vlab = 'epsvmc' ; xlab = 'VMCON error tolerance'
case (33, 34, 35, 36, 37, 47)
write(*,*) 'Kallenbach model has been removed, remove the kallenbach scan variables'
stop 1
case (38)
boundu(129) = swp(iscn)
vlab = 'boundu(129)' ; xlab = ' Neon upper limit'
case (39)
boundu(131) = swp(iscn)
vlab = 'boundu(131)' ; xlab = ' Argon upper limit'
case (40)
boundu(135) = swp(iscn)
vlab = 'boundu(135)' ; xlab = ' Xenon upper limit'
case (41)
blnkoth = swp(iscn)
vlab = 'blnkoth' ; xlab = 'Outboard blanket thick.'
case (42)
fimp(9) = swp(iscn)
impurity_arr_frac(9) = fimp(9)
vlab = 'fimp(9)' ; xlab = 'Argon fraction'
case (43)
! rho_ecrh = swp(iscn)
vlab = 'obsolete' ; xlab = 'obsolete'
case (44)
sig_tf_case_max = swp(iscn)
vlab = 'sig_tf_case_max' ; xlab = 'Allowable_stress_in_tf_coil_case_Tresca_(pa)'
case (45)
tmargmin_tf = swp(iscn)
vlab = 'tmargmin_tf' ; xlab = 'Minimum_allowable_temperature_margin'
case (46)
boundu(152) = swp(iscn)
vlab = 'boundu(152)' ; xlab = 'Max allowable fgwsep'
case (48)
n_pancake = int(swp(iscn))
vlab = 'n_pancake' ; xlab = 'TF Coil - n_pancake'
case (49)
n_layer = int(swp(iscn))
vlab = 'n_layer' ; xlab = 'TF Coil - n_layer'
case (50)
fimp(13) = swp(iscn)
impurity_arr_frac(13) = fimp(13)
vlab = 'fimp(13)' ; xlab = 'Xenon fraction'
case (51)
ftar = swp(iscn)
vlab = 'ftar' ; xlab = 'lower_divertor_power_fraction'
case (52)
rad_fraction_sol = swp(iscn)
vlab = 'rad_fraction_sol' ; xlab = 'SoL radiation fraction'
case (53)
boundu(157) = swp(iscn)
vlab = 'boundu(157)' ; xlab = 'Max allowable fvssu'
case (54)
b_crit_upper_nbti = swp(iscn)
vlab = 'Bc2(0K)' ; xlab = 'GL_NbTi Bc2(0K)'
case(55)
shldith = swp(iscn)
vlab = 'shldith' ; xlab = 'Inboard neutronic shield'
case(56)
crypmw_max = swp(iscn)
vlab = 'crypmw_max' ; xlab = 'max allowable crypmw'
case(57)
boundl(2) = swp(iscn)
vlab = 'boundl(2)' ; xlab = 'bt minimum'
case(58)
scrapli = swp(iscn)
vlab = 'scrapli' ; xlab = 'Inboard FW-plasma sep gap'
case(59)
scraplo = swp(iscn)
vlab = 'scraplo' ; xlab = 'Outboard FW-plasma sep gap'
case (60)
sig_tf_wp_max = swp(iscn)
vlab = 'sig_tf_wp_max' ; xlab = 'Allowable_stress_in_tf_coil_conduit_Tresca_(pa)'
case (61)
copperaoh_m2_max = swp(iscn)
vlab = 'copperaoh_m2_max' ; xlab = 'Max CS coil current / copper area'
case (62)
coheof = swp(iscn)
vlab = 'coheof' ; xlab = 'CS coil current density at EOF (A/m2)'
case (63)
ohcth = swp(iscn)
vlab = 'ohcth' ; xlab = 'CS coil thickness (m)'
case (64)
ohhghf = swp(iscn)
vlab = 'ohhghf' ; xlab = 'CS height (m)'
case (65)
n_cycle_min = swp(iscn)
vlab = 'n_cycle_min' ; xlab = 'CS stress cycles min'
case (66)
oh_steel_frac = swp(iscn)
vlab = 'oh_steel_frac' ; xlab = 'CS steel fraction'
case (67)
t_crack_vertical = swp(iscn)
vlab = 't_crack_vertical' ; xlab = 'Initial crack vertical size (m)'
case (68)
inlet_temp_liq = swp(iscn)
vlab = 'inlet_temp_liq' ; xlab = 'Inlet Temperature Liquid Metal Breeder/Coolant (K)'
case (69)
outlet_temp_liq = swp(iscn)
vlab = 'outlet_temp_liq' ; xlab = 'Outlet Temperature Liquid Metal Breeder/Coolant (K)'
case(70)
blpressure_liq = swp(iscn)
vlab = 'blpressure_liq' ; xlab = 'Blanket liquid metal breeder/coolant pressure (Pa)'
case(71)
n_liq_recirc = swp(iscn)
vlab = 'n_liq_recirc' ; xlab = 'Selected number of liquid metal breeder recirculations per day'
case(72)
bz_channel_conduct_liq = swp(iscn)
vlab = 'bz_channel_conduct_liq' ; xlab = 'Conductance of liquid metal breeder duct walls (A V-1 m-1)'
case(73)
pnuc_fw_ratio_dcll = swp(iscn)
vlab = 'pnuc_fw_ratio_dcll' ; xlab = 'Ratio of FW nuclear power as fraction of total (FW+BB)'
case(74)
f_nuc_pow_bz_struct = swp(iscn)
vlab = 'f_nuc_pow_bz_struct' ; xlab = 'Fraction of BZ power cooled by primary coolant for dual-coolant balnket'
case(75)
pitch = swp(iscn)
vlab = 'pitch' ; xlab = 'pitch of first wall cooling channels (m)'
case (76)
etath = swp(iscn)
vlab = 'etath' ; xlab = 'Thermal conversion eff.'
case (77)
startupratio = swp(iscn)
vlab = 'startupratio' ; xlab = 'Gyrotron redundancy'
case (78)
fkind = swp(iscn)
vlab = 'fkind' ; xlab = 'Multiplier for Nth of a kind costs'
case (79)
etaech = swp(iscn)
vlab = 'etaech' ; xlab = 'ECH wall plug to injector efficiency'
case (80)
fcoolcp = swp(iscn)
vlab = 'fcoolcp' ; xlab = 'Coolant fraction of TF'
case (81)
n_tf_turn = swp(iscn)
vlab = 'n_tf_turn' ; xlab = 'Number of turns in TF'
case default
idiags(1) = nwp ; call report_error(96)
end select
end subroutine scan_select