ifebdg Subroutine

public subroutine ifebdg(outfile, iprint)

Routine to calculate the volumes of the buildings required for an Inertial Fusion Energy power plant 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 calculates the volumes of the buildings required for an Inertial Fusion Energy power plant. The method is based closely on that for tokamaks etc. in routine BLDGS. F/MI/PJK/LOGBOOK12, p.87

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: outfile
integer, intent(in) :: iprint

Contents

Source Code


Source Code

  subroutine ifebdg(outfile,iprint)
    !! Routine to calculate the volumes of the buildings required for
    !! an Inertial Fusion Energy power plant
    !! 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 calculates the volumes of the buildings required for
    !! an Inertial Fusion Energy power plant. The method is based
    !! closely on that for tokamaks etc. in routine
    !! <A HREF="bldgs.html">BLDGS</A>.
    !! F/MI/PJK/LOGBOOK12, p.87
    !!
    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

    use buildings_variables, only: wrbi, rbwt, rbrt, fndt, trcl, hcwt, hccl, &
      wgt2, stcl, pibv, efloor, triv, conv, admv, shov, admvol, convol, elevol, &
      rbvol, rmbvol, shovol, volrci, wsvol, volnucb
    use process_output, only: oheadr, ovarre
    use fwbs_variables, only: whtshld
    use heat_transport_variables, only: helpow
    use ife_variables, only: zl7, zu7, r7, zl6, zu6, r6

    implicit none

    !  Arguments
    integer, intent(in) :: outfile,iprint

    !  Local variables

    real(dp) :: cran,dcl,dcw,fac2,fac3,hcl,hcw,hrbi,rbh,rbl,rbw, &
         rmbh,rmbl,rmbw,rwl,rww,shh,tch,tcl,tcw,wgts,wsa
    real(dp), save :: cryv,elev,rbv,rmbv,vrci,wsv

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

    !  Reactor building
    !  ================

    !  Total internal height
    hrbi = zl7 + zu7

    !  Distance from centre of device to wall
    wrbi = r7

    !  Internal volume (square floor)
    vrci = (2.0D0 * wrbi)**2 * hrbi

    !  External dimensions
    !  RBWT = wall thickness
    !  RBRT = roof thickness
    !  FNDT = foundation thickness
    rbw = 2.0D0 * (r7 + rbwt)
    rbl = rbw
    rbh = hrbi + rbrt + fndt

    !  External volume
    rbv = rbw * rbl * rbh

    !  Maintenance building
    !  ====================

    !  The reactor maintenance building includes the hot cells, the
    !  decontamination chamber, the transfer corridors, and the waste
    !  treatment building.  The dimensions of these areas are scaled
    !  from a reference (tokamak) design based on the shield sector size.

    !  Shield height

    shh = zl6 + zu6

    !  Transport corridor size

    tcw = r6 + 4.0D0*trcl
    tcl = 5.0D0*tcw + 2.0D0*hcwt

    !  Decontamination cell size

    dcw = 2.0D0*tcw + 1.0D0
    dcl = 2.0D0*tcw + 1.0D0

    !  Hot cell size

    hcw = r6 + 3.0D0*hccl + 2.0D0
    hcl = 3.0D0*r6 + 4.0D0*hccl + tcw

    !  Radioactive waste treatment

    rww = dcw
    rwl = hcl-dcl-hcwt

    !  Maintenance building dimensions

    rmbw = hcw + dcw + 3.0D0*hcwt
    rmbl = hcl + 2.0D0*hcwt

    !  Height

    if (wgt2 > 1.0D0) then
       wgts = wgt2
    else
       wgts = whtshld
    end if
    cran = 9.41D-6*wgts + 5.1D0
    rmbh = 10.0D0 + (zl6+zu6) + trcl + cran + 5.1D0 + stcl + fndt
    tch = shh + stcl + fndt

    !  Volume

    fac2 = 2.8D0
    rmbv = fac2*rmbw*rmbl*rmbh + tcw*tcl*tch

    !  Warm shop and hot cell gallery

    wsa = (rmbw+7.0D0)*20.0D0 + rmbl*7.0D0
    fac3 = 1.9D0
    wsv = fac3 * wsa*rmbh

    !  Cryogenic building volume

    cryv = 55.0D0 * sqrt(helpow)

    !  Electrical building volume
    !  (set equal to power injection (i.e. driver) building volume)

    elev = pibv

    !  Calculate effective floor area for ac power module

    efloor = (rbv+rmbv+wsv+triv+elev+conv+cryv+admv+shov)/6.0D0

    !  Convert local into global variables

    admvol = admv
    convol = conv
    elevol = elev
    rbvol  = rbv
    rmbvol = rmbv
    shovol = shov
    volrci = vrci
    wsvol  = wsv

    !  Total volume of nuclear buildings

    volnucb = ( vrci + rmbv + wsv + triv + cryv )

    if (iprint == 0) return

    !  Output section

    call oheadr(outfile,'Plant Buildings System')
    call ovarre(outfile,'Internal volume of reactor building (m3)', &
         '(vrci)',vrci)
    call ovarre(outfile,'Dist from device centre to bldg wall (m)', &
         '(wrbi)',wrbi)
    call ovarre(outfile,'Effective floor area (m2)','(efloor)',efloor)
    call ovarre(outfile,'Reactor building volume (m3)','(rbv)',rbv)
    call ovarre(outfile,'Reactor maintenance building volume (m3)', &
         '(rmbv)',rmbv)
    call ovarre(outfile,'Warmshop volume (m3)','(wsv)',wsv)
    call ovarre(outfile,'Tritium building volume (m3)','(triv)',triv)
    call ovarre(outfile,'Electrical building volume (m3)','(elev)',elev)
    call ovarre(outfile,'Control building volume (m3)','(conv)',conv)
    call ovarre(outfile,'Cryogenics building volume (m3)','(cryv)',cryv)
    call ovarre(outfile,'Administration building volume (m3)','(admv)', &
         admv)
    call ovarre(outfile,'Shops volume (m3)','(shov)',shov)
    call ovarre(outfile,'Total volume of nuclear buildings (m3)', &
         '(volnucb)',volnucb)

  end subroutine ifebdg