GADGET-4
driftfac.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#include <gsl/gsl_integration.h>
15#include <gsl/gsl_math.h>
16#include <math.h>
17#include <mpi.h>
18#include <stdio.h>
19#include <stdlib.h>
20#include <string.h>
21
22#include "../data/allvars.h"
23#include "../data/dtypes.h"
24#include "../time_integration/driftfac.h"
25
27{
28#define WORKSIZE 100000
29
30 gsl_function F;
31 gsl_integration_workspace *workspace;
32
33 logTimeBegin = log(All.TimeBegin);
34 logTimeMax = log(All.TimeMax);
35
36 workspace = gsl_integration_workspace_alloc(WORKSIZE);
37
38 for(int i = 0; i < DRIFT_TABLE_LENGTH; i++)
39 {
40 double result, abserr;
41
42 F.function = &driftfac::drift_integ;
43 gsl_integration_qag(&F, exp(logTimeBegin), exp(logTimeBegin + ((logTimeMax - logTimeBegin) / DRIFT_TABLE_LENGTH) * (i + 1)), 0,
44 1.0e-8, WORKSIZE, GSL_INTEG_GAUSS41, workspace, &result, &abserr);
45 DriftTable[i] = result;
46
47 F.function = &gravkick_integ;
48 gsl_integration_qag(&F, exp(logTimeBegin), exp(logTimeBegin + ((logTimeMax - logTimeBegin) / DRIFT_TABLE_LENGTH) * (i + 1)), 0,
49 1.0e-8, WORKSIZE, GSL_INTEG_GAUSS41, workspace, &result, &abserr);
50 GravKickTable[i] = result;
51
52 F.function = &hydrokick_integ;
53 gsl_integration_qag(&F, exp(logTimeBegin), exp(logTimeBegin + ((logTimeMax - logTimeBegin) / DRIFT_TABLE_LENGTH) * (i + 1)), 0,
54 1.0e-8, WORKSIZE, GSL_INTEG_GAUSS41, workspace, &result, &abserr);
55 HydroKickTable[i] = result;
56 }
57
58 gsl_integration_workspace_free(workspace);
59}
60
69{
70 static integertime last_time0 = -1, last_time1 = -1;
71 static double last_value;
72
73 if(time0 == last_time0 && time1 == last_time1)
74 return last_value;
75
76 /* note: will only be called for cosmological integration */
77
78 double a1 = logTimeBegin + time0 * All.Timebase_interval;
79 double a2 = logTimeBegin + time1 * All.Timebase_interval;
80
81 double u1 = (a1 - logTimeBegin) / (logTimeMax - logTimeBegin) * DRIFT_TABLE_LENGTH;
82 int i1 = (int)u1;
83 if(i1 >= DRIFT_TABLE_LENGTH)
84 i1 = DRIFT_TABLE_LENGTH - 1;
85
86 double df1;
87 if(i1 <= 1)
88 df1 = u1 * DriftTable[0];
89 else
90 df1 = DriftTable[i1 - 1] + (DriftTable[i1] - DriftTable[i1 - 1]) * (u1 - i1);
91
92 double u2 = (a2 - logTimeBegin) / (logTimeMax - logTimeBegin) * DRIFT_TABLE_LENGTH;
93 int i2 = (int)u2;
94 if(i2 >= DRIFT_TABLE_LENGTH)
95 i2 = DRIFT_TABLE_LENGTH - 1;
96
97 double df2;
98 if(i2 <= 1)
99 df2 = u2 * DriftTable[0];
100 else
101 df2 = DriftTable[i2 - 1] + (DriftTable[i2] - DriftTable[i2 - 1]) * (u2 - i2);
102
103 last_time0 = time0;
104 last_time1 = time1;
105
106 return last_value = (df2 - df1);
107}
108
110{
111 static integertime last_time0 = -1, last_time1 = -1;
112 static double last_value;
113
114 if(time0 == last_time0 && time1 == last_time1)
115 return last_value;
116
117 /* note: will only be called for cosmological integration */
118
119 double a1 = logTimeBegin + time0 * All.Timebase_interval;
120 double a2 = logTimeBegin + time1 * All.Timebase_interval;
121
122 double u1 = (a1 - logTimeBegin) / (logTimeMax - logTimeBegin) * DRIFT_TABLE_LENGTH;
123 int i1 = (int)u1;
124 if(i1 >= DRIFT_TABLE_LENGTH)
125 i1 = DRIFT_TABLE_LENGTH - 1;
126
127 double df1;
128 if(i1 <= 1)
129 df1 = u1 * GravKickTable[0];
130 else
131 df1 = GravKickTable[i1 - 1] + (GravKickTable[i1] - GravKickTable[i1 - 1]) * (u1 - i1);
132
133 double u2 = (a2 - logTimeBegin) / (logTimeMax - logTimeBegin) * DRIFT_TABLE_LENGTH;
134 int i2 = (int)u2;
135 if(i2 >= DRIFT_TABLE_LENGTH)
136 i2 = DRIFT_TABLE_LENGTH - 1;
137
138 double df2;
139 if(i2 <= 1)
140 df2 = u2 * GravKickTable[0];
141 else
142 df2 = GravKickTable[i2 - 1] + (GravKickTable[i2] - GravKickTable[i2 - 1]) * (u2 - i2);
143
144 last_time0 = time0;
145 last_time1 = time1;
146
147 return last_value = (df2 - df1);
148}
149
151{
152 static integertime last_time0 = -1, last_time1 = -1;
153 static double last_value;
154
155 if(time0 == last_time0 && time1 == last_time1)
156 return last_value;
157
158 /* note: will only be called for cosmological integration */
159
160 double a1 = logTimeBegin + time0 * All.Timebase_interval;
161 double a2 = logTimeBegin + time1 * All.Timebase_interval;
162
163 double u1 = (a1 - logTimeBegin) / (logTimeMax - logTimeBegin) * DRIFT_TABLE_LENGTH;
164 int i1 = (int)u1;
165 if(i1 >= DRIFT_TABLE_LENGTH)
166 i1 = DRIFT_TABLE_LENGTH - 1;
167
168 double df1;
169 if(i1 <= 1)
170 df1 = u1 * HydroKickTable[0];
171 else
172 df1 = HydroKickTable[i1 - 1] + (HydroKickTable[i1] - HydroKickTable[i1 - 1]) * (u1 - i1);
173
174 double u2 = (a2 - logTimeBegin) / (logTimeMax - logTimeBegin) * DRIFT_TABLE_LENGTH;
175 int i2 = (int)u2;
176 if(i2 >= DRIFT_TABLE_LENGTH)
177 i2 = DRIFT_TABLE_LENGTH - 1;
178
179 double df2;
180 if(i2 <= 1)
181 df2 = u2 * HydroKickTable[0];
182 else
183 df2 = HydroKickTable[i2 - 1] + (HydroKickTable[i2] - HydroKickTable[i2 - 1]) * (u2 - i2);
184
185 last_time0 = time0;
186 last_time1 = time1;
187
188 return last_value = (df2 - df1);
189}
190
192{
193 /* we just need to multiply this with the speed of light to get the correct cosmological factor */
194 double fac = get_gravkick_factor(time0, TIMEBASE);
195
196 return fac * (CLIGHT / All.UnitVelocity_in_cm_per_s);
197}
198
200{
201 integertime time0 = log(ascale / All.TimeBegin) / All.Timebase_interval;
202
203 /* we just need to multiply this with the speed of light to get the correct cosmological factor */
204 double fac = get_gravkick_factor(time0, TIMEBASE);
205
206 return fac * (CLIGHT / All.UnitVelocity_in_cm_per_s);
207}
208
210{
211 double fac = dist / (CLIGHT / All.UnitVelocity_in_cm_per_s);
212
214
215 double ascale = All.TimeBegin * exp(time0 * All.Timebase_interval);
216
217 return ascale;
218}
219
221{
222 integertime time1 = TIMEBASE;
223
224 double a2 = logTimeBegin + time1 * All.Timebase_interval;
225
226 double u2 = (a2 - logTimeBegin) / (logTimeMax - logTimeBegin) * DRIFT_TABLE_LENGTH;
227 int i2 = (int)u2;
228 if(i2 >= DRIFT_TABLE_LENGTH)
229 i2 = DRIFT_TABLE_LENGTH - 1;
230
231 double df2;
232 if(i2 <= 1)
233 df2 = u2 * GravKickTable[0];
234 else
235 df2 = GravKickTable[i2 - 1] + (GravKickTable[i2] - GravKickTable[i2 - 1]) * (u2 - i2);
236
237 double df1 = df2 - fac;
238
239 int i0 = 0;
240 int i1 = DRIFT_TABLE_LENGTH - 1;
241
242 if(df1 < 0 || df1 > GravKickTable[i1])
243 Terminate("out of range: df1=%g GravKickTable[i0]=%g GravKickTable[i1]=%g\n", df1, GravKickTable[i0], GravKickTable[i1]);
244
245 double u1;
246
247 if(df1 <= GravKickTable[0])
248 u1 = df1 / GravKickTable[0];
249 else
250 {
251 while(i1 - i0 > 1)
252 {
253 int im = (i0 + i1) / 2;
254 if(df1 < GravKickTable[im])
255 i1 = im;
256 else
257 i0 = im;
258 }
259
260 u1 = (df1 - GravKickTable[i0]) / (GravKickTable[i1] - GravKickTable[i0]) + i1;
261 }
262
263 double a1 = u1 * (logTimeMax - logTimeBegin) / DRIFT_TABLE_LENGTH + logTimeBegin;
264
265 integertime time0 = (a1 - logTimeBegin) / All.Timebase_interval;
266
267 return time0;
268}
global_data_all_processes All
Definition: main.cc:40
double get_drift_factor(integertime time0, integertime time1)
Definition: driftfac.cc:68
double get_comoving_distance_for_scalefactor(double ascale)
Definition: driftfac.cc:199
double get_gravkick_factor(integertime time0, integertime time1)
Definition: driftfac.cc:109
double get_hydrokick_factor(integertime time0, integertime time1)
Definition: driftfac.cc:150
double get_scalefactor_for_comoving_distance(double dist)
Definition: driftfac.cc:209
double get_comoving_distance(integertime time0)
Definition: driftfac.cc:191
integertime get_gravkick_factor_inverse(double fac)
Definition: driftfac.cc:220
void init_drift_table(void)
Definition: driftfac.cc:26
#define CLIGHT
Definition: constants.h:121
int integertime
Definition: constants.h:331
#define TIMEBASE
Definition: constants.h:333
#define WORKSIZE
#define DRIFT_TABLE_LENGTH
Definition: driftfac.h:49
#define Terminate(...)
Definition: macros.h:15
expr exp(half arg)
Definition: half.hpp:2724
expr log(half arg)
Definition: half.hpp:2745
double UnitVelocity_in_cm_per_s
Definition: allvars.h:119