CMS 3D CMS Logo

Pi0FixedMassWindowCalibration Class Reference

#include <Calibration/EcalCalibAlgos/interface/Pi0FixedMassWindowCalibration.h>

Inheritance diagram for Pi0FixedMassWindowCalibration:

edm::ESProducerLooper edm::ESProducer edm::EventSetupRecordIntervalFinder edm::EDLooper edm::ESProxyFactoryProducer edm::eventsetup::DataProxyProvider

List of all members.

Public Member Functions

virtual void beginOfJob (const edm::EventSetup &)
 Called at beginning of job.
virtual Status duringLoop (const edm::Event &, const edm::EventSetup &)
 Called at each event.
virtual void endOfJob ()
 Called at end of job.
virtual Status endOfLoop (const edm::EventSetup &, unsigned int iLoop)
 Called at end of loop.
 Pi0FixedMassWindowCalibration (const edm::ParameterSet &iConfig)
 Constructor.
virtual void produce (edm::Event &, const edm::EventSetup &)
 Dummy implementation (job done in duringLoop).
virtual void startingNewLoop (unsigned int iLoop)
 Called at beginning of loop.
 ~Pi0FixedMassWindowCalibration ()
 Destructor.

Private Attributes

std::vector< DetIdbarrelCells
std::string barrelClusterCollection_
std::string barrelClusterShapeAssociation_
std::string barrelHitCollection_
std::string barrelHitProducer_
std::string barrelHits_
std::string clustershapecollectionEB_
std::string ecalHitsProducer_
const EcalRecHitCollectionecalRecHitBarrelCollection
IslandClusterAlgoisland_p
double mwxtals [85][360][2]
int nevent
double newCalibs_barl [85][360][2]
double oldCalibs_barl [85][360][2]
PositionCalc posCalculator_
const EcalRecHitCollectionrecalibEcalRecHitCollection
std::map< DetId, EcalRecHit > * recHitsEB_map
double selePi0DetaBelt_
double selePi0DRBelt_
double selePi0EtBeltIsoRatioMax_
double selePi0MinvMeanFixed_
double selePi0MinvSigmaFixed_
double selePi0PtGammaOneMin_
double selePi0PtGammaTwoMin_
double selePi0PtPi0Min_
double selePi0S4S9GammaOneMin_
double selePi0S4S9GammaTwoMin_
double selePi0S9S25GammaOneMin_
double selePi0S9S25GammaTwoMin_
ClusterShapeAlgo shapeAlgo_
TFile * theFile
unsigned int theMaxLoops
edm::ParameterSet theParameterSet
IslandClusterAlgo::VerbosityLevel verbosity
double wxtals [85][360][2]


Detailed Description

Definition at line 58 of file Pi0FixedMassWindowCalibration.h.


Constructor & Destructor Documentation

Pi0FixedMassWindowCalibration::Pi0FixedMassWindowCalibration ( const edm::ParameterSet iConfig  ) 

Constructor.

Definition at line 34 of file Pi0FixedMassWindowCalibration.cc.

References barrelClusterCollection_, barrelClusterShapeAssociation_, clustershapecollectionEB_, GenMuonPlsPt100GeV_cfg::cout, lat::endl(), edm::ParameterSet::getParameter(), island_p, IslandClusterAlgo::pDEBUG, IslandClusterAlgo::pERROR, IslandClusterAlgo::pINFO, posCalculator_, IslandClusterAlgo::pWARNING, selePi0DetaBelt_, selePi0DRBelt_, selePi0EtBeltIsoRatioMax_, selePi0MinvMeanFixed_, selePi0MinvSigmaFixed_, selePi0PtGammaOneMin_, selePi0PtGammaTwoMin_, selePi0PtPi0Min_, selePi0S4S9GammaOneMin_, selePi0S4S9GammaTwoMin_, selePi0S9S25GammaOneMin_, selePi0S9S25GammaTwoMin_, shapeAlgo_, theParameterSet, and verbosity.

00034                                                                                            :
00035   theMaxLoops( iConfig.getUntrackedParameter<unsigned int>("maxLoops",0) ),
00036   ecalHitsProducer_( iConfig.getParameter< std::string > ("ecalRecHitsProducer") ),
00037   barrelHits_( iConfig.getParameter< std::string > ("barrelHitCollection") )
00038 {
00039 
00040   std::cout << "[Pi0FixedMassWindowCalibration] Constructor "<<std::endl;
00041   // The verbosity level
00042   std::string verbosityString = iConfig.getParameter<std::string>("VerbosityLevel");
00043   if      (verbosityString == "DEBUG")   verbosity = IslandClusterAlgo::pDEBUG;
00044   else if (verbosityString == "WARNING") verbosity = IslandClusterAlgo::pWARNING;
00045   else if (verbosityString == "INFO")    verbosity = IslandClusterAlgo::pINFO;
00046   else                                   verbosity = IslandClusterAlgo::pERROR;
00047 
00048   // The names of the produced cluster collections
00049   barrelClusterCollection_  = iConfig.getParameter<std::string>("barrelClusterCollection");
00050 
00051   // Island algorithm parameters
00052   double barrelSeedThreshold = iConfig.getParameter<double>("IslandBarrelSeedThr");
00053   double endcapSeedThreshold = iConfig.getParameter<double>("IslandEndcapSeedThr");
00054 
00055 
00056   // Selection algorithm parameters
00057   selePi0PtGammaOneMin_ = iConfig.getParameter<double>("selePi0PtGammaOneMin");
00058   selePi0PtGammaTwoMin_ = iConfig.getParameter<double>("selePi0PtGammaTwoMin");
00059 
00060   selePi0DRBelt_ = iConfig.getParameter<double>("selePi0DRBelt");
00061   selePi0DetaBelt_ = iConfig.getParameter<double>("selePi0DetaBelt");
00062 
00063   selePi0PtPi0Min_ = iConfig.getParameter<double>("selePi0PtPi0Min");
00064 
00065   selePi0S4S9GammaOneMin_ = iConfig.getParameter<double>("selePi0S4S9GammaOneMin");
00066   selePi0S4S9GammaTwoMin_ = iConfig.getParameter<double>("selePi0S4S9GammaTwoMin");
00067   selePi0S9S25GammaOneMin_ = iConfig.getParameter<double>("selePi0S9S25GammaOneMin");
00068   selePi0S9S25GammaTwoMin_ = iConfig.getParameter<double>("selePi0S9S25GammaTwoMin");
00069  
00070   selePi0EtBeltIsoRatioMax_ = iConfig.getParameter<double>("selePi0EtBeltIsoRatioMax");
00071 
00072   selePi0MinvMeanFixed_ = iConfig.getParameter<double>("selePi0MinvMeanFixed");
00073   selePi0MinvSigmaFixed_ = iConfig.getParameter<double>("selePi0MinvSigmaFixed");
00074 
00075 
00076 
00077   // Parameters for the position calculation:
00078   std::map<std::string,double> providedParameters;
00079   providedParameters.insert(std::make_pair("LogWeighted",iConfig.getParameter<bool>("posCalc_logweight")));
00080   providedParameters.insert(std::make_pair("T0_barl",iConfig.getParameter<double>("posCalc_t0_barl")));
00081   providedParameters.insert(std::make_pair("T0_endc",iConfig.getParameter<double>("posCalc_t0_endc")));
00082   providedParameters.insert(std::make_pair("T0_endcPresh",iConfig.getParameter<double>("posCalc_t0_endcPresh")));
00083   providedParameters.insert(std::make_pair("W0",iConfig.getParameter<double>("posCalc_w0")));
00084   providedParameters.insert(std::make_pair("X0",iConfig.getParameter<double>("posCalc_x0")));
00085   posCalculator_ = PositionCalc(providedParameters);
00086   //  shapeAlgo_ = ClusterShapeAlgo(posCalculator_);
00087   shapeAlgo_ = ClusterShapeAlgo(providedParameters);
00088   clustershapecollectionEB_ = iConfig.getParameter<std::string>("clustershapecollectionEB");
00089 
00090   //AssociationMap
00091   barrelClusterShapeAssociation_ = iConfig.getParameter<std::string>("barrelShapeAssociation");
00092 
00093   island_p = new IslandClusterAlgo(barrelSeedThreshold, endcapSeedThreshold, posCalculator_,verbosity);
00094 
00095   theParameterSet=iConfig;
00096 
00097 
00098 }

Pi0FixedMassWindowCalibration::~Pi0FixedMassWindowCalibration (  ) 

Destructor.

Definition at line 104 of file Pi0FixedMassWindowCalibration.cc.

00105 {
00106 
00107 }


Member Function Documentation

void Pi0FixedMassWindowCalibration::beginOfJob ( const edm::EventSetup iSetup  )  [virtual]

Called at beginning of job.

Reimplemented from edm::EDLooper.

Definition at line 112 of file Pi0FixedMassWindowCalibration.cc.

References funct::abs(), barrelCells, calib, TestMuL1L2Filter_cff::cerr, GenMuonPlsPt100GeV_cfg::cout, DetId::Ecal, EcalBarrel, EcalCondObjectContainer< T >::end(), lat::endl(), exception, EcalCondObjectContainer< T >::find(), edm::EventSetup::get(), EcalCondObjectContainer< T >::getMap(), CaloGeometry::getValidDetIds(), EBDetId::ieta(), EBDetId::iphi(), mwxtals, newCalibs_barl, oldCalibs_barl, edm::ESHandle< T >::product(), DetId::rawId(), EcalCondObjectContainer< T >::size(), wxtals, and EBDetId::zside().

00113 {
00114 
00115   std::cout << "[Pi0FixedMassWindowCalibration] beginOfJob "<<std::endl;
00116 
00117   // initialize arrays
00118 
00119   for (int sign=0; sign<2; sign++) {
00120     for (int ieta=0; ieta<85; ieta++) {
00121       for (int iphi=0; iphi<360; iphi++) {
00122         oldCalibs_barl[ieta][iphi][sign]=0.;
00123         newCalibs_barl[ieta][iphi][sign]=0.;
00124         wxtals[ieta][iphi][sign]=0.;
00125         mwxtals[ieta][iphi][sign]=0.;
00126       }
00127     }
00128   }
00129 
00130   // get initial constants out of DB
00131 
00132   edm::ESHandle<EcalIntercalibConstants> pIcal;
00133   EcalIntercalibConstantMap imap;
00134 
00135   try {
00136     iSetup.get<EcalIntercalibConstantsRcd>().get(pIcal);
00137     std::cout << "Taken EcalIntercalibConstants" << std::endl;
00138     imap = pIcal.product()->getMap();
00139     std::cout << "imap.size() = " << imap.size() << std::endl;
00140   } catch ( std::exception& ex ) {     
00141     std::cerr << "Error! can't get EcalIntercalibConstants " << std::endl;
00142   } 
00143 
00144   // get the ecal geometry:
00145   edm::ESHandle<CaloGeometry> geoHandle;
00146   iSetup.get<CaloGeometryRecord>().get(geoHandle);
00147   const CaloGeometry& geometry = *geoHandle;
00148   //const CaloSubdetectorGeometry *barrelGeometry = geometry.getSubdetectorGeometry(DetId::Ecal, EcalBarrel);
00149 
00150   // loop over all barrel crystals
00151   barrelCells = geometry.getValidDetIds(DetId::Ecal, EcalBarrel);
00152   std::vector<DetId>::const_iterator barrelIt;
00153   for (barrelIt=barrelCells.begin(); barrelIt!=barrelCells.end(); barrelIt++) {
00154     EBDetId eb(*barrelIt);
00155 
00156     // get the initial calibration constants
00157     EcalIntercalibConstantMap::const_iterator itcalib = imap.find(eb.rawId());
00158     if ( itcalib == imap.end() ) {
00159             // FIXME -- throw error
00160     }
00161     EcalIntercalibConstant calib = (*itcalib);
00162     int sign = eb.zside()>0 ? 1 : 0;
00163     oldCalibs_barl[abs(eb.ieta())-1][eb.iphi()-1][sign] = calib;
00164     if (eb.iphi()==1) std::cout << "Read old constant for crystal "
00165                            << " (" << eb.ieta() << "," << eb.iphi()
00166                            << ") : " << calib << std::endl;
00167 
00168 
00169 
00170   }
00171 
00172 }

edm::EDLooper::Status Pi0FixedMassWindowCalibration::duringLoop ( const edm::Event event,
const edm::EventSetup setup 
) [virtual]

Called at each event.

Implements edm::EDLooper.

Definition at line 274 of file Pi0FixedMassWindowCalibration.cc.

References funct::abs(), IslandClusterAlgo::barrel, barrelHits_, edm::SortedCollection< T, SORT >::begin(), ClusterShapeAlgo::Calculate(), funct::cos(), GenMuonPlsPt100GeV_cfg::cout, DetId::Ecal, EcalBarrel, ecalHitsProducer_, EcalPreshower, ecalRecHitBarrelCollection, edm::SortedCollection< T, SORT >::end(), lat::endl(), funct::exp(), edm::EventSetup::get(), i, EBDetId::ieta(), int, EBDetId::iphi(), island_p, j, k, edm::EDLooper::kContinue, kk, funct::log(), IslandClusterAlgo::makeClusters(), mwxtals, nevent, oldCalibs_barl, recalibEcalRecHitCollection, recHitsEB_map, selePi0DetaBelt_, selePi0DRBelt_, selePi0EtBeltIsoRatioMax_, selePi0MinvMeanFixed_, selePi0MinvSigmaFixed_, selePi0PtGammaOneMin_, selePi0PtPi0Min_, selePi0S4S9GammaOneMin_, selePi0S4S9GammaTwoMin_, selePi0S9S25GammaOneMin_, selePi0S9S25GammaTwoMin_, shapeAlgo_, funct::sin(), edm::SortedCollection< T, SORT >::size(), funct::sqrt(), std, funct::tan(), wxtals, and EBDetId::zside().

00276 {
00277   using namespace edm;
00278   using namespace std;
00279 
00280   nevent++;
00281 
00282   if ((nevent<100 && nevent%10==0) 
00283       ||(nevent<1000 && nevent%100==0) 
00284       ||(nevent<10000 && nevent%100==0) 
00285       ||(nevent<100000 && nevent%1000==0) 
00286       ||(nevent<10000000 && nevent%1000==0))
00287     std::cout << "[Pi0FixedMassWindowCalibration] Events processed: "<<nevent<<std::endl;
00288 
00289   recHitsEB_map = new std::map<DetId, EcalRecHit>();
00290 
00291   EcalRecHitCollection* recalibEcalRecHitCollection( new EcalRecHitCollection );
00292 
00293   int nRecHitsEB=0;
00294   Handle<EcalRecHitCollection> pEcalRecHitBarrelCollection;
00295   event.getByLabel(ecalHitsProducer_, barrelHits_, pEcalRecHitBarrelCollection);
00296   const EcalRecHitCollection* ecalRecHitBarrelCollection = pEcalRecHitBarrelCollection.product();
00297   cout << " ECAL Barrel RecHits # "<< ecalRecHitBarrelCollection->size() <<endl;
00298   for(EcalRecHitCollection::const_iterator aRecHitEB = ecalRecHitBarrelCollection->begin(); aRecHitEB != ecalRecHitBarrelCollection->end(); aRecHitEB++) {
00299     //cout << " ECAL Barrel RecHit #,E,time,det,subdetid: "<<nRecHitsEB<<" "<<aRecHitEB->energy()<<" "<<aRecHitEB->time()<<" "<<aRecHitEB->detid().det()<<" "<<aRecHitEB->detid().subdetId()<<endl;
00300 
00301 
00302     EBDetId ebrhdetid = aRecHitEB->detid();
00303     //cout << " EBDETID: z,ieta,iphi "<<ebrhdetid.zside()<<" "<<ebrhdetid.ieta()<<" "<<ebrhdetid.iphi()<<endl;
00304     //cout << " EBDETID: tower_ieta,tower_iphi "<<ebrhdetid.tower_ieta()<<" "<<ebrhdetid.tower_iphi()<<endl;
00305     //cout << " EBDETID: iSM, ic "<<ebrhdetid.ism()<<" "<<ebrhdetid.ic()<<endl;
00306 
00307     int sign = ebrhdetid.zside()>0 ? 1 : 0;
00308     EcalRecHit aHit(aRecHitEB->id(),aRecHitEB->energy()*oldCalibs_barl[abs(ebrhdetid.ieta())-1][ebrhdetid.iphi()-1][sign],aRecHitEB->time());
00309     recalibEcalRecHitCollection->push_back(aHit);
00310 
00311     nRecHitsEB++;
00312 
00313   }
00314 
00315   //  cout<<" Recalib size: "<<recalibEcalRecHitCollection->size()<<endl;
00316   int irecalib=0;
00317   for(EcalRecHitCollection::const_iterator aRecHitEB = recalibEcalRecHitCollection->begin(); aRecHitEB != recalibEcalRecHitCollection->end(); aRecHitEB++) {
00318     //cout << " [recalibrated] ECAL Barrel RecHit #,E,time,det,subdetid: "<<irecalib<<" "<<aRecHitEB->energy()<<" "<<aRecHitEB->time()<<" "<<aRecHitEB->detid().det()<<" "<<aRecHitEB->detid().subdetId()<<endl;
00319 
00320     //    EBDetId ebrhdetid = aRecHitEB->detid();
00321     //cout << " [recalibrated] EBDETID: z,ieta,iphi "<<ebrhdetid.zside()<<" "<<ebrhdetid.ieta()<<" "<<ebrhdetid.iphi()<<endl;
00322     //cout << " [recalibrated] EBDETID: tower_ieta,tower_iphi "<<ebrhdetid.tower_ieta()<<" "<<ebrhdetid.tower_iphi()<<endl;
00323     //cout << " [recalibrated] EBDETID: iSM, ic "<<ebrhdetid.ism()<<" "<<ebrhdetid.ic()<<endl;
00324 
00325     std::pair<DetId, EcalRecHit> map_entry(aRecHitEB->id(), *aRecHitEB);
00326     recHitsEB_map->insert(map_entry);
00327 
00328 
00329     irecalib++;
00330 
00331   }
00332 
00333   
00334   // get the geometry and topology from the event setup:
00335   edm::ESHandle<CaloGeometry> geoHandle;
00336   setup.get<CaloGeometryRecord>().get(geoHandle);
00337 
00338   const CaloSubdetectorGeometry *geometry_p;
00339   CaloSubdetectorTopology *topology_p;
00340 
00341   std::string clustershapetag;
00342   geometry_p = geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalBarrel);
00343   topology_p = new EcalBarrelTopology(geoHandle);
00344 
00345   const CaloSubdetectorGeometry *geometryES_p;
00346   geometryES_p = geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalPreshower);
00347 
00348   /*
00349   reco::BasicClusterCollection clusters;
00350   clusters = island_p->makeClusters(ecalRecHitBarrelCollection,geometry_p,topology_p,geometryES_p,IslandClusterAlgo::barrel);
00351   
00352   //Create associated ClusterShape objects.
00353   std::vector <reco::ClusterShape> ClusVec;
00354   for (int erg=0;erg<int(clusters.size());++erg){
00355     reco::ClusterShape TestShape = shapeAlgo_.Calculate(clusters[erg],ecalRecHitBarrelCollection,geometry_p,topology_p);
00356     ClusVec.push_back(TestShape);
00357   }
00358 
00359   //Put clustershapes in event, but retain a Handle on them.
00360   std::auto_ptr<reco::ClusterShapeCollection> clustersshapes_p(new reco::ClusterShapeCollection);
00361   clustersshapes_p->assign(ClusVec.begin(), ClusVec.end());
00362 
00363   cout<<"[Pi0Calibration] Basic Cluster collection size: "<<clusters.size()<<endl;
00364   cout<<"[Pi0Calibration] Basic Cluster Shape Collection size: "<<clustersshapes_p->size()<<endl;
00365 
00366   int iClus=0;
00367   for(reco::BasicClusterCollection::const_iterator aClus = clusters.begin(); aClus != clusters.end(); aClus++) {
00368     cout<<" CLUSTER : #,NHits,e,et,eta,phi,e2x2,e3x3,e5x5: "<<iClus<<" "<<aClus->getHitsByDetId().size()<<" "<<aClus->energy()<<" "<<aClus->energy()*sin(aClus->position().theta())<<" "<<aClus->position().eta()<<" "<<aClus->position().phi()<<" "<<(*clustersshapes_p)[iClus].e2x2()<<" "<<(*clustersshapes_p)[iClus].e3x3()<<" "<<(*clustersshapes_p)[iClus].e5x5()<<endl; 
00369     iClus++;
00370   }
00371   */
00372 
00373   // recalibrated clusters
00374   reco::BasicClusterCollection clusters_recalib;
00375   clusters_recalib = island_p->makeClusters(recalibEcalRecHitCollection,geometry_p,topology_p,geometryES_p,IslandClusterAlgo::barrel);
00376   
00377   //Create associated ClusterShape objects.
00378   std::vector <reco::ClusterShape> ClusVec_recalib;
00379   for (int erg=0;erg<int(clusters_recalib.size());++erg){
00380     reco::ClusterShape TestShape_recalib = shapeAlgo_.Calculate(clusters_recalib[erg],recalibEcalRecHitCollection,geometry_p,topology_p);
00381     ClusVec_recalib.push_back(TestShape_recalib);
00382   }
00383 
00384   //Put clustershapes in event, but retain a Handle on them.
00385   std::auto_ptr<reco::ClusterShapeCollection> clustersshapes_p_recalib(new reco::ClusterShapeCollection);
00386   clustersshapes_p_recalib->assign(ClusVec_recalib.begin(), ClusVec_recalib.end());
00387 
00388   cout<<"[Pi0FixedMassWindowCalibration][recalibration] Basic Cluster collection size: "<<clusters_recalib.size()<<endl;
00389   cout<<"[Pi0FixedMassWindowCalibration][recalibration] Basic Cluster Shape Collection size: "<<clustersshapes_p_recalib->size()<<endl;
00390 
00391 
00392   // pizero selection
00393 
00394   // Get ECAL Barrel Island Basic Clusters collection
00395   // ECAL Barrel Island Basic Clusters 
00396   static const int MAXBCEB = 200;
00397   static const int MAXBCEBRH = 200;
00398   int nIslandBCEB;
00399   float eIslandBCEB[MAXBCEB];
00400   float etIslandBCEB[MAXBCEB];
00401   float etaIslandBCEB[MAXBCEB];
00402   float phiIslandBCEB[MAXBCEB];
00403   int ietaIslandBCEB[MAXBCEB];
00404   int iphiIslandBCEB[MAXBCEB];
00405   float e2x2IslandBCEB[MAXBCEB];
00406   float e3x3IslandBCEB[MAXBCEB];
00407   float e5x5IslandBCEB[MAXBCEB];
00408   // indexes to the RecHits assiciated with 
00409   // ECAL Barrel Island Basic Clusters
00410   int nIslandBCEBRecHits[MAXBCEB];
00411   //  int indexIslandBCEBRecHits[MAXBCEB][MAXBCEBRH];
00412   int ietaIslandBCEBRecHits[MAXBCEB][MAXBCEBRH];
00413   int iphiIslandBCEBRecHits[MAXBCEB][MAXBCEBRH];
00414   int zsideIslandBCEBRecHits[MAXBCEB][MAXBCEBRH];
00415   float eIslandBCEBRecHits[MAXBCEB][MAXBCEBRH];
00416 
00417   nIslandBCEB=0;
00418   for(int i=0; i<MAXBCEB; i++){
00419     eIslandBCEB[i] = 0;
00420     etIslandBCEB[i] = 0;
00421     etaIslandBCEB[i] = 0;
00422     phiIslandBCEB[i] = 0;
00423     ietaIslandBCEB[i] = 0;
00424     iphiIslandBCEB[i] = 0;
00425     e2x2IslandBCEB[i] = 0;
00426     e3x3IslandBCEB[i] = 0;
00427     e5x5IslandBCEB[i] = 0;
00428      nIslandBCEBRecHits[i] = 0;
00429      for(int j=0;j<MAXBCEBRH;j++){
00430        //       indexIslandBCEBRecHits[i][j] = 0;
00431        ietaIslandBCEBRecHits[i][j] = 0;
00432        iphiIslandBCEBRecHits[i][j] = 0;
00433        zsideIslandBCEBRecHits[i][j] = 0;
00434        eIslandBCEBRecHits[i][j] = 0;
00435      }
00436   }
00437 
00438 
00439   int iClus_recalib=0;
00440   for(reco::BasicClusterCollection::const_iterator aClus = clusters_recalib.begin(); aClus != clusters_recalib.end(); aClus++) {
00441     cout<<" CLUSTER [recalibration] : #,NHits,e,et,eta,phi,e2x2,e3x3,e5x5: "<<iClus_recalib<<" "<<aClus->getHitsByDetId().size()<<" "<<aClus->energy()<<" "<<aClus->energy()*sin(aClus->position().theta())<<" "<<aClus->position().eta()<<" "<<aClus->position().phi()<<" "<<(*clustersshapes_p_recalib)[iClus_recalib].e2x2()<<" "<<(*clustersshapes_p_recalib)[iClus_recalib].e3x3()<<" "<<(*clustersshapes_p_recalib)[iClus_recalib].e5x5()<<endl; 
00442 
00443     eIslandBCEB[nIslandBCEB] = aClus->energy();
00444     etIslandBCEB[nIslandBCEB] = aClus->energy()*sin(aClus->position().theta());
00445     etaIslandBCEB[nIslandBCEB] = aClus->position().eta();
00446     phiIslandBCEB[nIslandBCEB] = aClus->position().phi();
00447     
00448     e2x2IslandBCEB[nIslandBCEB] = (*clustersshapes_p_recalib)[nIslandBCEB].e2x2();    
00449     e3x3IslandBCEB[nIslandBCEB] = (*clustersshapes_p_recalib)[nIslandBCEB].e3x3();    
00450     e5x5IslandBCEB[nIslandBCEB] = (*clustersshapes_p_recalib)[nIslandBCEB].e5x5();    
00451 
00452     nIslandBCEBRecHits[nIslandBCEB] = aClus->getHitsByDetId().size();
00453 
00454     std::vector<DetId> hits = aClus->getHitsByDetId();
00455     std::vector<DetId>::iterator hit;
00456     std::map<DetId, EcalRecHit>::iterator aHit;
00457 
00458     int irhcount=0;
00459     for(hit = hits.begin(); hit != hits.end(); hit++)
00460       {
00461         // need to get hit by DetID in order to get energy
00462         aHit = recHitsEB_map->find(*hit);
00463         //cout << "       RecHit #: "<<irhcount<<" from Basic Cluster with Energy: "<<aHit->second.energy()<<endl; 
00464 
00465         EBDetId sel_rh = aHit->second.detid();
00466         //cout << "       RecHit: z,ieta,iphi "<<sel_rh.zside()<<" "<<sel_rh.ieta()<<" "<<sel_rh.iphi()<<endl;
00467         //cout << "       RecHit: tower_ieta,tower_iphi "<<sel_rh.tower_ieta()<<" "<<sel_rh.tower_iphi()<<endl;
00468         //cout << "       RecHit: iSM, ic "<<sel_rh.ism()<<" "<<sel_rh.ic()<<endl;
00469 
00470         ietaIslandBCEBRecHits[nIslandBCEB][irhcount] = sel_rh.ieta();
00471         iphiIslandBCEBRecHits[nIslandBCEB][irhcount] = sel_rh.iphi();
00472         zsideIslandBCEBRecHits[nIslandBCEB][irhcount] = sel_rh.zside();
00473         eIslandBCEBRecHits[nIslandBCEB][irhcount] = aHit->second.energy();
00474 
00475         irhcount++;
00476       }
00477     nIslandBCEB++;
00478     iClus_recalib++;
00479 
00480 
00481   }
00482 
00483   // Selection, based on ECAL Barrel Basic Clusters
00484 
00485   if (nIslandBCEB > 1) 
00486     {
00487       for(int i=0 ; i<nIslandBCEB ; i++)
00488         {
00489           for(int j=i+1 ; j<nIslandBCEB ; j++)
00490             {
00491 
00492               if( etIslandBCEB[i]>selePi0PtGammaOneMin_ && etIslandBCEB[j]>selePi0PtGammaOneMin_) 
00493                 {
00494                 
00495                   float theta_0 = 2. * atan(exp(-etaIslandBCEB[i]));
00496                   float theta_1 = 2. * atan(exp(-etaIslandBCEB[j]));
00497                 
00498                   float p0x = eIslandBCEB[i] * sin(theta_0) * cos(phiIslandBCEB[i]);
00499                   float p1x = eIslandBCEB[j] * sin(theta_1) * cos(phiIslandBCEB[j]);
00500                 
00501                   float p0y = eIslandBCEB[i] * sin(theta_0) * sin(phiIslandBCEB[i]);
00502                   float p1y = eIslandBCEB[j] * sin(theta_1) * sin(phiIslandBCEB[j]);
00503                 
00504                   float p0z = eIslandBCEB[i] * cos(theta_0);
00505                   float p1z = eIslandBCEB[j] * cos(theta_1);
00506 
00507                   float pi0_px = p0x + p1x;
00508                   float pi0_py = p0y + p1y;
00509                   float pi0_pz = p0z + p1z;
00510                 
00511                   float pi0_ptot = sqrt (pi0_px*pi0_px + pi0_py*pi0_py + pi0_pz*pi0_pz); 
00512 
00513                   float pi0_theta = acos(pi0_pz/pi0_ptot);
00514                   float pi0_eta = -log(tan(pi0_theta/2));
00515                   float pi0_phi = atan(pi0_py/pi0_px);
00516                   //cout << " pi0_theta, pi0_eta, pi0_phi "<<pi0_theta<<" "<<pi0_eta<<" "<<pi0_phi<<endl;
00517 
00518                   // belt isolation
00519                 
00520                   float et_belt=0;
00521                   for(Int_t k=0 ; k<nIslandBCEB ; k++)
00522                     {
00523                       if( (k != i) && (k != j) ) 
00524                         {
00525                           float dr_pi0_k = sqrt( (etaIslandBCEB[k]-pi0_eta)*(etaIslandBCEB[k]-pi0_eta) + (phiIslandBCEB[k]-pi0_phi)*(phiIslandBCEB[k]-pi0_phi) );
00526                           float deta_pi0_k = fabs(etaIslandBCEB[k]-pi0_eta);
00527                           if ( (dr_pi0_k<selePi0DRBelt_) && (deta_pi0_k<selePi0DetaBelt_) ) et_belt = et_belt + etIslandBCEB[k];
00528                         }
00529                     }
00530 
00531                 
00532                   float pt_pi0 = sqrt( (p0x+p1x)*(p0x+p1x) + (p0y+p1y)*(p0y+p1y));
00533                   //float dr_pi0 = sqrt ( (etaIslandBCEB[i]-etaIslandBCEB[j])*(etaIslandBCEB[i]-etaIslandBCEB[j]) + (phiIslandBCEB[i]-phiIslandBCEB[j])*(phiIslandBCEB[i]-phiIslandBCEB[j]) );
00534 
00535                   //cout <<" pi0 pt,dr:  "<<pt_pi0<<" "<<dr_pi0<<endl; 
00536                   if (pt_pi0 > selePi0PtPi0Min_) 
00537                     {
00538                       float m_inv = sqrt ( (eIslandBCEB[i] + eIslandBCEB[j])*(eIslandBCEB[i] + eIslandBCEB[j]) - (p0x+p1x)*(p0x+p1x) - (p0y+p1y)*(p0y+p1y) - (p0z+p1z)*(p0z+p1z) );  
00539                       cout <<" pi0 (pt>2.5 GeV) m_inv = "<<m_inv<<endl;  
00540 
00541                       float s4s9_1 = e2x2IslandBCEB[i]/e3x3IslandBCEB[i];
00542                       float s4s9_2 = e2x2IslandBCEB[j]/e3x3IslandBCEB[j];
00543 
00544                       float s9s25_1 = e3x3IslandBCEB[i]/e5x5IslandBCEB[i];
00545                       float s9s25_2 = e3x3IslandBCEB[j]/e5x5IslandBCEB[j];
00546 
00547                       //float s9Esc_1 = e3x3IslandBCEB[i]/eIslandBCEB[i];
00548                       //float s9Esc_2 = e3x3IslandBCEB[j]/eIslandBCEB[j];
00549 
00550                       if (s4s9_1>selePi0S4S9GammaOneMin_ && s4s9_2>selePi0S4S9GammaTwoMin_ && s9s25_1>selePi0S9S25GammaOneMin_ && s9s25_2>selePi0S9S25GammaTwoMin_ && (et_belt/pt_pi0 < selePi0EtBeltIsoRatioMax_) ) 
00551                         {
00552                           //good pizero candidate
00553                           if ( m_inv>(selePi0MinvMeanFixed_ - 2*selePi0MinvSigmaFixed_) && m_inv<(selePi0MinvMeanFixed_ + 2*selePi0MinvSigmaFixed_) ) 
00554                             {
00555                               //fill wxtals and mwxtals weights
00556                               cout<<" Pi0 Good candidate : minv = "<<m_inv<<endl;
00557                               for(int kk=0 ; kk<nIslandBCEBRecHits[i] ; kk++)
00558                                 {
00559                                   int ieta_xtal = ietaIslandBCEBRecHits[i][kk];
00560                                   int iphi_xtal = iphiIslandBCEBRecHits[i][kk];
00561                                   int sign = zsideIslandBCEBRecHits[i][kk]>0 ? 1 : 0;
00562                                   wxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] = wxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] + eIslandBCEBRecHits[i][kk]/e3x3IslandBCEB[i];
00563                                   mwxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] = mwxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] + (selePi0MinvMeanFixed_/m_inv)*(selePi0MinvMeanFixed_/m_inv)*(eIslandBCEBRecHits[i][kk]/e3x3IslandBCEB[i]);
00564                                   cout<< "[Pi0FixedMassWindowCalibration] eta, phi, sign, e, e3x3, wxtals and mwxtals: "<<ieta_xtal<<" "<<iphi_xtal<<" "<<sign<<" "<<eIslandBCEBRecHits[i][kk]<<" "<<e3x3IslandBCEB[i]<<" "<<wxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign]<<" "<<mwxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign]<<endl;
00565                                 }
00566 
00567                               for(int kk=0 ; kk<nIslandBCEBRecHits[j] ; kk++)
00568                                 {
00569                                   int ieta_xtal = ietaIslandBCEBRecHits[j][kk];
00570                                   int iphi_xtal = iphiIslandBCEBRecHits[j][kk];
00571                                   int sign = zsideIslandBCEBRecHits[j][kk]>0 ? 1 : 0;
00572                                   wxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] = wxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] + eIslandBCEBRecHits[j][kk]/e3x3IslandBCEB[j];
00573                                   mwxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] = mwxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] + (selePi0MinvMeanFixed_/m_inv)*(selePi0MinvMeanFixed_/m_inv)*(eIslandBCEBRecHits[j][kk]/e3x3IslandBCEB[j]);
00574                                   cout<< "[Pi0FixedMassWindowCalibration] eta, phi, sign, e, e3x3, wxtals and mwxtals: "<<ieta_xtal<<" "<<iphi_xtal<<" "<<sign<<" "<<eIslandBCEBRecHits[j][kk]<<" "<<e3x3IslandBCEB[j]<<" "<<wxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign]<<" "<<mwxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign]<<endl;
00575                                 }
00576 
00577                             }
00578                         }
00579                     }
00580                 }
00581                 
00582               
00583             } // End of the "j" loop over BCEB
00584         } // End of the "i" loop over BCEB
00585 
00586 
00587     } else {
00588       cout<< " Not enough ECAL Barrel Basic Clusters: "<<nIslandBCEB<<endl;
00589     }
00590 
00591 
00592 
00593 
00594 
00595 
00596   return kContinue;
00597 }

void Pi0FixedMassWindowCalibration::endOfJob (  )  [virtual]

Called at end of job.

Reimplemented from edm::EDLooper.

Definition at line 175 of file Pi0FixedMassWindowCalibration.cc.

References funct::abs(), barrelCells, GenMuonPlsPt100GeV_cfg::cout, EcalBarrel, lat::endl(), EBDetId::ieta(), EBDetId::iphi(), newCalibs_barl, Pi0CalibXMLwriter::writeLine(), and EBDetId::zside().

00176 {
00177   
00178   std::cout << "[Pi0FixedMassWindowCalibration] endOfJob"<<endl;
00179 
00180   // Write new calibration constants
00181 
00182   Pi0CalibXMLwriter barrelWriter(EcalBarrel,99);
00183 
00184   std::vector<DetId>::const_iterator barrelIt=barrelCells.begin();
00185   for (; barrelIt!=barrelCells.end(); barrelIt++) {
00186     EBDetId eb(*barrelIt);
00187     int ieta = eb.ieta();
00188     int iphi = eb.iphi();
00189     int sign = eb.zside()>0 ? 1 : 0;
00190     barrelWriter.writeLine(eb,newCalibs_barl[abs(ieta)-1][iphi-1][sign]);
00191   }
00192 
00193 }

edm::EDLooper::Status Pi0FixedMassWindowCalibration::endOfLoop ( const edm::EventSetup iSetup,
unsigned int  iLoop 
) [virtual]

Called at end of loop.

Implements edm::EDLooper.

Definition at line 217 of file Pi0FixedMassWindowCalibration.cc.

References funct::abs(), barrelCells, GenMuonPlsPt100GeV_cfg::cout, EcalBarrel, lat::endl(), EBDetId::ieta(), EBDetId::iphi(), edm::EDLooper::kContinue, edm::EDLooper::kStop, mwxtals, newCalibs_barl, oldCalibs_barl, theMaxLoops, Pi0CalibXMLwriter::writeLine(), wxtals, and EBDetId::zside().

00218 {
00219   std::cout << "[Pi0FixedMassWindowCalibration] Ending loop " << iLoop<<std::endl;
00220 
00221   for (int sign=0; sign<2; sign++) {
00222     for (int ieta=0; ieta<85; ieta++) {
00223       for (int iphi=0; iphi<360; iphi++) {
00224 
00225         if (wxtals[ieta][iphi][sign] == 0 ) 
00226           {
00227             newCalibs_barl[ieta][iphi][sign]=oldCalibs_barl[ieta][iphi][sign];
00228           } else {
00229             newCalibs_barl[ieta][iphi][sign]=oldCalibs_barl[ieta][iphi][sign]*(mwxtals[ieta][iphi][sign]/wxtals[ieta][iphi][sign]);
00230           }
00231         cout<< " New calibration constant: ieta iphi sign - old,mwxtals,wxtals,new: "<<ieta<<" "<<iphi<<" "<<sign<<" - "<<oldCalibs_barl[ieta][iphi][sign]<<" "<<mwxtals[ieta][iphi][sign]<<" "<<wxtals[ieta][iphi][sign]<<" "<<newCalibs_barl[ieta][iphi][sign]<<endl;
00232 
00233       }
00234     }
00235   }
00236     
00237   Pi0CalibXMLwriter barrelWriter(EcalBarrel,iLoop+1);
00238 
00239   std::vector<DetId>::const_iterator barrelIt=barrelCells.begin();
00240   for (; barrelIt!=barrelCells.end(); barrelIt++) {
00241     EBDetId eb(*barrelIt);
00242     int ieta = eb.ieta();
00243     int iphi = eb.iphi();
00244     int sign = eb.zside()>0 ? 1 : 0;
00245     barrelWriter.writeLine(eb,newCalibs_barl[abs(ieta)-1][iphi-1][sign]);
00246     if (iphi==1) {
00247       std::cout << "Calib constant for barrel crystal "
00248                 << " (" << ieta << "," << iphi << ") changed from "
00249                 << oldCalibs_barl[abs(ieta)-1][iphi-1][sign] << " to "
00250                 << newCalibs_barl[abs(ieta)-1][iphi-1][sign] << std::endl;
00251     }
00252   }
00253 
00254   // replace old calibration constants with new one
00255 
00256   for (int sign=0; sign<2; sign++) {
00257     for (int ieta=0; ieta<85; ieta++) {
00258       for (int iphi=0; iphi<360; iphi++) {
00259         oldCalibs_barl[ieta][iphi][sign]=newCalibs_barl[ieta][iphi][sign];
00260         newCalibs_barl[ieta][iphi][sign]=0;
00261       }
00262     }
00263   }
00264 
00265 
00266   if ( iLoop == theMaxLoops-1 || iLoop >= theMaxLoops ) return kStop;
00267   else return kContinue;
00268 }

virtual void Pi0FixedMassWindowCalibration::produce ( edm::Event ,
const edm::EventSetup  
) [inline, virtual]

Dummy implementation (job done in duringLoop).

Definition at line 70 of file Pi0FixedMassWindowCalibration.h.

00070 {};

void Pi0FixedMassWindowCalibration::startingNewLoop ( unsigned int  iLoop  )  [virtual]

Called at beginning of loop.

Implements edm::EDLooper.

Definition at line 197 of file Pi0FixedMassWindowCalibration.cc.

References GenMuonPlsPt100GeV_cfg::cout, lat::endl(), mwxtals, and wxtals.

00198 {
00199 
00200   for (int sign=0; sign<2; sign++) {
00201     for (int ieta=0; ieta<85; ieta++) {
00202       for (int iphi=0; iphi<360; iphi++) {
00203         wxtals[ieta][iphi][sign]=0.;
00204         mwxtals[ieta][iphi][sign]=0.;
00205       }
00206     }
00207   }
00208   std::cout << "[Pi0FixedMassWindowCalibration] Starting loop number " << iLoop<<std::endl;
00209 
00210 }


Member Data Documentation

std::vector<DetId> Pi0FixedMassWindowCalibration::barrelCells [private]

Definition at line 131 of file Pi0FixedMassWindowCalibration.h.

Referenced by beginOfJob(), endOfJob(), and endOfLoop().

std::string Pi0FixedMassWindowCalibration::barrelClusterCollection_ [private]

Definition at line 102 of file Pi0FixedMassWindowCalibration.h.

Referenced by Pi0FixedMassWindowCalibration().

std::string Pi0FixedMassWindowCalibration::barrelClusterShapeAssociation_ [private]

Definition at line 104 of file Pi0FixedMassWindowCalibration.h.

Referenced by Pi0FixedMassWindowCalibration().

std::string Pi0FixedMassWindowCalibration::barrelHitCollection_ [private]

Definition at line 101 of file Pi0FixedMassWindowCalibration.h.

std::string Pi0FixedMassWindowCalibration::barrelHitProducer_ [private]

Definition at line 100 of file Pi0FixedMassWindowCalibration.h.

std::string Pi0FixedMassWindowCalibration::barrelHits_ [private]

Definition at line 95 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop().

std::string Pi0FixedMassWindowCalibration::clustershapecollectionEB_ [private]

Definition at line 103 of file Pi0FixedMassWindowCalibration.h.

Referenced by Pi0FixedMassWindowCalibration().

std::string Pi0FixedMassWindowCalibration::ecalHitsProducer_ [private]

Definition at line 94 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop().

const EcalRecHitCollection* Pi0FixedMassWindowCalibration::ecalRecHitBarrelCollection [private]

Definition at line 149 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop().

IslandClusterAlgo* Pi0FixedMassWindowCalibration::island_p [private]

Definition at line 108 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop(), and Pi0FixedMassWindowCalibration().

double Pi0FixedMassWindowCalibration::mwxtals[85][360][2] [private]

Definition at line 137 of file Pi0FixedMassWindowCalibration.h.

Referenced by beginOfJob(), duringLoop(), endOfLoop(), and startingNewLoop().

int Pi0FixedMassWindowCalibration::nevent [private]

Definition at line 91 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop().

double Pi0FixedMassWindowCalibration::newCalibs_barl[85][360][2] [private]

Definition at line 135 of file Pi0FixedMassWindowCalibration.h.

Referenced by beginOfJob(), endOfJob(), and endOfLoop().

double Pi0FixedMassWindowCalibration::oldCalibs_barl[85][360][2] [private]

Definition at line 134 of file Pi0FixedMassWindowCalibration.h.

Referenced by beginOfJob(), duringLoop(), and endOfLoop().

PositionCalc Pi0FixedMassWindowCalibration::posCalculator_ [private]

Definition at line 106 of file Pi0FixedMassWindowCalibration.h.

Referenced by Pi0FixedMassWindowCalibration().

const EcalRecHitCollection* Pi0FixedMassWindowCalibration::recalibEcalRecHitCollection [private]

Definition at line 150 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop().

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

Definition at line 147 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop().

double Pi0FixedMassWindowCalibration::selePi0DetaBelt_ [private]

Definition at line 115 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop(), and Pi0FixedMassWindowCalibration().

double Pi0FixedMassWindowCalibration::selePi0DRBelt_ [private]

Definition at line 114 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop(), and Pi0FixedMassWindowCalibration().

double Pi0FixedMassWindowCalibration::selePi0EtBeltIsoRatioMax_ [private]

Definition at line 124 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop(), and Pi0FixedMassWindowCalibration().

double Pi0FixedMassWindowCalibration::selePi0MinvMeanFixed_ [private]

Definition at line 126 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop(), and Pi0FixedMassWindowCalibration().

double Pi0FixedMassWindowCalibration::selePi0MinvSigmaFixed_ [private]

Definition at line 127 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop(), and Pi0FixedMassWindowCalibration().

double Pi0FixedMassWindowCalibration::selePi0PtGammaOneMin_ [private]

Definition at line 111 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop(), and Pi0FixedMassWindowCalibration().

double Pi0FixedMassWindowCalibration::selePi0PtGammaTwoMin_ [private]

Definition at line 112 of file Pi0FixedMassWindowCalibration.h.

Referenced by Pi0FixedMassWindowCalibration().

double Pi0FixedMassWindowCalibration::selePi0PtPi0Min_ [private]

Definition at line 117 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop(), and Pi0FixedMassWindowCalibration().

double Pi0FixedMassWindowCalibration::selePi0S4S9GammaOneMin_ [private]

Definition at line 119 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop(), and Pi0FixedMassWindowCalibration().

double Pi0FixedMassWindowCalibration::selePi0S4S9GammaTwoMin_ [private]

Definition at line 120 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop(), and Pi0FixedMassWindowCalibration().

double Pi0FixedMassWindowCalibration::selePi0S9S25GammaOneMin_ [private]

Definition at line 121 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop(), and Pi0FixedMassWindowCalibration().

double Pi0FixedMassWindowCalibration::selePi0S9S25GammaTwoMin_ [private]

Definition at line 122 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop(), and Pi0FixedMassWindowCalibration().

ClusterShapeAlgo Pi0FixedMassWindowCalibration::shapeAlgo_ [private]

Definition at line 107 of file Pi0FixedMassWindowCalibration.h.

Referenced by duringLoop(), and Pi0FixedMassWindowCalibration().

TFile* Pi0FixedMassWindowCalibration::theFile [private]

Definition at line 154 of file Pi0FixedMassWindowCalibration.h.

unsigned int Pi0FixedMassWindowCalibration::theMaxLoops [private]

Definition at line 93 of file Pi0FixedMassWindowCalibration.h.

Referenced by endOfLoop().

edm::ParameterSet Pi0FixedMassWindowCalibration::theParameterSet [private]

Definition at line 141 of file Pi0FixedMassWindowCalibration.h.

Referenced by Pi0FixedMassWindowCalibration().

IslandClusterAlgo::VerbosityLevel Pi0FixedMassWindowCalibration::verbosity [private]

Definition at line 98 of file Pi0FixedMassWindowCalibration.h.

Referenced by Pi0FixedMassWindowCalibration().

double Pi0FixedMassWindowCalibration::wxtals[85][360][2] [private]

Definition at line 136 of file Pi0FixedMassWindowCalibration.h.

Referenced by beginOfJob(), duringLoop(), endOfLoop(), and startingNewLoop().


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