/* ======================================================================== *\ ! ! * ! * 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="/mnt/Data/rootdata/Miscellaneous/2003_12_19/20031218_03522_P_Park_E.root", TString calname="/mnt/Data/rootdata/Miscellaneous/2003_12_19/20031218_03527_C_Park_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")); tlist2.AddToList(&geomapl); // // Now setup the new tasks and tasklist for the calibration // --------------------------------------------------- // MReadMarsFile read2("Events", calname); read2.DisableAutoScheme(); MExtractSignal sigsig; MCalibrationCalc calcalc; MExtractedSignalCam sigcam; plist2.AddToList(&geomcam); plist2.AddToList(&sigcam); // // As long, as we don't have digital modules, // we have to set the color of the pulser LED by hand // calcalc.SetPulserColor(MCalibrationCalc::kECT1); tlist2.AddToList(&read2); tlist2.AddToList(&sigsig); tlist2.AddToList(&calcalc); // // 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(); MHCamEvent camevt; MHCamera disp1 (geomcam, "MCalibrationPix;Charge", "Fitted Mean Charges"); MHCamera disp3 (geomcam, "MCalibrationPix;SigmaCharge", "Sigma of Fitted Charges"); MHCamera disp5 (geomcam, "MCalibrationPix;ChargeProb", "Probability of Fit"); MHCamera disp6 (geomcam, "MCalibrationPix;Time", "Arrival Times"); MHCamera disp7 (geomcam, "MCalibrationPix;SigmaTime", "Sigma of Arrival Times"); MHCamera disp8 (geomcam, "MCalibrationPix;TimeChiSquare", "Chi Square of Time Fit"); MHCamera disp9 (geomcam, "MCalibrationPix;Ped", "Pedestals"); MHCamera disp10 (geomcam, "MCalibrationPix;PedRms", "Pedestal RMS"); MHCamera disp11 (geomcam, "MCalibrationPix;RSigma", "Reduced Sigmas"); MHCamera disp12 (geomcam, "MCalibrationPix;PheFFactorMethod", "Nr. of Phe's (F-Factor Method)"); MHCamera disp13 (geomcam, "MCalibrationPix;MeanConversionFFactorMethod", "Conversion Factor (F-Factor Method)"); MHCamera disp14 (geomcam, "MCalibrationPix;MeanPhotInsidePlexiglass", "Nr. of Photons (Blind Pixel Method)"); MHCamera disp15 (geomcam, "MCalibrationPix;MeanConversionBlindPixelMethod", "Conversion Factor (Blind Pixel Method)"); MHCamera disp16 (geomcam, "MCalibrationPix;RSigma/Charge", "Reduced Sigma per Charge"); disp1.SetCamContent(*cam, 0); disp1.SetCamError(*cam,1); disp3.SetCamContent(*cam, 2); disp3.SetCamError(*cam,3); disp5.SetCamContent(*cam, 4); disp6.SetCamContent(*cam, 5); disp6.SetCamError(*cam, 6); disp7.SetCamContent(*cam, 6); disp8.SetCamContent(*cam, 7); disp9.SetCamContent(*cam, 8); disp9.SetCamError(*cam, 9); disp10.SetCamContent(*cam, 9); disp11.SetCamContent(*cam, 10); disp12.SetCamContent(*cam, 11); disp12.SetCamError(*cam, 12); disp13.SetCamContent(*cam, 13); disp13.SetCamError(*cam, 14); disp14.SetCamContent(*cam, 15); disp15.SetCamContent(*cam, 16); disp16.SetCamContent(*cam, 17); disp1.SetYTitle("Charge [FADC counts]"); // disp2.SetYTitle("\\Delta Q [FADC counts]"); disp3.SetYTitle("\\sigma_{Charge} [FADC counts]"); // disp4.SetYTitle("\\Delta {\\sigma_{Q}} [FADC counts]"); disp5.SetYTitle("P_{Charge} [1]"); disp6.SetYTitle("Arr. Time [Time Slice Nr.]"); disp7.SetYTitle("\\sigma_{Time} [Time Slices]"); disp8.SetYTitle("\\chi^{2}_{Time} [1]"); disp9.SetYTitle("Ped [FADC Counts ]"); disp10.SetYTitle("RMS_{Ped} [FADC Counts ]"); disp11.SetYTitle("\\sqrt{\\sigma^{2}_{Charge} - RMS^{2}_{Ped}} [FADC Counts]"); disp12.SetYTitle("Nr. Photo-Electrons [1]"); disp13.SetYTitle("Conversion Factor [PhE/FADC Count]"); disp14.SetYTitle("Nr. Photons [1]"); disp15.SetYTitle("Conversion Factor [Phot/FADC Count]"); disp16.SetYTitle("Reduced Sigma / Charge [1]"); MStatusDisplay *d3 = new MStatusDisplay; d3->SetUpdateTime(3000); d3->Resize(1180,900); gStyle->SetOptStat(1111); gStyle->SetOptFit(); // Charges TCanvas *c1 = &d3->AddTab("Fitted Charges"); c1->Divide(2,3); TObject *obj1; TObject *obj2; c1->cd(1); obj1=disp1.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c1->cd(3); gPad->SetBorderMode(0); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c1->cd(5); gPad->SetBorderMode(0); obj2 = obj1->DrawClone("proj"); gPad->Update(); ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q"); ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow); c1->cd(2); obj1=disp3.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c1->cd(4); gPad->SetBorderMode(0); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c1->cd(6); gPad->SetBorderMode(0); obj2 = obj1->DrawClone("proj"); gPad->Update(); ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q"); ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow); // Fit Probability TCanvas *c12 = &d3->AddTab("Fit Prob."); c12->Divide(1, 3); c12->cd(1); obj1=disp5.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c12->cd(2); gPad->SetBorderMode(0); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c12->cd(3); gPad->SetBorderMode(0); obj2 = obj1->DrawClone("proj"); gPad->Update(); ((MHCamera*)obj2)->GetYProj()->Fit("pol0","Q"); ((MHCamera*)obj2)->GetYProj()->GetFunction("pol0")->SetLineColor(kYellow); // Times TCanvas *c2 = &d3->AddTab("Fitted Times"); c2->Divide(3, 3); c2->cd(1); obj1=disp6.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c2->cd(4); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c2->cd(7); gPad->SetBorderMode(0); obj2 = obj1->DrawClone("proj"); gPad->Update(); ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q"); ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow); c2->cd(2); obj1=disp7.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c2->cd(5); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c2->cd(8); gPad->SetBorderMode(0); obj2 = obj1->DrawClone("proj"); gPad->Update(); ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q"); ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow); c2->cd(3); obj1=disp8.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c2->cd(6); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c2->cd(9); gPad->SetBorderMode(0); obj2 = obj1->DrawClone("proj"); gPad->Update(); ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q"); ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow); // Pedestals TCanvas *c3 = &d3->AddTab("Pedestals"); c3->Divide(2, 3); c3->cd(1); obj1=disp9.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c3->cd(3); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c3->cd(5); gPad->SetBorderMode(0); obj2 = obj1->DrawClone("proj"); gPad->Update(); ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q"); ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow); c3->cd(2); obj1=disp10.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c3->cd(4); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c3->cd(6); gPad->SetBorderMode(0); obj2 = obj1->DrawClone("proj"); gPad->Update(); ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q"); ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow); // Reduced Sigmas TCanvas *c4 = &d3->AddTab("Reduced Sigmas"); c4->Divide(2,3); c4->cd(1); obj1=disp11.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c4->cd(3); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c4->cd(5); gPad->SetBorderMode(0); obj2 = obj1->DrawClone("proj"); gPad->Update(); ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q"); ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow); c4->cd(2); obj1=disp16.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c4->cd(4); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c4->cd(6); gPad->SetBorderMode(0); obj2 = obj1->DrawClone("proj"); gPad->Update(); ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q"); ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow); // F-Factor Method TCanvas *c5 = &d3->AddTab("F-Factor Method"); c5->Divide(2, 3); c5->cd(1); obj1=disp12.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c5->cd(3); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c5->cd(5); gPad->SetBorderMode(0); obj2 = obj1->DrawClone("proj"); gPad->Update(); ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q"); ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow); c5->cd(2); obj1=disp13.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c5->cd(4); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c5->cd(6); gPad->SetBorderMode(0); obj2 = obj1->DrawClone("proj"); gPad->Update(); ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q"); ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow); // Blind Pixel Method TCanvas *c6 = &d3->AddTab("Blind Pixel Method"); c6->Divide(2, 3); c6->cd(1); obj1=disp14.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c6->cd(3); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c6->cd(2); obj1=disp15.DrawCopy("hist"); ((MHCamera*)obj1)->AddNotify(*cam); c6->cd(4); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c6->cd(6); gPad->SetBorderMode(0); obj2 = obj1->DrawClone("proj"); gPad->Update(); ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q"); ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow); }