CMS 3D CMS Logo

/data/git/CMSSW_5_3_11_patch5/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 maxHFEnergy = 0.;
00673   
00674   Double_t maxHBPhi = -1000.;
00675   Double_t maxHEPhi = -1000.;
00676   Double_t maxHOPhi = -1000.;
00677   Double_t maxHFPhi = -1000.;
00678   
00679   
00680   Double_t PI = 3.141592653589;
00681   
00683   // get HBHE information
00685   std::vector<edm::Handle<HBHERecHitCollection> > hbhe;
00686   iEvent.getManyByType(hbhe);
00687   bool validHBHE = true;
00688   if (!hbhe[0].isValid()) {
00689     LogDebug(MsgLoggerCat)
00690       << "Unable to find any HBHERecHitCollections in event!";
00691     validHBHE = false;
00692   } 
00693 
00694   if (validHBHE) {
00695     std::vector<edm::Handle<HBHERecHitCollection> >::iterator ihbhe;
00696     
00697     int iHB = 0;
00698     int iHE = 0; 
00699     for (ihbhe = hbhe.begin(); ihbhe != hbhe.end(); ++ihbhe) {
00700       
00701       // find max values
00702       for (HBHERecHitCollection::const_iterator jhbhe = (*ihbhe)->begin();
00703            jhbhe != (*ihbhe)->end(); ++jhbhe) {
00704         
00705         HcalDetId cell(jhbhe->id());
00706         
00707         if (cell.subdet() == sdHcalBrl) {
00708           
00709           const CaloCellGeometry* cellGeometry =
00710             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00711           double fPhi = cellGeometry->getPosition().phi () ;
00712           if ( (jhbhe->energy()) > maxHBEnergy ) {
00713             maxHBEnergy = jhbhe->energy();
00714             maxHBPhi = fPhi;
00715             maxHOPhi = maxHBPhi;
00716           }       
00717         }
00718       
00719         if (cell.subdet() == sdHcalEC) {
00720           
00721           const CaloCellGeometry* cellGeometry =
00722             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00723           double fPhi = cellGeometry->getPosition().phi () ;
00724           if ( (jhbhe->energy()) > maxHEEnergy ) {
00725             maxHEEnergy = jhbhe->energy();
00726             maxHEPhi = fPhi;
00727           }       
00728         }
00729       } // end find max values
00730       
00731       for (HBHERecHitCollection::const_iterator jhbhe = (*ihbhe)->begin();
00732            jhbhe != (*ihbhe)->end(); ++jhbhe) {
00733         
00734         HcalDetId cell(jhbhe->id());
00735         
00736         if (cell.subdet() == sdHcalBrl) {
00737           
00738           ++iHB;
00739           
00740           const CaloCellGeometry* cellGeometry =
00741             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00742           double fPhi = cellGeometry->getPosition().phi () ;
00743           
00744           float deltaphi = maxHBPhi - fPhi;
00745           if (fPhi > maxHBPhi) { deltaphi = fPhi - maxHBPhi;}
00746           if (deltaphi > PI) { deltaphi = 2.0 * PI - deltaphi;}
00747           
00748           mehHcalRes[0]->Fill(jhbhe->energy() - 
00749                               fHBEnergySimHits[cell.rawId()]);
00750         }
00751         
00752         if (cell.subdet() == sdHcalEC) {
00753           
00754           ++iHE;
00755           
00756           const CaloCellGeometry* cellGeometry =
00757             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00758           double fPhi = cellGeometry->getPosition().phi () ;
00759           
00760           float deltaphi = maxHEPhi - fPhi;
00761           if (fPhi > maxHEPhi) { deltaphi = fPhi - maxHEPhi;}
00762           if (deltaphi > PI) { deltaphi = 2.0 * PI - deltaphi;}
00763           mehHcalRes[1]->Fill(jhbhe->energy() - 
00764                               fHEEnergySimHits[cell.rawId()]);
00765         }
00766       }
00767     } // end loop through collection
00768   
00769     
00770     if (verbosity > 1) {
00771       eventout += "\n          Number of HBRecHits collected:............ ";
00772       eventout += iHB;
00773     }
00774     
00775     if (verbosity > 1) {
00776       eventout += "\n          Number of HERecHits collected:............ ";
00777       eventout += iHE;
00778     }
00779     mehHcaln[0]->Fill((float)iHB);
00780     mehHcaln[1]->Fill((float)iHE);
00781   }
00782 
00784   // get HF information
00786   std::vector<edm::Handle<HFRecHitCollection> > hf;
00787   iEvent.getManyByType(hf);
00788   bool validHF = true;
00789   if (!hf[0].isValid()) {
00790     LogDebug(MsgLoggerCat)
00791       << "Unable to find any HFRecHitCollections in event!";
00792     validHF = false;
00793   } 
00794   if (validHF) {
00795     std::vector<edm::Handle<HFRecHitCollection> >::iterator ihf;
00796     
00797     int iHF = 0; 
00798     for (ihf = hf.begin(); ihf != hf.end(); ++ihf) {
00799       
00800       // find max values
00801       for (HFRecHitCollection::const_iterator jhf = (*ihf)->begin();
00802            jhf != (*ihf)->end(); ++jhf) {
00803         
00804         HcalDetId cell(jhf->id());
00805         
00806         if (cell.subdet() == sdHcalFwd) {
00807           
00808           const CaloCellGeometry* cellGeometry =
00809             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00810           double fPhi = cellGeometry->getPosition().phi () ;
00811           if ( (jhf->energy()) > maxHFEnergy ) {
00812             maxHFEnergy = jhf->energy();
00813             maxHFPhi = fPhi;
00814           }       
00815         }
00816       } // end find max values
00817       
00818       for (HFRecHitCollection::const_iterator jhf = (*ihf)->begin();
00819            jhf != (*ihf)->end(); ++jhf) {
00820         
00821         HcalDetId cell(jhf->id());
00822         
00823         if (cell.subdet() == sdHcalFwd) {
00824           
00825           ++iHF;
00826           
00827           const CaloCellGeometry* cellGeometry =
00828             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00829           double fPhi = cellGeometry->getPosition().phi () ;
00830           
00831           float deltaphi = maxHBPhi - fPhi;
00832           if (fPhi > maxHFPhi) { deltaphi = fPhi - maxHFPhi;}
00833           if (deltaphi > PI) { deltaphi = 2.0 * PI - deltaphi;}
00834           
00835           mehHcalRes[2]->Fill(jhf->energy()-fHFEnergySimHits[cell.rawId()]);
00836         }
00837       }
00838     } // end loop through collection
00839     
00840     if (verbosity > 1) {
00841       eventout += "\n          Number of HFDigis collected:.............. ";
00842       eventout += iHF;
00843     }
00844     mehHcaln[2]->Fill((float)iHF);
00845   }    
00846 
00848   // get HO information
00850   std::vector<edm::Handle<HORecHitCollection> > ho;
00851   iEvent.getManyByType(ho);
00852   bool validHO = true;
00853   if (!ho[0].isValid()) {
00854     LogDebug(MsgLoggerCat)
00855       << "Unable to find any HORecHitCollections in event!";
00856     validHO = false;
00857   } 
00858 
00859   if (validHO) {
00860     std::vector<edm::Handle<HORecHitCollection> >::iterator iho;
00861     
00862     int iHO = 0; 
00863     for (iho = ho.begin(); iho != ho.end(); ++iho) {
00864       
00865       for (HORecHitCollection::const_iterator jho = (*iho)->begin();
00866            jho != (*iho)->end(); ++jho) {
00867         
00868         HcalDetId cell(jho->id());
00869         
00870         if (cell.subdet() == sdHcalOut) {
00871           
00872           ++iHO;
00873           
00874           const CaloCellGeometry* cellGeometry =
00875             geometry->getSubdetectorGeometry (cell)->getGeometry (cell) ;
00876           double fPhi = cellGeometry->getPosition().phi () ;
00877           
00878           float deltaphi = maxHOPhi - fPhi;
00879           if (fPhi > maxHOPhi) { deltaphi = fPhi - maxHOPhi;}
00880           if (deltaphi > PI) { deltaphi = 2.0 * PI - deltaphi;}
00881           mehHcalRes[3]->Fill(jho->energy()-fHOEnergySimHits[cell.rawId()]);
00882         }
00883       }
00884     } // end loop through collection
00885     
00886     if (verbosity > 1) {
00887       eventout += "\n          Number of HODigis collected:.............. ";
00888       eventout += iHO;
00889     }
00890     mehHcaln[3]->Fill((float)iHO);
00891   }
00892 
00893   if (verbosity > 0)
00894     edm::LogInfo(MsgLoggerCat) << eventout << "\n";
00895   
00896   return;
00897 }
00898 
00899 void GlobalRecHitsAnalyzer::fillTrk(const edm::Event& iEvent, 
00900                                     const edm::EventSetup& iSetup)
00901 {
00902   std::string MsgLoggerCat = "GlobalRecHitsAnalyzer_fillTrk";
00903   
00904   TString eventout;
00905   if (verbosity > 0)
00906     eventout = "\nGathering info:";  
00907   
00908   // get strip information
00909   edm::Handle<SiStripMatchedRecHit2DCollection> rechitsmatched;
00910   iEvent.getByLabel(SiStripSrc_, rechitsmatched);
00911   bool validstrip = true;
00912   if (!rechitsmatched.isValid()) {
00913     LogDebug(MsgLoggerCat)
00914       << "Unable to find stripmatchedrechits in event!";
00915     validstrip = false;
00916   }  
00917   
00918   TrackerHitAssociator associate(iEvent,conf_);
00919   
00920   edm::ESHandle<TrackerGeometry> pDD;
00921   iSetup.get<TrackerDigiGeometryRecord>().get(pDD);
00922   if (!pDD.isValid()) {
00923     edm::LogWarning(MsgLoggerCat)
00924       << "Unable to find TrackerDigiGeometry in event!";
00925     return;
00926   }
00927   const TrackerGeometry &tracker(*pDD);
00928   
00929   if (validstrip) {
00930     int nStripBrl = 0, nStripFwd = 0;
00931     
00932     // loop over det units
00933     for (TrackerGeometry::DetContainer::const_iterator it = 
00934            pDD->dets().begin();
00935          it != pDD->dets().end(); ++it) {
00936       
00937       uint32_t myid = ((*it)->geographicalId()).rawId();
00938       DetId detid = ((*it)->geographicalId());
00939       
00940       //loop over rechits-matched in the same subdetector
00941       SiStripMatchedRecHit2DCollection::const_iterator rechitmatchedMatch = rechitsmatched->find(detid);
00942       
00943       if (rechitmatchedMatch != rechitsmatched->end()) {
00944         SiStripMatchedRecHit2DCollection::DetSet rechitmatchedRange = *rechitmatchedMatch;
00945         SiStripMatchedRecHit2DCollection::DetSet::const_iterator rechitmatchedRangeIteratorBegin = rechitmatchedRange.begin();
00946         SiStripMatchedRecHit2DCollection::DetSet::const_iterator rechitmatchedRangeIteratorEnd   = rechitmatchedRange.end();
00947         SiStripMatchedRecHit2DCollection::DetSet::const_iterator itermatched = rechitmatchedRangeIteratorBegin;
00948         
00949         for ( itermatched = rechitmatchedRangeIteratorBegin; 
00950               itermatched != rechitmatchedRangeIteratorEnd;
00951               ++itermatched) {
00952           
00953           SiStripMatchedRecHit2D const rechit = *itermatched;
00954           LocalPoint position = rechit.localPosition();
00955           
00956           float mindist = 999999.;
00957           float distx = 999999.;
00958           float disty = 999999.;
00959           float dist = 999999.;
00960           std::pair<LocalPoint,LocalVector> closestPair;
00961           matched.clear();
00962           
00963           float rechitmatchedx = position.x();
00964           float rechitmatchedy = position.y();
00965           
00966           matched = associate.associateHit(rechit);
00967           
00968           if (!matched.empty()) {
00969             //project simhit;
00970             const GluedGeomDet* gluedDet = 
00971               (const GluedGeomDet*)tracker.idToDet(rechit.geographicalId());
00972             const StripGeomDetUnit* partnerstripdet =
00973               (StripGeomDetUnit*) gluedDet->stereoDet();
00974             std::pair<LocalPoint,LocalVector> hitPair;
00975             
00976             for(std::vector<PSimHit>::const_iterator m = matched.begin(); 
00977                 m != matched.end(); m++){
00978               //project simhit;
00979               hitPair = projectHit((*m),partnerstripdet,gluedDet->surface());
00980               distx = fabs(rechitmatchedx - hitPair.first.x());
00981               disty = fabs(rechitmatchedy - hitPair.first.y());
00982               dist = sqrt(distx*distx+disty*disty);
00983               
00984               if(dist < mindist){
00985                 mindist = dist;
00986                 closestPair = hitPair;
00987               }
00988             }
00989             
00990             // get TIB
00991             if (detid.subdetId() == sdSiTIB) {
00992               
00993               TIBDetId tibid(myid);
00994               ++nStripBrl;
00995               
00996               if (tibid.layer() == 1) {
00997                 mehSiStripResX[8]->Fill(rechitmatchedx-closestPair.first.x());
00998                 mehSiStripResY[8]->Fill(rechitmatchedy-closestPair.first.y());
00999               }
01000               if (tibid.layer() == 2) {
01001                 mehSiStripResX[9]->Fill(rechitmatchedx-closestPair.first.x());
01002                 mehSiStripResY[9]->Fill(rechitmatchedy-closestPair.first.y());
01003               } 
01004               if (tibid.layer() == 3) {
01005                 mehSiStripResX[10]->Fill(rechitmatchedx-closestPair.first.x());
01006                 mehSiStripResY[10]->Fill(rechitmatchedy-closestPair.first.y());
01007                 
01008               }
01009               if (tibid.layer() == 4) {
01010                 mehSiStripResX[11]->Fill(rechitmatchedx-closestPair.first.x());
01011                 mehSiStripResY[11]->Fill(rechitmatchedy-closestPair.first.y());
01012               }
01013             }
01014             
01015             // get TOB
01016             if (detid.subdetId() == sdSiTOB) {
01017               
01018               TOBDetId tobid(myid);
01019               ++nStripBrl;
01020               
01021               if (tobid.layer() == 1) {
01022                 mehSiStripResX[15]->Fill(rechitmatchedx-closestPair.first.x());
01023                 mehSiStripResY[15]->Fill(rechitmatchedy-closestPair.first.y());
01024               }
01025               if (tobid.layer() == 2) {
01026                 mehSiStripResX[16]->Fill(rechitmatchedx-closestPair.first.x());
01027                 mehSiStripResY[16]->Fill(rechitmatchedy-closestPair.first.y());
01028               } 
01029               if (tobid.layer() == 3) {
01030                 mehSiStripResX[17]->Fill(rechitmatchedx-closestPair.first.x());
01031                 mehSiStripResY[17]->Fill(rechitmatchedy-closestPair.first.y());
01032               }
01033               if (tobid.layer() == 4) {
01034                 mehSiStripResX[18]->Fill(rechitmatchedx-closestPair.first.x());
01035                 mehSiStripResY[18]->Fill(rechitmatchedy-closestPair.first.y());
01036               }
01037             }
01038             
01039             // get TID
01040             if (detid.subdetId() == sdSiTID) {
01041               
01042               TIDDetId tidid(myid);
01043               ++nStripFwd;
01044               
01045               if (tidid.wheel() == 1) {
01046                 mehSiStripResX[12]->Fill(rechitmatchedx-closestPair.first.x());
01047                 mehSiStripResY[12]->Fill(rechitmatchedy-closestPair.first.y());
01048               }
01049               if (tidid.wheel() == 2) {
01050                 mehSiStripResX[13]->Fill(rechitmatchedx-closestPair.first.x());
01051                 mehSiStripResY[13]->Fill(rechitmatchedy-closestPair.first.y());
01052               } 
01053               if (tidid.wheel() == 3) {
01054                 mehSiStripResX[14]->Fill(rechitmatchedx-closestPair.first.x());
01055                 mehSiStripResY[14]->Fill(rechitmatchedy-closestPair.first.y());
01056               }
01057             }
01058             
01059             // get TEC
01060             if (detid.subdetId() == sdSiTEC) {
01061               
01062               TECDetId tecid(myid);
01063               ++nStripFwd;
01064               
01065               if (tecid.wheel() == 1) {
01066                 mehSiStripResX[0]->Fill(rechitmatchedx-closestPair.first.x());
01067                 mehSiStripResY[0]->Fill(rechitmatchedy-closestPair.first.y());
01068               }
01069               if (tecid.wheel() == 2) {
01070                 mehSiStripResX[1]->Fill(rechitmatchedx-closestPair.first.x());
01071                 mehSiStripResY[1]->Fill(rechitmatchedy-closestPair.first.y());
01072               } 
01073               if (tecid.wheel() == 3) {
01074                 mehSiStripResX[2]->Fill(rechitmatchedx-closestPair.first.x());
01075                 mehSiStripResY[2]->Fill(rechitmatchedy-closestPair.first.y());
01076               }
01077               if (tecid.wheel() == 4) {
01078                 mehSiStripResX[3]->Fill(rechitmatchedx-closestPair.first.x());
01079                 mehSiStripResY[3]->Fill(rechitmatchedy-closestPair.first.y());
01080                 
01081               }
01082               if (tecid.wheel() == 5) {
01083                 mehSiStripResX[4]->Fill(rechitmatchedx-closestPair.first.x());
01084                 mehSiStripResY[4]->Fill(rechitmatchedy-closestPair.first.y());
01085               } 
01086               if (tecid.wheel() == 6) {
01087                 mehSiStripResX[5]->Fill(rechitmatchedx-closestPair.first.x());
01088                 mehSiStripResY[5]->Fill(rechitmatchedy-closestPair.first.y());
01089               }
01090               if (tecid.wheel() == 7) {
01091                 mehSiStripResX[6]->Fill(rechitmatchedx-closestPair.first.x());
01092                 mehSiStripResY[6]->Fill(rechitmatchedy-closestPair.first.y());
01093               } 
01094               if (tecid.wheel() == 8) {
01095                 mehSiStripResX[7]->Fill(rechitmatchedx-closestPair.first.x());
01096                 mehSiStripResY[7]->Fill(rechitmatchedy-closestPair.first.y()); 
01097               }
01098             }
01099             
01100           } // end if matched empty
01101         }
01102       }
01103     } // end loop over det units
01104                                                                       
01105     if (verbosity > 1) {
01106       eventout += "\n          Number of BrlStripRecHits collected:...... ";
01107       eventout += nStripBrl;
01108     }
01109     
01110     for(int i =8; i<12; ++i)
01111       {mehSiStripn[i]->Fill((float)nStripBrl);}
01112     for(int i =16; i<19; ++i)
01113       {mehSiStripn[i]->Fill((float)nStripBrl);}
01114     
01115     if (verbosity > 1) {
01116       eventout += "\n          Number of FrwdStripRecHits collected:..... ";
01117       eventout += nStripFwd;
01118     }
01119     for(int i =0; i<8; ++i)
01120       {mehSiStripn[i]->Fill((float)nStripFwd);}
01121     for(int i =12; i<16; ++i)
01122       {mehSiStripn[i]->Fill((float)nStripFwd);}
01123   }
01124 
01125   // get pixel information
01126   //Get RecHits
01127   edm::Handle<SiPixelRecHitCollection> recHitColl;
01128   iEvent.getByLabel(SiPxlSrc_, recHitColl);
01129   bool validpixel = true;
01130   if (!recHitColl.isValid()) {
01131     LogDebug(MsgLoggerCat)
01132       << "Unable to find SiPixelRecHitCollection in event!";
01133     validpixel = false;
01134   }  
01135   
01136   //Get event setup
01137   edm::ESHandle<TrackerGeometry> geom;
01138   iSetup.get<TrackerDigiGeometryRecord>().get(geom); 
01139   if (!geom.isValid()) {
01140     edm::LogWarning(MsgLoggerCat)
01141       << "Unable to find TrackerDigiGeometry in event!";
01142     return;
01143   }
01144 
01145   if (validpixel) {
01146     int nPxlBrl = 0, nPxlFwd = 0;    
01147     //iterate over detunits
01148     for (TrackerGeometry::DetContainer::const_iterator it = 
01149            geom->dets().begin();
01150          it != geom->dets().end(); ++it) {
01151       
01152       uint32_t myid = ((*it)->geographicalId()).rawId();
01153       DetId detId = ((*it)->geographicalId());
01154       int subid = detId.subdetId();
01155       
01156       if (! ((subid == sdPxlBrl) || (subid == sdPxlFwd))) continue;
01157       
01158       SiPixelRecHitCollection::const_iterator pixeldet = recHitColl->find(detId);
01159       if (pixeldet == recHitColl->end()) continue;
01160       SiPixelRecHitCollection::DetSet pixelrechitRange = *pixeldet;
01161       SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorBegin = pixelrechitRange.begin();
01162       SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorEnd   = pixelrechitRange.end();
01163       SiPixelRecHitCollection::DetSet::const_iterator pixeliter = pixelrechitRangeIteratorBegin;
01164 
01165       
01166       std::vector<PSimHit> matched;
01167       
01168       //----Loop over rechits for this detId
01169       for ( ; pixeliter != pixelrechitRangeIteratorEnd; ++pixeliter) {
01170         
01171         matched.clear();
01172         matched = associate.associateHit(*pixeliter);
01173         
01174         if ( !matched.empty() ) {
01175           
01176           float closest = 9999.9;
01177           LocalPoint lp = pixeliter->localPosition();
01178           float rechit_x = lp.x();
01179           float rechit_y = lp.y();
01180           
01181           float sim_x = 0.;
01182           float sim_y = 0.;
01183           
01184           //loop over sim hits and fill closet
01185           for (std::vector<PSimHit>::const_iterator m = matched.begin(); 
01186                m != matched.end(); ++m) {
01187             
01188             float sim_x1 = (*m).entryPoint().x();
01189             float sim_x2 = (*m).exitPoint().x();
01190             float sim_xpos = 0.5*(sim_x1+sim_x2);
01191             
01192             float sim_y1 = (*m).entryPoint().y();
01193             float sim_y2 = (*m).exitPoint().y();
01194             float sim_ypos = 0.5*(sim_y1+sim_y2);
01195             
01196             float x_res = fabs(sim_xpos - rechit_x);
01197             float y_res = fabs(sim_ypos - rechit_y);
01198             
01199             float dist = sqrt(x_res*x_res + y_res*y_res);
01200             
01201             if ( dist < closest ) {
01202               closest = dist;
01203               sim_x = sim_xpos;
01204               sim_y = sim_ypos;
01205             }
01206           } // end sim hit loop
01207           
01208           // get Barrel pixels ***************Pixel STuff******************
01209           if (subid == sdPxlBrl) {
01210             PXBDetId bdetid(myid);
01211             ++nPxlBrl;
01212             
01213             if (bdetid.layer() == 1) {
01214               mehSiPixelResX[0]->Fill(rechit_x-sim_x);
01215               mehSiPixelResY[0]->Fill(rechit_y-sim_y); 
01216               
01217             }
01218             if (bdetid.layer() == 2) {
01219               mehSiPixelResX[1]->Fill(rechit_x-sim_x);
01220               mehSiPixelResY[1]->Fill(rechit_y-sim_y); 
01221             }
01222             if (bdetid.layer() == 3) {
01223               mehSiPixelResX[2]->Fill(rechit_x-sim_x);
01224               mehSiPixelResY[2]->Fill(rechit_y-sim_y); 
01225             }
01226           }
01227           
01228           // get Forward pixels
01229           if (subid == sdPxlFwd) {
01230             PXFDetId fdetid(myid);
01231             ++nPxlFwd;
01232             
01233             if (fdetid.disk() == 1) {
01234               if (fdetid.side() == 1) {
01235                 mehSiPixelResX[3]->Fill(rechit_x-sim_x);
01236                 mehSiPixelResY[3]->Fill(rechit_y-sim_y);
01237               }
01238               if (fdetid.side() == 2) {
01239                 mehSiPixelResX[4]->Fill(rechit_x-sim_x);
01240                 mehSiPixelResY[4]->Fill(rechit_y-sim_y); 
01241               }
01242             }
01243             if (fdetid.disk() == 2) {
01244               if (fdetid.side() == 1) {
01245                 mehSiPixelResX[5]->Fill(rechit_x-sim_x);
01246                 mehSiPixelResY[5]->Fill(rechit_y-sim_y);
01247               }
01248               if (fdetid.side() == 2) {
01249                 mehSiPixelResX[6]->Fill(rechit_x-sim_x);
01250                 mehSiPixelResY[6]->Fill(rechit_y-sim_y); 
01251               }
01252             }
01253           }      
01254         } // end matched emtpy
01255       } // <-----end rechit loop 
01256     } // <------ end detunit loop  
01257     
01258     
01259     if (verbosity > 1) {
01260       eventout += "\n          Number of BrlPixelRecHits collected:...... ";
01261       eventout += nPxlBrl;
01262     }
01263     for(int i=0; i<3; ++i) {
01264       mehSiPixeln[i]->Fill((float)nPxlBrl);
01265     }
01266     
01267     if (verbosity > 1) {
01268       eventout += "\n          Number of FrwdPixelRecHits collected:..... ";
01269       eventout += nPxlFwd;
01270     }
01271     
01272     for(int i=3; i<7; ++i) {
01273       mehSiPixeln[i]->Fill((float)nPxlFwd);
01274     }
01275   }
01276    
01277   if (verbosity > 0)
01278     edm::LogInfo(MsgLoggerCat) << eventout << "\n";
01279   
01280   return;
01281 }
01282 
01283 void GlobalRecHitsAnalyzer::fillMuon(const edm::Event& iEvent, 
01284                                      const edm::EventSetup& iSetup)
01285 {
01286   std::string MsgLoggerCat = "GlobalRecHitsAnalyzer_fillMuon";
01287   
01288   TString eventout;
01289   if (verbosity > 0)
01290     eventout = "\nGathering info:";  
01291 
01292   // get DT information
01293   edm::ESHandle<DTGeometry> dtGeom;
01294   iSetup.get<MuonGeometryRecord>().get(dtGeom);
01295   if (!dtGeom.isValid()) {
01296     edm::LogWarning(MsgLoggerCat)
01297       << "Unable to find DTMuonGeometryRecord in event!";
01298     return;
01299   }  
01300 
01301   edm::Handle<edm::PSimHitContainer> dtsimHits;
01302   iEvent.getByLabel(MuDTSimSrc_, dtsimHits);
01303   bool validdtsim = true;
01304   if (!dtsimHits.isValid()) {
01305     LogDebug(MsgLoggerCat)
01306       << "Unable to find dtsimHits in event!";
01307     validdtsim = false;
01308   } 
01309 
01310   edm::Handle<DTRecHitCollection> dtRecHits;
01311   iEvent.getByLabel(MuDTSrc_, dtRecHits);
01312   bool validdtrec = true;
01313   if (!dtRecHits.isValid()) {
01314     LogDebug(MsgLoggerCat)
01315       << "Unable to find dtRecHits in event!";
01316     validdtrec = false;
01317   }   
01318     
01319   if (validdtsim && validdtrec) {
01320 
01321     std::map<DTWireId, edm::PSimHitContainer> simHitsPerWire =
01322       DTHitQualityUtils::mapSimHitsPerWire(*(dtsimHits.product()));
01323 
01324     std::map<DTWireId, std::vector<DTRecHit1DPair> > recHitsPerWire =
01325       map1DRecHitsPerWire(dtRecHits.product());
01326 
01327     int nDt = compute(dtGeom.product(), simHitsPerWire, recHitsPerWire, 1);
01328     
01329     if (verbosity > 1) {
01330       eventout += "\n          Number of DtMuonRecHits collected:........ ";
01331       eventout += nDt;
01332     }
01333     mehDtMuonn->Fill(float(nDt));
01334   }
01335 
01336   // get CSC Strip information
01337   // get map of sim hits
01338   theMap.clear();
01339   edm::Handle<CrossingFrame<PSimHit> > cf;
01340 
01341   iEvent.getByLabel("mix",hitsProducer+"MuonCSCHits",cf);
01342   bool validXFrame = true;
01343   if (!cf.isValid()) {
01344     LogDebug(MsgLoggerCat)
01345       << "Unable to find muo CSC crossingFrame in event!";
01346     validXFrame = false;
01347   }
01348   if (validXFrame) {
01349     MixCollection<PSimHit> simHits(cf.product());
01350     
01351     // arrange the hits by detUnit
01352     for(MixCollection<PSimHit>::MixItr hitItr = simHits.begin();
01353         hitItr != simHits.end(); ++hitItr) {
01354       theMap[hitItr->detUnitId()].push_back(*hitItr);
01355     }  
01356   }
01357 
01358   // get geometry
01359   edm::ESHandle<CSCGeometry> hGeom;
01360   iSetup.get<MuonGeometryRecord>().get(hGeom);
01361   if (!hGeom.isValid()) {
01362     edm::LogWarning(MsgLoggerCat)
01363       << "Unable to find CSCMuonGeometryRecord in event!";
01364     return;
01365   }    
01366   const CSCGeometry *theCSCGeometry = &*hGeom;
01367 
01368   // get rechits
01369   edm::Handle<CSCRecHit2DCollection> hRecHits;
01370   iEvent.getByLabel(MuCSCSrc_, hRecHits);
01371   bool validCSC = true;
01372   if (!hRecHits.isValid()) {
01373     LogDebug(MsgLoggerCat)
01374       << "Unable to find CSC RecHits in event!";
01375     validCSC = false;
01376   }    
01377 
01378   if (validCSC) {
01379     const CSCRecHit2DCollection *cscRecHits = hRecHits.product();
01380     
01381     int nCSC = 0;
01382     for (CSCRecHit2DCollection::const_iterator recHitItr = cscRecHits->begin();
01383          recHitItr != cscRecHits->end(); ++recHitItr) {
01384       
01385       int detId = (*recHitItr).cscDetId().rawId();
01386       
01387       edm::PSimHitContainer simHits;   
01388       std::map<int, edm::PSimHitContainer>::const_iterator mapItr = 
01389         theMap.find(detId);
01390       if (mapItr != theMap.end()) {
01391         simHits = mapItr->second;
01392       }
01393       
01394       if (simHits.size() == 1) {
01395         ++nCSC;
01396         
01397         const GeomDetUnit* detUnit = 
01398           theCSCGeometry->idToDetUnit(CSCDetId(detId));
01399         const CSCLayer *layer = dynamic_cast<const CSCLayer *>(detUnit); 
01400         
01401         int chamberType = layer->chamber()->specs()->chamberType();
01402         plotResolution(simHits[0], *recHitItr, layer, chamberType);
01403       }
01404     }
01405     
01406     if (verbosity > 1) {
01407       eventout += "\n          Number of CSCRecHits collected:........... ";
01408       eventout += nCSC;
01409     }
01410     mehCSCn->Fill((float)nCSC);
01411   }
01412 
01413   // get RPC information
01414   std::map<double, int> mapsim, maprec;
01415   std::map<int, double> nmapsim, nmaprec;
01416 
01417   edm::ESHandle<RPCGeometry> rpcGeom;
01418   iSetup.get<MuonGeometryRecord>().get(rpcGeom);
01419   if (!rpcGeom.isValid()) {
01420     edm::LogWarning(MsgLoggerCat)
01421       << "Unable to find RPCMuonGeometryRecord in event!";
01422     return;
01423   }  
01424 
01425   edm::Handle<edm::PSimHitContainer> simHit;
01426   iEvent.getByLabel(MuRPCSimSrc_, simHit);
01427   bool validrpcsim = true;
01428   if (!simHit.isValid()) {
01429     LogDebug(MsgLoggerCat)
01430       << "Unable to find RPCSimHit in event!";
01431     validrpcsim = false;
01432   }    
01433 
01434   edm::Handle<RPCRecHitCollection> recHit;
01435   iEvent.getByLabel(MuRPCSrc_, recHit);
01436   bool validrpc = true;
01437   if (!simHit.isValid()) {
01438     LogDebug(MsgLoggerCat)
01439       << "Unable to find RPCRecHit in event!";
01440     validrpc = false;
01441   } 
01442 
01443   if (validrpc) {
01444     int nRPC = 0;
01445     RPCRecHitCollection::const_iterator recIt;
01446     int nrec = 0;
01447     for (recIt = recHit->begin(); recIt != recHit->end(); ++recIt) {
01448       RPCDetId Rid = (RPCDetId)(*recIt).rpcId();
01449       const RPCRoll *roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(Rid));
01450       if (roll->isForward()) {
01451         
01452         if (verbosity > 1) {
01453           eventout += 
01454             "\n          Number of RPCRecHits collected:........... ";
01455           eventout += nRPC;
01456         }
01457         
01458         if (verbosity > 0)
01459           edm::LogInfo(MsgLoggerCat) << eventout << "\n";
01460         return;
01461       }
01462       nrec = nrec + 1;
01463       LocalPoint rhitlocal = (*recIt).localPosition();
01464       double rhitlocalx = rhitlocal.x();
01465       maprec[rhitlocalx] = nrec; 
01466     }
01467     
01468     int i = 0;
01469     for (std::map<double,int>::iterator iter = maprec.begin();
01470          iter != maprec.end(); ++iter) {
01471       i = i + 1;
01472       nmaprec[i] = (*iter).first;
01473     }
01474     
01475     int nsim = 0;
01476     if (validrpcsim) {
01477       edm::PSimHitContainer::const_iterator simIt;
01478       for (simIt = simHit->begin(); simIt != simHit->end(); simIt++) {
01479         int ptype = (*simIt).particleType();
01480         if (ptype == 13 || ptype == -13) {
01481           nsim = nsim + 1;
01482           LocalPoint shitlocal = (*simIt).localPosition();
01483           double shitlocalx = shitlocal.x();
01484           mapsim[shitlocalx] = nsim;
01485         }
01486       }
01487       
01488       i = 0;
01489       for (std::map<double,int>::iterator iter = mapsim.begin();
01490            iter != mapsim.end(); ++iter) {
01491         i = i + 1;
01492         nmapsim[i] = (*iter).first;
01493       }
01494     }
01495 
01496     if (nsim == nrec) {
01497       for (int r = 0; r < nsim; r++) {
01498         ++nRPC;
01499         mehRPCResX->Fill(nmaprec[r+1]-nmapsim[r+1]);
01500       }
01501     }
01502                                                                   
01503     if (verbosity > 1) {
01504       eventout += "\n          Number of RPCRecHits collected:........... ";
01505       eventout += nRPC;
01506     }
01507     mehRPCn->Fill((float)nRPC);
01508   }
01509 
01510   if (verbosity > 0)
01511     edm::LogInfo(MsgLoggerCat) << eventout << "\n";
01512   
01513   return;
01514 }
01515 
01516 //needed by to do the residual for matched hits in SiStrip
01517 std::pair<LocalPoint,LocalVector> 
01518 GlobalRecHitsAnalyzer::projectHit(const PSimHit& hit, 
01519                                   const StripGeomDetUnit* stripDet,
01520                                   const BoundPlane& plane) 
01521 {
01522   
01523   const StripTopology& topol = stripDet->specificTopology();
01524   GlobalPoint globalpos= stripDet->surface().toGlobal(hit.localPosition());
01525   LocalPoint localHit = plane.toLocal(globalpos);
01526   //track direction
01527   LocalVector locdir=hit.localDirection();
01528   //rotate track in new frame
01529   
01530   GlobalVector globaldir= stripDet->surface().toGlobal(locdir);
01531   LocalVector dir=plane.toLocal(globaldir);
01532   float scale = -localHit.z() / dir.z();
01533   
01534   LocalPoint projectedPos = localHit + scale*dir;
01535     
01536   float selfAngle = topol.stripAngle( topol.strip( hit.localPosition()));
01537 
01538   // vector along strip in hit frame 
01539   LocalVector stripDir( sin(selfAngle), cos(selfAngle), 0); 
01540   
01541   LocalVector 
01542     localStripDir(plane.toLocal(stripDet->surface().toGlobal(stripDir)));
01543   
01544   return std::pair<LocalPoint,LocalVector>( projectedPos, localStripDir);
01545 }
01546 
01547 // Return a map between DTRecHit1DPair and wireId
01548 std::map<DTWireId, std::vector<DTRecHit1DPair> >
01549 GlobalRecHitsAnalyzer::map1DRecHitsPerWire(const DTRecHitCollection* 
01550                                            dt1DRecHitPairs) {
01551   std::map<DTWireId, std::vector<DTRecHit1DPair> > ret;
01552   
01553   for(DTRecHitCollection::const_iterator rechit = dt1DRecHitPairs->begin();
01554       rechit != dt1DRecHitPairs->end(); rechit++) {
01555     ret[(*rechit).wireId()].push_back(*rechit);
01556   }
01557   
01558   return ret;
01559 }
01560 
01561 // Compute SimHit distance from wire (cm)
01562 float GlobalRecHitsAnalyzer::simHitDistFromWire(const DTLayer* layer,
01563                                                 DTWireId wireId,
01564                                                 const PSimHit& hit) {
01565   float xwire = layer->specificTopology().wirePosition(wireId.wire());
01566   LocalPoint entryP = hit.entryPoint();
01567   LocalPoint exitP = hit.exitPoint();
01568   float xEntry = entryP.x()-xwire;
01569   float xExit  = exitP.x()-xwire;
01570 
01571   //FIXME: check...  
01572   return fabs(xEntry - (entryP.z()*(xExit-xEntry))/(exitP.z()-entryP.z()));
01573 }
01574 
01575 // Find the RecHit closest to the muon SimHit
01576 template  <typename type>
01577 const type* 
01578 GlobalRecHitsAnalyzer::findBestRecHit(const DTLayer* layer,
01579                                       DTWireId wireId,
01580                                       const std::vector<type>& recHits,
01581                                       const float simHitDist) {
01582   float res = 99999;
01583   const type* theBestRecHit = 0;
01584   // Loop over RecHits within the cell
01585   for(typename std::vector<type>::const_iterator recHit = recHits.begin();
01586       recHit != recHits.end();
01587       recHit++) {
01588     float distTmp = recHitDistFromWire(*recHit, layer);
01589     if(fabs(distTmp-simHitDist) < res) {
01590       res = fabs(distTmp-simHitDist);
01591       theBestRecHit = &(*recHit);
01592     }
01593   } // End of loop over RecHits within the cell
01594   
01595   return theBestRecHit;
01596 }
01597 
01598 // Compute the distance from wire (cm) of a hits in a DTRecHit1DPair
01599 float 
01600 GlobalRecHitsAnalyzer::recHitDistFromWire(const DTRecHit1DPair& hitPair, 
01601                                           const DTLayer* layer) {
01602   // Compute the rechit distance from wire
01603   return fabs(hitPair.localPosition(DTEnums::Left).x() -
01604               hitPair.localPosition(DTEnums::Right).x())/2.;
01605 }
01606 
01607 // Compute the distance from wire (cm) of a hits in a DTRecHit1D
01608 float 
01609 GlobalRecHitsAnalyzer::recHitDistFromWire(const DTRecHit1D& recHit, 
01610                                           const DTLayer* layer) {
01611   return fabs(recHit.localPosition().x() - 
01612               layer->specificTopology().wirePosition(recHit.wireId().wire()));
01613 }
01614 
01615 template  <typename type>
01616 int GlobalRecHitsAnalyzer::compute(const DTGeometry *dtGeom,
01617                                    std::map<DTWireId, std::vector<PSimHit> > 
01618                                    simHitsPerWire,
01619                                    std::map<DTWireId, std::vector<type> > 
01620                                    recHitsPerWire,
01621                                    int step) {
01622 
01623   int nDt = 0;
01624   // Loop over cells with a muon SimHit
01625   for(std::map<DTWireId, std::vector<PSimHit> >::const_iterator wireAndSHits = 
01626         simHitsPerWire.begin();
01627       wireAndSHits != simHitsPerWire.end();
01628       wireAndSHits++) {
01629     DTWireId wireId = (*wireAndSHits).first;
01630     std::vector<PSimHit> simHitsInCell = (*wireAndSHits).second;
01631     
01632     // Get the layer
01633     const DTLayer* layer = dtGeom->layer(wireId);
01634     
01635     // Look for a mu hit in the cell
01636     const PSimHit* muSimHit = DTHitQualityUtils::findMuSimHit(simHitsInCell);
01637     if (muSimHit==0) {
01638       continue; // Skip this cell
01639     }
01640 
01641     // Find the distance of the simhit from the wire
01642     float simHitWireDist = simHitDistFromWire(layer, wireId, *muSimHit);
01643     // Skip simhits out of the cell
01644     if(simHitWireDist>2.1) {
01645       continue; // Skip this cell
01646     }
01647 
01648     // Look for RecHits in the same cell
01649     if(recHitsPerWire.find(wireId) == recHitsPerWire.end()) {
01650       continue; // No RecHit found in this cell
01651     } else {
01652 
01653       std::vector<type> recHits = recHitsPerWire[wireId];
01654          
01655       // Find the best RecHit
01656       const type* theBestRecHit = 
01657         findBestRecHit(layer, wireId, recHits, simHitWireDist);
01658  
01659       float recHitWireDist =  recHitDistFromWire(*theBestRecHit, layer);
01660       
01661       ++nDt;
01662 
01663       mehDtMuonRes->Fill(recHitWireDist-simHitWireDist);
01664       
01665     } // find rechits
01666   } // loop over simhits
01667 
01668   return nDt;
01669 }
01670 
01671 void 
01672 GlobalRecHitsAnalyzer::plotResolution(const PSimHit & simHit, 
01673                                       const CSCRecHit2D & recHit,
01674                                       const CSCLayer * layer, 
01675                                       int chamberType) {
01676   GlobalPoint simHitPos = layer->toGlobal(simHit.localPosition());
01677   GlobalPoint recHitPos = layer->toGlobal(recHit.localPosition());
01678   
01679   mehCSCResRDPhi->Fill(recHitPos.phi()-simHitPos.phi());
01680 }
01681