simsopt.geo package

Submodules

simsopt.geo.biotsavart module

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

Bases: object

A(compute_derivatives=0)
B(compute_derivatives=0)
B_and_dB_vjp(v, vgrad)
B_vjp(v)
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)
d2B_by_dXdcoilcurrents(compute_derivatives=1)
d3B_by_dXdXdcoilcurrents(compute_derivatives=2)
dA_by_dX(compute_derivatives=1)
dB_by_dX(compute_derivatives=1)
dB_by_dcoilcurrents(compute_derivatives=0)
set_points(points)

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.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.jit module

simsopt.geo.jit.jit(fun)

simsopt.geo.objectives module

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)

Boozer coordinates are (varphi, theta)

For a given surface, this function computes its aspect ratio using the VMEC definition: AR = R_major / R_minor where R_major = (volume enclosed by surface) / (2 pi^2 * R_minor^2)

R_minor = sqrt[ (mean cross sectional area) / pi ]

The main difficult part of this calculation is the (mean cross sectional area). This is given by the integral (mean cross sectional area) = 1/(2*pi)int^{2 pi}_{0} int_{S_phi} dS dphi 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 using the divergence theorem int_{S_phi}dS = int_{S_phi} dR dZ

= int_{partial S_phi} nabla_{R,Z} cdot [R,0] cdot n dl where n = [n_R, n_Z] is the outward pointing normal = int_{partial S_phi} R * n_R dl

Consider the surface written in cylindrical coordinates in terms of its Boozer angles [R(varphi,theta), phi(varphi,theta), Z(varphi,theta)]. partial S_phi is given by the points theta->[R(varphi(phi,theta),theta), Z(varphi(phi,theta),theta)] for fixed phi. The cross sectional area of S_phi becomes

= int^{2pi}_{0} R(varphi(phi,theta),theta)

d/dtheta[Z(varphi(phi,theta),theta)] dtheta

Now, substituting this into the formula for the mean cross sectional area, we have 1/(2*pi)int^{2 pi}_{0}int^{2pi}_{0} R(varphi(phi,theta),theta)

d/dtheta[Z(varphi(phi,theta),theta)] dtheta dphi

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

[phi,theta] <- [atan2(y(varphi,theta), x(varphi,theta)), theta] After the change of variables, the integral becomes: 1/(2*pi)int^{2 pi}_{0}int^{2pi}_{0} R(varphi,theta) [dZ_dvarphi dvarphi_dtheta

  • dZ_dtheta ] detJ dtheta dvarphi

where detJ is the determinant of the mapping’s Jacobian.

cross_section(phi, varphi_resolution=None, theta_resolution=None)

This function takes in a cylindrical angle phi and returns the cross section of the surface in that plane. 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.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:

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

x_{c,m,n} cos(m theta - n nfp phi)

  • x_{s,m,n} sin(m theta - n nfp phi)

]

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

y_{c,m,n} cos(m theta - n nfp phi)

  • y_{s,m,n} sin(m theta - n nfp phi)

]

x = hat x * cos(phi) - hat y * sin(phi) y = hat x * sin(phi) + hat y * cos(phi)

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

z_{c,m,n} cos(m theta - n nfp phi) + z_{s,m,n} sin(m theta - n nfp phi)

]

Note that for m=0 we skip the n<0 term for the cos terms, and the n<=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.

Module contents