2021-12-19 16:39:19 +01:00
! This modules contains the core GRAY routines
2021-12-15 02:31:14 +01:00
module gray_core
2021-12-19 16:39:19 +01:00
2015-11-18 17:34:33 +01:00
use const_and_precisions , only : wp_
2021-12-19 16:39:19 +01:00
2015-11-18 17:34:33 +01:00
implicit none
contains
2021-12-15 02:31:09 +01:00
subroutine gray_main ( params , data , results , error , rhout )
2019-03-26 15:21:22 +01:00
use const_and_precisions , only : zero , one , degree , comp_tiny
2021-12-19 16:39:19 +01:00
use coreprofiles , only : temp , fzeff
2021-12-15 02:31:09 +01:00
use dispersion , only : expinit
2021-12-19 16:39:19 +01:00
use gray_params , only : gray_parameters , gray_data , gray_results , &
2022-05-10 12:36:37 +02:00
print_parameters
2021-12-15 02:31:09 +01:00
use beams , only : xgygcoeff , launchangles2n
use beamdata , only : pweight , rayi2jk
2022-07-14 00:38:34 +02:00
use gray_errors , only : is_critical , print_err_raytracing , print_err_ecrh_cd
2016-06-01 15:49:35 +02:00
use magsurf_data , only : flux_average , dealloc_surfvec
2022-05-10 12:36:37 +02:00
use beamdata , only : init_btr , dealloc_beam
2021-12-15 02:31:09 +01:00
use pec , only : pec_init , spec , postproc_profiles , dealloc_pec , &
rhop_tab , rhot_tab
use utils , only : vmaxmin
use reflections , only : inside
2021-12-19 16:39:19 +01:00
use multipass , only : alloc_multipass , dealloc_multipass , initbeam , &
initmultipass , turnoffray , plasma_in , plasma_out , &
wall_out
2021-12-18 18:57:38 +01:00
use logger , only : log_info , log_debug
2019-03-26 15:21:22 +01:00
2015-11-18 17:34:33 +01:00
implicit none
2021-12-15 02:31:09 +01:00
2021-12-19 16:39:19 +01:00
! subroutine arguments
2022-05-10 12:36:37 +02:00
type ( gray_parameters ) , intent ( inout ) :: params
type ( gray_data ) , intent ( in ) :: data
type ( gray_results ) , intent ( out ) :: results
2021-12-15 02:31:09 +01:00
! Predefined grid for the output profiles (optional)
2015-11-24 17:36:20 +01:00
real ( wp_ ) , dimension ( : ) , intent ( in ) , optional :: rhout
2015-11-18 17:34:33 +01:00
2021-12-15 02:31:09 +01:00
! Exit code
integer , intent ( out ) :: error
! local variables
2019-03-26 15:21:22 +01:00
real ( wp_ ) , parameter :: taucr = 1 2._wp_ , etaucr = exp ( - taucr )
2022-05-10 12:36:37 +02:00
character , dimension ( 2 ) , parameter :: mode = [ 'O' , 'X' ]
2021-12-18 18:57:38 +01:00
2022-05-20 13:10:49 +02:00
integer :: sox
real ( wp_ ) :: ak0 , bres , xgcn , xg , yg , rrm , zzm , alpha , didp , anpl , anpr , anprim , anprre
2019-03-26 15:21:22 +01:00
real ( wp_ ) :: chipol , psipol , btot , psinv , dens , tekev , dersdst , derdnm
2015-11-19 19:20:58 +01:00
real ( wp_ ) :: tau , pow , dids , ddr , ddi , taumn , taumx
2015-11-23 18:55:27 +01:00
real ( wp_ ) :: rhotpav , drhotpav , rhotjava , drhotjava , dpdvp , jphip
real ( wp_ ) :: rhotp , drhotp , rhotj , drhotj , dpdvmx , jphimx , ratjamx , ratjbmx
2019-03-26 15:21:22 +01:00
real ( wp_ ) :: pabs_beam , icd_beam , cpl_beam1 , cpl_beam2 , cpl_cbeam1 , cpl_cbeam2
2015-11-23 18:55:27 +01:00
2019-03-26 15:21:22 +01:00
real ( wp_ ) , dimension ( 2 ) :: pabs_pass , icd_pass , cpl , cpl0
2017-09-12 21:37:06 +02:00
real ( wp_ ) , dimension ( 3 ) :: xv , anv0 , anv , bv , derxg
2021-12-15 02:31:09 +01:00
! Ray variables
real ( wp_ ) , dimension ( : , : ) , pointer :: yw = > null ( ) , ypw = > null ( ) , gri = > null ( )
2016-06-01 15:49:35 +02:00
real ( wp_ ) , dimension ( : , : , : ) , pointer :: xc = > null ( ) , du1 = > null ( ) , ggri = > null ( )
2021-12-15 02:31:09 +01:00
! i: integration step, jk: global ray index
integer :: i , jk
2022-07-14 00:38:34 +02:00
integer :: iox , nharm , nhf , nnd , iokhawa , ierrn , ierrhcd , index_rt
2019-03-26 15:21:22 +01:00
integer :: ip , ib , iopmin , ipar , iO
2021-12-15 02:31:09 +01:00
integer :: igrad_b , istop_pass , nbeam_pass , nlim
2019-03-26 15:21:22 +01:00
logical :: ins_pl , ins_wl , ent_pl , ext_pl , ent_wl , ext_wl , iboff
real ( wp_ ) , dimension ( : , : , : ) , pointer :: yynext = > null ( ) , yypnext = > null ( )
real ( wp_ ) , dimension ( : , : ) , pointer :: psjki = > null ( ) , ppabs = > null ( ) , ccci = > null ( )
real ( wp_ ) , dimension ( : , : ) , pointer :: taus = > null ( ) , stnext = > null ( ) , &
yw0 = > null ( ) , ypw0 = > null ( ) , cpls = > null ( )
real ( wp_ ) , dimension ( : ) , pointer :: p0ray = > null ( ) , tau0 = > null ( ) , alphaabs0 = > null ( ) , &
dids0 = > null ( ) , ccci0 = > null ( ) , tau1 = > null ( ) , etau1 = > null ( ) , cpl1 = > null ( ) , lgcpl1 = > null ( )
real ( wp_ ) , dimension ( : ) , pointer :: p0jk = > null ( )
real ( wp_ ) , dimension ( : ) , pointer :: jphi_beam = > null ( ) , pins_beam = > null ( ) , &
currins_beam = > null ( ) , dpdv_beam = > null ( ) , jcd_beam = > null ( ) , stv = > null ( ) , &
psipv = > null ( ) , chipv = > null ( )
complex ( wp_ ) , dimension ( : ) , pointer :: ext = > null ( ) , eyt = > null ( )
integer , dimension ( : ) , pointer :: iiv = > null ( ) , iop = > null ( ) , iow = > null ( )
logical , dimension ( : ) , pointer :: iwait = > null ( )
logical , dimension ( : , : ) , pointer :: iroff = > null ( )
2015-11-18 17:34:33 +01:00
2021-12-18 18:57:38 +01:00
! buffer for formatting log messages
character ( 256 ) :: msg
2016-04-27 16:37:57 +02:00
2021-12-15 02:31:09 +01:00
! ======== set environment BEGIN ========
! Number of limiter contourn points
nlim = size ( data % equilibrium % zlim )
2015-11-18 17:34:33 +01:00
2021-12-15 02:31:09 +01:00
! Compute X=ω/ω_ce and Y=(ω/ω_pe)² (with B=1)
call xgygcoeff ( params % antenna % fghz , ak0 , bres , xgcn )
2015-11-18 17:34:33 +01:00
2021-12-15 02:31:09 +01:00
! Compute the initial cartesian wavevector (anv0)
2023-03-26 16:03:57 +02:00
! from launch angles α ,β and the position
2021-12-15 02:31:09 +01:00
call launchangles2n ( params % antenna , anv0 )
2015-11-18 17:34:33 +01:00
2021-12-15 02:31:09 +01:00
! Initialise the ray variables (beamtracing)
call init_btr ( params % raytracing , yw , ypw , xc , du1 , gri , ggri , psjki , ppabs , ccci , &
tau0 , alphaabs0 , dids0 , ccci0 , p0jk , ext , eyt , iiv )
2015-11-18 17:34:33 +01:00
2021-12-15 02:31:09 +01:00
! Initialise the dispersion module
2022-05-10 12:36:37 +02:00
if ( params % ecrh_cd % iwarm > 1 ) call expinit
2021-12-15 02:31:09 +01:00
! Initialise the magsurf_data module
call flux_average ! requires frhotor for dadrhot,dvdrhot
2022-05-10 12:36:37 +02:00
2021-12-15 02:31:09 +01:00
! Initialise the output profiles
2022-05-10 12:36:37 +02:00
call pec_init ( params % output % ipec , rhout )
2021-12-15 02:31:09 +01:00
nnd = size ( rhop_tab ) ! number of radial profile points
call alloc_multipass ( nnd , iwait , iroff , iop , iow , yynext , yypnext , yw0 , ypw0 , stnext , &
stv , p0ray , taus , tau1 , etau1 , cpls , cpl1 , lgcpl1 , jphi_beam , &
pins_beam , currins_beam , dpdv_beam , jcd_beam , psipv , chipv )
! Allocate memory for the results...
allocate ( results % dpdv ( params % output % nrho ) )
allocate ( results % jcd ( params % output % nrho ) )
! ...and initialise them
results % pabs = zero
results % icd = zero
results % dpdv = zero
results % jcd = zero
! ========= set environment END =========
! ======== pre-proc prints BEGIN ========
2022-05-03 23:16:21 +02:00
call print_headers ( params )
2021-12-15 02:31:09 +01:00
2021-12-18 18:57:38 +01:00
! print ψ surface for q=1.5 and q=2 on file and log psi,rhot,rhop
2021-12-15 02:31:09 +01:00
call print_surfq ( [ 1.5_wp_ , 2.0_wp_ ] )
! print initial position
2021-12-18 18:57:38 +01:00
write ( msg , '("initial position:",3(x,g0.3))' ) params % antenna % pos
call log_info ( msg , mod = 'gray_core' , proc = 'gray_main' )
write ( msg , '("initial direction:",2(x,a,"=",g0.2))' ) &
'α ' , params % antenna % alpha , 'β' , params % antenna % beta
call log_info ( msg , mod = 'gray_core' , proc = 'gray_main' )
2021-12-15 02:31:09 +01:00
! print Btot=Bres
! print ne, Te, q, Jphi versus psi, rhop, rhot
2015-11-23 18:55:27 +01:00
call print_bres ( bres )
2023-05-23 17:49:51 +02:00
call print_prof ( params % profiles )
2021-12-19 02:35:24 +01:00
call print_maps ( bres , xgcn , &
norm2 ( params % antenna % pos ( 1 : 2 ) ) * 0.01_wp_ , &
2021-12-15 02:31:09 +01:00
sin ( params % antenna % beta * degree ) )
! ========= pre-proc prints END =========
! =========== main loop BEGIN ===========
call initmultipass ( params % raytracing % ipol , params % antenna % iox , &
iroff , yynext , yypnext , yw0 , ypw0 , &
stnext , p0ray , taus , tau1 , etau1 , cpls , cpl1 , lgcpl1 , psipv , chipv )
2022-05-10 12:36:37 +02:00
if ( params % raytracing % ipol == 0 ) then
if ( params % antenna % iox == 2 ) then ! only X mode on 1st pass
cpl0 = [ zero , one ]
2019-03-26 15:21:22 +01:00
else ! only O mode on 1st pass
2022-05-10 12:36:37 +02:00
cpl0 = [ one , zero ]
2019-03-26 15:21:22 +01:00
end if
end if
2021-12-15 02:31:09 +01:00
2022-05-20 13:10:49 +02:00
sox = 1 ! mode inverted for each beam
2019-03-26 15:21:22 +01:00
iox = 2 ! start with O: sox=-1, iox=1
2022-05-10 12:36:37 +02:00
2021-12-15 02:31:09 +01:00
psipol = params % antenna % psi
chipol = params % antenna % chi
call pweight ( params % antenna % power , p0jk )
2015-11-18 17:34:33 +01:00
2019-03-26 15:21:22 +01:00
nbeam_pass = 1 ! max n of beam per pass
index_rt = 0 ! global beam index: 1,O 2,X 1st pass
2019-03-28 10:50:28 +01:00
! | | | |
2021-12-18 18:57:38 +01:00
call log_debug ( 'pass loop start' , mod = 'gray_core' , proc = 'gray_main' ) ! 3,O 4,X 5,O 6,X 2nd pass
2022-05-10 12:36:37 +02:00
main_loop : do ip = 1 , params % raytracing % ipass
2021-12-18 18:57:38 +01:00
write ( msg , '("pass: ",g0)' ) ip
call log_info ( msg , mod = 'gray_core' , proc = 'gray_main' )
2019-03-26 15:21:22 +01:00
pabs_pass = zero
icd_pass = zero
istop_pass = 0 ! stop flag for current pass
nbeam_pass = 2 * nbeam_pass ! max n of beams in current pass
2022-05-10 12:36:37 +02:00
if ( ip > 1 ) then
2019-03-28 10:50:28 +01:00
du1 = zero
2019-03-26 15:21:22 +01:00
gri = zero
ggri = zero
2022-05-10 12:36:37 +02:00
if ( ip == params % raytracing % ipass ) cpl = [ zero , zero ] ! no successive passes
2019-03-26 15:21:22 +01:00
end if
2022-05-10 12:36:37 +02:00
! =========== beam loop BEGIN ===========
2021-12-18 18:57:38 +01:00
call log_debug ( 'beam loop start' , mod = 'gray_core' , proc = 'gray_main' )
2022-05-10 12:36:37 +02:00
beam_loop : do ib = 1 , nbeam_pass
2019-03-26 15:21:22 +01:00
sox = - sox ! invert mode
iox = 3 - iox ! O-mode at ip=1,ib=1
2022-05-10 12:36:37 +02:00
2019-03-26 15:21:22 +01:00
index_rt = index_rt + 1
iO = 2 * index_rt + 1 ! * index_rt of O-mode derived ray (iX=iO+1)
2015-11-18 17:34:33 +01:00
2019-03-26 15:21:22 +01:00
call initbeam ( index_rt , iroff , iboff , iwait , stv , jphi_beam , &
pins_beam , currins_beam , dpdv_beam , jcd_beam )
2021-12-18 18:57:38 +01:00
write ( msg , '(" beam: ",g0," (",a1," mode)")' ) index_rt , mode ( iox )
call log_info ( msg , mod = 'gray_core' , proc = 'gray_main' )
2019-03-26 15:21:22 +01:00
if ( iboff ) then ! no propagation for current beam
istop_pass = istop_pass + 1 ! * +1 non propagating beam
2021-12-18 18:57:38 +01:00
call log_info ( " beam is off" , mod = 'gray_core' , proc = 'gray_main' )
2019-03-26 15:21:22 +01:00
cycle
2015-11-18 17:34:33 +01:00
end if
2022-05-10 12:36:37 +02:00
2019-03-26 15:21:22 +01:00
call vectinit ( psjki , ppabs , ccci , tau0 , alphaabs0 , dids0 , ccci0 , iiv )
2022-05-10 12:36:37 +02:00
if ( ip == 1 ) then ! 1st pass
igrad_b = params % raytracing % igrad ! * input value, igrad_b=0 from 2nd pass
2019-03-28 10:50:28 +01:00
tau1 = zero ! * tau from previous passes
2019-03-26 15:21:22 +01:00
etau1 = one
2019-03-28 10:50:28 +01:00
cpl1 = one ! * coupling from previous passes
2019-03-26 15:21:22 +01:00
lgcpl1 = zero
p0ray = p0jk ! * initial beam power
2022-05-10 12:36:37 +02:00
call ic_gb ( params , anv0 , ak0 , yw , ypw , & ! * initial conditions
xc , du1 , gri , ggri , index_rt )
call set_pol ( yw , bres , sox , params % raytracing % ipol , & ! * initial polarization
psipol , chipol , ext , eyt )
do jk = 1 , params % raytracing % nray
2019-03-26 15:21:22 +01:00
zzm = yw ( 3 , jk ) * 0.01_wp_
rrm = sqrt ( yw ( 1 , jk ) * yw ( 1 , jk ) + yw ( 2 , jk ) * yw ( 2 , jk ) ) * 0.01_wp_
2022-05-10 12:36:37 +02:00
2021-12-15 02:31:09 +01:00
if ( inside ( data % equilibrium % rlim , data % equilibrium % zlim , &
nlim , rrm , zzm ) ) then ! * start propagation in/outside vessel?
2019-03-26 15:21:22 +01:00
iow ( jk ) = 1 ! + inside
else
iow ( jk ) = 0 ! + outside
end if
end do
2022-05-10 12:36:37 +02:00
2019-03-28 10:50:28 +01:00
else ! 2nd+ passes
2019-03-26 15:21:22 +01:00
ipar = ( index_rt + 1 ) / 2 - 1 ! * parent beam index
yw = yynext ( : , : , ipar ) ! * starting coordinates from
ypw = yypnext ( : , : , ipar ) ! parent beam last step
stv = stnext ( : , ipar ) ! * starting step from parent beam last step
iow = 1 ! * start propagation inside vessel
tau1 = taus ( : , index_rt ) ! * tau from previous passes
etau1 = exp ( - tau1 )
cpl1 = cpls ( : , index_rt ) ! * coupling from previous passes
lgcpl1 = - log ( cpl1 )
p0ray = p0jk * etau1 * cpl1 ! * initial beam power
end if
iop = 0 ! start propagation outside plasma
2022-05-10 12:36:37 +02:00
if ( params % raytracing % nray > 1 . and . all ( . not . iwait ) ) & ! iproj=0 ==> nfilp=8
call print_projxyzt ( stv , yw , 0 )
! ======= propagation loop BEGIN =======
2021-12-18 18:57:38 +01:00
call log_debug ( ' propagation loop start' , mod = 'gray_core' , proc = 'gray_main' )
2022-05-10 12:36:37 +02:00
propagation_loop : do i = 1 , params % raytracing % nstep
! advance one step with "frozen" grad(S_I)
do jk = 1 , params % raytracing % nray
2019-03-26 15:21:22 +01:00
if ( iwait ( jk ) ) cycle ! jk ray is waiting for next pass
2022-05-10 12:36:37 +02:00
stv ( jk ) = stv ( jk ) + params % raytracing % dst ! current ray step
2019-03-26 15:21:22 +01:00
call rkstep ( sox , bres , xgcn , yw ( : , jk ) , ypw ( : , jk ) , gri ( : , jk ) , ggri ( : , : , jk ) , igrad_b )
end do
2022-05-10 12:36:37 +02:00
! update position and grad
2019-03-26 15:21:22 +01:00
if ( igrad_b == 1 ) call gradi_upd ( yw , ak0 , xc , du1 , gri , ggri )
2022-05-10 12:36:37 +02:00
2021-12-15 02:31:09 +01:00
error = 0
2019-03-26 15:21:22 +01:00
iopmin = 10
2022-05-10 12:36:37 +02:00
! =========== rays loop BEGIN ===========
rays_loop : do jk = 1 , params % raytracing % nray
2019-03-26 15:21:22 +01:00
if ( iwait ( jk ) ) cycle ! jk ray is waiting for next pass
2022-05-10 12:36:37 +02:00
! compute derivatives with updated gradient and local plasma values
2019-03-26 15:21:22 +01:00
xv = yw ( 1 : 3 , jk )
anv = yw ( 4 : 6 , jk )
call ywppla_upd ( xv , anv , gri ( : , jk ) , ggri ( : , : , jk ) , sox , bres , xgcn , ypw ( : , jk ) , &
psinv , dens , btot , bv , xg , yg , derxg , anpl , anpr , ddr , ddi , dersdst , derdnm , &
ierrn , igrad_b )
2022-05-10 12:36:37 +02:00
! update global error code and print message
2019-03-26 15:21:22 +01:00
if ( ierrn / = 0 ) then
2021-12-15 02:31:09 +01:00
error = ior ( error , ierrn )
2022-07-14 00:38:34 +02:00
call print_err_raytracing ( ierrn , i , anpl )
2019-03-26 15:21:22 +01:00
end if
2022-05-10 12:36:37 +02:00
! check entrance/exit plasma/wall
2019-03-26 15:21:22 +01:00
zzm = xv ( 3 ) * 0.01_wp_
rrm = sqrt ( xv ( 1 ) * xv ( 1 ) + xv ( 2 ) * xv ( 2 ) ) * 0.01_wp_
2022-05-10 12:36:37 +02:00
ins_pl = ( psinv > = zero . and . psinv < params % profiles % psnbnd ) ! in/out plasma?
2021-12-15 02:31:09 +01:00
ins_wl = inside ( data % equilibrium % rlim , data % equilibrium % zlim , &
nlim , rrm , zzm ) ! in/out vessel?
2022-05-10 12:36:37 +02:00
ent_pl = ( mod ( iop ( jk ) , 2 ) == 0 . and . ins_pl ) ! enter plasma
ext_pl = ( mod ( iop ( jk ) , 2 ) == 1 . and . . not . ins_pl ) ! exit plasma
ent_wl = ( mod ( iow ( jk ) , 2 ) == 0 . and . ins_wl ) ! enter vessel
ext_wl = ( mod ( iow ( jk ) , 2 ) == 1 . and . . not . ins_wl ) ! exit vessel
2019-03-26 15:21:22 +01:00
if ( ent_pl ) then ! ray enters plasma
2022-05-26 02:33:55 +02:00
write ( msg , '(" ray ",g0," entered plasma (",g0," steps)")' ) jk , i
call log_debug ( msg , mod = 'gray_core' , proc = 'gray_main' )
2019-03-26 15:21:22 +01:00
call plasma_in ( jk , xv , anv , bres , sox , cpl , psipol , chipol , iop , ext , eyt )
2022-05-10 12:36:37 +02:00
if ( iop ( jk ) == 1 . and . ip == 1 ) then ! * 1st entrance on 1st pass (ray hasn't entered in plasma yet) => continue current pass
if ( params % raytracing % ipol == 0 ) then ! + IF single mode propagation
2019-03-26 15:21:22 +01:00
cpl = cpl0
p0ray ( jk ) = p0ray ( jk ) * cpl ( iox )
2022-05-10 12:36:37 +02:00
else if ( cpl ( iox ) < etaucr ) then ! + ELSE IF low coupled power for current mode => de-activate derived rays
2019-03-26 15:21:22 +01:00
call turnoffray ( jk , ip + 1 , 2 * ib + 2 - iox , iroff )
iwait ( jk ) = . true . ! . stop advancement and H&CD computation for current ray
if ( cpl ( iox ) . le . comp_tiny ) cpl ( iox ) = etaucr
else ! + ELSE assign coupled power to current ray
p0ray ( jk ) = p0ray ( jk ) * cpl ( iox )
end if
cpls ( jk , index_rt ) = cpl ( iox )
2022-05-10 12:36:37 +02:00
if ( jk == 1 ) then
2021-12-18 18:57:38 +01:00
write ( msg , '(" 1st pass - central ray (",a1,"-mode) c=",g0.4)' ) &
mode ( iox ) , cpl ( iox )
call log_info ( msg , mod = 'gray_core' , proc = 'gray_main' )
2019-03-26 15:21:22 +01:00
psipv ( index_rt ) = psipol ! + polarization angles at plasma boundary for central ray
chipv ( index_rt ) = chipol
end if
2022-05-10 12:36:37 +02:00
else if ( iop ( jk ) > 2 ) then ! * 2nd entrance on 1st pass / entrance on 2nd+ pass => end of current pass for ray jk
2019-03-26 15:21:22 +01:00
igrad_b = 0 ! + switch to ray-tracing
iwait ( jk ) = . true . ! + stop advancement and H&CD computation for current ray
2022-05-10 12:36:37 +02:00
if ( ip < params % raytracing % ipass ) then ! + not last pass
2019-03-28 10:50:28 +01:00
yynext ( : , jk , index_rt ) = yw0 ( : , jk ) ! . copy starting coordinates
yypnext ( : , jk , index_rt ) = ypw0 ( : , jk ) ! for next pass from last step
2022-05-10 12:36:37 +02:00
stnext ( jk , index_rt ) = stv ( jk ) - params % raytracing % dst ! . starting step for next pass = last step
if ( cpl ( 1 ) < etaucr ) then ! . low coupled power for O-mode => de-activate derived rays
2019-03-26 15:21:22 +01:00
call turnoffray ( jk , ip + 1 , 2 * ib - 1 , iroff )
if ( cpl ( 1 ) . le . comp_tiny ) cpl ( 1 ) = etaucr
end if
2022-05-10 12:36:37 +02:00
if ( cpl ( 2 ) < etaucr ) then ! . low coupled power for X-mode => de-activate derived rays
2019-03-26 15:21:22 +01:00
call turnoffray ( jk , ip + 1 , 2 * ib , iroff )
if ( cpl ( 2 ) . le . comp_tiny ) cpl ( 2 ) = etaucr
end if
2022-05-10 12:36:37 +02:00
2019-03-28 10:50:28 +01:00
taus ( jk , iO : iO + 1 ) = tau1 ( jk ) + tau0 ( jk ) ! . starting tau for next O-mode pass
cpls ( jk , iO ) = cpl1 ( jk ) * cpl ( 1 ) ! . cumulative coupling for next O-mode pass
cpls ( jk , iO + 1 ) = cpl1 ( jk ) * cpl ( 2 ) ! . cumulative coupling for next X-mode pass
2022-05-10 12:36:37 +02:00
if ( jk == 1 ) then ! . polarization angles at plasma boundary for central ray
2019-03-26 15:21:22 +01:00
psipv ( iO : iO + 1 ) = psipol
chipv ( iO : iO + 1 ) = chipol
end if
2019-03-28 10:50:28 +01:00
else ! * 1st entrance on 2nd+ pass (ray hasn't entered in plasma since end of previous pass) => continue current pass
2022-05-10 12:36:37 +02:00
cpl = [ zero , zero ]
2019-03-26 15:21:22 +01:00
end if
end if
2022-05-10 12:36:37 +02:00
2019-03-26 15:21:22 +01:00
else if ( ext_pl ) then ! ray exits plasma
2022-05-26 02:33:55 +02:00
write ( msg , '(" ray ", g0, " left plasma")' ) jk
call log_debug ( msg , mod = 'gray_core' , proc = 'gray_main' )
2019-03-26 15:21:22 +01:00
call plasma_out ( jk , xv , anv , bres , sox , iop , ext , eyt )
end if
2022-05-10 12:36:37 +02:00
2019-03-26 15:21:22 +01:00
if ( ent_wl ) then ! ray enters vessel
2019-03-28 10:50:28 +01:00
iow ( jk ) = iow ( jk ) + 1 ! * out->in
2022-05-10 12:36:37 +02:00
2019-03-26 15:21:22 +01:00
else if ( ext_wl ) then ! ray exit vessel
call wall_out ( jk , ins_pl , xv , anv , bres , sox , psipol , chipol , iow , iop , ext , eyt )
2022-05-10 12:36:37 +02:00
yw ( : , jk ) = [ xv , anv ] ! * updated coordinates (reflected)
2019-03-26 15:21:22 +01:00
igrad_b = 0 ! * switch to ray-tracing
2022-05-10 12:36:37 +02:00
2019-03-26 15:21:22 +01:00
call ywppla_upd ( xv , anv , gri ( : , jk ) , ggri ( : , : , jk ) , sox , bres , xgcn , ypw ( : , jk ) , &
psinv , dens , btot , bv , xg , yg , derxg , anpl , anpr , ddr , ddi , dersdst , derdnm , &
ierrn , igrad_b ) ! * update derivatives after reflection
if ( ierrn / = 0 ) then ! * update global error code and print message
2021-12-15 02:31:09 +01:00
error = ior ( error , ierrn )
2022-07-14 00:38:34 +02:00
call print_err_raytracing ( ierrn , i , anpl )
2019-03-26 15:21:22 +01:00
end if
2022-05-10 12:36:37 +02:00
if ( jk == 1 . and . ip == 1 ) then ! * 1st pass, polarization angles at reflection for central ray
2019-03-26 15:21:22 +01:00
psipv ( index_rt ) = psipol
chipv ( index_rt ) = chipol
end if
2022-05-10 12:36:37 +02:00
2019-03-28 10:50:28 +01:00
if ( ins_pl ) then ! * plasma-wall overlapping => wall+plasma crossing => end of current pass
2019-03-26 15:21:22 +01:00
iwait ( jk ) = . true . ! + stop advancement and H&CD computation for current ray
2022-05-10 12:36:37 +02:00
if ( ip < params % raytracing % ipass ) then ! + not last pass
yynext ( : , jk , index_rt ) = [ xv , anv ] ! . starting coordinates
2019-03-26 15:21:22 +01:00
yypnext ( : , jk , index_rt ) = ypw ( : , jk ) ! for next pass = reflection point
stnext ( jk , index_rt ) = stv ( jk ) ! . starting step for next pass = step after reflection
2022-05-10 12:36:37 +02:00
2019-03-26 15:21:22 +01:00
call plasma_in ( jk , xv , anv , bres , sox , cpl , psipol , chipol , iop , ext , eyt ) ! . ray re-enters plasma after reflection
2022-05-10 12:36:37 +02:00
if ( cpl ( 1 ) < etaucr ) then ! . low coupled power for O-mode? => de-activate derived rays
2019-03-26 15:21:22 +01:00
call turnoffray ( jk , ip + 1 , 2 * ib - 1 , iroff )
if ( cpl ( 1 ) . le . comp_tiny ) cpl ( 1 ) = etaucr
end if
2022-05-10 12:36:37 +02:00
if ( cpl ( 2 ) < etaucr ) then ! . low coupled power for X-mode? => de-activate derived rays
2019-03-26 15:21:22 +01:00
call turnoffray ( jk , ip + 1 , 2 * ib , iroff )
if ( cpl ( 2 ) . le . comp_tiny ) cpl ( 2 ) = etaucr
end if
2022-05-10 12:36:37 +02:00
2019-03-26 15:21:22 +01:00
taus ( jk , iO : iO + 1 ) = tau1 ( jk ) + tau0 ( jk ) ! . starting tau for next O-mode pass
cpls ( jk , iO ) = cpl1 ( jk ) * cpl ( 1 ) ! . cumulative coupling for next O-mode pass
cpls ( jk , iO + 1 ) = cpl1 ( jk ) * cpl ( 2 ) ! . cumulative coupling for next X-mode pass
2022-05-10 12:36:37 +02:00
if ( jk == 1 ) then ! + polarization angles at plasma boundary for central ray
2019-03-26 15:21:22 +01:00
psipv ( iO : iO + 1 ) = psipol
chipv ( iO : iO + 1 ) = chipol
end if
end if
end if
end if
2022-05-10 12:36:37 +02:00
2019-03-26 15:21:22 +01:00
iopmin = min ( iopmin , iop ( jk ) )
2022-05-10 12:36:37 +02:00
if ( ip < params % raytracing % ipass ) then ! not last pass
2019-03-28 10:50:28 +01:00
yw0 ( : , jk ) = yw ( : , jk ) ! * store current coordinates in case
ypw0 ( : , jk ) = ypw ( : , jk ) ! current pass ends on next step
2019-03-26 15:21:22 +01:00
end if
2022-05-10 12:36:37 +02:00
! Compute ECRH&CD if (inside plasma & power available>0 & ray still active)
! Note: power check is τ + τ₀ + log(coupling O/X) < τ_critic
if ( ierrn == 0 . and . params % ecrh_cd % iwarm > 0 . and . ins_pl . and . &
2019-03-26 15:21:22 +01:00
( tau1 ( jk ) + tau0 ( jk ) + lgcpl1 ( jk ) ) < = taucr . and . . not . iwait ( jk ) ) then ! H&CD computation check
tekev = temp ( psinv )
2022-05-20 13:10:49 +02:00
block
complex ( wp_ ) :: Npr_warm
call alpha_effj ( params % ecrh_cd , psinv , xg , yg , dens , tekev , &
ak0 , bres , derdnm , anpl , anpr , sox , Npr_warm , &
alpha , didp , nharm , nhf , iokhawa , ierrhcd )
anprre = Npr_warm % re
anprim = Npr_warm % im
2022-07-14 00:38:34 +02:00
if ( ierrhcd / = 0 ) then
error = ior ( error , ierrhcd )
call print_err_ecrh_cd ( ierrhcd , i , Npr_warm , alpha )
end if
2022-05-20 13:10:49 +02:00
end block
2019-03-26 15:21:22 +01:00
else
tekev = zero
alpha = zero
didp = zero
anprim = zero
anprre = anpr
nharm = 0
nhf = 0
iokhawa = 0
end if
if ( nharm > 0 ) iiv ( jk ) = i
2022-05-10 12:36:37 +02:00
2019-03-26 15:21:22 +01:00
psjki ( jk , i ) = psinv
2022-05-10 12:36:37 +02:00
! Computation of the ray τ, dP/ds, P(s), dI/ds, I(s)
! optical depth: τ = ∫α(s)ds using the trapezoid rule
tau = tau0 ( jk ) + 0.5_wp_ * ( alphaabs0 ( jk ) + alpha ) * dersdst * params % raytracing % dst
pow = p0ray ( jk ) * exp ( - tau ) ! residual power: P = P₀exp(-τ)
ppabs ( jk , i ) = p0ray ( jk ) - pow ! absorbed power: P_abs = P₀ - P
2023-03-29 11:19:00 +02:00
dids = didp * pow * alpha ! current driven: dI/ds = dI/dP⋅dP/ds = dI/dP⋅P⋅α
2022-05-10 12:36:37 +02:00
! current: I = ∫dI/ds⋅ds using the trapezoid rule
ccci ( jk , i ) = ccci0 ( jk ) + 0.5_wp_ * ( dids0 ( jk ) + dids ) * dersdst * params % raytracing % dst
tau0 ( jk ) = tau
alphaabs0 ( jk ) = alpha
dids0 ( jk ) = dids
ccci0 ( jk ) = ccci ( jk , i )
2019-03-28 10:50:28 +01:00
if ( iwait ( jk ) ) then ! copy values from last pass for inactive ray
2022-05-10 12:36:37 +02:00
ppabs ( jk , i : params % raytracing % nstep ) = ppabs ( jk , i - 1 )
ccci ( jk , i : params % raytracing % nstep ) = ccci ( jk , i - 1 )
psjki ( jk , i : params % raytracing % nstep ) = psjki ( jk , i - 1 )
2019-03-26 15:21:22 +01:00
else
2019-12-09 15:49:37 +01:00
call print_output ( i , jk , stv ( jk ) , p0ray ( jk ) , xv , psinv , &
2019-03-26 15:21:22 +01:00
btot , bv , ak0 , anpl , anpr , anv , anprim , dens , tekev , alpha , tau , dids , &
2019-03-28 10:50:28 +01:00
nharm , nhf , iokhawa , index_rt , ddr , ddi , xg , yg , derxg ) ! p0ray/etau1 [dids normalization] = fraction of p0 coupled to this ray (not including absorption from previous passes)
2019-03-26 15:21:22 +01:00
end if
2022-05-10 12:36:37 +02:00
end do rays_loop
! ============ rays loop END ============
if ( i == params % raytracing % nstep ) then ! step limit reached?
do jk = 1 , params % raytracing % nray
2019-03-26 15:21:22 +01:00
if ( iop ( jk ) < 3 ) call turnoffray ( jk , ip , ib , iroff ) ! * ray hasn't exited+reentered the plasma by last step => stop ray
end do
end if
2022-05-10 12:36:37 +02:00
! print ray positions for j=nrayr in local reference system
if ( mod ( i , params % output % istpr ) == 0 ) then
if ( params % raytracing % nray > 1 . and . all ( . not . iwait ) ) &
call print_projxyzt ( stv , yw , 0 )
2019-03-26 15:21:22 +01:00
end if
2022-05-10 12:36:37 +02:00
! test whether further trajectory integration is unnecessary
call vmaxmin ( tau1 + tau0 + lgcpl1 , params % raytracing % nray , taumn , taumx ) ! test on tau + coupling
2022-07-14 00:38:34 +02:00
if ( is_critical ( error ) ) then ! stop propagation for current beam
2019-03-26 15:21:22 +01:00
istop_pass = istop_pass + 1 ! * +1 non propagating beam
2022-05-10 12:36:37 +02:00
if ( ip < params % raytracing % ipass ) call turnoffray ( 0 , ip , ib , iroff ) ! * de-activate derived beams
2019-03-26 15:21:22 +01:00
exit
2019-03-28 10:50:28 +01:00
else if ( all ( iwait ) ) then ! all rays in current beam are waiting for next pass => do not increase istop_pass
2019-03-26 15:21:22 +01:00
exit
2015-11-19 19:20:58 +01:00
end if
2022-05-10 12:36:37 +02:00
end do propagation_loop
2021-12-18 18:57:38 +01:00
call log_debug ( ' propagation loop end' , mod = 'gray_core' , proc = 'gray_main' )
2022-05-10 12:36:37 +02:00
! ======== propagation loop END ========
! print all ray positions in local reference system
if ( params % raytracing % nray > 1 . and . all ( . not . iwait ) ) call print_projxyzt ( stv , yw , 1 )
! =========== post-proc BEGIN ===========
! compute total absorbed power and driven current for current beam
if ( i > params % raytracing % nstep ) i = params % raytracing % nstep
2019-03-26 15:21:22 +01:00
pabs_beam = sum ( ppabs ( : , i ) )
icd_beam = sum ( ccci ( : , i ) )
2022-05-10 12:36:37 +02:00
call vmaxmin ( tau0 , params % raytracing % nray , taumn , taumx ) ! taumn,taumx for print
! compute power and current density profiles for all rays
2019-03-26 15:21:22 +01:00
call spec ( psjki , ppabs , ccci , iiv , pabs_beam , icd_beam , dpdv_beam , jphi_beam , jcd_beam , &
pins_beam , currins_beam )
2022-05-10 12:36:37 +02:00
2019-03-26 15:21:22 +01:00
pabs_pass ( iox ) = pabs_pass ( iox ) + pabs_beam ! 0D results for current pass, sum on O/X mode beams
icd_pass ( iox ) = icd_pass ( iox ) + icd_beam
2022-05-10 12:36:37 +02:00
if ( ip < params % raytracing % ipass . and . iopmin > 2 ) then ! not last pass AND at least one ray re-entered plasma
cpl_beam1 = sum ( p0ray * exp ( - tau0 ) * cpls ( : , iO ) / cpl1 , MASK = iop > 2 ) / &
sum ( p0ray * exp ( - tau0 ) , MASK = iop > 2 ) ! * average O-mode coupling for next beam (on active rays)
2019-03-28 12:33:43 +01:00
cpl_beam2 = one - cpl_beam1 ! * average X-mode coupling for next beam
2022-05-10 12:36:37 +02:00
if ( iop ( 1 ) > 2 ) then ! * central ray O/X-mode coupling for next beam
2019-03-28 12:33:43 +01:00
cpl_cbeam1 = cpls ( 1 , iO ) / cpl1 ( 1 )
cpl_cbeam2 = one - cpl_cbeam1
end if
else ! last pass OR no ray re-entered plasma
2019-03-26 15:21:22 +01:00
cpl_beam1 = zero
cpl_beam2 = zero
2015-11-18 17:34:33 +01:00
end if
2021-12-18 18:57:38 +01:00
! print final results for pass on screen
call log_info ( ' partial results:' , mod = 'gray_core' , proc = 'gray_main' )
2022-05-26 02:33:55 +02:00
write ( msg , '(" final step:", x,a,g0, x,a,g0.4)' ) 'n=' , i , '(s, ct, Sr)=' , stv ( 1 )
2021-12-18 18:57:38 +01:00
call log_info ( msg , mod = 'gray_core' , proc = 'gray_main' )
write ( msg , '(3x,a,2(x,a,"=",g0.4))' ) 'optical depth:' , 'τ_min' , taumn , 'τ_max' , taumx
call log_info ( msg , mod = 'gray_core' , proc = 'gray_main' )
write ( msg , '(3x,a,g0.3," MW")' ) 'absoption: P=' , pabs_beam
call log_info ( msg , mod = 'gray_core' , proc = 'gray_main' )
2023-10-04 16:56:05 +02:00
write ( msg , '(3x,a,g0.3," kA")' ) 'current drive: I=' , icd_beam * 1.0e3_wp_
2021-12-18 18:57:38 +01:00
call log_info ( msg , mod = 'gray_core' , proc = 'gray_main' )
2022-05-10 12:36:37 +02:00
if ( ip < params % raytracing % ipass ) then
2021-12-18 18:57:38 +01:00
write ( msg , '(3x,a,(g0.4,", ",g0.4))' ) & ! average coupling for next O/X beams (=0 if no ray re-entered plasma)
'next couplings [O,X mode]: c=' , cpl_beam1 , cpl_beam2
call log_info ( msg , mod = 'gray_core' , proc = 'gray_main' )
if ( iop ( 1 ) > 2 ) then
write ( msg , '(3x,a,(g0.4,", ",g0.4))' ) &
'coupling [ctr ray, O/X]:' , cpl_cbeam1 , cpl_cbeam2 ! central ray coupling for next O/X beams
end if
2019-03-26 15:21:22 +01:00
end if
2021-12-18 18:57:38 +01:00
2019-03-26 15:21:22 +01:00
call print_pec ( rhop_tab , rhot_tab , jphi_beam , jcd_beam , dpdv_beam , currins_beam , &
pins_beam , ip ) ! *print power and current density profiles for current beam
2022-05-10 12:36:37 +02:00
2019-03-26 15:21:22 +01:00
call postproc_profiles ( pabs_beam , icd_beam , rhot_tab , dpdv_beam , jphi_beam , &
rhotpav , drhotpav , rhotjava , drhotjava , dpdvp , jphip , rhotp , drhotp , rhotj , &
drhotj , dpdvmx , jphimx , ratjamx , ratjbmx ) ! *compute profiles width for current beam
2022-05-10 12:36:37 +02:00
2019-03-26 15:21:22 +01:00
call print_finals ( pabs_beam , icd_beam , dpdvp , jphip , rhotpav , rhotjava , &
drhotpav , drhotjava , dpdvmx , jphimx , rhotp , rhotj , drhotp , drhotj , ratjamx , &
ratjbmx , stv ( 1 ) , psipv ( index_rt ) , chipv ( index_rt ) , index_rt , sum ( p0ray ) , &
cpl_beam1 , cpl_beam2 ) ! *print 0D results for current beam
2022-05-10 12:36:37 +02:00
! ============ post-proc END ============
end do beam_loop
2021-12-18 18:57:38 +01:00
call log_debug ( 'beam loop end' , mod = 'gray_core' , proc = 'gray_main' )
2022-05-10 12:36:37 +02:00
! ============ beam loop END ============
! ======= cumulative prints BEGIN =======
2021-12-15 02:31:09 +01:00
results % pabs = results % pabs + sum ( pabs_pass ) ! *final results (O+X) [gray_main output]
results % icd = results % icd + sum ( icd_pass )
2021-12-18 18:57:38 +01:00
! print final results for pass on screen
call log_info ( ' comulative results:' , mod = 'gray_core' , proc = 'gray_main' )
write ( msg , '(" absoption [O,X mode] P=",g0.4,", ",g0.4," MW")' ) &
pabs_pass ( 1 ) , pabs_pass ( 2 )
call log_info ( msg , mod = 'gray_core' , proc = 'gray_main' )
write ( msg , '(" current drive [O,X mode] I=",g0.4,", ",g0.4," kA")' ) &
icd_pass ( 1 ) * 1.0e3_wp_ , icd_pass ( 2 ) * 1.0e3_wp_
call log_info ( msg , mod = 'gray_core' , proc = 'gray_main' )
2022-05-10 12:36:37 +02:00
! ======== cumulative prints END ========
2019-03-26 15:21:22 +01:00
if ( istop_pass == nbeam_pass ) exit ! no active beams
2022-05-10 12:36:37 +02:00
end do main_loop
2021-12-18 18:57:38 +01:00
call log_debug ( 'pass loop end' , mod = 'gray_core' , proc = 'gray_main' )
2022-05-10 12:36:37 +02:00
! ============ main loop END ============
2019-03-26 15:21:22 +01:00
2022-05-10 12:36:37 +02:00
! ========== free memory BEGIN ==========
2016-06-01 15:49:35 +02:00
call dealloc_surfvec
2019-03-26 15:21:22 +01:00
call dealloc_beam ( yw , ypw , xc , du1 , gri , ggri , psjki , ppabs , ccci , tau0 , &
alphaabs0 , dids0 , ccci0 , p0jk , ext , eyt , iiv )
2016-06-01 15:49:35 +02:00
call dealloc_pec
2019-03-26 15:21:22 +01:00
call dealloc_multipass ( iwait , iroff , iop , iow , yynext , yypnext , yw0 , ypw0 , &
stnext , stv , p0ray , taus , tau1 , etau1 , cpls , cpl1 , lgcpl1 , jphi_beam , &
pins_beam , currins_beam , dpdv_beam , jcd_beam , psipv , chipv )
2022-05-10 12:36:37 +02:00
! =========== free memory END ===========
2015-11-23 18:55:27 +01:00
end subroutine gray_main
2015-11-18 17:34:33 +01:00
2021-12-15 02:30:55 +01:00
2015-11-19 19:20:58 +01:00
subroutine vectinit ( psjki , ppabs , ccci , tau0 , alphaabs0 , dids0 , ccci0 , iiv )
2021-12-19 16:39:19 +01:00
use const_and_precisions , only : zero
2015-11-18 17:34:33 +01:00
implicit none
! arguments
2015-11-19 19:20:58 +01:00
real ( wp_ ) , dimension ( : , : ) , intent ( out ) :: psjki , ppabs , ccci
real ( wp_ ) , dimension ( : ) , intent ( out ) :: tau0 , alphaabs0 , dids0 , ccci0
2015-11-18 17:34:33 +01:00
integer , dimension ( : ) , intent ( out ) :: iiv
!! common/external functions/variables
! integer :: jclosest
! real(wp_), dimension(3) :: anwcl,xwcl
!
! common/refln/anwcl,xwcl,jclosest
!
! jclosest=nrayr+1
! anwcl(1:3)=0.0_wp_
! xwcl(1:3)=0.0_wp_
2015-11-19 19:20:58 +01:00
psjki = zero
ppabs = zero
ccci = zero
tau0 = zero
alphaabs0 = zero
dids0 = zero
ccci0 = zero
iiv = 1
2015-11-18 17:34:33 +01:00
end subroutine vectinit
2015-11-23 18:55:27 +01:00
2022-05-10 12:36:37 +02:00
subroutine ic_gb ( params , anv0c , ak0 , ywrk0 , ypwrk0 , &
xc0 , du10 , gri , ggri , index_rt )
! beam tracing initial conditions igrad=1
! !!!!!! check ray tracing initial conditions igrad=0 !!!!!!
2021-12-19 16:39:19 +01:00
use const_and_precisions , only : zero , one , pi , half , two , degree , ui = > im
2022-05-10 12:36:37 +02:00
use gray_params , only : gray_parameters
use beamdata , only : nray , nrayr , nrayth , rwmax
2015-11-18 17:34:33 +01:00
implicit none
2022-05-10 12:36:37 +02:00
! subroutine arguments
type ( gray_parameters ) , intent ( in ) :: params
real ( wp_ ) , dimension ( 3 ) , intent ( in ) :: anv0c
real ( wp_ ) , intent ( in ) :: ak0
real ( wp_ ) , dimension ( 6 , nray ) , intent ( out ) :: ywrk0 , ypwrk0
real ( wp_ ) , dimension ( 3 , nrayth , nrayr ) , intent ( out ) :: xc0 , du10
real ( wp_ ) , dimension ( 3 , nray ) , intent ( out ) :: gri
real ( wp_ ) , dimension ( 3 , 3 , nray ) , intent ( out ) :: ggri
integer , intent ( in ) :: index_rt
! local variables
real ( wp_ ) , dimension ( 3 ) :: xv0c
real ( wp_ ) :: wcsi , weta , rcicsi , rcieta , phiw , phir
2015-11-18 17:34:33 +01:00
integer :: j , k , jk
real ( wp_ ) :: csth , snth , csps , snps , phiwrad , phirrad , csphiw , snphiw , alfak
real ( wp_ ) :: wwcsi , wweta , sk , sw , dk , dw , rci1 , ww1 , rci2 , ww2 , wwxx , wwyy , wwxy
real ( wp_ ) :: rcixx , rciyy , rcixy , dwwxx , dwwyy , dwwxy , d2wwxx , d2wwyy , d2wwxy
real ( wp_ ) :: drcixx , drciyy , drcixy , dr , da , ddfu , dcsiw , detaw , dx0t , dy0t
real ( wp_ ) :: x0t , y0t , z0t , dx0 , dy0 , dz0 , x0 , y0 , z0 , gxt , gyt , gzt , gr2
real ( wp_ ) :: gxxt , gyyt , gzzt , gxyt , gxzt , gyzt , dgr2xt , dgr2yt , dgr2zt
real ( wp_ ) :: dgr2x , dgr2y , dgr2z , pppx , pppy , denpp , ppx , ppy
real ( wp_ ) :: anzt , anxt , anyt , anx , any , anz , an20 , an0
real ( wp_ ) :: du1tx , du1ty , du1tz , denom , ddr , ddi
real ( wp_ ) , dimension ( nrayr ) :: uj
real ( wp_ ) , dimension ( nrayth ) :: sna , csa
complex ( wp_ ) :: sss , ddd , phic , qi1 , qi2 , tc , ts , qqxx , qqxy , qqyy , dqi1 , dqi2
complex ( wp_ ) :: dqqxx , dqqyy , dqqxy , d2qi1 , d2qi2 , d2qqxx , d2qqyy , d2qqxy
2022-05-10 12:36:37 +02:00
! initial beam parameters
xv0c = params % antenna % pos
wcsi = params % antenna % w ( 1 )
weta = params % antenna % w ( 2 )
rcicsi = params % antenna % ri ( 1 )
rcieta = params % antenna % ri ( 2 )
phiw = params % antenna % phi ( 1 )
phir = params % antenna % phi ( 2 )
2015-11-18 17:34:33 +01:00
csth = anv0c ( 3 )
snth = sqrt ( one - csth ** 2 )
if ( snth > zero ) then
csps = anv0c ( 2 ) / snth
snps = anv0c ( 1 ) / snth
else
csps = one
snps = zero
end if
2022-05-10 12:36:37 +02:00
! Gaussian beam: exp[-ik0 zt] exp[-i k0/2 S(xt,yt,zt)]
! xt,yt,zt, cartesian coordinate system with zt along the beamline and xt in the z = 0 plane
! S(xt,yt,zt) = S_real +i S_imag = Qxx(zt) xt^2 + Qyy(zt) yt^2 + 2 Qxy(zt) xt yt
! (csiw, etaw) and (csiR, etaR) intensity and phase ellipse, rotated by angle phiw and phiR
! S(xt,yt,zt) = csiR^2 / Rccsi +etaR^2 /Rceta - i (csiw^2 Wcsi +etaw^2 Weta)
! Rccsi,eta curvature radius at the launching point
! Wcsi,eta =2/(k0 wcsi,eta^2) with wcsi,eta^2 beam size at the launching point
2015-11-19 18:44:17 +01:00
2015-11-18 17:34:33 +01:00
phiwrad = phiw * degree
phirrad = phir * degree
csphiw = cos ( phiwrad )
snphiw = sin ( phiwrad )
2022-05-10 12:36:37 +02:00
!csphir = cos(phirrad)
!snphir = sin(phirrad)
2015-11-18 17:34:33 +01:00
wwcsi = two / ( ak0 * wcsi ** 2 )
wweta = two / ( ak0 * weta ** 2 )
if ( phir / = phiw ) then
2015-11-19 18:44:17 +01:00
sk = rcicsi + rcieta
sw = wwcsi + wweta
dk = rcicsi - rcieta
dw = wwcsi - wweta
ts = - ( dk * sin ( 2 * phirrad ) - ui * dw * sin ( 2 * phiwrad ) )
tc = ( dk * cos ( 2 * phirrad ) - ui * dw * cos ( 2 * phiwrad ) )
2022-07-15 00:09:36 +02:00
phic = half * atan ( ts / tc )
2015-11-19 18:44:17 +01:00
ddd = dk * cos ( 2 * ( phirrad + phic ) ) - ui * dw * cos ( 2 * ( phiwrad + phic ) )
sss = sk - ui * sw
qi1 = half * ( sss + ddd )
qi2 = half * ( sss - ddd )
2023-09-20 16:03:08 +02:00
rci1 = real ( qi1 )
rci2 = real ( qi2 )
ww1 = - imag ( qi1 )
ww2 = - imag ( qi2 )
2015-11-18 17:34:33 +01:00
else
rci1 = rcicsi
rci2 = rcieta
ww1 = wwcsi
ww2 = wweta
phic = - phiwrad
qi1 = rci1 - ui * ww1
qi2 = rci2 - ui * ww2
end if
2022-05-10 12:36:37 +02:00
!w01=sqrt(2.0_wp_/(ak0*ww1))
!d01=-rci1/(rci1**2+ww1**2)
!w02=sqrt(2.0_wp_/(ak0*ww2))
!d02=-rci2/(rci2**2+ww2**2)
2015-11-18 17:34:33 +01:00
qqxx = qi1 * cos ( phic ) ** 2 + qi2 * sin ( phic ) ** 2
qqyy = qi1 * sin ( phic ) ** 2 + qi2 * cos ( phic ) ** 2
2015-11-19 18:44:17 +01:00
qqxy = - ( qi1 - qi2 ) * sin ( phic ) * cos ( phic )
2023-09-20 16:03:08 +02:00
wwxx = - imag ( qqxx )
wwyy = - imag ( qqyy )
wwxy = - imag ( qqxy )
rcixx = real ( qqxx )
rciyy = real ( qqyy )
rcixy = real ( qqxy )
2015-11-18 17:34:33 +01:00
dqi1 = - qi1 ** 2
dqi2 = - qi2 ** 2
d2qi1 = 2 * qi1 ** 3
d2qi2 = 2 * qi2 ** 3
dqqxx = dqi1 * cos ( phic ) ** 2 + dqi2 * sin ( phic ) ** 2
dqqyy = dqi1 * sin ( phic ) ** 2 + dqi2 * cos ( phic ) ** 2
2015-11-19 18:44:17 +01:00
dqqxy = - ( dqi1 - dqi2 ) * sin ( phic ) * cos ( phic )
2015-11-18 17:34:33 +01:00
d2qqxx = d2qi1 * cos ( phic ) ** 2 + d2qi2 * sin ( phic ) ** 2
d2qqyy = d2qi1 * sin ( phic ) ** 2 + d2qi2 * cos ( phic ) ** 2
2015-11-19 18:44:17 +01:00
d2qqxy = - ( d2qi1 - d2qi2 ) * sin ( phic ) * cos ( phic )
2023-09-20 16:03:08 +02:00
dwwxx = - imag ( dqqxx )
dwwyy = - imag ( dqqyy )
dwwxy = - imag ( dqqxy )
d2wwxx = - imag ( d2qqxx )
d2wwyy = - imag ( d2qqyy )
d2wwxy = - imag ( d2qqxy )
drcixx = real ( dqqxx )
drciyy = real ( dqqyy )
drcixy = real ( dqqxy )
2015-11-18 17:34:33 +01:00
if ( nrayr > 1 ) then
2023-09-20 16:03:08 +02:00
dr = rwmax / ( nrayr - 1 )
2015-11-18 17:34:33 +01:00
else
dr = one
end if
2023-09-20 16:03:08 +02:00
ddfu = two * dr ** 2 / ak0 ! twodr2 = 2*dr**2 = 2*rwmax/real(nrayr-1)
2015-11-18 17:34:33 +01:00
do j = 1 , nrayr
2023-09-20 16:03:08 +02:00
uj ( j ) = real ( j - 1 , wp_ )
2022-05-10 12:36:37 +02:00
end do
2015-11-18 17:34:33 +01:00
2023-09-20 16:03:08 +02:00
da = 2 * pi / nrayth
2015-11-18 17:34:33 +01:00
do k = 1 , nrayth
alfak = ( k - 1 ) * da
sna ( k ) = sin ( alfak )
csa ( k ) = cos ( alfak )
end do
2022-05-10 12:36:37 +02:00
! central ray
2015-11-18 17:34:33 +01:00
jk = 1
gri ( : , 1 ) = zero
ggri ( : , : , 1 ) = zero
ywrk0 ( 1 : 3 , 1 ) = xv0c
ywrk0 ( 4 : 6 , 1 ) = anv0c
ypwrk0 ( 1 : 3 , 1 ) = anv0c
ypwrk0 ( 4 : 6 , 1 ) = zero
2022-05-10 12:36:37 +02:00
2015-11-18 17:34:33 +01:00
do k = 1 , nrayth
dcsiw = dr * csa ( k ) * wcsi
detaw = dr * sna ( k ) * weta
dx0t = dcsiw * csphiw - detaw * snphiw
dy0t = dcsiw * snphiw + detaw * csphiw
du1tx = ( dx0t * wwxx + dy0t * wwxy ) / ddfu
du1ty = ( dx0t * wwxy + dy0t * wwyy ) / ddfu
xc0 ( : , k , 1 ) = xv0c
du10 ( 1 , k , 1 ) = du1tx * csps + snps * du1ty * csth
du10 ( 2 , k , 1 ) = - du1tx * snps + csps * du1ty * csth
du10 ( 3 , k , 1 ) = - du1ty * snth
end do
ddr = zero
ddi = zero
2022-05-10 12:36:37 +02:00
! loop on rays jk>1
2015-11-18 17:34:33 +01:00
j = 2
k = 0
do jk = 2 , nray
k = k + 1
if ( k > nrayth ) then
j = j + 1
k = 1
end if
2022-05-10 12:36:37 +02:00
!csiw=u*dcsiw
!etaw=u*detaw
!csir=x0t*csphir+y0t*snphir
!etar=-x0t*snphir+y0t*csphir
2015-11-18 17:34:33 +01:00
dcsiw = dr * csa ( k ) * wcsi
detaw = dr * sna ( k ) * weta
dx0t = dcsiw * csphiw - detaw * snphiw
dy0t = dcsiw * snphiw + detaw * csphiw
x0t = uj ( j ) * dx0t
y0t = uj ( j ) * dy0t
2015-11-23 18:55:27 +01:00
z0t = - ( half * ( rcixx * x0t ** 2 + rciyy * y0t ** 2 ) + rcixy * x0t * y0t )
2015-11-18 17:34:33 +01:00
dx0 = x0t * csps + snps * ( y0t * csth + z0t * snth )
dy0 = - x0t * snps + csps * ( y0t * csth + z0t * snth )
dz0 = z0t * csth - y0t * snth
x0 = xv0c ( 1 ) + dx0
y0 = xv0c ( 2 ) + dy0
z0 = xv0c ( 3 ) + dz0
gxt = x0t * wwxx + y0t * wwxy
gyt = x0t * wwxy + y0t * wwyy
gzt = half * ( x0t ** 2 * dwwxx + y0t ** 2 * dwwyy ) + x0t * y0t * dwwxy
gr2 = gxt * gxt + gyt * gyt + gzt * gzt
gxxt = wwxx
gyyt = wwyy
gzzt = half * ( x0t ** 2 * d2wwxx + y0t ** 2 * d2wwyy ) + x0t * y0t * d2wwxy
gxyt = wwxy
gxzt = x0t * dwwxx + y0t * dwwxy
gyzt = x0t * dwwxy + y0t * dwwyy
dgr2xt = 2 * ( gxt * gxxt + gyt * gxyt + gzt * gxzt )
dgr2yt = 2 * ( gxt * gxyt + gyt * gyyt + gzt * gyzt )
dgr2zt = 2 * ( gxt * gxzt + gyt * gyzt + gzt * gzzt )
dgr2x = dgr2xt * csps + snps * ( dgr2yt * csth + dgr2zt * snth )
dgr2y = - dgr2xt * snps + csps * ( dgr2yt * csth + dgr2zt * snth )
dgr2z = dgr2zt * csth - dgr2yt * snth
2022-05-10 12:36:37 +02:00
2015-11-18 17:34:33 +01:00
gri ( 1 , jk ) = gxt * csps + snps * ( gyt * csth + gzt * snth )
gri ( 2 , jk ) = - gxt * snps + csps * ( gyt * csth + gzt * snth )
gri ( 3 , jk ) = gzt * csth - gyt * snth
ggri ( 1 , 1 , jk ) = gxxt * csps ** 2 &
+ snps ** 2 * ( gyyt * csth ** 2 + gzzt * snth ** 2 + 2 * snth * csth * gyzt ) &
+ 2 * snps * csps * ( gxyt * csth + gxzt * snth )
ggri ( 2 , 1 , jk ) = csps * snps &
* ( - gxxt + csth ** 2 * gyyt + snth ** 2 * gzzt + 2 * csth * snth * gyzt ) &
+ ( csps ** 2 - snps ** 2 ) * ( snth * gxzt + csth * gxyt )
ggri ( 3 , 1 , jk ) = csth * snth * snps * ( gzzt - gyyt ) + ( csth ** 2 - snth ** 2 ) &
* snps * gyzt + csps * ( csth * gxzt - snth * gxyt )
ggri ( 1 , 2 , jk ) = ggri ( 2 , 1 , jk )
ggri ( 2 , 2 , jk ) = gxxt * snps ** 2 &
+ csps ** 2 * ( gyyt * csth ** 2 + gzzt * snth ** 2 + 2 * snth * csth * gyzt ) &
- 2 * snps * csps * ( gxyt * csth + gxzt * snth )
ggri ( 3 , 2 , jk ) = csth * snth * csps * ( gzzt - gyyt ) + ( csth ** 2 - snth ** 2 ) &
* csps * gyzt + snps * ( snth * gxyt - csth * gxzt )
ggri ( 1 , 3 , jk ) = ggri ( 3 , 1 , jk )
ggri ( 2 , 3 , jk ) = ggri ( 3 , 2 , jk )
ggri ( 3 , 3 , jk ) = gzzt * csth ** 2 + gyyt * snth ** 2 - 2 * csth * snth * gyzt
du1tx = ( dx0t * wwxx + dy0t * wwxy ) / ddfu
du1ty = ( dx0t * wwxy + dy0t * wwyy ) / ddfu
du1tz = half * uj ( j ) * ( dx0t ** 2 * dwwxx + dy0t ** 2 * dwwyy + 2 * dx0t * dy0t * dwwxy ) / ddfu
du10 ( 1 , k , j ) = du1tx * csps + snps * ( du1ty * csth + du1tz * snth )
du10 ( 2 , k , j ) = - du1tx * snps + csps * ( du1ty * csth + du1tz * snth )
du10 ( 3 , k , j ) = du1tz * csth - du1ty * snth
pppx = x0t * rcixx + y0t * rcixy
pppy = x0t * rcixy + y0t * rciyy
denpp = pppx * gxt + pppy * gyt
if ( denpp / = zero ) then
ppx = - pppx * gzt / denpp
ppy = - pppy * gzt / denpp
else
ppx = zero
ppy = zero
end if
anzt = sqrt ( ( one + gr2 ) / ( one + ppx ** 2 + ppy ** 2 ) )
anxt = ppx * anzt
anyt = ppy * anzt
anx = anxt * csps + snps * ( anyt * csth + anzt * snth )
any = - anxt * snps + csps * ( anyt * csth + anzt * snth )
anz = anzt * csth - anyt * snth
an20 = one + gr2
an0 = sqrt ( an20 )
xc0 ( 1 , k , j ) = x0
xc0 ( 2 , k , j ) = y0
xc0 ( 3 , k , j ) = z0
ywrk0 ( 1 , jk ) = x0
ywrk0 ( 2 , jk ) = y0
ywrk0 ( 3 , jk ) = z0
ywrk0 ( 4 , jk ) = anx
ywrk0 ( 5 , jk ) = any
ywrk0 ( 6 , jk ) = anz
2022-05-10 12:36:37 +02:00
! integration variable
select case ( params % raytracing % idst )
case ( 0 )
! optical path: s
denom = an0
case ( 1 )
! "time": c⋅t
denom = one
case ( 2 )
! real eikonal: S_R
denom = an20
2015-11-18 17:34:33 +01:00
end select
ypwrk0 ( 1 , jk ) = anx / denom
ypwrk0 ( 2 , jk ) = any / denom
ypwrk0 ( 3 , jk ) = anz / denom
ypwrk0 ( 4 , jk ) = dgr2x / ( 2 * denom )
ypwrk0 ( 5 , jk ) = dgr2y / ( 2 * denom )
ypwrk0 ( 6 , jk ) = dgr2z / ( 2 * denom )
ddr = anx ** 2 + any ** 2 + anz ** 2 - an20
ddi = 2 * ( anxt * gxt + anyt * gyt + anzt * gzt )
2022-05-10 12:36:37 +02:00
! st=0, index_rt=1, B=0, N=0, psin=-1, Xg=0, Yg=0, gradXg=0
call print_output ( 0 , jk , zero , one , xc0 ( : , k , j ) , - one , zero , [ zero , zero , zero ] , &
ak0 , zero , zero , [ zero , zero , zero ] , zero , zero , zero , zero , zero , zero , &
0 , 0 , 0 , index_rt , ddr , ddi , zero , zero , [ zero , zero , zero ] )
2015-11-18 17:34:33 +01:00
end do
end subroutine ic_gb
2015-11-23 18:55:27 +01:00
2019-03-26 15:21:22 +01:00
subroutine rkstep ( sox , bres , xgcn , y , yp , dgr , ddgr , igrad )
2015-11-18 17:34:33 +01:00
! Runge-Kutta integrator
! use gray_params, only : igrad
use beamdata , only : h , hh , h6
implicit none
2022-05-20 13:10:49 +02:00
real ( wp_ ) , intent ( in ) :: bres , xgcn
2015-11-18 17:34:33 +01:00
real ( wp_ ) , dimension ( 6 ) , intent ( inout ) :: y
real ( wp_ ) , dimension ( 6 ) , intent ( in ) :: yp
real ( wp_ ) , dimension ( 3 ) , intent ( in ) :: dgr
real ( wp_ ) , dimension ( 3 , 3 ) , intent ( in ) :: ddgr
2022-05-20 13:10:49 +02:00
integer , intent ( in ) :: igrad , sox
2015-11-18 17:34:33 +01:00
real ( wp_ ) , dimension ( 6 ) :: yy , fk1 , fk2 , fk3 , fk4
2023-09-21 13:59:21 +02:00
fk1 = yp
2015-11-18 17:34:33 +01:00
yy = y + fk1 * hh
2023-09-21 13:59:21 +02:00
call rhs ( sox , bres , xgcn , yy , dgr , ddgr , fk2 , igrad )
2015-11-18 17:34:33 +01:00
yy = y + fk2 * hh
2023-09-21 13:59:21 +02:00
call rhs ( sox , bres , xgcn , yy , dgr , ddgr , fk3 , igrad )
2015-11-18 17:34:33 +01:00
yy = y + fk3 * h
2023-09-21 13:59:21 +02:00
call rhs ( sox , bres , xgcn , yy , dgr , ddgr , fk4 , igrad )
2015-11-18 17:34:33 +01:00
y = y + h6 * ( fk1 + 2 * fk2 + 2 * fk3 + fk4 )
end subroutine rkstep
2015-11-23 18:55:27 +01:00
2023-09-21 13:59:21 +02:00
subroutine rhs ( sox , bres , xgcn , y , dgr , ddgr , dery , igrad )
2015-11-18 17:34:33 +01:00
! Compute right-hand side terms of the ray equations (dery)
! used in R-K integrator
implicit none
! arguments
real ( wp_ ) , dimension ( 6 ) , intent ( in ) :: y
2023-09-21 13:59:21 +02:00
real ( wp_ ) , intent ( in ) :: bres , xgcn
real ( wp_ ) , dimension ( 3 ) , intent ( in ) :: dgr
2015-11-18 17:34:33 +01:00
real ( wp_ ) , dimension ( 3 , 3 ) , intent ( in ) :: ddgr
real ( wp_ ) , dimension ( 6 ) , intent ( out ) :: dery
2022-05-20 13:10:49 +02:00
integer , intent ( in ) :: igrad , sox
2015-11-18 17:34:33 +01:00
! local variables
2023-09-21 13:59:21 +02:00
real ( wp_ ) :: dens , btot , xg , yg
2015-11-18 17:34:33 +01:00
real ( wp_ ) , dimension ( 3 ) :: xv , anv , bv , derxg , deryg
real ( wp_ ) , dimension ( 3 , 3 ) :: derbv
xv = y ( 1 : 3 )
2023-03-29 22:13:43 +02:00
call plas_deriv ( xv , bres , xgcn , dens , btot , bv , derbv , xg , yg , derxg , deryg )
2015-11-18 17:34:33 +01:00
anv = y ( 4 : 6 )
2023-03-29 22:13:43 +02:00
call disp_deriv ( anv , sox , xg , yg , derxg , deryg , bv , derbv , dgr , ddgr , igrad , dery )
2015-11-18 17:34:33 +01:00
end subroutine rhs
2015-11-23 18:55:27 +01:00
2015-11-18 17:34:33 +01:00
subroutine ywppla_upd ( xv , anv , dgr , ddgr , sox , bres , xgcn , dery , psinv , dens , btot , &
2021-12-15 02:31:09 +01:00
bv , xg , yg , derxg , anpl , anpr , ddr , ddi , dersdst , derdnm , error , igrad )
2015-11-18 17:34:33 +01:00
! Compute right-hand side terms of the ray equations (dery)
! used after full R-K step and grad(S_I) update
2022-07-14 00:38:34 +02:00
use gray_errors , only : raise_error , large_npl
2015-11-18 17:34:33 +01:00
implicit none
! arguments
real ( wp_ ) , dimension ( 3 ) , intent ( in ) :: xv , anv
real ( wp_ ) , dimension ( 3 ) , intent ( in ) :: dgr
real ( wp_ ) , dimension ( 3 , 3 ) , intent ( in ) :: ddgr
2022-05-20 13:10:49 +02:00
integer , intent ( in ) :: sox
real ( wp_ ) , intent ( in ) :: bres , xgcn
2015-11-18 17:34:33 +01:00
real ( wp_ ) , dimension ( 6 ) , intent ( out ) :: dery
real ( wp_ ) , intent ( out ) :: psinv , dens , btot , xg , yg , anpl , anpr
real ( wp_ ) , intent ( out ) :: ddr , ddi , dersdst , derdnm
2016-02-12 17:49:00 +01:00
real ( wp_ ) , dimension ( 3 ) , intent ( out ) :: bv
2021-12-15 02:31:09 +01:00
integer , intent ( out ) :: error
2017-09-12 21:37:06 +02:00
real ( wp_ ) , dimension ( 3 ) , intent ( out ) :: derxg
2019-03-26 15:21:22 +01:00
integer , intent ( in ) :: igrad
2015-11-18 17:34:33 +01:00
! local variables
2023-03-29 22:13:43 +02:00
real ( wp_ ) , dimension ( 3 ) :: deryg
2015-11-18 17:34:33 +01:00
real ( wp_ ) , dimension ( 3 , 3 ) :: derbv
2015-11-19 19:20:58 +01:00
real ( wp_ ) , parameter :: anplth1 = 0.99_wp_ , anplth2 = 1.05_wp_
2015-11-18 17:34:33 +01:00
2023-03-29 22:13:43 +02:00
call plas_deriv ( xv , bres , xgcn , dens , btot , bv , derbv , xg , yg , derxg , deryg , psinv )
call disp_deriv ( anv , sox , xg , yg , derxg , deryg , bv , derbv , dgr , ddgr , igrad , &
dery , anpl , anpr , ddr , ddi , dersdst , derdnm )
2015-11-19 19:20:58 +01:00
2022-07-14 00:38:34 +02:00
if ( abs ( anpl ) > anplth2 ) then
error = raise_error ( error , large_npl , 1 )
else if ( abs ( anpl ) > anplth1 ) then
error = raise_error ( error , large_npl , 0 )
else
error = 0
2015-11-19 19:20:58 +01:00
end if
2015-11-18 17:34:33 +01:00
end subroutine ywppla_upd
2015-11-23 18:55:27 +01:00
2015-11-18 17:34:33 +01:00
subroutine gradi_upd ( ywrk , ak0 , xc , du1 , gri , ggri )
2021-12-19 16:39:19 +01:00
use const_and_precisions , only : zero , half
2015-11-18 17:34:33 +01:00
use beamdata , only : nray , nrayr , nrayth , twodr2
implicit none
real ( wp_ ) , intent ( in ) :: ak0
real ( wp_ ) , dimension ( 6 , nray ) , intent ( in ) :: ywrk
real ( wp_ ) , dimension ( 3 , nrayth , nrayr ) , intent ( inout ) :: xc , du1
real ( wp_ ) , dimension ( 3 , nray ) , intent ( out ) :: gri
real ( wp_ ) , dimension ( 3 , 3 , nray ) , intent ( out ) :: ggri
! local variables
real ( wp_ ) , dimension ( 3 , nrayth , nrayr ) :: xco , du1o
integer :: jk , j , jm , jp , k , km , kp
real ( wp_ ) :: ux , uxx , uxy , uxz , uy , uyy , uyz , uz , uzz
real ( wp_ ) :: dfuu , dffiu , gx , gxx , gxy , gxz , gy , gyy , gyz , gz , gzz
real ( wp_ ) , dimension ( 3 ) :: dxv1 , dxv2 , dxv3 , dgu
real ( wp_ ) , dimension ( 3 , 3 ) :: dgg , dff
! update position and du1 vectors
xco = xc
du1o = du1
jk = 1
do j = 1 , nrayr
do k = 1 , nrayth
if ( j > 1 ) jk = jk + 1
xc ( 1 : 3 , k , j ) = ywrk ( 1 : 3 , jk )
end do
end do
! compute grad u1 for central ray
j = 1
jp = 2
do k = 1 , nrayth
if ( k == 1 ) then
km = nrayth
else
km = k - 1
end if
if ( k == nrayth ) then
kp = 1
else
kp = k + 1
end if
dxv1 = xc ( : , k , jp ) - xc ( : , k , j )
dxv2 = xc ( : , kp , jp ) - xc ( : , km , jp )
dxv3 = xc ( : , k , j ) - xco ( : , k , j )
call solg0 ( dxv1 , dxv2 , dxv3 , dgu )
du1 ( : , k , j ) = dgu
end do
gri ( : , 1 ) = zero
2022-05-10 12:36:37 +02:00
2015-11-18 17:34:33 +01:00
! compute grad u1 and grad(S_I) for all the other rays
dfuu = twodr2 / ak0 ! twodr2 = 2*dr**2 = 2*(rwmax/(nrayr-1))**2
jm = 1
j = 2
k = 0
dffiu = dfuu
do jk = 2 , nray
k = k + 1
if ( k > nrayth ) then
jm = j
j = j + 1
k = 1
dffiu = dfuu * jm
end if
kp = k + 1
km = k - 1
if ( k == 1 ) then
km = nrayth
else if ( k == nrayth ) then
kp = 1
end if
dxv1 = xc ( : , k , j ) - xc ( : , k , jm )
dxv2 = xc ( : , kp , j ) - xc ( : , km , j )
dxv3 = xc ( : , k , j ) - xco ( : , k , j )
call solg0 ( dxv1 , dxv2 , dxv3 , dgu )
du1 ( : , k , j ) = dgu
gri ( : , jk ) = dgu ( : ) * dffiu
end do
! compute derivatives of grad u and grad(S_I) for rays jk>1
ggri ( : , : , 1 ) = zero
jm = 1
j = 2
k = 0
dffiu = dfuu
do jk = 2 , nray
k = k + 1
if ( k > nrayth ) then
jm = j
j = j + 1
k = 1
dffiu = dfuu * jm
end if
kp = k + 1
km = k - 1
if ( k == 1 ) then
km = nrayth
else if ( k == nrayth ) then
kp = 1
end if
dxv1 = xc ( : , k , j ) - xc ( : , k , jm )
dxv2 = xc ( : , kp , j ) - xc ( : , km , j )
dxv3 = xc ( : , k , j ) - xco ( : , k , j )
dff ( : , 1 ) = du1 ( : , k , j ) - du1 ( : , k , jm )
dff ( : , 2 ) = du1 ( : , kp , j ) - du1 ( : , km , j )
dff ( : , 3 ) = du1 ( : , k , j ) - du1o ( : , k , j )
call solg3 ( dxv1 , dxv2 , dxv3 , dff , dgg )
! derivatives of u
ux = du1 ( 1 , k , j )
uy = du1 ( 2 , k , j )
uz = du1 ( 3 , k , j )
uxx = dgg ( 1 , 1 )
uyy = dgg ( 2 , 2 )
uzz = dgg ( 3 , 3 )
uxy = ( dgg ( 1 , 2 ) + dgg ( 2 , 1 ) ) * half
uxz = ( dgg ( 1 , 3 ) + dgg ( 3 , 1 ) ) * half
uyz = ( dgg ( 2 , 3 ) + dgg ( 3 , 2 ) ) * half
! derivatives of S_I and Grad(S_I)
gx = ux * dffiu
gy = uy * dffiu
gz = uz * dffiu
gxx = dfuu * ux * ux + dffiu * uxx
gyy = dfuu * uy * uy + dffiu * uyy
gzz = dfuu * uz * uz + dffiu * uzz
gxy = dfuu * ux * uy + dffiu * uxy
gxz = dfuu * ux * uz + dffiu * uxz
gyz = dfuu * uy * uz + dffiu * uyz
ggri ( 1 , 1 , jk ) = gxx
ggri ( 2 , 1 , jk ) = gxy
ggri ( 3 , 1 , jk ) = gxz
ggri ( 1 , 2 , jk ) = gxy
ggri ( 2 , 2 , jk ) = gyy
ggri ( 3 , 2 , jk ) = gyz
ggri ( 1 , 3 , jk ) = gxz
ggri ( 2 , 3 , jk ) = gyz
ggri ( 3 , 3 , jk ) = gzz
end do
2022-05-10 12:36:37 +02:00
2015-11-18 17:34:33 +01:00
end subroutine gradi_upd
2015-11-23 18:55:27 +01:00
2015-11-18 17:34:33 +01:00
subroutine solg0 ( dxv1 , dxv2 , dxv3 , dgg )
! solution of the linear system of 3 eqs : dgg . dxv = dff
! input vectors : dxv1, dxv2, dxv3, dff
! output vector : dgg
! dff=(1,0,0)
implicit none
! arguments
real ( wp_ ) , dimension ( 3 ) , intent ( in ) :: dxv1 , dxv2 , dxv3
real ( wp_ ) , dimension ( 3 ) , intent ( out ) :: dgg
! local variables
real ( wp_ ) :: denom , aa1 , aa2 , aa3
aa1 = ( dxv2 ( 2 ) * dxv3 ( 3 ) - dxv3 ( 2 ) * dxv2 ( 3 ) )
aa2 = ( dxv1 ( 2 ) * dxv3 ( 3 ) - dxv3 ( 2 ) * dxv1 ( 3 ) )
aa3 = ( dxv1 ( 2 ) * dxv2 ( 3 ) - dxv2 ( 2 ) * dxv1 ( 3 ) )
denom = dxv1 ( 1 ) * aa1 - dxv2 ( 1 ) * aa2 + dxv3 ( 1 ) * aa3
dgg ( 1 ) = aa1 / denom
dgg ( 2 ) = - ( dxv2 ( 1 ) * dxv3 ( 3 ) - dxv3 ( 1 ) * dxv2 ( 3 ) ) / denom
dgg ( 3 ) = ( dxv2 ( 1 ) * dxv3 ( 2 ) - dxv3 ( 1 ) * dxv2 ( 2 ) ) / denom
end subroutine solg0
subroutine solg3 ( dxv1 , dxv2 , dxv3 , dff , dgg )
! rhs "matrix" dff, result in dgg
implicit none
! arguments
real ( wp_ ) , dimension ( 3 ) , intent ( in ) :: dxv1 , dxv2 , dxv3
real ( wp_ ) , dimension ( 3 , 3 ) , intent ( in ) :: dff
real ( wp_ ) , dimension ( 3 , 3 ) , intent ( out ) :: dgg
! local variables
real ( wp_ ) denom , a11 , a21 , a31 , a12 , a22 , a32 , a13 , a23 , a33
a11 = ( dxv2 ( 2 ) * dxv3 ( 3 ) - dxv3 ( 2 ) * dxv2 ( 3 ) )
a21 = ( dxv1 ( 2 ) * dxv3 ( 3 ) - dxv3 ( 2 ) * dxv1 ( 3 ) )
a31 = ( dxv1 ( 2 ) * dxv2 ( 3 ) - dxv2 ( 2 ) * dxv1 ( 3 ) )
a12 = ( dxv2 ( 1 ) * dxv3 ( 3 ) - dxv3 ( 1 ) * dxv2 ( 3 ) )
a22 = ( dxv1 ( 1 ) * dxv3 ( 3 ) - dxv3 ( 1 ) * dxv1 ( 3 ) )
a32 = ( dxv1 ( 1 ) * dxv2 ( 3 ) - dxv2 ( 1 ) * dxv1 ( 3 ) )
a13 = ( dxv2 ( 1 ) * dxv3 ( 2 ) - dxv3 ( 1 ) * dxv2 ( 2 ) )
a23 = ( dxv1 ( 1 ) * dxv3 ( 2 ) - dxv3 ( 1 ) * dxv1 ( 2 ) )
a33 = ( dxv1 ( 1 ) * dxv2 ( 2 ) - dxv2 ( 1 ) * dxv1 ( 2 ) )
denom = dxv1 ( 1 ) * a11 - dxv2 ( 1 ) * a21 + dxv3 ( 1 ) * a31
dgg ( : , 1 ) = ( dff ( : , 1 ) * a11 - dff ( : , 2 ) * a21 + dff ( : , 3 ) * a31 ) / denom
dgg ( : , 2 ) = ( - dff ( : , 1 ) * a12 + dff ( : , 2 ) * a22 - dff ( : , 3 ) * a32 ) / denom
dgg ( : , 3 ) = ( dff ( : , 1 ) * a13 - dff ( : , 2 ) * a23 + dff ( : , 3 ) * a33 ) / denom
end subroutine solg3
2015-11-23 18:55:27 +01:00
2023-03-29 22:13:43 +02:00
subroutine plas_deriv ( xv , bres , xgcn , dens , btot , bv , derbv , &
xg , yg , derxg , deryg , psinv_opt )
use const_and_precisions , only : zero
2022-05-10 12:36:37 +02:00
use gray_params , only : iequil
use equilibrium , only : psia , equinum_fpol , equinum_psi , equian , sgnbphi
use coreprofiles , only : density
2015-11-18 17:34:33 +01:00
implicit none
2022-05-10 12:36:37 +02:00
! subroutine arguments
2023-03-29 22:13:43 +02:00
real ( wp_ ) , dimension ( 3 ) , intent ( in ) :: xv
real ( wp_ ) , intent ( in ) :: xgcn , bres
real ( wp_ ) , intent ( out ) :: dens , btot , xg , yg
real ( wp_ ) , dimension ( 3 ) , intent ( out ) :: bv , derxg , deryg
2015-11-18 17:34:33 +01:00
real ( wp_ ) , dimension ( 3 , 3 ) , intent ( out ) :: derbv
2023-03-29 22:13:43 +02:00
real ( wp_ ) , optional , intent ( out ) :: psinv_opt
2022-05-10 12:36:37 +02:00
! local variables
2015-11-18 17:34:33 +01:00
integer :: jv
2023-03-29 22:13:43 +02:00
real ( wp_ ) :: xx , yy , zz
real ( wp_ ) :: psinv
2022-05-06 00:31:18 +02:00
real ( wp_ ) :: csphi , drrdx , drrdy , dphidx , dphidy , rr , rr2 , rrm , snphi , zzm
2015-11-18 17:34:33 +01:00
real ( wp_ ) , dimension ( 3 ) :: dbtot , bvc
real ( wp_ ) , dimension ( 3 , 3 ) :: dbvcdc , dbvdc , dbv
2023-03-29 22:13:43 +02:00
real ( wp_ ) :: brr , bphi , bzz , dxgdpsi
2015-11-18 17:34:33 +01:00
real ( wp_ ) :: dpsidr , dpsidz , ddpsidrr , ddpsidzz , ddpsidrz , fpolv , dfpv , ddenspsin
xg = zero
yg = 9 9._wp_
psinv = - 1._wp_
dens = zero
btot = zero
derxg = zero
deryg = zero
bv = zero
derbv = zero
2023-03-29 22:13:43 +02:00
if ( iequil == 0 ) then
! copy optional output
if ( present ( psinv_opt ) ) psinv_opt = psinv
return
end if
2015-11-18 17:34:33 +01:00
dbtot = zero
dbv = zero
dbvcdc = zero
dbvcdc = zero
dbvdc = zero
xx = xv ( 1 )
yy = xv ( 2 )
zz = xv ( 3 )
2022-05-10 12:36:37 +02:00
! cylindrical coordinates
2015-11-18 17:34:33 +01:00
rr2 = xx ** 2 + yy ** 2
rr = sqrt ( rr2 )
csphi = xx / rr
snphi = yy / rr
bv ( 1 ) = - snphi * sgnbphi
bv ( 2 ) = csphi * sgnbphi
2022-05-10 12:36:37 +02:00
! convert from cm to meters
2015-11-18 17:34:33 +01:00
zzm = 1.0e-2_wp_ * zz
rrm = 1.0e-2_wp_ * rr
if ( iequil == 1 ) then
call equian ( rrm , zzm , psinv , fpolv , dfpv , dpsidr , dpsidz , &
ddpsidrr , ddpsidzz , ddpsidrz )
else
call equinum_psi ( rrm , zzm , psinv , dpsidr , dpsidz , ddpsidrr , ddpsidzz , ddpsidrz )
call equinum_fpol ( psinv , fpolv , dfpv )
end if
2023-03-29 22:13:43 +02:00
! copy optional output
if ( present ( psinv_opt ) ) psinv_opt = psinv
2022-05-10 12:36:37 +02:00
! compute yg and derivative
2015-11-18 17:34:33 +01:00
if ( psinv < zero ) then
bphi = fpolv / rrm
btot = abs ( bphi )
yg = btot / bres
return
end if
2022-05-10 12:36:37 +02:00
! compute xg and derivative
2015-11-18 17:34:33 +01:00
call density ( psinv , dens , ddenspsin )
xg = xgcn * dens
dxgdpsi = xgcn * ddenspsin / psia
2022-05-10 12:36:37 +02:00
! B = f(psi)/R e_phi+ grad psi x e_phi/R
2015-11-18 17:34:33 +01:00
bphi = fpolv / rrm
brr = - dpsidz / rrm
bzz = dpsidr / rrm
2022-05-10 12:36:37 +02:00
! bvc(i) = B_i in cylindrical coordinates
bvc = [ brr , bphi , bzz ]
2015-11-18 17:34:33 +01:00
2022-05-10 12:36:37 +02:00
! bv(i) = B_i in cartesian coordinates
2015-11-18 17:34:33 +01:00
bv ( 1 ) = bvc ( 1 ) * csphi - bvc ( 2 ) * snphi
bv ( 2 ) = bvc ( 1 ) * snphi + bvc ( 2 ) * csphi
bv ( 3 ) = bvc ( 3 )
2022-05-10 12:36:37 +02:00
! dbvcdc(iv,jv) = d Bcil(iv) / dxvcil(jv)
2015-11-18 17:34:33 +01:00
dbvcdc ( 1 , 1 ) = - ddpsidrz / rrm - brr / rrm
dbvcdc ( 2 , 1 ) = dfpv * dpsidr / rrm - bphi / rrm
dbvcdc ( 3 , 1 ) = ddpsidrr / rrm - bzz / rrm
dbvcdc ( 1 , 3 ) = - ddpsidzz / rrm
dbvcdc ( 2 , 3 ) = dfpv * dpsidz / rrm
dbvcdc ( 3 , 3 ) = ddpsidrz / rrm
2022-05-10 12:36:37 +02:00
! dbvdc(iv,jv) = d Bcart(iv) / dxvcil(jv)
2015-11-18 17:34:33 +01:00
dbvdc ( 1 , 1 ) = dbvcdc ( 1 , 1 ) * csphi - dbvcdc ( 2 , 1 ) * snphi
dbvdc ( 2 , 1 ) = dbvcdc ( 1 , 1 ) * snphi + dbvcdc ( 2 , 1 ) * csphi
dbvdc ( 3 , 1 ) = dbvcdc ( 3 , 1 )
dbvdc ( 1 , 2 ) = - bv ( 2 )
dbvdc ( 2 , 2 ) = bv ( 1 )
2023-10-04 16:56:05 +02:00
dbvdc ( 3 , 2 ) = dbvcdc ( 3 , 2 ) ! = 0
2015-11-18 17:34:33 +01:00
dbvdc ( 1 , 3 ) = dbvcdc ( 1 , 3 ) * csphi - dbvcdc ( 2 , 3 ) * snphi
dbvdc ( 2 , 3 ) = dbvcdc ( 1 , 3 ) * snphi + dbvcdc ( 2 , 3 ) * csphi
dbvdc ( 3 , 3 ) = dbvcdc ( 3 , 3 )
drrdx = csphi
drrdy = snphi
dphidx = - snphi / rrm
dphidy = csphi / rrm
2022-05-10 12:36:37 +02:00
! dbv(iv,jv) = d Bcart(iv) / dxvcart(jv)
2015-11-18 17:34:33 +01:00
dbv ( : , 1 ) = drrdx * dbvdc ( : , 1 ) + dphidx * dbvdc ( : , 2 )
dbv ( : , 2 ) = drrdy * dbvdc ( : , 1 ) + dphidy * dbvdc ( : , 2 )
dbv ( : , 3 ) = dbvdc ( : , 3 )
2022-05-10 12:36:37 +02:00
! B magnitude and derivatives
2022-05-06 00:31:18 +02:00
btot = norm2 ( bv )
2015-11-18 17:34:33 +01:00
2023-10-04 16:56:05 +02:00
! dbtot(i) = d |B| / dxvcart(i)
2022-05-06 00:31:18 +02:00
dbtot = matmul ( bv , dbv ) / btot
2015-11-18 17:34:33 +01:00
yg = btot / Bres
2022-05-10 12:36:37 +02:00
! convert spatial derivatives from dummy/m -> dummy/cm
! to be used in rhs
2015-11-18 17:34:33 +01:00
2022-05-10 12:36:37 +02:00
! bv(i) = B_i / B ; derbv(i,j) = d (B_i / B) /d x,y,z
2015-11-18 17:34:33 +01:00
deryg = 1.0e-2_wp_ * dbtot / Bres
bv = bv / btot
do jv = 1 , 3
derbv ( : , jv ) = 1.0e-2_wp_ * ( dbv ( : , jv ) - bv ( : ) * dbtot ( jv ) ) / btot
end do
derxg ( 1 ) = 1.0e-2_wp_ * drrdx * dpsidr * dxgdpsi
derxg ( 2 ) = 1.0e-2_wp_ * drrdy * dpsidr * dxgdpsi
derxg ( 3 ) = 1.0e-2_wp_ * dpsidz * dxgdpsi
end subroutine plas_deriv
2015-11-23 18:55:27 +01:00
2023-03-29 22:13:43 +02:00
subroutine disp_deriv ( anv , sox , xg , yg , derxg , deryg , bv , derbv , &
dgr , ddgr , igrad , dery , anpl_ , anpr , ddr , ddi , &
dersdst , derdnm_ )
! Computes the dispersion relation, derivatives and other
! related quantities
!
! The mandatory outputs are used for both integrating the ray
! trajectory (`rhs` subroutine) and updating the ray state
! (`ywppla_upd` suborutine); while the optional ones are used for
! computing the absoprtion and current drive.
2022-05-10 12:36:37 +02:00
use const_and_precisions , only : zero , one , half , two
use gray_params , only : idst
2015-11-18 17:34:33 +01:00
implicit none
2022-05-10 12:36:37 +02:00
! subroutine arguments
2023-03-29 22:13:43 +02:00
! Inputs
! refractive index N̅ vector, b̅ = B̅/B magnetic field direction
real ( wp_ ) , dimension ( 3 ) , intent ( in ) :: anv , bv
! sign of polarisation mode: -1 ⇒ O, +1 ⇒ X
2022-05-20 13:10:49 +02:00
integer , intent ( in ) :: sox
2023-03-29 22:13:43 +02:00
! CMA diagram variables: X=(ω_pe/ω)², Y=ω_ce/ω
real ( wp_ ) , intent ( in ) :: xg , yg
! gradients of X, Y
real ( wp_ ) , dimension ( 3 ) , intent ( in ) :: derxg , deryg
! gradients of the complex eikonal
real ( wp_ ) , dimension ( 3 ) , intent ( in ) :: dgr ! ∇S_I
real ( wp_ ) , dimension ( 3 , 3 ) , intent ( in ) :: ddgr ! ∇∇S_I
! gradient of the magnetic field direction
real ( wp_ ) , dimension ( 3 , 3 ) , intent ( in ) :: derbv ! ∇b̅
! raytracing/beamtracing switch
integer , intent ( in ) :: igrad
! Outputs
! the actual derivatives: (∂Λ/∂x̅, -∂Λ/∂N̅)
real ( wp_ ) , dimension ( 6 ) , intent ( out ) :: dery
! additional quantities:
! refractive index
real ( wp_ ) , optional , intent ( out ) :: anpl_ , anpr ! N∥, N⊥
! real and imaginary part of the dispersion
real ( wp_ ) , optional , intent ( out ) :: ddr , ddi ! Λ, ∂Λ/∂N̅⋅∇S_I
! Jacobian ds/dσ , where s arclength, σ integration variable
real ( wp_ ) , optional , intent ( out ) :: dersdst
! |∂Λ/∂N̅|
real ( wp_ ) , optional , intent ( out ) :: derdnm_
! Note: assign values to missing optional arguments results in a segfault.
! Since some are always needed anyway, we store them here and copy
! them later, if needed:
real ( wp_ ) :: anpl , derdnm
2022-05-10 12:36:37 +02:00
! local variables
2023-03-29 22:13:43 +02:00
real ( wp_ ) :: gr2 , yg2 , anpl2 , del , dnl , duh , dan2sdnpl , an2 , an2s
real ( wp_ ) :: dan2sdxg , dan2sdyg , denom
2022-05-10 12:36:37 +02:00
real ( wp_ ) :: derdom , dfdiadnpl , dfdiadxg , dfdiadyg , fdia , bdotgr
real ( wp_ ) , dimension ( 3 ) :: derdxv , danpldxv , derdnv , dbgr
2015-11-18 17:34:33 +01:00
2023-03-29 22:13:43 +02:00
an2 = dot_product ( anv , anv ) ! N²
anpl = dot_product ( anv , bv ) ! N∥ = N̅⋅B̅
2015-11-18 17:34:33 +01:00
2023-03-29 22:13:43 +02:00
! Shorthands used in the expressions below
yg2 = yg ** 2 ! Y²
anpl2 = anpl ** 2 ! N∥²
dnl = one - anpl2 ! 1 - N∥²
duh = one - xg - yg2 ! UH denom (duh=0 on the upper-hybrid resonance)
! Compute/copy optional outputs
if ( present ( anpr ) ) anpr = sqrt ( max ( an2 - anpl2 , zero ) ) ! N⊥
if ( present ( anpl_ ) ) anpl_ = anpl ! N∥
2015-11-18 17:34:33 +01:00
an2s = one
dan2sdxg = zero
dan2sdyg = zero
dan2sdnpl = zero
del = zero
fdia = zero
dfdiadnpl = zero
dfdiadxg = zero
dfdiadyg = zero
if ( xg > zero ) then
2023-03-29 22:13:43 +02:00
! Derivatives of the cold plasma refractive index
2022-05-10 12:36:37 +02:00
!
2023-03-29 22:13:43 +02:00
! N²s = 1 - X - XY²⋅(1 + N∥² ± √Δ)/[2(1 - X - Y²)]
2022-05-10 12:36:37 +02:00
!
2023-03-29 22:13:43 +02:00
! where Δ = (1 - N∥²)² + 4N∥²⋅(1 - X)/Y²
! + for the X mode, - for the O mode
! √Δ
del = sqrt ( dnl ** 2 + 4.0_wp_ * anpl2 * ( one - xg ) / yg2 )
2015-11-18 17:34:33 +01:00
2022-05-10 12:36:37 +02:00
! ∂(N²s)/∂X
2023-03-29 11:19:00 +02:00
! Note: this term is nonzero for X=0, but it multiplies terms
! proportional to X or ∂X/∂ψ which are zero outside the plasma.
2015-11-18 17:34:33 +01:00
dan2sdxg = - half * yg2 * ( one - yg2 ) * ( one + anpl2 + sox * del ) / duh ** 2 &
+ sox * xg * anpl2 / ( del * duh ) - one
2022-05-10 12:36:37 +02:00
! ∂(N²s)/∂Y
2015-11-18 17:34:33 +01:00
dan2sdyg = - xg * yg * ( one - xg ) * ( one + anpl2 + sox * del ) / duh ** 2 &
+ two * sox * xg * ( one - xg ) * anpl2 / ( yg * del * duh )
2022-05-10 12:36:37 +02:00
! ∂(N²s)/∂N∥
2015-11-18 17:34:33 +01:00
dan2sdnpl = - xg * yg2 * anpl / duh &
- sox * xg * anpl * ( two * ( one - xg ) - yg2 * dnl ) / ( del * duh )
if ( igrad > 0 ) then
2023-03-29 22:13:43 +02:00
! Derivatives used in the complex eikonal terms (beamtracing only)
block
real ( wp_ ) :: ddelnpl2 , ddelnpl2x , ddelnpl2y , derdel
! ∂²Δ/∂N∥²
ddelnpl2 = two * ( two * ( one - xg ) * ( one + 3.0_wp_ * anpl2 ** 2 ) &
- yg2 * dnl ** 3 ) / yg2 / del ** 3
! ∂²(N²s)/∂N∥²
fdia = - xg * yg2 * ( one + half * sox * ddelnpl2 ) / duh
! Intermediates results used right below
derdel = two * ( one - xg ) * anpl2 * ( one + 3.0_wp_ * anpl2 ** 2 ) &
- dnl ** 2 * ( one + 3.0_wp_ * anpl2 ) * yg2
derdel = 4.0_wp_ * derdel / ( yg * del ) ** 5
ddelnpl2y = two * ( one - xg ) * derdel
ddelnpl2x = yg * derdel
! ∂³(N²s)/∂N∥³
dfdiadnpl = 2 4.0_wp_ * sox * xg * ( one - xg ) * anpl * ( one - anpl2 ** 2 ) &
/ ( yg2 * del ** 5 )
! ∂³(N²s)/∂N∥²∂X
dfdiadxg = - yg2 * ( one - yg2 ) / duh ** 2 - sox * yg2 * ( ( one - yg2 ) &
* ddelnpl2 + xg * duh * ddelnpl2x ) / ( two * duh ** 2 )
! ∂³(N²s)/∂N∥²∂Y
dfdiadyg = - two * yg * xg * ( one - xg ) / duh ** 2 &
- sox * xg * yg * ( two * ( one - xg ) * ddelnpl2 &
+ yg * duh * ddelnpl2y ) / ( two * duh ** 2 )
end block
2015-11-18 17:34:33 +01:00
end if
end if
2023-03-29 22:13:43 +02:00
! ∇(N∥) = ∇(N̅⋅b̅) = N̅⋅∇b̅
2022-05-10 12:36:37 +02:00
danpldxv = matmul ( anv , derbv )
2015-11-18 17:34:33 +01:00
2023-03-29 22:13:43 +02:00
! ∂Λ/∂x̅ = - ∇(N²s) = -∂Λ/∂X ∇X - ∂Λ/∂Y ∇Y - ∂Λ/∂N∥ ∇(N∥)
derdxv = - ( derxg * dan2sdxg + deryg * dan2sdyg + danpldxv * dan2sdnpl )
! ∂Λ/∂N̅ = 2N̅ - ∂(N²s)/∂N̅ = 2N̅ - ∂(N²s)/∂N∥ b̅
! Note: we used the identity ∇f(v̅⋅b̅) = f' ∇(v̅⋅b̅) = f'b̅.
derdnv = two * anv - dan2sdnpl * bv
! ∂Λ/∂ω = ∂N²/∂ω - ∂N²s/∂X⋅∂X/∂ω - ∂N²s/∂Y⋅∂Y/∂ω - ∂N²s/∂N∥⋅∂N∥/∂ω
! Notes: 1. N depends on ω: N²=c²k²/ω² ⇒ ∂N²/∂ω = -2N²/ω
! N∥=ck∥/ω ⇒ ∂N∥/∂ω = -N∥/ω
! 2. derdom is actually ω⋅∂Λ/∂ω, see below for the reason.
! 3. N gains a dependency on ω because Λ(∇S, ω) is computed
! on the constrains Λ=0.
derdom = - two * an2 + two * xg * dan2sdxg + yg * dan2sdyg + anpl * dan2sdnpl
if ( igrad > 0 ) then
! Complex eikonal terms added to the above expressions
block
real ( wp_ ) :: dgr2 ( 3 )
bdotgr = dot_product ( bv , dgr ) ! b̅⋅∇S_I
gr2 = dot_product ( dgr , dgr ) ! |∇S_I|²
dgr2 = 2 * matmul ( ddgr , dgr ) ! ∇|∇S_I|² = 2 ∇∇S_I⋅∇S_I
! ∇(b̅⋅∇S_I) = ∇b̅ᵗ ∇S_I + b̅ᵗ ∇∇S_I
! Notes: 1. We are using the convention ∇b̅_ij = ∂b_i/∂x_j.
! 2. Fortran doesn't distinguish between column/row vectors,
! so matmul(A, v̅) ≅ Av̅ and matmul(v̅, A) ≅ v̅ᵗA ≅ Aᵗv̅.
dbgr = matmul ( dgr , derbv ) + matmul ( bv , ddgr )
! ∂Λ/∂x̅ += - ∇(|∇S_I|²) + ½ ∇[(b̅⋅∇S_I)² ∂²N²s/∂N∥²]
derdxv = derdxv - dgr2 + fdia * bdotgr * dbgr + half * bdotgr ** 2 &
* ( derxg * dfdiadxg + deryg * dfdiadyg + danpldxv * dfdiadnpl )
! ∂Λ/∂N̅ += ½(b̅⋅∇S_I)² ∂³(N²s)/∂N∥³ b̅
! Note: we used again the identity ∇f(v̅⋅b̅) = f'b̅.
derdnv = derdnv + half * bdotgr ** 2 * dfdiadnpl * bv
! ∂Λ/∂ω += ∂|∇S_I|²/∂ω + ½∂(b⋅∇S_I)²/∂ω + ½(b̅⋅∇S_I)² ∂/∂ω (∂²N²s/∂N∥²)
! Note: as above ∇S_I gains a dependency on ω
derdom = derdom + two * gr2 - bdotgr ** 2 &
* ( fdia + xg * dfdiadxg + half * yg * dfdiadyg &
+ half * anpl * dfdiadnpl )
end block
end if
2015-11-18 17:34:33 +01:00
2023-03-29 22:13:43 +02:00
derdnm = norm2 ( derdnv )
! Denominator of the r.h.s, depending on the
! choice of the integration variable σ :
!
! dx̅/dσ = + ∂Λ/∂N̅ / denom
! dN̅/dσ = - ∂Λ/∂x̅ / denom
!
select case ( idst )
case ( 0 ) ! σ =s, arclength
! denom = |∂Λ/∂N̅|
! Proof: Normalising ∂Λ/∂N̅ (∥ to the group velocity)
! simply makes σ the arclength parameter s.
2022-05-10 12:36:37 +02:00
denom = derdnm
2023-03-29 22:13:43 +02:00
case ( 1 ) ! σ =ct, time
! denom = -ω⋅∂Λ/∂ω
! Proof: The Hamilton equations are
!
! dx̅/dt = + ∂H/∂k̅ (H=ℏΩ, p=ℏk̅)
! dp̅/dt = - ∂H/∂x̅
!
! where Ω(x̅, k̅) is implicitely defined by the dispersion
! D(k̅, ω, x̅) = 0. By differentiating the latter we get:
!
! ∂Ω/∂x̅ = - ∂D/∂x̅ / ∂D/∂ω = - ∂Λ/∂x̅ / ∂Λ/∂ω
! ∂Ω/∂k̅ = - ∂D/∂k̅ / ∂D/∂ω = - ∂Λ/∂k̅ / ∂Λ/∂ω
!
! with Λ=D/k₀, k₀=ω/c. Finally, substituting k̅=k₀N̅:
!
! dx̅/d(ct) = + ∂Λ/∂N̅ / (-ω⋅∂Λ/∂ω)
! dN̅/d(ct) = - ∂Λ/∂x̅ / (-ω⋅∂Λ/∂ω)
!
2022-05-10 12:36:37 +02:00
denom = - derdom
2023-03-29 22:13:43 +02:00
case ( 2 ) ! s=S_R, phase
! denom = N̅⋅∂Λ/∂N̅
! Note: This parametrises the curve by the value
! of the (real) phase, with the wave frozen in time.
! Proof: By definition N̅ = k̅/k₀ = ∇S. Using the gradient
! theorem on the ray curve parametrised by the arclength
!
! S(s) = ∫₀^x̅(s) N̅⋅dl̅ = ∫₀^s N̅(σ )⋅dx̅/ds ds
!
! Differentiating gives dS = N̅⋅dx̅/ds ds, so
!
! dS = N̅⋅∂Λ/∂N̅ / |∂Λ/∂N̅| ds ⇒
!
! dx̅/dS = + ∂Λ/∂N̅ / N̅⋅∂Λ/∂N̅
! dN̅/dS = - ∂Λ/∂x̅ / N̅⋅∂Λ/∂N̅
!
2022-05-10 12:36:37 +02:00
denom = dot_product ( anv , derdnv )
end select
2015-11-18 17:34:33 +01:00
2023-03-29 22:13:43 +02:00
! Jacobian ds/dσ , where s is the arclength,
! for computing integrals in dσ , like:
!
! τ = ∫α(s)ds = ∫α(σ )(ds/dσ )dσ
!
if ( present ( dersdst ) ) dersdst = derdnm / denom
! |∂Λ/∂N̅|: used in the α computation (see alpha_effj)
if ( present ( derdnm_ ) ) derdnm_ = derdnm
2015-11-18 17:34:33 +01:00
2022-05-10 12:36:37 +02:00
! r.h.s. vector
2023-03-29 22:13:43 +02:00
dery ( 1 : 3 ) = derdnv ( : ) / denom ! +∂Λ/∂N̅
2022-05-10 12:36:37 +02:00
dery ( 4 : 6 ) = - derdxv ( : ) / denom ! -∂Λ/∂x̅
2015-11-18 17:34:33 +01:00
2023-03-29 22:13:43 +02:00
if ( present ( ddr ) . or . present ( ddi ) ) then
! Dispersion relation (geometric optics)
! ddr ← Λ = N² - N²s(X,Y,N∥) = 0
an2s = one - xg - half * xg * yg2 * ( one + anpl2 + sox * del ) / duh
ddr = an2 - an2s
end if
if ( present ( ddr ) . and . igrad > 0 ) then
! Dispersion relation (complex eikonal)
! ddr ← Λ = N² - N²s - |∇S_I|² + ½ (b̅⋅∇S_I)² ∂²N²s/∂N∥²
ddr = ddr - gr2 - half * bdotgr ** 2 * fdia ! real part
end if
! Note: we have to return ddi even for igrad=0 because
! it's printed to udisp unconditionally
if ( present ( ddi ) ) then
! Dispersion relation (complex eikonal)
! ddi ← ∂Λ/∂N̅⋅∇S_I
ddi = merge ( dot_product ( derdnv , dgr ) , zero , igrad > 0 ) ! imaginary part
end if
2015-11-18 17:34:33 +01:00
end subroutine disp_deriv
2015-11-23 18:55:27 +01:00
2022-05-20 13:10:49 +02:00
subroutine alpha_effj ( params , psinv , X , Y , density , temperature , &
k0 , Bres , derdnm , Npl , Npr_cold , sox , &
Npr , alpha , dIdp , nhmin , nhmax , iokhawa , error )
2022-05-10 12:36:37 +02:00
! Computes the absorption coefficient and effective current
2022-05-20 13:10:49 +02:00
use const_and_precisions , only : pi , mc2 = > mc2_
2022-05-10 12:36:37 +02:00
use gray_params , only : ecrh_cd_parameters
use coreprofiles , only : fzeff
use equilibrium , only : sgnbphi
use dispersion , only : harmnumber , warmdisp
use eccd , only : setcdcoeff , eccdeff , fjch0 , fjch , fjncl
2022-07-14 00:38:34 +02:00
use gray_errors , only : negative_absorption , raise_error
2022-05-10 12:36:37 +02:00
use magsurf_data , only : fluxval
2015-11-18 17:34:33 +01:00
implicit none
2022-05-10 12:36:37 +02:00
! subroutine arguments
! Inputs
! ECRH & CD parameters
type ( ecrh_cd_parameters ) :: params
! poloidal flux ψ
real ( wp_ ) , intent ( in ) :: psinv
! CMA diagram variables: X=(ω_pe/ω)², Y=ω_ce/ω
2022-05-20 13:10:49 +02:00
real ( wp_ ) , intent ( in ) :: X , Y
! densityity [10¹⁹ m⁻³], temperature [keV]
real ( wp_ ) , intent ( in ) :: density , temperature
2022-05-10 12:36:37 +02:00
! vacuum wavenumber k₀=ω/c, resonant B field
2022-05-20 13:10:49 +02:00
real ( wp_ ) , intent ( in ) :: k0 , Bres
2022-05-10 12:36:37 +02:00
! group velocity: |∂Λ/∂N̅| where Λ=c²/ω²⋅D_R
real ( wp_ ) , intent ( in ) :: derdnm
! refractive index: N∥, N⊥ (cold dispersion)
2022-05-20 13:10:49 +02:00
real ( wp_ ) , intent ( in ) :: Npl , Npr_cold
2022-05-10 12:36:37 +02:00
! sign of polarisation mode: -1 ⇒ O, +1 ⇒ X
2022-05-20 13:10:49 +02:00
integer , intent ( in ) :: sox
2022-05-10 12:36:37 +02:00
! Outputs
2022-05-20 13:10:49 +02:00
! orthogonal refractive index N⊥ (solution of the warm dispersion)
complex ( wp_ ) , intent ( out ) :: Npr
2022-05-10 12:36:37 +02:00
! absorption coefficient, current density
2022-05-20 13:10:49 +02:00
real ( wp_ ) , intent ( out ) :: alpha , dIdp
2022-05-10 12:36:37 +02:00
! lowest/highest resonant harmonic numbers
integer , intent ( out ) :: nhmin , nhmax
! ECCD/overall error codes
integer , intent ( out ) :: iokhawa , error
! local variables
real ( wp_ ) :: rbavi , rrii , rhop
2022-05-20 13:10:49 +02:00
integer :: nlarmor , ithn , ierrcd
real ( wp_ ) :: mu , rbn , rbx
2022-05-10 12:36:37 +02:00
real ( wp_ ) :: zeff , cst2 , bmxi , bmni , fci
real ( wp_ ) , dimension ( : ) , pointer :: eccdpar = > null ( )
2022-05-20 13:10:49 +02:00
real ( wp_ ) :: effjcd , effjcdav , Btot
complex ( wp_ ) :: e ( 3 )
alpha = 0
Npr = 0
dIdp = 0
2022-05-10 12:36:37 +02:00
nhmin = 0
nhmax = 0
iokhawa = 0
error = 0
! Absorption computation
! Skip when the temperature is zero (no plasma)
2022-05-20 13:10:49 +02:00
if ( temperature < = 0 ) return
2022-05-10 12:36:37 +02:00
! Skip when the lowest resonant harmonic number is zero
2022-05-20 13:10:49 +02:00
mu = mc2 / temperature ! μ=(mc²/kT)
call harmnumber ( Y , mu , Npl ** 2 , params % iwarm == 1 , nhmin , nhmax )
2022-05-10 12:36:37 +02:00
if ( nhmin < = 0 ) return
2022-05-20 13:10:49 +02:00
! Solve the full dispersion only when needed
if ( params % iwarm / = 4 . or . params % ieccd / = 0 ) then
nlarmor = max ( params % ilarm , nhmax )
if ( params % ieccd / = 0 ) then
! Compute the polarisation vector only when current drive is on
call warmdisp ( X , Y , mu , Npl , Npr_cold , sox , error , Npr , e , &
model = params % iwarm , nlarmor = nlarmor , &
max_iters = abs ( params % imx ) , fallback = params % imx < 0 )
else
call warmdisp ( X , Y , mu , Npl , Npr_cold , sox , error , Npr , &
model = params % iwarm , nlarmor = nlarmor , &
max_iters = abs ( params % imx ) , fallback = params % imx < 0 )
end if
end if
2022-05-10 12:36:37 +02:00
2022-05-20 13:10:49 +02:00
! Compute α from the solution of the dispersion relation
2022-05-10 12:36:37 +02:00
! The absoption coefficient is defined as
!
! α = 2 Im(k̅)⋅s̅
!
! where s̅ = v̅_g/|v_g|, the direction of the energy flow.
2022-05-20 13:10:49 +02:00
! Since v̅_g = - ∂Λ/∂N̅ / ∂Λ/∂ω, using the cold dispersion
! relation, we have that
2022-05-10 12:36:37 +02:00
!
! s̅ = ∂Λ/∂N̅ / |∂Λ/∂N̅|
2022-05-20 13:10:49 +02:00
! = [2N̅ - ∂(N²s)/∂N∥ b̅
! + ½(b̅⋅∇S_I)² ∂³(N²s)/∂N∥³ b̅] / |∂Λ/∂N̅|
2022-05-10 12:36:37 +02:00
!
! Assuming Im(k∥)=0:
!
! α = 4 Im(k⊥)⋅N⊥ / |∂Λ/∂N̅|
!
2022-05-20 13:10:49 +02:00
block
real ( wp_ ) :: k_im
k_im = k0 * Npr % im ! imaginary part of k⊥
alpha = 4 * k_im * Npr_cold / derdnm
end block
2022-05-10 12:36:37 +02:00
2022-05-20 13:10:49 +02:00
if ( alpha < 0 ) then
2022-07-14 00:38:34 +02:00
error = raise_error ( error , negative_absorption )
2022-05-10 12:36:37 +02:00
return
2015-11-18 17:34:33 +01:00
end if
2022-05-10 12:36:37 +02:00
! Current drive computation
if ( params % ieccd < = 0 ) return
zeff = fzeff ( psinv )
ithn = 1
2022-05-20 13:10:49 +02:00
if ( nlarmor > nhmin ) ithn = 2
2022-05-10 12:36:37 +02:00
rhop = sqrt ( psinv )
call fluxval ( rhop , rri = rrii , rbav = rbavi , bmn = bmni , bmx = bmxi , fc = fci )
2022-05-20 13:10:49 +02:00
Btot = Y * Bres
rbn = Btot / bmni
rbx = Btot / bmxi
2022-05-10 12:36:37 +02:00
select case ( params % ieccd )
case ( 1 )
! Cohen model
2022-05-20 13:10:49 +02:00
call setcdcoeff ( zeff , rbn , rbx , cst2 , eccdpar )
call eccdeff ( Y , Npl , Npr % re , density , mu , e , nhmin , nhmax , &
ithn , cst2 , fjch , eccdpar , effjcd , iokhawa , ierrcd )
2022-05-10 12:36:37 +02:00
case ( 2 )
! No trapping
2022-05-20 13:10:49 +02:00
call setcdcoeff ( zeff , cst2 , eccdpar )
call eccdeff ( Y , Npl , Npr % re , density , mu , e , nhmin , nhmax , &
ithn , cst2 , fjch0 , eccdpar , effjcd , iokhawa , ierrcd )
2022-05-10 12:36:37 +02:00
case default
! Neoclassical model
2022-05-20 13:10:49 +02:00
call setcdcoeff ( zeff , rbx , fci , mu , rhop , cst2 , eccdpar )
call eccdeff ( Y , Npl , Npr % re , density , mu , e , nhmin , nhmax , &
ithn , cst2 , fjncl , eccdpar , effjcd , iokhawa , ierrcd )
2022-05-10 12:36:37 +02:00
end select
error = error + ierrcd
if ( associated ( eccdpar ) ) deallocate ( eccdpar )
! current drive efficiency <j/p> [A⋅m/W]
effjcdav = rbavi * effjcd
2022-05-20 13:10:49 +02:00
dIdp = sgnbphi * effjcdav / ( 2 * pi * rrii )
2022-05-10 12:36:37 +02:00
2015-11-18 17:34:33 +01:00
end subroutine alpha_effj
2022-05-10 12:36:37 +02:00
subroutine set_pol ( ywrk0 , bres , sox , ipol , psipol0 , chipol0 , ext0 , eyt0 )
use const_and_precisions , only : degree , zero , one , half , im
use beamdata , only : nray , nrayth
use equilibrium , only : bfield
use polarization , only : pol_limit , polellipse , &
stokes_ce , stokes_ell
2015-11-18 17:34:33 +01:00
implicit none
2022-05-10 12:36:37 +02:00
! subroutine arguments
real ( wp_ ) , dimension ( 6 , nray ) , intent ( in ) :: ywrk0
2022-05-20 13:10:49 +02:00
real ( wp_ ) , intent ( in ) :: bres
integer , intent ( in ) :: sox , ipol
2022-05-10 12:36:37 +02:00
real ( wp_ ) , intent ( inout ) :: psipol0 , chipol0
complex ( wp_ ) , dimension ( nray ) , intent ( out ) :: ext0 , eyt0
! local variables
2015-11-18 17:34:33 +01:00
integer :: j , k , jk
real ( wp_ ) , dimension ( 3 ) :: xmv , anv , bv
real ( wp_ ) :: rm , csphi , snphi , bphi , br , bz , qq , uu , vv , deltapol
j = 1
k = 0
do jk = 1 , nray
k = k + 1
if ( jk == 2 . or . k > nrayth ) then
j = j + 1
k = 1
end if
if ( ipol == 0 ) then
xmv = ywrk0 ( 1 : 3 , jk ) * 0.01_wp_ ! convert from cm to m
anv = ywrk0 ( 4 : 6 , jk )
rm = sqrt ( xmv ( 1 ) ** 2 + xmv ( 2 ) ** 2 )
csphi = xmv ( 1 ) / rm
snphi = xmv ( 2 ) / rm
call bfield ( rm , xmv ( 3 ) , bphi , br , bz )
2022-05-10 12:36:37 +02:00
! bv(i) = B_i in cartesian coordinates
2015-11-18 17:34:33 +01:00
bv ( 1 ) = br * csphi - bphi * snphi
bv ( 2 ) = br * snphi + bphi * csphi
bv ( 3 ) = bz
call pol_limit ( anv , bv , bres , sox , ext0 ( jk ) , eyt0 ( jk ) )
if ( jk == 1 ) then
call stokes_ce ( ext0 ( jk ) , eyt0 ( jk ) , qq , uu , vv )
call polellipse ( qq , uu , vv , psipol0 , chipol0 )
psipol0 = psipol0 / degree ! convert from rad to degree
chipol0 = chipol0 / degree
end if
else
call stokes_ell ( chipol0 * degree , psipol0 * degree , qq , uu , vv )
if ( qq ** 2 < one ) then
deltapol = asin ( vv / sqrt ( one - qq ** 2 ) )
ext0 ( jk ) = sqrt ( half * ( one + qq ) )
eyt0 ( jk ) = sqrt ( half * ( one - qq ) ) * exp ( - im * deltapol )
else
ext0 ( jk ) = one
eyt0 ( jk ) = zero
end if
endif
end do
end subroutine set_pol
2015-11-23 18:55:27 +01:00
subroutine cniteq ( rqgrid , zqgrid , matr2dgrid , nr , nz , h , ncon , npts , icount , rcon , zcon )
! v2.01 12/07/95 -- written by d v bartlett, jet joint undertaking.
! (based on an older code)
implicit none
! arguments
integer , intent ( in ) :: nr , nz
real ( wp_ ) , dimension ( nr ) , intent ( in ) :: rqgrid
real ( wp_ ) , dimension ( nz ) , intent ( in ) :: zqgrid
real ( wp_ ) , dimension ( nr , nz ) , intent ( in ) :: matr2dgrid
real ( wp_ ) , intent ( in ) :: h
integer , intent ( inout ) :: ncon , icount
integer , dimension ( ncon ) , intent ( out ) :: npts
real ( wp_ ) , dimension ( icount ) , intent ( out ) :: rcon , zcon
! local variables
integer :: i , j , k , l , nrqmax , iclast , mpl , ix , jx , mxr , n1 , jm , jfor , lda , ldb
integer :: jabs , jnb , kx , ikx , itm , inext , in
integer , dimension ( 3 , 2 ) :: ja
integer , dimension ( icount / 2 - 1 ) :: lx
real ( wp_ ) :: drgrd , dzgrd , ah , adn , px , x , y
real ( wp_ ) , dimension ( nr * nz ) :: a
logical :: flag1
px = 0.5_wp_
2022-05-10 12:36:37 +02:00
a = reshape ( matr2dgrid , [ nr * nz ] )
2015-11-23 18:55:27 +01:00
rcon = 0.0_wp_
zcon = 0.0_wp_
nrqmax = nr
drgrd = rqgrid ( 2 ) - rqgrid ( 1 )
dzgrd = zqgrid ( 2 ) - zqgrid ( 1 )
ncon = 0
npts = 0
iclast = 0
icount = 0
mpl = 0
ix = 0
mxr = nrqmax * ( nz - 1 )
n1 = nr - 1
do jx = 2 , n1
do jm = jx , mxr , nrqmax
j = jm + nrqmax
ah = a ( j ) - h
if ( ah < = 0.0_wp_ . and . a ( jm ) > h . or . &
ah > 0.0_wp_ . and . a ( jm ) < = h ) then
ix = ix + 1
lx ( ix ) = - j
end if
if ( ah < = 0.0_wp_ . and . a ( j - 1 ) > h . or . &
ah > 0.0_wp_ . and . a ( j - 1 ) < = h ) then
ix = ix + 1
lx ( ix ) = j
end if
end do
end do
do jm = nr , mxr , nrqmax
j = jm + nrqmax
ah = a ( j ) - h
if ( ah < = 0.0_wp_ . and . a ( j - 1 ) > h . or . &
ah > 0.0_wp_ . and . a ( j - 1 ) < = h ) then
ix = ix + 1
lx ( ix ) = j
end if
if ( ah < = 0.0_wp_ . and . a ( jm ) > h . or . &
ah > 0.0_wp_ . and . a ( jm ) < = h ) then
ix = ix + 1
lx ( ix ) = - j
end if
end do
do jm = 1 , mxr , nrqmax
j = jm + nrqmax
if ( a ( j ) < = h . and . a ( jm ) > h . or . &
a ( j ) > h . and . a ( jm ) < = h ) then
ix = ix + 1
lx ( ix ) = - j
end if
end do
do j = 2 , nr
if ( a ( j ) < = h . and . a ( j - 1 ) > h . or . &
a ( j ) > h . and . a ( j - 1 ) < = h ) then
ix = ix + 1
lx ( ix ) = j
end if
end do
if ( ix < = 0 ) return
2022-05-10 12:36:37 +02:00
2015-11-23 18:55:27 +01:00
bb : do
in = ix
jx = lx ( in )
jfor = 0
lda = 1
ldb = 2
do
if ( jx < 0 ) then
jabs = - jx
jnb = jabs - nrqmax
else
jabs = jx
jnb = jabs - 1
end if
adn = a ( jabs ) - a ( jnb )
if ( adn / = 0 ) px = ( a ( jabs ) - h ) / adn
kx = ( jabs - 1 ) / nrqmax
ikx = jabs - nrqmax * kx - 1
if ( jx < 0 ) then
x = drgrd * ikx
y = dzgrd * ( kx - px )
else
x = drgrd * ( ikx - px )
y = dzgrd * kx
end if
icount = icount + 1
rcon ( icount ) = x + rqgrid ( 1 )
zcon ( icount ) = y + zqgrid ( 1 )
mpl = icount
itm = 1
ja ( 1 , 1 ) = jabs + nrqmax
j = 1
if ( jx < = 0 ) then
ja ( 1 , 1 ) = - jabs - 1
j = 2
end if
ja ( 2 , 1 ) = - ja ( 1 , 1 )
ja ( 3 , 1 ) = - jx + 1 - nrqmax
ja ( 3 , 2 ) = - jx
ja ( j , 2 ) = jabs - nrqmax
k = 3 - j
ja ( k , 2 ) = 1 - jabs
if ( kx < = 0 . or . ikx < = 0 ) then
lda = 1
ldb = lda
else if ( ikx + 1 - nr > = 0 . and . jx < = 0 ) then
lda = 2
ldb = lda
else if ( jfor / = 0 ) then
lda = 2
do i = 1 , 3
if ( jfor == ja ( i , 2 ) ) then
lda = 1
exit
end if
end do
ldb = lda
end if
flag1 = . false .
aa : do k = 1 , 3
do l = lda , ldb
do i = 1 , ix
if ( lx ( i ) == ja ( k , l ) ) then
itm = itm + 1
inext = i
if ( jfor / = 0 ) exit aa
2022-05-10 12:36:37 +02:00
if ( itm > 3 ) then
2015-11-23 18:55:27 +01:00
flag1 = . true .
exit aa
end if
end if
end do
end do
end do aa
if ( . not . flag1 ) then
lx ( in ) = 0
2022-05-10 12:36:37 +02:00
if ( itm == 1 ) exit
end if
2015-11-23 18:55:27 +01:00
jfor = jx
jx = lx ( inext )
in = inext
end do
do
if ( lx ( ix ) / = 0 ) then
if ( mpl > = 4 ) then
ncon = ncon + 1
npts ( ncon ) = icount - iclast
iclast = icount
end if
exit
end if
ix = ix - 1
if ( ix < = 0 ) exit bb
end do
end do bb
if ( mpl > = 4 ) then
ncon = ncon + 1
npts ( ncon ) = icount - iclast
iclast = icount
end if
end subroutine cniteq
2022-05-03 23:16:21 +02:00
subroutine print_headers ( params )
2021-12-19 16:39:19 +01:00
! Prints the headers of all gray_main output tables
2022-05-03 23:16:21 +02:00
use units , only : uprj0 , uwbm , udisp , ucenr , uoutr , upec , usumm , &
unit_active , active_units
use gray_params , only : gray_parameters , headw , headl , print_parameters
2021-12-19 16:39:19 +01:00
2015-11-23 18:55:27 +01:00
implicit none
2021-12-19 16:39:19 +01:00
2021-12-15 02:31:09 +01:00
! subroutine arguments
2022-05-03 23:16:21 +02:00
type ( gray_parameters ) , intent ( in ) :: params
2021-12-19 16:39:19 +01:00
2021-12-15 02:31:09 +01:00
! local variables
2021-12-19 16:39:19 +01:00
integer :: i , j
integer :: main_units ( 8 )
character ( 256 ) :: main_headers ( 8 )
2022-05-03 23:16:21 +02:00
character ( len = headw ) , dimension ( headl ) :: header
2021-12-19 16:39:19 +01:00
main_units = [ uprj0 , uprj0 + 1 , uwbm , udisp , ucenr , uoutr , upec , usumm ]
main_headers ( 1 ) = '#sst j k xt yt zt rt'
main_headers ( 2 ) = '#sst j k xt yt zt rt'
main_headers ( 3 ) = '#sst w1 w2'
main_headers ( 4 ) = '#sst Dr_Nr Di_Nr'
main_headers ( 5 ) = '#sst R z phi psin rhot ne Te Btot Bx By Bz Nperp ' &
/ / 'Npl Nx Ny Nz ki alpha tau Pt dIds nhmin nhmax ' &
/ / 'iohkw index_rt ddr Xg Yg dXgdx dXgdy dXgdz'
main_headers ( 6 ) = '#i k sst x y R z psin tau Npl alpha index_rt'
main_headers ( 7 ) = '#rhop rhot Jphi Jcdb dPdV Icdins Pins index_rt'
main_headers ( 8 ) = '#Icd Pa Jphip dPdVp rhotj rhotjava rhotp rhotpav ' &
/ / 'drhotjava drhotpav ratjamx ratjbmx stmx psipol ' &
/ / 'chipol index_rt Jphimx dPdVmx drhotj drhotp P0 ' &
/ / 'cplO cplX'
2022-05-03 23:16:21 +02:00
if ( all ( active_units == 0 ) ) return
call print_parameters ( params , header )
2021-12-19 16:39:19 +01:00
do i = 1 , size ( main_units )
if ( unit_active ( main_units ( i ) ) ) then
2022-05-03 23:16:21 +02:00
do j = 1 , size ( header )
write ( main_units ( i ) , '(1x,a)' ) header ( j )
2021-12-19 16:39:19 +01:00
end do
write ( main_units ( i ) , '(1x,a)' ) trim ( main_headers ( i ) )
end if
2016-04-27 16:37:57 +02:00
end do
2015-11-23 18:55:27 +01:00
2021-12-19 16:39:19 +01:00
end subroutine print_headers
2015-11-23 18:55:27 +01:00
2023-05-23 17:49:51 +02:00
subroutine print_prof ( params )
! Prints the (input) plasma kinetic profiles (unit 55)
2021-12-19 16:39:19 +01:00
2023-05-23 17:49:51 +02:00
use gray_params , only : profiles_parameters
2022-12-18 14:09:40 +01:00
use equilibrium , only : q_spline , fq , frhotor , tor_curr_psi
use coreprofiles , only : density , temp
use units , only : uprfin , unit_active
2021-12-19 16:39:19 +01:00
2015-11-23 18:55:27 +01:00
implicit none
2021-12-19 16:39:19 +01:00
2023-05-23 17:49:51 +02:00
! suborutine arguments
type ( profiles_parameters ) , intent ( in ) :: params
2021-12-19 16:39:19 +01:00
! local constants
real ( wp_ ) , parameter :: eps = 1.e-4_wp_
! local variables
2023-05-23 17:49:51 +02:00
integer :: i , N_data , N_tail
real ( wp_ ) :: rhot , jphi , dens , ddens
real ( wp_ ) :: psin , dpsin , psin_last
2021-12-19 16:39:19 +01:00
if ( . not . unit_active ( uprfin ) ) return
2015-11-23 18:55:27 +01:00
2023-05-23 17:49:51 +02:00
! N of profiles data points
N_data = q_spline % ndata
! parameters for the density tail (numerical profiles only)
if ( params % iprof == 1 ) then
N_tail = N_data / 10 ! N of density tail points
psin_last = q_spline % data ( N_data ) ! last data point
dpsin = ( params % psnbnd - psin_last ) / N_tail ! Δψ for uniform grid
else
N_tail = 0
end if
2021-12-19 16:39:19 +01:00
write ( uprfin , * ) '#psi rhot ne Te q Jphi'
2023-05-23 17:49:51 +02:00
do i = 1 , N_data + N_tail
if ( i > N_data ) then
psin = psin_last + dpsin * ( i - N_data )
else
psin = q_spline % data ( i )
end if
2021-12-19 16:39:19 +01:00
rhot = frhotor ( sqrt ( psin ) )
call density ( psin , dens , ddens )
2022-12-18 14:09:40 +01:00
jphi = tor_curr_psi ( max ( eps , psin ) )
2015-11-23 18:55:27 +01:00
2021-12-19 16:39:19 +01:00
write ( uprfin , '(12(1x,e12.5))' ) &
2022-12-18 14:09:40 +01:00
psin , rhot , dens , temp ( psin ) , fq ( psin ) , jphi * 1.e-6_wp_
2015-11-23 18:55:27 +01:00
end do
end subroutine print_prof
subroutine print_bres ( bres )
2021-12-19 16:39:19 +01:00
! Prints the EC resonance surface table (unit 70)
2022-12-18 14:09:40 +01:00
use equilibrium , only : rmnm , rmxm , zmnm , zmxm , bfield , q_spline
use units , only : ubres , unit_active
2021-12-19 16:39:19 +01:00
2015-11-23 18:55:27 +01:00
implicit none
2021-12-19 16:39:19 +01:00
! subroutine arguments
real ( wp_ ) , intent ( in ) :: bres
! local constants
integer , parameter :: icmx = 2002
! local variables
2015-11-23 18:55:27 +01:00
integer :: j , k , n , nconts , nctot
integer , dimension ( 10 ) :: ncpts
real ( wp_ ) :: dr , dz , btmx , btmn , zzk , rrj , bbphi , bbr , bbz , bbb
real ( wp_ ) , dimension ( icmx ) :: rrcb , zzcb
2022-12-18 14:09:40 +01:00
real ( wp_ ) :: rv ( q_spline % ndata ) , zv ( q_spline % ndata )
real ( wp_ ) , dimension ( q_spline % ndata , q_spline % ndata ) :: btotal
2015-11-23 18:55:27 +01:00
2021-12-19 16:39:19 +01:00
if ( . not . unit_active ( ubres ) ) return
2022-05-10 12:36:37 +02:00
2022-12-18 14:09:40 +01:00
dr = ( rmxm - rmnm ) / ( q_spline % ndata - 1 )
dz = ( zmxm - zmnm ) / ( q_spline % ndata - 1 )
do j = 1 , q_spline % ndata
2021-12-19 16:39:19 +01:00
rv ( j ) = rmnm + dr * ( j - 1 )
zv ( j ) = zmnm + dz * ( j - 1 )
2022-05-10 12:36:37 +02:00
end do
2015-11-23 18:55:27 +01:00
2021-12-19 16:39:19 +01:00
! Btotal on psi grid
btmx = - 1.0e30_wp_
btmn = 1.0e30_wp_
2022-12-18 14:09:40 +01:00
do k = 1 , q_spline % ndata
2021-12-19 16:39:19 +01:00
zzk = zv ( k )
2022-12-18 14:09:40 +01:00
do j = 1 , q_spline % ndata
2021-12-19 16:39:19 +01:00
rrj = rv ( j )
call bfield ( rrj , zzk , bbphi , bbr , bbz )
btotal ( j , k ) = sqrt ( bbr ** 2 + bbz ** 2 + bbphi ** 2 )
if ( btotal ( j , k ) > = btmx ) btmx = btotal ( j , k )
if ( btotal ( j , k ) < = btmn ) btmn = btotal ( j , k )
end do
end do
! compute Btot=Bres/n with n=1,5
write ( ubres , * ) '#i Btot R z'
2022-12-18 14:09:40 +01:00
do n = 1 , 5
2023-09-20 16:03:08 +02:00
bbb = bres / n
2021-12-19 16:39:19 +01:00
if ( bbb > = btmn . and . bbb < = btmx ) then
nconts = size ( ncpts )
nctot = size ( rrcb )
2022-12-18 14:09:40 +01:00
call cniteq ( rv , zv , btotal , q_spline % ndata , q_spline % ndata , bbb , &
2021-12-19 16:39:19 +01:00
nconts , ncpts , nctot , rrcb , zzcb )
2022-12-18 14:09:40 +01:00
do j = 1 , nctot
2021-12-19 16:39:19 +01:00
write ( ubres , '(i6,12(1x,e12.5))' ) j , bbb , rrcb ( j ) , zzcb ( j )
2015-11-23 18:55:27 +01:00
end do
end if
2021-12-19 16:39:19 +01:00
write ( ubres , * )
2015-11-23 18:55:27 +01:00
end do
2021-12-19 16:39:19 +01:00
2015-11-23 18:55:27 +01:00
end subroutine print_bres
2021-12-19 16:39:19 +01:00
subroutine print_maps ( bres , xgcn , r0 , anpl0 )
! Prints several input quantities on a regular grid (unit 72)
2015-11-23 18:55:27 +01:00
2021-12-19 16:39:19 +01:00
use gray_params , only : iequil
use equilibrium , only : rmnm , rmxm , zmnm , zmxm , equian , equinum_psi , &
2022-12-18 14:09:40 +01:00
equinum_fpol , q_spline
2016-02-12 17:49:00 +01:00
use coreprofiles , only : density , temp
2021-12-19 16:39:19 +01:00
use units , only : umaps , unit_active
2016-02-12 17:49:00 +01:00
implicit none
2021-12-19 16:39:19 +01:00
! subroutine arguments
2016-02-12 17:49:00 +01:00
real ( wp_ ) , intent ( in ) :: bres , xgcn , r0 , anpl0
2021-12-19 16:39:19 +01:00
! local variables
2016-02-12 17:49:00 +01:00
integer :: j , k
2021-12-19 16:39:19 +01:00
real ( wp_ ) :: dr , dz , zk , rj , bphi , br , bz , btot , &
psin , ne , dne , te , xg , yg , anpl
2022-12-18 14:09:40 +01:00
real ( wp_ ) , dimension ( q_spline % ndata ) :: r , z
2016-02-12 17:49:00 +01:00
2021-12-19 16:39:19 +01:00
if ( . not . unit_active ( umaps ) ) return
2022-12-18 14:09:40 +01:00
dr = ( rmxm - rmnm ) / ( q_spline % ndata - 1 )
dz = ( zmxm - zmnm ) / ( q_spline % ndata - 1 )
do j = 1 , q_spline % ndata
2021-12-19 16:39:19 +01:00
r ( j ) = rmnm + dr * ( j - 1 )
z ( j ) = zmnm + dz * ( j - 1 )
2016-02-12 17:49:00 +01:00
end do
2022-05-10 12:36:37 +02:00
2021-12-19 16:39:19 +01:00
write ( umaps , * ) '#R z psin Br Bphi Bz Btot ne Te X Y Npl'
2022-12-18 14:09:40 +01:00
do j = 1 , q_spline % ndata
2021-12-19 16:39:19 +01:00
rj = r ( j )
anpl = anpl0 * r0 / rj
2022-12-18 14:09:40 +01:00
do k = 1 , q_spline % ndata
2021-12-19 16:39:19 +01:00
zk = z ( k )
2016-02-12 17:49:00 +01:00
if ( iequil < 2 ) then
rework the analytical model
This change modifies the analytical equilibrium in order to simplify the
computation of the poloidal flux normalization and the derivatives.
In the power law parametrisation of the safety factor, ρ_t is replaced
with ρ_p and, similarly, the normalised poloidal radius is now
identified with ρ_p, instead of ρ_t.
With the same parameters (q₀,q₁,α...), this choice slightly changes the
plasma current distribution, but enables us to obtain a closed form for
ψ_a = ψ(r=a) and the relation ρ_t(ρ_p). In fact, both expressions are
now obtained by integrating the q(ρ_p), instead of 1/q(ρ_t), which has
no elementary antiderivative.
As the normalisation is now computed exactly, the values of the
normalised flux ψ_n = ψ/ψ_a and the gradient ∇ψ (entering the raytracing
equations in X and ∇X, respectively) are computed to the same precision.
Previously, ψ_n was computed to a lower precision due to the use of a
simple trapezoid integration of 1/q(ρ_p) for ψ_a, while ∇ψ was computed
up to machine precision using an exact formula.
This error effectively caused a very slight decoupling between X=ω_p²/ω²
and ∇X that introduced a systematic error in the numerical solution of
the raytracing equations.
The error manifests itself as a bias with a weak dependency on X in the
values taken by the dispersion function Λ(r̅, n̅) on the phase-space
points generated by the integrator. More specifically,
lim h→0 Λ(r̅_i, n̅_i) = -kX(r̅_i)
where h is the integrator step size;
r̅_i is the position at the i-th step;
k ≈ -3.258⋅10⁻⁵ and depends only on the number of points used to
perform the trapedoid integral for ψ_a (as ~ 1/n²).
After this change Λ behaves consistently with being a conserved quantity
(zero) up to the cumulative integration error of the 4° order
Runge-Kutta method. In fact we now have that:
Λ(r̅_i, n̅_i) ∝ - h⁴ ‖∂⁴X(r̅_i)/∂r̅⁴‖
It must be said that within this model the relation ρ_p(ρ_t) can't be
computed analytically (inverting ρ_t(ρ_p) produces a trascendental
equation of the form b = x + c x^α). However, this relation is not
necessary for raytracing and is easily solved, up to machine
precision, using minpack.
In addition, this change also makes the model consistetly use the
cocos=3 and fully implements the ability to force the signs of I_p, B_φ
(via equilibrium.sgni,sgnb) and rescaling the field (via
equilibrium.factb).
2023-10-11 17:29:24 +02:00
call equian ( rj , zk , psin = psin , fpolv = bphi , dpsidr = bz , dpsidz = br )
2016-02-12 17:49:00 +01:00
else
2022-12-18 14:09:40 +01:00
call equinum_psi ( rj , zk , psi = psin , dpsidr = bz , dpsidz = br )
call equinum_fpol ( psin , fpol = bphi )
2016-02-12 17:49:00 +01:00
end if
br = - br / rj
bphi = bphi / rj
bz = bz / rj
2021-12-19 16:39:19 +01:00
btot = sqrt ( br ** 2 + bphi ** 2 + bz ** 2 )
2016-02-12 17:49:00 +01:00
yg = btot / bres
te = temp ( psin )
2021-12-19 16:39:19 +01:00
call density ( psin , ne , dne )
2016-02-12 17:49:00 +01:00
xg = xgcn * ne
2021-12-19 16:39:19 +01:00
write ( umaps , '(12(x,e12.5))' ) &
rj , zk , psin , br , bphi , bz , btot , ne , te , xg , yg , anpl
end do
write ( umaps , * )
end do
2016-02-12 17:49:00 +01:00
end subroutine print_maps
2015-11-23 18:55:27 +01:00
subroutine print_surfq ( qval )
2021-12-19 16:39:19 +01:00
! Print constant ψ surfaces for a given `q` value
2022-12-18 14:09:40 +01:00
use equilibrium , only : q_spline , fq , frhotor , &
2021-12-18 18:57:38 +01:00
rmaxis , zmaxis , zbsup , zbinf
use magsurf_data , only : contours_psi , npoints , print_contour
use utils , only : locate , intlin
use logger , only : log_info
2015-11-23 18:55:27 +01:00
implicit none
2021-12-18 18:57:38 +01:00
! subroutine arguments
2015-11-23 18:55:27 +01:00
real ( wp_ ) , dimension ( : ) , intent ( in ) :: qval
2021-12-18 18:57:38 +01:00
! local variables
2016-06-01 15:49:35 +02:00
integer :: i1 , i
2015-11-23 18:55:27 +01:00
real ( wp_ ) :: rup , zup , rlw , zlw , rhot , psival
real ( wp_ ) , dimension ( npoints ) :: rcn , zcn
2022-12-18 14:09:40 +01:00
real ( wp_ ) , dimension ( q_spline % ndata ) :: qpsi
2021-12-18 18:57:38 +01:00
character ( 256 ) :: msg ! for log messages formatting
2015-11-23 18:55:27 +01:00
2022-12-18 14:09:40 +01:00
! build the q profile on the ψ grid
do i = 1 , q_spline % ndata
qpsi ( i ) = fq ( q_spline % data ( i ) )
2015-11-23 18:55:27 +01:00
end do
2021-12-18 18:57:38 +01:00
! locate ψ surface for q=qval
call log_info ( 'constant ψ surfaces for:' , &
mod = 'gray_core' , proc = 'print_surfq' )
2022-12-18 14:09:40 +01:00
do i = 1 , size ( qval )
2021-12-18 18:57:38 +01:00
! FIXME: check for non monotonous q profile
2022-12-18 14:09:40 +01:00
call locate ( abs ( qpsi ) , q_spline % ndata , qval ( i ) , i1 )
if ( i1 > 0 . and . i1 < q_spline % ndata ) then
call intlin ( abs ( qpsi ( i1 ) ) , q_spline % data ( i1 ) , abs ( qpsi ( i1 + 1 ) ) , &
q_spline % data ( i1 + 1 ) , qval ( i ) , psival )
2021-12-19 16:39:19 +01:00
rup = rmaxis
rlw = rmaxis
zup = ( zbsup + zmaxis ) / 2.0_wp_
zlw = ( zmaxis + zbinf ) / 2.0_wp_
call contours_psi ( psival , rcn , zcn , rup , zup , rlw , zlw )
call print_contour ( psival , rcn , zcn )
rhot = frhotor ( sqrt ( psival ) )
2021-12-18 18:57:38 +01:00
write ( msg , '(4(x,a,"=",g0.3))' ) &
'q' , qval ( i ) , 'ψ' , psival , 'rhop' , sqrt ( psival ) , 'rhot' , rhot
call log_info ( msg , mod = 'gray_core' , proc = 'print_surfq' )
2015-11-23 18:55:27 +01:00
end if
end do
2021-12-18 18:57:38 +01:00
end subroutine print_surfq
2015-11-23 18:55:27 +01:00
2021-12-19 16:39:19 +01:00
subroutine print_projxyzt ( stv , ywrk , iproj )
! Prints the beam shape tables (unit 8, 12)
use const_and_precisions , only : comp_huge , zero , one
use beamdata , only : nray , nrayr , nrayth , rayi2jk
use units , only : uprj0 , uwbm , unit_active
2015-11-23 18:55:27 +01:00
implicit none
2021-12-19 16:39:19 +01:00
! subroutine arguments
real ( wp_ ) , dimension ( : ) , intent ( in ) :: stv
2015-11-23 18:55:27 +01:00
real ( wp_ ) , dimension ( : , : ) , intent ( in ) :: ywrk
2021-12-19 16:39:19 +01:00
integer , intent ( in ) :: iproj
! local variables
integer :: jk , jkz , uprj
2015-11-23 18:55:27 +01:00
integer , dimension ( 2 ) :: jkv
2021-12-19 16:39:19 +01:00
real ( wp_ ) , dimension ( 3 ) :: xv1 , dir , dxv
real ( wp_ ) :: rtimn , rtimx , csth1 , snth1 , &
csps1 , snps1 , xti , yti , zti , rti
if ( . not . ( unit_active ( uprj0 ) . or . unit_active ( uwbm ) ) ) return
! common/external functions/variables
2015-11-23 18:55:27 +01:00
uprj = uprj0 + iproj
2021-12-19 16:39:19 +01:00
xv1 = ywrk ( 1 : 3 , 1 )
dir = ywrk ( 4 : 6 , 1 )
dir = dir / norm2 ( dir )
2015-11-23 18:55:27 +01:00
csth1 = dir ( 3 )
snth1 = sqrt ( one - csth1 ** 2 )
if ( snth1 > zero ) then
2021-12-19 16:39:19 +01:00
csps1 = dir ( 2 ) / snth1
snps1 = dir ( 1 ) / snth1
2015-11-23 18:55:27 +01:00
else
2021-12-19 16:39:19 +01:00
csps1 = one
snps1 = zero
2015-11-23 18:55:27 +01:00
end if
if ( iproj == 0 ) then
jkz = nray - nrayth + 1
else
jkz = 1
end if
rtimn = comp_huge
rtimx = zero
do jk = jkz , nray
dxv = ywrk ( 1 : 3 , jk ) - xv1
xti = dxv ( 1 ) * csps1 - dxv ( 2 ) * snps1
yti = ( dxv ( 1 ) * snps1 + dxv ( 2 ) * csps1 ) * csth1 - dxv ( 3 ) * snth1
zti = ( dxv ( 1 ) * snps1 + dxv ( 2 ) * csps1 ) * snth1 + dxv ( 3 ) * csth1
2021-12-19 16:39:19 +01:00
rti = hypot ( xti , yti )
2015-11-23 18:55:27 +01:00
jkv = rayi2jk ( jk )
2021-12-19 16:39:19 +01:00
if ( . not . ( iproj == 0 . and . jk == 1 ) ) &
write ( uprj , '(1x,e16.8e3,2i5,4(1x,e16.8e3))' ) &
stv ( jk ) , jkv , xti , yti , zti , rti
if ( iproj == 1 . and . jkv ( 2 ) == nrayth ) write ( uprj , * )
2015-11-23 18:55:27 +01:00
if ( rti > = rtimx . and . jkv ( 1 ) == nrayr ) rtimx = rti
if ( rti < = rtimn . and . jkv ( 1 ) == nrayr ) rtimn = rti
end do
2021-12-19 16:39:19 +01:00
write ( uprj , * )
write ( uwbm , '(3(1x,e16.8e3))' ) stv ( 1 ) , rtimn , rtimx
2015-11-23 18:55:27 +01:00
end subroutine print_projxyzt
2021-12-19 16:39:19 +01:00
subroutine print_output ( i , jk , st , qj , xv , psinv , btot , bv , ak0 , &
anpl , anpr , anv , anprim , dens , tekev , &
alpha , tau , dids , nhm , nhf , iokhawa , &
index_rt , ddr , ddi , xg , yg , derxg )
! Prints the ray variables and Hamiltonian error (units 4, 33, 17)
2015-11-23 18:55:27 +01:00
2015-11-18 17:34:33 +01:00
use const_and_precisions , only : degree , zero , one
2021-12-19 16:39:19 +01:00
use equilibrium , only : frhotor
use gray_params , only : istpl0
use beamdata , only : nray , nrayth , jkray1
use units , only : ucenr , uoutr , udisp , unit_active
2015-11-18 17:34:33 +01:00
implicit none
2021-12-19 16:39:19 +01:00
! subroutine arguments
integer , intent ( in ) :: i , jk , nhm , nhf , iokhawa , index_rt
real ( wp_ ) , dimension ( 3 ) , intent ( in ) :: xv , bv , anv
real ( wp_ ) , intent ( in ) :: st , qj , psinv , btot , ak0 , anpl , anpr , anprim
real ( wp_ ) , intent ( in ) :: dens , tekev , alpha , tau , dids , ddr , ddi
real ( wp_ ) , intent ( in ) :: xg , yg
2017-09-12 21:37:06 +02:00
real ( wp_ ) , dimension ( 3 ) , intent ( in ) :: derxg
2021-12-19 16:39:19 +01:00
! local variables
real ( wp_ ) :: stm , xxm , yym , zzm , rrm , phideg , rhot , akim , pt , didsn
2015-11-18 17:34:33 +01:00
integer :: k
2021-12-19 16:39:19 +01:00
stm = st * 1.0e-2_wp_
xxm = xv ( 1 ) * 1.0e-2_wp_
yym = xv ( 2 ) * 1.0e-2_wp_
zzm = xv ( 3 ) * 1.0e-2_wp_
rrm = hypot ( xxm , yym )
2015-11-18 17:34:33 +01:00
2021-12-19 16:39:19 +01:00
! print central ray trajectory. dIds in A/m/W, ki in m⁻¹
if ( jk == 1 ) then
phideg = atan2 ( yym , xxm ) / degree
2015-11-18 17:34:33 +01:00
if ( psinv > = zero . and . psinv < = one ) then
2021-12-19 16:39:19 +01:00
rhot = frhotor ( sqrt ( psinv ) )
2015-11-18 17:34:33 +01:00
else
2021-12-19 16:39:19 +01:00
rhot = 1.0_wp_
end if
akim = anprim * ak0 * 1.0e2_wp_
pt = qj * exp ( - tau )
didsn = dids * 1.0e2_wp_ / qj
! print central ray data
if ( unit_active ( ucenr ) ) then
write ( ucenr , '(22(1x,e16.8e3),4i5,6(1x,e16.8e3))' ) &
stm , rrm , zzm , phideg , psinv , rhot , dens , tekev , &
btot , bv , anpr , anpl , anv , akim , alpha , tau , pt , &
didsn , nhm , nhf , iokhawa , index_rt , ddr , xg , yg , &
derxg
2015-11-18 17:34:33 +01:00
end if
end if
2021-12-19 16:39:19 +01:00
! print conservation of dispersion relation (Hamiltonian error)
if ( unit_active ( udisp ) . and . jk == nray ) then
write ( udisp , '(30(1x,e16.8e3))' ) st , ddr , ddi
end if
2015-11-18 17:34:33 +01:00
2021-12-19 16:39:19 +01:00
! print outer ray trajectories
if ( mod ( i , istpl0 ) == 0 ) then
2015-11-23 18:55:27 +01:00
k = jk - jkray1 + 1
2021-12-19 16:39:19 +01:00
if ( unit_active ( uoutr ) . and . k > 0 . and . k < = nrayth ) then
write ( uoutr , '(2i5,9(1x,e16.8e3),i5)' ) &
i , k , stm , xxm , yym , rrm , zzm , &
psinv , tau , anpl , alpha , index_rt
2015-11-18 17:34:33 +01:00
end if
end if
end subroutine print_output
2015-11-23 18:55:27 +01:00
2021-12-19 16:39:19 +01:00
subroutine print_pec ( rhop_tab , rhot_tab , jphi , jcd , &
dpdv , currins , pins , index_rt )
! Prints the ECRH and CD profiles (unit 48)
use units , only : upec , unit_active
2015-11-23 18:55:27 +01:00
implicit none
2021-12-19 16:39:19 +01:00
! subroutine arguments
real ( wp_ ) , dimension ( : ) , intent ( in ) :: rhop_tab , rhot_tab , jphi , jcd , &
dpdv , currins , pins
integer , intent ( in ) :: index_rt
! local variables
2015-11-23 18:55:27 +01:00
integer :: i
2021-12-19 16:39:19 +01:00
if ( . not . unit_active ( upec ) ) return
do i = 1 , size ( rhop_tab )
write ( upec , '(7(1x, e16.8e3), i5)' ) rhop_tab ( i ) , rhot_tab ( i ) , &
jphi ( i ) , jcd ( i ) , dpdv ( i ) , currins ( i ) , pins ( i ) , index_rt
2022-05-10 12:36:37 +02:00
end do
2021-12-19 16:39:19 +01:00
write ( upec , * ) ''
2015-11-23 18:55:27 +01:00
end subroutine print_pec
2021-12-19 16:39:19 +01:00
subroutine print_finals ( pabs , icd , dpdvp , jphip , rhotpav , rhotjava , drhotpav , &
drhotjava , dpdvmx , jphimx , rhotp , rhotj , drhotp , &
drhotj , ratjamx , ratjbmx , stmx , psipol , chipol , &
index_rt , p0 , cpl1 , cpl2 )
! Prints global data and final results (unit 7)
use units , only : ucenr , usumm , unit_active
2015-11-23 18:55:27 +01:00
implicit none
2021-12-19 16:39:19 +01:00
! subroutine arguments
real ( wp_ ) , intent ( in ) :: pabs , icd , dpdvp , jphip , rhotpav , rhotjava , &
drhotpav , drhotjava , dpdvmx , jphimx , rhotp , &
rhotj , drhotp , drhotj , ratjamx , ratjbmx , &
stmx , psipol , chipol , p0 , cpl1 , cpl2
2015-11-23 18:55:27 +01:00
integer , intent ( in ) :: index_rt
2021-12-19 16:39:19 +01:00
if ( unit_active ( ucenr ) ) then
write ( ucenr , * ) ''
end if
if ( unit_active ( usumm ) ) then
write ( usumm , '(15(1x,e12.5),i5,7(1x,e12.5))' ) &
icd , pabs , jphip , dpdvp , rhotj , rhotjava , rhotp , rhotpav , &
drhotjava , drhotpav , ratjamx , ratjbmx , stmx , psipol , &
chipol , index_rt , jphimx , dpdvmx , drhotj , drhotp , &
p0 , cpl1 , cpl2
end if
2015-11-23 18:55:27 +01:00
end subroutine print_finals
2021-12-15 02:31:14 +01:00
end module gray_core