aboutsummaryrefslogtreecommitdiff
path: root/src/electron.c
blob: 2f6902ef55372609f8574206fb7f483f45f3ca62 (plain)
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
#include <stdio.h>
#include <errno.h>
#include <string.h>
#include <stdlib.h>
#include <gsl/gsl_errno.h>
#include <gsl/gsl_spline.h>
#include <math.h>
#include "optics.h"
#include "quantum_efficiency.h"
#include "solid_angle.h"
#include "pdg.h"
#include "vector.h"
#include "electron.h"
#include "sno.h"
#include "scattering.h"
#include "pmt_response.h"
#include "misc.h"
#include <gsl/gsl_integration.h>
#include <math.h> /* for fmax() */
#include <gsl/gsl_sf_gamma.h>
#include "util.h"

static int initialized = 0;

static double *x, *dEdx_rad, *dEdx, *csda_range;
static size_t size;

static gsl_interp_accel *acc_dEdx_rad;
static gsl_spline *spline_dEdx_rad;

static gsl_interp_accel *acc_dEdx;
static gsl_spline *spline_dEdx;

static gsl_interp_accel *acc_range;
static gsl_spline *spline_range;

/* Electron critical energy in H2O and D2O.
 *
 * These values come from
 * http://pdgprod.lbl.gov/~deg/AtomicNuclearProperties/HTML/deuterium_oxide_liquid.html for D2O and
 * http://pdg.lbl.gov/2018/AtomicNuclearProperties/HTML/water_liquid.html for H2O.
 */
static const double ELECTRON_CRITICAL_ENERGY_H2O = 78.33;
static const double ELECTRON_CRITICAL_ENERGY_D2O = 78.33;

void electron_get_position_distribution_parameters(double T0, double *a, double *b)
{
    /* Computes the gamma distribution parameters describing the longitudinal
     * profile of the photons generated from an electromagnetic shower.
     *
     * From the PDG "Passage of Particles" section 32.5:
     *
     * "The mean longitudinal profile of the energy deposition in an
     * electromagnetic cascade is reasonably well described by a gamma
     * distribution."
     *
     * Here we use a slightly different form of the gamma distribution:
     *
     *     f(x) = x**(a-1)*exp(-x/b)/(Gamma(a)*b**a)
     *
     * I determined the b parameter by simulating high energy electrons using
     * RAT-PAC and determined that it's roughly equal to the radiation length.
     * To calculate the a parameter we use the formula from the PDG, i.e.
     *
     *     tmax = (a-1)/b = ln(E/E_C) - 0.5
     *
     * Therefore, we calculate a as:
     *
     *     a = tmax*b+1.
     *
     * `T` should be in units of MeV.
     *
     * Example:
     *
     *     double a, b;
     *     electron_get_position_distribution_parameters(1000.0, &a, &b);
     *     double pdf = gamma_pdf(x, a, b);
     *
     * See http://pdg.lbl.gov/2014/reviews/rpp2014-rev-passage-particles-matter.pdf. */
    double tmax;

    *b = RADIATION_LENGTH;
    tmax = log(T0/ELECTRON_CRITICAL_ENERGY_D2O) - 0.5;
    *a = fmax(1.1,tmax*(*b)/RADIATION_LENGTH + 1);
}

double electron_get_angular_distribution_alpha(double T0)
{
    /* To describe the angular distribution of photons in an electromagnetic
     * shower I came up with the heuristic form:
     *
     *     f(cos_theta) ~ exp(-abs(cos_theta-mu)^alpha/beta)
     *
     * I simulated high energy electrons using RAT-PAC in heavy water to fit
     * for the alpha and beta parameters as a function of energy and determined
     * that a reasonably good fit to the data was:
     *
     *     alpha = c0 + c1/log(c2*T0 + c3)
     *
     * where T0 is the initial energy of the electron in MeV and c0, c1, c2,
     * and c3 are constants which I fit out. */
    return fmin(100.0,3.141318e-1 + 2.08198e-01/log(6.33331e-03*T0 + 1.19213e+00));
}

double electron_get_angular_distribution_beta(double T0)
{
    /* To describe the angular distribution of photons in an electromagnetic
     * shower I came up with the heuristic form:
     *
     * f(cos_theta) ~ exp(-abs(cos_theta-mu)^alpha/beta)
     *
     * I simulated high energy electrons using RAT-PAC in heavy water to fit
     * for the alpha and beta parameters as a function of energy and determined
     * that a reasonably good fit to the data was:
     *
     * beta = c0 + c1/log(c2*T0 + c3)
     *
     * where T0 is the initial energy of the electron in MeV and c0, c1, c2,
     * and c3 are constants which I fit out. */
    return 1.35372e-01 + 2.22344e-01/log(1.96249e-02*T0 + 1.23912e+00);
}

double electron_get_angular_pdf_no_norm(double cos_theta, double alpha, double beta, double mu)
{
    return exp(-pow(fabs(cos_theta-mu),alpha)/beta);
}

double electron_get_angular_pdf_norm(double alpha, double beta, double mu)
{
    return pow(beta,1/alpha)*gsl_sf_gamma(1/alpha)*(gsl_sf_gamma_inc_P(1/alpha,pow(1-mu,alpha)/beta)+gsl_sf_gamma_inc_P(1/alpha,pow(1+mu,alpha)/beta))/alpha;
}

double electron_get_angular_pdf_delta_ray(double cos_theta, double alpha, double beta, double mu)
{
    /* Returns the probability density that a photon in the wavelength range
     * 200 nm - 800 nm is emitted at an angle cos_theta.
     *
     * The angular distribution is modelled by the pdf:
     *
     *     f(cos_theta) ~ exp(-abs(cos_theta-mu)^alpha/beta)
     *
     * where alpha and beta are constants which depend on the initial energy of
     * the particle.
     *
     * There is no nice closed form solution for the integral of this function,
     * so we just compute it on the fly. To make things faster, we keep track
     * of the last alpha, beta, and mu parameters that were passed in so we
     * avoid recomputing the normalization factor. */
    size_t i;
    static double last_alpha, last_beta, last_mu, norm;
    static int first = 1;
    static double x[10000], y[10000];

    if (first || alpha != last_alpha || beta != last_beta || mu != last_mu) {
        norm = electron_get_angular_pdf_norm(alpha, beta, mu);

        last_alpha = alpha;
        last_beta = beta;
        last_mu = mu;

        for (i = 0; i < LEN(x); i++) {
            x[i] = -1.0 + 2.0*i/(LEN(x)-1);
            y[i] = electron_get_angular_pdf_no_norm(x[i], alpha, beta, mu)/norm;
        }

        first = 0;
    }

    return interp1d(cos_theta,x,y,LEN(x));
}

double electron_get_angular_pdf(double cos_theta, double alpha, double beta, double mu)
{
    /* Returns the probability density that a photon in the wavelength range
     * 200 nm - 800 nm is emitted at an angle cos_theta.
     *
     * The angular distribution is modelled by the pdf:
     *
     *     f(cos_theta) ~ exp(-abs(cos_theta-mu)^alpha/beta)
     *
     * where alpha and beta are constants which depend on the initial energy of
     * the particle.
     *
     * There is no nice closed form solution for the integral of this function,
     * so we just compute it on the fly. To make things faster, we keep track
     * of the last alpha, beta, and mu parameters that were passed in so we
     * avoid recomputing the normalization factor. */
    size_t i;
    static double last_alpha, last_beta, last_mu, norm;
    static int first = 1;
    static double x[10000], y[10000];

    if (first || alpha != last_alpha || beta != last_beta || mu != last_mu) {
        norm = electron_get_angular_pdf_norm(alpha, beta, mu);

        last_alpha = alpha;
        last_beta = beta;
        last_mu = mu;

        for (i = 0; i < LEN(x); i++) {
            x[i] = -1.0 + 2.0*i/(LEN(x)-1);
            y[i] = electron_get_angular_pdf_no_norm(x[i], alpha, beta, mu)/norm;
        }

        first = 0;
    }

    return interp1d(cos_theta,x,y,LEN(x));
}

static int init()
{
    int i, j;
    char line[256];
    char *str;
    double value;
    int n;

    FILE *f = open_file("e_water_liquid.txt", "r");

    if (!f) {
        fprintf(stderr, "failed to open e_water_liquid.txt: %s\n", strerror(errno));
        return -1;
    }

    i = 0;
    n = 0;
    /* For the first pass, we just count how many values there are. */
    while (fgets(line, sizeof(line), f)) {
        size_t len = strlen(line);
        if (len && (line[len-1] != '\n')) {
            fprintf(stderr, "got incomplete line on line %i: '%s'\n", i, line);
            goto err;
        }

        i += 1;

        /* Skip the first 8 lines since it's just a header. */
        if (i <= 8) continue;

        if (!len) continue;
        else if (line[0] == '#') continue;

        str = strtok(line," \n");

        while (str) {
            value = strtod(str, NULL);
            str = strtok(NULL," \n");
        }

        n += 1;
    }

    x = malloc(sizeof(double)*n);
    dEdx_rad = malloc(sizeof(double)*n);
    dEdx = malloc(sizeof(double)*n);
    csda_range = malloc(sizeof(double)*n);
    size = n;

    i = 0;
    n = 0;
    /* Now, we actually store the values. */
    rewind(f);
    while (fgets(line, sizeof(line), f)) {
        size_t len = strlen(line);
        if (len && (line[len-1] != '\n')) {
            fprintf(stderr, "got incomplete line on line %i: '%s'\n", i, line);
            goto err;
        }

        i += 1;

        /* Skip the first 8 lines since it's just a header. */
        if (i <= 8) continue;

        if (!len) continue;
        else if (line[0] == '#') continue;

        str = strtok(line," \n");

        j = 0;
        while (str) {
            value = strtod(str, NULL);
            switch (j) {
            case 0:
                x[n] = value;
                break;
            case 2:
                dEdx_rad[n] = value;
                break;
            case 3:
                dEdx[n] = value;
                break;
            case 4:
                csda_range[n] = value;
                break;
            }
            j += 1;
            str = strtok(NULL," \n");
        }

        n += 1;
    }

    fclose(f);

    acc_dEdx_rad = gsl_interp_accel_alloc();
    spline_dEdx_rad = gsl_spline_alloc(gsl_interp_linear, size);
    gsl_spline_init(spline_dEdx_rad, x, dEdx_rad, size);

    acc_dEdx = gsl_interp_accel_alloc();
    spline_dEdx = gsl_spline_alloc(gsl_interp_linear, size);
    gsl_spline_init(spline_dEdx, x, dEdx, size);

    acc_range = gsl_interp_accel_alloc();
    spline_range = gsl_spline_alloc(gsl_interp_linear, size);
    gsl_spline_init(spline_range, x, csda_range, size);

    initialized = 1;

    return 0;

err:
    fclose(f);

    return -1;
}

double electron_get_range(double T, double rho)
{
    /* Returns the approximate range a electron with kinetic energy `T` will travel
     * in water before losing all of its energy. This range is interpolated
     * based on data from the PDG which uses the continuous slowing down
     * approximation.
     *
     * `T` should be in MeV, and `rho` should be in g/cm^3.
     *
     * Return value is in cm.
     *
     * See http://pdg.lbl.gov/2018/AtomicNuclearProperties/adndt.pdf. */
    if (!initialized) {
        if (init()) {
            exit(1);
        }
    }

    return gsl_spline_eval(spline_range, T, acc_range)/rho;
}

double electron_get_dEdx_rad(double T, double rho)
{
    /* Returns the approximate radiative dE/dx for a electron in water with
     * kinetic energy `T`.
     *
     * `T` should be in MeV and `rho` in g/cm^3.
     *
     * Return value is in MeV/cm.
     *
     * See http://pdg.lbl.gov/2018/AtomicNuclearProperties/adndt.pdf. */
    if (!initialized) {
        if (init()) {
            exit(1);
        }
    }

    if (T < spline_dEdx_rad->x[0]) return spline_dEdx_rad->y[0];

    return gsl_spline_eval(spline_dEdx_rad, T, acc_dEdx_rad)*rho;
}

double electron_get_dEdx(double T, double rho)
{
    /* Returns the approximate dE/dx for a electron in water with kinetic energy
     * `T`.
     *
     * `T` should be in MeV and `rho` in g/cm^3.
     *
     * Return value is in MeV/cm.
     *
     * See http://pdg.lbl.gov/2018/AtomicNuclearProperties/adndt.pdf. */
    if (!initialized) {
        if (init()) {
            exit(1);
        }
    }

    if (T < spline_dEdx->x[0]) return spline_dEdx->y[0];

    return gsl_spline_eval(spline_dEdx, T, acc_dEdx)*rho;
}