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.

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.

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.

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.

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.

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.

getUseMF()

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

getVIInactiveSet()

Return the index set for the inactive set.

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.

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.

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.

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.

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.

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.

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_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:227

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:1171

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:911

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:930

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:1017

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:953

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:1201

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:145

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:1997

destroy()#

Destroy the solver.

Collective.

See also

SNESDestroy

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

Return type:

Self

getApplicationContext()#

Return the application context.

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

Return type:

Any

getCompositeNumber()#

Return the number of solvers in the composite.

Not collective.

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

Return type:

int

getCompositeSNES(n)#

Return the n-th solver in the composite.

Not collective.

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

Parameters:

n (int)

Return type:

SNES

getConvergedReason()#

Return the termination flag.

Not collective.

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

Return type:

ConvergedReason

getConvergenceHistory()#

Return the convergence history.

Not collective.

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

Return type:

tuple[ArrayReal, ArrayInt]

getConvergenceTest()#

Return the callback to used as convergence test.

Not collective.

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

Return type:

SNESConvergedFunction

getDM()#

Return the DM associated with the solver.

Not collective.

See also

setDM, SNESGetDM

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

Return type:

DM

getErrorIfNotConverged()#

Return the flag indicating error on divergence.

Not collective.

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

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:478

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:461

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:384

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:324

Parameters:

level (int)

Return type:

Mat

getFASLevels()#

Return the number of levels used.

Not collective.

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

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:354

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:493

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:510

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:527

Parameters:

level (int)

Return type:

SNES

getFunction()#

Return the callback to compute the nonlinear function.

Not collective.

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

Return type:

SNESFunction

getFunctionEvaluations()#

Return the current number of function evaluations.

Not collective.

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

Return type:

int

getFunctionNorm()#

Return the function norm.

Not collective.

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

Return type:

float

getInitialGuess()#

Return the callback to compute the initial guess.

Not collective.

See also

setInitialGuess

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

Return type:

SNESGuessFunction

getIterationNumber()#

Return the current iteration number.

Not collective.

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

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:839

getKSP()#

Return the linear solver used by the nonlinear solver.

Not collective.

See also

setKSP, SNESGetKSP

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

Return type:

KSP

getKSPFailures()#

Return the current number of linear solve failures.

Not collective.

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

Return type:

int

getLinearSolveIterations()#

Return the total number of linear iterations.

Not collective.

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

Return type:

int

getMaxFunctionEvaluations()#

Return the maximum allowed number of function evaluations.

Not collective.

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

Return type:

int

getMaxKSPFailures()#

Return the maximum allowed number of linear solve failures.

Not collective.

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

Return type:

int

getMaxStepFailures()#

Return the maximum allowed number of step failures.

Not collective.

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

Return type:

int

getMonitor()#

Return the callback used to monitor solver convergence.

Not collective.

See also

setMonitor

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

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:2130

Return type:

int

getNASMSNES(n)#

Return the n-th solver in NASM.

Not collective.

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

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:1003

Return type:

SNESNGSFunction

getNPC()#

Return the nonlinear preconditioner associated with the solver.

Not collective.

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

Return type:

SNES

getNPCSide()#

Return the nonlinear preconditioning side.

Not collective.

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

Return type:

Side

getNormSchedule()#

Return the norm schedule.

Not collective.

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

Return type:

NormSchedule

getObjective()#

Return the objective callback tuple.

Not collective.

See also

setObjective

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

Return type:

SNESObjFunction

getOptionsPrefix()#

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

Not collective.

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

Return type:

str

getParamsEW()#

Get the parameters of the Eisenstat and Walker trick.

Not collective.

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

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:2034

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:2064

Return type:

str

getRhs()#

Return the vector holding the right-hand side.

Not collective.

See also

SNESGetRhs

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

Return type:

Vec

getSolution()#

Return the vector holding the solution.

Not collective.

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

Return type:

Vec

getSolutionUpdate()#

Return the vector holding the solution update.

Not collective.

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

Return type:

Vec

getStepFailures()#

Return the current number of step failures.

Not collective.

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

Return type:

int

getTolerances()#

Return the tolerance parameters used in the solver convergence tests.

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:1071

getType()#

Return the type of the solver.

Not collective.

See also

setType, SNESGetType

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

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:788

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:1819

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:1952

Return type:

False

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:1925

Return type:

bool

getVIInactiveSet()#

Return the index set for the inactive set.

Not collective.

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

Return type:

IS

hasNPC()#

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

Not collective.

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

Return type:

bool

logConvergenceHistory(norm, linear_its=0)#

Log residual norm and linear iterations.

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

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:1351

monitorCancel()#

Cancel all the monitors of the solver.

Logically collective.

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

Return type:

None

reset()#

Reset the solver.

Collective.

See also

SNESReset

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

Return type:

None

setApplicationContext(appctx)#

Set the application context.

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

Parameters:

appctx (Any)

Return type:

None

setConvergedReason(reason)#

Set the termination flag.

Collective.

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

Parameters:

reason (ConvergedReason)

Return type:

None

setConvergenceHistory(length=None, reset=False)#

Set the convergence history.

Logically collective.

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

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:1124

setDM(dm)#

Associate a DM with the solver.

Not collective.

See also

getDM, SNESSetDM

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

Parameters:

dm (DM)

Return type:

None

setErrorIfNotConverged(flag)#

Immediately generate an error if the solver has not converged.

Collective.

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

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:370

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:310

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:413

setFASRScale(level, vec)#

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

Collective.

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

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:340

Parameters:
Return type:

None

setForceIteration(force)#

Force solve to take at least one iteration.

Collective.

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

Parameters:

force (bool)

Return type:

None

setFromOptions()#

Configure the solver from the options database.

Collective.

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

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:692

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:1667

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:649

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:1625

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:800

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:1769

Parameters:

ksp (KSP)

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:615

setMaxFunctionEvaluations(max_funcs)#

Set the maximum allowed number of function evaluations.

Collective.

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

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:1460

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:1419

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:1291

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:972

setNPC(snes)#

Set the nonlinear preconditioner.

Logically collective.

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

Parameters:

snes (SNES)

Return type:

None

setNPCSide(side)#

Set the nonlinear preconditioning side.

Collective.

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

Parameters:

side (Side)

Return type:

None

setNormSchedule(normsched)#

Set the norm schedule.

Collective.

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

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:866

setOptionsPrefix(prefix)#

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

Logically collective.

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

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:
  • version (int) – The version of the algorithm. Defaults to DEFAULT.

  • rtol_0 (float) – The initial relative residual norm. Defaults to DEFAULT.

  • rtol_max (float) – The maximum relative residual norm. Defaults to DEFAULT.

  • gamma (float) – Parameter. Defaults to DEFAULT.

  • alpha (float) – Parameter. Defaults to DEFAULT.

  • alpha2 (float) – Parameter. Defaults to DEFAULT.

  • threshold (float) – Parameter. Defaults to DEFAULT.

Return type:

None

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

setPatchCellNumbering(sec)#

Set cell patch numbering.

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

Parameters:

sec (Section)

Return type:

None

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

Set patch compute function.

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

Return type:

None

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

Set patch compute operator.

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

Return type:

None

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

Set patch construct type.

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

Return type:

None

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

Set patch discretisation information.

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

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:2022

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:2049

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:1276

Parameters:

reset (bool)

Return type:

None

setSolution(vec)#

Set the vector used to store the solution.

Collective.

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

Parameters:

vec (Vec)

Return type:

None

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

Set the tolerance parameters used in the solver convergence tests.

Collective.

Parameters:
  • rtol (float) – The relative norm of the residual. Defaults to DEFAULT.

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

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

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

Return type:

None

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

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:166

setUp()#

Set up the internal data structures for using the solver.

Collective.

See also

SNESSetUp

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

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:1522

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:757

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:1796

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:1939

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:1912

Return type:

None

setVariableBounds(xl, xu)#

Set the vector for the variable bounds.

Collective.

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

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:1548

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:113

Attributes Documentation

appctx#

Application context.

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

atol#

Absolute residual tolerance.

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

dm#

DM.

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

history#

Convergence history.

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

is_converged#

Boolean indicating if the solver has converged.

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

is_diverged#

Boolean indicating if the solver has failed.

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

is_iterating#

Boolean indicating if the solver has not converged yet.

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

its#

Number of iterations.

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

ksp#

Linear solver.

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

max_funcs#

Maximum number of function evaluations.

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

max_it#

Maximum number of iterations.

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

norm#

Function norm.

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

npc#

Nonlinear preconditioner.

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

reason#

Converged reason.

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

rtol#

Relative residual tolerance.

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

stol#

Solution update tolerance.

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

use_ew#

Use the Eisenstat-Walker trick.

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

use_fd#

Boolean indicating if the solver uses coloring finite-differencing.

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

use_mf#

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

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

vec_rhs#

Right-hand side vector.

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

vec_sol#

Solution vector.

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

vec_upd#

Update vector.

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