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
|
<!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/KSPHPDDM.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>KSPHPDDM</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/KSPHPDDM.html "><small>Report Typos and Errors</small></a></div>
<A NAME="KSPHPDDM"><H1>KSPHPDDM</H1></A>
Interface with the HPDDM library. This <A HREF="../KSP/KSP.html#KSP">KSP</A> may be used to further select methods that are currently not implemented natively in PETSc, e.g., GCRODR [2006], a recycled Krylov method which is similar to <A HREF="../KSP/KSPLGMRES.html#KSPLGMRES">KSPLGMRES</A>, see [2016] for a comparison. ex75.c shows how to reproduce the results from the aforementioned paper [2006]. A chronological bibliography of relevant publications linked with <A HREF="../KSP/KSP.html#KSP">KSP</A> available in HPDDM through <A HREF="../KSP/KSPHPDDM.html#KSPHPDDM">KSPHPDDM</A>, and not available directly in PETSc, may be found below.
<P>
<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_gmres_restart <restart, default=30> </B></TD><TD>- see <A HREF="../KSP/KSPGMRES.html#KSPGMRES">KSPGMRES</A>
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_hpddm_type <type, default=gmres> </B></TD><TD>- any of gmres, bgmres, cg, bcg, gcrodr, bgcrodr, bfbcg, or preonly, see <A HREF="../KSP/KSPHPDDMType.html#KSPHPDDMType">KSPHPDDMType</A>
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_hpddm_deflation_tol <eps, default=-1.0> </B></TD><TD>- tolerance when deflating right-hand sides inside block methods (no deflation by default, only relevant with block methods)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_hpddm_enlarge_krylov_subspace <p, default=1> </B></TD><TD>- split the initial right-hand side into multiple vectors (only relevant with nonblock methods)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_hpddm_orthogonalization <type, default=cgs> </B></TD><TD>- any of cgs or mgs, see <A HREF="../KSP/KSPGMRES.html#KSPGMRES">KSPGMRES</A>
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_hpddm_qr <type, default=cholqr> </B></TD><TD>- distributed QR factorizations with any of cholqr, cgs, or mgs (only relevant with block methods)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_hpddm_variant <type, default=left> </B></TD><TD>- any of left, right, or flexible (this option is superseded by <A HREF="../KSP/KSPSetPCSide.html#KSPSetPCSide">KSPSetPCSide</A>())
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_hpddm_recycle <n, default=0> </B></TD><TD>- number of harmonic Ritz vectors to compute (only relevant with GCRODR or BGCRODR)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_hpddm_recycle_target <type, default=SM> </B></TD><TD>- criterion to select harmonic Ritz vectors using either SM, LM, SR, LR, SI, or LI (only relevant with GCRODR or BGCRODR). For BGCRODR, if PETSc is compiled with SLEPc, this option is not relevant, since SLEPc is used instead. Options are set with the prefix -ksp_hpddm_recycle_eps_
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_hpddm_recycle_strategy <type, default=A> </B></TD><TD>- generalized eigenvalue problem A or B to solve for recycling (only relevant with flexible GCRODR or BGCRODR)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_hpddm_recycle_symmetric <true, default=false> </B></TD><TD>- symmetric generalized eigenproblems in BGCRODR, useful to switch to distributed solvers like EPSELEMENTAL (only relevant when PETSc is compiled with SLEPc)
</TD></TR></TABLE>
<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>1980 </B></TD><TD>- The Block Conjugate Gradient Algorithm and Related Methods. O'Leary. Linear Algebra and its Applications.
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>2006 </B></TD><TD>- Recycling Krylov Subspaces for Sequences of Linear Systems. Parks, de Sturler, Mackey, Johnson, and Maiti. SIAM Journal on Scientific Computing
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>2013 </B></TD><TD>- A Modified Block Flexible GMRES Method with Deflation at Each Iteration for the Solution of Non-Hermitian Linear Systems with Multiple Right-Hand Sides. Calandra, Gratton, Lago, Vasseur, and Carvalho. SIAM Journal on Scientific Computing.
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>2016 </B></TD><TD>- Block Iterative Methods and Recycling for Improved Scalability of Linear Solvers. Jolivet and Tournier. SC16.
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>2017 </B></TD><TD>- A breakdown-free block conjugate gradient method. Ji and Li. BIT Numerical Mathematics.
</TD></TR></TABLE>
<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>, <A HREF="../KSP/KSPGMRES.html#KSPGMRES">KSPGMRES</A>, <A HREF="../KSP/KSPCG.html#KSPCG">KSPCG</A>, <A HREF="../KSP/KSPLGMRES.html#KSPLGMRES">KSPLGMRES</A>, <A HREF="../KSP/KSPDGMRES.html#KSPDGMRES">KSPDGMRES</A>
<BR><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/hpddm/hpddm.cxx.html#KSPHPDDM">src/ksp/ksp/impls/hpddm/hpddm.cxx</A>
<P><H3><FONT COLOR="#CC3333">Examples</FONT></H3>
<A HREF="../../../src/ksp/ksp/tutorials/ex75.c.html">src/ksp/ksp/tutorials/ex75.c.html</A><BR>
<A HREF="../../../src/ksp/ksp/tutorials/ex76.c.html">src/ksp/ksp/tutorials/ex76.c.html</A><BR>
<A HREF="../../../src/ksp/ksp/tutorials/ex77.c.html">src/ksp/ksp/tutorials/ex77.c.html</A><BR>
<A HREF="../../../src/ksp/ksp/tutorials/ex78.c.html">src/ksp/ksp/tutorials/ex78.c.html</A><BR>
<A HREF="../../../src/ksp/ksp/tutorials/ex75f.F90.html">src/ksp/ksp/tutorials/ex75f.F90.html</A><BR>
<A HREF="../../../src/ksp/ksp/tutorials/ex77f.F90.html">src/ksp/ksp/tutorials/ex77f.F90.html</A><BR>
<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>
|