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