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
|
<!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/PCASMSetOverlap.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>PCASMSetOverlap</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/PC/PCASMSetOverlap.html "><small>Report Typos and Errors</small></a></div>
<A NAME="PCASMSetOverlap"><H1>PCASMSetOverlap</H1></A>
Sets the overlap between a pair of subdomains for the additive Schwarz preconditioner. Either all or no processors in the <A HREF="../PC/PC.html#PC">PC</A> communicator must call this routine.
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscpc.h"
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../PC/PCASMSetOverlap.html#PCASMSetOverlap">PCASMSetOverlap</A>(<A HREF="../PC/PC.html#PC">PC</A> pc,<A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> ovl)
</PRE>
Logically Collective on pc
<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>ovl </B></TD><TD>- the amount of overlap between subdomains (ovl >= 0, default value = 1)
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Options Database Key</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_asm_overlap <ovl> </B></TD><TD>- Sets overlap
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
By default the ASM preconditioner uses 1 block per processor. To use
multiple blocks per perocessor, see <A HREF="../PC/PCASMSetTotalSubdomains.html#PCASMSetTotalSubdomains">PCASMSetTotalSubdomains</A>() and
<A HREF="../PC/PCASMSetLocalSubdomains.html#PCASMSetLocalSubdomains">PCASMSetLocalSubdomains</A>() (and the option -pc_asm_blocks <blks>).
<P>
The overlap defaults to 1, so if one desires that no additional
overlap be computed beyond what may have been set with a call to
<A HREF="../PC/PCASMSetTotalSubdomains.html#PCASMSetTotalSubdomains">PCASMSetTotalSubdomains</A>() or <A HREF="../PC/PCASMSetLocalSubdomains.html#PCASMSetLocalSubdomains">PCASMSetLocalSubdomains</A>(), then ovl
must be set to be 0. In particular, if one does not explicitly set
the subdomains an application code, then all overlap would be computed
internally by PETSc, and using an overlap of 0 would result in an ASM
variant that is equivalent to the block Jacobi preconditioner.
<P>
The default algorithm used by PETSc to increase overlap is fast, but not scalable,
use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
<P>
Note that one can define initial index sets with any overlap via
<A HREF="../PC/PCASMSetLocalSubdomains.html#PCASMSetLocalSubdomains">PCASMSetLocalSubdomains</A>(); the routine
<A HREF="../PC/PCASMSetOverlap.html#PCASMSetOverlap">PCASMSetOverlap</A>() merely allows PETSc to extend that overlap further
if desired.
<P>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../PC/PCASMSetTotalSubdomains.html#PCASMSetTotalSubdomains">PCASMSetTotalSubdomains</A>(), <A HREF="../PC/PCASMSetLocalSubdomains.html#PCASMSetLocalSubdomains">PCASMSetLocalSubdomains</A>(), <A HREF="../PC/PCASMGetSubKSP.html#PCASMGetSubKSP">PCASMGetSubKSP</A>(),
<BR><A HREF="../PC/PCASMCreateSubdomains2D.html#PCASMCreateSubdomains2D">PCASMCreateSubdomains2D</A>(), <A HREF="../PC/PCASMGetLocalSubdomains.html#PCASMGetLocalSubdomains">PCASMGetLocalSubdomains</A>(), <A HREF="../Mat/MatIncreaseOverlap.html#MatIncreaseOverlap">MatIncreaseOverlap</A>()
<P><B></B><H3><FONT COLOR="#CC3333">Level</FONT></H3>intermediate<BR>
<H3><FONT COLOR="#CC3333">Location</FONT></H3>
</B><A HREF="../../../src/ksp/pc/impls/asm/asm.c.html#PCASMSetOverlap">src/ksp/pc/impls/asm/asm.c</A>
<P><H3><FONT COLOR="#CC3333">Examples</FONT></H3>
<A HREF="../../../src/ksp/ksp/tutorials/ex8.c.html">src/ksp/ksp/tutorials/ex8.c.html</A><BR>
<P><H3><FONT COLOR="CC3333">Implementations</FONT></H3><A HREF="../../../src/ksp/pc/impls/asm/asm.c.html#PCASMSetOverlap_ASM">PCASMSetOverlap_ASM in src/ksp/pc/impls/asm/asm.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>
|