add_sphere Subroutine

private subroutine add_sphere(me, center, radius, num_lat_points, num_lon_points)

Add a sphere to an STL file.

Type Bound

stl_file

Arguments

Type IntentOptional Attributes Name
class(stl_file), intent(inout) :: me
real(kind=wp), intent(in), dimension(3) :: center

coordinates of sphere center [x,y,z]

real(kind=wp), intent(in) :: radius

radius of the sphere

integer, intent(in) :: num_lat_points

number of latitude points (not counting poles)

integer, intent(in) :: num_lon_points

number of longitude points


Calls

proc~~add_sphere~~CallsGraph proc~add_sphere stl_file%add_sphere proc~add_plate stl_file%add_plate proc~add_sphere->proc~add_plate proc~spherical_to_cartesian spherical_to_cartesian proc~add_sphere->proc~spherical_to_cartesian

Source Code

    subroutine add_sphere(me,center,radius,num_lat_points,num_lon_points)

    implicit none

    class(stl_file),intent(inout)    :: me
    real(wp),dimension(3),intent(in) :: center         !! coordinates of sphere center [x,y,z]
    real(wp),intent(in)              :: radius         !! radius of the sphere
    integer,intent(in)               :: num_lat_points !! number of latitude points (not counting poles)
    integer,intent(in)               :: num_lon_points !! number of longitude points

    integer :: i  !! counter
    integer :: j  !! counter
    real(wp) :: delta_lat  !! step in latitude (deg)
    real(wp) :: delta_lon  !! step in longitude (deg)
    real(wp),dimension(:),allocatable :: lat !! array of latitude values (deg)
    real(wp),dimension(:),allocatable :: lon !! array of longitude value (deg)
    real(wp),dimension(3) :: v1,v2,v3,v4 !! vertices

    ! Example:
    !
    ! num_lat_points = 3
    ! num_lon_points = 5
    !
    !  90 -------------  North pole
    !     |  *  *  *  |
    !     |  *  *  *  |
    !     |  *  *  *  |
    ! -90 -------------  South pole
    !     0          360

    delta_lat = 180.0_wp / (1+num_lat_points)
    delta_lon = 360.0_wp / (1+num_lon_points)

    lat = -90.0_wp + [(delta_lat*(i-1), i = 1,num_lat_points+2)]
    lon =            [(delta_lon*(i-1), i = 1,num_lon_points+2)]

    ! generate all the plates on the sphere.
    ! start at bottom left and go right then up.
    ! each box is two triangular plates.
    do i = 1, num_lat_points+1
        do j = 1, num_lon_points+1

            !   3----2
            !   |  / |
            !   | /  |
            ! i 1----4
            !   j

            v1 = spherical_to_cartesian(radius,lon(j),  lat(i)  ) + center
            v2 = spherical_to_cartesian(radius,lon(j+1),lat(i+1)) + center
            v3 = spherical_to_cartesian(radius,lon(j),  lat(i+1)) + center
            v4 = spherical_to_cartesian(radius,lon(j+1),lat(i)  ) + center
            call me%add_plate(v1,v2,v3)
            call me%add_plate(v1,v4,v2)

        end do
    end do

    end subroutine add_sphere