hydro.h 19.9 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
/*******************************************************************************
 * This file is part of SWIFT.
 * Copyright (c) 2016 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/>.
 *
 ******************************************************************************/
#ifndef SWIFT_PRESSURE_ENTROPY_HYDRO_H
#define SWIFT_PRESSURE_ENTROPY_HYDRO_H

/**
23
24
25
 * @file PressureEntropy/hydro.h
 * @brief Pressure-Entropy implementation of SPH (Non-neighbour loop
 * equations)
26
27
28
29
 *
 * The thermal variable is the entropy (S) and the entropy is smoothed over
 * contact discontinuities to prevent spurious surface tension.
 *
30
31
 * Follows eqautions (19), (21) and (22) of Hopkins, P., MNRAS, 2013,
 * Volume 428, Issue 4, pp. 2840-2856 with a simple Balsara viscosity term.
32
33
34
 */

#include "adiabatic_index.h"
35
#include "approx_math.h"
36
#include "cosmology.h"
37
38
39
#include "dimension.h"
#include "equation_of_state.h"
#include "hydro_properties.h"
40
#include "hydro_space.h"
41
#include "kernel_hydro.h"
42
#include "minmax.h"
43
44

/**
45
 * @brief Returns the comoving internal energy of a particle
46
47
48
 *
 * @param p The particle of interest
 */
49
50
__attribute__((always_inline)) INLINE static float
hydro_get_comoving_internal_energy(const struct part *restrict p) {
51

52
  return gas_internal_energy_from_entropy(p->rho_bar, p->entropy);
53
54
55
}

/**
56
 * @brief Returns the physical internal energy of a particle
57
58
 *
 * @param p The particle of interest
59
 * @param cosmo The cosmological model.
60
 */
61
62
63
64
65
66
67
68
69
70
71
72
73
74
__attribute__((always_inline)) INLINE static float
hydro_get_physical_internal_energy(const struct part *restrict p,
                                   const struct cosmology *cosmo) {

  return gas_internal_energy_from_entropy(p->rho_bar * cosmo->a3_inv,
                                          p->entropy);
}

/**
 * @brief Returns the comoving pressure of a particle
 *
 * @param p The particle of interest
 */
__attribute__((always_inline)) INLINE static float hydro_get_comoving_pressure(
75
    const struct part *restrict p) {
76

77
  return gas_pressure_from_entropy(p->rho_bar, p->entropy);
78
79
80
}

/**
81
82
83
84
85
86
87
88
89
90
91
92
 * @brief Returns the physical pressure of a particle
 *
 * @param p The particle of interest
 */
__attribute__((always_inline)) INLINE static float hydro_get_physical_pressure(
    const struct part *restrict p, const struct cosmology *cosmo) {

  return gas_pressure_from_entropy(p->rho_bar * cosmo->a3_inv, p->entropy);
}

/**
 * @brief Returns the comoving entropy of a particle
93
94
95
 *
 * @param p The particle of interest
 */
96
__attribute__((always_inline)) INLINE static float hydro_get_comoving_entropy(
97
    const struct part *restrict p) {
98

99
  return p->entropy;
100
101
102
}

/**
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
 * @brief Returns the physical entropy of a particle
 *
 * @param p The particle of interest.
 * @param cosmo The cosmological model.
 */
__attribute__((always_inline)) INLINE static float hydro_get_physical_entropy(
    const struct part *restrict p, const struct cosmology *cosmo) {

  /* Note: no cosmological conversion required here with our choice of
   * coordinates. */
  return p->entropy;
}

/**
 * @brief Returns the comoving sound speed of a particle
118
119
120
 *
 * @param p The particle of interest
 */
121
122
__attribute__((always_inline)) INLINE static float
hydro_get_comoving_soundspeed(const struct part *restrict p) {
123
124
125
126

  return p->force.soundspeed;
}

127
/**
128
 * @brief Returns the physical sound speed of a particle
129
130
 *
 * @param p The particle of interest
131
 * @param cosmo The cosmological model.
132
 */
133
134
135
136
137
138
139
140
141
142
143
144
145
__attribute__((always_inline)) INLINE static float
hydro_get_physical_soundspeed(const struct part *restrict p,
                              const struct cosmology *cosmo) {

  return cosmo->a_factor_sound_speed * p->force.soundspeed;
}

/**
 * @brief Returns the comoving density of a particle
 *
 * @param p The particle of interest
 */
__attribute__((always_inline)) INLINE static float hydro_get_comoving_density(
146
147
148
149
150
    const struct part *restrict p) {

  return p->rho;
}

151
152
153
154
155
156
157
158
159
160
161
/**
 * @brief Returns the physical density of a particle
 *
 * @param p The particle of interest
 */
__attribute__((always_inline)) INLINE static float hydro_get_physical_density(
    const struct part *restrict p, const struct cosmology *cosmo) {

  return p->rho * cosmo->a3_inv;
}

162
163
164
165
166
167
168
169
170
171
172
/**
 * @brief Returns the mass of a particle
 *
 * @param p The particle of interest
 */
__attribute__((always_inline)) INLINE static float hydro_get_mass(
    const struct part *restrict p) {

  return p->mass;
}

173
174
175
176
177
178
179
180
181
182
183
184
/**
 * @brief Sets the mass of a particle
 *
 * @param p The particle of interest
 * @param m The mass to set.
 */
__attribute__((always_inline)) INLINE static void hydro_set_mass(
    struct part *restrict p, float m) {

  p->mass = m;
}

185
186
187
188
189
/**
 * @brief Returns the velocities drifted to the current time of a particle.
 *
 * @param p The particle of interest
 * @param xp The extended data of the particle.
190
191
 * @param dt_kick_hydro The time (for hydro accelerations) since the last kick.
 * @param dt_kick_grav The time (for gravity accelerations) since the last kick.
192
193
194
 * @param v (return) The velocities at the current time.
 */
__attribute__((always_inline)) INLINE static void hydro_get_drifted_velocities(
195
196
197
198
199
200
201
202
203
    const struct part *restrict p, const struct xpart *xp, float dt_kick_hydro,
    float dt_kick_grav, float v[3]) {

  v[0] = xp->v_full[0] + p->a_hydro[0] * dt_kick_hydro +
         xp->a_grav[0] * dt_kick_grav;
  v[1] = xp->v_full[1] + p->a_hydro[1] * dt_kick_hydro +
         xp->a_grav[1] * dt_kick_grav;
  v[2] = xp->v_full[2] + p->a_hydro[2] * dt_kick_hydro +
         xp->a_grav[2] * dt_kick_grav;
204
205
}

206
/**
207
 * @brief Returns the time derivative of internal energy of a particle
208
 *
209
 * We assume a constant density.
210
 *
211
 * @param p The particle of interest
212
 */
213
214
__attribute__((always_inline)) INLINE static float hydro_get_internal_energy_dt(
    const struct part *restrict p) {
215

216
  return gas_internal_energy_from_entropy(p->rho_bar, p->entropy_dt);
217
218
219
}

/**
220
 * @brief Returns the time derivative of internal energy of a particle
221
 *
222
 * We assume a constant density.
223
 *
224
225
 * @param p The particle of interest.
 * @param du_dt The new time derivative of the internal energy.
226
 */
227
228
__attribute__((always_inline)) INLINE static void hydro_set_internal_energy_dt(
    struct part *restrict p, float du_dt) {
229

230
  p->entropy_dt = gas_entropy_from_internal_energy(p->rho_bar, du_dt);
231
232
233
234
235
}

/**
 * @brief Computes the hydro time-step of a given particle
 *
236
237
238
239
 * @param p Pointer to the particle data.
 * @param xp Pointer to the extended particle data.
 * @param hydro_properties The constants used in the scheme.
 * @param cosmo The current cosmological model.
240
241
242
 */
__attribute__((always_inline)) INLINE static float hydro_compute_timestep(
    const struct part *restrict p, const struct xpart *restrict xp,
243
244
    const struct hydro_props *restrict hydro_properties,
    const struct cosmology *restrict cosmo) {
245

246
  const float CFL = hydro_properties->CFL_condition;
247
248

  /* CFL condition */
249
250
  const float dt_cfl = 2.f * kernel_gamma * CFL * cosmo->a * p->h /
                       (cosmo->a_factor_sound_speed * p->force.v_sig);
251
252
253
254

  return dt_cfl;
}

255
256
257
258
259
260
261
262
263
264
/**
 * @brief Does some extra hydro operations once the actual physical time step
 * for the particle is known.
 *
 * @param p The particle to act upon.
 * @param dt Physical time step of the particle during the next step.
 */
__attribute__((always_inline)) INLINE static void hydro_timestep_extra(
    struct part *p, float dt) {}

265
266
267
268
269
270
271
/**
 * @brief Prepares a particle for the density calculation.
 *
 * Zeroes all the relevant arrays in preparation for the sums taking place in
 * the variaous density tasks
 *
 * @param p The particle to act upon
272
 * @param hs #hydro_space containing hydro specific space information.
273
274
 */
__attribute__((always_inline)) INLINE static void hydro_init_part(
275
    struct part *restrict p, const struct hydro_space *hs) {
276

277
  p->rho = 0.f;
278
  p->rho_bar = 0.f;
279
280
  p->density.wcount = 0.f;
  p->density.wcount_dh = 0.f;
281
282
  p->density.rho_dh = 0.f;
  p->density.pressure_dh = 0.f;
283

284
285
286
287
  p->density.div_v = 0.f;
  p->density.rot_v[0] = 0.f;
  p->density.rot_v[1] = 0.f;
  p->density.rot_v[2] = 0.f;
288
289
290
291
292
293
294
295
296
}

/**
 * @brief Finishes the density calculation.
 *
 * Multiplies the density and number of neighbours by the appropiate constants
 * and add the self-contribution term.
 *
 * @param p The particle to act upon
297
 * @param cosmo The current cosmological model.
298
299
 */
__attribute__((always_inline)) INLINE static void hydro_end_density(
300
    struct part *restrict p, const struct cosmology *cosmo) {
301
302
303

  /* Some smoothing length multiples. */
  const float h = p->h;
304
305
306
  const float h_inv = 1.0f / h;                       /* 1/h */
  const float h_inv_dim = pow_dimension(h_inv);       /* 1/h^d */
  const float h_inv_dim_plus_one = h_inv_dim * h_inv; /* 1/h^(d+1) */
307
308
309

  /* Final operation on the density (add self-contribution). */
  p->rho += p->mass * kernel_root;
310
  p->rho_bar += p->mass * p->entropy_one_over_gamma * kernel_root;
311
312
  p->density.rho_dh -= hydro_dimension * p->mass * kernel_root;
  p->density.pressure_dh -=
313
      hydro_dimension * p->mass * p->entropy_one_over_gamma * kernel_root;
314
  p->density.wcount += kernel_root;
315
  p->density.wcount_dh -= hydro_dimension * kernel_root;
316
317
318

  /* Finish the calculation by inserting the missing h-factors */
  p->rho *= h_inv_dim;
319
  p->rho_bar *= h_inv_dim;
320
321
  p->density.rho_dh *= h_inv_dim_plus_one;
  p->density.pressure_dh *= h_inv_dim_plus_one;
322
323
  p->density.wcount *= h_inv_dim;
  p->density.wcount_dh *= h_inv_dim_plus_one;
324

325
  const float rho_inv = 1.f / p->rho;
326
  const float a_inv2 = cosmo->a2_inv;
327
  const float entropy_minus_one_over_gamma = 1.f / p->entropy_one_over_gamma;
328

329
330
  /* Final operation on the weighted density */
  p->rho_bar *= entropy_minus_one_over_gamma;
331
332

  /* Finish calculation of the velocity curl components */
333
334
335
  p->density.rot_v[0] *= h_inv_dim_plus_one * a_inv2 * rho_inv;
  p->density.rot_v[1] *= h_inv_dim_plus_one * a_inv2 * rho_inv;
  p->density.rot_v[2] *= h_inv_dim_plus_one * a_inv2 * rho_inv;
336
337

  /* Finish calculation of the velocity divergence */
338
  p->density.div_v *= h_inv_dim_plus_one * rho_inv;
339
340
}

341
342
343
344
345
/**
 * @brief Sets all particle fields to sensible values when the #part has 0 ngbs.
 *
 * @param p The particle to act upon
 * @param xp The extended particle data to act upon
346
 * @param cosmo The current cosmological model.
347
348
 */
__attribute__((always_inline)) INLINE static void hydro_part_has_no_neighbours(
349
350
    struct part *restrict p, struct xpart *restrict xp,
    const struct cosmology *cosmo) {
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369

  /* Some smoothing length multiples. */
  const float h = p->h;
  const float h_inv = 1.0f / h;                 /* 1/h */
  const float h_inv_dim = pow_dimension(h_inv); /* 1/h^d */

  /* Re-set problematic values */
  p->rho = p->mass * kernel_root * h_inv_dim;
  p->rho_bar = p->mass * kernel_root * h_inv_dim;
  p->density.wcount = kernel_root * kernel_norm * h_inv_dim;
  p->density.rho_dh = 0.f;
  p->density.wcount_dh = 0.f;
  p->density.pressure_dh = 0.f;
  p->density.div_v = 0.f;
  p->density.rot_v[0] = 0.f;
  p->density.rot_v[1] = 0.f;
  p->density.rot_v[2] = 0.f;
}

370
371
372
373
374
375
376
/**
 * @brief Prepare a particle for the force calculation.
 *
 * Computes viscosity term, conduction term and smoothing length gradient terms.
 *
 * @param p The particle to act upon
 * @param xp The extended particle data to act upon
377
 * @param cosmo The current cosmological model.
378
379
 */
__attribute__((always_inline)) INLINE static void hydro_prepare_force(
380
381
    struct part *restrict p, struct xpart *restrict xp,
    const struct cosmology *cosmo) {
382

383
  const float fac_mu = cosmo->a_factor_mu;
384
385
386
387
388
389
390
391
392

  /* Compute the norm of the curl */
  const float curl_v = sqrtf(p->density.rot_v[0] * p->density.rot_v[0] +
                             p->density.rot_v[1] * p->density.rot_v[1] +
                             p->density.rot_v[2] * p->density.rot_v[2]);

  /* Compute the norm of div v */
  const float abs_div_v = fabsf(p->density.div_v);

393
  /* Compute the pressure */
394
  const float pressure = gas_pressure_from_entropy(p->rho_bar, p->entropy);
395

396
  /* Compute the sound speed from the pressure*/
397
398
  const float soundspeed = gas_soundspeed_from_pressure(p->rho_bar, pressure);

399
400
  /* Compute the Balsara switch */
  const float balsara =
401
      abs_div_v / (abs_div_v + curl_v + 0.0001f * soundspeed * fac_mu / p->h);
402

Matthieu Schaller's avatar
Matthieu Schaller committed
403
404
405
406
407
  /* Divide the pressure by the density squared to get the SPH term */
  const float rho_bar_inv = 1.f / p->rho_bar;
  const float P_over_rho2 = pressure * rho_bar_inv * rho_bar_inv;

  /* Compute "grad h" term (note we use rho here and not rho_bar !)*/
408
409
410
  const float rho_inv = 1.f / p->rho;
  const float rho_dh =
      1.f / (1.f + hydro_dimension_inv * p->h * p->density.rho_dh * rho_inv);
411
412
  const float pressure_dh =
      p->density.pressure_dh * rho_inv * p->h * hydro_dimension_inv;
413
414

  const float grad_h_term = rho_dh * pressure_dh;
415

416
417
  /* Update variables. */
  p->force.soundspeed = soundspeed;
418
  p->force.P_over_rho2 = P_over_rho2;
419
420
  p->force.balsara = balsara;
  p->force.f = grad_h_term;
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
}

/**
 * @brief Reset acceleration fields of a particle
 *
 * Resets all hydro acceleration and time derivative fields in preparation
 * for the sums taking place in the variaous force tasks
 *
 * @param p The particle to act upon
 */
__attribute__((always_inline)) INLINE static void hydro_reset_acceleration(
    struct part *restrict p) {

  /* Reset the acceleration. */
  p->a_hydro[0] = 0.0f;
  p->a_hydro[1] = 0.0f;
  p->a_hydro[2] = 0.0f;

  /* Reset the time derivatives. */
  p->entropy_dt = 0.0f;
  p->force.h_dt = 0.0f;

  /* Reset maximal signal velocity */
  p->force.v_sig = 0.0f;
}

447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
/**
 * @brief Sets the values to be predicted in the drifts to their values at a
 * kick time
 *
 * @param p The particle.
 * @param xp The extended data of this particle.
 */
__attribute__((always_inline)) INLINE static void hydro_reset_predicted_values(
    struct part *restrict p, const struct xpart *restrict xp) {

  /* Re-set the predicted velocities */
  p->v[0] = xp->v_full[0];
  p->v[1] = xp->v_full[1];
  p->v[2] = xp->v_full[2];

  /* Re-set the entropy */
  p->entropy = xp->entropy_full;
}

466
467
468
469
470
/**
 * @brief Predict additional particle fields forward in time when drifting
 *
 * @param p The particle
 * @param xp The extended data of the particle
471
472
 * @param dt_drift The drift time-step for positions.
 * @param dt_therm The drift time-step for thermal quantities.
473
474
 */
__attribute__((always_inline)) INLINE static void hydro_predict_extra(
475
476
    struct part *restrict p, const struct xpart *restrict xp, float dt_drift,
    float dt_therm) {
477

478
  const float h_inv = 1.f / p->h;
479

480
  /* Predict smoothing length */
481
  const float w1 = p->force.h_dt * h_inv * dt_drift;
482
  if (fabsf(w1) < 0.2f)
483
    p->h *= approx_expf(w1); /* 4th order expansion of exp(w) */
484
  else
485
486
487
488
489
490
491
492
493
494
495
496
    p->h *= expf(w1);

  /* Predict density */
  const float w2 = -hydro_dimension * w1;
  if (fabsf(w2) < 0.2f) {
    p->rho *= approx_expf(w2); /* 4th order expansion of exp(w) */
    p->rho_bar *= approx_expf(w2);
  } else {
    p->rho *= expf(w2);
    p->rho_bar *= expf(w2);
  }

497
  /* Predict the entropy */
498
  p->entropy += p->entropy_dt * dt_therm;
499
500

  /* Compute the pressure */
501
  const float pressure = gas_pressure_from_entropy(p->rho_bar, p->entropy);
502

503
504
505
506
  /* Compute the new sound speed */
  const float soundspeed = gas_soundspeed_from_pressure(p->rho_bar, pressure);

  /* Divide the pressure by the density squared to get the SPH term */
507
  const float rho_bar_inv = 1.f / p->rho_bar;
508
  const float P_over_rho2 = pressure * rho_bar_inv * rho_bar_inv;
509

510
  /* Update the variables */
511
  p->entropy_one_over_gamma = pow_one_over_gamma(p->entropy);
512
  p->force.soundspeed = soundspeed;
513
  p->force.P_over_rho2 = P_over_rho2;
514
515
516
517
518
519
520
521
}

/**
 * @brief Finishes the force calculation.
 *
 * Multiplies the forces and accelerationsby the appropiate constants
 *
 * @param p The particle to act upon
522
 * @param cosmo The current cosmological model.
523
524
 */
__attribute__((always_inline)) INLINE static void hydro_end_force(
525
    struct part *restrict p, const struct cosmology *cosmo) {
526
527
528

  p->force.h_dt *= p->h * hydro_dimension_inv;

529
530
  p->entropy_dt = 0.5f * cosmo->a2_inv *
                  gas_entropy_from_internal_energy(p->rho_bar, p->entropy_dt);
531
532
533
534
535
536
537
}

/**
 * @brief Kick the additional variables
 *
 * @param p The particle to act upon
 * @param xp The particle extended data to act upon
538
 * @param dt_therm The time-step for this kick (for thermodynamic quantities)
539
540
 */
__attribute__((always_inline)) INLINE static void hydro_kick_extra(
541
    struct part *restrict p, struct xpart *restrict xp, float dt_therm) {
542
543

  /* Do not decrease the entropy (temperature) by more than a factor of 2*/
544
545
  if (dt_therm > 0. && p->entropy_dt * dt_therm < -0.5f * xp->entropy_full) {
    p->entropy_dt = -0.5f * xp->entropy_full / dt_therm;
546
  }
547
  xp->entropy_full += p->entropy_dt * dt_therm;
548
549

  /* Compute the pressure */
550
551
  const float pressure =
      gas_pressure_from_entropy(p->rho_bar, xp->entropy_full);
552

553
554
555
556
  /* Compute the new sound speed */
  const float soundspeed = gas_soundspeed_from_pressure(p->rho_bar, pressure);

  /* Divide the pressure by the density squared to get the SPH term */
557
  const float rho_bar_inv = 1.f / p->rho_bar;
558
  const float P_over_rho2 = pressure * rho_bar_inv * rho_bar_inv;
559
560
561

  p->entropy_one_over_gamma = pow_one_over_gamma(p->entropy);
  p->force.soundspeed = soundspeed;
562
  p->force.P_over_rho2 = P_over_rho2;
563
564
565
566
567
568
569
570
}

/**
 *  @brief Converts hydro quantity of a particle at the start of a run
 *
 * Requires the density to be known
 *
 * @param p The particle to act upon
571
572
 * @param xp The extended data.
 * @param cosmo The cosmological model.
573
574
 */
__attribute__((always_inline)) INLINE static void hydro_convert_quantities(
575
576
    struct part *restrict p, struct xpart *restrict xp,
    const struct cosmology *cosmo) {
577
578

  /* We read u in the entropy field. We now get S from u */
579
580
  xp->entropy_full =
      gas_entropy_from_internal_energy(p->rho_bar * cosmo->a3_inv, p->entropy);
581
  p->entropy = xp->entropy_full;
582
583
584
  p->entropy_one_over_gamma = pow_one_over_gamma(p->entropy);

  /* Compute the pressure */
585
  const float pressure = gas_pressure_from_entropy(p->rho_bar, p->entropy);
586

587
588
589
590
  /* Compute the sound speed */
  const float soundspeed = gas_soundspeed_from_pressure(p->rho_bar, pressure);

  /* Divide the pressure by the density squared to get the SPH term */
591
  const float rho_bar_inv = 1.f / p->rho_bar;
592
593
  const float P_over_rho2 = pressure * rho_bar_inv * rho_bar_inv;

594
  p->force.soundspeed = soundspeed;
595
  p->force.P_over_rho2 = P_over_rho2;
596
597
}

598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
/**
 * @brief Initialises the particles for the first time
 *
 * This function is called only once just after the ICs have been
 * read in to do some conversions.
 *
 * @param p The particle to act upon
 * @param xp The extended particle data to act upon
 */
__attribute__((always_inline)) INLINE static void hydro_first_init_part(
    struct part *restrict p, struct xpart *restrict xp) {

  p->time_bin = 0;
  p->rho_bar = 0.f;
  p->entropy_one_over_gamma = pow_one_over_gamma(p->entropy);
  xp->v_full[0] = p->v[0];
  xp->v_full[1] = p->v[1];
  xp->v_full[2] = p->v[2];
616
617
618
  xp->a_grav[0] = 0.f;
  xp->a_grav[1] = 0.f;
  xp->a_grav[2] = 0.f;
619
620

  hydro_reset_acceleration(p);
621
  hydro_init_part(p, NULL);
622
623
}

624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
/**
 * @brief Overwrite the initial internal energy of a particle.
 *
 * Note that in the cases where the thermodynamic variable is not
 * internal energy but gets converted later, we must overwrite that
 * field. The conversion to the actual variable happens later after
 * the initial fake time-step.
 *
 * @param p The #part to write to.
 * @param u_init The new initial internal energy.
 */
__attribute__((always_inline)) INLINE static void
hydro_set_init_internal_energy(struct part *p, float u_init) {

  p->entropy = u_init;
}

641
#endif /* SWIFT_PRESSURE_ENTROPY_HYDRO_H */