GitLab maintenance scheduled for Tomorrow, 2020-08-11, from 17:00 to 18:00 CT - Services will be unavailable during this time.

miniHACC-SoA.cpp 9.13 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15
#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[])
{
16
  int world_numtasks, world_myrank, mycolor, mykey, sub_numtasks, sub_myrank, i, file_id;
17
  int64_t num_particles;
18 19 20 21 22 23 24 25 26 27 28 29 30 31 32
  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_Status status;
  char output[100];
  Tapioca tp;
  int64_t chunkCount[9], chunkOffset[9];
  int chunkSize[9];
  
  MPI_Init(&argc, &argv);
  MPI_Comm_size(MPI_COMM_WORLD, &world_numtasks);
  MPI_Comm_rank(MPI_COMM_WORLD, &world_myrank);

33
  mycolor = tp.topology.BridgeNodeId ();
34 35 36 37 38 39
  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);

40 41
  num_particles = atoi ( argv[1] );

42
  file_id = 0;
43 44
  if ( argv[2] != NULL )
    file_id = atoi ( argv[2] );
45 46
       
#ifdef BGQ
47
  snprintf (output, 100, "/projects/visualization/ftessier/debug/HACC-SOA-%08d-%d.dat", mycolor, file_id);
48
#elif XC40
49
  snprintf (output, 100, "/lus/theta-fs0/projects/Performance/ftessier/HACC/HACC-AOS-%08d-%d.dat", mycolor, file_id);
50
#else
51
  snprintf (output, 100, "./HACC-SOA-%08d-%d.dat", mycolor, file_id);
52 53
#endif

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 90 91 92 93 94

  /*****************/
  /*     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;

  MPI_Exscan (&num_particles, &scan_size, 1, MPI_LONG_LONG, MPI_SUM, sub_comm);
    
  if (0 == sub_myrank) {
95
    fprintf (stdout, GREEN "[INFO]" RESET " [%08d] miniHACC-SoA\n", mycolor);
96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125
    fprintf (stdout, GREEN "[INFO]" RESET " [%08d] Write output file (SoA 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 * chunkSize[0];
  for ( i = 1; i < 9; i++ ) {
    chunkOffset[i] = chunkOffset[i - 1];
    chunkOffset[i] += (sub_particles - scan_size) * chunkSize[i - 1];
    chunkOffset[i] += scan_size * chunkSize[i];
  }

126
  tp.Init (chunkCount, chunkSize, chunkOffset, 9, hdr, sub_comm);
127
  tp.setAggregationTier (2, NVR, "/scratch/tmp");
128
  tp.setTargetTier (HDD, file_size, output);
129 130 131 132 133 134
  /*****************/

  start_time = MPI_Wtime();

  offset = scan_size * sizeof(float);

135
  tp.Write (offset, xx, num_particles, MPI_FLOAT, &status);
136 137
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(float);

138
  tp.Write (offset, yy, num_particles, MPI_FLOAT, &status);
139 140
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(float);

141
  tp.Write (offset, zz, num_particles, MPI_FLOAT, &status);
142 143
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(float);

144
  tp.Write (offset, vx, num_particles, MPI_FLOAT, &status);
145 146
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(float);

147
  tp.Write (offset, vy, num_particles, MPI_FLOAT, &status);
148 149
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(float);

150
  tp.Write (offset, vz, num_particles, MPI_FLOAT, &status);
151 152
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(float);

153
  tp.Write (offset, phi, num_particles, MPI_FLOAT, &status);
154 155
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(int64_t);

156
  tp.Write (offset, pid, num_particles, MPI_LONG_LONG, &status);
157 158
  offset += (sub_particles - scan_size) * sizeof(int64_t) + scan_size * sizeof(uint16_t);

159
  tp.Write (offset, mask, num_particles, MPI_UNSIGNED_SHORT, &status);
160 161 162 163 164 165 166 167 168 169 170 171 172
  
  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);

173 174
  tp.Finalize ();

175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191
  /*****************/
  /*     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];

192 193 194
  /*****************/
  /* INIT TAPIOCA  */
  /*****************/
195
  tp.Init (chunkCount, chunkSize, chunkOffset, 9, hdr, sub_comm);
196
  tp.setAggregationTier (2, NVR, "/scratch/tmp");
197
  tp.setTargetTier (HDD, file_size, output);
198 199
  /*****************/

200 201 202 203 204 205 206
  start_time = MPI_Wtime();

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

  offset = scan_size * sizeof(float);

207
  tp.Read (offset, xx_r, num_particles, MPI_FLOAT, &status);
208 209
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(float);

210
  tp.Read (offset, yy_r, num_particles, MPI_FLOAT, &status);
211 212
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(float);

213
  tp.Read (offset, zz_r, num_particles, MPI_FLOAT, &status);
214 215
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(float);

216
  tp.Read (offset, vx_r, num_particles, MPI_FLOAT, &status);
217 218
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(float);

219
  tp.Read (offset, vy_r, num_particles, MPI_FLOAT, &status);
220 221
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(float);

222
  tp.Read (offset, vz_r, num_particles, MPI_FLOAT, &status);
223 224
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(float);

225
  tp.Read (offset, phi_r, num_particles, MPI_FLOAT, &status);
226 227
  offset += (sub_particles - scan_size) * sizeof(float) + scan_size * sizeof(int64_t);

228
  tp.Read (offset, pid_r, num_particles, MPI_LONG_LONG, &status);
229 230
  offset += (sub_particles - scan_size) * sizeof(int64_t) + scan_size * sizeof(uint16_t);

231
  tp.Read (offset, mask_r, num_particles, MPI_UNSIGNED_SHORT, &status);
232 233 234 235 236 237 238 239 240 241 242
  
  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);
  }

243 244
  tp.Finalize ();

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 275 276 277 278 279 280 281 282 283 284 285
  /*****************/
  /* 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 ();
}