Skip to content

Modules

Grids

Chmy.Grids.AbstractAxis Type
julia
abstract type AbstractAxis{T}

Abstract type representing an axis in a grid, where the axis is parameterized by the type T of the coordinates.

source
Chmy.Grids.Center Type
julia
struct Center <: Location

The Center struct represents a location at the center along a dimension of a grid cell.

source
Chmy.Grids.Connectivity Type
julia
abstract type Connectivity

Abstract type representing the connectivity of grid elements.

source
Chmy.Grids.Location Type
julia
abstract type Location

Abstract type representing a location in a grid cell.

source
Chmy.Grids.StructuredGrid Type
julia
StructuredGrid

Represents a structured grid with orthogonal axes.

source
Chmy.Grids.UniformGrid Method
julia
UniformGrid(arch; origin, extent, dims, topology=nothing)

Constructs a uniform grid with specified origin, extent, dimensions, and topology.

Arguments

  • arch::Architecture: The associated architecture.

  • origin::NTuple{N,Number}: The origin of the grid.

  • extent::NTuple{N,Number}: The extent of the grid.

  • dims::NTuple{N,Integer}: The dimensions of the grid.

  • topology=nothing: The topology of the grid. If not provided, a default Bounded topology is used.

source
Chmy.Grids.Vertex Type
julia
struct Vertex <: Location

The Vertex struct represents a location at the vertex along a dimension of a grid cell.

source
Chmy.Grids.axes_names Method
julia
axes_names(::SG{1})

Returns the names of the axes for a 1-dimensional structured grid.

source
Chmy.Grids.axes_names Method
julia
axes_names(::SG{2})

Returns the names of the axes for a 2-dimensional structured grid.

source
Chmy.Grids.axes_names Method
julia
axes_names(::SG{3})

Returns the names of the axes for a 3-dimensional structured grid.

source
Chmy.Grids.axis Method
julia
axis(grid, dim::Dim)

Return the axis corresponding to the spatial dimension dim.

source
Chmy.Grids.bounds Method
julia
bounds(grid, loc, [dim::Dim])

Return the bounds of a structured grid at the specified location(s).

source
Chmy.Grids.connectivity Method
julia
connectivity(grid, dim::Dim, side::Side)

Return the connectivity of the structured grid grid for the given dimension dim and side side.

source
Chmy.Grids.coord Method
julia
coord(grid, loc, I...)

Return a tuple of spatial coordinates of a grid point at location loc and indices I.

For vertex locations, first grid point is at the origin. For center locations, first grid point at half-spacing distance from the origin.

source
Chmy.Grids.extent Method
julia
extent(grid, loc, [dim::Dim])

Return the extent of a structured grid at the specified location(s).

source
Chmy.Grids.inv_spacing Method
julia
inv_spacing(grid, loc, I...)

Return a tuple of inverse grid spacings at location loc and indices I.

source
Chmy.Grids.inv_spacing Method
julia
inv_spacing(grid::UniformGrid)

Return a tuple of inverse grid spacing for a uniform grid grid.

source
Chmy.Grids.inv_volume Method
julia
inv_volume(grid, loc, I...)

Return the inverse of control volume at location loc and indices I.

source
Chmy.Grids.iΔ Function
julia

Alias for the inv_spacing method that returns the reciprocal of the spacing between grid points.

source
Chmy.Grids.origin Method
julia
origin(grid, loc, [dim::Dim])

Return the origin of a structured grid at the specified location(s).

source
Chmy.Grids.spacing Method
julia
spacing(grid, loc, I...)

Return a tuple of grid spacings at location loc and indices I.

source
Chmy.Grids.spacing Method
julia
spacing(grid::UniformGrid)

Return a tuple of grid spacing for a uniform grid grid.

source
Chmy.Grids.volume Method
julia
volume(grid, loc, I...)

Return the control volume at location loc and indices I.

source
Chmy.Grids.Δ Function
julia
Δ

Alias for the spacing method that returns the spacing between grid points.

source

Architectures

Chmy.Architectures.Architecture Type
julia
Architecture

Abstract type representing an architecture.

source
Chmy.Architectures.SingleDeviceArchitecture Type
julia
SingleDeviceArchitecture <: Architecture

A struct representing an architecture that operates on a single CPU or GPU device.

source
Chmy.Architectures.SingleDeviceArchitecture Method
julia
SingleDeviceArchitecture(arch::Architecture)

Create a SingleDeviceArchitecture object retrieving backend and device from arch.

source
Chmy.Architectures.Arch Method
julia
Arch(backend::Backend; device_id::Integer=1)

Create an architecture object for the specified backend and device.

Arguments

  • backend: The backend to use for computation.

  • device_id=1: The ID of the device to use.

source
Chmy.Architectures.activate! Method
julia
activate!(arch::SingleDeviceArchitecture; priority=:normal)

Activate the given architecture on the specified device and set the priority of the backend. For the priority accepted values are :normal, :low and :high.

source
Chmy.Architectures.get_device Method
julia
get_device(arch::SingleDeviceArchitecture)

Get the device associated with a SingleDeviceArchitecture.

source
KernelAbstractions.get_backend Method
julia
get_backend(arch::SingleDeviceArchitecture)

Get the backend associated with a SingleDeviceArchitecture.

source

Fields

Chmy.Fields.AbstractField Type
julia
abstract type AbstractField{T,N,L} <: AbstractArray{T,N}

Abstract type representing a field with data type T, number of dimensions N, location L where the field should be defined on.

See also: abstract type ConstantField

source
Chmy.Fields.ConstantField Type

Scalar field with a constant value

source
Chmy.Fields.Field Type
julia
struct Field{T,N,L,H,A} <: AbstractField{T,N,L}

Field represents a discrete scalar field with specified type, number of dimensions, location, and halo size.

source
Chmy.Fields.Field Method
julia
Field(arch::Architecture, args...; kwargs...)

Create a Field object on the specified architecture.

Arguments:

  • arch::Architecture: The architecture for which to create the Field.

  • args...: Additional positional arguments to pass to the Field constructor.

  • kwargs...: Additional keyword arguments to pass to the Field constructor.

source
Chmy.Fields.Field Method
julia
Field(backend, grid, loc, type=eltype(grid); halo=1)

Constructs a field on a structured grid at the specified location.

Arguments:

  • backend: The backend to use for memory allocation.

  • grid: The structured grid on which the field is constructed.

  • loc: The location or locations on the grid where the field is constructed.

  • type: The element type of the field. Defaults to the element type of the grid.

  • halo: The halo size for the field. Defaults to 1.

source
Chmy.Fields.FunctionField Type
julia
FunctionField <: AbstractField

Continuous or discrete field with values computed at runtime.

Constructors

  • FunctionField(func, grid, loc; [discrete], [parameters]): Create a new FunctionField object.
source
Chmy.Fields.FunctionField Method
julia
FunctionField(func::F, grid::StructuredGrid{N}, loc; discrete=false, parameters=nothing) where {F,N}

Create a FunctionField on the given grid using the specified function func.

Arguments:

  • func::F: The function used to generate the field values.

  • grid::StructuredGrid{N}: The structured grid defining the computational domain.

  • loc: The nodal location on the grid grid where the function field is defined on.

  • discrete::Bool=false: A flag indicating whether the field should be discrete. Defaults to false.

  • parameters=nothing: Additional parameters to be used by the function. Defaults to nothing.

source
Chmy.Fields.OneField Type

Constant field with values equal to one(T)

source
Chmy.Fields.ValueField Type

Field with a constant value

source
Chmy.Fields.ZeroField Type

Constant field with values equal to zero(T)

source
Chmy.Fields.TensorField Method
julia
TensorField(backend::Backend, grid::StructuredGrid{2}, args...; kwargs...)

Create a 2D tensor field in the form of a named tuple on the given grid using the specified backend, with components xx, yy, and xy each being a Field.

Arguments:

  • backend::Backend: The backend to be used for computation.

  • grid::StructuredGrid{2}: The 2D structured grid defining the computational domain.

  • args...: Additional positional arguments to pass to the Field constructor.

  • kwargs...: Additional keyword arguments to pass to the Field constructor.

source
Chmy.Fields.TensorField Method
julia
TensorField(backend::Backend, grid::StructuredGrid{3}, args...; kwargs...)

Create a 3D tensor field in the form of a named tuple on the given grid using the specified backend, with components xx, yy, zz, xy, xz, and yz each being a Field.

Arguments:

  • backend::Backend: The backend to be used for computation.

  • grid::StructuredGrid{3}: The 3D structured grid defining the computational domain.

  • args...: Additional positional arguments to pass to the Field constructor.

  • kwargs...: Additional keyword arguments to pass to the Field constructor.

source
Chmy.Fields.VectorField Method
julia
VectorField(backend::Backend, grid::StructuredGrid{N}, args...; kwargs...) where {N}

Create a vector field in the form of a NamedTuple on the given grid using the specified backend. With each component being a Field.

Arguments:

  • backend::Backend: The backend to be used for computation.

  • grid::StructuredGrid{N}: The structured grid defining the computational domain.

  • args...: Additional positional arguments to pass to the Field constructor.

  • kwargs...: Additional keyword arguments to pass to the Field constructor.

source
Chmy.Fields.interior Method
julia
interior(f::Field; with_halo=false)

Displays the field on the interior of the grid on which it is defined on. One could optionally specify to display the halo regions on the grid with with_halo=true.

source
Chmy.Fields.set! Method
julia
set!(f::Field, A::AbstractArray)

Set the elements of the Field f using the values from the AbstractArray A.

Arguments:

  • f::Field: The Field object to be modified.

  • A::AbstractArray: The array whose values are to be copied to the Field.

source
Chmy.Fields.set! Method
julia
set!(f::Field, other::AbstractField)

Set the elements of the Field f using the values from another AbstractField other.

Arguments:

  • f::Field: The destination Field object to be modified.

  • other::AbstractField: The source AbstractField whose values are to be copied to f.

source
Chmy.Fields.set! Method
julia
set!(f::Field, val::Number)

Set all elements of the Field f to the specified numeric value val.

Arguments:

  • f::Field: The Field object to be modified.

  • val::Number: The numeric value to set in the Field.

source

Grid Operators

Chmy.GridOperators.AbstractMask Type
julia
abstract type AbstractMask{T,N}

Abstract type representing the data transformation to be performed on elements in a field of dimension N, where each element is of typeT.

source
Chmy.GridOperators.InterpolationRule Type
julia
abstract type InterpolationRule

A type representing an interpolation rule that specifies how the interpolant f should be reconstructed using a data set on a given grid.

source
Chmy.GridOperators.hlerp Method
julia
hlerp(f, to, grid, I...)

Interpolate a field f to location to using harmonic linear interpolation rule.

rule(t, v0, v1) = 1/(1/v0 + t * (1/v1 - 1/v0))

source
Chmy.GridOperators.itp Method
julia
itp(f, to, r, grid, I...)

Interpolates the field f from its current location to the specified location(s) to using the given interpolation rule r. The indices specify the position within the grid at location(s) to.

source
Chmy.GridOperators.leftx Method
julia
leftx(f, ω, I)

"left side" of a field ([1:end-1]) in x direction, masked with ω.

source
Chmy.GridOperators.leftx Method
julia
leftx(f, I)

"left side" of a field ([1:end-1]) in x direction.

source
Chmy.GridOperators.lefty Method
julia
lefty(f, ω, I)

"left side" of a field ([1:end-1]) in y direction, masked with ω.

source
Chmy.GridOperators.lefty Method
julia
lefty(f, I)

"left side" of a field ([1:end-1]) in y direction.

source
Chmy.GridOperators.leftz Method
julia
leftz(f, ω, I)

"left side" of a field ([1:end-1]) in z direction, masked with ω.

source
Chmy.GridOperators.leftz Method
julia
leftz(f, I)

"left side" of a field ([1:end-1]) in z direction.

source
Chmy.GridOperators.lerp Method
julia
lerp(f, to, grid, I...)

Linearly interpolate values of a field f to location to.

source
Chmy.GridOperators.rightx Method
julia
rightx(f, ω, I)

"right side" of a field ([2:end]) in x direction, masked with ω.

source
Chmy.GridOperators.rightx Method
julia
rightx(f, I)

"right side" of a field ([2:end]) in x direction.

source
Chmy.GridOperators.righty Method
julia
righty(f, ω, I)

"right side" of a field ([2:end]) in y direction, masked with ω.

source
Chmy.GridOperators.righty Method
julia
righty(f, I)

"right side" of a field ([2:end]) in y direction.

source
Chmy.GridOperators.rightz Method
julia
rightz(f, ω, I)

"right side" of a field ([2:end]) in z direction, masked with ω.

source
Chmy.GridOperators.rightz Method
julia
rightz(f, I)

"right side" of a field ([2:end]) in z direction.

source
Chmy.GridOperators.δx Method
julia
δx(f, ω, I)

Finite difference in x direction, masked with ω.

source
Chmy.GridOperators.δx Method
julia
δx(f, I)

Finite difference in x direction.

source
Chmy.GridOperators.δy Method
julia
δy(f, ω, I)

Finite difference in y direction, masked with ω.

source
Chmy.GridOperators.δy Method
julia
δy(f, I)

Finite difference in y direction.

source
Chmy.GridOperators.δz Method
julia
δz(f, ω, I)

Finite difference in z direction, masked with ω.

source
Chmy.GridOperators.δz Method
julia
δz(f, I)

Finite difference in z direction.

source
Chmy.GridOperators.∂x Method
julia
∂x(f, ω, grid, I)

Directional partial derivative in x direction, masked with ω.

source
Chmy.GridOperators.∂x Method
julia
∂x(f, grid, I)

Directional partial derivative in x direction.

source
Chmy.GridOperators.∂y Method
julia
∂y(f, ω, grid, I)

Directional partial derivative in y direction, masked with ω.

source
Chmy.GridOperators.∂y Method
julia
∂y(f, grid, I)

Directional partial derivative in y direction.

source
Chmy.GridOperators.∂z Method
julia
∂z(f, ω, grid, I)

Directional partial derivative in z direction, masked with ω.

source
Chmy.GridOperators.∂z Method
julia
∂z(f, grid, I)

Directional partial derivative in z direction.

source

Boundary Conditions

Chmy.BoundaryConditions.AbstractBatch Type
julia
AbstractBatch

Abstract type representing a batch of boundary conditions.

source
Chmy.BoundaryConditions.BoundaryFunction Type
julia
abstract type BoundaryFunction{F}

Abstract type for boundary condition functions with function type F.

source
Chmy.BoundaryConditions.Dirichlet Type
julia
Dirichlet(value=nothing)

Create a Dirichlet object representing the Dirichlet boundary condition with the specified value.

source
Chmy.BoundaryConditions.EmptyBatch Type
julia
EmptyBatch <: AbstractBatch

EmptyBatch represents no boundary conditions.

source
Chmy.BoundaryConditions.ExchangeBatch Type
julia
ExchangeBatch <: AbstractBatch

ExchangeBatch represents a batch used for MPI communication.

source
Chmy.BoundaryConditions.FieldBatch Type
julia
FieldBatch <: AbstractBatch

FieldBatch is a batch of boundary conditions, where each field has one boundary condition.

source
Chmy.BoundaryConditions.FieldBoundaryCondition Type
julia
FieldBoundaryCondition

Abstract supertype for all boundary conditions that are specified per-field.

source
Chmy.BoundaryConditions.FirstOrderBC Type
julia
struct FirstOrderBC{T,Kind} <: FieldBoundaryCondition

A struct representing a boundary condition of first-order accuracy.

source
Chmy.BoundaryConditions.Neumann Type
julia
Neumann(value=nothing)

Create a Neumann object representing the Neumann boundary condition with the specified value.

source
Chmy.BoundaryConditions.bc! Method
julia
bc!(arch::Architecture, grid::StructuredGrid, batch::BatchSet)

Apply boundary conditions using a batch set batch containing an AbstractBatch per dimension per side of grid.

Arguments

  • arch: The architecture.

  • grid: The grid.

  • batch:: The batch set to apply boundary conditions to.

source

Kernel launcher

Chmy.KernelLaunch.Launcher Type
julia
struct Launcher{Worksize,OuterWidth,Workers}

A struct representing a launcher for asynchronous kernel execution.

source
Chmy.KernelLaunch.Launcher Method
julia
Launcher(arch, grid; outer_width=nothing)

Constructs a Launcher object configured based on the input parameters.

Arguments:

  • arch: The associated architecture.

  • grid: The grid defining the computational domain.

  • outer_width: Optional parameter specifying outer width.

Warning

worksize for the last dimension N takes into account only last outer width W[N], N-1 uses W[N] and W[N-1], N-2 uses W[N], W[N-1], and W[N-2].

source
Chmy.KernelLaunch.Launcher Method
julia
(launcher::Launcher)(arch::Architecture, grid, kernel_and_args::Pair{F,Args}; bc=nothing, async=false) where {F,Args}

Launches a computational kernel using the specified arch, grid, kernel_and_args, and optional boundary conditions (bc).

Arguments:

  • arch::Architecture: The architecture on which to execute the computation.

  • grid: The grid defining the computational domain.

  • kernel_and_args::Pair{F,Args}: A pair consisting of the computational kernel F and its arguments Args.

  • bc=nothing: Optional boundary conditions for the computation.

  • async=false: If true, launches the kernel asynchronously.

Warning

  • arch should be compatible with the Launcher's architecture.

  • If bc is nothing, the kernel is launched without boundary conditions.

  • If async is false (default), the function waits for the computation to complete before returning.

source

Distributed

Chmy.Distributed.CartesianTopology Type
julia
CartesianTopology

Represents N-dimensional Cartesian topology of distributed processes.

source
Chmy.Distributed.CartesianTopology Method
julia
CartesianTopology(comm, dims)

Create an N-dimensional Cartesian topology using base MPI communicator comm with dimensions dims. If all entries in dims are not equal to 0, the product of dims should be equal to the total number of MPI processes MPI.Comm_size(comm). If any (or all) entries of dims are 0, the dimensions in the corresponding spatial directions will be picked automatically.

source
Chmy.Distributed.DistributedArchitecture Type
julia
DistributedArchitecture <: Architecture

A struct representing a distributed architecture.

source
Chmy.Distributed.StackAllocator Type
julia
mutable struct StackAllocator

Simple stack (a.k.a. bump/arena) allocator. Maintains an internal buffer that grows dynamically if the requested allocation exceeds current buffer size.

source
Chmy.Distributed.StackAllocator Method
julia
StackAllocator(backend::Backend)

Create a stack allocator using the specified backend to store allocations.

source
Base.resize! Method
julia
resize!(sa::StackAllocator, sz::Integer)

Resize the StackAllocator's buffer to capacity of sz bytes. This method will throw an error if any arrays were already allocated using this allocator.

source
Chmy.Architectures.Arch Method
julia
Architectures.Arch(backend::Backend, comm::MPI.Comm, dims; device_id=nothing)

Create a distributed Architecture using backend backend and comm. For GPU backends, device will be selected automatically based on a process id within a node, unless specified by device_id.

Arguments

  • backend::Backend: The backend to use for the architecture.

  • comm::MPI.Comm: The MPI communicator to use for the architecture.

  • dims: The dimensions of the architecture.

Keyword Arguments

  • device_id: The ID of the device to use. If not provided, the shared rank of the topology plus one is used.
source
Chmy.Architectures.activate! Method
julia
activate!(arch::DistributedArchitecture; kwargs...)

Activate the given DistributedArchitecture by delegating to the child architecture, and pass through any keyword arguments. For example, the priority can be set with accepted values being :normal, :low, and :high.

source
Chmy.Architectures.get_device Method
julia
get_device(arch::DistributedArchitecture)

Get the device associated with a DistributedArchitecture by delegating to the child architecture.

source
Chmy.BoundaryConditions.bc! Method
julia
BoundaryConditions.bc!(side::Side, dim::Dim,
                            arch::DistributedArchitecture,
                            grid::StructuredGrid,
                            batch::ExchangeBatch; kwargs...)

Apply boundary conditions on a distributed grid with halo exchange performed internally.

Arguments

  • side: The side of the grid where the halo exchange is performed.

  • dim: The dimension along which the halo exchange is performed.

  • arch: The distributed architecture used for communication.

  • grid: The structured grid on which the halo exchange is performed.

  • batch: The batch set to apply boundary conditions to.

source
Chmy.Distributed.allocate Function
julia
allocate(sa::StackAllocator, T::DataType, dims, [align=sizeof(T)])

Allocate a buffer of type T with dimensions dims using a stack allocator. The align parameter specifies the alignment of the buffer elements.

Arguments

  • sa::StackAllocator: The stack allocator object.

  • T::DataType: The data type of the requested allocation.

  • dims: The dimensions of the requested allocation.

  • align::Integer: The alignment of the allocated buffer in bytes.

Warning

Arrays allocated with StackAllocator are not managed by Julia runtime. User is responsible for ensuring correct lifetimes, i.e., that the reference to allocator outlives all arrays allocated using this allocator.

source
Chmy.Distributed.cart_comm Method
julia
cart_comm(topo)

MPI Cartesian communicator for the topology.

source
Chmy.Distributed.cart_coords Method
julia
cart_coords(topo)

Coordinates of a current process within a Cartesian topology.

source
Chmy.Distributed.dims Method
julia
dims(topo)

Dimensions of the topology as NTuple.

source
Chmy.Distributed.exchange_halo! Method
julia
exchange_halo!(side::Side, dim::Dim, arch, grid, fields...; async=false)

Perform halo exchange communication between neighboring processes in a distributed architecture.

Arguments

  • side: The side of the grid where the halo exchange is performed.

  • dim: The dimension along which the halo exchange is performed.

  • arch: The distributed architecture used for communication.

  • grid: The structured grid on which the halo exchange is performed.

  • fields...: The fields to be exchanged.

Optional Arguments

  • async=false: Whether to perform the halo exchange asynchronously.
source
Chmy.Distributed.exchange_halo! Method
julia
exchange_halo!(arch, grid, fields...)

Perform halo exchange for the given architecture, grid, and fields.

Arguments

  • arch: The distributed architecture to perform halo exchange on.

  • grid: The structured grid on which halo exchange is performed.

  • fields: The fields on which halo exchange is performed.

source
Chmy.Distributed.gather! Method
julia
gather!(arch, dst, src::Field; kwargs...)

Gather the interior of a field src into a global array dst on the CPU.

source
Chmy.Distributed.gather! Method
julia
gather!(dst, src, comm::MPI.Comm; root=0)

Gather local array src into a global array dst. Size of the global array size(dst) should be equal to the product of the size of a local array size(src) and the dimensions of a Cartesian communicator comm. The array will be gathered on the process with id root (root=0 by default). Note that the memory for a global array should be allocated only on the process with id root, on other processes dst can be set to nothing.

source
Chmy.Distributed.global_rank Method
julia
global_rank(topo)

Global id of a process in a Cartesian topology.

source
Chmy.Distributed.global_size Method
julia
global_size(topo)

Total number of processes within the topology.

source
Chmy.Distributed.has_neighbor Method
julia
has_neighbor(topo, dim, side)

Returns true if there a neighbor process in spatial direction dim on the side side, or false otherwise.

source
Chmy.Distributed.nallocs Method
julia
nallocs(sa::StackAllocator)

Get the number of allocations made by the given StackAllocator.

source
Chmy.Distributed.neighbor Method
julia
neighbor(topo, dim, side)

Returns id of a neighbor process in spatial direction dim on the side side, if this neighbor exists, or MPI.PROC_NULL otherwise.

source
Chmy.Distributed.neighbors Method
julia
neighbors(topo)

Neighbors of a current process.

Returns tuple of ranks of the two immediate neighbors in each spatial direction, or MPI.PROC_NULL if there is no neighbor on a corresponding side.

source
Chmy.Distributed.node_name Method
julia
node_name(topo)

Name of a node according to MPI.Get_processor_name().

source
Chmy.Distributed.node_size Method
julia
node_size(topo)

Number of processes sharing the same node.

source
Chmy.Distributed.reset! Method
julia
reset!(sa::StackAllocator)

Reset the stack allocator by resetting the pointer. Doesn't free the internal memory buffer.

source
Chmy.Distributed.shared_comm Method
julia
shared_comm(topo)

MPI communicator for the processes sharing the same node.

source
Chmy.Distributed.shared_rank Method
julia
shared_rank(topo)

Local id of a process within a single node. Can be used to set the GPU device.

source
Chmy.Distributed.topology Method
julia
topology(arch::DistributedArchitecture)

Get the virtual MPI topology of a distributed architecture

source
KernelAbstractions.get_backend Method
julia
get_backend(arch::DistributedArchitecture)

Get the backend associated with a DistributedArchitecture by delegating to the child architecture.

source

Workers

Chmy.Workers.Worker Type
julia
Worker

A worker that performs tasks asynchronously.

Constructor

Worker{T}(; [setup], [teardown]) where {T}

Constructs a new Worker object.

Arguments

  • setup: A function to be executed before the worker starts processing tasks. (optional)

  • teardown: A function to be executed after the worker finishes processing tasks. (optional)

source