PetscSection - Low-level Interface

The PetscSection component provides a flexible mechanism for describing the layout of degrees of freedom (DOFs) in a discretization, particularly for finite element methods and other multi-field problems where different points have different numbers of DOFs.

Overview

PetscSection is used to:

  • Define DOF layouts: Specify how many DOFs exist at each point (vertex, edge, face, cell)
  • Multi-field problems: Handle systems with multiple variables per point
  • Local-to-global mapping: Convert between local and global numbering
  • FEM discretizations: Manage DOF distributions for finite elements
  • DM integration: Work with DM objects for mesh-based computations

A section maps from points (mesh entities like vertices, cells) to DOF ranges, answering:

  • How many DOFs are at point p?
  • What is the offset of the DOFs for point p?
  • Which field do specific DOFs belong to?

Basic Usage

using PETSc, MPI

# Initialize MPI and PETSc
MPI.Init()
petsclib = PETSc.getlib()
PETSc.initialize(petsclib)
PetscInt = petsclib.PetscInt

# Create a section
section = Ref{LibPETSc.PetscSection}()
LibPETSc.PetscSectionCreate(petsclib, LibPETSc.PETSC_COMM_SELF, section)

# Set chart: range of valid point indices [pStart, pEnd)
LibPETSc.PetscSectionSetChart(petsclib, section[], 0, 10)

# Set DOF count for each point
for p in 0:9
    num_dofs = (p < 4) ? 1 : 2  # Different DOFs per point
    LibPETSc.PetscSectionSetDof(petsclib, section[], p, num_dofs)
end

# Setup: compute offsets
LibPETSc.PetscSectionSetUp(petsclib, section[])

# Query the section
dof = LibPETSc.PetscSectionGetDof(petsclib, section[], 5)
println("DOF count for point 5: ", dof)

offset = LibPETSc.PetscSectionGetOffset(petsclib, section[], 5)
println("Offset for point 5: ", offset)

# Get total storage size
storage_size = LibPETSc.PetscSectionGetStorageSize(petsclib, section[])
println("Total storage size: ", storage_size)

# Cleanup
LibPETSc.PetscSectionDestroy(petsclib, section)  # pass the Ref directly

PETSc.finalize(petsclib)
MPI.Finalize()

Multi-Field Sections

For problems with multiple fields (e.g., velocity + pressure):

# Create section with 2 fields
section = Ref{LibPETSc.PetscSection}()
LibPETSc.PetscSectionCreate(petsclib, LibPETSc.PETSC_COMM_SELF, section)
@assert section[] != C_NULL

LibPETSc.PetscSectionSetNumFields(petsclib, section[], 2)

# Set field names
LibPETSc.PetscSectionSetFieldName(petsclib, section[], 0, "velocity")
LibPETSc.PetscSectionSetFieldName(petsclib, section[], 1, "pressure")

# Set chart
LibPETSc.PetscSectionSetChart(petsclib, section[], 0, 10)

# Set field components: velocity has 3 components (vx, vy, vz), pressure has 1
LibPETSc.PetscSectionSetFieldComponents(petsclib, section[], 0, 3)
LibPETSc.PetscSectionSetFieldComponents(petsclib, section[], 1, 1)

# Set DOFs per field per point
for p in 0:9
    LibPETSc.PetscSectionSetFieldDof(petsclib, section[], p, 0, 3)  # 3 velocity DOFs
    LibPETSc.PetscSectionSetFieldDof(petsclib, section[], p, 1, 1)  # 1 pressure DOF
    LibPETSc.PetscSectionSetDof(petsclib, section[], p, 4)          # Total: 4 DOFs
end

LibPETSc.PetscSectionSetUp(petsclib, section[])

Constrained DOFs

Mark certain DOFs as constrained (e.g., for boundary conditions):

# Set chart and DOFs...
# Note: Constraint-related functions may have wrapper issues

# Set constraint DOF count
# LibPETSc.PetscSectionSetConstraintDof(petsclib, section[], point, num_constrained)

# Specify which DOFs are constrained
# constrained_indices = PetscInt[0, 2]
# LibPETSc.PetscSectionSetConstraintIndices(petsclib, section[], point, constrained_indices)

# Query constrained storage size
# LibPETSc.PetscSectionGetConstrainedStorageSize(petsclib, section[])

Integration with DM

Sections are commonly used with DM objects:

# Get section from DM
dm_section = Ref{LibPETSc.PetscSection}()
# LibPETSc.DMGetSection(petsclib, dm, dm_section)

# Set section on DM
# LibPETSc.DMSetSection(petsclib, dm, section[])

# Get local section (ghosted)
# local_section = Ref{LibPETSc.PetscSection}()
# LibPETSc.DMGetLocalSection(petsclib, dm, local_section)

Common Workflows

1. FEM Discretization

# Create section for Q1 finite elements on structured grid
# - Vertices: 1 DOF each
# - Cells: 0 DOFs

section = Ref{LibPETSc.PetscSection}()
LibPETSc.PetscSectionCreate(petsclib, comm, section)
@assert section[] != C_NULL

LibPETSc.PetscSectionSetChart(petsclib, section[], vStart, cEnd)

# Set DOFs (vStart to vEnd are vertices, vEnd to cEnd are cells)
for v in vStart:vEnd-1
    LibPETSc.PetscSectionSetDof(petsclib, section[], v, 1)
end

LibPETSc.PetscSectionSetUp(petsclib, section[])

2. Point Closure

Get all DOFs in the closure of a point (point + its boundary):

# Get closure DOFs for a cell
closure_size = Ref{PetscInt}()
closure = Ref{Ptr{PetscInt}}()
# LibPETSc.DMPlexGetTransitiveClosure(petsclib, dm, cell, PETSC_TRUE, closure_size, closure)

# Map closure points to DOF offsets using section
# ... use section to get DOF offsets for each point in closure ...

# LibPETSc.DMPlexRestoreTransitiveClosure(petsclib, dm, cell, PETSC_TRUE, closure_size, closure)

Querying Section Properties

# Get total number of fields
num_fields = LibPETSc.PetscSectionGetNumFields(petsclib, section[])

# Get field components
components = LibPETSc.PetscSectionGetFieldComponents(petsclib, section[], field)

# Get field name (returns String directly)
name = LibPETSc.PetscSectionGetFieldName(petsclib, section[], field)

# Get maximum DOF count across all points
max_dof = LibPETSc.PetscSectionGetMaxDof(petsclib, section[])

Cloning and Permutation

# Clone a section (wrapper may have issues, use with caution)
new_section = Ref{LibPETSc.PetscSection}()
# LibPETSc.PetscSectionClone(petsclib, section[], new_section)

# Note: Clone/Permute functions may require direct ccall if wrapper signatures are incorrect

Function Reference

PETSc.LibPETSc.PetscSectionAddConstraintDofMethod
PetscSectionAddConstraintDof(petsclib::PetscLibType,s::PetscSection, point::PetscInt, numDof::PetscInt)

Increment the number of constrained degrees of freedom associated with a given point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • numDof - the number of additional dof which are fixed by constraints

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionAddDof(), PetscSectionGetConstraintDof(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionAddDofMethod
PetscSectionAddDof(petsclib::PetscLibType,s::PetscSection, point::PetscInt, numDof::PetscInt)

Adds to the total number of degrees of freedom associated with a given point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • numDof - the number of additional dof

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetDof(), PetscSectionSetDof(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionAddFieldConstraintDofMethod
PetscSectionAddFieldConstraintDof(petsclib::PetscLibType,s::PetscSection, point::PetscInt, field::PetscInt, numDof::PetscInt)

Increment the number of constrained degrees of freedom associated with a given field on a point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • field - the field
  • numDof - the number of additional dof which are fixed by constraints

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionAddDof(), PetscSectionGetFieldConstraintDof(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionAddFieldDofMethod
PetscSectionAddFieldDof(petsclib::PetscLibType,s::PetscSection, point::PetscInt, field::PetscInt, numDof::PetscInt)

Adds a number of degrees of freedom associated with a field on a given point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • field - the field
  • numDof - the number of dof

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionSetFieldDof(), PetscSectionGetFieldDof(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionArrayViewMethod
PetscSectionArrayView(petsclib::PetscLibType,s::PetscSection, array::Cvoid, data_type::PetscDataType, viewer::PetscViewer)

View an array, using the section to structure the values

Collective

Input Parameters:

  • s - the organizing PetscSection
  • array - the array of values
  • data_type - the PetscDataType of the array
  • viewer - the PetscViewer

Level: developer

-seealso: PetscSection, PetscViewer, PetscSectionCreate(), VecSetValuesSection(), PetscSectionVecView()

External Links

source
PETSc.LibPETSc.PetscSectionCloneMethod
PetscSectionClone(petsclib::PetscLibType,section::PetscSection, newSection::PetscSection)

Creates a shallow (if possible) copy of the PetscSection

Collective

Input Parameter:

  • section - the PetscSection

Output Parameter:

  • newSection - the copy

Level: beginner

-seealso: PetscSection, PetscSection, PetscSectionCreate(), PetscSectionDestroy(), PetscSectionCopy()

External Links

source
PETSc.LibPETSc.PetscSectionCompareMethod
congruent::PetscBool = PetscSectionCompare(petsclib::PetscLibType,s1::PetscSection, s2::PetscSection)

Compares two sections

Collective

Input Parameters:

  • s1 - the first PetscSection
  • s2 - the second PetscSection

Output Parameter:

  • congruent - PETSC_TRUE if the two sections are congruent, PETSC_FALSE otherwise

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionCreate(), PetscSectionCopy(), PetscSectionClone()

External Links

source
PETSc.LibPETSc.PetscSectionCopyMethod
PetscSectionCopy(petsclib::PetscLibType,section::PetscSection, newSection::PetscSection)

Creates a shallow (if possible) copy of the PetscSection

Collective

Input Parameter:

  • section - the PetscSection

Output Parameter:

  • newSection - the copy

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionCreate(), PetscSectionDestroy()

External Links

source
PETSc.LibPETSc.PetscSectionCreateMethod
PetscSectionCreate(petsclib::PetscLibType,comm::MPI_Comm, s::PetscSection)

Allocates a PetscSection and sets the map contents to the default.

Collective

Input Parameters:

  • comm - the MPI communicator
  • s - pointer to the section

Level: beginner

-seealso: PetscSection, PetscSection, PetscSectionSetChart(), PetscSectionDestroy(), PetscSectionCreateGlobalSection()

External Links

source
PETSc.LibPETSc.PetscSectionCreateComponentSubsectionMethod
subs::PetscSection = PetscSectionCreateComponentSubsection(petsclib::PetscLibType,s::PetscSection, len::PetscInt, comps::Vector{PetscInt})

Create a new, smaller PetscSection composed of only selected components

Collective

Input Parameters:

  • s - the PetscSection
  • len - the number of components
  • comps - the component numbers

Output Parameter:

  • subs - the subsection

Level: advanced

-seealso: PetscSection, PetscSection, PetscSectionCreateSupersection(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionCreateGlobalSectionMethod
gsection::PetscSection = PetscSectionCreateGlobalSection(petsclib::PetscLibType,s::PetscSection, sf::PetscSF, usePermutation::PetscBool, includeConstraints::PetscBool, locOffsets::PetscBool)

Create a parallel section describing the global layout using a local (sequential) PetscSection on each MPI process and a PetscSF describing the section point overlap.

Input Parameters:

  • s - The PetscSection for the local field layout
  • sf - The PetscSF describing parallel layout of the section points (leaves are unowned local points)
  • usePermutation - By default this is PETSC_TRUE, meaning any permutation of the local section is transferred to the global section
  • includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
  • localOffsets - If PETSC_TRUE, use local rather than global offsets for the points

Output Parameter:

  • gsection - The PetscSection for the global field layout

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionCreate(), PetscSectionCreateGlobalSectionCensored()

External Links

source
PETSc.LibPETSc.PetscSectionCreateGlobalSectionCensoredMethod
gsection::PetscSection = PetscSectionCreateGlobalSectionCensored(petsclib::PetscLibType,s::PetscSection, sf::PetscSF, includeConstraints::PetscBool, numExcludes::PetscInt, excludes::Vector{PetscInt})

Create a PetscSection describing the globallayout using a local (sequential) PetscSection on each MPI process and an PetscSF describing the section point overlap.

Input Parameters:

  • s - The PetscSection for the local field layout
  • sf - The PetscSF describing parallel layout of the section points
  • includeConstraints - By default this is PETSC_FALSE, meaning that the global vector will not possess constrained dofs
  • numExcludes - The number of exclusion ranges, this must have the same value on all MPI processes
  • excludes - An array [start0, end0, start1, end1, ...] where there are numExcludes pairs and must have the same values on all MPI processes

Output Parameter:

  • gsection - The PetscSection for the global field layout

Level: advanced

-seealso: PetscSection, PetscSection, PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionCreateGlobalSectionLabelMethod
gsection::PetscSection = PetscSectionCreateGlobalSectionLabel(petsclib::PetscLibType,s::PetscSection, sf::PetscSF, includeConstraints::PetscBool, label::DMLabel, labelValue::PetscInt)

Create a section describing the global field layout using the local section and an PetscSF describing the section point overlap.

Collective

Input Parameters:

  • s - The PetscSection for the local field layout
  • sf - The PetscSF describing parallel layout of the section points
  • includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
  • label - The label specifying the points
  • labelValue - The label stratum specifying the points

Output Parameter:

  • gsection - The PetscSection for the global field layout

Level: developer

-seealso: DMLabel, DM, PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionCreateSubdomainSectionMethod
subs::PetscSection = PetscSectionCreateSubdomainSection(petsclib::PetscLibType,s::PetscSection, subpointMap::IS)

Create a new, smaller section with support on a subdomain of the mesh

Collective

Input Parameters:

  • s - the PetscSection
  • subpointMap - a sorted list of points in the original mesh which are in the subdomain

Output Parameter:

  • subs - the subsection

Level: advanced

-seealso: PetscSection, PetscSection, PetscSectionCreateSubmeshSection(), PetscSectionCreateSubsection(), DMPlexGetSubpointMap(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionCreateSubmeshSectionMethod
subs::PetscSection = PetscSectionCreateSubmeshSection(petsclib::PetscLibType,s::PetscSection, subpointIS::IS)

Create a new, smaller section with support on the submesh

Collective

Input Parameters:

  • s - the PetscSection
  • subpointIS - a sorted list of points in the original mesh which are in the submesh

Output Parameter:

  • subs - the subsection

Level: advanced

-seealso: PetscSection, PetscSection, PetscSectionCreateSubdomainSection(), PetscSectionCreateSubsection(), DMPlexGetSubpointMap(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionCreateSubsectionMethod
subs::PetscSection = PetscSectionCreateSubsection(petsclib::PetscLibType,s::PetscSection, len::PetscInt, fields::Vector{PetscInt})

Create a new, smaller PetscSection composed of only selected fields

Collective

Input Parameters:

  • s - the PetscSection
  • len - the number of subfields
  • fields - the subfield numbers

Output Parameter:

  • subs - the subsection

Level: advanced

-seealso: PetscSection, PetscSection, PetscSectionCreateSupersection(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionCreateSupersectionMethod
supers::PetscSection = PetscSectionCreateSupersection(petsclib::PetscLibType,s::Vector{PetscSection}, len::PetscInt)

Create a new, larger section composed of multiple PetscSections

Collective

Input Parameters:

  • s - the input sections
  • len - the number of input sections

Output Parameter:

  • supers - the supersection

Level: advanced

-seealso: PetscSection, PetscSection, PetscSectionCreateSubsection(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetBlockStartsMethod
PetscSectionGetBlockStarts(petsclib::PetscLibType,s::PetscSection, blockStarts::PetscBT)

Returns a table indicating which points start new blocks

Not Collective, No Fortran Support

Input Parameter:

  • s - the PetscSection

Output Parameter:

  • blockStarts - The PetscBT with a 1 for each point that begins a block

-seealso: , IS, PetscSection, PetscSectionSetBlockStarts(), PetscSectionCreate(), DMCreateMatrix(), MatSetVariableBlockSizes()

External Links

source
PETSc.LibPETSc.PetscSectionGetChartMethod
pStart::PetscInt,pEnd::PetscInt = PetscSectionGetChart(petsclib::PetscLibType,s::PetscSection)

Returns the range [pStart, pEnd) in which points (indices) lie for this PetscSection on this MPI process

Not Collective

Input Parameter:

  • s - the PetscSection

Output Parameters:

  • pStart - the first point
  • pEnd - one past the last point

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionSetChart(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetClosureIndexMethod
PetscSectionGetClosureIndex(petsclib::PetscLibType,section::PetscSection, obj::PetscObject, clSection::PetscSection, clPoints::IS)

Get the cache of points in the closure of each point in the section set with PetscSectionSetClosureIndex()

Collective

Input Parameters:

  • section - The PetscSection
  • obj - A PetscObject which serves as the key for this index

Output Parameters:

  • clSection - PetscSection giving the size of the closure of each point
  • clPoints - IS giving the points in each closure

Level: advanced

-seealso: PetscSection, PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex()

External Links

source
PETSc.LibPETSc.PetscSectionGetClosureInversePermutationMethod
PetscSectionGetClosureInversePermutation(petsclib::PetscLibType,section::PetscSection, obj::PetscObject, depth::PetscInt, clSize::PetscInt, perm::IS)

Get the inverse dof permutation for the closure of each cell in the section, meaning clPerm[oldIndex] = newIndex.

Not Collective

Input Parameters:

  • section - The PetscSection
  • obj - A PetscObject which serves as the key for this index (usually a DM)
  • depth - Depth stratum on which to obtain closure permutation
  • clSize - Closure size to be permuted (e.g., may vary with element topology and degree)

Output Parameter:

  • perm - The dof closure permutation

Level: intermediate

-seealso: PetscSection, PetscSection, IS, PetscSectionSetClosurePermutation(), PetscSectionGetClosureIndex(), PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex()

External Links

source
PETSc.LibPETSc.PetscSectionGetClosurePermutationMethod
PetscSectionGetClosurePermutation(petsclib::PetscLibType,section::PetscSection, obj::PetscObject, depth::PetscInt, clSize::PetscInt, perm::IS)

Get the dof permutation for the closure of each cell in the section, meaning clPerm[newIndex] = oldIndex.

Not Collective

Input Parameters:

  • section - The PetscSection
  • obj - A PetscObject which serves as the key for this index (usually a DM)
  • depth - Depth stratum on which to obtain closure permutation
  • clSize - Closure size to be permuted (e.g., may vary with element topology and degree)

Output Parameter:

  • perm - The dof closure permutation

Level: intermediate

-seealso: PetscSection, PetscSection, IS, PetscSectionSetClosurePermutation(), PetscSectionGetClosureInversePermutation(), PetscSectionGetClosureIndex(), PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex()

External Links

source
PETSc.LibPETSc.PetscSectionGetComponentNameMethod
PetscSectionGetComponentName(petsclib::PetscLibType,s::PetscSection, field::PetscInt, comp::PetscInt, compName::String)

Gets the name of a field component in the PetscSection

Not Collective

Input Parameters:

  • s - the PetscSection
  • field - the field number
  • comp - the component number

Output Parameter:

  • compName - the component name

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetFieldName(), PetscSectionSetNumFields(), PetscSectionGetNumFields(), PetscSectionSetComponentName(), PetscSectionSetFieldName(), PetscSectionGetFieldComponents(), PetscSectionSetFieldComponents()

External Links

source
PETSc.LibPETSc.PetscSectionGetConstrainedStorageSizeMethod
size::PetscInt = PetscSectionGetConstrainedStorageSize(petsclib::PetscLibType,s::PetscSection)

Return the size of an array or local Vec capable of holding all unconstrained degrees of freedom in a PetscSection

Not Collective

Input Parameter:

  • s - the PetscSection

Output Parameter:

  • size - the size of an array which can hold all unconstrained dofs

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetStorageSize(), PetscSectionGetOffset(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetConstraintDofMethod
numDof::PetscInt = PetscSectionGetConstraintDof(petsclib::PetscLibType,s::PetscSection, point::PetscInt)

Return the number of constrained degrees of freedom associated with a given point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point

Output Parameter:

  • numDof - the number of dof which are fixed by constraints

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetDof(), PetscSectionSetConstraintDof(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetConstraintIndicesMethod
indices::Vector{PetscInt} = PetscSectionGetConstraintIndices(petsclib::PetscLibType,s::PetscSection, point::PetscInt)

Get the point dof numbers, in [0, dof), which are constrained for a given point

Not Collective

Input Parameters:

  • s - The PetscSection
  • point - The point

Output Parameter:

  • indices - The constrained dofs

Level: intermediate

-seealso: PetscSection, PetscSectionSetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection

External Links

source
PETSc.LibPETSc.PetscSectionGetDofMethod
numDof::PetscInt = PetscSectionGetDof(petsclib::PetscLibType,s::PetscSection, point::PetscInt)

Return the total number of degrees of freedom associated with a given point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point

Output Parameter:

  • numDof - the number of dof

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionSetDof(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetFieldMethod
PetscSectionGetField(petsclib::PetscLibType,s::PetscSection, field::PetscInt, subs::PetscSection)

Get the PetscSection associated with a single field

Input Parameters:

  • s - The PetscSection
  • field - The field number

Output Parameter:

  • subs - The PetscSection for the given field, note the chart of subs is not set

Level: intermediate

-seealso: PetscSection, PetscSection, IS, PetscSectionSetNumFields()

External Links

source
PETSc.LibPETSc.PetscSectionGetFieldComponentsMethod
numComp::PetscInt = PetscSectionGetFieldComponents(petsclib::PetscLibType,s::PetscSection, field::PetscInt)

Returns the number of field components for the given field.

Not Collective

Input Parameters:

  • s - the PetscSection
  • field - the field number

Output Parameter:

  • numComp - the number of field components

Level: advanced

-seealso: PetscSection, PetscSection, PetscSectionSetFieldComponents(), PetscSectionGetNumFields(), PetscSectionSetComponentName(), PetscSectionGetComponentName()

External Links

source
PETSc.LibPETSc.PetscSectionGetFieldConstraintDofMethod
numDof::PetscInt = PetscSectionGetFieldConstraintDof(petsclib::PetscLibType,s::PetscSection, point::PetscInt, field::PetscInt)

Return the number of constrained degrees of freedom associated with a given field on a point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • field - the field

Output Parameter:

  • numDof - the number of dof which are fixed by constraints

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetDof(), PetscSectionSetFieldConstraintDof(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetFieldConstraintIndicesMethod
indices::Vector{PetscInt} = PetscSectionGetFieldConstraintIndices(petsclib::PetscLibType,s::PetscSection, point::PetscInt, field::PetscInt)

Get the field dof numbers, in [0, fdof), which are constrained

Not Collective

Input Parameters:

  • s - The PetscSection
  • field - The field number
  • point - The point

Output Parameter:

  • indices - The constrained dofs sorted in ascending order, the length is returned by PetscSectionGetConstraintDof().

Level: intermediate

-seealso: PetscSection, PetscSectionSetFieldConstraintIndices(), PetscSectionGetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection

External Links

source
PETSc.LibPETSc.PetscSectionGetFieldDofMethod
numDof::PetscInt = PetscSectionGetFieldDof(petsclib::PetscLibType,s::PetscSection, point::PetscInt, field::PetscInt)

Return the number of degrees of freedom associated with a field on a given point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • field - the field

Output Parameter:

  • numDof - the number of dof

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionSetFieldDof(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetFieldNameMethod
PetscSectionGetFieldName(petsclib::PetscLibType,s::PetscSection, field::PetscInt, fieldName::String)

Returns the name of a field in the PetscSection

Not Collective

Input Parameters:

  • s - the PetscSection
  • field - the field number

Output Parameter:

  • fieldName - the field name

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionSetFieldName(), PetscSectionSetNumFields(), PetscSectionGetNumFields()

External Links

source
PETSc.LibPETSc.PetscSectionGetFieldOffsetMethod
offset::PetscInt = PetscSectionGetFieldOffset(petsclib::PetscLibType,s::PetscSection, point::PetscInt, field::PetscInt)

Return the offset into an array or Vec for the field dof associated with the given point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • field - the field

Output Parameter:

  • offset - the offset

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetOffset(), PetscSectionCreate(), PetscSectionGetFieldPointOffset()

External Links

source
PETSc.LibPETSc.PetscSectionGetFieldPointOffsetMethod
offset::PetscInt = PetscSectionGetFieldPointOffset(petsclib::PetscLibType,s::PetscSection, point::PetscInt, field::PetscInt)

Return the offset for the first field dof associated with the given point relative to the offset for that point for the unnamed default field's first dof

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • field - the field

Output Parameter:

  • offset - the offset

Level: advanced

-seealso: PetscSection, PetscSection, PetscSectionGetOffset(), PetscSectionCreate(), PetscSectionGetFieldOffset()

External Links

source
PETSc.LibPETSc.PetscSectionGetFieldPointSymsMethod
perms::PetscInt,rots::PetscScalar = PetscSectionGetFieldPointSyms(petsclib::PetscLibType,section::PetscSection, field::PetscInt, numPoints::PetscInt, points::PetscInt)

Get the symmetries for a set of points in a field of a PetscSection under specific orientations.

Not Collective

Input Parameters:

  • section - the section
  • field - the field of the section
  • numPoints - the number of points
  • points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an

arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that context, see DMPlexGetConeOrientation()).

Output Parameters:

  • perms - The permutations for the given orientations (or NULL if there is no symmetry or the permutation is the identity).
  • rots - The field rotations symmetries for the given orientations (or NULL if there is no symmetry or the rotations are all

identity).

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionGetPointSyms(), PetscSectionRestoreFieldPointSyms()

External Links

source
PETSc.LibPETSc.PetscSectionGetFieldSymMethod
PetscSectionGetFieldSym(petsclib::PetscLibType,section::PetscSection, field::PetscInt, sym::PetscSectionSym)

Get the symmetries for the data referred to by a field of the section

Collective

Input Parameters:

  • section - the section describing data layout
  • field - the field number

Output Parameter:

  • sym - the symmetry describing the affect of orientation on the access of the data

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionSetFieldSym(), PetscSectionSymCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetIncludesConstraintsMethod
includesConstraints::PetscBool = PetscSectionGetIncludesConstraints(petsclib::PetscLibType,s::PetscSection)

Returns the flag indicating if constrained dofs were included when computing offsets in the PetscSection. The value is set with PetscSectionSetIncludesConstraints()

Not Collective

Input Parameter:

  • s - the PetscSection

Output Parameter:

  • includesConstraints - the flag indicating if constrained dofs were included when computing offsets

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionSetIncludesConstraints()

External Links

source
PETSc.LibPETSc.PetscSectionGetMaxDofMethod
maxDof::PetscInt = PetscSectionGetMaxDof(petsclib::PetscLibType,s::PetscSection)

Return the maximum number of degrees of freedom on any point in the PetscSection

Not Collective

Input Parameter:

  • s - the PetscSection

Output Parameter:

  • maxDof - the maximum dof

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetDof(), PetscSectionSetDof(), PetscSectionAddDof(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetNumFieldsMethod
numFields::PetscInt = PetscSectionGetNumFields(petsclib::PetscLibType,s::PetscSection)

Returns the number of fields in a PetscSection, or 0 if no fields were defined.

Not Collective

Input Parameter:

  • s - the PetscSection

Output Parameter:

  • numFields - the number of fields defined, or 0 if none were defined

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionSetNumFields()

External Links

source
PETSc.LibPETSc.PetscSectionGetOffsetMethod
offset::PetscInt = PetscSectionGetOffset(petsclib::PetscLibType,s::PetscSection, point::PetscInt)

Return the offset into an array or Vec for the dof associated with the given point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point

Output Parameter:

  • offset - the offset

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetFieldOffset(), PetscSectionCreate(), PetscSectionSetPointMajor()

External Links

source
PETSc.LibPETSc.PetscSectionGetOffsetRangeMethod
start::PetscInt,end_::PetscInt = PetscSectionGetOffsetRange(petsclib::PetscLibType,s::PetscSection)

Return the full range of offsets [start, end) for a PetscSection

Not Collective

Input Parameter:

  • s - the PetscSection

Output Parameters:

  • start - the minimum offset
  • end - one more than the maximum offset

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetOffset(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetPermutationMethod
PetscSectionGetPermutation(petsclib::PetscLibType,s::PetscSection, perm::IS)

Returns the permutation of [0, pEnd

Not Collective

Input Parameter:

  • s - the PetscSection

Output Parameter:

  • perm - The permutation as an IS

Level: intermediate

-seealso: , IS, PetscSection, PetscSectionSetPermutation(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetPointLayoutMethod
PetscSectionGetPointLayout(petsclib::PetscLibType,comm::MPI_Comm, s::PetscSection, layout::PetscLayout)

Get a PetscLayout for the points with nonzero dof counts of the unnamed default field within this PetscSections local chart

Collective

Input Parameters:

  • comm - The MPI_Comm
  • s - The PetscSection

Output Parameter:

  • layout - The point layout for the data that defines the section

Level: advanced

-seealso: PetscSection, PetscSection, PetscSectionGetValueLayout(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetPointMajorMethod
pm::PetscBool = PetscSectionGetPointMajor(petsclib::PetscLibType,s::PetscSection)

Returns the flag for dof ordering, PETSC_TRUE if it is point major, PETSC_FALSE if it is field major

Not Collective

Input Parameter:

  • s - the PetscSection

Output Parameter:

  • pm - the flag for point major ordering

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionSetPointMajor()

External Links

source
PETSc.LibPETSc.PetscSectionGetPointSymsMethod
perms::PetscInt,rots::PetscScalar = PetscSectionGetPointSyms(petsclib::PetscLibType,section::PetscSection, numPoints::PetscInt, points::PetscInt)

Get the symmetries for a set of points in a PetscSection under specific orientations.

Not Collective

Input Parameters:

  • section - the section
  • numPoints - the number of points
  • points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an

arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that context, see DMPlexGetConeOrientation()).

Output Parameters:

  • perms - The permutations for the given orientations (or NULL if there is no symmetry or the permutation is the identity).
  • rots - The field rotations symmetries for the given orientations (or NULL if there is no symmetry or the rotations are all

identity).

Example of usage, gathering dofs into a local array (lArray) from a section array (sArray): -seealso: PetscSection, PetscSectionSym, PetscSectionRestorePointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym()

External Links

source
PETSc.LibPETSc.PetscSectionGetStorageSizeMethod
size::PetscInt = PetscSectionGetStorageSize(petsclib::PetscLibType,s::PetscSection)

Return the size of an array or local Vec capable of holding all the degrees of freedom defined in a PetscSection

Not Collective

Input Parameter:

  • s - the PetscSection

Output Parameter:

  • size - the size of an array which can hold all the dofs

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetOffset(), PetscSectionGetConstrainedStorageSize(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetSymMethod
PetscSectionGetSym(petsclib::PetscLibType,section::PetscSection, sym::PetscSectionSym)

Get the symmetries for the data referred to by the section

Not Collective

Input Parameter:

  • section - the section describing data layout

Output Parameter:

  • sym - the symmetry describing the affect of orientation on the access of the data, provided previously by PetscSectionSetSym()

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionSetSym(), PetscSectionSymCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetUseFieldOffsetsMethod
flg::PetscBool = PetscSectionGetUseFieldOffsets(petsclib::PetscLibType,s::PetscSection)

Get the flag indicating if field offsets are used directly in a global section, rather than just the point offset

Not Collective

Input Parameter:

  • s - the global PetscSection

Output Parameter:

  • flg - the flag

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionSetChart(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionGetValueLayoutMethod
PetscSectionGetValueLayout(petsclib::PetscLibType,comm::MPI_Comm, s::PetscSection, layout::PetscLayout)

Get the PetscLayout associated with the section dofs of a PetscSection

Collective

Input Parameters:

  • comm - The MPI_Comm
  • s - The PetscSection

Output Parameter:

  • layout - The dof layout for the section

Level: advanced

-seealso: PetscSection, PetscSection, PetscSectionGetPointLayout(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionHasConstraintsMethod
hasConstraints::PetscBool = PetscSectionHasConstraints(petsclib::PetscLibType,s::PetscSection)

Determine whether a PetscSection has constrained dofs

Not Collective

Input Parameter:

  • s - The PetscSection

Output Parameter:

  • hasConstraints - flag indicating that the section has constrained dofs

Level: intermediate

-seealso: PetscSection, PetscSectionSetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection

External Links

source
PETSc.LibPETSc.PetscSectionLoadMethod
PetscSectionLoad(petsclib::PetscLibType,s::PetscSection, viewer::PetscViewer)

Loads a PetscSection

Collective

Input Parameters:

  • s - the PetscSection object to load
  • viewer - the viewer

Level: beginner

-seealso: PetscSection, PetscSection, PetscSectionCreate(), PetscSectionDestroy(), PetscSectionView()

External Links

source
PETSc.LibPETSc.PetscSectionPermuteMethod
PetscSectionPermute(petsclib::PetscLibType,section::PetscSection, permutation::IS, sectionNew::PetscSection)

Reorder the section according to the input point permutation

Collective

Input Parameters:

  • section - The PetscSection object
  • permutation - The point permutation, old point p becomes new point perm[p]

Output Parameter:

  • sectionNew - The permuted PetscSection

Level: intermediate

-seealso: PetscSection, IS, PetscSection, MatPermute(), PetscSectionSetPermutation()

External Links

source
PETSc.LibPETSc.PetscSectionResetMethod
PetscSectionReset(petsclib::PetscLibType,s::PetscSection)

Frees all section data, the section is then as if PetscSectionCreate() had just been called.

Not Collective

Input Parameter:

  • s - the PetscSection

Level: beginner

-seealso: PetscSection, PetscSection, PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionRestoreFieldPointSymsMethod
PetscSectionRestoreFieldPointSyms(petsclib::PetscLibType,section::PetscSection, field::PetscInt, numPoints::PetscInt, points::PetscInt, perms::PetscInt, rots::PetscScalar)

Restore the symmetries returned by PetscSectionGetFieldPointSyms()

Not Collective

Input Parameters:

  • section - the section
  • field - the field number
  • numPoints - the number of points
  • points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an

arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that context, see DMPlexGetConeOrientation()).

  • perms - The permutations for the given orientations: set to NULL at conclusion
  • rots - The field rotations symmetries for the given orientations: set to NULL at conclusion

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionRestorePointSyms(), petscSectionGetFieldPointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym()

External Links

source
PETSc.LibPETSc.PetscSectionRestorePointSymsMethod
PetscSectionRestorePointSyms(petsclib::PetscLibType,section::PetscSection, numPoints::PetscInt, points::PetscInt, perms::PetscInt, rots::PetscScalar)

Restore the symmetries returned by PetscSectionGetPointSyms()

Not Collective

Input Parameters:

  • section - the section
  • numPoints - the number of points
  • points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an

arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that context, see DMPlexGetConeOrientation()).

  • perms - The permutations for the given orientations: set to NULL at conclusion
  • rots - The field rotations symmetries for the given orientations: set to NULL at conclusion

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionGetPointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym()

External Links

source
PETSc.LibPETSc.PetscSectionSetBlockStartsMethod
PetscSectionSetBlockStarts(petsclib::PetscLibType,s::PetscSection, blockStarts::PetscBT)

Sets a table indicating which points start new blocks

Not Collective, No Fortran Support

Input Parameters:

  • s - the PetscSection
  • blockStarts - The PetscBT with a 1 for each point that begins a block

Level: intermediate

-seealso: , IS, PetscSection, PetscSectionGetBlockStarts(), PetscSectionCreate(), DMCreateMatrix(), MatSetVariableBlockSizes()

External Links

source
PETSc.LibPETSc.PetscSectionSetChartMethod
PetscSectionSetChart(petsclib::PetscLibType,s::PetscSection, pStart::PetscInt, pEnd::PetscInt)

Sets the range [pStart, pEnd) in which points (indices) lie for this PetscSection on this MPI process

Not Collective

Input Parameters:

  • s - the PetscSection
  • pStart - the first point
  • pEnd - one past the last point, pStartpEnd

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetChart(), PetscSectionCreate(), PetscSectionSetNumFields()

External Links

source
PETSc.LibPETSc.PetscSectionSetClosureIndexMethod
PetscSectionSetClosureIndex(petsclib::PetscLibType,section::PetscSection, obj::PetscObject, clSection::PetscSection, clPoints::IS)

Create an internal data structure to speed up closure queries.

Collective

Input Parameters:

  • section - The PetscSection
  • obj - A PetscObject which serves as the key for this index
  • clSection - PetscSection giving the size of the closure of each point
  • clPoints - IS giving the points in each closure

Level: advanced

-seealso: PetscSection, PetscSection, PetscSectionGetClosureIndex(), DMPlexCreateClosureIndex()

External Links

source
PETSc.LibPETSc.PetscSectionSetClosurePermutationMethod
PetscSectionSetClosurePermutation(petsclib::PetscLibType,section::PetscSection, obj::PetscObject, depth::PetscInt, perm::IS)

Set the dof permutation for the closure of each cell in the section, meaning clPerm[newIndex] = oldIndex.

Not Collective

Input Parameters:

  • section - The PetscSection
  • obj - A PetscObject which serves as the key for this index (usually a DM)
  • depth - Depth of points on which to apply the given permutation
  • perm - Permutation of the cell dof closure

Level: intermediate

-seealso: PetscSection, PetscSection, IS, PetscSectionGetClosurePermutation(), PetscSectionGetClosureIndex(), DMPlexCreateClosureIndex(), PetscCopyMode

External Links

source
PETSc.LibPETSc.PetscSectionSetComponentNameMethod
PetscSectionSetComponentName(petsclib::PetscLibType,s::PetscSection, field::PetscInt, comp::PetscInt, compName::String)

Sets the name of a field component in the PetscSection

Not Collective

Input Parameters:

  • s - the PetscSection
  • field - the field number
  • comp - the component number
  • compName - the component name

Level: advanced

-seealso: PetscSection, PetscSection, PetscSectionGetComponentName(), PetscSectionSetNumFields(), PetscSectionGetNumFields(), PetscSectionSetFieldName(), PetscSectionGetFieldComponents(), PetscSectionSetFieldComponents()

External Links

source
PETSc.LibPETSc.PetscSectionSetConstraintDofMethod
PetscSectionSetConstraintDof(petsclib::PetscLibType,s::PetscSection, point::PetscInt, numDof::PetscInt)

Set the number of constrained degrees of freedom associated with a given point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • numDof - the number of dof which are fixed by constraints

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionSetDof(), PetscSectionGetConstraintDof(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionSetConstraintIndicesMethod
PetscSectionSetConstraintIndices(petsclib::PetscLibType,s::PetscSection, point::PetscInt, indices::Vector{PetscInt})

Set the point dof numbers, in [0, dof), which are constrained

Not Collective

Input Parameters:

  • s - The PetscSection
  • point - The point
  • indices - The constrained dofs

Level: intermediate

-seealso: PetscSection, PetscSectionGetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection

External Links

source
PETSc.LibPETSc.PetscSectionSetDofMethod
PetscSectionSetDof(petsclib::PetscLibType,s::PetscSection, point::PetscInt, numDof::PetscInt)

Sets the total number of degrees of freedom associated with a given point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • numDof - the number of dof, these values may be negative -(dof+1) to indicate they are off process

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetDof(), PetscSectionAddDof(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionSetFieldComponentsMethod
PetscSectionSetFieldComponents(petsclib::PetscLibType,s::PetscSection, field::PetscInt, numComp::PetscInt)

Sets the number of field components for the given field.

Not Collective

Input Parameters:

  • s - the PetscSection
  • field - the field number
  • numComp - the number of field components

Level: advanced

-seealso: PetscSection, PetscSection, PetscSectionGetFieldComponents(), PetscSectionSetComponentName(), PetscSectionGetComponentName(), PetscSectionGetNumFields()

External Links

source
PETSc.LibPETSc.PetscSectionSetFieldConstraintDofMethod
PetscSectionSetFieldConstraintDof(petsclib::PetscLibType,s::PetscSection, point::PetscInt, field::PetscInt, numDof::PetscInt)

Set the number of constrained degrees of freedom associated with a given field on a point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • field - the field
  • numDof - the number of dof which are fixed by constraints

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionSetDof(), PetscSectionGetFieldConstraintDof(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionSetFieldConstraintIndicesMethod
PetscSectionSetFieldConstraintIndices(petsclib::PetscLibType,s::PetscSection, point::PetscInt, field::PetscInt, indices::Vector{PetscInt})

Set the field dof numbers, in [0, fdof), which are constrained

Not Collective

Input Parameters:

  • s - The PetscSection
  • point - The point
  • field - The field number
  • indices - The constrained dofs

Level: intermediate

-seealso: PetscSection, PetscSectionSetConstraintIndices(), PetscSectionGetFieldConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection

External Links

source
PETSc.LibPETSc.PetscSectionSetFieldDofMethod
PetscSectionSetFieldDof(petsclib::PetscLibType,s::PetscSection, point::PetscInt, field::PetscInt, numDof::PetscInt)

Sets the number of degrees of freedom associated with a field on a given point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • field - the field
  • numDof - the number of dof, these values may be negative -(dof+1) to indicate they are off process

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetFieldDof(), PetscSectionCreate(), PetscSectionAddDof(), PetscSectionSetDof()

External Links

source
PETSc.LibPETSc.PetscSectionSetFieldNameMethod
PetscSectionSetFieldName(petsclib::PetscLibType,s::PetscSection, field::PetscInt, fieldName::String)

Sets the name of a field in the PetscSection

Not Collective

Input Parameters:

  • s - the PetscSection
  • field - the field number
  • fieldName - the field name

Level: intermediate

-seealso: PetscSection, PetscSectionGetFieldName(), PetscSectionSetNumFields(), PetscSectionGetNumFields()

External Links

source
PETSc.LibPETSc.PetscSectionSetFieldOffsetMethod
PetscSectionSetFieldOffset(petsclib::PetscLibType,s::PetscSection, point::PetscInt, field::PetscInt, offset::PetscInt)

Set the offset into an array or Vec for the dof associated with the given field at a point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • field - the field
  • offset - the offset, these values may be negative indicating the values are off process

Level: developer

-seealso: PetscSection, PetscSection, PetscSectionGetFieldOffset(), PetscSectionSetOffset(), PetscSectionCreate(), PetscSectionSetUp()

External Links

source
PETSc.LibPETSc.PetscSectionSetFieldSymMethod
PetscSectionSetFieldSym(petsclib::PetscLibType,section::PetscSection, field::PetscInt, sym::PetscSectionSym)

Set the symmetries for the data referred to by a field of the section

Collective

Input Parameters:

  • section - the section describing data layout
  • field - the field number
  • sym - the symmetry describing the affect of orientation on the access of the data

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionGetFieldSym(), PetscSectionSymCreate()

External Links

source
PETSc.LibPETSc.PetscSectionSetFromOptionsMethod
PetscSectionSetFromOptions(petsclib::PetscLibType,s::PetscSection)

sets parameters in a PetscSection from the options database

Collective

Input Parameter:

  • s - the PetscSection

Options Database Key:

  • -petscsection_point_major - PETSC_TRUE for point-major order

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionCreate(), PetscSectionDestroy()

External Links

source
PETSc.LibPETSc.PetscSectionSetIncludesConstraintsMethod
PetscSectionSetIncludesConstraints(petsclib::PetscLibType,s::PetscSection, includesConstraints::PetscBool)

Sets the flag indicating if constrained dofs are to be included when computing offsets

Not Collective

Input Parameters:

  • s - the PetscSection
  • includesConstraints - the flag indicating if constrained dofs are to be included when computing offsets

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetIncludesConstraints()

External Links

source
PETSc.LibPETSc.PetscSectionSetNumFieldsMethod
PetscSectionSetNumFields(petsclib::PetscLibType,s::PetscSection, numFields::PetscInt)

Sets the number of fields in a PetscSection

Not Collective

Input Parameters:

  • s - the PetscSection
  • numFields - the number of fields

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetNumFields(), PetscSectionSetChart(), PetscSectionReset()

External Links

source
PETSc.LibPETSc.PetscSectionSetOffsetMethod
PetscSectionSetOffset(petsclib::PetscLibType,s::PetscSection, point::PetscInt, offset::PetscInt)

Set the offset into an array or Vec for the dof associated with the given point.

Not Collective

Input Parameters:

  • s - the PetscSection
  • point - the point
  • offset - the offset, these values may be negative indicating the values are off process

Level: developer

-seealso: PetscSection, PetscSection, PetscSectionGetFieldOffset(), PetscSectionCreate(), PetscSectionSetUp()

External Links

source
PETSc.LibPETSc.PetscSectionSetPermutationMethod
PetscSectionSetPermutation(petsclib::PetscLibType,s::PetscSection, perm::IS)

Sets a permutation of the chart for this section, [0, pEnd

Not Collective

Input Parameters:

  • s - the PetscSection
  • perm - the permutation of points

Level: intermediate

-seealso: , IS, PetscSection, PetscSectionSetUp(), PetscSectionGetPermutation(), PetscSectionPermute(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionSetPointMajorMethod
PetscSectionSetPointMajor(petsclib::PetscLibType,s::PetscSection, pm::PetscBool)

Sets the flag for dof ordering, PETSC_TRUE for point major, otherwise it will be field major

Not Collective

Input Parameters:

  • s - the PetscSection
  • pm - the flag for point major ordering

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionGetPointMajor(), PetscSectionSetPermutation()

External Links

source
PETSc.LibPETSc.PetscSectionSetSymMethod
PetscSectionSetSym(petsclib::PetscLibType,section::PetscSection, sym::PetscSectionSym)

Set the symmetries for the data referred to by the section

Collective

Input Parameters:

  • section - the section describing data layout
  • sym - the symmetry describing the affect of orientation on the access of the data

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionGetSym(), PetscSectionSymCreate()

External Links

source
PETSc.LibPETSc.PetscSectionSetUpMethod
PetscSectionSetUp(petsclib::PetscLibType,s::PetscSection)

Calculate offsets based upon the number of degrees of freedom for each point in preparation for use of the PetscSection

Not Collective

Input Parameter:

  • s - the PetscSection

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionCreate(), PetscSectionSetPermutation()

External Links

source
PETSc.LibPETSc.PetscSectionSetUseFieldOffsetsMethod
PetscSectionSetUseFieldOffsets(petsclib::PetscLibType,s::PetscSection, flg::PetscBool)

Set the flag to use field offsets directly in a global section, rather than just the point offset

Not Collective

Input Parameters:

  • s - the global PetscSection
  • flg - the flag

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionGetUseFieldOffsets(), PetscSectionSetChart(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionSymCopyMethod
PetscSectionSymCopy(petsclib::PetscLibType,sym::PetscSectionSym, nsym::PetscSectionSym)

Copy the symmetries, assuming that the point structure is compatible

Not Collective

Input Parameter:

  • sym - the PetscSectionSym

Output Parameter:

  • nsym - the equivalent symmetries

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()

External Links

source
PETSc.LibPETSc.PetscSectionSymCreateMethod
sym::PetscSectionSym = PetscSectionSymCreate(petsclib::PetscLibType,comm::MPI_Comm)

Creates an empty PetscSectionSym object.

Collective

Input Parameter:

  • comm - the MPI communicator

Output Parameter:

  • sym - pointer to the new set of symmetries

Level: developer

-seealso: PetscSection, PetscSection, PetscSectionSym, PetscSectionSymDestroy()

External Links

source
PETSc.LibPETSc.PetscSectionSymCreateLabelMethod
sym::PetscSectionSym = PetscSectionSymCreateLabel(petsclib::PetscLibType,comm::MPI_Comm, label::DMLabel)

Create a section symmetry that assigns one symmetry to each stratum of a label

Collective

Input Parameters:

  • comm - the MPI communicator for the new symmetry
  • label - the label defining the strata

Output Parameter:

  • sym - the section symmetries

Level: developer

-seealso: DMLabel, DM, PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()

External Links

source
PETSc.LibPETSc.PetscSectionSymDistributeMethod
PetscSectionSymDistribute(petsclib::PetscLibType,sym::PetscSectionSym, migrationSF::PetscSF, dsym::PetscSectionSym)

Distribute the symmetries in accordance with the input PetscSF

Collective

Input Parameters:

  • sym - the PetscSectionSym
  • migrationSF - the distribution map from roots to leaves

Output Parameter:

  • dsym - the redistributed symmetries

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms()

External Links

source
PETSc.LibPETSc.PetscSectionSymGetTypeMethod
type::PetscSectionSymType = PetscSectionSymGetType(petsclib::PetscLibType,sym::PetscSectionSym)

Gets the section symmetry type name (as a string) from the PetscSectionSym.

Not Collective

Input Parameter:

  • sym - The section symmetry

Output Parameter:

  • type - The index set type name

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionSymType, PetscSectionSymSetType(), PetscSectionSymCreate()

External Links

source
PETSc.LibPETSc.PetscSectionSymLabelGetStratumMethod
size::PetscInt,minOrient::PetscInt,maxOrient::PetscInt,perms::PetscInt,rots::PetscScalar = PetscSectionSymLabelGetStratum(petsclib::PetscLibType,sym::PetscSectionSym, stratum::PetscInt)

get the symmetries for the orientations of a stratum

Logically Collective

Input Parameters:

  • sym - the section symmetries
  • stratum - the stratum value in the label that we are assigning symmetries for

Output Parameters:

  • size - the number of dofs for points in the stratum of the label
  • minOrient - the smallest orientation for a point in this stratum
  • maxOrient - one greater than the largest orientation for a ppoint in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
  • perms - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation. A NULL permutation is the identity
  • rots - NULL if there are no rotations, or (maxOrient - minOrient) sets of rotations, one for each orientation. A NULL set of orientations is the identity

Level: developer

-seealso: DMLabel, DM, PetscSectionSymLabelSetStratum(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()

External Links

source
PETSc.LibPETSc.PetscSectionSymLabelSetLabelMethod
PetscSectionSymLabelSetLabel(petsclib::PetscLibType,sym::PetscSectionSym, label::DMLabel)

set the label whose strata will define the points that receive symmetries

Logically

Input Parameters:

  • sym - the section symmetries
  • label - the DMLabel describing the types of points

Level: developer:

-seealso: DMLabel, DM, PetscSectionSymLabelSetStratum(), PetscSectionSymCreateLabel(), PetscSectionGetPointSyms()

External Links

source
PETSc.LibPETSc.PetscSectionSymLabelSetStratumMethod
PetscSectionSymLabelSetStratum(petsclib::PetscLibType,sym::PetscSectionSym, stratum::PetscInt, size::PetscInt, minOrient::PetscInt, maxOrient::PetscInt, mode::PetscCopyMode, perms::PetscInt, rots::PetscScalar)

set the symmetries for the orientations of a stratum

Logically

Input Parameters:

  • sym - the section symmetries
  • stratum - the stratum value in the label that we are assigning symmetries for
  • size - the number of dofs for points in the stratum of the label
  • minOrient - the smallest orientation for a point in this stratum
  • maxOrient - one greater than the largest orientation for a point in this stratum (i.e., orientations are in the range [minOrient, maxOrient))
  • mode - how sym should copy the perms and rots arrays
  • perms - NULL if there are no permutations, or (maxOrient - minOrient) permutations, one for each orientation. A NULL permutation is the identity
  • rots - NULL if there are no rotations, or (maxOrient - minOrient) sets of rotations, one for each orientation. A NULL set of orientations is the identity

Level: developer

-seealso: DMLabel, DM, PetscSectionSymLabelGetStratum(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetPointSyms(), PetscSectionSymCreateLabel()

External Links

source
PETSc.LibPETSc.PetscSectionSymRegisterMethod
PetscSectionSymRegister(petsclib::PetscLibType,sname::String, fnc::external)

Registers a new section symmetry implementation

Not Collective, No Fortran Support

Input Parameters:

  • sname - The name of a new user-defined creation routine
  • function - The creation routine itself

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionSymType, PetscSectionSymCreate(), PetscSectionSymSetType()

External Links

source
PETSc.LibPETSc.PetscSectionSymSetTypeMethod
PetscSectionSymSetType(petsclib::PetscLibType,sym::PetscSectionSym, method::PetscSectionSymType)

Builds a PetscSectionSym, for a particular implementation.

Collective

Input Parameters:

  • sym - The section symmetry object
  • method - The name of the section symmetry type

Level: developer

-seealso: PetscSection, PetscSectionSym, PetscSectionSymType, PetscSectionSymGetType(), PetscSectionSymCreate()

External Links

source
PETSc.LibPETSc.PetscSectionSymViewMethod
PetscSectionSymView(petsclib::PetscLibType,sym::PetscSectionSym, viewer::PetscViewer)

Displays a section symmetry

Collective

Input Parameters:

  • sym - the index set
  • viewer - viewer used to display the set, for example PETSC_VIEWER_STDOUT_SELF.

Level: developer

-seealso: PetscSectionSym, PetscViewer, PetscViewerASCIIOpen()

External Links

source
PETSc.LibPETSc.PetscSectionVecNormMethod
val::Vector{PetscReal} = PetscSectionVecNorm(petsclib::PetscLibType,s::PetscSection, gs::PetscSection, x::PetscVec, type::NormType)

Computes the vector norm of each field

Input Parameters:

  • s - the local Section
  • gs - the global section
  • x - the vector
  • type - one of NORM_1, NORM_2, NORM_INFINITY.

Output Parameter:

  • val - the array of norms

Level: intermediate

-seealso: VecNorm(), PetscSectionCreate()

External Links

source
PETSc.LibPETSc.PetscSectionVecViewMethod
PetscSectionVecView(petsclib::PetscLibType,s::PetscSection, v::PetscVec, viewer::PetscViewer)

View a vector, using the section to structure the values

Collective

Input Parameters:

  • s - the organizing PetscSection
  • v - the Vec
  • viewer - the PetscViewer

Level: developer

-seealso: PetscSection, PetscViewer, PetscSectionCreate(), VecSetValuesSection(), PetscSectionArrayView()

External Links

source
PETSc.LibPETSc.PetscSectionViewMethod
PetscSectionView(petsclib::PetscLibType,s::PetscSection, viewer::PetscViewer)

Views a PetscSection

Collective

Input Parameters:

  • s - the PetscSection object to view
  • viewer - the viewer

Level: beginner

-seealso: PetscSection, PetscSection, PetscSectionCreate(), PetscSectionDestroy(), PetscSectionLoad(), PetscViewer

External Links

source
PETSc.LibPETSc.PetscSectionViewFromOptionsMethod
PetscSectionViewFromOptions(petsclib::PetscLibType,A::PetscSection, obj::PetscObject, name::String)

View the PetscSection based on values in the options database

Collective

Input Parameters:

  • A - the PetscSection object to view
  • obj - Optional object that provides the options prefix used for the options
  • name - command line option

Level: intermediate

-seealso: PetscSection, PetscSection, PetscSectionView, PetscObjectViewFromOptions(), PetscSectionCreate(), PetscSectionView()

External Links

source