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 0.5*|| f ||^2_2, where f is the residual computed by boozer_surface_residual (see surfaceobjectives.py). This objective is zero when (phi,theta) that parametrize the surface correspond to Boozer angles.
The surface label can be area, volume, or toroidal flux. The surface is constrained by the user-provided targetlabel.
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
minimize_boozer_penalty_constraints_LBFGS minimize_boozer_penalty_constraints_newton 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
minimize_boozer_exact_constraints_newton
where Newton is used to solve the first order optimality condition.
-
boozer_exact_constraints
(xl, derivatives=0, optimize_G=True)¶ This function returns the optimality conditions corresponding to the minimisation problem
\[\text{min} \frac{1}{2} \| f(x) \|^2_2\]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. 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
\[r(x) = [f_1(x),...,f_n(x), \sqrt{w_c} (l-l_0), \sqrt{w_c} (z(\varphi=0, \theta=0) - 0)]\]where \(w_c\) is the constraint weight, \({f_i}_i\) are the Boozer residuals at quadrature points 1,…,n, \(l\) is the surface label, and \(l_0\) is the target surface label.
For
scalarized=False
, this function returns \(r(x)\) and optionally the Jacobian of \(r(x)\).for
scalarized=True
, this function returns\[g(x) = \frac{1}{2}r(x)^Tr(x),\]i.e. the least squares residual and optionally the gradient and the Hessian of \(g(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} \frac{1}{2} \| f(x) \|^2_2\]subject to
\[\begin{split}l - l_0 &= 0 \\ z(\varphi=0,\theta=0) - 0 &= 0\end{split}\]using Lagrange multipliers and Newton’s method. The final constraint is not necessary for stellarator symmetric surfaces as it is automatically satisfied by stellarator symmetric surfaces.
-
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} \frac{1}{2} \| f(x) \|^2_2 + \frac{1}{2} w_c (l - l_0)^2 + \frac{1}{2} c_w (z(\varphi=0, \theta=0) - 0)^2\]where \(|| f(x)||^2_2\) is the sum of squares of the Boozer residual at the quadrature points. 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 the above, but instead of LBFGS it uses a nonlinear least squares algorithm. Options for method are the same as for scipy.optimize.least_squares.
-
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 the above, 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
Given some input coils and currents, this performs the reflection and rotation to generate a full set of stellarator coils.
simsopt.geo.config module¶
simsopt.geo.curve module¶
-
class
simsopt.geo.curve.
Curve
¶ Bases:
simsopt._core.optimizable.Optimizable
-
dfrenet_frame_by_dcoeff
()¶
-
dincremental_arclength_by_dcoeff_vjp
(v)¶
-
dkappa_by_dcoeff_impl
(dkappa_by_dcoeff)¶
-
dkappa_by_dcoeff_vjp
(v)¶
-
dkappadash_by_dcoeff
()¶
-
dtorsion_by_dcoeff_impl
(dtorsion_by_dcoeff)¶
-
dtorsion_by_dcoeff_vjp
(v)¶
-
frenet_frame
()¶
-
kappa_impl
(kappa)¶
-
kappadash
()¶
-
plot
(ax=None, show=True, plot_derivative=False, closed_loop=True, color=None, linestyle=None)¶
-
plot_mayavi
(show=True)¶
-
torsion_impl
(torsion)¶
-
-
class
simsopt.geo.curve.
JaxCurve
(quadpoints, gamma_pure)¶ Bases:
simsgeopp.Curve
,simsopt.geo.curve.Curve
-
dgamma_by_dcoeff_impl
(dgamma_by_dcoeff)¶
-
dgamma_by_dcoeff_vjp
(self: simsgeopp.Curve, arg0: numpy.ndarray[numpy.float64]) → numpy.ndarray[numpy.float64]¶
-
dgammadash_by_dcoeff_impl
(dgammadash_by_dcoeff)¶
-
dgammadash_by_dcoeff_vjp
(self: simsgeopp.Curve, arg0: numpy.ndarray[numpy.float64]) → numpy.ndarray[numpy.float64]¶
-
dgammadashdash_by_dcoeff_impl
(dgammadashdash_by_dcoeff)¶
-
dgammadashdash_by_dcoeff_vjp
(self: simsgeopp.Curve, arg0: numpy.ndarray[numpy.float64]) → numpy.ndarray[numpy.float64]¶
-
dgammadashdashdash_by_dcoeff_impl
(dgammadashdashdash_by_dcoeff)¶
-
dgammadashdashdash_by_dcoeff_vjp
(self: simsgeopp.Curve, arg0: numpy.ndarray[numpy.float64]) → numpy.ndarray[numpy.float64]¶
-
dkappa_by_dcoeff_vjp
(v)¶
-
dtorsion_by_dcoeff_vjp
(v)¶
-
gamma_impl
(self: simsgeopp.Curve, arg0: numpy.ndarray[numpy.float64], arg1: numpy.ndarray[numpy.float64]) → None¶
-
gammadash_impl
(gammadash)¶
-
gammadashdash_impl
(gammadashdash)¶
-
gammadashdashdash_impl
(gammadashdashdash)¶
-
-
class
simsopt.geo.curve.
RotatedCurve
(curve, theta, flip)¶ Bases:
simsgeopp.Curve
,simsopt.geo.curve.Curve
-
dgamma_by_dcoeff_impl
(dgamma_by_dcoeff)¶
-
dgamma_by_dcoeff_vjp
(self: simsgeopp.Curve, arg0: numpy.ndarray[numpy.float64]) → numpy.ndarray[numpy.float64]¶
-
dgammadash_by_dcoeff_impl
(dgammadash_by_dcoeff)¶
-
dgammadash_by_dcoeff_vjp
(self: simsgeopp.Curve, arg0: numpy.ndarray[numpy.float64]) → numpy.ndarray[numpy.float64]¶
-
dgammadashdash_by_dcoeff_impl
(dgammadashdash_by_dcoeff)¶
-
dgammadashdash_by_dcoeff_vjp
(self: simsgeopp.Curve, arg0: numpy.ndarray[numpy.float64]) → numpy.ndarray[numpy.float64]¶
-
dgammadashdashdash_by_dcoeff_impl
(dgammadashdashdash_by_dcoeff)¶
-
dgammadashdashdash_by_dcoeff_vjp
(self: simsgeopp.Curve, arg0: numpy.ndarray[numpy.float64]) → numpy.ndarray[numpy.float64]¶
-
gamma_impl
(self: simsgeopp.Curve, arg0: numpy.ndarray[numpy.float64], arg1: numpy.ndarray[numpy.float64]) → None¶
-
gammadash_impl
(gammadash)¶
-
gammadashdash_impl
(gammadashdash)¶
-
gammadashdashdash_impl
(gammadashdashdash)¶
-
get_dofs
(self: simsgeopp.Curve) → List[float]¶
-
num_dofs
(self: simsgeopp.Curve) → int¶
-
set_dofs_impl
(d)¶
-
-
simsopt.geo.curve.
incremental_arclength_pure
(d1gamma)¶
-
simsopt.geo.curve.
incremental_arclength_vjp
(d1gamma, v)¶
-
simsopt.geo.curve.
kappa_pure
(d1gamma, d2gamma)¶
-
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)¶
-
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.curveobjectives module¶
-
class
simsopt.geo.curveobjectives.
CurveLength
(curve)¶ Bases:
simsopt._core.optimizable.Optimizable
-
J
()¶
-
dJ
()¶
-
-
class
simsopt.geo.curveobjectives.
LpCurveCurvature
(curve, p, desired_length=None)¶ Bases:
simsopt._core.optimizable.Optimizable
-
J
()¶
-
dJ
()¶
-
-
class
simsopt.geo.curveobjectives.
LpCurveTorsion
(curve, p)¶ Bases:
simsopt._core.optimizable.Optimizable
-
J
()¶
-
dJ
()¶
-
-
simsopt.geo.curveobjectives.
Lp_curvature_pure
(kappa, gammadash, p, desired_kappa)¶
-
simsopt.geo.curveobjectives.
Lp_torsion_pure
(torsion, gammadash, p)¶
-
class
simsopt.geo.curveobjectives.
MinimumDistance
(curves, minimum_distance)¶ Bases:
simsopt._core.optimizable.Optimizable
-
J
()¶
-
dJ
()¶
-
-
simsopt.geo.curveobjectives.
curve_length_pure
(l)¶
-
simsopt.geo.curveobjectives.
distance_pure
(gamma1, l1, gamma2, l2, minimum_distance)¶
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:
r(phi) = sum_{n=0}^{order} x_{c,n}cos(n*nfp*phi) + sum_{n=1}^order x_{s,n}sin(n*nfp*phi) z(phi) = sum_{n=0}^{order} z_{c,n}cos(n*nfp*phi) + sum_{n=1}^order z_{s,n}sin(n*nfp*phi)
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},…,r_{c,order},r_{s,1},…,r_{s,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
(self: simsgeopp.CurveRZFourier) → List[float]¶
-
set_dofs
(self: simsgeopp.CurveRZFourier, arg0: List[float]) → None¶
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:
x(phi) = sum_{m=0}^{order} x_{c,m}cos(m*phi) + sum_{m=1}^order x_{s,m}sin(m*phi) y(phi) = sum_{m=0}^{order} y_{c,m}cos(m*phi) + sum_{m=1}^order y_{s,m}sin(m*phi) z(phi) = sum_{m=0}^{order} z_{c,m}cos(m*phi) + sum_{m=1}^order z_{s,m}sin(m*phi)
The dofs are stored in the order
[x_{c,0},…,x_{c,order},x_{s,1},…,x_{s,order},y_{c,0},….]
-
get_dofs
(self: simsgeopp.CurveXYZFourier) → List[float]¶
-
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
(self: simsgeopp.CurveXYZFourier, arg0: List[float]) → None¶
-
-
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
(self: simsgeopp.Curve) → List[float]¶
-
num_dofs
(self: simsgeopp.Curve) → int¶
-
set_dofs_impl
(dofs)¶
-
-
simsopt.geo.curvexyzfourier.
jaxfouriercurve_pure
(dofs, quadpoints, order)¶
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.
-
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)¶
-
to_RZFourier
()¶ Return a 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 computation
-
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\).
-
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 computation
-
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 with points x on it, 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.
\(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}^{mpol} sum_{n=-ntor}^{ntor} [
r_{c,m,n} cos(m theta - n nfp phi) + r_{s,m,n} sin(m theta - n nfp 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<=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.
-
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
()¶
-
dvolume
()¶
-
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
(self: simsgeopp.SurfaceRZFourier) → List[float]¶
-
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{ nfp } \phi) + x_{s,m,n} \sin(m \theta - n_\text{ nfp } \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{nfp} \phi) + y_{s,m,n} \sin(m \theta - n_\text{nfp} \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{ nfp} \phi) + z_{s,m,n} \sin(m \theta - n_\text{ nfp} \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,*,*} and z_{c,*,*}
terms to zero.
-
get_dofs
(self: simsgeopp.SurfaceXYZFourier) → List[float]¶
-
set_dofs
(self: simsgeopp.SurfaceXYZFourier, arg0: List[float]) → None¶
-
to_RZFourier
()¶ Return a SurfaceRZFourier instance corresponding to the shape of this surface. All subclasses should implement this abstract method.
-
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
-
get_dofs
(self: simsgeopp.SurfaceXYZTensorFourier) → List[float]¶
-
set_dofs
(self: simsgeopp.SurfaceXYZTensorFourier, arg0: List[float]) → None¶
-
to_RZFourier
()¶ Return a SurfaceRZFourier instance corresponding to the shape of this surface. All subclasses should implement this abstract method.
-