1*c73416fdSAidan Hamilton static char help[] = "Test query functions for DMNetwork \n\n"; 2*c73416fdSAidan Hamilton 3*c73416fdSAidan Hamilton #include <petscdmnetwork.h> 4*c73416fdSAidan Hamilton 5*c73416fdSAidan Hamilton /* 6*c73416fdSAidan Hamilton CreateStarGraphEdgeList - Create a k-Star Graph Edgelist on current processor 7*c73416fdSAidan Hamilton Not Collective 8*c73416fdSAidan Hamilton 9*c73416fdSAidan Hamilton Input Parameters: 10*c73416fdSAidan Hamilton . k - order of the star graph (number of edges) 11*c73416fdSAidan Hamilton . directin - if true direction of edges is towards the center vertex, otherwise they are directed out of the center vertex. 12*c73416fdSAidan Hamilton 13*c73416fdSAidan Hamilton Output Parameters: 14*c73416fdSAidan Hamilton . ne - number of edges of this star graph 15*c73416fdSAidan 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*c73416fdSAidan Hamilton [first vertex of first edge, second vertex of first edge, first vertex of second edge, second vertex of second edge, etc]. 17*c73416fdSAidan Hamilton 18*c73416fdSAidan Hamilton User is responsible for deallocating this memory. 19*c73416fdSAidan Hamilton */ 20*c73416fdSAidan Hamilton PetscErrorCode StarGraphCreateEdgeList(PetscInt k, PetscBool directin, PetscInt *ne, PetscInt *edgelist[]) 21*c73416fdSAidan Hamilton { 22*c73416fdSAidan Hamilton PetscInt i; 23*c73416fdSAidan Hamilton 24*c73416fdSAidan Hamilton PetscFunctionBegin; 25*c73416fdSAidan Hamilton *ne = k; 26*c73416fdSAidan Hamilton PetscCall(PetscCalloc1(2 * k, edgelist)); 27*c73416fdSAidan Hamilton 28*c73416fdSAidan Hamilton if (directin) { 29*c73416fdSAidan Hamilton for (i = 0; i < k; i++) { 30*c73416fdSAidan Hamilton (*edgelist)[2 * i] = i + 1; 31*c73416fdSAidan Hamilton (*edgelist)[2 * i + 1] = 0; 32*c73416fdSAidan Hamilton } 33*c73416fdSAidan Hamilton } else { 34*c73416fdSAidan Hamilton for (i = 0; i < k; i++) { 35*c73416fdSAidan Hamilton (*edgelist)[2 * i] = 0; 36*c73416fdSAidan Hamilton (*edgelist)[2 * i + 1] = i + 1; 37*c73416fdSAidan Hamilton } 38*c73416fdSAidan Hamilton } 39*c73416fdSAidan Hamilton PetscFunctionReturn(0); 40*c73416fdSAidan Hamilton } 41*c73416fdSAidan Hamilton 42*c73416fdSAidan Hamilton /* 43*c73416fdSAidan Hamilton CreateSimpleStarGraph - Create a Distributed k-Star Graph DMNetwork with a single PetscInt component on 44*c73416fdSAidan Hamilton all edges and vertices, a selectable number of dofs on vertices and edges. Intended mostly to be used for testing purposes. 45*c73416fdSAidan Hamilton 46*c73416fdSAidan Hamilton Input Parameters: 47*c73416fdSAidan Hamilton . comm - the communicator of the dm 48*c73416fdSAidan Hamilton . numdofvert - number of degrees of freedom (dofs) on vertices 49*c73416fdSAidan Hamilton . numdofedge - number of degrees of freedom (dofs) on edges 50*c73416fdSAidan Hamilton . k - order of the star graph (number of edges) 51*c73416fdSAidan Hamilton . directin - if true direction of edges is towards the center vertex, otherwise they are directed out of the center vertex 52*c73416fdSAidan Hamilton 53*c73416fdSAidan Hamilton Output Parameters: 54*c73416fdSAidan Hamilton . newdm - The created and distributed simple Star Graph 55*c73416fdSAidan Hamilton */ 56*c73416fdSAidan Hamilton PetscErrorCode StarGraphCreate(MPI_Comm comm, PetscInt numdofvert, PetscInt numdofedge, PetscInt k, PetscBool directin, DM *newdm) 57*c73416fdSAidan Hamilton { 58*c73416fdSAidan Hamilton DM dm; 59*c73416fdSAidan Hamilton PetscMPIInt rank; 60*c73416fdSAidan Hamilton PetscInt ne = 0, compkey, eStart, eEnd, vStart, vEnd, e, v; 61*c73416fdSAidan Hamilton PetscInt *edgelist = NULL, *compedge, *compvert; 62*c73416fdSAidan Hamilton 63*c73416fdSAidan Hamilton PetscFunctionBegin; 64*c73416fdSAidan Hamilton PetscCall(DMNetworkCreate(comm, &dm)); 65*c73416fdSAidan Hamilton PetscCall(DMNetworkSetNumSubNetworks(dm, PETSC_DECIDE, 1)); 66*c73416fdSAidan Hamilton PetscCallMPI(MPI_Comm_rank(comm, &rank)); 67*c73416fdSAidan Hamilton if (rank == 0) PetscCall(StarGraphCreateEdgeList(k, directin, &ne, &edgelist)); 68*c73416fdSAidan Hamilton PetscCall(DMNetworkAddSubnetwork(dm, "Main", ne, edgelist, NULL)); 69*c73416fdSAidan Hamilton PetscCall(DMNetworkRegisterComponent(dm, "dummy", sizeof(PetscInt), &compkey)); 70*c73416fdSAidan Hamilton PetscCall(DMNetworkLayoutSetUp(dm)); 71*c73416fdSAidan Hamilton PetscCall(PetscFree(edgelist)); 72*c73416fdSAidan Hamilton PetscCall(DMNetworkGetEdgeRange(dm, &eStart, &eEnd)); 73*c73416fdSAidan Hamilton PetscCall(DMNetworkGetVertexRange(dm, &vStart, &vEnd)); 74*c73416fdSAidan Hamilton PetscCall(PetscMalloc2(eEnd - eStart, &compedge, vEnd - vStart, &compvert)); 75*c73416fdSAidan Hamilton for (e = eStart; e < eEnd; e++) { 76*c73416fdSAidan Hamilton compedge[e - eStart] = e; 77*c73416fdSAidan Hamilton PetscCall(DMNetworkAddComponent(dm, e, compkey, &compedge[e - eStart], numdofedge)); 78*c73416fdSAidan Hamilton } 79*c73416fdSAidan Hamilton for (v = vStart; v < vEnd; v++) { 80*c73416fdSAidan Hamilton compvert[v - vStart] = v; 81*c73416fdSAidan Hamilton PetscCall(DMNetworkAddComponent(dm, v, compkey, &compvert[v - vStart], numdofvert)); 82*c73416fdSAidan Hamilton } 83*c73416fdSAidan Hamilton PetscCall(DMSetFromOptions(dm)); 84*c73416fdSAidan Hamilton PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 85*c73416fdSAidan Hamilton PetscCall(DMSetUp(dm)); 86*c73416fdSAidan Hamilton PetscCall(PetscFree2(compedge, compvert)); 87*c73416fdSAidan Hamilton *newdm = dm; 88*c73416fdSAidan Hamilton PetscFunctionReturn(0); 89*c73416fdSAidan Hamilton } 90*c73416fdSAidan Hamilton 91*c73416fdSAidan Hamilton PetscErrorCode StarGraphTestQuery(DM dm, PetscInt ne) 92*c73416fdSAidan Hamilton { 93*c73416fdSAidan Hamilton PetscInt globalnumvert, localnumvert, globalnumedge, localnumedge; 94*c73416fdSAidan Hamilton 95*c73416fdSAidan Hamilton PetscFunctionBegin; 96*c73416fdSAidan Hamilton PetscCall(DMNetworkGetNumEdges(dm, &localnumedge, &globalnumedge)); 97*c73416fdSAidan Hamilton PetscCall(DMNetworkGetNumVertices(dm, &localnumvert, &globalnumvert)); 98*c73416fdSAidan Hamilton 99*c73416fdSAidan Hamilton PetscCheck(globalnumedge == ne, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Global number of edges should be %" PetscInt_FMT "instead was %" PetscInt_FMT, ne, globalnumedge); 100*c73416fdSAidan Hamilton PetscCheck(globalnumvert == ne + 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Global number of vertices should be %" PetscInt_FMT "instead was %" PetscInt_FMT, ne + 1, globalnumvert); 101*c73416fdSAidan Hamilton PetscFunctionReturn(0); 102*c73416fdSAidan Hamilton } 103*c73416fdSAidan Hamilton 104*c73416fdSAidan Hamilton int main(int argc, char **argv) 105*c73416fdSAidan Hamilton { 106*c73416fdSAidan Hamilton DM dm; 107*c73416fdSAidan Hamilton PetscInt ne = 1; 108*c73416fdSAidan Hamilton PetscMPIInt rank; 109*c73416fdSAidan Hamilton 110*c73416fdSAidan Hamilton PetscFunctionBeginUser; 111*c73416fdSAidan Hamilton PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 112*c73416fdSAidan Hamilton PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank)); 113*c73416fdSAidan Hamilton /* create a distributed k-Star graph DMNetwork */ 114*c73416fdSAidan Hamilton PetscCall(PetscOptionsGetInt(NULL, NULL, "-ne", &ne, NULL)); 115*c73416fdSAidan Hamilton PetscCall(StarGraphCreate(PETSC_COMM_WORLD, 1, 0, ne, PETSC_TRUE, &dm)); 116*c73416fdSAidan Hamilton PetscCall(DMNetworkDistribute(&dm, 0)); 117*c73416fdSAidan Hamilton /* Test if query functions for DMNetwork run sucessfully */ 118*c73416fdSAidan Hamilton PetscCall(StarGraphTestQuery(dm, ne)); 119*c73416fdSAidan Hamilton PetscCall(DMDestroy(&dm)); 120*c73416fdSAidan Hamilton PetscCall(PetscFinalize()); 121*c73416fdSAidan Hamilton } 122*c73416fdSAidan Hamilton 123*c73416fdSAidan Hamilton /*TEST 124*c73416fdSAidan Hamilton test: 125*c73416fdSAidan Hamilton suffix: 0 126*c73416fdSAidan Hamilton args: -ne 5 127*c73416fdSAidan Hamilton test: 128*c73416fdSAidan Hamilton suffix: 1 129*c73416fdSAidan Hamilton nsize: 2 130*c73416fdSAidan Hamilton args: -ne 5 -petscpartitioner_type simple 131*c73416fdSAidan Hamilton TEST*/ 132