diff options
author | tlatorre <tlatorre@uchicago.edu> | 2018-08-14 09:53:09 -0500 |
---|---|---|
committer | tlatorre <tlatorre@uchicago.edu> | 2018-08-14 09:53:09 -0500 |
commit | 0b7f199c0d93074484ea580504485a32dc29f5e2 (patch) | |
tree | e167b6d102b87b7a5eca4558e7f39265d5edc502 /scattering.c | |
parent | 636595905c9f63e6bfcb6d331312090ac2075377 (diff) | |
download | sddm-0b7f199c0d93074484ea580504485a32dc29f5e2.tar.gz sddm-0b7f199c0d93074484ea580504485a32dc29f5e2.tar.bz2 sddm-0b7f199c0d93074484ea580504485a32dc29f5e2.zip |
initial commit of likelihood fit for muons
This commit contains code to fit for the energy, position, and direction of
muons in the SNO detector. Currently, we read events from SNOMAN zebra files
and fill an event struct containing the PMT hits and fit it with the Nelder
Mead simplex algorithm from GSL.
I've also added code to read in ZEBRA title bank files to read in the DQXX
files for a specific run. Any problems with channels in the DQCH and DQCR banks
are flagged in the event struct by masking in a bit in the flags variable and
these PMT hits are not included in the likelihood calculation.
The likelihood for an event is calculated by integrating along the particle
track for each PMT and computing the expected number of PE. The charge
likelihood is then calculated by looping over all possible number of PE and
computing:
P(q|n)*P(n|mu)
where q is the calibrated QHS charge, n is the number of PE, and mu is the
expected number of photoelectrons. The latter is calculated assuming the
distribution of PE at a given PMT follows a Poisson distribution (which I think
should be correct given the track, but is probably not perfect for tracks which
scatter a lot).
The time part of the likelihood is calculated by integrating over the track for
each PMT and calculating the average time at which the PMT is hit. We then
assume the PDF for the photons to arrive is approximately a delta function and
compute the first order statistic for a given time to compute the probability
that the first photon arrived at a given time. So far I've only tested this
with single tracks but the method was designed to be easy to use when you are
fitting for multiple particles.
Diffstat (limited to 'scattering.c')
-rw-r--r-- | scattering.c | 121 |
1 files changed, 121 insertions, 0 deletions
diff --git a/scattering.c b/scattering.c new file mode 100644 index 0000000..a7fa717 --- /dev/null +++ b/scattering.c @@ -0,0 +1,121 @@ +#include "scattering.h" +#include "quantum_efficiency.h" +#include <gsl/gsl_sf_bessel.h> +#include <gsl/gsl_integration.h> +#include <gsl/gsl_interp2d.h> +#include <gsl/gsl_spline2d.h> +#include "optics.h" +#include "sno.h" +#include <stddef.h> /* for size_t */ +#include <stdlib.h> /* for exit() */ +#include <stdio.h> /* for fprintf() */ +#include <gsl/gsl_errno.h> /* for gsl_strerror() */ + +static double xlo = -1.0; +static double xhi = 1.0; +static size_t nx = 200; +static double ylo = 0.0; +static double yhi = 1.0; +static size_t ny = 1000; + +static double *x; +static double *y; +static double *z; + +static gsl_spline2d *spline; +static gsl_interp_accel *xacc; +static gsl_interp_accel *yacc; + +static double prob_scatter(double wavelength, void *params) +{ + /* Calculate the number of photons emitted per unit solid angle per cm at + * an angle `theta` for a particle travelling with velocity `beta` with an + * angular distribution of width `theta0`. */ + double qe, delta, index; + + double beta_cos_theta = ((double *) params)[0]; + double beta_sin_theta_theta0 = ((double *) params)[1]; + + qe = get_quantum_efficiency(wavelength); + + index = get_index(HEAVY_WATER_DENSITY, wavelength, 10.0); + + delta = (1.0/index - beta_cos_theta)/(2*beta_sin_theta_theta0); + + /* FIXME: ignore GSL error for underflow here. */ + if (delta*delta > 500) return 0.0; + else if (delta*delta == 0.0) return INFINITY; + return qe*exp(-delta*delta)*gsl_sf_bessel_K0(delta*delta)/pow(wavelength,2)*1e7/(4*M_PI*M_PI); +} + +void init_interpolation(void) +{ + size_t i, j; + double params[2]; + double result, error; + size_t nevals; + int status; + + x = malloc(nx*sizeof(double)); + y = malloc(ny*sizeof(double)); + z = malloc(nx*ny*sizeof(double)); + + spline = gsl_spline2d_alloc(gsl_interp2d_bilinear, nx, ny); + xacc = gsl_interp_accel_alloc(); + yacc = gsl_interp_accel_alloc(); + + for (i = 0; i < nx; i++) { + x[i] = xlo + (xhi-xlo)*i/(nx-1); + } + + for (i = 0; i < ny; i++) { + y[i] = ylo + (yhi-ylo)*i/(ny-1); + } + + gsl_integration_cquad_workspace *w = gsl_integration_cquad_workspace_alloc(100); + + gsl_function F; + F.function = &prob_scatter; + F.params = params; + + for (i = 0; i < nx; i++) { + for (j = 0; j < ny; j++) { + params[0] = x[i]; + params[1] = y[j]; + + status = gsl_integration_cquad(&F, 200, 800, 0, 1e-2, w, &result, &error, &nevals); + + if (status) { + fprintf(stderr, "error integrating photon angular distribution: %s\n", gsl_strerror(status)); + exit(1); + } + gsl_spline2d_set(spline, z, i, j, result); + } + } + + gsl_spline2d_init(spline, x, y, z, nx, ny); + + gsl_integration_cquad_workspace_free(w); +} + +double get_probability(double beta, double cos_theta, double theta0) +{ + double sin_theta; + + /* Technically this isn't defined up to a sign, but it doesn't matter since + * we are going to square it everywhere. */ + sin_theta = fabs(sin(acos(cos_theta))); + + return gsl_spline2d_eval(spline, beta*cos_theta, beta*sin_theta*theta0, xacc, yacc)/sin_theta; +} + +void free_interpolation(void) +{ + free(x); + free(y); + free(z); + + gsl_spline2d_free(spline); + gsl_interp_accel_free(xacc); + gsl_interp_accel_free(yacc); +} |