MP-Gadget  5.0.1.dev1-76bc7d4726-dirty
Classes | Enumerations | Functions
power.h File Reference
#include <libgadget/cosmology.h>
#include <bigfile-mpi.h>
Include dependency graph for power.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

struct  power_params
 

Enumerations

enum  TransferType {
  DELTA_BAR = 0 , DELTA_CDM = 1 , DELTA_NU = 2 , DELTA_CB = 3 ,
  VEL_BAR = 4 , VEL_CDM = 5 , VEL_NU = 6 , VEL_CB = 7 ,
  VEL_TOT = 8 , DELTA_TOT = 9
}
 

Functions

double DeltaSpec (double kmag, enum TransferType Type)
 
double dlogGrowth (double kmag, enum TransferType Type)
 
int init_powerspectrum (int ThisTask, double InitTime, double UnitLength_in_cm_in, Cosmology *CPin, struct power_params *ppar)
 
void save_all_transfer_tables (BigFile *bf, int ThisTask)
 

Enumeration Type Documentation

◆ TransferType

Enumerator
DELTA_BAR 
DELTA_CDM 
DELTA_NU 
DELTA_CB 
VEL_BAR 
VEL_CDM 
VEL_NU 
VEL_CB 
VEL_TOT 
DELTA_TOT 

Definition at line 40 of file power.h.

41 {
42  DELTA_BAR = 0,
43  DELTA_CDM = 1,
44  DELTA_NU = 2,
45  DELTA_CB = 3,
46  VEL_BAR = 4,
47  VEL_CDM = 5,
48  VEL_NU = 6,
49  VEL_CB = 7,
50  VEL_TOT = 8,
51  /*Always unity, so there is no memory for it*/
52  DELTA_TOT = 9,
53 };
@ DELTA_CB
Definition: power.h:45
@ DELTA_BAR
Definition: power.h:42
@ DELTA_TOT
Definition: power.h:52
@ VEL_TOT
Definition: power.h:50
@ VEL_BAR
Definition: power.h:46
@ VEL_CB
Definition: power.h:49
@ DELTA_CDM
Definition: power.h:43
@ VEL_CDM
Definition: power.h:47
@ DELTA_NU
Definition: power.h:44
@ VEL_NU
Definition: power.h:48

Function Documentation

◆ DeltaSpec()

double DeltaSpec ( double  kmag,
enum TransferType  Type 
)

Definition at line 52 of file power.c.

53 {
54  double power;
55 
56  if(WhichSpectrum == 2)
57  power = Delta_Tabulated(k, Type);
58  else
59  power = Delta_EH(k);
60 
61  /*Normalise the power spectrum*/
62  power *= Norm;
63 
64  return power;
65 }
static int WhichSpectrum
Definition: power.c:23
static double Delta_Tabulated(double k, enum TransferType Type)
Definition: power.c:93
static double Delta_EH(double k)
Definition: power.c:432
static double Norm
Definition: power.c:22

References Delta_EH(), Delta_Tabulated(), Norm, and WhichSpectrum.

Referenced by density_transfer(), disp_transfer(), print_spec(), sigma2_int(), test_growth_numerical(), test_read_no_rescale(), and test_read_rescale_sigma8().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ dlogGrowth()

double dlogGrowth ( double  kmag,
enum TransferType  Type 
)

Definition at line 101 of file power.c.

102 {
103  /*Default to total growth: type 3 is cdm + baryons.*/
104  if(Type < DELTA_BAR || Type > DELTA_CB)
105  Type = VEL_TOT;
106  else
107  /*Type should be an offset from the first velocity*/
108  Type = (enum TransferType) ((int) VEL_BAR + ((int) Type - (int) DELTA_BAR));
109  return get_Tabulated(kmag, Type, 1);
110 }
static double get_Tabulated(double k, enum TransferType Type, double oobval)
Definition: power.c:69
TransferType
Definition: power.h:41

References DELTA_BAR, DELTA_CB, get_Tabulated(), VEL_BAR, and VEL_TOT.

Referenced by disp_transfer(), and test_growth_numerical().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ init_powerspectrum()

int init_powerspectrum ( int  ThisTask,
double  InitTime,
double  UnitLength_in_cm_in,
Cosmology CPin,
struct power_params ppar 
)

Definition at line 394 of file power.c.

395 {
397  /*Used only for tk_eh*/
399  UnitLength_in_cm = UnitLength_in_cm_in;
400  CP = CPin;
401 
402  if(ppar->WhichSpectrum == 2) {
404  /*Initialise the interpolation*/
407  if(ppar->DifferentTransferFunctions || ppar->ScaleDepVelocity) {
408  init_transfer_table(ThisTask, InitTime, ppar);
409  }
410  }
411 
412  Norm = 1.0;
413  if(ppar->InputPowerRedshift >= 0 || ppar->Sigma8 > 0) {
414  const double R8 = 8 * (CM_PER_MPC / UnitLength_in_cm); /* 8 Mpc/h */
415  if(ppar->Sigma8 > 0) {
416  double res = TopHatSigma2(R8);
417  if(isfinite(res) && res > 0)
418  Norm = ppar->Sigma8 / sqrt(res);
419  else
420  endrun(1, "Could not normalize P(k) to Sigma8=%g! Measured Sigma8^2 is %g\n", ppar->Sigma8, res);
421  }
422  double Dplus = GrowthFactor(CP, InitTime, 1/(1+ppar->InputPowerRedshift));
423  if(ppar->InputPowerRedshift >= 0) {
424  Norm *= Dplus;
425  message(0,"Growth factor from z=%g (InputPowerRedshift) to z=%g (Init): %g \n", ppar->InputPowerRedshift, 1. / InitTime - 1, Dplus);
426  }
427  message(0, "Normalization adjusted to Sigma8=%g (at z=0) (Normfac=%g). \n", sqrt(TopHatSigma2(R8))/Dplus, Norm);
428  }
429  return power_table.Nentry;
430 }
double GrowthFactor(Cosmology *CP, double astart, double aend)
Definition: cosmology.c:82
void message(int where, const char *fmt,...)
Definition: endrun.c:175
void endrun(int where, const char *fmt,...)
Definition: endrun.c:147
#define CM_PER_MPC
Definition: physconst.h:20
int init_transfer_table(int ThisTask, double InitTime, const struct power_params *const ppar)
Definition: power.c:322
static double PrimordialIndex
Definition: power.c:25
static struct table power_table
Definition: power.c:46
static struct table transfer_table
Definition: power.c:48
void read_power_table(int ThisTask, const char *inputfile, const int ncols, struct table *out_tab, const double InitTime, _parse_fn parse_line)
Definition: power.c:238
static Cosmology * CP
Definition: power.c:27
void parse_power(int i, double k, char *line, struct table *out_tab, int *InputInLog10, const double InitTime, int NumCol)
Definition: power.c:158
static double UnitLength_in_cm
Definition: power.c:26
static double TopHatSigma2(double R)
Definition: power.c:470
int ScaleDepVelocity
Definition: power.h:11
int DifferentTransferFunctions
Definition: power.h:10
double PrimordialIndex
Definition: power.h:16
int WhichSpectrum
Definition: power.h:9
double Sigma8
Definition: power.h:14
char * FileWithInputSpectrum
Definition: power.h:13
double InputPowerRedshift
Definition: power.h:15
int Nentry
Definition: power.c:36
double * logD[MAXCOLS]
Definition: power.c:38
double * logk
Definition: power.c:37
gsl_interp * mat_intp[MAXCOLS]
Definition: power.c:39
int ThisTask
Definition: test_exchange.c:23

References CM_PER_MPC, CP, power_params::DifferentTransferFunctions, endrun(), power_params::FileWithInputSpectrum, GrowthFactor(), init_transfer_table(), power_params::InputPowerRedshift, table::logD, table::logk, table::mat_intp, message(), table::Nentry, Norm, parse_power(), power_table, PrimordialIndex, power_params::PrimordialIndex, read_power_table(), power_params::ScaleDepVelocity, power_params::Sigma8, ThisTask, TopHatSigma2(), transfer_table, UnitLength_in_cm, WhichSpectrum, and power_params::WhichSpectrum.

Referenced by main(), test_growth_numerical(), test_read_no_rescale(), and test_read_rescale_sigma8().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ save_all_transfer_tables()

void save_all_transfer_tables ( BigFile *  bf,
int  ThisTask 
)

Definition at line 150 of file power.c.

151 {
152  const char * pname = "DELTA_MAT";
153  save_transfer(bf, 1, &power_table, "ICPower", ThisTask, &pname);
154  save_transfer(bf, MAXCOLS, &transfer_table, "ICTransfers", ThisTask, tnames);
155 }
#define MAXCOLS
Definition: power.c:32
static const char * tnames[MAXCOLS]
Definition: power.c:50
static void save_transfer(BigFile *bf, int ncol, struct table *ttable, const char *bname, int ThisTask, const char *colnames[])
Definition: power.c:113

References MAXCOLS, power_table, save_transfer(), ThisTask, tnames, and transfer_table.

Referenced by main().

Here is the call graph for this function:
Here is the caller graph for this function: