liquid_breeder_pressure_drop_mhd Function

public function liquid_breeder_pressure_drop_mhd(ip, ofile, vel, vsc, conduct_liq, l_channel, num_pol, label)

Calculates the pressure drop in a liquid metal flow channel due to MHD effects. The total pressure drop is the sum of contributions. This is only used for secondary coolant/breeder so rectangular flow channels are assumed.

author: G Graham, CCFE

References:

  [Miy1986]   Miyazaki et al. (1986), Magneto-Hydro-Dynamic Pressure Drop of Lithium
              Flow in Rectangular Ducts, Fusion Technology, 10:3P2A, 830-836, DOI: 10.13182/FST10-830

  [Mal1995]   Malang and Mattas (1995), Comparison of lithium and the eutectic
              lead-lithium alloy, two candidate liquid metal breeder materials
              for self-cooled blankets, Fusion Engineering and Design 27, 399-406

  [Iba2013]   Ibano et al (2013), Nutronics and pumping power analysis on the
              Tokamak reactor for the fusion-biomass hybrid concept,
              Fusion Engineering and Design, 88

  [Sho2018]   Shoki et al (2018), MHD pressure drop measurement of PbLi flow
              in double-bended pipe, Fusion Engineering and Design, 136, 17-23

  [Klu2019]   Kluber et al. (2019), Numerical simulations of 3D magnetohydrodynamic
              flows in dual-coolant lead lithium blankets, Fusion Engineering and Design,
              146, 684-687

  [Sua2021]   MHD effects in geometrical sigularities on high velocity breeding
              blanket designs. Part II, ENR-PRD.BB-T007-D002, EFDA_D_2PDT9U.
              Also, see asssociated paper: Suarez et al. (2021), On the use of CFD
              to obtain head loss coefficients in hydraulic systems and it's appliaction
              to liquid metal flows in nuclear fusion reactor blankets, Plasma. Phys.
              Control fusion, 63, 124002

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

Function return parameter !!!!!

IB OB

Half-widths N.B. a_bz_liq (width in the toroidal direction) is in B direction If have thin conducting walls... Caculate resistances of fluid and walls Calculate pressure drop for conducting wall [Miy1986] If have perfcetly insulating FCIs... Calculate pressure drop for (perfectly) insulating FCI [Mal1995] Total (Pa)

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: ip

Liquid metal coolant/breeder flow velocity (m/s)

integer, intent(in) :: ofile

Liquid metal coolant/breeder flow velocity (m/s)

real(kind=dp), intent(in) :: vel

Liquid metal visosity

real(kind=dp), intent(in) :: vsc

Liquid metal conductivity

real(kind=dp), intent(in) :: conduct_liq

Length long poloidal sections of channel

real(kind=dp), intent(in) :: l_channel

Number long poloidal sections of channel

integer, intent(in) :: num_pol

Description of this calculation

character(len=*), intent(in) :: label

Local Variables !!!!!!!!!!!!!!!

Half-widths of channel (m)

Return Value real(kind=dp)

Arguments !!!!!!!!!!!!!!!!!!!!!


Contents


Source Code

    function liquid_breeder_pressure_drop_mhd(ip, ofile, vel, vsc, conduct_liq, l_channel, num_pol, label)

        !! Calculates the pressure drop in a liquid metal flow channel due to MHD effects. The total pressure
        !! drop is the sum of contributions. This is only used for secondary coolant/breeder so rectangular flow
        !! channels are assumed.
        !!
        !! author: G Graham, CCFE
        !!
        !! References:
        !!
        !!      [Miy1986]   Miyazaki et al. (1986), Magneto-Hydro-Dynamic Pressure Drop of Lithium
        !!                  Flow in Rectangular Ducts, Fusion Technology, 10:3P2A, 830-836, DOI: 10.13182/FST10-830
        !!
        !!      [Mal1995]   Malang and Mattas (1995), Comparison of lithium and the eutectic
        !!                  lead-lithium alloy, two candidate liquid metal breeder materials
        !!                  for self-cooled blankets, Fusion Engineering and Design 27, 399-406
        !!
        !!      [Iba2013]   Ibano et al (2013), Nutronics and pumping power analysis on the
        !!                  Tokamak reactor for the fusion-biomass hybrid concept,
        !!                  Fusion Engineering and Design, 88
        !!
        !!      [Sho2018]   Shoki et al (2018), MHD pressure drop measurement of PbLi flow
        !!                  in double-bended pipe, Fusion Engineering and Design, 136, 17-23
        !!
        !!      [Klu2019]   Kluber et al. (2019), Numerical simulations of 3D magnetohydrodynamic
        !!                  flows in dual-coolant lead lithium blankets, Fusion Engineering and Design,
        !!                  146, 684-687
        !!
        !!      [Sua2021]   MHD effects in geometrical sigularities on high velocity breeding
        !!                  blanket designs. Part II, ENR-PRD.BB-T007-D002, EFDA_D_2PDT9U.
        !!                  Also, see asssociated paper: Suarez et al. (2021), On the use of CFD
        !!                  to obtain head loss coefficients in hydraulic systems and it's appliaction
        !!                  to liquid metal flows in nuclear fusion reactor blankets, Plasma. Phys.
        !!                  Control fusion, 63, 124002
        !!
        !! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

        use fwbs_variables, only: ifci, a_bz_liq, b_bz_liq, hartmann_liq, b_mag_blkt, &
        bz_channel_conduct_liq, th_wall_secondary
        use constraint_variables, only: maxradwallload
        use physics_variables, only: btot

        implicit none

        !! Function return parameter !!!!!

        real(dp) :: liquid_breeder_pressure_drop_mhd

        !! Arguments !!!!!!!!!!!!!!!!!!!!!

        integer, intent(in) :: ip, ofile

        !! Liquid metal coolant/breeder  flow velocity (m/s)
        real(dp), intent(in) :: vel

        !! Liquid metal visosity
        real(dp), intent(in) :: vsc

        !! Liquid metal conductivity
        real(dp), intent(in) :: conduct_liq

        !! Length long poloidal sections of channel
        real(dp), intent(in) :: l_channel

        !! Number long poloidal sections of channel
        integer, intent(in) :: num_pol

        !! Description of this calculation
        character(len=*), intent(in) :: label

        !! Local Variables !!!!!!!!!!!!!!!

        !! Half-widths of channel (m)
        real(dp)  :: half_wth_a, half_wth_b

        !! MHD pressure drop for single channel
        real(dp)  :: mhd_pressure_drop

        !! Magnetic field strenght (T)
        real(dp) :: b_mag

        !! Internal resistance of fluid
        real(dp) :: r_i

        !! Resistance of side walls
        real(dp) :: r_w

        !! Resistance of electrode walls
        real(dp) :: r_e

        !! Wall-to-fluid conductance ratio
        real(dp) :: big_c

        !!  Term calculated for [Miy1986] pressure drop
        real(dp) :: kp

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

        !! Magnetic feild strength in IB or OB blanket
        if (label=='Inboard blanket breeder liquid') b_mag = b_mag_blkt(1) !! IB
        if (label=='Outboard blanket breeder liquid') b_mag = b_mag_blkt(2) !! OB

        !! Half-widths
        !! N.B. a_bz_liq (width in the toroidal direction) is in B direction
        half_wth_a = a_bz_liq * 0.5
        half_wth_b = b_bz_liq *0.5

        !! If have thin conducting walls...
        if (ifci/=1) then

            !! Caculate resistances of fluid and walls
            r_i = half_wth_b/(conduct_liq * half_wth_a)
            r_w = half_wth_b/(bz_channel_conduct_liq * th_wall_secondary)
            r_e = half_wth_a/(bz_channel_conduct_liq * th_wall_secondary)
            big_c = r_i/r_w
            !!  Calculate pressure drop for conducting wall [Miy1986]
            kp = big_c/(1 + half_wth_a/(3 * half_wth_b) + big_c)
            mhd_pressure_drop = kp * conduct_liq * vel * (b_mag**2) * l_channel

        !! If have perfcetly insulating FCIs...
        else

            !! Calculate pressure drop for (perfectly) insulating FCI [Mal1995]
            mhd_pressure_drop =  vel * b_mag * l_channel * sqrt(conduct_liq * vsc / half_wth_a)

        endif

        !! Total (Pa)
        liquid_breeder_pressure_drop_mhd = num_pol * mhd_pressure_drop

        if (ip==0) return
        call write_output_liquid_breeder_pressure_drop_mhd

        contains

        subroutine write_output_liquid_breeder_pressure_drop_mhd

            use process_output, only: oheadr, osubhd, ovarrf, ovarre, &
            ocmmnt, ovarin, ovarst

            implicit none

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

            call osubhd(ofile, 'Liquid metal breeder/coolant MHD pressure drop for ' // label)

            if (ifci==0) then
                call ocmmnt(ofile, 'Flow channels have thin conducting walls (ifci==0)')
                call ovarre(ofile, 'Wall conductance (A V-1 m-1)', '(bz_channel_conduct_liq)', bz_channel_conduct_liq, 'OP ')
            else if (ifci==2) then
                call ocmmnt(ofile, 'Flow Channel Inserts (FCIs) used (ifci==2)')
                call ovarre(ofile, 'FCI conductance (A V-1 m-1)', '(bz_channel_conduct_liq)', bz_channel_conduct_liq, 'OP ')
            else
                call ocmmnt(ofile, 'Flow Channel Inserts - assumed perfect insulator (ifci==1)')
            endif
            call ovarre(ofile, 'Length of long poloidal secion of channel (m)', '(l_channel)', l_channel, 'OP ')
            call ovarin(ofile, 'Number of long poloidal secions of channel', '(num_pol)', num_pol, 'OP ')
            call ovarre(ofile, 'MHD pressure drop (Pa)', '(liquid_breeder_pressure_drop_mhd)', liquid_breeder_pressure_drop_mhd, 'OP ')

        end subroutine write_output_liquid_breeder_pressure_drop_mhd

    end function liquid_breeder_pressure_drop_mhd