CMS 3D CMS Logo

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

Generated on Tue Jun 9 17:49:22 2009 for CMSSW by  doxygen 1.5.4