Reference
Base.getindex
— MethodGet symbolic representation for function index object
Base.isconst
— Methodisconst(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.
Base.ndims
— Methodndims(x::DiscreteFunction)
Return the number of dimensions corresponding to the discrete function x
.
Base.parent
— Methodparent(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)
`
Base.size
— Methodsize(x::DiscreteFunction)
Return the shape of the grid for the discrete function x
.
Devito.Derivative
— MethodDerivative(x::Union{Constant, Number}, args...; kwargs...)
Returns the derivative of a constant or number, which is zero.
Devito.Derivative
— MethodDerivative(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))
```
Devito.Max
— FunctionMax(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
Devito.Min
— FunctionMin(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
Devito.Mod
— MethodMod(x::AbstractDimension,y::Int)
Perform Modular division on a dimension
Devito.apply
— Methodapply( 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.
Devito.backward
— Methodbackward(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
Devito.ccode
— Methodccode(x::Operator; filename="")
Print the ccode associated with a devito operator. If filename is provided, writes ccode to disk using that filename
Devito.configuration!
— Methodconfiguration!(key, value)
Configure Devito. Examples include
configuration!("log-level", "DEBUG")
configuration!("language", "openmp")
configuration!("mpi", false)
Devito.coordinates
— Methodcoordinates(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.
Devito.coordinates_data
— Methodcoordinates_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.
Devito.data
— Methoddata(x::Constant{T})
Returns value(x::Constant{T})
. See value
documentation for more information.
Devito.data
— Methoddata(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.
Devito.data_allocated
— Methoddata_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.
Devito.data_with_halo
— Methoddata_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.
Devito.data_with_inhalo
— Methoddata_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.
Devito.dimensions
— Methoddimensions(x::Union{Grid,DiscreteFunction})
Returns a tuple with the dimensions associated with the Devito grid.
Devito.dt
— Functiondt(f::TimeFunction, args...; kwargs...)
Returns the symbol for the first time derivative of a time function
Devito.dt2
— Functiondt2(f::TimeFunction, args...; kwargs...)
Returns the symbol for the second time derivative of a time function
Devito.dx
— Functiondx(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.
Devito.dx2
— Functiondx2(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.
Devito.dxl
— Functiondxl(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.
Devito.dxr
— Functiondxr(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.
Devito.dy
— Functiondy(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.
Devito.dy2
— Functiondy2(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.
Devito.dyl
— Functiondyl(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.
Devito.dyr
— Functiondyr(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.
Devito.dz
— Functiondzr(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.
Devito.dzl
— Functiondzl(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.
Devito.evaluate
— Methodevaluate(x::PyObject)
Evaluate a PyCall expression
Devito.forward
— Methodforward(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
Devito.grid
— Methodgrid(f::DiscreteFunction)
Return the grid corresponding to the discrete function f
.
Devito.halo
— Methodhalo(x::DiscreteFunction)
Return the Devito "outer" halo size corresponding to the discrete function f
.
Devito.indexed
— MethodThe wrapped IndexedData object.
Devito.inhalo
— Methodinhalo(x::DiscreteFunction)
Return the Devito "inner" halo size used for domain decomposition, and corresponding to the discrete function f
.
Devito.inject
— Methodinject(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)
Devito.interior
— Methodinterior(x::grid)
returns the interior subdomain of a Devito grid
Devito.interpolate
— Methodinterpolate(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)
Devito.is_Derived
— Functionis_Derived(x::Dimension)
Returns true when dimension is derived, false when it is not
Devito.name
— Methodname(x::Union{SubDomain, DiscreteFunction, TimeFunction, Function, Constant, AbstractDimension, Operator})
returns the name of the Devito object
Devito.nsimplify
— Methodnsimplify(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
Devito.origin
— Methodorigin(grid)
returns the tuple corresponding to the grid's origin
Devito.size_with_halo
— Methodsize_with_halo(x::DiscreteFunction)
Return the size of the grid associated with x
, inclusive of the Devito "outer" halo.
Devito.size_with_inhalo
— Methodsize_with_inhalo(x::DiscreteFunction)
Return the size of the grid associated with z
, inclusive the the Devito "inner" and "outer" halos.
Devito.solve
— Methodsolve(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 aFunction
or any other symbolic object.
kwargs
- Symbolic optimizations applied while rearranging the equation. For more information. refer to
sympy.solve.__doc__
.
Devito.space_order
— Methodspace_order(x::Union{TimeFunction,Function})
Returns the space order for spatial derivatives defined on the associated TimeFunction or Function
Devito.spacing
— Functionspacing(x::Dimension)
Symbol representing the physical spacing along the Dimension.
Devito.stepping_dim
— Methodstepping_dim(x::Grid)
Returns the stepping dimension for the associated grid.
Devito.subdomains
— Methodsubdomains(grid)
returns subdomains associated with a Devito grid
Devito.subs
— Methodsubs(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)
Devito.symbolic_max
— Functionsymbolic_max(x::Dimension)
Symbol defining the maximum point of the Dimension
Devito.symbolic_min
— Functionsymbolic_min(x::Dimension)
Symbol defining the minimum point of the Dimension
Devito.symbolic_size
— Functionsymbolic_size(x::Dimension)
Symbol defining the size of the Dimension
Devito.thickness
— Methodthickness(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)
Devito.time_dim
— Methodtime_dim(x::Union{Grid,TimeFunction})
Returns the time dimension for the associated object.
Devito.value!
— Methodvalue!(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.
Devito.value
— Methodvalue(x::Constant)
Returns the value of a devito constant. Can not be used to change constant value, for that use value!(x,y)
Devito.Buffer
— MethodBuffer(value::Int)
Construct a devito buffer. This may be used as a save= keyword argument in the construction of TimeFunctions.
Devito.Byref
— TypeSymbolic representation of the C notation &expr
.
Devito.Cast
— TypeSymbolic representation of the C notation (type)expr
.
Devito.ConditionalDimension
— TypeConditionalDimension(;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")
Devito.Constant
— MethodConstant(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 fromFloat32
orFloat64
. Default isFloat32
Devito.Deref
— TypeSymbolic representation of the C notation *expr
.
Devito.Function
— MethodDevito.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)
Devito.Grid
— MethodGrid(; 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)
Devito.Operator
— TypeOperator(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")
Devito.Pointer
— TypeSymbolic representation of a pointer in C
Devito.SpaceDimension
— TypeSpaceDimension(;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))
`
Devito.SparseFunction
— MethodSparseFunction(; 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)
Devito.SparseTimeFunction
— MethodSparseTimeFunction(; 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))
Devito.SubDimensionLeft
— TypeSubDimensionLeft(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)
Devito.SubDimensionMiddle
— TypeSubDimensionMiddle(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)
Devito.SubDimensionRight
— TypeSubDimensionRight(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)
Devito.SubDomain
— MethodSubDomain(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))