CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/Validation/GlobalRecHits/src/GlobalRecHitsAnalyzer.cc

Go to the documentation of this file.
00001 
00010 using namespace std;
00011 #include "Validation/GlobalRecHits/interface/GlobalRecHitsAnalyzer.h"
00012 #include "DQMServices/Core/interface/DQMStore.h"
00013 #include "DQMServices/Core/interface/MonitorElement.h"
00014 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00015 
00016 GlobalRecHitsAnalyzer::GlobalRecHitsAnalyzer(const edm::ParameterSet& iPSet) :
00017   fName(""), verbosity(0), frequency(0), label(""), getAllProvenances(false),
00018   printProvenanceInfo(false), count(0)
00019 {
00020   std::string MsgLoggerCat = "GlobalRecHitsAnalyzer_GlobalRecHitsAnalyzer";
00021 
00022   // get information from parameter set
00023   fName = iPSet.getUntrackedParameter<std::string>("Name");
00024   verbosity = iPSet.getUntrackedParameter<int>("Verbosity");
00025   frequency = iPSet.getUntrackedParameter<int>("Frequency");
00026   edm::ParameterSet m_Prov =
00027     iPSet.getParameter<edm::ParameterSet>("ProvenanceLookup");
00028   getAllProvenances = 
00029     m_Prov.getUntrackedParameter<bool>("GetAllProvenances");
00030   printProvenanceInfo = 
00031     m_Prov.getUntrackedParameter<bool>("PrintProvenanceInfo");
00032   hitsProducer = iPSet.getParameter<std::string>("hitsProducer");
00033 
00034   //get Labels to use to extract information
00035   ECalEBSrc_ = iPSet.getParameter<edm::InputTag>("ECalEBSrc");
00036   ECalUncalEBSrc_ = iPSet.getParameter<edm::InputTag>("ECalUncalEBSrc");
00037   ECalEESrc_ = iPSet.getParameter<edm::InputTag>("ECalEESrc");
00038   ECalUncalEESrc_ = iPSet.getParameter<edm::InputTag>("ECalUncalEESrc");
00039   ECalESSrc_ = iPSet.getParameter<edm::InputTag>("ECalESSrc");
00040   HCalSrc_ = iPSet.getParameter<edm::InputTag>("HCalSrc");
00041   SiStripSrc_ = iPSet.getParameter<edm::InputTag>("SiStripSrc"); 
00042   SiPxlSrc_ = iPSet.getParameter<edm::InputTag>("SiPxlSrc");
00043   MuDTSrc_ = iPSet.getParameter<edm::InputTag>("MuDTSrc");
00044   MuDTSimSrc_ = iPSet.getParameter<edm::InputTag>("MuDTSimSrc");
00045   MuCSCSrc_ = iPSet.getParameter<edm::InputTag>("MuCSCSrc");
00046   MuRPCSrc_ = iPSet.getParameter<edm::InputTag>("MuRPCSrc");
00047   MuRPCSimSrc_ = iPSet.getParameter<edm::InputTag>("MuRPCSimSrc");
00048 
00049   conf_ = iPSet;
00050 
00051   // use value of first digit to determine default output level (inclusive)
00052   // 0 is none, 1 is basic, 2 is fill output, 3 is gather output
00053   verbosity %= 10;
00054 
00055   // create persistent object
00056   // produces<PGlobalRecHit>(label);
00057 
00058   // print out Parameter Set information being used
00059   if (verbosity >= 0) {
00060     edm::LogInfo(MsgLoggerCat) 
00061       << "\n===============================\n"
00062       << "Initialized as EDProducer with parameter values:\n"
00063       << "    Name           = " << fName << "\n"
00064       << "    Verbosity      = " << verbosity << "\n"
00065       << "    Frequency      = " << frequency << "\n"
00066       << "    GetProv        = " << getAllProvenances << "\n"
00067       << "    PrintProv      = " << printProvenanceInfo << "\n"
00068       << "    ECalEBSrc      = " << ECalEBSrc_.label() 
00069       << ":" << ECalEBSrc_.instance() << "\n"
00070       << "    ECalUncalEBSrc = " << ECalUncalEBSrc_.label() 
00071       << ":" << ECalUncalEBSrc_.instance() << "\n"
00072       << "    ECalEESrc      = " << ECalEESrc_.label() 
00073       << ":" << ECalUncalEESrc_.instance() << "\n"
00074       << "    ECalUncalEESrc = " << ECalUncalEESrc_.label() 
00075       << ":" << ECalEESrc_.instance() << "\n"
00076       << "    ECalESSrc      = " << ECalESSrc_.label() 
00077       << ":" << ECalESSrc_.instance() << "\n"
00078       << "    HCalSrc        = " << HCalSrc_.label() 
00079       << ":" << HCalSrc_.instance() << "\n"
00080       << "    SiStripSrc     = " << SiStripSrc_.label() 
00081       << ":" << SiStripSrc_.instance() << "\n" 
00082       << "    SiPixelSrc     = " << SiPxlSrc_.label()
00083       << ":" << SiPxlSrc_.instance() << "\n"
00084       << "    MuDTSrc        = " << MuDTSrc_.label()
00085       << ":" << MuDTSrc_.instance() << "\n"
00086       << "    MuDTSimSrc     = " << MuDTSimSrc_.label()
00087       << ":" << MuDTSimSrc_.instance() << "\n"
00088       << "    MuCSCSrc       = " << MuCSCSrc_.label()
00089       << ":" << MuCSCSrc_.instance() << "\n"
00090       << "    MuRPCSrc       = " << MuRPCSrc_.label()
00091       << ":" << MuRPCSrc_.instance() << "\n"
00092       << "    MuRPCSimSrc    = " << MuRPCSimSrc_.label()
00093       << ":" << MuRPCSimSrc_.instance() << "\n"
00094       << "===============================\n";
00095   }
00096   //Put in analyzer stuff here....
00097 
00098   dbe = 0;
00099   dbe = edm::Service<DQMStore>().operator->();
00100   if (dbe) {
00101     if (verbosity > 0 ) {
00102       dbe->setVerbose(1);
00103     } else {
00104       dbe->setVerbose(0);
00105     }
00106   }
00107   if (dbe) {
00108     if (verbosity > 0 ) dbe->showDirStructure();
00109   }
00110 
00111   //monitor elements 
00112   
00113   //Si Strip
00114   if(dbe) {
00115     string SiStripString[19] = {"TECW1", "TECW2", "TECW3", "TECW4", "TECW5", 
00116                                 "TECW6", "TECW7", "TECW8", "TIBL1", "TIBL2", 
00117                                 "TIBL3", "TIBL4", "TIDW1", "TIDW2", "TIDW3", 
00118                                 "TOBL1", "TOBL2", "TOBL3", "TOBL4"};
00119     for(int i = 0; i<19; ++i) {
00120       mehSiStripn[i]=0;
00121       mehSiStripResX[i]=0;
00122       mehSiStripResY[i]=0;
00123     }
00124     string hcharname, hchartitle;
00125     dbe->setCurrentFolder("GlobalRecHitsV/SiStrips");
00126     for(int amend = 0; amend < 19; ++amend) { 
00127       hcharname = "hSiStripn_"+SiStripString[amend];
00128       hchartitle= SiStripString[amend]+"  rechits";
00129       mehSiStripn[amend] = dbe->book1D(hcharname,hchartitle,200,0.,200.);
00130       mehSiStripn[amend]->setAxisTitle("Number of hits in "+
00131                                        SiStripString[amend],1);
00132       mehSiStripn[amend]->setAxisTitle("Count",2);
00133       hcharname = "hSiStripResX_"+SiStripString[amend];
00134       hchartitle= SiStripString[amend]+" rechit x resolution";
00135       mehSiStripResX[amend] = dbe->book1D(hcharname,hchartitle,200,-0.02,.02);
00136       mehSiStripResX[amend]->setAxisTitle("X-resolution in "
00137                                           +SiStripString[amend],1);
00138       mehSiStripResX[amend]->setAxisTitle("Count",2);
00139       hcharname = "hSiStripResY_"+SiStripString[amend];
00140       hchartitle= SiStripString[amend]+" rechit y resolution";
00141       mehSiStripResY[amend] = dbe->book1D(hcharname,hchartitle,200,-0.02,.02);
00142       mehSiStripResY[amend]->setAxisTitle("Y-resolution in "+
00143                                           SiStripString[amend],1);
00144       mehSiStripResY[amend]->setAxisTitle("Count",2);
00145     }
00146     
00147     
00148     //HCal
00149     //string hcharname, hchartitle;
00150     string HCalString[4]={"HB", "HE", "HF", "HO"};
00151     float HCalnUpper[4]={3000.,3000.,3000.,3000.}; 
00152     float HCalnLower[4]={0.,0.,0.,0.};
00153     for(int j =0; j <4; ++j) {
00154       mehHcaln[j]=0;
00155       mehHcalRes[j]=0;
00156     }
00157     
00158     dbe->setCurrentFolder("GlobalRecHitsV/HCals");
00159     for(int amend = 0; amend < 4; ++amend) {
00160       hcharname = "hHcaln_"+HCalString[amend];
00161       hchartitle= HCalString[amend]+"  rechits";
00162       mehHcaln[amend] = dbe->book1D(hcharname,hchartitle, 1000, HCalnLower[amend], 
00163                                     HCalnUpper[amend]);
00164       mehHcaln[amend]->setAxisTitle("Number of RecHits",1);
00165       mehHcaln[amend]->setAxisTitle("Count",2);
00166       hcharname = "hHcalRes_"+HCalString[amend];
00167       hchartitle= HCalString[amend]+"  rechit resolution";
00168       mehHcalRes[amend] = dbe->book1D(hcharname,hchartitle, 25, -2., 2.);
00169       mehHcalRes[amend]->setAxisTitle("RecHit E - SimHit E",1);
00170       mehHcalRes[amend]->setAxisTitle("Count",2);
00171     }
00172     
00173     
00174     //Ecal
00175     string ECalString[3] = {"EB","EE", "ES"}; 
00176     int ECalnBins[3] = {1000,3000,150};
00177     float ECalnUpper[3] = {20000., 62000., 3000.};
00178     float ECalnLower[3] = {0., 0., 0.};
00179     int ECalResBins[3] = {200,200,200};
00180     float ECalResUpper[3] = {1., 0.3, .0002};
00181     float ECalResLower[3] = {-1., -0.3, -.0002};
00182     for(int i =0; i<3; ++i) {
00183       mehEcaln[i]=0;
00184       mehEcalRes[i]=0;
00185     }
00186     dbe->setCurrentFolder("GlobalRecHitsV/ECals");
00187     
00188     for(int amend = 0; amend < 3; ++amend) {
00189       hcharname = "hEcaln_"+ECalString[amend];
00190       hchartitle= ECalString[amend]+"  rechits";
00191       mehEcaln[amend] = dbe->book1D(hcharname,hchartitle, ECalnBins[amend], 
00192                                     ECalnLower[amend], ECalnUpper[amend]);
00193       mehEcaln[amend]->setAxisTitle("Number of RecHits",1);
00194       mehEcaln[amend]->setAxisTitle("Count",2);
00195       hcharname = "hEcalRes_"+ECalString[amend];
00196       hchartitle= ECalString[amend]+"  rechit resolution";
00197       mehEcalRes[amend] = dbe->book1D(hcharname,hchartitle,ECalResBins[amend], 
00198                                       ECalResLower[amend], 
00199                                       ECalResUpper[amend]);
00200       mehEcalRes[amend]->setAxisTitle("RecHit E - SimHit E",1);
00201       mehEcalRes[amend]->setAxisTitle("Count",2);
00202     }
00203     
00204     //Si Pixels
00205     string SiPixelString[7] = {"BRL1", "BRL2", "BRL3", "FWD1n", "FWD1p", 
00206                                "FWD2n", "FWD2p"};
00207     for(int j =0; j<7; ++j) {
00208       mehSiPixeln[j]=0;
00209       mehSiPixelResX[j]=0;
00210       mehSiPixelResY[j]=0;
00211     }
00212     
00213     dbe->setCurrentFolder("GlobalRecHitsV/SiPixels");
00214     for(int amend = 0; amend < 7; ++amend) {
00215       hcharname = "hSiPixeln_"+SiPixelString[amend];
00216       hchartitle= SiPixelString[amend]+" rechits";
00217       mehSiPixeln[amend] = dbe->book1D(hcharname,hchartitle,200,0.,200.);
00218       mehSiPixeln[amend]->setAxisTitle("Number of hits in "+
00219                                        SiPixelString[amend],1);
00220       mehSiPixeln[amend]->setAxisTitle("Count",2);
00221       hcharname = "hSiPixelResX_"+SiPixelString[amend];
00222       hchartitle= SiPixelString[amend]+" rechit x resolution";
00223       mehSiPixelResX[amend] = dbe->book1D(hcharname,hchartitle,200,-0.02,.02);
00224       mehSiPixelResX[amend]->setAxisTitle("X-resolution in "+
00225                                           SiPixelString[amend],1);
00226       mehSiPixelResX[amend]->setAxisTitle("Count",2);
00227       hcharname = "hSiPixelResY_"+SiPixelString[amend];
00228       hchartitle= SiPixelString[amend]+" rechit y resolution";
00229       
00230       mehSiPixelResY[amend] = dbe->book1D(hcharname,hchartitle,200,-0.02,.02);
00231       mehSiPixelResY[amend]->setAxisTitle("Y-resolution in "+
00232                                           SiPixelString[amend],1);
00233       mehSiPixelResY[amend]->setAxisTitle("Count",2);
00234     }
00235  
00236     //Muons 
00237     dbe->setCurrentFolder("GlobalRecHitsV/Muons");
00238     
00239     mehDtMuonn = 0;
00240     mehCSCn = 0;
00241     mehRPCn = 0;
00242     
00243     string n_List[3] = {"hDtMuonn", "hCSCn", "hRPCn"};
00244     string hist_string[3] = {"Dt", "CSC", "RPC"};
00245     
00246     for(int amend=0; amend<3; ++amend) {
00247       hchartitle = hist_string[amend]+" rechits";
00248       if(amend==0) {
00249         mehDtMuonn=dbe->book1D(n_List[amend],hchartitle,50, 0., 500.);
00250         mehDtMuonn->setAxisTitle("Number of Rechits",1);
00251         mehDtMuonn->setAxisTitle("Count",2);
00252       }
00253       if(amend==1) {
00254         mehCSCn=dbe->book1D(n_List[amend],hchartitle,50, 0., 500.);
00255         mehCSCn->setAxisTitle("Number of Rechits",1);
00256         mehCSCn->setAxisTitle("Count",2);
00257       }
00258       if(amend==2){
00259         mehRPCn=dbe->book1D(n_List[amend],hchartitle,50, 0., 500.);
00260         mehRPCn->setAxisTitle("Number of Rechits",1);
00261         mehRPCn->setAxisTitle("Count",2);
00262       }
00263     }
00264     
00265     mehDtMuonRes=0;
00266     mehCSCResRDPhi=0;
00267     mehRPCResX=0;
00268     
00269     hcharname = "hDtMuonRes";
00270     hchartitle = "DT wire distance resolution";
00271     mehDtMuonRes = dbe->book1D(hcharname, hchartitle, 200, -0.2, 0.2);
00272     hcharname = "CSCResRDPhi";
00273     hchartitle = "CSC perp*dphi resolution";
00274     mehCSCResRDPhi = dbe->book1D(hcharname, hchartitle, 200, -0.2, 0.2);
00275     hcharname = "hRPCResX";
00276     hchartitle = "RPC rechits x resolution";
00277     mehRPCResX = dbe->book1D(hcharname, hchartitle, 50, -5., 5.);
00278   } 
00279 }
00280 
00281 GlobalRecHitsAnalyzer::~GlobalRecHitsAnalyzer() {}
00282 
00283 void GlobalRecHitsAnalyzer::beginJob()
00284 {
00285   return;
00286 }
00287 
00288 void GlobalRecHitsAnalyzer::endJob()
00289 {
00290   std::string MsgLoggerCat = "GlobalRecHitsAnalyzer_endJob";
00291   if (verbosity >= 0)
00292     edm::LogInfo(MsgLoggerCat) 
00293       << "Terminating having processed " << count << " events.";
00294   return;
00295 }
00296 
00297 void GlobalRecHitsAnalyzer::analyze(const edm::Event& iEvent, 
00298                                     const edm::EventSetup& iSetup)
00299 {
00300   std::string MsgLoggerCat = "GlobalRecHitsAnalyzer_analyze";
00301   
00302   // keep track of number of events processed
00303   ++count;
00304   
00305   // get event id information
00306   int nrun = iEvent.id().run();
00307   int nevt = iEvent.id().event();
00308   
00309   if (verbosity > 0) {
00310     edm::LogInfo(MsgLoggerCat)
00311       << "Processing run " << nrun << ", event " << nevt
00312       << " (" << count << " events total)";
00313   } else if (verbosity == 0) {
00314     if (nevt%frequency == 0 || nevt == 1) {
00315       edm::LogInfo(MsgLoggerCat)
00316         << "Processing run " << nrun << ", event " << nevt
00317         << " (" << count << " events total)";
00318     }
00319   }
00320   
00321   // look at information available in the event
00322   if (getAllProvenances) {
00323     
00324     std::vector<const edm::Provenance*> AllProv;
00325     iEvent.getAllProvenance(AllProv);
00326     
00327     if (verbosity >= 0)
00328       edm::LogInfo(MsgLoggerCat)
00329         << "Number of Provenances = " << AllProv.size();
00330     
00331     if (printProvenanceInfo && (verbosity >= 0)) {
00332       TString eventout("\nProvenance info:\n");      
00333       
00334       for (unsigned int i = 0; i < AllProv.size(); ++i) {
00335         eventout += "\n       ******************************";
00336         eventout += "\n       Module       : ";
00337         eventout += AllProv[i]->moduleLabel();
00338         eventout += "\n       ProductID    : ";
00339         eventout += AllProv[i]->productID().id();
00340         eventout += "\n       ClassName    : ";
00341         eventout += AllProv[i]->className();
00342         eventout += "\n       InstanceName : ";
00343         eventout += AllProv[i]->productInstanceName();
00344         eventout += "\n       BranchName   : ";
00345         eventout += AllProv[i]->branchName();
00346       }
00347       eventout += "\n       ******************************\n";
00348       edm::LogInfo(MsgLoggerCat) << eventout << "\n";
00349       printProvenanceInfo = false;
00350     }
00351     getAllProvenances = false;
00352   }
00353   
00354   // call fill functions
00355   // gather Ecal information from event
00356   fillECal(iEvent, iSetup);
00357   // gather Hcal information from event
00358   fillHCal(iEvent, iSetup);
00359   // gather Track information from event
00360   fillTrk(iEvent, iSetup);
00361   // gather Muon information from event
00362   fillMuon(iEvent, iSetup);
00363   
00364   if (verbosity > 0)
00365     edm::LogInfo (MsgLoggerCat)
00366       << "Done gathering data from event.";
00367   
00368   return;
00369 }
00370 
00371 void GlobalRecHitsAnalyzer::fillECal(const edm::Event& iEvent, 
00372                                      const edm::EventSetup& iSetup)
00373 {
00374   std::string MsgLoggerCat = "GlobalRecHitsAnalyzer_fillECal";
00375   
00376   TString eventout;
00377   if (verbosity > 0)
00378     eventout = "\nGathering info:";  
00379   
00380   // extract crossing frame from event
00381   edm::Handle<CrossingFrame<PCaloHit> > crossingFrame;
00382   
00384   //extract EB information
00386   edm::Handle<EBUncalibratedRecHitCollection> EcalUncalibRecHitEB;
00387   iEvent.getByLabel(ECalUncalEBSrc_, EcalUncalibRecHitEB);
00388   bool validUncalibRecHitEB = true;
00389   if (!EcalUncalibRecHitEB.isValid()) {
00390     LogDebug(MsgLoggerCat)
00391       << "Unable to find EcalUncalRecHitEB in event!";
00392     validUncalibRecHitEB = false;
00393   }  
00394   
00395   edm::Handle<EBRecHitCollection> EcalRecHitEB;
00396   iEvent.getByLabel(ECalEBSrc_, EcalRecHitEB);
00397   bool validRecHitEB = true;
00398   if (!EcalRecHitEB.isValid()) {
00399     LogDebug(MsgLoggerCat)
00400       << "Unable to find EcalRecHitEB in event!";
00401     validRecHitEB = false;
00402   }  
00403 
00404   // loop over simhits
00405   const std::string barrelHitsName(hitsProducer+"EcalHitsEB");
00406   iEvent.getByLabel("mix",barrelHitsName,crossingFrame);
00407   bool validXFrame = true;
00408   if (!crossingFrame.isValid()) {
00409     LogDebug(MsgLoggerCat)
00410       << "Unable to find cal barrel crossingFrame in event!";
00411     validXFrame = false;
00412   }
00413 
00414   MapType ebSimMap;
00415   if (validXFrame) {
00416     std::auto_ptr<MixCollection<PCaloHit> >
00417       barrelHits(new MixCollection<PCaloHit>(crossingFrame.product()));  
00418     
00419     // keep track of sum of simhit energy in each crystal
00420     for (MixCollection<PCaloHit>::MixItr hitItr 
00421            = barrelHits->begin();
00422          hitItr != barrelHits->end();
00423          ++hitItr) {
00424       
00425       EBDetId ebid = EBDetId(hitItr->id());
00426       
00427       uint32_t crystid = ebid.rawId();
00428       ebSimMap[crystid] += hitItr->energy();
00429     }
00430   }  
00431 
00432   int nEBRecHits = 0;
00433   // loop over RecHits
00434   if (validUncalibRecHitEB && validRecHitEB) {
00435     const EBUncalibratedRecHitCollection *EBUncalibRecHit = 
00436       EcalUncalibRecHitEB.product();
00437     const EBRecHitCollection *EBRecHit = EcalRecHitEB.product();
00438   
00439     for (EcalUncalibratedRecHitCollection::const_iterator uncalibRecHit =
00440            EBUncalibRecHit->begin();
00441          uncalibRecHit != EBUncalibRecHit->end();
00442          ++uncalibRecHit) {
00443       
00444       EBDetId EBid = EBDetId(uncalibRecHit->id());
00445       
00446       EcalRecHitCollection::const_iterator myRecHit = EBRecHit->find(EBid);
00447       
00448       if (myRecHit != EBRecHit->end()) {
00449         ++nEBRecHits;
00450         mehEcalRes[1]->Fill(myRecHit->energy()-ebSimMap[EBid.rawId()]);
00451       }
00452     }
00453     
00454     if (verbosity > 1) {
00455       eventout += "\n          Number of EBRecHits collected:............ ";
00456       eventout += nEBRecHits;
00457     }
00458     mehEcaln[1]->Fill((float)nEBRecHits); 
00459   }
00460 
00462   //extract EE information
00464   edm::Handle<EEUncalibratedRecHitCollection> EcalUncalibRecHitEE;
00465   iEvent.getByLabel(ECalUncalEESrc_, EcalUncalibRecHitEE);
00466   bool validuncalibRecHitEE = true;
00467   if (!EcalUncalibRecHitEE.isValid()) {
00468     LogDebug(MsgLoggerCat)
00469       << "Unable to find EcalUncalRecHitEE in event!";
00470     validuncalibRecHitEE = false;
00471   }  
00472   
00473   edm::Handle<EERecHitCollection> EcalRecHitEE;
00474   iEvent.getByLabel(ECalEESrc_, EcalRecHitEE);
00475   bool validRecHitEE = true;
00476   if (!EcalRecHitEE.isValid()) {
00477     LogDebug(MsgLoggerCat)
00478       << "Unable to find EcalRecHitEE in event!";
00479     validRecHitEE = false;
00480   }  
00481   
00482   // loop over simhits
00483   const std::string endcapHitsName(hitsProducer+"EcalHitsEE");
00484   iEvent.getByLabel("mix",endcapHitsName,crossingFrame);
00485   validXFrame = true;
00486   if (!crossingFrame.isValid()) {
00487     LogDebug(MsgLoggerCat)
00488       << "Unable to find cal endcap crossingFrame in event!";
00489     validXFrame = false;
00490   }
00491 
00492   MapType eeSimMap;
00493   if (validXFrame) {
00494     std::auto_ptr<MixCollection<PCaloHit> >
00495       endcapHits(new MixCollection<PCaloHit>(crossingFrame.product()));  
00496     
00497     // keep track of sum of simhit energy in each crystal
00498     for (MixCollection<PCaloHit>::MixItr hitItr 
00499            = endcapHits->begin();
00500          hitItr != endcapHits->end();
00501          ++hitItr) {
00502       
00503       EEDetId eeid = EEDetId(hitItr->id());
00504       
00505       uint32_t crystid = eeid.rawId();
00506       eeSimMap[crystid] += hitItr->energy();
00507     }
00508   }    
00509 
00510   int nEERecHits = 0;
00511   if (validuncalibRecHitEE && validRecHitEE) {
00512     // loop over RecHits
00513     const EEUncalibratedRecHitCollection *EEUncalibRecHit = 
00514       EcalUncalibRecHitEE.product();
00515     const EERecHitCollection *EERecHit = EcalRecHitEE.product();
00516     
00517     for (EcalUncalibratedRecHitCollection::const_iterator uncalibRecHit =
00518            EEUncalibRecHit->begin();
00519          uncalibRecHit != EEUncalibRecHit->end();
00520          ++uncalibRecHit) {
00521       
00522       EEDetId EEid = EEDetId(uncalibRecHit->id());
00523       
00524       EcalRecHitCollection::const_iterator myRecHit = EERecHit->find(EEid);
00525       
00526       if (myRecHit != EERecHit->end()) {
00527         ++nEERecHits;
00528         mehEcalRes[0]->Fill(myRecHit->energy()-eeSimMap[EEid.rawId()]);
00529       }
00530     }
00531     
00532     if (verbosity > 1) {
00533       eventout += "\n          Number of EERecHits collected:............ ";
00534       eventout += nEERecHits;
00535     }
00536     mehEcaln[0]->Fill((float)nEERecHits);
00537   }    
00538 
00540   //extract ES information
00542   edm::Handle<ESRecHitCollection> EcalRecHitES;
00543   iEvent.getByLabel(ECalESSrc_, EcalRecHitES);
00544   bool validRecHitES = true;
00545   if (!EcalRecHitES.isValid()) {
00546     LogDebug(MsgLoggerCat)
00547       << "Unable to find EcalRecHitES in event!";
00548     validRecHitES = false;
00549   }  
00550 
00551   // loop over simhits
00552   const std::string preshowerHitsName(hitsProducer+"EcalHitsES");
00553   iEvent.getByLabel("mix",preshowerHitsName,crossingFrame);
00554   validXFrame = true;
00555   if (!crossingFrame.isValid()) {
00556     LogDebug(MsgLoggerCat)
00557       << "Unable to find cal preshower crossingFrame in event!";
00558     validXFrame = false;
00559   }
00560 
00561   MapType esSimMap;
00562   if (validXFrame) {
00563     std::auto_ptr<MixCollection<PCaloHit> >
00564       preshowerHits(new MixCollection<PCaloHit>(crossingFrame.product()));  
00565     
00566     // keep track of sum of simhit energy in each crystal
00567     for (MixCollection<PCaloHit>::MixItr hitItr 
00568            = preshowerHits->begin();
00569          hitItr != preshowerHits->end();
00570          ++hitItr) {
00571       
00572       ESDetId esid = ESDetId(hitItr->id());
00573       
00574       uint32_t crystid = esid.rawId();
00575       esSimMap[crystid] += hitItr->energy();
00576     }
00577   }
00578 
00579   int nESRecHits = 0;
00580   if (validRecHitES) {
00581     // loop over RecHits
00582     const ESRecHitCollection *ESRecHit = EcalRecHitES.product();
00583     for (EcalRecHitCollection::const_iterator recHit =
00584            ESRecHit->begin();
00585          recHit != ESRecHit->end();
00586          ++recHit) {
00587       
00588       ESDetId ESid = ESDetId(recHit->id());
00589       
00590       ++nESRecHits;
00591       mehEcalRes[2]->Fill(recHit->energy()-esSimMap[ESid.rawId()]);
00592     }
00593     
00594     if (verbosity > 1) {
00595       eventout += "\n          Number of ESRecHits collected:............ ";
00596       eventout += nESRecHits;
00597     }
00598     mehEcaln[2]->Fill(float(nESRecHits));
00599   }
00600 
00601   if (verbosity > 0)
00602     edm::LogInfo(MsgLoggerCat) << eventout << "\n";
00603   
00604   return;
00605 }
00606 
00607 void GlobalRecHitsAnalyzer::fillHCal(const edm::Event& iEvent, 
00608                                      const edm::EventSetup& iSetup)
00609 {
00610   std::string MsgLoggerCat = "GlobalRecHitsAnalyzer_fillHCal";
00611   
00612   TString eventout;
00613   if (verbosity > 0)
00614     eventout = "\nGathering info:";  
00615   
00616   // get geometry
00617   edm::ESHandle<CaloGeometry> geometry;
00618   iSetup.get<CaloGeometryRecord>().get(geometry);
00619   if (!geometry.isValid()) {
00620     edm::LogWarning(MsgLoggerCat)
00621       << "Unable to find CaloGeometry in event!";
00622     return;
00623   }
00624   
00625   // iterator to access containers
00626   edm::PCaloHitContainer::const_iterator itHit;
00627   
00629   // extract simhit info
00631   edm::Handle<edm::PCaloHitContainer> hcalHits;
00632   iEvent.getByLabel(HCalSrc_,hcalHits);
00633   bool validhcalHits = true;
00634   if (!hcalHits.isValid()) {
00635     LogDebug(MsgLoggerCat)
00636       << "Unable to find hcalHits in event!";
00637     validhcalHits = false;
00638   }  
00639 
00640   MapType fHBEnergySimHits;
00641   MapType fHEEnergySimHits;
00642   MapType fHOEnergySimHits;
00643   MapType fHFEnergySimHits;
00644   if (validhcalHits) {
00645     const edm::PCaloHitContainer *simhitResult = hcalHits.product();
00646   
00647     for (std::vector<PCaloHit>::const_iterator simhits = simhitResult->begin();
00648          simhits != simhitResult->end();
00649          ++simhits) {
00650       
00651       HcalDetId detId(simhits->id());
00652       uint32_t cellid = detId.rawId();
00653       
00654       if (detId.subdet() == sdHcalBrl){  
00655         fHBEnergySimHits[cellid] += simhits->energy(); 
00656       }
00657       if (detId.subdet() == sdHcalEC){  
00658         fHEEnergySimHits[cellid] += simhits->energy(); 
00659       }    
00660       if (detId.subdet() == sdHcalOut){  
00661         fHOEnergySimHits[cellid] += simhits->energy(); 
00662       }    
00663       if (detId.subdet() == sdHcalFwd){  
00664         fHFEnergySimHits[cellid] += simhits->energy(); 
00665       }    
00666     }
00667   }
00668 
00669   // max values to be used (HO is found in HB)
00670   Double_t maxHBEnergy = 0.;
00671   Double_t maxHEEnergy = 0.;
00672   Double_t maxHOEnergy = 0.;
00673   Double_t maxHFEnergy = 0.;
00674   
00675   Double_t maxHBPhi = -1000.;
00676   Double_t maxHEPhi = -1000.;
00677   Double_t maxHOPhi = -1000.;
00678   Double_t maxHFPhi = -1000.;
00679   
00680   Double_t maxHBEta = -1000.;
00681   Double_t maxHEEta = -1000.;
00682   Double_t maxHOEta = -1000.;
00683   
00684   Double_t PI = 3.141592653589;
00685   
00687   // get HBHE information
00689   std::vector<edm::Handle<HBHERecHitCollection> > hbhe;
00690   iEvent.getManyByType(hbhe);
00691   bool validHBHE = true;
00692   if (!hbhe[0].isValid()) {
00693     LogDebug(MsgLoggerCat)
00694       << "Unable to find any HBHERecHitCollections in event!";
00695     validHBHE = false;
00696   } 
00697 
00698   if (validHBHE) {
00699     std::vector<edm::Handle<HBHERecHitCollection> >::iterator ihbhe;
00700     
00701     int iHB = 0;
00702     int iHE = 0; 
00703     for (ihbhe = hbhe.begin(); ihbhe != hbhe.end(); ++ihbhe) {
00704       
00705       // find max values
00706       for (HBHERecHitCollection::const_iterator jhbhe = (*ihbhe)->begin();
00707            jhbhe != (*ihbhe)->end(); ++jhbhe) {
00708         
00709         HcalDetId cell(jhbhe->id());
00710         
00711         if (cell.subdet() == sdHcalBrl) {
00712           
00713           const CaloCellGeometry* cellGeometry =
00714             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00715           double fEta = cellGeometry->getPosition().eta () ;
00716           double fPhi = cellGeometry->getPosition().phi () ;
00717           if ( (jhbhe->energy()) > maxHBEnergy ) {
00718             maxHBEnergy = jhbhe->energy();
00719             maxHOEnergy = maxHBEnergy;
00720             maxHBPhi = fPhi;
00721             maxHOPhi = maxHBPhi;
00722             maxHBEta = fEta;
00723             maxHOEta = maxHBEta;
00724           }       
00725         }
00726       
00727         if (cell.subdet() == sdHcalEC) {
00728           
00729           const CaloCellGeometry* cellGeometry =
00730             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00731           double fEta = cellGeometry->getPosition().eta () ;
00732           double fPhi = cellGeometry->getPosition().phi () ;
00733           if ( (jhbhe->energy()) > maxHEEnergy ) {
00734             maxHEEnergy = jhbhe->energy();
00735             maxHEPhi = fPhi;
00736             maxHEEta = fEta;
00737           }       
00738         }
00739       } // end find max values
00740       
00741       for (HBHERecHitCollection::const_iterator jhbhe = (*ihbhe)->begin();
00742            jhbhe != (*ihbhe)->end(); ++jhbhe) {
00743         
00744         HcalDetId cell(jhbhe->id());
00745         
00746         if (cell.subdet() == sdHcalBrl) {
00747           
00748           ++iHB;
00749           
00750           const CaloCellGeometry* cellGeometry =
00751             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00752           double fPhi = cellGeometry->getPosition().phi () ;
00753           
00754           float deltaphi = maxHBPhi - fPhi;
00755           if (fPhi > maxHBPhi) { deltaphi = fPhi - maxHBPhi;}
00756           if (deltaphi > PI) { deltaphi = 2.0 * PI - deltaphi;}
00757           
00758           mehHcalRes[0]->Fill(jhbhe->energy() - 
00759                               fHBEnergySimHits[cell.rawId()]);
00760         }
00761         
00762         if (cell.subdet() == sdHcalEC) {
00763           
00764           ++iHE;
00765           
00766           const CaloCellGeometry* cellGeometry =
00767             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00768           double fPhi = cellGeometry->getPosition().phi () ;
00769           
00770           float deltaphi = maxHEPhi - fPhi;
00771           if (fPhi > maxHEPhi) { deltaphi = fPhi - maxHEPhi;}
00772           if (deltaphi > PI) { deltaphi = 2.0 * PI - deltaphi;}
00773           mehHcalRes[1]->Fill(jhbhe->energy() - 
00774                               fHEEnergySimHits[cell.rawId()]);
00775         }
00776       }
00777     } // end loop through collection
00778   
00779     
00780     if (verbosity > 1) {
00781       eventout += "\n          Number of HBRecHits collected:............ ";
00782       eventout += iHB;
00783     }
00784     
00785     if (verbosity > 1) {
00786       eventout += "\n          Number of HERecHits collected:............ ";
00787       eventout += iHE;
00788     }
00789     mehHcaln[0]->Fill((float)iHB);
00790     mehHcaln[1]->Fill((float)iHE);
00791   }
00792 
00794   // get HF information
00796   std::vector<edm::Handle<HFRecHitCollection> > hf;
00797   iEvent.getManyByType(hf);
00798   bool validHF = true;
00799   if (!hf[0].isValid()) {
00800     LogDebug(MsgLoggerCat)
00801       << "Unable to find any HFRecHitCollections in event!";
00802     validHF = false;
00803   } 
00804   if (validHF) {
00805     std::vector<edm::Handle<HFRecHitCollection> >::iterator ihf;
00806     
00807     int iHF = 0; 
00808     for (ihf = hf.begin(); ihf != hf.end(); ++ihf) {
00809       
00810       // find max values
00811       for (HFRecHitCollection::const_iterator jhf = (*ihf)->begin();
00812            jhf != (*ihf)->end(); ++jhf) {
00813         
00814         HcalDetId cell(jhf->id());
00815         
00816         if (cell.subdet() == sdHcalFwd) {
00817           
00818           const CaloCellGeometry* cellGeometry =
00819             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00820           double fPhi = cellGeometry->getPosition().phi () ;
00821           if ( (jhf->energy()) > maxHFEnergy ) {
00822             maxHFEnergy = jhf->energy();
00823             maxHFPhi = fPhi;
00824           }       
00825         }
00826       } // end find max values
00827       
00828       for (HFRecHitCollection::const_iterator jhf = (*ihf)->begin();
00829            jhf != (*ihf)->end(); ++jhf) {
00830         
00831         HcalDetId cell(jhf->id());
00832         
00833         if (cell.subdet() == sdHcalFwd) {
00834           
00835           ++iHF;
00836           
00837           const CaloCellGeometry* cellGeometry =
00838             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00839           double fPhi = cellGeometry->getPosition().phi () ;
00840           
00841           float deltaphi = maxHBPhi - fPhi;
00842           if (fPhi > maxHFPhi) { deltaphi = fPhi - maxHFPhi;}
00843           if (deltaphi > PI) { deltaphi = 2.0 * PI - deltaphi;}
00844           
00845           mehHcalRes[2]->Fill(jhf->energy()-fHFEnergySimHits[cell.rawId()]);
00846         }
00847       }
00848     } // end loop through collection
00849     
00850     if (verbosity > 1) {
00851       eventout += "\n          Number of HFDigis collected:.............. ";
00852       eventout += iHF;
00853     }
00854     mehHcaln[2]->Fill((float)iHF);
00855   }    
00856 
00858   // get HO information
00860   std::vector<edm::Handle<HORecHitCollection> > ho;
00861   iEvent.getManyByType(ho);
00862   bool validHO = true;
00863   if (!ho[0].isValid()) {
00864     LogDebug(MsgLoggerCat)
00865       << "Unable to find any HORecHitCollections in event!";
00866     validHO = false;
00867   } 
00868 
00869   if (validHO) {
00870     std::vector<edm::Handle<HORecHitCollection> >::iterator iho;
00871     
00872     int iHO = 0; 
00873     for (iho = ho.begin(); iho != ho.end(); ++iho) {
00874       
00875       for (HORecHitCollection::const_iterator jho = (*iho)->begin();
00876            jho != (*iho)->end(); ++jho) {
00877         
00878         HcalDetId cell(jho->id());
00879         
00880         if (cell.subdet() == sdHcalOut) {
00881           
00882           ++iHO;
00883           
00884           const CaloCellGeometry* cellGeometry =
00885             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00886           double fPhi = cellGeometry->getPosition().phi () ;
00887           
00888           float deltaphi = maxHOPhi - fPhi;
00889           if (fPhi > maxHOPhi) { deltaphi = fPhi - maxHOPhi;}
00890           if (deltaphi > PI) { deltaphi = 2.0 * PI - deltaphi;}
00891           mehHcalRes[3]->Fill(jho->energy()-fHOEnergySimHits[cell.rawId()]);
00892         }
00893       }
00894     } // end loop through collection
00895     
00896     if (verbosity > 1) {
00897       eventout += "\n          Number of HODigis collected:.............. ";
00898       eventout += iHO;
00899     }
00900     mehHcaln[3]->Fill((float)iHO);
00901   }
00902 
00903   if (verbosity > 0)
00904     edm::LogInfo(MsgLoggerCat) << eventout << "\n";
00905   
00906   return;
00907 }
00908 
00909 void GlobalRecHitsAnalyzer::fillTrk(const edm::Event& iEvent, 
00910                                     const edm::EventSetup& iSetup)
00911 {
00912   std::string MsgLoggerCat = "GlobalRecHitsAnalyzer_fillTrk";
00913   
00914   TString eventout;
00915   if (verbosity > 0)
00916     eventout = "\nGathering info:";  
00917   
00918   // get strip information
00919   edm::Handle<SiStripMatchedRecHit2DCollection> rechitsmatched;
00920   iEvent.getByLabel(SiStripSrc_, rechitsmatched);
00921   bool validstrip = true;
00922   if (!rechitsmatched.isValid()) {
00923     LogDebug(MsgLoggerCat)
00924       << "Unable to find stripmatchedrechits in event!";
00925     validstrip = false;
00926   }  
00927   
00928   TrackerHitAssociator associate(iEvent,conf_);
00929   
00930   edm::ESHandle<TrackerGeometry> pDD;
00931   iSetup.get<TrackerDigiGeometryRecord>().get(pDD);
00932   if (!pDD.isValid()) {
00933     edm::LogWarning(MsgLoggerCat)
00934       << "Unable to find TrackerDigiGeometry in event!";
00935     return;
00936   }
00937   const TrackerGeometry &tracker(*pDD);
00938   
00939   if (validstrip) {
00940     int nStripBrl = 0, nStripFwd = 0;
00941     
00942     // loop over det units
00943     for (TrackerGeometry::DetContainer::const_iterator it = 
00944            pDD->dets().begin();
00945          it != pDD->dets().end(); ++it) {
00946       
00947       uint32_t myid = ((*it)->geographicalId()).rawId();
00948       DetId detid = ((*it)->geographicalId());
00949       
00950       //loop over rechits-matched in the same subdetector
00951       SiStripMatchedRecHit2DCollection::const_iterator rechitmatchedMatch = rechitsmatched->find(detid);
00952       
00953       if (rechitmatchedMatch != rechitsmatched->end()) {
00954         SiStripMatchedRecHit2DCollection::DetSet rechitmatchedRange = *rechitmatchedMatch;
00955         SiStripMatchedRecHit2DCollection::DetSet::const_iterator rechitmatchedRangeIteratorBegin = rechitmatchedRange.begin();
00956         SiStripMatchedRecHit2DCollection::DetSet::const_iterator rechitmatchedRangeIteratorEnd   = rechitmatchedRange.end();
00957         SiStripMatchedRecHit2DCollection::DetSet::const_iterator itermatched = rechitmatchedRangeIteratorBegin;
00958         
00959         for ( itermatched = rechitmatchedRangeIteratorBegin; 
00960               itermatched != rechitmatchedRangeIteratorEnd;
00961               ++itermatched) {
00962           
00963           SiStripMatchedRecHit2D const rechit = *itermatched;
00964           LocalPoint position = rechit.localPosition();
00965           
00966           float mindist = 999999.;
00967           float distx = 999999.;
00968           float disty = 999999.;
00969           float dist = 999999.;
00970           std::pair<LocalPoint,LocalVector> closestPair;
00971           matched.clear();
00972           
00973           float rechitmatchedx = position.x();
00974           float rechitmatchedy = position.y();
00975           
00976           matched = associate.associateHit(rechit);
00977           
00978           if (!matched.empty()) {
00979             //project simhit;
00980             const GluedGeomDet* gluedDet = 
00981               (const GluedGeomDet*)tracker.idToDet(rechit.geographicalId());
00982             const StripGeomDetUnit* partnerstripdet =
00983               (StripGeomDetUnit*) gluedDet->stereoDet();
00984             std::pair<LocalPoint,LocalVector> hitPair;
00985             
00986             for(std::vector<PSimHit>::const_iterator m = matched.begin(); 
00987                 m != matched.end(); m++){
00988               //project simhit;
00989               hitPair = projectHit((*m),partnerstripdet,gluedDet->surface());
00990               distx = fabs(rechitmatchedx - hitPair.first.x());
00991               disty = fabs(rechitmatchedy - hitPair.first.y());
00992               dist = sqrt(distx*distx+disty*disty);
00993               
00994               if(dist < mindist){
00995                 mindist = dist;
00996                 closestPair = hitPair;
00997               }
00998             }
00999             
01000             // get TIB
01001             if (detid.subdetId() == sdSiTIB) {
01002               
01003               TIBDetId tibid(myid);
01004               ++nStripBrl;
01005               
01006               if (tibid.layer() == 1) {
01007                 mehSiStripResX[8]->Fill(rechitmatchedx-closestPair.first.x());
01008                 mehSiStripResY[8]->Fill(rechitmatchedy-closestPair.first.y());
01009               }
01010               if (tibid.layer() == 2) {
01011                 mehSiStripResX[9]->Fill(rechitmatchedx-closestPair.first.x());
01012                 mehSiStripResY[9]->Fill(rechitmatchedy-closestPair.first.y());
01013               } 
01014               if (tibid.layer() == 3) {
01015                 mehSiStripResX[10]->Fill(rechitmatchedx-closestPair.first.x());
01016                 mehSiStripResY[10]->Fill(rechitmatchedy-closestPair.first.y());
01017                 
01018               }
01019               if (tibid.layer() == 4) {
01020                 mehSiStripResX[11]->Fill(rechitmatchedx-closestPair.first.x());
01021                 mehSiStripResY[11]->Fill(rechitmatchedy-closestPair.first.y());
01022               }
01023             }
01024             
01025             // get TOB
01026             if (detid.subdetId() == sdSiTOB) {
01027               
01028               TOBDetId tobid(myid);
01029               ++nStripBrl;
01030               
01031               if (tobid.layer() == 1) {
01032                 mehSiStripResX[15]->Fill(rechitmatchedx-closestPair.first.x());
01033                 mehSiStripResY[15]->Fill(rechitmatchedy-closestPair.first.y());
01034               }
01035               if (tobid.layer() == 2) {
01036                 mehSiStripResX[16]->Fill(rechitmatchedx-closestPair.first.x());
01037                 mehSiStripResY[16]->Fill(rechitmatchedy-closestPair.first.y());
01038               } 
01039               if (tobid.layer() == 3) {
01040                 mehSiStripResX[17]->Fill(rechitmatchedx-closestPair.first.x());
01041                 mehSiStripResY[17]->Fill(rechitmatchedy-closestPair.first.y());
01042               }
01043               if (tobid.layer() == 4) {
01044                 mehSiStripResX[18]->Fill(rechitmatchedx-closestPair.first.x());
01045                 mehSiStripResY[18]->Fill(rechitmatchedy-closestPair.first.y());
01046               }
01047             }
01048             
01049             // get TID
01050             if (detid.subdetId() == sdSiTID) {
01051               
01052               TIDDetId tidid(myid);
01053               ++nStripFwd;
01054               
01055               if (tidid.wheel() == 1) {
01056                 mehSiStripResX[12]->Fill(rechitmatchedx-closestPair.first.x());
01057                 mehSiStripResY[12]->Fill(rechitmatchedy-closestPair.first.y());
01058               }
01059               if (tidid.wheel() == 2) {
01060                 mehSiStripResX[13]->Fill(rechitmatchedx-closestPair.first.x());
01061                 mehSiStripResY[13]->Fill(rechitmatchedy-closestPair.first.y());
01062               } 
01063               if (tidid.wheel() == 3) {
01064                 mehSiStripResX[14]->Fill(rechitmatchedx-closestPair.first.x());
01065                 mehSiStripResY[14]->Fill(rechitmatchedy-closestPair.first.y());
01066               }
01067             }
01068             
01069             // get TEC
01070             if (detid.subdetId() == sdSiTEC) {
01071               
01072               TECDetId tecid(myid);
01073               ++nStripFwd;
01074               
01075               if (tecid.wheel() == 1) {
01076                 mehSiStripResX[0]->Fill(rechitmatchedx-closestPair.first.x());
01077                 mehSiStripResY[0]->Fill(rechitmatchedy-closestPair.first.y());
01078               }
01079               if (tecid.wheel() == 2) {
01080                 mehSiStripResX[1]->Fill(rechitmatchedx-closestPair.first.x());
01081                 mehSiStripResY[1]->Fill(rechitmatchedy-closestPair.first.y());
01082               } 
01083               if (tecid.wheel() == 3) {
01084                 mehSiStripResX[2]->Fill(rechitmatchedx-closestPair.first.x());
01085                 mehSiStripResY[2]->Fill(rechitmatchedy-closestPair.first.y());
01086               }
01087               if (tecid.wheel() == 4) {
01088                 mehSiStripResX[3]->Fill(rechitmatchedx-closestPair.first.x());
01089                 mehSiStripResY[3]->Fill(rechitmatchedy-closestPair.first.y());
01090                 
01091               }
01092               if (tecid.wheel() == 5) {
01093                 mehSiStripResX[4]->Fill(rechitmatchedx-closestPair.first.x());
01094                 mehSiStripResY[4]->Fill(rechitmatchedy-closestPair.first.y());
01095               } 
01096               if (tecid.wheel() == 6) {
01097                 mehSiStripResX[5]->Fill(rechitmatchedx-closestPair.first.x());
01098                 mehSiStripResY[5]->Fill(rechitmatchedy-closestPair.first.y());
01099               }
01100               if (tecid.wheel() == 7) {
01101                 mehSiStripResX[6]->Fill(rechitmatchedx-closestPair.first.x());
01102                 mehSiStripResY[6]->Fill(rechitmatchedy-closestPair.first.y());
01103               } 
01104               if (tecid.wheel() == 8) {
01105                 mehSiStripResX[7]->Fill(rechitmatchedx-closestPair.first.x());
01106                 mehSiStripResY[7]->Fill(rechitmatchedy-closestPair.first.y()); 
01107               }
01108             }
01109             
01110           } // end if matched empty
01111         }
01112       }
01113     } // end loop over det units
01114                                                                       
01115     if (verbosity > 1) {
01116       eventout += "\n          Number of BrlStripRecHits collected:...... ";
01117       eventout += nStripBrl;
01118     }
01119     
01120     for(int i =8; i<12; ++i)
01121       {mehSiStripn[i]->Fill((float)nStripBrl);}
01122     for(int i =16; i<19; ++i)
01123       {mehSiStripn[i]->Fill((float)nStripBrl);}
01124     
01125     if (verbosity > 1) {
01126       eventout += "\n          Number of FrwdStripRecHits collected:..... ";
01127       eventout += nStripFwd;
01128     }
01129     for(int i =0; i<8; ++i)
01130       {mehSiStripn[i]->Fill((float)nStripFwd);}
01131     for(int i =12; i<16; ++i)
01132       {mehSiStripn[i]->Fill((float)nStripFwd);}
01133   }
01134 
01135   // get pixel information
01136   //Get RecHits
01137   edm::Handle<SiPixelRecHitCollection> recHitColl;
01138   iEvent.getByLabel(SiPxlSrc_, recHitColl);
01139   bool validpixel = true;
01140   if (!recHitColl.isValid()) {
01141     LogDebug(MsgLoggerCat)
01142       << "Unable to find SiPixelRecHitCollection in event!";
01143     validpixel = false;
01144   }  
01145   
01146   //Get event setup
01147   edm::ESHandle<TrackerGeometry> geom;
01148   iSetup.get<TrackerDigiGeometryRecord>().get(geom); 
01149   if (!geom.isValid()) {
01150     edm::LogWarning(MsgLoggerCat)
01151       << "Unable to find TrackerDigiGeometry in event!";
01152     return;
01153   }
01154 
01155   if (validpixel) {
01156     int nPxlBrl = 0, nPxlFwd = 0;    
01157     //iterate over detunits
01158     for (TrackerGeometry::DetContainer::const_iterator it = 
01159            geom->dets().begin();
01160          it != geom->dets().end(); ++it) {
01161       
01162       uint32_t myid = ((*it)->geographicalId()).rawId();
01163       DetId detId = ((*it)->geographicalId());
01164       int subid = detId.subdetId();
01165       
01166       if (! ((subid == sdPxlBrl) || (subid == sdPxlFwd))) continue;
01167       
01168       SiPixelRecHitCollection::const_iterator pixeldet = recHitColl->find(detId);
01169       if (pixeldet == recHitColl->end()) continue;
01170       SiPixelRecHitCollection::DetSet pixelrechitRange = *pixeldet;
01171       SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorBegin = pixelrechitRange.begin();
01172       SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorEnd   = pixelrechitRange.end();
01173       SiPixelRecHitCollection::DetSet::const_iterator pixeliter = pixelrechitRangeIteratorBegin;
01174 
01175       
01176       std::vector<PSimHit> matched;
01177       
01178       //----Loop over rechits for this detId
01179       for ( ; pixeliter != pixelrechitRangeIteratorEnd; ++pixeliter) {
01180         
01181         matched.clear();
01182         matched = associate.associateHit(*pixeliter);
01183         
01184         if ( !matched.empty() ) {
01185           
01186           float closest = 9999.9;
01187           LocalPoint lp = pixeliter->localPosition();
01188           float rechit_x = lp.x();
01189           float rechit_y = lp.y();
01190           
01191           float sim_x = 0.;
01192           float sim_y = 0.;
01193           
01194           //loop over sim hits and fill closet
01195           for (std::vector<PSimHit>::const_iterator m = matched.begin(); 
01196                m != matched.end(); ++m) {
01197             
01198             float sim_x1 = (*m).entryPoint().x();
01199             float sim_x2 = (*m).exitPoint().x();
01200             float sim_xpos = 0.5*(sim_x1+sim_x2);
01201             
01202             float sim_y1 = (*m).entryPoint().y();
01203             float sim_y2 = (*m).exitPoint().y();
01204             float sim_ypos = 0.5*(sim_y1+sim_y2);
01205             
01206             float x_res = fabs(sim_xpos - rechit_x);
01207             float y_res = fabs(sim_ypos - rechit_y);
01208             
01209             float dist = sqrt(x_res*x_res + y_res*y_res);
01210             
01211             if ( dist < closest ) {
01212               closest = dist;
01213               sim_x = sim_xpos;
01214               sim_y = sim_ypos;
01215             }
01216           } // end sim hit loop
01217           
01218           // get Barrel pixels ***************Pixel STuff******************
01219           if (subid == sdPxlBrl) {
01220             PXBDetId bdetid(myid);
01221             ++nPxlBrl;
01222             
01223             if (bdetid.layer() == 1) {
01224               mehSiPixelResX[0]->Fill(rechit_x-sim_x);
01225               mehSiPixelResY[0]->Fill(rechit_y-sim_y); 
01226               
01227             }
01228             if (bdetid.layer() == 2) {
01229               mehSiPixelResX[1]->Fill(rechit_x-sim_x);
01230               mehSiPixelResY[1]->Fill(rechit_y-sim_y); 
01231             }
01232             if (bdetid.layer() == 3) {
01233               mehSiPixelResX[2]->Fill(rechit_x-sim_x);
01234               mehSiPixelResY[2]->Fill(rechit_y-sim_y); 
01235             }
01236           }
01237           
01238           // get Forward pixels
01239           if (subid == sdPxlFwd) {
01240             PXFDetId fdetid(myid);
01241             ++nPxlFwd;
01242             
01243             if (fdetid.disk() == 1) {
01244               if (fdetid.side() == 1) {
01245                 mehSiPixelResX[3]->Fill(rechit_x-sim_x);
01246                 mehSiPixelResY[3]->Fill(rechit_y-sim_y);
01247               }
01248               if (fdetid.side() == 2) {
01249                 mehSiPixelResX[4]->Fill(rechit_x-sim_x);
01250                 mehSiPixelResY[4]->Fill(rechit_y-sim_y); 
01251               }
01252             }
01253             if (fdetid.disk() == 2) {
01254               if (fdetid.side() == 1) {
01255                 mehSiPixelResX[5]->Fill(rechit_x-sim_x);
01256                 mehSiPixelResY[5]->Fill(rechit_y-sim_y);
01257               }
01258               if (fdetid.side() == 2) {
01259                 mehSiPixelResX[6]->Fill(rechit_x-sim_x);
01260                 mehSiPixelResY[6]->Fill(rechit_y-sim_y); 
01261               }
01262             }
01263           }      
01264         } // end matched emtpy
01265       } // <-----end rechit loop 
01266     } // <------ end detunit loop  
01267     
01268     
01269     if (verbosity > 1) {
01270       eventout += "\n          Number of BrlPixelRecHits collected:...... ";
01271       eventout += nPxlBrl;
01272     }
01273     for(int i=0; i<3; ++i) {
01274       mehSiPixeln[i]->Fill((float)nPxlBrl);
01275     }
01276     
01277     if (verbosity > 1) {
01278       eventout += "\n          Number of FrwdPixelRecHits collected:..... ";
01279       eventout += nPxlFwd;
01280     }
01281     
01282     for(int i=3; i<7; ++i) {
01283       mehSiPixeln[i]->Fill((float)nPxlFwd);
01284     }
01285   }
01286    
01287   if (verbosity > 0)
01288     edm::LogInfo(MsgLoggerCat) << eventout << "\n";
01289   
01290   return;
01291 }
01292 
01293 void GlobalRecHitsAnalyzer::fillMuon(const edm::Event& iEvent, 
01294                                      const edm::EventSetup& iSetup)
01295 {
01296   std::string MsgLoggerCat = "GlobalRecHitsAnalyzer_fillMuon";
01297   
01298   TString eventout;
01299   if (verbosity > 0)
01300     eventout = "\nGathering info:";  
01301 
01302   // get DT information
01303   edm::ESHandle<DTGeometry> dtGeom;
01304   iSetup.get<MuonGeometryRecord>().get(dtGeom);
01305   if (!dtGeom.isValid()) {
01306     edm::LogWarning(MsgLoggerCat)
01307       << "Unable to find DTMuonGeometryRecord in event!";
01308     return;
01309   }  
01310 
01311   edm::Handle<edm::PSimHitContainer> dtsimHits;
01312   iEvent.getByLabel(MuDTSimSrc_, dtsimHits);
01313   bool validdtsim = true;
01314   if (!dtsimHits.isValid()) {
01315     LogDebug(MsgLoggerCat)
01316       << "Unable to find dtsimHits in event!";
01317     validdtsim = false;
01318   } 
01319 
01320   edm::Handle<DTRecHitCollection> dtRecHits;
01321   iEvent.getByLabel(MuDTSrc_, dtRecHits);
01322   bool validdtrec = true;
01323   if (!dtRecHits.isValid()) {
01324     LogDebug(MsgLoggerCat)
01325       << "Unable to find dtRecHits in event!";
01326     validdtrec = false;
01327   }   
01328     
01329   if (validdtsim && validdtrec) {
01330 
01331     std::map<DTWireId, edm::PSimHitContainer> simHitsPerWire =
01332       DTHitQualityUtils::mapSimHitsPerWire(*(dtsimHits.product()));
01333 
01334     std::map<DTWireId, std::vector<DTRecHit1DPair> > recHitsPerWire =
01335       map1DRecHitsPerWire(dtRecHits.product());
01336 
01337     int nDt = compute(dtGeom.product(), simHitsPerWire, recHitsPerWire, 1);
01338     
01339     if (verbosity > 1) {
01340       eventout += "\n          Number of DtMuonRecHits collected:........ ";
01341       eventout += nDt;
01342     }
01343     mehDtMuonn->Fill(float(nDt));
01344   }
01345 
01346   // get CSC Strip information
01347   // get map of sim hits
01348   theMap.clear();
01349   edm::Handle<CrossingFrame<PSimHit> > cf;
01350 
01351   iEvent.getByLabel("mix",hitsProducer+"MuonCSCHits",cf);
01352   bool validXFrame = true;
01353   if (!cf.isValid()) {
01354     LogDebug(MsgLoggerCat)
01355       << "Unable to find muo CSC crossingFrame in event!";
01356     validXFrame = false;
01357   }
01358   if (validXFrame) {
01359     MixCollection<PSimHit> simHits(cf.product());
01360     
01361     // arrange the hits by detUnit
01362     for(MixCollection<PSimHit>::MixItr hitItr = simHits.begin();
01363         hitItr != simHits.end(); ++hitItr) {
01364       theMap[hitItr->detUnitId()].push_back(*hitItr);
01365     }  
01366   }
01367 
01368   // get geometry
01369   edm::ESHandle<CSCGeometry> hGeom;
01370   iSetup.get<MuonGeometryRecord>().get(hGeom);
01371   if (!hGeom.isValid()) {
01372     edm::LogWarning(MsgLoggerCat)
01373       << "Unable to find CSCMuonGeometryRecord in event!";
01374     return;
01375   }    
01376   const CSCGeometry *theCSCGeometry = &*hGeom;
01377 
01378   // get rechits
01379   edm::Handle<CSCRecHit2DCollection> hRecHits;
01380   iEvent.getByLabel(MuCSCSrc_, hRecHits);
01381   bool validCSC = true;
01382   if (!hRecHits.isValid()) {
01383     LogDebug(MsgLoggerCat)
01384       << "Unable to find CSC RecHits in event!";
01385     validCSC = false;
01386   }    
01387 
01388   if (validCSC) {
01389     const CSCRecHit2DCollection *cscRecHits = hRecHits.product();
01390     
01391     int nCSC = 0;
01392     for (CSCRecHit2DCollection::const_iterator recHitItr = cscRecHits->begin();
01393          recHitItr != cscRecHits->end(); ++recHitItr) {
01394       
01395       int detId = (*recHitItr).cscDetId().rawId();
01396       
01397       edm::PSimHitContainer simHits;   
01398       std::map<int, edm::PSimHitContainer>::const_iterator mapItr = 
01399         theMap.find(detId);
01400       if (mapItr != theMap.end()) {
01401         simHits = mapItr->second;
01402       }
01403       
01404       if (simHits.size() == 1) {
01405         ++nCSC;
01406         
01407         const GeomDetUnit* detUnit = 
01408           theCSCGeometry->idToDetUnit(CSCDetId(detId));
01409         const CSCLayer *layer = dynamic_cast<const CSCLayer *>(detUnit); 
01410         
01411         int chamberType = layer->chamber()->specs()->chamberType();
01412         plotResolution(simHits[0], *recHitItr, layer, chamberType);
01413       }
01414     }
01415     
01416     if (verbosity > 1) {
01417       eventout += "\n          Number of CSCRecHits collected:........... ";
01418       eventout += nCSC;
01419     }
01420     mehCSCn->Fill((float)nCSC);
01421   }
01422 
01423   // get RPC information
01424   std::map<double, int> mapsim, maprec;
01425   std::map<int, double> nmapsim, nmaprec;
01426 
01427   edm::ESHandle<RPCGeometry> rpcGeom;
01428   iSetup.get<MuonGeometryRecord>().get(rpcGeom);
01429   if (!rpcGeom.isValid()) {
01430     edm::LogWarning(MsgLoggerCat)
01431       << "Unable to find RPCMuonGeometryRecord in event!";
01432     return;
01433   }  
01434 
01435   edm::Handle<edm::PSimHitContainer> simHit;
01436   iEvent.getByLabel(MuRPCSimSrc_, simHit);
01437   bool validrpcsim = true;
01438   if (!simHit.isValid()) {
01439     LogDebug(MsgLoggerCat)
01440       << "Unable to find RPCSimHit in event!";
01441     validrpcsim = false;
01442   }    
01443 
01444   edm::Handle<RPCRecHitCollection> recHit;
01445   iEvent.getByLabel(MuRPCSrc_, recHit);
01446   bool validrpc = true;
01447   if (!simHit.isValid()) {
01448     LogDebug(MsgLoggerCat)
01449       << "Unable to find RPCRecHit in event!";
01450     validrpc = false;
01451   } 
01452 
01453   if (validrpc) {
01454     int nRPC = 0;
01455     RPCRecHitCollection::const_iterator recIt;
01456     int nrec = 0;
01457     for (recIt = recHit->begin(); recIt != recHit->end(); ++recIt) {
01458       RPCDetId Rid = (RPCDetId)(*recIt).rpcId();
01459       const RPCRoll *roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(Rid));
01460       if (roll->isForward()) {
01461         
01462         if (verbosity > 1) {
01463           eventout += 
01464             "\n          Number of RPCRecHits collected:........... ";
01465           eventout += nRPC;
01466         }
01467         
01468         if (verbosity > 0)
01469           edm::LogInfo(MsgLoggerCat) << eventout << "\n";
01470         return;
01471       }
01472       nrec = nrec + 1;
01473       LocalPoint rhitlocal = (*recIt).localPosition();
01474       double rhitlocalx = rhitlocal.x();
01475       maprec[rhitlocalx] = nrec; 
01476     }
01477     
01478     int i = 0;
01479     for (std::map<double,int>::iterator iter = maprec.begin();
01480          iter != maprec.end(); ++iter) {
01481       i = i + 1;
01482       nmaprec[i] = (*iter).first;
01483     }
01484     
01485     int nsim = 0;
01486     if (validrpcsim) {
01487       edm::PSimHitContainer::const_iterator simIt;
01488       for (simIt = simHit->begin(); simIt != simHit->end(); simIt++) {
01489         int ptype = (*simIt).particleType();
01490         if (ptype == 13 || ptype == -13) {
01491           nsim = nsim + 1;
01492           LocalPoint shitlocal = (*simIt).localPosition();
01493           double shitlocalx = shitlocal.x();
01494           mapsim[shitlocalx] = nsim;
01495         }
01496       }
01497       
01498       i = 0;
01499       for (std::map<double,int>::iterator iter = mapsim.begin();
01500            iter != mapsim.end(); ++iter) {
01501         i = i + 1;
01502         nmapsim[i] = (*iter).first;
01503       }
01504     }
01505 
01506     if (nsim == nrec) {
01507       for (int r = 0; r < nsim; r++) {
01508         ++nRPC;
01509         mehRPCResX->Fill(nmaprec[r+1]-nmapsim[r+1]);
01510       }
01511     }
01512                                                                   
01513     if (verbosity > 1) {
01514       eventout += "\n          Number of RPCRecHits collected:........... ";
01515       eventout += nRPC;
01516     }
01517     mehRPCn->Fill((float)nRPC);
01518   }
01519 
01520   if (verbosity > 0)
01521     edm::LogInfo(MsgLoggerCat) << eventout << "\n";
01522   
01523   return;
01524 }
01525 
01526 //needed by to do the residual for matched hits in SiStrip
01527 std::pair<LocalPoint,LocalVector> 
01528 GlobalRecHitsAnalyzer::projectHit(const PSimHit& hit, 
01529                                   const StripGeomDetUnit* stripDet,
01530                                   const BoundPlane& plane) 
01531 {
01532   
01533   const StripTopology& topol = stripDet->specificTopology();
01534   GlobalPoint globalpos= stripDet->surface().toGlobal(hit.localPosition());
01535   LocalPoint localHit = plane.toLocal(globalpos);
01536   //track direction
01537   LocalVector locdir=hit.localDirection();
01538   //rotate track in new frame
01539   
01540   GlobalVector globaldir= stripDet->surface().toGlobal(locdir);
01541   LocalVector dir=plane.toLocal(globaldir);
01542   float scale = -localHit.z() / dir.z();
01543   
01544   LocalPoint projectedPos = localHit + scale*dir;
01545     
01546   float selfAngle = topol.stripAngle( topol.strip( hit.localPosition()));
01547 
01548   // vector along strip in hit frame 
01549   LocalVector stripDir( sin(selfAngle), cos(selfAngle), 0); 
01550   
01551   LocalVector 
01552     localStripDir(plane.toLocal(stripDet->surface().toGlobal(stripDir)));
01553   
01554   return std::pair<LocalPoint,LocalVector>( projectedPos, localStripDir);
01555 }
01556 
01557 // Return a map between DTRecHit1DPair and wireId
01558 std::map<DTWireId, std::vector<DTRecHit1DPair> >
01559 GlobalRecHitsAnalyzer::map1DRecHitsPerWire(const DTRecHitCollection* 
01560                                            dt1DRecHitPairs) {
01561   std::map<DTWireId, std::vector<DTRecHit1DPair> > ret;
01562   
01563   for(DTRecHitCollection::const_iterator rechit = dt1DRecHitPairs->begin();
01564       rechit != dt1DRecHitPairs->end(); rechit++) {
01565     ret[(*rechit).wireId()].push_back(*rechit);
01566   }
01567   
01568   return ret;
01569 }
01570 
01571 // Compute SimHit distance from wire (cm)
01572 float GlobalRecHitsAnalyzer::simHitDistFromWire(const DTLayer* layer,
01573                                                 DTWireId wireId,
01574                                                 const PSimHit& hit) {
01575   float xwire = layer->specificTopology().wirePosition(wireId.wire());
01576   LocalPoint entryP = hit.entryPoint();
01577   LocalPoint exitP = hit.exitPoint();
01578   float xEntry = entryP.x()-xwire;
01579   float xExit  = exitP.x()-xwire;
01580 
01581   //FIXME: check...  
01582   return fabs(xEntry - (entryP.z()*(xExit-xEntry))/(exitP.z()-entryP.z()));
01583 }
01584 
01585 // Find the RecHit closest to the muon SimHit
01586 template  <typename type>
01587 const type* 
01588 GlobalRecHitsAnalyzer::findBestRecHit(const DTLayer* layer,
01589                                       DTWireId wireId,
01590                                       const std::vector<type>& recHits,
01591                                       const float simHitDist) {
01592   float res = 99999;
01593   const type* theBestRecHit = 0;
01594   // Loop over RecHits within the cell
01595   for(typename std::vector<type>::const_iterator recHit = recHits.begin();
01596       recHit != recHits.end();
01597       recHit++) {
01598     float distTmp = recHitDistFromWire(*recHit, layer);
01599     if(fabs(distTmp-simHitDist) < res) {
01600       res = fabs(distTmp-simHitDist);
01601       theBestRecHit = &(*recHit);
01602     }
01603   } // End of loop over RecHits within the cell
01604   
01605   return theBestRecHit;
01606 }
01607 
01608 // Compute the distance from wire (cm) of a hits in a DTRecHit1DPair
01609 float 
01610 GlobalRecHitsAnalyzer::recHitDistFromWire(const DTRecHit1DPair& hitPair, 
01611                                           const DTLayer* layer) {
01612   // Compute the rechit distance from wire
01613   return fabs(hitPair.localPosition(DTEnums::Left).x() -
01614               hitPair.localPosition(DTEnums::Right).x())/2.;
01615 }
01616 
01617 // Compute the distance from wire (cm) of a hits in a DTRecHit1D
01618 float 
01619 GlobalRecHitsAnalyzer::recHitDistFromWire(const DTRecHit1D& recHit, 
01620                                           const DTLayer* layer) {
01621   return fabs(recHit.localPosition().x() - 
01622               layer->specificTopology().wirePosition(recHit.wireId().wire()));
01623 }
01624 
01625 template  <typename type>
01626 int GlobalRecHitsAnalyzer::compute(const DTGeometry *dtGeom,
01627                                    std::map<DTWireId, std::vector<PSimHit> > 
01628                                    simHitsPerWire,
01629                                    std::map<DTWireId, std::vector<type> > 
01630                                    recHitsPerWire,
01631                                    int step) {
01632 
01633   int nDt = 0;
01634   // Loop over cells with a muon SimHit
01635   for(std::map<DTWireId, std::vector<PSimHit> >::const_iterator wireAndSHits = 
01636         simHitsPerWire.begin();
01637       wireAndSHits != simHitsPerWire.end();
01638       wireAndSHits++) {
01639     DTWireId wireId = (*wireAndSHits).first;
01640     std::vector<PSimHit> simHitsInCell = (*wireAndSHits).second;
01641     
01642     // Get the layer
01643     const DTLayer* layer = dtGeom->layer(wireId);
01644     
01645     // Look for a mu hit in the cell
01646     const PSimHit* muSimHit = DTHitQualityUtils::findMuSimHit(simHitsInCell);
01647     if (muSimHit==0) {
01648       continue; // Skip this cell
01649     }
01650 
01651     // Find the distance of the simhit from the wire
01652     float simHitWireDist = simHitDistFromWire(layer, wireId, *muSimHit);
01653     // Skip simhits out of the cell
01654     if(simHitWireDist>2.1) {
01655       continue; // Skip this cell
01656     }
01657 
01658     // Look for RecHits in the same cell
01659     if(recHitsPerWire.find(wireId) == recHitsPerWire.end()) {
01660       continue; // No RecHit found in this cell
01661     } else {
01662 
01663       std::vector<type> recHits = recHitsPerWire[wireId];
01664          
01665       // Find the best RecHit
01666       const type* theBestRecHit = 
01667         findBestRecHit(layer, wireId, recHits, simHitWireDist);
01668  
01669       float recHitWireDist =  recHitDistFromWire(*theBestRecHit, layer);
01670       
01671       ++nDt;
01672 
01673       mehDtMuonRes->Fill(recHitWireDist-simHitWireDist);
01674       
01675     } // find rechits
01676   } // loop over simhits
01677 
01678   return nDt;
01679 }
01680 
01681 void 
01682 GlobalRecHitsAnalyzer::plotResolution(const PSimHit & simHit, 
01683                                       const CSCRecHit2D & recHit,
01684                                       const CSCLayer * layer, 
01685                                       int chamberType) {
01686   GlobalPoint simHitPos = layer->toGlobal(simHit.localPosition());
01687   GlobalPoint recHitPos = layer->toGlobal(recHit.localPosition());
01688   
01689   mehCSCResRDPhi->Fill(recHitPos.phi()-simHitPos.phi());
01690 }
01691