1 static char help[] = "Test query functions for DMNetwork \n\n";
2
3 #include <petscdmnetwork.h>
4
5 /*
6 CreateStarGraphEdgeList - Create a k-Star Graph Edgelist on current processor
7 Not Collective
8
9 Input Parameters:
10 . k - order of the star graph (number of edges)
11 . directin - if true direction of edges is towards the center vertex, otherwise they are directed out of the center vertex.
12
13 Output Parameters:
14 . ne - number of edges of this star graph
15 . edgelist - list of edges for this star graph, this is a one dimensional array with pairs of entries being the two vertices (in global numbering of the vertices) of each edge,
16 [first vertex of first edge, second vertex of first edge, first vertex of second edge, second vertex of second edge, etc].
17
18 User is responsible for deallocating this memory.
19 */
StarGraphCreateEdgeList(PetscInt k,PetscBool directin,PetscInt * ne,PetscInt * edgelist[])20 PetscErrorCode StarGraphCreateEdgeList(PetscInt k, PetscBool directin, PetscInt *ne, PetscInt *edgelist[])
21 {
22 PetscInt i;
23
24 PetscFunctionBegin;
25 *ne = k;
26 PetscCall(PetscCalloc1(2 * k, edgelist));
27
28 if (directin) {
29 for (i = 0; i < k; i++) {
30 (*edgelist)[2 * i] = i + 1;
31 (*edgelist)[2 * i + 1] = 0;
32 }
33 } else {
34 for (i = 0; i < k; i++) {
35 (*edgelist)[2 * i] = 0;
36 (*edgelist)[2 * i + 1] = i + 1;
37 }
38 }
39 PetscFunctionReturn(PETSC_SUCCESS);
40 }
41
42 /*
43 CreateSimpleStarGraph - Create a Distributed k-Star Graph DMNetwork with a single PetscInt component on
44 all edges and vertices, a selectable number of dofs on vertices and edges. Intended mostly to be used for testing purposes.
45
46 Input Parameters:
47 . comm - the communicator of the dm
48 . numdofvert - number of degrees of freedom (dofs) on vertices
49 . numdofedge - number of degrees of freedom (dofs) on edges
50 . k - order of the star graph (number of edges)
51 . directin - if true direction of edges is towards the center vertex, otherwise they are directed out of the center vertex
52
53 Output Parameter:
54 . newdm - The created and distributed simple Star Graph
55 */
StarGraphCreate(MPI_Comm comm,PetscInt numdofvert,PetscInt numdofedge,PetscInt k,PetscBool directin,DM * newdm)56 PetscErrorCode StarGraphCreate(MPI_Comm comm, PetscInt numdofvert, PetscInt numdofedge, PetscInt k, PetscBool directin, DM *newdm)
57 {
58 DM dm;
59 PetscMPIInt rank;
60 PetscInt ne = 0, compkey, eStart, eEnd, vStart, vEnd, e, v;
61 PetscInt *edgelist = NULL, *compedge, *compvert;
62
63 PetscFunctionBegin;
64 PetscCall(DMNetworkCreate(comm, &dm));
65 PetscCall(DMNetworkSetNumSubNetworks(dm, PETSC_DECIDE, 1));
66 PetscCallMPI(MPI_Comm_rank(comm, &rank));
67 if (rank == 0) PetscCall(StarGraphCreateEdgeList(k, directin, &ne, &edgelist));
68 PetscCall(DMNetworkAddSubnetwork(dm, "Main", ne, edgelist, NULL));
69 PetscCall(DMNetworkRegisterComponent(dm, "unused", sizeof(PetscInt), &compkey));
70 PetscCall(DMNetworkLayoutSetUp(dm));
71 PetscCall(PetscFree(edgelist));
72 PetscCall(DMNetworkGetEdgeRange(dm, &eStart, &eEnd));
73 PetscCall(DMNetworkGetVertexRange(dm, &vStart, &vEnd));
74 PetscCall(PetscMalloc2(eEnd - eStart, &compedge, vEnd - vStart, &compvert));
75 for (e = eStart; e < eEnd; e++) {
76 compedge[e - eStart] = e;
77 PetscCall(DMNetworkAddComponent(dm, e, compkey, &compedge[e - eStart], numdofedge));
78 }
79 for (v = vStart; v < vEnd; v++) {
80 compvert[v - vStart] = v;
81 PetscCall(DMNetworkAddComponent(dm, v, compkey, &compvert[v - vStart], numdofvert));
82 }
83 PetscCall(DMSetFromOptions(dm));
84 PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
85 PetscCall(DMSetUp(dm));
86 PetscCall(PetscFree2(compedge, compvert));
87 *newdm = dm;
88 PetscFunctionReturn(PETSC_SUCCESS);
89 }
90
StarGraphTestQuery(DM dm,PetscInt ne)91 PetscErrorCode StarGraphTestQuery(DM dm, PetscInt ne)
92 {
93 PetscInt globalnumvert, localnumvert, globalnumedge, localnumedge;
94
95 PetscFunctionBegin;
96 PetscCall(DMNetworkGetNumEdges(dm, &localnumedge, &globalnumedge));
97 PetscCall(DMNetworkGetNumVertices(dm, &localnumvert, &globalnumvert));
98
99 PetscCheck(globalnumedge == ne, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Global number of edges should be %" PetscInt_FMT "instead was %" PetscInt_FMT, ne, globalnumedge);
100 PetscCheck(globalnumvert == ne + 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Global number of vertices should be %" PetscInt_FMT "instead was %" PetscInt_FMT, ne + 1, globalnumvert);
101 PetscFunctionReturn(PETSC_SUCCESS);
102 }
103
main(int argc,char ** argv)104 int main(int argc, char **argv)
105 {
106 DM dm;
107 PetscInt ne = 1;
108 PetscMPIInt rank;
109
110 PetscFunctionBeginUser;
111 PetscCall(PetscInitialize(&argc, &argv, NULL, help));
112 PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank));
113 /* create a distributed k-Star graph DMNetwork */
114 PetscCall(PetscOptionsGetInt(NULL, NULL, "-ne", &ne, NULL));
115 PetscCall(StarGraphCreate(PETSC_COMM_WORLD, 1, 0, ne, PETSC_TRUE, &dm));
116 PetscCall(DMNetworkDistribute(&dm, 0));
117 /* Test if query functions for DMNetwork run successfully */
118 PetscCall(StarGraphTestQuery(dm, ne));
119 PetscCall(DMDestroy(&dm));
120 PetscCall(PetscFinalize());
121 }
122
123 /*TEST
124 test:
125 suffix: 0
126 args: -ne 5
127 output_file: output/empty.out
128 test:
129 suffix: 1
130 nsize: 2
131 args: -ne 5 -petscpartitioner_type simple
132 output_file: output/empty.out
133 TEST*/
134