f82f91bc8d
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.
101 lines
2.5 KiB
Fortran
101 lines
2.5 KiB
Fortran
module polarization
|
|
|
|
implicit none
|
|
|
|
interface stokes
|
|
module procedure stokes_ce,stokes_ell
|
|
end interface
|
|
|
|
contains
|
|
subroutine stokes_ce(ext,eyt,qq,uu,vv)
|
|
use const_and_precisions, only : wp_
|
|
! arguments
|
|
complex(wp_), intent(in) :: ext,eyt
|
|
real(wp_), intent(out) :: qq,uu,vv
|
|
|
|
qq = abs(ext)**2 - abs(eyt)**2
|
|
uu = 2*real(ext*conjg(eyt))
|
|
vv = 2*imag(ext*conjg(eyt))
|
|
end subroutine stokes_ce
|
|
|
|
|
|
subroutine stokes_ell(chi,psi,qq,uu,vv)
|
|
use const_and_precisions, only : wp_,two
|
|
! arguments
|
|
real(wp_), intent(in) :: chi,psi
|
|
real(wp_), intent(out) :: qq,uu,vv
|
|
|
|
qq=cos(two*chi)*cos(two*psi)
|
|
uu=cos(two*chi)*sin(two*psi)
|
|
vv=sin(two*chi)
|
|
end subroutine stokes_ell
|
|
|
|
|
|
subroutine polellipse(qq,uu,vv,psi,chi)
|
|
use const_and_precisions, only : wp_,half
|
|
! arguments
|
|
real(wp_), intent(in) :: qq,uu,vv
|
|
real(wp_), intent(out) :: psi,chi
|
|
! real(wp_) :: ll,aa,bb,ell
|
|
|
|
! ll = sqrt(qq**2 + uu**2)
|
|
! aa = sqrt(half*(1 + ll))
|
|
! bb = sqrt(half*(1 - ll))
|
|
! ell = bb/aa
|
|
psi = half*atan2(uu,qq)
|
|
chi = half*asin(vv)
|
|
end subroutine polellipse
|
|
|
|
subroutine pol_limit(anv,bv,bres,sox,ext,eyt) !,gam)
|
|
use const_and_precisions, only : wp_,ui=>im,zero,one
|
|
! arguments
|
|
real(wp_), dimension(3), intent(in) :: anv,bv
|
|
real(wp_), intent(in) :: bres
|
|
integer, intent(in) :: sox
|
|
complex(wp_), intent(out) :: ext,eyt
|
|
! real(wp_), optional, intent(out) :: gam
|
|
! local variables
|
|
real(wp_), dimension(3) :: bnv
|
|
real(wp_) :: anx,any,anz,an2,an,anpl2,anpl,anpr,anxy, &
|
|
btot,yg,den,dnl,del0,ff,ff2,sngam,csgam
|
|
!
|
|
btot = sqrt(bv(1)**2+bv(2)**2+bv(3)**2)
|
|
bnv = bv/btot
|
|
yg = btot/bres
|
|
|
|
anx = anv(1)
|
|
any = anv(2)
|
|
anz = anv(3)
|
|
an2 = anx**2 + any**2 + anz**2
|
|
an = sqrt(an2)
|
|
anxy = sqrt(anx**2 + any**2)
|
|
|
|
anpl = (anv(1)*bnv(1) + anv(2)*bnv(2) + anv(3)*bnv(3))
|
|
anpl2= anpl**2
|
|
anpr = sqrt(an2 - anpl2)
|
|
|
|
dnl = one - anpl2
|
|
del0 = sqrt(dnl**2 + 4.0_wp_*anpl2/yg**2)
|
|
|
|
sngam = (anz*anpl - an2*bnv(3))/(an*anxy*anpr)
|
|
csgam = -(any*bnv(1) - anx*bnv(2))/ (anxy*anpr)
|
|
|
|
ff = 0.5_wp_*yg*(dnl - sox*del0)
|
|
ff2 = ff**2
|
|
den = ff2 + anpl2
|
|
if (den>zero) then
|
|
ext = (ff*csgam - ui*anpl*sngam)/sqrt(den)
|
|
eyt = (-ff*sngam - ui*anpl*csgam)/sqrt(den)
|
|
den = sqrt(abs(ext)**2+abs(eyt)**2)
|
|
ext = ext/den
|
|
eyt = eyt/den
|
|
else ! only for XM (sox=+1) when N//=0
|
|
ext = -ui*sngam
|
|
eyt = -ui*csgam
|
|
end if
|
|
|
|
! gam = atan2(sngam,csgam)/degree
|
|
end subroutine pol_limit
|
|
|
|
end module polarization
|