Particle Swarms

Swarm

class underworld3.swarm.Swarm[source]

Bases: Stateful, uw_object

A basic particle swarm implementation for Lagrangian particle tracking and data storage.

The UW Swarm class provides a simplified particle management system that uses PETSc’s DMSWARM_BASIC type. Unlike the standard Swarm class, this implementation does not rely on PETSc to determine ranks for particle migration but instead uses our own kdtree neighbour-domain computations.

This class is preferred for most operations except where particle / cell relationships are always required.

Parameters:
  • mesh (uw.discretisation.Mesh) – The mesh object that defines the computational domain for particle operations. Particles will be associated with this mesh for spatial queries and operations.

  • recycle_rate (int, optional) – Rate at which particles are recycled for streak management. If > 1, enables streak particle functionality where particles are duplicated and tracked across multiple cycles. Default is 0 (no recycling).

  • verbose (bool, optional) – Enable verbose output for debugging and monitoring particle operations. Default is False.

Examples

Create a basic swarm and populate with particles:

>>> import underworld3 as uw
>>> mesh = uw.meshing.UnstructuredSimplexBox(minCoords=(0,0), maxCoords=(1,1))
>>> swarm = uw.swarm.Swarm(mesh=mesh)
>>> swarm.populate(fill_param=2)

Create a streak swarm with recycling:

>>> streak_swarm = uw.swarm.Swarm(mesh=mesh, recycle_rate=5)
>>> streak_swarm.populate(fill_param=1)

Add custom particle data:

>>> temperature = swarm.add_variable("temperature", 1)
>>> velocity = swarm.add_variable("velocity", mesh.dim)

Manual particle migration after coordinate updates:

Note: particle migration is still called automatically when we access and update the particle_coordinates variables

Note: swarm.populate uses a the mesh point locations for discontinuous interpolants to determine the particle locations.

instances = 0
__init__(mesh, recycle_rate=0, verbose=False, clip_to_mesh=True)[source]
property mesh

The mesh this swarm operates on

property local_size

Number of particles on this MPI rank.

Returns:

Local particle count.

Return type:

int

See also

dm.getLocalSize

Underlying PETSc method.

property data

Particle coordinates (alias for points).

Deprecated since version 0.99.0: Use direct DM field access for particle coordinates.

Returns:

Particle coordinate array of shape (n_particles, dim).

Return type:

numpy.ndarray

property points

Swarm particle coordinates in physical units.

Deprecated since version 0.99.0: Use swarm variables or direct DM access instead. swarm.points is being deprecated.

When the mesh has coordinate scaling applied (via model units), this property automatically converts from internal model coordinates to physical coordinates for user access.

When the mesh has coordinate units specified, returns a unit-aware array.

Returns:

numpy.ndarray or UnitAwareArray: Particle coordinates (with units if mesh.units is set)

property coords

Swarm particle coordinates in physical units.

This is the primary public interface for accessing particle coordinates. Coordinates are automatically converted from internal model units to physical units based on the model’s reference quantities.

Returns:

Particle coordinates in physical units with shape (n_particles, dim). If model has reference quantities, returns UWQuantity with appropriate length units. Otherwise returns plain array.

Return type:

UWQuantity or numpy.ndarray

Notes

  • Coordinates are converted from model units to physical units automatically

  • For internal use with model units, access swarm._particle_coordinates.data

  • Setting coordinates accepts either physical units or plain numbers

Examples

>>> coords_physical = swarm.coords  # Get physical coordinates
>>> swarm.coords = new_coords_with_units  # Set from physical units

See also

swarm.units

Get the unit specification for coordinates

property units

Unit specification for swarm coordinates.

Returns the physical unit string for coordinates based on the model’s reference quantities. This indicates what units the coordinates are in when accessed via the coords property.

Returns:

Unit string for coordinates (e.g., ‘kilometer’, ‘meter’), or None if no reference quantities are set

Return type:

str or None

Examples

>>> print(swarm.units)  # 'kilometer' if length_scale was set in km
>>> coords = swarm.coords  # Coordinates in kilometers
property clip_to_mesh

Whether particles are clipped to remain within mesh boundaries.

When True (default), particles that move outside the mesh domain during advection or coordinate updates are removed or repositioned to stay within bounds. When False, particles can exist outside the mesh domain.

Returns:

Current clipping state.

Return type:

bool

See also

dont_clip_to_mesh

Context manager to temporarily disable clipping.

dont_clip_to_mesh()[source]

Context manager that temporarily disables mesh clipping for the swarm. swarm.migrate is called automatically when exiting the context.

Usage:
with swarm.dont_clip_to_mesh():

# swarm operations that should not be clipped to mesh swarm.data = new_positions

migration_disabled()[source]

Legacy context manager that completely disables migration. Use migration_control(disable=True) for new code.

Context manager that temporarily disables particle migration for the swarm. Migration is NOT called when exiting the context.

Usage:
with swarm.migration_disabled():

# swarm operations that should not trigger migration swarm.data = new_positions # … other operations … # migrate() will be skipped during these operations

migration_control(disable=False)[source]

Context manager to control particle migration behavior.

Parameters:

disable (bool) – If False (default), migration is deferred until context exit. If True, migration is completely disabled.

Examples

Defer migration until end (default):

with swarm.migration_control():
    swarm.points[mask1] += delta1
    swarm.points[mask2] *= scale
    # Migration happens HERE on exit

Completely disable migration:

with swarm.migration_control(disable=True):
    # Operations where migration should never happen
    # No migration on exit
populate(**kwargs)

Populate the swarm with particles throughout the domain.

Parameters:

fill_param – Parameter determining the particle count per cell (per dimension) for the given layout, using the mesh degree.

Raises:

RuntimeError – If the swarm has already been initialized with particles.

migrate(**kwargs)

Migrate swarm across processes after coordinates have been updated.

The algorithm uses a global kD-tree for the centroids of the domains to decide the particle mpi.rank (send to the closest) If the particles are mis-assigned to a particular mpi.rank, the next choice is the second-closest and so on.

A few particles are still not found after this distribution process which probably means they are just outside the mesh. If some points remain lost, they will be deleted if delete_lost_points is set.

Implementation note:

We retained (above) the name DMSwarmPIC_coor for the particle field to allow this routine to be inherited by a PIC swarm which has this field pre-defined. (We’d need to add a cellid field as well, and re-compute it upon landing)

Note: This is a COLLECTIVE operation - all MPI ranks must call it.

add_particles_with_coordinates(coordinatesArray)[source]

Add particles at given coordinates, keeping only locally-owned points.

Each rank filters the input array and adds only the points that fall within its local domain partition. Non-local points are silently ignored, so it is safe to pass the same global coordinate array to every rank — no duplicates will be created.

This is the recommended method for user code. For pre-partitioned data where each rank already holds only its own points, this method also works correctly.

Parameters:

coordinatesArray (numpy.ndarray) – Coordinates of new particles, shape (n, dim).

Returns:

npoints – Number of points actually added on this rank.

Return type:

int

add_particles_with_global_coordinates(globalCoordinatesArray, migrate=True, delete_lost_points=True)[source]

Add particles on every rank, then migrate to correct owners.

Every rank inserts all supplied points into the local swarm, then dm.migrate() moves each particle to the rank that owns its spatial location. If the same array is passed on every rank, this produces one copy of each point in the correct partition — but calling it with different arrays per rank will accumulate all of them.

This is primarily an internal method used by global-evaluation and mesh-transfer utilities. For general use, prefer add_particles_with_coordinates(), which filters non-local points automatically and never creates duplicates.

Parameters:
  • globalCoordinatesArray (numpy.ndarray) – Coordinates of new particles, shape (n, dim).

  • migrate (bool) – Run PETSc swarm migration after insertion (default True).

  • delete_lost_points (bool) – Remove particles that fall outside any rank’s domain during migration (default True).

Returns:

npoints – Number of points added on this rank before migration.

Return type:

int

save(filename, compression=False, compressionType='gzip', force_sequential=False)[source]

Save the swarm coordinates to a h5 file.

Parameters:
  • filename (int) – The filename of the swarm checkpoint file to save to disk.

  • compression (bool | None) – Add compression to the h5 files (saves space but increases write times with increasing no. of processors)

  • compressionType (str | None) – Type of compression to use, ‘gzip’ and ‘lzf’ supported. ‘gzip’ is default. Compression also needs to be set to ‘True’.

read_timestep(base_filename, swarm_id, index, outputPath='', migrate=True)[source]
Parameters:
  • base_filename (str)

  • swarm_id (str)

  • index (int)

  • outputPath (str | None)

add_variable(name, size=1, dtype=<class 'float'>, proxy_degree=2, _nn_proxy=False, units=None)[source]

Add a variable to the swarm.

Variables must be created before the swarm is populated with particles. Once swarm.populate() or similar methods are called, PETSc finalizes field registration and no new variables can be added.

Parameters:
  • name (str) – Variable name

  • size (int, default 1) – Number of components (1 for scalar, 2-3 for vector, etc.)

  • dtype (type, default float) – Data type (float or int)

  • proxy_degree (int, default 2) – Degree for mesh proxy variable interpolation

  • _nn_proxy (bool, default False) – Internal parameter for nearest-neighbor proxy

  • units (str, optional) – Physical units for this variable (e.g., “kg/m^3”, “m/s”)

Returns:

The created swarm variable

Return type:

SwarmVariable

Raises:

RuntimeError – If swarm is already populated with particles

Examples

Correct usage: >>> swarm = uw.swarm.Swarm(mesh) >>> material = swarm.add_variable(“material”, 1, dtype=int) >>> temperature = swarm.add_variable(“temperature”, 1) >>> swarm.populate(fill_param=3) # Populate after creating variables

Incorrect usage (will raise error): >>> swarm = uw.swarm.Swarm(mesh) >>> swarm.populate(fill_param=3) >>> material = swarm.add_variable(“material”, 1) # ERROR!

petsc_save_checkpoint(swarmName, index, outputPath='')[source]

Use PETSc to save the swarm and attached data to a .pbin and xdmf file.

Parameters:
  • swarmName (str) – Name of the swarm to save.

  • index (int) – An index which might correspond to the timestep or output number (for example).

  • outputPath (str | None) – Path to save the data. If left empty it will save the data in the current working directory.

write_timestep(filename, swarmname, index, swarmVars=None, outputPath='', time=None, compression=False, compressionType='gzip', force_sequential=False)[source]

Save data to h5 and a corresponding xdmf for visualisation using h5py.

Parameters:
  • swarmName – Name of the swarm to save.

  • swarmVars (list | None) – List of swarm objects to save.

  • index (int) – An index which might correspond to the timestep or output number (for example).

  • outputPath (str | None) – Path to save the data. If left empty it will save the data in the current working directory.

  • time (int | None) – Attach the time to the generated xdmf.

  • compression (bool | None) – Whether to compress the h5 files [bool].

  • compressionType (str | None) – The type of compression to use. ‘gzip’ and ‘lzf’ are the supported types, with ‘gzip’ as the default.

  • filename (str)

  • swarmname (str)

  • force_sequential (bool | None)

property vars

List of SwarmVariables attached to this swarm.

Returns:

List of SwarmVariable objects defined on this swarm.

Return type:

list

snapshot_payload()[source]

Return a self-contained dict describing this swarm’s state.

Captured: per-rank particle coordinates (from DMSwarmPIC_coor) and every user swarm-variable’s data array. PETSc-internal variables (DMSwarmPIC_coor, DMSwarm_X0, DMSwarm_remeshed) are excluded — their contents either come from the captured coords or are regenerated on the next solve.

Return type:

dict

apply_snapshot_payload(payload)[source]

Rebuild this swarm’s local particle population from a payload.

Algorithm:

  1. Drop every current local particle (dm.removePoint from the end is O(1) per call, O(N) total).

  2. Add the captured-rank’s particles back via the raw PETSc primitives — addNPoints then writing the coord field directly. We deliberately bypass add_particles_with_coordinates() because that method filters via points_in_domain (slow) and triggers dm.migrate (unnecessary — saved coords were already local at capture time, and the mesh hasn’t changed).

  3. Write captured per-variable data back. The local particle count matches the captured count because we just put the same particles back in the same order.

This bumps _population_generation once (from the addNPoints step in restore), which is correct: the population did just change. Downstream consumers that care can compare against the captured value in payload['captured_population_generation'].

Parameters:

payload (dict)

Return type:

None

access(*writeable_vars)[source]

Dummy access manager that provides deferred sync for backward compatibility. Uses NDArray_With_Callback.delay_callbacks_global() internally.

This is a compatibility wrapper that allows existing code using the access() context manager to work with the new direct-access variable interfaces. All variable modifications are deferred and synchronized at context exit.

Parameters:

writeable_vars (SwarmVariable) – Variables that will be modified (ignored - all variables are writable with the new interface, this parameter is kept for API compatibility)

Return type:

Context manager that defers variable synchronization until exit

Notes

This method is deprecated. New code should access variable.data or variable.array directly without requiring an access context.

advection(V_fn, delta_t, order=2, corrector=False, restore_points_to_domain_func=None, evalf=False, step_limit=False)[source]
estimate_dt(V_fn)[source]

Calculates an appropriate advective timestep for the given mesh and velocity configuration.

SwarmVariable

class underworld3.swarm.SwarmVariable[source]

Bases: DimensionalityMixin, MathematicalMixin, Stateful, uw_object

Variable supported by a particle swarm (point cloud).

A SwarmVariable stores values at discrete particle locations and provides a mesh-based proxy representation for use in symbolic expressions. This enables Lagrangian tracking of material properties through deformation.

Parameters:
  • name (str) – Identifier for this variable (must be unique within the swarm).

  • swarm (Swarm) – The supporting particle swarm.

  • size (int or tuple, optional) – Shape specification: int for vectors, tuple for matrices. If None, inferred from vtype.

  • vtype (VarType, optional) – Variable type (SCALAR, VECTOR, TENSOR, SYM_TENSOR, MATRIX). If None, inferred from size.

  • dtype (type, default=float) – Data type for storage (float or int).

  • proxy_degree (int, default=1) – Polynomial degree for the mesh proxy variable.

  • proxy_continuous (bool, default=True) – Whether the proxy uses continuous (True) or discontinuous (False) interpolation.

  • varsymbol (str, optional) – LaTeX symbol for display. Defaults to name.

  • rebuild_on_cycle (bool, default=True) – If True, rebuild the proxy when particles cycle through periodic boundaries. Recommended for continuous fields.

  • units (str or pint.Unit, optional) – Physical units for this variable (e.g., ‘kelvin’, ‘Pa’). Requires reference quantities to be set on the model.

See also

MeshVariable

Variable supported by mesh nodes.

Swarm

Container for particle locations.

Examples

Create a temperature field on a swarm:

>>> swarm = uw.swarm.Swarm(mesh)
>>> T = swarm.add_variable("T", size=1, vtype=uw.VarType.SCALAR)
>>> T.data[:] = 1600.0  # Set initial temperature

Create a velocity field:

>>> v = swarm.add_variable("v", size=mesh.dim, vtype=uw.VarType.VECTOR)

Notes

SwarmVariables are essential for tracking material properties that advect with the flow. The mesh proxy enables their use in finite element formulations while particle storage preserves Lagrangian history.

__init__(name, swarm, size=None, vtype=None, dtype=<class 'float'>, proxy_degree=1, proxy_continuous=True, _register=True, _proxy=True, _nn_proxy=False, varsymbol=None, rebuild_on_cycle=True, units=None, units_backend=None)[source]

Initialize dimensionality tracking attributes

property units

Return the units associated with this variable.

property has_units

Check if this variable has units.

use_legacy_array()[source]

Deprecated: Array interface is now unified using NDArray_With_Callback

use_enhanced_array()[source]

Deprecated: Array interface is now unified using NDArray_With_Callback

sync_disabled(description='batch operation')[source]

Context manager to disable automatic synchronization for batch operations. Now uses NDArray_With_Callback’s delay_callback mechanism.

Parameters:

description (str) – Description of the batch operation for debugging

pack_uw_data_to_petsc(data_array, sync=True)[source]

Enhanced pack method that directly accesses PETSc field without access() context. Designed for the new swarmVariable.array interface.

Parameters:
  • data_array (numpy.ndarray) – Array data to pack into PETSc field

  • sync (bool) – Whether to sync parallel operations (default True)

unpack_uw_data_from_petsc(squeeze=True, sync=True)[source]

Enhanced unpack method that directly accesses PETSc field without access() context. Designed for the new swarmVariable.array interface.

Parameters:
  • squeeze (bool) – Whether to squeeze singleton dimensions (default True)

  • sync (bool) – Whether to sync parallel operations (default True)

pack_raw_data_to_petsc(data_array, sync=True)[source]

Pack data array to PETSc using traditional data shape (-1, num_components). Direct PETSc access without access() context for backward compatibility.

Parameters:
  • data_array (numpy.ndarray) – Array data in traditional flat format (-1, num_components)

  • sync (bool) – Whether to sync parallel operations (default True)

unpack_raw_data_from_petsc(squeeze=True, sync=True)[source]

Unpack data from PETSc in traditional data shape (-1, num_components). Direct PETSc access without access() context for backward compatibility.

Parameters:
  • squeeze (bool) – Whether to remove singleton dimensions (default True)

  • sync (bool) – Whether to sync parallel operations (default True)

Returns:

Array data in traditional flat format (-1, num_components)

Return type:

numpy.ndarray

rbf_interpolate(new_coords, verbose=False, nnn=None)[source]

Radial basis function interpolation of particle data to arbitrary points.

Uses inverse-distance weighting to interpolate particle values to new coordinate locations.

Parameters:
  • new_coords (numpy.ndarray) – Target coordinates of shape (N, dim) to interpolate to.

  • verbose (bool, default=False) – Print diagnostic information during interpolation.

  • nnn (int, optional) – Number of nearest neighbors to use. Defaults to mesh.dim + 1.

Returns:

Interpolated values at the target coordinates.

Return type:

numpy.ndarray

property swarm

The swarm this variable belongs to (accessed via weak reference). Raises RuntimeError if the swarm has been garbage collected.

property old_data

Original data property implementation.

Type:

TESTING

property data

Canonical data storage in flat format for internal operations.

Returns particle data in shape (-1, num_components) regardless of variable type. Values are always non-dimensional (no unit conversion).

This property handles PETSc synchronization. The .array property is a view of this with shape conversion.

When to Use

  • Variable-to-variable transfers: Copying data between swarm variables avoids redundant unit conversions

  • Low-level operations: Direct access to particle data

  • Backward compatibility: Existing code using flat format

For general user access, prefer .array which provides a structured shape.

returns:

Array with shape (-1, num_components) with automatic PETSc sync.

rtype:

NDArray_With_Callback

Examples

>>> # Efficient variable-to-variable copy
>>> new_material.data[...] = old_material.data[...]
>>> # Check data shape
>>> scalar_var.data.shape  # (n_particles, 1)
>>> vector_var.data.shape  # (n_particles, dim)

See also

array

Structured format (N, a, b) for general access.

sym

Symbolic representation for equations.

property array

Primary interface for reading and writing particle data.

Returns a structured array view with shape (N, a, b) where N is the number of particles and (a, b) depends on the variable type:

  • Scalar: (N, 1, 1)

  • Vector: (N, 1, dim)

  • Tensor: (N, dim, dim)

Returns:

Array view that delegates changes back to canonical storage.

Return type:

NDArray

Examples

>>> # Scalar field initialization
>>> material_property.array[:, 0, 0] = 1e21  # Set viscosity
>>> # Vector field
>>> velocity.array[:, 0, 0] = vx_values  # x-component
>>> velocity.array[:, 0, 1] = vy_values  # y-component
>>> # Reading values
>>> max_visc = material_property.array[:, 0, 0].max()

Notes

This property is a view of the canonical .data property with automatic shape conversion. All modifications trigger proxy variable updates for mesh interpolation.

See also

data

Flat format (-1, components) for variable-to-variable transfers.

sym

Symbolic representation for use in equations.

property sym

Symbolic representation for use in equations.

Returns the symbolic expression from the proxy mesh variable, which can be used in SymPy expressions for constitutive models, boundary conditions, etc.

Returns:

Symbolic matrix expression.

Return type:

sympy.Matrix

Notes

The proxy is automatically updated if particle data has changed.

property sym_1d

Flattened symbolic representation.

Returns the symbolic expression as a 1D (column) vector form, useful for Voigt notation in tensor calculations.

Returns:

Flattened symbolic expression.

Return type:

sympy.Matrix

global_max(**kwargs)

Maximum value across all MPI ranks.

Finds the maximum value of the particle property across all processors. Useful for finding extreme values in particle swarm data.

Parameters:
  • axis (None, int, or tuple of ints, optional) – Axis or axes along which to operate. By default, flattened input is used.

  • out (None, optional) – Alternative output array (not supported, kept for API compatibility).

  • keepdims (bool, optional) – If True, reduced axes are left as dimensions with size one.

Returns:

Maximum value with units preserved (if variable has units).

Return type:

UWQuantity or scalar

Examples

>>> max_temp = temperature_swarm.global_max()
>>> print(f"Maximum temperature: {max_temp}")

Notes

This is a collective operation - all ranks must call it. The result is identical on all ranks.

global_min(**kwargs)

Minimum value across all MPI ranks.

Finds the minimum value of the particle property across all processors. Useful for finding extreme values in particle swarm data.

Parameters:
  • axis (None, int, or tuple of ints, optional) – Axis or axes along which to operate. By default, flattened input is used.

  • out (None, optional) – Alternative output array (not supported, kept for API compatibility).

  • keepdims (bool, optional) – If True, reduced axes are left as dimensions with size one.

Returns:

Minimum value with units preserved (if variable has units).

Return type:

UWQuantity or scalar

Examples

>>> min_pressure = pressure_swarm.global_min()
>>> print(f"Minimum pressure: {min_pressure}")

Notes

This is a collective operation - all ranks must call it. The result is identical on all ranks.

global_sum(**kwargs)

Sum of values across all MPI ranks.

Computes the sum of particle property values across all processors.

Parameters:
  • axis (None, int, or tuple of ints, optional) – Axis or axes along which to operate. By default, flattened input is used.

  • out (None, optional) – Alternative output array (not supported, kept for API compatibility).

  • keepdims (bool, optional) – If True, reduced axes are left as dimensions with size one.

Returns:

Sum with units preserved (if variable has units).

Return type:

UWQuantity or scalar

Notes

This is a collective operation - all ranks must call it. The result is identical on all ranks.

Warning: This sum is NOT a physical domain-integrated quantity because particles are non-uniformly distributed. For domain integration, use the proxy mesh variable with uw.maths.Integral().

global_norm(**kwargs)

L2 norm (Frobenius norm) across all MPI ranks.

Computes the L2 norm of particle property values: sqrt(sum(x**2)) across all processors.

Parameters:

ord ({non-zero int, inf, -inf, 'fro', 'nuc'}, optional) – Order of the norm (default: None = 2-norm)

Returns:

L2 norm with units preserved (if variable has units).

Return type:

UWQuantity or scalar

Notes

This is a collective operation - all ranks must call it. The result is identical on all ranks.

For vectors, computes the Frobenius norm treating the array as flattened.

Warning: This norm is NOT a physical domain-integrated quantity because particles are non-uniformly distributed.

global_size(**kwargs)

Total particle count across all MPI ranks.

Returns the total number of particles across all processors. Useful for population monitoring and load balancing diagnostics.

Returns:

Total number of particles across all ranks.

Return type:

int

Examples

>>> total_particles = swarm_var.global_size()
>>> local_particles = swarm_var.data.shape[0]
>>> print(f"Rank has {local_particles} of {total_particles} particles")

Notes

This is a collective operation - all ranks must call it. The result is identical on all ranks.

save(filename, compression=False, compressionType='gzip', force_sequential=False)[source]

Save the swarm variable to a h5 file.

Parameters:
  • filename (int) – The filename of the swarm variable to save to disk.

  • compression (bool | None) – Add compression to the h5 files (saves space but increases write times with increasing no. of processors)

  • compressionType (str | None) – Type of compression to use, ‘gzip’ and ‘lzf’ supported. ‘gzip’ is default. Compression also needs to be set to ‘True’.

  • force_sequential (activate the serial version of hdf5)

write_proxy(filename)[source]
Parameters:

filename (str)

read_timestep(data_filename, swarmID, data_name, index, outputPath='')[source]
Parameters:
  • data_filename (str)

  • swarmID (str)

  • data_name (str)

  • index (int)

Index Swarm Variable

class underworld3.swarm.IndexSwarmVariable[source]

Bases: SwarmVariable

Integer-valued swarm variable for material tracking.

IndexSwarmVariable stores integer indices at particle locations, typically used for tracking distinct material types. It automatically generates symbolic mask expressions for each material index, enabling material- dependent properties in constitutive models.

Parameters:
  • name (str) – Variable name for identification and I/O.

  • swarm (Swarm) – Parent swarm object.

  • indices (int) – Number of distinct material indices (default 1).

  • proxy_degree (int) – Polynomial degree for mesh projection (default 1).

  • proxy_continuous (bool) – Whether mesh proxy is continuous (default True).

Examples

>>> material = IndexSwarmVariable("M", swarm, indices=3)
>>> material.data[:] = 0  # Set all particles to material 0
>>> # Use sym[i] as multiplier for material i properties
>>> viscosity = material.sym[0] * 1e20 + material.sym[1] * 1e21

See also

SwarmVariable

Base class for particle-supported variables.

__init__(name, swarm, indices=1, proxy_degree=1, proxy_continuous=True, update_type=0, npoints=5, radius=0.5, npoints_bc=2, ind_bc=None, varsymbol=None)[source]

Initialize dimensionality tracking attributes

property sym

Lazy evaluation of symbolic mask array.

Only updates proxy variables when they’re actually needed (when sym is accessed) and only if the proxy variables are marked as stale due to data changes. This avoids expensive RBF interpolation during data assignment operations.

property sym_1d

One-dimensional symbolic mask array (alias for sym).

Returns the same symbolic mask array as sym, provided for API compatibility with other variable types that distinguish between multi-dimensional and flattened representations.

Returns:

Symbolic mask array of shape (indices, 1).

Return type:

sympy.Matrix

See also

sym

Primary symbolic mask array access.

createMask(funcsList)[source]

Create a material-weighted symbolic expression from per-material values.

This method creates a SymPy expression that combines multiple material properties using the index variable’s symbolic masks. The result can be used directly in Underworld’s solver equations.

Parameters:

funcsList (list or tuple) – List of values or symbolic expressions, one per material index. Length must equal indices.

Returns:

Symbolic expression: sum(funcsList[i] * mask[i] for i in indices).

Return type:

sympy.Basic

Raises:

RuntimeError – If funcsList is not a list/tuple or has wrong length.

Examples

>>> # Define viscosity per material
>>> viscosity = material.createMask([1e21, 1e20, 1e22])  # 3 materials
>>> # Use in solver
>>> stokes.constitutive_model.viscosity = viscosity

See also

visMask

Create visualization mask showing material indices.

viewMask(expr)[source]

Decompose a masked expression into per-material components.

Note

This method is not yet implemented. Currently returns None.

Takes a symbolic expression created by createMask() and extracts the individual material-specific components.

Parameters:

expr (sympy.Basic) – A masked symbolic expression created by createMask().

Returns:

List of symbolic expressions, one per material index. Currently returns None (not implemented).

Return type:

list or None

See also

createMask

Create a masked expression from per-material values.

visMask()[source]

Create a visualization mask showing material indices.

Returns a symbolic expression where each material region shows its index value (0, 1, 2, …). Useful for visualization and debugging of material distributions.

Returns:

Symbolic expression evaluating to material index at each point.

Return type:

sympy.Basic

Examples

>>> vis_field = material.visMask()
>>> values = uw.function.evaluate(vis_field, swarm.data)
>>> # values[i] gives material index at particle i

See also

createMask

Create arbitrary material-weighted expressions.

view()[source]

Show information on IndexSwarmVariable

Nodal Point Swarm

A swarm type where particles are located at mesh nodes.

class underworld3.swarm.NodalPointSwarm[source]

Bases: Swarm

BASIC_Swarm with particles located at the coordinate points of a meshVariable

The swarmVariable X0 is defined so that the particles can “snap back” to their original locations after they have been moved.

The purpose of this Swarm is to manage sample points for advection schemes based on upstream sampling (method of characteristics etc)

__init__(trackedVariable, verbose=False)[source]
Parameters:

trackedVariable (EnhancedMeshVariable)

advection(V_fn, delta_t, order=2, corrector=False, restore_points_to_domain_func=None, evalf=False, step_limit=True)[source]