petsc4py.PETSc.IS#
- class petsc4py.PETSc.IS#
Bases:
Object
A collection of indices.
IS objects are used to index into vectors and matrices and to set up vector scatters.
See also
Enumerations
Methods Summary
Concatenate index sets stored across processors.
buildTwoSided
([toindx])Create an index set describing a global mapping.
complement
(nmin, nmax)Create a complement index set.
copy
([result])Copy the contents of the index set into another.
create
([comm])Create an IS.
createBlock
(bsize, indices[, comm])Create a blocked index set.
createGeneral
(indices[, comm])Create an IS with indices.
createStride
(size[, first, step, comm])Create an index set consisting of evenly spaced values.
destroy
()Destroy the index set.
difference
(iset)Return the difference between two index sets.
Create a copy of the index set.
embed
(iset, drop)Embed
self
intoiset
.equal
(iset)Return whether the index sets have the same set of indices or not.
expand
(iset)Return the union of two (possibly unsorted) index sets.
Return the indices of an index set with type
IS.Type.BLOCK
.Return the number of elements in a block.
Return the indices of the index set.
getInfo
()Return stride information for an index set with type
IS.Type.STRIDE
.Return the process-local length of the index set.
getSize
()Return the global length of an index set.
getSizes
()Return the local and global sizes of the index set.
Return size and stride information.
getType
()Return the index set type associated with the IS.
invertPermutation
([nlocal])Invert the index set.
Return whether the index set has been declared as an identity.
Return whether an index set has been declared to be a permutation.
isSorted
()Return whether the indices have been sorted.
load
(viewer)Load a stored index set.
renumber
([mult])Renumber the non-negative entries of an index set, starting from 0.
setBlockIndices
(bsize, indices)Set the indices for an index set with type
IS.Type.BLOCK
.setBlockSize
(bs)Set the block size of the index set.
Mark the index set as being an identity.
setIndices
(indices)Set the indices of an index set.
Mark the index set as being a permutation.
setStride
(size[, first, step])Set the stride information for an index set with type
IS.Type.STRIDE
.setType
(is_type)Set the type of the index set.
sort
()Sort the indices of an index set.
sum
(iset)Return the union of two (sorted) index sets.
Convert the index set type to
IS.Type.GENERAL
.union
(iset)Return the union of two (possibly unsorted) index sets.
view
([viewer])Display the index set.
Attributes Summary
View of the index set as an array of integers.
The number of elements in a block.
The indices of the index set.
The local size of the index set.
The global size of the index set.
The local and global sizes of the index set.
Methods Documentation
- allGather()#
Concatenate index sets stored across processors.
Collective.
The returned index set will be the same on every processor.
See also
Source code at petsc4py/PETSc/IS.pyx:306
- Return type:
- buildTwoSided(toindx=None)#
Create an index set describing a global mapping.
Collective.
This function generates an index set that contains new numbers from remote or local on the index set.
- Parameters:
toindx (IS | None) – Index set describing which indices to send, default is to send natural numbering.
- Returns:
New index set containing the new numbers from remote or local.
- Return type:
See also
- complement(nmin, nmax)#
Create a complement index set.
Collective.
The complement set of indices is all indices that are not in the provided set (and within the provided bounds).
- Parameters:
- Return type:
Notes
For a parallel index set, this will generate the local part of the complement on each process.
To generate the entire complement (on each process) of a parallel index set, first call
IS.allGather
and then call this method.See also
- copy(result=None)#
Copy the contents of the index set into another.
Collective.
- Parameters:
result (IS | None) – The target index set. If
None
thenIS.duplicate
is called first.- Returns:
The copied index set. If
result
is notNone
then this is returned here.- Return type:
See also
- create(comm=None)#
Create an IS.
Collective.
- Parameters:
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.- Return type:
See also
- createBlock(bsize, indices, comm=None)#
Create a blocked index set.
Collective.
- Parameters:
- Return type:
See also
- createGeneral(indices, comm=None)#
Create an IS with indices.
Collective.
- Parameters:
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.
- Return type:
See also
- createStride(size, first=0, step=0, comm=None)#
Create an index set consisting of evenly spaced values.
Collective.
- Parameters:
- Return type:
See also
- destroy()#
Destroy the index set.
Collective.
See also
Source code at petsc4py/PETSc/IS.pyx:74
- Return type:
- difference(iset)#
Return the difference between two index sets.
Collective.
- Parameters:
iset (IS) – Index set to compute the difference with.
- Returns:
Index set representing the difference between
self
andiset
.- Return type:
See also
- duplicate()#
Create a copy of the index set.
Collective.
See also
Source code at petsc4py/PETSc/IS.pyx:241
- Return type:
- embed(iset, drop)#
Embed
self
intoiset
.Not collective.
The embedding is performed by finding the locations in
iset
that have the same indices asself
.- Parameters:
- Returns:
The embedded index set.
- Return type:
See also
- equal(iset)#
Return whether the index sets have the same set of indices or not.
Collective.
See also
- expand(iset)#
Return the union of two (possibly unsorted) index sets.
Collective.
To compute the union,
expand
concatenates the two index sets and removes any duplicates.- Parameters:
iset (IS) – Index set to compute the union with.
- Returns:
The new, combined, index set.
- Return type:
See also
- getBlockIndices()#
Return the indices of an index set with type
IS.Type.BLOCK
.Not collective.
See also
Source code at petsc4py/PETSc/IS.pyx:835
- Return type:
- getBlockSize()#
Return the number of elements in a block.
Not collective.
See also
Source code at petsc4py/PETSc/IS.pyx:439
- Return type:
- getIndices()#
Return the indices of the index set.
Not collective.
See also
Source code at petsc4py/PETSc/IS.pyx:791
- Return type:
- getInfo()#
Return stride information for an index set with type
IS.Type.STRIDE
.Not collective.
- Returns:
- Return type:
See also
- getLocalSize()#
Return the process-local length of the index set.
Not collective.
See also
Source code at petsc4py/PETSc/IS.pyx:403
- Return type:
- getSize()#
Return the global length of an index set.
Not collective.
See also
Source code at petsc4py/PETSc/IS.pyx:389
- Return type:
- getSizes()#
Return the local and global sizes of the index set.
Not collective.
- Returns:
- Return type:
See also
- getStride()#
Return size and stride information.
Not collective.
- Returns:
- Return type:
See also
- getType()#
Return the index set type associated with the IS.
Not collective.
See also
Source code at petsc4py/PETSc/IS.pyx:127
- Return type:
- invertPermutation(nlocal=None)#
Invert the index set.
Collective.
For this to be correct the index set must be a permutation.
- Parameters:
nlocal (int | None) – The number of indices on this processor in the resulting index set, defaults to
PETSC_DECIDE
.- Return type:
See also
- isIdentity()#
Return whether the index set has been declared as an identity.
Collective.
See also
Source code at petsc4py/PETSc/IS.pyx:538
- Return type:
- isPermutation()#
Return whether an index set has been declared to be a permutation.
Logically collective.
See also
Source code at petsc4py/PETSc/IS.pyx:511
- Return type:
- isSorted()#
Return whether the indices have been sorted.
Collective.
See also
Source code at petsc4py/PETSc/IS.pyx:484
- Return type:
- load(viewer)#
Load a stored index set.
Collective.
- Parameters:
viewer (Viewer) – Binary file viewer, either
Viewer.Type.BINARY
orViewer.Type.HDF5
.- Return type:
See also
- renumber(mult=None)#
Renumber the non-negative entries of an index set, starting from 0.
Collective.
- Parameters:
mult (IS | None) – The multiplicity of each entry in
self
, default implies a multiplicity of 1.- Returns:
- Return type:
See also
- setBlockIndices(bsize, indices)#
Set the indices for an index set with type
IS.Type.BLOCK
.Collective.
- Parameters:
- Return type:
See also
- setBlockSize(bs)#
Set the block size of the index set.
Logically collective.
See also
- setIdentity()#
Mark the index set as being an identity.
Logically collective.
See also
Source code at petsc4py/PETSc/IS.pyx:525
- Return type:
- setIndices(indices)#
Set the indices of an index set.
Logically collective.
The index set is assumed to be of type
IS.Type.GENERAL
.See also
- setPermutation()#
Mark the index set as being a permutation.
Logically collective.
See also
Source code at petsc4py/PETSc/IS.pyx:498
- Return type:
- setStride(size, first=0, step=1)#
Set the stride information for an index set with type
IS.Type.STRIDE
.Logically collective.
- Parameters:
- Return type:
See also
- setType(is_type)#
Set the type of the index set.
Collective.
See also
- sort()#
Sort the indices of an index set.
Collective.
See also
Source code at petsc4py/PETSc/IS.pyx:471
- Return type:
- sum(iset)#
Return the union of two (sorted) index sets.
Collective.
See also
- toGeneral()#
Convert the index set type to
IS.Type.GENERAL
.Collective.
See also
Source code at petsc4py/PETSc/IS.pyx:322
- Return type:
- union(iset)#
Return the union of two (possibly unsorted) index sets.
Collective.
This function will call either
ISSum
orISExpand
depending on whether or not the input sets are already sorted.Sequential only (as
ISSum
is sequential only).- Parameters:
iset (IS) – Index set to compute the union with.
- Returns:
The new, combined, index set.
- Return type:
- view(viewer=None)#
Display the index set.
Collective.
See also
Attributes Documentation
- array#
View of the index set as an array of integers.
Not collective.
- block_size#
The number of elements in a block.
Not collective.
See also
- indices#
The indices of the index set.
Not collective.
See also
- local_size#
The local size of the index set.
Not collective.
See also
- size#
The global size of the index set.
Not collective.
See also
- sizes#
The local and global sizes of the index set.
Not collective.
See also