Python API Reference¶
This section provides the API reference for the µGrid Python module. The Python bindings provide access to the core µGrid functionality with a Pythonic interface.
Field Collections¶
Field collections manage groups of fields on structured grids.
GlobalFieldCollection¶
- class muGrid.GlobalFieldCollection(nb_grid_pts, nb_sub_pts=None, nb_ghosts_left=None, nb_ghosts_right=None, device=None)¶
A GlobalFieldCollection manages a set of fields that share the same global grid structure. It can allocate fields in either host (CPU) or device (GPU) memory.
- Parameters:
nb_grid_pts (sequence of int) – Grid dimensions, e.g.,
[64, 64]for 2D or[32, 32, 32]for 3D.nb_sub_pts (dict, optional) – Number of sub-points per pixel for each sub-point type. Default is
{}.nb_ghosts_left (sequence of int, optional) – Ghost cells on low-index side. Default is no ghosts.
nb_ghosts_right (sequence of int, optional) – Ghost cells on high-index side. Default is no ghosts.
device (str or Device, optional) – Device for field allocation. Can be a string (
"cpu","cuda","cuda:N","rocm","rocm:N") or aDeviceinstance. Default is"cpu".
Example:
>>> fc = muGrid.GlobalFieldCollection([64, 64]) >>> field = fc.real_field("temperature") >>> field.p[:] = 300.0 # Set temperature to 300K >>> # GPU field collection >>> fc_gpu = muGrid.GlobalFieldCollection([64, 64], device="cuda")
- real_field(name, components=(), sub_pt='pixel')¶
Create a real-valued field.
- complex_field(name, components=(), sub_pt='pixel')¶
Create a complex-valued field.
- int_field(name, components=(), sub_pt='pixel')¶
Create an integer field.
- uint_field(name, components=(), sub_pt='pixel')¶
Create an unsigned integer field.
LocalFieldCollection¶
- class muGrid.LocalFieldCollection(spatial_dim, name='', nb_sub_pts=None, device=None)¶
A LocalFieldCollection manages fields on a subset of pixels, typically used for material-specific data in heterogeneous simulations.
- Parameters:
spatial_dim (int) – Spatial dimension (2 or 3).
name (str, optional) – Name for the collection. Default is
"".nb_sub_pts (dict, optional) – Number of sub-points per pixel for each sub-point type.
device (str or Device, optional) – Device for field allocation. Can be a string (
"cpu","cuda","cuda:N","rocm","rocm:N") or aDeviceinstance. Default is"cpu".
The field creation methods (
real_field,complex_field,int_field,uint_field) are the same as forGlobalFieldCollection.
Fields¶
The Field class wraps C++ fields with convenient numpy array access.
Field¶
- class muGrid.Field(cpp_field)¶
Python wrapper for muGrid fields providing numpy/cupy array views.
This class wraps a C++ muGrid field and provides the following properties for accessing the underlying data as arrays:
s: SubPt layout, excluding ghost regionssg: SubPt layout, including ghost regionsp: Pixel layout, excluding ghost regionspg: Pixel layout, including ghost regions
For CPU fields, the arrays are numpy arrays. For GPU fields (CUDA/ROCm), the arrays are CuPy arrays. Both are views into the underlying C++ data, so modifications to the arrays will modify the field data directly (zero-copy).
- Parameters:
cpp_field – The underlying C++ field object.
- property s¶
SubPt layout array excluding ghost regions.
Shape:
(*components_shape, nb_sub_pts, *spatial_dims_without_ghosts)
- property sg¶
SubPt layout array including ghost regions.
Shape:
(*components_shape, nb_sub_pts, *spatial_dims_with_ghosts)
- property p¶
Pixel layout array excluding ghost regions.
Shape:
(nb_components * nb_sub_pts, *spatial_dims_without_ghosts)
- property pg¶
Pixel layout array including ghost regions.
Shape:
(nb_components * nb_sub_pts, *spatial_dims_with_ghosts)
wrap_field¶
Domain Decomposition¶
Classes for parallel domain decomposition with MPI.
CartesianDecomposition¶
- class muGrid.CartesianDecomposition(communicator, nb_domain_grid_pts, nb_subdivisions=None, nb_ghosts_left=None, nb_ghosts_right=None, nb_sub_pts=None, device=None)¶
CartesianDecomposition manages domain decomposition for MPI-parallel computations on structured grids, including ghost buffer regions for stencil operations.
- Parameters:
communicator (Communicator) – MPI communicator for parallel execution.
nb_domain_grid_pts (sequence of int) – Global domain grid dimensions.
nb_subdivisions (sequence of int, optional) – Number of subdivisions in each dimension. Default is automatic.
nb_ghosts_left (sequence of int, optional) – Ghost cells on low-index side. Default is no ghosts.
nb_ghosts_right (sequence of int, optional) – Ghost cells on high-index side. Default is no ghosts.
nb_sub_pts (dict, optional) – Number of sub-points per pixel for each sub-point type.
device (str or Device, optional) – Device for field allocation. Can be a string (
"cpu","cuda","cuda:N","rocm","rocm:N") or aDeviceinstance. Default is"cpu".
Example:
>>> from muGrid import Communicator, CartesianDecomposition >>> comm = Communicator() >>> decomp = CartesianDecomposition( ... comm, ... nb_domain_grid_pts=[128, 128], ... nb_subdivisions=[1, 1], ... nb_ghosts_left=[1, 1], ... nb_ghosts_right=[1, 1] ... ) >>> field = decomp.real_field("displacement", components=(3,))
- property nb_grid_pts¶
Local subdomain grid dimensions (alias for
nb_subdomain_grid_pts).
- set_nb_sub_pts(sub_pt_type, nb_sub_pts)¶
Set the number of sub-points for a given sub-point type.
- communicate_ghosts(field)¶
Exchange ghost buffer data for a field.
- Parameters:
field (Field) – The field whose ghost buffers should be filled from neighbors.
- reduce_ghosts(field)¶
Accumulate ghost buffer contributions back to the interior domain.
This is the adjoint operation of
communicate_ghostsand is needed for transpose operations (e.g., divergence) with periodic BCs. After the operation, ghost buffers are zeroed.- Parameters:
field (Field) – The field whose ghost buffers should be reduced to interior.
Communicator¶
- muGrid.Communicator(communicator=None)¶
Factory function for the communicator class.
- Parameters:
communicator – The bare MPI communicator (mpi4py or muGrid communicator object). Default is a serial communicator containing just the present process.
- Returns:
A muGrid communicator object.
Example:
>>> from muGrid import Communicator >>> comm = Communicator() # Serial communicator >>> # Or with MPI: >>> from mpi4py import MPI >>> comm = Communicator(MPI.COMM_WORLD)
Operators¶
Discrete operators for stencil-based computations.
ConvolutionOperator¶
- class muGrid.GenericLinearOperator(offset, stencil)¶
Applies convolution (stencil) operations to fields. Useful for computing gradients, Laplacians, and other discrete differential operators.
- Parameters:
offset (sequence of int) – Offset of the stencil origin relative to the current pixel.
stencil (array_like) – Stencil coefficients. Shape determines the stencil size.
Example:
>>> # Create a 2D Laplacian stencil >>> import numpy as np >>> stencil = np.array([[0, 1, 0], [1, -4, 1], [0, 1, 0]]) >>> laplace = muGrid.GenericLinearOperator([-1, -1], stencil) >>> laplace.apply(input_field, output_field)
- apply(nodal_field, quadrature_point_field)¶
Apply convolution to fields.
LaplaceOperator¶
- class muGrid.LaplaceOperator(spatial_dim, scale=1.0)¶
A hard-coded, optimized Laplacian stencil operator using the standard 5-point (2D) or 7-point (3D) finite difference stencil.
- Parameters:
Example:
>>> laplace = muGrid.LaplaceOperator(2, scale=-1.0) >>> laplace.apply(input_field, output_field)
- apply(input_field, output_field)¶
Apply the Laplacian operator.
- apply_increment(input_field, alpha, output_field)¶
Apply Laplacian and add scaled result to output:
output += alpha * L(input).
- transpose(input_field, output_field, weights=None)¶
Apply transpose operator. For Laplacian, this is the same as apply.
FEMGradientOperator¶
- class muGrid.FEMGradientOperator(spatial_dim, grid_spacing=None)¶
A hard-coded, optimized gradient operator using linear finite element shape functions on triangles (2D) or tetrahedra (3D).
- Parameters:
Example:
>>> grad = muGrid.FEMGradientOperator(2) >>> grad.apply(nodal_field, quadrature_point_gradient_field)
- apply(nodal_field, quadrature_point_field)¶
Apply the gradient operator (nodal values → quadrature point gradients).
- apply_increment(nodal_field, alpha, quadrature_point_field)¶
Apply gradient and add scaled result to output:
output += alpha * grad(input).
- transpose(quadrature_point_field, nodal_field, weights=None)¶
Apply transpose (divergence) operator (quadrature points → nodal values).
- transpose_increment(quadrature_point_field, alpha, nodal_field, weights=None)¶
Apply transpose and add scaled result:
output += alpha * div(input).
IsotropicStiffnessOperator¶
- class muGrid.IsotropicStiffnessOperator2D(grid_spacing)¶
- class muGrid.IsotropicStiffnessOperator3D(grid_spacing)¶
Fused stiffness operators for isotropic linear elastic materials.
These operators compute
K @ u = B^T C B @ ufor linear finite elements without storing the full stiffness matrix. Instead, they exploit the isotropic material structure:\[K = 2\mu G + \lambda V\]where G and V are geometry-only matrices (same for all voxels) and λ, μ are the Lamé parameters which can vary spatially.
This reduces memory from O(N × 24²) for full K storage to O(N × 2) for spatially-varying isotropic materials, plus O(1) for the shared geometry matrices.
- Parameters:
grid_spacing (sequence of float) – Grid spacing in each direction
[hx, hy](2D) or[hx, hy, hz](3D).
Example:
>>> import numpy as np >>> import muGrid >>> >>> # Create 3D operator >>> stiffness = muGrid.IsotropicStiffnessOperator3D([0.1, 0.1, 0.1]) >>> >>> # Setup fields (displacement at nodes, materials at elements) >>> fc = muGrid.CartesianDecomposition( ... muGrid.Communicator(), ... nb_domain_grid_pts=[32, 32, 32], ... nb_ghosts_left=[1, 1, 1], ... nb_ghosts_right=[1, 1, 1], ... ) >>> u = fc.real_field("displacement", (3,)) >>> f = fc.real_field("force", (3,)) >>> >>> # Material parameters (one fewer element in each direction) >>> lam = fc.real_field("lambda") # Lamé first parameter >>> mu = fc.real_field("mu") # Shear modulus >>> lam.p[...] = 1.0 # Uniform material >>> mu.p[...] = 1.0 >>> >>> # Apply stiffness: f = K @ u >>> fc.communicate_ghosts(u) >>> stiffness.apply(u, lam, mu, f)
- apply(displacement, lam, mu, force)¶
Apply the stiffness operator:
force = K @ displacement.
- apply_increment(displacement, lam, mu, alpha, force)¶
Apply stiffness and add scaled result:
force += alpha * K @ displacement.
FFT Engine¶
Distributed FFT operations using pencil decomposition.
FFTEngine¶
- class muGrid.FFTEngine(nb_domain_grid_pts, communicator=None, nb_ghosts_left=None, nb_ghosts_right=None, nb_sub_pts=None)¶
The FFTEngine provides distributed FFT operations on structured grids with MPI parallelization using pencil (2D) decomposition.
- Parameters:
nb_domain_grid_pts (sequence of int) – Global grid dimensions
[Nx, Ny]or[Nx, Ny, Nz].communicator (Communicator, optional) – MPI communicator. Default is serial execution.
nb_ghosts_left (sequence of int, optional) – Ghost cells on low-index side of each dimension.
nb_ghosts_right (sequence of int, optional) – Ghost cells on high-index side of each dimension.
nb_sub_pts (dict, optional) – Number of sub-points per pixel.
Example:
>>> engine = muGrid.FFTEngine([64, 64]) >>> real_field = engine.real_space_field("displacement", nb_components=3) >>> fourier_field = engine.fourier_space_field("displacement_k", nb_components=3) >>> engine.fft(real_field, fourier_field) >>> engine.ifft(fourier_field, real_field) >>> real_field.s[:] *= engine.normalisation
- fft(input_field, output_field)¶
Forward FFT: real space → Fourier space.
The transform is unnormalized. To recover original data after
ifft(fft(x)), multiply bynormalisation.
- ifft(input_field, output_field)¶
Inverse FFT: Fourier space → real space.
The transform is unnormalized. To recover original data after
ifft(fft(x)), multiply bynormalisation.
- register_real_space_field(name, nb_components=1)¶
Register a new real-space field.
Raises an error if a field with the given name already exists.
- Parameters:
- Returns:
Wrapped real-valued field with array accessors.
- Return type:
- Raises:
RuntimeError – If a field with the given name already exists.
- register_fourier_space_field(name, nb_components=1)¶
Register a new Fourier-space field.
Raises an error if a field with the given name already exists.
- Parameters:
- Returns:
Wrapped complex-valued field with array accessors.
- Return type:
- Raises:
RuntimeError – If a field with the given name already exists.
- real_space_field(name, nb_components=1)¶
Get or create a real-space field for FFT operations.
If a field with the given name already exists, returns it. Otherwise creates a new field with the specified number of components.
- fourier_space_field(name, nb_components=1)¶
Get or create a Fourier-space field for FFT operations.
If a field with the given name already exists, returns it. Otherwise creates a new field with the specified number of components.
- normalisation¶
Normalization factor for FFT. Multiply by this after
ifft(fft(x))to recover the original data.- Type:
- fftfreq¶
Normalized FFT frequencies for the local Fourier subdomain.
Returns an array of shape
[dim, local_fx, local_fy, ...]where each element is the normalized frequency in range[-0.5, 0.5).For MPI-parallel runs, this returns only the frequencies for the local subdomain owned by this rank.
- Type:
- ifftfreq¶
Integer FFT frequency indices for the local Fourier subdomain.
Returns an array of shape
[dim, local_fx, local_fy, ...]where each element is the integer frequency index.- Type:
- coords¶
Normalized real-space coordinates for the local subdomain (excluding ghost cells).
Returns an array of shape
[dim, local_nx, local_ny, ...]where each element is the normalized coordinate in range[0, 1).- Type:
- icoords¶
Integer real-space coordinate indices for the local subdomain (excluding ghost cells).
Returns an array of shape
[dim, local_nx, local_ny, ...]where each element is the integer coordinate index.- Type:
- coordsg¶
Normalized real-space coordinates including ghost cells.
Same as
coordsbut includes ghost cells if configured.- Type:
- icoordsg¶
Integer real-space coordinate indices including ghost cells.
Same as
icoordsbut includes ghost cells if configured.- Type:
FFT Utilities¶
- muGrid.fft_normalization(nb_grid_pts)¶
Return FFT normalization factor.
File I/O¶
Classes for reading and writing fields to NetCDF files.
FileIONetCDF¶
- class muGrid.FileIONetCDF(file_name, open_mode='read', communicator=None)¶
Provides NetCDF file I/O for muGrid fields with optional MPI support.
- Parameters:
Example:
>>> file = muGrid.FileIONetCDF("output.nc", open_mode="overwrite") >>> file.register_field_collection(field_collection) >>> file.append_frame().write()
- register_field_collection(collection)¶
Register a field collection for I/O.
- Parameters:
collection (GlobalFieldCollection, LocalFieldCollection, or CartesianDecomposition) – The field collection to register. If a
CartesianDecompositionis passed, its underlying field collection is used.
- append_frame()¶
Append a new frame to the file.
- Returns:
Frame object for writing.
- OpenMode¶
Enum for file open modes:
Read,Write,Overwrite,Append.
Utilities¶
Timer¶
- class muGrid.Timer¶
Hierarchical timer utility for performance measurement with context manager support.
The Timer supports nested timing regions and can be used as a context manager for convenient timing of code blocks.
Example:
>>> timer = muGrid.Timer() >>> with timer("outer"): ... with timer("inner"): ... # ... some computation ... ... pass >>> timer.print_summary() Timer Summary: outer: 0.1234s (1 calls) inner: 0.0567s (1 calls)
- __call__(name)¶
Context manager for timing a named code block.
- Parameters:
name (str) – Name of the timing region.
- Returns:
Context manager that times the enclosed block.
- get_time(name)¶
Get total elapsed time for a timing region.
- get_calls(name)¶
Get the number of times a timing region was entered.
- print_summary()¶
Print a formatted summary of all timing regions.
Solvers¶
The Solvers module provides simple parallel iterative solvers.
- muGrid.Solvers.conjugate_gradients(comm, fc, b, x, hessp, tol=1e-6, maxiter=1000, callback=None, timer=None)¶
Conjugate gradient method for matrix-free solution of the linear problem
Ax = b, whereAis represented by the functionhessp(which computes the product ofAwith a vector). The method iteratively refines the solutionxuntil the residual||Ax - b||is less thantolor untilmaxiteriterations are reached.- Parameters:
comm (muGrid.Communicator) – Communicator for parallel processing.
fc (muGrid.GlobalFieldCollection, muGrid.LocalFieldCollection, or muGrid.CartesianDecomposition) – Collection for creating temporary fields used by the CG algorithm.
b (muGrid.Field) – Right-hand side vector.
x (muGrid.Field) – Initial guess for the solution (modified in place).
hessp (callable) – Function that computes the product of the Hessian matrix
Awith a vector. Signature:hessp(input_field, output_field)where both aremuGrid.Field.tol (float, optional) – Tolerance for convergence. Default is
1e-6.maxiter (int, optional) – Maximum number of iterations. Default is
1000.callback (callable, optional) – Function to call after each iteration with signature:
callback(iteration, state_dict)wherestate_dictcontains keys"x","r","p", and"rr"(squared residual norm).timer (muGrid.Timer, optional) – Timer object for performance profiling. If provided, the solver will record timing for various operations (hessp, dot products, updates).
- Returns:
Solution to the system
Ax = b(same as input fieldx).- Return type:
- Raises:
RuntimeError – If the algorithm does not converge within
maxiteriterations, or if the Hessian is not positive definite.
Device Selection¶
The Device class and DeviceType enum provide a Pythonic way to specify
where field data should be allocated.
Device¶
- class muGrid.Device¶
Represents a compute device (CPU or GPU) for field allocation.
Factory methods provide convenient construction:
- static cuda(device_id=0)¶
Create a CUDA GPU device.
- static rocm(device_id=0)¶
Create a ROCm GPU device.
- static gpu(device_id=0)¶
Create a GPU device using the default backend.
Automatically selects the available GPU backend:
Returns CUDA device if CUDA is available
Returns ROCm device if ROCm is available (and CUDA is not)
Returns CPU device if no GPU backend is available
This is the recommended way to request GPU execution without hard-coding a specific backend.
- is_host()¶
Check if this is a host (CPU) device.
- Returns:
Trueif CPU device,Falseotherwise.- Return type:
- is_device()¶
Check if this is a device (GPU) memory location.
- Returns:
Trueif GPU device,Falseotherwise.- Return type:
- get_type()¶
Get the device type.
- Returns:
The device type enum value.
- Return type:
- get_device_id()¶
Get the device ID for multi-GPU systems.
- Returns:
Device ID (0 for single-GPU or CPU).
- Return type:
Example:
>>> import muGrid >>> # Create devices >>> cpu = muGrid.Device.cpu() >>> gpu0 = muGrid.Device.cuda() >>> gpu1 = muGrid.Device.cuda(1) >>> # Check device type >>> cpu.is_host() True >>> gpu0.is_device() True >>> # Use with field collections >>> fc = muGrid.GlobalFieldCollection([64, 64], device=muGrid.Device.cuda())
DeviceType¶
Enumerations¶
Module Constants¶
The following constants indicate compile-time configuration:
- muGrid.has_mpi¶
Trueif MPI support is enabled.
- muGrid.has_cuda¶
Trueif CUDA GPU support is compiled in.
- muGrid.has_rocm¶
Trueif ROCm/HIP GPU support is compiled in.
- muGrid.has_gpu¶
Trueif any GPU support is available.
- muGrid.has_netcdf¶
Trueif NetCDF I/O support is available.