xref: /petsc/src/dm/impls/network/tests/ex1.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
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 {
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(0);
40 }
41 
42 /*
43 CreateSimpleStarGraph - Create a Distributed k-Star Graph DMNetwork with a single PetscInt component on
44 all edges and vertices, aselectable 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 CreateSimpleStarGraph(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) {
68     PetscCall(CreateStarGraphEdgeList(k,directin,&ne,&edgelist));
69   }
70   PetscCall(DMNetworkAddSubnetwork(dm,"Main",ne,edgelist,NULL));
71   PetscCall(DMNetworkRegisterComponent(dm,"dummy",sizeof(PetscInt),&compkey));
72   PetscCall(DMNetworkLayoutSetUp(dm));
73   PetscCall(PetscFree(edgelist));
74   PetscCall(DMNetworkGetEdgeRange(dm,&eStart,&eEnd));
75   PetscCall(DMNetworkGetVertexRange(dm,&vStart,&vEnd));
76   PetscCall(PetscMalloc2(eEnd-eStart,&compedge,vEnd-vStart,&compvert));
77   for (e=eStart; e<eEnd; e++) {
78     compedge[e-eStart] = e;
79     PetscCall(DMNetworkAddComponent(dm,e,compkey,&compedge[e-eStart],numdofedge));
80   }
81   for (v=vStart; v<vEnd; v++) {
82     compvert[v-vStart] = v;
83     PetscCall(DMNetworkAddComponent(dm,v,compkey,&compvert[v-vStart],numdofvert));
84   }
85   PetscCall(DMSetFromOptions(dm));
86   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
87   PetscCall(DMSetUp(dm));
88   PetscCall(PetscFree2(compedge,compvert));
89   PetscCall(DMNetworkDistribute(&dm,0));
90   *newdm = dm;
91   PetscFunctionReturn(0);
92 }
93 
94 int main(int argc, char **argv)
95 {
96   DM             dm, dmclone,plex;
97   PetscInt       e,eStart,eEnd,ndofs,ndofsprev;
98   PetscInt       *compprev,*comp,compkey;
99   PetscInt       dofv=1,dofe=1,ne=1;
100   PetscSection   sec;
101   Vec            vec;
102 
103   PetscFunctionBeginUser;
104   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
105   /* create a distributed k-Star graph DMNetwork */
106   PetscCall(PetscOptionsGetInt(NULL,NULL,"-dofv",&dofv,NULL));
107   PetscCall(PetscOptionsGetInt(NULL,NULL,"-dofe",&dofe,NULL));
108   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ne",&ne,NULL));
109   PetscCall(CreateSimpleStarGraph(PETSC_COMM_WORLD,dofv,dofe,ne,PETSC_TRUE,&dm));
110   PetscCall(DMNetworkGetEdgeRange(dm,&eStart,&eEnd));
111 
112   /* check if cloning changed any componenent */
113   if (eStart < eEnd) {
114     PetscCall(DMNetworkGetComponent(dm,eStart,0,NULL,(void**)&compprev,&ndofsprev));
115   }
116   PetscCall(DMClone(dm,&dmclone));
117   if (eStart < eEnd) {
118     PetscCall(DMNetworkGetComponent(dm,eStart,0,NULL,(void**)&comp,&ndofs));
119     PetscCheck(*comp == *compprev,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Cloning changed the Original, comp (previous) : %" PetscInt_FMT " comp (now) : %" PetscInt_FMT,*compprev,*comp);
120     PetscCheck(ndofsprev == ndofs,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Cloning changed the Original, ndofs (previous) : %" PetscInt_FMT " ndofs (now) : %" PetscInt_FMT,ndofsprev,ndofs);
121   }
122 
123   /* register new components to the clone and add a dummy component to every point */
124   PetscCall(DMNetworkRegisterComponent(dmclone,"dummyclone",sizeof(PetscInt),&compkey));
125   PetscCall(DMNetworkGetEdgeRange(dmclone,&eStart,&eEnd));
126   PetscCall(PetscMalloc1(eEnd-eStart,&comp));
127   for (e=eStart; e<eEnd; e++) {
128     comp[e-eStart] = e;
129     PetscCall(DMNetworkAddComponent(dmclone,e,compkey,&comp[e-eStart],2));
130   }
131   PetscCall(DMNetworkFinalizeComponents(dmclone));
132   PetscCall(PetscFree(comp));
133 
134   PetscCall(PetscPrintf(PETSC_COMM_WORLD," dm: \n"));
135   PetscCall(DMView(dm,PETSC_VIEWER_STDOUT_WORLD));
136   PetscCall(DMNetworkGetPlex(dm,&plex));
137   PetscCall(DMGetLocalSection(plex,&sec));
138   PetscCall(PetscSectionView(sec,PETSC_VIEWER_STDOUT_WORLD));
139 
140   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n dmclone: \n"));
141   PetscCall(DMView(dmclone,PETSC_VIEWER_STDOUT_WORLD));
142   PetscCall(DMNetworkGetPlex(dmclone,&plex));
143   PetscCall(DMGetLocalSection(plex,&sec));
144   PetscCall(PetscSectionView(sec,PETSC_VIEWER_STDOUT_WORLD));
145 
146   /* create Vectors */
147   PetscCall(DMCreateGlobalVector(dm,&vec));
148   PetscCall(VecSet(vec,1.0));
149   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n dm vec:\n"));
150   PetscCall(VecView(vec,PETSC_VIEWER_STDOUT_WORLD));
151   PetscCall(VecDestroy(&vec));
152 
153   PetscCall(DMCreateGlobalVector(dmclone,&vec));
154   PetscCall(VecSet(vec,2.0));
155   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n dmclone vec:\n"));
156   PetscCall(VecView(vec,PETSC_VIEWER_STDOUT_WORLD));
157   PetscCall(VecDestroy(&vec));
158 
159   PetscCall(DMDestroy(&dm));
160   PetscCall(DMDestroy(&dmclone));
161   PetscCall(PetscFinalize());
162 }
163 
164 /*TEST
165 
166   test:
167     suffix: 0
168     args:
169 
170   test:
171     suffix: 1
172     nsize: 2
173     args: -dofv 2 -dofe 2 -ne 2
174 
175  TEST*/
176