////////////////////////////////////////////////////////////////////
/// \file PlotHitTimes.cc
///
/// \brief Functions to plot hit times
///
/// \author P G Jones
///
/// REVISION HISTORY:\n
/// 2014-03-27 : P G Jones - First Revision.\n
///
/// \details Hit times of the photoelectron, MCHit, uncalibrated and
/// calibrated PMTs can be plotted. Uncalibrated hit times are ADC
/// counts and thus not easily compared.
///
////////////////////////////////////////////////////////////////////
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
/// Plot the MC photoelectron hit times
///
/// @param[in] fileName of the RAT::DS root file to analyse
/// @return the histogram plot
TH1D* PlotMCPEHitTimes( const std::string& fileName )
{
TH1D* hHitTimes = new TH1D( "hHitTimes", "Hit times of the of photoelectrons", 500, 0.0, 500.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 );
const RAT::DS::MC& rMC = rDS.GetMC();
for( size_t iPMT = 0; iPMT < rMC.GetMCPMTCount(); iPMT++ )
{
const RAT::DS::MCPMT& rMCPMT = rMC.GetMCPMT( iPMT );
for( size_t iPE = 0; iPE < rMCPMT.GetMCPECount(); iPE++ )
hHitTimes->Fill( rMCPMT.GetMCPE( iPE ).GetCreationTime() );
}
}
hHitTimes->GetYaxis()->SetTitle( "Count per 1 ns bin" );
hHitTimes->GetXaxis()->SetTitle( "Hit times [ns]" );
hHitTimes->Draw();
return hHitTimes;
}
/// Plot the MCHit hit times per event
///
/// @param[in] fileName of the RAT::DS root file to analyse
/// @return the histogram plot
TH1D* PlotMCHitTimes( const std::string& fileName )
{
TH1D* hHitTimes = new TH1D( "hHitTimes", "Hit times of the MCHit per event", 500, 0.0, 500.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 iMCEvent = 0; iMCEvent < rDS.GetMCEVCount(); iMCEvent++ ) // Loop over triggered events
{
const RAT::DS::MCHits& mcHits = rDS.GetMCEV( iMCEvent ).GetMCHits();
for( size_t iMCHit = 0; iMCHit < mcHits.GetAllCount(); iMCHit++ )
hHitTimes->Fill( mcHits.GetAllPMT( iMCHit ).GetTime() );
}
}
hHitTimes->GetYaxis()->SetTitle( "Count per 1 ns bin" );
hHitTimes->GetXaxis()->SetTitle( "Hit times [ns]" );
hHitTimes->Draw();
return hHitTimes;
}
/// Plot the calibrated hit times
///
/// @param[in] fileName of the RAT::DS root file to analyse
/// @return the histogram plot
TH1D* PlotCalibratedHitTimes( const std::string& fileName )
{
TH1D* hHitTimes = new TH1D( "hHitTimes", "Hit times for the calibrated PMTs", 500, 0.0, 500.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++ )
{
const RAT::DS::CalPMTs& calibratedPMTs = rDS.GetEV( iEV ).GetCalPMTs();
for( size_t iPMT = 0; iPMT < calibratedPMTs.GetAllCount(); iPMT++ )
hHitTimes->Fill( calibratedPMTs.GetAllPMT( iPMT ).GetTime() );
}
}
hHitTimes->GetYaxis()->SetTitle( "Count per 1 ns bin" );
hHitTimes->GetXaxis()->SetTitle( "Hit times [ns]" );
hHitTimes->Draw();
return hHitTimes;
}
/// Plot all NHit numbers
///
/// @param[in] fileName of the RAT::DS root file to analyse
void PlotHitTimes( const std::string& fileName )
{
TCanvas* c1 = new TCanvas();
TH1D* numPE = PlotMCPEHitTimes( fileName );
TH1D* mcHits = PlotMCHitTimes( fileName );
TH1D* calibrated = PlotCalibratedHitTimes( fileName );
numPE->Draw();
mcHits->SetLineColor( kGreen + 2 );
mcHits->Draw("SAME");
calibrated->SetLineColor( kRed );
calibrated->Draw("SAME");
TLegend* t1 = new TLegend( 0.7, 0.7, 0.88, 0.88 );
t1->AddEntry( numPE, "Photoelectrons", "l" );
t1->AddEntry( mcHits, "MC Hits", "l" );
t1->AddEntry( calibrated, "Calibrated", "l" );
t1->SetLineColor( kWhite );
t1->SetFillColor( kWhite );
t1->Draw();
c1->Update();
}