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