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
|
<!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/PCApplyRichardson.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>PCApplyRichardson</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.10.3 2018-12-18</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.10.3 v3.10.3 docs/manualpages/PC/PCApplyRichardson.html "><small>Report Typos and Errors</small></a></div>
<A NAME="PCApplyRichardson"><H1>PCApplyRichardson</H1></A>
Applies several steps of Richardson iteration with the particular preconditioner. This routine is usually used by the Krylov solvers and not the application code directly.
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscksp.h"
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../PC/PCApplyRichardson.html#PCApplyRichardson">PCApplyRichardson</A>(<A HREF="../PC/PC.html#PC">PC</A> pc,<A HREF="../Vec/Vec.html#Vec">Vec</A> b,<A HREF="../Vec/Vec.html#Vec">Vec</A> y,<A HREF="../Vec/Vec.html#Vec">Vec</A> w,<A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> rtol,<A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> abstol, <A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> dtol,<A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> its,<A HREF="../Sys/PetscBool.html#PetscBool">PetscBool</A> guesszero,<A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> *outits,<A HREF="../PC/PCRichardsonConvergedReason.html#PCRichardsonConvergedReason">PCRichardsonConvergedReason</A> *reason)
</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>b </B></TD><TD>- the right hand side
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>w </B></TD><TD>- one work vector
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>rtol </B></TD><TD>- relative decrease in residual norm convergence criteria
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>abstol </B></TD><TD>- absolute residual norm convergence criteria
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>dtol </B></TD><TD>- divergence residual norm increase criteria
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>its </B></TD><TD>- the number of iterations to apply.
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>guesszero </B></TD><TD>- if the input x contains nonzero initial guess
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Output Parameter</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>outits </B></TD><TD>- number of iterations actually used (for SOR this always equals its)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>reason </B></TD><TD>- the reason the apply terminated
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>y </B></TD><TD>- the solution (also contains initial guess if guesszero is <A HREF="../Sys/PETSC_FALSE.html#PETSC_FALSE">PETSC_FALSE</A>
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
Most preconditioners do not support this function. Use the command
<A HREF="../PC/PCApplyRichardsonExists.html#PCApplyRichardsonExists">PCApplyRichardsonExists</A>() to determine if one does.
<P>
Except for the multigrid <A HREF="../PC/PC.html#PC">PC</A> this routine ignores the convergence tolerances
and always runs for the number of iterations
<P>
<P>
<H3><FONT COLOR="#CC3333">Keywords</FONT></H3>
<A HREF="../PC/PC.html#PC">PC</A>, apply, Richardson
<BR>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../PC/PCApplyRichardsonExists.html#PCApplyRichardsonExists">PCApplyRichardsonExists</A>()
<BR><P><B></B><H3><FONT COLOR="#CC3333">Level</FONT></H3>developer<BR>
<H3><FONT COLOR="#CC3333">Location</FONT></H3>
</B><A HREF="../../../src/ksp/pc/interface/precon.c.html#PCApplyRichardson">src/ksp/pc/interface/precon.c</A>
<P><H3><FONT COLOR="CC3333">Implementations</FONT></H3><A HREF="../../../src/ksp/pc/impls/hypre/hypre.c.html#PCApplyRichardson_HYPRE_BoomerAMG">PCApplyRichardson_HYPRE_BoomerAMG in src/ksp/pc/impls/hypre/hypre.c</A><BR>
<A HREF="../../../src/ksp/pc/impls/hypre/hypre.c.html#PCApplyRichardson_PFMG">PCApplyRichardson_PFMG in src/ksp/pc/impls/hypre/hypre.c</A><BR>
<A HREF="../../../src/ksp/pc/impls/hypre/hypre.c.html#PCApplyRichardson_SysPFMG">PCApplyRichardson_SysPFMG in src/ksp/pc/impls/hypre/hypre.c</A><BR>
<A HREF="../../../src/ksp/pc/impls/mg/mg.c.html#PCApplyRichardson_MG">PCApplyRichardson_MG in src/ksp/pc/impls/mg/mg.c</A><BR>
<A HREF="../../../src/ksp/pc/impls/shell/shellpc.c.html#PCApplyRichardson_Shell">PCApplyRichardson_Shell in src/ksp/pc/impls/shell/shellpc.c</A><BR>
<A HREF="../../../src/ksp/pc/impls/sor/sor.c.html#PCApplyRichardson_SOR">PCApplyRichardson_SOR in src/ksp/pc/impls/sor/sor.c</A><BR>
<A HREF="../../../src/ksp/pc/impls/telescope/telescope.c.html#PCApplyRichardson_Telescope">PCApplyRichardson_Telescope in src/ksp/pc/impls/telescope/telescope.c</A><BR>
<A HREF="../../../src/ksp/pc/impls/telescope/telescope_dmda.c.html#PCApplyRichardson_Telescope_dmda">PCApplyRichardson_Telescope_dmda in src/ksp/pc/impls/telescope/telescope_dmda.c</A><BR>
<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>
</BODY></HTML>
|