xref: /petsc/src/dm/impls/network/tests/ex2.c (revision cc13d4123af1c7f72a268c2d2542044f452cd6ac)
18afb7921SAidan Hamilton static char help[] = "Test DMCreateCoordinateDM_Network, and related functions \n\n";
28afb7921SAidan Hamilton 
38afb7921SAidan Hamilton #include <petscdmnetwork.h>
48afb7921SAidan Hamilton 
58afb7921SAidan Hamilton /*
68afb7921SAidan Hamilton CreateStarGraphEdgeList - Create a k-Star Graph Edgelist on current processor
78afb7921SAidan Hamilton   Not Collective
88afb7921SAidan Hamilton 
98afb7921SAidan Hamilton   Input Parameters:
108afb7921SAidan Hamilton . k    - order of the star graph (number of edges)
118afb7921SAidan Hamilton . directin - if true direction of edges is towards the center vertex, otherwise they are directed out of the center vertex.
128afb7921SAidan Hamilton 
138afb7921SAidan Hamilton   Output Parameters:
148afb7921SAidan Hamilton .  ne - number of edges of this star graph
158afb7921SAidan 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,
168afb7921SAidan Hamilton               [first vertex of first edge, second vertex of first edge, first vertex of second edge, second vertex of second edge, etc].
178afb7921SAidan Hamilton 
188afb7921SAidan Hamilton               User is responsible for deallocating this memory.
198afb7921SAidan Hamilton */
208afb7921SAidan Hamilton PetscErrorCode StarGraphCreateEdgeList(PetscInt k, PetscBool directin, PetscInt *ne, PetscInt *edgelist[])
218afb7921SAidan Hamilton {
228afb7921SAidan Hamilton   PetscInt i;
238afb7921SAidan Hamilton 
248afb7921SAidan Hamilton   PetscFunctionBegin;
258afb7921SAidan Hamilton   *ne = k;
268afb7921SAidan Hamilton   PetscCall(PetscCalloc1(2 * k, edgelist));
278afb7921SAidan Hamilton 
288afb7921SAidan Hamilton   if (directin) {
298afb7921SAidan Hamilton     for (i = 0; i < k; i++) {
308afb7921SAidan Hamilton       (*edgelist)[2 * i]     = i + 1;
318afb7921SAidan Hamilton       (*edgelist)[2 * i + 1] = 0;
328afb7921SAidan Hamilton     }
338afb7921SAidan Hamilton   } else {
348afb7921SAidan Hamilton     for (i = 0; i < k; i++) {
358afb7921SAidan Hamilton       (*edgelist)[2 * i]     = 0;
368afb7921SAidan Hamilton       (*edgelist)[2 * i + 1] = i + 1;
378afb7921SAidan Hamilton     }
388afb7921SAidan Hamilton   }
393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
408afb7921SAidan Hamilton }
418afb7921SAidan Hamilton 
428afb7921SAidan Hamilton /*
438afb7921SAidan Hamilton CreateSimpleStarGraph - Create a Distributed k-Star Graph DMNetwork with a single PetscInt component on
448afb7921SAidan Hamilton all edges and vertices, a selectable number of dofs on vertices and edges. Intended mostly to be used for testing purposes.
458afb7921SAidan Hamilton 
468afb7921SAidan Hamilton   Input Parameters:
478afb7921SAidan Hamilton . comm       - the communicator of the dm
488afb7921SAidan Hamilton . numdofvert - number of degrees of freedom (dofs) on vertices
498afb7921SAidan Hamilton . numdofedge - number of degrees of freedom (dofs) on edges
508afb7921SAidan Hamilton . k          - order of the star graph (number of edges)
518afb7921SAidan Hamilton . directin   - if true direction of edges is towards the center vertex, otherwise they are directed out of the center vertex
528afb7921SAidan Hamilton 
538afb7921SAidan Hamilton   Output Parameters:
548afb7921SAidan Hamilton . newdm       - The created and distributed simple Star Graph
558afb7921SAidan Hamilton */
568afb7921SAidan Hamilton PetscErrorCode StarGraphCreate(MPI_Comm comm, PetscInt numdofvert, PetscInt numdofedge, PetscInt k, PetscBool directin, DM *newdm)
578afb7921SAidan Hamilton {
588afb7921SAidan Hamilton   DM          dm;
598afb7921SAidan Hamilton   PetscMPIInt rank;
608afb7921SAidan Hamilton   PetscInt    ne       = 0, compkey, eStart, eEnd, vStart, vEnd, e, v;
618afb7921SAidan Hamilton   PetscInt   *edgelist = NULL, *compedge, *compvert;
628afb7921SAidan Hamilton 
638afb7921SAidan Hamilton   PetscFunctionBegin;
648afb7921SAidan Hamilton   PetscCall(DMNetworkCreate(comm, &dm));
658afb7921SAidan Hamilton   PetscCall(DMNetworkSetNumSubNetworks(dm, PETSC_DECIDE, 1));
668afb7921SAidan Hamilton   PetscCallMPI(MPI_Comm_rank(comm, &rank));
678afb7921SAidan Hamilton   if (rank == 0) PetscCall(StarGraphCreateEdgeList(k, directin, &ne, &edgelist));
688afb7921SAidan Hamilton   PetscCall(DMNetworkAddSubnetwork(dm, "Main", ne, edgelist, NULL));
698afb7921SAidan Hamilton   PetscCall(DMNetworkRegisterComponent(dm, "dummy", sizeof(PetscInt), &compkey));
708afb7921SAidan Hamilton   PetscCall(DMNetworkLayoutSetUp(dm));
718afb7921SAidan Hamilton   PetscCall(PetscFree(edgelist));
728afb7921SAidan Hamilton   PetscCall(DMNetworkGetEdgeRange(dm, &eStart, &eEnd));
738afb7921SAidan Hamilton   PetscCall(DMNetworkGetVertexRange(dm, &vStart, &vEnd));
748afb7921SAidan Hamilton   PetscCall(PetscMalloc2(eEnd - eStart, &compedge, vEnd - vStart, &compvert));
758afb7921SAidan Hamilton   for (e = eStart; e < eEnd; e++) {
768afb7921SAidan Hamilton     compedge[e - eStart] = e;
778afb7921SAidan Hamilton     PetscCall(DMNetworkAddComponent(dm, e, compkey, &compedge[e - eStart], numdofedge));
788afb7921SAidan Hamilton   }
798afb7921SAidan Hamilton   for (v = vStart; v < vEnd; v++) {
808afb7921SAidan Hamilton     compvert[v - vStart] = v;
818afb7921SAidan Hamilton     PetscCall(DMNetworkAddComponent(dm, v, compkey, &compvert[v - vStart], numdofvert));
828afb7921SAidan Hamilton   }
838afb7921SAidan Hamilton   PetscCall(DMSetFromOptions(dm));
848afb7921SAidan Hamilton   PetscCall(DMSetUp(dm));
858afb7921SAidan Hamilton   PetscCall(PetscFree2(compedge, compvert));
868afb7921SAidan Hamilton   *newdm = dm;
873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
888afb7921SAidan Hamilton }
898afb7921SAidan Hamilton 
908afb7921SAidan Hamilton /* Simple Circular embedding of the star graph */
918afb7921SAidan Hamilton PetscErrorCode StarGraphSetCoordinates(DM dm)
928afb7921SAidan Hamilton {
938afb7921SAidan Hamilton   DM           cdm;
948afb7921SAidan Hamilton   Vec          Coord;
958afb7921SAidan Hamilton   PetscScalar *coord;
968afb7921SAidan Hamilton   PetscInt     vStart, vEnd, v, vglobal, compkey, off, NVert;
978afb7921SAidan Hamilton   PetscReal    theta;
988afb7921SAidan Hamilton 
998afb7921SAidan Hamilton   PetscFunctionBegin;
1008afb7921SAidan Hamilton   PetscCall(DMGetCoordinateDM(dm, &cdm));
1018afb7921SAidan Hamilton   PetscCall(DMSetCoordinateDim(dm, 2));
1028afb7921SAidan Hamilton   PetscCall(DMNetworkGetVertexRange(cdm, &vStart, &vEnd));
1038afb7921SAidan Hamilton   PetscCall(DMNetworkRegisterComponent(cdm, "coordinates", 0, &compkey));
1048afb7921SAidan Hamilton   for (v = vStart; v < vEnd; v++) { PetscCall(DMNetworkAddComponent(cdm, v, compkey, NULL, 2)); }
1058afb7921SAidan Hamilton   PetscCall(DMNetworkFinalizeComponents(cdm));
1068afb7921SAidan Hamilton 
1078afb7921SAidan Hamilton   PetscCall(DMCreateLocalVector(cdm, &Coord));
1088afb7921SAidan Hamilton   PetscCall(VecGetArray(Coord, &coord));
1098afb7921SAidan Hamilton   PetscCall(DMNetworkGetNumVertices(cdm, NULL, &NVert));
1108afb7921SAidan Hamilton   theta = 2 * PETSC_PI / (NVert - 1);
1118afb7921SAidan Hamilton   for (v = vStart; v < vEnd; v++) {
1128afb7921SAidan Hamilton     PetscCall(DMNetworkGetGlobalVertexIndex(cdm, v, &vglobal));
1138afb7921SAidan Hamilton     PetscCall(DMNetworkGetLocalVecOffset(cdm, v, 0, &off));
1148afb7921SAidan Hamilton     if (vglobal == 0) {
1158afb7921SAidan Hamilton       coord[off]     = 0.0;
1168afb7921SAidan Hamilton       coord[off + 1] = 0.0;
1178afb7921SAidan Hamilton     } else {
1188afb7921SAidan Hamilton       /* embed on the unit circle */
1198afb7921SAidan Hamilton       coord[off]     = PetscCosReal(theta * (vglobal - 1));
1208afb7921SAidan Hamilton       coord[off + 1] = PetscSinReal(theta * (vglobal - 1));
1218afb7921SAidan Hamilton     }
1228afb7921SAidan Hamilton   }
1238afb7921SAidan Hamilton   PetscCall(VecRestoreArray(Coord, &coord));
1248afb7921SAidan Hamilton   PetscCall(DMSetCoordinatesLocal(dm, Coord));
1258afb7921SAidan Hamilton   PetscCall(VecDestroy(&Coord));
1263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1278afb7921SAidan Hamilton }
1288afb7921SAidan Hamilton 
129*cc13d412SHong Zhang /* This subroutine is used in petsc/src/snes/tutorials/network/ex1.c */
130*cc13d412SHong Zhang static PetscErrorCode CoordinatePrint(DM dm)
131*cc13d412SHong Zhang {
132*cc13d412SHong Zhang   DM                 dmclone;
133*cc13d412SHong Zhang   PetscInt           cdim, v, off, vglobal, vStart, vEnd;
134*cc13d412SHong Zhang   const PetscScalar *carray;
135*cc13d412SHong Zhang   Vec                coords;
136*cc13d412SHong Zhang   MPI_Comm           comm;
137*cc13d412SHong Zhang   PetscMPIInt        rank;
138*cc13d412SHong Zhang 
139*cc13d412SHong Zhang   PetscFunctionBegin;
140*cc13d412SHong Zhang   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
141*cc13d412SHong Zhang   PetscCallMPI(MPI_Comm_rank(comm, &rank));
142*cc13d412SHong Zhang 
143*cc13d412SHong Zhang   PetscCall(DMGetCoordinateDM(dm, &dmclone));
144*cc13d412SHong Zhang   PetscCall(DMNetworkGetVertexRange(dm, &vStart, &vEnd));
145*cc13d412SHong Zhang   PetscCall(DMGetCoordinatesLocal(dm, &coords));
146*cc13d412SHong Zhang 
147*cc13d412SHong Zhang   PetscCall(DMGetCoordinateDim(dm, &cdim));
148*cc13d412SHong Zhang   PetscCall(VecGetArrayRead(coords, &carray));
149*cc13d412SHong Zhang 
150*cc13d412SHong Zhang   PetscCall(PetscPrintf(MPI_COMM_WORLD, "\nCoordinatePrint, cdim %" PetscInt_FMT ":\n", cdim));
151*cc13d412SHong Zhang   PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "[%i]\n", rank));
152*cc13d412SHong Zhang   for (v = vStart; v < vEnd; v++) {
153*cc13d412SHong Zhang     PetscCall(DMNetworkGetLocalVecOffset(dmclone, v, 0, &off));
154*cc13d412SHong Zhang     PetscCall(DMNetworkGetGlobalVertexIndex(dmclone, v, &vglobal));
155*cc13d412SHong Zhang     switch (cdim) {
156*cc13d412SHong Zhang     case 2:
157*cc13d412SHong Zhang       PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "Vertex: %" PetscInt_FMT ", x =  %f y = %f \n", vglobal, (double)PetscRealPart(carray[off]), (double)PetscRealPart(carray[off + 1])));
158*cc13d412SHong Zhang       break;
159*cc13d412SHong Zhang     default:
160*cc13d412SHong Zhang       PetscCheck(cdim == 2, MPI_COMM_WORLD, PETSC_ERR_SUP, "Only supports Network embedding dimension of 2, not supplied  %" PetscInt_FMT, cdim);
161*cc13d412SHong Zhang       break;
162*cc13d412SHong Zhang     }
163*cc13d412SHong Zhang   }
164*cc13d412SHong Zhang   PetscCall(PetscSynchronizedFlush(MPI_COMM_WORLD, NULL));
165*cc13d412SHong Zhang   PetscCall(VecRestoreArrayRead(coords, &carray));
166*cc13d412SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
167*cc13d412SHong Zhang }
168*cc13d412SHong Zhang 
1698afb7921SAidan Hamilton int main(int argc, char **argv)
1708afb7921SAidan Hamilton {
171*cc13d412SHong Zhang   DM          dm;
172*cc13d412SHong Zhang   PetscInt    dofv = 1, dofe = 1, ne = 1;
1738afb7921SAidan Hamilton   PetscMPIInt rank;
1748afb7921SAidan Hamilton   PetscBool   testdistribute = PETSC_FALSE;
1758afb7921SAidan Hamilton 
1768afb7921SAidan Hamilton   PetscFunctionBeginUser;
1778afb7921SAidan Hamilton   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1788afb7921SAidan Hamilton   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank));
179*cc13d412SHong Zhang 
1808afb7921SAidan Hamilton   /* create a distributed k-Star graph DMNetwork */
1818afb7921SAidan Hamilton   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dofv", &dofv, NULL));
1828afb7921SAidan Hamilton   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dofe", &dofe, NULL));
1838afb7921SAidan Hamilton   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ne", &ne, NULL));
1848afb7921SAidan Hamilton   PetscCall(PetscOptionsGetBool(NULL, NULL, "-testdistribute", &testdistribute, NULL));
1858afb7921SAidan Hamilton   PetscCall(StarGraphCreate(PETSC_COMM_WORLD, dofv, dofe, ne, PETSC_TRUE, &dm));
1868afb7921SAidan Hamilton 
1878afb7921SAidan Hamilton   /* setup a quick R^2 embedding of the star graph */
1888afb7921SAidan Hamilton   PetscCall(StarGraphSetCoordinates(dm));
1898afb7921SAidan Hamilton 
1908afb7921SAidan Hamilton   if (testdistribute) {
1918afb7921SAidan Hamilton     PetscCall(DMNetworkDistribute(&dm, 0));
1928afb7921SAidan Hamilton     PetscCall(DMView(dm, PETSC_VIEWER_STDOUT_WORLD));
1938afb7921SAidan Hamilton   }
1948afb7921SAidan Hamilton 
195*cc13d412SHong Zhang   /* print or view the coordinates of each vertex */
196*cc13d412SHong Zhang   PetscCall(CoordinatePrint(dm));
1978afb7921SAidan Hamilton 
1988afb7921SAidan Hamilton   PetscCall(DMDestroy(&dm));
1998afb7921SAidan Hamilton   PetscCall(PetscFinalize());
2008afb7921SAidan Hamilton }
2018afb7921SAidan Hamilton 
2028afb7921SAidan Hamilton /*TEST
2038afb7921SAidan Hamilton 
2048afb7921SAidan Hamilton   test:
2058afb7921SAidan Hamilton     suffix: 0
2068afb7921SAidan Hamilton     args: -ne 4 -testdistribute
2078afb7921SAidan Hamilton 
2088afb7921SAidan Hamilton   test:
2098afb7921SAidan Hamilton     suffix: 1
2108afb7921SAidan Hamilton     nsize: 2
2118afb7921SAidan Hamilton     args: -ne 4 -testdistribute -petscpartitioner_type simple
2128afb7921SAidan Hamilton 
2138afb7921SAidan Hamilton  TEST*/
214