CMS 3D CMS Logo

AlCaPi0RecHitsProducer Class Reference

#include <Calibration/EcalAlCaRecoProducers/interface/AlCaPi0RecHitsProducer.h>

Inheritance diagram for AlCaPi0RecHitsProducer:

edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 AlCaPi0RecHitsProducer (const edm::ParameterSet &)
virtual void produce (edm::Event &, const edm::EventSetup &)
 ~AlCaPi0RecHitsProducer ()

Private Attributes

std::string barrelHits_
int clusEtaSize_
int clusPhiSize_
double clusSeedThr_
std::string ecalHitsProducer_
int gammaCandEtaSize_
int gammaCandPhiSize_
std::string pi0BarrelHits_
std::map< DetId, EcalRecHit > * recHitsEB_map
double seleMinvMaxPi0_
double seleMinvMinPi0_
int seleNRHMax_
double selePtGammaOne_
double selePtGammaTwo_
double selePtPi0_
double seleXtalMinEnergy_


Detailed Description

Definition at line 47 of file AlCaPi0RecHitsProducer.h.


Constructor & Destructor Documentation

AlCaPi0RecHitsProducer::AlCaPi0RecHitsProducer ( const edm::ParameterSet iConfig  )  [explicit]

Definition at line 18 of file AlCaPi0RecHitsProducer.cc.

References barrelHits_, clusEtaSize_, clusPhiSize_, clusSeedThr_, ecalHitsProducer_, gammaCandEtaSize_, gammaCandPhiSize_, edm::ParameterSet::getParameter(), pi0BarrelHits_, seleMinvMaxPi0_, seleMinvMinPi0_, seleNRHMax_, selePtGammaOne_, selePtGammaTwo_, selePtPi0_, and seleXtalMinEnergy_.

00019 {
00020   ecalHitsProducer_ = iConfig.getParameter< std::string > ("ecalRecHitsProducer");
00021   barrelHits_ = iConfig.getParameter< std::string > ("barrelHitCollection");
00022   pi0BarrelHits_ = iConfig.getParameter< std::string > ("pi0BarrelHitCollection");
00023 
00024   gammaCandEtaSize_ = iConfig.getParameter<int> ("gammaCandEtaSize");
00025   gammaCandPhiSize_ = iConfig.getParameter<int> ("gammaCandPhiSize");
00026   if ( gammaCandPhiSize_ % 2 == 0 ||  gammaCandEtaSize_ % 2 == 0)
00027     edm::LogError("AlCaPi0RecHitsProducerError") << "Size of eta/phi for sliding window should be odd numbers";
00028 
00029   clusSeedThr_ = iConfig.getParameter<double> ("clusSeedThr");
00030   clusEtaSize_ = iConfig.getParameter<int> ("clusEtaSize");
00031   clusPhiSize_ = iConfig.getParameter<int> ("clusPhiSize");
00032   if ( clusPhiSize_ % 2 == 0 ||  clusEtaSize_ % 2 == 0)
00033     edm::LogError("AlCaPi0RecHitsProducerError") << "Size of eta/phi for simple clustering should be odd numbers";
00034 
00035   selePtGammaOne_ = iConfig.getParameter<double> ("selePtGammaOne");  
00036   selePtGammaTwo_ = iConfig.getParameter<double> ("selePtGammaTwo");  
00037   selePtPi0_ = iConfig.getParameter<double> ("selePtPi0");  
00038   seleMinvMaxPi0_ = iConfig.getParameter<double> ("seleMinvMaxPi0");  
00039   seleMinvMinPi0_ = iConfig.getParameter<double> ("seleMinvMinPi0");  
00040   seleXtalMinEnergy_ = iConfig.getParameter<double> ("seleXtalMinEnergy");
00041   seleNRHMax_ = iConfig.getParameter<int> ("seleNRHMax");
00042 
00043   //register your products
00044   produces< EBRecHitCollection >(pi0BarrelHits_);
00045 }

AlCaPi0RecHitsProducer::~AlCaPi0RecHitsProducer (  ) 

Definition at line 48 of file AlCaPi0RecHitsProducer.cc.

References GenMuonPlsPt100GeV_cfg::cout, TimingReport::current(), and TimingReport::dump().

00049 {
00050  
00051   TimingReport::current()->dump(std::cout);
00052 
00053 }


Member Function Documentation

void AlCaPi0RecHitsProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements edm::EDProducer.

Definition at line 58 of file AlCaPi0RecHitsProducer.cc.

References barrelHits_, PositionCalc::Calculate_Location(), category, TimerStack::clear_stack(), clusEtaSize_, clusPhiSize_, clusSeedThr_, funct::cos(), GenMuonPlsPt100GeV_cfg::cout, DetId::Ecal, EcalBarrel, ecalHitsProducer_, EcalPreshower, lat::endl(), relval_parameters_module::energy, EBDetId::ETAPHIMODE, funct::exp(), gammaCandEtaSize_, gammaCandPhiSize_, edm::EventSetup::get(), edm::Event::getByLabel(), i, EBDetId::ieta(), EBDetId::iphi(), j, LogDebug, pi0BarrelHits_, TimerStack::pop_and_push(), TimerStack::push(), edm::Event::put(), recHitsEB_map, row, seleMinvMaxPi0_, seleMinvMinPi0_, seleNRHMax_, selePtGammaOne_, selePtGammaTwo_, selePtPi0_, seleXtalMinEnergy_, funct::sin(), python::multivaluedict::sort(), funct::sqrt(), and std.

00059 {
00060   using namespace edm;
00061   using namespace std;
00062 
00063   // Timer
00064   const std::string category = "AlCaPi0RecHitsProducer";
00065   TimerStack timers;
00066   string timerName = category + "::Total";
00067   timers.push(timerName);
00068 
00069 
00070   Handle<EBRecHitCollection> barrelRecHitsHandle;
00071 
00072   iEvent.getByLabel(ecalHitsProducer_,barrelHits_,barrelRecHitsHandle);
00073   if (!barrelRecHitsHandle.isValid()) {
00074     LogDebug("") << "AlCaPi0RecHitsProducer: Error! can't get product!" << std::endl;
00075   }
00076 
00077   recHitsEB_map = new std::map<DetId, EcalRecHit>();
00078 
00079   std::vector<EcalRecHit> seeds;
00080   seeds.clear();
00081 
00082   vector<EBDetId> usedXtals;
00083   usedXtals.clear();
00084 
00085 
00086   //Create empty output collections
00087   std::auto_ptr< EBRecHitCollection > pi0EBRecHitCollection( new EBRecHitCollection );
00088   std::auto_ptr< EBRecHitCollection > pi0EBDummyRecHitCollection( new EBRecHitCollection );
00089 
00090   //Select interesting EcalRecHits (barrel)
00091   EBRecHitCollection::const_iterator itb;
00092   cout<< "   EB RecHits #: "<<barrelRecHitsHandle->size()<<endl;
00093   for (itb=barrelRecHitsHandle->begin(); itb!=barrelRecHitsHandle->end(); itb++) {
00094 
00095     std::pair<DetId, EcalRecHit> map_entry(itb->id(), *itb);
00096     recHitsEB_map->insert(map_entry);
00097 
00098     double energy = itb->energy();
00099     if (energy > clusSeedThr_) seeds.push_back(*itb);
00100   }
00101 
00102   timerName = category + "::readEBRecHitsCollection";
00103   timers.push(timerName);
00104 
00105 
00106   // Initialize the Position Calc
00107   const CaloSubdetectorGeometry *geometry_p;    
00108   const CaloSubdetectorTopology *topology_p;
00109   const CaloSubdetectorGeometry *geometryES_p;
00110   const EcalRecHitCollection *hitCollection_p = barrelRecHitsHandle.product();
00111 
00112   edm::ESHandle<CaloGeometry> geoHandle;
00113   iSetup.get<CaloGeometryRecord>().get(geoHandle);     
00114   geometry_p = geoHandle->getSubdetectorGeometry(DetId::Ecal,EcalBarrel);
00115   topology_p = new EcalBarrelTopology(geoHandle);
00116   geometryES_p = geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalPreshower);
00117 
00118   std::map<std::string,double> providedParameters;  
00119   providedParameters.insert(std::make_pair("LogWeighted",1));
00120   providedParameters.insert(std::make_pair("X0",0.89));
00121   providedParameters.insert(std::make_pair("T0_barl",5.7));
00122   providedParameters.insert(std::make_pair("W0",4.2));
00123 
00124   PositionCalc posCalculator_ = PositionCalc(providedParameters);
00125   //PositionCalc::Initialize(providedParameters, &recHitsEB_map, &(*geometry_p));
00126 
00127   static const int MAXCLUS = 200;
00128   int nClus;
00129   float eClus[MAXCLUS];
00130   float etClus[MAXCLUS];
00131   float etaClus[MAXCLUS];
00132   float phiClus[MAXCLUS];
00133   EBDetId max_hit[MAXCLUS];  
00134 
00135   nClus=0;
00136   for(int i=0; i<MAXCLUS; i++){
00137     eClus[i] = 0;
00138     etClus[i] = 0;
00139     etaClus[i] = 0;
00140     phiClus[i] = 0;
00141     max_hit[i] = EBDetId(0);
00142   }
00143 
00144   // Make own simple clusters (3x3, 5x5 or clusPhiSize_ x clusEtaSize_)
00145   sort(seeds.begin(), seeds.end(), ecalRecHitLess());
00146 
00147 
00148   for (std::vector<EcalRecHit>::iterator itseed=seeds.begin(); itseed!=seeds.end(); itseed++) {
00149     EBDetId seed_id = itseed->id();
00150     std::vector<EBDetId>::const_iterator usedIds;
00151     
00152     //cout<< " Start: Seed with energy "<<itseed->energy()<<endl;
00153     //cout<< " Start: Seed with z,ieta,iphi : "<<seed_id.zside()<<" "<<seed_id.ieta()<<" " <<seed_id.iphi()<<endl;
00154     bool seedAlreadyUsed=false;
00155     for(usedIds=usedXtals.begin(); usedIds!=usedXtals.end(); usedIds++)
00156       if(*usedIds==seed_id)
00157         {
00158           seedAlreadyUsed=true;
00159           //cout<< " Seed with energy "<<itseed->energy()<<" was used !"<<endl;
00160           break;
00161         }
00162     if(!seedAlreadyUsed)            
00163       {
00164         double simple_energy = 0; 
00165         std::vector<DetId> clus_v;
00166         clus_v.clear();
00167     for (int icry=0;icry< clusEtaSize_*clusPhiSize_;icry++)
00168       {
00169         
00170         int row = icry / clusEtaSize_;
00171         int column= icry % clusEtaSize_;
00172         int curr_eta=seed_id.ieta() + column - (clusEtaSize_/2);
00173         int curr_phi=seed_id.iphi() + row - (clusPhiSize_/2);
00174         
00175         if (curr_eta * seed_id.ieta() <= 0) {if (seed_id.ieta() > 0) curr_eta--; else curr_eta++; }  // JUMP over 0
00176         if (curr_phi < 1) curr_phi += 360;
00177         if (curr_phi > 360) curr_phi -= 360;
00178         
00179         try
00180           {
00181             EBDetId det = EBDetId(curr_eta,curr_phi,EBDetId::ETAPHIMODE);
00182             std::vector<EBDetId>::const_iterator usedIds;
00183                           
00184             std::map<DetId, EcalRecHit>::iterator aHit;
00185             bool HitAlreadyUsed=false;
00186             for(usedIds=usedXtals.begin(); usedIds!=usedXtals.end(); usedIds++)
00187               if(*usedIds==det)
00188                 {
00189                   HitAlreadyUsed=true;
00190                   break;
00191                 }
00192             
00193             if(!HitAlreadyUsed)
00194               if (recHitsEB_map->find(det) != recHitsEB_map->end())
00195                 {
00196                   aHit = recHitsEB_map->find(det);
00197                   usedXtals.push_back(det);
00198                   clus_v.push_back(det);
00199 
00200                   simple_energy = simple_energy + aHit->second.energy();
00201                   
00202                   //EBDetId sel_rh = aHit->second.detid();
00203                   //cout << "       Simple Clustering: RecHit Ok 3x3 matrix inside cluster : z,ieta,iphi "<<sel_rh.zside()<<" "<<sel_rh.ieta()<<" "<<sel_rh.iphi()<<endl;    
00204                   //cout << "       Simple Clustering: RecHit Ok 3x3 matrix inside cluster : tower_ieta,tower_iphi "<<sel_rh.tower_ieta()<<" "<<sel_rh.tower_iphi()<<endl;   
00205                   //cout << "       Simple Clustering: RecHit Ok 3x3 matrix inside cluster : iSM, ic "<<sel_rh.ism()<<" "<<sel_rh.ic()<<endl;
00206                   
00207                 }
00208           }
00209         catch (...)
00210           {
00211           }
00212       }
00213     math::XYZPoint clus_pos = posCalculator_.Calculate_Location(clus_v,hitCollection_p,geometry_p,geometryES_p);
00214     //cout<< "       Simple Clustering: Total energy for this simple cluster : "<<simple_energy<<endl; 
00215     //cout<< "       Simple Clustering: eta phi : "<<clus_pos.eta()<<" "<<clus_pos.phi()<<endl; 
00216     //cout<< "       Simple Clustering: x y z : "<<clus_pos.x()<<" "<<clus_pos.y()<<" "<<clus_pos.z()<<endl; 
00217 
00218             float theta_s = 2. * atan(exp(-clus_pos.eta()));
00219             float p0x_s = simple_energy * sin(theta_s) * cos(clus_pos.phi());
00220             float p0y_s = simple_energy * sin(theta_s) * sin(clus_pos.phi());
00221             float p0z_s = simple_energy * cos(theta_s);
00222             float et_s = sqrt( p0x_s*p0x_s + p0y_s*p0y_s);
00223 
00224             //cout << "       Simple Clustering: E,Et,px,py,pz: "<<simple_energy<<" "<<et_s<<" "<<p0x_s<<" "<<p0y_s<<" "<<p0z_s<<endl;
00225     
00226             eClus[nClus] = simple_energy;
00227             etClus[nClus] = et_s;
00228             etaClus[nClus] = clus_pos.eta();
00229             phiClus[nClus] = clus_pos.phi();
00230             max_hit[nClus] = seed_id;
00231             
00232             nClus++;
00233       }
00234   }
00235   
00236   timerName = category + "::makeSimpleClusters";
00237   timers.pop_and_push(timerName);
00238 
00239 
00240   // Selection, based on Simple clustering
00241   //pi0 candidates
00242   static const int MAXPI0S = 20;
00243   int npi0_s=0;
00244   int sClus_1[MAXPI0S],sClus_2[MAXPI0S];
00245   for(int i=0; i<MAXPI0S; i++){
00246     sClus_1[i]=0;
00247     sClus_2[i]=0;
00248   }
00249 
00250   if (nClus > 1) 
00251     {
00252       for(Int_t i=0 ; i<nClus ; i++)
00253         {
00254           for(Int_t j=i+1 ; j<nClus ; j++)
00255             {
00256 
00257               if( etClus[i]>selePtGammaOne_ && etClus[j]>selePtGammaTwo_) 
00258                 {
00259                 
00260                   float theta_0 = 2. * atan(exp(-etaClus[i]));
00261                   float theta_1 = 2. * atan(exp(-etaClus[j]));
00262                 
00263                   float p0x = eClus[i] * sin(theta_0) * cos(phiClus[i]);
00264                   float p1x = eClus[j] * sin(theta_1) * cos(phiClus[j]);
00265                   float p0y = eClus[i] * sin(theta_0) * sin(phiClus[i]);
00266                   float p1y = eClus[j] * sin(theta_1) * sin(phiClus[j]);
00267                   float p0z = eClus[i] * cos(theta_0);
00268                   float p1z = eClus[j] * cos(theta_1);
00269                 
00270                   //
00271                 
00272                   float pt_pi0 = sqrt( (p0x+p1x)*(p0x+p1x) + (p0y+p1y)*(p0y+p1y));
00273                   //float dr_pi0 = sqrt ( (etaIslandBCEB[i]-etaIslandBCEB[j])*(etaIslandBCEB[i]-etaIslandBCEB[j]) + (phiIslandBCEB[i]-phiIslandBCEB[j])*(phiIslandBCEB[i]-phiIslandBCEB[j]) );
00274                   if (pt_pi0 > selePtPi0_) 
00275                     {
00276                       float m_inv = sqrt ( (eClus[i] + eClus[j])*(eClus[i] + eClus[j]) - (p0x+p1x)*(p0x+p1x) - (p0y+p1y)*(p0y+p1y) - (p0z+p1z)*(p0z+p1z) );  
00277                       if ( (m_inv<seleMinvMaxPi0_) && (m_inv>seleMinvMinPi0_) )
00278                         {
00279                           cout <<"  Simple Clustering: pi0 Candidate (pt>2.5 GeV, m_inv<0.2) pt,m_inv,i,j :   "<<pt_pi0<<" "<<m_inv<<" "<<i<<" "<<j<<" "<<endl;  
00280                           sClus_1[npi0_s]=i;
00281                           sClus_2[npi0_s]=j;
00282                           npi0_s++;
00283                         }
00284 
00285                     }
00286                 
00287                 }
00288             } // End of the "j" loop over Simple Clusters
00289         } // End of the "i" loop over Simple Clusters
00290 
00291     }
00292 
00293   timerName = category + "::makePi0Cand";
00294   timers.pop_and_push(timerName);
00295 
00296 
00297       cout<<" "<<endl;
00298       cout<<"  (Simple Clustering) Pi0 candidates #: "<<npi0_s<<endl;
00299 
00300       for(Int_t jj=0;jj<nClus;jj++)
00301         {
00302           EBDetId maxHit = max_hit[jj];
00303           //cout << "        Simple Clustering: maxHit jj eta phi : "<<jj<<" "<<maxHit.ieta()<<" "<<maxHit.iphi()<<endl;
00304         }
00305 
00306       vector<EBDetId> scXtals;
00307       scXtals.clear();
00308       for(Int_t i=0 ; i<npi0_s ; i++)
00309         {
00310           //cout<<"     Pi0: i, Bc1, Bc2 "<<i<<" "<<sClus_1[i]<<" "<<sClus_2[i]<<endl; 
00311 
00312           for(Int_t j=0; j<nClus;j++)
00313             {
00314 
00315               EBDetId maxHit = max_hit[j];
00316               if( (sClus_1[i] == j) || (sClus_2[i] == j) )
00317                 {
00318 
00319                   for (int icry=0;icry< gammaCandEtaSize_*gammaCandPhiSize_;icry++)
00320                     {
00321               
00322                       int row = icry / gammaCandEtaSize_;
00323                       int column= icry % gammaCandEtaSize_;
00324                       int curr_eta=maxHit.ieta() + column - (gammaCandEtaSize_/2);
00325                       int curr_phi=maxHit.iphi() + row - (gammaCandPhiSize_/2);
00326               
00327                       if (curr_eta * maxHit.ieta() <= 0) {if (maxHit.ieta() > 0) curr_eta--; else curr_eta++; }  // JUMP over 0
00328                       if (curr_phi < 1) curr_phi += 360;
00329                       if (curr_phi > 360) curr_phi -= 360;
00330               
00331                       try
00332                         {
00333                           EBDetId det = EBDetId(curr_eta,curr_phi,EBDetId::ETAPHIMODE);
00334                           std::vector<EBDetId>::const_iterator usedIds;
00335                           
00336                           bool HitAlreadyUsed=false;
00337                           for(usedIds=scXtals.begin(); usedIds!=scXtals.end(); usedIds++)
00338                             if(*usedIds==det)
00339                               {
00340                                 HitAlreadyUsed=true;
00341                                 break;
00342                               }
00343                           
00344                           if(!HitAlreadyUsed)
00345                             if (recHitsEB_map->find(det) != recHitsEB_map->end())
00346                               {
00347                                 std::map<DetId, EcalRecHit>::iterator aHit;
00348                                 aHit = recHitsEB_map->find(det);
00349                                 if (aHit->second.energy() > seleXtalMinEnergy_)
00350                                   {
00351                                     pi0EBRecHitCollection->push_back(aHit->second);
00352                                     scXtals.push_back(det);
00353                                   }
00354 
00355                                 //EBDetId sel_rh = aHit->second.detid();
00356                                 //cout << "       RecHit Ok 20x20 matrix outside cluster : z,ieta,iphi "<<sel_rh.zside()<<" "<<sel_rh.ieta()<<" "<<sel_rh.iphi()<<endl;    
00357                                 //cout << "       RecHit Ok 20x20 matrix outside cluster : tower_ieta,tower_iphi "<<sel_rh.tower_ieta()<<" "<<sel_rh.tower_iphi()<<endl;   
00358                                 //cout << "       RecHit Ok 20x20 matrix outside cluster : iSM, ic "<<sel_rh.ism()<<" "<<sel_rh.ic()<<endl;
00359                                 
00360                               }
00361                         }
00362                       catch (...)
00363                         {
00364                         }
00365                     }
00366 
00367                 }
00368         
00369             }
00370           //        intbc++;
00371           
00372         }
00373 
00374 
00375 
00376       timerName = category + "::preparePi0RecHitsCollection";
00377       timers.pop_and_push(timerName);
00378 
00379 
00380       //Put selected information in the event
00381       //      if (npi0>0) iEvent.put( pi0EBRecHitCollection, pi0BarrelHits_);
00382       if ( pi0EBRecHitCollection->size() > seleNRHMax_ )
00383         {
00384           //      pi0EBRecHitCollection->clear();
00385           cout<< "   Max RH limit exceeded: "<<pi0EBRecHitCollection->size()<<" Max RH: "<<seleNRHMax_<<endl;
00386           iEvent.put( pi0EBDummyRecHitCollection, pi0BarrelHits_);
00387           cout<< "   EB RecHits # in Collection: "<<0<<endl;
00388         } else {
00389 
00390           cout<< "   EB RecHits # in Collection: "<<pi0EBRecHitCollection->size()<<endl;
00391           iEvent.put( pi0EBRecHitCollection, pi0BarrelHits_);
00392         }
00393   
00394       timerName = category + "::storePi0RecHitsCollection";
00395       timers.pop_and_push(timerName);
00396 
00397       timers.clear_stack();
00398       
00399       delete recHitsEB_map;
00400       delete topology_p;
00401 
00402 }


Member Data Documentation

std::string AlCaPi0RecHitsProducer::barrelHits_ [private]

Definition at line 58 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

int AlCaPi0RecHitsProducer::clusEtaSize_ [private]

Definition at line 65 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

int AlCaPi0RecHitsProducer::clusPhiSize_ [private]

Definition at line 66 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

double AlCaPi0RecHitsProducer::clusSeedThr_ [private]

Definition at line 64 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

std::string AlCaPi0RecHitsProducer::ecalHitsProducer_ [private]

Definition at line 57 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

int AlCaPi0RecHitsProducer::gammaCandEtaSize_ [private]

Definition at line 61 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

int AlCaPi0RecHitsProducer::gammaCandPhiSize_ [private]

Definition at line 62 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

std::string AlCaPi0RecHitsProducer::pi0BarrelHits_ [private]

Definition at line 59 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

std::map<DetId, EcalRecHit>* AlCaPi0RecHitsProducer::recHitsEB_map [private]

Definition at line 76 of file AlCaPi0RecHitsProducer.h.

Referenced by produce().

double AlCaPi0RecHitsProducer::seleMinvMaxPi0_ [private]

Definition at line 71 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

double AlCaPi0RecHitsProducer::seleMinvMinPi0_ [private]

Definition at line 72 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

int AlCaPi0RecHitsProducer::seleNRHMax_ [private]

Definition at line 74 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

double AlCaPi0RecHitsProducer::selePtGammaOne_ [private]

Definition at line 68 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

double AlCaPi0RecHitsProducer::selePtGammaTwo_ [private]

Definition at line 69 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

double AlCaPi0RecHitsProducer::selePtPi0_ [private]

Definition at line 70 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().

double AlCaPi0RecHitsProducer::seleXtalMinEnergy_ [private]

Definition at line 73 of file AlCaPi0RecHitsProducer.h.

Referenced by AlCaPi0RecHitsProducer(), and produce().


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