1 | /* ======================================================================== *\
|
---|
2 | ! $Name: not supported by cvs2svn $:$Id: MHSingleMuon.cc,v 1.17 2007-12-03 17:44:59 tbretz Exp $
|
---|
3 | ! --------------------------------------------------------------------------
|
---|
4 | !
|
---|
5 | ! *
|
---|
6 | ! * This file is part of MARS, the MAGIC Analysis and Reconstruction
|
---|
7 | ! * Software. It is distributed to you in the hope that it can be a useful
|
---|
8 | ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
|
---|
9 | ! * It is distributed WITHOUT ANY WARRANTY.
|
---|
10 | ! *
|
---|
11 | ! * Permission to use, copy, modify and distribute this software and its
|
---|
12 | ! * documentation for any purpose is hereby granted without fee,
|
---|
13 | ! * provided that the above copyright notice appear in all copies and
|
---|
14 | ! * that both that copyright notice and this permission notice appear
|
---|
15 | ! * in supporting documentation. It is provided "as is" without express
|
---|
16 | ! * or implied warranty.
|
---|
17 | ! *
|
---|
18 | !
|
---|
19 | !
|
---|
20 | ! Author(s): Keiichi Mase, 10/2004
|
---|
21 | ! Author(s): Markus Meyer, 02/2005 <mailto:meyer@astro.uni-wuerzburg.de>
|
---|
22 | ! Author(s): Thomas Bretz, 04/2005 <mailto:tbretz@astro.uni-wuerzburg.de>
|
---|
23 | !
|
---|
24 | ! Copyright: MAGIC Software Development, 2000-2005
|
---|
25 | !
|
---|
26 | !
|
---|
27 | \* ======================================================================== */
|
---|
28 |
|
---|
29 | /////////////////////////////////////////////////////////////////////////////
|
---|
30 | //
|
---|
31 | // MHSingleMuon
|
---|
32 | //
|
---|
33 | // This class is a histogram class for displaying the radial (fHistWidth)
|
---|
34 | // and the azimuthal (fHistPhi) intensity distribution for one muon.
|
---|
35 | // You can retrieve the histogram (TH1F) using the function GetHistPhi()
|
---|
36 | // or GetHistWidth().
|
---|
37 | // From these histograms the fraction of the ring segment (ArcPhi) and the
|
---|
38 | // Width of the muon ring (ArcWidth) is calculated.
|
---|
39 | //
|
---|
40 | // First, the radius and center of the ring has to be calculted by
|
---|
41 | // MMuonSearchParCalc
|
---|
42 | // After that the histograms has to be filled in the following way:
|
---|
43 | //
|
---|
44 | // MFillH fillmuon("MHSingleMuon", "", "FillMuon");
|
---|
45 | //
|
---|
46 | // The allowed region to estimate ArcPhi is a certain margin around the
|
---|
47 | // radius. The default value is 0.2 deg (60mm). If the estimated radius
|
---|
48 | // of the arc is 1.0 deg, the pixel contents in the radius range from
|
---|
49 | // 0.8 deg to 1.2 deg are fill in the histogram.
|
---|
50 | //
|
---|
51 | // For ArcPhi only bins over a certain threshold are supposed to be part
|
---|
52 | // of the ring.
|
---|
53 | // For ArcWidth, the same algorithm is used to determine the fit region
|
---|
54 | // for a gaussian fit to the radial intensity distribution. The ArcWidth
|
---|
55 | // is defined as the sigma value of the gaussian fit.
|
---|
56 | //
|
---|
57 | // The binning of the histograms can be changed in the following way:
|
---|
58 | //
|
---|
59 | // MBinning bins1("BinningMuonWidth");
|
---|
60 | // MBinning bins2("BinningArcPhi");
|
---|
61 | // bins1.SetEdges(28, 0.3, 1.7);
|
---|
62 | // bins2.SetEdges(20, -180,180);
|
---|
63 | // plist.AddToList(&bins1);
|
---|
64 | // plist.AddToList(&bins2);
|
---|
65 | //
|
---|
66 | // The values for the thresholds and the margin are saved in MMuonSetup.
|
---|
67 | // They can be easily changed in star.rc.
|
---|
68 | //
|
---|
69 | // Please have in mind, that changes in this basic parameters will change
|
---|
70 | // your results!!
|
---|
71 | //
|
---|
72 | // InputContainer:
|
---|
73 | // - MGeomCam
|
---|
74 | // - MMuonSearchPar
|
---|
75 | //
|
---|
76 | //
|
---|
77 | // Class Version 2:
|
---|
78 | // ----------------
|
---|
79 | // + Double_t fRelTimeMean; // Result of the gaus fit to the arrival time
|
---|
80 | // + Double_t fRelTimeSigma; // Result of the gaus fit to the arrival time
|
---|
81 | //
|
---|
82 | ////////////////////////////////////////////////////////////////////////////
|
---|
83 | #include "MHSingleMuon.h"
|
---|
84 |
|
---|
85 | #include <TF1.h>
|
---|
86 | #include <TMinuit.h>
|
---|
87 | #include <TPad.h>
|
---|
88 | #include <TCanvas.h>
|
---|
89 |
|
---|
90 | #include "MLog.h"
|
---|
91 | #include "MLogManip.h"
|
---|
92 |
|
---|
93 | #include "MBinning.h"
|
---|
94 | #include "MParList.h"
|
---|
95 |
|
---|
96 | #include "MGeomCam.h"
|
---|
97 | #include "MGeomPix.h"
|
---|
98 |
|
---|
99 | #include "MSignalCam.h"
|
---|
100 | #include "MSignalPix.h"
|
---|
101 |
|
---|
102 | #include "MMuonSetup.h"
|
---|
103 | #include "MMuonCalibPar.h"
|
---|
104 | #include "MMuonSearchPar.h"
|
---|
105 |
|
---|
106 | ClassImp(MHSingleMuon);
|
---|
107 |
|
---|
108 | using namespace std;
|
---|
109 |
|
---|
110 | // --------------------------------------------------------------------------
|
---|
111 | //
|
---|
112 | // Setup histograms
|
---|
113 | //
|
---|
114 | MHSingleMuon::MHSingleMuon(const char *name, const char *title) :
|
---|
115 | fSignalCam(0), fMuonSearchPar(0), fGeomCam(0), fMargin(0)
|
---|
116 | {
|
---|
117 | fName = name ? name : "MHSingleMuon";
|
---|
118 | fTitle = title ? title : "Histograms of muon parameters";
|
---|
119 |
|
---|
120 | fHistPhi.SetName("HistPhi");
|
---|
121 | fHistPhi.SetTitle("HistPhi");
|
---|
122 | fHistPhi.SetXTitle("\\phi [\\circ]");
|
---|
123 | fHistPhi.SetYTitle("sum of ADC");
|
---|
124 | fHistPhi.SetDirectory(NULL);
|
---|
125 | fHistPhi.SetFillStyle(4000);
|
---|
126 | fHistPhi.UseCurrentStyle();
|
---|
127 |
|
---|
128 | fHistWidth.SetName("HistWidth");
|
---|
129 | fHistWidth.SetTitle("HistWidth");
|
---|
130 | fHistWidth.SetXTitle("distance from the ring center [\\circ]");
|
---|
131 | fHistWidth.SetYTitle("sum of ADC");
|
---|
132 | fHistWidth.SetDirectory(NULL);
|
---|
133 | fHistWidth.SetFillStyle(4000);
|
---|
134 | fHistWidth.UseCurrentStyle();
|
---|
135 |
|
---|
136 | fHistTime.SetName("HistTime");
|
---|
137 | fHistTime.SetTitle("HistTime");
|
---|
138 | fHistTime.SetXTitle("timing difference");
|
---|
139 | fHistTime.SetYTitle("Counts");
|
---|
140 | fHistTime.SetDirectory(NULL);
|
---|
141 | fHistTime.SetFillStyle(4000);
|
---|
142 | fHistTime.UseCurrentStyle();
|
---|
143 |
|
---|
144 | MBinning bins;
|
---|
145 | bins.SetEdges(20, -180, 180);
|
---|
146 | bins.Apply(fHistPhi);
|
---|
147 |
|
---|
148 | bins.SetEdges(28, 0.3, 1.7);
|
---|
149 | bins.Apply(fHistWidth);
|
---|
150 |
|
---|
151 | bins.SetEdges(101, -33, 33); // +/- 33ns
|
---|
152 | bins.Apply(fHistTime);
|
---|
153 | }
|
---|
154 |
|
---|
155 | // --------------------------------------------------------------------------
|
---|
156 | //
|
---|
157 | // Setup the Binning for the histograms automatically if the correct
|
---|
158 | // instances of MBinning
|
---|
159 | //
|
---|
160 | Bool_t MHSingleMuon::SetupFill(const MParList *plist)
|
---|
161 | {
|
---|
162 | fGeomCam = (MGeomCam*)plist->FindObject("MGeomCam");
|
---|
163 | if (!fGeomCam)
|
---|
164 | {
|
---|
165 | *fLog << warn << "MGeomCam not found... abort." << endl;
|
---|
166 | return kFALSE;
|
---|
167 | }
|
---|
168 | fMuonSearchPar = (MMuonSearchPar*)plist->FindObject("MMuonSearchPar");
|
---|
169 | if (!fMuonSearchPar)
|
---|
170 | {
|
---|
171 | *fLog << warn << "MMuonSearchPar not found... abort." << endl;
|
---|
172 | return kFALSE;
|
---|
173 | }
|
---|
174 | fSignalCam = (MSignalCam*)plist->FindObject("MSignalCam");
|
---|
175 | if (!fSignalCam)
|
---|
176 | {
|
---|
177 | *fLog << warn << "MSignalCam not found... abort." << endl;
|
---|
178 | return kFALSE;
|
---|
179 | }
|
---|
180 |
|
---|
181 | MMuonSetup *setup = (MMuonSetup*)const_cast<MParList*>(plist)->FindCreateObj("MMuonSetup");
|
---|
182 | if (!setup)
|
---|
183 | return kFALSE;
|
---|
184 |
|
---|
185 | fMargin = setup->GetMargin()/fGeomCam->GetConvMm2Deg();
|
---|
186 |
|
---|
187 | ApplyBinning(*plist, "ArcPhi", &fHistPhi);
|
---|
188 | ApplyBinning(*plist, "MuonWidth", &fHistWidth);
|
---|
189 | ApplyBinning(*plist, "MuonTime", &fHistTime);
|
---|
190 |
|
---|
191 | return kTRUE;
|
---|
192 | }
|
---|
193 |
|
---|
194 | // --------------------------------------------------------------------------
|
---|
195 | //
|
---|
196 | // Fill the histograms with data from a MMuonCalibPar and
|
---|
197 | // MMuonSearchPar container.
|
---|
198 | //
|
---|
199 | Bool_t MHSingleMuon::Fill(const MParContainer *par, const Stat_t w)
|
---|
200 | {
|
---|
201 | fRelTimeMean = 0;
|
---|
202 | fRelTimeSigma = -1;
|
---|
203 |
|
---|
204 | fHistPhi.Reset();
|
---|
205 | fHistWidth.Reset();
|
---|
206 | fHistTime.Reset();
|
---|
207 |
|
---|
208 | const Int_t entries = fSignalCam->GetNumPixels();
|
---|
209 |
|
---|
210 | // the position of the center of a muon ring
|
---|
211 | const Float_t cenx = fMuonSearchPar->GetCenterX();
|
---|
212 | const Float_t ceny = fMuonSearchPar->GetCenterY();
|
---|
213 |
|
---|
214 | for (Int_t i=0; i<entries; i++)
|
---|
215 | {
|
---|
216 | const MSignalPix &pix = (*fSignalCam)[i];
|
---|
217 | const MGeomPix &gpix = (*fGeomCam)[i];
|
---|
218 |
|
---|
219 | const Float_t dx = gpix.GetX() - cenx;
|
---|
220 | const Float_t dy = gpix.GetY() - ceny;
|
---|
221 |
|
---|
222 | const Float_t dist = TMath::Hypot(dx, dy);
|
---|
223 |
|
---|
224 | // if the signal is not near the estimated circle, it is ignored.
|
---|
225 | if (TMath::Abs(dist-fMuonSearchPar->GetRadius())<fMargin)
|
---|
226 | {
|
---|
227 | // The arrival time is aligned around 0 for smaller
|
---|
228 | // and more stable histogram range
|
---|
229 | fHistTime.Fill(pix.GetArrivalTime()-fMuonSearchPar->GetTime());
|
---|
230 | }
|
---|
231 |
|
---|
232 | // use only the inner pixles. FIXME: This is geometry dependent
|
---|
233 | if(gpix.GetAidx()>0)
|
---|
234 | continue;
|
---|
235 |
|
---|
236 | fHistWidth.Fill(dist*fGeomCam->GetConvMm2Deg(), pix.GetNumPhotons());
|
---|
237 | }
|
---|
238 | // Setup the function and perform the fit
|
---|
239 | TF1 g1("g1", "gaus");//, -fHistTime.GetXmin(), fHistTime.GetXmax());
|
---|
240 |
|
---|
241 | // Choose starting values as accurate as possible
|
---|
242 | g1.SetParameter(0, fHistTime.GetMaximum());
|
---|
243 | g1.SetParameter(1, 0);
|
---|
244 | g1.SetParameter(2, 0.7); // FIXME! GetRMS instead???
|
---|
245 |
|
---|
246 | // According to fMuonSearchPar->GetTimeRMS() identified muons
|
---|
247 | // do not have an arrival time rms>3
|
---|
248 | g1.SetParLimits(1, -1.7, 1.7);
|
---|
249 | g1.SetParLimits(2, 0, 3.4);
|
---|
250 |
|
---|
251 | // options : N do not store the function, do not draw
|
---|
252 | // I use integral of function in bin rather than value at bin center
|
---|
253 | // R use the range specified in the function range
|
---|
254 | // Q quiet mode
|
---|
255 | fHistTime.Fit(&g1, "QNB");
|
---|
256 |
|
---|
257 | Double_t dummy;
|
---|
258 | gMinuit->GetParameter(1, fRelTimeMean, dummy); // get the mean value
|
---|
259 | gMinuit->GetParameter(2, fRelTimeSigma, dummy); // get the sigma value
|
---|
260 |
|
---|
261 | // The mean arrival time which was subtracted before will
|
---|
262 | // be added again, now
|
---|
263 | const Double_t tm0 = fMuonSearchPar->GetTime()+fRelTimeMean;
|
---|
264 |
|
---|
265 | for (Int_t i=0; i<entries; i++)
|
---|
266 | {
|
---|
267 | const MSignalPix &pix = (*fSignalCam)[i];
|
---|
268 | const MGeomPix &gpix = (*fGeomCam)[i];
|
---|
269 |
|
---|
270 | const Float_t dx = gpix.GetX() - cenx;
|
---|
271 | const Float_t dy = gpix.GetY() - ceny;
|
---|
272 |
|
---|
273 | const Float_t dist = TMath::Hypot(dx, dy);
|
---|
274 |
|
---|
275 | // if the signal is not near the estimated circle, it is ignored.
|
---|
276 | if (TMath::Abs(dist-fMuonSearchPar->GetRadius())<fMargin &&
|
---|
277 | TMath::Abs(pix.GetArrivalTime()-tm0) < 2*fRelTimeSigma)
|
---|
278 | {
|
---|
279 | fHistPhi.Fill(TMath::ATan2(dx, dy)*TMath::RadToDeg(), pix.GetNumPhotons());
|
---|
280 | }
|
---|
281 | }
|
---|
282 |
|
---|
283 | return kTRUE;
|
---|
284 |
|
---|
285 | /*
|
---|
286 | // Because the errors (sqrt(content)) are only scaled by a fixed
|
---|
287 | // factor, and the absolute value of the error is nowhere
|
---|
288 | // needed we skip this step
|
---|
289 |
|
---|
290 | // error estimation (temporarily)
|
---|
291 | // The error is estimated from the signal. In order to do so, we have to
|
---|
292 | // once convert the signal from ADC to photo-electron. Then we can get
|
---|
293 | // the fluctuation such as F-factor*sqrt(phe).
|
---|
294 | // Up to now, the error of pedestal is not taken into accout. This is not
|
---|
295 | // of course correct. We will include this soon.
|
---|
296 | const Double_t Ffactor = 1.4;
|
---|
297 | for (Int_t i=0; i<fHistPhi.GetNbinsX()+1; i++)
|
---|
298 | fHistPhi.SetBinError(i, fHistPhi.GetBinError(i)*Ffactor);
|
---|
299 |
|
---|
300 | for (Int_t i=0; i<fHistWidth.GetNbinsX()+1; i++)
|
---|
301 | fHistWidth.SetBinError(i, fHistWidth.GetBinError(i)*Ffactor);
|
---|
302 |
|
---|
303 | return kTRUE;
|
---|
304 | */
|
---|
305 | }
|
---|
306 |
|
---|
307 | // --------------------------------------------------------------------------
|
---|
308 | //
|
---|
309 | // Find the first bins starting at the bin with maximum content in both
|
---|
310 | // directions which are below threshold.
|
---|
311 | // If in a range of half the histogram size in both directions no bin
|
---|
312 | // below the threshold is found, kFALSE is returned.
|
---|
313 | //
|
---|
314 | Bool_t MHSingleMuon::FindRangeAboveThreshold(const TProfile &h, Float_t thres, Int_t &first, Int_t &last) const
|
---|
315 | {
|
---|
316 | const Int_t n = h.GetNbinsX();
|
---|
317 | const Int_t maxbin = h.GetMaximumBin();
|
---|
318 | const Int_t edge = maxbin+n/2;
|
---|
319 |
|
---|
320 | // Search from the peak to the right
|
---|
321 | last = -1;
|
---|
322 | for (Int_t i=maxbin; i<edge; i++)
|
---|
323 | {
|
---|
324 | const Float_t val = h.GetBinContent(i%n + 1);
|
---|
325 | if (val<thres)
|
---|
326 | {
|
---|
327 | last = i%n+1;
|
---|
328 | break;
|
---|
329 | }
|
---|
330 | }
|
---|
331 |
|
---|
332 | // Search from the peak to the left
|
---|
333 | first = -1;
|
---|
334 | for (Int_t i=maxbin+n-1; i>=edge; i--)
|
---|
335 | {
|
---|
336 | const Float_t val = h.GetBinContent(i%n + 1);
|
---|
337 | if (val<thres)
|
---|
338 | {
|
---|
339 | first = i%n+1;
|
---|
340 | break;
|
---|
341 | }
|
---|
342 | }
|
---|
343 |
|
---|
344 | return first>=0 && last>=0;
|
---|
345 | }
|
---|
346 |
|
---|
347 | // --------------------------------------------------------------------------
|
---|
348 | //
|
---|
349 | // Photon distribution along the estimated circle is fitted with theoritical
|
---|
350 | // function in order to get some more information such as Arc Phi and Arc
|
---|
351 | // Length.
|
---|
352 | //
|
---|
353 | Bool_t MHSingleMuon::CalcPhi(Double_t thres, Double_t &peakphi, Double_t &arcphi) const
|
---|
354 | {
|
---|
355 | if (fHistPhi.GetMaximum()<thres)
|
---|
356 | return kFALSE;
|
---|
357 |
|
---|
358 | peakphi = 180.-fHistPhi.GetBinCenter(fHistPhi.GetMaximumBin());
|
---|
359 |
|
---|
360 | // Now find the position at which the peak edges crosses the threshold
|
---|
361 | Int_t first, last;
|
---|
362 |
|
---|
363 | FindRangeAboveThreshold(fHistPhi, thres, first, last);
|
---|
364 |
|
---|
365 | const Int_t n = fHistPhi.GetNbinsX();
|
---|
366 | const Int_t edge = fHistPhi.GetMaximumBin()+n/2;
|
---|
367 | if (first<0)
|
---|
368 | first = (edge-1)%n+1;
|
---|
369 | if (last<0)
|
---|
370 | last = edge%n+1;;
|
---|
371 |
|
---|
372 | const Float_t startfitval = fHistPhi.GetBinLowEdge(first+1);
|
---|
373 | const Float_t endfitval = fHistPhi.GetBinLowEdge(last);
|
---|
374 |
|
---|
375 | arcphi = last-1<first ? 360+(endfitval-startfitval) : endfitval-startfitval;
|
---|
376 |
|
---|
377 | //if (fEnableImpactCalc)
|
---|
378 | // CalcImpact(effbinnum, startfitval, endfitval);
|
---|
379 |
|
---|
380 | return kTRUE;
|
---|
381 | }
|
---|
382 |
|
---|
383 | // --------------------------------------------------------------------------
|
---|
384 | //
|
---|
385 | // Photon distribution of distance from the center of estimated ring is
|
---|
386 | // fitted in order to get some more information such as ARC WIDTH which
|
---|
387 | // can represent to the PSF of our reflector.
|
---|
388 | //
|
---|
389 | // thres: Threshold above zero to determin the edges of the peak which
|
---|
390 | // is used as fit range
|
---|
391 | // width: ArcWidth returned in deg
|
---|
392 | // chi: Chi^2/NDF of the fit
|
---|
393 | //
|
---|
394 | Bool_t MHSingleMuon::CalcWidth(Double_t thres, Double_t &width, Double_t &chi)
|
---|
395 | {
|
---|
396 | Int_t first, last;
|
---|
397 |
|
---|
398 | if (!FindRangeAboveThreshold(fHistWidth, thres, first, last))
|
---|
399 | return kFALSE;
|
---|
400 |
|
---|
401 | // This happens in some cases
|
---|
402 | const Int_t n = fHistWidth.GetNbinsX()/2;
|
---|
403 | const Int_t m = fHistWidth.GetMaximumBin();
|
---|
404 | if (first>last)
|
---|
405 | if (m>n) // If maximum is on the right side of histogram
|
---|
406 | last = n;
|
---|
407 | else
|
---|
408 | first = 0; // If maximum is on the left side of histogram
|
---|
409 |
|
---|
410 | if (last-first<=3)
|
---|
411 | return kFALSE;
|
---|
412 |
|
---|
413 | // Now get the fit range
|
---|
414 | const Float_t startfitval = fHistWidth.GetBinLowEdge(first+1);
|
---|
415 | const Float_t endfitval = fHistWidth.GetBinLowEdge(last);
|
---|
416 |
|
---|
417 | // Setup the function and perform the fit
|
---|
418 | TF1 f1("f1", "gaus + [3]", startfitval, endfitval);
|
---|
419 | f1.SetLineColor(kBlue);
|
---|
420 |
|
---|
421 | // Choose starting values as accurate as possible
|
---|
422 | f1.SetParameter(0, fHistWidth.GetMaximum());
|
---|
423 | f1.SetParameter(1, fHistWidth.GetBinCenter(m));
|
---|
424 | // f1.SetParameter(2, (endfitval-startfitval)/2);
|
---|
425 | f1.SetParameter(2, 0.1);
|
---|
426 | f1.SetParameter(3, 1.8);
|
---|
427 |
|
---|
428 | // options : N do not store the function, do not draw
|
---|
429 | // I use integral of function in bin rather than value at bin center
|
---|
430 | // R use the range specified in the function range
|
---|
431 | // Q quiet mode
|
---|
432 | // fHistWidth.Fit(&f1, "QRO");
|
---|
433 | fHistWidth.Fit(&f1, "QRN");
|
---|
434 |
|
---|
435 | chi = f1.GetChisquare()/f1.GetNDF();
|
---|
436 |
|
---|
437 | Double_t ferr;
|
---|
438 | gMinuit->GetParameter(2, width, ferr); // get the sigma value
|
---|
439 |
|
---|
440 | return kTRUE;
|
---|
441 | }
|
---|
442 |
|
---|
443 | /*
|
---|
444 | // --------------------------------------------------------------------------
|
---|
445 | //
|
---|
446 | // An impact parameter is calculated by fitting the histogram of photon
|
---|
447 | // distribution along the circle with a theoritical model.
|
---|
448 | // (See G. Vacanti et. al., Astroparticle Physics 2, 1994, 1-11.
|
---|
449 | // The function (6) is used here.)
|
---|
450 | //
|
---|
451 | // By default this calculation is suppressed because this calculation is
|
---|
452 | // very time consuming. If you want to calculate an impact parameter,
|
---|
453 | // you can call the function of EnableImpactCalc().
|
---|
454 | //
|
---|
455 | void MMuonCalibParCalc::CalcImpact(Int_t effbinnum, Float_t startfitval, Float_t endfitval)
|
---|
456 | {
|
---|
457 | // Fit the distribution with Vacanti function. The function is different
|
---|
458 | // for the impact parameter of inside or outside of our reflector.
|
---|
459 | // Then two different functions are applied to the photon distribution,
|
---|
460 | // and the one which give us smaller chisquare value is taken as a
|
---|
461 | // proper one.
|
---|
462 |
|
---|
463 | Double_t val1,err1,val2,err2;
|
---|
464 | // impact parameter inside mirror radius (8.5m)
|
---|
465 | TString func1;
|
---|
466 | Float_t tmpval = (*fMuonSearchPar).GetRadius()*(*fGeomCam).GetConvMm2Deg()*TMath::DegToRad();
|
---|
467 | tmpval = sin(2.*tmpval)*8.5;
|
---|
468 | func1 += "[0]*";
|
---|
469 | func1 += tmpval;
|
---|
470 | func1 += "*(sqrt(1.-([1]/8.5)**2*sin((x-[2])*3.1415926/180.)**2)+([1]/8.5)*cos((x-[2])*3.1415926/180.))";
|
---|
471 |
|
---|
472 | TF1 f1("f1",func1,startfitval,endfitval);
|
---|
473 | f1.SetParameters(2000,3,0);
|
---|
474 | f1.SetParLimits(1,0,8.5);
|
---|
475 | f1.SetParLimits(2,-180.,180.);
|
---|
476 |
|
---|
477 | fMuonCalibPar->fHistPhi->Fit("f1","RQEM");
|
---|
478 |
|
---|
479 | Float_t chi1 = -1;
|
---|
480 | Float_t chi2 = -1;
|
---|
481 | if(effbinnum>3)
|
---|
482 | chi1 = f1.GetChisquare()/((Float_t)(effbinnum-3));
|
---|
483 |
|
---|
484 | gMinuit->GetParameter(1,val1,err1); // get the estimated IP
|
---|
485 | Float_t estip1 = val1;
|
---|
486 |
|
---|
487 | // impact parameter beyond mirror area (8.5m)
|
---|
488 | TString func2;
|
---|
489 | Float_t tmpval2 = (*fMuonSearchPar).GetRadius()*(*fGeomCam).GetConvMm2Deg()*TMath::DegToRad();
|
---|
490 | tmpval2 = sin(2.*tmpval2)*8.5*2.;
|
---|
491 | func2 += "[0]*";
|
---|
492 | func2 += tmpval2;
|
---|
493 | func2 += "*sqrt(1.-(([1]/8.5)*sin((x-[2])*3.1415926/180.))**2)";
|
---|
494 |
|
---|
495 | TF1 f2("f2",func2,startfitval,endfitval);
|
---|
496 | f2.SetParameters(2000,20,0);
|
---|
497 | f2.SetParLimits(1,8.5,300.);
|
---|
498 | f2.SetParLimits(2,-180.,180.);
|
---|
499 |
|
---|
500 | fMuonCalibPar->fHistPhi->Fit("f2","RQEM+");
|
---|
501 |
|
---|
502 | if(effbinnum>3)
|
---|
503 | chi2 = f2.GetChisquare()/((Float_t)(effbinnum-3));
|
---|
504 |
|
---|
505 | gMinuit->GetParameter(1,val2,err2); // get the estimated IP
|
---|
506 | Float_t estip2 = val2;
|
---|
507 |
|
---|
508 | if(effbinnum<=3)
|
---|
509 | {
|
---|
510 | fMuonCalibPar->SetEstImpact(-1.);
|
---|
511 | }
|
---|
512 | if(chi2 > chi1)
|
---|
513 | {
|
---|
514 | fMuonCalibPar->SetEstImpact(estip1);
|
---|
515 | fMuonCalibPar->SetChiArcPhi(chi1);
|
---|
516 | }
|
---|
517 | else
|
---|
518 | {
|
---|
519 | fMuonCalibPar->SetEstImpact(estip2);
|
---|
520 | fMuonCalibPar->SetChiArcPhi(chi2);
|
---|
521 | }
|
---|
522 | }
|
---|
523 | */
|
---|
524 |
|
---|
525 | Float_t MHSingleMuon::CalcSize() const
|
---|
526 | {
|
---|
527 | const Int_t n = fHistPhi.GetNbinsX();
|
---|
528 |
|
---|
529 | Double_t sz=0;
|
---|
530 | for (Int_t i=1; i<=n; i++)
|
---|
531 | sz += fHistPhi.GetBinContent(i)*fHistPhi.GetBinEntries(i);
|
---|
532 |
|
---|
533 | return sz;
|
---|
534 | }
|
---|
535 |
|
---|
536 | void MHSingleMuon::Paint(Option_t *o)
|
---|
537 | {
|
---|
538 | TF1 *f = fHistWidth.GetFunction("f1");
|
---|
539 | if (f)
|
---|
540 | f->ResetBit(1<<9);
|
---|
541 | }
|
---|
542 |
|
---|
543 | void MHSingleMuon::Draw(Option_t *o)
|
---|
544 | {
|
---|
545 | TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
|
---|
546 | pad->SetBorderMode(0);
|
---|
547 |
|
---|
548 | AppendPad("");
|
---|
549 |
|
---|
550 | pad->Divide(1,2);
|
---|
551 |
|
---|
552 | pad->cd(1);
|
---|
553 | gPad->SetBorderMode(0);
|
---|
554 | fHistPhi.Draw();
|
---|
555 |
|
---|
556 | pad->cd(2);
|
---|
557 | gPad->SetBorderMode(0);
|
---|
558 | fHistWidth.Draw();
|
---|
559 | }
|
---|