cooling_rates.h 28.4 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
/*******************************************************************************
 * This file is part of SWIFT.
 * Copyright (c) 2020 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/>.
 *
 ******************************************************************************/
#ifndef SWIFT_QLA_COOLING_RATES_H
#define SWIFT_QLA_COOLING_RATES_H

22
/* Config parameters. */
23
24
25
#include "../config.h"

/* Local includes. */
26
#include "chemistry_struct.h"
27
#include "cooling_tables.h"
28
#include "error.h"
29
30
31
32
33
34
35
36
37
#include "exp10.h"
#include "interpolate.h"

/**
 * @brief Compute ratio of mass fraction to solar mass fraction
 * for each element carried by a given particle.
 *
 * The solar abundances are taken from the tables themselves.
 *
38
39
40
41
42
43
44
45
 * The COLIBRE chemistry model does not track S and Ca. We assume
 * that their abundance with respect to solar is the same as
 * the ratio for Si.
 *
 * The other un-tracked elements are scaled with the total metallicity.
 *
 * We optionally apply a correction if the user asked for a different
 * ratio.
46
47
 *
 * We also re-order the elements such that they match the order of the
48
 * tables. This is [H, He, C, N, O, Ne, Mg, Si, S, Ca, Fe, OA].
49
 *
50
51
52
53
54
55
 * The solar abundances table (from the cooling struct) is arranged as
 * [H, He, C, N, O, Ne, Mg, Si, S, Ca, Fe].
 *
 * @param p Pointer to #part struct.
 * @param cooling #cooling_function_data struct.
 * @param ratio_solar (return) Array of ratios to solar abundances.
56
57
58
 *
 * @return The log10 of the total metallicity with respect to solar, i.e.
 * log10(Z / Z_sun).
59
 */
60
__attribute__((always_inline)) INLINE static float abundance_ratio_to_solar(
61
62
    const struct part *p, const struct cooling_function_data *cooling,
    const struct phys_const *phys_const,
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
    float ratio_solar[qla_cooling_N_elementtypes]) {

  /* Convert mass fractions to abundances (nx/nH) and compute metal mass */
  float totmass = 0., metalmass = 0.;
  for (enum qla_cooling_element elem = element_H; elem < element_OA; elem++) {

    double Z_mass_frac = -1.f;
    double Z_mass_frac_H = 1. - phys_const->const_primordial_He_fraction;

    /* Assume primordial abundances */
    if (elem == element_H) {
      Z_mass_frac = 1. - phys_const->const_primordial_He_fraction;
    } else if (elem == element_He) {
      Z_mass_frac = phys_const->const_primordial_He_fraction;
    } else {
      Z_mass_frac = 0.f;
    }

    const int indx1d =
        row_major_index_2d(cooling->indxZsol, elem, qla_cooling_N_metallicity,
                           qla_cooling_N_elementtypes);

    const float Mfrac = Z_mass_frac;

    /* ratio_X = ((M_x/M) / (M_H/M)) * (m_H / m_X) * (1 / Z_sun_X) */
    ratio_solar[elem] =
        (Mfrac / Z_mass_frac_H) * cooling->atomicmass[element_H] *
        cooling->atomicmass_inv[elem] * cooling->Abundances_inv[indx1d];

    totmass += Mfrac;
    if (elem > element_He) metalmass += Mfrac;
  }

  /* Compute metallicity (Z / Z_sun) of the elements we explicitly track. */
  /* float ZZsol = (metalmass / totmass) * cooling->Zsol_inv[0];
  if (ZZsol <= 0.0f) ZZsol = FLT_MIN;
  const float logZZsol = log10f(ZZsol); */

  /* Get total metallicity [Z/Z_sun] from the particle data */
  const float Z_total =
      (float)chemistry_get_total_metal_mass_fraction_for_cooling(p);
  float ZZsol = Z_total * cooling->Zsol_inv[0];
  if (ZZsol <= 0.0f) ZZsol = FLT_MIN;
  const float logZZsol = log10f(ZZsol);

  /* All other elements (element_OA): scale with metallicity */
  ratio_solar[element_OA] = ZZsol;

  /* Get index and offset from the metallicity table conresponding to this
   * metallicity */
  int met_index;
  float d_met;

  get_index_1d(cooling->Metallicity, qla_cooling_N_metallicity, logZZsol,
               &met_index, &d_met);

  /* At this point ratio_solar is (nx/nH) / (nx/nH)_sol.
   * To multiply with the tables, we want the individual abundance ratio
   * relative to what is used in the tables for each metallicity
   */

  /* For example: for a metallicity of 1 per cent solar, the metallicity bin
   * for logZZsol = -2 has already the reduced cooling rates for each element;
   * it should therefore NOT be multiplied by 0.01 again.
   *
   * BUT: if e.g. Carbon is twice as abundant as the solar abundance ratio,
   * i.e. nC / nH = 0.02 * (nC/nH)_sol for the overall metallicity of 0.01,
   * the Carbon cooling rate is multiplied by 2
   *
   * We only do this if we are not in the primodial metallicity bin in which
   * case the ratio to solar should be 0.
   */

  for (int i = 0; i < qla_cooling_N_elementtypes; i++) {

    /* Are we considering a metal and are *not* in the primodial metallicity
     * bin? Or are we looking at H or He? */
    if ((met_index > 0) || (i == element_H) || (i == element_He)) {

      /* Get min/max abundances */
      const float log_nx_nH_min = cooling->LogAbundances[row_major_index_2d(
          met_index, i, qla_cooling_N_metallicity, qla_cooling_N_elementtypes)];

      const float log_nx_nH_max = cooling->LogAbundances[row_major_index_2d(
          met_index + 1, i, qla_cooling_N_metallicity,
          qla_cooling_N_elementtypes)];

      /* Get solar abundances */
      const float log_nx_nH_sol = cooling->LogAbundances[row_major_index_2d(
          cooling->indxZsol, i, qla_cooling_N_metallicity,
          qla_cooling_N_elementtypes)];

      /* Interpolate ! (linearly in log-space) */
      const float log_nx_nH =
          (log_nx_nH_min * (1.f - d_met) + log_nx_nH_max * d_met) -
          log_nx_nH_sol;

      ratio_solar[i] *= exp10f(-log_nx_nH);

    } else {

      /* Primordial bin --> Z/Z_sun is 0 for that element */
      ratio_solar[i] = 0.;
    }
  }

  /* at this point ratio_solar is (nx/nH) / (nx/nH)_table [Z],
   * the metallicity dependent abundance ratio for solar abundances.
   * We also return the total metallicity */

  return logZZsol;
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
}

/**
 * @brief Computes the extra heat from Helium reionisation at a given redshift.
 *
 * We follow the implementation of Wiersma et al. 2009, MNRAS, 399, 574-600,
 * section. 2. The calculation returns energy in CGS.
 *
 * Note that delta_z is negative.
 *
 * @param z The current redshift.
 * @param delta_z The change in redhsift over the course of this time-step.
 * @param cooling The #cooling_function_data used in the run.
 * @return Helium reionization energy in CGS units.
 */
189
190
191
__attribute__((always_inline)) INLINE static double
eagle_helium_reionization_extraheat(
    double z, double delta_z, const struct cooling_function_data *cooling) {
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215

#ifdef SWIFT_DEBUG_CHECKS
  if (delta_z > 0.f) error("Invalid value for delta_z. Should be negative.");
#endif

  /* Recover the values we need */
  const double z_centre = cooling->He_reion_z_centre;
  const double z_sigma = cooling->He_reion_z_sigma;
  const double heat_cgs = cooling->He_reion_heat_cgs;

  double extra_heat = 0.;

  /* Integral of the Gaussian between z and z - delta_z */
  extra_heat += erf((z - delta_z - z_centre) / (M_SQRT2 * z_sigma));
  extra_heat -= erf((z - z_centre) / (M_SQRT2 * z_sigma));

  /* Multiply by the normalisation factor */
  extra_heat *= heat_cgs * 0.5;

  return extra_heat;
}

/**
 * @brief Computes the log_10 of the temperature corresponding to a given
216
 * internal energy, hydrogen number density, metallicity and redshift
217
218
219
220
221
 *
 * @param log_10_u_cgs Log base 10 of internal energy in cgs.
 * @param redshift Current redshift.
 * @param n_H_index Index along the Hydrogen density dimension.
 * @param d_n_H Offset between Hydrogen density and table[n_H_index].
222
223
224
225
 * @param met_index Index along the metallicity dimension.
 * @param d_met Offset between metallicity and table[met_index].
 * @param red_index Index along the redshift dimension.
 * @param d_red Offset between redshift and table[red_index].
226
227
228
 * @param cooling #cooling_function_data structure.
 *
 * @return log_10 of the temperature.
229
230
 *
 * TO DO: outside table ranges, it uses at the moment the minimum, maximu value
231
 */
232
__attribute__((always_inline)) INLINE static float qla_convert_u_to_temp(
233
234
235
    const double log_10_u_cgs, const float redshift, const int n_H_index,
    const float d_n_H, const int met_index, const float d_met,
    const int red_index, const float d_red,
236
237
238
239
240
    const struct cooling_function_data *cooling) {

  /* Get index of u along the internal energy axis */
  int u_index;
  float d_u;
241
242

  get_index_1d(cooling->Therm, qla_cooling_N_internalenergy, log_10_u_cgs,
243
244
245
               &u_index, &d_u);

  /* Interpolate temperature table to return temperature for current
246
   * internal energy */
247
  float log_10_T;
248
249
250
251
252
253

  /* Temperature from internal energy */
  log_10_T = interpolation_4d(
      cooling->table.T_from_U, red_index, u_index, met_index, n_H_index, d_red,
      d_u, d_met, d_n_H, qla_cooling_N_redshifts, qla_cooling_N_internalenergy,
      qla_cooling_N_metallicity, qla_cooling_N_density);
254
255
256
257
258
259
260
261
262
263
264
265
266

  /* Special case for temperatures below the start of the table */
  if (u_index == 0 && d_u == 0.f) {

    /* The temperature is multiplied by u / 10^T[0]
     * where T[0] is the first entry in the table */
    log_10_T += log_10_u_cgs - cooling->Temp[0];
  }

  return log_10_T;
}

/**
267
268
269
270
271
272
273
274
275
276
277
278
 * @brief Computes the log_10 of the internal energy corresponding to a given
 * temperature, hydrogen number density, metallicity and redshift
 *
 * @param log_10_T Log base 10 of temperature in K
 * @param redshift Current redshift.
 * @param n_H_index Index along the Hydrogen density dimension.
 * @param d_n_H Offset between Hydrogen density and table[n_H_index].
 * @param met_index Index along the metallicity dimension.
 * @param d_met Offset between metallicity and table[met_index].
 * @param red_index Index along the redshift dimension.
 * @param d_red Offset between redshift and table[red_index].
 * @param cooling #cooling_function_data structure.
279
 *
280
 * @return log_10 of the internal energy in cgs
281
 *
282
 * TO DO: outside table ranges, it uses at the moment the minimum, maximu value
283
 */
284
__attribute__((always_inline)) INLINE static float qla_convert_temp_to_u(
285
286
287
    const double log_10_T, const float redshift, const int n_H_index,
    const float d_n_H, const int met_index, const float d_met,
    const int red_index, const float d_red,
288
289
290
291
292
    const struct cooling_function_data *cooling) {

  /* Get index of u along the internal energy axis */
  int T_index;
  float d_T;
293

294
295
  get_index_1d(cooling->Temp, qla_cooling_N_temperature, log_10_T, &T_index,
               &d_T);
296

297
298
299
  /* Interpolate internal energy table to return internal energy for current
   * temperature */
  float log_10_U;
300

301
302
303
304
305
306
307
  /* Internal energy from temperature*/
  log_10_U = interpolation_4d(
      cooling->table.U_from_T, red_index, T_index, met_index, n_H_index, d_red,
      d_T, d_met, d_n_H, qla_cooling_N_redshifts, qla_cooling_N_temperature,
      qla_cooling_N_metallicity, qla_cooling_N_density);

  return log_10_U;
308
309
310
}

/**
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
 * @brief Computes the mean particle mass for a given
 * metallicity, temperature, redshift, and density.
 *
 * @param log_T_cgs Log base 10 of temperature in K
 * @param redshift Current redshift
 * @param n_H_cgs Hydrogen number density in cgs
 * @param ZZsol Metallicity relative to the solar value from the tables
 * @param n_H_index Index along the Hydrogen number density dimension
 * @param d_n_H Offset between Hydrogen density and table[n_H_index]
 * @param met_index Index along the metallicity dimension
 * @param d_met Offset between metallicity and table[met_index]
 * @param red_index Index along redshift dimension
 * @param d_red Offset between redshift and table[red_index]
 * @param cooling #cooling_function_data structure
 *
 * @return linear electron density in cm-3 (NOT the electron fraction)
327
 */
328
INLINE static float qla_meanparticlemass_temperature(
329
330
331
332
    const double log_T_cgs, const double redshift, const double n_H_cgs,
    const float ZZsol, const int n_H_index, const float d_n_H,
    const int met_index, const float d_met, const int red_index,
    const float d_red, const struct cooling_function_data *cooling) {
333

334
  /* Get index of T along the temperature axis */
335
336
  int T_index;
  float d_T;
337
338

  get_index_1d(cooling->Temp, qla_cooling_N_temperature, log_T_cgs, &T_index,
339
340
               &d_T);

341
342
343
344
  const float mu = interpolation_4d(
      cooling->table.Tmu, red_index, T_index, met_index, n_H_index, d_red, d_T,
      d_met, d_n_H, qla_cooling_N_redshifts, qla_cooling_N_temperature,
      qla_cooling_N_metallicity, qla_cooling_N_density);
345

346
347
  return mu;
}
348

349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
/**
 * @brief Computes the electron density for a given element
 * abundance ratio, internal energy, redshift, and density.
 *
 * @param log_u_cgs Log base 10 of internal energy in cgs [erg g-1]
 * @param redshift Current redshift
 * @param n_H_cgs Hydrogen number density in cgs
 * @param ZZsol Metallicity relative to the solar value from the tables
 * @param abundance_ratio Abundance ratio for each element x relative to solar
 * @param n_H_index Index along the Hydrogen number density dimension
 * @param d_n_H Offset between Hydrogen density and table[n_H_index]
 * @param met_index Index along the metallicity dimension
 * @param d_met Offset between metallicity and table[met_index]
 * @param red_index Index along redshift dimension
 * @param d_red Offset between redshift and table[red_index]
 * @param cooling #cooling_function_data structure
 *
 * @return linear electron density in cm-3 (NOT the electron fraction)
 */
INLINE static float qla_electron_density(
369
370
371
372
    const double log_u_cgs, const double redshift, const double n_H_cgs,
    const float ZZsol, const float abundance_ratio[qla_cooling_N_elementtypes],
    const int n_H_index, const float d_n_H, const int met_index,
    const float d_met, const int red_index, const float d_red,
373
    const struct cooling_function_data *cooling) {
374

375
376
377
  /* Get index of u along the internal energy axis */
  int U_index;
  float d_U;
378

379
380
  get_index_1d(cooling->Therm, qla_cooling_N_internalenergy, log_u_cgs,
               &U_index, &d_U);
381

382
383
384
385
386
387
388
  /* n_e / n_H */
  const float electron_fraction = interpolation4d_plus_summation(
      cooling->table.Uelectron_fraction, abundance_ratio, element_H,
      qla_cooling_N_electrontypes - 4, red_index, U_index, met_index, n_H_index,
      d_red, d_U, d_met, d_n_H, qla_cooling_N_redshifts,
      qla_cooling_N_internalenergy, qla_cooling_N_metallicity,
      qla_cooling_N_density, qla_cooling_N_electrontypes);
389

390
391
  return electron_fraction * n_H_cgs;
}
392

393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
/**
 * @brief Computes the electron density for a given element
 * abundance ratio, temperature, redshift, and density.
 *
 * @param log_T_cgs Log base 10 of temperature
 * @param redshift Current redshift
 * @param n_H_cgs Hydrogen number density in cgs
 * @param ZZsol Metallicity relative to the solar value from the tables
 * @param abundance_ratio Abundance ratio for each element x relative to solar
 * @param n_H_index Index along the Hydrogen number density dimension
 * @param d_n_H Offset between Hydrogen density and table[n_H_index]
 * @param met_index Index along the metallicity dimension
 * @param d_met Offset between metallicity and table[met_index]
 * @param red_index Index along redshift dimension
 * @param d_red Offset between redshift and table[red_index]
 * @param cooling #cooling_function_data structure
 *
 * @return linear electron density in cm-3 (NOT the electron fraction)
 */
INLINE static float qla_electron_density_temperature(
413
414
415
416
    const double log_T_cgs, const double redshift, const double n_H_cgs,
    const float ZZsol, const float abundance_ratio[qla_cooling_N_elementtypes],
    const int n_H_index, const float d_n_H, const int met_index,
    const float d_met, const int red_index, const float d_red,
417
    const struct cooling_function_data *cooling) {
418

419
420
421
  /* Get index of u along the internal energy axis */
  int T_index;
  float d_T;
422

423
424
  get_index_1d(cooling->Temp, qla_cooling_N_temperature, log_T_cgs, &T_index,
               &d_T);
425

426
427
428
429
430
431
432
  /* n_e / n_H */
  const float electron_fraction = interpolation4d_plus_summation(
      cooling->table.Telectron_fraction, abundance_ratio, element_H,
      qla_cooling_N_electrontypes - 4, red_index, T_index, met_index, n_H_index,
      d_red, d_T, d_met, d_n_H, qla_cooling_N_redshifts,
      qla_cooling_N_internalenergy, qla_cooling_N_metallicity,
      qla_cooling_N_density, qla_cooling_N_electrontypes);
433

434
435
  return electron_fraction * n_H_cgs;
}
436

437
/**
438
 * @brief Computes the net cooling rate (heating - cooling) for a given element
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
 * abundance ratio, internal energy, redshift, and density. The unit of the net
 * cooling rate is Lambda / nH**2 [erg cm^3 s-1] and all input values are in
 * cgs. The Compton cooling is not taken from the tables but calculated
 * analytically and added separately
 *
 * @param log_u_cgs Log base 10 of internal energy in cgs [erg g-1]
 * @param redshift Current redshift
 * @param n_H_cgs Hydrogen number density in cgs
 * @param abundance_ratio Abundance ratio for each element x relative to solar
 * @param n_H_index Index along the Hydrogen number density dimension
 * @param d_n_H Offset between Hydrogen density and table[n_H_index]
 * @param met_index Index along the metallicity dimension
 * @param d_met Offset between metallicity and table[met_index]
 * @param red_index Index along redshift dimension
 * @param d_red Offset between redshift and table[red_index]
 * @param cooling #cooling_function_data structure
 *
 * @param onlyicool if true / 1 only plot cooling channel icool
 * @param onlyiheat if true / 1 only plot cooling channel iheat
 * @param icool cooling channel to be used
 * @param iheat heating channel to be used
 *
 * Throughout the code: onlyicool = onlyiheat = icool = iheat = 0
 * These are only used for testing.
 */
INLINE static double qla_cooling_rate(
465
466
467
468
469
470
    const double log_u_cgs, const double redshift, const double n_H_cgs,
    const float abundance_ratio[qla_cooling_N_elementtypes],
    const int n_H_index, const float d_n_H, const int met_index,
    const float d_met, const int red_index, const float d_red,
    const struct cooling_function_data *cooling, const int onlyicool,
    const int onlyiheat, const int icool, const int iheat) {
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485

  /* Set weights for cooling rates */
  float weights_cooling[qla_cooling_N_cooltypes - 2];
  for (int i = 0; i < qla_cooling_N_cooltypes - 2; i++) {

    if (i < qla_cooling_N_elementtypes) {
      /* Use abundance ratios */
      weights_cooling[i] = abundance_ratio[i];
    } else if (i == cooltype_Compton) {
      /* added analytically later, do not use value from table*/
      weights_cooling[i] = 0.f;
    } else {
      /* use same abundances as in the tables */
      weights_cooling[i] = 1.f;
    }
486
487
  }

488
489
490
491
492
  /* If we care only about one channel, cancel all the other ones */
  if (onlyicool != 0) {
    for (int i = 0; i < qla_cooling_N_cooltypes - 2; i++) {
      if (i != icool) weights_cooling[i] = 0.f;
    }
493
494
  }

495
496
497
498
499
500
501
502
  /* Set weights for heating rates */
  float weights_heating[qla_cooling_N_heattypes - 2];
  for (int i = 0; i < qla_cooling_N_heattypes - 2; i++) {
    if (i < qla_cooling_N_elementtypes) {
      weights_heating[i] = abundance_ratio[i];
    } else {
      weights_heating[i] = 1.f; /* use same abundances as in the tables */
    }
503
504
  }

505
506
507
508
  /* If we care only about one channel, cancel all the other ones */
  if (onlyiheat != 0) {
    for (int i = 0; i < qla_cooling_N_heattypes - 2; i++) {
      if (i != iheat) weights_heating[i] = 0.f;
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
  /* Get index of u along the internal energy axis */
  int U_index;
  float d_U;
  get_index_1d(cooling->Therm, qla_cooling_N_internalenergy, log_u_cgs,
               &U_index, &d_U);

  /* n_e / n_H */
  const double electron_fraction = interpolation4d_plus_summation(
      cooling->table.Uelectron_fraction, abundance_ratio, /* */
      element_H, qla_cooling_N_electrontypes - 4,         /* */
      red_index, U_index, met_index, n_H_index,           /* */
      d_red, d_U, d_met, d_n_H,                           /* */
      qla_cooling_N_redshifts,                            /* */
      qla_cooling_N_internalenergy,                       /* */
      qla_cooling_N_metallicity,                          /* */
      qla_cooling_N_density,                              /* */
      qla_cooling_N_electrontypes);                       /* */

  /* Lambda / n_H**2 */
  const double cooling_rate = interpolation4d_plus_summation(
      cooling->table.Ucooling, weights_cooling, /* */
      element_H, qla_cooling_N_cooltypes - 3,   /* */
      red_index, U_index, met_index, n_H_index, /* */
      d_red, d_U, d_met, d_n_H,                 /* */
      qla_cooling_N_redshifts,                  /* */
      qla_cooling_N_internalenergy,             /* */
      qla_cooling_N_metallicity,                /* */
      qla_cooling_N_density,                    /* */
      qla_cooling_N_cooltypes);                 /* */

  /* Gamma / n_H**2 */
  const double heating_rate = interpolation4d_plus_summation(
      cooling->table.Uheating, weights_heating, /* */
      element_H, qla_cooling_N_heattypes - 3,   /* */
      red_index, U_index, met_index, n_H_index, /* */
      d_red, d_U, d_met, d_n_H,                 /* */
      qla_cooling_N_redshifts,                  /* */
      qla_cooling_N_internalenergy,             /* */
      qla_cooling_N_metallicity,                /* */
      qla_cooling_N_density,                    /* */
      qla_cooling_N_heattypes);                 /* */

  /* Temperature from internal energy */
  const double logtemp =
      interpolation_4d(cooling->table.T_from_U,                  /* */
                       red_index, U_index, met_index, n_H_index, /* */
                       d_red, d_U, d_met, d_n_H,                 /* */
                       qla_cooling_N_redshifts,                  /* */
                       qla_cooling_N_internalenergy,             /* */
                       qla_cooling_N_metallicity,                /* */
                       qla_cooling_N_density);                   /* */
                                                                 /* */
  const double temp = exp10(logtemp);

  /* Compton cooling/heating */
  double Compton_cooling_rate = 0.;
  if (onlyicool == 0 || (onlyicool == 1 && icool == cooltype_Compton)) {

    const double zp1 = 1. + redshift;
    const double zp1p2 = zp1 * zp1;
    const double zp1p4 = zp1p2 * zp1p2;

    /* CMB temperature at this redshift */
    const double T_CMB = cooling->T_CMB_0 * zp1;

    /* Analytic Compton cooling rate: Lambda_Compton / n_H**2 */
    Compton_cooling_rate = cooling->compton_rate_cgs * (temp - T_CMB) * zp1p4 *
                           electron_fraction / n_H_cgs;
  }
581

582
583
584
  /* Return the net heating rate (Lambda_heat - Lambda_cool) */
  return heating_rate - cooling_rate - Compton_cooling_rate;
}
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
/**
 * @brief Computes the net cooling rate (cooling - heating) for a given element
 * abundance ratio, temperature, redshift, and density. The unit of the net
 * cooling rate is Lambda / nH**2 [erg cm^3 s-1] and all input values are in
 * cgs. The Compton cooling is not taken from the tables but calculated
 * analytically and added separately
 *
 * @param log_T_cgs Log base 10 of temperature in K
 * @param redshift Current redshift
 * @param n_H_cgs Hydrogen number density in cgs
 * @param abundance_ratio Abundance ratio for each element x relative to solar
 * @param n_H_index Index along the Hydrogen number density dimension
 * @param d_n_H Offset between Hydrogen density and table[n_H_index]
 * @param met_index Index along the metallicity dimension
 * @param d_met Offset between metallicity and table[met_index]
 * @param red_index Index along redshift dimension
 * @param d_red Offset between redshift and table[red_index]
 * @param cooling #cooling_function_data structure
 *
 * @param onlyicool if true / 1 only plot cooling channel icool
 * @param onlyiheat if true / 1 only plot cooling channel iheat
 * @param icool cooling channel to be used
 * @param iheat heating channel to be used
 *
 * Throughout the code: onlyicool = onlyiheat = icool = iheat = 0
 * These are only used for testing.
 */
INLINE static double qla_cooling_rate_temperature(
614
615
616
617
618
619
    const double log_T_cgs, const double redshift, const double n_H_cgs,
    const float abundance_ratio[qla_cooling_N_elementtypes],
    const int n_H_index, const float d_n_H, const int met_index,
    const float d_met, const int red_index, const float d_red,
    const struct cooling_function_data *cooling, const int onlyicool,
    const int onlyiheat, const int icool, const int iheat) {
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635

  /* Set weights for cooling rates */
  float weights_cooling[qla_cooling_N_cooltypes - 2];
  for (int i = 0; i < qla_cooling_N_cooltypes - 2; i++) {

    if (i < qla_cooling_N_elementtypes) {
      /* Use abundance ratios */
      weights_cooling[i] = abundance_ratio[i];
    } else if (i == cooltype_Compton) {
      /* added analytically later, do not use value from table*/
      weights_cooling[i] = 0.f;
    } else {
      /* use same abundances as in the tables */
      weights_cooling[i] = 1.f;
    }
  }
636

637
638
639
640
  /* If we care only about one channel, cancel all the other ones */
  if (onlyicool != 0) {
    for (int i = 0; i < qla_cooling_N_cooltypes - 2; i++) {
      if (i != icool) weights_cooling[i] = 0.f;
641
642
643
    }
  }

644
645
646
647
648
649
650
  /* Set weights for heating rates */
  float weights_heating[qla_cooling_N_heattypes - 2];
  for (int i = 0; i < qla_cooling_N_heattypes - 2; i++) {
    if (i < qla_cooling_N_elementtypes) {
      weights_heating[i] = abundance_ratio[i];
    } else {
      weights_heating[i] = 1.f; /* use same abundances as in the tables */
651
652
653
    }
  }

654
655
656
657
658
  /* If we care only about one channel, cancel all the other ones */
  if (onlyiheat != 0) {
    for (int i = 0; i < qla_cooling_N_heattypes - 2; i++) {
      if (i != iheat) weights_heating[i] = 0.f;
    }
659
660
  }

661
662
663
664
665
  /* Get index of T along the internal energy axis */
  int T_index;
  float d_T;
  get_index_1d(cooling->Temp, qla_cooling_N_temperature, log_T_cgs, &T_index,
               &d_T);
666

667
668
669
670
671
672
673
674
675
676
677
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
  /* n_e / n_H */
  const double electron_fraction = interpolation4d_plus_summation(
      cooling->table.Telectron_fraction, abundance_ratio, /* */
      element_H, qla_cooling_N_electrontypes - 4,         /* */
      red_index, T_index, met_index, n_H_index,           /* */
      d_red, d_T, d_met, d_n_H,                           /* */
      qla_cooling_N_redshifts,                            /* */
      qla_cooling_N_temperature,                          /* */
      qla_cooling_N_metallicity,                          /* */
      qla_cooling_N_density,                              /* */
      qla_cooling_N_electrontypes);                       /* */

  /* Lambda / n_H**2 */
  const double cooling_rate = interpolation4d_plus_summation(
      cooling->table.Tcooling, weights_cooling, /* */
      element_H, qla_cooling_N_cooltypes - 3,   /* */
      red_index, T_index, met_index, n_H_index, /* */
      d_red, d_T, d_met, d_n_H,                 /* */
      qla_cooling_N_redshifts,                  /* */
      qla_cooling_N_temperature,                /* */
      qla_cooling_N_metallicity,                /* */
      qla_cooling_N_density,                    /* */
      qla_cooling_N_cooltypes);                 /* */

  /* Gamma / n_H**2 */
  const double heating_rate = interpolation4d_plus_summation(
      cooling->table.Theating, weights_heating, /* */
      element_H, qla_cooling_N_heattypes - 3,   /* */
      red_index, T_index, met_index, n_H_index, /* */
      d_red, d_T, d_met, d_n_H,                 /* */
      qla_cooling_N_redshifts,                  /* */
      qla_cooling_N_temperature,                /* */
      qla_cooling_N_metallicity,                /* */
      qla_cooling_N_density,                    /* */
      qla_cooling_N_heattypes);                 /* */

  const double temp = exp10(log_T_cgs);

  double Compton_cooling_rate = 0.;
  if (onlyicool == 0 || (onlyicool == 1 && icool == cooltype_Compton)) {

    const double zp1 = 1. + redshift;
    const double zp1p2 = zp1 * zp1;
    const double zp1p4 = zp1p2 * zp1p2;

    /* CMB temperature at this redshift */
    const double T_CMB = cooling->T_CMB_0 * zp1;

    /* Analytic Compton cooling rate: Lambda_Compton / n_H**2 */
    Compton_cooling_rate = cooling->compton_rate_cgs * (temp - T_CMB) * zp1p4 *
                           electron_fraction / n_H_cgs;
  }
719

720
721
  /* Return the net heating rate (Lambda_heat - Lambda_cool) */
  return heating_rate - cooling_rate - Compton_cooling_rate;
722
723
724
}

#endif /* SWIFT_QLA_COOLING_RATES_H */