1560e081fSCameron Smith #include <mpi.h> 2560e081fSCameron Smith #include <stdio.h> 3560e081fSCameron Smith #include <stdlib.h> 4560e081fSCameron Smith #include <unistd.h> 5*4894d94eSMichel Rasquin #include <set> 6560e081fSCameron Smith #include "phIO.h" 7560e081fSCameron Smith #include "posixio.h" 8560e081fSCameron Smith #include "phio_posix.h" 9560e081fSCameron Smith 10560e081fSCameron Smith int main(int argc, char* argv[]) { 11560e081fSCameron Smith MPI_Init(&argc,&argv); 12560e081fSCameron Smith if( argc != 2 ) { 13560e081fSCameron Smith fprintf(stderr, "Usage: %s <geombc posix file>\n",argv[0]); 14560e081fSCameron Smith MPI_Finalize(); 15560e081fSCameron Smith return 1; 16560e081fSCameron Smith } 17560e081fSCameron Smith const char* filename = argv[1]; 18560e081fSCameron Smith const char* phrase = "ilwork"; 19560e081fSCameron Smith const char* type = "integer"; 20560e081fSCameron Smith const char* iotype = "binary"; 21560e081fSCameron Smith int* ilwork = NULL; 22560e081fSCameron Smith int len = 0; 23560e081fSCameron Smith int one = 2; 24560e081fSCameron Smith 25560e081fSCameron Smith phio_fp file; 26560e081fSCameron Smith posixio_setup(&(file), 'r'); 27560e081fSCameron Smith posix_openfile_single(filename, file); 28560e081fSCameron Smith phio_readheader(file, phrase, &len, &one, type, iotype); 29560e081fSCameron Smith fprintf(stderr, "len %d\n", len); 30560e081fSCameron Smith ilwork = (int*) malloc(len*sizeof(int)); 31560e081fSCameron Smith phio_readdatablock(file, phrase, ilwork, &len, type, iotype); 32560e081fSCameron Smith phio_closefile(file); 33*4894d94eSMichel Rasquin 34*4894d94eSMichel Rasquin // Initialization 35*4894d94eSMichel Rasquin int itkbeg = 0; 36*4894d94eSMichel Rasquin int m = 0; 37*4894d94eSMichel Rasquin int idl = 0; 38*4894d94eSMichel Rasquin std::set<int> neighbors; 39*4894d94eSMichel Rasquin 40*4894d94eSMichel Rasquin // Compute summary info first such as number of communication tasks, neighboring parts, onwer parts, etc 41*4894d94eSMichel Rasquin int numtask = ilwork[0]; 42*4894d94eSMichel Rasquin int numowner = 0; 43*4894d94eSMichel Rasquin for(int itask=0;itask<numtask;itask++) { 44*4894d94eSMichel Rasquin int iacc = ilwork [itkbeg + 2]; 45*4894d94eSMichel Rasquin int iother = ilwork [itkbeg + 3]; 46*4894d94eSMichel Rasquin int numseg = ilwork [itkbeg + 4]; 47*4894d94eSMichel Rasquin if(iacc == 1) numowner++; 48*4894d94eSMichel Rasquin neighbors.insert(iother); 49*4894d94eSMichel Rasquin itkbeg = itkbeg + 4 + 2*numseg; 50*4894d94eSMichel Rasquin } 51*4894d94eSMichel Rasquin printf("Number of neighboring parts: %d\n",neighbors.size()); 52*4894d94eSMichel Rasquin printf("Number of communication tasks: %d\n",numtask); 53*4894d94eSMichel Rasquin printf("Number of owner communications: %d\n",numowner); 54*4894d94eSMichel Rasquin printf("Number of non-owner communications: %d\n",numtask-numowner); 55*4894d94eSMichel Rasquin 56*4894d94eSMichel Rasquin // Print now communication info 57*4894d94eSMichel Rasquin printf("\n"); 58*4894d94eSMichel Rasquin printf("Communication info for this part:\n"); 59*4894d94eSMichel Rasquin itkbeg = 0; 60*4894d94eSMichel Rasquin for(int itask=0;itask<numtask;itask++) { 61*4894d94eSMichel Rasquin int itag = ilwork [itkbeg + 1]; 62*4894d94eSMichel Rasquin int iacc = ilwork [itkbeg + 2]; 63*4894d94eSMichel Rasquin int iother = ilwork [itkbeg + 3]; 64*4894d94eSMichel Rasquin int numseg = ilwork [itkbeg + 4]; 65*4894d94eSMichel Rasquin int isgbeg = ilwork [itkbeg + 5]; 66*4894d94eSMichel Rasquin 67*4894d94eSMichel Rasquin // Toal number of nodes involved in this communication (lfront), including all segments 68*4894d94eSMichel Rasquin int lfront = 0; 69*4894d94eSMichel Rasquin for(int is=1;is<=numseg;is++) { 70*4894d94eSMichel Rasquin int lenseg = ilwork [itkbeg + 4 + 2*is]; 71*4894d94eSMichel Rasquin lfront = lfront + lenseg; 72*4894d94eSMichel Rasquin } 73*4894d94eSMichel Rasquin printf("Communication %d:\ttag %d\townership %d\trank %d\tnumseg %d\tnumvtx %d\n",itask,itag,iacc,iother,numseg,lfront); 74*4894d94eSMichel Rasquin itkbeg = itkbeg + 4 + 2*numseg; 75*4894d94eSMichel Rasquin } 76*4894d94eSMichel Rasquin 77*4894d94eSMichel Rasquin // Print now the raw ilwork array 78*4894d94eSMichel Rasquin printf("\n"); 79*4894d94eSMichel Rasquin printf("ilwork array:\n"); 80*4894d94eSMichel Rasquin for(int i=0;i<len;i++) { 81*4894d94eSMichel Rasquin printf("%d\n",ilwork[i]); 82*4894d94eSMichel Rasquin } 83*4894d94eSMichel Rasquin 84560e081fSCameron Smith free(ilwork); 85560e081fSCameron Smith 86560e081fSCameron Smith MPI_Finalize(); 87560e081fSCameron Smith return 0; 88560e081fSCameron Smith } 89