simsopt.util package

Submodules

simsopt.util.constants module

simsopt.util.famus_helpers module

This file contains a class for reading in FOCUS-style data files.

This file was copied over from the MAGPIE code used for generating permanent magnet structures. All credit for this code is thanks to the PM4Stell team and Ken Hammond for his consent to use this file and work with the permanent magnet branch of SIMSOPT.

class simsopt.util.famus_helpers.FocusData(filename, downsample=1)

Bases: object

Class object for reading in FOCUS-style data files.

Parameters:

filename – a FOCUS file

adjust_rho(q_new)

Adjust rho according to the desired momentq

Parameters:

q_new – New exponent to use for FAMUS representation of dipole magnitudes.

flip_negative_magnets()

For magnets with negative rho values, adjust the orientation angles (mt and mp) so that the rho value can be negated to a positive value

float_inds = [3, 4, 5, 7, 8, 10, 11]
init_pol_vecs(n_pol)

Initializes arrays for the x, y, and z components of allowable polarization vectors, as an array giving the ID of the polarization vector actually used for the respective magnet (0 means magnet is off)

Parameters:

n_pol – int Number of allowed polarization vectors, typically 3 or 12.

perp_vector(inds)

Returns the x, y, and z components of a Cartesian unit vector perpendicular to the polarizaton direction of the nth magnet

print_to_file(filename)

Write the FocusData class object information to a FOCUS style file.

Parameters:

filename – string denoting the name of the output file.

propNames = ['type', 'symm', 'coilname', 'ox', 'oy', 'oz', 'Ic', 'M_0', 'pho', 'Lc', 'mp', 'mt', 'op']
read_from_file(filename, downsample)
repeat_hp_to_fp(nfp, magnet_sector=1)

Duplicates the magnets to the adjacent half-period, implementing the appropriate transformations to uphold stellarator symmetry.

Notes

  • At the present time, duplicate magnets will have the same pol_id as their corresponding originals

  • At the present time, duplicate magnets will have the same name their corresponding originals.

Parameters:
  • nfp (integer) – Number of field periods in the configuration

  • magnet_sector (integer (optional)) – Sector (half-period) of the torus in which the magnets are assumed to be located. Must be between 1 and 2*nfp, inclusive. Sector 1 starts at toroidal angle 0.

unit_vector(inds)

Returns the x, y, and z components of a Cartesian unit vector in the polarizaton direction of the nth magnet

class simsopt.util.famus_helpers.FocusPlasmaBnormal(fname)

Bases: object

bnormal_grid(nphi, ntheta, range)

Calculates the normal component of the magnetic field on a regularly-spaced 2D grid of points on the surface with a specified resolution.

Parameters:
  • nphi (int) – Number of grid points in the toroidal dimension per field period.

  • ntheta (int) – Number of grid points in the poloidal dimension.

  • range (str) – Extent of the grid in the toroidal dimension. Choices include half period, field period, or full torus.

Returns:

bnormal – Normal component of the magnetic field evaluated at each of the grid points.

Return type:

2D array

simsopt.util.famus_helpers.stell_point_transform(mode, phi, x_in, y_in, z_in)
Transforms a point in one of two ways, depending on the mode selected:
reflect: Reflects a point in a given poloidal plane presumed to form

the boundary between two stellarator-symmetryc half-periods. The output point will have the equivalent location associated with the adjacent half-period.

translate: Translates a point in the toroidal direction (i.e., along

a circle with a fixed radius about the z axis) by a given angle.

Parameters:
  • mode (string) – ‘translate’ or ‘reflect’ (see description above)

  • phi (double) – ‘reflect’ mode: toroidal angle (rad.) of the symmetry plane ‘translate’ mode: toroidal interval (rad.) along which to translate

  • x_in (double (possibly arrays of equal dimensions)) – x, y, z coordinates of the point(s) to be transformed

  • y_in (double (possibly arrays of equal dimensions)) – x, y, z coordinates of the point(s) to be transformed

  • z_in (double (possibly arrays of equal dimensions)) – x, y, z coordinates of the point(s) to be transformed

Returns:

x_out, y_out, z_out – x, y, z coordinates of the transformed point(s)

Return type:

double (possibly arrays)

simsopt.util.famus_helpers.stell_vector_transform(mode, phi, vx_in, vy_in, vz_in)
Transforms a vector in one of two ways, depending on the mode selected:
reflect: Reflects a vector, with a defined origin point, in a given

poloidal plane that is presumed to form the boundary between two stellarator-symmetryc half-periods. The output vector will have the equivalent origin point and direction associated with the target half-period. Vector magnitude is preserved.

translate: Translates a vector, with a defined origin point, in the

toroidal direction. The origin thus moves along a circle of fixed radius about the z axis by a given angle. The azimuthal components of the vector change in order to preserve the radial and toroidal components. The vertical component of the vector, as well as its length, are held fixed.

Parameters:
  • mode (string) – ‘translate’ or ‘reflect’ (see description above)

  • phi (double) – ‘reflect’ mode: toroidal angle (rad.) of the symmetry plane ‘translate’ mode: toroidal interval (rad.) along which to translate

  • vx_in (double (possibly arrays of equal dimensions)) – x, y, z components of the input vectors

  • vy_in (double (possibly arrays of equal dimensions)) – x, y, z components of the input vectors

  • vz_in (double (possibly arrays of equal dimensions)) – x, y, z components of the input vectors

Returns:

x, y, z components of the output transformed vector(s)

Return type:

vx_out, vy_out, vz_out

simsopt.util.fourier_interpolation module

This module contains a subroutine for spectrally accurate interpolation of data that is known on a uniform grid in a periodic domain.

simsopt.util.fourier_interpolation.fourier_interpolation(fk, x)

Interpolate data that is known on a uniform grid in [0, 2pi).

This routine is based on the matlab routine fourint.m in the DMSuite package by S.C. Reddy and J.A.C. Weideman, available at http://www.mathworks.com/matlabcentral/fileexchange/29 or here: http://dip.sun.ac.za/~weideman/research/differ.html

Parameters:
  • fk – Vector of y-coordinates of data, at equidistant points x(k) = (k-1)*2*pi/N, k = 1…N

  • x – Vector of x-values where interpolant is to be evaluated.

Returns:

Array of length len(x) with the interpolated values.

simsopt.util.logger module

simsopt.util.logger.initialize_logging(filename: str | None = None, level: str | None = None, mpi: bool = False) None

Initializes logging in a simple way for both serial and MPI jobs. The MPI logging uses MPILogger package.

Parameters:
  • filename – Name of file to store the logging info

  • level – Logging level. Could be ‘INFO’, ‘DEBUG’, ‘WARNING’, etc.

  • mpi – If True MPI logging is used provided mpi4py is installed.

simsopt.util.mpi module

This module contains the MpiPartition class.

This module should be completely self-contained, depending only on mpi4py and numpy, not on any other simsopt components.

simsopt.util.mpi.log(level: int = 20)

Turn on logging. If MPI is available, the processor number will be added to all logging entries.

Parameters:

level – Typically logging.INFO for regular output, or logging.DEBUG for more extensive output.

simsopt.util.mpi.proc0_print(*args, **kwargs)

simsopt.util.mpi_logger module

simsopt.util.permanent_magnet_helper_functions module

This module contains the a number of useful functions for using the permanent magnets functionality in the SIMSOPT code.

simsopt.util.permanent_magnet_helper_functions.calculate_on_axis_B(bs, s)

Check the average, approximate, on-axis magnetic field strength to make sure the configuration is scaled correctly.

Parameters:
  • bs – MagneticField or BiotSavart class object.

  • s – plasma boundary surface.

Returns:

Average magnetic field strength along

the major radius of the device.

Return type:

B0avg

simsopt.util.permanent_magnet_helper_functions.coil_optimization(s, bs, base_curves, curves, out_dir='')

Optimize the coils for the QA, QH, or other configurations.

Parameters:
  • s – plasma boundary.

  • bs – Biot Savart class object, presumably representing the magnetic fields generated by the coils.

  • base_curves – List of CurveXYZFourier class objects.

  • curves – List of Curve class objects.

  • out_dir – Path or string for the output directory for saved files.

Returns:

Biot Savart class object, presumably representing the

OPTIMIZED magnetic fields generated by the coils.

Return type:

bs

simsopt.util.permanent_magnet_helper_functions.initialize_coils(config_flag, TEST_DIR, s, out_dir='')

Initializes coils for each of the target configurations that are used for permanent magnet optimization.

Parameters:
  • config_flag – String denoting the stellarator configuration being initialized.

  • TEST_DIR – String denoting where to find the input files.

  • out_dir – Path or string for the output directory for saved files.

  • s – plasma boundary surface.

Returns:

List of CurveXYZ class objects. curves: List of Curve class objects. coils: List of Coil class objects.

Return type:

base_curves

simsopt.util.permanent_magnet_helper_functions.initialize_default_kwargs(algorithm='RS')

Keywords to the permanent magnet optimizers are now passed by the kwargs dictionary. Default dictionaries are initialized here.

Parameters:

algorithm – String denoting which algorithm is being used for permanent magnet optimization. Options are ‘RS’ (relax and split), ‘GPMO’ (greedy placement), ‘backtracking’ (GPMO with backtracking), ‘multi’ (GPMO, placing multiple magnets each iteration), ‘ArbVec’ (GPMO with arbitrary dipole orientiations), ‘ArbVec_backtracking’ (GPMO with backtracking and arbitrary dipole orientations).

Returns:

Dictionary of keywords to pass to a permanent magnet

optimizer.

Return type:

kwargs

simsopt.util.permanent_magnet_helper_functions.make_Bnormal_plots(bs, s_plot, out_dir='', bs_filename='Bnormal')

Plot Bnormal on plasma surface from a MagneticField object. Do this quite a bit in the permanent magnet optimization and initialization so this is a wrapper function to reduce the amount of code.

Parameters:
  • bs – MagneticField class object.

  • s_plot – plasma boundary surface with range = ‘full torus’.

  • out_dir – Path or string for the output directory for saved files.

  • bs_filename – String denoting the name of the output file.

simsopt.util.permanent_magnet_helper_functions.make_optimization_plots(RS_history, m_history, m_proxy_history, pm_opt, out_dir='')

Make line plots of the algorithm convergence and make histograms of m, m_proxy, and if available, the FAMUS solution.

Parameters:
  • RS_history – List of the relax-and-split optimization progress.

  • m_history – List of the permanent magnet solutions generated as the relax-and-split optimization progresses.

  • m_proxy_history – Same but for the ‘proxy’ variable in the relax-and-split optimization.

  • pm_opt – PermanentMagnetGrid class object that was optimized.

  • out_dir – Path or string for the output directory for saved files.

simsopt.util.permanent_magnet_helper_functions.make_qfm(s, Bfield)

Given some Bfield generated by dipoles AND a set of TF coils, compute a quadratic flux-minimizing surface (QFMS) for the total field configuration on the surface s.

Parameters:
  • s – plasma boundary surface.

  • Bfield – MagneticField or inherited MagneticField class object.

Returns:

The identified QfmSurface class object.

Return type:

qfm_surface

simsopt.util.permanent_magnet_helper_functions.read_focus_coils(filename)

Reads in the coils from a FOCUS file. For instance, this is used for loading in the MUSE phased TF coils.

Parameters:

filename – String denoting the name of the coils file.

Returns:

List of CurveXYZFourier class objects. base_currents: List of Current class objects. ncoils: Integer representing the number of coils.

Return type:

coils

simsopt.util.permanent_magnet_helper_functions.run_Poincare_plots(s_plot, bs, b_dipole, comm, filename_poincare, out_dir='')

Wrapper function for making Poincare plots.

Parameters:
  • s_plot – plasma boundary surface with range = ‘full torus’.

  • bs – MagneticField class object.

  • b_dipole – DipoleField class object.

  • comm – MPI COMM_WORLD object for using MPI for tracing.

  • filename_poincare – Filename for the output poincare picture.

  • out_dir – Path or string for the output directory for saved files.

simsopt.util.permanent_magnet_helper_functions.trace_fieldlines(bfield, label, s, comm, out_dir='')

Make Poincare plots on a surface as in the trace_fieldlines example in the examples/1_Simple/ directory.

Parameters:
  • bfield – MagneticField or InterpolatedField class object.

  • label – Name of the file to write to.

  • s – plasma boundary surface.

  • comm – MPI COMM_WORLD object for using MPI for tracing.

  • out_dir – Path or string for the output directory for saved files.

simsopt.util.polarization_project module

This file contains definitions and functions for generating the PM4Stell magnet orientations and other magnet functionality.

This file was copied over from the MAGPIE code used for generating permanent magnet structures. All credit for this code is thanks to the PM4Stell team and Ken Hammond for his consent to use this file and work with the permanent magnet branch of SIMSOPT.

simsopt.util.polarization_project.discretize_polarizations(mag_data, orientation_phi, pol_axes, pol_type)

Adjusts the polarization axes of a set of dipole moments to the closest from a set of allowable axes.

Also populates the pol_x, pol_y, pol_z, and pol_id arrays within the input instance of the mag_data class, providing all of the allowable polarization vectors for each magnet in the lab frame.

Parameters:
  • mag_data – mag_data object Instance of the mag_data class characterizing the dipole moments of a set of magnets, whose orientations are to be updated.

  • orientation_phi – double array Orientation angle defining a reference axes for the allowable axes for each respective magnet (as described in the documentation for orientation_phi().) Number of elements must agree with the number of magnets in mag_data.

  • pol_axes – double array with 3 columns Set of Cartesian unit vectors defining the allowable polarization vectors for each magnet. The first, second, and third columns contain the x, y, and z components, respectively.

  • pol_type – integer array Type IDs corresponding to each row in pol_axes.

simsopt.util.polarization_project.edge_triplet(theta_fe, theta_fc)

Generates a set of vectors corresponding to the “edge-triplet” set of polarization types: edge, face-edge, and face-corner.

Parameters:
  • theta_fe – double Angle with respect to the nearest face-normal for the face-edge type

  • theta_fc – double Angle with respect to the nearest face-normal for face-corner type

Returns:

double array

3-column array in which each row is one unit vector in the full set

Return type:

vectors

simsopt.util.polarization_project.face_triplet(theta_fe, theta_fc)

Generates a set of vectors corresponding to the “face-triplet” set of polarization types: face, face-edge, and face-corner.

Parameters:
  • theta_fe – double Angle with respect to the nearest face-normal for the face-edge type

  • theta_fc – double Angle with respect to the nearest face-normal for face-corner type

Returns:

double array

3-column array in which each row is one unit vector in the full set

Return type:

vectors

simsopt.util.polarization_project.orientation_phi(corners_fname)

Determines the azimuthal (phi) angle that sets the orientations of rectangular prisms from an arrangement based on data from the corners file

Phi is the azimuthal angle of the normal vector of the face whose normal vector is in the “radial-like” direction with a polar angle of pi/2 radians (i.e., parallel to the x-y plane). The remaining faces will be assumed to have normal vectors in either:

  1. the z-direction, or

  2. parallel to the x-y plane with an azimuthal angle equal to the orientation phi plus pi/2 radians (the “phi-like” direction).

Parameters:

corners_fname (string) – Corners file

Returns:

Orientation phi angle (radians) of each prism in the arrangement

Return type:

orientation_phi (double array)

simsopt.util.polarization_project.polarization_axes(polarizations)

Populates an array of allowable polarization vectors based on an input list of names. Also outputs a 1D array with an integer corresponding to the polarization type in the order it appears in the list of names.

Parameters:

polarizations – list of strings List of names of polarization types. Each type should not be listed more than once

Returns:

2d double array

Allowable axes for the magnets expressed as unit vectors in local cylindrical coordinates. The first, second, and third columns contain the r, phi, and z components, respectively.

pol_type: 1d integer array

ID number for the polarization type associated with each row in pol_axes. The ID corresponds to the order in which the type appears in the input list (polarizations). The first type is given an ID of 1.

Return type:

pol_axes