constraint_eqn_045 Subroutine

public subroutine constraint_eqn_045(tmp_cc, tmp_con, tmp_err, tmp_symbol, tmp_units)

Equation for edge safety factor lower limit (TART) author: P B Lloyd, CCFE, Culham Science Centre args : output structure : residual error; constraint value; residual error in physical units; output string; units string Equation for edge safety factor lower limit (TART) #=# tfcoil #=#=# fq, qlim and hence also optional here. Logic change during pre-factoring: err, symbol, units will be assigned only if present. fq : input real : f-value for edge safety factor q : safety factor 'near' plasma edge: equal to q95 (unless i_plasma_current = 2 (ST current scaling), in which case q = mean edge safety factor qbar) qlim : input real : lower limit for edge safety factor itart : input integer : switch for spherical tokamak (ST) models:

  • = 0 use conventional aspect ratio models;
  • = 1 use spherical tokamak models

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(out) :: tmp_cc
real(kind=dp), intent(out) :: tmp_con
real(kind=dp), intent(out) :: tmp_err
character(len=1), intent(out) :: tmp_symbol
character(len=10), intent(out) :: tmp_units

Contents

Source Code


Source Code

   subroutine constraint_eqn_045(tmp_cc, tmp_con, tmp_err, tmp_symbol, tmp_units)
      !! Equation for edge safety factor lower limit (TART)
      !! author: P B Lloyd, CCFE, Culham Science Centre
      !! args : output structure : residual error; constraint value;
      !! residual error in physical units; output string; units string
      !! Equation for edge safety factor lower limit (TART)
      !! #=# tfcoil
      !! #=#=# fq, qlim
      !! and hence also optional here.
      !! Logic change during pre-factoring: err, symbol, units will be assigned only if present.
      !! fq : input real : f-value for edge safety factor
      !! q : safety factor 'near' plasma edge: equal to q95
      !! (unless i_plasma_current = 2 (ST current scaling), in which case q = mean edge safety factor qbar)
      !! qlim : input real :  lower limit for edge safety factor
      !! itart : input integer : switch for spherical tokamak (ST) models:<UL>
      !! <LI> = 0 use conventional aspect ratio models;
      !! <LI> = 1 use spherical tokamak models</UL>
      use constraint_variables, only: fq
      use physics_variables, only: q, qlim, itart
      implicit none
            real(dp), intent(out) :: tmp_cc
      real(dp), intent(out) :: tmp_con
      real(dp), intent(out) :: tmp_err
      character(len=1), intent(out) :: tmp_symbol
      character(len=10), intent(out) :: tmp_units

      ! if the machine isn't a ST then report error
      if (itart == 0) call report_error(9)
      tmp_cc =   1.0D0 - fq * q/qlim
      tmp_con = qlim * (1.0D0 - tmp_cc)
      tmp_err = qlim * tmp_cc
      tmp_symbol = '<'
      tmp_units = ''

   end subroutine constraint_eqn_045