gray/src/gray_params.f90

374 lines
15 KiB
Fortran
Raw Normal View History

2015-11-18 17:34:33 +01:00
module gray_params
2015-11-18 17:34:33 +01:00
use const_and_precisions, only : wp_
implicit none
integer, parameter :: lenfnm = 256
integer, parameter :: headw = 132, headl = 21
! Antenna/wave launcher parameters
type antenna_parameters
! From gray_params.data:
real(wp_) :: alpha, beta ! Launching angles
real(wp_) :: power ! Initial power
real(wp_) :: psi, chi ! Initial polarisation angles
integer :: iox ! Initial wave mode
integer :: ibeam ! Beam kind
character(len=lenfnm) :: filenm ! beamdata.txt filename
! From beamdata.txt:
real(wp_) :: fghz ! EC frequency
real(wp_) :: pos(3) ! Launcher position (tokamak frame)
real(wp_) :: w(2) ! Beam waists w(z) for the amplitude (local frame)
real(wp_) :: ri(2) ! Beam inverse radii 1/R(z) for the phase (local frame)
real(wp_) :: phi(2) ! Axes orientation angles for amplitude, phase (local frame)
end type
! MHD equilibrium parameters
type equilibrium_parameters
2015-11-18 17:34:33 +01:00
real(wp_) :: ssplps, ssplf, factb
integer :: sgnb, sgni, ixp
integer :: iequil, icocos, ipsinorm, idesc, ifreefmt
character(len=lenfnm) :: filenm
end type
2015-11-18 17:34:33 +01:00
! Kinetic plasma profiles
type profiles_parameters
real(wp_) :: psnbnd ! plasma boundary (ψ: ne(ψ)=0)
real(wp_) :: sspld, factne, factte
integer :: iscal, irho
2015-11-18 17:34:33 +01:00
integer :: iprof
character(len=lenfnm) :: filenm
end type
2015-11-18 17:34:33 +01:00
! Raytracing parameters
type raytracing_parameters
2015-11-18 17:34:33 +01:00
real(wp_) :: rwmax, dst
integer :: nrayr, nrayth, nstep
integer :: igrad, idst, ipass, ipol
end type
2015-11-18 17:34:33 +01:00
! EC resonant heating & Current Drive parameters
type ecrh_cd_parameters
2015-11-18 17:34:33 +01:00
integer :: iwarm, ilarm, imx, ieccd
end type
! Output data parameters
type output_parameters
integer :: ipec ! Grid spacing for profiles (profili EC)
integer :: nrho ! Number of grid steps for EC profiles + 1
integer :: istpr ! Subsampling factor for beam cross-section (units 8, 12)
integer :: istpl ! Subsampling factor for outer rays (unit 33)
end type
! Other parameters
type misc_parameters
real(wp_) :: rwall ! R of the limiter (fallback)
end type
! MHD equilibrium data
type equilibrium_data
real(wp_), allocatable :: rv(:) ! R of the uniform grid
real(wp_), allocatable :: zv(:) ! Z of the uniform grid
real(wp_), allocatable :: rlim(:) ! R of the limiter contour (wall)
real(wp_), allocatable :: zlim(:) ! Z of the limiter contour
real(wp_), allocatable :: rbnd(:) ! R of the boundary contour (plasma)
real(wp_), allocatable :: zbnd(:) ! Z of the boundary contour
real(wp_), allocatable :: fpol(:) ! Poloidal current function
real(wp_), allocatable :: qpsi(:) ! Safety factor on the flux grid
real(wp_), allocatable :: psin(:,:) ! Poloidal flux on a uniform grid
real(wp_), allocatable :: psinr(:) ! Poloidal flux
real(wp_) :: psia ! Poloidal flux at edge - flux at magnetic axis
real(wp_) :: rvac ! Reference R₀ (B = B₀R₀ without the plasma)
real(wp_) :: rax ! R of the magnetic axis
real(wp_) :: zax ! Z of the magnetic axis
end type
! Kinetic plasma profiles data
type profiles_data
real(wp_), allocatable :: psrad(:) ! Radial coordinate
real(wp_), allocatable :: terad(:) ! Electron temperature profile
real(wp_), allocatable :: derad(:) ! Electron density profile
real(wp_), allocatable :: zfc(:) ! Effective charge profile
end type
! All GRAY parameters
type gray_parameters
type(antenna_parameters) :: antenna
type(equilibrium_parameters) :: equilibrium
type(profiles_parameters) :: profiles
type(raytracing_parameters) :: raytracing
type(ecrh_cd_parameters) :: ecrh_cd
type(output_parameters) :: output
type(misc_parameters) :: misc
end type
! All GRAY input data
type gray_data
type(equilibrium_data) :: equilibrium
type(profiles_data) :: profiles
end type
! GRAY final results (only some)
type gray_results
real(wp_) :: pabs ! Total absorbed power
real(wp_) :: icd ! Total driven current
real(wp_), allocatable :: dpdv(:) ! Absorbed power density
real(wp_), allocatable :: jcd(:) ! Driven current density
end type
! Global variables exposed for graycore
integer, save :: iequil, iprof, ipol
integer, save :: iwarm, ilarm, imx, ieccd
integer, save :: igrad, idst, ipass
integer, save :: istpr0, istpl0
integer, save :: ipec, nnd
2015-11-18 17:34:33 +01:00
contains
2021-12-15 02:30:52 +01:00
subroutine print_parameters(params, strout)
implicit none
2021-12-15 02:30:52 +01:00
! subroutine arguments
type(gray_parameters), intent(in) :: params
character(len=*), dimension(:), intent(out) :: strout ! min len=110, dimension(21)
2021-12-15 02:30:52 +01:00
! local variables
character(len=8) :: date
character(len=10) :: time
#ifndef REVISION
character(len=*), parameter :: REVISION="unknown"
#endif
! Date and time
call date_and_time(date, time)
write(strout(1), '("# Run date/time: ",a4,2("/",a2),1x,2(a2,":"),a6)') &
date(1:4), date(5:6), date(7:8), &
time(1:2), time(3:4), time(5:10)
2021-12-15 02:30:52 +01:00
! Git revision
write(strout(2), '("# GRAY Git revision: ",a)') REVISION
2021-12-15 02:30:52 +01:00
! MHD equilibrium parameters
if (params%equilibrium%iequil > 0) then
write(strout(3), '("# EQL input: ",a)') trim(params%equilibrium%filenm)
! TODO add missing values
write(strout(7), '("# EQL B0 R0 aminor Rax zax:",5(1x,e12.5))') &
0._wp_, 0._wp_, 0._wp_, 0._wp_, 0._wp_
else
write(strout(3), '("# EQL input: N/A (vacuum)")')
write(strout(7), '("# EQL B0 R0 aminor Rax zax: N/A (vacuum)")')
end if
2021-12-15 02:30:52 +01:00
write(strout(4), '("# EQL iequil sgnb sgni factb:",3(1x,i4),1x,e12.5)') &
params%equilibrium%iequil, params%equilibrium%sgnb, params%equilibrium%sgni, params%equilibrium%factb
if (params%equilibrium%iequil > 1) then
write(strout(5), '("# EQL icocos ipsinorm idesc ifreefmt:",4(1x,i4))') &
params%equilibrium%icocos, params%equilibrium%ipsinorm, params%equilibrium%idesc, params%equilibrium%ifreefmt
write(strout(6), '("# EQL ssplps ssplf ixp:",2(1x,e12.5),1x,i4)') &
params%equilibrium%ssplps, params%equilibrium%ssplf, params%equilibrium%ixp
else
write(strout(5), '("# EQL icocos ipsinorm idesc ifreefmt: N/A (analytical)")')
write(strout(6), '("# EQL ssplps ssplf ixp: N/A (analytical)")')
end if
2021-12-15 02:30:52 +01:00
! Profiles parameters
if (params%equilibrium%iequil > 0) then
write(strout(8), '("# PRF input: ",a)') trim(params%profiles%filenm)
write(strout(9), '("# PRF iprof iscal factne factte:",2(1x,i4),2(1x,e12.5))') &
params%profiles%iprof, params%profiles%iscal, params%profiles%factne, params%profiles%factte
if (params%profiles%iprof > 0) then
write(strout(10), '("# PRF irho psnbnd sspld:",1x,i4,2(1x,e12.5))') &
params%profiles%irho,params%profiles%psnbnd,params%profiles%sspld
else
write(strout(10), '("# PRF irho psnbnd sspld: N/A (analytical)")')
end if
! TODO: add missing values
write(strout(11), '("# PRF Te0 ne0 Zeff0:",3(1x,e12.5))') &
0._wp_, 0._wp_, 0._wp_
else
write(strout(8), '("# PRF input: N/A (vacuum)")')
write(strout(9), '("# PRF iprof iscal factne factte: N/A (vacuum)")')
write(strout(10), '("# PRF irho psnbnd sspld: N/A (vacuum)")')
write(strout(11), '("# PRF Te0 ne0 Zeff0: N/A (vacuum)")')
end if
2021-12-15 02:30:52 +01:00
! Antenna parameters
write(strout(12), '("# ANT input: ",a)') trim(params%antenna%filenm)
write(strout(13), '("# ANT ibeam iox psi chi:",2(1x,i4),2(1x,e12.5))') &
params%antenna%ibeam, params%antenna%iox, params%antenna%psi, params%antenna%chi
write(strout(14), '("# ANT alpha beta power:",3(1x,e12.5))') &
params%antenna%alpha, params%antenna%beta, params%antenna%power
! TODO: add missing values
write(strout(15), '("# ANT x0 y0 z0:",3(1x,e12.5))') &
0._wp_, 0._wp_, 0._wp_
! TODO: add missing values
write(strout(16), '("# ANT wx wy Rcix Rciy psiw psir:",6(1x,e12.5))') &
0._wp_, 0._wp_, 0._wp_, 0._wp_, 0._wp_, 0._wp_
2021-12-15 02:30:52 +01:00
! Other parameters
write(strout(17), '("# RFL rwall:",1x,e12.5)') params%misc%rwall
2021-12-15 02:30:52 +01:00
! code parameters
write(strout(18), '("# COD igrad idst ipass ipol:",4(1x,i4))') &
params%raytracing%igrad, params%raytracing%idst, params%raytracing%ipass, params%raytracing%ipol
write(strout(19), '("# COD nrayr nrayth nstep rwmax dst:",3(1x,i4),2(1x,e12.5))') &
params%raytracing%nrayr, params%raytracing%nrayth, params%raytracing%nstep, params%raytracing%rwmax, params%raytracing%dst
write(strout(20), '("# COD iwarm ilarm imx ieccd:",4(1x,i4))') &
params%ecrh_cd%iwarm, params%ecrh_cd%ilarm, params%ecrh_cd%imx, params%ecrh_cd%ieccd
write(strout(21), '("# COD ipec nrho istpr istpl:",4(1x,i4))') &
params%output%ipec, params%output%nrho, params%output%istpr, params%output%istpl
end subroutine print_parameters
subroutine read_parameters(filename, params, unit)
2015-11-18 17:34:33 +01:00
use utils, only : get_free_unit
implicit none
2021-12-15 02:30:52 +01:00
! subrouting arguments
character(len=*), intent(in) :: filename
type(gray_parameters), intent(out) :: params
2015-11-18 17:34:33 +01:00
integer, intent(in), optional :: unit
2021-12-15 02:30:52 +01:00
! local variables
integer :: u, iostat
2015-11-18 17:34:33 +01:00
if (present(unit)) then
u = unit
2015-11-18 17:34:33 +01:00
else
u = get_free_unit()
end if
open(u, file=filename, status='old', action='read', iostat=iostat)
if (iostat > 0) then
print '(a)', 'gray_params file not found!'
call EXIT(1)
end if
2015-11-18 17:34:33 +01:00
2021-12-15 02:30:52 +01:00
! Raytracing
! ========================================================================
! nrayr :number of rays in radial direction
! nrayth :number of rays in angular direction
! rwmax :normalized maximum radius of beam power
! rwmax=1 -> :last ray at radius = waist
read(u, *) params%raytracing%nrayr, params%raytracing%nrayth, params%raytracing%rwmax
2021-12-15 02:30:52 +01:00
! igrad=0 :optical ray-tracing, initial conditions as for beam
! igrad=1 :quasi-optical ray-tracing
! igrad=-1 :ray-tracing, init. condit.
! from center of mirror and with angular spread
! ipass=1/2 :1 or 2 passes into plasma
! ipol=0 :compute mode polarization at antenna, ipol=1 use polariz angles
read(u, *) params%raytracing%igrad, params%raytracing%ipass, params%raytracing%ipol
2021-12-15 02:30:52 +01:00
! dst :integration step
! nstep :maximum number of integration steps
! idst=0/1/2 :0 integration in s, 1 integr. in ct, 2 integr. in Sr
read(u, *) params%raytracing%dst, params%raytracing%nstep, params%raytracing%idst
2021-12-15 02:30:52 +01:00
! Heating & Current drive
! ========================================================================
! iwarm=0 :no absorption and cd
! iwarm=1 :weakly relativistic absorption
! iwarm=2 :relativistic absorption, n<1 asymptotic expansion
! iwarm=3 :relativistic absorption, numerical integration
! ilarm :order of larmor expansion
! imx :max n of iterations in dispersion, imx<0 uses 1st
! iteration in case of failure after |imx| iterations
read(u, *) params%ecrh_cd%iwarm,params%ecrh_cd%ilarm,params%ecrh_cd%imx
2021-12-15 02:30:52 +01:00
! ieccd 0/1 NO/YES ECCD calculation ieccd>0 different CD models
read(u, *) params%ecrh_cd%ieccd
2021-12-15 02:30:52 +01:00
! Wave launcher
! ========================================================================
! alpha0, beta0 (cartesian) launching angles
read(u, *) params%antenna%alpha, params%antenna%beta
2021-12-15 02:30:52 +01:00
! p0mw injected power (MW)
read(u, *) params%antenna%power
2021-12-15 02:30:52 +01:00
! abs(iox)=1/2 OM/XM
! psipol0,chipol0 polarization angles at the antenna (if iox<0)
read(u, *) params%antenna%iox, params%antenna%psi, params%antenna%chi
read(u, *) params%antenna%ibeam
read(u, *) params%antenna%filenm
! MHD equilibrium
2021-12-15 02:30:52 +01:00
! ========================================================================
! iequil=0 :vacuum (no plasma)
2021-12-15 02:30:52 +01:00
! iequil=1 :analytical equilibrium
! iequil=2 :read eqdsk
! iequil=2 :read eqdsk, data only valid inside last closed flux surface
read(u, *) params%equilibrium%iequil
read(u, *) params%equilibrium%filenm
2021-12-15 02:30:52 +01:00
! icocos :index for equilibrium from COCOS - O. Sauter Feb 2012
! ipsinorm :0 standard EQDSK format, 1 format Portone summer 2004
read(u, *) params%equilibrium%icocos, params%equilibrium%ipsinorm, params%equilibrium%idesc, params%equilibrium%ifreefmt
2021-12-15 02:30:52 +01:00
! ixp=0,-1,+1 : no X point , bottom/up X point
! ssplps : spline parameter for psi interpolation
read(u, *) params%equilibrium%ixp, params%equilibrium%ssplps !, params%equilibrium%ssplf
params%equilibrium%ssplf=0.01_wp_
2021-12-15 02:30:52 +01:00
! signum of toroidal B and I
! factb factor for magnetic field (only for numerical equil)
! scaling adopted: beta=const, qpsi=const, nustar=const
read(u, *) params%equilibrium%sgnb, params%equilibrium%sgni, params%equilibrium%factb
2021-12-15 02:30:52 +01:00
! Wall
! ========================================================================
read(u, *) params%misc%rwall
2021-12-15 02:30:52 +01:00
! Profiles
! ========================================================================
! iprof=0 :analytical density and temp. profiles
! iprof>0 :numerical density and temp. profiles
read(u, *) params%profiles%iprof, params%profiles%irho ! irho=0,1,2 -> num profiles vs rhot,rhop,psin
read(u, *) params%profiles%filenm
2021-12-15 02:30:52 +01:00
! psbnd value of psi ( > 1 ) of density boundary
read(u, *) params%profiles%psnbnd, params%profiles%sspld
2021-12-15 02:30:52 +01:00
! prfparam%sspld=0.001_wp_
! iscal :ne Te scaling 0: nustar=const, 1: n_greenw=const; 2 no rescaling
! factT factn :factor for Te&ne scaling
read(u, *) params%profiles%factte, params%profiles%factne, params%profiles%iscal
2021-12-15 02:30:52 +01:00
! Output
! ========================================================================
! ipec=0/1 :pec profiles grid in psi/rhop
! nrho :number of grid steps for pec profiles +1
read(u, *) params%output%ipec, params%output%nrho
2021-12-15 02:30:52 +01:00
! istpr0 :projection step = dsdt*istprj
! istpl0 :plot step = dsdt*istpl
read(u, *) params%output%istpr, params%output%istpl
2015-11-18 17:34:33 +01:00
close(u)
end subroutine read_parameters
2015-11-18 17:34:33 +01:00
2021-12-15 02:30:52 +01:00
subroutine set_globals(params)
2015-11-18 17:34:33 +01:00
implicit none
! subroutine arguments
type(gray_parameters), intent(in) :: params
iequil = params%equilibrium%iequil
iprof = params%profiles%iprof
ipec = params%output%ipec
nnd = params%output%nrho
istpr0 = params%output%istpr
istpl0 = params%output%istpl
ipol = params%raytracing%ipol
igrad = params%raytracing%igrad
idst = params%raytracing%idst
ipass = params%raytracing%ipass
if (params%raytracing%nrayr < 5) then
igrad = 0
print *, ' nrayr < 5 ! => OPTICAL CASE ONLY'
print *, ' '
2015-11-18 17:34:33 +01:00
end if
iwarm = params%ecrh_cd%iwarm
ilarm = params%ecrh_cd%ilarm
imx = params%ecrh_cd%imx
ieccd = params%ecrh_cd%ieccd
2015-11-18 17:34:33 +01:00
end subroutine set_globals
2015-11-18 17:34:33 +01:00
end module gray_params