2015-11-18 17:34:33 +01:00
|
|
|
module polarization
|
2024-01-27 12:09:56 +01:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2015-11-18 17:34:33 +01:00
|
|
|
interface stokes
|
|
|
|
module procedure stokes_ce,stokes_ell
|
|
|
|
end interface
|
|
|
|
|
|
|
|
contains
|
|
|
|
subroutine stokes_ce(ext,eyt,qq,uu,vv)
|
2023-09-20 16:03:08 +02:00
|
|
|
use const_and_precisions, only : wp_
|
2015-11-18 17:34:33 +01:00
|
|
|
! arguments
|
|
|
|
complex(wp_), intent(in) :: ext,eyt
|
|
|
|
real(wp_), intent(out) :: qq,uu,vv
|
|
|
|
|
|
|
|
qq = abs(ext)**2 - abs(eyt)**2
|
2023-09-20 16:03:08 +02:00
|
|
|
uu = 2*real(ext*conjg(eyt))
|
|
|
|
vv = 2*imag(ext*conjg(eyt))
|
2015-11-18 17:34:33 +01:00
|
|
|
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)
|
2021-12-15 02:31:12 +01:00
|
|
|
use const_and_precisions, only : wp_,ui=>im,zero,one
|
2015-11-18 17:34:33 +01:00
|
|
|
! arguments
|
|
|
|
real(wp_), dimension(3), intent(in) :: anv,bv
|
2022-05-20 13:10:49 +02:00
|
|
|
real(wp_), intent(in) :: bres
|
|
|
|
integer, intent(in) :: sox
|
2015-11-18 17:34:33 +01:00
|
|
|
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)
|
2019-03-26 15:21:22 +01:00
|
|
|
den = sqrt(abs(ext)**2+abs(eyt)**2)
|
|
|
|
ext = ext/den
|
|
|
|
eyt = eyt/den
|
2015-11-18 17:34:33 +01:00
|
|
|
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
|