diff --git a/source/Utils/include/FilterJetConeHits.h b/source/Utils/include/FilterJetConeHits.h new file mode 100644 index 0000000..1c8f05f --- /dev/null +++ b/source/Utils/include/FilterJetConeHits.h @@ -0,0 +1,92 @@ +#ifndef FilterJetConeHits_h +#define FilterJetConeHits_h 1 + +#include "marlin/Processor.h" +#include "lcio.h" +#include +#include + +#include + +#include +#include +#include "IMPL/ReconstructedParticleImpl.h" + +using namespace lcio ; +using namespace marlin ; + + +/** Utility processor that selects and saves the tracker hits that are included in + * a DeltaR cone around the jet direction direction along with the corresponding + * sim hits and the reco-sim relations. + * + * @parameter MCParticleCollection name of the MCParticle collection + * @parameter TrackerHitInputCollections name of the tracker hit input collections + * @parameter TrackerSimHitInputCollections name of the tracker simhit input collections + * @parameter TrackerHitInputRelations name of the tracker hit relation input collections + * @parameter TrackerHitOutputCollections name of the tracker hit output collections + * @parameter TrackerSimHitOutputCollections name of the tracker simhit output collections + * @parameter TrackerHitOutputRelations name of the tracker hit relation output collections + * @parameter DeltaRCut maximum angular distance between the hits and the particle direction + * @parameter FillHistograms flag to fill the diagnostic histograms + * + * @author L. Sestini, INFN Padova + * @date 18 March 2021 + * @version $Id: $ + */ + +class FilterJetConeHits : public Processor { + + public: + + virtual Processor* newProcessor() { return new FilterJetConeHits ; } + + + FilterJetConeHits() ; + + virtual void init() ; + + virtual void processRunHeader( LCRunHeader* run ) ; + + virtual void processEvent( LCEvent * evt ) ; + + virtual void check( LCEvent * evt ) ; + + virtual void end() ; + + bool filterJetBib(ReconstructedParticle* jet); + + + protected: + + // --- Input/output collection names: + std::string m_inputJetCaloCollName{} ; + std::vector m_inputTrackerHitsCollNames{} ; + std::vector m_inputTrackerSimHitsCollNames{} ; + std::vector m_inputTrackerHitRelNames{} ; + std::vector m_outputTrackerHitsCollNames{} ; + std::vector m_outputTrackerSimHitsCollNames{} ; + std::vector m_outputTrackerHitRelNames{} ; + + + // --- Processor parameters: + bool m_fillHistos{} ; + double m_deltaRCut{} ; + + // Jet filter parameters with BIB: + double m_minDaughterMaxPt{} ; + int m_minNTracks{} ; + + // --- Diagnostic histograms: + TH1F* m_dist = nullptr ; + + // --- Run and event counters: + int _nRun{} ; + int _nEvt{} ; + +} ; + +#endif + + + diff --git a/source/Utils/src/FilterJetConeHits.cc b/source/Utils/src/FilterJetConeHits.cc new file mode 100644 index 0000000..bc54335 --- /dev/null +++ b/source/Utils/src/FilterJetConeHits.cc @@ -0,0 +1,434 @@ +#include "FilterJetConeHits.h" +#include +#include +#include + + +#include + +#include +#include +#include +#include +#include + +#include +#include + +#include "HelixClass_double.h" + +using namespace lcio ; +using namespace marlin ; + + +FilterJetConeHits aFilterJetConeHits ; + + +FilterJetConeHits::FilterJetConeHits() : Processor("FilterJetConeHits") { + + // --- Processor description: + + _description = "FilterJetConeHits selects tracker hits in a cone opened around a jet cone direction"; + + + // --- Processor parameters: + + registerProcessorParameter("JetCaloCollection", + "Name of the JetCalo collection", + m_inputJetCaloCollName, + std::string("JetCalo") ); + + registerProcessorParameter("TrackerHitInputCollections", + "Name of the tracker hit input collections", + m_inputTrackerHitsCollNames, + {} ); + + registerProcessorParameter("TrackerSimHitInputCollections", + "Name of the tracker simhit input collections", + m_inputTrackerSimHitsCollNames, + {} ); + + registerProcessorParameter("TrackerHitInputRelations", + "Name of the tracker hit relation collections", + m_inputTrackerHitRelNames, + {} ); + + registerProcessorParameter("TrackerHitOutputCollections", + "Name of the tracker hit output collections", + m_outputTrackerHitsCollNames, + {} ); + + registerProcessorParameter("TrackerSimHitOutputCollections", + "Name of the tracker simhit output collections", + m_outputTrackerSimHitsCollNames, + {} ); + + registerProcessorParameter("TrackerHitOutputRelations", + "Name of the tracker hit relation collections", + m_outputTrackerHitRelNames, + {} ); + + registerProcessorParameter( "DeltaRCut" , + "Maximum angular distance between the hits and the particle direction" , + m_deltaRCut, + double(1.) ); + + registerProcessorParameter( "FillHistograms", + "Flag to fill the diagnostic histograms", + m_fillHistos, + false ); + + registerProcessorParameter( "MinJetTracks", + "Min number of tracks in jet to use it as filter", + m_minNTracks, + int(1) ); + + registerProcessorParameter( "MinDaughterMaxPt", + "Min pT of the highest-pT track in jet to use it as filter", + m_minDaughterMaxPt, + double(2.) ); + + +} + + + +void FilterJetConeHits::init() { + + streamlog_out(DEBUG) << " init called " << std::endl ; + + // --- Print the processor parameters: + + printParameters() ; + + + // --- Initialize the run and event counters: + + _nRun = 0 ; + _nEvt = 0 ; + + + // --- Initialize the AIDAProcessor and book the diagnostic histograms: + + AIDAProcessor::histogramFactory(this); + + m_dist = new TH1F("m_dist", "deltaR distance between hit and jet axis", 1000, 0., 6.); + +} + + +void FilterJetConeHits::processRunHeader( LCRunHeader* run) { + + _nRun++ ; + +} + + + +void FilterJetConeHits::processEvent( LCEvent * evt ) { + + // --- Check whether the number of input and output collections match + bool noSimColl = false, noRelColl = false; + if ( m_inputTrackerSimHitsCollNames.size() == 0 ) noSimColl = true; + if ( m_inputTrackerHitRelNames.size() == 0 ) noRelColl = true; + if ( noSimColl != noRelColl ) { + std::stringstream err_msg; + err_msg << "Error: RelColl are necessary to save SimColl!" + << std::endl ; + + throw EVENT::Exception( err_msg.str() ) ; + } + + if ( m_inputTrackerHitsCollNames.size() != m_outputTrackerHitsCollNames.size() || + ( !noSimColl && m_inputTrackerSimHitsCollNames.size() != m_outputTrackerSimHitsCollNames.size() ) || + ( !noRelColl && m_inputTrackerHitRelNames.size() != m_outputTrackerHitRelNames.size() ) ) { + + std::stringstream err_msg; + err_msg << "Mismatch between input and output collections" + << std::endl ; + + throw EVENT::Exception( err_msg.str() ) ; + } + + if ( !noSimColl && ( m_inputTrackerHitsCollNames.size() != m_inputTrackerSimHitsCollNames.size() || + m_inputTrackerHitsCollNames.size() != m_inputTrackerHitRelNames.size() ) ){ + + std::stringstream err_msg; + err_msg << "Mismatch between the reco and sim hits input collections" + << std::endl ; + + throw EVENT::Exception( err_msg.str() ) ; + + } + + if ( !noSimColl && ( m_outputTrackerHitsCollNames.size() != m_outputTrackerSimHitsCollNames.size() || + m_outputTrackerHitsCollNames.size() != m_outputTrackerHitRelNames.size() ) ){ + + std::stringstream err_msg; + err_msg << "Mismatch between the reco and sim hits output collections" + << std::endl ; + + throw EVENT::Exception( err_msg.str() ) ; + + } + + + // --- Get the JetCalo collection: + + LCCollection* m_inputJetCalo = NULL; + try { + m_inputJetCalo = evt->getCollection( m_inputJetCaloCollName ); + } + catch( lcio::DataNotAvailableException& e ) { + streamlog_out(WARNING) << m_inputJetCaloCollName << " collection not available" << std::endl; + return; + } + + + // --- Get the input hit collections and create the corresponding output collections: + + const unsigned int nTrackerHitCol = m_inputTrackerHitsCollNames.size(); + std::vector inputHitColls(nTrackerHitCol); + std::vector inputSimHitColls(nTrackerHitCol); + std::vector inputHitRels(nTrackerHitCol); + + std::vector outputTrackerHitColls(nTrackerHitCol); + std::vector outputTrackerSimHitColls(nTrackerHitCol); + std::vector outputTrackerHitRels(nTrackerHitCol); + + for (unsigned int icol=0; icolgetCollection(m_inputTrackerHitsCollNames[icol]); + } + catch( lcio::DataNotAvailableException& e ) { + streamlog_out(WARNING) << m_inputTrackerHitsCollNames[icol] + << " collection not available" << std::endl; + continue; + } + + // get the sim hits + if(!noSimColl){ + + try { + inputSimHitColls[icol] = evt->getCollection(m_inputTrackerSimHitsCollNames[icol]); + } + catch( lcio::DataNotAvailableException& e ) { + streamlog_out(WARNING) << m_inputTrackerSimHitsCollNames[icol] + << " collection not available" << std::endl; + continue; + } + + try { + inputHitRels[icol] = evt->getCollection(m_inputTrackerHitRelNames[icol]); + } + catch( lcio::DataNotAvailableException& e ) { + streamlog_out(WARNING) << m_inputTrackerHitRelNames[icol] + << " collection not available" << std::endl; + continue; + } + } + + // reco hit output collections + std::string encoderString = inputHitColls[icol]->getParameters().getStringVal( "CellIDEncoding" ); + outputTrackerHitColls[icol] = new LCCollectionVec( inputHitColls[icol]->getTypeName() ); + outputTrackerHitColls[icol]->parameters().setValue( "CellIDEncoding", encoderString ); + LCFlagImpl lcFlag(inputHitColls[icol]->getFlag()); + outputTrackerHitColls[icol]->setFlag(lcFlag.getFlag()); + + // sim hit output collections + if(!noSimColl){ + outputTrackerSimHitColls[icol] = new LCCollectionVec( inputSimHitColls[icol]->getTypeName() ); + outputTrackerSimHitColls[icol]->parameters().setValue( "CellIDEncoding", encoderString ); + LCFlagImpl lcFlag_sim(inputSimHitColls[icol]->getFlag()); + outputTrackerSimHitColls[icol]->setFlag(lcFlag_sim.getFlag()); + + outputTrackerHitRels[icol] = new LCCollectionVec( inputHitRels[icol]->getTypeName() ); + LCFlagImpl lcFlag_rel(inputHitRels[icol]->getFlag()) ; + outputTrackerHitRels[icol]->setFlag( lcFlag_rel.getFlag() ) ; + } + + } + + + // --- Loop over the JetCalo: + + std::vector > hits_to_save(nTrackerHitCol); + + for (int ipart=0; ipartgetNumberOfElements(); ++ipart){ + + ReconstructedParticle* part = dynamic_cast( m_inputJetCalo->getElementAt(ipart) ); + + if( !FilterJetConeHits::filterJetBib(part) ) continue; //quality selection to reject fake jets + + // double part_p = sqrt( part->getMomentum()[0]*part->getMomentum()[0] +part->getMomentum()[1]*part->getMomentum()[1] + part->getMomentum()[2]*part->getMomentum()[2] ); + + // --- Loop over the tracker hits and select hits inside a cone around the jet axis: + + for (unsigned int icol=0; icolgetNumberOfElements(); ++ihit){ + + TrackerHitPlane* hit = dynamic_cast(hit_col->getElementAt(ihit)); + + // --- Skip hits that are in the opposite hemisphere w.r.t. the jet axis: + if ( ( hit->getPosition()[0]*part->getMomentum()[0] + + hit->getPosition()[1]*part->getMomentum()[1] + + hit->getPosition()[2]*part->getMomentum()[2] ) < 0. ) continue; + + + // --- Get the distance between the hit and the jet axis + + double jet_p = sqrt( pow(part->getMomentum()[0],2) + pow(part->getMomentum()[1],2) + pow(part->getMomentum()[2],2) ); + double jet_theta = acos(part->getMomentum()[2]/jet_p); + double jet_eta = -std::log(tan(jet_theta/2.)); + + double hit_d = sqrt( pow(hit->getPosition()[0],2) + pow(hit->getPosition()[1],2) + pow(hit->getPosition()[2],2) ); + double hit_theta = acos(hit->getPosition()[2]/hit_d); + double hit_eta = -std::log(tan(hit_theta/2.)); + + double jet_pxy = sqrt( pow(part->getMomentum()[0],2) + pow(part->getMomentum()[1],2) ); + double hit_dxy = sqrt( pow(hit->getPosition()[0],2) + pow(hit->getPosition()[1],2) ); + double deltaPhi = acos( (part->getMomentum()[0]*hit->getPosition()[0] + part->getMomentum()[1]*hit->getPosition()[1] )/jet_pxy/hit_dxy); + + double deltaR = sqrt( pow(deltaPhi,2) + pow(jet_eta-hit_eta,2) ); + + if ( m_fillHistos ){ + + m_dist->Fill(deltaR); + + } + + if ( deltaR < m_deltaRCut ) + hits_to_save[icol].insert(ihit); + + } // ihit loop + + } // icol loop + + } // ipart loop + + + // --- Add the filtered hits to the output collections: + + for (unsigned int icol=0; icol(inputHitColls[icol]->getElementAt(ihit)); + TrackerHitPlaneImpl* hit_new = new TrackerHitPlaneImpl(); + + hit_new->setCellID0(hit->getCellID0()); + hit_new->setCellID1(hit->getCellID1()); + hit_new->setType(hit->getType()); + hit_new->setPosition(hit->getPosition()); + hit_new->setU(hit->getU()); + hit_new->setV(hit->getV()); + hit_new->setdU(hit->getdU()); + hit_new->setdV(hit->getdV()); + hit_new->setEDep(hit->getEDep()); + hit_new->setEDepError(hit->getEDepError()); + hit_new->setTime(hit->getTime()); + hit_new->setQuality(hit->getQuality()); + + outputTrackerHitColls[icol]->addElement( hit_new ); + + if( noSimColl ) continue; + + LCRelation* rel = dynamic_cast(inputHitRels[icol]->getElementAt(ihit)); + + + SimTrackerHit* simhit = dynamic_cast(rel->getTo()); + SimTrackerHitImpl* simhit_new = new SimTrackerHitImpl(); + + simhit_new->setCellID0(simhit->getCellID0()); + simhit_new->setCellID1(simhit->getCellID1()); + simhit_new->setPosition(simhit->getPosition()); + simhit_new->setEDep(simhit->getEDep()); + simhit_new->setTime(simhit->getTime()); + simhit_new->setMCParticle(simhit->getMCParticle()); + simhit_new->setMomentum(simhit->getMomentum()); + simhit_new->setPathLength(simhit->getPathLength()); + simhit_new->setQuality(simhit->getQuality()); + simhit_new->setOverlay(simhit->isOverlay()); + simhit_new->setProducedBySecondary(simhit->isProducedBySecondary()); + + outputTrackerSimHitColls[icol]->addElement( simhit_new ); + + + LCRelationImpl* rel_new = new LCRelationImpl(); + + rel_new->setFrom(hit_new); + rel_new->setTo(simhit_new); + rel_new->setWeight(rel->getWeight()); + + outputTrackerHitRels[icol]->addElement( rel_new ); + + } // ihit loop + + streamlog_out( MESSAGE ) << " " << hits_to_save[icol].size() << " hits added to the collections: " + << m_outputTrackerHitsCollNames[icol] << std::endl; + + evt->addCollection( outputTrackerHitColls[icol], m_outputTrackerHitsCollNames[icol] ) ; + if(!noSimColl){ + evt->addCollection( outputTrackerSimHitColls[icol], m_outputTrackerSimHitsCollNames[icol] ) ; + evt->addCollection( outputTrackerHitRels[icol], m_outputTrackerHitRelNames[icol] ) ; + } + + streamlog_out( DEBUG5 ) << " output collection " << m_outputTrackerHitsCollNames[icol] << " of type " + << outputTrackerHitColls[icol]->getTypeName() << " added to the event" << std::endl; + if(!noSimColl) streamlog_out( DEBUG5 ) << " output collection " << m_outputTrackerSimHitsCollNames[icol] << " of type " + << outputTrackerSimHitColls[icol]->getTypeName() << " added to the event \n" + << " output collection " << m_outputTrackerHitRelNames[icol] << " of type " + << outputTrackerHitRels[icol]->getTypeName() << " added to the event " + << std::endl ; + + } // icol loop + + streamlog_out(DEBUG) << " processing event: " << evt->getEventNumber() + << " in run: " << evt->getRunNumber() << std::endl ; + + _nEvt ++ ; + +} + + + +void FilterJetConeHits::check( LCEvent * evt ) { +} + + + +void FilterJetConeHits::end(){ + + std::cout << "FilterJetConeHits::end() " << name() + << " processed " << _nEvt << " events in " << _nRun << " runs " + << std::endl ; + +} + +// jet quality selection +bool FilterJetConeHits::filterJetBib(ReconstructedParticle* jet){ + + auto daughters = jet->getParticles(); + + double maxPt = 0., pt = 0.; + int ntracks = 0; + + for(unsigned int i = 0; i < daughters.size(); i++){ + + if(daughters[i]->getCharge() == 0) continue; //only care about tracks + ntracks++; + pt = sqrt( pow( daughters[i]->getMomentum()[0] ,2) + pow( daughters[i]->getMomentum()[1] ,2) ); + if(pt > maxPt) maxPt = pt; + } + + if(ntracks < m_minNTracks || maxPt < m_minDaughterMaxPt) return false; //check quality of jet, if bad reject + + streamlog_out(DEBUG) << "Accepted jet with ntracks = " << ntracks << " and max daughter pT = " << maxPt << std::endl; + return true; +}