XGCa
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Data Types | Public Member Functions | Public Attributes | List of all members
grid_class Module Reference
Collaboration diagram for grid_class:
Collaboration graph
[legend]

Data Types

type  decomp_type
 
interface  gradparx2_cpp
 
type  grid_type
 

Public Member Functions

subroutine init_fortran_grid (grid_nnode_in, grid_ntriangle_in, grid_nwall_in, grid_x_in, grid_rgn_in, grid_nd_in, grid_mapping_in, grid_wall_nodes_in, grid_node_to_wall_in)
 
subroutine init_fortran_grid_surf (n_xpt, i_x1, i_x2, nsurf1, nsurf2, nsurf3a, nsurf3b, isurf_sep1, isurf_sep2, npsi_surf, surf_maxlen, num_non_aligned, surf_len, surf_idx, non_aligned_vert, non_aligned_nsurf, non_aligned_surf_idx)
 
subroutine init_grid (grid)
 Initializes the grid data structure from input files containing i) The (R,Z) coordinates of all vertices including whether they are wall vertices (.node file), ii) the connectivity of the vertices defining the triangles of the unstructured mesh (.ele file), iii) the flux-surface connectivity, i.e. which vertices are on the same flux-surface, and for vertices not aligned on flux-surfaces between which flux-surfaces those vertices lie. More...
 
subroutine init_triangle (grid)
 
subroutine t_coeff (grid, itr, x, p)
 
subroutine t_coeff_mod (grid, xy, psiin, itr, p)
 
subroutine search_tr (grid, x, itr, p)
 
subroutine guess (grid, x, init)
 
subroutine init_guess_table (grid)
 
subroutine search_tr_with_guess (grid, x, init, itr, p, count)
 
subroutine calc_gen_theta_psi (grid)
 Calculates the average psi on the complete flux-surfaces and the generalized poloidal angle (=normalized poloidal arc length) More...
 
subroutine search_ptl_1d_init (grid)
 Initializes a guess list for the 1D search on a non-uniform psi-grid. Uses psi00 grid. (Assumes that flux-surfaces on the non-uniform grid are in ascending order) More...
 
subroutine search_ptl_1d (grid, psi, wp, ip)
 1D binary search on a non-uniform psi grid with initial guess based on uniform psi-grid (Assumes that flux-surfaces on the non-uniform grid are in ascending order) More...
 
subroutine calc_mag_drift_flux_avg (grid)
 Computes the flux-surface average of 1/B^3 Bxgrad(B) and 1/B^2 curl(B). Those are zero analytically, but not numerically. This numerical cancellation error affects the calculation of radial fluxes a lot –> can be removed from the calculation. More...
 
subroutine init_wall (grid)
 
subroutine sendright_recvleft (grid, input)
 Sends the iphi=1 plane to the toroidal domain to the right and receives the iphi=0 plane from the toroidal domain to the left. It is used for exchanging field information that is known only on the iphi=1 plane on all processes. More...
 
subroutine sendright_recvleft_sendleft_recvright (grid, input, output)
 Sends the (R,Z) plane data to the toroidal domain to the right and receives the (R,Z) plane data from the toroidal domain to the left and vice versa to collect the data from 3 adjacent plances It is typically used for exchanging field information that is known only on the iphi=1 plane on all processes. More...
 
subroutine send_recv_potential (dpot, recvr, sendl, nnode)
 Collect (R,Z) plane data on 4 adjacent planes iphi=(-1:2) iphi=1 is assumed to be local data on the MPI rank. More...
 
subroutine gradparx (grid, tr, p, dx, bd, input, output)
 Calculates the parallel gradient like for the right plane (iphi=1). The input is supposed to contain the data of four adjacent planes (0:2) More...
 
subroutine gradparx2 (grid, tr, p, dx, bd, input, output)
 Calculates the parallel gradient like GradParX, but for the left and the right plane at the same time. The input is supposed to contain the data of four adjacent planes (-1:2 !!!) The routine calculates the parallel gradient for both planes in the sector (0 and 1) More...
 

Public Attributes

integer, parameter grid_rgn_wall = 100
 
real(8), parameter grid_search_tr2_psi_null = -1D10
 
integer, parameter grid_min_surf_size = 5
 
type(grid_type), target grid_global
 

Member Function/Subroutine Documentation

subroutine grid_class::calc_gen_theta_psi ( type(grid_type), intent(inout)  grid)

Calculates the average psi on the complete flux-surfaces and the generalized poloidal angle (=normalized poloidal arc length)

Parameters
grid(inout) type(grid_type), grid information

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine grid_class::calc_mag_drift_flux_avg ( type(grid_type), intent(inout)  grid)

Computes the flux-surface average of 1/B^3 Bxgrad(B) and 1/B^2 curl(B). Those are zero analytically, but not numerically. This numerical cancellation error affects the calculation of radial fluxes a lot –> can be removed from the calculation.

Parameters
[in,out]gridgrid data, type(grid_type)

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine grid_class::gradparx ( type(grid_type), intent(in)  grid,
integer, dimension(grid%nnode,0:1), intent(in)  tr,
real (kind=8), dimension(3,grid%nnode,0:1), intent(in)  p,
real (kind=8), dimension(grid%nnode,0:1), intent(in)  dx,
type(boundary2_type), intent(in)  bd,
real (kind=8), dimension(grid%nnode,3), intent(in)  input,
real (kind=8), dimension(grid%nnode), intent(out)  output 
)

Calculates the parallel gradient like for the right plane (iphi=1). The input is supposed to contain the data of four adjacent planes (0:2)

Parameters
[in]gridXGC grid data structure, type(grid_type)
[in]trparallel triangle connectivity, integer
[in]pinterpolation weights corresponding to tr, real(8)
[in]dxdistances along the field line, real(8)
[in]bdboundary conditions, type(boundary2_type)
[in]inputinput data for all planes from 0 to 2, real(8)
[out]outputparallel derivative on planes 0 and 1, real(8)

Here is the call graph for this function:

subroutine grid_class::gradparx2 ( type(grid_type), intent(in)  grid,
integer, dimension(grid%nnode,0:1), intent(in)  tr,
real (kind=8), dimension(3,grid%nnode,0:1), intent(in)  p,
real (kind=8), dimension(grid%nnode,0:1), intent(in)  dx,
type(boundary2_type), intent(in)  bd,
real (kind=8), dimension(grid%nnode,-1:2), intent(in)  input,
real (kind=8), dimension(grid%nnode,0:1), intent(out)  output 
)

Calculates the parallel gradient like GradParX, but for the left and the right plane at the same time. The input is supposed to contain the data of four adjacent planes (-1:2 !!!) The routine calculates the parallel gradient for both planes in the sector (0 and 1)

Parameters
[in]gridXGC grid data structure, type(grid_type)
[in]trparallel triangle connectivity, integer
[in]pinterpolation weights corresponding to tr, real(8)
[in]dxdistances along the field line, real(8)
[in]bdboundary conditions, type(boundary2_type)
[in]inputinput data for all planes from -1 to 2, real(8)
[out]outputparallel derivative on planes 0 and 1, real(8)

Here is the call graph for this function:

subroutine grid_class::guess ( type(grid_type), intent(in)  grid,
real (kind=8), dimension(2), intent(in)  x,
integer, intent(out)  init 
)

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine grid_class::init_fortran_grid ( integer(c_int), intent(in)  grid_nnode_in,
integer(c_int), intent(in)  grid_ntriangle_in,
integer(c_int), intent(in)  grid_nwall_in,
real (kind=8), dimension(2,grid_nnode_in), target  grid_x_in,
integer(c_int), dimension(grid_nnode_in), target  grid_rgn_in,
integer(c_int), dimension(3,grid_ntriangle_in), target  grid_nd_in,
real (kind=8), dimension(2,3,grid_ntriangle_in), target  grid_mapping_in,
integer(c_int), dimension(grid_nwall_in), target  grid_wall_nodes_in,
integer(c_int), dimension(grid_nnode_in), target  grid_node_to_wall_in 
)
subroutine grid_class::init_fortran_grid_surf ( integer(c_int), intent(in)  n_xpt,
integer(c_int), intent(in)  i_x1,
integer(c_int), intent(in)  i_x2,
integer(c_int), intent(in)  nsurf1,
integer(c_int), intent(in)  nsurf2,
integer(c_int), intent(in)  nsurf3a,
integer(c_int), intent(in)  nsurf3b,
integer(c_int), intent(in)  isurf_sep1,
integer(c_int), intent(in)  isurf_sep2,
integer(c_int), intent(in)  npsi_surf,
integer(c_int), intent(in)  surf_maxlen,
integer(c_int), intent(in)  num_non_aligned,
integer(c_int), dimension(npsi_surf), target  surf_len,
integer(c_int), dimension(surf_maxlen, npsi_surf), target  surf_idx,
integer(c_int), dimension(num_non_aligned), target  non_aligned_vert,
integer(c_int), dimension(num_non_aligned), target  non_aligned_nsurf,
integer(c_int), dimension(4,num_non_aligned), target  non_aligned_surf_idx 
)
subroutine grid_class::init_grid ( type(grid_type), intent(inout)  grid)

Initializes the grid data structure from input files containing i) The (R,Z) coordinates of all vertices including whether they are wall vertices (.node file), ii) the connectivity of the vertices defining the triangles of the unstructured mesh (.ele file), iii) the flux-surface connectivity, i.e. which vertices are on the same flux-surface, and for vertices not aligned on flux-surfaces between which flux-surfaces those vertices lie.

Parameters
[in,out]gridmemory for grid data, type(grid_type)
[in]filename3File containing the flux-surface connectivity

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine grid_class::init_guess_table ( type(grid_type grid)

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine grid_class::init_triangle ( type(grid_type grid)

Here is the caller graph for this function:

subroutine grid_class::init_wall ( type(grid_type grid)

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine grid_class::search_ptl_1d ( type(grid_type), intent(in)  grid,
real (kind=8), intent(in)  psi,
real (kind=8), intent(out)  wp,
integer, intent(out)  ip 
)

1D binary search on a non-uniform psi grid with initial guess based on uniform psi-grid (Assumes that flux-surfaces on the non-uniform grid are in ascending order)

Parameters
[in]grid(grid_type), grid-data
[in]psi(real (8)), psi-value to be located
[out]wp(real (8)), linear interpolation weight for the lower bound
[out]ip(integer), index of the lower boundary surface

Here is the caller graph for this function:

subroutine grid_class::search_ptl_1d_init ( type(grid_type), intent(inout)  grid)

Initializes a guess list for the 1D search on a non-uniform psi-grid. Uses psi00 grid. (Assumes that flux-surfaces on the non-uniform grid are in ascending order)

Parameters
[in,out]grid(grid_type), grid data

Here is the caller graph for this function:

subroutine grid_class::search_tr ( type (grid_type), intent(in)  grid,
real (kind=8), dimension(2), intent(in)  x,
integer, intent(out)  itr,
real (kind=8), dimension(3), intent(out)  p 
)

Here is the call graph for this function:

subroutine grid_class::search_tr_with_guess ( type(grid_type), target  grid,
real (kind=8), dimension(2), intent(in)  x,
integer, intent(in)  init,
integer, intent(out)  itr,
real (kind=8), dimension(3), intent(out)  p,
integer  count 
)

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine grid_class::send_recv_potential ( real (kind=8), dimension(nnode,-1:2)  dpot,
real (kind=8), dimension(nnode)  recvr,
real (kind=8), dimension(nnode)  sendl,
integer, intent(in)  nnode 
)

Collect (R,Z) plane data on 4 adjacent planes iphi=(-1:2) iphi=1 is assumed to be local data on the MPI rank.

Here is the caller graph for this function:

subroutine grid_class::sendright_recvleft ( type(grid_type), intent(in)  grid,
real (kind=8), dimension(grid%nnode,0:1), intent(inout)  input 
)

Sends the iphi=1 plane to the toroidal domain to the right and receives the iphi=0 plane from the toroidal domain to the left. It is used for exchanging field information that is known only on the iphi=1 plane on all processes.

Parameters
[in]gridXGC grid data structure, type(grid_type)
[in,out]inputData that is to be exchanged, real(8)
subroutine grid_class::sendright_recvleft_sendleft_recvright ( type(grid_type), intent(in)  grid,
real (kind=8), dimension(grid%nnode), intent(in)  input,
real (kind=8), dimension(grid%nnode,3), intent(out)  output 
)

Sends the (R,Z) plane data to the toroidal domain to the right and receives the (R,Z) plane data from the toroidal domain to the left and vice versa to collect the data from 3 adjacent plances It is typically used for exchanging field information that is known only on the iphi=1 plane on all processes.

Parameters
[in]gridgrid information,type(grid_type)
[in]inputThe scalar that is to be exchanged, real(8)
[out]outputThe collected data of 3 planes, real(8)
subroutine grid_class::t_coeff ( type(grid_type), intent(in)  grid,
integer, intent(in)  itr,
real (kind=8), dimension(2), intent(in)  x,
real (kind=8), dimension(3), intent(out)  p 
)

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine grid_class::t_coeff_mod ( type (grid_type), intent(in)  grid,
real (kind=8), dimension(2), intent(in)  xy,
real (kind=8), intent(in)  psiin,
integer, intent(in)  itr,
real (kind=8), dimension(3), intent(out)  p 
)

Here is the caller graph for this function:

Member Data Documentation

type(grid_type), target grid_class::grid_global
integer, parameter grid_class::grid_min_surf_size = 5
integer, parameter grid_class::grid_rgn_wall = 100
real (8), parameter grid_class::grid_search_tr2_psi_null = -1D10

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