MP-Gadget  5.0.1.dev1-76bc7d4726-dirty
density.h
Go to the documentation of this file.
1 #ifndef DENSITY_H
2 #define DENSITY_H
3 
4 #include "forcetree.h"
5 #include "timestep.h"
6 #include "densitykernel.h"
7 #include "utils/paramset.h"
8 
10 {
14  /* These are for black hole neighbour finding and so belong in the density module, not the black hole module.*/
17 
18  enum DensityKernelType DensityKernelType; /* 0 for Cubic Spline, (recmd NumNgb = 33)
19  1 for Quintic spline (recmd NumNgb = 97) */
20 
23 };
24 
26 {
29  /* VelPred can always be derived from the current time and acceleration.
30  * However, doing so makes the SPH and hydro code much (a factor of two)
31  * slower. It requires computing get_gravkick_factor twice with different arguments,
32  * which defeats the lookup cache in timefac.c. Because VelPred is used multiple times,
33  * it is much quicker to compute it once and re-use this*/
35 };
36 
37 /*Set the parameters of the density module*/
39 /*Set cooling module parameters from a density_params struct for the tests*/
40 void set_densitypar(struct density_params dp);
41 
42 /* This routine computes the particle densities. If update_hsml is true
43  * it runs multiple times, changing the smoothing length until
44  * there are enough neighbours. If update_hsml is false (when initializing the EgyWtDensity)
45  * it just computes densities.
46  * If DoEgyDensity is true it also computes the entropy-weighted density for
47  * pressure-entropy SPH. */
48 void density(const ActiveParticles * act, int update_hsml, int DoEgyDensity, int BlackHoleOn, double MinEgySpec, const DriftKickTimes times, Cosmology * CP, struct sph_pred_data * SPH_predicted, MyFloat * GradRho, const ForceTree * const tree);
49 
50 /* Get the desired nuber of neighbours for the supplied kernel*/
51 double GetNumNgb(enum DensityKernelType KernelType);
52 
53 /* Get the current density kernel type*/
55 
57 void slots_free_sph_pred_data(struct sph_pred_data * sph_pred);
58 
59 /* Predicted quantity computation used in hydro*/
60 MyFloat SPH_EntVarPred(int PI, double MinEgySpec, double a3inv, double dloga);
61 void SPH_VelPred(int i, MyFloat * VelPred, const double FgravkickB, double gravkick, double hydrokick);
62 
63 
64 #endif
void SPH_VelPred(int i, MyFloat *VelPred, const double FgravkickB, double gravkick, double hydrokick)
Definition: density.c:89
double GetNumNgb(enum DensityKernelType KernelType)
Definition: density.c:55
MyFloat SPH_EntVarPred(int PI, double MinEgySpec, double a3inv, double dloga)
Definition: density.c:71
void slots_free_sph_pred_data(struct sph_pred_data *sph_pred)
Definition: density.c:665
void set_densitypar(struct density_params dp)
Definition: density.c:25
struct sph_pred_data slots_allocate_sph_pred_data(int nsph)
Definition: density.c:654
void set_density_params(ParameterSet *ps)
Definition: density.c:32
enum DensityKernelType GetDensityKernelType(void)
Definition: density.c:63
void density(const ActiveParticles *act, int update_hsml, int DoEgyDensity, int BlackHoleOn, double MinEgySpec, const DriftKickTimes times, Cosmology *CP, struct sph_pred_data *SPH_predicted, MyFloat *GradRho, const ForceTree *const tree)
Definition: density.c:203
DensityKernelType
Definition: densitykernel.h:17
static double dloga
Definition: lightcone.c:21
static Cosmology * CP
Definition: power.c:27
double MaxNumNgbDeviation
Definition: density.h:12
double DensityResolutionEta
Definition: density.h:11
double MinGasHsmlFractional
Definition: density.h:22
double BlackHoleNgbFactor
Definition: density.h:15
enum DensityKernelType DensityKernelType
Definition: density.h:18
double BlackHoleMaxAccretionRadius
Definition: density.h:16
MyFloat * EntVarPred
Definition: density.h:28
MyFloat * VelPred
Definition: density.h:34
LOW_PRECISION MyFloat
Definition: types.h:19