/* ======================================================================== *\ ! ! * ! * 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/DarkPatch4/2003_11_27/20031126_03030_P_DarkPatch4_E.root"; // const TString calfile = "/mnt/Data/rootdata/DarkPatch4/2004_01_17/20040116_*_C_DarkPatch2_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; MPedCalcPedRun pedcalc; MGeomCamMagic geomcam; MPedestalCam pedcam; tlist.AddToList(&read); tlist.AddToList(&geomapl); tlist.AddToList(&pedcalc); plist.AddToList(&pedcam); MHCamEvent hist("Pedestal"); hist.SetType(1); 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); MExtractedSignalCam sigcam; MCalibrationCam calcam; // // Get the previously created MPedestalCam into the new Parameter List // plist2.AddToList(&pedcam); plist2.AddToList(&sigcam); plist2.AddToList(&calcam); // // Get the MAGIC geometry // tlist2.AddToList(&geomapl); plist2.AddToList(&geomcam); // // Now setup the new tasks and tasklist for the calibration // --------------------------------------------------- // MReadMarsFile read2("Events", calname); read2.DisableAutoScheme(); MExtractSignal sigsig; MCalibrationCalc calcalc; // // In case, we want to exclude a pre-defined list of bad pixels: // (This is a preliminary feature) // // calcalc.ExcludePixelsFromAsciiFile("badpixels.dat"); // // 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 apply another fit function to the // blind pixel // MCalibrationBlindPix *bp = calcam.GetBlindPixel(); bp->ChangeFitFunc(MHCalibrationBlindPixel::kEPoisson5); 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]"); 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,cam,1,2,1); CamDraw(c1,disp3,cam,2,2,2); // Fit Probability TCanvas &c2 = d3->AddTab("Fit Prob."); c2.Divide(1,3); CamDraw(c2,disp5,cam,1,1,3); // Times TCanvas &c3 = d3->AddTab("Fitted Times"); c3.Divide(3,3); CamDraw(c3,disp6,cam,1,3,1); CamDraw(c3,disp7,cam,2,3,0); CamDraw(c3,disp8,cam,3,3,0); // Pedestals TCanvas &c4 = d3->AddTab("Pedestals"); c4.Divide(2,3); CamDraw(c4,disp9,cam,1,2,0); CamDraw(c4,disp10,cam,2,2,1); // Reduced Sigmas TCanvas &c5 = d3->AddTab("Reduced Sigmas"); c5.Divide(2,3); CamDraw(c5,disp11,cam,1,2,2); CamDraw(c5,disp16,cam,2,2,2); // F-Factor Method TCanvas &c6 = d3->AddTab("F-Factor Method"); c6.Divide(2,3); CamDraw(c6,disp12,cam,1,2,1); CamDraw(c6,disp13,cam,2,2,1); // Blind Pixel Method TCanvas &c7 = d3->AddTab("Blind Pixel Method"); c7.Divide(2, 3); CamDraw(c7,disp14,cam,1,2,9); CamDraw(c7,disp15,cam,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); Float_t he = gStyle->GetStatH(); Float_t wi = gStyle->GetStatH(); gStyle->SetStatH(0.4); gStyle->SetStatW(0.25); 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.5; const Double_t mean = obj2->GetMean(); const Double_t rms = obj2->GetRMS(); const Double_t width = max-min; switch (fit) { case 0: TF1 *sgaus = new TF1("sgaus","gaus(0)",min,max); 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: TF1 *dgaus = new TF1("dgaus","gaus(0)+gaus(3)",min,max); dgaus->SetParNames("A1","#mu1","#sigma1","A2","#mu2","#sigma2"); dgaus->SetParameters(integ/width,max-width/6.,width/4., integ/width,min+width/6.,width/4.); dgaus->SetParLimits(0,0,integ); dgaus->SetParLimits(1,min,max); dgaus->SetParLimits(2,0,width/3.); dgaus->SetParLimits(3,0,integ); dgaus->SetParLimits(4,min,max); dgaus->SetParLimits(5,0,width/3.); obj2->Fit("dgaus","QLR"); obj2->GetFunction("dgaus")->SetLineColor(kYellow); break; case 2: TF1 *tgaus = new TF1("tgaus","gaus(0)+gaus(3)+gaus(6)",min,max); tgaus->SetParNames("A1","#mu1","#sigma1","A2","#mu2","#sigma2","A3","#mu3","#sigma3"); tgaus->SetParameters(integ/width,max-width/6.,width/4., integ/width,min+width/6.,width/4., integ/width,min+width/6.,width/2.); tgaus->SetParLimits(0,0,integ); tgaus->SetParLimits(1,min,max); tgaus->SetParLimits(2,0,width/4.); tgaus->SetParLimits(3,0,integ); tgaus->SetParLimits(4,min,max); tgaus->SetParLimits(5,0,width/4.); tgaus->SetParLimits(6,0,integ); tgaus->SetParLimits(7,min,max); tgaus->SetParLimits(8,0,width/2.); obj2->Fit("tgaus","QLR"); 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; } gStyle->SetStatH(he); gStyle->SetStatW(wi); }