Index

CartesianGrids.CircularConvolutionType
CircularConvolution{M, N, T}

A preplanned, circular convolution operator on an M × N matrix of data of type T

Fields

  • : DFT coefficients of the convolution kernel
  • F: preplanned rFFT operator
  • F⁻¹: preplanned irFFT operator
  • nthreads : optimized number of threads to use, if appropriate
  • paddedSpace: scratch space to zero-pad the input matrix
  • : scratch space to store the DFT coefficients of the zero-padded input matrix

Constructors:

  • CircularConvolution(G::Matrix{T})

Example:

julia> G = repeat(1.0:3,1,4)
3×4 Array{Float64,2}:
 1.0  1.0  1.0  1.0
 2.0  2.0  2.0  2.0
 3.0  3.0  3.0  3.0

julia> C = CircularConvolution(G)
Circular convolution on a 3 × 4 matrix of data type Float64

julia> C*reshape(1:12, 3, 4)
3×4 Array{Int64,2}:
 164  164  164  164
 130  130  130  130
 148  148  148  148
source
CartesianGrids.DDFMethod
DDF([ddftype=Roma],[dx=1.0])

Construct a discrete delta function operator. This is generally only needed internally by the Regularize operator, so the user doesn't have much need for accessing this directly. The default DDF is the Roma function, which has a support of 3 grid cells. Other choices are the Goza operator, which is a truncated Gaussian with 28 cells support, and the Witchhat, which has 2 cells support. The resulting operator is evaluated with one, two or three coordinate arguments, producing, respectively, 1-d, 2-d, or 3-d smeared delta functions. It can also be called with the usual Julia vectorized dot notation with arrays of arguments. The optional cell spacing argument dx rescales the coordinates by this spacing, and the result is also rescaled by this spacing (raised to the number of dimensions). This spacing argument defaults to 1.0.

julia> ddf = DDF(ddftype=Roma)
Discrete delta function operator of type CartesianGrids.Roma, with spacing 1.0

julia> ddf(1)
0.16666666666666666

julia> ddf(-1)
0.16666666666666666

julia> ddf.([-1,0,1])
3-element Array{Float64,1}:
 0.16666666666666666
 0.6666666666666666
 0.16666666666666666
source
CartesianGrids.EdgeGradientType
EdgeGradient

EdgeGradient is a wrapper for tensor-valued data that lie partly at the nodes of dual cells and primary cells. EdgeGradient type data have fields dudx, dudy, dvdx, dvdy for the components of the tensor field. The diagonal components lie at one set of nodes (e.g. Primal), and the offdiagonal at the other set (e.g. Dual).

Constructors

  • EdgeGradient(C,dims) creates a tensor field of zeros in cells of type C (where C is either Dual or Primal), on a grid of dimensions dims. Note that dims represent the number of dual cells on the grid.
  • EdgeGradient(C,w) performs the same construction, but uses existing field data w of GridData type to determine the size of the grid.
  • Adding the dtype= keyword allows the data type of the field data to be changed. The default is Float64, but can be changed to, e.g., ComplexF64
source
CartesianGrids.EdgesType
Edges

Edges is a wrapper for vector-valued data that lie at the faces of either dual cells or primary cells. Edges type data have fields u and v for the components of the vector field. These are the normal components of the vector field on the vertical and horizontal faces of the corresponding cell.

Constructors

  • Edges(C,dims) creates a vector field of zeros in cells of type C (where C is either Dual or Primal), on a grid of dimensions dims. Note that dims represent the number of dual cells on the grid.
  • Edges(C,w) performs the same construction, but uses existing field data w of GridData type to determine the size of the grid.
  • Adding the dtype= keyword allows the data type of the field data to be changed. The default is Float64, but can be changed to, e.g., ComplexF64
source
CartesianGrids.GeneratedFieldType
GeneratedField(d::GridData,field::AbstractSpatialField...,grid::PhysicalGrid)

Create an instance of a spatial field function field on scalar grid data d, based on a grid grid. After creating the instance g = GeneratedField(d,field,grid), then the resulting grid data can be accessed by typing g(). For vector grid data, a separate field must be supplied for each component.

If the fields are time dependent, then you can also evaluate g(t) at the desired time. The time argument is ignored if the fields are static.

source
CartesianGrids.InterpolationMatrixType
InterpolationMatrix(H::Regularize,u::CellData,f::PointData) -> Emat

Construct and store a matrix representation of interpolation associated with H for data of type u to data of type f. The resulting matrix Emat can then be used to apply on grid data of type u to interpolate it to point data of type f, using mul!(f,Emat,u). It can also be used as just Emat*u.

source
CartesianGrids.ModulatedFieldType
ModulatedField(g::GeneratedField,modfcn::Abstract1DProfile)

Create a time-modulated form of a generated spatial field, useful for introducing a forcing field onto the grid. The supplied field g is modulated by a function modfcn with a specified profle shape. The resulting object can be evaluated with a single argument (time) and returns a GridData type object of the same type contained in g.

source
CartesianGrids.NodePairType
NodePair

NodePair is a wrapper for vector-valued data that lie at the nodes of dual cells and primal cells. NodePair type data have fields u and v for the components of the vector field. These are the normal components of a vector field on nodes that form the faces of a virtual cell centered at one of the faces of the primal cell.

Constructors

  • NodePair(C,dims) creates a vector field of zeros in cells of type C (where C is either Dual or Primal), on a grid of dimensions dims. Note that dims represent the number of dual cells on the grid.
  • NodePair(C,w) performs the same construction, but uses existing field data w of GridData type to determine the size of the grid.
  • Adding the dtype= keyword allows the data type of the field data to be changed. The default is Float64, but can be changed to, e.g., ComplexF64
source
CartesianGrids.NodesType
Nodes

Nodes is a wrapper for scalar-valued data that lie at the centers of either dual cells or primary cells. A Nodes type can be accessed by indexing like any other array, and allows the use of size, similar, zero functions.

Constructors

  • Nodes(C,dims) creates a field of zeros in cells of type C (where C is either Dual or Primal), on a grid of dimensions dims (a tuple). Note that dims represent the number of dual cells on the grid, even if C is Primal.
  • Nodes(C,w) performs the same construction, but uses existing field data w of GridData type to determine the size of the grid.
  • Adding the dtype= keyword allows the data type of the field data to be changed. The default is Float64, but can be changed to, e.g., ComplexF64
source
CartesianGrids.PhysicalGridMethod
PhysicalGrid(xlim::Tuple{Real,Real},ylim::Tuple{Real,Real},Δx::Float64;[optimize=true,nthreads_max=length(Sys.cpu_info())])

Constructor to set up a grid connected to physical space. The region to be discretized by the grid is defined by the limits xlim and ylim, and the cell spacing (uniform and indentical in each direction) is specified by Δx. The constructor uses this information to determine the number of cells in each direction, expanding the given range if necessary to accommodate an integer number. It also pads each side with a ghost cell. It also determines the indices corresponding to the corner of the cell to which the physical origin corresponds. Note that the corner corresponding to the lowest limit in each direction has indices (1,1).

source
CartesianGrids.RegularizationMatrixType
RegularizationMatrix(H::Regularize,f::PointData,u::CellData) -> Hmat

Construct and store a matrix representation of regularization associated with H for data of type f to data of type u. The resulting matrix Hmat can then be used to apply on point data of type f to regularize it to grid data of type u, using mul!(u,Hmat,f). It can also be used as just Hmat*f.

If H is a symmetric regularization and interpolation operator, then this actually returns a tuple Hmat, Emat, where Emat is the interpolation matrix.

source
CartesianGrids.RegularizeMethod
Regularize(x,y,dx,[ddftype=Yang3],[graddir=0],[I0=(1,1)], [weights=1.0], [filter=false],
                   [issymmetric=false])

Constructor to set up an operator for regularizing and interpolating data from/to points immersed in the grid to/from fields on the grid itself. The supplied x and y represent physical coordinates of the immersed points, and dx denotes a uniform physical cell size of the grid. The separate arguments x and y can be replaced by a single argument X of type VectorData holding the coordinates.

The operations of regularization and interpolation are carried out with a discrete delta function (ddf), which defaults to the type Yang3. Others are also possible, such as Roma, Goza or M3. The optional argument graddir, if set to 1 or 2, will generate an interpolation operator that evaluates the negative of the respective component of the gradient of a grid field at the immersed points. The default value of this argument is 0, which simply interpolates. Note that the regularization form of this gradient type is also possible.

The optional tuple I0 represents the indices of the primary node that coincides with (x,y) = (0,0). This defaults to (1,1), which leaves one layer of ghost (dual) cells and sets the physical origin in the lower left corner of the grid of interior dual cells.

Another optional parameter, weights, sets the weight of each point in the regularization. This would generally be set with, say, the differential arc length for regularization of data on a curve. It can be a vector (of the same length as x and y) or a scalar if uniform. It defaults to 1.0.

The optional Boolean parameter filter can be set to true if it is desired to apply filtering (see Goza et al, J Comput Phys 2016) to the grid data before interpolating. This is generally only used in the context of preconditioning the solution for forces on the immersed points.

If the optional Boolean parameter issymmetric is set to true, then the regularization and interpolation are constructed to be transposes of each other. Note that this option overrides any supplied weights. The default of this parameter is false.

The resulting operator can be used in either direction, regularization and interpolation, with the first argument representing the target (the entity to regularize/interpolate to), and the second argument the source (the entity to regularize/interpolate from). The regularization does not use the filtering option.

Example

In the example below, we set up a 12 x 12 grid. Using the default value for I0 and setting dx = 0.1, the physical dimensions of the non-ghost part of the grid are 1.0 x 1.0. Three points are set up in the interior, and a vector field is assigned to them, with the x component of each of them set to 1.0. These data are regularized to a field of primal edges on the grid, using the Roma DDF kernel.

julia> x = [0.25,0.75,0.25]; y = [0.75,0.25,0.25];

julia> X = VectorData(x,y);

julia> q = Edges(Primal,(12,12));

julia> dx = 0.1;

julia> H = Regularize(x,y,dx;ddftype=Roma)
Regularization/interpolation operator with non-filtered interpolation
  DDF type CartesianGrids.Roma
  3 points in grid with cell area 0.01

julia> f = VectorData(X);

julia> fill!(f.u,1.0);

julia> H(q,f)
Edges{Primal,12,12,Float64} data
u (in grid orientation)
11×12 Array{Float64,2}:
 0.0  0.0  0.0       0.0     0.0      …  0.0       0.0     0.0      0.0  0.0
 0.0  0.0  0.0       0.0     0.0         0.0       0.0     0.0      0.0  0.0
 0.0  0.0  8.33333  33.3333  8.33333     0.0       0.0     0.0      0.0  0.0
 0.0  0.0  8.33333  33.3333  8.33333     0.0       0.0     0.0      0.0  0.0
 0.0  0.0  0.0       0.0     0.0         0.0       0.0     0.0      0.0  0.0
 0.0  0.0  0.0       0.0     0.0      …  0.0       0.0     0.0      0.0  0.0
 0.0  0.0  0.0       0.0     0.0         0.0       0.0     0.0      0.0  0.0
 0.0  0.0  8.33333  33.3333  8.33333     8.33333  33.3333  8.33333  0.0  0.0
 0.0  0.0  8.33333  33.3333  8.33333     8.33333  33.3333  8.33333  0.0  0.0
 0.0  0.0  0.0       0.0     0.0         0.0       0.0     0.0      0.0  0.0
 0.0  0.0  0.0       0.0     0.0      …  0.0       0.0     0.0      0.0  0.0
v (in grid orientation)
12×11 Array{Float64,2}:
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
source
CartesianGrids.ScalarDataType
ScalarData <: PointData

A wrapper for a one-dimensional array of scalar-valued data. The resulting wrapper can be indexed in the same way as the array itself.

Constructors

  • ScalarData(d::AbstractVector[,dtype=Float64]) constructs a wrapper for the one-dimensional array of data d
  • ScalarData(n::Int) constructs a wrapper for an array of zeros of length n.
  • ScalarData(x::PointData) constructs a wrapper for an array of zeros of the same length as that wrapped by x.
  • ScalarData(n::Int,dtype=ComplexF64) constructs a wrapper for complex-valued data.
  • ScalarData(x::PointData,dtype=ComplexF64) constructs a wrapper for an array of complex zeros of the same length as that wrapped by x.

Example

julia> f = ScalarData(10);

julia> f[5] = 1.0;

julia> f
10 points of scalar-valued Float64 data
10-element Array{Float64,1}:
 0.0
 0.0
 0.0
 0.0
 1.0
 0.0
 0.0
 0.0
 0.0
 0.0
source
CartesianGrids.TensorDataType
TensorData <: PointData

A wrapper for a one-dimensional array of 2x2 tensor-valued data, with fields dudx, dudy, dvdx, dvdy. The resulting wrapper can be indexed as though these four components are stacked on top of each other.

Constructors

  • TensorData(d::AbstractVector[,dtype=Float64]) constructs a wrapper for the one-dimensional array of data d, splitting d into the four components evenly.
  • TensorData(dudx,dudy,dvdx,dvdy) constructs a wrapper for the tensor components data, each of type AbstractVector
  • TensorData(n::Int) constructs a wrapper with zeros of length n for all components.
  • TensorData(x::PointData[,dtype=Float64]) constructs a wrapper for zero components of the same length as that wrapped by x.

Example

source
CartesianGrids.VectorDataType
VectorData <: PointData

A wrapper for a one-dimensional array of two-component vector-valued data. The resulting wrapper can be indexed as though the first component and second component are stacked on top of each other.

Constructors

  • VectorData(d::AbstractVector[,dtype=Float64]) constructs a wrapper for the one-dimensional array of data d, splitting d into the u and v components evenly.
  • VectorData(u::AbstractVector,v::AbstractVector) constructs a wrapper for the vector components data u and v.
  • VectorData(n::Int) constructs a wrapper with zeros of length n for both components.
  • VectorData(x::PointData) constructs a wrapper for zero components of the same length as that wrapped by x.
  • VectorData(n::Int,dtype=ComplexF64) constructs a wrapper with complex-valued zeros of length n for both components.

Example

julia> f = VectorData(10,dtype=ComplexF64);

julia> f.v[1:5] = 1:5;

julia> f
10 points of vector-valued Complex{Float64} data
10×2 Array{Complex{Float64},2}:
 0.0+0.0im  1.0+0.0im
 0.0+0.0im  2.0+0.0im
 0.0+0.0im  3.0+0.0im
 0.0+0.0im  4.0+0.0im
 0.0+0.0im  5.0+0.0im
 0.0+0.0im  0.0+0.0im
 0.0+0.0im  0.0+0.0im
 0.0+0.0im  0.0+0.0im
 0.0+0.0im  0.0+0.0im
 0.0+0.0im  0.0+0.0im

julia> f[7] = 1im; f[18] = 0.2;

julia> f
10 points of vector-valued Complex{Float64} data
10×2 Array{Complex{Float64},2}:
 0.0+0.0im  1.0+0.0im
 0.0+0.0im  2.0+0.0im
 0.0+0.0im  3.0+0.0im
 0.0+0.0im  4.0+0.0im
 0.0+0.0im  5.0+0.0im
 0.0+0.0im  0.0+0.0im
 0.0+1.0im  0.0+0.0im
 0.0+0.0im  0.2+0.0im
 0.0+0.0im  0.0+0.0im
 0.0+0.0im  0.0+0.0im
source
CartesianGrids.XEdgesType
XEdges

XEdges is a wrapper for scalar-valued data that lie at the centers of either dual cells or primary cells. A XEdges type can be accessed by indexing like any other array, and allows the use of size, similar, zero functions.

Constructors

  • XEdges(C,dims) creates a field of zeros in cells of type C (where C is either Dual or Primal), on a grid of dimensions dims (a tuple). Note that dims represent the number of dual cells on the grid, even if C is Primal.
  • XEdges(C,w) performs the same construction, but uses existing field data w of GridData type to determine the size of the grid.
  • Adding the dtype= keyword allows the data type of the field data to be changed. The default is Float64, but can be changed to, e.g., ComplexF64
source
CartesianGrids.YEdgesType
YEdges

YEdges is a wrapper for scalar-valued data that lie at the centers of either dual cells or primary cells. A YEdges type can be accessed by indexing like any other array, and allows the use of size, similar, zero functions.

Constructors

  • YEdges(C,dims) creates a field of zeros in cells of type C (where C is either Dual or Primal), on a grid of dimensions dims (a tuple). Note that dims represent the number of dual cells on the grid, even if C is Primal.
  • YEdges(C,w) performs the same construction, but uses existing field data w of GridData type to determine the size of the grid.
  • Adding the dtype= keyword allows the data type of the field data to be changed. The default is Float64, but can be changed to, e.g., ComplexF64
source
Base.:*Method
(*)(p::VectorData,q::VectorData) -> TensorData

Calculate the element by element tensor product between vectors p and q, returning TensorData type.

source
Base.:+Method
(+)(X::VectorData,a::Tuple{T,T}) where {T<:Number} -> VectorData
(-)(X::VectorData,a::Tuple{T,T}) where {T<:Number} -> VectorData

Adds or subtracts the tuple a component by component to each element of X. All data in a are converted to Float64. Can also switch the arguments.

Example

julia> f = VectorData(5);

julia> f + (2,3)
5 points of vector-valued Float64 data
5×2 Array{Float64,2}:
 2.0  3.0
 2.0  3.0
 2.0  3.0
 2.0  3.0
 2.0  3.0
source
Base.expMethod
exp(L::Laplacian,a[,Nodes(Dual)][;nthreads=L.conv.nthreads])

Create the integrating factor exp(L*a). The default size of the operator is the one appropriate for dual nodes; another size can be specified by supplying grid data in the optional third argument. Note that, if L contains a factor, it scales the exponent with this factor. The number of threads used by the resulting operator can be set by the nthreads optional argument; by default, it takes this number from L.

source
Base.lengthMethod
length(g::PhysicalGrid,d::Int) -> Int

Return the total number of cells in grid g.

source
Base.sizeMethod
size(g::PhysicalGrid,d::Int) -> Int

Return the number of cells in direction d in grid g.

source
Base.sizeMethod
size(g::PhysicalGrid) -> Tuple

Return a tuple of the number of cells in all directions in grid g.

source
Base.transposeMethod
transpose(p::TensorData) -> TensorData

Element-by-element transpose of TensorData p

source
Base.transposeMethod
transpose(p::EdgeGradient) -> EdgeGradient

Element-by-element transpose of EdgeGradient p

source
Base.zeroMethod
zero(p::GridData)

Return data of the same type as p, filled with zeros.

source
Base.zeroMethod
zero(p::PointData)

Return data of the same type as p, filled with zeros.

source
CartesianGrids.convective_derivative!Method
convective_derivative!(out,q)

Compute the convective derivative of q in the form $u\cdot\nabla u$ and put the result into out. Note that the result is not scaled by any grid spacing.

source
CartesianGrids.convective_derivative_rot!Method
convective_derivative_rot!(out,q)

Compute the rotational form of the convective derivative of q in the form $\frac{1}{2}\nabla|u|^2-u\times(\nabla\times u)$ and put the result into out. Note that the result is not scaled by any grid spacing.

source
CartesianGrids.coordinatesFunction
coordinates(w::GridData;[dx=1.0],[I0=(1,1)])

Return a tuple of the ranges of the physical coordinates in each direction for grid data w. If w is of Nodes type, then it returns a tuple of the form xg,yg. If w is of Edges or NodePair type, then it returns a tuple of the form xgu,ygu,xgv,ygv.

The optional keyword argument dx sets the grid spacing; its default is 1.0. The optional keyword I0 accepts a tuple of integers to set the index pair of the primal nodes that coincide with the origin. The default is (1,1).

Example

julia> w = Nodes(Dual,(12,22));

julia> xg, yg = coordinates(w,dx=0.1)
(-0.05:0.1:1.05, -0.05:0.1:2.0500000000000003)
source
CartesianGrids.curl!Method
curl!(q::Edges{Primal},w::Nodes{Dual})

Evaluate the discrete curl of w and return it as q.

Example

julia> w = Nodes(Dual,(8,6));

julia> w[3,4] = 1.0;

julia> q = Edges(Primal,w);

julia> curl!(q,w)
Edges{Primal,8,6,Float64} data
u (in grid orientation)
5×8 Array{Float64,2}:
 0.0  0.0   0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  -1.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0   1.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0  0.0
v (in grid orientation)
6×7 Array{Float64,2}:
 0.0   0.0  0.0  0.0  0.0  0.0  0.0
 0.0   0.0  0.0  0.0  0.0  0.0  0.0
 0.0  -1.0  1.0  0.0  0.0  0.0  0.0
 0.0   0.0  0.0  0.0  0.0  0.0  0.0
 0.0   0.0  0.0  0.0  0.0  0.0  0.0
 0.0   0.0  0.0  0.0  0.0  0.0  0.0
source
CartesianGrids.curlMethod
curl(w::Nodes{Dual}) --> Edges{Primal}

Evaluate the discrete curl of w. Another way to perform this operation is to construct a Curl object and apply it with *.

Example

julia> C = Curl();

julia> w = Nodes(Dual,(8,6));

julia> w[3,4] = 1.0;

julia> C*w
Edges{Primal,8,6,Float64} data
u (in grid orientation)
5×8 Array{Float64,2}:
 0.0  0.0   0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  -1.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0   1.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0  0.0
v (in grid orientation)
6×7 Array{Float64,2}:
 0.0   0.0  0.0  0.0  0.0  0.0  0.0
 0.0   0.0  0.0  0.0  0.0  0.0  0.0
 0.0  -1.0  1.0  0.0  0.0  0.0  0.0
 0.0   0.0  0.0  0.0  0.0  0.0  0.0
 0.0   0.0  0.0  0.0  0.0  0.0  0.0
 0.0   0.0  0.0  0.0  0.0  0.0  0.0
source
CartesianGrids.curl_cross!Method
curl_cross!(out,a,b)

Compute the curl of the cross product of a and b, $\nabla\times(a\times b)$, and put the result into out. Note that the result is not scaled by any grid spacing.

source
CartesianGrids.diff!Function
diff!(out:ScalarGridData,in::ScalarGridData) -> ScalarGridData

Return the 1-d central finite difference of scalar grid data in in scalar grid data out. Either in or out must be edge component data and the other must be node data. The direction of differencing is determined by the relationship of in and out. For example, if in is dual nodes (cell by cell) and out is primal x-edge components (cell by edge), then the differencing takes place in the y direction, since they are different types in this direction.

source
CartesianGrids.directional_derivative!Method
directional_derivative!(out,f,q)

Compute the directional derivative of f in the direction of q, $q\cdot\nabla f$, and put the result into out. Note that the result is not scaled by any grid spacing.

source
CartesianGrids.directional_derivative_conserve!Method
directional_derivative_conserve!(out,f,q)

Compute the conservative form of the directional derivative of f in the direction of q, $\nabla\cdot(qf)$, and put the result into out. Note that the result is not scaled by any grid spacing. This form is only appropriate if q is divergence free.

source
CartesianGrids.divergence!Method
divergence!(w::Edges,q::EdgeGradient)

Evaluate the discrete divergence of edge gradient tensor data q and return it as data w. Note that q can be either primal/dual or dual/primal tensor data, and w must be, respectively, primal or edges type.

source
CartesianGrids.divergence!Method
divergence!(w::Nodes,q::Edges)

Evaluate the discrete divergence of edge data q and return it as nodal data w. Note that q can be either primal or dual edge data, but w must be of the same cell type.

Example

julia> q = Edges(Primal,(8,6));

julia> q.u[3,2] = 1.0;

julia> w = Nodes(Primal,(8,6));

julia> divergence!(w,q)
Nodes{Primal,8,6,Float64} data
Printing in grid orientation (lower left is (1,1))
5×7 Array{Float64,2}:
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  1.0  -1.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
source
CartesianGrids.divergence!Method
divergence!(w::XEdges/YEdges,q::NodePair)

Evaluate the discrete divergence of node pair data q and return it as data w. Note that q can be either primal/dual or dual/primal node data, and w must be, respectively, primal x-edges/dual y-edges or primal y-edges/dual x-edges type.

source
CartesianGrids.divergenceMethod
divergence(q::Edges) --> Nodes

Evaluate the discrete divergence of edge data q. Can also perform this operation by creating an object of Divergence type and applying it with *.

Example

julia> D = Divergence();

julia> q = Edges(Primal,(8,6));

julia> q.u[3,2] = 1.0;

julia> D*q
Nodes{Primal,8,6,Float64} data
Printing in grid orientation (lower left is (1,1))
5×7 Array{Float64,2}:
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  1.0  -1.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
source
CartesianGrids.exp!Method
exp!(L::Laplacian,a[,Nodes(Dual)][;nthreads=L.conv.nthreads])

Create the in-place version of the integrating factor exp(L*a).

source
CartesianGrids.grad!Method
grad!(d::EdgeGradient{Dual,Primal},q::Edges{Dual})

Evaluate the discrete gradient of dual edge data q and return it as edge gradient data d, where the diagonal entries of the gradient lie on dual nodes and the off-diagonal entries lie at primal nodes.

source
CartesianGrids.grad!Method
grad!(d::EdgeGradient{Primal,Dual},q::Edges{Primal})

Evaluate the discrete gradient of primal edge data q and return it as edge gradient data d, where the diagonal entries of the gradient lie on primal nodes and the off-diagonal entries lie at dual nodes.

source
CartesianGrids.grad!Method
grad!(q::Edges{Dual},w::Nodes{Dual})

Evaluate the discrete gradient of dual nodal data w and return it as dual edge data q.

source
CartesianGrids.grad!Method
grad!(q::Edges{Primal},w::Nodes{Primal})

Evaluate the discrete gradient of primal nodal data w and return it as primal edge data q.

Example

julia> w = Nodes(Primal,(8,6));

julia> w[3,4] = 1.0;

julia> q = Edges(Primal,(8,6));

julia> grad!(q,w)
Edges{Primal,8,6,Float64} data
u (in grid orientation)
5×8 Array{Float64,2}:
 0.0  0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0  1.0  -1.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0  0.0
v (in grid orientation)
6×7 Array{Float64,2}:
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0  -1.0  0.0  0.0  0.0  0.0
 0.0  0.0   1.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
source
CartesianGrids.gradMethod
grad(q::Edges{Primal/Dual}) --> EdgeGradient{Dual/Primal,Primal/Dual}

Evaluate the discrete gradient of primal or dual edge data q. Can also perform this operation by creating an object of Grad type and applying it with *.

source
CartesianGrids.gradMethod
grad(w::Nodes{Dual}) --> Edges{Dual}

Evaluate the discrete gradient of dual nodal data w. Can also perform this operation by creating an object of Grad type and applying it with *.

source
CartesianGrids.gradMethod
grad(w::Nodes{Primal}) --> Edges{Primal}

Evaluate the discrete gradient of primal nodal data w. Can also perform this operation by creating an object of Grad type and applying it with *.

Example

julia> w = Nodes(Primal,(8,6));

julia> w[3,4] = 1.0;

julia> G = Grad();

julia> G*w
Edges{Primal,8,6,Float64} data
u (in grid orientation)
5×8 Array{Float64,2}:
 0.0  0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0  1.0  -1.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0  0.0
v (in grid orientation)
6×7 Array{Float64,2}:
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0  -1.0  0.0  0.0  0.0  0.0
 0.0  0.0   1.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0   0.0  0.0  0.0  0.0  0.0
source
CartesianGrids.grid_interpolate!Function
grid_interpolate!(out:ScalarGridData,in::ScalarGridData) -> ScalarGridData

Return the 1-d symmetric interpolation of scalar grid data in onto scalar grid data out. Either in or out must be edge component data and the other must be node data. The direction of interpolation is determined by the relationship of in and out. For example, if in is dual nodes (cell by cell) and out is primal x-edge components (cell by edge), then the interpolation takes place in the y direction, since they are different types in this direction. (In the x direction, they are of the same type (cell), so there is no interpolation in that direction.)

source
CartesianGrids.grid_interpolate!Method
grid_interpolate!((wx::Nodes,wy::Nodes),q::Edges)

Interpolate the edge data q (of either dual or primal type) to the dual or primal nodes, and return the result in wx and wy. wx holds the shifted q.u data and wy the shifted q.v data.

Example

julia> q = Edges(Primal,(8,6));

julia> q.u[3,2] = 1.0;

julia> wx = Nodes(Dual,(8,6)); wy = deepcopy(wx);

julia> grid_interpolate!((wx,wy),q);

julia> wx
Nodes{Dual,8,6,Float64} data
Printing in grid orientation (lower left is (1,1))
6×8 Array{Float64,2}:
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.5  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.5  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0

julia> wy
Nodes{Dual,8,6,Float64} data
Printing in grid orientation (lower left is (1,1))
6×8 Array{Float64,2}:
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
source
CartesianGrids.grid_interpolate!Method
grid_interpolate!(dq::EdgeGradient{Primal/Dual},q::Edges{Primal/Dual})

Interpolate the primal (dual) edge data q to primal (dual) tensor positions and hold it in dq.

source
CartesianGrids.grid_interpolate!Method
grid_interpolate!(q::Edges{Dual},w::Nodes{Dual})

Interpolate the dual nodal data w to the edges of the dual cells, and return the result in q.

Example

julia> w = Nodes(Dual,(8,6));

julia> w[3,4] = 1.0;

julia> q = Edges(Dual,w);

julia> grid_interpolate!(q,w)
Edges{Dual,8,6,Float64} data
u (in grid orientation)
6×7 Array{Float64,2}:
 0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.5  0.5  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0
v (in grid orientation)
5×8 Array{Float64,2}:
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.5  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.5  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
source
CartesianGrids.grid_interpolate!Method
grid_interpolate!(q::Edges{Primal/Dual},dq::EdgeGradient{Primal/Dual})

Interpolate the primal (dual) tensor data dq to primal (dual) edge positions and hold it in q.

source
CartesianGrids.grid_interpolate!Method
grid_interpolate!(q::Edges{Primal},w::Nodes{Primal})

Interpolate the primal nodal data w to the edges of the primal cells, and return the result in q.

source
CartesianGrids.grid_interpolate!Method
grid_interpolate!(v::Nodes{Dual/Primal},q::Nodes{Primal/Dual})

Interpolate the primal (resp. dual) node data q to the edges of the dual (resp. primal) nodes, and return the result in v.

source
CartesianGrids.gridwise_dot!Method

gridwise_dot!(udotv::Nodes{Primal/Dual},u::Edges{Primal/Dual},v::Edges{Primal/Dual})

Calculate the in-placed dot product of vector grid data u and v, placing the result on the cell centers or nodes in udotv.

source
CartesianGrids.gridwise_dotMethod

gridwise_dot(u::Edges{Primal/Dual},v::Edges{Primal/Dual})

Calculate the dot product of vector grid data u and v, placing the result on the cell centers (if data are Edges{Primal}) or nodes (if data are Edges{Dual}).

source
CartesianGrids.helmholtz!Method
helmholtz!(v,w,α)

Evaluate the discrete Helmholtz operator (iα - L) on w and return it as v. The data w can be of type dual/primary nodes or edges; v must be of the same type. However, both have to be of complex data type.

source
CartesianGrids.helmholtzMethod
helmholtz(w,α)

Evaluate the discrete Helmholtz operator (iα - L) on w. The data w can be of complex type dual/primary nodes or edges. The returned result is of the same type as w.

Example

source
CartesianGrids.interpolatable_fieldMethod
interpolatable_field(x,y,f::ScalarGridData)

Generates an interpolatable version of grid data f, based on coordinates in x and y (which should be in range form). The output can be called as a function with coordinate pairs as arguments.

source
CartesianGrids.interpolatable_fieldMethod
interpolatable_field(xu,yu,xv,yv,f::EdgeData/NodePair)

Generates an interpolatable version of grid data f, based on coordinates in xu, yu, xv, yv (which should be in range form). The output can be called as a function with coordinate pairs as arguments.

source
CartesianGrids.interpolatable_fieldMethod
interpolatable_field(f::GridData,g::PhysicalGrid)

Generates an interpolatable version of grid data f, based on grid g. The output can be called as a function with coordinate pairs as arguments.

source
CartesianGrids.laplacian!Method
laplacian!(v,w)

Evaluate the discrete Laplacian of w and return it as v. The data w can be of type dual/primary nodes or edge components or edges; v must be of the same type.

Example

julia> w = Nodes(Dual,(8,6));

julia> v = deepcopy(w);

julia> w[4,3] = 1.0;

julia> laplacian!(v,w)
Nodes{Dual,8,6,Float64} data
Printing in grid orientation (lower left is (1,1))
6×8 Array{Float64,2}:
 0.0  0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0   1.0  0.0  0.0  0.0  0.0
 0.0  0.0  1.0  -4.0  1.0  0.0  0.0  0.0
 0.0  0.0  0.0   1.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0  0.0
source
CartesianGrids.laplacianMethod
laplacian(w)

Evaluate the discrete Laplacian of w. The data w can be of type dual/primary nodes or edges. The returned result is of the same type as w.

Example

julia> q = Edges(Primal,(8,6));

julia> q.u[2,2] = 1.0;

julia> laplacian(q)
Edges{Primal,8,6,Float64} data
u (in grid orientation)
5×8 Array{Float64,2}:
 0.0   0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0   0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0   1.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  -4.0  1.0  0.0  0.0  0.0  0.0  0.0
 0.0   0.0  0.0  0.0  0.0  0.0  0.0  0.0
v (in grid orientation)
6×7 Array{Float64,2}:
 0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0
source
CartesianGrids.laplacian_symm!Method
laplacian_symm!(v,w)

Evaluate the symmetric 5-point discrete Laplacian of w and return it as v. The data w can be of type dual nodes only for now. This symmetric Laplacian also evaluates the partial Laplacians (using only available stencil data) on the ghost nodes.

source
CartesianGrids.limitsMethod
limits(g::PhysicalGrid,d::Int) -> Tuple

Return the minimum and maximum physical dimensions in direction d for grid g.

source
CartesianGrids.mag!Method
mag!(magu::Nodes{Primal/Dual},u::Edges{Primal/Dual})

Calculate the in-placed magnitude of vector grid data u, placing the result on the cell centers in magu.

source
CartesianGrids.magMethod

mag(u::Edges{Primal/Dual}) -> Nodes{Primal/Dual}

Calculate the magnitude of vector grid data u, placing the result on the cell centers.

source
CartesianGrids.magsq!Method

magsq!(magusq::Nodes{Primal/Dual},u::Edges{Primal/Dual})

Calculate the in-placed squared magnitude of vector grid data u, placing the result on the cell centers in magusq.

source
CartesianGrids.magsqMethod

magsq(u::Edges{Primal/Dual}) -> Nodes{Primal/Dual}

Calculate the squared magnitude of vector grid data u, placing the result on the cell centers.

source
CartesianGrids.optimize_gridsizeMethod
optimize_gridsize(nx0,ny0[;region_size=4,nthreads_max=length(cpu_info()),nsamp=1])

Given a nominal grid size (nx0 x ny0), determine the optimal grid size that minimizes the compute time. Optional arguments are the maximum number of threads (if multithreading is allowed) and the number of samples to take of the cpu time for each trial. Returns optimal nx, ny, and the corresponding CPU time.

source
CartesianGrids.originMethod
origin(g::PhysicalGrid) -> Tuple{Int,Int}

Return a tuple of the indices of the primal node that corresponds to the physical origin of the coordinate system used by g. Note that these indices need not lie inside the range of indices occupied by the grid. For example, if the range of physical coordinates occupied by the grid is (1.0,3.0) x (2.0,4.0), then the origin is not inside the grid.

source
CartesianGrids.plan_helmholtzFunction
plan_helmholtz(dims::Tuple,α::Number,[with_inverse=false],[fftw_flags=FFTW.ESTIMATE],
                      [factor=1.0],[dx=1.0])

Constructor to set up an operator for evaluating the discrete Helmholtz operator on complex dual or primal nodal data of dimension dims. If the optional keyword with_inverse is set to true, then it also sets up the inverse Helmholtz operator (the lattice Green's function, LGF). These can then be applied, respectively, with * and \ operations on data of the appropriate size. The optional parameter factor is a scalar used to multiply the result of the operator and divide the inverse. The optional parameter dx is used in adjusting the uniform value of the LGF to match the behavior of the continuous analog at large distances; this is set to 1.0 by default.

Instead of the first argument, one can also supply w::Nodes to specify the size of the domain.

Example

source
CartesianGrids.plan_helmholtz!Function
plan_helmholtz!(dims::Tuple,α::Number,[with_inverse=false],[fftw_flags=FFTW.ESTIMATE],
                      [factor=1.0],[dx=1.0])

Same as plan_helmholtz, but forms an operator that works in-place on the data it operates on.

source
CartesianGrids.plan_implicit_diffusionFunction
plan_implicit_diffusion(a::Real,dims::Tuple,[fftw_flags=FFTW.ESTIMATE])

Constructor to set up the forward and inverse operators of I - a*L where L is the discrete Laplacian (not scaled by grid spacing) and a is a real-valued parameter. This can then be applied with the * or `` operation on data of the appropriate size.

The dims argument can be replaced with data of type ScalarGridData to specify the size of the domain.

source
CartesianGrids.plan_implicit_diffusion!Function
plan_implicit_diffusion!(a::Real,dims::Tuple,[fftw_flags=FFTW.ESTIMATE][,nthreads=length(Sys.cpu_info())])

Same as plan_implicit_diffusion, but the resulting operator performs an in-place operation on data. The number of threads threads defaults to the number of logical CPU cores on the system.

source
CartesianGrids.plan_intfactFunction
plan_intfact(a::Real,dims::Tuple,[fftw_flags=FFTW.ESTIMATE])

Constructor to set up an operator for evaluating the integrating factor with real-valued parameter a. This can then be applied with the * or `` operation on data of the appropriate size.

Note that a can be positive, negative, or zero. However, if a is negative, then only the $operation is actually correct; the `*` operation merely returns the identity to avoid excessive (and noisy) calculation. Similarly, if `a` is positive, the$ operation returns the identity. Thus, these operations are not inverses of one another. If a is zero, both operations return the identity.

The dims argument can be replaced with data of type ScalarGridData to specify the size of the domain.

Example

julia> w = Nodes(Dual,(6,6));

julia> w[4,4] = 1.0;

julia> E = plan_intfact(1.0,(6,6))
Integrating factor with parameter 1.0 on a (nx = 6, ny = 6) grid

julia> E*w
Nodes{Dual,6,6,Float64} data
Printing in grid orientation (lower left is (1,1))
6×6 Array{Float64,2}:
 0.00268447   0.00869352  0.0200715   0.028765    0.0200715   0.00869352
 0.00619787   0.0200715   0.0463409   0.0664124   0.0463409   0.0200715
 0.00888233   0.028765    0.0664124   0.0951774   0.0664124   0.028765
 0.00619787   0.0200715   0.0463409   0.0664124   0.0463409   0.0200715
 0.00268447   0.00869352  0.0200715   0.028765    0.0200715   0.00869352
 0.000828935  0.00268447  0.00619787  0.00888233  0.00619787  0.00268447
source
CartesianGrids.plan_intfact!Function
plan_intfact!(a::Real,dims::Tuple,[fftw_flags=FFTW.ESTIMATE][,nthreads=length(Sys.cpu_info())])

Same as plan_intfact, but the resulting operator performs an in-place operation on data. The number of threads threads defaults to the number of logical CPU cores on the system.

source
CartesianGrids.plan_laplacianFunction
plan_laplacian(dims::Tuple,[with_inverse=false],[fftw_flags=FFTW.ESTIMATE],
                      [factor=1.0],[dx=1.0],[dtype=Float64])

Constructor to set up an operator for evaluating the discrete Laplacian on dual or primal nodal data of dimension dims. If the optional keyword with_inverse is set to true, then it also sets up the inverse Laplacian (the lattice Green's function, LGF). These can then be applied, respectively, with * and \ operations on data of the appropriate size. The optional parameter factor is a scalar used to multiply the result of the operator and divide the inverse. The optional parameter dx is used in adjusting the uniform value of the LGF to match the behavior of the continuous analog at large distances; this is set to 1.0 by default. The type of data on which to act is floating point by default, but can also be ComplexF64. This is specified with the optional parameter dtype

Instead of the first argument, one can also supply w::Nodes to specify the size of the domain.

Example

julia> w = Nodes(Dual,(5,5));

julia> w[3,3] = 1.0;

julia> L = plan_laplacian(5,5;with_inverse=true)
Discrete Laplacian (and inverse) on a (nx = 5, ny = 5) grid acting on Float64 data with spacing 1.0

julia> s = L\w
Nodes{Dual,5,5,Float64} data
Printing in grid orientation (lower left is (1,1))
5×5 Array{Float64,2}:
 0.16707    0.129276     0.106037     0.129276    0.16707
 0.129276   0.0609665   -0.00734343   0.0609665   0.129276
 0.106037  -0.00734343  -0.257343    -0.00734343  0.106037
 0.129276   0.0609665   -0.00734343   0.0609665   0.129276
 0.16707    0.129276     0.106037     0.129276    0.16707

julia> L*s ≈ w
true
source
CartesianGrids.plan_laplacian!Function
plan_laplacian!(dims::Tuple,[with_inverse=false],[fftw_flags=FFTW.ESTIMATE],
                      [factor=1.0][,nthreads=length(Sys.cpu_info())])

Same as plan_laplacian, but operates in-place on data. The number of threads threads defaults to the number of logical CPU cores on the system.

source
CartesianGrids.pointwise_cross!Method
pointwise_cross!(C::ScalarData,A::VectorData,B::VectorData) -> ScalarData

Compute the cross product between the vector point data A and B and return the result as scalar data C (treated as an out-of-plane component of a vector).

source
CartesianGrids.pointwise_cross!Method
pointwise_cross!(C::VectorData,A::ScalarData/VectorData,B::VectorData/ScalarData) -> VectorData

Compute the cross product between the point data A and B, one of which is scalar data and treated as an out-of-plane component of a vector, while the other is in-plane vector data, and return the result as vector data C.

source
CartesianGrids.pointwise_crossMethod
pointwise_cross(a::Number/ScalarData,A::VectorData) -> VectorData

Compute the cross product between the scalar a (treated as an out-of-plane component of a vector) and the planar vector data A.

source
CartesianGrids.pointwise_crossMethod
pointwise_cross(A::VectorData/ScalarData,B::ScalarData/VectorData) -> VectorData

Compute the cross product between the point data A and B, one of which is scalar data and treated as an out-of-plane component of a vector, while the other is in-plane vector data, and return the result as vector data

source
CartesianGrids.pointwise_crossMethod
pointwise_cross(A::VectorData,B::VectorData) -> ScalarData

Compute the cross product between the vector point data A and B and return the result as scalar data (treated as an out-of-plane component of a vector).

source
CartesianGrids.pointwise_dot!Method
pointwise_dot!(C::ScalarData,A::VectorData,B::VectorData) -> ScalarData

Compute the element by element dot product between A and B and return the result in C.

source
CartesianGrids.pointwise_dot!Method
pointwise_dot!(C::VectorData,A::TensorData/VectorData,B::VectorData/TensorData) -> VectorData

Compute the element by element dot product between A and B, where one is TensorData and the other is VectorData and return the result as VectorData in C.

source
CartesianGrids.pointwise_dotMethod
pointwise_dot(A::TensorData/VectorData,B::VectorData/TensorData) -> VectorData

Compute the element by element dot product between A and B, where one is TensorData and the other is VectorData and return the result as VectorData.

source
CartesianGrids.pointwise_dotMethod
pointwise_dot(A::VectorData,B::VectorData) -> ScalarData

Compute the element by element dot product between A and B and return the result as ScalarData.

source
CartesianGrids.pointwise_tensorproduct!Method
pointwise_tensorproduct!(pq::TensorData,p::VectorData,q::VectorData) -> TensorData

Calculate the element by element tensor product between vectors p and q, returning data in pq

source
CartesianGrids.productFunction
product(p::PointData,q::PointData) -> PointData
(∘)(p::PointData,q::PointData) -> PointData

Compute the Hadamard (i.e. element by element) product of point data data p and q. Works similarly to product!.

source
CartesianGrids.product!Method
product!(out::PointData,p::PointData,q::PointData)

Compute the Hadamard (i.e. element by element) product of point data data p and q and return the result in out. Note that p and q can be of mixed type (scalar, vector, tensor), as long as one of them is a scalar. Also, out must have element type that is consistent with the promoted type of p and q.

Example

julia> fcs = ScalarData(5,dtype=ComplexF64);

julia> fill!(fcs,2im)
5 points of scalar-valued Complex{Float64} data
5-element Array{Complex{Float64},1}:
 0.0 + 2.0im
 0.0 + 2.0im
 0.0 + 2.0im
 0.0 + 2.0im
 0.0 + 2.0im

julia> frt = TensorData(fcs,dtype=Float64);

julia> fill!(frt,1.0)
5 points of tensor-valued Float64 data dudx, dudy, dvdx, dvdy
5×4 Array{Float64,2}:
 1.0  1.0  1.0  1.0
 1.0  1.0  1.0  1.0
 1.0  1.0  1.0  1.0
 1.0  1.0  1.0  1.0
 1.0  1.0  1.0  1.0

julia> out = similar(frt,element_type=ComplexF64);

julia> product!(out,frt,fcs)
5 points of tensor-valued Complex{Float64} data dudx, dudy, dvdx, dvdy
5×4 Array{Complex{Float64},2}:
 0.0+2.0im  0.0+2.0im  0.0+2.0im  0.0+2.0im
 0.0+2.0im  0.0+2.0im  0.0+2.0im  0.0+2.0im
 0.0+2.0im  0.0+2.0im  0.0+2.0im  0.0+2.0im
 0.0+2.0im  0.0+2.0im  0.0+2.0im  0.0+2.0im
 0.0+2.0im  0.0+2.0im  0.0+2.0im  0.0+2.0im
source
CartesianGrids.product!Method
product!(out::GridData,p::GridData,q::GridData)

Compute the Hadamard (i.e. element by element) product of grid data p and q (of the same type) and return the result in out.

Example

julia> q = Edges(Dual,(8,6));

julia> out = p = deepcopy(q);

julia> q.u[3,2] = 0.3;

julia> p.u[3,2] = 0.2;

julia> product!(out,p,q)
Edges{Dual,8,6,Float64} data
u (in grid orientation)
6×7 Array{Float64,2}:
 0.0  0.0  0.0   0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0
 0.0  0.0  0.06  0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0
v (in grid orientation)
5×8 Array{Float64,2}:
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
source
CartesianGrids.productMethod
product(p::Edges/Nodes,q::Edges/Nodes) --> Edges/Nodes

Compute the Hadamard product of edge or nodal (primal or dual) data p and q and return the result. This operation can also be carried out with the operator:

Example

julia> q = Edges(Dual,(8,6));

julia> p = deepcopy(q);

julia> q.u[3,2] = 0.3;

julia> p.u[3,2] = 0.2;

julia> p∘q
Edges{Dual,8,6,Float64} data
u (in grid orientation)
6×7 Array{Float64,2}:
 0.0  0.0  0.0   0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0
 0.0  0.0  0.06  0.0  0.0  0.0  0.0
 0.0  0.0  0.0   0.0  0.0  0.0  0.0
v (in grid orientation)
5×8 Array{Float64,2}:
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
 0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0
source
CartesianGrids.tensorproduct!Method
tensorproduct!(q::EdgeGradient,u::Edges,v::Edges,ut::EdgeGradient,vt::EdgeGradient)

In-place tensor product of u and v, with the result returned in q. The ut and vt are supplied as temporary storage.

source
CartesianGrids.test_cputimeMethod
test_cputime(nx,ny,nthreads_max;[nsamp=1][,testtype=:laplacian][,kwargs]) -> Float64, Float64

Evaluate a sample FFT-based problem with the given size of grid nx x ny and the provided maximum number of threads nthreads_max. Returns the mean and standard deviation of the computational time. The optional argument nsamp can be used to perform an average timing over multiple samples. The test type is specified with the testtype optional argument. The default test is inversion of a Laplacian (:laplacian). Other options are :intfact and :helmholtz.

source
LinearAlgebra.dotMethod
dot(p1::Edges{Dual/Primal},p2::Edges{Dual/Primal}) -> Real

Computes the inner product between two sets of edge data on the same grid.

source
LinearAlgebra.dotMethod
dot(p1::EdgeGradient{Dual/Primal},p2::EdgeGradient{Dual/Primal}) -> Real

Computes the inner product between two sets of edge gradient data on the same grid.

source
LinearAlgebra.dotMethod
dot(p1::Nodes{Dual},p2::Nodes{Dual}) -> Real

Computes the inner product between two sets of dual node data on the same grid.

source
LinearAlgebra.dotMethod
dot(p1::Nodes{Primal},p2::Nodes{Primal}) -> Real

Computes the inner product between two sets of primal node data on the same grid.

source
LinearAlgebra.dotMethod
dot(p1::XEdges{Dual},p2::XEdges{Dual}) -> Real

Computes the inner product between two sets of dual x-edge component data on the same grid.

source
LinearAlgebra.dotMethod
dot(p1::XEdges{Primal},p2::XEdges{Primal}) -> Real

Computes the inner product between two sets of primal x-edge component data on the same grid.

source
LinearAlgebra.dotMethod
dot(p1::YEdges{Dual},p2::YEdges{Dual}) -> Real

Computes the inner product between two sets of dual y-edge component data on the same grid.

source
LinearAlgebra.dotMethod
dot(p1::YEdges{Primal},p2::YEdges{Primal}) -> Real

Computes the inner product between two sets of primal y-edge component data on the same grid.

source
LinearAlgebra.dotMethod
dot(A::Tuple{T,T},B::TensorData) where {T<:Number} -> VectorData
⋅(A::Tuple{T,T},B::TensorData) where {T<:Number} -> VectorData

Computes the dot product between the tuple A and the elements of a tensor B on a set of points and returns vector data on the same set of points.

source
LinearAlgebra.dotMethod
dot(A::Tuple{T,T},B::VectorData) where {T<:Number} -> ScalarData
⋅(A::Tuple{T,T},B::VectorData) where {T<:Number} -> ScalarData

Computes the dot product between the tuple v and the elements of a tensor B on a set of points and returns scalar data on the same set of points.

source