CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/Alignment/TrackerAlignment/plugins/TkAlCaSkimTreeMerger.cc

Go to the documentation of this file.
00001 
00002 #ifndef TrackerAlignment_TkAlCaSkimTreeMerger_H
00003 #define TrackerAlignment_TkAlCaSkimTreeMerger_H
00004 
00005 #include <Riostream.h>
00006 #include <string>
00007 #include <fstream>
00008 #include <map>
00009 
00010 #include "FWCore/Framework/interface/EDAnalyzer.h"
00011 #include "FWCore/Framework/interface/EventPrincipal.h" 
00012 #include "FWCore/Framework/interface/Event.h"
00013 #include "FWCore/Framework/interface/EventSetup.h"
00014 #include "FWCore/Framework/interface/ESHandle.h"
00015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "FWCore/Utilities/interface/InputTag.h"
00018 
00019 
00020 #include "TFile.h"
00021 #include "TString.h"
00022 #include "TChain.h"
00023 #include "TStopwatch.h"
00024 
00025 class TkAlCaSkimTreeMerger : public edm::EDAnalyzer{
00026 
00027  public:
00028   TkAlCaSkimTreeMerger(const edm::ParameterSet &iConfig);
00029   ~TkAlCaSkimTreeMerger();
00030   void beginJob();
00031   void endJob();
00032   void analyze(const edm::Event&, const edm::EventSetup&);
00033 
00034  private:
00035   TTree *out_;//TTree containing the merged result
00036   TTree *firsttree_;//first tree of the list; this gives the structure to all the others 
00037   TChain *ch_;//chain containing all the tree you want to merge
00038   std::string filelist_;//text file containing the list of input files
00039   std::string firstfilename_;
00040   std::string treename_;//name of the tree you want to merge (contained in the file)
00041   std::string outfilename_;//name of the file where you want to save the output
00042  
00043   //Hit Population
00044   typedef map<uint32_t,uint32_t>DetHitMap;
00045   DetHitMap hitmap_;
00046   DetHitMap overlapmap_;
00047   int maxhits_;//above this number, the hit population is prescaled. Configurable for each subdet 
00048   edm::ParameterSet maxhitsSet_;
00049   int maxPXBhits_, maxPXFhits_, maxTIBhits_, maxTIDhits_, maxTOBhits_, maxTEChits_;
00050  
00051 
00052   TStopwatch myclock;
00053 
00054 };
00055 
00056 
00057 
00058 
00059 #endif
00060 
00061 
00062 TkAlCaSkimTreeMerger::TkAlCaSkimTreeMerger(const edm::ParameterSet &iConfig) :
00063     filelist_(iConfig.getParameter<string>("FileList")), 
00064     treename_(iConfig.getParameter<string>("TreeName")),
00065     outfilename_(iConfig.getParameter<string>("OutputFile")),
00066     // maxhits_(iConfig.getParameter<vint>("NhitsMaxLimit"))
00067     maxhits_(iConfig.getParameter<int32_t>("NhitsMaxLimit")),
00068     maxhitsSet_(iConfig.getParameter<edm::ParameterSet>("NhitsMaxSet"))
00069 {
00070   maxPXBhits_=maxhitsSet_.getParameter<int32_t>("PXBmaxhits");
00071   maxPXFhits_=maxhitsSet_.getParameter<int32_t>("PXFmaxhits");
00072   maxTIBhits_=maxhitsSet_.getParameter<int32_t>("TIBmaxhits");
00073   maxTIDhits_=maxhitsSet_.getParameter<int32_t>("TIDmaxhits");
00074   maxTOBhits_=maxhitsSet_.getParameter<int32_t>("TOBmaxhits");
00075   maxTEChits_=maxhitsSet_.getParameter<int32_t>("TECmaxhits");
00076   //anything you want to do for initializing
00077   cout<<"\n\n*** MAX N HITS = "<<maxhits_<<endl<<endl;
00078   out_=0;
00079   firsttree_=0;
00080   ch_=0;
00081 }
00082 
00083 
00084 TkAlCaSkimTreeMerger::~TkAlCaSkimTreeMerger(){
00085   //default destructor
00086   // delete out_;
00087   // delete firsttree_;
00088 
00089   delete ch_;
00090   cout<<"finished."<<endl;
00091 }
00092 
00093 // ------------ method called before analyzing the first event  ------------
00094 void TkAlCaSkimTreeMerger::beginJob(){
00095  
00096   myclock.Start();
00097 
00098   //prepare the chain
00099   ch_=new TChain(treename_.c_str());
00100   cout<<"The chain contains "<<ch_->GetNtrees()<<" trees"<<endl;
00101   
00102   //load the trees into the chain
00103   ifstream flist(filelist_.c_str(),ios::in);
00104   std::string filename;
00105   std::string firstfilename;
00106   bool first=true;
00107   while(!flist.eof()){
00108     filename="";
00109     flist>>filename;
00110     if(filename.empty())continue;
00111     //cout<<"Adding "<<filename<<endl;
00112     ch_->Add(filename.c_str());
00113     if(first){
00114       firstfilename_=filename;
00115       first=false;
00116     }
00117    
00118   }
00119   cout<<"Now the chain contains "<<ch_->GetNtrees()<<" trees ("<<ch_->GetEntries()<<" entries)"<<endl;
00120 
00121  
00122   unsigned int id_ch=0;
00123   uint32_t nhits_ch=0,noverlaps_ch=0;
00124   ch_->SetBranchAddress("DetId",    &id_ch);
00125   ch_->SetBranchAddress("Nhits",    &nhits_ch);
00126   ch_->SetBranchAddress("Noverlaps",&noverlaps_ch);
00127 
00128   ch_->SetBranchStatus("SubDet",0);
00129   ch_->SetBranchStatus("Layer",0);
00130   ch_->SetBranchStatus("is2D",0);
00131   ch_->SetBranchStatus("isStereo",0);
00132   ch_->SetBranchStatus("posX",0);
00133   ch_->SetBranchStatus("posY",0);
00134   ch_->SetBranchStatus("posZ",0);
00135   ch_->SetBranchStatus("posR",0);
00136   ch_->SetBranchStatus("posEta",0);
00137   ch_->SetBranchStatus("posPhi",0);//now only id, nhits and noverlaps are on...
00138 
00139 
00140   int totnhits(0),totnoverlaps(0);
00141 
00142 
00143   //look if you find this detid in the map
00144   DetHitMap::iterator mapiter;
00145   DetHitMap::iterator overlapiter;
00146   
00147   for(int ent=0;ent<ch_->GetEntries();++ent){
00148   //  for(int ent=0;ent<100;++ent){
00149     ch_->GetEntry(ent);
00150     totnhits+=nhits_ch;
00151     totnoverlaps+=noverlaps_ch;
00152 
00153     mapiter= hitmap_.find(id_ch);
00154     if(mapiter!=hitmap_.end() ){//present, increase its value
00155       hitmap_[id_ch]=hitmap_[id_ch]+nhits_ch;
00156     }
00157     else{//not present, let's add this key to the map with value=1
00158       hitmap_.insert(pair<uint32_t, uint32_t>(id_ch, nhits_ch));
00159     }
00160     //do the same for overlaps
00161     overlapiter= overlapmap_.find(id_ch);
00162     if(overlapiter!=overlapmap_.end() ){
00163       overlapmap_[id_ch]=overlapmap_[id_ch]+noverlaps_ch;
00164     }
00165     else{
00166       overlapmap_.insert(pair<uint32_t, uint32_t>(id_ch, noverlaps_ch));
00167     }
00168     
00169   }//end loop on ent - entries in the chain
00170 
00171 
00172   cout<<"Nhits in the chain: "<<totnhits<<endl;
00173   cout<<"NOverlaps in the chain: "<<totnoverlaps<<endl;
00174 
00175 
00176   myclock.Stop();
00177   cout<<"Finished beginJob after "<<myclock.RealTime()<<" s (real time) / "<<myclock.CpuTime()<<" s (cpu time)"<<endl;
00178   myclock.Continue();
00179 }//end beginJob
00180 
00181 
00182 // ------------ method called to for each event  ------------
00183 void TkAlCaSkimTreeMerger::analyze(const edm::Event&, const edm::EventSetup&){
00184     // cout<<firsttree_->GetEntries()<<endl;
00185 }//end analyze
00186 
00187 // ------------ method called after having analyzed all the events  ------------
00188 void TkAlCaSkimTreeMerger::endJob(){
00189 
00190 
00191   //address variables in the first tree and in the chain
00192   TFile *firstfile=new TFile(firstfilename_.c_str(),"READ");
00193   firsttree_=(TTree*)firstfile->Get(treename_.c_str());
00194   cout<<"the first tree has "<<firsttree_->GetEntries() <<" entries"<<endl; 
00195   unsigned int id=0;
00196   uint32_t nhits=0,noverlaps=0;
00197   float posX(-99999.0),posY(-77777.0),posZ(-88888.0);
00198   float posEta(-6666.0),posPhi(-5555.0),posR(-4444.0);
00199   int subdet=0;
00200   unsigned int layer=0; 
00201   // bool is2D=false,isStereo=false;
00202   firsttree_->SetBranchAddress("DetId",    &id);
00203   firsttree_->SetBranchAddress("Nhits",    &nhits);
00204   firsttree_->SetBranchAddress("Noverlaps",&noverlaps);
00205   firsttree_->SetBranchAddress("SubDet",   &subdet);
00206   firsttree_->SetBranchAddress("Layer",    &layer);
00207   //  firsttree_->SetBranchAddress("is2D" ,    &is2D);
00208   // firsttree_->SetBranchAddress("isStereo", &isStereo);
00209   firsttree_->SetBranchAddress("posX",     &posX);
00210   firsttree_->SetBranchAddress("posY",     &posY);
00211   firsttree_->SetBranchAddress("posZ",     &posZ);
00212   firsttree_->SetBranchAddress("posR",     &posR);
00213   firsttree_->SetBranchAddress("posEta",   &posEta);
00214   firsttree_->SetBranchAddress("posPhi",   &posPhi);
00215 
00216 
00217   //create and book the output
00218  
00219  
00220   TFile *outfile=new TFile(outfilename_.c_str(),"RECREATE");
00221   out_=new TTree(treename_.c_str(),"AlignmentHitMapsTOTAL");
00222   unsigned int id_out=0;
00223   uint32_t nhits_out=0,noverlaps_out=0;
00224   float posX_out(-99999.0),posY_out(-77777.0),posZ_out(-88888.0);
00225   float posEta_out(-6666.0),posPhi_out(-5555.0),posR_out(-4444.0);
00226   int subdet_out=0;
00227   unsigned int layer_out=0; 
00228   bool is2D_out=false,isStereo_out=false;
00229   float prescfact_out=1.0;
00230   float prescfact_overlap_out=1.0;
00231 
00232   out_->Branch("DetId",    &id_out ,      "DetId/i");
00233   out_->Branch("Nhits",    &nhits_out ,   "Nhits/i");
00234   out_->Branch("Noverlaps",&noverlaps_out,"Noverlaps/i");
00235   out_->Branch("SubDet",   &subdet_out,   "SubDet/I");
00236   out_->Branch("Layer",    &layer_out,    "Layer/i");
00237   out_->Branch("is2D" ,    &is2D_out,     "is2D/B");
00238   out_->Branch("isStereo", &isStereo_out, "isStereo/B");
00239   out_->Branch("posX",     &posX_out,     "posX/F");
00240   out_->Branch("posY",     &posY_out,     "posY/F");
00241   out_->Branch("posZ",     &posZ_out,     "posZ/F");
00242   out_->Branch("posR",     &posR_out,     "posR/F");
00243   out_->Branch("posEta",   &posEta_out,   "posEta/F");
00244   out_->Branch("posPhi",   &posPhi_out,   "posPhi/F");
00245   out_->Branch("PrescaleFactor",&prescfact_out,"PrescaleFact/F");  
00246   out_->Branch("PrescaleFactorOverlap",&prescfact_overlap_out,"PrescaleFactOverlap/F"); 
00247 
00248 
00249   //look if you find this detid in the map
00250   DetHitMap::iterator mapiter;
00251   
00252    for(int mod=0;mod<firsttree_->GetEntries();mod++){
00253      //for(int mod=0;mod<100;++mod){
00254     //   nhits_out=0;
00255     // noverlaps_out=0;
00256  
00257     firsttree_->GetEntry(mod);
00258     nhits_out=hitmap_[id];
00259     noverlaps_out=overlapmap_[id];
00260     // if(mod<25)cout<<"Nhits 1st tree: "<<nhits<<"\tTotal nhits chain: "<<nhits_out<<endl;
00261     id_out=id;
00262     subdet_out=subdet;
00263     layer_out=layer;
00264     posX_out=posX;
00265     posY_out=posY;
00266     posZ_out=posZ;
00267     posR_out=posR;
00268     posEta_out=posEta;
00269     posPhi_out=posPhi;
00270 
00271     //calculate prescaling factor
00272     int subdetmax=-1;
00273     if(subdet_out==1)subdetmax=maxPXBhits_;
00274     else if(subdet_out==2)subdetmax=maxPXFhits_;
00275     else if(subdet_out==3)subdetmax=maxTIBhits_;
00276     else if(subdet_out==4)subdetmax=maxTIDhits_;
00277     else if(subdet_out==5)subdetmax=maxTOBhits_;
00278     else if(subdet_out==6)subdetmax=maxTEChits_;
00279     else subdetmax=-9999;
00280 
00281     if(maxhits_>-1){
00282       if(int(nhits_out)>maxhits_){
00283         prescfact_out=float(maxhits_)/float(nhits_out);
00284       }
00285       if(int(noverlaps_out)>maxhits_){
00286         prescfact_overlap_out=float(maxhits_)/float(noverlaps_out);
00287       }
00288     }
00289     else if(subdetmax>0){//calculate different prescaling factors for each subdet
00290       if(int(nhits_out)>subdetmax){
00291         prescfact_out=float(subdetmax/nhits_out);
00292       }
00293      if(int(noverlaps_out)>subdetmax){
00294         prescfact_overlap_out=float(subdetmax)/float(noverlaps_out);
00295       }
00296     }
00297     else{
00298       prescfact_out=1.0;
00299       prescfact_overlap_out=1.0;
00300     }
00301     out_->Fill();
00302 
00303   }//end loop on mod - first tree modules
00304 
00305 
00306   myclock.Stop();
00307   cout<<"Finished endJob after "<<myclock.RealTime()<<" s (real time) / "<<myclock.CpuTime()<<" s (cpu time)"<<endl;
00308   cout<<"Ending the tree merging."<<endl;
00309   out_->Write();
00310   cout<<"Deleting..."<<flush;
00311   delete firstfile;
00312   delete outfile;
00313  
00314   
00315 }//end endJob
00316 
00317 
00318 // ========= MODULE DEF ==============
00319 #include "FWCore/PluginManager/interface/ModuleDef.h"
00320 #include "FWCore/Framework/interface/MakerMacros.h"
00321 DEFINE_FWK_MODULE(TkAlCaSkimTreeMerger);
00322