File: KSPPIPEGCR.html

package info (click to toggle)
petsc 3.14.5%2Bdfsg1-4
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 266,472 kB
  • sloc: ansic: 680,898; python: 33,303; cpp: 16,324; makefile: 14,022; f90: 13,731; javascript: 10,713; fortran: 9,581; sh: 1,373; xml: 619; objc: 445; csh: 192; pascal: 148; java: 13
file content (56 lines) | stat: -rw-r--r-- 4,368 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
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPPIPEGCR.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>KSPPIPEGCR</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/KSP/KSPPIPEGCR.html "><small>Report Typos and Errors</small></a></div>
<A NAME="KSPPIPEGCR"><H1>KSPPIPEGCR</H1></A>
Implements a Pipelined Generalized Conjugate Residual method. 
<H3><FONT COLOR="#CC3333">Options Database Keys</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_pipegcr_mmax &lt;N&gt;  </B></TD><TD>- the max number of Krylov directions to orthogonalize against
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_pipegcr_unroll_w </B></TD><TD>- unroll w at the storage cost of a maximum of (mmax+1) extra vectors with the benefit of better pipelining (default: <A HREF="../Sys/PETSC_TRUE.html#PETSC_TRUE">PETSC_TRUE</A>)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_pipegcr_nprealloc &lt;N&gt; </B></TD><TD>- the number of vectors to preallocated for storing Krylov directions. Once exhausted new directions are allocated blockwise (default: 5)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_pipegcr_truncation_type &lt;standard,notay&gt; </B></TD><TD>- which previous search directions to orthogonalize against
</TD></TR></TABLE>
<P>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
The PIPEGCR Krylov method supports non-symmetric matrices and permits the use of a preconditioner
which may vary from one iteration to the next. Users can can define a method to vary the
preconditioner between iterates via <A HREF="../KSP/KSPPIPEGCRSetModifyPC.html#KSPPIPEGCRSetModifyPC">KSPPIPEGCRSetModifyPC</A>().
Restarts are solves with x0 not equal to zero. When a restart occurs, the initial starting
solution is given by the current estimate for x which was obtained by the last restart
iterations of the PIPEGCR algorithm.
The method implemented requires at most the storage of 4 x mmax + 5 vectors, roughly twice as much as GCR.
<P>
Only supports left preconditioning.
<P>
The natural "norm" for this method is (u,Au), where u is the preconditioned residual. This norm is available at no additional computational cost, as with standard CG. Choosing preconditioned or unpreconditioned norm types involves a blocking reduction which prevents any benefit from pipelining.
<P>
<H3><FONT COLOR="#CC3333">Reference</FONT></H3>
P. Sanan, S.M. Schnepp, and D.A. May,
"Pipelined, Flexible Krylov Subspace Methods,"
SIAM Journal on Scientific Computing 2016 38:5, C441-C470,
DOI: 10.1137/15M1049130
<P>

<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
  <A HREF="../KSP/KSPCreate.html#KSPCreate">KSPCreate</A>(), <A HREF="../KSP/KSPSetType.html#KSPSetType">KSPSetType</A>(), <A HREF="../KSP/KSPType.html#KSPType">KSPType</A> (for list of available types), <A HREF="../KSP/KSP.html#KSP">KSP</A>,
<BR><A HREF="../KSP/KSPPIPEFGMRES.html#KSPPIPEFGMRES">KSPPIPEFGMRES</A>, <A HREF="../KSP/KSPPIPECG.html#KSPPIPECG">KSPPIPECG</A>, <A HREF="../KSP/KSPPIPECR.html#KSPPIPECR">KSPPIPECR</A>, <A HREF="../KSP/KSPPIPEFCG.html#KSPPIPEFCG">KSPPIPEFCG</A>,<A HREF="../KSP/KSPPIPEGCRSetTruncationType.html#KSPPIPEGCRSetTruncationType">KSPPIPEGCRSetTruncationType</A>(),<A HREF="../KSP/KSPPIPEGCRSetNprealloc.html#KSPPIPEGCRSetNprealloc">KSPPIPEGCRSetNprealloc</A>(),<A HREF="../KSP/KSPPIPEGCRSetUnrollW.html#KSPPIPEGCRSetUnrollW">KSPPIPEGCRSetUnrollW</A>(),<A HREF="../KSP/KSPPIPEGCRSetMmax.html#KSPPIPEGCRSetMmax">KSPPIPEGCRSetMmax</A>()
<P>
<P>
<P><B></B><H3><FONT COLOR="#CC3333">Level</FONT></H3>intermediate<BR>
<H3><FONT COLOR="#CC3333">Location</FONT></H3>
</B><A HREF="../../../src/ksp/ksp/impls/gcr/pipegcr/pipegcr.c.html#KSPPIPEGCR">src/ksp/ksp/impls/gcr/pipegcr/pipegcr.c</A>
<BR><A HREF="./index.html">Index of all KSP 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>