libfida Module

Main FIDASIM library


Uses

  • module~~libfida~~UsesGraph module~libfida libfida HDF5 HDF5 module~libfida->HDF5 H5LT H5LT module~libfida->H5LT module~utilities utilities module~libfida->module~utilities module~hdf5_extra hdf5_extra module~libfida->module~hdf5_extra module~eigensystem eigensystem module~libfida->module~eigensystem omp_lib omp_lib module~utilities->omp_lib module~hdf5_extra->HDF5 module~hdf5_extra->H5LT

Used by

  • module~~libfida~~UsedByGraph module~libfida libfida program~fidasim fidasim program~fidasim->module~libfida

Contents


Variables

TypeVisibility AttributesNameInitial
character(len=30), public :: version =''

FIDASIM version number

integer, public, parameter:: charlim =150

Defines character limit for files and directories

character(len=charlim), public :: namelist_file

Input namelist file

integer, public, parameter:: nbif_type =1

Identifier for full energy NBI neutral interaction

integer, public, parameter:: nbih_type =2

Identifier for half energy NBI neutral interaction

integer, public, parameter:: nbit_type =3

Identifier for third energy NBI neutral interaction

integer, public, parameter:: halo_type =4

Identifier for halo neutral interaction

integer, public, parameter:: fida_type =5

Identifier for fida neutral interaction

integer, public, parameter:: brems_type =6

Identifier for bremsstrahlung interaction. Acts as dummy type

integer, public, parameter:: ntypes =6

Number of different types of neutrals

integer, public, parameter:: beam_ion =1

Identifier for a beam ion

integer, public, parameter:: thermal_ion =2

Identifier for a thermal ion

real(kind=Float64), public, parameter:: e_amu =5.485799093287202d-4

Atomic mass of an electron [amu]

real(kind=Float64), public, parameter:: H_1_amu =1.00782504d0

Atomic mass of Hydrogen-1 [amu]

real(kind=Float64), public, parameter:: H_2_amu =2.0141017778d0

Atomic mass of Hydrogen-2 [amu]

real(kind=Float64), public, parameter:: B5_amu =10.81d0

Atomic mass of Boron [amu]

real(kind=Float64), public, parameter:: C6_amu =12.011d0

Atomic mass of Carbon [amu]

real(kind=Float64), public, parameter:: mass_u =1.6605402d-27

Atomic mass unit [kg]

real(kind=Float64), public, parameter:: e0 =1.60217733d-19

Electron charge [C]

real(kind=Float64), public, parameter:: pi =3.14159265358979323846264d0

Pi

real(kind=Float64), public, parameter:: c0 =2.99792458d+08

Speed of light [m/s]

real(kind=Float64), public, parameter:: h_planck =4.135667516d-15

Planck's constant [eV*s]

real(kind=Float64), public, parameter:: lambda0 =656.1d0

D-alpha emission line [nm]

real(kind=Float64), public, parameter:: v2_to_E_per_amu =mass_u/(2.*e0*1.d3)*1.d-4

to keV conversion factor

integer, public, parameter:: n_stark =15

Number of Stark lines

real(kind=Float64), public, parameter, dimension(n_stark):: stark_wavel =[-2.20200d-07, -1.65200d-07, -1.37700d-07, -1.10200d-07, -8.26400d-08, -5.51000d-08, -2.75600d-08, 0.00000d0, 2.75700d-08, 5.51500d-08, 8.27400d-08, 1.10300d-07, 1.38000d-07, 1.65600d-07, 2.20900d-07]

Stark wavelengths [nm*m/V]

real(kind=Float64), public, parameter, dimension(n_stark):: stark_intens =[1.000d0, 18.00d0, 16.00d0, 1681.d0, 2304.d0, 729.0d0, 1936.d0, 5490.d0, 1936.d0, 729.0d0, 2304.d0, 1681.d0, 16.00d0, 18.00d0, 1.000d0]

Stark Intensities

integer, public, parameter, dimension(n_stark):: stark_pi =[1, 0, 0, 1, 1, 1, 0, 0, 0, 1, 1, 1, 0, 0, 1]

Pi line indicators

integer, public, parameter, dimension(n_stark):: stark_sigma =1-stark_pi

Sigma line indicators

integer, public, parameter:: nlevs =6

Number of atomic energy levels

real(kind=Float64), public :: colrad_threshold =1.d6

colrad threshold

real(kind=Float64), public, dimension(ntypes):: halo_iter_dens =0.d0

Keeps track of how of each generations halo density

integer, public :: nbi_outside =0

Keeps track of how many beam neutrals do not hit the beam_grid

type(BeamGrid), public, save:: beam_grid

Variable containing beam grid definition

type(InterpolationGrid), public, save:: inter_grid

Variable containing interpolation grid definition

type(FastIonDistribution), public, save:: fbm

Variable containing the fast-ion distribution function

type(FastIonParticles), public, save:: particles

Variable containing a MC fast-ion distribution

type(Equilibrium), public, save:: equil

Variable containing the plasma parameters and fields

type(NeutralBeam), public, save:: nbi

Variable containing the neutral beam geometry and settings

type(AtomicTables), public, save:: tables

Variable containing the atomic tables

type(NPAResults), public, save:: npa

Variable for storing the calculated NPA results

type(SpectralChords), public, save:: spec_chords

Variable containing the spectral system definition

type(NPAChords), public, save:: npa_chords

Variable containing the NPA system definition

type(SimulationInputs), public, save:: inputs

Variable containing the simulation inputs

type(BirthProfile), public, save:: birth

Variable for storing the calculated birth profile

type(NeutralDensity), public, save:: neut

Variable for storing the calculated beam density

type(Spectra), public, save:: spec

Variable for storing the calculated spectra

type(NeutronRate), public, save:: neutron

Variable for storing the neutron rate

type(FIDAWeights), public, save:: fweight

Variable for storing the calculated FIDA weights

type(NPAWeights), public, save:: nweight

Variable for storing the calculated NPA weights


Interfaces

public interface assignment(=)

public interface operator(+)

Allows for adding Profiles,LocalProfiles, EMFields, and LocalEMFields

public interface operator(-)

Allows for subtracting Profiles,LocalProfiles, EMFields, and LocalEMFields

public interface operator(*)

Allows for multiplying Profiles,LocalProfiles, EMFields, and LocalEMFields by scalars

  • public function sp_multiply(real_scalar, p1) result(p3)

    Defines how to multiply Profiles types by a scalar

    Arguments

    Type IntentOptional AttributesName
    real(kind=Float64), intent(in) :: real_scalar
    type(Profiles), intent(in) :: p1

    Return Value type(Profiles)

  • public function ps_multiply(p1, real_scalar) result(p3)

    Defines how to multiply Profiles types by a scalar

    Arguments

    Type IntentOptional AttributesName
    type(Profiles), intent(in) :: p1
    real(kind=Float64), intent(in) :: real_scalar

    Return Value type(Profiles)

  • public function lps_multiply(p1, real_scalar) result(p3)

    Defines how to multiply LocalProfiles types by a scalar

    Arguments

    Type IntentOptional AttributesName
    type(LocalProfiles), intent(in) :: p1
    real(kind=Float64), intent(in) :: real_scalar

    Return Value type(LocalProfiles)

  • public function slp_multiply(real_scalar, p1) result(p3)

    Defines how to multiply LocalProfiles types by a scalar

    Arguments

    Type IntentOptional AttributesName
    real(kind=Float64), intent(in) :: real_scalar
    type(LocalProfiles), intent(in) :: p1

    Return Value type(LocalProfiles)

  • public function sf_multiply(real_scalar, p1) result(p3)

    Defines how to multiply EMFields types by a scalar

    Arguments

    Type IntentOptional AttributesName
    real(kind=Float64), intent(in) :: real_scalar
    type(EMFields), intent(in) :: p1

    Return Value type(EMFields)

  • public function fs_multiply(p1, real_scalar) result(p3)

    Defines how to multiply EMFields types by a scalar

    Arguments

    Type IntentOptional AttributesName
    type(EMFields), intent(in) :: p1
    real(kind=Float64), intent(in) :: real_scalar

    Return Value type(EMFields)

  • public function lfs_multiply(p1, real_scalar) result(p3)

    Defines how to multiply LocalEMFields types by a scalar

    Arguments

    Type IntentOptional AttributesName
    type(LocalEMFields), intent(in) :: p1
    real(kind=Float64), intent(in) :: real_scalar

    Return Value type(LocalEMFields)

  • public function slf_multiply(real_scalar, p1) result(p3)

    Defines how to multiply LocalEMFields types by a scalar

    Arguments

    Type IntentOptional AttributesName
    real(kind=Float64), intent(in) :: real_scalar
    type(LocalEMFields), intent(in) :: p1

    Return Value type(LocalEMFields)

public interface operator(/)

Allows for dividing Profiles,LocalProfiles, EMFields, and LocalEMFields by scalars

  • public function ps_divide(p1, real_scalar) result(p3)

    Defines how to divide Profiles types by a scalar

    Arguments

    Type IntentOptional AttributesName
    type(Profiles), intent(in) :: p1
    real(kind=Float64), intent(in) :: real_scalar

    Return Value type(Profiles)

  • public function lps_divide(p1, real_scalar) result(p3)

    Defines how to divide LocalProfiles types by a scalar

    Arguments

    Type IntentOptional AttributesName
    type(LocalProfiles), intent(in) :: p1
    real(kind=Float64), intent(in) :: real_scalar

    Return Value type(LocalProfiles)

  • public function fs_divide(p1, real_scalar) result(p3)

    Defines how to divide EMFields types by a scalar

    Arguments

    Type IntentOptional AttributesName
    type(EMFields), intent(in) :: p1
    real(kind=Float64), intent(in) :: real_scalar

    Return Value type(EMFields)

  • public function lfs_divide(p1, real_scalar) result(p3)

    Defines how to divide LocalEMFields types by a scalar

    Arguments

    Type IntentOptional AttributesName
    type(LocalEMFields), intent(in) :: p1
    real(kind=Float64), intent(in) :: real_scalar

    Return Value type(LocalEMFields)

public interface interpol_coeff

Calculates linear interpolation coefficients

  • public subroutine interpol1D_coeff(xmin, dx, nx, xout, c, err)

    Linear interpolation coefficients and index for a 1D grid y(x)

    Arguments

    Type IntentOptional AttributesName
    real(kind=Float64), intent(in) :: xmin

    Minimum abscissa value

    real(kind=Float64), intent(in) :: dx

    Absissa spacing

    integer, intent(in) :: nx

    Number of abscissa

    real(kind=Float64), intent(in) :: xout

    Abscissa value to interpolate

    type(InterpolCoeffs1D), intent(out) :: c

    Interpolation Coefficients

    integer, intent(out), optional :: err

    Error code

  • public subroutine interpol1D_coeff_arr(x, xout, c, err)

    Linear interpolation coefficients and index for a 1D grid y(x)

    Arguments

    Type IntentOptional AttributesName
    real(kind=Float64), intent(in), dimension(:):: x

    Abscissa values

    real(kind=Float64), intent(in) :: xout

    Abscissa value to interpolate

    type(InterpolCoeffs1D), intent(out) :: c

    Interpolation Coefficients

    integer, intent(out), optional :: err

    Error code

  • public subroutine interpol2D_coeff(xmin, dx, nx, ymin, dy, ny, xout, yout, c, err)

    Bilinear interpolation coefficients and indicies for a 2D grid z(x,y)

    Arguments

    Type IntentOptional AttributesName
    real(kind=Float64), intent(in) :: xmin

    Minimum abscissa

    real(kind=Float64), intent(in) :: dx

    Abscissa spacing

    integer, intent(in) :: nx

    Number of abscissa

    real(kind=Float64), intent(in) :: ymin

    Minimum ordinate

    real(kind=Float64), intent(in) :: dy

    Ordinate spacing

    integer, intent(in) :: ny

    Number of ordinates points

    real(kind=Float64), intent(in) :: xout

    Abscissa value to interpolate

    real(kind=Float64), intent(in) :: yout

    Ordinate value to interpolate

    type(InterpolCoeffs2D), intent(out) :: c

    Interpolation Coefficients

    integer, intent(out), optional :: err

    Error code

  • public subroutine interpol2D_coeff_arr(x, y, xout, yout, c, err)

    Arguments

    Type IntentOptional AttributesName
    real(kind=Float64), intent(in), dimension(:):: x

    Abscissa values

    real(kind=Float64), intent(in), dimension(:):: y

    Ordinate values

    real(kind=Float64), intent(in) :: xout

    Abscissa value to interpolate

    real(kind=Float64), intent(in) :: yout

    Ordinate value to interpolate

    type(InterpolCoeffs2D), intent(out) :: c

    Interpolation Coefficients

    integer, intent(out), optional :: err

    Error code

public interface interpol

Performs linear/bilinear interpolation

  • public subroutine interpol1D_arr(x, y, xout, yout, err, coeffs)

    Performs linear interpolation on a uniform 1D grid y(x)

    Arguments

    Type IntentOptional AttributesName
    real(kind=Float64), intent(in), dimension(:):: x

    The abscissa values of y

    real(kind=Float64), intent(in), dimension(:):: y

    Values at abscissa values x: y(x)

    real(kind=Float64), intent(in) :: xout

    Abscissa value to interpolate

    real(kind=Float64), intent(out) :: yout

    Interpolant: y(xout)

    integer, intent(out), optional :: err

    Error code

    type(InterpolCoeffs1D), intent(in), optional :: coeffs

    Precomputed Linear Interpolation Coefficients

  • public subroutine interpol2D_arr(x, y, z, xout, yout, zout, err, coeffs)

    Performs bilinear interpolation on a 2D grid z(x,y)

    Arguments

    Type IntentOptional AttributesName
    real(kind=Float64), intent(in), dimension(:):: x

    The abscissa values of z

    real(kind=Float64), intent(in), dimension(:):: y

    The ordinate values of z

    real(kind=Float64), intent(in), dimension(:,:):: z

    Values at the abscissa/ordinates: z(x,y)

    real(kind=Float64), intent(in) :: xout

    The abscissa value to interpolate

    real(kind=Float64), intent(in) :: yout

    The ordinate value to interpolate

    real(kind=Float64), intent(out) :: zout

    Interpolant: z(xout,yout)

    integer, intent(out), optional :: err

    Error code

    type(InterpolCoeffs2D), intent(in), optional :: coeffs

    Precomputed Linear Interpolation Coefficients

  • public subroutine interpol2D_2D_arr(x, y, z, xout, yout, zout, err, coeffs)

    Performs bilinear interpolation on a 2D grid of 2D arrays z(:,:,x,y)

    Arguments

    Type IntentOptional AttributesName
    real(kind=Float64), intent(in), dimension(:):: x

    The abscissa values of z

    real(kind=Float64), intent(in), dimension(:):: y

    The ordinate values of z

    real(kind=Float64), intent(in), dimension(:,:,:,:):: z

    Values at the abscissa/ordinates: z(:,:,x,y)

    real(kind=Float64), intent(in) :: xout

    The abscissa value to interpolate

    real(kind=Float64), intent(in) :: yout

    The ordinate value to interpolate

    real(kind=Float64), intent(out), dimension(:,:):: zout

    Interpolant: z(:,:,xout,yout)

    integer, intent(out), optional :: err

    Error code

    type(InterpolCoeffs2D), intent(in), optional :: coeffs

    Precomputed Linear Interpolation Coefficients

public interface store_neutrals

  • public subroutine store_neutrals_cell(ind, neut_type, dens, store_iter)

    Arguments

    Type IntentOptional AttributesName
    integer(kind=Int32), intent(in), dimension(3):: ind

    beam_grid indices

    integer, value:: neut_type

    Neutral type

    real(kind=Float64), intent(in), dimension(:):: dens

    Neutral density [neutrals/cm^3]

    logical, intent(in), optional :: store_iter

    Store DCX/Halo iteration density in halo_iter_dens

  • public subroutine store_neutrals_track(tracks, ncell, neut_type)

    Arguments

    Type IntentOptional AttributesName
    type(ParticleTrack), intent(in), dimension(:):: tracks

    Neutral Particle Track

    integer, intent(in) :: ncell

    Number of cell in the particle track

    integer, value:: neut_type

    Neutral type


Derived Types

type, public :: InterpolCoeffs1D

Linear Interpolation Coefficients and indices

Components

TypeVisibility AttributesNameInitial
integer, public :: i =0

Index of position right before xout

real(kind=Float64), public :: b1 =0.d0

Coefficient for y(i) term

real(kind=Float64), public :: b2 =0.d0

Coefficient for y(i+1) term

type, public :: InterpolCoeffs2D

2D Linear Interpolation Coefficients and indices

Components

TypeVisibility AttributesNameInitial
integer, public :: i =0

Index of abscissa before xout

integer, public :: j =0

Index of ordinate before yout

real(kind=Float64), public :: b11 =0.d0

Coefficient for z(i,j) term

real(kind=Float64), public :: b12 =0.d0

Coefficient for z(i,j+1) term

real(kind=Float64), public :: b21 =0.d0

Coefficient for z(i+1,j) term

real(kind=Float64), public :: b22 =0.d0

Coefficient for z(i+1,j+1) term

type, public :: BeamGrid

Defines a 3D grid for neutral beam calculations

Components

TypeVisibility AttributesNameInitial
integer(kind=Int32), public :: nx

Number of cells in the x direction

integer(kind=Int32), public :: ny

Number of cells in the y direction

integer(kind=Int32), public :: nz

Number of cells in the z direction

real(kind=Float64), public :: xmin

Minimum x value

real(kind=Float64), public :: xmax

Maximum x value

real(kind=Float64), public :: ymin

Minimum y value

real(kind=Float64), public :: ymax

Maximum y value

real(kind=Float64), public :: zmin

Minimum z value

real(kind=Float64), public :: zmax

Maximum z value

real(kind=Float64), public :: alpha

Tait-Bryan angle for a rotation about z [radians]

real(kind=Float64), public :: beta

Tait-Bryan angle for a rotation about y' [radians]

real(kind=Float64), public :: gamma

Tait-Bryan angle for a rotation about x" [radians]

real(kind=Float64), public :: drmin

Minimum cell spacing: min(dx,dy,dz)

real(kind=Float64), public :: dv

Cell volume []

real(kind=Float64), public :: volume

Grid volume []

integer(kind=Int32), public :: ntrack

Maximum number of cell for particle tracking

integer(kind=Int32), public :: ngrid

Number of cells

real(kind=Float64), public, dimension(3):: origin

Origin of beam grid in machine coordinates

real(kind=Float64), public, dimension(3):: center

Center of beam grid in beam coordinates

real(kind=Float64), public, dimension(3):: dr

Cell spacings [dx, dy, dz]

real(kind=Float64), public, dimension(3):: lwh

Grid [length(x), width(y), height(z)]

real(kind=Float64), public, dimension(3,3):: basis

Beam grid basis for converting from beam coordinates(xyz) to machine coordinates(uvw): (\uvw = B*xyz + origin)

real(kind=Float64), public, dimension(3,3):: inv_basis

Inverse basis for reverse transformation: (\xyz = B^{-1}*(uvw - origin))

real(kind=Float64), public, dimension(:), allocatable:: xc

x positions of cell centers

real(kind=Float64), public, dimension(:), allocatable:: yc

y positions of cell centers

real(kind=Float64), public, dimension(:), allocatable:: zc

z positions of cell centers

type, public :: InterpolationGrid

Defines a 2D R-Z grid for interpolating plasma parameters and fields

Components

TypeVisibility AttributesNameInitial
integer(kind=Int32), public :: nr

Number of Radii

integer(kind=Int32), public :: nz

Number of Z values

real(kind=Float64), public :: dr

Radial spacing [cm]

real(kind=Float64), public :: dz

Vertical spacing [cm]

real(kind=Float64), public :: da

Grid element area []

real(kind=Float64), public, dimension(:), allocatable:: r

Radii values [cm]

real(kind=Float64), public, dimension(:), allocatable:: z

Z values [cm]

real(kind=Float64), public, dimension(:,:), allocatable:: r2d

2D R grid [cm]

real(kind=Float64), public, dimension(:,:), allocatable:: z2d

2D Z grid [cm]

type, public :: Profiles

Torodial symmetric plasma parameters at a given R-Z

Components

TypeVisibility AttributesNameInitial
real(kind=Float64), public :: dene =0.d0

Electron density []

real(kind=Float64), public :: denp =0.d0

Ion density []

real(kind=Float64), public :: denimp =0.d0

Impurity density []

real(kind=Float64), public :: denf =0.d0

Fast-ion density []

real(kind=Float64), public :: te =0.d0

Electron temperature [kev]

real(kind=Float64), public :: ti =0.d0

Ion temperature [kev]

real(kind=Float64), public :: zeff =0.d0

Effective Nuclear Charge

real(kind=Float64), public :: vr =0.d0

Plasma rotation in radial direction

real(kind=Float64), public :: vt =0.d0

Plasma rotation in torodial/phi direction

real(kind=Float64), public :: vz =0.d0

Plasma rotation in z direction

type, public, extends(Profiles) :: LocalProfiles

Plasma parameters at given position

Components

TypeVisibility AttributesNameInitial
real(kind=Float64), public :: dene =0.d0

Electron density []

real(kind=Float64), public :: denp =0.d0

Ion density []

real(kind=Float64), public :: denimp =0.d0

Impurity density []

real(kind=Float64), public :: denf =0.d0

Fast-ion density []

real(kind=Float64), public :: te =0.d0

Electron temperature [kev]

real(kind=Float64), public :: ti =0.d0

Ion temperature [kev]

real(kind=Float64), public :: zeff =0.d0

Effective Nuclear Charge

real(kind=Float64), public :: vr =0.d0

Plasma rotation in radial direction

real(kind=Float64), public :: vt =0.d0

Plasma rotation in torodial/phi direction

real(kind=Float64), public :: vz =0.d0

Plasma rotation in z direction

logical, public :: in_plasma =.False.

Indicates whether plasma parameters are valid/known

logical, public :: machine_coords =.False.

Indicates whether vectors are in machine coordinates

real(kind=Float64), public, dimension(3):: pos =0.d0

Position in beam grid coordinates

real(kind=Float64), public, dimension(3):: uvw =0.d0

Position in machine coordinates

real(kind=Float64), public, dimension(3):: vrot =0.d0

Plasma rotation in beam grid coordinates

type(InterpolCoeffs2D), public :: c

Linear Interpolation Coefficients and indicies for interpolation at pos

type, public :: EMFields

Torodial symmetric electro-magnetic fields at given R-Z

Components

TypeVisibility AttributesNameInitial
real(kind=Float64), public :: br =0.d0

Radial magnetic field [T]

real(kind=Float64), public :: bt =0.d0

Torodial magnetic field [T]

real(kind=Float64), public :: bz =0.d0

Vertical magnetic field [T]

real(kind=Float64), public :: er =0.d0

Radial electric field [V/m]

real(kind=Float64), public :: et =0.d0

Torodial electric field [V/m]

real(kind=Float64), public :: ez =0.d0

Vertical electric field [V/m]

real(kind=Float64), public :: dbr_dr =0.d0

Radial derivative of the radial magnetic field [T/m]

real(kind=Float64), public :: dbr_dz =0.d0

Vertical derivative of the radial magnetic field [T/m]

real(kind=Float64), public :: dbt_dr =0.d0

Radial derivative of the torodial magnetic field [T/m]

real(kind=Float64), public :: dbt_dz =0.d0

Vertical derivative of the torodial magnetic field [T/m]

real(kind=Float64), public :: dbz_dr =0.d0

Radial derivative of the radial magnetic field [T/m]

real(kind=Float64), public :: dbz_dz =0.d0

Vertical derivative of the vertical magnetic field [T/m]

type, public, extends(EMFields) :: LocalEMFields

Electro-magnetic fields at given position

Components

TypeVisibility AttributesNameInitial
real(kind=Float64), public :: br =0.d0

Radial magnetic field [T]

real(kind=Float64), public :: bt =0.d0

Torodial magnetic field [T]

real(kind=Float64), public :: bz =0.d0

Vertical magnetic field [T]

real(kind=Float64), public :: er =0.d0

Radial electric field [V/m]

real(kind=Float64), public :: et =0.d0

Torodial electric field [V/m]

real(kind=Float64), public :: ez =0.d0

Vertical electric field [V/m]

real(kind=Float64), public :: dbr_dr =0.d0

Radial derivative of the radial magnetic field [T/m]

real(kind=Float64), public :: dbr_dz =0.d0

Vertical derivative of the radial magnetic field [T/m]

real(kind=Float64), public :: dbt_dr =0.d0

Radial derivative of the torodial magnetic field [T/m]

real(kind=Float64), public :: dbt_dz =0.d0

Vertical derivative of the torodial magnetic field [T/m]

real(kind=Float64), public :: dbz_dr =0.d0

Radial derivative of the radial magnetic field [T/m]

real(kind=Float64), public :: dbz_dz =0.d0

Vertical derivative of the vertical magnetic field [T/m]

logical, public :: in_plasma =.False.

Indicates whether fields are valid/known

logical, public :: machine_coords =.False.

Indicates whether vectors are in machine coordinates

real(kind=Float64), public :: b_abs =0.d0

Magnitude of magnetic field

real(kind=Float64), public :: e_abs =0.d0

Magnitude of electrin field

real(kind=Float64), public, dimension(3):: pos =0.d0

Position in beam grid coordinates

real(kind=Float64), public, dimension(3):: uvw =0.d0

Position in machine coordinates

real(kind=Float64), public, dimension(3):: b_norm =0.d0

Direction of magnetic field in beam grid coordinates

real(kind=Float64), public, dimension(3):: a_norm =0.d0

Vector perpendicular to b_norm and c_norm

real(kind=Float64), public, dimension(3):: c_norm =0.d0

Vector perpendicular to b_norm and a_norm

real(kind=Float64), public, dimension(3):: e_norm =0.d0

Direction of electric field in beam grid coordinates

type(InterpolCoeffs2D), public :: c

Linear Interpolation Coefficients and indicies for interpolation at pos

type, public :: Equilibrium

MHD Equilbrium

Components

TypeVisibility AttributesNameInitial
type(EMFields), public, dimension(:,:), allocatable:: fields

Electro-magnetic fields at points defined in inter_grid

type(Profiles), public, dimension(:,:), allocatable:: plasma

Plasma parameters at points defined in inter_grid

real(kind=Float64), public, dimension(:,:), allocatable:: mask

Indicates whether fields and plasma are well-defined at points defined in inter_grid

type, public :: FastIonDistribution

Defines a Guiding Center Fast-ion Distribution Function: F(E,p,R,Z)

Components

TypeVisibility AttributesNameInitial
integer(kind=Int32), public :: nenergy

Number of energies

integer(kind=Int32), public :: npitch

Number of pitches

integer(kind=Int32), public :: nr

Number of radii

integer(kind=Int32), public :: nz

Number of z values

real(kind=Float64), public :: dE

Energy spacing [keV]

real(kind=Float64), public :: dp

Pitch spacing

real(kind=Float64), public :: dr

Radial spacing [cm]

real(kind=Float64), public :: dz

Z spacing [cm]

real(kind=Float64), public :: emin

Minimum energy [keV]

real(kind=Float64), public :: emax

Maximum energy [keV]

real(kind=Float64), public :: e_range

Energy interval length [keV]

real(kind=Float64), public :: pmin

Minimum pitch

real(kind=Float64), public :: pmax

Maximum pitch

real(kind=Float64), public :: p_range

Pitch interval length

real(kind=Float64), public :: n_tot =0.d0

Total Number of fast-ions

real(kind=Float64), public, dimension(:), allocatable:: energy

Energy values [keV]

real(kind=Float64), public, dimension(:), allocatable:: pitch

Pitch w.r.t. the magnetic field

real(kind=Float64), public, dimension(:), allocatable:: r

Radius [cm]

real(kind=Float64), public, dimension(:), allocatable:: z

Z [cm]

real(kind=Float64), public, dimension(:,:), allocatable:: denf

Fast-ion density defined on the inter_grid: denf(R,Z)

real(kind=Float64), public, dimension(:,:,:,:), allocatable:: f

Fast-ion distribution function defined on the inter_grid: F(E,p,R,Z)

type, public :: FastIon

Defines a fast-ion

Components

TypeVisibility AttributesNameInitial
logical, public :: cross_grid =.False.

Indicates whether the fast-ion crosses the beam_grid

real(kind=Float64), public :: r =0.d0

Radial position of fast-ion [cm]

real(kind=Float64), public :: z =0.d0

Vertical position of fast-ion [cm]

real(kind=Float64), public :: phi_enter =0.d0

Torodial/phi position where fast-ion enters the beam_grid [radians]

real(kind=Float64), public :: delta_phi =2*pi

Angle subtended by the beam_grid at (r,z)

real(kind=Float64), public :: energy =0.d0

Energy [keV]

real(kind=Float64), public :: pitch =0.d0

Pitch w.r.t. the magnetic field

real(kind=Float64), public :: vabs =0.d0

Speed [cm/s]

real(kind=Float64), public :: vr =0.d0

Radial velocity [cm/s]

real(kind=Float64), public :: vt =0.d0

Torodial velocity [cm/s]

real(kind=Float64), public :: vz =0.d0

Z velocity [cm/s]

real(kind=Float64), public :: weight =0.d0

Particle weight: How many fast-ions does particle represent.

integer(kind=Int32), public :: class =0

Orbit class id

type, public :: FastIonParticles

Collection of fast-ion particles

Components

TypeVisibility AttributesNameInitial
integer(kind=Int32), public :: nparticle =0

Number of particles

integer(kind=Int32), public :: nclass =1

Number of orbit classes

type(FastIon), public, dimension(:), allocatable:: fast_ion

Fast-ion particles

type, public :: NeutralBeam

Defines a neutral beam with +x defined to be into the plasma

Components

TypeVisibility AttributesNameInitial
character(len=25), public :: name =''

Beam name

integer, public :: shape

Beam source shape 1="rectangular", 2="circular"

real(kind=Float64), public :: widy

Half width of source in y direction

real(kind=Float64), public :: widz

Half height of source in z direction

real(kind=Float64), public :: focy

Focal length in y direction

real(kind=Float64), public :: focz

Focal length in z direction

real(kind=Float64), public :: einj

NBI voltage [kV]

real(kind=Float64), public :: pinj

NBI power [MW]

real(kind=Float64), public :: vinj

NBI velocity [cm/s]

real(kind=Float64), public :: alpha

Z rotation not same as beam_grid alpha

real(kind=Float64), public :: beta

Tilt rotation not same as beam_grid beta

real(kind=Float64), public, dimension(3):: divy

Energy dependent divergence in y direction

real(kind=Float64), public, dimension(3):: divz

Energy dependent divergence in z direction

real(kind=Float64), public, dimension(3):: current_fractions

Fractions of full, half, and third energy neutrals

real(kind=Float64), public, dimension(3):: src

Position of source in beam grid coordinates [cm]

real(kind=Float64), public, dimension(3):: axis

Beam centerline

integer, public :: naperture

Number of beam apertures

integer, public, dimension(:), allocatable:: ashape

Aperture shape 1="rectangular", 2="circular"

real(kind=Float64), public, dimension(:), allocatable:: awidy

Half width of the aperture(s) in y direction

real(kind=Float64), public, dimension(:), allocatable:: awidz

Half height of the aperture(s) in z direction

real(kind=Float64), public, dimension(:), allocatable:: aoffy

Horizontal (y) offset of the aperture(s) relative to the beam centerline [cm]

real(kind=Float64), public, dimension(:), allocatable:: aoffz

Vertical (z) offset of the aperture(s) relative to the beam centerline [cm]

real(kind=Float64), public, dimension(:), allocatable:: adist

Distance from the center of the beam source grid to the aperture(s) plane [cm]

real(kind=Float64), public, dimension(3,3):: basis

Beam basis for converting from centerline coordinates to beam grid coordinates

real(kind=Float64), public, dimension(3,3):: inv_basis

Inverse basis for reverse transfomation

type, public :: AtomicCrossSection

Defines a n/m-resolved atomic cross section table

Components

TypeVisibility AttributesNameInitial
integer, public :: nenergy =1

Number of beam energies

real(kind=Float64), public :: logemin =0.d0

Log-10 minimum energy

real(kind=Float64), public :: logemax =0.d0

Log-10 maximum energy

integer, public :: n_max =nlevs

Number of initial atomic energy levels

integer, public :: m_max =nlevs

Number of final atomic energy levels

real(kind=Float64), public :: dlogE =0.d0

Log-10 energy spacing

real(kind=Float64), public :: minlog_cross

Log-10 minimum cross section

real(kind=Float64), public, dimension(:,:,:), allocatable:: log_cross

Log-10 cross sections

type, public :: AtomicRates

Defines a n/m-resolved atomic cross section table

Components

TypeVisibility AttributesNameInitial
integer, public :: nenergy =1

Number of beam energies

real(kind=Float64), public :: logemin =0.d0

Log-10 minimum energy

real(kind=Float64), public :: logemax =0.d0

Log-10 maximum energy

integer, public :: ntemp =1

Number of target temperatures

real(kind=Float64), public :: logtmin =0.d0

Log-10 minimum temperature

real(kind=Float64), public :: logtmax =0.d0

Log-10 maximum temperature

integer, public :: n_max =nlevs

Number of initial atomic energy levels

integer, public :: m_max =nlevs

Number of final atomic energy levels

real(kind=Float64), public :: dlogE =0.d0

Log-10 energy spacing

real(kind=Float64), public :: dlogT =0.d0

Log-10 temperature spacing

real(kind=Float64), public :: minlog_rate =0.d0

Log-10 minimum reaction rate

real(kind=Float64), public, dimension(2):: ab =0.d0

Atomic mass of beam and thermal ions respectively [amu]

real(kind=Float64), public, dimension(:,:,:,:,:), allocatable:: log_rate

Log-10 beam-target rates

type, public :: AtomicTransitions

Defines an atomic table for populating and de-populating reaction rates

Components

TypeVisibility AttributesNameInitial
integer, public :: nenergy =1

Number of beam energies

real(kind=Float64), public :: logemin =0.d0

Log-10 minimum energy

real(kind=Float64), public :: logemax =0.d0

Log-10 maximum energy

integer, public :: ntemp =1

Number of target temperatures

real(kind=Float64), public :: logtmin =0.d0

Log-10 minimum temperature

real(kind=Float64), public :: logtmax =0.d0

Log-10 maximum temperature

integer, public :: n_max =nlevs

Number of initial atomic energy levels

integer, public :: m_max =nlevs

Number of final atomic energy levels

real(kind=Float64), public :: dlogE =0.d0

Log-10 energy spacing

real(kind=Float64), public :: dlogT =0.d0

Log-10 temperature spacing

real(kind=Float64), public :: minlog_pop =0.d0

Log-10 minimum reaction rates for populating transistions

real(kind=Float64), public :: minlog_depop =0.d0

Log-10 minimum reaction rates for de-populating transistions

real(kind=Float64), public, dimension(2):: ab =0.d0

Atomic mass of beam and thermal ions respectively [amu]

real(kind=Float64), public, dimension(:,:,:,:,:), allocatable:: log_pop

Log-10 reaction rates for populating transistions

real(kind=Float64), public, dimension(:,:,:,:), allocatable:: log_depop

Log-10 reaction rates for de-populating transistions

type, public :: NuclearRates

Nuclear reaction rates

Components

TypeVisibility AttributesNameInitial
integer, public :: nbranch =1

Number of reaction branches

integer, public :: nenergy =1

Number of beam energies

real(kind=Float64), public :: logemin =0.d0

Log-10 minimum energy

real(kind=Float64), public :: logemax =0.d0

Log-10 maximum energy

integer, public :: ntemp =1

Number of target temperatures

real(kind=Float64), public :: logtmin =0.d0

Log-10 minimum temperature

real(kind=Float64), public :: logtmax =0.d0

Log-10 maximum temperature

real(kind=Float64), public :: dlogE =0.d0

Log-10 energy spacing

real(kind=Float64), public :: dlogT =0.d0

Log-10 temperature spacing

real(kind=Float64), public :: minlog_rate =0.d0

Log-10 minimum reaction rate

real(kind=Float64), public, dimension(2):: bt_amu =0.d0

Isotope mass of beam and thermal ions respectively [amu]

real(kind=Float64), public, dimension(:,:,:), allocatable:: log_rate

Log-10 reaction rates: log_rate(energy, temperature, branch)

type, public :: AtomicTables

Atomic tables for various types of interactions

Components

TypeVisibility AttributesNameInitial
type(AtomicCrossSection), public :: H_H_cx_cross

Hydrogen-Hydrogen charge exchange n/m-resolved cross sections

type(AtomicRates), public :: H_H_cx_rate

Hydrogen-Hydrogen charge exchange n/m-resolved beam-target rates

type(AtomicTransitions), public :: H_H

Hydrogen-Hydrogen atomic transitions

type(AtomicTransitions), public :: H_e

Hydrogen-Electron atomic transitions

type(AtomicTransitions), public :: H_Aq

Hydrogen-Impurity atomic transitions

real(kind=Float64), public, dimension(nlevs,nlevs):: einstein

Einstein coefficients for spontaneous emission

type(NuclearRates), public :: D_D

Deuterium-Deuterium reaction rates

type, public :: LineOfSight

Defines a line of sight

Components

TypeVisibility AttributesNameInitial
real(kind=Float64), public :: sigma_pi =1.d0

Ratio of sigma to pi line intensity

real(kind=Float64), public :: spot_size =0.d0

Radius of spot size [cm]

real(kind=Float64), public, dimension(3):: lens =0.d0

Lens location in beam grid coordinates

real(kind=Float64), public, dimension(3):: axis =0.d0

Optical axis in beam grid coordinates

type, public :: LOSElement

Defines a element of a line of sight and cell intersection

Components

TypeVisibility AttributesNameInitial
integer, public :: id

Line of sight index

real(kind=Float64), public :: length

Length of crossing

type, public :: LOSInters

Defines the channels that intersect a cell

Components

TypeVisibility AttributesNameInitial
integer, public :: nchan =0

Number of channels that intersect

type(LOSElement), public, dimension(:), allocatable:: los_elem

Array of crossing

type, public :: SpectralChords

Defines an spectral diagnostic system

Components

TypeVisibility AttributesNameInitial
integer, public :: nchan =0

Number of channels

type(LineOfSight), public, dimension(:), allocatable:: los

Line of sight array

real(kind=Float64), public, dimension(:), allocatable:: radius

Radius of each line of sight

type(LOSInters), public, dimension(:,:,:), allocatable:: inter

Array of LOS intersections with beam_grid

type, public :: BoundedPlane

Defines a plane with a circular or rectangular boundary

Components

TypeVisibility AttributesNameInitial
integer, public :: shape =0

Boundary shape 1="Rectangular", 2="circular"

real(kind=Float64), public :: hh =0.d0

Half height of boundary [cm]

real(kind=Float64), public :: hw =0.d0

Half width of boundary [cm]

real(kind=Float64), public, dimension(3):: origin =0.d0

Origin of plane in machine coordinates

real(kind=Float64), public, dimension(3,3):: basis =0.d0

Basis vectors basis(:,1) = u_1 is plane normal

real(kind=Float64), public, dimension(3,3):: inv_basis =0.d0

Inverse basis

type, public :: NPADetector

Defines a NPA detector

Components

TypeVisibility AttributesNameInitial
type(BoundedPlane), public :: detector

Detecting plane of NPA detector

type(BoundedPlane), public :: aperture

Aperture plane of NPA detector

type, public :: NPAProbability

Type to contain the probability of hitting a NPA detector

Components

TypeVisibility AttributesNameInitial
real(kind=Float64), public :: p =0.d0

Hit probability

real(kind=Float64), public :: pitch =-2.d0

Pitch

real(kind=Float64), public, dimension(3):: eff_rd =0.d0

Effective position of detector

real(kind=Float64), public, dimension(3):: dir =0.d0

Trajectory direction

type, public :: NPAChords

Defines a NPA system

Components

TypeVisibility AttributesNameInitial
integer, public :: nchan =0

Number of channels

type(NPADetector), public, dimension(:), allocatable:: det

NPA detector array

real(kind=Float64), public, dimension(:), allocatable:: radius

Radius [cm]

logical, public, dimension(:,:,:), allocatable:: hit

Indicates whether a particle can hit any NPA detector from a grid cell: hit(x,y,z)

type(NPAProbability), public, dimension(:,:,:,:), allocatable:: phit

Probability of hitting a detector from a grid cell: phit(x,y,z,chan)

type, public :: NPAParticle

Defines a NPA particle

Components

TypeVisibility AttributesNameInitial
integer, public :: detector =0

Detector NPA particle hit

real(kind=Float64), public :: xi =0.d0

Initial x position

real(kind=Float64), public :: yi =0.d0

Initial y position

real(kind=Float64), public :: zi =0.d0

Initial z position

real(kind=Float64), public :: xf =0.d0

Final x position

real(kind=Float64), public :: yf =0.d0

Final y position

real(kind=Float64), public :: zf =0.d0

Final z position

real(kind=Float64), public :: weight =0.d0

NPA particle weight

real(kind=Float64), public :: energy =0.d0

Birth Energy [keV]

real(kind=Float64), public :: pitch =0.d0

Birth Pitch

type, public :: NPAResults

MC NPA result structure

Components

TypeVisibility AttributesNameInitial
integer(kind=Int32), public :: nchan =0

Number of NPA channels

integer(kind=Int32), public :: npart =0

Number of particles that hit a detector

integer(kind=Int32), public :: nmax =1000000

Maximum allowed number of particles grows if necessary

integer(kind=Int32), public :: nenergy =100

Number of energy values

type(NPAParticle), public, dimension(:), allocatable:: part

Array of NPA particles

real(kind=Float64), public, dimension(:), allocatable:: energy

Energy array [keV]

real(kind=Float64), public, dimension(:,:,:), allocatable:: flux

Neutral particle flux: flux(energy,chan, orbit_type) [neutrals/(s*dE)]

type, public :: BirthProfile

Birth profile structure

Components

TypeVisibility AttributesNameInitial
integer, public :: cnt =1

Particle counter

integer, public, dimension(:), allocatable:: neut_type

Particle birth type (1=Full, 2=Half, 3=Third)

real(kind=Float64), public, dimension(:,:), allocatable:: ri

Particle birth position [cm]

real(kind=Float64), public, dimension(:,:), allocatable:: vi

Particle birth velocity [cm/s]

integer, public, dimension(:,:), allocatable:: ind

Particle beam_grid indices

real(kind=Float64), public, dimension(:,:,:,:), allocatable:: dens

Birth density: dens(neutral_type,x,y,z) [fast-ions/(s*cm^3)]

type, public :: Spectra

Spectra storage structure

Components

TypeVisibility AttributesNameInitial
real(kind=Float64), public, dimension(:,:), allocatable:: brems

Bremsstruhlung: brems(lambda,chan)

real(kind=Float64), public, dimension(:,:,:), allocatable:: bes

Beam emission: bes(lambda,chan,neutral_type)

real(kind=Float64), public, dimension(:,:,:), allocatable:: fida

FIDA emission: fida(lambda,chan,orbit_type)

type, public :: NeutronRate

Neutron storage structure

Components

TypeVisibility AttributesNameInitial
real(kind=Float64), public, dimension(:), allocatable:: rate

Neutron rate: rate(orbit_type) [neutrons/sec]

real(kind=Float64), public, dimension(:,:,:,:), allocatable:: weight

Neutron rate weight: weight(E,p,R,Z)

type, public :: NeutralDensity

Neutral density structure

Components

TypeVisibility AttributesNameInitial
real(kind=Float64), public, dimension(:,:,:,:,:), allocatable:: dens

Neutral density: dens(lev,neutral_type,x,y,z)

type, public :: FIDAWeights

FIDA weights structure

Components

TypeVisibility AttributesNameInitial
real(kind=Float64), public, dimension(:,:,:), allocatable:: mean_f

Estimate of mean fast-ion distribution function "seen" by LOS: mean_f(E,p,chan)

real(kind=Float64), public, dimension(:,:,:,:), allocatable:: weight

FIDA weight function: weight(lambda,E,p,chan)

type, public :: NPAWeights

NPA weights structure

Components

TypeVisibility AttributesNameInitial
real(kind=Float64), public, dimension(:,:,:,:,:), allocatable:: attenuation

Attenuation fraction: attenuation(E,x,y,z,chan)

real(kind=Float64), public, dimension(:,:,:,:,:), allocatable:: cx

Charge Exchange reaction rates: cx(E,x,y,z,chan)

real(kind=Float64), public, dimension(:,:,:,:), allocatable:: emissivity

Emissivity: emissivity(x,y,z,chan) [neutrals/(s*dV)]

real(kind=Float64), public, dimension(:,:,:), allocatable:: weight

NPA weight function: weight(E,p,chan) [neutrals/(sfast-iondE*dP)]

real(kind=Float64), public, dimension(:,:), allocatable:: flux

Neutral particle flux: flux(E,chan) [neutrals/(s*dE)]

type, public :: SimulationInputs

Simulation settings structure

Components

TypeVisibility AttributesNameInitial
integer(kind=Int32), public :: shot_number

Shot Number

real(kind=Float64), public :: time

Shot time [s]

character(len=charlim), public :: runid =''

FIDASIM run ID

character(len=charlim), public :: result_dir =''

Result directory

character(len=charlim), public :: tables_file =''

Atomic tables file

character(len=charlim), public :: geometry_file =''

FIDASIM input file containing geometric quantities

character(len=charlim), public :: equilibrium_file =''

FIDASIM input file containing the plasma parameters and fields

character(len=charlim), public :: distribution_file =''

FIDASIM input file containing the fast-ion distribution

character(len=charlim), public :: neutrals_file =''

FIDASIM output/input file containing beam neutral density. Used when load_neutrals is set.

integer(kind=Int64), public :: n_fida

Number of FIDA mc markers

integer(kind=Int64), public :: n_npa

Number of NPA mc markers

integer(kind=Int64), public :: n_nbi

Number of neutral beam mc markers

integer(kind=Int64), public :: n_dcx

Number of direct charge exchange (DCX) mc markers

integer(kind=Int64), public :: n_halo

Number of halo mc markers

integer(kind=Int64), public :: n_birth

Number of birth particles per n_nbi

integer(kind=Int32), public :: calc_spec

Calculate spectra: 0 = off, 1=on

integer(kind=Int32), public :: calc_brems

Calculate bremmstruhlung: 0 = off, 1=on

integer(kind=Int32), public :: calc_bes

Calculate BES: 0 = off, 1=on

integer(kind=Int32), public :: calc_fida

Calculate FIDA: 0 = off, 1=on

integer(kind=Int32), public :: load_neutrals

Load neutrals from file: 0 = off, 1=on

integer(kind=Int32), public :: calc_npa

Calculate NPA: 0 = off, 1=on, 2=on++

integer(kind=Int32), public :: calc_fida_wght

Calculate FIDA weight: 0 = off, 1=on, 2=on++

integer(kind=Int32), public :: calc_npa_wght

Calculate NPA weights: 0 = off, 1=on, 2=on++

integer(kind=Int32), public :: calc_birth

Calculate birth profile: 0 = off, 1=on

integer(kind=Int32), public :: calc_neutron

Calculate neutron flux: 0 = off, 1=on

integer(kind=Int32), public :: no_flr

Turns off Finite Larmor Radius effects: 0=off, 1=on

integer(kind=Int32), public :: dump_dcx

Output DCX density and spectra: 0 = off, 1=on

integer(kind=Int32), public :: verbose
real(kind=Float64), public :: ab

Atomic mass of beam neutrals

integer(kind=Int32), public :: impurity_charge

Impurity proton number

real(kind=Float64), public :: ai

Atomic mass of thermal ions

integer(kind=Int32), public :: dist_type

Type of fast-ion distribution

integer(kind=Int32), public :: nlambda

Number of wavelength to calculate

real(kind=Float64), public :: dlambda

Wavelength spacing [nm]

real(kind=Float64), public :: lambdamin

Minimum wavelength [nm]

real(kind=Float64), public :: lambdamax

Maximum wavelength [nm]

integer(kind=Int32), public :: ne_wght

Number of energies in weight functions

integer(kind=Int32), public :: np_wght

Number of pitches in weight functions

integer(kind=Int32), public :: nphi_wght

Number of gyro-angles to average over in weight functions

integer(kind=Int32), public :: nlambda_wght

Number of wavelength to calculate in weight functions

real(kind=Float64), public :: emax_wght

Maximum energy in weight functions [keV]

real(kind=Float64), public :: lambdamin_wght

Minimum wavelength in weight functions [nm]

real(kind=Float64), public :: lambdamax_wght

Maximum wavelength in weight functions [nm]

type, public :: ParticleTrack

Stores properties seen when traveling through a 3D grid

Components

TypeVisibility AttributesNameInitial
real(kind=Float64), public :: time =0.d0

Time/distance/... in cell

real(kind=Float64), public :: flux =0.d0

Flux/density/... in cell

real(kind=Float64), public, dimension(nlevs):: dens =0.d0

Density [cm^-3]

integer(kind=Int32), public, dimension(3):: ind =0

Indices of cell

real(kind=Float64), public, dimension(3):: pos =0.d0

Midpoint of track in cell [cm]

logical, public :: in_plasma =.False.

Indicates whether we are in the plasma

type, public :: GyroSurface

Surface containing the fast-ion velocity vectors for all values of the gyro-angle. It takes the form of a hyperboloid where is the gyro-angle, is the ion gyro-frequency and

Components

TypeVisibility AttributesNameInitial
real(kind=Float64), public :: v =0.d0

Particle speed

real(kind=Float64), public :: omega =0.d0

Ion gyro-frequency

real(kind=Float64), public, dimension(3):: axes

Semi-axes of the hyperboloid, i.e. a, b, c coefficients

real(kind=Float64), public, dimension(3):: center =0.d0

Center of the gyrosurface

real(kind=Float64), public, dimension(3,3):: A =0.d0

Coefficients of quartic surface i.e. basis*diagm(1/a^2,1/b^2,1/c^2)*basis'

real(kind=Float64), public, dimension(3,3):: basis =0.d0

Basis of coordinate system of gyrosurface


Functions

public function pp_add(p1, p2) result(p3)

Defines how to add two Profiles types

Arguments

Type IntentOptional AttributesName
type(Profiles), intent(in) :: p1
type(Profiles), intent(in) :: p2

Return Value type(Profiles)

public function pp_subtract(p1, p2) result(p3)

Defines how to subtract two Profiles types

Arguments

Type IntentOptional AttributesName
type(Profiles), intent(in) :: p1
type(Profiles), intent(in) :: p2

Return Value type(Profiles)

public function lplp_add(p1, p2) result(p3)

Defines how to add two LocalProfiles types

Arguments

Type IntentOptional AttributesName
type(LocalProfiles), intent(in) :: p1
type(LocalProfiles), intent(in) :: p2

Return Value type(LocalProfiles)

public function lplp_subtract(p1, p2) result(p3)

Defines how to subtract two LocalProfiles types

Arguments

Type IntentOptional AttributesName
type(LocalProfiles), intent(in) :: p1
type(LocalProfiles), intent(in) :: p2

Return Value type(LocalProfiles)

public function ps_multiply(p1, real_scalar) result(p3)

Defines how to multiply Profiles types by a scalar

Arguments

Type IntentOptional AttributesName
type(Profiles), intent(in) :: p1
real(kind=Float64), intent(in) :: real_scalar

Return Value type(Profiles)

public function sp_multiply(real_scalar, p1) result(p3)

Defines how to multiply Profiles types by a scalar

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: real_scalar
type(Profiles), intent(in) :: p1

Return Value type(Profiles)

public function ps_divide(p1, real_scalar) result(p3)

Defines how to divide Profiles types by a scalar

Arguments

Type IntentOptional AttributesName
type(Profiles), intent(in) :: p1
real(kind=Float64), intent(in) :: real_scalar

Return Value type(Profiles)

public function lps_multiply(p1, real_scalar) result(p3)

Defines how to multiply LocalProfiles types by a scalar

Arguments

Type IntentOptional AttributesName
type(LocalProfiles), intent(in) :: p1
real(kind=Float64), intent(in) :: real_scalar

Return Value type(LocalProfiles)

public function slp_multiply(real_scalar, p1) result(p3)

Defines how to multiply LocalProfiles types by a scalar

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: real_scalar
type(LocalProfiles), intent(in) :: p1

Return Value type(LocalProfiles)

public function lps_divide(p1, real_scalar) result(p3)

Defines how to divide LocalProfiles types by a scalar

Arguments

Type IntentOptional AttributesName
type(LocalProfiles), intent(in) :: p1
real(kind=Float64), intent(in) :: real_scalar

Return Value type(LocalProfiles)

public function ff_add(p1, p2) result(p3)

Defines how to add two EMFields types

Arguments

Type IntentOptional AttributesName
type(EMFields), intent(in) :: p1
type(EMFields), intent(in) :: p2

Return Value type(EMFields)

public function ff_subtract(p1, p2) result(p3)

Defines how to subtract two EMFields types

Arguments

Type IntentOptional AttributesName
type(EMFields), intent(in) :: p1
type(EMFields), intent(in) :: p2

Return Value type(EMFields)

public function fs_multiply(p1, real_scalar) result(p3)

Defines how to multiply EMFields types by a scalar

Arguments

Type IntentOptional AttributesName
type(EMFields), intent(in) :: p1
real(kind=Float64), intent(in) :: real_scalar

Return Value type(EMFields)

public function sf_multiply(real_scalar, p1) result(p3)

Defines how to multiply EMFields types by a scalar

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: real_scalar
type(EMFields), intent(in) :: p1

Return Value type(EMFields)

public function fs_divide(p1, real_scalar) result(p3)

Defines how to divide EMFields types by a scalar

Arguments

Type IntentOptional AttributesName
type(EMFields), intent(in) :: p1
real(kind=Float64), intent(in) :: real_scalar

Return Value type(EMFields)

public function lflf_add(p1, p2) result(p3)

Defines how to add two LocalEMFields types

Arguments

Type IntentOptional AttributesName
type(LocalEMFields), intent(in) :: p1
type(LocalEMFields), intent(in) :: p2

Return Value type(LocalEMFields)

public function lflf_subtract(p1, p2) result(p3)

Defines how to subtract two LocalEMFields types

Arguments

Type IntentOptional AttributesName
type(LocalEMFields), intent(in) :: p1
type(LocalEMFields), intent(in) :: p2

Return Value type(LocalEMFields)

public function lfs_multiply(p1, real_scalar) result(p3)

Defines how to multiply LocalEMFields types by a scalar

Arguments

Type IntentOptional AttributesName
type(LocalEMFields), intent(in) :: p1
real(kind=Float64), intent(in) :: real_scalar

Return Value type(LocalEMFields)

public function slf_multiply(real_scalar, p1) result(p3)

Defines how to multiply LocalEMFields types by a scalar

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: real_scalar
type(LocalEMFields), intent(in) :: p1

Return Value type(LocalEMFields)

public function lfs_divide(p1, real_scalar) result(p3)

Defines how to divide LocalEMFields types by a scalar

Arguments

Type IntentOptional AttributesName
type(LocalEMFields), intent(in) :: p1
real(kind=Float64), intent(in) :: real_scalar

Return Value type(LocalEMFields)

public function approx_eq(x, y, tol) result(a)

Inexact equality comparison: x ~= y true if abs(x-y) <= tol else false

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: x

First value in comparison

real(kind=Float64), intent(in) :: y

Second value in comparison

real(kind=Float64), intent(in) :: tol

Equality tolerance

Return Value logical

public function approx_ge(x, y, tol) result(a)

Inexact greater than or equal to comparison: x >~= y

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: x

First value in comparison

real(kind=Float64), intent(in) :: y

Second value in comparison

real(kind=Float64), intent(in) :: tol

Equality tolerance

Return Value logical

public function approx_le(x, y, tol) result(a)

Inexact less then or equal to comparison: x <~= y

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: x

First value in comparison

real(kind=Float64), intent(in) :: y

Second value in comparison

real(kind=Float64), intent(in) :: tol

Equality tolerance

Return Value logical

public function cross_product(u, v) result(s)

Calculates the cross product of two vectors: uxv

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: u
real(kind=Float64), intent(in), dimension(3):: v

Return Value real(kind=Float64), dimension(3)

public function in_boundary(bplane, p) result(in_b)

Indicator function for determining if a point on a plane is within the plane boundary

Arguments

Type IntentOptional AttributesName
type(BoundedPlane), intent(in) :: bplane

Plane with boundary

real(kind=Float64), intent(in), dimension(3):: p

Point on plane

Return Value logical

public function in_gyro_surface(gs, p) result(in_gs)

Indicator function for determining if a point is inside the gyro_surface

Arguments

Type IntentOptional AttributesName
type(GyroSurface), intent(in) :: gs

Gyro-surface

real(kind=Float64), intent(in), dimension(3):: p

Point

Return Value logical

public function in_grid(xyz) result(ing)

Determines if a position pos is in the beam_grid

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: xyz

Position in beam grid coordinates [cm]

Return Value logical

Indicates whether the position is in the beam grid

public function gyro_radius(fields, energy, pitch) result(gyro_rad)

Calculates mean gyro-radius

Arguments

Type IntentOptional AttributesName
type(LocalEMFields), intent(in) :: fields

Electromagnetic fields at guiding center

real(kind=Float64), intent(in) :: energy

Energy of particle

real(kind=Float64), intent(in) :: pitch

Particle pitch w.r.t the magnetic field

Return Value real(kind=Float64)

Mean gyro-radius


Subroutines

public subroutine print_banner()

Prints FIDASIM banner

Arguments

None

public subroutine fast_ion_assign(p1, p2)

Defines how to assign FastIon types to eachother

Arguments

Type IntentOptional AttributesName
type(FastIon), intent(out) :: p1
type(FastIon), intent(in) :: p2

public subroutine npa_part_assign(p1, p2)

Defines how to assign NPAParticle types to eachother

Arguments

Type IntentOptional AttributesName
type(NPAParticle), intent(out) :: p1
type(NPAParticle), intent(in) :: p2

public subroutine pp_assign(p1, p2)

Defines how to assign Profiles types to eachother

Arguments

Type IntentOptional AttributesName
type(Profiles), intent(inout) :: p1
type(Profiles), intent(in) :: p2

public subroutine lpp_assign(p1, p2)

Defines how to assign a Profiles type to a LocalProfiles type

Arguments

Type IntentOptional AttributesName
type(LocalProfiles), intent(inout) :: p1
type(Profiles), intent(in) :: p2

public subroutine plp_assign(p1, p2)

Defines how to assign a LocalProfiles type to a Profiles type

Arguments

Type IntentOptional AttributesName
type(Profiles), intent(inout) :: p1
type(LocalProfiles), intent(in) :: p2

public subroutine lplp_assign(p1, p2)

Defines how to assign LocalProfiles types to eachother

Arguments

Type IntentOptional AttributesName
type(LocalProfiles), intent(inout) :: p1
type(LocalProfiles), intent(in) :: p2

public subroutine ff_assign(p1, p2)

Defines how to assign EMFields types to eachother

Arguments

Type IntentOptional AttributesName
type(EMFields), intent(inout) :: p1
type(EMFields), intent(in) :: p2

public subroutine lff_assign(p1, p2)

Defines how to assign a EMFields type to a LocalEMFields type

Arguments

Type IntentOptional AttributesName
type(LocalEMFields), intent(inout) :: p1
type(EMFields), intent(in) :: p2

public subroutine flf_assign(p1, p2)

Defines how to assign a LocalEMFields type to a EMFields type

Arguments

Type IntentOptional AttributesName
type(EMFields), intent(inout) :: p1
type(LocalEMFields), intent(in) :: p2

public subroutine lflf_assign(p1, p2)

Defines how to assign LocalEMFields types to eachother

Arguments

Type IntentOptional AttributesName
type(LocalEMFields), intent(inout) :: p1
type(LocalEMFields), intent(in) :: p2

public subroutine read_inputs()

Reads input namelist file and stores the results into inputs, nbi, and beam_grid 20 for suffixes and seperators e.g. /, _npa.h5, ...

Arguments

None

public subroutine make_beam_grid()

Makes [[libfida:beam_grid] from user defined inputs

Arguments

None

public subroutine read_beam()

Reads neutral beam geometry and stores the quantities in nbi

Arguments

None

public subroutine read_chords()

Reads the spectral geometry and stores the quantities in spec_chords

Arguments

None

public subroutine read_npa()

Reads the NPA geometry and stores the quantities in npa_chords

Arguments

None

public subroutine read_equilibrium()

Reads in the interpolation grid, plasma parameters, and fields and stores the quantities in inter_grid and equil

Arguments

None

public subroutine read_f(fid, error)

Reads in the fast-ion distribution function and stores the quantities in fbm

Arguments

Type IntentOptional AttributesName
integer(kind=HID_T), intent(inout) :: fid

HDF5 file ID

integer, intent(out) :: error

Error code

public subroutine read_mc(fid, error)

Reads in a MC particle fast-ion distribution and puts them in particles

Arguments

Type IntentOptional AttributesName
integer(kind=HID_T), intent(inout) :: fid

HDF5 file ID

integer, intent(out) :: error

Error code

public subroutine read_distribution()

Reads in the fast-ion distribution

Arguments

None

public subroutine read_atomic_cross(fid, grp, cross)

Reads in a cross section table from file and puts it into a AtomicCrossSection type

Arguments

Type IntentOptional AttributesName
integer(kind=HID_T), intent(in) :: fid

HDF5 file ID

character(len=*), intent(in) :: grp

HDF5 group to read from

type(AtomicCrossSection), intent(inout) :: cross

Atomic cross section

public subroutine read_atomic_rate(fid, grp, b_amu, t_amu, rates)

Reads in a atomic rate table from file and puts it into a AtomicRates type

Arguments

Type IntentOptional AttributesName
integer(kind=HID_T), intent(in) :: fid

HDF5 file ID

character(len=*), intent(in) :: grp

HDF5 group to read from

real(kind=Float64), intent(in), dimension(2):: b_amu

Atomic masses of "beam" species (beam ion and thermal ion)

real(kind=Float64), intent(in) :: t_amu

Atomic mass of "target" species (thermal ion)

type(AtomicRates), intent(inout) :: rates

Atomic reaction rates

public subroutine read_atomic_transitions(fid, grp, b_amu, t_amu, rates)

Reads in a atomic transitions table from file and puts it into a AtomicTransitions type

Arguments

Type IntentOptional AttributesName
integer(kind=HID_T), intent(in) :: fid

HDF5 file ID

character(len=*), intent(in) :: grp

HDF5 group to read from

real(kind=Float64), intent(in), dimension(2):: b_amu

Atomic masses of "beam" species (beam ion and thermal ion)

real(kind=Float64), intent(in) :: t_amu

Atomic mass of "target" species (thermal ion)

type(AtomicTransitions), intent(inout) :: rates

Atomic transitions

public subroutine read_nuclear_rates(fid, grp, rates)

Reads in a nuclear reaction rates table from file and puts it into a NuclearRates type

Arguments

Type IntentOptional AttributesName
integer(kind=HID_T), intent(in) :: fid

HDF5 file ID

character(len=*), intent(in) :: grp

HDF5 group to read from

type(NuclearRates), intent(inout) :: rates

Atomic reaction rates

public subroutine read_tables()

Reads in atomic tables from file and stores them in tables

Arguments

None

public subroutine write_beam_grid(id, error)

Write beam_grid to an HDF5 file

Arguments

Type IntentOptional AttributesName
integer(kind=HID_T), intent(inout) :: id

HDF5 file ID

integer, intent(out) :: error

Error code

public subroutine write_birth_profile()

Writes birth to a HDF5 file

Arguments

None

public subroutine write_dcx()

Writes the direct charge exchange (DCX) neutrals and spectra to a HDF5 file

Arguments

None

public subroutine write_neutrals()

Writes neut to a HDF5 file

Arguments

None

public subroutine write_npa()

Writes npa to a HDF5 file

Arguments

None

public subroutine write_spectra()

Writes Spectra to a HDF5 file

Arguments

None

public subroutine write_neutrons()

Writes neutron to a HDF5 file

Arguments

None

public subroutine write_fida_weights()

Writes fweight to a HDF5 file

Arguments

None

public subroutine write_npa_weights()

Writes nweight to a HDF5 file

Arguments

None

public subroutine read_neutrals()

Reads neutral density from file and puts it in neut

Arguments

None

public subroutine tb_zyx(alpha, beta, gamma, basis, inv_basis)

Creates active rotation matrix for z-y'-x" rotation given Tait-Bryan angles

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: alpha

Angle of rotation about z

real(kind=Float64), intent(in) :: beta

Angle of rotation about y'

real(kind=Float64), intent(in) :: gamma

Angle of rotation about x"

real(kind=Float64), intent(out), dimension(3,3):: basis

Rotation matrix/basis for transforming from rotated to non-rotated coordinates

real(kind=Float64), intent(out), optional dimension(3,3):: inv_basis

Inverse basis for reverse transformation

public subroutine line_basis(r0, v0, basis, inv_basis)

Calculates basis from a line with +x in the direction of line

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: r0

Starting point of line [cm]

real(kind=Float64), intent(in), dimension(3):: v0

Direction of line

real(kind=Float64), intent(out), dimension(3,3):: basis

Basis for transforming from line coordinates to cartesian

real(kind=Float64), intent(out), optional dimension(3,3):: inv_basis

Inverse basis for the reverse transformation cartesian to line

public subroutine plane_basis(center, redge, tedge, basis, inv_basis)

Calculates basis from 3 points on a plane with +z being the plane normal

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: center

Plane origin

real(kind=Float64), intent(in), dimension(3):: redge

Right edge of plane

real(kind=Float64), intent(in), dimension(3):: tedge

Top edge of plane

real(kind=Float64), intent(out), dimension(3,3):: basis

Basis for transforming from plane to cartesian coordinates

real(kind=Float64), intent(out), optional dimension(3,3):: inv_basis

Inverse basis for the reverse transformation cartesian to plane

public subroutine line_plane_intersect(l0, l, p0, n, p, t)

Calculates the intersection of a line and a plane

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: l0

Point on line

real(kind=Float64), intent(in), dimension(3):: l

Ray of line

real(kind=Float64), intent(in), dimension(3):: p0

Point on plane

real(kind=Float64), intent(in), dimension(3):: n

Normal vector of plane

real(kind=Float64), intent(out), dimension(3):: p

Line-plane intersect point

real(kind=Float64), intent(out) :: t

"time" to intersect

public subroutine boundary_edge(bplane, bedge, nb)

Returns 3 x nb array containing points along the BoundedPlane's boundary edge

Arguments

Type IntentOptional AttributesName
type(BoundedPlane), intent(in) :: bplane

Bounded plane

real(kind=Float64), intent(out), dimension(:,:):: bedge

Boundary edge points of bounded plane

integer, intent(out) :: nb

Number of points in boundary edge

public subroutine gyro_surface(fields, energy, pitch, gs)

Calculates the surface of all possible trajectories

Arguments

Type IntentOptional AttributesName
type(LocalEMFields), intent(in) :: fields

Electromagnetic fields at guiding center

real(kind=Float64), intent(in) :: energy

Energy of particle

real(kind=Float64), intent(in) :: pitch

Particle pitch w.r.t the magnetic field

type(GyroSurface), intent(out) :: gs

Gyro-surface

public subroutine line_gyro_surface_intersect(r0, v0, gs, t)

Calculates the times of intersection of a line and a gyro-surface

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: r0

Point on line

real(kind=Float64), intent(in), dimension(3):: v0

Direction of line

type(GyroSurface), intent(in) :: gs

Gyro-surface

real(kind=Float64), intent(out), dimension(2):: t

"time" to intersect

public subroutine gyro_surface_coordinates(gs, p, u)

Calculates the parametric coordinates, u, of point p on the gyro_surface

Arguments

Type IntentOptional AttributesName
type(GyroSurface), intent(in) :: gs

Gyro_surface

real(kind=Float64), intent(in), dimension(3):: p

Point on gyro_surface

real(kind=Float64), intent(out), dimension(2):: u

Parametric coordinates (gyro-angle, t)

public subroutine gyro_trajectory(gs, theta, ri, vi)

Calculate particle trajectory for a given gyro-angle and gyro-surface

Arguments

Type IntentOptional AttributesName
type(GyroSurface), intent(in) :: gs

Gyro-Surface

real(kind=Float64), intent(in) :: theta

Gyro-angle

real(kind=Float64), dimension(3):: ri

Particle position

real(kind=Float64), dimension(3):: vi

Particle Velocity

public subroutine gyro_range(b, gs, gyrange, nrange)

Calculates the range(s) of gyro-angles that would land within a bounded plane

Arguments

Type IntentOptional AttributesName
type(BoundedPlane), intent(in) :: b

Bounded Plane

type(GyroSurface), intent(in) :: gs

Gyro-surface

real(kind=Float64), intent(out), dimension(2,4):: gyrange

(theta, dtheta) values

integer, intent(out) :: nrange

Number of ranges. 1 <= nrange <= 4

public subroutine npa_gyro_range(ichan, gs, gyrange, nrange)

Calculates range of gyro-angles that would hit the NPA detector

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: ichan

Index of NPA detector

type(GyroSurface), intent(in) :: gs
real(kind=Float64), intent(out), dimension(2,4):: gyrange
integer, intent(out) :: nrange

public subroutine hit_npa_detector(r0, v0, d_index, rd, det)

Routine to check if a particle will hit a NPA detector

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: r0

Starting point of particle

real(kind=Float64), intent(in), dimension(3):: v0

Particle velocity

integer, intent(out) :: d_index

Index of NPA detector. Zero if particle doesn't hit

real(kind=Float64), intent(out), optional dimension(3):: rd

Point where particle hit detector

integer, intent(in), optional :: det

Index of NPA detector to check

public subroutine xyz_to_uvw(xyz, uvw)

Convert beam coordinate xyz to machine coordinate uvw

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: xyz
real(kind=Float64), intent(out), dimension(3):: uvw

public subroutine uvw_to_xyz(uvw, xyz)

Convert machine coordinate uvw to beam coordinate xyz

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: uvw
real(kind=Float64), intent(out), dimension(3):: xyz

public subroutine grid_intersect(r0, v0, length, r_enter, r_exit, center_in, lwh_in)

Calculates a particles intersection length with the beam_grid

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: r0

Initial position of particle [cm]

real(kind=Float64), intent(in), dimension(3):: v0

Velocity of particle [cm/s]

real(kind=Float64), intent(out) :: length

Intersection length [cm]

real(kind=Float64), intent(out), dimension(3):: r_enter

Point where particle enters beam_grid

real(kind=Float64), intent(out), dimension(3):: r_exit

Point where particle exits beam_grid

real(kind=Float64), intent(in), optional dimension(3):: center_in

Alternative grid center

real(kind=Float64), intent(in), optional dimension(3):: lwh_in

Alternative grid [length,width,height]

public subroutine circle_grid_intersect(r0, e1, e2, radius, phi_enter, phi_exit)

Calculates the intersection arclength of a circle with the beam_grid

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: r0

Position of center enter of the circle in beam grid coordinates [cm]

real(kind=Float64), intent(in), dimension(3):: e1

Unit vector pointing towards (R, 0) (r,phi) position of the circle in beam grid coordinates

real(kind=Float64), intent(in), dimension(3):: e2

Unit vector pointing towards (R, pi/2) (r,phi) position of the circle in beam grid coordinates

real(kind=Float64), intent(in) :: radius

Radius of circle [cm]

real(kind=Float64), intent(out) :: phi_enter

Phi value where the circle entered the beam_grid [rad]

real(kind=Float64), intent(out) :: phi_exit

Phi value where the circle exits the beam_grid [rad]

public subroutine get_indices(pos, ind)

Find closests beam_grid indices ind to position pos

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: pos

Position [cm]

integer(kind=Int32), intent(out), dimension(3):: ind

Closest indices to position

public subroutine get_position(ind, pos)

Get position pos given beam_grid indices ind

Arguments

Type IntentOptional AttributesName
integer(kind=Int32), intent(in), dimension(3):: ind

beam_grid indices

real(kind=Float64), intent(out), dimension(3):: pos

Position [cm]

public subroutine track(rin, vin, tracks, ncell, los_intersect)

Computes the path of a neutral through the beam_grid

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: rin

Initial position of particle

real(kind=Float64), intent(in), dimension(3):: vin

Initial velocity of particle

type(ParticleTrack), intent(inout), dimension(:):: tracks

Array of ParticleTrack type

integer(kind=Int32), intent(out) :: ncell

Number of cells that a particle crosses

logical, intent(out), optional :: los_intersect

Indicator whether particle intersects a LOS in spec_chords

public subroutine interpol1D_coeff(xmin, dx, nx, xout, c, err)

Linear interpolation coefficients and index for a 1D grid y(x)

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: xmin

Minimum abscissa value

real(kind=Float64), intent(in) :: dx

Absissa spacing

integer, intent(in) :: nx

Number of abscissa

real(kind=Float64), intent(in) :: xout

Abscissa value to interpolate

type(InterpolCoeffs1D), intent(out) :: c

Interpolation Coefficients

integer, intent(out), optional :: err

Error code

public subroutine interpol1D_coeff_arr(x, xout, c, err)

Linear interpolation coefficients and index for a 1D grid y(x)

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(:):: x

Abscissa values

real(kind=Float64), intent(in) :: xout

Abscissa value to interpolate

type(InterpolCoeffs1D), intent(out) :: c

Interpolation Coefficients

integer, intent(out), optional :: err

Error code

public subroutine interpol2D_coeff(xmin, dx, nx, ymin, dy, ny, xout, yout, c, err)

Bilinear interpolation coefficients and indicies for a 2D grid z(x,y)

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: xmin

Minimum abscissa

real(kind=Float64), intent(in) :: dx

Abscissa spacing

integer, intent(in) :: nx

Number of abscissa

real(kind=Float64), intent(in) :: ymin

Minimum ordinate

real(kind=Float64), intent(in) :: dy

Ordinate spacing

integer, intent(in) :: ny

Number of ordinates points

real(kind=Float64), intent(in) :: xout

Abscissa value to interpolate

real(kind=Float64), intent(in) :: yout

Ordinate value to interpolate

type(InterpolCoeffs2D), intent(out) :: c

Interpolation Coefficients

integer, intent(out), optional :: err

Error code

public subroutine interpol2D_coeff_arr(x, y, xout, yout, c, err)

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(:):: x

Abscissa values

real(kind=Float64), intent(in), dimension(:):: y

Ordinate values

real(kind=Float64), intent(in) :: xout

Abscissa value to interpolate

real(kind=Float64), intent(in) :: yout

Ordinate value to interpolate

type(InterpolCoeffs2D), intent(out) :: c

Interpolation Coefficients

integer, intent(out), optional :: err

Error code

public subroutine interpol1D_arr(x, y, xout, yout, err, coeffs)

Performs linear interpolation on a uniform 1D grid y(x)

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(:):: x

The abscissa values of y

real(kind=Float64), intent(in), dimension(:):: y

Values at abscissa values x: y(x)

real(kind=Float64), intent(in) :: xout

Abscissa value to interpolate

real(kind=Float64), intent(out) :: yout

Interpolant: y(xout)

integer, intent(out), optional :: err

Error code

type(InterpolCoeffs1D), intent(in), optional :: coeffs

Precomputed Linear Interpolation Coefficients

public subroutine interpol2D_arr(x, y, z, xout, yout, zout, err, coeffs)

Performs bilinear interpolation on a 2D grid z(x,y)

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(:):: x

The abscissa values of z

real(kind=Float64), intent(in), dimension(:):: y

The ordinate values of z

real(kind=Float64), intent(in), dimension(:,:):: z

Values at the abscissa/ordinates: z(x,y)

real(kind=Float64), intent(in) :: xout

The abscissa value to interpolate

real(kind=Float64), intent(in) :: yout

The ordinate value to interpolate

real(kind=Float64), intent(out) :: zout

Interpolant: z(xout,yout)

integer, intent(out), optional :: err

Error code

type(InterpolCoeffs2D), intent(in), optional :: coeffs

Precomputed Linear Interpolation Coefficients

public subroutine interpol2D_2D_arr(x, y, z, xout, yout, zout, err, coeffs)

Performs bilinear interpolation on a 2D grid of 2D arrays z(:,:,x,y)

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(:):: x

The abscissa values of z

real(kind=Float64), intent(in), dimension(:):: y

The ordinate values of z

real(kind=Float64), intent(in), dimension(:,:,:,:):: z

Values at the abscissa/ordinates: z(:,:,x,y)

real(kind=Float64), intent(in) :: xout

The abscissa value to interpolate

real(kind=Float64), intent(in) :: yout

The ordinate value to interpolate

real(kind=Float64), intent(out), dimension(:,:):: zout

Interpolant: z(:,:,xout,yout)

integer, intent(out), optional :: err

Error code

type(InterpolCoeffs2D), intent(in), optional :: coeffs

Precomputed Linear Interpolation Coefficients

public subroutine in_plasma(xyz, inp, machine_coords, coeffs, uvw_out)

Indicator subroutine to determine if a position is in a region where the plasma parameter and fields are valid/known

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: xyz

Position in beam coordinates

logical, intent(out) :: inp

Indicates whether plasma parameters and fields are valid/known

logical, intent(in), optional :: machine_coords

Indicates that xyz is in machine coordinates

type(InterpolCoeffs2D), intent(out), optional :: coeffs

Linear Interpolation coefficients used in calculation

real(kind=Float64), intent(out), optional dimension(3):: uvw_out

Position in machine coordinates

public subroutine get_plasma(plasma, pos, ind)

Gets plasma parameters at position pos or beam_grid indices ind

Arguments

Type IntentOptional AttributesName
type(LocalProfiles), intent(out) :: plasma

Plasma parameters at pos/ind

real(kind=Float64), intent(in), optional dimension(3):: pos

Position in beam grid coordinates

integer(kind=Int32), intent(in), optional dimension(3):: ind

beam_grid indices

public subroutine calc_perp_vectors(b, a, c)

Calculates normalized vectors that are perpendicular to b such that a x c = b_norm

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: b
real(kind=Float64), intent(out), dimension(3):: a
real(kind=Float64), intent(out), dimension(3):: c

public subroutine get_fields(fields, pos, ind, machine_coords)

Gets electro-magnetic fields at position pos or beam_grid indices ind

Arguments

Type IntentOptional AttributesName
type(LocalEMFields), intent(out) :: fields

Electro-magnetic fields at pos/ind

real(kind=Float64), intent(in), optional dimension(3):: pos

Position in beam grid coordinates

integer(kind=Int32), intent(in), optional dimension(3):: ind

beam_grid indices

logical, intent(in), optional :: machine_coords

Indicates that pos is machine coordinates

public subroutine get_distribution(fbeam, denf, pos, ind, coeffs)

Gets Guiding Center distribution at position pos or beam_grid indices ind

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(out), dimension(:,:):: fbeam

Guiding Center Fast-ion distribution at pos/ind: F(E,p)

real(kind=Float64), intent(out) :: denf

Guiding Center Fast-ion density at pos/ind [fast-ions/cm^3]

real(kind=Float64), intent(in), optional dimension(3):: pos

Position in beam grid coordinates

integer(kind=Int32), intent(in), optional dimension(3):: ind

beam_grid indices

type(InterpolCoeffs2D), intent(in), optional :: coeffs

Precomputed Linear Interpolation Coefficients

public subroutine get_ep_denf(energy, pitch, denf, pos, ind, coeffs)

Get guiding center fast-ion density at given energy and pitch at position pos or beam_grid indices ind

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: energy

Energy [keV]

real(kind=Float64), intent(in) :: pitch

Pitch

real(kind=Float64), intent(out) :: denf

Fast-ion density [fast-ions/(cm^3dEdp)]

real(kind=Float64), intent(in), optional dimension(3):: pos

Position in beam grid coordinates

integer(kind=Int32), intent(in), optional dimension(3):: ind

beam_grid indices

type(InterpolCoeffs2D), intent(in), optional :: coeffs

Precomputed Linear Interpolation Coefficients

public subroutine store_neutrals_cell(ind, neut_type, dens, store_iter)

Arguments

Type IntentOptional AttributesName
integer(kind=Int32), intent(in), dimension(3):: ind

beam_grid indices

integer, value:: neut_type

Neutral type

real(kind=Float64), intent(in), dimension(:):: dens

Neutral density [neutrals/cm^3]

logical, intent(in), optional :: store_iter

Store DCX/Halo iteration density in halo_iter_dens

public subroutine store_neutrals_track(tracks, ncell, neut_type)

Arguments

Type IntentOptional AttributesName
type(ParticleTrack), intent(in), dimension(:):: tracks

Neutral Particle Track

integer, intent(in) :: ncell

Number of cell in the particle track

integer, value:: neut_type

Neutral type

public subroutine store_births(ind, neut_type, dflux)

Store birth particles/density in birth

Arguments

Type IntentOptional AttributesName
integer(kind=Int32), intent(in), dimension(3):: ind

beam_grid indices

integer(kind=Int32), intent(in) :: neut_type

Neutral type

real(kind=Float64), intent(in) :: dflux

Deposited flux

public subroutine store_npa(det, ri, rf, vn, flux, orbit_class)

Store NPA particles in npa

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: det

Detector/Channel Number

real(kind=Float64), intent(in), dimension(3):: ri

Birth position in beam coordinates [cm]

real(kind=Float64), intent(in), dimension(3):: rf

Detector position in beam coordinates [cm]

real(kind=Float64), intent(in), dimension(3):: vn

Particle velocity [cm/s]

real(kind=Float64), intent(in) :: flux

Neutral flux [neutrals/s]

integer, intent(in), optional :: orbit_class

Orbit class ID

public subroutine bb_cx_rates(denn, vi, vn, rates)

Get beam-beam neutralization/cx rates

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(nlevs):: denn

Neutral density [cm^-3]

real(kind=Float64), intent(in), dimension(3):: vi

Ion velocity [cm/s]

real(kind=Float64), intent(in), dimension(3):: vn

Neutral velocity [cm/s]

real(kind=Float64), intent(out), dimension(nlevs):: rates

Reaction rates [1/s]

public subroutine bt_cx_rates(plasma, denn, vi, i_type, rates)

Get beam-target neutralization/cx rates

Arguments

Type IntentOptional AttributesName
type(LocalProfiles), intent(in) :: plasma

Plasma parameters

real(kind=Float64), intent(in), dimension(nlevs):: denn

Neutral density [cm^-3]

real(kind=Float64), intent(in), dimension(3):: vi

Ion velocity [cm/s]

integer, intent(in) :: i_type

Ion type

real(kind=Float64), intent(out), dimension(nlevs):: rates

Reaction rates [1/s]

public subroutine get_neutron_rate(plasma, eb, rate)

Gets neutron rate for a beam with energy eb interacting with a target plasma

Arguments

Type IntentOptional AttributesName
type(LocalProfiles), intent(in) :: plasma

Plasma Paramters

real(kind=Float64), intent(in) :: eb

Beam energy [keV]

real(kind=Float64), intent(out) :: rate

Neutron reaction rate [1/s]

public subroutine get_beam_cx_rate(ind, pos, v_ion, i_type, types, prob)

Get probability of a thermal ion charge exchanging with types neutrals

Arguments

Type IntentOptional AttributesName
integer(kind=Int32), intent(in), dimension(3):: ind

beam_grid indices

real(kind=Float64), intent(in), dimension(3):: pos

Interaction position in beam grid coordinates

real(kind=Float64), intent(in), dimension(3):: v_ion

Ion velocity [cm/s]

integer, intent(in) :: i_type

Ion type

integer(kind=Int32), intent(in), dimension(:):: types

Neutral types

real(kind=Float64), intent(out), dimension(nlevs):: prob

Charge exchange rate/probability [1/s]

public subroutine get_rate_matrix(plasma, i_type, eb, rmat)

Gets rate matrix for use in colrad

Arguments

Type IntentOptional AttributesName
type(LocalProfiles), intent(in) :: plasma

Plasma parameters

integer, intent(in) :: i_type

Ion type

real(kind=Float64), intent(in) :: eb

Ion energy [keV]

real(kind=Float64), intent(out), dimension(nlevs,nlevs):: rmat

Rate matrix

public subroutine colrad(plasma, i_type, vn, dt, states, dens, photons)

Evolve density of states in time dt via collisional radiative model

Arguments

Type IntentOptional AttributesName
type(LocalProfiles), intent(in) :: plasma

Plasma parameters

integer, intent(in) :: i_type

Ion/Neutral type (beam,thermal)

real(kind=Float64), intent(in), dimension(:):: vn

Neutral velocitiy [cm/s]

real(kind=Float64), intent(in) :: dt

Time interval [s]

real(kind=Float64), intent(inout), dimension(:):: states

Density of states

real(kind=Float64), intent(out), dimension(nlevs):: dens

Density of neutrals

real(kind=Float64), intent(out) :: photons

Emitted photons(3->2)

public subroutine attenuate(ri, rf, vi, states, dstep_in)

Attenuate states along a trajectory

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: ri

Initial position in beam grid coordinates

real(kind=Float64), intent(in), dimension(3):: rf

Final position in beam grid coordinates

real(kind=Float64), intent(in), dimension(3):: vi

Initial velocity of neutral

real(kind=Float64), intent(inout), dimension(nlevs):: states

Density of states

real(kind=Float64), intent(in), optional :: dstep_in

Step length [cm]

public subroutine spectrum(vecp, vi, fields, sigma_pi, photons, dlength, lambda, intensity)

Calculates doppler shift and stark splitting

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: vecp

Vector directing towards optical head

real(kind=Float64), intent(in), dimension(3):: vi

Particle velocity

type(LocalEMFields), intent(in) :: fields

Electro-magnetic fields

real(kind=Float64), intent(in) :: sigma_pi

Sigma-pi ratio

real(kind=Float64), intent(in) :: photons

Photon density from colrad

real(kind=Float64), intent(in) :: dlength

LOS intersection length with beam_grid cell particle is in

real(kind=Float64), intent(out), dimension(n_stark):: lambda

Wavelengths [nm]

real(kind=Float64), intent(out), dimension(n_stark):: intensity

Spectra intensities [Ph/(s cm^2 starkline)]

public subroutine store_bes_photons(pos, vi, photons, neut_type)

Store BES photons in Spectra

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: pos

Position of neutral in beam grid coordinates

real(kind=Float64), intent(in), dimension(3):: vi

Velocitiy of neutral [cm/s]

real(kind=Float64), intent(in) :: photons

Photons from colrad [Ph/(s*cm^3)]

integer, intent(in) :: neut_type

Neutral type (full,half,third,halo)

public subroutine store_fida_photons(pos, vi, photons, orbit_class)

Store fida photons in Spectra

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: pos

Position of neutral in beam grid coordinates

real(kind=Float64), intent(in), dimension(3):: vi

Velocitiy of neutral [cm/s]

real(kind=Float64), intent(in) :: photons

Photons from colrad [Ph/(s*cm^3)]

integer, intent(in), optional :: orbit_class

Orbit class ID

public subroutine store_neutrons(rate, orbit_class)

Store neutron rate in neutron

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: rate

Neutron rate [neutrons/sec]

integer, intent(in), optional :: orbit_class

Orbit class ID

public subroutine store_fw_photons_at_chan(ichan, eind, pind, vp, vi, fields, dlength, sigma_pi, denf, photons)

Store FIDA weight photons in fweight for a specific channel

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: ichan

Channel index

integer, intent(in) :: eind

Energy index

integer, intent(in) :: pind

Pitch index

real(kind=Float64), intent(in), dimension(3):: vp

Vector pointing toward optical head

real(kind=Float64), intent(in), dimension(3):: vi

Velocity of neutral [cm/s]

type(LocalEMFields), intent(in) :: fields

Electro-magnetic fields

real(kind=Float64), intent(in) :: dlength

LOS intersection length with beam_grid cell particle is in

real(kind=Float64), intent(in) :: sigma_pi

Sigma-pi ratio for channel

real(kind=Float64), intent(in) :: denf

Fast-ion density [cm^-3]

real(kind=Float64), intent(in) :: photons

Photons from colrad [Ph/(s*cm^3)]

public subroutine store_fw_photons(eind, pind, pos, vi, denf, photons)

Store FIDA weight photons in fweight

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: eind

Energy index

integer, intent(in) :: pind

Pitch index

real(kind=Float64), intent(in), dimension(3):: pos

Position of neutral

real(kind=Float64), intent(in), dimension(3):: vi

Velocity of neutral [cm/s]

real(kind=Float64), intent(in) :: denf

Fast-ion density [cm^-3]

real(kind=Float64), intent(in) :: photons

Photons from colrad [Ph/(s*cm^3)]

public subroutine get_nlaunch(nr_markers, papprox, papprox_tot, nlaunch)

Sets the number of MC markers launched from each beam_grid cell

Arguments

Type IntentOptional AttributesName
integer(kind=Int64), intent(in) :: nr_markers

Approximate total number of markers to launch

real(kind=Float64), intent(in), dimension(:,:,:):: papprox

beam_grid cell weights

real(kind=Float64), intent(in) :: papprox_tot

Total cell weights

real(kind=Float64), intent(out), dimension(:,:,:):: nlaunch

Number of mc markers to launch for each cell: nlaunch(x,y,z)

public subroutine pitch_to_vec(pitch, gyroangle, fields, vi_norm)

Calculates velocity vector from pitch, gyroangle and fields

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: pitch

Pitch

real(kind=Float64), intent(in) :: gyroangle

Gyroangle [radians]

type(LocalEMFields), intent(in) :: fields

Electromagnetic fields

real(kind=Float64), intent(out), dimension(3):: vi_norm

Normalized velocity vector

public subroutine gyro_step(vi, fields, r_gyro)

Calculates gyro-step

Read more…

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(3):: vi

Ion velocity

type(LocalEMFields), intent(in) :: fields

Electro-magnetic fields

real(kind=Float64), intent(out), dimension(3):: r_gyro

Gyro-step Gyro-radius vector from particle position to guiding center

public subroutine gyro_correction(fields, energy, pitch, rp, vp, phi_in)

Calculates gyro correction for Guiding Center MC distribution calculation

Arguments

Type IntentOptional AttributesName
type(LocalEMFields), intent(in) :: fields

Electromagnetic fields at guiding center

real(kind=Float64), intent(in) :: energy

Energy of particle

real(kind=Float64), intent(in) :: pitch

Particle pitch w.r.t the magnetic field

real(kind=Float64), intent(out), dimension(3):: rp

Particle position

real(kind=Float64), intent(out), dimension(3):: vp

Particle velocity

real(kind=Float64), intent(in), optional :: phi_in

Gyro-angle

public subroutine mc_fastion(ind, fields, eb, ptch, denf)

Samples a Guiding Center Fast-ion distribution function at a given beam_grid index

Arguments

Type IntentOptional AttributesName
integer, intent(in), dimension(3):: ind

beam_grid index

type(LocalEMFields), intent(out) :: fields

Electromagnetic fields at the guiding center

real(kind=Float64), intent(out) :: eb

Energy of the fast ion

real(kind=Float64), intent(out) :: ptch

Pitch of the fast ion

real(kind=Float64), intent(out) :: denf

Fast-ion density at guiding center

public subroutine mc_halo(ind, vhalo, ri, plasma_in)

Sample thermal Maxwellian distribution at beam_grid indices ind

Arguments

Type IntentOptional AttributesName
integer, intent(in), dimension(3):: ind

beam_grid indices

real(kind=Float64), intent(out), dimension(3):: vhalo

Velocity [cm/s]

real(kind=Float64), intent(out), optional dimension(3):: ri

Position in beam_grid cell

type(LocalProfiles), intent(in), optional :: plasma_in

Plasma parameters

public subroutine mc_nbi(vnbi, efrac, rnbi, err)

Generates a neutral beam particle trajectory

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(out), dimension(3):: vnbi

Velocity [cm/s]

integer, intent(in) :: efrac

Beam neutral type (1,2,3)

real(kind=Float64), intent(out), dimension(3):: rnbi

Starting position on beam_grid

logical, intent(out) :: err

Error Code

public subroutine ndmc()

Calculates neutral beam deposition and spectra

Arguments

None

public subroutine bremsstrahlung()

Calculates bremsstrahlung

Arguments

None

public subroutine dcx()

Calculates Direct Charge Exchange (DCX) neutral density and spectra

Arguments

None

public subroutine halo()

Calculates halo neutral density and spectra

Arguments

None

public subroutine fida_f()

Calculate FIDA emission using a Fast-ion distribution function F(E,p,r,z)

Arguments

None

public subroutine fida_mc()

Calculate FIDA emission using a Monte Carlo Fast-ion distribution

Arguments

None

public subroutine npa_f()

Calculate NPA flux using a fast-ion distribution function F(E,p,r,z)

Arguments

None

public subroutine npa_mc()

Calculate NPA flux using a Monte Carlo fast-ion distribution

Arguments

None

public subroutine neutron_f()

Calculate neutron emission rate using a fast-ion distribution function F(E,p,r,z)

Arguments

None

public subroutine neutron_mc()

Calculate neutron flux using a Monte Carlo Fast-ion distribution

Arguments

None

public subroutine fida_weights_mc()

Calculates FIDA weights

Arguments

None

public subroutine fida_weights_los()

Calculates LOS averaged FIDA weights

Arguments

None

public subroutine npa_weights()

Calculates NPA weights

Arguments

None