#include <AlignmentPrescaler.h>
Public Member Functions | |
AlignmentPrescaler (const edm::ParameterSet &iConfig) | |
void | beginJob () |
void | endJob () |
virtual void | produce (edm::Event &iEvent, const edm::EventSetup &iSetup) |
~AlignmentPrescaler () | |
Private Member Functions | |
int | layerFromId (const DetId &id) const |
Private Attributes | |
unsigned int | detid_ |
TFile * | fpresc_ |
float | hitPrescFactor_ |
TRandom3 * | myrand_ |
float | overlapPrescFactor_ |
std::string | prescfilename_ |
std::string | presctreename_ |
edm::InputTag | src_ |
edm::InputTag | srcQualityMap_ |
int | totnhitspxl_ |
TTree * | tpresc_ |
Definition at line 23 of file AlignmentPrescaler.h.
AlignmentPrescaler::AlignmentPrescaler | ( | const edm::ParameterSet & | iConfig | ) |
Definition at line 29 of file AlignmentPrescaler.cc.
: src_(iConfig.getParameter<edm::InputTag>("src")), srcQualityMap_(iConfig.getParameter<edm::InputTag>("assomap")), prescfilename_(iConfig.getParameter<std::string>("PrescFileName")), presctreename_(iConfig.getParameter<std::string>("PrescTreeName")) { // issue the produce<> produces<AliClusterValueMap>(); produces<AliTrackTakenClusterValueMap>(); }
AlignmentPrescaler::~AlignmentPrescaler | ( | ) |
Definition at line 41 of file AlignmentPrescaler.cc.
{
//
}
void AlignmentPrescaler::beginJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDProducer.
Definition at line 45 of file AlignmentPrescaler.cc.
References gather_cfg::cout, detid_, fpresc_, hitPrescFactor_, myrand_, overlapPrescFactor_, prescfilename_, presctreename_, and tpresc_.
{ // std::cout<<"in AlignmentPrescaler::beginJob"<<std::flush; fpresc_=new TFile(prescfilename_.c_str(),"READ"); tpresc_=(TTree*)fpresc_->Get(presctreename_.c_str()); tpresc_->BuildIndex("DetId"); tpresc_->SetBranchStatus("*",0); tpresc_->SetBranchStatus("DetId",1); tpresc_->SetBranchStatus("PrescaleFactor",1); tpresc_->SetBranchStatus("PrescaleFactorOverlap",1); cout<<" Branches activated "<<std::flush; detid_=0; hitPrescFactor_=99.0; overlapPrescFactor_=88.0; tpresc_->SetBranchAddress("DetId",&detid_); tpresc_->SetBranchAddress("PrescaleFactor",&hitPrescFactor_); tpresc_->SetBranchAddress("PrescaleFactorOverlap",&overlapPrescFactor_); cout<<" addressed "<<std::flush; myrand_=new TRandom3(); // myrand_->SetSeed(); cout<<" ok "<<std::endl; }
void AlignmentPrescaler::endJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDProducer.
Definition at line 70 of file AlignmentPrescaler.cc.
int AlignmentPrescaler::layerFromId | ( | const DetId & | id | ) | const [private] |
Definition at line 243 of file AlignmentPrescaler.cc.
References PXFDetId::disk(), TIBDetId::layer(), TOBDetId::layer(), PXBDetId::layer(), PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, PXFDetId::side(), TIDDetId::side(), TECDetId::side(), StripSubdetector::TEC, sistripsummary::TIB, sistripsummary::TID, StripSubdetector::TOB, TIDDetId::wheel(), and TECDetId::wheel().
{ if ( uint32_t(id.subdetId())==PixelSubdetector::PixelBarrel ) { PXBDetId tobId(id); return tobId.layer(); } else if ( uint32_t(id.subdetId())==PixelSubdetector::PixelEndcap ) { PXFDetId tobId(id); return tobId.disk() + (3*(tobId.side()-1)); } else if ( id.subdetId()==StripSubdetector::TIB ) { TIBDetId tibId(id); return tibId.layer(); } else if ( id.subdetId()==StripSubdetector::TOB ) { TOBDetId tobId(id); return tobId.layer(); } else if ( id.subdetId()==StripSubdetector::TEC ) { TECDetId tobId(id); return tobId.wheel() + (9*(tobId.side()-1)); } else if ( id.subdetId()==StripSubdetector::TID ) { TIDDetId tobId(id); return tobId.wheel() + (3*(tobId.side()-1)); } return -1; }//end layerfromId
void AlignmentPrescaler::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 78 of file AlignmentPrescaler.cc.
References SiStripRecHit1D::cluster(), SiStripRecHit2D::cluster(), SiPixelRecHit::cluster(), gather_cfg::cout, edm::helper::Filler< Map >::fill(), TrackingRecHit::geographicalId(), edm::Event::getByLabel(), hitPrescFactor_, edm::helper::Filler< Map >::insert(), TrackingRecHit::isValid(), myrand_, overlapPrescFactor_, edm::Event::put(), DetId::rawId(), src_, srcQualityMap_, DetId::subdetId(), and tpresc_.
{ // std::cout<<"\n\n#################\n### Starting the AlignmentPrescaler::produce ; Event: "<<iEvent.id().run() <<", "<<iEvent.id().event()<<std::endl; edm::Handle<reco::TrackCollection> Tracks; iEvent.getByLabel(src_, Tracks); //take HitAssomap edm::Handle<AliClusterValueMap> hMap; iEvent.getByLabel(srcQualityMap_, hMap); AliClusterValueMap InValMap=*hMap; //prepare the output of the ValueMap flagging tracks std::vector<int> trackflags(Tracks->size(),0); //int npxlhits=0; //loop on tracks for(std::vector<reco::Track>::const_iterator ittrk = Tracks->begin(), edtrk = Tracks->end(); ittrk != edtrk; ++ittrk){ //loop on tracking rechits // std::cout << "Loop on hits of track #" << (ittrk - Tracks->begin()) << std::endl; int nhit=0; int ntakenhits=0; bool firstTakenHit=false; for (trackingRecHit_iterator ith = ittrk->recHitsBegin(), edh = ittrk->recHitsEnd(); ith != edh; ++ith) { const TrackingRecHit *hit = ith->get(); // ith is an iterator on edm::Ref to rechit if(! hit->isValid()){ nhit++; continue; } uint32_t tmpdetid = hit->geographicalId().rawId(); tpresc_->GetEntryWithIndex(tmpdetid); //------------- //decide whether to take this hit or not bool takeit=false; int subdetId=hit->geographicalId().subdetId(); //check first if the cluster is also in the overlap asso map bool isOverlapHit=false; // bool first=true; //ugly... const SiPixelRecHit* pixelhit= dynamic_cast<const SiPixelRecHit*>(hit); const SiStripRecHit1D* stripHit1D = dynamic_cast<const SiStripRecHit1D*>(hit); const SiStripRecHit2D* stripHit2D = dynamic_cast<const SiStripRecHit2D*>(hit); AlignmentClusterFlag tmpflag(hit->geographicalId()); int stripType=0; if(subdetId>2){// SST case const std::type_info &type = typeid(*hit); if (type == typeid(SiStripRecHit1D)) stripType=1; else if (type == typeid(SiStripRecHit2D)) stripType=2; else stripType=3; if(stripType==1) { // const SiStripRecHit1D* stripHit1D = dynamic_cast<const SiStripRecHit1D*>(hit); if(stripHit1D!=0){ SiStripRecHit1D::ClusterRef stripclust(stripHit1D->cluster()); tmpflag=InValMap[stripclust]; tmpflag.SetDetId(hit->geographicalId()); if(tmpflag.isOverlap())isOverlapHit=true; // std::cout<<"~*~*~* Prescale (1D) for module "<<tmpflag.detId().rawId()<<"("<<InValMap[stripclust].detId().rawId() <<") is "<<hitPrescFactor_<<std::flush; // if(tmpflag.isOverlap())cout<<" (it is Overlap)"<<endl; // else cout<<endl; }//end if striphit1D!=0 } else if (stripType==2) { //const SiStripRecHit2D* stripHit2D = dynamic_cast<const SiStripRecHit2D*>(hit); if(stripHit2D!=0){ SiStripRecHit2D::ClusterRef stripclust(stripHit2D->cluster()); tmpflag=InValMap[stripclust]; tmpflag.SetDetId(hit->geographicalId()); if(tmpflag.isOverlap())isOverlapHit=true; // std::cout<<"~*~*~* Prescale (2D) for module "<<tmpflag.detId().rawId()<<"("<<InValMap[stripclust].detId().rawId() <<") is "<<hitPrescFactor_<<std::flush; // if(tmpflag.isOverlap())cout<<" (it is Overlap)"<<endl; // else cout<<endl; }//end if striphit2D!=0 } }//end if is a strip hit else{ // const SiPixelRecHit* pixelhit= dynamic_cast<const SiPixelRecHit*>(hit); if(pixelhit!=0){ //npxlhits++; SiPixelClusterRefNew pixclust(pixelhit->cluster()); tmpflag=InValMap[pixclust]; tmpflag.SetDetId(hit->geographicalId()); if(tmpflag.isOverlap())isOverlapHit=true; } }//end else is a pixel hit // tmpflag.SetDetId(hit->geographicalId()); if( isOverlapHit ){ //cout<<" DetId="<<tmpdetid<<" is Overlap! "<<flush; takeit=(float(myrand_->Rndm())<=overlapPrescFactor_); } if( !takeit ){ float rr=float(myrand_->Rndm()); takeit=(rr<=hitPrescFactor_); } if(takeit){//HIT TAKEN ! //cout<<" DetId="<<tmpdetid<<" taken!"<<flush; tmpflag.SetTakenFlag(); if(subdetId>2){ if(stripType==1){ SiStripRecHit1D::ClusterRef stripclust(stripHit1D->cluster()); InValMap[stripclust]=tmpflag;//.SetTakenFlag(); } else if(stripType==2){ SiStripRecHit1D::ClusterRef stripclust(stripHit2D->cluster()); InValMap[stripclust]=tmpflag;//.SetTakenFlag(); } else std::cout<<"Unknown type of strip hit"<<std::endl; } else{ SiPixelClusterRefNew pixclust(pixelhit->cluster()); InValMap[pixclust]=tmpflag;//.SetTakenFlag(); } if(!firstTakenHit){ firstTakenHit=true; //std::cout<<"Index of the track iterator is "<< ittrk-Tracks->begin() <<endl; } ntakenhits++; }//end if take this hit //cout<<endl; nhit++; //cout<<endl; }//end loop on RecHits trackflags[ittrk-Tracks->begin()]=ntakenhits; }//end loop on tracks // totnhitspxl_+=ntakenhits; //cout<<"AlignmentPrescaler::produce says that in this event "<<ntakenhits<<" pixel clusters were taken (out of "<<npxlhits<<" total pixel hits."<<endl; //save the asso map, tracks... // prepare output std::auto_ptr<AliClusterValueMap> OutVM( new AliClusterValueMap); *OutVM=InValMap; iEvent.put(OutVM); std::auto_ptr<AliTrackTakenClusterValueMap> trkVM( new AliTrackTakenClusterValueMap); AliTrackTakenClusterValueMap::Filler trkmapfiller(*trkVM); trkmapfiller.insert(Tracks,trackflags.begin(),trackflags.end() ); trkmapfiller.fill(); iEvent.put(trkVM); }//end produce
unsigned int AlignmentPrescaler::detid_ [private] |
Definition at line 46 of file AlignmentPrescaler.h.
Referenced by beginJob().
TFile* AlignmentPrescaler::fpresc_ [private] |
Definition at line 39 of file AlignmentPrescaler.h.
Referenced by beginJob(), and endJob().
float AlignmentPrescaler::hitPrescFactor_ [private] |
Definition at line 47 of file AlignmentPrescaler.h.
Referenced by beginJob(), and produce().
TRandom3* AlignmentPrescaler::myrand_ [private] |
Definition at line 41 of file AlignmentPrescaler.h.
Referenced by beginJob(), endJob(), and produce().
float AlignmentPrescaler::overlapPrescFactor_ [private] |
Definition at line 47 of file AlignmentPrescaler.h.
Referenced by beginJob(), and produce().
std::string AlignmentPrescaler::prescfilename_ [private] |
Definition at line 36 of file AlignmentPrescaler.h.
Referenced by beginJob().
std::string AlignmentPrescaler::presctreename_ [private] |
Definition at line 37 of file AlignmentPrescaler.h.
Referenced by beginJob().
edm::InputTag AlignmentPrescaler::src_ [private] |
Definition at line 33 of file AlignmentPrescaler.h.
Referenced by produce().
Definition at line 34 of file AlignmentPrescaler.h.
Referenced by produce().
int AlignmentPrescaler::totnhitspxl_ [private] |
Definition at line 48 of file AlignmentPrescaler.h.
TTree* AlignmentPrescaler::tpresc_ [private] |
Definition at line 40 of file AlignmentPrescaler.h.
Referenced by beginJob(), endJob(), and produce().