Forked from
faser / calypso
1248 commits behind the upstream repository.
-
Savannah Rose Shively authoredSavannah Rose Shively authored
Code owners
Assign users and groups as approvers for specific file changes. Learn more.
SimHitAlg.cxx 4.60 KiB
#include "SimHitAlg.h"
SimHitAlg::SimHitAlg(const std::string& name, ISvcLocator* pSvcLocator)
: AthHistogramAlgorithm(name, pSvcLocator) { m_hist = nullptr; }
SimHitAlg::~SimHitAlg() { }
StatusCode SimHitAlg::initialize()
{
// initialize a histogram
// letter at end of TH1 indicated variable type (D double, F float etc)
m_hist = new TH1D("eLoss", "SCT Hit Energy Loss", 100, 0, 1); //first string is root object name, second is histogram title
m_module = new TH2D("module", "SCT Hit Module", 3, -1.5, 1.5, 4, -0.5, 3.5 );
m_moduleSide1 = new TH2D("moduleSide1", "SCT Hit Module", 3, -1.5, 1.5, 4, -0.5, 3.5 );
m_moduleSide2 = new TH2D("moduleSide2", "SCT Hit Module", 3, -1.5, 1.5, 4, -0.5, 3.5 );
ATH_CHECK(histSvc()->regHist("/HIST/eloss", m_hist));
ATH_CHECK(histSvc()->regHist("/HIST/modules", m_module));
ATH_CHECK(histSvc()->regHist("/HIST/modulesSide1", m_moduleSide1));
ATH_CHECK(histSvc()->regHist("/HIST/modulesSide2", m_moduleSide2));
m_plate = new TH2D("plate", "Scint Hit Plate", 3, -1, 1, 4, 0, 1 );
ATH_CHECK(histSvc()->regHist("/HIST/plates", m_plate));
// initialize data handle keys
ATH_CHECK( m_mcEventKey.initialize() );
ATH_CHECK( m_faserSiHitKey.initialize() );
//Sav things
ATH_CHECK( m_preshowerHitKey.initialize() );
ATH_CHECK( m_triggerHitKey.initialize() );
ATH_CHECK( m_vetoHitKey.initialize() );
ATH_MSG_INFO( "Using GenEvent collection with key " << m_mcEventKey.key());
ATH_MSG_INFO( "Using Faser SiHit collection with key " << m_faserSiHitKey.key());
ATH_MSG_INFO( "Using ScintHit collection with key " << m_preshowerHitKey.key());
ATH_MSG_INFO( "Using ScintHit collection with key " << m_triggerHitKey.key());
ATH_MSG_INFO( "Using ScintHit collection with key " << m_vetoHitKey.key());
return StatusCode::SUCCESS;
}
StatusCode SimHitAlg::execute()
{
// Handles created from handle keys behave like pointers to the corresponding container
SG::ReadHandle<McEventCollection> h_mcEvents(m_mcEventKey);
ATH_MSG_INFO("Read McEventContainer with " << h_mcEvents->size() << " events");
if (h_mcEvents->size() == 0) return StatusCode::FAILURE;
SG::ReadHandle<FaserSiHitCollection> h_siHits(m_faserSiHitKey);
ATH_MSG_INFO("Read FaserSiHitCollection with " << h_siHits->size() << " hits");
SG::ReadHandle<ScintHitCollection> h_preshowerHits(m_preshowerHitKey);
ATH_MSG_INFO("Read ScintHitCollection/Preshower with " << h_preshowerHits->size() << " hits");
SG::ReadHandle<ScintHitCollection> h_triggerHits(m_triggerHitKey);
ATH_MSG_INFO("Read ScintHitCollection/Trigger with " << h_triggerHits->size() << " hits");
SG::ReadHandle<ScintHitCollection> h_vetoHits(m_vetoHitKey);
ATH_MSG_INFO("Read ScintHitCollectionVeto with " << h_vetoHits->size() << " hits");
// Since we have no pile-up, there should always be a single GenEvent in the container
const HepMC::GenEvent* ev = (*h_mcEvents)[0];
if (ev == nullptr)
{
ATH_MSG_FATAL("GenEvent pointer is null");
return StatusCode::FAILURE;
}
ATH_MSG_INFO("Event contains " << ev->particles_size() << " truth particles" );
// The hit container might be empty because particles missed the wafers
//if (h_siHits->size() == 0) return StatusCode::SUCCESS;
if (h_siHits->size()!=0){
// Loop over all hits; print and fill histogram
for (const FaserSiHit& hit : *h_siHits)
{
hit.print();
m_hist->Fill( hit.energyLoss() );
m_module->Fill( hit.getModule(), hit.getRow() );
if (hit.getSensor() == 0)
{
m_moduleSide1->Fill( hit.getModule(), hit.getRow());
}
else
{
m_moduleSide2->Fill( hit.getModule(), hit.getRow());
}
}
}
if (h_preshowerHits->size()!=0){
for (const ScintHit& hit : *h_preshowerHits)
{
hit.print();
m_hist->Fill(hit.energyLoss());
m_plate->Fill(hit.getStation(),hit.getPlate());
}
}
if (h_triggerHits->size()!=0){
for (const ScintHit& hit : *h_triggerHits)
{
hit.print();
m_hist->Fill(hit.energyLoss());
m_plate->Fill(hit.getStation(),hit.getPlate());
}
}
if (h_vetoHits->size()!=0){
for (const ScintHit& hit : *h_vetoHits)
{
hit.print();
m_hist->Fill(hit.energyLoss());
m_plate->Fill(hit.getStation(),hit.getPlate());
}
}
return StatusCode::SUCCESS;
}
StatusCode SimHitAlg::finalize()
{
return StatusCode::SUCCESS;
}