List of all functions

Here a summary of all functions:

PETSc.DMGlobalVecType
DMGlobalVec(v::CVec, dm::AbstractDM)

Container for an PETSc vector we know is "global"

External Links

source
PETSc.DMLocalVecType
DMLocalVec(v::CVec, dm::AbstractDM)

Container for an PETSc vector we know is "local"

External Links

source
PETSc.MatSeqAIJType
MatSeqAIJ{T}

PETSc sparse array using AIJ format (also known as a compressed sparse row or CSR format).

Memory allocation is handled by PETSc.

source
PETSc.MatShellType
MatShell{T}(obj, m, n)

Create a m×n PETSc shell matrix object wrapping obj.

If obj is a Function, then the multiply action obj(y,x); otherwise it calls mul!(y, obj, x). This can be changed by defining PETSc._mul!.

source
PETSc.OptionsType
Options{PetscLib <: PetscLibType}(kw -> arg, ...)
Options(petsclib, kw -> arg, ...)

Create a PETSc options data structure for the petsclib.

For construction a set of keyword argment pairs should be given. If the option has no value it should be set to nothing or true. Setting an option to false will cause the option not to be set on the PETSc options table.

Examples

julia> using PETSc

julia> petsclib = PETSc.petsclibs[1];

julia> PETSc.initialize(petsclib)

julia> opt = PETSc.Options(
                         petsclib,
                         ksp_monitor = nothing,
                         ksp_view = true,
                         pc_type = "mg",
                         pc_mg_levels = 1,
                         false_opt = false,
                     )
#PETSc Option Table entries:
-ksp_monitor
-ksp_view
-pc_mg_levels 1
-pc_type mg
#End of PETSc Option Table entries


julia> opt["ksp_monitor"]
""

julia> opt["pc_type"]
"mg"

julia> opt["pc_type"] = "ilu"
"ilu"

julia> opt["pc_type"]
"ilu"

julia> opt["false_opt"]
ERROR: KeyError: key "bad_key" not found

julia> opt["bad_key"]
ERROR: KeyError: key "bad_key" not found

External Links

source
PETSc.PetscLibTypeType
PetscLibType{PetscScalar, PetscInt}(petsc_library)

A container for specific PETSc libraries.

All other containers for PETSc objects should be typed on this to ensure that dispatch is correct.

source
PETSc.SNESMethod
SNES{PetscScalar}(
    ::UnionPetscLib,
    comm::MPI.Comm; 
    snessetfromoptions = true,
    options...)

Initializes a SNES nonlinear solver object

source
PETSc.VecSeqType
VecSeq(v::Vector)

A standard, sequentially-stored serial PETSc vector, wrapping the Julia vector v.

This reuses the array v as storage, and so v should not be resize!-ed or otherwise have its length modified while the PETSc object exists.

This should only be need to be called for more advanced uses, for most simple usecases, users should be able to pass Vectors directly and have the wrapping performed automatically

External Links

source
PETSc.DMDACreate1dFunction
DMDACreate1d(
    ::PetscLib
    comm::MPI.Comm,
    boundary_type::DMBoundaryType,
    global_dim,
    dof_per_node,
    stencil_width,
    points_per_proc::Union{Nothing, Vector{PetscInt}};
    dmsetfromoptions=true,
    dmsetup=true,
    options...
)

Creates a 1-D distributed array with the options specified using keyword arguments.

If keyword argument dmsetfromoptions == true then setfromoptions! called. If keyword argument dmsetup == true then setup! is called.

External Links

source
PETSc.DMDACreate2dFunction
DMDACreate2d(
    ::PetscLib
    comm::MPI.Comm,
    boundary_type_x::DMBoundaryType,
    boundary_type_y::DMBoundaryType,
    stencil_type::DMDAStencilType,
    global_dim_x,
    global_dim_y,
    procs_x,
    procs_y,
    dof_per_node,
    stencil_width,
    points_per_proc_x::Union{Nothing, Vector{PetscInt}};
    points_per_proc_y::Union{Nothing, Vector{PetscInt}};
    dmsetfromoptions=true,
    dmsetup=true,
    options...
)

Creates a 2-D distributed array with the options specified using keyword arguments.

If keyword argument dmsetfromoptions == true then setfromoptions! called. If keyword argument dmsetup == true then setup! is called.

External Links

source
PETSc.DMDACreate3dFunction
DMDACreate3d(
    ::PetscLib
    comm::MPI.Comm,
    boundary_type_x::DMBoundaryType,
    boundary_type_y::DMBoundaryType,
    boundary_type_z::DMBoundaryType,
    stencil_type::DMDAStencilType,
    global_dim_x,
    global_dim_y,
    global_dim_z,
    procs_x,
    procs_y,
    procs_z,
    global_dim_z,
    dof_per_node,
    stencil_width,
    points_per_proc_x::Union{Nothing, Vector{PetscInt}};
    points_per_proc_y::Union{Nothing, Vector{PetscInt}};
    points_per_proc_z::Union{Nothing, Vector{PetscInt}};
    dmsetfromoptions=true,
    dmsetup=true,
    options...
)

Creates a 3-D distributed array with the options specified using keyword arguments.

If keyword argument dmsetfromoptions == true then setfromoptions! called. If keyword argument dmsetup == true then setup! is called.

External Links

source
PETSc.MatSetValuesStencil!Function
MatSetValuesStencil!(mat::AbstractMat{PetscScalar},
    rows::Vector{MatStencil{PetscInt}}, 
    cols::Vector{MatStencil{PetscInt}}, 
    vals::Vector{PetscScalar},
    mode;
    num_cols = length(col),
    num_rows = length(row)
  )

Insert the vals specified by rows and cols stencil indices into the mat. The optional arguments num_cosl and num_rows allow the limiting of the elements of the rows and cols vectors.

see PETSc manual

source
PETSc.decrefMethod
decref(obj)

Decrement the reference counter for obj.

In general we don't need to use this, as we can call destroy instead.

source
PETSc.finalizeFunction

finalize(petsclib)

Finalize the petsclib, if no petsclib is given then all PETSc.petsclibs will be finalized.

External Links

source
PETSc.getcornersMethod
getcorners(da::DMDA)

Returns a NamedTuple with the global indices (excluding ghost points) of the lower and upper corners as well as the size.

External Links

source
PETSc.getvalues!Method
getvalues!(
    vector::AbstractVec{PetscScalar},
    indices::Vector{PetscInt},
    vals::Vector{PetscScalar};
    num_vals = length(inds)
)

Get a set of values from the vector. Equivalent to one of the following

vals[1:num_vals] .= vector[indices[1:num_vals]]
Warning

indices should use 0-based indexing!

External Links

source
PETSc.increfMethod
incref(obj)

Increment the reference counter fo obj. This usually only needs to be called when accessing objects owned by other objects, e.g. via KSPGetPC.

source
PETSc.inttypeMethod
inttype(petsclib::PetscLibType)

return the int type for the associated petsclib

source
PETSc.map_unsafe_localarray!Method
map_unsafe_localarray!(f!, x::AbstractVec{T}; read=true, write=true)

Convert x to an Array{T} and apply the function f!.

Use read=false if the array is write-only; write=false if read-only.

Examples

```julia-repl julia> mapunsafelocalarray(x; write=true) do x @. x .*= 2 end

Note

Base.finalize should is automatically called on the array.

source
PETSc.ownershiprangeFunction
ownershiprange(vec::AbstractVec)

The range of indices owned by this processor, assuming that the vectors are laid out with the first n1 elements on the first processor, next n2 elements on the second, etc. For certain parallel layouts this range may not be well defined.

Note: unlike the C function, the range returned is inclusive (idx_first:idx_last)

External Links

source
PETSc.parse_optionsMethod
parse_options(args::Vector{String})

Parse the args vector into a NamedTuple that can be used as the options for the PETSc solvers.

julia --project file.jl -ksp_monitor -pc_type mg -ksp_view
source
PETSc.realtypeMethod
realtype(petsclib::PetscLibType)

return the real type for the associated petsclib

source
PETSc.scalartypeMethod
scalartype(petsclib::PetscLibType)

return the scalar type for the associated petsclib

source
PETSc.setuniformcoordinates!Function
setuniformcoordinates!(
    da::DMDA
    xyzmin::NTuple{N, Real},
    xyzmax::NTuple{N, Real},
) where {N}

Set uniform coordinates for the da using the lower and upper corners defined by the NTuples xyzmin and xyzmax. If N is less than the dimension of the da then the value of the trailing coordinates is set to 0.

External Links

source
PETSc.setvalues!Method
setvalues!(
    vector::AbstractVec{PetscScalar},
    indices::Vector{PetscInt},
    vals::Vector{PetscScalar},
    mode::InsertMode;
    num_vals = length(ind)
)

Insert a set of values into the vector. Equivalent to one of the following depending on the mode

vector[indices[1:num_vals]] .= vals[1:num_vals]
vector[indices[1:num_vals]] .+= vals[1:num_vals]
Warning

indices should use 0-based indexing!

External Links

source
PETSc.unsafe_localarrayFunction
unsafe_localarray(PetscScalar, ptr:CVec; read=true, write=true)
unsafe_localarray(ptr:AbstractVec; read=true, write=true)

Return an Array{PetscScalar} containing local portion of the PETSc data.

Use read=false if the array is write-only; write=false if read-only.

Note

Base.finalize should be called on the Array before the data can be used.

source
PETSc.update!Method
update!(
    global_vec::DMGlobalVec,
    local_vec::DMLocalVec,
    mode::InsertMode,
)

Updates global_vec from local_vec with insert mode

External Links

source
PETSc.update!Method
update!(
    global_ptr::CVec,
    local_vec::DMLocalVec,
    mode::InsertMode,
)

Updates pointer to global vec global_ptr from local_vec with insert mode

External Links

source
PETSc.update!Method
update!(
    local_vec::DMLocalVec,
    global_vec::DMGlobalVec,
    mode::InsertMode
)

Updates local_vec from a pointer to a global vec global_ptr with insert mode

External Links

source
PETSc.update!Method
update!(
    local_vec::DMLocalVec,
    global_ptr::CVec,
    mode::InsertMode,
)

Updates local_vec from pointer to global vec global_ptr with insert mode

External Links

source
PETSc.update!Method
update!(
    local_ptr::CVec,
    global_vec::DMGlobalVec
    mode::InsertMode,
)

Updates pointer to local vec local_ptr from global_vec with insert mode

External Links

source
PETSc.update_global2local!Function
update_global2local!(
    local_ptr::CVec,
    global_ptr::CVec,
    mode::InsertMode,
    dm::AbstractDM
)

Updates pointer to local vector from pointer to global vector with insert mode, assuming that both belong to the same dm

This is a low-level routine that is typically called by update!

source
PETSc.update_local2global!Function
update_local2global!(
    global_ptr::CVec,
    local_ptr::CVec,
    mode::InsertMode,
    dm::AbstractDM
)

Updates pointer of global_vec from pointer of local_vec with insert mode. Both vectors should belong to the same dm

This is a low-level routine that is typically called by update!

source
PETSc.viewMethod
view(dm::AbstractDM, viewer::Viewer=ViewerStdout(petsclib, getcomm(dm)))

view a dm with viewer

External Links

source
PETSc.withMethod
with(f, opts::Options)

Call f() with the Options opts set temporarily (in addition to any global options).

source