1560e081fSCameron Smith #include <mpi.h> 2560e081fSCameron Smith #include <stdio.h> 3560e081fSCameron Smith #include <stdlib.h> 4560e081fSCameron Smith #include <unistd.h> 54894d94eSMichel 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); 12*25b215f0SCameron Smith if( argc != 3 ) { 13*25b215f0SCameron Smith fprintf(stderr, "Usage: %s <geombc posix file> <write ilwork=0|1>\n",argv[0]); 14560e081fSCameron Smith MPI_Finalize(); 15560e081fSCameron Smith return 1; 16560e081fSCameron Smith } 17560e081fSCameron Smith const char* filename = argv[1]; 18*25b215f0SCameron Smith const int writeilwork = atoi(argv[2]); 19560e081fSCameron Smith const char* phrase = "ilwork"; 20560e081fSCameron Smith const char* type = "integer"; 21560e081fSCameron Smith const char* iotype = "binary"; 22560e081fSCameron Smith int* ilwork = NULL; 23560e081fSCameron Smith int len = 0; 24560e081fSCameron Smith int one = 2; 25560e081fSCameron Smith 26560e081fSCameron Smith phio_fp file; 27560e081fSCameron Smith posixio_setup(&(file), 'r'); 28560e081fSCameron Smith posix_openfile_single(filename, file); 29560e081fSCameron Smith phio_readheader(file, phrase, &len, &one, type, iotype); 30560e081fSCameron Smith fprintf(stderr, "len %d\n", len); 31560e081fSCameron Smith ilwork = (int*) malloc(len*sizeof(int)); 32560e081fSCameron Smith phio_readdatablock(file, phrase, ilwork, &len, type, iotype); 33560e081fSCameron Smith phio_closefile(file); 344894d94eSMichel Rasquin 354894d94eSMichel Rasquin // Initialization 364894d94eSMichel Rasquin int itkbeg = 0; 374894d94eSMichel Rasquin int m = 0; 384894d94eSMichel Rasquin int idl = 0; 394894d94eSMichel Rasquin std::set<int> neighbors; 404894d94eSMichel Rasquin 414894d94eSMichel Rasquin // Compute summary info first such as number of communication tasks, neighboring parts, onwer parts, etc 424894d94eSMichel Rasquin int numtask = ilwork[0]; 434894d94eSMichel Rasquin int numowner = 0; 444894d94eSMichel Rasquin for(int itask=0;itask<numtask;itask++) { 454894d94eSMichel Rasquin int iacc = ilwork [itkbeg + 2]; 464894d94eSMichel Rasquin int iother = ilwork [itkbeg + 3]; 474894d94eSMichel Rasquin int numseg = ilwork [itkbeg + 4]; 484894d94eSMichel Rasquin if(iacc == 1) numowner++; 494894d94eSMichel Rasquin neighbors.insert(iother); 504894d94eSMichel Rasquin itkbeg = itkbeg + 4 + 2*numseg; 514894d94eSMichel Rasquin } 524894d94eSMichel Rasquin printf("Number of neighboring parts: %d\n",neighbors.size()); 534894d94eSMichel Rasquin printf("Number of communication tasks: %d\n",numtask); 544894d94eSMichel Rasquin printf("Number of owner communications: %d\n",numowner); 554894d94eSMichel Rasquin printf("Number of non-owner communications: %d\n",numtask-numowner); 564894d94eSMichel Rasquin 574894d94eSMichel Rasquin // Print now communication info 584894d94eSMichel Rasquin printf("\n"); 594894d94eSMichel Rasquin printf("Communication info for this part:\n"); 604894d94eSMichel Rasquin itkbeg = 0; 614894d94eSMichel Rasquin for(int itask=0;itask<numtask;itask++) { 624894d94eSMichel Rasquin int itag = ilwork [itkbeg + 1]; 634894d94eSMichel Rasquin int iacc = ilwork [itkbeg + 2]; 644894d94eSMichel Rasquin int iother = ilwork [itkbeg + 3]; 654894d94eSMichel Rasquin int numseg = ilwork [itkbeg + 4]; 664894d94eSMichel Rasquin int isgbeg = ilwork [itkbeg + 5]; 674894d94eSMichel Rasquin 684894d94eSMichel Rasquin // Toal number of nodes involved in this communication (lfront), including all segments 694894d94eSMichel Rasquin int lfront = 0; 704894d94eSMichel Rasquin for(int is=1;is<=numseg;is++) { 714894d94eSMichel Rasquin int lenseg = ilwork [itkbeg + 4 + 2*is]; 724894d94eSMichel Rasquin lfront = lfront + lenseg; 734894d94eSMichel Rasquin } 744894d94eSMichel Rasquin printf("Communication %d:\ttag %d\townership %d\trank %d\tnumseg %d\tnumvtx %d\n",itask,itag,iacc,iother,numseg,lfront); 754894d94eSMichel Rasquin itkbeg = itkbeg + 4 + 2*numseg; 764894d94eSMichel Rasquin } 774894d94eSMichel Rasquin 784894d94eSMichel Rasquin // Print now the raw ilwork array 794894d94eSMichel Rasquin printf("\n"); 80*25b215f0SCameron Smith if( writeilwork ) { 814894d94eSMichel Rasquin printf("ilwork array:\n"); 824894d94eSMichel Rasquin for(int i=0;i<len;i++) { 834894d94eSMichel Rasquin printf("%d\n",ilwork[i]); 844894d94eSMichel Rasquin } 85*25b215f0SCameron Smith } 864894d94eSMichel Rasquin 87560e081fSCameron Smith free(ilwork); 88560e081fSCameron Smith 89560e081fSCameron Smith MPI_Finalize(); 90560e081fSCameron Smith return 0; 91560e081fSCameron Smith } 92