xref: /petsc/src/dm/impls/network/tests/ex2.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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   }
39*3ba16761SJacob 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(DMViewFromOptions(dm, NULL, "-dm_view"));
858afb7921SAidan Hamilton   PetscCall(DMSetUp(dm));
868afb7921SAidan Hamilton   PetscCall(PetscFree2(compedge, compvert));
878afb7921SAidan Hamilton   *newdm = dm;
88*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
898afb7921SAidan Hamilton }
908afb7921SAidan Hamilton 
918afb7921SAidan Hamilton /* Simple Circular embedding of the star graph */
928afb7921SAidan Hamilton PetscErrorCode StarGraphSetCoordinates(DM dm)
938afb7921SAidan Hamilton {
948afb7921SAidan Hamilton   DM           cdm;
958afb7921SAidan Hamilton   Vec          Coord;
968afb7921SAidan Hamilton   PetscScalar *coord;
978afb7921SAidan Hamilton   PetscInt     vStart, vEnd, v, vglobal, compkey, off, NVert;
988afb7921SAidan Hamilton   PetscReal    theta;
998afb7921SAidan Hamilton 
1008afb7921SAidan Hamilton   PetscFunctionBegin;
1018afb7921SAidan Hamilton   PetscCall(DMGetCoordinateDM(dm, &cdm));
1028afb7921SAidan Hamilton   PetscCall(DMSetCoordinateDim(dm, 2));
1038afb7921SAidan Hamilton   PetscCall(DMNetworkGetVertexRange(cdm, &vStart, &vEnd));
1048afb7921SAidan Hamilton   PetscCall(DMNetworkRegisterComponent(cdm, "coordinates", 0, &compkey));
1058afb7921SAidan Hamilton   for (v = vStart; v < vEnd; v++) { PetscCall(DMNetworkAddComponent(cdm, v, compkey, NULL, 2)); }
1068afb7921SAidan Hamilton   PetscCall(DMNetworkFinalizeComponents(cdm));
1078afb7921SAidan Hamilton 
1088afb7921SAidan Hamilton   PetscCall(DMCreateLocalVector(cdm, &Coord));
1098afb7921SAidan Hamilton   PetscCall(VecGetArray(Coord, &coord));
1108afb7921SAidan Hamilton   PetscCall(DMNetworkGetNumVertices(cdm, NULL, &NVert));
1118afb7921SAidan Hamilton   theta = 2 * PETSC_PI / (NVert - 1);
1128afb7921SAidan Hamilton   for (v = vStart; v < vEnd; v++) {
1138afb7921SAidan Hamilton     PetscCall(DMNetworkGetGlobalVertexIndex(cdm, v, &vglobal));
1148afb7921SAidan Hamilton     PetscCall(DMNetworkGetLocalVecOffset(cdm, v, 0, &off));
1158afb7921SAidan Hamilton     if (vglobal == 0) {
1168afb7921SAidan Hamilton       coord[off]     = 0.0;
1178afb7921SAidan Hamilton       coord[off + 1] = 0.0;
1188afb7921SAidan Hamilton     } else {
1198afb7921SAidan Hamilton       /* embed on the unit circle */
1208afb7921SAidan Hamilton       coord[off]     = PetscCosReal(theta * (vglobal - 1));
1218afb7921SAidan Hamilton       coord[off + 1] = PetscSinReal(theta * (vglobal - 1));
1228afb7921SAidan Hamilton     }
1238afb7921SAidan Hamilton   }
1248afb7921SAidan Hamilton   PetscCall(VecRestoreArray(Coord, &coord));
1258afb7921SAidan Hamilton   PetscCall(DMSetCoordinatesLocal(dm, Coord));
1268afb7921SAidan Hamilton   PetscCall(VecDestroy(&Coord));
127*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1288afb7921SAidan Hamilton }
1298afb7921SAidan Hamilton 
1308afb7921SAidan Hamilton int main(int argc, char **argv)
1318afb7921SAidan Hamilton {
1328afb7921SAidan Hamilton   DM           dm, cdm;
1338afb7921SAidan Hamilton   PetscInt     dofv = 1, dofe = 1, ne = 1, cdim, v, vStart, vEnd, off, vglobal;
1348afb7921SAidan Hamilton   Vec          Coord;
1358afb7921SAidan Hamilton   PetscScalar *coord;
1368afb7921SAidan Hamilton   PetscMPIInt  rank;
1378afb7921SAidan Hamilton   PetscBool    testdistribute = PETSC_FALSE;
1388afb7921SAidan Hamilton 
1398afb7921SAidan Hamilton   PetscFunctionBeginUser;
1408afb7921SAidan Hamilton   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1418afb7921SAidan Hamilton   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank));
1428afb7921SAidan Hamilton   /* create a distributed k-Star graph DMNetwork */
1438afb7921SAidan Hamilton   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dofv", &dofv, NULL));
1448afb7921SAidan Hamilton   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dofe", &dofe, NULL));
1458afb7921SAidan Hamilton   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ne", &ne, NULL));
1468afb7921SAidan Hamilton   PetscCall(PetscOptionsGetBool(NULL, NULL, "-testdistribute", &testdistribute, NULL));
1478afb7921SAidan Hamilton   PetscCall(StarGraphCreate(PETSC_COMM_WORLD, dofv, dofe, ne, PETSC_TRUE, &dm));
1488afb7921SAidan Hamilton 
1498afb7921SAidan Hamilton   /* setup a quick R^2 embedding of the star graph */
1508afb7921SAidan Hamilton   PetscCall(StarGraphSetCoordinates(dm));
1518afb7921SAidan Hamilton 
1528afb7921SAidan Hamilton   if (testdistribute) {
1538afb7921SAidan Hamilton     PetscCall(DMNetworkDistribute(&dm, 0));
1548afb7921SAidan Hamilton     PetscCall(DMView(dm, PETSC_VIEWER_STDOUT_WORLD));
1558afb7921SAidan Hamilton   }
1568afb7921SAidan Hamilton 
1578afb7921SAidan Hamilton   /* print the coordinates of each vertex */
1588afb7921SAidan Hamilton   PetscCall(DMGetCoordinateDim(dm, &cdim));
1598afb7921SAidan Hamilton   PetscCall(DMGetCoordinateDM(dm, &cdm));
1608afb7921SAidan Hamilton   PetscCall(DMNetworkGetVertexRange(cdm, &vStart, &vEnd));
1618afb7921SAidan Hamilton   PetscCall(DMGetCoordinatesLocal(dm, &Coord));
1628afb7921SAidan Hamilton   PetscCall(VecGetArray(Coord, &coord));
1638afb7921SAidan Hamilton   PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "\n Rank %i \n\n", rank));
1648afb7921SAidan Hamilton   for (v = vStart; v < vEnd; v++) {
1658afb7921SAidan Hamilton     PetscCall(DMNetworkGetLocalVecOffset(cdm, v, 0, &off));
1668afb7921SAidan Hamilton     PetscCall(DMNetworkGetGlobalVertexIndex(cdm, v, &vglobal));
1678afb7921SAidan Hamilton     switch (cdim) {
1688afb7921SAidan Hamilton     case 2:
1698afb7921SAidan Hamilton       PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "Vertex: %" PetscInt_FMT ", x =  %f y = %f \n", vglobal, (double)PetscRealPart(coord[off]), (double)PetscRealPart(coord[off + 1])));
1708afb7921SAidan Hamilton       break;
1718afb7921SAidan Hamilton     default:
1728afb7921SAidan Hamilton       SETERRQ(MPI_COMM_WORLD, PETSC_ERR_SUP, "Only supports Network embedding dimension of 2, not supplied  %" PetscInt_FMT, cdim);
1738afb7921SAidan Hamilton       break;
1748afb7921SAidan Hamilton     }
1758afb7921SAidan Hamilton   }
1768afb7921SAidan Hamilton   PetscCall(PetscSynchronizedFlush(MPI_COMM_WORLD, NULL));
1778afb7921SAidan Hamilton 
1788afb7921SAidan Hamilton   PetscCall(DMDestroy(&dm));
1798afb7921SAidan Hamilton   PetscCall(PetscFinalize());
1808afb7921SAidan Hamilton }
1818afb7921SAidan Hamilton 
1828afb7921SAidan Hamilton /*TEST
1838afb7921SAidan Hamilton 
1848afb7921SAidan Hamilton   test:
1858afb7921SAidan Hamilton     suffix: 0
1868afb7921SAidan Hamilton     args: -ne 4 -testdistribute
1878afb7921SAidan Hamilton 
1888afb7921SAidan Hamilton   test:
1898afb7921SAidan Hamilton     suffix: 1
1908afb7921SAidan Hamilton     nsize: 2
1918afb7921SAidan Hamilton     args: -ne 4 -testdistribute -petscpartitioner_type simple
1928afb7921SAidan Hamilton 
1938afb7921SAidan Hamilton  TEST*/
194