API guide#

This guide takes you through the classes used in Struphy launch files generated by:

struphy params MODEL

for a valid MODEL from the list of Models. All information for running simulations are contained and can be altered in such launch (parameter) files. It is recommended to generate such a file and inspect it with an editor like VScode (using Pylance).

Simulation class#

class struphy.Simulation(model: StruphyModel, name: str = '', description: str = '', params_path: str | None = None, env: EnvironmentOptions = EnvironmentOptions(out_folders='/home/runner/work/struphy/struphy/doc', sim_folder='sim_1', restart=False, max_runtime=300, save_step=1, sort_step=0, num_clones=1, profiling_activated=False, profiling_trace=False), base_units: BaseUnits = BaseUnits(x=1.0, B=1.0, n=1.0, kBT=None), time_opts: Time = Time(dt=0.01, Tend=0.03, split_algo='LieTrotter'), domain: Domain = Cuboid(l1=0.0, r1=1.0, l2=0.0, r2=1.0, l3=0.0, r3=1.0), equil: FluidEquilibrium | None = None, grid: TensorProductGrid = TensorProductGrid(Nel=(24, 10, 1), mpi_dims_mask=(True, True, True)), derham_opts: DerhamOptions = DerhamOptions(p=(1, 1, 1), spl_kind=(True, True, True), dirichlet_bc=((False, False), (False, False), (False, False)), nquads=None, nq_pr=None, polar_ck=-1, local_projectors=False), verbose: bool = False)[source]#

Bases: SimulationBase

Top-level class to configure and run a Struphy simulation.

The Simulation class wraps model setup, MPI configuration, output management, normalization (units), memory allocation and time stepping. It initializes the model’s variables and propagators, prepares runtime metadata and output folders, and provides the main run() entry point to execute the simulation.

Parameters:
  • model (StruphyModel) – Physics model that provides species, propagators and variables.

  • name (str, optional) – Name of the simulation.

  • description (str, optional) – Description of the simulation.

  • params_path (str, optional) – Path to a Python parameter file to save alongside outputs.

  • env (EnvironmentOptions) – Runtime and output environment options.

  • base_units (BaseUnits) – Units used for normalization.

  • time_opts (Time) – Time-stepping options (dt, Tend, split algorithm, …).

  • domain (Domain) – Computational domain description.

  • equil (FluidEquilibrium, optional) – Initial fluid equilibrium (may be None).

  • grid (TensorProductGrid) – Spatial grid used for FEEC variables.

  • derham_opts (DerhamOptions) – Options for discrete differential operators.

  • verbose (bool, optional) – If True, print additional setup information.

Variables:
  • meta (dict) – Metadata about the run (platform, python version, model name, etc.).

  • units (Units) – Unit/normalization helper created from base_units.

  • data (DataContainer) – Output container used to store simulation data.

  • start_time (float) – Wall-clock time when the simulation object was created.

Simulation parameters#

For available models see Models.

For domains see Geometry.

For fluid backgrounds see Available fluid equilibria.

class struphy.EnvironmentOptions(out_folders: str = '/home/runner/work/struphy/struphy/doc', sim_folder: str = 'sim_1', restart: bool = False, max_runtime: int = 300, save_step: int = 1, sort_step: int = 0, num_clones: int = 1, profiling_activated: bool = False, profiling_trace: bool = False)[source]#

Set environment options for launching run on current architecture (these options do not influence the simulation result).

Parameters:
  • out_folders (str) – Absolute path to directory for sim_folder.

  • sim_folder (str) – Folder in out_folders/ for the current simulation (default= sim_1/ ). Will create the folder if it does not exist OR cleans the folder for new runs.

  • restart (bool) – Whether to restart a run (default=False).

  • max_runtime (int) – Maximum run time of simulation in minutes. Will finish the time integration once this limit is reached (default=300).

  • save_step (int) – When to save data output: every time step (save_step=1), every second time step (save_step=2), etc (default=1).

  • sort_step (int, optional) – Sort markers in memory every N time steps (default=0, which means markers are sorted only at the start of simulation)

  • num_clones (int, optional) – Number of domain clones (default=1)

  • profiling_activated (bool, optional) – Activate profiling with scope-profiler (default=False)

  • profiling_trace (bool, optional) – Save time-trace of each profiling region (default=False)

class struphy.BaseUnits(x: float = 1.0, B: float = 1.0, n: float = 1.0, kBT: float | None = None)[source]#

Set base units in parameter/launch files from which other units are derived. See Normalization.

Parameters:
  • x (float) – Unit of length in meters.

  • B (float) – Unit of magnetic field in Tesla.

  • n (float) – Unit of particle number density in 1e20/m^3.

  • kBT (float, optional) – Unit of internal energy in keV. Only in effect if the velocity scale is set to ‘thermal’.

class struphy.Time(dt: float = 0.01, Tend: float = 0.03, split_algo: Literal['LieTrotter', 'Strang'] = 'LieTrotter')[source]#

Set options for time stepping in parameter/launch files.

Parameters:
  • dt (float) – Time step.

  • Tend (float) – End time.

  • split_algo (LiteralOptions.SplitAlgos) – Splitting algorithm (the order of the propagators is defined in the model).

class struphy.grids.TensorProductGrid(Nel: tuple = (24, 10, 1), mpi_dims_mask: tuple = (True, True, True))#

Bases: object

Grid as a tensor product of 1d grids.

Parameters:
  • Nel (tuple[int]) – Number of elements in each direction.

  • mpi_dims_mask (Tuple of bool) – True if the dimension is to be used in the domain decomposition (=default for each dimension). If mpi_dims_mask[i]=False, the i-th dimension will not be decomposed.

class struphy.DerhamOptions(p: tuple = (1, 1, 1), spl_kind: tuple = (True, True, True), dirichlet_bc: tuple = ((False, False), (False, False), (False, False)), nquads: tuple | None = None, nq_pr: tuple | None = None, polar_ck: Literal[-1, 1] = -1, local_projectors: bool = False)[source]#

Set options for the Derham spaces in parameter/launch files. See Geometric finite elements (FEEC).

Parameters:
  • p (tuple[int]) – Spline degree in each direction.

  • spl_kind (tuple[bool]) – Kind of spline in each direction (True=periodic, False=clamped).

  • dirichlet_bc (tuple[tuple[bool]]) – Whether to apply homogeneous Dirichlet boundary conditions (at left or right boundary in each direction).

  • nquads (tuple[int]) – Number of Gauss-Legendre quadrature points in each direction (default = p, leads to exact integration of degree 2p-1 polynomials).

  • nq_pr (tuple[int]) – Number of Gauss-Legendre quadrature points in each direction for geometric projectors (default = p+1, leads to exact integration of degree 2p+1 polynomials).

  • polar_ck (PolarRegularity) – Smoothness at a polar singularity at eta_1=0 (default -1 : standard tensor product splines, OR 1 : C1 polar splines)

  • local_projectors (bool) – Whether to build the local commuting projectors based on quasi-inter-/histopolation.

Species types#

Each Struphy model is a collection of species of one of the following types:

class struphy.models.species.FieldSpecies[source]#

Represents a field species with zero mass and charge.

Field species are used to represent electromagnetic or other non-particle fields in a plasma model. They have no direct physical mass or charge properties (charge_number = 0, mass_number = 0), but may have associated equation parameters for scaled formulations.

Examples

>>> E_field = FieldSpecies()
>>> E_field.set_species_properties(alpha=0.5, epsilon=0.1, kappa=0.2)
class struphy.models.species.FluidSpecies[source]#

Represents a single fluid species evolving in 3D configuration space.

FluidSpecies describes macroscopic plasma dynamics using fluid or moment-based equations (e.g., Euler equations, MHD equations). Each fluid species has a specific charge and mass number and evolves according to fluid propagators.

Examples

>>> ions = FluidSpecies()
>>> ions.set_species_properties(charge_number=-1, mass_number=1/1836)  # electrons
class struphy.models.species.ParticleSpecies[source]#

Represents a single kinetic species in 3D configuration space plus velocity space.

ParticleSpecies describes plasma dynamics using kinetic theory via particles or markers in 3D + vdim (where vdim is 2 or 3) phase space. This class manages particle initialization, sorting, and diagnostic output.

set_markers(loading_params, weights_params, boundary_params, bufsize)[source]#

Configure marker initialization and weight parameters.

set_sorting_boxes(do_sort, sorting_frequency, boxes_per_dim, box_bufsize, dims_mask)[source]#

Configure spatial sorting for memory efficiency and kernel evaluation.

set_save_data(n_markers, binning_plots, kernel_density_plots)[source]#

Configure diagnostic output for particles and distribution functions.

Examples

>>> electrons = ParticleSpecies()
>>> electrons.set_species_properties(charge_number=-1, mass_number=1/1836)
>>> load_params = LoadingParameters(Np=100000)
>>> electrons.set_markers(loading_params=load_params)
>>> electrons.set_sorting_boxes(do_sort=True, sorting_frequency=10)
Species.set_species_properties(charge_number: int = 1, mass_number: int = 1, alpha: float | None = None, epsilon: float | None = None, kappa: float | None = None)[source]#

Set physical and equation parameters for a plasma species.

Sets the charge and mass numbers, and optionally overrides normalized equation parameters (alpha, epsilon, kappa) that would otherwise be computed from physical units.

Parameters:
  • charge_number (int, optional) – Charge number in units of elementary charge (default = 1).

  • mass_number (int, optional) – Mass number in units of proton mass (default = 1).

  • alpha (float, optional) – Dimensionless parameter: plasma frequency / cyclotron frequency. If None, computed from units and charge/mass numbers (default = None).

  • epsilon (float, optional) – Normalized cyclotron period: 1 / (cyclotron frequency × time unit). If None, computed from units and charge/mass numbers (default = None).

  • kappa (float, optional) – Normalized plasma frequency: plasma frequency × time unit. If None, computed from units and charge/mass numbers (default = None).

Notes

This method should be called BEFORE instantiating a Simulation object. For existing simulation objects, call Simulation.normalize_model() to apply changes. A warning will be issued if this requirement is not followed.

Variable types#

Each species can contain multiple variables.

class struphy.models.variables.PICVariable(space: Literal['Particles6D', 'DeltaFParticles6D', 'Particles5D', 'Particles3D'] = 'Particles6D')[source]#

Kinetic variable using Particle-In-Cell (PIC) discretization in phase space.

PICVariable represents kinetic distribution functions discretized via marker/particle methods in 3D configuration space plus velocity space (3D+2D or 3D+3D). Supports various phase-space decompositions for efficient kinetic simulations.

Variables:
  • space (str) – Phase space decomposition. Options include, among others: - ‘Particles6D’: Full 3D config + 3D velocity space - ‘Particles5D’: 3D config + 2D gyro-velocity space (for strong magnetic fields)

  • particles (Particles) – The marker/particle object managing kinetic data. Only available after allocate().

  • particles_class (type) – Reference to the Particles class corresponding to the space.

  • species (ParticleSpecies) – Parent kinetic species containing this variable.

  • backgrounds (KineticBackground) – Mandatory equilibrium distribution function (e.g., Maxwellian).

  • initial_condition (KineticBackground) – The initial kinetic distribution. If not explicitly set, uses the background.

  • perturbations (Perturbation | list[Perturbation] | None) – Perturbations to moments of the distribution function.

  • n_as_volume_form (bool) – Whether number density is represented as a differential form (volume-weighted) or scalar.

  • n_to_save (int) – Number of markers for which trajectories are saved.

  • saved_markers (ndarray) – Array storing saved marker data.

add_background(background, n_as_volume_form)[source]#

Set mandatory equilibrium distribution (e.g., Maxwellian).

add_initial_condition(init)[source]#

Set initial kinetic distribution (must be consistent with background).

show_initial_condition()[source]#

Display current initial condition information.

allocate(clone_config, derham, domain, equil, projected_equil, verbose)[source]#

Initialize particles and allocate marker arrays.

Notes

The background is mandatory and often used for noise reduction. The initial condition should be consistent with the background to avoid numerical artifacts.

If no initial condition is specified, the background is used as the initial condition. If both are specified, they should be the same base distribution with perturbations on top.

Examples

>>> electrons = PICVariable(space='Particles6D')
>>> maxwellian = Maxwellian3D(n=(1.0, None), vth2=(0.1, None))
>>> electrons.add_background(maxwellian)
>>> pert = TorusModesCos(amps=(0.1,))
>>> electrons.add_perturbation(pert)
class struphy.models.variables.FEECVariable(space: Literal['H1', 'Hcurl', 'Hdiv', 'L2', 'H1vec'] = 'H1')[source]#

Grid-based finite element variable using FEEC (Finite Element Exterior Calculus) discretization.

FEECVariable represents field quantities discretized on a computational grid using finite element methods. It supports arbitrary FEEC spaces (H1, Hdiv, Hcurl, L2) and is used for electromagnetic fields, fluid moments, or other spatially distributed quantities.

Variables:
  • space (str) – The FEEC function space. Options: ‘H1’, ‘HDiv’, ‘HCurl’, ‘L2’. Determines the continuity and smoothness properties of the discretization.

  • spline (SplineFunction) – The underlying spline-based representation of the field on the computational mesh. Only available after calling allocate().

  • species (FieldSpecies | FluidSpecies) – Parent species containing this variable.

  • backgrounds (FieldsBackground | None) – Equilibrium fields added to initial conditions.

  • perturbations (Perturbation | list[Perturbation] | None) – Initial perturbations added to backgrounds.

add_background(background)[source]#

Add an equilibrium field background.

add_perturbation(perturbation)[source]#

Add initial perturbations to the field.

allocate(derham, domain, equil, verbose)[source]#

Allocate spline function and initialize on the mesh.

Notes

Initial conditions are constructed by combining backgrounds and perturbations: initial_field = background(s) + perturbation(s)

If neither background nor perturbation is specified, the variable initializes to zero.

Examples

>>> E_field = FEECVariable(space='Hdiv')  # For electric field
>>> E_field.add_background(FieldsBackground(type='uniform', magnitude=1.0))
class struphy.models.variables.SPHVariable[source]#

Fluid variable using Smoothed Particle Hydrodynamics (SPH) discretization.

SPHVariable represents macroscopic fluid quantities (density, velocity, temperature) using SPH methods, where fields are reconstructed from marker positions and properties. This is a particle-based Lagrangian approach to fluid dynamics.

Variables:
  • space (str) – Always ‘ParticlesSPH’ for SPH method.

  • particles (ParticlesSPH) – The SPH marker/particle object. Only available after allocate().

  • particles_class (type) – Reference to the ParticlesSPH class.

  • particle_data (dict) – Dictionary storing particle-associated data fields.

  • species (ParticleSpecies) – Parent SPH species containing this variable.

  • backgrounds (FluidEquilibrium) – Background fluid state (density, velocity, pressure profiles).

  • perturbations (dict[str, Perturbation]) – Perturbations to density (‘n’) and velocity components (‘u1’, ‘u2’, ‘u3’). Each component can have independent perturbations or None.

  • n_as_volume_form (bool) – Always True for SPH; number density represented as volume-weighted quantity.

  • n_to_save (int) – Number of SPH particles for which data is saved.

  • saved_markers (ndarray) – Array storing saved particle data.

add_background(background)[source]#

Set the background fluid equilibrium state.

add_perturbation(del_n, del_u1, del_u2, del_u3)[source]#

Add perturbations to density and/or velocity components.

show_perturbations()[source]#

Display detailed information about density and velocity perturbations.

allocate(derham, domain, equil, projected_equil, verbose)[source]#

Initialize SPH particles and allocate marker arrays.

Notes

Initial conditions combine background and perturbations: - density: background + del_n - velocity: background + (del_u1, del_u2, del_u3)

If neither background nor perturbations are specified, fields initialize to zero.

Examples

>>> fluid = SPHVariable()
>>> equil = HomogenSlab()
>>> fluid.add_background(equil)
>>> pert = TorusModesCos(amps=(0.1,))
>>> fluid.add_perturbation(del_n=pert, del_u1=pert)

Particle parameters#

class struphy.particles.parameters.LoadingParameters(Np: int | None = None, ppc: int | None = None, ppb: int = 10, loading: Literal['pseudo_random', 'sobol_standard', 'sobol_antithetic', 'external', 'restart', 'tesselation'] = 'pseudo_random', seed: int | None = None, moments: tuple | None = None, spatial: Literal['uniform', 'disc'] = 'uniform', specific_markers: tuple[tuple] | None = None, set_zero_velocity: tuple[bool] = (False, False, False), n_quad: int = 1, dir_exrernal: str | None = None, dir_particles: str | None = None, dir_particles_abs: str | None = None, restart_key: str | None = None)[source]#

Configuration for particle (marker) loading strategies and data sources.

This class encapsulates all parameters needed to initialize particles in simulations, including population size, spatial and velocity distributions, loading algorithms, and restart/external data sources. Supports multiple loading strategies: Monte-Carlo based distributions with customizable moments, regular grid tesselation, specific manually-defined markers, and loading from restart or external HDF5 files.

Parameters:
  • Np (int, optional) – Total number of particles to load into the simulation.

  • ppc (int, optional) – Particles per cell to load if a grid is defined. Cell divisions follow domain_array.

  • ppb (int, default=10) – Particles per sorting box. Sorting boxes are defined by boxes_per_dim.

  • loading (LiteralOptions.OptsLoading, default="pseudo_random") – Loading algorithm strategy. Options include various Monte-Carlo methods or ‘tesselation’ for regular grid positioning.

  • seed (int, optional) – Seed for the random number generator for reproducible results. If None, no seed is applied.

  • moments (tuple, optional) – Mean velocities and temperatures defining the Gaussian velocity distribution. If None, automatically computed from the background distribution.

  • spatial (LiteralOptions.OptsSpatialLoading, default="uniform") – Spatial sampling method: ‘uniform’ samples uniformly in (eta1, eta2) coordinates, while ‘disc’ samples uniformly on the disc image of the coordinate space.

  • specific_markers (tuple[tuple], optional) – Manually-defined markers as tuples of phase space coordinates (floats). Each tuple represents a single particle’s initial state.

  • n_quad (int, default=1) – Number of quadrature points used for tesselation-based particle loading.

  • set_zero_velocity (tuple) – Initialize velocity of Maxwellain along selected axis to be zero.

  • dir_external (str, optional) – Absolute path to HDF5 file from which to load external marker data.

  • dir_particles_abs (str, optional) – Absolute path to HDF5 file from which to load restart marker data.

  • dir_particles (str, optional) – Relative path (relative to output folder) to HDF5 restart file for loading markers.

  • restart_key (str, optional) – HDF5 dataset key within the ‘restart/’ folder containing marker array data.

class struphy.particles.parameters.WeightsParameters(control_variate: bool = False, reject_weights: bool = False, threshold: float = 0.0)[source]#

Configuration for particle weight handling and variance reduction.

Manages particle weighting strategies used in Monte-Carlo type simulations, including control variate techniques for variance reduction and weight thresholding to eliminate negligibly-weighted particles.

Parameters:
  • control_variate (bool, default=False) – Whether to apply control variate variance reduction technique to particle weights.

  • reject_weights (bool, default=False) – Whether to filter out particles with weights below the specified threshold.

  • threshold (float, default=0.0) – Minimum weight threshold. Particles with weights below this value are rejected when reject_weights is True.

class struphy.particles.parameters.BoundaryParameters(bc: tuple[Literal['periodic', 'reflect']] = ('periodic', 'periodic', 'periodic'), bc_refill=None, bc_sph: tuple[Literal['periodic', 'mirror', 'fixed']] = ('periodic', 'periodic', 'periodic'))[source]#

Configuration for boundary conditions applied to particles and kernel reconstructions.

Defines how particles behave at domain boundaries (particle boundary conditions) and how smoothed particle hydrodynamics (SPH) kernel reconstructions are handled at domain edges. Supports independent boundary conditions per spatial dimension.

Parameters:
  • bc (tuple[LiteralOptions.OptsMarkerBC], default=("periodic", "periodic", "periodic")) – Particle boundary conditions for each spatial direction (3D). Options per direction: ‘remove’ (delete particles), ‘reflect’ (specular reflection), ‘periodic’ (wrap around), or ‘refill’ (reload particles).

  • bc_refill (list, optional) – Refill strategy when ‘refill’ boundary condition is active. Either ‘inner’ or ‘outer’.

  • bc_sph (tuple[LiteralOptions.OptsRecontructBC], default=("periodic", "periodic", "periodic")) – Boundary conditions for SPH kernel reconstruction in each spatial direction. Typically matches or differs from bc depending on reconstruction needs.

class struphy.particles.parameters.BinningPlot(slice: str = 'e1', n_bins: int | tuple[int] = 128, ranges: tuple[float] | tuple[tuple[float]] = (0.0, 1.0), divide_by_jac: bool = True, output_quantity: Literal['density', 'current_1', 'current_2', 'current_3', 'energy_tensor_11', 'energy_tensor_22', 'energy_tensor_33', 'energy_tensor_12', 'energy_tensor_13', 'energy_tensor_23', 'heat_flux_1', 'heat_flux_2', 'heat_flux_3'] = 'density')[source]#

Configuration for particle phase-space binning and histogram generation.

Produces binned distributions from particle data across specified phase-space coordinates. Supports arbitrary combinations of spatial (eta) and velocity (v) coordinates with flexible binning resolution and coordinate ranges. Automatically computes bin edges and allocates storage for full and delta-f distributions.

Parameters:
  • slice (str, default="e1") – Phase-space coordinates to bin, specified as underscore-separated coordinate names (e.g., ‘e1’, ‘e1_e2’, ‘e1_v1’). Valid coordinates: ‘e1’, ‘e2’, ‘e3’ (spatial), ‘v1’, ‘v2’, ‘v3’ (velocity). Example: ‘e1’ produces 1D binning over eta1; ‘e1_v1’ produces 2D binning over eta1 and v1.

  • n_bins (int | tuple[int], default=128) – Number of bins per coordinate. If int, applies to all coordinates; if tuple, specifies bins for each coordinate separately.

  • ranges (tuple[float] | tuple[tuple[float]], default=(0.0, 1.0)) – Binning ranges as intervals [min, max] for each coordinate. If a single tuple, applies to all coordinates; if nested tuples, specifies range for each coordinate.

  • divide_by_jac (bool, default=True) – Whether to normalize distributions by the Jacobian determinant (volume-to-0-form conversion). Set False to use unnormalized bin counts.

  • output_quantity (LiteralOptions.BinningQuantity, default="density") – Quantity to compute in binning: determines weighting scheme and output format.

class struphy.particles.parameters.KernelDensityPlot(pts_e1: int = 16, pts_e2: int = 16, pts_e3: int = 1)[source]#

Configuration for smoothed particle hydrodynamics (SPH) density reconstructions.

Evaluates particle density fields at structured grid points using SPH kernel interpolation. Supports 1D, 2D, and 3D evaluations with independent resolution control per dimension.

Parameters:
  • pts_e1 (int, default=16) – Number of evaluation grid points in the first spatial direction (eta1).

  • pts_e2 (int, default=16) – Number of evaluation grid points in the second spatial direction (eta2).

  • pts_e3 (int, default=1) – Number of evaluation grid points in the third spatial direction (eta3). Set to 1 for 2D density plots.

Fields background#

class struphy.io.options.FieldsBackground(type: Literal['LogicalConst', 'FluidEquilibrium'] = 'LogicalConst', values: tuple = (1.5, 0.7, 2.4), variable: str | None = None)[source]#

Set options for static fluid backgrounds/equilibria in parameter/launch files.

Parameters:
  • type (BackgroundTypes) – Type of background.

  • values (tuple[float]) – Values for LogicalConst on the unit cube. Can be length 1 for scalar functions; must be length 3 for vector-valued functions.

  • variable (str) – Name of the method in FluidEquilibrium that should be the background.

Perturbation functions#

class struphy.initial.base.Perturbation[source]#

Abstract base class for perturbation functions used as initial conditions in simulations.

This class provides the interface and common functionality for defining perturbation fields in logical (eta) or physical coordinate spaces. Subclasses must implement the __call__ method to define the perturbation as a callable function of spatial coordinates.

The class supports flexible representation bases (p-forms, vector fields, physical coordinates) and allows specification of which component is perturbed for vector-valued quantities.

Variables:
  • given_in_basis (str) –

    Specifies the basis representation of the perturbation. Options:

    • ’0’, ‘1’, ‘2’, ‘3’ : Differential form basis (0-form=scalar, 1-form, etc.)

    • ’v’ : Vector field basis

    • ’physical’ : Physical (mapped) domain coordinates

    • ’physical_at_eta’ : Physical components evaluated in logical (eta) domain, u(F(eta))

    • ’norm’ : Normalized co-variant basis (\(delta_i / |delta_i|\))

  • comp (int) – Component index for vector-valued perturbations (0-2 for vector components, 0 for scalar-valued functions). Default is 0.

Examples

Subclasses should override __call__ to implement specific perturbation fields:

>>> class CustomPerturbation(Perturbation):
...     def __init__(self):
...         self.given_in_basis = 'physical'
...     def __call__(self, eta1, eta2, eta3, flat_eval=False):
...         return eta1 * eta2  # Example perturbation field

See available Perturbation functions.