Pressure drops are calculated for a pipe with a number of 90 and 180 degree bends. The pressure drop due to frictional forces along the total straight length of the pipe is calculated, then the pressure drop due to the bends is calculated. The total pressure drop is the sum of all contributions.
original author: P. J. Knight, CCFE moved from previous version of pumppower function by: G Graham, CCFE
N.B Darcy-Weisbach Equation (straight pipe):
kstrght = lambda * L/D
pressure drop = kstrght * (rho*V^2)/2
lambda - Darcy friction factor, L - pipe length, D - hydraulic diameter, rho - fluid density, V - fluid flow average velocity
This function also calculates pressure drop equations for elbow bends, with modified coefficients.
N.B. Darcy friction factor is estimated from the Haaland approximation.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
Function return parameter !!!!!
Pressure drop along the pipe (Pa) Reynolds number Calculate Darcy friction factor N.B. friction function Uses Haaland approx. which assumes a filled circular pipe. Use dh which allows us to do fluid calculations for non-cicular tubes (dh is estimate appropriate for fully developed flow). Pressure drop coefficient !!!!!
Straight section In preveious version of pumppower: - elbow radius assumed = 0.018m for 90 degree elbow, from WCLL - elbow radius assumed half that of 90 deg case for 180 deg elbow Intialised value for afw is 0.006m, so elbow radius = 3 * afw, aka 1.5 * pipe diameter, which seems to be engineering standard for a steel pipe long-radius elbow (short-radius elbow = 2 * afw).
If primary coolant... If secondary coolant... See DCLL 90 degree elbow pressure drop coefficient 180 degree elbow pressure drop coefficient Total (Pa)
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | ip | Swicth for primary or secondary coolant - =1 primary - =2 secondary |
||
integer, | intent(in) | :: | ofile | Swicth for primary or secondary coolant - =1 primary - =2 secondary |
||
integer, | intent(in) | :: | i_ps | Number of 90 and 180 degree bends in pipe |
||
integer, | intent(in) | :: | num_90 | Total flow length along pipe (m) |
||
integer, | intent(in) | :: | num_180 | Total flow length along pipe (m) |
||
real(kind=dp), | intent(in) | :: | l_pipe | Coolant density (kg/m3) |
||
real(kind=dp), | intent(in) | :: | den | Coolant viscosity (Pa s) |
||
real(kind=dp), | intent(in) | :: | vsc | Coolant flow velocity (m/s) |
||
real(kind=dp), | intent(in) | :: | vv | Description of this calculation |
||
character(len=*), | intent(in) | :: | label | Local variables !!!!!!!!!!!!!!! Hydraulic diameter of coolant flow channels (m) |
Arguments !!!!!!!!!!!!!!!!!!!!!
function pressure_drop(ip, ofile, i_ps, num_90, num_180, l_pipe, den, vsc, vv, label)
!! Pressure drops are calculated for a pipe with a number of 90
!! and 180 degree bends. The pressure drop due to frictional forces along
!! the total straight length of the pipe is calculated, then the pressure
!! drop due to the bends is calculated. The total pressure drop is the sum
!! of all contributions.
!!
!! original author: P. J. Knight, CCFE
!! moved from previous version of pumppower function by: G Graham, CCFE
!!
!! N.B Darcy-Weisbach Equation (straight pipe):
!!
!! kstrght = lambda * L/D
!!
!! pressure drop = kstrght * (rho*V^2)/2
!!
!! lambda - Darcy friction factor, L - pipe length, D - hydraulic diameter,
!! rho - fluid density, V - fluid flow average velocity
!!
!! This function also calculates pressure drop equations for elbow bends,
!! with modified coefficients.
!!
!! N.B. Darcy friction factor is estimated from the Haaland approximation.
!!
!! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
use fwbs_variables, only: afw, a_bz_liq, b_bz_liq
use fw_module, only: friction
implicit none
!! Function return parameter !!!!!
!! Pressure drop along the pipe (Pa)
real(dp) :: pressure_drop
!! Arguments !!!!!!!!!!!!!!!!!!!!!
integer, intent(in) :: ip, ofile
!! Swicth for primary or secondary coolant
!! - =1 primary
!! - =2 secondary
integer, intent(in) :: i_ps
!! Number of 90 and 180 degree bends in pipe
integer, intent(in) :: num_90, num_180
!! Total flow length along pipe (m)
real(dp), intent(in) :: l_pipe
!! Coolant density (kg/m3)
real(dp), intent(in) :: den
!! Coolant viscosity (Pa s)
real(dp), intent(in) :: vsc
!! Coolant flow velocity (m/s)
real(dp), intent(in) :: vv
!! Description of this calculation
character(len=*), intent(in) :: label
!! Local variables !!!!!!!!!!!!!!!
!! Hydraulic diameter of coolant flow channels (m)
real(dp) :: dh
!! Reynolds number
real(dp) :: reyn
!! Darcy friction factor
real(dp) :: lambda
!! Pressure drops for straight setions, 90 bends and 180 bends (Pa)
real(dp) :: pdropstraight, pdrop90, pdrop180
!! Elbow radius for 90 and 180 deg brends (m)
real(dp) :: elbow_radius
!! Pressure drop coefficients
real(dp) :: kelbwn, kelbwt, kstrght
!! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!! Calculate hydraulic dimater for round or retancular pipe (m)
dh = hydraulic_diameter(i_ps)
!! Reynolds number
reyn = den * vv * dh / vsc
!! Calculate Darcy friction factor
!! N.B. friction function Uses Haaland approx. which assumes a filled circular pipe.
!! Use dh which allows us to do fluid calculations for non-cicular tubes
!! (dh is estimate appropriate for fully developed flow).
call friction(reyn,lambda)
!! Pressure drop coefficient !!!!!
!! Straight section
kstrght = lambda * l_pipe/dh
!! In preveious version of pumppower:
!! - elbow radius assumed = 0.018m for 90 degree elbow, from WCLL
!! - elbow radius assumed half that of 90 deg case for 180 deg elbow
!! Intialised value for afw is 0.006m, so elbow radius = 3 * afw,
!! aka 1.5 * pipe diameter, which seems to be engineering standard for
!! a steel pipe long-radius elbow (short-radius elbow = 2 * afw).
!! If primary coolant...
if (i_ps==1) then
elbow_radius = 3 * afw
!! If secondary coolant...
else
!! See DCLL
elbow_radius = b_bz_liq
endif
!! 90 degree elbow pressure drop coefficient
kelbwn = elbow_coeff(elbow_radius, 90.0D0, lambda, dh)
!! 180 degree elbow pressure drop coefficient
kelbwt = elbow_coeff(elbow_radius/2, 180.0D0, lambda, dh)
!! Total (Pa)
pdropstraight = kstrght * 0.5D0*den*vv*vv
pdrop90 = num_90*kelbwn * 0.5D0*den*vv*vv
pdrop180 = num_180*kelbwt * 0.5D0*den*vv*vv
pressure_drop = pdropstraight + pdrop90 + pdrop180
if (ip==0) return
call write_output_pressure_drop
contains
subroutine write_output_pressure_drop
use process_output, only: oheadr, osubhd, ovarrf, ovarre, &
ocmmnt, ovarin, ovarst, oblnkl
use global_variables, only: verbose
implicit none
!! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
call osubhd(ofile, 'Pressure drop (friction) for ' // label)
call ovarre(ofile, 'Reynolds number', '(reyn)', reyn, 'OP ')
call ovarre(ofile, 'Darcy friction factor', '(lambda)', lambda, 'OP ')
call ovarre(ofile, 'Pressure drop (Pa)', '(pressure_drop)', pressure_drop, 'OP ')
call ocmmnt(ofile, 'This is the sum of the following:')
call ovarre(ofile, ' Straight sections (Pa)', '(pdropstraight)', &
pdropstraight, 'OP ')
call ovarre(ofile, ' 90 degree bends (Pa)', '(pdrop90)', pdrop90, 'OP ')
call ovarre(ofile, ' 180 degree bends (Pa)', '(pdrop180)', pdrop180, 'OP ')
call ocmmnt(ofile, 'Additional information is printed when verbose = 1')
if (verbose==1) then
call oblnkl(ofile)
call ovarre(ofile, 'Straight section pressure drop coefficient', &
'(kstrght)', kstrght, 'OP ')
call ovarre(ofile, '90 degree elbow coefficient', '(kelbwn)', kelbwn, 'OP ')
call ovarre(ofile, '180 degree elbow coefficient coefficient', '(kelbwt)', kelbwt, 'OP ')
end if
end subroutine write_output_pressure_drop
end function pressure_drop