File: PCLSC.html

package info (click to toggle)
petsc 3.4.2.dfsg1-8.1
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 129,104 kB
  • ctags: 516,422
  • sloc: ansic: 395,939; cpp: 47,201; python: 34,788; makefile: 17,193; fortran: 16,251; f90: 1,592; objc: 954; sh: 822; xml: 621; java: 381; lisp: 293; csh: 241
file content (84 lines) | stat: -rw-r--r-- 5,054 bytes parent folder | download
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
<!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/PCLSC.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>PCLSC</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
   <div id="version" align=right><b>petsc-3.4.2 2013-07-02</b></div>
<A NAME="PCLSC"><H1>PCLSC</H1></A>
Preconditioning for Schur complements, based on Least Squares Commutators 
<H3><FONT COLOR="#CC3333">Options Database Key</FONT></H3>
<DT><B>-pc_lsc_scale_diag </B> -Use the diagonal of A for scaling
<br>
<P>

<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
This preconditioner will normally be used with PCFieldSplit to precondition the Schur complement, but
it can be used for any Schur complement system.  Consider the Schur complement
<P>
<PRE>
   S = A11 - A10 inv(A00) A01
</PRE>

<P>
<A HREF="../PC/PCLSC.html#PCLSC">PCLSC</A> currently doesn't do anything with A11, so let's assume it is 0.  The idea is that a good approximation to
inv(S) is given by
<P>
<PRE>
   inv(A10 A01) A10 A00 A01 inv(A10 A01)
</PRE>

<P>
The product A10 A01 can be computed for you, but you can provide it (this is
usually more efficient anyway).  In the case of incompressible flow, A10 A10 is a Laplacian, call it L.  The current
interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix.
<P>
If you had called <A HREF="../KSP/KSPSetOperators.html#KSPSetOperators">KSPSetOperators</A>(ksp,S,Sp,flg), S should have type MATSCHURCOMPLEMENT and Sp can be any type you
like (<A HREF="../PC/PCLSC.html#PCLSC">PCLSC</A> doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp".
For example, you might have setup code like this
<P>
<PRE>
   <A HREF="../Sys/PetscObjectCompose.html#PetscObjectCompose">PetscObjectCompose</A>((<A HREF="../Sys/PetscObject.html#PetscObject">PetscObject</A>)Sp,"LSC_L",(<A HREF="../Sys/PetscObject.html#PetscObject">PetscObject</A>)L);
   <A HREF="../Sys/PetscObjectCompose.html#PetscObjectCompose">PetscObjectCompose</A>((<A HREF="../Sys/PetscObject.html#PetscObject">PetscObject</A>)Sp,"LSC_Lp",(<A HREF="../Sys/PetscObject.html#PetscObject">PetscObject</A>)Lp);
</PRE>

<P>
And then your Jacobian assembly would look like
<P>
<PRE>
   <A HREF="../Sys/PetscObjectQuery.html#PetscObjectQuery">PetscObjectQuery</A>((<A HREF="../Sys/PetscObject.html#PetscObject">PetscObject</A>)Sp,"LSC_L",(<A HREF="../Sys/PetscObject.html#PetscObject">PetscObject</A>*)&amp;L);
   <A HREF="../Sys/PetscObjectQuery.html#PetscObjectQuery">PetscObjectQuery</A>((<A HREF="../Sys/PetscObject.html#PetscObject">PetscObject</A>)Sp,"LSC_Lp",(<A HREF="../Sys/PetscObject.html#PetscObject">PetscObject</A>*)&amp;Lp);
   if (L) { assembly L }
   if (Lp) { assemble Lp }
</PRE>

<P>
With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L
<P>
<PRE>
   -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml
</PRE>

<P>
Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners.
<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>Elman, Howle, Shadid, Shuttleworth, and Tuminaro, Block preconditioners based on approximate commutators, 2006.</B></TD><TD>- -  Silvester, Elman, Kay, Wathen, Efficient preconditioning of the linearized Navier-Stokes equations for incompressible flow, 2001.
</TD></TR>
<P>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
  <A HREF="../PC/PCCreate.html#PCCreate">PCCreate</A>(), <A HREF="../PC/PCSetType.html#PCSetType">PCSetType</A>(), <A HREF="../PC/PCType.html#PCType">PCType</A> (for list of available types), <A HREF="../PC/PC.html#PC">PC</A>, Block_Preconditioners, <A HREF="../PC/PCFIELDSPLIT.html#PCFIELDSPLIT">PCFIELDSPLIT</A>,
<BR><A HREF="../PC/PCFieldSplitGetSubKSP.html#PCFieldSplitGetSubKSP">PCFieldSplitGetSubKSP</A>(), <A HREF="../PC/PCFieldSplitSetFields.html#PCFieldSplitSetFields">PCFieldSplitSetFields</A>(), <A HREF="../PC/PCFieldSplitSetType.html#PCFieldSplitSetType">PCFieldSplitSetType</A>(), <A HREF="../PC/PCFieldSplitSetIS.html#PCFieldSplitSetIS">PCFieldSplitSetIS</A>(), <A HREF="../PC/PCFieldSplitSchurPrecondition.html#PCFieldSplitSchurPrecondition">PCFieldSplitSchurPrecondition</A>(),
<A HREF="../KSP/MatCreateSchurComplement.html#MatCreateSchurComplement">MatCreateSchurComplement</A>()
<P><B><P><B><FONT COLOR="#CC3333">Level:</FONT></B>intermediate
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/ksp/pc/impls/lsc/lsc.c.html#PCLSC">src/ksp/pc/impls/lsc/lsc.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/ex55.c.html">src/snes/examples/tutorials/ex55.c.html</A><BR>
</BODY></HTML>