XGC1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Public Attributes | List of all members
sheath_module Module Reference

Module for variables and routines related to the logical sheath boundary condition. More...

Public Member Functions

subroutine init_sheath_mod (grid, 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 finalize_sheath_mod (sheath_mode)
 Finalize sheath module: deallocation of sheath variables. More...
 
subroutine distribute_leftover_charge (iref, iref_max, solved, loss_shift, q_minmax, 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 sheath_adjust_mode2_p2 (ipc)
 
subroutine find_sheath_pot (dep, loss_shift, iep0, loss, q_minmax, 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_diag_output (grid, i, ipc, sheath_dep, sheath_elost_en)
 Outputs some diagnostic data that helps to check the correctness of the algorithm that determines the sheath potential when sheath_mode==2. More...
 
subroutine sheath_adjust_mode1 ()
 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...
 

Public Attributes

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...
 

Detailed Description

Module for variables and routines related to the logical sheath boundary condition.

Member Function/Subroutine Documentation

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.

Here is the call graph for this function:

subroutine sheath_module::distribute_leftover_charge ( integer  iref,
integer  iref_max,
integer, dimension(grid_global%nwall,sheath_nphi), intent(in)  solved,
real(8), dimension(grid_global%nwall,sheath_nphi), intent(inout)  loss_shift,
real(8), dimension(grid_global%nwall,sheath_nphi,2), intent(inout)  q_minmax,
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.

Parameters
[in]gridUnstructured mesh data, type(grid_type)
[in]irefLoop iteration
[in]iref_maxTotal number of iterations
[in]solvedWhether this bin has been solved
[in,out]loss_shiftThe loss shift
[in,out]q_minmaxq_minmax
[in,out]num_solvednumber of bins solved
[in,out]num_failednumber of bins failed
[in,out]q_availq_avail
[in,out]delta_qdelta_q
[in,out]num_rebalancednum_rebalanced

Here is the call graph for this function:

subroutine sheath_module::finalize_sheath_mod ( integer, intent(in)  sheath_mode)

Finalize sheath module: deallocation of sheath variables.

Parameters
[in]sheath_modeOperating mode of the sheath calculation
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,
real (kind=8), dimension(grid_global%nwall,sheath_nphi,2), intent(out)  q_minmax,
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)

Parameters
[in]gridUnstructured mesh data, type(grid_type)
[in]depResolution of energy grid, real(8)
[in]loss_shiftShift in lost charge used for global loss balance, real(8)
[in,out]iep0Used for refinement of sheath potential on fine energy grid, integer
[out]lossApprox. charge loss determined on energy grid, real(8)
[out]q_minmaxMin. and max. of possible electron loss, real(8)
[out]solvedWhether local charge balance is successful on a node, logical

Here is the call graph for this function:

subroutine sheath_module::init_sheath_mod ( type(grid_type), intent(in)  grid,
integer, intent(in)  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.

Parameters
[in]gridUnstructured mesh data, type(grid_type)
[in]sheath_modeOperating mode of the sheath calculation: 1) 1) feedback control adjustment 2) Instantaneous sheath adjustment based on energy distribution of incoming electrons

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine sheath_module::sheath_adjust_mode1 ( )

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.

Parameters
[in]gridUnstructured mesh data, type(grid_type)

Here is the call graph for this function:

subroutine sheath_module::sheath_adjust_mode2_p2 ( integer  ipc)

Here is the call graph for this function:

subroutine sheath_module::sheath_diag_output ( type(grid_type), intent(in)  grid,
integer, intent(in)  i,
integer, intent(in)  ipc,
real(8), dimension(grid%nwall,sheath_nphi sheath_dep,
real (kind=8), dimension(0:sheath_nep,grid_global%nwall,sheath_nphi,sml_nthreads)  sheath_elost_en 
)

Outputs some diagnostic data that helps to check the correctness of the algorithm that determines the sheath potential when sheath_mode==2.

Parameters
[in]gridUnstructured mesh data, type(grid_type)
[in]iFlag indicating whether data should be written, integer
[in]ipcRunge-Kutta step

Member Data Documentation

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.


The documentation for this module was generated from the following file: