/* ======================================================================== *\ ! ! * ! * 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): Markus Gaug, 11/2003 ! ! Copyright: MAGIC Software Development, 2000-2003 ! ! \* ======================================================================== */ void calibration(TString pedname="../../Mars-0.8.2/20031102_02399_P_Unavailable_E.root", TString calname="../../Mars-0.8.2/20031102_02400_D_Flip500Hz_E.root") { // // Create a empty Parameter List and an empty Task List // The tasklist is identified in the eventloop by its name // MParList plist; MTaskList tlist; plist.AddToList(&tlist); // // Now setup the tasks and tasklist for the pedestals: // --------------------------------------------------- // MReadMarsFile read("Events", pedname); read.DisableAutoScheme(); MGeomApply geomapl; MPedCalcPedRun pedcalc; MGeomCamMagic geomcam; MPedestalCam pedcam; tlist.AddToList(&read); tlist.AddToList(&geomapl); tlist.AddToList(&pedcalc); plist.AddToList(&pedcam); MHCamEvent hist("Pedestal"); hist.SetType(0); plist.AddToList(&hist); MFillH fill(&hist, "MPedestalCam"); tlist.AddToList(&fill); // MStatusDisplay *d1 = new MStatusDisplay; // Set update time to 3s // d1->SetUpdateTime(3000); // // Create and setup the eventloop // MEvtLoop evtloop; evtloop.SetParList(&plist); // evtloop.SetDisplay(d1); // // Execute first analysis // if (!evtloop.Eventloop()) return; tlist.PrintStatistics(); // // Create a empty Parameter List and an empty Task List // MParList plist2; MTaskList tlist2; plist2.AddToList(&tlist2); plist2.AddToList((MPedestalCam*)plist.FindObject("MPedestalCam")); // MGeomApply geomapl2; tlist2.AddToList(&geomapl); // // Now setup the new tasks and tasklist for the calibration // --------------------------------------------------- // MReadMarsFile read2("Events", calname); read2.DisableAutoScheme(); MCalibrationCalc calcalc; // calcalc.SetSkipTFits(); plist2.AddToList(&geomcam); // // As long, as we don't have digital modules, // we have to set the color of the pulser LED by hand // calcalc.SetPulserColor(MCalibrationCalc::kEBlue); tlist2.AddToList(&read2); tlist2.AddToList(&calcalc); MHCamEvent hist2; hist2.SetType(8); plist2.AddToList(&hist2); MFillH fill2("MHCamEvent", "MCalibrationCam"); tlist2.AddToList(&fill2); // // Create and setup the eventloop // MEvtLoop evtloop2; evtloop2.SetParList(&plist2); // // Execute second analysis // if (!evtloop2.Eventloop()) return; tlist2.PrintStatistics(); // // just one example how to get the plots of individual pixels // MCalibrationCam *cam = plist2.FindObject("MCalibrationCam"); cam.Print(); Int_t pixnr; while (1) { cout << "Which pixel number do you want to display? (Press 0 to exit)" << endl; cin >> pixnr; if (pixnr == 0) break; if (pixnr >= 577) break; MCalibrationPix *pix = cam->GetCalibrationPix(pixnr); pix->Draw(); } // // Here we are confronted to a serious bug in ROOT: // If we do not apply the next command, gPad will get // screwed up completely: (Thanks to tbretz for finding out // the reason during several hours!!!) // gROOT->GetListOfCanvases()->Delete(); MHCamEvent &h = *(MHCamEvent*)plist2->FindObject("MHCamEvent"); MHCamera &disp0 = *h.GetHistByName(); MHCamera disp1 (geomcam, "MCalibrationCam;q", "Fitted Mean Charges"); MHCamera disp2 (geomcam, "MCalibrationCam;errq", "Error of Fitted Mean Charges"); MHCamera disp3 (geomcam, "MCalibrationCam;sigmaq", "Sigma of Fitted Mean Charges"); MHCamera disp4 (geomcam, "MCalibrationCam;errsigmaq", "Error of Sigma of Fitted Mean Charges"); MHCamera disp5 (geomcam, "MCalibrationCam;probq", "Probability of Fit"); MHCamera disp6 (geomcam, "MCalibrationCam;t", "Arrival Times"); MHCamera disp7 (geomcam, "MCalibrationCam;sigmat", "Sigma of Arrival Times"); MHCamera disp8 (geomcam, "MCalibrationCam;probt", "Probability of Time Fit"); MHCamera disp9 (geomcam, "MCalibrationCam;ped", "Pedestals"); MHCamera disp10 (geomcam, "MCalibrationCam;pedrms", "Pedestal RMS"); MHCamera disp11 (geomcam, "MCalibrationCam;rq", "Reduced Charges"); MHCamera disp12 (geomcam, "MCalibrationCam;rsigma", "Reduced Sigmas"); MHCamera disp13 (geomcam, "MCalibrationCam;phe", "Nr. of Phe's (F-Factor Method)"); MHCamera disp14 (geomcam, "MCalibrationCam;convphe", "Conversion Factor (F-Factor Method)"); MHCamera disp15 (geomcam, "MCalibrationCam;photons", "Nr. of Photons (Blind Pixel Method)"); MHCamera disp16 (geomcam, "MCalibrationCam;convphot", "Conversion Factor (Blind Pixel Method)"); disp1.SetCamContent(*cam, 0); disp2.SetCamContent(*cam, 1); disp3.SetCamContent(*cam, 2); disp4.SetCamContent(*cam, 3); disp5.SetCamContent(*cam, 4); disp6.SetCamContent(*cam, 5); disp7.SetCamContent(*cam, 6); disp8.SetCamContent(*cam, 7); disp9.SetCamContent(*cam, 8); disp10.SetCamContent(*cam, 9); disp11.SetCamContent(*cam, 10); disp12.SetCamContent(*cam, 11); disp13.SetCamContent(*cam, 12); disp14.SetCamContent(*cam, 13); disp15.SetCamContent(*cam, 14); disp16.SetCamContent(*cam, 15); disp1.SetYTitle("Q [FADC counts]"); disp2.SetYTitle("\\Delta Q [FADC counts]"); disp3.SetYTitle("\\sigma_{Q} [FADC counts]"); disp4.SetYTitle("\\Delta {\\sigma_{Q}} [FADC counts]"); disp5.SetYTitle("P [au]"); disp6.SetYTitle("T [FADC slices]"); disp7.SetYTitle("\\Delta T [FADC slices]"); disp8.SetYTitle("P [au]"); disp9.SetYTitle("P [FADC counts/ slice ]"); disp10.SetYTitle("RMS_{P} [FADC counts / slice ]"); disp11.SetYTitle("Q [FADC counts]"); disp12.SetYTitle("\\sigma^2_{Q} - RMS^2_{P} [FADC counts^2]"); disp13.SetYTitle("Nr Phe's"); disp14.SetYTitle("Conversion Factor [Phe/FADC count]"); disp15.SetYTitle("Nr Photons"); disp16.SetYTitle("Conversion Factor [Ph/FADC count]"); MStatusDisplay *d2 = new MStatusDisplay; // Set update time to 1s d2->SetUpdateTime(1000); TCanvas *c1 = &d2->AddTab("Charges Mean"); c1->Divide(2, 2); TObject *obj; c1->cd(1); gStyle->SetOptStat(1111); obj=disp1.DrawCopy("hist"); c1->cd(3); gPad->SetBorderMode(0); obj->Draw(); c1->cd(2); gStyle->SetOptStat(1101); obj=disp2.DrawCopy("hist"); c1->cd(4); gPad->SetBorderMode(0); obj->Draw(); TCanvas *c11 = &d2->AddTab("Charges Sigma"); c11->Divide(2, 2); c11->cd(1); gStyle->SetOptStat(1101); obj=disp3.DrawCopy("hist"); c11->cd(3); gPad->SetBorderMode(0); obj->Draw(); c11->cd(2); gStyle->SetOptStat(1101); obj=disp4.DrawCopy("hist"); c11->cd(4); gPad->SetBorderMode(0); obj->Draw(); TCanvas *c12 = &d2->AddTab("Fit Prob."); c12->Divide(1, 2); c12->cd(1); gStyle->SetOptStat(1101); obj=disp5.DrawCopy("hist"); c12->cd(2); gPad->SetBorderMode(0); obj->Draw(); TCanvas *c2 = &d2->AddTab("Fitted Times"); c2->Divide(3, 2); c2->cd(1); gStyle->SetOptStat(1111); obj=disp6.DrawCopy("hist"); c2->cd(4); obj->Draw(); c2->cd(2); gStyle->SetOptStat(1101); obj=disp7.DrawCopy("hist"); c2->cd(5); obj->Draw(); c2->cd(3); gStyle->SetOptStat(1101); obj=disp8.DrawCopy("hist"); c2->cd(6); obj->Draw(); TCanvas *c3 = &d2->AddTab("Pedestals"); c3->Divide(2, 2); c3->cd(1); gStyle->SetOptStat(1111); obj=disp9.DrawCopy("hist"); c3->cd(3); obj->Draw(); c3->cd(2); gStyle->SetOptStat(1111); obj=disp10.DrawCopy("hist"); c3->cd(4); obj->Draw(); TCanvas *c4 = &d2->AddTab("Reduced Charges"); c4->Divide(2, 2); c4->cd(1); gStyle->SetOptStat(1111); obj=disp11.DrawCopy("hist"); c4->cd(3); obj->Draw(); c4->cd(2); gStyle->SetOptStat(1101); obj=disp12.DrawCopy("hist"); c4->cd(4); obj->Draw(); TCanvas *c5 = &d2->AddTab("F-Factor Method"); c5->Divide(2, 2); c5->cd(1); gStyle->SetOptStat(1111); obj=disp13.DrawCopy("hist"); c5->cd(3); obj->Draw(); c5->cd(2); gStyle->SetOptStat(1101); obj=disp14.DrawCopy("hist"); c5->cd(4); obj->Draw(); TCanvas *c6 = &d2->AddTab("Blind Pixel Method"); c6->Divide(2, 2); c6->cd(1); gStyle->SetOptStat(1111); obj=disp15.DrawCopy("hist"); c6->cd(3); obj->Draw(); c6->cd(2); gStyle->SetOptStat(1101); obj=disp16.DrawCopy("hist"); c6->cd(4); obj->Draw(); #endif }