MP-Gadget  5.0.1.dev1-76bc7d4726-dirty
Classes | Macros | Functions
test_timefac.c File Reference
#include <stdarg.h>
#include <stddef.h>
#include <setjmp.h>
#include <cmocka.h>
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include <gsl/gsl_integration.h>
#include <stdint.h>
#include "stub.h"
#include <libgadget/timefac.h>
Include dependency graph for test_timefac.c:

Go to the source code of this file.

Classes

struct  fac_params
 

Macros

#define AMIN   0.005
 
#define AMAX   1.0
 
#define LTIMEBINS   29
 

Functions

double hubble_function (const Cosmology *CP, double a)
 
double fac_integ (double a, void *param)
 
double loga_from_ti (int ti)
 
static int get_ti (double aa)
 
double exact_drift_factor (Cosmology *CP, double a1, double a2, int exp)
 
void test_drift_factor (void **state)
 
int main (void)
 

Macro Definition Documentation

◆ AMAX

#define AMAX   1.0

Definition at line 17 of file test_timefac.c.

◆ AMIN

#define AMIN   0.005

Definition at line 16 of file test_timefac.c.

◆ LTIMEBINS

#define LTIMEBINS   29

Definition at line 18 of file test_timefac.c.

Function Documentation

◆ exact_drift_factor()

double exact_drift_factor ( Cosmology CP,
double  a1,
double  a2,
int  exp 
)

Definition at line 64 of file test_timefac.c.

65 {
66  double result, abserr;
67  gsl_function F;
68  gsl_integration_workspace *workspace;
69  workspace = gsl_integration_workspace_alloc(10000);
70  F.function = &fac_integ;
71  struct fac_params ff = {CP, exp};
72  F.params = &ff;
73  gsl_integration_qag(&F, a1,a2, 0, 1.0e-8, 10000, GSL_INTEG_GAUSS61, workspace, &result, &abserr);
74  gsl_integration_workspace_free(workspace);
75  return result;
76 }
static Cosmology * CP
Definition: power.c:27
double fac_integ(double a, void *param)
Definition: test_timefac.c:40

References CP, fac_params::exp, and fac_integ().

Referenced by test_drift_factor().

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

◆ fac_integ()

double fac_integ ( double  a,
void *  param 
)

Definition at line 40 of file test_timefac.c.

41 {
42  double h;
43  struct fac_params * ff = (struct fac_params *) param;
44 
45  h = hubble_function(ff->CP, a);
46 
47  return 1 / (h * pow(a,ff->exp));
48 }
Cosmology * CP
Definition: test_timefac.c:36
double hubble_function(const Cosmology *CP, double a)
Definition: test_timefac.c:21

References fac_params::CP, fac_params::exp, and hubble_function().

Referenced by exact_drift_factor().

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

◆ get_ti()

static int get_ti ( double  aa)
inlinestatic

Definition at line 58 of file test_timefac.c.

59 {
60  double logDTime = (log(AMAX) - log(AMIN)) / (1 << LTIMEBINS);
61  return (log(aa) - log(AMIN))/logDTime;
62 }
#define AMIN
Definition: test_timefac.c:16
#define AMAX
Definition: test_timefac.c:17
#define LTIMEBINS
Definition: test_timefac.c:18

References AMAX, AMIN, and LTIMEBINS.

Referenced by test_drift_factor().

Here is the caller graph for this function:

◆ hubble_function()

double hubble_function ( const Cosmology CP,
double  a 
)

Definition at line 21 of file test_timefac.c.

22 {
23  double hubble_a;
24  /* lambda and matter */
25  hubble_a = 1 - CP->Omega0;
26  hubble_a += CP->Omega0 / (a * a * a);
27  /*Radiation*/
28  hubble_a += 5.045e-5*(1+7./8.*pow(pow(4/11.,1/3.)*1.00328,4)*3) / (a * a * a * a);
29  /* Hard-code default Gadget unit system. */
30  hubble_a = 0.1 * sqrt(hubble_a);
31  return (hubble_a);
32 }
double Omega0
Definition: cosmology.h:10

References CP, and Cosmology::Omega0.

Referenced by fac_integ().

Here is the caller graph for this function:

◆ loga_from_ti()

double loga_from_ti ( int  ti)

Definition at line 51 of file test_timefac.c.

52 {
53  double logDTime = (log(AMAX) - log(AMIN)) / (1 << LTIMEBINS);
54  return log(AMIN) + ti * logDTime;
55 }

References AMAX, AMIN, and LTIMEBINS.

Referenced by cooling_direct(), density(), get_atime(), get_exact_factor(), and test_conversions().

Here is the caller graph for this function:

◆ main()

int main ( void  )

Definition at line 109 of file test_timefac.c.

109  {
110  const struct CMUnitTest tests[] = {
111  cmocka_unit_test(test_drift_factor),
112  };
113  return cmocka_run_group_tests_mpi(tests, NULL, NULL);
114 }
void test_drift_factor(void **state)
Definition: test_timefac.c:78

References test_drift_factor().

Here is the call graph for this function:

◆ test_drift_factor()

void test_drift_factor ( void **  state)

Definition at line 78 of file test_timefac.c.

79 {
80  /*Initialise the table: default values from z=200 to z=0*/
81  Cosmology CP;
82  CP.Omega0 = 1.;
83  /* Check default scaling: for total matter domination
84  * we should have a drift factor like 1/sqrt(a)*/
85  assert_true(fabs(get_exact_drift_factor(&CP, get_ti(0.8), get_ti(0.85)) + 2/0.1*(1/sqrt(0.85) - 1/sqrt(0.8))) < 5e-5);
86  /*Test the kick table*/
87  assert_true(fabs(get_exact_gravkick_factor(&CP, get_ti(0.8), get_ti(0.85)) - 2/0.1*(sqrt(0.85) - sqrt(0.8))) < 5e-5);
88 
89  //Chosen so we get the same bin
90  assert_true(fabs(get_exact_drift_factor(&CP, get_ti(0.8), get_ti(0.8003)) + 2/0.1*(1/sqrt(0.8003) - 1/sqrt(0.8))) < 5e-6);
91  //Now choose a more realistic cosmology
92  CP.Omega0 = 0.25;
93  /*Check late and early times*/
94  assert_true(fabs(get_exact_drift_factor(&CP, get_ti(0.95), get_ti(0.98)) - exact_drift_factor(&CP, 0.95, 0.98,3)) < 5e-5);
95  assert_true(fabs(get_exact_drift_factor(&CP, get_ti(0.05), get_ti(0.06)) - exact_drift_factor(&CP, 0.05, 0.06,3)) < 5e-5);
96  /*Check boundary conditions*/
97  double logDtime = (log(AMAX)-log(AMIN))/(1<<LTIMEBINS);
98  assert_true(fabs(get_exact_drift_factor(&CP, ((1<<LTIMEBINS)-1), 1<<LTIMEBINS) - exact_drift_factor(&CP, AMAX-logDtime, AMAX,3)) < 5e-5);
99  assert_true(fabs(get_exact_drift_factor(&CP, 0, 1) - exact_drift_factor(&CP, 1.0 - exp(log(AMAX)-log(AMIN))/(1<<LTIMEBINS), 1.0,3)) < 5e-5);
100  /*Gravkick*/
101  assert_true(fabs(get_exact_gravkick_factor(&CP, get_ti(0.8), get_ti(0.85)) - exact_drift_factor(&CP, 0.8, 0.85, 2)) < 5e-5);
102  assert_true(fabs(get_exact_gravkick_factor(&CP, get_ti(0.05), get_ti(0.06)) - exact_drift_factor(&CP, 0.05, 0.06, 2)) < 5e-5);
103 
104  /*Test the hydrokick table: always the same as drift*/
105  assert_true(fabs(get_exact_hydrokick_factor(&CP, get_ti(0.8), get_ti(0.85)) - get_exact_drift_factor(&CP, get_ti(0.8), get_ti(0.85))) < 5e-5);
106 
107 }
double exact_drift_factor(Cosmology *CP, double a1, double a2, int exp)
Definition: test_timefac.c:64
static int get_ti(double aa)
Definition: test_timefac.c:58
double get_exact_hydrokick_factor(Cosmology *CP, inttime_t ti0, inttime_t ti1)
Definition: timefac.c:75
double get_exact_drift_factor(Cosmology *CP, inttime_t ti0, inttime_t ti1)
Definition: timefac.c:64
double get_exact_gravkick_factor(Cosmology *CP, inttime_t ti0, inttime_t ti1)
Definition: timefac.c:70

References AMAX, AMIN, CP, exact_drift_factor(), fac_params::exp, get_exact_drift_factor(), get_exact_gravkick_factor(), get_exact_hydrokick_factor(), get_ti(), LTIMEBINS, and Cosmology::Omega0.

Referenced by main().

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