cooling_tables.c 21.2 KB
Newer Older
Matthieu Schaller's avatar
Matthieu Schaller 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
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
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
90
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
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
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
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
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
343
344
345
346
347
348
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
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
/*******************************************************************************
 * This file is part of SWIFT.
 * Copyright (c) 2017 Matthieu Schaller (matthieu.schaller@durham.ac.uk)
 *
 * 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/>.
 *
 ******************************************************************************/

/**
 * @file src/cooling/COLIBRE/cooling_tables.c
 * @brief Functions to read COLIBRE tables
 */

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

/* This file's header */
#include "cooling_tables.h"

/* Standard includes */
#include <hdf5.h>
#include <math.h>
#include <stdlib.h>
#include <string.h>

/* Local includes. */
#include "chemistry_struct.h"
#include "cooling_struct.h"
#include "error.h"
#include "exp10.h"
#include "interpolate.h"

/**
 * @brief Reads in COLIBRE cooling table header. Consists of tables
 * of values for temperature, hydrogen number density, metallicity,
 * abundance ratios, and elements used to index the cooling tables.
 *
 * @param cooling Cooling data structure
 */
void read_cooling_header(struct cooling_function_data *cooling) {

#ifdef HAVE_HDF5

  hid_t dataset;
  herr_t status;

  /* read sizes of array dimensions */
  hid_t tempfile_id =
      H5Fopen(cooling->cooling_table_path, H5F_ACC_RDONLY, H5P_DEFAULT);
  if (tempfile_id < 0)
    error("unable to open file %s\n", cooling->cooling_table_path);

  /* allocate arrays of bins */
  if (posix_memalign((void **)&cooling->Temp, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_temperature * sizeof(float)) != 0)
    error("Failed to allocate temperature table\n");

  if (posix_memalign((void **)&cooling->Redshifts, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_redshifts * sizeof(float)) != 0)
    error("Failed to allocate redshift table\n");

  if (posix_memalign((void **)&cooling->nH, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_density * sizeof(float)) != 0)
    error("Failed to allocate density table\n");

  if (posix_memalign((void **)&cooling->Metallicity, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_metallicity * sizeof(float)) != 0)
    error("Failed to allocate metallicity table\n");

  if (posix_memalign((void **)&cooling->Therm, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_internalenergy * sizeof(float)) != 0)
    error("Failed to allocate internal energy table\n");

  if (posix_memalign((void **)&cooling->LogAbundances, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_metallicity *
                         colibre_cooling_N_elementtypes * sizeof(float)) != 0)
    error("Failed to allocate abundance array\n");

  if (posix_memalign((void **)&cooling->Abundances, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_metallicity *
                         colibre_cooling_N_elementtypes * sizeof(float)) != 0)
    error("Failed to allocate abundance array\n");

  if (posix_memalign((void **)&cooling->Abundances_inv, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_metallicity *
                         colibre_cooling_N_elementtypes * sizeof(float)) != 0)
    error("Failed to allocate abundance array\n");

  if (posix_memalign((void **)&cooling->atomicmass, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_elementtypes * sizeof(float)) != 0)
    error("Failed to allocate atomic masses array\n");

  if (posix_memalign((void **)&cooling->atomicmass_inv, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_elementtypes * sizeof(float)) != 0)
    error("Failed to allocate inverse atomic masses array\n");

  if (posix_memalign((void **)&cooling->Zsol, SWIFT_STRUCT_ALIGNMENT,
                     1 * sizeof(float)) != 0)
    error("Failed to allocate solar metallicity array\n");

  if (posix_memalign((void **)&cooling->Zsol_inv, SWIFT_STRUCT_ALIGNMENT,
                     1 * sizeof(float)) != 0)
    error("Failed to allocate inverse solar metallicity array\n");

  if (posix_memalign((void **)&cooling->LogMassFractions,
                     SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_metallicity *
                         colibre_cooling_N_elementtypes * sizeof(float)) != 0)
    error("Failed to allocate log mass fraction array\n");

  if (posix_memalign((void **)&cooling->MassFractions, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_metallicity *
                         colibre_cooling_N_elementtypes * sizeof(float)) != 0)
    error("Failed to allocate mass fraction array\n");

  /* read in bins and misc information */
  dataset = H5Dopen(tempfile_id, "/TableBins/TemperatureBins", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->Temp);
  if (status < 0) error("error reading temperature bins\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  dataset = H5Dopen(tempfile_id, "/TableBins/RedshiftBins", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->Redshifts);
  if (status < 0) error("error reading redshift bins\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  dataset = H5Dopen(tempfile_id, "/TableBins/DensityBins", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->nH);
  if (status < 0) error("error reading density bins\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  dataset = H5Dopen(tempfile_id, "/TableBins/MetallicityBins", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->Metallicity);
  if (status < 0) error("error reading metallicity bins\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  dataset = H5Dopen(tempfile_id, "/TableBins/InternalEnergyBins", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->Therm);
  if (status < 0) error("error reading internal energy bins\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  dataset = H5Dopen(tempfile_id, "/TotalAbundances", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->LogAbundances);
  if (status < 0) error("error reading total abundances\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  dataset = H5Dopen(tempfile_id, "/TotalMassFractions", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->LogMassFractions);
  if (status < 0) error("error reading total mass fractions\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  dataset = H5Dopen(tempfile_id, "/ElementMasses", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->atomicmass);
  if (status < 0) error("error reading element masses\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  dataset = H5Dopen(tempfile_id, "/SolarMetallicity", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->Zsol);
  if (status < 0) error("error reading solar metallicity \n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  /* Close the file */
  H5Fclose(tempfile_id);

  cooling->Zsol_inv[0] = 1.f / cooling->Zsol[0];

  /* find the metallicity bin that refers to solar metallicity */
  const float tol = 1.e-3;
  for (int i = 0; i < colibre_cooling_N_metallicity; i++) {
    if (fabsf(cooling->Metallicity[i]) < tol) {
      cooling->indxZsol = i;
    }
  }

#if defined(__ICC)
#pragma novector
#endif
  for (int i = 0; i < colibre_cooling_N_elementtypes; i++) {
    cooling->atomicmass_inv[i] = 1.f / cooling->atomicmass[i];
  }

  /* set some additional useful abundance arrays */
  for (int i = 0; i < colibre_cooling_N_metallicity; i++) {

#if defined(__ICC)
#pragma novector
#endif
    for (int j = 0; j < colibre_cooling_N_elementtypes; j++) {
      const int indx1d = row_major_index_2d(i, j, colibre_cooling_N_metallicity,
                                            colibre_cooling_N_elementtypes);
      cooling->Abundances[indx1d] = exp10f(cooling->LogAbundances[indx1d]);
      cooling->Abundances_inv[indx1d] = 1.f / cooling->Abundances[indx1d];
      cooling->MassFractions[indx1d] =
          exp10f(cooling->LogMassFractions[indx1d]);
    }
  }

#else
  error("Need HDF5 to read cooling tables");
#endif
}

/**
 * @brief Allocate space for cooling tables and read them
 *
 * @param cooling #cooling_function_data structure
 */
void read_cooling_tables(struct cooling_function_data *restrict cooling) {

#ifdef HAVE_HDF5
  hid_t dataset;
  herr_t status;

  /* open hdf5 file */
  hid_t tempfile_id =
      H5Fopen(cooling->cooling_table_path, H5F_ACC_RDONLY, H5P_DEFAULT);
  if (tempfile_id < 0)
    error("unable to open file %s\n", cooling->cooling_table_path);

  /* Allocate and read arrays to store cooling tables. */

  /* Mean particle mass (temperature) */
  if (posix_memalign((void **)&cooling->table.Tmu, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_redshifts *
                         colibre_cooling_N_temperature *
                         colibre_cooling_N_metallicity *
                         colibre_cooling_N_density * sizeof(float)) != 0)
    error("Failed to allocate Tmu array\n");

  dataset = H5Dopen(tempfile_id, "/Tdep/MeanParticleMass", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.Tmu);
  if (status < 0) error("error reading Tmu\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing mean particle mass dataset");

  /* Mean particle mass (internal energy) */
  if (posix_memalign((void **)&cooling->table.Umu, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_redshifts *
                         colibre_cooling_N_internalenergy *
                         colibre_cooling_N_metallicity *
                         colibre_cooling_N_density * sizeof(float)) != 0)
    error("Failed to allocate Umu array\n");

  dataset = H5Dopen(tempfile_id, "/Udep/MeanParticleMass", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.Umu);
  if (status < 0) error("error reading Umu\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing mean particle mass dataset");

  /* Cooling (temperature) */
  if (posix_memalign(
          (void **)&cooling->table.Tcooling, SWIFT_STRUCT_ALIGNMENT,
          colibre_cooling_N_redshifts * colibre_cooling_N_temperature *
              colibre_cooling_N_metallicity * colibre_cooling_N_density *
              colibre_cooling_N_cooltypes * sizeof(float)) != 0)
    error("Failed to allocate Tcooling array\n");

  dataset = H5Dopen(tempfile_id, "/Tdep/Cooling", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.Tcooling);
  if (status < 0) error("error reading Tcooling\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  /* Cooling (internal energy) */
  if (posix_memalign(
          (void **)&cooling->table.Ucooling, SWIFT_STRUCT_ALIGNMENT,
          colibre_cooling_N_redshifts * colibre_cooling_N_internalenergy *
              colibre_cooling_N_metallicity * colibre_cooling_N_density *
              colibre_cooling_N_cooltypes * sizeof(float)) != 0)
    error("Failed to allocate Ucooling array\n");

  dataset = H5Dopen(tempfile_id, "/Udep/Cooling", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.Ucooling);
  if (status < 0) error("error reading Ucooling\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  /* Heating (temperature) */
  if (posix_memalign(
          (void **)&cooling->table.Theating, SWIFT_STRUCT_ALIGNMENT,
          colibre_cooling_N_redshifts * colibre_cooling_N_temperature *
              colibre_cooling_N_metallicity * colibre_cooling_N_density *
              colibre_cooling_N_heattypes * sizeof(float)) != 0)
    error("Failed to allocate Theating array\n");

  dataset = H5Dopen(tempfile_id, "/Tdep/Heating", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.Theating);
  if (status < 0) error("error reading Theating\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  /* Heating (internal energy) */
  if (posix_memalign(
          (void **)&cooling->table.Uheating, SWIFT_STRUCT_ALIGNMENT,
          colibre_cooling_N_redshifts * colibre_cooling_N_internalenergy *
              colibre_cooling_N_metallicity * colibre_cooling_N_density *
              colibre_cooling_N_heattypes * sizeof(float)) != 0)
    error("Failed to allocate Uheating array\n");

  dataset = H5Dopen(tempfile_id, "/Udep/Heating", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.Uheating);
  if (status < 0) error("error reading Uheating\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  /* Electron fraction (temperature) */
  if (posix_memalign(
          (void **)&cooling->table.Telectron_fraction, SWIFT_STRUCT_ALIGNMENT,
          colibre_cooling_N_redshifts * colibre_cooling_N_temperature *
              colibre_cooling_N_metallicity * colibre_cooling_N_density *
              colibre_cooling_N_electrontypes * sizeof(float)) != 0)
    error("Failed to allocate Telectron_fraction array\n");

  dataset = H5Dopen(tempfile_id, "/Tdep/ElectronFractionsVol", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.Telectron_fraction);
  if (status < 0) error("error reading electron_fraction (temperature)\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  /* Electron fraction (internal energy) */
  if (posix_memalign(
          (void **)&cooling->table.Uelectron_fraction, SWIFT_STRUCT_ALIGNMENT,
          colibre_cooling_N_redshifts * colibre_cooling_N_internalenergy *
              colibre_cooling_N_metallicity * colibre_cooling_N_density *
              colibre_cooling_N_electrontypes * sizeof(float)) != 0)
    error("Failed to allocate Uelectron_fraction array\n");

  dataset = H5Dopen(tempfile_id, "/Udep/ElectronFractionsVol", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.Uelectron_fraction);
  if (status < 0) error("error reading electron_fraction (internal energy)\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  /* Internal energy from temperature */
  if (posix_memalign((void **)&cooling->table.U_from_T, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_redshifts *
                         colibre_cooling_N_temperature *
                         colibre_cooling_N_metallicity *
                         colibre_cooling_N_density * sizeof(float)) != 0)
    error("Failed to allocate U_from_T array\n");

  dataset = H5Dopen(tempfile_id, "/Tdep/U_from_T", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.U_from_T);
  if (status < 0) error("error reading U_from_T array\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  /* Temperature from interal energy */
  if (posix_memalign((void **)&cooling->table.T_from_U, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_redshifts *
                         colibre_cooling_N_internalenergy *
                         colibre_cooling_N_metallicity *
                         colibre_cooling_N_density * sizeof(float)) != 0)
    error("Failed to allocate T_from_U array\n");

  dataset = H5Dopen(tempfile_id, "/Udep/T_from_U", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.T_from_U);
  if (status < 0) error("error reading T_from_U array\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing cooling dataset");

  /* Thermal equilibrium temperature */
  if (posix_memalign((void **)&cooling->table.logTeq, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_redshifts *
                         colibre_cooling_N_metallicity *
                         colibre_cooling_N_density * sizeof(float)) != 0)
    error("Failed to allocate logTeq array\n");

  dataset = H5Dopen(tempfile_id, "/ThermEq/Temperature", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.logTeq);
  if (status < 0) error("error reading Teq array\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing logTeq dataset");

  /* Mean particle mass at thermal equilibrium temperature */
  if (posix_memalign(
          (void **)&cooling->table.meanpartmass_Teq, SWIFT_STRUCT_ALIGNMENT,
          colibre_cooling_N_redshifts * colibre_cooling_N_metallicity *
              colibre_cooling_N_density * sizeof(float)) != 0)
    error("Failed to allocate mu array\n");

  dataset = H5Dopen(tempfile_id, "/ThermEq/MeanParticleMass", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.meanpartmass_Teq);
  if (status < 0) error("error reading mu array\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing mu dataset");

  /* Hydrogen fractions at thermal equilibirum temperature */
  if (posix_memalign(
          (void **)&cooling->table.logHfracs_Teq, SWIFT_STRUCT_ALIGNMENT,
          colibre_cooling_N_redshifts * colibre_cooling_N_metallicity *
              colibre_cooling_N_density * 3 * sizeof(float)) != 0)
    error("Failed to allocate hydrogen fractions array\n");

  dataset = H5Dopen(tempfile_id, "/ThermEq/HydrogenFractionsVol", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.logHfracs_Teq);
  if (status < 0) error("error reading hydrogen fractions array\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing hydrogen fractions dataset");

  /* All hydrogen fractions */
  if (posix_memalign(
          (void **)&cooling->table.logHfracs_all, SWIFT_STRUCT_ALIGNMENT,
          colibre_cooling_N_redshifts * colibre_cooling_N_temperature *
              colibre_cooling_N_metallicity * colibre_cooling_N_density * 3 *
              sizeof(float)) != 0)
    error("Failed to allocate big hydrogen fractions array\n");

  dataset = H5Dopen(tempfile_id, "/Tdep/HydrogenFractionsVol", H5P_DEFAULT);
  status = H5Dread(dataset, H5T_NATIVE_FLOAT, H5S_ALL, H5S_ALL, H5P_DEFAULT,
                   cooling->table.logHfracs_all);
  if (status < 0) error("error reading big hydrogen fractions array\n");
  status = H5Dclose(dataset);
  if (status < 0) error("error closing big hydrogen fractions dataset");

  /* Close the file */
  H5Fclose(tempfile_id);

  /* Pressure at thermal equilibrium temperature */
  if (posix_memalign((void **)&cooling->table.logPeq, SWIFT_STRUCT_ALIGNMENT,
                     colibre_cooling_N_redshifts *
                         colibre_cooling_N_metallicity *
                         colibre_cooling_N_density * sizeof(float)) != 0)
    error("Failed to allocate logPeq array\n");

  const float log10_kB_cgs = cooling->log10_kB_cgs;

  /* Compute the pressures at thermal eq. */
  for (int ired = 0; ired < colibre_cooling_N_redshifts; ired++) {
    for (int imet = 0; imet < colibre_cooling_N_metallicity; imet++) {

      const int index_XH =
          row_major_index_2d(imet, 0, colibre_cooling_N_metallicity,
                             colibre_cooling_N_elementtypes);

      const float log10_XH = cooling->LogMassFractions[index_XH];

      for (int iden = 0; iden < colibre_cooling_N_density; iden++) {

        const int index_Peq = row_major_index_3d(
            ired, imet, iden, colibre_cooling_N_redshifts,
            colibre_cooling_N_metallicity, colibre_cooling_N_density);

        cooling->table.logPeq[index_Peq] =
            cooling->nH[iden] + cooling->table.logTeq[index_Peq] - log10_XH -
            log10(cooling->table.meanpartmass_Teq[index_Peq]) + log10_kB_cgs;
      }
    }
  }

#ifdef SWIFT_DEBUG_CHECKS
  message("Done reading in general cooling table");
#endif

#else
  error("Need HDF5 to read cooling tables");
#endif
}