1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESASPIN.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>SNESASPIN</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.14.5 2021-03-03</b></div>
<div id="bugreport" align=right><a href="mailto:petsc-maint@mcs.anl.gov?subject=Typo or Error in Documentation &body=Please describe the typo or error in the documentation: petsc-3.14.5 v3.14.5 docs/manualpages/SNES/SNESASPIN.html "><small>Report Typos and Errors</small></a></div>
<A NAME="SNESASPIN"><H1>SNESASPIN</H1></A>
Helper <A HREF="../SNES/SNES.html#SNES">SNES</A> type for Additive-Schwarz Preconditioned Inexact Newton
<H3><FONT COLOR="#CC3333">Options Database</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-npc_snes_ </B></TD><TD>- options prefix of the nonlinear subdomain solver (must be of type NASM)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-npc_sub_snes_ </B></TD><TD>- options prefix of the subdomain nonlinear solves
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-npc_sub_ksp_ </B></TD><TD>- options prefix of the subdomain Krylov solver
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-npc_sub_pc_ </B></TD><TD>- options prefix of the subdomain preconditioner
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
This routine sets up an instance of NETWONLS with nonlinear left preconditioning. It differs from other
similar functionality in <A HREF="../SNES/SNES.html#SNES">SNES</A> as it creates a linear shell matrix that corresponds to the product
<P>
\sum_{i=0}^{N_b}J_b({X^b_{converged}})^{-1}J(X + \sum_{i=0}^{N_b}(X^b_{converged} - X^b))
<P>
which is the ASPIN preconditioned matrix. Similar solvers may be constructed by having matrix-free differencing of
nonlinear solves per linear iteration, but this is far more efficient when subdomain sparse-direct preconditioner
factorizations are reused on each application of J_b^{-1}.
<P>
The Krylov method used in this nonlinear solver is run with NO preconditioner, because the preconditioning is done
at the nonlinear level, but the Jacobian for the original function must be provided (or calculated via coloring and
finite differences automatically) in the Pmat location of <A HREF="../SNES/SNESSetJacobian.html#SNESSetJacobian">SNESSetJacobian</A>() because the action of the original Jacobian
is needed by the shell matrix used to apply the Jacobian of the nonlinear preconditioned problem (see above).
Note that since the Pmat is not used to construct a preconditioner it could be provided in a matrix-free form.
The code for this implementation is a bit confusing because the Amat of <A HREF="../SNES/SNESSetJacobian.html#SNESSetJacobian">SNESSetJacobian</A>() applies the Jacobian of the
nonlinearly preconditioned function Jacobian while the Pmat provides the Jacobian of the original user provided function.
Note that the original <A HREF="../SNES/SNES.html#SNES">SNES</A> and nonlinear preconditioner preconditioner (see <A HREF="../SNES/SNESGetNPC.html#SNESGetNPC">SNESGetNPC</A>()), in this case NASM, share
the same Jacobian matrices. <A HREF="../SNES/SNESNASM.html#SNESNASM">SNESNASM</A> computes the needed Jacobian in SNESNASMComputeFinalJacobian_Private().
<P>
<P>
<H3><FONT COLOR="#CC3333">References</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>1. </B></TD><TD>- X. C. Cai and D. E. Keyes, "Nonlinearly preconditioned inexact Newton algorithms", SIAM J. Sci. Comput., 24, 2002.
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>2. </B></TD><TD>- Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
SIAM Review, 57(4), 2015
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../SNES/SNESCreate.html#SNESCreate">SNESCreate</A>(), <A HREF="../SNES/SNES.html#SNES">SNES</A>, <A HREF="../SNES/SNESSetType.html#SNESSetType">SNESSetType</A>(), <A HREF="../SNES/SNESNEWTONLS.html#SNESNEWTONLS">SNESNEWTONLS</A>, <A HREF="../SNES/SNESNASM.html#SNESNASM">SNESNASM</A>, <A HREF="../SNES/SNESGetNPC.html#SNESGetNPC">SNESGetNPC</A>(), <A HREF="../SNES/SNESGetNPCSide.html#SNESGetNPCSide">SNESGetNPCSide</A>()
<BR>
<P>
<P><B></B><H3><FONT COLOR="#CC3333">Level</FONT></H3>intermediate<BR>
<H3><FONT COLOR="#CC3333">Location</FONT></H3>
</B><A HREF="../../../src/snes/impls/nasm/aspin.c.html#SNESASPIN">src/snes/impls/nasm/aspin.c</A>
<BR><A HREF="./index.html">Index of all SNES routines</A>
<BR><A HREF="../../index.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
</BODY></HTML>
|