/* ======================================================================== *\ ! ! * ! * 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 ! ! \* ======================================================================== */ const TString pedfile = "/mnt/Data/rootdata/CrabNebula/2004_01_27/20040126_12149_P_Cab-On_E.root"; const TString calfile = "/mnt/Data/rootdata/CrabNebula/2004_01_27/20040126_1252*_C_Cab-On_E.root"; //const TString pedfile = "/mnt/Data/rootdata/Miscellaneous/2003_12_19/20031218_03522_P_Park_E.root"; //const TString calfile = "/mnt/Data/rootdata/Miscellaneous/2003_12_19/20031218_03527_C_Park_E.root"; void calibration(TString pedname=pedfile, TString calname=calfile) { // // 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; MExtractSignal sigcalc; sigcalc.SetRange(1,14,1,14); MPedCalcPedRun pedcalc; pedcalc.SetUseHists(); tlist.AddToList(&read); tlist.AddToList(&geomapl); tlist.AddToList(&sigcalc); tlist.AddToList(&pedcalc); MGeomCamMagic geomcam; MPedestalCam pedcam; plist.AddToList(&geomcam); plist.AddToList(&pedcam); // // Create and setup the eventloop // MEvtLoop evtloop; evtloop.SetParList(&plist); // // Execute first analysis // if (!evtloop.Eventloop()) return; tlist.PrintStatistics(); pedcam[17].DrawClone(); MStatusDisplay *d1 = new MStatusDisplay; d1->SetUpdateTime(3000); d1->Resize(850,700); MHCamera dispped0 (geomcam, "MPedestalPix;Pedestal", "Mean per Slice"); MHCamera dispped1 (geomcam, "MPedestalPix;PedestalErr", "Mean Error per Slice"); MHCamera dispped2 (geomcam, "MPedestalPix;PedestalRms", "RMS per Slice"); MHCamera dispped3 (geomcam, "MPedestalPix;PedestalRmsErr", "RMS Error per Slice"); MHCamera dispped4 (geomcam, "MPedestalPix;Mean", "Fitted Mean per Slice"); MHCamera dispped5 (geomcam, "MPedestalPix;MeanErr", "Fitted Error of Mean per Slice"); MHCamera dispped6 (geomcam, "MPedestalPix;Sigma", "Fitted Sigma per Slice"); MHCamera dispped7 (geomcam, "MPedestalPix;SigmaErr", "Fitted Error of Sigma per Slice"); MHCamera dispped8 (geomcam, "MPedestalPix;Prob", "Probability of Fit"); MHCamera dispped9 (geomcam, "MPedestalPix;DeltaPedestalMean", "Rel. Diff. Mean per Slice (Calc.-Fitte)"); MHCamera dispped11 (geomcam, "MPedestalPix;DeltaPedestalMeanError", "Rel. Diff. Mean Error per Slice (Calc.-Fitted)"); MHCamera dispped12 (geomcam, "MPedestalPix;DeltaRmsSigma", "Rel. Diff. RMS per Slice (Calc.-Fitted)"); MHCamera dispped14 (geomcam, "MPedestalPix;DeltaRmsSigmaError", "Rel. Diff. RMS Error per Slice (Calc.-Fitted)"); dispped0.SetCamContent(pedcam, 0); dispped0.SetCamError(pedcam, 1); dispped1.SetCamContent(pedcam, 1); dispped2.SetCamContent(pedcam, 2); dispped2.SetCamError(pedcam,3); dispped3.SetCamContent(pedcam, 3); dispped4.SetCamContent(pedcam, 4); dispped4.SetCamError(pedcam, 5); dispped5.SetCamContent(pedcam, 5); dispped6.SetCamContent(pedcam, 6); dispped6.SetCamError(pedcam, 7); dispped7.SetCamContent(pedcam, 7); dispped8.SetCamContent(pedcam, 8); dispped9.SetCamContent(pedcam, 9); dispped9.SetCamError(pedcam, 10); dispped11.SetCamContent(pedcam, 11); dispped12.SetCamContent(pedcam, 12); dispped12.SetCamError(pedcam, 13); dispped14.SetCamContent(pedcam, 14); dispped0.SetYTitle("Calc. Pedestal per slice [FADC counts]"); dispped1.SetYTitle("Calc. Pedestal Error per slice [FADC counts]"); dispped2.SetYTitle("Calc. Pedestal RMS per slice [FADC counts]"); dispped3.SetYTitle("Calc. Pedestal RMS Error per slice [FADC counts]"); dispped4.SetYTitle("Fitted Mean per slice [FADC counts]"); dispped5.SetYTitle("Error of Fitted Mean per slice [FADC counts]"); dispped6.SetYTitle("Fitted Sigma per slice [FADC counts]"); dispped7.SetYTitle("Error of Fitted Sigma per slice [FADC counts]"); dispped8.SetYTitle("Fit Probability [1]"); dispped9.SetYTitle("Rel. Diff. Pedestal Calc.-Fitted per slice [1]"); dispped11.SetYTitle("Rel. Diff. Pedestal Error Calc.-Fitted per slice [1]"); dispped12.SetYTitle("Rel. Diff. Pedestal RMS Calc.-Fitted per slice [1]"); dispped14.SetYTitle("Rel. Diff. Pedestal RMS Error Calc.-Fitted per slice [1]"); gStyle->SetOptStat(1111); gStyle->SetOptFit(); // Histogram values TCanvas &b1 = d1->AddTab("Direct Calculation"); b1.Divide(4,3); CamDraw(b1,dispped0,pedcam,1,4,0); CamDraw(b1,dispped1,pedcam,2,4,2); CamDraw(b1,dispped2,pedcam,3,4,2); CamDraw(b1,dispped3,pedcam,4,4,2); // Fitted values TCanvas &b2 = d1->AddTab("Fits"); b2.Divide(4,3); CamDraw(b2,dispped4,pedcam,1,4,0); CamDraw(b2,dispped5,pedcam,2,4,1); CamDraw(b2,dispped6,pedcam,3,4,0); CamDraw(b2,dispped7,pedcam,4,4,1); // Fits Probability TCanvas &b3 = d1->AddTab("Fit Probabilities"); b3.Divide(1,3); CamDraw(b3,dispped8,pedcam,1,1,3); // Differences TCanvas &c4 = d1->AddTab("Relative Difference Calculation-Fits"); c4.Divide(4,3); CamDraw(c4,dispped9,pedcam,1,4,1); CamDraw(c4,dispped11,pedcam,2,4,1); CamDraw(c4,dispped12,pedcam,3,4,1); CamDraw(c4,dispped14,pedcam,4,4,1); // // Create a empty Parameter List and an empty Task List // MParList plist2; MTaskList tlist2; plist2.AddToList(&tlist2); MExtractedSignalCam sigcam; MCalibrationCam calcam; // // Get the previously created MPedestalCam into the new Parameter List // plist2.AddToList(&geomcam); plist2.AddToList(&pedcam); plist2.AddToList(&sigcam); plist2.AddToList(&calcam); // // Get the MAGIC geometry // tlist2.AddToList(&geomapl); // // Now setup the new tasks and tasklist for the calibration // --------------------------------------------------- // MReadMarsFile read2("Events", calname); read2.DisableAutoScheme(); // MExtractSignal sigcalc; MArrivalTimeCalc timecalc; MCalibrationCalc calcalc; // // Making the step size a bit bigger, gives us // faster results // timecalc.SetStepSize(0.5); // // 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); // // In case, we want to exclude a pre-defined list of bad pixels: // (This is a preliminary feature) // //calcalc.ExcludePixelsFromAsciiFile("badpixels_all.dat"); // // In case, you want to skip the blind pixel method: // (NOT RECOMMENDED!!!) // // calcalc.SkipBlindPixelFit(); // // In case, you want to skip the cosmics rejection // (NOT RECOMMENDED!!!) // // calcalc.SkipCosmicsRejection(); // // In case, you want to skip the quality checks // (NOT RECOMMENDED!!!) // // calcalc.SkipQualityChecks(); // // In case, we want to apply another fit function to the // blind pixel // MCalibrationBlindPix *bp = calcam.GetBlindPixel(); // bp->ChangeFitFunc(MHCalibrationBlindPixel::kEPolya); // bp->ChangeFitFunc(MHCalibrationBlindPixel::kEPoisson4); tlist2.AddToList(&read2); tlist2.AddToList(&sigcalc); // // In case, you want to skip the somewhat lengthy calculation // of the arrival times using a spline, uncomment the next line // tlist2.AddToList(&timecalc); tlist2.AddToList(&calcalc); // // Create and setup the eventloop // MEvtLoop evtloop2; evtloop2.SetParList(&plist2); // // Execute second analysis // if (!evtloop2.Eventloop()) return; tlist2.PrintStatistics(); // // print the most important results of all pixels // calcam.Print(); // // just one example how to get the plots of individual pixels // calcam[17].DrawClone(); 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(calcam, 0); disp1.SetCamError(calcam,1); disp3.SetCamContent(calcam, 2); disp3.SetCamError(calcam,3); disp5.SetCamContent(calcam, 4); disp6.SetCamContent(calcam, 5); disp6.SetCamError(calcam, 6); disp7.SetCamContent(calcam, 6); disp8.SetCamContent(calcam, 7); disp9.SetCamContent(calcam, 8); disp9.SetCamError(calcam, 9); disp10.SetCamContent(calcam, 9); disp11.SetCamContent(calcam, 10); disp12.SetCamContent(calcam, 11); disp12.SetCamError(calcam, 12); disp13.SetCamContent(calcam, 13); disp13.SetCamError(calcam, 14); disp14.SetCamContent(calcam, 15); disp15.SetCamContent(calcam, 16); disp16.SetCamContent(calcam, 17); disp1.SetYTitle("Charge [FADC counts]"); disp3.SetYTitle("\\sigma_{Charge} [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(850,700); gStyle->SetOptStat(1111); gStyle->SetOptFit(); // Charges TCanvas &c1 = d3->AddTab("Fitted Charges"); c1.Divide(2,3); CamDraw(c1,disp1,calcam,1,2,1); CamDraw(c1,disp3,calcam,2,2,1); // Fit Probability TCanvas &c2 = d3->AddTab("Fit Prob."); c2.Divide(1,3); CamDraw(c2,disp5,calcam,1,1,3); // Times TCanvas &c3 = d3->AddTab("Fitted Times"); c3.Divide(3,3); CamDraw(c3,disp6,calcam,1,3,1); CamDraw(c3,disp7,calcam,2,3,0); CamDraw(c3,disp8,calcam,3,3,0); // Pedestals TCanvas &c4 = d3->AddTab("Pedestals"); c4.Divide(2,3); CamDraw(c4,disp9,calcam,1,2,0); CamDraw(c4,disp10,calcam,2,2,1); // Reduced Sigmas TCanvas &c5 = d3->AddTab("Reduced Sigmas"); c5.Divide(2,3); // CamDraw(c5,disp11,calcam,1,2,1); CamDraw(c5,disp11,calcam,1,2,2); CamDraw(c5,disp16,calcam,2,2,1); // F-Factor Method TCanvas &c6 = d3->AddTab("F-Factor Method"); c6.Divide(2,3); CamDraw(c6,disp12,calcam,1,2,1); CamDraw(c6,disp13,calcam,2,2,1); // Blind Pixel Method TCanvas &c7 = d3->AddTab("Blind Pixel Method"); c7.Divide(2, 3); CamDraw(c7,disp14,calcam,1,2,9); CamDraw(c7,disp15,calcam,2,2,1); } void CamDraw(TCanvas &c, MHCamera &cam, MCamEvent &evt, Int_t i, Int_t j, Int_t fit) { c.cd(i); gPad->SetBorderMode(0); MHCamera *obj1=(MHCamera*)cam.DrawCopy("hist"); obj1->AddNotify(evt); c.cd(i+j); gPad->SetBorderMode(0); obj1->Draw(); ((MHCamera*)obj1)->SetPrettyPalette(); c.cd(i+2*j); gPad->SetBorderMode(0); TH1D *obj2 = (TH1D*)obj1->Projection(); obj2->Draw(); obj2->SetBit(kCanDelete); const Double_t min = obj2->GetBinCenter(obj2->GetXaxis()->GetFirst()); const Double_t max = obj2->GetBinCenter(obj2->GetXaxis()->GetLast()); const Double_t integ = obj2->Integral("width")/2.5066283; const Double_t mean = obj2->GetMean(); const Double_t rms = obj2->GetRMS(); const Double_t width = max-min; if (rms == 0. || width == 0. ) return; switch (fit) { case 0: TF1 *sgaus = new TF1("sgaus","gaus(0)",min,max); sgaus->SetBit(kCanDelete); sgaus->SetParNames("Area","#mu","#sigma"); sgaus->SetParameters(integ/rms,mean,rms); sgaus->SetParLimits(0,0.,integ); sgaus->SetParLimits(1,min,max); sgaus->SetParLimits(2,0,width/1.5); obj2->Fit("sgaus","QLR"); obj2->GetFunction("sgaus")->SetLineColor(kYellow); break; case 1: TString dgausform = "([0]-[3])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"; dgausform += "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])"; TF1 *dgaus = new TF1("dgaus",dgausform.Data(),min,max); dgaus->SetBit(kCanDelete); dgaus->SetParNames("A_{tot}","#mu_{1}","#sigma_{1}","A_{2}","#mu_{2}","#sigma_{2}"); dgaus->SetParameters(integ,(min+mean)/2.,width/4., integ/width/2.,(max+mean)/2.,width/4.); // The left-sided Gauss dgaus->SetParLimits(0,integ-1.5,integ+1.5); dgaus->SetParLimits(1,min+(width/10.),mean); dgaus->SetParLimits(2,0,width/2.); // The right-sided Gauss dgaus->SetParLimits(3,0,integ); dgaus->SetParLimits(4,mean,max-(width/10.)); dgaus->SetParLimits(5,0,width/2.); obj2->Fit("dgaus","QLRM"); obj2->GetFunction("dgaus")->SetLineColor(kYellow); break; case 2: TString tgausform = "([0]-[3]-[6])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"; tgausform += "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])"; tgausform += "+[6]/[8]*exp(-0.5*(x-[7])*(x-[7])/[8]/[8])"; TF1 *tgaus = new TF1("tgaus",tgausform.Data(),min,max); tgaus->SetBit(kCanDelete); tgaus->SetParNames("A_{tot}","#mu_{1}","#sigma_{1}", "A_{2}","#mu_{2}","#sigma_{2}", "A_{3}","#mu_{3}","#sigma_{3}"); tgaus->SetParameters(integ,(min+mean)/2,width/4., integ/width/3.,(max+mean)/2.,width/4., integ/width/3.,mean,width/2.); // The left-sided Gauss tgaus->SetParLimits(0,integ-1.5,integ+1.5); tgaus->SetParLimits(1,min+(width/10.),mean); tgaus->SetParLimits(2,width/15.,width/2.); // The right-sided Gauss tgaus->SetParLimits(3,0.,integ); tgaus->SetParLimits(4,mean,max-(width/10.)); tgaus->SetParLimits(5,width/15.,width/2.); // The Gauss describing the outliers tgaus->SetParLimits(6,0.,integ); tgaus->SetParLimits(7,min,max); tgaus->SetParLimits(8,width/4.,width/1.5); obj2->Fit("tgaus","QLRM"); obj2->GetFunction("tgaus")->SetLineColor(kYellow); break; case 3: obj2->Fit("pol0","Q"); obj2->GetFunction("pol0")->SetLineColor(kYellow); break; case 9: break; default: obj2->Fit("gaus","Q"); obj2->GetFunction("gaus")->SetLineColor(kYellow); break; } }