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, 5/2005 <mailto:tbretz@astro.uni-wuerzburg.de>
|
---|
19 | !
|
---|
20 | ! Copyright: MAGIC Software Development, 2000-2005
|
---|
21 | !
|
---|
22 | !
|
---|
23 | \* ======================================================================== */
|
---|
24 |
|
---|
25 | //////////////////////////////////////////////////////////////////////////////
|
---|
26 | //
|
---|
27 | // MHDisp
|
---|
28 | //
|
---|
29 | // Create a false source plot using disp.
|
---|
30 | //
|
---|
31 | //////////////////////////////////////////////////////////////////////////////
|
---|
32 | #include "MHDisp.h"
|
---|
33 |
|
---|
34 | #include <TStyle.h>
|
---|
35 | #include <TCanvas.h>
|
---|
36 |
|
---|
37 | #include <TF1.h>
|
---|
38 | #include <TF2.h>
|
---|
39 | #include <TProfile.h>
|
---|
40 |
|
---|
41 | #include "MParList.h"
|
---|
42 | #include "MParameters.h"
|
---|
43 |
|
---|
44 | #include "MHillasExt.h"
|
---|
45 | #include "MHillasSrc.h"
|
---|
46 | #include "MSrcPosCam.h"
|
---|
47 | #include "MPointingPos.h"
|
---|
48 |
|
---|
49 | #include "MLog.h"
|
---|
50 | #include "MLogManip.h"
|
---|
51 |
|
---|
52 | ClassImp(MHDisp);
|
---|
53 |
|
---|
54 | using namespace std;
|
---|
55 |
|
---|
56 | // --------------------------------------------------------------------------
|
---|
57 | //
|
---|
58 | // Default Constructor
|
---|
59 | //
|
---|
60 | MHDisp::MHDisp(const char *name, const char *title)
|
---|
61 | : fHilExt(0), fDisp(0)
|
---|
62 | {
|
---|
63 | //
|
---|
64 | // set the name and title of this object
|
---|
65 | //
|
---|
66 | fName = name ? name : "MHDisp";
|
---|
67 | fTitle = title ? title : "3D-plot using Disp vs x, y";
|
---|
68 |
|
---|
69 | fHist.SetDirectory(NULL);
|
---|
70 |
|
---|
71 | fHist.SetName("Alpha");
|
---|
72 | fHist.SetTitle("3D-plot of ThetaSq vs x, y");
|
---|
73 | fHist.SetXTitle("x [\\circ]");
|
---|
74 | fHist.SetYTitle("y [\\circ]");
|
---|
75 | fHist.SetZTitle("\\vartheta [deg^{2}]");
|
---|
76 | }
|
---|
77 |
|
---|
78 | // --------------------------------------------------------------------------
|
---|
79 | //
|
---|
80 | // Set binnings (takes BinningFalseSource) and prepare filling of the
|
---|
81 | // histogram.
|
---|
82 | //
|
---|
83 | // Also search for MTime, MObservatory and MPointingPos
|
---|
84 | //
|
---|
85 | //MHillasExt *ext=0;
|
---|
86 | Bool_t MHDisp::SetupFill(const MParList *plist)
|
---|
87 | {
|
---|
88 | if (!MHFalseSource::SetupFill(plist))
|
---|
89 | return kFALSE;
|
---|
90 |
|
---|
91 | fDisp = (MParameterD*)plist->FindObject("Disp", "MParameterD");
|
---|
92 | if (!fDisp)
|
---|
93 | {
|
---|
94 | *fLog << err << "Disp [MParameterD] not found... abort." << endl;
|
---|
95 | return kFALSE;
|
---|
96 | }
|
---|
97 |
|
---|
98 | MParameterD *p = (MParameterD*)plist->FindObject("M3lCut", "MParameterD");
|
---|
99 | if (!p)
|
---|
100 | {
|
---|
101 | *fLog << err << "M3lCut [MParameterD] not found... abort." << endl;
|
---|
102 | return kFALSE;
|
---|
103 | }
|
---|
104 | fM3lCut = TMath::Abs(p->GetVal());
|
---|
105 |
|
---|
106 | p = (MParameterD*)plist->FindObject("DispXi", "MParameterD");
|
---|
107 | if (!p)
|
---|
108 | {
|
---|
109 | *fLog << err << "DispXi [MParameterD] not found... abort." << endl;
|
---|
110 | return kFALSE;
|
---|
111 | }
|
---|
112 | fXi = p->GetVal();
|
---|
113 |
|
---|
114 | p = (MParameterD*)plist->FindObject("DispXiTheta", "MParameterD");
|
---|
115 | if (!p)
|
---|
116 | {
|
---|
117 | *fLog << err << "DispXiTheta [MParameterD] not found... abort." << endl;
|
---|
118 | return kFALSE;
|
---|
119 | }
|
---|
120 | fXiTheta = p->GetVal();
|
---|
121 |
|
---|
122 | fHilExt = (MHillasExt*)plist->FindObject("MHillasExt");
|
---|
123 | if (!fHilExt)
|
---|
124 | {
|
---|
125 | *fLog << err << "MHillasExt not found... aborting." << endl;
|
---|
126 | return kFALSE;
|
---|
127 | }
|
---|
128 |
|
---|
129 | // Initialize all bins with a small (=0) value otherwise
|
---|
130 | // these bins are not displayed
|
---|
131 | for (int x=0; x<fHist.GetNbinsX(); x++)
|
---|
132 | for (int y=0; y<fHist.GetNbinsX(); y++)
|
---|
133 | fHist.Fill(fHist.GetXaxis()->GetBinCenter(x+1),
|
---|
134 | fHist.GetYaxis()->GetBinCenter(y+1),
|
---|
135 | 0.0, 1e-10);
|
---|
136 |
|
---|
137 | return kTRUE;
|
---|
138 | }
|
---|
139 |
|
---|
140 | // --------------------------------------------------------------------------
|
---|
141 | //
|
---|
142 | // Fill the histogram. For details see the code or the class description
|
---|
143 | //
|
---|
144 | Bool_t MHDisp::Fill(const MParContainer *par, const Stat_t w)
|
---|
145 | {
|
---|
146 | const MHillas *hil = dynamic_cast<const MHillas*>(par);
|
---|
147 | if (!hil)
|
---|
148 | {
|
---|
149 | *fLog << err << "MHDisp::Fill: No container specified!" << endl;
|
---|
150 | return kFALSE;
|
---|
151 | }
|
---|
152 |
|
---|
153 | // Get camera rotation angle
|
---|
154 | Double_t rho = 0;
|
---|
155 | if (fTime && fObservatory && fPointPos)
|
---|
156 | rho = fPointPos->RotationAngle(*fObservatory, *fTime);
|
---|
157 |
|
---|
158 | // Get Disp from Parlist
|
---|
159 | const Double_t disp = fDisp->GetVal();
|
---|
160 |
|
---|
161 | // Calculate where disp is pointing
|
---|
162 | TVector2 pos1(hil->GetCosDelta(), hil->GetSinDelta());
|
---|
163 | TVector2 pos2(hil->GetCosDelta(), hil->GetSinDelta());
|
---|
164 | pos1 *= disp;
|
---|
165 | pos2 *= -disp;
|
---|
166 |
|
---|
167 | pos1 += hil->GetMean()*fMm2Deg;
|
---|
168 | pos2 += hil->GetMean()*fMm2Deg;
|
---|
169 |
|
---|
170 | // gammaweight: If we couldn't decide which position makes the
|
---|
171 | // event a gamma, both position are assigned 'half of a gamma'
|
---|
172 | Double_t gweight = 0.5;
|
---|
173 |
|
---|
174 | // Check whether our g/h-separation allows to asign definitly
|
---|
175 | // to one unique position. Therefor we requeire that the event
|
---|
176 | // would be considered a gamma for one, but not for the other
|
---|
177 | // position. This can only be the case if the third moment
|
---|
178 | // has a value higher than the absolute cut value.
|
---|
179 | if (TMath::Abs(fHilExt->GetM3Long()) > fM3lCut)
|
---|
180 | {
|
---|
181 | // Because at one position the event is considered a gamma
|
---|
182 | // we have to find out which position it is...
|
---|
183 | MSrcPosCam src;
|
---|
184 | MHillasSrc hsrc;
|
---|
185 | hsrc.SetSrcPos(&src);
|
---|
186 |
|
---|
187 | // Calculate the sign for one of the desired source positions
|
---|
188 | // The other position must have the opposite sign
|
---|
189 | src.SetXY(pos1/fMm2Deg);
|
---|
190 | if (hsrc.Calc(*hil)>0)
|
---|
191 | {
|
---|
192 | *fLog << warn << "Calculation of MHillasSrc failed." << endl;
|
---|
193 | return kFALSE;
|
---|
194 | }
|
---|
195 | const Double_t m3l = fHilExt->GetM3Long()*TMath::Sign(1.0f, hsrc.GetCosDeltaAlpha())*fMm2Deg;
|
---|
196 |
|
---|
197 | gweight = m3l>fM3lCut ? 1 : 0;
|
---|
198 | }
|
---|
199 |
|
---|
200 | // Now we can safly derotate both position...
|
---|
201 | //TVector2 srcp(fSrcPos->GetXY());
|
---|
202 | if (rho!=0)
|
---|
203 | {
|
---|
204 | pos1=pos1.Rotate(-rho);
|
---|
205 | pos2=pos2.Rotate(-rho);
|
---|
206 | //srcp=srcp.Rotate(-rho);
|
---|
207 | }
|
---|
208 |
|
---|
209 | /*if (fSrcPos)
|
---|
210 | {
|
---|
211 | pos1 -= srcp*fMm2Deg;
|
---|
212 | pos2 -= srcp*fMm2Deg;
|
---|
213 | }*/
|
---|
214 |
|
---|
215 | // Workaround: Number-of-entries
|
---|
216 | if (gweight>0.25)
|
---|
217 | fHist.Fill(pos1.X(), pos1.Y(), 0.0, w*gweight);
|
---|
218 | if (gweight<0.75)
|
---|
219 | fHist.Fill(pos2.X(), pos2.Y(), 0.0, w*(1-gweight));
|
---|
220 |
|
---|
221 | return kTRUE;
|
---|
222 | }
|
---|
223 | /*
|
---|
224 | static Double_t FcnGauss2d(Double_t *x, Double_t *par)
|
---|
225 | {
|
---|
226 | TVector2 v = TVector2(x[0], x[1]).Rotate(par[5]*TMath::DegToRad());
|
---|
227 |
|
---|
228 | const Double_t g0 = TMath::Gaus(v.X(), par[1], par[2]);
|
---|
229 | const Double_t g1 = TMath::Gaus(v.Y(), par[3], par[4]);
|
---|
230 |
|
---|
231 | //Gaus(Double_t x, Double_t mean=0, Double_t sigma=1, Bool_t norm=kFALSE);
|
---|
232 | return par[0]*g0*g1 + par[6]*(v.X()*v.X() + v.Y()*v.Y()) + par[7];
|
---|
233 | }*/
|
---|
234 |
|
---|
235 | void MHDisp::Paint(Option_t *o)
|
---|
236 | {
|
---|
237 | TVirtualPad *pad = gPad;
|
---|
238 |
|
---|
239 | pad->cd(1);
|
---|
240 |
|
---|
241 | fHist.GetZaxis()->SetRange(0,9999);
|
---|
242 | TH1 *h1=fHist.Project3D("yx_on");
|
---|
243 | gStyle->SetPalette(1, 0);
|
---|
244 |
|
---|
245 | Int_t ix, iy, iz;
|
---|
246 | TF1 func("fcn", "gaus + [3]*x*x + [4]");
|
---|
247 |
|
---|
248 | pad->cd(3);
|
---|
249 | TH1 *h2 = (TH1*)gPad->FindObject("Radial");
|
---|
250 | if (h1 && h2)
|
---|
251 | {
|
---|
252 | h2->Reset();
|
---|
253 |
|
---|
254 | h1->GetMaximumBin(ix, iy, iz);
|
---|
255 |
|
---|
256 | const Double_t x0 = h1->GetXaxis()->GetBinCenter(ix);
|
---|
257 | const Double_t y0 = h1->GetYaxis()->GetBinCenter(iy);
|
---|
258 |
|
---|
259 | h2->SetTitle(Form("Profile @ (%.2f\\circ, %.2f\\circ)", x0, y0));
|
---|
260 |
|
---|
261 | for (int x=0; x<h1->GetNbinsX(); x++)
|
---|
262 | for (int y=0; y<h1->GetNbinsY(); y++)
|
---|
263 | {
|
---|
264 | const Double_t r = TMath::Hypot(h1->GetXaxis()->GetBinCenter(x+1)-x0,
|
---|
265 | h1->GetXaxis()->GetBinCenter(y+1)-y0);
|
---|
266 | h2->Fill(r, h1->GetBinContent(x+1, y+1));
|
---|
267 | }
|
---|
268 |
|
---|
269 | func.SetLineWidth(1);
|
---|
270 | func.SetLineColor(kBlue);
|
---|
271 | func.SetParameter(0, h2->GetBinContent(1));
|
---|
272 | func.FixParameter(1, 0);
|
---|
273 | func.SetParameter(2, 0.05);
|
---|
274 | func.SetParameter(4, h2->GetBinContent(10));
|
---|
275 | h2->Fit(&func, "IMQ");
|
---|
276 | }
|
---|
277 | /*
|
---|
278 | if (h1)
|
---|
279 | {
|
---|
280 | const Double_t maxr = 0.9*TMath::Abs(fHist.GetXaxis()->GetXmax());
|
---|
281 |
|
---|
282 | TF2 f2d("Gaus2D", FcnGauss2d, -maxr, maxr, -maxr, maxr, 8);
|
---|
283 | f2d.SetLineWidth(1);
|
---|
284 |
|
---|
285 | f2d.SetParameter(0, h1->GetMaximum()*5); // A
|
---|
286 | f2d.SetParLimits(1, h1->GetXaxis()->GetBinCenter(ix)-h1->GetXaxis()->GetBinWidth(ix)*5, h1->GetXaxis()->GetBinCenter(ix)+h1->GetXaxis()->GetBinWidth(ix)); // mu_1
|
---|
287 | f2d.SetParLimits(3, h1->GetYaxis()->GetBinCenter(iy)-h1->GetYaxis()->GetBinWidth(iy)*5, h1->GetYaxis()->GetBinCenter(iy)+h1->GetYaxis()->GetBinWidth(iy)); // mu_2
|
---|
288 | f2d.SetParLimits(2, 0, func.GetParameter(2)*5); // sigma_1
|
---|
289 | f2d.SetParLimits(4, 0, func.GetParameter(2)*5); // sigma_2
|
---|
290 | f2d.SetParLimits(5, 0, 45); // phi
|
---|
291 | f2d.SetParLimits(6, 0, func.GetParameter(3)*5);
|
---|
292 | f2d.SetParLimits(7, 0, func.GetParameter(4)*5);
|
---|
293 |
|
---|
294 | f2d.SetParameter(0, h1->GetMaximum()); // A
|
---|
295 | f2d.SetParameter(1, h1->GetXaxis()->GetBinCenter(ix)); // mu_1
|
---|
296 | f2d.SetParameter(2, func.GetParameter(2)); // sigma_1
|
---|
297 | f2d.SetParameter(3, h1->GetYaxis()->GetBinCenter(iy)); // mu_2
|
---|
298 | f2d.SetParameter(4, func.GetParameter(2)); // sigma_2
|
---|
299 | f2d.FixParameter(5, 0); // phi
|
---|
300 | f2d.SetParameter(6, func.GetParameter(3));
|
---|
301 | f2d.SetParameter(7, func.GetParameter(4));
|
---|
302 | h1->Fit(&f2d, "Q", "cont2");
|
---|
303 | //f2d.DrawCopy("cont2same");
|
---|
304 | }*/
|
---|
305 | }
|
---|
306 |
|
---|
307 | void MHDisp::Draw(Option_t *o)
|
---|
308 | {
|
---|
309 | TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
|
---|
310 | const Int_t col = pad->GetFillColor();
|
---|
311 | pad->SetBorderMode(0);
|
---|
312 |
|
---|
313 | AppendPad("");
|
---|
314 |
|
---|
315 | TString name = Form("%s_1", pad->GetName());
|
---|
316 | TPad *p = new TPad(name,name, 0.005, 0.005, 0.65, 0.995, col, 0, 0);
|
---|
317 | p->SetNumber(1);
|
---|
318 | p->Draw();
|
---|
319 | p->cd();
|
---|
320 |
|
---|
321 | TH1 *h3 = fHist.Project3D("yx_on");
|
---|
322 | h3->SetTitle("Distribution of equivalent events");
|
---|
323 | h3->SetDirectory(NULL);
|
---|
324 | h3->SetXTitle(fHist.GetXaxis()->GetTitle());
|
---|
325 | h3->SetYTitle(fHist.GetYaxis()->GetTitle());
|
---|
326 | h3->SetMinimum(0);
|
---|
327 | h3->Draw("colz");
|
---|
328 | h3->SetBit(kCanDelete);
|
---|
329 | // catalog->Draw("mirror same *");
|
---|
330 |
|
---|
331 | pad->cd();
|
---|
332 | name = Form("%s_2", pad->GetName());
|
---|
333 | p = new TPad(name,name, 0.66, 0.005, 0.995, 0.5, col, 0, 0);
|
---|
334 | p->SetNumber(2);
|
---|
335 | p->Draw();
|
---|
336 | p->cd();
|
---|
337 | h3->Draw("surf3");
|
---|
338 |
|
---|
339 | pad->cd();
|
---|
340 | name = Form("%s_3", pad->GetName());
|
---|
341 | p = new TPad(name,name, 0.66, 0.5, 0.995, 0.995, col, 0, 0);
|
---|
342 | p->SetNumber(3);
|
---|
343 | p->Draw();
|
---|
344 | p->cd();
|
---|
345 |
|
---|
346 | const Double_t maxr = TMath::Hypot(h3->GetXaxis()->GetXmax(), h3->GetYaxis()->GetXmax());
|
---|
347 | const Int_t nbin = (h3->GetNbinsX()+h3->GetNbinsY())/2;
|
---|
348 | TProfile *h = new TProfile("Radial", "Radial Profile", nbin, 0, maxr, "s");
|
---|
349 | h->SetBit(kCanDelete);
|
---|
350 | h->Draw();
|
---|
351 | }
|
---|