/* ======================================================================== *\ ! ! * ! * 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, 02/2004 ! ! Copyright: MAGIC Software Development, 2000-2004 ! ! \* ======================================================================== */ ///////////////////////////////////////////////////////////////////////////// // // MJCalibration // // Do one calibration loop over serious of runs with the same pulser // colour and the same intensity. The following containers (rectangular) and // tasks (ellipses) are called to produce an MCalibrationChargeCam and to // update the MCalibrationQECam: (MCalibrate is not called from this class) // //Begin_Html /* */ //End_Html // // Different signal extractors can be chosen via the command SetExtractorLevel(UInt_t i) // Up to now, the following extractors are available: // i=1: Use MExtractSignal (fixed window) // i=2: Use MExtractSignal2 (sliding window: default) // i=3: Use MExtractSignal3 (coherent sliding window for all pixels) // // At the end of the eventloop, plots and results are displayed, depending on // the flags set (see DisplayResult()) // // If the flag SetFullDisplay() is set, all MHCameras will be displayed. // if the flag SetDataCheck() is set, only the most important ones are displayed // Otherwise, (default: SetNormalDisplay()), a good selection of plots is given // // See also: MHCalibrationChargePix, MHCalibrationChargeCam, MHGausEvents // MHCalibrationChargeBlindPix, MHCalibrationChargePINDiode // MCalibrationChargePix, MCalibrationChargeCam, MCalibrationChargeCalc // MCalibrationChargeBlindPix, MCalibrationChargePINDiode, // MCalibrationQECam, MBadPixelsPix, MBadPixelsCam // // If the flag RelTimeCalibration() is set, a calibration of the relative arrival // times is also performed. The following containers (rectangular) and // tasks (ellipses) are called to produce an MCalibrationRelTimeCam used by // MCalibrateTime to correct timing offset between pixels: (MCalibrateTime is not // called from this class) // //Begin_Html /* */ //End_Html // // Different arrival time extractors can be chosen via the command SetArrivalTimeLevel(UInt_t i) // Up to now, the following extractors are available: // i=1: Use MArrivalTimeCalc (using MCubicSpline, arrival time == position at half maximum) // i=2: Use MArrivalTimeCalc2 (mean time of fWindowSize time slices with the highest integral content: default) // // See also: MHCalibrationRelTimePix, MHCalibrationRelTimeCam, MHGausEvents // MCalibrationRelTimePix, MCalibrationRelTimeCam // MBadPixelsPix, MBadPixelsCam // ///////////////////////////////////////////////////////////////////////////// #include "MJCalibration.h" #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 "MHCamera.h" #include "MGeomCam.h" #include "MPedestalCam.h" #include "MCalibrationCam.h" #include "MCalibrationQECam.h" #include "MCalibrationChargeCam.h" #include "MCalibrationChargePINDiode.h" #include "MCalibrationChargeBlindPix.h" #include "MCalibrationChargeCalc.h" #include "MHGausEvents.h" #include "MHCalibrationCam.h" #include "MHCalibrationChargeCam.h" #include "MHCalibrationRelTimeCam.h" #include "MCalibrationRelTimeCam.h" #include "MReadMarsFile.h" #include "MGeomApply.h" #include "MBadPixelsMerge.h" #include "MBadPixelsCam.h" #include "MExtractSignal.h" #include "MExtractPINDiode.h" #include "MExtractBlindPixel.h" #include "MExtractSignal2.h" #include "MExtractSignal3.h" #include "MFCosmics.h" #include "MContinue.h" #include "MFillH.h" #include "MArrivalTimeCam.h" #include "MArrivalTimeCalc.h" #include "MArrivalTimeCalc2.h" #include "MJCalibration.h" #include "MStatusDisplay.h" ClassImp(MJCalibration); using namespace std; // -------------------------------------------------------------------------- // // Default constructor. // // Sets fRuns to 0, fColor to kNONE, fDisplay to kNormalDisplay, // fRelTime to kFALSE, fExtractorLevel to 2, fArrivalTimeLevel to 2 // MJCalibration::MJCalibration(const char *name, const char *title) : fRuns(0), fColor(MCalibrationCam::kNONE), fDisplayType(kNormalDisplay), fRelTimes(kFALSE), fExtractorLevel(2), fArrivalTimeLevel(2) { fName = name ? name : "MJCalibration"; fTitle = title ? title : "Tool to create a pedestal file (MPedestalCam)"; } // -------------------------------------------------------------------------- // // Draw the MHCamera into the MStatusDisplay: // // 1) Draw it as histogram (MHCamera::DrawCopy("hist") // 2) Draw it as a camera, with MHCamera::SetPrettyPalette() set. // 3) If "rad" is not zero, draw its values vs. the radius from the camera center. // (DrawRadialProfile()) // 4) Depending on the variable "fit", draw the values projection on the y-axis // (DrawProjection()): // 0: don't draw // 1: Draw fit to Single Gauss (for distributions flat-fielded over the whole camera) // 2: Draw and fit to Double Gauss (for distributions different for inner and outer pixels) // 3: Draw and fit to Triple Gauss (for distributions with inner, outer pixels and outliers) // 4: Draw and fit to Polynomial grade 0: (for the probability distributions) // >4: Draw and don;t fit. // void MJCalibration::CamDraw(TCanvas &c, const Int_t x, const Int_t y, const MHCamera &cam1, const Int_t fit, const Int_t rad) { c.cd(x); gPad->SetBorderMode(0); gPad->SetTicks(); MHCamera *obj1=(MHCamera*)cam1.DrawCopy("hist"); obj1->SetDirectory(NULL); c.cd(x+y); gPad->SetBorderMode(0); obj1->SetPrettyPalette(); obj1->AddNotify(&fCalibrationCam); obj1->Draw(); if (rad) { c.cd(x+2*y); gPad->SetBorderMode(0); gPad->SetTicks(); DrawRadialProfile(obj1); } if (!fit) return; c.cd(rad ? x+3*y : x+2*y); gPad->SetBorderMode(0); gPad->SetTicks(); DrawProjection(obj1, fit); } // -------------------------------------------------------------------------- // // Display the results in MStatusDisplay: // // - Add "Calibration" to the MStatusDisplay title // - Retrieve the MGeomCam from MParList // - Initialize the following MHCamera's: // 1) MCalibrationPix::GetMean() // 2) MCalibrationPix::Sigma() // 3) MCalibrationChargePix::GetRSigma() // 4) MCalibrationChargePix::GetRSigmaPerCharge() // 5) MCalibrationChargePix::GetPheFFactorMethod() // 6) MCalibrationChargePix::GetMeanConvFADC2Phe() // 7) MCalibrationChargePix::GetMeanFFactorFADC2Phot() // 8) MCalibrationQEPix::GetQECascadesFFactor() // 9) MCalibrationQEPix::GetQECascadesBlindPixel() // 10) MCalibrationQEPix::GetQECascadesPINDiode() // 11) MCalibrationQEPix::GetQECascadesCombined() // 12) MCalibrationQEPix::IsAverageQEFFactorAvailable() // 13) MCalibrationQEPix::IsAverageQEBlindPixelAvailable() // 14) MCalibrationQEPix::IsAverageQEPINDiodeAvailable() // 15) MCalibrationQEPix::IsAverageQECombinedAvailable() // 16) MCalibrationChargePix::IsHiGainSaturation() // 17) MCalibrationPix::GetHiLoMeansDivided() // 18) MCalibrationPix::GetHiLoSigmasDivided() // 19) MCalibrationChargePix::GetHiGainPickup() // 20) MCalibrationChargePix::GetLoGainPickup() // 21) MCalibrationChargePix::GetHiGainBlackout() // 22) MCalibrationChargePix::GetLoGainBlackout() // 23) MCalibrationPix::IsExcluded() // 24) MBadPixelsPix::IsUnsuitable(MBadPixelsPix::kUnsuitableRun) // 25) MBadPixelsPix::IsUnsuitable(MBadPixelsPix::kUnreliableRun) // 26) MBadPixelsPix::IsUncalibrated(MBadPixelsPix::kHiGainOscillating) // 27) MBadPixelsPix::IsUncalibrated(MBadPixelsPix::kLoGainOscillating) // 28) MCalibrationChargePix::GetAbsTimeMean() // 29) MCalibrationChargePix::GetAbsTimeRms() // // If the flag SetFullDisplay() is set, all MHCameras will be displayed. // if the flag SetDataCheck() is set, only the most important ones are displayed // and otherwise, (default: SetNormalDisplay()), a good selection of plots is given // void MJCalibration::DisplayResult(MParList &plist) { if (!fDisplay) return; // // Update display // TString title = fDisplay->GetTitle(); title += "-- Calibration "; title += fRuns->GetRunsAsString(); title += " --"; fDisplay->SetTitle(title); // // Get container from list // MGeomCam &geomcam = *(MGeomCam*)plist.FindObject("MGeomCam"); // Create histograms to display MHCamera disp1 (geomcam, "Cal;Charge", "Fitted Mean Charges"); MHCamera disp2 (geomcam, "Cal;SigmaCharge", "Sigma of Fitted Charges"); MHCamera disp3 (geomcam, "Cal;RSigma", "Reduced Sigmas"); MHCamera disp4 (geomcam, "Cal;RSigma/Charge", "Reduced Sigma per Charge"); MHCamera disp5 (geomcam, "Cal;FFactorPhe", "Nr. of Phe's (F-Factor Method)"); MHCamera disp6 (geomcam, "Cal;FFactorConv", "Conversion Factor (F-Factor Method)"); MHCamera disp7 (geomcam, "Cal;FFactorFFactor", "Total F-Factor (F-Factor Method)"); MHCamera disp8 (geomcam, "Cal;CascadesQEFFactor", "Cascades QE (F-Factor Method)"); MHCamera disp9 (geomcam, "Cal;CascadesQEBlindPix","Cascades QE (Blind Pixel Method)"); MHCamera disp10(geomcam, "Cal;CascadesQEPINDiode","Cascades QE (PIN Diode Method)"); MHCamera disp11(geomcam, "Cal;CascadesQECombined","Cascades QE (Combined Method)"); MHCamera disp12(geomcam, "Cal;FFactorValid", "Pixels with valid F-Factor calibration"); MHCamera disp13(geomcam, "Cal;BlindPixelValid", "Pixels with valid BlindPixel calibration"); MHCamera disp14(geomcam, "Cal;PINdiodeValid", "Pixels with valid PINDiode calibration"); MHCamera disp15(geomcam, "Cal;CombinedValid", "Pixels with valid Combined calibration"); MHCamera disp16(geomcam, "Cal;Saturation", "Pixels with saturated Hi Gain"); MHCamera disp17(geomcam, "Cal;ConversionMeans", "Conversion HiGain.vs.LoGain Means"); MHCamera disp18(geomcam, "Cal;ConversionSigmas", "Conversion HiGain.vs.LoGain Sigmas"); MHCamera disp19(geomcam, "Cal;HiGainPickup", "Number Pickup events Hi Gain"); MHCamera disp20(geomcam, "Cal;LoGainPickup", "Number Pickup events Lo Gain"); MHCamera disp21(geomcam, "Cal;HiGainBlackout", "Number Blackout events Hi Gain"); MHCamera disp22(geomcam, "Cal;LoGainBlackout", "Number Blackout events Lo Gain"); MHCamera disp23(geomcam, "Cal;Excluded", "Pixels previously excluded"); MHCamera disp24(geomcam, "Bad;UnSuitable", "Pixels not suited for further analysis"); MHCamera disp25(geomcam, "Bad;UnReliable", "Pixels not reliable for further analysis"); MHCamera disp26(geomcam, "Bad;HiGainOscillating", "Oscillating Pixels High Gain"); MHCamera disp27(geomcam, "Bad;LoGainOscillating", "Oscillating Pixels Low Gain"); MHCamera disp28(geomcam, "Cal;AbsTimeMean", "Abs. Arrival Times"); MHCamera disp29(geomcam, "Cal;AbsTimeRms", "RMS of Arrival Times"); MHCamera disp30(geomcam, "time;MeanTime", "Mean Rel. Arrival Times"); MHCamera disp31(geomcam, "time;SigmaTime", "Sigma Rel. Arrival Times"); MHCamera disp32(geomcam, "time;TimeProb", "Probability of Time Fit"); MHCamera disp33(geomcam, "time;NotFitValid", "Pixels with not valid fit results"); MHCamera disp34(geomcam, "time;Oscillating", "Oscillating Pixels"); // Fitted charge means and sigmas disp1.SetCamContent(fCalibrationCam, 0); disp1.SetCamError( fCalibrationCam, 1); disp2.SetCamContent(fCalibrationCam, 2); disp2.SetCamError( fCalibrationCam, 3); // Reduced Sigmas and reduced sigmas per charge disp3.SetCamContent(fCalibrationCam, 5); disp3.SetCamError( fCalibrationCam, 6); disp4.SetCamContent(fCalibrationCam, 7); disp4.SetCamError( fCalibrationCam, 8); // F-Factor Method disp5.SetCamContent(fCalibrationCam, 9); disp5.SetCamError( fCalibrationCam, 10); disp6.SetCamContent(fCalibrationCam, 11); disp6.SetCamError( fCalibrationCam, 12); disp7.SetCamContent(fCalibrationCam, 13); disp7.SetCamError( fCalibrationCam, 14); // Quantum Efficiencies disp8.SetCamContent (fQECam, 0 ); disp8.SetCamError (fQECam, 1 ); disp9.SetCamContent (fQECam, 2 ); disp9.SetCamError (fQECam, 3 ); disp10.SetCamContent(fQECam, 4 ); disp10.SetCamError (fQECam, 5 ); disp11.SetCamContent(fQECam, 6 ); disp11.SetCamError (fQECam, 7 ); // Valid flags disp12.SetCamContent(fQECam, 8 ); disp13.SetCamContent(fQECam, 9 ); disp14.SetCamContent(fQECam, 10); disp15.SetCamContent(fQECam, 11); // Conversion Hi-Lo disp16.SetCamContent(fCalibrationCam, 25); disp17.SetCamContent(fCalibrationCam, 16); disp17.SetCamError (fCalibrationCam, 17); disp18.SetCamContent(fCalibrationCam, 18); disp18.SetCamError (fCalibrationCam, 19); // Pickup and Blackout disp19.SetCamContent(fCalibrationCam, 21); disp20.SetCamContent(fCalibrationCam, 22); disp21.SetCamContent(fCalibrationCam, 23); disp22.SetCamContent(fCalibrationCam, 24); // Pixels with defects disp23.SetCamContent(fCalibrationCam, 20); disp24.SetCamContent(fBadPixels, 1); disp25.SetCamContent(fBadPixels, 3); // Oscillations disp26.SetCamContent(fBadPixels, 10); disp27.SetCamContent(fBadPixels, 11); // Arrival Times disp28.SetCamContent(fCalibrationCam, 26); disp28.SetCamError( fCalibrationCam, 27); disp29.SetCamContent(fCalibrationCam, 27); disp1.SetYTitle("Q [FADC counts]"); disp2.SetYTitle("\\sigma_{Q} [FADC counts]"); disp3.SetYTitle("\\sqrt{\\sigma^{2}_{Q} - RMS^{2}_{Ped}} [FADC Counts]"); disp4.SetYTitle("Red.Sigma/ [1]"); disp5.SetYTitle("Nr. Phe's [1]"); disp6.SetYTitle("Conv.Factor [PhE/FADC counts]"); disp7.SetYTitle("Total F-Factor [1]"); disp8.SetYTitle("QE [1]"); disp9.SetYTitle("QE [1]"); disp10.SetYTitle("QE [1]"); disp11.SetYTitle("QE [1]"); disp12.SetYTitle("[1]"); disp13.SetYTitle("[1]"); disp14.SetYTitle("[1]"); disp15.SetYTitle("[1]"); disp16.SetYTitle("[1]"); disp17.SetYTitle("(High)/(Low) [1]"); disp18.SetYTitle("\\sigma_{Q}(High)/\\sigma_{Q}(Low) [1]"); disp19.SetYTitle("[1]"); disp20.SetYTitle("[1]"); disp21.SetYTitle("[1]"); disp22.SetYTitle("[1]"); disp23.SetYTitle("[1]"); disp24.SetYTitle("[1]"); disp25.SetYTitle("[1]"); disp26.SetYTitle("[1]"); disp27.SetYTitle("[1]"); disp28.SetYTitle("Mean Abs. Time [FADC slice]"); disp29.SetYTitle("RMS Abs. Time [FADC slices]"); if (fRelTimes) { disp30.SetCamContent(fRelTimeCam,0); disp30.SetCamError( fRelTimeCam,1); disp31.SetCamContent(fRelTimeCam,2); disp31.SetCamError( fRelTimeCam,3); disp32.SetCamContent(fRelTimeCam,4); disp33.SetCamContent(fBadPixels,20); disp34.SetCamContent(fBadPixels,21); disp30.SetYTitle("Time Offset [ns]"); disp31.SetYTitle("Timing resolution [ns]"); disp32.SetYTitle("P_{Time} [1]"); disp33.SetYTitle("[1]"); disp34.SetYTitle("[1]"); } gStyle->SetOptStat(1111); gStyle->SetOptFit(); if (fDisplayType == kDataCheck) { TCanvas &c1 = fDisplay->AddTab("Fit.Charge"); c1.Divide(3, 3); CamDraw(c1, 1, 3, disp1, 2); CamDraw(c1, 2, 3, disp4, 2); CamDraw(c1, 3, 3, disp28, 2); // F-Factor TCanvas &c2 = fDisplay->AddTab("Phe's"); c2.Divide(3,4); CamDraw(c2, 1, 3, disp6, 2, 1); CamDraw(c2, 2, 3, disp7, 2, 1); CamDraw(c2, 3, 3, disp8, 2, 1); // QE's TCanvas &c3 = fDisplay->AddTab("QE's"); c3.Divide(3,4); CamDraw(c3, 1, 3, disp8, 2, 1); CamDraw(c3, 2, 3, disp9, 2, 1); CamDraw(c3, 3, 3, disp10, 2, 1); // Defects TCanvas &c4 = fDisplay->AddTab("Defect"); c4.Divide(3,2); CamDraw(c4, 1, 3, disp23, 0); CamDraw(c4, 2, 3, disp24, 0); CamDraw(c4, 3, 3, disp25, 0); if (fRelTimes) { // Rel. Times TCanvas &c5 = fDisplay->AddTab("Rel. Times"); c5.Divide(2,4); CamDraw(c5, 1, 2, disp30, 2); CamDraw(c5, 2, 2, disp31, 2); } return; } if (fDisplayType == kNormalDisplay) { // Charges TCanvas &c11 = fDisplay->AddTab("Fit.Charge"); c11.Divide(2, 4); CamDraw(c11, 1, 2, disp1, 2, 1); CamDraw(c11, 2, 2, disp2, 2, 1); // Reduced Sigmas TCanvas &c12 = fDisplay->AddTab("Red.Sigma"); c12.Divide(2,4); CamDraw(c12, 1, 2, disp3, 2, 1); CamDraw(c12, 2, 2, disp4, 2, 1); // F-Factor TCanvas &c13 = fDisplay->AddTab("Phe's"); c13.Divide(3,4); CamDraw(c13, 1, 3, disp5, 2, 1); CamDraw(c13, 2, 3, disp6, 2, 1); CamDraw(c13, 3, 3, disp7, 2, 1); // QE's TCanvas &c14 = fDisplay->AddTab("QE's"); c14.Divide(4,4); CamDraw(c14, 1, 4, disp8, 2, 1); CamDraw(c14, 2, 4, disp9, 2, 1); CamDraw(c14, 3, 4, disp10, 2, 1); CamDraw(c14, 4, 4, disp11, 2, 1); // Defects TCanvas &c15 = fDisplay->AddTab("Defect"); c15.Divide(5,2); CamDraw(c15, 1, 5, disp23, 0); CamDraw(c15, 2, 5, disp24, 0); CamDraw(c15, 3, 5, disp25, 0); CamDraw(c15, 4, 5, disp26, 0); CamDraw(c15, 5, 5, disp27, 0); // Abs. Times TCanvas &c16 = fDisplay->AddTab("Abs. Times"); c16.Divide(2,3); CamDraw(c16, 1, 2, disp28, 2); CamDraw(c16, 2, 2, disp29, 1); if (fRelTimes) { // Rel. Times TCanvas &c17 = fDisplay->AddTab("Rel. Times"); c17.Divide(2,4); CamDraw(c17, 1, 2, disp30, 2, 1); CamDraw(c17, 2, 2, disp31, 2, 1); } return; } if (fDisplayType == kFullDisplay) { MHCalibrationCam *cam = (MHCalibrationCam*)plist.FindObject("MHCalibrationChargeCam"); for (Int_t aidx=0;aidxGetAverageAreas();aidx++) { cam->GetAverageHiGainArea(aidx).DrawClone("all"); cam->GetAverageLoGainArea(aidx).DrawClone("all"); } for (Int_t sector=1;sectorGetAverageSectors();sector++) { cam->GetAverageHiGainSector(sector).DrawClone("all"); cam->GetAverageLoGainSector(sector).DrawClone("all"); } // Charges TCanvas &c21 = fDisplay->AddTab("Fit.Charge"); c21.Divide(2, 4); CamDraw(c21, 1, 2, disp1, 2, 1); CamDraw(c21, 2, 2, disp2, 2, 1); // Reduced Sigmas TCanvas &c23 = fDisplay->AddTab("Red.Sigma"); c23.Divide(2,4); CamDraw(c23, 1, 2, disp3, 2, 1); CamDraw(c23, 2, 2, disp4, 2, 1); // F-Factor TCanvas &c24 = fDisplay->AddTab("Phe's"); c24.Divide(3,4); CamDraw(c24, 1, 3, disp5, 2, 1); CamDraw(c24, 2, 3, disp6, 2, 1); CamDraw(c24, 3, 3, disp7, 2, 1); // QE's TCanvas &c25 = fDisplay->AddTab("QE's"); c25.Divide(4,4); CamDraw(c25, 1, 4, disp8, 2, 1); CamDraw(c25, 2, 4, disp9, 2, 1); CamDraw(c25, 3, 4, disp10, 2, 1); CamDraw(c25, 4, 4, disp11, 2, 1); // Validity TCanvas &c26 = fDisplay->AddTab("Valid"); c26.Divide(4,2); CamDraw(c26, 1, 4, disp12, 0); CamDraw(c26, 2, 4, disp13, 0); CamDraw(c26, 3, 4, disp14, 0); CamDraw(c26, 4, 4, disp15, 0); // Other info TCanvas &c27 = fDisplay->AddTab("HiLoGain"); c27.Divide(3,3); CamDraw(c27, 1, 3, disp16, 0); CamDraw(c27, 2, 3, disp17, 1); CamDraw(c27, 3, 3, disp18, 1); // Pickup TCanvas &c28 = fDisplay->AddTab("Pickup"); c28.Divide(4,2); CamDraw(c28, 1, 4, disp19, 0); CamDraw(c28, 2, 4, disp20, 0); CamDraw(c28, 3, 4, disp21, 0); CamDraw(c28, 4, 4, disp22, 0); // Defects TCanvas &c29 = fDisplay->AddTab("Defect"); c29.Divide(5,2); CamDraw(c29, 1, 5, disp23, 0); CamDraw(c29, 2, 5, disp24, 0); CamDraw(c29, 3, 5, disp25, 0); CamDraw(c29, 4, 5, disp26, 0); CamDraw(c29, 5, 5, disp27, 0); // Abs. Times TCanvas &c30 = fDisplay->AddTab("Abs. Times"); c30.Divide(2,3); CamDraw(c30, 1, 2, disp28, 2); CamDraw(c30, 2, 2, disp29, 1); if (fRelTimes) { // Rel. Times TCanvas &c31 = fDisplay->AddTab("Rel. Times"); c31.Divide(3,4); CamDraw(c31, 1, 3, disp30, 2, 1); CamDraw(c31, 2, 3, disp31, 2, 1); CamDraw(c31, 3, 3, disp32, 4, 1); // Time Defects TCanvas &c32 = fDisplay->AddTab("Time Def."); c32.Divide(2,2); CamDraw(c32, 1, 2, disp33,0); CamDraw(c32, 2, 2, disp34,0); MHCalibrationCam *cam = (MHCalibrationCam*)plist.FindObject("MHCalibrationRelTimeCam"); for (Int_t aidx=0;aidxGetAverageAreas();aidx++) { cam->GetAverageHiGainArea(aidx).DrawClone("fourierevents"); cam->GetAverageLoGainArea(aidx).DrawClone("fourierevents"); } for (Int_t sector=1;sectorGetAverageSectors();sector++) { cam->GetAverageHiGainSector(sector).DrawClone("fourierevents"); cam->GetAverageLoGainSector(sector).DrawClone("fourierevents"); } } return; } } // -------------------------------------------------------------------------- // // Draw a projection of MHCamera onto the y-axis values. Depending on the // variable fit, the following fits are performed: // // 1: Single Gauss (for distributions flat-fielded over the whole camera) // 2: Double Gauss (for distributions different for inner and outer pixels) // 3: Triple Gauss (for distributions with inner, outer pixels and outliers) // 4: flat (for the probability distributions) // // Moreover, sectors 6,1 and 2 of the camera and sectors 3,4 and 5 are // drawn separately, for inner and outer pixels. // void MJCalibration::DrawProjection(MHCamera *obj, Int_t fit) const { TH1D *obj2 = (TH1D*)obj->Projection(obj->GetName()); obj2->SetDirectory(0); obj2->Draw(); obj2->SetBit(kCanDelete); if (obj->GetGeomCam().InheritsFrom("MGeomCamMagic")) { TArrayI s0(3); s0[0] = 6; s0[1] = 1; s0[2] = 2; TArrayI s1(3); s1[0] = 3; s1[1] = 4; s1[2] = 5; TArrayI inner(1); inner[0] = 0; TArrayI outer(1); outer[0] = 1; // Just to get the right (maximum) binning TH1D *half[4]; half[0] = obj->ProjectionS(s0, inner, "Sector 6-1-2 Inner"); half[1] = obj->ProjectionS(s1, inner, "Sector 3-4-5 Inner"); half[2] = obj->ProjectionS(s0, outer, "Sector 6-1-2 Outer"); half[3] = obj->ProjectionS(s1, outer, "Sector 3-4-5 Outer"); for (int i=0; i<4; i++) { half[i]->SetLineColor(kRed+i); half[i]->SetDirectory(0); half[i]->SetBit(kCanDelete); half[i]->Draw("same"); } } 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; const TString dgausformula = "([0]-[3])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])" "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])"; const TString tgausformula = "([0]-[3]-[6])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])" "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])" "+[6]/[8]*exp(-0.5*(x-[7])*(x-[7])/[8]/[8])"; TF1 *f=0; switch (fit) { case 1: f = new TF1("sgaus", "gaus(0)", min, max); f->SetLineColor(kYellow); f->SetBit(kCanDelete); f->SetParNames("Area", "#mu", "#sigma"); f->SetParameters(integ/rms, mean, rms); f->SetParLimits(0, 0, integ); f->SetParLimits(1, min, max); f->SetParLimits(2, 0, width/1.5); obj2->Fit(f, "QLR"); break; case 2: f = new TF1("dgaus",dgausformula.Data(),min,max); f->SetLineColor(kYellow); f->SetBit(kCanDelete); f->SetParNames("A_{tot}", "#mu1", "#sigma1", "A2", "#mu2", "#sigma2"); f->SetParameters(integ,(min+mean)/2.,width/4., integ/width/2.,(max+mean)/2.,width/4.); // The left-sided Gauss f->SetParLimits(0,integ-1.5 , integ+1.5); f->SetParLimits(1,min+(width/10.), mean); f->SetParLimits(2,0 , width/2.); // The right-sided Gauss f->SetParLimits(3,0 , integ); f->SetParLimits(4,mean, max-(width/10.)); f->SetParLimits(5,0 , width/2.); obj2->Fit(f,"QLRM"); break; case 3: f = new TF1("tgaus",tgausformula.Data(),min,max); f->SetLineColor(kYellow); f->SetBit(kCanDelete); f->SetParNames("A_{tot}","#mu_{1}","#sigma_{1}", "A_{2}","#mu_{2}","#sigma_{2}", "A_{3}","#mu_{3}","#sigma_{3}"); f->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 f->SetParLimits(0,integ-1.5,integ+1.5); f->SetParLimits(1,min+(width/10.),mean); f->SetParLimits(2,width/15.,width/2.); // The right-sided Gauss f->SetParLimits(3,0.,integ); f->SetParLimits(4,mean,max-(width/10.)); f->SetParLimits(5,width/15.,width/2.); // The Gauss describing the outliers f->SetParLimits(6,0.,integ); f->SetParLimits(7,min,max); f->SetParLimits(8,width/4.,width/1.5); obj2->Fit(f,"QLRM"); break; case 4: obj2->Fit("pol0", "Q"); obj2->GetFunction("pol0")->SetLineColor(kYellow); break; case 9: break; default: obj2->Fit("gaus", "Q"); obj2->GetFunction("gaus")->SetLineColor(kYellow); break; } } // -------------------------------------------------------------------------- // // Draw a projection of MHCamera vs. the radius from the central pixel. // // The inner and outer pixels are drawn separately, both fitted by a polynomial // of grade 1. // void MJCalibration::DrawRadialProfile(MHCamera *obj) const { TProfile *obj2 = (TProfile*)obj->RadialProfile(obj->GetName()); obj2->SetDirectory(0); obj2->Draw(); obj2->SetBit(kCanDelete); if (obj->GetGeomCam().InheritsFrom("MGeomCamMagic")) { TArrayI s0(6); s0[0] = 1; s0[1] = 2; s0[2] = 3; s0[3] = 4; s0[4] = 5; s0[5] = 6; TArrayI inner(1); inner[0] = 0; TArrayI outer(1); outer[0] = 1; // Just to get the right (maximum) binning TProfile *half[2]; half[0] = obj->RadialProfileS(s0, inner,Form("%s%s",obj->GetName(),"Inner")); half[1] = obj->RadialProfileS(s0, outer,Form("%s%s",obj->GetName(),"Outer")); for (Int_t i=0; i<2; i++) { Double_t min = obj->GetGeomCam().GetMinRadius(i); Double_t max = obj->GetGeomCam().GetMaxRadius(i); half[i]->SetLineColor(kRed+i); half[i]->SetDirectory(0); half[i]->SetBit(kCanDelete); half[i]->Draw("same"); half[i]->Fit("pol1","Q","",min,max); half[i]->GetFunction("pol1")->SetLineColor(kRed+i); half[i]->GetFunction("pol1")->SetLineWidth(1); } } } // -------------------------------------------------------------------------- // // Retrieve the output file written by WriteResult() // TString MJCalibration::GetOutputFile() const { if (!fRuns) return ""; return Form("%s/%s-F1.root", (const char*)fOutputPath, (const char*)fRuns->GetRunsAsFileName()); } // -------------------------------------------------------------------------- // // Call the ProcessFile(MPedestalCam) // Bool_t MJCalibration::Process(MPedestalCam &pedcam) { if (!ReadCalibrationCam()) return ProcessFile(pedcam); return kTRUE; } // -------------------------------------------------------------------------- // // Execute the task list and the eventloop: // // - Check if there are fRuns, otherwise return // - Check for consistency between run numbers and number of files // - Add fRuns to MReadMarsFile // - Put into MParList: // 1) MPedestalCam (pedcam) // 2) MCalibrationQECam (fQECam) // 3) MCalibrationChargeCam (fCalibrationCam) // 4) MCalibrationRelTimeCam (fRelTimeCam) (only if flag fRelTimes is chosen) // 5) MBadPixelsCam (fBadPixels) // 6) MCalibrationChargePINDiode // 7) MCalibrationChargeBlindPix // - Put into the MTaskList: // 1) MReadMarsFile // 2) MBadPixelsMerge // 3) MGeomApply // 4) MExtractSignal, MExtractSignal2 or MExtractSignal3, depending on fExtractorLevel // 5) MExtractPINDiode // 6) MExtractBlindPixel // 7) MArrivalTimeCalc, MArrivalTimeCalc2, depending on fArrivalTimeLevel (only if flag fRelTimes is chosen) // 8) MContinue(MFCosmics) // 9) MFillH("MHCalibrationChargePINDiode", "MExtractedSignalPINDiode") // 10) MFillH("MHCalibrationChargeBlindPix", "MExtractedSignalBlindPixel") // 11) MFillH("MHCalibrationChargeCam", "MExtractedSignalCam") // 12) MFillH("MHCalibrationChargeCam", "MExtractedSignalCam") // 13) MCalibrationChargeCalc // 14) MFillH("MHCalibrationRelTimeCam", "MArrivalTimeCam") (only if flag fRelTimes is chosen) // - Execute MEvtLoop // - DisplayResult() // - WriteResult() // Bool_t MJCalibration::ProcessFile(MPedestalCam &pedcam) { if (!fRuns) { *fLog << err << "No Runs choosen... abort." << endl; return kFALSE; } if (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 MCalibrationCam from Runs " << fRuns->GetRunsAsString() << endl; *fLog << endl; MReadMarsFile read("Events"); read.DisableAutoScheme(); static_cast(read).AddFiles(*fRuns); MCalibrationChargePINDiode pindiode; MCalibrationChargeBlindPix blindpix; // Setup Tasklist MParList plist; plist.AddToList(&pedcam); plist.AddToList(&fBadPixels); plist.AddToList(&fQECam); plist.AddToList(&fCalibrationCam); plist.AddToList(&pindiode); plist.AddToList(&blindpix); MTaskList tlist; plist.AddToList(&tlist); MGeomApply apply; // MBadPixelsMerge merge(&fBadPixels); MExtractPINDiode pinext; MExtractBlindPixel blindext; MExtractSignal extract1; MExtractSignal2 extract2; MExtractSignal3 extract3; MArrivalTimeCalc tmecalc1; MArrivalTimeCalc2 tmecalc2; MCalibrationChargeCalc calcalc; // // As long as there are no DM's, have to colour by hand // calcalc.SetPulserColor(fColor); MFillH fillpin("MHCalibrationChargePINDiode", "MExtractedSignalPINDiode"); MFillH fillbnd("MHCalibrationChargeBlindPix", "MExtractedSignalBlindPixel"); MFillH fillcam("MHCalibrationChargeCam", "MExtractedSignalCam"); MFillH filltme("MHCalibrationRelTimeCam", "MArrivalTimeCam"); fillpin.SetNameTab("PINDiode"); fillbnd.SetNameTab("BlindPix"); fillcam.SetNameTab("Charge"); filltme.SetNameTab("RelTimes"); // // Apply a filter against cosmics // (will have to be needed in the future // when the calibration hardware-trigger is working) // MFCosmics cosmics; MContinue cont(&cosmics); tlist.AddToList(&read); // tlist.AddToList(&merge); tlist.AddToList(&apply); if (fExtractorLevel <= 1) tlist.AddToList(&extract1); else if (fExtractorLevel == 2) tlist.AddToList(&extract2); else if (fExtractorLevel == 3) tlist.AddToList(&extract3); else { *fLog << err << GetDescriptor() << ": No valid Signal extractor has been chosen, have only: " << fExtractorLevel << " aborting..." << endl; return kFALSE; } tlist.AddToList(&pinext); tlist.AddToList(&blindext); if (fRelTimes) { plist.AddToList(&fRelTimeCam); if (fArrivalTimeLevel <= 1) tlist.AddToList(&tmecalc1); else if (fArrivalTimeLevel == 2) tlist.AddToList(&tmecalc2); else { *fLog << err << GetDescriptor() << ": No valid Signal ArrivalTime Level has been chosen, have only: " << fArrivalTimeLevel << " aborting..." << endl; return kFALSE; } } tlist.AddToList(&cont); tlist.AddToList(&fillcam); tlist.AddToList(&fillpin); tlist.AddToList(&fillbnd); tlist.AddToList(&calcalc); if (fRelTimes) tlist.AddToList(&filltme); // Create and setup the eventloop MEvtLoop evtloop(fName); evtloop.SetParList(&plist); evtloop.SetDisplay(fDisplay); evtloop.SetLogStream(fLog); // 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; } // -------------------------------------------------------------------------- // // Read the following containers from GetOutputFile() // - MCalibrationChargeCam // - MCalibrationQECam // - MBadPixelsCam // Bool_t MJCalibration::ReadCalibrationCam() { const TString fname = GetOutputFile(); if (gSystem->AccessPathName(fname, kFileExists)) { *fLog << err << "Input file " << fname << " doesn't exist." << endl; return kFALSE; } *fLog << inf << "Reading from file: " << fname << endl; TFile file(fname, "READ"); if (fCalibrationCam.Read()<=0) { *fLog << err << "Unable to read MCalibrationChargeCam from " << fname << endl; return kFALSE; } if (fQECam.Read()<=0) { *fLog << err << "Unable to read MCalibrationQECam 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")*/) // FIXME! fDisplay->Read(); return kTRUE; } // -------------------------------------------------------------------------- // // Set the path for output files, written by WriteResult() // void MJCalibration::SetOutputPath(const char *path) { fOutputPath = path; if (fOutputPath.EndsWith("/")) fOutputPath = fOutputPath(0, fOutputPath.Length()-1); } // -------------------------------------------------------------------------- // // Write the result into the output file GetOutputFile(), if fOutputPath exists. // // The following containers are written: // - MStatusDisplay // - MCalibrationChargeCam // - MCalibrationQECam // - MBadPixelsCam // Bool_t MJCalibration::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 (fCalibrationCam.Write()<=0) { *fLog << err << "Unable to write MCalibrationChargeCam to " << oname << endl; return kFALSE; } if (fQECam.Write()<=0) { *fLog << err << "Unable to write MCalibrationQECam to " << oname << endl; return kFALSE; } if (fBadPixels.Write()<=0) { *fLog << err << "Unable to write MBadPixelsCam to " << oname << endl; return kFALSE; } return kTRUE; }