Evaluate magnetizing field (angular-component) for Mie-Scattering of electromagnetic wave at dielectric cylinder.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_miescatter_field_type), | intent(in) | :: | me |
The function to evaluate |
||
real(kind=rk), | intent(in) | :: | coord(n,3) |
barycentric Ids of an elements. 1st index goes over number of elements and 2nd index goes over x,y,z coordinates |
||
real(kind=rk), | intent(in) | :: | time |
The time to evaluate the function at. |
||
integer, | intent(in) | :: | n |
Number of points to evaluate the function for. |
return value of the function
function tem_eval_miescatter_magnangular(me, coord, time, n) result(res) ! --------------------------------------------------------------------------- !> The function to evaluate type(tem_miescatter_field_type), intent(in) :: me !> Number of points to evaluate the function for. integer, intent(in) :: n !> barycentric Ids of an elements. !! 1st index goes over number of elements and !! 2nd index goes over x,y,z coordinates real(kind=rk), intent( in ) :: coord(n,3) !> The time to evaluate the function at. real(kind=rk), intent( in ) :: time !> return value of the function real(kind=rk) :: res(n) ! --------------------------------------------------------------------------- integer :: iPoint, iCoeff ! Polar coordinate vector in x-y plane. ! First entry: radius \n ! Second entry: angle real(kind=rk) :: polar(2) complex(kind=rk) :: tmp ! --------------------------------------------------------------------------- do iPoint = 1, n ! Convert to polar coordinates (relative to the center of the ! cylinder. polar = cart2polar( coord(iPoint,1)-me%miescatter%center(1), & & coord(iPoint,2)-me%miescatter%center(2) ) ! Inside the cylinder if(polar(1).le.me%miescatter%radius) then tmp = (0.0_rk, 0.0_rk) do iCoeff = -me%mieexpansion%nCoeffs , me%mieexpansion%nCoeffs tmp = tmp & & + me%mieexpansion%c_tot(iCoeff) & & * bessel_jn_derivative(iCoeff, & & me%miescatter%wavenumber_cylinder*polar(1)) & & * exp( (0.0_rk,1.0_rk) * iCoeff * polar(2) ) end do ! Twiddle by phase and get the real part res(iPoint) = real( tmp & ! JZ: the paper has an additional factor of -1 here. However, I think ! that this factor is wrong at this place. So, I commented out ! the factor. ! & * (-1.0_rk) & & * (exp( (0.0_rk,1.0_rk)*me%miescatter%omega*time)) & & * (0.0_rk,-1.0_rk)*me%miescatter%wavenumber_cylinder & & / (me%miescatter%omega*me%miescatter%permeability_cylinder) ) ! Outside the cylinder else ! Add up incoming and scattered field tmp = (0.0_rk, 0.0_rk) do iCoeff = -me%mieexpansion%nCoeffs , me%mieexpansion%nCoeffs tmp = tmp + & & ( (0.0_rk,1.0_rk)**(-iCoeff) & & * bessel_jn_derivative(iCoeff, & & me%miescatter%wavenumber_background*polar(1)) & & + me%mieexpansion%c_scat(iCoeff) & & * hankel2_n_derivative(iCoeff, & & me%miescatter%wavenumber_background*polar(1)) & & ) * exp( (0.0_rk,1.0_rk) * iCoeff * polar(2) ) end do ! Twiddle by phase and get the real part res(iPoint) = real( tmp & ! JZ: the paper has an additional factor of -1 here. However, I think ! that this factor is wrong at this place. So, I commented out ! the factor. ! & * (-1.0_rk) & & * (exp( (0.0_rk,1.0_rk) * me%miescatter%omega * time )) & & * (0.0_rk,-1.0_rk)*me%miescatter%wavenumber_background & & / (me%miescatter%omega*me%miescatter%permeability_background) ) end if end do end function tem_eval_miescatter_magnangular