1 | /* ======================================================================== *\
|
---|
2 | !
|
---|
3 | ! *
|
---|
4 | ! * This file is part of MARS, the MAGIC Analysis and Reconstruction
|
---|
5 | ! * Software. It is distributed to you in the hope that it can be a useful
|
---|
6 | ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
|
---|
7 | ! * It is distributed WITHOUT ANY WARRANTY.
|
---|
8 | ! *
|
---|
9 | ! * Permission to use, copy, modify and distribute this software and its
|
---|
10 | ! * documentation for any purpose is hereby granted without fee,
|
---|
11 | ! * provided that the above copyright notice appear in all copies and
|
---|
12 | ! * that both that copyright notice and this permission notice appear
|
---|
13 | ! * in supporting documentation. It is provided "as is" without express
|
---|
14 | ! * or implied warranty.
|
---|
15 | ! *
|
---|
16 | !
|
---|
17 | !
|
---|
18 | ! Author(s): Thomas Bretz, 3/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
|
---|
19 | !
|
---|
20 | ! Copyright: MAGIC Software Development, 2000-2004
|
---|
21 | !
|
---|
22 | !
|
---|
23 | \* ======================================================================== */
|
---|
24 |
|
---|
25 | //////////////////////////////////////////////////////////////////////////////
|
---|
26 | //
|
---|
27 | // MHFalseSource
|
---|
28 | //
|
---|
29 | // Create a false source plot. For the Binning in x,y the object MBinning
|
---|
30 | // "BinningFalseSource" is taken from the paremeter list.
|
---|
31 | //
|
---|
32 | // The binning in alpha is currently fixed as 18bins from 0 to 90deg.
|
---|
33 | //
|
---|
34 | // If MTime, MObservatory and MPointingPos is available in the paremeter
|
---|
35 | // list each image is derotated.
|
---|
36 | //
|
---|
37 | // MHFalseSource fills a 3D histogram with alpha distribution for
|
---|
38 | // each (derotated) x and y position.
|
---|
39 | //
|
---|
40 | // Only a radius of 1.2deg around the camera center is taken into account.
|
---|
41 | //
|
---|
42 | // The displayed histogram is the projection of alpha<fAlphaCut into
|
---|
43 | // the x,y plain and the projection of alpha>90-fAlphaCut
|
---|
44 | //
|
---|
45 | // By using the context menu (find it in a small region between the single
|
---|
46 | // pads) you can change the AlphaCut 'online'
|
---|
47 | //
|
---|
48 | // Each Z-Projection (Alpha-histogram) is scaled such, that the number
|
---|
49 | // of entries fits the maximum number of entries in all Z-Projections.
|
---|
50 | // This should correct for the different acceptance in the center
|
---|
51 | // and at the vorder of the camera. This, however, produces more noise
|
---|
52 | // at the border.
|
---|
53 | //
|
---|
54 | // Here is a slightly simplified version of the algorithm:
|
---|
55 | // ------------------------------------------------------
|
---|
56 | // MHillas hil; // Taken as second argument in MFillH
|
---|
57 | //
|
---|
58 | // MSrcPosCam src;
|
---|
59 | // MHillasSrc hsrc;
|
---|
60 | // hsrc.SetSrcPos(&src);
|
---|
61 | //
|
---|
62 | // for (int ix=0; ix<nx; ix++)
|
---|
63 | // for (int iy=0; iy<ny; iy++)
|
---|
64 | // {
|
---|
65 | // TVector2 v(cx[ix], cy[iy]); //cx center of x-bin ix
|
---|
66 | // if (rho!=0) //cy center of y-bin iy
|
---|
67 | // v=v.Rotate(rho); //image rotation angle
|
---|
68 | //
|
---|
69 | // src.SetXY(v); //source position in the camera
|
---|
70 | //
|
---|
71 | // if (!hsrc.Calc(hil)) //calc source dependant hillas
|
---|
72 | // return;
|
---|
73 | //
|
---|
74 | // //fill absolute alpha into histogram
|
---|
75 | // const Double_t alpha = hsrc.GetAlpha();
|
---|
76 | // fHist.Fill(cx[ix], cy[iy], TMath::Abs(alpha), w);
|
---|
77 | // }
|
---|
78 | // }
|
---|
79 | //
|
---|
80 | // Use MHFalse Source like this:
|
---|
81 | // -----------------------------
|
---|
82 | // MFillH fill("MHFalseSource", "MHillas");
|
---|
83 | // or
|
---|
84 | // MHFalseSource hist;
|
---|
85 | // hist.SetAlphaCut(12.5); // The default value
|
---|
86 | // hist.SetBgmean(55); // The default value
|
---|
87 | // MFillH fill(&hist, "MHillas");
|
---|
88 | //
|
---|
89 | // To be implemented:
|
---|
90 | // ------------------
|
---|
91 | // - a switch to use alpha or |alpha|
|
---|
92 | // - taking the binning for alpha from the parlist (binning is
|
---|
93 | // currently fixed)
|
---|
94 | // - a possibility to change the fit-function (eg using a different TF1)
|
---|
95 | // - a possibility to change the fit algorithm (eg which paremeters
|
---|
96 | // are fixed at which time)
|
---|
97 | // - use a different varaible than alpha
|
---|
98 | // - a possiblity to change the algorithm which is used to calculate
|
---|
99 | // alpha (or another parameter) is missing (this could be done using
|
---|
100 | // a tasklist somewhere...)
|
---|
101 | // - a more clever (and faster) algorithm to fill the histogram, eg by
|
---|
102 | // calculating alpha once and fill the two coils around the mean
|
---|
103 | // - more drawing options...
|
---|
104 | // - Move Significance() to a more 'general' place and implement
|
---|
105 | // also different algorithms like (Li/Ma)
|
---|
106 | // - implement fit for best alpha distribution -- online (Paint)
|
---|
107 | // - currently a constant pointing position is assumed in Fill()
|
---|
108 | // - the center of rotation need not to be the camera center
|
---|
109 | // - ERRORS on each alpha bin to estimate the chi^2 correctly!
|
---|
110 | // (sqrt(N)/binwidth) needed for WOlfgangs proposed caluclation
|
---|
111 | // of alpha(Li/Ma)
|
---|
112 | //
|
---|
113 | //////////////////////////////////////////////////////////////////////////////
|
---|
114 | #include "MHFalseSource.h"
|
---|
115 |
|
---|
116 | #include <TF1.h>
|
---|
117 | #include <TH2.h>
|
---|
118 | #include <TGraph.h>
|
---|
119 | #include <TStyle.h>
|
---|
120 | #include <TCanvas.h>
|
---|
121 | #include <TPaveText.h>
|
---|
122 | #include <TStopwatch.h>
|
---|
123 |
|
---|
124 | #include "MGeomCam.h"
|
---|
125 | #include "MSrcPosCam.h"
|
---|
126 | #include "MHillasSrc.h"
|
---|
127 | #include "MTime.h"
|
---|
128 | #include "MObservatory.h"
|
---|
129 | #include "MPointingPos.h"
|
---|
130 | #include "MAstroCatalog.h"
|
---|
131 | #include "MAstroSky2Local.h"
|
---|
132 | #include "MStatusDisplay.h"
|
---|
133 | #include "MMath.h"
|
---|
134 |
|
---|
135 | #include "MBinning.h"
|
---|
136 | #include "MParList.h"
|
---|
137 |
|
---|
138 | #include "MLog.h"
|
---|
139 | #include "MLogManip.h"
|
---|
140 |
|
---|
141 | ClassImp(MHFalseSource);
|
---|
142 |
|
---|
143 | using namespace std;
|
---|
144 |
|
---|
145 | // --------------------------------------------------------------------------
|
---|
146 | //
|
---|
147 | // Default Constructor
|
---|
148 | //
|
---|
149 | MHFalseSource::MHFalseSource(const char *name, const char *title)
|
---|
150 | : fTime(0), fPointPos(0), fObservatory(0), fMm2Deg(-1), fAlphaCut(12.5),
|
---|
151 | fBgMean(55), fMinDist(-1), fMaxDist(-1), fMinLD(-1), fMaxLD(-1)
|
---|
152 | {
|
---|
153 | //
|
---|
154 | // set the name and title of this object
|
---|
155 | //
|
---|
156 | fName = name ? name : "MHFalseSource";
|
---|
157 | fTitle = title ? title : "3D-plot of Alpha vs x, y";
|
---|
158 |
|
---|
159 | fHist.SetDirectory(NULL);
|
---|
160 |
|
---|
161 | fHist.SetName("Alpha");
|
---|
162 | fHist.SetTitle("3D-plot of Alpha vs x, y");
|
---|
163 | fHist.SetXTitle("x [\\circ]");
|
---|
164 | fHist.SetYTitle("y [\\circ]");
|
---|
165 | fHist.SetZTitle("\\alpha [\\circ]");
|
---|
166 | }
|
---|
167 |
|
---|
168 | void MHFalseSource::MakeSymmetric(TH1 *h)
|
---|
169 | {
|
---|
170 | const Float_t min = TMath::Abs(h->GetMinimum());
|
---|
171 | const Float_t max = TMath::Abs(h->GetMaximum());
|
---|
172 |
|
---|
173 | const Float_t absmax = TMath::Max(min, max)*1.002;
|
---|
174 |
|
---|
175 | h->SetMaximum( absmax);
|
---|
176 | h->SetMinimum(-absmax);
|
---|
177 | }
|
---|
178 |
|
---|
179 | // --------------------------------------------------------------------------
|
---|
180 | //
|
---|
181 | // Set the alpha cut (|alpha|<fAlphaCut) which is use to estimate the
|
---|
182 | // number of excess events
|
---|
183 | //
|
---|
184 | void MHFalseSource::SetAlphaCut(Float_t alpha)
|
---|
185 | {
|
---|
186 | if (alpha<0)
|
---|
187 | *fLog << warn << "Alpha<0... taking absolute value." << endl;
|
---|
188 |
|
---|
189 | fAlphaCut = TMath::Abs(alpha);
|
---|
190 |
|
---|
191 | Modified();
|
---|
192 | }
|
---|
193 |
|
---|
194 | // --------------------------------------------------------------------------
|
---|
195 | //
|
---|
196 | // Set mean alpha around which the off sample is determined
|
---|
197 | // (fBgMean-fAlphaCut/2<|fAlpha|<fBgMean+fAlphaCut/2) which is use
|
---|
198 | // to estimate the number of off events
|
---|
199 | //
|
---|
200 | void MHFalseSource::SetBgMean(Float_t alpha)
|
---|
201 | {
|
---|
202 | if (alpha<0)
|
---|
203 | *fLog << warn << "Alpha<0... taking absolute value." << endl;
|
---|
204 |
|
---|
205 | fBgMean = TMath::Abs(alpha);
|
---|
206 |
|
---|
207 | Modified();
|
---|
208 | }
|
---|
209 |
|
---|
210 | // --------------------------------------------------------------------------
|
---|
211 | //
|
---|
212 | // Calculate Significance as
|
---|
213 | // significance = (s-b)/sqrt(s+k*k*b) mit k=s/b
|
---|
214 | //
|
---|
215 | // s: total number of events in signal region
|
---|
216 | // b: number of background events in signal region
|
---|
217 | //
|
---|
218 | Double_t MHFalseSource::Significance(Double_t s, Double_t b)
|
---|
219 | {
|
---|
220 | return MMath::SignificanceSym(s, b);
|
---|
221 | }
|
---|
222 |
|
---|
223 | // --------------------------------------------------------------------------
|
---|
224 | //
|
---|
225 | // calculates the significance according to Li & Ma
|
---|
226 | // ApJ 272 (1983) 317, Formula 17
|
---|
227 | //
|
---|
228 | // s // s: number of on events
|
---|
229 | // b // b: number of off events
|
---|
230 | // alpha = t_on/t_off; // t: observation time
|
---|
231 | //
|
---|
232 | Double_t MHFalseSource::SignificanceLiMa(Double_t s, Double_t b, Double_t alpha)
|
---|
233 | {
|
---|
234 | return MMath::SignificanceLiMaSigned(s, b);
|
---|
235 | }
|
---|
236 |
|
---|
237 | // --------------------------------------------------------------------------
|
---|
238 | //
|
---|
239 | // Set binnings (takes BinningFalseSource) and prepare filling of the
|
---|
240 | // histogram.
|
---|
241 | //
|
---|
242 | // Also search for MTime, MObservatory and MPointingPos
|
---|
243 | //
|
---|
244 | Bool_t MHFalseSource::SetupFill(const MParList *plist)
|
---|
245 | {
|
---|
246 | const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
|
---|
247 | if (!geom)
|
---|
248 | {
|
---|
249 | *fLog << err << "MGeomCam not found... aborting." << endl;
|
---|
250 | return kFALSE;
|
---|
251 | }
|
---|
252 |
|
---|
253 | fMm2Deg = geom->GetConvMm2Deg();
|
---|
254 |
|
---|
255 | MBinning binsa;
|
---|
256 | binsa.SetEdges(18, 0, 90);
|
---|
257 |
|
---|
258 | const MBinning *bins = (MBinning*)plist->FindObject("BinningFalseSource");
|
---|
259 | if (!bins)
|
---|
260 | {
|
---|
261 | const Float_t r = (geom ? geom->GetMaxRadius()/3 : 200)*fMm2Deg;
|
---|
262 |
|
---|
263 | MBinning b;
|
---|
264 | b.SetEdges(20, -r, r);
|
---|
265 | SetBinning(&fHist, &b, &b, &binsa);
|
---|
266 | }
|
---|
267 | else
|
---|
268 | SetBinning(&fHist, bins, bins, &binsa);
|
---|
269 |
|
---|
270 | fPointPos = (MPointingPos*)plist->FindObject(AddSerialNumber("MPointingPos"));
|
---|
271 | if (!fPointPos)
|
---|
272 | *fLog << warn << "MPointingPos not found... no derotation." << endl;
|
---|
273 |
|
---|
274 | fTime = (MTime*)plist->FindObject(AddSerialNumber("MTime"));
|
---|
275 | if (!fTime)
|
---|
276 | *fLog << warn << "MTime not found... no derotation." << endl;
|
---|
277 |
|
---|
278 | fSrcPos = (MSrcPosCam*)plist->FindObject(AddSerialNumber("MSrcPosCam"));
|
---|
279 | if (!fSrcPos)
|
---|
280 | *fLog << warn << "MSrcPosCam not found... no translation." << endl;
|
---|
281 |
|
---|
282 | fObservatory = (MObservatory*)plist->FindObject(AddSerialNumber("MObservatory"));
|
---|
283 | if (!fObservatory)
|
---|
284 | *fLog << warn << "MObservatory not found... no derotation." << endl;
|
---|
285 |
|
---|
286 | // FIXME: Because the pointing position could change we must check
|
---|
287 | // for the current pointing position and add a offset in the
|
---|
288 | // Fill function!
|
---|
289 | fRa = fPointPos->GetRa();
|
---|
290 | fDec = fPointPos->GetDec();
|
---|
291 |
|
---|
292 | return kTRUE;
|
---|
293 | }
|
---|
294 |
|
---|
295 | // --------------------------------------------------------------------------
|
---|
296 | //
|
---|
297 | // Fill the histogram. For details see the code or the class description
|
---|
298 | //
|
---|
299 | Bool_t MHFalseSource::Fill(const MParContainer *par, const Stat_t w)
|
---|
300 | {
|
---|
301 | const MHillas *hil = dynamic_cast<const MHillas*>(par);
|
---|
302 | if (!hil)
|
---|
303 | {
|
---|
304 | *fLog << err << "MHFalseSource::Fill: No container specified!" << endl;
|
---|
305 | return kFALSE;
|
---|
306 | }
|
---|
307 |
|
---|
308 | // Get max radius...
|
---|
309 | const Double_t maxr = 0.98*TMath::Abs(fHist.GetBinCenter(1));
|
---|
310 |
|
---|
311 | // Get camera rotation angle
|
---|
312 | Double_t rho = 0;
|
---|
313 | if (fTime && fObservatory && fPointPos)
|
---|
314 | rho = fPointPos->RotationAngle(*fObservatory, *fTime);
|
---|
315 | //if (fPointPos)
|
---|
316 | // rho = fPointPos->RotationAngle(*fObservatory);
|
---|
317 |
|
---|
318 | // Create necessary containers for calculation
|
---|
319 | MSrcPosCam src;
|
---|
320 | MHillasSrc hsrc;
|
---|
321 | hsrc.SetSrcPos(&src);
|
---|
322 |
|
---|
323 | // Get number of bins and bin-centers
|
---|
324 | const Int_t nx = fHist.GetNbinsX();
|
---|
325 | const Int_t ny = fHist.GetNbinsY();
|
---|
326 | Axis_t cx[nx];
|
---|
327 | Axis_t cy[ny];
|
---|
328 | fHist.GetXaxis()->GetCenter(cx);
|
---|
329 | fHist.GetYaxis()->GetCenter(cy);
|
---|
330 |
|
---|
331 | for (int ix=0; ix<nx; ix++)
|
---|
332 | {
|
---|
333 | for (int iy=0; iy<ny; iy++)
|
---|
334 | {
|
---|
335 | // check distance... to get a circle plot
|
---|
336 | if (TMath::Hypot(cx[ix], cy[iy])>maxr)
|
---|
337 | continue;
|
---|
338 |
|
---|
339 | // rotate center of bin
|
---|
340 | // precalculation of sin/cos doesn't accelerate
|
---|
341 | TVector2 v(cx[ix], cy[iy]);
|
---|
342 | if (rho!=0)
|
---|
343 | v=v.Rotate(rho);
|
---|
344 |
|
---|
345 | // convert degrees to millimeters
|
---|
346 | v *= 1./fMm2Deg;
|
---|
347 |
|
---|
348 | if (fSrcPos)
|
---|
349 | v += fSrcPos->GetXY();
|
---|
350 |
|
---|
351 | src.SetXY(v);
|
---|
352 |
|
---|
353 | // Source dependant hillas parameters
|
---|
354 | if (hsrc.Calc(*hil)>0)
|
---|
355 | {
|
---|
356 | *fLog << warn << "Calculation of MHillasSrc failed for x=" << cx[ix] << " y=" << cy[iy] << endl;
|
---|
357 | return kFALSE;
|
---|
358 | }
|
---|
359 |
|
---|
360 | // FIXME: This should be replaced by an external MFilter
|
---|
361 | // and/or MTaskList
|
---|
362 | // Source dependant distance cut
|
---|
363 | if (fMinDist>0 && hsrc.GetDist()*fMm2Deg<fMinDist)
|
---|
364 | continue;
|
---|
365 | if (fMaxDist>0 && hsrc.GetDist()*fMm2Deg>fMaxDist)
|
---|
366 | continue;
|
---|
367 |
|
---|
368 | if (fMaxLD>0 && hil->GetLength()>fMaxLD*hsrc.GetDist())
|
---|
369 | continue;
|
---|
370 | if (fMinLD>0 && hil->GetLength()<fMinLD*hsrc.GetDist())
|
---|
371 | continue;
|
---|
372 |
|
---|
373 | // Fill histogram
|
---|
374 | const Double_t alpha = hsrc.GetAlpha();
|
---|
375 | fHist.Fill(cx[ix], cy[iy], TMath::Abs(alpha), w);
|
---|
376 | }
|
---|
377 | }
|
---|
378 |
|
---|
379 | return kTRUE;
|
---|
380 | }
|
---|
381 |
|
---|
382 | // --------------------------------------------------------------------------
|
---|
383 | //
|
---|
384 | // Create projection for off data, taking sum of bin contents of
|
---|
385 | // range (fBgMean-fAlphaCut/2, fBgMean+fAlphaCut) Making sure to take
|
---|
386 | // the same number of bins than for on-data
|
---|
387 | //
|
---|
388 | void MHFalseSource::ProjectOff(TH2D *h2, TH2D *all)
|
---|
389 | {
|
---|
390 | TAxis &axe = *fHist.GetZaxis();
|
---|
391 |
|
---|
392 | // Find range to cut (left edge and width)
|
---|
393 | const Int_t f = axe.FindBin(fBgMean-fAlphaCut/2);
|
---|
394 | const Int_t l = axe.FindBin(fAlphaCut)+f-1;
|
---|
395 |
|
---|
396 | axe.SetRange(f, l);
|
---|
397 | const Float_t cut1 = axe.GetBinLowEdge(f);
|
---|
398 | const Float_t cut2 = axe.GetBinUpEdge(l);
|
---|
399 | h2->SetTitle(Form("Distribution of %.1f\\circ<|\\alpha|<%.1f\\circ in x,y", cut1, cut2));
|
---|
400 |
|
---|
401 | // Get projection for range
|
---|
402 | TH2D *p = (TH2D*)fHist.Project3D("yx_off");
|
---|
403 |
|
---|
404 | // Reset range
|
---|
405 | axe.SetRange(0,9999);
|
---|
406 |
|
---|
407 | // Move contents from projection to h2
|
---|
408 | h2->Reset();
|
---|
409 | h2->Add(p, all->GetMaximum());
|
---|
410 | h2->Divide(all);
|
---|
411 |
|
---|
412 | // Delete p
|
---|
413 | delete p;
|
---|
414 |
|
---|
415 | // Set Minimum as minimum value Greater Than 0
|
---|
416 | h2->SetMinimum(GetMinimumGT(*h2));
|
---|
417 | }
|
---|
418 |
|
---|
419 | // --------------------------------------------------------------------------
|
---|
420 | //
|
---|
421 | // Create projection for on data, taking sum of bin contents of
|
---|
422 | // range (0, fAlphaCut)
|
---|
423 | //
|
---|
424 | void MHFalseSource::ProjectOn(TH2D *h3, TH2D *all)
|
---|
425 | {
|
---|
426 | TAxis &axe = *fHist.GetZaxis();
|
---|
427 |
|
---|
428 | // Find range to cut
|
---|
429 | axe.SetRangeUser(0, fAlphaCut);
|
---|
430 | const Float_t cut = axe.GetBinUpEdge(axe.GetLast());
|
---|
431 | h3->SetTitle(Form("Distribution of |\\alpha|<%.1f\\circ in x,y", cut));
|
---|
432 |
|
---|
433 | // Get projection for range
|
---|
434 | TH2D *p = (TH2D*)fHist.Project3D("yx_on");
|
---|
435 |
|
---|
436 | // Reset range
|
---|
437 | axe.SetRange(0,9999);
|
---|
438 |
|
---|
439 | // Move contents from projection to h3
|
---|
440 | h3->Reset();
|
---|
441 | h3->Add(p, all->GetMaximum());
|
---|
442 | h3->Divide(all);
|
---|
443 |
|
---|
444 | // Delete p
|
---|
445 | delete p;
|
---|
446 |
|
---|
447 | // Set Minimum as minimum value Greater Than 0
|
---|
448 | h3->SetMinimum(GetMinimumGT(*h3));
|
---|
449 | }
|
---|
450 |
|
---|
451 | // --------------------------------------------------------------------------
|
---|
452 | //
|
---|
453 | // Create projection for all data, taking sum of bin contents of
|
---|
454 | // range (0, 90) - corresponding to the number of entries in this slice.
|
---|
455 | //
|
---|
456 | void MHFalseSource::ProjectAll(TH2D *h3)
|
---|
457 | {
|
---|
458 | h3->SetTitle("Number of entries");
|
---|
459 |
|
---|
460 | // Get projection for range
|
---|
461 | TH2D *p = (TH2D*)fHist.Project3D("yx_all");
|
---|
462 |
|
---|
463 | // Move contents from projection to h3
|
---|
464 | h3->Reset();
|
---|
465 | h3->Add(p);
|
---|
466 | delete p;
|
---|
467 |
|
---|
468 | // Set Minimum as minimum value Greater Than 0
|
---|
469 | h3->SetMinimum(GetMinimumGT(*h3));
|
---|
470 | }
|
---|
471 |
|
---|
472 | // --------------------------------------------------------------------------
|
---|
473 | //
|
---|
474 | // Update the projections and paint them
|
---|
475 | //
|
---|
476 | void MHFalseSource::Paint(Option_t *opt)
|
---|
477 | {
|
---|
478 | // Set pretty color palette
|
---|
479 | gStyle->SetPalette(1, 0);
|
---|
480 |
|
---|
481 | TVirtualPad *padsave = gPad;
|
---|
482 |
|
---|
483 | TH1D* h1;
|
---|
484 | TH2D* h0;
|
---|
485 | TH2D* h2;
|
---|
486 | TH2D* h3;
|
---|
487 | TH2D* h4;
|
---|
488 | TH2D* h5;
|
---|
489 |
|
---|
490 | // Update projection of all-events
|
---|
491 | padsave->GetPad(2)->cd(3);
|
---|
492 | if ((h0 = (TH2D*)gPad->FindObject("Alpha_yx_all")))
|
---|
493 | ProjectAll(h0);
|
---|
494 |
|
---|
495 | // Update projection of on-events
|
---|
496 | padsave->GetPad(1)->cd(1);
|
---|
497 | if ((h3 = (TH2D*)gPad->FindObject("Alpha_yx_on")))
|
---|
498 | ProjectOn(h3, h0);
|
---|
499 |
|
---|
500 | // Update projection of off-events
|
---|
501 | padsave->GetPad(1)->cd(3);
|
---|
502 | if ((h2 = (TH2D*)gPad->FindObject("Alpha_yx_off")))
|
---|
503 | ProjectOff(h2, h0);
|
---|
504 |
|
---|
505 | padsave->GetPad(2)->cd(2);
|
---|
506 | if ((h5 = (TH2D*)gPad->FindObject("Alpha_yx_diff")))
|
---|
507 | {
|
---|
508 | h5->Add(h2, h3, -1);
|
---|
509 | MakeSymmetric(h5);
|
---|
510 | }
|
---|
511 |
|
---|
512 | // Update projection of significance
|
---|
513 | padsave->GetPad(1)->cd(2);
|
---|
514 | if (h2 && h3 && (h4 = (TH2D*)gPad->FindObject("Alpha_yx_sig")))
|
---|
515 | {
|
---|
516 | const Int_t nx = h4->GetXaxis()->GetNbins();
|
---|
517 | const Int_t ny = h4->GetYaxis()->GetNbins();
|
---|
518 | //const Int_t nr = nx*nx + ny*ny;
|
---|
519 |
|
---|
520 | Int_t maxx=nx/2;
|
---|
521 | Int_t maxy=ny/2;
|
---|
522 |
|
---|
523 | Int_t max = h4->GetBin(nx, ny);
|
---|
524 |
|
---|
525 | for (int ix=1; ix<=nx; ix++)
|
---|
526 | for (int iy=1; iy<=ny; iy++)
|
---|
527 | {
|
---|
528 | const Int_t n = h4->GetBin(ix, iy);
|
---|
529 |
|
---|
530 | const Double_t s = h3->GetBinContent(n);
|
---|
531 | const Double_t b = h2->GetBinContent(n);
|
---|
532 |
|
---|
533 | const Double_t sig = SignificanceLiMa(s, b);
|
---|
534 |
|
---|
535 | h4->SetBinContent(n, sig);
|
---|
536 |
|
---|
537 | if (sig>h4->GetBinContent(max) && sig>0/* && (ix-nx/2)*(ix-nx/2)+(iy-ny/2)*(iy-ny/2)<nr*nr/9*/)
|
---|
538 | {
|
---|
539 | max = n;
|
---|
540 | maxx=ix;
|
---|
541 | maxy=iy;
|
---|
542 | }
|
---|
543 | }
|
---|
544 |
|
---|
545 | MakeSymmetric(h4);
|
---|
546 |
|
---|
547 | // Update projection of 'the best alpha-plot'
|
---|
548 | padsave->GetPad(2)->cd(1);
|
---|
549 | if ((h1 = (TH1D*)gPad->FindObject("Alpha")) && max>0)
|
---|
550 | {
|
---|
551 | const Double_t x = h4->GetXaxis()->GetBinCenter(maxx);
|
---|
552 | const Double_t y = h4->GetYaxis()->GetBinCenter(maxy);
|
---|
553 | const Double_t s = h4->GetBinContent(max);
|
---|
554 |
|
---|
555 | // This is because a 3D histogram x and y are vice versa
|
---|
556 | // Than for their projections
|
---|
557 | TH1 *h = fHist.ProjectionZ("Alpha", maxx, maxx, maxy, maxy);
|
---|
558 | h->SetTitle(Form("Distribution of \\alpha for x=%.2f y=%.2f (\\sigma_{max}=%.1f)", x, y, s));
|
---|
559 | }
|
---|
560 | }
|
---|
561 |
|
---|
562 | gPad = padsave;
|
---|
563 | }
|
---|
564 |
|
---|
565 | // --------------------------------------------------------------------------
|
---|
566 | //
|
---|
567 | // Get the MAstroCatalog corresponding to fRa, fDec. The limiting magnitude
|
---|
568 | // is set to 9, while the fov is equal to the current fov of the false
|
---|
569 | // source plot.
|
---|
570 | //
|
---|
571 | TObject *MHFalseSource::GetCatalog()
|
---|
572 | {
|
---|
573 | const Double_t maxr = 0.98*TMath::Abs(fHist.GetBinCenter(1));
|
---|
574 |
|
---|
575 | // Create catalog...
|
---|
576 | MAstroCatalog *stars = new MAstroCatalog;
|
---|
577 | stars->SetLimMag(9);
|
---|
578 | stars->SetGuiActive(kFALSE);
|
---|
579 | stars->SetRadiusFOV(maxr);
|
---|
580 | stars->SetRaDec(fRa*TMath::DegToRad()*15, fDec*TMath::DegToRad());
|
---|
581 | stars->ReadBSC("bsc5.dat");
|
---|
582 |
|
---|
583 | *fLog << err << "FIXME - The catalog will never be deleted, because this crashes!" << endl;
|
---|
584 |
|
---|
585 | // stars->SetBit(kCanDelete);
|
---|
586 |
|
---|
587 | return stars;
|
---|
588 | }
|
---|
589 |
|
---|
590 | // --------------------------------------------------------------------------
|
---|
591 | //
|
---|
592 | // Draw the histogram
|
---|
593 | //
|
---|
594 | void MHFalseSource::Draw(Option_t *opt)
|
---|
595 | {
|
---|
596 | TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
|
---|
597 | pad->SetBorderMode(0);
|
---|
598 |
|
---|
599 | AppendPad("");
|
---|
600 |
|
---|
601 | pad->Divide(1, 2, 0, 0.03);
|
---|
602 |
|
---|
603 | *fLog << err << "FIXME - Plotting the catalog is broken!" << endl;
|
---|
604 |
|
---|
605 | TObject *catalog = GetCatalog();
|
---|
606 |
|
---|
607 | // Initialize upper part
|
---|
608 | pad->cd(1);
|
---|
609 | gPad->SetBorderMode(0);
|
---|
610 | gPad->Divide(3, 1);
|
---|
611 |
|
---|
612 | // PAD #1
|
---|
613 | pad->GetPad(1)->cd(1);
|
---|
614 | gPad->SetBorderMode(0);
|
---|
615 | fHist.GetZaxis()->SetRangeUser(0,fAlphaCut);
|
---|
616 | TH1 *h3 = fHist.Project3D("yx_on");
|
---|
617 | fHist.GetZaxis()->SetRange(0,9999);
|
---|
618 | h3->SetDirectory(NULL);
|
---|
619 | h3->SetXTitle(fHist.GetXaxis()->GetTitle());
|
---|
620 | h3->SetYTitle(fHist.GetYaxis()->GetTitle());
|
---|
621 | h3->Draw("colz");
|
---|
622 | h3->SetBit(kCanDelete);
|
---|
623 | catalog->Draw("mirror same *");
|
---|
624 |
|
---|
625 | // PAD #2
|
---|
626 | pad->GetPad(1)->cd(2);
|
---|
627 | gPad->SetBorderMode(0);
|
---|
628 | fHist.GetZaxis()->SetRange(0,0);
|
---|
629 | TH1 *h4 = fHist.Project3D("yx_sig"); // Do this to get the correct binning....
|
---|
630 | fHist.GetZaxis()->SetRange(0,9999);
|
---|
631 | h4->SetTitle("Significance");
|
---|
632 | h4->SetDirectory(NULL);
|
---|
633 | h4->SetXTitle(fHist.GetXaxis()->GetTitle());
|
---|
634 | h4->SetYTitle(fHist.GetYaxis()->GetTitle());
|
---|
635 | h4->Draw("colz");
|
---|
636 | h4->SetBit(kCanDelete);
|
---|
637 | catalog->Draw("mirror same *");
|
---|
638 |
|
---|
639 | // PAD #3
|
---|
640 | pad->GetPad(1)->cd(3);
|
---|
641 | gPad->SetBorderMode(0);
|
---|
642 | fHist.GetZaxis()->SetRangeUser(fBgMean-fAlphaCut/2, fBgMean+fAlphaCut/2);
|
---|
643 | TH1 *h2 = fHist.Project3D("yx_off");
|
---|
644 | h2->SetDirectory(NULL);
|
---|
645 | h2->SetXTitle(fHist.GetXaxis()->GetTitle());
|
---|
646 | h2->SetYTitle(fHist.GetYaxis()->GetTitle());
|
---|
647 | h2->Draw("colz");
|
---|
648 | h2->SetBit(kCanDelete);
|
---|
649 | catalog->Draw("mirror same *");
|
---|
650 |
|
---|
651 | // Initialize lower part
|
---|
652 | pad->cd(2);
|
---|
653 | gPad->SetBorderMode(0);
|
---|
654 | gPad->Divide(3, 1);
|
---|
655 |
|
---|
656 | // PAD #4
|
---|
657 | pad->GetPad(2)->cd(1);
|
---|
658 | gPad->SetBorderMode(0);
|
---|
659 | TH1 *h1 = fHist.ProjectionZ("Alpha");
|
---|
660 | h1->SetDirectory(NULL);
|
---|
661 | h1->SetTitle("Distribution of \\alpha");
|
---|
662 | h1->SetXTitle(fHist.GetZaxis()->GetTitle());
|
---|
663 | h1->SetYTitle("Counts");
|
---|
664 | h1->Draw(opt);
|
---|
665 | h1->SetBit(kCanDelete);
|
---|
666 |
|
---|
667 | // PAD #5
|
---|
668 | pad->GetPad(2)->cd(2);
|
---|
669 | gPad->SetBorderMode(0);
|
---|
670 | TH1 *h5 = (TH1*)h3->Clone("Alpha_yx_diff");
|
---|
671 | h5->Add(h2, -1);
|
---|
672 | h5->SetTitle("Difference of on- and off-distribution");
|
---|
673 | h5->SetDirectory(NULL);
|
---|
674 | h5->Draw("colz");
|
---|
675 | h5->SetBit(kCanDelete);
|
---|
676 | catalog->Draw("mirror same *");
|
---|
677 |
|
---|
678 | // PAD #6
|
---|
679 | pad->GetPad(2)->cd(3);
|
---|
680 | gPad->SetBorderMode(0);
|
---|
681 | TH1 *h0 = fHist.Project3D("yx_all");
|
---|
682 | h0->SetDirectory(NULL);
|
---|
683 | h0->SetXTitle(fHist.GetXaxis()->GetTitle());
|
---|
684 | h0->SetYTitle(fHist.GetYaxis()->GetTitle());
|
---|
685 | h0->Draw("colz");
|
---|
686 | h0->SetBit(kCanDelete);
|
---|
687 | catalog->Draw("mirror same *");
|
---|
688 | }
|
---|
689 |
|
---|
690 | // --------------------------------------------------------------------------
|
---|
691 | //
|
---|
692 | // Everything which is in the main pad belongs to this class!
|
---|
693 | //
|
---|
694 | Int_t MHFalseSource::DistancetoPrimitive(Int_t px, Int_t py)
|
---|
695 | {
|
---|
696 | return 0;
|
---|
697 | }
|
---|
698 |
|
---|
699 | // --------------------------------------------------------------------------
|
---|
700 | //
|
---|
701 | // Set all sub-pads to Modified()
|
---|
702 | //
|
---|
703 | void MHFalseSource::Modified()
|
---|
704 | {
|
---|
705 | if (!gPad)
|
---|
706 | return;
|
---|
707 |
|
---|
708 | TVirtualPad *padsave = gPad;
|
---|
709 | padsave->Modified();
|
---|
710 | padsave->GetPad(1)->cd(1);
|
---|
711 | gPad->Modified();
|
---|
712 | padsave->GetPad(1)->cd(3);
|
---|
713 | gPad->Modified();
|
---|
714 | padsave->GetPad(2)->cd(1);
|
---|
715 | gPad->Modified();
|
---|
716 | padsave->GetPad(2)->cd(2);
|
---|
717 | gPad->Modified();
|
---|
718 | padsave->GetPad(2)->cd(3);
|
---|
719 | gPad->Modified();
|
---|
720 | gPad->cd();
|
---|
721 | }
|
---|
722 |
|
---|
723 | // --------------------------------------------------------------------------
|
---|
724 | //
|
---|
725 | // This is a preliminary implementation of a alpha-fit procedure for
|
---|
726 | // all possible source positions. It will be moved into its own
|
---|
727 | // more powerfull class soon.
|
---|
728 | //
|
---|
729 | // The fit function is "gaus(0)+pol2(3)" which is equivalent to:
|
---|
730 | // [0]*exp(-0.5*((x-[1])/[2])^2) + [3] + [4]*x + [5]*x^2
|
---|
731 | // or
|
---|
732 | // A*exp(-0.5*((x-mu)/sigma)^2) + a + b*x + c*x^2
|
---|
733 | //
|
---|
734 | // Parameter [1] is fixed to 0 while the alpha peak should be
|
---|
735 | // symmetric around alpha=0.
|
---|
736 | //
|
---|
737 | // Parameter [4] is fixed to 0 because the first derivative at
|
---|
738 | // alpha=0 should be 0, too.
|
---|
739 | //
|
---|
740 | // In a first step the background is fitted between bgmin and bgmax,
|
---|
741 | // while the parameters [0]=0 and [2]=1 are fixed.
|
---|
742 | //
|
---|
743 | // In a second step the signal region (alpha<sigmax) is fittet using
|
---|
744 | // the whole function with parameters [1], [3], [4] and [5] fixed.
|
---|
745 | //
|
---|
746 | // The number of excess and background events are calculated as
|
---|
747 | // s = int(0, sigint, gaus(0)+pol2(3))
|
---|
748 | // b = int(0, sigint, pol2(3))
|
---|
749 | //
|
---|
750 | // The Significance is calculated using the Significance() member
|
---|
751 | // function.
|
---|
752 | //
|
---|
753 | void MHFalseSource::FitSignificance(Float_t sigint, Float_t sigmax, Float_t bgmin, Float_t bgmax, Byte_t polynom)
|
---|
754 | {
|
---|
755 | TObject *catalog = GetCatalog();
|
---|
756 |
|
---|
757 | TH1D h0a("A", "", 50, 0, 4000);
|
---|
758 | TH1D h4a("chisq1", "", 50, 0, 35);
|
---|
759 | //TH1D h5a("prob1", "", 50, 0, 1.1);
|
---|
760 | TH1D h6("signifcance", "", 50, -20, 20);
|
---|
761 |
|
---|
762 | TH1D h1("mu", "Parameter \\mu", 50, -1, 1);
|
---|
763 | TH1D h2("sigma", "Parameter \\sigma", 50, 0, 90);
|
---|
764 | TH1D h3("b", "Parameter b", 50, -0.1, 0.1);
|
---|
765 |
|
---|
766 | TH1D h0b("a", "Parameter a (red), A (blue)", 50, 0, 4000);
|
---|
767 | TH1D h4b("\\chi^{2}", "\\chi^{2} (red, green) / significance (black)", 50, 0, 35);
|
---|
768 | //TH1D h5b("prob", "Fit probability: Bg(red), F(blue)", 50, 0, 1.1);
|
---|
769 |
|
---|
770 | h0a.SetLineColor(kBlue);
|
---|
771 | h4a.SetLineColor(kBlue);
|
---|
772 | //h5a.SetLineColor(kBlue);
|
---|
773 | h0b.SetLineColor(kRed);
|
---|
774 | h4b.SetLineColor(kRed);
|
---|
775 | //h5b.SetLineColor(kRed);
|
---|
776 |
|
---|
777 | TH1 *hist = fHist.Project3D("yx_fit");
|
---|
778 | hist->SetDirectory(0);
|
---|
779 | TH1 *hists = fHist.Project3D("yx_fit");
|
---|
780 | hists->SetDirectory(0);
|
---|
781 | TH1 *histb = fHist.Project3D("yx_fit");
|
---|
782 | histb->SetDirectory(0);
|
---|
783 | hist->Reset();
|
---|
784 | hists->Reset();
|
---|
785 | histb->Reset();
|
---|
786 | hist->SetNameTitle("Significance",
|
---|
787 | Form("Fit Region: Signal<%.1f\\circ, %.1f\\circ<Bg<%.1f\\circ",
|
---|
788 | sigmax, bgmin, bgmax));
|
---|
789 | hists->SetName("Excess");
|
---|
790 | histb->SetName("Background");
|
---|
791 | hist->SetXTitle(fHist.GetXaxis()->GetTitle());
|
---|
792 | hists->SetXTitle(fHist.GetXaxis()->GetTitle());
|
---|
793 | histb->SetXTitle(fHist.GetXaxis()->GetTitle());
|
---|
794 | hist->SetYTitle(fHist.GetYaxis()->GetTitle());
|
---|
795 | hists->SetYTitle(fHist.GetYaxis()->GetTitle());
|
---|
796 | histb->SetYTitle(fHist.GetYaxis()->GetTitle());
|
---|
797 |
|
---|
798 | const Double_t w = fHist.GetZaxis()->GetBinWidth(1);
|
---|
799 |
|
---|
800 | // xmin, xmax, npar
|
---|
801 | //TF1 func("MyFunc", fcn, 0, 90, 6);
|
---|
802 | // Implementing the function yourself is only about 5% faster
|
---|
803 | TF1 func("", Form("gaus(0) + pol%d(3)", polynom), 0, 90);
|
---|
804 | //TF1 func("", Form("[0]*(TMath::Gaus(x, [1], [2])+TMath::Gaus(x, -[1], [2]))+pol%d(3)", polynom), 0, 90);
|
---|
805 | TArrayD maxpar(func.GetNpar());
|
---|
806 |
|
---|
807 | /* func.SetParName(0, "A");
|
---|
808 | * func.SetParName(1, "mu");
|
---|
809 | * func.SetParName(2, "sigma");
|
---|
810 | */
|
---|
811 |
|
---|
812 | func.FixParameter(1, 0);
|
---|
813 | func.FixParameter(4, 0);
|
---|
814 | func.SetParLimits(2, 0, 90);
|
---|
815 | func.SetParLimits(3, -1, 1);
|
---|
816 |
|
---|
817 | const Int_t nx = hist->GetXaxis()->GetNbins();
|
---|
818 | const Int_t ny = hist->GetYaxis()->GetNbins();
|
---|
819 | //const Int_t nr = nx*nx+ny*ny;
|
---|
820 |
|
---|
821 | Double_t maxalpha0=0;
|
---|
822 | Double_t maxs=3;
|
---|
823 |
|
---|
824 | Int_t maxx=0;
|
---|
825 | Int_t maxy=0;
|
---|
826 |
|
---|
827 | TStopwatch clk;
|
---|
828 | clk.Start();
|
---|
829 |
|
---|
830 | *fLog << inf;
|
---|
831 | *fLog << "Signal fit: alpha < " << sigmax << endl;
|
---|
832 | *fLog << "Integration: alpha < " << sigint << endl;
|
---|
833 | *fLog << "Background fit: " << bgmin << " < alpha < " << bgmax << endl;
|
---|
834 | *fLog << "Polynom order: " << (int)polynom << endl;
|
---|
835 | *fLog << "Fitting False Source Plot..." << flush;
|
---|
836 |
|
---|
837 | TH1 *h0 = fHist.Project3D("yx_entries");
|
---|
838 | Float_t entries = h0->GetMaximum();
|
---|
839 | delete h0;
|
---|
840 |
|
---|
841 | TH1 *h=0;
|
---|
842 | for (int ix=1; ix<=nx; ix++)
|
---|
843 | for (int iy=1; iy<=ny; iy++)
|
---|
844 | {
|
---|
845 | // This is because a 3D histogram x and y are vice versa
|
---|
846 | // Than for their projections
|
---|
847 | h = fHist.ProjectionZ("AlphaFit", ix, ix, iy, iy);
|
---|
848 |
|
---|
849 | const Double_t alpha0 = h->GetBinContent(1);
|
---|
850 |
|
---|
851 | // Check for the regios which is not filled...
|
---|
852 | if (alpha0==0)
|
---|
853 | continue;
|
---|
854 |
|
---|
855 | h->Scale(entries/h->GetEntries());
|
---|
856 |
|
---|
857 | if (alpha0>maxalpha0)
|
---|
858 | maxalpha0=alpha0;
|
---|
859 |
|
---|
860 | // First fit a polynom in the off region
|
---|
861 | func.FixParameter(0, 0);
|
---|
862 | func.FixParameter(2, 1);
|
---|
863 | func.ReleaseParameter(3);
|
---|
864 |
|
---|
865 | for (int i=5; i<func.GetNpar(); i++)
|
---|
866 | func.ReleaseParameter(i);
|
---|
867 |
|
---|
868 | h->Fit(&func, "N0Q", "", bgmin, bgmax);
|
---|
869 |
|
---|
870 | h4a.Fill(func.GetChisquare());
|
---|
871 | //h5a.Fill(func.GetProb());
|
---|
872 |
|
---|
873 | //func.SetParLimits(0, 0.5*h->GetBinContent(1), 1.5*h->GetBinContent(1));
|
---|
874 | //func.SetParLimits(2, 5, 90);
|
---|
875 |
|
---|
876 | func.ReleaseParameter(0);
|
---|
877 | //func.ReleaseParameter(1);
|
---|
878 | func.ReleaseParameter(2);
|
---|
879 | func.FixParameter(3, func.GetParameter(3));
|
---|
880 | for (int i=5; i<func.GetNpar(); i++)
|
---|
881 | func.FixParameter(i, func.GetParameter(i));
|
---|
882 |
|
---|
883 | // Do not allow signals smaller than the background
|
---|
884 | const Double_t A = alpha0-func.GetParameter(3);
|
---|
885 | const Double_t dA = TMath::Abs(A);
|
---|
886 | func.SetParLimits(0, -dA*4, dA*4);
|
---|
887 | func.SetParLimits(2, 0, 90);
|
---|
888 |
|
---|
889 | // Now fit a gaus in the on region on top of the polynom
|
---|
890 | func.SetParameter(0, A);
|
---|
891 | func.SetParameter(2, sigmax*0.75);
|
---|
892 |
|
---|
893 | h->Fit(&func, "N0Q", "", 0, sigmax);
|
---|
894 |
|
---|
895 | TArrayD p(func.GetNpar(), func.GetParameters());
|
---|
896 |
|
---|
897 | // Fill results into some histograms
|
---|
898 | h0a.Fill(p[0]);
|
---|
899 | h0b.Fill(p[3]);
|
---|
900 | h1.Fill(p[1]);
|
---|
901 | h2.Fill(p[2]);
|
---|
902 | if (polynom>1)
|
---|
903 | h3.Fill(p[5]);
|
---|
904 | h4b.Fill(func.GetChisquare());
|
---|
905 | //h5b.Fill(func.GetProb());
|
---|
906 |
|
---|
907 | // Implementing the integral as analytical function
|
---|
908 | // gives the same result in the order of 10e-5
|
---|
909 | // and it is not faster at all...
|
---|
910 | //const Bool_t ok = /*p[0]>=0 && /*p[0]<alpha0*2 &&*/ p[2]>1.75;// && p[2]<88.5;
|
---|
911 | /*
|
---|
912 | if (p[0]<-fabs(alpha0-p[3])*7 && p[2]<alphaw/3)
|
---|
913 | {
|
---|
914 | func.SetParameter(0, alpha0-p[3]);
|
---|
915 | cout << p[2] << " " << p[0] << " " << alpha0-p[3] << endl;
|
---|
916 | }
|
---|
917 | */
|
---|
918 |
|
---|
919 | // The fitted function returned units of
|
---|
920 | // counts bin binwidth. To get the correct number
|
---|
921 | // of events we must adapt the functions by dividing
|
---|
922 | // the result of the integration by the bin-width
|
---|
923 | const Double_t s = func.Integral(0, sigint)/w;
|
---|
924 |
|
---|
925 | func.SetParameter(0, 0);
|
---|
926 | func.SetParameter(2, 1);
|
---|
927 |
|
---|
928 | const Double_t b = func.Integral(0, sigint)/w;
|
---|
929 | const Double_t sig = SignificanceLiMa(s, b);
|
---|
930 |
|
---|
931 | const Int_t n = hist->GetBin(ix, iy);
|
---|
932 | hists->SetBinContent(n, s-b);
|
---|
933 | histb->SetBinContent(n, b);
|
---|
934 |
|
---|
935 | hist->SetBinContent(n, sig);
|
---|
936 | if (sig!=0)
|
---|
937 | h6.Fill(sig);
|
---|
938 |
|
---|
939 | if (sig>maxs/* && (ix-nx/2)*(ix-nx/2)+(iy-ny/2)*(iy-ny/2)<nr*nr/9*/)
|
---|
940 | {
|
---|
941 | maxs = sig;
|
---|
942 | maxx = ix;
|
---|
943 | maxy = iy;
|
---|
944 | maxpar = p;
|
---|
945 | }
|
---|
946 | }
|
---|
947 |
|
---|
948 | *fLog << inf << "done." << endl;
|
---|
949 |
|
---|
950 | h0a.GetXaxis()->SetRangeUser(0, maxalpha0*1.5);
|
---|
951 | h0b.GetXaxis()->SetRangeUser(0, maxalpha0*1.5);
|
---|
952 |
|
---|
953 | hists->SetTitle(Form("Excess events for \\alpha<%.0f\\circ (N_{max}=%d)", sigint, (int)hists->GetMaximum()));
|
---|
954 | histb->SetTitle(Form("Background events for \\alpha<%.0f\\circ", sigint));
|
---|
955 |
|
---|
956 | //hists->SetMinimum(GetMinimumGT(*hists));
|
---|
957 | histb->SetMinimum(GetMinimumGT(*histb));
|
---|
958 |
|
---|
959 | MakeSymmetric(hists);
|
---|
960 | MakeSymmetric(hist);
|
---|
961 |
|
---|
962 | clk.Stop();
|
---|
963 | clk.Print("m");
|
---|
964 |
|
---|
965 | TCanvas *c=new TCanvas;
|
---|
966 |
|
---|
967 | gStyle->SetPalette(1, 0);
|
---|
968 |
|
---|
969 | c->Divide(3,2, 0, 0);
|
---|
970 | c->cd(1);
|
---|
971 | gPad->SetBorderMode(0);
|
---|
972 | hists->Draw("colz");
|
---|
973 | hists->SetBit(kCanDelete);
|
---|
974 | catalog->Draw("mirror same *");
|
---|
975 | c->cd(2);
|
---|
976 | gPad->SetBorderMode(0);
|
---|
977 | hist->Draw("colz");
|
---|
978 | hist->SetBit(kCanDelete);
|
---|
979 | catalog->Draw("mirror same *");
|
---|
980 | c->cd(3);
|
---|
981 | gPad->SetBorderMode(0);
|
---|
982 | histb->Draw("colz");
|
---|
983 | histb->SetBit(kCanDelete);
|
---|
984 | catalog->Draw("mirror same *");
|
---|
985 | c->cd(4);
|
---|
986 | gPad->Divide(1,3, 0, 0);
|
---|
987 | TVirtualPad *p=gPad;
|
---|
988 | p->SetBorderMode(0);
|
---|
989 | p->cd(1);
|
---|
990 | gPad->SetBorderMode(0);
|
---|
991 | h0b.DrawCopy();
|
---|
992 | h0a.DrawCopy("same");
|
---|
993 | p->cd(2);
|
---|
994 | gPad->SetBorderMode(0);
|
---|
995 | h3.DrawCopy();
|
---|
996 | p->cd(3);
|
---|
997 | gPad->SetBorderMode(0);
|
---|
998 | h2.DrawCopy();
|
---|
999 | c->cd(6);
|
---|
1000 | gPad->Divide(1,2, 0, 0);
|
---|
1001 | TVirtualPad *q=gPad;
|
---|
1002 | q->SetBorderMode(0);
|
---|
1003 | q->cd(1);
|
---|
1004 | gPad->SetBorderMode(0);
|
---|
1005 | gPad->SetBorderMode(0);
|
---|
1006 | h4b.DrawCopy();
|
---|
1007 | h4a.DrawCopy("same");
|
---|
1008 | h6.DrawCopy("same");
|
---|
1009 | q->cd(2);
|
---|
1010 | gPad->SetBorderMode(0);
|
---|
1011 | //h5b.DrawCopy();
|
---|
1012 | //h5a.DrawCopy("same");
|
---|
1013 | c->cd(5);
|
---|
1014 | gPad->SetBorderMode(0);
|
---|
1015 | if (maxx>0 && maxy>0)
|
---|
1016 | {
|
---|
1017 | const char *title = Form(" \\alpha for x=%.2f y=%.2f (\\sigma_{max}=%.1f) ",
|
---|
1018 | hist->GetXaxis()->GetBinCenter(maxx),
|
---|
1019 | hist->GetYaxis()->GetBinCenter(maxy), maxs);
|
---|
1020 |
|
---|
1021 | TH1 *result = fHist.ProjectionZ("AlphaFit", maxx, maxx, maxy, maxy);
|
---|
1022 | result->Scale(entries/h->GetEntries());
|
---|
1023 |
|
---|
1024 | result->SetDirectory(NULL);
|
---|
1025 | result->SetNameTitle("Result \\alpha", title);
|
---|
1026 | result->SetBit(kCanDelete);
|
---|
1027 | result->SetXTitle("\\alpha [\\circ]");
|
---|
1028 | result->SetYTitle("Counts");
|
---|
1029 | result->Draw();
|
---|
1030 |
|
---|
1031 | TF1 f1("", func.GetExpFormula(), 0, 90);
|
---|
1032 | TF1 f2("", func.GetExpFormula(), 0, 90);
|
---|
1033 | f1.SetParameters(maxpar.GetArray());
|
---|
1034 | f2.SetParameters(maxpar.GetArray());
|
---|
1035 | f2.FixParameter(0, 0);
|
---|
1036 | f2.FixParameter(1, 0);
|
---|
1037 | f2.FixParameter(2, 1);
|
---|
1038 | f1.SetLineColor(kGreen);
|
---|
1039 | f2.SetLineColor(kRed);
|
---|
1040 |
|
---|
1041 | f2.DrawCopy("same");
|
---|
1042 | f1.DrawCopy("same");
|
---|
1043 |
|
---|
1044 | TPaveText *leg = new TPaveText(0.35, 0.10, 0.90, 0.35, "brNDC");
|
---|
1045 | leg->SetBorderSize(1);
|
---|
1046 | leg->SetTextSize(0.04);
|
---|
1047 | leg->AddText(0.5, 0.82, Form("A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + pol%d", polynom))->SetTextAlign(22);
|
---|
1048 | //leg->AddText(0.5, 0.82, "A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + b*x^{2} + a")->SetTextAlign(22);
|
---|
1049 | leg->AddLine(0, 0.65, 0, 0.65);
|
---|
1050 | leg->AddText(0.06, 0.54, Form("A=%.2f", maxpar[0]))->SetTextAlign(11);
|
---|
1051 | leg->AddText(0.06, 0.34, Form("\\sigma=%.2f", maxpar[2]))->SetTextAlign(11);
|
---|
1052 | leg->AddText(0.06, 0.14, Form("\\mu=%.2f (fix)", maxpar[1]))->SetTextAlign(11);
|
---|
1053 | leg->AddText(0.60, 0.54, Form("a=%.2f", maxpar[3]))->SetTextAlign(11);
|
---|
1054 | leg->AddText(0.60, 0.34, Form("b=%.2f (fix)", maxpar[4]))->SetTextAlign(11);
|
---|
1055 | if (polynom>1)
|
---|
1056 | leg->AddText(0.60, 0.14, Form("c=%.2f", maxpar[5]))->SetTextAlign(11);
|
---|
1057 | leg->SetBit(kCanDelete);
|
---|
1058 | leg->Draw();
|
---|
1059 |
|
---|
1060 | q->cd(2);
|
---|
1061 |
|
---|
1062 | TGraph *g = new TGraph;
|
---|
1063 | g->SetPoint(0, 0, 0);
|
---|
1064 |
|
---|
1065 | Int_t max=0;
|
---|
1066 | Float_t maxsig=0;
|
---|
1067 | for (int i=1; i<89; i++)
|
---|
1068 | {
|
---|
1069 | const Double_t s = f1.Integral(0, (float)i)/w;
|
---|
1070 | const Double_t b = f2.Integral(0, (float)i)/w;
|
---|
1071 |
|
---|
1072 | const Double_t sig = SignificanceLiMa(s, b);
|
---|
1073 |
|
---|
1074 | g->SetPoint(g->GetN(), i, sig);
|
---|
1075 |
|
---|
1076 | if (sig>maxsig)
|
---|
1077 | {
|
---|
1078 | max = i;
|
---|
1079 | maxsig = sig;
|
---|
1080 | }
|
---|
1081 | }
|
---|
1082 |
|
---|
1083 | g->SetNameTitle("SigVs\\alpha", "Significance vs \\alpha");
|
---|
1084 | g->GetHistogram()->SetXTitle("\\alpha_{0} [\\circ]");
|
---|
1085 | g->GetHistogram()->SetYTitle("Significance");
|
---|
1086 | g->SetBit(kCanDelete);
|
---|
1087 | g->Draw("AP");
|
---|
1088 |
|
---|
1089 | leg = new TPaveText(0.35, 0.10, 0.90, 0.25, "brNDC");
|
---|
1090 | leg->SetBorderSize(1);
|
---|
1091 | leg->SetTextSize(0.1);
|
---|
1092 | leg->AddText(Form("\\sigma_{max}=%.1f at \\alpha_{max}=%d\\circ", maxsig, max));
|
---|
1093 | leg->SetBit(kCanDelete);
|
---|
1094 | leg->Draw();
|
---|
1095 | }
|
---|
1096 | }
|
---|