|
XGCa
|
Module for variables and routines related to the logical sheath boundary condition. More...
Functions/Subroutines | |
| subroutine | init_sheath_mod (sheath_mode) |
| Initializes sheath module, i.e. memory allocation initialization of variables, and for sheath_mode==2 identification of vertex groups among the wall vertices. More... | |
| subroutine | distribute_leftover_charge (iref, iref_max, solved, use_global_balance, sheath_elost_en, loss_shift, num_solved, num_failed, q_avail, delta_q, num_rebalanced) |
| Distributes the leftover charge. More... | |
| subroutine | calculate_sheath_dep (sheath_dep, sheath_dep_inv) |
| This new version of the logical sheath adustment does three things: 1) compute use energy distribution of incoming electrons to find the new sheath potential and 2) loop over electrons and decide which ones penetrate the sheath. More... | |
| subroutine | find_sheath_pot (dep, loss_shift, iep0, loss, solved, sheath_dep, sheath_elost_en) |
| Finds sheath potential at which the electron loss is (approx.) equal to the ion loss plus a given shift (loss_shift) More... | |
| subroutine | sheath_adjust_mode1_fort () |
| Old sheath calculation (sheath_mode==1) Based on the electron and ion losses at the wall nodes a new value of the sheath potential is determined by adjusting the previous value by a certain amount --> feedback control loop. More... | |
Variables | |
| integer | sheath_nep =100 |
| number of energy grid bins for sheath calculation More... | |
| real(kind=8), dimension(:,:,:), allocatable, target | sheath_en_max |
| Parallel velo cutoff (normalized to sqrt(T_e)) for sheath potential calculation. More... | |
| real(kind=8), dimension(:), allocatable, target | sheath_te_ev |
| Initial wall temperature for normalization. More... | |
| integer | sheath_nphi |
| Number of poloidal planes for sheath potential. More... | |
| real(8), dimension(:,:), allocatable, target | sheath_pot |
| Sheath potential on the wall nodes. More... | |
| real(8), dimension(:,:), allocatable, target | sheath_lost |
| Total charge lost to the wall. More... | |
| real(8), dimension(:,:), allocatable, target | sheath_ilost |
| Ion charge lost to the wall. More... | |
| integer | sheath_ngroups |
| Number of vertex groups for sheath calculation. More... | |
| integer | sheath_group_maxlen |
| Max. number of members of a sheath vertex group. More... | |
| integer, dimension(:), allocatable | sheath_igroup |
| Index of group to which a wall vertex belongs. More... | |
| integer, dimension(:), allocatable | sheath_group_len |
| Number of vertices that belong to a group. More... | |
| integer, dimension(:,:), allocatable | sheath_group_idx |
| Indices of wall vertices that belong to a group. More... | |
Module for variables and routines related to the logical sheath boundary condition.
| subroutine sheath_module::calculate_sheath_dep | ( | real(8), dimension(grid_global%nwall,sheath_nphi) | sheath_dep, |
| real(8), dimension(grid_global%nwall,sheath_nphi) | sheath_dep_inv | ||
| ) |
This new version of the logical sheath adustment does three things: 1) compute use energy distribution of incoming electrons to find the new sheath potential and 2) loop over electrons and decide which ones penetrate the sheath.

| subroutine sheath_module::distribute_leftover_charge | ( | integer, intent(in), value | iref, |
| integer, intent(in), value | iref_max, | ||
| integer, dimension(grid_global%nwall,sheath_nphi), intent(in) | solved, | ||
| integer, intent(in), value | use_global_balance, | ||
| real (kind=8), dimension(0:sheath_nep,grid_global%nwall,sheath_nphi), intent(in) | sheath_elost_en, | ||
| real(8), dimension(grid_global%nwall,sheath_nphi), intent(inout) | loss_shift, | ||
| integer, dimension(grid_global%nwall,sheath_nphi), intent(inout) | num_solved, | ||
| integer, dimension(grid_global%nwall,sheath_nphi), intent(inout) | num_failed, | ||
| real(8), dimension(grid_global%nwall), intent(inout) | q_avail, | ||
| real(8), dimension(grid_global%nwall,sheath_nphi), intent(inout) | delta_q, | ||
| integer, intent(inout) | num_rebalanced | ||
| ) |
Distributes the leftover charge.
| [in] | grid | Unstructured mesh data, type(grid_type) |
| [in] | iref | Loop iteration |
| [in] | iref_max | Total number of iterations |
| [in] | solved | Whether this bin has been solved |
| [in,out] | loss_shift | The loss shift |
| [in,out] | num_solved | number of bins solved |
| [in,out] | num_failed | number of bins failed |
| [in,out] | q_avail | q_avail |
| [in,out] | delta_q | delta_q |
| [in,out] | num_rebalanced | num_rebalanced |

| subroutine sheath_module::find_sheath_pot | ( | real (kind=8), dimension(grid_global%nwall,sheath_nphi), intent(in) | dep, |
| real (kind=8), dimension(grid_global%nwall,sheath_nphi), intent(in) | loss_shift, | ||
| integer, dimension(grid_global%nwall,sheath_nphi), intent(inout) | iep0, | ||
| real (kind=8), dimension(grid_global%nwall,sheath_nphi), intent(out) | loss, | ||
| integer, dimension(grid_global%nwall,sheath_nphi), intent(out) | solved, | ||
| real (kind=8), dimension(grid_global%nwall,sheath_nphi), intent(in) | sheath_dep, | ||
| real (kind=8), dimension(0:sheath_nep,grid_global%nwall,sheath_nphi) | sheath_elost_en | ||
| ) |
Finds sheath potential at which the electron loss is (approx.) equal to the ion loss plus a given shift (loss_shift)
| [in] | grid | Unstructured mesh data, type(grid_type) |
| [in] | dep | Resolution of energy grid, real(8) |
| [in] | loss_shift | Shift in lost charge used for global loss balance, real(8) |
| [in,out] | iep0 | Used for refinement of sheath potential on fine energy grid, integer |
| [out] | loss | Approx. charge loss determined on energy grid, real(8) |
| [out] | solved | Whether local charge balance is successful on a node, logical |

| subroutine sheath_module::init_sheath_mod | ( | integer, intent(in), value | sheath_mode | ) |
Initializes sheath module, i.e. memory allocation initialization of variables, and for sheath_mode==2 identification of vertex groups among the wall vertices.
| [in] | grid | Unstructured mesh data, type(grid_type) |
| [in] | sheath_mode | Operating mode of the sheath calculation: 1) 1) feedback control adjustment 2) Instantaneous sheath adjustment based on energy distribution of incoming electrons |

| subroutine sheath_module::sheath_adjust_mode1_fort |
Old sheath calculation (sheath_mode==1) Based on the electron and ion losses at the wall nodes a new value of the sheath potential is determined by adjusting the previous value by a certain amount --> feedback control loop.
| [in] | grid | Unstructured mesh data, type(grid_type) |

| real (kind=8), dimension(:,:,:), allocatable, target sheath_module::sheath_en_max |
Parallel velo cutoff (normalized to sqrt(T_e)) for sheath potential calculation.
| integer, dimension(:,:), allocatable sheath_module::sheath_group_idx |
Indices of wall vertices that belong to a group.
| integer, dimension(:), allocatable sheath_module::sheath_group_len |
Number of vertices that belong to a group.
| integer sheath_module::sheath_group_maxlen |
Max. number of members of a sheath vertex group.
| integer, dimension(:), allocatable sheath_module::sheath_igroup |
Index of group to which a wall vertex belongs.
| real (8), dimension(:,:), allocatable, target sheath_module::sheath_ilost |
Ion charge lost to the wall.
| real (8), dimension(:,:), allocatable, target sheath_module::sheath_lost |
Total charge lost to the wall.
| integer sheath_module::sheath_nep =100 |
number of energy grid bins for sheath calculation
| integer sheath_module::sheath_ngroups |
Number of vertex groups for sheath calculation.
| integer sheath_module::sheath_nphi |
Number of poloidal planes for sheath potential.
| real (8), dimension(:,:), allocatable, target sheath_module::sheath_pot |
Sheath potential on the wall nodes.
| real (kind=8), dimension(:), allocatable, target sheath_module::sheath_te_ev |
Initial wall temperature for normalization.