CMS 3D CMS Logo

DQMHLTSourcePi0 Class Reference

#include <HLTriggerOffline/special/src/DQMHLTSourcePi0.h>

Inheritance diagram for DQMHLTSourcePi0:

edm::EDAnalyzer

List of all members.

Public Member Functions

 DQMHLTSourcePi0 (const edm::ParameterSet &)
 ~DQMHLTSourcePi0 ()

Protected Member Functions

void analyze (const edm::Event &e, const edm::EventSetup &c)
void beginJob (const edm::EventSetup &c)
void beginLuminosityBlock (const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &context)
void beginRun (const edm::Run &r, const edm::EventSetup &c)
void endJob ()
void endLuminosityBlock (const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &c)
void endRun (const edm::Run &r, const edm::EventSetup &c)

Private Attributes

int clusEtaSize_
int clusPhiSize_
double clusSeedThr_
DQMStoredbe_
int eventCounter_
std::string fileName_
 Output file name if required.
std::string folderName_
 DQM folder name.
int gammaCandEtaSize_
int gammaCandPhiSize_
MonitorElementhEventEnergyEB_
 Distribution of total event energy.
MonitorElementhEventEnergyEE_
 Distribution of total event energy.
MonitorElementhiEtaDistrEB_
 Distribution of rechits in iEta.
MonitorElementhiPhiDistrEB_
 Distribution of rechits in iPhi.
MonitorElementhIsoPi0EB_
 Pi0 Iso.
MonitorElementhMeanRecHitEnergyEB_
 Distribution of Mean energy per rechit.
MonitorElementhMeanRecHitEnergyEE_
 Distribution of Mean energy per rechit.
MonitorElementhMinvPi0EB_
 Pi0 invariant mass in EB.
MonitorElementhNRecHitsEB_
 Distribution of number of RecHits.
MonitorElementhNRecHitsEE_
 Distribution of number of RecHits.
MonitorElementhPt1Pi0EB_
 Pt of the 1st most energetic Pi0 photon in EB.
MonitorElementhPt2Pi0EB_
 Pt of the 2nd most energetic Pi0 photon in EB.
MonitorElementhPtPi0EB_
 Pi0 Pt in EB.
MonitorElementhRechitEnergyEB_
 Energy Distribution of rechits.
MonitorElementhRechitEnergyEE_
 Energy Distribution of rechits.
MonitorElementhS4S91EB_
 S4S9 of the 1st most energetic pi0 photon.
MonitorElementhS4S92EB_
 S4S9 of the 2nd most energetic pi0 photon.
bool isMonEB_
 which subdet will be monitored
bool isMonEE_
bool ParameterLogWeighted_
double ParameterT0_barl_
double ParameterW0_
double ParameterX0_
unsigned int prescaleFactor_
 Monitor every prescaleFactor_ events.
edm::InputTag productMonitoredEB_
 object to monitor
edm::InputTag productMonitoredEE_
 object to monitor
bool saveToFile_
 Write to file.
double seleMinvMaxPi0_
double seleMinvMinPi0_
int seleNRHMax_
double selePi0BeltDeta_
double selePi0BeltDR_
double selePi0Iso_
double selePtGammaOne_
double selePtGammaTwo_
double selePtPi0_
double seleS4S9GammaOne_
double seleS4S9GammaTwo_
double seleXtalMinEnergy_


Detailed Description

Definition at line 40 of file DQMHLTSourcePi0.h.


Constructor & Destructor Documentation

DQMHLTSourcePi0::DQMHLTSourcePi0 ( const edm::ParameterSet ps  ) 

Definition at line 52 of file DQMHLTSourcePi0.cc.

References clusEtaSize_, clusPhiSize_, clusSeedThr_, dbe_, fileName_, folderName_, gammaCandEtaSize_, gammaCandPhiSize_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), isMonEB_, isMonEE_, ParameterLogWeighted_, ParameterT0_barl_, ParameterW0_, ParameterX0_, prescaleFactor_, productMonitoredEB_, productMonitoredEE_, saveToFile_, seleMinvMaxPi0_, seleMinvMinPi0_, seleNRHMax_, selePi0BeltDeta_, selePi0BeltDR_, selePi0Iso_, selePtGammaOne_, selePtGammaTwo_, selePtPi0_, seleS4S9GammaOne_, seleS4S9GammaTwo_, and seleXtalMinEnergy_.

00052                                                             :
00053 eventCounter_(0)
00054 {
00055   dbe_ = Service<DQMStore>().operator->();
00056   folderName_ = ps.getUntrackedParameter<string>("FolderName","HLT/AlCaEcalPi0");
00057   prescaleFactor_ = ps.getUntrackedParameter<int>("prescaleFactor",1);
00058   productMonitoredEB_= ps.getUntrackedParameter<edm::InputTag>("AlCaStreamEBTag");
00059   productMonitoredEE_= ps.getUntrackedParameter<edm::InputTag>("AlCaStreamEETag");
00060   isMonEB_ = ps.getUntrackedParameter<bool>("isMonEB",false);
00061   isMonEE_ = ps.getUntrackedParameter<bool>("isMonEE",false);
00062 
00063   saveToFile_=ps.getUntrackedParameter<bool>("SaveToFile",false);
00064   fileName_=  ps.getUntrackedParameter<string>("FileName","MonitorAlCaEcalPi0.root");
00065 
00066   gammaCandEtaSize_ = ps.getParameter<int> ("gammaCandEtaSize");
00067   gammaCandPhiSize_ = ps.getParameter<int> ("gammaCandPhiSize");
00068   if ( gammaCandPhiSize_ % 2 == 0 ||  gammaCandEtaSize_ % 2 == 0)
00069     edm::LogError("AlCaPi0RecHitsProducerError") << "Size of eta/phi for sliding window should be odd numbers";
00070 
00071   clusSeedThr_ = ps.getParameter<double> ("clusSeedThr");
00072   clusEtaSize_ = ps.getParameter<int> ("clusEtaSize");
00073   clusPhiSize_ = ps.getParameter<int> ("clusPhiSize");
00074   if ( clusPhiSize_ % 2 == 0 ||  clusEtaSize_ % 2 == 0)
00075     edm::LogError("AlCaPi0RecHitsProducerError") << "Size of eta/phi for simple clustering should be odd numbers";
00076 
00077   selePtGammaOne_ = ps.getParameter<double> ("selePtGammaOne");  
00078   selePtGammaTwo_ = ps.getParameter<double> ("selePtGammaTwo");  
00079   selePtPi0_ = ps.getParameter<double> ("selePtPi0");  
00080   seleMinvMaxPi0_ = ps.getParameter<double> ("seleMinvMaxPi0");  
00081   seleMinvMinPi0_ = ps.getParameter<double> ("seleMinvMinPi0");  
00082   seleXtalMinEnergy_ = ps.getParameter<double> ("seleXtalMinEnergy");
00083   seleNRHMax_ = ps.getParameter<int> ("seleNRHMax");
00084   //New Selection
00085   seleS4S9GammaOne_ = ps.getParameter<double> ("seleS4S9GammaOne");  
00086   seleS4S9GammaTwo_ = ps.getParameter<double> ("seleS4S9GammaTwo");  
00087   selePi0Iso_ = ps.getParameter<double> ("selePi0Iso");  
00088   selePi0BeltDR_ = ps.getParameter<double> ("selePi0BeltDR");  
00089   selePi0BeltDeta_ = ps.getParameter<double> ("selePi0BeltDeta");  
00090 
00091   ParameterLogWeighted_ = ps.getParameter<bool> ("ParameterLogWeighted");
00092   ParameterX0_ = ps.getParameter<double> ("ParameterX0");
00093   ParameterT0_barl_ = ps.getParameter<double> ("ParameterT0_barl");
00094   ParameterW0_ = ps.getParameter<double> ("ParameterW0");
00095 
00096 
00097 
00098 }

DQMHLTSourcePi0::~DQMHLTSourcePi0 (  ) 

Definition at line 101 of file DQMHLTSourcePi0.cc.

00102 {}


Member Function Documentation

void DQMHLTSourcePi0::analyze ( const edm::Event e,
const edm::EventSetup c 
) [protected, virtual]

Implements edm::EDAnalyzer.

Definition at line 192 of file DQMHLTSourcePi0.cc.

References PositionCalc::Calculate_Location(), clusEtaSize_, clusPhiSize_, clusSeedThr_, funct::cos(), GenMuonPlsPt100GeV_cfg::cout, DetId::Ecal, EcalBarrel, EcalPreshower, lat::endl(), relval_parameters_module::energy, eventCounter_, funct::exp(), MonitorElement::Fill(), edm::EventSetup::get(), edm::Event::getByLabel(), CaloSubdetectorTopology::getWindow(), hEventEnergyEB_, hEventEnergyEE_, hiEtaDistrEB_, hiPhiDistrEB_, hIsoPi0EB_, hMeanRecHitEnergyEB_, hMeanRecHitEnergyEE_, hMinvPi0EB_, hNRecHitsEB_, hNRecHitsEE_, hPt1Pi0EB_, hPt2Pi0EB_, hPtPi0EB_, hRechitEnergyEB_, hRechitEnergyEE_, hS4S91EB_, hS4S92EB_, i, id, EBDetId::ieta(), EBDetId::iphi(), isMonEB_, isMonEE_, edm::Handle< T >::isValid(), j, k, ParameterLogWeighted_, ParameterT0_barl_, ParameterW0_, ParameterX0_, prescaleFactor_, edm::Handle< T >::product(), productMonitoredEB_, productMonitoredEE_, seleMinvMaxPi0_, seleMinvMinPi0_, selePi0BeltDeta_, selePi0BeltDR_, selePi0Iso_, selePtGammaOne_, selePtGammaTwo_, selePtPi0_, seleS4S9GammaOne_, seleS4S9GammaTwo_, seleXtalMinEnergy_, funct::sin(), python::multivaluedict::sort(), funct::sqrt(), and muonGeometry::TVector3().

00193                                                          {  
00194  
00195   if (eventCounter_% prescaleFactor_ ) return; 
00196   eventCounter_++;
00197 
00198   edm::ESHandle<CaloTopology> theCaloTopology;
00199   iSetup.get<CaloTopologyRecord>().get(theCaloTopology);
00200   //  bool accept=false;
00201 
00202   std::map<DetId, EcalRecHit>  recHitsEB_map;
00203   //  recHitsEB_map= new std::map<DetId, EcalRecHit>();
00204 
00205   std::vector<EcalRecHit> seeds;
00206   seeds.clear();
00207 
00208   vector<EBDetId> usedXtals;
00209   usedXtals.clear();
00210 
00211 
00212 
00213 
00214 
00215 
00216     
00217   edm::Handle<EcalRecHitCollection> rhEB;
00218   edm::Handle<EcalRecHitCollection> rhEE;
00219  
00220   if(isMonEB_) iEvent.getByLabel(productMonitoredEB_, rhEB); 
00221   if(isMonEE_) iEvent.getByLabel(productMonitoredEE_, rhEE);
00222 
00223   EcalRecHitCollection::const_iterator itb;
00224 
00225   // fill EB histos
00226   if(isMonEB_){
00227     if (rhEB.isValid()){
00228       float etot =0;
00229       for(itb=rhEB->begin(); itb!=rhEB->end(); ++itb){
00230         
00231         EBDetId id(itb->id());
00232         double energy = itb->energy();
00233         if (energy > seleXtalMinEnergy_) {
00234           std::pair<DetId, EcalRecHit> map_entry(itb->id(), *itb);
00235           recHitsEB_map.insert(map_entry);
00236         }
00237         if (energy > clusSeedThr_) seeds.push_back(*itb);
00238 
00239         hiPhiDistrEB_->Fill(id.iphi());
00240         hiEtaDistrEB_->Fill(id.ieta());
00241         hRechitEnergyEB_->Fill(itb->energy());
00242         
00243         etot+= itb->energy();    
00244       } // Eb rechits
00245       
00246       hNRecHitsEB_->Fill(rhEB->size());
00247       hMeanRecHitEnergyEB_->Fill(etot/rhEB->size());
00248       hEventEnergyEB_->Fill(etot);
00249 
00250       // Pi0 maker
00251 
00252       //cout<< " RH coll size: "<<rhEB->size()<<endl;
00253       //cout<< " Pi0 seeds: "<<seeds.size()<<endl;
00254 
00255       // Initialize the Position Calc
00256       const CaloSubdetectorGeometry *geometry_p;    
00257       const CaloSubdetectorTopology *topology_p;
00258       const CaloSubdetectorGeometry *geometryES_p;
00259       const EcalRecHitCollection *hitCollection_p = rhEB.product();
00260 
00261       edm::ESHandle<CaloGeometry> geoHandle;
00262       iSetup.get<CaloGeometryRecord>().get(geoHandle);     
00263       geometry_p = geoHandle->getSubdetectorGeometry(DetId::Ecal,EcalBarrel);
00264       geometryES_p = geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalPreshower);
00265 
00266       std::map<std::string,double> providedParameters;  
00267       providedParameters.insert(std::make_pair("LogWeighted",ParameterLogWeighted_));
00268       providedParameters.insert(std::make_pair("X0",ParameterX0_));
00269       providedParameters.insert(std::make_pair("T0_barl",ParameterT0_barl_));
00270       providedParameters.insert(std::make_pair("W0",ParameterW0_));
00271 
00272       PositionCalc posCalculator_ = PositionCalc(providedParameters);
00273 
00274 
00275       static const int MAXCLUS = 2000;
00276       int nClus;
00277       vector<float> eClus;
00278       vector<float> etClus;
00279       vector<float> etaClus;
00280       vector<float> phiClus;
00281       vector<EBDetId> max_hit;
00282       vector< vector<EcalRecHit> > RecHitsCluster;
00283       vector<float> s4s9Clus;
00284 
00285       nClus=0;
00286 
00287       // Make own simple clusters (3x3, 5x5 or clusPhiSize_ x clusEtaSize_)
00288       sort(seeds.begin(), seeds.end(), ecalRecHitLess());
00289 
00290       for (std::vector<EcalRecHit>::iterator itseed=seeds.begin(); itseed!=seeds.end(); itseed++) {
00291         EBDetId seed_id = itseed->id();
00292         std::vector<EBDetId>::const_iterator usedIds;
00293 
00294         bool seedAlreadyUsed=false;
00295         for(usedIds=usedXtals.begin(); usedIds!=usedXtals.end(); usedIds++){
00296           if(*usedIds==seed_id){
00297             seedAlreadyUsed=true;
00298             //cout<< " Seed with energy "<<itseed->energy()<<" was used !"<<endl;
00299             break; 
00300           }
00301         }
00302         if(seedAlreadyUsed)continue;
00303         topology_p = theCaloTopology->getSubdetectorTopology(DetId::Ecal,EcalBarrel);
00304         std::vector<DetId> clus_v = topology_p->getWindow(seed_id,clusEtaSize_,clusPhiSize_);       
00305         std::vector<DetId> clus_used;
00306         //Reject the seed if not able to build the cluster around it correctly
00307         //if(clus_v.size() < clusEtaSize_*clusPhiSize_){cout<<" Not enough RecHits "<<endl; continue;}
00308         vector<EcalRecHit> RecHitsInWindow;
00309 
00310         double simple_energy = 0; 
00311 
00312         for (std::vector<DetId>::iterator det=clus_v.begin(); det!=clus_v.end(); det++) {
00313           EBDetId EBdet = *det;
00314           //      cout<<" det "<< EBdet<<" ieta "<<EBdet.ieta()<<" iphi "<<EBdet.iphi()<<endl;
00315           bool  HitAlreadyUsed=false;
00316           for(usedIds=usedXtals.begin(); usedIds!=usedXtals.end(); usedIds++){
00317             if(*usedIds==*det){
00318               HitAlreadyUsed=true;
00319               break;
00320             }
00321           }
00322           if(HitAlreadyUsed)continue;
00323           if (recHitsEB_map.find(*det) != recHitsEB_map.end()){
00324             //      cout<<" Used det "<< EBdet<<endl;
00325             std::map<DetId, EcalRecHit>::iterator aHit;
00326             aHit = recHitsEB_map.find(*det);
00327             usedXtals.push_back(*det);
00328             RecHitsInWindow.push_back(aHit->second);
00329             clus_used.push_back(*det);
00330             simple_energy = simple_energy + aHit->second.energy();
00331           }
00332         }
00333    
00334         math::XYZPoint clus_pos = posCalculator_.Calculate_Location(clus_used,hitCollection_p,geometry_p,geometryES_p);
00335         //cout<< "       Simple Clustering: Total energy for this simple cluster : "<<simple_energy<<endl; 
00336         //cout<< "       Simple Clustering: eta phi : "<<clus_pos.eta()<<" "<<clus_pos.phi()<<endl; 
00337         //cout<< "       Simple Clustering: x y z : "<<clus_pos.x()<<" "<<clus_pos.y()<<" "<<clus_pos.z()<<endl; 
00338 
00339         float theta_s = 2. * atan(exp(-clus_pos.eta()));
00340         float p0x_s = simple_energy * sin(theta_s) * cos(clus_pos.phi());
00341         float p0y_s = simple_energy * sin(theta_s) * sin(clus_pos.phi());
00342         //      float p0z_s = simple_energy * cos(theta_s);
00343         float et_s = sqrt( p0x_s*p0x_s + p0y_s*p0y_s);
00344     
00345         //cout << "       Simple Clustering: E,Et,px,py,pz: "<<simple_energy<<" "<<et_s<<" "<<p0x_s<<" "<<p0y_s<<" "<<endl;
00346 
00347         eClus.push_back(simple_energy);
00348         etClus.push_back(et_s);
00349         etaClus.push_back(clus_pos.eta());
00350         phiClus.push_back(clus_pos.phi());
00351         max_hit.push_back(seed_id);
00352         RecHitsCluster.push_back(RecHitsInWindow);
00353         //Compute S4/S9 variable
00354         //We are not sure to have 9 RecHits so need to check eta and phi:
00355         float s4s9_[4];
00356         for(int i=0;i<4;i++)s4s9_[i]= itseed->energy();
00357         for(unsigned int j=0; j<RecHitsInWindow.size();j++){
00358           //cout << " Simple cluster rh, ieta, iphi : "<<((EBDetId)RecHitsInWindow[j].id()).ieta()<<" "<<((EBDetId)RecHitsInWindow[j].id()).iphi()<<endl;
00359         if((((EBDetId)RecHitsInWindow[j].id()).ieta() == seed_id.ieta()-1 && seed_id.ieta()!=1 ) || ( seed_id.ieta()==1 && (((EBDetId)RecHitsInWindow[j].id()).ieta() == seed_id.ieta()-2))){
00360           if(((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi()-1 ||((EBDetId)RecHitsInWindow[j].id()).iphi()-360 == seed_id.iphi()-1 ){
00361             s4s9_[0]+=RecHitsInWindow[j].energy();
00362           }else{
00363             if(((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi()){
00364               s4s9_[0]+=RecHitsInWindow[j].energy();
00365               s4s9_[1]+=RecHitsInWindow[j].energy();
00366             }else{
00367               if(((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi()+1 ||((EBDetId)RecHitsInWindow[j].id()).iphi()-360 == seed_id.iphi()+1 ){
00368                 s4s9_[1]+=RecHitsInWindow[j].energy(); 
00369               }
00370             }
00371           }
00372         }else{
00373           if(((EBDetId)RecHitsInWindow[j].id()).ieta() == seed_id.ieta()){
00374             if(((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi()-1 ||((EBDetId)RecHitsInWindow[j].id()).iphi()-360 == seed_id.iphi()-1 ){
00375               s4s9_[0]+=RecHitsInWindow[j].energy();
00376               s4s9_[3]+=RecHitsInWindow[j].energy();
00377             }else{
00378               if(((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi()+1 ||((EBDetId)RecHitsInWindow[j].id()).iphi()-360 == seed_id.iphi()+1 ){
00379                 s4s9_[1]+=RecHitsInWindow[j].energy(); 
00380                 s4s9_[2]+=RecHitsInWindow[j].energy(); 
00381               }
00382             }
00383           }else{
00384             if((((EBDetId)RecHitsInWindow[j].id()).ieta() == seed_id.ieta()+1 && seed_id.ieta()!=-1 ) || ( seed_id.ieta()==-1 && (((EBDetId)RecHitsInWindow[j].id()).ieta() == seed_id.ieta()+2))){
00385               if(((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi()-1 ||((EBDetId)RecHitsInWindow[j].id()).iphi()-360 == seed_id.iphi()-1 ){
00386                 s4s9_[3]+=RecHitsInWindow[j].energy();
00387               }else{
00388                 if(((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi()){
00389                   s4s9_[2]+=RecHitsInWindow[j].energy();
00390                   s4s9_[3]+=RecHitsInWindow[j].energy();
00391                 }else{
00392                   if(((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi()+1 ||((EBDetId)RecHitsInWindow[j].id()).iphi()-360 == seed_id.iphi()+1 ){
00393                     s4s9_[2]+=RecHitsInWindow[j].energy(); 
00394                   }
00395                 }
00396               }
00397             }else{
00398               cout<<" (EBDetId)RecHitsInWindow[j].id()).ieta() "<<((EBDetId)RecHitsInWindow[j].id()).ieta()<<" seed_id.ieta() "<<seed_id.ieta()<<endl;
00399               cout<<" Problem with S4 calculation "<<endl;return;
00400             }
00401           }
00402         }
00403         }
00404         s4s9Clus.push_back(*max_element( s4s9_,s4s9_+4)/simple_energy);
00405         //    cout<<" s4s9Clus[0] "<<s4s9_[0]/simple_energy<<" s4s9Clus[1] "<<s4s9_[1]/simple_energy<<" s4s9Clus[2] "<<s4s9_[2]/simple_energy<<" s4s9Clus[3] "<<s4s9_[3]/simple_energy<<endl;
00406     //    cout<<" Max "<<*max_element( s4s9_,s4s9_+4)/simple_energy<<endl;
00407     nClus++;
00408     if (nClus == MAXCLUS) return;
00409   }
00410  
00411       // cout<< " Pi0 clusters: "<<nClus<<endl;
00412 
00413       // Selection, based on Simple clustering
00414       //pi0 candidates
00415       static const int MAXPI0S = 200;
00416       int npi0_s=0;
00417 
00418       vector<EBDetId> scXtals;
00419       scXtals.clear();
00420 
00421       if (nClus <= 1) return;
00422       for(Int_t i=0 ; i<nClus ; i++){
00423         for(Int_t j=i+1 ; j<nClus ; j++){
00424           //      cout<<" i "<<i<<"  etClus[i] "<<etClus[i]<<" j "<<j<<"  etClus[j] "<<etClus[j]<<endl;
00425           if( etClus[i]>selePtGammaOne_ && etClus[j]>selePtGammaTwo_ && s4s9Clus[i]>seleS4S9GammaOne_ && s4s9Clus[j]>seleS4S9GammaTwo_){
00426             float theta_0 = 2. * atan(exp(-etaClus[i]));
00427             float theta_1 = 2. * atan(exp(-etaClus[j]));
00428         
00429             float p0x = eClus[i] * sin(theta_0) * cos(phiClus[i]);
00430             float p1x = eClus[j] * sin(theta_1) * cos(phiClus[j]);
00431             float p0y = eClus[i] * sin(theta_0) * sin(phiClus[i]);
00432             float p1y = eClus[j] * sin(theta_1) * sin(phiClus[j]);
00433             float p0z = eClus[i] * cos(theta_0);
00434             float p1z = eClus[j] * cos(theta_1);
00435         
00436             float pt_pi0 = sqrt( (p0x+p1x)*(p0x+p1x) + (p0y+p1y)*(p0y+p1y));
00437             //      cout<<" pt_pi0 "<<pt_pi0<<endl;
00438             if (pt_pi0 < selePtPi0_)continue;
00439             float m_inv = sqrt ( (eClus[i] + eClus[j])*(eClus[i] + eClus[j]) - (p0x+p1x)*(p0x+p1x) - (p0y+p1y)*(p0y+p1y) - (p0z+p1z)*(p0z+p1z) );  
00440             if ( (m_inv<seleMinvMaxPi0_) && (m_inv>seleMinvMinPi0_) ){
00441 
00442               //New Loop on cluster to measure isolation:
00443               vector<int> IsoClus;
00444               IsoClus.clear();
00445               float Iso = 0;
00446               TVector3 pi0vect = TVector3((p0x+p1x), (p0y+p1y), (p0z+p1z));
00447               for(Int_t k=0 ; k<nClus ; k++){
00448                 if(k==i || k==j)continue;
00449                 TVector3 Clusvect = TVector3(eClus[k] * sin(2. * atan(exp(-etaClus[k]))) * cos(phiClus[k]), eClus[k] * sin(2. * atan(exp(-etaClus[k]))) * sin(phiClus[k]) , eClus[k] * cos(2. * atan(exp(-etaClus[k]))));
00450                 float dretaclpi0 = fabs(etaClus[k] - pi0vect.Eta());
00451                 float drclpi0 = Clusvect.DeltaR(pi0vect);
00452                 //      cout<< "   Iso: k, E, drclpi0, detaclpi0, dphiclpi0 "<<k<<" "<<eClus[k]<<" "<<drclpi0<<" "<<dretaclpi0<<endl;
00453                 if((drclpi0<selePi0BeltDR_) && (dretaclpi0<selePi0BeltDeta_) ){
00454                   //              cout<< "   ... good iso cluster #: "<<k<<" etClus[k] "<<etClus[k] <<endl;
00455                   Iso = Iso + etClus[k];
00456                   IsoClus.push_back(k);
00457                 }
00458               }
00459 
00460               //      cout<<"  Iso/pt_pi0 "<<Iso/pt_pi0<<endl;
00461               if(Iso/pt_pi0<selePi0Iso_){
00462                 //for(unsigned int Rec=0;Rec<RecHitsCluster[i].size();Rec++)pi0EBRecHitCollection->push_back(RecHitsCluster[i][Rec]);
00463                 //for(unsigned int Rec2=0;Rec2<RecHitsCluster[j].size();Rec2++)pi0EBRecHitCollection->push_back(RecHitsCluster[j][Rec2]);
00464                 
00465                 
00466                 hMinvPi0EB_->Fill(m_inv);
00467                 hPt1Pi0EB_->Fill(etClus[i]);
00468                 hPt2Pi0EB_->Fill(etClus[j]);
00469                 hPtPi0EB_->Fill(pt_pi0);
00470                 hIsoPi0EB_->Fill(Iso/pt_pi0);
00471                 hS4S91EB_->Fill(s4s9Clus[i]);
00472                 hS4S92EB_->Fill(s4s9Clus[j]);
00473                 
00474                 npi0_s++;
00475               }
00476               //                for(unsigned int iii=0 ; iii<IsoClus.size() ; iii++){   
00477                 //cout<< "    Iso cluster: # "<<iii<<" "<<IsoClus[iii]<<endl;  
00478               //                  for(unsigned int Rec3=0;Rec3<RecHitsCluster[IsoClus[iii]].size();Rec3++)pi0EBRecHitCollection->push_back(RecHitsCluster[IsoClus[iii]][Rec3]);
00479               //}   
00480             
00481                 //cout <<"  Simple Clustering: pi0 Candidate pt,m_inv,i,j :   "<<pt_pi0<<" "<<m_inv<<" "<<i<<" "<<j<<" "<<endl;  
00482 
00483           
00484               if(npi0_s == MAXPI0S) return;
00485             }
00486           }
00487         } // End of the "j" loop over Simple Clusters
00488       } // End of the "i" loop over Simple Clusters
00489 
00490       //cout<<"  (Simple Clustering) Pi0 candidates #: "<<npi0_s<<endl;
00491 
00492       //      delete recHitsEB_map;
00493 
00494 
00495 
00496 
00497 
00498 
00499 
00500 
00501 
00502       
00503     } // if valid
00504 
00505   } // if isMonEB
00506 
00507   // fill EE histos
00508 
00509   if(isMonEE_){  
00510     EcalRecHitCollection::const_iterator ite;
00511     
00512     if (rhEE.isValid()){
00513       
00514       float etot =0;
00515       for(ite=rhEE->begin(); ite!=rhEE->end(); ++ite){
00516         
00517         EEDetId id(ite->id());
00518         hRechitEnergyEE_->Fill(ite->energy());
00519         etot+= ite->energy();    
00520       } // EE rechits
00521       
00522       hNRecHitsEE_->Fill(rhEE->size());
00523       hMeanRecHitEnergyEE_->Fill(etot/rhEE->size());
00524       hEventEnergyEE_->Fill(etot);
00525     }
00526     
00527   }//isMonEE
00528 } //analyze

void DQMHLTSourcePi0::beginJob ( const edm::EventSetup c  )  [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 106 of file DQMHLTSourcePi0.cc.

References DQMStore::book1D(), dbe_, folderName_, hEventEnergyEB_, hEventEnergyEE_, hiEtaDistrEB_, hiPhiDistrEB_, hIsoPi0EB_, hMeanRecHitEnergyEB_, hMeanRecHitEnergyEE_, hMinvPi0EB_, hNRecHitsEB_, hNRecHitsEE_, hPt1Pi0EB_, hPt2Pi0EB_, hPtPi0EB_, hRechitEnergyEB_, hRechitEnergyEE_, hS4S91EB_, hS4S92EB_, MonitorElement::setAxisTitle(), and DQMStore::setCurrentFolder().

00106                                                        {
00107 
00108 
00109   // create and cd into new folder
00110   dbe_->setCurrentFolder(folderName_);
00111 
00112   // book some histograms 1D
00113 
00114   hiPhiDistrEB_ = dbe_->book1D("iphiDistributionEB", "RechitEB iphi", 361, 1,361);
00115 
00116   hiPhiDistrEB_->setAxisTitle("i#phi ", 1);
00117   hiPhiDistrEB_->setAxisTitle("# rechits", 2);
00118 
00119 
00120   hiEtaDistrEB_ = dbe_->book1D("iEtaDistributionEB", "RechitEB ieta", 171, -85, 86);
00121   hiEtaDistrEB_->setAxisTitle("eta", 1);
00122   hiEtaDistrEB_->setAxisTitle("#rechits", 2);
00123 
00124 
00125   hRechitEnergyEB_ = dbe_->book1D("rhEnergyEB","rechits energy EB",160,0.,2.0);
00126   hRechitEnergyEB_->setAxisTitle("energy (GeV) ",1);
00127   hRechitEnergyEB_->setAxisTitle("#rechits",2);
00128 
00129   hEventEnergyEB_ = dbe_->book1D("eventEnergyEB","event energy EB",100,0.,20.0);
00130   hEventEnergyEB_->setAxisTitle("energy (GeV) ",1);
00131 
00132   hNRecHitsEB_ = dbe_->book1D("nRechitsEB","#rechits in event EB",100,0.,250.);
00133   hNRecHitsEB_->setAxisTitle("rechits ",1);
00134   
00135   hMeanRecHitEnergyEB_ = dbe_->book1D("meanEnergyEB","Mean rechit energy EB",50,0.,2.);
00136   hMeanRecHitEnergyEB_->setAxisTitle("Mean Energy [GeV] ",1);
00137   
00138   hMinvPi0EB_ = dbe_->book1D("Pi0InvmassEB","Pi0 Invariant Mass in EB",100,0.,0.5);
00139   hMinvPi0EB_->setAxisTitle("Inv Mass [GeV] ",1);
00140 
00141   
00142   hPt1Pi0EB_ = dbe_->book1D("Pt1Pi0EB","Pt 1st most energetic Pi0 photon in EB",100,0.,20.);
00143   hPt1Pi0EB_->setAxisTitle("1st photon Pt [GeV] ",1);
00144   
00145   hPt2Pi0EB_ = dbe_->book1D("Pt2Pi0EB","Pt 2nd most energetic Pi0 photon in EB",100,0.,20.);
00146   hPt2Pi0EB_->setAxisTitle("2nd photon Pt [GeV] ",1);
00147 
00148   
00149   hPtPi0EB_ = dbe_->book1D("PtPi0EB","Pi0 Pt in EB",100,0.,20.);
00150   hPtPi0EB_->setAxisTitle("Pi0 Pt [GeV] ",1);
00151 
00152   hIsoPi0EB_ = dbe_->book1D("IsoPi0EB","Pi0 Iso in EB",50,0.,1.);
00153   hIsoPi0EB_->setAxisTitle("Pi0 Iso",1);
00154 
00155   hS4S91EB_ = dbe_->book1D("S4S91EB","S4S9 1st most energetic Pi0 photon in EB",50,0.,1.);
00156   hS4S91EB_->setAxisTitle("S4S9 of the 1st Pi0 Photon ",1);
00157 
00158   hS4S92EB_ = dbe_->book1D("S4S92EB","S4S9 2nd most energetic Pi0 photon in EB",50,0.,1.);
00159   hS4S92EB_->setAxisTitle("S4S9 of the 2nd Pi0 Photon",1);
00160 
00161   
00162 
00163   hRechitEnergyEE_ = dbe_->book1D("rhEnergyEE","rechits energy EE",160,0.,3.0);
00164   hRechitEnergyEE_->setAxisTitle("energy (GeV) ",1);
00165   hRechitEnergyEE_->setAxisTitle("#rechits",2);
00166 
00167   hEventEnergyEE_ = dbe_->book1D("eventEnergyEE","event energy EE",100,0.,20.0);
00168   hEventEnergyEE_->setAxisTitle("energy (GeV) ",1);
00169 
00170   hNRecHitsEE_ = dbe_->book1D("nRechitsEE","#rechits in event EE" ,100,0.,250.);
00171   hNRecHitsEE_->setAxisTitle("rechits ",1);
00172  
00173   hMeanRecHitEnergyEE_ = dbe_->book1D("meanEnergyEE","Mean rechit energy EE",50,0.,5.);
00174   hMeanRecHitEnergyEE_-> setAxisTitle("Mean Energy [GeV] ",1);
00175   
00176 
00177 }

void DQMHLTSourcePi0::beginLuminosityBlock ( const edm::LuminosityBlock lumiSeg,
const edm::EventSetup context 
) [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 185 of file DQMHLTSourcePi0.cc.

00186                                 {
00187   
00188 }

void DQMHLTSourcePi0::beginRun ( const edm::Run r,
const edm::EventSetup c 
) [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 180 of file DQMHLTSourcePi0.cc.

00180                                                                          {
00181 
00182 }

void DQMHLTSourcePi0::endJob ( void   )  [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 542 of file DQMHLTSourcePi0.cc.

References dbe_, fileName_, DQMStore::save(), and saveToFile_.

00542                             {
00543 
00544   if(dbe_) {  
00545     if (saveToFile_) {
00546       dbe_->save(fileName_);
00547     }
00548   }
00549 }

void DQMHLTSourcePi0::endLuminosityBlock ( const edm::LuminosityBlock lumiSeg,
const edm::EventSetup c 
) [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 534 of file DQMHLTSourcePi0.cc.

00535                                                                      {
00536 }

void DQMHLTSourcePi0::endRun ( const edm::Run r,
const edm::EventSetup c 
) [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 538 of file DQMHLTSourcePi0.cc.

00538                                                                    {
00539 
00540 }


Member Data Documentation

int DQMHLTSourcePi0::clusEtaSize_ [private]

Definition at line 137 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

int DQMHLTSourcePi0::clusPhiSize_ [private]

Definition at line 138 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

double DQMHLTSourcePi0::clusSeedThr_ [private]

Definition at line 136 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

DQMStore* DQMHLTSourcePi0::dbe_ [private]

Definition at line 68 of file DQMHLTSourcePi0.h.

Referenced by beginJob(), DQMHLTSourcePi0(), and endJob().

int DQMHLTSourcePi0::eventCounter_ [private]

Definition at line 69 of file DQMHLTSourcePi0.h.

Referenced by analyze().

std::string DQMHLTSourcePi0::fileName_ [private]

Output file name if required.

Definition at line 175 of file DQMHLTSourcePi0.h.

Referenced by DQMHLTSourcePi0(), and endJob().

std::string DQMHLTSourcePi0::folderName_ [private]

DQM folder name.

Definition at line 165 of file DQMHLTSourcePi0.h.

Referenced by beginJob(), and DQMHLTSourcePi0().

int DQMHLTSourcePi0::gammaCandEtaSize_ [private]

Definition at line 133 of file DQMHLTSourcePi0.h.

Referenced by DQMHLTSourcePi0().

int DQMHLTSourcePi0::gammaCandPhiSize_ [private]

Definition at line 134 of file DQMHLTSourcePi0.h.

Referenced by DQMHLTSourcePi0().

MonitorElement* DQMHLTSourcePi0::hEventEnergyEB_ [private]

Distribution of total event energy.

Definition at line 82 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hEventEnergyEE_ [private]

Distribution of total event energy.

Definition at line 119 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hiEtaDistrEB_ [private]

Distribution of rechits in iEta.

Definition at line 76 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hiPhiDistrEB_ [private]

Distribution of rechits in iPhi.

Definition at line 73 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hIsoPi0EB_ [private]

Pi0 Iso.

Definition at line 105 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hMeanRecHitEnergyEB_ [private]

Distribution of Mean energy per rechit.

Definition at line 88 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hMeanRecHitEnergyEE_ [private]

Distribution of Mean energy per rechit.

Definition at line 125 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hMinvPi0EB_ [private]

Pi0 invariant mass in EB.

Definition at line 91 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hNRecHitsEB_ [private]

Distribution of number of RecHits.

Definition at line 85 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hNRecHitsEE_ [private]

Distribution of number of RecHits.

Definition at line 122 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hPt1Pi0EB_ [private]

Pt of the 1st most energetic Pi0 photon in EB.

Definition at line 94 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hPt2Pi0EB_ [private]

Pt of the 2nd most energetic Pi0 photon in EB.

Definition at line 98 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hPtPi0EB_ [private]

Pi0 Pt in EB.

Definition at line 102 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hRechitEnergyEB_ [private]

Energy Distribution of rechits.

Definition at line 79 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hRechitEnergyEE_ [private]

Energy Distribution of rechits.

Definition at line 116 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hS4S91EB_ [private]

S4S9 of the 1st most energetic pi0 photon.

Definition at line 108 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

MonitorElement* DQMHLTSourcePi0::hS4S92EB_ [private]

S4S9 of the 2nd most energetic pi0 photon.

Definition at line 111 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and beginJob().

bool DQMHLTSourcePi0::isMonEB_ [private]

which subdet will be monitored

Definition at line 171 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

bool DQMHLTSourcePi0::isMonEE_ [private]

Definition at line 172 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

bool DQMHLTSourcePi0::ParameterLogWeighted_ [private]

Definition at line 153 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

double DQMHLTSourcePi0::ParameterT0_barl_ [private]

Definition at line 155 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

double DQMHLTSourcePi0::ParameterW0_ [private]

Definition at line 156 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

double DQMHLTSourcePi0::ParameterX0_ [private]

Definition at line 154 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

unsigned int DQMHLTSourcePi0::prescaleFactor_ [private]

Monitor every prescaleFactor_ events.

Definition at line 162 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

edm::InputTag DQMHLTSourcePi0::productMonitoredEB_ [private]

object to monitor

Definition at line 128 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

edm::InputTag DQMHLTSourcePi0::productMonitoredEE_ [private]

object to monitor

Definition at line 131 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

bool DQMHLTSourcePi0::saveToFile_ [private]

Write to file.

Definition at line 168 of file DQMHLTSourcePi0.h.

Referenced by DQMHLTSourcePi0(), and endJob().

double DQMHLTSourcePi0::seleMinvMaxPi0_ [private]

Definition at line 143 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

double DQMHLTSourcePi0::seleMinvMinPi0_ [private]

Definition at line 144 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

int DQMHLTSourcePi0::seleNRHMax_ [private]

Definition at line 146 of file DQMHLTSourcePi0.h.

Referenced by DQMHLTSourcePi0().

double DQMHLTSourcePi0::selePi0BeltDeta_ [private]

Definition at line 151 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

double DQMHLTSourcePi0::selePi0BeltDR_ [private]

Definition at line 150 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

double DQMHLTSourcePi0::selePi0Iso_ [private]

Definition at line 152 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

double DQMHLTSourcePi0::selePtGammaOne_ [private]

Definition at line 140 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

double DQMHLTSourcePi0::selePtGammaTwo_ [private]

Definition at line 141 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

double DQMHLTSourcePi0::selePtPi0_ [private]

Definition at line 142 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

double DQMHLTSourcePi0::seleS4S9GammaOne_ [private]

Definition at line 148 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

double DQMHLTSourcePi0::seleS4S9GammaTwo_ [private]

Definition at line 149 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().

double DQMHLTSourcePi0::seleXtalMinEnergy_ [private]

Definition at line 145 of file DQMHLTSourcePi0.h.

Referenced by analyze(), and DQMHLTSourcePi0().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:18:35 2009 for CMSSW by  doxygen 1.5.4