Modules
Grids
Chmy.Grids.AbstractAxis Type
abstract type AbstractAxis{T}
Abstract type representing an axis in a grid, where the axis is parameterized by the type T
of the coordinates.
Chmy.Grids.Center Type
struct Center <: Location
The Center
struct represents a location at the center along a dimension of a grid cell.
Chmy.Grids.Connectivity Type
abstract type Connectivity
Abstract type representing the connectivity of grid elements.
sourceChmy.Grids.Location Type
abstract type Location
Abstract type representing a location in a grid cell.
sourceChmy.Grids.StructuredGrid Type
StructuredGrid
Represents a structured grid with orthogonal axes.
sourceChmy.Grids.UniformGrid Method
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 defaultBounded
topology is used.
Chmy.Grids.Vertex Type
struct Vertex <: Location
The Vertex
struct represents a location at the vertex along a dimension of a grid cell.
Chmy.Grids.axes_names Method
axes_names(::SG{1})
Returns the names of the axes for a 1-dimensional structured grid.
sourceChmy.Grids.axes_names Method
axes_names(::SG{2})
Returns the names of the axes for a 2-dimensional structured grid.
sourceChmy.Grids.axes_names Method
axes_names(::SG{3})
Returns the names of the axes for a 3-dimensional structured grid.
sourceChmy.Grids.axis Method
axis(grid, dim::Dim)
Return the axis corresponding to the spatial dimension dim
.
Chmy.Grids.bounds Method
bounds(grid, loc, [dim::Dim])
Return the bounds of a structured grid at the specified location(s).
sourceChmy.Grids.connectivity Method
connectivity(grid, dim::Dim, side::Side)
Return the connectivity of the structured grid grid
for the given dimension dim
and side side
.
Chmy.Grids.coord Method
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.
sourceChmy.Grids.extent Method
extent(grid, loc, [dim::Dim])
Return the extent of a structured grid at the specified location(s).
sourceChmy.Grids.inv_spacing Method
inv_spacing(grid, loc, I...)
Return a tuple of inverse grid spacings at location loc
and indices I
.
Chmy.Grids.inv_spacing Method
inv_spacing(grid::UniformGrid)
Return a tuple of inverse grid spacing for a uniform grid grid
.
Chmy.Grids.inv_volume Method
inv_volume(grid, loc, I...)
Return the inverse of control volume at location loc
and indices I
.
Chmy.Grids.iΔ Function
iΔ
Alias for the inv_spacing
method that returns the reciprocal of the spacing between grid points.
Chmy.Grids.origin Method
origin(grid, loc, [dim::Dim])
Return the origin of a structured grid at the specified location(s).
sourceChmy.Grids.spacing Method
spacing(grid, loc, I...)
Return a tuple of grid spacings at location loc
and indices I
.
Chmy.Grids.spacing Method
spacing(grid::UniformGrid)
Return a tuple of grid spacing for a uniform grid grid
.
Chmy.Grids.volume Method
volume(grid, loc, I...)
Return the control volume at location loc
and indices I
.
Chmy.Grids.Δ Function
Δ
Alias for the spacing
method that returns the spacing between grid points.
Architectures
Chmy.Architectures.Architecture Type
Architecture
Abstract type representing an architecture.
sourceChmy.Architectures.SingleDeviceArchitecture Type
SingleDeviceArchitecture <: Architecture
A struct representing an architecture that operates on a single CPU or GPU device.
sourceChmy.Architectures.SingleDeviceArchitecture Method
SingleDeviceArchitecture(arch::Architecture)
Create a SingleDeviceArchitecture
object retrieving backend and device from arch
.
Chmy.Architectures.Arch Method
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.
Chmy.Architectures.activate! Method
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
.
Chmy.Architectures.get_device Method
get_device(arch::SingleDeviceArchitecture)
Get the device associated with a SingleDeviceArchitecture.
sourceKernelAbstractions.get_backend Method
get_backend(arch::SingleDeviceArchitecture)
Get the backend associated with a SingleDeviceArchitecture.
sourceFields
Chmy.Fields.AbstractField Type
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
Chmy.Fields.Field Type
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.
sourceChmy.Fields.Field Method
Field(arch::Architecture, args...; kwargs...)
Create a Field
object on the specified architecture.
Arguments:
arch::Architecture
: The architecture for which to create theField
.args...
: Additional positional arguments to pass to theField
constructor.kwargs...
: Additional keyword arguments to pass to theField
constructor.
Chmy.Fields.Field Method
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.
Chmy.Fields.FunctionField Type
FunctionField <: AbstractField
Continuous or discrete field with values computed at runtime.
Constructors
FunctionField(func, grid, loc; [discrete], [parameters])
: Create a newFunctionField
object.
Chmy.Fields.FunctionField Method
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 tofalse
.parameters=nothing
: Additional parameters to be used by the function. Defaults tonothing
.
Chmy.Fields.TensorField Method
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 theField
constructor.kwargs...
: Additional keyword arguments to pass to theField
constructor.
Chmy.Fields.TensorField Method
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 theField
constructor.kwargs...
: Additional keyword arguments to pass to theField
constructor.
Chmy.Fields.VectorField Method
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 theField
constructor.kwargs...
: Additional keyword arguments to pass to theField
constructor.
Chmy.Fields.interior Method
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
.
Chmy.Fields.set! Method
set!(f::Field, A::AbstractArray)
Set the elements of the Field
f
using the values from the AbstractArray
A
.
Arguments:
f::Field
: TheField
object to be modified.A::AbstractArray
: The array whose values are to be copied to theField
.
Chmy.Fields.set! Method
set!(f::Field, other::AbstractField)
Set the elements of the Field
f
using the values from another AbstractField
other
.
Arguments:
f::Field
: The destinationField
object to be modified.other::AbstractField
: The sourceAbstractField
whose values are to be copied tof
.
Chmy.Fields.set! Method
set!(f::Field, val::Number)
Set all elements of the Field
f
to the specified numeric value val
.
Arguments:
f::Field
: TheField
object to be modified.val::Number
: The numeric value to set in theField
.
Grid Operators
Chmy.GridOperators.AbstractMask Type
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
.
Chmy.GridOperators.InterpolationRule Type
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.
Chmy.GridOperators.hlerp Method
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))
Chmy.GridOperators.itp Method
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
.
Chmy.GridOperators.leftx Method
leftx(f, ω, I)
"left side" of a field ([1:end-1]
) in x direction, masked with ω
.
Chmy.GridOperators.leftx Method
leftx(f, I)
"left side" of a field ([1:end-1]
) in x direction.
Chmy.GridOperators.lefty Method
lefty(f, ω, I)
"left side" of a field ([1:end-1]
) in y direction, masked with ω
.
Chmy.GridOperators.lefty Method
lefty(f, I)
"left side" of a field ([1:end-1]
) in y direction.
Chmy.GridOperators.leftz Method
leftz(f, ω, I)
"left side" of a field ([1:end-1]
) in z direction, masked with ω
.
Chmy.GridOperators.leftz Method
leftz(f, I)
"left side" of a field ([1:end-1]
) in z direction.
Chmy.GridOperators.lerp Method
lerp(f, to, grid, I...)
Linearly interpolate values of a field f
to location to
.
Chmy.GridOperators.rightx Method
rightx(f, ω, I)
"right side" of a field ([2:end]
) in x direction, masked with ω
.
Chmy.GridOperators.rightx Method
rightx(f, I)
"right side" of a field ([2:end]
) in x direction.
Chmy.GridOperators.righty Method
righty(f, ω, I)
"right side" of a field ([2:end]
) in y direction, masked with ω
.
Chmy.GridOperators.righty Method
righty(f, I)
"right side" of a field ([2:end]
) in y direction.
Chmy.GridOperators.rightz Method
rightz(f, ω, I)
"right side" of a field ([2:end]
) in z direction, masked with ω
.
Chmy.GridOperators.rightz Method
rightz(f, I)
"right side" of a field ([2:end]
) in z direction.
Chmy.GridOperators.∂x Method
∂x(f, ω, grid, I)
Directional partial derivative in x direction, masked with ω
.
Chmy.GridOperators.∂y Method
∂y(f, ω, grid, I)
Directional partial derivative in y direction, masked with ω
.
Chmy.GridOperators.∂z Method
∂z(f, ω, grid, I)
Directional partial derivative in z direction, masked with ω
.
Boundary Conditions
Chmy.BoundaryConditions.AbstractBatch Type
AbstractBatch
Abstract type representing a batch of boundary conditions.
sourceChmy.BoundaryConditions.BoundaryFunction Type
abstract type BoundaryFunction{F}
Abstract type for boundary condition functions with function type F
.
Chmy.BoundaryConditions.Dirichlet Type
Dirichlet(value=nothing)
Create a Dirichlet
object representing the Dirichlet boundary condition with the specified value.
Chmy.BoundaryConditions.EmptyBatch Type
EmptyBatch <: AbstractBatch
EmptyBatch represents no boundary conditions.
sourceChmy.BoundaryConditions.ExchangeBatch Type
ExchangeBatch <: AbstractBatch
ExchangeBatch represents a batch used for MPI communication.
sourceChmy.BoundaryConditions.FieldBatch Type
FieldBatch <: AbstractBatch
FieldBatch is a batch of boundary conditions, where each field has one boundary condition.
sourceChmy.BoundaryConditions.FieldBoundaryCondition Type
FieldBoundaryCondition
Abstract supertype for all boundary conditions that are specified per-field.
sourceChmy.BoundaryConditions.FirstOrderBC Type
struct FirstOrderBC{T,Kind} <: FieldBoundaryCondition
A struct representing a boundary condition of first-order accuracy.
sourceChmy.BoundaryConditions.Neumann Type
Neumann(value=nothing)
Create a Neumann
object representing the Neumann boundary condition with the specified value.
Chmy.BoundaryConditions.bc! Method
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.
Kernel launcher
Chmy.KernelLaunch.Launcher Type
struct Launcher{Worksize,OuterWidth,Workers}
A struct representing a launcher for asynchronous kernel execution.
sourceChmy.KernelLaunch.Launcher Method
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].
Chmy.KernelLaunch.Launcher Method
(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 kernelF
and its argumentsArgs
.bc=nothing
: Optional boundary conditions for the computation.async=false
: Iftrue
, launches the kernel asynchronously.
Warning
arch
should be compatible with theLauncher
's architecture.If
bc
isnothing
, the kernel is launched without boundary conditions.If
async
isfalse
(default), the function waits for the computation to complete before returning.
Distributed
Chmy.Distributed.CartesianTopology Type
CartesianTopology
Represents N-dimensional Cartesian topology of distributed processes.
sourceChmy.Distributed.CartesianTopology Method
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.
Chmy.Distributed.DistributedArchitecture Type
DistributedArchitecture <: Architecture
A struct representing a distributed architecture.
sourceChmy.Distributed.StackAllocator Type
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.
sourceChmy.Distributed.StackAllocator Method
StackAllocator(backend::Backend)
Create a stack allocator using the specified backend to store allocations.
sourceBase.resize! Method
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.
Chmy.Architectures.Arch Method
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.
Chmy.Architectures.activate! Method
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
.
Chmy.Architectures.get_device Method
get_device(arch::DistributedArchitecture)
Get the device associated with a DistributedArchitecture by delegating to the child architecture.
sourceChmy.BoundaryConditions.bc! Method
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.
Chmy.Distributed.allocate Function
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.
Chmy.Distributed.cart_comm Method
cart_comm(topo)
MPI Cartesian communicator for the topology.
sourceChmy.Distributed.cart_coords Method
cart_coords(topo)
Coordinates of a current process within a Cartesian topology.
sourceChmy.Distributed.exchange_halo! Method
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.
Chmy.Distributed.exchange_halo! Method
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.
Chmy.Distributed.gather! Method
gather!(arch, dst, src::Field; kwargs...)
Gather the interior of a field src
into a global array dst
on the CPU.
Chmy.Distributed.gather! Method
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
.
Chmy.Distributed.global_rank Method
global_rank(topo)
Global id of a process in a Cartesian topology.
sourceChmy.Distributed.global_size Method
global_size(topo)
Total number of processes within the topology.
sourceChmy.Distributed.has_neighbor Method
has_neighbor(topo, dim, side)
Returns true if there a neighbor process in spatial direction dim
on the side side
, or false otherwise.
Chmy.Distributed.nallocs Method
nallocs(sa::StackAllocator)
Get the number of allocations made by the given StackAllocator
.
Chmy.Distributed.neighbor Method
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.
Chmy.Distributed.neighbors Method
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.
sourceChmy.Distributed.node_name Method
node_name(topo)
Name of a node according to MPI.Get_processor_name()
.
Chmy.Distributed.node_size Method
node_size(topo)
Number of processes sharing the same node.
sourceChmy.Distributed.reset! Method
reset!(sa::StackAllocator)
Reset the stack allocator by resetting the pointer. Doesn't free the internal memory buffer.
sourceChmy.Distributed.shared_comm Method
shared_comm(topo)
MPI communicator for the processes sharing the same node.
sourceChmy.Distributed.shared_rank Method
shared_rank(topo)
Local id of a process within a single node. Can be used to set the GPU device.
sourceChmy.Distributed.topology Method
topology(arch::DistributedArchitecture)
Get the virtual MPI topology of a distributed architecture
sourceKernelAbstractions.get_backend Method
get_backend(arch::DistributedArchitecture)
Get the backend associated with a DistributedArchitecture by delegating to the child architecture.
sourceWorkers
Chmy.Workers.Worker Type
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)