12#include "gadgetconfig.h"
16#include <gsl/gsl_integration.h>
17#include <gsl/gsl_rng.h>
22#include "../data/allvars.h"
23#include "../data/dtypes.h"
24#include "../data/mymalloc.h"
25#include "../logs/timer.h"
26#include "../main/simulation.h"
27#include "../mpi_utils/mpi_utils.h"
28#include "../ngenic/ngenic.h"
29#include "../pm/pm_mpi_fft.h"
30#include "../system/system.h"
32double ngenic::ngenic_power_spec(
double k)
36#if defined(MULTICOMPONENTGLASSFILE) && defined(DIFFERENT_TRANSFER_FUNC)
39 switch(
All.PowerSpectrumType)
42 power = ngenic_powerspec_eh(k);
46 power = ngenic_powerspec_tabulated(k);
50 power = ngenic_powerspec_efstathiou(k);
54#if defined(MULTICOMPONENTGLASSFILE) && defined(DIFFERENT_TRANSFER_FUNC)
57 power = PowerSpec_DM_2ndSpecies(k);
61 power *=
pow(k,
All.PrimordialIndex - 1.0);
66void ngenic::free_power_table(
void) {
Mem.myfree(PowerTable); }
68void ngenic::read_power_table(
void)
74 sprintf(buf,
All.PowerSpectrumFile);
76 if(!(fd = fopen(buf,
"r")))
78 Terminate(
"can't read input spectrum in file '%s' on task %d\n", buf, ThisTask);
84 if(fscanf(fd,
" %lg %lg ", &k, &p) == 2)
93 mpi_printf(
"found %d rows in input spectrum table\n", NPowerTable);
95 PowerTable = (pow_table *)
Mem.mymalloc(
"PowerTable", NPowerTable *
sizeof(pow_table));
97 sprintf(buf,
All.PowerSpectrumFile);
99 if(!(fd = fopen(buf,
"r")))
101 Terminate(
"can't read input spectrum in file '%s' on task %d\n", buf, ThisTask);
109 if(fscanf(fd,
" %lg %lg ", &k, &p) == 2)
111 PowerTable[NPowerTable].logk = k;
112 PowerTable[NPowerTable].logD = p;
122 std::sort(PowerTable, PowerTable + NPowerTable);
125void ngenic::ngenic_initialize_powerspectrum(
void)
136 if(
All.PowerSpectrumType == 2)
139 if(
All.ReNormalizeInputSpectrum == 0 &&
All.PowerSpectrumType == 2)
147#ifdef DIFFERENT_TRANSFER_FUNC
151 res = ngenic_tophat_sigma2(R8);
153 if(ThisTask == 0 &&
All.PowerSpectrumType == 2)
154 printf(
"\nNormalization of spectrum in file: Sigma8 = %g\n",
sqrt(res));
156 Norm =
All.Sigma8 *
All.Sigma8 / res;
158 if(ThisTask == 0 &&
All.PowerSpectrumType == 2)
159 printf(
"Normalization adjusted to Sigma8=%g (Normfac=%g)\n\n",
All.Sigma8, Norm);
163 mpi_printf(
"NGENIC: Dplus=%g\n", Dplus);
166double ngenic::ngenic_powerspec_tabulated(
double k)
172 double logk =
log10(k);
174 if(logk < PowerTable[0].logk || logk > PowerTable[NPowerTable - 1].logk)
178 int binhigh = NPowerTable - 1;
180 while(binhigh - binlow > 1)
182 int binmid = (binhigh + binlow) / 2;
183 if(logk < PowerTable[binmid].logk)
189 double dlogk = PowerTable[binhigh].logk - PowerTable[binlow].logk;
194 double u = (logk - PowerTable[binlow].logk) / dlogk;
196 double logD = (1 - u) * PowerTable[binlow].logD + u * PowerTable[binhigh].logD;
198 double Delta2 =
pow(10.0, logD);
200 double P = Norm * Delta2 / (4 *
M_PI * kold * kold * kold);
205double ngenic::ngenic_powerspec_efstathiou(
double k)
207 return Norm * k /
pow(1 +
pow(AA * k +
pow(BB * k, 1.5) + CC * CC * k * k, nu), 2 / nu);
210double ngenic::ngenic_powerspec_eh(
double k) {
return Norm * k *
pow(ngenic_tk_eh(k), 2); }
212double ngenic::ngenic_tk_eh(
double k)
214 double q, theta, ommh2, a, s, gamma, L0, C0;
216 double omegam, ombh2;
231 a = 1. - 0.328 *
log(431. * ommh2) * ombh2 / ommh2 + 0.380 *
log(22.3 * ommh2) * (ombh2 / ommh2) * (ombh2 / ommh2);
232 gamma = a + (1. - a) / (1. +
exp(4 *
log(0.43 * k * s)));
234 q = k * theta * theta / gamma;
235 L0 =
log(2. *
exp(1.) + 1.8 * q);
236 C0 = 14.2 + 731. / (1. + 62.5 * q);
237 tmp = L0 / (L0 + C0 * q * q);
241double ngenic::ngenic_tophat_sigma2(
double R)
243 const int worksize = 1000000;
245 double result, abserr, kmin, kmax;
248 myparams par = {R,
this};
250 gsl_integration_workspace *workspace = gsl_integration_workspace_alloc(worksize);
251 F.function = &sigma2_int;
254 if(
All.PowerSpectrumType == 2)
265 if(
All.PowerSpectrumType == 2)
272 gsl_integration_qag(&F,
log(kmin),
log(kmax), 0, 0.1, worksize, GSL_INTEG_GAUSS15, workspace, &result, &abserr);
275 double errbound = 1.0e-8 / NPowerTable * result;
278 for(
int i = 0; i < NPowerTable - 2; i++)
284 gsl_integration_qag(&F,
log(k0),
log(k1), errbound, 0, worksize, GSL_INTEG_GAUSS15, workspace, &x, &abserr);
292 gsl_integration_qag(&F,
log(kmin),
log(kmax), 0, 1.0e-8, worksize, GSL_INTEG_GAUSS15, workspace, &result, &abserr);
295 gsl_integration_workspace_free(workspace);
300double ngenic::ngenic_growth_factor(
double astart,
double aend) {
return ngenic_growth(aend) / ngenic_growth(astart); }
302double ngenic::ngenic_growth(
double a)
308 const int worksize = 100000;
310 double result, abserr;
313 gsl_integration_workspace *workspace = gsl_integration_workspace_alloc(worksize);
314 F.function = &ngenic_growth_int;
316 gsl_integration_qag(&F, 0, a, 0, 1.0e-8, worksize, GSL_INTEG_GAUSS41, workspace, &result, &abserr);
318 gsl_integration_workspace_free(workspace);
320 return hubble_a * result;
323double ngenic::ngenic_f1_omega(
double a)
329 return pow(omega_a, 0.6);
332double ngenic::ngenic_f2_omega(
double a)
338 return 2 *
pow(omega_a, 4.0 / 7);
global_data_all_processes All
expr pow(half base, half exp)