File: PCASMSetLocalSubdomains.html

package info (click to toggle)
petsc 3.10.3%2Bdfsg1-5
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 209,064 kB
  • sloc: ansic: 587,333; python: 29,696; makefile: 12,445; fortran: 11,626; f90: 9,677; cpp: 8,768; sh: 1,027; xml: 621; objc: 445; csh: 194; java: 13
file content (70 lines) | stat: -rw-r--r-- 5,075 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
<!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/PCASMSetLocalSubdomains.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>PCASMSetLocalSubdomains</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/PCASMSetLocalSubdomains.html "><small>Report Typos and Errors</small></a></div>
<A NAME="PCASMSetLocalSubdomains"><H1>PCASMSetLocalSubdomains</H1></A>
Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner. 
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscpc.h" 
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A>  <A HREF="../PC/PCASMSetLocalSubdomains.html#PCASMSetLocalSubdomains">PCASMSetLocalSubdomains</A>(<A HREF="../PC/PC.html#PC">PC</A> pc,<A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> n,<A HREF="../IS/IS.html#IS">IS</A> is[],<A HREF="../IS/IS.html#IS">IS</A> is_local[])
</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>n </B></TD><TD>- the number of subdomains for this processor (default value = 1)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>is </B></TD><TD>- the index set that defines the subdomains for this processor
(or NULL for PETSc to determine subdomains)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>is_local </B></TD><TD>- the index sets that define the local part of the subdomains for this processor, not used unless <A HREF="../PC/PCASMType.html#PCASMType">PCASMType</A> is <A HREF="../PC/PCASMType.html#PCASMType">PC_ASM_RESTRICT</A>
(or NULL to not provide these)
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Options Database Key</FONT></H3>
To set the total number of subdomain blocks rather than specify the
index sets, use the option
<DT><B>-pc_asm_local_blocks &lt;blks&gt; </B> -Sets local blocks
<br>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
The <A HREF="../IS/IS.html#IS">IS</A> numbering is in the parallel, global numbering of the vector for both is and is_local
<P>
By default the ASM preconditioner uses 1 block per processor.
<P>
Use <A HREF="../PC/PCASMSetTotalSubdomains.html#PCASMSetTotalSubdomains">PCASMSetTotalSubdomains</A>() to set the subdomains for all processors.
<P>
If is_local is provided and <A HREF="../PC/PCASMType.html#PCASMType">PCASMType</A> is <A HREF="../PC/PCASMType.html#PCASMType">PC_ASM_RESTRICT</A> then the solution only over the is_local region is interpolated
back to form the global solution (this is the standard restricted additive Schwarz method)
<P>
If the is_local is provided and <A HREF="../PC/PCASMType.html#PCASMType">PCASMType</A> is <A HREF="../PC/PCASMType.html#PCASMType">PC_ASM_INTERPOLATE</A> or <A HREF="../PC/PCASMType.html#PCASMType">PC_ASM_NONE</A> then an error is generated since there is
no code to handle that case.
<P>

<P>
<H3><FONT COLOR="#CC3333">Keywords</FONT></H3>
 <A HREF="../PC/PC.html#PC">PC</A>, ASM, set, local, subdomains, additive Schwarz
<BR>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
 <A HREF="../PC/PCASMSetTotalSubdomains.html#PCASMSetTotalSubdomains">PCASMSetTotalSubdomains</A>(), <A HREF="../PC/PCASMSetOverlap.html#PCASMSetOverlap">PCASMSetOverlap</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="../PC/PCASMType.html#PCASMType">PCASMType</A>, <A HREF="../PC/PCASMSetType.html#PCASMSetType">PCASMSetType</A>()
<P><B></B><H3><FONT COLOR="#CC3333">Level</FONT></H3>advanced<BR>
<H3><FONT COLOR="#CC3333">Location</FONT></H3>
</B><A HREF="../../../src/ksp/pc/impls/asm/asm.c.html#PCASMSetLocalSubdomains">src/ksp/pc/impls/asm/asm.c</A>
<P><H3><FONT COLOR="#CC3333">Examples</FONT></H3>
<A HREF="../../../src/ksp/ksp/examples/tutorials/ex6.c.html">src/ksp/ksp/examples/tutorials/ex6.c.html</A><BR>
<A HREF="../../../src/ksp/ksp/examples/tutorials/ex8.c.html">src/ksp/ksp/examples/tutorials/ex8.c.html</A><BR>
<P><H3><FONT COLOR="CC3333">Implementations</FONT></H3><A HREF="../../../src/ksp/pc/impls/asm/asm.c.html#PCASMSetLocalSubdomains_ASM">PCASMSetLocalSubdomains_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>