gs2_io_grids Module

A module for writing geometry to netcdf.


Uses


Contents


Variables

Type Visibility Attributes Name Initial
integer(kind=kind_nf), private :: ncid

Subroutines

public subroutine nc_grid_file_open(filename, mode)

Open the netCDF file If the file already exists, remove it.

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename
character(len=*), intent(in) :: mode

public subroutine nc_grid_file_close()

Close the netCDF file

Arguments

None

public subroutine nc_get_grid_sizes(ntheta, ntgrid, nperiod)

Read sizes of geometry data from the input netCDF file

Arguments

Type IntentOptional Attributes Name
integer, intent(out), optional :: ntheta
integer, intent(out), optional :: ntgrid
integer, intent(out), optional :: nperiod

public subroutine nc_write_grid_sizes(ntheta, ntgrid, nperiod)

Writes sizes of geometry data to the output netCDF file

Arguments

Type IntentOptional Attributes Name
integer, intent(in), optional :: ntheta
integer, intent(in), optional :: ntgrid
integer, intent(in), optional :: nperiod

public subroutine nc_get_lambda_grid_size(nlambda)

Read sizes of lambda grid data from the input netCDF file

Arguments

Type IntentOptional Attributes Name
integer, intent(out) :: nlambda

public subroutine nc_get_grid_scalars(shat, drhodpsi, kxfac, qval, rmaj, B_T, aminor, grhoavg, surfarea)

Read scalar geometry data from the input netCDF file

Arguments

Type IntentOptional Attributes Name
real, intent(out), optional :: shat
real, intent(out), optional :: drhodpsi
real, intent(out), optional :: kxfac
real, intent(out), optional :: qval
real, intent(out), optional :: rmaj
real, intent(out), optional :: B_T
real, intent(out), optional :: aminor
real, intent(out), optional :: grhoavg
real, intent(out), optional :: surfarea

public subroutine nc_write_grid_scalars(shat, drhodpsi, kxfac, qval, rmaj, B_T, aminor, grhoavg, surfarea)

Write scalar geometry data to the output netCDF file

Arguments

Type IntentOptional Attributes Name
real, intent(in), optional :: shat
real, intent(in), optional :: drhodpsi
real, intent(in), optional :: kxfac
real, intent(in), optional :: qval
real, intent(in), optional :: rmaj
real, intent(in), optional :: B_T
real, intent(in), optional :: aminor
real, intent(in), optional :: grhoavg
real, intent(in), optional :: surfarea

public subroutine nc_get_grids(ntgrid, bmag, gradpar, gbdrift, gbdrift0, cvdrift, cvdrift0, gds2, gds21, gds22, grho, theta, cdrift, cdrift0, Rplot, Zplot, Rprime, Zprime, aplot, aprime, Bpol, jacob)

Read geometry data from the input netCDF file The no_end_point flag in the input netCDF file tells if the grid data have the duplicate end point data at ntgrid. The default no_end_point behavior: The grid file generated by this module should always store the no_end_point flag. The grid file came from outside (GX) may not have this flag. (GX grid file does not have the end point data.) read_and_pad handles the end point data

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: ntgrid
real, intent(out), dimension(-ntgrid:ntgrid) :: bmag
real, intent(out), dimension(-ntgrid:ntgrid) :: gradpar
real, intent(out), dimension(-ntgrid:ntgrid) :: gbdrift
real, intent(out), dimension(-ntgrid:ntgrid) :: gbdrift0
real, intent(out), dimension(-ntgrid:ntgrid) :: cvdrift
real, intent(out), dimension(-ntgrid:ntgrid) :: cvdrift0
real, intent(out), dimension(-ntgrid:ntgrid) :: gds2
real, intent(out), dimension(-ntgrid:ntgrid) :: gds21
real, intent(out), dimension(-ntgrid:ntgrid) :: gds22
real, intent(out), dimension(-ntgrid:ntgrid) :: grho
real, intent(out), dimension(-ntgrid:ntgrid) :: theta
real, intent(out), optional, dimension(-ntgrid:ntgrid) :: cdrift
real, intent(out), optional, dimension(-ntgrid:ntgrid) :: cdrift0
real, intent(out), optional, dimension(-ntgrid:ntgrid) :: Rplot
real, intent(out), optional, dimension(-ntgrid:ntgrid) :: Zplot
real, intent(out), optional, dimension(-ntgrid:ntgrid) :: Rprime
real, intent(out), optional, dimension(-ntgrid:ntgrid) :: Zprime
real, intent(out), optional, dimension(-ntgrid:ntgrid) :: aplot
real, intent(out), optional, dimension(-ntgrid:ntgrid) :: aprime
real, intent(out), optional, dimension(-ntgrid:ntgrid) :: Bpol
real, intent(out), optional, dimension(-ntgrid:ntgrid) :: jacob

public subroutine nc_write_grids(ntgrid, bmag, gradpar, gbdrift, gbdrift0, cvdrift, cvdrift0, gds2, gds21, gds22, grho, theta, cdrift, cdrift0, Rplot, Rprime, Zplot, Zprime, aplot, aprime, Bpol, jacob, no_end_point_in)

Write geometry data to the output netCDF file Get the no_end_point flag from the input and store it in the output netCDF file. If no_end_point=.true., data at ntgrid are not stored. Data at ntgrid may be generated from the data at -ntgrid with some symmetry assumption. This option is meant to provide compatibility with GX.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: ntgrid
real, intent(in), dimension(-ntgrid:ntgrid) :: bmag
real, intent(in), dimension(-ntgrid:ntgrid) :: gradpar
real, intent(in), dimension(-ntgrid:ntgrid) :: gbdrift
real, intent(in), dimension(-ntgrid:ntgrid) :: gbdrift0
real, intent(in), dimension(-ntgrid:ntgrid) :: cvdrift
real, intent(in), dimension(-ntgrid:ntgrid) :: cvdrift0
real, intent(in), dimension(-ntgrid:ntgrid) :: gds2
real, intent(in), dimension(-ntgrid:ntgrid) :: gds21
real, intent(in), dimension(-ntgrid:ntgrid) :: gds22
real, intent(in), dimension(-ntgrid:ntgrid) :: grho
real, intent(in), dimension(-ntgrid:ntgrid) :: theta
real, intent(in), optional, dimension(-ntgrid:ntgrid) :: cdrift
real, intent(in), optional, dimension(-ntgrid:ntgrid) :: cdrift0
real, intent(in), optional, dimension(-ntgrid:ntgrid) :: Rplot
real, intent(in), optional, dimension(-ntgrid:ntgrid) :: Rprime
real, intent(in), optional, dimension(-ntgrid:ntgrid) :: Zplot
real, intent(in), optional, dimension(-ntgrid:ntgrid) :: Zprime
real, intent(in), optional, dimension(-ntgrid:ntgrid) :: aplot
real, intent(in), optional, dimension(-ntgrid:ntgrid) :: aprime
real, intent(in), optional, dimension(-ntgrid:ntgrid) :: Bpol
real, intent(in), optional, dimension(-ntgrid:ntgrid) :: jacob
logical, intent(in), optional :: no_end_point_in

public subroutine nc_get_lambda_grid(nlambda, lambda)

Read lambda grid data from the input netCDF file

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nlambda
real, intent(out), dimension(nlambda) :: lambda

public subroutine nc_write_lambda_grid(nlambda, lambda)

Write lambda grid data to the output netCDF file

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nlambda
real, intent(in), dimension(nlambda) :: lambda

private subroutine read_and_pad(ncid, array, name, no_end_point)

Helper subroutine to appropriately get geometry data depending on the no_end_point flag. Read geometry data from the input netCDF file with an appropriate size and pad ntgrid data if no_end_point=.true. The array corresponding to the dummy argument array in the upstream subroutine may be non-present optional argument, whose size will be 1 /= 2*ntgrid+1.

Arguments

Type IntentOptional Attributes Name
integer(kind=kind_nf), intent(in) :: ncid
real, intent(inout), dimension(:) :: array
character(len=*), intent(in) :: name
logical, intent(in) :: no_end_point

private subroutine ntgrid_padding(array)

Generate data at ntgrid point using some symmetry assumption

Arguments

Type IntentOptional Attributes Name
real, intent(inout), dimension (:) :: array

If grid file does not have the end point data, they are generated using information on the other end by assuming (anti-)symmetry w.r.t theta=0. If symmetric, array(ntgrid) = array(-ntgrid) If anti-symmetric, array(ntgrid) = -array(-ntgrid)

Read more…