ifefbs Subroutine

public subroutine ifefbs(outfile, iprint)

Routine to calculate the first wall, blanket and shield volumes, masses and other parameters, for an Inertial Fusion Energy device 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 first wall, blanket and shield volumes, masses and other parameters, for an Inertial Fusion Energy device. F/MI/PJK/LOGBOOK12, p.86 Moir et al., Fusion Technology, vol.25 (1994) p.5

Arguments

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

Contents

Source Code


Source Code

  subroutine ifefbs(outfile,iprint)
    !! Routine to calculate the first wall, blanket and shield volumes,
    !! masses and other parameters, for an Inertial Fusion Energy device
    !! 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 first wall, blanket and shield volumes,
    !! masses and other parameters, for an Inertial Fusion Energy device.
    !! F/MI/PJK/LOGBOOK12, p.86
    !! Moir et al., Fusion Technology, vol.25 (1994) p.5
    !!
    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

    use build_variables, only: fwarea
    use error_handling, only: report_error, fdiags
    use process_output, only: ovarre, oheadr
    use cost_variables, only: tlife, abktflnc, cfactr
    use fwbs_variables, only: denstl, fwmass, whtblkt, whtshld, whtblbe, whtblvd, &
      whtblss, wtblli2o, whtblli, bktlife, fwlife
    use ife_variables, only: maxmat, chmatm, chmatv, fwmatm, fwmatv, v1matm, &
      v1matv, blmatm, blmatv, v2matm, v2matv, shmatm, shmatv, v3matm, v3matv, &
      mflibe, fbreed, ifetyp
    use physics_variables, only: wallmw

    implicit none

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

    !  Local variables

    real(dp) :: den,life
    real(dp), dimension(0:maxmat) :: matden
    integer :: i,j

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

    !  Material densities
    !  0 = void
    !  1 = steel
    !  2 = carbon
    !  3 = FLiBe (inferred from Moir et al)
    !  4 = Li2O
    !  5 = concrete
    !  6 = helium (at typical coolant temperatures)
    !  7 = xenon (taken as ten times the normal tabulated value)
    !  8 = lithium (liquid, source Wikipedia)

    matden(0) = 0.0D0
    matden(1) = denstl
    matden(2) = 2300.0D0
    matden(3) = 2020.0D0
    matden(4) = 2010.0D0
    matden(5) = 2400.0D0
    matden(6) = 1.517D0
    matden(7) = 55.0D0
    matden(8) = 512.0D0

    !  Material masses
    do i = 0,maxmat
      den = matden(i)
      chmatm(i) = chmatv(i) * den
      do j = 1,3
        fwmatm(j,i) = fwmatv(j,i) * den
        v1matm(j,i) = v1matv(j,i) * den
        blmatm(j,i) = blmatv(j,i) * den
        v2matm(j,i) = v2matv(j,i) * den
        shmatm(j,i) = shmatv(j,i) * den
        v3matm(j,i) = v3matv(j,i) * den
      end do
    end do

    !  Total masses of components (excluding coolant)
    fwmass = 0.0D0
    whtblkt = 0.0D0
    whtshld = 0.0D0
    do i = 1,5
      do j = 1,3
        fwmass = fwmass + fwmatm(j,i)
        whtblkt = whtblkt + blmatm(j,i)
        whtshld = whtshld + shmatm(j,i)
      end do
    end do

    !  Other masses
    whtblbe = 0.0D0
    whtblvd = 0.0D0
    whtblss = 0.0D0
    wtblli2o = 0.0D0
    whtblli = 0.0D0
    do j = 1,3
      whtblss = whtblss + blmatm(j,1)
      wtblli2o = wtblli2o + blmatm(j,4)
      whtblli = whtblli + blmatm(j,8)
    end do

    !  Total mass of FLiBe
    mflibe = chmatm(3)
    do j = 1,3
      mflibe = mflibe + fwmatm(j,3) + v1matm(j,3) + blmatm(j,3) + &
        v2matm(j,3) + shmatm(j,3) + v3matm(j,3)
    end do

    !  A fraction FBREED of the total breeder inventory is outside the
    !  core region, i.e. is in the rest of the heat transport system
    if ((fbreed < 0.0D0).or.(fbreed > 0.999D0)) then
      fdiags(1) = fbreed ; call report_error(26)
    end if

    !  Following assumes that use of FLiBe and Li2O are
    !  mutually exclusive
    mflibe = mflibe / (1.0D0 - fbreed)
    wtblli2o = wtblli2o / (1.0D0 - fbreed)
    whtblli = whtblli / (1.0D0 - fbreed)

    !  Blanket and first wall lifetimes (HYLIFE-II: = plant life)
    if ((ifetyp == 3).or.(ifetyp == 4)) THEN
      life = tlife
    else
      life = min( tlife, abktflnc/(wallmw*cfactr) )
    end if

    bktlife = life
    fwlife = life

    !  Cryostat mass (=zero)
    !cryomass = 0.0D0

    if (iprint == 0) return

    !  Output section
    call oheadr(outfile,'First Wall, Blanket, Shield')
    call ovarre(outfile,'First wall area (m2)','(fwarea)',fwarea)
    call ovarre(outfile,'First wall mass (kg)','(fwmass)',fwmass)
    call ovarre(outfile,'Blanket mass (kg)','(whtblkt)',whtblkt)
    call ovarre(outfile,'Blanket lithium mass (kg)','(whtblli)',whtblli)
    call ovarre(outfile,'Total mass of FLiBe (kg)','(mflibe)',mflibe)
    call ovarre(outfile,'Shield mass (kg)','(whtshld)',whtshld)

  end subroutine ifefbs