CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/Validation/HcalHits/src/HcalSimHitsValidation.cc

Go to the documentation of this file.
00001 #include "Validation/HcalHits/interface/HcalSimHitsValidation.h"
00002 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00003 
00004 
00005 HcalSimHitsValidation::HcalSimHitsValidation(edm::ParameterSet const& conf) {
00006   // DQM ROOT output
00007   outputFile_ = conf.getUntrackedParameter<std::string>("outputFile", "myfile.root");
00008   
00009   if ( outputFile_.size() != 0 ) {    edm::LogInfo("OutputInfo") << " Hcal SimHit Task histograms will be saved to '" << outputFile_.c_str() << "'";
00010   } else {
00011     edm::LogInfo("OutputInfo") << " Hcal SimHit Task histograms will NOT be saved";
00012   }
00013   
00014   nevtot = 0;
00015   
00016   dbe_ = 0;
00017   // get hold of back-end interface
00018   dbe_ = edm::Service<DQMStore>().operator->();
00019    
00020   Char_t histo[200];
00021 
00022   if ( dbe_ ) {
00023     dbe_->setCurrentFolder("HcalSimHitsV/HcalSimHitTask");
00024 
00025     // General counters
00026     sprintf  (histo, "N_HB" );
00027     Nhb = dbe_->book1D(histo, histo, 2600,0.,2600.);
00028     sprintf  (histo, "N_HE" );
00029     Nhe = dbe_->book1D(histo, histo, 2600,0.,2600.);
00030     sprintf  (histo, "N_HO" );
00031     Nho = dbe_->book1D(histo, histo, 2200,0.,2200.);
00032     sprintf  (histo, "N_HF" );
00033     Nhf = dbe_->book1D(histo, histo, 1800,0., 1800.);
00034 
00035     //Mean energy vs iEta TProfiles
00036     sprintf  (histo, "emean_vs_ieta_HB1" );
00037     emean_vs_ieta_HB1 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s");
00038     sprintf  (histo, "emean_vs_ieta_HB2" );
00039     emean_vs_ieta_HB2 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s");
00040     sprintf  (histo, "emean_vs_ieta_HE1" );
00041     emean_vs_ieta_HE1 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10. ,2000., "s" );
00042     sprintf  (histo, "emean_vs_ieta_HE2" );
00043     emean_vs_ieta_HE2 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s");
00044     sprintf  (histo, "emean_vs_ieta_HE3" );
00045     emean_vs_ieta_HE3 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s" );
00046     sprintf  (histo, "emean_vs_ieta_HO" );
00047     emean_vs_ieta_HO = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s" );
00048     sprintf  (histo, "emean_vs_ieta_HF1" );
00049     emean_vs_ieta_HF1 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s" );
00050     sprintf  (histo, "emean_vs_ieta_HF2" );
00051     emean_vs_ieta_HF2 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s" );
00052 
00053     //Occupancy vs. iEta TH1Fs
00054     sprintf  (histo, "occupancy_vs_ieta_HB1" );
00055     occupancy_vs_ieta_HB1 = dbe_->book1D(histo, histo, 82, -41., 41.);
00056     sprintf  (histo, "occupancy_vs_ieta_HB2" );
00057     occupancy_vs_ieta_HB2 = dbe_->book1D(histo, histo, 82, -41., 41.);
00058     sprintf  (histo, "occupancy_vs_ieta_HE1" );
00059     occupancy_vs_ieta_HE1 = dbe_->book1D(histo, histo, 82, -41., 41.);
00060     sprintf  (histo, "occupancy_vs_ieta_HE2" );
00061     occupancy_vs_ieta_HE2 = dbe_->book1D(histo, histo, 82, -41., 41.);
00062     sprintf  (histo, "occupancy_vs_ieta_HE3" );
00063     occupancy_vs_ieta_HE3 = dbe_->book1D(histo, histo, 82, -41., 41.);
00064     sprintf  (histo, "occupancy_vs_ieta_HO" );
00065     occupancy_vs_ieta_HO = dbe_->book1D(histo, histo, 82, -41., 41.);
00066     sprintf  (histo, "occupancy_vs_ieta_HF1" );
00067     occupancy_vs_ieta_HF1 = dbe_->book1D(histo, histo, 82, -41., 41.);
00068     sprintf  (histo, "occupancy_vs_ieta_HF2" );
00069     occupancy_vs_ieta_HF2 = dbe_->book1D(histo, histo, 82, -41., 41.);
00070 
00071     //Energy spectra
00072     sprintf (histo, "HcalSimHitTask_energy_of_simhits_HB" ) ;
00073     meSimHitsEnergyHB = dbe_->book1D(histo, histo, 510 , -0.1 , 5.); 
00074 
00075     sprintf (histo, "HcalSimHitTask_energy_of_simhits_HE" ) ;
00076     meSimHitsEnergyHE = dbe_->book1D(histo, histo, 510, -0.02, 2.); 
00077 
00078     sprintf (histo, "HcalSimHitTask_energy_of_simhits_HO" ) ;
00079     meSimHitsEnergyHO = dbe_->book1D(histo, histo, 510 , -0.1 , 5.); 
00080     
00081     sprintf (histo, "HcalSimHitTask_energy_of_simhits_HF" ) ;
00082     meSimHitsEnergyHF = dbe_->book1D(histo, histo, 1010 , -5. , 500.); 
00083 
00084     //Energy in Cone
00085     sprintf (histo, "HcalSimHitTask_En_simhits_cone_profile_vs_ieta_all_depths");
00086     meEnConeEtaProfile = dbe_->bookProfile(histo, histo, 82, -41., 41., 210, -10., 200.);  
00087     
00088     sprintf (histo, "HcalSimHitTask_En_simhits_cone_profile_vs_ieta_all_depths_E");
00089     meEnConeEtaProfile_E = dbe_->bookProfile(histo, histo, 82, -41., 41., 210, -10., 200.);  
00090       
00091     sprintf (histo, "HcalSimHitTask_En_simhits_cone_profile_vs_ieta_all_depths_EH");
00092     meEnConeEtaProfile_EH = dbe_->bookProfile(histo, histo, 82, -41., 41., 210, -10., 200.);  
00093     
00094    }  //end-of if(_dbe) 
00095 
00096 }
00097 
00098 
00099 HcalSimHitsValidation::~HcalSimHitsValidation() { }
00100 
00101 void HcalSimHitsValidation::endJob() { 
00102   //before check that histos are there....
00103 
00104   // check if ME still there (and not killed by MEtoEDM for memory saving)
00105   if( dbe_ )
00106     {
00107       // check existence of first histo in the list
00108       if (! dbe_->get("HcalSimHitsV/HcalSimHitTask/N_HB")) return;
00109     }
00110   else
00111     return;
00112   
00113   //======================================
00114 
00115   for (int i = 1; i <= occupancy_vs_ieta_HB1->getNbinsX(); i++){
00116 
00117     int ieta = i - 42;        // -41 -1, 0 40 
00118     if(ieta >=0 ) ieta +=1;   // -41 -1, 1 41  - to make it detector-like
00119 
00120     float phi_factor;
00121 
00122     if      (fabs(ieta) <= 20) phi_factor = 72.;
00123     else if (fabs(ieta) <  40) phi_factor = 36.;
00124     else                       phi_factor = 18.;
00125     
00126     float cnorm;
00127 
00128     //Occupancy vs. iEta TH1Fs
00129     cnorm = occupancy_vs_ieta_HB1->getBinContent(i) / (phi_factor * nevtot); 
00130     occupancy_vs_ieta_HB1->setBinContent(i, cnorm);
00131     cnorm = occupancy_vs_ieta_HB2->getBinContent(i) / (phi_factor * nevtot); 
00132     occupancy_vs_ieta_HB2->setBinContent(i, cnorm);
00133 
00134     cnorm = occupancy_vs_ieta_HE1->getBinContent(i) / (phi_factor * nevtot); 
00135     occupancy_vs_ieta_HE1->setBinContent(i, cnorm);
00136     cnorm = occupancy_vs_ieta_HE2->getBinContent(i) / (phi_factor * nevtot); 
00137     occupancy_vs_ieta_HE2->setBinContent(i, cnorm);
00138     cnorm = occupancy_vs_ieta_HE3->getBinContent(i) / (phi_factor * nevtot); 
00139     occupancy_vs_ieta_HE3->setBinContent(i, cnorm);
00140 
00141     cnorm = occupancy_vs_ieta_HO->getBinContent(i) / (phi_factor * nevtot); 
00142     occupancy_vs_ieta_HO->setBinContent(i, cnorm);
00143 
00144     cnorm = occupancy_vs_ieta_HF1->getBinContent(i) / (phi_factor * nevtot); 
00145     occupancy_vs_ieta_HF1->setBinContent(i, cnorm);
00146     cnorm = occupancy_vs_ieta_HF2->getBinContent(i) / (phi_factor * nevtot); 
00147     occupancy_vs_ieta_HF2->setBinContent(i, cnorm);
00148   }
00149 
00150 
00151   if ( outputFile_.size() != 0 && dbe_ ) dbe_->save(outputFile_);
00152 }
00153 
00154 
00155 void HcalSimHitsValidation::beginJob(){ }
00156 
00157 void HcalSimHitsValidation::analyze(edm::Event const& ev, edm::EventSetup const& c) {
00158 
00159   using namespace edm;
00160   using namespace std;
00161 
00162   //===========================================================================
00163   // Getting SimHits
00164   //===========================================================================
00165 
00166   double phi_MC = -999.;  // phi of initial particle from HepMC
00167   double eta_MC = -999.;  // eta of initial particle from HepMC
00168 
00169   edm::Handle<edm::HepMCProduct> evtMC;
00170   //  ev.getByLabel("VtxSmeared",evtMC);
00171   ev.getByLabel("generator",evtMC);  // generator in late 310_preX
00172   if (!evtMC.isValid()) {
00173     std::cout << "no HepMCProduct found" << std::endl;    
00174   }
00175 
00176   // MC particle with highest pt is taken as a direction reference  
00177   double maxPt = -99999.;
00178   int npart    = 0;
00179 
00180   const HepMC::GenEvent * myGenEvent = evtMC->GetEvent();
00181   for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
00182         p != myGenEvent->particles_end(); ++p ) {
00183     double phip = (*p)->momentum().phi();
00184     double etap = (*p)->momentum().eta();
00185     double pt   = (*p)->momentum().perp();
00186     if(pt > maxPt) {npart++; maxPt = pt; phi_MC = phip; eta_MC = etap; }
00187   }
00188 
00189   double partR   = 0.3;
00190 
00191 
00192   //Hcal SimHits
00193 
00194   //Approximate calibration constants
00195   const float calib_HB = 120.;
00196   const float calib_HE = 190.;
00197   const float calib_HF1 = 1.0/0.383;
00198   const float calib_HF2 = 1.0/0.368;
00199   
00200   edm::Handle<PCaloHitContainer> hcalHits;
00201   ev.getByLabel("g4SimHits","HcalHits",hcalHits);
00202   const PCaloHitContainer * SimHitResult = hcalHits.product () ;
00203     
00204   float eta_diff;
00205   float etaMax  = 9999;
00206   int ietaMax = 0;
00207 
00208   double HcalCone = 0;
00209 
00210   c.get<CaloGeometryRecord>().get (geometry);
00211     
00212   for (std::vector<PCaloHit>::const_iterator SimHits = SimHitResult->begin () ; SimHits != SimHitResult->end(); ++SimHits) {
00213     HcalDetId cell(SimHits->id());
00214     const CaloCellGeometry* cellGeometry = geometry->getSubdetectorGeometry (cell)->getGeometry (cell);
00215     double etaS = cellGeometry->getPosition().eta () ;
00216     double phiS = cellGeometry->getPosition().phi () ;
00217     double en   = SimHits->energy();    
00218     
00219     int sub     = cell.subdet();
00220     int depth   = cell.depth();
00221     double ieta = cell.ieta();
00222 
00223     //Energy in Cone 
00224     double r  = dR(eta_MC, phi_MC, etaS, phiS);
00225     
00226     if (r < partR){      
00227       eta_diff = fabs(eta_MC - etaS);
00228       if(eta_diff < etaMax) {
00229         etaMax  = eta_diff; 
00230         ietaMax = cell.ieta(); 
00231       }
00232       //Approximation of calibration
00233       if      (sub == 1)               HcalCone += en*calib_HB;
00234       else if (sub == 2)               HcalCone += en*calib_HE;
00235       else if (sub == 4 && depth == 1) HcalCone += en*calib_HF1;
00236       else if (sub == 4 && depth == 2) HcalCone += en*calib_HF2;
00237     }
00238     
00239     //Account for lack of ieta = 0
00240     if (ieta > 0) ieta--;
00241 
00242     //HB
00243     if (sub == 1){
00244       meSimHitsEnergyHB->Fill(en);
00245       if (depth == 1){
00246         emean_vs_ieta_HB1->Fill(double(ieta), en);
00247         occupancy_vs_ieta_HB1->Fill(double(ieta));
00248       }
00249       if (depth == 2){
00250         emean_vs_ieta_HB2->Fill(double(ieta), en);
00251         occupancy_vs_ieta_HB2->Fill(double(ieta));
00252       }
00253     }
00254     //HE
00255     if (sub == 2){
00256       meSimHitsEnergyHE->Fill(en);
00257       if (depth == 1){
00258         emean_vs_ieta_HE1->Fill(double(ieta), en);
00259         occupancy_vs_ieta_HE1->Fill(double(ieta));
00260       }
00261       if (depth == 2){
00262         emean_vs_ieta_HE2->Fill(double(ieta), en);
00263         occupancy_vs_ieta_HE2->Fill(double(ieta));
00264       }
00265       if (depth == 3){
00266         emean_vs_ieta_HE3->Fill(double(ieta), en);
00267         occupancy_vs_ieta_HE3->Fill(double(ieta));
00268       }
00269     }
00270     //HO
00271     if (sub == 3){
00272       meSimHitsEnergyHO->Fill(en);
00273       emean_vs_ieta_HO->Fill(double(ieta), en);
00274       occupancy_vs_ieta_HO->Fill(double(ieta));
00275     }
00276     //HF
00277     if (sub == 4){
00278       meSimHitsEnergyHF->Fill(en);
00279       if (depth == 1){
00280         emean_vs_ieta_HF1->Fill(double(ieta), en);
00281         occupancy_vs_ieta_HF1->Fill(double(ieta));
00282       }
00283       if (depth == 2){
00284         emean_vs_ieta_HF2->Fill(double(ieta), en);
00285         occupancy_vs_ieta_HF2->Fill(double(ieta));
00286       }
00287     }     
00288   }
00289 
00290   //Ecal EB SimHits
00291   edm::Handle<PCaloHitContainer> ecalEBHits;
00292   ev.getByLabel("g4SimHits","EcalHitsEB",ecalEBHits);
00293   const PCaloHitContainer * SimHitResultEB = ecalEBHits.product () ;
00294 
00295   double EcalCone = 0;
00296 
00297   for (std::vector<PCaloHit>::const_iterator SimHits = SimHitResultEB->begin () ; SimHits != SimHitResultEB->end(); ++SimHits) {
00298 
00299     EBDetId EBid = EBDetId(SimHits->id());
00300 
00301     const CaloCellGeometry* cellGeometry = geometry->getSubdetectorGeometry (EBid)->getGeometry (EBid) ;
00302     double etaS = cellGeometry->getPosition().eta () ;
00303     double phiS = cellGeometry->getPosition().phi () ;
00304     double en   = SimHits->energy();    
00305   
00306     double r  = dR(eta_MC, phi_MC, etaS, phiS);
00307     
00308     if (r < partR) EcalCone += en;   
00309   }
00310 
00311   //Ecal EE SimHits
00312   edm::Handle<PCaloHitContainer> ecalEEHits;
00313   ev.getByLabel("g4SimHits","EcalHitsEE",ecalEEHits);
00314   const PCaloHitContainer * SimHitResultEE = ecalEEHits.product () ;
00315 
00316   for (std::vector<PCaloHit>::const_iterator SimHits = SimHitResultEE->begin () ; SimHits != SimHitResultEE->end(); ++SimHits) {
00317 
00318     EEDetId EEid = EEDetId(SimHits->id());
00319 
00320     const CaloCellGeometry* cellGeometry = geometry->getSubdetectorGeometry (EEid)->getGeometry (EEid) ;
00321     double etaS = cellGeometry->getPosition().eta () ;
00322     double phiS = cellGeometry->getPosition().phi () ;
00323     double en   = SimHits->energy();    
00324     
00325     double r  = dR(eta_MC, phi_MC, etaS, phiS);
00326     
00327     if (r < partR) EcalCone += en;   
00328   }
00329 
00330   if (ietaMax != 0){            //If ietaMax == 0, there were no good HCAL SimHits 
00331     if (ietaMax > 0) ietaMax--; //Account for lack of ieta = 0
00332     
00333     meEnConeEtaProfile       ->Fill(double(ietaMax), HcalCone);    
00334     meEnConeEtaProfile_E     ->Fill(double(ietaMax), EcalCone);
00335     meEnConeEtaProfile_EH    ->Fill(double(ietaMax), HcalCone+EcalCone); 
00336   }
00337   
00338   nevtot++;
00339 }
00340 
00341 
00342 double HcalSimHitsValidation::dR(double eta1, double phi1, double eta2, double phi2) { 
00343   double PI = 3.1415926535898;
00344   double deltaphi= phi1 - phi2;
00345   if( phi2 > phi1 ) { deltaphi= phi2 - phi1;}
00346   if(deltaphi > PI) { deltaphi = 2.*PI - deltaphi;}
00347   double deltaeta = eta2 - eta1;
00348   double tmp = sqrt(deltaeta* deltaeta + deltaphi*deltaphi);
00349   return tmp;
00350 }
00351 
00352 double HcalSimHitsValidation::phi12(double phi1, double en1, double phi2, double en2) {
00353   // weighted mean value of phi1 and phi2
00354   
00355   double tmp;
00356   double PI = 3.1415926535898;
00357   double a1 = phi1; double a2  = phi2;
00358 
00359   if( a1 > 0.5*PI  && a2 < 0.) a2 += 2*PI; 
00360   if( a2 > 0.5*PI  && a1 < 0.) a1 += 2*PI; 
00361   tmp = (a1 * en1 + a2 * en2)/(en1 + en2);
00362   if(tmp > PI) tmp -= 2.*PI; 
00363  
00364   return tmp;
00365 
00366 }
00367 
00368 double HcalSimHitsValidation::dPhiWsign(double phi1, double phi2) {
00369   // clockwise      phi2 w.r.t phi1 means "+" phi distance
00370   // anti-clockwise phi2 w.r.t phi1 means "-" phi distance 
00371 
00372   double PI = 3.1415926535898;
00373   double a1 = phi1; double a2  = phi2;
00374   double tmp =  a2 - a1;
00375   if( a1*a2 < 0.) {
00376     if(a1 > 0.5 * PI)  tmp += 2.*PI;
00377     if(a2 > 0.5 * PI)  tmp -= 2.*PI;
00378   }
00379   return tmp;
00380 
00381 }
00382 
00383 
00384 
00385 DEFINE_FWK_MODULE(HcalSimHitsValidation);
00386