tprofile Function

public function tprofile(rho, rhopedt, t0, tped, tsep, alphat, tbeta)

Implementation of HELIOS-type temperature pedestal profile author: R Kemp, CCFE, Culham Science Centre author: H Lux, CCFE, Culham Science Centre author: P J Knight, CCFE, Culham Science Centre rho : input real : normalised minor radius rhopedt : input real : normalised minor radius pedestal position t0 : input real : central temperature (keV) tped : input real : pedestal temperature (keV) tsep : input real : separatrix temperature (keV) alphat : input real : temperature peaking parameter tbeta : input real : second temperature exponent trho : output real : T(rho) (keV) This routine calculates the temperature at a normalised minor radius position rho for a pedestalised profile.

If ipedestal = 0 the original parabolic profile form is used instead. J.Johner, Fusion Science and Technology 59 (2011), pp 308-349

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in) :: rho
real(kind=dp), intent(in) :: rhopedt
real(kind=dp), intent(in) :: t0
real(kind=dp), intent(in) :: tped
real(kind=dp), intent(in) :: tsep
real(kind=dp), intent(in) :: alphat
real(kind=dp), intent(in) :: tbeta

Return Value real(kind=dp)


Contents

Source Code


Source Code

  function tprofile(rho, rhopedt, t0, tped, tsep, alphat, tbeta)

    !! Implementation of HELIOS-type temperature pedestal profile
    !! author: R Kemp, CCFE, Culham Science Centre
    !! author: H Lux, CCFE, Culham Science Centre
    !! author: P J Knight, CCFE, Culham Science Centre
    !! rho     : input real : normalised minor radius
    !! rhopedt : input real : normalised minor radius pedestal position
    !! t0      : input real : central temperature (keV)
    !! tped    : input real : pedestal temperature (keV)
    !! tsep    : input real : separatrix temperature (keV)
    !! alphat  : input real : temperature peaking parameter
    !! tbeta   : input real : second temperature exponent
    !! trho    : output real : T(rho) (keV)
    !! This routine calculates the temperature at a normalised minor
    !! radius position rho for a pedestalised profile.
    !! <P>If <CODE>ipedestal = 0</CODE> the original parabolic
    !! profile form is used instead.
    !! J.Johner, Fusion Science and Technology 59 (2011), pp 308-349
    !
    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

		use error_handling, only: fdiags, report_error
		use physics_variables, only: ipedestal
    implicit none

    real(dp) :: tprofile

    !  Arguments

    real(dp), intent(in) :: rho, rhopedt, t0, tped, tsep, alphat, tbeta

    !  Local variables

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

    if (ipedestal == 0) then
       tprofile = t0 * (1.0D0 - rho**2)**alphat
       return
    end if

    !  Error trap; shouldn't happen unless volume-averaged temperature has
    !  been allowed to drop below tped. This may happen during a HYBRD case,
    !  but should have been prevented for optimisation runs.

    if (t0 < tped) then
       fdiags(1) = tped ; fdiags(2) = t0
       call report_error(148)
    end if

    if (rho <= rhopedt) then
       tprofile = tped + (t0 - tped) * (1.0D0 - (rho/rhopedT)**tbeta)**alphat
    else
       tprofile = tsep + (tped - tsep) * (1.0D0 - rho)/(1.0D0 - rhopedt)
    end if

  end function tprofile