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
|
<!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/PCGASMCreateSubdomains.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>PCGASMCreateSubdomains</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.7.5 2017-01-01</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.7.5 v3.7.5 docs/manualpages/PC/PCGASMCreateSubdomains.html "><small>Report Typos and Errors</small></a></div>
<A NAME="PCGASMCreateSubdomains"><H1>PCGASMCreateSubdomains</H1></A>
Creates n index sets defining n nonoverlapping subdomains for the additive Schwarz preconditioner for a any problem based on its matrix.
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscpc.h"
PetscErrorCode PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[])
</PRE>
Collective
<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>A </B></TD><TD>- The global matrix operator
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>N </B></TD><TD>- the number of global subdomains requested
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Output Parameters</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>n </B></TD><TD>- the number of subdomains created on this processor
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>iis </B></TD><TD>- the array of index sets defining the local inner subdomains (on which the correction is applied)
</TD></TR></TABLE>
<P>
<P>
Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor.
When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local.
The resulting subdomains can be use in <A HREF="../PC/PCGASMSetSubdomains.html#PCGASMSetSubdomains">PCGASMSetSubdomains</A>(pc,n,iss,NULL). The overlapping
outer subdomains will be automatically generated from these according to the requested amount of
overlap; this is currently supported only with local subdomains.
<P>
<P>
<H3><FONT COLOR="#CC3333">Keywords</FONT></H3>
<A HREF="../PC/PC.html#PC">PC</A>, GASM, additive Schwarz, create, subdomains, unstructured grid
<BR>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../PC/PCGASMSetSubdomains.html#PCGASMSetSubdomains">PCGASMSetSubdomains</A>(), <A HREF="../PC/PCGASMDestroySubdomains.html#PCGASMDestroySubdomains">PCGASMDestroySubdomains</A>()
<BR><P><B><P><B><FONT COLOR="#CC3333">Level:</FONT></B>advanced
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/ksp/pc/impls/gasm/gasm.c.html#PCGASMCreateSubdomains">src/ksp/pc/impls/gasm/gasm.c</A>
<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>
|