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 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCFieldSplitSetSchurPre.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>PCFieldSplitSetSchurPre</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.7.5 2017-01-01</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.7.5 v3.7.5 docs/manualpages/PC/PCFieldSplitSetSchurPre.html "><small>Report Typos and Errors</small></a></div>
<A NAME="PCFieldSplitSetSchurPre"><H1>PCFieldSplitSetSchurPre</H1></A>
Indicates if the Schur complement is preconditioned by a preconditioner constructed by the A11 matrix. Otherwise no preconditioner is used.
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscpc.h"
PetscErrorCode PCFieldSplitSetSchurPre(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
</PRE>
Collective on <A HREF="../PC/PC.html#PC">PC</A>
<P>
<H3><FONT COLOR="#CC3333">Input Parameters</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>pc </B></TD><TD>- the preconditioner context
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ptype </B></TD><TD>- which matrix to use for preconditioning the Schur complement: PC_FIELDSPLIT_SCHUR_PRE_A11 (default), PC_FIELDSPLIT_SCHUR_PRE_SELF, PC_FIELDSPLIT_SCHUR_PRE_USER
PC_FIELDSPLIT_SCHUR_PRE_SELFP, and PC_FIELDSPLIT_SCHUR_PRE_FULL
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>userpre </B></TD><TD>- matrix to use for preconditioning, or NULL
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Options Database</FONT></H3>
<DT><B>-pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> </B> -default is a11. See notes for meaning of various arguments
<br>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
<pre>
If ptype is
</pre>
<pre>
a11 then the preconditioner for the Schur complement is generated from the block diagonal part of the preconditioner
</pre>
<pre>
matrix associated with the Schur complement (i.e. A11), not he Schur complement matrix
</pre>
<pre>
self the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix:
</pre>
<pre>
The only preconditioner that currently works with this symbolic respresentation matrix object is the <A HREF="../PC/PCLSC.html#PCLSC">PCLSC</A>
</pre>
<pre>
preconditioner
</pre>
<pre>
user then the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
</pre>
<pre>
to this function).
</pre>
<pre>
selfp then the preconditioning for the Schur complement is generated from an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
</pre>
<pre>
This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
</pre>
<pre>
lumped before extracting the diagonal using the additional option -fieldsplit_1_mat_schur_complement_ainv_type lump
</pre>
<pre>
full then the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation computed internally by PFIELDSPLIT (this is expensive)
</pre>
<pre>
useful mostly as a test that the Schur complement approach can work for your problem
</pre>
<P>
When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense
with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
-fieldsplit_1_pc_type lsc which uses the least squares commutator to compute a preconditioner for the Schur complement.
<P>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../PC/PCFieldSplitGetSchurPre.html#PCFieldSplitGetSchurPre">PCFieldSplitGetSchurPre</A>(), <A HREF="../PC/PCFieldSplitGetSubKSP.html#PCFieldSplitGetSubKSP">PCFieldSplitGetSubKSP</A>(), <A HREF="../PC/PCFIELDSPLIT.html#PCFIELDSPLIT">PCFIELDSPLIT</A>, <A HREF="../PC/PCFieldSplitSetFields.html#PCFieldSplitSetFields">PCFieldSplitSetFields</A>(), <A HREF="../PC/PCFieldSplitSchurPreType.html#PCFieldSplitSchurPreType">PCFieldSplitSchurPreType</A>,
<BR><A HREF="../KSP/MatSchurComplementSetAinvType.html#MatSchurComplementSetAinvType">MatSchurComplementSetAinvType</A>(), <A HREF="../PC/PCLSC.html#PCLSC">PCLSC</A>
<P>
<P><B><P><B><FONT COLOR="#CC3333">Level:</FONT></B>intermediate
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/ksp/pc/impls/fieldsplit/fieldsplit.c.html#PCFieldSplitSetSchurPre">src/ksp/pc/impls/fieldsplit/fieldsplit.c</A>
<BR><A HREF="./index.html">Index of all PC 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>
<P><H3><FONT COLOR="#CC3333">Examples</FONT></H3>
<A HREF="../../../src/snes/examples/tutorials/ex70.c.html">src/snes/examples/tutorials/ex70.c.html</A><BR>
</BODY></HTML>
|