12#include "gadgetconfig.h"
14#include <gsl/gsl_integration.h>
15#include <gsl/gsl_math.h>
22#include "../data/allvars.h"
23#include "../data/dtypes.h"
24#include "../time_integration/driftfac.h"
28#define WORKSIZE 100000
31 gsl_integration_workspace *workspace;
36 workspace = gsl_integration_workspace_alloc(
WORKSIZE);
40 double result, abserr;
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;
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;
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;
58 gsl_integration_workspace_free(workspace);
70 static integertime last_time0 = -1, last_time1 = -1;
71 static double last_value;
73 if(time0 == last_time0 && time1 == last_time1)
88 df1 = u1 * DriftTable[0];
90 df1 = DriftTable[i1 - 1] + (DriftTable[i1] - DriftTable[i1 - 1]) * (u1 - i1);
99 df2 = u2 * DriftTable[0];
101 df2 = DriftTable[i2 - 1] + (DriftTable[i2] - DriftTable[i2 - 1]) * (u2 - i2);
106 return last_value = (df2 - df1);
111 static integertime last_time0 = -1, last_time1 = -1;
112 static double last_value;
114 if(time0 == last_time0 && time1 == last_time1)
129 df1 = u1 * GravKickTable[0];
131 df1 = GravKickTable[i1 - 1] + (GravKickTable[i1] - GravKickTable[i1 - 1]) * (u1 - i1);
140 df2 = u2 * GravKickTable[0];
142 df2 = GravKickTable[i2 - 1] + (GravKickTable[i2] - GravKickTable[i2 - 1]) * (u2 - i2);
147 return last_value = (df2 - df1);
152 static integertime last_time0 = -1, last_time1 = -1;
153 static double last_value;
155 if(time0 == last_time0 && time1 == last_time1)
170 df1 = u1 * HydroKickTable[0];
172 df1 = HydroKickTable[i1 - 1] + (HydroKickTable[i1] - HydroKickTable[i1 - 1]) * (u1 - i1);
181 df2 = u2 * HydroKickTable[0];
183 df2 = HydroKickTable[i2 - 1] + (HydroKickTable[i2] - HydroKickTable[i2 - 1]) * (u2 - i2);
188 return last_value = (df2 - df1);
233 df2 = u2 * GravKickTable[0];
235 df2 = GravKickTable[i2 - 1] + (GravKickTable[i2] - GravKickTable[i2 - 1]) * (u2 - i2);
237 double df1 = df2 - fac;
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]);
247 if(df1 <= GravKickTable[0])
248 u1 = df1 / GravKickTable[0];
253 int im = (i0 + i1) / 2;
254 if(df1 < GravKickTable[im])
260 u1 = (df1 - GravKickTable[i0]) / (GravKickTable[i1] - GravKickTable[i0]) + i1;
global_data_all_processes All
double get_drift_factor(integertime time0, integertime time1)
double get_comoving_distance_for_scalefactor(double ascale)
double get_gravkick_factor(integertime time0, integertime time1)
double get_hydrokick_factor(integertime time0, integertime time1)
double get_scalefactor_for_comoving_distance(double dist)
double get_comoving_distance(integertime time0)
integertime get_gravkick_factor_inverse(double fac)
void init_drift_table(void)
#define DRIFT_TABLE_LENGTH
double UnitVelocity_in_cm_per_s