//////////////////////////////////////////////////////////////////// /// \file PlotNhit.cc /// /// \brief Functions to plot nhit (various definitions) /// /// \author P G Jones /// /// REVISION HISTORY:\n /// 2014-03-26 : P G Jones - First Revision.\n /// /// \details Nhit could be the MC photoelectrons, the MC hits, the /// EV uncalibrated or the EV calibrated. All examples are included. /// //////////////////////////////////////////////////////////////////// #include #include #include #include #include #include #include #include #include /// Plot the number of NhitsCleaned per event /// /// @param[in] fileName of the RAT::DS root file to analyse /// @return the histogram plot TH1D* PlotNhitsCleaned(const std::string& fileName) { TH1D* hNHitsCleaned = new TH1D( "hNHitsCleaned", "Number of cleaned Nhits", 2000, 0.0, 2000.0); RAT::DU::DSReader dsReader(fileName); for( size_t iEntry = 0; iEntry < dsReader.GetEntryCount(); iEntry++ ) { const RAT::DS::Entry& rDS = dsReader.GetEntry( iEntry ); for(size_t iev = 0; iev < rDS.GetEVCount(); iev++) // Loop over triggered events hNHitsCleaned->Fill(rDS.GetEV(iev).GetNhitsCleaned()); } hNHitsCleaned->GetYaxis()->SetTitle( "Count per bin" ); hNHitsCleaned->GetXaxis()->SetTitle( "Number of NhitsCleaned per event" ); hNHitsCleaned->Draw(); return hNHitsCleaned; } /// Plot the number of NhitsCleaned per event re-normalised for detector state /// This is useful for like-to-like comparison between runs. Uses ideal /// scintillator state as reference by default. /// /// @param[in] fileName of the RAT::DS root file to analyse /// @return the histogram plot TH1D* PlotNhitsCleanedCorrected(const std::string& fileName) { TH1D* hNHitsCleanedCorrected = new TH1D( "hNHitsCleanedCorrected", "Number of cleaned Nhits corrected for detector state", 2000, 0.0, 2000.0); RAT::DU::DSReader dsReader(fileName); RAT::DU::DetectorStateCorrection stateCorrUtility = RAT::DU::Utility::Get()->GetDetectorStateCorrection(); for( size_t iEntry = 0; iEntry < dsReader.GetEntryCount(); iEntry++ ) { const RAT::DS::Entry& rDS = dsReader.GetEntry( iEntry ); for(size_t iev = 0; iev < rDS.GetEVCount(); iev++){ // Loop over triggered events double rawNhitsCleaned = rDS.GetEV(iev).GetNhitsCleaned(); double corrected = stateCorrUtility.ApplyCorrectionNhits(rawNhitsCleaned); hNHitsCleanedCorrected->Fill(corrected); } } hNHitsCleanedCorrected->GetYaxis()->SetTitle( "Count per bin" ); hNHitsCleanedCorrected->GetXaxis()->SetTitle( "Number of corrected NhitsCleaned per event" ); hNHitsCleanedCorrected->Draw(); return hNHitsCleanedCorrected; } /// Plot the number of MC photoelectrons per event (or NumPE) /// /// @param[in] fileName of the RAT::DS root file to analyse /// @return the histogram plot TH1D* PlotMCPhotoelectronNhit( const std::string& fileName ) { TH1D* hNumPE = new TH1D( "hNumPE", "Number of photoelectrons per event", 2000, 0.0, 2000.0 ); // If this is being done on data that does not require remote database connection // eg.: a simple simulation with default run number (0) // We can disable the remote connections: // // NOTE: Don't do this if you are using real data!!! RAT::DB::Get()->SetAirplaneModeStatus(true); RAT::DU::DSReader dsReader( fileName ); for( size_t iEntry = 0; iEntry < dsReader.GetEntryCount(); iEntry++ ) { const RAT::DS::Entry& rDS = dsReader.GetEntry( iEntry ); hNumPE->Fill( rDS.GetMC().GetMCPECount() ); } hNumPE->GetYaxis()->SetTitle( "Count per 1 pe bin" ); hNumPE->GetXaxis()->SetTitle( "Number of photoelectrons per event" ); hNumPE->Draw(); return hNumPE; } /// Plot the number of MC hits per event /// /// @param[in] fileName of the RAT::DS root file to analyse /// @return the histogram plot TH1D* PlotMCHitsNhit( const std::string& fileName ) { TH1D* hNumHits = new TH1D( "hNumHits", "Number of MC hits per event", 2000, 0.0, 2000.0 ); RAT::DU::DSReader dsReader( fileName ); for( size_t iEntry = 0; iEntry < dsReader.GetEntryCount(); iEntry++ ) { const RAT::DS::Entry& rDS = dsReader.GetEntry( iEntry ); for( size_t iMCEvent = 0; iMCEvent < rDS.GetMCEVCount(); iMCEvent++ ) // Loop over triggered events hNumHits->Fill( rDS.GetMCEV( iMCEvent ).GetMCHits().GetAllCount() ); // All can be changed to Normal, OWL etc... } hNumHits->GetYaxis()->SetTitle( "Count per 1 hit bin" ); hNumHits->GetXaxis()->SetTitle( "Number of MC hits per event" ); hNumHits->Draw(); return hNumHits; } /// Plot the number of Uncalibrated hits per event /// /// @param[in] fileName of the RAT::DS root file to analyse /// @return the histogram plot TH1D* PlotUncalibratedNhit( const std::string& fileName ) { TH1D* hUncalHits = new TH1D( "hUncalHits", "Number of Uncalibrated hits per event", 2000, 0.0, 2000.0 ); // If this is being done on data that does not require remote database connection // eg.: a simple simulation with default run number (0) // We can disable the remote connections: // // NOTE: Don't do this if you are using real data!!! RAT::DB::Get()->SetAirplaneModeStatus(true); RAT::DU::DSReader dsReader( fileName ); for( size_t iEntry = 0; iEntry < dsReader.GetEntryCount(); iEntry++ ) { const RAT::DS::Entry& rDS = dsReader.GetEntry( iEntry ); for( size_t iEV = 0; iEV < rDS.GetEVCount(); iEV++ ) hUncalHits->Fill( rDS.GetEV( iEV ).GetUncalPMTs().GetAllCount() ); } hUncalHits->GetYaxis()->SetTitle( "Count per 1 hit bin" ); hUncalHits->GetXaxis()->SetTitle( "Number of Uncalibrated hits per event" ); hUncalHits->Draw(); return hUncalHits; } /// Plot the number of Calibrated hits per event /// /// @param[in] fileName of the RAT::DS root file to analyse /// @return the histogram plot TH1D* PlotCalibratedNhit( const std::string& fileName ) { TH1D* hCalHits = new TH1D( "hCalHits", "Number of Calibrated hits per event", 2000, 0.0, 2000.0 ); // If this is being done on data that does not require remote database connection // eg.: a simple simulation with default run number (0) // We can disable the remote connections: // // NOTE: Don't do this if you are using real data!!! RAT::DB::Get()->SetAirplaneModeStatus(true); RAT::DU::DSReader dsReader( fileName ); for( size_t iEntry = 0; iEntry < dsReader.GetEntryCount(); iEntry++ ) { const RAT::DS::Entry& rDS = dsReader.GetEntry( iEntry ); for( size_t iEV = 0; iEV < rDS.GetEVCount(); iEV++ ) hCalHits->Fill( rDS.GetEV( iEV ).GetCalPMTs().GetAllCount() ); } hCalHits->GetYaxis()->SetTitle( "Count per 1 hit bin" ); hCalHits->GetXaxis()->SetTitle( "Number of calibrated hits per event" ); hCalHits->Draw(); return hCalHits; } /// Plot all NHit numbers /// /// @param[in] fileName of the RAT::DS root file to analyse void PlotNhit( const std::string& fileName ) { TCanvas* c1 = new TCanvas(); TH1D* numPE = PlotMCPhotoelectronNhit( fileName ); TH1D* mcHits = PlotMCHitsNhit( fileName ); TH1D* uncalibrated = PlotUncalibratedNhit( fileName ); TH1D* calibrated = PlotCalibratedNhit( fileName ); TH1D* cleaned = PlotNhitsCleaned( fileName ); TH1D* corrected = PlotNhitsCleanedCorrected( fileName ); numPE->Draw(); mcHits->SetLineColor( kGreen + 2 ); mcHits->Draw("SAME"); uncalibrated->SetLineColor( kBlue ); uncalibrated->Draw("SAME"); calibrated->SetLineColor( kRed ); calibrated->Draw("SAME"); cleaned->SetLineColor( kMagenta+2 ); cleaned->Draw("SAME"); corrected->SetLineColor( kOrange-3 ); corrected->Draw("SAME"); TLegend* t1 = new TLegend( 0.65, 0.65, 0.9, 0.9 ); t1->AddEntry( numPE, "Photoelectrons", "l" ); t1->AddEntry( mcHits, "MC Hits", "l" ); t1->AddEntry( uncalibrated, "Uncalibrated", "l" ); t1->AddEntry( calibrated, "Calibrated", "l" ); t1->AddEntry( cleaned, "Cleaned", "l" ); t1->AddEntry( corrected, "Cleaned and Corrected", "l" ); t1->Draw(); c1->Update(); }