File: Omega_D_effWEC.c

package info (click to toggle)
inhomog 0.1.9.2-1
  • links: PTS, VCS
  • area: main
  • in suites: bookworm, bullseye, buster, sid, trixie
  • size: 2,476 kB
  • sloc: ansic: 10,454; sh: 4,380; makefile: 173
file content (434 lines) | stat: -rw-r--r-- 16,205 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
/*
   Omega_D_effWEC for weak energy condition - RZA2 arXiv:1303.6193v2 - near clone of test_Omega_D

   Copyright (C) 2014-2015 Boud Roukema, Jan Ostrowski

   This program is free software; you can redistribute it and/or modify
   it under the terms of the GNU General Public License as published by
   the Free Software Foundation; either version 2, 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 General Public License
   along with this program; if not, write to the Free Software Foundation,
   Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.

   See also http://www.gnu.org/licenses/gpl.html

*/

#include <stdio.h>
#include <sys/types.h>
#include "config.h"
#include <math.h>

#include <gsl/gsl_rng.h>
#include <gsl/gsl_errno.h>
#include <gsl/gsl_statistics.h>

/* for malloc_usable_size if available */
#ifdef __GNUC__
#include <malloc.h>
#endif

#include "lib/inhomog.h"

/* int test_Omega_D(
                 struct rza_integrand_params_s rza_integrand_params,
*/
int main(void)
{
  struct rza_integrand_params_s rza_integrand_params;
  struct rza_integrand_params_s rza_integrand_params_sigma8;
  struct background_cosm_params_s background_cosm_params;

  int want_verbose = 0;

#define N_PLOT_T 3
#define N_PLOT_R 256
#define N_COSM_I 0
#define N_COSM 2
  double R_domain_mod[N_PLOT_R];

  /*  double t_i = 0.004574998; */ /* for z=200, H_0=50, EdS */
  /* double t_0 = 13.04;       */ /* for z=200, H_0=50, EdS */
  double t_i[N_COSM];
  double t_0[N_COSM];
  double t_background[N_COSM][N_PLOT_T];
  double rza_Q_D[N_COSM][N_PLOT_R][N_PLOT_T];
  double rza_R_D[N_COSM][N_PLOT_R][N_PLOT_T];
  double a_D[N_COSM][N_PLOT_R][N_PLOT_T];
  double dot_a_D[N_COSM][N_PLOT_R][N_PLOT_T];
  double unphysical[N_COSM][N_PLOT_R][N_PLOT_T];
  double H_D[N_COSM][N_PLOT_R][N_PLOT_T];
  double Omm_D[N_COSM][N_PLOT_R][N_PLOT_T];
  double OmQ_D[N_COSM][N_PLOT_R][N_PLOT_T];
  double OmLam_D[N_COSM][N_PLOT_R][N_PLOT_T];
  double OmR_D[N_COSM][N_PLOT_R][N_PLOT_T];
  double delta_t;
  double delta_ln_R;
  double a_FLRW_expected;
  double a_dot_FLRW_expected;
  double H_FLRW_expected;

  /* ordinary sigma^2 (e.g. \sigma_8^2) */
  double sigmaIsq[N_COSM][N_PLOT_R][N_PLOT_T];
  double sigmaIsq_err[N_COSM][N_PLOT_R][N_PLOT_T];
  /* double sigmaIsq_grown[N_COSM][N_PLOT_R][N_PLOT_T]; */
  double sigmaI_normalised[N_COSM][N_PLOT_R][N_PLOT_T];

  double sigma8sq[N_COSM][N_PLOT_T];
  double sigma8sq_err[N_COSM][N_PLOT_T];
  /* double sigma8sq_grown[N_COSM][N_PLOT_T]; */

  /* test mode: compare with RZA2 arXiv:1303.6193v2 */
  int want_planar =0; /* cf RZA2 V.A */
  int want_spherical = 0; /* cf RZA2 V.B.3 */
  /* TODO: choose pseudo-tophat window function from RZA2 */

  /* TODO: choose power spectrum from RZA2 */


  const gsl_rng_type * T_gsl;
  gsl_rng * r_gsl;
  static unsigned long int local_gsl_seed=0;

  /* long   n_calls_invariants = 40000; */
  long   n_calls_invariants = 100000;
  /* long   n_calls_invariants = 500000; */
  double n_sigma[3] =
  /*    {2.0, 2.0, 2.0};  */
  /*  {1.0, 1.0, 1.0}; */
      {-1.0, -1.0, -1.0};
  /* {0.0, -1.0, 0.0};  */
  /* {-1.0, -1.0, -1.0};   */
  /* {1.0, 0.0, 1.0};   */
  /*    {-1.0, 0.3333333333333, -0.037037037};  */

  /* values for several different cases */
#define N_RZA2 27
  double n_sigma_RZA2[N_RZA2][3]; /* =
    { {-1.0, 0.0, 0.0},
    {1.0, 0.0, 0.0} }; */
  int i_sigma,ii_sigma,iii_sigma; /* initialisation of n_sigma_RZA2 values */
  int i_sigma_min = 0;
  int i_sigma_max = 2;
  int ii_sigma_min = -1;
  int ii_sigma_max = 1;
  int iii_sigma_min = -1;
  int iii_sigma_max = 1;
  int i_RZA2; /* iterate over RZA2 figure requirements */
  int rza2_figures=1; /* Enable calculation for all RZA2 figures? plotting test?
                       * In this clone: use this for the BAO scale for + and - fluctuations
                       */
  int n_rza2; /* either N_RZA2 or an override */

  /* BAO scale: sqrt(90*123) \approx 105; /0.5 = h^{-1} */
  double R_domain_lower_BKS00 =
    0.5 * 1.0/0.67*INHOMOG_A_SCALE_FACTOR_INITIAL / INHOMOG_A_SCALE_FACTOR_NOW;
  double R_domain_upper_BKS00 =
    0.5 * 100.0/0.67*INHOMOG_A_SCALE_FACTOR_INITIAL / INHOMOG_A_SCALE_FACTOR_NOW;
#define EIGHT_MPC_H100 8.0

  /* n_calls_invariants = 300000; */
  /*
    { { {1.1178, 0.1758, 0.000103667, 0.663858},
        {1.02825, 0.241843, -6.35389e-05, 0.710868} },
      { {1.31723, 0.290207, -0.0105453, 0},
        {1.06082, 0.744882, -0.004969, 0} } };
  */

  /* n_calls_invariants = 10000;
    { { {1.11369, 0.175526, 3.59408e-05, 0.663589},
        {1.02291, 0.245114, -8.68882e-05, 0.713223} },
      { {1.30339, 0.288883, -0.0113041, 0  },
        {1.05623, 0.737999, -0.00639778, 0  } } };
  */



  /* int    I_invariant; */

  int i_t, i_R; /* counters in time and length scale */
  int i_EdS;

  gsl_rng_env_setup();
  T_gsl = gsl_rng_default;
  gsl_rng_default_seed += 2467 + local_gsl_seed;
  local_gsl_seed = gsl_rng_default_seed;
  r_gsl = gsl_rng_alloc (T_gsl); /* this gets reallocated each time the
                                    function gets called */


  /* cf BKS00 Table 1, Fig 1 : final FLRW comoving units  16 Mpc, 100 Mpc */
  /*
    R_domain_lower_BKS00 = 0.08;
    R_domain_upper_BKS00 = 0.5;
  */

  rza_integrand_params.w_type = 1; /* test standard spherical domains */

  rza_integrand_params.delta_R = DELTA_R_DEFAULT; /*   if(2==rza_integrand_params.w_type) */

  /* 'B' = BBKS (BKS version); 'e' = EisHu98 short formula set ;
     'E' = EisHu98 full formula set */

  rza_integrand_params.pow_spec_type = 'E';

  /* possible overrides: this may make   expected_a_D_Omegas  cause failure */
  /* RZA2  Fig 2 :
     physical radii approx 0.125 Mpc, 0.5 Mpc;
     final FLRW comoving diameters (2R) approx 50 Mpc, 200 Mpc */
  /* R_domain_lower_BKS00 = 0.5 * 20.0*INHOMOG_A_SCALE_FACTOR_INITIAL / INHOMOG_A_SCALE_FACTOR_NOW;  */
  /*   R_domain_lower_BKS00 = 0.5 * 50.0*INHOMOG_A_SCALE_FACTOR_INITIAL / INHOMOG_A_SCALE_FACTOR_NOW;  */
  /* R_domain_lower_BKS00 =  33.8 *INHOMOG_A_SCALE_FACTOR_INITIAL / INHOMOG_A_SCALE_FACTOR_NOW;  */
  /* TODO: use the constants in inhomog.h */
  /* R_domain_lower_BKS00 = 0.5 * 25.0/(1.0+z_initial);  */
  /* R_domain_upper_BKS00 = 0.5 * 100.0*INHOMOG_A_SCALE_FACTOR_INITIAL / INHOMOG_A_SCALE_FACTOR_NOW;   */
  /*   R_domain_upper_BKS00 = 0.5 * 200.0*INHOMOG_A_SCALE_FACTOR_INITIAL / INHOMOG_A_SCALE_FACTOR_NOW;  */

  if(rza2_figures){
    n_rza2 = N_RZA2;
    i_RZA2=0;
    /* assign initial invariants I, II, III in terms of sigma's */
    for(i_sigma=i_sigma_min; i_sigma <= i_sigma_max; i_sigma++){
      for(ii_sigma=ii_sigma_min; ii_sigma <= ii_sigma_max; ii_sigma++){
        for(iii_sigma=iii_sigma_min; iii_sigma <= iii_sigma_max; iii_sigma++){
          n_sigma_RZA2[i_RZA2][0] = (double)i_sigma;
          n_sigma_RZA2[i_RZA2][1] = (double)ii_sigma;
          n_sigma_RZA2[i_RZA2][2] = (double)iii_sigma;
          i_RZA2++;
        };
      };
    }; /* for(i_sigma=i_sigma_min; i_sigma < i_sigma_max; i_sigma++) */
    if(i_RZA2 != n_rza2){
      printf("Error in assigning initial I, II, III values. Correct the code, please.\n");
      exit(1);
    };
  }else{
    n_rza2 = 1;
  };

  for(i_RZA2=0; i_RZA2<n_rza2; i_RZA2++){

    if(rza2_figures){
      n_sigma[0]= n_sigma_RZA2[i_RZA2][0];
      n_sigma[1]= n_sigma_RZA2[i_RZA2][1];
      n_sigma[2]= n_sigma_RZA2[i_RZA2][2];
    };

    background_cosm_params.flatFLRW = 1;
    for(i_EdS=N_COSM_I; i_EdS<N_COSM; i_EdS++){
      background_cosm_params.EdS = i_EdS;
      printf("background_cosm_params.EdS  = %d\n",
             background_cosm_params.EdS);

      background_cosm_params.recalculate_t_0 = 1; /* initially recalculate for this test */
      background_cosm_params.Theta2p7 = 2.726/2.7; /* e.g. arXiv:0911.1955 */
      if(1 == background_cosm_params.EdS){
        background_cosm_params.H_0 = 50.0;
        background_cosm_params.Omm_0 = 1.0;
        background_cosm_params.OmLam_0 = 0.0;
        background_cosm_params.Ombary_0 = 0.088;
        /* Planck: Ade et al 2013 arXiv:1303.5076v3, Table 2 */
        background_cosm_params.sigma8 = 0.83;

        t_i[i_EdS] = t_EdS(&background_cosm_params,
                           (double)INHOMOG_A_SCALE_FACTOR_INITIAL,
                           want_verbose);
      }else if(1 == background_cosm_params.flatFLRW){
        /*
           WMAP:
           background_cosm_params.H_0 = 70.0;
           background_cosm_params.OmLam_0 = 0.73;
        */
        /*
           Planck:
        */
        background_cosm_params.H_0 = 67.0;
        background_cosm_params.OmLam_0 = 0.68;
        background_cosm_params.Ombary_0 = 0.044;
        /* Planck: Ade et al 2013 arXiv:1303.5076v3, Table 2 */
        background_cosm_params.sigma8 = 0.83;
        /*
          background_cosm_params.H_0 = 50.0;
          background_cosm_params.OmLam_0 = 0.001;
        */
        background_cosm_params.Omm_0 = 1.0 - background_cosm_params.OmLam_0;
        t_i[i_EdS] = t_flatFLRW(&background_cosm_params,
                                (double)INHOMOG_A_SCALE_FACTOR_INITIAL,
                                want_verbose);
      }else{
        printf("No other options for background_cosm_params so far in program.\n");
        exit(1);
      };

      background_cosm_params.recalculate_t_0 = 1;



      printf("\nOmega_D_effWEC:\n");

      t_0[i_EdS] = background_cosm_params.t_0;
      delta_t = (t_0[i_EdS]-t_i[i_EdS])/((double)(N_PLOT_T-1));
      delta_ln_R = (log(R_domain_upper_BKS00) - log(R_domain_lower_BKS00))/(double)(N_PLOT_R-1);

      /* random (50:50) choice of + or - fluctuations, independently for
         each of the three invariants  */

      /* not used in RZA2 article:
         for(I_invariant=0; I_invariant<3; I_invariant++){
         if(gsl_rng_uniform(r_gsl) > 0.5)
         n_sigma[I_invariant] = -n_sigma[I_invariant];
         };
      */

      for(i_t=0; i_t<N_PLOT_T; i_t++){
        t_background[i_EdS][i_t] = t_i[i_EdS] + (double)i_t * delta_t;
      };

      /* allow nested openmp if possible - so that the invariants can be calculated in parallel */
      /* #ifdef _OPENMP
         omp_set_nested(1);
         #endif
      */

#pragma omp parallel                            \
  default(shared)                               \
  private(i_R,                                                          \
          rza_integrand_params_sigma8)                                  \
  firstprivate(rza_integrand_params,                                    \
               background_cosm_params)
      {
#pragma omp for schedule(dynamic)
        for(i_R=0; i_R<N_PLOT_R; i_R++){
          switch(rza_integrand_params.w_type)
            {
            case 1:
            default:
              R_domain_mod[i_R] = exp( log(R_domain_lower_BKS00) + ((double)i_R)* delta_ln_R );
              rza_integrand_params.R_domain = R_domain_mod[i_R];
              break;

            case 2: /* shell-shaped domains */
              R_domain_mod[i_R] = exp( log(R_domain_lower_BKS00) + ((double)i_R)* delta_ln_R );
              rza_integrand_params.R_domain_1 = R_domain_mod[i_R] -
                0.5*rza_integrand_params.delta_R;
              rza_integrand_params.R_domain_2 = R_domain_mod[i_R] +
                0.5*rza_integrand_params.delta_R;
              break;
            };

          /* calculate the first invariant, i.e. same as for FLRW linear theory */
          rza_integrand_params.background_cosm_params =
            background_cosm_params;
          sigma_sq_invariant_I(rza_integrand_params,
                               n_calls_invariants,
                               want_verbose,
                               (double*)&(sigmaIsq[i_EdS][i_R][0]),
                               (double*)&(sigmaIsq_err[i_EdS][i_R][0]));

          /* calculate sigma_8^2 using the first invariant, i.e. same
             as for FLRW linear theory */
          if(0==i_R){
            rza_integrand_params_sigma8 = rza_integrand_params;
            rza_integrand_params_sigma8.R_domain = 0.5 * EIGHT_MPC_H100/
              (background_cosm_params.H_0/100.0) *
              INHOMOG_A_SCALE_FACTOR_INITIAL / INHOMOG_A_SCALE_FACTOR_NOW;
            sigma_sq_invariant_I(rza_integrand_params_sigma8,
                                 n_calls_invariants,
                                 want_verbose,
                                 (double*)&(sigma8sq[i_EdS][0]),
                                 (double*)&(sigma8sq_err[i_EdS][0]));
          };

          Omega_D(&rza_integrand_params,
                  background_cosm_params,
                  (double*)&(t_background[i_EdS][0]), (int)N_PLOT_T,
                  n_sigma,
                  n_calls_invariants,
                  want_planar, /* cf RZA2 V.A */
                  want_spherical,
                  want_verbose,
                  (double*)&(rza_Q_D[i_EdS][i_R][0]),
                  (double*)&(rza_R_D[i_EdS][i_R][0]),
                  (double*)&(a_D[i_EdS][i_R][0]),
                  (double*)&(dot_a_D[i_EdS][i_R][0]),
                  (int*)&(unphysical[i_EdS][i_R][0]),
                  (double*)&(H_D[i_EdS][i_R][0]),
                  (double*)&(Omm_D[i_EdS][i_R][0]),
                  (double*)&(OmQ_D[i_EdS][i_R][0]),
                  (double*)&(OmLam_D[i_EdS][i_R][0]),
                  (double*)&(OmR_D[i_EdS][i_R][0])
                  );
        };  /*    for(i_R=0; i_R<N_PLOT_R; i_R++) */
      }; /* #pragma omp parallel                            */


      /* normalise the sigma_I^2 estimates according to FLRW linear
         theory */
      for(i_R=0; i_R<N_PLOT_R; i_R++){
        for(i_t=0; i_t<N_PLOT_T; i_t++){
          sigmaI_normalised[i_EdS][i_R][i_t] =
            sqrt(sigmaIsq[i_EdS][i_R][0])
            / sqrt(sigma8sq[i_EdS][0])
            * background_cosm_params.sigma8;
            /* * growth_FLRW(&background_cosm_params,
                              t_background[i_EdS][i_t],
                              want_verbose)/
                growth_FLRW(&background_cosm_params,
                            t_background[i_EdS][0],
                            want_verbose)
            */
        };
      };


      for(i_t=0; i_t<N_PLOT_T; i_t++){
        for(i_R=0; i_R<N_PLOT_R; i_R++){
          a_FLRW_expected = exp(log(t_background[i_EdS][i_t]/t_i[i_EdS])*(2.0/3.0))
            * INHOMOG_A_SCALE_FACTOR_INITIAL;
          H_FLRW_expected = 2.0/(3.0*t_background[i_EdS][i_t]);
          a_dot_FLRW_expected = H_FLRW_expected * a_FLRW_expected;
          (void)a_dot_FLRW_expected; /* avoid compile warning */

          /* printf("test_scale_factor_D: Warning: curvature backreaction inferred\n");
             printf("from Hamilton constraint; not calculated directly.\n"); */
          /*
          OmR_D[i_EdS][i_R][i_t] =
            1.0  - Omm_D[i_EdS][i_R][i_t]
            - OmQ_D[i_EdS][i_R][i_t]
            - OmLam_D[i_EdS][i_R][i_t];
          */

          printf(/*nsig=*/" %4.1f %4.1f %4.1f: tB2 R : a_D H_D : %g %g : %g %g : %g %g %g : %g : %g %g\n",
                 n_sigma[0],
                 n_sigma[1],
                 n_sigma[2],
                 t_background[i_EdS][i_t],
                 R_domain_mod[i_R],
                 a_D[i_EdS][i_R][i_t],
                 H_D[i_EdS][i_R][i_t],
                 Omm_D[i_EdS][i_R][i_t],
                 OmQ_D[i_EdS][i_R][i_t],
                 /* OmLam_D[i_EdS][i_R][i_t], */
                 OmR_D[i_EdS][i_R][i_t],
                 sigmaI_normalised[i_EdS][i_R][i_t],
                 Omm_D[i_EdS][i_R][i_t] + OmQ_D[i_EdS][i_R][i_t] + OmR_D[i_EdS][i_R][i_t],
                 Omm_D[i_EdS][i_R][i_t] + 2.0* OmQ_D[i_EdS][i_R][i_t] + (2.0/3.0)*OmR_D[i_EdS][i_R][i_t]
                 );
        }; /* for(i_R=0; i_R<N_PLOT_R; i_R++) */
      }; /* for(i_t=0; i_t<N_PLOT_T; i_t++) */

    }; /*   for(i_EdS=0; i_EdS<N_COSM; i_EdS++) */

  }; /*   for(i_RZA2=0; i_RZA2<n_rza2; i_RZA2++) */

  gsl_rng_free(r_gsl);

  return 0;
}