simsopt.geo package

Submodules

simsopt.geo.biotsavart module

class simsopt.geo.biotsavart.BiotSavart(coils, coil_currents)

Bases: simsopt.geo.magneticfield.MagneticField

B_and_dB_vjp(v, vgrad)
B_vjp(v)
compute(points, compute_derivatives=0)
compute_A(points, compute_derivatives=0)
d2B_by_dXdcoilcurrents(compute_derivatives=1)
d3B_by_dXdXdcoilcurrents(compute_derivatives=2)
dB_by_dcoilcurrents(compute_derivatives=0)

simsopt.geo.boozersurface module

class simsopt.geo.boozersurface.BoozerSurface(biotsavart, surface, label, targetlabel)

Bases: object

BoozerSurface and its associated methods can be used to compute the Boozer angles on a surface. It takes a Surface representation (e.g. SurfaceXYZFourier, or SurfaceXYZTensorFourier), a magnetic field evaluator, surface label evaluator, and a target surface label.

The Boozer angles are computed by solving a constrained least squares problem. The least squares objective is given by \(J(x) = \frac{1}{2} \mathbf r^T(x) \mathbf r(x)\), where \(\mathbf r\) is a vector of residuals computed by boozer_surface_residual (see surfaceobjectives.py), and some constraints. This objective is zero when the surface corresponds to a magnetic surface of the field and \((\phi,\theta)\) that parametrize the surface correspond to Boozer angles, and the constraints are satisfied.

The surface label can be area, volume, or toroidal flux. The label on the computed surface will be equal or close to the user-provided targetlabel, depending on how the label constraint is imposed. This constrained least squares problem can be solved by scalarizing and adding the constraint as an additional penalty term to the objective. This is done in

  1. minimize_boozer_penalty_constraints_LBFGS

  2. minimize_boozer_penalty_constraints_newton

  3. minimize_boozer_penalty_constraints_ls

where LBFGS, Newton, or scipy.optimize.least_squares optimizers are used, respectively. Alternatively, the exactly constrained least squares optimization problem can be solved. This is done in

  1. minimize_boozer_exact_constraints_newton

where Newton is used to solve the first order necessary conditions for optimality.

boozer_exact_constraints(xl, derivatives=0, optimize_G=True)

This function returns the optimality conditions corresponding to the minimization problem

\[\text{min}_x ~J(x)\]

subject to

\[\begin{split}l - l_0 &= 0 \\ z(\varphi=0,\theta=0) - 0 &= 0\end{split}\]

where \(l\) is the surface label and \(l_0\) is the target surface label, \(J(x) = \frac{1}{2}\mathbf r(x)^T \mathbf r(x)\), and \(\mathbf r(x)\) contains the Boozer residuals at quadrature points \(1,\dots,n\). We can also optionally return the first derivatives of these optimality conditions.

boozer_penalty_constraints(x, derivatives=0, constraint_weight=1.0, scalarize=True, optimize_G=False)

Define the residual

\[\mathbf r(x) = [r_1(x),...,r_n(x), \sqrt{w_c} (l-l_0), \sqrt{w_c} (z(\varphi=0, \theta=0) - 0)]\]

where \(w_c\) is the constraint weight, \(r_i\) are the Boozer residuals at quadrature points \(1,\dots,n\), \(l\) is the surface label, and \(l_0\) is the target surface label.

For scalarized=False, this function returns \(\mathbf r(x)\) and optionally the Jacobian of \(\mathbf r(x)\).

for scalarized=True, this function returns

\[J(x) = \frac{1}{2}\mathbf r(x)^T \mathbf r(x),\]

i.e. the least squares residual and optionally the gradient and the Hessian of \(J(x)\).

minimize_boozer_exact_constraints_newton(tol=1e-12, maxiter=10, iota=0.0, G=None, lm=[0.0, 0.0])

This function solves the constrained optimization problem

\[\text{min}_x ~ J(x)\]

subject to

\[\begin{split}l - l_0 &= 0 \\ z(\varphi=0,\theta=0) - 0 &= 0\end{split}\]

using Lagrange multipliers and Newton’s method. In the above, \(J(x) = \frac{1}{2}\mathbf r(x)^T \mathbf r(x)\), and \(\mathbf r(x)\) contains the Boozer residuals at quadrature points \(1,\dots,n\).

The final constraint is not necessary for stellarator symmetric surfaces as it is automatically satisfied by the stellarator symmetric surface parametrization.

minimize_boozer_penalty_constraints_LBFGS(tol=0.001, maxiter=1000, constraint_weight=1.0, iota=0.0, G=None)

This function tries to find the surface that approximately solves

\[\text{min}_x ~J(x) + \frac{1}{2} w_c (l - l_0)^2 + \frac{1}{2} w_c (z(\varphi=0, \theta=0) - 0)^2\]

where \(J(x) = \frac{1}{2}\mathbf r(x)^T \mathbf r(x)\), and \(\mathbf r(x)\) contains the Boozer residuals at quadrature points \(1,\dots,n\). This is done using LBFGS.

minimize_boozer_penalty_constraints_ls(tol=1e-12, maxiter=10, constraint_weight=1.0, iota=0.0, G=None, method='lm')

This function does the same as minimize_boozer_penalty_constraints_LBFGS, but instead of LBFGS it uses a nonlinear least squares algorithm when method='lm'. Options for the method are the same as for scipy.optimize.least_squares. If method='manual', then a damped Gauss-Newton method is used.

minimize_boozer_penalty_constraints_newton(tol=1e-12, maxiter=10, constraint_weight=1.0, iota=0.0, G=None, stab=0.0)

This function does the same as minimize_boozer_penalty_constraints_LBFGS, but instead of LBFGS it uses Newton’s method.

simsopt.geo.coilcollection module

class simsopt.geo.coilcollection.CoilCollection(coils, currents, nfp, stellarator_symmetry)

Bases: object

This class represents a collection of coils and currents. For stellarators with symmetries (rotational, stellarator symmetry), this class performs reflections and rotations to generate a full set of stellarator coils from a base set of modular coils.

simsopt.geo.config module

simsopt.geo.curve module

class simsopt.geo.curve.Curve

Bases: simsopt._core.optimizable.Optimizable

Curve is a base class for various representations of curves in SIMSOPT.

dfrenet_frame_by_dcoeff()

This function returns the derivative of the curve’s Frenet frame,

\[\left(\frac{\partial \mathbf{t}}{\partial \mathbf{c}}, \frac{\partial \mathbf{n}}{\partial \mathbf{c}}, \frac{\partial \mathbf{b}}{\partial \mathbf{c}}\right),\]

with respect to the curve dofs, where \((\mathbf t, \mathbf n, \mathbf b)\) correspond to the Frenet frame, and \(\mathbf c\) are the curve dofs.

dincremental_arclength_by_dcoeff_vjp(v)

This function returns the vector Jacobian product

\[v^T \frac{\partial \|\Gamma'\|}{\partial \mathbf{c}}\]

where \(\|\Gamma'\|\) is the incremental arclength, \(\Gamma'\) is the tangent to the curve and \(\mathbf{c}\) are the curve dofs.

dkappa_by_dcoeff_impl(dkappa_by_dcoeff)

This function computes the derivative of the curvature with respect to the curve coefficients.

\[\frac{\partial \kappa}{\partial \mathbf c}\]

where \(\mathbf c\) are the curve dofs, \(\kappa\) is the curvature.

dkappa_by_dcoeff_vjp(v)

This function returns the vector Jacobian product

\[v^T \frac{\partial \kappa}{\partial \mathbf{c}}\]

where \(\mathbf c\) are the curve dofs and \(\kappa\) is the curvature.

dkappadash_by_dcoeff()

This function returns

\[\frac{\partial \kappa'(\phi)}{\partial \mathbf{c}}.\]

where \(\mathbf{c}\) are the curve dofs, and \(\kappa\) is the curvature.

dtorsion_by_dcoeff_impl(dtorsion_by_dcoeff)

This function returns the derivative of torsion with respect to the curve dofs.

\[\frac{\partial \tau}{\partial \mathbf c}\]

where \(\mathbf c\) are the curve dofs, and \(\tau\) is the torsion.

dtorsion_by_dcoeff_vjp(v)

This function returns the vector Jacobian product

\[v^T \frac{\partial \tau}{\partial \mathbf{c}}\]

where \(\mathbf c\) are the curve dofs, and \(\tau\) is the torsion.

frenet_frame()

This function returns the Frenet frame, \((\mathbf{t}, \mathbf{n}, \mathbf{b})\), associated to the curve.

kappa_impl(kappa)

This function implements the curvature, \(\kappa(\varphi)\).

kappadash()

This function returns \(\kappa'(\phi)\), where \(\kappa\) is the curvature.

plot(ax=None, show=True, plot_derivative=False, closed_loop=True, color=None, linestyle=None)

Plot the curve using matplotlib.pyplot, along with optionally its tangent when plot_derivative=True. When close_loop=False the first and final point on the surface will not be connected, and when it is True, they will be connected by a line segment and a closed curve will be plotted.

plot_mayavi(show=True)

Plot the curve using mayavi.mlab rather than matplotlib.pyplot.

torsion_impl(torsion)

This function returns the torsion, \(\tau\), of a curve.

class simsopt.geo.curve.JaxCurve(quadpoints, gamma_pure)

Bases: simsgeopp.Curve, simsopt.geo.curve.Curve

dgamma_by_dcoeff_impl(dgamma_by_dcoeff)

This function returns

\[\frac{\partial \Gamma}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgamma_by_dcoeff_vjp(v)

This function returns the vector Jacobian product

\[v^T \frac{\partial \Gamma}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgammadash_by_dcoeff_impl(dgammadash_by_dcoeff)

This function returns

\[\frac{\partial \Gamma'}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgammadash_by_dcoeff_vjp(v)

This function returns

\[\mathbf v^T \frac{\partial \Gamma'}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgammadashdash_by_dcoeff_impl(dgammadashdash_by_dcoeff)

This function returns

\[\frac{\partial \Gamma''}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgammadashdash_by_dcoeff_vjp(v)

This function returns the vector Jacobian product

\[v^T \frac{\partial \Gamma''}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgammadashdashdash_by_dcoeff_impl(dgammadashdashdash_by_dcoeff)

This function returns

\[\frac{\partial \Gamma'''}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgammadashdashdash_by_dcoeff_vjp(v)

This function returns the vector Jacobian product

\[v^T \frac{\partial \Gamma'''}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dkappa_by_dcoeff_vjp(v)

This function returns the vector Jacobian product

\[v^T \frac{\partial \kappa}{\partial \mathbf{c}}\]

where \(\mathbf{c}\) are the curve dofs and \(\kappa\) is the curvature.

dtorsion_by_dcoeff_vjp(v)

This function returns the vector Jacobian product

\[v^T \frac{\partial \tau}{\partial \mathbf{c}}\]

where \(\mathbf{c}\) are the curve dofs, and \(\tau\) is the torsion.

gamma_impl(gamma, quadpoints)

This function returns the x,y,z coordinates of the curve \(\Gamma\).

gammadash_impl(gammadash)

This function returns \(\Gamma'(\varphi)\), where \(\Gamma\) are the x, y, z coordinates of the curve.

gammadashdash_impl(gammadashdash)

This function returns \(\Gamma''(\varphi)\), and \(\Gamma\) are the x, y, z coordinates of the curve.

gammadashdashdash_impl(gammadashdashdash)

This function returns \(\Gamma'''(\varphi)\), and \(\Gamma\) are the x, y, z coordinates of the curve.

class simsopt.geo.curve.RotatedCurve(curve, theta, flip)

Bases: simsgeopp.Curve, simsopt.geo.curve.Curve

RotatedCurve inherits from the Curve base class. It takes an input a Curve, rotates it by theta, and optionally completes a reflection when flip=True.

dgamma_by_dcoeff_impl(dgamma_by_dcoeff)

This function returns

\[\frac{\partial \Gamma}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgamma_by_dcoeff_vjp(v)

This function returns the vector Jacobian product

\[v^T \frac{\partial \Gamma}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgammadash_by_dcoeff_impl(dgammadash_by_dcoeff)

This function returns

\[\frac{\partial \Gamma'}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgammadash_by_dcoeff_vjp(v)

This function returns the vector Jacobian product

\[v^T \frac{\partial \Gamma'}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgammadashdash_by_dcoeff_impl(dgammadashdash_by_dcoeff)

This function returns

\[\frac{\partial \Gamma''}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgammadashdash_by_dcoeff_vjp(v)

This function returns the vector Jacobian product

\[v^T \frac{\partial \Gamma''}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgammadashdashdash_by_dcoeff_impl(dgammadashdashdash_by_dcoeff)

This function returns

\[\frac{\partial \Gamma'''}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

dgammadashdashdash_by_dcoeff_vjp(v)

This function returns the vector Jacobian product

\[v^T \frac{\partial \Gamma'''}{\partial \mathbf c}\]

where \(\mathbf{c}\) are the curve dofs, and \(\Gamma\) are the x, y, z coordinates of the curve.

gamma_impl(gamma, quadpoints)

This function returns the x,y,z coordinates of the curve, \(\Gamma\), where \(\Gamma\) are the x, y, z coordinates of the curve.

gammadash_impl(gammadash)

This function returns \(\Gamma'(\varphi)\), where \(\Gamma\) are the x, y, z coordinates of the curve.

gammadashdash_impl(gammadashdash)

This function returns \(\Gamma''(\varphi)\), where \(\Gamma\) are the x, y, z coordinates of the curve.

gammadashdashdash_impl(gammadashdashdash)

This function returns \(\Gamma'''(\varphi)\), where \(\Gamma\) are the x, y, z coordinates of the curve.

get_dofs()

This function returns the curve dofs.

num_dofs()

This function returns the number of dofs associated to the curve.

set_dofs_impl(d)

This function sets the curve dofs.

simsopt.geo.curve.incremental_arclength_pure(d1gamma)

This function is used in a Python+Jax implementation of the curve arc length formula.

simsopt.geo.curve.incremental_arclength_vjp(d1gamma, v)
simsopt.geo.curve.kappa_pure(d1gamma, d2gamma)

This function is used in a Python+Jax implementation of formula for curvature.

simsopt.geo.curve.kappagrad0(d1gamma, d2gamma)
simsopt.geo.curve.kappagrad1(d1gamma, d2gamma)
simsopt.geo.curve.kappavjp0(d1gamma, d2gamma, v)
simsopt.geo.curve.kappavjp1(d1gamma, d2gamma, v)
simsopt.geo.curve.torsion_pure(d1gamma, d2gamma, d3gamma)

This function is used in a Python+Jax implementation of formula for torsion.

simsopt.geo.curve.torsionvjp0(d1gamma, d2gamma, d3gamma, v)
simsopt.geo.curve.torsionvjp1(d1gamma, d2gamma, d3gamma, v)
simsopt.geo.curve.torsionvjp2(d1gamma, d2gamma, d3gamma, v)

simsopt.geo.curvehelical module

class simsopt.geo.curvehelical.CurveHelical(quadpoints, order, n0, l0, R0, r0)

Bases: simsopt.geo.curve.JaxCurve

Curve representation of a helical coil. The helical coil positions are specified by a poloidal angle eta that is a function of the toroidal angle phi with Fourier coefficients \(A_k\) and \(B_k\). The poloidal angle is represented as

\[\eta = m_0 \phi/l_0 + \sum_k A_k \cos(n_0 \phi k/l_0) + B_k \sin(n_0 \phi k/l_0)\]
Parameters
  • quadpoints – number of grid points/resolution along the curve;

  • order – number of fourier coefficients, i.e., the length of the array A (or B);

  • n0 – toroidal periodicity/number of field periods

  • l0 – number of \(2\pi\) turns in \(\phi\)

  • R0 – major radius

  • r0 – minor radius

get_dofs(self: simsgeopp.Curve)List[float]
num_dofs(self: simsgeopp.Curve)int
set_dofs_impl(dofs)
simsopt.geo.curvehelical.jaxHelicalfouriercurve_pure(dofs, quadpoints, order, n0, l0, R0, r0)

simsopt.geo.curveobjectives module

class simsopt.geo.curveobjectives.CurveLength(curve)

Bases: simsopt._core.optimizable.Optimizable

CurveLength is a class that computes the length of a curve, i.e.

\[J = \int_{\text{curve}}~dl.\]
J()

This returns the value of the quantity.

dJ()

This returns the derivative of the quantity with respect to the curve dofs.

class simsopt.geo.curveobjectives.LpCurveCurvature(curve, p, desired_length=None)

Bases: simsopt._core.optimizable.Optimizable

This class computes a penalty term based on the \(L_p\) norm of the curve’s curvature, and penalizes where the local curve curvature exceeds a threshold

\[J = \frac{1}{p} \int_{\text{curve}} \text{max}(\kappa - \kappa_0, 0)^p ~dl\]

where \(\kappa_0\) is a threshold curvature. When desired_length is provided, \(\kappa_0 = 2 \pi / \text{desired_length}\), otherwise \(\kappa_0=0\).

J()

This returns the value of the quantity.

dJ()

This returns the derivative of the quantity with respect to the curve dofs.

class simsopt.geo.curveobjectives.LpCurveTorsion(curve, p)

Bases: simsopt._core.optimizable.Optimizable

LpCurveTorsion is a class that computes a penalty term based on the \(L_p\) norm of the curve’s torsion:

\[J = \frac{1}{p} \int_{\text{curve}} \tau^p ~dl.\]
J()

This returns the value of the quantity.

dJ()

This returns the derivative of the quantity with respect to the curve dofs.

simsopt.geo.curveobjectives.Lp_curvature_pure(kappa, gammadash, p, desired_kappa)

This function is used in a Python+Jax implementation of the curvature penalty term.

simsopt.geo.curveobjectives.Lp_torsion_pure(torsion, gammadash, p)

This function is used in a Python+Jax implementation of the formula for the torsion penalty term.

class simsopt.geo.curveobjectives.MinimumDistance(curves, minimum_distance)

Bases: simsopt._core.optimizable.Optimizable

MinimumDistance is a class that computes

\[J = \sum_{i = 1}^{\text{num_coils}} \sum_{j = 1}^{i-1} d_{i,j}\]

where

\[\begin{split}d_{i,j} = \int_{\text{curve}_i} \int_{\text{curve}_j} \max(0, d_{\min} - \| \mathbf{r}_i - \mathbf{r}_j \|_2)^2 ~dl_j ~dl_i\\\end{split}\]

and \(\mathbf{r}_i\), \(\mathbf{r}_j\) are points on coils \(i\) and \(j\), respectively. \(d_\min\) is a desired threshold minimum intercoil distance. This penalty term is zero when the points on coil \(i\) and coil \(j\) lie more than \(d_\min\) away from one another, for \(i, j \in \{1, \cdots, \text{num_coils}\}\)

J()

This returns the value of the quantity.

dJ()

This returns the derivative of the quantity with respect to the curve dofs.

simsopt.geo.curveobjectives.curve_length_pure(l)

This function is used in a Python+Jax implementation of the curve length formula.

simsopt.geo.curveobjectives.distance_pure(gamma1, l1, gamma2, l2, minimum_distance)

This function is used in a Python+Jax implementation of the distance formula.

simsopt.geo.curverzfourier module

class simsopt.geo.curverzfourier.CurveRZFourier(quadpoints, order, nfp, stellsym)

Bases: simsgeopp.CurveRZFourier, simsopt.geo.curve.Curve

CurveRZFourier is a curve that is represented in cylindrical

coordinates using the following Fourier series:

\[\begin{split}r(\phi) &= \sum_{m=0}^{\text{order}} x_{c,m}\cos(n_{\text{fp}} m \phi) + \sum_{m=1}^{\text{order}} x_{s,m}\sin(n_{\text{fp}} m \phi) \\ z(\phi) &= \sum_{m=0}^{\text{order}} z_{c,m}\cos(n_{\text{fp}} m \phi) + \sum_{m=1}^{\text{order}} z_{s,m}\sin(n_{\text{fp}} m \phi)\end{split}\]

If stellsym = True, then the \(\sin\) terms for \(r\) and the \(\cos\) terms for \(z\) are zero.

For the stellsym = False case, the dofs are stored in the order

\[[r_{c,0}, \cdots, r_{c,\text{order}}, r_{s,1}, \cdots, r_{s,\text{order}}, z_{c,0},....]\]

or in the stellsym = True case they are stored

\[[r_{c,0},...,r_{c,order},z_{s,1},...,z_{s,order}]\]
get_dofs()

This function returns the dofs associated to this object.

set_dofs(dofs)

This function sets the dofs associated to this object.

simsopt.geo.curvexyzfourier module

class simsopt.geo.curvexyzfourier.CurveXYZFourier(quadpoints, order)

Bases: simsgeopp.CurveXYZFourier, simsopt.geo.curve.Curve

CurveXYZFourier is a curve that is represented in Cartesian coordinates using the following Fourier series:

\[\begin{split}x(\phi) &= \sum_{m=0}^{\text{order}} x_{c,m}\cos(m\phi) + \sum_{m=1}^{\text{order}} x_{s,m}\sin(m\phi) \\ y(\phi) &= \sum_{m=0}^{\text{order}} y_{c,m}\cos(m\phi) + \sum_{m=1}^{\text{order}} y_{s,m}\sin(m\phi) \\ z(\phi) &= \sum_{m=0}^{\text{order}} z_{c,m}\cos(m\phi) + \sum_{m=1}^{\text{order}} z_{s,m}\sin(m\phi)\end{split}\]

The dofs are stored in the order

\[[x_{c,0}, x_{c,1}, \cdots, x_{c,\text{order}}, x_{s,1}, \cdots, x_{s,\text{order}}, y_{c,0}, \cdots]\]
get_dofs()

This function returns the dofs associated to this object.

static load_curves_from_file(filename, order=None, ppp=20, delimiter=',')

This function loads a file containing Fourier coefficients for several coils. The file is expected to have 6*num_coils many columns, and order+1 many rows. The columns are in the following order,

sin_x_coil1, cos_x_coil1, sin_y_coil1, cos_y_coil1, sin_z_coil1, cos_z_coil1, sin_x_coil2, cos_x_coil2, sin_y_coil2, cos_y_coil2, sin_z_coil2, cos_z_coil2, …

set_dofs(dofs)

This function sets the dofs associated to this object.

class simsopt.geo.curvexyzfourier.JaxCurveXYZFourier(quadpoints, order)

Bases: simsopt.geo.curve.JaxCurve

A Python+Jax implementation of the CurveXYZFourier class. There is actually no reason why one should use this over the C++ implementation in simsgeopp, but the point of this class is to illustrate how jax can be used to define a geometric object class and calculate all the derivatives (both with respect to dofs and with respect to the angle \(\phi\)) automatically.

get_dofs()

This function returns the dofs associated to this object.

num_dofs()

This function returns the number of dofs associated to this object.

set_dofs_impl(dofs)

This function sets the dofs associated to this object.

simsopt.geo.curvexyzfourier.jaxfouriercurve_pure(dofs, quadpoints, order)

simsopt.geo.graph_surface module

This module provides several classes for representing toroidal surfaces. There is a base class Surface, and several child classes corresponding to different discrete representations.

class simsopt.geo.graph_surface.Surface(nfp=1, stellsym=True)

Bases: abc.ABC

Surface is a base class for various representations of toroidal surfaces in simsopt.

_abc_impl = <_abc_data object>
abstract to_RZFourier()

Return a SurfaceRZFourier instance corresponding to the shape of this surface. All subclasses should implement this abstract method.

class simsopt.geo.graph_surface.SurfaceGarabedian(nfp=1, mmax=1, mmin=0, nmax=0, nmin=None)

Bases: simsopt.geo.graph_surface.Surface, simsopt._core.graph_optimizable.Optimizable

SurfaceGarabedian represents a toroidal surface for which the shape is parameterized using Garabedian’s Delta_{m,n} coefficients.

The present implementation assumes stellarator symmetry. Note that non-stellarator-symmetric surfaces require that the Delta_{m,n} coefficients be imaginary.

property Delta
_abc_impl = <_abc_data object>
_ids = count(1)
area()

Return the area of the surface.

area_volume()

Compute the surface area and the volume enclosed by the surface.

f()
get_Delta(m, n)

Return a particular Delta_{m,n} coefficient.

return_fn_map: Dict[str, Callable] = {'area': <function SurfaceGarabedian.area>, 'volume': <function SurfaceGarabedian.volume>}
set_Delta(m, n, val)

Set a particular Delta_{m,n} coefficient.

to_RZFourier()

Return a SurfaceRZFourier object with the identical shape.

For a derivation of the transformation here, see https://terpconnect.umd.edu/~mattland/assets/notes/toroidal_surface_parameterizations.pdf

volume()

Return the volume of the surface.

class simsopt.geo.graph_surface.SurfaceRZFourier(nfp=1, stellsym=True, mpol=1, ntor=0)

Bases: simsopt.geo.graph_surface.Surface, simsopt._core.graph_optimizable.Optimizable

SurfaceRZFourier is a surface that is represented in cylindrical coordinates using the following Fourier series:

r(theta, phi) = sum_{m=0}^{mpol} sum_{n=-ntor}^{ntor} [

r_{c,m,n} cos(m heta - n nfp phi) + r_{s,m,n} sin(m heta - n nfp phi) ]

and the same for z(theta, phi).

Here, (r, phi, z) are standard cylindrical coordinates, and theta is any poloidal angle.

_abc_impl = <_abc_data object>
_ids = count(1)
_validate_mn(m, n)

Check whether m and n are in the allowed range.

area()

Return the area of the surface.

area_volume()

Compute the surface area and the volume enclosed by the surface.

darea()

Return the area of the surface.

darea_volume()

Compute the derivative of the surface area and the volume enclosed by the surface.

dvolume()

Return the volume of the surface.

f()
classmethod from_focus(filename)

Read in a surface from a FOCUS-format file.

get_rc(m, n)

Return a particular rc Parameter.

get_rs(m, n)

Return a particular rs Parameter.

get_zc(m, n)

Return a particular zc Parameter.

get_zs(m, n)

Return a particular zs Parameter.

make_names(prefix)

Create the arrays for the rc, rs, zc, and zs coefficients.

property rc
return_fn_map: Dict[str, Callable] = {'area': <function SurfaceRZFourier.area>, 'volume': <function SurfaceRZFourier.volume>}
property rs
set_rc(m, n, val)

Set a particular rc Parameter.

set_rs(m, n, val)

Set a particular rs Parameter.

set_zc(m, n, val)

Set a particular zc Parameter.

set_zs(m, n, val)

Set a particular zs Parameter.

to_Garabedian()

Return a SurfaceGarabedian object with the identical shape.

For a derivation of the transformation here, see https://terpconnect.umd.edu/~mattland/assets/notes/toroidal_surface_parameterizations.pdf

to_RZFourier()

No conversion necessary.

volume()

Return the volume of the surface.

property zc
property zs
simsopt.geo.graph_surface.area_volume_pure(rc, rs, zc, zs, stellsym, nfp, mpol, ntor, ntheta, nphi)

Compute the area and volume of a surface. This pure function is designed for automatic differentiation.

simsopt.geo.graph_surface.dataset_area_volume(surface_rz_fourier, stellsym, nfp, mpol, ntor, ntheta, nphi)

Compute the area and volume of a surface. This pure function is designed for automatic differentiation.

simsopt.geo.graph_surface.jit_area_volume_pure(rc, rs, zc, zs, stellsym, nfp, mpol, ntor, ntheta, nphi)

Compute the area and volume of a surface. This pure function is designed for automatic differentiation.

simsopt.geo.graph_surface.jit_dataset_area_volume(surface_rz_fourier, stellsym, nfp, mpol, ntor, ntheta, nphi)

Compute the area and volume of a surface. This pure function is designed for automatic differentiation.

simsopt.geo.jit module

simsopt.geo.jit.jit(fun)

simsopt.geo.magneticfield module

class simsopt.geo.magneticfield.MagneticField

Bases: object

Generic class that represents any magnetic field for which each magnetic field class inherits.

A(compute_derivatives=0)
B(compute_derivatives=0)
clear_cached_properties()
compute(points, compute_derivatives=0)
compute_A(points, compute_derivatives=0)
d2A_by_dXdX(compute_derivatives=2)
d2B_by_dXdX(compute_derivatives=2)
dA_by_dX(compute_derivatives=1)
dB_by_dX(compute_derivatives=1)
set_points(points)
class simsopt.geo.magneticfield.MagneticFieldMultiply(scalar, Bfield)

Bases: simsopt.geo.magneticfield.MagneticField

Class used to multiply a magnetic field by a scalar. It takes as input a MagneticField class and a scalar and multiplies B, A and their derivatives by that value

compute(points, compute_derivatives=0)
compute_A(points, compute_derivatives=0)
set_points(points)
class simsopt.geo.magneticfield.MagneticFieldSum(Bfields)

Bases: simsopt.geo.magneticfield.MagneticField

Class used to sum two or more magnetic field together. It can either be called directly with a list of magnetic fields given as input and outputing another magnetic field with B, A and its derivatives added together or it can be called by summing magnetic fields classes as Bfield1 + Bfield1

compute(points, compute_derivatives=0)
compute_A(points, compute_derivatives=0)
set_points(points)

simsopt.geo.magneticfieldclasses module

class simsopt.geo.magneticfieldclasses.CircularCoil(r0=0.1, center=[0, 0, 0], I=159154.94309189534, normal=[0, 0])

Bases: simsopt.geo.magneticfield.MagneticField

Magnetic field created by a single circular coil evaluated using analytical functions, including complete elliptic integrals of the first and second kind. As inputs, it takes the radius of the coil (r0), its center, current (I) and its normal vector [either spherical angle components (normal=[theta,phi]) or (x,y,z) components of a vector (normal=[x,y,z])]).

Parameters
  • r0 – radius of the coil

  • center – point at the coil center

  • I – current of the coil in Ampere’s

  • normal – if list with two values treats it as spherical angles theta and phi of the normal vector to the plane of the coil centered at the coil center, if list with three values treats it a vector

compute(points, compute_derivatives=0)
compute_A(points, compute_derivatives=0)
class simsopt.geo.magneticfieldclasses.Dommaschk(mn=[[0, 0]], coeffs=[[0, 0]])

Bases: simsopt.geo.magneticfield.MagneticField

Vacuum magnetic field created by an explicit representation of the magnetic field scalar potential as proposed by W. Dommaschk (1986), Computer Physics Communications 40, 203-218. As inputs, it takes the arrays for the harmonics m, n and its corresponding coefficients.

Parameters
  • m – first harmonic array

  • n – second harmonic array

  • coeffs – coefficient for Vml for each of the ith index of the harmonics m and n

compute(points, compute_derivatives=0)
class simsopt.geo.magneticfieldclasses.Reiman(iota0=0.15, iota1=0.38, k=[6], epsilonk=[0.01], m0=1)

Bases: simsopt.geo.magneticfield.MagneticField

Magnetic field model in section 5 of Reiman and Greenside, Computer Physics Communications 43 (1986) 157—167. This field allows for an analytical expression of the magnetic island width that can be used for island optimization. However, the field is not completely physical as it does not have nested flux surfaces.

Parameters
  • iota0 – unperturbed rotational transform

  • iota1 – unperturbed global magnetic shear

  • k – integer array specifying the Fourier modes used

  • epsilonk – coefficient of the Fourier modes

  • m0 – toroidal symmetry parameter (normally m0=1)

compute(points, compute_derivatives=0)
class simsopt.geo.magneticfieldclasses.ScalarPotentialRZMagneticField(phi_str)

Bases: simsopt.geo.magneticfield.MagneticField

Vacuum magnetic field as a solution of B = grad(phi) where phi is the magnetic field scalar potential. It takes phi as an input string, which should contain an expression involving the standard cylindrical coordinates (R, phi, Z) Example: ScalarPotentialRZMagneticField(“2*phi”) yields a magnetic field B = grad(2*phi) = (0,2/R,0). Note: this function needs sympy.

Parameters

phi_str – string containing vacuum scalar potential expression as a function of R, Z and phi

compute(points, compute_derivatives=0)
class simsopt.geo.magneticfieldclasses.ToroidalField(R0, B0)

Bases: simsopt.geo.magneticfield.MagneticField

Magnetic field purely in the toroidal direction, that is, in the phi direction with (R,phi,Z) the standard cylindrical coordinates. Its modulus is given by B = B0*R0/R where R0 is the first input and B0 the second input to the function.

Parameters
  • B0 – modulus of the magnetic field at R0

  • R0 – radius of normalization

compute(points, compute_derivatives=0)
compute_A(points, compute_derivatives=0)

simsopt.geo.surface module

class simsopt.geo.surface.Surface

Bases: simsopt._core.optimizable.Optimizable

Surface is a base class for various representations of toroidal surfaces in simsopt.

A Surface is modelled as a function \(\Gamma:[0, 1] \times [0, 1] \to R^3\) and is evaluated at quadrature points \(\{\phi_1, \ldots, \phi_{n_\phi}\}\times\{\theta_1, \ldots, \theta_{n_\theta}\}\).

aspect_ratio()

Note: cylindrical coordinates are \((R, \phi, Z)\), where \(\phi \in [-\pi,\pi)\) and the angles that parametrize the surface are \((\varphi, \theta) \in [0,1)^2\) For a given surface, this function computes its aspect ratio using the VMEC definition:

\[AR = R_{\text{major}} / R_{\text{minor}}\]

where

\[\begin{split}R_{\text{minor}} &= \sqrt{ \overline{A} / \pi } \\ R_{\text{major}} &= \frac{V}{2 \pi^2 R_{\text{minor}}^2}\end{split}\]

and \(V\) is the volume enclosed by the surface, and \(\overline{A}\) is the average cross sectional area. The main difficult part of this calculation is the mean cross sectional area. This is given by the integral

\[\overline{A} = \frac{1}{2\pi} \int_{S_{\phi}} ~dS ~d\phi\]

where \(S_\phi\) is the cross section of the surface at the cylindrical angle \(\phi\). Note that \(\int_{S_\phi} ~dS\) can be rewritten as a line integral

\[\begin{split}\int_{S_\phi}~dS &= \int_{S_\phi} ~dR dZ \\ &= \int_{\partial S_\phi} [R,0] \cdot \mathbf n/\|\mathbf n\| ~dl \\ &= \int^1_{0} R \frac{\partial Z}{\partial \theta}~d\theta\end{split}\]

where \(\mathbf n = [n_R, n_Z] = [\partial Z/\partial \theta, -\partial R/\partial \theta]\) is the outward pointing normal.

Consider the surface in cylindrical coordinates terms of its angles \([R(\varphi,\theta), \phi(\varphi,\theta), Z(\varphi,\theta)]\). The boundary of the cross section \(\partial S_\phi\) is given by the points \(\theta\rightarrow[R(\varphi(\phi,\theta),\theta),\phi, Z(\varphi(\phi,\theta),\theta)]\) for fixed \(\phi\). The cross sectional area of \(S_\phi\) becomes

\[\int^{1}_{0} R(\varphi(\phi,\theta),\theta) \frac{\partial}{\partial \theta}[Z(\varphi(\phi,\theta),\theta)] ~d\theta\]

Now, substituting this into the formula for the mean cross sectional area, we have

\[\overline{A} = \frac{1}{2\pi}\int^{\pi}_{-\pi}\int^{1}_{0} R(\varphi(\phi,\theta),\theta) \frac{\partial}{\partial \theta}[Z(\varphi(\phi,\theta),\theta)] ~d\theta ~d\phi\]

Instead of integrating over cylindrical \(\phi\), let’s complete the change of variables and integrate over \(\varphi\) using the mapping:

\[[\phi,\theta] \leftarrow [\text{atan2}(y(\varphi,\theta), x(\varphi,\theta)), \theta]\]

After the change of variables, the integral becomes:

\[\overline{A} = \frac{1}{2\pi}\int^{1}_{0}\int^{1}_{0} R(\varphi,\theta) \left[\frac{\partial Z}{\partial \varphi} \frac{\partial \varphi}{d \theta} + \frac{\partial Z}{\partial \theta} \right] \text{det} J ~d\theta ~d\varphi\]

where \(\text{det}J\) is the determinant of the mapping’s Jacobian.

cross_section(phi, thetas=None)

This function takes in a cylindrical angle \(\phi\) and returns the cross section of the surface in that plane evaluated at thetas. This is done using the method of bisection.

This function assumes that the surface intersection with the plane is a single curve.

plot(ax=None, show=True, plot_normal=False, plot_derivative=False, scalars=None, wireframe=True)

Plot the surface using mayavi. Note: the ax and show parameter can be used to plot more than one surface:

ax = surface1.plot(show=False)
ax = surface2.plot(ax=ax, show=False)
surface3.plot(ax=ax, show=True)
to_RZFourier()

Return a simsopt.geo.surfacerzfourier.SurfaceRZFourier instance corresponding to the shape of this surface. All subclasses should implement this abstract method.

simsopt.geo.surfacegarabedian module

class simsopt.geo.surfacegarabedian.SurfaceGarabedian(nfp=1, mmax=1, mmin=0, nmax=0, nmin=None)

Bases: simsopt.geo.surface.Surface

SurfaceGarabedian represents a toroidal surface for which the shape is parameterized using Garabedian’s \(\Delta_{m,n}\) coefficients.

The present implementation assumes stellarator symmetry. Note that non-stellarator-symmetric surfaces require that the \(\Delta_{m,n}\) coefficients be imaginary.

allocate()

Create the array for the \(\Delta_{m,n}\) coefficients.

area()

Return the area of the surface.

area_volume()

Compute the surface area and the volume enclosed by the surface.

fixed_range(mmin, mmax, nmin, nmax, fixed=True)

Set the ‘fixed’ property for a range of m and n values.

All modes with m in the interval [mmin, mmax] and n in the interval [nmin, nmax] will have their fixed property set to the value of the ‘fixed’ parameter. Note that mmax and nmax are included (unlike the upper bound in python’s range(min, max).)

get_Delta(m, n)

Return a particular \(\Delta_{m,n}\) coefficient.

get_dofs()

Return a 1D numpy array with all the degrees of freedom.

set_Delta(m, n, val)

Set a particular \(\Delta_{m,n}\) coefficient.

set_dofs(v)

Set the shape coefficients from a 1D list/array

to_RZFourier()

Return a SurfaceRZFourier object with the identical shape.

For a derivation of the transformation here, see https://terpconnect.umd.edu/~mattland/assets/notes/toroidal_surface_parameterizations.pdf

volume()

Return the volume of the surface.

simsopt.geo.surfaceobjectives module

class simsopt.geo.surfaceobjectives.Area(surface)

Bases: object

Wrapper class for surface area label.

J()

Compute the area of a surface.

d2J_by_dsurfacecoefficientsdsurfacecoefficients()

Calculate the second derivatives with respect to the surface coefficients.

dJ_by_dsurfacecoefficients()

Calculate the derivatives with respect to the surface coefficients.

class simsopt.geo.surfaceobjectives.ToroidalFlux(surface, biotsavart, idx=0)

Bases: object

Given a surface and Biot Savart kernel, this objective calculates

\[\begin{split}J &= \int_{S_{\varphi}} \mathbf{B} \cdot \mathbf{n} ~ds, \\ &= \int_{S_{\varphi}} \text{curl} \mathbf{A} \cdot \mathbf{n} ~ds, \\ &= \int_{\partial S_{\varphi}} \mathbf{A} \cdot \mathbf{t}~dl,\end{split}\]

where \(S_{\varphi}\) is a surface of constant \(\varphi\), and \(\mathbf A\) is the magnetic vector potential.

J()

Compute the toroidal flux on the surface where \(\varphi = \texttt{quadpoints_varphi}[\texttt{idx}]\).

d2J_by_dsurfacecoefficientsdsurfacecoefficients()

Calculate the second derivatives with respect to the surface coefficients.

dJ_by_dsurfacecoefficients()

Calculate the derivatives with respect to the surface coefficients.

invalidate_cache()
class simsopt.geo.surfaceobjectives.Volume(surface)

Bases: object

Wrapper class for volume label.

J()

Compute the volume enclosed by the surface.

d2J_by_dsurfacecoefficientsdsurfacecoefficients()

Calculate the second derivatives with respect to the surface coefficients.

dJ_by_dsurfacecoefficients()

Calculate the derivatives with respect to the surface coefficients.

simsopt.geo.surfaceobjectives.boozer_surface_residual(surface, iota, G, biotsavart, derivatives=0)

For a given surface, this function computes the residual

\[G\mathbf B_\text{BS}(\mathbf x) - ||\mathbf B_\text{BS}(\mathbf x)||^2 (\mathbf x_\varphi + \iota \mathbf x_\theta)\]

as well as the derivatives of this residual with respect to surface dofs, iota, and G. In the above, \(\mathbf x\) are points on the surface, \(\iota\) is the rotational transform on that surface, and \(\mathbf B_{\text{BS}}\) is the magnetic field computed using the Biot-Savart law.

\(G\) is known for exact boozer surfaces, so if G=None is passed, then that value is used instead.

simsopt.geo.surfacerzfourier module

class simsopt.geo.surfacerzfourier.SurfaceRZFourier(nfp=1, stellsym=True, mpol=1, ntor=0, quadpoints_phi=63, quadpoints_theta=62)

Bases: simsgeopp.SurfaceRZFourier, simsopt.geo.surface.Surface

SurfaceRZFourier is a surface that is represented in cylindrical coordinates using the following Fourier series:

\[r(\theta, \phi) = \sum_{m=0}^{m_{\text{pol}}} \sum_{n=-n_{\text{tor}}}^{n_\text{tor}} [ r_{c,m,n} \cos(m \theta - n_{\text{fp}} n \phi) + r_{s,m,n} \sin(m \theta - n_{\text{fp}} n \phi) ]\]

and the same for \(z(\theta, \phi)\).

Here, \((r,\phi, z)\) are standard cylindrical coordinates, and theta is any poloidal angle.

Note that for \(m=0\) we skip the \(n<0\) term for the cos terms, and the \(n \leq 0\) for the sin terms.

In addition, in the stellsym=True case, we skip the sin terms for \(r\), and the cos terms for \(z\).

_validate_mn(m, n)

Check whether m and n are in the allowed range.

change_resolution(mpol, ntor)

Change the values of mpol and ntor. Any new Fourier amplitudes will have a magnitude of zero. Any previous nonzero Fourier amplitudes that are not within the new range will be discarded.

darea()

Short hand for Surface.darea_by_dcoeff()

dvolume()

Short hand for Surface.dvolume_by_dcoeff()

fixed_range(mmin, mmax, nmin, nmax, fixed=True)

Set the ‘fixed’ property for a range of m and n values.

All modes with m in the interval [mmin, mmax] and n in the interval [nmin, nmax] will have their fixed property set to the value of the fixed parameter. Note that mmax and nmax are included (unlike the upper bound in python’s range(min, max).)

classmethod from_focus(filename, nphi=32, ntheta=32)

Read in a surface from a FOCUS-format file.

get_dofs()

Return the dofs associated to this surface.

get_rc(m, n)

Return a particular rc Parameter.

get_rs(m, n)

Return a particular rs Parameter.

get_zc(m, n)

Return a particular zc Parameter.

get_zs(m, n)

Return a particular zs Parameter.

make_names()

Form a list of names of the rc, zs, rs, or zc array elements.

make_names_helper(prefix, include0)
set_dofs(self: simsgeopp.SurfaceRZFourier, arg0: List[float])None
set_rc(m, n, val)

Set a particular rc Parameter.

set_rs(m, n, val)

Set a particular rs Parameter.

set_zc(m, n, val)

Set a particular zc Parameter.

set_zs(m, n, val)

Set a particular zs Parameter.

to_Garabedian()

Return a SurfaceGarabedian object with the identical shape.

For a derivation of the transformation here, see https://terpconnect.umd.edu/~mattland/assets/notes/toroidal_surface_parameterizations.pdf

to_RZFourier()

No conversion necessary.

simsopt.geo.surfacexyzfourier module

class simsopt.geo.surfacexyzfourier.SurfaceXYZFourier(nfp=1, stellsym=True, mpol=1, ntor=0, quadpoints_phi=32, quadpoints_theta=32)

Bases: simsgeopp.SurfaceXYZFourier, simsopt.geo.surface.Surface

SurfaceXYZFourier is a surface that is represented in Cartesian coordinates using the following Fourier series:

\[\begin{split}\hat x(\phi,\theta) &= \sum_{m=0}^{m_\text{pol}} \sum_{n=-n_{\text{tor}}}^{n_{tor}} [ x_{c,m,n} \cos(m \theta - n_\text{ fp } n \phi) + x_{s,m,n} \sin(m \theta - n_\text{ fp } n \phi)]\\ \hat y(\phi,\theta) &= \sum_{m=0}^{m_\text{pol}} \sum_{n=-n_\text{tor}}^{n_\text{tor}} [ y_{c,m,n} \cos(m \theta - n_\text{fp} n \phi) + y_{s,m,n} \sin(m \theta - n_\text{fp} n \phi)]\\ z(\phi,\theta) &= \sum_{m=0}^{m_\text{pol}} \sum_{n=-n_\text{tor}}^{n_\text{tor}} [ z_{c,m,n} \cos(m \theta - n_\text{fp}n \phi) + z_{s,m,n} \sin(m \theta - n_\text{fp}n \phi)]\end{split}\]

where

\[\begin{split}x &= \hat x \cos(\phi) - \hat y \sin(\phi)\\ y &= \hat x \sin(\phi) + \hat y \cos(\phi)\end{split}\]

Note that for \(m=0\) we skip the \(n<0\) term for the cos terms, and the \(n \leq 0\) for the sin terms.

When enforcing stellarator symmetry, we set the

\[x_{s,*,*}, ~y_{c,*,*}, \text{and} ~z_{c,*,*}\]

terms to zero.

get_dofs()

Return the dofs associated to this surface.

set_dofs(dofs)

Set the dofs associated to this surface.

to_RZFourier()

Return a SurfaceRZFourier instance corresponding to the shape of this surface.

simsopt.geo.surfacexyztensorfourier module

class simsopt.geo.surfacexyztensorfourier.SurfaceXYZTensorFourier(nfp=1, stellsym=True, mpol=1, ntor=1, clamped_dims=[False, False, False], quadpoints_phi=32, quadpoints_theta=32)

Bases: simsgeopp.SurfaceXYZTensorFourier, simsopt.geo.surface.Surface

SurfaceXYZTensorFourier is a surface that is represented in cartesian coordinates using the following Fourier series:

\[\begin{split}\hat x(\theta, \phi) &= \sum_{i=0}^{2m_\text{pol}} \sum_{j=0}^{2n_\text{tor}} x_{ij} w_i(\theta)v_j(\phi)\\ \hat y(\theta, \phi) &= \sum_{i=0}^{2m_\text{pol}} \sum_{j=0}^{2n_\text{tor}} y_{ij} w_i(\theta)v_j(\phi)\\ x(\phi, \theta) &= \hat x(\phi, \theta) \cos(\phi) - \hat y(\phi, \theta) \sin(\phi)\\ y(\phi, \theta) &= \hat x(\phi, \theta) \sin(\phi) + \hat y(\phi, \theta) \cos(\phi)\\ z(\theta, \phi) &= \sum_{i=0}^{2m_\text{pol}} \sum_{j=0}^{2n_\text{tor}} z_{ij} w_i(\theta)v_j(\phi)\end{split}\]

where the basis functions \({v_j}\) are given by

\[\{1, \cos(1\,\mathrm{nfp}\,\phi), \ldots, \cos(n_\text{tor}\,\mathrm{nfp}\,\phi), \sin(1\,\mathrm{nfp}\,\phi), \ldots, \sin(n_\text{tor}\,\mathrm{nfp}\,\phi)\}\]

and \({w_i}\) are given by

\[\{1, \cos(1\theta), \ldots, \cos(m_\text{pol}\theta), \sin(1\theta), \ldots, \sin(m_\text{pol}\theta)\}\]

When stellsym=True the sums above change as follows:

\[\begin{split}\hat x(\theta, \phi) &= \sum_{i=0}^{m_\text{pol}} \sum_{j=0}^{n_\text{tor}} x_{ij} w_i(\theta)v_j(\phi) + \sum_{i=m_\text{pol}+1}^{2m_\text{pol}} \sum_{j=n_\text{tor}+1}^{2n_\text{tor}} x_{ij} w_i(\theta)v_j(\phi)\\ \hat y(\theta, \phi) &= \sum_{i=0}^{m_\text{pol}} \sum_{j=n_\text{tor}+1}^{2n_\text{tor}} y_{ij} w_i(\theta)v_j(\phi) + \sum_{i=m_\text{pol}+1}^{2m_\text{pol}} \sum_{j=0}^{n_\text{tor}} y_{ij} w_i(\theta)v_j(\phi)\\\\ z(\theta, \phi) &= \sum_{i=0}^{m_\text{pol}} \sum_{j=n_\text{tor}+1}^{2n_\text{tor}} z_{ij} w_i(\theta)v_j(\phi) + \sum_{i=m_\text{pol}+1}^{2m_\text{pol}} \sum_{j=0}^{n_\text{tor}} z_{ij} w_i(\theta)v_j(\phi)\end{split}\]
get_dofs()

Return the dofs associated to this surface.

set_dofs(dofs)

Set the dofs associated to this surface.

to_RZFourier()

Return a SurfaceRZFourier instance corresponding to the shape of this surface.

Module contents