serial_io.c 23.7 KB
Newer Older
1
2
3
4
/*******************************************************************************
 * This file is part of SWIFT.
 * Coypright (c) 2012 Pedro Gonnet (pedro.gonnet@durham.ac.uk),
 *                    Matthieu Schaller (matthieu.schaller@durham.ac.uk).
5
 *
6
7
8
9
 * This program is free software: you can redistribute it and/or modify
 * it under the terms of the GNU Lesser General Public License as published
 * by the Free Software Foundation, either version 3 of the License, or
 * (at your option) any later version.
10
 *
11
12
13
14
 * This program is distributed in the hope that it will be useful,
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU General Public License for more details.
15
 *
16
17
 * You should have received a copy of the GNU Lesser General Public License
 * along with this program.  If not, see <http://www.gnu.org/licenses/>.
18
 *
19
20
21
22
23
 ******************************************************************************/

/* Config parameters. */
#include "../config.h"

24
#if defined(HAVE_HDF5) && defined(WITH_MPI) && !defined(HAVE_PARALLEL_HDF5)
25
26
27
28
29
30
31

/* Some standard headers. */
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <stddef.h>
#include <hdf5.h>
32
#include <math.h>
33

34
35
#include "mpi.h"

36
#include "const.h"
Pedro Gonnet's avatar
Pedro Gonnet committed
37
#include "cycle.h"
38
39
40
#include "lock.h"
#include "task.h"
#include "space.h"
41
#include "scheduler.h"
42
#include "engine.h"
43
#include "error.h"
44
#include "kernel.h"
45
#include "common_io.h"
46
47
48
49
50
51
52
53
54
55
56
57
58

/*-----------------------------------------------------------------------------
 * Routines reading an IC file
 *-----------------------------------------------------------------------------*/

/**
 * @brief Reads a data array from a given HDF5 group.
 *
 * @param grp The group from which to read.
 * @param name The name of the array to read.
 * @param type The #DATA_TYPE of the attribute.
 * @param N The number of particles.
 * @param dim The dimension of the data (1 for scalar, 3 for vector)
59
60
61
62
 * @param part_c A (char*) pointer on the first occurence of the field of
 *interest in the parts array
 * @param importance If COMPULSORY, the data must be present in the IC file. If
 *OPTIONAL, the array will be zeroed when the data is not present.
63
 *
64
65
 * @todo A better version using HDF5 hyperslabs to read the file directly into
 *the part array
66
 * will be written once the strucutres have been stabilized.
67
 *
68
69
 * Calls #error() if an error occurs.
 */
70
71
72
73
void readArrayBackEnd(hid_t grp, char* name, enum DATA_TYPE type, int N,
                      int dim, long long N_total, long long offset,
                      char* part_c, enum DATA_IMPORTANCE importance) {
  hid_t h_data = 0, h_err = 0, h_type = 0, h_memspace = 0, h_filespace = 0;
74
  hsize_t shape[2], offsets[2];
75
  htri_t exist = 0;
76
  void* temp;
77
  int i = 0, rank = 0;
78
79
80
81
82
83
84
  const size_t typeSize = sizeOfType(type);
  const size_t copySize = typeSize * dim;
  const size_t partSize = sizeof(struct part);
  char* temp_c = 0;

  /* Check whether the dataspace exists or not */
  exist = H5Lexists(grp, name, 0);
85
86
87
88
89
90
91
92
  if (exist < 0) {
    error("Error while checking the existence of data set '%s'.", name);
  } else if (exist == 0) {
    if (importance == COMPULSORY) {
      error("Compulsory data set '%s' not present in the file.", name);
    } else {
      for (i = 0; i < N; ++i) memset(part_c + i * partSize, 0, copySize);
      return;
93
    }
94
  }
95

96
97
  /* message( "Reading %s '%s' array...", importance == COMPULSORY ?
   * "compulsory": "optional  ", name); */
98
99
100

  /* Open data space */
  h_data = H5Dopen1(grp, name);
101
  if (h_data < 0) error("Error while opening data space '%s'.", name);
102
103
104

  /* Check data type */
  h_type = H5Dget_type(h_data);
105
106
  if (h_type < 0) error("Unable to retrieve data type from the file");
  if (!H5Tequal(h_type, hdf5Type(type)))
107
    error("Non-matching types between the code and the file");
108

109
110
  /* Allocate temporary buffer */
  temp = malloc(N * dim * sizeOfType(type));
111
  if (temp == NULL) error("Unable to allocate memory for temporary buffer");
112

113
  /* Prepare information for hyperslab */
114
115
116
117
118
119
120
121
122
123
124
125
126
  if (dim > 1) {
    rank = 2;
    shape[0] = N;
    shape[1] = dim;
    offsets[0] = offset;
    offsets[1] = 0;
  } else {
    rank = 1;
    shape[0] = N;
    shape[1] = 0;
    offsets[0] = offset;
    offsets[1] = 0;
  }
127
128
129
130
131
132
133
134

  /* Create data space in memory */
  h_memspace = H5Screate_simple(rank, shape, NULL);

  /* Select hyperslab in file */
  h_filespace = H5Dget_space(h_data);
  H5Sselect_hyperslab(h_filespace, H5S_SELECT_SET, offsets, NULL, shape, NULL);

135
136
137
  /* Read HDF5 dataspace in temporary buffer */
  /* Dirty version that happens to work for vectors but should be improved */
  /* Using HDF5 dataspaces would be better */
138
139
140
141
142
  h_err = H5Dread(h_data, hdf5Type(type), h_memspace, h_filespace, H5P_DEFAULT,
                  temp);
  if (h_err < 0) {
    error("Error while reading data array '%s'.", name);
  }
143
144
145

  /* Copy temporary buffer to particle data */
  temp_c = temp;
146
147
148
  for (i = 0; i < N; ++i)
    memcpy(part_c + i * partSize, &temp_c[i * copySize], copySize);

149
150
  /* Free and close everything */
  free(temp);
151
152
  H5Sclose(h_filespace);
  H5Sclose(h_memspace);
153
154
155
156
157
158
159
160
161
162
163
164
165
  H5Tclose(h_type);
  H5Dclose(h_data);
}

/**
 * @brief A helper macro to call the readArrayBackEnd function more easily.
 *
 * @param grp The group from which to read.
 * @param name The name of the array to read.
 * @param type The #DATA_TYPE of the attribute.
 * @param N The number of particles.
 * @param dim The dimension of the data (1 for scalar, 3 for vector)
 * @param part The array of particles to fill
166
167
 * @param N_total Total number of particles
 * @param offset Offset in the array where this task starts reading
168
169
170
171
 * @param field The name of the field (C code name as defined in part.h) to fill
 * @param importance Is the data compulsory or not
 *
 */
172
173
174
175
#define readArray(grp, name, type, N, dim, part, N_total, offset, field, \
                  importance)                                            \
  readArrayBackEnd(grp, name, type, N, dim, N_total, offset,             \
                   (char*)(&(part[0]).field), importance)
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195

/**
 * @brief Reads an HDF5 initial condition file (GADGET-3 type)
 *
 * @param fileName The file to read.
 * @param dim (output) The dimension of the volume read from the file.
 * @param parts (output) The array of #part read from the file.
 * @param N (output) The number of particles read from the file.
 * @param periodic (output) 1 if the volume is periodic, 0 if not.
 *
 * Opens the HDF5 file fileName and reads the particles contained
 * in the parts array. N is the returned number of particles found
 * in the file.
 *
 * @warning Can not read snapshot distributed over more than 1 file !!!
 * @todo Read snaphsots distributed in more than one file.
 *
 * Calls #error() if an error occurs.
 *
 */
196
197
198
199
200
201
202
203
204
void read_ic_serial(char* fileName, double dim[3], struct part** parts, int* N,
                    int* periodic, int mpi_rank, int mpi_size, MPI_Comm comm,
                    MPI_Info info) {
  hid_t h_file = 0, h_grp = 0;
  double boxSize[3] = {
      0.0, -1.0, -1.0}; /* GADGET has only cubic boxes (in cosmological mode) */
  int numParticles[6] = {
      0}; /* GADGET has 6 particle types. We only keep the type 0*/
  int numParticles_highWord[6] = {0};
205
206
207
208
209
210
211
212
213
214
  long long offset = 0;
  long long N_total = 0;
  int rank;

  /* First read some information about the content */
  if (mpi_rank == 0) {

    /* Open file */
    /* message("Opening file '%s' as IC.", fileName); */
    h_file = H5Fopen(fileName, H5F_ACC_RDONLY, H5P_DEFAULT);
215
216
217
    if (h_file < 0)
      error("Error while opening file '%s' for inital read.", fileName);

218
219
220
    /* Open header to read simulation properties */
    /* message("Reading runtime parameters..."); */
    h_grp = H5Gopen1(h_file, "/RuntimePars");
221
222
    if (h_grp < 0) error("Error while opening runtime parameters\n");

223
224
    /* Read the relevant information */
    readAttribute(h_grp, "PeriodicBoundariesOn", INT, periodic);
225

226
227
    /* Close runtime parameters */
    H5Gclose(h_grp);
228

229
230
231
    /* Open header to read simulation properties */
    /* message("Reading file header..."); */
    h_grp = H5Gopen1(h_file, "/Header");
232
233
    if (h_grp < 0) error("Error while opening file header\n");

234
235
236
237
238
    /* Read the relevant information and print status */
    readAttribute(h_grp, "BoxSize", DOUBLE, boxSize);
    readAttribute(h_grp, "NumPart_Total", UINT, numParticles);
    readAttribute(h_grp, "NumPart_Total_HighWord", UINT, numParticles_highWord);

239
240
    N_total = ((long long)numParticles[0]) +
              ((long long)numParticles_highWord[0] << 32);
241
    dim[0] = boxSize[0];
242
243
    dim[1] = (boxSize[1] < 0) ? boxSize[0] : boxSize[1];
    dim[2] = (boxSize[2] < 0) ? boxSize[0] : boxSize[2];
244
245
246

    /* message("Found %d particles in a %speriodic box of size [%f %f %f].",  */
    /* 	 *N, (periodic ? "": "non-"), dim[0], dim[1], dim[2]); */
247

248
249
    /* Close header */
    H5Gclose(h_grp);
250

251
252
253
    /* Close file */
    H5Fclose(h_file);
  }
254

255
256
257
258
  /* Now need to broadcast that information to all ranks. */
  MPI_Bcast(periodic, 1, MPI_INT, 0, comm);
  MPI_Bcast(&N_total, 1, MPI_LONG_LONG, 0, comm);
  MPI_Bcast(dim, 3, MPI_DOUBLE, 0, comm);
259

260
261
262
  /* Divide the particles among the tasks. */
  offset = mpi_rank * N_total / mpi_size;
  *N = (mpi_rank + 1) * N_total / mpi_size - offset;
263
264

  /* Allocate memory to store particles */
265
  if (posix_memalign((void*)parts, part_align, (*N) * sizeof(struct part)) != 0)
266
    error("Error while allocating memory for particles");
267
268
269
  bzero(*parts, *N * sizeof(struct part));
  /* message("Allocated %8.2f MB for particles.", *N * sizeof(struct part) /
   * (1024.*1024.)); */
270
271

  /* Now loop over ranks and read the data */
272
  for (rank = 0; rank < mpi_size; ++rank) {
273
274

    /* Is it this rank's turn to read ? */
275
    if (rank == mpi_rank) {
276
277

      h_file = H5Fopen(fileName, H5F_ACC_RDONLY, H5P_DEFAULT);
278
279
280
      if (h_file < 0)
        error("Error while opening file '%s' on rank %d.", fileName, mpi_rank);

281
282
283
      /* Open SPH particles group */
      /* message("Reading particle arrays..."); */
      h_grp = H5Gopen1(h_file, "/PartType0");
284
285
286
      if (h_grp < 0)
        error("Error while opening particle group on rank %d.\n", mpi_rank);

287
      /* Read arrays */
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
      readArray(h_grp, "Coordinates", DOUBLE, *N, 3, *parts, N_total, offset, x,
                COMPULSORY);
      readArray(h_grp, "Velocities", FLOAT, *N, 3, *parts, N_total, offset, v,
                COMPULSORY);
      readArray(h_grp, "Masses", FLOAT, *N, 1, *parts, N_total, offset, mass,
                COMPULSORY);
      readArray(h_grp, "SmoothingLength", FLOAT, *N, 1, *parts, N_total, offset,
                h, COMPULSORY);
      readArray(h_grp, "InternalEnergy", FLOAT, *N, 1, *parts, N_total, offset,
                u, COMPULSORY);
      readArray(h_grp, "ParticleIDs", ULONGLONG, *N, 1, *parts, N_total, offset,
                id, COMPULSORY);
      readArray(h_grp, "TimeStep", FLOAT, *N, 1, *parts, N_total, offset, dt,
                OPTIONAL);
      readArray(h_grp, "Acceleration", FLOAT, *N, 3, *parts, N_total, offset, a,
                OPTIONAL);
      readArray(h_grp, "Density", FLOAT, *N, 1, *parts, N_total, offset, rho,
                OPTIONAL);

307
308
309
310
311
312
313
314
315
316
      /* Close particle group */
      H5Gclose(h_grp);

      /* Close file */
      H5Fclose(h_file);
    }

    /* Wait for the read of the reading to complete */
    MPI_Barrier(comm);
  }
317

318
  /* message("Done Reading particles..."); */
319
320
321
322
323
}

/*-----------------------------------------------------------------------------
 * Routines writing an output file
 *-----------------------------------------------------------------------------*/
324

325
326
327
328
void prepareArray(hid_t grp, char* fileName, FILE* xmfFile, char* name,
                  enum DATA_TYPE type, long long N_total, int dim,
                  struct UnitSystem* us, enum UnitConversionFactor convFactor) {
  hid_t h_data = 0, h_err = 0, h_space = 0;
329
  void* temp = 0;
330
  int i = 0, rank = 0;
331
332
333
334
335
336
337
338
339
  const size_t typeSize = sizeOfType(type);
  const size_t copySize = typeSize * dim;
  const size_t partSize = sizeof(struct part);
  char* temp_c = 0;
  hsize_t shape[2];
  char buffer[150];

  /* Create data space */
  h_space = H5Screate(H5S_SIMPLE);
340
341
342
343
344
345
346
347
348
349
350
351
352
353
  if (h_space < 0) {
    error("Error while creating data space for field '%s'.", name);
  }

  if (dim > 1) {
    rank = 2;
    shape[0] = N_total;
    shape[1] = dim;
  } else {
    rank = 1;
    shape[0] = N_total;
    shape[1] = 0;
  }

354
355
  /* Change shape of data space */
  h_err = H5Sset_extent_simple(h_space, rank, shape, NULL);
356
357
358
359
  if (h_err < 0) {
    error("Error while changing data space shape for field '%s'.", name);
  }

360
361
  /* Create dataset */
  h_data = H5Dcreate1(grp, name, hdf5Type(type), h_space, H5P_DEFAULT);
362
363
364
  if (h_data < 0) {
    error("Error while creating dataspace '%s'.", name);
  }
365
366
367
368
369

  /* Write XMF description for this data set */
  writeXMFline(xmfFile, fileName, name, N_total, dim, type);

  /* Write unit conversion factors for this data set */
370
371
372
373
374
375
  conversionString(buffer, us, convFactor);
  writeAttribute_d(h_data, "CGS conversion factor",
                   conversionFactor(us, convFactor));
  writeAttribute_f(h_data, "h-scale exponant", hFactor(us, convFactor));
  writeAttribute_f(h_data, "a-scale exponant", aFactor(us, convFactor));
  writeAttribute_s(h_data, "Conversion factor", buffer);
376
377
378
379
380

  H5Dclose(h_data);
  H5Sclose(h_space);
}

381
382
383
384
/**
 * @brief Writes a data array in given HDF5 group.
 *
 * @param grp The group in which to write.
385
386
 * @param fileName The name of the file in which the data is written
 * @param xmfFile The FILE used to write the XMF description
387
388
389
390
 * @param name The name of the array to write.
 * @param type The #DATA_TYPE of the array.
 * @param N The number of particles to write.
 * @param dim The dimension of the data (1 for scalar, 3 for vector)
391
392
 * @param part_c A (char*) pointer on the first occurence of the field of
 *interest in the parts array
393
394
 * @param us The UnitSystem currently in use
 * @param convFactor The UnitConversionFactor for this array
395
396
397
398
 *
 *
 * Calls #error() if an error occurs.
 */
399
400
401
402
void writeArrayBackEnd(hid_t grp, char* name, enum DATA_TYPE type, int N,
                       int dim, long long N_total, long long offset,
                       char* part_c) {
  hid_t h_data = 0, h_err = 0, h_memspace = 0, h_filespace = 0;
403
  hsize_t shape[2], shape_total[2], offsets[2];
404
  void* temp = 0;
405
  int i = 0, rank = 0;
406
407
408
409
410
  const size_t typeSize = sizeOfType(type);
  const size_t copySize = typeSize * dim;
  const size_t partSize = sizeof(struct part);
  char* temp_c = 0;

411
  /* message("Writing '%s' array...", name); */
412
413
414

  /* Allocate temporary buffer */
  temp = malloc(N * dim * sizeOfType(type));
415
  if (temp == NULL) error("Unable to allocate memory for temporary buffer");
416

417
  /* Copy particle data to temporary buffer */
418
  temp_c = temp;
419
420
  for (i = 0; i < N; ++i)
    memcpy(&temp_c[i * copySize], part_c + i * partSize, copySize);
421

422
  /* Construct information for the hyperslab */
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
  if (dim > 1) {
    rank = 2;
    shape[0] = N;
    shape[1] = dim;
    shape_total[0] = N_total;
    shape_total[1] = dim;
    offsets[0] = offset;
    offsets[1] = 0;
  } else {
    rank = 1;
    shape[0] = N;
    shape[1] = 0;
    shape_total[0] = N_total;
    shape_total[1] = 0;
    offsets[0] = offset;
    offsets[1] = 0;
  }
440
441
442

  /* Create data space in memory */
  h_memspace = H5Screate(H5S_SIMPLE);
443
444
  if (h_memspace < 0)
    error("Error while creating data space (memory) for field '%s'.", name);
445
446
447

  /* Change shape of memory data space */
  h_err = H5Sset_extent_simple(h_memspace, rank, shape, NULL);
448
449
450
451
  if (h_err < 0)
    error("Error while changing data space (memory) shape for field '%s'.",
          name);

452
453
  /* Open pre-existing data set */
  h_data = H5Dopen(grp, name, H5P_DEFAULT);
454
  if (h_data < 0) error("Error while opening dataset '%s'.", name);
455

456
457
458
  /* Select data space in that data set */
  h_filespace = H5Dget_space(h_data);
  H5Sselect_hyperslab(h_filespace, H5S_SELECT_SET, offsets, NULL, shape, NULL);
459

460
  /* Write temporary buffer to HDF5 dataspace */
461
462
463
464
  h_err = H5Dwrite(h_data, hdf5Type(type), h_memspace, h_filespace, H5P_DEFAULT,
                   temp);
  if (h_err < 0) error("Error while writing data array '%s'.", name);

465
466
467
  /* Free and close everything */
  free(temp);
  H5Dclose(h_data);
468
469
  H5Sclose(h_memspace);
  H5Sclose(h_filespace);
470
471
472
473
474
475
}

/**
 * @brief A helper macro to call the readArrayBackEnd function more easily.
 *
 * @param grp The group in which to write.
476
477
 * @param fileName The name of the file in which the data is written
 * @param xmfFile The FILE used to write the XMF description
478
479
480
481
 * @param name The name of the array to write.
 * @param type The #DATA_TYPE of the array.
 * @param N The number of particles to write.
 * @param dim The dimension of the data (1 for scalar, 3 for vector)
482
483
 * @param part A (char*) pointer on the first occurence of the field of interest
 *in the parts array
484
 * @param field The name (code name) of the field to read from.
485
486
 * @param us The UnitSystem currently in use
 * @param convFactor The UnitConversionFactor for this array
487
488
 *
 */
489
490
491
#define writeArray(grp, name, type, N, dim, N_total, offset, part, field) \
  writeArrayBackEnd(grp, name, type, N, dim, N_total, offset,             \
                    (char*)(&(part[0]).field))
492
493

/**
494
 * @brief Writes an HDF5 output file (GADGET-3 type) with its XMF descriptor
495
 *
496
 * @param e The engine containing all the system.
497
 * @param us The UnitSystem used for the conversion of units in the output
498
 *
499
 * Creates an HDF5 output file and writes the particles contained
500
 * in the engine. If such a file already exists, it is erased and replaced
501
 * by the new one.
502
 * The companion XMF file is also updated accordingly.
503
504
505
506
 *
 * Calls #error() if an error occurs.
 *
 */
507
508
509
void write_output_serial(struct engine* e, struct UnitSystem* us, int mpi_rank,
                         int mpi_size, MPI_Comm comm, MPI_Info info) {
  hid_t h_file = 0, h_grp = 0;
510
511
  int N = e->s->nr_parts;
  int periodic = e->s->periodic;
512
513
514
  int numParticles[6] = {N, 0};
  int numParticlesHighWord[6] = {0};
  unsigned int flagEntropy[6] = {0};
515
516
  long long N_total = 0, offset = 0;
  double offset_d = 0., N_d = 0., N_total_d = 0.;
517
  int numFiles = 1;
518
  int rank = 0;
519
  struct part* parts = e->s->parts;
520
521
  FILE* xmfFile = 0;
  static int outputCount = 0;
522

523
524
525
526
  /* File name */
  char fileName[200];
  sprintf(fileName, "output_%03i.hdf5", outputCount);

527
528
529
530
531
  /* Compute offset in the file and total number of particles */
  /* Done using double to allow for up to 2^50=10^15 particles */
  N_d = (double)N;
  MPI_Exscan(&N_d, &offset_d, 1, MPI_DOUBLE, MPI_SUM, comm);
  N_total_d = offset_d + N_d;
532
533
534
535
536
537
538
  MPI_Bcast(&N_total_d, 1, MPI_DOUBLE, mpi_size - 1, comm);
  if (N_total_d > 1.e15)
    error(
        "Error while computing the offest for parallel output: Simulation has "
        "more than 10^15 particles.\n");
  N_total = (long long)N_total_d;
  offset = (long long)offset_d;
539

540
  /* Do common stuff first */
541
  if (mpi_rank == 0) {
542

543
    /* First time, we need to create the XMF file */
544
545
    if (outputCount == 0) createXMFfile();

546
547
    /* Prepare the XMF file for the new entry */
    xmfFile = prepareXMFfile();
548

549
550
    /* Write the part corresponding to this specific output */
    writeXMFheader(xmfFile, N_total, fileName, e->time);
551

552
553
    /* Open file */
    /* message("Opening file '%s'.", fileName); */
554
555
556
557
    h_file = H5Fcreate(fileName, H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
    if (h_file < 0) {
      error("Error while opening file '%s'.", fileName);
    }
558

559
560
561
    /* Open header to write simulation properties */
    /* message("Writing runtime parameters..."); */
    h_grp = H5Gcreate1(h_file, "/RuntimePars", 0);
562
    if (h_grp < 0) error("Error while creating runtime parameters group\n");
563

564
565
566
567
568
    /* Write the relevant information */
    writeAttribute(h_grp, "PeriodicBoundariesOn", INT, &periodic, 1);

    /* Close runtime parameters */
    H5Gclose(h_grp);
569

570
571
572
    /* Open header to write simulation properties */
    /* message("Writing file header..."); */
    h_grp = H5Gcreate1(h_file, "/Header", 0);
573
574
    if (h_grp < 0) error("Error while creating file header\n");

575
576
577
578
579
    /* Print the relevant information and print status */
    writeAttribute(h_grp, "BoxSize", DOUBLE, e->s->dim, 3);
    writeAttribute(h_grp, "Time", DOUBLE, &e->time, 1);

    /* GADGET-2 legacy values */
580
    numParticles[0] = (unsigned int)N_total;
581
582
    writeAttribute(h_grp, "NumPart_ThisFile", UINT, numParticles, 6);
    writeAttribute(h_grp, "NumPart_Total", UINT, numParticles, 6);
583
584
585
    numParticlesHighWord[0] = (unsigned int)(N_total >> 32);
    writeAttribute(h_grp, "NumPart_Total_HighWord", UINT, numParticlesHighWord,
                   6);
586
587
588
589
590
591
592
593
594
595
596
597
598
    double MassTable[6] = {0., 0., 0., 0., 0., 0.};
    writeAttribute(h_grp, "MassTable", DOUBLE, MassTable, 6);
    writeAttribute(h_grp, "Flag_Entropy_ICs", UINT, flagEntropy, 6);
    writeAttribute(h_grp, "NumFilesPerSnapshot", INT, &numFiles, 1);

    /* Close header */
    H5Gclose(h_grp);

    /* Print the SPH parameters */
    writeSPHflavour(h_file);

    /* Print the system of Units */
    writeUnitSystem(h_file, us);
599

600
601
602
    /* Create SPH particles group */
    /* message("Writing particle arrays..."); */
    h_grp = H5Gcreate1(h_file, "/PartType0", 0);
603
    if (h_grp < 0) error("Error while creating particle group.\n");
604
605

    /* Prepare the arrays in the file */
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
    prepareArray(h_grp, fileName, xmfFile, "Coordinates", DOUBLE, N_total, 3,
                 us, UNIT_CONV_LENGTH);
    prepareArray(h_grp, fileName, xmfFile, "Velocities", FLOAT, N_total, 3, us,
                 UNIT_CONV_SPEED);
    prepareArray(h_grp, fileName, xmfFile, "Masses", FLOAT, N_total, 1, us,
                 UNIT_CONV_MASS);
    prepareArray(h_grp, fileName, xmfFile, "SmoothingLength", FLOAT, N_total, 1,
                 us, UNIT_CONV_LENGTH);
    prepareArray(h_grp, fileName, xmfFile, "InternalEnergy", FLOAT, N_total, 1,
                 us, UNIT_CONV_ENERGY_PER_UNIT_MASS);
    prepareArray(h_grp, fileName, xmfFile, "ParticleIDs", ULONGLONG, N_total, 1,
                 us, UNIT_CONV_NO_UNITS);
    prepareArray(h_grp, fileName, xmfFile, "TimeStep", FLOAT, N_total, 1, us,
                 UNIT_CONV_TIME);
    prepareArray(h_grp, fileName, xmfFile, "Acceleration", FLOAT, N_total, 3,
                 us, UNIT_CONV_ACCELERATION);
    prepareArray(h_grp, fileName, xmfFile, "Density", FLOAT, N_total, 1, us,
                 UNIT_CONV_DENSITY);

625
626
    /* Close particle group */
    H5Gclose(h_grp);
627

628
629
630
631
632
633
    /* Close file */
    H5Fclose(h_file);

    /* Write footer of LXMF file descriptor */
    writeXMFfooter(xmfFile);
  }
634

635
  /* Now loop over ranks and write the data */
636
  for (rank = 0; rank < mpi_size; ++rank) {
637
638

    /* Is it this rank's turn to write ? */
639
    if (rank == mpi_rank) {
640
641

      h_file = H5Fopen(fileName, H5F_ACC_RDWR, H5P_DEFAULT);
642
643
      if (h_file < 0)
        error("Error while opening file '%s' on rank %d.", fileName, mpi_rank);
644
645
646
647

      /* Open SPH particles group */
      /* message("Reading particle arrays..."); */
      h_grp = H5Gopen1(h_file, "/PartType0");
648
649
      if (h_grp < 0)
        error("Error while opening particle group on rank %d.\n", mpi_rank);
650
651
652
653
654

      /* Write arrays */
      writeArray(h_grp, "Coordinates", DOUBLE, N, 3, N_total, offset, parts, x);
      writeArray(h_grp, "Velocities", FLOAT, N, 3, N_total, offset, parts, v);
      writeArray(h_grp, "Masses", FLOAT, N, 1, N_total, offset, parts, mass);
655
656
657
658
659
660
      writeArray(h_grp, "SmoothingLength", FLOAT, N, 1, N_total, offset, parts,
                 h);
      writeArray(h_grp, "InternalEnergy", FLOAT, N, 1, N_total, offset, parts,
                 u);
      writeArray(h_grp, "ParticleIDs", ULONGLONG, N, 1, N_total, offset, parts,
                 id);
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
      writeArray(h_grp, "TimeStep", FLOAT, N, 1, N_total, offset, parts, dt);
      writeArray(h_grp, "Acceleration", FLOAT, N, 3, N_total, offset, parts, a);
      writeArray(h_grp, "Density", FLOAT, N, 1, N_total, offset, parts, rho);

      /* Close particle group */
      H5Gclose(h_grp);

      /* Close file */
      H5Fclose(h_file);
    }

    /* Wait for the read of the reading to complete */
    MPI_Barrier(comm);
  }

  /* message("Done writing particles..."); */
677
  ++outputCount;
678
679
}

680
#endif /* HAVE_HDF5 */