physics_functions_module Module

Module containing physics subfunctions author: K Ellis, CCFE, Culham Science Centre N/A This module contains physics routines which can be called by physics or other modules !



Contents


Variables

TypeVisibility AttributesNameInitial
real(kind=dp), public :: vcritx

Functions

public function bosch_hale(t, reaction)

Routine to calculate the fusion reaction rate author: R Kemp, CCFE, Culham Science Centre author: P J Knight, CCFE, Culham Science Centre t : input real : Maxwellian density-weighted ion temperature (keV) reaction : input integer : flag for fusion reaction to use: 1 : D-T reaction 2 : D-3He reaction 3 : D-D 1st reaction (50% probability) 4 : D-D 2nd reaction (50% probability) This routine calculates the volumetric fusion reaction rate <sigma v> in m3/s for one of four nuclear reactions, using the Bosch-Hale parametrization.

The valid range of the fit is 0.2 keV < t < 100 keV Bosch and Hale, Nuclear Fusion 32 (1992) 611-631

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in) :: t
integer, intent(in) :: reaction

Return Value real(kind=dp)

public function fsv(u)

Integrand function for the hot beam fusion reaction rate author: P J Knight, CCFE, Culham Science Centre u : input real : abscissa of integration, = ratio of beam velocity to the critical velocity This is the integrand function for the hot beam fusion reaction rate. !

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in) :: u

Return Value real(kind=dp)

public function plasma_elongation_IPB()

Author
H Lux (UKAEA)

Volume measure of plasma elongation using the IPB definition

Read more…

Arguments

None

Return Value real(kind=dp)

Plasma elongation (IPB)

public function total_mag_field()

Author
J. Morris (UKAEA)

Calculates the total magnetic field

Read more…

Arguments

None

Return Value real(kind=dp)

public function beta_poloidal()

Author
J. Morris (UKAEA)

Calculates total poloidal beta

Read more…

Arguments

None

Return Value real(kind=dp)

public function res_diff_time()

Author
J. Morris (UKAEA)

Calculates resistive diffusion time

Read more…

Arguments

None

Return Value real(kind=dp)


Subroutines

public subroutine init_physics_functions()

Initialise module variables

Arguments

None

public subroutine pthresh(dene, dnla, bt, rmajor, kappa, sarea, aion, aspect, pthrmw)

Author
P J Knight, CCFE, Culham Science Centre

L-mode to H-mode power threshold calculation dene : input real : volume-averaged electron density (/m3) dnla : input real : line-averaged electron density (/m3) bt : input real : toroidal field on axis (T) rmajor : input real : plasma major radius (m) kappa : input real : plasma elongation sarea : input real : plasma surface area (m**2) aion : input real : average mass of all ions (amu) aspect : input real : aspect ratio pthrmw(17) : output real array : power threshold (different scalings) This routine calculates the power threshold for the L-mode to H-mode transition. ITER Physics Design Description Document, p.2-2 ITER-FDR Plasma Performance Assessments, p.III-9 Snipes, 24th EPS Conference, Berchtesgaden 1997, p.961 Martin et al, 11th IAEA Tech. Meeting on H-mode Physics and Transport Barriers, Journal of Physics: Conference Series 123 (2008) 012033 J A Snipes and the International H-mode Threshold Database Working Group, 2000, Plasma Phys. Control. Fusion, 42, A299 !

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in) :: dene
real(kind=dp), intent(in) :: dnla
real(kind=dp), intent(in) :: bt
real(kind=dp), intent(in) :: rmajor
real(kind=dp), intent(in) :: kappa
real(kind=dp), intent(in) :: sarea
real(kind=dp), intent(in) :: aion
real(kind=dp), intent(in) :: aspect
real(kind=dp), intent(out), dimension(21):: pthrmw

public subroutine palph(alphan, alphat, deni, fdeut, fhe3, ftrit, ti, palppv, pchargepv, pneutpv, sigvdt, fusionrate, alpharate, protonrate, pdtpv, pdhe3pv, pddpv)

(Initial part of) fusion power and fast alpha pressure calculations author: P J Knight, CCFE, Culham Science Centre alphan : input real : density profile index alphat : input real : temperature profile index deni : input real : fuel ion density (/m3) fdeut : input real : deuterium fuel fraction fhe3 : input real : helium-3 fuel fraction ftrit : input real : tritium fuel fraction ti : input real : ion temperature (keV) palppv : output real : alpha particle fusion power per volume (MW/m3) pchargepv : output real : other charged particle fusion power/volume (MW/m3) pneutpv : output real : neutron fusion power per volume (MW/m3) sigvdt : output real : profile averaged (m3/s) fusionrate : output real : fusion reaction rate (reactions/m3/s) alpharate : output real : alpha particle production rate (/m3/s) protonrate : output real : proton production rate (/m3/s) pdtpv : output real : D-T fusion power (MW/m3) pdhe3pv : output real : D-He3 fusion power (MW/m3) pddpv : output real : D-D fusion power (MW/m3) This subroutine numerically integrates over plasma cross-section to find the core plasma fusion power. T&M/PKNIGHT/LOGBOOK24, p.6

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in) :: alphan
real(kind=dp), intent(in) :: alphat
real(kind=dp), intent(in) :: deni
real(kind=dp), intent(in) :: fdeut
real(kind=dp), intent(in) :: fhe3
real(kind=dp), intent(in) :: ftrit
real(kind=dp), intent(in) :: ti
real(kind=dp), intent(out) :: palppv
real(kind=dp), intent(out) :: pchargepv
real(kind=dp), intent(out) :: pneutpv
real(kind=dp), intent(out) :: sigvdt
real(kind=dp), intent(out) :: fusionrate
real(kind=dp), intent(out) :: alpharate
real(kind=dp), intent(out) :: protonrate
real(kind=dp), intent(out) :: pdtpv
real(kind=dp), intent(out) :: pdhe3pv
real(kind=dp), intent(out) :: pddpv

public subroutine palph2(bt, bp, dene, deni, dnitot, falpe, falpi, palpnb, ifalphap, pchargepv, pneutpv, ten, tin, vol, palpmw, pneutmw, pchargemw, betaft, palppv, palpipv, palpepv, pfuscmw, powfmw)

(Concluding part of) fusion power and fast alpha pressure calculations author: P J Knight, CCFE, Culham Science Centre bp : input real : poloidal field (T) bt : input real : toroidal field on axis (T) dene : input real : electron density (/m3) deni : input real : fuel ion density (/m3) dnitot : input real : total ion density (/m3) falpe : input real : fraction of alpha energy to electrons falpi : input real : fraction of alpha energy to ions ifalphap : input integer : switch for fast alpha pressure method palpnb : input real : alpha power from hot neutral beam ions (MW) pchargepv : input real : other charged particle fusion power/volume (MW/m3) pneutpv : input/output real : neutron fusion power per volume (MW/m3) ten : input real : density-weighted electron temperature (keV) tin : input real : density-weighted ion temperature (keV) vol : input real : plasma volume (m3) palpmw : output real : alpha power (MW) pneutmw : output real : neutron fusion power (MW) pchargemw : output real : other charged particle fusion power (MW) betaft : output real : fast alpha beta component palppv : input/output real : alpha power per volume (MW/m3) palpepv : output real : alpha power per volume to electrons (MW/m3) palpipv : output real : alpha power per volume to ions (MW/m3) pfuscmw : output real : charged particle fusion power (MW) powfmw : output real : fusion power (MW) This subroutine completes the calculation of the fusion power fast alpha pressure, and determines other alpha particle quantities. ITER Physics Design Guidelines: 1989 [IPDG89], N. A. Uckan et al, ITER Documentation Series No.10, IAEA/ITER/DS/10, IAEA, Vienna, 1990 D J Ward, UKAEA Fusion: F/PL/PJK/PROCESS/CODE/050

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in) :: bt
real(kind=dp), intent(in) :: bp
real(kind=dp), intent(in) :: dene
real(kind=dp), intent(in) :: deni
real(kind=dp), intent(in) :: dnitot
real(kind=dp), intent(in) :: falpe
real(kind=dp), intent(in) :: falpi
real(kind=dp), intent(in) :: palpnb
integer, intent(in) :: ifalphap
real(kind=dp), intent(in) :: pchargepv
real(kind=dp), intent(inout) :: pneutpv
real(kind=dp), intent(in) :: ten
real(kind=dp), intent(in) :: tin
real(kind=dp), intent(in) :: vol
real(kind=dp), intent(out) :: palpmw
real(kind=dp), intent(out) :: pneutmw
real(kind=dp), intent(out) :: pchargemw
real(kind=dp), intent(out) :: betaft
real(kind=dp), intent(inout) :: palppv
real(kind=dp), intent(out) :: palpipv
real(kind=dp), intent(out) :: palpepv
real(kind=dp), intent(out) :: pfuscmw
real(kind=dp), intent(out) :: powfmw

public subroutine beamfus(beamfus0, betbm0, bp, bt, cnbeam, dene, deni, dlamie, ealphadt, enbeam, fdeut, ftrit, ftritbm, sigvdt, ten, tin, vol, zeffai, betanb, dnbeam2, palpnb)

Routine to calculate beam slowing down properties author: P J Knight, CCFE, Culham Science Centre beamfus0: input real : multiplier for beam-background fusion calculation betbm0 : input real : leading coefficient for neutral beam beta fraction bp : input real : poloidal field (T) bt : input real : toroidal field on axis (T) cnbeam : input real : neutral beam current (A) dene : input real : electron density (/m3) deni : input real : fuel ion density (/m3) dlamie : input real : ion-electron coulomb logarithm ealphadt : input real : alpha particle birth energy (D-T) (keV) enbeam : input real : neutral beam energy (keV) fdeut : input real : deuterium fraction of main plasma ftrit : input real : tritium fraction of main plasma ftritbm: input real : tritium fraction of neutral beam sigvdt : input real : profile averaged for D-T (m3/s) ten : input real : density weighted average electron temperature (keV) tin : input real : density weighted average ion temperature (keV) vol : input real : plasma volume (m3) zeffai : input real : mass weighted plasma effective charge betanb : output real : neutral beam beta component dnbeam2: output real : hot beam ion density (/m3) palpnb : output real : alpha power from hot neutral beam ions (MW) This routine calculates the beam slowing down properties. !

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in) :: beamfus0
real(kind=dp), intent(in) :: betbm0
real(kind=dp), intent(in) :: bp
real(kind=dp), intent(in) :: bt
real(kind=dp), intent(in) :: cnbeam
real(kind=dp), intent(in) :: dene
real(kind=dp), intent(in) :: deni
real(kind=dp), intent(in) :: dlamie
real(kind=dp), intent(in) :: ealphadt
real(kind=dp), intent(in) :: enbeam
real(kind=dp), intent(in) :: fdeut
real(kind=dp), intent(in) :: ftrit
real(kind=dp), intent(in) :: ftritbm
real(kind=dp), intent(in) :: sigvdt
real(kind=dp), intent(in) :: ten
real(kind=dp), intent(in) :: tin
real(kind=dp), intent(in) :: vol
real(kind=dp), intent(in) :: zeffai
real(kind=dp), intent(out) :: betanb
real(kind=dp), intent(out) :: dnbeam2
real(kind=dp), intent(out) :: palpnb

public subroutine beamcalc(nd, nt, ealphadt, ebeam, ecritd, ecritt, tausbme, ftritbm, ibeam, ti, vol, svdt, palfdb, palftb, nhot, ehot)

Neutral beam alpha power and ion energy author: P J Knight, CCFE, Culham Science Centre ealphadt : input real : alpha particle birth energy (D-T) (keV) ebeam : input real : beam energy (keV) ecritd : input real : critical energy for electron/ion slowing down of the beam ion (deuterium neutral beam) (keV) ecritt : input real : critical energy for beam slowing down (tritium neutral beam) (keV) ftritbm: input real : beam tritium fraction (0.0 = deuterium beam) ibeam : input real : beam current (A) nd : input real : thermal deuterium density (/m3) nt : input real : thermal tritium density (/m3) svdt : input real : profile averaged for D-T (m3/s) tausbme: input real : beam ion slowing down time on electrons (s) ti : input real : thermal ion temperature (keV) vol : input real : plasma volume (m3) (95% flux surface) ehot : output real : average hot beam ion energy (keV) nhot : output real : hot beam ion density (/m3) palfdb : output real : alpha power from deut. beam-background fusion (MW) palftb : output real : alpha power from trit. beam-background fusion (MW) This routine calculates the neutral beam alpha power and ion energy. !

Arguments

Type IntentOptional AttributesName
real(kind=kind(1.0D0)), intent(in) :: nd
real(kind=kind(1.0D0)), intent(in) :: nt
real(kind=kind(1.0D0)), intent(in) :: ealphadt
real(kind=kind(1.0D0)), intent(in) :: ebeam
real(kind=kind(1.0D0)), intent(in) :: ecritd
real(kind=kind(1.0D0)), intent(in) :: ecritt
real(kind=kind(1.0D0)), intent(in) :: tausbme
real(kind=kind(1.0D0)), intent(in) :: ftritbm
real(kind=kind(1.0D0)), intent(in) :: ibeam
real(kind=kind(1.0D0)), intent(in) :: ti
real(kind=kind(1.0D0)), intent(in) :: vol
real(kind=kind(1.0D0)), intent(in) :: svdt
real(kind=kind(1.0D0)), intent(out) :: palfdb
real(kind=kind(1.0D0)), intent(out) :: palftb
real(kind=kind(1.0D0)), intent(out) :: nhot
real(kind=kind(1.0D0)), intent(out) :: ehot