Routine to create the build of an inertial fusion energy device and to calculate the material volumes for the device core author: P J Knight, CCFE, Culham Science Centre outfile : input integer : output file unit iprint : input integer : switch for writing to output file (1=yes) This routine constructs the build of an inertial fusion energy device and calculates the material volumes for the device core. F/MI/PJK/LOGBOOK12, p.52
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | outfile | |||
integer, | intent(in) | :: | iprint |
subroutine ifebld(outfile,iprint)
!! Routine to create the build of an inertial fusion energy device
!! and to calculate the material volumes for the device core
!! author: P J Knight, CCFE, Culham Science Centre
!! outfile : input integer : output file unit
!! iprint : input integer : switch for writing to output file (1=yes)
!! This routine constructs the build of an inertial fusion energy device
!! and calculates the material volumes for the device core.
!! F/MI/PJK/LOGBOOK12, p.52
!!
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
use process_output, only: oheadr, obuild, ovarre
use ife_variables, only: ifetyp, chrad, r1, fwdr, r2, v1dr, r3, bldr, r4, &
v2dr, r5, shdr, r6, v3dr, r7, zl7, v3dzl, zl6, shdzl, zl5, v2dzl, zl4, &
bldzl, zl3, v1dzl, zl2, fwdzl, zl1, chdzl, chdzu, zu1, fwdzu, zu2, &
v1dzu, zu3, bldzu, zu4, v2dzu, zu5, shdzu, zu6, v3dzu, zu7, fwmatv, &
v1matv, blmatv, v2matv, shmatv, v3matv, chmatv
use constants, only: mfile
implicit none
! Arguments
integer, intent(in) :: outfile,iprint
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
select case (ifetyp)
case (1) ! OSIRIS-type device
call osibld
case (2) ! SOMBRERO-type device
call sombld
case (3) ! HYLIFE-II-type device
call hylbld
case (4) ! 2019 device
call bld2019
case default ! Generic device
call genbld
end select
if (iprint == 0) return
! Output section
call oheadr(outfile,'Radial Build')
write(outfile,20)
20 format(T43,'Thickness (m)',T60,'Radius (m)')
call obuild(outfile,'Device centreline',0.0D0,0.0D0)
call obuild(outfile,'Chamber',chrad,r1)
call ovarre(mfile,'Chamber (m)','(chrad)',chrad)
call obuild(outfile,'First Wall',fwdr,r2)
call ovarre(mfile,'First Wall (m)','(fwdr)',fwdr)
call obuild(outfile,'Void 1',v1dr,r3)
call ovarre(mfile,'Void 1 (m)','(v1dr)',v1dr)
call obuild(outfile,'Blanket',bldr,r4)
call ovarre(mfile,'Blanket (m)','(bldr)',bldr)
call obuild(outfile,'Void 2',v2dr,r5)
call ovarre(mfile,'Void 2 (m)','(v2dr)',v2dr)
call obuild(outfile,'Shield',shdr,r6)
call ovarre(mfile,'Shield (m)','(shdr)',shdr)
call obuild(outfile,'Void 3',v3dr,r7)
call ovarre(mfile,'Void 3 (m)','(v3dr)',v3dr)
30 format(T43,'Thickness (m)',T60,'Height (m)')
if (ifetyp /= 4) then
call oheadr(outfile,'Vertical Build')
write(outfile,30)
call obuild(outfile,'Base of device',0.0D0,-zl7)
call obuild(outfile,'Void 3',v3dzl,-zl6)
call obuild(outfile,'Shield',shdzl,-zl5)
call obuild(outfile,'Void 2',v2dzl,-zl4)
call obuild(outfile,'Blanket',bldzl,-zl3)
call obuild(outfile,'Void 1',v1dzl,-zl2)
call obuild(outfile,'First Wall',fwdzl,-zl1)
call obuild(outfile,'Chamber',chdzl,0.0D0)
call obuild(outfile,'Chamber',chdzu,zu1)
call obuild(outfile,'First Wall',fwdzu,zu2)
call obuild(outfile,'Void 1',v1dzu,zu3)
call obuild(outfile,'Blanket',bldzu,zu4)
call obuild(outfile,'Void 2',v2dzu,zu5)
call obuild(outfile,'Shield',shdzu,zu6)
call obuild(outfile,'Void 3',v3dzu,zu7)
else
call oheadr(outfile,'Vertical Build - Midplane')
write(outfile,30)
call obuild(outfile,'Base of device',0.0D0,-zl7)
call obuild(outfile,'Void 3',v3dzl,-zl6)
call obuild(outfile,'Shield',shdzl,-zl5)
call obuild(outfile,'Void 2',v2dzl,-zl4)
call obuild(outfile,'Blanket',bldzl,-zl3)
call obuild(outfile,'Void 1',v1dzl,-zl2)
call obuild(outfile,'First Wall',fwdzl,-zl1)
call obuild(outfile,'Chamber',chdzl,0.0D0)
call obuild(outfile,'Chamber',chdzu,zu1)
call obuild(outfile,'First Wall',fwdzu,zu2)
call obuild(outfile,'Void 1',v1dzu,zu3)
call obuild(outfile,'Blanket',bldzu-bldzu,zu4-bldzu)
call obuild(outfile,'Void 2',v2dzu,zu5-bldzu)
call obuild(outfile,'Shield',shdzu,zu6-bldzu)
call obuild(outfile,'Void 3',v3dzu+bldzu,zu7)
call oheadr(outfile,'Vertical Build - Edge')
write(outfile,30)
call obuild(outfile,'Base of device',0.0D0,-zl7)
call obuild(outfile,'Void 3',v3dzl,-zl6)
call ovarre(mfile,'Void 3 lower (m)','(v3dzl)',v3dzl)
call obuild(outfile,'Shield',shdzl,-zl5)
call ovarre(mfile,'Shield lower (m)','(shdzl)',shdzl)
call obuild(outfile,'Void 2',v2dzl,-zl4)
call ovarre(mfile,'Void 2 lower (m)','(v2dzl)',v2dzl)
call obuild(outfile,'Blanket',bldzl,-zl3)
call ovarre(mfile,'Blanket lower (m)','(bldzl)',bldzl)
call obuild(outfile,'Void 1',v1dzl,-zl2)
call ovarre(mfile,'Void 1 lower (m)','(v1dzl)',v1dzl)
call obuild(outfile,'First Wall',fwdzl,-zl1)
call ovarre(mfile,'First Wall lower (m)','(fwdzl)',fwdzl)
call obuild(outfile,'Chamber',chdzl,0.0D0)
call ovarre(mfile,'Chamber lower (m)','(chdzl)',chdzl)
call obuild(outfile,'Chamber',chdzu,zu1)
call ovarre(mfile,'Chamber upper (m)','(chdzu)',chdzu)
call obuild(outfile,'First Wall',fwdzu,zu2)
call ovarre(mfile,'First Wall upper (m)','(fwdzu)',fwdzu)
call obuild(outfile,'Void 1',v1dzu,zu3)
call ovarre(mfile,'Void 1 upper (m)','(v1dzu)',v1dzu)
call obuild(outfile,'Blanket',bldzu,zu4)
call ovarre(mfile,'Blanket upper (m)','(bldzu)',bldzu)
call obuild(outfile,'Void 2',v2dzu,zu5)
call ovarre(mfile,'Void 2 upper (m)','(v2dzu)',v2dzu)
call obuild(outfile,'Shield',shdzu,zu6)
call ovarre(mfile,'Shield upper (m)','(shdzu)',shdzu)
call obuild(outfile,'Void 3',v3dzu,zu7)
call ovarre(mfile,'Void 3 upper (m)','(v3dzu)',v3dzu)
end if
! Print matrix of material volumes
call oheadr(outfile,'Material volumes')
write(outfile,*) ' Chamber 1st wall Void 1 Blanket ' &
//' Void 2 Shield Void 3'
write(outfile,'(A9,7(1pe9.2))') 'void ', &
chmatv(0), &
(fwmatv(1,0)+fwmatv(2,0)+fwmatv(3,0)), &
(v1matv(1,0)+v1matv(2,0)+v1matv(3,0)), &
(blmatv(1,0)+blmatv(2,0)+blmatv(3,0)), &
(v2matv(1,0)+v2matv(2,0)+v2matv(3,0)), &
(shmatv(1,0)+shmatv(2,0)+shmatv(3,0)), &
(v3matv(1,0)+v3matv(2,0)+v3matv(3,0))
write(outfile,'(A9,7(1pe9.2))') 'steel ', &
chmatv(1), &
(fwmatv(1,1)+fwmatv(2,1)+fwmatv(3,1)), &
(v1matv(1,1)+v1matv(2,1)+v1matv(3,1)), &
(blmatv(1,1)+blmatv(2,1)+blmatv(3,1)), &
(v2matv(1,1)+v2matv(2,1)+v2matv(3,1)), &
(shmatv(1,1)+shmatv(2,1)+shmatv(3,1)), &
(v3matv(1,1)+v3matv(2,1)+v3matv(3,1))
write(outfile,'(A9,7(1pe9.2))') 'carbon ', &
chmatv(2), &
(fwmatv(1,2)+fwmatv(2,2)+fwmatv(3,2)), &
(v1matv(1,2)+v1matv(2,2)+v1matv(3,2)), &
(blmatv(1,2)+blmatv(2,2)+blmatv(3,2)), &
(v2matv(1,2)+v2matv(2,2)+v2matv(3,2)), &
(shmatv(1,2)+shmatv(2,2)+shmatv(3,2)), &
(v3matv(1,2)+v3matv(2,2)+v3matv(3,2))
write(outfile,'(A9,7(1pe9.2))') 'FLiBe ', &
chmatv(3), &
(fwmatv(1,3)+fwmatv(2,3)+fwmatv(3,3)), &
(v1matv(1,3)+v1matv(2,3)+v1matv(3,3)), &
(blmatv(1,3)+blmatv(2,3)+blmatv(3,3)), &
(v2matv(1,3)+v2matv(2,3)+v2matv(3,3)), &
(shmatv(1,3)+shmatv(2,3)+shmatv(3,3)), &
(v3matv(1,3)+v3matv(2,3)+v3matv(3,3))
write(outfile,'(A9,7(1pe9.2))') 'Li2O ', &
chmatv(4), &
(fwmatv(1,4)+fwmatv(2,4)+fwmatv(3,4)), &
(v1matv(1,4)+v1matv(2,4)+v1matv(3,4)), &
(blmatv(1,4)+blmatv(2,4)+blmatv(3,4)), &
(v2matv(1,4)+v2matv(2,4)+v2matv(3,4)), &
(shmatv(1,4)+shmatv(2,4)+shmatv(3,4)), &
(v3matv(1,4)+v3matv(2,4)+v3matv(3,4))
write(outfile,'(A9,7(1pe9.2))') 'concrete ', &
chmatv(5), &
(fwmatv(1,5)+fwmatv(2,5)+fwmatv(3,5)), &
(v1matv(1,5)+v1matv(2,5)+v1matv(3,5)), &
(blmatv(1,5)+blmatv(2,5)+blmatv(3,5)), &
(v2matv(1,5)+v2matv(2,5)+v2matv(3,5)), &
(shmatv(1,5)+shmatv(2,5)+shmatv(3,5)), &
(v3matv(1,5)+v3matv(2,5)+v3matv(3,5))
write(outfile,'(A9,7(1pe9.2))') 'helium ', &
chmatv(6), &
(fwmatv(1,6)+fwmatv(2,6)+fwmatv(3,6)), &
(v1matv(1,6)+v1matv(2,6)+v1matv(3,6)), &
(blmatv(1,6)+blmatv(2,6)+blmatv(3,6)), &
(v2matv(1,6)+v2matv(2,6)+v2matv(3,6)), &
(shmatv(1,6)+shmatv(2,6)+shmatv(3,6)), &
(v3matv(1,6)+v3matv(2,6)+v3matv(3,6))
write(outfile,'(A9,7(1pe9.2))') 'xenon ', &
chmatv(7), &
(fwmatv(1,7)+fwmatv(2,7)+fwmatv(3,7)), &
(v1matv(1,7)+v1matv(2,7)+v1matv(3,7)), &
(blmatv(1,7)+blmatv(2,7)+blmatv(3,7)), &
(v2matv(1,7)+v2matv(2,7)+v2matv(3,7)), &
(shmatv(1,7)+shmatv(2,7)+shmatv(3,7)), &
(v3matv(1,7)+v3matv(2,7)+v3matv(3,7))
write(outfile,'(A9,7(1pe9.2))') 'lithium ', &
chmatv(8), &
(fwmatv(1,8)+fwmatv(2,8)+fwmatv(3,8)), &
(v1matv(1,8)+v1matv(2,8)+v1matv(3,8)), &
(blmatv(1,8)+blmatv(2,8)+blmatv(3,8)), &
(v2matv(1,8)+v2matv(2,8)+v2matv(3,8)), &
(shmatv(1,8)+shmatv(2,8)+shmatv(3,8)), &
(v3matv(1,8)+v3matv(2,8)+v3matv(3,8))
contains
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine bld2019
!! Routine to create the build of a 2019 inertial fusion energy
!! device, and to calculate the material volumes for the device core
!! author: S I Muldrew, CCFE, Culham Science Centre
!! None
!! This routine constructs the build of a modern inertial fusion energy
!! device, assumed to be cylindrically-symmetric, with a pool at bottom
!! and top corners and with a lower shield at the centre. See diagram
!! attached to Issue #907.
!! Issue #907
!! !
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
use constants, only: pi
use build_variables, only: fwarea
use error_handling, only: report_error
use buildings_variables, only: trcl, stcl
use fwbs_variables, only: tbr, emult
use ife_variables, only: blmatf, bldrc, lipmw, etali, taufall, rrmax, chvol, &
fwvol, v1vol, blvol, v2vol, shvol, flirad, v3vol, maxmat, chmatf, fwmatf, &
v1matf, v2matf, shmatf, v3matf
implicit none
! Arguments
! Local variables
integer :: i,j
real(dp), save :: g, vel, acurt, mdot, phi, sang, li_frac
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Check input
if ((fwdr.gt.1.0D-9).or.(v1dr.gt.1.0d-9)) then
call report_error(230)
else if ((fwdzu.gt.1.0D-9).or.(v1dzu.gt.1.0D-9).or.(v2dzu.gt.1.0D-9)) then
call report_error(231)
else if ((fwdzl.gt.1.0D-9).or.(v1dzl.gt.1.0d-9).or.(v2dzu.gt.1.0D-9)) then
call report_error(232)
end if
! Lithium Pump
! Acceleration due to gravity (m/s^2)
g = 9.81D0
! Velocity
vel = sqrt(2.0D0*g*(chdzu+bldzu))
! Lithium Fraction
blmatf(1,8) = 0.91*sqrt(bldzu/(chdzu+bldzu))
blmatf(1,0) = 1.0D0 - blmatf(1,8)
! Spatial Thickness
bldr = bldrc / blmatf(1,8)
! Area
acurt = pi * ((chrad+bldr)**2.0D0-chrad**2.0D0)
! Mass Flow
mdot = 512.0D0 * vel * blmatf(1,8) * acurt
! Pump Power (MW)
lipmw = 1.0D-6 * mdot * g * (chdzl+chdzu+bldzu+bldzl) / etali
! Fall Time
taufall = 2.0D0 * (chdzl+chdzu+bldzu) / vel
rrmax = 1.0D0 / taufall
! TBR and Emult model was for spherical lithium
! Remove reactor head
phi = atan(r1/zu1)
sang = 1.0D0 - cos(phi)
li_frac = 1.0D0 - 0.5D0 * sang
! TBR
tbr = 3.7418D0 * (1.0D0 / (1.0D0 + exp(-2.6366D0 * bldrc)) - 0.5D0)
tbr = tbr * li_frac
! Energy Multiplication
emult = 2.2414D0 * (1.0D0 / (1.0D0 + exp(-3.0038D0 * bldrc)) - 0.5D0)
emult = emult * li_frac
! Radial build
r1 = chrad
r2 = r1 + fwdr
r3 = r2 + v1dr
r4 = r3 + bldr
r5 = r4 + v2dr
r6 = r5 + shdr
r7 = r6 + v3dr
! Vertical build (below midplane)
zl1 = chdzl
zl2 = zl1 + fwdzl
zl3 = zl2 + v1dzl
zl4 = zl3 + bldzl
zl5 = zl4 + v2dzl
zl6 = zl5 + shdzl
zl7 = zl6 + v3dzl
! Vertical build (above midplane)
zu1 = chdzu
zu2 = zu1 + fwdzu
zu3 = zu2 + v1dzu
zu4 = zu3 + bldzu
zu5 = zu4 + v2dzu
zu6 = zu5 + shdzu
v3dzu = (zu6+zl6) + trcl + stcl + 5.1D0 + 9.41D-6 * 1.0D5
zu7 = zu6 + v3dzu
! Component volumes
! The following notation applies below:
! J=1 : side part
! J=2 : top part
! J=3 : bottom part
! Chamber
chvol = pi * r1*r1 * (zu1 + zl1)
! First wall
fwvol(1) = pi * (r2*r2 - r1*r1) * (zu1 + zl1)
fwvol(2) = pi * r2*r2 * (zu2 - zu1)
fwvol(3) = pi * r2*r2 * (zl2 - zl1)
! First void
v1vol(1) = pi * (r3*r3 - r2*r2) * (zu2 + zl2)
v1vol(2) = pi * r3*r3 * (zu3 - zu2)
v1vol(3) = pi * r3*r3 * (zl3 - zl2)
! Blanket
! Radial Blanket - between void 2 and chamber
blvol(1) = pi * (r4*r4 - r3*r3) * (zu3 + zl3)
! Upper Blanket - Pool radially between shield and
! chamber of input height.
blvol(2) = pi * (r5*r5 - r3*r3) * bldzu
! Lower Blanket - Pool filling base of device
blvol(3) = pi * r5*r5 * (zl4 - zl3)
! Second void
v2vol(1) = pi * (r5*r5 - r4*r4) * (chdzl+chdzu)
v2vol(2) = 0.0D0
v2vol(3) = 0.0D0
! Shield
shvol(1) = pi * (r6*r6 - r5*r5) * (zu5 + zl5)
! Top Section is in three parts to account for the dip at
! the centre. The first is the horizontal top, the second is the
! horizontal
shvol(2) = pi * (((r6*r6 - (chrad-shdr)*(chrad-shdr)) * shdzu) &
+ ((r1*r1 - flirad*flirad) * shdzu) &
+ ( (r1*r1 - (r1-shdzu)*(r1-shdzu)) *(bldzu-shdzu) ))
shvol(3) = pi * r6*r6 * (zl6 - zl5)
! Third void
v3vol(1) = pi * (r7*r7 - r6*r6) * (zu6 + zl6)
v3vol(2) = pi * r7*r7 * (zu7 - zu6) + pi * ((r1-shdzu)*(r1-shdzu) - flirad*flirad) * bldzu
v3vol(3) = pi * r7*r7 * (zl7 - zl6)
! Material volumes
do i = 0,maxmat
chmatv(i) = max(0.0D0, chvol * chmatf(i))
do j = 1,3
fwmatv(j,i) = max(0.0D0, fwvol(j) * fwmatf(j,i))
v1matv(j,i) = max(0.0D0, v1vol(j) * v1matf(j,i))
blmatv(j,i) = max(0.0D0, blvol(j) * blmatf(j,i))
v2matv(j,i) = max(0.0D0, v2vol(j) * v2matf(j,i))
shmatv(j,i) = max(0.0D0, shvol(j) * shmatf(j,i))
v3matv(j,i) = max(0.0D0, v3vol(j) * v3matf(j,i))
end do
end do
! First wall area
! The chamber is surrounded by liquid on three sides
! with only the top being solid. This is considered part
! of the shield. There is a target injector tube at the
! centre of this area.
fwarea = pi * (r1*r1 - flirad*flirad)
end subroutine bld2019
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine genbld
!! Routine to create the build of a generic inertial fusion energy
!! device, and to calculate the material volumes for the device core
!! author: P J Knight, CCFE, Culham Science Centre
!! None
!! This routine constructs the build of a generic inertial fusion energy
!! device, assumed to be cylindrically-symmetric, and to calculate
!! the material volumes for the device core.
!! F/MI/PJK/LOGBOOK12, p.52
!! !
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
use constants, only: pi
use build_variables, only: fwarea
use ife_variables, only: chvol, fwvol, v1vol, blvol, v2vol, shvol, v3vol, &
maxmat, chmatf, fwmatf, v1matf, blmatf, v2matf, shmatf, v3matf
implicit none
! Local variables
integer :: i,j
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Radial build
r1 = chrad
r2 = r1 + fwdr
r3 = r2 + v1dr
r4 = r3 + bldr
r5 = r4 + v2dr
r6 = r5 + shdr
r7 = r6 + v3dr
! Vertical build (below midplane)
zl1 = chdzl
zl2 = zl1 + fwdzl
zl3 = zl2 + v1dzl
zl4 = zl3 + bldzl
zl5 = zl4 + v2dzl
zl6 = zl5 + shdzl
zl7 = zl6 + v3dzl
! Vertical build (above midplane)
zu1 = chdzu
zu2 = zu1 + fwdzu
zu3 = zu2 + v1dzu
zu4 = zu3 + bldzu
zu5 = zu4 + v2dzu
zu6 = zu5 + shdzu
zu7 = zu6 + v3dzu
! Component volumes
! The following notation applies below:
! J=1 : side part
! J=2 : top part
! J=3 : bottom part
! Chamber
chvol = pi * r1*r1 * (zu1 + zl1)
! First wall
fwvol(1) = pi * (r2*r2 - r1*r1) * (zu1 + zl1)
fwvol(2) = pi * r2*r2 * (zu2 - zu1)
fwvol(3) = pi * r2*r2 * (zl2 - zl1)
! First void
v1vol(1) = pi * (r3*r3 - r2*r2) * (zu2 + zl2)
v1vol(2) = pi * r3*r3 * (zu3 - zu2)
v1vol(3) = pi * r3*r3 * (zl3 - zl2)
! Blanket
blvol(1) = pi * (r4*r4 - r3*r3) * (zu3 + zl3)
blvol(2) = pi * r4*r4 * (zu4 - zu3)
blvol(3) = pi * r4*r4 * (zl4 - zl3)
! Second void
v2vol(1) = pi * (r5*r5 - r4*r4) * (zu4 + zl4)
v2vol(2) = pi * r5*r5 * (zu5 - zu4)
v2vol(3) = pi * r5*r5 * (zl5 - zl4)
! Shield
shvol(1) = pi * (r6*r6 - r5*r5) * (zu5 + zl5)
shvol(2) = pi * r6*r6 * (zu6 - zu5)
shvol(3) = pi * r6*r6 * (zl6 - zl5)
! Third void
v3vol(1) = pi * (r7*r7 - r6*r6) * (zu6 + zl6)
v3vol(2) = pi * r7*r7 * (zu7 - zu6)
v3vol(3) = pi * r7*r7 * (zl7 - zl6)
! Material volumes
do i = 0,maxmat
chmatv(i) = max(0.0D0, chvol * chmatf(i))
do j = 1,3
fwmatv(j,i) = max(0.0D0, fwvol(j) * fwmatf(j,i))
v1matv(j,i) = max(0.0D0, v1vol(j) * v1matf(j,i))
blmatv(j,i) = max(0.0D0, blvol(j) * blmatf(j,i))
v2matv(j,i) = max(0.0D0, v2vol(j) * v2matf(j,i))
shmatv(j,i) = max(0.0D0, shvol(j) * shmatf(j,i))
v3matv(j,i) = max(0.0D0, v3vol(j) * v3matf(j,i))
end do
end do
! First wall area
fwarea = 2.0D0*pi*r1 * ((zu1 + zl1) + r1)
end subroutine genbld
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine osibld
!! Routine to create the build of an inertial fusion energy
!! device, based on the design of the OSIRIS study,
!! and to calculate the material volumes for the device core
!! author: P J Knight, CCFE, Culham Science Centre
!! None
!! This routine constructs the build of an inertial fusion energy
!! device, based on the design of the OSIRIS study, and to calculate
!! the material volumes for the device core.
!! F/MI/PJK/LOGBOOK12, p.56
!! !
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
use constants, only: pi
use build_variables, only: fwarea
implicit none
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Careful choice of thicknesses, and assuming that the FLiBe
! inlet radius is small, allows the generic build calculation
! to be roughly applicable.
call genbld
! First wall area: no true first wall at bottom of chamber
fwarea = 2.0D0*pi*r1*(zu1 + zl1) + pi*r1*r1
end subroutine osibld
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine sombld
!! Routine to create the build of an inertial fusion energy
!! device, based on the design of the SOMBRERO study,
!! and to calculate the material volumes for the device core
!! author: P J Knight, CCFE, Culham Science Centre
!! None
!! This routine constructs the build of an inertial fusion energy
!! device, based on the design of the SOMBRERO study, and to calculate
!! the material volumes for the device core.
!! Sviatoslavsky et al, Fusion Technology vol.21 (1992) 1470
!! F/MI/PJK/LOGBOOK12, p.53
!! !
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
use constants, only: pi
use build_variables, only: fwarea
use ife_variables, only: chvol, fwvol, v1vol, blvol, v2vol, somtdr, sombdr, &
shvol, v3vol, maxmat, chmatf, fwmatf, blmatf, v1matf, v2matf, shmatf, v3matf
implicit none
! Local variables
real(dp), parameter :: third = 1.0D0/3.0D0
real(dp) :: chcylh,ddz,dvol
integer :: i,j
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Radial build
r1 = chrad
r2 = r1 + fwdr
r3 = r2 + v1dr
r4 = r3 + bldr
r5 = r4 + v2dr
r6 = r5 + shdr
r7 = r6 + v3dr
! Vertical build (below midplane)
zl1 = chdzl
zl2 = zl1 + fwdzl
zl3 = zl2 + v1dzl
zl4 = zl3 + bldzl
zl5 = zl4 + v2dzl
zl6 = zl5 + shdzl
zl7 = zl6 + v3dzl
! Vertical build (above midplane)
zu1 = chdzu
zu2 = zu1 + fwdzu
zu3 = zu2 + v1dzu
zu4 = zu3 + bldzu
zu5 = zu4 + v2dzu
zu6 = zu5 + shdzu
zu7 = zu6 + v3dzu
! The SOMBRERO chamber is made up of a cylindrical first wall/
! blanket, with conical regions above and below. Outside this is
! a cylindrical shield.
! Component volumes
! The following notation applies below:
! J=1 : side part
! J=2 : top part
! J=3 : bottom part
! Chamber : CHCYLH is the height of the cylindrical part
chcylh = chdzu + chdzl - 2.0D0*chrad
chvol = pi * r1*r1 * (chcylh + 2.0D0*third * chrad)
! First wall
fwvol(1) = pi * (r2*r2 - r1*r1) * chcylh
fwvol(2) = third * pi * ( r2*r2*(chrad+fwdzu) - r1*r1*chrad )
fwvol(3) = third * pi * ( r2*r2*(chrad+fwdzl) - r1*r1*chrad )
! First void
v1vol(1) = pi * (r3*r3 - r2*r2) * chcylh
v1vol(2) = third * pi * ( r3*r3*(chrad+fwdzu+v1dzu) - &
r2*r2*(chrad+fwdzu) )
v1vol(3) = third * pi * ( r3*r3*(chrad+fwdzl+v1dzl) - &
r2*r2*(chrad+fwdzl) )
! Blanket: SOMTDR and SOMBDR are the radii of the cylindrical
! sections at the top/bottom of the blanket
! DDZ = Height of top cylindrical section (by similar triangles)
! DVOL = Volume of top cylindrical section, less the internal cone
blvol(1) = pi * (r4*r4 - r3*r3) * chcylh
blvol(2) = third * pi * ( r4*r4*(chrad+fwdzu+v1dzu+bldzu) - &
r3*r3*(chrad+fwdzu+v1dzu) )
ddz = (chrad+fwdzu+v1dzu+bldzu)/(chrad+fwdr+v1dr+bldr)*somtdr
dvol = 2.0D0*third * pi * somtdr*somtdr * ddz
blvol(2) = blvol(2) + dvol
! Ditto for bottom region...
blvol(3) = third * pi * ( r4*r4*(chrad+fwdzl+v1dzl+bldzl) - &
r3*r3*(chrad+fwdzl+v1dzl) )
ddz = (chrad+fwdzl+v1dzl+bldzl)/(chrad+fwdr+v1dr+bldr)*sombdr
dvol = 2.0D0*third * pi * sombdr*sombdr * ddz
blvol(3) = blvol(3) + dvol
! Second void
v2vol(1) = pi * (r5*r5 - r4*r4) * chcylh
v2vol(2) = pi * r5*r5 * (zu5 - chdzu + chrad) - ( &
fwvol(2) + v1vol(2) + blvol(2) + (third*pi*r1*r1*chrad) )
v2vol(3) = pi * r5*r5 * (zl5 - chdzl + chrad) - ( &
fwvol(3) + v1vol(3) + blvol(3) + (third*pi*r1*r1*chrad) )
! Shield
shvol(1) = pi * (r6*r6 - r5*r5) * (zu6 + zl6)
shvol(2) = pi * r5*r5 * (zu6 - zu5)
shvol(3) = pi * r5*r5 * (zl6 - zl5)
! Third void
v3vol(1) = pi * (r7*r7 - r6*r6) * (zu7 + zl7)
v3vol(2) = pi * r6*r6 * (zu7 - zu6)
v3vol(3) = pi * r6*r6 * (zl7 - zl6)
! Material volumes
do i = 0,maxmat
chmatv(i) = max(0.0D0, chvol * chmatf(i))
do j = 1,3
fwmatv(j,i) = max(0.0D0, fwvol(j) * fwmatf(j,i))
v1matv(j,i) = max(0.0D0, v1vol(j) * v1matf(j,i))
blmatv(j,i) = max(0.0D0, blvol(j) * blmatf(j,i))
v2matv(j,i) = max(0.0D0, v2vol(j) * v2matf(j,i))
shmatv(j,i) = max(0.0D0, shvol(j) * shmatf(j,i))
v3matv(j,i) = max(0.0D0, v3vol(j) * v3matf(j,i))
end do
end do
! First wall area
fwarea = 2.0D0*pi*r1*( (zu1 + zl1) + r1*sqrt(2.0D0) )
end subroutine sombld
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine hylbld
!! Routine to create the build of an inertial fusion energy
!! device, based on the design of the HYLIFE-II study,
!! and to calculate the material volumes for the device core
!! author: P J Knight, CCFE, Culham Science Centre
!! None
!! This routine constructs the build of an inertial fusion energy
!! device, based on the design of the HYLIFE-II study, and to calculate
!! the material volumes for the device core.
!! Moir, Fusion Technology vol.21 (1992) 1475
!! F/MI/PJK/LOGBOOK12, p.57
!! !
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
use constants, only: pi
use build_variables, only: fwarea
use ife_variables, only: chvol, fwvol, v1vol, blvol, v2vol, shvol, v3vol, &
maxmat, chmatf, fwmatf, blmatf, v1matf, v2matf, shmatf, v3matf, flirad
implicit none
! Local variables
real(dp), parameter :: third = 1.0D0/3.0D0
real(dp) :: chcylh,ddz,dvol
integer :: i,j
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Radial build
r1 = chrad
r2 = r1 + fwdr
r3 = r2 + v1dr
r4 = r3 + bldr
r5 = r4 + v2dr
r6 = r5 + shdr
r7 = r6 + v3dr
! Vertical build (below midplane)
zl1 = chdzl
zl2 = zl1 + fwdzl
zl3 = zl2 + v1dzl
zl4 = zl3 + bldzl
zl5 = zl4 + v2dzl
zl6 = zl5 + shdzl
zl7 = zl6 + v3dzl
! Vertical build (above midplane)
zu1 = chdzu
zu2 = zu1 + fwdzu
zu3 = zu2 + v1dzu
zu4 = zu3 + bldzu
zu5 = zu4 + v2dzu
zu6 = zu5 + shdzu
zu7 = zu6 + v3dzu
! The HYLIFE-II chamber is assumed to be mostly cylindrical, but
! with a conical region below the midplane that causes the Flibe
! to flow downwards and outwards towards the outlet.
! Component volumes
! The following notation applies below:
! J=1 : side part
! J=2 : top part
! J=3 : bottom part
chvol = pi * r1*r1 * ( (zu1 + zl5) - third*(zl5 - zl1) )
! First wall
! FLIRAD is the radius of the Flibe inlet
fwvol(1) = pi * (r2*r2 - r1*r1) * (zu2 + zl5)
fwvol(2) = pi * (r1*r1 - flirad*flirad) * (zu2 - zu1)
fwvol(3) = third * pi * ( r2*r2*(zl5-zl1) - r1*r1*(zl5-zl2) )
! First void
v1vol(1) = pi * (r3*r3 - r2*r2) * (zu2 + zl3)
v1vol(2) = pi * (r4*r4 - flirad*flirad) * (zu3 - zu2)
v1vol(3) = third * pi * r1*r1 * (zl3 - zl2)
! Blanket
blvol(1) = pi * (r4*r4 - r3*r3) * (zu2 + zl3)
blvol(2) = pi * (r4*r4 - flirad*flirad) * (zu4 - zu3)
blvol(3) = pi * r4*r4 * (zl4 - zl3)
! Second void
v2vol(1) = pi * (r5*r5 - r4*r4) * (zu4 + zl4)
v2vol(2) = pi * (r5*r5 - flirad*flirad) * (zu5 - zu4)
v2vol(3) = pi * r5*r5 * (zl5 - zl4)
! Shield
shvol(1) = pi * (r6*r6 - r5*r5) * (zu5 + zl5)
shvol(2) = pi * r6*r6 * (zu6 - zu5)
shvol(3) = pi * r6*r6 * (zl6 - zl5)
! Third void
v3vol(1) = pi * (r7*r7 - r6*r6) * (zu6 + zl6)
v3vol(2) = pi * r7*r7 * (zu7 - zu6)
v3vol(3) = pi * r7*r7 * (zl7 - zl6)
! Material volumes
do i = 0,maxmat
chmatv(i) = max(0.0D0, chvol * chmatf(i))
do j = 1,3
fwmatv(j,i) = max(0.0D0, fwvol(j) * fwmatf(j,i))
v1matv(j,i) = max(0.0D0, v1vol(j) * v1matf(j,i))
blmatv(j,i) = max(0.0D0, blvol(j) * blmatf(j,i))
v2matv(j,i) = max(0.0D0, v2vol(j) * v2matf(j,i))
shmatv(j,i) = max(0.0D0, shvol(j) * shmatf(j,i))
v3matv(j,i) = max(0.0D0, v3vol(j) * v3matf(j,i))
end do
end do
! First wall area
fwarea = 2.0D0 * pi * r1 * (zu1 + zl5)
fwarea = fwarea + pi * (r1*r1 - flirad*flirad)
fwarea = fwarea + pi * r1 * sqrt(r1*r1 + (zl3-zl1)**2)
end subroutine hylbld
end subroutine ifebld