constraint_eqn_084 Subroutine

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

Equation for the lower limit of beta author: J Lion, IPP Greifswald args : output structure : residual error; constraint value; residual error in physical units; output string; units string (beta-betaft) > betalim_lower #=# physics #=#=# betaft, beta, fbetatry_lower Logic change during pre-factoring: err, symbol, units will be assigned only if present. fbetatry_lower : input real : f-value for constraint beta-betaft > betalim_lower betalim_lower : input real : Lower limit for beta beta : input real : plasma beta betaft : input real : Alpha particle beta

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_084(tmp_cc, tmp_con, tmp_err, tmp_symbol, tmp_units)
      !! Equation for the lower limit of beta
      !! author: J Lion, IPP Greifswald
      !! args : output structure : residual error; constraint value;
      !! residual error in physical units; output string; units string
      !!  (beta-betaft) > betalim_lower
      !! #=# physics
      !! #=#=# betaft, beta, fbetatry_lower
      !! Logic change during pre-factoring: err, symbol, units will be assigned only if present.
      !! fbetatry_lower : input real : f-value for constraint beta-betaft > betalim_lower
      !! betalim_lower : input real :  Lower limit for beta
      !! beta : input real :  plasma beta
      !! betaft : input real : Alpha particle beta

      use physics_variables, only: betalim_lower, beta, betaft
      use constraint_variables, only: fbetatry_lower
      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


      tmp_cc = 1.0D0 - fbetatry_lower * (beta-betaft)/betalim_lower
      tmp_con = betalim_lower * (1.0D0 - tmp_cc)
      tmp_err = (beta-betaft) * tmp_cc
      tmp_symbol = '>'
      tmp_units = ''


   end subroutine constraint_eqn_084