GADGET-4
power.cc
Go to the documentation of this file.
1/*******************************************************************************
2 * \copyright This file is part of the GADGET4 N-body/SPH code developed
3 * \copyright by Volker Springel. Copyright (C) 2014-2020 by Volker Springel
4 * \copyright (vspringel@mpa-garching.mpg.de) and all contributing authors.
5 *******************************************************************************/
6
12#include "gadgetconfig.h"
13
14#ifdef NGENIC
15
16#include <gsl/gsl_integration.h>
17#include <gsl/gsl_rng.h>
18#include <math.h>
19#include <mpi.h>
20#include <stdlib.h>
21
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"
31
32double ngenic::ngenic_power_spec(double k)
33{
34 double power = 0;
35
36#if defined(MULTICOMPONENTGLASSFILE) && defined(DIFFERENT_TRANSFER_FUNC)
37 if(Type == 1)
38#endif
39 switch(All.PowerSpectrumType)
40 {
41 case 1:
42 power = ngenic_powerspec_eh(k);
43 break;
44
45 case 2:
46 power = ngenic_powerspec_tabulated(k);
47 break;
48
49 default:
50 power = ngenic_powerspec_efstathiou(k);
51 break;
52 }
53
54#if defined(MULTICOMPONENTGLASSFILE) && defined(DIFFERENT_TRANSFER_FUNC)
55 if(Type == 2)
56 {
57 power = PowerSpec_DM_2ndSpecies(k);
58 }
59#endif
60
61 power *= pow(k, All.PrimordialIndex - 1.0);
62
63 return power;
64}
65
66void ngenic::free_power_table(void) { Mem.myfree(PowerTable); }
67
68void ngenic::read_power_table(void)
69{
70 FILE *fd;
71 char buf[MAXLEN_PATH];
72 double k, p;
73
74 sprintf(buf, All.PowerSpectrumFile);
75
76 if(!(fd = fopen(buf, "r")))
77 {
78 Terminate("can't read input spectrum in file '%s' on task %d\n", buf, ThisTask);
79 }
80
81 NPowerTable = 0;
82 do
83 {
84 if(fscanf(fd, " %lg %lg ", &k, &p) == 2)
85 NPowerTable++;
86 else
87 break;
88 }
89 while(1);
90
91 fclose(fd);
92
93 mpi_printf("found %d rows in input spectrum table\n", NPowerTable);
94
95 PowerTable = (pow_table *)Mem.mymalloc("PowerTable", NPowerTable * sizeof(pow_table));
96
97 sprintf(buf, All.PowerSpectrumFile);
98
99 if(!(fd = fopen(buf, "r")))
100 {
101 Terminate("can't read input spectrum in file '%s' on task %d\n", buf, ThisTask);
102 }
103
104 NPowerTable = 0;
105 do
106 {
107 double p;
108
109 if(fscanf(fd, " %lg %lg ", &k, &p) == 2)
110 {
111 PowerTable[NPowerTable].logk = k;
112 PowerTable[NPowerTable].logD = p;
113 NPowerTable++;
114 }
115 else
116 break;
117 }
118 while(1);
119
120 fclose(fd);
121
122 std::sort(PowerTable, PowerTable + NPowerTable);
123}
124
125void ngenic::ngenic_initialize_powerspectrum(void)
126{
127 double res;
128
129 AA = 6.4 / All.ShapeGamma * (3.085678e24 / All.UnitLength_in_cm);
130 BB = 3.0 / All.ShapeGamma * (3.085678e24 / All.UnitLength_in_cm);
131 CC = 1.7 / All.ShapeGamma * (3.085678e24 / All.UnitLength_in_cm);
132 nu = 1.13;
133
134 R8 = 8 * (3.085678e24 / All.UnitLength_in_cm); /* 8 Mpc/h */
135
136 if(All.PowerSpectrumType == 2)
137 read_power_table();
138
139 if(All.ReNormalizeInputSpectrum == 0 && All.PowerSpectrumType == 2)
140 {
141 Norm = 1.0;
142 /* tabulated file is already at the initial redshift */
143 Dplus = 1.0;
144 }
145 else
146 {
147#ifdef DIFFERENT_TRANSFER_FUNC
148 Type = 1;
149#endif
150 Norm = 1.0;
151 res = ngenic_tophat_sigma2(R8);
152
153 if(ThisTask == 0 && All.PowerSpectrumType == 2)
154 printf("\nNormalization of spectrum in file: Sigma8 = %g\n", sqrt(res));
155
156 Norm = All.Sigma8 * All.Sigma8 / res;
157
158 if(ThisTask == 0 && All.PowerSpectrumType == 2)
159 printf("Normalization adjusted to Sigma8=%g (Normfac=%g)\n\n", All.Sigma8, Norm);
160
161 Dplus = ngenic_growth_factor(All.cf_atime, 1.0);
162 }
163 mpi_printf("NGENIC: Dplus=%g\n", Dplus);
164}
165
166double ngenic::ngenic_powerspec_tabulated(double k)
167{
168 double kold = k;
169
170 k *= (All.InputSpectrum_UnitLength_in_cm / All.UnitLength_in_cm); // convert to h/Mpc
171
172 double logk = log10(k);
173
174 if(logk < PowerTable[0].logk || logk > PowerTable[NPowerTable - 1].logk)
175 return 0;
176
177 int binlow = 0;
178 int binhigh = NPowerTable - 1;
179
180 while(binhigh - binlow > 1)
181 {
182 int binmid = (binhigh + binlow) / 2;
183 if(logk < PowerTable[binmid].logk)
184 binhigh = binmid;
185 else
186 binlow = binmid;
187 }
188
189 double dlogk = PowerTable[binhigh].logk - PowerTable[binlow].logk;
190
191 if(dlogk == 0)
192 Terminate("dlogk == 0");
193
194 double u = (logk - PowerTable[binlow].logk) / dlogk;
195
196 double logD = (1 - u) * PowerTable[binlow].logD + u * PowerTable[binhigh].logD;
197
198 double Delta2 = pow(10.0, logD);
199
200 double P = Norm * Delta2 / (4 * M_PI * kold * kold * kold);
201
202 return P;
203}
204
205double ngenic::ngenic_powerspec_efstathiou(double k)
206{
207 return Norm * k / pow(1 + pow(AA * k + pow(BB * k, 1.5) + CC * CC * k * k, nu), 2 / nu);
208}
209
210double ngenic::ngenic_powerspec_eh(double k) /* Eisenstein & Hu */ { return Norm * k * pow(ngenic_tk_eh(k), 2); }
211
212double ngenic::ngenic_tk_eh(double k) /* from Martin White */
213{
214 double q, theta, ommh2, a, s, gamma, L0, C0;
215 double tmp;
216 double omegam, ombh2;
217
218 /* other input parameters */
219
220 omegam = All.Omega0;
222
223 if(All.OmegaBaryon == 0)
224 ombh2 = 0.04 * All.HubbleParam * All.HubbleParam;
225
226 k *= (3.085678e24 / All.UnitLength_in_cm); /* convert to h/Mpc */
227
228 theta = 2.728 / 2.7;
229 ommh2 = omegam * All.HubbleParam * All.HubbleParam;
230 s = 44.5 * log(9.83 / ommh2) / sqrt(1. + 10. * exp(0.75 * log(ombh2))) * All.HubbleParam;
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)));
233 gamma *= omegam * All.HubbleParam;
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);
238 return (tmp);
239}
240
241double ngenic::ngenic_tophat_sigma2(double R)
242{
243 const int worksize = 1000000;
244
245 double result, abserr, kmin, kmax;
246 gsl_function F;
247
248 myparams par = {R, this};
249
250 gsl_integration_workspace *workspace = gsl_integration_workspace_alloc(worksize);
251 F.function = &sigma2_int;
252 F.params = &par;
253
254 if(All.PowerSpectrumType == 2)
255 {
256 kmin = pow(10.0, PowerTable[0].logk) * (All.UnitLength_in_cm / All.InputSpectrum_UnitLength_in_cm);
257 kmax = pow(10.0, PowerTable[NPowerTable - 1].logk) * (All.UnitLength_in_cm / All.InputSpectrum_UnitLength_in_cm);
258 }
259 else
260 {
261 kmin = 1.0e-15 / R;
262 kmax = 1.0e3 / R;
263 }
264
265 if(All.PowerSpectrumType == 2)
266 {
267 /* because of the oscillatory behaviour of the integrand, the gsl_integration_qag() has trouble with its error estimates
268 * when the function is piece-wise interpolated. That's why we integrate the tabulated function segment by segment.
269 */
270
271 /* first get a rough result with up to 10% relative error */
272 gsl_integration_qag(&F, log(kmin), log(kmax), 0, 0.1, worksize, GSL_INTEG_GAUSS15, workspace, &result, &abserr);
273
274 /* now set a low absolute error bound for each segment */
275 double errbound = 1.0e-8 / NPowerTable * result;
276 result = 0.0;
277
278 for(int i = 0; i < NPowerTable - 2; i++)
279 {
280 double k0 = pow(10.0, PowerTable[i].logk) * (All.UnitLength_in_cm / All.InputSpectrum_UnitLength_in_cm);
281 double k1 = pow(10.0, PowerTable[i + 1].logk) * (All.UnitLength_in_cm / All.InputSpectrum_UnitLength_in_cm);
282 double x;
283
284 gsl_integration_qag(&F, log(k0), log(k1), errbound, 0, worksize, GSL_INTEG_GAUSS15, workspace, &x, &abserr);
285
286 result += x;
287 }
288 }
289 else
290 {
291 /* for the smooth analytic function, we integrate directly with a relative error estimate */
292 gsl_integration_qag(&F, log(kmin), log(kmax), 0, 1.0e-8, worksize, GSL_INTEG_GAUSS15, workspace, &result, &abserr);
293 }
294
295 gsl_integration_workspace_free(workspace);
296
297 return result;
298}
299
300double ngenic::ngenic_growth_factor(double astart, double aend) { return ngenic_growth(aend) / ngenic_growth(astart); }
301
302double ngenic::ngenic_growth(double a)
303{
304 double hubble_a;
305
306 hubble_a = sqrt(All.Omega0 / (a * a * a) + (1 - All.Omega0 - All.OmegaLambda) / (a * a) + All.OmegaLambda);
307
308 const int worksize = 100000;
309
310 double result, abserr;
311 gsl_function F;
312
313 gsl_integration_workspace *workspace = gsl_integration_workspace_alloc(worksize);
314 F.function = &ngenic_growth_int;
315
316 gsl_integration_qag(&F, 0, a, 0, 1.0e-8, worksize, GSL_INTEG_GAUSS41, workspace, &result, &abserr);
317
318 gsl_integration_workspace_free(workspace);
319
320 return hubble_a * result;
321}
322
323double ngenic::ngenic_f1_omega(double a)
324{
325 double omega_a;
326
327 omega_a = All.Omega0 / (All.Omega0 + a * (1 - All.Omega0 - All.OmegaLambda) + a * a * a * All.OmegaLambda);
328
329 return pow(omega_a, 0.6);
330}
331
332double ngenic::ngenic_f2_omega(double a)
333{
334 double omega_a;
335
336 omega_a = All.Omega0 / (All.Omega0 + a * (1 - All.Omega0 - All.OmegaLambda) + a * a * a * All.OmegaLambda);
337
338 return 2 * pow(omega_a, 4.0 / 7);
339}
340
341#endif
global_data_all_processes All
Definition: main.cc:40
#define MAXLEN_PATH
Definition: constants.h:300
#define M_PI
Definition: constants.h:56
#define Terminate(...)
Definition: macros.h:15
memory Mem
Definition: main.cc:44
expr exp(half arg)
Definition: half.hpp:2724
expr log(half arg)
Definition: half.hpp:2745
expr log10(half arg)
Definition: half.hpp:2752
expr pow(half base, half exp)
Definition: half.hpp:2803
expr sqrt(half arg)
Definition: half.hpp:2777