xref: /phasta/phSolver/common/test/phIOreadIlwork.cc (revision eeef5fd685651e13b19a5148edf0b6c69dc87f53)
1560e081fSCameron Smith #include <mpi.h>
2560e081fSCameron Smith #include <stdio.h>
3560e081fSCameron Smith #include <stdlib.h>
4560e081fSCameron Smith #include <unistd.h>
54894d94eSMichel Rasquin #include <set>
6*eeef5fd6SCameron Smith #include <string>
7*eeef5fd6SCameron Smith #include <sstream>
8560e081fSCameron Smith #include "phIO.h"
9560e081fSCameron Smith #include "posixio.h"
10560e081fSCameron Smith #include "phio_posix.h"
11560e081fSCameron Smith 
12560e081fSCameron Smith int main(int argc, char* argv[]) {
13560e081fSCameron Smith   MPI_Init(&argc,&argv);
14302b84e8SCameron Smith   int commrank,commsize;
15302b84e8SCameron Smith   MPI_Comm_rank(MPI_COMM_WORLD,&commrank);
16302b84e8SCameron Smith   MPI_Comm_size(MPI_COMM_WORLD,&commsize);
17*eeef5fd6SCameron Smith   if( argc != 5 ) {
18302b84e8SCameron Smith     if( !commrank )
19*eeef5fd6SCameron Smith       fprintf(stderr, "Usage: %s <geombc posix file> <verbosity=0|1|2> <rankoffset> <outfile>\n",argv[0]);
20*eeef5fd6SCameron Smith       fprintf(stderr, "verbosity=0 will only write to the specified \'outfile\'\n",argv[0]);
21*eeef5fd6SCameron Smith       fprintf(stderr, "verbosity>0 will write to the specified \'outfile\' and to stdout\n",argv[0]);
22560e081fSCameron Smith     MPI_Finalize();
23560e081fSCameron Smith     return 1;
24560e081fSCameron Smith   }
25560e081fSCameron Smith   const char* filename = argv[1];
26302b84e8SCameron Smith   const int verbosity = atoi(argv[2]);
27302b84e8SCameron Smith   const int rankoffset = atoi(argv[3]);
28*eeef5fd6SCameron Smith   char* outfilename = argv[4];
29560e081fSCameron Smith   const char* phrase = "ilwork";
30560e081fSCameron Smith   const char* type = "integer";
31560e081fSCameron Smith   const char* iotype = "binary";
32560e081fSCameron Smith   int* ilwork = NULL;
33560e081fSCameron Smith   int len = 0;
34560e081fSCameron Smith   int one = 2;
35560e081fSCameron Smith 
36302b84e8SCameron Smith 
37302b84e8SCameron Smith   int phastarank = commrank+1+rankoffset;
38302b84e8SCameron Smith 
39302b84e8SCameron Smith   char geomfilename[4096];
40302b84e8SCameron Smith   sprintf(geomfilename, "%s/geombc.dat.%d",filename,phastarank);
41302b84e8SCameron Smith 
42560e081fSCameron Smith   phio_fp file;
43560e081fSCameron Smith   posixio_setup(&(file), 'r');
44302b84e8SCameron Smith   posix_openfile_single(geomfilename, file);
45560e081fSCameron Smith   phio_readheader(file, phrase, &len, &one, type, iotype);
46560e081fSCameron Smith   ilwork = (int*) malloc(len*sizeof(int));
47560e081fSCameron Smith   phio_readdatablock(file, phrase, ilwork, &len, type, iotype);
48560e081fSCameron Smith   phio_closefile(file);
494894d94eSMichel Rasquin 
504894d94eSMichel Rasquin   // Initialization
514894d94eSMichel Rasquin   int itkbeg = 0;
524894d94eSMichel Rasquin   int m = 0;
534894d94eSMichel Rasquin   int idl = 0;
544894d94eSMichel Rasquin   std::set<int> neighbors;
554894d94eSMichel Rasquin 
564894d94eSMichel Rasquin   // Compute summary info first such as number of communication tasks, neighboring parts, onwer parts, etc
574894d94eSMichel Rasquin   int numtask = ilwork[0];
584894d94eSMichel Rasquin   int numowner = 0;
594894d94eSMichel Rasquin   for(int itask=0;itask<numtask;itask++) {
604894d94eSMichel Rasquin     int iacc   = ilwork [itkbeg + 2];
614894d94eSMichel Rasquin     int iother = ilwork [itkbeg + 3];
624894d94eSMichel Rasquin     int numseg = ilwork [itkbeg + 4];
634894d94eSMichel Rasquin     if(iacc == 1) numowner++;
644894d94eSMichel Rasquin     neighbors.insert(iother);
654894d94eSMichel Rasquin     itkbeg = itkbeg + 4 + 2*numseg;
664894d94eSMichel Rasquin   }
67*eeef5fd6SCameron Smith   const int numints = 5;
68*eeef5fd6SCameron Smith   MPI_Status status;
69*eeef5fd6SCameron Smith   MPI_File outfile;
70*eeef5fd6SCameron Smith   MPI_File_open(MPI_COMM_WORLD,outfilename,
71*eeef5fd6SCameron Smith       MPI_MODE_CREATE|MPI_MODE_WRONLY,MPI_INFO_NULL,&outfile);
72*eeef5fd6SCameron Smith   std::string header("rank peers tasks owned notowned\n");
73*eeef5fd6SCameron Smith   if( !commrank ) //write header
74*eeef5fd6SCameron Smith     MPI_File_write_at(outfile,0,(void*)header.c_str(),header.size(),MPI_CHAR,&status);
75*eeef5fd6SCameron Smith   std::stringstream ss;
76*eeef5fd6SCameron Smith   ss << phastarank << ","
77*eeef5fd6SCameron Smith      << neighbors.size() << ","
78*eeef5fd6SCameron Smith      << numtask << ","
79*eeef5fd6SCameron Smith      << numowner << ","
80*eeef5fd6SCameron Smith      << numtask-numowner << "\n";
81*eeef5fd6SCameron Smith   std::string s = ss.str();
82*eeef5fd6SCameron Smith   int size = s.size();
83*eeef5fd6SCameron Smith   int offset = 0;
84*eeef5fd6SCameron Smith   MPI_Exscan(&size,&offset,1,MPI_INT,MPI_SUM,MPI_COMM_WORLD);
85*eeef5fd6SCameron Smith   offset += header.size();
86*eeef5fd6SCameron Smith   MPI_File_write_at(outfile,offset,(void*)s.c_str(),s.size(),MPI_CHAR,&status);
87302b84e8SCameron Smith   if( verbosity > 0 ) {
884894d94eSMichel Rasquin     // Print now communication info
894894d94eSMichel Rasquin     printf("\n");
904894d94eSMichel Rasquin     printf("Communication info for this part:\n");
914894d94eSMichel Rasquin     itkbeg = 0;
924894d94eSMichel Rasquin     for(int itask=0;itask<numtask;itask++) {
934894d94eSMichel Rasquin       int itag   = ilwork [itkbeg + 1];
944894d94eSMichel Rasquin       int iacc   = ilwork [itkbeg + 2];
954894d94eSMichel Rasquin       int iother = ilwork [itkbeg + 3];
964894d94eSMichel Rasquin       int numseg = ilwork [itkbeg + 4];
974894d94eSMichel Rasquin       int isgbeg = ilwork [itkbeg + 5];
984894d94eSMichel Rasquin 
994894d94eSMichel Rasquin       // Toal number of nodes involved in this communication (lfront), including all segments
1004894d94eSMichel Rasquin       int lfront = 0;
1014894d94eSMichel Rasquin       for(int is=1;is<=numseg;is++) {
1024894d94eSMichel Rasquin         int lenseg = ilwork [itkbeg + 4 + 2*is];
1034894d94eSMichel Rasquin         lfront = lfront + lenseg;
1044894d94eSMichel Rasquin       }
1054894d94eSMichel Rasquin       printf("Communication %d:\ttag %d\townership %d\trank %d\tnumseg %d\tnumvtx %d\n",itask,itag,iacc,iother,numseg,lfront);
1064894d94eSMichel Rasquin       itkbeg = itkbeg + 4 + 2*numseg;
1074894d94eSMichel Rasquin     }
108302b84e8SCameron Smith     printf("\n");
109302b84e8SCameron Smith   }
1104894d94eSMichel Rasquin 
1114894d94eSMichel Rasquin   // Print now the raw ilwork array
112302b84e8SCameron Smith   if( verbosity > 1 ) {
1134894d94eSMichel Rasquin     printf("ilwork array:\n");
1144894d94eSMichel Rasquin     for(int i=0;i<len;i++) {
1154894d94eSMichel Rasquin       printf("%d\n",ilwork[i]);
1164894d94eSMichel Rasquin     }
11725b215f0SCameron Smith   }
1184894d94eSMichel Rasquin 
119560e081fSCameron Smith   free(ilwork);
120*eeef5fd6SCameron Smith   MPI_File_close(&outfile);
121560e081fSCameron Smith   MPI_Finalize();
122560e081fSCameron Smith   return 0;
123560e081fSCameron Smith }
124