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): Wolfgang Wittek, 04/2003 <mailto:wittek@mppmu.mpg.de>
|
---|
19 | !
|
---|
20 | ! Copyright: MAGIC Software Development, 2000-2003
|
---|
21 | !
|
---|
22 | !
|
---|
23 | \* ======================================================================== */
|
---|
24 |
|
---|
25 | /////////////////////////////////////////////////////////////////////////////
|
---|
26 | // //
|
---|
27 | // MCT1SupercutsCalc //
|
---|
28 | // //
|
---|
29 | /////////////////////////////////////////////////////////////////////////////
|
---|
30 | #include "MCT1SupercutsCalc.h"
|
---|
31 |
|
---|
32 | #include <math.h>
|
---|
33 | #include <fstream>
|
---|
34 |
|
---|
35 | #include "TFile.h"
|
---|
36 | #include "TArrayD.h"
|
---|
37 |
|
---|
38 | #include "MParList.h"
|
---|
39 | #include "MHillasExt.h"
|
---|
40 | #include "MHillasSrc.h"
|
---|
41 | #include "MMcEvt.hxx"
|
---|
42 | #include "MCerPhotEvt.h"
|
---|
43 | #include "MGeomCam.h"
|
---|
44 | #include "MHadronness.h"
|
---|
45 | #include "MHMatrix.h"
|
---|
46 | #include "MCT1Supercuts.h"
|
---|
47 |
|
---|
48 | #include "MLog.h"
|
---|
49 | #include "MLogManip.h"
|
---|
50 |
|
---|
51 | ClassImp(MCT1SupercutsCalc);
|
---|
52 |
|
---|
53 | using namespace std;
|
---|
54 |
|
---|
55 |
|
---|
56 | // --------------------------------------------------------------------------
|
---|
57 | //
|
---|
58 | // constructor
|
---|
59 | //
|
---|
60 |
|
---|
61 | MCT1SupercutsCalc::MCT1SupercutsCalc(const char *hilname,
|
---|
62 | const char *hilsrcname, const char *name, const char *title)
|
---|
63 | : fHadronnessName("MHadronness"), fHilName(hilname), fHilSrcName(hilsrcname),
|
---|
64 | fSuperName("MCT1Supercuts")
|
---|
65 | {
|
---|
66 | fName = name ? name : "MCT1SupercutsCalc";
|
---|
67 | fTitle = title ? title : "Class to evaluate the Supercuts";
|
---|
68 |
|
---|
69 | fMatrix = NULL;
|
---|
70 | }
|
---|
71 |
|
---|
72 | // --------------------------------------------------------------------------
|
---|
73 | //
|
---|
74 | Int_t MCT1SupercutsCalc::PreProcess(MParList *pList)
|
---|
75 | {
|
---|
76 | MGeomCam *cam = (MGeomCam*)pList->FindObject("MGeomCam");
|
---|
77 | if (!cam)
|
---|
78 | {
|
---|
79 | *fLog << err << "MGeomCam (Camera Geometry) not found... aborting." << endl;
|
---|
80 | return kFALSE;
|
---|
81 | }
|
---|
82 |
|
---|
83 | fMm2Deg = cam->GetConvMm2Deg();
|
---|
84 |
|
---|
85 | fHadronness = (MHadronness*)pList->FindCreateObj("MHadronness", fHadronnessName);
|
---|
86 | if (!fHadronness)
|
---|
87 | {
|
---|
88 | *fLog << err << fHadronnessName << " [MHadronness] not found... aborting." << endl;
|
---|
89 | return kFALSE;
|
---|
90 | }
|
---|
91 |
|
---|
92 | fSuper = (MCT1Supercuts*)pList->FindObject(fSuperName, "MCT1Supercuts");
|
---|
93 | if (!fSuper)
|
---|
94 | {
|
---|
95 | *fLog << err << fSuperName << " [MCT1Supercuts] not found... aborting." << endl;
|
---|
96 | return kFALSE;
|
---|
97 | }
|
---|
98 |
|
---|
99 |
|
---|
100 | if (fMatrix)
|
---|
101 | return kTRUE;
|
---|
102 |
|
---|
103 | //-----------------------------------------------------------
|
---|
104 | fHil = (MHillas*)pList->FindObject(fHilName, "MHillas");
|
---|
105 | if (!fHil)
|
---|
106 | {
|
---|
107 | *fLog << err << fHilName << " [MHillas] not found... aborting." << endl;
|
---|
108 | return kFALSE;
|
---|
109 | }
|
---|
110 |
|
---|
111 | fHilSrc = (MHillasSrc*)pList->FindObject(fHilSrcName, "MHillasSrc");
|
---|
112 | if (!fHilSrc)
|
---|
113 | {
|
---|
114 | *fLog << err << fHilSrcName << " [MHillasSrc] not found... aborting." << endl;
|
---|
115 | return kFALSE;
|
---|
116 | }
|
---|
117 |
|
---|
118 | fMcEvt = (MMcEvt*)pList->FindObject("MMcEvt");
|
---|
119 | if (!fMcEvt)
|
---|
120 | {
|
---|
121 | *fLog << err << "MMcEvt not found... aborting." << endl;
|
---|
122 | return kFALSE;
|
---|
123 | }
|
---|
124 |
|
---|
125 |
|
---|
126 | return kTRUE;
|
---|
127 | }
|
---|
128 |
|
---|
129 | // --------------------------------------------------------------------------
|
---|
130 | //
|
---|
131 | // Calculation of upper and lower limits
|
---|
132 | //
|
---|
133 | Double_t MCT1SupercutsCalc::CtsMCut(TArrayD &a, Double_t ls, Double_t ct,
|
---|
134 | Double_t ls2, Double_t dd2)
|
---|
135 | {
|
---|
136 | // define cut-function
|
---|
137 | //
|
---|
138 | // dNOMLOGSIZE = 4.1 (=log(60.0)
|
---|
139 | // dNOMCOSZA = 1.0
|
---|
140 | //
|
---|
141 | // a: array of cut parameters
|
---|
142 | // ls: log(SIZE) - dNOMLOGSIZE
|
---|
143 | // ls2: ls^2
|
---|
144 | // ct: Cos(ZA.) - dNOMCOSZA
|
---|
145 | // dd2: DIST^2
|
---|
146 | const Double_t limit =
|
---|
147 | a[0] + a[1] * dd2 + a[2] * ct +
|
---|
148 | ls * (a[3] + a[4] * dd2 + a[5] * ct) +
|
---|
149 | ls2 * (a[6] + a[7] * dd2);
|
---|
150 |
|
---|
151 | //*fLog << "MCT1SupercutsCalc::CtsMCut; *a = "
|
---|
152 | // << *a << ", " << *(a+1) << ", " << *(a+2) << ", "
|
---|
153 | // << *(a+3) << ", " << *(a+4) << ", " << *(a+5) << ", "
|
---|
154 | // << *(a+6) << ", " << *(a+7) << endl;
|
---|
155 |
|
---|
156 | //*fLog << "MCT1SupercutsCalc::CtsMCut; ls, ls2, ct, dd2, limit = " << ls
|
---|
157 | // << ", " << ls2 << ", " << ct << ", " << dd2 << ", "
|
---|
158 | // << limit << endl;
|
---|
159 |
|
---|
160 | return limit;
|
---|
161 | }
|
---|
162 |
|
---|
163 | // --------------------------------------------------------------------------
|
---|
164 | //
|
---|
165 | // Returns the mapped value from the Matrix
|
---|
166 | //
|
---|
167 | Double_t MCT1SupercutsCalc::GetVal(Int_t i) const
|
---|
168 | {
|
---|
169 | return (*fMatrix)[fMap[i]];
|
---|
170 | }
|
---|
171 |
|
---|
172 | // --------------------------------------------------------------------------
|
---|
173 | //
|
---|
174 | // You can use this function if you want to use a MHMatrix instead of the
|
---|
175 | // given containers. This function adds all necessary columns to the
|
---|
176 | // given matrix. Afterward you should fill the matrix with the corresponding
|
---|
177 | // data (eg from a file by using MHMatrix::Fill). If you now loop
|
---|
178 | // through the matrix (eg using MMatrixLoop) MEnergyEstParam::Process
|
---|
179 | // will take the values from the matrix instead of the containers.
|
---|
180 | //
|
---|
181 | void MCT1SupercutsCalc::InitMapping(MHMatrix *mat)
|
---|
182 | {
|
---|
183 | if (fMatrix)
|
---|
184 | return;
|
---|
185 |
|
---|
186 | fMatrix = mat;
|
---|
187 |
|
---|
188 | fMap[0] = fMatrix->AddColumn("MMcEvt.fTelescopeTheta");
|
---|
189 | fMap[1] = fMatrix->AddColumn("MHillas.fWidth");
|
---|
190 | fMap[2] = fMatrix->AddColumn("MHillas.fLength");
|
---|
191 | fMap[3] = fMatrix->AddColumn("MHillas.fSize");
|
---|
192 | fMap[4] = fMatrix->AddColumn("MHillas.fMeanX");
|
---|
193 | fMap[5] = fMatrix->AddColumn("MHillas.fMeanY");
|
---|
194 | fMap[6] = fMatrix->AddColumn("MHillasSrc.fDist");
|
---|
195 | fMap[7] = fMatrix->AddColumn("fabs(MHillasSrc.fAlpha)");
|
---|
196 | }
|
---|
197 |
|
---|
198 | // ---------------------------------------------------------------------------
|
---|
199 | //
|
---|
200 | // Evaluate dynamical supercuts for CT1 Mkn421 2001 (Daniel Kranich)
|
---|
201 | // optimized for mkn 421 2001 data
|
---|
202 | //
|
---|
203 | // set hadronness to 0.25 if cuts are fullfilled
|
---|
204 | // 0.75 otherwise
|
---|
205 | //
|
---|
206 | Int_t MCT1SupercutsCalc::Process()
|
---|
207 | {
|
---|
208 | const Double_t kNomLogSize = 4.1;
|
---|
209 | const Double_t kNomCosZA = 1.0;
|
---|
210 |
|
---|
211 | const Double_t theta = fMatrix ? GetVal(0) : fMcEvt->GetTelescopeTheta();
|
---|
212 | const Double_t width0 = fMatrix ? GetVal(1) : fHil->GetWidth();
|
---|
213 | const Double_t length0 = fMatrix ? GetVal(2) : fHil->GetLength();
|
---|
214 | const Double_t size = fMatrix ? GetVal(3) : fHil->GetSize();
|
---|
215 | const Double_t meanx = fMatrix ? GetVal(4) : fHil->GetMeanX();
|
---|
216 | const Double_t meany = fMatrix ? GetVal(5) : fHil->GetMeanY();
|
---|
217 | const Double_t dist0 = fMatrix ? GetVal(6) : fHilSrc->GetDist();
|
---|
218 |
|
---|
219 | const Double_t newdist = dist0 * fMm2Deg;
|
---|
220 |
|
---|
221 | const Double_t dist2 = meanx*meanx + meany*meany;
|
---|
222 | const Double_t dist = sqrt(dist2) * fMm2Deg;
|
---|
223 | const Double_t dd2 = dist*dist;
|
---|
224 |
|
---|
225 |
|
---|
226 | const Double_t dmls = log(size) - kNomLogSize;
|
---|
227 | const Double_t dmls2 = dmls * dmls;
|
---|
228 |
|
---|
229 | const Double_t dmcza = cos(theta) - kNomCosZA;
|
---|
230 |
|
---|
231 | const Double_t length = length0 * fMm2Deg;
|
---|
232 | const Double_t width = width0 * fMm2Deg;
|
---|
233 |
|
---|
234 | if (newdist < 1.05 &&
|
---|
235 | newdist < CtsMCut (fSuper->GetDistUp(), dmls, dmcza, dmls2, dd2) &&
|
---|
236 | newdist > CtsMCut (fSuper->GetDistLo(), dmls, dmcza, dmls2, dd2) &&
|
---|
237 | dist < 1.05 &&
|
---|
238 | length < CtsMCut (fSuper->GetLengthUp(), dmls, dmcza, dmls2, dd2) &&
|
---|
239 | length > CtsMCut (fSuper->GetLengthLo(), dmls, dmcza, dmls2, dd2) &&
|
---|
240 | width < CtsMCut (fSuper->GetWidthUp(), dmls, dmcza, dmls2, dd2) &&
|
---|
241 | width > CtsMCut (fSuper->GetWidthLo(), dmls, dmcza, dmls2, dd2) &&
|
---|
242 | //asym < CtsMCut (fSuper->GetAsymUp(), dmls, dmcza, dmls2, dd2) &&
|
---|
243 | //asym > CtsMCut (fSuper->GetAsymLo(), dmls, dmcza, dmls2, dd2) &&
|
---|
244 | dist < CtsMCut (fSuper->GetDistUp(), dmls, dmcza, dmls2, dd2) &&
|
---|
245 | dist > CtsMCut (fSuper->GetDistLo(), dmls, dmcza, dmls2, dd2) )
|
---|
246 | fHadronness->SetHadronness(0.25);
|
---|
247 | else
|
---|
248 | fHadronness->SetHadronness(0.75);
|
---|
249 |
|
---|
250 | fHadronness->SetReadyToSave();
|
---|
251 |
|
---|
252 | return kTRUE;
|
---|
253 | }
|
---|
254 | //==========================================================================
|
---|
255 |
|
---|
256 |
|
---|
257 |
|
---|
258 |
|
---|