cart2polar Function

public function cart2polar(coord, n) result(polar)

Convert from cartesian coordinates (in the x-y plane) to polar coordinates (radius,angle)

Arguments

Type IntentOptional Attributes Name
real(kind=rk), intent(in) :: coord(n,3)
integer, intent(in) :: n

Return Value real(kind=rk), (n,2)

Polar coordinates, radius (first entry) and angle (second entry)


Called by

proc~~cart2polar~~CalledByGraph proc~cart2polar cart2polar proc~tem_eval_cylindricalwave tem_eval_cylindricalWave proc~tem_eval_cylindricalwave->proc~cart2polar proc~tem_spatial_for_coord tem_spatial_for_coord proc~tem_spatial_for_coord->proc~tem_eval_cylindricalwave proc~tem_spacetime_for_coord tem_spacetime_for_coord proc~tem_spacetime_for_coord->proc~tem_eval_cylindricalwave proc~tem_spatial_scalar_for_index tem_spatial_scalar_for_index proc~tem_spatial_scalar_for_index->proc~tem_spatial_for_coord proc~tem_spacetime_for_stcoord tem_spacetime_for_stcoord proc~tem_spacetime_for_stcoord->proc~tem_spacetime_for_coord interface~tem_spatial_for tem_spatial_for interface~tem_spatial_for->proc~tem_spatial_for_coord interface~tem_spacetime_for tem_spacetime_for interface~tem_spacetime_for->proc~tem_spacetime_for_coord proc~tem_spacetime_scalar_for_index tem_spacetime_scalar_for_index proc~tem_spacetime_scalar_for_index->proc~tem_spacetime_for_coord

Contents

Source Code


Source Code

  function cart2polar(coord, n) result(polar)
    ! --------------------------------------------------------------------------
    integer, intent(in) :: n
    real(kind=rk),intent(in) :: coord(n,3)
    !> Polar coordinates, radius (first entry) and angle (second entry)
    real(kind=rk) :: polar(n,2)
    ! --------------------------------------------------------------------------
    !> X coordinate
    real(kind=rk) :: x
    !> Y coordinate
    real(kind=rk) :: y
    integer :: i
    ! --------------------------------------------------------------------------

    do i = 1,n

      x = coord(i,1)
      y = coord(i,2)

      polar(i,1) = sqrt(x*x + y*y)

      ! Atan2 is not defined when both coordinates are zero. To cover
      ! this situation correctly, we define the angle to be 0.
      if(polar(i,1) > epsilon(1.0_rk) ) then
        polar(i,2) = atan2(y,x)
      else
        polar(i,2) = 0.0_rk
      end if

    end do

  end function cart2polar