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
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
|
/* Copyright (c) 2019, Anthony Latorre <tlatorre at uchicago>
*
* This program is free software: you can redistribute it and/or modify it
* under the terms of the GNU General Public License as published by the Free
* Software Foundation, either version 3 of the License, or (at your option)
* any later version.
* This program is distributed in the hope that it will be useful, but WITHOUT
* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for
* more details.
* You should have received a copy of the GNU General Public License along with
* this program. If not, see <https://www.gnu.org/licenses/>.
*/
#include <math.h>
#include "optics.h"
#include "misc.h"
#include "quantum_efficiency.h"
#include <stdio.h> /* for fprintf() */
#include <gsl/gsl_integration.h>
#include <gsl/gsl_errno.h> /* for gsl_strerror() */
#include <gsl/gsl_spline.h>
#include "misc.h"
#include "db.h"
#include "sno.h"
char optics_err[256];
/* Index of refraction in water and heavy water averaged over the Cerenkov
* spectrum and the PMT quantum efficiency. */
static double avg_index_h2o = 0.0, avg_index_d2o = 0.0;
/* Absorption coefficients for H2O, D2O, and acrylic as a function of
* wavelength from SNOMAN. */
static double absorption_coefficient_d2o_wavelengths[6];
static double absorption_coefficient_d2o[6];
static double rayleigh_scl_fac_d2o;
static double isothermal_comp_d2o;
static double absorption_coefficient_h2o_wavelengths[6];
static double absorption_coefficient_h2o[6];
static double rayleigh_scl_fac_acrylic;
static double isothermal_comp_acrylic;
static double absorption_coefficient_acrylic_wavelengths[6];
static double absorption_coefficient_acrylic[6];
static double rayleigh_scl_fac_h2o;
static double isothermal_comp_h2o;
typedef struct optics_constants {
double absorption_coefficient_d2o_wavelengths[6];
double absorption_coefficient_d2o[6];
double rayleigh_scl_fac_d2o;
double isothermal_comp_d2o;
double absorption_coefficient_h2o_wavelengths[6];
double absorption_coefficient_h2o[6];
double rayleigh_scl_fac_acrylic;
double isothermal_comp_acrylic;
double absorption_coefficient_acrylic_wavelengths[6];
double absorption_coefficient_acrylic[6];
double rayleigh_scl_fac_h2o;
double isothermal_comp_h2o;
} optics_constants;
/* Note: These numbers come from mcprod/media_qoca_d2o_20060110.cmd. */
optics_constants d2o_constants = {
.absorption_coefficient_d2o_wavelengths = {337.0, 365.0, 386.0, 420.0, 500.0, 620.0},
.absorption_coefficient_d2o = {6.057e-05, 2.680e-05, 3.333e-05, 3.006e-05, 2.773e-05, 2.736e-05},
.rayleigh_scl_fac_d2o = 1.000,
.isothermal_comp_d2o = 4.92e-10,
.absorption_coefficient_h2o_wavelengths = {337.0, 365.0, 386.0, 420.0, 500.0, 620.0},
.absorption_coefficient_h2o = {8.410e-05, 2.880e-06, 1.431e-04, 1.034e-04, 5.239e-04, 2.482e-03},
.rayleigh_scl_fac_acrylic = 0.950,
.isothermal_comp_acrylic = 3.55e-10,
.absorption_coefficient_acrylic_wavelengths = {337.0, 365.0, 386.0, 420.0, 500.0, 620.0},
.absorption_coefficient_acrylic = {5.610e-02, 2.279e-02, 1.204e-02, 7.587e-03, 7.036e-03, 7.068e-03},
.rayleigh_scl_fac_h2o = 0.870,
.isothermal_comp_h2o = 4.78e-10,
};
/* Note: These numbers come from mcprod/media_qoca_salt_20060420.cmd. */
optics_constants salt_constants = {
.absorption_coefficient_d2o_wavelengths = {337.0, 365.0, 386.0, 420.0, 500.0, 620.0},
.absorption_coefficient_d2o = {1.024e-04, 8.922e-05, 8.763e-05, 9.653e-05, 1.207e-05, 1.177e-05},
.rayleigh_scl_fac_d2o = 1.289,
.isothermal_comp_d2o = 4.92e-10,
.absorption_coefficient_h2o_wavelengths = {337.0, 365.0, 386.0, 420.0, 500.0, 620.0},
.absorption_coefficient_h2o = {1.068e-06, 1.022e-05, 4.123e-05, 2.193e-04, 4.344e-04, 2.498e-03},
.rayleigh_scl_fac_h2o = 0.870,
.isothermal_comp_h2o = 4.78e-10,
.absorption_coefficient_acrylic_wavelengths = {337.0, 365.0, 386.0, 420.0, 500.0, 620.0},
.absorption_coefficient_acrylic = {5.610e-02, 2.279e-02, 1.204e-02, 7.587e-03, 7.036e-03, 7.068e-03},
.rayleigh_scl_fac_acrylic = 0.950,
.isothermal_comp_acrylic = 3.55e-10,
};
static int initialized = 0;
/* Number of points to sample when computing the average absorption and
* scattering probabilities.
*
* Note that we use this array to sample for very small distances when
* computing the absorption in the acrylic so (xhi-xlo)/N should be less than
* 10 cm. */
#define N 10000
/* Lower and upper bounds for the distances used to calculate the average
* absorption and scattering probabilities. */
static const double xlo = 0.0;
static const double xhi = 2*PSUP_RADIUS;
/* Array of lengths used in the lookup tables for the average absorption and
* scattering probabilities. */
static double fabs_x[N];
/* Average fraction of light which is not absorbed by D2O as a function of
* length weighted by the Cerenkov spectrum and the quantum efficiency.
*
* We calculate the fraction not absorbed instead of absorbed */
static double fabs_d2o[N];
static double fabs_h2o[N];
static double fabs_acrylic[N];
static double fsct_d2o[N];
static double fsct_h2o[N];
static double fsct_acrylic[N];
/* From Table 4 in the paper. */
static double A0 = 0.243905091;
static double A1 = 9.53518094e-3;
static double A2 = -3.64358110e-3;
static double A3 = 2.65666426e-4;
static double A4 = 1.59189325e-3;
static double A5 = 2.45733798e-3;
static double A6 = 0.897478251;
static double A7 = -1.63066183e-2;
static double UV = 0.2292020;
static double IR = 5.432937;
static double RIND_H2O_C1 = 1.302;
static double RIND_H2O_C2 = 0.01562;
static double RIND_H2O_C3 = 0.32;
static double RIND_D2O_C1 = 1.302;
static double RIND_D2O_C2 = 0.01333;
static double RIND_D2O_C3 = 0.32;
static double RIND_ACRYLIC_C1 = 1.452;
static double RIND_ACRYLIC_C2 = 0.02;
static double RIND_ACRYLIC_C3 = 0.32;
/* Wavelength of 1 eV particle in nm.
*
* From http://pdg.lbl.gov/2018/reviews/rpp2018-rev-phys-constants.pdf. */
static double HC = 1.239841973976e3;
double get_avg_index_d2o(void)
{
if (!initialized) {
fprintf(stderr, "average index in d2o hasn't been initialized!\n");
exit(1);
}
return avg_index_d2o;
}
double get_avg_index_h2o(void)
{
if (!initialized) {
fprintf(stderr, "average index in h2o hasn't been initialized!\n");
exit(1);
}
return avg_index_h2o;
}
double get_index(double p, double wavelength, double T)
{
/* Returns the index of refraction of pure water for a given density,
* wavelength, and temperature. The density should be in units of g/cm^3,
* the wavelength in nm, and the temperature in Celsius.
*
* See "Refractive Index of Water and Steam as a function of Wavelength,
* Temperature, and Density" by Schiebener et al. 1990. */
/* normalize the temperature and pressure */
wavelength = wavelength/589.0;
T = (T+273.15)/273.15;
/* first we compute the right hand side of Equation 7 */
double c = A0 + A1*p + A2*T + A3*pow(wavelength,2)*T + A4/pow(wavelength,2) + A5/(pow(wavelength,2)-pow(UV,2)) + A6/(pow(wavelength,2)-pow(IR,2)) + A7*pow(p,2);
c *= p;
return sqrt((2*c+1)/(1-c));
}
double get_index_snoman_h2o(double wavelength)
{
/* Returns the index of refraction of light water that SNOMAN uses for a
* given wavelength. The wavelength should be in units of nm.
*
* The coefficients come from media.dat in SNOMAN version 5.0294 and the
* formula used to compute the index of refraction comes from the SNOMAN
* companion page for titles MEDA. */
double E;
/* Calculate the energy of the photon in eV. */
E = HC/wavelength;
return RIND_H2O_C1 + RIND_H2O_C2*exp(RIND_H2O_C3*E);
}
double get_index_snoman_d2o(double wavelength)
{
/* Returns the index of refraction of heavy water that SNOMAN uses for a
* given wavelength. The wavelength should be in units of nm.
*
* The coefficients come from media.dat in SNOMAN version 5.0294 and the
* formula used to compute the index of refraction comes from the SNOMAN
* companion page for titles MEDA. */
double E;
/* Calculate the energy of the photon in eV. */
E = HC/wavelength;
return RIND_D2O_C1 + RIND_D2O_C2*exp(RIND_D2O_C3*E);
}
double get_index_snoman_acrylic(double wavelength)
{
/* Returns the index of refraction of SNO acrylic for a given wavelength.
* The wavelength should be in units of nm.
*
* The coefficients come from media.dat in SNOMAN version 5.0294 and the
* formula used to compute the index of refraction comes from the SNOMAN
* companion page for titles MEDA. */
double E;
/* Calculate the energy of the photon in eV. */
E = HC/wavelength;
return RIND_ACRYLIC_C1 + RIND_ACRYLIC_C2*exp(RIND_ACRYLIC_C3*E);
}
/* Returns the probability per unit length (1/cm) for Rayleigh scattering in a
* material with an isothermal compressibility of `isothermal_comp` at a given
* wavelength in nm.
*
* The isothermal compressibility should be in units of N^-1 m^2.
*
* This function is mostly copied from snoman's rayint_prob.for. */
double rayint_prob(double wavelength, double n, double isothermal_comp)
{
double E;
double confac = 1.53E26;
/* Calculate the energy of the photon in MeV. */
E = 1e-6*HC/wavelength;
return isothermal_comp*confac*pow(E,4)*pow(n*n - 1.0,2)*pow(n*n + 2.0,2);
}
/* Returns the probability per unit length (1/cm) for Rayleigh scattering in
* D2O at a given wavelength in nm. */
double rayint_prob_d2o(double wavelength)
{
return rayleigh_scl_fac_d2o*rayint_prob(wavelength, get_index_snoman_d2o(wavelength), isothermal_comp_d2o);
}
/* Returns the probability per unit length (1/cm) for Rayleigh scattering in
* H2O at a given wavelength in nm. */
double rayint_prob_h2o(double wavelength)
{
return rayleigh_scl_fac_h2o*rayint_prob(wavelength, get_index_snoman_h2o(wavelength), isothermal_comp_h2o);
}
/* Returns the probability per unit length (1/cm) for Rayleigh scattering in
* acrylic at a given wavelength in nm. */
double rayint_prob_acrylic(double wavelength)
{
return rayleigh_scl_fac_acrylic*rayint_prob(wavelength, get_index_snoman_acrylic(wavelength), isothermal_comp_acrylic);
}
/* Returns the average probability that a photon is not absorbed after a
* distance `x` in cm in D2O.
*
* This average probability is found by averaging
*
* e^(-x/l(wavelength))
*
* over the wavelength range 200-800 nm weighted by the Cerenkov wavelenght
* distribution and the quantum efficiency. */
double get_fabs_d2o(double x)
{
if (!initialized) {
fprintf(stderr, "weighted absorption length hasn't been initialized!\n");
exit(1);
}
return interp1d(x,fabs_x,fabs_d2o,LEN(fabs_x));
}
/* Returns the average probability that a photon is not absorbed after a
* distance `x` in cm in H2O. See get_fabs_d2o() for how this is calculated. */
double get_fabs_h2o(double x)
{
if (!initialized) {
fprintf(stderr, "average probability hasn't been initialized!\n");
exit(1);
}
return interp1d(x,fabs_x,fabs_h2o,LEN(fabs_x));
}
/* Returns the average probability that a photon is not absorbed after a
* distance `x` in cm in the SNO acrylic. See get_fabs_d2o() for how this is
* calculated.
*
* FIXME: Should include a term for the "dark" acrylic. */
double get_fabs_acrylic(double x)
{
if (!initialized) {
fprintf(stderr, "weighted absorption length hasn't been initialized!\n");
exit(1);
}
return interp1d(x,fabs_x,fabs_acrylic,LEN(fabs_x));
}
/* Returns the average probability that a photon is not scattered after a
* distance `x` in cm in D2O. See get_fabs_d2o() for how this is calculated. */
double get_fsct_d2o(double x)
{
if (!initialized) {
fprintf(stderr, "weighted rayleigh scattering length hasn't been initialized!\n");
exit(1);
}
return interp1d(x,fabs_x,fsct_d2o,LEN(fabs_x));
}
/* Returns the average probability that a photon is not scattered after a
* distance `x` in cm in H2O. See get_fabs_d2o() for how this is calculated. */
double get_fsct_h2o(double x)
{
if (!initialized) {
fprintf(stderr, "weighted rayleigh scattering length hasn't been initialized!\n");
exit(1);
}
return interp1d(x,fabs_x,fsct_h2o,LEN(fabs_x));
}
/* Returns the average probability that a photon is not scattered after a
* distance `x` in cm in SNO acrylic. See get_fabs_d2o() for how this is
* calculated. */
double get_fsct_acrylic(double x)
{
if (!initialized) {
fprintf(stderr, "weighted rayleigh scattering length hasn't been initialized!\n");
exit(1);
}
return interp1d(x,fabs_x,fsct_acrylic,LEN(fabs_x));
}
/* Returns the absorption length as a function of wavelength in D2O.
*
* The wavelength should be in nm and the absorption length is in cm. */
double get_absorption_length_snoman_d2o(double wavelength)
{
/* Note: We use GSL splines here instead of interp1d() since we aren't
* guaranteed that the values in the lookup table are evenly spaced. */
static gsl_spline *spline;
static gsl_interp_accel *acc;
if (!spline) {
spline = gsl_spline_alloc(gsl_interp_linear, LEN(absorption_coefficient_d2o_wavelengths));
gsl_spline_init(spline, absorption_coefficient_d2o_wavelengths, absorption_coefficient_d2o, LEN(absorption_coefficient_d2o_wavelengths));
acc = gsl_interp_accel_alloc();
}
/* If the wavelength is outside the range of the lookup table we just
* return the value at the beginning or end. We do this because by default
* gsl_spline_eval will exit if the wavelength is out of bounds. */
if (wavelength < absorption_coefficient_d2o_wavelengths[0])
wavelength = absorption_coefficient_d2o_wavelengths[0];
if (wavelength > absorption_coefficient_d2o_wavelengths[LEN(absorption_coefficient_d2o_wavelengths)-1])
wavelength = absorption_coefficient_d2o_wavelengths[LEN(absorption_coefficient_d2o_wavelengths)-1];
return 1.0/gsl_spline_eval(spline, wavelength, acc);
}
/* Returns the absorption length as a function of wavelength in H2O.
*
* The wavelength should be in nm and the absorption length is in cm. */
double get_absorption_length_snoman_h2o(double wavelength)
{
/* Note: We use GSL splines here instead of interp1d() since we aren't
* guaranteed that the values in the lookup table are evenly spaced. */
static gsl_spline *spline;
static gsl_interp_accel *acc;
if (!spline) {
spline = gsl_spline_alloc(gsl_interp_linear, LEN(absorption_coefficient_h2o_wavelengths));
gsl_spline_init(spline, absorption_coefficient_h2o_wavelengths, absorption_coefficient_h2o, LEN(absorption_coefficient_h2o_wavelengths));
acc = gsl_interp_accel_alloc();
}
/* If the wavelength is outside the range of the lookup table we just
* return the value at the beginning or end. We do this because by default
* gsl_spline_eval will exit if the wavelength is out of bounds. */
if (wavelength < absorption_coefficient_h2o_wavelengths[0])
wavelength = absorption_coefficient_h2o_wavelengths[0];
if (wavelength > absorption_coefficient_h2o_wavelengths[LEN(absorption_coefficient_h2o_wavelengths)-1])
wavelength = absorption_coefficient_h2o_wavelengths[LEN(absorption_coefficient_h2o_wavelengths)-1];
return 1.0/gsl_spline_eval(spline, wavelength, acc);
}
/* Returns the absorption length as a function of wavelength in acrylic.
*
* The wavelength should be in nm and the absorption length is in cm. */
double get_absorption_length_snoman_acrylic(double wavelength)
{
/* Note: We use GSL splines here instead of interp1d() since we aren't
* guaranteed that the values in the lookup table are evenly spaced. */
static gsl_spline *spline;
static gsl_interp_accel *acc;
if (!spline) {
spline = gsl_spline_alloc(gsl_interp_linear, LEN(absorption_coefficient_acrylic_wavelengths));
gsl_spline_init(spline, absorption_coefficient_acrylic_wavelengths, absorption_coefficient_acrylic, LEN(absorption_coefficient_acrylic_wavelengths));
acc = gsl_interp_accel_alloc();
}
/* If the wavelength is outside the range of the lookup table we just
* return the value at the beginning or end. We do this because by default
* gsl_spline_eval will exit if the wavelength is out of bounds. */
if (wavelength < absorption_coefficient_acrylic_wavelengths[0])
wavelength = absorption_coefficient_acrylic_wavelengths[0];
if (wavelength > absorption_coefficient_acrylic_wavelengths[LEN(absorption_coefficient_acrylic_wavelengths)-1])
wavelength = absorption_coefficient_acrylic_wavelengths[LEN(absorption_coefficient_acrylic_wavelengths)-1];
return 1.0/gsl_spline_eval(spline, wavelength, acc);
}
static double gsl_rayleigh_scattering_length_snoman_d2o(double wavelength, void *params)
{
/* Returns the Rayleigh scattering length in D2O weighted by the quantum
* efficiency and the Cerenkov spectrum. */
double qe, x;
x = *((double *) params);
qe = get_quantum_efficiency(wavelength);
return qe*exp(-x*rayint_prob_d2o(wavelength))/pow(wavelength,2);
}
static double gsl_rayleigh_scattering_length_snoman_h2o(double wavelength, void *params)
{
/* Returns the Rayleigh scattering length in H2O weighted by the quantum
* efficiency and the Cerenkov spectrum. */
double qe, x;
x = *((double *) params);
qe = get_quantum_efficiency(wavelength);
return qe*exp(-x*rayint_prob_h2o(wavelength))/pow(wavelength,2);
}
static double gsl_rayleigh_scattering_length_snoman_acrylic(double wavelength, void *params)
{
/* Returns the Rayleigh scattering length in acrylic weighted by the
* quantum efficiency and the Cerenkov spectrum. */
double qe, x;
x = *((double *) params);
qe = get_quantum_efficiency(wavelength);
return qe*exp(-x*rayint_prob_acrylic(wavelength))/pow(wavelength,2);
}
static double gsl_absorption_length_snoman_d2o(double wavelength, void *params)
{
/* Returns the absorption length in D2O weighted by the quantum efficiency
* and the Cerenkov spectrum. */
double qe, x;
x = *((double *) params);
qe = get_quantum_efficiency(wavelength);
return qe*exp(-x/get_absorption_length_snoman_d2o(wavelength))/pow(wavelength,2);
}
static double gsl_absorption_length_snoman_h2o(double wavelength, void *params)
{
/* Returns the absorption length in H2O weighted by the quantum efficiency
* and the Cerenkov spectrum. */
double qe, x;
x = *((double *) params);
qe = get_quantum_efficiency(wavelength);
return qe*exp(-x/get_absorption_length_snoman_h2o(wavelength))/pow(wavelength,2);
}
static double gsl_absorption_length_snoman_acrylic(double wavelength, void *params)
{
/* Returns the absorption length in acrylic weighted by the quantum
* efficiency and the Cerenkov spectrum. */
double qe, x;
x = *((double *) params);
qe = get_quantum_efficiency(wavelength);
return qe*exp(-x/get_absorption_length_snoman_acrylic(wavelength))/pow(wavelength,2);
}
/* Returns the product of the index of refraction in D2O multiplied by the
* quantum efficiency and the Cerenkov spectrum. */
static double gsl_index_d2o(double wavelength, void *params)
{
double qe;
qe = get_quantum_efficiency(wavelength);
return qe*get_index_snoman_d2o(wavelength)/pow(wavelength,2);
}
/* Returns the product of the index of refraction in H2O multiplied by the
* quantum efficiency and the Cerenkov spectrum. */
static double gsl_index_h2o(double wavelength, void *params)
{
double qe;
qe = get_quantum_efficiency(wavelength);
return qe*get_index_snoman_h2o(wavelength)/pow(wavelength,2);
}
static double gsl_cerenkov(double wavelength, void *params)
{
/* Returns the quantum efficiency multiplied by the Cerenkov spectrum. */
double qe;
qe = get_quantum_efficiency(wavelength);
return qe/pow(wavelength,2);
}
int optics_init(int run)
{
int i;
double norm;
double result, error;
size_t nevals;
int status;
gsl_integration_cquad_workspace *w;
gsl_function F;
if (run < 19999) {
for (i = 0; i < LEN(absorption_coefficient_d2o_wavelengths); i++) {
absorption_coefficient_d2o_wavelengths[i] = d2o_constants.absorption_coefficient_d2o_wavelengths[i];
absorption_coefficient_d2o[i] = d2o_constants.absorption_coefficient_d2o[i];
}
rayleigh_scl_fac_d2o = d2o_constants.rayleigh_scl_fac_d2o;
isothermal_comp_d2o = d2o_constants.isothermal_comp_d2o;
for (i = 0; i < LEN(absorption_coefficient_h2o_wavelengths); i++) {
absorption_coefficient_h2o_wavelengths[i] = d2o_constants.absorption_coefficient_h2o_wavelengths[i];
absorption_coefficient_h2o[i] = d2o_constants.absorption_coefficient_h2o[i];
}
rayleigh_scl_fac_h2o = d2o_constants.rayleigh_scl_fac_h2o;
isothermal_comp_h2o = d2o_constants.isothermal_comp_h2o;
for (i = 0; i < LEN(absorption_coefficient_acrylic_wavelengths); i++) {
absorption_coefficient_acrylic_wavelengths[i] = d2o_constants.absorption_coefficient_acrylic_wavelengths[i];
absorption_coefficient_acrylic[i] = d2o_constants.absorption_coefficient_acrylic[i];
}
rayleigh_scl_fac_acrylic = d2o_constants.rayleigh_scl_fac_acrylic;
isothermal_comp_acrylic = d2o_constants.isothermal_comp_acrylic;
} else if (run < 33907) {
for (i = 0; i < LEN(absorption_coefficient_d2o_wavelengths); i++) {
absorption_coefficient_d2o_wavelengths[i] = salt_constants.absorption_coefficient_d2o_wavelengths[i];
absorption_coefficient_d2o[i] = salt_constants.absorption_coefficient_d2o[i];
}
rayleigh_scl_fac_d2o = salt_constants.rayleigh_scl_fac_d2o;
isothermal_comp_d2o = salt_constants.isothermal_comp_d2o;
for (i = 0; i < LEN(absorption_coefficient_h2o_wavelengths); i++) {
absorption_coefficient_h2o_wavelengths[i] = salt_constants.absorption_coefficient_h2o_wavelengths[i];
absorption_coefficient_h2o[i] = salt_constants.absorption_coefficient_h2o[i];
}
rayleigh_scl_fac_h2o = salt_constants.rayleigh_scl_fac_h2o;
isothermal_comp_h2o = salt_constants.isothermal_comp_h2o;
for (i = 0; i < LEN(absorption_coefficient_acrylic_wavelengths); i++) {
absorption_coefficient_acrylic_wavelengths[i] = salt_constants.absorption_coefficient_acrylic_wavelengths[i];
absorption_coefficient_acrylic[i] = salt_constants.absorption_coefficient_acrylic[i];
}
rayleigh_scl_fac_acrylic = salt_constants.rayleigh_scl_fac_acrylic;
isothermal_comp_acrylic = salt_constants.isothermal_comp_acrylic;
} else {
sprintf(optics_err, "no optics info for ncd or h2o phase!");
return -1;
}
w = gsl_integration_cquad_workspace_alloc(100);
F.function = &gsl_cerenkov;
F.params = NULL;
status = gsl_integration_cquad(&F, 200, 800, 0, 1e-2, w, &norm, &error, &nevals);
if (status) {
sprintf(optics_err, "error integrating cerenkov distribution: %s\n", gsl_strerror(status));
return -1;
}
F.function = &gsl_index_d2o;
F.params = NULL;
status = gsl_integration_cquad(&F, 200, 800, 0, 1e-2, w, &result, &error, &nevals);
if (status) {
sprintf(optics_err, "error integrating cerenkov distribution: %s\n", gsl_strerror(status));
return -1;
}
avg_index_d2o = result/norm;
F.function = &gsl_index_h2o;
F.params = NULL;
status = gsl_integration_cquad(&F, 200, 800, 0, 1e-2, w, &result, &error, &nevals);
if (status) {
sprintf(optics_err, "error integrating cerenkov distribution: %s\n", gsl_strerror(status));
return -1;
}
avg_index_h2o = result/norm;
for (i = 0; i < LEN(fabs_x); i++) {
fabs_x[i] = xlo + (xhi-xlo)*i/(LEN(fabs_x)-1);
}
F.function = &gsl_absorption_length_snoman_d2o;
for (i = 0; i < LEN(fabs_x); i++) {
F.params = fabs_x+i;
status = gsl_integration_cquad(&F, 200, 800, 0, 1e-2, w, &result, &error, &nevals);
if (status) {
sprintf(optics_err, "error integrating cerenkov distribution: %s\n", gsl_strerror(status));
return -1;
}
fabs_d2o[i] = result/norm;
}
F.function = &gsl_absorption_length_snoman_h2o;
for (i = 0; i < LEN(fabs_x); i++) {
F.params = fabs_x+i;
status = gsl_integration_cquad(&F, 200, 800, 0, 1e-2, w, &result, &error, &nevals);
if (status) {
sprintf(optics_err, "error integrating cerenkov distribution: %s\n", gsl_strerror(status));
return -1;
}
fabs_h2o[i] = result/norm;
}
F.function = &gsl_absorption_length_snoman_acrylic;
for (i = 0; i < LEN(fabs_x); i++) {
F.params = fabs_x+i;
status = gsl_integration_cquad(&F, 200, 800, 0, 1e-2, w, &result, &error, &nevals);
if (status) {
sprintf(optics_err, "error integrating cerenkov distribution: %s\n", gsl_strerror(status));
return -1;
}
fabs_acrylic[i] = result/norm;
}
F.function = &gsl_rayleigh_scattering_length_snoman_d2o;
for (i = 0; i < LEN(fabs_x); i++) {
F.params = fabs_x+i;
status = gsl_integration_cquad(&F, 200, 800, 0, 1e-2, w, &result, &error, &nevals);
if (status) {
sprintf(optics_err, "error integrating cerenkov distribution: %s\n", gsl_strerror(status));
return -1;
}
fsct_d2o[i] = result/norm;
}
F.function = &gsl_rayleigh_scattering_length_snoman_h2o;
for (i = 0; i < LEN(fabs_x); i++) {
F.params = fabs_x+i;
status = gsl_integration_cquad(&F, 200, 800, 0, 1e-2, w, &result, &error, &nevals);
if (status) {
sprintf(optics_err, "error integrating cerenkov distribution: %s\n", gsl_strerror(status));
return -1;
}
fsct_h2o[i] = result/norm;
}
F.function = &gsl_rayleigh_scattering_length_snoman_acrylic;
for (i = 0; i < LEN(fabs_x); i++) {
F.params = fabs_x+i;
status = gsl_integration_cquad(&F, 200, 800, 0, 1e-2, w, &result, &error, &nevals);
if (status) {
sprintf(optics_err, "error integrating cerenkov distribution: %s\n", gsl_strerror(status));
return -1;
}
fsct_acrylic[i] = result/norm;
}
gsl_integration_cquad_workspace_free(w);
initialized = 1;
return 0;
}
|