aboutsummaryrefslogtreecommitdiff
path: root/utils/dc-closure-test
blob: c09e5f93616c9e680aca613e0a21fd66691ce4d1 (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
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
#!/usr/bin/env python3
# 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/>.

from __future__ import print_function, division
import numpy as np
from scipy.stats import iqr
import nlopt
from scipy.stats import poisson, norm
import contextlib
import sys
from math import exp
import emcee
from scipy.optimize import brentq
from scipy.stats import truncnorm
from matplotlib.lines import Line2D
from sddm.plot import despine
from sddm.dc import *
from sddm.plot_energy import *

try:
    from emcee import moves
except ImportError:
    print("emcee version 2.2.1 is required",file=sys.stderr)
    sys.exit(1)

# from https://stackoverflow.com/questions/2891790/how-to-pretty-print-a-numpy-array-without-scientific-notation-and-with-given-pre
@contextlib.contextmanager
def printoptions(*args, **kwargs):
    original = np.get_printoptions()
    np.set_printoptions(*args, **kwargs)
    try:
        yield
    finally: 
        np.set_printoptions(**original)

def radius_cut(ev):
    ev['radius_cut'] = np.digitize((ev.r/PSUP_RADIUS)**3,(0.9,))
    return ev

def udotr_cut(ev):
    ev['udotr_cut'] = np.digitize(ev.udotr,(-0.5,))
    return ev

def psi_cut(ev):
    ev['psi_cut'] = np.digitize(ev.psi,(6.0,))
    return ev

def cos_theta_cut(ev):
    ev['cos_theta_cut'] = np.digitize(ev.cos_theta,(-0.5,))
    return ev

def z_cut(ev):
    ev['z_cut'] = np.digitize(ev.z,(0.0,))
    return ev

# Constraint to enforce the fact that P(r,psi,z,udotr|muon) all add up to 1.0.
# In the likelihood function we set the last possibility for r and udotr equal
# to 1.0 minus the others. Therefore, we need to enforce the fact that the
# others must add up to less than 1.
muon_r_psi_z_udotr = Constraint(range(11,26))

# Constraint to enforce the fact that P(z,udotr|noise) all add up to 1.0. In
# the likelihood function we set the last possibility for r and udotr equal to
# 1.0 minus the others. Therefore, we need to enforce the fact that the others
# must add up to less than 1.
noise_z_udotr = Constraint(range(28,31))

# Constraint to enforce the fact that P(r,z,udotr|neck) all add up to 1.0. In
# the likelihood function we set the last possibility for r and udotr equal to
# 1.0 minus the others. Therefore, we need to enforce the fact that the others
# must add up to less than 1.
neck_r_z_udotr = Constraint(range(31,38))

# Constraint to enforce the fact that P(r,udotr|flasher) all add up to 1.0. In
# the likelihood function we set the last possibility for r and udotr equal to
# 1.0 minus the others. Therefore, we need to enforce the fact that the others
# must add up to less than 1
flasher_r_udotr = Constraint(range(39,42))

# Constraint to enforce the fact that P(r,udotr|breakdown) all add up to 1.0.
# In the likelihood function we set the last possibility for r and udotr equal
# to 1.0 minus the others. Therefore, we need to enforce the fact that the
# others must add up to less than 1.
breakdown_r_udotr = Constraint(range(44,47))

def make_nll(data, sacrifice, constraints, fitted_fraction):
    def nll(x, grad=None, fill_value=1e9):
        if grad is not None and grad.size > 0:
            raise Exception("nll got passed grad!")

        nll = 0.0
        # Here we explicitly return a crazy high value if one of the
        # constraints is violated. When using nlopt it should respect all the
        # constraints, *but* later when we do the Metropolis Hastings algorithm
        # we don't have any way to add the constraints explicitly.
        for constraint in constraints:
            if constraint(x) > 0:
                nll += fill_value + 1e4*constraint(x)**2

        if (x <= 0).any() or (x[6:] >= 1).any():
            nll += fill_value + 1e4*np.sum((x[x < 0])**2) + 1e4*np.sum((x[6:][x[6:] > 1]-1)**2)

        if nll:
            return nll

        (mu_signal, mu_muon, mu_noise, mu_neck, mu_flasher, mu_breakdown,
         contamination_muon, contamination_noise, contamination_neck, contamination_flasher, contamination_breakdown,
         p_r_psi_z_udotr_muon_lolololo, # 11
         p_r_psi_z_udotr_muon_lololohi,
         p_r_psi_z_udotr_muon_lolohilo,
         p_r_psi_z_udotr_muon_lolohihi,
         p_r_psi_z_udotr_muon_lohilolo,
         p_r_psi_z_udotr_muon_lohilohi,
         p_r_psi_z_udotr_muon_lohihilo,
         p_r_psi_z_udotr_muon_lohihihi,
         p_r_psi_z_udotr_muon_hilololo,
         p_r_psi_z_udotr_muon_hilolohi,
         p_r_psi_z_udotr_muon_hilohilo,
         p_r_psi_z_udotr_muon_hilohihi,
         p_r_psi_z_udotr_muon_hihilolo,
         p_r_psi_z_udotr_muon_hihilohi,
         p_r_psi_z_udotr_muon_hihihilo,
         p_r_noise_lo, p_psi_noise_lo, # 26, 27
         p_z_udotr_noise_lolo, # 28
         p_z_udotr_noise_lohi,
         p_z_udotr_noise_hilo,
         p_r_z_udotr_neck_lololo, # 31
         p_r_z_udotr_neck_lolohi,
         p_r_z_udotr_neck_lohilo,
         p_r_z_udotr_neck_lohihi,
         p_r_z_udotr_neck_hilolo,
         p_r_z_udotr_neck_hilohi,
         p_r_z_udotr_neck_hihilo,
         p_psi_neck_lo, # 38
         p_r_udotr_flasher_lolo, p_r_udotr_flasher_lohi, p_r_udotr_flasher_hilo, # 39, ..., 41
         p_psi_flasher_lo, p_z_flasher_lo,
         p_r_udotr_breakdown_lolo, p_r_udotr_breakdown_lohi, p_r_udotr_breakdown_hilo, # 44, ..., 46
         p_psi_breakdown_lo, p_z_breakdown_lo,
         p_neck_given_muon) = x

        p_r_udotr_flasher_hihi = 1-p_r_udotr_flasher_lolo-p_r_udotr_flasher_lohi-p_r_udotr_flasher_hilo
        p_r_udotr_breakdown_hihi = 1-p_r_udotr_breakdown_lolo-p_r_udotr_breakdown_lohi-p_r_udotr_breakdown_hilo
        p_r_psi_z_udotr_muon_hihihihi = 1 - \
            p_r_psi_z_udotr_muon_lolololo - \
            p_r_psi_z_udotr_muon_lololohi - \
            p_r_psi_z_udotr_muon_lolohilo - \
            p_r_psi_z_udotr_muon_lolohihi - \
            p_r_psi_z_udotr_muon_lohilolo - \
            p_r_psi_z_udotr_muon_lohilohi - \
            p_r_psi_z_udotr_muon_lohihilo - \
            p_r_psi_z_udotr_muon_lohihihi - \
            p_r_psi_z_udotr_muon_hilololo - \
            p_r_psi_z_udotr_muon_hilolohi - \
            p_r_psi_z_udotr_muon_hilohilo - \
            p_r_psi_z_udotr_muon_hilohihi - \
            p_r_psi_z_udotr_muon_hihilolo - \
            p_r_psi_z_udotr_muon_hihilohi - \
            p_r_psi_z_udotr_muon_hihihilo
        p_r_z_udotr_neck_hihihi = 1 - p_r_z_udotr_neck_lololo - p_r_z_udotr_neck_lolohi - p_r_z_udotr_neck_lohilo - p_r_z_udotr_neck_lohihi - p_r_z_udotr_neck_hilolo - p_r_z_udotr_neck_hilohi - p_r_z_udotr_neck_hihilo
        p_z_udotr_noise_hihi = 1 - p_z_udotr_noise_lolo - p_z_udotr_noise_lohi - p_z_udotr_noise_hilo

        # Muon events
        # first 6 parameters are the mean number of signal and bgs
        p_muon = np.array([\
            [[[p_r_psi_z_udotr_muon_lolololo, p_r_psi_z_udotr_muon_lololohi], \
              [p_r_psi_z_udotr_muon_lolohilo, p_r_psi_z_udotr_muon_lolohihi]], \
             [[p_r_psi_z_udotr_muon_lohilolo, p_r_psi_z_udotr_muon_lohilohi], \
              [p_r_psi_z_udotr_muon_lohihilo, p_r_psi_z_udotr_muon_lohihihi]]], \
            [[[p_r_psi_z_udotr_muon_hilololo, p_r_psi_z_udotr_muon_hilolohi], \
              [p_r_psi_z_udotr_muon_hilohilo, p_r_psi_z_udotr_muon_hilohihi]], \
             [[p_r_psi_z_udotr_muon_hihilolo, p_r_psi_z_udotr_muon_hihilohi], \
              [p_r_psi_z_udotr_muon_hihihilo, p_r_psi_z_udotr_muon_hihihihi]]]])
        expected_muon = p_muon*contamination_muon*mu_muon*fitted_fraction['muon'] + sacrifice['muon']*mu_signal

        nll -= fast_poisson_logpmf(data['muon'],expected_muon).sum()

        # Noise events
        p_r_noise = np.array([p_r_noise_lo,1-p_r_noise_lo])
        p_psi_noise = np.array([p_psi_noise_lo,1-p_psi_noise_lo])
        p_z_udotr_noise = np.array([\
            [p_z_udotr_noise_lolo,p_z_udotr_noise_lohi],
            [p_z_udotr_noise_hilo,p_z_udotr_noise_hihi]])
        p_noise = p_r_noise[:,np.newaxis,np.newaxis,np.newaxis]*p_psi_noise[:,np.newaxis,np.newaxis]*p_z_udotr_noise
        expected_noise = p_noise*contamination_noise*mu_noise*fitted_fraction['noise'] + sacrifice['noise']*mu_signal

        nll -= fast_poisson_logpmf(data['noise'],expected_noise).sum()

        # Neck events
        # FIXME: for now assume parameterized same as muon
        p_r_z_udotr_neck = np.array([\
            [[p_r_z_udotr_neck_lololo, p_r_z_udotr_neck_lolohi], \
             [p_r_z_udotr_neck_lohilo, p_r_z_udotr_neck_lohihi]], \
            [[p_r_z_udotr_neck_hilolo, p_r_z_udotr_neck_hilohi], \
             [p_r_z_udotr_neck_hihilo, p_r_z_udotr_neck_hihihi]]])
        p_psi_neck = np.array([p_psi_neck_lo,1-p_psi_neck_lo])
        p_neck = p_r_z_udotr_neck[:,np.newaxis,:,:]*p_psi_neck[:,np.newaxis,np.newaxis]
        expected_neck = p_neck*contamination_neck*mu_neck*fitted_fraction['neck'] + sacrifice['neck']*mu_signal
        # FIXME: pdf should be different for muon given neck
        expected_neck += p_muon*p_neck_given_muon*mu_muon*fitted_fraction['neck']

        nll -= fast_poisson_logpmf(data['neck'],expected_neck).sum()

        # Flasher events
        p_r_udotr_flasher = np.array([\
            [p_r_udotr_flasher_lolo,p_r_udotr_flasher_lohi], \
            [p_r_udotr_flasher_hilo,p_r_udotr_flasher_hihi]])
        p_psi_flasher = np.array([p_psi_flasher_lo,1-p_psi_flasher_lo])
        p_z_flasher = np.array([p_z_flasher_lo,1-p_z_flasher_lo])
        p_flasher = p_r_udotr_flasher[:,np.newaxis,np.newaxis,:]*p_psi_flasher[:,np.newaxis,np.newaxis]*p_z_flasher[:,np.newaxis]
        expected_flasher = p_flasher*contamination_flasher*mu_flasher*fitted_fraction['flasher'] + sacrifice['flasher']*mu_signal

        nll -= fast_poisson_logpmf(data['flasher'],expected_flasher).sum()

        # Breakdown events
        p_r_udotr_breakdown = np.array([\
            [p_r_udotr_breakdown_lolo,p_r_udotr_breakdown_lohi], \
            [p_r_udotr_breakdown_hilo,p_r_udotr_breakdown_hihi]])
        p_psi_breakdown = np.array([p_psi_breakdown_lo,1-p_psi_breakdown_lo])
        p_z_breakdown = np.array([p_z_breakdown_lo,1-p_z_breakdown_lo])
        p_breakdown = p_r_udotr_breakdown[:,np.newaxis,np.newaxis,:]*p_psi_breakdown[:,np.newaxis,np.newaxis]*p_z_breakdown[:,np.newaxis]
        expected_breakdown = p_breakdown*contamination_breakdown*mu_breakdown*fitted_fraction['breakdown'] + sacrifice['breakdown']*mu_signal

        nll -= fast_poisson_logpmf(data['breakdown'],expected_breakdown).sum()

        # Signal like events
        expected_signal = np.zeros_like(expected_muon)
        expected_signal += mu_signal*sacrifice['signal']
        expected_signal += p_muon*(1-contamination_muon)*mu_muon
        expected_signal += p_neck*(1-contamination_neck)*mu_neck
        expected_signal += p_noise*(1-contamination_noise)*mu_noise
        expected_signal += p_flasher*(1-contamination_flasher)*mu_flasher
        expected_signal += p_breakdown*(1-contamination_breakdown)*mu_breakdown

        nll -= fast_poisson_logpmf(data['signal'],expected_signal).sum()

        if not np.isfinite(nll):
            print("x = ", x)
            print("p_r_z_udotr_neck = ", p_r_z_udotr_neck)
            print("expected_muon = ", expected_muon)
            print("expected_noise = ", expected_noise)
            print("expected_neck = ", expected_neck)
            print("expected_flasher = ", expected_flasher)
            print("expected_breakdown = ", expected_breakdown)
            print("nll is not finite!")
            sys.exit(0)

        return nll
    return nll

def fit(data, sacrifice, steps):
    constraints = [flasher_r_udotr, breakdown_r_udotr,muon_r_psi_z_udotr,neck_r_z_udotr,noise_z_udotr]
    nll = make_nll(data,sacrifice,constraints,fitted_fraction)

    x0 = []
    for bg in ['signal','muon','noise','neck','flasher','breakdown']:
        x0.append(data[bg].sum())

    # contamination
    x0 += [0.99]*5

    if data['muon'].sum() > 0:
        # P(r,psi,z,udotr|muon)
        x0 += [data['muon'][0,0,0,0].sum()/data['muon'].sum()]
        x0 += [data['muon'][0,0,0,1].sum()/data['muon'].sum()]
        x0 += [data['muon'][0,0,1,0].sum()/data['muon'].sum()]
        x0 += [data['muon'][0,0,1,1].sum()/data['muon'].sum()]
        x0 += [data['muon'][0,1,0,0].sum()/data['muon'].sum()]
        x0 += [data['muon'][0,1,0,1].sum()/data['muon'].sum()]
        x0 += [data['muon'][0,1,1,0].sum()/data['muon'].sum()]
        x0 += [data['muon'][0,1,1,1].sum()/data['muon'].sum()]
        x0 += [data['muon'][1,0,0,0].sum()/data['muon'].sum()]
        x0 += [data['muon'][1,0,0,1].sum()/data['muon'].sum()]
        x0 += [data['muon'][1,0,1,0].sum()/data['muon'].sum()]
        x0 += [data['muon'][1,0,1,1].sum()/data['muon'].sum()]
        x0 += [data['muon'][1,1,0,0].sum()/data['muon'].sum()]
        x0 += [data['muon'][1,1,0,1].sum()/data['muon'].sum()]
        x0 += [data['muon'][1,1,1,0].sum()/data['muon'].sum()]
    else:
        x0 += [0.1]*15

    if data['noise'].sum() > 0:
        # P(r|noise)
        x0 += [data['noise'][0].sum()/data['noise'].sum()]
        # P(psi|noise)
        x0 += [data['noise'][:,0].sum()/data['noise'].sum()]
        # P(z,udotr|noise)
        x0 += [data['noise'][:,:,0,0].sum()/data['noise'].sum()]
        x0 += [data['noise'][:,:,0,1].sum()/data['noise'].sum()]
        x0 += [data['noise'][:,:,1,0].sum()/data['noise'].sum()]
    else:
        x0 += [0.1]*5

    if data['neck'].sum() > 0:
        # P(r,z,udotr|neck)
        x0 += [data['neck'][0,:,0,0].sum()/data['neck'].sum()]
        x0 += [data['neck'][0,:,0,1].sum()/data['neck'].sum()]
        x0 += [data['neck'][0,:,1,0].sum()/data['neck'].sum()]
        x0 += [data['neck'][0,:,1,1].sum()/data['neck'].sum()]
        x0 += [data['neck'][1,:,0,0].sum()/data['neck'].sum()]
        x0 += [data['neck'][1,:,0,1].sum()/data['neck'].sum()]
        x0 += [data['neck'][1,:,1,0].sum()/data['neck'].sum()]
        # P(psi|neck)
        x0 += [data['neck'][:,0].sum()/data['neck'].sum()]
    else:
        x0 += [0.1]*8

    if data['flasher'].sum() > 0:
        # P(r,udotr|flasher)
        x0 += [data['flasher'][0,:,:,0].sum()/data['flasher'].sum()]
        x0 += [data['flasher'][0,:,:,1].sum()/data['flasher'].sum()]
        x0 += [data['flasher'][1,:,:,0].sum()/data['flasher'].sum()]
        # P(psi|flasher)
        x0 += [data['flasher'][:,0].sum()/data['flasher'].sum()]
        # P(z|flasher)
        x0 += [data['flasher'][:,:,0].sum()/data['flasher'].sum()]
    else:
        x0 += [0.1]*5

    if data['breakdown'].sum() > 0:
        # P(r,udotr|breakdown)
        x0 += [data['breakdown'][0,:,:,0].sum()/data['breakdown'].sum()]
        x0 += [data['breakdown'][0,:,:,1].sum()/data['breakdown'].sum()]
        x0 += [data['breakdown'][1,:,:,0].sum()/data['breakdown'].sum()]
        # P(psi|breakdown)
        x0 += [data['breakdown'][:,0].sum()/data['breakdown'].sum()]
        # P(z|breakdown)
        x0 += [data['breakdown'][:,:,0].sum()/data['breakdown'].sum()]
    else:
        x0 += [0.1]*5

    # P(neck|muon)
    x0 += [EPSILON]

    x0 = np.array(x0)

    # Use the COBYLA algorithm here because it is the only derivative free
    # minimization routine which honors inequality constraints
    # Edit: SBPLX seems to work better
    opt = nlopt.opt(nlopt.LN_SBPLX, len(x0))
    opt.set_min_objective(nll)
    # set lower bounds to 1e-10 to prevent nans if we predict something should
    # be 0 but observe an event.
    low = np.ones_like(x0)*EPSILON
    high = np.array([1e9]*6 + [1-EPSILON]*(len(x0)-6))
    x0[x0 < low] = low[x0 < low]
    x0[x0 > high] = high[x0 > high]
    opt.set_lower_bounds(low)
    opt.set_upper_bounds(high)
    opt.set_ftol_abs(1e-10)
    opt.set_initial_step([1]*6 + [0.01]*(len(x0)-6))
    #for constraint in constraints:
        #opt.add_inequality_constraint(constraint,0)

    xopt = opt.optimize(x0)
    nll_xopt = nll(xopt)
    print("nll(xopt) = ", nll(xopt))

    while True:
        xopt = opt.optimize(xopt)
        if not nll(xopt) < nll_xopt - 1e-10:
            break
        nll_xopt = nll(xopt)
        print("nll(xopt) = ", nll(xopt))
        #print("n = ", opt.get_numevals())

    stepsizes = estimate_errors(nll,xopt,low,high,constraints)
    with printoptions(precision=3, suppress=True):
        print("Errors: ", stepsizes)

    #samples = metropolis_hastings(nll,xopt,stepsizes,100000)
    #print("nll(xopt) = %.2g" % nll(xopt))

    pos = np.empty((10, len(x0)),dtype=np.double)
    for i in range(pos.shape[0]):
        pos[i] = xopt + np.random.randn(len(x0))*stepsizes
        pos[i,:6] = np.clip(pos[i,:6],EPSILON,1e9)
        pos[i,6:] = np.clip(pos[i,6:],EPSILON,1-EPSILON)

        for constraint in constraints:
            if constraint(pos[i]) >= 0:
                pos[i] = constraint.renormalize_no_fix(pos[i])

    nwalkers, ndim = pos.shape

    proposal = get_proposal_func(stepsizes*0.5,low,high)
    sampler = emcee.EnsembleSampler(nwalkers, ndim, lambda x, grad, fill_value: -nll(x,grad,fill_value), moves=emcee.moves.MHMove(proposal),args=[None,np.inf])
    with np.errstate(invalid='ignore'):
        sampler.run_mcmc(pos, steps)

    print("Mean acceptance fraction: {0:.3f}".format(np.mean(sampler.acceptance_fraction)))

    try:
        print("autocorrelation time: ", sampler.get_autocorr_time(quiet=True))
    except Exception as e:
        print(e)

    samples = sampler.chain.reshape((-1,len(x0)))

    return samples

if __name__ == '__main__':
    import argparse
    import numpy as np
    import pandas as pd
    import sys
    import h5py
    from sddm import setup_matplotlib

    parser = argparse.ArgumentParser("plot fit results")
    parser.add_argument("filenames", nargs='+', help="input files")
    parser.add_argument("--steps", type=int, default=100000, help="number of steps in the MCMC chain")
    parser.add_argument("--save", action="store_true", default=False, help="save plots")
    parser.add_argument("--mc", nargs='+', required=True, help="atmospheric MC files")
    parser.add_argument("-n", type=int, default=10, help="number of fits to run")
    parser.add_argument("--nhit-thresh", type=int, default=None, help="nhit threshold to apply to events before processing (should only be used for testing to speed things up)")
    args = parser.parse_args()

    setup_matplotlib(args.save)

    import matplotlib.pyplot as plt

    # Loop over runs to prevent using too much memory
    evs = []
    rhdr = pd.concat([read_hdf(filename, "rhdr").assign(filename=filename) for filename in args.filenames],ignore_index=True)
    for run, df in rhdr.groupby('run'):
        evs.append(get_events(df.filename.values, merge_fits=True, nhit_thresh=args.nhit_thresh))
    ev = pd.concat(evs)

    # Hack to get rid of flasher and muon events in the breakdown sample.
    ev.breakdown &= ev.r_psup < 0.9

    ev = ev[ev.prompt]
    ev = ev[ev.nhit_cal > 100]

    # Note: Technically we want to know the fitted fraction only for events
    # which *would* reconstruct above 20 MeV. However, there is no way to know
    # if the energy is above 20 MeV without fitting it. However, since I only
    # skip fitting events based on the gtid, there shouldn't be any correlation
    # with energy and so the fitted fraction here should be correct.
    fitted_fraction = {}
    for bg in ['signal','muon','noise','neck','flasher','breakdown']:
        if np.count_nonzero(ev[bg]):
            fitted_fraction[bg] = np.count_nonzero(ev[bg] & ~np.isnan(ev.fmin))/np.count_nonzero(ev[bg])
            print("Fitted fraction for %s: %.0f %%" % (bg,fitted_fraction[bg]*100))
        else:
            print_warning("Warning: No %s events in sample!" % bg)
            sys.exit(1)

    ev = ev[~np.isnan(ev.fmin)]
    ev = ev[ev.ke > 20]

    # figure out bins for high level variables
    ev = radius_cut(ev)
    ev = psi_cut(ev)
    ev = cos_theta_cut(ev)
    ev = z_cut(ev)
    ev = udotr_cut(ev)

    # Note: We loop over the MC filenames here instead of just passing the
    # whole list to get_events() because I had to rerun some of the MC events
    # using SNOMAN and so most of the runs actually have two different files
    # and otherwise the GTIDs will clash
    ev_mcs = []
    for filename in args.mc:
        ev_mcs.append(get_events([filename], merge_fits=True, nhit_thresh=args.nhit_thresh, mc=True))
    ev_mc = pd.concat(ev_mcs)

    ev_mc = ev_mc[ev_mc.prompt]
    ev_mc = ev_mc[ev_mc.nhit_cal > 100]
    ev_mc = ev_mc[~np.isnan(ev_mc.fmin)]
    ev_mc = ev_mc[ev_mc.ke > 20]

    # figure out bins for high level variables
    ev_mc = radius_cut(ev_mc)
    ev_mc = psi_cut(ev_mc)
    ev_mc = cos_theta_cut(ev_mc)
    ev_mc = z_cut(ev_mc)
    ev_mc = udotr_cut(ev_mc)

    contamination_pull = {}

    nbg = {}
    for bg in ['signal','muon','noise','neck','flasher','breakdown']:
        nbg[bg] = min(100,len(ev[ev[bg]]))
        contamination_pull[bg] = []

    for i in range(args.n):
        data = {}
        for bg in ['signal','muon','noise','neck','flasher','breakdown']:
            data[bg] = np.zeros((2,2,2,2),dtype=int)
            if bg == 'signal':
                for bg2 in ['signal','muon','noise','neck','flasher','breakdown']:
                    if bg2 == 'signal':
                        for _, row in ev_mc[ev_mc[bg2]].sample(n=nbg[bg2],replace=True).iterrows():
                            data[bg][row.radius_cut][row.psi_cut][row.z_cut][row.udotr_cut] += 1
                    else:
                        for _, row in ev[ev[bg2]].sample(n=nbg[bg2],replace=True).iterrows():
                            data[bg][row.radius_cut][row.psi_cut][row.z_cut][row.udotr_cut] += 1
            else:
                for _, row in ev[ev[bg]].iterrows():
                    data[bg][row.radius_cut][row.psi_cut][row.z_cut][row.udotr_cut] += 1

        # FIXME: Double check that what I'm calculating here matches with what I
        # expect
        sacrifice = {}
        for bg in ['signal','muon','noise','neck','flasher','breakdown']:
            sacrifice[bg] = np.zeros((2,2,2,2),dtype=float)
            for _, row in ev_mc[ev_mc[bg]].iterrows():
                sacrifice[bg][row.radius_cut][row.psi_cut][row.z_cut][row.udotr_cut] += 1

            sacrifice[bg] /= len(ev_mc)

        samples = fit(data, sacrifice, args.steps)

        (mu_signal, mu_muon, mu_noise, mu_neck, mu_flasher, mu_breakdown,
         contamination_muon, contamination_noise, contamination_neck, contamination_flasher, contamination_breakdown,
         p_r_psi_z_udotr_muon_lolololo, # 11
         p_r_psi_z_udotr_muon_lololohi,
         p_r_psi_z_udotr_muon_lolohilo,
         p_r_psi_z_udotr_muon_lolohihi,
         p_r_psi_z_udotr_muon_lohilolo,
         p_r_psi_z_udotr_muon_lohilohi,
         p_r_psi_z_udotr_muon_lohihilo,
         p_r_psi_z_udotr_muon_lohihihi,
         p_r_psi_z_udotr_muon_hilololo,
         p_r_psi_z_udotr_muon_hilolohi,
         p_r_psi_z_udotr_muon_hilohilo,
         p_r_psi_z_udotr_muon_hilohihi,
         p_r_psi_z_udotr_muon_hihilolo,
         p_r_psi_z_udotr_muon_hihilohi,
         p_r_psi_z_udotr_muon_hihihilo,
         p_r_noise_lo, p_psi_noise_lo, # 26, 27
         p_z_udotr_noise_lolo, # 28
         p_z_udotr_noise_lohi,
         p_z_udotr_noise_hilo,
         p_r_z_udotr_neck_lololo, # 31
         p_r_z_udotr_neck_lolohi,
         p_r_z_udotr_neck_lohilo,
         p_r_z_udotr_neck_lohihi,
         p_r_z_udotr_neck_hilolo,
         p_r_z_udotr_neck_hilohi,
         p_r_z_udotr_neck_hihilo,
         p_psi_neck_lo, # 38
         p_r_udotr_flasher_lolo, p_r_udotr_flasher_lohi, p_r_udotr_flasher_hilo, # 39, ..., 41
         p_psi_flasher_lo, p_z_flasher_lo,
         p_r_udotr_breakdown_lolo, p_r_udotr_breakdown_lohi, p_r_udotr_breakdown_hilo, # 44, ..., 46
         p_psi_breakdown_lo, p_z_breakdown_lo,
         p_neck_given_muon) = samples.T

        for i, bg in enumerate(['signal','muon','noise','neck','flasher','breakdown']):
            if i == 0:
                contamination = samples[:,i]
            else:
                contamination = samples[:,i]*(1-samples[:,5+i])
            mean = np.mean(contamination)
            std = np.std(contamination)
            contamination_pull[bg].append((mean - nbg[bg])/std)

    bins = np.linspace(-10,10,101)
    bincenters = (bins[1:] + bins[:-1])/2

    fig = plt.figure()
    axes = []
    for i, bg in enumerate(['signal','muon','noise','neck','flasher','breakdown']):
        axes.append(plt.subplot(3,2,i+1))
        plt.hist(contamination_pull[bg],bins=bins,histtype='step',normed=True)
        plt.plot(bincenters,norm.pdf(bincenters))
        plt.title(bg.capitalize())
    for ax in axes:
        ax.set_xlim((-10,10))
        despine(ax=ax,left=True,trim=True)
        ax.get_yaxis().set_visible(False)
    plt.tight_layout()

    if args.save:
        fig.savefig("contamination_pull_plot.pdf")
        fig.savefig("contamination_pull_plot.eps")
    else:
        plt.show()