petsc4py.PETSc.SNES
- class petsc4py.PETSc.SNES
Bases:
Object
Nonlinear equations solver.
SNES is described in the
PETSc manual
.See also
petsc.SNES
Enumerations
SNES solver termination reason.
SNES norm schedule.
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.
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.
Return the application context.
Return the number of solvers in the composite.
Return the n-th solver in the composite.
Return the termination flag.
Return the convergence history.
Return the callback to used as convergence test.
getDM
()Return the
DM
associated with the solver.Return the flag indicating error on divergence.
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.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.
Return the callback to compute the nonlinear function.
Return the current number of function evaluations.
Return the function norm.
Return the callback to compute the initial guess.
Return the current iteration number.
Return the matrices used to compute the Jacobian and the callback tuple.
getKSP
()Return the linear solver used by the nonlinear solver.
Return the current number of linear solve failures.
Return the total number of linear iterations.
Return the maximum allowed number of function evaluations.
Return the maximum allowed number of linear solve failures.
Return the maximum allowed number of step failures.
Return the callback used to monitor solver convergence.
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.
Return the nonlinear preconditioning side.
Return the norm schedule.
Return the objective callback tuple.
Return the prefix used for searching for options in the database.
Get the parameters of the Eisenstat and Walker trick.
Return the instance of the class implementing the required Python methods.
Return the fully qualified Python name of the class used by the solver.
getRhs
()Return the vector holding the right-hand side.
Return the vector holding the solution.
Return the vector holding the solution update.
Return the current number of step failures.
Return the tolerance parameters used in the solver convergence tests.
getType
()Return the type of the solver.
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.
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.
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.
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.
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.
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
Application context.
Absolute residual tolerance.
DM
.Convergence history.
Boolean indicating if the solver has converged.
Boolean indicating if the solver has failed.
Boolean indicating if the solver has not converged yet.
Number of iterations.
Linear solver.
Maximum number of function evaluations.
Maximum number of iterations.
Function norm.
Nonlinear preconditioner.
Converged reason.
Relative residual tolerance.
Solution update tolerance.
Boolean indicating if the solver uses coloring finite-differencing.
Boolean indicating if the solver uses matrix-free finite-differencing.
Right-hand side vector.
Solution vector.
Update vector.
Methods Documentation
- appendOptionsPrefix(prefix)
Append to the prefix used for searching for options in the database.
Logically collective.
See also
Working with PETSc options (TODO),
setOptionsPrefix
,petsc.SNESAppendOptionsPrefix
- callConvergenceTest(its, xnorm, ynorm, fnorm)
Compute the convergence test.
Collective.
- Parameters:
- Return type:
See also
- computeFunction(x, f)
Compute the function.
Collective.
See also
setFunction
,petsc.SNESComputeFunction
- computeJacobian(x, J, P=None)
Compute the Jacobian.
Collective.
- Parameters:
- Return type:
See also
setJacobian
,petsc.SNESComputeJacobian
- computeNGS(x, b=None)
Compute a nonlinear Gauss-Seidel step.
Collective.
- Parameters:
- Return type:
- computeObjective(x)
Compute the value of the objective function.
Collective.
See also
setObjective
,petsc.SNESComputeObjective
- converged(its, xnorm, ynorm, fnorm)
Compute the convergence test and update the solver converged reason.
Collective.
- Parameters:
- Return type:
See also
setConvergenceTest
,getConvergenceTest
,petsc.SNESConverged
- create(comm=None)
Create a SNES solver.
Collective.
- Parameters:
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.- Return type:
See also
Sys.getDefaultComm
,petsc.SNESCreate
- 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:
- destroy()
Destroy the solver.
Collective.
See also
petsc.SNESDestroy
Source code at petsc4py/PETSc/SNES.pyx:129
- Return type:
- getApplicationContext()
Return the application context.
Source code at petsc4py/PETSc/SNES.pyx:261
- Return type:
- getCompositeNumber()
Return the number of solvers in the composite.
Not collective.
See also
getCompositeSNES
,petsc.SNESCompositeGetNumber
,petsc.SNESCOMPOSITE
Source code at petsc4py/PETSc/SNES.pyx:2084
- Return type:
- getCompositeSNES(n)
Return the n-th solver in the composite.
Not collective.
See also
getCompositeNumber
,petsc.SNESCompositeGetSNES
,petsc.SNESCOMPOSITE
- getConvergedReason()
Return the termination flag.
Not collective.
See also
setConvergedReason
,petsc.SNESGetConvergedReason
Source code at petsc4py/PETSc/SNES.pyx:1570
- Return type:
- getConvergenceHistory()
Return the convergence history.
Not collective.
See also
petsc.SNESGetConvergenceHistory
- getConvergenceTest()
Return the callback to used as convergence test.
Not collective.
See also
Source code at petsc4py/PETSc/SNES.pyx:1147
- Return type:
- getDM()
Return the
DM
associated with the solver.Not collective.
See also
setDM
,petsc.SNESGetDM
Source code at petsc4py/PETSc/SNES.pyx:276
- Return type:
- getErrorIfNotConverged()
Return the flag indicating error on divergence.
Not collective.
See also
setErrorIfNotConverged
,petsc.SNESGetErrorIfNotConverged
Source code at petsc4py/PETSc/SNES.pyx:1597
- Return type:
- getFASCoarseSolve()
Return the
SNES
used at the coarsest level of the FAS hierarchy.Not collective.
See also
getFASSmoother
,petsc.SNESFASGetCoarseSolve
,petsc.SNESFAS
Source code at petsc4py/PETSc/SNES.pyx:475
- Return type:
- getFASCycleSNES(level)
Return the
SNES
corresponding to a particular level of the FAS hierarchy.Not collective.
See also
setFASLevels
,getFASCoarseSolve
,getFASSmoother
,petsc.SNESFASGetCycleSNES
,petsc.SNESFAS
- getFASInjection(level)
Return the
Mat
used to apply the injection from level-1 to level.Not collective.
See also
setFASInjection
,petsc.SNESFASGetInjection
,petsc.SNESFAS
- getFASInterpolation(level)
Return the
Mat
used to apply the interpolation from level-1 to level.Not collective.
See also
setFASInterpolation
,petsc.SNESFASGetInterpolation
,petsc.SNESFAS
- getFASLevels()
Return the number of levels used.
Not collective.
See also
setFASLevels
,petsc.SNESFASGetLevels
,petsc.SNESFAS
Source code at petsc4py/PETSc/SNES.pyx:444
- Return type:
- getFASRestriction(level)
Return the
Mat
used to apply the restriction from level-1 to level.Not collective.
See also
setFASRestriction
,petsc.SNESFASGetRestriction
,petsc.SNESFAS
- getFASSmoother(level)
Return the smoother used at a given level of the FAS hierarchy.
Not collective.
See also
setFASLevels
,getFASCoarseSolve
,getFASSmootherDown
,getFASSmootherUp
,petsc.SNESFASGetSmoother
,petsc.SNESFAS
- getFASSmootherDown(level)
Return the downsmoother used at a given level of the FAS hierarchy.
Not collective.
See also
setFASLevels
,getFASCoarseSolve
,getFASSmoother
,getFASSmootherUp
,petsc.SNESFASGetSmootherDown
,petsc.SNESFAS
- getFASSmootherUp(level)
Return the upsmoother used at a given level of the FAS hierarchy.
Not collective.
See also
setFASLevels
,getFASCoarseSolve
,getFASSmoother
,getFASSmootherDown
,petsc.SNESFASGetSmootherUp
,petsc.SNESFAS
- getFunction()
Return the callback to compute the nonlinear function.
Not collective.
See also
setFunction
,petsc.SNESGetFunction
Source code at petsc4py/PETSc/SNES.pyx:723
- Return type:
- getFunctionEvaluations()
Return the current number of function evaluations.
Not collective.
See also
setMaxFunctionEvaluations
,petsc.SNESGetNumberFunctionEvals
Source code at petsc4py/PETSc/SNES.pyx:1391
- Return type:
- getFunctionNorm()
Return the function norm.
Not collective.
See also
setFunctionNorm
,petsc.SNESGetFunctionNorm
Source code at petsc4py/PETSc/SNES.pyx:1668
- Return type:
- getInitialGuess()
Return the callback to compute the initial guess.
Not collective.
See also
Source code at petsc4py/PETSc/SNES.pyx:677
- Return type:
- getIterationNumber()
Return the current iteration number.
Not collective.
See also
setIterationNumber
,petsc.SNESGetIterationNumber
Source code at petsc4py/PETSc/SNES.pyx:1626
- Return type:
- 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:
See also
setJacobian
,petsc.SNESGetJacobian
- getKSP()
Return the linear solver used by the nonlinear solver.
Not collective.
See also
setKSP
,petsc.SNESGetKSP
Source code at petsc4py/PETSc/SNES.pyx:1767
- Return type:
- getKSPFailures()
Return the current number of linear solve failures.
Not collective.
See also
getMaxKSPFailures
,petsc.SNESGetLinearSolveFailures
Source code at petsc4py/PETSc/SNES.pyx:1473
- Return type:
- getLinearSolveIterations()
Return the total number of linear iterations.
Not collective.
See also
petsc.SNESGetLinearSolveIterations
Source code at petsc4py/PETSc/SNES.pyx:1682
- Return type:
- getMaxFunctionEvaluations()
Return the maximum allowed number of function evaluations.
Not collective.
See also
setMaxFunctionEvaluations
,petsc.SNESSetTolerances
Source code at petsc4py/PETSc/SNES.pyx:1375
- Return type:
- getMaxKSPFailures()
Return the maximum allowed number of linear solve failures.
Not collective.
See also
setMaxKSPFailures
,petsc.SNESGetMaxLinearSolveFailures
Source code at petsc4py/PETSc/SNES.pyx:1459
- Return type:
- getMaxStepFailures()
Return the maximum allowed number of step failures.
Not collective.
See also
setMaxStepFailures
,petsc.SNESGetMaxNonlinearStepFailures
Source code at petsc4py/PETSc/SNES.pyx:1418
- Return type:
- getMonitor()
Return the callback used to monitor solver convergence.
Not collective.
See also
- getNASMNumber()
Return the number of solvers in NASM.
Not collective.
See also
getNASMSNES
,petsc.SNESNASMGetNumber
,petsc.SNESNASM
Source code at petsc4py/PETSc/SNES.pyx:2116
- Return type:
- getNASMSNES(n)
Return the n-th solver in NASM.
Not collective.
See also
getNASMNumber
,petsc.SNESNASMGetSNES
,petsc.SNESNASM
- getNGS()
Return the nonlinear Gauss-Seidel callback tuple.
Not collective.
See also
Source code at petsc4py/PETSc/SNES.pyx:991
- Return type:
- getNPC()
Return the nonlinear preconditioner associated with the solver.
Not collective.
See also
setNPC
,hasNPC
,setNPCSide
,getNPCSide
,petsc.SNESGetNPC
Source code at petsc4py/PETSc/SNES.pyx:543
- Return type:
- getNPCSide()
Return the nonlinear preconditioning side.
Not collective.
See also
setNPC
,getNPC
,hasNPC
,setNPCSide
,petsc.SNESGetNPCSide
Source code at petsc4py/PETSc/SNES.pyx:596
- Return type:
- getNormSchedule()
Return the norm schedule.
Not collective.
See also
setNormSchedule
,petsc.SNESGetNormSchedule
Source code at petsc4py/PETSc/SNES.pyx:1098
- Return type:
- getObjective()
Return the objective callback tuple.
Not collective.
See also
Source code at petsc4py/PETSc/SNES.pyx:889
- Return type:
- getOptionsPrefix()
Return the prefix used for searching for options in the database.
Not collective.
See also
Working with PETSc options (TODO),
setOptionsPrefix
,petsc.SNESGetOptionsPrefix
Source code at petsc4py/PETSc/SNES.pyx:210
- Return type:
- getParamsEW()
Get the parameters of the Eisenstat and Walker trick.
Not collective.
See also
setUseEW
,setParamsEW
,petsc.SNESKSPGetParametersEW
- getPythonContext()
Return the instance of the class implementing the required Python methods.
Not collective.
Source code at petsc4py/PETSc/SNES.pyx:2020
- Return type:
- getPythonType()
Return the fully qualified Python name of the class used by the solver.
Not collective.
See also
PETSc Python nonlinear solver type (TODO),
setPythonContext
,setPythonType
,petsc.SNESPythonGetType
Source code at petsc4py/PETSc/SNES.pyx:2050
- Return type:
- getRhs()
Return the vector holding the right-hand side.
Not collective.
See also
petsc.SNESGetRhs
Source code at petsc4py/PETSc/SNES.pyx:1696
- Return type:
- getSolution()
Return the vector holding the solution.
Not collective.
See also
setSolution
,petsc.SNESGetSolution
Source code at petsc4py/PETSc/SNES.pyx:1711
- Return type:
- getSolutionUpdate()
Return the vector holding the solution update.
Not collective.
See also
petsc.SNESGetSolutionUpdate
Source code at petsc4py/PETSc/SNES.pyx:1738
- Return type:
- getStepFailures()
Return the current number of step failures.
Not collective.
See also
getMaxStepFailures
,petsc.SNESGetNonlinearStepFailures
Source code at petsc4py/PETSc/SNES.pyx:1432
- Return type:
- getTolerances()
Return the tolerance parameters used in the solver convergence tests.
Collective.
- Returns:
- Return type:
See also
setTolerances
,petsc.SNESGetTolerances
- getType()
Return the type of the solver.
Not collective.
See also
setType
,petsc.SNESGetType
Source code at petsc4py/PETSc/SNES.pyx:182
- Return type:
- getUpdate()
Return the callback to compute the update at the beginning of each step.
Not collective.
See also
Source code at petsc4py/PETSc/SNES.pyx:783
- Return type:
- getUseEW()
Return the flag indicating if the solver uses the Eisenstat-Walker trick.
Not Collective.
See also
setUseEW
,setParamsEW
,petsc.SNESKSPGetUseEW
Source code at petsc4py/PETSc/SNES.pyx:1805
- Return type:
- getUseFD()
Return
true
if the solver uses color finite-differencing for the Jacobian.Not collective.
See also
Source code at petsc4py/PETSc/SNES.pyx:1938
- Return type:
- getUseMF()
Return the flag indicating if the solver uses matrix-free finite-differencing.
Not collective.
See also
Source code at petsc4py/PETSc/SNES.pyx:1911
- Return type:
- getVIInactiveSet()
Return the index set for the inactive set.
Not collective.
See also
petsc.SNESVIGetInactiveSet
Source code at petsc4py/PETSc/SNES.pyx:1966
- Return type:
- hasNPC()
Return a boolean indicating whether the solver has a nonlinear preconditioner.
Not collective.
See also
setNPC
,getNPC
,setNPCSide
,getNPCSide
,petsc.SNESHasNPC
Source code at petsc4py/PETSc/SNES.pyx:558
- Return type:
- logConvergenceHistory(norm, linear_its=0)
Log residual norm and linear iterations.
- monitor(its, rnorm)
Monitor the solver.
Collective.
- Parameters:
its – Current number of iterations.
rnorm – Current value of the residual norm.
- Return type:
See also
setMonitor
,petsc.SNESMonitor
- monitorCancel()
Cancel all the monitors of the solver.
Logically collective.
See also
setMonitor
,petsc.SNESMonitorCancel
Source code at petsc4py/PETSc/SNES.pyx:1322
- Return type:
- reset()
Reset the solver.
Collective.
See also
petsc.SNESReset
Source code at petsc4py/PETSc/SNES.pyx:1522
- Return type:
- setApplicationContext(appctx)
Set the application context.
- setConvergedReason(reason)
Set the termination flag.
Collective.
See also
getConvergedReason
,petsc.SNESSetConvergedReason
Source code at petsc4py/PETSc/SNES.pyx:1557
- Parameters:
reason (ConvergedReason) –
- Return type:
- setConvergenceHistory(length=None, reset=False)
Set the convergence history.
Logically collective.
See also
petsc.SNESSetConvergenceHistory
Source code at petsc4py/PETSc/SNES.pyx:1217
- Return type:
- setConvergenceTest(converged, args=None, kargs=None)
Set the callback to use as convergence test.
Logically collective.
- Parameters:
- Return type:
See also
getConvergenceTest
,callConvergenceTest
,petsc.SNESSetConvergenceTest
- setErrorIfNotConverged(flag)
Immediately generate an error if the solver has not converged.
Collective.
See also
getErrorIfNotConverged
,petsc.SNESSetErrorIfNotConverged
- setFASInjection(level, mat)
Set the
Mat
to be used to apply the injection from level-1 to level.Collective.
See also
getFASInjection
,setFASInterpolation
,setFASRestriction
,petsc.SNESFASSetInjection
,petsc.SNESFAS
- setFASInterpolation(level, mat)
Set the
Mat
to be used to apply the interpolation from level-1 to level.Collective.
See also
getFASInterpolation
,setFASRestriction
,setFASInjection
,petsc.SNESFASSetInterpolation
,petsc.SNESFAS
- setFASLevels(levels, comms=None)
Set the number of levels to use with FAS.
Collective.
- Parameters:
levels (int) – The number of levels
comms (Sequence[Comm]) – An optional sequence of communicators of length
Logging Levels
, orNone
for the default communicatorSys.getDefaultComm
.
- Return type:
See also
getFASLevels
,petsc.SNESFASSetLevels
,petsc.SNESFAS
- setFASRScale(level, vec)
Set the scaling factor of the restriction operator from level to level-1.
Collective.
See also
setFASRestriction
,petsc.SNESFASSetRScale
,petsc.SNESFAS
- setFASRestriction(level, mat)
Set the
Mat
to be used to apply the restriction from level-1 to level.Collective.
See also
setFASRScale
,getFASRestriction
,setFASInterpolation
,setFASInjection
,petsc.SNESFASSetRestriction
,petsc.SNESFAS
- setForceIteration(force)
Force solve to take at least one iteration.
Collective.
See also
petsc.SNESSetForceIteration
- setFromOptions()
Configure the solver from the options database.
Collective.
See also
Working with PETSc options (TODO),
petsc.SNESSetFromOptions
Source code at petsc4py/PETSc/SNES.pyx:238
- Return type:
- setFunction(function, f=None, args=None, kargs=None)
Set the callback to compute the nonlinear function.
Logically collective.
- Parameters:
- Return type:
See also
getFunction
,petsc.SNESSetFunction
- setFunctionNorm(norm)
Set the function norm value.
Collective.
This is only of use to implementers of custom SNES types.
See also
getFunctionNorm
,petsc.SNESSetFunctionNorm
- setInitialGuess(initialguess, args=None, kargs=None)
Set the callback to compute the initial guess.
Logically collective.
- Parameters:
- Return type:
See also
getInitialGuess
,petsc.SNESSetComputeInitialGuess
- setIterationNumber(its)
Set the current iteration number.
Collective.
This is only of use to implementers of custom SNES types.
See also
getIterationNumber
,petsc.SNESSetIterationNumber
- setJacobian(jacobian, J=None, P=None, args=None, kargs=None)
Set the callback to compute the Jacobian.
Logically collective.
- Parameters:
- Return type:
See also
getJacobian
,petsc.SNESSetJacobian
- setKSP(ksp)
Set the linear solver that will be used by the nonlinear solver.
Logically collective.
See also
getKSP
,petsc.SNESSetKSP
- setLineSearchPreCheck(precheck, args=None, kargs=None)
Set the callback that will be called before applying the linesearch.
Logically collective.
- Parameters:
- Return type:
See also
petsc.SNESLineSearchSetPreCheck
- setMaxFunctionEvaluations(max_funcs)
Set the maximum allowed number of function evaluations.
Collective.
See also
getMaxFunctionEvaluations
,petsc.SNESSetTolerances
- setMaxKSPFailures(max_fails)
Set the maximum allowed number of linear solve failures.
Collective.
See also
getMaxKSPFailures
,petsc.SNESSetMaxLinearSolveFailures
- setMaxStepFailures(max_fails)
Set the maximum allowed number of step failures.
Collective.
See also
getMaxStepFailures
,petsc.SNESSetMaxNonlinearStepFailures
- setMonitor(monitor, args=None, kargs=None)
Set the callback used to monitor solver convergence.
Logically collective.
- Parameters:
- Return type:
See also
getMonitor
,petsc.SNESMonitorSet
- setNGS(ngs, args=None, kargs=None)
Set the callback to compute nonlinear Gauss-Seidel.
Logically collective.
- Parameters:
- Return type:
See also
getNGS
,computeNGS
,petsc.SNESSetNGS
- setNPC(snes)
Set the nonlinear preconditioner.
Logically collective.
See also
getNPC
,hasNPC
,setNPCSide
,getNPCSide
,petsc.SNESSetNPC
- setNPCSide(side)
Set the nonlinear preconditioning side.
Collective.
See also
setNPC
,getNPC
,hasNPC
,getNPCSide
,petsc.SNESSetNPCSide
- setNormSchedule(normsched)
Set the norm schedule.
Collective.
See also
getNormSchedule
,petsc.SNESSetNormSchedule
Source code at petsc4py/PETSc/SNES.pyx:1086
- Parameters:
normsched (NormSchedule) –
- Return type:
- setObjective(objective, args=None, kargs=None)
Set the callback to compute the objective function.
Logically collective.
- Parameters:
- Return type:
See also
getObjective
,petsc.SNESSetObjective
- setOptionsPrefix(prefix)
Set the prefix used for searching for options in the database.
Logically collective.
See also
Working with PETSc options (TODO),
petsc.SNESSetOptionsPrefix
- 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:
See also
setUseEW
,getParamsEW
,petsc.SNESKSPSetParametersEW
- setPatchCellNumbering(sec)
Set cell patch numbering.
- setPatchComputeFunction(function, args=None, kargs=None)
Set patch compute function.
Source code at petsc4py/PETSc/SNES.pyx:2185
- Return type:
- setPatchComputeOperator(operator, args=None, kargs=None)
Set patch compute operator.
Source code at petsc4py/PETSc/SNES.pyx:2177
- Return type:
- setPatchConstructType(typ, operator=None, args=None, kargs=None)
Set patch construct type.
Source code at petsc4py/PETSc/SNES.pyx:2193
- Return type:
- setPatchDiscretisationInfo(dms, bs, cellNodeMaps, subspaceOffsets, ghostBcNodes, globalBcNodes)
Set patch discretisation information.
Source code at petsc4py/PETSc/SNES.pyx:2136
- Return type:
- setPythonContext(context)
Set the instance of the class implementing the required Python methods.
Not collective.
- setPythonType(py_type)
Set the fully qualified Python name of the class to be used.
Collective.
See also
PETSc Python nonlinear solver type (TODO),
setPythonContext
,getPythonType
,petsc.SNESPythonSetType
- setResetCounters(reset=True)
Set the flag to reset the counters.
Collective.
See also
petsc.SNESSetCountersReset
- setSolution(vec)
Set the vector used to store the solution.
Collective.
See also
getSolution
,petsc.SNESSetSolution
- setTolerances(rtol=None, atol=None, stol=None, max_it=None)
Set the tolerance parameters used in the solver convergence tests.
Collective.
- Parameters:
- Return type:
See also
getTolerances
,petsc.SNESSetTolerances
- setType(snes_type)
Set the type of the solver.
Logically collective.
See also
getType
,petsc.SNESSetType
- setUp()
Set up the internal data structures for using the solver.
Collective.
See also
petsc.SNESSetUp
Source code at petsc4py/PETSc/SNES.pyx:1496
- Return type:
- setUpMatrices()
Ensures that matrices are available for Newton-like methods.
Collective.
This is only of use to implementers of custom SNES types.
See also
setUp
,petsc.SNESSetUpMatrices
Source code at petsc4py/PETSc/SNES.pyx:1508
- Return type:
- setUpdate(update, args=None, kargs=None)
Set the callback to compute update at the beginning of each step.
Logically collective.
- Parameters:
- Return type:
See also
getUpdate
,petsc.SNESSetUpdate
- 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:
See also
getUseEW
,setParamsEW
,petsc.SNESKSPSetUseEW
- setUseFD(flag=True)
Set the boolean flag to use coloring finite-differencing for Jacobian assembly.
Logically collective.
See also
Source code at petsc4py/PETSc/SNES.pyx:1925
- Return type:
- setUseMF(flag=True)
Set the boolean flag indicating to use matrix-free finite-differencing.
Logically collective.
See also
Source code at petsc4py/PETSc/SNES.pyx:1898
- Return type:
- setVariableBounds(xl, xu)
Set the vector for the variable bounds.
Collective.
See also
petsc.SNESVISetVariableBounds
- solve(b=None, x=None)
Solve the nonlinear equations.
Collective.
- Parameters:
- Return type:
See also
setSolution
,getSolution
,petsc.SNESSolve
- view(viewer=None)
View the solver.
Collective.
See also
petsc.SNESView
Attributes Documentation
- appctx
Application context.
- atol
Absolute residual tolerance.
- 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_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.