petsc4py.PETSc.Vec#

class petsc4py.PETSc.Vec#

Bases: Object

A vector object.

See also

Vec

Enumerations

Option

Vector assembly option.

Type

The vector type.

Methods Summary

abs()

Replace each entry (xₙ) in the vector by abs|xₙ|.

appendOptionsPrefix(prefix)

Append to the prefix used for searching for options in the database.

assemble()

Assemble the vector.

assemblyBegin()

Begin an assembling stage of the vector.

assemblyEnd()

Finish the assembling stage initiated with assemblyBegin.

attachDLPackInfo([vec, dltensor])

Attach tensor information from another vector or DLPack tensor.

axpby(alpha, beta, x)

Compute and store y = ɑ·x + β·y.

axpy(alpha, x)

Compute and store y = ɑ·x + y.

aypx(alpha, x)

Compute and store y = x + ɑ·y.

bindToCPU(flg)

Bind vector operations execution on the CPU.

boundToCPU()

Return whether the vector has been bound to the CPU.

chop(tol)

Set all vector entries less than some absolute tolerance to zero.

clearDLPackInfo()

Clear tensor information.

conjugate()

Conjugate the vector.

copy([result])

Return a copy of the vector.

create([comm])

Create a vector object.

createCUDAWithArrays([cpuarray, cudahandle, ...])

Create a Type.CUDA vector with optional arrays.

createGhost(ghosts, size[, bsize, comm])

Create a parallel vector with ghost padding on each processor.

createGhostWithArray(ghosts, array[, size, ...])

Create a parallel vector with ghost padding and provided arrays.

createHIPWithArrays([cpuarray, hiphandle, ...])

Create a Type.HIP vector with optional arrays.

createLocalVector()

Create a local vector.

createMPI(size[, bsize, comm])

Create a parallel Type.MPI vector.

createNest(vecs[, isets, comm])

Create a Type.NEST vector containing multiple nested subvectors.

createSeq(size[, bsize, comm])

Create a sequential Type.SEQ vector.

createShared(size[, bsize, comm])

Create a Type.SHARED vector that uses shared memory.

createViennaCLWithArrays([cpuarray, ...])

Create a Type.VIENNACL vector with optional arrays.

createWithArray(array[, size, bsize, comm])

Create a vector using a provided array.

createWithDLPack(dltensor[, size, bsize, comm])

Create a vector wrapping a DLPack object, sharing the same memory.

destroy()

Destroy the vector.

dot(vec)

Return the dot product with vec.

dotBegin(vec)

Begin computing the dot product.

dotEnd(vec)

Finish computing the dot product initiated with dotBegin.

dotNorm2(vec)

Return the dot product with vec and its squared norm.

duplicate([array])

Create a new vector with the same type, optionally with data.

equal(vec)

Return whether the vector is equal to another.

exp()

Replace each entry (xₙ) in the vector by exp(xₙ).

getArray([readonly])

Return local portion of the vector as an ndarray.

getBlockSize()

Return the block size of the vector.

getBuffer([readonly])

Return a buffered view of the local portion of the vector.

getCLContextHandle()

Return the OpenCL context associated with the vector.

getCLMemHandle([mode])

Return the OpenCL buffer associated with the vector.

getCLQueueHandle()

Return the OpenCL command queue associated with the vector.

getCUDAHandle([mode])

Return a pointer to the device buffer.

getDM()

Return the DM associated to the vector.

getHIPHandle([mode])

Return a pointer to the device buffer.

getLGMap()

Return the local-to-global mapping.

getLocalSize()

Return the local size of the vector.

getLocalVector(lvec[, readonly])

Maps the local portion of the vector into a local vector.

getNestSubVecs()

Return all the vectors contained in the nested vector.

getOffloadMask()

Return the offloading status of the vector.

getOptionsPrefix()

Return the prefix used for searching for options in the database.

getOwnershipRange()

Return the locally owned range of indices (start, end).

getOwnershipRanges()

Return the range of indices owned by each process.

getSize()

Return the global size of the vector.

getSizes()

Return the vector sizes.

getSubVector(iset[, subvec])

Return a subvector from given indices.

getType()

Return the type of the vector.

getValue(index)

Return a single value from the vector.

getValues(indices[, values])

Return values from certain locations in the vector.

getValuesStagStencil(indices[, values])

Not implemented.

ghostUpdate([addv, mode])

Update ghosted vector entries.

ghostUpdateBegin([addv, mode])

Begin updating ghosted vector entries.

ghostUpdateEnd([addv, mode])

Finish updating ghosted vector entries initiated with ghostUpdateBegin.

isaxpy(idx, alpha, x)

Add a scaled reduced-space vector to a subset of the vector.

isset(idx, alpha)

Set specific elements of the vector to the same value.

load(viewer)

Load a vector.

localForm()

Return a context manager for viewing ghost vectors in local form.

log()

Replace each entry in the vector by its natural logarithm.

mDot(vecs[, out])

Compute Xᴴ·y with X an array of vectors.

mDotBegin(vecs, out)

Starts a split phase multiple dot product computation.

mDotEnd(vecs, out)

Ends a split phase multiple dot product computation.

max()

Return the vector entry with maximum real part and its location.

maxPointwiseDivide(vec)

Return the maximum of the component-wise absolute value division.

maxpy(alphas, vecs)

Compute and store y = Σₙ(ɑₙ·Xₙ) + y with X an array of vectors.

min()

Return the vector entry with minimum real part and its location.

mtDot(vecs[, out])

Compute Xᵀ·y with X an array of vectors.

mtDotBegin(vecs, out)

Starts a split phase transpose multiple dot product computation.

mtDotEnd(vecs, out)

Ends a split phase transpose multiple dot product computation.

norm([norm_type])

Compute the vector norm.

normBegin([norm_type])

Begin computing the vector norm.

normEnd([norm_type])

Finish computations initiated with normBegin.

normalize()

Normalize the vector by its 2-norm.

permute(order[, invert])

Permute the vector in-place with a provided ordering.

placeArray(array)

Set the local portion of the vector to a provided array.

pointwiseDivide(x, y)

Compute and store the component-wise division of two vectors.

pointwiseMax(x, y)

Compute and store the component-wise maximum of two vectors.

pointwiseMaxAbs(x, y)

Compute and store the component-wise maximum absolute values.

pointwiseMin(x, y)

Compute and store the component-wise minimum of two vectors.

pointwiseMult(x, y)

Compute and store the component-wise multiplication of two vectors.

reciprocal()

Replace each entry in the vector by its reciprocal.

resetArray([force])

Reset the vector to use its default array.

restoreCLMemHandle()

Restore a pointer to the OpenCL buffer obtained with getCLMemHandle.

restoreCUDAHandle(handle[, mode])

Restore a pointer to the device buffer obtained with getCUDAHandle.

restoreHIPHandle(handle[, mode])

Restore a pointer to the device buffer obtained with getHIPHandle.

restoreLocalVector(lvec[, readonly])

Unmap a local access obtained with getLocalVector.

restoreSubVector(iset, subvec)

Restore a subvector extracted using getSubVector.

scale(alpha)

Scale all entries of the vector.

set(alpha)

Set all components of the vector to the same value.

setArray(array)

Set values for the local portion of the vector.

setBlockSize(bsize)

Set the block size of the vector.

setDM(dm)

Associate a DM to the vector.

setFromOptions()

Configure the vector from the options database.

setLGMap(lgmap)

Set the local-to-global mapping.

setMPIGhost(ghosts)

Set the ghost points for a ghosted vector.

setNestSubVecs(sx[, idxm])

Set the component vectors at specified indices in the nested vector.

setOption(option, flag)

Set option.

setOptionsPrefix(prefix)

Set the prefix used for searching for options in the database.

setRandom([random])

Set all components of the vector to random numbers.

setSizes(size[, bsize])

Set the local and global sizes of the vector.

setType(vec_type)

Set the vector type.

setUp()

Set up the internal data structures for using the vector.

setValue(index, value[, addv])

Insert or add a single value in the vector.

setValueLocal(index, value[, addv])

Insert or add a single value in the vector using a local numbering.

setValues(indices, values[, addv])

Insert or add multiple values in the vector.

setValuesBlocked(indices, values[, addv])

Insert or add blocks of values in the vector.

setValuesBlockedLocal(indices, values[, addv])

Insert or add blocks of values in the vector with a local numbering.

setValuesLocal(indices, values[, addv])

Insert or add multiple values in the vector with a local numbering.

setValuesStagStencil(indices, values[, addv])

Not implemented.

shift(alpha)

Shift all entries in the vector.

sqrtabs()

Replace each entry (xₙ) in the vector by √|xₙ|.

strideGather(field, vec[, addv])

Insert component values into a single-component vector.

strideMax(field)

Return the maximum of entries in a subvector.

strideMin(field)

Return the minimum of entries in a subvector.

strideNorm(field[, norm_type])

Return the norm of entries in a subvector.

strideScale(field, alpha)

Scale a component of the vector.

strideScatter(field, vec[, addv])

Scatter entries into a component of another vector.

strideSum(field)

Sum subvector entries.

sum()

Return the sum of all the entries of the vector.

swap(vec)

Swap the content of two vectors.

tDot(vec)

Return the indefinite dot product with vec.

tDotBegin(vec)

Begin computing the indefinite dot product.

tDotEnd(vec)

Finish computing the indefinite dot product initiated with tDotBegin.

toDLPack([mode])

Return a DLPack PyCapsule wrapping the vector data.

view([viewer])

Display the vector.

waxpy(alpha, x, y)

Compute and store w = ɑ·x + y.

zeroEntries()

Set all entries in the vector to zero.

Attributes Summary

array

Alias for array_w.

array_r

Read-only ndarray containing the local portion of the vector.

array_w

Writeable ndarray containing the local portion of the vector.

block_size

The block size.

buffer

Alias for buffer_w.

buffer_r

Read-only buffered view of the local portion of the vector.

buffer_w

Writeable buffered view of the local portion of the vector.

local_size

The local vector size.

owner_range

The locally owned range of indices in the form [low, high).

owner_ranges

The range of indices owned by each process.

size

The global vector size.

sizes

The local and global vector sizes.

Methods Documentation

abs()#

Replace each entry (xₙ) in the vector by abs|xₙ|.

Logically collective.

See also

VecAbs

Source code at petsc4py/PETSc/Vec.pyx:2293

Return type:

None

appendOptionsPrefix(prefix)#

Append to the prefix used for searching for options in the database.

Logically collective.

Source code at petsc4py/PETSc/Vec.pyx:1028

Parameters:

prefix (str | None)

Return type:

None

assemble()#

Assemble the vector.

Collective.

Source code at petsc4py/PETSc/Vec.pyx:3047

Return type:

None

assemblyBegin()#

Begin an assembling stage of the vector.

Collective.

Source code at petsc4py/PETSc/Vec.pyx:3023

Return type:

None

assemblyEnd()#

Finish the assembling stage initiated with assemblyBegin.

Collective.

Source code at petsc4py/PETSc/Vec.pyx:3035

Return type:

None

attachDLPackInfo(vec=None, dltensor=None)#

Attach tensor information from another vector or DLPack tensor.

Logically collective.

This tensor information is required when converting a Vec to a DLPack object.

Parameters:
  • vec (Vec | None) – Vector with attached tensor information. This is typically created by calling createWithDLPack.

  • dltensor – DLPack tensor. This will only be used if vec is None.

Return type:

Self

Notes

This operation does not copy any data from vec or dltensor.

Source code at petsc4py/PETSc/Vec.pyx:643

axpby(alpha, beta, x)#

Compute and store y = ɑ·x + β·y.

Logically collective.

Parameters:
  • alpha (Scalar) – First scale factor.

  • beta (Scalar) – Second scale factor.

  • x (Vec) – Input vector, must not be the current vector.

Return type:

None

See also

axpy, aypx, waxpy, VecAXPBY

Source code at petsc4py/PETSc/Vec.pyx:2524

axpy(alpha, x)#

Compute and store y = ɑ·x + y.

Logically collective.

Parameters:
  • alpha (Scalar) – Scale factor.

  • x (Vec) – Input vector.

Return type:

None

See also

isaxpy, VecAXPY

Source code at petsc4py/PETSc/Vec.pyx:2460

aypx(alpha, x)#

Compute and store y = x + ɑ·y.

Logically collective.

Parameters:
  • alpha (Scalar) – Scale factor.

  • x (Vec) – Input vector, must not be the current vector.

Return type:

None

See also

axpy, axpby, VecAYPX

Source code at petsc4py/PETSc/Vec.pyx:2504

bindToCPU(flg)#

Bind vector operations execution on the CPU.

Logically collective.

Source code at petsc4py/PETSc/Vec.pyx:1396

Parameters:

flg (bool)

Return type:

None

boundToCPU()#

Return whether the vector has been bound to the CPU.

Not collective.

Source code at petsc4py/PETSc/Vec.pyx:1409

Return type:

bool

chop(tol)#

Set all vector entries less than some absolute tolerance to zero.

Collective.

Parameters:

tol (float) – The absolute tolerance below which entries are set to zero.

Return type:

None

See also

VecFilter

Source code at petsc4py/PETSc/Vec.pyx:1732

clearDLPackInfo()#

Clear tensor information.

Logically collective.

Source code at petsc4py/PETSc/Vec.pyx:705

Return type:

Self

conjugate()#

Conjugate the vector.

Logically collective.

See also

VecConjugate

Source code at petsc4py/PETSc/Vec.pyx:2305

Return type:

None

copy(result=None)#

Return a copy of the vector.

Logically collective.

This operation copies vector entries to the new vector.

Parameters:

result (Vec | None) – Target vector for the copy. If None then a new vector is created internally.

Return type:

Vec

See also

duplicate, VecCopy

Source code at petsc4py/PETSc/Vec.pyx:1707

create(comm=None)#

Create a vector object.

Collective.

After creation the vector type can then be set with setType.

Parameters:

comm (Comm | None) – MPI communicator, defaults to Sys.getDefaultComm.

Return type:

Self

See also

destroy, VecCreate

Source code at petsc4py/PETSc/Vec.pyx:176

createCUDAWithArrays(cpuarray=None, cudahandle=None, size=None, bsize=None, comm=None)#

Create a Type.CUDA vector with optional arrays.

Collective.

Parameters:
  • cpuarray (Sequence[Scalar] | None) – Host array. Will be lazily allocated if not provided.

  • cudahandle (Any | None) – Address of the array on the GPU. Will be lazily allocated if not provided.

  • size (LayoutSizeSpec | None) – Vector size.

  • bsize (int | None) – Vector block size. If None, bsize = 1.

  • comm (Comm | None) – MPI communicator, defaults to Sys.getDefaultComm.

Return type:

Self

Source code at petsc4py/PETSc/Vec.pyx:370

createGhost(ghosts, size, bsize=None, comm=None)#

Create a parallel vector with ghost padding on each processor.

Collective.

Parameters:
Return type:

Self

Source code at petsc4py/PETSc/Vec.pyx:819

createGhostWithArray(ghosts, array, size=None, bsize=None, comm=None)#

Create a parallel vector with ghost padding and provided arrays.

Collective.

Parameters:
  • ghosts (Sequence[int]) – Global indices of ghost points.

  • array (Sequence[Scalar]) – Array to store the vector values. Must be at least as large as the local size of the vector (including ghost points).

  • size (LayoutSizeSpec | None) – Vector size.

  • bsize (int | None) – Vector block size. If None, bsize = 1.

  • comm (Comm | None) – MPI communicator, defaults to Sys.getDefaultComm.

Return type:

Self

Source code at petsc4py/PETSc/Vec.pyx:861

createHIPWithArrays(cpuarray=None, hiphandle=None, size=None, bsize=None, comm=None)#

Create a Type.HIP vector with optional arrays.

Collective.

Parameters:
  • cpuarray (Sequence[Scalar] | None) – Host array. Will be lazily allocated if not provided.

  • hiphandle (Any | None) – Address of the array on the GPU. Will be lazily allocated if not provided.

  • size (LayoutSizeSpec | None) – Vector size.

  • bsize (int | None) – Vector block size. If None, bsize = 1.

  • comm (Comm | None) – MPI communicator, defaults to Sys.getDefaultComm.

Return type:

Self

Source code at petsc4py/PETSc/Vec.pyx:428

createLocalVector()#

Create a local vector.

Not collective.

Returns:

The local vector.

Return type:

Vec

Source code at petsc4py/PETSc/Vec.pyx:1204

createMPI(size, bsize=None, comm=None)#

Create a parallel Type.MPI vector.

Collective.

Parameters:
Return type:

Self

Source code at petsc4py/PETSc/Vec.pyx:283

createNest(vecs, isets=None, comm=None)#

Create a Type.NEST vector containing multiple nested subvectors.

Collective.

Parameters:
Return type:

Self

See also

VecCreateNest

Source code at petsc4py/PETSc/Vec.pyx:952

createSeq(size, bsize=None, comm=None)#

Create a sequential Type.SEQ vector.

Collective.

Parameters:
Return type:

Self

Source code at petsc4py/PETSc/Vec.pyx:247

createShared(size, bsize=None, comm=None)#

Create a Type.SHARED vector that uses shared memory.

Collective.

Parameters:
Return type:

Self

See also

VecCreateShared

Source code at petsc4py/PETSc/Vec.pyx:918

createViennaCLWithArrays(cpuarray=None, viennaclvechandle=None, size=None, bsize=None, comm=None)#

Create a Type.VIENNACL vector with optional arrays.

Collective.

Parameters:
  • cpuarray (Sequence[Scalar] | None) – Host array. Will be lazily allocated if not provided.

  • viennaclvechandle (Any | None) – Address of the array on the GPU. Will be lazily allocated if not provided.

  • size (LayoutSizeSpec | None) – Vector size.

  • bsize (int | None) – Vector block size. If None, bsize = 1.

  • comm (Comm | None) – MPI communicator, defaults to Sys.getDefaultComm.

Return type:

Self

Source code at petsc4py/PETSc/Vec.pyx:486

createWithArray(array, size=None, bsize=None, comm=None)#

Create a vector using a provided array.

Collective.

This method will create either a Type.SEQ or Type.MPI depending on the size of the communicator.

Parameters:
Return type:

Self

Source code at petsc4py/PETSc/Vec.pyx:319

createWithDLPack(dltensor, size=None, bsize=None, comm=None)#

Create a vector wrapping a DLPack object, sharing the same memory.

Collective.

This operation does not modify the storage of the original tensor and should be used with contiguous tensors only. If the tensor is stored in row-major order (e.g. PyTorch tensors), the resulting vector will look like an unrolled tensor using row-major order.

The resulting vector type will be one of Type.SEQ, Type.MPI, Type.SEQCUDA, Type.MPICUDA, Type.SEQHIP or Type.MPIHIP depending on the type of dltensor and the number of processes in the communicator.

Parameters:
Return type:

Self

Source code at petsc4py/PETSc/Vec.pyx:545

destroy()#

Destroy the vector.

Collective.

See also

create, VecDestroy

Source code at petsc4py/PETSc/Vec.pyx:163

Return type:

Self

dot(vec)#

Return the dot product with vec.

Collective.

For complex numbers this computes yᴴ·x with self as x, vec as y and where yᴴ denotes the conjugate transpose of y.

Use tDot for the indefinite form yᵀ·x where yᵀ denotes the transpose of y.

Parameters:

vec (Vec) – Vector to compute the dot product with.

Return type:

Scalar

See also

dotBegin, dotEnd, tDot, VecDot

Source code at petsc4py/PETSc/Vec.pyx:1787

dotBegin(vec)#

Begin computing the dot product.

Collective.

This should be paired with a call to dotEnd.

Parameters:

vec (Vec) – Vector to compute the dot product with.

Return type:

None

See also

dotEnd, dot, VecDotBegin

Source code at petsc4py/PETSc/Vec.pyx:1812

dotEnd(vec)#

Finish computing the dot product initiated with dotBegin.

Collective.

See also

dotBegin, dot, VecDotEnd

Source code at petsc4py/PETSc/Vec.pyx:1832

Parameters:

vec (Vec)

Return type:

Scalar

dotNorm2(vec)#

Return the dot product with vec and its squared norm.

Collective.

See also

dot, norm, VecDotNorm2

Source code at petsc4py/PETSc/Vec.pyx:2151

Parameters:

vec (Vec)

Return type:

tuple[Scalar, float]

duplicate(array=None)#

Create a new vector with the same type, optionally with data.

Collective.

Parameters:

array (Sequence[Scalar] | None) – Optional values to store in the new vector.

Return type:

Vec

See also

copy, VecDuplicate

Source code at petsc4py/PETSc/Vec.pyx:1682

equal(vec)#

Return whether the vector is equal to another.

Collective.

Parameters:

vec (Vec) – Vector to compare with.

Return type:

bool

See also

VecEqual

Source code at petsc4py/PETSc/Vec.pyx:1768

exp()#

Replace each entry (xₙ) in the vector by exp(xₙ).

Logically collective.

See also

log, VecExp

Source code at petsc4py/PETSc/Vec.pyx:2257

Return type:

None

getArray(readonly=False)#

Return local portion of the vector as an ndarray.

Logically collective.

Parameters:

readonly (bool) – Request read-only access.

Return type:

ArrayScalar

See also

setArray, getBuffer

Source code at petsc4py/PETSc/Vec.pyx:1313

getBlockSize()#

Return the block size of the vector.

Not collective.

See also

VecGetBlockSize

Source code at petsc4py/PETSc/Vec.pyx:1149

Return type:

int

getBuffer(readonly=False)#

Return a buffered view of the local portion of the vector.

Logically collective.

Parameters:

readonly (bool) – Request read-only access.

Returns:

Buffer object wrapping the local portion of the vector data. This can be used either as a context manager providing access as a numpy array or can be passed to array constructors accepting buffered objects such as numpy.asarray.

Return type:

typing.Any

Examples

Accessing the data with a context manager:

>>> vec = PETSc.Vec().createWithArray([1, 2, 3])
>>> with vec.getBuffer() as arr:
...     arr
array([1., 2., 3.])

Converting the buffer to an ndarray:

>>> buf = PETSc.Vec().createWithArray([1, 2, 3]).getBuffer()
>>> np.asarray(buf)
array([1., 2., 3.])

See also

getArray

Source code at petsc4py/PETSc/Vec.pyx:1270

getCLContextHandle()#

Return the OpenCL context associated with the vector.

Not collective.

Returns:

Pointer to underlying CL context. This can be used with pyopencl through pyopencl.Context.from_int_ptr.

Return type:

int

Source code at petsc4py/PETSc/Vec.pyx:1593

getCLMemHandle(mode='rw')#

Return the OpenCL buffer associated with the vector.

Not collective.

Returns:

Pointer to the device buffer. This can be used with pyopencl through pyopencl.Context.from_int_ptr.

Return type:

int

Parameters:

mode (AccessModeSpec)

Notes

This method may incur a host-to-device copy if the device data is out of date and mode is "r" or "rw".

Source code at petsc4py/PETSc/Vec.pyx:1633

getCLQueueHandle()#

Return the OpenCL command queue associated with the vector.

Not collective.

Returns:

Pointer to underlying CL command queue. This can be used with pyopencl through pyopencl.Context.from_int_ptr.

Return type:

int

Source code at petsc4py/PETSc/Vec.pyx:1613

getCUDAHandle(mode='rw')#

Return a pointer to the device buffer.

Not collective.

The returned pointer should be released using restoreCUDAHandle with the same access mode.

Returns:

CUDA device pointer.

Return type:

typing.Any

Parameters:

mode (AccessModeSpec)

Notes

This method may incur a host-to-device copy if the device data is out of date and mode is "r" or "rw".

Source code at petsc4py/PETSc/Vec.pyx:1423

getDM()#

Return the DM associated to the vector.

Not collective.

See also

setDM, VecGetDM

Source code at petsc4py/PETSc/Vec.pyx:3500

Return type:

DM

getHIPHandle(mode='rw')#

Return a pointer to the device buffer.

Not collective.

The returned pointer should be released using restoreHIPHandle with the same access mode.

Returns:

HIP device pointer.

Return type:

typing.Any

Parameters:

mode (AccessModeSpec)

Notes

This method may incur a host-to-device copy if the device data is out of date and mode is "r" or "rw".

Source code at petsc4py/PETSc/Vec.pyx:1495

getLGMap()#

Return the local-to-global mapping.

Not collective.

Source code at petsc4py/PETSc/Vec.pyx:2899

Return type:

LGMap

getLocalSize()#

Return the local size of the vector.

Not collective.

Source code at petsc4py/PETSc/Vec.pyx:1107

Return type:

int

getLocalVector(lvec, readonly=False)#

Maps the local portion of the vector into a local vector.

Logically collective.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:1223

getNestSubVecs()#

Return all the vectors contained in the nested vector.

Not collective.

Source code at petsc4py/PETSc/Vec.pyx:3432

Return type:

list[Vec]

getOffloadMask()#

Return the offloading status of the vector.

Not collective.

Common return values include:

  • 1: PETSC_OFFLOAD_CPU - CPU has valid entries

  • 2: PETSC_OFFLOAD_GPU - GPU has valid entries

  • 3: PETSC_OFFLOAD_BOTH - CPU and GPU are in sync

Returns:

Enum value from PetscOffloadMask describing the offloading status.

Return type:

int

Source code at petsc4py/PETSc/Vec.pyx:1567

getOptionsPrefix()#

Return the prefix used for searching for options in the database.

Not collective.

Source code at petsc4py/PETSc/Vec.pyx:1014

Return type:

str

getOwnershipRange()#

Return the locally owned range of indices (start, end).

Not collective.

Returns:

  • start (int) – The first local element.

  • end (int) – One more than the last local element.

Return type:

tuple[int, int]

Source code at petsc4py/PETSc/Vec.pyx:1163

getOwnershipRanges()#

Return the range of indices owned by each process.

Not collective.

The returned array is the result of exclusive scan of the local sizes.

Source code at petsc4py/PETSc/Vec.pyx:1184

Return type:

ArrayInt

getSize()#

Return the global size of the vector.

Not collective.

Source code at petsc4py/PETSc/Vec.pyx:1093

Return type:

int

getSizes()#

Return the vector sizes.

Not collective.

Source code at petsc4py/PETSc/Vec.pyx:1121

Return type:

LayoutSizeSpec

getSubVector(iset, subvec=None)#

Return a subvector from given indices.

Collective.

Once finished with the subvector it should be returned with restoreSubVector.

Parameters:
  • iset (IS) – Index set describing which indices to extract into the subvector.

  • subvec (Vec | None) – Subvector to copy entries into. If None then a new Vec will be created.

Return type:

Vec

Source code at petsc4py/PETSc/Vec.pyx:3387

getType()#

Return the type of the vector.

Not collective.

See also

setType, VecGetType

Source code at petsc4py/PETSc/Vec.pyx:1079

Return type:

str

getValue(index)#

Return a single value from the vector.

Not collective.

Only values locally stored may be accessed.

Parameters:

index (int) – Location of the value to read.

Return type:

Scalar

Source code at petsc4py/PETSc/Vec.pyx:2720

getValues(indices, values=None)#

Return values from certain locations in the vector.

Not collective.

Only values locally stored may be accessed.

Parameters:
  • indices (Sequence[int]) – Locations of the values to read.

  • values (Sequence[Scalar] | None) – Location to store the collected values. If not provided then a new array will be allocated.

Return type:

ArrayScalar

Source code at petsc4py/PETSc/Vec.pyx:2742

getValuesStagStencil(indices, values=None)#

Not implemented.

Source code at petsc4py/PETSc/Vec.pyx:2767

Return type:

None

ghostUpdate(addv=None, mode=None)#

Update ghosted vector entries.

Neighborwise collective.

Parameters:
Return type:

None

Examples

To accumulate ghost region values onto owning processes:

>>> vec.ghostUpdate(InsertMode.ADD_VALUES, ScatterMode.REVERSE)

Update ghost regions:

>>> vec.ghostUpdate(InsertMode.INSERT_VALUES, ScatterMode.FORWARD)

Source code at petsc4py/PETSc/Vec.pyx:3331

ghostUpdateBegin(addv=None, mode=None)#

Begin updating ghosted vector entries.

Neighborwise collective.

Source code at petsc4py/PETSc/Vec.pyx:3297

Parameters:
Return type:

None

ghostUpdateEnd(addv=None, mode=None)#

Finish updating ghosted vector entries initiated with ghostUpdateBegin.

Neighborwise collective.

Source code at petsc4py/PETSc/Vec.pyx:3314

Parameters:
Return type:

None

isaxpy(idx, alpha, x)#

Add a scaled reduced-space vector to a subset of the vector.

Logically collective.

This is equivalent to y[idx[i]] += alpha*x[i].

Parameters:
  • idx (IS) – Index set for the reduced space. Negative indices are skipped.

  • alpha (Scalar) – Scale factor.

  • x (Vec) – Reduced-space vector.

Return type:

None

See also

axpy, aypx, axpby, VecISAXPY

Source code at petsc4py/PETSc/Vec.pyx:2480

isset(idx, alpha)#

Set specific elements of the vector to the same value.

Not collective.

Parameters:
  • idx (IS) – Index set specifying the vector entries to set.

  • alpha (Scalar) – Value to set the selected entries to.

Return type:

None

See also

set, zeroEntries, VecISSet

Source code at petsc4py/PETSc/Vec.pyx:2383

load(viewer)#

Load a vector.

Collective.

See also

view, VecLoad

Source code at petsc4py/PETSc/Vec.pyx:1750

Parameters:

viewer (Viewer)

Return type:

Self

localForm()#

Return a context manager for viewing ghost vectors in local form.

Logically collective.

Returns:

Context manager yielding the vector in local (ghosted) form.

Return type:

typing.Any

Notes

This operation does not perform a copy. To obtain up-to-date ghost values ghostUpdateBegin and ghostUpdateEnd must be called first.

Non-ghost values can be found at values[0:nlocal] and ghost values at values[nlocal:nlocal+nghost].

Examples

>>> with vec.localForm() as lf:
...     # compute with lf

Source code at petsc4py/PETSc/Vec.pyx:3264

log()#

Replace each entry in the vector by its natural logarithm.

Logically collective.

See also

exp, VecLog

Source code at petsc4py/PETSc/Vec.pyx:2269

Return type:

None

mDot(vecs, out=None)#

Compute Xᴴ·y with X an array of vectors.

Collective.

Parameters:
  • vecs (Sequence[Vec]) – Array of vectors.

  • out (ArrayScalar | None) – Optional placeholder for the result.

Return type:

ArrayScalar

Source code at petsc4py/PETSc/Vec.pyx:1902

mDotBegin(vecs, out)#

Starts a split phase multiple dot product computation.

Collective.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:1935

mDotEnd(vecs, out)#

Ends a split phase multiple dot product computation.

Collective.

Parameters:
Return type:

ArrayScalar

Source code at petsc4py/PETSc/Vec.pyx:1965

max()#

Return the vector entry with maximum real part and its location.

Collective.

Returns:

  • p (int) – Location of the maximum value. If multiple entries exist with the same value then the smallest index will be returned.

  • val (Scalar) – Minimum value.

Return type:

tuple[int, float]

See also

min, VecMax

Source code at petsc4py/PETSc/Vec.pyx:2203

maxPointwiseDivide(vec)#

Return the maximum of the component-wise absolute value division.

Logically collective.

Equivalent to result = max_i abs(x[i] / y[i]).

Parameters:
  • x – Numerator vector.

  • y – Denominator vector.

  • vec (Vec)

Return type:

float

Source code at petsc4py/PETSc/Vec.pyx:2696

maxpy(alphas, vecs)#

Compute and store y = Σₙ(ɑₙ·Xₙ) + y with X an array of vectors.

Logically collective.

Equivalent to y[:] = alphas[i]*vecs[i, :] + y[:].

Parameters:
  • alphas (Sequence[Scalar]) – Array of scale factors, one for each vector in vecs.

  • vecs (Sequence[Vec]) – Array of vectors.

Return type:

None

See also

axpy, aypx, axpby, waxpy, VecMAXPY

Source code at petsc4py/PETSc/Vec.pyx:2569

min()#

Return the vector entry with minimum real part and its location.

Collective.

Returns:

  • p (int) – Location of the minimum value. If multiple entries exist with the same value then the smallest index will be returned.

  • val (Scalar) – Minimum value.

Return type:

tuple[int, float]

See also

max, VecMin

Source code at petsc4py/PETSc/Vec.pyx:2180

mtDot(vecs, out=None)#

Compute Xᵀ·y with X an array of vectors.

Collective.

Parameters:
  • vecs (Sequence[Vec]) – Array of vectors.

  • out (ArrayScalar | None) – Optional placeholder for the result.

Return type:

ArrayScalar

Source code at petsc4py/PETSc/Vec.pyx:1996

mtDotBegin(vecs, out)#

Starts a split phase transpose multiple dot product computation.

Collective.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:2029

mtDotEnd(vecs, out)#

Ends a split phase transpose multiple dot product computation.

Collective.

Parameters:
Return type:

ArrayScalar

Source code at petsc4py/PETSc/Vec.pyx:2059

norm(norm_type=None)#

Compute the vector norm.

Collective.

A 2-tuple is returned if NormType.NORM_1_AND_2 is specified.

See also

VecNorm, NormType

Source code at petsc4py/PETSc/Vec.pyx:2090

Parameters:

norm_type (NormTypeSpec)

Return type:

float | tuple[float, float]

normBegin(norm_type=None)#

Begin computing the vector norm.

Collective.

This should be paired with a call to normEnd.

Source code at petsc4py/PETSc/Vec.pyx:2112

Parameters:

norm_type (NormTypeSpec)

Return type:

None

normEnd(norm_type=None)#

Finish computations initiated with normBegin.

Collective.

Source code at petsc4py/PETSc/Vec.pyx:2131

Parameters:

norm_type (NormTypeSpec)

Return type:

float | tuple[float, float]

normalize()#

Normalize the vector by its 2-norm.

Collective.

Returns:

The vector norm before normalization.

Return type:

float

See also

norm, VecNormalize

Source code at petsc4py/PETSc/Vec.pyx:2226

permute(order, invert=False)#

Permute the vector in-place with a provided ordering.

Collective.

Parameters:
  • order (IS) – Ordering for the permutation.

  • invert (bool) – Whether to invert the permutation.

Return type:

None

See also

VecPermute

Source code at petsc4py/PETSc/Vec.pyx:2337

placeArray(array)#

Set the local portion of the vector to a provided array.

Not collective.

Source code at petsc4py/PETSc/Vec.pyx:1345

Parameters:

array (Sequence[Scalar])

Return type:

None

pointwiseDivide(x, y)#

Compute and store the component-wise division of two vectors.

Logically collective.

Equivalent to w[i] = x[i] / y[i].

Parameters:
  • x (Vec) – Numerator vector.

  • y (Vec) – Denominator vector.

Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:2618

pointwiseMax(x, y)#

Compute and store the component-wise maximum of two vectors.

Logically collective.

Equivalent to w[i] = max(x[i], y[i]).

Parameters:
  • x (Vec) – Input vectors to find the component-wise maxima.

  • y (Vec) – Input vectors to find the component-wise maxima.

Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:2658

pointwiseMaxAbs(x, y)#

Compute and store the component-wise maximum absolute values.

Logically collective.

Equivalent to w[i] = max(abs(x[i]), abs(y[i])).

Parameters:
  • x (Vec) – Input vectors to find the component-wise maxima.

  • y (Vec) – Input vectors to find the component-wise maxima.

Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:2677

pointwiseMin(x, y)#

Compute and store the component-wise minimum of two vectors.

Logically collective.

Equivalent to w[i] = min(x[i], y[i]).

Parameters:
  • x (Vec) – Input vectors to find the component-wise minima.

  • y (Vec) – Input vectors to find the component-wise minima.

Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:2639

pointwiseMult(x, y)#

Compute and store the component-wise multiplication of two vectors.

Logically collective.

Equivalent to w[i] = x[i] * y[i].

Parameters:
  • x (Vec) – Input vectors to multiply component-wise.

  • y (Vec) – Input vectors to multiply component-wise.

Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:2599

reciprocal()#

Replace each entry in the vector by its reciprocal.

Logically collective.

See also

VecReciprocal

Source code at petsc4py/PETSc/Vec.pyx:2245

Return type:

None

resetArray(force=False)#

Reset the vector to use its default array.

Not collective.

Parameters:

force (bool) – Force the calling of VecResetArray even if no user array has been placed with placeArray.

Returns:

The array previously provided by the user with placeArray. Can be None if force is True and no array was placed before.

Return type:

ArrayScalar

Source code at petsc4py/PETSc/Vec.pyx:1366

restoreCLMemHandle()#

Restore a pointer to the OpenCL buffer obtained with getCLMemHandle.

Not collective.

Source code at petsc4py/PETSc/Vec.pyx:1670

Return type:

None

restoreCUDAHandle(handle, mode='rw')#

Restore a pointer to the device buffer obtained with getCUDAHandle.

Not collective.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:1462

restoreHIPHandle(handle, mode='rw')#

Restore a pointer to the device buffer obtained with getHIPHandle.

Not collective.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:1534

restoreLocalVector(lvec, readonly=False)#

Unmap a local access obtained with getLocalVector.

Logically collective.

Parameters:
  • lvec (Vec) – The local vector.

  • readonly (bool) – Request read-only access.

Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:1246

restoreSubVector(iset, subvec)#

Restore a subvector extracted using getSubVector.

Collective.

Parameters:
  • iset (IS) – Index set describing the indices represented by the subvector.

  • subvec (Vec) – Subvector to be restored.

Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:3413

scale(alpha)#

Scale all entries of the vector.

Collective.

This method sets each entry (xₙ) in the vector to ɑ·xₙ.

Parameters:

alpha (Scalar) – The scaling factor.

Return type:

None

See also

shift, VecScale

Source code at petsc4py/PETSc/Vec.pyx:2403

set(alpha)#

Set all components of the vector to the same value.

Collective.

See also

zeroEntries, isset, VecSet

Source code at petsc4py/PETSc/Vec.pyx:2370

Parameters:

alpha (Scalar)

Return type:

None

setArray(array)#

Set values for the local portion of the vector.

Logically collective.

See also

placeArray

Source code at petsc4py/PETSc/Vec.pyx:1333

Parameters:

array (Sequence[Scalar])

Return type:

None

setBlockSize(bsize)#

Set the block size of the vector.

Logically collective.

See also

VecSetBlockSize

Source code at petsc4py/PETSc/Vec.pyx:1136

Parameters:

bsize (int)

Return type:

None

setDM(dm)#

Associate a DM to the vector.

Not collective.

See also

getDM, VecSetDM

Source code at petsc4py/PETSc/Vec.pyx:3488

Parameters:

dm (DM)

Return type:

None

setFromOptions()#

Configure the vector from the options database.

Collective.

Source code at petsc4py/PETSc/Vec.pyx:1042

Return type:

None

setLGMap(lgmap)#

Set the local-to-global mapping.

Logically collective.

This allows users to insert vector entries using a local numbering with setValuesLocal.

Source code at petsc4py/PETSc/Vec.pyx:2884

Parameters:

lgmap (LGMap)

Return type:

None

setMPIGhost(ghosts)#

Set the ghost points for a ghosted vector.

Collective.

Parameters:

ghosts (Sequence[int]) – Global indices of ghost points.

Return type:

None

See also

createGhost

Source code at petsc4py/PETSc/Vec.pyx:3366

setNestSubVecs(sx, idxm=None)#

Set the component vectors at specified indices in the nested vector.

Not collective.

Parameters:
  • sx (Sequence[Vec]) – Array of component vectors.

  • idxm (Sequence[int] | None) – Indices of the component vectors, defaults to range(len(sx)).

Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:3454

setOption(option, flag)#

Set option.

Collective.

See also

VecSetOption

Source code at petsc4py/PETSc/Vec.pyx:1067

Parameters:
Return type:

None

setOptionsPrefix(prefix)#

Set the prefix used for searching for options in the database.

Logically collective.

Source code at petsc4py/PETSc/Vec.pyx:1000

Parameters:

prefix (str | None)

Return type:

None

setRandom(random=None)#

Set all components of the vector to random numbers.

Collective.

Parameters:

random (Random | None) – Random number generator. If None then one will be created internally.

Return type:

None

See also

VecSetRandom

Source code at petsc4py/PETSc/Vec.pyx:2317

setSizes(size, bsize=None)#

Set the local and global sizes of the vector.

Collective.

Parameters:
Return type:

None

See also

getSizes, VecSetSizes

Source code at petsc4py/PETSc/Vec.pyx:218

setType(vec_type)#

Set the vector type.

Collective.

Parameters:

vec_type (Type | str) – The vector type.

Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:199

setUp()#

Set up the internal data structures for using the vector.

Collective.

See also

create, destroy, VecSetUp

Source code at petsc4py/PETSc/Vec.pyx:1054

Return type:

Self

setValue(index, value, addv=None)#

Insert or add a single value in the vector.

Not collective.

Parameters:
  • index (int) – Location to write to. Negative indices are ignored.

  • value (Scalar) – Value to insert at index.

  • addv (InsertModeSpec) – Insertion mode.

Return type:

None

Notes

The values may be cached so assemblyBegin and assemblyEnd must be called after all calls of this method are completed.

Multiple calls to setValue cannot be made with different values for addv without intermediate calls to assemblyBegin and assemblyEnd.

Source code at petsc4py/PETSc/Vec.pyx:2771

setValueLocal(index, value, addv=None)#

Insert or add a single value in the vector using a local numbering.

Not collective.

Parameters:
  • index (int) – Location to write to.

  • value (Scalar) – Value to insert at index.

  • addv (InsertModeSpec) – Insertion mode.

Return type:

None

Notes

The values may be cached so assemblyBegin and assemblyEnd must be called after all calls of this method are completed.

Multiple calls to setValueLocal cannot be made with different values for addv without intermediate calls to assemblyBegin and assemblyEnd.

Source code at petsc4py/PETSc/Vec.pyx:2914

setValues(indices, values, addv=None)#

Insert or add multiple values in the vector.

Not collective.

Parameters:
Return type:

None

Notes

The values may be cached so assemblyBegin and assemblyEnd must be called after all calls of this method are completed.

Multiple calls to setValues cannot be made with different values for addv without intermediate calls to assemblyBegin and assemblyEnd.

Source code at petsc4py/PETSc/Vec.pyx:2808

setValuesBlocked(indices, values, addv=None)#

Insert or add blocks of values in the vector.

Not collective.

Equivalent to x[bs*indices[i]+j] = y[bs*i+j] for 0 <= i < len(indices), 0 <= j < bs and bs block_size.

Parameters:
  • indices (Sequence[int]) – Block indices to write to. Negative indices are ignored.

  • values (Sequence[Scalar]) – Values to insert at indices. Should have length len(indices) * vec.block_size.

  • addv (InsertModeSpec) – Insertion mode.

Return type:

None

Notes

The values may be cached so assemblyBegin and assemblyEnd must be called after all calls of this method are completed.

Multiple calls to setValuesBlocked cannot be made with different values for addv without intermediate calls to assemblyBegin and assemblyEnd.

Source code at petsc4py/PETSc/Vec.pyx:2842

setValuesBlockedLocal(indices, values, addv=None)#

Insert or add blocks of values in the vector with a local numbering.

Not collective.

Equivalent to x[bs*indices[i]+j] = y[bs*i+j] for 0 <= i < len(indices), 0 <= j < bs and bs block_size.

Parameters:
  • indices (Sequence[int]) – Local block indices to write to.

  • values (Sequence[Scalar]) – Values to insert at indices. Should have length len(indices) * vec.block_size.

  • addv (InsertModeSpec) – Insertion mode.

Return type:

None

Notes

The values may be cached so assemblyBegin and assemblyEnd must be called after all calls of this method are completed.

Multiple calls to setValuesBlockedLocal cannot be made with different values for addv without intermediate calls to assemblyBegin and assemblyEnd.

Source code at petsc4py/PETSc/Vec.pyx:2985

setValuesLocal(indices, values, addv=None)#

Insert or add multiple values in the vector with a local numbering.

Not collective.

Parameters:
Return type:

None

Notes

The values may be cached so assemblyBegin and assemblyEnd must be called after all calls of this method are completed.

Multiple calls to setValuesLocal cannot be made with different values for addv without intermediate calls to assemblyBegin and assemblyEnd.

Source code at petsc4py/PETSc/Vec.pyx:2951

setValuesStagStencil(indices, values, addv=None)#

Not implemented.

Source code at petsc4py/PETSc/Vec.pyx:2880

Return type:

None

shift(alpha)#

Shift all entries in the vector.

Collective.

This method sets each entry (xₙ) in the vector to xₙ + ɑ.

Parameters:

alpha (Scalar) – The shift to apply to the vector values.

Return type:

None

See also

scale, VecShift

Source code at petsc4py/PETSc/Vec.pyx:2423

sqrtabs()#

Replace each entry (xₙ) in the vector by √|xₙ|.

Logically collective.

See also

VecSqrtAbs

Source code at petsc4py/PETSc/Vec.pyx:2281

Return type:

None

strideGather(field, vec, addv=None)#

Insert component values into a single-component vector.

Collective.

The current vector is expected to be multi-component (block_size greater than 1) and the target vector is expected to be single-component.

Parameters:
  • field (int) – Component index. Must be between 0 and vec.block_size.

  • vec (Vec) – Single-component vector to be inserted into.

  • addv (InsertModeSpec) – Insertion mode.

Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:3231

strideMax(field)#

Return the maximum of entries in a subvector.

Collective.

Equivalent to max(x[field], x[field+bs], x[field+2*bs], ...) where bs is block_size.

Parameters:

field (int) – Component index. Must be between 0 and vec.block_size.

Returns:

  • int – Location of maximum.

  • float – Maximum value.

Return type:

tuple[int, float]

Source code at petsc4py/PETSc/Vec.pyx:3137

strideMin(field)#

Return the minimum of entries in a subvector.

Collective.

Equivalent to min(x[field], x[field+bs], x[field+2*bs], ...) where bs is block_size.

Parameters:

field (int) – Component index. Must be between 0 and vec.block_size.

Returns:

  • int – Location of minimum.

  • float – Minimum value.

Return type:

tuple[int, float]

Source code at petsc4py/PETSc/Vec.pyx:3106

strideNorm(field, norm_type=None)#

Return the norm of entries in a subvector.

Collective.

Equivalent to norm(x[field], x[field+bs], x[field+2*bs], ...) where bs is block_size.

Parameters:
  • field (int) – Component index. Must be between 0 and vec.block_size.

  • norm_type (NormTypeSpec) – The norm type.

Return type:

float | tuple[float, float]

Source code at petsc4py/PETSc/Vec.pyx:3168

strideScale(field, alpha)#

Scale a component of the vector.

Logically collective.

Parameters:
  • field (int) – Component index. Must be between 0 and vec.block_size.

  • alpha (Scalar) – Factor to multiple the component entries by.

Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:3062

strideScatter(field, vec, addv=None)#

Scatter entries into a component of another vector.

Collective.

The current vector is expected to be single-component (block_size of 1) and the target vector is expected to be multi-component.

Parameters:
  • field (int) – Component index. Must be between 0 and vec.block_size.

  • vec (Vec) – Multi-component vector to be scattered into.

  • addv (InsertModeSpec) – Insertion mode.

Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:3200

strideSum(field)#

Sum subvector entries.

Collective.

Equivalent to sum(x[field], x[field+bs], x[field+2*bs], ...) where bs is block_size.

Parameters:

field (int) – Component index. Must be between 0 and vec.block_size.

Return type:

Scalar

Source code at petsc4py/PETSc/Vec.pyx:3083

sum()#

Return the sum of all the entries of the vector.

Collective.

See also

VecSum

Source code at petsc4py/PETSc/Vec.pyx:2166

Return type:

Scalar

swap(vec)#

Swap the content of two vectors.

Logically collective.

Parameters:

vec (Vec) – The vector to swap data with.

Return type:

None

See also

VecSwap

Source code at petsc4py/PETSc/Vec.pyx:2443

tDot(vec)#

Return the indefinite dot product with vec.

Collective.

This computes yᵀ·x with self as x, vec as y and where yᵀ denotes the transpose of y.

Parameters:

vec (Vec) – Vector to compute the indefinite dot product with.

Return type:

Scalar

Source code at petsc4py/PETSc/Vec.pyx:1846

tDotBegin(vec)#

Begin computing the indefinite dot product.

Collective.

This should be paired with a call to tDotEnd.

Parameters:

vec (Vec) – Vector to compute the indefinite dot product with.

Return type:

None

Source code at petsc4py/PETSc/Vec.pyx:1868

tDotEnd(vec)#

Finish computing the indefinite dot product initiated with tDotBegin.

Collective.

Source code at petsc4py/PETSc/Vec.pyx:1888

Parameters:

vec (Vec)

Return type:

Scalar

toDLPack(mode='rw')#

Return a DLPack PyCapsule wrapping the vector data.

Collective.

Parameters:

mode (AccessModeSpec) – Access mode for the vector.

Returns:

Capsule of a DLPack tensor wrapping a Vec.

Return type:

PyCapsule

Notes

It is important that the access mode is respected by the consumer as this is not enforced internally.

See also

createWithDLPack

Source code at petsc4py/PETSc/Vec.pyx:726

view(viewer=None)#

Display the vector.

Collective.

Parameters:

viewer (Viewer | None) – A Viewer instance or None for the default viewer.

Return type:

None

See also

load, VecView

Source code at petsc4py/PETSc/Vec.pyx:144

waxpy(alpha, x, y)#

Compute and store w = ɑ·x + y.

Logically collective.

Parameters:
  • alpha (Scalar) – Scale factor.

  • x (Vec) – First input vector.

  • y (Vec) – Second input vector.

Return type:

None

See also

axpy, aypx, axpby, maxpy, VecWAXPY

Source code at petsc4py/PETSc/Vec.pyx:2547

zeroEntries()#

Set all entries in the vector to zero.

Logically collective.

See also

set, VecZeroEntries

Source code at petsc4py/PETSc/Vec.pyx:2358

Return type:

None

Attributes Documentation

array#

Alias for array_w.

Source code at petsc4py/PETSc/Vec.pyx:3578

array_r#

Read-only ndarray containing the local portion of the vector.

Source code at petsc4py/PETSc/Vec.pyx:3568

array_w#

Writeable ndarray containing the local portion of the vector.

Source code at petsc4py/PETSc/Vec.pyx:3559

block_size#

The block size.

Source code at petsc4py/PETSc/Vec.pyx:3534

buffer#

Alias for buffer_w.

Source code at petsc4py/PETSc/Vec.pyx:3573

buffer_r#

Read-only buffered view of the local portion of the vector.

Source code at petsc4py/PETSc/Vec.pyx:3554

buffer_w#

Writeable buffered view of the local portion of the vector.

Source code at petsc4py/PETSc/Vec.pyx:3549

local_size#

The local vector size.

Source code at petsc4py/PETSc/Vec.pyx:3529

owner_range#

The locally owned range of indices in the form [low, high).

Source code at petsc4py/PETSc/Vec.pyx:3539

owner_ranges#

The range of indices owned by each process.

Source code at petsc4py/PETSc/Vec.pyx:3544

size#

The global vector size.

Source code at petsc4py/PETSc/Vec.pyx:3524

sizes#

The local and global vector sizes.

Source code at petsc4py/PETSc/Vec.pyx:3516