blanket_mod_pol_height Subroutine

public subroutine blanket_mod_pol_height()

Calculations for blanket module poloidal height author: J. Morris, CCFE, Culham Science Centre Calculations for blanket module poloidal height for D shaped and elliptical machines

Arguments

None

Contents


Source Code

    subroutine blanket_mod_pol_height

        !! Calculations for blanket module poloidal height
        !! author: J. Morris, CCFE, Culham Science Centre
        !! Calculations for blanket module poloidal height for D shaped and elliptical machines
        !
        ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

        use constants, only: pi
        use build_variables, only: scrapli, scraplo
        use fwbs_variables, only: fwbsshape, nblktmodpi, fdiv, nblktmodpo, &
        iblanket
        use physics_variables, only: itart, rminor, idivrt, rmajor, triang

        implicit none

        ! Local variables

        ! Mid-plane distance from inboard to outboard side (m)
        real(dp) :: a

        ! Internal half-height of blanket (m)
        real(dp) :: b

        ! Calculate ellipse circumference using Ramanujan approximation (m)
        real(dp) :: ptor

        ! Major radius where half-ellipses 'meet' (m)
        real(dp) :: r1

        ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

        if ((itart == 1).or.(fwbsshape == 1)) then  ! D-shaped machine

            ! Segment vertical inboard surface (m)
            bllengi = (2.0D0*hblnkt) / nblktmodpi

            ! Calculate perimeter of ellipse that defines the internal
            ! surface of the outboard first wall / blanket

            ! Mid-plane distance from inboard to outboard side (m)
            a = scrapli + 2.0D0*rminor + scraplo

            ! Internal half-height of blanket (m)
            b = hblnkt

            ! Calculate ellipse circumference using Ramanujan approximation (m)
            ptor = pi * ( 3.0D0*(a+b) - sqrt( (3.0D0*a + b)*(a + 3.0D0*b) ) )

            ! Calculate blanket poloidal length and segment, subtracting divertor length (m)
            ! kit hcll version only had the single null option
            if (idivrt == 2) then
                ! Double null configuration
                bllengo = 0.5D0*ptor * (1.0D0 - 2.0D0*fdiv) / nblktmodpo
            else
                ! single null configuration
                bllengo = 0.5D0*ptor * (1.0D0 - fdiv) / nblktmodpo
            end if

        ! shape defined by two half-ellipses
        else

            ! Major radius where half-ellipses 'meet' (m)
            r1 = rmajor - rminor*triang

            ! Internal half-height of blanket (m)
            b = hblnkt

            ! Distance between r1 and nearest edge of inboard first wall / blanket (m)
            a = r1 - (rmajor - rminor - scrapli)

            ! Calculate ellipse circumference using Ramanujan approximation (m)
            ptor = pi * ( 3.0D0*(a+b) - sqrt( (3.0D0*a + b)*(a + 3.0D0*b) ) )

            ! Calculate inboard blanket poloidal length and segment, subtracting divertor length (m)
            ! Assume divertor lies between the two ellipses, so fraction fdiv still applies

            ! kit hcll version only had the single null option
            if (idivrt == 2) then
                ! Double null configuration
                bllengi = 0.5D0*ptor * (1.0D0 - 2.0D0*fdiv) / nblktmodpi
            else
                ! single null configuration
                bllengi = 0.5D0*ptor * (1.0D0 - fdiv) / nblktmodpi
            end if

            ! Distance between r1 and inner edge of outboard first wall / blanket (m)
            a = rmajor + rminor + scraplo - r1

            ! Calculate ellipse circumference using Ramanujan approximation (m)
            ptor = pi * ( 3.0D0*(a+b) - sqrt( (3.0D0*a + b)*(a + 3.0D0*b) ) )

            ! kit hcll version only had the single null option
            ! Calculate outboard blanket poloidal length and segment, subtracting divertor length (m)
            if (idivrt == 2) then
                ! Double null configuration
                bllengo = 0.5D0*ptor * (1.0D0 - 2.0D0*fdiv) / nblktmodpo
            else
                ! single null configuration
                bllengo = 0.5D0*ptor * (1.0D0 - fdiv) / nblktmodpo
            end if

        end if

    end subroutine blanket_mod_pol_height