distributed_io.c 31.7 KB
Newer Older
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 26 27 28 29 30 31
/*******************************************************************************
 * This file is part of SWIFT.
 * Copyright (c) 2019 Matthieu Schaller (schaller@strw.leidenuniv.nl)
 *
 * 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.
 *
 * 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.
 *
 * 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/>.
 *
 ******************************************************************************/

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

#if defined(HAVE_HDF5) && defined(WITH_MPI)

/* Some standard headers. */
#include <hdf5.h>
#include <math.h>
#include <stddef.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
32
#include <sys/stat.h>
33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52
#include <time.h>

/* This object's header. */
#include "distributed_io.h"

/* Local includes. */
#include "black_holes_io.h"
#include "chemistry_io.h"
#include "common_io.h"
#include "cooling_io.h"
#include "dimension.h"
#include "engine.h"
#include "error.h"
#include "fof_io.h"
#include "gravity_io.h"
#include "gravity_properties.h"
#include "hydro_io.h"
#include "hydro_properties.h"
#include "io_properties.h"
#include "memuse.h"
53 54
#include "output_list.h"
#include "output_options.h"
55 56
#include "part.h"
#include "part_type.h"
Loic Hausammann's avatar
Loic Hausammann committed
57
#include "sink_io.h"
58 59
#include "star_formation_io.h"
#include "stars_io.h"
60
#include "tools.h"
61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81
#include "tracers_io.h"
#include "units.h"
#include "velociraptor_io.h"
#include "xmf.h"

/**
 * @brief Writes a data array in given HDF5 group.
 *
 * @param e The #engine we are writing from.
 * @param grp The group in which to write.
 * @param fileName The name of the file in which the data is written
 * @param partTypeGroupName The name of the group containing the particles in
 * the HDF5 file.
 * @param props The #io_props of the field to read
 * @param N The number of particles to write.
 * @param internal_units The #unit_system used internally
 * @param snapshot_units The #unit_system used in the snapshots
 *
 * @todo A better version using HDF5 hyper-slabs to write the file directly from
 * the part array will be written once the structures have been stabilized.
 */
Matthieu Schaller's avatar
Matthieu Schaller committed
82 83 84 85 86 87
void write_distributed_array(
    const struct engine* e, hid_t grp, const char* fileName,
    const char* partTypeGroupName, const struct io_props props, const size_t N,
    const enum lossy_compression_schemes lossy_compression,
    const struct unit_system* internal_units,
    const struct unit_system* snapshot_units) {
88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103

  const size_t typeSize = io_sizeof_type(props.type);
  const size_t num_elements = N * props.dimension;

  /* message("Writing '%s' array...", props.name); */

  /* Allocate temporary buffer */
  void* temp = NULL;
  if (swift_memalign("writebuff", (void**)&temp, IO_BUFFER_ALIGNMENT,
                     num_elements * typeSize) != 0)
    error("Unable to allocate temporary i/o buffer");

  /* Copy the particle data to the temporary buffer */
  io_copy_temp_buffer(temp, e, props, N, internal_units, snapshot_units);

  /* Create data space */
104 105 106 107 108 109
  hid_t h_space;
  if (N > 0)
    h_space = H5Screate(H5S_SIMPLE);
  else
    h_space = H5Screate(H5S_NULL);

110 111 112
  if (h_space < 0)
    error("Error while creating data space for field '%s'.", props.name);

113
  /* Decide what chunk size to use based on compression */
Matthieu Schaller's avatar
Matthieu Schaller committed
114
  int log2_chunk_size = 20;
115

116 117 118 119 120 121 122 123
  int rank;
  hsize_t shape[2];
  hsize_t chunk_shape[2];

  if (props.dimension > 1) {
    rank = 2;
    shape[0] = N;
    shape[1] = props.dimension;
124
    chunk_shape[0] = 1 << log2_chunk_size;
125 126 127 128 129
    chunk_shape[1] = props.dimension;
  } else {
    rank = 1;
    shape[0] = N;
    shape[1] = 0;
130
    chunk_shape[0] = 1 << log2_chunk_size;
131 132 133 134 135 136 137 138 139 140 141
    chunk_shape[1] = 0;
  }

  /* Make sure the chunks are not larger than the dataset */
  if (chunk_shape[0] > N) chunk_shape[0] = N;

  /* Change shape of data space */
  hid_t h_err = H5Sset_extent_simple(h_space, rank, shape, shape);
  if (h_err < 0)
    error("Error while changing data space shape for field '%s'.", props.name);

Matthieu Schaller's avatar
Matthieu Schaller committed
142 143 144
  /* Dataset type */
  hid_t h_type = H5Tcopy(io_hdf5_type(props.type));

145
  /* Dataset properties */
Matthieu Schaller's avatar
Matthieu Schaller committed
146
  hid_t h_prop = H5Pcreate(H5P_DATASET_CREATE);
147

148 149
  /* Create filters and set compression level if we have something to write */
  if (N > 0) {
150

151 152
    /* Set chunk size */
    h_err = H5Pset_chunk(h_prop, rank, chunk_shape);
153
    if (h_err < 0)
154 155
      error("Error while setting chunk size (%llu, %llu) for field '%s'.",
            chunk_shape[0], chunk_shape[1], props.name);
156

Matthieu Schaller's avatar
Matthieu Schaller committed
157 158 159 160
    /* Are we imposing some form of lossy compression filter? */
    if (lossy_compression != compression_write_lossless)
      set_hdf5_lossy_compression(&h_prop, &h_type, lossy_compression,
                                 props.name);
161

Matthieu Schaller's avatar
Matthieu Schaller committed
162
    /* Impose GZIP data compression */
163 164 165 166 167 168 169 170 171 172 173
    if (e->snapshot_compression > 0) {
      h_err = H5Pset_shuffle(h_prop);
      if (h_err < 0)
        error("Error while setting shuffling options for field '%s'.",
              props.name);

      h_err = H5Pset_deflate(h_prop, e->snapshot_compression);
      if (h_err < 0)
        error("Error while setting compression options for field '%s'.",
              props.name);
    }
Matthieu Schaller's avatar
Matthieu Schaller committed
174 175 176 177 178

    /* Impose check-sum to verify data corruption */
    h_err = H5Pset_fletcher32(h_prop);
    if (h_err < 0)
      error("Error while setting checksum options for field '%s'.", props.name);
179 180 181
  }

  /* Create dataset */
Matthieu Schaller's avatar
Matthieu Schaller committed
182 183
  const hid_t h_data = H5Dcreate(grp, props.name, h_type, h_space, H5P_DEFAULT,
                                 h_prop, H5P_DEFAULT);
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
  if (h_data < 0) error("Error while creating dataspace '%s'.", props.name);

  /* Write temporary buffer to HDF5 dataspace */
  h_err = H5Dwrite(h_data, io_hdf5_type(props.type), h_space, H5S_ALL,
                   H5P_DEFAULT, temp);
  if (h_err < 0) error("Error while writing data array '%s'.", props.name);

  /* Write unit conversion factors for this data set */
  char buffer[FIELD_BUFFER_SIZE] = {0};
  units_cgs_conversion_string(buffer, snapshot_units, props.units,
                              props.scale_factor_exponent);
  float baseUnitsExp[5];
  units_get_base_unit_exponents_array(baseUnitsExp, props.units);
  io_write_attribute_f(h_data, "U_M exponent", baseUnitsExp[UNIT_MASS]);
  io_write_attribute_f(h_data, "U_L exponent", baseUnitsExp[UNIT_LENGTH]);
  io_write_attribute_f(h_data, "U_t exponent", baseUnitsExp[UNIT_TIME]);
  io_write_attribute_f(h_data, "U_I exponent", baseUnitsExp[UNIT_CURRENT]);
  io_write_attribute_f(h_data, "U_T exponent", baseUnitsExp[UNIT_TEMPERATURE]);
  io_write_attribute_f(h_data, "h-scale exponent", 0.f);
  io_write_attribute_f(h_data, "a-scale exponent", props.scale_factor_exponent);
  io_write_attribute_s(h_data, "Expression for physical CGS units", buffer);

  /* Write the actual number this conversion factor corresponds to */
  const double factor =
      units_cgs_conversion_factor(snapshot_units, props.units);
  io_write_attribute_d(
      h_data,
      "Conversion factor to CGS (not including cosmological corrections)",
      factor);
  io_write_attribute_d(
      h_data,
      "Conversion factor to physical CGS (including cosmological corrections)",
      factor * pow(e->cosmology->a, props.scale_factor_exponent));

#ifdef SWIFT_DEBUG_CHECKS
  if (strlen(props.description) == 0)
    error("Invalid (empty) description of the field '%s'", props.name);
#endif

  /* Write the full description */
  io_write_attribute_s(h_data, "Description", props.description);

  /* Free and close everything */
  swift_free("writebuff", temp);
Matthieu Schaller's avatar
Matthieu Schaller committed
228
  H5Tclose(h_type);
229 230 231 232 233 234 235 236 237 238 239
  H5Pclose(h_prop);
  H5Dclose(h_data);
  H5Sclose(h_space);
}

/**
 * @brief Writes a snapshot distributed into multiple files.
 *
 * @param e The engine containing all the system.
 * @param internal_units The #unit_system used internally
 * @param snapshot_units The #unit_system used in the snapshots
240 241 242 243
 * @param mpi_rank The rank number of the calling MPI rank.
 * @param mpi_size the number of MPI ranks.
 * @param comm The communicator used by the MPI ranks.
 * @param info The MPI information object.
244 245 246 247 248 249
 *
 * Creates a series of HDF5 output files (1 per MPI node) as a snapshot.
 * Writes the particles contained in the engine.
 * If such files already exist, it is erased and replaced by the new one.
 * The companion XMF file is also updated accordingly.
 */
250
void write_output_distributed(struct engine* e,
251 252 253 254 255 256 257 258 259 260
                              const struct unit_system* internal_units,
                              const struct unit_system* snapshot_units,
                              const int mpi_rank, const int mpi_size,
                              MPI_Comm comm, MPI_Info info) {

  hid_t h_file = 0, h_grp = 0;
  int numFiles = mpi_size;
  const struct part* parts = e->s->parts;
  const struct xpart* xparts = e->s->xparts;
  const struct gpart* gparts = e->s->gparts;
Loic Hausammann's avatar
Loic Hausammann committed
261
  const struct sink* sinks = e->s->sinks;
262 263
  const struct spart* sparts = e->s->sparts;
  const struct bpart* bparts = e->s->bparts;
264 265
  struct output_options* output_options = e->output_options;
  struct output_list* output_list = e->output_list_snapshots;
266 267 268 269 270 271 272 273 274 275 276 277 278 279 280
  const int with_cosmology = e->policy & engine_policy_cosmology;
  const int with_cooling = e->policy & engine_policy_cooling;
  const int with_temperature = e->policy & engine_policy_temperature;
  const int with_fof = e->policy & engine_policy_fof;
  const int with_DM_background = e->s->with_DM_background;
#ifdef HAVE_VELOCIRAPTOR
  const int with_stf = (e->policy & engine_policy_structure_finding) &&
                       (e->s->gpart_group_data != NULL);
#else
  const int with_stf = 0;
#endif

  /* Number of particles currently in the arrays */
  const size_t Ntot = e->s->nr_gparts;
  const size_t Ngas = e->s->nr_parts;
Loic Hausammann's avatar
Loic Hausammann committed
281
  const size_t Nsinks = e->s->nr_sinks;
282 283 284 285 286 287 288 289 290 291 292 293 294 295
  const size_t Nstars = e->s->nr_sparts;
  const size_t Nblackholes = e->s->nr_bparts;

  size_t Ndm_background = 0;
  if (with_DM_background) {
    Ndm_background = io_count_dm_background_gparts(gparts, Ntot);
  }

  /* Number of particles that we will write in this file.
   * Recall that background particles are never inhibited and have no extras */
  const size_t Ntot_written =
      e->s->nr_gparts - e->s->nr_inhibited_gparts - e->s->nr_extra_gparts;
  const size_t Ngas_written =
      e->s->nr_parts - e->s->nr_inhibited_parts - e->s->nr_extra_parts;
Loic Hausammann's avatar
Loic Hausammann committed
296 297
  const size_t Nsinks_written =
      e->s->nr_sinks - e->s->nr_inhibited_sinks - e->s->nr_extra_sinks;
298 299 300 301 302
  const size_t Nstars_written =
      e->s->nr_sparts - e->s->nr_inhibited_sparts - e->s->nr_extra_sparts;
  const size_t Nblackholes_written =
      e->s->nr_bparts - e->s->nr_inhibited_bparts - e->s->nr_extra_bparts;
  const size_t Nbaryons_written =
Loic Hausammann's avatar
Loic Hausammann committed
303
      Ngas_written + Nstars_written + Nblackholes_written + Nsinks_written;
304 305 306
  const size_t Ndm_written =
      Ntot_written > 0 ? Ntot_written - Nbaryons_written - Ndm_background : 0;

307
  int snap_count = -1;
308
  if (e->snapshot_int_time_label_on)
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342
    snap_count = (int)round(e->time);
  else if (e->snapshot_invoke_stf)
    snap_count = e->stf_output_count;
  else
    snap_count = e->snapshot_output_count;

  int number_digits = -1;
  if (e->snapshot_int_time_label_on)
    number_digits = 6;
  else
    number_digits = 4;

  /* Directory and file name */
  char dirName[1024];
  char fileName[1024];

  if (strnlen(e->snapshot_subdir, PARSER_MAX_LINE_SIZE) > 0) {
    sprintf(dirName, "%s/%s_%0*d", e->snapshot_subdir, e->snapshot_base_name,
            number_digits, snap_count);

    sprintf(fileName, "%s/%s_%0*d/%s_%0*d.%d.hdf5", e->snapshot_subdir,
            e->snapshot_base_name, number_digits, snap_count,
            e->snapshot_base_name, number_digits, snap_count, mpi_rank);

  } else {
    sprintf(dirName, "%s_%0*d", e->snapshot_base_name, number_digits,
            snap_count);

    sprintf(fileName, "%s_%0*d/%s_%0*d.%d.hdf5", e->snapshot_base_name,
            number_digits, snap_count, e->snapshot_base_name, number_digits,
            snap_count, mpi_rank);
  }

  /* Create the directory */
343 344
  if (mpi_rank == 0) safe_checkdir(dirName, /*create=*/1);
  MPI_Barrier(comm);
345 346 347

  /* Compute offset in the file and total number of particles */
  const long long N[swift_type_count] = {Ngas_written,   Ndm_written,
Loic Hausammann's avatar
Loic Hausammann committed
348
                                         Ndm_background, Nsinks_written,
349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374
                                         Nstars_written, Nblackholes_written};

  /* Gather the total number of particles to write */
  long long N_total[swift_type_count] = {0};
  MPI_Allreduce(N, N_total, swift_type_count, MPI_LONG_LONG_INT, MPI_SUM, comm);

  /* Open file */
  /* message("Opening file '%s'.", fileName); */
  h_file = H5Fcreate(fileName, H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
  if (h_file < 0) error("Error while opening file '%s'.", fileName);

  /* Open header to write simulation properties */
  /* message("Writing file header..."); */
  h_grp = H5Gcreate(h_file, "/Header", H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
  if (h_grp < 0) error("Error while creating file header\n");

  /* Convert basic output information to snapshot units */
  const double factor_time =
      units_conversion_factor(internal_units, snapshot_units, UNIT_CONV_TIME);
  const double factor_length =
      units_conversion_factor(internal_units, snapshot_units, UNIT_CONV_LENGTH);
  const double dblTime = e->time * factor_time;
  const double dim[3] = {e->s->dim[0] * factor_length,
                         e->s->dim[1] * factor_length,
                         e->s->dim[2] * factor_length};

375 376
  /* Determine if we are writing a reduced snapshot, and if so which
   * output selection type to use */
377 378
  char current_selection_name[FIELD_BUFFER_SIZE] =
      select_output_header_default_name;
379 380 381
  if (output_list) {
    /* Users could have specified a different Select Output scheme for each
     * snapshot. */
382
    output_list_get_current_select_output(output_list, current_selection_name);
383 384
  }

385 386 387 388 389 390 391 392 393 394
  /* Print the relevant information and print status */
  io_write_attribute(h_grp, "BoxSize", DOUBLE, dim, 3);
  io_write_attribute(h_grp, "Time", DOUBLE, &dblTime, 1);
  const int dimension = (int)hydro_dimension;
  io_write_attribute(h_grp, "Dimension", INT, &dimension, 1);
  io_write_attribute(h_grp, "Redshift", DOUBLE, &e->cosmology->z, 1);
  io_write_attribute(h_grp, "Scale-factor", DOUBLE, &e->cosmology->a, 1);
  io_write_attribute_s(h_grp, "Code", "SWIFT");
  io_write_attribute_s(h_grp, "RunName", e->run_name);

395 396 397 398 399 400 401
  /* Store the time at which the snapshot was written */
  time_t tm = time(NULL);
  struct tm* timeinfo = localtime(&tm);
  char snapshot_date[64];
  strftime(snapshot_date, 64, "%T %F %Z", timeinfo);
  io_write_attribute_s(h_grp, "Snapshot date", snapshot_date);

402
  /* GADGET-2 legacy values:  Number of particles of each type */
403 404
  unsigned int numParticles[swift_type_count] = {0};
  unsigned int numParticlesHighWord[swift_type_count] = {0};
405 406 407 408

  /* Total number of fields to write per ptype */
  int numFields[swift_type_count] = {0};

409 410 411
  for (int ptype = 0; ptype < swift_type_count; ++ptype) {
    numParticles[ptype] = (unsigned int)N_total[ptype];
    numParticlesHighWord[ptype] = (unsigned int)(N_total[ptype] >> 32);
412 413 414

    numFields[ptype] = output_options_get_num_fields_to_write(
        output_options, current_selection_name, ptype);
415
  }
416

417 418 419 420 421 422 423 424 425 426 427 428 429
  io_write_attribute(h_grp, "NumPart_ThisFile", LONGLONG, N, swift_type_count);
  io_write_attribute(h_grp, "NumPart_Total", UINT, numParticles,
                     swift_type_count);
  io_write_attribute(h_grp, "NumPart_Total_HighWord", UINT,
                     numParticlesHighWord, swift_type_count);
  double MassTable[swift_type_count] = {0};
  io_write_attribute(h_grp, "MassTable", DOUBLE, MassTable, swift_type_count);
  unsigned int flagEntropy[swift_type_count] = {0};
  flagEntropy[0] = writeEntropyFlag();
  io_write_attribute(h_grp, "Flag_Entropy_ICs", UINT, flagEntropy,
                     swift_type_count);
  io_write_attribute_i(h_grp, "NumFilesPerSnapshot", numFiles);
  io_write_attribute_i(h_grp, "ThisFile", mpi_rank);
430
  io_write_attribute_s(h_grp, "OutputType", "FullVolume");
431
  io_write_attribute_s(h_grp, "SelectOutput", current_selection_name);
432 433 434 435

  /* Close header */
  H5Gclose(h_grp);

436 437
  /* Write all the meta-data */
  io_write_meta_data(h_file, e, internal_units, snapshot_units);
438

439 440 441 442
  /* Now write the top-level cell structure
   * We use a global offset of 0 here. This means that the cells will write
   * their offset with respect to the start of the file they belong to and
   * not a global offset */
443 444 445 446 447 448
  long long global_offsets[swift_type_count] = {0};
  h_grp = H5Gcreate(h_file, "/Cells", H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
  if (h_grp < 0) error("Error while creating cells group");

  /* Write the location of the particles in the arrays */
  io_write_cell_offsets(h_grp, e->s->cdim, e->s->dim, e->s->pos_dithering,
449
                        e->s->cells_top, e->s->nr_cells, e->s->width, mpi_rank,
450
                        /*distributed=*/1, N_total, global_offsets, numFields,
451
                        internal_units, snapshot_units);
452 453 454 455 456
  H5Gclose(h_grp);

  /* Loop over all particle types */
  for (int ptype = 0; ptype < swift_type_count; ptype++) {

457 458 459
    /* Don't do anything if there are (a) no particles of this kind, or (b)
     * if we have disabled every field of this particle type. */
    if (numParticles[ptype] == 0 || numFields[ptype] == 0) continue;
460 461 462 463 464 465 466 467 468

    /* Open the particle group in the file */
    char partTypeGroupName[PARTICLE_GROUP_BUFFER_SIZE];
    snprintf(partTypeGroupName, PARTICLE_GROUP_BUFFER_SIZE, "/PartType%d",
             ptype);
    h_grp = H5Gcreate(h_file, partTypeGroupName, H5P_DEFAULT, H5P_DEFAULT,
                      H5P_DEFAULT);
    if (h_grp < 0) error("Error while creating particle group.\n");

469 470 471 472 473 474 475 476 477 478 479
    /* Add an alias name for convenience */
    char aliasName[PARTICLE_GROUP_BUFFER_SIZE];
    snprintf(aliasName, PARTICLE_GROUP_BUFFER_SIZE, "/%sParticles",
             part_type_names[ptype]);
    hid_t h_err = H5Lcreate_soft(partTypeGroupName, h_grp, aliasName,
                                 H5P_DEFAULT, H5P_DEFAULT);
    if (h_err < 0) error("Error while creating alias for particle group.\n");

    /* Write the number of particles as an attribute */
    io_write_attribute_l(h_grp, "NumberOfParticles", N[ptype]);

480 481 482 483 484 485 486 487
    int num_fields = 0;
    struct io_props list[100];
    size_t Nparticles = 0;

    struct part* parts_written = NULL;
    struct xpart* xparts_written = NULL;
    struct gpart* gparts_written = NULL;
    struct velociraptor_gpart_data* gpart_group_data_written = NULL;
Loic Hausammann's avatar
Loic Hausammann committed
488
    struct sink* sinks_written = NULL;
489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649
    struct spart* sparts_written = NULL;
    struct bpart* bparts_written = NULL;

    /* Write particle fields from the particle structure */
    switch (ptype) {

      case swift_type_gas: {
        if (Ngas == Ngas_written) {

          /* No inhibted particles: easy case */
          Nparticles = Ngas;
          hydro_write_particles(parts, xparts, list, &num_fields);
          num_fields += chemistry_write_particles(parts, list + num_fields);
          if (with_cooling || with_temperature) {
            num_fields += cooling_write_particles(
                parts, xparts, list + num_fields, e->cooling_func);
          }
          if (with_fof) {
            num_fields += fof_write_parts(parts, xparts, list + num_fields);
          }
          if (with_stf) {
            num_fields +=
                velociraptor_write_parts(parts, xparts, list + num_fields);
          }
          num_fields += tracers_write_particles(
              parts, xparts, list + num_fields, with_cosmology);
          num_fields +=
              star_formation_write_particles(parts, xparts, list + num_fields);

        } else {

          /* Ok, we need to fish out the particles we want */
          Nparticles = Ngas_written;

          /* Allocate temporary arrays */
          if (swift_memalign("parts_written", (void**)&parts_written,
                             part_align,
                             Ngas_written * sizeof(struct part)) != 0)
            error("Error while allocating temporary memory for parts");
          if (swift_memalign("xparts_written", (void**)&xparts_written,
                             xpart_align,
                             Ngas_written * sizeof(struct xpart)) != 0)
            error("Error while allocating temporary memory for xparts");

          /* Collect the particles we want to write */
          io_collect_parts_to_write(parts, xparts, parts_written,
                                    xparts_written, Ngas, Ngas_written);

          /* Select the fields to write */
          hydro_write_particles(parts_written, xparts_written, list,
                                &num_fields);
          num_fields +=
              chemistry_write_particles(parts_written, list + num_fields);
          if (with_cooling || with_temperature) {
            num_fields +=
                cooling_write_particles(parts_written, xparts_written,
                                        list + num_fields, e->cooling_func);
          }
          if (with_fof) {
            num_fields += fof_write_parts(parts_written, xparts_written,
                                          list + num_fields);
          }
          if (with_stf) {
            num_fields += velociraptor_write_parts(
                parts_written, xparts_written, list + num_fields);
          }
          num_fields += tracers_write_particles(
              parts_written, xparts_written, list + num_fields, with_cosmology);
          num_fields += star_formation_write_particles(
              parts_written, xparts_written, list + num_fields);
        }
      } break;

      case swift_type_dark_matter: {
        if (Ntot == Ndm_written) {

          /* This is a DM-only run without background or inhibited particles */
          Nparticles = Ntot;
          darkmatter_write_particles(gparts, list, &num_fields);
          if (with_fof) {
            num_fields += fof_write_gparts(gparts, list + num_fields);
          }
          if (with_stf) {
            num_fields += velociraptor_write_gparts(e->s->gpart_group_data,
                                                    list + num_fields);
          }
        } else {

          /* Ok, we need to fish out the particles we want */
          Nparticles = Ndm_written;

          /* Allocate temporary array */
          if (swift_memalign("gparts_written", (void**)&gparts_written,
                             gpart_align,
                             Ndm_written * sizeof(struct gpart)) != 0)
            error("Error while allocating temporary memory for gparts");

          if (with_stf) {
            if (swift_memalign(
                    "gpart_group_written", (void**)&gpart_group_data_written,
                    gpart_align,
                    Ndm_written * sizeof(struct velociraptor_gpart_data)) != 0)
              error(
                  "Error while allocating temporary memory for gparts STF "
                  "data");
          }

          /* Collect the non-inhibited DM particles from gpart */
          io_collect_gparts_to_write(gparts, e->s->gpart_group_data,
                                     gparts_written, gpart_group_data_written,
                                     Ntot, Ndm_written, with_stf);

          /* Select the fields to write */
          darkmatter_write_particles(gparts_written, list, &num_fields);
          if (with_fof) {
            num_fields += fof_write_gparts(gparts_written, list + num_fields);
          }
          if (with_stf) {
            num_fields += velociraptor_write_gparts(gpart_group_data_written,
                                                    list + num_fields);
          }
        }
      } break;

      case swift_type_dark_matter_background: {

        /* Ok, we need to fish out the particles we want */
        Nparticles = Ndm_background;

        /* Allocate temporary array */
        if (swift_memalign("gparts_written", (void**)&gparts_written,
                           gpart_align,
                           Ndm_background * sizeof(struct gpart)) != 0)
          error("Error while allocating temporart memory for gparts");

        if (with_stf) {
          if (swift_memalign(
                  "gpart_group_written", (void**)&gpart_group_data_written,
                  gpart_align,
                  Ndm_background * sizeof(struct velociraptor_gpart_data)) != 0)
            error(
                "Error while allocating temporart memory for gparts STF "
                "data");
        }

        /* Collect the non-inhibited DM particles from gpart */
        io_collect_gparts_background_to_write(
            gparts, e->s->gpart_group_data, gparts_written,
            gpart_group_data_written, Ntot, Ndm_background, with_stf);

        /* Select the fields to write */
        darkmatter_write_particles(gparts_written, list, &num_fields);
        if (with_fof) {
          num_fields += fof_write_gparts(gparts_written, list + num_fields);
        }
        if (with_stf) {
          num_fields += velociraptor_write_gparts(gpart_group_data_written,
                                                  list + num_fields);
        }
      } break;

Loic Hausammann's avatar
Loic Hausammann committed
650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677
      case swift_type_sink: {
        if (Nsinks == Nsinks_written) {

          /* No inhibted particles: easy case */
          Nparticles = Nsinks;
          sink_write_particles(sinks, list, &num_fields, with_cosmology);

        } else {

          /* Ok, we need to fish out the particles we want */
          Nparticles = Nsinks_written;

          /* Allocate temporary arrays */
          if (swift_memalign("sinks_written", (void**)&sinks_written,
                             sink_align,
                             Nsinks_written * sizeof(struct sink)) != 0)
            error("Error while allocating temporary memory for sinks");

          /* Collect the particles we want to write */
          io_collect_sinks_to_write(sinks, sinks_written, Nsinks,
                                    Nsinks_written);

          /* Select the fields to write */
          sink_write_particles(sinks_written, list, &num_fields,
                               with_cosmology);
        }
      } break;

678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772
      case swift_type_stars: {
        if (Nstars == Nstars_written) {

          /* No inhibted particles: easy case */
          Nparticles = Nstars;
          stars_write_particles(sparts, list, &num_fields, with_cosmology);
          num_fields += chemistry_write_sparticles(sparts, list + num_fields);
          num_fields += tracers_write_sparticles(sparts, list + num_fields,
                                                 with_cosmology);
          if (with_fof) {
            num_fields += fof_write_sparts(sparts, list + num_fields);
          }
          if (with_stf) {
            num_fields += velociraptor_write_sparts(sparts, list + num_fields);
          }
        } else {

          /* Ok, we need to fish out the particles we want */
          Nparticles = Nstars_written;

          /* Allocate temporary arrays */
          if (swift_memalign("sparts_written", (void**)&sparts_written,
                             spart_align,
                             Nstars_written * sizeof(struct spart)) != 0)
            error("Error while allocating temporary memory for sparts");

          /* Collect the particles we want to write */
          io_collect_sparts_to_write(sparts, sparts_written, Nstars,
                                     Nstars_written);

          /* Select the fields to write */
          stars_write_particles(sparts_written, list, &num_fields,
                                with_cosmology);
          num_fields +=
              chemistry_write_sparticles(sparts_written, list + num_fields);
          num_fields += tracers_write_sparticles(
              sparts_written, list + num_fields, with_cosmology);
          if (with_fof) {
            num_fields += fof_write_sparts(sparts_written, list + num_fields);
          }
          if (with_stf) {
            num_fields +=
                velociraptor_write_sparts(sparts_written, list + num_fields);
          }
        }
      } break;

      case swift_type_black_hole: {
        if (Nblackholes == Nblackholes_written) {

          /* No inhibted particles: easy case */
          Nparticles = Nblackholes;
          black_holes_write_particles(bparts, list, &num_fields,
                                      with_cosmology);
          num_fields += chemistry_write_bparticles(bparts, list + num_fields);
          if (with_fof) {
            num_fields += fof_write_bparts(bparts, list + num_fields);
          }
          if (with_stf) {
            num_fields += velociraptor_write_bparts(bparts, list + num_fields);
          }
        } else {

          /* Ok, we need to fish out the particles we want */
          Nparticles = Nblackholes_written;

          /* Allocate temporary arrays */
          if (swift_memalign("bparts_written", (void**)&bparts_written,
                             bpart_align,
                             Nblackholes_written * sizeof(struct bpart)) != 0)
            error("Error while allocating temporary memory for bparts");

          /* Collect the particles we want to write */
          io_collect_bparts_to_write(bparts, bparts_written, Nblackholes,
                                     Nblackholes_written);

          /* Select the fields to write */
          black_holes_write_particles(bparts_written, list, &num_fields,
                                      with_cosmology);
          num_fields +=
              chemistry_write_bparticles(bparts_written, list + num_fields);
          if (with_fof) {
            num_fields += fof_write_bparts(bparts_written, list + num_fields);
          }
          if (with_stf) {
            num_fields +=
                velociraptor_write_bparts(bparts_written, list + num_fields);
          }
        }
      } break;

      default:
        error("Particle Type %d not yet supported. Aborting", ptype);
    }

773 774
    /* Did the user specify a non-standard default for the entire particle
     * type? */
Matthieu Schaller's avatar
Matthieu Schaller committed
775 776 777 778
    const enum lossy_compression_schemes compression_level_current_default =
        output_options_get_ptype_default_compression(
            output_options->select_output, current_selection_name,
            (enum part_type)ptype);
779

780 781
    /* Write everything that is not cancelled */
    int num_fields_written = 0;
782 783 784
    for (int i = 0; i < num_fields; ++i) {

      /* Did the user cancel this field? */
Matthieu Schaller's avatar
Matthieu Schaller committed
785 786 787 788
      const enum lossy_compression_schemes compression_level =
          output_options_get_field_compression(
              output_options, current_selection_name, list[i].name,
              (enum part_type)ptype, compression_level_current_default);
789

Matthieu Schaller's avatar
Matthieu Schaller committed
790
      if (compression_level != compression_do_not_write) {
791
        write_distributed_array(e, h_grp, fileName, partTypeGroupName, list[i],
Matthieu Schaller's avatar
Matthieu Schaller committed
792 793
                                Nparticles, compression_level, internal_units,
                                snapshot_units);
794 795
        num_fields_written++;
      }
796
    }
797 798 799

    /* Only write this now that we know exactly how many fields there are. */
    io_write_attribute_i(h_grp, "NumberOfFields", num_fields_written);
800 801 802 803 804 805 806

    /* Free temporary arrays */
    if (parts_written) swift_free("parts_written", parts_written);
    if (xparts_written) swift_free("xparts_written", xparts_written);
    if (gparts_written) swift_free("gparts_written", gparts_written);
    if (gpart_group_data_written)
      swift_free("gpart_group_written", gpart_group_data_written);
Loic Hausammann's avatar
Loic Hausammann committed
807
    if (sinks_written) swift_free("sinks_written", sinks_written);
808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824
    if (sparts_written) swift_free("sparts_written", sparts_written);
    if (bparts_written) swift_free("bparts_written", bparts_written);

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

  /* message("Done writing particles..."); */

  /* Close file */
  H5Fclose(h_file);

  e->snapshot_output_count++;
  if (e->snapshot_invoke_stf) e->stf_output_count++;
}

#endif /* HAVE_HDF5 && WITH_MPI */