Python Library Documentation: prock in module __main__ object

 
class prock(builtins.object)
     Methods defined here:
deprecated_make_cold_spectrum(vn_z, mvn_mod, dim_freq, spectrum, im=0) -> 'int'
deprecated_make_cold_spectrum(vn_z: real[:], mvn_mod: real[:,_:,_:], dim_freq: int, spectrum: real[:,_:], im: int = 0) -> int
 
Parameters:
    vn_z:     Shape (dim_z,)
    mvn_mod:  Shape (dim_m, DIM_MOD, dim_z)
    dim_freq (int):   Number of spectral points
    spectrum: Shape (2, dim_freq) array whose vectors
                      [0, :] are the frequencies,
                      [1, :] are the amplitudes.
    im (int):         Index of mode
 
Returns
    Error code
errmsg(err) -> 'str'
errmsg(err: int) -> str
 
Get the message of an error code
 
Returns:
    str: The error message
find_all_cold_resonances( vn_z, vn_mod, tdm_c, tdm_zk, buf_cfreq, opt_iz_min_fcut=None ) -> 'tuple[int, int]'
find_all_cold_resonances(vn_z: real[:], vn_mod: real[:,_:], tdm_c: real[:], tdm_zk: real[:], buf_cfreq: double_complex[:], opt_iz_min_fcut=None) -> tuple[int, int]
 
Parameters:
    vn_z:       Shape (dim_z,)
    vn_mod:     Shape (DIM_MOD, dim_z)
    tdm_c:      Shape (3*stride_real(dim_z),)
    tdm_zk:     Shape (3*stride_real(dim_z),)
                The FEM zk matrix of a certain mode.
    buf_cfreq:  Shape (dim_freq,)
                where dim_freq is a larger number than the estimated
                dimension of spectrum from `guess_dim_cold_spectrum()`.
                The first `dim_res` element buffer also holds the found
                resonances.
    opt_iz_min_fcut (int):  The index of grid node, whose cutoff frequency
                            becomes lowest boundary of search interval.
                            Ignored if opt_iz_min_fcut >= dim_z.
 
Returns
    A tuple (dim_res, e)
    dim_res: Number of resounances found
    e:       Error code
find_cold_resonance( vn_z, mvn_mod, tdm_c, mv_tdm_zk, out_vm_cfreq, opt_vm_freq_init=None, opt_iz_min_fcut=None ) -> 'int'
find_cold_resonance(vn_z: real[:], mvn_mod: real[:,_:,_:], tdm_c: real[:], mv_tdm_zk: real[:,_:], out_vm_cfreq: double_complex[:], opt_vm_freq_init: real[:] = None, opt_iz_min_fcut=None) -> int
 
Parameters:
    vn_z:         Shape (dim_z,)
    mvn_mod:      Shape (dim_m, DIM_MOD, dim_z)
    tdm_c:        Shape (3*stride_real(dim_z),)
    mv_tdm_zk:    Shape (dim_m, 3*stride_real(dim_z))
    out_vm_cfreq: Shape (dim_m,)
                          Complex resonance frequency
    opt_vm_freq_init: Shape (dim_m,)
                          The initial frequencies for the searching.
                          If `opt_vm_freq_init` is `None` or
                          `opt_vm_freq_init[im] <= 0`, the resonances
                          with the highest Q are searched; otherwise,
                          for each mode, a resonance close to the given
                          frequency is to be found.
    opt_iz_min_fcut (int):  The index of grid node, whose cutoff frequency
                            becomes lowest boundary of search interval.
                            Ignored if opt_iz_min_fcut >= dim_z.
 
Returns
    Error code
get_crock_version() -> 'str'
get_crock_version() -> str
guess_dim_cold_spectrum(vn_z, mvn_mod, im=0) -> 'int'
guess_dim_cold_spectrum(vn_z: real[:], mvn_mod: real[:,_:,_:], im: int = 0) -> int
 
Parameters:
    vn_z: Shape (dim_z,)
    mvn_mod: Shape (dim_m, DIM_MOD, dim_z)
    im (int): Index of mode
 
Returns
    int: A suggested number of frequency points in the spectrum,
         or 0 when error occurs
init_bd(vm_bc, mvn_env, out_vm_bd) -> 'None'
init_bd(vm_bc: double_complex[:,_:], mvn_env: double_complex[:,_:], out_vm_bd: double_complex[:,_:]) -> None
 
Initialize boundary data.
 
It assumes a stationary state a priori to the start of the simulation.
 
Since this function is used in all system simulations, instead of letting
user manually initialize the vm_bd every time, CROCK provides a thin wrap
of the process.
 
Parameters:
    vm_bc:     Shape (DIM_BC, dim_m) boundary condition
    mvn_env:   Shape (dim_m, dim_z) Complex envelope profile
    out_vm_bd: Shape (DIM_BD, dim_m) boundary data
                       If None then new array will be allocated.
 
Returns:
    vm_bd: Shape (DIM_BD, dim_m) boundary data
init_bp(mvn_env, out_vm_bp) -> 'None'
init_bp(mvn_env: double_complex[:,_:], out_vm_bp: real[:,_:]) -> None
 
Parameters:
    out_vm_bp: Shape (2, dim_m) boundary phase
make_bc_simple(vm_freq, mvn_mod, out_vm_bc) -> 'int'
make_bc_simple(vm_freq: real[:], mvn_mod: real[:,_:,_:], out_vm_bc: double_complex[:,_:]) -> int
 
Parameters:
    mvn_mod:    Shape (dim_m, DIM_MOD, dim_z)
    out_vm_bc:  Shape (DIM_BC, dim_m)
                        Boundary condition
 
Returns
    Error code
make_bcc(r_center, r_span, vm_num, vn_mag_axis, mvn_mod, out_nvmh_bcc) -> 'int'
make_bcc(r_center: float, r_span: float, vm_num: int16_t[:,_:], vn_mag_axis: real[:], mvn_mod: real[:,_:,_:], out_nvmh_bcc: real[:,_:,_:,_:]) -> int
 
Make grid for beam coupling coefficient
 
This is not an equi-flux grid, but equidistant grid, in order to
make `dim_bcc=1` at `r_center` rather than `flux_center`
 
Parameters:
    r_center : Center radius (in meter) of BCC grid, at the
               **entrance** of the cavity. Typically one can choose
               the radius of the norminal beam guiding center.
    r_span   : Relative span, should be larger than the beam radius
               span, otherwise the range of BCC grid will be cropped.
 
Returns:
    Error code
make_chi_TE_corr_holl_imp(d, w, num_azi, num_rad, out_chi) -> 'bool'
make_chi_TE_corr_holl_imp(d: float, w: float, num_azi: int, num_rad: int, out_chi: real[:]) -> bool
 
List eigenvalues of hollow mode with longitudinal corrugated outer wall
from TE(m,1) to TE(m,n) assuming infinitely dense of corrugations
(Impedance Model)
 
Parameters:
    d:          Ratio corrugation_depth / r_out
    w:          Ratio of (corrugation width) / (width of a period)
    num_azi:    Azimuthal number (i.e. the `m`)
    num_rad:    The maximum order of solution (i.e. the `n`)
    out_chi:    Shape(num_rad,)
                The buffer to store the results, which a vector
                of eigen values for TE(m,1) to TE(m, num_rad).
                Attention: TE(m,1) has the [0] index!
 
Returns:
    False if error occurs, otherwise True
make_chi_TE_corr_ins(r, d, w, num_azi, num_rad, out_chi) -> 'bool'
make_chi_TE_corr_ins(r: float, d: float, w: float, num_azi: int, num_rad: int, out_chi: real[:]) -> bool
 
List eigenvalues of coaxial mode with corrugated insert
from TE(m,1) to TE(m,n).
 
Parameters:
    r:          Ratio r_in / r_out
    d:          Ratio corrugation_depth / r_out
    w:          Ratio of (corrugation width) / (width of a period)
    num_azi:    Azimuthal number (i.e. the `m`)
    num_rad:    The maximum order of solution (i.e. the `n`)
    out_chi:    Shape(num_rad,)
                The buffer to store the results, which a vector
                of eigen values for TE(m,1) to TE(m, num_rad).
                Attention: TE(m,1) has the [0] index!
 
Returns:
    False if error occurs, otherwise True
make_cold_envelope(vm_cfreq, vn_z, mvn_mod, out_mvn_env) -> 'int'
make_cold_envelope(vm_cfreq: double_complex[:], vn_z: real[:], mvn_mod: real[:,_:,_:], out_mvn_env: double_complex[:,_:]) -> int
 
Parameters:
    vm_cfreq: Shape (dim_z,) Complex resonance frequency of modes
    vn_z: Shape (dim_z,)
    mvn_mod: Shape (dim_m, DIM_MOD, dim_z)
    out_mvn_env: Shape (dim_m, dim_z) Complex envelope function
Returns
    Error code
make_cold_spectrum( vn_z, mvn_mod, tdm_c, mv_tdm_zk, dim_freq, spectrum, im=0, max_iter=0 ) -> 'int'
make_cold_spectrum(vn_z: real[:], mvn_mod: real[:,_:,_:], tdm_c: real[:], mv_tdm_zk: real[:,_:], dim_freq: int, spectrum: real[:,_:], im: int = 0, max_iter: int = 0) -> int
 
Parameters:
    vn_z:       Shape (dim_z,)
    mvn_mod:    Shape (dim_m, DIM_MOD, dim_z)
    dim_freq (int):     Number of spectral points
    tdm_c:      Shape (3*stride_real(dim_z),)
    mv_tdm_zk:  Shape (dim_m, 3*stride_real(dim_z))
    spectrum: Shape (2, dim_freq) array whose vectors
                      [0, :] are the frequencies,
                      [1, :] are the amplitudes.
    im (int):           Index of mode
 
Returns
    Error code
make_df(elapsed_time, mvn_new_env, mvn_old_env_phase, out_vm_df) -> 'None'
make_df(elapsed_time: float, mvn_new_env: double_complex[:,_:], mvn_old_env_phase: real[:,_:], out_vm_df: real[:]) -> None
 
Make delta frequency
 
Calculate the frequency shift and overwrite the
mvn_old_env_phase by the mvn_new_env
 
Parameters:
    mvn_new_env:        Shape (dim_m, dim_z)
                        Current complex envelope
    mvn_old_env_phase:  Shape (dim_m, dim_z)
                        Previous (real) phase
    out_vm_df:          Shape (dim_m,)
                        Buffer for the output
 
Return:
    Frequency shift (vm_df)
make_fem(vn_z, mvn_mod, out_tdm_c, out_mv_tdm_zk) -> 'int'
make_fem(vn_z: real[:], mvn_mod: real[:,_:,_:], out_tdm_c: real[:], out_mv_tdm_zk: real[:,_:]) -> int
 
Parameters:
    vn_z: Shape (dim_z,)
    mvn_mod : Shape (dim_m, DIM_MOD, dim_z)
    out_tdm_c : Shape (3*stride_real(dim_z),)
    out_mv_tdm_zk: Shape (dim_m, 3*stride_real(dim_z))
 
Returns
    Error code
make_mode_coax( vm_num, vn_geo_rout, vn_geo_rin, out_mvn_mod, vn_ins_corr_ratio=None, vn_ins_corr_depth=None ) -> 'int'
make_mode_coax(vm_num: int16_t[:,_:], vn_geo_rout: real[:], vn_geo_rin: real[:], out_mvn_mod: real[:,_:,_:], vn_ins_corr_ratio: real[:] = None, vn_ins_corr_depth: real[:] = None) -> int
 
Makes a hollow cavity
 
Parameters:
    vm_num:             Shape (2, dim_m)
                        Mode number, where vm_num[0,:] are the azimuthal
                        indices, vm_num[1,:] the radial indices
    vn_geo_rout:        Shape (dim_z,)
                        Outer radius
    out_mvn_mod:        Shape (dim_m, DIM_MOD, dim_z)
    vn_ins_corr_ratio:  Shape (dim_z,)
                        Corrugation ratio (dimensionless)
    vn_ins_corr_depth:  Shape (dim_z,)
                        Corrugation depth in meter
 
Returns:
    Error code
make_mode_corr_holl_imp( vm_num, vn_geo_rout, out_mvn_mod, vn_ins_corr_ratio=None, vn_ins_corr_depth=None ) -> 'int'
make_mode_corr_holl_imp(vm_num: int16_t[:,_:], vn_geo_rout: real[:], out_mvn_mod: real[:,_:,_:], vn_ins_corr_ratio: real[:] = None, vn_ins_corr_depth: real[:] = None) -> int
 
    Makes a hollow cavity
 
    Parameters:
        vm_num:             Shape (2, dim_m)
                            Mode number, where vm_num[0,:] are the azimuthal
                            indices, vm_num[1,:] the radial indices
        vn_geo_rout:        Shape (dim_z,)
                            Outer radius
        out_mvn_mod:        Shape (dim_m, DIM_MOD, dim_z)
        vn_ins_corr_ratio:  Shape (dim_z,)
                            Corrugation ratio (dimensionless)
        vn_ins_corr_depth:  Shape (dim_z,)
                            Corrugation depth in meter
Returns:
        Error code
make_mode_holl(vm_num, vn_geo_rout, out_mvn_mod) -> 'int'
make_mode_holl(vm_num: int16_t[:,_:], vn_geo_rout: real[:], out_mvn_mod: real[:,_:,_:]) -> int
 
Makes a hollow cavity
 
Parameters:
    vm_num:         Shape (2, dim_m)
                    Mode number, where vm_num[0,:] are the azimuthal
                    indices, vm_num[1,:] the radial indices
    vn_geo_rout:    Shape (dim_z,)
                    Outer radius
    out_mvn_mod:    Shape (dim_m, DIM_MOD, dim_z)
 
Returns:
    Error code
make_potential( pot_outer, beam_init, vn_mag_axis, vn_rhocd, vn_r_outer, out_vn_pot, opt_pot_inner=0.0, opt_vn_r_inner=None ) -> 'None'
make_potential(pot_outer: float, beam_init: real[:], vn_mag_axis: real[:], vn_rhocd: real[:], vn_r_outer: real[:], out_vn_pot: real[:], opt_pot_inner: float = 0, opt_vn_r_inner: real[:] = None) -> None
 
Make depression voltage from charge density and geometry
 
This is a simple mode for a thin beam, which supports also coaxial cavity.
However, the user is encouraged to use own advanced models.
 
Parameters:
    pot_outer       : Potential of outer boundary
    beam_init       : From which the initial guiding center is extracted
    vn_mag_axis     : B field on the axis
    vn_rhocd        : Charge density (which has negative values!)
    vn_r_outer      : Radius of outer boundary
    out_vn_pot      : Electric potential
    opt_pot_inner   : Potential of coaxial insert (unused for hollow cavity)
    opt_vn_r_inner  : Radius of coaxial insert (NULL for hollow cavity)
 
Returns:
    None
make_power_out( vm_cf, mvn_mod, mvn_env, opt_vm_pout_left=None, opt_vm_pout_right=None ) -> 'None'
make_power_out(vm_cf: real[:], mvn_mod: real[:,_:,_:], mvn_env: double_complex[:,_:], opt_vm_pout_left: real[:] = None, opt_vm_pout_right: real[:] = None) -> None
make_uload(mvn_env, mvn_uuload, out_vn_uload) -> 'None'
make_uload(mvn_env: double_complex[:,_:], mvn_uuload: real[:,_:], out_vn_uload: real[:]) -> None
 
Make unscaled ohmic load
 
Parameters:
    mvn_uuload:     The output from `make_uuload`
    out_mvn_uload:  This vector divdied by sqrt(sigma) at each nodes
                    results in is the ohmic load in W/m^2
make_uuload(vm_num, vm_freq, vn_radius, mvn_mod, out_mvn_uuload) -> 'None'
make_uuload(vm_num: int16_t[:,_:], vm_freq: real[:], vn_radius: real[:], mvn_mod: real[:,_:,_:], out_mvn_uuload: real[:,_:]) -> None
 
Make very unscaled ohmic load per mode
 
Parameters:
    vn_radius:      Shape (dim_z,)
                    Geometry of outer wall or coaxial insert
    out_mvn_uuload: Shape (dim_m, dim_z)
                    Will be the input of `make_uload`
numa_bless(array)
numa_bless(array)
 
Set the data of `array` to zero
 
It performs the NUMA first-touch onto array.shape[0] static parallel nodes,
but not assert the alignment of page boundaries on array.shape[0].
refine_z(weight, vn_z_old, vn_mrs, vn_z_new) -> 'int'
refine_z(weight: float, vn_z_old: real[:], vn_mrs: real[:], vn_z_new: real[:]) -> int
 
Remap meshes considering the mesh refinement score
 
Parameters:
    weight:     weight=0 means not using `vn_mrs`
                while keeping the original density;
                weight=1 is a moderate value which does not make
                         the coarse segments too much coarser;
                weight=infinite means only weighted by `vn_mrs`.
    vn_z_old:   Shape (dim_z_old,)
    vn_mrs:     Shape (dim_z_old,)
                Mesh Refinement Score
    vn_z_new:   Shape (dim_z_new,)
 
Return:
    Error code
search_cold_resonance( vn_z, mvn_mod, out_vm_freq, opt_ref_z: 'Optional[float]' = None ) -> 'int'
search_cold_resonance(vn_z: real[:], mvn_mod: real[:,_:,_:], out_vm_freq: double_complex[:], opt_ref_z: Optional[float] = None) -> int
 
Parameters:
    vn_z:        Shape (dim_z,)
    mvn_mod:     Shape (dim_m, DIM_MOD, dim_z)
    out_vm_freq: Shape (dim_m,) Complex resonance frequency
    opt_ref_z:   An optional z-coordinate whose cut-off frequency
                 will become the minimum frequency of the
                 searching range. This is designed for gyrotron
                 cavity simulation including up-taper, where the
                 low-frequency resonances with profile maximum
                 within the up-taper can be skipped. If None, the
                 minimum cutoff of the whole geometry is the
                 lower boundary of frequency range.
 
Returns
    Error code
shoot_beamlet( options, beam_init, vm_num, vm_cf, vm_cp, vn_z, vn_pot, vn_mag_axis, mvn_mod, mvn_env, dim_b, bcc_r_center, bcc_r_span, nvmh_bcc, io_seed, out_vn_rhocd, out_mvn_src, opt_out_vn_mrs=None, opt_out_vn_blt_avg=None, opt_out_bvn_blt_trj=None, opt_cache: 'Cache' = None ) -> 'int'
shoot_beamlet(options: int, beam_init: real[:], vm_num: int16_t[:,_:], vm_cf: real[:], vm_cp: real[:], vn_z: real[:], vn_pot: real[:], vn_mag_axis: real[:], mvn_mod: real[:,_:,_:], mvn_env: double_complex[:,_:], dim_b: int, bcc_r_center: float, bcc_r_span: float, nvmh_bcc: real[:,_:,_:,_:], io_seed: uint32_t[:], out_vn_rhocd: real[:], out_mvn_src: double_complex[:,_:], opt_out_vn_mrs: real[:] = None, opt_out_vn_blt_avg: real[:,_:] = None, opt_out_bvn_blt_trj: real[:,_:,_:] = None, opt_cache: Cache = None) -> int
sizeof_real() -> 'int'
sizeof_real() -> int
step_field( dt, vm_cf, vm_bc, tdm_c, mv_tdm_zk, mvn_src, io_vm_bd, io_mvn_env, opt_io_vm_cp=None, out_opt_vm_df=None, opt_cache: 'Cache' = None ) -> 'int'
step_field(dt: float, vm_cf: real[:], vm_bc: double_complex[:,_:], tdm_c: real[:], mv_tdm_zk: real[:,_:], mvn_src: double_complex[:,_:], io_vm_bd: double_complex[:,_:], io_mvn_env: double_complex[:,_:], opt_io_vm_cp: real[:] = None, out_opt_vm_df: real[:] = None, opt_cache: Cache = None) -> int
 
Step the field envelope profile
 
Parameters
    vm_cf:          Shape (dim_m,)
                    Carrier frequency
    vm_bc:          Shape (DIM_BC, dim_m)
                    Boundary condition
    tdm_c:          Shape (3*stride_real(dim_z),)
                    Matrix from make_fem()
    mv_tdm_zk:      Shape (dim_m, 3*stride_real(dim_z))
                    Matrix from make_fem()
    mvn_src:        Shape (dim_m, dim_z)
                    Non-weighted excitation per mode (unit: V/m^2)
    vm_bd:          Shape (DIM_BD, dim_m)
                    Boundary data
    io_mvn_env:     Shape (dim_m, dim_z)
                    Envelope profile
    opt_io_vm_cp:   Shape (dim_m,)
                    Carrier phase, can be None
    out_opt_vm_df:  Shape (dim_m,)
                    Frequency shift, set to None to reduce computational effort
    opt_cache:      Empty or previously used cache with same problem dimensions
 
Returns:
    Error code
stride_cplx(dim: 'int') -> 'int'
stride_cplx(dim: int) -> int
stride_real(dim: 'int') -> 'int'
stride_real(dim: int) -> int
update_bc_simple(dt, vm_cf, mvn_mod, mvn_env, io_vm_bp, out_vm_bc) -> 'None'
update_bc_simple(dt: float, vm_cf: real[:], mvn_mod: real[:,_:,_:], mvn_env: double_complex[:,_:], io_vm_bp: real[:,_:], out_vm_bc: double_complex[:,_:]) -> None
 
Parameters:
    vm_cf:     Shape (dim_m,) carrier frequency
    mvn_mod:   Shape (dim_m, DIM_MOD, dim_z)
    io_vm_bp:  Shape (2, dim_m) boundary phase
    out_vm_bc: Shape (DIM_BC, dim_m) boundary condition

Data descriptors defined here:
__dict__
dictionary for instance variables
__weakref__
list of weak references to the object

Data and other attributes defined here:
Cache = <class 'prock.Cache'>
DIM_BC = <ROW_OF_BC.DIM_BC: 4>
DIM_BD = <ROW_OF_BD.DIM_BD: 2>
DIM_BEAM_INIT = <ROW_OF_BEAM_INIT.DIM_BEAM_INIT: 7>
DIM_BLT_AVG = <ROW_OF_BLT_AVG.DIM_BLT_AVG: 3>
DIM_BLT_TRJ = <ROW_OF_BLT_TRJ.DIM_BLT_TRJ: 6>
DIM_HARM = 4
DIM_MOD = <ROW_OF_MOD.DIM_MOD: 3>
E_TOTAL = 45
MASK_BLT_ALL_FEATURES = 63
MASK_BLT_CALI_CYC_PHASE = 4
MASK_BLT_CALI_GC_AZI = 8
MASK_BLT_CALI_GC_R = 16
MASK_BLT_DEBUG = 2147483648
MASK_BLT_DOPPLER = 32
MASK_BLT_MONOGC = 2147483648
MASK_BLT_SRC_MONOHARM = 1073741824
MASK_BLT_TAYLOR = 2
MASK_BLT_VAR_UZ = 1
MAX_HARM = 3
MIN_HARM = 0
RID_BEAM_A = <ROW_OF_BEAM_INIT.RID_BEAM_A: 2>
RID_BEAM_AS = <ROW_OF_BEAM_INIT.RID_BEAM_AS: 4>
RID_BEAM_G = <ROW_OF_BEAM_INIT.RID_BEAM_G: 3>
RID_BEAM_GS = <ROW_OF_BEAM_INIT.RID_BEAM_GS: 5>
RID_BEAM_I = <ROW_OF_BEAM_INIT.RID_BEAM_I: 0>
RID_BEAM_R = <ROW_OF_BEAM_INIT.RID_BEAM_R: 1>
RID_BEAM_RS = <ROW_OF_BEAM_INIT.RID_BEAM_RS: 6>
RID_BLT_AVG_ALPHA = <ROW_OF_BLT_AVG.RID_BLT_AVG_ALPHA: 0>
RID_BLT_AVG_FREQ = <ROW_OF_BLT_AVG.RID_BLT_AVG_FREQ: 2>
RID_BLT_AVG_GAMMA = <ROW_OF_BLT_AVG.RID_BLT_AVG_GAMMA: 1>
RID_BLT_TRJ_ALPHA = <ROW_OF_BLT_TRJ.RID_BLT_TRJ_ALPHA: 4>
RID_BLT_TRJ_CYC_ANG = <ROW_OF_BLT_TRJ.RID_BLT_TRJ_CYC_ANG: 3>
RID_BLT_TRJ_CYC_R = <ROW_OF_BLT_TRJ.RID_BLT_TRJ_CYC_R: 2>
RID_BLT_TRJ_GAMMA = <ROW_OF_BLT_TRJ.RID_BLT_TRJ_GAMMA: 5>
RID_BLT_TRJ_GC_ANG = <ROW_OF_BLT_TRJ.RID_BLT_TRJ_GC_ANG: 1>
RID_BLT_TRJ_GC_R = <ROW_OF_BLT_TRJ.RID_BLT_TRJ_GC_R: 0>
RID_C_BESL_J = <ROW_OF_MOD.RID_C_BESL_J: 1>
RID_C_BESL_Y = <ROW_OF_MOD.RID_C_BESL_Y: 2>
RID_K_PERP = <ROW_OF_MOD.RID_K_PERP: 0>
RID_UZL = <ROW_OF_BD.RID_UZL: 0>
RID_UZR = <ROW_OF_BD.RID_UZR: 1>
RID_X1L = <ROW_OF_BC.RID_X1L: 0>
RID_X1R = <ROW_OF_BC.RID_X1R: 2>
RID_X2L = <ROW_OF_BC.RID_X2L: 1>
RID_X2R = <ROW_OF_BC.RID_X2R: 3>
ROW_OF_BC = <flag 'ROW_OF_BC'>
ROW_OF_BD = <flag 'ROW_OF_BD'>
ROW_OF_BEAM_INIT = <flag 'ROW_OF_BEAM_INIT'>
ROW_OF_BLT_AVG = <flag 'ROW_OF_BLT_AVG'>
ROW_OF_BLT_TRJ = <flag 'ROW_OF_BLT_TRJ'>
ROW_OF_MOD = <flag 'ROW_OF_MOD'>
W_TOTAL = 12
developer_mode = True