| 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, 1/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
|
|---|
| 19 | ! Author(s): Markus Gaug, 4/2004 <mailto:markus@ifae.es>
|
|---|
| 20 | !
|
|---|
| 21 | ! Copyright: MAGIC Software Development, 2000-2004
|
|---|
| 22 | !
|
|---|
| 23 | !
|
|---|
| 24 | \* ======================================================================== */
|
|---|
| 25 |
|
|---|
| 26 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 27 | //
|
|---|
| 28 | // MJPedestal
|
|---|
| 29 | //
|
|---|
| 30 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 31 | #include "MJPedestal.h"
|
|---|
| 32 |
|
|---|
| 33 | #include <TF1.h>
|
|---|
| 34 | #include <TEnv.h>
|
|---|
| 35 | #include <TFile.h>
|
|---|
| 36 | #include <TLine.h>
|
|---|
| 37 | #include <TLatex.h>
|
|---|
| 38 | #include <TString.h>
|
|---|
| 39 | #include <TCanvas.h>
|
|---|
| 40 | #include <TSystem.h>
|
|---|
| 41 | #include <TLegend.h>
|
|---|
| 42 |
|
|---|
| 43 | #include "MLog.h"
|
|---|
| 44 | #include "MLogManip.h"
|
|---|
| 45 |
|
|---|
| 46 | #include "MTaskEnv.h"
|
|---|
| 47 | #include "MSequence.h"
|
|---|
| 48 | #include "MRunIter.h"
|
|---|
| 49 | #include "MParList.h"
|
|---|
| 50 | #include "MTaskList.h"
|
|---|
| 51 | #include "MEvtLoop.h"
|
|---|
| 52 | #include "MExtractor.h"
|
|---|
| 53 |
|
|---|
| 54 | #include "MStatusDisplay.h"
|
|---|
| 55 |
|
|---|
| 56 | #include "MGeomCam.h"
|
|---|
| 57 | #include "MHCamera.h"
|
|---|
| 58 | #include "MPedestalCam.h"
|
|---|
| 59 | #include "MBadPixelsCam.h"
|
|---|
| 60 |
|
|---|
| 61 | #include "MReadMarsFile.h"
|
|---|
| 62 | #include "MRawFileRead.h"
|
|---|
| 63 | #include "MGeomApply.h"
|
|---|
| 64 | #include "MBadPixelsMerge.h"
|
|---|
| 65 | #include "MPedCalcPedRun.h"
|
|---|
| 66 | #include "MPedCalcFromLoGain.h"
|
|---|
| 67 |
|
|---|
| 68 | ClassImp(MJPedestal);
|
|---|
| 69 |
|
|---|
| 70 | using namespace std;
|
|---|
| 71 |
|
|---|
| 72 | const Double_t MJPedestal::fgPedestalMin = 4.;
|
|---|
| 73 | const Double_t MJPedestal::fgPedestalMax = 16.;
|
|---|
| 74 | const Double_t MJPedestal::fgPedRmsMin = 0.;
|
|---|
| 75 | const Double_t MJPedestal::fgPedRmsMax = 20.;
|
|---|
| 76 |
|
|---|
| 77 | const Float_t MJPedestal::fgRefPedClosedLids = 9.635;
|
|---|
| 78 | const Float_t MJPedestal::fgRefPedExtraGalactic = 9.93;
|
|---|
| 79 | const Float_t MJPedestal::fgRefPedGalactic = 10.03;
|
|---|
| 80 | const Float_t MJPedestal::fgRefPedRmsClosedLidsInner = 1.7;
|
|---|
| 81 | const Float_t MJPedestal::fgRefPedRmsExtraGalacticInner = 5.6;
|
|---|
| 82 | const Float_t MJPedestal::fgRefPedRmsGalacticInner = 6.92;
|
|---|
| 83 | const Float_t MJPedestal::fgRefPedRmsClosedLidsOuter = 1.7;
|
|---|
| 84 | const Float_t MJPedestal::fgRefPedRmsExtraGalacticOuter = 3.35;
|
|---|
| 85 | const Float_t MJPedestal::fgRefPedRmsGalacticOuter = 4.2;
|
|---|
| 86 |
|
|---|
| 87 | // --------------------------------------------------------------------------
|
|---|
| 88 | //
|
|---|
| 89 | // Default constructor.
|
|---|
| 90 | //
|
|---|
| 91 | // Sets:
|
|---|
| 92 | // - fRuns to 0,
|
|---|
| 93 | // - fExtractor to NULL,
|
|---|
| 94 | // - fDataCheck to kFALSE,
|
|---|
| 95 | // - fExtractType to kUsePedRun
|
|---|
| 96 | // - fStorage to Normal Storage
|
|---|
| 97 | //
|
|---|
| 98 | MJPedestal::MJPedestal(const char *name, const char *title)
|
|---|
| 99 | : fRuns(0), fExtractor(NULL), fDisplayType(kNormalDisplay),
|
|---|
| 100 | fDataCheck(kFALSE)
|
|---|
| 101 | {
|
|---|
| 102 | fName = name ? name : "MJPedestal";
|
|---|
| 103 | fTitle = title ? title : "Tool to create a pedestal file (MPedestalCam)";
|
|---|
| 104 |
|
|---|
| 105 | SetNormalStorage();
|
|---|
| 106 | SetUsePedRun();
|
|---|
| 107 | }
|
|---|
| 108 |
|
|---|
| 109 | const char* MJPedestal::GetOutputFile() const
|
|---|
| 110 | {
|
|---|
| 111 | if (fSequence.IsValid())
|
|---|
| 112 | return Form("%s/pedest%06d.root", (const char*)fPathOut, fSequence.GetSequence());
|
|---|
| 113 |
|
|---|
| 114 | if (!fRuns)
|
|---|
| 115 | return "";
|
|---|
| 116 |
|
|---|
| 117 | return Form("%s/%s-F0.root", (const char*)fPathOut, (const char*)fRuns->GetRunsAsFileName());
|
|---|
| 118 | }
|
|---|
| 119 |
|
|---|
| 120 | //---------------------------------------------------------------------------------
|
|---|
| 121 | //
|
|---|
| 122 | // Try to read an existing MPedestalCam from a previously created output file.
|
|---|
| 123 | // If found, also an MBadPixelsCam and the corresponding display is read.
|
|---|
| 124 | //
|
|---|
| 125 | // In case of Storage type "kNoStorage" or if the file is not found or the
|
|---|
| 126 | // MPedestalCam cannot be read, return kFALSE, otherwise kTRUE.
|
|---|
| 127 | //
|
|---|
| 128 | Bool_t MJPedestal::ReadPedestalCam()
|
|---|
| 129 | {
|
|---|
| 130 | if (IsNoStorage())
|
|---|
| 131 | return kFALSE;
|
|---|
| 132 |
|
|---|
| 133 | const TString fname = GetOutputFile();
|
|---|
| 134 |
|
|---|
| 135 | if (gSystem->AccessPathName(fname, kFileExists))
|
|---|
| 136 | {
|
|---|
| 137 | *fLog << warn << "Input file " << fname << " doesn't exist, will create it." << endl;
|
|---|
| 138 | return kFALSE;
|
|---|
| 139 | }
|
|---|
| 140 |
|
|---|
| 141 | *fLog << inf << "Reading from file: " << fname << endl;
|
|---|
| 142 |
|
|---|
| 143 | TFile file(fname, "READ");
|
|---|
| 144 | if (fPedestalCam.Read()<=0)
|
|---|
| 145 | {
|
|---|
| 146 | *fLog << err << "Unable to read MPedestalCam from " << fname << endl;
|
|---|
| 147 | return kFALSE;
|
|---|
| 148 | }
|
|---|
| 149 |
|
|---|
| 150 | if (file.FindKey("MBadPixelsCam"))
|
|---|
| 151 | {
|
|---|
| 152 | MBadPixelsCam bad;
|
|---|
| 153 | if (bad.Read()<=0)
|
|---|
| 154 | {
|
|---|
| 155 | *fLog << err << "Unable to read MBadPixelsCam from " << fname << endl;
|
|---|
| 156 | return kFALSE;
|
|---|
| 157 | }
|
|---|
| 158 | fBadPixels.Merge(bad);
|
|---|
| 159 | }
|
|---|
| 160 |
|
|---|
| 161 | if (fDisplay && !fDisplay->GetCanvas("Pedestals"))
|
|---|
| 162 | fDisplay->Read();
|
|---|
| 163 |
|
|---|
| 164 | return kTRUE;
|
|---|
| 165 | }
|
|---|
| 166 |
|
|---|
| 167 | //---------------------------------------------------------------------------------
|
|---|
| 168 | //
|
|---|
| 169 | // Display the results.
|
|---|
| 170 | // If Display type "kDataCheck" was chosen, also the reference lines are displayed.
|
|---|
| 171 | //
|
|---|
| 172 | void MJPedestal::DisplayResult(MParList &plist)
|
|---|
| 173 | {
|
|---|
| 174 | if (!fDisplay)
|
|---|
| 175 | return;
|
|---|
| 176 |
|
|---|
| 177 | //
|
|---|
| 178 | // Update display
|
|---|
| 179 | //
|
|---|
| 180 | TString title = fDisplay->GetTitle();
|
|---|
| 181 | title += "-- Pedestal ";
|
|---|
| 182 | if (fSequence.IsValid())
|
|---|
| 183 | title += fSequence.GetName();
|
|---|
| 184 | else
|
|---|
| 185 | if (fRuns) // FIXME: What to do if an environmentfile was used?
|
|---|
| 186 | title += fRuns->GetRunsAsString();
|
|---|
| 187 | title += " --";
|
|---|
| 188 | fDisplay->SetTitle(title);
|
|---|
| 189 |
|
|---|
| 190 | //
|
|---|
| 191 | // Get container from list
|
|---|
| 192 | //
|
|---|
| 193 | MGeomCam &geomcam = *(MGeomCam*)plist.FindObject("MGeomCam");
|
|---|
| 194 |
|
|---|
| 195 | //
|
|---|
| 196 | // Create container to display
|
|---|
| 197 | //
|
|---|
| 198 | MHCamera disp0(geomcam, "MPedestalCam;ped", "Mean Pedestal");
|
|---|
| 199 | MHCamera disp1(geomcam, "MPedestalCam;RMS", "Pedestal RMS");
|
|---|
| 200 |
|
|---|
| 201 | disp0.SetCamContent(fPedestalCam, 0);
|
|---|
| 202 | disp0.SetCamError (fPedestalCam, 1);
|
|---|
| 203 |
|
|---|
| 204 | disp1.SetCamContent(fPedestalCam, 2);
|
|---|
| 205 | disp1.SetCamError (fPedestalCam, 3);
|
|---|
| 206 |
|
|---|
| 207 | disp0.SetYTitle("Pedestal [counts/slice]");
|
|---|
| 208 | disp1.SetYTitle("RMS [counts/slice]");
|
|---|
| 209 |
|
|---|
| 210 | //
|
|---|
| 211 | // Display data
|
|---|
| 212 | //
|
|---|
| 213 | TCanvas &c3 = fDisplay->AddTab("Pedestals");
|
|---|
| 214 | c3.Divide(2,3);
|
|---|
| 215 |
|
|---|
| 216 | if (fDisplayType != kDataCheckDisplay)
|
|---|
| 217 | {
|
|---|
| 218 | disp0.CamDraw(c3, 1, 2, 1);
|
|---|
| 219 | disp1.CamDraw(c3, 2, 2, 6);
|
|---|
| 220 | return;
|
|---|
| 221 | }
|
|---|
| 222 |
|
|---|
| 223 | c3.cd(1);
|
|---|
| 224 | gPad->SetBorderMode(0);
|
|---|
| 225 | gPad->SetTicks();
|
|---|
| 226 | MHCamera *obj1=(MHCamera*)disp0.DrawCopy("hist");
|
|---|
| 227 | //
|
|---|
| 228 | // for the datacheck, fix the ranges!!
|
|---|
| 229 | //
|
|---|
| 230 | obj1->SetMinimum(fgPedestalMin);
|
|---|
| 231 | obj1->SetMaximum(fgPedestalMax);
|
|---|
| 232 |
|
|---|
| 233 | //
|
|---|
| 234 | // Set the datacheck sizes:
|
|---|
| 235 | //
|
|---|
| 236 | FixDataCheckHist((TH1D*)obj1);
|
|---|
| 237 | //
|
|---|
| 238 | // set reference lines
|
|---|
| 239 | //
|
|---|
| 240 | DisplayReferenceLines(obj1,0);
|
|---|
| 241 |
|
|---|
| 242 | //
|
|---|
| 243 | // end reference lines
|
|---|
| 244 | //
|
|---|
| 245 | c3.cd(3);
|
|---|
| 246 | gPad->SetBorderMode(0);
|
|---|
| 247 | obj1->SetPrettyPalette();
|
|---|
| 248 | obj1->Draw();
|
|---|
| 249 |
|
|---|
| 250 | c3.cd(5);
|
|---|
| 251 | gPad->SetBorderMode(0);
|
|---|
| 252 | gPad->SetTicks();
|
|---|
| 253 | TH1D *obj2 = (TH1D*)obj1->Projection(obj1->GetName());
|
|---|
| 254 | obj2->Draw();
|
|---|
| 255 | obj2->SetBit(kCanDelete);
|
|---|
| 256 | obj2->Fit("gaus","Q");
|
|---|
| 257 | obj2->GetFunction("gaus")->SetLineColor(kYellow);
|
|---|
| 258 | //
|
|---|
| 259 | // Set the datacheck sizes:
|
|---|
| 260 | //
|
|---|
| 261 | FixDataCheckHist(obj2);
|
|---|
| 262 | obj2->SetStats(1);
|
|---|
| 263 |
|
|---|
| 264 | c3.cd(2);
|
|---|
| 265 | gPad->SetBorderMode(0);
|
|---|
| 266 | gPad->SetTicks();
|
|---|
| 267 | MHCamera *obj3=(MHCamera*)disp1.DrawCopy("hist");
|
|---|
| 268 | //
|
|---|
| 269 | // for the datacheck, fix the ranges!!
|
|---|
| 270 | //
|
|---|
| 271 | obj3->SetMinimum(fgPedRmsMin);
|
|---|
| 272 | obj3->SetMaximum(fgPedRmsMax);
|
|---|
| 273 | //
|
|---|
| 274 | // Set the datacheck sizes:
|
|---|
| 275 | //
|
|---|
| 276 | FixDataCheckHist((TH1D*)obj3);
|
|---|
| 277 | //
|
|---|
| 278 | // set reference lines
|
|---|
| 279 | //
|
|---|
| 280 | DisplayReferenceLines(obj1,1);
|
|---|
| 281 |
|
|---|
| 282 | c3.cd(4);
|
|---|
| 283 | gPad->SetBorderMode(0);
|
|---|
| 284 | obj3->SetPrettyPalette();
|
|---|
| 285 | obj3->Draw();
|
|---|
| 286 |
|
|---|
| 287 | c3.cd(6);
|
|---|
| 288 | gPad->SetBorderMode(0);
|
|---|
| 289 |
|
|---|
| 290 | if (geomcam.InheritsFrom("MGeomCamMagic"))
|
|---|
| 291 | {
|
|---|
| 292 | TArrayI inner(1);
|
|---|
| 293 | inner[0] = 0;
|
|---|
| 294 |
|
|---|
| 295 | TArrayI outer(1);
|
|---|
| 296 | outer[0] = 1;
|
|---|
| 297 |
|
|---|
| 298 | TArrayI s0(6);
|
|---|
| 299 | s0[0] = 6;
|
|---|
| 300 | s0[1] = 1;
|
|---|
| 301 | s0[2] = 2;
|
|---|
| 302 | s0[3] = 3;
|
|---|
| 303 | s0[4] = 4;
|
|---|
| 304 | s0[5] = 5;
|
|---|
| 305 |
|
|---|
| 306 | TArrayI s1(3);
|
|---|
| 307 | s1[0] = 6;
|
|---|
| 308 | s1[1] = 1;
|
|---|
| 309 | s1[2] = 2;
|
|---|
| 310 |
|
|---|
| 311 | TArrayI s2(3);
|
|---|
| 312 | s2[0] = 3;
|
|---|
| 313 | s2[1] = 4;
|
|---|
| 314 | s2[2] = 5;
|
|---|
| 315 |
|
|---|
| 316 | TVirtualPad *pad = gPad;
|
|---|
| 317 | pad->Divide(2,1);
|
|---|
| 318 |
|
|---|
| 319 | TH1D *inout[2];
|
|---|
| 320 | inout[0] = disp1.ProjectionS(s0, inner, "Inner");
|
|---|
| 321 | inout[1] = disp1.ProjectionS(s0, outer, "Outer");
|
|---|
| 322 | FixDataCheckHist(inout[0]);
|
|---|
| 323 | FixDataCheckHist(inout[1]);
|
|---|
| 324 |
|
|---|
| 325 | inout[0]->SetTitle(Form("%s %s",disp1.GetTitle(),"Inner"));
|
|---|
| 326 | inout[1]->SetTitle(Form("%s %s",disp1.GetTitle(),"Outer"));
|
|---|
| 327 |
|
|---|
| 328 |
|
|---|
| 329 | for (int i=0; i<2; i++)
|
|---|
| 330 | {
|
|---|
| 331 | pad->cd(i+1);
|
|---|
| 332 | gPad->SetBorderMode(0);
|
|---|
| 333 | gPad->SetTicks();
|
|---|
| 334 |
|
|---|
| 335 | inout[i]->SetDirectory(NULL);
|
|---|
| 336 | inout[i]->SetLineColor(kRed+i);
|
|---|
| 337 | inout[i]->SetBit(kCanDelete);
|
|---|
| 338 | inout[i]->Draw();
|
|---|
| 339 | inout[i]->Fit("gaus", "Q");
|
|---|
| 340 |
|
|---|
| 341 | TLegend *leg2 = new TLegend(0.6,0.2,0.9,0.55);
|
|---|
| 342 | leg2->SetHeader(inout[i]->GetName());
|
|---|
| 343 | leg2->AddEntry(inout[i], inout[i]->GetName(), "l");
|
|---|
| 344 |
|
|---|
| 345 | //
|
|---|
| 346 | // Display the outliers as dead and noisy pixels
|
|---|
| 347 | //
|
|---|
| 348 | DisplayOutliers(inout[i]);
|
|---|
| 349 |
|
|---|
| 350 | //
|
|---|
| 351 | // Display the two half of the camera separately
|
|---|
| 352 | //
|
|---|
| 353 | TH1D *half[2];
|
|---|
| 354 | half[0] = disp1.ProjectionS(s1, i==0 ? inner : outer , "Sector 6-1-2");
|
|---|
| 355 | half[1] = disp1.ProjectionS(s2, i==0 ? inner : outer , "Sector 3-4-5");
|
|---|
| 356 |
|
|---|
| 357 | for (int j=0; j<2; j++)
|
|---|
| 358 | {
|
|---|
| 359 | half[j]->SetLineColor(kRed+i+2*j+1);
|
|---|
| 360 | half[j]->SetDirectory(NULL);
|
|---|
| 361 | half[j]->SetBit(kCanDelete);
|
|---|
| 362 | half[j]->Draw("same");
|
|---|
| 363 | leg2->AddEntry(half[j], half[j]->GetName(), "l");
|
|---|
| 364 | }
|
|---|
| 365 | leg2->Draw();
|
|---|
| 366 | }
|
|---|
| 367 | }
|
|---|
| 368 | }
|
|---|
| 369 |
|
|---|
| 370 | void MJPedestal::DisplayReferenceLines(MHCamera *cam, const Int_t what) const
|
|---|
| 371 | {
|
|---|
| 372 | Double_t x = cam->GetNbinsX();
|
|---|
| 373 |
|
|---|
| 374 | const MGeomCam *geom = cam->GetGeometry();
|
|---|
| 375 |
|
|---|
| 376 | if (geom->InheritsFrom("MGeomCamMagic"))
|
|---|
| 377 | x = what ? 397 : cam->GetNbinsX();
|
|---|
| 378 |
|
|---|
| 379 | TLine line;
|
|---|
| 380 | line.SetLineStyle(kDashed);
|
|---|
| 381 | line.SetLineWidth(3);
|
|---|
| 382 |
|
|---|
| 383 | line.SetLineColor(kBlue);
|
|---|
| 384 | TLine *l1 = line.DrawLine(0, what ? fgRefPedRmsGalacticInner : fgRefPedGalactic,
|
|---|
| 385 | x, what ? fgRefPedRmsGalacticInner : fgRefPedGalactic);
|
|---|
| 386 |
|
|---|
| 387 | line.SetLineColor(kYellow);
|
|---|
| 388 | TLine *l2 = line.DrawLine(0, what ? fgRefPedRmsExtraGalacticInner : fgRefPedExtraGalactic,
|
|---|
| 389 | x, what ? fgRefPedRmsExtraGalacticInner : fgRefPedExtraGalactic);
|
|---|
| 390 |
|
|---|
| 391 | line.SetLineColor(kMagenta);
|
|---|
| 392 | TLine *l3 = line.DrawLine(0, what ? fgRefPedRmsClosedLidsInner : fgRefPedClosedLids,
|
|---|
| 393 | x, what ? fgRefPedRmsClosedLidsInner : fgRefPedClosedLids);
|
|---|
| 394 |
|
|---|
| 395 | if (geom->InheritsFrom("MGeomCamMagic"))
|
|---|
| 396 | if (what)
|
|---|
| 397 | {
|
|---|
| 398 | const Double_t x2 = cam->GetNbinsX();
|
|---|
| 399 |
|
|---|
| 400 | line.SetLineColor(kBlue);
|
|---|
| 401 | line.DrawLine(398, fgRefPedRmsGalacticOuter,
|
|---|
| 402 | x2, fgRefPedRmsGalacticOuter);
|
|---|
| 403 |
|
|---|
| 404 | line.SetLineColor(kYellow);
|
|---|
| 405 | line.DrawLine(398, fgRefPedRmsExtraGalacticOuter,
|
|---|
| 406 | x2, fgRefPedRmsExtraGalacticOuter);
|
|---|
| 407 |
|
|---|
| 408 | line.SetLineColor(kMagenta);
|
|---|
| 409 | line.DrawLine(398, fgRefPedRmsClosedLidsOuter,
|
|---|
| 410 | x2, fgRefPedRmsClosedLidsOuter);
|
|---|
| 411 | }
|
|---|
| 412 |
|
|---|
| 413 |
|
|---|
| 414 | TLegend *leg = new TLegend(0.4,0.75,0.7,0.99);
|
|---|
| 415 | leg->SetBit(kCanDelete);
|
|---|
| 416 | leg->AddEntry(l1, "Galactic Source","l");
|
|---|
| 417 | leg->AddEntry(l2, "Extra-Galactic Source","l");
|
|---|
| 418 | leg->AddEntry(l3, "Closed Lids","l");
|
|---|
| 419 | leg->Draw();
|
|---|
| 420 | }
|
|---|
| 421 |
|
|---|
| 422 | void MJPedestal::DisplayOutliers(TH1D *hist) const
|
|---|
| 423 | {
|
|---|
| 424 | const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
|
|---|
| 425 | const Float_t lolim = mean - 3.5*hist->GetFunction("gaus")->GetParameter(2);
|
|---|
| 426 | const Float_t uplim = mean + 3.5*hist->GetFunction("gaus")->GetParameter(2);
|
|---|
| 427 | const Stat_t dead = hist->Integral(0,hist->FindBin(lolim)-1);
|
|---|
| 428 | const Stat_t noisy = hist->Integral(hist->FindBin(uplim)+1,hist->GetNbinsX()+1);
|
|---|
| 429 |
|
|---|
| 430 | TLatex deadtex;
|
|---|
| 431 | deadtex.SetTextSize(0.06);
|
|---|
| 432 | deadtex.DrawLatex(0.1,hist->GetBinContent(hist->GetMaximumBin())/1.1,Form("%3i dead pixels",(Int_t)dead));
|
|---|
| 433 |
|
|---|
| 434 | TLatex noisytex;
|
|---|
| 435 | noisytex.SetTextSize(0.06);
|
|---|
| 436 | noisytex.DrawLatex(0.1,hist->GetBinContent(hist->GetMaximumBin())/1.2,Form("%3i noisy pixels",(Int_t)noisy));
|
|---|
| 437 | }
|
|---|
| 438 |
|
|---|
| 439 | void MJPedestal::FixDataCheckHist(TH1D *hist) const
|
|---|
| 440 | {
|
|---|
| 441 | hist->SetDirectory(NULL);
|
|---|
| 442 | hist->SetStats(0);
|
|---|
| 443 |
|
|---|
| 444 | //
|
|---|
| 445 | // set the labels bigger
|
|---|
| 446 | //
|
|---|
| 447 | TAxis *xaxe = hist->GetXaxis();
|
|---|
| 448 | TAxis *yaxe = hist->GetYaxis();
|
|---|
| 449 |
|
|---|
| 450 | xaxe->CenterTitle();
|
|---|
| 451 | yaxe->CenterTitle();
|
|---|
| 452 | xaxe->SetTitleSize(0.06);
|
|---|
| 453 | yaxe->SetTitleSize(0.06);
|
|---|
| 454 | xaxe->SetTitleOffset(0.8);
|
|---|
| 455 | yaxe->SetTitleOffset(0.5);
|
|---|
| 456 | xaxe->SetLabelSize(0.05);
|
|---|
| 457 | yaxe->SetLabelSize(0.05);
|
|---|
| 458 | }
|
|---|
| 459 |
|
|---|
| 460 | /*
|
|---|
| 461 | Bool_t MJPedestal::WriteEventloop(MEvtLoop &evtloop) const
|
|---|
| 462 | {
|
|---|
| 463 | if (fOutputPath.IsNull())
|
|---|
| 464 | return kTRUE;
|
|---|
| 465 |
|
|---|
| 466 | const TString oname(GetOutputFile());
|
|---|
| 467 |
|
|---|
| 468 | *fLog << inf << "Writing to file: " << oname << endl;
|
|---|
| 469 |
|
|---|
| 470 | TFile file(oname, fOverwrite?"RECREATE":"NEW", "File created by MJPedestal", 9);
|
|---|
| 471 | if (!file.IsOpen())
|
|---|
| 472 | {
|
|---|
| 473 | *fLog << err << "ERROR - Couldn't open file " << oname << " for writing..." << endl;
|
|---|
| 474 | return kFALSE;
|
|---|
| 475 | }
|
|---|
| 476 |
|
|---|
| 477 | if (evtloop.Write(fName)<=0)
|
|---|
| 478 | {
|
|---|
| 479 | *fLog << err << "Unable to write MEvtloop to " << oname << endl;
|
|---|
| 480 | return kFALSE;
|
|---|
| 481 | }
|
|---|
| 482 |
|
|---|
| 483 | return kTRUE;
|
|---|
| 484 | }
|
|---|
| 485 | */
|
|---|
| 486 |
|
|---|
| 487 | // --------------------------------------------------------------------------
|
|---|
| 488 | //
|
|---|
| 489 | // The following resource options are available:
|
|---|
| 490 | //
|
|---|
| 491 | // Do a datacheck run (read raw-data and enable display)
|
|---|
| 492 | // Prefix.DataCheck: Yes, No <default>
|
|---|
| 493 | //
|
|---|
| 494 | // Show data check display
|
|---|
| 495 | // Prefix.DataCheckDisplay: Yes, No <default>
|
|---|
| 496 | //
|
|---|
| 497 | // Use cosmic data instead of pedestal data (DatRuns)
|
|---|
| 498 | // Prefix.UseData: Yes, No <default>
|
|---|
| 499 | //
|
|---|
| 500 | // Write an output file with pedestals and status-display
|
|---|
| 501 | // Prefix.DisableOutput: Yes, No <default>
|
|---|
| 502 | //
|
|---|
| 503 | Bool_t MJPedestal::CheckEnvLocal()
|
|---|
| 504 | {
|
|---|
| 505 |
|
|---|
| 506 | SetDataCheck(GetEnv("DataCheck", fDataCheck));
|
|---|
| 507 |
|
|---|
| 508 | if (HasEnv("DataCheckDisplay"))
|
|---|
| 509 | fDisplayType = GetEnv("DataCheckDisplay", kFALSE) ? kDataCheckDisplay : kNormalDisplay;
|
|---|
| 510 |
|
|---|
| 511 | if (HasEnv("UseData"))
|
|---|
| 512 | fExtractType = GetEnv("UseData",kFALSE) ? kUseData : kUsePedRun;
|
|---|
| 513 |
|
|---|
| 514 | if (HasEnv("UseExtractor"))
|
|---|
| 515 | if (GetEnv("UseExtractor",kFALSE))
|
|---|
| 516 | fExtractType = kUseExtractor;
|
|---|
| 517 |
|
|---|
| 518 | SetNoStorage ( GetEnv("DisableOutput", IsNoStorage() ));
|
|---|
| 519 |
|
|---|
| 520 | return kTRUE;
|
|---|
| 521 | }
|
|---|
| 522 |
|
|---|
| 523 | //---------------------------------------------------------------------------------
|
|---|
| 524 | //
|
|---|
| 525 | // Try to write the created MPedestalCam in the output file.
|
|---|
| 526 | // If possible, also an MBadPixelsCam and the corresponding display is written.
|
|---|
| 527 | //
|
|---|
| 528 | // In case of Storage type "kNoStorage" or if any of the containers
|
|---|
| 529 | // cannot be written, return kFALSE, otherwise kTRUE.
|
|---|
| 530 | //
|
|---|
| 531 | Bool_t MJPedestal::WriteResult()
|
|---|
| 532 | {
|
|---|
| 533 | if (IsNoStorage())
|
|---|
| 534 | return kTRUE;
|
|---|
| 535 |
|
|---|
| 536 | if (fPathOut.IsNull())
|
|---|
| 537 | return kTRUE;
|
|---|
| 538 |
|
|---|
| 539 | const TString oname(GetOutputFile());
|
|---|
| 540 |
|
|---|
| 541 | *fLog << inf << "Writing to file: " << oname << endl;
|
|---|
| 542 |
|
|---|
| 543 | TFile file(oname, "UPDATE", "File created by MJPedestal", 9);
|
|---|
| 544 | if (!file.IsOpen())
|
|---|
| 545 | {
|
|---|
| 546 | *fLog << err << "ERROR - Couldn't open file " << oname << " for writing..." << endl;
|
|---|
| 547 | return kFALSE;
|
|---|
| 548 | }
|
|---|
| 549 |
|
|---|
| 550 | if (fDisplay && fDisplay->Write()<=0)
|
|---|
| 551 | {
|
|---|
| 552 | *fLog << err << "Unable to write MStatusDisplay to " << oname << endl;
|
|---|
| 553 | return kFALSE;
|
|---|
| 554 | }
|
|---|
| 555 |
|
|---|
| 556 | if (fPedestalCam.Write()<=0)
|
|---|
| 557 | {
|
|---|
| 558 | *fLog << err << "Unable to write MPedestalCam to " << oname << endl;
|
|---|
| 559 | return kFALSE;
|
|---|
| 560 | }
|
|---|
| 561 |
|
|---|
| 562 | if (fBadPixels.Write()<=0)
|
|---|
| 563 | {
|
|---|
| 564 | *fLog << err << "Unable to write MBadPixelsCam to " << oname << endl;
|
|---|
| 565 | return kFALSE;
|
|---|
| 566 | }
|
|---|
| 567 |
|
|---|
| 568 | return kTRUE;
|
|---|
| 569 | }
|
|---|
| 570 |
|
|---|
| 571 | Bool_t MJPedestal::Process()
|
|---|
| 572 | {
|
|---|
| 573 | if (!ReadPedestalCam())
|
|---|
| 574 | return ProcessFile();
|
|---|
| 575 |
|
|---|
| 576 | return kTRUE;
|
|---|
| 577 | }
|
|---|
| 578 |
|
|---|
| 579 | Bool_t MJPedestal::ProcessFile()
|
|---|
| 580 | {
|
|---|
| 581 | if (!fSequence.IsValid())
|
|---|
| 582 | {
|
|---|
| 583 | if (!fRuns)
|
|---|
| 584 | {
|
|---|
| 585 | *fLog << err << "Neither AddRuns nor SetSequence nor SetEnv was called... abort." << endl;
|
|---|
| 586 | return kFALSE;
|
|---|
| 587 | }
|
|---|
| 588 | if (fRuns && fRuns->GetNumRuns() != fRuns->GetNumEntries())
|
|---|
| 589 | {
|
|---|
| 590 | *fLog << err << "Number of files found doesn't match number of runs... abort." << endl;
|
|---|
| 591 | return kFALSE;
|
|---|
| 592 | }
|
|---|
| 593 | }
|
|---|
| 594 |
|
|---|
| 595 | //if (!CheckEnv())
|
|---|
| 596 | // return kFALSE;
|
|---|
| 597 |
|
|---|
| 598 | CheckEnv();
|
|---|
| 599 |
|
|---|
| 600 | // --------------------------------------------------------------------------------
|
|---|
| 601 |
|
|---|
| 602 | const TString type = IsUseData() ? "data" : "pedestal";
|
|---|
| 603 |
|
|---|
| 604 | *fLog << inf;
|
|---|
| 605 | fLog->Separator(GetDescriptor());
|
|---|
| 606 | *fLog << "Calculate MPedestalCam from " << type << "-runs ";
|
|---|
| 607 | if (fSequence.IsValid())
|
|---|
| 608 | *fLog << fSequence.GetName() << endl;
|
|---|
| 609 | else
|
|---|
| 610 | if (fRuns)
|
|---|
| 611 | *fLog << fRuns->GetRunsAsString() << endl;
|
|---|
| 612 | else
|
|---|
| 613 | *fLog << "in Resource File." << endl;
|
|---|
| 614 | *fLog << endl;
|
|---|
| 615 |
|
|---|
| 616 | // --------------------------------------------------------------------------------
|
|---|
| 617 |
|
|---|
| 618 | MParList plist;
|
|---|
| 619 | MTaskList tlist;
|
|---|
| 620 | plist.AddToList(&tlist);
|
|---|
| 621 | plist.AddToList(this); // take care of fDisplay!
|
|---|
| 622 |
|
|---|
| 623 | MReadMarsFile read("Events");
|
|---|
| 624 | MRawFileRead rawread(NULL);
|
|---|
| 625 |
|
|---|
| 626 | MDirIter iter;
|
|---|
| 627 | if (fSequence.IsValid())
|
|---|
| 628 | {
|
|---|
| 629 | const Int_t n0 = IsUseData() ? fSequence.SetupDatRuns(iter, fPathData, "D", fDataCheck) : fSequence.SetupPedRuns(iter, fPathData, "P", fDataCheck);
|
|---|
| 630 | const Int_t n1 = IsUseData() ? fSequence.GetNumDatRuns() : fSequence.GetNumPedRuns();
|
|---|
| 631 | if (n0==0)
|
|---|
| 632 | {
|
|---|
| 633 | *fLog << err << "ERROR - No " << type << " input files of sequence found in " << (fPathData.IsNull()?"<defaul>":fPathData.Data()) << endl;
|
|---|
| 634 | return kFALSE;
|
|---|
| 635 | }
|
|---|
| 636 | if (n0!=n1)
|
|---|
| 637 | {
|
|---|
| 638 | *fLog << err << "ERROR - Number of " << type << " files found (" << n0 << ") in " << (fPathData.IsNull()?"<defaul>":fPathData.Data()) << " doesn't match number of files in sequence (" << n1 << ")" << endl;
|
|---|
| 639 | return kFALSE;
|
|---|
| 640 | }
|
|---|
| 641 | }
|
|---|
| 642 |
|
|---|
| 643 | if (fDataCheck)
|
|---|
| 644 | {
|
|---|
| 645 | if (fRuns || fSequence.IsValid())
|
|---|
| 646 | rawread.AddFiles(fSequence.IsValid() ? iter : *fRuns);
|
|---|
| 647 | tlist.AddToList(&rawread);
|
|---|
| 648 | }
|
|---|
| 649 | else
|
|---|
| 650 | {
|
|---|
| 651 | read.DisableAutoScheme();
|
|---|
| 652 | if (fRuns || fSequence.IsValid())
|
|---|
| 653 | read.AddFiles(fSequence.IsValid() ? iter : *fRuns);
|
|---|
| 654 | tlist.AddToList(&read);
|
|---|
| 655 | }
|
|---|
| 656 |
|
|---|
| 657 | // Setup Tasklist
|
|---|
| 658 | plist.AddToList(&fPedestalCam);
|
|---|
| 659 | plist.AddToList(&fBadPixels);
|
|---|
| 660 |
|
|---|
| 661 | MGeomApply geomapl;
|
|---|
| 662 | MBadPixelsMerge merge(&fBadPixels);
|
|---|
| 663 |
|
|---|
| 664 | MPedCalcPedRun pedcalc;
|
|---|
| 665 | MPedCalcFromLoGain pedlogain;
|
|---|
| 666 |
|
|---|
| 667 | MTaskEnv taskenv("ExtractPedestal");
|
|---|
| 668 | if (IsUseData())
|
|---|
| 669 | taskenv.SetDefault(&pedlogain);
|
|---|
| 670 | else if (IsUseExtractor())
|
|---|
| 671 | taskenv.SetDefault(fExtractor);
|
|---|
| 672 | else
|
|---|
| 673 | taskenv.SetDefault(&pedcalc);
|
|---|
| 674 |
|
|---|
| 675 | if (fExtractor)
|
|---|
| 676 | {
|
|---|
| 677 | if (IsUseData())
|
|---|
| 678 | pedlogain.SetExtractWindow(15,(Int_t)fExtractor->GetNumHiGainSamples());
|
|---|
| 679 | else if (IsUsePedRun())
|
|---|
| 680 | {
|
|---|
| 681 | pedcalc.SetWindowSize((Int_t)fExtractor->GetNumHiGainSamples());
|
|---|
| 682 | pedcalc.SetRange(fExtractor->GetHiGainFirst(), fExtractor->GetHiGainLast());
|
|---|
| 683 | }
|
|---|
| 684 | else if (IsUseExtractor())
|
|---|
| 685 | if (!fExtractor->IsNoiseCalculation())
|
|---|
| 686 | {
|
|---|
| 687 | *fLog << err << GetDescriptor()
|
|---|
| 688 | << "Extraction Type kUseExtractor is chosen, but the extractor has no flag kNoiseCalculation set! " << endl;
|
|---|
| 689 | return kFALSE;
|
|---|
| 690 | }
|
|---|
| 691 | }
|
|---|
| 692 | else if (IsUseExtractor())
|
|---|
| 693 | {
|
|---|
| 694 | *fLog << err << GetDescriptor()
|
|---|
| 695 | << "Extraction Type kUseExtractor is chosen, but no extractor has been handed over" << endl;
|
|---|
| 696 | return kFALSE;
|
|---|
| 697 | }
|
|---|
| 698 |
|
|---|
| 699 | /*
|
|---|
| 700 | else
|
|---|
| 701 | {
|
|---|
| 702 | *fLog << warn << GetDescriptor();
|
|---|
| 703 | *fLog << ": No extractor has been chosen, take default number of FADC samples " << endl;
|
|---|
| 704 | }
|
|---|
| 705 | */
|
|---|
| 706 |
|
|---|
| 707 | tlist.AddToList(&geomapl);
|
|---|
| 708 | tlist.AddToList(&merge);
|
|---|
| 709 | tlist.AddToList(&taskenv);
|
|---|
| 710 |
|
|---|
| 711 | //
|
|---|
| 712 | // Execute the eventloop
|
|---|
| 713 | //
|
|---|
| 714 | MEvtLoop evtloop(fName);
|
|---|
| 715 | evtloop.SetParList(&plist);
|
|---|
| 716 | evtloop.SetDisplay(fDisplay);
|
|---|
| 717 | evtloop.SetLogStream(fLog);
|
|---|
| 718 | if (!SetupEnv(evtloop))
|
|---|
| 719 | return kFALSE;
|
|---|
| 720 |
|
|---|
| 721 | // if (!WriteEventloop(evtloop))
|
|---|
| 722 | // return kFALSE;
|
|---|
| 723 |
|
|---|
| 724 | // Execute first analysis
|
|---|
| 725 | if (!evtloop.Eventloop(fMaxEvents))
|
|---|
| 726 | {
|
|---|
| 727 | *fLog << err << GetDescriptor() << ": Failed." << endl;
|
|---|
| 728 | return kFALSE;
|
|---|
| 729 | }
|
|---|
| 730 |
|
|---|
| 731 | tlist.PrintStatistics();
|
|---|
| 732 |
|
|---|
| 733 | DisplayResult(plist);
|
|---|
| 734 |
|
|---|
| 735 | if (!WriteResult())
|
|---|
| 736 | return kFALSE;
|
|---|
| 737 |
|
|---|
| 738 | *fLog << all << GetDescriptor() << ": Done." << endl;
|
|---|
| 739 | *fLog << endl << endl;
|
|---|
| 740 |
|
|---|
| 741 | return kTRUE;
|
|---|
| 742 | }
|
|---|