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 */ 20*7cd471e7SHong Zhang static 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 /* Simple Circular embedding of the star graph */ 43*7cd471e7SHong Zhang static PetscErrorCode StarGraphSetCoordinates(DM dm, PetscReal *vcolor) 448afb7921SAidan Hamilton { 458afb7921SAidan Hamilton DM cdm; 468afb7921SAidan Hamilton Vec Coord; 478afb7921SAidan Hamilton PetscScalar *coord; 48*7cd471e7SHong Zhang PetscInt vStart, vEnd, v, vglobal, compkey = 0, off, NVert; 498afb7921SAidan Hamilton PetscReal theta; 508afb7921SAidan Hamilton 518afb7921SAidan Hamilton PetscFunctionBegin; 528afb7921SAidan Hamilton PetscCall(DMSetCoordinateDim(dm, 2)); 53*7cd471e7SHong Zhang PetscCall(DMGetCoordinateDM(dm, &cdm)); 54*7cd471e7SHong Zhang 558afb7921SAidan Hamilton PetscCall(DMNetworkGetVertexRange(cdm, &vStart, &vEnd)); 56*7cd471e7SHong Zhang PetscCall(DMNetworkRegisterComponent(cdm, "coordinate", sizeof(PetscReal), &compkey)); 57*7cd471e7SHong Zhang for (v = vStart; v < vEnd; v++) { 58*7cd471e7SHong Zhang PetscCall(DMNetworkGetGlobalVertexIndex(cdm, v, &vglobal)); 59*7cd471e7SHong Zhang vcolor[v - vStart] = vglobal; 60*7cd471e7SHong Zhang PetscCall(DMNetworkAddComponent(cdm, v, compkey, &vcolor[v - vStart], 2)); 61*7cd471e7SHong Zhang } 628afb7921SAidan Hamilton PetscCall(DMNetworkFinalizeComponents(cdm)); 638afb7921SAidan Hamilton 648afb7921SAidan Hamilton PetscCall(DMCreateLocalVector(cdm, &Coord)); 658afb7921SAidan Hamilton PetscCall(VecGetArray(Coord, &coord)); 668afb7921SAidan Hamilton PetscCall(DMNetworkGetNumVertices(cdm, NULL, &NVert)); 678afb7921SAidan Hamilton theta = 2 * PETSC_PI / (NVert - 1); 688afb7921SAidan Hamilton for (v = vStart; v < vEnd; v++) { 698afb7921SAidan Hamilton PetscCall(DMNetworkGetGlobalVertexIndex(cdm, v, &vglobal)); 708afb7921SAidan Hamilton PetscCall(DMNetworkGetLocalVecOffset(cdm, v, 0, &off)); 718afb7921SAidan Hamilton if (vglobal == 0) { 728afb7921SAidan Hamilton coord[off] = 0.0; 738afb7921SAidan Hamilton coord[off + 1] = 0.0; 748afb7921SAidan Hamilton } else { 758afb7921SAidan Hamilton /* embed on the unit circle */ 768afb7921SAidan Hamilton coord[off] = PetscCosReal(theta * (vglobal - 1)); 778afb7921SAidan Hamilton coord[off + 1] = PetscSinReal(theta * (vglobal - 1)); 788afb7921SAidan Hamilton } 798afb7921SAidan Hamilton } 808afb7921SAidan Hamilton PetscCall(VecRestoreArray(Coord, &coord)); 818afb7921SAidan Hamilton PetscCall(DMSetCoordinatesLocal(dm, Coord)); 828afb7921SAidan Hamilton PetscCall(VecDestroy(&Coord)); 833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 848afb7921SAidan Hamilton } 858afb7921SAidan Hamilton 86*7cd471e7SHong Zhang /* 87*7cd471e7SHong Zhang CreateSimpleStarGraph - Create a Distributed k-Star Graph DMNetwork with a single PetscInt component on 88*7cd471e7SHong Zhang all edges and vertices, a selectable number of dofs on vertices and edges. Intended mostly to be used for testing purposes. 89*7cd471e7SHong Zhang 90*7cd471e7SHong Zhang Input Parameters: 91*7cd471e7SHong Zhang . comm - the communicator of the dm 92*7cd471e7SHong Zhang . numdofvert - number of degrees of freedom (dofs) on vertices 93*7cd471e7SHong Zhang . numdofedge - number of degrees of freedom (dofs) on edges 94*7cd471e7SHong Zhang . k - order of the star graph (number of edges) 95*7cd471e7SHong Zhang . directin - if true direction of edges is towards the center vertex, otherwise they are directed out of the center vertex 96*7cd471e7SHong Zhang 97*7cd471e7SHong Zhang Output Parameters: 98*7cd471e7SHong Zhang . newdm - The created and distributed simple Star Graph 99*7cd471e7SHong Zhang */ 100*7cd471e7SHong Zhang static PetscErrorCode StarGraphCreate(MPI_Comm comm, PetscInt numdofvert, PetscInt numdofedge, PetscInt k, PetscBool directin, DM *newdm) 101*7cd471e7SHong Zhang { 102*7cd471e7SHong Zhang DM dm; 103*7cd471e7SHong Zhang PetscMPIInt rank; 104*7cd471e7SHong Zhang PetscInt ne = 0, compkey, eStart, eEnd, vStart, vEnd, e, v; 105*7cd471e7SHong Zhang PetscInt *edgelist = NULL, *compedge, *compvert; 106*7cd471e7SHong Zhang PetscReal *vcolor; 107*7cd471e7SHong Zhang 108*7cd471e7SHong Zhang PetscFunctionBegin; 109*7cd471e7SHong Zhang PetscCall(DMNetworkCreate(comm, &dm)); 110*7cd471e7SHong Zhang PetscCall(DMNetworkSetNumSubNetworks(dm, PETSC_DECIDE, 1)); 111*7cd471e7SHong Zhang PetscCallMPI(MPI_Comm_rank(comm, &rank)); 112*7cd471e7SHong Zhang if (rank == 0) PetscCall(StarGraphCreateEdgeList(k, directin, &ne, &edgelist)); 113*7cd471e7SHong Zhang PetscCall(DMNetworkAddSubnetwork(dm, "Main", ne, edgelist, NULL)); 114*7cd471e7SHong Zhang PetscCall(DMNetworkRegisterComponent(dm, "dummy", sizeof(PetscInt), &compkey)); 115*7cd471e7SHong Zhang PetscCall(DMNetworkLayoutSetUp(dm)); 116*7cd471e7SHong Zhang PetscCall(PetscFree(edgelist)); 117*7cd471e7SHong Zhang PetscCall(DMNetworkGetEdgeRange(dm, &eStart, &eEnd)); 118*7cd471e7SHong Zhang PetscCall(DMNetworkGetVertexRange(dm, &vStart, &vEnd)); 119*7cd471e7SHong Zhang PetscCall(PetscCalloc3(eEnd - eStart, &compedge, vEnd - vStart, &compvert, vEnd - vStart, &vcolor)); 120*7cd471e7SHong Zhang for (e = eStart; e < eEnd; e++) { 121*7cd471e7SHong Zhang compedge[e - eStart] = e; 122*7cd471e7SHong Zhang PetscCall(DMNetworkAddComponent(dm, e, compkey, &compedge[e - eStart], numdofedge)); 123*7cd471e7SHong Zhang } 124*7cd471e7SHong Zhang for (v = vStart; v < vEnd; v++) { 125*7cd471e7SHong Zhang compvert[v - vStart] = v; 126*7cd471e7SHong Zhang PetscCall(DMNetworkAddComponent(dm, v, compkey, &compvert[v - vStart], numdofvert)); 127*7cd471e7SHong Zhang } 128*7cd471e7SHong Zhang 129*7cd471e7SHong Zhang PetscCall(StarGraphSetCoordinates(dm, vcolor)); 130*7cd471e7SHong Zhang 131*7cd471e7SHong Zhang PetscCall(DMSetFromOptions(dm)); 132*7cd471e7SHong Zhang PetscCall(DMSetUp(dm)); 133*7cd471e7SHong Zhang PetscCall(PetscFree3(compedge, compvert, vcolor)); 134*7cd471e7SHong Zhang 135*7cd471e7SHong Zhang *newdm = dm; 136*7cd471e7SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 137*7cd471e7SHong Zhang } 138*7cd471e7SHong Zhang 139cc13d412SHong Zhang /* This subroutine is used in petsc/src/snes/tutorials/network/ex1.c */ 140cc13d412SHong Zhang static PetscErrorCode CoordinatePrint(DM dm) 141cc13d412SHong Zhang { 142cc13d412SHong Zhang DM dmclone; 143cc13d412SHong Zhang PetscInt cdim, v, off, vglobal, vStart, vEnd; 144cc13d412SHong Zhang const PetscScalar *carray; 145cc13d412SHong Zhang Vec coords; 146cc13d412SHong Zhang MPI_Comm comm; 147cc13d412SHong Zhang PetscMPIInt rank; 148cc13d412SHong Zhang 149cc13d412SHong Zhang PetscFunctionBegin; 150cc13d412SHong Zhang PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 151cc13d412SHong Zhang PetscCallMPI(MPI_Comm_rank(comm, &rank)); 152cc13d412SHong Zhang 153cc13d412SHong Zhang PetscCall(DMGetCoordinateDM(dm, &dmclone)); 154cc13d412SHong Zhang PetscCall(DMNetworkGetVertexRange(dm, &vStart, &vEnd)); 155cc13d412SHong Zhang PetscCall(DMGetCoordinatesLocal(dm, &coords)); 156cc13d412SHong Zhang 157cc13d412SHong Zhang PetscCall(DMGetCoordinateDim(dm, &cdim)); 158cc13d412SHong Zhang PetscCall(VecGetArrayRead(coords, &carray)); 159cc13d412SHong Zhang 160cc13d412SHong Zhang PetscCall(PetscPrintf(MPI_COMM_WORLD, "\nCoordinatePrint, cdim %" PetscInt_FMT ":\n", cdim)); 161*7cd471e7SHong Zhang PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "[%d]\n", rank)); 162cc13d412SHong Zhang for (v = vStart; v < vEnd; v++) { 163cc13d412SHong Zhang PetscCall(DMNetworkGetLocalVecOffset(dmclone, v, 0, &off)); 164cc13d412SHong Zhang PetscCall(DMNetworkGetGlobalVertexIndex(dmclone, v, &vglobal)); 165cc13d412SHong Zhang switch (cdim) { 166cc13d412SHong Zhang case 2: 167cc13d412SHong Zhang PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "Vertex: %" PetscInt_FMT ", x = %f y = %f \n", vglobal, (double)PetscRealPart(carray[off]), (double)PetscRealPart(carray[off + 1]))); 168cc13d412SHong Zhang break; 169cc13d412SHong Zhang default: 170cc13d412SHong Zhang PetscCheck(cdim == 2, MPI_COMM_WORLD, PETSC_ERR_SUP, "Only supports Network embedding dimension of 2, not supplied %" PetscInt_FMT, cdim); 171cc13d412SHong Zhang break; 172cc13d412SHong Zhang } 173cc13d412SHong Zhang } 174cc13d412SHong Zhang PetscCall(PetscSynchronizedFlush(MPI_COMM_WORLD, NULL)); 175cc13d412SHong Zhang PetscCall(VecRestoreArrayRead(coords, &carray)); 176cc13d412SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 177cc13d412SHong Zhang } 178cc13d412SHong Zhang 1798afb7921SAidan Hamilton int main(int argc, char **argv) 1808afb7921SAidan Hamilton { 181cc13d412SHong Zhang DM dm; 182cc13d412SHong Zhang PetscInt dofv = 1, dofe = 1, ne = 1; 1838afb7921SAidan Hamilton PetscMPIInt rank; 184*7cd471e7SHong Zhang PetscBool viewCSV = PETSC_FALSE; 1858afb7921SAidan Hamilton 1868afb7921SAidan Hamilton PetscFunctionBeginUser; 1878afb7921SAidan Hamilton PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1888afb7921SAidan Hamilton PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank)); 189cc13d412SHong Zhang 1908afb7921SAidan Hamilton /* create a distributed k-Star graph DMNetwork */ 1918afb7921SAidan Hamilton PetscCall(PetscOptionsGetInt(NULL, NULL, "-dofv", &dofv, NULL)); 1928afb7921SAidan Hamilton PetscCall(PetscOptionsGetInt(NULL, NULL, "-dofe", &dofe, NULL)); 1938afb7921SAidan Hamilton PetscCall(PetscOptionsGetInt(NULL, NULL, "-ne", &ne, NULL)); 194*7cd471e7SHong Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewCSV", &viewCSV, NULL)); 195*7cd471e7SHong Zhang 196*7cd471e7SHong Zhang /* create and setup a quick R^2 embedding of the star graph */ 1978afb7921SAidan Hamilton PetscCall(StarGraphCreate(PETSC_COMM_WORLD, dofv, dofe, ne, PETSC_TRUE, &dm)); 1988afb7921SAidan Hamilton 1998afb7921SAidan Hamilton PetscCall(DMNetworkDistribute(&dm, 0)); 200*7cd471e7SHong Zhang if (viewCSV) { /* CSV View of network with coordinates */ 201*7cd471e7SHong Zhang PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_CSV)); 202*7cd471e7SHong Zhang PetscCall(DMView(dm, PETSC_VIEWER_STDOUT_WORLD)); 203*7cd471e7SHong Zhang PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 204*7cd471e7SHong Zhang } else { 2058afb7921SAidan Hamilton PetscCall(DMView(dm, PETSC_VIEWER_STDOUT_WORLD)); 2068afb7921SAidan Hamilton } 2078afb7921SAidan Hamilton 208cc13d412SHong Zhang /* print or view the coordinates of each vertex */ 209cc13d412SHong Zhang PetscCall(CoordinatePrint(dm)); 2108afb7921SAidan Hamilton 2118afb7921SAidan Hamilton PetscCall(DMDestroy(&dm)); 2128afb7921SAidan Hamilton PetscCall(PetscFinalize()); 2138afb7921SAidan Hamilton } 2148afb7921SAidan Hamilton 2158afb7921SAidan Hamilton /*TEST 2168afb7921SAidan Hamilton 2178afb7921SAidan Hamilton test: 2188afb7921SAidan Hamilton suffix: 0 219*7cd471e7SHong Zhang args: -ne 4 2208afb7921SAidan Hamilton 2218afb7921SAidan Hamilton test: 2228afb7921SAidan Hamilton suffix: 1 2238afb7921SAidan Hamilton nsize: 2 224*7cd471e7SHong Zhang args: -ne 4 -petscpartitioner_type simple 2258afb7921SAidan Hamilton 2268afb7921SAidan Hamilton TEST*/ 227