gray/src/multipass.f90

302 lines
16 KiB
Fortran
Raw Normal View History

module multipass
use const_and_precisions, only : wp_, zero, half, one, degree, czero
use beamdata, only : dst, nray
use gray_params, only : ipass
use polarization, only : pol_limit, stokes_ce, polellipse
use reflections, only : wall_refl
use equilibrium, only : bfield
implicit none
2019-03-28 10:50:28 +01:00
integer, save :: nbeam_max ! max n of beams active at a time
integer, save :: nbeam_tot ! total n of beams
contains
fix coupling for subsequent beams In situations when multiple beams are traced, either when allowing multiple plasma crossings (raytracing.ipass > 0) or the initial polarisation is mixed (raytracing.ipol == .true.), the couplings of all but the first beam (with least index_rt) were invalid. The bug is due to the re-use of the psipol,chipol variables as the beams are traced sequentially over the beam_loop. For the first beam being traced the psipol,chipol are correctly initialised to the user-defined value and the resulting coupling is correct. However, in each subsequent beam the values were not set to those of the parent beam (or to the user-defined value in the case of the first X mode beam), but to those of the previous beams (current index_rt - 1). This change repurposes the psipv,chipv arrays to store the polarisation of the parent beams, including the initial user-defined value and makes plasma_in always use these to compute the coupling. In addition, in the case the polarisation is not immediately known (i.e. if raytracing.ipol == .false.), this change postpones the computation of the Jones vector (ext, eyt) from the launch point, if the magnetic equilibrium is available, to when the ray actually crosses the plasma boundary. The original code, besides being strictly incorrect, can lead to non-negligible alterations to the coupling. This change also mean: 1. most of the functionality of `set_pol` has been merged with `plasma_in` 2. the polarisation is undefined and the Jones vector is set to the placeholder value [1, 0] till `plasma_im` is called Finally, `polarcold` is removed because it's unused.
2024-03-02 17:32:03 +01:00
subroutine plasma_in(i, x, N, Bres, sox, cpl, psi, chi, iop, ext, eyt, perfect)
! Computes the ray polarisation and power couplings when it enteres the plasma
use const_and_precisions, only: cm
! subroutine arguments
integer, intent(in) :: i ! ray index
real(wp_), intent(in) :: x(3), N(3) ! position, refactive index
real(wp_), intent(in) :: Bres ! resonant B field
integer, intent(in) :: sox ! sign of polarisation mode: -1 ⇒ O, +1 ⇒ X
real(wp_), intent(out) :: cpl(2) ! power coupling vector (O, X)
real(wp_), intent(out) :: psi, chi ! polarisation ellipse angles
integer, intent(inout), pointer :: iop(:) ! inside/outside plasma flag
complex(wp_), intent(inout), pointer :: ext(:), eyt(:) ! ray polarisation vector (e_x, e_y)
logical, intent(in) :: perfect ! whether to assume perfect coupling
! local variables
real(wp_) :: R, z, cosphi, sinphi, B_phi, B_R, B_z
real(wp_) :: B(3)
real(wp_) :: c
complex(wp_) :: e_mode(2), e_ray(2)
! Update the inside/outside flag
iop(i) = iop(i) + 1
! Compute B in cartesian coordinates
R = norm2(x(1:2)) * cm
z = x(3) * cm
cosphi = x(1)/R * cm
sinphi = x(2)/R * cm
call bfield(R, z, B_R, B_z, B_phi)
B(1) = B_R*cosphi - B_phi*sinphi
B(2) = B_R*sinphi + B_phi*cosphi
B(3) = B_z
! Get the polarisation vector of the given mode
call pol_limit(N, B, Bres, sox, e_mode(1), e_mode(2))
if(i == 1) then
! For the central ray, compute the polarization ellipse
chi = asin(2*imag(e_mode(1) * conjg(e_mode(2)))) / 2
psi = atan2(real(2 * e_mode(1) * conjg(e_mode(2))), &
real(dot_product(e_mode, conjg(e_mode)))) / 2
! and convert from radians to degree
psi = psi / degree
chi = chi / degree
else
fix coupling for subsequent beams In situations when multiple beams are traced, either when allowing multiple plasma crossings (raytracing.ipass > 0) or the initial polarisation is mixed (raytracing.ipol == .true.), the couplings of all but the first beam (with least index_rt) were invalid. The bug is due to the re-use of the psipol,chipol variables as the beams are traced sequentially over the beam_loop. For the first beam being traced the psipol,chipol are correctly initialised to the user-defined value and the resulting coupling is correct. However, in each subsequent beam the values were not set to those of the parent beam (or to the user-defined value in the case of the first X mode beam), but to those of the previous beams (current index_rt - 1). This change repurposes the psipv,chipv arrays to store the polarisation of the parent beams, including the initial user-defined value and makes plasma_in always use these to compute the coupling. In addition, in the case the polarisation is not immediately known (i.e. if raytracing.ipol == .false.), this change postpones the computation of the Jones vector (ext, eyt) from the launch point, if the magnetic equilibrium is available, to when the ray actually crosses the plasma boundary. The original code, besides being strictly incorrect, can lead to non-negligible alterations to the coupling. This change also mean: 1. most of the functionality of `set_pol` has been merged with `plasma_in` 2. the polarisation is undefined and the Jones vector is set to the placeholder value [1, 0] till `plasma_im` is called Finally, `polarcold` is removed because it's unused.
2024-03-02 17:32:03 +01:00
psi = 0
chi = 0
end if
if (perfect) then
! Ignore the given vector and use the expected one
! Note: this will give 100% coupling to the current mode
ext(i) = e_mode(1)
eyt(i) = e_mode(2)
end if
fix coupling for subsequent beams In situations when multiple beams are traced, either when allowing multiple plasma crossings (raytracing.ipass > 0) or the initial polarisation is mixed (raytracing.ipol == .true.), the couplings of all but the first beam (with least index_rt) were invalid. The bug is due to the re-use of the psipol,chipol variables as the beams are traced sequentially over the beam_loop. For the first beam being traced the psipol,chipol are correctly initialised to the user-defined value and the resulting coupling is correct. However, in each subsequent beam the values were not set to those of the parent beam (or to the user-defined value in the case of the first X mode beam), but to those of the previous beams (current index_rt - 1). This change repurposes the psipv,chipv arrays to store the polarisation of the parent beams, including the initial user-defined value and makes plasma_in always use these to compute the coupling. In addition, in the case the polarisation is not immediately known (i.e. if raytracing.ipol == .false.), this change postpones the computation of the Jones vector (ext, eyt) from the launch point, if the magnetic equilibrium is available, to when the ray actually crosses the plasma boundary. The original code, besides being strictly incorrect, can lead to non-negligible alterations to the coupling. This change also mean: 1. most of the functionality of `set_pol` has been merged with `plasma_in` 2. the polarisation is undefined and the Jones vector is set to the placeholder value [1, 0] till `plasma_im` is called Finally, `polarcold` is removed because it's unused.
2024-03-02 17:32:03 +01:00
! Compute the power coupling with the current mode
e_ray = [ext(i), eyt(i)]
c = abs(dot_product(e_mode, e_ray))**2
! Store both O and X couplings, in this order
c = merge(c, 1-c, sox == -1)
cpl = [c, 1-c]
end subroutine plasma_in
fix coupling for subsequent beams In situations when multiple beams are traced, either when allowing multiple plasma crossings (raytracing.ipass > 0) or the initial polarisation is mixed (raytracing.ipol == .true.), the couplings of all but the first beam (with least index_rt) were invalid. The bug is due to the re-use of the psipol,chipol variables as the beams are traced sequentially over the beam_loop. For the first beam being traced the psipol,chipol are correctly initialised to the user-defined value and the resulting coupling is correct. However, in each subsequent beam the values were not set to those of the parent beam (or to the user-defined value in the case of the first X mode beam), but to those of the previous beams (current index_rt - 1). This change repurposes the psipv,chipv arrays to store the polarisation of the parent beams, including the initial user-defined value and makes plasma_in always use these to compute the coupling. In addition, in the case the polarisation is not immediately known (i.e. if raytracing.ipol == .false.), this change postpones the computation of the Jones vector (ext, eyt) from the launch point, if the magnetic equilibrium is available, to when the ray actually crosses the plasma boundary. The original code, besides being strictly incorrect, can lead to non-negligible alterations to the coupling. This change also mean: 1. most of the functionality of `set_pol` has been merged with `plasma_in` 2. the polarisation is undefined and the Jones vector is set to the placeholder value [1, 0] till `plasma_im` is called Finally, `polarcold` is removed because it's unused.
2024-03-02 17:32:03 +01:00
subroutine plasma_out(i,xv,anv,bres,sox,iop,ext,eyt) ! ray exits plasma
2019-03-28 10:50:28 +01:00
! arguments
integer, intent(in) :: i ! ray index
real(wp_), dimension(3), intent(in) :: xv,anv
real(wp_), intent(in) :: bres
integer, intent(in) :: sox
2019-03-28 10:50:28 +01:00
integer, dimension(:), intent(inout), pointer :: iop ! in/out plasma flag
complex(wp_), dimension(:), intent(out), pointer :: ext,eyt
2019-03-28 10:50:28 +01:00
! local variables
real(wp_) :: rm,csphi,snphi,bphi,br,bz
real(wp_), dimension(3) :: bv,xmv
!
2019-03-28 10:50:28 +01:00
iop(i)=iop(i)+1 ! in->out
xmv=xv*0.01_wp_ ! convert from cm to m
rm=sqrt(xmv(1)**2+xmv(2)**2)
csphi=xmv(1)/rm
snphi=xmv(2)/rm
call bfield(rm,xmv(3),br,bz,bphi)
bv(1)=br*csphi-bphi*snphi
bv(2)=br*snphi+bphi*csphi
bv(3)=bz
call pol_limit(anv,bv,bres,sox,ext(i),eyt(i)) ! polarization at plasma exit
end subroutine plasma_out
! ------------------------------
! subroutine wall_in(i) ! ray enters vessel
! integer, intent(in) :: i ! ray index
!
! iow(i)=iow(i)+1
! end subroutine wall_in
! ------------------------------
subroutine wall_out(i,ins,xv,anv,bres,sox,psipol1,chipol1,iow,iop,ext,eyt) ! ray exits vessel
2019-03-28 10:50:28 +01:00
! arguments
integer, intent(in) :: i ! ray index
logical, intent(in) :: ins ! inside plasma? (ins=1 plasma/wall overlap)
real(wp_), dimension(3), intent(inout) :: xv,anv
real(wp_), intent(in) :: bres
integer, intent(in) :: sox
real(wp_), intent(out) :: psipol1,chipol1
2019-03-28 10:50:28 +01:00
integer, dimension(:), intent(inout), pointer :: iow,iop ! in/out vessel and plasma flags
complex(wp_), dimension(:), intent(inout), pointer :: ext,eyt
2019-03-28 10:50:28 +01:00
! local variables
integer :: irfl
real(wp_) :: qq1,uu1,vv1
real(wp_), dimension(3) :: xvrfl,anvrfl,walln
complex(wp_) :: ext1,eyt1
2019-03-28 10:50:28 +01:00
!
iow(i)=iow(i)+1 ! out->in
if(ins) call plasma_out(i,xv,anv,bres,sox,iop,ext,eyt) ! plasma-wall overlapping
call wall_refl(xv-dst*anv,anv,ext(i),eyt(i),xvrfl,anvrfl,ext1,eyt1,walln,irfl) ! ray reflects at wall
ext(i) = ext1 ! save parameters at wall reflection
eyt(i) = eyt1
xv = xvrfl
anv = anvrfl
if(i.eq.1) then ! polarization angles at wall reflection for central ray
call stokes_ce(ext1,eyt1,qq1,uu1,vv1)
call polellipse(qq1,uu1,vv1,psipol1,chipol1)
psipol1=psipol1/degree ! convert from rad to degree
chipol1=chipol1/degree
else
psipol1 = zero
chipol1 = zero
end if
end subroutine wall_out
! ------------------------------
subroutine initbeam(i,iroff,iboff,iwait,stv,jphi_beam,pins_beam,currins_beam, &
dpdv_beam,jcd_beam) ! initialization at beam propagation start
2021-12-18 18:57:38 +01:00
use logger, only : log_info, log_warning
2019-03-28 10:50:28 +01:00
! arguments
integer, intent(in) :: i ! beam index
2019-03-28 10:50:28 +01:00
logical, dimension(:,:), intent(in), pointer :: iroff ! global ray status (F = active, T = inactive)
logical, intent(out) :: iboff
logical, dimension(:), intent(out), pointer :: iwait
real(wp_), dimension(:), intent(out), pointer :: jphi_beam,pins_beam, &
currins_beam,dpdv_beam,jcd_beam,stv
2021-12-18 18:57:38 +01:00
character(256) :: msg ! buffer for formatting log messages
iboff = .false. ! beam status (F = active, T = inactive)
iwait = iroff(:,i) ! copy ray status for current beam from global ray status
if(all(iwait)) then ! no rays active => stop beam
iboff = .true.
else if (any(iwait)) then
! only some rays active
write (msg,'(" beam ",g0,": only some rays are active!")') i
2021-12-18 18:57:38 +01:00
call log_warning(msg, mod='multipass', proc='initbeam')
end if
stv = zero ! starting step
jphi_beam = zero ! 1D beam profiles
pins_beam = zero
currins_beam = zero
dpdv_beam = zero
jcd_beam = zero
end subroutine initbeam
! ------------------------------
subroutine initmultipass(i,iox,iroff,yynext,yypnext,yw0,ypw0,stnext,p0ray, &
taus,tau1,etau1,cpls,cpl1,lgcpl1,psipv,chipv) ! initialization before pass loop
2019-03-28 10:50:28 +01:00
! arguments
logical, intent(in) :: i ! ipol
integer, intent(in) :: iox ! mode active on 1st pass
2019-03-28 10:50:28 +01:00
logical, dimension(:,:), intent(out), pointer :: iroff ! global ray status (F = active, T = inactive)
real(wp_), dimension(:), intent(out), pointer :: p0ray,tau1,etau1,cpl1,lgcpl1, &
psipv,chipv
real(wp_), dimension(:,:), intent(out), pointer :: yw0,ypw0,stnext,taus,cpls
real(wp_), dimension(:,:,:), intent(out), pointer :: yynext,yypnext
!
iroff = .false. ! global ray status (F = active, T = inactive)
if(.not. i) call turnoffray(0,1,3-iox,iroff) ! !ipol => stop other mode (iox=1/2 -> stop ib=2/1 at first pass)
yynext = zero ! starting beam coordinates (1)
yypnext = zero ! starting beam coordinates (2)
yw0 = zero ! temporary beam coordinates (1)
ypw0 = zero ! temporary beam coordinates (2)
stnext = zero ! starting beam step
p0ray = zero ! starting beam power
2019-03-28 10:50:28 +01:00
taus = zero ! beam tau from previous passes
tau1 = zero
etau1 = one
2019-03-28 10:50:28 +01:00
cpls = one ! beam coupling from previous passes
cpl1 = one
lgcpl1 = zero
psipv = zero ! psi polarization angle at vacuum-plasma boundary
chipv = zero ! chi polarization angle at vacuum-plasma boundary
end subroutine initmultipass
! ------------------------------
subroutine turnoffray(jk,ip,ib,iroff) ! turn off ray propagation
2019-03-28 10:50:28 +01:00
! arguments
integer, intent(in) :: jk, ip, ib ! ray (0=all rays), pass, beam indexes
2019-03-28 10:50:28 +01:00
logical, dimension(:,:), intent(out), pointer :: iroff ! global ray status (F = active, T = inactive)
! local variables
integer :: ipx, i1, i2
!
if(jk==0) then ! stop all rays
2019-03-28 10:50:28 +01:00
do ipx=ip,ipass ! from pass ip to last pass
i1 = 2**ipx-2+2**(ipx-ip)*(ib-1)+1 ! first derived beam at pass ipx
i2 = 2**ipx-2+2**(ipx-ip)*ib ! last derived beam at pass ipx (i1=i2 for ipx=ip)
iroff(:,i1:i2) = .true.
end do
else ! only stop ray jk
do ipx=ip,ipass
i1 = 2**ipx-2+2**(ipx-ip)*(ib-1)+1
i2 = 2**ipx-2+2**(ipx-ip)*ib
iroff(jk,i1:i2) = .true.
end do
end if
end subroutine turnoffray
! ------------------------------
subroutine alloc_multipass(dim,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)
integer :: dim
logical, dimension(:), intent(out), pointer :: iwait
logical, dimension(:,:), intent(out), pointer :: iroff
integer, dimension(:), intent(out), pointer :: iop,iow
real(wp_), dimension(:), intent(out), pointer :: jphi_beam,pins_beam,currins_beam, &
dpdv_beam,jcd_beam,stv,tau1,etau1,cpl1,lgcpl1,p0ray,psipv,chipv
real(wp_), dimension(:,:), intent(out), pointer :: taus,cpls,stnext,yw0,ypw0
real(wp_), dimension(:,:,:), intent(out), pointer :: yynext,yypnext
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)
nbeam_max = 2**ipass ! max n of beams active at a time
nbeam_tot = 2**(ipass+1)-2 ! total n of beams
allocate(iwait(nray),iroff(nray,nbeam_tot),iop(nray),iow(nray), &
yynext(6,nray,nbeam_max-2),yypnext(6,nray,nbeam_max-2), &
yw0(6,nray),ypw0(6,nray),stnext(nray,nbeam_tot),stv(nray), &
p0ray(nray),taus(nray,nbeam_tot),tau1(nray),etau1(nray), &
cpls(nray,nbeam_tot),cpl1(nray),lgcpl1(nray),jphi_beam(dim), &
pins_beam(dim),currins_beam(dim),dpdv_beam(dim),jcd_beam(dim), &
fix coupling for subsequent beams In situations when multiple beams are traced, either when allowing multiple plasma crossings (raytracing.ipass > 0) or the initial polarisation is mixed (raytracing.ipol == .true.), the couplings of all but the first beam (with least index_rt) were invalid. The bug is due to the re-use of the psipol,chipol variables as the beams are traced sequentially over the beam_loop. For the first beam being traced the psipol,chipol are correctly initialised to the user-defined value and the resulting coupling is correct. However, in each subsequent beam the values were not set to those of the parent beam (or to the user-defined value in the case of the first X mode beam), but to those of the previous beams (current index_rt - 1). This change repurposes the psipv,chipv arrays to store the polarisation of the parent beams, including the initial user-defined value and makes plasma_in always use these to compute the coupling. In addition, in the case the polarisation is not immediately known (i.e. if raytracing.ipol == .false.), this change postpones the computation of the Jones vector (ext, eyt) from the launch point, if the magnetic equilibrium is available, to when the ray actually crosses the plasma boundary. The original code, besides being strictly incorrect, can lead to non-negligible alterations to the coupling. This change also mean: 1. most of the functionality of `set_pol` has been merged with `plasma_in` 2. the polarisation is undefined and the Jones vector is set to the placeholder value [1, 0] till `plasma_im` is called Finally, `polarcold` is removed because it's unused.
2024-03-02 17:32:03 +01:00
psipv(0:nbeam_tot),chipv(0:nbeam_tot))
end subroutine alloc_multipass
! ------------------------------
subroutine 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)
logical, dimension(:), intent(out), pointer :: iwait
logical, dimension(:,:), intent(out), pointer :: iroff
integer, dimension(:), intent(out), pointer :: iop,iow
real(wp_), dimension(:), intent(out), pointer :: stv,p0ray,tau1,etau1,cpl1,lgcpl1, &
jphi_beam,pins_beam,currins_beam,dpdv_beam,jcd_beam,psipv,chipv
real(wp_), dimension(:,:), intent(out), pointer :: yw0,ypw0,stnext,taus,cpls
real(wp_), dimension(:,:,:), intent(out), pointer :: yynext,yypnext
if (associated(iwait)) deallocate(iwait)
if (associated(iroff)) deallocate(iroff)
if (associated(iop)) deallocate(iop)
if (associated(iow)) deallocate(iow)
if (associated(yynext)) deallocate(yynext)
if (associated(yypnext)) deallocate(yypnext)
if (associated(yw0)) deallocate(yw0)
if (associated(ypw0)) deallocate(ypw0)
if (associated(stnext)) deallocate(stnext)
if (associated(stv)) deallocate(stv)
if (associated(p0ray)) deallocate(p0ray)
if (associated(taus)) deallocate(taus)
if (associated(tau1)) deallocate(tau1)
if (associated(etau1)) deallocate(etau1)
if (associated(cpls)) deallocate(cpls)
if (associated(cpl1)) deallocate(cpl1)
if (associated(lgcpl1)) deallocate(lgcpl1)
if (associated(jphi_beam)) deallocate(jphi_beam)
if (associated(pins_beam)) deallocate(pins_beam)
if (associated(currins_beam)) deallocate(currins_beam)
if (associated(dpdv_beam)) deallocate(dpdv_beam)
if (associated(jcd_beam)) deallocate(jcd_beam)
if (associated(psipv)) deallocate(psipv)
if (associated(chipv)) deallocate(chipv)
end subroutine dealloc_multipass
end module multipass