miniHACC-AoS.cpp 8.54 KB
Newer Older
Francois Tessier's avatar
Francois Tessier committed
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25
#include <stdio.h>
#include <stdlib.h>
#include <assert.h>
#include <stdint.h>
#include <mpi.h>

#include "tapioca.hpp"

#define RED     "\x1b[31m"
#define GREEN   "\x1b[32m"
#define BLUE    "\x1b[34m"
#define RESET   "\x1b[0m"

int main (int argc, char * argv[])
{
  int world_numtasks, world_myrank, mycolor, mykey, sub_numtasks, sub_myrank, i;
  int64_t num_particles = 25000;
  int64_t sub_particles, tot_particles, particle_size, file_size, tot_size;
  int64_t scan_size = 0, offset, hdr = 0;
  double start_time, end_time, tot_time, max_time;
  double io_bw;
  MPI_Comm sub_comm;
  MPI_File file_handle;
  MPI_Status status;
  char output[100];
26 27 28
  Tapioca tp;
  int64_t chunkCount[9], chunkOffset[9];
  int chunkSize[9];
Francois Tessier's avatar
Francois Tessier committed
29 30 31 32 33
  
  MPI_Init(&argc, &argv);
  MPI_Comm_size(MPI_COMM_WORLD, &world_numtasks);
  MPI_Comm_rank(MPI_COMM_WORLD, &world_myrank);

34
  mycolor = tp.topology.BridgeNodeId ();
Francois Tessier's avatar
Francois Tessier committed
35 36 37 38 39 40
  mykey = world_myrank;

  MPI_Comm_split (MPI_COMM_WORLD, mycolor, mykey, &sub_comm);
  MPI_Comm_size(sub_comm, &sub_numtasks);
  MPI_Comm_rank(sub_comm, &sub_myrank);

41
  snprintf (output, 100, "/projects/visualization/ftessier/debug/HACC-AOS-%08d-%d.dat", mycolor, atoi(argv[1]));
Francois Tessier's avatar
Francois Tessier committed
42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89

  /*****************/
  /*     WRITE     */
  /*****************/
  float *xx, *yy, *zz, *vx, *vy, *vz, *phi;
  int64_t* pid;
  uint16_t* mask;

  xx = new float[num_particles];
  yy = new float[num_particles];
  zz = new float[num_particles];
  vx = new float[num_particles];
  vy = new float[num_particles];
  vz = new float[num_particles];
  phi = new float[num_particles];
  pid = new int64_t[num_particles];
  mask = new uint16_t[num_particles];

  for (uint64_t i = 0; i< num_particles; i++)
    {
      xx[i] = (float)i;
      yy[i] = (float)i;
      zz[i] = (float)i;
      vx[i] = (float)i;
      vy[i] = (float)i;
      vz[i] = (float)i;
      phi[i] = (float)i;
      pid[i] =  (int64_t)i;
      mask[i] = (uint16_t)world_myrank;
    }

  MPI_Allreduce(&num_particles, &sub_particles, 1, MPI_LONG_LONG, MPI_SUM, sub_comm);
  MPI_Allreduce(&num_particles, &tot_particles, 1, MPI_LONG_LONG, MPI_SUM, MPI_COMM_WORLD);
  
  particle_size = (7 * sizeof(float)) + sizeof(int64_t) + sizeof(uint16_t);
  file_size = particle_size * sub_particles;
  tot_size = particle_size * tot_particles;

  if (sub_myrank == 0) {
    MPI_File_open(MPI_COMM_SELF, output,
		  MPI_MODE_WRONLY | MPI_MODE_CREATE, MPI_INFO_NULL, &file_handle);
    MPI_File_set_size(file_handle, file_size);
    MPI_File_close (&file_handle);
  }

  MPI_Exscan (&num_particles, &scan_size, 1, MPI_LONG_LONG, MPI_SUM, sub_comm);
    
  if (0 == sub_myrank) {
90
    fprintf (stdout, GREEN "[INFO]" RESET " [%08d] miniHACC-AoS\n", mycolor);
Francois Tessier's avatar
Francois Tessier committed
91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118
    fprintf (stdout, GREEN "[INFO]" RESET " [%08d] Write output file (AoS data layout)\n", mycolor);
    fprintf (stdout, GREEN "[INFO]" RESET " [%08d] --> %lld particles per rank\n", mycolor, num_particles);
    fprintf (stdout, GREEN "[INFO]" RESET " [%08d] --> File size: %.2f MB (%lld particles)\n", 
	     mycolor, (double)file_size/(1024*1024), sub_particles);
  }

  /*****************/
  /* INIT TAPIOCA  */
  /*****************/
  for ( i = 0; i < 9; i++ ) {
    chunkCount[i] = num_particles;
  }
  
  chunkSize[0] = sizeof(float);
  chunkSize[1] = sizeof(float);
  chunkSize[2] = sizeof(float);
  chunkSize[3] = sizeof(float);
  chunkSize[4] = sizeof(float);
  chunkSize[5] = sizeof(float);
  chunkSize[6] = sizeof(float);
  chunkSize[7] = sizeof(int64_t);
  chunkSize[8] = sizeof(uint16_t);
  
  chunkOffset[0] = hdr + scan_size * particle_size;
  for ( i = 1; i < 9; i++ ) {
    chunkOffset[i] = chunkOffset[i - 1] + chunkCount[i - 1] * chunkSize[i - 1];
  }
  
119
  tp.WriteInitialize (chunkCount, chunkSize, chunkOffset, 9, hdr, ARRAY_OF_STRUCTURES, sub_comm);
Francois Tessier's avatar
Francois Tessier committed
120 121 122 123 124 125 126 127 128
  /*****************/
  
  start_time = MPI_Wtime();

  MPI_File_open(sub_comm, output,
		MPI_MODE_WRONLY, MPI_INFO_NULL, &file_handle);
  
  offset = scan_size * particle_size;

129
  tp.Write (file_handle, offset, xx, num_particles, MPI_FLOAT, &status);
Francois Tessier's avatar
Francois Tessier committed
130 131
  offset += num_particles * sizeof(float);

132
  tp.Write (file_handle, offset, yy, num_particles, MPI_FLOAT, &status);
Francois Tessier's avatar
Francois Tessier committed
133 134
  offset += num_particles * sizeof(float);

135
  tp.Write (file_handle, offset, zz, num_particles, MPI_FLOAT, &status);
Francois Tessier's avatar
Francois Tessier committed
136 137
  offset += num_particles * sizeof(float);

138
  tp.Write (file_handle, offset, vx, num_particles, MPI_FLOAT, &status);
Francois Tessier's avatar
Francois Tessier committed
139 140
  offset += num_particles * sizeof(float);

141
  tp.Write (file_handle, offset, vy, num_particles, MPI_FLOAT, &status);
Francois Tessier's avatar
Francois Tessier committed
142 143
  offset += num_particles * sizeof(float);

144
  tp.Write (file_handle, offset, vz, num_particles, MPI_FLOAT, &status);
Francois Tessier's avatar
Francois Tessier committed
145 146
  offset += num_particles * sizeof(float);

147
  tp.Write (file_handle, offset, phi, num_particles, MPI_FLOAT, &status);
Francois Tessier's avatar
Francois Tessier committed
148 149
  offset += num_particles * sizeof(float);

150
  tp.Write (file_handle, offset, pid, num_particles, MPI_LONG_LONG, &status);
Francois Tessier's avatar
Francois Tessier committed
151 152
  offset += num_particles * sizeof(int64_t);

153
  tp.Write (file_handle, offset, mask, num_particles, MPI_UNSIGNED_SHORT, &status);
Francois Tessier's avatar
Francois Tessier committed
154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274
  
  MPI_File_close (&file_handle);

  end_time = MPI_Wtime();
  tot_time = end_time - start_time;
  MPI_Reduce (&tot_time, &max_time, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD);
  
  if (0 == world_myrank) {
    io_bw = (double)tot_size / max_time / (1024 * 1024);
    fprintf (stdout, BLUE "[TIMING]" RESET " Write I/O bandwidth: %.2f MBps (%.2f MB in %.2f ms)\n",
	     io_bw, (double)tot_size/(1024*1024), max_time * 1000);
  }

  MPI_Barrier (MPI_COMM_WORLD);

  /*****************/
  /*     READ      */
  /*****************/
  float *xx_r, *yy_r, *zz_r, *vx_r, *vy_r, *vz_r, *phi_r;
  int64_t* pid_r;
  uint16_t* mask_r;

  xx_r = new float[num_particles];
  yy_r = new float[num_particles];
  zz_r = new float[num_particles];
  vx_r = new float[num_particles];
  vy_r = new float[num_particles];
  vz_r = new float[num_particles];
  phi_r = new float[num_particles];
  pid_r = new int64_t[num_particles];
  mask_r = new uint16_t[num_particles];

  start_time = MPI_Wtime();

  MPI_File_open(sub_comm, output,
		MPI_MODE_RDONLY, MPI_INFO_NULL, &file_handle);

  if (0 == sub_myrank)
    fprintf (stdout, GREEN "[INFO]" RESET " [%08d] Read output file\n", mycolor);

  offset = scan_size * particle_size;

  MPI_File_read_at_all (file_handle, offset, xx_r, num_particles, MPI_FLOAT, &status);
  offset += num_particles * sizeof(float);

  MPI_File_read_at_all (file_handle, offset, yy_r, num_particles, MPI_FLOAT, &status);
  offset += num_particles * sizeof(float);

  MPI_File_read_at_all (file_handle, offset, zz_r, num_particles, MPI_FLOAT, &status);
  offset += num_particles * sizeof(float);

  MPI_File_read_at_all (file_handle, offset, vx_r, num_particles, MPI_FLOAT, &status);
  offset += num_particles * sizeof(float);

  MPI_File_read_at_all (file_handle, offset, vy_r, num_particles, MPI_FLOAT, &status);
  offset += num_particles * sizeof(float);

  MPI_File_read_at_all (file_handle, offset, vz_r, num_particles, MPI_FLOAT, &status);
  offset += num_particles * sizeof(float);

  MPI_File_read_at_all (file_handle, offset, phi_r, num_particles, MPI_FLOAT, &status);
  offset += num_particles * sizeof(float);

  MPI_File_read_at_all (file_handle, offset, pid_r, num_particles, MPI_LONG_LONG, &status);
  offset += num_particles * sizeof(int64_t);

  MPI_File_read_at_all (file_handle, offset, mask_r, num_particles, MPI_UNSIGNED_SHORT, &status);
  
  MPI_File_close (&file_handle);

  end_time = MPI_Wtime();
  tot_time = end_time - start_time;
  MPI_Reduce (&tot_time, &max_time, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD);
  
  if (0 == world_myrank) {
    io_bw = (double)tot_size / max_time / (1024 * 1024);
    fprintf (stdout, BLUE "[TIMING]" RESET " Read I/O bandwidth: %.2f MBps (%.2f MB in %.2f ms)\n",
	     io_bw, (double)tot_size/(1024*1024), max_time * 1000);
  }

  /*****************/
  /* VERIFICATION  */
  /*****************/
  for (uint64_t i = 0; i< num_particles; i++) {
    if ((xx[i] != xx_r[i]) || (yy[i] != yy_r[i]) || (zz[i] != zz_r[i])
	|| (vx[i] != vx_r[i]) || (vy[i] != vy_r[i]) || (vz[i] != vz_r[i])
	|| (phi[i] != phi_r[i])|| (pid[i] != pid_r[i]) || (mask[i] != mask_r[i]))
      {
	fprintf (stdout, RED "[ERROR]" RESET " Wrong value for particle %d\n", i);
	MPI_Abort (MPI_COMM_WORLD, -1);
      }
  }

  if (0 == sub_myrank)
    fprintf (stdout, GREEN "[INFO]" RESET " [%08d] Content verified and consistent\n", mycolor);

  /*****************/
  /*      FREE     */
  /*****************/
  delete [] xx;
  delete [] xx_r;
  delete [] yy;
  delete [] yy_r;
  delete [] zz;
  delete [] zz_r;
  delete [] vx;
  delete [] vx_r;
  delete [] vy;
  delete [] vy_r;
  delete [] vz;
  delete [] vz_r;
  delete [] phi;
  delete [] phi_r;
  delete [] pid;
  delete [] pid_r;
  delete [] mask;
  delete [] mask_r;
  
  MPI_Finalize ();
}