/* ======================================================================== *\ ! ! * ! * This file is part of MARS, the MAGIC Analysis and Reconstruction ! * Software. It is distributed to you in the hope that it can be a useful ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes. ! * It is distributed WITHOUT ANY WARRANTY. ! * ! * Permission to use, copy, modify and distribute this software and its ! * documentation for any purpose is hereby granted without fee, ! * provided that the above copyright notice appear in all copies and ! * that both that copyright notice and this permission notice appear ! * in supporting documentation. It is provided "as is" without express ! * or implied warranty. ! * ! ! ! Author(s): Thomas Bretz, 1/2004 ! Markus Gaug ,04/2004 ! ! Copyright: MAGIC Software Development, 2000-2004 ! ! \* ======================================================================== */ ///////////////////////////////////////////////////////////////////////////// // // MJPedestal // ///////////////////////////////////////////////////////////////////////////// #include "MJPedestal.h" #include #include #include #include #include #include #include #include #include #include "MLog.h" #include "MLogManip.h" #include "MRunIter.h" #include "MParList.h" #include "MTaskList.h" #include "MEvtLoop.h" #include "MExtractor.h" #include "MStatusDisplay.h" #include "MGeomCam.h" #include "MHCamera.h" #include "MPedestalCam.h" #include "MBadPixelsCam.h" #include "MReadMarsFile.h" #include "MRawFileRead.h" #include "MGeomApply.h" #include "MBadPixelsMerge.h" #include "MPedCalcPedRun.h" ClassImp(MJPedestal); using namespace std; const Double_t MJPedestal::fgPedestalMin = 4.; const Double_t MJPedestal::fgPedestalMax = 16.; const Double_t MJPedestal::fgPedRmsMin = 0.; const Double_t MJPedestal::fgPedRmsMax = 20.; const Float_t MJPedestal::fgRefPedClosedLids = 9.635; const Float_t MJPedestal::fgRefPedExtraGalactic = 9.93; const Float_t MJPedestal::fgRefPedGalactic = 10.03; const Float_t MJPedestal::fgRefPedRmsClosedLids = 1.7; const Float_t MJPedestal::fgRefPedRmsExtraGalactic = 5.6; const Float_t MJPedestal::fgRefPedRmsGalactic = 6.92; // -------------------------------------------------------------------------- // // Default constructor. // // Sets fRuns to 0, fExtractor to NULL, fDataCheck to kFALSE // MJPedestal::MJPedestal(const char *name, const char *title) : fEnv(0), fRuns(0), fExtractor(NULL), fDataCheck(kFALSE) { fName = name ? name : "MJPedestal"; fTitle = title ? title : "Tool to create a pedestal file (MPedestalCam)"; } MJPedestal::~MJPedestal() { if (fEnv) delete fEnv; } const char* MJPedestal::GetOutputFile() const { if (!fRuns) return ""; return Form("%s/%s-F0.root", (const char*)fOutputPath, (const char*)fRuns->GetRunsAsFileName()); } Bool_t MJPedestal::ReadPedestalCam() { const TString fname = GetOutputFile(); if (gSystem->AccessPathName(fname, kFileExists)) { *fLog << warn << "Input file " << fname << " doesn't exist, will create it." << endl; return kFALSE; } *fLog << inf << "Reading from file: " << fname << endl; TFile file(fname, "READ"); if (fPedestalCam.Read()<=0) { *fLog << err << "Unable to read MPedestalCam from " << fname << endl; return kFALSE; } if (file.FindKey("MBadPixelsCam")) { MBadPixelsCam bad; if (bad.Read()<=0) { *fLog << err << "Unable to read MBadPixelsCam from " << fname << endl; return kFALSE; } fBadPixels.Merge(bad); } if (fDisplay && !fDisplay->GetCanvas("Pedestals")) fDisplay->Read(); return kTRUE; } void MJPedestal::DisplayResult(MParList &plist) { if (!fDisplay) return; // // Update display // TString title = fDisplay->GetTitle(); title += "-- Pedestal "; if (fRuns) // FIXME: What to do if an environmentfile was used? title += fRuns->GetRunsAsString(); title += " --"; fDisplay->SetTitle(title); // // Get container from list // MGeomCam &geomcam = *(MGeomCam*)plist.FindObject("MGeomCam"); // // Create container to display // MHCamera disp0(geomcam, "MPedestalCam;ped", "Mean Pedestal"); MHCamera disp1(geomcam, "MPedestalCam;RMS", "Pedestal RMS"); disp0.SetCamContent(fPedestalCam, 0); disp0.SetCamError (fPedestalCam, 1); disp1.SetCamContent(fPedestalCam, 2); disp1.SetCamError (fPedestalCam, 3); disp0.SetYTitle("P [fadc/slice]"); disp1.SetYTitle("\\sigma_{P} [fadc/slice]"); // // Display data // TCanvas &c3 = fDisplay->AddTab("Pedestals"); c3.Divide(2,3); if (!fDataCheck) { disp0.CamDraw(c3, 1, 2, 1); disp1.CamDraw(c3, 2, 2, 6); return; } c3.cd(1); gPad->SetBorderMode(0); gPad->SetTicks(); MHCamera *obj1=(MHCamera*)disp0.DrawCopy("hist"); obj1->SetDirectory(NULL); // // for the datacheck, fix the ranges!! // obj1->SetMinimum(fgPedestalMin); obj1->SetMaximum(fgPedestalMax); // // set reference lines // DisplayReferenceLines(obj1,0); // // end reference lines // c3.cd(3); gPad->SetBorderMode(0); obj1->SetPrettyPalette(); obj1->Draw(); c3.cd(5); gPad->SetBorderMode(0); gPad->SetTicks(); obj1->DrawProjection(7); c3.cd(2); gPad->SetBorderMode(0); gPad->SetTicks(); MHCamera *obj2=(MHCamera*)disp1.DrawCopy("hist"); obj2->SetDirectory(NULL); obj2->SetMinimum(fgPedRmsMin); obj2->SetMaximum(fgPedRmsMax); // // set reference lines // DisplayReferenceLines(obj1,1); c3.cd(4); gPad->SetBorderMode(0); obj2->SetPrettyPalette(); obj2->Draw(); c3.cd(6); gPad->SetBorderMode(0); gPad->SetTicks(); TArrayI inner(1); inner[0] = 0; TArrayI outer(1); outer[0] = 1; if (geomcam.InheritsFrom("MGeomCamMagic")) { TArrayI s0(6); s0[0] = 6; s0[1] = 1; s0[2] = 2; s0[3] = 3; s0[4] = 4; s0[5] = 5; TArrayI s1(3); s1[0] = 6; s1[1] = 1; s1[2] = 2; TArrayI s2(3); s2[0] = 3; s2[1] = 4; s2[2] = 5; gPad->Clear(); TVirtualPad *pad = gPad; pad->Divide(2,1); TH1D *inout[2]; inout[0] = disp1.ProjectionS(s0, inner, "Inner"); inout[1] = disp1.ProjectionS(s0, outer, "Outer"); inout[0]->SetDirectory(NULL); inout[1]->SetDirectory(NULL); for (int i=0; i<2; i++) { TLegend *leg2 = new TLegend(0.6,0.2,0.9,0.55); leg2->SetHeader(inout[i]->GetName()); pad->cd(i+1); inout[i]->SetLineColor(kRed+i); inout[i]->SetBit(kCanDelete); inout[i]->Draw(); inout[i]->Fit("gaus","Q"); leg2->AddEntry(inout[i],inout[i]->GetName(),"l"); // // Display the outliers as dead and noisy pixels // DisplayOutliers(inout[i]); // // Display the two half of the camera separately // TH1D *half[2]; half[0] = disp1.ProjectionS(s1, i==0 ? inner : outer , "Sector 6-1-2"); half[1] = disp1.ProjectionS(s2, i==0 ? inner : outer , "Sector 3-4-5"); for (int j=0; j<2; j++) { half[j]->SetLineColor(kRed+i+2*j+1); half[j]->SetDirectory(NULL); half[j]->SetBit(kCanDelete); half[j]->Draw("same"); leg2->AddEntry(half[j],half[j]->GetName(),"l"); } leg2->Draw(); } } } void MJPedestal::DisplayReferenceLines(MHCamera *cam, const Int_t what) const { const Double_t x = cam->GetNbinsX(); TLine line; line.SetLineStyle(kDashed); line.SetLineWidth(3); line.SetLineColor(kBlue); TLine *l1 = line.DrawLine(0, what ? fgRefPedRmsGalactic : fgRefPedGalactic, x, what ? fgRefPedRmsGalactic : fgRefPedGalactic); line.SetLineColor(kYellow); TLine *l2 = line.DrawLine(0, what ? fgRefPedRmsExtraGalactic : fgRefPedExtraGalactic, x, what ? fgRefPedRmsExtraGalactic : fgRefPedExtraGalactic); line.SetLineColor(kMagenta); TLine *l3 = line.DrawLine(0, what ? fgRefPedRmsClosedLids : fgRefPedClosedLids, x, what ? fgRefPedRmsClosedLids : fgRefPedClosedLids); TLegend *leg = new TLegend(0.4,0.75,0.7,0.99); leg->SetBit(kCanDelete); leg->AddEntry(l1, "Galactic Source","l"); leg->AddEntry(l2, "Extra-Galactic Source","l"); leg->AddEntry(l3, "Closed Lids","l"); leg->Draw(); } void MJPedestal::DisplayOutliers(TH1D *hist) const { const Float_t mean = hist->GetFunction("gaus")->GetParameter(1); const Float_t lolim = mean - 3.5*hist->GetFunction("gaus")->GetParameter(2); const Float_t uplim = mean + 3.5*hist->GetFunction("gaus")->GetParameter(2); const Stat_t dead = hist->Integral(0,hist->FindBin(lolim)-1); const Stat_t noisy = hist->Integral(hist->FindBin(uplim)+1,hist->GetNbinsX()+1); TLatex deadtex; deadtex.SetTextSize(0.06); deadtex.DrawLatex(0.1,hist->GetBinContent(hist->GetMaximumBin())/1.1,Form("%3i dead pixels",(Int_t)dead)); TLatex noisytex; noisytex.SetTextSize(0.06); noisytex.DrawLatex(0.1,hist->GetBinContent(hist->GetMaximumBin())/1.2,Form("%3i noisy pixels",(Int_t)noisy)); } Bool_t MJPedestal::WriteResult() { if (fOutputPath.IsNull()) return kTRUE; const TString oname(GetOutputFile()); *fLog << inf << "Writing to file: " << oname << endl; TFile file(oname, "UPDATE"); if (fDisplay && fDisplay->Write()<=0) { *fLog << err << "Unable to write MStatusDisplay to " << oname << endl; return kFALSE; } if (fPedestalCam.Write()<=0) { *fLog << err << "Unable to write MPedestalCam to " << oname << endl; return kFALSE; } if (fBadPixels.Write()<=0) { *fLog << err << "Unable to write MBadPixelsCam to " << oname << endl; return kFALSE; } return kTRUE; } void MJPedestal::SetOutputPath(const char *path) { fOutputPath = path; if (fOutputPath.EndsWith("/")) fOutputPath = fOutputPath(0, fOutputPath.Length()-1); } void MJPedestal::SetEnv(const char *env) { if (fEnv) delete fEnv; fEnv = new TEnv(env); } Bool_t MJPedestal::Process() { if (!ReadPedestalCam()) return ProcessFile(); return kTRUE; } Bool_t MJPedestal::ProcessFile() { if (!fRuns && !fEnv) { *fLog << err << "Neither AddRuns was called nor SetEnv was used... abort." << endl; return kFALSE; } if (fRuns && fRuns->GetNumRuns() != fRuns->GetNumEntries()) { *fLog << err << "Number of files found doesn't match number of runs... abort." << endl; return kFALSE; } *fLog << inf; fLog->Separator(GetDescriptor()); *fLog << "Calculate MPedestalCam from Runs "; if (fRuns) *fLog << fRuns->GetRunsAsString() << endl; else *fLog << "in " << fEnv->GetName() << endl; *fLog << endl; MParList plist; MTaskList tlist; plist.AddToList(&tlist); MReadMarsFile read("Events"); MRawFileRead rawread(NULL); if (fDataCheck) { if (fRuns) rawread.AddFiles(*fRuns); tlist.AddToList(&rawread); } else { read.DisableAutoScheme(); if (fRuns) static_cast(read).AddFiles(*fRuns); tlist.AddToList(&read); } // Enable logging to file //*fLog.SetOutputFile(lname, kTRUE); // Setup Tasklist plist.AddToList(&fPedestalCam); plist.AddToList(&fBadPixels); MGeomApply geomapl; MBadPixelsMerge merge(&fBadPixels); MPedCalcPedRun pedcalc; if (fExtractor) { pedcalc.SetWindowSize((Int_t)fExtractor->GetNumHiGainSamples()); pedcalc.SetRange(fExtractor->GetHiGainFirst(), fExtractor->GetHiGainLast()); } else { *fLog << warn << GetDescriptor(); *fLog << ": No extractor has been chosen, take default number of FADC samples " << endl; } tlist.AddToList(&geomapl); tlist.AddToList(&merge); tlist.AddToList(&pedcalc); // // Execute the eventloop // MEvtLoop evtloop(fName); evtloop.SetParList(&plist); evtloop.SetDisplay(fDisplay); evtloop.SetLogStream(fLog); if (fEnv) evtloop.ReadEnv(*fEnv); // Execute first analysis if (!evtloop.Eventloop()) { *fLog << err << GetDescriptor() << ": Failed." << endl; return kFALSE; } tlist.PrintStatistics(); DisplayResult(plist); if (!WriteResult()) return kFALSE; *fLog << inf << GetDescriptor() << ": Done." << endl; return kTRUE; }