MP-Gadget  5.0.1.dev1-76bc7d4726-dirty
Macros | Functions
timefac.c File Reference
#include <mpi.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <math.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_integration.h>
#include "physconst.h"
#include "timefac.h"
#include "timebinmgr.h"
#include "utils.h"
Include dependency graph for timefac.c:

Go to the source code of this file.

Macros

#define WORKSIZE   10000
 

Functions

static double drift_integ (double a, void *param)
 
static double gravkick_integ (double a, void *param)
 
static double hydrokick_integ (double a, void *param)
 
static double get_exact_factor (Cosmology *CP, inttime_t t0, inttime_t t1, double(*factor)(double, void *))
 
double get_exact_drift_factor (Cosmology *CP, inttime_t ti0, inttime_t ti1)
 
double get_exact_gravkick_factor (Cosmology *CP, inttime_t ti0, inttime_t ti1)
 
double get_exact_hydrokick_factor (Cosmology *CP, inttime_t ti0, inttime_t ti1)
 

Macro Definition Documentation

◆ WORKSIZE

#define WORKSIZE   10000

Definition at line 14 of file timefac.c.

Function Documentation

◆ drift_integ()

static double drift_integ ( double  a,
void *  param 
)
static

Definition at line 17 of file timefac.c.

18 {
19  Cosmology * CP = (Cosmology *) param;
20  double h = hubble_function(CP, a);
21  return 1 / (h * a * a * a);
22 }
double hubble_function(const Cosmology *CP, double a)
Definition: cosmology.c:58
static Cosmology * CP
Definition: power.c:27

References CP, and hubble_function().

Referenced by get_exact_drift_factor().

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

◆ get_exact_drift_factor()

double get_exact_drift_factor ( Cosmology CP,
inttime_t  ti0,
inttime_t  ti1 
)

Definition at line 64 of file timefac.c.

65 {
66  return get_exact_factor(CP, ti0, ti1, &drift_integ);
67 }
static double drift_integ(double a, void *param)
Definition: timefac.c:17
static double get_exact_factor(Cosmology *CP, inttime_t t0, inttime_t t1, double(*factor)(double, void *))
Definition: timefac.c:46

References CP, drift_integ(), and get_exact_factor().

Referenced by domain_build_exchange_list(), drift_all_particles(), hydro_force(), lightcone_compute(), and test_drift_factor().

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

◆ get_exact_factor()

static double get_exact_factor ( Cosmology CP,
inttime_t  t0,
inttime_t  t1,
double(*)(double, void *)  factor 
)
static

Definition at line 46 of file timefac.c.

47 {
48  double result, abserr;
49  if(t0 == t1)
50  return 0;
51  double a0 = exp(loga_from_ti(t0));
52  double a1 = exp(loga_from_ti(t1));
53  gsl_function F;
54  gsl_integration_workspace *workspace;
55  workspace = gsl_integration_workspace_alloc(WORKSIZE);
56  F.function = factor;
57  F.params = CP;
58  gsl_integration_qag(&F, a0, a1, 0, 1.0e-8, WORKSIZE, GSL_INTEG_GAUSS61, workspace, &result, &abserr);
59  gsl_integration_workspace_free(workspace);
60  return result;
61 }
double loga_from_ti(int ti)
Definition: test_timefac.c:51
#define WORKSIZE
Definition: timefac.c:14

References CP, loga_from_ti(), and WORKSIZE.

Referenced by get_exact_drift_factor(), get_exact_gravkick_factor(), and get_exact_hydrokick_factor().

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

◆ get_exact_gravkick_factor()

double get_exact_gravkick_factor ( Cosmology CP,
inttime_t  ti0,
inttime_t  ti1 
)

Definition at line 70 of file timefac.c.

71 {
72  return get_exact_factor(CP, ti0, ti1, &gravkick_integ);
73 }
static double gravkick_integ(double a, void *param)
Definition: timefac.c:25

References CP, get_exact_factor(), and gravkick_integ().

Referenced by apply_half_kick(), apply_PM_half_kick(), density(), hydro_force(), and test_drift_factor().

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

◆ get_exact_hydrokick_factor()

double get_exact_hydrokick_factor ( Cosmology CP,
inttime_t  ti0,
inttime_t  ti1 
)

Definition at line 75 of file timefac.c.

76 {
77  return get_exact_factor(CP, ti0, ti1, &hydrokick_integ);
78 }
static double hydrokick_integ(double a, void *param)
Definition: timefac.c:35

References CP, get_exact_factor(), and hydrokick_integ().

Referenced by apply_half_kick(), density(), hydro_force(), and test_drift_factor().

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

◆ gravkick_integ()

static double gravkick_integ ( double  a,
void *  param 
)
static

Definition at line 25 of file timefac.c.

26 {
27  Cosmology * CP = (Cosmology *) param;
28  double h = hubble_function(CP, a);
29 
30  return 1 / (h * a * a);
31 }

References CP, and hubble_function().

Referenced by get_exact_gravkick_factor().

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

◆ hydrokick_integ()

static double hydrokick_integ ( double  a,
void *  param 
)
static

Definition at line 35 of file timefac.c.

36 {
37  double h;
38 
39  Cosmology * CP = (Cosmology *) param;
40  h = hubble_function(CP, a);
41 
42  return 1 / (h * pow(a, 3 * GAMMA_MINUS1) * a);
43 }
#define GAMMA_MINUS1
Definition: physconst.h:35

References CP, GAMMA_MINUS1, and hubble_function().

Referenced by get_exact_hydrokick_factor().

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