CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/Alignment/TrackerAlignment/plugins/AlignmentPrescaler.cc

Go to the documentation of this file.
00001 #include "Alignment/TrackerAlignment/plugins/AlignmentPrescaler.h"
00002 
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "DataFormats/Common/interface/View.h"
00005 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00006 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00007 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00008 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00009 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00010 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00011 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00012 #include "Alignment/TrackerAlignment/interface/AlignableTracker.h"
00013 #include "Alignment/CommonAlignment/interface/Alignable.h"
00014 #include "Alignment/CommonAlignment/interface/Utilities.h"
00015 #include "Utilities/General/interface/ClassName.h"
00016 
00017 #include "DataFormats/TrackReco/interface/Track.h"
00018 #include "AnalysisDataFormats/TrackInfo/interface/TrackInfo.h"
00019 
00020 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00021 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00022 #include "RecoTracker/TransientTrackingRecHit/interface/TSiStripRecHit2DLocalPos.h"
00023 #include "RecoTracker/TransientTrackingRecHit/interface/TSiPixelRecHit.h"
00024 #include "DataFormats/Alignment/interface/AlignmentClusterFlag.h"
00025 #include "DataFormats/Alignment/interface/AliClusterValueMap.h"
00026 
00027 #include "Riostream.h"
00028 
00029 AlignmentPrescaler::AlignmentPrescaler(const edm::ParameterSet &iConfig):
00030   src_(iConfig.getParameter<edm::InputTag>("src")),
00031   srcQualityMap_(iConfig.getParameter<edm::InputTag>("assomap")),
00032   prescfilename_(iConfig.getParameter<std::string>("PrescFileName")),
00033   presctreename_(iConfig.getParameter<std::string>("PrescTreeName"))
00034 {
00035   // issue the produce<>
00036   produces<AliClusterValueMap>();
00037   produces<AliTrackTakenClusterValueMap>();
00038 
00039 }
00040 
00041 AlignmentPrescaler::~AlignmentPrescaler(){
00042   //  
00043 }
00044 
00045 void AlignmentPrescaler::beginJob(){
00046   //
00047   std::cout<<"in AlignmentPrescaler::beginJob"<<std::flush;
00048    fpresc_=new TFile(prescfilename_.c_str(),"READ");
00049    tpresc_=(TTree*)fpresc_->Get(presctreename_.c_str());
00050    tpresc_->BuildIndex("DetId");
00051    tpresc_->SetBranchStatus("*",0);
00052    tpresc_->SetBranchStatus("DetId",1);
00053    tpresc_->SetBranchStatus("PrescaleFactor",1);
00054    tpresc_->SetBranchStatus("PrescaleFactorOverlap",1);
00055    cout<<" Branches activated "<<std::flush;
00056    detid_=0;
00057    hitPrescFactor_=99.0;
00058    overlapPrescFactor_=88.0;
00059    
00060    tpresc_->SetBranchAddress("DetId",&detid_);
00061    tpresc_->SetBranchAddress("PrescaleFactor",&hitPrescFactor_);
00062    tpresc_->SetBranchAddress("PrescaleFactorOverlap",&overlapPrescFactor_);
00063    cout<<" addressed "<<std::flush;   
00064    myrand_=new TRandom3();
00065    //   myrand_->SetSeed();
00066    cout<<" ok "<<std::endl;
00067 
00068 }
00069 
00070 void AlignmentPrescaler::endJob( ){
00071 
00072   delete tpresc_;
00073   fpresc_->Close();
00074   delete fpresc_;
00075   delete myrand_;
00076 }
00077 
00078 void AlignmentPrescaler::produce(edm::Event &iEvent, const edm::EventSetup &iSetup){
00079   // std::cout<<"\n\n#################\n### Starting the AlignmentPrescaler::produce ; Event: "<<iEvent.id().run() <<", "<<iEvent.id().event()<<std::endl;
00080   edm::Handle<reco::TrackCollection> Tracks;
00081   iEvent.getByLabel(src_, Tracks);
00082  
00083   //take  HitAssomap
00084   edm::Handle<AliClusterValueMap> hMap;
00085   iEvent.getByLabel(srcQualityMap_, hMap);
00086   AliClusterValueMap InValMap=*hMap;
00087 
00088   //prepare the output of the ValueMap flagging tracks
00089   std::vector<int> trackflags(Tracks->size(),0);
00090 
00091 
00092   //int npxlhits=0;
00093   
00094     //loop on tracks
00095   for(std::vector<reco::Track>::const_iterator ittrk = Tracks->begin(), edtrk = Tracks->end(); ittrk != edtrk; ++ittrk){
00096     //loop on tracking rechits
00097     // std::cout << "Loop on hits of track #" << (ittrk - Tracks->begin()) << std::endl;
00098     int nhit=0;
00099     int ntakenhits=0;
00100     bool firstTakenHit=false;
00101 
00102     for (trackingRecHit_iterator ith = ittrk->recHitsBegin(), edh = ittrk->recHitsEnd(); ith != edh; ++ith) {
00103       const TrackingRecHit *hit = ith->get(); // ith is an iterator on edm::Ref to rechit
00104       if(! hit->isValid()){
00105         nhit++;
00106         continue;
00107       }
00108       uint32_t tmpdetid = hit->geographicalId().rawId();
00109       tpresc_->GetEntryWithIndex(tmpdetid);
00110       
00111 
00112       //-------------
00113       //decide whether to take this hit or not
00114       bool takeit=false;  
00115       int subdetId=hit->geographicalId().subdetId();   
00116  
00117 
00118       //check first if the cluster is also in the overlap asso map
00119       bool isOverlapHit=false;
00120       //  bool first=true;
00121       //ugly...
00122       const SiPixelRecHit*   pixelhit= dynamic_cast<const SiPixelRecHit*>(hit);
00123       const SiStripRecHit1D* stripHit1D = dynamic_cast<const SiStripRecHit1D*>(hit);
00124       const SiStripRecHit2D* stripHit2D = dynamic_cast<const SiStripRecHit2D*>(hit);
00125 
00126       AlignmentClusterFlag tmpflag(hit->geographicalId());
00127       int stripType=0;
00128       if(subdetId>2){// SST case
00129         const std::type_info &type = typeid(*hit);       
00130         if (type == typeid(SiStripRecHit1D))    stripType=1;
00131         else  if (type == typeid(SiStripRecHit2D))      stripType=2;
00132         else    stripType=3;
00133 
00134         if(stripType==1) { 
00135           //      const SiStripRecHit1D* stripHit1D = dynamic_cast<const SiStripRecHit1D*>(hit);
00136           
00137           if(stripHit1D!=0){
00138             SiStripRecHit1D::ClusterRef stripclust(stripHit1D->cluster());
00139             tmpflag=InValMap[stripclust];
00140             tmpflag.SetDetId(hit->geographicalId());
00141             if(tmpflag.isOverlap())isOverlapHit=true;
00142             // std::cout<<"~*~*~* Prescale (1D) for module "<<tmpflag.detId().rawId()<<"("<<InValMap[stripclust].detId().rawId() <<") is "<<hitPrescFactor_<<std::flush;
00143             //  if(tmpflag.isOverlap())cout<<" (it is Overlap)"<<endl;
00144             // else cout<<endl;
00145             
00146           }//end if striphit1D!=0
00147         }
00148         else if (stripType==2) {
00149           //const SiStripRecHit2D* stripHit2D = dynamic_cast<const SiStripRecHit2D*>(hit);
00150           if(stripHit2D!=0){
00151             SiStripRecHit2D::ClusterRef stripclust(stripHit2D->cluster());
00152             tmpflag=InValMap[stripclust];
00153             tmpflag.SetDetId(hit->geographicalId());
00154             if(tmpflag.isOverlap())isOverlapHit=true;
00155             // std::cout<<"~*~*~* Prescale (2D) for module "<<tmpflag.detId().rawId()<<"("<<InValMap[stripclust].detId().rawId() <<") is "<<hitPrescFactor_<<std::flush;
00156             //  if(tmpflag.isOverlap())cout<<" (it is Overlap)"<<endl;
00157             // else cout<<endl;
00158           
00159           }//end if striphit2D!=0
00160         }
00161       }//end if is a strip hit
00162       else{
00163         //      const SiPixelRecHit*   pixelhit= dynamic_cast<const SiPixelRecHit*>(hit);
00164         if(pixelhit!=0){
00165           //npxlhits++;
00166           SiPixelClusterRefNew pixclust(pixelhit->cluster());
00167           tmpflag=InValMap[pixclust];
00168           tmpflag.SetDetId(hit->geographicalId());
00169           if(tmpflag.isOverlap())isOverlapHit=true;
00170         }
00171       }//end else is a pixel hit
00172       //      tmpflag.SetDetId(hit->geographicalId());
00173 
00174       if( isOverlapHit ){
00175         //cout<<"  DetId="<<tmpdetid<<" is Overlap! "<<flush;
00176         takeit=(float(myrand_->Rndm())<=overlapPrescFactor_);
00177       }
00178       if( !takeit ){
00179         float rr=float(myrand_->Rndm());
00180         takeit=(rr<=hitPrescFactor_);
00181       }
00182       if(takeit){//HIT TAKEN !
00183         //cout<<"  DetId="<<tmpdetid<<" taken!"<<flush;
00184         tmpflag.SetTakenFlag();
00185 
00186         if(subdetId>2){
00187           if(stripType==1){
00188             SiStripRecHit1D::ClusterRef stripclust(stripHit1D->cluster());
00189             InValMap[stripclust]=tmpflag;//.SetTakenFlag();
00190           }
00191           else if(stripType==2){
00192             SiStripRecHit1D::ClusterRef stripclust(stripHit2D->cluster());
00193             InValMap[stripclust]=tmpflag;//.SetTakenFlag();
00194           }
00195           else std::cout<<"Unknown type of strip hit"<<std::endl;
00196         }
00197         else{
00198           SiPixelClusterRefNew pixclust(pixelhit->cluster());
00199           InValMap[pixclust]=tmpflag;//.SetTakenFlag();
00200         }
00201         
00202         if(!firstTakenHit){
00203           firstTakenHit=true;
00204           //std::cout<<"Index of the track iterator is "<< ittrk-Tracks->begin() <<endl;
00205           
00206         }
00207         ntakenhits++;
00208       }//end if take this hit
00209       //cout<<endl;
00210 
00211         nhit++;
00212       //cout<<endl;
00213     }//end loop on RecHits
00214     trackflags[ittrk-Tracks->begin()]=ntakenhits;
00215   
00216   }//end loop on tracks
00217   
00218 
00219 
00220   // totnhitspxl_+=ntakenhits;
00221   //cout<<"AlignmentPrescaler::produce says that in this event "<<ntakenhits<<" pixel clusters were taken (out of "<<npxlhits<<" total pixel hits."<<endl;
00222 
00223 
00224 
00225   //save the asso map, tracks...
00226   // prepare output 
00227   std::auto_ptr<AliClusterValueMap> OutVM( new AliClusterValueMap);
00228   *OutVM=InValMap;
00229 
00230   iEvent.put(OutVM);
00231   
00232   
00233   std::auto_ptr<AliTrackTakenClusterValueMap> trkVM( new AliTrackTakenClusterValueMap);
00234   AliTrackTakenClusterValueMap::Filler trkmapfiller(*trkVM);
00235   trkmapfiller.insert(Tracks,trackflags.begin(),trackflags.end() );
00236   trkmapfiller.fill();
00237   iEvent.put(trkVM);
00238 
00239 
00240 }//end produce
00241 
00242 
00243 int AlignmentPrescaler::layerFromId (const DetId& id) const
00244 {
00245  if ( uint32_t(id.subdetId())==PixelSubdetector::PixelBarrel ) {
00246     PXBDetId tobId(id);
00247     return tobId.layer();
00248   }
00249   else if ( uint32_t(id.subdetId())==PixelSubdetector::PixelEndcap ) {
00250     PXFDetId tobId(id);
00251     return tobId.disk() + (3*(tobId.side()-1));
00252   }
00253   else if ( id.subdetId()==StripSubdetector::TIB ) {
00254     TIBDetId tibId(id);
00255     return tibId.layer();
00256   }
00257   else if ( id.subdetId()==StripSubdetector::TOB ) {
00258     TOBDetId tobId(id);
00259     return tobId.layer();
00260   }
00261   else if ( id.subdetId()==StripSubdetector::TEC ) {
00262     TECDetId tobId(id);
00263     return tobId.wheel() + (9*(tobId.side()-1));
00264   }
00265   else if ( id.subdetId()==StripSubdetector::TID ) {
00266     TIDDetId tobId(id);
00267     return tobId.wheel() + (3*(tobId.side()-1));
00268   }
00269   return -1;
00270 
00271 }//end layerfromId
00272 
00273 // ========= MODULE DEF ==============
00274 #include "FWCore/PluginManager/interface/ModuleDef.h"
00275 #include "FWCore/Framework/interface/MakerMacros.h"
00276 DEFINE_FWK_MODULE(AlignmentPrescaler);