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