petsc4py.PETSc.SNES#

class petsc4py.PETSc.SNES#

Bases: Object

Nonlinear equations solver.

SNES is described in the PETSc manual.

See also

SNES

Enumerations

ConvergedReason

SNES solver termination reason.

NewtonALCorrectionType

SNESNEWTONAL correction type.

NormSchedule

SNES norm schedule.

Type

SNES solver type.

Methods Summary

appendOptionsPrefix(prefix)

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

callConvergenceTest(its, xnorm, ynorm, fnorm)

Compute the convergence test.

computeFunction(x, f)

Compute the function.

computeJacobian(x, J[, P])

Compute the Jacobian.

computeNGS(x[, b])

Compute a nonlinear Gauss-Seidel step.

computeObjective(x)

Compute the value of the objective function.

converged(its, xnorm, ynorm, fnorm)

Compute the convergence test and update the solver converged reason.

create([comm])

Create a SNES solver.

createPython([context, comm])

Create a nonlinear solver of Python type.

destroy()

Destroy the solver.

getApplicationContext()

Return the application context.

getCompositeNumber()

Return the number of solvers in the composite.

getCompositeSNES(n)

Return the n-th solver in the composite.

getConvergedReason()

Return the termination flag.

getConvergenceHistory()

Return the convergence history.

getConvergenceTest()

Return the callback to used as convergence test.

getDM()

Return the DM associated with the solver.

getDivergenceTolerance()

Get the divergence tolerance parameter used in the convergence tests.

getErrorIfNotConverged()

Return the flag indicating error on divergence.

getFASCoarseSolve()

Return the SNES used at the coarsest level of the FAS hierarchy.

getFASCycleSNES(level)

Return the SNES corresponding to a particular level of the FAS hierarchy.

getFASInjection(level)

Return the Mat used to apply the injection from level-1 to level.

getFASInterpolation(level)

Return the Mat used to apply the interpolation from level-1 to level.

getFASLevels()

Return the number of levels used.

getFASRestriction(level)

Return the Mat used to apply the restriction from level-1 to level.

getFASSmoother(level)

Return the smoother used at a given level of the FAS hierarchy.

getFASSmootherDown(level)

Return the downsmoother used at a given level of the FAS hierarchy.

getFASSmootherUp(level)

Return the upsmoother used at a given level of the FAS hierarchy.

getFunction()

Return the callback to compute the nonlinear function.

getFunctionEvaluations()

Return the current number of function evaluations.

getFunctionNorm()

Return the function norm.

getInitialGuess()

Return the callback to compute the initial guess.

getIterationNumber()

Return the current iteration number.

getJacobian()

Return the matrices used to compute the Jacobian and the callback tuple.

getKSP()

Return the linear solver used by the nonlinear solver.

getKSPFailures()

Return the current number of linear solve failures.

getLineSearch()

Return the linesearch object associated with this SNES.

getLinearSolveIterations()

Return the total number of linear iterations.

getMaxFunctionEvaluations()

Return the maximum allowed number of function evaluations.

getMaxKSPFailures()

Return the maximum allowed number of linear solve failures.

getMaxStepFailures()

Return the maximum allowed number of step failures.

getMonitor()

Return the callback used to monitor solver convergence.

getNASMNumber()

Return the number of solvers in NASM.

getNASMSNES(n)

Return the n-th solver in NASM.

getNGS()

Return the nonlinear Gauss-Seidel callback tuple.

getNPC()

Return the nonlinear preconditioner associated with the solver.

getNPCSide()

Return the nonlinear preconditioning side.

getNewtonALLoadParameter()

Return the load parameter for SNESNEWTONAL.

getNormSchedule()

Return the norm schedule.

getObjective()

Return the objective callback tuple.

getOptionsPrefix()

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

getParamsEW()

Get the parameters of the Eisenstat and Walker trick.

getPythonContext()

Return the instance of the class implementing the required Python methods.

getPythonType()

Return the fully qualified Python name of the class used by the solver.

getRhs()

Return the vector holding the right-hand side.

getSolution()

Return the vector holding the solution.

getSolutionUpdate()

Return the vector holding the solution update.

getStepFailures()

Return the current number of step failures.

getTRTolerances()

Return the tolerance parameters used for the trust region.

getTRUpdateParameters()

Return the update parameters used for the trust region.

getTolerances()

Return the tolerance parameters used in the solver convergence tests.

getType()

Return the type of the solver.

getUpdate()

Return the callback to compute the update at the beginning of each step.

getUseEW()

Return the flag indicating if the solver uses the Eisenstat-Walker trick.

getUseFD()

Return true if the solver uses color finite-differencing for the Jacobian.

getUseKSP()

Return the flag indicating if the solver uses a linear solver.

getUseMF()

Return the flag indicating if the solver uses matrix-free finite-differencing.

getVIInactiveSet()

Return the index set for the inactive set.

getVariableBounds()

Get the vectors for the variable bounds.

hasNPC()

Return a boolean indicating whether the solver has a nonlinear preconditioner.

logConvergenceHistory(norm[, linear_its])

Log residual norm and linear iterations.

monitor(its, rnorm)

Monitor the solver.

monitorCancel()

Cancel all the monitors of the solver.

reset()

Reset the solver.

setApplicationContext(appctx)

Set the application context.

setConvergedReason(reason)

Set the termination flag.

setConvergenceHistory([length, reset])

Set the convergence history.

setConvergenceTest(converged[, args, kargs])

Set the callback to use as convergence test.

setDM(dm)

Associate a DM with the solver.

setDivergenceTolerance(dtol)

Set the divergence tolerance parameter used in the convergence tests.

setErrorIfNotConverged(flag)

Immediately generate an error if the solver has not converged.

setFASInjection(level, mat)

Set the Mat to be used to apply the injection from level-1 to level.

setFASInterpolation(level, mat)

Set the Mat to be used to apply the interpolation from level-1 to level.

setFASLevels(levels[, comms])

Set the number of levels to use with FAS.

setFASRScale(level, vec)

Set the scaling factor of the restriction operator from level to level-1.

setFASRestriction(level, mat)

Set the Mat to be used to apply the restriction from level-1 to level.

setForceIteration(force)

Force solve to take at least one iteration.

setFromOptions()

Configure the solver from the options database.

setFunction(function[, f, args, kargs])

Set the callback to compute the nonlinear function.

setFunctionNorm(norm)

Set the function norm value.

setInitialGuess(initialguess[, args, kargs])

Set the callback to compute the initial guess.

setIterationNumber(its)

Set the current iteration number.

setJacobian(jacobian[, J, P, args, kargs])

Set the callback to compute the Jacobian.

setKSP(ksp)

Set the linear solver that will be used by the nonlinear solver.

setLineSearch(linesearch)

Set the linesearch object to be used by this SNES.

setLineSearchPreCheck(precheck[, args, kargs])

Set the callback that will be called before applying the linesearch.

setMaxFunctionEvaluations(max_funcs)

Set the maximum allowed number of function evaluations.

setMaxKSPFailures(max_fails)

Set the maximum allowed number of linear solve failures.

setMaxStepFailures(max_fails)

Set the maximum allowed number of step failures.

setMonitor(monitor[, args, kargs])

Set the callback used to monitor solver convergence.

setNGS(ngs[, args, kargs])

Set the callback to compute nonlinear Gauss-Seidel.

setNPC(snes)

Set the nonlinear preconditioner.

setNPCSide(side)

Set the nonlinear preconditioning side.

setNewtonALCorrectionType(corrtype)

Set the correction type for SNESNEWTONAL.

setNewtonALFunction(function[, args, kargs])

Set the callback to compute the tangent load vector for SNESNEWTONAL.

setNormSchedule(normsched)

Set the norm schedule.

setObjective(objective[, args, kargs])

Set the callback to compute the objective function.

setOptionsPrefix(prefix)

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

setParamsEW([version, rtol_0, rtol_max, ...])

Set the parameters for the Eisenstat and Walker trick.

setPatchCellNumbering(sec)

Set cell patch numbering.

setPatchComputeFunction(function[, args, kargs])

Set patch compute function.

setPatchComputeOperator(operator[, args, kargs])

Set patch compute operator.

setPatchConstructType(typ[, operator, args, ...])

Set patch construct type.

setPatchDiscretisationInfo(dms, bs, ...)

Set patch discretisation information.

setPythonContext(context)

Set the instance of the class implementing the required Python methods.

setPythonType(py_type)

Set the fully qualified Python name of the class to be used.

setResetCounters([reset])

Set the flag to reset the counters.

setSolution(vec)

Set the vector used to store the solution.

setTRTolerances([delta_min, delta_max, delta_0])

Set the tolerance parameters used for the trust region.

setTRUpdateParameters([eta1, eta2, eta3, t1, t2])

Set the update parameters used for the trust region.

setTolerances([rtol, atol, stol, max_it])

Set the tolerance parameters used in the solver convergence tests.

setType(snes_type)

Set the type of the solver.

setUp()

Set up the internal data structures for using the solver.

setUpMatrices()

Ensures that matrices are available for Newton-like methods.

setUpdate(update[, args, kargs])

Set the callback to compute update at the beginning of each step.

setUseEW([flag])

Tell the solver to use the Eisenstat-Walker trick.

setUseFD([flag])

Set the boolean flag to use coloring finite-differencing for Jacobian assembly.

setUseKSP([flag])

Set the boolean flag indicating to use a linear solver.

setUseMF([flag])

Set the boolean flag indicating to use matrix-free finite-differencing.

setVariableBounds(xl, xu)

Set the vector for the variable bounds.

solve([b, x])

Solve the nonlinear equations.

view([viewer])

View the solver.

Attributes Summary

appctx

Application context.

atol

Absolute residual tolerance.

dm

DM.

history

Convergence history.

is_converged

Boolean indicating if the solver has converged.

is_diverged

Boolean indicating if the solver has failed.

is_iterating

Boolean indicating if the solver has not converged yet.

its

Number of iterations.

ksp

Linear solver.

linesearch

The linesearch object associated with this SNES.

max_funcs

Maximum number of function evaluations.

max_it

Maximum number of iterations.

norm

Function norm.

npc

Nonlinear preconditioner.

reason

Converged reason.

rtol

Relative residual tolerance.

stol

Solution update tolerance.

use_ew

Use the Eisenstat-Walker trick.

use_fd

Boolean indicating if the solver uses coloring finite-differencing.

use_ksp

Boolean indicating if the solver uses a linear solver.

use_mf

Boolean indicating if the solver uses matrix-free finite-differencing.

vec_rhs

Right-hand side vector.

vec_sol

Solution vector.

vec_upd

Update vector.

Methods Documentation

appendOptionsPrefix(prefix)#

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

Logically collective.

Source code at petsc4py/PETSc/SNES.pyx:243

Parameters:

prefix (str | None)

Return type:

None

callConvergenceTest(its, xnorm, ynorm, fnorm)#

Compute the convergence test.

Collective.

Parameters:
  • its (int) – Iteration number.

  • xnorm (float) – Solution norm.

  • ynorm (float) – Update norm.

  • fnorm (float) – Function norm.

Return type:

ConvergedReason

Source code at petsc4py/PETSc/SNES.pyx:1342

computeFunction(x, f)#

Compute the function.

Collective.

Parameters:
  • x (Vec) – The input state vector.

  • f (Vec) – The output vector.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:1043

computeJacobian(x, J, P=None)#

Compute the Jacobian.

Collective.

Parameters:
  • x (Vec) – The input state vector.

  • J (Mat) – The output Jacobian matrix.

  • P (Mat | None) – The output Jacobian matrix used to construct the preconditioner.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:1062

computeNGS(x, b=None)#

Compute a nonlinear Gauss-Seidel step.

Collective.

Parameters:
  • x (Vec) – The input/output state vector.

  • b (Vec | None) – The input right-hand side vector.

Return type:

None

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

computeObjective(x)#

Compute the value of the objective function.

Collective.

Parameters:

x (Vec) – The input state vector.

Return type:

float

Source code at petsc4py/PETSc/SNES.pyx:1085

converged(its, xnorm, ynorm, fnorm)#

Compute the convergence test and update the solver converged reason.

Collective.

Parameters:
  • its (int) – Iteration number.

  • xnorm (float) – Solution norm.

  • ynorm (float) – Update norm.

  • fnorm (float) – Function norm.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:1372

create(comm=None)#

Create a SNES solver.

Collective.

Parameters:

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

Return type:

Self

Source code at petsc4py/PETSc/SNES.pyx:161

createPython(context=None, comm=None)#

Create a nonlinear solver of Python type.

Collective.

Parameters:
  • context (Any) – An instance of the Python class implementing the required methods.

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

Return type:

Self

Source code at petsc4py/PETSc/SNES.pyx:2210

destroy()#

Destroy the solver.

Collective.

See also

SNESDestroy

Source code at petsc4py/PETSc/SNES.pyx:148

Return type:

Self

getApplicationContext()#

Return the application context.

Source code at petsc4py/PETSc/SNES.pyx:280

Return type:

Any

getCompositeNumber()#

Return the number of solvers in the composite.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:2311

Return type:

int

getCompositeSNES(n)#

Return the n-th solver in the composite.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:2294

Parameters:

n (int)

Return type:

SNES

getConvergedReason()#

Return the termination flag.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1755

Return type:

ConvergedReason

getConvergenceHistory()#

Return the convergence history.

Not collective.

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

Return type:

tuple[ArrayReal, ArrayInt]

getConvergenceTest()#

Return the callback to used as convergence test.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1330

Return type:

SNESConvergedFunction

getDM()#

Return the DM associated with the solver.

Not collective.

See also

setDM, SNESGetDM

Source code at petsc4py/PETSc/SNES.pyx:295

Return type:

DM

getDivergenceTolerance()#

Get the divergence tolerance parameter used in the convergence tests.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1255

Return type:

float

getErrorIfNotConverged()#

Return the flag indicating error on divergence.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1782

Return type:

bool

getFASCoarseSolve()#

Return the SNES used at the coarsest level of the FAS hierarchy.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:610

Return type:

SNES

getFASCycleSNES(level)#

Return the SNES corresponding to a particular level of the FAS hierarchy.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:593

Parameters:

level (int)

Return type:

SNES

getFASInjection(level)#

Return the Mat used to apply the injection from level-1 to level.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:516

Parameters:

level (int)

Return type:

Mat

getFASInterpolation(level)#

Return the Mat used to apply the interpolation from level-1 to level.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:456

Parameters:

level (int)

Return type:

Mat

getFASLevels()#

Return the number of levels used.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:579

Return type:

int

getFASRestriction(level)#

Return the Mat used to apply the restriction from level-1 to level.

Not collective.

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

Parameters:

level (int)

Return type:

Mat

getFASSmoother(level)#

Return the smoother used at a given level of the FAS hierarchy.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:625

Parameters:

level (int)

Return type:

SNES

getFASSmootherDown(level)#

Return the downsmoother used at a given level of the FAS hierarchy.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:642

Parameters:

level (int)

Return type:

SNES

getFASSmootherUp(level)#

Return the upsmoother used at a given level of the FAS hierarchy.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:659

Parameters:

level (int)

Return type:

SNES

getFunction()#

Return the callback to compute the nonlinear function.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:859

Return type:

SNESFunction

getFunctionEvaluations()#

Return the current number of function evaluations.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1576

Return type:

int

getFunctionNorm()#

Return the function norm.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1853

Return type:

float

getInitialGuess()#

Return the callback to compute the initial guess.

Not collective.

See also

setInitialGuess

Source code at petsc4py/PETSc/SNES.pyx:812

Return type:

SNESGuessFunction

getIterationNumber()#

Return the current iteration number.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1811

Return type:

int

getJacobian()#

Return the matrices used to compute the Jacobian and the callback tuple.

Not collective.

Returns:

  • J (Mat) – The matrix to store the Jacobian.

  • P (Mat) – The matrix to construct the preconditioner.

  • callback (SNESJacobianFunction) – callback, positional and keyword arguments.

Return type:

tuple[Mat, Mat, SNESJacobianFunction]

Source code at petsc4py/PETSc/SNES.pyx:971

getKSP()#

Return the linear solver used by the nonlinear solver.

Not collective.

See also

setKSP, SNESGetKSP

Source code at petsc4py/PETSc/SNES.pyx:1952

Return type:

KSP

getKSPFailures()#

Return the current number of linear solve failures.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1658

Return type:

int

getLineSearch()#

Return the linesearch object associated with this SNES.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:2619

Return type:

SNESLineSearch

getLinearSolveIterations()#

Return the total number of linear iterations.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1867

Return type:

int

getMaxFunctionEvaluations()#

Return the maximum allowed number of function evaluations.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1560

Return type:

int

getMaxKSPFailures()#

Return the maximum allowed number of linear solve failures.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1644

Return type:

int

getMaxStepFailures()#

Return the maximum allowed number of step failures.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1603

Return type:

int

getMonitor()#

Return the callback used to monitor solver convergence.

Not collective.

See also

setMonitor

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

Return type:

list[tuple[SNESMonitorFunction, tuple[Any, …], dict[str, Any]]]

getNASMNumber()#

Return the number of solvers in NASM.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:2343

Return type:

int

getNASMSNES(n)#

Return the n-th solver in NASM.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:2327

Parameters:

n (int)

Return type:

SNES

getNGS()#

Return the nonlinear Gauss-Seidel callback tuple.

Not collective.

See also

setNGS, computeNGS

Source code at petsc4py/PETSc/SNES.pyx:1135

Return type:

SNESNGSFunction

getNPC()#

Return the nonlinear preconditioner associated with the solver.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:678

Return type:

SNES

getNPCSide()#

Return the nonlinear preconditioning side.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:731

Return type:

Side

getNewtonALLoadParameter()#

Return the load parameter for SNESNEWTONAL.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:2687

Return type:

float

getNormSchedule()#

Return the norm schedule.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1281

Return type:

NormSchedule

getObjective()#

Return the objective callback tuple.

Not collective.

See also

setObjective

Source code at petsc4py/PETSc/SNES.pyx:1029

Return type:

SNESObjFunction

getOptionsPrefix()#

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

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:229

Return type:

str

getParamsEW()#

Get the parameters of the Eisenstat and Walker trick.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:2056

Return type:

dict[str, int | float]

getPythonContext()#

Return the instance of the class implementing the required Python methods.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:2247

Return type:

Any

getPythonType()#

Return the fully qualified Python name of the class used by the solver.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:2277

Return type:

str

getRhs()#

Return the vector holding the right-hand side.

Not collective.

See also

SNESGetRhs

Source code at petsc4py/PETSc/SNES.pyx:1881

Return type:

Vec

getSolution()#

Return the vector holding the solution.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1896

Return type:

Vec

getSolutionUpdate()#

Return the vector holding the solution update.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1923

Return type:

Vec

getStepFailures()#

Return the current number of step failures.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1617

Return type:

int

getTRTolerances()#

Return the tolerance parameters used for the trust region.

Not collective.

Returns:

  • delta_min (float) – The minimum allowed trust region size.

  • delta_max (float) – The maximum allowed trust region size.

  • delta_0 (float) – The initial trust region size.

Return type:

tuple[float, float, float]

Source code at petsc4py/PETSc/SNES.pyx:353

getTRUpdateParameters()#

Return the update parameters used for the trust region.

Not collective.

Returns:

  • eta1 (float) – The step acceptance tolerance.

  • eta2 (float) – The shrinking tolerance.

  • eta3 (float) – The enlarging tolerance.

  • t1 (float) – The shrinking factor.

  • t2 (float) – The enlarging factor.

Return type:

tuple[float, float, float, float, float]

Source code at petsc4py/PETSc/SNES.pyx:412

getTolerances()#

Return the tolerance parameters used in the solver convergence tests.

Not collective.

Returns:

  • rtol (float) – The relative norm of the residual.

  • atol (float) – The absolute norm of the residual.

  • stol (float) – The absolute norm of the step.

  • max_it (int) – The maximum allowed number of iterations.

Return type:

tuple[float, float, float, int]

Source code at petsc4py/PETSc/SNES.pyx:1211

getType()#

Return the type of the solver.

Not collective.

See also

setType, SNESGetType

Source code at petsc4py/PETSc/SNES.pyx:201

Return type:

str

getUpdate()#

Return the callback to compute the update at the beginning of each step.

Not collective.

See also

setUpdate

Source code at petsc4py/PETSc/SNES.pyx:920

Return type:

SNESUpdateFunction

getUseEW()#

Return the flag indicating if the solver uses the Eisenstat-Walker trick.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:1990

Return type:

bool

getUseFD()#

Return true if the solver uses color finite-differencing for the Jacobian.

Not collective.

See also

setUseFD

Source code at petsc4py/PETSc/SNES.pyx:2150

Return type:

bool

getUseKSP()#

Return the flag indicating if the solver uses a linear solver.

Not collective.

See also

setUseKSP

Source code at petsc4py/PETSc/SNES.pyx:2094

Return type:

bool

getUseMF()#

Return the flag indicating if the solver uses matrix-free finite-differencing.

Not collective.

See also

setUseMF

Source code at petsc4py/PETSc/SNES.pyx:2123

Return type:

bool

getVIInactiveSet()#

Return the index set for the inactive set.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:2193

Return type:

IS

getVariableBounds()#

Get the vectors for the variable bounds.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:2178

Return type:

tuple[Vec, Vec]

hasNPC()#

Return a boolean indicating whether the solver has a nonlinear preconditioner.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:693

Return type:

bool

logConvergenceHistory(norm, linear_its=0)#

Log residual norm and linear iterations.

Source code at petsc4py/PETSc/SNES.pyx:1441

Parameters:
Return type:

None

monitor(its, rnorm)#

Monitor the solver.

Collective.

Parameters:
  • its – Current number of iterations.

  • rnorm – Current value of the residual norm.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:1522

monitorCancel()#

Cancel all the monitors of the solver.

Logically collective.

Source code at petsc4py/PETSc/SNES.pyx:1507

Return type:

None

reset()#

Reset the solver.

Collective.

See also

SNESReset

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

Return type:

None

setApplicationContext(appctx)#

Set the application context.

Source code at petsc4py/PETSc/SNES.pyx:271

Parameters:

appctx (Any)

Return type:

None

setConvergedReason(reason)#

Set the termination flag.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:1742

Parameters:

reason (ConvergedReason)

Return type:

None

setConvergenceHistory(length=None, reset=False)#

Set the convergence history.

Logically collective.

Source code at petsc4py/PETSc/SNES.pyx:1399

Return type:

None

setConvergenceTest(converged, args=None, kargs=None)#

Set the callback to use as convergence test.

Logically collective.

Parameters:
  • converged (SNESConvergedFunction | Literal['skip', 'default']) – The convergence testing callback.

  • args (tuple[Any, ...] | None) – Positional arguments for the callback.

  • kargs (dict[str, Any] | None) – Keyword arguments for the callback.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:1295

setDM(dm)#

Associate a DM with the solver.

Not collective.

See also

getDM, SNESSetDM

Source code at petsc4py/PETSc/SNES.pyx:312

Parameters:

dm (DM)

Return type:

None

setDivergenceTolerance(dtol)#

Set the divergence tolerance parameter used in the convergence tests.

Logically collective.

Parameters:

dtol (float) – The divergence tolerance.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:1238

setErrorIfNotConverged(flag)#

Immediately generate an error if the solver has not converged.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:1769

Parameters:

flag (bool)

Return type:

None

setFASInjection(level, mat)#

Set the Mat to be used to apply the injection from level-1 to level.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:502

Parameters:
Return type:

None

setFASInterpolation(level, mat)#

Set the Mat to be used to apply the interpolation from level-1 to level.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:442

Parameters:
Return type:

None

setFASLevels(levels, comms=None)#

Set the number of levels to use with FAS.

Collective.

Parameters:
Return type:

None

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

setFASRScale(level, vec)#

Set the scaling factor of the restriction operator from level to level-1.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:532

Parameters:
Return type:

None

setFASRestriction(level, mat)#

Set the Mat to be used to apply the restriction from level-1 to level.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:472

Parameters:
Return type:

None

setForceIteration(force)#

Force solve to take at least one iteration.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:1825

Parameters:

force (bool)

Return type:

None

setFromOptions()#

Configure the solver from the options database.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:257

Return type:

None

setFunction(function, f=None, args=None, kargs=None)#

Set the callback to compute the nonlinear function.

Logically collective.

Parameters:
  • function (SNESFunction | None) – The callback.

  • f (Vec | None) – An optional vector to store the result.

  • args (tuple[Any, ...] | None) – Positional arguments for the callback.

  • kargs (dict[str, Any] | None) – Keyword arguments for the callback.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:824

setFunctionNorm(norm)#

Set the function norm value.

Collective.

This is only of use to implementers of custom SNES types.

Source code at petsc4py/PETSc/SNES.pyx:1838

Parameters:

norm (float)

Return type:

None

setInitialGuess(initialguess, args=None, kargs=None)#

Set the callback to compute the initial guess.

Logically collective.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:781

setIterationNumber(its)#

Set the current iteration number.

Collective.

This is only of use to implementers of custom SNES types.

Source code at petsc4py/PETSc/SNES.pyx:1796

Parameters:

its (int)

Return type:

None

setJacobian(jacobian, J=None, P=None, args=None, kargs=None)#

Set the callback to compute the Jacobian.

Logically collective.

Parameters:
  • jacobian (SNESJacobianFunction | None) – The Jacobian callback.

  • J (Mat | None) – The matrix to store the Jacobian.

  • P (Mat | None) – The matrix to construct the preconditioner.

  • args (tuple[Any, ...] | None) – Positional arguments for the callback.

  • kargs (dict[str, Any] | None) – Keyword arguments for the callback.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:932

setKSP(ksp)#

Set the linear solver that will be used by the nonlinear solver.

Logically collective.

See also

getKSP, SNESSetKSP

Source code at petsc4py/PETSc/SNES.pyx:1940

Parameters:

ksp (KSP)

Return type:

None

setLineSearch(linesearch)#

Set the linesearch object to be used by this SNES.

Logically collective.

Source code at petsc4py/PETSc/SNES.pyx:2634

Parameters:

linesearch (SNESLineSearch)

Return type:

None

setLineSearchPreCheck(precheck, args=None, kargs=None)#

Set the callback that will be called before applying the linesearch.

Logically collective.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:747

setMaxFunctionEvaluations(max_funcs)#

Set the maximum allowed number of function evaluations.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:1545

Parameters:

max_funcs (int)

Return type:

None

setMaxKSPFailures(max_fails)#

Set the maximum allowed number of linear solve failures.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:1631

Parameters:

max_fails (int)

Return type:

None

setMaxStepFailures(max_fails)#

Set the maximum allowed number of step failures.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:1590

Parameters:

max_fails (int)

Return type:

None

setMonitor(monitor, args=None, kargs=None)#

Set the callback used to monitor solver convergence.

Logically collective.

Parameters:
Return type:

None

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

setNGS(ngs, args=None, kargs=None)#

Set the callback to compute nonlinear Gauss-Seidel.

Logically collective.

Parameters:
  • ngs (SNESNGSFunction | None) – The nonlinear Gauss-Seidel callback.

  • args (tuple[Any, ...] | None) – Positional arguments for the callback.

  • kargs (dict[str, Any] | None) – Keyword arguments for the callback.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:1104

setNPC(snes)#

Set the nonlinear preconditioner.

Logically collective.

Source code at petsc4py/PETSc/SNES.pyx:707

Parameters:

snes (SNES)

Return type:

None

setNPCSide(side)#

Set the nonlinear preconditioning side.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:719

Parameters:

side (Side)

Return type:

None

setNewtonALCorrectionType(corrtype)#

Set the correction type for SNESNEWTONAL.

Logically collective.

Source code at petsc4py/PETSc/SNES.pyx:2702

Parameters:

corrtype (NewtonALCorrectionType)

Return type:

None

setNewtonALFunction(function, args=None, kargs=None)#

Set the callback to compute the tangent load vector for SNESNEWTONAL.

Logically collective.

Parameters:
  • function (SNESFunction | None) – The callback.

  • args (tuple[Any, ...] | None) – Positional arguments for the callback.

  • kargs (dict[str, Any] | None) – Keyword arguments for the callback.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:2656

setNormSchedule(normsched)#

Set the norm schedule.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:1269

Parameters:

normsched (NormSchedule)

Return type:

None

setObjective(objective, args=None, kargs=None)#

Set the callback to compute the objective function.

Logically collective.

Parameters:
  • objective (SNESObjFunction | None) – The Jacobian callback.

  • args (tuple[Any, ...] | None) – Positional arguments for the callback.

  • kargs (dict[str, Any] | None) – Keyword arguments for the callback.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:998

setOptionsPrefix(prefix)#

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

Logically collective.

Source code at petsc4py/PETSc/SNES.pyx:215

Parameters:

prefix (str | None)

Return type:

None

setParamsEW(version=None, rtol_0=None, rtol_max=None, gamma=None, alpha=None, alpha2=None, threshold=None)#

Set the parameters for the Eisenstat and Walker trick.

Logically collective.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:2004

setPatchCellNumbering(sec)#

Set cell patch numbering.

Source code at petsc4py/PETSc/SNES.pyx:2359

Parameters:

sec (Section)

Return type:

None

setPatchComputeFunction(function, args=None, kargs=None)#

Set patch compute function.

Source code at petsc4py/PETSc/SNES.pyx:2412

Return type:

None

setPatchComputeOperator(operator, args=None, kargs=None)#

Set patch compute operator.

Source code at petsc4py/PETSc/SNES.pyx:2404

Return type:

None

setPatchConstructType(typ, operator=None, args=None, kargs=None)#

Set patch construct type.

Source code at petsc4py/PETSc/SNES.pyx:2420

Return type:

None

setPatchDiscretisationInfo(dms, bs, cellNodeMaps, subspaceOffsets, ghostBcNodes, globalBcNodes)#

Set patch discretisation information.

Source code at petsc4py/PETSc/SNES.pyx:2363

Return type:

None

setPythonContext(context)#

Set the instance of the class implementing the required Python methods.

Not collective.

Source code at petsc4py/PETSc/SNES.pyx:2235

Parameters:

context (Any)

Return type:

None

setPythonType(py_type)#

Set the fully qualified Python name of the class to be used.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:2262

Parameters:

py_type (str)

Return type:

None

setResetCounters(reset=True)#

Set the flag to reset the counters.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:1447

Parameters:

reset (bool)

Return type:

None

setSolution(vec)#

Set the vector used to store the solution.

Collective.

Source code at petsc4py/PETSc/SNES.pyx:1911

Parameters:

vec (Vec)

Return type:

None

setTRTolerances(delta_min=None, delta_max=None, delta_0=None)#

Set the tolerance parameters used for the trust region.

Logically collective.

Parameters:
  • delta_min (float | None) – The minimum allowed trust region size. Defaults to CURRENT.

  • delta_max (float | None) – The maximum allowed trust region size. Defaults to CURRENT.

  • delta_0 (float | None) – The initial trust region size. Defaults to CURRENT.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:326

setTRUpdateParameters(eta1=None, eta2=None, eta3=None, t1=None, t2=None)#

Set the update parameters used for the trust region.

Logically collective.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:377

setTolerances(rtol=None, atol=None, stol=None, max_it=None)#

Set the tolerance parameters used in the solver convergence tests.

Logically collective.

Parameters:
  • rtol (float | None) – The relative norm of the residual, or DETERMINE to use the value when the object’s type was set. Defaults to CURRENT.

  • atol (float | None) – The absolute norm of the residual, or DETERMINE to use the value when the object’s type was set. Defaults to CURRENT.

  • stol (float | None) – The absolute norm of the step, or DETERMINE to use the value when the object’s type was set. Defaults to CURRENT.

  • max_it (int | None) – The maximum allowed number of iterations, or DETERMINE to use the value when the object’s type was set. Defaults to CURRENT. Use UNLIMITED to have no maximum.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:1172

setType(snes_type)#

Set the type of the solver.

Logically collective.

Parameters:

snes_type (Type | str) – The type of the solver.

Return type:

None

See also

getType, SNESSetType

Source code at petsc4py/PETSc/SNES.pyx:182

setUp()#

Set up the internal data structures for using the solver.

Collective.

See also

SNESSetUp

Source code at petsc4py/PETSc/SNES.pyx:1681

Return type:

None

setUpMatrices()#

Ensures that matrices are available for Newton-like methods.

Collective.

This is only of use to implementers of custom SNES types.

Source code at petsc4py/PETSc/SNES.pyx:1693

Return type:

None

setUpdate(update, args=None, kargs=None)#

Set the callback to compute update at the beginning of each step.

Logically collective.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:889

setUseEW(flag=True, *targs, **kargs)#

Tell the solver to use the Eisenstat-Walker trick.

Logically collective.

Parameters:
  • flag (bool) – Whether or not to use the Eisenstat-Walker trick.

  • *targs (Any) – Positional arguments for setParamsEW.

  • **kargs (Any) – Keyword arguments for setParamsEW.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:1967

setUseFD(flag=True)#

Set the boolean flag to use coloring finite-differencing for Jacobian assembly.

Logically collective.

See also

getUseFD

Source code at petsc4py/PETSc/SNES.pyx:2137

Return type:

None

setUseKSP(flag=True)#

Set the boolean flag indicating to use a linear solver.

Logically collective.

See also

getUseKSP

Source code at petsc4py/PETSc/SNES.pyx:2081

Return type:

None

setUseMF(flag=True)#

Set the boolean flag indicating to use matrix-free finite-differencing.

Logically collective.

See also

getUseMF

Source code at petsc4py/PETSc/SNES.pyx:2110

Return type:

None

setVariableBounds(xl, xu)#

Set the vector for the variable bounds.

Collective.

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

Parameters:
Return type:

None

solve(b=None, x=None)#

Solve the nonlinear equations.

Collective.

Parameters:
  • b (Vec | None) – The affine right-hand side or None to use zero.

  • x (Vec | None) – The starting vector or None to use the vector stored internally.

Return type:

None

Source code at petsc4py/PETSc/SNES.pyx:1719

view(viewer=None)#

View the solver.

Collective.

Parameters:

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

Return type:

None

See also

SNESView

Source code at petsc4py/PETSc/SNES.pyx:129

Attributes Documentation

appctx#

Application context.

Source code at petsc4py/PETSc/SNES.pyx:2436

atol#

Absolute residual tolerance.

Source code at petsc4py/PETSc/SNES.pyx:2517

dm#

DM.

Source code at petsc4py/PETSc/SNES.pyx:2446

history#

Convergence history.

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

is_converged#

Boolean indicating if the solver has converged.

Source code at petsc4py/PETSc/SNES.pyx:2589

is_diverged#

Boolean indicating if the solver has failed.

Source code at petsc4py/PETSc/SNES.pyx:2594

is_iterating#

Boolean indicating if the solver has not converged yet.

Source code at petsc4py/PETSc/SNES.pyx:2584

its#

Number of iterations.

Source code at petsc4py/PETSc/SNES.pyx:2553

ksp#

Linear solver.

Source code at petsc4py/PETSc/SNES.pyx:2483

linesearch#

The linesearch object associated with this SNES.

Source code at petsc4py/PETSc/SNES.pyx:2646

max_funcs#

Maximum number of function evaluations.

Source code at petsc4py/PETSc/SNES.pyx:2543

max_it#

Maximum number of iterations.

Source code at petsc4py/PETSc/SNES.pyx:2533

norm#

Function norm.

Source code at petsc4py/PETSc/SNES.pyx:2561

npc#

Nonlinear preconditioner.

Source code at petsc4py/PETSc/SNES.pyx:2456

reason#

Converged reason.

Source code at petsc4py/PETSc/SNES.pyx:2576

rtol#

Relative residual tolerance.

Source code at petsc4py/PETSc/SNES.pyx:2509

stol#

Solution update tolerance.

Source code at petsc4py/PETSc/SNES.pyx:2525

use_ew#

Use the Eisenstat-Walker trick.

Source code at petsc4py/PETSc/SNES.pyx:2499

use_fd#

Boolean indicating if the solver uses coloring finite-differencing.

Source code at petsc4py/PETSc/SNES.pyx:2609

use_ksp#

Boolean indicating if the solver uses a linear solver.

Source code at petsc4py/PETSc/SNES.pyx:2491

use_mf#

Boolean indicating if the solver uses matrix-free finite-differencing.

Source code at petsc4py/PETSc/SNES.pyx:2601

vec_rhs#

Right-hand side vector.

Source code at petsc4py/PETSc/SNES.pyx:2476

vec_sol#

Solution vector.

Source code at petsc4py/PETSc/SNES.pyx:2466

vec_upd#

Update vector.

Source code at petsc4py/PETSc/SNES.pyx:2471