XGC1
|
#include "t_coeff_mod_macro.h"
Functions/Subroutines | |
subroutine | neutral_totalf (istep, grid, sp, f0_f, f0_df0g, f0_node_cost, axisym_parallelization) |
Calculates kinetic source due to interaction of plasma with neutral species via a test particle Monte Carlo method. Presently only hydrogen atom ionization and charge exchange processes are considered. The neutral background, characterized by a Maxwellian density and temperature, is computed and updated by the Monte Carlo algorithm implemented in subroutine neutral_totalf_step. More... | |
subroutine | neutral_totalf_setup (grid) |
Sets up neutral_totalf variables, and birth locations. More... | |
subroutine | neutral_totalf_step (grid, mass) |
2D Monte Carlo neutral transport calculation to establish the background used in subroutine neutral_totalf. The neutral birth locations are now at the material boundaries, an improvement over past versions. More... | |
subroutine | maxwell_neut (tn, vr, vz, vangle) |
Return a particle with velocities sampled from a Maxwellian distribution. More... | |
subroutine | neut_combine_tallies |
Sums the neutral particle weights over all MPI processes, distributes to all compute nodes, then calculates the actual neutral fluid quantities (neutral density, temperature, velocity), also calculating relative standard deviation). More... | |
subroutine | startrz (iwall, iwall_next, wall_frac, wall_length, r, z, wall_angle, Lwalli) |
Determines the starting (R,Z) position to launch a MC neutral particle. More... | |
subroutine | clear_neutral_tallies |
subroutine | accumulate_tally (grid, r, z, vr, vz, weight, vparan) |
Adds neutral particle weight contribution, to eventually calculate neutral fluid quantities (density, temperature, velocity). Calls search_tr2 to determine triangle. More... | |
subroutine | neut_pos_test (grid, r, z, ptest, itr, p) |
Determine whether neutral particle at a point (R,Z) is in the simulation domain or not. More... | |
subroutine | plasma (grid, itr, p, dene_out, deni_out, Te_out, Ti_out, Vparai_out) |
Calculate the plasma density, temperature, and parallel velocity for a point in triangle itr using plasma moments from mom_module subroutine mom_module::calc_moments. More... | |
subroutine | advance_neutral_sample (weight, r, z, vr, vz, vangle, grid, mass) |
subroutine accumulate_tally | ( | type(grid_type) | grid, |
real (kind=8), intent(in) | r, | ||
real (kind=8), intent(in) | z, | ||
real (kind=8), intent(in) | vr, | ||
real (kind=8), intent(in) | vz, | ||
real (kind=8), intent(in) | weight, | ||
real (kind=8), intent(in) | vparan | ||
) |
Adds neutral particle weight contribution, to eventually calculate neutral fluid quantities (density, temperature, velocity). Calls search_tr2 to determine triangle.
grid | (in) real(8), Neutral temperature, units eV |
r | (in) real(8), Neutral particle R position |
z | (in) real(8), Neutral particle Z position |
vr | (in) real(8), Neutral velocity R-direction |
vz | (in) real(8), Neutral velocity Z-direction |
weight | (in) real(8), Neutral Monte-Carlo weight |
vparan | (in) real(8), Neutral parallel velocity |
subroutine advance_neutral_sample | ( | real (kind=8), intent(in) | weight, |
real (kind=8), intent(inout) | r, | ||
real (kind=8), intent(inout) | z, | ||
real (kind=8), intent(inout) | vr, | ||
real (kind=8), intent(inout) | vz, | ||
real (kind=8), intent(inout) | vangle, | ||
type(grid_type) | grid, | ||
real (kind=8), intent(in) | mass | ||
) |
subroutine clear_neutral_tallies | ( | ) |
subroutine neutral_totalf_step::maxwell_neut | ( | real(kind=8), intent(in) | tn, |
real(kind=8), intent(out) | vr, | ||
real(kind=8), intent(out) | vz, | ||
real(kind=8), intent(out) | vangle | ||
) |
Return a particle with velocities sampled from a Maxwellian distribution.
tn | (in) real(8), Neutral temperature, units eV |
vr | (out) real(8), Neutral velocity R-direction |
vz | (out) real(8), Neutral velocity Z-direction |
vangle | (out) real(8), Neutral velocity angle |
subroutine neutral_totalf_step::neut_combine_tallies | ( | ) |
Sums the neutral particle weights over all MPI processes, distributes to all compute nodes, then calculates the actual neutral fluid quantities (neutral density, temperature, velocity), also calculating relative standard deviation).
subroutine neut_pos_test | ( | type(grid_type) | grid, |
real (kind=8), intent(in) | r, | ||
real (kind=8), intent(in) | z, | ||
integer, intent(out) | ptest, | ||
integer, intent(out) | itr, | ||
real (kind=8), dimension(3), intent(out) | p | ||
) |
Determine whether neutral particle at a point (R,Z) is in the simulation domain or not.
grid | (in) type(grid_type), XGC grid information |
r | (in) real(8), Neutral particle R position |
z | (in) real(8), Neutral particle Z position |
ptest | (out) integer, 1 if neutral in simulation domain, -1 otherwise |
itr | (out) integer, Triangle containing point |
p | (out) real(8) [3], Triangle weights for interpolating node quantities onto point |
subroutine neutral_totalf | ( | integer, intent(in) | istep, |
type(grid_type) | grid, | ||
type(species_type) | sp, | ||
real (8), dimension(-f0_nvp:f0_nvp, f0_inode1:f0_inode2, f0_imu1:f0_imu2, ptl_isp:ptl_nsp), intent(in) | f0_f, | ||
real (8), dimension(-f0_nvp:f0_nvp, f0_inode1:f0_inode2, f0_imu1:f0_imu2, ptl_isp:ptl_nsp), intent(inout) | f0_df0g, | ||
real (8), dimension(f0_inode1:f0_inode2) | f0_node_cost, | ||
logical, intent(in) | axisym_parallelization | ||
) |
Calculates kinetic source due to interaction of plasma with neutral species via a test particle Monte Carlo method. Presently only hydrogen atom ionization and charge exchange processes are considered. The neutral background, characterized by a Maxwellian density and temperature, is computed and updated by the Monte Carlo algorithm implemented in subroutine neutral_totalf_step.
istep | (in) integer, Global XGC step |
grid | (in) type(grid_type), XGC grid information |
sp(in) | type(species_type), XGC species information |
subroutine neutral_totalf_setup | ( | type(grid_type) | grid | ) |
Sets up neutral_totalf variables, and birth locations.
grid | (in) type(grid_type), XGC grid information |
subroutine neutral_totalf_step | ( | type(grid_type) | grid, |
real (kind=8), intent(in) | mass | ||
) |
2D Monte Carlo neutral transport calculation to establish the background used in subroutine neutral_totalf. The neutral birth locations are now at the material boundaries, an improvement over past versions.
grid | (in) type(grid_type), XGC grid information |
mass | (in) real(8), neutral mass |
subroutine plasma | ( | type(grid_type) | grid, |
integer, intent(in) | itr, | ||
real (kind=8), dimension(3), intent(in) | p, | ||
real (kind=8), intent(out) | dene_out, | ||
real (kind=8), intent(out) | deni_out, | ||
real (kind=8), intent(out) | Te_out, | ||
real (kind=8), intent(out) | Ti_out, | ||
real (kind=8), intent(out) | Vparai_out | ||
) |
Calculate the plasma density, temperature, and parallel velocity for a point in triangle itr using plasma moments from mom_module subroutine mom_module::calc_moments.
grid | (in) type(grid_type), XGC grid information |
itr | (in) integer, Triangle for point to calculate plasma quantities |
p | (in) real(8) [3], Triangle weights for interpolating node quantities onto point |
dene_out | (out) real(8), Electron density at point |
deni_out | (out) real(8), Ion density at point |
Te_out | (out) real(8), Electron temperature at point |
Ti_out | (out) real(8), Ion temperature at point |
Vparai_out | (out) real(8), Ion parallel velocity at point |
subroutine neutral_totalf_step::startrz | ( | integer, intent(in) | iwall, |
integer, intent(in) | iwall_next, | ||
real (kind=8), intent(in) | wall_frac, | ||
real (kind=8), dimension(neu_nbirth), intent(in) | wall_length, | ||
real (kind=8), intent(out) | r, | ||
real (kind=8), intent(out) | z, | ||
real (kind=8), intent(out) | wall_angle, | ||
real (kind=8), intent(out) | Lwalli | ||
) |
Determines the starting (R,Z) position to launch a MC neutral particle.
iwall | (in) integer, Start index of wall node bouding position to launch from |
iwall_next | (in) integer, End index of wall node bounding position to launch from |
wall_frac | (in) real(8), Fraction of the distance between iwall and iwall_next to launch neutral at |
wall_length | (in) real(8) [neu_nbirth], Array of wall lengths for each node point (taken as sum of half distance to each neighboring wall node) |
r | (out) real(8), Neutral particle R position of launch |
z | (out) real(8), Neutral particle Z position of launch |
wall_angle | (out) real(8), Neutral particle launch angle w.r.t. wall |
Lwalli | (out) real(8), Wall length for (R,Z) position |