File: PetscSFGetSubSF.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 (111 lines) | stat: -rw-r--r-- 8,532 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
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMNetwork/PetscSFGetSubSF.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>PetscSFGetSubSF</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/DMNetwork/PetscSFGetSubSF.html "><small>Report Typos and Errors</small></a></div>
<A NAME="PetscSFGetSubSF"><H1>PetscSFGetSubSF</H1></A>
Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering. 
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscdmnetwork.h"  
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> DMNetworkGetSupportingEdges(<A HREF="../DM/DM.html#DM">DM</A> dm,<A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> vertex,<A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> *nedges,const <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> *edges[])
</PRE>
<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>masterSF </B></TD><TD>- the original SF structure
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>map      </B></TD><TD>- a ISLocalToGlobal mapping that contains the subset of points
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Output Parameters</FONT></H3>
<DT><B>subSF    </B> -a subset of the masterSF for the desired subset.
*/
<br>
<P>
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../DMNetwork/PetscSFGetSubSF.html#PetscSFGetSubSF">PetscSFGetSubSF</A>(<A HREF="../PetscSF/PetscSF.html#PetscSF">PetscSF</A> mastersf, <A HREF="../IS/ISLocalToGlobalMapping.html#ISLocalToGlobalMapping">ISLocalToGlobalMapping</A> map, <A HREF="../PetscSF/PetscSF.html#PetscSF">PetscSF</A> *subSF) {
<P>
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A>        ierr;
<A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A>              nroots, nleaves, *ilocal_sub;
<A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A>              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
<A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A>              *local_points, *remote_points;
<A HREF="../PetscSF/PetscSFNode.html#PetscSFNode">PetscSFNode</A>           *iremote_sub;
const <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A>        *ilocal;
const <A HREF="../PetscSF/PetscSFNode.html#PetscSFNode">PetscSFNode</A>     *iremote;
<P>
<A HREF="../Sys/PetscFunctionBegin.html#PetscFunctionBegin">PetscFunctionBegin</A>;
ierr = <A HREF="../PetscSF/PetscSFGetGraph.html#PetscSFGetGraph">PetscSFGetGraph</A>(mastersf,&amp;nroots,&amp;nleaves,&amp;ilocal,&amp;iremote);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
<P>
/* Look for leaves that pertain to the subset of points. Get the local ordering */
ierr = <A HREF="../Sys/PetscMalloc1.html#PetscMalloc1">PetscMalloc1</A>(nleaves,&amp;ilocal_map);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr = <A HREF="../IS/ISGlobalToLocalMappingApply.html#ISGlobalToLocalMappingApply">ISGlobalToLocalMappingApply</A>(map,<A HREF="../IS/ISGlobalToLocalMappingMode.html#ISGlobalToLocalMappingMode">IS_GTOLM_MASK</A>,nleaves,ilocal,NULL,ilocal_map);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
for (i = 0; i &lt; nleaves; i++) {
if (ilocal_map[i] != -1) nleaves_sub += 1;
}
/* Re-number ilocal with subset numbering. Need information from roots */
ierr = <A HREF="../Sys/PetscMalloc2.html#PetscMalloc2">PetscMalloc2</A>(nroots,&amp;local_points,nroots,&amp;remote_points);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
for (i = 0; i &lt; nroots; i++) local_points[i] = i;
ierr = <A HREF="../IS/ISGlobalToLocalMappingApply.html#ISGlobalToLocalMappingApply">ISGlobalToLocalMappingApply</A>(map,<A HREF="../IS/ISGlobalToLocalMappingMode.html#ISGlobalToLocalMappingMode">IS_GTOLM_MASK</A>,nroots,local_points,NULL,local_points);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr = <A HREF="../PetscSF/PetscSFBcastBegin.html#PetscSFBcastBegin">PetscSFBcastBegin</A>(mastersf, <A HREF="../Sys/MPIU_INT.html#MPIU_INT">MPIU_INT</A>, local_points, remote_points);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr = <A HREF="../PetscSF/PetscSFBcastEnd.html#PetscSFBcastEnd">PetscSFBcastEnd</A>(mastersf, <A HREF="../Sys/MPIU_INT.html#MPIU_INT">MPIU_INT</A>, local_points, remote_points);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
/* Fill up graph using local (that is, local to the subset) numbering. */
ierr = <A HREF="../Sys/PetscMalloc1.html#PetscMalloc1">PetscMalloc1</A>(nleaves_sub,&amp;ilocal_sub);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr = <A HREF="../Sys/PetscMalloc1.html#PetscMalloc1">PetscMalloc1</A>(nleaves_sub,&amp;iremote_sub);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
nleaves_sub = 0;
for (i = 0; i &lt; nleaves; i++) {
if (ilocal_map[i] != -1) {
ilocal_sub[nleaves_sub] = ilocal_map[i];
iremote_sub[nleaves_sub].rank = iremote[i].rank;
iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
nleaves_sub += 1;
}
}
ierr = <A HREF="../Sys/PetscFree2.html#PetscFree2">PetscFree2</A>(local_points,remote_points);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr = <A HREF="../IS/ISLocalToGlobalMappingGetSize.html#ISLocalToGlobalMappingGetSize">ISLocalToGlobalMappingGetSize</A>(map,&amp;nroots_sub);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
<P>
/* Create new subSF */
ierr = <A HREF="../PetscSF/PetscSFCreate.html#PetscSFCreate">PetscSFCreate</A>(<A HREF="../Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</A>,subSF);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr = <A HREF="../PetscSF/PetscSFSetFromOptions.html#PetscSFSetFromOptions">PetscSFSetFromOptions</A>(*subSF);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr = <A HREF="../PetscSF/PetscSFSetGraph.html#PetscSFSetGraph">PetscSFSetGraph</A>(*subSF,nroots_sub,nleaves_sub,ilocal_sub,<A HREF="../Sys/PetscCopyMode.html#PetscCopyMode">PETSC_OWN_POINTER</A>,iremote_sub,<A HREF="../Sys/PetscCopyMode.html#PetscCopyMode">PETSC_COPY_VALUES</A>);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr = <A HREF="../Sys/PetscFree.html#PetscFree">PetscFree</A>(ilocal_map);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr = <A HREF="../Sys/PetscFree.html#PetscFree">PetscFree</A>(iremote_sub);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
<A HREF="../Sys/PetscFunctionReturn.html#PetscFunctionReturn">PetscFunctionReturn</A>(0);
}
<P>
/*@C
DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
<P>
Not 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>dm </B></TD><TD>- The DMNetwork object
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>p  </B></TD><TD>- the vertex point
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Output Paramters</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>nedges </B></TD><TD>- number of edges connected to this vertex point
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>edges  </B></TD><TD>- List of edge points
</TD></TR></TABLE>
<P>

<P>
<H3><FONT COLOR="#CC3333">Fortran Notes</FONT></H3>
Since it returns an array, this routine is only available in Fortran 90, and you must
include petsc.h90 in your code.
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
 <A HREF="../DMNetwork/DMNetworkCreate.html#DMNetworkCreate">DMNetworkCreate</A>, <A HREF="../DMNetwork/DMNetworkGetConnectedVertices.html#DMNetworkGetConnectedVertices">DMNetworkGetConnectedVertices</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/dm/impls/network/network.c.html#PetscSFGetSubSF">src/dm/impls/network/network.c</A>
<BR><A HREF="./index.html">Index of all DMNetwork 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>