! This modules contains the core GRAY routines module gray_core use const_and_precisions, only : wp_ implicit none contains subroutine gray_main(params, data, plasma, results, error, rhout) use const_and_precisions, only : zero, one, comp_tiny use polarization, only : ellipse_to_field use types, only : table, wrap use gray_params, only : gray_parameters, gray_data, gray_results, EQ_VACUUM use gray_plasma, only : abstract_plasma use gray_tables, only : init_tables, store_ec_profiles, store_ray_data, & store_beam_shape, find_flux_surfaces, & kinetic_profiles, ec_resonance, inputs_maps use gray_errors, only : is_critical, print_err_raytracing, print_err_ecrh_cd use beams, only : xgygcoeff, launchangles2n use beamdata, only : pweight use magsurf_data, only : compute_flux_averages, dealloc_surfvec use pec, only : pec_init, spec, postproc_profiles, dealloc_pec, & rhop_tab, rhot_tab use utils, only : vmaxmin, inside use multipass, only : initbeam, initmultipass, turnoffray, & plasma_in, plasma_out, wall_out use logger, only : log_info, log_debug ! subroutine arguments type(gray_parameters), intent(inout) :: params type(gray_data), intent(in) :: data class(abstract_plasma), intent(in) :: plasma type(gray_results), intent(out) :: results ! Predefined grid for the output profiles (optional) real(wp_), dimension(:), intent(in), optional :: rhout ! Exit code integer, intent(out) :: error ! local variables real(wp_), parameter :: taucr = 12._wp_, etaucr = exp(-taucr) character, dimension(2), parameter :: mode=['O','X'] integer :: sox real(wp_) :: ak0,bres,xgcn,xg,yg,rrm,zzm,alpha,didp,anpl,anpr,anprim,anprre real(wp_) :: chipol,psipol,btot,psinv,dens,tekev,dersdst,derdnm real(wp_) :: tau,pow,dids,ddr,ddi,taumn,taumx real(wp_) :: rhotpav,drhotpav,rhotjava,drhotjava,dpdvp,jphip real(wp_) :: rhotp,drhotp,rhotj,drhotj,dpdvmx,jphimx,ratjamx,ratjbmx real(wp_) :: pabs_beam,icd_beam,cpl_beam1,cpl_beam2,cpl_cbeam1,cpl_cbeam2 integer :: iox,nharm,nhf,iokhawa,ierrn,ierrhcd, index_rt, parent_index_rt integer :: ip,ib,iopmin,ipar,child_index_rt integer :: igrad_b,istop_pass,nbeam_pass logical :: ins_pl,ins_wl,ent_pl,ext_pl,ent_wl,ext_wl,iboff ! i: integration step, jk: global ray index integer :: i, jk ! buffer for formatting log messages character(256) :: msg real(wp_), dimension(2) :: pabs_pass,icd_pass,cpl real(wp_), dimension(3) :: xv,anv0,anv,bv,derxg associate (nray => params%raytracing%nray, & nrayr => params%raytracing%nrayr, & nrayth => params%raytracing%nrayth, & nstep => params%raytracing%nstep, & npass => params%raytracing%ipass, & nbeam_tot => 2**(params%raytracing%ipass+1)-2, & nbeam_max => 2**(params%raytracing%ipass)) block ! ray variables real(wp_), dimension(6, nray) :: yw, ypw real(wp_), dimension(3, nray) :: gri real(wp_), dimension(3, 3, nray) :: ggri real(wp_), dimension(3, nrayth, nrayr) :: xc, du1 real(wp_), dimension(nray) :: ccci0, p0jk real(wp_), dimension(nray) :: tau0, alphaabs0 real(wp_), dimension(nray) :: dids0 integer, dimension(nray) :: iiv real(wp_), dimension(nray, nstep) :: psjki, ppabs, ccci complex(wp_), dimension(nray) :: ext, eyt ! multipass variables logical, dimension(nray) :: iwait integer, dimension(nray) :: iow, iop logical, dimension(nray, nbeam_tot) :: iroff real(wp_), dimension(6, nray, nbeam_max-2) :: yynext, yypnext real(wp_), dimension(6, nray) :: yw0, ypw0 real(wp_), dimension(nray, nbeam_tot) :: stnext real(wp_), dimension(nray) :: stv, p0ray real(wp_), dimension(nray, nbeam_tot) :: taus, cpls real(wp_), dimension(nray) :: tau1, etau1, cpl1, lgcpl1 real(wp_), dimension(0:nbeam_tot) :: psipv, chipv ! beam variables real(wp_), dimension(params%output%nrho) :: jphi_beam, pins_beam real(wp_), dimension(params%output%nrho) :: currins_beam, dpdv_beam, jcd_beam ! ======== set environment BEGIN ======== ! Compute X=ω/ω_ce and Y=(ω/ω_pe)² (with B=1) call xgygcoeff(params%antenna%fghz, ak0, bres, xgcn) ! Initialise the output tables call init_tables(params, results%tables) ! Compute the initial cartesian wavevector (anv0) ! from launch angles α,β and the position call launchangles2n(params%antenna, anv0) if (params%equilibrium%iequil /= EQ_VACUUM) then ! Initialise the magsurf_data module call compute_flux_averages(params, results%tables) ! Initialise the output profiles call pec_init(params%output, rhout) end if ! 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-determinted tables results%tables%kinetic_profiles = kinetic_profiles(params, plasma) results%tables%ec_resonance = ec_resonance(params, bres) results%tables%inputs_maps = inputs_maps(params, plasma, bres, xgcn) ! print ψ surface for q=3/2 and q=2/1 call find_flux_surfaces( & qvals=[1.5_wp_, 2.0_wp_], params=params, & tbl=results%tables%flux_surfaces) ! print initial position 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') ! =========== main loop BEGIN =========== call initmultipass(params, iroff, yynext, yypnext, yw0, ypw0, & stnext, p0ray, taus, tau1, etau1, cpls, & cpl1, lgcpl1, psipv, chipv) sox=1 ! mode inverted for each beam iox=2 ! start with O: sox=-1, iox=1 call pweight(params%raytracing, params%antenna%power, p0jk) ! Set original polarisation psipv(0) = params%antenna%psi chipv(0) = params%antenna%chi nbeam_pass=1 ! max n of beam per pass index_rt=0 ! global beam index: 1,O 2,X 1st pass ! | | | | call log_debug('pass loop start', mod='gray_core', proc='gray_main') ! 3,O 4,X 5,O 6,X 2nd pass main_loop: do ip=1,params%raytracing%ipass write (msg, '("pass: ",g0)') ip call log_info(msg, mod='gray_core', proc='gray_main') 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 if(ip > 1) then du1 = zero gri = zero ggri = zero if(ip == params%raytracing%ipass) cpl = [zero, zero] ! no successive passes end if ! =========== beam loop BEGIN =========== call log_debug('beam loop start', mod='gray_core', proc='gray_main') beam_loop: do ib=1,nbeam_pass sox = -sox ! invert mode iox = 3-iox ! O-mode at ip=1,ib=1 index_rt = index_rt +1 child_index_rt = 2*index_rt + 1 ! * index_rt of O-mode child beam parent_index_rt = ceiling(index_rt / 2.0_wp_) - 1 ! * index_rt of parent beam call initbeam(index_rt,iroff,iboff,iwait,stv,jphi_beam, & pins_beam,currins_beam,dpdv_beam,jcd_beam) write(msg, '(" beam: ",g0," (",a1," mode)")') index_rt, mode(iox) call log_info(msg, mod='gray_core', proc='gray_main') if(iboff) then ! no propagation for current beam istop_pass = istop_pass +1 ! * +1 non propagating beam call log_info(" beam is off", mod='gray_core', proc='gray_main') cycle end if call vectinit(psjki,ppabs,ccci,tau0,alphaabs0,dids0,ccci0,iiv) if(ip == 1) then ! 1st pass igrad_b = params%raytracing%igrad ! * input value, igrad_b=0 from 2nd pass tau1 = zero ! * tau from previous passes etau1 = one cpl1 = one ! * coupling from previous passes lgcpl1 = zero p0ray = p0jk ! * initial beam power call ic_gb(params, anv0, ak0, yw, ypw, stv, xc, du1, & gri, ggri, index_rt, results%tables) ! * initial conditions do jk=1,params%raytracing%nray zzm = yw(3,jk)*0.01_wp_ rrm = sqrt(yw(1,jk)*yw(1,jk)+yw(2,jk)*yw(2,jk))*0.01_wp_ if (data%equilibrium%limiter%contains(rrm, zzm)) then ! * start propagation in/outside vessel? iow(jk) = 1 ! + inside else iow(jk) = 0 ! + outside end if end do else ! 2nd+ passes 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 if(params%raytracing%nray>1 .and. all(.not.iwait)) & ! iproj=0 ==> nfilp=8 call store_beam_shape(params%raytracing, results%tables, stv, yw) ! ======= propagation loop BEGIN ======= call log_debug(' propagation loop start', mod='gray_core', proc='gray_main') propagation_loop: do i=1,params%raytracing%nstep ! advance one step with "frozen" grad(S_I) do jk=1,params%raytracing%nray if(iwait(jk)) cycle ! jk ray is waiting for next pass stv(jk) = stv(jk) + params%raytracing%dst ! current ray step call rkstep(params, plasma, sox, bres, xgcn, yw(:,jk), & ypw(:,jk), gri(:,jk), ggri(:,:,jk), igrad_b) end do ! update position and grad if(igrad_b == 1) call gradi_upd(params%raytracing, yw, ak0, xc, du1, gri, ggri) error = 0 iopmin = 10 ! =========== rays loop BEGIN =========== rays_loop: do jk=1,params%raytracing%nray if(iwait(jk)) cycle ! jk ray is waiting for next pass ! compute derivatives with updated gradient and local plasma values xv = yw(1:3,jk) anv = yw(4:6,jk) call ywppla_upd(params, plasma, 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 global error code and print message if(ierrn/=0) then error = ior(error,ierrn) call print_err_raytracing(ierrn, i, anpl) end if ! check entrance/exit plasma/wall zzm = xv(3)*0.01_wp_ rrm = sqrt(xv(1)*xv(1)+xv(2)*xv(2))*0.01_wp_ ins_pl = (psinv>=zero .and. psinv continue current pass if(cpl(iox) < etaucr) then ! + IF low coupled power for current mode => de-activate derived rays call turnoffray(jk,ip+1,npass,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) if(jk == 1) then 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') write (msg,'(" polarisation: ψ=",g0.5,"°, χ=",g0.5,"°")') & psipol, chipol call log_debug(msg, mod='gray_core', proc='gray_main') psipv(index_rt) = psipol ! + polarization angles at plasma boundary for central ray chipv(index_rt) = chipol end if else if(iop(jk) > 2) then ! * 2nd entrance on 1st pass / entrance on 2nd+ pass => end of current pass for ray jk igrad_b = 0 ! + switch to ray-tracing iwait(jk) = .true. ! + stop advancement and H&CD computation for current ray if(ip < params%raytracing%ipass) then ! + not last pass yynext(:,jk,index_rt) = yw0(:,jk) ! . copy starting coordinates yypnext(:,jk,index_rt) = ypw0(:,jk) ! for next pass from last step 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 call turnoffray(jk,ip+1,npass,2*ib-1,iroff) if(cpl(1).le.comp_tiny) cpl(1)=etaucr end if if(cpl(2) < etaucr) then ! . low coupled power for X-mode => de-activate derived rays call turnoffray(jk,ip+1,npass,2*ib,iroff) if(cpl(2).le.comp_tiny) cpl(2)=etaucr end if taus(jk,child_index_rt:child_index_rt+1) = tau1(jk) + tau0(jk) ! . starting tau for next O-mode pass cpls(jk,child_index_rt) = cpl1(jk) * cpl(1) ! . cumulative coupling for next O-mode pass cpls(jk,child_index_rt+1) = cpl1(jk) * cpl(2) ! . cumulative coupling for next X-mode pass if(jk == 1) then ! . polarization angles at plasma boundary for central ray psipv(child_index_rt:child_index_rt+1) = psipol chipv(child_index_rt:child_index_rt+1) = chipol end if else ! * 1st entrance on 2nd+ pass (ray hasn't entered in plasma since end of previous pass) => continue current pass cpl = [zero, zero] end if end if else if(ext_pl) then ! ray exits plasma write (msg, '(" ray ", g0, " left plasma")') jk call log_debug(msg, mod='gray_core', proc='gray_main') call plasma_out(jk,xv,anv,bres,sox,iop,ext,eyt) end if if(ent_wl) then ! ray enters vessel iow(jk)=iow(jk)+1 ! * out->in else if(ext_wl) then ! ray exit vessel call wall_out(jk, data%equilibrium%limiter, ins_pl, xv, anv, & params%raytracing%dst, bres, sox, psipol, chipol, & iow, iop, ext, eyt) yw(:,jk) = [xv, anv] ! * updated coordinates (reflected) igrad_b = 0 ! * switch to ray-tracing call ywppla_upd(params, plasma, 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 error = ior(error,ierrn) call print_err_raytracing(ierrn, i, anpl) end if if(jk == 1 .and. ip == 1) then ! * 1st pass, polarization angles at reflection for central ray psipv(index_rt) = psipol chipv(index_rt) = chipol end if if(ins_pl) then ! * plasma-wall overlapping => wall+plasma crossing => end of current pass iwait(jk) = .true. ! + stop advancement and H&CD computation for current ray if(ip < params%raytracing%ipass) then ! + not last pass yynext(:,jk,index_rt) = [xv, anv] ! . starting coordinates 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 call plasma_in(jk, xv, anv, bres, sox, cpl, psipol, chipol, & ! . ray re-enters plasma after reflection iop, ext, eyt, perfect=.false.) if(cpl(1) < etaucr) then ! . low coupled power for O-mode? => de-activate derived rays call turnoffray(jk,ip+1,npass,2*ib-1,iroff) if(cpl(1).le.comp_tiny) cpl(1)=etaucr end if if(cpl(2) < etaucr) then ! . low coupled power for X-mode? => de-activate derived rays call turnoffray(jk,ip+1,npass,2*ib,iroff) if(cpl(2).le.comp_tiny) cpl(2)=etaucr end if taus(jk,child_index_rt:child_index_rt+1) = tau1(jk) + tau0(jk) ! . starting tau for next O-mode pass cpls(jk,child_index_rt) = cpl1(jk) * cpl(1) ! . cumulative coupling for next O-mode pass cpls(jk,child_index_rt+1) = cpl1(jk) * cpl(2) ! . cumulative coupling for next X-mode pass if(jk == 1) then ! + polarization angles at plasma boundary for central ray psipv(child_index_rt:child_index_rt+1) = psipol chipv(child_index_rt:child_index_rt+1) = chipol end if end if end if end if iopmin = min(iopmin,iop(jk)) if(ip < params%raytracing%ipass) then ! not last pass yw0(:,jk) = yw(:,jk) ! * store current coordinates in case ypw0(:,jk) = ypw(:,jk) ! current pass ends on next step end if ! 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. & (tau1(jk)+tau0(jk)+lgcpl1(jk))<=taucr .and. .not.iwait(jk)) then ! H&CD computation check tekev = plasma%temp(psinv) block complex(wp_) :: Npr_warm call alpha_effj(params%ecrh_cd, plasma, 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 if(ierrhcd /= 0) then error = ior(error, ierrhcd) call print_err_ecrh_cd(ierrhcd, i, Npr_warm, alpha) end if end block 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 psjki(jk,i) = psinv ! 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 dids = didp * pow * alpha ! current driven: dI/ds = dI/dP⋅dP/ds = dI/dP⋅P⋅α ! 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) if(iwait(jk)) then ! copy values from last pass for inactive ray 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) else call store_ray_data(params, results%tables, & i, jk, stv(jk), p0ray(jk), xv, psinv, btot, bv, ak0, & anpl, anpr, anv, anprim, dens, tekev, alpha, tau, dids, & 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) end if end do rays_loop ! ============ rays loop END ============ if(i==params%raytracing%nstep) then ! step limit reached? do jk=1,params%raytracing%nray if(iop(jk)<3) call turnoffray(jk,ip,npass,ib,iroff) ! * ray hasn't exited+reentered the plasma by last step => stop ray end do end if ! 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 store_beam_shape(params%raytracing, results%tables, stv, yw) end if ! test whether further trajectory integration is unnecessary call vmaxmin(tau1+tau0+lgcpl1,params%raytracing%nray,taumn,taumx) ! test on tau + coupling if(is_critical(error)) then ! stop propagation for current beam istop_pass = istop_pass +1 ! * +1 non propagating beam if(ip < params%raytracing%ipass) call turnoffray(0,ip,npass,ib,iroff) ! * de-activate derived beams exit else if(all(iwait)) then ! all rays in current beam are waiting for next pass => do not increase istop_pass exit end if end do propagation_loop call log_debug(' propagation loop end', mod='gray_core', proc='gray_main') ! ======== propagation loop END ======== ! print all ray positions in local reference system if(params%raytracing%nray > 1 .and. all(.not.iwait)) & call store_beam_shape(params%raytracing, results%tables, & stv, yw, full=.true.) ! =========== post-proc BEGIN =========== ! compute total absorbed power and driven current for current beam if(i>params%raytracing%nstep) i=params%raytracing%nstep pabs_beam = sum(ppabs(:,i)) icd_beam = sum(ccci(:,i)) call vmaxmin(tau0,params%raytracing%nray,taumn,taumx) ! taumn,taumx for print if (params%equilibrium%iequil /= EQ_VACUUM) then ! compute power and current density profiles for all rays call spec(psjki, ppabs, ccci, iiv, pabs_beam, & icd_beam, dpdv_beam, jphi_beam, jcd_beam, & pins_beam, currins_beam) end if 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 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(:,child_index_rt)/cpl1, MASK=iop > 2) / & sum(p0ray * exp(-tau0), MASK=iop > 2) ! * average O-mode coupling for next beam (on active rays) cpl_beam2 = one - cpl_beam1 ! * average X-mode coupling for next beam if(iop(1) > 2) then ! * central ray O/X-mode coupling for next beam cpl_cbeam1 = cpls(1,child_index_rt)/cpl1(1) cpl_cbeam2 = one - cpl_cbeam1 end if else ! last pass OR no ray re-entered plasma cpl_beam1 = zero cpl_beam2 = zero end if ! print final results for pass on screen call log_info(' partial results:', mod='gray_core', proc='gray_main') write(msg, '(" final step:", x,a,g0, x,a,g0.4)') 'n=', i, '(s, ct, Sr)=', stv(1) 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') write(msg, '(3x,a,g0.3," kA")') 'current drive: I=', icd_beam*1.0e3_wp_ call log_info(msg, mod='gray_core', proc='gray_main') if(ip < params%raytracing%ipass) then 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 call log_info(msg, mod='gray_core', proc='gray_main') end if end if if (params%equilibrium%iequil /= EQ_VACUUM) then call store_ec_profiles( & results%tables%ec_profiles, rhop_tab, rhot_tab, jphi_beam, & jcd_beam, dpdv_beam, currins_beam, pins_beam, ip) 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 if (results%tables%summary%active) then call results%tables%summary%append([ & wrap(icd_beam), wrap(pabs_beam), wrap(jphip), wrap(dpdvp), & wrap(rhotj), wrap(rhotjava), wrap(rhotp), wrap(rhotpav), & wrap(drhotjava), wrap(drhotpav), wrap(ratjamx), wrap(ratjbmx), & wrap(stv(1)), wrap(psipv(index_rt)), wrap(chipv(index_rt)), & wrap(index_rt), wrap(jphimx), wrap(dpdvmx), wrap(drhotj), & wrap(drhotp), wrap(sum(p0ray)), wrap(cpl_beam1), wrap(cpl_beam2)]) ! *print 0D results for current beam end if end if ! ============ post-proc END ============ end do beam_loop call log_debug('beam loop end', mod='gray_core', proc='gray_main') ! ============ beam loop END ============ ! ======= cumulative prints BEGIN ======= results%pabs = results%pabs + sum(pabs_pass) ! *final results (O+X) [gray_main output] results%icd = results%icd + sum(icd_pass) ! 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') ! ======== cumulative prints END ======== if(istop_pass == nbeam_pass) exit ! no active beams end do main_loop call log_debug('pass loop end', mod='gray_core', proc='gray_main') ! ============ main loop END ============ ! Free memory call dealloc_surfvec call dealloc_pec end block end associate end subroutine gray_main subroutine vectinit(psjki,ppabs,ccci,tau0,alphaabs0,dids0,ccci0,iiv) use const_and_precisions, only : zero ! arguments real(wp_), dimension(:,:), intent(out) :: psjki,ppabs,ccci real(wp_), dimension(:), intent(out) :: tau0,alphaabs0,dids0,ccci0 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_ psjki = zero ppabs = zero ccci = zero tau0 = zero alphaabs0 = zero dids0 = zero ccci0 = zero iiv = 1 end subroutine vectinit subroutine ic_gb(params, anv0c, ak0, ywrk0, ypwrk0, & stv, xc0, du10, gri, ggri, index_rt, & tables) ! beam tracing initial conditions igrad=1 ! !!!!!! check ray tracing initial conditions igrad=0 !!!!!! use const_and_precisions, only : zero,one,pi,half,two,degree,ui=>im use gray_params, only : gray_parameters, gray_tables use types, only : table use gray_tables, only : store_ray_data ! subroutine arguments type(gray_parameters), intent(in) :: params real(wp_), dimension(3), intent(in) :: anv0c real(wp_), intent(in) :: ak0 real(wp_), dimension(6, params%raytracing%nray), intent(out) :: ywrk0, ypwrk0 real(wp_), dimension(params%raytracing%nrayth * params%raytracing%nrayr), intent(out) :: stv real(wp_), dimension(3, params%raytracing%nrayth, params%raytracing%nrayr), intent(out) :: xc0, du10 real(wp_), dimension(3, params%raytracing%nray), intent(out) :: gri real(wp_), dimension(3, 3, params%raytracing%nray), intent(out) :: ggri integer, intent(in) :: index_rt ! tables for storing initial rays data type(gray_tables), intent(inout), optional :: tables ! local variables real(wp_), dimension(3) :: xv0c real(wp_) :: wcsi,weta,rcicsi,rcieta,phiw,phir 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(params%raytracing%nrayr) :: uj real(wp_), dimension(params%raytracing%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 ! 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) 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 ! 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 phiwrad = phiw*degree phirrad = phir*degree csphiw = cos(phiwrad) snphiw = sin(phiwrad) !csphir = cos(phirrad) !snphir = sin(phirrad) wwcsi = two/(ak0*wcsi**2) wweta = two/(ak0*weta**2) if(phir/=phiw) then 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)) phic = half*atan(ts/tc) ddd = dk*cos(2*(phirrad+phic)) - ui*dw*cos(2*(phiwrad+phic)) sss = sk - ui*sw qi1 = half*(sss + ddd) qi2 = half*(sss - ddd) rci1 = real(qi1) rci2 = real(qi2) ww1 = -imag(qi1) ww2 = -imag(qi2) else rci1 = rcicsi rci2 = rcieta ww1 = wwcsi ww2 = wweta phic = -phiwrad qi1 = rci1 - ui*ww1 qi2 = rci2 - ui*ww2 end if !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) qqxx = qi1*cos(phic)**2 + qi2*sin(phic)**2 qqyy = qi1*sin(phic)**2 + qi2*cos(phic)**2 qqxy = -(qi1 - qi2)*sin(phic)*cos(phic) wwxx = -imag(qqxx) wwyy = -imag(qqyy) wwxy = -imag(qqxy) rcixx = real(qqxx) rciyy = real(qqyy) rcixy = real(qqxy) 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 dqqxy = -(dqi1 - dqi2)*sin(phic)*cos(phic) d2qqxx = d2qi1*cos(phic)**2 + d2qi2*sin(phic)**2 d2qqyy = d2qi1*sin(phic)**2 + d2qi2*cos(phic)**2 d2qqxy = -(d2qi1 - d2qi2)*sin(phic)*cos(phic) 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) if(params%raytracing%nrayr > 1) then dr = params%raytracing%rwmax / (params%raytracing%nrayr - 1) else dr = 1 end if ddfu = 2*dr**2/ak0 do j = 1, params%raytracing%nrayr uj(j) = real(j-1, wp_) end do da=2*pi/params%raytracing%nrayth do k=1,params%raytracing%nrayth alfak = (k-1)*da sna(k) = sin(alfak) csa(k) = cos(alfak) end do ! central ray 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 do k=1,params%raytracing%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 ! loop on rays jk>1 j=2 k=0 do jk = 2, params%raytracing%nray k=k+1 if(k > params%raytracing%nrayth) then j=j+1 k=1 end if !csiw=u*dcsiw !etaw=u*detaw !csir=x0t*csphir+y0t*snphir !etar=-x0t*snphir+y0t*csphir 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 z0t = 0 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 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 ! integration variable select case(params%raytracing%idst) case(0) ! optical path: s stv(jk) = 0 denom = an0 case(1) ! "time": c⋅t stv(jk) = 0 denom = one case(2) ! real eikonal: S_R stv(jk) = half*(rcixx*x0t**2 + rciyy*y0t**2) + rcixy*x0t*y0t denom = an20 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) ! save step "zero" data if (present(tables)) & call store_ray_data(params, tables, & i=0, jk=jk, s=stv(jk), P0=one, pos=xc0(:,k,j), & psi_n=-one, B=zero, b_n=[zero,zero,zero], k0=ak0, & N_pl=zero, N_pr=zero, N=ywrk0(:, jk), N_pr_im=zero, & n_e=zero, T_e=zero, alpha=zero, tau=zero, dIds=zero, & nhm=0, nhf=0, iokhawa=0, index_rt=index_rt, & lambda_r=ddr, lambda_i=ddi, X=zero, Y=zero, & grad_X=[zero,zero,zero]) end do end subroutine ic_gb subroutine rkstep(params, plasma, sox, bres, xgcn, y, yp, dgr, ddgr, igrad) ! Runge-Kutta integrator use gray_params, only : gray_parameters use gray_plasma, only : abstract_plasma ! subroutine arguments type(gray_parameters), intent(in) :: params class(abstract_plasma), intent(in) :: plasma real(wp_), intent(in) :: bres, xgcn real(wp_), intent(inout) :: y(6) real(wp_), intent(in) :: yp(6) real(wp_), intent(in) :: dgr(3) real(wp_), intent(in) :: ddgr(3,3) integer, intent(in) :: igrad,sox ! local variables real(wp_), dimension(6) :: k1, k2, k3, k4 associate (h => params%raytracing%dst) k1 = yp k2 = f(y + k1*h/2) k3 = f(y + k2*h/2) k4 = f(y + k3*h) y = y + h * (k1/6 + k2/3 + k3/3 + k4/6) end associate contains function f(y) real(wp_), intent(in) :: y(6) real(wp_) :: f(6) call rhs(params, plasma, sox, bres, xgcn, y, dgr, ddgr, f, igrad) end function end subroutine rkstep subroutine rhs(params, plasma, sox, bres, xgcn, y, dgr, ddgr, dery, igrad) ! Compute right-hand side terms of the ray equations (dery) ! used in R-K integrator use gray_params, only : gray_parameters use gray_plasma, only : abstract_plasma ! subroutine arguments type(gray_parameters), intent(in) :: params class(abstract_plasma), intent(in) :: plasma real(wp_), intent(in) :: y(6) real(wp_), intent(in) :: bres, xgcn real(wp_), intent(in) :: dgr(3) real(wp_), intent(in) :: ddgr(3,3) real(wp_), intent(out) :: dery(6) integer, intent(in) :: igrad, sox ! local variables real(wp_) :: dens,btot,xg,yg real(wp_), dimension(3) :: xv,anv,bv,derxg,deryg real(wp_), dimension(3,3) :: derbv xv = y(1:3) call plas_deriv(params, plasma, xv, bres, xgcn, dens, btot, & bv, derbv, xg, yg, derxg, deryg) anv = y(4:6) call disp_deriv(params, anv, sox, xg, yg, derxg, deryg, & bv, derbv, dgr, ddgr, igrad, dery) end subroutine rhs subroutine ywppla_upd(params, plasma, xv, anv, dgr, ddgr, sox, bres, xgcn, dery, & psinv, dens, btot, bv, xg, yg, derxg, anpl, anpr, & ddr, ddi, dersdst, derdnm, error, igrad) ! Compute right-hand side terms of the ray equations (dery) ! used after full R-K step and grad(S_I) update use gray_errors, only : raise_error, large_npl use gray_params, only : gray_parameters use gray_plasma, only : abstract_plasma ! subroutine rguments type(gray_parameters), intent(in) :: params class(abstract_plasma), intent(in) :: plasma real(wp_), intent(in) :: xv(3), anv(3) real(wp_), intent(in) :: dgr(3) real(wp_), intent(in) :: ddgr(3,3) integer, intent(in) :: sox real(wp_), intent(in) :: bres,xgcn real(wp_), intent(out) :: dery(6) real(wp_), intent(out) :: psinv, dens, btot, xg, yg, anpl, anpr real(wp_), intent(out) :: ddr, ddi, dersdst, derdnm real(wp_), intent(out) :: bv(3) integer, intent(out) :: error real(wp_), intent(out) :: derxg(3) integer, intent(in) :: igrad ! local variables real(wp_), dimension(3) :: deryg real(wp_), dimension(3,3) :: derbv real(wp_), parameter :: anplth1 = 0.99_wp_, anplth2 = 1.05_wp_ call plas_deriv(params, plasma, xv,bres,xgcn,dens,btot,bv,derbv,xg,yg,derxg,deryg,psinv) call disp_deriv(params, anv,sox,xg,yg,derxg,deryg,bv,derbv,dgr,ddgr,igrad, & dery,anpl,anpr,ddr,ddi,dersdst,derdnm) 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 end if end subroutine ywppla_upd subroutine gradi_upd(params, ywrk,ak0,xc,du1,gri,ggri) use const_and_precisions, only : zero, half use gray_params, only : raytracing_parameters ! subroutine parameters type(raytracing_parameters), intent(in) :: params real(wp_), intent(in) :: ak0 real(wp_), dimension(6,params%nray), intent(in) :: ywrk real(wp_), dimension(3,params%nrayth,params%nrayr), intent(inout) :: xc, du1 real(wp_), dimension(3,params%nray), intent(out) :: gri real(wp_), dimension(3,3,params%nray), intent(out) :: ggri ! local variables real(wp_), dimension(3,params%nrayth,params%nrayr) :: xco, du1o integer :: jk, j, jm, jp, k, km, kp real(wp_) :: ux, uxx, uxy, uxz, uy, uyy, uyz, uz, uzz real(wp_) :: dr, 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,params%nrayr do k=1,params%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,params%nrayth if(k == 1) then km = params%nrayth else km = k-1 end if if(k == params%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 ! compute grad u1 and grad(S_I) for all the other rays if (params%nrayr > 1) then dr = params%rwmax / (params%nrayr - 1) else dr = 1 end if dfuu=2*dr**2/ak0 jm=1 j=2 k=0 dffiu = dfuu do jk=2,params%nray k=k+1 if(k > params%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=params%nrayth else if (k == params%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,params%nray k=k+1 if(k > params%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=params%nrayth else if (k == params%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 end subroutine gradi_upd 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) ! 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 ! 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 subroutine plas_deriv(params, plasma, xv, bres, xgcn, dens, btot, bv, & derbv, xg, yg, derxg, deryg, psinv_opt) use const_and_precisions, only : zero, cm use gray_params, only : gray_parameters, EQ_VACUUM use gray_plasma, only : abstract_plasma use equilibrium, only : psia, pol_flux, pol_curr, sgnbphi ! subroutine arguments type(gray_parameters), intent(in) :: params class(abstract_plasma), intent(in) :: plasma 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 real(wp_), dimension(3,3), intent(out) :: derbv real(wp_), optional, intent(out) :: psinv_opt ! local variables integer :: jv real(wp_) :: xx,yy,zz real(wp_) :: psinv real(wp_) :: csphi,drrdx,drrdy,dphidx,dphidy,rr,rr2,rrm,snphi,zzm real(wp_), dimension(3) :: dbtot,bvc real(wp_), dimension(3,3) :: dbvcdc,dbvdc,dbv real(wp_) :: brr,bphi,bzz,dxgdpsi real(wp_) :: dpsidr,dpsidz,ddpsidrr,ddpsidzz,ddpsidrz,fpolv,dfpv,ddensdpsi xg = zero yg = 99._wp_ psinv = -1._wp_ dens = zero btot = zero derxg = zero deryg = zero bv = zero derbv = zero if (params%equilibrium%iequil == EQ_VACUUM) then ! copy optional output if (present(psinv_opt)) psinv_opt = psinv return end if dbtot = zero dbv = zero dbvcdc = zero dbvcdc = zero dbvdc = zero xx = xv(1) yy = xv(2) zz = xv(3) ! cylindrical coordinates rr2 = xx**2 + yy**2 rr = sqrt(rr2) csphi = xx/rr snphi = yy/rr bv(1) = -snphi*sgnbphi bv(2) = csphi*sgnbphi ! convert from cm to meters zzm = 1.0e-2_wp_*zz rrm = 1.0e-2_wp_*rr call pol_flux(rrm, zzm, psinv, dpsidr, dpsidz, ddpsidrr, ddpsidzz, ddpsidrz) call pol_curr(psinv, fpolv, dfpv) ! copy optional output if (present(psinv_opt)) psinv_opt = psinv ! compute yg and derivative if(psinv < zero) then bphi = fpolv/rrm btot = abs(bphi) yg = btot/bres return end if ! compute xg and derivative call plasma%density(psinv, dens, ddensdpsi) xg = xgcn*dens dxgdpsi = xgcn*ddensdpsi ! B = f(psi)/R e_phi+ grad psi x e_phi/R bphi = fpolv/rrm brr = -dpsidz*psia/rrm bzz = +dpsidr*psia/rrm ! bvc(i) = B_i in cylindrical coordinates bvc = [brr, bphi, bzz] ! bv(i) = B_i in cartesian coordinates bv(1)=bvc(1)*csphi - bvc(2)*snphi bv(2)=bvc(1)*snphi + bvc(2)*csphi bv(3)=bvc(3) ! dbvcdc(iv,jv) = d Bcil(iv) / dxvcil(jv) dbvcdc(1,1) = -ddpsidrz*psia/rrm - brr/rrm dbvcdc(1,3) = -ddpsidzz*psia/rrm dbvcdc(2,1) = dfpv*dpsidr/rrm - bphi/rrm dbvcdc(2,3) = dfpv*dpsidz/rrm dbvcdc(3,1) = +ddpsidrr*psia/rrm - bzz/rrm dbvcdc(3,3) = +ddpsidrz*psia/rrm ! dbvdc(iv,jv) = d Bcart(iv) / dxvcil(jv) 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) dbvdc(3,2) = dbvcdc(3,2) ! = 0 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 ! dbv(iv,jv) = d Bcart(iv) / dxvcart(jv) dbv(:,1) = drrdx*dbvdc(:,1) + dphidx*dbvdc(:,2) dbv(:,2) = drrdy*dbvdc(:,1) + dphidy*dbvdc(:,2) dbv(:,3) = dbvdc(:,3) ! B magnitude and derivatives btot = norm2(bv) ! dbtot(i) = d |B| / dxvcart(i) dbtot = matmul(bv, dbv)/btot yg = btot/Bres ! convert spatial derivatives from dummy/m -> dummy/cm ! to be used in rhs ! bv(i) = B_i / B ; derbv(i,j) = d (B_i / B) /d x,y,z deryg = 1.0e-2_wp_*dbtot/Bres bv = bv/btot do jv=1,3 derbv(:,jv) = cm * (dbv(:,jv) - bv(:)*dbtot(jv))/btot end do derxg(1) = cm * drrdx*dpsidr*dxgdpsi derxg(2) = cm * drrdy*dpsidr*dxgdpsi derxg(3) = cm * dpsidz *dxgdpsi end subroutine plas_deriv subroutine disp_deriv(params, 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. use const_and_precisions, only : zero, one, half, two use gray_params, only : gray_parameters, STEP_ARCLEN, & STEP_TIME, STEP_PHASE ! subroutine arguments ! Inputs type(gray_parameters), intent(in) :: params ! 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 integer, intent(in) :: sox ! 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 ! local variables real(wp_) :: gr2, yg2, anpl2, del, dnl, duh, dan2sdnpl, an2, an2s real(wp_) :: dan2sdxg, dan2sdyg, denom real(wp_) :: derdom, dfdiadnpl, dfdiadxg, dfdiadyg, fdia, bdotgr real(wp_), dimension(3) :: derdxv, danpldxv, derdnv, dbgr an2 = dot_product(anv, anv) ! N² anpl = dot_product(anv, bv) ! N∥ = N̅⋅B̅ ! 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∥ an2s = one dan2sdxg = zero dan2sdyg = zero dan2sdnpl = zero del = zero fdia = zero dfdiadnpl = zero dfdiadxg = zero dfdiadyg = zero if(xg > zero) then ! Derivatives of the cold plasma refractive index ! ! N²s = 1 - X - XY²⋅(1 + N∥² ± √Δ)/[2(1 - X - Y²)] ! ! 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) ! ∂(N²s)/∂X ! Note: this term is nonzero for X=0, but it multiplies terms ! proportional to X or ∂X/∂ψ which are zero outside the plasma. dan2sdxg = - half*yg2*(one - yg2)*(one + anpl2 + sox*del)/duh**2 & + sox*xg*anpl2/(del*duh) - one ! ∂(N²s)/∂Y dan2sdyg = - xg*yg*(one - xg)*(one + anpl2 + sox*del)/duh**2 & + two*sox*xg*(one - xg)*anpl2/(yg*del*duh) ! ∂(N²s)/∂N∥ dan2sdnpl = - xg*yg2*anpl/duh & - sox*xg*anpl*(two*(one - xg) - yg2*dnl)/(del*duh) if(igrad > 0) then ! 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 = 24.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 end if end if ! ∇(N∥) = ∇(N̅⋅b̅) = N̅⋅∇b̅ danpldxv = matmul(anv, derbv) ! ∂Λ/∂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 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 (params%raytracing%idst) case (STEP_ARCLEN) ! σ=s, arclength ! denom = |∂Λ/∂N̅| ! Proof: Normalising ∂Λ/∂N̅ (∥ to the group velocity) ! simply makes σ the arclength parameter s. denom = derdnm case (STEP_TIME) ! σ=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̅ / (-ω⋅∂Λ/∂ω) ! denom = -derdom case (STEP_PHASE) ! 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̅ ! denom = dot_product(anv, derdnv) end select ! 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 ! r.h.s. vector dery(1:3) = derdnv(:)/denom ! +∂Λ/∂N̅ dery(4:6) = -derdxv(:)/denom ! -∂Λ/∂x̅ 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 end subroutine disp_deriv subroutine alpha_effj(params, plasma, psinv, X, Y, density, temperature, & k0, Bres, derdnm, Npl, Npr_cold, sox, Npr, & alpha, dIdp, nhmin, nhmax, iokhawa, error) ! Computes the absorption coefficient and effective current use const_and_precisions, only : pi, mc2=>mc2_ use gray_params, only : ecrh_cd_parameters use gray_plasma, only : abstract_plasma use equilibrium, only : sgnbphi use dispersion, only : harmnumber, warmdisp use eccd, only : setcdcoeff, eccdeff, fjch0, fjch, fjncl use gray_errors, only : negative_absorption, raise_error use magsurf_data, only : fluxval ! subroutine arguments ! Inputs ! ECRH & CD parameters type(ecrh_cd_parameters), intent(in) :: params ! plasma object class(abstract_plasma), intent(in) :: plasma ! poloidal flux ψ real(wp_), intent(in) :: psinv ! CMA diagram variables: X=(ω_pe/ω)², Y=ω_ce/ω real(wp_), intent(in) :: X, Y ! densityity [10¹⁹ m⁻³], temperature [keV] real(wp_), intent(in) :: density, temperature ! vacuum wavenumber k₀=ω/c, resonant B field real(wp_), intent(in) :: k0, Bres ! group velocity: |∂Λ/∂N̅| where Λ=c²/ω²⋅D_R real(wp_), intent(in) :: derdnm ! refractive index: N∥, N⊥ (cold dispersion) real(wp_), intent(in) :: Npl, Npr_cold ! sign of polarisation mode: -1 ⇒ O, +1 ⇒ X integer, intent(in) :: sox ! Outputs ! orthogonal refractive index N⊥ (solution of the warm dispersion) complex(wp_), intent(out) :: Npr ! absorption coefficient, current density real(wp_), intent(out) :: alpha, dIdp ! 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 integer :: nlarmor, ithn, ierrcd real(wp_) :: mu, rbn, rbx real(wp_) :: zeff, cst2, bmxi, bmni, fci real(wp_), dimension(:), allocatable :: eccdpar real(wp_) :: effjcd, effjcdav, Btot complex(wp_) :: e(3) alpha = 0 Npr = 0 dIdp = 0 nhmin = 0 nhmax = 0 iokhawa = 0 error = 0 ! Absorption computation ! Skip when the temperature is zero (no plasma) if (temperature <= 0) return ! Skip when the lowest resonant harmonic number is zero mu = mc2/temperature ! μ=(mc²/kT) call harmnumber(Y, mu, Npl**2, params%iwarm == 1, nhmin, nhmax) if (nhmin <= 0) return ! 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 ! Compute α from the solution of the dispersion relation ! The absoption coefficient is defined as ! ! α = 2 Im(k̅)⋅s̅ ! ! where s̅ = v̅_g/|v_g|, the direction of the energy flow. ! Since v̅_g = - ∂Λ/∂N̅ / ∂Λ/∂ω, using the cold dispersion ! relation, we have that ! ! s̅ = ∂Λ/∂N̅ / |∂Λ/∂N̅| ! = [2N̅ - ∂(N²s)/∂N∥ b̅ ! + ½(b̅⋅∇S_I)² ∂³(N²s)/∂N∥³ b̅] / |∂Λ/∂N̅| ! ! Assuming Im(k∥)=0: ! ! α = 4 Im(k⊥)⋅N⊥ / |∂Λ/∂N̅| ! block real(wp_) :: k_im k_im = k0 * Npr%im ! imaginary part of k⊥ alpha = 4 * k_im*Npr_cold / derdnm end block if (alpha < 0) then error = raise_error(error, negative_absorption) return end if ! Current drive computation if (params%ieccd <= 0) return zeff = plasma%zeff(psinv) ithn = 1 if (nlarmor > nhmin) ithn = 2 rhop = sqrt(psinv) call fluxval(rhop, rri=rrii, rbav=rbavi, bmn=bmni, bmx=bmxi, fc=fci) Btot = Y*Bres rbn = Btot/bmni rbx = Btot/bmxi select case(params%ieccd) case(1) ! Cohen model 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) case(2) ! No trapping call setcdcoeff(zeff, cst2, eccdpar) call eccdeff(Y, Npl, Npr%re, density, mu, e, nhmin, nhmax, & ithn, cst2, fjch0, eccdpar, effjcd, iokhawa, ierrcd) case default ! Neoclassical model 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) end select error = error + ierrcd if (allocated(eccdpar)) deallocate(eccdpar) ! current drive efficiency [A⋅m/W] effjcdav = rbavi*effjcd dIdp = sgnbphi * effjcdav / (2*pi*rrii) end subroutine alpha_effj end module gray_core