desc.grid.LinearGrid

class desc.grid.LinearGrid(L=None, M=None, N=None, NFP=1, sym=False, axis=True, endpoint=False, rho=None, theta=None, zeta=None)Source

Grid in which the nodes are linearly spaced in each coordinate.

Useful for plotting and other analysis, though not very efficient for using as the solution grid.

Parameters:
  • L (int, optional) – Radial grid resolution.

  • M (int, optional) – Poloidal grid resolution.

  • N (int, optional) – Toroidal grid resolution.

  • NFP (int) – Number of field periods (Default = 1). Change this only if your nodes are placed within one field period or should be interpreted as spanning one field period.

  • sym (bool) – True for poloidal up/down symmetry, False otherwise. Default is False. Whether to truncate the poloidal domain to [0, π] ⊂ [0, 2π) to take advantage of poloidal up/down symmetry, which is a stronger condition than stellarator symmetry. Still, when stellarator symmetry exists, flux surface integrals and volume integrals are invariant to this truncation.

  • axis (bool) – True to include a point at rho=0 (default), False for rho[0] = rho[1]/2.

  • endpoint (bool) – If True, theta=0 and zeta=0 are duplicated after a full period. Should be False for use with FFT. (Default = False). This boolean is ignored if an array is given for theta or zeta.

  • rho (int or ndarray of float, optional) – Radial coordinates (Default = 1.0). Alternatively, the number of radial coordinates (if an integer). Note that if supplied the values may be reordered in the resulting grid.

  • theta (int or ndarray of float, optional) – Poloidal coordinates (Default = 0.0). Alternatively, the number of poloidal coordinates (if an integer). Note that if supplied the values may be reordered in the resulting grid.

  • zeta (int or ndarray of float, optional) – Toroidal coordinates (Default = 0.0). Alternatively, the number of toroidal coordinates (if an integer). Note that if supplied the values may be reordered in the resulting grid.

Methods

change_resolution(L, M, N[, NFP])

Change the resolution of the grid.

compress(x[, surface_label])

Return elements of x at indices of unique surface label values.

copy([deepcopy])

Return a (deep)copy of this object.

copy_data_from_other(x, other_grid[, ...])

Copy data x from other_grid to this grid at matching surface label.

equiv(other)

Compare equivalence between DESC objects.

expand(x[, surface_label])

Expand x by duplicating elements to match the grid's pattern.

get_label(label)

Get general label that specifies direction given label.

load(load_from[, file_format])

Initialize from file.

meshgrid_flatten(x, order)

Flatten data to match standard ordering.

meshgrid_reshape(x, order)

Reshape data to match grid coordinates.

replace_at_axis(x, y[, copy])

Replace elements of x with elements of y at the axis of grid.

save(file_name[, file_format, file_mode])

Save the object.

Attributes

L

Radial grid resolution.

M

Poloidal grid resolution.

N

Toroidal grid resolution.

NFP

Number of (toroidal) field periods.

axis

Indices of nodes at magnetic axis.

can_fft2

Whether this grid is compatible with 2D FFT.

coordinates

Coordinates specified by the nodes.

endpoint

Whether the grid is made of open or closed intervals.

fft_poloidal

whether this grid is compatible with fft in the poloidal direction.

fft_toroidal

whether this grid is compatible with fft in the toroidal direction.

inverse_alpha_idx

Indices that recover field line poloidal angles.

inverse_poloidal_idx

Indices that recover the unique poloidal coordinates.

inverse_rho_idx

Indices of unique_rho_idx that recover the rho coordinates.

inverse_theta_PEST_idx

Indices that recover unique straight field line poloidal angles.

inverse_theta_idx

Indices that recover unique theta coordinates.

inverse_zeta_idx

Indices of unique_zeta_idx that recover the zeta coordinates.

is_meshgrid

Whether this grid is a tensor-product grid.

node_pattern

Pattern for placement of nodes in (rho,theta,zeta).

nodes

Node coordinates, in (rho,theta,zeta).

num_alpha

Number of unique field line poloidal angles.

num_nodes

Total number of nodes.

num_poloidal

Number of unique poloidal angle coordinates.

num_rho

Number of unique rho coordinates.

num_theta

Number of unique theta coordinates.

num_theta_PEST

Number of unique straight field line poloidal angles.

num_zeta

Number of unique zeta coordinates.

period

Periodicity of coordinates.

spacing

Quadrature weights for integration over surfaces.

sym

True for poloidal up/down symmetry, False otherwise.

unique_alpha_idx

Indices of unique field line poloidal angles.

unique_poloidal_idx

Indices of unique poloidal angle coordinates.

unique_rho_idx

Indices of unique rho coordinates.

unique_theta_PEST_idx

Indices of unique straight field line poloidal angles.

unique_theta_idx

Indices of unique theta coordinates.

unique_zeta_idx

Indices of unique zeta coordinates.

weights

Weight for each node, either exact quadrature or volume based.