vmecpp package¶
- vmecpp.run(input, magnetic_field=None, *, max_threads=None, verbose=True, restart_from=None)¶
Run VMEC++ using the provided input. This is the main entrypoint for both fixed- and free-boundary calculations.
- Parameters:
input (
VmecInput
) – a VmecInput instance, corresponding to the contents of a classic VMEC input filemagnetic_field (
Optional
[MagneticFieldResponseTable
]) – if present, VMEC++ will pass the magnetic field object in memory instead of reading it from an mgrid file (only relevant in free-boundary runs).max_threads (
Optional
[int
]) – maximum number of threads that VMEC++ should spawn. The actual number might still be lower that this in case there are too few flux surfaces to keep these many threads busy. If None, a number of threads equal to the number of logical cores is used.verbose (
bool
) – if True, VMEC++ logs its progress to standard output.restart_from (
Optional
[VmecOutput
]) – if present, VMEC++ is initialized using the converged equilibrium from the provided VmecOutput. This can dramatically decrease the number of iterations to convergence when running VMEC++ on a configuration that is very similar to the restart_from equilibrium.
- Return type:
Example
>>> import vmecpp >>> path = "examples/data/solovev.json" >>> vmec_input = vmecpp.VmecInput.from_file(path) >>> output = vmecpp.run(vmec_input, verbose=False, max_threads=1) >>> round(output.wout.b0, 14) # Exact value may differ by C library 0.20333137113443
- class vmecpp.VmecInput(**data)¶
Bases:
BaseModelWithNumpy
The input to a VMEC++ run. Contains settings as well as the definition of the plasma boundary.
Python equivalent of a VMEC++ JSON input file or a classic INDATA file (e.g. “input.best”).
Deserialize from JSON and serialize to JSON using the usual pydantic methods:
model_validate_json
andmodel_dump_json
.- Parameters:
data (
Any
)
- model_config: ClassVar[ConfigDict] = {'arbitrary_types_allowed': True, 'ser_json_inf_nan': 'strings'}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- lasym: bool¶
Flag to indicate non-stellarator-symmetry.
False, assumes stellarator symmetry (only cosine/sine coefficients used).
True, (currently unsupported) allows for non-stellarator-symmetric terms.
- nfp: int¶
Number of toroidal field periods (=1 for Tokamak)
- mpol: int¶
Number of poloidal Fourier harmonics; m = 0, 1, …, (mpol-1)
- ntor: int¶
Number of toroidal Fourier harmonics; n = -ntor, -ntor+1, …, -1, 0, 1, …, ntor-1, ntor.
- ntheta: int¶
Number of poloidal grid points (ntheta >= 0).
Controls the poloidal resolution in real space. If 0, chosen automatically as minimally allowed. Must be at least 2*mpol + 6.
- nzeta: int¶
Number of toroidal grid points (nzeta >= 0).
Controls the toroidal resolution in real space. If 0, chosen automatically as minimally allowed. Must be at least 2*ntor + 4. We typically use use phi as the convention for the toroidal angle, the name nzeta is due to beckwards compatibility.
- ns_array: jt.Int[np.ndarray, 'num_grids']¶
Number of flux surfaces per multigrid step.
Each entry >= 3 and >= previous entry.
- ftol_array: jt.Float[np.ndarray, 'num_grids']¶
Requested force tolerance for convergence per multigrid step.
- niter_array: jt.Int[np.ndarray, 'num_grids']¶
Maximum number of iterations per multigrid step.
- phiedge: float¶
Total enclosed toroidal magnetic flux in Vs == Wb.
In fixed-boundary, this determines the magnetic field strength.
In free-boundary, the magnetic field strength is given externally, so this determines cross-section area and volume of the plasma.
- ncurr: typing.Literal[0, 1]¶
Select constraint on iota or enclosed toroidal current profiles.
0: constrained-iota (rotational transform profile specified)
1: constrained-current (toroidal current profile specified)
- pmass_type: ProfileType¶
Parametrization of mass/pressure profile.
- am: jt.Float[np.ndarray, 'am_len']¶
Mass/pressure profile coefficients.
Units: Pascals for pressure.
- am_aux_s: jt.Float[np.ndarray, 'am_aux_len']¶
knot locations in s
- Type:
Spline mass/pressure profile
- am_aux_f: jt.Float[np.ndarray, 'am_aux_len']¶
values at knots
- Type:
Spline mass/pressure profile
- pres_scale: float¶
Global scaling factor for mass/pressure profile.
- gamma: float¶
Adiabatic index \(\gamma\) (ratio of specific heats).
Specifying 0 implies that the pressure profile is specified. For all other values, the mass profile is specified.
- spres_ped: float¶
Location of pressure pedestal in s.
Outside this radial location, pressure is constant.
- piota_type: ProfileType¶
Parametrization of iota (rotational transform) profile.
- ai: jt.Float[np.ndarray, 'ai_len']¶
Iota profile coefficients.
- ai_aux_s: jt.Float[np.ndarray, 'ai_aux_len']¶
knot locations in s
- Type:
Spline iota profile
- ai_aux_f: jt.Float[np.ndarray, 'ai_aux_len']¶
values at knots
- Type:
Spline iota profile
- pcurr_type: ProfileType¶
Parametrization of toroidal current profile.
- ac: jt.Float[np.ndarray, 'ac_len']¶
Enclosed toroidal current profile coefficients.
- ac_aux_s: jt.Float[np.ndarray, 'ac_aux_len']¶
knot locations in s
- Type:
Spline toroidal current profile
- ac_aux_f: jt.Float[np.ndarray, 'ac_aux_len']¶
values at knots
- Type:
Spline toroidal current profile
- curtor: float¶
Net toroidal current in A.
The toroidal current profile is scaled to yield this total.
- bloat: float¶
Bloating factor (for constrained toroidal current)
- lfreeb: bool¶
Flag to indicate free-boundary.
If True, run in free-boundary mode; if False, fixed-boundary.
- mgrid_file: typing.Annotated[str, pydantic.Field(max_length=200)]¶
Full path for vacuum Green’s function data.
NetCDF MGRID file with magnetic field response factors for external coils.
- extcur: jt.Float[np.ndarray, 'ext_current']¶
Coil currents in A.
- nvacskip: int¶
Number of iterations between full vacuum calculations.
- nstep: int¶
Printout interval at which convergence progress is logged.
- aphi: jt.Float[np.ndarray, 'aphi_len']¶
Radial flux zoning profile coefficients.
- delt: float¶
Initial value for artificial time step in iterative solver.
- tcon0: float¶
Constraint force scaling factor for ns –> 0.
- lforbal: bool¶
directly compute innermost flux surface geometry from radial force balance
- Type:
Hack
- return_outputs_even_if_not_converged: bool¶
If true, return the outputs even if VMEC++ did not converge.
Otherwise a RuntimeError will be raised.
- raxis_c: jt.Float[np.ndarray, 'ntor_plus_1']¶
Magnetic axis coefficients for R ~ cos(n*v); stellarator-symmetric.
At least 1 value required, up to n=ntor considered.
- zaxis_s: jt.Float[np.ndarray, 'ntor_plus_1']¶
Magnetic axis coefficients for Z ~ sin(n*v); stellarator-symmetric.
Up to n=ntor considered; first entry (n=0) is ignored.
- raxis_s: jt.Float[np.ndarray, 'ntor_plus_1'] | None¶
Magnetic axis coefficients for R ~ sin(n*v); non-stellarator-symmetric.
Up to n=ntor considered; first entry (n=0) is ignored. Only used if lasym=True.
- zaxis_c: jt.Float[np.ndarray, 'ntor_plus_1'] | None¶
Magnetic axis coefficients for Z ~ cos(n*v); non-stellarator-symmetric.
Only used if lasym=True.
- rbc: SerializableSparseCoefficientArray[jt.Float[np.ndarray, 'mpol two_ntor_plus_one']]¶
Boundary coefficients for R ~ cos(m*u - n*v); stellarator-symmetric
- zbs: SerializableSparseCoefficientArray[jt.Float[np.ndarray, 'mpol two_ntor_plus_one']]¶
Boundary coefficients for Z ~ sin(m*u - n*v); stellarator-symmetric
- rbs: SerializableSparseCoefficientArray[jt.Float[np.ndarray, 'mpol two_ntor_plus_one']] | None¶
Boundary coefficients for R ~ sin(m*u - n*v); non-stellarator-symmetric.
Only used if lasym=True.
- zbc: SerializableSparseCoefficientArray[jt.Float[np.ndarray, 'mpol two_ntor_plus_one']] | None¶
Boundary coefficients for Z ~ cos(m*u - n*v); non-stellarator-symmetric.
Only used if lasym=True.
- static resize_2d_coeff(coeff, mpol_new, ntor_new)¶
Resizes a 2D NumPy array representing Fourier coefficients, padding with zeros or truncating as needed.
- Parameters:
- Return type:
Float[ndarray, 'mpol_new two_ntor_new_plus_one']
Examples
>>> coeff = np.array([[1, 2, 3, 4, 5], [6, 7, 8, 9, 10]]) >>> VmecInput.resize_2d_coeff(coeff, 3, 3) array([[ 0., 1., 2., 3., 4., 5., 0.], [ 0., 6., 7., 8., 9., 10., 0.], [ 0., 0., 0., 0., 0., 0., 0.]])
>>> VmecInput.resize_2d_coeff(coeff, 1, 1) array([[2., 3., 4.]])
>>> VmecInput.resize_2d_coeff(coeff, 4, 1) array([[2., 3., 4.], [7., 8., 9.], [0., 0., 0.], [0., 0., 0.]])
- static from_file(input_file)¶
Build a VmecInput from either a VMEC++ JSON input file or a classic INDATA file.
- static default()¶
Return a
VmecInput
with VMEC++ default values.
- to_json(**kwargs)¶
Serialize the object to JSON.
- Keyword Arguments:
**kwargs – Additional keyword arguments to forward to the model_dump_json method.
- Return type:
- class vmecpp.VmecOutput(**data)¶
Bases:
BaseModelWithNumpy
Container for the full output of a VMEC run.
- Parameters:
data (
Any
)
- model_config: ClassVar[ConfigDict] = {'arbitrary_types_allowed': True, 'ser_json_inf_nan': 'strings'}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- input: VmecInput¶
The input to the VMEC run that produced this output.
- jxbout: JxBOut¶
Python equivalent of VMEC’s “jxbout” file.
- mercier: Mercier¶
Python equivalent of VMEC’s “mercier” file.
Contains radial profiles and stability criteria relevant for Mercier stability analysis, including the Mercier criterion and its decomposition into shear, well, current, and geodesic contributions. Also includes profiles of rotational transform, toroidal flux, pressure, and their derivatives.
- threed1_volumetrics: Threed1Volumetrics¶
Python equivalent of VMEC’s volumetrics section in the “threed1” file.
Contains global and flux-surface-averaged quantities such as total and average pressure, poloidal and toroidal magnetic field energies, kinetic energy, and related integrals. Useful for postprocessing and global equilibrium characterization.
- wout: VmecWOut¶
Python equivalent of VMEC’s “wout” file.
- class vmecpp.VmecWOut(**data)¶
Bases:
BaseModelWithNumpy
Python equivalent of a VMEC “wout file”.
VmecWOut exposes the layout that SIMSOPT expects. The
save
method produces a NetCDF file compatible with SIMSOPT/Fortran VMECwout.nc
.- Parameters:
data (
Any
)
- model_config: ClassVar[ConfigDict] = {'arbitrary_types_allowed': True, 'extra': 'allow', 'ser_json_inf_nan': 'strings', 'serialize_by_alias': False, 'validate_by_alias': False, 'validate_by_name': True}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- input_extension: typing.Annotated[str, pydantic.Field(max_length=100)]¶
File extension of the input file.
- ier_flag: int¶
Status code indicating success or problems during the VMEC++ run.
See the
reason
property for a human-readable description.
- nfp: int¶
Number of toroidal field periods.
- ns: int¶
Number of radial grid points (=number of flux surfaces).
- mpol: int¶
Number of poloidal Fourier modes.
- ntor: int¶
Number of toroidal Fourier modes.
- mnmax: int¶
Number of Fourier coefficients for the state vector.
- mnmax_nyq: int¶
Number of Fourier coefficients for the Nyquist-quantities.
- lasym: typing.Annotated[bool, pydantic.PlainSerializer(lambda x: int(x)), pydantic.BeforeValidator(lambda x: bool(x)), pydantic.Field(alias='lasym__logical__')]¶
Flag indicating non-stellarator-symmetry.
Non-stellarator symmetric fields are only populated if this is True.
- lfreeb: typing.Annotated[bool, pydantic.PlainSerializer(lambda x: int(x)), pydantic.BeforeValidator(lambda x: bool(x)), pydantic.Field(alias='lfreeb__logical__')]¶
Flag indicating free-boundary computation.
- wb: float¶
volume integral of |B|^2/(2 mu0).
- Type:
Magnetic energy
- wp: float¶
volume integral of p.
- Type:
Kinetic energy
- rmax_surf: float¶
Maximum
R
on the plasma boundary over all grid points.
- rmin_surf: float¶
Minimum
R
on the plasma boundary over all grid points.
- zmax_surf: float¶
Maximum
Z
on the plasma boundary over all grid points.
- aspect: float¶
Aspect ratio (major radius over minor radius) of the plasma boundary.
- betapol: float¶
Poloidal plasma beta.
The ratio of the total thermal energy of the plasma to the total poloidal magnetic energy. \(\beta = W_{th} / W_{B_\theta} = \int p\, dV / \left( \int B_\theta^2 / (2 \mu_0)\, dV \right )\)
- betator: float¶
Toroidal plasma beta.
The ratio of the total thermal energy of the plasma to the total toroidal magnetic energy. \(\beta = W_{th} / W_{B_\phi} = \int p\, dV / \left( \int B_\phi^2 / (2 \mu_0)\, dV \right )\)
- betaxis: float¶
Plasma beta on the magnetic axis.
- b0: float¶
Toroidal magnetic flux density from poloidal current and magnetic axis position at
phi=0
.
- rbtor0: float¶
Poloidal ribbon current at the axis.
- rbtor: float¶
Poloidal ribbon current at the plasma boundary.
- IonLarmor: float¶
Larmor radius of plasma ions.
- ctor: float¶
Net toroidal plasma current.
- Aminor_p: float¶
Minor radius of the plasma.
- Rmajor_p: float¶
Major radius of the plasma.
- volume: typing.Annotated[float, pydantic.Field(alias='volume_p')]¶
Plasma volume.
- fsqr: float¶
Invariant force residual of the force on
R
at end of the run.
- fsqz: float¶
Invariant force residual of the force on
Z
at end of the run.
- fsql: float¶
Invariant force residual of the force on
lambda
at end of the run.
- ftolv: float¶
Force tolerance value used to determine convergence.
- itfsq: int¶
Number of force-balance iterations after which the run terminated.
- phipf: jt.Float[np.ndarray, 'n_surfaces']¶
Radial derivative of enclosed toroidal magnetic flux
phi'
on the full- grid.
- chipf: jt.Float[np.ndarray, 'n_surfaces']¶
Radial derivative of enclosed poloidal magnetic flux
chi'
on the full- grid.
- jcuru: jt.Float[np.ndarray, 'n_surfaces']¶
Radial derivative of enclosed poloidal current on full-grid.
- jcurv: jt.Float[np.ndarray, 'n_surfaces']¶
Radial derivative of enclosed toroidal current on full-grid.
- fsqt: jt.Float[np.ndarray, 'time']¶
Evolution of the total force residual along the run.
This is the sum of
force_residual_r
,force_residual_z
, andforce_residual_lambda
.
- force_residual_r: jt.Float[np.ndarray, 'time']¶
Evolution of the r radial force residual along the run.
- force_residual_z: jt.Float[np.ndarray, 'time']¶
Evolution of the z vertical force residual along the run.
- force_residual_lambda: jt.Float[np.ndarray, 'time']¶
Evolution of the lambda force residual along the run.
- delbsq: jt.Float[np.ndarray, 'time']¶
Evolution of the force residual at the vacuum boundary along the run.
- restart_reason_timetrace: typing.Annotated[jt.Int[np.ndarray, 'time'], pydantic.Field(alias='restart_reasons'), pydantic.BeforeValidator(lambda x: np.array(x).astype(np.int64))]¶
Internal restart reasons at each step along the run. (debugging quantity).
Use the
restart_reasons
field to access a more readable enum version of this instead of integer status codes.
- wdot: jt.Float[np.ndarray, 'time']¶
Evolution of the MHD energy decay along the run.
- jdotb: jt.Float[np.ndarray, 'n_surfaces']¶
Flux-surface-averaged \(\langle j \cdot B \rangle\) on full-grid.
- bdotb: jt.Float[np.ndarray, 'n_surfaces']¶
Flux-surface-averaged \(\langle B \cdot B \rangle\) on full-grid.
- bdotgradv: jt.Float[np.ndarray, 'n_surfaces']¶
Flux-surface-averaged toroidal magnetic field component \(B \cdot \nabla v\) on full-grid.
- DMerc: jt.Float[np.ndarray, 'n_surfaces']¶
Full Mercier stability criterion on the full-grid.
- equif: jt.Float[np.ndarray, 'n_surfaces']¶
Radial force balance residual on full-grid.
- xm: SerializeIntAsFloat[jt.Int[np.ndarray, 'mn_mode']]¶
Poloidal mode numbers
m
for the Fourier coefficients in the state vector.
- xn: SerializeIntAsFloat[jt.Int[np.ndarray, 'mn_mode']]¶
Toroidal mode numbers times number of toroidal field periods
n * nfp
for the Fourier coefficients in the state vector.
- xm_nyq: SerializeIntAsFloat[jt.Int[np.ndarray, 'mn_mode_nyq']]¶
Poloidal mode numbers
m
for the Fourier coefficients in the Nyquist- quantities.
- xn_nyq: SerializeIntAsFloat[jt.Int[np.ndarray, 'mn_mode_nyq']]¶
Toroidal mode numbers times number of toroidal field periods
n * nfp
for the Fourier coefficients in the Nyquist-quantities.
- mass: jt.Float[np.ndarray, 'n_surfaces']¶
Plasma mass profile
m
on half-grid.
- buco: jt.Float[np.ndarray, 'n_surfaces']¶
Profile of enclosed toroidal current
I
on half-grid.
- bvco: jt.Float[np.ndarray, 'n_surfaces']¶
Profile of enclosed poloidal ribbon current
G
on half-grid.
- phips: jt.Float[np.ndarray, 'n_surfaces']¶
Radial derivative of enclosed toroidal magnetic flux
phi'
on the half- grid.
- bmnc: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces']¶
Fourier coefficients (cos) of the magnetic field strength
|B|
on the half- grid.
- gmnc: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces']¶
Fourier coefficients (cos) of the Jacobian \(\sqrt{g}\) on the half-grid.
- bsubumnc: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces']¶
Fourier coefficients (cos) of the covariant magnetic field component \(B_{\theta}\) on the half-grid.
- bsubvmnc: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces']¶
Fourier coefficients (cos) of the covariant magnetic field component \(B_{\phi}\) on the half-grid.
- bsubsmns: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces']¶
Fourier coefficients (sin) of the covariant magnetic field component \(B_{s}\) on the full- grid.
- bsupumnc: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces']¶
Fourier coefficients (cos) of the contravariant magnetic field component \(B^{\theta}\) on the half-grid.
- bsupvmnc: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces']¶
Fourier coefficients (cos) of the contravariant magnetic field component \(B^{\phi}\) on the half-grid.
- rmnc: jt.Float[np.ndarray, 'mn_mode n_surfaces']¶
Fourier coefficients (cos) for
R
of the geometry of the flux surfaces on the full- grid.
- zmns: jt.Float[np.ndarray, 'mn_mode n_surfaces']¶
Fourier coefficients (sin) for
Z
of the geometry of the flux surfaces on the full- grid.
- lmns: jt.Float[np.ndarray, 'mn_mode n_surfaces']¶
Fourier coefficients (sin) for
lambda
stream function on the half-grid.
- lmns_full: jt.Float[np.ndarray, 'mn_mode n_surfaces']¶
Fourier coefficients (sin) for
lambda
stream function on the full-grid.This quantity is VMEC++ specific and required for hot-restart to work properly. We store it with the Fortran convention for the order of the dimensions for consistency with lmns.
- rmns: jt.Float[np.ndarray, 'mn_mode n_surfaces'] | None¶
Fourier coefficients (sin) for R of the geometry of the flux surfaces on the full-grid; non-stellarator-symmetric.
- zmnc: jt.Float[np.ndarray, 'mn_mode n_surfaces'] | None¶
Fourier coefficients (cos) for Z of the geometry of the flux surfaces on the full-grid; non-stellarator-symmetric.
- lmnc: jt.Float[np.ndarray, 'mn_mode n_surfaces'] | None¶
Fourier coefficients (cos) for lambda stream function on the half-grid; non- stellarator-symmetric.
- lmnc_full: jt.Float[np.ndarray, 'mn_mode n_surfaces'] | None¶
Fourier coefficients (cos) for lambda stream function on the full-grid; non- stellarator-symmetric.
This quantity is VMEC++ specific and required for hot-restart to work properly. We store it with the Fortran convention for the order of the dimensions for consistency with lmnc. Only populated when lasym=True (non-stellarator-symmetric configurations).
- gmns: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces'] | None¶
Fourier coefficients (sin) of the Jacobian \(\sqrt{g}\) on the half-grid; non-stellarator-symmetric.
- bmns: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces'] | None¶
Fourier coefficients (sin) of the magnetic field strength
|B|
on the half- grid; non-stellarator-symmetric.
- bsubumns: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces'] | None¶
Fourier coefficients (sin) of the covariant magnetic field component \(B_{\theta}\) on the half-grid; non-stellarator-symmetric.
- bsubvmns: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces'] | None¶
Fourier coefficients (sin) of the covariant magnetic field component \(B_{\phi}\) on the half-grid; non-stellarator-symmetric.
- bsubsmnc: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces'] | None¶
Fourier coefficients (cos) of the covariant magnetic field component \(B_{s}\) on the full- grid; non-stellarator-symmetric.
- bsupumns: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces'] | None¶
Fourier coefficients (sin) of the contravariant magnetic field component \(B^{\theta}\) on the half-grid; non-stellarator-symmetric.
- bsupvmns: jt.Float[np.ndarray, 'mn_mode_nyq n_surfaces'] | None¶
Fourier coefficients (sin) of the contravariant magnetic field component \(B^{\phi}\) on the half-grid; non-stellarator-symmetric.
- pcurr_type: ProfileType¶
Parametrization of toroidal current profile (copied from input).
- pmass_type: ProfileType¶
Parametrization of mass/pressure profile (copied from input).
- piota_type: ProfileType¶
Parametrization of iota profile (copied from input).
- am: jt.Float[np.ndarray, 'preset']¶
Mass/pressure profile coefficients (copied from input).
- ac: jt.Float[np.ndarray, 'preset']¶
Enclosed toroidal current profile coefficients (copied from input).
- ai: jt.Float[np.ndarray, 'preset']¶
Iota profile coefficients (copied from input).
- am_aux_s: AuxSType[jt.Float[np.ndarray, 'ndfmax']]¶
knot locations in
s
(copied from input).- Type:
Spline mass/pressure profile
- am_aux_f: AuxFType[jt.Float[np.ndarray, 'ndfmax']]¶
values at knots (copied from input).
- Type:
Spline mass/pressure profile
- ac_aux_s: AuxSType[jt.Float[np.ndarray, 'ndfmax']]¶
knot locations in
s
(copied from input).- Type:
Spline toroidal current profile
- ac_aux_f: AuxFType[jt.Float[np.ndarray, 'ndfmax']]¶
values at knots (copied from input).
- Type:
Spline toroidal current profile
- ai_aux_s: AuxSType[jt.Float[np.ndarray, 'ndfmax']]¶
knot locations in
s
(copied from input).- Type:
Spline iota profile
- ai_aux_f: AuxFType[jt.Float[np.ndarray, 'ndfmax']]¶
values at knots (copied from input).
- Type:
Spline iota profile
- gamma: float¶
Adiabatic index \(\gamma\) (copied from input).
- mgrid_file: typing.Annotated[str, pydantic.Field(max_length=200)]¶
Full path for vacuum Green’s function data (copied from input).
- nextcur: int¶
Number of external coil currents.
- extcur: typing.Annotated[jt.Float[np.ndarray, 'ext_current'], pydantic.BeforeValidator(lambda x: x if np.shape(x) != () else np.array([])), pydantic.PlainSerializer(lambda x: x if np.shape(x) != (0,) else netCDF4.default_fillvals['f8'])]¶
Coil currents in A.
for free-boundary runs,
extcur
has shape (nextcur,) for fixed-boundary it is a scalar float extcur=nan
- mgrid_mode: MgridModeType¶
Indicates if the mgrid file was normalized to unit currents (“S”) or not (“R”).
- iotas: jt.Float[np.ndarray, 'n_surfaces']¶
Rotational transform \(\iota\) on the half-grid.
- iotaf: jt.Float[np.ndarray, 'n_surfaces']¶
Rotational transform \(\iota\) on the full-grid.
- betatotal: float¶
Total plasma beta.
The ratio of the total thermal energy of the plasma to the total magnetic energy.
\(\beta = W_{th} / W_B = \int p\, dV / \left( \int B^2 / (2 \mu_0)\, dV \right )\)
- raxis_cc: jt.Float[np.ndarray, 'ntor_plus_1']¶
Fourier coefficients of \(R(phi)\) of the magnetic axis geometry.
- zaxis_cs: jt.Float[np.ndarray, 'ntor_plus_1']¶
Fourier coefficients of \(Z(phi)\) of the magnetic axis geometry.
- raxis_cs: jt.Float[np.ndarray, 'ntor_plus_1'] | None¶
Fourier coefficients of \(R(phi)\) of the magnetic axis geometry; non- stellarator-symmetric.
- zaxis_cc: jt.Float[np.ndarray, 'ntor_plus_1'] | None¶
Fourier coefficients of \(Z(phi)\) of the magnetic axis geometry; non- stellarator-symmetric.
- vp: jt.Float[np.ndarray, 'n_surfaces']¶
Differential volume \(V' = \frac{\partial V}{\partial s}\) on half-grid.
Note: called
dVds
in cpp
- presf: jt.Float[np.ndarray, 'n_surfaces']¶
Kinetic pressure
p
on the full-grid.
- pres: jt.Float[np.ndarray, 'n_surfaces']¶
Kinetic pressure
p
on the half-grid.
- phi: jt.Float[np.ndarray, 'n_surfaces']¶
Enclosed toroidal magnetic flux \(\phi\) on the full-grid.
- signgs: int¶
Sign of the Jacobian of the coordinate transform between flux coordinates and cylindrical coordinates.
- volavgB: float¶
Volume-averaged magnetic field strength.
- q_factor: jt.Float[np.ndarray, 'n_surfaces']¶
Safety factor \(q = 1/\iota\) on the full-grid.
- chi: jt.Float[np.ndarray, 'n_surfaces']¶
Enclosed poloidal magnetic flux \(\chi\) on the full-grid.
- specw: jt.Float[np.ndarray, 'n_surfaces']¶
Spectral width
M
on the full-grid.
- over_r: jt.Float[np.ndarray, 'n_surfaces']¶
<\tau / R> / V'
on half-grid.\(\left\langle \frac{\tau}{R} \right\rangle / V'\)
- DShear: jt.Float[np.ndarray, 'n_surfaces']¶
Mercier stability criterion contribution due to magnetic shear.
- DWell: jt.Float[np.ndarray, 'n_surfaces']¶
Mercier stability criterion contribution due to magnetic well.
- DCurr: jt.Float[np.ndarray, 'n_surfaces']¶
Mercier stability criterion contribution due to plasma currents.
- DGeod: jt.Float[np.ndarray, 'n_surfaces']¶
Mercier stability criterion contribution due to geodesic curvature.
- niter: int¶
Maximum number of force-balance iterations allowed.
- beta_vol: jt.Float[np.ndarray, 'n_surfaces']¶
Flux-surface averaged plasma beta on half-grid.
- version_: float¶
Version number of VMEC, that this VMEC++ wout file is compatible with.
Some codes change how they interpret values in the wout file depending on this number. (E.g. COBRAVMEC checks if >6 or not)
- property volume_p¶
The attribute is called volume_p in the Fortran wout file, while simsopt.mhd.Vmec.wout uses volume.
We expose both.
- property lasym__logical__¶
This is how the attribute is called in the Fortran wout file.
- property lfreeb__logical__¶
This is how the attribute is called in the Fortran wout file.
- property restart_reasons: list[tuple[int, RestartReason]]¶
Get the restart reasons as a list of tuples.
Each tuple contains the iteration number and the reason for the restart.
- save(out_path)¶
Save contents in NetCDF3 format, e.g.
wout.nc
.This is the format used by Fortran VMEC implementations and the one expected by SIMSOPT.
- static from_wout_file(wout_filename)¶
Load wout contents in NetCDF format.
This is the format used by Fortran VMEC implementations and the one expected by SIMSOPT. We allow for additional attributes to be present in the file, for compatibility with wouf files from other VMEC versions, but require at least the fields produced by VMEC++.
- class vmecpp.JxBOut(**data)¶
Bases:
BaseModelWithNumpy
- Parameters:
data (
Any
)
- model_config: ClassVar[ConfigDict] = {'arbitrary_types_allowed': True, 'extra': 'forbid', 'ser_json_inf_nan': 'strings'}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- itheta: jt.Float[np.ndarray, 'num_half nZnT']¶
Poloidal surface current.
\(itheta = (\frac{\partial B_s}{\partial \Phi} - \frac{\partial B_\phi}{\partial s}) / \mu_0\)
- izeta: jt.Float[np.ndarray, 'num_half nZnT']¶
Toroidal surface current.
\(izeta = (-\frac{\partial B_s}{\partial \Theta} + \frac{\partial B_\theta}{\partial s}) / \mu_0\)
- bdotk: jt.Float[np.ndarray, 'num_full nZnT']¶
- amaxfor: jt.Float[np.ndarray, 'n_surfaces']¶
100 times the maximum value of the real space force residual on each radial surface.
- aminfor: jt.Float[np.ndarray, 'n_surfaces']¶
100 times the minimum value of the real space force residual on each radial surface.
- avforce: jt.Float[np.ndarray, 'n_surfaces']¶
Average force residual on each radial surface.
- pprim: jt.Float[np.ndarray, 'n_surfaces']¶
Radial derivative of the pressure profile.
- jdotb: jt.Float[np.ndarray, 'n_surfaces']¶
Flux-surface-averaged \(\langle j \cdot B \rangle\) on full-grid.
- bdotb: jt.Float[np.ndarray, 'n_surfaces']¶
Flux-surface-averaged \(\langle B \cdot B \rangle\) on full-grid.
- bdotgradv: jt.Float[np.ndarray, 'n_surfaces']¶
- jpar2: jt.Float[np.ndarray, 'n_surfaces']¶
langle j_{||}^2 rangle.
- Type:
Flux-surface-averaged squared parallel current density
- Type:
math
- jperp2: jt.Float[np.ndarray, 'n_surfaces']¶
langle j_{perp}^2 rangle.
- Type:
Flux-surface-averaged squared perpendicular current density
- Type:
math
- phin: jt.Float[np.ndarray, 'n_surfaces']¶
Normalized, enclosed toroidal flux at each radial surface.
phin = toroidal_flux/toroidal_flux[-1]
- jsupu3: jt.Float[np.ndarray, 'num_full nZnT']¶
Contravariant current density component j^u on the full grid.
\(j^u = itheta/V'\)
- jsupv3: jt.Float[np.ndarray, 'num_full nZnT']¶
Contravariant current density component j^v on the full grid.
\(j^u = izeta/V'\)
- jsups3: jt.Float[np.ndarray, 'num_half nZnT']¶
Contravariant current density component \(j^s\) on the half grid.
\(j^s = \frac{\partial B_\theta}{\partial \Phi} - \frac{\partial B_\phi}{\partial \Theta}{\mu_0 V'}\)
- bsupu3: jt.Float[np.ndarray, 'num_full nZnT']¶
- bsupv3: jt.Float[np.ndarray, 'num_full nZnT']¶
- jcrossb: jt.Float[np.ndarray, 'num_full nZnT']¶
Magnitude of \(j \times B\) at each grid point.
- jxb_gradp: jt.Float[np.ndarray, 'num_full nZnT']¶
Dot product of \(j \times B\) and \(\nabla p\) at each grid point.
- jdotb_sqrtg: jt.Float[np.ndarray, 'num_full nZnT']¶
Product of \(j \cdot B\) and \(\sqrt{g}\) at each grid point.
- sqrtg3: jt.Float[np.ndarray, 'num_full nZnT']¶
Jacobian determinant \(\sqrt{g}\) at each grid point.
- bsubu3: jt.Float[np.ndarray, 'num_half nZnT']¶
- bsubv3: jt.Float[np.ndarray, 'num_half nZnT']¶
- bsubs3: jt.Float[np.ndarray, 'num_full nZnT']¶
- class vmecpp.Mercier(**data)¶
Bases:
BaseModelWithNumpy
- Parameters:
data (
Any
)
- model_config: ClassVar[ConfigDict] = {'arbitrary_types_allowed': True, 'extra': 'forbid', 'ser_json_inf_nan': 'strings'}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- s: jt.Float[np.ndarray, 'n_surfaces']¶
Normalized toroidal flux coordinate s.
- toroidal_flux: jt.Float[np.ndarray, 'n_surfaces']¶
Enclosed toroidal magnetic flux phi.
- iota: jt.Float[np.ndarray, 'n_surfaces']¶
Rotational transform iota.
- shear: jt.Float[np.ndarray, 'n_surfaces']¶
Magnetic shear profile.
- d_volume_d_s: jt.Float[np.ndarray, 'n_surfaces']¶
Radial derivative of plasma volume with respect to s.
- well: jt.Float[np.ndarray, 'n_surfaces']¶
Magnetic well profile.
- toroidal_current: jt.Float[np.ndarray, 'n_surfaces']¶
Enclosed toroidal current profile.
- d_toroidal_current_d_s: jt.Float[np.ndarray, 'n_surfaces']¶
Radial derivative of enclosed toroidal current.
- pressure: jt.Float[np.ndarray, 'n_surfaces']¶
Pressure profile p.
- d_pressure_d_s: jt.Float[np.ndarray, 'n_surfaces']¶
Radial derivative of pressure profile.
- DMerc: jt.Float[np.ndarray, 'n_surfaces']¶
Full Mercier stability criterion.
- Dshear: jt.Float[np.ndarray, 'n_surfaces']¶
Mercier criterion contribution due to magnetic shear.
- Dwell: jt.Float[np.ndarray, 'n_surfaces']¶
Mercier criterion contribution due to magnetic well.
- Dcurr: jt.Float[np.ndarray, 'n_surfaces']¶
Mercier criterion contribution due to plasma currents.
- Dgeod: jt.Float[np.ndarray, 'n_surfaces']¶
Mercier criterion contribution due to geodesic curvature.
- class vmecpp.Threed1Volumetrics(**data)¶
Bases:
BaseModelWithNumpy
- Parameters:
data (
Any
)
- model_config: ClassVar[ConfigDict] = {'arbitrary_types_allowed': True, 'extra': 'forbid', 'ser_json_inf_nan': 'strings'}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- int_p: float¶
Total plasma pressure integrated over the plasma volume.
- avg_p: float¶
Volume-averaged plasma pressure.
- int_bpol: float¶
Total poloidal magnetic field energy B_phi^2/(2 mu0) integrated over the plasma volume.
- avg_bpol: float¶
Volume-averaged poloidal magnetic field energy.
- int_btor: float¶
Total toroidal magnetic field energy integrated over the plasma volume.
- avg_btor: float¶
Volume-averaged toroidal magnetic field energy.
- int_modb: float¶
Total |B| integrated over the plasma volume.
- avg_modb: float¶
Volume-averaged |B|.
- int_ekin: float¶
Total kinetic energy integrated over the plasma volume.
- avg_ekin: float¶
Volume-averaged kinetic energy.
- class vmecpp.MakegridParameters(**data)¶
Bases:
BaseModelWithNumpy
Pydantic model mirroring the C++ makegrid::MakegridParameters struct.
Parameters defining the grid for external magnetic field calculations (mgrid file).
- Parameters:
data (
Any
)
- model_config: ClassVar[ConfigDict] = {'arbitrary_types_allowed': True, 'extra': 'forbid', 'ser_json_inf_nan': 'strings'}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- normalize_by_currents: bool¶
If true, normalize the magnetic field by coil currents and windings.
- assume_stellarator_symmetry: bool¶
If true, compute on half-period and mirror.
- number_of_field_periods: int¶
Number of toroidal field periods.
- r_grid_minimum: float¶
Radial coordinate of the first grid point.
- r_grid_maximum: float¶
Radial coordinate of the last grid point.
- number_of_r_grid_points: int¶
Number of radial grid points.
- z_grid_minimum: float¶
Vertical coordinate of the first grid point.
- z_grid_maximum: float¶
Vertical coordinate of the last grid point.
- number_of_z_grid_points: int¶
Number of vertical grid points.
- number_of_phi_grid_points: int¶
Number of toroidal grid points per field period.
- static from_file(input_file)¶
- Parameters:
- Return type:
- class vmecpp.MagneticFieldResponseTable(**data)¶
Bases:
BaseModelWithNumpy
Pydantic model mirroring the C++ makegrid::MagneticFieldResponseTable struct.
Holds the precomputed magnetic field response on a grid, separated by coil circuit. Each field component (b_r, b_p, b_z) is a list where each element corresponds to a circuit and contains the flattened 1D array of field values on the grid for that circuit carrying unit current.
- Parameters:
data (
Any
)
- model_config: ClassVar[ConfigDict] = {'arbitrary_types_allowed': True, 'extra': 'forbid', 'ser_json_inf_nan': 'strings'}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- parameters: MakegridParameters¶
- b_r: jt.Float[np.ndarray, 'num_coils num_mgrid_cells']¶
Cylindrical R components of magnetic field per circuit.
- b_p: jt.Float[np.ndarray, 'num_coils num_mgrid_cells']¶
Cylindrical Phi components of magnetic field per circuit.
- b_z: jt.Float[np.ndarray, 'num_coils num_mgrid_cells']¶
Cylindrical Z components of magnetic field per circuit.
- static from_coils_file(coils_path, makegrid_parameters)¶
- Parameters:
makegrid_parameters (
MakegridParameters
)
- Return type:
- vmecpp.populate_raw_profile(vmec_input, field, f)¶
Populate a line segment profile using callable
f
.The callable is evaluated on all unique
s
values required for the multi-grid steps (full and half grids). The resulting knots and values are stored in the auxiliary arrays for the chosen profile.
Submodules¶
- vmecpp.simsopt_compat module
Vmec
Vmec.verbose
Vmec.mpi
Vmec.n_pressure
Vmec.n_current
Vmec.n_iota
Vmec.wout
Vmec.s_full_grid
Vmec.ds
Vmec.s_half_grid
Vmec.indata
Vmec.input_file
Vmec.iter
Vmec.free_boundary
Vmec.runnable
Vmec.output_file
Vmec.need_to_run_code
Vmec.recompute_bell()
Vmec.run()
Vmec.load_wout_from_outfile()
Vmec.aspect()
Vmec.volume()
Vmec.iota_axis()
Vmec.iota_edge()
Vmec.mean_iota()
Vmec.mean_shear()
Vmec.get_dofs()
Vmec.set_dofs()
Vmec.vacuum_well()
Vmec.external_current()
Vmec.boundary
Vmec.set_indata()
Vmec.get_input()
Vmec.write_input()
Vmec.set_mpol_ntor()