vmecpp package¶
- class vmecpp.JxBOut(**data)¶
Bases:
BaseModel
- Parameters:
data (
Any
)
- amaxfor: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- aminfor: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- avforce: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- bdotb: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- bdotgradv: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- bdotk: npyd.NDArray[npyd.Shape['* num_full, * nZnT'], float]¶
- bsubs3: npyd.NDArray[npyd.Shape['* num_full, * nZnT'], float]¶
- bsubu3: npyd.NDArray[npyd.Shape['* num_half, * nZnT'], float]¶
- bsubv3: npyd.NDArray[npyd.Shape['* num_half, * nZnT'], float]¶
- bsupu3: npyd.NDArray[npyd.Shape['* num_full, * nZnT'], float]¶
- bsupv3: npyd.NDArray[npyd.Shape['* num_full, * nZnT'], float]¶
- itheta: npyd.NDArray[npyd.Shape['* num_full, * nZnT'], float]¶
- izeta: npyd.NDArray[npyd.Shape['* num_full, * nZnT'], float]¶
- jcrossb: npyd.NDArray[npyd.Shape['* num_full, * nZnT'], float]¶
- jdotb: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- jdotb_sqrtg: npyd.NDArray[npyd.Shape['* num_full, * nZnT'], float]¶
- jpar2: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- jperp2: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- jsups3: npyd.NDArray[npyd.Shape['* num_half, * nZnT'], float]¶
- jsupu3: npyd.NDArray[npyd.Shape['* num_full, * nZnT'], float]¶
- jsupv3: npyd.NDArray[npyd.Shape['* num_full, * nZnT'], float]¶
- jxb_gradp: npyd.NDArray[npyd.Shape['* num_full, * nZnT'], float]¶
- model_config: ClassVar[ConfigDict] = {'extra': 'forbid'}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- phin: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- pprim: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- sqrtg3: npyd.NDArray[npyd.Shape['* num_full, * nZnT'], float]¶
- class vmecpp.Mercier(**data)¶
Bases:
BaseModel
- Parameters:
data (
Any
)
- DMerc: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- Dcurr: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- Dgeod: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- Dshear: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- Dwell: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- d_pressure_d_s: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- d_toroidal_current_d_s: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- d_volume_d_s: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- iota: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- model_config: ClassVar[ConfigDict] = {'extra': 'forbid'}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- pressure: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- s: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- shear: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- toroidal_current: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- toroidal_flux: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- well: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- class vmecpp.Threed1Volumetrics(**data)¶
Bases:
BaseModel
- Parameters:
data (
Any
)
- avg_bpol: float¶
- avg_btor: float¶
- avg_ekin: float¶
- avg_modb: float¶
- avg_p: float¶
- int_bpol: float¶
- int_btor: float¶
- int_ekin: float¶
- int_modb: float¶
- int_p: float¶
- model_config: ClassVar[ConfigDict] = {'extra': 'forbid'}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- class vmecpp.VmecInput(**data)¶
Bases:
BaseModel
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 and model_dump_json.
- Parameters:
data (
Any
)
- ac: npyd.NDArray[npyd.Shape['* ac_len'], float]¶
Enclosed toroidal current profile coefficients.
- ac_aux_f: npyd.NDArray[npyd.Shape['* ac_aux_len'], float]¶
values at knots
- Type:
Spline toroidal current profile
- ac_aux_s: npyd.NDArray[npyd.Shape['* ac_aux_len'], float]¶
knot locations in s
- Type:
Spline toroidal current profile
- ai: npyd.NDArray[npyd.Shape['* ai_len'], float]¶
Iota profile coefficients.
- ai_aux_f: npyd.NDArray[npyd.Shape['* ai_aux_len'], float]¶
values at knots
- Type:
Spline iota profile
- ai_aux_s: npyd.NDArray[npyd.Shape['* ai_aux_len'], float]¶
knot locations in s
- Type:
Spline iota profile
- am: npyd.NDArray[npyd.Shape['* am_len'], float]¶
Mass/pressure profile coefficients.
- am_aux_f: npyd.NDArray[npyd.Shape['* am_aux_len'], float]¶
values at knots
- Type:
Spline mass/pressure profile
- am_aux_s: npyd.NDArray[npyd.Shape['* am_aux_len'], float]¶
knot locations in s
- Type:
Spline mass/pressure profile
- aphi: npyd.NDArray[npyd.Shape['* aphi_len'], float]¶
Radial flux zoning profile coefficients.
- bloat: float¶
Bloating factor (for constrained toroidal current)
- curtor: float¶
Toroidal current in A.
- delt: float¶
Initial value for artificial time step in iterative solver.
- extcur: npyd.NDArray[npyd.Shape['* extcur_len'], float]¶
Coil currents in A.
- static from_file(input_file)¶
Build a VmecInput from either a VMEC++ JSON input file or a classic INDATA file.
- ftol_array: npyd.NDArray[npyd.Shape['* num_grids'], float]¶
Requested force tolerance for convergence per multigrid step.
- gamma: float¶
Adiabatic index.
- lasym: bool¶
Flag to indicate non-stellarator-symmetry.
Note: this flag is False if stellarator symmetry is present, True if not.
- lforbal: bool¶
directly compute innermost flux surface geometry from radial force balance
- Type:
Hack
- lfreeb: bool¶
Flag to indicate free-boundary.
- mgrid_file: str¶
Full path for vacuum Green’s function data.
- model_config: ClassVar[ConfigDict] = {'ser_json_inf_nan': 'strings'}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- mpol: int¶
Number of poloidal Fourier harmonics; m = 0, 1, …, (mpol-1)
- ncurr: int¶
constrained-current
- Type:
Select constraint on iota or enclosed toroidal current profiles 0
- Type:
constrained-iota; 1
- nfp: int¶
Number of toroidal field periods (=1 for Tokamak)
- niter_array: npyd.NDArray[npyd.Shape['* num_grids'], int]¶
Maximum number of iterations per multigrid step.
- ns_array: npyd.NDArray[npyd.Shape['* num_grids'], int]¶
Number of flux surfaces per multigrid step.
- nstep: int¶
Printout interval.
- ntheta: int¶
is rounded to next smaller even number.
- Type:
Number of poloidal grid points; if odd
- ntor: int¶
Number of toroidal Fourier harmonics; n = -ntor, -ntor+1, …, -1, 0, 1, …, ntor-1, ntor.
- nvacskip: int¶
Number of iterations between full vacuum calculations.
- nzeta: int¶
Number of toroidal grid points; must match nzeta of mgrid file if using free- boundary.
- pcurr_type: str¶
Parametrization of toroidal current profile.
- phiedge: float¶
Total enclosed toroidal magnetic flux in Vs == Wb.
- piota_type: str¶
Parametrization of iota profile.
- pmass_type: str¶
Parametrization of mass/pressure profile.
- pres_scale: float¶
Global scaling factor for mass/pressure profile.
- raxis_c: npyd.NDArray[npyd.Shape['* ntor_plus_1'], float]¶
Magnetic axis coefficients for R ~ cos(n*v); stellarator-symmetric.
- raxis_s: npyd.NDArray[npyd.Shape['* ntor_plus_1'], float] | None¶
Magnetic axis coefficients for R ~ sin(n*v); non-stellarator-symmetric.
- rbc: SerializableSparseCoefficientArray¶
Boundary coefficients for R ~ cos(m*u - n*v); stellarator-symmetric
- rbs: SerializableSparseCoefficientArray | None¶
Boundary coefficients for R ~ sin(m*u - n*v); non-stellarator-symmetric
- 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:
NDArray
[Shape
[* mpol_new, * two_ntor_new_plus_one], (float16
,float32
,float64
,float32
,float64
)]
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.]])
- return_outputs_even_if_not_converged: bool¶
If true, return the outputs even if VMEC++ did not converge.
Otherwise a RuntimeError will be raised.
- spres_ped: float¶
Location of pressure pedestal in s.
- tcon0: float¶
Constraint force scaling factor for ns –> 0.
- zaxis_c: npyd.NDArray[npyd.Shape['* ntor_plus_1'], float] | None¶
Magnetic axis coefficients for Z ~ cos(n*v); non-stellarator-symmetric.
- zaxis_s: npyd.NDArray[npyd.Shape['* ntor_plus_1'], float]¶
Magnetic axis coefficients for Z ~ sin(n*v); stellarator-symmetric.
- zbc: SerializableSparseCoefficientArray | None¶
Boundary coefficients for Z ~ cos(m*u - n*v); non-stellarator-symmetric
- zbs: SerializableSparseCoefficientArray¶
Boundary coefficients for Z ~ sin(m*u - n*v); stellarator-symmetric
- class vmecpp.VmecOutput(**data)¶
Bases:
BaseModel
Container for the full output of a VMEC run.
- Parameters:
data (
Any
)
- 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.
- model_config: ClassVar[ConfigDict] = {}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- threed1_volumetrics: Threed1Volumetrics¶
Python equivalent of VMEC’s volumetrics section in the “threed1” file.
- wout: VmecWOut¶
Python equivalent of VMEC’s “wout” file.
- class vmecpp.VmecWOut(**data)¶
Bases:
BaseModel
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 VMEC.
- Parameters:
data (
Any
)
- Aminor_p: float¶
- DCurr: npyd.NDArray[npyd.Shape['*'], float]¶
- DGeod: npyd.NDArray[npyd.Shape['*'], float]¶
- DMerc: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- DShear: npyd.NDArray[npyd.Shape['*'], float]¶
- DWell: npyd.NDArray[npyd.Shape['*'], float]¶
- IonLarmor: float¶
- Rmajor_p: float¶
- ac: npyd.NDArray[npyd.Shape[str(_PRESET_DIM)], float]¶
- ac_aux_f: npyd.NDArray[npyd.Shape[str(_NDF_MAX_DIM)], float]¶
- ac_aux_s: npyd.NDArray[npyd.Shape[str(_NDF_MAX_DIM)], float]¶
- ai: npyd.NDArray[npyd.Shape[str(_PRESET_DIM)], float]¶
- ai_aux_f: npyd.NDArray[npyd.Shape[str(_NDF_MAX_DIM)], float]¶
- ai_aux_s: npyd.NDArray[npyd.Shape[str(_NDF_MAX_DIM)], float]¶
- am: npyd.NDArray[npyd.Shape[str(_PRESET_DIM)], float]¶
- am_aux_f: npyd.NDArray[npyd.Shape[str(_NDF_MAX_DIM)], float]¶
- am_aux_s: npyd.NDArray[npyd.Shape[str(_NDF_MAX_DIM)], float]¶
- aspect: float¶
- b0: float¶
- bdotb: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- bdotgradv: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- beta_vol: npyd.NDArray[npyd.Shape['*'], float]¶
- betapol: float¶
- betator: float¶
- betatotal: float¶
- betaxis: float¶
- bmnc: npyd.NDArray[npyd.Shape['* dim1, * dim2'], float]¶
- bsubsmns: npyd.NDArray[npyd.Shape['* dim1, * dim2'], float]¶
- bsubumnc: npyd.NDArray[npyd.Shape['* dim1, * dim2'], float]¶
- bsubvmnc: npyd.NDArray[npyd.Shape['* dim1, * dim2'], float]¶
- bsupumnc: npyd.NDArray[npyd.Shape['* dim1, * dim2'], float]¶
- bsupvmnc: npyd.NDArray[npyd.Shape['* dim1, * dim2'], float]¶
- buco: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- bvco: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- chi: npyd.NDArray[npyd.Shape['*'], float]¶
- chipf: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- ctor: float¶
- equif: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- extcur: npyd.NDArray[npyd.Shape['* extcur_len'], float] | float¶
Coil currents in A.
for free-boundary runs, extcur has shape (nextcur,) for fixed-boundary it is a scalar float extcur=nan
- 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.
- fsql: float¶
- fsqr: float¶
- fsqz: float¶
- ftolv: float¶
- gamma: float¶
- gmnc: npyd.NDArray[npyd.Shape['* dim1, * dim2'], float]¶
- ier_flag: int¶
- input_extension: str¶
- iotaf: npyd.NDArray[npyd.Shape['*'], float]¶
- iotas: npyd.NDArray[npyd.Shape['*'], float]¶
- jcuru: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- jcurv: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- jdotb: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- lasym: bool¶
- property lasym__logical__¶
This is how the attribute is called in the Fortran wout file.
- lfreeb: bool¶
- property lfreeb__logical__¶
This is how the attribute is called in the Fortran wout file.
- lmns: npyd.NDArray[npyd.Shape['* mnmax, * n_surfaces'], float]¶
- lmns_full: npyd.NDArray[npyd.Shape['* mnmax, * n_surfaces'], float]¶
- mass: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- mgrid_file: str¶
- mgrid_mode: MgridModeType¶
- mnmax: int¶
- mnmax_nyq: int¶
- model_config: ClassVar[ConfigDict] = {'extra': 'forbid'}¶
Configuration for the model, should be a dictionary conforming to [ConfigDict][pydantic.config.ConfigDict].
- mpol: int¶
- nextcur: int¶
- nfp: int¶
- niter: int¶
- ns: int¶
- ntor: int¶
- over_r: npyd.NDArray[npyd.Shape['*'], float]¶
- pcurr_type: str¶
- phi: npyd.NDArray[npyd.Shape['*'], float]¶
- phipf: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- phips: npyd.NDArray[npyd.Shape['* dim1'], float]¶
- piota_type: str¶
- pmass_type: str¶
- pres: npyd.NDArray[npyd.Shape['*'], float]¶
- presf: npyd.NDArray[npyd.Shape['*'], float]¶
- q_factor: npyd.NDArray[npyd.Shape['*'], float]¶
- raxis_cc: npyd.NDArray[npyd.Shape['*'], float]¶
- rbtor: float¶
- rbtor0: float¶
- rmax_surf: float¶
- rmin_surf: float¶
- rmnc: npyd.NDArray[npyd.Shape['* dim1, * dim2'], float]¶
- save(out_path)¶
Save contents in NetCDF3 format.
This is the format used by Fortran VMEC implementations and the one expected by SIMSOPT.
- signgs: int¶
- specw: npyd.NDArray[npyd.Shape['*'], float]¶
- version_: float¶
- volavgB: float¶
- volume: float¶
- property volume_p¶
The attribute is called volume_p in the Fortran wout file, while simsopt.mhd.Vmec.wout uses volume.
We expose both.
- vp: npyd.NDArray[npyd.Shape['*'], float]¶
- wb: float¶
- wp: float¶
- xm: npyd.NDArray[npyd.Shape['* dim1'], int]¶
- xm_nyq: npyd.NDArray[npyd.Shape['* dim1'], int]¶
- xn: npyd.NDArray[npyd.Shape['* dim1'], int]¶
- xn_nyq: npyd.NDArray[npyd.Shape['* dim1'], int]¶
- zaxis_cs: npyd.NDArray[npyd.Shape['*'], float]¶
- zmax_surf: float¶
- zmns: npyd.NDArray[npyd.Shape['* dim1, * dim2'], float]¶
- 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).
- Keyword Arguments:
max_threads – 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 – if True, VMEC++ logs its progress to standard output.
restart_from – 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
- Parameters:
verbose (
bool
)restart_from (
Optional
[VmecOutput
])
Submodules¶
- vmecpp.simsopt_compat module
Vmec
Vmec.aspect()
Vmec.boundary
Vmec.ds
Vmec.external_current()
Vmec.free_boundary
Vmec.get_dofs()
Vmec.get_input()
Vmec.indata
Vmec.input_file
Vmec.iota_axis()
Vmec.iota_edge()
Vmec.iter
Vmec.load_wout_from_outfile()
Vmec.mean_iota()
Vmec.mean_shear()
Vmec.mpi
Vmec.n_current
Vmec.n_iota
Vmec.n_pressure
Vmec.need_to_run_code
Vmec.output_file
Vmec.recompute_bell()
Vmec.run()
Vmec.runnable
Vmec.s_full_grid
Vmec.s_half_grid
Vmec.set_dofs()
Vmec.set_indata()
Vmec.set_mpol_ntor()
Vmec.vacuum_well()
Vmec.verbose
Vmec.volume()
Vmec.wout
Vmec.write_input()