aboutsummaryrefslogtreecommitdiff
path: root/src/test.c
blob: 2ea73d1a56dbc1d079c2d0b5debcb51ffdb318bb (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
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
#include "solid_angle.h"
#include <math.h>
#include <stdio.h>
#include "optics.h"
#include "muon.h"
#include "mt19937ar.h"
#include <gsl/gsl_sf_gamma.h>
#include "misc.h"
#include "sno_charge.h"
#include <gsl/gsl_integration.h>
#include "path.h"
#include "random.h"
#include "pdg.h"
#include <gsl/gsl_spline.h>

typedef int testFunction(char *err);

/* Table of some of the tabulated values of the refractive index of water as a
 * function of wavelength and temperature. In all cases, I just used the values
 * for standard atmospheric pressure and assume this corresponds approximately
 * to a density of 1000 kg/m^3.
 *
 * See Table 7 in https://aip.scitation.org/doi/pdf/10.1063/1.555859. */

struct refractive_index_results {
    double p;
    double T;
    double wavelength;
    double n;
} refractive_index_results[] = {
    {1.0, 0 ,  226.50, 1.39468},
    {1.0, 10,  226.50, 1.39439},
    {1.0, 20,  226.50, 1.39353},
    {1.0, 30,  226.50, 1.39224},
    {1.0, 0 ,  404.41, 1.34431},
    {1.0, 10,  404.41, 1.34404},
    {1.0, 20,  404.41, 1.34329},
    {1.0, 30,  404.41, 1.34218},
    {1.0, 0 ,  589.00, 1.33447},
    {1.0, 10,  589.00, 1.33422},
    {1.0, 20,  589.00, 1.33350},
    {1.0, 30,  589.00, 1.33243},
    {1.0, 0 ,  632.80, 1.33321},
    {1.0, 10,  632.80, 1.33296},
    {1.0, 20,  632.80, 1.33224},
    {1.0, 30,  632.80, 1.33118},
    {1.0, 0 , 1013.98, 1.32626},
    {1.0, 10, 1013.98, 1.32604},
    {1.0, 20, 1013.98, 1.32537},
    {1.0, 30, 1013.98, 1.32437},
    {1.0, 0 , 2325.42, 1.27663},
    {1.0, 10, 2325.42, 1.27663},
    {1.0, 20, 2325.42, 1.27627},
    {1.0, 30, 2325.42, 1.27563},
};

/* Table of the values of solid angle for various values of r0/r and L/r.
 *
 * See Table 1 in http://www.umich.edu/~ners312/CourseLibrary/SolidAngleOfADiskOffAxis.pdf. */

struct solid_angle_results {
    double L;
    double r0;
    double omega;
} solid_angle_results[] = {
    {0.5,0.0,3.4732594},
    {0.5,0.2,3.4184435},
    {0.5,0.4,3.2435434},
    {0.5,0.6,2.9185178},
    {0.5,0.8,2.4122535},
    {0.5,1.0,1.7687239},
    {0.5,1.2,1.1661307},
    {0.5,1.4,0.7428889},
    {0.5,1.6,0.4841273},
    {0.5,1.8,0.3287007},
    {0.5,2.0,0.2324189},
    {1.0,0.0,1.8403024},
    {1.0,0.2,1.8070933},
    {1.0,0.4,1.7089486},
    {1.0,0.6,1.5517370},
    {1.0,0.8,1.3488367},
    {1.0,1.0,1.1226876},
    {1.0,1.2,0.9003572},
    {1.0,1.4,0.7039130},
    {1.0,1.6,0.5436956},
    {1.0,1.8,0.4195415},
    {1.0,2.0,0.3257993},
    {1.5,0.0,1.0552591},
    {1.5,0.2,1.0405177},
    {1.5,0.4,0.9975504},
    {1.5,0.6,0.9301028},
    {1.5,0.8,0.8441578},
    {1.5,1.0,0.7472299},
    {1.5,1.2,0.6472056},
    {1.5,1.4,0.5509617},
    {1.5,1.6,0.4632819},
    {1.5,1.8,0.3866757},
    {1.5,2.0,0.3217142},
    {2.0,0.0,0.6633335},
    {2.0,0.2,0.6566352},
    {2.0,0.4,0.6370508},
    {2.0,0.6,0.6060694},
    {2.0,0.8,0.5659755},
    {2.0,1.0,0.5195359},
    {2.0,1.2,0.4696858},
    {2.0,1.4,0.4191714},
    {2.0,1.6,0.3702014},
    {2.0,1.8,0.3243908},
    {2.0,2.0,0.282707}
};

int test_muon_get_T(char *err)
{
    /* A very simple test to make sure the energy as a function of distance
     * along the track makes sense. Should include more detailed tests later. */
    double T, E, range;

    /* Assume initial kinetic energy is 1 GeV. */
    T = 1000.0;
    E = get_T(T,1e-9,1.0);

    /* At the beginning of the track we should have roughly the same energy. */
    if (!isclose(E, T, 1e-5, 0)) {
        sprintf(err, "KE = %.5f, but expected %.5f", E, T);
        return 1;
    }

    /* Assume initial kinetic energy is 1 GeV. */
    T = 1000.0;
    range = get_range(T,1.0);
    E = get_T(T,range,1.0);

    /* At the end of the track we should have roughly the same energy. */
    if (!isclose(E, 0, 1e-5, 1e-5)) {
        sprintf(err, "KE = %.5f, but expected %.5f", E, 0.0);
        return 1;
    }

    return 0;
}

int test_muon_get_range(char *err)
{
    /* A very simple test to make sure we read in the PDG table correctly. */
    double value;

    value = get_range(1.0,1.0);

    if (!isclose(value, 2.071e-3, 1e-5, 0)) {
        sprintf(err, "range = %.5f, but expected %.5f", value, 1.863e-3);
        return 1;
    }

    return 0;
}

int test_muon_get_dEdx(char *err)
{
    /* A very simple test to make sure we read in the PDG table correctly. */
    double value;

    value = get_dEdx(1.0,1.0);

    if (!isclose(value, 5.471, 1e-5, 0)) {
        sprintf(err, "dE/dx = %.5f, but expected %.5f", value, 6.097);
        return 1;
    }

    return 0;
}

int test_refractive_index(char *err)
{
    /* Tests the get_index() function. */
    int i;
    double n;
    struct refractive_index_results result;

    for (i = 0; i < sizeof(refractive_index_results)/sizeof(struct refractive_index_results); i++) {
        result = refractive_index_results[i];
        n = get_index(result.p, result.wavelength, result.T);

        if (!isclose(n, result.n, 1e-2, 0)) {
            sprintf(err, "n = %.5f, but expected %.5f", n, result.n);
            return 1;
        }
    }

    return 0;
}

int test_solid_angle_approx(char *err)
{
    /* Tests the get_solid_angle_approx() function. */
    int i;
    double pmt[3] = {0,0,0};
    double pos[3] = {0,0,1};
    double n[3] = {0,0,1};
    double r = 1.0;
    double solid_angle;

    solid_angle = get_solid_angle_approx(pos,pmt,n,r);

    if (!isclose(solid_angle, 2*M_PI*(1-1/sqrt(2)), 1e-2, 0)) {
        sprintf(err, "solid angle = %.5f, but expected %.5f", solid_angle, 2*M_PI*(1-1/sqrt(2)));
        return 1;
    }

    for (i = 0; i < sizeof(solid_angle_results)/sizeof(struct solid_angle_results); i++) {
        pos[0] = solid_angle_results[i].r0*r;
        pos[2] = solid_angle_results[i].L*r;

        solid_angle = get_solid_angle_approx(pos,pmt,n,r);

        if (!isclose(solid_angle, solid_angle_results[i].omega, 1e-2, 0)) {
            sprintf(err, "solid angle = %.5f, but expected %.5f", solid_angle, solid_angle_results[i].omega);
            return 1;
        }
    }

    return 0;
}

int test_solid_angle(char *err)
{
    /* Tests the get_solid_angle() function. */
    int i;
    double pmt[3] = {0,0,0};
    double pos[3] = {0,0,1};
    double n[3] = {0,0,1};
    double r = 1.0;
    double solid_angle;

    solid_angle = get_solid_angle(pos,pmt,n,r);

    if (!isclose(solid_angle, 2*M_PI*(1-1/sqrt(2)), 1e-9, 0)) {
        sprintf(err, "solid angle = %.5f, but expected %.5f", solid_angle, 2*M_PI*(1-1/sqrt(2)));
        return 1;
    }

    for (i = 0; i < sizeof(solid_angle_results)/sizeof(struct solid_angle_results); i++) {
        pos[0] = solid_angle_results[i].r0*r;
        pos[2] = solid_angle_results[i].L*r;

        solid_angle = get_solid_angle(pos,pmt,n,r);

        if (!isclose(solid_angle, solid_angle_results[i].omega, 1e-4, 0)) {
            sprintf(err, "solid angle = %.5f, but expected %.5f", solid_angle, solid_angle_results[i].omega);
            return 1;
        }
    }

    return 0;
}

static double sno_charge(double q, void *params)
{
    int n = ((int *) params)[0];
    return pq(q,n);
}

int test_sno_charge_integral(char *err)
{
    /* Test that the single PE charge distribution integrates to one. */
    double result, error;
    gsl_function F;
    size_t nevals;
    gsl_integration_cquad_workspace *w;
    int params[1];

    w = gsl_integration_cquad_workspace_alloc(100);

    F.function = &sno_charge;
    params[0] = 1;
    F.params = params;

    init_charge();

    gsl_integration_cquad(&F, -10.0, 1000.0, 0, 1e-9, w, &result, &error, &nevals);

    if (!isclose(result, 1.0, 1e-9, 1e-9)) {
        sprintf(err, "integral of single PE charge distribution is %.2f", result);
        goto err;
    }

    gsl_integration_cquad_workspace_free(w);

    return 0;

err:
    gsl_integration_cquad_workspace_free(w);
    return 1;
}

int test_logsumexp(char *err)
{
    /* Tests the logsumexp() function. */
    int i, j;
    double logp[100], mu, result, expected;

    init_genrand(0);

    for (i = 0; i < 100; i++) {
        mu = genrand_real2();

        for (j = 0; j < sizeof(logp)/sizeof(double); j++) {
            logp[j] = -mu + j*log(mu) - gsl_sf_lnfact(j);
        }
        result = logsumexp(logp, sizeof(logp)/sizeof(double));

        expected = 0.0;

        if (!isclose(result, expected, 1e-9, 1e-9)) {
            sprintf(err, "logsumexp(logp) = %.5g, but expected %.5g", result, expected);
            return 1;
        }
    }

    return 0;
}

double getKineticEnergy(double x, double T0)
{
    return 1.0;
}

int test_path(char *err)
{
    /* Tests the logsumexp() function. */
    int i, j, k;
    double pos0[3], dir0[3], T0, range, m;
    double T, t;
    double pos_expected[3], t_expected;
    double pos[3], dir[3];
    double E, mom, beta;
    double s;
    double theta0;

    T0 = 1.0;
    range = 1.0;
    m = 1.0;

    init_genrand(0);

    for (i = 0; i < 100; i++) {
        pos0[0] = genrand_real2();
        pos0[1] = genrand_real2();
        pos0[2] = genrand_real2();

        rand_sphere(dir0);

        path *p = path_init(pos0,dir0,T0,range,0.1,getKineticEnergy,NULL,NULL,0,m);

        for (j = 0; j < 100; j++) {
            s = range*j/99;
            pos_expected[0] = pos0[0] + dir0[0]*s;
            pos_expected[1] = pos0[1] + dir0[1]*s;
            pos_expected[2] = pos0[2] + dir0[2]*s;

            /* Calculate total energy */
            E = T0 + m;
            mom = sqrt(E*E - m*m);
            beta = mom/E;

            t_expected = s/(beta*SPEED_OF_LIGHT);

            path_eval(p,s,pos,dir,&T,&t,&theta0);

            for (k = 0; k < 3; k++) {
                if (!isclose(pos[k], pos_expected[k], 1e-9, 1e-9)) {
                    sprintf(err, "path_eval(%.2g) returned pos[%i] = %.5g, but expected %.5g", s, k, pos[k], pos_expected[k]);
                    return 1;
                }
            }

            for (k = 0; k < 3; k++) {
                if (!isclose(dir[k], dir0[k], 1e-9, 1e-9)) {
                    sprintf(err, "path_eval(%.2g) returned dir[%i] = %.5g, but expected %.5g", s, k, dir[k], dir0[k]);
                    return 1;
                }
            }

            if (!isclose(T, 1.0, 1e-9, 1e-9)) {
                sprintf(err, "path_eval(%.2g) returned T = %.5g, but expected %.5g", s, T, 1.0);
                return 1;
            }

            if (!isclose(t, t_expected, 1e-9, 1e-9)) {
                sprintf(err, "path_eval(%.2g) returned t = %.5g, but expected %.5g", s, t, t_expected);
                return 1;
            }
        }

        path_free(p);
    }

    return 0;
}

int test_interp1d(char *err)
{
    /* Tests the interp1d() function. */
    size_t i;
    double xp[100], yp[100], range, x, y, expected;
    gsl_interp_accel *acc;
    gsl_spline *spline;

    range = 1.0;

    init_genrand(0);

    for (i = 0; i < sizeof(xp)/sizeof(xp[0]); i++) {
        xp[i] = range*i/(100-1);
        yp[i] = genrand_real2();
    }

    acc = gsl_interp_accel_alloc();
    spline = gsl_spline_alloc(gsl_interp_linear,100);
    gsl_spline_init(spline,xp,yp,100);

    for (i = 0; i < 1000; i++) {
        x = genrand_real2()*range;
        expected = gsl_spline_eval(spline,x,acc);
        y = interp1d(x,xp,yp,100);
        if (!isclose(y, expected, 1e-9, 1e-9)) {
            sprintf(err, "interp1d(%.2g) returned %.5g, but expected %.5g", x, y, expected);
            goto err;
        }
    }

    gsl_interp_accel_free(acc);
    gsl_spline_free(spline);

    return 0;

err:
    gsl_interp_accel_free(acc);
    gsl_spline_free(spline);

    return 1;
}

struct tests {
    testFunction *test;
    char *name;
} tests[] = {
    {test_solid_angle, "test_solid_angle"},
    {test_solid_angle_approx, "test_solid_angle_approx"},
    {test_refractive_index, "test_refractive_index"},
    {test_muon_get_dEdx, "test_muon_get_dEdx"},
    {test_muon_get_range, "test_muon_get_range"},
    {test_muon_get_T, "test_muon_get_T"},
    {test_logsumexp, "test_logsumexp"},
    {test_sno_charge_integral, "test_sno_charge_integral"},
    {test_path, "test_path"},
    {test_interp1d, "test_interp1d"}
};

int main(int argc, char **argv)
{
    int i;
    char err[256];
    int retval = 0;
    struct tests test;

    for (i = 0; i < sizeof(tests)/sizeof(struct tests); i++) {
        test = tests[i];

        if (!test.test(err)) {
            printf("[\033[92mok\033[0m] %s\n", test.name);
        } else {
            printf("[\033[91mfail\033[0m] %s: %s\n", test.name, err);
            retval = 1;
        }
    }

    return retval;
}