1*8afb7921SAidan Hamilton static char help[] = "Test DMCreateCoordinateDM_Network, and related functions \n\n"; 2*8afb7921SAidan Hamilton 3*8afb7921SAidan Hamilton #include <petscdmnetwork.h> 4*8afb7921SAidan Hamilton 5*8afb7921SAidan Hamilton /* 6*8afb7921SAidan Hamilton CreateStarGraphEdgeList - Create a k-Star Graph Edgelist on current processor 7*8afb7921SAidan Hamilton Not Collective 8*8afb7921SAidan Hamilton 9*8afb7921SAidan Hamilton Input Parameters: 10*8afb7921SAidan Hamilton . k - order of the star graph (number of edges) 11*8afb7921SAidan Hamilton . directin - if true direction of edges is towards the center vertex, otherwise they are directed out of the center vertex. 12*8afb7921SAidan Hamilton 13*8afb7921SAidan Hamilton Output Parameters: 14*8afb7921SAidan Hamilton . ne - number of edges of this star graph 15*8afb7921SAidan 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*8afb7921SAidan Hamilton [first vertex of first edge, second vertex of first edge, first vertex of second edge, second vertex of second edge, etc]. 17*8afb7921SAidan Hamilton 18*8afb7921SAidan Hamilton User is responsible for deallocating this memory. 19*8afb7921SAidan Hamilton */ 20*8afb7921SAidan Hamilton PetscErrorCode StarGraphCreateEdgeList(PetscInt k, PetscBool directin, PetscInt *ne, PetscInt *edgelist[]) 21*8afb7921SAidan Hamilton { 22*8afb7921SAidan Hamilton PetscInt i; 23*8afb7921SAidan Hamilton 24*8afb7921SAidan Hamilton PetscFunctionBegin; 25*8afb7921SAidan Hamilton *ne = k; 26*8afb7921SAidan Hamilton PetscCall(PetscCalloc1(2 * k, edgelist)); 27*8afb7921SAidan Hamilton 28*8afb7921SAidan Hamilton if (directin) { 29*8afb7921SAidan Hamilton for (i = 0; i < k; i++) { 30*8afb7921SAidan Hamilton (*edgelist)[2 * i] = i + 1; 31*8afb7921SAidan Hamilton (*edgelist)[2 * i + 1] = 0; 32*8afb7921SAidan Hamilton } 33*8afb7921SAidan Hamilton } else { 34*8afb7921SAidan Hamilton for (i = 0; i < k; i++) { 35*8afb7921SAidan Hamilton (*edgelist)[2 * i] = 0; 36*8afb7921SAidan Hamilton (*edgelist)[2 * i + 1] = i + 1; 37*8afb7921SAidan Hamilton } 38*8afb7921SAidan Hamilton } 39*8afb7921SAidan Hamilton PetscFunctionReturn(0); 40*8afb7921SAidan Hamilton } 41*8afb7921SAidan Hamilton 42*8afb7921SAidan Hamilton /* 43*8afb7921SAidan Hamilton CreateSimpleStarGraph - Create a Distributed k-Star Graph DMNetwork with a single PetscInt component on 44*8afb7921SAidan Hamilton all edges and vertices, a selectable number of dofs on vertices and edges. Intended mostly to be used for testing purposes. 45*8afb7921SAidan Hamilton 46*8afb7921SAidan Hamilton Input Parameters: 47*8afb7921SAidan Hamilton . comm - the communicator of the dm 48*8afb7921SAidan Hamilton . numdofvert - number of degrees of freedom (dofs) on vertices 49*8afb7921SAidan Hamilton . numdofedge - number of degrees of freedom (dofs) on edges 50*8afb7921SAidan Hamilton . k - order of the star graph (number of edges) 51*8afb7921SAidan Hamilton . directin - if true direction of edges is towards the center vertex, otherwise they are directed out of the center vertex 52*8afb7921SAidan Hamilton 53*8afb7921SAidan Hamilton Output Parameters: 54*8afb7921SAidan Hamilton . newdm - The created and distributed simple Star Graph 55*8afb7921SAidan Hamilton */ 56*8afb7921SAidan Hamilton PetscErrorCode StarGraphCreate(MPI_Comm comm, PetscInt numdofvert, PetscInt numdofedge, PetscInt k, PetscBool directin, DM *newdm) 57*8afb7921SAidan Hamilton { 58*8afb7921SAidan Hamilton DM dm; 59*8afb7921SAidan Hamilton PetscMPIInt rank; 60*8afb7921SAidan Hamilton PetscInt ne = 0, compkey, eStart, eEnd, vStart, vEnd, e, v; 61*8afb7921SAidan Hamilton PetscInt *edgelist = NULL, *compedge, *compvert; 62*8afb7921SAidan Hamilton 63*8afb7921SAidan Hamilton PetscFunctionBegin; 64*8afb7921SAidan Hamilton PetscCall(DMNetworkCreate(comm, &dm)); 65*8afb7921SAidan Hamilton PetscCall(DMNetworkSetNumSubNetworks(dm, PETSC_DECIDE, 1)); 66*8afb7921SAidan Hamilton PetscCallMPI(MPI_Comm_rank(comm, &rank)); 67*8afb7921SAidan Hamilton if (rank == 0) PetscCall(StarGraphCreateEdgeList(k, directin, &ne, &edgelist)); 68*8afb7921SAidan Hamilton PetscCall(DMNetworkAddSubnetwork(dm, "Main", ne, edgelist, NULL)); 69*8afb7921SAidan Hamilton PetscCall(DMNetworkRegisterComponent(dm, "dummy", sizeof(PetscInt), &compkey)); 70*8afb7921SAidan Hamilton PetscCall(DMNetworkLayoutSetUp(dm)); 71*8afb7921SAidan Hamilton PetscCall(PetscFree(edgelist)); 72*8afb7921SAidan Hamilton PetscCall(DMNetworkGetEdgeRange(dm, &eStart, &eEnd)); 73*8afb7921SAidan Hamilton PetscCall(DMNetworkGetVertexRange(dm, &vStart, &vEnd)); 74*8afb7921SAidan Hamilton PetscCall(PetscMalloc2(eEnd - eStart, &compedge, vEnd - vStart, &compvert)); 75*8afb7921SAidan Hamilton for (e = eStart; e < eEnd; e++) { 76*8afb7921SAidan Hamilton compedge[e - eStart] = e; 77*8afb7921SAidan Hamilton PetscCall(DMNetworkAddComponent(dm, e, compkey, &compedge[e - eStart], numdofedge)); 78*8afb7921SAidan Hamilton } 79*8afb7921SAidan Hamilton for (v = vStart; v < vEnd; v++) { 80*8afb7921SAidan Hamilton compvert[v - vStart] = v; 81*8afb7921SAidan Hamilton PetscCall(DMNetworkAddComponent(dm, v, compkey, &compvert[v - vStart], numdofvert)); 82*8afb7921SAidan Hamilton } 83*8afb7921SAidan Hamilton PetscCall(DMSetFromOptions(dm)); 84*8afb7921SAidan Hamilton PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 85*8afb7921SAidan Hamilton PetscCall(DMSetUp(dm)); 86*8afb7921SAidan Hamilton PetscCall(PetscFree2(compedge, compvert)); 87*8afb7921SAidan Hamilton *newdm = dm; 88*8afb7921SAidan Hamilton PetscFunctionReturn(0); 89*8afb7921SAidan Hamilton } 90*8afb7921SAidan Hamilton 91*8afb7921SAidan Hamilton /* Simple Circular embedding of the star graph */ 92*8afb7921SAidan Hamilton PetscErrorCode StarGraphSetCoordinates(DM dm) 93*8afb7921SAidan Hamilton { 94*8afb7921SAidan Hamilton DM cdm; 95*8afb7921SAidan Hamilton Vec Coord; 96*8afb7921SAidan Hamilton PetscScalar *coord; 97*8afb7921SAidan Hamilton PetscInt vStart, vEnd, v, vglobal, compkey, off, NVert; 98*8afb7921SAidan Hamilton PetscReal theta; 99*8afb7921SAidan Hamilton 100*8afb7921SAidan Hamilton PetscFunctionBegin; 101*8afb7921SAidan Hamilton PetscCall(DMGetCoordinateDM(dm, &cdm)); 102*8afb7921SAidan Hamilton PetscCall(DMSetCoordinateDim(dm, 2)); 103*8afb7921SAidan Hamilton PetscCall(DMNetworkGetVertexRange(cdm, &vStart, &vEnd)); 104*8afb7921SAidan Hamilton PetscCall(DMNetworkRegisterComponent(cdm, "coordinates", 0, &compkey)); 105*8afb7921SAidan Hamilton for (v = vStart; v < vEnd; v++) { PetscCall(DMNetworkAddComponent(cdm, v, compkey, NULL, 2)); } 106*8afb7921SAidan Hamilton PetscCall(DMNetworkFinalizeComponents(cdm)); 107*8afb7921SAidan Hamilton 108*8afb7921SAidan Hamilton PetscCall(DMCreateLocalVector(cdm, &Coord)); 109*8afb7921SAidan Hamilton PetscCall(VecGetArray(Coord, &coord)); 110*8afb7921SAidan Hamilton PetscCall(DMNetworkGetNumVertices(cdm, NULL, &NVert)); 111*8afb7921SAidan Hamilton theta = 2 * PETSC_PI / (NVert - 1); 112*8afb7921SAidan Hamilton for (v = vStart; v < vEnd; v++) { 113*8afb7921SAidan Hamilton PetscCall(DMNetworkGetGlobalVertexIndex(cdm, v, &vglobal)); 114*8afb7921SAidan Hamilton PetscCall(DMNetworkGetLocalVecOffset(cdm, v, 0, &off)); 115*8afb7921SAidan Hamilton if (vglobal == 0) { 116*8afb7921SAidan Hamilton coord[off] = 0.0; 117*8afb7921SAidan Hamilton coord[off + 1] = 0.0; 118*8afb7921SAidan Hamilton } else { 119*8afb7921SAidan Hamilton /* embed on the unit circle */ 120*8afb7921SAidan Hamilton coord[off] = PetscCosReal(theta * (vglobal - 1)); 121*8afb7921SAidan Hamilton coord[off + 1] = PetscSinReal(theta * (vglobal - 1)); 122*8afb7921SAidan Hamilton } 123*8afb7921SAidan Hamilton } 124*8afb7921SAidan Hamilton PetscCall(VecRestoreArray(Coord, &coord)); 125*8afb7921SAidan Hamilton PetscCall(DMSetCoordinatesLocal(dm, Coord)); 126*8afb7921SAidan Hamilton PetscCall(VecDestroy(&Coord)); 127*8afb7921SAidan Hamilton PetscFunctionReturn(0); 128*8afb7921SAidan Hamilton } 129*8afb7921SAidan Hamilton 130*8afb7921SAidan Hamilton int main(int argc, char **argv) 131*8afb7921SAidan Hamilton { 132*8afb7921SAidan Hamilton DM dm, cdm; 133*8afb7921SAidan Hamilton PetscInt dofv = 1, dofe = 1, ne = 1, cdim, v, vStart, vEnd, off, vglobal; 134*8afb7921SAidan Hamilton Vec Coord; 135*8afb7921SAidan Hamilton PetscScalar *coord; 136*8afb7921SAidan Hamilton PetscMPIInt rank; 137*8afb7921SAidan Hamilton PetscBool testdistribute = PETSC_FALSE; 138*8afb7921SAidan Hamilton 139*8afb7921SAidan Hamilton PetscFunctionBeginUser; 140*8afb7921SAidan Hamilton PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 141*8afb7921SAidan Hamilton PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank)); 142*8afb7921SAidan Hamilton /* create a distributed k-Star graph DMNetwork */ 143*8afb7921SAidan Hamilton PetscCall(PetscOptionsGetInt(NULL, NULL, "-dofv", &dofv, NULL)); 144*8afb7921SAidan Hamilton PetscCall(PetscOptionsGetInt(NULL, NULL, "-dofe", &dofe, NULL)); 145*8afb7921SAidan Hamilton PetscCall(PetscOptionsGetInt(NULL, NULL, "-ne", &ne, NULL)); 146*8afb7921SAidan Hamilton PetscCall(PetscOptionsGetBool(NULL, NULL, "-testdistribute", &testdistribute, NULL)); 147*8afb7921SAidan Hamilton PetscCall(StarGraphCreate(PETSC_COMM_WORLD, dofv, dofe, ne, PETSC_TRUE, &dm)); 148*8afb7921SAidan Hamilton 149*8afb7921SAidan Hamilton /* setup a quick R^2 embedding of the star graph */ 150*8afb7921SAidan Hamilton PetscCall(StarGraphSetCoordinates(dm)); 151*8afb7921SAidan Hamilton 152*8afb7921SAidan Hamilton if (testdistribute) { 153*8afb7921SAidan Hamilton PetscCall(DMNetworkDistribute(&dm, 0)); 154*8afb7921SAidan Hamilton PetscCall(DMView(dm, PETSC_VIEWER_STDOUT_WORLD)); 155*8afb7921SAidan Hamilton } 156*8afb7921SAidan Hamilton 157*8afb7921SAidan Hamilton /* print the coordinates of each vertex */ 158*8afb7921SAidan Hamilton PetscCall(DMGetCoordinateDim(dm, &cdim)); 159*8afb7921SAidan Hamilton PetscCall(DMGetCoordinateDM(dm, &cdm)); 160*8afb7921SAidan Hamilton PetscCall(DMNetworkGetVertexRange(cdm, &vStart, &vEnd)); 161*8afb7921SAidan Hamilton PetscCall(DMGetCoordinatesLocal(dm, &Coord)); 162*8afb7921SAidan Hamilton PetscCall(VecGetArray(Coord, &coord)); 163*8afb7921SAidan Hamilton PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "\n Rank %i \n\n", rank)); 164*8afb7921SAidan Hamilton for (v = vStart; v < vEnd; v++) { 165*8afb7921SAidan Hamilton PetscCall(DMNetworkGetLocalVecOffset(cdm, v, 0, &off)); 166*8afb7921SAidan Hamilton PetscCall(DMNetworkGetGlobalVertexIndex(cdm, v, &vglobal)); 167*8afb7921SAidan Hamilton switch (cdim) { 168*8afb7921SAidan Hamilton case 2: 169*8afb7921SAidan Hamilton PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "Vertex: %" PetscInt_FMT ", x = %f y = %f \n", vglobal, (double)PetscRealPart(coord[off]), (double)PetscRealPart(coord[off + 1]))); 170*8afb7921SAidan Hamilton break; 171*8afb7921SAidan Hamilton default: 172*8afb7921SAidan Hamilton SETERRQ(MPI_COMM_WORLD, PETSC_ERR_SUP, "Only supports Network embedding dimension of 2, not supplied %" PetscInt_FMT, cdim); 173*8afb7921SAidan Hamilton break; 174*8afb7921SAidan Hamilton } 175*8afb7921SAidan Hamilton } 176*8afb7921SAidan Hamilton PetscCall(PetscSynchronizedFlush(MPI_COMM_WORLD, NULL)); 177*8afb7921SAidan Hamilton 178*8afb7921SAidan Hamilton PetscCall(DMDestroy(&dm)); 179*8afb7921SAidan Hamilton PetscCall(PetscFinalize()); 180*8afb7921SAidan Hamilton } 181*8afb7921SAidan Hamilton 182*8afb7921SAidan Hamilton /*TEST 183*8afb7921SAidan Hamilton 184*8afb7921SAidan Hamilton test: 185*8afb7921SAidan Hamilton suffix: 0 186*8afb7921SAidan Hamilton args: -ne 4 -testdistribute 187*8afb7921SAidan Hamilton 188*8afb7921SAidan Hamilton test: 189*8afb7921SAidan Hamilton suffix: 1 190*8afb7921SAidan Hamilton nsize: 2 191*8afb7921SAidan Hamilton args: -ne 4 -testdistribute -petscpartitioner_type simple 192*8afb7921SAidan Hamilton 193*8afb7921SAidan Hamilton TEST*/ 194