xref: /petsc/src/dm/impls/network/tests/ex2.c (revision 37d05b0256c1e9ba4bc423c4eccb3df226931ef0)
1 static char help[] = "Test DMCreateCoordinateDM_Network, and related functions \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 */
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 Parameters:
54 . newdm       - The created and distributed simple Star Graph
55 */
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, "dummy", 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 
91 /* Simple Circular embedding of the star graph */
92 PetscErrorCode StarGraphSetCoordinates(DM dm)
93 {
94   DM           cdm;
95   Vec          Coord;
96   PetscScalar *coord;
97   PetscInt     vStart, vEnd, v, vglobal, compkey, off, NVert;
98   PetscReal    theta;
99 
100   PetscFunctionBegin;
101   PetscCall(DMGetCoordinateDM(dm, &cdm));
102   PetscCall(DMSetCoordinateDim(dm, 2));
103   PetscCall(DMNetworkGetVertexRange(cdm, &vStart, &vEnd));
104   PetscCall(DMNetworkRegisterComponent(cdm, "coordinates", 0, &compkey));
105   for (v = vStart; v < vEnd; v++) { PetscCall(DMNetworkAddComponent(cdm, v, compkey, NULL, 2)); }
106   PetscCall(DMNetworkFinalizeComponents(cdm));
107 
108   PetscCall(DMCreateLocalVector(cdm, &Coord));
109   PetscCall(VecGetArray(Coord, &coord));
110   PetscCall(DMNetworkGetNumVertices(cdm, NULL, &NVert));
111   theta = 2 * PETSC_PI / (NVert - 1);
112   for (v = vStart; v < vEnd; v++) {
113     PetscCall(DMNetworkGetGlobalVertexIndex(cdm, v, &vglobal));
114     PetscCall(DMNetworkGetLocalVecOffset(cdm, v, 0, &off));
115     if (vglobal == 0) {
116       coord[off]     = 0.0;
117       coord[off + 1] = 0.0;
118     } else {
119       /* embed on the unit circle */
120       coord[off]     = PetscCosReal(theta * (vglobal - 1));
121       coord[off + 1] = PetscSinReal(theta * (vglobal - 1));
122     }
123   }
124   PetscCall(VecRestoreArray(Coord, &coord));
125   PetscCall(DMSetCoordinatesLocal(dm, Coord));
126   PetscCall(VecDestroy(&Coord));
127   PetscFunctionReturn(PETSC_SUCCESS);
128 }
129 
130 int main(int argc, char **argv)
131 {
132   DM           dm, cdm;
133   PetscInt     dofv = 1, dofe = 1, ne = 1, cdim, v, vStart, vEnd, off, vglobal;
134   Vec          Coord;
135   PetscScalar *coord;
136   PetscMPIInt  rank;
137   PetscBool    testdistribute = PETSC_FALSE;
138 
139   PetscFunctionBeginUser;
140   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
141   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank));
142   /* create a distributed k-Star graph DMNetwork */
143   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dofv", &dofv, NULL));
144   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dofe", &dofe, NULL));
145   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ne", &ne, NULL));
146   PetscCall(PetscOptionsGetBool(NULL, NULL, "-testdistribute", &testdistribute, NULL));
147   PetscCall(StarGraphCreate(PETSC_COMM_WORLD, dofv, dofe, ne, PETSC_TRUE, &dm));
148 
149   /* setup a quick R^2 embedding of the star graph */
150   PetscCall(StarGraphSetCoordinates(dm));
151 
152   if (testdistribute) {
153     PetscCall(DMNetworkDistribute(&dm, 0));
154     PetscCall(DMView(dm, PETSC_VIEWER_STDOUT_WORLD));
155   }
156 
157   /* print the coordinates of each vertex */
158   PetscCall(DMGetCoordinateDim(dm, &cdim));
159   PetscCall(DMGetCoordinateDM(dm, &cdm));
160   PetscCall(DMNetworkGetVertexRange(cdm, &vStart, &vEnd));
161   PetscCall(DMGetCoordinatesLocal(dm, &Coord));
162   PetscCall(VecGetArray(Coord, &coord));
163   PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "\n Rank %i \n\n", rank));
164   for (v = vStart; v < vEnd; v++) {
165     PetscCall(DMNetworkGetLocalVecOffset(cdm, v, 0, &off));
166     PetscCall(DMNetworkGetGlobalVertexIndex(cdm, v, &vglobal));
167     switch (cdim) {
168     case 2:
169       PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "Vertex: %" PetscInt_FMT ", x =  %f y = %f \n", vglobal, (double)PetscRealPart(coord[off]), (double)PetscRealPart(coord[off + 1])));
170       break;
171     default:
172       SETERRQ(MPI_COMM_WORLD, PETSC_ERR_SUP, "Only supports Network embedding dimension of 2, not supplied  %" PetscInt_FMT, cdim);
173       break;
174     }
175   }
176   PetscCall(PetscSynchronizedFlush(MPI_COMM_WORLD, NULL));
177 
178   PetscCall(DMDestroy(&dm));
179   PetscCall(PetscFinalize());
180 }
181 
182 /*TEST
183 
184   test:
185     suffix: 0
186     args: -ne 4 -testdistribute
187 
188   test:
189     suffix: 1
190     nsize: 2
191     args: -ne 4 -testdistribute -petscpartitioner_type simple
192 
193  TEST*/
194