xref: /petsc/src/dm/impls/network/tests/ex1.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 static char help[] = "Check if DMClone for DMNetwork Correctly Shallow Clones Topology Only \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 CreateStarGraphEdgeList(PetscInt k, PetscBool directin, PetscInt *ne, PetscInt *edgelist[]) {
21   PetscInt i;
22 
23   PetscFunctionBegin;
24   *ne = k;
25   PetscCall(PetscCalloc1(2 * k, edgelist));
26 
27   if (directin) {
28     for (i = 0; i < k; i++) {
29       (*edgelist)[2 * i]     = i + 1;
30       (*edgelist)[2 * i + 1] = 0;
31     }
32   } else {
33     for (i = 0; i < k; i++) {
34       (*edgelist)[2 * i]     = 0;
35       (*edgelist)[2 * i + 1] = i + 1;
36     }
37   }
38   PetscFunctionReturn(0);
39 }
40 
41 /*
42 CreateSimpleStarGraph - Create a Distributed k-Star Graph DMNetwork with a single PetscInt component on
43 all edges and vertices, aselectable number of dofs on vertices and edges. Intended mostly to be used for testing purposes.
44 
45   Input Parameters:
46 . comm       - the communicator of the dm
47 . numdofvert - number of degrees of freedom (dofs) on vertices
48 . numdofedge - number of degrees of freedom (dofs) on edges
49 . k          - order of the star graph (number of edges)
50 . directin   - if true direction of edges is towards the center vertex, otherwise they are directed out of the center vertex
51 
52   Output Parameters:
53 . newdm       - The created and distributed simple Star Graph
54 */
55 PetscErrorCode CreateSimpleStarGraph(MPI_Comm comm, PetscInt numdofvert, PetscInt numdofedge, PetscInt k, PetscBool directin, DM *newdm) {
56   DM          dm;
57   PetscMPIInt rank;
58   PetscInt    ne       = 0, compkey, eStart, eEnd, vStart, vEnd, e, v;
59   PetscInt   *edgelist = NULL, *compedge, *compvert;
60 
61   PetscFunctionBegin;
62   PetscCall(DMNetworkCreate(comm, &dm));
63   PetscCall(DMNetworkSetNumSubNetworks(dm, PETSC_DECIDE, 1));
64   PetscCallMPI(MPI_Comm_rank(comm, &rank));
65   if (rank == 0) { PetscCall(CreateStarGraphEdgeList(k, directin, &ne, &edgelist)); }
66   PetscCall(DMNetworkAddSubnetwork(dm, "Main", ne, edgelist, NULL));
67   PetscCall(DMNetworkRegisterComponent(dm, "dummy", sizeof(PetscInt), &compkey));
68   PetscCall(DMNetworkLayoutSetUp(dm));
69   PetscCall(PetscFree(edgelist));
70   PetscCall(DMNetworkGetEdgeRange(dm, &eStart, &eEnd));
71   PetscCall(DMNetworkGetVertexRange(dm, &vStart, &vEnd));
72   PetscCall(PetscMalloc2(eEnd - eStart, &compedge, vEnd - vStart, &compvert));
73   for (e = eStart; e < eEnd; e++) {
74     compedge[e - eStart] = e;
75     PetscCall(DMNetworkAddComponent(dm, e, compkey, &compedge[e - eStart], numdofedge));
76   }
77   for (v = vStart; v < vEnd; v++) {
78     compvert[v - vStart] = v;
79     PetscCall(DMNetworkAddComponent(dm, v, compkey, &compvert[v - vStart], numdofvert));
80   }
81   PetscCall(DMSetFromOptions(dm));
82   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
83   PetscCall(DMSetUp(dm));
84   PetscCall(PetscFree2(compedge, compvert));
85   PetscCall(DMNetworkDistribute(&dm, 0));
86   *newdm = dm;
87   PetscFunctionReturn(0);
88 }
89 
90 int main(int argc, char **argv) {
91   DM           dm, dmclone, plex;
92   PetscInt     e, eStart, eEnd, ndofs, ndofsprev;
93   PetscInt    *compprev, *comp, compkey;
94   PetscInt     dofv = 1, dofe = 1, ne = 1;
95   PetscSection sec;
96   Vec          vec;
97 
98   PetscFunctionBeginUser;
99   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
100   /* create a distributed k-Star graph DMNetwork */
101   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dofv", &dofv, NULL));
102   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dofe", &dofe, NULL));
103   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ne", &ne, NULL));
104   PetscCall(CreateSimpleStarGraph(PETSC_COMM_WORLD, dofv, dofe, ne, PETSC_TRUE, &dm));
105   PetscCall(DMNetworkGetEdgeRange(dm, &eStart, &eEnd));
106 
107   /* check if cloning changed any componenent */
108   if (eStart < eEnd) { PetscCall(DMNetworkGetComponent(dm, eStart, 0, NULL, (void **)&compprev, &ndofsprev)); }
109   PetscCall(DMClone(dm, &dmclone));
110   if (eStart < eEnd) {
111     PetscCall(DMNetworkGetComponent(dm, eStart, 0, NULL, (void **)&comp, &ndofs));
112     PetscCheck(*comp == *compprev, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Cloning changed the Original, comp (previous) : %" PetscInt_FMT " comp (now) : %" PetscInt_FMT, *compprev, *comp);
113     PetscCheck(ndofsprev == ndofs, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Cloning changed the Original, ndofs (previous) : %" PetscInt_FMT " ndofs (now) : %" PetscInt_FMT, ndofsprev, ndofs);
114   }
115 
116   /* register new components to the clone and add a dummy component to every point */
117   PetscCall(DMNetworkRegisterComponent(dmclone, "dummyclone", sizeof(PetscInt), &compkey));
118   PetscCall(DMNetworkGetEdgeRange(dmclone, &eStart, &eEnd));
119   PetscCall(PetscMalloc1(eEnd - eStart, &comp));
120   for (e = eStart; e < eEnd; e++) {
121     comp[e - eStart] = e;
122     PetscCall(DMNetworkAddComponent(dmclone, e, compkey, &comp[e - eStart], 2));
123   }
124   PetscCall(DMNetworkFinalizeComponents(dmclone));
125   PetscCall(PetscFree(comp));
126 
127   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " dm: \n"));
128   PetscCall(DMView(dm, PETSC_VIEWER_STDOUT_WORLD));
129   PetscCall(DMNetworkGetPlex(dm, &plex));
130   PetscCall(DMGetLocalSection(plex, &sec));
131   PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
132 
133   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n dmclone: \n"));
134   PetscCall(DMView(dmclone, PETSC_VIEWER_STDOUT_WORLD));
135   PetscCall(DMNetworkGetPlex(dmclone, &plex));
136   PetscCall(DMGetLocalSection(plex, &sec));
137   PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
138 
139   /* create Vectors */
140   PetscCall(DMCreateGlobalVector(dm, &vec));
141   PetscCall(VecSet(vec, 1.0));
142   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n dm vec:\n"));
143   PetscCall(VecView(vec, PETSC_VIEWER_STDOUT_WORLD));
144   PetscCall(VecDestroy(&vec));
145 
146   PetscCall(DMCreateGlobalVector(dmclone, &vec));
147   PetscCall(VecSet(vec, 2.0));
148   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n dmclone vec:\n"));
149   PetscCall(VecView(vec, PETSC_VIEWER_STDOUT_WORLD));
150   PetscCall(VecDestroy(&vec));
151 
152   PetscCall(DMDestroy(&dm));
153   PetscCall(DMDestroy(&dmclone));
154   PetscCall(PetscFinalize());
155 }
156 
157 /*TEST
158 
159   test:
160     suffix: 0
161     args:
162 
163   test:
164     suffix: 1
165     nsize: 2
166     args: -dofv 2 -dofe 2 -ne 2
167 
168  TEST*/
169