constraint_eqn_024 Subroutine

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

Equation for beta upper limit 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 beta upper limit #=# physics #=#=# fbetatry, betalim and hence also optional here. Logic change during pre-factoring: err, symbol, units will be assigned only if present. iculbl : input integer : switch for beta limit scaling (constraint equation 24):

  • = 0 apply limit to total beta;
  • = 1 apply limit to thermal beta;
  • = 2 apply limit to thermal + neutral beam beta
  • = 3 apply limit to toroidal beta
istell : input integer : switch for stellarator option (set via device.dat):
  • = 0 use tokamak model;
  • = 1 use stellarator model
fbetatry : input real : f-value for beta limit betalim : input real : allowable beta beta : input real : total plasma beta (calculated if ipedestal =3) betaft : input real : fast alpha beta component beta_beam : input real : neutral beam beta component bt : input real : toroidal field btot : input real : total field

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_024(tmp_cc, tmp_con, tmp_err, tmp_symbol, tmp_units)
      !! Equation for beta upper limit
      !! 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 beta upper limit
      !! #=# physics
      !! #=#=# fbetatry, betalim
      !! and hence also optional here.
      !! Logic change during pre-factoring: err, symbol, units will be assigned only if present.
      !! iculbl : input integer : switch for beta limit scaling (constraint equation  24):<UL>
      !! <LI> = 0 apply limit to total beta;
      !! <LI> = 1 apply limit to thermal beta;
      !! <LI> = 2 apply limit to thermal + neutral beam beta
      !! <LI> = 3 apply limit to toroidal beta </UL>
      !! istell : input integer : switch for stellarator option (set via <CODE>device.dat</CODE>):<UL>
      !! <LI> = 0 use tokamak model;
      !! <LI> = 1 use stellarator model</UL>
      !! fbetatry : input real : f-value for beta limit
      !! betalim : input real : allowable beta
      !! beta : input real : total plasma beta (calculated if ipedestal =3)
      !! betaft : input real : fast alpha beta component
      !! beta_beam : input real : neutral beam beta component
      !! bt : input real : toroidal field
      !! btot : input real : total field
      use physics_variables, only: iculbl, betalim, beta, beta_beam, betaft, bt, btot
      use stellarator_variables, only: istell
      use constraint_variables, only: fbetatry
      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

      ! Include all beta components: relevant for both tokamaks and stellarators
      if ((iculbl == 0).or.(istell /= 0)) then
         tmp_cc =  1.0D0 - fbetatry * betalim/beta
         tmp_con = betalim
         tmp_err = betalim - beta / fbetatry
         tmp_symbol = '<'
         tmp_units = ''
      ! Here, the beta limit applies to only the thermal component, not the fast alpha or neutral beam parts
      else if (iculbl == 1) then
         tmp_cc = 1.0D0 - fbetatry * betalim/(beta-betaft-beta_beam)
         tmp_con = betalim
         tmp_err = betalim - (beta-betaft-beta_beam) / fbetatry
         tmp_symbol = '<'
         tmp_units = ''
      ! Beta limit applies to thermal + neutral beam: components of the total beta, i.e. excludes alphas
      else if (iculbl == 2) then
         tmp_cc = 1.0D0 - fbetatry * betalim/(beta-betaft)
         tmp_con = betalim * (1.0D0 - tmp_cc)
         tmp_err = (beta-betaft) * tmp_cc
         tmp_symbol = '<'
         tmp_units = ''
      ! Beta limit applies to toroidal beta
      else if (iculbl == 3) then
         tmp_cc =  1.0D0 - fbetatry * betalim/(beta*(btot/bt)**2)
         tmp_con = betalim
         tmp_err = betalim - (beta*(btot/bt)**2) / fbetatry
         tmp_symbol = '<'
         tmp_units = ''
      end if

   end subroutine constraint_eqn_024