vmecpp.simsopt_compat module

SIMSOPT compatibility layer for VMEC++.

class vmecpp.simsopt_compat.Vmec(filename, verbose=True, ntheta=50, nphi=50, range_surface='full torus', mpi=None, keep_all_files=False)

Bases: Optimizable

A SIMSOPT-compatible Python wrapper for VMEC++.

Based on the original SIMSOPT wrapper for VMEC, see https://github.com/hiddenSymmetries/simsopt/blob/master/src/simsopt/mhd/vmec.py.

Parameters:
aspect()

Return the plasma aspect ratio.

Return type:

float

property boundary: SurfaceRZFourier
ds: float | None
external_current()

Return the total electric current associated with external currents, i.e. the current through the “doughnut hole”. This number is useful for coil optimization, to know what the sum of the coil currents must be.

Return type:

float

Returns:

float with the total external electric current in Amperes.

free_boundary: bool
get_dofs()
Return type:

NDArray[Shape[*], (float16, float32, float64, float32, float64)]

get_input()

Generate a VMEC++ input file (in JSON format).

The JSON data will be returned as a string. To save a file, see the write_input() function.

Return type:

str

indata: VmecInput | None
input_file: str | None
iota_axis()

Return the rotational transform on axis.

Return type:

float

iota_edge()

Return the rotational transform at the boundary.

Return type:

float

iter: int
load_wout_from_outfile()

Load data from self.output_file into self.wout.

Return type:

None

mean_iota()

Return the mean rotational transform.

The average is taken over the normalized toroidal flux s.

Return type:

float

mean_shear()

Return an average magnetic shear, d(iota)/ds, where s is the normalized toroidal flux.

This is computed by fitting the rotational transform to a linear (plus constant) function in s. The slope of this fit function is returned.

Return type:

float

mpi: Optional[MpiPartition]
n_current: int
n_iota: int
n_pressure: int
need_to_run_code: bool
output_file: str | None
recompute_bell(parent=None)

Function to be called whenever new DOFs input is given or if the parent Optimizable’s data changed, so the output from the current Optimizable object is invalid.

This method gets called by various DOF setters. If only the local DOFs of an object are being set, the recompute_bell method is called in that object and also in the descendent objects that have a dependency on the object, whose local DOFs are being changed. If gloabl DOFs of an object are being set, the recompute_bell method is called in the object, ancestors of the object, as well as the descendents of the object.

Need to be implemented by classes that provide a dof_setter for external handling of DOFs.

Return type:

None

run(restart_from=None, max_threads=1)

Run VMEC if need_to_run_code is True.

The max_threads argument is not present in SIMSOPT’s original implementation as it is specific to VMEC++, which will spawn the corresponding number of OpenMP threads to parallelize execution. By default max_threads=1, so VMEC++ runs on a single thread. To automatically enable all threads, pass max_threads=None

Most optimization frameworks use multi-process parallelization for finite differencing, and we do not want to end up overcommitting the machine with NCPU processes running NCPU threads each – especially when OpenMP is involved, as OpenMP threads are generally bad at resource-sharing.

Parameters:
Return type:

None

runnable: bool
s_full_grid: NDArray[Shape[* ns], (float16, float32, float64, float32, float64)] | None
s_half_grid: NDArray[Shape[* nshalf], (float16, float32, float64, float32, float64)] | None
set_dofs(x)
Parameters:

x (list[float])

Return type:

None

set_indata()

Transfer data from simsopt objects to vmec.indata.

Presently, this function sets the boundary shape and magnetic axis shape. In the future, the input profiles will be set here as well. This data transfer is performed before writing a Vmec input file or running Vmec. The boundary surface object converted to SurfaceRZFourier is returned.

Return type:

None

set_mpol_ntor(new_mpol, new_ntor)
Parameters:
  • new_mpol (int)

  • new_ntor (int)

vacuum_well()

Compute a single number W that summarizes the vacuum magnetic well, given by the formula.

W = (dV/ds(s=0) - dV/ds(s=1)) / (dV/ds(s=0)

where dVds is the derivative of the flux surface volume with respect to the radial coordinate s. Positive values of W are favorable for stability to interchange modes. This formula for W is motivated by the fact that

d^2 V / d s^2 < 0

is favorable for stability. Integrating over s from 0 to 1 and normalizing gives the above formula for W. Notice that W is dimensionless, and it scales as the square of the minor radius. To compute dV/ds, we use

dV/ds = 4 * pi**2 * abs(sqrt(g)_{0,0})

where sqrt(g) is the Jacobian of (s, theta, phi) coordinates, computed by VMEC in the gmnc array, and _{0,0} indicates the m=n=0 Fourier component. Since gmnc is reported by VMEC on the half mesh, we extrapolate by half of a radial grid point to s = 0 and 1.

Return type:

float

verbose: bool
volume()

Return the volume inside the VMEC last closed flux surface.

Return type:

float

wout: VmecWOut | None
write_input(filename)

Write a VMEC++ input file (in JSON format).

To just get the result as a string without saving a file, see the get_input() function.

Parameters:

filename (str | Path)

Return type:

None