1 |
|
---|
2 | /* ======================================================================== *\
|
---|
3 | !
|
---|
4 | ! *
|
---|
5 | ! * This file is part of MARS, the MAGIC Analysis and Reconstruction
|
---|
6 | ! * Software. It is distributed to you in the hope that it can be a useful
|
---|
7 | ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
|
---|
8 | ! * It is distributed WITHOUT ANY WARRANTY.
|
---|
9 | ! *
|
---|
10 | ! * Permission to use, copy, modify and distribute this software and its
|
---|
11 | ! * documentation for any purpose is hereby granted without fee,
|
---|
12 | ! * provided that the above copyright notice appear in all copies and
|
---|
13 | ! * that both that copyright notice and this permission notice appear
|
---|
14 | ! * in supporting documentation. It is provided "as is" without express
|
---|
15 | ! * or implied warranty.
|
---|
16 | ! *
|
---|
17 | !
|
---|
18 | !
|
---|
19 | ! Author(s): Thomas Bretz 04/2003 <mailto:tbretz@astro.uni-wuerzburg.de>
|
---|
20 | !
|
---|
21 | ! Copyright: MAGIC Software Development, 2000-2003
|
---|
22 | !
|
---|
23 | !
|
---|
24 | \* ======================================================================== */
|
---|
25 |
|
---|
26 | /////////////////////////////////////////////////////////////////////////////
|
---|
27 | //
|
---|
28 | // MHBlindPixels
|
---|
29 | //
|
---|
30 | ////////////////////////////////////////////////////////////////////////////
|
---|
31 | #include "MHBlindPixels.h"
|
---|
32 |
|
---|
33 | #include <TCanvas.h>
|
---|
34 |
|
---|
35 | #include "MMcEvt.hxx"
|
---|
36 | #include "MBlindPixels.h"
|
---|
37 | #include "MGeomCam.h"
|
---|
38 | #include "MPedPhotCam.h"
|
---|
39 | #include "MParList.h"
|
---|
40 | #include "MBinning.h"
|
---|
41 |
|
---|
42 | #include "MLog.h"
|
---|
43 | #include "MLogManip.h"
|
---|
44 |
|
---|
45 | ClassImp(MHBlindPixels);
|
---|
46 |
|
---|
47 | using namespace std;
|
---|
48 |
|
---|
49 | // -------------------------------------------------------------------------
|
---|
50 | //
|
---|
51 | // Default Constructor.
|
---|
52 | //
|
---|
53 | MHBlindPixels::MHBlindPixels(const char *name, const char *title)
|
---|
54 | {
|
---|
55 | fName = name ? name : "MHBlindPixels";
|
---|
56 | fTitle = title ? title : "Histogram for Blind Pixels vs. Theta";
|
---|
57 |
|
---|
58 | // - we initialize the histogram
|
---|
59 | // - we have to give different names for the different histograms because
|
---|
60 | // root don't allow us to have diferent histograms with the same name
|
---|
61 |
|
---|
62 | fBlindId.SetName("2D-IdBlindPixels");
|
---|
63 | fBlindId.SetTitle("2D-IdBlindPixels");
|
---|
64 | fBlindId.SetDirectory(NULL);
|
---|
65 | fBlindId.SetXTitle("\\Theta [\\circ]");
|
---|
66 | fBlindId.SetYTitle("pixel Id");
|
---|
67 |
|
---|
68 | fBlindN.SetName("2D-NBlindPixels");
|
---|
69 | fBlindN.SetTitle("2D-NBlindPixels");
|
---|
70 | fBlindN.SetDirectory(NULL);
|
---|
71 | fBlindN.SetXTitle("\\Theta [\\circ]");
|
---|
72 | fBlindN.SetYTitle("number of blind pixels");
|
---|
73 | }
|
---|
74 |
|
---|
75 | // --------------------------------------------------------------------------
|
---|
76 | //
|
---|
77 | // Set the binnings and prepare the filling of the histogram
|
---|
78 | //
|
---|
79 | Bool_t MHBlindPixels::SetupFill(const MParList *plist)
|
---|
80 | {
|
---|
81 | MGeomCam *fCam = (MGeomCam*)plist->FindObject(AddSerialNumber("MGeomCam"));
|
---|
82 | if (!fCam)
|
---|
83 | {
|
---|
84 | *fLog << err << "MHBlindPixels::SetupFill; MGeomCam not found... aborting." << endl;
|
---|
85 | return kFALSE;
|
---|
86 | }
|
---|
87 |
|
---|
88 | fMcEvt = (MMcEvt*)plist->FindObject("MMcEvt");
|
---|
89 | if (!fMcEvt)
|
---|
90 | {
|
---|
91 | *fLog << err << "MMcEvt not found... aborting." << endl;
|
---|
92 | return kFALSE;
|
---|
93 | }
|
---|
94 |
|
---|
95 |
|
---|
96 | fPedPhot = (MPedPhotCam*)plist->FindObject("MPedPhotCam");
|
---|
97 | if (!fPedPhot)
|
---|
98 | {
|
---|
99 | *fLog << err << "MPedPhotCam not found... aborting." << endl;
|
---|
100 | return kFALSE;
|
---|
101 | }
|
---|
102 | fPedPhot->InitSize(fCam->GetNumPixels());
|
---|
103 |
|
---|
104 |
|
---|
105 | // Get Theta Binning
|
---|
106 | MBinning* binstheta = (MBinning*)plist->FindObject("BinningTheta", "MBinning");
|
---|
107 | if (!binstheta)
|
---|
108 | {
|
---|
109 | *fLog << err << "Object 'BinningTheta' [MBinning] not found... aborting" << endl;
|
---|
110 | return kFALSE;
|
---|
111 | }
|
---|
112 |
|
---|
113 | // Get binning for pixel number
|
---|
114 | const UInt_t npix1 = fPedPhot->GetSize()+1;
|
---|
115 |
|
---|
116 | MBinning binspix("BinningPixel");
|
---|
117 | binspix.SetEdges(npix1, -0.5, npix1-0.5);
|
---|
118 |
|
---|
119 | // Set binnings in histograms
|
---|
120 | SetBinning(&fBlindId, binstheta, &binspix);
|
---|
121 | SetBinning(&fBlindN, binstheta, &binspix);
|
---|
122 |
|
---|
123 | return kTRUE;
|
---|
124 | }
|
---|
125 |
|
---|
126 | // ------------------------------------------------------------------------
|
---|
127 | //
|
---|
128 | // Drawing function. It creates its own canvas.
|
---|
129 | //
|
---|
130 | void MHBlindPixels::Draw(Option_t *option)
|
---|
131 | {
|
---|
132 | TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
|
---|
133 | pad->SetBorderMode(0);
|
---|
134 |
|
---|
135 | pad->Divide(2,2);
|
---|
136 |
|
---|
137 | TH1D *h;
|
---|
138 |
|
---|
139 | pad->cd(1);
|
---|
140 | fBlindId.Draw(option);
|
---|
141 |
|
---|
142 | pad->cd(2);
|
---|
143 | fBlindN.Draw(option);
|
---|
144 |
|
---|
145 | pad->cd(3);
|
---|
146 | gPad->SetBorderMode(0);
|
---|
147 | h = ((TH2*)&fBlindId)->ProjectionY("ProjY-pixId", -1, 9999, "");
|
---|
148 | h->SetDirectory(NULL);
|
---|
149 | h->SetTitle("Distribution of blind pixel Id");
|
---|
150 | h->SetXTitle("Id of blind pixel");
|
---|
151 | h->SetYTitle("No. of events");
|
---|
152 | h->Draw(option);
|
---|
153 | h->SetBit(kCanDelete);
|
---|
154 |
|
---|
155 | pad->cd(4);
|
---|
156 | gPad->SetBorderMode(0);
|
---|
157 | h = ((TH2*)&fBlindN)->ProjectionY("ProjY-pixN", -1, 9999, "");
|
---|
158 | h->SetDirectory(NULL);
|
---|
159 | h->SetTitle("Distribution of no.of blind pixels");
|
---|
160 | h->SetXTitle("No. of blind pixels");
|
---|
161 | h->SetYTitle("No. of events");
|
---|
162 | h->Draw(option);
|
---|
163 | h->SetBit(kCanDelete);
|
---|
164 |
|
---|
165 | pad->Modified();
|
---|
166 | pad->Update();
|
---|
167 | }
|
---|
168 |
|
---|
169 | Bool_t MHBlindPixels::Fill(const MParContainer *par, const Stat_t w)
|
---|
170 | {
|
---|
171 | if (!par)
|
---|
172 | return kFALSE;
|
---|
173 |
|
---|
174 | Double_t theta = fMcEvt->GetTelescopeTheta()*kRad2Deg;
|
---|
175 | const MBlindPixels &bp = *(MBlindPixels*)par;
|
---|
176 |
|
---|
177 | // FIXME: Slow.
|
---|
178 | const UInt_t npix = fPedPhot->GetSize();
|
---|
179 |
|
---|
180 | UInt_t nb = 0;
|
---|
181 | for (UInt_t i=0; i<npix; i++)
|
---|
182 | {
|
---|
183 | if (bp.IsBlind(i))
|
---|
184 | {
|
---|
185 | fBlindId.Fill(theta, i, w);
|
---|
186 | nb++;
|
---|
187 | }
|
---|
188 | }
|
---|
189 | fBlindN.Fill(theta, nb, w);
|
---|
190 |
|
---|
191 | return kTRUE;
|
---|
192 | }
|
---|
193 |
|
---|
194 |
|
---|
195 |
|
---|
196 |
|
---|
197 |
|
---|
198 |
|
---|
199 |
|
---|
200 |
|
---|
201 |
|
---|
202 |
|
---|
203 |
|
---|
204 |
|
---|