CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/SimG4CMS/HcalTestBeam/src/HcalTB02Analysis.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:     HcalTestBeam
00004 // Class  :     HcalTB02Analysis
00005 //
00006 // Implementation:
00007 //     Main analysis class for Hcal Test Beam 2002 Analysis
00008 //
00009 // Original Author:
00010 //         Created:  Sun May 21 10:14:34 CEST 2006
00011 // $Id: HcalTB02Analysis.cc,v 1.5 2009/05/24 15:31:15 fabiocos Exp $
00012 //
00013   
00014 // system include files
00015 #include <cmath>
00016 #include <iostream>
00017 #include <iomanip>
00018   
00019 // user include files
00020 #include "SimG4CMS/HcalTestBeam/interface/HcalTB02Analysis.h"
00021 
00022 // to retreive hits
00023 #include "SimG4CMS/Calo/interface/CaloG4Hit.h"
00024 #include "SimG4CMS/Calo/interface/CaloG4HitCollection.h"
00025 #include "SimG4CMS/HcalTestBeam/interface/HcalTB02HcalNumberingScheme.h"
00026 #include "SimG4CMS/HcalTestBeam/interface/HcalTB02XtalNumberingScheme.h"
00027 
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00030 #include "FWCore/ServiceRegistry/interface/Service.h"
00031 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00032 #include "CLHEP/Random/RandGaussQ.h"
00033 #include "FWCore/Utilities/interface/Exception.h"
00034 
00035 #include "G4SDManager.hh"
00036 #include "G4VProcess.hh"
00037 #include "G4HCofThisEvent.hh"
00038 #include "CLHEP/Units/GlobalSystemOfUnits.h"
00039 #include "CLHEP/Units/GlobalPhysicalConstants.h"
00040 
00041 //
00042 // constructors and destructor
00043 //
00044 
00045 HcalTB02Analysis::HcalTB02Analysis(const edm::ParameterSet &p): histo(0) {
00046 
00047   edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("HcalTB02Analysis");
00048   hcalOnly      = m_Anal.getUntrackedParameter<bool>("HcalClusterOnly",true);
00049   names         = m_Anal.getParameter<std::vector<std::string> >("Names");
00050   
00051   produces<HcalTB02HistoClass>();
00052 
00053   edm::LogInfo("HcalTBSim") << "HcalTB02Analysis:: Initialised as observer of "
00054                             << "BeginOfJob/BeginOfEvent/EndOfEvent with "
00055                             << "Parameter values:\n \thcalOnly = " << hcalOnly;
00056 
00057   histo  = new HcalTB02Histo(m_Anal);
00058 } 
00059    
00060 HcalTB02Analysis::~HcalTB02Analysis() {
00061 
00062   finish();
00063 
00064   if (histo)   {
00065     delete histo;
00066     histo  = 0;
00067   }
00068   edm::LogInfo("HcalTBSim") << "HcalTB02Analysis is deleting";
00069 }
00070 
00071 //
00072 // member functions
00073 //
00074 
00075 void HcalTB02Analysis::produce(edm::Event& e, const edm::EventSetup&) {
00076 
00077   std::auto_ptr<HcalTB02HistoClass> product(new HcalTB02HistoClass);
00078   fillEvent(*product);
00079   e.put(product);
00080 }
00081 
00082 void HcalTB02Analysis::update(const BeginOfEvent * evt) {
00083  
00084   edm::LogInfo("HcalTBSim") << "HcalTB02Analysis: =====> Begin of event = "
00085                             << (*evt) ()->GetEventID();
00086   clear();
00087 }
00088 
00089 void HcalTB02Analysis::update(const EndOfEvent * evt) {
00090 
00091   edm::Service<edm::RandomNumberGenerator> rng;
00092   if ( ! rng.isAvailable()) {
00093     throw cms::Exception("Configuration")
00094       << "HcalTB02Analysis requires the RandomNumberGeneratorService\n"
00095       << "which is not present in the configuration file. "
00096       << "You must add the service\n in the configuration file or "
00097       << "remove the modules that require it.";
00098   }
00099   CLHEP::RandGaussQ  randGauss(rng->getEngine());
00100 
00101   // Look for the Hit Collection
00102   LogDebug("HcalTBSim") << "HcalTB02Analysis::Fill event " 
00103                         << (*evt)()->GetEventID();
00104 
00105   // access to the G4 hit collections 
00106   G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
00107   int ihit = 0;    
00108  
00109   // HCAL
00110   std::string sd = names[0];
00111   int HCHCid = G4SDManager::GetSDMpointer()->GetCollectionID(sd);
00112   CaloG4HitCollection* theHCHC = (CaloG4HitCollection*) allHC->GetHC(HCHCid);
00113   HcalTB02HcalNumberingScheme *org = new HcalTB02HcalNumberingScheme();   
00114   LogDebug("HcalTBSim") << "HcalTB02Analysis :: Hit Collection for " << sd 
00115                         << " of ID " << HCHCid << " is obtained at " <<theHCHC;
00116 
00117   int nentries = 0;
00118   nentries = theHCHC->entries();
00119   if (nentries==0) return;
00120 
00121   int xentries = 0;  
00122   int XTALSid=0;
00123   CaloG4HitCollection* theXTHC=0;
00124 
00125   if (!hcalOnly) {
00126     // XTALS
00127     sd      = names[1];
00128     XTALSid = G4SDManager::GetSDMpointer()->GetCollectionID(sd);
00129     //    assert (XTALSid != 0);
00130     theXTHC = (CaloG4HitCollection*) allHC->GetHC(XTALSid);
00131     //    assert (theXTHC != 0);
00132     //HcalTB02XtalNumberingScheme *xorg = new HcalTB02XtalNumberingScheme();
00133     LogDebug("HcalTBSim") << "HcalTB02Analysis :: Hit Collection for " << sd
00134                           << " of ID " << XTALSid << " is obtained at " 
00135                           << theXTHC;
00136     xentries = theXTHC->entries();
00137     if (xentries==0) return;
00138   }
00139 
00140   LogDebug("HcalTBSim") << "HcalTB02Analysis :: There are " << nentries 
00141                         << " HCal hits, and" << xentries  << " xtal hits";
00142 
00143   float ETot=0., xETot=0.;
00144   //float maxE = 0.; 
00145   //int maxI=0, 
00146   int scintID=0, xtalID=0;
00147 
00148   // HCAL
00149 
00150   if (HCHCid >= 0 && theHCHC > 0) {
00151     for ( ihit = 0; ihit < nentries; ihit++) {
00152 
00153       CaloG4Hit* aHit = (*theHCHC)[ihit]; 
00154       scintID     = aHit->getUnitID();
00155       int layer   = org->getlayerID(scintID);
00156       float enEm  = aHit->getEM();
00157       float enhad = aHit->getHadr();
00158       
00159       if (layer==0) {
00160         enEm =enEm/2.;
00161         enhad=enhad/2.;
00162       }
00163 
00164       energyInScints[scintID]+= enEm + enhad;
00165       primaries[aHit->getTrackID()]+= enEm + enhad;
00166       float time = aHit->getTimeSliceID();
00167       if (time > maxTime) maxTime=(int)time;
00168       histo->fillAllTime(time);   
00169       
00170     }           
00171  
00172     //
00173     // Profile
00174     //
00175 
00176     float TowerEne[8][18], TowerEneCF[8][18], LayerEne[19], EnRing[100];
00177     for (int iphi=0 ; iphi<8; iphi++) {
00178       for (int jeta=0 ; jeta<18; jeta++) {
00179         TowerEne[iphi][jeta]=0.;
00180         TowerEneCF[iphi][jeta]=0.;
00181       }
00182     }
00183     
00184     for (int ilayer=0; ilayer<19; ilayer++) LayerEne[ilayer]=0.;
00185     for (int iring=0; iring<100; iring++) EnRing[iring]=0.;
00186     
00187     for (std::map<int,float>::iterator is = energyInScints.begin();
00188          is!= energyInScints.end(); is++) {
00189 
00190       ETot = (*is).second;
00191 
00192       int layer = org->getlayerID((*is).first);
00193 
00194       if ( (layer!=17) && (layer!=18) )  {
00195             
00196         float eta = org->getetaID((*is).first);
00197         float phi = org->getphiID((*is).first);
00198             
00199         SEnergy += ETot;
00200         TowerEne[(int)phi][(int)eta] += ETot;
00201 
00202         TowerEneCF[(int)phi][(int)eta] += ETot*(1.+ 0.1*randGauss.fire() );
00203         double dR=0.08727*std::sqrt( (eta-8.)*(eta-8.) + (phi-3.)*(phi-3.) );
00204         EnRing[(int)(dR/0.01)] += ETot;
00205       }
00206 
00207       LayerEne[layer] += ETot;
00208         
00209     }
00210     
00211     for (int ilayer=0 ; ilayer<19 ; ilayer++) {
00212       histo->fillProfile(ilayer,LayerEne[ilayer]/GeV);
00213     }
00214     
00215     for (int iring=0; iring<100; iring++)
00216       histo->fillTransProf(0.01*iring+0.005,EnRing[iring]/GeV);
00217     
00218     for (int iphi=0 ; iphi<8; iphi++) {
00219       for (int jeta=0 ; jeta<18; jeta++) {
00220         
00221         //SEnergyN += TowerEneCF[iphi][jeta] + 3.2*randGauss.fire(); // LHEP
00222         SEnergyN += TowerEneCF[iphi][jeta] + 3.*randGauss.fire(); // QGSP
00223 
00224         //double dR=0.08727*sqrt( (jeta-8.)*(jeta-8.)+(iphi-3.)*(iphi-3.) );
00225         //cout.testOut << " phi= " << iphi << " eta= " << jeta 
00226         //           << " TowerEne[iphi,jeta]= " << TowerEne[iphi][jeta] 
00227         //           << "dR= "  << dR << endl;
00228         
00229         //double Rand = 3.2*randGauss.fire(); // LHEP
00230         double Rand = 3.*randGauss.fire(); // QGSP
00231         
00232         if ( (iphi>=0) && (iphi<7) ) {
00233           if ( (jeta>=5) && (jeta<12) ) {
00234                 
00235             E7x7Matrix += TowerEne[iphi][jeta];
00236             E7x7MatrixN += TowerEneCF[iphi][jeta] + Rand;
00237                 
00238             if ( (iphi>=1) && (iphi<6) ) {
00239               if ( (jeta>=6) && (jeta<11) ) {
00240                         
00241                 E5x5Matrix += TowerEne[iphi][jeta];
00242                 E5x5MatrixN += TowerEneCF[iphi][jeta] + Rand;
00243                         
00244               }
00245             }
00246                 
00247           }
00248         }
00249         
00250       }
00251     }
00252     
00253     //
00254     // Find Primary info:
00255     //  
00256     int trackID = 0;
00257     G4PrimaryParticle* thePrim=0;
00258     G4int nvertex = (*evt)()->GetNumberOfPrimaryVertex();
00259     LogDebug("HcalTBSim") << "HcalTB02Analysis :: Event has " << nvertex 
00260                           << " vertex";
00261     if (nvertex==0)
00262       edm::LogWarning("HcalTBSim") << "HcalTB02Analysis:: End Of Event  "
00263                                    << "ERROR: no vertex";
00264 
00265     for (int i = 0 ; i<nvertex; i++) {
00266         
00267       G4PrimaryVertex* avertex = (*evt)()->GetPrimaryVertex(i);
00268       if (avertex == 0) {
00269         edm::LogWarning("HcalTBSim") << "HcalTB02Analysis:: End Of Event "
00270                                      << "ERROR: pointer to vertex = 0";
00271       } else {
00272         int npart = avertex->GetNumberOfParticle();
00273         LogDebug("HcalTBSim") << "HcalTB02Analysis::Vertex number :" << i 
00274                               << " with " << npart << " particles";
00275         if (thePrim==0) thePrim=avertex->GetPrimary(trackID);
00276       }
00277     }
00278     
00279     double px=0.,py=0.,pz=0.;
00280     
00281     if (thePrim != 0) {
00282       px = thePrim->GetPx();
00283       py = thePrim->GetPy();
00284       pz = thePrim->GetPz();
00285       pInit = std::sqrt(pow(px,2.)+pow(py,2.)+pow(pz,2.));
00286       if (pInit==0) {
00287         edm::LogWarning("HcalTBSim") << "HcalTB02Analysis:: End Of Event "
00288                                      << " ERROR: primary has p=0 ";
00289       } else {   
00290         float costheta = pz/pInit;
00291         float theta = acos(std::min(std::max(costheta,float(-1.)),float(1.)));
00292         eta = -log(tan(theta/2));
00293         if (px != 0) phi = atan(py/px);  
00294       }
00295       particleType      = thePrim->GetPDGcode();
00296     } else {
00297       LogDebug("HcalTBSim") << "HcalTB02Analysis:: End Of Event ERROR: could"
00298                             << " not find primary ";
00299     }
00300     
00301     CaloG4Hit* firstHit =(*theHCHC)[0];
00302     incidentEnergy = (firstHit->getIncidentEnergy()/GeV);
00303     
00304   }// number of Hits > 0
00305 
00306   if (!hcalOnly) {
00307 
00308     // XTALS
00309 
00310     if (XTALSid >= 0 && theXTHC > 0) {
00311       for (int xihit = 0; xihit < xentries; xihit++) {
00312 
00313         CaloG4Hit* xaHit = (*theXTHC)[xihit]; 
00314 
00315         float xenEm = xaHit->getEM();
00316         float xenhad = xaHit->getHadr();
00317         xtalID = xaHit->getUnitID();
00318           
00319         energyInCrystals[xtalID]+= xenEm + xenhad;
00320       }
00321       
00322       float xCrysEne[7][7];
00323       for (int irow=0 ; irow<7; irow++) {
00324         for (int jcol=0 ; jcol<7; jcol++) {
00325           xCrysEne[irow][jcol]=0.;
00326         }
00327       }
00328       
00329       for (std::map<int,float>::iterator is = energyInCrystals.begin();
00330            is!= energyInCrystals.end(); is++) {
00331         int xtalID = (*is).first;
00332         xETot = (*is).second;
00333             
00334         int irow = (int)(xtalID/100.);
00335         int jcol = (int)(xtalID-100.*irow);
00336         
00337         xSEnergy += xETot;
00338         xCrysEne[irow][jcol] = xETot;
00339             
00340         float dR=std::sqrt( 0.01619*0.01619*(jcol-3)*(jcol-3) + 
00341                             0.01606*0.01606*(irow-3)*(irow-3) );
00342         histo->fillTransProf(dR,xETot*1.05);
00343             
00344         if ( (irow>0) && (irow<6) ) {
00345           if ( (jcol>0) && (jcol<6) ) {
00346                     
00347             xE5x5Matrix += xCrysEne[irow][jcol];
00348             xE5x5MatrixN += xCrysEne[irow][jcol] + 108.5*randGauss.fire();
00349                     
00350             if ( (irow>1) && (irow<5) ) {
00351               if ( (jcol>1) && (jcol<5) ) {             
00352                 xE3x3Matrix += xCrysEne[irow][jcol];
00353                 xE3x3MatrixN += xCrysEne[irow][jcol] +108.5*randGauss.fire();
00354               }    
00355             }
00356           }
00357         }  
00358         
00359       }      
00360 
00361       if (!hcalOnly) {
00362         //      assert(theXTHC);
00363         if ( theXTHC != 0 ) {
00364           CaloG4Hit* xfirstHit =(*theXTHC)[0];
00365           xIncidentEnergy = xfirstHit->getIncidentEnergy()/GeV;
00366         }
00367       }
00368       
00369     }// number of Hits > 0
00370 
00371   }
00372 
00373   int iEvt = (*evt)()->GetEventID();
00374   if (iEvt < 10) 
00375     std::cout << " Event " << iEvt << std::endl;
00376   else if ((iEvt < 100) && (iEvt%10 == 0)) 
00377     std::cout << " Event " << iEvt << std::endl;
00378   else if ((iEvt < 1000) && (iEvt%100 == 0)) 
00379     std::cout << " Event " << iEvt << std::endl;
00380   else if ((iEvt < 10000) && (iEvt%1000 == 0)) 
00381     std::cout << " Event " << iEvt << std::endl;
00382 }
00383 
00384 void HcalTB02Analysis::fillEvent(HcalTB02HistoClass& product) {
00385 
00386   //Beam information
00387   product.set_Nprim(float(primaries.size()));
00388   product.set_partType(particleType);
00389   product.set_Einit(pInit/GeV);
00390   product.set_eta(eta);
00391   product.set_phi(phi);
00392   product.set_Eentry(incidentEnergy);
00393 
00394   //Calorimeter energy
00395   product.set_ETot(SEnergy/GeV );
00396   product.set_E7x7(E7x7Matrix/GeV );
00397   product.set_E5x5(E5x5Matrix/GeV );
00398   product.set_ETotN(SEnergyN/GeV);
00399   product.set_E7x7N(E7x7MatrixN/GeV );
00400   product.set_E5x5N(E5x5MatrixN/GeV );
00401   product.set_NUnit(float(energyInScints.size()));
00402   product.set_Ntimesli(float(maxTime));
00403 
00404   //crystal information
00405   product.set_xEentry(xIncidentEnergy);
00406   product.set_xNUnit(float(energyInCrystals.size()));
00407   product.set_xETot(xSEnergy/GeV);
00408   product.set_xETotN(xSEnergyN/GeV);
00409   product.set_xE5x5(xE5x5Matrix/GeV);
00410   product.set_xE3x3(xE3x3Matrix/GeV);
00411   product.set_xE5x5N(xE5x5MatrixN/GeV);
00412   product.set_xE3x3N(xE3x3MatrixN/GeV);
00413 }
00414 
00415 void HcalTB02Analysis::clear() {
00416 
00417   primaries.clear();
00418   particleType = 0;
00419   pInit = eta = phi = incidentEnergy = 0;
00420 
00421   SEnergy = E7x7Matrix = E5x5Matrix = SEnergyN = 0;
00422   E7x7MatrixN = E5x5MatrixN = 0;
00423   energyInScints.clear();
00424   maxTime = 0;
00425 
00426   xIncidentEnergy = 0;
00427   energyInCrystals.clear();
00428   xSEnergy = xSEnergyN = xE5x5Matrix = xE3x3Matrix = 0;
00429   xE5x5MatrixN = xE3x3MatrixN = 0;
00430 }
00431 
00432 void HcalTB02Analysis::finish() {
00433 
00434   /*
00435   //Profile 
00436   std::ofstream   oFile;
00437   oFile.open("profile.dat");
00438   float st[19] = {0.8,
00439                   0.4,  0.4,  0.4,  0.4,  0.4,   
00440                   0.4,  0.4,  0.4,  0.4,  0.4,
00441                   0.4,  0.4,  0.4,  0.4,  0.4, 
00442                   0.8,  1.0,  1.0};
00443                  
00444   //cm of material (brass) in front of scintillator layer i:
00445  
00446   float w[19] = {7.45,                         //iron !
00447                  6.00, 6.00, 6.00, 6.00, 6.00, //brass
00448                  6.00, 6.00, 6.00, 6.00, 6.60, //brass
00449                  6.60, 6.60, 6.60, 6.60, 6.60, //brass
00450                  8.90, 20.65, 19.5};            //brass,iron !
00451 
00452   for (int ilayer = 0; ilayer<19; ilayer++) {
00453 
00454     // Histogram mean and sigma calculated from the ROOT histos
00455     edm::LogInfo("HcalTBSim") << "Layer number: " << ilayer << " Mean = " 
00456                               << histo->getMean(ilayer) << " sigma = "   
00457                               << histo->getRMS(ilayer) << " LThick= "   
00458                               << w[ilayer] << " SThick= "   << st[ilayer];
00459       
00460     oFile << ilayer << " "  << histo->getMean(ilayer) << " " 
00461           << histo->getRMS(ilayer)  << " " << w[ilayer] << " " << st[ilayer] 
00462           << std::endl;
00463 
00464   } 
00465   oFile.close(); 
00466   */
00467 }