Reference

Base.isconstMethod
isconst(x::Constant)

True if the symbol value cannot be modified within an Operator (and thus its value is provided by the user directly from Python-land), False otherwise.

source
Base.ndimsMethod
ndims(x::DiscreteFunction)

Return the number of dimensions corresponding to the discrete function x.

source
Base.parentMethod
parent(x::AbstractSubDimension)

Returns the parent dimension of a subdimension.

Example

julia x = SpaceDimension(name="x") xr = SubDimensionRight(name="xr", parent=x, thickness=2) parent(xr)`

source
Base.sizeMethod
size(x::DiscreteFunction)

Return the shape of the grid for the discrete function x.

source
Devito.DerivativeMethod
Derivative(x::Union{Constant, Number}, args...; kwargs...)

Returns the derivative of a constant or number, which is zero.

source
Devito.DerivativeMethod
Derivative(x::Union{DiscreteFunction,PyObject}, args...; kwargs...)


An unevaluated Derivative, which carries metadata (Dimensions,
derivative order, etc) describing how the derivative will be expanded
upon evaluation.

Parameters
----------
expr : expr-like
    Expression for which the Derivative is produced.
dims : Dimension or tuple of Dimension
    Dimenions w.r.t. which to differentiate.
fd_order : int or tuple of int, optional
    Coefficient discretization order. Note: this impacts the width of
    the resulting stencil. Defaults to 1.
deriv_order: int or tuple of int, optional
    Derivative order. Defaults to 1.
side : Side or tuple of Side, optional
    Side of the finite difference location, centered (at x), left (at x - 1)
    or right (at x +1). Defaults to ``centered``.
transpose : Transpose, optional
    Forward (matvec=direct) or transpose (matvec=transpose) mode of the
    finite difference. Defaults to ``direct``.
subs : dict, optional
    Substitutions to apply to the finite-difference expression after evaluation.
x0 : dict, optional
    Origin (where the finite-difference is evaluated at) for the finite-difference
    scheme, e.g. Dict(x=> x, y => y + spacing(y)/2).

Examples
--------
Creation

```julia
using Devito
grid = Grid((10, 10))
y, x = dimensions(grid)
u = Devito.Function(name="u", grid=grid, space_order=2)
Derivative(u, x)

# You can also specify the order as a keyword argument

Derivative(u, x, deriv_order=2)

# Or as a tuple

Derivative(u, (x, 2))
```
source
Devito.MaxFunction
Max(args...)

Can be used in a Devito.Eq to return the minimum of a collection of arguments Example:

    eqmax = Eq(f,Max(g,1))

Is equivalent to f = Max(g,1) for Devito functions f,g

source
Devito.MinFunction
Min(args...)

Can be used in a Devito.Eq to return the minimum of a collection of arguments Example:

    eqmin = Eq(f,Min(g,1))

Is equivalent to f = Min(g,1) for Devito functions f,g

source
Devito.ModMethod
Mod(x::AbstractDimension,y::Int)

Perform Modular division on a dimension

source
Devito.applyMethod

apply( operator::Operator; kwargs...)

Execute the Devito operator, Operator.

See: https://www.devitoproject.org/devito/operator.html?highlight=apply#devito.operator.operator.Operator.apply

Note that this returns a summary::Dict of the action of applying the operator. This contains information such as the number of floating point operations executed per second.

source
Devito.backwardMethod
backward(x::TimeFunction)

Returns the symbol for the time-backward state of the TimeFunction.

See: https://www.devitoproject.org/devito/timefunction.html?highlight=forward#devito.types.TimeFunction.backward

source
Devito.ccodeMethod
ccode(x::Operator; filename="")

Print the ccode associated with a devito operator. If filename is provided, writes ccode to disk using that filename

source
Devito.configuration!Method
configuration!(key, value)

Configure Devito. Examples include

configuration!("log-level", "DEBUG")
configuration!("language", "openmp")
configuration!("mpi", false)
source
Devito.coordinatesMethod
coordinates(x::SparseDiscreteFunction)

Returns a Devito function associated with the coordinates of a sparse time function. Note that contrary to typical Julia convention, coordinate order is from slow-to-fast (Python ordering). Thus, for a 3D grid, the sparse time function coordinates would be ordered x,y,z.

source
Devito.coordinates_dataMethod
coordinates_data(x::SparseDiscreteFunction)

Returns a Devito array associated with the coordinates of a sparse time function. Note that contrary to typical Julia convention, coordinate order is from slow-to-fast (Python ordering). Thus, for a 3D grid, the sparse time function coordinates would be ordered x,y,z.

source
Devito.dataMethod
data(x::Constant{T})

Returns value(x::Constant{T}). See value documentation for more information.

source
Devito.dataMethod
data(x::DiscreteFunction)

Return the data associated with the grid that corresponds to the discrete function x. This is the portion of the grid that excludes the halo. In the case of non-MPI Devito, this returns an array of type DevitoArray. In the case of the MPI Devito, this returns an array of type DevitoMPIArray.

The data can be converted to an Array via convert(Array, data(x)). In the case where data(x)::DevitoMPIArray, this also collects the data onto MPI rank 0.

source
Devito.data_allocatedMethod
data_allocated(x::DiscreteFunction)

Return the data associated with the grid that corresponds to the discrete function x. This is the portion of the grid that includes the inner halo and includes the outer halo. We expect this to be equivalent to data_with_inhalo.

The data can be converted to an Array via convert(Array, data(x)). In the case where data(x)::DevitoMPIArray, this also collects the data onto MPI rank 0.

source
Devito.data_with_haloMethod
data_with_halo(x::DiscreteFunction)

Return the data associated with the grid that corresponds to the discrete function x. This is the portion of the grid that excludes the inner halo and includes the outer halo. In the case of non-MPI Devito, this returns an array of type DevitoArray. In the case of the MPI Devito, this returns an array of type DevitoMPIArray.

The data can be converted to an Array via convert(Array, data(x)). In the case where data(x)::DevitoMPIArray, this also collects the data onto MPI rank 0.

source
Devito.data_with_inhaloMethod
data_with_inhalo(x::DiscreteFunction)

Return the data associated with the grid that corresponds to the discrete function x. This is the portion of the grid that includes the inner halo and includes the outer halo. In the case of non-MPI Devito, this returns an array of type DevitoArray. In the case of the MPI Devito, this returns an array of type DevitoMPIArray.

The data can be converted to an Array via convert(Array, data(x)). In the case where data(x)::DevitoMPIArray, this also collects the data onto MPI rank 0.

source
Devito.dimensionsMethod
dimensions(x::Union{Grid,DiscreteFunction})

Returns a tuple with the dimensions associated with the Devito grid.

source
Devito.dtFunction
dt(f::TimeFunction, args...; kwargs...)

Returns the symbol for the first time derivative of a time function

source
Devito.dt2Function
dt2(f::TimeFunction, args...; kwargs...)

Returns the symbol for the second time derivative of a time function

source
Devito.dxFunction
dx(f::Union{DiscreteFunction,PyObject,Constant,Number}, args...; kwargs...)

Returns the symbol for the first derivative with respect to x if f is a Function with dimension x. Otherwise returns 0. Thus, the derivative of a function with respect to a dimension it doesn't have is zero, as is the derivative of a constant.

source
Devito.dx2Function
dx2(f::Union{DiscreteFunction,PyObject,Constant,Number}, args...; kwargs...)

Returns the symbol for the second derivative with respect to x if f is a Function with dimension x. Otherwise returns 0. Thus, the derivative of a function with respect to a dimension it doesn't have is zero, as is the derivative of a constant.

source
Devito.dxlFunction
dxl(f::DiscreteFunction, args...; kwargs...)

Returns the symbol for the first backward one-sided derivative with respect to x if f is a Function with dimension x. Otherwise returns 0. Thus, the derivative of a function with respect to a dimension it doesn't have is zero, as is the derivative of a constant.

source
Devito.dxrFunction
dxr(f::DiscreteFunction, args...; kwargs...)

Returns the symbol for the first forward one-sided derivative with respect to x if f is a Function with dimension x. Otherwise returns 0. Thus, the derivative of a function with respect to a dimension it doesn't have is zero, as is the derivative of a constant.

source
Devito.dyFunction
dy(f::DiscreteFunction, args...; kwargs...)

Returns the symbol for the first derivative with respect to y if f is a Function with dimension y. Otherwise returns 0. Thus, the derivative of a function with respect to a dimension it doesn't have is zero, as is the derivative of a constant.

source
Devito.dy2Function
dy2(f::DiscreteFunction, args...; kwargs...)

Returns the symbol for the second derivative with respect to y if f is a Function with dimension y. Otherwise returns 0. Thus, the derivative of a function with respect to a dimension it doesn't have is zero, as is the derivative of a constant.

source
Devito.dylFunction
dyl(f::DiscreteFunction, args...; kwargs...)

Returns the symbol for the first backward one-sided derivative with respect to y if f is a Function with dimension y. Otherwise returns 0. Thus, the derivative of a function with respect to a dimension it doesn't have is zero, as is the derivative of a constant.

source
Devito.dyrFunction
dyr(f::DiscreteFunction, args...; kwargs...)

Returns the symbol for the first forward one-sided derivative with respect to y if f is a Function with dimension y. Otherwise returns 0. Thus, the derivative of a function with respect to a dimension it doesn't have is zero, as is the derivative of a constant.

source
Devito.dzFunction
dzr(f::DiscreteFunction, args...; kwargs...)

Returns the symbol for the first forward one-sided derivative with respect to z if f is a Function with dimension z. Otherwise returns 0. Thus, the derivative of a function with respect to a dimension it doesn't have is zero, as is the derivative of a constant.

source
Devito.dzlFunction
dzl(f::DiscreteFunction, args...; kwargs...)

Returns the symbol for the first backward one-sided derivative with respect to z if f is a Function with dimension y. Otherwise returns 0. Thus, the derivative of a function with respect to a dimension it doesn't have is zero, as is the derivative of a constant.

source
Devito.forwardMethod
forward(x::TimeFunction)

Returns the symbol for the time-forward state of the TimeFunction.

See: https://www.devitoproject.org/devito/timefunction.html?highlight=forward#devito.types.TimeFunction.forward

source
Devito.gridMethod
grid(f::DiscreteFunction)

Return the grid corresponding to the discrete function f.

source
Devito.haloMethod
halo(x::DiscreteFunction)

Return the Devito "outer" halo size corresponding to the discrete function f.

source
Devito.inhaloMethod
inhalo(x::DiscreteFunction)

Return the Devito "inner" halo size used for domain decomposition, and corresponding to the discrete function f.

source
Devito.injectMethod
inject(x::SparseDiscreteFunction; kwargs...)

Generate equations injecting an arbitrary expression into a field.

See: Generate equations injecting an arbitrary expression into a field.

Example

x = SpaceDimension(name="x", spacing=Constant(name="h_x", value=5.0))
z = SpaceDimension(name="z", spacing=Constant(name="h_z", value=5.0))
grid = Grid(
    dimensions = (x,z),
    shape = (251,501), # assume x is first, z is second (i.e. z is fast in python)
    origin = (0.0,0.0),
    extent = (1250.0,2500.0),
    dtype = Float32)

time_range = 0.0f0:0.5f0:1000.0f0
src = SparseTimeFunction(name="src", grid=grid, npoint=1, nt=length(time_range))
src_term = inject(src; field=forward(p), expr=2*src)
source
Devito.interpolateMethod
interpolate(x::SparseDiscreteFunction; kwargs...)

Generate equations interpolating an arbitrary expression into self.

See: https://www.devitoproject.org/devito/sparsetimefunction.html#devito.types.SparseTimeFunction.interpolate

Example

x = SpaceDimension(name="x", spacing=Constant(name="h_x", value=5.0))
z = SpaceDimension(name="z", spacing=Constant(name="h_z", value=5.0))
grid = Grid(
    dimensions = (z,x),
    shape = (501,251), # assume z is first, x is second
    origin = (0.0,0.0),
    extent = (2500.0,1250.0),
    dtype = Float32)

p = TimeFunction(name="p", grid=grid, time_order=2, space_order=8)

time_range = 0.0f0:0.5f0:1000.0f0
nz,nx,δz,δx = size(grid)...,spacing(grid)...
rec = SparseTimeFunction(name="rec", grid=grid, npoint=nx, nt=length(time_range))
rec_coords = coordinates_data(rec)
rec_coords[1,:] .= 10.0
rec_coords[2,:] .= δx*(0:nx-1)


rec_term = interpolate(rec, expr=p)
source
Devito.is_DerivedFunction
is_Derived(x::Dimension)

Returns true when dimension is derived, false when it is not

source
Devito.nameMethod
name(x::Union{SubDomain, DiscreteFunction, TimeFunction, Function, Constant, AbstractDimension, Operator})

returns the name of the Devito object

source
Devito.nsimplifyMethod
nsimplify(expr::PyObject; constants=(), tolerance=none, full=false, rational=none, rational_conversion="base10")

Wrapper around sympy.nsimplify. Find a simple representation for a number or, if there are free symbols or if $rational=True$, then replace Floats with their Rational equivalents. If no change is made and rational is not False then Floats will at least be converted to Rationals.

Explanation

For numerical expressions, a simple formula that numerically matches the given numerical expression is sought (and the input should be possible to evalf to a precision of at least 30 digits). Optionally, a list of (rationally independent) constants to include in the formula may be given. A lower tolerance may be set to find less exact matches. If no tolerance is given then the least precise value will set the tolerance (e.g. Floats default to 15 digits of precision, so would be tolerance=10**-15). With $full=True$, a more extensive search is performed (this is useful to find simpler numbers when the tolerance is set low). When converting to rational, if rationalconversion='base10' (the default), then convert floats to rationals using their base-10 (string) representation. When rationalconversion='exact' it uses the exact, base-2 representation.

See https://github.com/sympy/sympy/blob/52f606a503cea5e9588de14150ccb9f7f9ed4752/sympy/simplify/simplify.py .

Examples:

nsimplify(π) # PyObject 314159265358979/100000000000000
nsimplify(π; tolerance=0.1) # PyObject 22/7
source
Devito.size_with_haloMethod
size_with_halo(x::DiscreteFunction)

Return the size of the grid associated with x, inclusive of the Devito "outer" halo.

source
Devito.size_with_inhaloMethod
size_with_inhalo(x::DiscreteFunction)

Return the size of the grid associated with z, inclusive the the Devito "inner" and "outer" halos.

source
Devito.solveMethod
solve(eq::PyObject, target::PyObject; kwargs...)

Algebraically rearrange an Eq w.r.t. a given symbol. This is a wrapper around $devito.solve$, which in turn is a wrapper around $sympy.solve$.

Parameters

  • eq::PyObject expr-like. The equation to be rearranged.
  • target::PyObject The symbol w.r.t. which the equation is rearranged. May be a Function or any other symbolic object.

kwargs

  • Symbolic optimizations applied while rearranging the equation. For more information. refer to sympy.solve.__doc__.
source
Devito.space_orderMethod
space_order(x::Union{TimeFunction,Function})

Returns the space order for spatial derivatives defined on the associated TimeFunction or Function

source
Devito.spacingFunction
spacing(x::Dimension)

Symbol representing the physical spacing along the Dimension.

source
Devito.subsMethod
subs(f::DiscreteFunction{T,N,M},dict::Dict)

Perform substitution on the dimensions of Devito Discrete Function f based on a dictionary.

Example

    grid = Grid(shape=(10,10,10))
    z,y,x = dimensions(grid)
    f = Devito.Function(grid=grid, name="f", staggered=x)
    subsdict = Dict(x=>x-spacing(x)/2)
    g = subs(f,subsdict)
source
Devito.thicknessMethod
thickness(x::AbstractSubDimension)

Returns a tuple of a tuple containing information about the left and right thickness of a SubDimension and a symbol corresponding to each side's thickness.

Example

x = SpaceDimension(name="x")
xr = SubDimensionRight(name="xr", parent=x, thickness=2)
thickness(xr)
source
Devito.time_dimMethod
time_dim(x::Union{Grid,TimeFunction})

Returns the time dimension for the associated object.

source
Devito.value!Method
value!(x::Constant{T},y::T)

Change the numerical value of a constant, x, after creation to y, after converting y to datatype T of constant x.

source
Devito.valueMethod
value(x::Constant)

Returns the value of a devito constant. Can not be used to change constant value, for that use value!(x,y)

source
Devito.BufferMethod
Buffer(value::Int)

Construct a devito buffer. This may be used as a save= keyword argument in the construction of TimeFunctions.

source
Devito.ConditionalDimensionType
ConditionalDimension(;kwargs)

Symbol defining a non-convex iteration sub-space derived from a parent Dimension, implemented by the compiler generating conditional “if-then” code within the parent Dimension’s iteration space.

See: https://www.devitoproject.org/devito/dimension.html?highlight=conditional#devito.types.dimension.ConditionalDimension

Example

size, factor = 16, 4
i  = SpaceDimension(name="i")
grid = Grid(shape=(size,),dimensions=(i,))
ci = ConditionalDimension(name="ci", parent=i, factor=factor)
g  = Devito.Function(name="g", grid=grid, shape=(size,), dimensions=(i,))
f  = Devito.Function(name="f", grid=grid, shape=(div(size,factor),), dimensions=(ci,))
op = Operator([Eq(g, 1), Eq(f, g)],name="Cond")
source
Devito.ConstantMethod
Constant(args...; kwargs...)

Symbol representing a constant, scalar value in symbolic equations. A Constant carries a scalar value.

kwargs

  • name::String Name of the symbol.
  • value::Real Value associated with the symbol. Defaults to 0.
  • dtype::Type{AbstractFloat} choose from Float32 or Float64. Default is Float32
source
Devito.FunctionMethod
Devito.Function(; kwargs...)

Tensor symbol representing a discrete function in symbolic equations.

See: https://www.devitoproject.org/devito/function.html?highlight=function#devito.types.Function

Example

x = SpaceDimension(name="x", spacing=Constant(name="h_x", value=5.0))
z = SpaceDimension(name="z", spacing=Constant(name="h_z", value=5.0))
grid = Grid(
    dimensions = (x,z),
    shape = (251,501), # assume x is first, z is second (i.e. z is fast in python)
    origin = (0.0,0.0),
    extent = (1250.0,2500.0),
    dtype = Float32)

b = Devito.Function(name="b", grid=grid, space_order=8)
source
Devito.GridMethod
Grid(; shape[, optional key-word arguments...])

Construct a grid that can be used in the construction of a Devito Functions.

See: https://www.devitoproject.org/devito/grid.html?highlight=grid#devito.types.Grid

Example

x = SpaceDimension(name="x", spacing=Constant(name="h_x", value=5.0))
z = SpaceDimension(name="z", spacing=Constant(name="h_z", value=5.0))
grid = Grid(
    dimensions = (x,z), # z is fast (row-major)
    shape = (251,501),
    origin = (0.0,0.0),
    extent = (1250.0,2500.0),
    dtype = Float32)
source
Devito.OperatorType
Operator(expressions...[; optional named arguments])

Generate, JIT-compile and run C code starting from an ordered sequence of symbolic expressions, and where you provide a list of expressions defining the computation.

See: https://www.devitoproject.org/devito/operator.html?highlight=operator#devito.operator.operator.Operator""

Optional named arguments

  • name::String Name of the Operator, defaults to “Kernel”.
  • subs::Dict Symbolic substitutions to be applied to expressions.
  • opt::String The performance optimization level. Defaults to configuration["opt"].
  • language::String The target language for shared-memory parallelism. Defaults to configuration["language"].
  • platform::String The architecture the code is generated for. Defaults to configuration["platform"].
  • compiler::String The backend compiler used to jit-compile the generated code. Defaults to configuration["compiler"].

Example

Assuming that one has constructed the following Devito expressions: stencil_p, src_term and rec_term,

op = Operator([stencil_p, src_term, rec_term]; name="opIso")
source
Devito.SpaceDimensionType
SpaceDimension(;kwargs...)

Construct a space dimension that defines the extend of a physical grid.

See https://www.devitoproject.org/devito/dimension.html?highlight=spacedimension#devito.types.dimension.SpaceDimension.

Example

julia x = SpaceDimension(name="x", spacing=Constant(name="h_x", value=5.0))`

source
Devito.SparseFunctionMethod
SparseFunction(; kwargs...)

Tensor symbol representing a sparse array in symbolic equations.

See: https://www.devitoproject.org/devito/sparsefunction.html

Example

x = SpaceDimension(name="x", spacing=Constant(name="h_x", value=5.0))
z = SpaceDimension(name="z", spacing=Constant(name="h_z", value=5.0))
grid = Grid(
    dimensions = (x,z),
    shape = (251,501), # assume x is first, z is second (i.e. z is fast in python)
    origin = (0.0,0.0),
    extent = (1250.0,2500.0),
    dtype = Float32)

src = SparseFunction(name="src", grid=grid, npoint=1)
source
Devito.SparseTimeFunctionMethod
SparseTimeFunction(; kwargs...)

Tensor symbol representing a space- and time-varying sparse array in symbolic equations.

See: https://www.devitoproject.org/devito/sparsetimefunction.html?highlight=sparsetimefunction

Example

x = SpaceDimension(name="x", spacing=Constant(name="h_x", value=5.0))
z = SpaceDimension(name="z", spacing=Constant(name="h_z", value=5.0))
grid = Grid(
    dimensions = (x,z),
    shape = (251,501), # assume x is first, z is second (i.e. z is fast in python)
    origin = (0.0,0.0),
    extent = (1250.0,2500.0),
    dtype = Float32)

time_range = 0.0f0:0.5f0:1000.0f0
src = SparseTimeFunction(name="src", grid=grid, npoint=1, nt=length(time_range))
source
Devito.SubDimensionLeftType
SubDimensionLeft(args...; kwargs...)

Creates middle a SubDimension. Equivalent to devito.SubDimension.left helper function.

Example

x = SpaceDimension(name="x")
xm = SubDimensionLeft(name="xl", parent=x, thickness=2)
source
Devito.SubDimensionMiddleType
SubDimensionMiddle(args...; kwargs...)

Creates middle a SubDimension. Equivalent to devito.SubDimension.middle helper function.

Example

x = SpaceDimension(name="x")
xm = SubDimensionMiddle(name="xm", parent=x, thickness_left=2, thickness_right=3)
source
Devito.SubDimensionRightType
SubDimensionRight(args...; kwargs...)

Creates right a SubDimension. Equivalent to devito.SubDimension.right helper function.

Example

x = SpaceDimension(name="x")
xr = SubDimensionRight(name="xr", parent=x, thickness=3)
source
Devito.SubDomainMethod
SubDomain(name, instructions)

Create a subdomain by passing a list of instructions for each dimension. Using an instruction with (nothing,) implies that the whole dimension should be used for that subdomain, as will ("middle",0,0)

Examples:

instructions = ("left",2),("middle",3,3)
SubDomain("subdomain_name",instructions)

or

instructions = [("right",4),("middle",1,2)]
SubDomain("subdomain_name",instructions)

or

SubDomain("subdomain_name",("right",2),("left",1))
source
Devito.TimeFunctionMethod
TimeFunction(; kwargs...)

Tensor symbol representing a discrete function in symbolic equations.

See https://www.devitoproject.org/devito/timefunction.html?highlight=timefunction#devito.types.TimeFunction. Note that setting lazy=false will force all of the time steps to be serialized to disk and disable lazy streaming.

Example

x = SpaceDimension(name="x", spacing=Constant(name="h_x", value=5.0))
z = SpaceDimension(name="z", spacing=Constant(name="h_z", value=5.0))
grid = Grid(
    dimensions = (x,z),
    shape = (251,501), # assume x is first, z is second (i.e. z is fast in python)
    origin = (0.0,0.0),
    extent = (1250.0,2500.0),
    dtype = Float32)

p = TimeFunction(name="p", grid=grid, time_order=2, space_order=8)
source