API Reference

Domain

Macchiato.DomainType
Domain{M, C}

Central container that ties together a point cloud, boundary conditions, and physics models.

Fields

  • cloud::PointCloud{M, C}: The discretized geometry (boundary + interior points)
  • boundaries::Dict{Symbol, Tuple{UnitRange, AbstractBoundaryCondition}}: Mapping from surface names to (index_range, bc) pairs, where index_range gives the global indices of that surface's points in the assembled system
  • models::AbstractVector{<:AbstractModel}: Physics models (e.g., SolidEnergy, LinearElasticity)
  • name::Symbol: Domain identifier (defaults to :domain1)

Constructors

Domain(cloud, boundaries, model)    # cloud + BCs + model(s)
Domain(cloud, model)                # cloud + model(s), no BCs

At construction, the Domain validates that:

  1. Every boundary condition key matches a surface name in the point cloud
  2. Every surface in the point cloud has a corresponding boundary condition entry
source
Macchiato.add!Function
add!(domain::Domain, model::AbstractModel)

Append a physics model to the domain's model list.

source
add!(domain::Domain, boundary::AbstractBoundaryCondition, name::Symbol)

Attach a boundary condition to the named surface on the domain.

source
Macchiato.delete!Function
delete!(domain::Domain, model::AbstractModel)

Remove a physics model from the domain.

source

Models

Energy

Macchiato.SolidEnergyType
SolidEnergy(; k, ρ, cₚ, source=nothing)

Solid-body energy (heat) transport model.

Solves the heat equation in a solid medium:

  • Steady-state: k ∇²T = -f (Poisson equation)
  • Transient: ρ cₚ ∂T/∂t = k ∇²T + f

Fields

  • k: Thermal conductivity
  • ρ: Density
  • cₚ: Specific heat capacity
  • source: Optional volumetric source term f(x, t) -> value (default: nothing)

Example

model = SolidEnergy(k=50.0, ρ=7800.0, cₚ=500.0)
model = SolidEnergy(k=1.0, ρ=1.0, cₚ=1.0, source=(x, t) -> -4.0)
source

Mechanics

Macchiato.LinearElasticityType
LinearElasticity{E, Nu, Rho, F} <: Solid

Linear isotropic elasticity model for solid mechanics (Navier-Cauchy equations).

Supports 2D plane stress formulation. The governing equations in displacement form:

(λ*+2μ) ∂²u/∂x² + μ ∂²u/∂y² + (λ*+μ) ∂²v/∂x∂y + fₓ = 0
(λ*+μ) ∂²u/∂x∂y + μ ∂²v/∂x² + (λ*+2μ) ∂²v/∂y² + fᵧ = 0

Fields

  • E: Young's modulus
  • ν: Poisson's ratio
  • ρ: Density (optional, for body forces / dynamics)
  • body_force: Body force function f(x) -> (fx, fy) (optional)

Example

model = LinearElasticity(E=200e3, ν=0.3)
model = LinearElasticity(E=200e3, ν=0.3, body_force=x -> (0.0, -9.81))
source
Macchiato.lame_parametersFunction
lame_parameters(model::LinearElasticity)

Compute Lamé parameters for plane stress:

  • μ = E / (2(1+ν))
  • λ = Eν / ((1+ν)(1-2ν))
  • λ* = 2μλ / (λ+2μ) (plane stress modification)

Returns (μ, λstar).

source

Fluids

Macchiato.IncompressibleNavierStokesType
IncompressibleNavierStokes(; μ, ρ)
IncompressibleNavierStokes(μ::Real, ρ)

Incompressible Navier-Stokes fluid model.

When μ is a plain number, it is automatically wrapped in NewtonianViscosity. Pass an AbstractViscosity subtype (e.g., CarreauYasudaViscosity) for non-Newtonian behavior.

Warning

The fluid solver is under active development and not yet fully functional.

Fields

  • μ::AbstractViscosity: Dynamic viscosity model
  • ρ: Fluid density

Example

model = IncompressibleNavierStokes(μ=0.001, ρ=1000.0)
model = IncompressibleNavierStokes(μ=CarreauYasudaViscosity(μ_inf=0.0035, μ_0=0.056), ρ=1060.0)
source
Macchiato.NewtonianViscosityType
NewtonianViscosity(μ)

Constant (Newtonian) viscosity model. Returns the same viscosity μ regardless of shear rate.

Example

visc = NewtonianViscosity(0.001)  # water-like viscosity
visc(100.0)  # => 0.001
source
Macchiato.CarreauYasudaViscosityType
CarreauYasudaViscosity(; μ_inf, μ_0, n=0.333, λ=0.31, a=2)

Generalized Newtonian viscosity model for shear-thinning (or shear-thickening) fluids.

The Carreau-Yasuda model:

μ(γ̇) = μ_inf + (μ_0 - μ_inf) * (1 + (λ γ̇)^a)^((n - 1) / a)

Fields

  • μ_inf: Infinite-shear-rate viscosity
  • μ_0: Zero-shear-rate viscosity
  • n: Power-law index (< 1 for shear-thinning)
  • λ: Relaxation time
  • a: Yasuda parameter (a = 2 recovers the Carreau model)

Example

visc = CarreauYasudaViscosity(μ_inf=0.0035, μ_0=0.056, n=0.333, λ=0.31)
source

Time

Macchiato.TimeType
Time <: AbstractModel

Abstract type for time-mode models that indicate whether a simulation is steady-state or transient. See Steady and Unsteady.

source
Macchiato.SteadyType
Steady(max_time)

Marker for steady-state simulations. max_time may be used as a convergence limit.

source
Macchiato.UnsteadyType
Unsteady(max_time)

Marker for transient (time-dependent) simulations. max_time sets the end time.

source

Model Interface

These are the functions to implement when defining a custom PDE.

Macchiato._num_varsFunction
_num_vars(model::AbstractModel, dim) -> Int

Return the number of solution variables per point for model in dim dimensions.

Examples: 1 for scalar PDEs, dim for vector PDEs, dim + 1 for velocity + pressure.

source
Macchiato.make_systemFunction
make_system(model::AbstractModel, domain; kwargs...) -> (A, b)

Assemble the system matrix A and right-hand side b for steady-state solving.

Required for steady-state simulations. Macchiato applies boundary conditions and solves Ax = b.

source
Macchiato.make_fFunction
make_f(model::AbstractModel, domain; kwargs...) -> f

Return an in-place ODE function f(du, u, p, t) for transient integration.

Required for transient simulations. Macchiato passes the returned function to OrdinaryDiffEq.jl.

source

Boundary Conditions

Core Types

Macchiato.DirichletType
Dirichlet <: AbstractBoundaryCondition

Essential boundary conditions that prescribe values at the boundary. Examples: Temperature, VelocityInlet, Displacement.

source
Macchiato.NeumannType
Neumann <: DerivativeBoundaryCondition

Natural boundary conditions that prescribe normal derivatives: ∂u/∂n = g. Examples: HeatFlux, Adiabatic, VelocityOutlet.

source
Macchiato.RobinType
Robin <: DerivativeBoundaryCondition

Mixed boundary conditions combining value and derivative: β·∂u/∂n + α·u = g. Example: Convection.

source

Generic BC Types

Macchiato.PrescribedValueType
PrescribedValue{F<:Function} <: Dirichlet

Generic Dirichlet BC that prescribes a value via a function.

The function has signature f(x, t) -> value where:

  • x: spatial coordinate of the boundary point
  • t: time
  • Returns the prescribed value at that location and time

Fields

  • f: Function with signature (x, t) -> value
  • name: Display name (e.g., :Temperature, :PrescribedValue)

For built-in physics, use named constructors (e.g., Temperature, VelocityInlet). For custom PDEs, use the unparameterized constructors directly: PrescribedValue(0.0).

source
Macchiato.PrescribedFluxType
PrescribedFlux{F<:Function} <: Neumann

Generic Neumann BC that prescribes a flux (normal derivative) via a function.

The flux condition is: ∂u/∂n = f(x, t)

The function has signature f(x, t) -> flux_value where:

  • x: spatial coordinate of the boundary point
  • t: time
  • Returns the prescribed flux value at that location and time

Fields

  • f: Function with signature (x, t) -> flux
  • name: Display name (e.g., :HeatFlux, :PrescribedFlux)

For built-in physics, use named constructors (e.g., HeatFlux, Traction). For custom PDEs, use the unparameterized constructors directly: PrescribedFlux(1.0).

source
Macchiato.ZeroFluxType
ZeroFlux <: Neumann

Generic Neumann BC with zero flux: ∂u/∂n = 0.

Represents symmetry, insulation, or fully-developed flow depending on context.

Fields

  • name: Display name (e.g., :Adiabatic, :VelocityOutlet, :ZeroFlux)

For built-in physics, use named constructors (e.g., Adiabatic, VelocityOutlet). For custom PDEs, use the unparameterized constructor directly: ZeroFlux().

source
Note

These generic types work with any physics model. For custom PDEs, use them directly: PrescribedValue(0.0), PrescribedFlux(1.0), ZeroFlux(). See Custom PDEs for a complete example.

Energy BCs

Macchiato.TemperatureFunction
Temperature(value)

Prescribed temperature BC. Value can be a Number or Function (x, t) -> value.

source
Macchiato.HeatFluxFunction
HeatFlux(flux)

Prescribed heat flux BC: ∂T/∂n = q. Flux can be a Number or Function (x, t) -> flux.

source
Macchiato.ConvectionType
Convection(h, k, T∞)

Convective heat transfer: h·T + k·∂T/∂n = h·T∞(x,t). T∞ can be a Number or Function (x, t) -> ambient_temp.

source

Mechanics BCs

Macchiato.DisplacementType
Displacement{F<:Function} <: Dirichlet

Prescribed displacement BC for solid mechanics. The function returns a tuple of displacement components: f(x, t) -> (ux, uy) for 2D.

Constructors

Displacement((x, t) -> (0.0, 0.0))          # Function returning tuple
Displacement(0.0, 0.0)                       # Constant displacement
Displacement(ux::Function, uy::Function)     # Per-component functions
source
Macchiato.TractionType
Traction{F<:Function} <: Neumann

Prescribed traction BC for solid mechanics. The function returns a tuple of traction components: f(x, t) -> (tx, ty) for 2D.

In terms of stress: t = σ·n where n is the outward normal.

Constructors

Traction((x, t) -> (0.0, -1000.0))          # Function returning tuple
Traction(tx::Number, ty::Number)             # Constant traction
Traction(tx::Function, ty::Function)         # Per-component functions
source

Fluid BCs

Macchiato.VelocityInletFunction
VelocityInlet(velocity)

Prescribed velocity at inlet. Value can be a Number or Function (x, t) -> velocity.

source
Macchiato.PressureOutletFunction
PressureOutlet(pressure)

Prescribed pressure at outlet. Value can be a Number or Function (x, t) -> pressure.

source

Wall BC

Macchiato.WallFunction
Wall(velocity)
Wall()

No-slip wall or moving wall BC. Value can be a Number or Function (x, t) -> velocity. No arguments creates stationary wall (v=0).

source

Simulation

Macchiato.SimulationType
Simulation{D<:Domain}

High-level simulation container that manages time stepping, callbacks, and output.

Fields

  • domain::D: The computational domain with models and boundary conditions
  • mode::SimulationMode: Whether this is a transient or steady-state simulation
  • Δt::Union{Nothing, Float64}: Time step (transient only)
  • stop_time::Union{Nothing, Float64}: End time (transient only)
  • stop_iteration::Int: Maximum iterations (0 = unlimited)
  • solver::Symbol: ODE solver to use (transient only)
  • solver_options::NamedTuple: Additional solver options
  • callbacks::Dict{Symbol, Callback}: Named callbacks
  • output_writers::Dict{Symbol, AbstractOutputWriter}: Named output writers
  • u0::Union{Nothing, Vector{Float64}}: Initial condition vector
  • iteration::Int: Current iteration count
  • time::Float64: Current simulation time
  • running::Bool: Whether simulation is currently running
  • _integrator: ODE integrator (transient only)
  • _solution: Solution vector (steady-state or final transient)

Constructors

Transient simulation

sim = Simulation(domain; Δt=0.001, stop_time=1.0)
sim = Simulation(domain; Δt=0.001, stop_iteration=1000)

SteadyState-state simulation

sim = Simulation(domain)  # No Δt → steady-state mode

Mode Detection

  • If Δt and/or stop_time provided → transient mode
  • If neither provided → steady-state mode
source
Macchiato.run!Function
run!(sim::Simulation)

Execute the simulation. Dispatches to _run_transient! or _run_steady! based on mode.

Returns the simulation object for chaining.

source
Macchiato.set!Function
set!(sim::Simulation; kwargs...)

Set initial conditions for simulation fields.

Arguments

  • sim: Simulation to set initial conditions for
  • kwargs: Field name/value pairs

Supported value types

  • Number: Uniform value for entire field
  • Function: Called as f(x) where x is coordinate vector [x, y] or [x, y, z]
  • Vector: Direct assignment (must match field length)

Examples

set!(sim, T=300.0)                           # Uniform temperature
set!(sim, T=x -> 300 + 10*x[1])              # Temperature function of position
set!(sim, u=0.0, v=0.0, p=0.0)               # Multiple fields
source

Field Extraction

Macchiato.velocityFunction
velocity(sim) -> Tuple{Vector{Float64}, ...}

Extract velocity components from simulation. Returns (u, v) for 2D or (u, v, w) for 3D.

source
Macchiato.displacementFunction
displacement(sim) -> Tuple{Vector{Float64}, ...}

Extract displacement components from simulation. Returns (ux, uy) for 2D or (ux, uy, uz) for 3D.

source

Callbacks and Schedules

Macchiato.CallbackType
Callback{F, S<:AbstractSchedule, P}

A callback that executes a function according to a schedule.

Fields

  • func::F: Function to call. Signature: func(sim) or func(sim, parameters)
  • schedule::S: When to execute the callback
  • parameters::P: Optional parameters passed to the function
  • _state::ScheduleState: Internal state for schedule tracking

Example

function print_progress(sim)
    println("Iteration $(sim.iteration), t = $(sim.time)")
end
callback = Callback(print_progress, IterationInterval(100))
source
Macchiato.IterationIntervalType
IterationInterval(N::Int)

Schedule that triggers every N iterations.

Example

callback = Callback(print_progress, IterationInterval(100))  # Every 100 iterations
source
Macchiato.TimeIntervalType
TimeInterval(Δt::Float64)

Schedule that triggers every Δt simulation time units.

Example

writer = VTKOutputWriter("results/", schedule=TimeInterval(0.1))  # Every 0.1 time units
source
Macchiato.WallTimeIntervalType
WallTimeInterval(seconds::Float64)

Schedule that triggers every seconds wall-clock seconds.

Example

callback = Callback(checkpoint, WallTimeInterval(300.0))  # Every 5 minutes
source
Macchiato.SpecifiedTimesType
SpecifiedTimes(times::Vector{Float64})

Schedule that triggers at specific simulation times.

Example

callback = Callback(snapshot, SpecifiedTimes([0.0, 0.5, 1.0, 2.0]))
source

Output Writers

Macchiato.AbstractOutputWriterType
AbstractOutputWriter

Abstract type for output writers. Implementations must define:

  • initialize!(writer, sim) - called before simulation starts
  • write_output!(writer, sim) - called according to schedule
  • finalize!(writer, sim) - called after simulation ends
source
Macchiato.VTKOutputWriterType
VTKOutputWriter{S<:AbstractSchedule}

Writes simulation output to VTK format for visualization in ParaView.

Fields

  • prefix::String: Output file prefix (directory/filename without extension)
  • schedule::S: When to write output
  • fields::Vector{Symbol}: Which fields to output (empty = all fields)
  • _state::ScheduleState: Internal state for schedule tracking
  • _pvd::Any: ParaView collection file handle
  • _output_count::Int: Number of outputs written

Example

writer = VTKOutputWriter("results/sim", schedule=TimeInterval(0.1))
writer = VTKOutputWriter("results/sim", schedule=IterationInterval(100), fields=[:T])
source
Macchiato.JLD2OutputWriterType
JLD2OutputWriter{S<:AbstractSchedule}

Writes simulation checkpoints to JLD2 format for restart capability.

Fields

  • prefix::String: Output file prefix
  • schedule::S: When to write checkpoints
  • _state::ScheduleState: Internal state for schedule tracking
  • _output_count::Int: Number of checkpoints written

Example

writer = JLD2OutputWriter("checkpoints/sim", schedule=WallTimeInterval(300.0))
source

Solvers

Macchiato.MultiphysicsProblemFunction
MultiphysicsProblem(domain, u0, tspan; kwargs...)

Construct an ODEProblem for transient simulation from a Domain, an initial condition vector u0, and a time span tspan = (t_start, t_stop).

Assembles all model physics functions and boundary condition functions from the domain, then composes them into a single ODE right-hand side f(du, u, p, t).

source
SciMLBase.LinearProblemType
LinearSolve.LinearProblem(domain::Domain; scheme=nothing, kwargs...)

Construct a LinearProblem for steady-state simulation from a Domain.

Assembles the system matrix A and right-hand side b from the physics model, then applies boundary conditions by modifying the appropriate rows of A and b. The resulting system Ax = b is solved with LinearSolve.solve.

source

I/O

Macchiato.exportvtkFunction
exportvtk(filename, points, data, names)

Export point-based simulation results to a VTK file.

Arguments

  • filename::String: Output file path (without .vtu extension)
  • points::AbstractVector: Point cloud or coordinate vector
  • data::AbstractVector{<:AbstractVector}: Field data arrays to export
  • names::AbstractVector: Corresponding field names (e.g., ["T", "u"])

Example

exportvtk("results/temperature", points(cloud), [T_values], ["T"])
source

Operators

Macchiato.upwindFunction
upwind(data, eval_points, dim[, basis]; Δ=nothing, k=autoselect_k(data, basis))
upwind(data, dim[, basis]; Δ=nothing, k=autoselect_k(data, basis))

Build an upwind finite-difference-style operator using RBF interpolation.

Computes backward, forward, and centered partial derivatives with respect to dimension dim, then returns a function (ϕ, v, θ) that blends them based on flow direction v and upwind parameter θ ∈ [0, 1] (1 = full upwind, 0 = centered).

The single-argument form upwind(data, dim) evaluates at the data points themselves.

Arguments

  • data: Stencil points for RBF approximation
  • eval_points: Points where the derivative is evaluated
  • dim: Spatial dimension (1 = x, 2 = y, …)
  • basis: Radial basis function (default: PHS(3; poly_deg=2))
  • Δ: Virtual node offset distance (auto-detected if nothing)
  • k: Number of nearest neighbors for stencil
source